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