xref: /petsc/src/ts/impls/glee/glee.c (revision 9657682daf79c7b3d0ddbe4fbfb7ddac10208b65)
12302aa1bSDebojyoti Ghosh /*
22302aa1bSDebojyoti Ghosh   Code for time stepping with the General Linear with Error Estimation method
32302aa1bSDebojyoti Ghosh 
42302aa1bSDebojyoti Ghosh   Notes:
52302aa1bSDebojyoti Ghosh   The general system is written as
62302aa1bSDebojyoti Ghosh 
72302aa1bSDebojyoti Ghosh   Udot = F(t,U)
82302aa1bSDebojyoti Ghosh 
92302aa1bSDebojyoti Ghosh */
102302aa1bSDebojyoti Ghosh #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
112302aa1bSDebojyoti Ghosh #include <petscdm.h>
122302aa1bSDebojyoti Ghosh 
132302aa1bSDebojyoti Ghosh static TSGLEEType   TSGLEEDefaultType = TSGLEE35;
142302aa1bSDebojyoti Ghosh static PetscBool    TSGLEERegisterAllCalled;
152302aa1bSDebojyoti Ghosh static PetscBool    TSGLEEPackageInitialized;
162302aa1bSDebojyoti Ghosh static PetscInt     explicit_stage_time_id;
172302aa1bSDebojyoti Ghosh 
182302aa1bSDebojyoti Ghosh typedef struct _GLEETableau *GLEETableau;
192302aa1bSDebojyoti Ghosh struct _GLEETableau {
202302aa1bSDebojyoti Ghosh   char      *name;
212302aa1bSDebojyoti Ghosh   PetscInt   order;               /* Classical approximation order of the method i*/
222302aa1bSDebojyoti Ghosh   PetscInt   s;                   /* Number of stages */
232302aa1bSDebojyoti Ghosh   PetscInt   r;                   /* Number of steps */
242302aa1bSDebojyoti Ghosh   PetscReal  gamma;               /* LTE ratio */
252302aa1bSDebojyoti Ghosh   PetscReal *A,*B,*U,*V,*S,*F,*c; /* Tableau */
262302aa1bSDebojyoti Ghosh   PetscReal *Fembed;              /* Embedded final method coefficients */
27fd96ebc2SDebojyoti Ghosh   PetscReal *Ferror;              /* Coefficients for computing error   */
2857df6a1bSDebojyoti Ghosh   PetscReal *Serror;              /* Coefficients for initializing the error   */
292302aa1bSDebojyoti Ghosh   PetscInt   pinterp;             /* Interpolation order */
302302aa1bSDebojyoti Ghosh   PetscReal *binterp;             /* Interpolation coefficients */
312302aa1bSDebojyoti Ghosh   PetscReal  ccfl;                /* Placeholder for CFL coefficient relative to forward Euler  */
322302aa1bSDebojyoti Ghosh };
332302aa1bSDebojyoti Ghosh typedef struct _GLEETableauLink *GLEETableauLink;
342302aa1bSDebojyoti Ghosh struct _GLEETableauLink {
352302aa1bSDebojyoti Ghosh   struct _GLEETableau tab;
362302aa1bSDebojyoti Ghosh   GLEETableauLink     next;
372302aa1bSDebojyoti Ghosh };
382302aa1bSDebojyoti Ghosh static GLEETableauLink GLEETableauList;
392302aa1bSDebojyoti Ghosh 
402302aa1bSDebojyoti Ghosh typedef struct {
412302aa1bSDebojyoti Ghosh   GLEETableau  tableau;
422302aa1bSDebojyoti Ghosh   Vec          *Y;         /* Solution vector (along with auxiliary solution y~ or eps) */
432302aa1bSDebojyoti Ghosh   Vec          *X;         /* Temporary solution vector */
442302aa1bSDebojyoti Ghosh   Vec          *YStage;    /* Stage values */
452302aa1bSDebojyoti Ghosh   Vec          *YdotStage; /* Stage right hand side */
462302aa1bSDebojyoti Ghosh   Vec          W;          /* Right-hand-side for implicit stage solve */
472302aa1bSDebojyoti Ghosh   Vec          Ydot;       /* Work vector holding Ydot during residual evaluation */
482302aa1bSDebojyoti Ghosh   PetscScalar  *work;      /* Scalar work */
492302aa1bSDebojyoti Ghosh   PetscReal    scoeff;     /* shift = scoeff/dt */
502302aa1bSDebojyoti Ghosh   PetscReal    stage_time;
512302aa1bSDebojyoti Ghosh   TSStepStatus status;
522302aa1bSDebojyoti Ghosh } TS_GLEE;
532302aa1bSDebojyoti Ghosh 
542302aa1bSDebojyoti Ghosh /*MC
552302aa1bSDebojyoti Ghosh      TSGLEE23 - Second order three stage GLEE method
562302aa1bSDebojyoti Ghosh 
572302aa1bSDebojyoti Ghosh      This method has three stages.
582302aa1bSDebojyoti Ghosh      s = 3, r = 2
592302aa1bSDebojyoti Ghosh 
602302aa1bSDebojyoti Ghosh      Level: advanced
612302aa1bSDebojyoti Ghosh 
622302aa1bSDebojyoti Ghosh .seealso: TSGLEE
632302aa1bSDebojyoti Ghosh */
642302aa1bSDebojyoti Ghosh /*MC
652302aa1bSDebojyoti Ghosh      TSGLEE24 - Second order four stage GLEE method
662302aa1bSDebojyoti Ghosh 
672302aa1bSDebojyoti Ghosh      This method has four stages.
682302aa1bSDebojyoti Ghosh      s = 4, r = 2
692302aa1bSDebojyoti Ghosh 
702302aa1bSDebojyoti Ghosh      Level: advanced
712302aa1bSDebojyoti Ghosh 
722302aa1bSDebojyoti Ghosh .seealso: TSGLEE
732302aa1bSDebojyoti Ghosh */
742302aa1bSDebojyoti Ghosh /*MC
752302aa1bSDebojyoti Ghosh      TSGLEE25i - Second order five stage GLEE method
762302aa1bSDebojyoti Ghosh 
772302aa1bSDebojyoti Ghosh      This method has five stages.
782302aa1bSDebojyoti Ghosh      s = 5, r = 2
792302aa1bSDebojyoti Ghosh 
802302aa1bSDebojyoti Ghosh      Level: advanced
812302aa1bSDebojyoti Ghosh 
822302aa1bSDebojyoti Ghosh .seealso: TSGLEE
832302aa1bSDebojyoti Ghosh */
842302aa1bSDebojyoti Ghosh /*MC
852302aa1bSDebojyoti Ghosh      TSGLEE35  - Third order five stage GLEE method
862302aa1bSDebojyoti Ghosh 
872302aa1bSDebojyoti Ghosh      This method has five stages.
882302aa1bSDebojyoti Ghosh      s = 5, r = 2
892302aa1bSDebojyoti Ghosh 
902302aa1bSDebojyoti Ghosh      Level: advanced
912302aa1bSDebojyoti Ghosh 
922302aa1bSDebojyoti Ghosh .seealso: TSGLEE
932302aa1bSDebojyoti Ghosh */
942302aa1bSDebojyoti Ghosh /*MC
952302aa1bSDebojyoti Ghosh      TSGLEEEXRK2A  - Second order six stage GLEE method
962302aa1bSDebojyoti Ghosh 
972302aa1bSDebojyoti Ghosh      This method has six stages.
982302aa1bSDebojyoti Ghosh      s = 6, r = 2
992302aa1bSDebojyoti Ghosh 
1002302aa1bSDebojyoti Ghosh      Level: advanced
1012302aa1bSDebojyoti Ghosh 
1022302aa1bSDebojyoti Ghosh .seealso: TSGLEE
1032302aa1bSDebojyoti Ghosh */
1042302aa1bSDebojyoti Ghosh /*MC
1052302aa1bSDebojyoti Ghosh      TSGLEERK32G1  - Third order eight stage GLEE method
1062302aa1bSDebojyoti Ghosh 
1072302aa1bSDebojyoti Ghosh      This method has eight stages.
1082302aa1bSDebojyoti Ghosh      s = 8, r = 2
1092302aa1bSDebojyoti Ghosh 
1102302aa1bSDebojyoti Ghosh      Level: advanced
1112302aa1bSDebojyoti Ghosh 
1122302aa1bSDebojyoti Ghosh .seealso: TSGLEE
1132302aa1bSDebojyoti Ghosh */
1142302aa1bSDebojyoti Ghosh /*MC
1152302aa1bSDebojyoti Ghosh      TSGLEERK285EX  - Second order nine stage GLEE method
1162302aa1bSDebojyoti Ghosh 
1172302aa1bSDebojyoti Ghosh      This method has nine stages.
1182302aa1bSDebojyoti Ghosh      s = 9, r = 2
1192302aa1bSDebojyoti Ghosh 
1202302aa1bSDebojyoti Ghosh      Level: advanced
1212302aa1bSDebojyoti Ghosh 
1222302aa1bSDebojyoti Ghosh .seealso: TSGLEE
1232302aa1bSDebojyoti Ghosh */
1242302aa1bSDebojyoti Ghosh 
1252302aa1bSDebojyoti Ghosh #undef __FUNCT__
1262302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEERegisterAll"
1272302aa1bSDebojyoti Ghosh /*@C
1282302aa1bSDebojyoti Ghosh   TSGLEERegisterAll - Registers all of the General Linear with Error Estimation methods in TSGLEE
1292302aa1bSDebojyoti Ghosh 
1302302aa1bSDebojyoti Ghosh   Not Collective, but should be called by all processes which will need the schemes to be registered
1312302aa1bSDebojyoti Ghosh 
1322302aa1bSDebojyoti Ghosh   Level: advanced
1332302aa1bSDebojyoti Ghosh 
1342302aa1bSDebojyoti Ghosh .keywords: TS, TSGLEE, register, all
1352302aa1bSDebojyoti Ghosh 
1362302aa1bSDebojyoti Ghosh .seealso:  TSGLEERegisterDestroy()
1372302aa1bSDebojyoti Ghosh @*/
1382302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEERegisterAll(void)
1392302aa1bSDebojyoti Ghosh {
1402302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
1412302aa1bSDebojyoti Ghosh 
1422302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
1432302aa1bSDebojyoti Ghosh   if (TSGLEERegisterAllCalled) PetscFunctionReturn(0);
1442302aa1bSDebojyoti Ghosh   TSGLEERegisterAllCalled = PETSC_TRUE;
1452302aa1bSDebojyoti Ghosh 
1462302aa1bSDebojyoti Ghosh   {
1472302aa1bSDebojyoti Ghosh     /* y-eps form */
1482302aa1bSDebojyoti Ghosh     const PetscInt
149ab8c5912SEmil Constantinescu       p = 1,
150ab8c5912SEmil Constantinescu       s = 3,
151ab8c5912SEmil Constantinescu       r = 2;
152ab8c5912SEmil Constantinescu     const PetscReal
153ab8c5912SEmil Constantinescu       gamma     = 0.5,
154ab8c5912SEmil Constantinescu       A[3][3]   = {{1.0,0,0},{0,0.5,0},{0,0.5,0.5}},
155ab8c5912SEmil Constantinescu       B[2][3]   = {{1.0,0,0},{-2.0,1.0,1.0}},
156ab8c5912SEmil Constantinescu       U[3][2]   = {{1.0,0},{1.0,0.5},{1.0,0.5}},
157ab8c5912SEmil Constantinescu       V[2][2]   = {{1,0},{0,1}},
158ab8c5912SEmil Constantinescu       S[2]      = {1,0},
159ab8c5912SEmil Constantinescu       F[2]      = {1,0},
160fd96ebc2SDebojyoti Ghosh       Fembed[2] = {1,1-gamma},
16157df6a1bSDebojyoti Ghosh       Ferror[2] = {0,1},
16257df6a1bSDebojyoti Ghosh       Serror[2] = {1,0};
16357df6a1bSDebojyoti Ghosh     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);
164ab8c5912SEmil Constantinescu   }
165ab8c5912SEmil Constantinescu   {
166ab8c5912SEmil Constantinescu     /* y-eps form */
167ab8c5912SEmil Constantinescu     const PetscInt
1682302aa1bSDebojyoti Ghosh       p = 2,
1692302aa1bSDebojyoti Ghosh       s = 3,
1702302aa1bSDebojyoti Ghosh       r = 2;
1712302aa1bSDebojyoti Ghosh     const PetscReal
1722302aa1bSDebojyoti Ghosh       gamma     = 0,
1732302aa1bSDebojyoti Ghosh       A[3][3]   = {{0,0,0},{1,0,0},{0.25,0.25,0}},
1742302aa1bSDebojyoti 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}},
1752302aa1bSDebojyoti Ghosh       U[3][2]   = {{1,0},{1,10},{1,-1}},
1762302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
1772302aa1bSDebojyoti Ghosh       S[2]      = {1,0},
1782302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
179fd96ebc2SDebojyoti Ghosh       Fembed[2] = {1,1-gamma},
18057df6a1bSDebojyoti Ghosh       Ferror[2] = {0,1},
18157df6a1bSDebojyoti Ghosh       Serror[2] = {1,0};
18257df6a1bSDebojyoti Ghosh     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);
1832302aa1bSDebojyoti Ghosh   }
1842302aa1bSDebojyoti Ghosh   {
1852302aa1bSDebojyoti Ghosh     /* y-y~ form */
1862302aa1bSDebojyoti Ghosh     const PetscInt
1872302aa1bSDebojyoti Ghosh       p = 2,
1882302aa1bSDebojyoti Ghosh       s = 4,
1892302aa1bSDebojyoti Ghosh       r = 2;
1902302aa1bSDebojyoti Ghosh     const PetscReal
1912302aa1bSDebojyoti Ghosh       gamma     = 0,
1922302aa1bSDebojyoti 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}},
1932302aa1bSDebojyoti 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}},
1942302aa1bSDebojyoti Ghosh       U[4][2]   = {{0,1},{75.0/58.0,-17.0/58.0},{0,1},{0,1}},
1952302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
1962302aa1bSDebojyoti Ghosh       S[2]      = {1,1},
1972302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
198fd96ebc2SDebojyoti Ghosh       Fembed[2] = {0,1},
19957df6a1bSDebojyoti Ghosh       Ferror[2] = {-1.0/(1.0-gamma),1.0/(1.0-gamma)},
20057df6a1bSDebojyoti Ghosh       Serror[2] = {1.0-gamma,1.0};
20157df6a1bSDebojyoti Ghosh       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);
2022302aa1bSDebojyoti Ghosh   }
2032302aa1bSDebojyoti Ghosh   {
2042302aa1bSDebojyoti Ghosh     /* y-y~ form */
2052302aa1bSDebojyoti Ghosh     const PetscInt
2062302aa1bSDebojyoti Ghosh       p = 2,
2072302aa1bSDebojyoti Ghosh       s = 5,
2082302aa1bSDebojyoti Ghosh       r = 2;
2092302aa1bSDebojyoti Ghosh     const PetscReal
2102302aa1bSDebojyoti Ghosh       gamma     = 0,
2112302aa1bSDebojyoti Ghosh       A[5][5]   = {{0,0,0,0,0},
2122302aa1bSDebojyoti Ghosh                    {-0.94079244066783383269,0,0,0,0},
2132302aa1bSDebojyoti Ghosh                    {0.64228187778301907108,0.10915356933958500042,0,0,0},
2142302aa1bSDebojyoti Ghosh                    {-0.51764297742287450812,0.74414270351096040738,-0.71404164927824538121,0,0},
2152302aa1bSDebojyoti Ghosh                    {-0.44696561556825969206,-0.76768425657590196518,0.20111608138142987881,0.93828186737840469796,0}},
2162302aa1bSDebojyoti Ghosh       B[2][5]   = {{-0.029309178948150356153,-0.49671981884013874923,0.34275801517650053274,0.32941112623949194988,0.85385985637229662276},
2172302aa1bSDebojyoti Ghosh                    {0.78133219686062535272,0.074238691892675897635,0.57957363498384957966,-0.24638502829674959968,-0.18875949544040123033}},
2182302aa1bSDebojyoti Ghosh       U[5][2]   = {{0.16911424754448327735,0.83088575245551672265},
2192302aa1bSDebojyoti Ghosh                    {0.53638465733199574340,0.46361534266800425660},
2202302aa1bSDebojyoti Ghosh                    {0.39901579167169582526,0.60098420832830417474},
2212302aa1bSDebojyoti Ghosh                    {0.87689005530618575480,0.12310994469381424520},
2222302aa1bSDebojyoti Ghosh                    {0.99056100455550913009,0.0094389954444908699092}},
2232302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
2242302aa1bSDebojyoti Ghosh       S[2]      = {1,1},
2252302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
226fd96ebc2SDebojyoti Ghosh       Fembed[2] = {0,1},
22757df6a1bSDebojyoti Ghosh       Ferror[2] = {-1.0/(1.0-gamma),1.0/(1.0-gamma)},
22857df6a1bSDebojyoti Ghosh       Serror[2] = {1.0-gamma,1.0};
22957df6a1bSDebojyoti Ghosh     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);
2302302aa1bSDebojyoti Ghosh   }
2312302aa1bSDebojyoti Ghosh   {
2322302aa1bSDebojyoti Ghosh     /* y-y~ form */
2332302aa1bSDebojyoti Ghosh     const PetscInt
2342302aa1bSDebojyoti Ghosh       p = 3,
2352302aa1bSDebojyoti Ghosh       s = 5,
2362302aa1bSDebojyoti Ghosh       r = 2;
2372302aa1bSDebojyoti Ghosh     const PetscReal
2382302aa1bSDebojyoti Ghosh       gamma     = 0,
2392302aa1bSDebojyoti Ghosh       A[5][5]   = {{0,0,0,0,0},
2402302aa1bSDebojyoti Ghosh                    {- 2169604947363702313.0 /  24313474998937147335.0,0,0,0,0},
2412302aa1bSDebojyoti Ghosh                    {46526746497697123895.0 /  94116917485856474137.0,-10297879244026594958.0 /  49199457603717988219.0,0,0,0},
2422302aa1bSDebojyoti Ghosh                    {23364788935845982499.0 /  87425311444725389446.0,-79205144337496116638.0 / 148994349441340815519.0,40051189859317443782.0 /  36487615018004984309.0,0,0},
2432302aa1bSDebojyoti Ghosh                    {42089522664062539205.0 / 124911313006412840286.0,-15074384760342762939.0 / 137927286865289746282.0,-62274678522253371016.0 / 125918573676298591413.0,13755475729852471739.0 /  79257927066651693390.0}},
2442302aa1bSDebojyoti Ghosh       B[2][5]   = {{61546696837458703723.0 /  56982519523786160813.0,-55810892792806293355.0 / 206957624151308356511.0,24061048952676379087.0 / 158739347956038723465.0,3577972206874351339.0 /   7599733370677197135.0,-59449832954780563947.0 / 137360038685338563670.0},
2452302aa1bSDebojyoti Ghosh                    {- 9738262186984159168.0 /  99299082461487742983.0,-32797097931948613195.0 /  61521565616362163366.0,42895514606418420631.0 /  71714201188501437336.0,22608567633166065068.0 /  55371917805607957003.0,94655809487476459565.0 / 151517167160302729021.0}},
2462302aa1bSDebojyoti Ghosh       U[5][2]   = {{70820309139834661559.0 /  80863923579509469826.0,10043614439674808267.0 /  80863923579509469826.0},
2472302aa1bSDebojyoti Ghosh                    {161694774978034105510.0 / 106187653640211060371.0,-55507121337823045139.0 / 106187653640211060371.0},
2482302aa1bSDebojyoti Ghosh                    {78486094644566264568.0 /  88171030896733822981.0,9684936252167558413.0 /  88171030896733822981.0},
2492302aa1bSDebojyoti Ghosh                    {65394922146334854435.0 /  84570853840405479554.0,19175931694070625119.0 /  84570853840405479554.0},
2502302aa1bSDebojyoti Ghosh                    {8607282770183754108.0 / 108658046436496925911.0,100050763666313171803.0 / 108658046436496925911.0}},
2512302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
2522302aa1bSDebojyoti Ghosh       S[2]      = {1,1},
2532302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
254fd96ebc2SDebojyoti Ghosh       Fembed[2] = {0,1},
25557df6a1bSDebojyoti Ghosh       Ferror[2] = {-1.0/(1.0-gamma),1.0/(1.0-gamma)},
25657df6a1bSDebojyoti Ghosh       Serror[2] = {1.0-gamma,1.0};
25757df6a1bSDebojyoti Ghosh     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);
2582302aa1bSDebojyoti Ghosh   }
2592302aa1bSDebojyoti Ghosh   {
2602302aa1bSDebojyoti Ghosh     /* y-eps form */
2612302aa1bSDebojyoti Ghosh     const PetscInt
2622302aa1bSDebojyoti Ghosh       p = 2,
2632302aa1bSDebojyoti Ghosh       s = 6,
2642302aa1bSDebojyoti Ghosh       r = 2;
2652302aa1bSDebojyoti Ghosh     const PetscReal
2662302aa1bSDebojyoti Ghosh       gamma     = 0.25,
2672302aa1bSDebojyoti Ghosh       A[6][6]   = {{0,0,0,0,0,0},
2682302aa1bSDebojyoti Ghosh                    {1,0,0,0,0,0},
2692302aa1bSDebojyoti Ghosh                    {0,0,0,0,0,0},
2702302aa1bSDebojyoti Ghosh                    {0,0,0.5,0,0,0},
2712302aa1bSDebojyoti Ghosh                    {0,0,0.25,0.25,0,0},
2722302aa1bSDebojyoti Ghosh                    {0,0,0.25,0.25,0.5,0}},
2732302aa1bSDebojyoti Ghosh       B[2][6]   = {{0.5,0.5,0,0,0,0},
2742302aa1bSDebojyoti 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}},
2752302aa1bSDebojyoti Ghosh       U[6][2]   = {{1,0},{1,0},{1,0.75},{1,0.75},{1,0.75},{1,0.75}},
2762302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
2772302aa1bSDebojyoti Ghosh       S[2]      = {1,0},
2782302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
279fd96ebc2SDebojyoti Ghosh       Fembed[2] = {1,1-gamma},
28057df6a1bSDebojyoti Ghosh       Ferror[2] = {0,1},
28157df6a1bSDebojyoti Ghosh       Serror[2] = {1,0};
28257df6a1bSDebojyoti Ghosh     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);
2832302aa1bSDebojyoti Ghosh   }
2842302aa1bSDebojyoti Ghosh   {
2852302aa1bSDebojyoti Ghosh     /* y-eps form */
2862302aa1bSDebojyoti Ghosh     const PetscInt
2872302aa1bSDebojyoti Ghosh       p = 3,
2882302aa1bSDebojyoti Ghosh       s = 8,
2892302aa1bSDebojyoti Ghosh       r = 2;
2902302aa1bSDebojyoti Ghosh     const PetscReal
2912302aa1bSDebojyoti Ghosh       gamma     = 0,
2922302aa1bSDebojyoti Ghosh       A[8][8]   = {{0,0,0,0,0,0,0,0},
2932302aa1bSDebojyoti Ghosh                    {0.5,0,0,0,0,0,0,0},
2942302aa1bSDebojyoti Ghosh                    {-1,2,0,0,0,0,0,0},
2952302aa1bSDebojyoti Ghosh                    {1.0/6.0,2.0/3.0,1.0/6.0,0,0,0,0,0},
2962302aa1bSDebojyoti Ghosh                    {0,0,0,0,0,0,0,0},
2972302aa1bSDebojyoti Ghosh                    {-7.0/24.0,1.0/3.0,1.0/12.0,-1.0/8.0,0.5,0,0,0},
2982302aa1bSDebojyoti Ghosh                    {7.0/6.0,-4.0/3.0,-1.0/3.0,0.5,-1.0,2.0,0,0},
2992302aa1bSDebojyoti Ghosh                    {0,0,0,0,1.0/6.0,2.0/3.0,1.0/6.0,0}},
3002302aa1bSDebojyoti Ghosh       B[2][8]   = {{1.0/6.0,2.0/3.0,1.0/6.0,0,0,0,0,0},
3012302aa1bSDebojyoti 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}},
3022302aa1bSDebojyoti Ghosh       U[8][2]   = {{1,0},{1,0},{1,0},{1,0},{1,1},{1,1},{1,1},{1,1}},
3032302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
3042302aa1bSDebojyoti Ghosh       S[2]      = {1,0},
3052302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
306fd96ebc2SDebojyoti Ghosh       Fembed[2] = {1,1-gamma},
30757df6a1bSDebojyoti Ghosh       Ferror[2] = {0,1},
30857df6a1bSDebojyoti Ghosh       Serror[2] = {1,0};
30957df6a1bSDebojyoti Ghosh     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);
3102302aa1bSDebojyoti Ghosh   }
3112302aa1bSDebojyoti Ghosh   {
3122302aa1bSDebojyoti Ghosh     /* y-eps form */
3132302aa1bSDebojyoti Ghosh     const PetscInt
3142302aa1bSDebojyoti Ghosh       p = 2,
3152302aa1bSDebojyoti Ghosh       s = 9,
3162302aa1bSDebojyoti Ghosh       r = 2;
3172302aa1bSDebojyoti Ghosh     const PetscReal
3182302aa1bSDebojyoti Ghosh       gamma     = 0.25,
3192302aa1bSDebojyoti Ghosh       A[9][9]   = {{0,0,0,0,0,0,0,0,0},
3202302aa1bSDebojyoti Ghosh                    {0.585786437626904966,0,0,0,0,0,0,0,0},
3212302aa1bSDebojyoti Ghosh                    {0.149999999999999994,0.849999999999999978,0,0,0,0,0,0,0},
3222302aa1bSDebojyoti Ghosh                    {0,0,0,0,0,0,0,0,0},
3232302aa1bSDebojyoti Ghosh                    {0,0,0,0.292893218813452483,0,0,0,0,0},
3242302aa1bSDebojyoti Ghosh                    {0,0,0,0.0749999999999999972,0.424999999999999989,0,0,0,0},
3252302aa1bSDebojyoti Ghosh                    {0,0,0,0.176776695296636893,0.176776695296636893,0.146446609406726241,0,0,0},
3262302aa1bSDebojyoti Ghosh                    {0,0,0,0.176776695296636893,0.176776695296636893,0.146446609406726241,0.292893218813452483,0,0},
3272302aa1bSDebojyoti Ghosh                    {0,0,0,0.176776695296636893,0.176776695296636893,0.146446609406726241,0.0749999999999999972,0.424999999999999989,0}},
3282302aa1bSDebojyoti Ghosh       B[2][9]   = {{0.353553390593273786,0.353553390593273786,0.292893218813452483,0,0,0,0,0,0},
3292302aa1bSDebojyoti Ghosh                    {-0.471404520791031678,-0.471404520791031678,-0.390524291751269959,0.235702260395515839,0.235702260395515839,0.195262145875634979,0.235702260395515839,0.235702260395515839,0.195262145875634979}},
3302302aa1bSDebojyoti 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}},
3312302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
3322302aa1bSDebojyoti Ghosh       S[2]      = {1,0},
3332302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
334fd96ebc2SDebojyoti Ghosh       Fembed[2] = {1,1-gamma},
33557df6a1bSDebojyoti Ghosh       Ferror[2] = {0,1},
33657df6a1bSDebojyoti Ghosh       Serror[2] = {1,0};
33757df6a1bSDebojyoti Ghosh     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);
3382302aa1bSDebojyoti Ghosh   }
3392302aa1bSDebojyoti Ghosh 
3402302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
3412302aa1bSDebojyoti Ghosh }
3422302aa1bSDebojyoti Ghosh 
3432302aa1bSDebojyoti Ghosh #undef __FUNCT__
3442302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEERegisterDestroy"
3452302aa1bSDebojyoti Ghosh /*@C
3462302aa1bSDebojyoti Ghosh    TSGLEERegisterDestroy - Frees the list of schemes that were registered by TSGLEERegister().
3472302aa1bSDebojyoti Ghosh 
3482302aa1bSDebojyoti Ghosh    Not Collective
3492302aa1bSDebojyoti Ghosh 
3502302aa1bSDebojyoti Ghosh    Level: advanced
3512302aa1bSDebojyoti Ghosh 
3522302aa1bSDebojyoti Ghosh .keywords: TSGLEE, register, destroy
3532302aa1bSDebojyoti Ghosh .seealso: TSGLEERegister(), TSGLEERegisterAll()
3542302aa1bSDebojyoti Ghosh @*/
3552302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEERegisterDestroy(void)
3562302aa1bSDebojyoti Ghosh {
3572302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
3582302aa1bSDebojyoti Ghosh   GLEETableauLink link;
3592302aa1bSDebojyoti Ghosh 
3602302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
3612302aa1bSDebojyoti Ghosh   while ((link = GLEETableauList)) {
3622302aa1bSDebojyoti Ghosh     GLEETableau t = &link->tab;
3632302aa1bSDebojyoti Ghosh     GLEETableauList = link->next;
3642302aa1bSDebojyoti Ghosh     ierr = PetscFree5(t->A,t->B,t->U,t->V,t->c);  CHKERRQ(ierr);
3652302aa1bSDebojyoti Ghosh     ierr = PetscFree2(t->S,t->F);                 CHKERRQ(ierr);
3662302aa1bSDebojyoti Ghosh     ierr = PetscFree (t->Fembed);                 CHKERRQ(ierr);
367fd96ebc2SDebojyoti Ghosh     ierr = PetscFree (t->Ferror);                 CHKERRQ(ierr);
36857df6a1bSDebojyoti Ghosh     ierr = PetscFree (t->Serror);                 CHKERRQ(ierr);
3692302aa1bSDebojyoti Ghosh     ierr = PetscFree (t->binterp);                CHKERRQ(ierr);
3702302aa1bSDebojyoti Ghosh     ierr = PetscFree (t->name);                   CHKERRQ(ierr);
3712302aa1bSDebojyoti Ghosh     ierr = PetscFree (link);                      CHKERRQ(ierr);
3722302aa1bSDebojyoti Ghosh   }
3732302aa1bSDebojyoti Ghosh   TSGLEERegisterAllCalled = PETSC_FALSE;
3742302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
3752302aa1bSDebojyoti Ghosh }
3762302aa1bSDebojyoti Ghosh 
3772302aa1bSDebojyoti Ghosh #undef __FUNCT__
3782302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEEInitializePackage"
3792302aa1bSDebojyoti Ghosh /*@C
3802302aa1bSDebojyoti Ghosh   TSGLEEInitializePackage - This function initializes everything in the TSGLEE package. It is called
3812302aa1bSDebojyoti Ghosh   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_GLEE()
3822302aa1bSDebojyoti Ghosh   when using static libraries.
3832302aa1bSDebojyoti Ghosh 
3842302aa1bSDebojyoti Ghosh   Level: developer
3852302aa1bSDebojyoti Ghosh 
3862302aa1bSDebojyoti Ghosh .keywords: TS, TSGLEE, initialize, package
3872302aa1bSDebojyoti Ghosh .seealso: PetscInitialize()
3882302aa1bSDebojyoti Ghosh @*/
3892302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEEInitializePackage(void)
3902302aa1bSDebojyoti Ghosh {
3912302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
3922302aa1bSDebojyoti Ghosh 
3932302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
3942302aa1bSDebojyoti Ghosh   if (TSGLEEPackageInitialized) PetscFunctionReturn(0);
3952302aa1bSDebojyoti Ghosh   TSGLEEPackageInitialized = PETSC_TRUE;
3962302aa1bSDebojyoti Ghosh   ierr = TSGLEERegisterAll();CHKERRQ(ierr);
3972302aa1bSDebojyoti Ghosh   ierr = PetscObjectComposedDataRegister(&explicit_stage_time_id);CHKERRQ(ierr);
3982302aa1bSDebojyoti Ghosh   ierr = PetscRegisterFinalize(TSGLEEFinalizePackage);CHKERRQ(ierr);
3992302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
4002302aa1bSDebojyoti Ghosh }
4012302aa1bSDebojyoti Ghosh 
4022302aa1bSDebojyoti Ghosh #undef __FUNCT__
4032302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEEFinalizePackage"
4042302aa1bSDebojyoti Ghosh /*@C
4052302aa1bSDebojyoti Ghosh   TSGLEEFinalizePackage - This function destroys everything in the TSGLEE package. It is
4062302aa1bSDebojyoti Ghosh   called from PetscFinalize().
4072302aa1bSDebojyoti Ghosh 
4082302aa1bSDebojyoti Ghosh   Level: developer
4092302aa1bSDebojyoti Ghosh 
4102302aa1bSDebojyoti Ghosh .keywords: Petsc, destroy, package
4112302aa1bSDebojyoti Ghosh .seealso: PetscFinalize()
4122302aa1bSDebojyoti Ghosh @*/
4132302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEEFinalizePackage(void)
4142302aa1bSDebojyoti Ghosh {
4152302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
4162302aa1bSDebojyoti Ghosh 
4172302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
4182302aa1bSDebojyoti Ghosh   TSGLEEPackageInitialized = PETSC_FALSE;
4192302aa1bSDebojyoti Ghosh   ierr = TSGLEERegisterDestroy();CHKERRQ(ierr);
4202302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
4212302aa1bSDebojyoti Ghosh }
4222302aa1bSDebojyoti Ghosh 
4232302aa1bSDebojyoti Ghosh #undef __FUNCT__
4242302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEERegister"
4252302aa1bSDebojyoti Ghosh /*@C
4262302aa1bSDebojyoti Ghosh    TSGLEERegister - register an GLEE scheme by providing the entries in the Butcher tableau
4272302aa1bSDebojyoti Ghosh 
4282302aa1bSDebojyoti Ghosh    Not Collective, but the same schemes should be registered on all processes on which they will be used
4292302aa1bSDebojyoti Ghosh 
4302302aa1bSDebojyoti Ghosh    Input Parameters:
4312302aa1bSDebojyoti Ghosh +  name   - identifier for method
4322302aa1bSDebojyoti Ghosh .  order  - order of method
4332302aa1bSDebojyoti Ghosh .  s - number of stages
4342302aa1bSDebojyoti Ghosh .  r - number of steps
4352302aa1bSDebojyoti Ghosh .  gamma - LTE ratio
4362302aa1bSDebojyoti Ghosh .  A - stage coefficients (dimension s*s, row-major)
4372302aa1bSDebojyoti Ghosh .  B - step completion coefficients (dimension r*s, row-major)
4382302aa1bSDebojyoti Ghosh .  U - method coefficients (dimension s*r, row-major)
4392302aa1bSDebojyoti Ghosh .  V - method coefficients (dimension r*r, row-major)
4402302aa1bSDebojyoti Ghosh .  S - starting coefficients
4412302aa1bSDebojyoti Ghosh .  F - finishing coefficients
4422302aa1bSDebojyoti Ghosh .  c - abscissa (dimension s; NULL to use row sums of A)
4432302aa1bSDebojyoti Ghosh .  Fembed - step completion coefficients for embedded method
444fd96ebc2SDebojyoti Ghosh .  Ferror - error computation coefficients
44557df6a1bSDebojyoti Ghosh .  Serror - error initialization coefficients
4462302aa1bSDebojyoti Ghosh .  pinterp - order of interpolation (0 if unavailable)
4472302aa1bSDebojyoti Ghosh -  binterp - array of interpolation coefficients (NULL if unavailable)
4482302aa1bSDebojyoti Ghosh 
4492302aa1bSDebojyoti Ghosh    Notes:
4502302aa1bSDebojyoti Ghosh    Several GLEE methods are provided, this function is only needed to create new methods.
4512302aa1bSDebojyoti Ghosh 
4522302aa1bSDebojyoti Ghosh    Level: advanced
4532302aa1bSDebojyoti Ghosh 
4542302aa1bSDebojyoti Ghosh .keywords: TS, register
4552302aa1bSDebojyoti Ghosh 
4562302aa1bSDebojyoti Ghosh .seealso: TSGLEE
4572302aa1bSDebojyoti Ghosh @*/
4582302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEERegister(TSGLEEType name,PetscInt order,PetscInt s, PetscInt r,
4592302aa1bSDebojyoti Ghosh                               PetscReal gamma,
4602302aa1bSDebojyoti Ghosh                               const PetscReal A[],const PetscReal B[],
4612302aa1bSDebojyoti Ghosh                               const PetscReal U[],const PetscReal V[],
4622302aa1bSDebojyoti Ghosh                               const PetscReal S[],const PetscReal F[],
463fd96ebc2SDebojyoti Ghosh                               const PetscReal c[],
464fd96ebc2SDebojyoti Ghosh                               const PetscReal Fembed[],const PetscReal Ferror[],
46557df6a1bSDebojyoti Ghosh                               const PetscReal Serror[],
4662302aa1bSDebojyoti Ghosh                               PetscInt pinterp, const PetscReal binterp[])
4672302aa1bSDebojyoti Ghosh {
4682302aa1bSDebojyoti Ghosh   PetscErrorCode    ierr;
4692302aa1bSDebojyoti Ghosh   GLEETableauLink   link;
4702302aa1bSDebojyoti Ghosh   GLEETableau       t;
4712302aa1bSDebojyoti Ghosh   PetscInt          i,j;
4722302aa1bSDebojyoti Ghosh 
4732302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
4742302aa1bSDebojyoti Ghosh   ierr     = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
4752302aa1bSDebojyoti Ghosh   ierr     = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
4762302aa1bSDebojyoti Ghosh   t        = &link->tab;
4772302aa1bSDebojyoti Ghosh   ierr     = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
4782302aa1bSDebojyoti Ghosh   t->order = order;
4792302aa1bSDebojyoti Ghosh   t->s     = s;
4802302aa1bSDebojyoti Ghosh   t->r     = r;
4812302aa1bSDebojyoti Ghosh   t->gamma = gamma;
4822302aa1bSDebojyoti Ghosh   ierr     = PetscMalloc5(s*s,&t->A,r*r,&t->V,s,&t->c,r*s,&t->B,s*r,&t->U); CHKERRQ(ierr);
4832302aa1bSDebojyoti Ghosh   ierr     = PetscMalloc2(r,&t->S,r,&t->F); CHKERRQ(ierr);
4842302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
4852302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->B,B,r*s*sizeof(B[0]));CHKERRQ(ierr);
4862302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->U,U,s*r*sizeof(U[0]));CHKERRQ(ierr);
4872302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->V,V,r*r*sizeof(V[0]));CHKERRQ(ierr);
4882302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->S,S,r  *sizeof(S[0]));CHKERRQ(ierr);
4892302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->F,F,r  *sizeof(F[0]));CHKERRQ(ierr);
4902302aa1bSDebojyoti Ghosh   if (c) {
4912302aa1bSDebojyoti Ghosh     ierr   = PetscMemcpy(t->c,c,s*sizeof(c[0])); CHKERRQ(ierr);
4922302aa1bSDebojyoti Ghosh   } else {
4932302aa1bSDebojyoti Ghosh     for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j];
4942302aa1bSDebojyoti Ghosh   }
4952302aa1bSDebojyoti Ghosh   ierr = PetscMalloc1(r,&t->Fembed);CHKERRQ(ierr);
496fd96ebc2SDebojyoti Ghosh   ierr = PetscMalloc1(r,&t->Ferror);CHKERRQ(ierr);
49757df6a1bSDebojyoti Ghosh   ierr = PetscMalloc1(r,&t->Serror);CHKERRQ(ierr);
4982302aa1bSDebojyoti Ghosh   ierr = PetscMemcpy(t->Fembed,Fembed,r*sizeof(Fembed[0]));CHKERRQ(ierr);
499fd96ebc2SDebojyoti Ghosh   ierr = PetscMemcpy(t->Ferror,Ferror,r*sizeof(Ferror[0]));CHKERRQ(ierr);
50057df6a1bSDebojyoti Ghosh   ierr = PetscMemcpy(t->Serror,Serror,r*sizeof(Serror[0]));CHKERRQ(ierr);
5012302aa1bSDebojyoti Ghosh   t->pinterp = pinterp;
5022302aa1bSDebojyoti Ghosh   ierr       = PetscMalloc1(s*pinterp,&t->binterp);CHKERRQ(ierr);
5032302aa1bSDebojyoti Ghosh   ierr       = PetscMemcpy(t->binterp,binterp,s*pinterp*sizeof(binterp[0]));CHKERRQ(ierr);
5042302aa1bSDebojyoti Ghosh 
5052302aa1bSDebojyoti Ghosh   link->next      = GLEETableauList;
5062302aa1bSDebojyoti Ghosh   GLEETableauList = link;
5072302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
5082302aa1bSDebojyoti Ghosh }
5092302aa1bSDebojyoti Ghosh 
5102302aa1bSDebojyoti Ghosh #undef __FUNCT__
5112302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSEvaluateStep_GLEE"
5122302aa1bSDebojyoti Ghosh static PetscErrorCode TSEvaluateStep_GLEE(TS ts,PetscInt order,Vec X,PetscBool *done)
5132302aa1bSDebojyoti Ghosh {
5142302aa1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*) ts->data;
5152302aa1bSDebojyoti Ghosh   GLEETableau     tab = glee->tableau;
5162302aa1bSDebojyoti Ghosh   PetscReal       h, *B = tab->B, *V = tab->V,
517fd96ebc2SDebojyoti Ghosh                   *F = tab->F,
518fd96ebc2SDebojyoti Ghosh                   *Fembed = tab->Fembed;
5192302aa1bSDebojyoti Ghosh   PetscInt        s = tab->s, r = tab->r, i, j;
5202302aa1bSDebojyoti Ghosh   Vec             *Y = glee->Y, *YdotStage = glee->YdotStage;
5212302aa1bSDebojyoti Ghosh   PetscScalar     *w = glee->work;
5222302aa1bSDebojyoti Ghosh   PetscErrorCode  ierr;
5232302aa1bSDebojyoti Ghosh 
5242302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
5252302aa1bSDebojyoti Ghosh 
5262302aa1bSDebojyoti Ghosh   switch (glee->status) {
5272302aa1bSDebojyoti Ghosh     case TS_STEP_INCOMPLETE:
5282302aa1bSDebojyoti Ghosh     case TS_STEP_PENDING:
5292302aa1bSDebojyoti Ghosh       h = ts->time_step; break;
5302302aa1bSDebojyoti Ghosh     case TS_STEP_COMPLETE:
5312302aa1bSDebojyoti Ghosh       h = ts->time_step_prev; break;
5322302aa1bSDebojyoti Ghosh     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
5332302aa1bSDebojyoti Ghosh   }
5342302aa1bSDebojyoti Ghosh 
5352302aa1bSDebojyoti Ghosh   if (order == tab->order) {
5362302aa1bSDebojyoti Ghosh 
5372302aa1bSDebojyoti Ghosh     /* Note: Irrespective of whether status is TS_STEP_INCOMPLETE
538fd96ebc2SDebojyoti Ghosh              or TS_STEP_COMPLETE, glee->X has the solution at the
5392302aa1bSDebojyoti Ghosh              beginning of the time step. So no need to roll-back.
5402302aa1bSDebojyoti Ghosh     */
5412302aa1bSDebojyoti Ghosh     if (glee->status == TS_STEP_INCOMPLETE) {
5422302aa1bSDebojyoti Ghosh       for (i=0; i<r; i++) {
5432302aa1bSDebojyoti Ghosh         ierr = VecZeroEntries(Y[i]); CHKERRQ(ierr);
5442302aa1bSDebojyoti Ghosh         ierr = VecMAXPY(Y[i],r,(V+i*r),glee->X); CHKERRQ(ierr);
5452302aa1bSDebojyoti Ghosh         for (j=0; j<s; j++) w[j] = h*B[i*s+j];
5462302aa1bSDebojyoti Ghosh         ierr = VecMAXPY(Y[i],s,w,YdotStage); CHKERRQ(ierr);
5472302aa1bSDebojyoti Ghosh       }
5482302aa1bSDebojyoti Ghosh       ierr = VecZeroEntries(X);CHKERRQ(ierr);
5492302aa1bSDebojyoti Ghosh       ierr = VecMAXPY(X,r,F,Y);CHKERRQ(ierr);
5502302aa1bSDebojyoti Ghosh     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
5512302aa1bSDebojyoti Ghosh     PetscFunctionReturn(0);
5522302aa1bSDebojyoti Ghosh 
5532302aa1bSDebojyoti Ghosh   } else if (order == tab->order-1) {
5542302aa1bSDebojyoti Ghosh 
5552302aa1bSDebojyoti Ghosh     /* Complete with the embedded method (Fembed) */
5562302aa1bSDebojyoti Ghosh     for (i=0; i<r; i++) {
5572302aa1bSDebojyoti Ghosh       ierr = VecZeroEntries(Y[i]); CHKERRQ(ierr);
5582302aa1bSDebojyoti Ghosh       ierr = VecMAXPY(Y[i],r,(V+i*r),glee->X); CHKERRQ(ierr);
5592302aa1bSDebojyoti Ghosh       for (j=0; j<s; j++) w[j] = h*B[i*s+j];
5602302aa1bSDebojyoti Ghosh       ierr = VecMAXPY(Y[i],s,w,YdotStage); CHKERRQ(ierr);
5612302aa1bSDebojyoti Ghosh     }
5622302aa1bSDebojyoti Ghosh     ierr = VecZeroEntries(X);CHKERRQ(ierr);
5632302aa1bSDebojyoti Ghosh     ierr = VecMAXPY(X,r,Fembed,Y);CHKERRQ(ierr);
564fd96ebc2SDebojyoti Ghosh 
5652302aa1bSDebojyoti Ghosh     if (done) *done = PETSC_TRUE;
5662302aa1bSDebojyoti Ghosh     PetscFunctionReturn(0);
5672302aa1bSDebojyoti Ghosh   }
5682302aa1bSDebojyoti Ghosh   if (done) *done = PETSC_FALSE;
5692302aa1bSDebojyoti 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);
5702302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
5712302aa1bSDebojyoti Ghosh }
5722302aa1bSDebojyoti Ghosh 
5732302aa1bSDebojyoti Ghosh #undef __FUNCT__
5742302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSStep_GLEE"
5752302aa1bSDebojyoti Ghosh static PetscErrorCode TSStep_GLEE(TS ts)
5762302aa1bSDebojyoti Ghosh {
5772302aa1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
5782302aa1bSDebojyoti Ghosh   GLEETableau     tab = glee->tableau;
5792302aa1bSDebojyoti Ghosh   const PetscInt  s = tab->s, r = tab->r;
5802302aa1bSDebojyoti Ghosh   PetscReal       *A = tab->A, *U = tab->U,
5812302aa1bSDebojyoti Ghosh                   *F = tab->F,
5822302aa1bSDebojyoti Ghosh                   *c = tab->c;
5832302aa1bSDebojyoti Ghosh   Vec             *Y = glee->Y, *X = glee->X,
5842302aa1bSDebojyoti Ghosh                   *YStage = glee->YStage,
5852302aa1bSDebojyoti Ghosh                   *YdotStage = glee->YdotStage,
5862302aa1bSDebojyoti Ghosh                   W = glee->W;
5872302aa1bSDebojyoti Ghosh   SNES            snes;
5882302aa1bSDebojyoti Ghosh   PetscScalar     *w = glee->work;
5892302aa1bSDebojyoti Ghosh   TSAdapt         adapt;
5902302aa1bSDebojyoti Ghosh   PetscInt        i,j,reject,next_scheme,its,lits;
5912302aa1bSDebojyoti Ghosh   PetscReal       next_time_step;
5922302aa1bSDebojyoti Ghosh   PetscReal       t;
5932302aa1bSDebojyoti Ghosh   PetscBool       accept;
5942302aa1bSDebojyoti Ghosh   PetscErrorCode  ierr;
5952302aa1bSDebojyoti Ghosh 
5962302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
5972302aa1bSDebojyoti Ghosh   for (i=0; i<r; i++) { ierr = VecCopy(Y[i],X[i]); CHKERRQ(ierr); }
5982302aa1bSDebojyoti Ghosh 
5992302aa1bSDebojyoti Ghosh   ierr           = TSGetSNES(ts,&snes);CHKERRQ(ierr);
6002302aa1bSDebojyoti Ghosh   next_time_step = ts->time_step;
6012302aa1bSDebojyoti Ghosh   t              = ts->ptime;
6022302aa1bSDebojyoti Ghosh   accept         = PETSC_TRUE;
6032302aa1bSDebojyoti Ghosh   glee->status   = TS_STEP_INCOMPLETE;
6042302aa1bSDebojyoti Ghosh 
6052302aa1bSDebojyoti Ghosh   for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) {
6062302aa1bSDebojyoti Ghosh 
6072302aa1bSDebojyoti Ghosh     PetscReal h = ts->time_step;
6082302aa1bSDebojyoti Ghosh     ierr = TSPreStep(ts);CHKERRQ(ierr);
6092302aa1bSDebojyoti Ghosh 
6102302aa1bSDebojyoti Ghosh     for (i=0; i<s; i++) {
6112302aa1bSDebojyoti Ghosh 
6122302aa1bSDebojyoti Ghosh       glee->stage_time = t + h*c[i];
6132302aa1bSDebojyoti Ghosh       ierr = TSPreStage(ts,glee->stage_time); CHKERRQ(ierr);
6142302aa1bSDebojyoti Ghosh 
6152302aa1bSDebojyoti Ghosh       if (A[i*s+i] == 0) {  /* Explicit stage */
6162302aa1bSDebojyoti Ghosh         ierr = VecZeroEntries(YStage[i]);CHKERRQ(ierr);
6172302aa1bSDebojyoti Ghosh         ierr = VecMAXPY(YStage[i],r,(U+i*r),X);CHKERRQ(ierr);
6182302aa1bSDebojyoti Ghosh         for (j=0; j<i; j++) w[j] = h*A[i*s+j];
6192302aa1bSDebojyoti Ghosh         ierr = VecMAXPY(YStage[i],i,w,YdotStage);CHKERRQ(ierr);
6202302aa1bSDebojyoti Ghosh       } else {              /* Implicit stage */
6212302aa1bSDebojyoti Ghosh         glee->scoeff = 1.0/A[i*s+i];
6222302aa1bSDebojyoti Ghosh         /* compute right-hand-side */
6232302aa1bSDebojyoti Ghosh         ierr = VecZeroEntries(W);CHKERRQ(ierr);
6242302aa1bSDebojyoti Ghosh         ierr = VecMAXPY(W,r,(U+i*r),X);CHKERRQ(ierr);
6252302aa1bSDebojyoti Ghosh         for (j=0; j<i; j++) w[j] = h*A[i*s+j];
6262302aa1bSDebojyoti Ghosh         ierr = VecMAXPY(W,i,w,YdotStage);CHKERRQ(ierr);
6272302aa1bSDebojyoti Ghosh         ierr = VecScale(W,glee->scoeff/h);CHKERRQ(ierr);
6282302aa1bSDebojyoti Ghosh         /* set initial guess */
6292302aa1bSDebojyoti Ghosh         ierr = VecCopy(i>0 ? YStage[i-1] : ts->vec_sol,YStage[i]);CHKERRQ(ierr);
6302302aa1bSDebojyoti Ghosh         /* solve for this stage */
6312302aa1bSDebojyoti Ghosh         ierr = SNESSolve(snes,W,YStage[i]);CHKERRQ(ierr);
6322302aa1bSDebojyoti Ghosh         ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
6332302aa1bSDebojyoti Ghosh         ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
6342302aa1bSDebojyoti Ghosh         ts->snes_its += its; ts->ksp_its += lits;
6352302aa1bSDebojyoti Ghosh       }
6362302aa1bSDebojyoti Ghosh       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
6372302aa1bSDebojyoti Ghosh       ierr = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr);
6382302aa1bSDebojyoti Ghosh       if (!accept) goto reject_step;
6392302aa1bSDebojyoti Ghosh       ierr = TSPostStage(ts,glee->stage_time,i,YStage); CHKERRQ(ierr);
6402302aa1bSDebojyoti Ghosh       ierr = TSComputeRHSFunction(ts,t+h*c[i],YStage[i],YdotStage[i]);CHKERRQ(ierr);
6412302aa1bSDebojyoti Ghosh     }
6422302aa1bSDebojyoti Ghosh     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
6432302aa1bSDebojyoti Ghosh     glee->status = TS_STEP_PENDING;
6442302aa1bSDebojyoti Ghosh 
6452302aa1bSDebojyoti Ghosh     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
6462302aa1bSDebojyoti Ghosh     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
6472302aa1bSDebojyoti Ghosh     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
6482302aa1bSDebojyoti Ghosh     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
6492302aa1bSDebojyoti Ghosh     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
6502302aa1bSDebojyoti Ghosh     if (accept) {
6512302aa1bSDebojyoti Ghosh       /* ignore next_scheme for now */
6522302aa1bSDebojyoti Ghosh       ts->ptime     += ts->time_step;
6532302aa1bSDebojyoti Ghosh       ts->time_step  = next_time_step;
6542302aa1bSDebojyoti Ghosh       ts->steps++;
6552302aa1bSDebojyoti Ghosh       glee->status = TS_STEP_COMPLETE;
6562302aa1bSDebojyoti Ghosh       ierr = PetscObjectComposedDataSetReal((PetscObject)ts->vec_sol,explicit_stage_time_id,ts->ptime);CHKERRQ(ierr);
6572302aa1bSDebojyoti Ghosh       break;
6582302aa1bSDebojyoti Ghosh     } else {                    /* Roll back the current step */
6592302aa1bSDebojyoti Ghosh       ierr = VecMAXPY(ts->vec_sol,r,F,X);CHKERRQ(ierr);
6602302aa1bSDebojyoti Ghosh       ts->time_step = next_time_step;
6612302aa1bSDebojyoti Ghosh       glee->status  = TS_STEP_INCOMPLETE;
6622302aa1bSDebojyoti Ghosh     }
6632302aa1bSDebojyoti Ghosh reject_step: continue;
6642302aa1bSDebojyoti Ghosh   }
6652302aa1bSDebojyoti Ghosh   if (glee->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
6662302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
6672302aa1bSDebojyoti Ghosh }
6682302aa1bSDebojyoti Ghosh 
6692302aa1bSDebojyoti Ghosh #undef __FUNCT__
6702302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSInterpolate_GLEE"
6712302aa1bSDebojyoti Ghosh static PetscErrorCode TSInterpolate_GLEE(TS ts,PetscReal itime,Vec X)
6722302aa1bSDebojyoti Ghosh {
6732302aa1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
6742302aa1bSDebojyoti Ghosh   PetscInt        s=glee->tableau->s, pinterp=glee->tableau->pinterp,i,j;
6752302aa1bSDebojyoti Ghosh   PetscReal       h,tt,t;
6762302aa1bSDebojyoti Ghosh   PetscScalar     *b;
6772302aa1bSDebojyoti Ghosh   const PetscReal *B = glee->tableau->binterp;
6782302aa1bSDebojyoti Ghosh   PetscErrorCode  ierr;
6792302aa1bSDebojyoti Ghosh 
6802302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
6812302aa1bSDebojyoti Ghosh   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSGLEE %s does not have an interpolation formula",glee->tableau->name);
6822302aa1bSDebojyoti Ghosh   switch (glee->status) {
6832302aa1bSDebojyoti Ghosh     case TS_STEP_INCOMPLETE:
6842302aa1bSDebojyoti Ghosh     case TS_STEP_PENDING:
6852302aa1bSDebojyoti Ghosh       h = ts->time_step;
6862302aa1bSDebojyoti Ghosh       t = (itime - ts->ptime)/h;
6872302aa1bSDebojyoti Ghosh       break;
6882302aa1bSDebojyoti Ghosh     case TS_STEP_COMPLETE:
6892302aa1bSDebojyoti Ghosh       h = ts->time_step_prev;
6902302aa1bSDebojyoti Ghosh       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
6912302aa1bSDebojyoti Ghosh       break;
6922302aa1bSDebojyoti Ghosh     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
6932302aa1bSDebojyoti Ghosh   }
6942302aa1bSDebojyoti Ghosh   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
6952302aa1bSDebojyoti Ghosh   for (i=0; i<s; i++) b[i] = 0;
6962302aa1bSDebojyoti Ghosh   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
6972302aa1bSDebojyoti Ghosh     for (i=0; i<s; i++) {
6982302aa1bSDebojyoti Ghosh       b[i]  += h * B[i*pinterp+j] * tt;
6992302aa1bSDebojyoti Ghosh     }
7002302aa1bSDebojyoti Ghosh   }
7012302aa1bSDebojyoti Ghosh   ierr = VecCopy(glee->YStage[0],X);CHKERRQ(ierr);
7022302aa1bSDebojyoti Ghosh   ierr = VecMAXPY(X,s,b,glee->YdotStage);CHKERRQ(ierr);
7032302aa1bSDebojyoti Ghosh   ierr = PetscFree(b);CHKERRQ(ierr);
7042302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
7052302aa1bSDebojyoti Ghosh }
7062302aa1bSDebojyoti Ghosh 
7072302aa1bSDebojyoti Ghosh /*------------------------------------------------------------*/
7082302aa1bSDebojyoti Ghosh #undef __FUNCT__
7092302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSReset_GLEE"
7102302aa1bSDebojyoti Ghosh static PetscErrorCode TSReset_GLEE(TS ts)
7112302aa1bSDebojyoti Ghosh {
7122302aa1bSDebojyoti Ghosh   TS_GLEE        *glee = (TS_GLEE*)ts->data;
7132302aa1bSDebojyoti Ghosh   PetscInt       s, r;
7142302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
7152302aa1bSDebojyoti Ghosh 
7162302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7172302aa1bSDebojyoti Ghosh   if (!glee->tableau) PetscFunctionReturn(0);
7182302aa1bSDebojyoti Ghosh   s    = glee->tableau->s;
7192302aa1bSDebojyoti Ghosh   r    = glee->tableau->r;
7202302aa1bSDebojyoti Ghosh   ierr = VecDestroyVecs(r,&glee->Y);CHKERRQ(ierr);
7212302aa1bSDebojyoti Ghosh   ierr = VecDestroyVecs(r,&glee->X);CHKERRQ(ierr);
7222302aa1bSDebojyoti Ghosh   ierr = VecDestroyVecs(s,&glee->YStage);CHKERRQ(ierr);
7232302aa1bSDebojyoti Ghosh   ierr = VecDestroyVecs(s,&glee->YdotStage);CHKERRQ(ierr);
7242302aa1bSDebojyoti Ghosh   ierr = VecDestroy(&glee->Ydot);CHKERRQ(ierr);
7252302aa1bSDebojyoti Ghosh   ierr = VecDestroy(&glee->W);CHKERRQ(ierr);
7262302aa1bSDebojyoti Ghosh   ierr = PetscFree(glee->work);CHKERRQ(ierr);
7272302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
7282302aa1bSDebojyoti Ghosh }
7292302aa1bSDebojyoti Ghosh 
7302302aa1bSDebojyoti Ghosh #undef __FUNCT__
7312302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSDestroy_GLEE"
7322302aa1bSDebojyoti Ghosh static PetscErrorCode TSDestroy_GLEE(TS ts)
7332302aa1bSDebojyoti Ghosh {
7342302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
7352302aa1bSDebojyoti Ghosh 
7362302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7372302aa1bSDebojyoti Ghosh   ierr = TSReset_GLEE(ts);CHKERRQ(ierr);
7382302aa1bSDebojyoti Ghosh   ierr = PetscFree(ts->data);CHKERRQ(ierr);
7392302aa1bSDebojyoti Ghosh   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLEEGetType_C",NULL);CHKERRQ(ierr);
7402302aa1bSDebojyoti Ghosh   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLEESetType_C",NULL);CHKERRQ(ierr);
7412302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
7422302aa1bSDebojyoti Ghosh }
7432302aa1bSDebojyoti Ghosh 
7442302aa1bSDebojyoti Ghosh #undef __FUNCT__
7452302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEEGetVecs"
7462302aa1bSDebojyoti Ghosh static PetscErrorCode TSGLEEGetVecs(TS ts,DM dm,Vec *Ydot)
7472302aa1bSDebojyoti Ghosh {
7482302aa1bSDebojyoti Ghosh   TS_GLEE     *glee = (TS_GLEE*)ts->data;
7492302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
7502302aa1bSDebojyoti Ghosh 
7512302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7522302aa1bSDebojyoti Ghosh   if (Ydot) {
7532302aa1bSDebojyoti Ghosh     if (dm && dm != ts->dm) {
7542302aa1bSDebojyoti Ghosh       ierr = DMGetNamedGlobalVector(dm,"TSGLEE_Ydot",Ydot);CHKERRQ(ierr);
7552302aa1bSDebojyoti Ghosh     } else *Ydot = glee->Ydot;
7562302aa1bSDebojyoti Ghosh   }
7572302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
7582302aa1bSDebojyoti Ghosh }
7592302aa1bSDebojyoti Ghosh 
7602302aa1bSDebojyoti Ghosh 
7612302aa1bSDebojyoti Ghosh #undef __FUNCT__
7622302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEERestoreVecs"
7632302aa1bSDebojyoti Ghosh static PetscErrorCode TSGLEERestoreVecs(TS ts,DM dm,Vec *Ydot)
7642302aa1bSDebojyoti Ghosh {
7652302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
7662302aa1bSDebojyoti Ghosh 
7672302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7682302aa1bSDebojyoti Ghosh   if (Ydot) {
7692302aa1bSDebojyoti Ghosh     if (dm && dm != ts->dm) {
7702302aa1bSDebojyoti Ghosh       ierr = DMRestoreNamedGlobalVector(dm,"TSGLEE_Ydot",Ydot);CHKERRQ(ierr);
7712302aa1bSDebojyoti Ghosh     }
7722302aa1bSDebojyoti Ghosh   }
7732302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
7742302aa1bSDebojyoti Ghosh }
7752302aa1bSDebojyoti Ghosh 
7762302aa1bSDebojyoti Ghosh /*
7772302aa1bSDebojyoti Ghosh   This defines the nonlinear equation that is to be solved with SNES
7782302aa1bSDebojyoti Ghosh */
7792302aa1bSDebojyoti Ghosh #undef __FUNCT__
7802302aa1bSDebojyoti Ghosh #define __FUNCT__ "SNESTSFormFunction_GLEE"
7812302aa1bSDebojyoti Ghosh static PetscErrorCode SNESTSFormFunction_GLEE(SNES snes,Vec X,Vec F,TS ts)
7822302aa1bSDebojyoti Ghosh {
7832302aa1bSDebojyoti Ghosh   TS_GLEE       *glee = (TS_GLEE*)ts->data;
7842302aa1bSDebojyoti Ghosh   DM             dm,dmsave;
7852302aa1bSDebojyoti Ghosh   Vec            Ydot;
7862302aa1bSDebojyoti Ghosh   PetscReal      shift = glee->scoeff / ts->time_step;
7872302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
7882302aa1bSDebojyoti Ghosh 
7892302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
7902302aa1bSDebojyoti Ghosh   ierr   = SNESGetDM(snes,&dm);CHKERRQ(ierr);
7912302aa1bSDebojyoti Ghosh   ierr   = TSGLEEGetVecs(ts,dm,&Ydot);CHKERRQ(ierr);
7922302aa1bSDebojyoti Ghosh   /* Set Ydot = shift*X */
7932302aa1bSDebojyoti Ghosh   ierr   = VecCopy(X,Ydot);CHKERRQ(ierr);
7942302aa1bSDebojyoti Ghosh   ierr   = VecScale(Ydot,shift);CHKERRQ(ierr);
7952302aa1bSDebojyoti Ghosh   dmsave = ts->dm;
7962302aa1bSDebojyoti Ghosh   ts->dm = dm;
7972302aa1bSDebojyoti Ghosh 
7982302aa1bSDebojyoti Ghosh   ierr = TSComputeIFunction(ts,glee->stage_time,X,Ydot,F,PETSC_FALSE);CHKERRQ(ierr);
7992302aa1bSDebojyoti Ghosh 
8002302aa1bSDebojyoti Ghosh   ts->dm = dmsave;
8012302aa1bSDebojyoti Ghosh   ierr   = TSGLEERestoreVecs(ts,dm,&Ydot);CHKERRQ(ierr);
8022302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8032302aa1bSDebojyoti Ghosh }
8042302aa1bSDebojyoti Ghosh 
8052302aa1bSDebojyoti Ghosh #undef __FUNCT__
8062302aa1bSDebojyoti Ghosh #define __FUNCT__ "SNESTSFormJacobian_GLEE"
8072302aa1bSDebojyoti Ghosh static PetscErrorCode SNESTSFormJacobian_GLEE(SNES snes,Vec X,Mat A,Mat B,TS ts)
8082302aa1bSDebojyoti Ghosh {
8092302aa1bSDebojyoti Ghosh   TS_GLEE        *glee = (TS_GLEE*)ts->data;
8102302aa1bSDebojyoti Ghosh   DM             dm,dmsave;
8112302aa1bSDebojyoti Ghosh   Vec            Ydot;
8122302aa1bSDebojyoti Ghosh   PetscReal      shift = glee->scoeff / ts->time_step;
8132302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
8142302aa1bSDebojyoti Ghosh 
8152302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8162302aa1bSDebojyoti Ghosh   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
8172302aa1bSDebojyoti Ghosh   ierr = TSGLEEGetVecs(ts,dm,&Ydot);CHKERRQ(ierr);
8182302aa1bSDebojyoti Ghosh   /* glee->Ydot has already been computed in SNESTSFormFunction_GLEE (SNES guarantees this) */
8192302aa1bSDebojyoti Ghosh   dmsave = ts->dm;
8202302aa1bSDebojyoti Ghosh   ts->dm = dm;
8212302aa1bSDebojyoti Ghosh 
8222302aa1bSDebojyoti Ghosh   ierr = TSComputeIJacobian(ts,glee->stage_time,X,Ydot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
8232302aa1bSDebojyoti Ghosh 
8242302aa1bSDebojyoti Ghosh   ts->dm = dmsave;
8252302aa1bSDebojyoti Ghosh   ierr   = TSGLEERestoreVecs(ts,dm,&Ydot);CHKERRQ(ierr);
8262302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8272302aa1bSDebojyoti Ghosh }
8282302aa1bSDebojyoti Ghosh 
8292302aa1bSDebojyoti Ghosh #undef __FUNCT__
8302302aa1bSDebojyoti Ghosh #define __FUNCT__ "DMCoarsenHook_TSGLEE"
8312302aa1bSDebojyoti Ghosh static PetscErrorCode DMCoarsenHook_TSGLEE(DM fine,DM coarse,void *ctx)
8322302aa1bSDebojyoti Ghosh {
8332302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8342302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8352302aa1bSDebojyoti Ghosh }
8362302aa1bSDebojyoti Ghosh 
8372302aa1bSDebojyoti Ghosh #undef __FUNCT__
8382302aa1bSDebojyoti Ghosh #define __FUNCT__ "DMRestrictHook_TSGLEE"
8392302aa1bSDebojyoti Ghosh static PetscErrorCode DMRestrictHook_TSGLEE(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
8402302aa1bSDebojyoti Ghosh {
8412302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8422302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8432302aa1bSDebojyoti Ghosh }
8442302aa1bSDebojyoti Ghosh 
8452302aa1bSDebojyoti Ghosh 
8462302aa1bSDebojyoti Ghosh #undef __FUNCT__
8472302aa1bSDebojyoti Ghosh #define __FUNCT__ "DMSubDomainHook_TSGLEE"
8482302aa1bSDebojyoti Ghosh static PetscErrorCode DMSubDomainHook_TSGLEE(DM dm,DM subdm,void *ctx)
8492302aa1bSDebojyoti Ghosh {
8502302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8512302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8522302aa1bSDebojyoti Ghosh }
8532302aa1bSDebojyoti Ghosh 
8542302aa1bSDebojyoti Ghosh #undef __FUNCT__
8552302aa1bSDebojyoti Ghosh #define __FUNCT__ "DMSubDomainRestrictHook_TSGLEE"
8562302aa1bSDebojyoti Ghosh static PetscErrorCode DMSubDomainRestrictHook_TSGLEE(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
8572302aa1bSDebojyoti Ghosh {
8582302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8592302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8602302aa1bSDebojyoti Ghosh }
8612302aa1bSDebojyoti Ghosh 
8622302aa1bSDebojyoti Ghosh #undef __FUNCT__
8632302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSSetUp_GLEE"
8642302aa1bSDebojyoti Ghosh static PetscErrorCode TSSetUp_GLEE(TS ts)
8652302aa1bSDebojyoti Ghosh {
8662302aa1bSDebojyoti Ghosh   TS_GLEE        *glee = (TS_GLEE*)ts->data;
8672302aa1bSDebojyoti Ghosh   GLEETableau    tab;
8682302aa1bSDebojyoti Ghosh   PetscInt       s,r;
8692302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
8702302aa1bSDebojyoti Ghosh   DM             dm;
8712302aa1bSDebojyoti Ghosh 
8722302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
8732302aa1bSDebojyoti Ghosh   if (!glee->tableau) {
8742302aa1bSDebojyoti Ghosh     ierr = TSGLEESetType(ts,TSGLEEDefaultType);CHKERRQ(ierr);
8752302aa1bSDebojyoti Ghosh   }
8762302aa1bSDebojyoti Ghosh   tab  = glee->tableau;
8772302aa1bSDebojyoti Ghosh   s    = tab->s;
8782302aa1bSDebojyoti Ghosh   r    = tab->r;
8792302aa1bSDebojyoti Ghosh   ierr = VecDuplicateVecs(ts->vec_sol,r,&glee->Y);CHKERRQ(ierr);
8802302aa1bSDebojyoti Ghosh   ierr = VecDuplicateVecs(ts->vec_sol,r,&glee->X);CHKERRQ(ierr);
8812302aa1bSDebojyoti Ghosh   ierr = VecDuplicateVecs(ts->vec_sol,s,&glee->YStage);CHKERRQ(ierr);
8822302aa1bSDebojyoti Ghosh   ierr = VecDuplicateVecs(ts->vec_sol,s,&glee->YdotStage);CHKERRQ(ierr);
8832302aa1bSDebojyoti Ghosh   ierr = VecDuplicate(ts->vec_sol,&glee->Ydot);CHKERRQ(ierr);
8842302aa1bSDebojyoti Ghosh   ierr = VecDuplicate(ts->vec_sol,&glee->W);CHKERRQ(ierr);
8852302aa1bSDebojyoti Ghosh   ierr = PetscMalloc1(s,&glee->work);CHKERRQ(ierr);
8862302aa1bSDebojyoti Ghosh   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
8872302aa1bSDebojyoti Ghosh   if (dm) {
8882302aa1bSDebojyoti Ghosh     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSGLEE,DMRestrictHook_TSGLEE,ts);CHKERRQ(ierr);
8892302aa1bSDebojyoti Ghosh     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSGLEE,DMSubDomainRestrictHook_TSGLEE,ts);CHKERRQ(ierr);
8902302aa1bSDebojyoti Ghosh   }
8912302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
8922302aa1bSDebojyoti Ghosh }
8932302aa1bSDebojyoti Ghosh 
8942302aa1bSDebojyoti Ghosh #undef __FUNCT__
8952302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSStartingMethod_GLEE"
8962302aa1bSDebojyoti Ghosh PetscErrorCode TSStartingMethod_GLEE(TS ts)
8972302aa1bSDebojyoti Ghosh {
8982302aa1bSDebojyoti Ghosh   TS_GLEE        *glee = (TS_GLEE*)ts->data;
8996e2e232bSDebojyoti Ghosh   GLEETableau    tab  = glee->tableau;
9006e2e232bSDebojyoti Ghosh   PetscInt       r=tab->r,i;
9016e2e232bSDebojyoti Ghosh   PetscReal      *S=tab->S;
9022302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
9032302aa1bSDebojyoti Ghosh 
9042302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
9052302aa1bSDebojyoti Ghosh   for (i=0; i<r; i++) {
9062302aa1bSDebojyoti Ghosh     ierr = VecZeroEntries(glee->Y[i]);CHKERRQ(ierr);
9072302aa1bSDebojyoti Ghosh     ierr = VecAXPY(glee->Y[i],S[i],ts->vec_sol);CHKERRQ(ierr);
9082302aa1bSDebojyoti Ghosh   }
9092302aa1bSDebojyoti Ghosh 
9102302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
9112302aa1bSDebojyoti Ghosh }
9122302aa1bSDebojyoti Ghosh 
9132302aa1bSDebojyoti Ghosh 
9142302aa1bSDebojyoti Ghosh /*------------------------------------------------------------*/
9152302aa1bSDebojyoti Ghosh 
9162302aa1bSDebojyoti Ghosh #undef __FUNCT__
9172302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSSetFromOptions_GLEE"
9182302aa1bSDebojyoti Ghosh static PetscErrorCode TSSetFromOptions_GLEE(PetscOptions *PetscOptionsObject,TS ts)
9192302aa1bSDebojyoti Ghosh {
9202302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
9212302aa1bSDebojyoti Ghosh   char           gleetype[256];
9222302aa1bSDebojyoti Ghosh 
9232302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
9242302aa1bSDebojyoti Ghosh   ierr = PetscOptionsHead(PetscOptionsObject,"GLEE ODE solver options");CHKERRQ(ierr);
9252302aa1bSDebojyoti Ghosh   {
9262302aa1bSDebojyoti Ghosh     GLEETableauLink link;
9272302aa1bSDebojyoti Ghosh     PetscInt        count,choice;
9282302aa1bSDebojyoti Ghosh     PetscBool       flg;
9292302aa1bSDebojyoti Ghosh     const char      **namelist;
9302302aa1bSDebojyoti Ghosh 
9312302aa1bSDebojyoti Ghosh     ierr = PetscStrncpy(gleetype,TSGLEEDefaultType,sizeof(gleetype));CHKERRQ(ierr);
9322302aa1bSDebojyoti Ghosh     for (link=GLEETableauList,count=0; link; link=link->next,count++) ;
9332302aa1bSDebojyoti Ghosh     ierr = PetscMalloc1(count,&namelist);CHKERRQ(ierr);
9342302aa1bSDebojyoti Ghosh     for (link=GLEETableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
9352302aa1bSDebojyoti Ghosh     ierr = PetscOptionsEList("-ts_glee_type","Family of GLEE method","TSGLEESetType",(const char*const*)namelist,count,gleetype,&choice,&flg);CHKERRQ(ierr);
9362302aa1bSDebojyoti Ghosh     ierr = TSGLEESetType(ts,flg ? namelist[choice] : gleetype);CHKERRQ(ierr);
9372302aa1bSDebojyoti Ghosh     ierr = PetscFree(namelist);CHKERRQ(ierr);
9382302aa1bSDebojyoti Ghosh   }
9392302aa1bSDebojyoti Ghosh   ierr = PetscOptionsTail();CHKERRQ(ierr);
9402302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
9412302aa1bSDebojyoti Ghosh }
9422302aa1bSDebojyoti Ghosh 
9432302aa1bSDebojyoti Ghosh #undef __FUNCT__
9442302aa1bSDebojyoti Ghosh #define __FUNCT__ "PetscFormatRealArray"
9452302aa1bSDebojyoti Ghosh static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
9462302aa1bSDebojyoti Ghosh {
9472302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
9482302aa1bSDebojyoti Ghosh   PetscInt       i;
9492302aa1bSDebojyoti Ghosh   size_t         left,count;
9502302aa1bSDebojyoti Ghosh   char           *p;
9512302aa1bSDebojyoti Ghosh 
9522302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
9532302aa1bSDebojyoti Ghosh   for (i=0,p=buf,left=len; i<n; i++) {
9542302aa1bSDebojyoti Ghosh     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
9552302aa1bSDebojyoti Ghosh     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
9562302aa1bSDebojyoti Ghosh     left -= count;
9572302aa1bSDebojyoti Ghosh     p    += count;
9582302aa1bSDebojyoti Ghosh     *p++  = ' ';
9592302aa1bSDebojyoti Ghosh   }
9602302aa1bSDebojyoti Ghosh   p[i ? 0 : -1] = 0;
9612302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
9622302aa1bSDebojyoti Ghosh }
9632302aa1bSDebojyoti Ghosh 
9642302aa1bSDebojyoti Ghosh #undef __FUNCT__
9652302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSView_GLEE"
9662302aa1bSDebojyoti Ghosh static PetscErrorCode TSView_GLEE(TS ts,PetscViewer viewer)
9672302aa1bSDebojyoti Ghosh {
9682302aa1bSDebojyoti Ghosh   TS_GLEE        *glee   = (TS_GLEE*)ts->data;
9692302aa1bSDebojyoti Ghosh   GLEETableau    tab  = glee->tableau;
9702302aa1bSDebojyoti Ghosh   PetscBool      iascii;
9712302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
9722302aa1bSDebojyoti Ghosh   TSAdapt        adapt;
9732302aa1bSDebojyoti Ghosh 
9742302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
9752302aa1bSDebojyoti Ghosh   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
9762302aa1bSDebojyoti Ghosh   if (iascii) {
9772302aa1bSDebojyoti Ghosh     TSGLEEType gleetype;
9782302aa1bSDebojyoti Ghosh     char   buf[512];
9792302aa1bSDebojyoti Ghosh     ierr = TSGLEEGetType(ts,&gleetype);CHKERRQ(ierr);
9802302aa1bSDebojyoti Ghosh     ierr = PetscViewerASCIIPrintf(viewer,"  GLEE %s\n",gleetype);CHKERRQ(ierr);
9812302aa1bSDebojyoti Ghosh     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
9822302aa1bSDebojyoti Ghosh     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa     c = %s\n",buf);CHKERRQ(ierr);
9832302aa1bSDebojyoti Ghosh     /* Note: print out r as well */
9842302aa1bSDebojyoti Ghosh   }
9852302aa1bSDebojyoti Ghosh   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
9862302aa1bSDebojyoti Ghosh   ierr = TSAdaptView(adapt,viewer);CHKERRQ(ierr);
9872302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
9882302aa1bSDebojyoti Ghosh }
9892302aa1bSDebojyoti Ghosh 
9902302aa1bSDebojyoti Ghosh #undef __FUNCT__
9912302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSLoad_GLEE"
9922302aa1bSDebojyoti Ghosh static PetscErrorCode TSLoad_GLEE(TS ts,PetscViewer viewer)
9932302aa1bSDebojyoti Ghosh {
9942302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
9952302aa1bSDebojyoti Ghosh   SNES           snes;
9962302aa1bSDebojyoti Ghosh   TSAdapt        tsadapt;
9972302aa1bSDebojyoti Ghosh 
9982302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
9992302aa1bSDebojyoti Ghosh   ierr = TSGetAdapt(ts,&tsadapt);CHKERRQ(ierr);
10002302aa1bSDebojyoti Ghosh   ierr = TSAdaptLoad(tsadapt,viewer);CHKERRQ(ierr);
10012302aa1bSDebojyoti Ghosh   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
10022302aa1bSDebojyoti Ghosh   ierr = SNESLoad(snes,viewer);CHKERRQ(ierr);
10032302aa1bSDebojyoti Ghosh   /* function and Jacobian context for SNES when used with TS is always ts object */
10042302aa1bSDebojyoti Ghosh   ierr = SNESSetFunction(snes,NULL,NULL,ts);CHKERRQ(ierr);
10052302aa1bSDebojyoti Ghosh   ierr = SNESSetJacobian(snes,NULL,NULL,NULL,ts);CHKERRQ(ierr);
10062302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
10072302aa1bSDebojyoti Ghosh }
10082302aa1bSDebojyoti Ghosh 
10092302aa1bSDebojyoti Ghosh #undef __FUNCT__
10102302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEESetType"
10112302aa1bSDebojyoti Ghosh /*@C
10122302aa1bSDebojyoti Ghosh   TSGLEESetType - Set the type of GLEE scheme
10132302aa1bSDebojyoti Ghosh 
10142302aa1bSDebojyoti Ghosh   Logically collective
10152302aa1bSDebojyoti Ghosh 
10162302aa1bSDebojyoti Ghosh   Input Parameter:
10172302aa1bSDebojyoti Ghosh +  ts - timestepping context
10182302aa1bSDebojyoti Ghosh -  gleetype - type of GLEE-scheme
10192302aa1bSDebojyoti Ghosh 
10202302aa1bSDebojyoti Ghosh   Level: intermediate
10212302aa1bSDebojyoti Ghosh 
10222302aa1bSDebojyoti Ghosh .seealso: TSGLEEGetType(), TSGLEE
10232302aa1bSDebojyoti Ghosh @*/
10242302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEESetType(TS ts,TSGLEEType gleetype)
10252302aa1bSDebojyoti Ghosh {
10262302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
10272302aa1bSDebojyoti Ghosh 
10282302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
10292302aa1bSDebojyoti Ghosh   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
10302302aa1bSDebojyoti Ghosh   ierr = PetscTryMethod(ts,"TSGLEESetType_C",(TS,TSGLEEType),(ts,gleetype));CHKERRQ(ierr);
10312302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
10322302aa1bSDebojyoti Ghosh }
10332302aa1bSDebojyoti Ghosh 
10342302aa1bSDebojyoti Ghosh #undef __FUNCT__
10352302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEEGetType"
10362302aa1bSDebojyoti Ghosh /*@C
10372302aa1bSDebojyoti Ghosh   TSGLEEGetType - Get the type of GLEE scheme
10382302aa1bSDebojyoti Ghosh 
10392302aa1bSDebojyoti Ghosh   Logically collective
10402302aa1bSDebojyoti Ghosh 
10412302aa1bSDebojyoti Ghosh   Input Parameter:
10422302aa1bSDebojyoti Ghosh .  ts - timestepping context
10432302aa1bSDebojyoti Ghosh 
10442302aa1bSDebojyoti Ghosh   Output Parameter:
10452302aa1bSDebojyoti Ghosh .  gleetype - type of GLEE-scheme
10462302aa1bSDebojyoti Ghosh 
10472302aa1bSDebojyoti Ghosh   Level: intermediate
10482302aa1bSDebojyoti Ghosh 
10492302aa1bSDebojyoti Ghosh .seealso: TSGLEESetType()
10502302aa1bSDebojyoti Ghosh @*/
10512302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEEGetType(TS ts,TSGLEEType *gleetype)
10522302aa1bSDebojyoti Ghosh {
10532302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
10542302aa1bSDebojyoti Ghosh 
10552302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
10562302aa1bSDebojyoti Ghosh   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
10572302aa1bSDebojyoti Ghosh   ierr = PetscUseMethod(ts,"TSGLEEGetType_C",(TS,TSGLEEType*),(ts,gleetype));CHKERRQ(ierr);
10582302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
10592302aa1bSDebojyoti Ghosh }
10602302aa1bSDebojyoti Ghosh 
10612302aa1bSDebojyoti Ghosh #undef __FUNCT__
10622302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEEGetType_GLEE"
10632302aa1bSDebojyoti Ghosh PetscErrorCode  TSGLEEGetType_GLEE(TS ts,TSGLEEType *gleetype)
10642302aa1bSDebojyoti Ghosh {
10652302aa1bSDebojyoti Ghosh   TS_GLEE     *glee = (TS_GLEE*)ts->data;
10662302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
10672302aa1bSDebojyoti Ghosh 
10682302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
10692302aa1bSDebojyoti Ghosh   if (!glee->tableau) {
10702302aa1bSDebojyoti Ghosh     ierr = TSGLEESetType(ts,TSGLEEDefaultType);CHKERRQ(ierr);
10712302aa1bSDebojyoti Ghosh   }
10722302aa1bSDebojyoti Ghosh   *gleetype = glee->tableau->name;
10732302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
10742302aa1bSDebojyoti Ghosh }
10752302aa1bSDebojyoti Ghosh #undef __FUNCT__
10762302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEESetType_GLEE"
10772302aa1bSDebojyoti Ghosh PetscErrorCode  TSGLEESetType_GLEE(TS ts,TSGLEEType gleetype)
10782302aa1bSDebojyoti Ghosh {
10792302aa1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
10802302aa1bSDebojyoti Ghosh   PetscErrorCode  ierr;
10812302aa1bSDebojyoti Ghosh   PetscBool       match;
10822302aa1bSDebojyoti Ghosh   GLEETableauLink link;
10832302aa1bSDebojyoti Ghosh 
10842302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
10852302aa1bSDebojyoti Ghosh   if (glee->tableau) {
10862302aa1bSDebojyoti Ghosh     ierr = PetscStrcmp(glee->tableau->name,gleetype,&match);CHKERRQ(ierr);
10872302aa1bSDebojyoti Ghosh     if (match) PetscFunctionReturn(0);
10882302aa1bSDebojyoti Ghosh   }
10892302aa1bSDebojyoti Ghosh   for (link = GLEETableauList; link; link=link->next) {
10902302aa1bSDebojyoti Ghosh     ierr = PetscStrcmp(link->tab.name,gleetype,&match);CHKERRQ(ierr);
10912302aa1bSDebojyoti Ghosh     if (match) {
10922302aa1bSDebojyoti Ghosh       ierr = TSReset_GLEE(ts);CHKERRQ(ierr);
10932302aa1bSDebojyoti Ghosh       glee->tableau = &link->tab;
10942302aa1bSDebojyoti Ghosh       PetscFunctionReturn(0);
10952302aa1bSDebojyoti Ghosh     }
10962302aa1bSDebojyoti Ghosh   }
10972302aa1bSDebojyoti Ghosh   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",gleetype);
10982302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
10992302aa1bSDebojyoti Ghosh }
11002302aa1bSDebojyoti Ghosh 
11012302aa1bSDebojyoti Ghosh #undef __FUNCT__
11022302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGetStages_GLEE"
11032302aa1bSDebojyoti Ghosh static PetscErrorCode  TSGetStages_GLEE(TS ts,PetscInt *ns,Vec **Y)
11042302aa1bSDebojyoti Ghosh {
11052302aa1bSDebojyoti Ghosh   TS_GLEE *glee = (TS_GLEE*)ts->data;
11062302aa1bSDebojyoti Ghosh 
11072302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
11082302aa1bSDebojyoti Ghosh   *ns = glee->tableau->s;
11092302aa1bSDebojyoti Ghosh   if(Y) *Y  = glee->Y;
11102302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
11112302aa1bSDebojyoti Ghosh }
11122302aa1bSDebojyoti Ghosh 
11132302aa1bSDebojyoti Ghosh #undef __FUNCT__
1114b2bf4f3aSDebojyoti Ghosh #define __FUNCT__ "TSGetSolutionComponents_GLEE"
1115b2bf4f3aSDebojyoti Ghosh PetscErrorCode TSGetSolutionComponents_GLEE(TS ts,PetscInt *n,Vec *Y)
11162302aa1bSDebojyoti Ghosh {
11172302aa1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
11182302aa1bSDebojyoti Ghosh   GLEETableau     tab   = glee->tableau;
11196e2e232bSDebojyoti Ghosh   PetscErrorCode  ierr;
11202302aa1bSDebojyoti Ghosh 
11212302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
11224cdd57e5SDebojyoti Ghosh   if (!Y) *n = tab->r;
11232302aa1bSDebojyoti Ghosh   else {
11244cdd57e5SDebojyoti Ghosh     if ((*n >= 0) && (*n < tab->r)) {
11254cdd57e5SDebojyoti Ghosh       ierr = VecCopy(glee->Y[*n],*Y); CHKERRQ(ierr);
1126*9657682dSDebojyoti Ghosh     } else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Second argument (%d) out of range[%d,%d].",*n,0,tab->r-1);
11272302aa1bSDebojyoti Ghosh   }
11282302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
11292302aa1bSDebojyoti Ghosh }
11302302aa1bSDebojyoti Ghosh 
11314cdd57e5SDebojyoti Ghosh #undef __FUNCT__
11324cdd57e5SDebojyoti Ghosh #define __FUNCT__ "TSGetAuxSolution_GLEE"
11334cdd57e5SDebojyoti Ghosh PetscErrorCode TSGetAuxSolution_GLEE(TS ts,Vec *X)
11344cdd57e5SDebojyoti Ghosh {
11354cdd57e5SDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
11364cdd57e5SDebojyoti Ghosh   GLEETableau     tab   = glee->tableau;
11374cdd57e5SDebojyoti Ghosh   PetscReal       *F    = tab->Fembed;
11384cdd57e5SDebojyoti Ghosh   PetscInt        r     = tab->r;
11394cdd57e5SDebojyoti Ghosh   Vec             *Y    = glee->Y;
11404cdd57e5SDebojyoti Ghosh   PetscErrorCode  ierr;
11414cdd57e5SDebojyoti Ghosh 
11424cdd57e5SDebojyoti Ghosh   PetscFunctionBegin;
11434cdd57e5SDebojyoti Ghosh   ierr = VecZeroEntries(*X);CHKERRQ(ierr);
11444cdd57e5SDebojyoti Ghosh   ierr = VecMAXPY((*X),r,F,Y);CHKERRQ(ierr);
11454cdd57e5SDebojyoti Ghosh   PetscFunctionReturn(0);
11464cdd57e5SDebojyoti Ghosh }
11474cdd57e5SDebojyoti Ghosh 
11484cdd57e5SDebojyoti Ghosh #undef __FUNCT__
11494cdd57e5SDebojyoti Ghosh #define __FUNCT__ "TSGetTimeError_GLEE"
11504cdd57e5SDebojyoti Ghosh PetscErrorCode TSGetTimeError_GLEE(TS ts,Vec *X)
11514cdd57e5SDebojyoti Ghosh {
11524cdd57e5SDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
11534cdd57e5SDebojyoti Ghosh   GLEETableau     tab   = glee->tableau;
11544cdd57e5SDebojyoti Ghosh   PetscReal       *F    = tab->Ferror;
11554cdd57e5SDebojyoti Ghosh   PetscInt        r     = tab->r;
11564cdd57e5SDebojyoti Ghosh   Vec             *Y    = glee->Y;
11574cdd57e5SDebojyoti Ghosh   PetscErrorCode  ierr;
11584cdd57e5SDebojyoti Ghosh 
11594cdd57e5SDebojyoti Ghosh   PetscFunctionBegin;
11604cdd57e5SDebojyoti Ghosh   ierr = VecZeroEntries(*X);CHKERRQ(ierr);
11614cdd57e5SDebojyoti Ghosh   ierr = VecMAXPY((*X),r,F,Y);CHKERRQ(ierr);
11624cdd57e5SDebojyoti Ghosh   PetscFunctionReturn(0);
11634cdd57e5SDebojyoti Ghosh }
11644cdd57e5SDebojyoti Ghosh 
116557df6a1bSDebojyoti Ghosh #undef __FUNCT__
116657df6a1bSDebojyoti Ghosh #define __FUNCT__ "TSSetTimeError_GLEE"
116757df6a1bSDebojyoti Ghosh PetscErrorCode TSSetTimeError_GLEE(TS ts,Vec X)
116857df6a1bSDebojyoti Ghosh {
116957df6a1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
117057df6a1bSDebojyoti Ghosh   GLEETableau     tab   = glee->tableau;
117157df6a1bSDebojyoti Ghosh   PetscReal       *S    = tab->Serror;
117257df6a1bSDebojyoti Ghosh   PetscInt        r     = tab->r,i;
117357df6a1bSDebojyoti Ghosh   Vec             *Y    = glee->Y;
117457df6a1bSDebojyoti Ghosh   PetscErrorCode  ierr;
117557df6a1bSDebojyoti Ghosh 
117657df6a1bSDebojyoti Ghosh   PetscFunctionBegin;
1177*9657682dSDebojyoti Ghosh   if (r != 2) SETERRQ2(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSSetTimeError_GLEE not supported for '%s' with r=%D.",tab->name,tab->r);
117857df6a1bSDebojyoti Ghosh   for (i=1; i<r; i++) {
117957df6a1bSDebojyoti Ghosh     ierr = VecCopy(ts->vec_sol,Y[i]); CHKERRQ(ierr);
118057df6a1bSDebojyoti Ghosh     ierr = VecAXPBY(Y[i],S[0],S[1],X); CHKERRQ(ierr);
118157df6a1bSDebojyoti Ghosh   }
118257df6a1bSDebojyoti Ghosh   PetscFunctionReturn(0);
118357df6a1bSDebojyoti Ghosh }
118457df6a1bSDebojyoti Ghosh 
11852302aa1bSDebojyoti Ghosh /* ------------------------------------------------------------ */
11862302aa1bSDebojyoti Ghosh /*MC
11872302aa1bSDebojyoti Ghosh       TSGLEE - ODE and DAE solver using General Linear with Error Estimation schemes
11882302aa1bSDebojyoti Ghosh 
11892302aa1bSDebojyoti Ghosh   The user should provide the right hand side of the equation
11902302aa1bSDebojyoti Ghosh   using TSSetRHSFunction().
11912302aa1bSDebojyoti Ghosh 
11922302aa1bSDebojyoti Ghosh   Notes:
11932302aa1bSDebojyoti Ghosh   The default is TSGLEE35, it can be changed with TSGLEESetType() or -ts_glee_type
11942302aa1bSDebojyoti Ghosh 
11952302aa1bSDebojyoti Ghosh   Level: beginner
11962302aa1bSDebojyoti Ghosh 
11972302aa1bSDebojyoti Ghosh .seealso:  TSCreate(), TS, TSSetType(), TSGLEESetType(), TSGLEEGetType(),
11982302aa1bSDebojyoti Ghosh            TSGLEE23, TTSGLEE24, TSGLEE35, TSGLEE25I, TSGLEEEXRK2A,
11992302aa1bSDebojyoti Ghosh            TSGLEERK32G1, TSGLEERK285EX, TSGLEEType, TSGLEERegister()
12002302aa1bSDebojyoti Ghosh 
12012302aa1bSDebojyoti Ghosh M*/
12022302aa1bSDebojyoti Ghosh #undef __FUNCT__
12032302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSCreate_GLEE"
12042302aa1bSDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSCreate_GLEE(TS ts)
12052302aa1bSDebojyoti Ghosh {
12062302aa1bSDebojyoti Ghosh   TS_GLEE         *th;
12072302aa1bSDebojyoti Ghosh   PetscErrorCode  ierr;
12082302aa1bSDebojyoti Ghosh 
12092302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
12102302aa1bSDebojyoti Ghosh #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
12112302aa1bSDebojyoti Ghosh   ierr = TSGLEEInitializePackage();CHKERRQ(ierr);
12122302aa1bSDebojyoti Ghosh #endif
12132302aa1bSDebojyoti Ghosh 
12142302aa1bSDebojyoti Ghosh   ts->ops->reset                  = TSReset_GLEE;
12152302aa1bSDebojyoti Ghosh   ts->ops->destroy                = TSDestroy_GLEE;
12162302aa1bSDebojyoti Ghosh   ts->ops->view                   = TSView_GLEE;
12172302aa1bSDebojyoti Ghosh   ts->ops->load                   = TSLoad_GLEE;
12182302aa1bSDebojyoti Ghosh   ts->ops->setup                  = TSSetUp_GLEE;
12192302aa1bSDebojyoti Ghosh   ts->ops->step                   = TSStep_GLEE;
12202302aa1bSDebojyoti Ghosh   ts->ops->interpolate            = TSInterpolate_GLEE;
12212302aa1bSDebojyoti Ghosh   ts->ops->evaluatestep           = TSEvaluateStep_GLEE;
12222302aa1bSDebojyoti Ghosh   ts->ops->setfromoptions         = TSSetFromOptions_GLEE;
12232302aa1bSDebojyoti Ghosh   ts->ops->getstages              = TSGetStages_GLEE;
12242302aa1bSDebojyoti Ghosh   ts->ops->snesfunction           = SNESTSFormFunction_GLEE;
12252302aa1bSDebojyoti Ghosh   ts->ops->snesjacobian           = SNESTSFormJacobian_GLEE;
1226b2bf4f3aSDebojyoti Ghosh   ts->ops->getsolutioncomponents  = TSGetSolutionComponents_GLEE;
12274cdd57e5SDebojyoti Ghosh   ts->ops->getauxsolution         = TSGetAuxSolution_GLEE;
12284cdd57e5SDebojyoti Ghosh   ts->ops->gettimeerror           = TSGetTimeError_GLEE;
122957df6a1bSDebojyoti Ghosh   ts->ops->settimeerror           = TSSetTimeError_GLEE;
12302302aa1bSDebojyoti Ghosh   ts->ops->startingmethod         = TSStartingMethod_GLEE;
12312302aa1bSDebojyoti Ghosh 
12322302aa1bSDebojyoti Ghosh   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
12332302aa1bSDebojyoti Ghosh   ts->data = (void*)th;
12342302aa1bSDebojyoti Ghosh 
12352302aa1bSDebojyoti Ghosh   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLEEGetType_C",TSGLEEGetType_GLEE);CHKERRQ(ierr);
12362302aa1bSDebojyoti Ghosh   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLEESetType_C",TSGLEESetType_GLEE);CHKERRQ(ierr);
12372302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
12382302aa1bSDebojyoti Ghosh }
1239