xref: /petsc/src/ts/impls/glee/glee.c (revision b1fbfd333c6a3a8f7cd29f5fb271e978ee11a1cb)
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
752302aa1bSDebojyoti Ghosh */
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
852302aa1bSDebojyoti Ghosh */
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
952302aa1bSDebojyoti Ghosh */
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
1052302aa1bSDebojyoti Ghosh */
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
1152302aa1bSDebojyoti Ghosh */
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
1252302aa1bSDebojyoti Ghosh */
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
1352302aa1bSDebojyoti Ghosh */
1362302aa1bSDebojyoti Ghosh 
1372302aa1bSDebojyoti Ghosh #undef __FUNCT__
1382302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEERegisterAll"
1392302aa1bSDebojyoti Ghosh /*@C
1402302aa1bSDebojyoti Ghosh   TSGLEERegisterAll - Registers all of the General Linear with Error Estimation methods in TSGLEE
1412302aa1bSDebojyoti Ghosh 
1422302aa1bSDebojyoti Ghosh   Not Collective, but should be called by all processes which will need the schemes to be registered
1432302aa1bSDebojyoti Ghosh 
1442302aa1bSDebojyoti Ghosh   Level: advanced
1452302aa1bSDebojyoti Ghosh 
1462302aa1bSDebojyoti Ghosh .keywords: TS, TSGLEE, register, all
1472302aa1bSDebojyoti Ghosh 
1482302aa1bSDebojyoti Ghosh .seealso:  TSGLEERegisterDestroy()
1492302aa1bSDebojyoti Ghosh @*/
1502302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEERegisterAll(void)
1512302aa1bSDebojyoti Ghosh {
1522302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
1532302aa1bSDebojyoti Ghosh 
1542302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
1552302aa1bSDebojyoti Ghosh   if (TSGLEERegisterAllCalled) PetscFunctionReturn(0);
1562302aa1bSDebojyoti Ghosh   TSGLEERegisterAllCalled = PETSC_TRUE;
1572302aa1bSDebojyoti Ghosh 
1582302aa1bSDebojyoti Ghosh   {
159*b1fbfd33SSatish Balay #define GAMMA 0.5
1602302aa1bSDebojyoti Ghosh     /* y-eps form */
1612302aa1bSDebojyoti Ghosh     const PetscInt
162ab8c5912SEmil Constantinescu       p = 1,
163ab8c5912SEmil Constantinescu       s = 3,
164ab8c5912SEmil Constantinescu       r = 2;
165ab8c5912SEmil Constantinescu     const PetscReal
166ab8c5912SEmil Constantinescu       A[3][3]   = {{1.0,0,0},{0,0.5,0},{0,0.5,0.5}},
167ab8c5912SEmil Constantinescu       B[2][3]   = {{1.0,0,0},{-2.0,1.0,1.0}},
168ab8c5912SEmil Constantinescu       U[3][2]   = {{1.0,0},{1.0,0.5},{1.0,0.5}},
169ab8c5912SEmil Constantinescu       V[2][2]   = {{1,0},{0,1}},
170ab8c5912SEmil Constantinescu       S[2]      = {1,0},
171ab8c5912SEmil Constantinescu       F[2]      = {1,0},
172*b1fbfd33SSatish Balay       Fembed[2] = {1,1-GAMMA},
17357df6a1bSDebojyoti Ghosh       Ferror[2] = {0,1},
17457df6a1bSDebojyoti Ghosh       Serror[2] = {1,0};
175*b1fbfd33SSatish 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);
176ab8c5912SEmil Constantinescu   }
177ab8c5912SEmil Constantinescu   {
178*b1fbfd33SSatish Balay #undef GAMMA
179*b1fbfd33SSatish Balay #define GAMMA 0.0
180ab8c5912SEmil Constantinescu     /* y-eps form */
181ab8c5912SEmil Constantinescu     const PetscInt
1822302aa1bSDebojyoti Ghosh       p = 2,
1832302aa1bSDebojyoti Ghosh       s = 3,
1842302aa1bSDebojyoti Ghosh       r = 2;
1852302aa1bSDebojyoti Ghosh     const PetscReal
1862302aa1bSDebojyoti Ghosh       A[3][3]   = {{0,0,0},{1,0,0},{0.25,0.25,0}},
1872302aa1bSDebojyoti 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}},
1882302aa1bSDebojyoti Ghosh       U[3][2]   = {{1,0},{1,10},{1,-1}},
1892302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
1902302aa1bSDebojyoti Ghosh       S[2]      = {1,0},
1912302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
192*b1fbfd33SSatish Balay       Fembed[2] = {1,1-GAMMA},
19357df6a1bSDebojyoti Ghosh       Ferror[2] = {0,1},
19457df6a1bSDebojyoti Ghosh       Serror[2] = {1,0};
195*b1fbfd33SSatish 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);
1962302aa1bSDebojyoti Ghosh   }
1972302aa1bSDebojyoti Ghosh   {
198*b1fbfd33SSatish Balay #undef GAMMA
199*b1fbfd33SSatish Balay #define GAMMA 0.0
2002302aa1bSDebojyoti Ghosh     /* y-y~ form */
2012302aa1bSDebojyoti Ghosh     const PetscInt
2022302aa1bSDebojyoti Ghosh       p = 2,
2032302aa1bSDebojyoti Ghosh       s = 4,
2042302aa1bSDebojyoti Ghosh       r = 2;
2052302aa1bSDebojyoti Ghosh     const PetscReal
2062302aa1bSDebojyoti 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}},
2072302aa1bSDebojyoti 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}},
2082302aa1bSDebojyoti Ghosh       U[4][2]   = {{0,1},{75.0/58.0,-17.0/58.0},{0,1},{0,1}},
2092302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
2102302aa1bSDebojyoti Ghosh       S[2]      = {1,1},
2112302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
212fd96ebc2SDebojyoti Ghosh       Fembed[2] = {0,1},
213*b1fbfd33SSatish Balay       Ferror[2] = {-1.0/(1.0-GAMMA),1.0/(1.0-GAMMA)},
214*b1fbfd33SSatish Balay       Serror[2] = {1.0-GAMMA,1.0};
215*b1fbfd33SSatish 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);
2162302aa1bSDebojyoti Ghosh   }
2172302aa1bSDebojyoti Ghosh   {
218*b1fbfd33SSatish Balay #undef GAMMA
219*b1fbfd33SSatish Balay #define GAMMA 0.0
2202302aa1bSDebojyoti Ghosh     /* y-y~ form */
2212302aa1bSDebojyoti Ghosh     const PetscInt
2222302aa1bSDebojyoti Ghosh       p = 2,
2232302aa1bSDebojyoti Ghosh       s = 5,
2242302aa1bSDebojyoti Ghosh       r = 2;
2252302aa1bSDebojyoti Ghosh     const PetscReal
2262302aa1bSDebojyoti Ghosh       A[5][5]   = {{0,0,0,0,0},
2272302aa1bSDebojyoti Ghosh                    {-0.94079244066783383269,0,0,0,0},
2282302aa1bSDebojyoti Ghosh                    {0.64228187778301907108,0.10915356933958500042,0,0,0},
2292302aa1bSDebojyoti Ghosh                    {-0.51764297742287450812,0.74414270351096040738,-0.71404164927824538121,0,0},
2302302aa1bSDebojyoti Ghosh                    {-0.44696561556825969206,-0.76768425657590196518,0.20111608138142987881,0.93828186737840469796,0}},
2312302aa1bSDebojyoti Ghosh       B[2][5]   = {{-0.029309178948150356153,-0.49671981884013874923,0.34275801517650053274,0.32941112623949194988,0.85385985637229662276},
2322302aa1bSDebojyoti Ghosh                    {0.78133219686062535272,0.074238691892675897635,0.57957363498384957966,-0.24638502829674959968,-0.18875949544040123033}},
2332302aa1bSDebojyoti Ghosh       U[5][2]   = {{0.16911424754448327735,0.83088575245551672265},
2342302aa1bSDebojyoti Ghosh                    {0.53638465733199574340,0.46361534266800425660},
2352302aa1bSDebojyoti Ghosh                    {0.39901579167169582526,0.60098420832830417474},
2362302aa1bSDebojyoti Ghosh                    {0.87689005530618575480,0.12310994469381424520},
2372302aa1bSDebojyoti Ghosh                    {0.99056100455550913009,0.0094389954444908699092}},
2382302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
2392302aa1bSDebojyoti Ghosh       S[2]      = {1,1},
2402302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
241fd96ebc2SDebojyoti Ghosh       Fembed[2] = {0,1},
242*b1fbfd33SSatish Balay       Ferror[2] = {-1.0/(1.0-GAMMA),1.0/(1.0-GAMMA)},
243*b1fbfd33SSatish Balay       Serror[2] = {1.0-GAMMA,1.0};
244*b1fbfd33SSatish 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);
2452302aa1bSDebojyoti Ghosh   }
2462302aa1bSDebojyoti Ghosh   {
247*b1fbfd33SSatish Balay #undef GAMMA
248*b1fbfd33SSatish Balay #define GAMMA 0.0
2492302aa1bSDebojyoti Ghosh     /* y-y~ form */
2502302aa1bSDebojyoti Ghosh     const PetscInt
2512302aa1bSDebojyoti Ghosh       p = 3,
2522302aa1bSDebojyoti Ghosh       s = 5,
2532302aa1bSDebojyoti Ghosh       r = 2;
2542302aa1bSDebojyoti Ghosh     const PetscReal
2552302aa1bSDebojyoti Ghosh       A[5][5]   = {{0,0,0,0,0},
2562302aa1bSDebojyoti Ghosh                    {- 2169604947363702313.0 /  24313474998937147335.0,0,0,0,0},
2572302aa1bSDebojyoti Ghosh                    {46526746497697123895.0 /  94116917485856474137.0,-10297879244026594958.0 /  49199457603717988219.0,0,0,0},
2582302aa1bSDebojyoti Ghosh                    {23364788935845982499.0 /  87425311444725389446.0,-79205144337496116638.0 / 148994349441340815519.0,40051189859317443782.0 /  36487615018004984309.0,0,0},
25917d8b14cSSatish Balay                    {42089522664062539205.0 / 124911313006412840286.0,-15074384760342762939.0 / 137927286865289746282.0,-62274678522253371016.0 / 125918573676298591413.0,13755475729852471739.0 /  79257927066651693390.0,0}},
2602302aa1bSDebojyoti Ghosh       B[2][5]   = {{61546696837458703723.0 /  56982519523786160813.0,-55810892792806293355.0 / 206957624151308356511.0,24061048952676379087.0 / 158739347956038723465.0,3577972206874351339.0 /   7599733370677197135.0,-59449832954780563947.0 / 137360038685338563670.0},
2612302aa1bSDebojyoti Ghosh                    {- 9738262186984159168.0 /  99299082461487742983.0,-32797097931948613195.0 /  61521565616362163366.0,42895514606418420631.0 /  71714201188501437336.0,22608567633166065068.0 /  55371917805607957003.0,94655809487476459565.0 / 151517167160302729021.0}},
2622302aa1bSDebojyoti Ghosh       U[5][2]   = {{70820309139834661559.0 /  80863923579509469826.0,10043614439674808267.0 /  80863923579509469826.0},
2632302aa1bSDebojyoti Ghosh                    {161694774978034105510.0 / 106187653640211060371.0,-55507121337823045139.0 / 106187653640211060371.0},
2642302aa1bSDebojyoti Ghosh                    {78486094644566264568.0 /  88171030896733822981.0,9684936252167558413.0 /  88171030896733822981.0},
2652302aa1bSDebojyoti Ghosh                    {65394922146334854435.0 /  84570853840405479554.0,19175931694070625119.0 /  84570853840405479554.0},
2662302aa1bSDebojyoti Ghosh                    {8607282770183754108.0 / 108658046436496925911.0,100050763666313171803.0 / 108658046436496925911.0}},
2672302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
2682302aa1bSDebojyoti Ghosh       S[2]      = {1,1},
2692302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
270fd96ebc2SDebojyoti Ghosh       Fembed[2] = {0,1},
271*b1fbfd33SSatish Balay       Ferror[2] = {-1.0/(1.0-GAMMA),1.0/(1.0-GAMMA)},
272*b1fbfd33SSatish Balay       Serror[2] = {1.0-GAMMA,1.0};
273*b1fbfd33SSatish 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);
2742302aa1bSDebojyoti Ghosh   }
2752302aa1bSDebojyoti Ghosh   {
276*b1fbfd33SSatish Balay #undef GAMMA
277*b1fbfd33SSatish Balay #define GAMMA 0.25
2782302aa1bSDebojyoti Ghosh     /* y-eps form */
2792302aa1bSDebojyoti Ghosh     const PetscInt
2802302aa1bSDebojyoti Ghosh       p = 2,
2812302aa1bSDebojyoti Ghosh       s = 6,
2822302aa1bSDebojyoti Ghosh       r = 2;
2832302aa1bSDebojyoti Ghosh     const PetscReal
2842302aa1bSDebojyoti Ghosh       A[6][6]   = {{0,0,0,0,0,0},
2852302aa1bSDebojyoti Ghosh                    {1,0,0,0,0,0},
2862302aa1bSDebojyoti Ghosh                    {0,0,0,0,0,0},
2872302aa1bSDebojyoti Ghosh                    {0,0,0.5,0,0,0},
2882302aa1bSDebojyoti Ghosh                    {0,0,0.25,0.25,0,0},
2892302aa1bSDebojyoti Ghosh                    {0,0,0.25,0.25,0.5,0}},
2902302aa1bSDebojyoti Ghosh       B[2][6]   = {{0.5,0.5,0,0,0,0},
2912302aa1bSDebojyoti 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}},
2922302aa1bSDebojyoti Ghosh       U[6][2]   = {{1,0},{1,0},{1,0.75},{1,0.75},{1,0.75},{1,0.75}},
2932302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
2942302aa1bSDebojyoti Ghosh       S[2]      = {1,0},
2952302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
296*b1fbfd33SSatish Balay       Fembed[2] = {1,1-GAMMA},
29757df6a1bSDebojyoti Ghosh       Ferror[2] = {0,1},
29857df6a1bSDebojyoti Ghosh       Serror[2] = {1,0};
299*b1fbfd33SSatish 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);
3002302aa1bSDebojyoti Ghosh   }
3012302aa1bSDebojyoti Ghosh   {
302*b1fbfd33SSatish Balay #undef GAMMA
303*b1fbfd33SSatish Balay #define GAMMA 0.0
3042302aa1bSDebojyoti Ghosh     /* y-eps form */
3052302aa1bSDebojyoti Ghosh     const PetscInt
3062302aa1bSDebojyoti Ghosh       p = 3,
3072302aa1bSDebojyoti Ghosh       s = 8,
3082302aa1bSDebojyoti Ghosh       r = 2;
3092302aa1bSDebojyoti Ghosh     const PetscReal
3102302aa1bSDebojyoti Ghosh       A[8][8]   = {{0,0,0,0,0,0,0,0},
3112302aa1bSDebojyoti Ghosh                    {0.5,0,0,0,0,0,0,0},
3122302aa1bSDebojyoti Ghosh                    {-1,2,0,0,0,0,0,0},
3132302aa1bSDebojyoti Ghosh                    {1.0/6.0,2.0/3.0,1.0/6.0,0,0,0,0,0},
3142302aa1bSDebojyoti Ghosh                    {0,0,0,0,0,0,0,0},
3152302aa1bSDebojyoti Ghosh                    {-7.0/24.0,1.0/3.0,1.0/12.0,-1.0/8.0,0.5,0,0,0},
3162302aa1bSDebojyoti Ghosh                    {7.0/6.0,-4.0/3.0,-1.0/3.0,0.5,-1.0,2.0,0,0},
3172302aa1bSDebojyoti Ghosh                    {0,0,0,0,1.0/6.0,2.0/3.0,1.0/6.0,0}},
3182302aa1bSDebojyoti Ghosh       B[2][8]   = {{1.0/6.0,2.0/3.0,1.0/6.0,0,0,0,0,0},
3192302aa1bSDebojyoti 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}},
3202302aa1bSDebojyoti Ghosh       U[8][2]   = {{1,0},{1,0},{1,0},{1,0},{1,1},{1,1},{1,1},{1,1}},
3212302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
3222302aa1bSDebojyoti Ghosh       S[2]      = {1,0},
3232302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
324*b1fbfd33SSatish Balay       Fembed[2] = {1,1-GAMMA},
32557df6a1bSDebojyoti Ghosh       Ferror[2] = {0,1},
32657df6a1bSDebojyoti Ghosh       Serror[2] = {1,0};
327*b1fbfd33SSatish 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);
3282302aa1bSDebojyoti Ghosh   }
3292302aa1bSDebojyoti Ghosh   {
330*b1fbfd33SSatish Balay #undef GAMMA
331*b1fbfd33SSatish Balay #define GAMMA 0.25
3322302aa1bSDebojyoti Ghosh     /* y-eps form */
3332302aa1bSDebojyoti Ghosh     const PetscInt
3342302aa1bSDebojyoti Ghosh       p = 2,
3352302aa1bSDebojyoti Ghosh       s = 9,
3362302aa1bSDebojyoti Ghosh       r = 2;
3372302aa1bSDebojyoti Ghosh     const PetscReal
3382302aa1bSDebojyoti Ghosh       A[9][9]   = {{0,0,0,0,0,0,0,0,0},
3392302aa1bSDebojyoti Ghosh                    {0.585786437626904966,0,0,0,0,0,0,0,0},
3402302aa1bSDebojyoti Ghosh                    {0.149999999999999994,0.849999999999999978,0,0,0,0,0,0,0},
3412302aa1bSDebojyoti Ghosh                    {0,0,0,0,0,0,0,0,0},
3422302aa1bSDebojyoti Ghosh                    {0,0,0,0.292893218813452483,0,0,0,0,0},
3432302aa1bSDebojyoti Ghosh                    {0,0,0,0.0749999999999999972,0.424999999999999989,0,0,0,0},
3442302aa1bSDebojyoti Ghosh                    {0,0,0,0.176776695296636893,0.176776695296636893,0.146446609406726241,0,0,0},
3452302aa1bSDebojyoti Ghosh                    {0,0,0,0.176776695296636893,0.176776695296636893,0.146446609406726241,0.292893218813452483,0,0},
3462302aa1bSDebojyoti Ghosh                    {0,0,0,0.176776695296636893,0.176776695296636893,0.146446609406726241,0.0749999999999999972,0.424999999999999989,0}},
3472302aa1bSDebojyoti Ghosh       B[2][9]   = {{0.353553390593273786,0.353553390593273786,0.292893218813452483,0,0,0,0,0,0},
3482302aa1bSDebojyoti Ghosh                    {-0.471404520791031678,-0.471404520791031678,-0.390524291751269959,0.235702260395515839,0.235702260395515839,0.195262145875634979,0.235702260395515839,0.235702260395515839,0.195262145875634979}},
3492302aa1bSDebojyoti 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}},
3502302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
3512302aa1bSDebojyoti Ghosh       S[2]      = {1,0},
3522302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
353*b1fbfd33SSatish Balay       Fembed[2] = {1,1-GAMMA},
35457df6a1bSDebojyoti Ghosh       Ferror[2] = {0,1},
35557df6a1bSDebojyoti Ghosh       Serror[2] = {1,0};
356*b1fbfd33SSatish 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);
3572302aa1bSDebojyoti Ghosh   }
3582302aa1bSDebojyoti Ghosh 
3592302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
3602302aa1bSDebojyoti Ghosh }
3612302aa1bSDebojyoti Ghosh 
3622302aa1bSDebojyoti Ghosh #undef __FUNCT__
3632302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEERegisterDestroy"
3642302aa1bSDebojyoti Ghosh /*@C
3652302aa1bSDebojyoti Ghosh    TSGLEERegisterDestroy - Frees the list of schemes that were registered by TSGLEERegister().
3662302aa1bSDebojyoti Ghosh 
3672302aa1bSDebojyoti Ghosh    Not Collective
3682302aa1bSDebojyoti Ghosh 
3692302aa1bSDebojyoti Ghosh    Level: advanced
3702302aa1bSDebojyoti Ghosh 
3712302aa1bSDebojyoti Ghosh .keywords: TSGLEE, register, destroy
3722302aa1bSDebojyoti Ghosh .seealso: TSGLEERegister(), TSGLEERegisterAll()
3732302aa1bSDebojyoti Ghosh @*/
3742302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEERegisterDestroy(void)
3752302aa1bSDebojyoti Ghosh {
3762302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
3772302aa1bSDebojyoti Ghosh   GLEETableauLink link;
3782302aa1bSDebojyoti Ghosh 
3792302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
3802302aa1bSDebojyoti Ghosh   while ((link = GLEETableauList)) {
3812302aa1bSDebojyoti Ghosh     GLEETableau t = &link->tab;
3822302aa1bSDebojyoti Ghosh     GLEETableauList = link->next;
3832302aa1bSDebojyoti Ghosh     ierr = PetscFree5(t->A,t->B,t->U,t->V,t->c);  CHKERRQ(ierr);
3842302aa1bSDebojyoti Ghosh     ierr = PetscFree2(t->S,t->F);                 CHKERRQ(ierr);
3852302aa1bSDebojyoti Ghosh     ierr = PetscFree (t->Fembed);                 CHKERRQ(ierr);
386fd96ebc2SDebojyoti Ghosh     ierr = PetscFree (t->Ferror);                 CHKERRQ(ierr);
38757df6a1bSDebojyoti Ghosh     ierr = PetscFree (t->Serror);                 CHKERRQ(ierr);
3882302aa1bSDebojyoti Ghosh     ierr = PetscFree (t->binterp);                CHKERRQ(ierr);
3892302aa1bSDebojyoti Ghosh     ierr = PetscFree (t->name);                   CHKERRQ(ierr);
3902302aa1bSDebojyoti Ghosh     ierr = PetscFree (link);                      CHKERRQ(ierr);
3912302aa1bSDebojyoti Ghosh   }
3922302aa1bSDebojyoti Ghosh   TSGLEERegisterAllCalled = PETSC_FALSE;
3932302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
3942302aa1bSDebojyoti Ghosh }
3952302aa1bSDebojyoti Ghosh 
3962302aa1bSDebojyoti Ghosh #undef __FUNCT__
3972302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEEInitializePackage"
3982302aa1bSDebojyoti Ghosh /*@C
3992302aa1bSDebojyoti Ghosh   TSGLEEInitializePackage - This function initializes everything in the TSGLEE package. It is called
4002302aa1bSDebojyoti Ghosh   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_GLEE()
4012302aa1bSDebojyoti Ghosh   when using static libraries.
4022302aa1bSDebojyoti Ghosh 
4032302aa1bSDebojyoti Ghosh   Level: developer
4042302aa1bSDebojyoti Ghosh 
4052302aa1bSDebojyoti Ghosh .keywords: TS, TSGLEE, initialize, package
4062302aa1bSDebojyoti Ghosh .seealso: PetscInitialize()
4072302aa1bSDebojyoti Ghosh @*/
4082302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEEInitializePackage(void)
4092302aa1bSDebojyoti Ghosh {
4102302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
4112302aa1bSDebojyoti Ghosh 
4122302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
4132302aa1bSDebojyoti Ghosh   if (TSGLEEPackageInitialized) PetscFunctionReturn(0);
4142302aa1bSDebojyoti Ghosh   TSGLEEPackageInitialized = PETSC_TRUE;
4152302aa1bSDebojyoti Ghosh   ierr = TSGLEERegisterAll();CHKERRQ(ierr);
4162302aa1bSDebojyoti Ghosh   ierr = PetscObjectComposedDataRegister(&explicit_stage_time_id);CHKERRQ(ierr);
4172302aa1bSDebojyoti Ghosh   ierr = PetscRegisterFinalize(TSGLEEFinalizePackage);CHKERRQ(ierr);
4182302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
4192302aa1bSDebojyoti Ghosh }
4202302aa1bSDebojyoti Ghosh 
4212302aa1bSDebojyoti Ghosh #undef __FUNCT__
4222302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEEFinalizePackage"
4232302aa1bSDebojyoti Ghosh /*@C
4242302aa1bSDebojyoti Ghosh   TSGLEEFinalizePackage - This function destroys everything in the TSGLEE package. It is
4252302aa1bSDebojyoti Ghosh   called from PetscFinalize().
4262302aa1bSDebojyoti Ghosh 
4272302aa1bSDebojyoti Ghosh   Level: developer
4282302aa1bSDebojyoti Ghosh 
4292302aa1bSDebojyoti Ghosh .keywords: Petsc, destroy, package
4302302aa1bSDebojyoti Ghosh .seealso: PetscFinalize()
4312302aa1bSDebojyoti Ghosh @*/
4322302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEEFinalizePackage(void)
4332302aa1bSDebojyoti Ghosh {
4342302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
4352302aa1bSDebojyoti Ghosh 
4362302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
4372302aa1bSDebojyoti Ghosh   TSGLEEPackageInitialized = PETSC_FALSE;
4382302aa1bSDebojyoti Ghosh   ierr = TSGLEERegisterDestroy();CHKERRQ(ierr);
4392302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
4402302aa1bSDebojyoti Ghosh }
4412302aa1bSDebojyoti Ghosh 
4422302aa1bSDebojyoti Ghosh #undef __FUNCT__
4432302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEERegister"
4442302aa1bSDebojyoti Ghosh /*@C
4452302aa1bSDebojyoti Ghosh    TSGLEERegister - register an GLEE scheme by providing the entries in the Butcher tableau
4462302aa1bSDebojyoti Ghosh 
4472302aa1bSDebojyoti Ghosh    Not Collective, but the same schemes should be registered on all processes on which they will be used
4482302aa1bSDebojyoti Ghosh 
4492302aa1bSDebojyoti Ghosh    Input Parameters:
4502302aa1bSDebojyoti Ghosh +  name   - identifier for method
4512302aa1bSDebojyoti Ghosh .  order  - order of method
4522302aa1bSDebojyoti Ghosh .  s - number of stages
4532302aa1bSDebojyoti Ghosh .  r - number of steps
4542302aa1bSDebojyoti Ghosh .  gamma - LTE ratio
4552302aa1bSDebojyoti Ghosh .  A - stage coefficients (dimension s*s, row-major)
4562302aa1bSDebojyoti Ghosh .  B - step completion coefficients (dimension r*s, row-major)
4572302aa1bSDebojyoti Ghosh .  U - method coefficients (dimension s*r, row-major)
4582302aa1bSDebojyoti Ghosh .  V - method coefficients (dimension r*r, row-major)
4592302aa1bSDebojyoti Ghosh .  S - starting coefficients
4602302aa1bSDebojyoti Ghosh .  F - finishing coefficients
4612302aa1bSDebojyoti Ghosh .  c - abscissa (dimension s; NULL to use row sums of A)
4622302aa1bSDebojyoti Ghosh .  Fembed - step completion coefficients for embedded method
463fd96ebc2SDebojyoti Ghosh .  Ferror - error computation coefficients
46457df6a1bSDebojyoti Ghosh .  Serror - error initialization coefficients
4652302aa1bSDebojyoti Ghosh .  pinterp - order of interpolation (0 if unavailable)
4662302aa1bSDebojyoti Ghosh -  binterp - array of interpolation coefficients (NULL if unavailable)
4672302aa1bSDebojyoti Ghosh 
4682302aa1bSDebojyoti Ghosh    Notes:
4692302aa1bSDebojyoti Ghosh    Several GLEE methods are provided, this function is only needed to create new methods.
4702302aa1bSDebojyoti Ghosh 
4712302aa1bSDebojyoti Ghosh    Level: advanced
4722302aa1bSDebojyoti Ghosh 
4732302aa1bSDebojyoti Ghosh .keywords: TS, register
4742302aa1bSDebojyoti Ghosh 
4752302aa1bSDebojyoti Ghosh .seealso: TSGLEE
4762302aa1bSDebojyoti Ghosh @*/
4772302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEERegister(TSGLEEType name,PetscInt order,PetscInt s, PetscInt r,
4782302aa1bSDebojyoti Ghosh                               PetscReal gamma,
4792302aa1bSDebojyoti Ghosh                               const PetscReal A[],const PetscReal B[],
4802302aa1bSDebojyoti Ghosh                               const PetscReal U[],const PetscReal V[],
4812302aa1bSDebojyoti Ghosh                               const PetscReal S[],const PetscReal F[],
482fd96ebc2SDebojyoti Ghosh                               const PetscReal c[],
483fd96ebc2SDebojyoti Ghosh                               const PetscReal Fembed[],const PetscReal Ferror[],
48457df6a1bSDebojyoti Ghosh                               const PetscReal Serror[],
4852302aa1bSDebojyoti Ghosh                               PetscInt pinterp, const PetscReal binterp[])
4862302aa1bSDebojyoti Ghosh {
4872302aa1bSDebojyoti Ghosh   PetscErrorCode    ierr;
4882302aa1bSDebojyoti Ghosh   GLEETableauLink   link;
4892302aa1bSDebojyoti Ghosh   GLEETableau       t;
4902302aa1bSDebojyoti Ghosh   PetscInt          i,j;
4912302aa1bSDebojyoti Ghosh 
4922302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
4932302aa1bSDebojyoti Ghosh   ierr     = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
4942302aa1bSDebojyoti Ghosh   ierr     = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
4952302aa1bSDebojyoti Ghosh   t        = &link->tab;
4962302aa1bSDebojyoti Ghosh   ierr     = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
4972302aa1bSDebojyoti Ghosh   t->order = order;
4982302aa1bSDebojyoti Ghosh   t->s     = s;
4992302aa1bSDebojyoti Ghosh   t->r     = r;
5002302aa1bSDebojyoti Ghosh   t->gamma = gamma;
5012302aa1bSDebojyoti Ghosh   ierr     = PetscMalloc5(s*s,&t->A,r*r,&t->V,s,&t->c,r*s,&t->B,s*r,&t->U); CHKERRQ(ierr);
5022302aa1bSDebojyoti Ghosh   ierr     = PetscMalloc2(r,&t->S,r,&t->F); CHKERRQ(ierr);
5032302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
5042302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->B,B,r*s*sizeof(B[0]));CHKERRQ(ierr);
5052302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->U,U,s*r*sizeof(U[0]));CHKERRQ(ierr);
5062302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->V,V,r*r*sizeof(V[0]));CHKERRQ(ierr);
5072302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->S,S,r  *sizeof(S[0]));CHKERRQ(ierr);
5082302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->F,F,r  *sizeof(F[0]));CHKERRQ(ierr);
5092302aa1bSDebojyoti Ghosh   if (c) {
5102302aa1bSDebojyoti Ghosh     ierr   = PetscMemcpy(t->c,c,s*sizeof(c[0])); CHKERRQ(ierr);
5112302aa1bSDebojyoti Ghosh   } else {
5122302aa1bSDebojyoti Ghosh     for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j];
5132302aa1bSDebojyoti Ghosh   }
5142302aa1bSDebojyoti Ghosh   ierr = PetscMalloc1(r,&t->Fembed);CHKERRQ(ierr);
515fd96ebc2SDebojyoti Ghosh   ierr = PetscMalloc1(r,&t->Ferror);CHKERRQ(ierr);
51657df6a1bSDebojyoti Ghosh   ierr = PetscMalloc1(r,&t->Serror);CHKERRQ(ierr);
5172302aa1bSDebojyoti Ghosh   ierr = PetscMemcpy(t->Fembed,Fembed,r*sizeof(Fembed[0]));CHKERRQ(ierr);
518fd96ebc2SDebojyoti Ghosh   ierr = PetscMemcpy(t->Ferror,Ferror,r*sizeof(Ferror[0]));CHKERRQ(ierr);
51957df6a1bSDebojyoti Ghosh   ierr = PetscMemcpy(t->Serror,Serror,r*sizeof(Serror[0]));CHKERRQ(ierr);
5202302aa1bSDebojyoti Ghosh   t->pinterp = pinterp;
5212302aa1bSDebojyoti Ghosh   ierr       = PetscMalloc1(s*pinterp,&t->binterp);CHKERRQ(ierr);
5222302aa1bSDebojyoti Ghosh   ierr       = PetscMemcpy(t->binterp,binterp,s*pinterp*sizeof(binterp[0]));CHKERRQ(ierr);
5232302aa1bSDebojyoti Ghosh 
5242302aa1bSDebojyoti Ghosh   link->next      = GLEETableauList;
5252302aa1bSDebojyoti Ghosh   GLEETableauList = link;
5262302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
5272302aa1bSDebojyoti Ghosh }
5282302aa1bSDebojyoti Ghosh 
5292302aa1bSDebojyoti Ghosh #undef __FUNCT__
5302302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSEvaluateStep_GLEE"
5312302aa1bSDebojyoti Ghosh static PetscErrorCode TSEvaluateStep_GLEE(TS ts,PetscInt order,Vec X,PetscBool *done)
5322302aa1bSDebojyoti Ghosh {
5332302aa1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*) ts->data;
5342302aa1bSDebojyoti Ghosh   GLEETableau     tab = glee->tableau;
5352302aa1bSDebojyoti Ghosh   PetscReal       h, *B = tab->B, *V = tab->V,
536fd96ebc2SDebojyoti Ghosh                   *F = tab->F,
537fd96ebc2SDebojyoti Ghosh                   *Fembed = tab->Fembed;
5382302aa1bSDebojyoti Ghosh   PetscInt        s = tab->s, r = tab->r, i, j;
5392302aa1bSDebojyoti Ghosh   Vec             *Y = glee->Y, *YdotStage = glee->YdotStage;
540ec0e3ee2SEmil Constantinescu   PetscScalar     *ws = glee->swork, *wr = glee->rwork;
5412302aa1bSDebojyoti Ghosh   PetscErrorCode  ierr;
5422302aa1bSDebojyoti Ghosh 
5432302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
5442302aa1bSDebojyoti Ghosh 
5452302aa1bSDebojyoti Ghosh   switch (glee->status) {
5462302aa1bSDebojyoti Ghosh     case TS_STEP_INCOMPLETE:
5472302aa1bSDebojyoti Ghosh     case TS_STEP_PENDING:
5482302aa1bSDebojyoti Ghosh       h = ts->time_step; break;
5492302aa1bSDebojyoti Ghosh     case TS_STEP_COMPLETE:
550ed450d42SEmil Constantinescu       h = ts->ptime - ts->ptime_prev; break;
5512302aa1bSDebojyoti Ghosh     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
5522302aa1bSDebojyoti Ghosh   }
5532302aa1bSDebojyoti Ghosh 
5542302aa1bSDebojyoti Ghosh   if (order == tab->order) {
5552302aa1bSDebojyoti Ghosh 
5562302aa1bSDebojyoti Ghosh     /* Note: Irrespective of whether status is TS_STEP_INCOMPLETE
557fd96ebc2SDebojyoti Ghosh              or TS_STEP_COMPLETE, glee->X has the solution at the
5582302aa1bSDebojyoti Ghosh              beginning of the time step. So no need to roll-back.
5592302aa1bSDebojyoti Ghosh     */
5602302aa1bSDebojyoti Ghosh     if (glee->status == TS_STEP_INCOMPLETE) {
5612302aa1bSDebojyoti Ghosh       for (i=0; i<r; i++) {
5622302aa1bSDebojyoti Ghosh         ierr = VecZeroEntries(Y[i]); CHKERRQ(ierr);
563ec0e3ee2SEmil Constantinescu         for (j=0; j<r; j++) wr[j] = V[i*r+j];
564ec0e3ee2SEmil Constantinescu         ierr = VecMAXPY(Y[i],r,wr,glee->X); CHKERRQ(ierr);
565ec0e3ee2SEmil Constantinescu         for (j=0; j<s; j++) ws[j] = h*B[i*s+j];
566ec0e3ee2SEmil Constantinescu         ierr = VecMAXPY(Y[i],s,ws,YdotStage); CHKERRQ(ierr);
5672302aa1bSDebojyoti Ghosh       }
5682302aa1bSDebojyoti Ghosh       ierr = VecZeroEntries(X);CHKERRQ(ierr);
569ec0e3ee2SEmil Constantinescu       for (j=0; j<r; j++) wr[j] = F[j];
570ec0e3ee2SEmil Constantinescu       ierr = VecMAXPY(X,r,wr,Y);CHKERRQ(ierr);
5712302aa1bSDebojyoti Ghosh     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
5722302aa1bSDebojyoti Ghosh     PetscFunctionReturn(0);
5732302aa1bSDebojyoti Ghosh 
5742302aa1bSDebojyoti Ghosh   } else if (order == tab->order-1) {
5752302aa1bSDebojyoti Ghosh 
5762302aa1bSDebojyoti Ghosh     /* Complete with the embedded method (Fembed) */
5772302aa1bSDebojyoti Ghosh     for (i=0; i<r; i++) {
5782302aa1bSDebojyoti Ghosh       ierr = VecZeroEntries(Y[i]); CHKERRQ(ierr);
579ec0e3ee2SEmil Constantinescu       for (j=0; j<r; j++) wr[j] = V[i*r+j];
580ec0e3ee2SEmil Constantinescu       ierr = VecMAXPY(Y[i],r,wr,glee->X); CHKERRQ(ierr);
581ec0e3ee2SEmil Constantinescu       for (j=0; j<s; j++) ws[j] = h*B[i*s+j];
582ec0e3ee2SEmil Constantinescu       ierr = VecMAXPY(Y[i],s,ws,YdotStage); CHKERRQ(ierr);
5832302aa1bSDebojyoti Ghosh     }
5842302aa1bSDebojyoti Ghosh     ierr = VecZeroEntries(X);CHKERRQ(ierr);
585ec0e3ee2SEmil Constantinescu     for (j=0; j<r; j++) wr[j] = Fembed[j];
586ec0e3ee2SEmil Constantinescu     ierr = VecMAXPY(X,r,wr,Y);CHKERRQ(ierr);
587fd96ebc2SDebojyoti Ghosh 
5882302aa1bSDebojyoti Ghosh     if (done) *done = PETSC_TRUE;
5892302aa1bSDebojyoti Ghosh     PetscFunctionReturn(0);
5902302aa1bSDebojyoti Ghosh   }
5912302aa1bSDebojyoti Ghosh   if (done) *done = PETSC_FALSE;
5922302aa1bSDebojyoti 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);
5932302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
5942302aa1bSDebojyoti Ghosh }
5952302aa1bSDebojyoti Ghosh 
5962302aa1bSDebojyoti Ghosh #undef __FUNCT__
5972302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSStep_GLEE"
5982302aa1bSDebojyoti Ghosh static PetscErrorCode TSStep_GLEE(TS ts)
5992302aa1bSDebojyoti Ghosh {
6002302aa1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
6012302aa1bSDebojyoti Ghosh   GLEETableau     tab = glee->tableau;
6022302aa1bSDebojyoti Ghosh   const PetscInt  s = tab->s, r = tab->r;
6032302aa1bSDebojyoti Ghosh   PetscReal       *A = tab->A, *U = tab->U,
6042302aa1bSDebojyoti Ghosh                   *F = tab->F,
6052302aa1bSDebojyoti Ghosh                   *c = tab->c;
6062302aa1bSDebojyoti Ghosh   Vec             *Y = glee->Y, *X = glee->X,
6072302aa1bSDebojyoti Ghosh                   *YStage = glee->YStage,
6082302aa1bSDebojyoti Ghosh                   *YdotStage = glee->YdotStage,
6092302aa1bSDebojyoti Ghosh                   W = glee->W;
6102302aa1bSDebojyoti Ghosh   SNES            snes;
611ec0e3ee2SEmil Constantinescu   PetscScalar     *ws = glee->swork, *wr = glee->rwork;
6122302aa1bSDebojyoti Ghosh   TSAdapt         adapt;
6132302aa1bSDebojyoti Ghosh   PetscInt        i,j,reject,next_scheme,its,lits;
6142302aa1bSDebojyoti Ghosh   PetscReal       next_time_step;
6152302aa1bSDebojyoti Ghosh   PetscReal       t;
6162302aa1bSDebojyoti Ghosh   PetscBool       accept;
6172302aa1bSDebojyoti Ghosh   PetscErrorCode  ierr;
6182302aa1bSDebojyoti Ghosh 
6192302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
620642f3702SEmil Constantinescu   ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr);
621642f3702SEmil Constantinescu 
6222302aa1bSDebojyoti Ghosh   for (i=0; i<r; i++) { ierr = VecCopy(Y[i],X[i]); CHKERRQ(ierr); }
6232302aa1bSDebojyoti Ghosh 
6242302aa1bSDebojyoti Ghosh   ierr           = TSGetSNES(ts,&snes);CHKERRQ(ierr);
6252302aa1bSDebojyoti Ghosh   next_time_step = ts->time_step;
6262302aa1bSDebojyoti Ghosh   t              = ts->ptime;
6272302aa1bSDebojyoti Ghosh   accept         = PETSC_TRUE;
6282302aa1bSDebojyoti Ghosh   glee->status   = TS_STEP_INCOMPLETE;
6292302aa1bSDebojyoti Ghosh 
6302302aa1bSDebojyoti Ghosh   for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) {
6312302aa1bSDebojyoti Ghosh 
6322302aa1bSDebojyoti Ghosh     PetscReal h = ts->time_step;
6332302aa1bSDebojyoti Ghosh     ierr = TSPreStep(ts);CHKERRQ(ierr);
6342302aa1bSDebojyoti Ghosh 
6352302aa1bSDebojyoti Ghosh     for (i=0; i<s; i++) {
6362302aa1bSDebojyoti Ghosh 
6372302aa1bSDebojyoti Ghosh       glee->stage_time = t + h*c[i];
6382302aa1bSDebojyoti Ghosh       ierr = TSPreStage(ts,glee->stage_time); CHKERRQ(ierr);
6392302aa1bSDebojyoti Ghosh 
6402302aa1bSDebojyoti Ghosh       if (A[i*s+i] == 0) {  /* Explicit stage */
6412302aa1bSDebojyoti Ghosh         ierr = VecZeroEntries(YStage[i]);CHKERRQ(ierr);
642ec0e3ee2SEmil Constantinescu         for (j=0; j<r; j++) wr[j] = U[i*r+j];
643ec0e3ee2SEmil Constantinescu         ierr = VecMAXPY(YStage[i],r,wr,X);CHKERRQ(ierr);
644ec0e3ee2SEmil Constantinescu         for (j=0; j<i; j++) ws[j] = h*A[i*s+j];
645ec0e3ee2SEmil Constantinescu         ierr = VecMAXPY(YStage[i],i,ws,YdotStage);CHKERRQ(ierr);
6462302aa1bSDebojyoti Ghosh       } else {              /* Implicit stage */
6472302aa1bSDebojyoti Ghosh         glee->scoeff = 1.0/A[i*s+i];
6482302aa1bSDebojyoti Ghosh         /* compute right-hand-side */
6492302aa1bSDebojyoti Ghosh         ierr = VecZeroEntries(W);CHKERRQ(ierr);
650ec0e3ee2SEmil Constantinescu         for (j=0; j<r; j++) wr[j] = U[i*r+j];
651ec0e3ee2SEmil Constantinescu         ierr = VecMAXPY(W,r,wr,X);CHKERRQ(ierr);
652ec0e3ee2SEmil Constantinescu         for (j=0; j<i; j++) ws[j] = h*A[i*s+j];
653ec0e3ee2SEmil Constantinescu         ierr = VecMAXPY(W,i,ws,YdotStage);CHKERRQ(ierr);
6542302aa1bSDebojyoti Ghosh         ierr = VecScale(W,glee->scoeff/h);CHKERRQ(ierr);
6552302aa1bSDebojyoti Ghosh         /* set initial guess */
6562302aa1bSDebojyoti Ghosh         ierr = VecCopy(i>0 ? YStage[i-1] : ts->vec_sol,YStage[i]);CHKERRQ(ierr);
6572302aa1bSDebojyoti Ghosh         /* solve for this stage */
6582302aa1bSDebojyoti Ghosh         ierr = SNESSolve(snes,W,YStage[i]);CHKERRQ(ierr);
6592302aa1bSDebojyoti Ghosh         ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
6602302aa1bSDebojyoti Ghosh         ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
6612302aa1bSDebojyoti Ghosh         ts->snes_its += its; ts->ksp_its += lits;
6622302aa1bSDebojyoti Ghosh       }
6632302aa1bSDebojyoti Ghosh       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
664d6629dc5SEmil Constantinescu       ierr = TSAdaptCheckStage(adapt,ts,glee->stage_time,YStage[i],&accept);CHKERRQ(ierr);
6652302aa1bSDebojyoti Ghosh       if (!accept) goto reject_step;
6662302aa1bSDebojyoti Ghosh       ierr = TSPostStage(ts,glee->stage_time,i,YStage); CHKERRQ(ierr);
6672302aa1bSDebojyoti Ghosh       ierr = TSComputeRHSFunction(ts,t+h*c[i],YStage[i],YdotStage[i]);CHKERRQ(ierr);
6682302aa1bSDebojyoti Ghosh     }
6692302aa1bSDebojyoti Ghosh     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
6702302aa1bSDebojyoti Ghosh     glee->status = TS_STEP_PENDING;
6712302aa1bSDebojyoti Ghosh 
6722302aa1bSDebojyoti Ghosh     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
6732302aa1bSDebojyoti Ghosh     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
6742302aa1bSDebojyoti Ghosh     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
6752302aa1bSDebojyoti Ghosh     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
6762302aa1bSDebojyoti Ghosh     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
6772302aa1bSDebojyoti Ghosh     if (accept) {
6782302aa1bSDebojyoti Ghosh       /* ignore next_scheme for now */
6792302aa1bSDebojyoti Ghosh       ts->ptime     += ts->time_step;
6802302aa1bSDebojyoti Ghosh       ts->time_step  = next_time_step;
6812302aa1bSDebojyoti Ghosh       glee->status = TS_STEP_COMPLETE;
6820a01e1b2SEmil Constantinescu       /* compute and store the global error */
6830a01e1b2SEmil Constantinescu       /* Note: this is not needed if TSAdaptGLEE is not used */
6840a01e1b2SEmil Constantinescu       ierr = TSGetTimeError(ts,0,&(glee->yGErr));CHKERRQ(ierr);
6852302aa1bSDebojyoti Ghosh       ierr = PetscObjectComposedDataSetReal((PetscObject)ts->vec_sol,explicit_stage_time_id,ts->ptime);CHKERRQ(ierr);
6862302aa1bSDebojyoti Ghosh       break;
6872302aa1bSDebojyoti Ghosh     } else {                    /* Roll back the current step */
688ec0e3ee2SEmil Constantinescu       for (j=0; j<r; j++) wr[j] = F[j];
689ec0e3ee2SEmil Constantinescu       ierr = VecMAXPY(ts->vec_sol,r,wr,X);CHKERRQ(ierr);
6902302aa1bSDebojyoti Ghosh       ts->time_step = next_time_step;
6912302aa1bSDebojyoti Ghosh       glee->status  = TS_STEP_INCOMPLETE;
6922302aa1bSDebojyoti Ghosh     }
6932302aa1bSDebojyoti Ghosh reject_step: continue;
6942302aa1bSDebojyoti Ghosh   }
6952302aa1bSDebojyoti Ghosh   if (glee->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
6962302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
6972302aa1bSDebojyoti Ghosh }
6982302aa1bSDebojyoti Ghosh 
6992302aa1bSDebojyoti Ghosh #undef __FUNCT__
7002302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSInterpolate_GLEE"
7012302aa1bSDebojyoti Ghosh static PetscErrorCode TSInterpolate_GLEE(TS ts,PetscReal itime,Vec X)
7022302aa1bSDebojyoti Ghosh {
7032302aa1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
7042302aa1bSDebojyoti Ghosh   PetscInt        s=glee->tableau->s, pinterp=glee->tableau->pinterp,i,j;
7052302aa1bSDebojyoti Ghosh   PetscReal       h,tt,t;
7062302aa1bSDebojyoti Ghosh   PetscScalar     *b;
7072302aa1bSDebojyoti Ghosh   const PetscReal *B = glee->tableau->binterp;
7082302aa1bSDebojyoti Ghosh   PetscErrorCode  ierr;
7092302aa1bSDebojyoti Ghosh 
7102302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7112302aa1bSDebojyoti Ghosh   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSGLEE %s does not have an interpolation formula",glee->tableau->name);
7122302aa1bSDebojyoti Ghosh   switch (glee->status) {
7132302aa1bSDebojyoti Ghosh     case TS_STEP_INCOMPLETE:
7142302aa1bSDebojyoti Ghosh     case TS_STEP_PENDING:
7152302aa1bSDebojyoti Ghosh       h = ts->time_step;
7162302aa1bSDebojyoti Ghosh       t = (itime - ts->ptime)/h;
7172302aa1bSDebojyoti Ghosh       break;
7182302aa1bSDebojyoti Ghosh     case TS_STEP_COMPLETE:
719ed450d42SEmil Constantinescu       h = ts->ptime - ts->ptime_prev;
7202302aa1bSDebojyoti Ghosh       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
7212302aa1bSDebojyoti Ghosh       break;
7222302aa1bSDebojyoti Ghosh     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
7232302aa1bSDebojyoti Ghosh   }
7242302aa1bSDebojyoti Ghosh   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
7252302aa1bSDebojyoti Ghosh   for (i=0; i<s; i++) b[i] = 0;
7262302aa1bSDebojyoti Ghosh   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
7272302aa1bSDebojyoti Ghosh     for (i=0; i<s; i++) {
7282302aa1bSDebojyoti Ghosh       b[i]  += h * B[i*pinterp+j] * tt;
7292302aa1bSDebojyoti Ghosh     }
7302302aa1bSDebojyoti Ghosh   }
7312302aa1bSDebojyoti Ghosh   ierr = VecCopy(glee->YStage[0],X);CHKERRQ(ierr);
7322302aa1bSDebojyoti Ghosh   ierr = VecMAXPY(X,s,b,glee->YdotStage);CHKERRQ(ierr);
7332302aa1bSDebojyoti Ghosh   ierr = PetscFree(b);CHKERRQ(ierr);
7342302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
7352302aa1bSDebojyoti Ghosh }
7362302aa1bSDebojyoti Ghosh 
7372302aa1bSDebojyoti Ghosh /*------------------------------------------------------------*/
7382302aa1bSDebojyoti Ghosh #undef __FUNCT__
7392302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSReset_GLEE"
7402302aa1bSDebojyoti Ghosh static PetscErrorCode TSReset_GLEE(TS ts)
7412302aa1bSDebojyoti Ghosh {
7422302aa1bSDebojyoti Ghosh   TS_GLEE        *glee = (TS_GLEE*)ts->data;
7432302aa1bSDebojyoti Ghosh   PetscInt       s, r;
7442302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
7452302aa1bSDebojyoti Ghosh 
7462302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7472302aa1bSDebojyoti Ghosh   if (!glee->tableau) PetscFunctionReturn(0);
7482302aa1bSDebojyoti Ghosh   s    = glee->tableau->s;
7492302aa1bSDebojyoti Ghosh   r    = glee->tableau->r;
7502302aa1bSDebojyoti Ghosh   ierr = VecDestroyVecs(r,&glee->Y);CHKERRQ(ierr);
7512302aa1bSDebojyoti Ghosh   ierr = VecDestroyVecs(r,&glee->X);CHKERRQ(ierr);
7522302aa1bSDebojyoti Ghosh   ierr = VecDestroyVecs(s,&glee->YStage);CHKERRQ(ierr);
7532302aa1bSDebojyoti Ghosh   ierr = VecDestroyVecs(s,&glee->YdotStage);CHKERRQ(ierr);
7542302aa1bSDebojyoti Ghosh   ierr = VecDestroy(&glee->Ydot);CHKERRQ(ierr);
7550a01e1b2SEmil Constantinescu   ierr = VecDestroy(&glee->yGErr);CHKERRQ(ierr);
7562302aa1bSDebojyoti Ghosh   ierr = VecDestroy(&glee->W);CHKERRQ(ierr);
757ec0e3ee2SEmil Constantinescu   ierr = PetscFree2(glee->swork,glee->rwork);CHKERRQ(ierr);
7582302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
7592302aa1bSDebojyoti Ghosh }
7602302aa1bSDebojyoti Ghosh 
7612302aa1bSDebojyoti Ghosh #undef __FUNCT__
7622302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSDestroy_GLEE"
7632302aa1bSDebojyoti Ghosh static PetscErrorCode TSDestroy_GLEE(TS ts)
7642302aa1bSDebojyoti Ghosh {
7652302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
7662302aa1bSDebojyoti Ghosh 
7672302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7682302aa1bSDebojyoti Ghosh   ierr = TSReset_GLEE(ts);CHKERRQ(ierr);
7692302aa1bSDebojyoti Ghosh   ierr = PetscFree(ts->data);CHKERRQ(ierr);
7702302aa1bSDebojyoti Ghosh   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLEEGetType_C",NULL);CHKERRQ(ierr);
7712302aa1bSDebojyoti Ghosh   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLEESetType_C",NULL);CHKERRQ(ierr);
7722302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
7732302aa1bSDebojyoti Ghosh }
7742302aa1bSDebojyoti Ghosh 
7752302aa1bSDebojyoti Ghosh #undef __FUNCT__
7762302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEEGetVecs"
7772302aa1bSDebojyoti Ghosh static PetscErrorCode TSGLEEGetVecs(TS ts,DM dm,Vec *Ydot)
7782302aa1bSDebojyoti Ghosh {
7792302aa1bSDebojyoti Ghosh   TS_GLEE     *glee = (TS_GLEE*)ts->data;
7802302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
7812302aa1bSDebojyoti Ghosh 
7822302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7832302aa1bSDebojyoti Ghosh   if (Ydot) {
7842302aa1bSDebojyoti Ghosh     if (dm && dm != ts->dm) {
7852302aa1bSDebojyoti Ghosh       ierr = DMGetNamedGlobalVector(dm,"TSGLEE_Ydot",Ydot);CHKERRQ(ierr);
7862302aa1bSDebojyoti Ghosh     } else *Ydot = glee->Ydot;
7872302aa1bSDebojyoti Ghosh   }
7882302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
7892302aa1bSDebojyoti Ghosh }
7902302aa1bSDebojyoti Ghosh 
7912302aa1bSDebojyoti Ghosh 
7922302aa1bSDebojyoti Ghosh #undef __FUNCT__
7932302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEERestoreVecs"
7942302aa1bSDebojyoti Ghosh static PetscErrorCode TSGLEERestoreVecs(TS ts,DM dm,Vec *Ydot)
7952302aa1bSDebojyoti Ghosh {
7962302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
7972302aa1bSDebojyoti Ghosh 
7982302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7992302aa1bSDebojyoti Ghosh   if (Ydot) {
8002302aa1bSDebojyoti Ghosh     if (dm && dm != ts->dm) {
8012302aa1bSDebojyoti Ghosh       ierr = DMRestoreNamedGlobalVector(dm,"TSGLEE_Ydot",Ydot);CHKERRQ(ierr);
8022302aa1bSDebojyoti Ghosh     }
8032302aa1bSDebojyoti Ghosh   }
8042302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8052302aa1bSDebojyoti Ghosh }
8062302aa1bSDebojyoti Ghosh 
8072302aa1bSDebojyoti Ghosh /*
8082302aa1bSDebojyoti Ghosh   This defines the nonlinear equation that is to be solved with SNES
8092302aa1bSDebojyoti Ghosh */
8102302aa1bSDebojyoti Ghosh #undef __FUNCT__
8112302aa1bSDebojyoti Ghosh #define __FUNCT__ "SNESTSFormFunction_GLEE"
8122302aa1bSDebojyoti Ghosh static PetscErrorCode SNESTSFormFunction_GLEE(SNES snes,Vec X,Vec F,TS ts)
8132302aa1bSDebojyoti Ghosh {
8142302aa1bSDebojyoti Ghosh   TS_GLEE       *glee = (TS_GLEE*)ts->data;
8152302aa1bSDebojyoti Ghosh   DM             dm,dmsave;
8162302aa1bSDebojyoti Ghosh   Vec            Ydot;
8172302aa1bSDebojyoti Ghosh   PetscReal      shift = glee->scoeff / ts->time_step;
8182302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
8192302aa1bSDebojyoti Ghosh 
8202302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8212302aa1bSDebojyoti Ghosh   ierr   = SNESGetDM(snes,&dm);CHKERRQ(ierr);
8222302aa1bSDebojyoti Ghosh   ierr   = TSGLEEGetVecs(ts,dm,&Ydot);CHKERRQ(ierr);
8232302aa1bSDebojyoti Ghosh   /* Set Ydot = shift*X */
8242302aa1bSDebojyoti Ghosh   ierr   = VecCopy(X,Ydot);CHKERRQ(ierr);
8252302aa1bSDebojyoti Ghosh   ierr   = VecScale(Ydot,shift);CHKERRQ(ierr);
8262302aa1bSDebojyoti Ghosh   dmsave = ts->dm;
8272302aa1bSDebojyoti Ghosh   ts->dm = dm;
8282302aa1bSDebojyoti Ghosh 
8292302aa1bSDebojyoti Ghosh   ierr = TSComputeIFunction(ts,glee->stage_time,X,Ydot,F,PETSC_FALSE);CHKERRQ(ierr);
8302302aa1bSDebojyoti Ghosh 
8312302aa1bSDebojyoti Ghosh   ts->dm = dmsave;
8322302aa1bSDebojyoti Ghosh   ierr   = TSGLEERestoreVecs(ts,dm,&Ydot);CHKERRQ(ierr);
8332302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8342302aa1bSDebojyoti Ghosh }
8352302aa1bSDebojyoti Ghosh 
8362302aa1bSDebojyoti Ghosh #undef __FUNCT__
8372302aa1bSDebojyoti Ghosh #define __FUNCT__ "SNESTSFormJacobian_GLEE"
8382302aa1bSDebojyoti Ghosh static PetscErrorCode SNESTSFormJacobian_GLEE(SNES snes,Vec X,Mat A,Mat B,TS ts)
8392302aa1bSDebojyoti Ghosh {
8402302aa1bSDebojyoti Ghosh   TS_GLEE        *glee = (TS_GLEE*)ts->data;
8412302aa1bSDebojyoti Ghosh   DM             dm,dmsave;
8422302aa1bSDebojyoti Ghosh   Vec            Ydot;
8432302aa1bSDebojyoti Ghosh   PetscReal      shift = glee->scoeff / ts->time_step;
8442302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
8452302aa1bSDebojyoti Ghosh 
8462302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8472302aa1bSDebojyoti Ghosh   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
8482302aa1bSDebojyoti Ghosh   ierr = TSGLEEGetVecs(ts,dm,&Ydot);CHKERRQ(ierr);
8492302aa1bSDebojyoti Ghosh   /* glee->Ydot has already been computed in SNESTSFormFunction_GLEE (SNES guarantees this) */
8502302aa1bSDebojyoti Ghosh   dmsave = ts->dm;
8512302aa1bSDebojyoti Ghosh   ts->dm = dm;
8522302aa1bSDebojyoti Ghosh 
8532302aa1bSDebojyoti Ghosh   ierr = TSComputeIJacobian(ts,glee->stage_time,X,Ydot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
8542302aa1bSDebojyoti Ghosh 
8552302aa1bSDebojyoti Ghosh   ts->dm = dmsave;
8562302aa1bSDebojyoti Ghosh   ierr   = TSGLEERestoreVecs(ts,dm,&Ydot);CHKERRQ(ierr);
8572302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8582302aa1bSDebojyoti Ghosh }
8592302aa1bSDebojyoti Ghosh 
8602302aa1bSDebojyoti Ghosh #undef __FUNCT__
8612302aa1bSDebojyoti Ghosh #define __FUNCT__ "DMCoarsenHook_TSGLEE"
8622302aa1bSDebojyoti Ghosh static PetscErrorCode DMCoarsenHook_TSGLEE(DM fine,DM coarse,void *ctx)
8632302aa1bSDebojyoti Ghosh {
8642302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8652302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8662302aa1bSDebojyoti Ghosh }
8672302aa1bSDebojyoti Ghosh 
8682302aa1bSDebojyoti Ghosh #undef __FUNCT__
8692302aa1bSDebojyoti Ghosh #define __FUNCT__ "DMRestrictHook_TSGLEE"
8702302aa1bSDebojyoti Ghosh static PetscErrorCode DMRestrictHook_TSGLEE(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
8712302aa1bSDebojyoti Ghosh {
8722302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8732302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8742302aa1bSDebojyoti Ghosh }
8752302aa1bSDebojyoti Ghosh 
8762302aa1bSDebojyoti Ghosh 
8772302aa1bSDebojyoti Ghosh #undef __FUNCT__
8782302aa1bSDebojyoti Ghosh #define __FUNCT__ "DMSubDomainHook_TSGLEE"
8792302aa1bSDebojyoti Ghosh static PetscErrorCode DMSubDomainHook_TSGLEE(DM dm,DM subdm,void *ctx)
8802302aa1bSDebojyoti Ghosh {
8812302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8822302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8832302aa1bSDebojyoti Ghosh }
8842302aa1bSDebojyoti Ghosh 
8852302aa1bSDebojyoti Ghosh #undef __FUNCT__
8862302aa1bSDebojyoti Ghosh #define __FUNCT__ "DMSubDomainRestrictHook_TSGLEE"
8872302aa1bSDebojyoti Ghosh static PetscErrorCode DMSubDomainRestrictHook_TSGLEE(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
8882302aa1bSDebojyoti Ghosh {
8892302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8902302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8912302aa1bSDebojyoti Ghosh }
8922302aa1bSDebojyoti Ghosh 
8932302aa1bSDebojyoti Ghosh #undef __FUNCT__
8942302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSSetUp_GLEE"
8952302aa1bSDebojyoti Ghosh static PetscErrorCode TSSetUp_GLEE(TS ts)
8962302aa1bSDebojyoti Ghosh {
8972302aa1bSDebojyoti Ghosh   TS_GLEE        *glee = (TS_GLEE*)ts->data;
8982302aa1bSDebojyoti Ghosh   GLEETableau    tab;
8992302aa1bSDebojyoti Ghosh   PetscInt       s,r;
9002302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
9012302aa1bSDebojyoti Ghosh   DM             dm;
9022302aa1bSDebojyoti Ghosh 
9032302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
9042302aa1bSDebojyoti Ghosh   if (!glee->tableau) {
9052302aa1bSDebojyoti Ghosh     ierr = TSGLEESetType(ts,TSGLEEDefaultType);CHKERRQ(ierr);
9062302aa1bSDebojyoti Ghosh   }
9072302aa1bSDebojyoti Ghosh   tab  = glee->tableau;
9082302aa1bSDebojyoti Ghosh   s    = tab->s;
9092302aa1bSDebojyoti Ghosh   r    = tab->r;
9102302aa1bSDebojyoti Ghosh   ierr = VecDuplicateVecs(ts->vec_sol,r,&glee->Y);CHKERRQ(ierr);
9112302aa1bSDebojyoti Ghosh   ierr = VecDuplicateVecs(ts->vec_sol,r,&glee->X);CHKERRQ(ierr);
9122302aa1bSDebojyoti Ghosh   ierr = VecDuplicateVecs(ts->vec_sol,s,&glee->YStage);CHKERRQ(ierr);
9132302aa1bSDebojyoti Ghosh   ierr = VecDuplicateVecs(ts->vec_sol,s,&glee->YdotStage);CHKERRQ(ierr);
9142302aa1bSDebojyoti Ghosh   ierr = VecDuplicate(ts->vec_sol,&glee->Ydot);CHKERRQ(ierr);
9150a01e1b2SEmil Constantinescu   ierr = VecDuplicate(ts->vec_sol,&glee->yGErr);CHKERRQ(ierr);
916657c1e31SEmil Constantinescu   ierr = VecZeroEntries(glee->yGErr);CHKERRQ(ierr);
9172302aa1bSDebojyoti Ghosh   ierr = VecDuplicate(ts->vec_sol,&glee->W);CHKERRQ(ierr);
918ec0e3ee2SEmil Constantinescu   ierr = PetscMalloc2(s,&glee->swork,r,&glee->rwork); CHKERRQ(ierr);
9192302aa1bSDebojyoti Ghosh   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
9202302aa1bSDebojyoti Ghosh   if (dm) {
9212302aa1bSDebojyoti Ghosh     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSGLEE,DMRestrictHook_TSGLEE,ts);CHKERRQ(ierr);
9222302aa1bSDebojyoti Ghosh     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSGLEE,DMSubDomainRestrictHook_TSGLEE,ts);CHKERRQ(ierr);
9232302aa1bSDebojyoti Ghosh   }
9242302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
9252302aa1bSDebojyoti Ghosh }
9262302aa1bSDebojyoti Ghosh 
9272302aa1bSDebojyoti Ghosh #undef __FUNCT__
9282302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSStartingMethod_GLEE"
9292302aa1bSDebojyoti Ghosh PetscErrorCode TSStartingMethod_GLEE(TS ts)
9302302aa1bSDebojyoti Ghosh {
9312302aa1bSDebojyoti Ghosh   TS_GLEE        *glee = (TS_GLEE*)ts->data;
9326e2e232bSDebojyoti Ghosh   GLEETableau    tab  = glee->tableau;
9336e2e232bSDebojyoti Ghosh   PetscInt       r=tab->r,i;
9346e2e232bSDebojyoti Ghosh   PetscReal      *S=tab->S;
9352302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
9362302aa1bSDebojyoti Ghosh 
9372302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
9382302aa1bSDebojyoti Ghosh   for (i=0; i<r; i++) {
9392302aa1bSDebojyoti Ghosh     ierr = VecZeroEntries(glee->Y[i]);CHKERRQ(ierr);
9402302aa1bSDebojyoti Ghosh     ierr = VecAXPY(glee->Y[i],S[i],ts->vec_sol);CHKERRQ(ierr);
9412302aa1bSDebojyoti Ghosh   }
9422302aa1bSDebojyoti Ghosh 
9432302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
9442302aa1bSDebojyoti Ghosh }
9452302aa1bSDebojyoti Ghosh 
9462302aa1bSDebojyoti Ghosh 
9472302aa1bSDebojyoti Ghosh /*------------------------------------------------------------*/
9482302aa1bSDebojyoti Ghosh 
9492302aa1bSDebojyoti Ghosh #undef __FUNCT__
9502302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSSetFromOptions_GLEE"
9516b234c65SEmil Constantinescu static PetscErrorCode TSSetFromOptions_GLEE(PetscOptionItems *PetscOptionsObject,TS ts)
9522302aa1bSDebojyoti Ghosh {
9532302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
9542302aa1bSDebojyoti Ghosh   char           gleetype[256];
9552302aa1bSDebojyoti Ghosh 
9562302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
9572302aa1bSDebojyoti Ghosh   ierr = PetscOptionsHead(PetscOptionsObject,"GLEE ODE solver options");CHKERRQ(ierr);
9582302aa1bSDebojyoti Ghosh   {
9592302aa1bSDebojyoti Ghosh     GLEETableauLink link;
9602302aa1bSDebojyoti Ghosh     PetscInt        count,choice;
9612302aa1bSDebojyoti Ghosh     PetscBool       flg;
9622302aa1bSDebojyoti Ghosh     const char      **namelist;
9632302aa1bSDebojyoti Ghosh 
9642302aa1bSDebojyoti Ghosh     ierr = PetscStrncpy(gleetype,TSGLEEDefaultType,sizeof(gleetype));CHKERRQ(ierr);
9652302aa1bSDebojyoti Ghosh     for (link=GLEETableauList,count=0; link; link=link->next,count++) ;
9662302aa1bSDebojyoti Ghosh     ierr = PetscMalloc1(count,&namelist);CHKERRQ(ierr);
9672302aa1bSDebojyoti Ghosh     for (link=GLEETableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
9682302aa1bSDebojyoti Ghosh     ierr = PetscOptionsEList("-ts_glee_type","Family of GLEE method","TSGLEESetType",(const char*const*)namelist,count,gleetype,&choice,&flg);CHKERRQ(ierr);
9692302aa1bSDebojyoti Ghosh     ierr = TSGLEESetType(ts,flg ? namelist[choice] : gleetype);CHKERRQ(ierr);
9702302aa1bSDebojyoti Ghosh     ierr = PetscFree(namelist);CHKERRQ(ierr);
9712302aa1bSDebojyoti Ghosh   }
9722302aa1bSDebojyoti Ghosh   ierr = PetscOptionsTail();CHKERRQ(ierr);
9732302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
9742302aa1bSDebojyoti Ghosh }
9752302aa1bSDebojyoti Ghosh 
9762302aa1bSDebojyoti Ghosh #undef __FUNCT__
9772302aa1bSDebojyoti Ghosh #define __FUNCT__ "PetscFormatRealArray"
9782302aa1bSDebojyoti Ghosh static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
9792302aa1bSDebojyoti Ghosh {
9802302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
9812302aa1bSDebojyoti Ghosh   PetscInt       i;
9822302aa1bSDebojyoti Ghosh   size_t         left,count;
9832302aa1bSDebojyoti Ghosh   char           *p;
9842302aa1bSDebojyoti Ghosh 
9852302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
9862302aa1bSDebojyoti Ghosh   for (i=0,p=buf,left=len; i<n; i++) {
9872302aa1bSDebojyoti Ghosh     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
9882302aa1bSDebojyoti Ghosh     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
9892302aa1bSDebojyoti Ghosh     left -= count;
9902302aa1bSDebojyoti Ghosh     p    += count;
9912302aa1bSDebojyoti Ghosh     *p++  = ' ';
9922302aa1bSDebojyoti Ghosh   }
9932302aa1bSDebojyoti Ghosh   p[i ? 0 : -1] = 0;
9942302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
9952302aa1bSDebojyoti Ghosh }
9962302aa1bSDebojyoti Ghosh 
9972302aa1bSDebojyoti Ghosh #undef __FUNCT__
9982302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSView_GLEE"
9992302aa1bSDebojyoti Ghosh static PetscErrorCode TSView_GLEE(TS ts,PetscViewer viewer)
10002302aa1bSDebojyoti Ghosh {
10012302aa1bSDebojyoti Ghosh   TS_GLEE        *glee   = (TS_GLEE*)ts->data;
10022302aa1bSDebojyoti Ghosh   GLEETableau    tab  = glee->tableau;
10032302aa1bSDebojyoti Ghosh   PetscBool      iascii;
10042302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
10052302aa1bSDebojyoti Ghosh   TSAdapt        adapt;
10062302aa1bSDebojyoti Ghosh 
10072302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
10082302aa1bSDebojyoti Ghosh   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
10092302aa1bSDebojyoti Ghosh   if (iascii) {
10102302aa1bSDebojyoti Ghosh     TSGLEEType gleetype;
10112302aa1bSDebojyoti Ghosh     char   buf[512];
10122302aa1bSDebojyoti Ghosh     ierr = TSGLEEGetType(ts,&gleetype);CHKERRQ(ierr);
10132302aa1bSDebojyoti Ghosh     ierr = PetscViewerASCIIPrintf(viewer,"  GLEE %s\n",gleetype);CHKERRQ(ierr);
10142302aa1bSDebojyoti Ghosh     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
10152302aa1bSDebojyoti Ghosh     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa     c = %s\n",buf);CHKERRQ(ierr);
10162302aa1bSDebojyoti Ghosh     /* Note: print out r as well */
10172302aa1bSDebojyoti Ghosh   }
10182302aa1bSDebojyoti Ghosh   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
10192302aa1bSDebojyoti Ghosh   ierr = TSAdaptView(adapt,viewer);CHKERRQ(ierr);
10202302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
10212302aa1bSDebojyoti Ghosh }
10222302aa1bSDebojyoti Ghosh 
10232302aa1bSDebojyoti Ghosh #undef __FUNCT__
10242302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSLoad_GLEE"
10252302aa1bSDebojyoti Ghosh static PetscErrorCode TSLoad_GLEE(TS ts,PetscViewer viewer)
10262302aa1bSDebojyoti Ghosh {
10272302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
10282302aa1bSDebojyoti Ghosh   SNES           snes;
10292302aa1bSDebojyoti Ghosh   TSAdapt        tsadapt;
10302302aa1bSDebojyoti Ghosh 
10312302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
10322302aa1bSDebojyoti Ghosh   ierr = TSGetAdapt(ts,&tsadapt);CHKERRQ(ierr);
10332302aa1bSDebojyoti Ghosh   ierr = TSAdaptLoad(tsadapt,viewer);CHKERRQ(ierr);
10342302aa1bSDebojyoti Ghosh   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
10352302aa1bSDebojyoti Ghosh   ierr = SNESLoad(snes,viewer);CHKERRQ(ierr);
10362302aa1bSDebojyoti Ghosh   /* function and Jacobian context for SNES when used with TS is always ts object */
10372302aa1bSDebojyoti Ghosh   ierr = SNESSetFunction(snes,NULL,NULL,ts);CHKERRQ(ierr);
10382302aa1bSDebojyoti Ghosh   ierr = SNESSetJacobian(snes,NULL,NULL,NULL,ts);CHKERRQ(ierr);
10392302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
10402302aa1bSDebojyoti Ghosh }
10412302aa1bSDebojyoti Ghosh 
10422302aa1bSDebojyoti Ghosh #undef __FUNCT__
10432302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEESetType"
10442302aa1bSDebojyoti Ghosh /*@C
10452302aa1bSDebojyoti Ghosh   TSGLEESetType - Set the type of GLEE scheme
10462302aa1bSDebojyoti Ghosh 
10472302aa1bSDebojyoti Ghosh   Logically collective
10482302aa1bSDebojyoti Ghosh 
10492302aa1bSDebojyoti Ghosh   Input Parameter:
10502302aa1bSDebojyoti Ghosh +  ts - timestepping context
10512302aa1bSDebojyoti Ghosh -  gleetype - type of GLEE-scheme
10522302aa1bSDebojyoti Ghosh 
10532302aa1bSDebojyoti Ghosh   Level: intermediate
10542302aa1bSDebojyoti Ghosh 
10552302aa1bSDebojyoti Ghosh .seealso: TSGLEEGetType(), TSGLEE
10562302aa1bSDebojyoti Ghosh @*/
10572302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEESetType(TS ts,TSGLEEType gleetype)
10582302aa1bSDebojyoti Ghosh {
10592302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
10602302aa1bSDebojyoti Ghosh 
10612302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
10622302aa1bSDebojyoti Ghosh   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
10632302aa1bSDebojyoti Ghosh   ierr = PetscTryMethod(ts,"TSGLEESetType_C",(TS,TSGLEEType),(ts,gleetype));CHKERRQ(ierr);
10642302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
10652302aa1bSDebojyoti Ghosh }
10662302aa1bSDebojyoti Ghosh 
10672302aa1bSDebojyoti Ghosh #undef __FUNCT__
10682302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEEGetType"
10692302aa1bSDebojyoti Ghosh /*@C
10702302aa1bSDebojyoti Ghosh   TSGLEEGetType - Get the type of GLEE scheme
10712302aa1bSDebojyoti Ghosh 
10722302aa1bSDebojyoti Ghosh   Logically collective
10732302aa1bSDebojyoti Ghosh 
10742302aa1bSDebojyoti Ghosh   Input Parameter:
10752302aa1bSDebojyoti Ghosh .  ts - timestepping context
10762302aa1bSDebojyoti Ghosh 
10772302aa1bSDebojyoti Ghosh   Output Parameter:
10782302aa1bSDebojyoti Ghosh .  gleetype - type of GLEE-scheme
10792302aa1bSDebojyoti Ghosh 
10802302aa1bSDebojyoti Ghosh   Level: intermediate
10812302aa1bSDebojyoti Ghosh 
10822302aa1bSDebojyoti Ghosh .seealso: TSGLEESetType()
10832302aa1bSDebojyoti Ghosh @*/
10842302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEEGetType(TS ts,TSGLEEType *gleetype)
10852302aa1bSDebojyoti Ghosh {
10862302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
10872302aa1bSDebojyoti Ghosh 
10882302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
10892302aa1bSDebojyoti Ghosh   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
10902302aa1bSDebojyoti Ghosh   ierr = PetscUseMethod(ts,"TSGLEEGetType_C",(TS,TSGLEEType*),(ts,gleetype));CHKERRQ(ierr);
10912302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
10922302aa1bSDebojyoti Ghosh }
10932302aa1bSDebojyoti Ghosh 
10942302aa1bSDebojyoti Ghosh #undef __FUNCT__
10952302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEEGetType_GLEE"
10962302aa1bSDebojyoti Ghosh PetscErrorCode  TSGLEEGetType_GLEE(TS ts,TSGLEEType *gleetype)
10972302aa1bSDebojyoti Ghosh {
10982302aa1bSDebojyoti Ghosh   TS_GLEE     *glee = (TS_GLEE*)ts->data;
10992302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
11002302aa1bSDebojyoti Ghosh 
11012302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
11022302aa1bSDebojyoti Ghosh   if (!glee->tableau) {
11032302aa1bSDebojyoti Ghosh     ierr = TSGLEESetType(ts,TSGLEEDefaultType);CHKERRQ(ierr);
11042302aa1bSDebojyoti Ghosh   }
11052302aa1bSDebojyoti Ghosh   *gleetype = glee->tableau->name;
11062302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
11072302aa1bSDebojyoti Ghosh }
11082302aa1bSDebojyoti Ghosh #undef __FUNCT__
11092302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEESetType_GLEE"
11102302aa1bSDebojyoti Ghosh PetscErrorCode  TSGLEESetType_GLEE(TS ts,TSGLEEType gleetype)
11112302aa1bSDebojyoti Ghosh {
11122302aa1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
11132302aa1bSDebojyoti Ghosh   PetscErrorCode  ierr;
11142302aa1bSDebojyoti Ghosh   PetscBool       match;
11152302aa1bSDebojyoti Ghosh   GLEETableauLink link;
11162302aa1bSDebojyoti Ghosh 
11172302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
11182302aa1bSDebojyoti Ghosh   if (glee->tableau) {
11192302aa1bSDebojyoti Ghosh     ierr = PetscStrcmp(glee->tableau->name,gleetype,&match);CHKERRQ(ierr);
11202302aa1bSDebojyoti Ghosh     if (match) PetscFunctionReturn(0);
11212302aa1bSDebojyoti Ghosh   }
11222302aa1bSDebojyoti Ghosh   for (link = GLEETableauList; link; link=link->next) {
11232302aa1bSDebojyoti Ghosh     ierr = PetscStrcmp(link->tab.name,gleetype,&match);CHKERRQ(ierr);
11242302aa1bSDebojyoti Ghosh     if (match) {
11252302aa1bSDebojyoti Ghosh       ierr = TSReset_GLEE(ts);CHKERRQ(ierr);
11262302aa1bSDebojyoti Ghosh       glee->tableau = &link->tab;
11272302aa1bSDebojyoti Ghosh       PetscFunctionReturn(0);
11282302aa1bSDebojyoti Ghosh     }
11292302aa1bSDebojyoti Ghosh   }
11302302aa1bSDebojyoti Ghosh   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",gleetype);
11312302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
11322302aa1bSDebojyoti Ghosh }
11332302aa1bSDebojyoti Ghosh 
11342302aa1bSDebojyoti Ghosh #undef __FUNCT__
11352302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGetStages_GLEE"
11362302aa1bSDebojyoti Ghosh static PetscErrorCode  TSGetStages_GLEE(TS ts,PetscInt *ns,Vec **Y)
11372302aa1bSDebojyoti Ghosh {
11382302aa1bSDebojyoti Ghosh   TS_GLEE *glee = (TS_GLEE*)ts->data;
11392302aa1bSDebojyoti Ghosh 
11402302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
11412302aa1bSDebojyoti Ghosh   *ns = glee->tableau->s;
1142d5509560SEmil Constantinescu   if(Y) *Y  = glee->YStage;
11432302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
11442302aa1bSDebojyoti Ghosh }
11452302aa1bSDebojyoti Ghosh 
11462302aa1bSDebojyoti Ghosh #undef __FUNCT__
1147b2bf4f3aSDebojyoti Ghosh #define __FUNCT__ "TSGetSolutionComponents_GLEE"
1148b2bf4f3aSDebojyoti Ghosh PetscErrorCode TSGetSolutionComponents_GLEE(TS ts,PetscInt *n,Vec *Y)
11492302aa1bSDebojyoti Ghosh {
11502302aa1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
11512302aa1bSDebojyoti Ghosh   GLEETableau     tab   = glee->tableau;
11526e2e232bSDebojyoti Ghosh   PetscErrorCode  ierr;
11532302aa1bSDebojyoti Ghosh 
11542302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
11554cdd57e5SDebojyoti Ghosh   if (!Y) *n = tab->r;
11562302aa1bSDebojyoti Ghosh   else {
11574cdd57e5SDebojyoti Ghosh     if ((*n >= 0) && (*n < tab->r)) {
11584cdd57e5SDebojyoti Ghosh       ierr = VecCopy(glee->Y[*n],*Y); CHKERRQ(ierr);
11599657682dSDebojyoti Ghosh     } else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Second argument (%d) out of range[%d,%d].",*n,0,tab->r-1);
11602302aa1bSDebojyoti Ghosh   }
11612302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
11622302aa1bSDebojyoti Ghosh }
11632302aa1bSDebojyoti Ghosh 
11644cdd57e5SDebojyoti Ghosh #undef __FUNCT__
11654cdd57e5SDebojyoti Ghosh #define __FUNCT__ "TSGetAuxSolution_GLEE"
11664cdd57e5SDebojyoti Ghosh PetscErrorCode TSGetAuxSolution_GLEE(TS ts,Vec *X)
11674cdd57e5SDebojyoti Ghosh {
11684cdd57e5SDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
11694cdd57e5SDebojyoti Ghosh   GLEETableau     tab   = glee->tableau;
11704cdd57e5SDebojyoti Ghosh   PetscReal       *F    = tab->Fembed;
11714cdd57e5SDebojyoti Ghosh   PetscInt        r     = tab->r;
11724cdd57e5SDebojyoti Ghosh   Vec             *Y    = glee->Y;
1173ec0e3ee2SEmil Constantinescu   PetscScalar     *wr   = glee->rwork;
1174ec0e3ee2SEmil Constantinescu   PetscInt        i;
11754cdd57e5SDebojyoti Ghosh   PetscErrorCode  ierr;
11764cdd57e5SDebojyoti Ghosh 
11774cdd57e5SDebojyoti Ghosh   PetscFunctionBegin;
11784cdd57e5SDebojyoti Ghosh   ierr = VecZeroEntries(*X);CHKERRQ(ierr);
1179ec0e3ee2SEmil Constantinescu   for (i=0; i<r; i++) wr[i] = F[i];
1180ec0e3ee2SEmil Constantinescu   ierr = VecMAXPY((*X),r,wr,Y);CHKERRQ(ierr);
11814cdd57e5SDebojyoti Ghosh   PetscFunctionReturn(0);
11824cdd57e5SDebojyoti Ghosh }
11834cdd57e5SDebojyoti Ghosh 
11844cdd57e5SDebojyoti Ghosh #undef __FUNCT__
11854cdd57e5SDebojyoti Ghosh #define __FUNCT__ "TSGetTimeError_GLEE"
11860a01e1b2SEmil Constantinescu PetscErrorCode TSGetTimeError_GLEE(TS ts,PetscInt n,Vec *X)
11874cdd57e5SDebojyoti Ghosh {
11884cdd57e5SDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
11894cdd57e5SDebojyoti Ghosh   GLEETableau     tab   = glee->tableau;
11904cdd57e5SDebojyoti Ghosh   PetscReal       *F    = tab->Ferror;
11914cdd57e5SDebojyoti Ghosh   PetscInt        r     = tab->r;
11924cdd57e5SDebojyoti Ghosh   Vec             *Y    = glee->Y;
1193ec0e3ee2SEmil Constantinescu   PetscScalar     *wr   = glee->rwork;
1194ec0e3ee2SEmil Constantinescu   PetscInt        i;
11954cdd57e5SDebojyoti Ghosh   PetscErrorCode  ierr;
11964cdd57e5SDebojyoti Ghosh 
11974cdd57e5SDebojyoti Ghosh   PetscFunctionBegin;
11984cdd57e5SDebojyoti Ghosh   ierr = VecZeroEntries(*X);CHKERRQ(ierr);
1199657c1e31SEmil Constantinescu   if(n==0){
1200ec0e3ee2SEmil Constantinescu     for (i=0; i<r; i++) wr[i] = F[i];
1201ec0e3ee2SEmil Constantinescu     ierr = VecMAXPY((*X),r,wr,Y);CHKERRQ(ierr);
12020a01e1b2SEmil Constantinescu   } else if(n==-1) {
1203657c1e31SEmil Constantinescu     *X=glee->yGErr;
12040a01e1b2SEmil Constantinescu   }
12054cdd57e5SDebojyoti Ghosh   PetscFunctionReturn(0);
12064cdd57e5SDebojyoti Ghosh }
12074cdd57e5SDebojyoti Ghosh 
120857df6a1bSDebojyoti Ghosh #undef __FUNCT__
120957df6a1bSDebojyoti Ghosh #define __FUNCT__ "TSSetTimeError_GLEE"
121057df6a1bSDebojyoti Ghosh PetscErrorCode TSSetTimeError_GLEE(TS ts,Vec X)
121157df6a1bSDebojyoti Ghosh {
121257df6a1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
121357df6a1bSDebojyoti Ghosh   GLEETableau     tab   = glee->tableau;
121457df6a1bSDebojyoti Ghosh   PetscReal       *S    = tab->Serror;
121557df6a1bSDebojyoti Ghosh   PetscInt        r     = tab->r,i;
121657df6a1bSDebojyoti Ghosh   Vec             *Y    = glee->Y;
121757df6a1bSDebojyoti Ghosh   PetscErrorCode  ierr;
121857df6a1bSDebojyoti Ghosh 
121957df6a1bSDebojyoti Ghosh   PetscFunctionBegin;
12209657682dSDebojyoti Ghosh   if (r != 2) SETERRQ2(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSSetTimeError_GLEE not supported for '%s' with r=%D.",tab->name,tab->r);
122157df6a1bSDebojyoti Ghosh   for (i=1; i<r; i++) {
122257df6a1bSDebojyoti Ghosh     ierr = VecCopy(ts->vec_sol,Y[i]); CHKERRQ(ierr);
122357df6a1bSDebojyoti Ghosh     ierr = VecAXPBY(Y[i],S[0],S[1],X); CHKERRQ(ierr);
1224642ca4e8SEmil Constantinescu     ierr = VecCopy(X,glee->yGErr); CHKERRQ(ierr);
122557df6a1bSDebojyoti Ghosh   }
122657df6a1bSDebojyoti Ghosh   PetscFunctionReturn(0);
122757df6a1bSDebojyoti Ghosh }
122857df6a1bSDebojyoti Ghosh 
12292302aa1bSDebojyoti Ghosh /* ------------------------------------------------------------ */
12302302aa1bSDebojyoti Ghosh /*MC
12312302aa1bSDebojyoti Ghosh       TSGLEE - ODE and DAE solver using General Linear with Error Estimation schemes
12322302aa1bSDebojyoti Ghosh 
12332302aa1bSDebojyoti Ghosh   The user should provide the right hand side of the equation
12342302aa1bSDebojyoti Ghosh   using TSSetRHSFunction().
12352302aa1bSDebojyoti Ghosh 
12362302aa1bSDebojyoti Ghosh   Notes:
12372302aa1bSDebojyoti Ghosh   The default is TSGLEE35, it can be changed with TSGLEESetType() or -ts_glee_type
12382302aa1bSDebojyoti Ghosh 
12392302aa1bSDebojyoti Ghosh   Level: beginner
12402302aa1bSDebojyoti Ghosh 
12412302aa1bSDebojyoti Ghosh .seealso:  TSCreate(), TS, TSSetType(), TSGLEESetType(), TSGLEEGetType(),
12422302aa1bSDebojyoti Ghosh            TSGLEE23, TTSGLEE24, TSGLEE35, TSGLEE25I, TSGLEEEXRK2A,
12432302aa1bSDebojyoti Ghosh            TSGLEERK32G1, TSGLEERK285EX, TSGLEEType, TSGLEERegister()
12442302aa1bSDebojyoti Ghosh 
12452302aa1bSDebojyoti Ghosh M*/
12462302aa1bSDebojyoti Ghosh #undef __FUNCT__
12472302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSCreate_GLEE"
12482302aa1bSDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSCreate_GLEE(TS ts)
12492302aa1bSDebojyoti Ghosh {
12502302aa1bSDebojyoti Ghosh   TS_GLEE         *th;
12512302aa1bSDebojyoti Ghosh   PetscErrorCode  ierr;
12522302aa1bSDebojyoti Ghosh 
12532302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
12542302aa1bSDebojyoti Ghosh #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
12552302aa1bSDebojyoti Ghosh   ierr = TSGLEEInitializePackage();CHKERRQ(ierr);
12562302aa1bSDebojyoti Ghosh #endif
12572302aa1bSDebojyoti Ghosh 
12582302aa1bSDebojyoti Ghosh   ts->ops->reset                  = TSReset_GLEE;
12592302aa1bSDebojyoti Ghosh   ts->ops->destroy                = TSDestroy_GLEE;
12602302aa1bSDebojyoti Ghosh   ts->ops->view                   = TSView_GLEE;
12612302aa1bSDebojyoti Ghosh   ts->ops->load                   = TSLoad_GLEE;
12622302aa1bSDebojyoti Ghosh   ts->ops->setup                  = TSSetUp_GLEE;
12632302aa1bSDebojyoti Ghosh   ts->ops->step                   = TSStep_GLEE;
12642302aa1bSDebojyoti Ghosh   ts->ops->interpolate            = TSInterpolate_GLEE;
12652302aa1bSDebojyoti Ghosh   ts->ops->evaluatestep           = TSEvaluateStep_GLEE;
12662302aa1bSDebojyoti Ghosh   ts->ops->setfromoptions         = TSSetFromOptions_GLEE;
12672302aa1bSDebojyoti Ghosh   ts->ops->getstages              = TSGetStages_GLEE;
12682302aa1bSDebojyoti Ghosh   ts->ops->snesfunction           = SNESTSFormFunction_GLEE;
12692302aa1bSDebojyoti Ghosh   ts->ops->snesjacobian           = SNESTSFormJacobian_GLEE;
1270b2bf4f3aSDebojyoti Ghosh   ts->ops->getsolutioncomponents  = TSGetSolutionComponents_GLEE;
12714cdd57e5SDebojyoti Ghosh   ts->ops->getauxsolution         = TSGetAuxSolution_GLEE;
12724cdd57e5SDebojyoti Ghosh   ts->ops->gettimeerror           = TSGetTimeError_GLEE;
127357df6a1bSDebojyoti Ghosh   ts->ops->settimeerror           = TSSetTimeError_GLEE;
12742302aa1bSDebojyoti Ghosh   ts->ops->startingmethod         = TSStartingMethod_GLEE;
12752302aa1bSDebojyoti Ghosh 
12762302aa1bSDebojyoti Ghosh   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
12772302aa1bSDebojyoti Ghosh   ts->data = (void*)th;
12782302aa1bSDebojyoti Ghosh 
12792302aa1bSDebojyoti Ghosh   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLEEGetType_C",TSGLEEGetType_GLEE);CHKERRQ(ierr);
12802302aa1bSDebojyoti Ghosh   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLEESetType_C",TSGLEESetType_GLEE);CHKERRQ(ierr);
12812302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
12822302aa1bSDebojyoti Ghosh }
1283