xref: /petsc/src/ts/impls/glee/glee.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
12302aa1bSDebojyoti Ghosh /*
22302aa1bSDebojyoti Ghosh   Code for time stepping with the General Linear with Error Estimation method
32302aa1bSDebojyoti Ghosh 
42302aa1bSDebojyoti Ghosh   Notes:
52302aa1bSDebojyoti Ghosh   The general system is written as
62302aa1bSDebojyoti Ghosh 
72302aa1bSDebojyoti Ghosh   Udot = F(t,U)
82302aa1bSDebojyoti Ghosh 
92302aa1bSDebojyoti Ghosh */
102302aa1bSDebojyoti Ghosh #include <petsc/private/tsimpl.h> /*I   "petscts.h"   I*/
112302aa1bSDebojyoti Ghosh #include <petscdm.h>
122302aa1bSDebojyoti Ghosh 
13642f3702SEmil Constantinescu static PetscBool  cited      = PETSC_FALSE;
149371c9d4SSatish Balay static const char citation[] = "@ARTICLE{Constantinescu_TR2016b,\n"
15642f3702SEmil Constantinescu                                " author = {Constantinescu, E.M.},\n"
16642f3702SEmil Constantinescu                                " title = {Estimating Global Errors in Time Stepping},\n"
17642f3702SEmil Constantinescu                                " journal = {ArXiv e-prints},\n"
18642f3702SEmil Constantinescu                                " year = 2016,\n"
19642f3702SEmil Constantinescu                                " adsurl = {http://adsabs.harvard.edu/abs/2015arXiv150305166C}\n}\n";
20642f3702SEmil Constantinescu 
212302aa1bSDebojyoti Ghosh static TSGLEEType TSGLEEDefaultType = TSGLEE35;
222302aa1bSDebojyoti Ghosh static PetscBool  TSGLEERegisterAllCalled;
232302aa1bSDebojyoti Ghosh static PetscBool  TSGLEEPackageInitialized;
242302aa1bSDebojyoti Ghosh static PetscInt   explicit_stage_time_id;
252302aa1bSDebojyoti Ghosh 
262302aa1bSDebojyoti Ghosh typedef struct _GLEETableau *GLEETableau;
272302aa1bSDebojyoti Ghosh struct _GLEETableau {
282302aa1bSDebojyoti Ghosh   char      *name;
292302aa1bSDebojyoti Ghosh   PetscInt   order;                     /* Classical approximation order of the method i*/
302302aa1bSDebojyoti Ghosh   PetscInt   s;                         /* Number of stages */
312302aa1bSDebojyoti Ghosh   PetscInt   r;                         /* Number of steps */
322302aa1bSDebojyoti Ghosh   PetscReal  gamma;                     /* LTE ratio */
332302aa1bSDebojyoti Ghosh   PetscReal *A, *B, *U, *V, *S, *F, *c; /* Tableau */
342302aa1bSDebojyoti Ghosh   PetscReal *Fembed;                    /* Embedded final method coefficients */
35fd96ebc2SDebojyoti Ghosh   PetscReal *Ferror;                    /* Coefficients for computing error   */
3657df6a1bSDebojyoti Ghosh   PetscReal *Serror;                    /* Coefficients for initializing the error   */
372302aa1bSDebojyoti Ghosh   PetscInt   pinterp;                   /* Interpolation order */
382302aa1bSDebojyoti Ghosh   PetscReal *binterp;                   /* Interpolation coefficients */
392302aa1bSDebojyoti Ghosh   PetscReal  ccfl;                      /* Placeholder for CFL coefficient relative to forward Euler  */
402302aa1bSDebojyoti Ghosh };
412302aa1bSDebojyoti Ghosh typedef struct _GLEETableauLink *GLEETableauLink;
422302aa1bSDebojyoti Ghosh struct _GLEETableauLink {
432302aa1bSDebojyoti Ghosh   struct _GLEETableau tab;
442302aa1bSDebojyoti Ghosh   GLEETableauLink     next;
452302aa1bSDebojyoti Ghosh };
462302aa1bSDebojyoti Ghosh static GLEETableauLink GLEETableauList;
472302aa1bSDebojyoti Ghosh 
482302aa1bSDebojyoti Ghosh typedef struct {
492302aa1bSDebojyoti Ghosh   GLEETableau  tableau;
502302aa1bSDebojyoti Ghosh   Vec         *Y;         /* Solution vector (along with auxiliary solution y~ or eps) */
512302aa1bSDebojyoti Ghosh   Vec         *X;         /* Temporary solution vector */
522302aa1bSDebojyoti Ghosh   Vec         *YStage;    /* Stage values */
532302aa1bSDebojyoti Ghosh   Vec         *YdotStage; /* Stage right hand side */
542302aa1bSDebojyoti Ghosh   Vec          W;         /* Right-hand-side for implicit stage solve */
552302aa1bSDebojyoti Ghosh   Vec          Ydot;      /* Work vector holding Ydot during residual evaluation */
560a01e1b2SEmil Constantinescu   Vec          yGErr;     /* Vector holding the global error after a step is completed */
57ec0e3ee2SEmil Constantinescu   PetscScalar *swork;     /* Scalar work (size of the number of stages)*/
58ec0e3ee2SEmil Constantinescu   PetscScalar *rwork;     /* Scalar work (size of the number of steps)*/
592302aa1bSDebojyoti Ghosh   PetscReal    scoeff;    /* shift = scoeff/dt */
602302aa1bSDebojyoti Ghosh   PetscReal    stage_time;
612302aa1bSDebojyoti Ghosh   TSStepStatus status;
622302aa1bSDebojyoti Ghosh } TS_GLEE;
632302aa1bSDebojyoti Ghosh 
642302aa1bSDebojyoti Ghosh /*MC
652302aa1bSDebojyoti Ghosh      TSGLEE23 - Second order three stage GLEE method
662302aa1bSDebojyoti Ghosh 
672302aa1bSDebojyoti Ghosh      This method has three stages.
682302aa1bSDebojyoti Ghosh      s = 3, r = 2
692302aa1bSDebojyoti Ghosh 
702302aa1bSDebojyoti Ghosh      Level: advanced
712302aa1bSDebojyoti Ghosh 
72db781477SPatrick Sanan .seealso: `TSGLEE`
7336c34d14SEmil Constantinescu M*/
742302aa1bSDebojyoti Ghosh /*MC
752302aa1bSDebojyoti Ghosh      TSGLEE24 - Second order four stage GLEE method
762302aa1bSDebojyoti Ghosh 
772302aa1bSDebojyoti Ghosh      This method has four stages.
782302aa1bSDebojyoti Ghosh      s = 4, r = 2
792302aa1bSDebojyoti Ghosh 
802302aa1bSDebojyoti Ghosh      Level: advanced
812302aa1bSDebojyoti Ghosh 
82db781477SPatrick Sanan .seealso: `TSGLEE`
8336c34d14SEmil Constantinescu M*/
842302aa1bSDebojyoti Ghosh /*MC
852302aa1bSDebojyoti Ghosh      TSGLEE25i - Second order five stage GLEE method
862302aa1bSDebojyoti Ghosh 
872302aa1bSDebojyoti Ghosh      This method has five stages.
882302aa1bSDebojyoti Ghosh      s = 5, r = 2
892302aa1bSDebojyoti Ghosh 
902302aa1bSDebojyoti Ghosh      Level: advanced
912302aa1bSDebojyoti Ghosh 
92db781477SPatrick Sanan .seealso: `TSGLEE`
9336c34d14SEmil Constantinescu M*/
942302aa1bSDebojyoti Ghosh /*MC
952302aa1bSDebojyoti Ghosh      TSGLEE35  - Third order five stage GLEE method
962302aa1bSDebojyoti Ghosh 
972302aa1bSDebojyoti Ghosh      This method has five stages.
982302aa1bSDebojyoti Ghosh      s = 5, r = 2
992302aa1bSDebojyoti Ghosh 
1002302aa1bSDebojyoti Ghosh      Level: advanced
1012302aa1bSDebojyoti Ghosh 
102db781477SPatrick Sanan .seealso: `TSGLEE`
10336c34d14SEmil Constantinescu M*/
1042302aa1bSDebojyoti Ghosh /*MC
1052302aa1bSDebojyoti Ghosh      TSGLEEEXRK2A  - Second order six stage GLEE method
1062302aa1bSDebojyoti Ghosh 
1072302aa1bSDebojyoti Ghosh      This method has six stages.
1082302aa1bSDebojyoti Ghosh      s = 6, r = 2
1092302aa1bSDebojyoti Ghosh 
1102302aa1bSDebojyoti Ghosh      Level: advanced
1112302aa1bSDebojyoti Ghosh 
112db781477SPatrick Sanan .seealso: `TSGLEE`
11336c34d14SEmil Constantinescu M*/
1142302aa1bSDebojyoti Ghosh /*MC
1152302aa1bSDebojyoti Ghosh      TSGLEERK32G1  - Third order eight stage GLEE method
1162302aa1bSDebojyoti Ghosh 
1172302aa1bSDebojyoti Ghosh      This method has eight stages.
1182302aa1bSDebojyoti Ghosh      s = 8, r = 2
1192302aa1bSDebojyoti Ghosh 
1202302aa1bSDebojyoti Ghosh      Level: advanced
1212302aa1bSDebojyoti Ghosh 
122db781477SPatrick Sanan .seealso: `TSGLEE`
12336c34d14SEmil Constantinescu M*/
1242302aa1bSDebojyoti Ghosh /*MC
1252302aa1bSDebojyoti Ghosh      TSGLEERK285EX  - Second order nine stage GLEE method
1262302aa1bSDebojyoti Ghosh 
1272302aa1bSDebojyoti Ghosh      This method has nine stages.
1282302aa1bSDebojyoti Ghosh      s = 9, r = 2
1292302aa1bSDebojyoti Ghosh 
1302302aa1bSDebojyoti Ghosh      Level: advanced
1312302aa1bSDebojyoti Ghosh 
132db781477SPatrick Sanan .seealso: `TSGLEE`
13336c34d14SEmil Constantinescu M*/
1342302aa1bSDebojyoti Ghosh 
1352302aa1bSDebojyoti Ghosh /*@C
1362302aa1bSDebojyoti Ghosh   TSGLEERegisterAll - Registers all of the General Linear with Error Estimation methods in TSGLEE
1372302aa1bSDebojyoti Ghosh 
1382302aa1bSDebojyoti Ghosh   Not Collective, but should be called by all processes which will need the schemes to be registered
1392302aa1bSDebojyoti Ghosh 
1402302aa1bSDebojyoti Ghosh   Level: advanced
1412302aa1bSDebojyoti Ghosh 
142db781477SPatrick Sanan .seealso: `TSGLEERegisterDestroy()`
1432302aa1bSDebojyoti Ghosh @*/
1449371c9d4SSatish Balay PetscErrorCode TSGLEERegisterAll(void) {
1452302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
1462302aa1bSDebojyoti Ghosh   if (TSGLEERegisterAllCalled) PetscFunctionReturn(0);
1472302aa1bSDebojyoti Ghosh   TSGLEERegisterAllCalled = PETSC_TRUE;
1482302aa1bSDebojyoti Ghosh 
1492302aa1bSDebojyoti Ghosh   {
150b1fbfd33SSatish Balay #define GAMMA 0.5
1512302aa1bSDebojyoti Ghosh     /* y-eps form */
1529371c9d4SSatish Balay     const PetscInt  p = 1, s = 3, r = 2;
1539371c9d4SSatish Balay     const PetscReal A[3][3] =
1549371c9d4SSatish Balay       {
1559371c9d4SSatish Balay         {1.0, 0,   0  },
1569371c9d4SSatish Balay         {0,   0.5, 0  },
1579371c9d4SSatish Balay         {0,   0.5, 0.5}
1589371c9d4SSatish Balay     },
1599371c9d4SSatish Balay                     B[2][3] = {{1.0, 0, 0}, {-2.0, 1.0, 1.0}}, U[3][2] = {{1.0, 0}, {1.0, 0.5}, {1.0, 0.5}}, V[2][2] = {{1, 0}, {0, 1}}, S[2] = {1, 0}, F[2] = {1, 0}, Fembed[2] = {1, 1 - GAMMA}, Ferror[2] = {0, 1}, Serror[2] = {1, 0};
1609566063dSJacob Faibussowitsch     PetscCall(TSGLEERegister(TSGLEEi1, p, s, r, GAMMA, &A[0][0], &B[0][0], &U[0][0], &V[0][0], S, F, NULL, Fembed, Ferror, Serror, 0, NULL));
161ab8c5912SEmil Constantinescu   }
162ab8c5912SEmil Constantinescu   {
163b1fbfd33SSatish Balay #undef GAMMA
164b1fbfd33SSatish Balay #define GAMMA 0.0
165ab8c5912SEmil Constantinescu     /* y-eps form */
1669371c9d4SSatish Balay     const PetscInt  p = 2, s = 3, r = 2;
1679371c9d4SSatish Balay     const PetscReal A[3][3] =
1689371c9d4SSatish Balay       {
1699371c9d4SSatish Balay         {0,    0,    0},
1709371c9d4SSatish Balay         {1,    0,    0},
1719371c9d4SSatish Balay         {0.25, 0.25, 0}
1729371c9d4SSatish Balay     },
1739371c9d4SSatish Balay                     B[2][3] = {{1.0 / 12.0, 1.0 / 12.0, 5.0 / 6.0}, {1.0 / 12.0, 1.0 / 12.0, -1.0 / 6.0}}, U[3][2] = {{1, 0}, {1, 10}, {1, -1}}, V[2][2] = {{1, 0}, {0, 1}}, S[2] = {1, 0}, F[2] = {1, 0}, Fembed[2] = {1, 1 - GAMMA}, Ferror[2] = {0, 1}, Serror[2] = {1, 0};
1749566063dSJacob Faibussowitsch     PetscCall(TSGLEERegister(TSGLEE23, p, s, r, GAMMA, &A[0][0], &B[0][0], &U[0][0], &V[0][0], S, F, NULL, Fembed, Ferror, Serror, 0, NULL));
1752302aa1bSDebojyoti Ghosh   }
1762302aa1bSDebojyoti Ghosh   {
177b1fbfd33SSatish Balay #undef GAMMA
178b1fbfd33SSatish Balay #define GAMMA 0.0
1792302aa1bSDebojyoti Ghosh     /* y-y~ form */
1809371c9d4SSatish Balay     const PetscInt  p = 2, s = 4, r = 2;
1819371c9d4SSatish Balay     const PetscReal A[4][4] =
1829371c9d4SSatish Balay       {
1839371c9d4SSatish Balay         {0,            0,            0,            0},
1849371c9d4SSatish Balay         {0.75,         0,            0,            0},
1859371c9d4SSatish Balay         {0.25,         29.0 / 60.0,  0,            0},
1869371c9d4SSatish Balay         {-21.0 / 44.0, 145.0 / 44.0, -20.0 / 11.0, 0}
1879371c9d4SSatish Balay     },
1889371c9d4SSatish Balay                     B[2][4] = {{109.0 / 275.0, 58.0 / 75.0, -37.0 / 110.0, 1.0 / 6.0}, {3.0 / 11.0, 0, 75.0 / 88.0, -1.0 / 8.0}}, U[4][2] = {{0, 1}, {75.0 / 58.0, -17.0 / 58.0}, {0, 1}, {0, 1}}, V[2][2] = {{1, 0}, {0, 1}}, S[2] = {1, 1}, F[2] = {1, 0}, Fembed[2] = {0, 1}, Ferror[2] = {-1.0 / (1.0 - GAMMA), 1.0 / (1.0 - GAMMA)}, Serror[2] = {1.0 - GAMMA, 1.0};
1899566063dSJacob Faibussowitsch     PetscCall(TSGLEERegister(TSGLEE24, p, s, r, GAMMA, &A[0][0], &B[0][0], &U[0][0], &V[0][0], S, F, NULL, Fembed, Ferror, Serror, 0, NULL));
1902302aa1bSDebojyoti Ghosh   }
1912302aa1bSDebojyoti Ghosh   {
192b1fbfd33SSatish Balay #undef GAMMA
193b1fbfd33SSatish Balay #define GAMMA 0.0
1942302aa1bSDebojyoti Ghosh     /* y-y~ form */
1959371c9d4SSatish Balay     const PetscInt  p = 2, s = 5, r = 2;
1969371c9d4SSatish Balay     const PetscReal A[5][5] =
1979371c9d4SSatish Balay       {
1989371c9d4SSatish Balay         {0,                       0,                       0,                       0,                      0},
1992302aa1bSDebojyoti Ghosh         {-0.94079244066783383269, 0,                       0,                       0,                      0},
2002302aa1bSDebojyoti Ghosh         {0.64228187778301907108,  0.10915356933958500042,  0,                       0,                      0},
2012302aa1bSDebojyoti Ghosh         {-0.51764297742287450812, 0.74414270351096040738,  -0.71404164927824538121, 0,                      0},
2029371c9d4SSatish Balay         {-0.44696561556825969206, -0.76768425657590196518, 0.20111608138142987881,  0.93828186737840469796, 0}
2039371c9d4SSatish Balay     },
2049371c9d4SSatish Balay                     B[2][5] = {{-0.029309178948150356153, -0.49671981884013874923, 0.34275801517650053274, 0.32941112623949194988, 0.85385985637229662276}, {0.78133219686062535272, 0.074238691892675897635, 0.57957363498384957966, -0.24638502829674959968, -0.18875949544040123033}}, U[5][2] = {{0.16911424754448327735, 0.83088575245551672265}, {0.53638465733199574340, 0.46361534266800425660}, {0.39901579167169582526, 0.60098420832830417474}, {0.87689005530618575480, 0.12310994469381424520}, {0.99056100455550913009, 0.0094389954444908699092}}, V[2][2] = {{1, 0}, {0, 1}}, S[2] = {1, 1}, F[2] = {1, 0}, Fembed[2] = {0, 1}, Ferror[2] = {-1.0 / (1.0 - GAMMA), 1.0 / (1.0 - GAMMA)}, Serror[2] = {1.0 - GAMMA, 1.0};
2059566063dSJacob Faibussowitsch     PetscCall(TSGLEERegister(TSGLEE25I, p, s, r, GAMMA, &A[0][0], &B[0][0], &U[0][0], &V[0][0], S, F, NULL, Fembed, Ferror, Serror, 0, NULL));
2062302aa1bSDebojyoti Ghosh   }
2072302aa1bSDebojyoti Ghosh   {
208b1fbfd33SSatish Balay #undef GAMMA
209b1fbfd33SSatish Balay #define GAMMA 0.0
2102302aa1bSDebojyoti Ghosh     /* y-y~ form */
2119371c9d4SSatish Balay     const PetscInt  p = 3, s = 5, r = 2;
2129371c9d4SSatish Balay     const PetscReal A[5][5] =
2139371c9d4SSatish Balay       {
2149371c9d4SSatish Balay         {0,                                                0,                                                 0,                                                 0,                                               0},
2152302aa1bSDebojyoti Ghosh         {-2169604947363702313.0 / 24313474998937147335.0,  0,                                                 0,                                                 0,                                               0},
2162302aa1bSDebojyoti Ghosh         {46526746497697123895.0 / 94116917485856474137.0,  -10297879244026594958.0 / 49199457603717988219.0,  0,                                                 0,                                               0},
2172302aa1bSDebojyoti Ghosh         {23364788935845982499.0 / 87425311444725389446.0,  -79205144337496116638.0 / 148994349441340815519.0, 40051189859317443782.0 / 36487615018004984309.0,   0,                                               0},
2189371c9d4SSatish Balay         {42089522664062539205.0 / 124911313006412840286.0, -15074384760342762939.0 / 137927286865289746282.0, -62274678522253371016.0 / 125918573676298591413.0, 13755475729852471739.0 / 79257927066651693390.0, 0}
2199371c9d4SSatish Balay     },
2209371c9d4SSatish Balay                     B[2][5] = {{61546696837458703723.0 / 56982519523786160813.0, -55810892792806293355.0 / 206957624151308356511.0, 24061048952676379087.0 / 158739347956038723465.0, 3577972206874351339.0 / 7599733370677197135.0, -59449832954780563947.0 / 137360038685338563670.0}, {-9738262186984159168.0 / 99299082461487742983.0, -32797097931948613195.0 / 61521565616362163366.0, 42895514606418420631.0 / 71714201188501437336.0, 22608567633166065068.0 / 55371917805607957003.0, 94655809487476459565.0 / 151517167160302729021.0}}, U[5][2] = {{70820309139834661559.0 / 80863923579509469826.0, 10043614439674808267.0 / 80863923579509469826.0}, {161694774978034105510.0 / 106187653640211060371.0, -55507121337823045139.0 / 106187653640211060371.0}, {78486094644566264568.0 / 88171030896733822981.0, 9684936252167558413.0 / 88171030896733822981.0}, {65394922146334854435.0 / 84570853840405479554.0, 19175931694070625119.0 / 84570853840405479554.0}, {8607282770183754108.0 / 108658046436496925911.0, 100050763666313171803.0 / 108658046436496925911.0}}, V[2][2] = {{1, 0}, {0, 1}}, S[2] = {1, 1}, F[2] = {1, 0}, Fembed[2] = {0, 1}, Ferror[2] = {-1.0 / (1.0 - GAMMA), 1.0 / (1.0 - GAMMA)}, Serror[2] = {1.0 - GAMMA, 1.0};
2219566063dSJacob Faibussowitsch     PetscCall(TSGLEERegister(TSGLEE35, p, s, r, GAMMA, &A[0][0], &B[0][0], &U[0][0], &V[0][0], S, F, NULL, Fembed, Ferror, Serror, 0, NULL));
2222302aa1bSDebojyoti Ghosh   }
2232302aa1bSDebojyoti Ghosh   {
224b1fbfd33SSatish Balay #undef GAMMA
225b1fbfd33SSatish Balay #define GAMMA 0.25
2262302aa1bSDebojyoti Ghosh     /* y-eps form */
2279371c9d4SSatish Balay     const PetscInt  p = 2, s = 6, r = 2;
2289371c9d4SSatish Balay     const PetscReal A[6][6] =
2299371c9d4SSatish Balay       {
2309371c9d4SSatish Balay         {0, 0, 0,    0,    0,   0},
2312302aa1bSDebojyoti Ghosh         {1, 0, 0,    0,    0,   0},
2322302aa1bSDebojyoti Ghosh         {0, 0, 0,    0,    0,   0},
2332302aa1bSDebojyoti Ghosh         {0, 0, 0.5,  0,    0,   0},
2342302aa1bSDebojyoti Ghosh         {0, 0, 0.25, 0.25, 0,   0},
2359371c9d4SSatish Balay         {0, 0, 0.25, 0.25, 0.5, 0}
2369371c9d4SSatish Balay     },
2379371c9d4SSatish Balay                     B[2][6] = {{0.5, 0.5, 0, 0, 0, 0}, {-2.0 / 3.0, -2.0 / 3.0, 1.0 / 3.0, 1.0 / 3.0, 1.0 / 3.0, 1.0 / 3.0}}, U[6][2] = {{1, 0}, {1, 0}, {1, 0.75}, {1, 0.75}, {1, 0.75}, {1, 0.75}}, V[2][2] = {{1, 0}, {0, 1}}, S[2] = {1, 0}, F[2] = {1, 0}, Fembed[2] = {1, 1 - GAMMA}, Ferror[2] = {0, 1}, Serror[2] = {1, 0};
2389566063dSJacob Faibussowitsch     PetscCall(TSGLEERegister(TSGLEEEXRK2A, p, s, r, GAMMA, &A[0][0], &B[0][0], &U[0][0], &V[0][0], S, F, NULL, Fembed, Ferror, Serror, 0, NULL));
2392302aa1bSDebojyoti Ghosh   }
2402302aa1bSDebojyoti Ghosh   {
241b1fbfd33SSatish Balay #undef GAMMA
242b1fbfd33SSatish Balay #define GAMMA 0.0
2432302aa1bSDebojyoti Ghosh     /* y-eps form */
2449371c9d4SSatish Balay     const PetscInt  p = 3, s = 8, r = 2;
2459371c9d4SSatish Balay     const PetscReal A[8][8] =
2469371c9d4SSatish Balay       {
2479371c9d4SSatish Balay         {0,           0,          0,          0,          0,         0,         0,         0},
2482302aa1bSDebojyoti Ghosh         {0.5,         0,          0,          0,          0,         0,         0,         0},
2492302aa1bSDebojyoti Ghosh         {-1,          2,          0,          0,          0,         0,         0,         0},
2502302aa1bSDebojyoti Ghosh         {1.0 / 6.0,   2.0 / 3.0,  1.0 / 6.0,  0,          0,         0,         0,         0},
2512302aa1bSDebojyoti Ghosh         {0,           0,          0,          0,          0,         0,         0,         0},
2522302aa1bSDebojyoti Ghosh         {-7.0 / 24.0, 1.0 / 3.0,  1.0 / 12.0, -1.0 / 8.0, 0.5,       0,         0,         0},
2532302aa1bSDebojyoti Ghosh         {7.0 / 6.0,   -4.0 / 3.0, -1.0 / 3.0, 0.5,        -1.0,      2.0,       0,         0},
2549371c9d4SSatish Balay         {0,           0,          0,          0,          1.0 / 6.0, 2.0 / 3.0, 1.0 / 6.0, 0}
2559371c9d4SSatish Balay     },
2569371c9d4SSatish Balay                     B[2][8] = {{1.0 / 6.0, 2.0 / 3.0, 1.0 / 6.0, 0, 0, 0, 0, 0}, {-1.0 / 6.0, -2.0 / 3.0, -1.0 / 6.0, 0, 1.0 / 6.0, 2.0 / 3.0, 1.0 / 6.0, 0}}, U[8][2] = {{1, 0}, {1, 0}, {1, 0}, {1, 0}, {1, 1}, {1, 1}, {1, 1}, {1, 1}}, V[2][2] = {{1, 0}, {0, 1}}, S[2] = {1, 0}, F[2] = {1, 0}, Fembed[2] = {1, 1 - GAMMA}, Ferror[2] = {0, 1}, Serror[2] = {1, 0};
2579566063dSJacob Faibussowitsch     PetscCall(TSGLEERegister(TSGLEERK32G1, p, s, r, GAMMA, &A[0][0], &B[0][0], &U[0][0], &V[0][0], S, F, NULL, Fembed, Ferror, Serror, 0, NULL));
2582302aa1bSDebojyoti Ghosh   }
2592302aa1bSDebojyoti Ghosh   {
260b1fbfd33SSatish Balay #undef GAMMA
261b1fbfd33SSatish Balay #define GAMMA 0.25
2622302aa1bSDebojyoti Ghosh     /* y-eps form */
2639371c9d4SSatish Balay     const PetscInt  p = 2, s = 9, r = 2;
2649371c9d4SSatish Balay     const PetscReal A[9][9] =
2659371c9d4SSatish Balay       {
2669371c9d4SSatish Balay         {0,                    0,                    0, 0,                     0,                    0,                    0,                     0,                    0},
2672302aa1bSDebojyoti Ghosh         {0.585786437626904966, 0,                    0, 0,                     0,                    0,                    0,                     0,                    0},
2682302aa1bSDebojyoti Ghosh         {0.149999999999999994, 0.849999999999999978, 0, 0,                     0,                    0,                    0,                     0,                    0},
2692302aa1bSDebojyoti Ghosh         {0,                    0,                    0, 0,                     0,                    0,                    0,                     0,                    0},
2702302aa1bSDebojyoti Ghosh         {0,                    0,                    0, 0.292893218813452483,  0,                    0,                    0,                     0,                    0},
2712302aa1bSDebojyoti Ghosh         {0,                    0,                    0, 0.0749999999999999972, 0.424999999999999989, 0,                    0,                     0,                    0},
2722302aa1bSDebojyoti Ghosh         {0,                    0,                    0, 0.176776695296636893,  0.176776695296636893, 0.146446609406726241, 0,                     0,                    0},
2732302aa1bSDebojyoti Ghosh         {0,                    0,                    0, 0.176776695296636893,  0.176776695296636893, 0.146446609406726241, 0.292893218813452483,  0,                    0},
2749371c9d4SSatish Balay         {0,                    0,                    0, 0.176776695296636893,  0.176776695296636893, 0.146446609406726241, 0.0749999999999999972, 0.424999999999999989, 0}
2759371c9d4SSatish Balay     },
2769371c9d4SSatish Balay                     B[2][9] = {{0.353553390593273786, 0.353553390593273786, 0.292893218813452483, 0, 0, 0, 0, 0, 0}, {-0.471404520791031678, -0.471404520791031678, -0.390524291751269959, 0.235702260395515839, 0.235702260395515839, 0.195262145875634979, 0.235702260395515839, 0.235702260395515839, 0.195262145875634979}}, U[9][2] = {{1, 0}, {1, 0}, {1, 0}, {1, 0.75}, {1, 0.75}, {1, 0.75}, {1, 0.75}, {1, 0.75}, {1, 0.75}}, V[2][2] = {{1, 0}, {0, 1}}, S[2] = {1, 0}, F[2] = {1, 0}, Fembed[2] = {1, 1 - GAMMA}, Ferror[2] = {0, 1}, Serror[2] = {1, 0};
2779566063dSJacob Faibussowitsch     PetscCall(TSGLEERegister(TSGLEERK285EX, p, s, r, GAMMA, &A[0][0], &B[0][0], &U[0][0], &V[0][0], S, F, NULL, Fembed, Ferror, Serror, 0, NULL));
2782302aa1bSDebojyoti Ghosh   }
2792302aa1bSDebojyoti Ghosh 
2802302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
2812302aa1bSDebojyoti Ghosh }
2822302aa1bSDebojyoti Ghosh 
2832302aa1bSDebojyoti Ghosh /*@C
2842302aa1bSDebojyoti Ghosh    TSGLEERegisterDestroy - Frees the list of schemes that were registered by TSGLEERegister().
2852302aa1bSDebojyoti Ghosh 
2862302aa1bSDebojyoti Ghosh    Not Collective
2872302aa1bSDebojyoti Ghosh 
2882302aa1bSDebojyoti Ghosh    Level: advanced
2892302aa1bSDebojyoti Ghosh 
290db781477SPatrick Sanan .seealso: `TSGLEERegister()`, `TSGLEERegisterAll()`
2912302aa1bSDebojyoti Ghosh @*/
2929371c9d4SSatish Balay PetscErrorCode TSGLEERegisterDestroy(void) {
2932302aa1bSDebojyoti Ghosh   GLEETableauLink link;
2942302aa1bSDebojyoti Ghosh 
2952302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
2962302aa1bSDebojyoti Ghosh   while ((link = GLEETableauList)) {
2972302aa1bSDebojyoti Ghosh     GLEETableau t   = &link->tab;
2982302aa1bSDebojyoti Ghosh     GLEETableauList = link->next;
2999566063dSJacob Faibussowitsch     PetscCall(PetscFree5(t->A, t->B, t->U, t->V, t->c));
300d0609cedSBarry Smith     PetscCall(PetscFree2(t->S, t->F));
301d0609cedSBarry Smith     PetscCall(PetscFree(t->Fembed));
302d0609cedSBarry Smith     PetscCall(PetscFree(t->Ferror));
303d0609cedSBarry Smith     PetscCall(PetscFree(t->Serror));
304d0609cedSBarry Smith     PetscCall(PetscFree(t->binterp));
305d0609cedSBarry Smith     PetscCall(PetscFree(t->name));
306d0609cedSBarry Smith     PetscCall(PetscFree(link));
3072302aa1bSDebojyoti Ghosh   }
3082302aa1bSDebojyoti Ghosh   TSGLEERegisterAllCalled = PETSC_FALSE;
3092302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
3102302aa1bSDebojyoti Ghosh }
3112302aa1bSDebojyoti Ghosh 
3122302aa1bSDebojyoti Ghosh /*@C
3132302aa1bSDebojyoti Ghosh   TSGLEEInitializePackage - This function initializes everything in the TSGLEE package. It is called
3148a690491SBarry Smith   from TSInitializePackage().
3152302aa1bSDebojyoti Ghosh 
3162302aa1bSDebojyoti Ghosh   Level: developer
3172302aa1bSDebojyoti Ghosh 
318db781477SPatrick Sanan .seealso: `PetscInitialize()`
3192302aa1bSDebojyoti Ghosh @*/
3209371c9d4SSatish Balay PetscErrorCode TSGLEEInitializePackage(void) {
3212302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
3222302aa1bSDebojyoti Ghosh   if (TSGLEEPackageInitialized) PetscFunctionReturn(0);
3232302aa1bSDebojyoti Ghosh   TSGLEEPackageInitialized = PETSC_TRUE;
3249566063dSJacob Faibussowitsch   PetscCall(TSGLEERegisterAll());
3259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposedDataRegister(&explicit_stage_time_id));
3269566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(TSGLEEFinalizePackage));
3272302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
3282302aa1bSDebojyoti Ghosh }
3292302aa1bSDebojyoti Ghosh 
3302302aa1bSDebojyoti Ghosh /*@C
3312302aa1bSDebojyoti Ghosh   TSGLEEFinalizePackage - This function destroys everything in the TSGLEE package. It is
3322302aa1bSDebojyoti Ghosh   called from PetscFinalize().
3332302aa1bSDebojyoti Ghosh 
3342302aa1bSDebojyoti Ghosh   Level: developer
3352302aa1bSDebojyoti Ghosh 
336db781477SPatrick Sanan .seealso: `PetscFinalize()`
3372302aa1bSDebojyoti Ghosh @*/
3389371c9d4SSatish Balay PetscErrorCode TSGLEEFinalizePackage(void) {
3392302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
3402302aa1bSDebojyoti Ghosh   TSGLEEPackageInitialized = PETSC_FALSE;
3419566063dSJacob Faibussowitsch   PetscCall(TSGLEERegisterDestroy());
3422302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
3432302aa1bSDebojyoti Ghosh }
3442302aa1bSDebojyoti Ghosh 
3452302aa1bSDebojyoti Ghosh /*@C
3462302aa1bSDebojyoti Ghosh    TSGLEERegister - register an GLEE scheme by providing the entries in the Butcher tableau
3472302aa1bSDebojyoti Ghosh 
3482302aa1bSDebojyoti Ghosh    Not Collective, but the same schemes should be registered on all processes on which they will be used
3492302aa1bSDebojyoti Ghosh 
3502302aa1bSDebojyoti Ghosh    Input Parameters:
3512302aa1bSDebojyoti Ghosh +  name   - identifier for method
3522302aa1bSDebojyoti Ghosh .  order  - order of method
3532302aa1bSDebojyoti Ghosh .  s - number of stages
3542302aa1bSDebojyoti Ghosh .  r - number of steps
3552302aa1bSDebojyoti Ghosh .  gamma - LTE ratio
3562302aa1bSDebojyoti Ghosh .  A - stage coefficients (dimension s*s, row-major)
3572302aa1bSDebojyoti Ghosh .  B - step completion coefficients (dimension r*s, row-major)
3582302aa1bSDebojyoti Ghosh .  U - method coefficients (dimension s*r, row-major)
3592302aa1bSDebojyoti Ghosh .  V - method coefficients (dimension r*r, row-major)
3602302aa1bSDebojyoti Ghosh .  S - starting coefficients
3612302aa1bSDebojyoti Ghosh .  F - finishing coefficients
3622302aa1bSDebojyoti Ghosh .  c - abscissa (dimension s; NULL to use row sums of A)
3632302aa1bSDebojyoti Ghosh .  Fembed - step completion coefficients for embedded method
364fd96ebc2SDebojyoti Ghosh .  Ferror - error computation coefficients
36557df6a1bSDebojyoti Ghosh .  Serror - error initialization coefficients
3662302aa1bSDebojyoti Ghosh .  pinterp - order of interpolation (0 if unavailable)
3672302aa1bSDebojyoti Ghosh -  binterp - array of interpolation coefficients (NULL if unavailable)
3682302aa1bSDebojyoti Ghosh 
3692302aa1bSDebojyoti Ghosh    Notes:
3702302aa1bSDebojyoti Ghosh    Several GLEE methods are provided, this function is only needed to create new methods.
3712302aa1bSDebojyoti Ghosh 
3722302aa1bSDebojyoti Ghosh    Level: advanced
3732302aa1bSDebojyoti Ghosh 
374db781477SPatrick Sanan .seealso: `TSGLEE`
3752302aa1bSDebojyoti Ghosh @*/
3769371c9d4SSatish Balay PetscErrorCode TSGLEERegister(TSGLEEType name, PetscInt order, PetscInt s, PetscInt r, PetscReal gamma, const PetscReal A[], const PetscReal B[], const PetscReal U[], const PetscReal V[], const PetscReal S[], const PetscReal F[], const PetscReal c[], const PetscReal Fembed[], const PetscReal Ferror[], const PetscReal Serror[], PetscInt pinterp, const PetscReal binterp[]) {
3772302aa1bSDebojyoti Ghosh   GLEETableauLink link;
3782302aa1bSDebojyoti Ghosh   GLEETableau     t;
3792302aa1bSDebojyoti Ghosh   PetscInt        i, j;
3802302aa1bSDebojyoti Ghosh 
3812302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
3829566063dSJacob Faibussowitsch   PetscCall(TSGLEEInitializePackage());
3839566063dSJacob Faibussowitsch   PetscCall(PetscNew(&link));
3842302aa1bSDebojyoti Ghosh   t = &link->tab;
3859566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &t->name));
3862302aa1bSDebojyoti Ghosh   t->order = order;
3872302aa1bSDebojyoti Ghosh   t->s     = s;
3882302aa1bSDebojyoti Ghosh   t->r     = r;
3892302aa1bSDebojyoti Ghosh   t->gamma = gamma;
3909566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(s * s, &t->A, r * r, &t->V, s, &t->c, r * s, &t->B, s * r, &t->U));
3919566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(r, &t->S, r, &t->F));
3929566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->A, A, s * s));
3939566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->B, B, r * s));
3949566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->U, U, s * r));
3959566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->V, V, r * r));
3969566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->S, S, r));
3979566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->F, F, r));
3982302aa1bSDebojyoti Ghosh   if (c) {
3999566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(t->c, c, s));
4002302aa1bSDebojyoti Ghosh   } else {
4019371c9d4SSatish Balay     for (i = 0; i < s; i++)
4029371c9d4SSatish Balay       for (j = 0, t->c[i] = 0; j < s; j++) t->c[i] += A[i * s + j];
4032302aa1bSDebojyoti Ghosh   }
4049566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(r, &t->Fembed));
4059566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(r, &t->Ferror));
4069566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(r, &t->Serror));
4079566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->Fembed, Fembed, r));
4089566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->Ferror, Ferror, r));
4099566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->Serror, Serror, r));
4102302aa1bSDebojyoti Ghosh   t->pinterp = pinterp;
4119566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s * pinterp, &t->binterp));
4129566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->binterp, binterp, s * pinterp));
4132302aa1bSDebojyoti Ghosh 
4142302aa1bSDebojyoti Ghosh   link->next      = GLEETableauList;
4152302aa1bSDebojyoti Ghosh   GLEETableauList = link;
4162302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
4172302aa1bSDebojyoti Ghosh }
4182302aa1bSDebojyoti Ghosh 
4199371c9d4SSatish Balay static PetscErrorCode TSEvaluateStep_GLEE(TS ts, PetscInt order, Vec X, PetscBool *done) {
4202302aa1bSDebojyoti Ghosh   TS_GLEE     *glee = (TS_GLEE *)ts->data;
4212302aa1bSDebojyoti Ghosh   GLEETableau  tab  = glee->tableau;
4229371c9d4SSatish Balay   PetscReal    h, *B = tab->B, *V = tab->V, *F = tab->F, *Fembed = tab->Fembed;
4232302aa1bSDebojyoti Ghosh   PetscInt     s = tab->s, r = tab->r, i, j;
4242302aa1bSDebojyoti Ghosh   Vec         *Y = glee->Y, *YdotStage = glee->YdotStage;
425ec0e3ee2SEmil Constantinescu   PetscScalar *ws = glee->swork, *wr = glee->rwork;
4262302aa1bSDebojyoti Ghosh 
4272302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
4282302aa1bSDebojyoti Ghosh 
4292302aa1bSDebojyoti Ghosh   switch (glee->status) {
4302302aa1bSDebojyoti Ghosh   case TS_STEP_INCOMPLETE:
4319371c9d4SSatish Balay   case TS_STEP_PENDING: h = ts->time_step; break;
4329371c9d4SSatish Balay   case TS_STEP_COMPLETE: h = ts->ptime - ts->ptime_prev; break;
4332302aa1bSDebojyoti Ghosh   default: SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus");
4342302aa1bSDebojyoti Ghosh   }
4352302aa1bSDebojyoti Ghosh 
4362302aa1bSDebojyoti Ghosh   if (order == tab->order) {
4372302aa1bSDebojyoti Ghosh     /* Note: Irrespective of whether status is TS_STEP_INCOMPLETE
438fd96ebc2SDebojyoti Ghosh              or TS_STEP_COMPLETE, glee->X has the solution at the
4392302aa1bSDebojyoti Ghosh              beginning of the time step. So no need to roll-back.
4402302aa1bSDebojyoti Ghosh     */
4412302aa1bSDebojyoti Ghosh     if (glee->status == TS_STEP_INCOMPLETE) {
4422302aa1bSDebojyoti Ghosh       for (i = 0; i < r; i++) {
4439566063dSJacob Faibussowitsch         PetscCall(VecZeroEntries(Y[i]));
444ec0e3ee2SEmil Constantinescu         for (j = 0; j < r; j++) wr[j] = V[i * r + j];
4459566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(Y[i], r, wr, glee->X));
446ec0e3ee2SEmil Constantinescu         for (j = 0; j < s; j++) ws[j] = h * B[i * s + j];
4479566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(Y[i], s, ws, YdotStage));
4482302aa1bSDebojyoti Ghosh       }
4499566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(X));
450ec0e3ee2SEmil Constantinescu       for (j = 0; j < r; j++) wr[j] = F[j];
4519566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(X, r, wr, Y));
4529566063dSJacob Faibussowitsch     } else PetscCall(VecCopy(ts->vec_sol, X));
4532302aa1bSDebojyoti Ghosh     PetscFunctionReturn(0);
4542302aa1bSDebojyoti Ghosh 
4552302aa1bSDebojyoti Ghosh   } else if (order == tab->order - 1) {
4562302aa1bSDebojyoti Ghosh     /* Complete with the embedded method (Fembed) */
4572302aa1bSDebojyoti Ghosh     for (i = 0; i < r; i++) {
4589566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(Y[i]));
459ec0e3ee2SEmil Constantinescu       for (j = 0; j < r; j++) wr[j] = V[i * r + j];
4609566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(Y[i], r, wr, glee->X));
461ec0e3ee2SEmil Constantinescu       for (j = 0; j < s; j++) ws[j] = h * B[i * s + j];
4629566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(Y[i], s, ws, YdotStage));
4632302aa1bSDebojyoti Ghosh     }
4649566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(X));
465ec0e3ee2SEmil Constantinescu     for (j = 0; j < r; j++) wr[j] = Fembed[j];
4669566063dSJacob Faibussowitsch     PetscCall(VecMAXPY(X, r, wr, Y));
467fd96ebc2SDebojyoti Ghosh 
4682302aa1bSDebojyoti Ghosh     if (done) *done = PETSC_TRUE;
4692302aa1bSDebojyoti Ghosh     PetscFunctionReturn(0);
4702302aa1bSDebojyoti Ghosh   }
4712302aa1bSDebojyoti Ghosh   if (done) *done = PETSC_FALSE;
47263a3b9bcSJacob Faibussowitsch   else SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "GLEE '%s' of order %" PetscInt_FMT " cannot evaluate step at order %" PetscInt_FMT, tab->name, tab->order, order);
4732302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
4742302aa1bSDebojyoti Ghosh }
4752302aa1bSDebojyoti Ghosh 
4769371c9d4SSatish Balay static PetscErrorCode TSStep_GLEE(TS ts) {
4772302aa1bSDebojyoti Ghosh   TS_GLEE       *glee = (TS_GLEE *)ts->data;
4782302aa1bSDebojyoti Ghosh   GLEETableau    tab  = glee->tableau;
4792302aa1bSDebojyoti Ghosh   const PetscInt s = tab->s, r = tab->r;
4809371c9d4SSatish Balay   PetscReal     *A = tab->A, *U = tab->U, *F = tab->F, *c = tab->c;
4819371c9d4SSatish Balay   Vec           *Y = glee->Y, *X = glee->X, *YStage = glee->YStage, *YdotStage = glee->YdotStage, W = glee->W;
4822302aa1bSDebojyoti Ghosh   SNES           snes;
483ec0e3ee2SEmil Constantinescu   PetscScalar   *ws = glee->swork, *wr = glee->rwork;
4842302aa1bSDebojyoti Ghosh   TSAdapt        adapt;
4852302aa1bSDebojyoti Ghosh   PetscInt       i, j, reject, next_scheme, its, lits;
4862302aa1bSDebojyoti Ghosh   PetscReal      next_time_step;
4872302aa1bSDebojyoti Ghosh   PetscReal      t;
4882302aa1bSDebojyoti Ghosh   PetscBool      accept;
4892302aa1bSDebojyoti Ghosh 
4902302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
4919566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(citation, &cited));
492642f3702SEmil Constantinescu 
4939566063dSJacob Faibussowitsch   for (i = 0; i < r; i++) PetscCall(VecCopy(Y[i], X[i]));
4942302aa1bSDebojyoti Ghosh 
4959566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
4962302aa1bSDebojyoti Ghosh   next_time_step = ts->time_step;
4972302aa1bSDebojyoti Ghosh   t              = ts->ptime;
4982302aa1bSDebojyoti Ghosh   accept         = PETSC_TRUE;
4992302aa1bSDebojyoti Ghosh   glee->status   = TS_STEP_INCOMPLETE;
5002302aa1bSDebojyoti Ghosh 
5012302aa1bSDebojyoti Ghosh   for (reject = 0; reject < ts->max_reject && !ts->reason; reject++, ts->reject++) {
5022302aa1bSDebojyoti Ghosh     PetscReal h = ts->time_step;
5039566063dSJacob Faibussowitsch     PetscCall(TSPreStep(ts));
5042302aa1bSDebojyoti Ghosh 
5052302aa1bSDebojyoti Ghosh     for (i = 0; i < s; i++) {
5062302aa1bSDebojyoti Ghosh       glee->stage_time = t + h * c[i];
5079566063dSJacob Faibussowitsch       PetscCall(TSPreStage(ts, glee->stage_time));
5082302aa1bSDebojyoti Ghosh 
5092302aa1bSDebojyoti Ghosh       if (A[i * s + i] == 0) { /* Explicit stage */
5109566063dSJacob Faibussowitsch         PetscCall(VecZeroEntries(YStage[i]));
511ec0e3ee2SEmil Constantinescu         for (j = 0; j < r; j++) wr[j] = U[i * r + j];
5129566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(YStage[i], r, wr, X));
513ec0e3ee2SEmil Constantinescu         for (j = 0; j < i; j++) ws[j] = h * A[i * s + j];
5149566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(YStage[i], i, ws, YdotStage));
5152302aa1bSDebojyoti Ghosh       } else { /* Implicit stage */
5162302aa1bSDebojyoti Ghosh         glee->scoeff = 1.0 / A[i * s + i];
5172302aa1bSDebojyoti Ghosh         /* compute right-hand-side */
5189566063dSJacob Faibussowitsch         PetscCall(VecZeroEntries(W));
519ec0e3ee2SEmil Constantinescu         for (j = 0; j < r; j++) wr[j] = U[i * r + j];
5209566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(W, r, wr, X));
521ec0e3ee2SEmil Constantinescu         for (j = 0; j < i; j++) ws[j] = h * A[i * s + j];
5229566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(W, i, ws, YdotStage));
5239566063dSJacob Faibussowitsch         PetscCall(VecScale(W, glee->scoeff / h));
5242302aa1bSDebojyoti Ghosh         /* set initial guess */
5259566063dSJacob Faibussowitsch         PetscCall(VecCopy(i > 0 ? YStage[i - 1] : ts->vec_sol, YStage[i]));
5262302aa1bSDebojyoti Ghosh         /* solve for this stage */
5279566063dSJacob Faibussowitsch         PetscCall(SNESSolve(snes, W, YStage[i]));
5289566063dSJacob Faibussowitsch         PetscCall(SNESGetIterationNumber(snes, &its));
5299566063dSJacob Faibussowitsch         PetscCall(SNESGetLinearSolveIterations(snes, &lits));
5309371c9d4SSatish Balay         ts->snes_its += its;
5319371c9d4SSatish Balay         ts->ksp_its += lits;
5322302aa1bSDebojyoti Ghosh       }
5339566063dSJacob Faibussowitsch       PetscCall(TSGetAdapt(ts, &adapt));
5349566063dSJacob Faibussowitsch       PetscCall(TSAdaptCheckStage(adapt, ts, glee->stage_time, YStage[i], &accept));
5352302aa1bSDebojyoti Ghosh       if (!accept) goto reject_step;
5369566063dSJacob Faibussowitsch       PetscCall(TSPostStage(ts, glee->stage_time, i, YStage));
5379566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSFunction(ts, t + h * c[i], YStage[i], YdotStage[i]));
5382302aa1bSDebojyoti Ghosh     }
5399566063dSJacob Faibussowitsch     PetscCall(TSEvaluateStep(ts, tab->order, ts->vec_sol, NULL));
5402302aa1bSDebojyoti Ghosh     glee->status = TS_STEP_PENDING;
5412302aa1bSDebojyoti Ghosh 
5422302aa1bSDebojyoti Ghosh     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
5439566063dSJacob Faibussowitsch     PetscCall(TSGetAdapt(ts, &adapt));
5449566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidatesClear(adapt));
5459566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidateAdd(adapt, tab->name, tab->order, 1, tab->ccfl, (PetscReal)tab->s, PETSC_TRUE));
5469566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, &next_scheme, &next_time_step, &accept));
5472302aa1bSDebojyoti Ghosh     if (accept) {
5482302aa1bSDebojyoti Ghosh       /* ignore next_scheme for now */
5492302aa1bSDebojyoti Ghosh       ts->ptime += ts->time_step;
5502302aa1bSDebojyoti Ghosh       ts->time_step = next_time_step;
5512302aa1bSDebojyoti Ghosh       glee->status  = TS_STEP_COMPLETE;
5520a01e1b2SEmil Constantinescu       /* compute and store the global error */
5530a01e1b2SEmil Constantinescu       /* Note: this is not needed if TSAdaptGLEE is not used */
5549566063dSJacob Faibussowitsch       PetscCall(TSGetTimeError(ts, 0, &(glee->yGErr)));
5559566063dSJacob Faibussowitsch       PetscCall(PetscObjectComposedDataSetReal((PetscObject)ts->vec_sol, explicit_stage_time_id, ts->ptime));
5562302aa1bSDebojyoti Ghosh       break;
5572302aa1bSDebojyoti Ghosh     } else { /* Roll back the current step */
558ec0e3ee2SEmil Constantinescu       for (j = 0; j < r; j++) wr[j] = F[j];
5599566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(ts->vec_sol, r, wr, X));
5602302aa1bSDebojyoti Ghosh       ts->time_step = next_time_step;
5612302aa1bSDebojyoti Ghosh       glee->status  = TS_STEP_INCOMPLETE;
5622302aa1bSDebojyoti Ghosh     }
5639371c9d4SSatish Balay   reject_step:
5649371c9d4SSatish Balay     continue;
5652302aa1bSDebojyoti Ghosh   }
5662302aa1bSDebojyoti Ghosh   if (glee->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
5672302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
5682302aa1bSDebojyoti Ghosh }
5692302aa1bSDebojyoti Ghosh 
5709371c9d4SSatish Balay static PetscErrorCode TSInterpolate_GLEE(TS ts, PetscReal itime, Vec X) {
5712302aa1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE *)ts->data;
5722302aa1bSDebojyoti Ghosh   PetscInt         s = glee->tableau->s, pinterp = glee->tableau->pinterp, i, j;
5732302aa1bSDebojyoti Ghosh   PetscReal        h, tt, t;
5742302aa1bSDebojyoti Ghosh   PetscScalar     *b;
5752302aa1bSDebojyoti Ghosh   const PetscReal *B = glee->tableau->binterp;
5762302aa1bSDebojyoti Ghosh 
5772302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
5783c633725SBarry Smith   PetscCheck(B, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSGLEE %s does not have an interpolation formula", glee->tableau->name);
5792302aa1bSDebojyoti Ghosh   switch (glee->status) {
5802302aa1bSDebojyoti Ghosh   case TS_STEP_INCOMPLETE:
5812302aa1bSDebojyoti Ghosh   case TS_STEP_PENDING:
5822302aa1bSDebojyoti Ghosh     h = ts->time_step;
5832302aa1bSDebojyoti Ghosh     t = (itime - ts->ptime) / h;
5842302aa1bSDebojyoti Ghosh     break;
5852302aa1bSDebojyoti Ghosh   case TS_STEP_COMPLETE:
586ed450d42SEmil Constantinescu     h = ts->ptime - ts->ptime_prev;
5872302aa1bSDebojyoti Ghosh     t = (itime - ts->ptime) / h + 1; /* In the interval [0,1] */
5882302aa1bSDebojyoti Ghosh     break;
5892302aa1bSDebojyoti Ghosh   default: SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus");
5902302aa1bSDebojyoti Ghosh   }
5919566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s, &b));
5922302aa1bSDebojyoti Ghosh   for (i = 0; i < s; i++) b[i] = 0;
5932302aa1bSDebojyoti Ghosh   for (j = 0, tt = t; j < pinterp; j++, tt *= t) {
594ad540459SPierre Jolivet     for (i = 0; i < s; i++) b[i] += h * B[i * pinterp + j] * tt;
5952302aa1bSDebojyoti Ghosh   }
5969566063dSJacob Faibussowitsch   PetscCall(VecCopy(glee->YStage[0], X));
5979566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(X, s, b, glee->YdotStage));
5989566063dSJacob Faibussowitsch   PetscCall(PetscFree(b));
5992302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
6002302aa1bSDebojyoti Ghosh }
6012302aa1bSDebojyoti Ghosh 
6022302aa1bSDebojyoti Ghosh /*------------------------------------------------------------*/
6039371c9d4SSatish Balay static PetscErrorCode TSReset_GLEE(TS ts) {
6042302aa1bSDebojyoti Ghosh   TS_GLEE *glee = (TS_GLEE *)ts->data;
6052302aa1bSDebojyoti Ghosh   PetscInt s, r;
6062302aa1bSDebojyoti Ghosh 
6072302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
6082302aa1bSDebojyoti Ghosh   if (!glee->tableau) PetscFunctionReturn(0);
6092302aa1bSDebojyoti Ghosh   s = glee->tableau->s;
6102302aa1bSDebojyoti Ghosh   r = glee->tableau->r;
6119566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(r, &glee->Y));
6129566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(r, &glee->X));
6139566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(s, &glee->YStage));
6149566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(s, &glee->YdotStage));
6159566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&glee->Ydot));
6169566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&glee->yGErr));
6179566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&glee->W));
6189566063dSJacob Faibussowitsch   PetscCall(PetscFree2(glee->swork, glee->rwork));
6192302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
6202302aa1bSDebojyoti Ghosh }
6212302aa1bSDebojyoti Ghosh 
6229371c9d4SSatish Balay static PetscErrorCode TSGLEEGetVecs(TS ts, DM dm, Vec *Ydot) {
6232302aa1bSDebojyoti Ghosh   TS_GLEE *glee = (TS_GLEE *)ts->data;
6242302aa1bSDebojyoti Ghosh 
6252302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
6262302aa1bSDebojyoti Ghosh   if (Ydot) {
6272302aa1bSDebojyoti Ghosh     if (dm && dm != ts->dm) {
6289566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm, "TSGLEE_Ydot", Ydot));
6292302aa1bSDebojyoti Ghosh     } else *Ydot = glee->Ydot;
6302302aa1bSDebojyoti Ghosh   }
6312302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
6322302aa1bSDebojyoti Ghosh }
6332302aa1bSDebojyoti Ghosh 
6349371c9d4SSatish Balay static PetscErrorCode TSGLEERestoreVecs(TS ts, DM dm, Vec *Ydot) {
6352302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
6362302aa1bSDebojyoti Ghosh   if (Ydot) {
63748a46eb9SPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSGLEE_Ydot", Ydot));
6382302aa1bSDebojyoti Ghosh   }
6392302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
6402302aa1bSDebojyoti Ghosh }
6412302aa1bSDebojyoti Ghosh 
6422302aa1bSDebojyoti Ghosh /*
6432302aa1bSDebojyoti Ghosh   This defines the nonlinear equation that is to be solved with SNES
6442302aa1bSDebojyoti Ghosh */
6459371c9d4SSatish Balay static PetscErrorCode SNESTSFormFunction_GLEE(SNES snes, Vec X, Vec F, TS ts) {
6462302aa1bSDebojyoti Ghosh   TS_GLEE  *glee = (TS_GLEE *)ts->data;
6472302aa1bSDebojyoti Ghosh   DM        dm, dmsave;
6482302aa1bSDebojyoti Ghosh   Vec       Ydot;
6492302aa1bSDebojyoti Ghosh   PetscReal shift = glee->scoeff / ts->time_step;
6502302aa1bSDebojyoti Ghosh 
6512302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
6529566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
6539566063dSJacob Faibussowitsch   PetscCall(TSGLEEGetVecs(ts, dm, &Ydot));
6542302aa1bSDebojyoti Ghosh   /* Set Ydot = shift*X */
6559566063dSJacob Faibussowitsch   PetscCall(VecCopy(X, Ydot));
6569566063dSJacob Faibussowitsch   PetscCall(VecScale(Ydot, shift));
6572302aa1bSDebojyoti Ghosh   dmsave = ts->dm;
6582302aa1bSDebojyoti Ghosh   ts->dm = dm;
6592302aa1bSDebojyoti Ghosh 
6609566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts, glee->stage_time, X, Ydot, F, PETSC_FALSE));
6612302aa1bSDebojyoti Ghosh 
6622302aa1bSDebojyoti Ghosh   ts->dm = dmsave;
6639566063dSJacob Faibussowitsch   PetscCall(TSGLEERestoreVecs(ts, dm, &Ydot));
6642302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
6652302aa1bSDebojyoti Ghosh }
6662302aa1bSDebojyoti Ghosh 
6679371c9d4SSatish Balay static PetscErrorCode SNESTSFormJacobian_GLEE(SNES snes, Vec X, Mat A, Mat B, TS ts) {
6682302aa1bSDebojyoti Ghosh   TS_GLEE  *glee = (TS_GLEE *)ts->data;
6692302aa1bSDebojyoti Ghosh   DM        dm, dmsave;
6702302aa1bSDebojyoti Ghosh   Vec       Ydot;
6712302aa1bSDebojyoti Ghosh   PetscReal shift = glee->scoeff / ts->time_step;
6722302aa1bSDebojyoti Ghosh 
6732302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
6749566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
6759566063dSJacob Faibussowitsch   PetscCall(TSGLEEGetVecs(ts, dm, &Ydot));
6762302aa1bSDebojyoti Ghosh   /* glee->Ydot has already been computed in SNESTSFormFunction_GLEE (SNES guarantees this) */
6772302aa1bSDebojyoti Ghosh   dmsave = ts->dm;
6782302aa1bSDebojyoti Ghosh   ts->dm = dm;
6792302aa1bSDebojyoti Ghosh 
6809566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts, glee->stage_time, X, Ydot, shift, A, B, PETSC_FALSE));
6812302aa1bSDebojyoti Ghosh 
6822302aa1bSDebojyoti Ghosh   ts->dm = dmsave;
6839566063dSJacob Faibussowitsch   PetscCall(TSGLEERestoreVecs(ts, dm, &Ydot));
6842302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
6852302aa1bSDebojyoti Ghosh }
6862302aa1bSDebojyoti Ghosh 
6879371c9d4SSatish Balay static PetscErrorCode DMCoarsenHook_TSGLEE(DM fine, DM coarse, void *ctx) {
6882302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
6892302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
6902302aa1bSDebojyoti Ghosh }
6912302aa1bSDebojyoti Ghosh 
6929371c9d4SSatish Balay static PetscErrorCode DMRestrictHook_TSGLEE(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) {
6932302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
6942302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
6952302aa1bSDebojyoti Ghosh }
6962302aa1bSDebojyoti Ghosh 
6979371c9d4SSatish Balay static PetscErrorCode DMSubDomainHook_TSGLEE(DM dm, DM subdm, void *ctx) {
6982302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
6992302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
7002302aa1bSDebojyoti Ghosh }
7012302aa1bSDebojyoti Ghosh 
7029371c9d4SSatish Balay static PetscErrorCode DMSubDomainRestrictHook_TSGLEE(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) {
7032302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7042302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
7052302aa1bSDebojyoti Ghosh }
7062302aa1bSDebojyoti Ghosh 
7079371c9d4SSatish Balay static PetscErrorCode TSSetUp_GLEE(TS ts) {
7082302aa1bSDebojyoti Ghosh   TS_GLEE    *glee = (TS_GLEE *)ts->data;
7092302aa1bSDebojyoti Ghosh   GLEETableau tab;
7102302aa1bSDebojyoti Ghosh   PetscInt    s, r;
7112302aa1bSDebojyoti Ghosh   DM          dm;
7122302aa1bSDebojyoti Ghosh 
7132302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
71448a46eb9SPierre Jolivet   if (!glee->tableau) PetscCall(TSGLEESetType(ts, TSGLEEDefaultType));
7152302aa1bSDebojyoti Ghosh   tab = glee->tableau;
7162302aa1bSDebojyoti Ghosh   s   = tab->s;
7172302aa1bSDebojyoti Ghosh   r   = tab->r;
7189566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vec_sol, r, &glee->Y));
7199566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vec_sol, r, &glee->X));
7209566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vec_sol, s, &glee->YStage));
7219566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vec_sol, s, &glee->YdotStage));
7229566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &glee->Ydot));
7239566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &glee->yGErr));
7249566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(glee->yGErr));
7259566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &glee->W));
7269566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(s, &glee->swork, r, &glee->rwork));
7279566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
7289566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSGLEE, DMRestrictHook_TSGLEE, ts));
7299566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSGLEE, DMSubDomainRestrictHook_TSGLEE, ts));
7302302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
7312302aa1bSDebojyoti Ghosh }
7322302aa1bSDebojyoti Ghosh 
7339371c9d4SSatish Balay PetscErrorCode TSStartingMethod_GLEE(TS ts) {
7342302aa1bSDebojyoti Ghosh   TS_GLEE    *glee = (TS_GLEE *)ts->data;
7356e2e232bSDebojyoti Ghosh   GLEETableau tab  = glee->tableau;
7366e2e232bSDebojyoti Ghosh   PetscInt    r    = tab->r, i;
7376e2e232bSDebojyoti Ghosh   PetscReal  *S    = tab->S;
7382302aa1bSDebojyoti Ghosh 
7392302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7402302aa1bSDebojyoti Ghosh   for (i = 0; i < r; i++) {
7419566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(glee->Y[i]));
7429566063dSJacob Faibussowitsch     PetscCall(VecAXPY(glee->Y[i], S[i], ts->vec_sol));
7432302aa1bSDebojyoti Ghosh   }
7442302aa1bSDebojyoti Ghosh 
7452302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
7462302aa1bSDebojyoti Ghosh }
7472302aa1bSDebojyoti Ghosh 
7482302aa1bSDebojyoti Ghosh /*------------------------------------------------------------*/
7492302aa1bSDebojyoti Ghosh 
7509371c9d4SSatish Balay static PetscErrorCode TSSetFromOptions_GLEE(TS ts, PetscOptionItems *PetscOptionsObject) {
7512302aa1bSDebojyoti Ghosh   char gleetype[256];
7522302aa1bSDebojyoti Ghosh 
7532302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
754d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "GLEE ODE solver options");
7552302aa1bSDebojyoti Ghosh   {
7562302aa1bSDebojyoti Ghosh     GLEETableauLink link;
7572302aa1bSDebojyoti Ghosh     PetscInt        count, choice;
7582302aa1bSDebojyoti Ghosh     PetscBool       flg;
7592302aa1bSDebojyoti Ghosh     const char    **namelist;
7602302aa1bSDebojyoti Ghosh 
7619566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(gleetype, TSGLEEDefaultType, sizeof(gleetype)));
7629371c9d4SSatish Balay     for (link = GLEETableauList, count = 0; link; link = link->next, count++)
7639371c9d4SSatish Balay       ;
7649566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(count, (char ***)&namelist));
7652302aa1bSDebojyoti Ghosh     for (link = GLEETableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name;
7669566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEList("-ts_glee_type", "Family of GLEE method", "TSGLEESetType", (const char *const *)namelist, count, gleetype, &choice, &flg));
7679566063dSJacob Faibussowitsch     PetscCall(TSGLEESetType(ts, flg ? namelist[choice] : gleetype));
7689566063dSJacob Faibussowitsch     PetscCall(PetscFree(namelist));
7692302aa1bSDebojyoti Ghosh   }
770d0609cedSBarry Smith   PetscOptionsHeadEnd();
7712302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
7722302aa1bSDebojyoti Ghosh }
7732302aa1bSDebojyoti Ghosh 
7749371c9d4SSatish Balay static PetscErrorCode TSView_GLEE(TS ts, PetscViewer viewer) {
7752302aa1bSDebojyoti Ghosh   TS_GLEE    *glee = (TS_GLEE *)ts->data;
7762302aa1bSDebojyoti Ghosh   GLEETableau tab  = glee->tableau;
7772302aa1bSDebojyoti Ghosh   PetscBool   iascii;
7782302aa1bSDebojyoti Ghosh 
7792302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7809566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
7812302aa1bSDebojyoti Ghosh   if (iascii) {
7822302aa1bSDebojyoti Ghosh     TSGLEEType gleetype;
7832302aa1bSDebojyoti Ghosh     char       buf[512];
7849566063dSJacob Faibussowitsch     PetscCall(TSGLEEGetType(ts, &gleetype));
7859566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  GLEE type %s\n", gleetype));
7869566063dSJacob Faibussowitsch     PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, tab->c));
7879566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Abscissa     c = %s\n", buf));
7882302aa1bSDebojyoti Ghosh     /* Note: print out r as well */
7892302aa1bSDebojyoti Ghosh   }
7902302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
7912302aa1bSDebojyoti Ghosh }
7922302aa1bSDebojyoti Ghosh 
7939371c9d4SSatish Balay static PetscErrorCode TSLoad_GLEE(TS ts, PetscViewer viewer) {
7942302aa1bSDebojyoti Ghosh   SNES    snes;
7952302aa1bSDebojyoti Ghosh   TSAdapt tsadapt;
7962302aa1bSDebojyoti Ghosh 
7972302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7989566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &tsadapt));
7999566063dSJacob Faibussowitsch   PetscCall(TSAdaptLoad(tsadapt, viewer));
8009566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
8019566063dSJacob Faibussowitsch   PetscCall(SNESLoad(snes, viewer));
8022302aa1bSDebojyoti Ghosh   /* function and Jacobian context for SNES when used with TS is always ts object */
8039566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, NULL, NULL, ts));
8049566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, ts));
8052302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8062302aa1bSDebojyoti Ghosh }
8072302aa1bSDebojyoti Ghosh 
8082302aa1bSDebojyoti Ghosh /*@C
8092302aa1bSDebojyoti Ghosh   TSGLEESetType - Set the type of GLEE scheme
8102302aa1bSDebojyoti Ghosh 
8112302aa1bSDebojyoti Ghosh   Logically collective
8122302aa1bSDebojyoti Ghosh 
813d8d19677SJose E. Roman   Input Parameters:
8142302aa1bSDebojyoti Ghosh +  ts - timestepping context
8152302aa1bSDebojyoti Ghosh -  gleetype - type of GLEE-scheme
8162302aa1bSDebojyoti Ghosh 
8172302aa1bSDebojyoti Ghosh   Level: intermediate
8182302aa1bSDebojyoti Ghosh 
819db781477SPatrick Sanan .seealso: `TSGLEEGetType()`, `TSGLEE`
8202302aa1bSDebojyoti Ghosh @*/
8219371c9d4SSatish Balay PetscErrorCode TSGLEESetType(TS ts, TSGLEEType gleetype) {
8222302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8232302aa1bSDebojyoti Ghosh   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
824b92453a8SLisandro Dalcin   PetscValidCharPointer(gleetype, 2);
825cac4c232SBarry Smith   PetscTryMethod(ts, "TSGLEESetType_C", (TS, TSGLEEType), (ts, gleetype));
8262302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8272302aa1bSDebojyoti Ghosh }
8282302aa1bSDebojyoti Ghosh 
8292302aa1bSDebojyoti Ghosh /*@C
8302302aa1bSDebojyoti Ghosh   TSGLEEGetType - Get the type of GLEE scheme
8312302aa1bSDebojyoti Ghosh 
8322302aa1bSDebojyoti Ghosh   Logically collective
8332302aa1bSDebojyoti Ghosh 
8342302aa1bSDebojyoti Ghosh   Input Parameter:
8352302aa1bSDebojyoti Ghosh .  ts - timestepping context
8362302aa1bSDebojyoti Ghosh 
8372302aa1bSDebojyoti Ghosh   Output Parameter:
8382302aa1bSDebojyoti Ghosh .  gleetype - type of GLEE-scheme
8392302aa1bSDebojyoti Ghosh 
8402302aa1bSDebojyoti Ghosh   Level: intermediate
8412302aa1bSDebojyoti Ghosh 
842db781477SPatrick Sanan .seealso: `TSGLEESetType()`
8432302aa1bSDebojyoti Ghosh @*/
8449371c9d4SSatish Balay PetscErrorCode TSGLEEGetType(TS ts, TSGLEEType *gleetype) {
8452302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8462302aa1bSDebojyoti Ghosh   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
847cac4c232SBarry Smith   PetscUseMethod(ts, "TSGLEEGetType_C", (TS, TSGLEEType *), (ts, gleetype));
8482302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8492302aa1bSDebojyoti Ghosh }
8502302aa1bSDebojyoti Ghosh 
8519371c9d4SSatish Balay PetscErrorCode TSGLEEGetType_GLEE(TS ts, TSGLEEType *gleetype) {
8522302aa1bSDebojyoti Ghosh   TS_GLEE *glee = (TS_GLEE *)ts->data;
8532302aa1bSDebojyoti Ghosh 
8542302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
85548a46eb9SPierre Jolivet   if (!glee->tableau) PetscCall(TSGLEESetType(ts, TSGLEEDefaultType));
8562302aa1bSDebojyoti Ghosh   *gleetype = glee->tableau->name;
8572302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8582302aa1bSDebojyoti Ghosh }
8599371c9d4SSatish Balay PetscErrorCode TSGLEESetType_GLEE(TS ts, TSGLEEType gleetype) {
8602302aa1bSDebojyoti Ghosh   TS_GLEE        *glee = (TS_GLEE *)ts->data;
8612302aa1bSDebojyoti Ghosh   PetscBool       match;
8622302aa1bSDebojyoti Ghosh   GLEETableauLink link;
8632302aa1bSDebojyoti Ghosh 
8642302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8652302aa1bSDebojyoti Ghosh   if (glee->tableau) {
8669566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(glee->tableau->name, gleetype, &match));
8672302aa1bSDebojyoti Ghosh     if (match) PetscFunctionReturn(0);
8682302aa1bSDebojyoti Ghosh   }
8692302aa1bSDebojyoti Ghosh   for (link = GLEETableauList; link; link = link->next) {
8709566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(link->tab.name, gleetype, &match));
8712302aa1bSDebojyoti Ghosh     if (match) {
8729566063dSJacob Faibussowitsch       PetscCall(TSReset_GLEE(ts));
8732302aa1bSDebojyoti Ghosh       glee->tableau = &link->tab;
8742302aa1bSDebojyoti Ghosh       PetscFunctionReturn(0);
8752302aa1bSDebojyoti Ghosh     }
8762302aa1bSDebojyoti Ghosh   }
87798921bdaSJacob Faibussowitsch   SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", gleetype);
8782302aa1bSDebojyoti Ghosh }
8792302aa1bSDebojyoti Ghosh 
8809371c9d4SSatish Balay static PetscErrorCode TSGetStages_GLEE(TS ts, PetscInt *ns, Vec **Y) {
8812302aa1bSDebojyoti Ghosh   TS_GLEE *glee = (TS_GLEE *)ts->data;
8822302aa1bSDebojyoti Ghosh 
8832302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8840429704eSStefano Zampini   if (ns) *ns = glee->tableau->s;
885d5509560SEmil Constantinescu   if (Y) *Y = glee->YStage;
8862302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8872302aa1bSDebojyoti Ghosh }
8882302aa1bSDebojyoti Ghosh 
8899371c9d4SSatish Balay PetscErrorCode TSGetSolutionComponents_GLEE(TS ts, PetscInt *n, Vec *Y) {
8902302aa1bSDebojyoti Ghosh   TS_GLEE    *glee = (TS_GLEE *)ts->data;
8912302aa1bSDebojyoti Ghosh   GLEETableau tab  = glee->tableau;
8922302aa1bSDebojyoti Ghosh 
8932302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8944cdd57e5SDebojyoti Ghosh   if (!Y) *n = tab->r;
8952302aa1bSDebojyoti Ghosh   else {
8964cdd57e5SDebojyoti Ghosh     if ((*n >= 0) && (*n < tab->r)) {
8979566063dSJacob Faibussowitsch       PetscCall(VecCopy(glee->Y[*n], *Y));
89863a3b9bcSJacob Faibussowitsch     } else SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Second argument (%" PetscInt_FMT ") out of range[0,%" PetscInt_FMT "].", *n, tab->r - 1);
8992302aa1bSDebojyoti Ghosh   }
9002302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
9012302aa1bSDebojyoti Ghosh }
9022302aa1bSDebojyoti Ghosh 
9039371c9d4SSatish Balay PetscErrorCode TSGetAuxSolution_GLEE(TS ts, Vec *X) {
9044cdd57e5SDebojyoti Ghosh   TS_GLEE     *glee = (TS_GLEE *)ts->data;
9054cdd57e5SDebojyoti Ghosh   GLEETableau  tab  = glee->tableau;
9064cdd57e5SDebojyoti Ghosh   PetscReal   *F    = tab->Fembed;
9074cdd57e5SDebojyoti Ghosh   PetscInt     r    = tab->r;
9084cdd57e5SDebojyoti Ghosh   Vec         *Y    = glee->Y;
909ec0e3ee2SEmil Constantinescu   PetscScalar *wr   = glee->rwork;
910ec0e3ee2SEmil Constantinescu   PetscInt     i;
9114cdd57e5SDebojyoti Ghosh 
9124cdd57e5SDebojyoti Ghosh   PetscFunctionBegin;
9139566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(*X));
914ec0e3ee2SEmil Constantinescu   for (i = 0; i < r; i++) wr[i] = F[i];
9159566063dSJacob Faibussowitsch   PetscCall(VecMAXPY((*X), r, wr, Y));
9164cdd57e5SDebojyoti Ghosh   PetscFunctionReturn(0);
9174cdd57e5SDebojyoti Ghosh }
9184cdd57e5SDebojyoti Ghosh 
9199371c9d4SSatish Balay PetscErrorCode TSGetTimeError_GLEE(TS ts, PetscInt n, Vec *X) {
9204cdd57e5SDebojyoti Ghosh   TS_GLEE     *glee = (TS_GLEE *)ts->data;
9214cdd57e5SDebojyoti Ghosh   GLEETableau  tab  = glee->tableau;
9224cdd57e5SDebojyoti Ghosh   PetscReal   *F    = tab->Ferror;
9234cdd57e5SDebojyoti Ghosh   PetscInt     r    = tab->r;
9244cdd57e5SDebojyoti Ghosh   Vec         *Y    = glee->Y;
925ec0e3ee2SEmil Constantinescu   PetscScalar *wr   = glee->rwork;
926ec0e3ee2SEmil Constantinescu   PetscInt     i;
9274cdd57e5SDebojyoti Ghosh 
9284cdd57e5SDebojyoti Ghosh   PetscFunctionBegin;
9299566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(*X));
930657c1e31SEmil Constantinescu   if (n == 0) {
931ec0e3ee2SEmil Constantinescu     for (i = 0; i < r; i++) wr[i] = F[i];
9329566063dSJacob Faibussowitsch     PetscCall(VecMAXPY((*X), r, wr, Y));
9330a01e1b2SEmil Constantinescu   } else if (n == -1) {
934657c1e31SEmil Constantinescu     *X = glee->yGErr;
9350a01e1b2SEmil Constantinescu   }
9364cdd57e5SDebojyoti Ghosh   PetscFunctionReturn(0);
9374cdd57e5SDebojyoti Ghosh }
9384cdd57e5SDebojyoti Ghosh 
9399371c9d4SSatish Balay PetscErrorCode TSSetTimeError_GLEE(TS ts, Vec X) {
94057df6a1bSDebojyoti Ghosh   TS_GLEE    *glee = (TS_GLEE *)ts->data;
94157df6a1bSDebojyoti Ghosh   GLEETableau tab  = glee->tableau;
94257df6a1bSDebojyoti Ghosh   PetscReal  *S    = tab->Serror;
94357df6a1bSDebojyoti Ghosh   PetscInt    r    = tab->r, i;
94457df6a1bSDebojyoti Ghosh   Vec        *Y    = glee->Y;
94557df6a1bSDebojyoti Ghosh 
94657df6a1bSDebojyoti Ghosh   PetscFunctionBegin;
94763a3b9bcSJacob Faibussowitsch   PetscCheck(r == 2, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSSetTimeError_GLEE not supported for '%s' with r=%" PetscInt_FMT ".", tab->name, tab->r);
94857df6a1bSDebojyoti Ghosh   for (i = 1; i < r; i++) {
9499566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vec_sol, Y[i]));
9509566063dSJacob Faibussowitsch     PetscCall(VecAXPBY(Y[i], S[0], S[1], X));
9519566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, glee->yGErr));
95257df6a1bSDebojyoti Ghosh   }
95357df6a1bSDebojyoti Ghosh   PetscFunctionReturn(0);
95457df6a1bSDebojyoti Ghosh }
95557df6a1bSDebojyoti Ghosh 
9569371c9d4SSatish Balay static PetscErrorCode TSDestroy_GLEE(TS ts) {
957b3a6b972SJed Brown   PetscFunctionBegin;
9589566063dSJacob Faibussowitsch   PetscCall(TSReset_GLEE(ts));
959b3a6b972SJed Brown   if (ts->dm) {
9609566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSGLEE, DMRestrictHook_TSGLEE, ts));
9619566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSGLEE, DMSubDomainRestrictHook_TSGLEE, ts));
962b3a6b972SJed Brown   }
9639566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
9649566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLEEGetType_C", NULL));
9659566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLEESetType_C", NULL));
966b3a6b972SJed Brown   PetscFunctionReturn(0);
967b3a6b972SJed Brown }
968b3a6b972SJed Brown 
9692302aa1bSDebojyoti Ghosh /* ------------------------------------------------------------ */
9702302aa1bSDebojyoti Ghosh /*MC
9712302aa1bSDebojyoti Ghosh       TSGLEE - ODE and DAE solver using General Linear with Error Estimation schemes
9722302aa1bSDebojyoti Ghosh 
9732302aa1bSDebojyoti Ghosh   The user should provide the right hand side of the equation
9742302aa1bSDebojyoti Ghosh   using TSSetRHSFunction().
9752302aa1bSDebojyoti Ghosh 
9762302aa1bSDebojyoti Ghosh   Notes:
9772302aa1bSDebojyoti Ghosh   The default is TSGLEE35, it can be changed with TSGLEESetType() or -ts_glee_type
9782302aa1bSDebojyoti Ghosh 
9792302aa1bSDebojyoti Ghosh   Level: beginner
9802302aa1bSDebojyoti Ghosh 
981db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSGLEESetType()`, `TSGLEEGetType()`,
982db781477SPatrick Sanan           `TSGLEE23`, `TTSGLEE24`, `TSGLEE35`, `TSGLEE25I`, `TSGLEEEXRK2A`,
983db781477SPatrick Sanan           `TSGLEERK32G1`, `TSGLEERK285EX`, `TSGLEEType`, `TSGLEERegister()`
9842302aa1bSDebojyoti Ghosh 
9852302aa1bSDebojyoti Ghosh M*/
9869371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TSCreate_GLEE(TS ts) {
9872302aa1bSDebojyoti Ghosh   TS_GLEE *th;
9882302aa1bSDebojyoti Ghosh 
9892302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
9909566063dSJacob Faibussowitsch   PetscCall(TSGLEEInitializePackage());
9912302aa1bSDebojyoti Ghosh 
9922302aa1bSDebojyoti Ghosh   ts->ops->reset                 = TSReset_GLEE;
9932302aa1bSDebojyoti Ghosh   ts->ops->destroy               = TSDestroy_GLEE;
9942302aa1bSDebojyoti Ghosh   ts->ops->view                  = TSView_GLEE;
9952302aa1bSDebojyoti Ghosh   ts->ops->load                  = TSLoad_GLEE;
9962302aa1bSDebojyoti Ghosh   ts->ops->setup                 = TSSetUp_GLEE;
9972302aa1bSDebojyoti Ghosh   ts->ops->step                  = TSStep_GLEE;
9982302aa1bSDebojyoti Ghosh   ts->ops->interpolate           = TSInterpolate_GLEE;
9992302aa1bSDebojyoti Ghosh   ts->ops->evaluatestep          = TSEvaluateStep_GLEE;
10002302aa1bSDebojyoti Ghosh   ts->ops->setfromoptions        = TSSetFromOptions_GLEE;
10012302aa1bSDebojyoti Ghosh   ts->ops->getstages             = TSGetStages_GLEE;
10022302aa1bSDebojyoti Ghosh   ts->ops->snesfunction          = SNESTSFormFunction_GLEE;
10032302aa1bSDebojyoti Ghosh   ts->ops->snesjacobian          = SNESTSFormJacobian_GLEE;
1004b2bf4f3aSDebojyoti Ghosh   ts->ops->getsolutioncomponents = TSGetSolutionComponents_GLEE;
10054cdd57e5SDebojyoti Ghosh   ts->ops->getauxsolution        = TSGetAuxSolution_GLEE;
10064cdd57e5SDebojyoti Ghosh   ts->ops->gettimeerror          = TSGetTimeError_GLEE;
100757df6a1bSDebojyoti Ghosh   ts->ops->settimeerror          = TSSetTimeError_GLEE;
10082302aa1bSDebojyoti Ghosh   ts->ops->startingmethod        = TSStartingMethod_GLEE;
1009b80d5313SLisandro Dalcin   ts->default_adapt_type         = TSADAPTGLEE;
10102302aa1bSDebojyoti Ghosh 
1011efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
1012efd4aadfSBarry Smith 
1013*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&th));
10142302aa1bSDebojyoti Ghosh   ts->data = (void *)th;
10152302aa1bSDebojyoti Ghosh 
10169566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLEEGetType_C", TSGLEEGetType_GLEE));
10179566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLEESetType_C", TSGLEESetType_GLEE));
10182302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
10192302aa1bSDebojyoti Ghosh }
1020