Actual source code: arkimex.c
1: /*
2: Code for timestepping with additive Runge-Kutta IMEX method
4: Notes:
5: The general system is written as
7: F(t,X,Xdot) = G(t,X)
9: where F represents the stiff part of the physics and G represents the non-stiff part.
11: */
12: #include <private/tsimpl.h> /*I "petscts.h" I*/
14: static const TSARKIMEXType TSARKIMEXDefault = TSARKIMEX2E;
15: static PetscBool TSARKIMEXRegisterAllCalled;
16: static PetscBool TSARKIMEXPackageInitialized;
18: typedef struct _ARKTableau *ARKTableau;
19: struct _ARKTableau {
20: char *name;
21: PetscInt order; /* Classical approximation order of the method */
22: PetscInt s; /* Number of stages */
23: PetscInt pinterp; /* Interpolation order */
24: PetscReal *At,*bt,*ct; /* Stiff tableau */
25: PetscReal *A,*b,*c; /* Non-stiff tableau */
26: PetscReal *binterpt,*binterp; /* Dense output formula */
27: };
28: typedef struct _ARKTableauLink *ARKTableauLink;
29: struct _ARKTableauLink {
30: struct _ARKTableau tab;
31: ARKTableauLink next;
32: };
33: static ARKTableauLink ARKTableauList;
35: typedef struct {
36: ARKTableau tableau;
37: Vec *Y; /* States computed during the step */
38: Vec *YdotI; /* Time derivatives for the stiff part */
39: Vec *YdotRHS; /* Function evaluations for the non-stiff part */
40: Vec Ydot; /* Work vector holding Ydot during residual evaluation */
41: Vec Work; /* Generic work vector */
42: Vec Z; /* Ydot = shift(Y-Z) */
43: PetscScalar *work; /* Scalar work */
44: PetscReal shift;
45: PetscReal stage_time;
46: PetscBool imex;
47: } TS_ARKIMEX;
49: /*MC
50: TSARKIMEX2D - Second order ARK IMEX scheme with L-stable implicit part.
52: This method has one explicit stage and two implicit stages.
54: .seealso: TSARKIMEX
55: M*/
56: /*MC
57: TSARKIMEX2E - Second order ARK IMEX scheme with L-stable implicit part.
59: This method has one explicit stage and two implicit stages. It is is an optimal method developed by Emil Constantinescu.
61: .seealso: TSARKIMEX
62: M*/
63: /*MC
64: TSARKIMEX3 - Third order ARK IMEX scheme with L-stable implicit part.
66: This method has one explicit stage and three implicit stages.
68: References:
69: Kennedy and Carpenter 2003.
71: .seealso: TSARKIMEX
72: M*/
73: /*MC
74: TSARKIMEX4 - Fourth order ARK IMEX scheme with L-stable implicit part.
76: This method has one explicit stage and four implicit stages.
78: References:
79: Kennedy and Carpenter 2003.
81: .seealso: TSARKIMEX
82: M*/
83: /*MC
84: TSARKIMEX5 - Fifth order ARK IMEX scheme with L-stable implicit part.
86: This method has one explicit stage and five implicit stages.
88: References:
89: Kennedy and Carpenter 2003.
91: .seealso: TSARKIMEX
92: M*/
96: /*@C
97: TSARKIMEXRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSARKIMEX
99: Not Collective, but should be called by all processes which will need the schemes to be registered
101: Level: advanced
103: .keywords: TS, TSARKIMEX, register, all
105: .seealso: TSARKIMEXRegisterDestroy()
106: @*/
107: PetscErrorCode TSARKIMEXRegisterAll(void)
108: {
112: if (TSARKIMEXRegisterAllCalled) return(0);
113: TSARKIMEXRegisterAllCalled = PETSC_TRUE;
115: {
116: const PetscReal
117: A[3][3] = {{0,0,0},
118: {0.41421356237309504880,0,0},
119: {0.75,0.25,0}},
120: At[3][3] = {{0,0,0},
121: {0.12132034355964257320,0.29289321881345247560,0},
122: {0.20710678118654752440,0.50000000000000000000,0.29289321881345247560}},
123: binterpt[3][2] = {{1,-0.5},{0,0},{0,0.5}};
124: TSARKIMEXRegister(TSARKIMEX2D,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,2,binterpt[0],PETSC_NULL);
125: }
126: { /* Optimal for linear implicit part */
127: const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),
128: A[3][3] = {{0,0,0},
129: {2-s2,0,0},
130: {(3-2*s2)/6,(3+2*s2)/6,0}},
131: At[3][3] = {{0,0,0},
132: {1-1/s2,1-1/s2,0},
133: {1/(2*s2),1/(2*s2),1-1/s2}},
134: binterpt[3][2] = {{1,-0.5},{0,0},{0,0.5}};
135: TSARKIMEXRegister(TSARKIMEX2E,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,2,binterpt[0],PETSC_NULL);
136: }
137: {
138: const PetscReal
139: A[4][4] = {{0,0,0,0},
140: {1767732205903./2027836641118.,0,0,0},
141: {5535828885825./10492691773637.,788022342437./10882634858940.,0,0},
142: {6485989280629./16251701735622.,-4246266847089./9704473918619.,10755448449292./10357097424841.,0}},
143: At[4][4] = {{0,0,0,0},
144: {1767732205903./4055673282236.,1767732205903./4055673282236.,0,0},
145: {2746238789719./10658868560708.,-640167445237./6845629431997.,1767732205903./4055673282236.,0},
146: {1471266399579./7840856788654.,-4482444167858./7529755066697.,11266239266428./11593286722821.,1767732205903./4055673282236.}},
147: binterpt[4][2] = {{4655552711362./22874653954995., -215264564351./13552729205753.},
148: {-18682724506714./9892148508045.,17870216137069./13817060693119.},
149: {34259539580243./13192909600954.,-28141676662227./17317692491321.},
150: {584795268549./6622622206610., 2508943948391./7218656332882.}};
151: TSARKIMEXRegister(TSARKIMEX3,3,4,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,2,binterpt[0],PETSC_NULL);
152: }
153: {
154: const PetscReal
155: A[6][6] = {{0,0,0,0,0,0},
156: {1./2,0,0,0,0,0},
157: {13861./62500.,6889./62500.,0,0,0,0},
158: {-116923316275./2393684061468.,-2731218467317./15368042101831.,9408046702089./11113171139209.,0,0,0},
159: {-451086348788./2902428689909.,-2682348792572./7519795681897.,12662868775082./11960479115383.,3355817975965./11060851509271.,0,0},
160: {647845179188./3216320057751.,73281519250./8382639484533.,552539513391./3454668386233.,3354512671639./8306763924573.,4040./17871.,0}},
161: At[6][6] = {{0,0,0,0,0,0},
162: {1./4,1./4,0,0,0,0},
163: {8611./62500.,-1743./31250.,1./4,0,0,0},
164: {5012029./34652500.,-654441./2922500.,174375./388108.,1./4,0,0},
165: {15267082809./155376265600.,-71443401./120774400.,730878875./902184768.,2285395./8070912.,1./4,0},
166: {82889./524892.,0,15625./83664.,69875./102672.,-2260./8211,1./4}},
167: binterpt[6][3] = {{6943876665148./7220017795957.,-54480133./30881146.,6818779379841./7100303317025.},
168: {0,0,0},
169: {7640104374378./9702883013639.,-11436875./14766696.,2173542590792./12501825683035.},
170: {-20649996744609./7521556579894.,174696575./18121608.,-31592104683404./5083833661969.},
171: {8854892464581./2390941311638.,-12120380./966161.,61146701046299./7138195549469.},
172: {-11397109935349./6675773540249.,3843./706.,-17219254887155./4939391667607.}};
173: TSARKIMEXRegister(TSARKIMEX4,4,6,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,3,binterpt[0],PETSC_NULL);
174: }
175: {
176: const PetscReal
177: A[8][8] = {{0,0,0,0,0,0,0,0},
178: {41./100,0,0,0,0,0,0,0},
179: {367902744464./2072280473677.,677623207551./8224143866563.,0,0,0,0,0,0},
180: {1268023523408./10340822734521.,0,1029933939417./13636558850479.,0,0,0,0,0},
181: {14463281900351./6315353703477.,0,66114435211212./5879490589093.,-54053170152839./4284798021562.,0,0,0,0},
182: {14090043504691./34967701212078.,0,15191511035443./11219624916014.,-18461159152457./12425892160975.,-281667163811./9011619295870.,0,0,0},
183: {19230459214898./13134317526959.,0,21275331358303./2942455364971.,-38145345988419./4862620318723.,-1./8,-1./8,0,0},
184: {-19977161125411./11928030595625.,0,-40795976796054./6384907823539.,177454434618887./12078138498510.,782672205425./8267701900261.,-69563011059811./9646580694205.,7356628210526./4942186776405.,0}},
185: At[8][8] = {{0,0,0,0,0,0,0,0},
186: {41./200.,41./200.,0,0,0,0,0,0},
187: {41./400.,-567603406766./11931857230679.,41./200.,0,0,0,0,0},
188: {683785636431./9252920307686.,0,-110385047103./1367015193373.,41./200.,0,0,0,0},
189: {3016520224154./10081342136671.,0,30586259806659./12414158314087.,-22760509404356./11113319521817.,41./200.,0,0,0},
190: {218866479029./1489978393911.,0,638256894668./5436446318841.,-1179710474555./5321154724896.,-60928119172./8023461067671.,41./200.,0,0},
191: {1020004230633./5715676835656.,0,25762820946817./25263940353407.,-2161375909145./9755907335909.,-211217309593./5846859502534.,-4269925059573./7827059040749.,41./200,0},
192: {-872700587467./9133579230613.,0,0,22348218063261./9555858737531.,-1143369518992./8141816002931.,-39379526789629./19018526304540.,32727382324388./42900044865799.,41./200.}},
193: binterpt[8][3] = {{-17674230611817./10670229744614. , 43486358583215./12773830924787. , -9257016797708./5021505065439.},
194: {0 , 0 , 0 },
195: {0 , 0 , 0 },
196: {65168852399939./7868540260826. , -91478233927265./11067650958493., 26096422576131./11239449250142.},
197: {15494834004392./5936557850923. , -79368583304911./10890268929626., 92396832856987./20362823103730.},
198: {-99329723586156./26959484932159., -12239297817655./9152339842473. , 30029262896817./10175596800299.},
199: {-19024464361622./5461577185407. , 115839755401235./10719374521269., -26136350496073./3983972220547.},
200: {-6511271360970./6095937251113. , 5843115559534./2180450260947. , -5289405421727./3760307252460. }};
201: TSARKIMEXRegister(TSARKIMEX5,5,8,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,3,binterpt[0],PETSC_NULL);
202: }
204: return(0);
205: }
209: /*@C
210: TSARKIMEXRegisterDestroy - Frees the list of schemes that were registered by TSARKIMEXRegister().
212: Not Collective
214: Level: advanced
216: .keywords: TSARKIMEX, register, destroy
217: .seealso: TSARKIMEXRegister(), TSARKIMEXRegisterAll(), TSARKIMEXRegisterDynamic()
218: @*/
219: PetscErrorCode TSARKIMEXRegisterDestroy(void)
220: {
222: ARKTableauLink link;
225: while ((link = ARKTableauList)) {
226: ARKTableau t = &link->tab;
227: ARKTableauList = link->next;
228: PetscFree6(t->At,t->bt,t->ct,t->A,t->b,t->c);
229: PetscFree2(t->binterpt,t->binterp);
230: PetscFree(t->name);
231: PetscFree(link);
232: }
233: TSARKIMEXRegisterAllCalled = PETSC_FALSE;
234: return(0);
235: }
239: /*@C
240: TSARKIMEXInitializePackage - This function initializes everything in the TSARKIMEX package. It is called
241: from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_ARKIMEX()
242: when using static libraries.
244: Input Parameter:
245: path - The dynamic library path, or PETSC_NULL
247: Level: developer
249: .keywords: TS, TSARKIMEX, initialize, package
250: .seealso: PetscInitialize()
251: @*/
252: PetscErrorCode TSARKIMEXInitializePackage(const char path[])
253: {
257: if (TSARKIMEXPackageInitialized) return(0);
258: TSARKIMEXPackageInitialized = PETSC_TRUE;
259: TSARKIMEXRegisterAll();
260: PetscRegisterFinalize(TSARKIMEXFinalizePackage);
261: return(0);
262: }
266: /*@C
267: TSARKIMEXFinalizePackage - This function destroys everything in the TSARKIMEX package. It is
268: called from PetscFinalize().
270: Level: developer
272: .keywords: Petsc, destroy, package
273: .seealso: PetscFinalize()
274: @*/
275: PetscErrorCode TSARKIMEXFinalizePackage(void)
276: {
280: TSARKIMEXPackageInitialized = PETSC_FALSE;
281: TSARKIMEXRegisterDestroy();
282: return(0);
283: }
287: /*@C
288: TSARKIMEXRegister - register an ARK IMEX scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
290: Not Collective, but the same schemes should be registered on all processes on which they will be used
292: Input Parameters:
293: + name - identifier for method
294: . order - approximation order of method
295: . s - number of stages, this is the dimension of the matrices below
296: . At - Butcher table of stage coefficients for stiff part (dimension s*s, row-major)
297: . bt - Butcher table for completing the stiff part of the step (dimension s; PETSC_NULL to use the last row of At)
298: . ct - Abscissa of each stiff stage (dimension s, PETSC_NULL to use row sums of At)
299: . A - Non-stiff stage coefficients (dimension s*s, row-major)
300: . b - Non-stiff step completion table (dimension s; PETSC_NULL to use last row of At)
301: . c - Non-stiff abscissa (dimension s; PETSC_NULL to use row sums of A)
302: . pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt and binterp
303: . binterpt - Coefficients of the interpolation formula for the stiff part (dimension s*pinterp)
304: - binterp - Coefficients of the interpolation formula for the non-stiff part (dimension s*pinterp; PETSC_NULL to reuse binterpt)
306: Notes:
307: Several ARK IMEX methods are provided, this function is only needed to create new methods.
309: Level: advanced
311: .keywords: TS, register
313: .seealso: TSARKIMEX
314: @*/
315: PetscErrorCode TSARKIMEXRegister(const TSARKIMEXType name,PetscInt order,PetscInt s,
316: const PetscReal At[],const PetscReal bt[],const PetscReal ct[],
317: const PetscReal A[],const PetscReal b[],const PetscReal c[],
318: PetscInt pinterp,const PetscReal binterpt[],const PetscReal binterp[])
319: {
321: ARKTableauLink link;
322: ARKTableau t;
323: PetscInt i,j;
326: PetscMalloc(sizeof(*link),&link);
327: PetscMemzero(link,sizeof(*link));
328: t = &link->tab;
329: PetscStrallocpy(name,&t->name);
330: t->order = order;
331: t->s = s;
332: PetscMalloc6(s*s,PetscReal,&t->At,s,PetscReal,&t->bt,s,PetscReal,&t->ct,s*s,PetscReal,&t->A,s,PetscReal,&t->b,s,PetscReal,&t->c);
333: PetscMemcpy(t->At,At,s*s*sizeof(At[0]));
334: PetscMemcpy(t->A,A,s*s*sizeof(A[0]));
335: if (bt) {PetscMemcpy(t->bt,bt,s*sizeof(bt[0]));}
336: else for (i=0; i<s; i++) t->bt[i] = At[(s-1)*s+i];
337: if (b) {PetscMemcpy(t->b,b,s*sizeof(b[0]));}
338: else for (i=0; i<s; i++) t->b[i] = At[(s-1)*s+i];
339: if (ct) {PetscMemcpy(t->ct,ct,s*sizeof(ct[0]));}
340: else for (i=0; i<s; i++) for (j=0,t->ct[i]=0; j<s; j++) t->ct[i] += At[i*s+j];
341: if (c) {PetscMemcpy(t->c,c,s*sizeof(c[0]));}
342: else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j];
343: t->pinterp = pinterp;
344: PetscMalloc2(s*pinterp,PetscReal,&t->binterpt,s*pinterp,PetscReal,&t->binterp);
345: PetscMemcpy(t->binterpt,binterpt,s*pinterp*sizeof(binterpt[0]));
346: PetscMemcpy(t->binterp,binterp?binterp:binterpt,s*pinterp*sizeof(binterpt[0]));
347: link->next = ARKTableauList;
348: ARKTableauList = link;
349: return(0);
350: }
354: static PetscErrorCode TSStep_ARKIMEX(TS ts)
355: {
356: TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
357: ARKTableau tab = ark->tableau;
358: const PetscInt s = tab->s;
359: const PetscReal *At = tab->At,*A = tab->A,*bt = tab->bt,*b = tab->b,*ct = tab->ct,*c = tab->c;
360: PetscScalar *w = ark->work;
361: Vec *Y = ark->Y,*YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS,Ydot = ark->Ydot,W = ark->Work,Z = ark->Z;
362: SNES snes;
363: PetscInt i,j,its,lits;
364: PetscReal next_time_step;
365: PetscReal h,t;
366: PetscErrorCode ierr;
369: TSGetSNES(ts,&snes);
370: next_time_step = ts->time_step;
371: h = ts->time_step;
372: t = ts->ptime;
374: for (i=0; i<s; i++) {
375: if (At[i*s+i] == 0) { /* This stage is explicit */
376: VecCopy(ts->vec_sol,Y[i]);
377: for (j=0; j<i; j++) w[j] = -h*At[i*s+j];
378: VecMAXPY(Y[i],i,w,YdotI);
379: for (j=0; j<i; j++) w[j] = h*A[i*s+j];
380: VecMAXPY(Y[i],i,w,YdotRHS);
381: } else {
382: ark->stage_time = t + h*ct[i];
383: ark->shift = 1./(h*At[i*s+i]);
384: /* Affine part */
385: VecZeroEntries(W);
386: for (j=0; j<i; j++) w[j] = h*A[i*s+j];
387: VecMAXPY(W,i,w,YdotRHS);
388: /* Ydot = shift*(Y-Z) */
389: VecCopy(ts->vec_sol,Z);
390: for (j=0; j<i; j++) w[j] = -h*At[i*s+j];
391: VecMAXPY(Z,i,w,YdotI);
392: /* Initial guess taken from last stage */
393: VecCopy(i>0?Y[i-1]:ts->vec_sol,Y[i]);
394: SNESSolve(snes,W,Y[i]);
395: SNESGetIterationNumber(snes,&its);
396: SNESGetLinearSolveIterations(snes,&lits);
397: ts->nonlinear_its += its; ts->linear_its += lits;
398: }
399: VecZeroEntries(Ydot);
400: TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],ark->imex);
401: if (ark->imex) {
402: TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);
403: } else {
404: VecZeroEntries(YdotRHS[i]);
405: }
406: }
407: for (j=0; j<s; j++) w[j] = -h*bt[j];
408: VecMAXPY(ts->vec_sol,s,w,YdotI);
409: for (j=0; j<s; j++) w[j] = h*b[j];
410: VecMAXPY(ts->vec_sol,s,w,YdotRHS);
412: ts->ptime += ts->time_step;
413: ts->time_step = next_time_step;
414: ts->steps++;
415: return(0);
416: }
420: static PetscErrorCode TSInterpolate_ARKIMEX(TS ts,PetscReal itime,Vec X)
421: {
422: TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
423: PetscInt s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j;
424: PetscReal tt,t = (itime - ts->ptime)/ts->time_step + 1; /* In the interval [0,1] */
425: PetscScalar *bt,*b;
426: const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp;
430: if (!Bt || !B) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name);
431: PetscMalloc2(s,PetscScalar,&bt,s,PetscScalar,&b);
432: for (i=0; i<s; i++) bt[i] = b[i] = 0;
433: for (j=0,tt=t; j<pinterp; j++,tt*=t) {
434: for (i=0; i<s; i++) {
435: bt[i] += ts->time_step * Bt[i*pinterp+j] * tt * -1.0;
436: b[i] += ts->time_step * B[i*pinterp+j] * tt;
437: }
438: }
439: if (ark->tableau->At[0*s+0] != 0.0) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_SUP,"First stage not explicit so starting stage not saved");
440: VecCopy(ark->Y[0],X);
441: VecMAXPY(X,s,bt,ark->YdotI);
442: VecMAXPY(X,s,b,ark->YdotRHS);
443: PetscFree2(bt,b);
444: return(0);
445: }
447: /*------------------------------------------------------------*/
450: static PetscErrorCode TSReset_ARKIMEX(TS ts)
451: {
452: TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
453: PetscInt s;
454: PetscErrorCode ierr;
457: if (!ark->tableau) return(0);
458: s = ark->tableau->s;
459: VecDestroyVecs(s,&ark->Y);
460: VecDestroyVecs(s,&ark->YdotI);
461: VecDestroyVecs(s,&ark->YdotRHS);
462: VecDestroy(&ark->Ydot);
463: VecDestroy(&ark->Work);
464: VecDestroy(&ark->Z);
465: PetscFree(ark->work);
466: return(0);
467: }
471: static PetscErrorCode TSDestroy_ARKIMEX(TS ts)
472: {
473: PetscErrorCode ierr;
476: TSReset_ARKIMEX(ts);
477: PetscFree(ts->data);
478: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXGetType_C","",PETSC_NULL);
479: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetType_C","",PETSC_NULL);
480: return(0);
481: }
483: /*
484: This defines the nonlinear equation that is to be solved with SNES
485: G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
486: */
489: static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec F,TS ts)
490: {
491: TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
495: VecAXPBYPCZ(ark->Ydot,-ark->shift,ark->shift,0,ark->Z,X); /* Ydot = shift*(X-Z) */
496: TSComputeIFunction(ts,ark->stage_time,X,ark->Ydot,F,ark->imex);
497: return(0);
498: }
502: static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
503: {
504: TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
508: /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */
509: TSComputeIJacobian(ts,ark->stage_time,X,ark->Ydot,ark->shift,A,B,str,PETSC_TRUE);
510: return(0);
511: }
515: static PetscErrorCode TSSetUp_ARKIMEX(TS ts)
516: {
517: TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
518: ARKTableau tab = ark->tableau;
519: PetscInt s = tab->s;
523: if (!ark->tableau) {
524: TSARKIMEXSetType(ts,TSARKIMEXDefault);
525: }
526: VecDuplicateVecs(ts->vec_sol,s,&ark->Y);
527: VecDuplicateVecs(ts->vec_sol,s,&ark->YdotI);
528: VecDuplicateVecs(ts->vec_sol,s,&ark->YdotRHS);
529: VecDuplicate(ts->vec_sol,&ark->Ydot);
530: VecDuplicate(ts->vec_sol,&ark->Work);
531: VecDuplicate(ts->vec_sol,&ark->Z);
532: PetscMalloc(s*sizeof(ark->work[0]),&ark->work);
533: return(0);
534: }
535: /*------------------------------------------------------------*/
539: static PetscErrorCode TSSetFromOptions_ARKIMEX(TS ts)
540: {
541: TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
543: char arktype[256];
546: PetscOptionsHead("ARKIMEX ODE solver options");
547: {
548: ARKTableauLink link;
549: PetscInt count,choice;
550: PetscBool flg;
551: const char **namelist;
552: PetscStrncpy(arktype,TSARKIMEXDefault,sizeof arktype);
553: for (link=ARKTableauList,count=0; link; link=link->next,count++) ;
554: PetscMalloc(count*sizeof(char*),&namelist);
555: for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
556: PetscOptionsEList("-ts_arkimex_type","Family of ARK IMEX method","TSARKIMEXSetType",(const char*const*)namelist,count,arktype,&choice,&flg);
557: TSARKIMEXSetType(ts,flg ? namelist[choice] : arktype);
558: PetscFree(namelist);
559: flg = (PetscBool)!ark->imex;
560: PetscOptionsBool("-ts_arkimex_fully_implicit","Solve the problem fully implicitly","TSARKIMEXSetFullyImplicit",flg,&flg,PETSC_NULL);
561: ark->imex = (PetscBool)!flg;
562: SNESSetFromOptions(ts->snes);
563: }
564: PetscOptionsTail();
565: return(0);
566: }
570: static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
571: {
573: int i,left,count;
574: char *p;
577: for (i=0,p=buf,left=(int)len; i<n; i++) {
578: PetscSNPrintf(p,left,fmt,&count,x[i]);
579: if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
580: left -= count;
581: p += count;
582: *p++ = ' ';
583: }
584: p[i ? 0 : -1] = 0;
585: return(0);
586: }
590: static PetscErrorCode TSView_ARKIMEX(TS ts,PetscViewer viewer)
591: {
592: TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
593: ARKTableau tab = ark->tableau;
594: PetscBool iascii;
598: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
599: if (iascii) {
600: const TSARKIMEXType arktype;
601: char buf[512];
602: TSARKIMEXGetType(ts,&arktype);
603: PetscViewerASCIIPrintf(viewer," ARK IMEX %s\n",arktype);
604: PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ct);
605: PetscViewerASCIIPrintf(viewer," Stiff abscissa ct = %s\n",buf);
606: PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->c);
607: PetscViewerASCIIPrintf(viewer," Nonstiff abscissa c = %s\n",buf);
608: }
609: SNESView(ts->snes,viewer);
610: return(0);
611: }
615: /*@C
616: TSARKIMEXSetType - Set the type of ARK IMEX scheme
618: Logically collective
620: Input Parameter:
621: + ts - timestepping context
622: - arktype - type of ARK-IMEX scheme
624: Level: intermediate
626: .seealso: TSARKIMEXGetType()
627: @*/
628: PetscErrorCode TSARKIMEXSetType(TS ts,const TSARKIMEXType arktype)
629: {
634: PetscTryMethod(ts,"TSARKIMEXSetType_C",(TS,const TSARKIMEXType),(ts,arktype));
635: return(0);
636: }
640: /*@C
641: TSARKIMEXGetType - Get the type of ARK IMEX scheme
643: Logically collective
645: Input Parameter:
646: . ts - timestepping context
648: Output Parameter:
649: . arktype - type of ARK-IMEX scheme
651: Level: intermediate
653: .seealso: TSARKIMEXGetType()
654: @*/
655: PetscErrorCode TSARKIMEXGetType(TS ts,const TSARKIMEXType *arktype)
656: {
661: PetscUseMethod(ts,"TSARKIMEXGetType_C",(TS,const TSARKIMEXType*),(ts,arktype));
662: return(0);
663: }
667: /*@C
668: TSARKIMEXSetFullyImplicit - Solve both parts of the equation implicitly
670: Logically collective
672: Input Parameter:
673: + ts - timestepping context
674: - flg - PETSC_TRUE for fully implicit
676: Level: intermediate
678: .seealso: TSARKIMEXGetType()
679: @*/
680: PetscErrorCode TSARKIMEXSetFullyImplicit(TS ts,PetscBool flg)
681: {
686: PetscTryMethod(ts,"TSARKIMEXSetFullyImplicit_C",(TS,PetscBool),(ts,flg));
687: return(0);
688: }
693: PetscErrorCode TSARKIMEXGetType_ARKIMEX(TS ts,const TSARKIMEXType *arktype)
694: {
695: TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
699: if (!ark->tableau) {TSARKIMEXSetType(ts,TSARKIMEXDefault);}
700: *arktype = ark->tableau->name;
701: return(0);
702: }
705: PetscErrorCode TSARKIMEXSetType_ARKIMEX(TS ts,const TSARKIMEXType arktype)
706: {
707: TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
709: PetscBool match;
710: ARKTableauLink link;
713: if (ark->tableau) {
714: PetscStrcmp(ark->tableau->name,arktype,&match);
715: if (match) return(0);
716: }
717: for (link = ARKTableauList; link; link=link->next) {
718: PetscStrcmp(link->tab.name,arktype,&match);
719: if (match) {
720: TSReset_ARKIMEX(ts);
721: ark->tableau = &link->tab;
722: return(0);
723: }
724: }
725: SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype);
726: return(0);
727: }
730: PetscErrorCode TSARKIMEXSetFullyImplicit_ARKIMEX(TS ts,PetscBool flg)
731: {
732: TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
735: ark->imex = (PetscBool)!flg;
736: return(0);
737: }
740: /* ------------------------------------------------------------ */
741: /*MC
742: TSARKIMEX - ODE solver using Additive Runge-Kutta IMEX schemes
744: These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
745: nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
746: of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
748: Notes:
749: The default is TSARKIMEX2E, it can be changed with TSARKIMEXSetType() or -ts_arkimex_type
751: This method currently only works with ODE, for which the stiff part G(t,X,Xdot) has the form Xdot + Ghat(t,X).
753: Level: beginner
755: .seealso: TSCreate(), TS, TSSetType(), TSARKIMEXSetType(), TSARKIMEXGetType(), TSARKIMEXSetFullyImplicit(), TSARKIMEX2D, TTSARKIMEX2E, TSARKIMEX3,
756: TSARKIMEX4, TSARKIMEX5, TSARKIMEXType, SARKIMEXRegister()
758: M*/
762: PetscErrorCode TSCreate_ARKIMEX(TS ts)
763: {
764: TS_ARKIMEX *th;
768: #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
769: TSARKIMEXInitializePackage(PETSC_NULL);
770: #endif
772: ts->ops->reset = TSReset_ARKIMEX;
773: ts->ops->destroy = TSDestroy_ARKIMEX;
774: ts->ops->view = TSView_ARKIMEX;
775: ts->ops->setup = TSSetUp_ARKIMEX;
776: ts->ops->step = TSStep_ARKIMEX;
777: ts->ops->interpolate = TSInterpolate_ARKIMEX;
778: ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX;
779: ts->ops->snesfunction = SNESTSFormFunction_ARKIMEX;
780: ts->ops->snesjacobian = SNESTSFormJacobian_ARKIMEX;
782: PetscNewLog(ts,TS_ARKIMEX,&th);
783: ts->data = (void*)th;
784: th->imex = PETSC_TRUE;
786: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXGetType_C","TSARKIMEXGetType_ARKIMEX",TSARKIMEXGetType_ARKIMEX);
787: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetType_C","TSARKIMEXSetType_ARKIMEX",TSARKIMEXSetType_ARKIMEX);
788: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C","TSARKIMEXSetFullyImplicit_ARKIMEX",TSARKIMEXSetFullyImplicit_ARKIMEX);
789: return(0);
790: }