Actual source code: arkimex.c

petsc-3.9.4 2018-09-11
Report Typos and Errors
  1: /*
  2:   Code for timestepping with additive Runge-Kutta IMEX method

  4:   Notes:
  5:   The general system is written as

  7:   F(t,U,Udot) = G(t,U)

  9:   where F represents the stiff part of the physics and G represents the non-stiff part.

 11: */
 12:  #include <petsc/private/tsimpl.h>
 13:  #include <petscdm.h>

 15: static TSARKIMEXType  TSARKIMEXDefault = TSARKIMEX3;
 16: static PetscBool      TSARKIMEXRegisterAllCalled;
 17: static PetscBool      TSARKIMEXPackageInitialized;
 18: static PetscErrorCode TSExtrapolate_ARKIMEX(TS,PetscReal,Vec);

 20: typedef struct _ARKTableau *ARKTableau;
 21: struct _ARKTableau {
 22:   char      *name;
 23:   PetscInt  order;                /* Classical approximation order of the method */
 24:   PetscInt  s;                    /* Number of stages */
 25:   PetscBool stiffly_accurate;     /* The implicit part is stiffly accurate*/
 26:   PetscBool FSAL_implicit;        /* The implicit part is FSAL*/
 27:   PetscBool explicit_first_stage; /* The implicit part has an explicit first stage*/
 28:   PetscInt  pinterp;              /* Interpolation order */
 29:   PetscReal *At,*bt,*ct;          /* Stiff tableau */
 30:   PetscReal *A,*b,*c;             /* Non-stiff tableau */
 31:   PetscReal *bembedt,*bembed;     /* Embedded formula of order one less (order-1) */
 32:   PetscReal *binterpt,*binterp;   /* Dense output formula */
 33:   PetscReal ccfl;                 /* Placeholder for CFL coefficient relative to forward Euler */
 34: };
 35: typedef struct _ARKTableauLink *ARKTableauLink;
 36: struct _ARKTableauLink {
 37:   struct _ARKTableau tab;
 38:   ARKTableauLink     next;
 39: };
 40: static ARKTableauLink ARKTableauList;

 42: typedef struct {
 43:   ARKTableau   tableau;
 44:   Vec          *Y;               /* States computed during the step */
 45:   Vec          *YdotI;           /* Time derivatives for the stiff part */
 46:   Vec          *YdotRHS;         /* Function evaluations for the non-stiff part */
 47:   Vec          *Y_prev;          /* States computed during the previous time step */
 48:   Vec          *YdotI_prev;      /* Time derivatives for the stiff part for the previous time step*/
 49:   Vec          *YdotRHS_prev;    /* Function evaluations for the non-stiff part for the previous time step*/
 50:   Vec          Ydot0;            /* Holds the slope from the previous step in FSAL case */
 51:   Vec          Ydot;             /* Work vector holding Ydot during residual evaluation */
 52:   Vec          Z;                /* Ydot = shift(Y-Z) */
 53:   PetscScalar  *work;            /* Scalar work */
 54:   PetscReal    scoeff;           /* shift = scoeff/dt */
 55:   PetscReal    stage_time;
 56:   PetscBool    imex;
 57:   PetscBool    extrapolate;      /* Extrapolate initial guess from previous time-step stage values */
 58:   TSStepStatus status;
 59: } TS_ARKIMEX;
 60: /*MC
 61:      TSARKIMEXARS122 - Second order ARK IMEX scheme.

 63:      This method has one explicit stage and one implicit stage.

 65:      Options Database:
 66: .      -ts_arkimex_type ars122

 68:      References:
 69: .   1. -  U. Ascher, S. Ruuth, R. J. Spiteri, Implicit explicit Runge Kutta methods for time dependent Partial Differential Equations. Appl. Numer. Math. 25, (1997).

 71:      Level: advanced

 73: .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
 74: M*/
 75: /*MC
 76:      TSARKIMEXA2 - Second order ARK IMEX scheme with A-stable implicit part.

 78:      This method has an explicit stage and one implicit stage, and has an A-stable implicit scheme. This method was provided by Emil Constantinescu.

 80:      Options Database:
 81: .      -ts_arkimex_type a2

 83:      Level: advanced

 85: .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
 86: M*/
 87: /*MC
 88:      TSARKIMEXL2 - Second order ARK IMEX scheme with L-stable implicit part.

 90:      This method has two implicit stages, and L-stable implicit scheme.

 92:      Options Database:
 93: .      -ts_arkimex_type l2

 95:     References:
 96: .   1. -  L. Pareschi, G. Russo, Implicit Explicit Runge Kutta schemes and applications to hyperbolic systems with relaxations. Journal of Scientific Computing Volume: 25, Issue: 1, October, 2005.

 98:      Level: advanced

100: .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
101: M*/
102: /*MC
103:      TSARKIMEX1BEE - First order Backward Euler represented as an ARK IMEX scheme with extrapolation as error estimator. This is a 3-stage method.

105:      This method is aimed at starting the integration of implicit DAEs when explicit first-stage ARK methods are used.

107:      Options Database:
108: .      -ts_arkimex_type 1bee

110:      Level: advanced

112: .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
113: M*/
114: /*MC
115:      TSARKIMEX2C - Second order ARK IMEX scheme with L-stable implicit part.

117:      This method has one explicit stage and two implicit stages. The implicit part is the same as in TSARKIMEX2D and TSARKIMEX2E, but the explicit part has a larger stability region on the negative real axis. This method was provided by Emil Constantinescu.

119:      Options Database:
120: .      -ts_arkimex_type 2c

122:      Level: advanced

124: .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
125: M*/
126: /*MC
127:      TSARKIMEX2D - Second order ARK IMEX scheme with L-stable implicit part.

129:      This method has one explicit stage and two implicit stages. The stability function is independent of the explicit part in the infinity limit of the implict component. This method was provided by Emil Constantinescu.

131:      Options Database:
132: .      -ts_arkimex_type 2d

134:      Level: advanced

136: .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
137: M*/
138: /*MC
139:      TSARKIMEX2E - Second order ARK IMEX scheme with L-stable implicit part.

141:      This method has one explicit stage and two implicit stages. It is is an optimal method developed by Emil Constantinescu.

143:      Options Database:
144: .      -ts_arkimex_type 2e
145:  
146:     Level: advanced

148: .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
149: M*/
150: /*MC
151:      TSARKIMEXPRSSP2 - Second order SSP ARK IMEX scheme.

153:      This method has three implicit stages.

155:      References:
156: .   1. -  L. Pareschi, G. Russo, Implicit Explicit Runge Kutta schemes and applications to hyperbolic systems with relaxations. Journal of Scientific Computing Volume: 25, Issue: 1, October, 2005.

158:      This method is referred to as SSP2-(3,3,2) in http://arxiv.org/abs/1110.4375

160:      Options Database:
161: .      -ts_arkimex_type prssp2

163:      Level: advanced

165: .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
166: M*/
167: /*MC
168:      TSARKIMEX3 - Third order ARK IMEX scheme with L-stable implicit part.

170:      This method has one explicit stage and three implicit stages.

172:      Options Database:
173: .      -ts_arkimex_type 3

175:      References:
176: .   1. -  Kennedy and Carpenter 2003.

178:      Level: advanced

180: .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
181: M*/
182: /*MC
183:      TSARKIMEXARS443 - Third order ARK IMEX scheme.

185:      This method has one explicit stage and four implicit stages.

187:      Options Database:
188: .      -ts_arkimex_type ars443

190:      References:
191: +   1. -  U. Ascher, S. Ruuth, R. J. Spiteri, Implicit explicit Runge Kutta methods for time dependent Partial Differential Equations. Appl. Numer. Math. 25, (1997).
192: -   2. -  This method is referred to as ARS(4,4,3) in http://arxiv.org/abs/1110.4375

194:      Level: advanced

196: .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
197: M*/
198: /*MC
199:      TSARKIMEXBPR3 - Third order ARK IMEX scheme.

201:      This method has one explicit stage and four implicit stages.

203:      Options Database:
204: .      -ts_arkimex_type bpr3

206:      References:
207:  .    This method is referred to as ARK3 in http://arxiv.org/abs/1110.4375

209:      Level: advanced

211: .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
212: M*/
213: /*MC
214:      TSARKIMEX4 - Fourth order ARK IMEX scheme with L-stable implicit part.

216:      This method has one explicit stage and four implicit stages.

218:      Options Database:
219: .      -ts_arkimex_type 4

221:      References:
222: .     Kennedy and Carpenter 2003.

224:      Level: advanced

226: .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
227: M*/
228: /*MC
229:      TSARKIMEX5 - Fifth order ARK IMEX scheme with L-stable implicit part.

231:      This method has one explicit stage and five implicit stages.

233:      Options Database:
234: .      -ts_arkimex_type 5

236:      References:
237: .     Kennedy and Carpenter 2003.

239:      Level: advanced

241: .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
242: M*/

244: /*@C
245:   TSARKIMEXRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSARKIMEX

247:   Not Collective, but should be called by all processes which will need the schemes to be registered

249:   Level: advanced

251: .keywords: TS, TSARKIMEX, register, all

253: .seealso:  TSARKIMEXRegisterDestroy()
254: @*/
255: PetscErrorCode TSARKIMEXRegisterAll(void)
256: {

260:   if (TSARKIMEXRegisterAllCalled) return(0);
261:   TSARKIMEXRegisterAllCalled = PETSC_TRUE;

263:   {
264:     const PetscReal
265:       A[3][3] = {{0.0,0.0,0.0},
266:                  {0.0,0.0,0.0},
267:                  {0.0,0.5,0.0}},
268:       At[3][3] = {{1.0,0.0,0.0},
269:                   {0.0,0.5,0.0},
270:                   {0.0,0.5,0.5}},
271:       b[3]       = {0.0,0.5,0.5},
272:       bembedt[3] = {1.0,0.0,0.0};
273:     TSARKIMEXRegister(TSARKIMEX1BEE,2,3,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL);
274:   }
275:   {
276:     const PetscReal
277:       A[2][2] = {{0.0,0.0},
278:                  {0.5,0.0}},
279:       At[2][2] = {{0.0,0.0},
280:                   {0.0,0.5}},
281:       b[2]       = {0.0,1.0},
282:       bembedt[2] = {0.5,0.5};
283:     /* binterpt[2][2] = {{1.0,-1.0},{0.0,1.0}};  second order dense output has poor stability properties and hence it is not currently in use*/
284:     TSARKIMEXRegister(TSARKIMEXARS122,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL);
285:   }
286:   {
287:     const PetscReal
288:       A[2][2] = {{0.0,0.0},
289:                  {1.0,0.0}},
290:       At[2][2] = {{0.0,0.0},
291:                   {0.5,0.5}},
292:       b[2]       = {0.5,0.5},
293:       bembedt[2] = {0.0,1.0};
294:     /* binterpt[2][2] = {{1.0,-0.5},{0.0,0.5}}  second order dense output has poor stability properties and hence it is not currently in use*/
295:     TSARKIMEXRegister(TSARKIMEXA2,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL);
296:   }
297:   {
298:     /* const PetscReal us2 = 1.0-1.0/PetscSqrtReal((PetscReal)2.0);    Direct evaluation: 0.2928932188134524755992. Used below to ensure all values are available at compile time   */
299:     const PetscReal
300:       A[2][2] = {{0.0,0.0},
301:                  {1.0,0.0}},
302:       At[2][2] = {{0.2928932188134524755992,0.0},
303:                   {1.0-2.0*0.2928932188134524755992,0.2928932188134524755992}},
304:       b[2]       = {0.5,0.5},
305:       bembedt[2] = {0.0,1.0},
306:       binterpt[2][2] = {{  (0.2928932188134524755992-1.0)/(2.0*0.2928932188134524755992-1.0),-1/(2.0*(1.0-2.0*0.2928932188134524755992))},
307:                         {1-(0.2928932188134524755992-1.0)/(2.0*0.2928932188134524755992-1.0),-1/(2.0*(1.0-2.0*0.2928932188134524755992))}},
308:       binterp[2][2] = {{1.0,-0.5},{0.0,0.5}};
309:     TSARKIMEXRegister(TSARKIMEXL2,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,2,binterpt[0],binterp[0]);
310:   }
311:   {
312:     /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),  Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time   */
313:     const PetscReal
314:       A[3][3] = {{0,0,0},
315:                  {2-1.414213562373095048802,0,0},
316:                  {0.5,0.5,0}},
317:       At[3][3] = {{0,0,0},
318:                   {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0},
319:                   {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}},
320:       bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)},
321:       binterpt[3][2] = {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
322:                         {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
323:                         {1.0-1.414213562373095048802,1.0/1.414213562373095048802}};
324:     TSARKIMEXRegister(TSARKIMEX2C,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);
325:   }
326:   {
327:     /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),  Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time   */
328:     const PetscReal
329:       A[3][3] = {{0,0,0},
330:                  {2-1.414213562373095048802,0,0},
331:                  {0.75,0.25,0}},
332:       At[3][3] = {{0,0,0},
333:                   {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0},
334:                   {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}},
335:       bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)},
336:       binterpt[3][2] =  {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
337:                          {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
338:                          {1.0-1.414213562373095048802,1.0/1.414213562373095048802}};
339:     TSARKIMEXRegister(TSARKIMEX2D,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);
340:   }
341:   {                             /* Optimal for linear implicit part */
342:     /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),  Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time   */
343:     const PetscReal
344:       A[3][3] = {{0,0,0},
345:                  {2-1.414213562373095048802,0,0},
346:                  {(3-2*1.414213562373095048802)/6,(3+2*1.414213562373095048802)/6,0}},
347:       At[3][3] = {{0,0,0},
348:                   {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0},
349:                   {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}},
350:       bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)},
351:       binterpt[3][2] =  {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
352:                          {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
353:                          {1.0-1.414213562373095048802,1.0/1.414213562373095048802}};
354:     TSARKIMEXRegister(TSARKIMEX2E,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);
355:   }
356:   {                             /* Optimal for linear implicit part */
357:     const PetscReal
358:       A[3][3] = {{0,0,0},
359:                  {0.5,0,0},
360:                  {0.5,0.5,0}},
361:       At[3][3] = {{0.25,0,0},
362:                   {0,0.25,0},
363:                   {1./3,1./3,1./3}};
364:     TSARKIMEXRegister(TSARKIMEXPRSSP2,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,NULL,NULL,0,NULL,NULL);
365:   }
366:   {
367:     const PetscReal
368:       A[4][4] = {{0,0,0,0},
369:                  {1767732205903./2027836641118.,0,0,0},
370:                  {5535828885825./10492691773637.,788022342437./10882634858940.,0,0},
371:                  {6485989280629./16251701735622.,-4246266847089./9704473918619.,10755448449292./10357097424841.,0}},
372:       At[4][4] = {{0,0,0,0},
373:                   {1767732205903./4055673282236.,1767732205903./4055673282236.,0,0},
374:                   {2746238789719./10658868560708.,-640167445237./6845629431997.,1767732205903./4055673282236.,0},
375:                   {1471266399579./7840856788654.,-4482444167858./7529755066697.,11266239266428./11593286722821.,1767732205903./4055673282236.}},
376:       bembedt[4]     = {2756255671327./12835298489170.,-10771552573575./22201958757719.,9247589265047./10645013368117.,2193209047091./5459859503100.},
377:       binterpt[4][2] = {{4655552711362./22874653954995., -215264564351./13552729205753.},
378:                         {-18682724506714./9892148508045.,17870216137069./13817060693119.},
379:                         {34259539580243./13192909600954.,-28141676662227./17317692491321.},
380:                         {584795268549./6622622206610.,   2508943948391./7218656332882.}};
381:     TSARKIMEXRegister(TSARKIMEX3,3,4,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);
382:   }
383:   {
384:     const PetscReal
385:       A[5][5] = {{0,0,0,0,0},
386:                  {1./2,0,0,0,0},
387:                  {11./18,1./18,0,0,0},
388:                  {5./6,-5./6,.5,0,0},
389:                  {1./4,7./4,3./4,-7./4,0}},
390:       At[5][5] = {{0,0,0,0,0},
391:                   {0,1./2,0,0,0},
392:                   {0,1./6,1./2,0,0},
393:                   {0,-1./2,1./2,1./2,0},
394:                   {0,3./2,-3./2,1./2,1./2}},
395:     *bembedt = NULL;
396:     TSARKIMEXRegister(TSARKIMEXARS443,3,5,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,0,NULL,NULL);
397:   }
398:   {
399:     const PetscReal
400:       A[5][5] = {{0,0,0,0,0},
401:                  {1,0,0,0,0},
402:                  {4./9,2./9,0,0,0},
403:                  {1./4,0,3./4,0,0},
404:                  {1./4,0,3./5,0,0}},
405:       At[5][5] = {{0,0,0,0,0},
406:                   {.5,.5,0,0,0},
407:                   {5./18,-1./9,.5,0,0},
408:                   {.5,0,0,.5,0},
409:                   {.25,0,.75,-.5,.5}},
410:     *bembedt = NULL;
411:     TSARKIMEXRegister(TSARKIMEXBPR3,3,5,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,0,NULL,NULL);
412:   }
413:   {
414:     const PetscReal
415:       A[6][6] = {{0,0,0,0,0,0},
416:                  {1./2,0,0,0,0,0},
417:                  {13861./62500.,6889./62500.,0,0,0,0},
418:                  {-116923316275./2393684061468.,-2731218467317./15368042101831.,9408046702089./11113171139209.,0,0,0},
419:                  {-451086348788./2902428689909.,-2682348792572./7519795681897.,12662868775082./11960479115383.,3355817975965./11060851509271.,0,0},
420:                  {647845179188./3216320057751.,73281519250./8382639484533.,552539513391./3454668386233.,3354512671639./8306763924573.,4040./17871.,0}},
421:       At[6][6] = {{0,0,0,0,0,0},
422:                   {1./4,1./4,0,0,0,0},
423:                   {8611./62500.,-1743./31250.,1./4,0,0,0},
424:                   {5012029./34652500.,-654441./2922500.,174375./388108.,1./4,0,0},
425:                   {15267082809./155376265600.,-71443401./120774400.,730878875./902184768.,2285395./8070912.,1./4,0},
426:                   {82889./524892.,0,15625./83664.,69875./102672.,-2260./8211,1./4}},
427:       bembedt[6]     = {4586570599./29645900160.,0,178811875./945068544.,814220225./1159782912.,-3700637./11593932.,61727./225920.},
428:       binterpt[6][3] = {{6943876665148./7220017795957.,-54480133./30881146.,6818779379841./7100303317025.},
429:                         {0,0,0},
430:                         {7640104374378./9702883013639.,-11436875./14766696.,2173542590792./12501825683035.},
431:                         {-20649996744609./7521556579894.,174696575./18121608.,-31592104683404./5083833661969.},
432:                         {8854892464581./2390941311638.,-12120380./966161.,61146701046299./7138195549469.},
433:                         {-11397109935349./6675773540249.,3843./706.,-17219254887155./4939391667607.}};
434:     TSARKIMEXRegister(TSARKIMEX4,4,6,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,3,binterpt[0],NULL);
435:   }
436:   {
437:     const PetscReal
438:       A[8][8] = {{0,0,0,0,0,0,0,0},
439:                  {41./100,0,0,0,0,0,0,0},
440:                  {367902744464./2072280473677.,677623207551./8224143866563.,0,0,0,0,0,0},
441:                  {1268023523408./10340822734521.,0,1029933939417./13636558850479.,0,0,0,0,0},
442:                  {14463281900351./6315353703477.,0,66114435211212./5879490589093.,-54053170152839./4284798021562.,0,0,0,0},
443:                  {14090043504691./34967701212078.,0,15191511035443./11219624916014.,-18461159152457./12425892160975.,-281667163811./9011619295870.,0,0,0},
444:                  {19230459214898./13134317526959.,0,21275331358303./2942455364971.,-38145345988419./4862620318723.,-1./8,-1./8,0,0},
445:                  {-19977161125411./11928030595625.,0,-40795976796054./6384907823539.,177454434618887./12078138498510.,782672205425./8267701900261.,-69563011059811./9646580694205.,7356628210526./4942186776405.,0}},
446:       At[8][8] = {{0,0,0,0,0,0,0,0},
447:                   {41./200.,41./200.,0,0,0,0,0,0},
448:                   {41./400.,-567603406766./11931857230679.,41./200.,0,0,0,0,0},
449:                   {683785636431./9252920307686.,0,-110385047103./1367015193373.,41./200.,0,0,0,0},
450:                   {3016520224154./10081342136671.,0,30586259806659./12414158314087.,-22760509404356./11113319521817.,41./200.,0,0,0},
451:                   {218866479029./1489978393911.,0,638256894668./5436446318841.,-1179710474555./5321154724896.,-60928119172./8023461067671.,41./200.,0,0},
452:                   {1020004230633./5715676835656.,0,25762820946817./25263940353407.,-2161375909145./9755907335909.,-211217309593./5846859502534.,-4269925059573./7827059040749.,41./200,0},
453:                   {-872700587467./9133579230613.,0,0,22348218063261./9555858737531.,-1143369518992./8141816002931.,-39379526789629./19018526304540.,32727382324388./42900044865799.,41./200.}},
454:       bembedt[8]     = {-975461918565./9796059967033.,0,0,78070527104295./32432590147079.,-548382580838./3424219808633.,-33438840321285./15594753105479.,3629800801594./4656183773603.,4035322873751./18575991585200.},
455:       binterpt[8][3] = {{-17674230611817./10670229744614.,  43486358583215./12773830924787., -9257016797708./5021505065439.},
456:                         {0,  0, 0                            },
457:                         {0,  0, 0                            },
458:                         {65168852399939./7868540260826.,  -91478233927265./11067650958493., 26096422576131./11239449250142.},
459:                         {15494834004392./5936557850923.,  -79368583304911./10890268929626., 92396832856987./20362823103730.},
460:                         {-99329723586156./26959484932159.,  -12239297817655./9152339842473., 30029262896817./10175596800299.},
461:                         {-19024464361622./5461577185407.,  115839755401235./10719374521269., -26136350496073./3983972220547.},
462:                         {-6511271360970./6095937251113.,  5843115559534./2180450260947., -5289405421727./3760307252460. }};
463:     TSARKIMEXRegister(TSARKIMEX5,5,8,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,3,binterpt[0],NULL);
464:   }
465:   return(0);
466: }

468: /*@C
469:    TSARKIMEXRegisterDestroy - Frees the list of schemes that were registered by TSARKIMEXRegister().

471:    Not Collective

473:    Level: advanced

475: .keywords: TSARKIMEX, register, destroy
476: .seealso: TSARKIMEXRegister(), TSARKIMEXRegisterAll()
477: @*/
478: PetscErrorCode TSARKIMEXRegisterDestroy(void)
479: {
481:   ARKTableauLink link;

484:   while ((link = ARKTableauList)) {
485:     ARKTableau t = &link->tab;
486:     ARKTableauList = link->next;
487:     PetscFree6(t->At,t->bt,t->ct,t->A,t->b,t->c);
488:     PetscFree2(t->bembedt,t->bembed);
489:     PetscFree2(t->binterpt,t->binterp);
490:     PetscFree(t->name);
491:     PetscFree(link);
492:   }
493:   TSARKIMEXRegisterAllCalled = PETSC_FALSE;
494:   return(0);
495: }

497: /*@C
498:   TSARKIMEXInitializePackage - This function initializes everything in the TSARKIMEX package. It is called
499:   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_ARKIMEX()
500:   when using static libraries.

502:   Level: developer

504: .keywords: TS, TSARKIMEX, initialize, package
505: .seealso: PetscInitialize()
506: @*/
507: PetscErrorCode TSARKIMEXInitializePackage(void)
508: {

512:   if (TSARKIMEXPackageInitialized) return(0);
513:   TSARKIMEXPackageInitialized = PETSC_TRUE;
514:   TSARKIMEXRegisterAll();
515:   PetscRegisterFinalize(TSARKIMEXFinalizePackage);
516:   return(0);
517: }

519: /*@C
520:   TSARKIMEXFinalizePackage - This function destroys everything in the TSARKIMEX package. It is
521:   called from PetscFinalize().

523:   Level: developer

525: .keywords: Petsc, destroy, package
526: .seealso: PetscFinalize()
527: @*/
528: PetscErrorCode TSARKIMEXFinalizePackage(void)
529: {

533:   TSARKIMEXPackageInitialized = PETSC_FALSE;
534:   TSARKIMEXRegisterDestroy();
535:   return(0);
536: }

538: /*@C
539:    TSARKIMEXRegister - register an ARK IMEX scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation

541:    Not Collective, but the same schemes should be registered on all processes on which they will be used

543:    Input Parameters:
544: +  name - identifier for method
545: .  order - approximation order of method
546: .  s - number of stages, this is the dimension of the matrices below
547: .  At - Butcher table of stage coefficients for stiff part (dimension s*s, row-major)
548: .  bt - Butcher table for completing the stiff part of the step (dimension s; NULL to use the last row of At)
549: .  ct - Abscissa of each stiff stage (dimension s, NULL to use row sums of At)
550: .  A - Non-stiff stage coefficients (dimension s*s, row-major)
551: .  b - Non-stiff step completion table (dimension s; NULL to use last row of At)
552: .  c - Non-stiff abscissa (dimension s; NULL to use row sums of A)
553: .  bembedt - Stiff part of completion table for embedded method (dimension s; NULL if not available)
554: .  bembed - Non-stiff part of completion table for embedded method (dimension s; NULL to use bembedt if provided)
555: .  pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt and binterp
556: .  binterpt - Coefficients of the interpolation formula for the stiff part (dimension s*pinterp)
557: -  binterp - Coefficients of the interpolation formula for the non-stiff part (dimension s*pinterp; NULL to reuse binterpt)

559:    Notes:
560:    Several ARK IMEX methods are provided, this function is only needed to create new methods.

562:    Level: advanced

564: .keywords: TS, register

566: .seealso: TSARKIMEX
567: @*/
568: PetscErrorCode TSARKIMEXRegister(TSARKIMEXType name,PetscInt order,PetscInt s,
569:                                  const PetscReal At[],const PetscReal bt[],const PetscReal ct[],
570:                                  const PetscReal A[],const PetscReal b[],const PetscReal c[],
571:                                  const PetscReal bembedt[],const PetscReal bembed[],
572:                                  PetscInt pinterp,const PetscReal binterpt[],const PetscReal binterp[])
573: {
575:   ARKTableauLink link;
576:   ARKTableau     t;
577:   PetscInt       i,j;

580:   PetscCalloc1(1,&link);
581:   t        = &link->tab;
582:   PetscStrallocpy(name,&t->name);
583:   t->order = order;
584:   t->s     = s;
585:   PetscMalloc6(s*s,&t->At,s,&t->bt,s,&t->ct,s*s,&t->A,s,&t->b,s,&t->c);
586:   PetscMemcpy(t->At,At,s*s*sizeof(At[0]));
587:   PetscMemcpy(t->A,A,s*s*sizeof(A[0]));
588:   if (bt) { PetscMemcpy(t->bt,bt,s*sizeof(bt[0])); }
589:   else for (i=0; i<s; i++) t->bt[i] = At[(s-1)*s+i];
590:   if (b)  { PetscMemcpy(t->b,b,s*sizeof(b[0])); }
591:   else for (i=0; i<s; i++) t->b[i] = t->bt[i];
592:   if (ct) { PetscMemcpy(t->ct,ct,s*sizeof(ct[0])); }
593:   else for (i=0; i<s; i++) for (j=0,t->ct[i]=0; j<s; j++) t->ct[i] += At[i*s+j];
594:   if (c)  { PetscMemcpy(t->c,c,s*sizeof(c[0])); }
595:   else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j];
596:   t->stiffly_accurate = PETSC_TRUE;
597:   for (i=0; i<s; i++) if (t->At[(s-1)*s+i] != t->bt[i]) t->stiffly_accurate = PETSC_FALSE;
598:   t->explicit_first_stage = PETSC_TRUE;
599:   for (i=0; i<s; i++) if (t->At[i] != 0.0) t->explicit_first_stage = PETSC_FALSE;
600:   /*def of FSAL can be made more precise*/
601:   t->FSAL_implicit = (PetscBool)(t->explicit_first_stage && t->stiffly_accurate);
602:   if (bembedt) {
603:     PetscMalloc2(s,&t->bembedt,s,&t->bembed);
604:     PetscMemcpy(t->bembedt,bembedt,s*sizeof(bembedt[0]));
605:     PetscMemcpy(t->bembed,bembed ? bembed : bembedt,s*sizeof(bembed[0]));
606:   }

608:   t->pinterp     = pinterp;
609:   PetscMalloc2(s*pinterp,&t->binterpt,s*pinterp,&t->binterp);
610:   PetscMemcpy(t->binterpt,binterpt,s*pinterp*sizeof(binterpt[0]));
611:   PetscMemcpy(t->binterp,binterp ? binterp : binterpt,s*pinterp*sizeof(binterpt[0]));
612:   link->next     = ARKTableauList;
613:   ARKTableauList = link;
614:   return(0);
615: }

617: /*
618:  The step completion formula is

620:  x1 = x0 - h bt^T YdotI + h b^T YdotRHS

622:  This function can be called before or after ts->vec_sol has been updated.
623:  Suppose we have a completion formula (bt,b) and an embedded formula (bet,be) of different order.
624:  We can write

626:  x1e = x0 - h bet^T YdotI + h be^T YdotRHS
627:      = x1 + h bt^T YdotI - h b^T YdotRHS - h bet^T YdotI + h be^T YdotRHS
628:      = x1 - h (bet - bt)^T YdotI + h (be - b)^T YdotRHS

630:  so we can evaluate the method with different order even after the step has been optimistically completed.
631: */
632: static PetscErrorCode TSEvaluateStep_ARKIMEX(TS ts,PetscInt order,Vec X,PetscBool *done)
633: {
634:   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
635:   ARKTableau     tab  = ark->tableau;
636:   PetscScalar    *w   = ark->work;
637:   PetscReal      h;
638:   PetscInt       s = tab->s,j;

642:   switch (ark->status) {
643:   case TS_STEP_INCOMPLETE:
644:   case TS_STEP_PENDING:
645:     h = ts->time_step; break;
646:   case TS_STEP_COMPLETE:
647:     h = ts->ptime - ts->ptime_prev; break;
648:   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
649:   }
650:   if (order == tab->order) {
651:     if (ark->status == TS_STEP_INCOMPLETE) {
652:       if (!ark->imex && tab->stiffly_accurate) { /* Only the stiffly accurate implicit formula is used */
653:         VecCopy(ark->Y[s-1],X);
654:       } else { /* Use the standard completion formula (bt,b) */
655:         VecCopy(ts->vec_sol,X);
656:         for (j=0; j<s; j++) w[j] = h*tab->bt[j];
657:         VecMAXPY(X,s,w,ark->YdotI);
658:         if (ark->imex) { /* Method is IMEX, complete the explicit formula */
659:           for (j=0; j<s; j++) w[j] = h*tab->b[j];
660:           VecMAXPY(X,s,w,ark->YdotRHS);
661:         }
662:       }
663:     } else {VecCopy(ts->vec_sol,X);}
664:     if (done) *done = PETSC_TRUE;
665:     return(0);
666:   } else if (order == tab->order-1) {
667:     if (!tab->bembedt) goto unavailable;
668:     if (ark->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (bet,be) */
669:       VecCopy(ts->vec_sol,X);
670:       for (j=0; j<s; j++) w[j] = h*tab->bembedt[j];
671:       VecMAXPY(X,s,w,ark->YdotI);
672:       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
673:       VecMAXPY(X,s,w,ark->YdotRHS);
674:     } else { /* Rollback and re-complete using (bet-be,be-b) */
675:       VecCopy(ts->vec_sol,X);
676:       for (j=0; j<s; j++) w[j] = h*(tab->bembedt[j] - tab->bt[j]);
677:       VecMAXPY(X,tab->s,w,ark->YdotI);
678:       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
679:       VecMAXPY(X,s,w,ark->YdotRHS);
680:     }
681:     if (done) *done = PETSC_TRUE;
682:     return(0);
683:   }
684: unavailable:
685:   if (done) *done = PETSC_FALSE;
686:   else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"ARKIMEX '%s' of order %D cannot evaluate step at order %D. Consider using -ts_adapt_type none or a different method that has an embedded estimate.",tab->name,tab->order,order);
687:   return(0);
688: }

690: static PetscErrorCode TSRollBack_ARKIMEX(TS ts)
691: {
692:   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
693:   ARKTableau      tab  = ark->tableau;
694:   const PetscInt  s    = tab->s;
695:   const PetscReal *bt  = tab->bt,*b = tab->b;
696:   PetscScalar     *w   = ark->work;
697:   Vec             *YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS;
698:   PetscInt        j;
699:   PetscReal       h;
700:   PetscErrorCode  ierr;

703:   switch (ark->status) {
704:   case TS_STEP_INCOMPLETE:
705:   case TS_STEP_PENDING:
706:     h = ts->time_step; break;
707:   case TS_STEP_COMPLETE:
708:     h = ts->ptime - ts->ptime_prev; break;
709:   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
710:   }
711:   for (j=0; j<s; j++) w[j] = -h*bt[j];
712:   VecMAXPY(ts->vec_sol,s,w,YdotI);
713:   for (j=0; j<s; j++) w[j] = -h*b[j];
714:   VecMAXPY(ts->vec_sol,s,w,YdotRHS);
715:   return(0);
716: }

718: static PetscErrorCode TSStep_ARKIMEX(TS ts)
719: {
720:   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
721:   ARKTableau      tab  = ark->tableau;
722:   const PetscInt  s    = tab->s;
723:   const PetscReal *At  = tab->At,*A = tab->A,*ct = tab->ct,*c = tab->c;
724:   PetscScalar     *w   = ark->work;
725:   Vec             *Y   = ark->Y,*YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS,Ydot = ark->Ydot,Ydot0 = ark->Ydot0,Z = ark->Z;
726:   PetscBool       extrapolate = ark->extrapolate;
727:   TSAdapt         adapt;
728:   SNES            snes;
729:   PetscInt        i,j,its,lits;
730:   PetscInt        rejections = 0;
731:   PetscBool       stageok,accept = PETSC_TRUE;
732:   PetscReal       next_time_step = ts->time_step;
733:   PetscErrorCode  ierr;

736:   if (ark->extrapolate && !ark->Y_prev) {
737:     VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y_prev);
738:     VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI_prev);
739:     VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS_prev);
740:   }

742:   if (!ts->steprollback) {
743:     if (ts->equation_type >= TS_EQ_IMPLICIT) { /* Save the initial slope for the next step */
744:       VecCopy(YdotI[s-1],Ydot0);
745:     }
746:     if (ark->extrapolate && !ts->steprestart) { /* Save the Y, YdotI, YdotRHS for extrapolation initial guess */
747:       for (i = 0; i<s; i++) {
748:         VecCopy(Y[i],ark->Y_prev[i]);
749:         VecCopy(YdotRHS[i],ark->YdotRHS_prev[i]);
750:         VecCopy(YdotI[i],ark->YdotI_prev[i]);
751:       }
752:     }
753:   }

755:   if (ts->equation_type >= TS_EQ_IMPLICIT && tab->explicit_first_stage && ts->steprestart) {
756:     TS ts_start;
757:     TSClone(ts,&ts_start);
758:     TSSetSolution(ts_start,ts->vec_sol);
759:     TSSetTime(ts_start,ts->ptime);
760:     TSSetMaxSteps(ts_start,ts->steps+1);
761:     TSSetMaxTime(ts_start,ts->ptime+ts->time_step);
762:     TSSetExactFinalTime(ts_start,TS_EXACTFINALTIME_STEPOVER);
763:     TSSetTimeStep(ts_start,ts->time_step);
764:     TSSetType(ts_start,TSARKIMEX);
765:     TSARKIMEXSetFullyImplicit(ts_start,PETSC_TRUE);
766:     TSARKIMEXSetType(ts_start,TSARKIMEX1BEE);

768:     TSRestartStep(ts_start);
769:     TSSolve(ts_start,ts->vec_sol);
770:     TSGetTime(ts_start,&ts->ptime);
771:     TSGetTimeStep(ts_start,&ts->time_step);

773:     { /* Save the initial slope for the next step */
774:       TS_ARKIMEX *ark_start = (TS_ARKIMEX*)ts_start->data;
775:       VecCopy(ark_start->YdotI[ark_start->tableau->s-1],Ydot0);
776:     }
777:     ts->steps++;

779:     /* Set the correct TS in SNES */
780:     /* We'll try to bypass this by changing the method on the fly */
781:     {
782:       TSGetSNES(ts,&snes);
783:       TSSetSNES(ts,snes);
784:     }
785:     TSDestroy(&ts_start);
786:   }

788:   ark->status = TS_STEP_INCOMPLETE;
789:   while (!ts->reason && ark->status != TS_STEP_COMPLETE) {
790:     PetscReal t = ts->ptime;
791:     PetscReal h = ts->time_step;
792:     for (i=0; i<s; i++) {
793:       ark->stage_time = t + h*ct[i];
794:       TSPreStage(ts,ark->stage_time);
795:       if (At[i*s+i] == 0) { /* This stage is explicit */
796:         if (i!=0 && ts->equation_type >= TS_EQ_IMPLICIT) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Explicit stages other than the first one are not supported for implicit problems");
797:         VecCopy(ts->vec_sol,Y[i]);
798:         for (j=0; j<i; j++) w[j] = h*At[i*s+j];
799:         VecMAXPY(Y[i],i,w,YdotI);
800:         for (j=0; j<i; j++) w[j] = h*A[i*s+j];
801:         VecMAXPY(Y[i],i,w,YdotRHS);
802:       } else {
803:         ark->scoeff = 1./At[i*s+i];
804:         /* Ydot = shift*(Y-Z) */
805:         VecCopy(ts->vec_sol,Z);
806:         for (j=0; j<i; j++) w[j] = h*At[i*s+j];
807:         VecMAXPY(Z,i,w,YdotI);
808:         for (j=0; j<i; j++) w[j] = h*A[i*s+j];
809:         VecMAXPY(Z,i,w,YdotRHS);
810:         if (extrapolate && !ts->steprestart) {
811:           /* Initial guess extrapolated from previous time step stage values */
812:           TSExtrapolate_ARKIMEX(ts,c[i],Y[i]);
813:         } else {
814:           /* Initial guess taken from last stage */
815:           VecCopy(i>0 ? Y[i-1] : ts->vec_sol,Y[i]);
816:         }
817:         TSGetSNES(ts,&snes);
818:         SNESSolve(snes,NULL,Y[i]);
819:         SNESGetIterationNumber(snes,&its);
820:         SNESGetLinearSolveIterations(snes,&lits);
821:         ts->snes_its += its; ts->ksp_its += lits;
822:         TSGetAdapt(ts,&adapt);
823:         TSAdaptCheckStage(adapt,ts,ark->stage_time,Y[i],&stageok);
824:         if (!stageok) {
825:           /* We are likely rejecting the step because of solver or function domain problems so we should not attempt to
826:            * use extrapolation to initialize the solves on the next attempt. */
827:           extrapolate = PETSC_FALSE;
828:           goto reject_step;
829:         }
830:       }
831:       if (ts->equation_type >= TS_EQ_IMPLICIT) {
832:         if (i==0 && tab->explicit_first_stage) {
833:           if (!tab->stiffly_accurate ) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s is not stiffly accurate and therefore explicit-first stage methods cannot be used if the equation is implicit because the slope cannot be evaluated",ark->tableau->name);
834:           VecCopy(Ydot0,YdotI[0]);                                      /* YdotI = YdotI(tn-1) */
835:         } else {
836:           VecAXPBYPCZ(YdotI[i],-ark->scoeff/h,ark->scoeff/h,0,Z,Y[i]);  /* YdotI = shift*(X-Z) */
837:         }
838:       } else {
839:         if (i==0 && tab->explicit_first_stage) {
840:           VecZeroEntries(Ydot);
841:           TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],ark->imex);/* YdotI = -G(t,Y,0)   */
842:           VecScale(YdotI[i],-1.0);
843:         } else {
844:           VecAXPBYPCZ(YdotI[i],-ark->scoeff/h,ark->scoeff/h,0,Z,Y[i]);  /* YdotI = shift*(X-Z) */
845:         }
846:         if (ark->imex) {
847:           TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);
848:         } else {
849:           VecZeroEntries(YdotRHS[i]);
850:         }
851:       }
852:       TSPostStage(ts,ark->stage_time,i,Y);
853:     }

855:     ark->status = TS_STEP_INCOMPLETE;
856:     TSEvaluateStep_ARKIMEX(ts,tab->order,ts->vec_sol,NULL);
857:     ark->status = TS_STEP_PENDING;
858:     TSGetAdapt(ts,&adapt);
859:     TSAdaptCandidatesClear(adapt);
860:     TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);
861:     TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);
862:     ark->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
863:     if (!accept) { /* Roll back the current step */
864:       TSRollBack_ARKIMEX(ts);
865:       ts->time_step = next_time_step;
866:       goto reject_step;
867:     }

869:     ts->ptime += ts->time_step;
870:     ts->time_step = next_time_step;
871:     break;

873:   reject_step:
874:     ts->reject++; accept = PETSC_FALSE;
875:     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
876:       ts->reason = TS_DIVERGED_STEP_REJECTED;
877:       PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);
878:     }
879:   }
880:   return(0);
881: }

883: static PetscErrorCode TSInterpolate_ARKIMEX(TS ts,PetscReal itime,Vec X)
884: {
885:   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
886:   PetscInt        s    = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j;
887:   PetscReal       h;
888:   PetscReal       tt,t;
889:   PetscScalar     *bt,*b;
890:   const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp;
891:   PetscErrorCode  ierr;

894:   if (!Bt || !B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name);
895:   switch (ark->status) {
896:   case TS_STEP_INCOMPLETE:
897:   case TS_STEP_PENDING:
898:     h = ts->time_step;
899:     t = (itime - ts->ptime)/h;
900:     break;
901:   case TS_STEP_COMPLETE:
902:     h = ts->ptime - ts->ptime_prev;
903:     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
904:     break;
905:   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
906:   }
907:   PetscMalloc2(s,&bt,s,&b);
908:   for (i=0; i<s; i++) bt[i] = b[i] = 0;
909:   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
910:     for (i=0; i<s; i++) {
911:       bt[i] += h * Bt[i*pinterp+j] * tt;
912:       b[i]  += h * B[i*pinterp+j] * tt;
913:     }
914:   }
915:   VecCopy(ark->Y[0],X);
916:   VecMAXPY(X,s,bt,ark->YdotI);
917:   VecMAXPY(X,s,b,ark->YdotRHS);
918:   PetscFree2(bt,b);
919:   return(0);
920: }

922: static PetscErrorCode TSExtrapolate_ARKIMEX(TS ts,PetscReal c,Vec X)
923: {
924:   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
925:   PetscInt        s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j;
926:   PetscReal       h,h_prev,t,tt;
927:   PetscScalar     *bt,*b;
928:   const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp;
929:   PetscErrorCode  ierr;

932:   if (!Bt || !B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name);
933:   PetscCalloc2(s,&bt,s,&b);
934:   h = ts->time_step;
935:   h_prev = ts->ptime - ts->ptime_prev;
936:   t = 1 + h/h_prev*c;
937:   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
938:     for (i=0; i<s; i++) {
939:       bt[i] += h * Bt[i*pinterp+j] * tt;
940:       b[i]  += h * B[i*pinterp+j] * tt;
941:     }
942:   }
943:   if (!ark->Y_prev) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Stages from previous step have not been stored");
944:   VecCopy(ark->Y_prev[0],X);
945:   VecMAXPY(X,s,bt,ark->YdotI_prev);
946:   VecMAXPY(X,s,b,ark->YdotRHS_prev);
947:   PetscFree2(bt,b);
948:   return(0);
949: }

951: /*------------------------------------------------------------*/

953: static PetscErrorCode TSARKIMEXTableauReset(TS ts)
954: {
955:   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
956:   ARKTableau     tab  = ark->tableau;

960:   if (!tab) return(0);
961:   PetscFree(ark->work);
962:   VecDestroyVecs(tab->s,&ark->Y);
963:   VecDestroyVecs(tab->s,&ark->YdotI);
964:   VecDestroyVecs(tab->s,&ark->YdotRHS);
965:   VecDestroyVecs(tab->s,&ark->Y_prev);
966:   VecDestroyVecs(tab->s,&ark->YdotI_prev);
967:   VecDestroyVecs(tab->s,&ark->YdotRHS_prev);
968:   return(0);
969: }

971: static PetscErrorCode TSReset_ARKIMEX(TS ts)
972: {
973:   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;

977:   TSARKIMEXTableauReset(ts);
978:   VecDestroy(&ark->Ydot);
979:   VecDestroy(&ark->Ydot0);
980:   VecDestroy(&ark->Z);
981:   return(0);
982: }

984: static PetscErrorCode TSARKIMEXGetVecs(TS ts,DM dm,Vec *Z,Vec *Ydot)
985: {
986:   TS_ARKIMEX     *ax = (TS_ARKIMEX*)ts->data;

990:   if (Z) {
991:     if (dm && dm != ts->dm) {
992:       DMGetNamedGlobalVector(dm,"TSARKIMEX_Z",Z);
993:     } else *Z = ax->Z;
994:   }
995:   if (Ydot) {
996:     if (dm && dm != ts->dm) {
997:       DMGetNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);
998:     } else *Ydot = ax->Ydot;
999:   }
1000:   return(0);
1001: }


1004: static PetscErrorCode TSARKIMEXRestoreVecs(TS ts,DM dm,Vec *Z,Vec *Ydot)
1005: {

1009:   if (Z) {
1010:     if (dm && dm != ts->dm) {
1011:       DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Z",Z);
1012:     }
1013:   }
1014:   if (Ydot) {
1015:     if (dm && dm != ts->dm) {
1016:       DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);
1017:     }
1018:   }
1019:   return(0);
1020: }

1022: /*
1023:   This defines the nonlinear equation that is to be solved with SNES
1024:   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
1025: */
1026: static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec F,TS ts)
1027: {
1028:   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1029:   DM             dm,dmsave;
1030:   Vec            Z,Ydot;
1031:   PetscReal      shift = ark->scoeff / ts->time_step;

1035:   SNESGetDM(snes,&dm);
1036:   TSARKIMEXGetVecs(ts,dm,&Z,&Ydot);
1037:   VecAXPBYPCZ(Ydot,-shift,shift,0,Z,X); /* Ydot = shift*(X-Z) */
1038:   dmsave = ts->dm;
1039:   ts->dm = dm;

1041:   TSComputeIFunction(ts,ark->stage_time,X,Ydot,F,ark->imex);

1043:   ts->dm = dmsave;
1044:   TSARKIMEXRestoreVecs(ts,dm,&Z,&Ydot);
1045:   return(0);
1046: }

1048: static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes,Vec X,Mat A,Mat B,TS ts)
1049: {
1050:   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1051:   DM             dm,dmsave;
1052:   Vec            Ydot;
1053:   PetscReal      shift = ark->scoeff / ts->time_step;

1057:   SNESGetDM(snes,&dm);
1058:   TSARKIMEXGetVecs(ts,dm,NULL,&Ydot);
1059:   /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */
1060:   dmsave = ts->dm;
1061:   ts->dm = dm;

1063:   TSComputeIJacobian(ts,ark->stage_time,X,Ydot,shift,A,B,ark->imex);

1065:   ts->dm = dmsave;
1066:   TSARKIMEXRestoreVecs(ts,dm,NULL,&Ydot);
1067:   return(0);
1068: }

1070: static PetscErrorCode DMCoarsenHook_TSARKIMEX(DM fine,DM coarse,void *ctx)
1071: {
1073:   return(0);
1074: }

1076: static PetscErrorCode DMRestrictHook_TSARKIMEX(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
1077: {
1078:   TS             ts = (TS)ctx;
1080:   Vec            Z,Z_c;

1083:   TSARKIMEXGetVecs(ts,fine,&Z,NULL);
1084:   TSARKIMEXGetVecs(ts,coarse,&Z_c,NULL);
1085:   MatRestrict(restrct,Z,Z_c);
1086:   VecPointwiseMult(Z_c,rscale,Z_c);
1087:   TSARKIMEXRestoreVecs(ts,fine,&Z,NULL);
1088:   TSARKIMEXRestoreVecs(ts,coarse,&Z_c,NULL);
1089:   return(0);
1090: }


1093: static PetscErrorCode DMSubDomainHook_TSARKIMEX(DM dm,DM subdm,void *ctx)
1094: {
1096:   return(0);
1097: }

1099: static PetscErrorCode DMSubDomainRestrictHook_TSARKIMEX(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
1100: {
1101:   TS             ts = (TS)ctx;
1103:   Vec            Z,Z_c;

1106:   TSARKIMEXGetVecs(ts,dm,&Z,NULL);
1107:   TSARKIMEXGetVecs(ts,subdm,&Z_c,NULL);

1109:   VecScatterBegin(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);
1110:   VecScatterEnd(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);

1112:   TSARKIMEXRestoreVecs(ts,dm,&Z,NULL);
1113:   TSARKIMEXRestoreVecs(ts,subdm,&Z_c,NULL);
1114:   return(0);
1115: }

1117: static PetscErrorCode TSARKIMEXTableauSetUp(TS ts)
1118: {
1119:   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1120:   ARKTableau     tab  = ark->tableau;

1124:   PetscMalloc1(tab->s,&ark->work);
1125:   VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y);
1126:   VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI);
1127:   VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS);
1128:   if (ark->extrapolate) {
1129:     VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y_prev);
1130:     VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI_prev);
1131:     VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS_prev);
1132:   }
1133:   return(0);
1134: }

1136: static PetscErrorCode TSSetUp_ARKIMEX(TS ts)
1137: {
1138:   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1140:   DM             dm;
1141:   SNES           snes;

1144:   TSARKIMEXTableauSetUp(ts);
1145:   VecDuplicate(ts->vec_sol,&ark->Ydot);
1146:   VecDuplicate(ts->vec_sol,&ark->Ydot0);
1147:   VecDuplicate(ts->vec_sol,&ark->Z);
1148:   TSGetDM(ts,&dm);
1149:   DMCoarsenHookAdd(dm,DMCoarsenHook_TSARKIMEX,DMRestrictHook_TSARKIMEX,ts);
1150:   DMSubDomainHookAdd(dm,DMSubDomainHook_TSARKIMEX,DMSubDomainRestrictHook_TSARKIMEX,ts);
1151:   TSGetSNES(ts,&snes);
1152:   return(0);
1153: }
1154: /*------------------------------------------------------------*/

1156: static PetscErrorCode TSSetFromOptions_ARKIMEX(PetscOptionItems *PetscOptionsObject,TS ts)
1157: {
1158:   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;

1162:   PetscOptionsHead(PetscOptionsObject,"ARKIMEX ODE solver options");
1163:   {
1164:     ARKTableauLink link;
1165:     PetscInt       count,choice;
1166:     PetscBool      flg;
1167:     const char     **namelist;
1168:     for (link=ARKTableauList,count=0; link; link=link->next,count++) ;
1169:     PetscMalloc1(count,(char***)&namelist);
1170:     for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
1171:     PetscOptionsEList("-ts_arkimex_type","Family of ARK IMEX method","TSARKIMEXSetType",(const char*const*)namelist,count,ark->tableau->name,&choice,&flg);
1172:     if (flg) {TSARKIMEXSetType(ts,namelist[choice]);}
1173:     PetscFree(namelist);

1175:     flg  = (PetscBool) !ark->imex;
1176:     PetscOptionsBool("-ts_arkimex_fully_implicit","Solve the problem fully implicitly","TSARKIMEXSetFullyImplicit",flg,&flg,NULL);
1177:     ark->imex = (PetscBool) !flg;
1178:     PetscOptionsBool("-ts_arkimex_initial_guess_extrapolate","Extrapolate the initial guess for the stage solution from stage values of the previous time step","",ark->extrapolate,&ark->extrapolate,NULL);
1179:   }
1180:   PetscOptionsTail();
1181:   return(0);
1182: }

1184: static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
1185: {
1187:   PetscInt       i;
1188:   size_t         left,count;
1189:   char           *p;

1192:   for (i=0,p=buf,left=len; i<n; i++) {
1193:     PetscSNPrintfCount(p,left,fmt,&count,(double)x[i]);
1194:     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
1195:     left -= count;
1196:     p    += count;
1197:     *p++  = ' ';
1198:   }
1199:   p[i ? 0 : -1] = 0;
1200:   return(0);
1201: }

1203: static PetscErrorCode TSView_ARKIMEX(TS ts,PetscViewer viewer)
1204: {
1205:   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1206:   PetscBool      iascii;

1210:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1211:   if (iascii) {
1212:     ARKTableau    tab = ark->tableau;
1213:     TSARKIMEXType arktype;
1214:     char          buf[512];
1215:     TSARKIMEXGetType(ts,&arktype);
1216:     PetscViewerASCIIPrintf(viewer,"  ARK IMEX %s\n",arktype);
1217:     PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->ct);
1218:     PetscViewerASCIIPrintf(viewer,"  Stiff abscissa       ct = %s\n",buf);
1219:     PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);
1220:     PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s\n",tab->stiffly_accurate ? "yes" : "no");
1221:     PetscViewerASCIIPrintf(viewer,"Explicit first stage: %s\n",tab->explicit_first_stage ? "yes" : "no");
1222:     PetscViewerASCIIPrintf(viewer,"FSAL property: %s\n",tab->FSAL_implicit ? "yes" : "no");
1223:     PetscViewerASCIIPrintf(viewer,"  Nonstiff abscissa     c = %s\n",buf);
1224:   }
1225:   return(0);
1226: }

1228: static PetscErrorCode TSLoad_ARKIMEX(TS ts,PetscViewer viewer)
1229: {
1231:   SNES           snes;
1232:   TSAdapt        adapt;

1235:   TSGetAdapt(ts,&adapt);
1236:   TSAdaptLoad(adapt,viewer);
1237:   TSGetSNES(ts,&snes);
1238:   SNESLoad(snes,viewer);
1239:   /* function and Jacobian context for SNES when used with TS is always ts object */
1240:   SNESSetFunction(snes,NULL,NULL,ts);
1241:   SNESSetJacobian(snes,NULL,NULL,NULL,ts);
1242:   return(0);
1243: }

1245: /*@C
1246:   TSARKIMEXSetType - Set the type of ARK IMEX scheme

1248:   Logically collective

1250:   Input Parameter:
1251: +  ts - timestepping context
1252: -  arktype - type of ARK-IMEX scheme

1254:   Options Database:
1255: .  -ts_arkimex_type <1bee,a2,l2,ars122,2c,2d,2e,prssp2,3,bpr3,ars443,4,5>

1257:   Level: intermediate

1259: .seealso: TSARKIMEXGetType(), TSARKIMEX, TSARKIMEXType, TSARKIMEX1BEE, TSARKIMEXA2, TSARKIMEXL2, TSARKIMEXARS122, TSARKIMEX2C, TSARKIMEX2D, TSARKIMEX2E, TSARKIMEXPRSSP2,
1260:           TSARKIMEX3, TSARKIMEXBPR3, TSARKIMEXARS443, TSARKIMEX4, TSARKIMEX5
1261: @*/
1262: PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType arktype)
1263: {

1269:   PetscTryMethod(ts,"TSARKIMEXSetType_C",(TS,TSARKIMEXType),(ts,arktype));
1270:   return(0);
1271: }

1273: /*@C
1274:   TSARKIMEXGetType - Get the type of ARK IMEX scheme

1276:   Logically collective

1278:   Input Parameter:
1279: .  ts - timestepping context

1281:   Output Parameter:
1282: .  arktype - type of ARK-IMEX scheme

1284:   Level: intermediate

1286: .seealso: TSARKIMEXGetType()
1287: @*/
1288: PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType *arktype)
1289: {

1294:   PetscUseMethod(ts,"TSARKIMEXGetType_C",(TS,TSARKIMEXType*),(ts,arktype));
1295:   return(0);
1296: }

1298: /*@
1299:   TSARKIMEXSetFullyImplicit - Solve both parts of the equation implicitly

1301:   Logically collective

1303:   Input Parameter:
1304: +  ts - timestepping context
1305: -  flg - PETSC_TRUE for fully implicit

1307:   Level: intermediate

1309: .seealso: TSARKIMEXGetType()
1310: @*/
1311: PetscErrorCode TSARKIMEXSetFullyImplicit(TS ts,PetscBool flg)
1312: {

1317:   PetscTryMethod(ts,"TSARKIMEXSetFullyImplicit_C",(TS,PetscBool),(ts,flg));
1318:   return(0);
1319: }

1321: static PetscErrorCode  TSARKIMEXGetType_ARKIMEX(TS ts,TSARKIMEXType *arktype)
1322: {
1323:   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;

1326:   *arktype = ark->tableau->name;
1327:   return(0);
1328: }
1329: static PetscErrorCode  TSARKIMEXSetType_ARKIMEX(TS ts,TSARKIMEXType arktype)
1330: {
1331:   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1333:   PetscBool      match;
1334:   ARKTableauLink link;

1337:   if (ark->tableau) {
1338:     PetscStrcmp(ark->tableau->name,arktype,&match);
1339:     if (match) return(0);
1340:   }
1341:   for (link = ARKTableauList; link; link=link->next) {
1342:     PetscStrcmp(link->tab.name,arktype,&match);
1343:     if (match) {
1344:       if (ts->setupcalled) {TSARKIMEXTableauReset(ts);}
1345:       ark->tableau = &link->tab;
1346:       if (ts->setupcalled) {TSARKIMEXTableauSetUp(ts);}
1347:       ts->default_adapt_type = ark->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
1348:       return(0);
1349:     }
1350:   }
1351:   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype);
1352:   return(0);
1353: }

1355: static PetscErrorCode  TSARKIMEXSetFullyImplicit_ARKIMEX(TS ts,PetscBool flg)
1356: {
1357:   TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;

1360:   ark->imex = (PetscBool)!flg;
1361:   return(0);
1362: }

1364: static PetscErrorCode TSDestroy_ARKIMEX(TS ts)
1365: {

1369:   TSReset_ARKIMEX(ts);
1370:   if (ts->dm) {
1371:     DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSARKIMEX,DMRestrictHook_TSARKIMEX,ts);
1372:     DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSARKIMEX,DMSubDomainRestrictHook_TSARKIMEX,ts);
1373:   }
1374:   PetscFree(ts->data);
1375:   PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C",NULL);
1376:   PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C",NULL);
1377:   PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",NULL);
1378:   return(0);
1379: }

1381: /* ------------------------------------------------------------ */
1382: /*MC
1383:       TSARKIMEX - ODE and DAE solver using Additive Runge-Kutta IMEX schemes

1385:   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
1386:   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
1387:   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().

1389:   Notes:
1390:   The default is TSARKIMEX3, it can be changed with TSARKIMEXSetType() or -ts_arkimex_type

1392:   If the equation is implicit or a DAE, then TSSetEquationType() needs to be set accordingly. Refer to the manual for further information.

1394:   Methods with an explicit stage can only be used with ODE in which the stiff part G(t,X,Xdot) has the form Xdot + Ghat(t,X).

1396:   Consider trying TSROSW if the stiff part is linear or weakly nonlinear.

1398:   Level: beginner

1400: .seealso:  TSCreate(), TS, TSSetType(), TSARKIMEXSetType(), TSARKIMEXGetType(), TSARKIMEXSetFullyImplicit(), TSARKIMEX1BEE,
1401:            TSARKIMEX2C, TSARKIMEX2D, TSARKIMEX2E, TSARKIMEX3, TSARKIMEXL2, TSARKIMEXA2, TSARKIMEXARS122,
1402:            TSARKIMEX4, TSARKIMEX5, TSARKIMEXPRSSP2, TSARKIMEXARS443, TSARKIMEXBPR3, TSARKIMEXType, TSARKIMEXRegister()

1404: M*/
1405: PETSC_EXTERN PetscErrorCode TSCreate_ARKIMEX(TS ts)
1406: {
1407:   TS_ARKIMEX     *th;

1411:   TSARKIMEXInitializePackage();

1413:   ts->ops->reset          = TSReset_ARKIMEX;
1414:   ts->ops->destroy        = TSDestroy_ARKIMEX;
1415:   ts->ops->view           = TSView_ARKIMEX;
1416:   ts->ops->load           = TSLoad_ARKIMEX;
1417:   ts->ops->setup          = TSSetUp_ARKIMEX;
1418:   ts->ops->step           = TSStep_ARKIMEX;
1419:   ts->ops->interpolate    = TSInterpolate_ARKIMEX;
1420:   ts->ops->evaluatestep   = TSEvaluateStep_ARKIMEX;
1421:   ts->ops->rollback       = TSRollBack_ARKIMEX;
1422:   ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX;
1423:   ts->ops->snesfunction   = SNESTSFormFunction_ARKIMEX;
1424:   ts->ops->snesjacobian   = SNESTSFormJacobian_ARKIMEX;

1426:   ts->usessnes = PETSC_TRUE;

1428:   PetscNewLog(ts,&th);
1429:   ts->data = (void*)th;
1430:   th->imex = PETSC_TRUE;

1432:   PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C",TSARKIMEXGetType_ARKIMEX);
1433:   PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C",TSARKIMEXSetType_ARKIMEX);
1434:   PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",TSARKIMEXSetFullyImplicit_ARKIMEX);

1436:   TSARKIMEXSetType(ts,TSARKIMEXDefault);
1437:   return(0);
1438: }