xref: /petsc/src/ts/impls/glee/glee.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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 */
53dd8e379bSPierre Jolivet   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
65efa39862SBarry Smith      TSGLEEi1 - Second order three stage implicit GLEE method
66efa39862SBarry Smith 
67efa39862SBarry Smith      This method has two stages.
68efa39862SBarry Smith      s = 3, r = 2
69efa39862SBarry Smith 
70efa39862SBarry Smith      Level: advanced
71efa39862SBarry Smith 
72efa39862SBarry Smith .seealso: [](ch_ts), `TSGLEE`
73efa39862SBarry Smith M*/
74efa39862SBarry Smith /*MC
75efa39862SBarry Smith      TSGLEE23 - Second order three stage explicit GLEE method
762302aa1bSDebojyoti Ghosh 
772302aa1bSDebojyoti Ghosh      This method has three stages.
782302aa1bSDebojyoti Ghosh      s = 3, r = 2
792302aa1bSDebojyoti Ghosh 
802302aa1bSDebojyoti Ghosh      Level: advanced
812302aa1bSDebojyoti Ghosh 
821cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEE`
8336c34d14SEmil Constantinescu M*/
842302aa1bSDebojyoti Ghosh /*MC
85efa39862SBarry Smith      TSGLEE24 - Second order four stage explicit GLEE method
862302aa1bSDebojyoti Ghosh 
872302aa1bSDebojyoti Ghosh      This method has four stages.
882302aa1bSDebojyoti Ghosh      s = 4, r = 2
892302aa1bSDebojyoti Ghosh 
902302aa1bSDebojyoti Ghosh      Level: advanced
912302aa1bSDebojyoti Ghosh 
921cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEE`
9336c34d14SEmil Constantinescu M*/
942302aa1bSDebojyoti Ghosh /*MC
953297b539SPierre Jolivet      TSGLEE25i - Second order five stage explicit 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 
1021cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEE`
10336c34d14SEmil Constantinescu M*/
1042302aa1bSDebojyoti Ghosh /*MC
105efa39862SBarry Smith      TSGLEE35  - Third order five stage explicit GLEE method
1062302aa1bSDebojyoti Ghosh 
1072302aa1bSDebojyoti Ghosh      This method has five stages.
1082302aa1bSDebojyoti Ghosh      s = 5, r = 2
1092302aa1bSDebojyoti Ghosh 
1102302aa1bSDebojyoti Ghosh      Level: advanced
1112302aa1bSDebojyoti Ghosh 
1121cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEE`
11336c34d14SEmil Constantinescu M*/
1142302aa1bSDebojyoti Ghosh /*MC
1153297b539SPierre Jolivet      TSGLEEEXRK2A  - Second order six stage explicit GLEE method
1162302aa1bSDebojyoti Ghosh 
1172302aa1bSDebojyoti Ghosh      This method has six stages.
1182302aa1bSDebojyoti Ghosh      s = 6, r = 2
1192302aa1bSDebojyoti Ghosh 
1202302aa1bSDebojyoti Ghosh      Level: advanced
1212302aa1bSDebojyoti Ghosh 
1221cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEE`
12336c34d14SEmil Constantinescu M*/
1242302aa1bSDebojyoti Ghosh /*MC
125efa39862SBarry Smith      TSGLEERK32G1  - Third order eight stage explicit GLEE method
1262302aa1bSDebojyoti Ghosh 
1272302aa1bSDebojyoti Ghosh      This method has eight stages.
1282302aa1bSDebojyoti Ghosh      s = 8, r = 2
1292302aa1bSDebojyoti Ghosh 
1302302aa1bSDebojyoti Ghosh      Level: advanced
1312302aa1bSDebojyoti Ghosh 
1321cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEE`
13336c34d14SEmil Constantinescu M*/
1342302aa1bSDebojyoti Ghosh /*MC
135efa39862SBarry Smith      TSGLEERK285EX  - Second order nine stage explicit GLEE method
1362302aa1bSDebojyoti Ghosh 
1372302aa1bSDebojyoti Ghosh      This method has nine stages.
1382302aa1bSDebojyoti Ghosh      s = 9, r = 2
1392302aa1bSDebojyoti Ghosh 
1402302aa1bSDebojyoti Ghosh      Level: advanced
1412302aa1bSDebojyoti Ghosh 
1421cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEE`
14336c34d14SEmil Constantinescu M*/
1442302aa1bSDebojyoti Ghosh 
1452302aa1bSDebojyoti Ghosh /*@C
146bcf0153eSBarry Smith   TSGLEERegisterAll - Registers all of the General Linear with Error Estimation methods in `TSGLEE`
1472302aa1bSDebojyoti Ghosh 
1482302aa1bSDebojyoti Ghosh   Not Collective, but should be called by all processes which will need the schemes to be registered
1492302aa1bSDebojyoti Ghosh 
1502302aa1bSDebojyoti Ghosh   Level: advanced
1512302aa1bSDebojyoti Ghosh 
1521cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEERegisterDestroy()`
1532302aa1bSDebojyoti Ghosh @*/
TSGLEERegisterAll(void)154d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLEERegisterAll(void)
155d71ae5a4SJacob Faibussowitsch {
1562302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
1573ba16761SJacob Faibussowitsch   if (TSGLEERegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
1582302aa1bSDebojyoti Ghosh   TSGLEERegisterAllCalled = PETSC_TRUE;
1592302aa1bSDebojyoti Ghosh 
1602302aa1bSDebojyoti Ghosh   {
161b1fbfd33SSatish Balay #define GAMMA 0.5
1622302aa1bSDebojyoti Ghosh     /* y-eps form */
1639371c9d4SSatish Balay     const PetscInt  p = 1, s = 3, r = 2;
1649371c9d4SSatish Balay     const PetscReal A[3][3] =
1659371c9d4SSatish Balay       {
1669371c9d4SSatish Balay         {1.0, 0,   0  },
1679371c9d4SSatish Balay         {0,   0.5, 0  },
1689371c9d4SSatish Balay         {0,   0.5, 0.5}
1699371c9d4SSatish Balay     },
1709371c9d4SSatish 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};
1719566063dSJacob 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));
172ab8c5912SEmil Constantinescu   }
173ab8c5912SEmil Constantinescu   {
174b1fbfd33SSatish Balay #undef GAMMA
175b1fbfd33SSatish Balay #define GAMMA 0.0
176ab8c5912SEmil Constantinescu     /* y-eps form */
1779371c9d4SSatish Balay     const PetscInt  p = 2, s = 3, r = 2;
1789371c9d4SSatish Balay     const PetscReal A[3][3] =
1799371c9d4SSatish Balay       {
1809371c9d4SSatish Balay         {0,    0,    0},
1819371c9d4SSatish Balay         {1,    0,    0},
1829371c9d4SSatish Balay         {0.25, 0.25, 0}
1839371c9d4SSatish Balay     },
1849371c9d4SSatish 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};
1859566063dSJacob 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));
1862302aa1bSDebojyoti Ghosh   }
1872302aa1bSDebojyoti Ghosh   {
188b1fbfd33SSatish Balay #undef GAMMA
189b1fbfd33SSatish Balay #define GAMMA 0.0
1902302aa1bSDebojyoti Ghosh     /* y-y~ form */
1919371c9d4SSatish Balay     const PetscInt  p = 2, s = 4, r = 2;
1929371c9d4SSatish Balay     const PetscReal A[4][4] =
1939371c9d4SSatish Balay       {
1949371c9d4SSatish Balay         {0,            0,            0,            0},
1959371c9d4SSatish Balay         {0.75,         0,            0,            0},
1969371c9d4SSatish Balay         {0.25,         29.0 / 60.0,  0,            0},
1979371c9d4SSatish Balay         {-21.0 / 44.0, 145.0 / 44.0, -20.0 / 11.0, 0}
1989371c9d4SSatish Balay     },
1999371c9d4SSatish 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};
2009566063dSJacob 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));
2012302aa1bSDebojyoti Ghosh   }
2022302aa1bSDebojyoti Ghosh   {
203b1fbfd33SSatish Balay #undef GAMMA
204b1fbfd33SSatish Balay #define GAMMA 0.0
2052302aa1bSDebojyoti Ghosh     /* y-y~ form */
2069371c9d4SSatish Balay     const PetscInt  p = 2, s = 5, r = 2;
2079371c9d4SSatish Balay     const PetscReal A[5][5] =
2089371c9d4SSatish Balay       {
2099371c9d4SSatish Balay         {0,                       0,                       0,                       0,                      0},
2102302aa1bSDebojyoti Ghosh         {-0.94079244066783383269, 0,                       0,                       0,                      0},
2112302aa1bSDebojyoti Ghosh         {0.64228187778301907108,  0.10915356933958500042,  0,                       0,                      0},
2122302aa1bSDebojyoti Ghosh         {-0.51764297742287450812, 0.74414270351096040738,  -0.71404164927824538121, 0,                      0},
2139371c9d4SSatish Balay         {-0.44696561556825969206, -0.76768425657590196518, 0.20111608138142987881,  0.93828186737840469796, 0}
2149371c9d4SSatish Balay     },
2159371c9d4SSatish 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};
2169566063dSJacob 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));
2172302aa1bSDebojyoti Ghosh   }
2182302aa1bSDebojyoti Ghosh   {
219b1fbfd33SSatish Balay #undef GAMMA
220b1fbfd33SSatish Balay #define GAMMA 0.0
2212302aa1bSDebojyoti Ghosh     /* y-y~ form */
2229371c9d4SSatish Balay     const PetscInt  p = 3, s = 5, r = 2;
2239371c9d4SSatish Balay     const PetscReal A[5][5] =
2249371c9d4SSatish Balay       {
2259371c9d4SSatish Balay         {0,                                                0,                                                 0,                                                 0,                                               0},
2262302aa1bSDebojyoti Ghosh         {-2169604947363702313.0 / 24313474998937147335.0,  0,                                                 0,                                                 0,                                               0},
2272302aa1bSDebojyoti Ghosh         {46526746497697123895.0 / 94116917485856474137.0,  -10297879244026594958.0 / 49199457603717988219.0,  0,                                                 0,                                               0},
2282302aa1bSDebojyoti Ghosh         {23364788935845982499.0 / 87425311444725389446.0,  -79205144337496116638.0 / 148994349441340815519.0, 40051189859317443782.0 / 36487615018004984309.0,   0,                                               0},
2299371c9d4SSatish Balay         {42089522664062539205.0 / 124911313006412840286.0, -15074384760342762939.0 / 137927286865289746282.0, -62274678522253371016.0 / 125918573676298591413.0, 13755475729852471739.0 / 79257927066651693390.0, 0}
2309371c9d4SSatish Balay     },
2319371c9d4SSatish 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};
2329566063dSJacob 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));
2332302aa1bSDebojyoti Ghosh   }
2342302aa1bSDebojyoti Ghosh   {
235b1fbfd33SSatish Balay #undef GAMMA
236b1fbfd33SSatish Balay #define GAMMA 0.25
2372302aa1bSDebojyoti Ghosh     /* y-eps form */
2389371c9d4SSatish Balay     const PetscInt  p = 2, s = 6, r = 2;
2399371c9d4SSatish Balay     const PetscReal A[6][6] =
2409371c9d4SSatish Balay       {
2419371c9d4SSatish Balay         {0, 0, 0,    0,    0,   0},
2422302aa1bSDebojyoti Ghosh         {1, 0, 0,    0,    0,   0},
2432302aa1bSDebojyoti Ghosh         {0, 0, 0,    0,    0,   0},
2442302aa1bSDebojyoti Ghosh         {0, 0, 0.5,  0,    0,   0},
2452302aa1bSDebojyoti Ghosh         {0, 0, 0.25, 0.25, 0,   0},
2469371c9d4SSatish Balay         {0, 0, 0.25, 0.25, 0.5, 0}
2479371c9d4SSatish Balay     },
2489371c9d4SSatish 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};
2499566063dSJacob 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));
2502302aa1bSDebojyoti Ghosh   }
2512302aa1bSDebojyoti Ghosh   {
252b1fbfd33SSatish Balay #undef GAMMA
253b1fbfd33SSatish Balay #define GAMMA 0.0
2542302aa1bSDebojyoti Ghosh     /* y-eps form */
2559371c9d4SSatish Balay     const PetscInt  p = 3, s = 8, r = 2;
2569371c9d4SSatish Balay     const PetscReal A[8][8] =
2579371c9d4SSatish Balay       {
2589371c9d4SSatish Balay         {0,           0,          0,          0,          0,         0,         0,         0},
2592302aa1bSDebojyoti Ghosh         {0.5,         0,          0,          0,          0,         0,         0,         0},
2602302aa1bSDebojyoti Ghosh         {-1,          2,          0,          0,          0,         0,         0,         0},
2612302aa1bSDebojyoti Ghosh         {1.0 / 6.0,   2.0 / 3.0,  1.0 / 6.0,  0,          0,         0,         0,         0},
2622302aa1bSDebojyoti Ghosh         {0,           0,          0,          0,          0,         0,         0,         0},
2632302aa1bSDebojyoti Ghosh         {-7.0 / 24.0, 1.0 / 3.0,  1.0 / 12.0, -1.0 / 8.0, 0.5,       0,         0,         0},
2642302aa1bSDebojyoti Ghosh         {7.0 / 6.0,   -4.0 / 3.0, -1.0 / 3.0, 0.5,        -1.0,      2.0,       0,         0},
2659371c9d4SSatish Balay         {0,           0,          0,          0,          1.0 / 6.0, 2.0 / 3.0, 1.0 / 6.0, 0}
2669371c9d4SSatish Balay     },
2679371c9d4SSatish 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};
2689566063dSJacob 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));
2692302aa1bSDebojyoti Ghosh   }
2702302aa1bSDebojyoti Ghosh   {
271b1fbfd33SSatish Balay #undef GAMMA
272b1fbfd33SSatish Balay #define GAMMA 0.25
2732302aa1bSDebojyoti Ghosh     /* y-eps form */
2749371c9d4SSatish Balay     const PetscInt  p = 2, s = 9, r = 2;
2759371c9d4SSatish Balay     const PetscReal A[9][9] =
2769371c9d4SSatish Balay       {
2779371c9d4SSatish Balay         {0,                    0,                    0, 0,                     0,                    0,                    0,                     0,                    0},
2782302aa1bSDebojyoti Ghosh         {0.585786437626904966, 0,                    0, 0,                     0,                    0,                    0,                     0,                    0},
2792302aa1bSDebojyoti Ghosh         {0.149999999999999994, 0.849999999999999978, 0, 0,                     0,                    0,                    0,                     0,                    0},
2802302aa1bSDebojyoti Ghosh         {0,                    0,                    0, 0,                     0,                    0,                    0,                     0,                    0},
2812302aa1bSDebojyoti Ghosh         {0,                    0,                    0, 0.292893218813452483,  0,                    0,                    0,                     0,                    0},
2822302aa1bSDebojyoti Ghosh         {0,                    0,                    0, 0.0749999999999999972, 0.424999999999999989, 0,                    0,                     0,                    0},
2832302aa1bSDebojyoti Ghosh         {0,                    0,                    0, 0.176776695296636893,  0.176776695296636893, 0.146446609406726241, 0,                     0,                    0},
2842302aa1bSDebojyoti Ghosh         {0,                    0,                    0, 0.176776695296636893,  0.176776695296636893, 0.146446609406726241, 0.292893218813452483,  0,                    0},
2859371c9d4SSatish Balay         {0,                    0,                    0, 0.176776695296636893,  0.176776695296636893, 0.146446609406726241, 0.0749999999999999972, 0.424999999999999989, 0}
2869371c9d4SSatish Balay     },
2879371c9d4SSatish 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};
2889566063dSJacob 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));
2892302aa1bSDebojyoti Ghosh   }
2903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2912302aa1bSDebojyoti Ghosh }
2922302aa1bSDebojyoti Ghosh 
2932302aa1bSDebojyoti Ghosh /*@C
294bcf0153eSBarry Smith   TSGLEERegisterDestroy - Frees the list of schemes that were registered by `TSGLEERegister()`.
2952302aa1bSDebojyoti Ghosh 
2962302aa1bSDebojyoti Ghosh   Not Collective
2972302aa1bSDebojyoti Ghosh 
2982302aa1bSDebojyoti Ghosh   Level: advanced
2992302aa1bSDebojyoti Ghosh 
3001cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEERegister()`, `TSGLEERegisterAll()`
3012302aa1bSDebojyoti Ghosh @*/
TSGLEERegisterDestroy(void)302d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLEERegisterDestroy(void)
303d71ae5a4SJacob Faibussowitsch {
3042302aa1bSDebojyoti Ghosh   GLEETableauLink link;
3052302aa1bSDebojyoti Ghosh 
3062302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
3072302aa1bSDebojyoti Ghosh   while ((link = GLEETableauList)) {
3082302aa1bSDebojyoti Ghosh     GLEETableau t   = &link->tab;
3092302aa1bSDebojyoti Ghosh     GLEETableauList = link->next;
3109566063dSJacob Faibussowitsch     PetscCall(PetscFree5(t->A, t->B, t->U, t->V, t->c));
311d0609cedSBarry Smith     PetscCall(PetscFree2(t->S, t->F));
312d0609cedSBarry Smith     PetscCall(PetscFree(t->Fembed));
313d0609cedSBarry Smith     PetscCall(PetscFree(t->Ferror));
314d0609cedSBarry Smith     PetscCall(PetscFree(t->Serror));
315d0609cedSBarry Smith     PetscCall(PetscFree(t->binterp));
316d0609cedSBarry Smith     PetscCall(PetscFree(t->name));
317d0609cedSBarry Smith     PetscCall(PetscFree(link));
3182302aa1bSDebojyoti Ghosh   }
3192302aa1bSDebojyoti Ghosh   TSGLEERegisterAllCalled = PETSC_FALSE;
3203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3212302aa1bSDebojyoti Ghosh }
3222302aa1bSDebojyoti Ghosh 
3232302aa1bSDebojyoti Ghosh /*@C
324bcf0153eSBarry Smith   TSGLEEInitializePackage - This function initializes everything in the `TSGLEE` package. It is called
325bcf0153eSBarry Smith   from `TSInitializePackage()`.
3262302aa1bSDebojyoti Ghosh 
3272302aa1bSDebojyoti Ghosh   Level: developer
3282302aa1bSDebojyoti Ghosh 
3291cc06b55SBarry Smith .seealso: [](ch_ts), `PetscInitialize()`
3302302aa1bSDebojyoti Ghosh @*/
TSGLEEInitializePackage(void)331d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLEEInitializePackage(void)
332d71ae5a4SJacob Faibussowitsch {
3332302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
3343ba16761SJacob Faibussowitsch   if (TSGLEEPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
3352302aa1bSDebojyoti Ghosh   TSGLEEPackageInitialized = PETSC_TRUE;
3369566063dSJacob Faibussowitsch   PetscCall(TSGLEERegisterAll());
3379566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposedDataRegister(&explicit_stage_time_id));
3389566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(TSGLEEFinalizePackage));
3393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3402302aa1bSDebojyoti Ghosh }
3412302aa1bSDebojyoti Ghosh 
3422302aa1bSDebojyoti Ghosh /*@C
343bcf0153eSBarry Smith   TSGLEEFinalizePackage - This function destroys everything in the `TSGLEE` package. It is
344bcf0153eSBarry Smith   called from `PetscFinalize()`.
3452302aa1bSDebojyoti Ghosh 
3462302aa1bSDebojyoti Ghosh   Level: developer
3472302aa1bSDebojyoti Ghosh 
3481cc06b55SBarry Smith .seealso: [](ch_ts), `PetscFinalize()`
3492302aa1bSDebojyoti Ghosh @*/
TSGLEEFinalizePackage(void)350d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLEEFinalizePackage(void)
351d71ae5a4SJacob Faibussowitsch {
3522302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
3532302aa1bSDebojyoti Ghosh   TSGLEEPackageInitialized = PETSC_FALSE;
3549566063dSJacob Faibussowitsch   PetscCall(TSGLEERegisterDestroy());
3553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3562302aa1bSDebojyoti Ghosh }
3572302aa1bSDebojyoti Ghosh 
3582302aa1bSDebojyoti Ghosh /*@C
359bcf0153eSBarry Smith   TSGLEERegister - register a new `TSGLEE` scheme by providing the entries in the Butcher tableau
3602302aa1bSDebojyoti Ghosh 
361cc4c1da9SBarry Smith   Not Collective, but the same schemes should be registered on all processes on which they will be used, No Fortran Support
3622302aa1bSDebojyoti Ghosh 
3632302aa1bSDebojyoti Ghosh   Input Parameters:
3642302aa1bSDebojyoti Ghosh + name    - identifier for method
3652302aa1bSDebojyoti Ghosh . order   - order of method
3662302aa1bSDebojyoti Ghosh . s       - number of stages
3672302aa1bSDebojyoti Ghosh . r       - number of steps
3682302aa1bSDebojyoti Ghosh . gamma   - LTE ratio
3692302aa1bSDebojyoti Ghosh . A       - stage coefficients (dimension s*s, row-major)
3702302aa1bSDebojyoti Ghosh . B       - step completion coefficients (dimension r*s, row-major)
3712302aa1bSDebojyoti Ghosh . U       - method coefficients (dimension s*r, row-major)
3722302aa1bSDebojyoti Ghosh . V       - method coefficients (dimension r*r, row-major)
3732302aa1bSDebojyoti Ghosh . S       - starting coefficients
3742302aa1bSDebojyoti Ghosh . F       - finishing coefficients
3752302aa1bSDebojyoti Ghosh . c       - abscissa (dimension s; NULL to use row sums of A)
3762302aa1bSDebojyoti Ghosh . Fembed  - step completion coefficients for embedded method
377fd96ebc2SDebojyoti Ghosh . Ferror  - error computation coefficients
37857df6a1bSDebojyoti Ghosh . Serror  - error initialization coefficients
3792302aa1bSDebojyoti Ghosh . pinterp - order of interpolation (0 if unavailable)
3802302aa1bSDebojyoti Ghosh - binterp - array of interpolation coefficients (NULL if unavailable)
3812302aa1bSDebojyoti Ghosh 
3822302aa1bSDebojyoti Ghosh   Level: advanced
3832302aa1bSDebojyoti Ghosh 
384bcf0153eSBarry Smith   Note:
385bcf0153eSBarry Smith   Several `TSGLEE` methods are provided, this function is only needed to create new methods.
386bcf0153eSBarry Smith 
3871cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEE`
3882302aa1bSDebojyoti Ghosh @*/
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[])389d71ae5a4SJacob Faibussowitsch 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[])
390d71ae5a4SJacob Faibussowitsch {
3912302aa1bSDebojyoti Ghosh   GLEETableauLink link;
3922302aa1bSDebojyoti Ghosh   GLEETableau     t;
3932302aa1bSDebojyoti Ghosh   PetscInt        i, j;
3942302aa1bSDebojyoti Ghosh 
3952302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
3969566063dSJacob Faibussowitsch   PetscCall(TSGLEEInitializePackage());
3979566063dSJacob Faibussowitsch   PetscCall(PetscNew(&link));
3982302aa1bSDebojyoti Ghosh   t = &link->tab;
3999566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &t->name));
4002302aa1bSDebojyoti Ghosh   t->order = order;
4012302aa1bSDebojyoti Ghosh   t->s     = s;
4022302aa1bSDebojyoti Ghosh   t->r     = r;
4032302aa1bSDebojyoti Ghosh   t->gamma = gamma;
4049566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(s * s, &t->A, r * r, &t->V, s, &t->c, r * s, &t->B, s * r, &t->U));
4059566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(r, &t->S, r, &t->F));
4069566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->A, A, s * s));
4079566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->B, B, r * s));
4089566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->U, U, s * r));
4099566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->V, V, r * r));
4109566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->S, S, r));
4119566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->F, F, r));
4122302aa1bSDebojyoti Ghosh   if (c) {
4139566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(t->c, c, s));
4142302aa1bSDebojyoti Ghosh   } else {
4159371c9d4SSatish Balay     for (i = 0; i < s; i++)
4169371c9d4SSatish Balay       for (j = 0, t->c[i] = 0; j < s; j++) t->c[i] += A[i * s + j];
4172302aa1bSDebojyoti Ghosh   }
4189566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(r, &t->Fembed));
4199566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(r, &t->Ferror));
4209566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(r, &t->Serror));
4219566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->Fembed, Fembed, r));
4229566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->Ferror, Ferror, r));
4239566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->Serror, Serror, r));
4242302aa1bSDebojyoti Ghosh   t->pinterp = pinterp;
4259566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s * pinterp, &t->binterp));
4269566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->binterp, binterp, s * pinterp));
4272302aa1bSDebojyoti Ghosh 
4282302aa1bSDebojyoti Ghosh   link->next      = GLEETableauList;
4292302aa1bSDebojyoti Ghosh   GLEETableauList = link;
4303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4312302aa1bSDebojyoti Ghosh }
4322302aa1bSDebojyoti Ghosh 
TSEvaluateStep_GLEE(TS ts,PetscInt order,Vec X,PetscBool * done)433d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateStep_GLEE(TS ts, PetscInt order, Vec X, PetscBool *done)
434d71ae5a4SJacob Faibussowitsch {
4352302aa1bSDebojyoti Ghosh   TS_GLEE     *glee = (TS_GLEE *)ts->data;
4362302aa1bSDebojyoti Ghosh   GLEETableau  tab  = glee->tableau;
4379371c9d4SSatish Balay   PetscReal    h, *B = tab->B, *V = tab->V, *F = tab->F, *Fembed = tab->Fembed;
4382302aa1bSDebojyoti Ghosh   PetscInt     s = tab->s, r = tab->r, i, j;
4392302aa1bSDebojyoti Ghosh   Vec         *Y = glee->Y, *YdotStage = glee->YdotStage;
440ec0e3ee2SEmil Constantinescu   PetscScalar *ws = glee->swork, *wr = glee->rwork;
4412302aa1bSDebojyoti Ghosh 
4422302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
4432302aa1bSDebojyoti Ghosh   switch (glee->status) {
4442302aa1bSDebojyoti Ghosh   case TS_STEP_INCOMPLETE:
445d71ae5a4SJacob Faibussowitsch   case TS_STEP_PENDING:
446d71ae5a4SJacob Faibussowitsch     h = ts->time_step;
447d71ae5a4SJacob Faibussowitsch     break;
448d71ae5a4SJacob Faibussowitsch   case TS_STEP_COMPLETE:
449d71ae5a4SJacob Faibussowitsch     h = ts->ptime - ts->ptime_prev;
450d71ae5a4SJacob Faibussowitsch     break;
451d71ae5a4SJacob Faibussowitsch   default:
452d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus");
4532302aa1bSDebojyoti Ghosh   }
4542302aa1bSDebojyoti Ghosh 
4552302aa1bSDebojyoti Ghosh   if (order == tab->order) {
4562302aa1bSDebojyoti Ghosh     /* Note: Irrespective of whether status is TS_STEP_INCOMPLETE
457fd96ebc2SDebojyoti Ghosh              or TS_STEP_COMPLETE, glee->X has the solution at the
4582302aa1bSDebojyoti Ghosh              beginning of the time step. So no need to roll-back.
4592302aa1bSDebojyoti Ghosh     */
4602302aa1bSDebojyoti Ghosh     if (glee->status == TS_STEP_INCOMPLETE) {
4612302aa1bSDebojyoti Ghosh       for (i = 0; i < r; i++) {
4629566063dSJacob Faibussowitsch         PetscCall(VecZeroEntries(Y[i]));
463ec0e3ee2SEmil Constantinescu         for (j = 0; j < r; j++) wr[j] = V[i * r + j];
4649566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(Y[i], r, wr, glee->X));
465ec0e3ee2SEmil Constantinescu         for (j = 0; j < s; j++) ws[j] = h * B[i * s + j];
4669566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(Y[i], s, ws, YdotStage));
4672302aa1bSDebojyoti Ghosh       }
4689566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(X));
469ec0e3ee2SEmil Constantinescu       for (j = 0; j < r; j++) wr[j] = F[j];
4709566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(X, r, wr, Y));
4719566063dSJacob Faibussowitsch     } else PetscCall(VecCopy(ts->vec_sol, X));
4723ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
4732302aa1bSDebojyoti Ghosh 
4742302aa1bSDebojyoti Ghosh   } else if (order == tab->order - 1) {
4752302aa1bSDebojyoti Ghosh     /* Complete with the embedded method (Fembed) */
4762302aa1bSDebojyoti Ghosh     for (i = 0; i < r; i++) {
4779566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(Y[i]));
478ec0e3ee2SEmil Constantinescu       for (j = 0; j < r; j++) wr[j] = V[i * r + j];
4799566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(Y[i], r, wr, glee->X));
480ec0e3ee2SEmil Constantinescu       for (j = 0; j < s; j++) ws[j] = h * B[i * s + j];
4819566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(Y[i], s, ws, YdotStage));
4822302aa1bSDebojyoti Ghosh     }
4839566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(X));
484ec0e3ee2SEmil Constantinescu     for (j = 0; j < r; j++) wr[j] = Fembed[j];
4859566063dSJacob Faibussowitsch     PetscCall(VecMAXPY(X, r, wr, Y));
486fd96ebc2SDebojyoti Ghosh 
4872302aa1bSDebojyoti Ghosh     if (done) *done = PETSC_TRUE;
4883ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
4892302aa1bSDebojyoti Ghosh   }
490966bd95aSPierre Jolivet   PetscCheck(done, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "GLEE '%s' of order %" PetscInt_FMT " cannot evaluate step at order %" PetscInt_FMT, tab->name, tab->order, order);
491966bd95aSPierre Jolivet   *done = PETSC_FALSE;
4923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4932302aa1bSDebojyoti Ghosh }
4942302aa1bSDebojyoti Ghosh 
TSStep_GLEE(TS ts)495d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_GLEE(TS ts)
496d71ae5a4SJacob Faibussowitsch {
4972302aa1bSDebojyoti Ghosh   TS_GLEE       *glee = (TS_GLEE *)ts->data;
4982302aa1bSDebojyoti Ghosh   GLEETableau    tab  = glee->tableau;
4992302aa1bSDebojyoti Ghosh   const PetscInt s = tab->s, r = tab->r;
5009371c9d4SSatish Balay   PetscReal     *A = tab->A, *U = tab->U, *F = tab->F, *c = tab->c;
5019371c9d4SSatish Balay   Vec           *Y = glee->Y, *X = glee->X, *YStage = glee->YStage, *YdotStage = glee->YdotStage, W = glee->W;
5022302aa1bSDebojyoti Ghosh   SNES           snes;
503ec0e3ee2SEmil Constantinescu   PetscScalar   *ws = glee->swork, *wr = glee->rwork;
5042302aa1bSDebojyoti Ghosh   TSAdapt        adapt;
5052302aa1bSDebojyoti Ghosh   PetscInt       i, j, reject, next_scheme, its, lits;
5062302aa1bSDebojyoti Ghosh   PetscReal      next_time_step;
5072302aa1bSDebojyoti Ghosh   PetscReal      t;
5082302aa1bSDebojyoti Ghosh   PetscBool      accept;
5092302aa1bSDebojyoti Ghosh 
5102302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
5119566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(citation, &cited));
512642f3702SEmil Constantinescu 
5139566063dSJacob Faibussowitsch   for (i = 0; i < r; i++) PetscCall(VecCopy(Y[i], X[i]));
5142302aa1bSDebojyoti Ghosh 
5159566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
5162302aa1bSDebojyoti Ghosh   next_time_step = ts->time_step;
5172302aa1bSDebojyoti Ghosh   t              = ts->ptime;
5182302aa1bSDebojyoti Ghosh   accept         = PETSC_TRUE;
5192302aa1bSDebojyoti Ghosh   glee->status   = TS_STEP_INCOMPLETE;
5202302aa1bSDebojyoti Ghosh 
5212302aa1bSDebojyoti Ghosh   for (reject = 0; reject < ts->max_reject && !ts->reason; reject++, ts->reject++) {
5222302aa1bSDebojyoti Ghosh     PetscReal h = ts->time_step;
5239566063dSJacob Faibussowitsch     PetscCall(TSPreStep(ts));
5242302aa1bSDebojyoti Ghosh 
5252302aa1bSDebojyoti Ghosh     for (i = 0; i < s; i++) {
5262302aa1bSDebojyoti Ghosh       glee->stage_time = t + h * c[i];
5279566063dSJacob Faibussowitsch       PetscCall(TSPreStage(ts, glee->stage_time));
5282302aa1bSDebojyoti Ghosh 
5292302aa1bSDebojyoti Ghosh       if (A[i * s + i] == 0) { /* Explicit stage */
5309566063dSJacob Faibussowitsch         PetscCall(VecZeroEntries(YStage[i]));
531ec0e3ee2SEmil Constantinescu         for (j = 0; j < r; j++) wr[j] = U[i * r + j];
5329566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(YStage[i], r, wr, X));
533ec0e3ee2SEmil Constantinescu         for (j = 0; j < i; j++) ws[j] = h * A[i * s + j];
5349566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(YStage[i], i, ws, YdotStage));
5352302aa1bSDebojyoti Ghosh       } else { /* Implicit stage */
5362302aa1bSDebojyoti Ghosh         glee->scoeff = 1.0 / A[i * s + i];
537dd8e379bSPierre Jolivet         /* compute right-hand side */
5389566063dSJacob Faibussowitsch         PetscCall(VecZeroEntries(W));
539ec0e3ee2SEmil Constantinescu         for (j = 0; j < r; j++) wr[j] = U[i * r + j];
5409566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(W, r, wr, X));
541ec0e3ee2SEmil Constantinescu         for (j = 0; j < i; j++) ws[j] = h * A[i * s + j];
5429566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(W, i, ws, YdotStage));
5439566063dSJacob Faibussowitsch         PetscCall(VecScale(W, glee->scoeff / h));
5442302aa1bSDebojyoti Ghosh         /* set initial guess */
5459566063dSJacob Faibussowitsch         PetscCall(VecCopy(i > 0 ? YStage[i - 1] : ts->vec_sol, YStage[i]));
5462302aa1bSDebojyoti Ghosh         /* solve for this stage */
5479566063dSJacob Faibussowitsch         PetscCall(SNESSolve(snes, W, YStage[i]));
5489566063dSJacob Faibussowitsch         PetscCall(SNESGetIterationNumber(snes, &its));
5499566063dSJacob Faibussowitsch         PetscCall(SNESGetLinearSolveIterations(snes, &lits));
5509371c9d4SSatish Balay         ts->snes_its += its;
5519371c9d4SSatish Balay         ts->ksp_its += lits;
5522302aa1bSDebojyoti Ghosh       }
5539566063dSJacob Faibussowitsch       PetscCall(TSGetAdapt(ts, &adapt));
5549566063dSJacob Faibussowitsch       PetscCall(TSAdaptCheckStage(adapt, ts, glee->stage_time, YStage[i], &accept));
5552302aa1bSDebojyoti Ghosh       if (!accept) goto reject_step;
5569566063dSJacob Faibussowitsch       PetscCall(TSPostStage(ts, glee->stage_time, i, YStage));
5579566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSFunction(ts, t + h * c[i], YStage[i], YdotStage[i]));
5582302aa1bSDebojyoti Ghosh     }
5599566063dSJacob Faibussowitsch     PetscCall(TSEvaluateStep(ts, tab->order, ts->vec_sol, NULL));
5602302aa1bSDebojyoti Ghosh     glee->status = TS_STEP_PENDING;
5612302aa1bSDebojyoti Ghosh 
5622302aa1bSDebojyoti Ghosh     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
5639566063dSJacob Faibussowitsch     PetscCall(TSGetAdapt(ts, &adapt));
5649566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidatesClear(adapt));
5659566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidateAdd(adapt, tab->name, tab->order, 1, tab->ccfl, (PetscReal)tab->s, PETSC_TRUE));
5669566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, &next_scheme, &next_time_step, &accept));
5672302aa1bSDebojyoti Ghosh     if (accept) {
5682302aa1bSDebojyoti Ghosh       /* ignore next_scheme for now */
5692302aa1bSDebojyoti Ghosh       ts->ptime += ts->time_step;
5702302aa1bSDebojyoti Ghosh       ts->time_step = next_time_step;
5712302aa1bSDebojyoti Ghosh       glee->status  = TS_STEP_COMPLETE;
5720a01e1b2SEmil Constantinescu       /* compute and store the global error */
5730a01e1b2SEmil Constantinescu       /* Note: this is not needed if TSAdaptGLEE is not used */
574f4f49eeaSPierre Jolivet       PetscCall(TSGetTimeError(ts, 0, &glee->yGErr));
5759566063dSJacob Faibussowitsch       PetscCall(PetscObjectComposedDataSetReal((PetscObject)ts->vec_sol, explicit_stage_time_id, ts->ptime));
5762302aa1bSDebojyoti Ghosh       break;
5772302aa1bSDebojyoti Ghosh     } else { /* Roll back the current step */
578ec0e3ee2SEmil Constantinescu       for (j = 0; j < r; j++) wr[j] = F[j];
5799566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(ts->vec_sol, r, wr, X));
5802302aa1bSDebojyoti Ghosh       ts->time_step = next_time_step;
5812302aa1bSDebojyoti Ghosh       glee->status  = TS_STEP_INCOMPLETE;
5822302aa1bSDebojyoti Ghosh     }
5839371c9d4SSatish Balay   reject_step:
5849371c9d4SSatish Balay     continue;
5852302aa1bSDebojyoti Ghosh   }
5862302aa1bSDebojyoti Ghosh   if (glee->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
5873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5882302aa1bSDebojyoti Ghosh }
5892302aa1bSDebojyoti Ghosh 
TSInterpolate_GLEE(TS ts,PetscReal itime,Vec X)590d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_GLEE(TS ts, PetscReal itime, Vec X)
591d71ae5a4SJacob Faibussowitsch {
5922302aa1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE *)ts->data;
5932302aa1bSDebojyoti Ghosh   PetscInt         s = glee->tableau->s, pinterp = glee->tableau->pinterp, i, j;
5942302aa1bSDebojyoti Ghosh   PetscReal        h, tt, t;
5952302aa1bSDebojyoti Ghosh   PetscScalar     *b;
5962302aa1bSDebojyoti Ghosh   const PetscReal *B = glee->tableau->binterp;
5972302aa1bSDebojyoti Ghosh 
5982302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
5993c633725SBarry Smith   PetscCheck(B, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSGLEE %s does not have an interpolation formula", glee->tableau->name);
6002302aa1bSDebojyoti Ghosh   switch (glee->status) {
6012302aa1bSDebojyoti Ghosh   case TS_STEP_INCOMPLETE:
6022302aa1bSDebojyoti Ghosh   case TS_STEP_PENDING:
6032302aa1bSDebojyoti Ghosh     h = ts->time_step;
6042302aa1bSDebojyoti Ghosh     t = (itime - ts->ptime) / h;
6052302aa1bSDebojyoti Ghosh     break;
6062302aa1bSDebojyoti Ghosh   case TS_STEP_COMPLETE:
607ed450d42SEmil Constantinescu     h = ts->ptime - ts->ptime_prev;
6082302aa1bSDebojyoti Ghosh     t = (itime - ts->ptime) / h + 1; /* In the interval [0,1] */
6092302aa1bSDebojyoti Ghosh     break;
610d71ae5a4SJacob Faibussowitsch   default:
611d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus");
6122302aa1bSDebojyoti Ghosh   }
6139566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s, &b));
6142302aa1bSDebojyoti Ghosh   for (i = 0; i < s; i++) b[i] = 0;
6152302aa1bSDebojyoti Ghosh   for (j = 0, tt = t; j < pinterp; j++, tt *= t) {
616ad540459SPierre Jolivet     for (i = 0; i < s; i++) b[i] += h * B[i * pinterp + j] * tt;
6172302aa1bSDebojyoti Ghosh   }
6189566063dSJacob Faibussowitsch   PetscCall(VecCopy(glee->YStage[0], X));
6199566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(X, s, b, glee->YdotStage));
6209566063dSJacob Faibussowitsch   PetscCall(PetscFree(b));
6213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6222302aa1bSDebojyoti Ghosh }
6232302aa1bSDebojyoti Ghosh 
TSReset_GLEE(TS ts)624d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_GLEE(TS ts)
625d71ae5a4SJacob Faibussowitsch {
6262302aa1bSDebojyoti Ghosh   TS_GLEE *glee = (TS_GLEE *)ts->data;
6272302aa1bSDebojyoti Ghosh   PetscInt s, r;
6282302aa1bSDebojyoti Ghosh 
6292302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
6303ba16761SJacob Faibussowitsch   if (!glee->tableau) PetscFunctionReturn(PETSC_SUCCESS);
6312302aa1bSDebojyoti Ghosh   s = glee->tableau->s;
6322302aa1bSDebojyoti Ghosh   r = glee->tableau->r;
6339566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(r, &glee->Y));
6349566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(r, &glee->X));
6359566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(s, &glee->YStage));
6369566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(s, &glee->YdotStage));
6379566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&glee->Ydot));
6389566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&glee->yGErr));
6399566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&glee->W));
6409566063dSJacob Faibussowitsch   PetscCall(PetscFree2(glee->swork, glee->rwork));
6413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6422302aa1bSDebojyoti Ghosh }
6432302aa1bSDebojyoti Ghosh 
TSGLEEGetVecs(TS ts,DM dm,Vec * Ydot)644d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLEEGetVecs(TS ts, DM dm, Vec *Ydot)
645d71ae5a4SJacob Faibussowitsch {
6462302aa1bSDebojyoti Ghosh   TS_GLEE *glee = (TS_GLEE *)ts->data;
6472302aa1bSDebojyoti Ghosh 
6482302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
6492302aa1bSDebojyoti Ghosh   if (Ydot) {
650ac530a7eSPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMGetNamedGlobalVector(dm, "TSGLEE_Ydot", Ydot));
651ac530a7eSPierre Jolivet     else *Ydot = glee->Ydot;
6522302aa1bSDebojyoti Ghosh   }
6533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6542302aa1bSDebojyoti Ghosh }
6552302aa1bSDebojyoti Ghosh 
TSGLEERestoreVecs(TS ts,DM dm,Vec * Ydot)656d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLEERestoreVecs(TS ts, DM dm, Vec *Ydot)
657d71ae5a4SJacob Faibussowitsch {
6582302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
6592302aa1bSDebojyoti Ghosh   if (Ydot) {
66048a46eb9SPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSGLEE_Ydot", Ydot));
6612302aa1bSDebojyoti Ghosh   }
6623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6632302aa1bSDebojyoti Ghosh }
6642302aa1bSDebojyoti Ghosh 
6652302aa1bSDebojyoti Ghosh /*
6662302aa1bSDebojyoti Ghosh   This defines the nonlinear equation that is to be solved with SNES
6672302aa1bSDebojyoti Ghosh */
SNESTSFormFunction_GLEE(SNES snes,Vec X,Vec F,TS ts)668d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_GLEE(SNES snes, Vec X, Vec F, TS ts)
669d71ae5a4SJacob Faibussowitsch {
6702302aa1bSDebojyoti Ghosh   TS_GLEE  *glee = (TS_GLEE *)ts->data;
6712302aa1bSDebojyoti Ghosh   DM        dm, dmsave;
6722302aa1bSDebojyoti Ghosh   Vec       Ydot;
6732302aa1bSDebojyoti Ghosh   PetscReal shift = glee->scoeff / ts->time_step;
6742302aa1bSDebojyoti Ghosh 
6752302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
6769566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
6779566063dSJacob Faibussowitsch   PetscCall(TSGLEEGetVecs(ts, dm, &Ydot));
6782302aa1bSDebojyoti Ghosh   /* Set Ydot = shift*X */
6799566063dSJacob Faibussowitsch   PetscCall(VecCopy(X, Ydot));
6809566063dSJacob Faibussowitsch   PetscCall(VecScale(Ydot, shift));
6812302aa1bSDebojyoti Ghosh   dmsave = ts->dm;
6822302aa1bSDebojyoti Ghosh   ts->dm = dm;
6832302aa1bSDebojyoti Ghosh 
6849566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts, glee->stage_time, X, Ydot, F, PETSC_FALSE));
6852302aa1bSDebojyoti Ghosh 
6862302aa1bSDebojyoti Ghosh   ts->dm = dmsave;
6879566063dSJacob Faibussowitsch   PetscCall(TSGLEERestoreVecs(ts, dm, &Ydot));
6883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6892302aa1bSDebojyoti Ghosh }
6902302aa1bSDebojyoti Ghosh 
SNESTSFormJacobian_GLEE(SNES snes,Vec X,Mat A,Mat B,TS ts)691d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_GLEE(SNES snes, Vec X, Mat A, Mat B, TS ts)
692d71ae5a4SJacob Faibussowitsch {
6932302aa1bSDebojyoti Ghosh   TS_GLEE  *glee = (TS_GLEE *)ts->data;
6942302aa1bSDebojyoti Ghosh   DM        dm, dmsave;
6952302aa1bSDebojyoti Ghosh   Vec       Ydot;
6962302aa1bSDebojyoti Ghosh   PetscReal shift = glee->scoeff / ts->time_step;
6972302aa1bSDebojyoti Ghosh 
6982302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
6999566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
7009566063dSJacob Faibussowitsch   PetscCall(TSGLEEGetVecs(ts, dm, &Ydot));
7012302aa1bSDebojyoti Ghosh   /* glee->Ydot has already been computed in SNESTSFormFunction_GLEE (SNES guarantees this) */
7022302aa1bSDebojyoti Ghosh   dmsave = ts->dm;
7032302aa1bSDebojyoti Ghosh   ts->dm = dm;
7042302aa1bSDebojyoti Ghosh 
7059566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts, glee->stage_time, X, Ydot, shift, A, B, PETSC_FALSE));
7062302aa1bSDebojyoti Ghosh 
7072302aa1bSDebojyoti Ghosh   ts->dm = dmsave;
7089566063dSJacob Faibussowitsch   PetscCall(TSGLEERestoreVecs(ts, dm, &Ydot));
7093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7102302aa1bSDebojyoti Ghosh }
7112302aa1bSDebojyoti Ghosh 
DMCoarsenHook_TSGLEE(DM fine,DM coarse,PetscCtx ctx)712*2a8381b2SBarry Smith static PetscErrorCode DMCoarsenHook_TSGLEE(DM fine, DM coarse, PetscCtx ctx)
713d71ae5a4SJacob Faibussowitsch {
7142302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7162302aa1bSDebojyoti Ghosh }
7172302aa1bSDebojyoti Ghosh 
DMRestrictHook_TSGLEE(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,PetscCtx ctx)718*2a8381b2SBarry Smith static PetscErrorCode DMRestrictHook_TSGLEE(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, PetscCtx ctx)
719d71ae5a4SJacob Faibussowitsch {
7202302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7222302aa1bSDebojyoti Ghosh }
7232302aa1bSDebojyoti Ghosh 
DMSubDomainHook_TSGLEE(DM dm,DM subdm,PetscCtx ctx)724*2a8381b2SBarry Smith static PetscErrorCode DMSubDomainHook_TSGLEE(DM dm, DM subdm, PetscCtx ctx)
725d71ae5a4SJacob Faibussowitsch {
7262302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7282302aa1bSDebojyoti Ghosh }
7292302aa1bSDebojyoti Ghosh 
DMSubDomainRestrictHook_TSGLEE(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,PetscCtx ctx)730*2a8381b2SBarry Smith static PetscErrorCode DMSubDomainRestrictHook_TSGLEE(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, PetscCtx ctx)
731d71ae5a4SJacob Faibussowitsch {
7322302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7342302aa1bSDebojyoti Ghosh }
7352302aa1bSDebojyoti Ghosh 
TSSetUp_GLEE(TS ts)736d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_GLEE(TS ts)
737d71ae5a4SJacob Faibussowitsch {
7382302aa1bSDebojyoti Ghosh   TS_GLEE    *glee = (TS_GLEE *)ts->data;
7392302aa1bSDebojyoti Ghosh   GLEETableau tab;
7402302aa1bSDebojyoti Ghosh   PetscInt    s, r;
7412302aa1bSDebojyoti Ghosh   DM          dm;
7422302aa1bSDebojyoti Ghosh 
7432302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
74448a46eb9SPierre Jolivet   if (!glee->tableau) PetscCall(TSGLEESetType(ts, TSGLEEDefaultType));
7452302aa1bSDebojyoti Ghosh   tab = glee->tableau;
7462302aa1bSDebojyoti Ghosh   s   = tab->s;
7472302aa1bSDebojyoti Ghosh   r   = tab->r;
7489566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vec_sol, r, &glee->Y));
7499566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vec_sol, r, &glee->X));
7509566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vec_sol, s, &glee->YStage));
7519566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vec_sol, s, &glee->YdotStage));
7529566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &glee->Ydot));
7539566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &glee->yGErr));
7549566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(glee->yGErr));
7559566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &glee->W));
7569566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(s, &glee->swork, r, &glee->rwork));
7579566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
7589566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSGLEE, DMRestrictHook_TSGLEE, ts));
7599566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSGLEE, DMSubDomainRestrictHook_TSGLEE, ts));
7603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7612302aa1bSDebojyoti Ghosh }
7622302aa1bSDebojyoti Ghosh 
TSStartingMethod_GLEE(TS ts)76366976f2fSJacob Faibussowitsch static PetscErrorCode TSStartingMethod_GLEE(TS ts)
764d71ae5a4SJacob Faibussowitsch {
7652302aa1bSDebojyoti Ghosh   TS_GLEE    *glee = (TS_GLEE *)ts->data;
7666e2e232bSDebojyoti Ghosh   GLEETableau tab  = glee->tableau;
7676e2e232bSDebojyoti Ghosh   PetscInt    r    = tab->r, i;
7686e2e232bSDebojyoti Ghosh   PetscReal  *S    = tab->S;
7692302aa1bSDebojyoti Ghosh 
7702302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7712302aa1bSDebojyoti Ghosh   for (i = 0; i < r; i++) {
7729566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(glee->Y[i]));
7739566063dSJacob Faibussowitsch     PetscCall(VecAXPY(glee->Y[i], S[i], ts->vec_sol));
7742302aa1bSDebojyoti Ghosh   }
7753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7762302aa1bSDebojyoti Ghosh }
7772302aa1bSDebojyoti Ghosh 
TSSetFromOptions_GLEE(TS ts,PetscOptionItems PetscOptionsObject)778ce78bad3SBarry Smith static PetscErrorCode TSSetFromOptions_GLEE(TS ts, PetscOptionItems PetscOptionsObject)
779d71ae5a4SJacob Faibussowitsch {
7802302aa1bSDebojyoti Ghosh   char gleetype[256];
7812302aa1bSDebojyoti Ghosh 
7822302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
783d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "GLEE ODE solver options");
7842302aa1bSDebojyoti Ghosh   {
7852302aa1bSDebojyoti Ghosh     GLEETableauLink link;
7862302aa1bSDebojyoti Ghosh     PetscInt        count, choice;
7872302aa1bSDebojyoti Ghosh     PetscBool       flg;
7882302aa1bSDebojyoti Ghosh     const char    **namelist;
7892302aa1bSDebojyoti Ghosh 
7909566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(gleetype, TSGLEEDefaultType, sizeof(gleetype)));
791fbccb6d4SPierre Jolivet     for (link = GLEETableauList, count = 0; link; link = link->next, count++);
7929566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(count, (char ***)&namelist));
7932302aa1bSDebojyoti Ghosh     for (link = GLEETableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name;
7949566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEList("-ts_glee_type", "Family of GLEE method", "TSGLEESetType", (const char *const *)namelist, count, gleetype, &choice, &flg));
7959566063dSJacob Faibussowitsch     PetscCall(TSGLEESetType(ts, flg ? namelist[choice] : gleetype));
7969566063dSJacob Faibussowitsch     PetscCall(PetscFree(namelist));
7972302aa1bSDebojyoti Ghosh   }
798d0609cedSBarry Smith   PetscOptionsHeadEnd();
7993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8002302aa1bSDebojyoti Ghosh }
8012302aa1bSDebojyoti Ghosh 
TSView_GLEE(TS ts,PetscViewer viewer)802d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_GLEE(TS ts, PetscViewer viewer)
803d71ae5a4SJacob Faibussowitsch {
8042302aa1bSDebojyoti Ghosh   TS_GLEE    *glee = (TS_GLEE *)ts->data;
8052302aa1bSDebojyoti Ghosh   GLEETableau tab  = glee->tableau;
8069f196a02SMartin Diehl   PetscBool   isascii;
8072302aa1bSDebojyoti Ghosh 
8082302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8099f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
8109f196a02SMartin Diehl   if (isascii) {
8112302aa1bSDebojyoti Ghosh     TSGLEEType gleetype;
8122302aa1bSDebojyoti Ghosh     char       buf[512];
8139566063dSJacob Faibussowitsch     PetscCall(TSGLEEGetType(ts, &gleetype));
8149566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  GLEE type %s\n", gleetype));
8159566063dSJacob Faibussowitsch     PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, tab->c));
8169566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Abscissa     c = %s\n", buf));
8172302aa1bSDebojyoti Ghosh     /* Note: print out r as well */
8182302aa1bSDebojyoti Ghosh   }
8193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8202302aa1bSDebojyoti Ghosh }
8212302aa1bSDebojyoti Ghosh 
TSLoad_GLEE(TS ts,PetscViewer viewer)822d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSLoad_GLEE(TS ts, PetscViewer viewer)
823d71ae5a4SJacob Faibussowitsch {
8242302aa1bSDebojyoti Ghosh   SNES    snes;
8252302aa1bSDebojyoti Ghosh   TSAdapt tsadapt;
8262302aa1bSDebojyoti Ghosh 
8272302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8289566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &tsadapt));
8299566063dSJacob Faibussowitsch   PetscCall(TSAdaptLoad(tsadapt, viewer));
8309566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
8319566063dSJacob Faibussowitsch   PetscCall(SNESLoad(snes, viewer));
8322302aa1bSDebojyoti Ghosh   /* function and Jacobian context for SNES when used with TS is always ts object */
8339566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, NULL, NULL, ts));
8349566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, ts));
8353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8362302aa1bSDebojyoti Ghosh }
8372302aa1bSDebojyoti Ghosh 
838cc4c1da9SBarry Smith /*@
839bcf0153eSBarry Smith   TSGLEESetType - Set the type of `TSGLEE` scheme
8402302aa1bSDebojyoti Ghosh 
84120f4b53cSBarry Smith   Logically Collective
8422302aa1bSDebojyoti Ghosh 
843d8d19677SJose E. Roman   Input Parameters:
8442302aa1bSDebojyoti Ghosh + ts       - timestepping context
845bcf0153eSBarry Smith - gleetype - type of `TSGLEE` scheme
8462302aa1bSDebojyoti Ghosh 
8472302aa1bSDebojyoti Ghosh   Level: intermediate
8482302aa1bSDebojyoti Ghosh 
8491cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEEGetType()`, `TSGLEE`
8502302aa1bSDebojyoti Ghosh @*/
TSGLEESetType(TS ts,TSGLEEType gleetype)851d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLEESetType(TS ts, TSGLEEType gleetype)
852d71ae5a4SJacob Faibussowitsch {
8532302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8542302aa1bSDebojyoti Ghosh   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
8554f572ea9SToby Isaac   PetscAssertPointer(gleetype, 2);
856cac4c232SBarry Smith   PetscTryMethod(ts, "TSGLEESetType_C", (TS, TSGLEEType), (ts, gleetype));
8573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8582302aa1bSDebojyoti Ghosh }
8592302aa1bSDebojyoti Ghosh 
860cc4c1da9SBarry Smith /*@
861bcf0153eSBarry Smith   TSGLEEGetType - Get the type of `TSGLEE` scheme
8622302aa1bSDebojyoti Ghosh 
86320f4b53cSBarry Smith   Logically Collective
8642302aa1bSDebojyoti Ghosh 
8652302aa1bSDebojyoti Ghosh   Input Parameter:
8662302aa1bSDebojyoti Ghosh . ts - timestepping context
8672302aa1bSDebojyoti Ghosh 
8682302aa1bSDebojyoti Ghosh   Output Parameter:
869bcf0153eSBarry Smith . gleetype - type of `TSGLEE` scheme
8702302aa1bSDebojyoti Ghosh 
8712302aa1bSDebojyoti Ghosh   Level: intermediate
8722302aa1bSDebojyoti Ghosh 
8731cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEE`, `TSGLEESetType()`
8742302aa1bSDebojyoti Ghosh @*/
TSGLEEGetType(TS ts,TSGLEEType * gleetype)875d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLEEGetType(TS ts, TSGLEEType *gleetype)
876d71ae5a4SJacob Faibussowitsch {
8772302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8782302aa1bSDebojyoti Ghosh   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
879cac4c232SBarry Smith   PetscUseMethod(ts, "TSGLEEGetType_C", (TS, TSGLEEType *), (ts, gleetype));
8803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8812302aa1bSDebojyoti Ghosh }
8822302aa1bSDebojyoti Ghosh 
TSGLEEGetType_GLEE(TS ts,TSGLEEType * gleetype)88366976f2fSJacob Faibussowitsch static PetscErrorCode TSGLEEGetType_GLEE(TS ts, TSGLEEType *gleetype)
884d71ae5a4SJacob Faibussowitsch {
8852302aa1bSDebojyoti Ghosh   TS_GLEE *glee = (TS_GLEE *)ts->data;
8862302aa1bSDebojyoti Ghosh 
8872302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
88848a46eb9SPierre Jolivet   if (!glee->tableau) PetscCall(TSGLEESetType(ts, TSGLEEDefaultType));
8892302aa1bSDebojyoti Ghosh   *gleetype = glee->tableau->name;
8903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8912302aa1bSDebojyoti Ghosh }
TSGLEESetType_GLEE(TS ts,TSGLEEType gleetype)89266976f2fSJacob Faibussowitsch static PetscErrorCode TSGLEESetType_GLEE(TS ts, TSGLEEType gleetype)
893d71ae5a4SJacob Faibussowitsch {
8942302aa1bSDebojyoti Ghosh   TS_GLEE        *glee = (TS_GLEE *)ts->data;
8952302aa1bSDebojyoti Ghosh   PetscBool       match;
8962302aa1bSDebojyoti Ghosh   GLEETableauLink link;
8972302aa1bSDebojyoti Ghosh 
8982302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8992302aa1bSDebojyoti Ghosh   if (glee->tableau) {
9009566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(glee->tableau->name, gleetype, &match));
9013ba16761SJacob Faibussowitsch     if (match) PetscFunctionReturn(PETSC_SUCCESS);
9022302aa1bSDebojyoti Ghosh   }
9032302aa1bSDebojyoti Ghosh   for (link = GLEETableauList; link; link = link->next) {
9049566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(link->tab.name, gleetype, &match));
9052302aa1bSDebojyoti Ghosh     if (match) {
9069566063dSJacob Faibussowitsch       PetscCall(TSReset_GLEE(ts));
9072302aa1bSDebojyoti Ghosh       glee->tableau = &link->tab;
9083ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
9092302aa1bSDebojyoti Ghosh     }
9102302aa1bSDebojyoti Ghosh   }
91198921bdaSJacob Faibussowitsch   SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", gleetype);
9122302aa1bSDebojyoti Ghosh }
9132302aa1bSDebojyoti Ghosh 
TSGetStages_GLEE(TS ts,PetscInt * ns,Vec ** Y)914d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGetStages_GLEE(TS ts, PetscInt *ns, Vec **Y)
915d71ae5a4SJacob Faibussowitsch {
9162302aa1bSDebojyoti Ghosh   TS_GLEE *glee = (TS_GLEE *)ts->data;
9172302aa1bSDebojyoti Ghosh 
9182302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
9190429704eSStefano Zampini   if (ns) *ns = glee->tableau->s;
920d5509560SEmil Constantinescu   if (Y) *Y = glee->YStage;
9213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9222302aa1bSDebojyoti Ghosh }
9232302aa1bSDebojyoti Ghosh 
TSGetSolutionComponents_GLEE(TS ts,PetscInt * n,Vec * Y)92466976f2fSJacob Faibussowitsch static PetscErrorCode TSGetSolutionComponents_GLEE(TS ts, PetscInt *n, Vec *Y)
925d71ae5a4SJacob Faibussowitsch {
9262302aa1bSDebojyoti Ghosh   TS_GLEE    *glee = (TS_GLEE *)ts->data;
9272302aa1bSDebojyoti Ghosh   GLEETableau tab  = glee->tableau;
9282302aa1bSDebojyoti Ghosh 
9292302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
9304cdd57e5SDebojyoti Ghosh   if (!Y) *n = tab->r;
9312302aa1bSDebojyoti Ghosh   else {
932966bd95aSPierre Jolivet     PetscCheck(*n >= 0 && *n < tab->r, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Second argument (%" PetscInt_FMT ") out of range[0,%" PetscInt_FMT "].", *n, tab->r - 1);
9339566063dSJacob Faibussowitsch     PetscCall(VecCopy(glee->Y[*n], *Y));
9342302aa1bSDebojyoti Ghosh   }
9353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9362302aa1bSDebojyoti Ghosh }
9372302aa1bSDebojyoti Ghosh 
TSGetAuxSolution_GLEE(TS ts,Vec * X)93866976f2fSJacob Faibussowitsch static PetscErrorCode TSGetAuxSolution_GLEE(TS ts, Vec *X)
939d71ae5a4SJacob Faibussowitsch {
9404cdd57e5SDebojyoti Ghosh   TS_GLEE     *glee = (TS_GLEE *)ts->data;
9414cdd57e5SDebojyoti Ghosh   GLEETableau  tab  = glee->tableau;
9424cdd57e5SDebojyoti Ghosh   PetscReal   *F    = tab->Fembed;
9434cdd57e5SDebojyoti Ghosh   PetscInt     r    = tab->r;
9444cdd57e5SDebojyoti Ghosh   Vec         *Y    = glee->Y;
945ec0e3ee2SEmil Constantinescu   PetscScalar *wr   = glee->rwork;
946ec0e3ee2SEmil Constantinescu   PetscInt     i;
9474cdd57e5SDebojyoti Ghosh 
9484cdd57e5SDebojyoti Ghosh   PetscFunctionBegin;
9499566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(*X));
950ec0e3ee2SEmil Constantinescu   for (i = 0; i < r; i++) wr[i] = F[i];
951f4f49eeaSPierre Jolivet   PetscCall(VecMAXPY(*X, r, wr, Y));
9523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9534cdd57e5SDebojyoti Ghosh }
9544cdd57e5SDebojyoti Ghosh 
TSGetTimeError_GLEE(TS ts,PetscInt n,Vec * X)95566976f2fSJacob Faibussowitsch static PetscErrorCode TSGetTimeError_GLEE(TS ts, PetscInt n, Vec *X)
956d71ae5a4SJacob Faibussowitsch {
9574cdd57e5SDebojyoti Ghosh   TS_GLEE     *glee = (TS_GLEE *)ts->data;
9584cdd57e5SDebojyoti Ghosh   GLEETableau  tab  = glee->tableau;
9594cdd57e5SDebojyoti Ghosh   PetscReal   *F    = tab->Ferror;
9604cdd57e5SDebojyoti Ghosh   PetscInt     r    = tab->r;
9614cdd57e5SDebojyoti Ghosh   Vec         *Y    = glee->Y;
962ec0e3ee2SEmil Constantinescu   PetscScalar *wr   = glee->rwork;
963ec0e3ee2SEmil Constantinescu   PetscInt     i;
9644cdd57e5SDebojyoti Ghosh 
9654cdd57e5SDebojyoti Ghosh   PetscFunctionBegin;
9669566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(*X));
967657c1e31SEmil Constantinescu   if (n == 0) {
968ec0e3ee2SEmil Constantinescu     for (i = 0; i < r; i++) wr[i] = F[i];
969f4f49eeaSPierre Jolivet     PetscCall(VecMAXPY(*X, r, wr, Y));
9700a01e1b2SEmil Constantinescu   } else if (n == -1) {
971657c1e31SEmil Constantinescu     *X = glee->yGErr;
9720a01e1b2SEmil Constantinescu   }
9733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9744cdd57e5SDebojyoti Ghosh }
9754cdd57e5SDebojyoti Ghosh 
TSSetTimeError_GLEE(TS ts,Vec X)97666976f2fSJacob Faibussowitsch static PetscErrorCode TSSetTimeError_GLEE(TS ts, Vec X)
977d71ae5a4SJacob Faibussowitsch {
97857df6a1bSDebojyoti Ghosh   TS_GLEE    *glee = (TS_GLEE *)ts->data;
97957df6a1bSDebojyoti Ghosh   GLEETableau tab  = glee->tableau;
98057df6a1bSDebojyoti Ghosh   PetscReal  *S    = tab->Serror;
98157df6a1bSDebojyoti Ghosh   PetscInt    r    = tab->r, i;
98257df6a1bSDebojyoti Ghosh   Vec        *Y    = glee->Y;
98357df6a1bSDebojyoti Ghosh 
98457df6a1bSDebojyoti Ghosh   PetscFunctionBegin;
98563a3b9bcSJacob Faibussowitsch   PetscCheck(r == 2, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSSetTimeError_GLEE not supported for '%s' with r=%" PetscInt_FMT ".", tab->name, tab->r);
98657df6a1bSDebojyoti Ghosh   for (i = 1; i < r; i++) {
9879566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vec_sol, Y[i]));
9889566063dSJacob Faibussowitsch     PetscCall(VecAXPBY(Y[i], S[0], S[1], X));
9899566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, glee->yGErr));
99057df6a1bSDebojyoti Ghosh   }
9913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
99257df6a1bSDebojyoti Ghosh }
99357df6a1bSDebojyoti Ghosh 
TSDestroy_GLEE(TS ts)994d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_GLEE(TS ts)
995d71ae5a4SJacob Faibussowitsch {
996b3a6b972SJed Brown   PetscFunctionBegin;
9979566063dSJacob Faibussowitsch   PetscCall(TSReset_GLEE(ts));
998b3a6b972SJed Brown   if (ts->dm) {
9999566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSGLEE, DMRestrictHook_TSGLEE, ts));
10009566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSGLEE, DMSubDomainRestrictHook_TSGLEE, ts));
1001b3a6b972SJed Brown   }
10029566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
10039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLEEGetType_C", NULL));
10049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLEESetType_C", NULL));
10053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1006b3a6b972SJed Brown }
1007b3a6b972SJed Brown 
10082302aa1bSDebojyoti Ghosh /*MC
10092302aa1bSDebojyoti Ghosh   TSGLEE - ODE and DAE solver using General Linear with Error Estimation schemes
10102302aa1bSDebojyoti Ghosh 
1011efa39862SBarry Smith   The user should provide the right-hand side of the equation using `TSSetRHSFunction()`
1012efa39862SBarry Smith   and for `TSGLEEi1` the Jacobian of the right-hand side using `TSSetRHSJacobian()`
10132302aa1bSDebojyoti Ghosh 
10142302aa1bSDebojyoti Ghosh   Level: beginner
10152302aa1bSDebojyoti Ghosh 
1016efa39862SBarry Smith   Notes:
1017efa39862SBarry Smith   The default is `TSGLEE35`, it can be changed with `TSGLEESetType()` or `-ts_glee_type type`
10182302aa1bSDebojyoti Ghosh 
1019efa39862SBarry Smith   The only implicit scheme is `TSGLEEi1`
1020efa39862SBarry Smith 
1021efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_glee), `TSCreate()`, `TS`, `TSSetType()`, `TSGLEESetType()`, `TSGLEEGetType()`,
1022f8d70eaaSPierre Jolivet           `TSGLEE23`, `TSGLEE24`, `TSGLEE35`, `TSGLEE25I`, `TSGLEEEXRK2A`,
1023bcf0153eSBarry Smith           `TSGLEERK32G1`, `TSGLEERK285EX`, `TSGLEEType`, `TSGLEERegister()`, `TSType`
10242302aa1bSDebojyoti Ghosh M*/
TSCreate_GLEE(TS ts)1025d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_GLEE(TS ts)
1026d71ae5a4SJacob Faibussowitsch {
10272302aa1bSDebojyoti Ghosh   TS_GLEE *th;
10282302aa1bSDebojyoti Ghosh 
10292302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
10309566063dSJacob Faibussowitsch   PetscCall(TSGLEEInitializePackage());
10312302aa1bSDebojyoti Ghosh 
10322302aa1bSDebojyoti Ghosh   ts->ops->reset                 = TSReset_GLEE;
10332302aa1bSDebojyoti Ghosh   ts->ops->destroy               = TSDestroy_GLEE;
10342302aa1bSDebojyoti Ghosh   ts->ops->view                  = TSView_GLEE;
10352302aa1bSDebojyoti Ghosh   ts->ops->load                  = TSLoad_GLEE;
10362302aa1bSDebojyoti Ghosh   ts->ops->setup                 = TSSetUp_GLEE;
10372302aa1bSDebojyoti Ghosh   ts->ops->step                  = TSStep_GLEE;
10382302aa1bSDebojyoti Ghosh   ts->ops->interpolate           = TSInterpolate_GLEE;
10392302aa1bSDebojyoti Ghosh   ts->ops->evaluatestep          = TSEvaluateStep_GLEE;
10402302aa1bSDebojyoti Ghosh   ts->ops->setfromoptions        = TSSetFromOptions_GLEE;
10412302aa1bSDebojyoti Ghosh   ts->ops->getstages             = TSGetStages_GLEE;
10422302aa1bSDebojyoti Ghosh   ts->ops->snesfunction          = SNESTSFormFunction_GLEE;
10432302aa1bSDebojyoti Ghosh   ts->ops->snesjacobian          = SNESTSFormJacobian_GLEE;
1044b2bf4f3aSDebojyoti Ghosh   ts->ops->getsolutioncomponents = TSGetSolutionComponents_GLEE;
10454cdd57e5SDebojyoti Ghosh   ts->ops->getauxsolution        = TSGetAuxSolution_GLEE;
10464cdd57e5SDebojyoti Ghosh   ts->ops->gettimeerror          = TSGetTimeError_GLEE;
104757df6a1bSDebojyoti Ghosh   ts->ops->settimeerror          = TSSetTimeError_GLEE;
10482302aa1bSDebojyoti Ghosh   ts->ops->startingmethod        = TSStartingMethod_GLEE;
1049b80d5313SLisandro Dalcin   ts->default_adapt_type         = TSADAPTGLEE;
10502302aa1bSDebojyoti Ghosh 
1051efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
1052efd4aadfSBarry Smith 
10534dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&th));
10542302aa1bSDebojyoti Ghosh   ts->data = (void *)th;
10552302aa1bSDebojyoti Ghosh 
10569566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLEEGetType_C", TSGLEEGetType_GLEE));
10579566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLEESetType_C", TSGLEESetType_GLEE));
10583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10592302aa1bSDebojyoti Ghosh }
1060