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: Level: advanced
56: .seealso: TSARKIMEX
57: M*/
58: /*MC
59: TSARKIMEX2E - Second order ARK IMEX scheme with L-stable implicit part.
61: This method has one explicit stage and two implicit stages. It is is an optimal method developed by Emil Constantinescu.
63: Level: advanced
65: .seealso: TSARKIMEX
66: M*/
67: /*MC
68: TSARKIMEX3 - Third order ARK IMEX scheme with L-stable implicit part.
70: This method has one explicit stage and three implicit stages.
72: References:
73: Kennedy and Carpenter 2003.
75: Level: advanced
77: .seealso: TSARKIMEX
78: M*/
79: /*MC
80: TSARKIMEX4 - Fourth order ARK IMEX scheme with L-stable implicit part.
82: This method has one explicit stage and four implicit stages.
84: References:
85: Kennedy and Carpenter 2003.
87: Level: advanced
89: .seealso: TSARKIMEX
90: M*/
91: /*MC
92: TSARKIMEX5 - Fifth order ARK IMEX scheme with L-stable implicit part.
94: This method has one explicit stage and five implicit stages.
96: References:
97: Kennedy and Carpenter 2003.
99: Level: advanced
101: .seealso: TSARKIMEX
102: M*/
106: /*@C
107: TSARKIMEXRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSARKIMEX
109: Not Collective, but should be called by all processes which will need the schemes to be registered
111: Level: advanced
113: .keywords: TS, TSARKIMEX, register, all
115: .seealso: TSARKIMEXRegisterDestroy()
116: @*/
117: PetscErrorCode TSARKIMEXRegisterAll(void)
118: {
122: if (TSARKIMEXRegisterAllCalled) return(0);
123: TSARKIMEXRegisterAllCalled = PETSC_TRUE;
125: {
126: const PetscReal
127: A[3][3] = {{0,0,0},
128: {0.41421356237309504880,0,0},
129: {0.75,0.25,0}},
130: At[3][3] = {{0,0,0},
131: {0.12132034355964257320,0.29289321881345247560,0},
132: {0.20710678118654752440,0.50000000000000000000,0.29289321881345247560}},
133: binterpt[3][2] = {{1,-0.5},{0,0},{0,0.5}};
134: TSARKIMEXRegister(TSARKIMEX2D,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,2,binterpt[0],PETSC_NULL);
135: }
136: { /* Optimal for linear implicit part */
137: const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),
138: A[3][3] = {{0,0,0},
139: {2-s2,0,0},
140: {(3-2*s2)/6,(3+2*s2)/6,0}},
141: At[3][3] = {{0,0,0},
142: {1-1/s2,1-1/s2,0},
143: {1/(2*s2),1/(2*s2),1-1/s2}},
144: binterpt[3][2] = {{1,-0.5},{0,0},{0,0.5}};
145: TSARKIMEXRegister(TSARKIMEX2E,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,2,binterpt[0],PETSC_NULL);
146: }
147: {
148: const PetscReal
149: A[4][4] = {{0,0,0,0},
150: {1767732205903./2027836641118.,0,0,0},
151: {5535828885825./10492691773637.,788022342437./10882634858940.,0,0},
152: {6485989280629./16251701735622.,-4246266847089./9704473918619.,10755448449292./10357097424841.,0}},
153: At[4][4] = {{0,0,0,0},
154: {1767732205903./4055673282236.,1767732205903./4055673282236.,0,0},
155: {2746238789719./10658868560708.,-640167445237./6845629431997.,1767732205903./4055673282236.,0},
156: {1471266399579./7840856788654.,-4482444167858./7529755066697.,11266239266428./11593286722821.,1767732205903./4055673282236.}},
157: binterpt[4][2] = {{4655552711362./22874653954995., -215264564351./13552729205753.},
158: {-18682724506714./9892148508045.,17870216137069./13817060693119.},
159: {34259539580243./13192909600954.,-28141676662227./17317692491321.},
160: {584795268549./6622622206610., 2508943948391./7218656332882.}};
161: TSARKIMEXRegister(TSARKIMEX3,3,4,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,2,binterpt[0],PETSC_NULL);
162: }
163: {
164: const PetscReal
165: A[6][6] = {{0,0,0,0,0,0},
166: {1./2,0,0,0,0,0},
167: {13861./62500.,6889./62500.,0,0,0,0},
168: {-116923316275./2393684061468.,-2731218467317./15368042101831.,9408046702089./11113171139209.,0,0,0},
169: {-451086348788./2902428689909.,-2682348792572./7519795681897.,12662868775082./11960479115383.,3355817975965./11060851509271.,0,0},
170: {647845179188./3216320057751.,73281519250./8382639484533.,552539513391./3454668386233.,3354512671639./8306763924573.,4040./17871.,0}},
171: At[6][6] = {{0,0,0,0,0,0},
172: {1./4,1./4,0,0,0,0},
173: {8611./62500.,-1743./31250.,1./4,0,0,0},
174: {5012029./34652500.,-654441./2922500.,174375./388108.,1./4,0,0},
175: {15267082809./155376265600.,-71443401./120774400.,730878875./902184768.,2285395./8070912.,1./4,0},
176: {82889./524892.,0,15625./83664.,69875./102672.,-2260./8211,1./4}},
177: binterpt[6][3] = {{6943876665148./7220017795957.,-54480133./30881146.,6818779379841./7100303317025.},
178: {0,0,0},
179: {7640104374378./9702883013639.,-11436875./14766696.,2173542590792./12501825683035.},
180: {-20649996744609./7521556579894.,174696575./18121608.,-31592104683404./5083833661969.},
181: {8854892464581./2390941311638.,-12120380./966161.,61146701046299./7138195549469.},
182: {-11397109935349./6675773540249.,3843./706.,-17219254887155./4939391667607.}};
183: TSARKIMEXRegister(TSARKIMEX4,4,6,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,3,binterpt[0],PETSC_NULL);
184: }
185: {
186: const PetscReal
187: A[8][8] = {{0,0,0,0,0,0,0,0},
188: {41./100,0,0,0,0,0,0,0},
189: {367902744464./2072280473677.,677623207551./8224143866563.,0,0,0,0,0,0},
190: {1268023523408./10340822734521.,0,1029933939417./13636558850479.,0,0,0,0,0},
191: {14463281900351./6315353703477.,0,66114435211212./5879490589093.,-54053170152839./4284798021562.,0,0,0,0},
192: {14090043504691./34967701212078.,0,15191511035443./11219624916014.,-18461159152457./12425892160975.,-281667163811./9011619295870.,0,0,0},
193: {19230459214898./13134317526959.,0,21275331358303./2942455364971.,-38145345988419./4862620318723.,-1./8,-1./8,0,0},
194: {-19977161125411./11928030595625.,0,-40795976796054./6384907823539.,177454434618887./12078138498510.,782672205425./8267701900261.,-69563011059811./9646580694205.,7356628210526./4942186776405.,0}},
195: At[8][8] = {{0,0,0,0,0,0,0,0},
196: {41./200.,41./200.,0,0,0,0,0,0},
197: {41./400.,-567603406766./11931857230679.,41./200.,0,0,0,0,0},
198: {683785636431./9252920307686.,0,-110385047103./1367015193373.,41./200.,0,0,0,0},
199: {3016520224154./10081342136671.,0,30586259806659./12414158314087.,-22760509404356./11113319521817.,41./200.,0,0,0},
200: {218866479029./1489978393911.,0,638256894668./5436446318841.,-1179710474555./5321154724896.,-60928119172./8023461067671.,41./200.,0,0},
201: {1020004230633./5715676835656.,0,25762820946817./25263940353407.,-2161375909145./9755907335909.,-211217309593./5846859502534.,-4269925059573./7827059040749.,41./200,0},
202: {-872700587467./9133579230613.,0,0,22348218063261./9555858737531.,-1143369518992./8141816002931.,-39379526789629./19018526304540.,32727382324388./42900044865799.,41./200.}},
203: binterpt[8][3] = {{-17674230611817./10670229744614. , 43486358583215./12773830924787. , -9257016797708./5021505065439.},
204: {0 , 0 , 0 },
205: {0 , 0 , 0 },
206: {65168852399939./7868540260826. , -91478233927265./11067650958493., 26096422576131./11239449250142.},
207: {15494834004392./5936557850923. , -79368583304911./10890268929626., 92396832856987./20362823103730.},
208: {-99329723586156./26959484932159., -12239297817655./9152339842473. , 30029262896817./10175596800299.},
209: {-19024464361622./5461577185407. , 115839755401235./10719374521269., -26136350496073./3983972220547.},
210: {-6511271360970./6095937251113. , 5843115559534./2180450260947. , -5289405421727./3760307252460. }};
211: TSARKIMEXRegister(TSARKIMEX5,5,8,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,3,binterpt[0],PETSC_NULL);
212: }
214: return(0);
215: }
219: /*@C
220: TSARKIMEXRegisterDestroy - Frees the list of schemes that were registered by TSARKIMEXRegister().
222: Not Collective
224: Level: advanced
226: .keywords: TSARKIMEX, register, destroy
227: .seealso: TSARKIMEXRegister(), TSARKIMEXRegisterAll(), TSARKIMEXRegisterDynamic()
228: @*/
229: PetscErrorCode TSARKIMEXRegisterDestroy(void)
230: {
232: ARKTableauLink link;
235: while ((link = ARKTableauList)) {
236: ARKTableau t = &link->tab;
237: ARKTableauList = link->next;
238: PetscFree6(t->At,t->bt,t->ct,t->A,t->b,t->c);
239: PetscFree2(t->binterpt,t->binterp);
240: PetscFree(t->name);
241: PetscFree(link);
242: }
243: TSARKIMEXRegisterAllCalled = PETSC_FALSE;
244: return(0);
245: }
249: /*@C
250: TSARKIMEXInitializePackage - This function initializes everything in the TSARKIMEX package. It is called
251: from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_ARKIMEX()
252: when using static libraries.
254: Input Parameter:
255: path - The dynamic library path, or PETSC_NULL
257: Level: developer
259: .keywords: TS, TSARKIMEX, initialize, package
260: .seealso: PetscInitialize()
261: @*/
262: PetscErrorCode TSARKIMEXInitializePackage(const char path[])
263: {
267: if (TSARKIMEXPackageInitialized) return(0);
268: TSARKIMEXPackageInitialized = PETSC_TRUE;
269: TSARKIMEXRegisterAll();
270: PetscRegisterFinalize(TSARKIMEXFinalizePackage);
271: return(0);
272: }
276: /*@C
277: TSARKIMEXFinalizePackage - This function destroys everything in the TSARKIMEX package. It is
278: called from PetscFinalize().
280: Level: developer
282: .keywords: Petsc, destroy, package
283: .seealso: PetscFinalize()
284: @*/
285: PetscErrorCode TSARKIMEXFinalizePackage(void)
286: {
290: TSARKIMEXPackageInitialized = PETSC_FALSE;
291: TSARKIMEXRegisterDestroy();
292: return(0);
293: }
297: /*@C
298: TSARKIMEXRegister - register an ARK IMEX scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
300: Not Collective, but the same schemes should be registered on all processes on which they will be used
302: Input Parameters:
303: + name - identifier for method
304: . order - approximation order of method
305: . s - number of stages, this is the dimension of the matrices below
306: . At - Butcher table of stage coefficients for stiff part (dimension s*s, row-major)
307: . bt - Butcher table for completing the stiff part of the step (dimension s; PETSC_NULL to use the last row of At)
308: . ct - Abscissa of each stiff stage (dimension s, PETSC_NULL to use row sums of At)
309: . A - Non-stiff stage coefficients (dimension s*s, row-major)
310: . b - Non-stiff step completion table (dimension s; PETSC_NULL to use last row of At)
311: . c - Non-stiff abscissa (dimension s; PETSC_NULL to use row sums of A)
312: . pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt and binterp
313: . binterpt - Coefficients of the interpolation formula for the stiff part (dimension s*pinterp)
314: - binterp - Coefficients of the interpolation formula for the non-stiff part (dimension s*pinterp; PETSC_NULL to reuse binterpt)
316: Notes:
317: Several ARK IMEX methods are provided, this function is only needed to create new methods.
319: Level: advanced
321: .keywords: TS, register
323: .seealso: TSARKIMEX
324: @*/
325: PetscErrorCode TSARKIMEXRegister(const TSARKIMEXType name,PetscInt order,PetscInt s,
326: const PetscReal At[],const PetscReal bt[],const PetscReal ct[],
327: const PetscReal A[],const PetscReal b[],const PetscReal c[],
328: PetscInt pinterp,const PetscReal binterpt[],const PetscReal binterp[])
329: {
331: ARKTableauLink link;
332: ARKTableau t;
333: PetscInt i,j;
336: PetscMalloc(sizeof(*link),&link);
337: PetscMemzero(link,sizeof(*link));
338: t = &link->tab;
339: PetscStrallocpy(name,&t->name);
340: t->order = order;
341: t->s = s;
342: 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);
343: PetscMemcpy(t->At,At,s*s*sizeof(At[0]));
344: PetscMemcpy(t->A,A,s*s*sizeof(A[0]));
345: if (bt) {PetscMemcpy(t->bt,bt,s*sizeof(bt[0]));}
346: else for (i=0; i<s; i++) t->bt[i] = At[(s-1)*s+i];
347: if (b) {PetscMemcpy(t->b,b,s*sizeof(b[0]));}
348: else for (i=0; i<s; i++) t->b[i] = At[(s-1)*s+i];
349: if (ct) {PetscMemcpy(t->ct,ct,s*sizeof(ct[0]));}
350: else for (i=0; i<s; i++) for (j=0,t->ct[i]=0; j<s; j++) t->ct[i] += At[i*s+j];
351: if (c) {PetscMemcpy(t->c,c,s*sizeof(c[0]));}
352: else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j];
353: t->pinterp = pinterp;
354: PetscMalloc2(s*pinterp,PetscReal,&t->binterpt,s*pinterp,PetscReal,&t->binterp);
355: PetscMemcpy(t->binterpt,binterpt,s*pinterp*sizeof(binterpt[0]));
356: PetscMemcpy(t->binterp,binterp?binterp:binterpt,s*pinterp*sizeof(binterpt[0]));
357: link->next = ARKTableauList;
358: ARKTableauList = link;
359: return(0);
360: }
364: static PetscErrorCode TSStep_ARKIMEX(TS ts)
365: {
366: TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
367: ARKTableau tab = ark->tableau;
368: const PetscInt s = tab->s;
369: const PetscReal *At = tab->At,*A = tab->A,*bt = tab->bt,*b = tab->b,*ct = tab->ct,*c = tab->c;
370: PetscScalar *w = ark->work;
371: Vec *Y = ark->Y,*YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS,Ydot = ark->Ydot,W = ark->Work,Z = ark->Z;
372: SNES snes;
373: PetscInt i,j,its,lits;
374: PetscReal next_time_step;
375: PetscReal h,t;
376: PetscErrorCode ierr;
379: TSGetSNES(ts,&snes);
380: next_time_step = ts->time_step;
381: h = ts->time_step;
382: t = ts->ptime;
384: for (i=0; i<s; i++) {
385: if (At[i*s+i] == 0) { /* This stage is explicit */
386: VecCopy(ts->vec_sol,Y[i]);
387: for (j=0; j<i; j++) w[j] = -h*At[i*s+j];
388: VecMAXPY(Y[i],i,w,YdotI);
389: for (j=0; j<i; j++) w[j] = h*A[i*s+j];
390: VecMAXPY(Y[i],i,w,YdotRHS);
391: } else {
392: ark->stage_time = t + h*ct[i];
393: ark->shift = 1./(h*At[i*s+i]);
394: /* Affine part */
395: VecZeroEntries(W);
396: for (j=0; j<i; j++) w[j] = h*A[i*s+j];
397: VecMAXPY(W,i,w,YdotRHS);
398: /* Ydot = shift*(Y-Z) */
399: VecCopy(ts->vec_sol,Z);
400: for (j=0; j<i; j++) w[j] = -h*At[i*s+j];
401: VecMAXPY(Z,i,w,YdotI);
402: /* Initial guess taken from last stage */
403: VecCopy(i>0?Y[i-1]:ts->vec_sol,Y[i]);
404: SNESSolve(snes,W,Y[i]);
405: SNESGetIterationNumber(snes,&its);
406: SNESGetLinearSolveIterations(snes,&lits);
407: ts->nonlinear_its += its; ts->linear_its += lits;
408: }
409: VecZeroEntries(Ydot);
410: TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],ark->imex);
411: if (ark->imex) {
412: TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);
413: } else {
414: VecZeroEntries(YdotRHS[i]);
415: }
416: }
417: for (j=0; j<s; j++) w[j] = -h*bt[j];
418: VecMAXPY(ts->vec_sol,s,w,YdotI);
419: for (j=0; j<s; j++) w[j] = h*b[j];
420: VecMAXPY(ts->vec_sol,s,w,YdotRHS);
422: ts->ptime += ts->time_step;
423: ts->time_step = next_time_step;
424: ts->steps++;
425: return(0);
426: }
430: static PetscErrorCode TSInterpolate_ARKIMEX(TS ts,PetscReal itime,Vec X)
431: {
432: TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
433: PetscInt s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j;
434: PetscReal tt,t = (itime - ts->ptime)/ts->time_step + 1; /* In the interval [0,1] */
435: PetscScalar *bt,*b;
436: const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp;
440: if (!Bt || !B) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name);
441: PetscMalloc2(s,PetscScalar,&bt,s,PetscScalar,&b);
442: for (i=0; i<s; i++) bt[i] = b[i] = 0;
443: for (j=0,tt=t; j<pinterp; j++,tt*=t) {
444: for (i=0; i<s; i++) {
445: bt[i] += ts->time_step * Bt[i*pinterp+j] * tt * -1.0;
446: b[i] += ts->time_step * B[i*pinterp+j] * tt;
447: }
448: }
449: 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");
450: VecCopy(ark->Y[0],X);
451: VecMAXPY(X,s,bt,ark->YdotI);
452: VecMAXPY(X,s,b,ark->YdotRHS);
453: PetscFree2(bt,b);
454: return(0);
455: }
457: /*------------------------------------------------------------*/
460: static PetscErrorCode TSReset_ARKIMEX(TS ts)
461: {
462: TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
463: PetscInt s;
464: PetscErrorCode ierr;
467: if (!ark->tableau) return(0);
468: s = ark->tableau->s;
469: VecDestroyVecs(s,&ark->Y);
470: VecDestroyVecs(s,&ark->YdotI);
471: VecDestroyVecs(s,&ark->YdotRHS);
472: VecDestroy(&ark->Ydot);
473: VecDestroy(&ark->Work);
474: VecDestroy(&ark->Z);
475: PetscFree(ark->work);
476: return(0);
477: }
481: static PetscErrorCode TSDestroy_ARKIMEX(TS ts)
482: {
483: PetscErrorCode ierr;
486: TSReset_ARKIMEX(ts);
487: PetscFree(ts->data);
488: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXGetType_C","",PETSC_NULL);
489: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetType_C","",PETSC_NULL);
490: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C","",PETSC_NULL);
491: return(0);
492: }
494: /*
495: This defines the nonlinear equation that is to be solved with SNES
496: G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
497: */
500: static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec F,TS ts)
501: {
502: TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
506: VecAXPBYPCZ(ark->Ydot,-ark->shift,ark->shift,0,ark->Z,X); /* Ydot = shift*(X-Z) */
507: TSComputeIFunction(ts,ark->stage_time,X,ark->Ydot,F,ark->imex);
508: return(0);
509: }
513: static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
514: {
515: TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
519: /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */
520: TSComputeIJacobian(ts,ark->stage_time,X,ark->Ydot,ark->shift,A,B,str,PETSC_TRUE);
521: return(0);
522: }
526: static PetscErrorCode TSSetUp_ARKIMEX(TS ts)
527: {
528: TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
529: ARKTableau tab = ark->tableau;
530: PetscInt s = tab->s;
534: if (!ark->tableau) {
535: TSARKIMEXSetType(ts,TSARKIMEXDefault);
536: }
537: VecDuplicateVecs(ts->vec_sol,s,&ark->Y);
538: VecDuplicateVecs(ts->vec_sol,s,&ark->YdotI);
539: VecDuplicateVecs(ts->vec_sol,s,&ark->YdotRHS);
540: VecDuplicate(ts->vec_sol,&ark->Ydot);
541: VecDuplicate(ts->vec_sol,&ark->Work);
542: VecDuplicate(ts->vec_sol,&ark->Z);
543: PetscMalloc(s*sizeof(ark->work[0]),&ark->work);
544: return(0);
545: }
546: /*------------------------------------------------------------*/
550: static PetscErrorCode TSSetFromOptions_ARKIMEX(TS ts)
551: {
552: TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
554: char arktype[256];
557: PetscOptionsHead("ARKIMEX ODE solver options");
558: {
559: ARKTableauLink link;
560: PetscInt count,choice;
561: PetscBool flg;
562: const char **namelist;
563: PetscStrncpy(arktype,TSARKIMEXDefault,sizeof arktype);
564: for (link=ARKTableauList,count=0; link; link=link->next,count++) ;
565: PetscMalloc(count*sizeof(char*),&namelist);
566: for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
567: PetscOptionsEList("-ts_arkimex_type","Family of ARK IMEX method","TSARKIMEXSetType",(const char*const*)namelist,count,arktype,&choice,&flg);
568: TSARKIMEXSetType(ts,flg ? namelist[choice] : arktype);
569: PetscFree(namelist);
570: flg = (PetscBool)!ark->imex;
571: PetscOptionsBool("-ts_arkimex_fully_implicit","Solve the problem fully implicitly","TSARKIMEXSetFullyImplicit",flg,&flg,PETSC_NULL);
572: ark->imex = (PetscBool)!flg;
573: SNESSetFromOptions(ts->snes);
574: }
575: PetscOptionsTail();
576: return(0);
577: }
581: static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
582: {
584: PetscInt i;
585: size_t left,count;
586: char *p;
589: for (i=0,p=buf,left=len; i<n; i++) {
590: PetscSNPrintfCount(p,left,fmt,&count,x[i]);
591: if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
592: left -= count;
593: p += count;
594: *p++ = ' ';
595: }
596: p[i ? 0 : -1] = 0;
597: return(0);
598: }
602: static PetscErrorCode TSView_ARKIMEX(TS ts,PetscViewer viewer)
603: {
604: TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
605: ARKTableau tab = ark->tableau;
606: PetscBool iascii;
610: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
611: if (iascii) {
612: const TSARKIMEXType arktype;
613: char buf[512];
614: TSARKIMEXGetType(ts,&arktype);
615: PetscViewerASCIIPrintf(viewer," ARK IMEX %s\n",arktype);
616: PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ct);
617: PetscViewerASCIIPrintf(viewer," Stiff abscissa ct = %s\n",buf);
618: PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->c);
619: PetscViewerASCIIPrintf(viewer," Nonstiff abscissa c = %s\n",buf);
620: }
621: SNESView(ts->snes,viewer);
622: return(0);
623: }
627: /*@C
628: TSARKIMEXSetType - Set the type of ARK IMEX scheme
630: Logically collective
632: Input Parameter:
633: + ts - timestepping context
634: - arktype - type of ARK-IMEX scheme
636: Level: intermediate
638: .seealso: TSARKIMEXGetType()
639: @*/
640: PetscErrorCode TSARKIMEXSetType(TS ts,const TSARKIMEXType arktype)
641: {
646: PetscTryMethod(ts,"TSARKIMEXSetType_C",(TS,const TSARKIMEXType),(ts,arktype));
647: return(0);
648: }
652: /*@C
653: TSARKIMEXGetType - Get the type of ARK IMEX scheme
655: Logically collective
657: Input Parameter:
658: . ts - timestepping context
660: Output Parameter:
661: . arktype - type of ARK-IMEX scheme
663: Level: intermediate
665: .seealso: TSARKIMEXGetType()
666: @*/
667: PetscErrorCode TSARKIMEXGetType(TS ts,const TSARKIMEXType *arktype)
668: {
673: PetscUseMethod(ts,"TSARKIMEXGetType_C",(TS,const TSARKIMEXType*),(ts,arktype));
674: return(0);
675: }
679: /*@C
680: TSARKIMEXSetFullyImplicit - Solve both parts of the equation implicitly
682: Logically collective
684: Input Parameter:
685: + ts - timestepping context
686: - flg - PETSC_TRUE for fully implicit
688: Level: intermediate
690: .seealso: TSARKIMEXGetType()
691: @*/
692: PetscErrorCode TSARKIMEXSetFullyImplicit(TS ts,PetscBool flg)
693: {
698: PetscTryMethod(ts,"TSARKIMEXSetFullyImplicit_C",(TS,PetscBool),(ts,flg));
699: return(0);
700: }
705: PetscErrorCode TSARKIMEXGetType_ARKIMEX(TS ts,const TSARKIMEXType *arktype)
706: {
707: TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
711: if (!ark->tableau) {TSARKIMEXSetType(ts,TSARKIMEXDefault);}
712: *arktype = ark->tableau->name;
713: return(0);
714: }
717: PetscErrorCode TSARKIMEXSetType_ARKIMEX(TS ts,const TSARKIMEXType arktype)
718: {
719: TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
721: PetscBool match;
722: ARKTableauLink link;
725: if (ark->tableau) {
726: PetscStrcmp(ark->tableau->name,arktype,&match);
727: if (match) return(0);
728: }
729: for (link = ARKTableauList; link; link=link->next) {
730: PetscStrcmp(link->tab.name,arktype,&match);
731: if (match) {
732: TSReset_ARKIMEX(ts);
733: ark->tableau = &link->tab;
734: return(0);
735: }
736: }
737: SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype);
738: return(0);
739: }
742: PetscErrorCode TSARKIMEXSetFullyImplicit_ARKIMEX(TS ts,PetscBool flg)
743: {
744: TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
747: ark->imex = (PetscBool)!flg;
748: return(0);
749: }
752: /* ------------------------------------------------------------ */
753: /*MC
754: TSARKIMEX - ODE solver using Additive Runge-Kutta IMEX schemes
756: These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
757: nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
758: of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
760: Notes:
761: The default is TSARKIMEX2E, it can be changed with TSARKIMEXSetType() or -ts_arkimex_type
763: This method currently only works with ODE, for which the stiff part G(t,X,Xdot) has the form Xdot + Ghat(t,X).
765: Level: beginner
767: .seealso: TSCreate(), TS, TSSetType(), TSARKIMEXSetType(), TSARKIMEXGetType(), TSARKIMEXSetFullyImplicit(), TSARKIMEX2D, TTSARKIMEX2E, TSARKIMEX3,
768: TSARKIMEX4, TSARKIMEX5, TSARKIMEXType, TSARKIMEXRegister()
770: M*/
774: PetscErrorCode TSCreate_ARKIMEX(TS ts)
775: {
776: TS_ARKIMEX *th;
780: #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
781: TSARKIMEXInitializePackage(PETSC_NULL);
782: #endif
784: ts->ops->reset = TSReset_ARKIMEX;
785: ts->ops->destroy = TSDestroy_ARKIMEX;
786: ts->ops->view = TSView_ARKIMEX;
787: ts->ops->setup = TSSetUp_ARKIMEX;
788: ts->ops->step = TSStep_ARKIMEX;
789: ts->ops->interpolate = TSInterpolate_ARKIMEX;
790: ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX;
791: ts->ops->snesfunction = SNESTSFormFunction_ARKIMEX;
792: ts->ops->snesjacobian = SNESTSFormJacobian_ARKIMEX;
794: PetscNewLog(ts,TS_ARKIMEX,&th);
795: ts->data = (void*)th;
796: th->imex = PETSC_TRUE;
798: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXGetType_C","TSARKIMEXGetType_ARKIMEX",TSARKIMEXGetType_ARKIMEX);
799: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetType_C","TSARKIMEXSetType_ARKIMEX",TSARKIMEXSetType_ARKIMEX);
800: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C","TSARKIMEXSetFullyImplicit_ARKIMEX",TSARKIMEXSetFullyImplicit_ARKIMEX);
801: return(0);
802: }