xref: /petsc/src/ts/impls/glee/glee.c (revision 2302aa1b855ebd0b8f717ac1679b74fa80c21cfb)
1*2302aa1bSDebojyoti Ghosh /*
2*2302aa1bSDebojyoti Ghosh   Code for time stepping with the General Linear with Error Estimation method
3*2302aa1bSDebojyoti Ghosh 
4*2302aa1bSDebojyoti Ghosh   Notes:
5*2302aa1bSDebojyoti Ghosh   The general system is written as
6*2302aa1bSDebojyoti Ghosh 
7*2302aa1bSDebojyoti Ghosh   Udot = F(t,U)
8*2302aa1bSDebojyoti Ghosh 
9*2302aa1bSDebojyoti Ghosh */
10*2302aa1bSDebojyoti Ghosh #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
11*2302aa1bSDebojyoti Ghosh #include <petscdm.h>
12*2302aa1bSDebojyoti Ghosh 
13*2302aa1bSDebojyoti Ghosh static TSGLEEType   TSGLEEDefaultType = TSGLEE35;
14*2302aa1bSDebojyoti Ghosh static PetscBool    TSGLEERegisterAllCalled;
15*2302aa1bSDebojyoti Ghosh static PetscBool    TSGLEEPackageInitialized;
16*2302aa1bSDebojyoti Ghosh static PetscInt     explicit_stage_time_id;
17*2302aa1bSDebojyoti Ghosh 
18*2302aa1bSDebojyoti Ghosh typedef struct _GLEETableau *GLEETableau;
19*2302aa1bSDebojyoti Ghosh struct _GLEETableau {
20*2302aa1bSDebojyoti Ghosh   char      *name;
21*2302aa1bSDebojyoti Ghosh   PetscInt   order;               /* Classical approximation order of the method i*/
22*2302aa1bSDebojyoti Ghosh   PetscInt   s;                   /* Number of stages */
23*2302aa1bSDebojyoti Ghosh   PetscInt   r;                   /* Number of steps */
24*2302aa1bSDebojyoti Ghosh   PetscReal  gamma;               /* LTE ratio */
25*2302aa1bSDebojyoti Ghosh   PetscReal *A,*B,*U,*V,*S,*F,*c; /* Tableau */
26*2302aa1bSDebojyoti Ghosh   PetscReal *Fembed;              /* Embedded final method coefficients */
27*2302aa1bSDebojyoti Ghosh   PetscInt   pinterp;             /* Interpolation order */
28*2302aa1bSDebojyoti Ghosh   PetscReal *binterp;             /* Interpolation coefficients */
29*2302aa1bSDebojyoti Ghosh   PetscReal  ccfl;                /* Placeholder for CFL coefficient relative to forward Euler  */
30*2302aa1bSDebojyoti Ghosh };
31*2302aa1bSDebojyoti Ghosh typedef struct _GLEETableauLink *GLEETableauLink;
32*2302aa1bSDebojyoti Ghosh struct _GLEETableauLink {
33*2302aa1bSDebojyoti Ghosh   struct _GLEETableau tab;
34*2302aa1bSDebojyoti Ghosh   GLEETableauLink     next;
35*2302aa1bSDebojyoti Ghosh };
36*2302aa1bSDebojyoti Ghosh static GLEETableauLink GLEETableauList;
37*2302aa1bSDebojyoti Ghosh 
38*2302aa1bSDebojyoti Ghosh typedef struct {
39*2302aa1bSDebojyoti Ghosh   GLEETableau  tableau;
40*2302aa1bSDebojyoti Ghosh   Vec          *Y;         /* Solution vector (along with auxiliary solution y~ or eps) */
41*2302aa1bSDebojyoti Ghosh   Vec          *X;         /* Temporary solution vector */
42*2302aa1bSDebojyoti Ghosh   Vec          *YStage;    /* Stage values */
43*2302aa1bSDebojyoti Ghosh   Vec          *YdotStage; /* Stage right hand side */
44*2302aa1bSDebojyoti Ghosh   Vec          W;          /* Right-hand-side for implicit stage solve */
45*2302aa1bSDebojyoti Ghosh   Vec          Ydot;       /* Work vector holding Ydot during residual evaluation */
46*2302aa1bSDebojyoti Ghosh   PetscScalar  *work;      /* Scalar work */
47*2302aa1bSDebojyoti Ghosh   PetscReal    scoeff;     /* shift = scoeff/dt */
48*2302aa1bSDebojyoti Ghosh   PetscReal    stage_time;
49*2302aa1bSDebojyoti Ghosh   TSStepStatus status;
50*2302aa1bSDebojyoti Ghosh } TS_GLEE;
51*2302aa1bSDebojyoti Ghosh 
52*2302aa1bSDebojyoti Ghosh /*MC
53*2302aa1bSDebojyoti Ghosh      TSGLEE23 - Second order three stage GLEE method
54*2302aa1bSDebojyoti Ghosh 
55*2302aa1bSDebojyoti Ghosh      This method has three stages.
56*2302aa1bSDebojyoti Ghosh      s = 3, r = 2
57*2302aa1bSDebojyoti Ghosh 
58*2302aa1bSDebojyoti Ghosh      Level: advanced
59*2302aa1bSDebojyoti Ghosh 
60*2302aa1bSDebojyoti Ghosh .seealso: TSGLEE
61*2302aa1bSDebojyoti Ghosh */
62*2302aa1bSDebojyoti Ghosh /*MC
63*2302aa1bSDebojyoti Ghosh      TSGLEE24 - Second order four stage GLEE method
64*2302aa1bSDebojyoti Ghosh 
65*2302aa1bSDebojyoti Ghosh      This method has four stages.
66*2302aa1bSDebojyoti Ghosh      s = 4, r = 2
67*2302aa1bSDebojyoti Ghosh 
68*2302aa1bSDebojyoti Ghosh      Level: advanced
69*2302aa1bSDebojyoti Ghosh 
70*2302aa1bSDebojyoti Ghosh .seealso: TSGLEE
71*2302aa1bSDebojyoti Ghosh */
72*2302aa1bSDebojyoti Ghosh /*MC
73*2302aa1bSDebojyoti Ghosh      TSGLEE25i - Second order five stage GLEE method
74*2302aa1bSDebojyoti Ghosh 
75*2302aa1bSDebojyoti Ghosh      This method has five stages.
76*2302aa1bSDebojyoti Ghosh      s = 5, r = 2
77*2302aa1bSDebojyoti Ghosh 
78*2302aa1bSDebojyoti Ghosh      Level: advanced
79*2302aa1bSDebojyoti Ghosh 
80*2302aa1bSDebojyoti Ghosh .seealso: TSGLEE
81*2302aa1bSDebojyoti Ghosh */
82*2302aa1bSDebojyoti Ghosh /*MC
83*2302aa1bSDebojyoti Ghosh      TSGLEE35  - Third order five stage GLEE method
84*2302aa1bSDebojyoti Ghosh 
85*2302aa1bSDebojyoti Ghosh      This method has five stages.
86*2302aa1bSDebojyoti Ghosh      s = 5, r = 2
87*2302aa1bSDebojyoti Ghosh 
88*2302aa1bSDebojyoti Ghosh      Level: advanced
89*2302aa1bSDebojyoti Ghosh 
90*2302aa1bSDebojyoti Ghosh .seealso: TSGLEE
91*2302aa1bSDebojyoti Ghosh */
92*2302aa1bSDebojyoti Ghosh /*MC
93*2302aa1bSDebojyoti Ghosh      TSGLEEEXRK2A  - Second order six stage GLEE method
94*2302aa1bSDebojyoti Ghosh 
95*2302aa1bSDebojyoti Ghosh      This method has six stages.
96*2302aa1bSDebojyoti Ghosh      s = 6, r = 2
97*2302aa1bSDebojyoti Ghosh 
98*2302aa1bSDebojyoti Ghosh      Level: advanced
99*2302aa1bSDebojyoti Ghosh 
100*2302aa1bSDebojyoti Ghosh .seealso: TSGLEE
101*2302aa1bSDebojyoti Ghosh */
102*2302aa1bSDebojyoti Ghosh /*MC
103*2302aa1bSDebojyoti Ghosh      TSGLEERK32G1  - Third order eight stage GLEE method
104*2302aa1bSDebojyoti Ghosh 
105*2302aa1bSDebojyoti Ghosh      This method has eight stages.
106*2302aa1bSDebojyoti Ghosh      s = 8, r = 2
107*2302aa1bSDebojyoti Ghosh 
108*2302aa1bSDebojyoti Ghosh      Level: advanced
109*2302aa1bSDebojyoti Ghosh 
110*2302aa1bSDebojyoti Ghosh .seealso: TSGLEE
111*2302aa1bSDebojyoti Ghosh */
112*2302aa1bSDebojyoti Ghosh /*MC
113*2302aa1bSDebojyoti Ghosh      TSGLEERK285EX  - Second order nine stage GLEE method
114*2302aa1bSDebojyoti Ghosh 
115*2302aa1bSDebojyoti Ghosh      This method has nine stages.
116*2302aa1bSDebojyoti Ghosh      s = 9, r = 2
117*2302aa1bSDebojyoti Ghosh 
118*2302aa1bSDebojyoti Ghosh      Level: advanced
119*2302aa1bSDebojyoti Ghosh 
120*2302aa1bSDebojyoti Ghosh .seealso: TSGLEE
121*2302aa1bSDebojyoti Ghosh */
122*2302aa1bSDebojyoti Ghosh 
123*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
124*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEERegisterAll"
125*2302aa1bSDebojyoti Ghosh /*@C
126*2302aa1bSDebojyoti Ghosh   TSGLEERegisterAll - Registers all of the General Linear with Error Estimation methods in TSGLEE
127*2302aa1bSDebojyoti Ghosh 
128*2302aa1bSDebojyoti Ghosh   Not Collective, but should be called by all processes which will need the schemes to be registered
129*2302aa1bSDebojyoti Ghosh 
130*2302aa1bSDebojyoti Ghosh   Level: advanced
131*2302aa1bSDebojyoti Ghosh 
132*2302aa1bSDebojyoti Ghosh .keywords: TS, TSGLEE, register, all
133*2302aa1bSDebojyoti Ghosh 
134*2302aa1bSDebojyoti Ghosh .seealso:  TSGLEERegisterDestroy()
135*2302aa1bSDebojyoti Ghosh @*/
136*2302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEERegisterAll(void)
137*2302aa1bSDebojyoti Ghosh {
138*2302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
139*2302aa1bSDebojyoti Ghosh 
140*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
141*2302aa1bSDebojyoti Ghosh   if (TSGLEERegisterAllCalled) PetscFunctionReturn(0);
142*2302aa1bSDebojyoti Ghosh   TSGLEERegisterAllCalled = PETSC_TRUE;
143*2302aa1bSDebojyoti Ghosh 
144*2302aa1bSDebojyoti Ghosh   {
145*2302aa1bSDebojyoti Ghosh     /* y-eps form */
146*2302aa1bSDebojyoti Ghosh     const PetscInt
147*2302aa1bSDebojyoti Ghosh       p = 2,
148*2302aa1bSDebojyoti Ghosh       s = 3,
149*2302aa1bSDebojyoti Ghosh       r = 2;
150*2302aa1bSDebojyoti Ghosh     const PetscReal
151*2302aa1bSDebojyoti Ghosh       gamma     = 0,
152*2302aa1bSDebojyoti Ghosh       A[3][3]   = {{0,0,0},{1,0,0},{0.25,0.25,0}},
153*2302aa1bSDebojyoti 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}},
154*2302aa1bSDebojyoti Ghosh       U[3][2]   = {{1,0},{1,10},{1,-1}},
155*2302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
156*2302aa1bSDebojyoti Ghosh       S[2]      = {1,0},
157*2302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
158*2302aa1bSDebojyoti Ghosh       Fembed[2] = {1,1};
159*2302aa1bSDebojyoti Ghosh     ierr = TSGLEERegister(TSGLEE23,p,s,r,gamma,&A[0][0],&B[0][0],&U[0][0],&V[0][0],S,F,NULL,Fembed,0,NULL);CHKERRQ(ierr);
160*2302aa1bSDebojyoti Ghosh   }
161*2302aa1bSDebojyoti Ghosh   {
162*2302aa1bSDebojyoti Ghosh     /* y-y~ form */
163*2302aa1bSDebojyoti Ghosh     const PetscInt
164*2302aa1bSDebojyoti Ghosh       p = 2,
165*2302aa1bSDebojyoti Ghosh       s = 4,
166*2302aa1bSDebojyoti Ghosh       r = 2;
167*2302aa1bSDebojyoti Ghosh     const PetscReal
168*2302aa1bSDebojyoti Ghosh       gamma     = 0,
169*2302aa1bSDebojyoti 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}},
170*2302aa1bSDebojyoti 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}},
171*2302aa1bSDebojyoti Ghosh       U[4][2]   = {{0,1},{75.0/58.0,-17.0/58.0},{0,1},{0,1}},
172*2302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
173*2302aa1bSDebojyoti Ghosh       S[2]      = {1,1},
174*2302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
175*2302aa1bSDebojyoti Ghosh       Fembed[2] = {0,1};
176*2302aa1bSDebojyoti Ghosh       ierr = TSGLEERegister(TSGLEE24,p,s,r,gamma,&A[0][0],&B[0][0],&U[0][0],&V[0][0],S,F,NULL,Fembed,0,NULL);CHKERRQ(ierr);
177*2302aa1bSDebojyoti Ghosh   }
178*2302aa1bSDebojyoti Ghosh   {
179*2302aa1bSDebojyoti Ghosh     /* y-y~ form */
180*2302aa1bSDebojyoti Ghosh     const PetscInt
181*2302aa1bSDebojyoti Ghosh       p = 2,
182*2302aa1bSDebojyoti Ghosh       s = 5,
183*2302aa1bSDebojyoti Ghosh       r = 2;
184*2302aa1bSDebojyoti Ghosh     const PetscReal
185*2302aa1bSDebojyoti Ghosh       gamma     = 0,
186*2302aa1bSDebojyoti Ghosh       A[5][5]   = {{0,0,0,0,0},
187*2302aa1bSDebojyoti Ghosh                    {-0.94079244066783383269,0,0,0,0},
188*2302aa1bSDebojyoti Ghosh                    {0.64228187778301907108,0.10915356933958500042,0,0,0},
189*2302aa1bSDebojyoti Ghosh                    {-0.51764297742287450812,0.74414270351096040738,-0.71404164927824538121,0,0},
190*2302aa1bSDebojyoti Ghosh                    {-0.44696561556825969206,-0.76768425657590196518,0.20111608138142987881,0.93828186737840469796,0}},
191*2302aa1bSDebojyoti Ghosh       B[2][5]   = {{-0.029309178948150356153,-0.49671981884013874923,0.34275801517650053274,0.32941112623949194988,0.85385985637229662276},
192*2302aa1bSDebojyoti Ghosh                    {0.78133219686062535272,0.074238691892675897635,0.57957363498384957966,-0.24638502829674959968,-0.18875949544040123033}},
193*2302aa1bSDebojyoti Ghosh       U[5][2]   = {{0.16911424754448327735,0.83088575245551672265},
194*2302aa1bSDebojyoti Ghosh                    {0.53638465733199574340,0.46361534266800425660},
195*2302aa1bSDebojyoti Ghosh                    {0.39901579167169582526,0.60098420832830417474},
196*2302aa1bSDebojyoti Ghosh                    {0.87689005530618575480,0.12310994469381424520},
197*2302aa1bSDebojyoti Ghosh                    {0.99056100455550913009,0.0094389954444908699092}},
198*2302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
199*2302aa1bSDebojyoti Ghosh       S[2]      = {1,1},
200*2302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
201*2302aa1bSDebojyoti Ghosh       Fembed[2] = {0,1};
202*2302aa1bSDebojyoti Ghosh     ierr = TSGLEERegister(TSGLEE25I,p,s,r,gamma,&A[0][0],&B[0][0],&U[0][0],&V[0][0],S,F,NULL,Fembed,0,NULL);CHKERRQ(ierr);
203*2302aa1bSDebojyoti Ghosh   }
204*2302aa1bSDebojyoti Ghosh   {
205*2302aa1bSDebojyoti Ghosh     /* y-y~ form */
206*2302aa1bSDebojyoti Ghosh     const PetscInt
207*2302aa1bSDebojyoti Ghosh       p = 3,
208*2302aa1bSDebojyoti Ghosh       s = 5,
209*2302aa1bSDebojyoti Ghosh       r = 2;
210*2302aa1bSDebojyoti Ghosh     const PetscReal
211*2302aa1bSDebojyoti Ghosh       gamma     = 0,
212*2302aa1bSDebojyoti Ghosh       A[5][5]   = {{0,0,0,0,0},
213*2302aa1bSDebojyoti Ghosh                    {- 2169604947363702313.0 /  24313474998937147335.0,0,0,0,0},
214*2302aa1bSDebojyoti Ghosh                    {46526746497697123895.0 /  94116917485856474137.0,-10297879244026594958.0 /  49199457603717988219.0,0,0,0},
215*2302aa1bSDebojyoti Ghosh                    {23364788935845982499.0 /  87425311444725389446.0,-79205144337496116638.0 / 148994349441340815519.0,40051189859317443782.0 /  36487615018004984309.0,0,0},
216*2302aa1bSDebojyoti Ghosh                    {42089522664062539205.0 / 124911313006412840286.0,-15074384760342762939.0 / 137927286865289746282.0,-62274678522253371016.0 / 125918573676298591413.0,13755475729852471739.0 /  79257927066651693390.0}},
217*2302aa1bSDebojyoti Ghosh       B[2][5]   = {{61546696837458703723.0 /  56982519523786160813.0,-55810892792806293355.0 / 206957624151308356511.0,24061048952676379087.0 / 158739347956038723465.0,3577972206874351339.0 /   7599733370677197135.0,-59449832954780563947.0 / 137360038685338563670.0},
218*2302aa1bSDebojyoti Ghosh                    {- 9738262186984159168.0 /  99299082461487742983.0,-32797097931948613195.0 /  61521565616362163366.0,42895514606418420631.0 /  71714201188501437336.0,22608567633166065068.0 /  55371917805607957003.0,94655809487476459565.0 / 151517167160302729021.0}},
219*2302aa1bSDebojyoti Ghosh       U[5][2]   = {{70820309139834661559.0 /  80863923579509469826.0,10043614439674808267.0 /  80863923579509469826.0},
220*2302aa1bSDebojyoti Ghosh                    {161694774978034105510.0 / 106187653640211060371.0,-55507121337823045139.0 / 106187653640211060371.0},
221*2302aa1bSDebojyoti Ghosh                    {78486094644566264568.0 /  88171030896733822981.0,9684936252167558413.0 /  88171030896733822981.0},
222*2302aa1bSDebojyoti Ghosh                    {65394922146334854435.0 /  84570853840405479554.0,19175931694070625119.0 /  84570853840405479554.0},
223*2302aa1bSDebojyoti Ghosh                    {8607282770183754108.0 / 108658046436496925911.0,100050763666313171803.0 / 108658046436496925911.0}},
224*2302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
225*2302aa1bSDebojyoti Ghosh       S[2]      = {1,1},
226*2302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
227*2302aa1bSDebojyoti Ghosh       Fembed[2] = {0,1};
228*2302aa1bSDebojyoti Ghosh     ierr = TSGLEERegister(TSGLEE35,p,s,r,gamma,&A[0][0],&B[0][0],&U[0][0],&V[0][0],S,F,NULL,Fembed,0,NULL);CHKERRQ(ierr);
229*2302aa1bSDebojyoti Ghosh   }
230*2302aa1bSDebojyoti Ghosh   {
231*2302aa1bSDebojyoti Ghosh     /* y-eps form */
232*2302aa1bSDebojyoti Ghosh     const PetscInt
233*2302aa1bSDebojyoti Ghosh       p = 2,
234*2302aa1bSDebojyoti Ghosh       s = 6,
235*2302aa1bSDebojyoti Ghosh       r = 2;
236*2302aa1bSDebojyoti Ghosh     const PetscReal
237*2302aa1bSDebojyoti Ghosh       gamma     = 0.25,
238*2302aa1bSDebojyoti Ghosh       A[6][6]   = {{0,0,0,0,0,0},
239*2302aa1bSDebojyoti Ghosh                    {1,0,0,0,0,0},
240*2302aa1bSDebojyoti Ghosh                    {0,0,0,0,0,0},
241*2302aa1bSDebojyoti Ghosh                    {0,0,0.5,0,0,0},
242*2302aa1bSDebojyoti Ghosh                    {0,0,0.25,0.25,0,0},
243*2302aa1bSDebojyoti Ghosh                    {0,0,0.25,0.25,0.5,0}},
244*2302aa1bSDebojyoti Ghosh       B[2][6]   = {{0.5,0.5,0,0,0,0},
245*2302aa1bSDebojyoti 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}},
246*2302aa1bSDebojyoti Ghosh       U[6][2]   = {{1,0},{1,0},{1,0.75},{1,0.75},{1,0.75},{1,0.75}},
247*2302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
248*2302aa1bSDebojyoti Ghosh       S[2]      = {1,0},
249*2302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
250*2302aa1bSDebojyoti Ghosh       Fembed[2] = {1,0.75};
251*2302aa1bSDebojyoti Ghosh     ierr = TSGLEERegister(TSGLEEEXRK2A,p,s,r,gamma,&A[0][0],&B[0][0],&U[0][0],&V[0][0],S,F,NULL,Fembed,0,NULL);CHKERRQ(ierr);
252*2302aa1bSDebojyoti Ghosh   }
253*2302aa1bSDebojyoti Ghosh   {
254*2302aa1bSDebojyoti Ghosh     /* y-eps form */
255*2302aa1bSDebojyoti Ghosh     const PetscInt
256*2302aa1bSDebojyoti Ghosh       p = 3,
257*2302aa1bSDebojyoti Ghosh       s = 8,
258*2302aa1bSDebojyoti Ghosh       r = 2;
259*2302aa1bSDebojyoti Ghosh     const PetscReal
260*2302aa1bSDebojyoti Ghosh       gamma     = 0,
261*2302aa1bSDebojyoti Ghosh       A[8][8]   = {{0,0,0,0,0,0,0,0},
262*2302aa1bSDebojyoti Ghosh                    {0.5,0,0,0,0,0,0,0},
263*2302aa1bSDebojyoti Ghosh                    {-1,2,0,0,0,0,0,0},
264*2302aa1bSDebojyoti Ghosh                    {1.0/6.0,2.0/3.0,1.0/6.0,0,0,0,0,0},
265*2302aa1bSDebojyoti Ghosh                    {0,0,0,0,0,0,0,0},
266*2302aa1bSDebojyoti Ghosh                    {-7.0/24.0,1.0/3.0,1.0/12.0,-1.0/8.0,0.5,0,0,0},
267*2302aa1bSDebojyoti Ghosh                    {7.0/6.0,-4.0/3.0,-1.0/3.0,0.5,-1.0,2.0,0,0},
268*2302aa1bSDebojyoti Ghosh                    {0,0,0,0,1.0/6.0,2.0/3.0,1.0/6.0,0}},
269*2302aa1bSDebojyoti Ghosh       B[2][8]   = {{1.0/6.0,2.0/3.0,1.0/6.0,0,0,0,0,0},
270*2302aa1bSDebojyoti 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}},
271*2302aa1bSDebojyoti Ghosh       U[8][2]   = {{1,0},{1,0},{1,0},{1,0},{1,1},{1,1},{1,1},{1,1}},
272*2302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
273*2302aa1bSDebojyoti Ghosh       S[2]      = {1,0},
274*2302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
275*2302aa1bSDebojyoti Ghosh       Fembed[2] = {1,1};
276*2302aa1bSDebojyoti Ghosh     ierr = TSGLEERegister(TSGLEERK32G1,p,s,r,gamma,&A[0][0],&B[0][0],&U[0][0],&V[0][0],S,F,NULL,Fembed,0,NULL);CHKERRQ(ierr);
277*2302aa1bSDebojyoti Ghosh   }
278*2302aa1bSDebojyoti Ghosh   {
279*2302aa1bSDebojyoti Ghosh     /* y-eps form */
280*2302aa1bSDebojyoti Ghosh     const PetscInt
281*2302aa1bSDebojyoti Ghosh       p = 2,
282*2302aa1bSDebojyoti Ghosh       s = 9,
283*2302aa1bSDebojyoti Ghosh       r = 2;
284*2302aa1bSDebojyoti Ghosh     const PetscReal
285*2302aa1bSDebojyoti Ghosh       gamma     = 0.25,
286*2302aa1bSDebojyoti Ghosh       A[9][9]   = {{0,0,0,0,0,0,0,0,0},
287*2302aa1bSDebojyoti Ghosh                    {0.585786437626904966,0,0,0,0,0,0,0,0},
288*2302aa1bSDebojyoti Ghosh                    {0.149999999999999994,0.849999999999999978,0,0,0,0,0,0,0},
289*2302aa1bSDebojyoti Ghosh                    {0,0,0,0,0,0,0,0,0},
290*2302aa1bSDebojyoti Ghosh                    {0,0,0,0.292893218813452483,0,0,0,0,0},
291*2302aa1bSDebojyoti Ghosh                    {0,0,0,0.0749999999999999972,0.424999999999999989,0,0,0,0},
292*2302aa1bSDebojyoti Ghosh                    {0,0,0,0.176776695296636893,0.176776695296636893,0.146446609406726241,0,0,0},
293*2302aa1bSDebojyoti Ghosh                    {0,0,0,0.176776695296636893,0.176776695296636893,0.146446609406726241,0.292893218813452483,0,0},
294*2302aa1bSDebojyoti Ghosh                    {0,0,0,0.176776695296636893,0.176776695296636893,0.146446609406726241,0.0749999999999999972,0.424999999999999989,0}},
295*2302aa1bSDebojyoti Ghosh       B[2][9]   = {{0.353553390593273786,0.353553390593273786,0.292893218813452483,0,0,0,0,0,0},
296*2302aa1bSDebojyoti Ghosh                    {-0.471404520791031678,-0.471404520791031678,-0.390524291751269959,0.235702260395515839,0.235702260395515839,0.195262145875634979,0.235702260395515839,0.235702260395515839,0.195262145875634979}},
297*2302aa1bSDebojyoti 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}},
298*2302aa1bSDebojyoti Ghosh       V[2][2]   = {{1,0},{0,1}},
299*2302aa1bSDebojyoti Ghosh       S[2]      = {1,0},
300*2302aa1bSDebojyoti Ghosh       F[2]      = {1,0},
301*2302aa1bSDebojyoti Ghosh       Fembed[2] = {1,0.75};
302*2302aa1bSDebojyoti Ghosh     ierr = TSGLEERegister(TSGLEERK285EX,p,s,r,gamma,&A[0][0],&B[0][0],&U[0][0],&V[0][0],S,F,NULL,Fembed,0,NULL);CHKERRQ(ierr);
303*2302aa1bSDebojyoti Ghosh   }
304*2302aa1bSDebojyoti Ghosh 
305*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
306*2302aa1bSDebojyoti Ghosh }
307*2302aa1bSDebojyoti Ghosh 
308*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
309*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEERegisterDestroy"
310*2302aa1bSDebojyoti Ghosh /*@C
311*2302aa1bSDebojyoti Ghosh    TSGLEERegisterDestroy - Frees the list of schemes that were registered by TSGLEERegister().
312*2302aa1bSDebojyoti Ghosh 
313*2302aa1bSDebojyoti Ghosh    Not Collective
314*2302aa1bSDebojyoti Ghosh 
315*2302aa1bSDebojyoti Ghosh    Level: advanced
316*2302aa1bSDebojyoti Ghosh 
317*2302aa1bSDebojyoti Ghosh .keywords: TSGLEE, register, destroy
318*2302aa1bSDebojyoti Ghosh .seealso: TSGLEERegister(), TSGLEERegisterAll()
319*2302aa1bSDebojyoti Ghosh @*/
320*2302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEERegisterDestroy(void)
321*2302aa1bSDebojyoti Ghosh {
322*2302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
323*2302aa1bSDebojyoti Ghosh   GLEETableauLink link;
324*2302aa1bSDebojyoti Ghosh 
325*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
326*2302aa1bSDebojyoti Ghosh   while ((link = GLEETableauList)) {
327*2302aa1bSDebojyoti Ghosh     GLEETableau t = &link->tab;
328*2302aa1bSDebojyoti Ghosh     GLEETableauList = link->next;
329*2302aa1bSDebojyoti Ghosh     ierr = PetscFree5(t->A,t->B,t->U,t->V,t->c);  CHKERRQ(ierr);
330*2302aa1bSDebojyoti Ghosh     ierr = PetscFree2(t->S,t->F);                 CHKERRQ(ierr);
331*2302aa1bSDebojyoti Ghosh     ierr = PetscFree (t->Fembed);                 CHKERRQ(ierr);
332*2302aa1bSDebojyoti Ghosh     ierr = PetscFree (t->binterp);                CHKERRQ(ierr);
333*2302aa1bSDebojyoti Ghosh     ierr = PetscFree (t->name);                   CHKERRQ(ierr);
334*2302aa1bSDebojyoti Ghosh     ierr = PetscFree (link);                      CHKERRQ(ierr);
335*2302aa1bSDebojyoti Ghosh   }
336*2302aa1bSDebojyoti Ghosh   TSGLEERegisterAllCalled = PETSC_FALSE;
337*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
338*2302aa1bSDebojyoti Ghosh }
339*2302aa1bSDebojyoti Ghosh 
340*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
341*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEEInitializePackage"
342*2302aa1bSDebojyoti Ghosh /*@C
343*2302aa1bSDebojyoti Ghosh   TSGLEEInitializePackage - This function initializes everything in the TSGLEE package. It is called
344*2302aa1bSDebojyoti Ghosh   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_GLEE()
345*2302aa1bSDebojyoti Ghosh   when using static libraries.
346*2302aa1bSDebojyoti Ghosh 
347*2302aa1bSDebojyoti Ghosh   Level: developer
348*2302aa1bSDebojyoti Ghosh 
349*2302aa1bSDebojyoti Ghosh .keywords: TS, TSGLEE, initialize, package
350*2302aa1bSDebojyoti Ghosh .seealso: PetscInitialize()
351*2302aa1bSDebojyoti Ghosh @*/
352*2302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEEInitializePackage(void)
353*2302aa1bSDebojyoti Ghosh {
354*2302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
355*2302aa1bSDebojyoti Ghosh 
356*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
357*2302aa1bSDebojyoti Ghosh   if (TSGLEEPackageInitialized) PetscFunctionReturn(0);
358*2302aa1bSDebojyoti Ghosh   TSGLEEPackageInitialized = PETSC_TRUE;
359*2302aa1bSDebojyoti Ghosh   ierr = TSGLEERegisterAll();CHKERRQ(ierr);
360*2302aa1bSDebojyoti Ghosh   ierr = PetscObjectComposedDataRegister(&explicit_stage_time_id);CHKERRQ(ierr);
361*2302aa1bSDebojyoti Ghosh   ierr = PetscRegisterFinalize(TSGLEEFinalizePackage);CHKERRQ(ierr);
362*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
363*2302aa1bSDebojyoti Ghosh }
364*2302aa1bSDebojyoti Ghosh 
365*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
366*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEEFinalizePackage"
367*2302aa1bSDebojyoti Ghosh /*@C
368*2302aa1bSDebojyoti Ghosh   TSGLEEFinalizePackage - This function destroys everything in the TSGLEE package. It is
369*2302aa1bSDebojyoti Ghosh   called from PetscFinalize().
370*2302aa1bSDebojyoti Ghosh 
371*2302aa1bSDebojyoti Ghosh   Level: developer
372*2302aa1bSDebojyoti Ghosh 
373*2302aa1bSDebojyoti Ghosh .keywords: Petsc, destroy, package
374*2302aa1bSDebojyoti Ghosh .seealso: PetscFinalize()
375*2302aa1bSDebojyoti Ghosh @*/
376*2302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEEFinalizePackage(void)
377*2302aa1bSDebojyoti Ghosh {
378*2302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
379*2302aa1bSDebojyoti Ghosh 
380*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
381*2302aa1bSDebojyoti Ghosh   TSGLEEPackageInitialized = PETSC_FALSE;
382*2302aa1bSDebojyoti Ghosh   ierr = TSGLEERegisterDestroy();CHKERRQ(ierr);
383*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
384*2302aa1bSDebojyoti Ghosh }
385*2302aa1bSDebojyoti Ghosh 
386*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
387*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEERegister"
388*2302aa1bSDebojyoti Ghosh /*@C
389*2302aa1bSDebojyoti Ghosh    TSGLEERegister - register an GLEE scheme by providing the entries in the Butcher tableau
390*2302aa1bSDebojyoti Ghosh 
391*2302aa1bSDebojyoti Ghosh    Not Collective, but the same schemes should be registered on all processes on which they will be used
392*2302aa1bSDebojyoti Ghosh 
393*2302aa1bSDebojyoti Ghosh    Input Parameters:
394*2302aa1bSDebojyoti Ghosh +  name   - identifier for method
395*2302aa1bSDebojyoti Ghosh .  order  - order of method
396*2302aa1bSDebojyoti Ghosh .  s - number of stages
397*2302aa1bSDebojyoti Ghosh .  r - number of steps
398*2302aa1bSDebojyoti Ghosh .  gamma - LTE ratio
399*2302aa1bSDebojyoti Ghosh .  A - stage coefficients (dimension s*s, row-major)
400*2302aa1bSDebojyoti Ghosh .  B - step completion coefficients (dimension r*s, row-major)
401*2302aa1bSDebojyoti Ghosh .  U - method coefficients (dimension s*r, row-major)
402*2302aa1bSDebojyoti Ghosh .  V - method coefficients (dimension r*r, row-major)
403*2302aa1bSDebojyoti Ghosh .  S - starting coefficients
404*2302aa1bSDebojyoti Ghosh .  F - finishing coefficients
405*2302aa1bSDebojyoti Ghosh .  c - abscissa (dimension s; NULL to use row sums of A)
406*2302aa1bSDebojyoti Ghosh .  Fembed - step completion coefficients for embedded method
407*2302aa1bSDebojyoti Ghosh .  pinterp - order of interpolation (0 if unavailable)
408*2302aa1bSDebojyoti Ghosh -  binterp - array of interpolation coefficients (NULL if unavailable)
409*2302aa1bSDebojyoti Ghosh 
410*2302aa1bSDebojyoti Ghosh    Notes:
411*2302aa1bSDebojyoti Ghosh    Several GLEE methods are provided, this function is only needed to create new methods.
412*2302aa1bSDebojyoti Ghosh 
413*2302aa1bSDebojyoti Ghosh    Level: advanced
414*2302aa1bSDebojyoti Ghosh 
415*2302aa1bSDebojyoti Ghosh .keywords: TS, register
416*2302aa1bSDebojyoti Ghosh 
417*2302aa1bSDebojyoti Ghosh .seealso: TSGLEE
418*2302aa1bSDebojyoti Ghosh @*/
419*2302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEERegister(TSGLEEType name,PetscInt order,PetscInt s, PetscInt r,
420*2302aa1bSDebojyoti Ghosh                               PetscReal gamma,
421*2302aa1bSDebojyoti Ghosh                               const PetscReal A[],const PetscReal B[],
422*2302aa1bSDebojyoti Ghosh                               const PetscReal U[],const PetscReal V[],
423*2302aa1bSDebojyoti Ghosh                               const PetscReal S[],const PetscReal F[],
424*2302aa1bSDebojyoti Ghosh                               const PetscReal c[],const PetscReal Fembed[],
425*2302aa1bSDebojyoti Ghosh                               PetscInt pinterp, const PetscReal binterp[])
426*2302aa1bSDebojyoti Ghosh {
427*2302aa1bSDebojyoti Ghosh   PetscErrorCode    ierr;
428*2302aa1bSDebojyoti Ghosh   GLEETableauLink   link;
429*2302aa1bSDebojyoti Ghosh   GLEETableau       t;
430*2302aa1bSDebojyoti Ghosh   PetscInt          i,j;
431*2302aa1bSDebojyoti Ghosh 
432*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
433*2302aa1bSDebojyoti Ghosh   ierr     = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
434*2302aa1bSDebojyoti Ghosh   ierr     = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
435*2302aa1bSDebojyoti Ghosh   t        = &link->tab;
436*2302aa1bSDebojyoti Ghosh   ierr     = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
437*2302aa1bSDebojyoti Ghosh   t->order = order;
438*2302aa1bSDebojyoti Ghosh   t->s     = s;
439*2302aa1bSDebojyoti Ghosh   t->r     = r;
440*2302aa1bSDebojyoti Ghosh   t->gamma = gamma;
441*2302aa1bSDebojyoti Ghosh   ierr     = PetscMalloc5(s*s,&t->A,r*r,&t->V,s,&t->c,r*s,&t->B,s*r,&t->U); CHKERRQ(ierr);
442*2302aa1bSDebojyoti Ghosh   ierr     = PetscMalloc2(r,&t->S,r,&t->F); CHKERRQ(ierr);
443*2302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
444*2302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->B,B,r*s*sizeof(B[0]));CHKERRQ(ierr);
445*2302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->U,U,s*r*sizeof(U[0]));CHKERRQ(ierr);
446*2302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->V,V,r*r*sizeof(V[0]));CHKERRQ(ierr);
447*2302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->S,S,r  *sizeof(S[0]));CHKERRQ(ierr);
448*2302aa1bSDebojyoti Ghosh   ierr     = PetscMemcpy(t->F,F,r  *sizeof(F[0]));CHKERRQ(ierr);
449*2302aa1bSDebojyoti Ghosh   if (c) {
450*2302aa1bSDebojyoti Ghosh     ierr   = PetscMemcpy(t->c,c,s*sizeof(c[0])); CHKERRQ(ierr);
451*2302aa1bSDebojyoti Ghosh   } else {
452*2302aa1bSDebojyoti Ghosh     for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j];
453*2302aa1bSDebojyoti Ghosh   }
454*2302aa1bSDebojyoti Ghosh   ierr = PetscMalloc1(r,&t->Fembed);CHKERRQ(ierr);
455*2302aa1bSDebojyoti Ghosh   ierr = PetscMemcpy(t->Fembed,Fembed,r*sizeof(Fembed[0]));CHKERRQ(ierr);
456*2302aa1bSDebojyoti Ghosh   t->pinterp = pinterp;
457*2302aa1bSDebojyoti Ghosh   ierr       = PetscMalloc1(s*pinterp,&t->binterp);CHKERRQ(ierr);
458*2302aa1bSDebojyoti Ghosh   ierr       = PetscMemcpy(t->binterp,binterp,s*pinterp*sizeof(binterp[0]));CHKERRQ(ierr);
459*2302aa1bSDebojyoti Ghosh 
460*2302aa1bSDebojyoti Ghosh   link->next      = GLEETableauList;
461*2302aa1bSDebojyoti Ghosh   GLEETableauList = link;
462*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
463*2302aa1bSDebojyoti Ghosh }
464*2302aa1bSDebojyoti Ghosh 
465*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
466*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSEvaluateStep_GLEE"
467*2302aa1bSDebojyoti Ghosh static PetscErrorCode TSEvaluateStep_GLEE(TS ts,PetscInt order,Vec X,PetscBool *done)
468*2302aa1bSDebojyoti Ghosh {
469*2302aa1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*) ts->data;
470*2302aa1bSDebojyoti Ghosh   GLEETableau     tab = glee->tableau;
471*2302aa1bSDebojyoti Ghosh   PetscReal       h, *B = tab->B, *V = tab->V,
472*2302aa1bSDebojyoti Ghosh                   *F = tab->F, *Fembed = tab->Fembed;
473*2302aa1bSDebojyoti Ghosh   PetscInt        s = tab->s, r = tab->r, i, j;
474*2302aa1bSDebojyoti Ghosh   Vec             *Y = glee->Y, *YdotStage = glee->YdotStage;
475*2302aa1bSDebojyoti Ghosh   PetscScalar     *w = glee->work;
476*2302aa1bSDebojyoti Ghosh   PetscErrorCode  ierr;
477*2302aa1bSDebojyoti Ghosh 
478*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
479*2302aa1bSDebojyoti Ghosh 
480*2302aa1bSDebojyoti Ghosh   switch (glee->status) {
481*2302aa1bSDebojyoti Ghosh     case TS_STEP_INCOMPLETE:
482*2302aa1bSDebojyoti Ghosh     case TS_STEP_PENDING:
483*2302aa1bSDebojyoti Ghosh       h = ts->time_step; break;
484*2302aa1bSDebojyoti Ghosh     case TS_STEP_COMPLETE:
485*2302aa1bSDebojyoti Ghosh       h = ts->time_step_prev; break;
486*2302aa1bSDebojyoti Ghosh     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
487*2302aa1bSDebojyoti Ghosh   }
488*2302aa1bSDebojyoti Ghosh 
489*2302aa1bSDebojyoti Ghosh   if (order == tab->order) {
490*2302aa1bSDebojyoti Ghosh 
491*2302aa1bSDebojyoti Ghosh     /* Note: Irrespective of whether status is TS_STEP_INCOMPLETE
492*2302aa1bSDebojyoti Ghosh              or TS_STEP_COMPLETE, glee-X has the solution at the
493*2302aa1bSDebojyoti Ghosh              beginning of the time step. So no need to roll-back.
494*2302aa1bSDebojyoti Ghosh     */
495*2302aa1bSDebojyoti Ghosh 
496*2302aa1bSDebojyoti Ghosh     if (glee->status == TS_STEP_INCOMPLETE) {
497*2302aa1bSDebojyoti Ghosh       for (i=0; i<r; i++) {
498*2302aa1bSDebojyoti Ghosh         ierr = VecZeroEntries(Y[i]); CHKERRQ(ierr);
499*2302aa1bSDebojyoti Ghosh         ierr = VecMAXPY(Y[i],r,(V+i*r),glee->X); CHKERRQ(ierr);
500*2302aa1bSDebojyoti Ghosh         for (j=0; j<s; j++) w[j] = h*B[i*s+j];
501*2302aa1bSDebojyoti Ghosh         ierr = VecMAXPY(Y[i],s,w,YdotStage); CHKERRQ(ierr);
502*2302aa1bSDebojyoti Ghosh       }
503*2302aa1bSDebojyoti Ghosh       ierr = VecZeroEntries(X);CHKERRQ(ierr);
504*2302aa1bSDebojyoti Ghosh       ierr = VecMAXPY(X,r,F,Y);CHKERRQ(ierr);
505*2302aa1bSDebojyoti Ghosh     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
506*2302aa1bSDebojyoti Ghosh     PetscFunctionReturn(0);
507*2302aa1bSDebojyoti Ghosh 
508*2302aa1bSDebojyoti Ghosh   } else if (order == tab->order-1) {
509*2302aa1bSDebojyoti Ghosh 
510*2302aa1bSDebojyoti Ghosh     /* Complete with the embedded method (Fembed) */
511*2302aa1bSDebojyoti Ghosh     for (i=0; i<r; i++) {
512*2302aa1bSDebojyoti Ghosh       ierr = VecZeroEntries(Y[i]); CHKERRQ(ierr);
513*2302aa1bSDebojyoti Ghosh       ierr = VecMAXPY(Y[i],r,(V+i*r),glee->X); CHKERRQ(ierr);
514*2302aa1bSDebojyoti Ghosh       for (j=0; j<s; j++) w[j] = h*B[i*s+j];
515*2302aa1bSDebojyoti Ghosh       ierr = VecMAXPY(Y[i],s,w,YdotStage); CHKERRQ(ierr);
516*2302aa1bSDebojyoti Ghosh     }
517*2302aa1bSDebojyoti Ghosh     ierr = VecZeroEntries(X);CHKERRQ(ierr);
518*2302aa1bSDebojyoti Ghosh     ierr = VecMAXPY(X,r,Fembed,Y);CHKERRQ(ierr);
519*2302aa1bSDebojyoti Ghosh     if (done) *done = PETSC_TRUE;
520*2302aa1bSDebojyoti Ghosh     PetscFunctionReturn(0);
521*2302aa1bSDebojyoti Ghosh   }
522*2302aa1bSDebojyoti Ghosh   if (done) *done = PETSC_FALSE;
523*2302aa1bSDebojyoti 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);
524*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
525*2302aa1bSDebojyoti Ghosh }
526*2302aa1bSDebojyoti Ghosh 
527*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
528*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSStep_GLEE"
529*2302aa1bSDebojyoti Ghosh static PetscErrorCode TSStep_GLEE(TS ts)
530*2302aa1bSDebojyoti Ghosh {
531*2302aa1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
532*2302aa1bSDebojyoti Ghosh   GLEETableau     tab = glee->tableau;
533*2302aa1bSDebojyoti Ghosh   const PetscInt  s = tab->s, r = tab->r;
534*2302aa1bSDebojyoti Ghosh   PetscReal       *A = tab->A, *U = tab->U,
535*2302aa1bSDebojyoti Ghosh                   *F = tab->F,
536*2302aa1bSDebojyoti Ghosh                   *c = tab->c;
537*2302aa1bSDebojyoti Ghosh   Vec             *Y = glee->Y, *X = glee->X,
538*2302aa1bSDebojyoti Ghosh                   *YStage = glee->YStage,
539*2302aa1bSDebojyoti Ghosh                   *YdotStage = glee->YdotStage,
540*2302aa1bSDebojyoti Ghosh                   W = glee->W;
541*2302aa1bSDebojyoti Ghosh   SNES            snes;
542*2302aa1bSDebojyoti Ghosh   PetscScalar     *w = glee->work;
543*2302aa1bSDebojyoti Ghosh   TSAdapt         adapt;
544*2302aa1bSDebojyoti Ghosh   PetscInt        i,j,reject,next_scheme,its,lits;
545*2302aa1bSDebojyoti Ghosh   PetscReal       next_time_step;
546*2302aa1bSDebojyoti Ghosh   PetscReal       t;
547*2302aa1bSDebojyoti Ghosh   PetscBool       accept;
548*2302aa1bSDebojyoti Ghosh   PetscErrorCode  ierr;
549*2302aa1bSDebojyoti Ghosh 
550*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
551*2302aa1bSDebojyoti Ghosh   for (i=0; i<r; i++) { ierr = VecCopy(Y[i],X[i]); CHKERRQ(ierr); }
552*2302aa1bSDebojyoti Ghosh 
553*2302aa1bSDebojyoti Ghosh   ierr           = TSGetSNES(ts,&snes);CHKERRQ(ierr);
554*2302aa1bSDebojyoti Ghosh   next_time_step = ts->time_step;
555*2302aa1bSDebojyoti Ghosh   t              = ts->ptime;
556*2302aa1bSDebojyoti Ghosh   accept         = PETSC_TRUE;
557*2302aa1bSDebojyoti Ghosh   glee->status   = TS_STEP_INCOMPLETE;
558*2302aa1bSDebojyoti Ghosh 
559*2302aa1bSDebojyoti Ghosh   for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) {
560*2302aa1bSDebojyoti Ghosh 
561*2302aa1bSDebojyoti Ghosh     PetscReal h = ts->time_step;
562*2302aa1bSDebojyoti Ghosh     ierr = TSPreStep(ts);CHKERRQ(ierr);
563*2302aa1bSDebojyoti Ghosh 
564*2302aa1bSDebojyoti Ghosh     for (i=0; i<s; i++) {
565*2302aa1bSDebojyoti Ghosh 
566*2302aa1bSDebojyoti Ghosh       glee->stage_time = t + h*c[i];
567*2302aa1bSDebojyoti Ghosh       ierr = TSPreStage(ts,glee->stage_time); CHKERRQ(ierr);
568*2302aa1bSDebojyoti Ghosh 
569*2302aa1bSDebojyoti Ghosh       if (A[i*s+i] == 0) {  /* Explicit stage */
570*2302aa1bSDebojyoti Ghosh         ierr = VecZeroEntries(YStage[i]);CHKERRQ(ierr);
571*2302aa1bSDebojyoti Ghosh         ierr = VecMAXPY(YStage[i],r,(U+i*r),X);CHKERRQ(ierr);
572*2302aa1bSDebojyoti Ghosh         for (j=0; j<i; j++) w[j] = h*A[i*s+j];
573*2302aa1bSDebojyoti Ghosh         ierr = VecMAXPY(YStage[i],i,w,YdotStage);CHKERRQ(ierr);
574*2302aa1bSDebojyoti Ghosh       } else {              /* Implicit stage */
575*2302aa1bSDebojyoti Ghosh         glee->scoeff = 1.0/A[i*s+i];
576*2302aa1bSDebojyoti Ghosh         /* compute right-hand-side */
577*2302aa1bSDebojyoti Ghosh         ierr = VecZeroEntries(W);CHKERRQ(ierr);
578*2302aa1bSDebojyoti Ghosh         ierr = VecMAXPY(W,r,(U+i*r),X);CHKERRQ(ierr);
579*2302aa1bSDebojyoti Ghosh         for (j=0; j<i; j++) w[j] = h*A[i*s+j];
580*2302aa1bSDebojyoti Ghosh         ierr = VecMAXPY(W,i,w,YdotStage);CHKERRQ(ierr);
581*2302aa1bSDebojyoti Ghosh         ierr = VecScale(W,glee->scoeff/h);CHKERRQ(ierr);
582*2302aa1bSDebojyoti Ghosh         /* set initial guess */
583*2302aa1bSDebojyoti Ghosh         ierr = VecCopy(i>0 ? YStage[i-1] : ts->vec_sol,YStage[i]);CHKERRQ(ierr);
584*2302aa1bSDebojyoti Ghosh         /* solve for this stage */
585*2302aa1bSDebojyoti Ghosh         ierr = SNESSolve(snes,W,YStage[i]);CHKERRQ(ierr);
586*2302aa1bSDebojyoti Ghosh         ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
587*2302aa1bSDebojyoti Ghosh         ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
588*2302aa1bSDebojyoti Ghosh         ts->snes_its += its; ts->ksp_its += lits;
589*2302aa1bSDebojyoti Ghosh       }
590*2302aa1bSDebojyoti Ghosh       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
591*2302aa1bSDebojyoti Ghosh       ierr = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr);
592*2302aa1bSDebojyoti Ghosh       if (!accept) goto reject_step;
593*2302aa1bSDebojyoti Ghosh       ierr = TSPostStage(ts,glee->stage_time,i,YStage); CHKERRQ(ierr);
594*2302aa1bSDebojyoti Ghosh       ierr = TSComputeRHSFunction(ts,t+h*c[i],YStage[i],YdotStage[i]);CHKERRQ(ierr);
595*2302aa1bSDebojyoti Ghosh     }
596*2302aa1bSDebojyoti Ghosh     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
597*2302aa1bSDebojyoti Ghosh     glee->status = TS_STEP_PENDING;
598*2302aa1bSDebojyoti Ghosh 
599*2302aa1bSDebojyoti Ghosh     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
600*2302aa1bSDebojyoti Ghosh     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
601*2302aa1bSDebojyoti Ghosh     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
602*2302aa1bSDebojyoti Ghosh     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
603*2302aa1bSDebojyoti Ghosh     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
604*2302aa1bSDebojyoti Ghosh     if (accept) {
605*2302aa1bSDebojyoti Ghosh       /* ignore next_scheme for now */
606*2302aa1bSDebojyoti Ghosh       ts->ptime     += ts->time_step;
607*2302aa1bSDebojyoti Ghosh       ts->time_step  = next_time_step;
608*2302aa1bSDebojyoti Ghosh       ts->steps++;
609*2302aa1bSDebojyoti Ghosh       glee->status = TS_STEP_COMPLETE;
610*2302aa1bSDebojyoti Ghosh       ierr = PetscObjectComposedDataSetReal((PetscObject)ts->vec_sol,explicit_stage_time_id,ts->ptime);CHKERRQ(ierr);
611*2302aa1bSDebojyoti Ghosh       break;
612*2302aa1bSDebojyoti Ghosh     } else {                    /* Roll back the current step */
613*2302aa1bSDebojyoti Ghosh       ierr = VecMAXPY(ts->vec_sol,r,F,X);CHKERRQ(ierr);
614*2302aa1bSDebojyoti Ghosh       ts->time_step = next_time_step;
615*2302aa1bSDebojyoti Ghosh       glee->status  = TS_STEP_INCOMPLETE;
616*2302aa1bSDebojyoti Ghosh     }
617*2302aa1bSDebojyoti Ghosh reject_step: continue;
618*2302aa1bSDebojyoti Ghosh   }
619*2302aa1bSDebojyoti Ghosh   if (glee->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
620*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
621*2302aa1bSDebojyoti Ghosh }
622*2302aa1bSDebojyoti Ghosh 
623*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
624*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSInterpolate_GLEE"
625*2302aa1bSDebojyoti Ghosh static PetscErrorCode TSInterpolate_GLEE(TS ts,PetscReal itime,Vec X)
626*2302aa1bSDebojyoti Ghosh {
627*2302aa1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
628*2302aa1bSDebojyoti Ghosh   PetscInt        s=glee->tableau->s, pinterp=glee->tableau->pinterp,i,j;
629*2302aa1bSDebojyoti Ghosh   PetscReal       h,tt,t;
630*2302aa1bSDebojyoti Ghosh   PetscScalar     *b;
631*2302aa1bSDebojyoti Ghosh   const PetscReal *B = glee->tableau->binterp;
632*2302aa1bSDebojyoti Ghosh   PetscErrorCode  ierr;
633*2302aa1bSDebojyoti Ghosh 
634*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
635*2302aa1bSDebojyoti Ghosh   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSGLEE %s does not have an interpolation formula",glee->tableau->name);
636*2302aa1bSDebojyoti Ghosh   switch (glee->status) {
637*2302aa1bSDebojyoti Ghosh     case TS_STEP_INCOMPLETE:
638*2302aa1bSDebojyoti Ghosh     case TS_STEP_PENDING:
639*2302aa1bSDebojyoti Ghosh       h = ts->time_step;
640*2302aa1bSDebojyoti Ghosh       t = (itime - ts->ptime)/h;
641*2302aa1bSDebojyoti Ghosh       break;
642*2302aa1bSDebojyoti Ghosh     case TS_STEP_COMPLETE:
643*2302aa1bSDebojyoti Ghosh       h = ts->time_step_prev;
644*2302aa1bSDebojyoti Ghosh       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
645*2302aa1bSDebojyoti Ghosh       break;
646*2302aa1bSDebojyoti Ghosh     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
647*2302aa1bSDebojyoti Ghosh   }
648*2302aa1bSDebojyoti Ghosh   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
649*2302aa1bSDebojyoti Ghosh   for (i=0; i<s; i++) b[i] = 0;
650*2302aa1bSDebojyoti Ghosh   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
651*2302aa1bSDebojyoti Ghosh     for (i=0; i<s; i++) {
652*2302aa1bSDebojyoti Ghosh       b[i]  += h * B[i*pinterp+j] * tt;
653*2302aa1bSDebojyoti Ghosh     }
654*2302aa1bSDebojyoti Ghosh   }
655*2302aa1bSDebojyoti Ghosh   ierr = VecCopy(glee->YStage[0],X);CHKERRQ(ierr);
656*2302aa1bSDebojyoti Ghosh   ierr = VecMAXPY(X,s,b,glee->YdotStage);CHKERRQ(ierr);
657*2302aa1bSDebojyoti Ghosh   ierr = PetscFree(b);CHKERRQ(ierr);
658*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
659*2302aa1bSDebojyoti Ghosh }
660*2302aa1bSDebojyoti Ghosh 
661*2302aa1bSDebojyoti Ghosh /*------------------------------------------------------------*/
662*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
663*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSReset_GLEE"
664*2302aa1bSDebojyoti Ghosh static PetscErrorCode TSReset_GLEE(TS ts)
665*2302aa1bSDebojyoti Ghosh {
666*2302aa1bSDebojyoti Ghosh   TS_GLEE        *glee = (TS_GLEE*)ts->data;
667*2302aa1bSDebojyoti Ghosh   PetscInt       s, r;
668*2302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
669*2302aa1bSDebojyoti Ghosh 
670*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
671*2302aa1bSDebojyoti Ghosh   if (!glee->tableau) PetscFunctionReturn(0);
672*2302aa1bSDebojyoti Ghosh   s    = glee->tableau->s;
673*2302aa1bSDebojyoti Ghosh   r    = glee->tableau->r;
674*2302aa1bSDebojyoti Ghosh   ierr = VecDestroyVecs(r,&glee->Y);CHKERRQ(ierr);
675*2302aa1bSDebojyoti Ghosh   ierr = VecDestroyVecs(r,&glee->X);CHKERRQ(ierr);
676*2302aa1bSDebojyoti Ghosh   ierr = VecDestroyVecs(s,&glee->YStage);CHKERRQ(ierr);
677*2302aa1bSDebojyoti Ghosh   ierr = VecDestroyVecs(s,&glee->YdotStage);CHKERRQ(ierr);
678*2302aa1bSDebojyoti Ghosh   ierr = VecDestroy(&glee->Ydot);CHKERRQ(ierr);
679*2302aa1bSDebojyoti Ghosh   ierr = VecDestroy(&glee->W);CHKERRQ(ierr);
680*2302aa1bSDebojyoti Ghosh   ierr = PetscFree(glee->work);CHKERRQ(ierr);
681*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
682*2302aa1bSDebojyoti Ghosh }
683*2302aa1bSDebojyoti Ghosh 
684*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
685*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSDestroy_GLEE"
686*2302aa1bSDebojyoti Ghosh static PetscErrorCode TSDestroy_GLEE(TS ts)
687*2302aa1bSDebojyoti Ghosh {
688*2302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
689*2302aa1bSDebojyoti Ghosh 
690*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
691*2302aa1bSDebojyoti Ghosh   ierr = TSReset_GLEE(ts);CHKERRQ(ierr);
692*2302aa1bSDebojyoti Ghosh   ierr = PetscFree(ts->data);CHKERRQ(ierr);
693*2302aa1bSDebojyoti Ghosh   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLEEGetType_C",NULL);CHKERRQ(ierr);
694*2302aa1bSDebojyoti Ghosh   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLEESetType_C",NULL);CHKERRQ(ierr);
695*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
696*2302aa1bSDebojyoti Ghosh }
697*2302aa1bSDebojyoti Ghosh 
698*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
699*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEEGetVecs"
700*2302aa1bSDebojyoti Ghosh static PetscErrorCode TSGLEEGetVecs(TS ts,DM dm,Vec *Ydot)
701*2302aa1bSDebojyoti Ghosh {
702*2302aa1bSDebojyoti Ghosh   TS_GLEE     *glee = (TS_GLEE*)ts->data;
703*2302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
704*2302aa1bSDebojyoti Ghosh 
705*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
706*2302aa1bSDebojyoti Ghosh   if (Ydot) {
707*2302aa1bSDebojyoti Ghosh     if (dm && dm != ts->dm) {
708*2302aa1bSDebojyoti Ghosh       ierr = DMGetNamedGlobalVector(dm,"TSGLEE_Ydot",Ydot);CHKERRQ(ierr);
709*2302aa1bSDebojyoti Ghosh     } else *Ydot = glee->Ydot;
710*2302aa1bSDebojyoti Ghosh   }
711*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
712*2302aa1bSDebojyoti Ghosh }
713*2302aa1bSDebojyoti Ghosh 
714*2302aa1bSDebojyoti Ghosh 
715*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
716*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEERestoreVecs"
717*2302aa1bSDebojyoti Ghosh static PetscErrorCode TSGLEERestoreVecs(TS ts,DM dm,Vec *Ydot)
718*2302aa1bSDebojyoti Ghosh {
719*2302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
720*2302aa1bSDebojyoti Ghosh 
721*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
722*2302aa1bSDebojyoti Ghosh   if (Ydot) {
723*2302aa1bSDebojyoti Ghosh     if (dm && dm != ts->dm) {
724*2302aa1bSDebojyoti Ghosh       ierr = DMRestoreNamedGlobalVector(dm,"TSGLEE_Ydot",Ydot);CHKERRQ(ierr);
725*2302aa1bSDebojyoti Ghosh     }
726*2302aa1bSDebojyoti Ghosh   }
727*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
728*2302aa1bSDebojyoti Ghosh }
729*2302aa1bSDebojyoti Ghosh 
730*2302aa1bSDebojyoti Ghosh /*
731*2302aa1bSDebojyoti Ghosh   This defines the nonlinear equation that is to be solved with SNES
732*2302aa1bSDebojyoti Ghosh */
733*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
734*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "SNESTSFormFunction_GLEE"
735*2302aa1bSDebojyoti Ghosh static PetscErrorCode SNESTSFormFunction_GLEE(SNES snes,Vec X,Vec F,TS ts)
736*2302aa1bSDebojyoti Ghosh {
737*2302aa1bSDebojyoti Ghosh   TS_GLEE       *glee = (TS_GLEE*)ts->data;
738*2302aa1bSDebojyoti Ghosh   DM             dm,dmsave;
739*2302aa1bSDebojyoti Ghosh   Vec            Ydot;
740*2302aa1bSDebojyoti Ghosh   PetscReal      shift = glee->scoeff / ts->time_step;
741*2302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
742*2302aa1bSDebojyoti Ghosh 
743*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
744*2302aa1bSDebojyoti Ghosh   ierr   = SNESGetDM(snes,&dm);CHKERRQ(ierr);
745*2302aa1bSDebojyoti Ghosh   ierr   = TSGLEEGetVecs(ts,dm,&Ydot);CHKERRQ(ierr);
746*2302aa1bSDebojyoti Ghosh   /* Set Ydot = shift*X */
747*2302aa1bSDebojyoti Ghosh   ierr   = VecCopy(X,Ydot);CHKERRQ(ierr);
748*2302aa1bSDebojyoti Ghosh   ierr   = VecScale(Ydot,shift);CHKERRQ(ierr);
749*2302aa1bSDebojyoti Ghosh   dmsave = ts->dm;
750*2302aa1bSDebojyoti Ghosh   ts->dm = dm;
751*2302aa1bSDebojyoti Ghosh 
752*2302aa1bSDebojyoti Ghosh   ierr = TSComputeIFunction(ts,glee->stage_time,X,Ydot,F,PETSC_FALSE);CHKERRQ(ierr);
753*2302aa1bSDebojyoti Ghosh 
754*2302aa1bSDebojyoti Ghosh   ts->dm = dmsave;
755*2302aa1bSDebojyoti Ghosh   ierr   = TSGLEERestoreVecs(ts,dm,&Ydot);CHKERRQ(ierr);
756*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
757*2302aa1bSDebojyoti Ghosh }
758*2302aa1bSDebojyoti Ghosh 
759*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
760*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "SNESTSFormJacobian_GLEE"
761*2302aa1bSDebojyoti Ghosh static PetscErrorCode SNESTSFormJacobian_GLEE(SNES snes,Vec X,Mat A,Mat B,TS ts)
762*2302aa1bSDebojyoti Ghosh {
763*2302aa1bSDebojyoti Ghosh   TS_GLEE        *glee = (TS_GLEE*)ts->data;
764*2302aa1bSDebojyoti Ghosh   DM             dm,dmsave;
765*2302aa1bSDebojyoti Ghosh   Vec            Ydot;
766*2302aa1bSDebojyoti Ghosh   PetscReal      shift = glee->scoeff / ts->time_step;
767*2302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
768*2302aa1bSDebojyoti Ghosh 
769*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
770*2302aa1bSDebojyoti Ghosh   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
771*2302aa1bSDebojyoti Ghosh   ierr = TSGLEEGetVecs(ts,dm,&Ydot);CHKERRQ(ierr);
772*2302aa1bSDebojyoti Ghosh   /* glee->Ydot has already been computed in SNESTSFormFunction_GLEE (SNES guarantees this) */
773*2302aa1bSDebojyoti Ghosh   dmsave = ts->dm;
774*2302aa1bSDebojyoti Ghosh   ts->dm = dm;
775*2302aa1bSDebojyoti Ghosh 
776*2302aa1bSDebojyoti Ghosh   ierr = TSComputeIJacobian(ts,glee->stage_time,X,Ydot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
777*2302aa1bSDebojyoti Ghosh 
778*2302aa1bSDebojyoti Ghosh   ts->dm = dmsave;
779*2302aa1bSDebojyoti Ghosh   ierr   = TSGLEERestoreVecs(ts,dm,&Ydot);CHKERRQ(ierr);
780*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
781*2302aa1bSDebojyoti Ghosh }
782*2302aa1bSDebojyoti Ghosh 
783*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
784*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "DMCoarsenHook_TSGLEE"
785*2302aa1bSDebojyoti Ghosh static PetscErrorCode DMCoarsenHook_TSGLEE(DM fine,DM coarse,void *ctx)
786*2302aa1bSDebojyoti Ghosh {
787*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
788*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
789*2302aa1bSDebojyoti Ghosh }
790*2302aa1bSDebojyoti Ghosh 
791*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
792*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "DMRestrictHook_TSGLEE"
793*2302aa1bSDebojyoti Ghosh static PetscErrorCode DMRestrictHook_TSGLEE(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
794*2302aa1bSDebojyoti Ghosh {
795*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
796*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
797*2302aa1bSDebojyoti Ghosh }
798*2302aa1bSDebojyoti Ghosh 
799*2302aa1bSDebojyoti Ghosh 
800*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
801*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "DMSubDomainHook_TSGLEE"
802*2302aa1bSDebojyoti Ghosh static PetscErrorCode DMSubDomainHook_TSGLEE(DM dm,DM subdm,void *ctx)
803*2302aa1bSDebojyoti Ghosh {
804*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
805*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
806*2302aa1bSDebojyoti Ghosh }
807*2302aa1bSDebojyoti Ghosh 
808*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
809*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "DMSubDomainRestrictHook_TSGLEE"
810*2302aa1bSDebojyoti Ghosh static PetscErrorCode DMSubDomainRestrictHook_TSGLEE(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
811*2302aa1bSDebojyoti Ghosh {
812*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
813*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
814*2302aa1bSDebojyoti Ghosh }
815*2302aa1bSDebojyoti Ghosh 
816*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
817*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSSetUp_GLEE"
818*2302aa1bSDebojyoti Ghosh static PetscErrorCode TSSetUp_GLEE(TS ts)
819*2302aa1bSDebojyoti Ghosh {
820*2302aa1bSDebojyoti Ghosh   TS_GLEE        *glee = (TS_GLEE*)ts->data;
821*2302aa1bSDebojyoti Ghosh   GLEETableau    tab;
822*2302aa1bSDebojyoti Ghosh   PetscInt       s,r;
823*2302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
824*2302aa1bSDebojyoti Ghosh   DM             dm;
825*2302aa1bSDebojyoti Ghosh 
826*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
827*2302aa1bSDebojyoti Ghosh   if (!glee->tableau) {
828*2302aa1bSDebojyoti Ghosh     ierr = TSGLEESetType(ts,TSGLEEDefaultType);CHKERRQ(ierr);
829*2302aa1bSDebojyoti Ghosh   }
830*2302aa1bSDebojyoti Ghosh   tab  = glee->tableau;
831*2302aa1bSDebojyoti Ghosh   s    = tab->s;
832*2302aa1bSDebojyoti Ghosh   r    = tab->r;
833*2302aa1bSDebojyoti Ghosh   ierr = VecDuplicateVecs(ts->vec_sol,r,&glee->Y);CHKERRQ(ierr);
834*2302aa1bSDebojyoti Ghosh   ierr = VecDuplicateVecs(ts->vec_sol,r,&glee->X);CHKERRQ(ierr);
835*2302aa1bSDebojyoti Ghosh   ierr = VecDuplicateVecs(ts->vec_sol,s,&glee->YStage);CHKERRQ(ierr);
836*2302aa1bSDebojyoti Ghosh   ierr = VecDuplicateVecs(ts->vec_sol,s,&glee->YdotStage);CHKERRQ(ierr);
837*2302aa1bSDebojyoti Ghosh   ierr = VecDuplicate(ts->vec_sol,&glee->Ydot);CHKERRQ(ierr);
838*2302aa1bSDebojyoti Ghosh   ierr = VecDuplicate(ts->vec_sol,&glee->W);CHKERRQ(ierr);
839*2302aa1bSDebojyoti Ghosh   ierr = PetscMalloc1(s,&glee->work);CHKERRQ(ierr);
840*2302aa1bSDebojyoti Ghosh   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
841*2302aa1bSDebojyoti Ghosh   if (dm) {
842*2302aa1bSDebojyoti Ghosh     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSGLEE,DMRestrictHook_TSGLEE,ts);CHKERRQ(ierr);
843*2302aa1bSDebojyoti Ghosh     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSGLEE,DMSubDomainRestrictHook_TSGLEE,ts);CHKERRQ(ierr);
844*2302aa1bSDebojyoti Ghosh   }
845*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
846*2302aa1bSDebojyoti Ghosh }
847*2302aa1bSDebojyoti Ghosh 
848*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
849*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSStartingMethod_GLEE"
850*2302aa1bSDebojyoti Ghosh PetscErrorCode TSStartingMethod_GLEE(TS ts)
851*2302aa1bSDebojyoti Ghosh {
852*2302aa1bSDebojyoti Ghosh   TS_GLEE        *glee = (TS_GLEE*)ts->data;
853*2302aa1bSDebojyoti Ghosh   GLEETableau    tab;
854*2302aa1bSDebojyoti Ghosh   PetscInt       r,i;
855*2302aa1bSDebojyoti Ghosh   PetscReal      *S;
856*2302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
857*2302aa1bSDebojyoti Ghosh 
858*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
859*2302aa1bSDebojyoti Ghosh   if (!glee->tableau) {
860*2302aa1bSDebojyoti Ghosh     ierr = TSGLEESetType(ts,TSGLEEDefaultType);CHKERRQ(ierr);
861*2302aa1bSDebojyoti Ghosh   }
862*2302aa1bSDebojyoti Ghosh   tab  = glee->tableau;
863*2302aa1bSDebojyoti Ghosh   r    = tab->r;
864*2302aa1bSDebojyoti Ghosh   S    = tab->S;
865*2302aa1bSDebojyoti Ghosh   ierr = VecDuplicateVecs(ts->vec_sol,r,&glee->Y);CHKERRQ(ierr);
866*2302aa1bSDebojyoti Ghosh   for (i=0; i<r; i++) {
867*2302aa1bSDebojyoti Ghosh     ierr = VecZeroEntries(glee->Y[i]);CHKERRQ(ierr);
868*2302aa1bSDebojyoti Ghosh     ierr = VecAXPY(glee->Y[i],S[i],ts->vec_sol);CHKERRQ(ierr);
869*2302aa1bSDebojyoti Ghosh   }
870*2302aa1bSDebojyoti Ghosh 
871*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
872*2302aa1bSDebojyoti Ghosh }
873*2302aa1bSDebojyoti Ghosh 
874*2302aa1bSDebojyoti Ghosh 
875*2302aa1bSDebojyoti Ghosh /*------------------------------------------------------------*/
876*2302aa1bSDebojyoti Ghosh 
877*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
878*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSSetFromOptions_GLEE"
879*2302aa1bSDebojyoti Ghosh static PetscErrorCode TSSetFromOptions_GLEE(PetscOptions *PetscOptionsObject,TS ts)
880*2302aa1bSDebojyoti Ghosh {
881*2302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
882*2302aa1bSDebojyoti Ghosh   char           gleetype[256];
883*2302aa1bSDebojyoti Ghosh 
884*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
885*2302aa1bSDebojyoti Ghosh   ierr = PetscOptionsHead(PetscOptionsObject,"GLEE ODE solver options");CHKERRQ(ierr);
886*2302aa1bSDebojyoti Ghosh   {
887*2302aa1bSDebojyoti Ghosh     GLEETableauLink link;
888*2302aa1bSDebojyoti Ghosh     PetscInt        count,choice;
889*2302aa1bSDebojyoti Ghosh     PetscBool       flg;
890*2302aa1bSDebojyoti Ghosh     const char      **namelist;
891*2302aa1bSDebojyoti Ghosh 
892*2302aa1bSDebojyoti Ghosh     ierr = PetscStrncpy(gleetype,TSGLEEDefaultType,sizeof(gleetype));CHKERRQ(ierr);
893*2302aa1bSDebojyoti Ghosh     for (link=GLEETableauList,count=0; link; link=link->next,count++) ;
894*2302aa1bSDebojyoti Ghosh     ierr = PetscMalloc1(count,&namelist);CHKERRQ(ierr);
895*2302aa1bSDebojyoti Ghosh     for (link=GLEETableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
896*2302aa1bSDebojyoti Ghosh     ierr = PetscOptionsEList("-ts_glee_type","Family of GLEE method","TSGLEESetType",(const char*const*)namelist,count,gleetype,&choice,&flg);CHKERRQ(ierr);
897*2302aa1bSDebojyoti Ghosh     ierr = TSGLEESetType(ts,flg ? namelist[choice] : gleetype);CHKERRQ(ierr);
898*2302aa1bSDebojyoti Ghosh     ierr = PetscFree(namelist);CHKERRQ(ierr);
899*2302aa1bSDebojyoti Ghosh   }
900*2302aa1bSDebojyoti Ghosh   ierr = PetscOptionsTail();CHKERRQ(ierr);
901*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
902*2302aa1bSDebojyoti Ghosh }
903*2302aa1bSDebojyoti Ghosh 
904*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
905*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "PetscFormatRealArray"
906*2302aa1bSDebojyoti Ghosh static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
907*2302aa1bSDebojyoti Ghosh {
908*2302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
909*2302aa1bSDebojyoti Ghosh   PetscInt       i;
910*2302aa1bSDebojyoti Ghosh   size_t         left,count;
911*2302aa1bSDebojyoti Ghosh   char           *p;
912*2302aa1bSDebojyoti Ghosh 
913*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
914*2302aa1bSDebojyoti Ghosh   for (i=0,p=buf,left=len; i<n; i++) {
915*2302aa1bSDebojyoti Ghosh     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
916*2302aa1bSDebojyoti Ghosh     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
917*2302aa1bSDebojyoti Ghosh     left -= count;
918*2302aa1bSDebojyoti Ghosh     p    += count;
919*2302aa1bSDebojyoti Ghosh     *p++  = ' ';
920*2302aa1bSDebojyoti Ghosh   }
921*2302aa1bSDebojyoti Ghosh   p[i ? 0 : -1] = 0;
922*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
923*2302aa1bSDebojyoti Ghosh }
924*2302aa1bSDebojyoti Ghosh 
925*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
926*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSView_GLEE"
927*2302aa1bSDebojyoti Ghosh static PetscErrorCode TSView_GLEE(TS ts,PetscViewer viewer)
928*2302aa1bSDebojyoti Ghosh {
929*2302aa1bSDebojyoti Ghosh   TS_GLEE        *glee   = (TS_GLEE*)ts->data;
930*2302aa1bSDebojyoti Ghosh   GLEETableau    tab  = glee->tableau;
931*2302aa1bSDebojyoti Ghosh   PetscBool      iascii;
932*2302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
933*2302aa1bSDebojyoti Ghosh   TSAdapt        adapt;
934*2302aa1bSDebojyoti Ghosh 
935*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
936*2302aa1bSDebojyoti Ghosh   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
937*2302aa1bSDebojyoti Ghosh   if (iascii) {
938*2302aa1bSDebojyoti Ghosh     TSGLEEType gleetype;
939*2302aa1bSDebojyoti Ghosh     char   buf[512];
940*2302aa1bSDebojyoti Ghosh     ierr = TSGLEEGetType(ts,&gleetype);CHKERRQ(ierr);
941*2302aa1bSDebojyoti Ghosh     ierr = PetscViewerASCIIPrintf(viewer,"  GLEE %s\n",gleetype);CHKERRQ(ierr);
942*2302aa1bSDebojyoti Ghosh     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
943*2302aa1bSDebojyoti Ghosh     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa     c = %s\n",buf);CHKERRQ(ierr);
944*2302aa1bSDebojyoti Ghosh     /* Note: print out r as well */
945*2302aa1bSDebojyoti Ghosh   }
946*2302aa1bSDebojyoti Ghosh   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
947*2302aa1bSDebojyoti Ghosh   ierr = TSAdaptView(adapt,viewer);CHKERRQ(ierr);
948*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
949*2302aa1bSDebojyoti Ghosh }
950*2302aa1bSDebojyoti Ghosh 
951*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
952*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSLoad_GLEE"
953*2302aa1bSDebojyoti Ghosh static PetscErrorCode TSLoad_GLEE(TS ts,PetscViewer viewer)
954*2302aa1bSDebojyoti Ghosh {
955*2302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
956*2302aa1bSDebojyoti Ghosh   SNES           snes;
957*2302aa1bSDebojyoti Ghosh   TSAdapt        tsadapt;
958*2302aa1bSDebojyoti Ghosh 
959*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
960*2302aa1bSDebojyoti Ghosh   ierr = TSGetAdapt(ts,&tsadapt);CHKERRQ(ierr);
961*2302aa1bSDebojyoti Ghosh   ierr = TSAdaptLoad(tsadapt,viewer);CHKERRQ(ierr);
962*2302aa1bSDebojyoti Ghosh   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
963*2302aa1bSDebojyoti Ghosh   ierr = SNESLoad(snes,viewer);CHKERRQ(ierr);
964*2302aa1bSDebojyoti Ghosh   /* function and Jacobian context for SNES when used with TS is always ts object */
965*2302aa1bSDebojyoti Ghosh   ierr = SNESSetFunction(snes,NULL,NULL,ts);CHKERRQ(ierr);
966*2302aa1bSDebojyoti Ghosh   ierr = SNESSetJacobian(snes,NULL,NULL,NULL,ts);CHKERRQ(ierr);
967*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
968*2302aa1bSDebojyoti Ghosh }
969*2302aa1bSDebojyoti Ghosh 
970*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
971*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEESetType"
972*2302aa1bSDebojyoti Ghosh /*@C
973*2302aa1bSDebojyoti Ghosh   TSGLEESetType - Set the type of GLEE scheme
974*2302aa1bSDebojyoti Ghosh 
975*2302aa1bSDebojyoti Ghosh   Logically collective
976*2302aa1bSDebojyoti Ghosh 
977*2302aa1bSDebojyoti Ghosh   Input Parameter:
978*2302aa1bSDebojyoti Ghosh +  ts - timestepping context
979*2302aa1bSDebojyoti Ghosh -  gleetype - type of GLEE-scheme
980*2302aa1bSDebojyoti Ghosh 
981*2302aa1bSDebojyoti Ghosh   Level: intermediate
982*2302aa1bSDebojyoti Ghosh 
983*2302aa1bSDebojyoti Ghosh .seealso: TSGLEEGetType(), TSGLEE
984*2302aa1bSDebojyoti Ghosh @*/
985*2302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEESetType(TS ts,TSGLEEType gleetype)
986*2302aa1bSDebojyoti Ghosh {
987*2302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
988*2302aa1bSDebojyoti Ghosh 
989*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
990*2302aa1bSDebojyoti Ghosh   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
991*2302aa1bSDebojyoti Ghosh   ierr = PetscTryMethod(ts,"TSGLEESetType_C",(TS,TSGLEEType),(ts,gleetype));CHKERRQ(ierr);
992*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
993*2302aa1bSDebojyoti Ghosh }
994*2302aa1bSDebojyoti Ghosh 
995*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
996*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEEGetType"
997*2302aa1bSDebojyoti Ghosh /*@C
998*2302aa1bSDebojyoti Ghosh   TSGLEEGetType - Get the type of GLEE scheme
999*2302aa1bSDebojyoti Ghosh 
1000*2302aa1bSDebojyoti Ghosh   Logically collective
1001*2302aa1bSDebojyoti Ghosh 
1002*2302aa1bSDebojyoti Ghosh   Input Parameter:
1003*2302aa1bSDebojyoti Ghosh .  ts - timestepping context
1004*2302aa1bSDebojyoti Ghosh 
1005*2302aa1bSDebojyoti Ghosh   Output Parameter:
1006*2302aa1bSDebojyoti Ghosh .  gleetype - type of GLEE-scheme
1007*2302aa1bSDebojyoti Ghosh 
1008*2302aa1bSDebojyoti Ghosh   Level: intermediate
1009*2302aa1bSDebojyoti Ghosh 
1010*2302aa1bSDebojyoti Ghosh .seealso: TSGLEESetType()
1011*2302aa1bSDebojyoti Ghosh @*/
1012*2302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEEGetType(TS ts,TSGLEEType *gleetype)
1013*2302aa1bSDebojyoti Ghosh {
1014*2302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
1015*2302aa1bSDebojyoti Ghosh 
1016*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
1017*2302aa1bSDebojyoti Ghosh   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1018*2302aa1bSDebojyoti Ghosh   ierr = PetscUseMethod(ts,"TSGLEEGetType_C",(TS,TSGLEEType*),(ts,gleetype));CHKERRQ(ierr);
1019*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
1020*2302aa1bSDebojyoti Ghosh }
1021*2302aa1bSDebojyoti Ghosh 
1022*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
1023*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEEGetType_GLEE"
1024*2302aa1bSDebojyoti Ghosh PetscErrorCode  TSGLEEGetType_GLEE(TS ts,TSGLEEType *gleetype)
1025*2302aa1bSDebojyoti Ghosh {
1026*2302aa1bSDebojyoti Ghosh   TS_GLEE     *glee = (TS_GLEE*)ts->data;
1027*2302aa1bSDebojyoti Ghosh   PetscErrorCode ierr;
1028*2302aa1bSDebojyoti Ghosh 
1029*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
1030*2302aa1bSDebojyoti Ghosh   if (!glee->tableau) {
1031*2302aa1bSDebojyoti Ghosh     ierr = TSGLEESetType(ts,TSGLEEDefaultType);CHKERRQ(ierr);
1032*2302aa1bSDebojyoti Ghosh   }
1033*2302aa1bSDebojyoti Ghosh   *gleetype = glee->tableau->name;
1034*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
1035*2302aa1bSDebojyoti Ghosh }
1036*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
1037*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEESetType_GLEE"
1038*2302aa1bSDebojyoti Ghosh PetscErrorCode  TSGLEESetType_GLEE(TS ts,TSGLEEType gleetype)
1039*2302aa1bSDebojyoti Ghosh {
1040*2302aa1bSDebojyoti Ghosh   TS_GLEE         *glee = (TS_GLEE*)ts->data;
1041*2302aa1bSDebojyoti Ghosh   PetscErrorCode  ierr;
1042*2302aa1bSDebojyoti Ghosh   PetscBool       match;
1043*2302aa1bSDebojyoti Ghosh   GLEETableauLink link;
1044*2302aa1bSDebojyoti Ghosh 
1045*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
1046*2302aa1bSDebojyoti Ghosh   if (glee->tableau) {
1047*2302aa1bSDebojyoti Ghosh     ierr = PetscStrcmp(glee->tableau->name,gleetype,&match);CHKERRQ(ierr);
1048*2302aa1bSDebojyoti Ghosh     if (match) PetscFunctionReturn(0);
1049*2302aa1bSDebojyoti Ghosh   }
1050*2302aa1bSDebojyoti Ghosh   for (link = GLEETableauList; link; link=link->next) {
1051*2302aa1bSDebojyoti Ghosh     ierr = PetscStrcmp(link->tab.name,gleetype,&match);CHKERRQ(ierr);
1052*2302aa1bSDebojyoti Ghosh     if (match) {
1053*2302aa1bSDebojyoti Ghosh       ierr = TSReset_GLEE(ts);CHKERRQ(ierr);
1054*2302aa1bSDebojyoti Ghosh       glee->tableau = &link->tab;
1055*2302aa1bSDebojyoti Ghosh       PetscFunctionReturn(0);
1056*2302aa1bSDebojyoti Ghosh     }
1057*2302aa1bSDebojyoti Ghosh   }
1058*2302aa1bSDebojyoti Ghosh   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",gleetype);
1059*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
1060*2302aa1bSDebojyoti Ghosh }
1061*2302aa1bSDebojyoti Ghosh 
1062*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
1063*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGetStages_GLEE"
1064*2302aa1bSDebojyoti Ghosh static PetscErrorCode  TSGetStages_GLEE(TS ts,PetscInt *ns,Vec **Y)
1065*2302aa1bSDebojyoti Ghosh {
1066*2302aa1bSDebojyoti Ghosh   TS_GLEE *glee = (TS_GLEE*)ts->data;
1067*2302aa1bSDebojyoti Ghosh 
1068*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
1069*2302aa1bSDebojyoti Ghosh   *ns = glee->tableau->s;
1070*2302aa1bSDebojyoti Ghosh   if(Y) *Y  = glee->Y;
1071*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
1072*2302aa1bSDebojyoti Ghosh }
1073*2302aa1bSDebojyoti Ghosh 
1074*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
1075*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGetAuxSolution_GLEE"
1076*2302aa1bSDebojyoti Ghosh PetscErrorCode TSGetAuxSolution_GLEE(TS ts,PetscInt *n,Vec *Y)
1077*2302aa1bSDebojyoti Ghosh {
1078*2302aa1bSDebojyoti Ghosh   TS_GLEE     *glee = (TS_GLEE*)ts->data;
1079*2302aa1bSDebojyoti Ghosh   GLEETableau tab   = glee->tableau;
1080*2302aa1bSDebojyoti Ghosh 
1081*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
1082*2302aa1bSDebojyoti Ghosh   if (!Y) *n = tab->r - 1;
1083*2302aa1bSDebojyoti Ghosh   else {
1084*2302aa1bSDebojyoti Ghosh     if ((*n > 0) && (*n < tab->r)) *Y = glee->Y[*n];
1085*2302aa1bSDebojyoti Ghosh     else {
1086*2302aa1bSDebojyoti Ghosh       /* XXX: Error message for invalid value of n */
1087*2302aa1bSDebojyoti Ghosh     }
1088*2302aa1bSDebojyoti Ghosh   }
1089*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
1090*2302aa1bSDebojyoti Ghosh }
1091*2302aa1bSDebojyoti Ghosh 
1092*2302aa1bSDebojyoti Ghosh /* ------------------------------------------------------------ */
1093*2302aa1bSDebojyoti Ghosh /*MC
1094*2302aa1bSDebojyoti Ghosh       TSGLEE - ODE and DAE solver using General Linear with Error Estimation schemes
1095*2302aa1bSDebojyoti Ghosh 
1096*2302aa1bSDebojyoti Ghosh   The user should provide the right hand side of the equation
1097*2302aa1bSDebojyoti Ghosh   using TSSetRHSFunction().
1098*2302aa1bSDebojyoti Ghosh 
1099*2302aa1bSDebojyoti Ghosh   Notes:
1100*2302aa1bSDebojyoti Ghosh   The default is TSGLEE35, it can be changed with TSGLEESetType() or -ts_glee_type
1101*2302aa1bSDebojyoti Ghosh 
1102*2302aa1bSDebojyoti Ghosh   Level: beginner
1103*2302aa1bSDebojyoti Ghosh 
1104*2302aa1bSDebojyoti Ghosh .seealso:  TSCreate(), TS, TSSetType(), TSGLEESetType(), TSGLEEGetType(),
1105*2302aa1bSDebojyoti Ghosh            TSGLEE23, TTSGLEE24, TSGLEE35, TSGLEE25I, TSGLEEEXRK2A,
1106*2302aa1bSDebojyoti Ghosh            TSGLEERK32G1, TSGLEERK285EX, TSGLEEType, TSGLEERegister()
1107*2302aa1bSDebojyoti Ghosh 
1108*2302aa1bSDebojyoti Ghosh M*/
1109*2302aa1bSDebojyoti Ghosh #undef __FUNCT__
1110*2302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSCreate_GLEE"
1111*2302aa1bSDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSCreate_GLEE(TS ts)
1112*2302aa1bSDebojyoti Ghosh {
1113*2302aa1bSDebojyoti Ghosh   TS_GLEE         *th;
1114*2302aa1bSDebojyoti Ghosh   PetscErrorCode  ierr;
1115*2302aa1bSDebojyoti Ghosh 
1116*2302aa1bSDebojyoti Ghosh   PetscFunctionBegin;
1117*2302aa1bSDebojyoti Ghosh #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
1118*2302aa1bSDebojyoti Ghosh   ierr = TSGLEEInitializePackage();CHKERRQ(ierr);
1119*2302aa1bSDebojyoti Ghosh #endif
1120*2302aa1bSDebojyoti Ghosh 
1121*2302aa1bSDebojyoti Ghosh   ts->ops->reset          = TSReset_GLEE;
1122*2302aa1bSDebojyoti Ghosh   ts->ops->destroy        = TSDestroy_GLEE;
1123*2302aa1bSDebojyoti Ghosh   ts->ops->view           = TSView_GLEE;
1124*2302aa1bSDebojyoti Ghosh   ts->ops->load           = TSLoad_GLEE;
1125*2302aa1bSDebojyoti Ghosh   ts->ops->setup          = TSSetUp_GLEE;
1126*2302aa1bSDebojyoti Ghosh   ts->ops->step           = TSStep_GLEE;
1127*2302aa1bSDebojyoti Ghosh   ts->ops->interpolate    = TSInterpolate_GLEE;
1128*2302aa1bSDebojyoti Ghosh   ts->ops->evaluatestep   = TSEvaluateStep_GLEE;
1129*2302aa1bSDebojyoti Ghosh   ts->ops->setfromoptions = TSSetFromOptions_GLEE;
1130*2302aa1bSDebojyoti Ghosh   ts->ops->getstages      = TSGetStages_GLEE;
1131*2302aa1bSDebojyoti Ghosh   ts->ops->snesfunction   = SNESTSFormFunction_GLEE;
1132*2302aa1bSDebojyoti Ghosh   ts->ops->snesjacobian   = SNESTSFormJacobian_GLEE;
1133*2302aa1bSDebojyoti Ghosh   ts->ops->getauxsolution = TSGetAuxSolution_GLEE;
1134*2302aa1bSDebojyoti Ghosh   ts->ops->startingmethod = TSStartingMethod_GLEE;
1135*2302aa1bSDebojyoti Ghosh 
1136*2302aa1bSDebojyoti Ghosh   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
1137*2302aa1bSDebojyoti Ghosh   ts->data = (void*)th;
1138*2302aa1bSDebojyoti Ghosh 
1139*2302aa1bSDebojyoti Ghosh   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLEEGetType_C",TSGLEEGetType_GLEE);CHKERRQ(ierr);
1140*2302aa1bSDebojyoti Ghosh   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLEESetType_C",TSGLEESetType_GLEE);CHKERRQ(ierr);
1141*2302aa1bSDebojyoti Ghosh   PetscFunctionReturn(0);
1142*2302aa1bSDebojyoti Ghosh }
1143