xref: /petsc/src/ts/impls/glee/glee.c (revision 0429704e04bf0203988aae2c826167b2b299adfb)
12302aa1bSDebojyoti Ghosh /*
22302aa1bSDebojyoti Ghosh   Code for time stepping with the General Linear with Error Estimation method
32302aa1bSDebojyoti Ghosh 
4642f3702SEmil Constantinescu 
52302aa1bSDebojyoti Ghosh   Notes:
62302aa1bSDebojyoti Ghosh   The general system is written as
72302aa1bSDebojyoti Ghosh 
82302aa1bSDebojyoti Ghosh   Udot = F(t,U)
92302aa1bSDebojyoti Ghosh 
102302aa1bSDebojyoti Ghosh */
112302aa1bSDebojyoti Ghosh #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
122302aa1bSDebojyoti Ghosh #include <petscdm.h>
132302aa1bSDebojyoti Ghosh 
14642f3702SEmil Constantinescu static PetscBool  cited = PETSC_FALSE;
15642f3702SEmil Constantinescu static const char citation[] =
16642f3702SEmil Constantinescu   "@ARTICLE{Constantinescu_TR2016b,\n"
17642f3702SEmil Constantinescu   " author = {Constantinescu, E.M.},\n"
18642f3702SEmil Constantinescu   " title = {Estimating Global Errors in Time Stepping},\n"
19642f3702SEmil Constantinescu   " journal = {ArXiv e-prints},\n"
20642f3702SEmil Constantinescu   " year = 2016,\n"
21642f3702SEmil Constantinescu   " adsurl = {http://adsabs.harvard.edu/abs/2015arXiv150305166C}\n}\n";
22642f3702SEmil Constantinescu 
232302aa1bSDebojyoti Ghosh static TSGLEEType   TSGLEEDefaultType = TSGLEE35;
242302aa1bSDebojyoti Ghosh static PetscBool    TSGLEERegisterAllCalled;
252302aa1bSDebojyoti Ghosh static PetscBool    TSGLEEPackageInitialized;
262302aa1bSDebojyoti Ghosh static PetscInt     explicit_stage_time_id;
272302aa1bSDebojyoti Ghosh 
282302aa1bSDebojyoti Ghosh typedef struct _GLEETableau *GLEETableau;
292302aa1bSDebojyoti Ghosh struct _GLEETableau {
302302aa1bSDebojyoti Ghosh   char      *name;
312302aa1bSDebojyoti Ghosh   PetscInt   order;               /* Classical approximation order of the method i*/
322302aa1bSDebojyoti Ghosh   PetscInt   s;                   /* Number of stages */
332302aa1bSDebojyoti Ghosh   PetscInt   r;                   /* Number of steps */
342302aa1bSDebojyoti Ghosh   PetscReal  gamma;               /* LTE ratio */
352302aa1bSDebojyoti Ghosh   PetscReal *A,*B,*U,*V,*S,*F,*c; /* Tableau */
362302aa1bSDebojyoti Ghosh   PetscReal *Fembed;              /* Embedded final method coefficients */
37fd96ebc2SDebojyoti Ghosh   PetscReal *Ferror;              /* Coefficients for computing error   */
3857df6a1bSDebojyoti Ghosh   PetscReal *Serror;              /* Coefficients for initializing the error   */
392302aa1bSDebojyoti Ghosh   PetscInt   pinterp;             /* Interpolation order */
402302aa1bSDebojyoti Ghosh   PetscReal *binterp;             /* Interpolation coefficients */
412302aa1bSDebojyoti Ghosh   PetscReal  ccfl;                /* Placeholder for CFL coefficient relative to forward Euler  */
422302aa1bSDebojyoti Ghosh };
432302aa1bSDebojyoti Ghosh typedef struct _GLEETableauLink *GLEETableauLink;
442302aa1bSDebojyoti Ghosh struct _GLEETableauLink {
452302aa1bSDebojyoti Ghosh   struct _GLEETableau tab;
462302aa1bSDebojyoti Ghosh   GLEETableauLink     next;
472302aa1bSDebojyoti Ghosh };
482302aa1bSDebojyoti Ghosh static GLEETableauLink GLEETableauList;
492302aa1bSDebojyoti Ghosh 
502302aa1bSDebojyoti Ghosh typedef struct {
512302aa1bSDebojyoti Ghosh   GLEETableau  tableau;
522302aa1bSDebojyoti Ghosh   Vec          *Y;         /* Solution vector (along with auxiliary solution y~ or eps) */
532302aa1bSDebojyoti Ghosh   Vec          *X;         /* Temporary solution vector */
542302aa1bSDebojyoti Ghosh   Vec          *YStage;    /* Stage values */
552302aa1bSDebojyoti Ghosh   Vec          *YdotStage; /* Stage right hand side */
562302aa1bSDebojyoti Ghosh   Vec          W;          /* Right-hand-side for implicit stage solve */
572302aa1bSDebojyoti Ghosh   Vec          Ydot;       /* Work vector holding Ydot during residual evaluation */
580a01e1b2SEmil Constantinescu   Vec          yGErr;      /* Vector holding the global error after a step is completed */
59ec0e3ee2SEmil Constantinescu   PetscScalar  *swork;     /* Scalar work (size of the number of stages)*/
60ec0e3ee2SEmil Constantinescu   PetscScalar  *rwork;     /* Scalar work (size of the number of steps)*/
612302aa1bSDebojyoti Ghosh   PetscReal    scoeff;     /* shift = scoeff/dt */
622302aa1bSDebojyoti Ghosh   PetscReal    stage_time;
632302aa1bSDebojyoti Ghosh   TSStepStatus status;
642302aa1bSDebojyoti Ghosh } TS_GLEE;
652302aa1bSDebojyoti Ghosh 
662302aa1bSDebojyoti Ghosh /*MC
672302aa1bSDebojyoti Ghosh      TSGLEE23 - Second order three stage GLEE method
682302aa1bSDebojyoti Ghosh 
692302aa1bSDebojyoti Ghosh      This method has three stages.
702302aa1bSDebojyoti Ghosh      s = 3, r = 2
712302aa1bSDebojyoti Ghosh 
722302aa1bSDebojyoti Ghosh      Level: advanced
732302aa1bSDebojyoti Ghosh 
742302aa1bSDebojyoti Ghosh .seealso: TSGLEE
7536c34d14SEmil Constantinescu M*/
762302aa1bSDebojyoti Ghosh /*MC
772302aa1bSDebojyoti Ghosh      TSGLEE24 - Second order four stage GLEE method
782302aa1bSDebojyoti Ghosh 
792302aa1bSDebojyoti Ghosh      This method has four stages.
802302aa1bSDebojyoti Ghosh      s = 4, r = 2
812302aa1bSDebojyoti Ghosh 
822302aa1bSDebojyoti Ghosh      Level: advanced
832302aa1bSDebojyoti Ghosh 
842302aa1bSDebojyoti Ghosh .seealso: TSGLEE
8536c34d14SEmil Constantinescu M*/
862302aa1bSDebojyoti Ghosh /*MC
872302aa1bSDebojyoti Ghosh      TSGLEE25i - Second order five stage GLEE method
882302aa1bSDebojyoti Ghosh 
892302aa1bSDebojyoti Ghosh      This method has five stages.
902302aa1bSDebojyoti Ghosh      s = 5, r = 2
912302aa1bSDebojyoti Ghosh 
922302aa1bSDebojyoti Ghosh      Level: advanced
932302aa1bSDebojyoti Ghosh 
942302aa1bSDebojyoti Ghosh .seealso: TSGLEE
9536c34d14SEmil Constantinescu M*/
962302aa1bSDebojyoti Ghosh /*MC
972302aa1bSDebojyoti Ghosh      TSGLEE35  - Third order five stage GLEE method
982302aa1bSDebojyoti Ghosh 
992302aa1bSDebojyoti Ghosh      This method has five stages.
1002302aa1bSDebojyoti Ghosh      s = 5, r = 2
1012302aa1bSDebojyoti Ghosh 
1022302aa1bSDebojyoti Ghosh      Level: advanced
1032302aa1bSDebojyoti Ghosh 
1042302aa1bSDebojyoti Ghosh .seealso: TSGLEE
10536c34d14SEmil Constantinescu M*/
1062302aa1bSDebojyoti Ghosh /*MC
1072302aa1bSDebojyoti Ghosh      TSGLEEEXRK2A  - Second order six stage GLEE method
1082302aa1bSDebojyoti Ghosh 
1092302aa1bSDebojyoti Ghosh      This method has six stages.
1102302aa1bSDebojyoti Ghosh      s = 6, r = 2
1112302aa1bSDebojyoti Ghosh 
1122302aa1bSDebojyoti Ghosh      Level: advanced
1132302aa1bSDebojyoti Ghosh 
1142302aa1bSDebojyoti Ghosh .seealso: TSGLEE
11536c34d14SEmil Constantinescu M*/
1162302aa1bSDebojyoti Ghosh /*MC
1172302aa1bSDebojyoti Ghosh      TSGLEERK32G1  - Third order eight stage GLEE method
1182302aa1bSDebojyoti Ghosh 
1192302aa1bSDebojyoti Ghosh      This method has eight stages.
1202302aa1bSDebojyoti Ghosh      s = 8, r = 2
1212302aa1bSDebojyoti Ghosh 
1222302aa1bSDebojyoti Ghosh      Level: advanced
1232302aa1bSDebojyoti Ghosh 
1242302aa1bSDebojyoti Ghosh .seealso: TSGLEE
12536c34d14SEmil Constantinescu M*/
1262302aa1bSDebojyoti Ghosh /*MC
1272302aa1bSDebojyoti Ghosh      TSGLEERK285EX  - Second order nine stage GLEE method
1282302aa1bSDebojyoti Ghosh 
1292302aa1bSDebojyoti Ghosh      This method has nine stages.
1302302aa1bSDebojyoti Ghosh      s = 9, r = 2
1312302aa1bSDebojyoti Ghosh 
1322302aa1bSDebojyoti Ghosh      Level: advanced
1332302aa1bSDebojyoti Ghosh 
1342302aa1bSDebojyoti Ghosh .seealso: TSGLEE
13536c34d14SEmil Constantinescu M*/
1362302aa1bSDebojyoti Ghosh 
1372302aa1bSDebojyoti Ghosh /*@C
1382302aa1bSDebojyoti Ghosh   TSGLEERegisterAll - Registers all of the General Linear with Error Estimation methods in TSGLEE
1392302aa1bSDebojyoti Ghosh 
1402302aa1bSDebojyoti Ghosh   Not Collective, but should be called by all processes which will need the schemes to be registered
1412302aa1bSDebojyoti Ghosh 
1422302aa1bSDebojyoti Ghosh   Level: advanced
1432302aa1bSDebojyoti Ghosh 
1442302aa1bSDebojyoti Ghosh .keywords: TS, TSGLEE, register, all
1452302aa1bSDebojyoti Ghosh 
1462302aa1bSDebojyoti Ghosh .seealso:  TSGLEERegisterDestroy()
1472302aa1bSDebojyoti Ghosh @*/
1482302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEERegisterAll(void)
1492302aa1bSDebojyoti Ghosh {
1502302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
1512302aa1bSDebojyoti Ghosh 
1522302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
1532302aa1bSDebojyoti Ghosh   if (TSGLEERegisterAllCalled) PetscFunctionReturn(0);
1542302aa1bSDebojyoti Ghosh   TSGLEERegisterAllCalled = PETSC_TRUE;
1552302aa1bSDebojyoti Ghosh 
1562302aa1bSDebojyoti Ghosh   {
157b1fbfd33SSatish Balay #define GAMMA 0.5
1582302aa1bSDebojyoti Ghosh     /* y-eps form */
1592302aa1bSDebojyoti Ghosh     const PetscInt
160ab8c5912SEmil Constantinescu       p = 1,
161ab8c5912SEmil Constantinescu       s = 3,
162ab8c5912SEmil Constantinescu       r = 2;
163ab8c5912SEmil Constantinescu     const PetscReal
164ab8c5912SEmil Constantinescu       A[3][3]   = {{1.0,0,0},{0,0.5,0},{0,0.5,0.5}},
165ab8c5912SEmil Constantinescu       B[2][3]   = {{1.0,0,0},{-2.0,1.0,1.0}},
166ab8c5912SEmil Constantinescu       U[3][2]   = {{1.0,0},{1.0,0.5},{1.0,0.5}},
167ab8c5912SEmil Constantinescu       V[2][2]   = {{1,0},{0,1}},
168ab8c5912SEmil Constantinescu       S[2]      = {1,0},
169ab8c5912SEmil Constantinescu       F[2]      = {1,0},
170b1fbfd33SSatish Balay       Fembed[2] = {1,1-GAMMA},
17157df6a1bSDebojyoti Ghosh       Ferror[2] = {0,1},
17257df6a1bSDebojyoti Ghosh       Serror[2] = {1,0};
173b1fbfd33SSatish Balay     ierr = 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);CHKERRQ(ierr);
174ab8c5912SEmil Constantinescu   }
175ab8c5912SEmil Constantinescu   {
176b1fbfd33SSatish Balay #undef GAMMA
177b1fbfd33SSatish Balay #define GAMMA 0.0
178ab8c5912SEmil Constantinescu     /* y-eps form */
179ab8c5912SEmil Constantinescu     const PetscInt
1802302aa1bSDebojyoti Ghosh       p = 2,
1812302aa1bSDebojyoti Ghosh       s = 3,
1822302aa1bSDebojyoti Ghosh       r = 2;
1832302aa1bSDebojyoti Ghosh     const PetscReal
1842302aa1bSDebojyoti Ghosh       A[3][3]   = {{0,0,0},{1,0,0},{0.25,0.25,0}},
1852302aa1bSDebojyoti Ghosh       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}},
1862302aa1bSDebojyoti Ghosh       U[3][2]   = {{1,0},{1,10},{1,-1}},
1872302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
1882302aa1bSDebojyoti Ghosh       S[2]      = {1,0},
1892302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
190b1fbfd33SSatish Balay       Fembed[2] = {1,1-GAMMA},
19157df6a1bSDebojyoti Ghosh       Ferror[2] = {0,1},
19257df6a1bSDebojyoti Ghosh       Serror[2] = {1,0};
193b1fbfd33SSatish Balay     ierr = 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);CHKERRQ(ierr);
1942302aa1bSDebojyoti Ghosh   }
1952302aa1bSDebojyoti Ghosh   {
196b1fbfd33SSatish Balay #undef GAMMA
197b1fbfd33SSatish Balay #define GAMMA 0.0
1982302aa1bSDebojyoti Ghosh     /* y-y~ form */
1992302aa1bSDebojyoti Ghosh     const PetscInt
2002302aa1bSDebojyoti Ghosh       p = 2,
2012302aa1bSDebojyoti Ghosh       s = 4,
2022302aa1bSDebojyoti Ghosh       r = 2;
2032302aa1bSDebojyoti Ghosh     const PetscReal
2042302aa1bSDebojyoti Ghosh       A[4][4]   = {{0,0,0,0},{0.75,0,0,0},{0.25,29.0/60.0,0,0},{-21.0/44.0,145.0/44.0,-20.0/11.0,0}},
2052302aa1bSDebojyoti Ghosh       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}},
2062302aa1bSDebojyoti Ghosh       U[4][2]   = {{0,1},{75.0/58.0,-17.0/58.0},{0,1},{0,1}},
2072302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
2082302aa1bSDebojyoti Ghosh       S[2]      = {1,1},
2092302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
210fd96ebc2SDebojyoti Ghosh       Fembed[2] = {0,1},
211b1fbfd33SSatish Balay       Ferror[2] = {-1.0/(1.0-GAMMA),1.0/(1.0-GAMMA)},
212b1fbfd33SSatish Balay       Serror[2] = {1.0-GAMMA,1.0};
213b1fbfd33SSatish Balay       ierr = 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);CHKERRQ(ierr);
2142302aa1bSDebojyoti Ghosh   }
2152302aa1bSDebojyoti Ghosh   {
216b1fbfd33SSatish Balay #undef GAMMA
217b1fbfd33SSatish Balay #define GAMMA 0.0
2182302aa1bSDebojyoti Ghosh     /* y-y~ form */
2192302aa1bSDebojyoti Ghosh     const PetscInt
2202302aa1bSDebojyoti Ghosh       p = 2,
2212302aa1bSDebojyoti Ghosh       s = 5,
2222302aa1bSDebojyoti Ghosh       r = 2;
2232302aa1bSDebojyoti Ghosh     const PetscReal
2242302aa1bSDebojyoti Ghosh       A[5][5]   = {{0,0,0,0,0},
2252302aa1bSDebojyoti Ghosh                    {-0.94079244066783383269,0,0,0,0},
2262302aa1bSDebojyoti Ghosh                    {0.64228187778301907108,0.10915356933958500042,0,0,0},
2272302aa1bSDebojyoti Ghosh                    {-0.51764297742287450812,0.74414270351096040738,-0.71404164927824538121,0,0},
2282302aa1bSDebojyoti Ghosh                    {-0.44696561556825969206,-0.76768425657590196518,0.20111608138142987881,0.93828186737840469796,0}},
2292302aa1bSDebojyoti Ghosh       B[2][5]   = {{-0.029309178948150356153,-0.49671981884013874923,0.34275801517650053274,0.32941112623949194988,0.85385985637229662276},
2302302aa1bSDebojyoti Ghosh                    {0.78133219686062535272,0.074238691892675897635,0.57957363498384957966,-0.24638502829674959968,-0.18875949544040123033}},
2312302aa1bSDebojyoti Ghosh       U[5][2]   = {{0.16911424754448327735,0.83088575245551672265},
2322302aa1bSDebojyoti Ghosh                    {0.53638465733199574340,0.46361534266800425660},
2332302aa1bSDebojyoti Ghosh                    {0.39901579167169582526,0.60098420832830417474},
2342302aa1bSDebojyoti Ghosh                    {0.87689005530618575480,0.12310994469381424520},
2352302aa1bSDebojyoti Ghosh                    {0.99056100455550913009,0.0094389954444908699092}},
2362302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
2372302aa1bSDebojyoti Ghosh       S[2]      = {1,1},
2382302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
239fd96ebc2SDebojyoti Ghosh       Fembed[2] = {0,1},
240b1fbfd33SSatish Balay       Ferror[2] = {-1.0/(1.0-GAMMA),1.0/(1.0-GAMMA)},
241b1fbfd33SSatish Balay       Serror[2] = {1.0-GAMMA,1.0};
242b1fbfd33SSatish Balay     ierr = 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);CHKERRQ(ierr);
2432302aa1bSDebojyoti Ghosh   }
2442302aa1bSDebojyoti Ghosh   {
245b1fbfd33SSatish Balay #undef GAMMA
246b1fbfd33SSatish Balay #define GAMMA 0.0
2472302aa1bSDebojyoti Ghosh     /* y-y~ form */
2482302aa1bSDebojyoti Ghosh     const PetscInt
2492302aa1bSDebojyoti Ghosh       p = 3,
2502302aa1bSDebojyoti Ghosh       s = 5,
2512302aa1bSDebojyoti Ghosh       r = 2;
2522302aa1bSDebojyoti Ghosh     const PetscReal
2532302aa1bSDebojyoti Ghosh       A[5][5]   = {{0,0,0,0,0},
2542302aa1bSDebojyoti Ghosh                    {- 2169604947363702313.0 /  24313474998937147335.0,0,0,0,0},
2552302aa1bSDebojyoti Ghosh                    {46526746497697123895.0 /  94116917485856474137.0,-10297879244026594958.0 /  49199457603717988219.0,0,0,0},
2562302aa1bSDebojyoti Ghosh                    {23364788935845982499.0 /  87425311444725389446.0,-79205144337496116638.0 / 148994349441340815519.0,40051189859317443782.0 /  36487615018004984309.0,0,0},
25717d8b14cSSatish Balay                    {42089522664062539205.0 / 124911313006412840286.0,-15074384760342762939.0 / 137927286865289746282.0,-62274678522253371016.0 / 125918573676298591413.0,13755475729852471739.0 /  79257927066651693390.0,0}},
2582302aa1bSDebojyoti Ghosh       B[2][5]   = {{61546696837458703723.0 /  56982519523786160813.0,-55810892792806293355.0 / 206957624151308356511.0,24061048952676379087.0 / 158739347956038723465.0,3577972206874351339.0 /   7599733370677197135.0,-59449832954780563947.0 / 137360038685338563670.0},
2592302aa1bSDebojyoti Ghosh                    {- 9738262186984159168.0 /  99299082461487742983.0,-32797097931948613195.0 /  61521565616362163366.0,42895514606418420631.0 /  71714201188501437336.0,22608567633166065068.0 /  55371917805607957003.0,94655809487476459565.0 / 151517167160302729021.0}},
2602302aa1bSDebojyoti Ghosh       U[5][2]   = {{70820309139834661559.0 /  80863923579509469826.0,10043614439674808267.0 /  80863923579509469826.0},
2612302aa1bSDebojyoti Ghosh                    {161694774978034105510.0 / 106187653640211060371.0,-55507121337823045139.0 / 106187653640211060371.0},
2622302aa1bSDebojyoti Ghosh                    {78486094644566264568.0 /  88171030896733822981.0,9684936252167558413.0 /  88171030896733822981.0},
2632302aa1bSDebojyoti Ghosh                    {65394922146334854435.0 /  84570853840405479554.0,19175931694070625119.0 /  84570853840405479554.0},
2642302aa1bSDebojyoti Ghosh                    {8607282770183754108.0 / 108658046436496925911.0,100050763666313171803.0 / 108658046436496925911.0}},
2652302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
2662302aa1bSDebojyoti Ghosh       S[2]      = {1,1},
2672302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
268fd96ebc2SDebojyoti Ghosh       Fembed[2] = {0,1},
269b1fbfd33SSatish Balay       Ferror[2] = {-1.0/(1.0-GAMMA),1.0/(1.0-GAMMA)},
270b1fbfd33SSatish Balay       Serror[2] = {1.0-GAMMA,1.0};
271b1fbfd33SSatish Balay     ierr = 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);CHKERRQ(ierr);
2722302aa1bSDebojyoti Ghosh   }
2732302aa1bSDebojyoti Ghosh   {
274b1fbfd33SSatish Balay #undef GAMMA
275b1fbfd33SSatish Balay #define GAMMA 0.25
2762302aa1bSDebojyoti Ghosh     /* y-eps form */
2772302aa1bSDebojyoti Ghosh     const PetscInt
2782302aa1bSDebojyoti Ghosh       p = 2,
2792302aa1bSDebojyoti Ghosh       s = 6,
2802302aa1bSDebojyoti Ghosh       r = 2;
2812302aa1bSDebojyoti Ghosh     const PetscReal
2822302aa1bSDebojyoti Ghosh       A[6][6]   = {{0,0,0,0,0,0},
2832302aa1bSDebojyoti Ghosh                    {1,0,0,0,0,0},
2842302aa1bSDebojyoti Ghosh                    {0,0,0,0,0,0},
2852302aa1bSDebojyoti Ghosh                    {0,0,0.5,0,0,0},
2862302aa1bSDebojyoti Ghosh                    {0,0,0.25,0.25,0,0},
2872302aa1bSDebojyoti Ghosh                    {0,0,0.25,0.25,0.5,0}},
2882302aa1bSDebojyoti Ghosh       B[2][6]   = {{0.5,0.5,0,0,0,0},
2892302aa1bSDebojyoti Ghosh                    {-2.0/3.0,-2.0/3.0,1.0/3.0,1.0/3.0,1.0/3.0,1.0/3.0}},
2902302aa1bSDebojyoti Ghosh       U[6][2]   = {{1,0},{1,0},{1,0.75},{1,0.75},{1,0.75},{1,0.75}},
2912302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
2922302aa1bSDebojyoti Ghosh       S[2]      = {1,0},
2932302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
294b1fbfd33SSatish Balay       Fembed[2] = {1,1-GAMMA},
29557df6a1bSDebojyoti Ghosh       Ferror[2] = {0,1},
29657df6a1bSDebojyoti Ghosh       Serror[2] = {1,0};
297b1fbfd33SSatish Balay     ierr = 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);CHKERRQ(ierr);
2982302aa1bSDebojyoti Ghosh   }
2992302aa1bSDebojyoti Ghosh   {
300b1fbfd33SSatish Balay #undef GAMMA
301b1fbfd33SSatish Balay #define GAMMA 0.0
3022302aa1bSDebojyoti Ghosh     /* y-eps form */
3032302aa1bSDebojyoti Ghosh     const PetscInt
3042302aa1bSDebojyoti Ghosh       p = 3,
3052302aa1bSDebojyoti Ghosh       s = 8,
3062302aa1bSDebojyoti Ghosh       r = 2;
3072302aa1bSDebojyoti Ghosh     const PetscReal
3082302aa1bSDebojyoti Ghosh       A[8][8]   = {{0,0,0,0,0,0,0,0},
3092302aa1bSDebojyoti Ghosh                    {0.5,0,0,0,0,0,0,0},
3102302aa1bSDebojyoti Ghosh                    {-1,2,0,0,0,0,0,0},
3112302aa1bSDebojyoti Ghosh                    {1.0/6.0,2.0/3.0,1.0/6.0,0,0,0,0,0},
3122302aa1bSDebojyoti Ghosh                    {0,0,0,0,0,0,0,0},
3132302aa1bSDebojyoti Ghosh                    {-7.0/24.0,1.0/3.0,1.0/12.0,-1.0/8.0,0.5,0,0,0},
3142302aa1bSDebojyoti Ghosh                    {7.0/6.0,-4.0/3.0,-1.0/3.0,0.5,-1.0,2.0,0,0},
3152302aa1bSDebojyoti Ghosh                    {0,0,0,0,1.0/6.0,2.0/3.0,1.0/6.0,0}},
3162302aa1bSDebojyoti Ghosh       B[2][8]   = {{1.0/6.0,2.0/3.0,1.0/6.0,0,0,0,0,0},
3172302aa1bSDebojyoti Ghosh                    {-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}},
3182302aa1bSDebojyoti Ghosh       U[8][2]   = {{1,0},{1,0},{1,0},{1,0},{1,1},{1,1},{1,1},{1,1}},
3192302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
3202302aa1bSDebojyoti Ghosh       S[2]      = {1,0},
3212302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
322b1fbfd33SSatish Balay       Fembed[2] = {1,1-GAMMA},
32357df6a1bSDebojyoti Ghosh       Ferror[2] = {0,1},
32457df6a1bSDebojyoti Ghosh       Serror[2] = {1,0};
325b1fbfd33SSatish Balay     ierr = 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);CHKERRQ(ierr);
3262302aa1bSDebojyoti Ghosh   }
3272302aa1bSDebojyoti Ghosh   {
328b1fbfd33SSatish Balay #undef GAMMA
329b1fbfd33SSatish Balay #define GAMMA 0.25
3302302aa1bSDebojyoti Ghosh     /* y-eps form */
3312302aa1bSDebojyoti Ghosh     const PetscInt
3322302aa1bSDebojyoti Ghosh       p = 2,
3332302aa1bSDebojyoti Ghosh       s = 9,
3342302aa1bSDebojyoti Ghosh       r = 2;
3352302aa1bSDebojyoti Ghosh     const PetscReal
3362302aa1bSDebojyoti Ghosh       A[9][9]   = {{0,0,0,0,0,0,0,0,0},
3372302aa1bSDebojyoti Ghosh                    {0.585786437626904966,0,0,0,0,0,0,0,0},
3382302aa1bSDebojyoti Ghosh                    {0.149999999999999994,0.849999999999999978,0,0,0,0,0,0,0},
3392302aa1bSDebojyoti Ghosh                    {0,0,0,0,0,0,0,0,0},
3402302aa1bSDebojyoti Ghosh                    {0,0,0,0.292893218813452483,0,0,0,0,0},
3412302aa1bSDebojyoti Ghosh                    {0,0,0,0.0749999999999999972,0.424999999999999989,0,0,0,0},
3422302aa1bSDebojyoti Ghosh                    {0,0,0,0.176776695296636893,0.176776695296636893,0.146446609406726241,0,0,0},
3432302aa1bSDebojyoti Ghosh                    {0,0,0,0.176776695296636893,0.176776695296636893,0.146446609406726241,0.292893218813452483,0,0},
3442302aa1bSDebojyoti Ghosh                    {0,0,0,0.176776695296636893,0.176776695296636893,0.146446609406726241,0.0749999999999999972,0.424999999999999989,0}},
3452302aa1bSDebojyoti Ghosh       B[2][9]   = {{0.353553390593273786,0.353553390593273786,0.292893218813452483,0,0,0,0,0,0},
3462302aa1bSDebojyoti Ghosh                    {-0.471404520791031678,-0.471404520791031678,-0.390524291751269959,0.235702260395515839,0.235702260395515839,0.195262145875634979,0.235702260395515839,0.235702260395515839,0.195262145875634979}},
3472302aa1bSDebojyoti Ghosh       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}},
3482302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
3492302aa1bSDebojyoti Ghosh       S[2]      = {1,0},
3502302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
351b1fbfd33SSatish Balay       Fembed[2] = {1,1-GAMMA},
35257df6a1bSDebojyoti Ghosh       Ferror[2] = {0,1},
35357df6a1bSDebojyoti Ghosh       Serror[2] = {1,0};
354b1fbfd33SSatish Balay     ierr = 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);CHKERRQ(ierr);
3552302aa1bSDebojyoti Ghosh   }
3562302aa1bSDebojyoti Ghosh 
3572302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
3582302aa1bSDebojyoti Ghosh }
3592302aa1bSDebojyoti Ghosh 
3602302aa1bSDebojyoti Ghosh /*@C
3612302aa1bSDebojyoti Ghosh    TSGLEERegisterDestroy - Frees the list of schemes that were registered by TSGLEERegister().
3622302aa1bSDebojyoti Ghosh 
3632302aa1bSDebojyoti Ghosh    Not Collective
3642302aa1bSDebojyoti Ghosh 
3652302aa1bSDebojyoti Ghosh    Level: advanced
3662302aa1bSDebojyoti Ghosh 
3672302aa1bSDebojyoti Ghosh .keywords: TSGLEE, register, destroy
3682302aa1bSDebojyoti Ghosh .seealso: TSGLEERegister(), TSGLEERegisterAll()
3692302aa1bSDebojyoti Ghosh @*/
3702302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEERegisterDestroy(void)
3712302aa1bSDebojyoti Ghosh {
3722302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
3732302aa1bSDebojyoti Ghosh   GLEETableauLink link;
3742302aa1bSDebojyoti Ghosh 
3752302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
3762302aa1bSDebojyoti Ghosh   while ((link = GLEETableauList)) {
3772302aa1bSDebojyoti Ghosh     GLEETableau t = &link->tab;
3782302aa1bSDebojyoti Ghosh     GLEETableauList = link->next;
3792302aa1bSDebojyoti Ghosh     ierr = PetscFree5(t->A,t->B,t->U,t->V,t->c);  CHKERRQ(ierr);
3802302aa1bSDebojyoti Ghosh     ierr = PetscFree2(t->S,t->F);                 CHKERRQ(ierr);
3812302aa1bSDebojyoti Ghosh     ierr = PetscFree (t->Fembed);                 CHKERRQ(ierr);
382fd96ebc2SDebojyoti Ghosh     ierr = PetscFree (t->Ferror);                 CHKERRQ(ierr);
38357df6a1bSDebojyoti Ghosh     ierr = PetscFree (t->Serror);                 CHKERRQ(ierr);
3842302aa1bSDebojyoti Ghosh     ierr = PetscFree (t->binterp);                CHKERRQ(ierr);
3852302aa1bSDebojyoti Ghosh     ierr = PetscFree (t->name);                   CHKERRQ(ierr);
3862302aa1bSDebojyoti Ghosh     ierr = PetscFree (link);                      CHKERRQ(ierr);
3872302aa1bSDebojyoti Ghosh   }
3882302aa1bSDebojyoti Ghosh   TSGLEERegisterAllCalled = PETSC_FALSE;
3892302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
3902302aa1bSDebojyoti Ghosh }
3912302aa1bSDebojyoti Ghosh 
3922302aa1bSDebojyoti Ghosh /*@C
3932302aa1bSDebojyoti Ghosh   TSGLEEInitializePackage - This function initializes everything in the TSGLEE package. It is called
3942302aa1bSDebojyoti Ghosh   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_GLEE()
3952302aa1bSDebojyoti Ghosh   when using static libraries.
3962302aa1bSDebojyoti Ghosh 
3972302aa1bSDebojyoti Ghosh   Level: developer
3982302aa1bSDebojyoti Ghosh 
3992302aa1bSDebojyoti Ghosh .keywords: TS, TSGLEE, initialize, package
4002302aa1bSDebojyoti Ghosh .seealso: PetscInitialize()
4012302aa1bSDebojyoti Ghosh @*/
4022302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEEInitializePackage(void)
4032302aa1bSDebojyoti Ghosh {
4042302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
4052302aa1bSDebojyoti Ghosh 
4062302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
4072302aa1bSDebojyoti Ghosh   if (TSGLEEPackageInitialized) PetscFunctionReturn(0);
4082302aa1bSDebojyoti Ghosh   TSGLEEPackageInitialized = PETSC_TRUE;
4092302aa1bSDebojyoti Ghosh   ierr = TSGLEERegisterAll();CHKERRQ(ierr);
4102302aa1bSDebojyoti Ghosh   ierr = PetscObjectComposedDataRegister(&explicit_stage_time_id);CHKERRQ(ierr);
4112302aa1bSDebojyoti Ghosh   ierr = PetscRegisterFinalize(TSGLEEFinalizePackage);CHKERRQ(ierr);
4122302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
4132302aa1bSDebojyoti Ghosh }
4142302aa1bSDebojyoti Ghosh 
4152302aa1bSDebojyoti Ghosh /*@C
4162302aa1bSDebojyoti Ghosh   TSGLEEFinalizePackage - This function destroys everything in the TSGLEE package. It is
4172302aa1bSDebojyoti Ghosh   called from PetscFinalize().
4182302aa1bSDebojyoti Ghosh 
4192302aa1bSDebojyoti Ghosh   Level: developer
4202302aa1bSDebojyoti Ghosh 
4212302aa1bSDebojyoti Ghosh .keywords: Petsc, destroy, package
4222302aa1bSDebojyoti Ghosh .seealso: PetscFinalize()
4232302aa1bSDebojyoti Ghosh @*/
4242302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEEFinalizePackage(void)
4252302aa1bSDebojyoti Ghosh {
4262302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
4272302aa1bSDebojyoti Ghosh 
4282302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
4292302aa1bSDebojyoti Ghosh   TSGLEEPackageInitialized = PETSC_FALSE;
4302302aa1bSDebojyoti Ghosh   ierr = TSGLEERegisterDestroy();CHKERRQ(ierr);
4312302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
4322302aa1bSDebojyoti Ghosh }
4332302aa1bSDebojyoti Ghosh 
4342302aa1bSDebojyoti Ghosh /*@C
4352302aa1bSDebojyoti Ghosh    TSGLEERegister - register an GLEE scheme by providing the entries in the Butcher tableau
4362302aa1bSDebojyoti Ghosh 
4372302aa1bSDebojyoti Ghosh    Not Collective, but the same schemes should be registered on all processes on which they will be used
4382302aa1bSDebojyoti Ghosh 
4392302aa1bSDebojyoti Ghosh    Input Parameters:
4402302aa1bSDebojyoti Ghosh +  name   - identifier for method
4412302aa1bSDebojyoti Ghosh .  order  - order of method
4422302aa1bSDebojyoti Ghosh .  s - number of stages
4432302aa1bSDebojyoti Ghosh .  r - number of steps
4442302aa1bSDebojyoti Ghosh .  gamma - LTE ratio
4452302aa1bSDebojyoti Ghosh .  A - stage coefficients (dimension s*s, row-major)
4462302aa1bSDebojyoti Ghosh .  B - step completion coefficients (dimension r*s, row-major)
4472302aa1bSDebojyoti Ghosh .  U - method coefficients (dimension s*r, row-major)
4482302aa1bSDebojyoti Ghosh .  V - method coefficients (dimension r*r, row-major)
4492302aa1bSDebojyoti Ghosh .  S - starting coefficients
4502302aa1bSDebojyoti Ghosh .  F - finishing coefficients
4512302aa1bSDebojyoti Ghosh .  c - abscissa (dimension s; NULL to use row sums of A)
4522302aa1bSDebojyoti Ghosh .  Fembed - step completion coefficients for embedded method
453fd96ebc2SDebojyoti Ghosh .  Ferror - error computation coefficients
45457df6a1bSDebojyoti Ghosh .  Serror - error initialization coefficients
4552302aa1bSDebojyoti Ghosh .  pinterp - order of interpolation (0 if unavailable)
4562302aa1bSDebojyoti Ghosh -  binterp - array of interpolation coefficients (NULL if unavailable)
4572302aa1bSDebojyoti Ghosh 
4582302aa1bSDebojyoti Ghosh    Notes:
4592302aa1bSDebojyoti Ghosh    Several GLEE methods are provided, this function is only needed to create new methods.
4602302aa1bSDebojyoti Ghosh 
4612302aa1bSDebojyoti Ghosh    Level: advanced
4622302aa1bSDebojyoti Ghosh 
4632302aa1bSDebojyoti Ghosh .keywords: TS, register
4642302aa1bSDebojyoti Ghosh 
4652302aa1bSDebojyoti Ghosh .seealso: TSGLEE
4662302aa1bSDebojyoti Ghosh @*/
4672302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEERegister(TSGLEEType name,PetscInt order,PetscInt s, PetscInt r,
4682302aa1bSDebojyoti Ghosh                               PetscReal gamma,
4692302aa1bSDebojyoti Ghosh                               const PetscReal A[],const PetscReal B[],
4702302aa1bSDebojyoti Ghosh                               const PetscReal U[],const PetscReal V[],
4712302aa1bSDebojyoti Ghosh                               const PetscReal S[],const PetscReal F[],
472fd96ebc2SDebojyoti Ghosh                               const PetscReal c[],
473fd96ebc2SDebojyoti Ghosh                               const PetscReal Fembed[],const PetscReal Ferror[],
47457df6a1bSDebojyoti Ghosh                               const PetscReal Serror[],
4752302aa1bSDebojyoti Ghosh                               PetscInt pinterp, const PetscReal binterp[])
4762302aa1bSDebojyoti Ghosh {
4772302aa1bSDebojyoti Ghosh   PetscErrorCode    ierr;
4782302aa1bSDebojyoti Ghosh   GLEETableauLink   link;
4792302aa1bSDebojyoti Ghosh   GLEETableau       t;
4802302aa1bSDebojyoti Ghosh   PetscInt          i,j;
4812302aa1bSDebojyoti Ghosh 
4822302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
4831d36bdfdSBarry Smith   ierr     = TSGLEEInitializePackage();CHKERRQ(ierr);
4842302aa1bSDebojyoti Ghosh   ierr     = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
4852302aa1bSDebojyoti Ghosh   ierr     = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
4862302aa1bSDebojyoti Ghosh   t        = &link->tab;
4872302aa1bSDebojyoti Ghosh   ierr     = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
4882302aa1bSDebojyoti Ghosh   t->order = order;
4892302aa1bSDebojyoti Ghosh   t->s     = s;
4902302aa1bSDebojyoti Ghosh   t->r     = r;
4912302aa1bSDebojyoti Ghosh   t->gamma = gamma;
4922302aa1bSDebojyoti Ghosh   ierr     = PetscMalloc5(s*s,&t->A,r*r,&t->V,s,&t->c,r*s,&t->B,s*r,&t->U); CHKERRQ(ierr);
4932302aa1bSDebojyoti Ghosh   ierr     = PetscMalloc2(r,&t->S,r,&t->F); CHKERRQ(ierr);
4942302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
4952302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->B,B,r*s*sizeof(B[0]));CHKERRQ(ierr);
4962302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->U,U,s*r*sizeof(U[0]));CHKERRQ(ierr);
4972302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->V,V,r*r*sizeof(V[0]));CHKERRQ(ierr);
4982302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->S,S,r  *sizeof(S[0]));CHKERRQ(ierr);
4992302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->F,F,r  *sizeof(F[0]));CHKERRQ(ierr);
5002302aa1bSDebojyoti Ghosh   if (c) {
5012302aa1bSDebojyoti Ghosh     ierr   = PetscMemcpy(t->c,c,s*sizeof(c[0])); CHKERRQ(ierr);
5022302aa1bSDebojyoti Ghosh   } else {
5032302aa1bSDebojyoti Ghosh     for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j];
5042302aa1bSDebojyoti Ghosh   }
5052302aa1bSDebojyoti Ghosh   ierr = PetscMalloc1(r,&t->Fembed);CHKERRQ(ierr);
506fd96ebc2SDebojyoti Ghosh   ierr = PetscMalloc1(r,&t->Ferror);CHKERRQ(ierr);
50757df6a1bSDebojyoti Ghosh   ierr = PetscMalloc1(r,&t->Serror);CHKERRQ(ierr);
5082302aa1bSDebojyoti Ghosh   ierr = PetscMemcpy(t->Fembed,Fembed,r*sizeof(Fembed[0]));CHKERRQ(ierr);
509fd96ebc2SDebojyoti Ghosh   ierr = PetscMemcpy(t->Ferror,Ferror,r*sizeof(Ferror[0]));CHKERRQ(ierr);
51057df6a1bSDebojyoti Ghosh   ierr = PetscMemcpy(t->Serror,Serror,r*sizeof(Serror[0]));CHKERRQ(ierr);
5112302aa1bSDebojyoti Ghosh   t->pinterp = pinterp;
5122302aa1bSDebojyoti Ghosh   ierr       = PetscMalloc1(s*pinterp,&t->binterp);CHKERRQ(ierr);
5132302aa1bSDebojyoti Ghosh   ierr       = PetscMemcpy(t->binterp,binterp,s*pinterp*sizeof(binterp[0]));CHKERRQ(ierr);
5142302aa1bSDebojyoti Ghosh 
5152302aa1bSDebojyoti Ghosh   link->next      = GLEETableauList;
5162302aa1bSDebojyoti Ghosh   GLEETableauList = link;
5172302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
5182302aa1bSDebojyoti Ghosh }
5192302aa1bSDebojyoti Ghosh 
5202302aa1bSDebojyoti Ghosh static PetscErrorCode TSEvaluateStep_GLEE(TS ts,PetscInt order,Vec X,PetscBool *done)
5212302aa1bSDebojyoti Ghosh {
5222302aa1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*) ts->data;
5232302aa1bSDebojyoti Ghosh   GLEETableau     tab = glee->tableau;
5242302aa1bSDebojyoti Ghosh   PetscReal       h, *B = tab->B, *V = tab->V,
525fd96ebc2SDebojyoti Ghosh                   *F = tab->F,
526fd96ebc2SDebojyoti Ghosh                   *Fembed = tab->Fembed;
5272302aa1bSDebojyoti Ghosh   PetscInt        s = tab->s, r = tab->r, i, j;
5282302aa1bSDebojyoti Ghosh   Vec             *Y = glee->Y, *YdotStage = glee->YdotStage;
529ec0e3ee2SEmil Constantinescu   PetscScalar     *ws = glee->swork, *wr = glee->rwork;
5302302aa1bSDebojyoti Ghosh   PetscErrorCode  ierr;
5312302aa1bSDebojyoti Ghosh 
5322302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
5332302aa1bSDebojyoti Ghosh 
5342302aa1bSDebojyoti Ghosh   switch (glee->status) {
5352302aa1bSDebojyoti Ghosh     case TS_STEP_INCOMPLETE:
5362302aa1bSDebojyoti Ghosh     case TS_STEP_PENDING:
5372302aa1bSDebojyoti Ghosh       h = ts->time_step; break;
5382302aa1bSDebojyoti Ghosh     case TS_STEP_COMPLETE:
539ed450d42SEmil Constantinescu       h = ts->ptime - ts->ptime_prev; break;
5402302aa1bSDebojyoti Ghosh     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
5412302aa1bSDebojyoti Ghosh   }
5422302aa1bSDebojyoti Ghosh 
5432302aa1bSDebojyoti Ghosh   if (order == tab->order) {
5442302aa1bSDebojyoti Ghosh 
5452302aa1bSDebojyoti Ghosh     /* Note: Irrespective of whether status is TS_STEP_INCOMPLETE
546fd96ebc2SDebojyoti Ghosh              or TS_STEP_COMPLETE, glee->X has the solution at the
5472302aa1bSDebojyoti Ghosh              beginning of the time step. So no need to roll-back.
5482302aa1bSDebojyoti Ghosh     */
5492302aa1bSDebojyoti Ghosh     if (glee->status == TS_STEP_INCOMPLETE) {
5502302aa1bSDebojyoti Ghosh       for (i=0; i<r; i++) {
5512302aa1bSDebojyoti Ghosh         ierr = VecZeroEntries(Y[i]); CHKERRQ(ierr);
552ec0e3ee2SEmil Constantinescu         for (j=0; j<r; j++) wr[j] = V[i*r+j];
553ec0e3ee2SEmil Constantinescu         ierr = VecMAXPY(Y[i],r,wr,glee->X); CHKERRQ(ierr);
554ec0e3ee2SEmil Constantinescu         for (j=0; j<s; j++) ws[j] = h*B[i*s+j];
555ec0e3ee2SEmil Constantinescu         ierr = VecMAXPY(Y[i],s,ws,YdotStage); CHKERRQ(ierr);
5562302aa1bSDebojyoti Ghosh       }
5572302aa1bSDebojyoti Ghosh       ierr = VecZeroEntries(X);CHKERRQ(ierr);
558ec0e3ee2SEmil Constantinescu       for (j=0; j<r; j++) wr[j] = F[j];
559ec0e3ee2SEmil Constantinescu       ierr = VecMAXPY(X,r,wr,Y);CHKERRQ(ierr);
5602302aa1bSDebojyoti Ghosh     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
5612302aa1bSDebojyoti Ghosh     PetscFunctionReturn(0);
5622302aa1bSDebojyoti Ghosh 
5632302aa1bSDebojyoti Ghosh   } else if (order == tab->order-1) {
5642302aa1bSDebojyoti Ghosh 
5652302aa1bSDebojyoti Ghosh     /* Complete with the embedded method (Fembed) */
5662302aa1bSDebojyoti Ghosh     for (i=0; i<r; i++) {
5672302aa1bSDebojyoti Ghosh       ierr = VecZeroEntries(Y[i]); CHKERRQ(ierr);
568ec0e3ee2SEmil Constantinescu       for (j=0; j<r; j++) wr[j] = V[i*r+j];
569ec0e3ee2SEmil Constantinescu       ierr = VecMAXPY(Y[i],r,wr,glee->X); CHKERRQ(ierr);
570ec0e3ee2SEmil Constantinescu       for (j=0; j<s; j++) ws[j] = h*B[i*s+j];
571ec0e3ee2SEmil Constantinescu       ierr = VecMAXPY(Y[i],s,ws,YdotStage); CHKERRQ(ierr);
5722302aa1bSDebojyoti Ghosh     }
5732302aa1bSDebojyoti Ghosh     ierr = VecZeroEntries(X);CHKERRQ(ierr);
574ec0e3ee2SEmil Constantinescu     for (j=0; j<r; j++) wr[j] = Fembed[j];
575ec0e3ee2SEmil Constantinescu     ierr = VecMAXPY(X,r,wr,Y);CHKERRQ(ierr);
576fd96ebc2SDebojyoti Ghosh 
5772302aa1bSDebojyoti Ghosh     if (done) *done = PETSC_TRUE;
5782302aa1bSDebojyoti Ghosh     PetscFunctionReturn(0);
5792302aa1bSDebojyoti Ghosh   }
5802302aa1bSDebojyoti Ghosh   if (done) *done = PETSC_FALSE;
5812302aa1bSDebojyoti Ghosh   else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"GLEE '%s' of order %D cannot evaluate step at order %D",tab->name,tab->order,order);
5822302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
5832302aa1bSDebojyoti Ghosh }
5842302aa1bSDebojyoti Ghosh 
5852302aa1bSDebojyoti Ghosh static PetscErrorCode TSStep_GLEE(TS ts)
5862302aa1bSDebojyoti Ghosh {
5872302aa1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
5882302aa1bSDebojyoti Ghosh   GLEETableau     tab = glee->tableau;
5892302aa1bSDebojyoti Ghosh   const PetscInt  s = tab->s, r = tab->r;
5902302aa1bSDebojyoti Ghosh   PetscReal       *A = tab->A, *U = tab->U,
5912302aa1bSDebojyoti Ghosh                   *F = tab->F,
5922302aa1bSDebojyoti Ghosh                   *c = tab->c;
5932302aa1bSDebojyoti Ghosh   Vec             *Y = glee->Y, *X = glee->X,
5942302aa1bSDebojyoti Ghosh                   *YStage = glee->YStage,
5952302aa1bSDebojyoti Ghosh                   *YdotStage = glee->YdotStage,
5962302aa1bSDebojyoti Ghosh                   W = glee->W;
5972302aa1bSDebojyoti Ghosh   SNES            snes;
598ec0e3ee2SEmil Constantinescu   PetscScalar     *ws = glee->swork, *wr = glee->rwork;
5992302aa1bSDebojyoti Ghosh   TSAdapt         adapt;
6002302aa1bSDebojyoti Ghosh   PetscInt        i,j,reject,next_scheme,its,lits;
6012302aa1bSDebojyoti Ghosh   PetscReal       next_time_step;
6022302aa1bSDebojyoti Ghosh   PetscReal       t;
6032302aa1bSDebojyoti Ghosh   PetscBool       accept;
6042302aa1bSDebojyoti Ghosh   PetscErrorCode  ierr;
6052302aa1bSDebojyoti Ghosh 
6062302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
607642f3702SEmil Constantinescu   ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr);
608642f3702SEmil Constantinescu 
6092302aa1bSDebojyoti Ghosh   for (i=0; i<r; i++) { ierr = VecCopy(Y[i],X[i]); CHKERRQ(ierr); }
6102302aa1bSDebojyoti Ghosh 
6112302aa1bSDebojyoti Ghosh   ierr           = TSGetSNES(ts,&snes);CHKERRQ(ierr);
6122302aa1bSDebojyoti Ghosh   next_time_step = ts->time_step;
6132302aa1bSDebojyoti Ghosh   t              = ts->ptime;
6142302aa1bSDebojyoti Ghosh   accept         = PETSC_TRUE;
6152302aa1bSDebojyoti Ghosh   glee->status   = TS_STEP_INCOMPLETE;
6162302aa1bSDebojyoti Ghosh 
6172302aa1bSDebojyoti Ghosh   for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) {
6182302aa1bSDebojyoti Ghosh 
6192302aa1bSDebojyoti Ghosh     PetscReal h = ts->time_step;
6202302aa1bSDebojyoti Ghosh     ierr = TSPreStep(ts);CHKERRQ(ierr);
6212302aa1bSDebojyoti Ghosh 
6222302aa1bSDebojyoti Ghosh     for (i=0; i<s; i++) {
6232302aa1bSDebojyoti Ghosh 
6242302aa1bSDebojyoti Ghosh       glee->stage_time = t + h*c[i];
6252302aa1bSDebojyoti Ghosh       ierr = TSPreStage(ts,glee->stage_time); CHKERRQ(ierr);
6262302aa1bSDebojyoti Ghosh 
6272302aa1bSDebojyoti Ghosh       if (A[i*s+i] == 0) {  /* Explicit stage */
6282302aa1bSDebojyoti Ghosh         ierr = VecZeroEntries(YStage[i]);CHKERRQ(ierr);
629ec0e3ee2SEmil Constantinescu         for (j=0; j<r; j++) wr[j] = U[i*r+j];
630ec0e3ee2SEmil Constantinescu         ierr = VecMAXPY(YStage[i],r,wr,X);CHKERRQ(ierr);
631ec0e3ee2SEmil Constantinescu         for (j=0; j<i; j++) ws[j] = h*A[i*s+j];
632ec0e3ee2SEmil Constantinescu         ierr = VecMAXPY(YStage[i],i,ws,YdotStage);CHKERRQ(ierr);
6332302aa1bSDebojyoti Ghosh       } else {              /* Implicit stage */
6342302aa1bSDebojyoti Ghosh         glee->scoeff = 1.0/A[i*s+i];
6352302aa1bSDebojyoti Ghosh         /* compute right-hand-side */
6362302aa1bSDebojyoti Ghosh         ierr = VecZeroEntries(W);CHKERRQ(ierr);
637ec0e3ee2SEmil Constantinescu         for (j=0; j<r; j++) wr[j] = U[i*r+j];
638ec0e3ee2SEmil Constantinescu         ierr = VecMAXPY(W,r,wr,X);CHKERRQ(ierr);
639ec0e3ee2SEmil Constantinescu         for (j=0; j<i; j++) ws[j] = h*A[i*s+j];
640ec0e3ee2SEmil Constantinescu         ierr = VecMAXPY(W,i,ws,YdotStage);CHKERRQ(ierr);
6412302aa1bSDebojyoti Ghosh         ierr = VecScale(W,glee->scoeff/h);CHKERRQ(ierr);
6422302aa1bSDebojyoti Ghosh         /* set initial guess */
6432302aa1bSDebojyoti Ghosh         ierr = VecCopy(i>0 ? YStage[i-1] : ts->vec_sol,YStage[i]);CHKERRQ(ierr);
6442302aa1bSDebojyoti Ghosh         /* solve for this stage */
6452302aa1bSDebojyoti Ghosh         ierr = SNESSolve(snes,W,YStage[i]);CHKERRQ(ierr);
6462302aa1bSDebojyoti Ghosh         ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
6472302aa1bSDebojyoti Ghosh         ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
6482302aa1bSDebojyoti Ghosh         ts->snes_its += its; ts->ksp_its += lits;
6492302aa1bSDebojyoti Ghosh       }
6502302aa1bSDebojyoti Ghosh       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
651d6629dc5SEmil Constantinescu       ierr = TSAdaptCheckStage(adapt,ts,glee->stage_time,YStage[i],&accept);CHKERRQ(ierr);
6522302aa1bSDebojyoti Ghosh       if (!accept) goto reject_step;
6532302aa1bSDebojyoti Ghosh       ierr = TSPostStage(ts,glee->stage_time,i,YStage); CHKERRQ(ierr);
6542302aa1bSDebojyoti Ghosh       ierr = TSComputeRHSFunction(ts,t+h*c[i],YStage[i],YdotStage[i]);CHKERRQ(ierr);
6552302aa1bSDebojyoti Ghosh     }
6562302aa1bSDebojyoti Ghosh     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
6572302aa1bSDebojyoti Ghosh     glee->status = TS_STEP_PENDING;
6582302aa1bSDebojyoti Ghosh 
6592302aa1bSDebojyoti Ghosh     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
6602302aa1bSDebojyoti Ghosh     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
6612302aa1bSDebojyoti Ghosh     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
6621917a363SLisandro Dalcin     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr);
6632302aa1bSDebojyoti Ghosh     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
6642302aa1bSDebojyoti Ghosh     if (accept) {
6652302aa1bSDebojyoti Ghosh       /* ignore next_scheme for now */
6662302aa1bSDebojyoti Ghosh       ts->ptime     += ts->time_step;
6672302aa1bSDebojyoti Ghosh       ts->time_step  = next_time_step;
6682302aa1bSDebojyoti Ghosh       glee->status = TS_STEP_COMPLETE;
6690a01e1b2SEmil Constantinescu       /* compute and store the global error */
6700a01e1b2SEmil Constantinescu       /* Note: this is not needed if TSAdaptGLEE is not used */
6710a01e1b2SEmil Constantinescu       ierr = TSGetTimeError(ts,0,&(glee->yGErr));CHKERRQ(ierr);
6722302aa1bSDebojyoti Ghosh       ierr = PetscObjectComposedDataSetReal((PetscObject)ts->vec_sol,explicit_stage_time_id,ts->ptime);CHKERRQ(ierr);
6732302aa1bSDebojyoti Ghosh       break;
6742302aa1bSDebojyoti Ghosh     } else {                    /* Roll back the current step */
675ec0e3ee2SEmil Constantinescu       for (j=0; j<r; j++) wr[j] = F[j];
676ec0e3ee2SEmil Constantinescu       ierr = VecMAXPY(ts->vec_sol,r,wr,X);CHKERRQ(ierr);
6772302aa1bSDebojyoti Ghosh       ts->time_step = next_time_step;
6782302aa1bSDebojyoti Ghosh       glee->status  = TS_STEP_INCOMPLETE;
6792302aa1bSDebojyoti Ghosh     }
6802302aa1bSDebojyoti Ghosh reject_step: continue;
6812302aa1bSDebojyoti Ghosh   }
6822302aa1bSDebojyoti Ghosh   if (glee->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
6832302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
6842302aa1bSDebojyoti Ghosh }
6852302aa1bSDebojyoti Ghosh 
6862302aa1bSDebojyoti Ghosh static PetscErrorCode TSInterpolate_GLEE(TS ts,PetscReal itime,Vec X)
6872302aa1bSDebojyoti Ghosh {
6882302aa1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
6892302aa1bSDebojyoti Ghosh   PetscInt        s=glee->tableau->s, pinterp=glee->tableau->pinterp,i,j;
6902302aa1bSDebojyoti Ghosh   PetscReal       h,tt,t;
6912302aa1bSDebojyoti Ghosh   PetscScalar     *b;
6922302aa1bSDebojyoti Ghosh   const PetscReal *B = glee->tableau->binterp;
6932302aa1bSDebojyoti Ghosh   PetscErrorCode  ierr;
6942302aa1bSDebojyoti Ghosh 
6952302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
6962302aa1bSDebojyoti Ghosh   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSGLEE %s does not have an interpolation formula",glee->tableau->name);
6972302aa1bSDebojyoti Ghosh   switch (glee->status) {
6982302aa1bSDebojyoti Ghosh     case TS_STEP_INCOMPLETE:
6992302aa1bSDebojyoti Ghosh     case TS_STEP_PENDING:
7002302aa1bSDebojyoti Ghosh       h = ts->time_step;
7012302aa1bSDebojyoti Ghosh       t = (itime - ts->ptime)/h;
7022302aa1bSDebojyoti Ghosh       break;
7032302aa1bSDebojyoti Ghosh     case TS_STEP_COMPLETE:
704ed450d42SEmil Constantinescu       h = ts->ptime - ts->ptime_prev;
7052302aa1bSDebojyoti Ghosh       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
7062302aa1bSDebojyoti Ghosh       break;
7072302aa1bSDebojyoti Ghosh     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
7082302aa1bSDebojyoti Ghosh   }
7092302aa1bSDebojyoti Ghosh   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
7102302aa1bSDebojyoti Ghosh   for (i=0; i<s; i++) b[i] = 0;
7112302aa1bSDebojyoti Ghosh   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
7122302aa1bSDebojyoti Ghosh     for (i=0; i<s; i++) {
7132302aa1bSDebojyoti Ghosh       b[i]  += h * B[i*pinterp+j] * tt;
7142302aa1bSDebojyoti Ghosh     }
7152302aa1bSDebojyoti Ghosh   }
7162302aa1bSDebojyoti Ghosh   ierr = VecCopy(glee->YStage[0],X);CHKERRQ(ierr);
7172302aa1bSDebojyoti Ghosh   ierr = VecMAXPY(X,s,b,glee->YdotStage);CHKERRQ(ierr);
7182302aa1bSDebojyoti Ghosh   ierr = PetscFree(b);CHKERRQ(ierr);
7192302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
7202302aa1bSDebojyoti Ghosh }
7212302aa1bSDebojyoti Ghosh 
7222302aa1bSDebojyoti Ghosh /*------------------------------------------------------------*/
7232302aa1bSDebojyoti Ghosh static PetscErrorCode TSReset_GLEE(TS ts)
7242302aa1bSDebojyoti Ghosh {
7252302aa1bSDebojyoti Ghosh   TS_GLEE        *glee = (TS_GLEE*)ts->data;
7262302aa1bSDebojyoti Ghosh   PetscInt       s, r;
7272302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
7282302aa1bSDebojyoti Ghosh 
7292302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7302302aa1bSDebojyoti Ghosh   if (!glee->tableau) PetscFunctionReturn(0);
7312302aa1bSDebojyoti Ghosh   s    = glee->tableau->s;
7322302aa1bSDebojyoti Ghosh   r    = glee->tableau->r;
7332302aa1bSDebojyoti Ghosh   ierr = VecDestroyVecs(r,&glee->Y);CHKERRQ(ierr);
7342302aa1bSDebojyoti Ghosh   ierr = VecDestroyVecs(r,&glee->X);CHKERRQ(ierr);
7352302aa1bSDebojyoti Ghosh   ierr = VecDestroyVecs(s,&glee->YStage);CHKERRQ(ierr);
7362302aa1bSDebojyoti Ghosh   ierr = VecDestroyVecs(s,&glee->YdotStage);CHKERRQ(ierr);
7372302aa1bSDebojyoti Ghosh   ierr = VecDestroy(&glee->Ydot);CHKERRQ(ierr);
7380a01e1b2SEmil Constantinescu   ierr = VecDestroy(&glee->yGErr);CHKERRQ(ierr);
7392302aa1bSDebojyoti Ghosh   ierr = VecDestroy(&glee->W);CHKERRQ(ierr);
740ec0e3ee2SEmil Constantinescu   ierr = PetscFree2(glee->swork,glee->rwork);CHKERRQ(ierr);
7412302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
7422302aa1bSDebojyoti Ghosh }
7432302aa1bSDebojyoti Ghosh 
7442302aa1bSDebojyoti Ghosh static PetscErrorCode TSGLEEGetVecs(TS ts,DM dm,Vec *Ydot)
7452302aa1bSDebojyoti Ghosh {
7462302aa1bSDebojyoti Ghosh   TS_GLEE     *glee = (TS_GLEE*)ts->data;
7472302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
7482302aa1bSDebojyoti Ghosh 
7492302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7502302aa1bSDebojyoti Ghosh   if (Ydot) {
7512302aa1bSDebojyoti Ghosh     if (dm && dm != ts->dm) {
7522302aa1bSDebojyoti Ghosh       ierr = DMGetNamedGlobalVector(dm,"TSGLEE_Ydot",Ydot);CHKERRQ(ierr);
7532302aa1bSDebojyoti Ghosh     } else *Ydot = glee->Ydot;
7542302aa1bSDebojyoti Ghosh   }
7552302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
7562302aa1bSDebojyoti Ghosh }
7572302aa1bSDebojyoti Ghosh 
7582302aa1bSDebojyoti Ghosh 
7592302aa1bSDebojyoti Ghosh static PetscErrorCode TSGLEERestoreVecs(TS ts,DM dm,Vec *Ydot)
7602302aa1bSDebojyoti Ghosh {
7612302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
7622302aa1bSDebojyoti Ghosh 
7632302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7642302aa1bSDebojyoti Ghosh   if (Ydot) {
7652302aa1bSDebojyoti Ghosh     if (dm && dm != ts->dm) {
7662302aa1bSDebojyoti Ghosh       ierr = DMRestoreNamedGlobalVector(dm,"TSGLEE_Ydot",Ydot);CHKERRQ(ierr);
7672302aa1bSDebojyoti Ghosh     }
7682302aa1bSDebojyoti Ghosh   }
7692302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
7702302aa1bSDebojyoti Ghosh }
7712302aa1bSDebojyoti Ghosh 
7722302aa1bSDebojyoti Ghosh /*
7732302aa1bSDebojyoti Ghosh   This defines the nonlinear equation that is to be solved with SNES
7742302aa1bSDebojyoti Ghosh */
7752302aa1bSDebojyoti Ghosh static PetscErrorCode SNESTSFormFunction_GLEE(SNES snes,Vec X,Vec F,TS ts)
7762302aa1bSDebojyoti Ghosh {
7772302aa1bSDebojyoti Ghosh   TS_GLEE       *glee = (TS_GLEE*)ts->data;
7782302aa1bSDebojyoti Ghosh   DM             dm,dmsave;
7792302aa1bSDebojyoti Ghosh   Vec            Ydot;
7802302aa1bSDebojyoti Ghosh   PetscReal      shift = glee->scoeff / ts->time_step;
7812302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
7822302aa1bSDebojyoti Ghosh 
7832302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7842302aa1bSDebojyoti Ghosh   ierr   = SNESGetDM(snes,&dm);CHKERRQ(ierr);
7852302aa1bSDebojyoti Ghosh   ierr   = TSGLEEGetVecs(ts,dm,&Ydot);CHKERRQ(ierr);
7862302aa1bSDebojyoti Ghosh   /* Set Ydot = shift*X */
7872302aa1bSDebojyoti Ghosh   ierr   = VecCopy(X,Ydot);CHKERRQ(ierr);
7882302aa1bSDebojyoti Ghosh   ierr   = VecScale(Ydot,shift);CHKERRQ(ierr);
7892302aa1bSDebojyoti Ghosh   dmsave = ts->dm;
7902302aa1bSDebojyoti Ghosh   ts->dm = dm;
7912302aa1bSDebojyoti Ghosh 
7922302aa1bSDebojyoti Ghosh   ierr = TSComputeIFunction(ts,glee->stage_time,X,Ydot,F,PETSC_FALSE);CHKERRQ(ierr);
7932302aa1bSDebojyoti Ghosh 
7942302aa1bSDebojyoti Ghosh   ts->dm = dmsave;
7952302aa1bSDebojyoti Ghosh   ierr   = TSGLEERestoreVecs(ts,dm,&Ydot);CHKERRQ(ierr);
7962302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
7972302aa1bSDebojyoti Ghosh }
7982302aa1bSDebojyoti Ghosh 
7992302aa1bSDebojyoti Ghosh static PetscErrorCode SNESTSFormJacobian_GLEE(SNES snes,Vec X,Mat A,Mat B,TS ts)
8002302aa1bSDebojyoti Ghosh {
8012302aa1bSDebojyoti Ghosh   TS_GLEE        *glee = (TS_GLEE*)ts->data;
8022302aa1bSDebojyoti Ghosh   DM             dm,dmsave;
8032302aa1bSDebojyoti Ghosh   Vec            Ydot;
8042302aa1bSDebojyoti Ghosh   PetscReal      shift = glee->scoeff / ts->time_step;
8052302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
8062302aa1bSDebojyoti Ghosh 
8072302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8082302aa1bSDebojyoti Ghosh   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
8092302aa1bSDebojyoti Ghosh   ierr = TSGLEEGetVecs(ts,dm,&Ydot);CHKERRQ(ierr);
8102302aa1bSDebojyoti Ghosh   /* glee->Ydot has already been computed in SNESTSFormFunction_GLEE (SNES guarantees this) */
8112302aa1bSDebojyoti Ghosh   dmsave = ts->dm;
8122302aa1bSDebojyoti Ghosh   ts->dm = dm;
8132302aa1bSDebojyoti Ghosh 
8142302aa1bSDebojyoti Ghosh   ierr = TSComputeIJacobian(ts,glee->stage_time,X,Ydot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
8152302aa1bSDebojyoti Ghosh 
8162302aa1bSDebojyoti Ghosh   ts->dm = dmsave;
8172302aa1bSDebojyoti Ghosh   ierr   = TSGLEERestoreVecs(ts,dm,&Ydot);CHKERRQ(ierr);
8182302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8192302aa1bSDebojyoti Ghosh }
8202302aa1bSDebojyoti Ghosh 
8212302aa1bSDebojyoti Ghosh static PetscErrorCode DMCoarsenHook_TSGLEE(DM fine,DM coarse,void *ctx)
8222302aa1bSDebojyoti Ghosh {
8232302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8242302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8252302aa1bSDebojyoti Ghosh }
8262302aa1bSDebojyoti Ghosh 
8272302aa1bSDebojyoti Ghosh static PetscErrorCode DMRestrictHook_TSGLEE(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
8282302aa1bSDebojyoti Ghosh {
8292302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8302302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8312302aa1bSDebojyoti Ghosh }
8322302aa1bSDebojyoti Ghosh 
8332302aa1bSDebojyoti Ghosh 
8342302aa1bSDebojyoti Ghosh static PetscErrorCode DMSubDomainHook_TSGLEE(DM dm,DM subdm,void *ctx)
8352302aa1bSDebojyoti Ghosh {
8362302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8372302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8382302aa1bSDebojyoti Ghosh }
8392302aa1bSDebojyoti Ghosh 
8402302aa1bSDebojyoti Ghosh static PetscErrorCode DMSubDomainRestrictHook_TSGLEE(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
8412302aa1bSDebojyoti Ghosh {
8422302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8432302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8442302aa1bSDebojyoti Ghosh }
8452302aa1bSDebojyoti Ghosh 
8462302aa1bSDebojyoti Ghosh static PetscErrorCode TSSetUp_GLEE(TS ts)
8472302aa1bSDebojyoti Ghosh {
8482302aa1bSDebojyoti Ghosh   TS_GLEE        *glee = (TS_GLEE*)ts->data;
8492302aa1bSDebojyoti Ghosh   GLEETableau    tab;
8502302aa1bSDebojyoti Ghosh   PetscInt       s,r;
8512302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
8522302aa1bSDebojyoti Ghosh   DM             dm;
8532302aa1bSDebojyoti Ghosh 
8542302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8552302aa1bSDebojyoti Ghosh   if (!glee->tableau) {
8562302aa1bSDebojyoti Ghosh     ierr = TSGLEESetType(ts,TSGLEEDefaultType);CHKERRQ(ierr);
8572302aa1bSDebojyoti Ghosh   }
8582302aa1bSDebojyoti Ghosh   tab  = glee->tableau;
8592302aa1bSDebojyoti Ghosh   s    = tab->s;
8602302aa1bSDebojyoti Ghosh   r    = tab->r;
8612302aa1bSDebojyoti Ghosh   ierr = VecDuplicateVecs(ts->vec_sol,r,&glee->Y);CHKERRQ(ierr);
8622302aa1bSDebojyoti Ghosh   ierr = VecDuplicateVecs(ts->vec_sol,r,&glee->X);CHKERRQ(ierr);
8632302aa1bSDebojyoti Ghosh   ierr = VecDuplicateVecs(ts->vec_sol,s,&glee->YStage);CHKERRQ(ierr);
8642302aa1bSDebojyoti Ghosh   ierr = VecDuplicateVecs(ts->vec_sol,s,&glee->YdotStage);CHKERRQ(ierr);
8652302aa1bSDebojyoti Ghosh   ierr = VecDuplicate(ts->vec_sol,&glee->Ydot);CHKERRQ(ierr);
8660a01e1b2SEmil Constantinescu   ierr = VecDuplicate(ts->vec_sol,&glee->yGErr);CHKERRQ(ierr);
867657c1e31SEmil Constantinescu   ierr = VecZeroEntries(glee->yGErr);CHKERRQ(ierr);
8682302aa1bSDebojyoti Ghosh   ierr = VecDuplicate(ts->vec_sol,&glee->W);CHKERRQ(ierr);
869ec0e3ee2SEmil Constantinescu   ierr = PetscMalloc2(s,&glee->swork,r,&glee->rwork); CHKERRQ(ierr);
8702302aa1bSDebojyoti Ghosh   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
8712302aa1bSDebojyoti Ghosh   ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSGLEE,DMRestrictHook_TSGLEE,ts);CHKERRQ(ierr);
8722302aa1bSDebojyoti Ghosh   ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSGLEE,DMSubDomainRestrictHook_TSGLEE,ts);CHKERRQ(ierr);
8732302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8742302aa1bSDebojyoti Ghosh }
8752302aa1bSDebojyoti Ghosh 
8762302aa1bSDebojyoti Ghosh PetscErrorCode TSStartingMethod_GLEE(TS ts)
8772302aa1bSDebojyoti Ghosh {
8782302aa1bSDebojyoti Ghosh   TS_GLEE        *glee = (TS_GLEE*)ts->data;
8796e2e232bSDebojyoti Ghosh   GLEETableau    tab  = glee->tableau;
8806e2e232bSDebojyoti Ghosh   PetscInt       r=tab->r,i;
8816e2e232bSDebojyoti Ghosh   PetscReal      *S=tab->S;
8822302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
8832302aa1bSDebojyoti Ghosh 
8842302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8852302aa1bSDebojyoti Ghosh   for (i=0; i<r; i++) {
8862302aa1bSDebojyoti Ghosh     ierr = VecZeroEntries(glee->Y[i]);CHKERRQ(ierr);
8872302aa1bSDebojyoti Ghosh     ierr = VecAXPY(glee->Y[i],S[i],ts->vec_sol);CHKERRQ(ierr);
8882302aa1bSDebojyoti Ghosh   }
8892302aa1bSDebojyoti Ghosh 
8902302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8912302aa1bSDebojyoti Ghosh }
8922302aa1bSDebojyoti Ghosh 
8932302aa1bSDebojyoti Ghosh 
8942302aa1bSDebojyoti Ghosh /*------------------------------------------------------------*/
8952302aa1bSDebojyoti Ghosh 
8966b234c65SEmil Constantinescu static PetscErrorCode TSSetFromOptions_GLEE(PetscOptionItems *PetscOptionsObject,TS ts)
8972302aa1bSDebojyoti Ghosh {
8982302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
8992302aa1bSDebojyoti Ghosh   char           gleetype[256];
9002302aa1bSDebojyoti Ghosh 
9012302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
9022302aa1bSDebojyoti Ghosh   ierr = PetscOptionsHead(PetscOptionsObject,"GLEE ODE solver options");CHKERRQ(ierr);
9032302aa1bSDebojyoti Ghosh   {
9042302aa1bSDebojyoti Ghosh     GLEETableauLink link;
9052302aa1bSDebojyoti Ghosh     PetscInt        count,choice;
9062302aa1bSDebojyoti Ghosh     PetscBool       flg;
9072302aa1bSDebojyoti Ghosh     const char      **namelist;
9082302aa1bSDebojyoti Ghosh 
9092302aa1bSDebojyoti Ghosh     ierr = PetscStrncpy(gleetype,TSGLEEDefaultType,sizeof(gleetype));CHKERRQ(ierr);
9102302aa1bSDebojyoti Ghosh     for (link=GLEETableauList,count=0; link; link=link->next,count++) ;
911f489ac74SBarry Smith     ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr);
9122302aa1bSDebojyoti Ghosh     for (link=GLEETableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
9132302aa1bSDebojyoti Ghosh     ierr = PetscOptionsEList("-ts_glee_type","Family of GLEE method","TSGLEESetType",(const char*const*)namelist,count,gleetype,&choice,&flg);CHKERRQ(ierr);
9142302aa1bSDebojyoti Ghosh     ierr = TSGLEESetType(ts,flg ? namelist[choice] : gleetype);CHKERRQ(ierr);
9152302aa1bSDebojyoti Ghosh     ierr = PetscFree(namelist);CHKERRQ(ierr);
9162302aa1bSDebojyoti Ghosh   }
9172302aa1bSDebojyoti Ghosh   ierr = PetscOptionsTail();CHKERRQ(ierr);
9182302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
9192302aa1bSDebojyoti Ghosh }
9202302aa1bSDebojyoti Ghosh 
9212302aa1bSDebojyoti Ghosh static PetscErrorCode TSView_GLEE(TS ts,PetscViewer viewer)
9222302aa1bSDebojyoti Ghosh {
9232302aa1bSDebojyoti Ghosh   TS_GLEE        *glee   = (TS_GLEE*)ts->data;
9242302aa1bSDebojyoti Ghosh   GLEETableau    tab  = glee->tableau;
9252302aa1bSDebojyoti Ghosh   PetscBool      iascii;
9262302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
9272302aa1bSDebojyoti Ghosh 
9282302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
9292302aa1bSDebojyoti Ghosh   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
9302302aa1bSDebojyoti Ghosh   if (iascii) {
9312302aa1bSDebojyoti Ghosh     TSGLEEType gleetype;
9322302aa1bSDebojyoti Ghosh     char   buf[512];
9332302aa1bSDebojyoti Ghosh     ierr = TSGLEEGetType(ts,&gleetype);CHKERRQ(ierr);
934efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  GLEE type %s\n",gleetype);CHKERRQ(ierr);
9352302aa1bSDebojyoti Ghosh     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
9362302aa1bSDebojyoti Ghosh     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa     c = %s\n",buf);CHKERRQ(ierr);
9372302aa1bSDebojyoti Ghosh     /* Note: print out r as well */
9382302aa1bSDebojyoti Ghosh   }
9392302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
9402302aa1bSDebojyoti Ghosh }
9412302aa1bSDebojyoti Ghosh 
9422302aa1bSDebojyoti Ghosh static PetscErrorCode TSLoad_GLEE(TS ts,PetscViewer viewer)
9432302aa1bSDebojyoti Ghosh {
9442302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
9452302aa1bSDebojyoti Ghosh   SNES           snes;
9462302aa1bSDebojyoti Ghosh   TSAdapt        tsadapt;
9472302aa1bSDebojyoti Ghosh 
9482302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
9492302aa1bSDebojyoti Ghosh   ierr = TSGetAdapt(ts,&tsadapt);CHKERRQ(ierr);
9502302aa1bSDebojyoti Ghosh   ierr = TSAdaptLoad(tsadapt,viewer);CHKERRQ(ierr);
9512302aa1bSDebojyoti Ghosh   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
9522302aa1bSDebojyoti Ghosh   ierr = SNESLoad(snes,viewer);CHKERRQ(ierr);
9532302aa1bSDebojyoti Ghosh   /* function and Jacobian context for SNES when used with TS is always ts object */
9542302aa1bSDebojyoti Ghosh   ierr = SNESSetFunction(snes,NULL,NULL,ts);CHKERRQ(ierr);
9552302aa1bSDebojyoti Ghosh   ierr = SNESSetJacobian(snes,NULL,NULL,NULL,ts);CHKERRQ(ierr);
9562302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
9572302aa1bSDebojyoti Ghosh }
9582302aa1bSDebojyoti Ghosh 
9592302aa1bSDebojyoti Ghosh /*@C
9602302aa1bSDebojyoti Ghosh   TSGLEESetType - Set the type of GLEE scheme
9612302aa1bSDebojyoti Ghosh 
9622302aa1bSDebojyoti Ghosh   Logically collective
9632302aa1bSDebojyoti Ghosh 
9642302aa1bSDebojyoti Ghosh   Input Parameter:
9652302aa1bSDebojyoti Ghosh +  ts - timestepping context
9662302aa1bSDebojyoti Ghosh -  gleetype - type of GLEE-scheme
9672302aa1bSDebojyoti Ghosh 
9682302aa1bSDebojyoti Ghosh   Level: intermediate
9692302aa1bSDebojyoti Ghosh 
9702302aa1bSDebojyoti Ghosh .seealso: TSGLEEGetType(), TSGLEE
9712302aa1bSDebojyoti Ghosh @*/
9722302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEESetType(TS ts,TSGLEEType gleetype)
9732302aa1bSDebojyoti Ghosh {
9742302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
9752302aa1bSDebojyoti Ghosh 
9762302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
9772302aa1bSDebojyoti Ghosh   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
978b92453a8SLisandro Dalcin   PetscValidCharPointer(gleetype,2);
9792302aa1bSDebojyoti Ghosh   ierr = PetscTryMethod(ts,"TSGLEESetType_C",(TS,TSGLEEType),(ts,gleetype));CHKERRQ(ierr);
9802302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
9812302aa1bSDebojyoti Ghosh }
9822302aa1bSDebojyoti Ghosh 
9832302aa1bSDebojyoti Ghosh /*@C
9842302aa1bSDebojyoti Ghosh   TSGLEEGetType - Get the type of GLEE scheme
9852302aa1bSDebojyoti Ghosh 
9862302aa1bSDebojyoti Ghosh   Logically collective
9872302aa1bSDebojyoti Ghosh 
9882302aa1bSDebojyoti Ghosh   Input Parameter:
9892302aa1bSDebojyoti Ghosh .  ts - timestepping context
9902302aa1bSDebojyoti Ghosh 
9912302aa1bSDebojyoti Ghosh   Output Parameter:
9922302aa1bSDebojyoti Ghosh .  gleetype - type of GLEE-scheme
9932302aa1bSDebojyoti Ghosh 
9942302aa1bSDebojyoti Ghosh   Level: intermediate
9952302aa1bSDebojyoti Ghosh 
9962302aa1bSDebojyoti Ghosh .seealso: TSGLEESetType()
9972302aa1bSDebojyoti Ghosh @*/
9982302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEEGetType(TS ts,TSGLEEType *gleetype)
9992302aa1bSDebojyoti Ghosh {
10002302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
10012302aa1bSDebojyoti Ghosh 
10022302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
10032302aa1bSDebojyoti Ghosh   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
10042302aa1bSDebojyoti Ghosh   ierr = PetscUseMethod(ts,"TSGLEEGetType_C",(TS,TSGLEEType*),(ts,gleetype));CHKERRQ(ierr);
10052302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
10062302aa1bSDebojyoti Ghosh }
10072302aa1bSDebojyoti Ghosh 
10082302aa1bSDebojyoti Ghosh PetscErrorCode  TSGLEEGetType_GLEE(TS ts,TSGLEEType *gleetype)
10092302aa1bSDebojyoti Ghosh {
10102302aa1bSDebojyoti Ghosh   TS_GLEE     *glee = (TS_GLEE*)ts->data;
10112302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
10122302aa1bSDebojyoti Ghosh 
10132302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
10142302aa1bSDebojyoti Ghosh   if (!glee->tableau) {
10152302aa1bSDebojyoti Ghosh     ierr = TSGLEESetType(ts,TSGLEEDefaultType);CHKERRQ(ierr);
10162302aa1bSDebojyoti Ghosh   }
10172302aa1bSDebojyoti Ghosh   *gleetype = glee->tableau->name;
10182302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
10192302aa1bSDebojyoti Ghosh }
10202302aa1bSDebojyoti Ghosh PetscErrorCode  TSGLEESetType_GLEE(TS ts,TSGLEEType gleetype)
10212302aa1bSDebojyoti Ghosh {
10222302aa1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
10232302aa1bSDebojyoti Ghosh   PetscErrorCode  ierr;
10242302aa1bSDebojyoti Ghosh   PetscBool       match;
10252302aa1bSDebojyoti Ghosh   GLEETableauLink link;
10262302aa1bSDebojyoti Ghosh 
10272302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
10282302aa1bSDebojyoti Ghosh   if (glee->tableau) {
10292302aa1bSDebojyoti Ghosh     ierr = PetscStrcmp(glee->tableau->name,gleetype,&match);CHKERRQ(ierr);
10302302aa1bSDebojyoti Ghosh     if (match) PetscFunctionReturn(0);
10312302aa1bSDebojyoti Ghosh   }
10322302aa1bSDebojyoti Ghosh   for (link = GLEETableauList; link; link=link->next) {
10332302aa1bSDebojyoti Ghosh     ierr = PetscStrcmp(link->tab.name,gleetype,&match);CHKERRQ(ierr);
10342302aa1bSDebojyoti Ghosh     if (match) {
10352302aa1bSDebojyoti Ghosh       ierr = TSReset_GLEE(ts);CHKERRQ(ierr);
10362302aa1bSDebojyoti Ghosh       glee->tableau = &link->tab;
10372302aa1bSDebojyoti Ghosh       PetscFunctionReturn(0);
10382302aa1bSDebojyoti Ghosh     }
10392302aa1bSDebojyoti Ghosh   }
10402302aa1bSDebojyoti Ghosh   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",gleetype);
10412302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
10422302aa1bSDebojyoti Ghosh }
10432302aa1bSDebojyoti Ghosh 
10442302aa1bSDebojyoti Ghosh static PetscErrorCode  TSGetStages_GLEE(TS ts,PetscInt *ns,Vec **Y)
10452302aa1bSDebojyoti Ghosh {
10462302aa1bSDebojyoti Ghosh   TS_GLEE *glee = (TS_GLEE*)ts->data;
10472302aa1bSDebojyoti Ghosh 
10482302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
1049*0429704eSStefano Zampini   if (ns) *ns = glee->tableau->s;
1050d5509560SEmil Constantinescu   if (Y)  *Y  = glee->YStage;
10512302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
10522302aa1bSDebojyoti Ghosh }
10532302aa1bSDebojyoti Ghosh 
1054b2bf4f3aSDebojyoti Ghosh PetscErrorCode TSGetSolutionComponents_GLEE(TS ts,PetscInt *n,Vec *Y)
10552302aa1bSDebojyoti Ghosh {
10562302aa1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
10572302aa1bSDebojyoti Ghosh   GLEETableau     tab   = glee->tableau;
10586e2e232bSDebojyoti Ghosh   PetscErrorCode  ierr;
10592302aa1bSDebojyoti Ghosh 
10602302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
10614cdd57e5SDebojyoti Ghosh   if (!Y) *n = tab->r;
10622302aa1bSDebojyoti Ghosh   else {
10634cdd57e5SDebojyoti Ghosh     if ((*n >= 0) && (*n < tab->r)) {
10644cdd57e5SDebojyoti Ghosh       ierr = VecCopy(glee->Y[*n],*Y); CHKERRQ(ierr);
10659657682dSDebojyoti Ghosh     } else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Second argument (%d) out of range[%d,%d].",*n,0,tab->r-1);
10662302aa1bSDebojyoti Ghosh   }
10672302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
10682302aa1bSDebojyoti Ghosh }
10692302aa1bSDebojyoti Ghosh 
10704cdd57e5SDebojyoti Ghosh PetscErrorCode TSGetAuxSolution_GLEE(TS ts,Vec *X)
10714cdd57e5SDebojyoti Ghosh {
10724cdd57e5SDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
10734cdd57e5SDebojyoti Ghosh   GLEETableau     tab   = glee->tableau;
10744cdd57e5SDebojyoti Ghosh   PetscReal       *F    = tab->Fembed;
10754cdd57e5SDebojyoti Ghosh   PetscInt        r     = tab->r;
10764cdd57e5SDebojyoti Ghosh   Vec             *Y    = glee->Y;
1077ec0e3ee2SEmil Constantinescu   PetscScalar     *wr   = glee->rwork;
1078ec0e3ee2SEmil Constantinescu   PetscInt        i;
10794cdd57e5SDebojyoti Ghosh   PetscErrorCode  ierr;
10804cdd57e5SDebojyoti Ghosh 
10814cdd57e5SDebojyoti Ghosh   PetscFunctionBegin;
10824cdd57e5SDebojyoti Ghosh   ierr = VecZeroEntries(*X);CHKERRQ(ierr);
1083ec0e3ee2SEmil Constantinescu   for (i=0; i<r; i++) wr[i] = F[i];
1084ec0e3ee2SEmil Constantinescu   ierr = VecMAXPY((*X),r,wr,Y);CHKERRQ(ierr);
10854cdd57e5SDebojyoti Ghosh   PetscFunctionReturn(0);
10864cdd57e5SDebojyoti Ghosh }
10874cdd57e5SDebojyoti Ghosh 
10880a01e1b2SEmil Constantinescu PetscErrorCode TSGetTimeError_GLEE(TS ts,PetscInt n,Vec *X)
10894cdd57e5SDebojyoti Ghosh {
10904cdd57e5SDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
10914cdd57e5SDebojyoti Ghosh   GLEETableau     tab   = glee->tableau;
10924cdd57e5SDebojyoti Ghosh   PetscReal       *F    = tab->Ferror;
10934cdd57e5SDebojyoti Ghosh   PetscInt        r     = tab->r;
10944cdd57e5SDebojyoti Ghosh   Vec             *Y    = glee->Y;
1095ec0e3ee2SEmil Constantinescu   PetscScalar     *wr   = glee->rwork;
1096ec0e3ee2SEmil Constantinescu   PetscInt        i;
10974cdd57e5SDebojyoti Ghosh   PetscErrorCode  ierr;
10984cdd57e5SDebojyoti Ghosh 
10994cdd57e5SDebojyoti Ghosh   PetscFunctionBegin;
11004cdd57e5SDebojyoti Ghosh   ierr = VecZeroEntries(*X);CHKERRQ(ierr);
1101657c1e31SEmil Constantinescu   if(n==0){
1102ec0e3ee2SEmil Constantinescu     for (i=0; i<r; i++) wr[i] = F[i];
1103ec0e3ee2SEmil Constantinescu     ierr = VecMAXPY((*X),r,wr,Y);CHKERRQ(ierr);
11040a01e1b2SEmil Constantinescu   } else if(n==-1) {
1105657c1e31SEmil Constantinescu     *X=glee->yGErr;
11060a01e1b2SEmil Constantinescu   }
11074cdd57e5SDebojyoti Ghosh   PetscFunctionReturn(0);
11084cdd57e5SDebojyoti Ghosh }
11094cdd57e5SDebojyoti Ghosh 
111057df6a1bSDebojyoti Ghosh PetscErrorCode TSSetTimeError_GLEE(TS ts,Vec X)
111157df6a1bSDebojyoti Ghosh {
111257df6a1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
111357df6a1bSDebojyoti Ghosh   GLEETableau     tab   = glee->tableau;
111457df6a1bSDebojyoti Ghosh   PetscReal       *S    = tab->Serror;
111557df6a1bSDebojyoti Ghosh   PetscInt        r     = tab->r,i;
111657df6a1bSDebojyoti Ghosh   Vec             *Y    = glee->Y;
111757df6a1bSDebojyoti Ghosh   PetscErrorCode  ierr;
111857df6a1bSDebojyoti Ghosh 
111957df6a1bSDebojyoti Ghosh   PetscFunctionBegin;
11209657682dSDebojyoti Ghosh   if (r != 2) SETERRQ2(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSSetTimeError_GLEE not supported for '%s' with r=%D.",tab->name,tab->r);
112157df6a1bSDebojyoti Ghosh   for (i=1; i<r; i++) {
112257df6a1bSDebojyoti Ghosh     ierr = VecCopy(ts->vec_sol,Y[i]); CHKERRQ(ierr);
112357df6a1bSDebojyoti Ghosh     ierr = VecAXPBY(Y[i],S[0],S[1],X); CHKERRQ(ierr);
1124642ca4e8SEmil Constantinescu     ierr = VecCopy(X,glee->yGErr); CHKERRQ(ierr);
112557df6a1bSDebojyoti Ghosh   }
112657df6a1bSDebojyoti Ghosh   PetscFunctionReturn(0);
112757df6a1bSDebojyoti Ghosh }
112857df6a1bSDebojyoti Ghosh 
1129b3a6b972SJed Brown static PetscErrorCode TSDestroy_GLEE(TS ts)
1130b3a6b972SJed Brown {
1131b3a6b972SJed Brown   PetscErrorCode ierr;
1132b3a6b972SJed Brown 
1133b3a6b972SJed Brown   PetscFunctionBegin;
1134b3a6b972SJed Brown   ierr = TSReset_GLEE(ts);CHKERRQ(ierr);
1135b3a6b972SJed Brown   if (ts->dm) {
1136b3a6b972SJed Brown     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSGLEE,DMRestrictHook_TSGLEE,ts);CHKERRQ(ierr);
1137b3a6b972SJed Brown     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSGLEE,DMSubDomainRestrictHook_TSGLEE,ts);CHKERRQ(ierr);
1138b3a6b972SJed Brown   }
1139b3a6b972SJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
1140b3a6b972SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLEEGetType_C",NULL);CHKERRQ(ierr);
1141b3a6b972SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLEESetType_C",NULL);CHKERRQ(ierr);
1142b3a6b972SJed Brown   PetscFunctionReturn(0);
1143b3a6b972SJed Brown }
1144b3a6b972SJed Brown 
11452302aa1bSDebojyoti Ghosh /* ------------------------------------------------------------ */
11462302aa1bSDebojyoti Ghosh /*MC
11472302aa1bSDebojyoti Ghosh       TSGLEE - ODE and DAE solver using General Linear with Error Estimation schemes
11482302aa1bSDebojyoti Ghosh 
11492302aa1bSDebojyoti Ghosh   The user should provide the right hand side of the equation
11502302aa1bSDebojyoti Ghosh   using TSSetRHSFunction().
11512302aa1bSDebojyoti Ghosh 
11522302aa1bSDebojyoti Ghosh   Notes:
11532302aa1bSDebojyoti Ghosh   The default is TSGLEE35, it can be changed with TSGLEESetType() or -ts_glee_type
11542302aa1bSDebojyoti Ghosh 
11552302aa1bSDebojyoti Ghosh   Level: beginner
11562302aa1bSDebojyoti Ghosh 
11572302aa1bSDebojyoti Ghosh .seealso:  TSCreate(), TS, TSSetType(), TSGLEESetType(), TSGLEEGetType(),
11582302aa1bSDebojyoti Ghosh            TSGLEE23, TTSGLEE24, TSGLEE35, TSGLEE25I, TSGLEEEXRK2A,
11592302aa1bSDebojyoti Ghosh            TSGLEERK32G1, TSGLEERK285EX, TSGLEEType, TSGLEERegister()
11602302aa1bSDebojyoti Ghosh 
11612302aa1bSDebojyoti Ghosh M*/
11622302aa1bSDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSCreate_GLEE(TS ts)
11632302aa1bSDebojyoti Ghosh {
11642302aa1bSDebojyoti Ghosh   TS_GLEE         *th;
11652302aa1bSDebojyoti Ghosh   PetscErrorCode  ierr;
11662302aa1bSDebojyoti Ghosh 
11672302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
11682302aa1bSDebojyoti Ghosh #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
11692302aa1bSDebojyoti Ghosh   ierr = TSGLEEInitializePackage();CHKERRQ(ierr);
11702302aa1bSDebojyoti Ghosh #endif
11712302aa1bSDebojyoti Ghosh 
11722302aa1bSDebojyoti Ghosh   ts->ops->reset                  = TSReset_GLEE;
11732302aa1bSDebojyoti Ghosh   ts->ops->destroy                = TSDestroy_GLEE;
11742302aa1bSDebojyoti Ghosh   ts->ops->view                   = TSView_GLEE;
11752302aa1bSDebojyoti Ghosh   ts->ops->load                   = TSLoad_GLEE;
11762302aa1bSDebojyoti Ghosh   ts->ops->setup                  = TSSetUp_GLEE;
11772302aa1bSDebojyoti Ghosh   ts->ops->step                   = TSStep_GLEE;
11782302aa1bSDebojyoti Ghosh   ts->ops->interpolate            = TSInterpolate_GLEE;
11792302aa1bSDebojyoti Ghosh   ts->ops->evaluatestep           = TSEvaluateStep_GLEE;
11802302aa1bSDebojyoti Ghosh   ts->ops->setfromoptions         = TSSetFromOptions_GLEE;
11812302aa1bSDebojyoti Ghosh   ts->ops->getstages              = TSGetStages_GLEE;
11822302aa1bSDebojyoti Ghosh   ts->ops->snesfunction           = SNESTSFormFunction_GLEE;
11832302aa1bSDebojyoti Ghosh   ts->ops->snesjacobian           = SNESTSFormJacobian_GLEE;
1184b2bf4f3aSDebojyoti Ghosh   ts->ops->getsolutioncomponents  = TSGetSolutionComponents_GLEE;
11854cdd57e5SDebojyoti Ghosh   ts->ops->getauxsolution         = TSGetAuxSolution_GLEE;
11864cdd57e5SDebojyoti Ghosh   ts->ops->gettimeerror           = TSGetTimeError_GLEE;
118757df6a1bSDebojyoti Ghosh   ts->ops->settimeerror           = TSSetTimeError_GLEE;
11882302aa1bSDebojyoti Ghosh   ts->ops->startingmethod         = TSStartingMethod_GLEE;
1189b80d5313SLisandro Dalcin   ts->default_adapt_type          = TSADAPTGLEE;
11902302aa1bSDebojyoti Ghosh 
1191efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
1192efd4aadfSBarry Smith 
11932302aa1bSDebojyoti Ghosh   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
11942302aa1bSDebojyoti Ghosh   ts->data = (void*)th;
11952302aa1bSDebojyoti Ghosh 
11962302aa1bSDebojyoti Ghosh   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLEEGetType_C",TSGLEEGetType_GLEE);CHKERRQ(ierr);
11972302aa1bSDebojyoti Ghosh   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLEESetType_C",TSGLEESetType_GLEE);CHKERRQ(ierr);
11982302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
11992302aa1bSDebojyoti Ghosh }
1200