Actual source code: rk.c

petsc-3.9.4 2018-09-11
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:   PetscNew(&link);
370:   t = &link->tab;

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

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

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

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

399: /*
400:  The step completion formula is

402:  x1 = x0 + h b^T YdotRHS

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

732: /*------------------------------------------------------------*/

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

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

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

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

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

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


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

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

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

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

797:   for (i=0; i<s; i++)
798:     for (j=0; j<s; j++) {
799:       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];
800:       PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);
801:     }

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

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

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

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


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

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


845: /*------------------------------------------------------------*/

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

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

870: static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
871: {
873:   PetscInt       i;
874:   size_t         left,count;
875:   char           *p;

878:   for (i=0,p=buf,left=len; i<n; i++) {
879:     PetscSNPrintfCount(p,left,fmt,&count,x[i]);
880:     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
881:     left -= count;
882:     p    += count;
883:     *p++  = ' ';
884:   }
885:   p[i ? 0 : -1] = 0;
886:   return(0);
887: }

889: static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
890: {
891:   TS_RK          *rk = (TS_RK*)ts->data;
892:   PetscBool      iascii;

896:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
897:   if (iascii) {
898:     RKTableau tab  = rk->tableau;
899:     TSRKType  rktype;
900:     char      buf[512];
901:     TSRKGetType(ts,&rktype);
902:     PetscViewerASCIIPrintf(viewer,"  RK type %s\n",rktype);
903:     PetscViewerASCIIPrintf(viewer,"  Order: %D\n",tab->order);
904:     PetscViewerASCIIPrintf(viewer,"  FSAL property: %s\n",tab->FSAL ? "yes" : "no");
905:     PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);
906:     PetscViewerASCIIPrintf(viewer,"  Abscissa c = %s\n",buf);
907:   }
908:   return(0);
909: }

911: static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer)
912: {
914:   TSAdapt        adapt;

917:   TSGetAdapt(ts,&adapt);
918:   TSAdaptLoad(adapt,viewer);
919:   return(0);
920: }

922: /*@C
923:   TSRKSetType - Set the type of RK scheme

925:   Logically collective

927:   Input Parameter:
928: +  ts - timestepping context
929: -  rktype - type of RK-scheme

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

934:   Level: intermediate

936: .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS
937: @*/
938: PetscErrorCode TSRKSetType(TS ts,TSRKType rktype)
939: {

945:   PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));
946:   return(0);
947: }

949: /*@C
950:   TSRKGetType - Get the type of RK scheme

952:   Logically collective

954:   Input Parameter:
955: .  ts - timestepping context

957:   Output Parameter:
958: .  rktype - type of RK-scheme

960:   Level: intermediate

962: .seealso: TSRKGetType()
963: @*/
964: PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype)
965: {

970:   PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));
971:   return(0);
972: }

974: static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype)
975: {
976:   TS_RK *rk = (TS_RK*)ts->data;

979:   *rktype = rk->tableau->name;
980:   return(0);
981: }
982: static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype)
983: {
984:   TS_RK          *rk = (TS_RK*)ts->data;
986:   PetscBool      match;
987:   RKTableauLink  link;

990:   if (rk->tableau) {
991:     PetscStrcmp(rk->tableau->name,rktype,&match);
992:     if (match) return(0);
993:   }
994:   for (link = RKTableauList; link; link=link->next) {
995:     PetscStrcmp(link->tab.name,rktype,&match);
996:     if (match) {
997:       if (ts->setupcalled) {TSRKTableauReset(ts);}
998:       rk->tableau = &link->tab;
999:       if (ts->setupcalled) {TSRKTableauSetUp(ts);}
1000:       ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
1001:       return(0);
1002:     }
1003:   }
1004:   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype);
1005:   return(0);
1006: }

1008: static PetscErrorCode  TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y)
1009: {
1010:   TS_RK *rk = (TS_RK*)ts->data;

1013:   *ns = rk->tableau->s;
1014:   if (Y) *Y = rk->Y;
1015:   return(0);
1016: }

1018: static PetscErrorCode TSDestroy_RK(TS ts)
1019: {

1023:   TSReset_RK(ts);
1024:   if (ts->dm) {
1025:     DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);
1026:     DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);
1027:   }
1028:   PetscFree(ts->data);
1029:   PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);
1030:   PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);
1031:   return(0);
1032: }

1034: /* ------------------------------------------------------------ */
1035: /*MC
1036:       TSRK - ODE and DAE solver using Runge-Kutta schemes

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

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

1044:   Level: beginner

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

1049: M*/
1050: PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
1051: {
1052:   TS_RK          *rk;

1056:   TSRKInitializePackage();

1058:   ts->ops->reset          = TSReset_RK;
1059:   ts->ops->destroy        = TSDestroy_RK;
1060:   ts->ops->view           = TSView_RK;
1061:   ts->ops->load           = TSLoad_RK;
1062:   ts->ops->setup          = TSSetUp_RK;
1063:   ts->ops->adjointsetup   = TSAdjointSetUp_RK;
1064:   ts->ops->step           = TSStep_RK;
1065:   ts->ops->interpolate    = TSInterpolate_RK;
1066:   ts->ops->evaluatestep   = TSEvaluateStep_RK;
1067:   ts->ops->rollback       = TSRollBack_RK;
1068:   ts->ops->setfromoptions = TSSetFromOptions_RK;
1069:   ts->ops->getstages      = TSGetStages_RK;
1070:   ts->ops->adjointstep    = TSAdjointStep_RK;

1072:   ts->ops->adjointintegral = TSAdjointCostIntegral_RK;
1073:   ts->ops->forwardintegral = TSForwardCostIntegral_RK;

1075:   PetscNewLog(ts,&rk);
1076:   ts->data = (void*)rk;

1078:   PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);
1079:   PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);

1081:   TSRKSetType(ts,TSRKDefault);
1082:   return(0);
1083: }