xref: /petsc/src/ts/impls/glee/glee.c (revision cc4c1da905d89950b196b027190941013bd3e15a)
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
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 
721cc06b55SBarry Smith .seealso: [](ch_ts), `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 
821cc06b55SBarry Smith .seealso: [](ch_ts), `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 
921cc06b55SBarry Smith .seealso: [](ch_ts), `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 
1021cc06b55SBarry Smith .seealso: [](ch_ts), `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 
1121cc06b55SBarry Smith .seealso: [](ch_ts), `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 
1221cc06b55SBarry Smith .seealso: [](ch_ts), `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 
1321cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEE`
13336c34d14SEmil Constantinescu M*/
1342302aa1bSDebojyoti Ghosh 
1352302aa1bSDebojyoti Ghosh /*@C
136bcf0153eSBarry Smith   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 
1421cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEERegisterDestroy()`
1432302aa1bSDebojyoti Ghosh @*/
144d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLEERegisterAll(void)
145d71ae5a4SJacob Faibussowitsch {
1462302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
1473ba16761SJacob Faibussowitsch   if (TSGLEERegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
1482302aa1bSDebojyoti Ghosh   TSGLEERegisterAllCalled = PETSC_TRUE;
1492302aa1bSDebojyoti Ghosh 
1502302aa1bSDebojyoti Ghosh   {
151b1fbfd33SSatish Balay #define GAMMA 0.5
1522302aa1bSDebojyoti Ghosh     /* y-eps form */
1539371c9d4SSatish Balay     const PetscInt  p = 1, s = 3, r = 2;
1549371c9d4SSatish Balay     const PetscReal A[3][3] =
1559371c9d4SSatish Balay       {
1569371c9d4SSatish Balay         {1.0, 0,   0  },
1579371c9d4SSatish Balay         {0,   0.5, 0  },
1589371c9d4SSatish Balay         {0,   0.5, 0.5}
1599371c9d4SSatish Balay     },
1609371c9d4SSatish 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};
1619566063dSJacob 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));
162ab8c5912SEmil Constantinescu   }
163ab8c5912SEmil Constantinescu   {
164b1fbfd33SSatish Balay #undef GAMMA
165b1fbfd33SSatish Balay #define GAMMA 0.0
166ab8c5912SEmil Constantinescu     /* y-eps form */
1679371c9d4SSatish Balay     const PetscInt  p = 2, s = 3, r = 2;
1689371c9d4SSatish Balay     const PetscReal A[3][3] =
1699371c9d4SSatish Balay       {
1709371c9d4SSatish Balay         {0,    0,    0},
1719371c9d4SSatish Balay         {1,    0,    0},
1729371c9d4SSatish Balay         {0.25, 0.25, 0}
1739371c9d4SSatish Balay     },
1749371c9d4SSatish 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};
1759566063dSJacob 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));
1762302aa1bSDebojyoti Ghosh   }
1772302aa1bSDebojyoti Ghosh   {
178b1fbfd33SSatish Balay #undef GAMMA
179b1fbfd33SSatish Balay #define GAMMA 0.0
1802302aa1bSDebojyoti Ghosh     /* y-y~ form */
1819371c9d4SSatish Balay     const PetscInt  p = 2, s = 4, r = 2;
1829371c9d4SSatish Balay     const PetscReal A[4][4] =
1839371c9d4SSatish Balay       {
1849371c9d4SSatish Balay         {0,            0,            0,            0},
1859371c9d4SSatish Balay         {0.75,         0,            0,            0},
1869371c9d4SSatish Balay         {0.25,         29.0 / 60.0,  0,            0},
1879371c9d4SSatish Balay         {-21.0 / 44.0, 145.0 / 44.0, -20.0 / 11.0, 0}
1889371c9d4SSatish Balay     },
1899371c9d4SSatish 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};
1909566063dSJacob 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));
1912302aa1bSDebojyoti Ghosh   }
1922302aa1bSDebojyoti Ghosh   {
193b1fbfd33SSatish Balay #undef GAMMA
194b1fbfd33SSatish Balay #define GAMMA 0.0
1952302aa1bSDebojyoti Ghosh     /* y-y~ form */
1969371c9d4SSatish Balay     const PetscInt  p = 2, s = 5, r = 2;
1979371c9d4SSatish Balay     const PetscReal A[5][5] =
1989371c9d4SSatish Balay       {
1999371c9d4SSatish Balay         {0,                       0,                       0,                       0,                      0},
2002302aa1bSDebojyoti Ghosh         {-0.94079244066783383269, 0,                       0,                       0,                      0},
2012302aa1bSDebojyoti Ghosh         {0.64228187778301907108,  0.10915356933958500042,  0,                       0,                      0},
2022302aa1bSDebojyoti Ghosh         {-0.51764297742287450812, 0.74414270351096040738,  -0.71404164927824538121, 0,                      0},
2039371c9d4SSatish Balay         {-0.44696561556825969206, -0.76768425657590196518, 0.20111608138142987881,  0.93828186737840469796, 0}
2049371c9d4SSatish Balay     },
2059371c9d4SSatish 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};
2069566063dSJacob 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));
2072302aa1bSDebojyoti Ghosh   }
2082302aa1bSDebojyoti Ghosh   {
209b1fbfd33SSatish Balay #undef GAMMA
210b1fbfd33SSatish Balay #define GAMMA 0.0
2112302aa1bSDebojyoti Ghosh     /* y-y~ form */
2129371c9d4SSatish Balay     const PetscInt  p = 3, s = 5, r = 2;
2139371c9d4SSatish Balay     const PetscReal A[5][5] =
2149371c9d4SSatish Balay       {
2159371c9d4SSatish Balay         {0,                                                0,                                                 0,                                                 0,                                               0},
2162302aa1bSDebojyoti Ghosh         {-2169604947363702313.0 / 24313474998937147335.0,  0,                                                 0,                                                 0,                                               0},
2172302aa1bSDebojyoti Ghosh         {46526746497697123895.0 / 94116917485856474137.0,  -10297879244026594958.0 / 49199457603717988219.0,  0,                                                 0,                                               0},
2182302aa1bSDebojyoti Ghosh         {23364788935845982499.0 / 87425311444725389446.0,  -79205144337496116638.0 / 148994349441340815519.0, 40051189859317443782.0 / 36487615018004984309.0,   0,                                               0},
2199371c9d4SSatish Balay         {42089522664062539205.0 / 124911313006412840286.0, -15074384760342762939.0 / 137927286865289746282.0, -62274678522253371016.0 / 125918573676298591413.0, 13755475729852471739.0 / 79257927066651693390.0, 0}
2209371c9d4SSatish Balay     },
2219371c9d4SSatish 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};
2229566063dSJacob 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));
2232302aa1bSDebojyoti Ghosh   }
2242302aa1bSDebojyoti Ghosh   {
225b1fbfd33SSatish Balay #undef GAMMA
226b1fbfd33SSatish Balay #define GAMMA 0.25
2272302aa1bSDebojyoti Ghosh     /* y-eps form */
2289371c9d4SSatish Balay     const PetscInt  p = 2, s = 6, r = 2;
2299371c9d4SSatish Balay     const PetscReal A[6][6] =
2309371c9d4SSatish Balay       {
2319371c9d4SSatish Balay         {0, 0, 0,    0,    0,   0},
2322302aa1bSDebojyoti Ghosh         {1, 0, 0,    0,    0,   0},
2332302aa1bSDebojyoti Ghosh         {0, 0, 0,    0,    0,   0},
2342302aa1bSDebojyoti Ghosh         {0, 0, 0.5,  0,    0,   0},
2352302aa1bSDebojyoti Ghosh         {0, 0, 0.25, 0.25, 0,   0},
2369371c9d4SSatish Balay         {0, 0, 0.25, 0.25, 0.5, 0}
2379371c9d4SSatish Balay     },
2389371c9d4SSatish 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};
2399566063dSJacob 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));
2402302aa1bSDebojyoti Ghosh   }
2412302aa1bSDebojyoti Ghosh   {
242b1fbfd33SSatish Balay #undef GAMMA
243b1fbfd33SSatish Balay #define GAMMA 0.0
2442302aa1bSDebojyoti Ghosh     /* y-eps form */
2459371c9d4SSatish Balay     const PetscInt  p = 3, s = 8, r = 2;
2469371c9d4SSatish Balay     const PetscReal A[8][8] =
2479371c9d4SSatish Balay       {
2489371c9d4SSatish Balay         {0,           0,          0,          0,          0,         0,         0,         0},
2492302aa1bSDebojyoti Ghosh         {0.5,         0,          0,          0,          0,         0,         0,         0},
2502302aa1bSDebojyoti Ghosh         {-1,          2,          0,          0,          0,         0,         0,         0},
2512302aa1bSDebojyoti Ghosh         {1.0 / 6.0,   2.0 / 3.0,  1.0 / 6.0,  0,          0,         0,         0,         0},
2522302aa1bSDebojyoti Ghosh         {0,           0,          0,          0,          0,         0,         0,         0},
2532302aa1bSDebojyoti Ghosh         {-7.0 / 24.0, 1.0 / 3.0,  1.0 / 12.0, -1.0 / 8.0, 0.5,       0,         0,         0},
2542302aa1bSDebojyoti Ghosh         {7.0 / 6.0,   -4.0 / 3.0, -1.0 / 3.0, 0.5,        -1.0,      2.0,       0,         0},
2559371c9d4SSatish Balay         {0,           0,          0,          0,          1.0 / 6.0, 2.0 / 3.0, 1.0 / 6.0, 0}
2569371c9d4SSatish Balay     },
2579371c9d4SSatish 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};
2589566063dSJacob 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));
2592302aa1bSDebojyoti Ghosh   }
2602302aa1bSDebojyoti Ghosh   {
261b1fbfd33SSatish Balay #undef GAMMA
262b1fbfd33SSatish Balay #define GAMMA 0.25
2632302aa1bSDebojyoti Ghosh     /* y-eps form */
2649371c9d4SSatish Balay     const PetscInt  p = 2, s = 9, r = 2;
2659371c9d4SSatish Balay     const PetscReal A[9][9] =
2669371c9d4SSatish Balay       {
2679371c9d4SSatish Balay         {0,                    0,                    0, 0,                     0,                    0,                    0,                     0,                    0},
2682302aa1bSDebojyoti Ghosh         {0.585786437626904966, 0,                    0, 0,                     0,                    0,                    0,                     0,                    0},
2692302aa1bSDebojyoti Ghosh         {0.149999999999999994, 0.849999999999999978, 0, 0,                     0,                    0,                    0,                     0,                    0},
2702302aa1bSDebojyoti Ghosh         {0,                    0,                    0, 0,                     0,                    0,                    0,                     0,                    0},
2712302aa1bSDebojyoti Ghosh         {0,                    0,                    0, 0.292893218813452483,  0,                    0,                    0,                     0,                    0},
2722302aa1bSDebojyoti Ghosh         {0,                    0,                    0, 0.0749999999999999972, 0.424999999999999989, 0,                    0,                     0,                    0},
2732302aa1bSDebojyoti Ghosh         {0,                    0,                    0, 0.176776695296636893,  0.176776695296636893, 0.146446609406726241, 0,                     0,                    0},
2742302aa1bSDebojyoti Ghosh         {0,                    0,                    0, 0.176776695296636893,  0.176776695296636893, 0.146446609406726241, 0.292893218813452483,  0,                    0},
2759371c9d4SSatish Balay         {0,                    0,                    0, 0.176776695296636893,  0.176776695296636893, 0.146446609406726241, 0.0749999999999999972, 0.424999999999999989, 0}
2769371c9d4SSatish Balay     },
2779371c9d4SSatish 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};
2789566063dSJacob 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));
2792302aa1bSDebojyoti Ghosh   }
2803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2812302aa1bSDebojyoti Ghosh }
2822302aa1bSDebojyoti Ghosh 
2832302aa1bSDebojyoti Ghosh /*@C
284bcf0153eSBarry Smith   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 
2901cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEERegister()`, `TSGLEERegisterAll()`
2912302aa1bSDebojyoti Ghosh @*/
292d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLEERegisterDestroy(void)
293d71ae5a4SJacob Faibussowitsch {
2942302aa1bSDebojyoti Ghosh   GLEETableauLink link;
2952302aa1bSDebojyoti Ghosh 
2962302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
2972302aa1bSDebojyoti Ghosh   while ((link = GLEETableauList)) {
2982302aa1bSDebojyoti Ghosh     GLEETableau t   = &link->tab;
2992302aa1bSDebojyoti Ghosh     GLEETableauList = link->next;
3009566063dSJacob Faibussowitsch     PetscCall(PetscFree5(t->A, t->B, t->U, t->V, t->c));
301d0609cedSBarry Smith     PetscCall(PetscFree2(t->S, t->F));
302d0609cedSBarry Smith     PetscCall(PetscFree(t->Fembed));
303d0609cedSBarry Smith     PetscCall(PetscFree(t->Ferror));
304d0609cedSBarry Smith     PetscCall(PetscFree(t->Serror));
305d0609cedSBarry Smith     PetscCall(PetscFree(t->binterp));
306d0609cedSBarry Smith     PetscCall(PetscFree(t->name));
307d0609cedSBarry Smith     PetscCall(PetscFree(link));
3082302aa1bSDebojyoti Ghosh   }
3092302aa1bSDebojyoti Ghosh   TSGLEERegisterAllCalled = PETSC_FALSE;
3103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3112302aa1bSDebojyoti Ghosh }
3122302aa1bSDebojyoti Ghosh 
3132302aa1bSDebojyoti Ghosh /*@C
314bcf0153eSBarry Smith   TSGLEEInitializePackage - This function initializes everything in the `TSGLEE` package. It is called
315bcf0153eSBarry Smith   from `TSInitializePackage()`.
3162302aa1bSDebojyoti Ghosh 
3172302aa1bSDebojyoti Ghosh   Level: developer
3182302aa1bSDebojyoti Ghosh 
3191cc06b55SBarry Smith .seealso: [](ch_ts), `PetscInitialize()`
3202302aa1bSDebojyoti Ghosh @*/
321d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLEEInitializePackage(void)
322d71ae5a4SJacob Faibussowitsch {
3232302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
3243ba16761SJacob Faibussowitsch   if (TSGLEEPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
3252302aa1bSDebojyoti Ghosh   TSGLEEPackageInitialized = PETSC_TRUE;
3269566063dSJacob Faibussowitsch   PetscCall(TSGLEERegisterAll());
3279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposedDataRegister(&explicit_stage_time_id));
3289566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(TSGLEEFinalizePackage));
3293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3302302aa1bSDebojyoti Ghosh }
3312302aa1bSDebojyoti Ghosh 
3322302aa1bSDebojyoti Ghosh /*@C
333bcf0153eSBarry Smith   TSGLEEFinalizePackage - This function destroys everything in the `TSGLEE` package. It is
334bcf0153eSBarry Smith   called from `PetscFinalize()`.
3352302aa1bSDebojyoti Ghosh 
3362302aa1bSDebojyoti Ghosh   Level: developer
3372302aa1bSDebojyoti Ghosh 
3381cc06b55SBarry Smith .seealso: [](ch_ts), `PetscFinalize()`
3392302aa1bSDebojyoti Ghosh @*/
340d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLEEFinalizePackage(void)
341d71ae5a4SJacob Faibussowitsch {
3422302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
3432302aa1bSDebojyoti Ghosh   TSGLEEPackageInitialized = PETSC_FALSE;
3449566063dSJacob Faibussowitsch   PetscCall(TSGLEERegisterDestroy());
3453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3462302aa1bSDebojyoti Ghosh }
3472302aa1bSDebojyoti Ghosh 
3482302aa1bSDebojyoti Ghosh /*@C
349bcf0153eSBarry Smith   TSGLEERegister - register a new `TSGLEE` scheme by providing the entries in the Butcher tableau
3502302aa1bSDebojyoti Ghosh 
351*cc4c1da9SBarry Smith   Not Collective, but the same schemes should be registered on all processes on which they will be used, No Fortran Support
3522302aa1bSDebojyoti Ghosh 
3532302aa1bSDebojyoti Ghosh   Input Parameters:
3542302aa1bSDebojyoti Ghosh + name    - identifier for method
3552302aa1bSDebojyoti Ghosh . order   - order of method
3562302aa1bSDebojyoti Ghosh . s       - number of stages
3572302aa1bSDebojyoti Ghosh . r       - number of steps
3582302aa1bSDebojyoti Ghosh . gamma   - LTE ratio
3592302aa1bSDebojyoti Ghosh . A       - stage coefficients (dimension s*s, row-major)
3602302aa1bSDebojyoti Ghosh . B       - step completion coefficients (dimension r*s, row-major)
3612302aa1bSDebojyoti Ghosh . U       - method coefficients (dimension s*r, row-major)
3622302aa1bSDebojyoti Ghosh . V       - method coefficients (dimension r*r, row-major)
3632302aa1bSDebojyoti Ghosh . S       - starting coefficients
3642302aa1bSDebojyoti Ghosh . F       - finishing coefficients
3652302aa1bSDebojyoti Ghosh . c       - abscissa (dimension s; NULL to use row sums of A)
3662302aa1bSDebojyoti Ghosh . Fembed  - step completion coefficients for embedded method
367fd96ebc2SDebojyoti Ghosh . Ferror  - error computation coefficients
36857df6a1bSDebojyoti Ghosh . Serror  - error initialization coefficients
3692302aa1bSDebojyoti Ghosh . pinterp - order of interpolation (0 if unavailable)
3702302aa1bSDebojyoti Ghosh - binterp - array of interpolation coefficients (NULL if unavailable)
3712302aa1bSDebojyoti Ghosh 
3722302aa1bSDebojyoti Ghosh   Level: advanced
3732302aa1bSDebojyoti Ghosh 
374bcf0153eSBarry Smith   Note:
375bcf0153eSBarry Smith   Several `TSGLEE` methods are provided, this function is only needed to create new methods.
376bcf0153eSBarry Smith 
3771cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEE`
3782302aa1bSDebojyoti Ghosh @*/
379d71ae5a4SJacob 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[])
380d71ae5a4SJacob Faibussowitsch {
3812302aa1bSDebojyoti Ghosh   GLEETableauLink link;
3822302aa1bSDebojyoti Ghosh   GLEETableau     t;
3832302aa1bSDebojyoti Ghosh   PetscInt        i, j;
3842302aa1bSDebojyoti Ghosh 
3852302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
3869566063dSJacob Faibussowitsch   PetscCall(TSGLEEInitializePackage());
3879566063dSJacob Faibussowitsch   PetscCall(PetscNew(&link));
3882302aa1bSDebojyoti Ghosh   t = &link->tab;
3899566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &t->name));
3902302aa1bSDebojyoti Ghosh   t->order = order;
3912302aa1bSDebojyoti Ghosh   t->s     = s;
3922302aa1bSDebojyoti Ghosh   t->r     = r;
3932302aa1bSDebojyoti Ghosh   t->gamma = gamma;
3949566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(s * s, &t->A, r * r, &t->V, s, &t->c, r * s, &t->B, s * r, &t->U));
3959566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(r, &t->S, r, &t->F));
3969566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->A, A, s * s));
3979566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->B, B, r * s));
3989566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->U, U, s * r));
3999566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->V, V, r * r));
4009566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->S, S, r));
4019566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->F, F, r));
4022302aa1bSDebojyoti Ghosh   if (c) {
4039566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(t->c, c, s));
4042302aa1bSDebojyoti Ghosh   } else {
4059371c9d4SSatish Balay     for (i = 0; i < s; i++)
4069371c9d4SSatish Balay       for (j = 0, t->c[i] = 0; j < s; j++) t->c[i] += A[i * s + j];
4072302aa1bSDebojyoti Ghosh   }
4089566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(r, &t->Fembed));
4099566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(r, &t->Ferror));
4109566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(r, &t->Serror));
4119566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->Fembed, Fembed, r));
4129566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->Ferror, Ferror, r));
4139566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->Serror, Serror, r));
4142302aa1bSDebojyoti Ghosh   t->pinterp = pinterp;
4159566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s * pinterp, &t->binterp));
4169566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->binterp, binterp, s * pinterp));
4172302aa1bSDebojyoti Ghosh 
4182302aa1bSDebojyoti Ghosh   link->next      = GLEETableauList;
4192302aa1bSDebojyoti Ghosh   GLEETableauList = link;
4203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4212302aa1bSDebojyoti Ghosh }
4222302aa1bSDebojyoti Ghosh 
423d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateStep_GLEE(TS ts, PetscInt order, Vec X, PetscBool *done)
424d71ae5a4SJacob Faibussowitsch {
4252302aa1bSDebojyoti Ghosh   TS_GLEE     *glee = (TS_GLEE *)ts->data;
4262302aa1bSDebojyoti Ghosh   GLEETableau  tab  = glee->tableau;
4279371c9d4SSatish Balay   PetscReal    h, *B = tab->B, *V = tab->V, *F = tab->F, *Fembed = tab->Fembed;
4282302aa1bSDebojyoti Ghosh   PetscInt     s = tab->s, r = tab->r, i, j;
4292302aa1bSDebojyoti Ghosh   Vec         *Y = glee->Y, *YdotStage = glee->YdotStage;
430ec0e3ee2SEmil Constantinescu   PetscScalar *ws = glee->swork, *wr = glee->rwork;
4312302aa1bSDebojyoti Ghosh 
4322302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
4332302aa1bSDebojyoti Ghosh   switch (glee->status) {
4342302aa1bSDebojyoti Ghosh   case TS_STEP_INCOMPLETE:
435d71ae5a4SJacob Faibussowitsch   case TS_STEP_PENDING:
436d71ae5a4SJacob Faibussowitsch     h = ts->time_step;
437d71ae5a4SJacob Faibussowitsch     break;
438d71ae5a4SJacob Faibussowitsch   case TS_STEP_COMPLETE:
439d71ae5a4SJacob Faibussowitsch     h = ts->ptime - ts->ptime_prev;
440d71ae5a4SJacob Faibussowitsch     break;
441d71ae5a4SJacob Faibussowitsch   default:
442d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus");
4432302aa1bSDebojyoti Ghosh   }
4442302aa1bSDebojyoti Ghosh 
4452302aa1bSDebojyoti Ghosh   if (order == tab->order) {
4462302aa1bSDebojyoti Ghosh     /* Note: Irrespective of whether status is TS_STEP_INCOMPLETE
447fd96ebc2SDebojyoti Ghosh              or TS_STEP_COMPLETE, glee->X has the solution at the
4482302aa1bSDebojyoti Ghosh              beginning of the time step. So no need to roll-back.
4492302aa1bSDebojyoti Ghosh     */
4502302aa1bSDebojyoti Ghosh     if (glee->status == TS_STEP_INCOMPLETE) {
4512302aa1bSDebojyoti Ghosh       for (i = 0; i < r; i++) {
4529566063dSJacob Faibussowitsch         PetscCall(VecZeroEntries(Y[i]));
453ec0e3ee2SEmil Constantinescu         for (j = 0; j < r; j++) wr[j] = V[i * r + j];
4549566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(Y[i], r, wr, glee->X));
455ec0e3ee2SEmil Constantinescu         for (j = 0; j < s; j++) ws[j] = h * B[i * s + j];
4569566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(Y[i], s, ws, YdotStage));
4572302aa1bSDebojyoti Ghosh       }
4589566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(X));
459ec0e3ee2SEmil Constantinescu       for (j = 0; j < r; j++) wr[j] = F[j];
4609566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(X, r, wr, Y));
4619566063dSJacob Faibussowitsch     } else PetscCall(VecCopy(ts->vec_sol, X));
4623ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
4632302aa1bSDebojyoti Ghosh 
4642302aa1bSDebojyoti Ghosh   } else if (order == tab->order - 1) {
4652302aa1bSDebojyoti Ghosh     /* Complete with the embedded method (Fembed) */
4662302aa1bSDebojyoti Ghosh     for (i = 0; i < r; i++) {
4679566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(Y[i]));
468ec0e3ee2SEmil Constantinescu       for (j = 0; j < r; j++) wr[j] = V[i * r + j];
4699566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(Y[i], r, wr, glee->X));
470ec0e3ee2SEmil Constantinescu       for (j = 0; j < s; j++) ws[j] = h * B[i * s + j];
4719566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(Y[i], s, ws, YdotStage));
4722302aa1bSDebojyoti Ghosh     }
4739566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(X));
474ec0e3ee2SEmil Constantinescu     for (j = 0; j < r; j++) wr[j] = Fembed[j];
4759566063dSJacob Faibussowitsch     PetscCall(VecMAXPY(X, r, wr, Y));
476fd96ebc2SDebojyoti Ghosh 
4772302aa1bSDebojyoti Ghosh     if (done) *done = PETSC_TRUE;
4783ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
4792302aa1bSDebojyoti Ghosh   }
4802302aa1bSDebojyoti Ghosh   if (done) *done = PETSC_FALSE;
48163a3b9bcSJacob 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);
4823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4832302aa1bSDebojyoti Ghosh }
4842302aa1bSDebojyoti Ghosh 
485d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_GLEE(TS ts)
486d71ae5a4SJacob Faibussowitsch {
4872302aa1bSDebojyoti Ghosh   TS_GLEE       *glee = (TS_GLEE *)ts->data;
4882302aa1bSDebojyoti Ghosh   GLEETableau    tab  = glee->tableau;
4892302aa1bSDebojyoti Ghosh   const PetscInt s = tab->s, r = tab->r;
4909371c9d4SSatish Balay   PetscReal     *A = tab->A, *U = tab->U, *F = tab->F, *c = tab->c;
4919371c9d4SSatish Balay   Vec           *Y = glee->Y, *X = glee->X, *YStage = glee->YStage, *YdotStage = glee->YdotStage, W = glee->W;
4922302aa1bSDebojyoti Ghosh   SNES           snes;
493ec0e3ee2SEmil Constantinescu   PetscScalar   *ws = glee->swork, *wr = glee->rwork;
4942302aa1bSDebojyoti Ghosh   TSAdapt        adapt;
4952302aa1bSDebojyoti Ghosh   PetscInt       i, j, reject, next_scheme, its, lits;
4962302aa1bSDebojyoti Ghosh   PetscReal      next_time_step;
4972302aa1bSDebojyoti Ghosh   PetscReal      t;
4982302aa1bSDebojyoti Ghosh   PetscBool      accept;
4992302aa1bSDebojyoti Ghosh 
5002302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
5019566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(citation, &cited));
502642f3702SEmil Constantinescu 
5039566063dSJacob Faibussowitsch   for (i = 0; i < r; i++) PetscCall(VecCopy(Y[i], X[i]));
5042302aa1bSDebojyoti Ghosh 
5059566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
5062302aa1bSDebojyoti Ghosh   next_time_step = ts->time_step;
5072302aa1bSDebojyoti Ghosh   t              = ts->ptime;
5082302aa1bSDebojyoti Ghosh   accept         = PETSC_TRUE;
5092302aa1bSDebojyoti Ghosh   glee->status   = TS_STEP_INCOMPLETE;
5102302aa1bSDebojyoti Ghosh 
5112302aa1bSDebojyoti Ghosh   for (reject = 0; reject < ts->max_reject && !ts->reason; reject++, ts->reject++) {
5122302aa1bSDebojyoti Ghosh     PetscReal h = ts->time_step;
5139566063dSJacob Faibussowitsch     PetscCall(TSPreStep(ts));
5142302aa1bSDebojyoti Ghosh 
5152302aa1bSDebojyoti Ghosh     for (i = 0; i < s; i++) {
5162302aa1bSDebojyoti Ghosh       glee->stage_time = t + h * c[i];
5179566063dSJacob Faibussowitsch       PetscCall(TSPreStage(ts, glee->stage_time));
5182302aa1bSDebojyoti Ghosh 
5192302aa1bSDebojyoti Ghosh       if (A[i * s + i] == 0) { /* Explicit stage */
5209566063dSJacob Faibussowitsch         PetscCall(VecZeroEntries(YStage[i]));
521ec0e3ee2SEmil Constantinescu         for (j = 0; j < r; j++) wr[j] = U[i * r + j];
5229566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(YStage[i], r, wr, X));
523ec0e3ee2SEmil Constantinescu         for (j = 0; j < i; j++) ws[j] = h * A[i * s + j];
5249566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(YStage[i], i, ws, YdotStage));
5252302aa1bSDebojyoti Ghosh       } else { /* Implicit stage */
5262302aa1bSDebojyoti Ghosh         glee->scoeff = 1.0 / A[i * s + i];
527dd8e379bSPierre Jolivet         /* compute right-hand side */
5289566063dSJacob Faibussowitsch         PetscCall(VecZeroEntries(W));
529ec0e3ee2SEmil Constantinescu         for (j = 0; j < r; j++) wr[j] = U[i * r + j];
5309566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(W, r, wr, X));
531ec0e3ee2SEmil Constantinescu         for (j = 0; j < i; j++) ws[j] = h * A[i * s + j];
5329566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(W, i, ws, YdotStage));
5339566063dSJacob Faibussowitsch         PetscCall(VecScale(W, glee->scoeff / h));
5342302aa1bSDebojyoti Ghosh         /* set initial guess */
5359566063dSJacob Faibussowitsch         PetscCall(VecCopy(i > 0 ? YStage[i - 1] : ts->vec_sol, YStage[i]));
5362302aa1bSDebojyoti Ghosh         /* solve for this stage */
5379566063dSJacob Faibussowitsch         PetscCall(SNESSolve(snes, W, YStage[i]));
5389566063dSJacob Faibussowitsch         PetscCall(SNESGetIterationNumber(snes, &its));
5399566063dSJacob Faibussowitsch         PetscCall(SNESGetLinearSolveIterations(snes, &lits));
5409371c9d4SSatish Balay         ts->snes_its += its;
5419371c9d4SSatish Balay         ts->ksp_its += lits;
5422302aa1bSDebojyoti Ghosh       }
5439566063dSJacob Faibussowitsch       PetscCall(TSGetAdapt(ts, &adapt));
5449566063dSJacob Faibussowitsch       PetscCall(TSAdaptCheckStage(adapt, ts, glee->stage_time, YStage[i], &accept));
5452302aa1bSDebojyoti Ghosh       if (!accept) goto reject_step;
5469566063dSJacob Faibussowitsch       PetscCall(TSPostStage(ts, glee->stage_time, i, YStage));
5479566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSFunction(ts, t + h * c[i], YStage[i], YdotStage[i]));
5482302aa1bSDebojyoti Ghosh     }
5499566063dSJacob Faibussowitsch     PetscCall(TSEvaluateStep(ts, tab->order, ts->vec_sol, NULL));
5502302aa1bSDebojyoti Ghosh     glee->status = TS_STEP_PENDING;
5512302aa1bSDebojyoti Ghosh 
5522302aa1bSDebojyoti Ghosh     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
5539566063dSJacob Faibussowitsch     PetscCall(TSGetAdapt(ts, &adapt));
5549566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidatesClear(adapt));
5559566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidateAdd(adapt, tab->name, tab->order, 1, tab->ccfl, (PetscReal)tab->s, PETSC_TRUE));
5569566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, &next_scheme, &next_time_step, &accept));
5572302aa1bSDebojyoti Ghosh     if (accept) {
5582302aa1bSDebojyoti Ghosh       /* ignore next_scheme for now */
5592302aa1bSDebojyoti Ghosh       ts->ptime += ts->time_step;
5602302aa1bSDebojyoti Ghosh       ts->time_step = next_time_step;
5612302aa1bSDebojyoti Ghosh       glee->status  = TS_STEP_COMPLETE;
5620a01e1b2SEmil Constantinescu       /* compute and store the global error */
5630a01e1b2SEmil Constantinescu       /* Note: this is not needed if TSAdaptGLEE is not used */
564f4f49eeaSPierre Jolivet       PetscCall(TSGetTimeError(ts, 0, &glee->yGErr));
5659566063dSJacob Faibussowitsch       PetscCall(PetscObjectComposedDataSetReal((PetscObject)ts->vec_sol, explicit_stage_time_id, ts->ptime));
5662302aa1bSDebojyoti Ghosh       break;
5672302aa1bSDebojyoti Ghosh     } else { /* Roll back the current step */
568ec0e3ee2SEmil Constantinescu       for (j = 0; j < r; j++) wr[j] = F[j];
5699566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(ts->vec_sol, r, wr, X));
5702302aa1bSDebojyoti Ghosh       ts->time_step = next_time_step;
5712302aa1bSDebojyoti Ghosh       glee->status  = TS_STEP_INCOMPLETE;
5722302aa1bSDebojyoti Ghosh     }
5739371c9d4SSatish Balay   reject_step:
5749371c9d4SSatish Balay     continue;
5752302aa1bSDebojyoti Ghosh   }
5762302aa1bSDebojyoti Ghosh   if (glee->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
5773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5782302aa1bSDebojyoti Ghosh }
5792302aa1bSDebojyoti Ghosh 
580d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_GLEE(TS ts, PetscReal itime, Vec X)
581d71ae5a4SJacob Faibussowitsch {
5822302aa1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE *)ts->data;
5832302aa1bSDebojyoti Ghosh   PetscInt         s = glee->tableau->s, pinterp = glee->tableau->pinterp, i, j;
5842302aa1bSDebojyoti Ghosh   PetscReal        h, tt, t;
5852302aa1bSDebojyoti Ghosh   PetscScalar     *b;
5862302aa1bSDebojyoti Ghosh   const PetscReal *B = glee->tableau->binterp;
5872302aa1bSDebojyoti Ghosh 
5882302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
5893c633725SBarry Smith   PetscCheck(B, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSGLEE %s does not have an interpolation formula", glee->tableau->name);
5902302aa1bSDebojyoti Ghosh   switch (glee->status) {
5912302aa1bSDebojyoti Ghosh   case TS_STEP_INCOMPLETE:
5922302aa1bSDebojyoti Ghosh   case TS_STEP_PENDING:
5932302aa1bSDebojyoti Ghosh     h = ts->time_step;
5942302aa1bSDebojyoti Ghosh     t = (itime - ts->ptime) / h;
5952302aa1bSDebojyoti Ghosh     break;
5962302aa1bSDebojyoti Ghosh   case TS_STEP_COMPLETE:
597ed450d42SEmil Constantinescu     h = ts->ptime - ts->ptime_prev;
5982302aa1bSDebojyoti Ghosh     t = (itime - ts->ptime) / h + 1; /* In the interval [0,1] */
5992302aa1bSDebojyoti Ghosh     break;
600d71ae5a4SJacob Faibussowitsch   default:
601d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus");
6022302aa1bSDebojyoti Ghosh   }
6039566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s, &b));
6042302aa1bSDebojyoti Ghosh   for (i = 0; i < s; i++) b[i] = 0;
6052302aa1bSDebojyoti Ghosh   for (j = 0, tt = t; j < pinterp; j++, tt *= t) {
606ad540459SPierre Jolivet     for (i = 0; i < s; i++) b[i] += h * B[i * pinterp + j] * tt;
6072302aa1bSDebojyoti Ghosh   }
6089566063dSJacob Faibussowitsch   PetscCall(VecCopy(glee->YStage[0], X));
6099566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(X, s, b, glee->YdotStage));
6109566063dSJacob Faibussowitsch   PetscCall(PetscFree(b));
6113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6122302aa1bSDebojyoti Ghosh }
6132302aa1bSDebojyoti Ghosh 
6142302aa1bSDebojyoti Ghosh /*------------------------------------------------------------*/
615d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_GLEE(TS ts)
616d71ae5a4SJacob Faibussowitsch {
6172302aa1bSDebojyoti Ghosh   TS_GLEE *glee = (TS_GLEE *)ts->data;
6182302aa1bSDebojyoti Ghosh   PetscInt s, r;
6192302aa1bSDebojyoti Ghosh 
6202302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
6213ba16761SJacob Faibussowitsch   if (!glee->tableau) PetscFunctionReturn(PETSC_SUCCESS);
6222302aa1bSDebojyoti Ghosh   s = glee->tableau->s;
6232302aa1bSDebojyoti Ghosh   r = glee->tableau->r;
6249566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(r, &glee->Y));
6259566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(r, &glee->X));
6269566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(s, &glee->YStage));
6279566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(s, &glee->YdotStage));
6289566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&glee->Ydot));
6299566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&glee->yGErr));
6309566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&glee->W));
6319566063dSJacob Faibussowitsch   PetscCall(PetscFree2(glee->swork, glee->rwork));
6323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6332302aa1bSDebojyoti Ghosh }
6342302aa1bSDebojyoti Ghosh 
635d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLEEGetVecs(TS ts, DM dm, Vec *Ydot)
636d71ae5a4SJacob Faibussowitsch {
6372302aa1bSDebojyoti Ghosh   TS_GLEE *glee = (TS_GLEE *)ts->data;
6382302aa1bSDebojyoti Ghosh 
6392302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
6402302aa1bSDebojyoti Ghosh   if (Ydot) {
6412302aa1bSDebojyoti Ghosh     if (dm && dm != ts->dm) {
6429566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm, "TSGLEE_Ydot", Ydot));
6432302aa1bSDebojyoti Ghosh     } else *Ydot = glee->Ydot;
6442302aa1bSDebojyoti Ghosh   }
6453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6462302aa1bSDebojyoti Ghosh }
6472302aa1bSDebojyoti Ghosh 
648d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGLEERestoreVecs(TS ts, DM dm, Vec *Ydot)
649d71ae5a4SJacob Faibussowitsch {
6502302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
6512302aa1bSDebojyoti Ghosh   if (Ydot) {
65248a46eb9SPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSGLEE_Ydot", Ydot));
6532302aa1bSDebojyoti Ghosh   }
6543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6552302aa1bSDebojyoti Ghosh }
6562302aa1bSDebojyoti Ghosh 
6572302aa1bSDebojyoti Ghosh /*
6582302aa1bSDebojyoti Ghosh   This defines the nonlinear equation that is to be solved with SNES
6592302aa1bSDebojyoti Ghosh */
660d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_GLEE(SNES snes, Vec X, Vec F, TS ts)
661d71ae5a4SJacob Faibussowitsch {
6622302aa1bSDebojyoti Ghosh   TS_GLEE  *glee = (TS_GLEE *)ts->data;
6632302aa1bSDebojyoti Ghosh   DM        dm, dmsave;
6642302aa1bSDebojyoti Ghosh   Vec       Ydot;
6652302aa1bSDebojyoti Ghosh   PetscReal shift = glee->scoeff / ts->time_step;
6662302aa1bSDebojyoti Ghosh 
6672302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
6689566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
6699566063dSJacob Faibussowitsch   PetscCall(TSGLEEGetVecs(ts, dm, &Ydot));
6702302aa1bSDebojyoti Ghosh   /* Set Ydot = shift*X */
6719566063dSJacob Faibussowitsch   PetscCall(VecCopy(X, Ydot));
6729566063dSJacob Faibussowitsch   PetscCall(VecScale(Ydot, shift));
6732302aa1bSDebojyoti Ghosh   dmsave = ts->dm;
6742302aa1bSDebojyoti Ghosh   ts->dm = dm;
6752302aa1bSDebojyoti Ghosh 
6769566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts, glee->stage_time, X, Ydot, F, PETSC_FALSE));
6772302aa1bSDebojyoti Ghosh 
6782302aa1bSDebojyoti Ghosh   ts->dm = dmsave;
6799566063dSJacob Faibussowitsch   PetscCall(TSGLEERestoreVecs(ts, dm, &Ydot));
6803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6812302aa1bSDebojyoti Ghosh }
6822302aa1bSDebojyoti Ghosh 
683d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_GLEE(SNES snes, Vec X, Mat A, Mat B, TS ts)
684d71ae5a4SJacob Faibussowitsch {
6852302aa1bSDebojyoti Ghosh   TS_GLEE  *glee = (TS_GLEE *)ts->data;
6862302aa1bSDebojyoti Ghosh   DM        dm, dmsave;
6872302aa1bSDebojyoti Ghosh   Vec       Ydot;
6882302aa1bSDebojyoti Ghosh   PetscReal shift = glee->scoeff / ts->time_step;
6892302aa1bSDebojyoti Ghosh 
6902302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
6919566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
6929566063dSJacob Faibussowitsch   PetscCall(TSGLEEGetVecs(ts, dm, &Ydot));
6932302aa1bSDebojyoti Ghosh   /* glee->Ydot has already been computed in SNESTSFormFunction_GLEE (SNES guarantees this) */
6942302aa1bSDebojyoti Ghosh   dmsave = ts->dm;
6952302aa1bSDebojyoti Ghosh   ts->dm = dm;
6962302aa1bSDebojyoti Ghosh 
6979566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts, glee->stage_time, X, Ydot, shift, A, B, PETSC_FALSE));
6982302aa1bSDebojyoti Ghosh 
6992302aa1bSDebojyoti Ghosh   ts->dm = dmsave;
7009566063dSJacob Faibussowitsch   PetscCall(TSGLEERestoreVecs(ts, dm, &Ydot));
7013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7022302aa1bSDebojyoti Ghosh }
7032302aa1bSDebojyoti Ghosh 
704d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCoarsenHook_TSGLEE(DM fine, DM coarse, void *ctx)
705d71ae5a4SJacob Faibussowitsch {
7062302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7082302aa1bSDebojyoti Ghosh }
7092302aa1bSDebojyoti Ghosh 
710d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRestrictHook_TSGLEE(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx)
711d71ae5a4SJacob Faibussowitsch {
7122302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7142302aa1bSDebojyoti Ghosh }
7152302aa1bSDebojyoti Ghosh 
716d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainHook_TSGLEE(DM dm, DM subdm, void *ctx)
717d71ae5a4SJacob Faibussowitsch {
7182302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7202302aa1bSDebojyoti Ghosh }
7212302aa1bSDebojyoti Ghosh 
722d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainRestrictHook_TSGLEE(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx)
723d71ae5a4SJacob Faibussowitsch {
7242302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7262302aa1bSDebojyoti Ghosh }
7272302aa1bSDebojyoti Ghosh 
728d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_GLEE(TS ts)
729d71ae5a4SJacob Faibussowitsch {
7302302aa1bSDebojyoti Ghosh   TS_GLEE    *glee = (TS_GLEE *)ts->data;
7312302aa1bSDebojyoti Ghosh   GLEETableau tab;
7322302aa1bSDebojyoti Ghosh   PetscInt    s, r;
7332302aa1bSDebojyoti Ghosh   DM          dm;
7342302aa1bSDebojyoti Ghosh 
7352302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
73648a46eb9SPierre Jolivet   if (!glee->tableau) PetscCall(TSGLEESetType(ts, TSGLEEDefaultType));
7372302aa1bSDebojyoti Ghosh   tab = glee->tableau;
7382302aa1bSDebojyoti Ghosh   s   = tab->s;
7392302aa1bSDebojyoti Ghosh   r   = tab->r;
7409566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vec_sol, r, &glee->Y));
7419566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vec_sol, r, &glee->X));
7429566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vec_sol, s, &glee->YStage));
7439566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vec_sol, s, &glee->YdotStage));
7449566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &glee->Ydot));
7459566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &glee->yGErr));
7469566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(glee->yGErr));
7479566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &glee->W));
7489566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(s, &glee->swork, r, &glee->rwork));
7499566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
7509566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSGLEE, DMRestrictHook_TSGLEE, ts));
7519566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSGLEE, DMSubDomainRestrictHook_TSGLEE, ts));
7523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7532302aa1bSDebojyoti Ghosh }
7542302aa1bSDebojyoti Ghosh 
75566976f2fSJacob Faibussowitsch static PetscErrorCode TSStartingMethod_GLEE(TS ts)
756d71ae5a4SJacob Faibussowitsch {
7572302aa1bSDebojyoti Ghosh   TS_GLEE    *glee = (TS_GLEE *)ts->data;
7586e2e232bSDebojyoti Ghosh   GLEETableau tab  = glee->tableau;
7596e2e232bSDebojyoti Ghosh   PetscInt    r    = tab->r, i;
7606e2e232bSDebojyoti Ghosh   PetscReal  *S    = tab->S;
7612302aa1bSDebojyoti Ghosh 
7622302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7632302aa1bSDebojyoti Ghosh   for (i = 0; i < r; i++) {
7649566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(glee->Y[i]));
7659566063dSJacob Faibussowitsch     PetscCall(VecAXPY(glee->Y[i], S[i], ts->vec_sol));
7662302aa1bSDebojyoti Ghosh   }
7673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7682302aa1bSDebojyoti Ghosh }
7692302aa1bSDebojyoti Ghosh 
7702302aa1bSDebojyoti Ghosh /*------------------------------------------------------------*/
7712302aa1bSDebojyoti Ghosh 
772d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_GLEE(TS ts, PetscOptionItems *PetscOptionsObject)
773d71ae5a4SJacob Faibussowitsch {
7742302aa1bSDebojyoti Ghosh   char gleetype[256];
7752302aa1bSDebojyoti Ghosh 
7762302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
777d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "GLEE ODE solver options");
7782302aa1bSDebojyoti Ghosh   {
7792302aa1bSDebojyoti Ghosh     GLEETableauLink link;
7802302aa1bSDebojyoti Ghosh     PetscInt        count, choice;
7812302aa1bSDebojyoti Ghosh     PetscBool       flg;
7822302aa1bSDebojyoti Ghosh     const char    **namelist;
7832302aa1bSDebojyoti Ghosh 
7849566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(gleetype, TSGLEEDefaultType, sizeof(gleetype)));
785fbccb6d4SPierre Jolivet     for (link = GLEETableauList, count = 0; link; link = link->next, count++);
7869566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(count, (char ***)&namelist));
7872302aa1bSDebojyoti Ghosh     for (link = GLEETableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name;
7889566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEList("-ts_glee_type", "Family of GLEE method", "TSGLEESetType", (const char *const *)namelist, count, gleetype, &choice, &flg));
7899566063dSJacob Faibussowitsch     PetscCall(TSGLEESetType(ts, flg ? namelist[choice] : gleetype));
7909566063dSJacob Faibussowitsch     PetscCall(PetscFree(namelist));
7912302aa1bSDebojyoti Ghosh   }
792d0609cedSBarry Smith   PetscOptionsHeadEnd();
7933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7942302aa1bSDebojyoti Ghosh }
7952302aa1bSDebojyoti Ghosh 
796d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_GLEE(TS ts, PetscViewer viewer)
797d71ae5a4SJacob Faibussowitsch {
7982302aa1bSDebojyoti Ghosh   TS_GLEE    *glee = (TS_GLEE *)ts->data;
7992302aa1bSDebojyoti Ghosh   GLEETableau tab  = glee->tableau;
8002302aa1bSDebojyoti Ghosh   PetscBool   iascii;
8012302aa1bSDebojyoti Ghosh 
8022302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8039566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
8042302aa1bSDebojyoti Ghosh   if (iascii) {
8052302aa1bSDebojyoti Ghosh     TSGLEEType gleetype;
8062302aa1bSDebojyoti Ghosh     char       buf[512];
8079566063dSJacob Faibussowitsch     PetscCall(TSGLEEGetType(ts, &gleetype));
8089566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  GLEE type %s\n", gleetype));
8099566063dSJacob Faibussowitsch     PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, tab->c));
8109566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Abscissa     c = %s\n", buf));
8112302aa1bSDebojyoti Ghosh     /* Note: print out r as well */
8122302aa1bSDebojyoti Ghosh   }
8133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8142302aa1bSDebojyoti Ghosh }
8152302aa1bSDebojyoti Ghosh 
816d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSLoad_GLEE(TS ts, PetscViewer viewer)
817d71ae5a4SJacob Faibussowitsch {
8182302aa1bSDebojyoti Ghosh   SNES    snes;
8192302aa1bSDebojyoti Ghosh   TSAdapt tsadapt;
8202302aa1bSDebojyoti Ghosh 
8212302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8229566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &tsadapt));
8239566063dSJacob Faibussowitsch   PetscCall(TSAdaptLoad(tsadapt, viewer));
8249566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
8259566063dSJacob Faibussowitsch   PetscCall(SNESLoad(snes, viewer));
8262302aa1bSDebojyoti Ghosh   /* function and Jacobian context for SNES when used with TS is always ts object */
8279566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, NULL, NULL, ts));
8289566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, ts));
8293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8302302aa1bSDebojyoti Ghosh }
8312302aa1bSDebojyoti Ghosh 
832*cc4c1da9SBarry Smith /*@
833bcf0153eSBarry Smith   TSGLEESetType - Set the type of `TSGLEE` scheme
8342302aa1bSDebojyoti Ghosh 
83520f4b53cSBarry Smith   Logically Collective
8362302aa1bSDebojyoti Ghosh 
837d8d19677SJose E. Roman   Input Parameters:
8382302aa1bSDebojyoti Ghosh + ts       - timestepping context
839bcf0153eSBarry Smith - gleetype - type of `TSGLEE` scheme
8402302aa1bSDebojyoti Ghosh 
8412302aa1bSDebojyoti Ghosh   Level: intermediate
8422302aa1bSDebojyoti Ghosh 
8431cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEEGetType()`, `TSGLEE`
8442302aa1bSDebojyoti Ghosh @*/
845d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLEESetType(TS ts, TSGLEEType gleetype)
846d71ae5a4SJacob Faibussowitsch {
8472302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8482302aa1bSDebojyoti Ghosh   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
8494f572ea9SToby Isaac   PetscAssertPointer(gleetype, 2);
850cac4c232SBarry Smith   PetscTryMethod(ts, "TSGLEESetType_C", (TS, TSGLEEType), (ts, gleetype));
8513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8522302aa1bSDebojyoti Ghosh }
8532302aa1bSDebojyoti Ghosh 
854*cc4c1da9SBarry Smith /*@
855bcf0153eSBarry Smith   TSGLEEGetType - Get the type of `TSGLEE` scheme
8562302aa1bSDebojyoti Ghosh 
85720f4b53cSBarry Smith   Logically Collective
8582302aa1bSDebojyoti Ghosh 
8592302aa1bSDebojyoti Ghosh   Input Parameter:
8602302aa1bSDebojyoti Ghosh . ts - timestepping context
8612302aa1bSDebojyoti Ghosh 
8622302aa1bSDebojyoti Ghosh   Output Parameter:
863bcf0153eSBarry Smith . gleetype - type of `TSGLEE` scheme
8642302aa1bSDebojyoti Ghosh 
8652302aa1bSDebojyoti Ghosh   Level: intermediate
8662302aa1bSDebojyoti Ghosh 
8671cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEE`, `TSGLEESetType()`
8682302aa1bSDebojyoti Ghosh @*/
869d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGLEEGetType(TS ts, TSGLEEType *gleetype)
870d71ae5a4SJacob Faibussowitsch {
8712302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8722302aa1bSDebojyoti Ghosh   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
873cac4c232SBarry Smith   PetscUseMethod(ts, "TSGLEEGetType_C", (TS, TSGLEEType *), (ts, gleetype));
8743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8752302aa1bSDebojyoti Ghosh }
8762302aa1bSDebojyoti Ghosh 
87766976f2fSJacob Faibussowitsch static PetscErrorCode TSGLEEGetType_GLEE(TS ts, TSGLEEType *gleetype)
878d71ae5a4SJacob Faibussowitsch {
8792302aa1bSDebojyoti Ghosh   TS_GLEE *glee = (TS_GLEE *)ts->data;
8802302aa1bSDebojyoti Ghosh 
8812302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
88248a46eb9SPierre Jolivet   if (!glee->tableau) PetscCall(TSGLEESetType(ts, TSGLEEDefaultType));
8832302aa1bSDebojyoti Ghosh   *gleetype = glee->tableau->name;
8843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8852302aa1bSDebojyoti Ghosh }
88666976f2fSJacob Faibussowitsch static PetscErrorCode TSGLEESetType_GLEE(TS ts, TSGLEEType gleetype)
887d71ae5a4SJacob Faibussowitsch {
8882302aa1bSDebojyoti Ghosh   TS_GLEE        *glee = (TS_GLEE *)ts->data;
8892302aa1bSDebojyoti Ghosh   PetscBool       match;
8902302aa1bSDebojyoti Ghosh   GLEETableauLink link;
8912302aa1bSDebojyoti Ghosh 
8922302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8932302aa1bSDebojyoti Ghosh   if (glee->tableau) {
8949566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(glee->tableau->name, gleetype, &match));
8953ba16761SJacob Faibussowitsch     if (match) PetscFunctionReturn(PETSC_SUCCESS);
8962302aa1bSDebojyoti Ghosh   }
8972302aa1bSDebojyoti Ghosh   for (link = GLEETableauList; link; link = link->next) {
8989566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(link->tab.name, gleetype, &match));
8992302aa1bSDebojyoti Ghosh     if (match) {
9009566063dSJacob Faibussowitsch       PetscCall(TSReset_GLEE(ts));
9012302aa1bSDebojyoti Ghosh       glee->tableau = &link->tab;
9023ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
9032302aa1bSDebojyoti Ghosh     }
9042302aa1bSDebojyoti Ghosh   }
90598921bdaSJacob Faibussowitsch   SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", gleetype);
9062302aa1bSDebojyoti Ghosh }
9072302aa1bSDebojyoti Ghosh 
908d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGetStages_GLEE(TS ts, PetscInt *ns, Vec **Y)
909d71ae5a4SJacob Faibussowitsch {
9102302aa1bSDebojyoti Ghosh   TS_GLEE *glee = (TS_GLEE *)ts->data;
9112302aa1bSDebojyoti Ghosh 
9122302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
9130429704eSStefano Zampini   if (ns) *ns = glee->tableau->s;
914d5509560SEmil Constantinescu   if (Y) *Y = glee->YStage;
9153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9162302aa1bSDebojyoti Ghosh }
9172302aa1bSDebojyoti Ghosh 
91866976f2fSJacob Faibussowitsch static PetscErrorCode TSGetSolutionComponents_GLEE(TS ts, PetscInt *n, Vec *Y)
919d71ae5a4SJacob Faibussowitsch {
9202302aa1bSDebojyoti Ghosh   TS_GLEE    *glee = (TS_GLEE *)ts->data;
9212302aa1bSDebojyoti Ghosh   GLEETableau tab  = glee->tableau;
9222302aa1bSDebojyoti Ghosh 
9232302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
9244cdd57e5SDebojyoti Ghosh   if (!Y) *n = tab->r;
9252302aa1bSDebojyoti Ghosh   else {
9264cdd57e5SDebojyoti Ghosh     if ((*n >= 0) && (*n < tab->r)) {
9279566063dSJacob Faibussowitsch       PetscCall(VecCopy(glee->Y[*n], *Y));
92863a3b9bcSJacob Faibussowitsch     } else SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Second argument (%" PetscInt_FMT ") out of range[0,%" PetscInt_FMT "].", *n, tab->r - 1);
9292302aa1bSDebojyoti Ghosh   }
9303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9312302aa1bSDebojyoti Ghosh }
9322302aa1bSDebojyoti Ghosh 
93366976f2fSJacob Faibussowitsch static PetscErrorCode TSGetAuxSolution_GLEE(TS ts, Vec *X)
934d71ae5a4SJacob Faibussowitsch {
9354cdd57e5SDebojyoti Ghosh   TS_GLEE     *glee = (TS_GLEE *)ts->data;
9364cdd57e5SDebojyoti Ghosh   GLEETableau  tab  = glee->tableau;
9374cdd57e5SDebojyoti Ghosh   PetscReal   *F    = tab->Fembed;
9384cdd57e5SDebojyoti Ghosh   PetscInt     r    = tab->r;
9394cdd57e5SDebojyoti Ghosh   Vec         *Y    = glee->Y;
940ec0e3ee2SEmil Constantinescu   PetscScalar *wr   = glee->rwork;
941ec0e3ee2SEmil Constantinescu   PetscInt     i;
9424cdd57e5SDebojyoti Ghosh 
9434cdd57e5SDebojyoti Ghosh   PetscFunctionBegin;
9449566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(*X));
945ec0e3ee2SEmil Constantinescu   for (i = 0; i < r; i++) wr[i] = F[i];
946f4f49eeaSPierre Jolivet   PetscCall(VecMAXPY(*X, r, wr, Y));
9473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9484cdd57e5SDebojyoti Ghosh }
9494cdd57e5SDebojyoti Ghosh 
95066976f2fSJacob Faibussowitsch static PetscErrorCode TSGetTimeError_GLEE(TS ts, PetscInt n, Vec *X)
951d71ae5a4SJacob Faibussowitsch {
9524cdd57e5SDebojyoti Ghosh   TS_GLEE     *glee = (TS_GLEE *)ts->data;
9534cdd57e5SDebojyoti Ghosh   GLEETableau  tab  = glee->tableau;
9544cdd57e5SDebojyoti Ghosh   PetscReal   *F    = tab->Ferror;
9554cdd57e5SDebojyoti Ghosh   PetscInt     r    = tab->r;
9564cdd57e5SDebojyoti Ghosh   Vec         *Y    = glee->Y;
957ec0e3ee2SEmil Constantinescu   PetscScalar *wr   = glee->rwork;
958ec0e3ee2SEmil Constantinescu   PetscInt     i;
9594cdd57e5SDebojyoti Ghosh 
9604cdd57e5SDebojyoti Ghosh   PetscFunctionBegin;
9619566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(*X));
962657c1e31SEmil Constantinescu   if (n == 0) {
963ec0e3ee2SEmil Constantinescu     for (i = 0; i < r; i++) wr[i] = F[i];
964f4f49eeaSPierre Jolivet     PetscCall(VecMAXPY(*X, r, wr, Y));
9650a01e1b2SEmil Constantinescu   } else if (n == -1) {
966657c1e31SEmil Constantinescu     *X = glee->yGErr;
9670a01e1b2SEmil Constantinescu   }
9683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9694cdd57e5SDebojyoti Ghosh }
9704cdd57e5SDebojyoti Ghosh 
97166976f2fSJacob Faibussowitsch static PetscErrorCode TSSetTimeError_GLEE(TS ts, Vec X)
972d71ae5a4SJacob Faibussowitsch {
97357df6a1bSDebojyoti Ghosh   TS_GLEE    *glee = (TS_GLEE *)ts->data;
97457df6a1bSDebojyoti Ghosh   GLEETableau tab  = glee->tableau;
97557df6a1bSDebojyoti Ghosh   PetscReal  *S    = tab->Serror;
97657df6a1bSDebojyoti Ghosh   PetscInt    r    = tab->r, i;
97757df6a1bSDebojyoti Ghosh   Vec        *Y    = glee->Y;
97857df6a1bSDebojyoti Ghosh 
97957df6a1bSDebojyoti Ghosh   PetscFunctionBegin;
98063a3b9bcSJacob Faibussowitsch   PetscCheck(r == 2, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSSetTimeError_GLEE not supported for '%s' with r=%" PetscInt_FMT ".", tab->name, tab->r);
98157df6a1bSDebojyoti Ghosh   for (i = 1; i < r; i++) {
9829566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vec_sol, Y[i]));
9839566063dSJacob Faibussowitsch     PetscCall(VecAXPBY(Y[i], S[0], S[1], X));
9849566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, glee->yGErr));
98557df6a1bSDebojyoti Ghosh   }
9863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
98757df6a1bSDebojyoti Ghosh }
98857df6a1bSDebojyoti Ghosh 
989d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_GLEE(TS ts)
990d71ae5a4SJacob Faibussowitsch {
991b3a6b972SJed Brown   PetscFunctionBegin;
9929566063dSJacob Faibussowitsch   PetscCall(TSReset_GLEE(ts));
993b3a6b972SJed Brown   if (ts->dm) {
9949566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSGLEE, DMRestrictHook_TSGLEE, ts));
9959566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSGLEE, DMSubDomainRestrictHook_TSGLEE, ts));
996b3a6b972SJed Brown   }
9979566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
9989566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLEEGetType_C", NULL));
9999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLEESetType_C", NULL));
10003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1001b3a6b972SJed Brown }
1002b3a6b972SJed Brown 
10032302aa1bSDebojyoti Ghosh /* ------------------------------------------------------------ */
10042302aa1bSDebojyoti Ghosh /*MC
10052302aa1bSDebojyoti Ghosh       TSGLEE - ODE and DAE solver using General Linear with Error Estimation schemes
10062302aa1bSDebojyoti Ghosh 
1007dd8e379bSPierre Jolivet   The user should provide the right-hand side of the equation using `TSSetRHSFunction()`.
10082302aa1bSDebojyoti Ghosh 
10092302aa1bSDebojyoti Ghosh   Level: beginner
10102302aa1bSDebojyoti Ghosh 
1011bcf0153eSBarry Smith   Note:
1012bcf0153eSBarry Smith   The default is `TSGLEE35`, it can be changed with `TSGLEESetType()` or -ts_glee_type
10132302aa1bSDebojyoti Ghosh 
10141cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSGLEESetType()`, `TSGLEEGetType()`,
1015bcf0153eSBarry Smith           `TSGLEE23`, `TTSGLEE24`, `TSGLEE35`, `TSGLEE25I`, `TSGLEEEXRK2A`,
1016bcf0153eSBarry Smith           `TSGLEERK32G1`, `TSGLEERK285EX`, `TSGLEEType`, `TSGLEERegister()`, `TSType`
10172302aa1bSDebojyoti Ghosh M*/
1018d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_GLEE(TS ts)
1019d71ae5a4SJacob Faibussowitsch {
10202302aa1bSDebojyoti Ghosh   TS_GLEE *th;
10212302aa1bSDebojyoti Ghosh 
10222302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
10239566063dSJacob Faibussowitsch   PetscCall(TSGLEEInitializePackage());
10242302aa1bSDebojyoti Ghosh 
10252302aa1bSDebojyoti Ghosh   ts->ops->reset                 = TSReset_GLEE;
10262302aa1bSDebojyoti Ghosh   ts->ops->destroy               = TSDestroy_GLEE;
10272302aa1bSDebojyoti Ghosh   ts->ops->view                  = TSView_GLEE;
10282302aa1bSDebojyoti Ghosh   ts->ops->load                  = TSLoad_GLEE;
10292302aa1bSDebojyoti Ghosh   ts->ops->setup                 = TSSetUp_GLEE;
10302302aa1bSDebojyoti Ghosh   ts->ops->step                  = TSStep_GLEE;
10312302aa1bSDebojyoti Ghosh   ts->ops->interpolate           = TSInterpolate_GLEE;
10322302aa1bSDebojyoti Ghosh   ts->ops->evaluatestep          = TSEvaluateStep_GLEE;
10332302aa1bSDebojyoti Ghosh   ts->ops->setfromoptions        = TSSetFromOptions_GLEE;
10342302aa1bSDebojyoti Ghosh   ts->ops->getstages             = TSGetStages_GLEE;
10352302aa1bSDebojyoti Ghosh   ts->ops->snesfunction          = SNESTSFormFunction_GLEE;
10362302aa1bSDebojyoti Ghosh   ts->ops->snesjacobian          = SNESTSFormJacobian_GLEE;
1037b2bf4f3aSDebojyoti Ghosh   ts->ops->getsolutioncomponents = TSGetSolutionComponents_GLEE;
10384cdd57e5SDebojyoti Ghosh   ts->ops->getauxsolution        = TSGetAuxSolution_GLEE;
10394cdd57e5SDebojyoti Ghosh   ts->ops->gettimeerror          = TSGetTimeError_GLEE;
104057df6a1bSDebojyoti Ghosh   ts->ops->settimeerror          = TSSetTimeError_GLEE;
10412302aa1bSDebojyoti Ghosh   ts->ops->startingmethod        = TSStartingMethod_GLEE;
1042b80d5313SLisandro Dalcin   ts->default_adapt_type         = TSADAPTGLEE;
10432302aa1bSDebojyoti Ghosh 
1044efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
1045efd4aadfSBarry Smith 
10464dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&th));
10472302aa1bSDebojyoti Ghosh   ts->data = (void *)th;
10482302aa1bSDebojyoti Ghosh 
10499566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLEEGetType_C", TSGLEEGetType_GLEE));
10509566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLEESetType_C", TSGLEESetType_GLEE));
10513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10522302aa1bSDebojyoti Ghosh }
1053