xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision 7dbacdf2ee95ecc467f6bb6cac89d21bae40344e)
1 /*
2   Code for time stepping with the Runge-Kutta method
3 
4   Notes:
5   The general system is written as
6 
7   Udot = F(t,U)
8 
9 */
10 
11 #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
12 #include <petscdm.h>
13 #include <../src/ts/impls/explicit/rk/rk.h>
14 #include <../src/ts/impls/explicit/rk/mrk.h>
15 
16 static TSRKType  TSRKDefault = TSRK3BS;
17 static PetscBool TSRKRegisterAllCalled;
18 static PetscBool TSRKPackageInitialized;
19 
20 static RKTableauLink RKTableauList;
21 
22 /*MC
23      TSRK1FE - First order forward Euler scheme.
24 
25      This method has one stage.
26 
27      Options database:
28 .     -ts_rk_type 1fe
29 
30      Level: advanced
31 
32 .seealso: TSRK, TSRKType, TSRKSetType()
33 M*/
34 /*MC
35      TSRK2A - Second order RK scheme.
36 
37      This method has two stages.
38 
39      Options database:
40 .     -ts_rk_type 2a
41 
42      Level: advanced
43 
44 .seealso: TSRK, TSRKType, TSRKSetType()
45 M*/
46 /*MC
47      TSRK3 - Third order RK scheme.
48 
49      This method has three stages.
50 
51      Options database:
52 .     -ts_rk_type 3
53 
54      Level: advanced
55 
56 .seealso: TSRK, TSRKType, TSRKSetType()
57 M*/
58 /*MC
59      TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method.
60 
61      This method has four stages with the First Same As Last (FSAL) property.
62 
63      Options database:
64 .     -ts_rk_type 3bs
65 
66      Level: advanced
67 
68      References: https://doi.org/10.1016/0893-9659(89)90079-7
69 
70 .seealso: TSRK, TSRKType, TSRKSetType()
71 M*/
72 /*MC
73      TSRK4 - Fourth order RK scheme.
74 
75      This is the classical Runge-Kutta method with four stages.
76 
77      Options database:
78 .     -ts_rk_type 4
79 
80      Level: advanced
81 
82 .seealso: TSRK, TSRKType, TSRKSetType()
83 M*/
84 /*MC
85      TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method.
86 
87      This method has six stages.
88 
89      Options database:
90 .     -ts_rk_type 5f
91 
92      Level: advanced
93 
94 .seealso: TSRK, TSRKType, TSRKSetType()
95 M*/
96 /*MC
97      TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method.
98 
99      This method has seven stages with the First Same As Last (FSAL) property.
100 
101      Options database:
102 .     -ts_rk_type 5dp
103 
104      Level: advanced
105 
106      References: https://doi.org/10.1016/0771-050X(80)90013-3
107 
108 .seealso: TSRK, TSRKType, TSRKSetType()
109 M*/
110 /*MC
111      TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method.
112 
113      This method has eight stages with the First Same As Last (FSAL) property.
114 
115      Options database:
116 .     -ts_rk_type 5bs
117 
118      Level: advanced
119 
120      References: https://doi.org/10.1016/0898-1221(96)00141-1
121 
122 .seealso: TSRK, TSRKType, TSRKSetType()
123 M*/
124 
125 /*@C
126   TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK
127 
128   Not Collective, but should be called by all processes which will need the schemes to be registered
129 
130   Level: advanced
131 
132 .keywords: TS, TSRK, register, all
133 
134 .seealso:  TSRKRegisterDestroy()
135 @*/
136 PetscErrorCode TSRKRegisterAll(void)
137 {
138   PetscErrorCode ierr;
139 
140   PetscFunctionBegin;
141   if (TSRKRegisterAllCalled) PetscFunctionReturn(0);
142   TSRKRegisterAllCalled = PETSC_TRUE;
143 
144 #define RC PetscRealConstant
145   {
146     const PetscReal
147       A[1][1] = {{0}},
148       b[1]    = {RC(1.0)};
149     ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
150   }
151   {
152     const PetscReal
153       A[2][2]   = {{0,0},
154                    {RC(1.0),0}},
155       b[2]      =  {RC(0.5),RC(0.5)},
156       bembed[2] =  {RC(1.0),0};
157     ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
158   }
159   {
160     const PetscReal
161       A[3][3] = {{0,0,0},
162                  {RC(2.0)/RC(3.0),0,0},
163                  {RC(-1.0)/RC(3.0),RC(1.0),0}},
164       b[3]    =  {RC(0.25),RC(0.5),RC(0.25)};
165     ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
166   }
167   {
168     const PetscReal
169       A[4][4]   = {{0,0,0,0},
170                    {RC(1.0)/RC(2.0),0,0,0},
171                    {0,RC(3.0)/RC(4.0),0,0},
172                    {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}},
173       b[4]      =  {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0},
174       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)};
175     ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
176   }
177   {
178     const PetscReal
179       A[4][4] = {{0,0,0,0},
180                  {RC(0.5),0,0,0},
181                  {0,RC(0.5),0,0},
182                  {0,0,RC(1.0),0}},
183       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)};
184     ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
185   }
186   {
187     const PetscReal
188       A[6][6]   = {{0,0,0,0,0,0},
189                    {RC(0.25),0,0,0,0,0},
190                    {RC(3.0)/RC(32.0),RC(9.0)/RC(32.0),0,0,0,0},
191                    {RC(1932.0)/RC(2197.0),RC(-7200.0)/RC(2197.0),RC(7296.0)/RC(2197.0),0,0,0},
192                    {RC(439.0)/RC(216.0),RC(-8.0),RC(3680.0)/RC(513.0),RC(-845.0)/RC(4104.0),0,0},
193                    {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}},
194       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)},
195       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};
196     ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
197   }
198   {
199     const PetscReal
200       A[7][7]   = {{0,0,0,0,0,0,0},
201                    {RC(1.0)/RC(5.0),0,0,0,0,0,0},
202                    {RC(3.0)/RC(40.0),RC(9.0)/RC(40.0),0,0,0,0,0},
203                    {RC(44.0)/RC(45.0),RC(-56.0)/RC(15.0),RC(32.0)/RC(9.0),0,0,0,0},
204                    {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},
205                    {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},
206                    {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}},
207       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},
208         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)},
209         binterp[7][5] =  {{RC(1.0),RC(-4034104133.0)/RC(1410260304.0),RC(105330401.0)/RC(33982176.0),RC(-13107642775.0)/RC(11282082432.0),RC(6542295.0)/RC(470086768.0)},
210                     {0,0,0,0,0},
211                     {0,RC(132343189600.0)/RC(32700410799.0),RC(-833316000.0)/RC(131326951.0),RC(91412856700.0)/RC(32700410799.0),RC(-523383600.0)/RC(10900136933.0)},
212                     {0,RC(-115792950.0)/RC(29380423.0),RC(185270875.0)/RC(16991088.0),RC(-12653452475.0)/RC(1880347072.0),RC(98134425.0)/RC(235043384.0)},
213                     {0,RC(70805911779.0)/RC(24914598704.0),RC(-4531260609.0)/RC(600351776.0),RC(988140236175.0)/RC(199316789632.0),RC(-14307999165.0)/RC(24914598704.0)},
214                     {0,RC(-331320693.0)/RC(205662961.0),RC(31361737.0)/RC(7433601.0),RC(-2426908385.0)/RC(822651844.0),RC(97305120.0)/RC(205662961.0)},
215                     {0,RC(44764047.0)/RC(29380423.0),RC(-1532549.0)/RC(353981.0),RC(90730570.0)/RC(29380423.0),RC(-8293050.0)/RC(29380423.0)}};
216 
217         ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,5,binterp[0]);CHKERRQ(ierr);
218   }
219   {
220     const PetscReal
221       A[8][8]   = {{0,0,0,0,0,0,0,0},
222                    {RC(1.0)/RC(6.0),0,0,0,0,0,0,0},
223                    {RC(2.0)/RC(27.0),RC(4.0)/RC(27.0),0,0,0,0,0,0},
224                    {RC(183.0)/RC(1372.0),RC(-162.0)/RC(343.0),RC(1053.0)/RC(1372.0),0,0,0,0,0},
225                    {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},
226                    {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},
227                    {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},
228                    {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}},
229       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},
230       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)};
231     ierr = TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
232   }
233 #undef RC
234   PetscFunctionReturn(0);
235 }
236 
237 /*@C
238    TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister().
239 
240    Not Collective
241 
242    Level: advanced
243 
244 .keywords: TSRK, register, destroy
245 .seealso: TSRKRegister(), TSRKRegisterAll()
246 @*/
247 PetscErrorCode TSRKRegisterDestroy(void)
248 {
249   PetscErrorCode ierr;
250   RKTableauLink  link;
251 
252   PetscFunctionBegin;
253   while ((link = RKTableauList)) {
254     RKTableau t = &link->tab;
255     RKTableauList = link->next;
256     ierr = PetscFree3(t->A,t->b,t->c);  CHKERRQ(ierr);
257     ierr = PetscFree (t->bembed);       CHKERRQ(ierr);
258     ierr = PetscFree (t->binterp);      CHKERRQ(ierr);
259     ierr = PetscFree (t->name);         CHKERRQ(ierr);
260     ierr = PetscFree (link);            CHKERRQ(ierr);
261   }
262   TSRKRegisterAllCalled = PETSC_FALSE;
263   PetscFunctionReturn(0);
264 }
265 
266 /*@C
267   TSRKInitializePackage - This function initializes everything in the TSRK package. It is called
268   from TSInitializePackage().
269 
270   Level: developer
271 
272 .keywords: TS, TSRK, initialize, package
273 .seealso: PetscInitialize()
274 @*/
275 PetscErrorCode TSRKInitializePackage(void)
276 {
277   PetscErrorCode ierr;
278 
279   PetscFunctionBegin;
280   if (TSRKPackageInitialized) PetscFunctionReturn(0);
281   TSRKPackageInitialized = PETSC_TRUE;
282   ierr = TSRKRegisterAll();CHKERRQ(ierr);
283   ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr);
284   PetscFunctionReturn(0);
285 }
286 
287 /*@C
288   TSRKFinalizePackage - This function destroys everything in the TSRK package. It is
289   called from PetscFinalize().
290 
291   Level: developer
292 
293 .keywords: Petsc, destroy, package
294 .seealso: PetscFinalize()
295 @*/
296 PetscErrorCode TSRKFinalizePackage(void)
297 {
298   PetscErrorCode ierr;
299 
300   PetscFunctionBegin;
301   TSRKPackageInitialized = PETSC_FALSE;
302   ierr = TSRKRegisterDestroy();CHKERRQ(ierr);
303   PetscFunctionReturn(0);
304 }
305 
306 /*@C
307    TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
308 
309    Not Collective, but the same schemes should be registered on all processes on which they will be used
310 
311    Input Parameters:
312 +  name - identifier for method
313 .  order - approximation order of method
314 .  s - number of stages, this is the dimension of the matrices below
315 .  A - stage coefficients (dimension s*s, row-major)
316 .  b - step completion table (dimension s; NULL to use last row of A)
317 .  c - abscissa (dimension s; NULL to use row sums of A)
318 .  bembed - completion table for embedded method (dimension s; NULL if not available)
319 .  p - Order of the interpolation scheme, equal to the number of columns of binterp
320 -  binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1)
321 
322    Notes:
323    Several RK methods are provided, this function is only needed to create new methods.
324 
325    Level: advanced
326 
327 .keywords: TS, register
328 
329 .seealso: TSRK
330 @*/
331 PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s,
332                             const PetscReal A[],const PetscReal b[],const PetscReal c[],
333                             const PetscReal bembed[],PetscInt p,const PetscReal binterp[])
334 {
335   PetscErrorCode  ierr;
336   RKTableauLink   link;
337   RKTableau       t;
338   PetscInt        i,j;
339 
340   PetscFunctionBegin;
341   PetscValidCharPointer(name,1);
342   PetscValidRealPointer(A,4);
343   if (b) PetscValidRealPointer(b,5);
344   if (c) PetscValidRealPointer(c,6);
345   if (bembed) PetscValidRealPointer(bembed,7);
346   if (binterp || p > 1) PetscValidRealPointer(binterp,9);
347 
348   ierr = TSRKInitializePackage();CHKERRQ(ierr);
349   ierr = PetscNew(&link);CHKERRQ(ierr);
350   t = &link->tab;
351 
352   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
353   t->order = order;
354   t->s = s;
355   ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr);
356   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
357   if (b)  { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); }
358   else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i];
359   if (c)  { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); }
360   else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j];
361   t->FSAL = PETSC_TRUE;
362   for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE;
363 
364   if (bembed) {
365     ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr);
366     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
367   }
368 
369   if (!binterp) { p = 1; binterp = t->b; }
370   t->p = p;
371   ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr);
372   ierr = PetscMemcpy(t->binterp,binterp,s*p*sizeof(binterp[0]));CHKERRQ(ierr);
373 
374   link->next = RKTableauList;
375   RKTableauList = link;
376   PetscFunctionReturn(0);
377 }
378 
379 /*
380  This is for single-step RK method
381  The step completion formula is
382 
383  x1 = x0 + h b^T YdotRHS
384 
385  This function can be called before or after ts->vec_sol has been updated.
386  Suppose we have a completion formula (b) and an embedded formula (be) of different order.
387  We can write
388 
389  x1e = x0 + h be^T YdotRHS
390      = x1 - h b^T YdotRHS + h be^T YdotRHS
391      = x1 + h (be - b)^T YdotRHS
392 
393  so we can evaluate the method with different order even after the step has been optimistically completed.
394 */
395 static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done)
396 {
397   TS_RK          *rk   = (TS_RK*)ts->data;
398   RKTableau      tab  = rk->tableau;
399   PetscScalar    *w    = rk->work;
400   PetscReal      h;
401   PetscInt       s    = tab->s,j;
402   PetscErrorCode ierr;
403 
404   PetscFunctionBegin;
405   switch (rk->status) {
406   case TS_STEP_INCOMPLETE:
407   case TS_STEP_PENDING:
408     h = ts->time_step; break;
409   case TS_STEP_COMPLETE:
410     h = ts->ptime - ts->ptime_prev; break;
411   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
412   }
413   if (order == tab->order) {
414     if (rk->status == TS_STEP_INCOMPLETE) {
415       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
416       for (j=0; j<s; j++) w[j] = h*tab->b[j]/rk->dtratio;
417       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
418     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
419     PetscFunctionReturn(0);
420   } else if (order == tab->order-1) {
421     if (!tab->bembed) goto unavailable;
422     if (rk->status == TS_STEP_INCOMPLETE) { /*Complete with the embedded method (be)*/
423       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
424       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
425       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
426     } else {  /*Rollback and re-complete using (be-b) */
427       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
428       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
429       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
430     }
431     if (done) *done = PETSC_TRUE;
432     PetscFunctionReturn(0);
433   }
434 unavailable:
435   if (done) *done = PETSC_FALSE;
436   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);
437   PetscFunctionReturn(0);
438 }
439 
440 static PetscErrorCode TSForwardCostIntegral_RK(TS ts)
441 {
442   TS_RK           *rk = (TS_RK*)ts->data;
443   RKTableau       tab = rk->tableau;
444   const PetscInt  s = tab->s;
445   const PetscReal *b = tab->b,*c = tab->c;
446   Vec             *Y = rk->Y;
447   PetscInt        i;
448   PetscErrorCode  ierr;
449 
450   PetscFunctionBegin;
451   /* backup cost integral */
452   ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr);
453   for (i=s-1; i>=0; i--) {
454     /* Evolve ts->vec_costintegral to compute integrals */
455     ierr = TSComputeCostIntegrand(ts,rk->ptime+rk->time_step*c[i],Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
456     ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
457   }
458   PetscFunctionReturn(0);
459 }
460 
461 static PetscErrorCode TSAdjointCostIntegral_RK(TS ts)
462 {
463   TS_RK           *rk = (TS_RK*)ts->data;
464   RKTableau       tab = rk->tableau;
465   const PetscInt  s = tab->s;
466   const PetscReal *b = tab->b,*c = tab->c;
467   Vec             *Y = rk->Y;
468   PetscInt        i;
469   PetscErrorCode  ierr;
470 
471   PetscFunctionBegin;
472   for (i=s-1; i>=0; i--) {
473     /* Evolve ts->vec_costintegral to compute integrals */
474     ierr = TSComputeCostIntegrand(ts,ts->ptime+ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
475     ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
476   }
477   PetscFunctionReturn(0);
478 }
479 
480 static PetscErrorCode TSRollBack_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;
486   PetscScalar     *w = rk->work;
487   Vec             *YdotRHS = rk->YdotRHS;
488   PetscInt        j;
489   PetscReal       h;
490   PetscErrorCode  ierr;
491 
492   PetscFunctionBegin;
493   switch (rk->status) {
494   case TS_STEP_INCOMPLETE:
495   case TS_STEP_PENDING:
496     h = ts->time_step; break;
497   case TS_STEP_COMPLETE:
498     h = ts->ptime - ts->ptime_prev; break;
499   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
500   }
501   for (j=0; j<s; j++) w[j] = -h*b[j];
502   ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr);
503   PetscFunctionReturn(0);
504 }
505 
506 static PetscErrorCode TSStep_RK(TS ts)
507 {
508   TS_RK           *rk  = (TS_RK*)ts->data;
509   RKTableau       tab  = rk->tableau;
510   const PetscInt  s = tab->s;
511   const PetscReal *A = tab->A,*c = tab->c;
512   PetscScalar     *w = rk->work;
513   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
514   PetscBool       FSAL = tab->FSAL;
515   TSAdapt         adapt;
516   PetscInt        i,j;
517   PetscInt        rejections = 0;
518   PetscBool       stageok,accept = PETSC_TRUE;
519   PetscReal       next_time_step = ts->time_step;
520   PetscErrorCode  ierr;
521 
522   PetscFunctionBegin;
523   if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE;
524   if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); }
525 
526   rk->status = TS_STEP_INCOMPLETE;
527   while (!ts->reason && rk->status != TS_STEP_COMPLETE) {
528     PetscReal t = ts->ptime;
529     PetscReal h = ts->time_step;
530     for (i=0; i<s; i++) {
531       rk->stage_time = t + h*c[i];
532       ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr);
533       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
534       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
535       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
536       ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
537       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
538       ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr);
539       if (!stageok) goto reject_step;
540       if (FSAL && !i) continue;
541       ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
542     }
543 
544     rk->status = TS_STEP_INCOMPLETE;
545     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
546     rk->status = TS_STEP_PENDING;
547     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
548     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
549     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr);
550     ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
551     rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
552     if (!accept) { /* Roll back the current step */
553       ierr = TSRollBack_RK(ts);CHKERRQ(ierr);
554       ts->time_step = next_time_step;
555       goto reject_step;
556     }
557 
558     if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
559       rk->ptime     = ts->ptime;
560       rk->time_step = ts->time_step;
561     }
562 
563     ts->ptime += ts->time_step;
564     ts->time_step = next_time_step;
565     break;
566 
567     reject_step:
568     ts->reject++; accept = PETSC_FALSE;
569     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
570       ts->reason = TS_DIVERGED_STEP_REJECTED;
571       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
572     }
573   }
574   PetscFunctionReturn(0);
575 }
576 
577 static PetscErrorCode TSAdjointSetUp_RK(TS ts)
578 {
579   TS_RK          *rk  = (TS_RK*)ts->data;
580   RKTableau      tab = rk->tableau;
581   PetscInt       s   = tab->s;
582   PetscErrorCode ierr;
583 
584   PetscFunctionBegin;
585   if (ts->adjointsetupcalled++) PetscFunctionReturn(0);
586   ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr);
587   if(ts->vecs_sensip) {
588     ierr = VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr);
589   }
590   PetscFunctionReturn(0);
591 }
592 
593 static PetscErrorCode TSAdjointStep_RK(TS ts)
594 {
595   TS_RK            *rk  = (TS_RK*)ts->data;
596   RKTableau        tab  = rk->tableau;
597   const PetscInt   s    = tab->s;
598   const PetscReal  *A = tab->A,*b = tab->b,*c = tab->c;
599   PetscScalar      *w    = rk->work;
600   Vec              *Y    = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu;
601   PetscInt         i,j,nadj;
602   PetscReal        t = ts->ptime;
603   PetscErrorCode   ierr;
604   PetscReal        h = ts->time_step;
605 
606   PetscFunctionBegin;
607   rk->status = TS_STEP_INCOMPLETE;
608   for (i=s-1; i>=0; i--) {
609     Mat       J;
610     PetscReal stage_time = t + h*(1.0-c[i]);
611     PetscBool zero = PETSC_FALSE;
612 
613     ierr = TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr);
614     ierr = TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);CHKERRQ(ierr);
615     if (ts->vec_costintegral) {
616       ierr = TSAdjointComputeDRDYFunction(ts,stage_time,Y[i],ts->vecs_drdy);CHKERRQ(ierr);
617     }
618     /* Stage values of mu */
619     if (ts->vecs_sensip) {
620       ierr = TSAdjointComputeRHSJacobian(ts,stage_time,Y[i],ts->Jacp);CHKERRQ(ierr);
621       if (ts->vec_costintegral) {
622         ierr = TSAdjointComputeDRDPFunction(ts,stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr);
623       }
624     }
625 
626     if (b[i] == 0 && i == s-1) zero = PETSC_TRUE;
627     for (nadj=0; nadj<ts->numcost; nadj++) {
628       DM  dm;
629       Vec VecSensiTemp;
630 
631       ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
632       ierr = DMGetGlobalVector(dm,&VecSensiTemp);CHKERRQ(ierr);
633       /* Stage values of lambda */
634       if (!zero) {
635         ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp);CHKERRQ(ierr);
636         ierr = VecScale(VecSensiTemp,-h*b[i]);CHKERRQ(ierr);
637         for (j=i+1; j<s; j++) w[j-i-1]  = -h*A[j*s+i];
638         ierr = VecMAXPY(VecSensiTemp,s-i-1,w,&VecDeltaLam[nadj*s+i+1]);CHKERRQ(ierr);
639         ierr = MatMultTranspose(J,VecSensiTemp,VecDeltaLam[nadj*s+i]);CHKERRQ(ierr);
640       } else {
641         ierr = VecSet(VecDeltaLam[nadj*s+i],0);CHKERRQ(ierr);
642       }
643       if (ts->vec_costintegral) {
644         ierr = VecAXPY(VecDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdy[nadj]);CHKERRQ(ierr);
645       }
646 
647       /* Stage values of mu */
648       if (ts->vecs_sensip) {
649         if (!zero) {
650           ierr = MatMultTranspose(ts->Jacp,VecSensiTemp,VecDeltaMu[nadj*s+i]);CHKERRQ(ierr);
651         } else {
652           ierr = VecSet(VecDeltaMu[nadj*s+i],0);CHKERRQ(ierr);
653         }
654         if (ts->vec_costintegral) {
655           ierr = VecAXPY(VecDeltaMu[nadj*s+i],-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr);
656         }
657       }
658       ierr = DMRestoreGlobalVector(dm,&VecSensiTemp);CHKERRQ(ierr);
659     }
660   }
661 
662   for (j=0; j<s; j++) w[j] = 1.0;
663   for (nadj=0; nadj<ts->numcost; nadj++) {
664     ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr);
665     if (ts->vecs_sensip) {
666       ierr = VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);CHKERRQ(ierr);
667     }
668   }
669   rk->status = TS_STEP_COMPLETE;
670   PetscFunctionReturn(0);
671 }
672 
673 static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X)
674 {
675   TS_RK            *rk = (TS_RK*)ts->data;
676   PetscInt         s  = rk->tableau->s,p = rk->tableau->p,i,j;
677   PetscReal        h;
678   PetscReal        tt,t;
679   PetscScalar      *b;
680   const PetscReal  *B = rk->tableau->binterp;
681   PetscErrorCode   ierr;
682 
683   PetscFunctionBegin;
684   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
685 
686   switch (rk->status) {
687     case TS_STEP_INCOMPLETE:
688     case TS_STEP_PENDING:
689       h = ts->time_step;
690       t = (itime - ts->ptime)/h;
691       break;
692     case TS_STEP_COMPLETE:
693       h = ts->ptime - ts->ptime_prev;
694       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
695       break;
696     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
697   }
698   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
699   for (i=0; i<s; i++) b[i] = 0;
700   for (j=0,tt=t; j<p; j++,tt*=t) {
701     for (i=0; i<s; i++) {
702       b[i]  += h * B[i*p+j] * tt;
703     }
704   }
705   ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr);
706   ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr);
707   ierr = PetscFree(b);CHKERRQ(ierr);
708   PetscFunctionReturn(0);
709 }
710 
711 /*------------------------------------------------------------*/
712 
713 static PetscErrorCode TSRKTableauReset(TS ts)
714 {
715   TS_RK          *rk = (TS_RK*)ts->data;
716   RKTableau      tab = rk->tableau;
717   PetscErrorCode ierr;
718 
719   PetscFunctionBegin;
720   if (!tab) PetscFunctionReturn(0);
721   ierr = PetscFree(rk->work);CHKERRQ(ierr);
722   ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr);
723   ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr);
724   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr);
725   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr);
726   PetscFunctionReturn(0);
727 }
728 
729 static PetscErrorCode TSReset_RK(TS ts)
730 {
731   TS_RK         *rk = (TS_RK*)ts->data;
732   PetscErrorCode ierr;
733 
734   PetscFunctionBegin;
735   ierr = TSRKTableauReset(ts);CHKERRQ(ierr);
736   ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr);
737   if (ts->use_splitrhsfunction) {
738     ierr = PetscTryMethod(ts,"TSReset_RK_MultirateSplit_C",(TS),(ts));CHKERRQ(ierr);
739   } else {
740     ierr = PetscTryMethod(ts,"TSReset_RK_MultirateNonsplit_C",(TS),(ts));CHKERRQ(ierr);
741   }
742   PetscFunctionReturn(0);
743 }
744 
745 static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx)
746 {
747   PetscFunctionBegin;
748   PetscFunctionReturn(0);
749 }
750 
751 static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
752 {
753   PetscFunctionBegin;
754   PetscFunctionReturn(0);
755 }
756 
757 
758 static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx)
759 {
760   PetscFunctionBegin;
761   PetscFunctionReturn(0);
762 }
763 
764 static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
765 {
766 
767   PetscFunctionBegin;
768   PetscFunctionReturn(0);
769 }
770 /*
771 static PetscErrorCode RKSetAdjCoe(RKTableau tab)
772 {
773   PetscReal *A,*b;
774   PetscInt        s,i,j;
775   PetscErrorCode  ierr;
776 
777   PetscFunctionBegin;
778   s    = tab->s;
779   ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr);
780 
781   for (i=0; i<s; i++)
782     for (j=0; j<s; j++) {
783       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];
784       ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr);
785     }
786 
787   for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i];
788 
789   ierr  = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
790   ierr  = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
791   ierr  = PetscFree2(A,b);CHKERRQ(ierr);
792   PetscFunctionReturn(0);
793 }
794 */
795 
796 static PetscErrorCode TSRKTableauSetUp(TS ts)
797 {
798   TS_RK          *rk  = (TS_RK*)ts->data;
799   RKTableau      tab = rk->tableau;
800   PetscErrorCode ierr;
801 
802   PetscFunctionBegin;
803   ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr);
804   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr);
805   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr);
806   PetscFunctionReturn(0);
807 }
808 
809 static PetscErrorCode TSSetUp_RK(TS ts)
810 {
811   TS_RK          *rk = (TS_RK*)ts->data;
812   PetscErrorCode ierr;
813   DM             dm;
814 
815   PetscFunctionBegin;
816   ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr);
817   ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);
818   if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */
819     ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr);
820   }
821   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
822   ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
823   ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
824   if (ts->use_splitrhsfunction) {
825     ierr = PetscTryMethod(ts,"TSSetUp_RK_MultirateSplit_C",(TS),(ts));CHKERRQ(ierr);
826   } else {
827     ierr = PetscTryMethod(ts,"TSSetUp_RK_MultirateNonsplit_C",(TS),(ts));CHKERRQ(ierr);
828   }
829   PetscFunctionReturn(0);
830 }
831 
832 static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts)
833 {
834   TS_RK          *rk = (TS_RK*)ts->data;
835   PetscErrorCode ierr;
836 
837   PetscFunctionBegin;
838   ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr);
839   {
840     RKTableauLink link;
841     PetscInt      count,choice;
842     PetscBool     flg,use_multirate = PETSC_FALSE;
843     const char    **namelist;
844 
845     for (link=RKTableauList,count=0; link; link=link->next,count++) ;
846     ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr);
847     for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
848     ierr = PetscOptionsBool("-ts_rk_multirate","Use interpolation-based multirate RK method","TSRKSetMultirate",rk->use_multirate,&use_multirate,&flg);CHKERRQ(ierr);
849     if (flg) {
850       ierr = TSRKSetMultirate(ts,use_multirate);CHKERRQ(ierr);
851     }
852     ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr);
853     if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);}
854     ierr = PetscFree(namelist);CHKERRQ(ierr);
855   }
856   ierr = PetscOptionsTail();CHKERRQ(ierr);
857   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Multirate methods options","");CHKERRQ(ierr);
858   ierr = PetscOptionsInt("-ts_rk_dtratio","time step ratio between slow and fast","",rk->dtratio,&rk->dtratio,NULL);CHKERRQ(ierr);
859   ierr = PetscOptionsEnd();CHKERRQ(ierr);
860   PetscFunctionReturn(0);
861 }
862 
863 static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
864 {
865   TS_RK          *rk = (TS_RK*)ts->data;
866   PetscBool      iascii;
867   PetscErrorCode ierr;
868 
869   PetscFunctionBegin;
870   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
871   if (iascii) {
872     RKTableau tab  = rk->tableau;
873     TSRKType  rktype;
874     char      buf[512];
875     ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr);
876     ierr = PetscViewerASCIIPrintf(viewer,"  RK type %s\n",rktype);CHKERRQ(ierr);
877     ierr = PetscViewerASCIIPrintf(viewer,"  Order: %D\n",tab->order);CHKERRQ(ierr);
878     ierr = PetscViewerASCIIPrintf(viewer,"  FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr);
879     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
880     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa c = %s\n",buf);CHKERRQ(ierr);
881   }
882   PetscFunctionReturn(0);
883 }
884 
885 static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer)
886 {
887   PetscErrorCode ierr;
888   TSAdapt        adapt;
889 
890   PetscFunctionBegin;
891   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
892   ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr);
893   PetscFunctionReturn(0);
894 }
895 
896 /*@C
897   TSRKSetType - Set the type of RK scheme
898 
899   Logically collective
900 
901   Input Parameter:
902 +  ts - timestepping context
903 -  rktype - type of RK-scheme
904 
905   Options Database:
906 .   -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs>
907 
908   Level: intermediate
909 
910 .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS
911 @*/
912 PetscErrorCode TSRKSetType(TS ts,TSRKType rktype)
913 {
914   PetscErrorCode ierr;
915 
916   PetscFunctionBegin;
917   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
918   PetscValidCharPointer(rktype,2);
919   ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr);
920   PetscFunctionReturn(0);
921 }
922 
923 /*@C
924   TSRKGetType - Get the type of RK scheme
925 
926   Logically collective
927 
928   Input Parameter:
929 .  ts - timestepping context
930 
931   Output Parameter:
932 .  rktype - type of RK-scheme
933 
934   Level: intermediate
935 
936 .seealso: TSRKGetType()
937 @*/
938 PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype)
939 {
940   PetscErrorCode ierr;
941 
942   PetscFunctionBegin;
943   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
944   ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr);
945   PetscFunctionReturn(0);
946 }
947 
948 static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype)
949 {
950   TS_RK *rk = (TS_RK*)ts->data;
951 
952   PetscFunctionBegin;
953   *rktype = rk->tableau->name;
954   PetscFunctionReturn(0);
955 }
956 
957 static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype)
958 {
959   TS_RK          *rk = (TS_RK*)ts->data;
960   PetscErrorCode ierr;
961   PetscBool      match;
962   RKTableauLink  link;
963 
964   PetscFunctionBegin;
965   if (rk->tableau) {
966     ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr);
967     if (match) PetscFunctionReturn(0);
968   }
969   for (link = RKTableauList; link; link=link->next) {
970     ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr);
971     if (match) {
972       if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);}
973       rk->tableau = &link->tab;
974       if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);}
975       ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
976       PetscFunctionReturn(0);
977     }
978   }
979   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype);
980   PetscFunctionReturn(0);
981 }
982 
983 static PetscErrorCode  TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y)
984 {
985   TS_RK *rk = (TS_RK*)ts->data;
986 
987   PetscFunctionBegin;
988   if (ns) *ns = rk->tableau->s;
989   if (Y)   *Y = rk->Y;
990   PetscFunctionReturn(0);
991 }
992 
993 static PetscErrorCode TSDestroy_RK(TS ts)
994 {
995   PetscErrorCode ierr;
996 
997   PetscFunctionBegin;
998   ierr = TSReset_RK(ts);CHKERRQ(ierr);
999   if (ts->dm) {
1000     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
1001     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
1002   }
1003   ierr = PetscFree(ts->data);CHKERRQ(ierr);
1004   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr);
1005   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr);
1006   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",NULL);CHKERRQ(ierr);
1007   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",NULL);CHKERRQ(ierr);
1008   PetscFunctionReturn(0);
1009 }
1010 
1011 /*@C
1012   TSRKSetMultirate - Use the interpolation-based multirate RK method
1013 
1014   Logically collective
1015 
1016   Input Parameter:
1017 +  ts - timestepping context
1018 -  use_multirate - PETSC_TRUE enables the multirate RK method, sets the basic method to be RK2A and sets the ratio between slow stepsize and fast stepsize to be 2
1019 
1020   Options Database:
1021 .   -ts_rk_multirate - <true,false>
1022 
1023   Notes:
1024   The multirate method requires interpolation. The default interpolation works for 1st- and 2nd- order RK, but not for high-order RKs except TSRK5DP which comes with the interpolation coeffcients (binterp).
1025 
1026   Level: intermediate
1027 
1028 .seealso: TSRKGetMultirate()
1029 @*/
1030 PetscErrorCode TSRKSetMultirate(TS ts,PetscBool use_multirate)
1031 {
1032   PetscErrorCode ierr;
1033 
1034   PetscFunctionBegin
1035   ierr = PetscTryMethod(ts,"TSRKSetMultirate_C",(TS,PetscBool),(ts,use_multirate));CHKERRQ(ierr);
1036   PetscFunctionReturn(0);
1037 }
1038 
1039 /*@C
1040   TSRKGetMultirate - Gets whether to Use the interpolation-based multirate RK method
1041 
1042   Not collective
1043 
1044   Input Parameter:
1045 .  ts - timestepping context
1046 
1047   Output Parameter:
1048 .  use_multirate - PETSC_TRUE if the multirate RK method is enabled, PETSC_FALSE otherwise
1049 
1050   Level: intermediate
1051 
1052 .seealso: TSRKSetMultirate()
1053 @*/
1054 PetscErrorCode TSRKGetMultirate(TS ts,PetscBool *use_multirate)
1055 {
1056   PetscErrorCode ierr;
1057 
1058   PetscFunctionBegin;
1059   ierr = PetscUseMethod(ts,"TSRKGetMultirate_C",(TS,PetscBool*),(ts,use_multirate));CHKERRQ(ierr);
1060   PetscFunctionReturn(0);
1061 }
1062 
1063 /*MC
1064       TSRK - ODE and DAE solver using Runge-Kutta schemes
1065 
1066   The user should provide the right hand side of the equation
1067   using TSSetRHSFunction().
1068 
1069   Notes:
1070   The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type
1071 
1072   Level: beginner
1073 
1074 .seealso:  TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3,
1075            TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister(), TSRKSetMultirate(), TSRKGetMultirate()
1076 
1077 M*/
1078 PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
1079 {
1080   TS_RK          *rk;
1081   PetscErrorCode ierr;
1082 
1083   PetscFunctionBegin;
1084   ierr = TSRKInitializePackage();CHKERRQ(ierr);
1085 
1086   ts->ops->reset          = TSReset_RK;
1087   ts->ops->destroy        = TSDestroy_RK;
1088   ts->ops->view           = TSView_RK;
1089   ts->ops->load           = TSLoad_RK;
1090   ts->ops->setup          = TSSetUp_RK;
1091   ts->ops->adjointsetup   = TSAdjointSetUp_RK;
1092   ts->ops->interpolate    = TSInterpolate_RK;
1093   ts->ops->step           = TSStep_RK;
1094   ts->ops->evaluatestep   = TSEvaluateStep_RK;
1095   ts->ops->rollback       = TSRollBack_RK;
1096   ts->ops->setfromoptions = TSSetFromOptions_RK;
1097   ts->ops->getstages      = TSGetStages_RK;
1098   ts->ops->adjointstep    = TSAdjointStep_RK;
1099 
1100   ts->ops->adjointintegral = TSAdjointCostIntegral_RK;
1101   ts->ops->forwardintegral = TSForwardCostIntegral_RK;
1102 
1103   ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr);
1104   ts->data = (void*)rk;
1105 
1106   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr);
1107   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr);
1108   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",TSRKSetMultirate_RK);CHKERRQ(ierr);
1109   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",TSRKGetMultirate_RK);CHKERRQ(ierr);
1110 
1111   ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr);
1112   rk->dtratio = 1;
1113   PetscFunctionReturn(0);
1114 }
1115