Actual source code: rk.c

petsc-3.10.5 2019-03-28
Report Typos and Errors
  1: /*
  2:   Code for time stepping with the Runge-Kutta method

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

  7:   Udot = F(t,U)

  9: */
 10:  #include <petsc/private/tsimpl.h>
 11:  #include <petscdm.h>

 13: static TSRKType  TSRKDefault = TSRK3BS;
 14: static PetscBool TSRKRegisterAllCalled;
 15: static PetscBool TSRKPackageInitialized;

 17: typedef struct _RKTableau *RKTableau;
 18: struct _RKTableau {
 19:   char      *name;
 20:   PetscInt   order;               /* Classical approximation order of the method i              */
 21:   PetscInt   s;                   /* Number of stages                                           */
 22:   PetscInt   p;                   /* Interpolation order                                        */
 23:   PetscBool  FSAL;                /* flag to indicate if tableau is FSAL                        */
 24:   PetscReal *A,*b,*c;             /* Tableau                                                    */
 25:   PetscReal *bembed;              /* Embedded formula of order one less (order-1)               */
 26:   PetscReal *binterp;             /* Dense output formula                                       */
 27:   PetscReal  ccfl;                /* Placeholder for CFL coefficient relative to forward Euler  */
 28: };
 29: typedef struct _RKTableauLink *RKTableauLink;
 30: struct _RKTableauLink {
 31:   struct _RKTableau tab;
 32:   RKTableauLink     next;
 33: };
 34: static RKTableauLink RKTableauList;

 36: typedef struct {
 37:   RKTableau    tableau;
 38:   Vec          *Y;               /* States computed during the step */
 39:   Vec          *YdotRHS;         /* Function evaluations for the non-stiff part */
 40:   Vec          *VecDeltaLam;     /* Increment of the adjoint sensitivity w.r.t IC at stage */
 41:   Vec          *VecDeltaMu;      /* Increment of the adjoint sensitivity w.r.t P at stage */
 42:   Vec          VecCostIntegral0; /* backup for roll-backs due to events */
 43:   PetscScalar  *work;            /* Scalar work */
 44:   PetscReal    stage_time;
 45:   TSStepStatus status;
 46:   PetscReal    ptime;
 47:   PetscReal    time_step;
 48: } TS_RK;

 50: /*MC
 51:      TSRK1FE - First order forward Euler scheme.

 53:      This method has one stage.

 55:      Options database:
 56: .     -ts_rk_type 1fe

 58:      Level: advanced

 60: .seealso: TSRK, TSRKType, TSRKSetType()
 61: M*/
 62: /*MC
 63:      TSRK2A - Second order RK scheme.

 65:      This method has two stages.

 67:      Options database:
 68: .     -ts_rk_type 2a

 70:      Level: advanced

 72: .seealso: TSRK, TSRKType, TSRKSetType()
 73: M*/
 74: /*MC
 75:      TSRK3 - Third order RK scheme.

 77:      This method has three stages.

 79:      Options database:
 80: .     -ts_rk_type 3

 82:      Level: advanced

 84: .seealso: TSRK, TSRKType, TSRKSetType()
 85: M*/
 86: /*MC
 87:      TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method.

 89:      This method has four stages with the First Same As Last (FSAL) property.

 91:      Options database:
 92: .     -ts_rk_type 3bs

 94:      Level: advanced

 96:      References: https://doi.org/10.1016/0893-9659(89)90079-7

 98: .seealso: TSRK, TSRKType, TSRKSetType()
 99: M*/
100: /*MC
101:      TSRK4 - Fourth order RK scheme.

103:      This is the classical Runge-Kutta method with four stages.

105:      Options database:
106: .     -ts_rk_type 4

108:      Level: advanced

110: .seealso: TSRK, TSRKType, TSRKSetType()
111: M*/
112: /*MC
113:      TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method.

115:      This method has six stages.

117:      Options database:
118: .     -ts_rk_type 5f

120:      Level: advanced

122: .seealso: TSRK, TSRKType, TSRKSetType()
123: M*/
124: /*MC
125:      TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method.

127:      This method has seven stages with the First Same As Last (FSAL) property.

129:      Options database:
130: .     -ts_rk_type 5dp

132:      Level: advanced

134:      References: https://doi.org/10.1016/0771-050X(80)90013-3

136: .seealso: TSRK, TSRKType, TSRKSetType()
137: M*/
138: /*MC
139:      TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method.

141:      This method has eight stages with the First Same As Last (FSAL) property.

143:      Options database:
144: .     -ts_rk_type 5bs

146:      Level: advanced

148:      References: https://doi.org/10.1016/0898-1221(96)00141-1

150: .seealso: TSRK, TSRKType, TSRKSetType()
151: M*/

153: /*@C
154:   TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK

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

158:   Level: advanced

160: .keywords: TS, TSRK, register, all

162: .seealso:  TSRKRegisterDestroy()
163: @*/
164: PetscErrorCode TSRKRegisterAll(void)
165: {

169:   if (TSRKRegisterAllCalled) return(0);
170:   TSRKRegisterAllCalled = PETSC_TRUE;

172: #define RC PetscRealConstant
173:   {
174:     const PetscReal
175:       A[1][1] = {{0}},
176:       b[1]    = {RC(1.0)};
177:     TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL);
178:   }
179:   {
180:     const PetscReal
181:       A[2][2]   = {{0,0},
182:                    {RC(1.0),0}},
183:       b[2]      =  {RC(0.5),RC(0.5)},
184:       bembed[2] =  {RC(1.0),0};
185:     TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL);
186:   }
187:   {
188:     const PetscReal
189:       A[3][3] = {{0,0,0},
190:                  {RC(2.0)/RC(3.0),0,0},
191:                  {RC(-1.0)/RC(3.0),RC(1.0),0}},
192:       b[3]    =  {RC(0.25),RC(0.5),RC(0.25)};
193:     TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);
194:   }
195:   {
196:     const PetscReal
197:       A[4][4]   = {{0,0,0,0},
198:                    {RC(1.0)/RC(2.0),0,0,0},
199:                    {0,RC(3.0)/RC(4.0),0,0},
200:                    {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}},
201:       b[4]      =  {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0},
202:       bembed[4] =  {RC(7.0)/RC(24.0),RC(1.0)/RC(4.0),RC(1.0)/RC(3.0),RC(1.0)/RC(8.0)};
203:     TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);
204:   }
205:   {
206:     const PetscReal
207:       A[4][4] = {{0,0,0,0},
208:                  {RC(0.5),0,0,0},
209:                  {0,RC(0.5),0,0},
210:                  {0,0,RC(1.0),0}},
211:       b[4]    =  {RC(1.0)/RC(6.0),RC(1.0)/RC(3.0),RC(1.0)/RC(3.0),RC(1.0)/RC(6.0)};
212:     TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);
213:   }
214:   {
215:     const PetscReal
216:       A[6][6]   = {{0,0,0,0,0,0},
217:                    {RC(0.25),0,0,0,0,0},
218:                    {RC(3.0)/RC(32.0),RC(9.0)/RC(32.0),0,0,0,0},
219:                    {RC(1932.0)/RC(2197.0),RC(-7200.0)/RC(2197.0),RC(7296.0)/RC(2197.0),0,0,0},
220:                    {RC(439.0)/RC(216.0),RC(-8.0),RC(3680.0)/RC(513.0),RC(-845.0)/RC(4104.0),0,0},
221:                    {RC(-8.0)/RC(27.0),RC(2.0),RC(-3544.0)/RC(2565.0),RC(1859.0)/RC(4104.0),RC(-11.0)/RC(40.0),0}},
222:       b[6]      =  {RC(16.0)/RC(135.0),0,RC(6656.0)/RC(12825.0),RC(28561.0)/RC(56430.0),RC(-9.0)/RC(50.0),RC(2.0)/RC(55.0)},
223:       bembed[6] =  {RC(25.0)/RC(216.0),0,RC(1408.0)/RC(2565.0),RC(2197.0)/RC(4104.0),RC(-1.0)/RC(5.0),0};
224:     TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);
225:   }
226:   {
227:     const PetscReal
228:       A[7][7]   = {{0,0,0,0,0,0,0},
229:                    {RC(1.0)/RC(5.0),0,0,0,0,0,0},
230:                    {RC(3.0)/RC(40.0),RC(9.0)/RC(40.0),0,0,0,0,0},
231:                    {RC(44.0)/RC(45.0),RC(-56.0)/RC(15.0),RC(32.0)/RC(9.0),0,0,0,0},
232:                    {RC(19372.0)/RC(6561.0),RC(-25360.0)/RC(2187.0),RC(64448.0)/RC(6561.0),RC(-212.0)/RC(729.0),0,0,0},
233:                    {RC(9017.0)/RC(3168.0),RC(-355.0)/RC(33.0),RC(46732.0)/RC(5247.0),RC(49.0)/RC(176.0),RC(-5103.0)/RC(18656.0),0,0},
234:                    {RC(35.0)/RC(384.0),0,RC(500.0)/RC(1113.0),RC(125.0)/RC(192.0),RC(-2187.0)/RC(6784.0),RC(11.0)/RC(84.0),0}},
235:       b[7]      =  {RC(35.0)/RC(384.0),0,RC(500.0)/RC(1113.0),RC(125.0)/RC(192.0),RC(-2187.0)/RC(6784.0),RC(11.0)/RC(84.0),0},
236:       bembed[7] =  {RC(5179.0)/RC(57600.0),0,RC(7571.0)/RC(16695.0),RC(393.0)/RC(640.0),RC(-92097.0)/RC(339200.0),RC(187.0)/RC(2100.0),RC(1.0)/RC(40.0)};
237:     TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,0,NULL);
238:   }
239:   {
240:     const PetscReal
241:       A[8][8]   = {{0,0,0,0,0,0,0,0},
242:                    {RC(1.0)/RC(6.0),0,0,0,0,0,0,0},
243:                    {RC(2.0)/RC(27.0),RC(4.0)/RC(27.0),0,0,0,0,0,0},
244:                    {RC(183.0)/RC(1372.0),RC(-162.0)/RC(343.0),RC(1053.0)/RC(1372.0),0,0,0,0,0},
245:                    {RC(68.0)/RC(297.0),RC(-4.0)/RC(11.0),RC(42.0)/RC(143.0),RC(1960.0)/RC(3861.0),0,0,0,0},
246:                    {RC(597.0)/RC(22528.0),RC(81.0)/RC(352.0),RC(63099.0)/RC(585728.0),RC(58653.0)/RC(366080.0),RC(4617.0)/RC(20480.0),0,0,0},
247:                    {RC(174197.0)/RC(959244.0),RC(-30942.0)/RC(79937.0),RC(8152137.0)/RC(19744439.0),RC(666106.0)/RC(1039181.0),RC(-29421.0)/RC(29068.0),RC(482048.0)/RC(414219.0),0,0},
248:                    {RC(587.0)/RC(8064.0),0,RC(4440339.0)/RC(15491840.0),RC(24353.0)/RC(124800.0),RC(387.0)/RC(44800.0),RC(2152.0)/RC(5985.0),RC(7267.0)/RC(94080.0),0}},
249:       b[8]      =  {RC(587.0)/RC(8064.0),0,RC(4440339.0)/RC(15491840.0),RC(24353.0)/RC(124800.0),RC(387.0)/RC(44800.0),RC(2152.0)/RC(5985.0),RC(7267.0)/RC(94080.0),0},
250:       bembed[8] =  {RC(2479.0)/RC(34992.0),0,RC(123.0)/RC(416.0),RC(612941.0)/RC(3411720.0),RC(43.0)/RC(1440.0),RC(2272.0)/RC(6561.0),RC(79937.0)/RC(1113912.0),RC(3293.0)/RC(556956.0)};
251:     TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);
252:   }
253: #undef RC
254:   return(0);
255: }

257: /*@C
258:    TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister().

260:    Not Collective

262:    Level: advanced

264: .keywords: TSRK, register, destroy
265: .seealso: TSRKRegister(), TSRKRegisterAll()
266: @*/
267: PetscErrorCode TSRKRegisterDestroy(void)
268: {
270:   RKTableauLink link;

273:   while ((link = RKTableauList)) {
274:     RKTableau t = &link->tab;
275:     RKTableauList = link->next;
276:     PetscFree3(t->A,t->b,t->c);
277:     PetscFree (t->bembed);
278:     PetscFree (t->binterp);
279:     PetscFree (t->name);
280:     PetscFree (link);
281:   }
282:   TSRKRegisterAllCalled = PETSC_FALSE;
283:   return(0);
284: }

286: /*@C
287:   TSRKInitializePackage - This function initializes everything in the TSRK package. It is called
288:   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RK()
289:   when using static libraries.

291:   Level: developer

293: .keywords: TS, TSRK, initialize, package
294: .seealso: PetscInitialize()
295: @*/
296: PetscErrorCode TSRKInitializePackage(void)
297: {

301:   if (TSRKPackageInitialized) return(0);
302:   TSRKPackageInitialized = PETSC_TRUE;
303:   TSRKRegisterAll();
304:   PetscRegisterFinalize(TSRKFinalizePackage);
305:   return(0);
306: }

308: /*@C
309:   TSRKFinalizePackage - This function destroys everything in the TSRK package. It is
310:   called from PetscFinalize().

312:   Level: developer

314: .keywords: Petsc, destroy, package
315: .seealso: PetscFinalize()
316: @*/
317: PetscErrorCode TSRKFinalizePackage(void)
318: {

322:   TSRKPackageInitialized = PETSC_FALSE;
323:   TSRKRegisterDestroy();
324:   return(0);
325: }

327: /*@C
328:    TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation

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

332:    Input Parameters:
333: +  name - identifier for method
334: .  order - approximation order of method
335: .  s - number of stages, this is the dimension of the matrices below
336: .  A - stage coefficients (dimension s*s, row-major)
337: .  b - step completion table (dimension s; NULL to use last row of A)
338: .  c - abscissa (dimension s; NULL to use row sums of A)
339: .  bembed - completion table for embedded method (dimension s; NULL if not available)
340: .  p - Order of the interpolation scheme, equal to the number of columns of binterp
341: -  binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1)

343:    Notes:
344:    Several RK methods are provided, this function is only needed to create new methods.

346:    Level: advanced

348: .keywords: TS, register

350: .seealso: TSRK
351: @*/
352: PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s,
353:                             const PetscReal A[],const PetscReal b[],const PetscReal c[],
354:                             const PetscReal bembed[],PetscInt p,const PetscReal binterp[])
355: {
356:   PetscErrorCode  ierr;
357:   RKTableauLink   link;
358:   RKTableau       t;
359:   PetscInt        i,j;


369:   TSRKInitializePackage();
370:   PetscNew(&link);
371:   t = &link->tab;

373:   PetscStrallocpy(name,&t->name);
374:   t->order = order;
375:   t->s = s;
376:   PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);
377:   PetscMemcpy(t->A,A,s*s*sizeof(A[0]));
378:   if (b)  { PetscMemcpy(t->b,b,s*sizeof(b[0])); }
379:   else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i];
380:   if (c)  { PetscMemcpy(t->c,c,s*sizeof(c[0])); }
381:   else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j];
382:   t->FSAL = PETSC_TRUE;
383:   for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE;

385:   if (bembed) {
386:     PetscMalloc1(s,&t->bembed);
387:     PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));
388:   }

390:   if (!binterp) { p = 1; binterp = t->b; }
391:   t->p = p;
392:   PetscMalloc1(s*p,&t->binterp);
393:   PetscMemcpy(t->binterp,binterp,s*p*sizeof(binterp[0]));

395:   link->next = RKTableauList;
396:   RKTableauList = link;
397:   return(0);
398: }

400: /*
401:  The step completion formula is

403:  x1 = x0 + h b^T YdotRHS

405:  This function can be called before or after ts->vec_sol has been updated.
406:  Suppose we have a completion formula (b) and an embedded formula (be) of different order.
407:  We can write

409:  x1e = x0 + h be^T YdotRHS
410:      = x1 - h b^T YdotRHS + h be^T YdotRHS
411:      = x1 + h (be - b)^T YdotRHS

413:  so we can evaluate the method with different order even after the step has been optimistically completed.
414: */
415: static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done)
416: {
417:   TS_RK         *rk   = (TS_RK*)ts->data;
418:   RKTableau      tab  = rk->tableau;
419:   PetscScalar   *w    = rk->work;
420:   PetscReal      h;
421:   PetscInt       s    = tab->s,j;

425:   switch (rk->status) {
426:   case TS_STEP_INCOMPLETE:
427:   case TS_STEP_PENDING:
428:     h = ts->time_step; break;
429:   case TS_STEP_COMPLETE:
430:     h = ts->ptime - ts->ptime_prev; break;
431:   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
432:   }
433:   if (order == tab->order) {
434:     if (rk->status == TS_STEP_INCOMPLETE) {
435:       VecCopy(ts->vec_sol,X);
436:       for (j=0; j<s; j++) w[j] = h*tab->b[j];
437:       VecMAXPY(X,s,w,rk->YdotRHS);
438:     } else {VecCopy(ts->vec_sol,X);}
439:     return(0);
440:   } else if (order == tab->order-1) {
441:     if (!tab->bembed) goto unavailable;
442:     if (rk->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (be) */
443:       VecCopy(ts->vec_sol,X);
444:       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
445:       VecMAXPY(X,s,w,rk->YdotRHS);
446:     } else { /* Rollback and re-complete using (be-b) */
447:       VecCopy(ts->vec_sol,X);
448:       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
449:       VecMAXPY(X,s,w,rk->YdotRHS);
450:     }
451:     if (done) *done = PETSC_TRUE;
452:     return(0);
453:   }
454: unavailable:
455:   if (done) *done = PETSC_FALSE;
456:   else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"RK '%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);
457:   return(0);
458: }

460: static PetscErrorCode TSForwardCostIntegral_RK(TS ts)
461: {
462:   TS_RK           *rk = (TS_RK*)ts->data;
463:   RKTableau       tab = rk->tableau;
464:   const PetscInt  s = tab->s;
465:   const PetscReal *b = tab->b,*c = tab->c;
466:   Vec             *Y = rk->Y;
467:   PetscInt        i;
468:   PetscErrorCode  ierr;

471:   /* backup cost integral */
472:   VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);
473:   for (i=s-1; i>=0; i--) {
474:     /* Evolve ts->vec_costintegral to compute integrals */
475:     TSComputeCostIntegrand(ts,rk->ptime+rk->time_step*c[i],Y[i],ts->vec_costintegrand);
476:     VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);
477:   }
478:   return(0);
479: }

481: static PetscErrorCode TSAdjointCostIntegral_RK(TS ts)
482: {
483:   TS_RK           *rk = (TS_RK*)ts->data;
484:   RKTableau       tab = rk->tableau;
485:   const PetscInt  s = tab->s;
486:   const PetscReal *b = tab->b,*c = tab->c;
487:   Vec             *Y = rk->Y;
488:   PetscInt        i;
489:   PetscErrorCode  ierr;

492:   for (i=s-1; i>=0; i--) {
493:     /* Evolve ts->vec_costintegral to compute integrals */
494:     TSComputeCostIntegrand(ts,ts->ptime+ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);
495:     VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);
496:   }
497:   return(0);
498: }

500: static PetscErrorCode TSRollBack_RK(TS ts)
501: {
502:   TS_RK           *rk = (TS_RK*)ts->data;
503:   RKTableau       tab = rk->tableau;
504:   const PetscInt  s  = tab->s;
505:   const PetscReal *b = tab->b;
506:   PetscScalar     *w = rk->work;
507:   Vec             *YdotRHS = rk->YdotRHS;
508:   PetscInt        j;
509:   PetscReal       h;
510:   PetscErrorCode  ierr;

513:   switch (rk->status) {
514:   case TS_STEP_INCOMPLETE:
515:   case TS_STEP_PENDING:
516:     h = ts->time_step; break;
517:   case TS_STEP_COMPLETE:
518:     h = ts->ptime - ts->ptime_prev; break;
519:   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
520:   }
521:   for (j=0; j<s; j++) w[j] = -h*b[j];
522:   VecMAXPY(ts->vec_sol,s,w,YdotRHS);
523:   return(0);
524: }

526: static PetscErrorCode TSStep_RK(TS ts)
527: {
528:   TS_RK           *rk   = (TS_RK*)ts->data;
529:   RKTableau        tab  = rk->tableau;
530:   const PetscInt   s = tab->s;
531:   const PetscReal *A = tab->A,*c = tab->c;
532:   PetscScalar     *w = rk->work;
533:   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
534:   PetscBool        FSAL = tab->FSAL;
535:   TSAdapt          adapt;
536:   PetscInt         i,j;
537:   PetscInt         rejections = 0;
538:   PetscBool        stageok,accept = PETSC_TRUE;
539:   PetscReal        next_time_step = ts->time_step;
540:   PetscErrorCode   ierr;

543:   if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE;
544:   if (FSAL) { VecCopy(YdotRHS[s-1],YdotRHS[0]); }

546:   rk->status = TS_STEP_INCOMPLETE;
547:   while (!ts->reason && rk->status != TS_STEP_COMPLETE) {
548:     PetscReal t = ts->ptime;
549:     PetscReal h = ts->time_step;
550:     for (i=0; i<s; i++) {
551:       rk->stage_time = t + h*c[i];
552:       TSPreStage(ts,rk->stage_time);
553:       VecCopy(ts->vec_sol,Y[i]);
554:       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
555:       VecMAXPY(Y[i],i,w,YdotRHS);
556:       TSPostStage(ts,rk->stage_time,i,Y);
557:       TSGetAdapt(ts,&adapt);
558:       TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);
559:       if (!stageok) goto reject_step;
560:       if (FSAL && !i) continue;
561:       TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);
562:     }

564:     rk->status = TS_STEP_INCOMPLETE;
565:     TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);
566:     rk->status = TS_STEP_PENDING;
567:     TSGetAdapt(ts,&adapt);
568:     TSAdaptCandidatesClear(adapt);
569:     TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);
570:     TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);
571:     rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
572:     if (!accept) { /* Roll back the current step */
573:       TSRollBack_RK(ts);
574:       ts->time_step = next_time_step;
575:       goto reject_step;
576:     }

578:     if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation*/
579:       rk->ptime     = ts->ptime;
580:       rk->time_step = ts->time_step;
581:     }

583:     ts->ptime += ts->time_step;
584:     ts->time_step = next_time_step;
585:     break;

587:   reject_step:
588:     ts->reject++; accept = PETSC_FALSE;
589:     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
590:       ts->reason = TS_DIVERGED_STEP_REJECTED;
591:       PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);
592:     }
593:   }
594:   return(0);
595: }

597: static PetscErrorCode TSAdjointSetUp_RK(TS ts)
598: {
599:   TS_RK         *rk  = (TS_RK*)ts->data;
600:   RKTableau      tab = rk->tableau;
601:   PetscInt       s   = tab->s;

605:   if (ts->adjointsetupcalled++) return(0);
606:   VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecDeltaLam);
607:   if(ts->vecs_sensip) {
608:     VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numcost,&rk->VecDeltaMu);
609:   }
610:   return(0);
611: }

613: static PetscErrorCode TSAdjointStep_RK(TS ts)
614: {
615:   TS_RK           *rk   = (TS_RK*)ts->data;
616:   RKTableau        tab  = rk->tableau;
617:   const PetscInt   s    = tab->s;
618:   const PetscReal *A = tab->A,*b = tab->b,*c = tab->c;
619:   PetscScalar     *w    = rk->work;
620:   Vec             *Y    = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu;
621:   PetscInt         i,j,nadj;
622:   PetscReal        t = ts->ptime;
623:   PetscErrorCode   ierr;
624:   PetscReal        h = ts->time_step;

627:   rk->status = TS_STEP_INCOMPLETE;
628:   for (i=s-1; i>=0; i--) {
629:     Mat       J;
630:     PetscReal stage_time = t + h*(1.0-c[i]);
631:     PetscBool zero = PETSC_FALSE;

633:     TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);
634:     TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);
635:     if (ts->vec_costintegral) {
636:       TSAdjointComputeDRDYFunction(ts,stage_time,Y[i],ts->vecs_drdy);
637:     }
638:     /* Stage values of mu */
639:     if (ts->vecs_sensip) {
640:       TSAdjointComputeRHSJacobian(ts,stage_time,Y[i],ts->Jacp);
641:       if (ts->vec_costintegral) {
642:         TSAdjointComputeDRDPFunction(ts,stage_time,Y[i],ts->vecs_drdp);
643:       }
644:     }

646:     if (b[i] == 0 && i == s-1) zero = PETSC_TRUE;
647:     for (nadj=0; nadj<ts->numcost; nadj++) {
648:       DM  dm;
649:       Vec VecSensiTemp;

651:       TSGetDM(ts,&dm);
652:       DMGetGlobalVector(dm,&VecSensiTemp);
653:       /* Stage values of lambda */
654:       if (!zero) {
655:         VecCopy(ts->vecs_sensi[nadj],VecSensiTemp);
656:         VecScale(VecSensiTemp,-h*b[i]);
657:         for (j=i+1; j<s; j++) w[j-i-1]  = -h*A[j*s+i];
658:         VecMAXPY(VecSensiTemp,s-i-1,w,&VecDeltaLam[nadj*s+i+1]);
659:         MatMultTranspose(J,VecSensiTemp,VecDeltaLam[nadj*s+i]);
660:       } else {
661:         VecSet(VecDeltaLam[nadj*s+i],0);
662:       }
663:       if (ts->vec_costintegral) {
664:         VecAXPY(VecDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdy[nadj]);
665:       }

667:       /* Stage values of mu */
668:       if (ts->vecs_sensip) {
669:         if (!zero) {
670:           MatMultTranspose(ts->Jacp,VecSensiTemp,VecDeltaMu[nadj*s+i]);
671:         } else {
672:           VecSet(VecDeltaMu[nadj*s+i],0);
673:         }
674:         if (ts->vec_costintegral) {
675:           VecAXPY(VecDeltaMu[nadj*s+i],-h*b[i],ts->vecs_drdp[nadj]);
676:         }
677:       }
678:       DMRestoreGlobalVector(dm,&VecSensiTemp);
679:     }
680:   }

682:   for (j=0; j<s; j++) w[j] = 1.0;
683:   for (nadj=0; nadj<ts->numcost; nadj++) {
684:     VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);
685:     if (ts->vecs_sensip) {
686:       VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);
687:     }
688:   }
689:   rk->status = TS_STEP_COMPLETE;
690:   return(0);
691: }

693: static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X)
694: {
695:   TS_RK           *rk = (TS_RK*)ts->data;
696:   PetscInt         s  = rk->tableau->s,p = rk->tableau->p,i,j;
697:   PetscReal        h;
698:   PetscReal        tt,t;
699:   PetscScalar     *b;
700:   const PetscReal *B = rk->tableau->binterp;
701:   PetscErrorCode   ierr;

704:   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);

706:   switch (rk->status) {
707:   case TS_STEP_INCOMPLETE:
708:   case TS_STEP_PENDING:
709:     h = ts->time_step;
710:     t = (itime - ts->ptime)/h;
711:     break;
712:   case TS_STEP_COMPLETE:
713:     h = ts->ptime - ts->ptime_prev;
714:     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
715:     break;
716:   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
717:   }
718:   PetscMalloc1(s,&b);
719:   for (i=0; i<s; i++) b[i] = 0;
720:   for (j=0,tt=t; j<p; j++,tt*=t) {
721:     for (i=0; i<s; i++) {
722:       b[i]  += h * B[i*p+j] * tt;
723:     }
724:   }

726:   VecCopy(rk->Y[0],X);
727:   VecMAXPY(X,s,b,rk->YdotRHS);

729:   PetscFree(b);
730:   return(0);
731: }

733: /*------------------------------------------------------------*/

735: static PetscErrorCode TSRKTableauReset(TS ts)
736: {
737:   TS_RK          *rk = (TS_RK*)ts->data;
738:   RKTableau      tab = rk->tableau;

742:   if (!tab) return(0);
743:   PetscFree(rk->work);
744:   VecDestroyVecs(tab->s,&rk->Y);
745:   VecDestroyVecs(tab->s,&rk->YdotRHS);
746:   VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaLam);
747:   VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaMu);
748:   return(0);
749: }

751: static PetscErrorCode TSReset_RK(TS ts)
752: {
753:   TS_RK         *rk = (TS_RK*)ts->data;

757:   TSRKTableauReset(ts);
758:   VecDestroy(&rk->VecCostIntegral0);
759:   return(0);
760: }

762: static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx)
763: {
765:   return(0);
766: }

768: static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
769: {
771:   return(0);
772: }


775: static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx)
776: {
778:   return(0);
779: }

781: static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
782: {

785:   return(0);
786: }
787: /*
788: static PetscErrorCode RKSetAdjCoe(RKTableau tab)
789: {
790:   PetscReal *A,*b;
791:   PetscInt        s,i,j;
792:   PetscErrorCode  ierr;

795:   s    = tab->s;
796:   PetscMalloc2(s*s,&A,s,&b);

798:   for (i=0; i<s; i++)
799:     for (j=0; j<s; j++) {
800:       A[i*s+j] = (tab->b[s-1-i]==0) ? 0: -tab->A[s-1-i+(s-1-j)*s] * tab->b[s-1-j] / tab->b[s-1-i];
801:       PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);
802:     }

804:   for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i];

806:   PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));
807:   PetscMemcpy(tab->b,b,s*sizeof(b[0]));
808:   PetscFree2(A,b);
809:   return(0);
810: }
811: */

813: static PetscErrorCode TSRKTableauSetUp(TS ts)
814: {
815:   TS_RK         *rk  = (TS_RK*)ts->data;
816:   RKTableau      tab = rk->tableau;

820:   PetscMalloc1(tab->s,&rk->work);
821:   VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);
822:   VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);
823:   return(0);
824: }


827: static PetscErrorCode TSSetUp_RK(TS ts)
828: {
829:   TS_RK         *rk = (TS_RK*)ts->data;
831:   DM             dm;

834:   TSCheckImplicitTerm(ts);
835:   TSRKTableauSetUp(ts);
836:   if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */
837:     VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);
838:   }
839:   TSGetDM(ts,&dm);
840:   DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);
841:   DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);
842:   return(0);
843: }


846: /*------------------------------------------------------------*/

848: static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts)
849: {
850:   TS_RK         *rk = (TS_RK*)ts->data;

854:   PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");
855:   {
856:     RKTableauLink  link;
857:     PetscInt       count,choice;
858:     PetscBool      flg;
859:     const char   **namelist;
860:     for (link=RKTableauList,count=0; link; link=link->next,count++) ;
861:     PetscMalloc1(count,(char***)&namelist);
862:     for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
863:     PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);
864:     if (flg) {TSRKSetType(ts,namelist[choice]);}
865:     PetscFree(namelist);
866:   }
867:   PetscOptionsTail();
868:   return(0);
869: }

871: static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
872: {
873:   TS_RK          *rk = (TS_RK*)ts->data;
874:   PetscBool      iascii;

878:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
879:   if (iascii) {
880:     RKTableau tab  = rk->tableau;
881:     TSRKType  rktype;
882:     char      buf[512];
883:     TSRKGetType(ts,&rktype);
884:     PetscViewerASCIIPrintf(viewer,"  RK type %s\n",rktype);
885:     PetscViewerASCIIPrintf(viewer,"  Order: %D\n",tab->order);
886:     PetscViewerASCIIPrintf(viewer,"  FSAL property: %s\n",tab->FSAL ? "yes" : "no");
887:     PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);
888:     PetscViewerASCIIPrintf(viewer,"  Abscissa c = %s\n",buf);
889:   }
890:   return(0);
891: }

893: static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer)
894: {
896:   TSAdapt        adapt;

899:   TSGetAdapt(ts,&adapt);
900:   TSAdaptLoad(adapt,viewer);
901:   return(0);
902: }

904: /*@C
905:   TSRKSetType - Set the type of RK scheme

907:   Logically collective

909:   Input Parameter:
910: +  ts - timestepping context
911: -  rktype - type of RK-scheme

913:   Options Database:
914: .   -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs>

916:   Level: intermediate

918: .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS
919: @*/
920: PetscErrorCode TSRKSetType(TS ts,TSRKType rktype)
921: {

927:   PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));
928:   return(0);
929: }

931: /*@C
932:   TSRKGetType - Get the type of RK scheme

934:   Logically collective

936:   Input Parameter:
937: .  ts - timestepping context

939:   Output Parameter:
940: .  rktype - type of RK-scheme

942:   Level: intermediate

944: .seealso: TSRKGetType()
945: @*/
946: PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype)
947: {

952:   PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));
953:   return(0);
954: }

956: static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype)
957: {
958:   TS_RK *rk = (TS_RK*)ts->data;

961:   *rktype = rk->tableau->name;
962:   return(0);
963: }
964: static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype)
965: {
966:   TS_RK          *rk = (TS_RK*)ts->data;
968:   PetscBool      match;
969:   RKTableauLink  link;

972:   if (rk->tableau) {
973:     PetscStrcmp(rk->tableau->name,rktype,&match);
974:     if (match) return(0);
975:   }
976:   for (link = RKTableauList; link; link=link->next) {
977:     PetscStrcmp(link->tab.name,rktype,&match);
978:     if (match) {
979:       if (ts->setupcalled) {TSRKTableauReset(ts);}
980:       rk->tableau = &link->tab;
981:       if (ts->setupcalled) {TSRKTableauSetUp(ts);}
982:       ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
983:       return(0);
984:     }
985:   }
986:   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype);
987:   return(0);
988: }

990: static PetscErrorCode  TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y)
991: {
992:   TS_RK *rk = (TS_RK*)ts->data;

995:   *ns = rk->tableau->s;
996:   if (Y) *Y = rk->Y;
997:   return(0);
998: }

1000: static PetscErrorCode TSDestroy_RK(TS ts)
1001: {

1005:   TSReset_RK(ts);
1006:   if (ts->dm) {
1007:     DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);
1008:     DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);
1009:   }
1010:   PetscFree(ts->data);
1011:   PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);
1012:   PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);
1013:   return(0);
1014: }

1016: /* ------------------------------------------------------------ */
1017: /*MC
1018:       TSRK - ODE and DAE solver using Runge-Kutta schemes

1020:   The user should provide the right hand side of the equation
1021:   using TSSetRHSFunction().

1023:   Notes:
1024:   The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type

1026:   Level: beginner

1028: .seealso:  TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3,
1029:            TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister()

1031: M*/
1032: PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
1033: {
1034:   TS_RK          *rk;

1038:   TSRKInitializePackage();

1040:   ts->ops->reset          = TSReset_RK;
1041:   ts->ops->destroy        = TSDestroy_RK;
1042:   ts->ops->view           = TSView_RK;
1043:   ts->ops->load           = TSLoad_RK;
1044:   ts->ops->setup          = TSSetUp_RK;
1045:   ts->ops->adjointsetup   = TSAdjointSetUp_RK;
1046:   ts->ops->step           = TSStep_RK;
1047:   ts->ops->interpolate    = TSInterpolate_RK;
1048:   ts->ops->evaluatestep   = TSEvaluateStep_RK;
1049:   ts->ops->rollback       = TSRollBack_RK;
1050:   ts->ops->setfromoptions = TSSetFromOptions_RK;
1051:   ts->ops->getstages      = TSGetStages_RK;
1052:   ts->ops->adjointstep    = TSAdjointStep_RK;

1054:   ts->ops->adjointintegral = TSAdjointCostIntegral_RK;
1055:   ts->ops->forwardintegral = TSForwardCostIntegral_RK;

1057:   PetscNewLog(ts,&rk);
1058:   ts->data = (void*)rk;

1060:   PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);
1061:   PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);

1063:   TSRKSetType(ts,TSRKDefault);
1064:   return(0);
1065: }