xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision 922a638cd8109adf4f4c5d6be0a8a03ab3e43f13)
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 TSForwardStep_RK(TS ts)
507 {
508   TS_RK           *rk   = (TS_RK*)ts->data;
509   RKTableau       tab  = rk->tableau;
510   Mat             J,*MatsFwdSensipTemp = rk->MatsFwdSensipTemp;
511   const PetscInt  s = tab->s;
512   const PetscReal *A = tab->A,*c = tab->c,*b = tab->b;
513   Vec             *Y = rk->Y;
514   PetscInt        i,j;
515   PetscReal       stage_time,h = ts->time_step;
516   PetscBool       zero;
517   PetscErrorCode  ierr;
518 
519   PetscFunctionBegin;
520   ierr = MatCopy(ts->mat_sensip,rk->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
521   ierr = TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr);
522 
523   for (i=0; i<s; i++) {
524     stage_time = ts->ptime + h*c[i];
525     zero = PETSC_FALSE;
526     if (b[i] == 0 && i == s-1) zero = PETSC_TRUE;
527     /* TLM Stage values */
528     if(!i) {
529       ierr = MatCopy(ts->mat_sensip,rk->MatsFwdStageSensip[i],SAME_NONZERO_PATTERN);CHKERRQ(ierr);
530     } else if (!zero) {
531       ierr = MatZeroEntries(rk->MatsFwdStageSensip[i]);CHKERRQ(ierr);
532       for (j=0; j<i; j++) {
533         ierr = MatAXPY(rk->MatsFwdStageSensip[i],h*A[i*s+j],MatsFwdSensipTemp[j],SAME_NONZERO_PATTERN);CHKERRQ(ierr);
534       }
535       ierr = MatAXPY(rk->MatsFwdStageSensip[i],1.,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
536     } else {
537       ierr = MatZeroEntries(rk->MatsFwdStageSensip[i]);CHKERRQ(ierr);
538     }
539 
540     ierr = TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);CHKERRQ(ierr);
541     ierr = MatMatMult(J,rk->MatsFwdStageSensip[i],MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatsFwdSensipTemp[i]);CHKERRQ(ierr);
542     if (ts->Jacp) {
543       ierr = TSComputeRHSJacobianP(ts,stage_time,Y[i],ts->Jacp);CHKERRQ(ierr); /* get f_p */
544       ierr = MatAXPY(MatsFwdSensipTemp[i],1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
545     }
546   }
547 
548   for (i=0; i<s; i++) {
549     ierr = MatAXPY(ts->mat_sensip,h*b[i],rk->MatsFwdSensipTemp[i],SAME_NONZERO_PATTERN);CHKERRQ(ierr);
550   }
551   rk->status = TS_STEP_COMPLETE;
552   PetscFunctionReturn(0);
553 }
554 
555 static PetscErrorCode TSForwardGetStages_RK(TS ts,PetscInt *ns,Mat **stagesensip)
556 {
557   TS_RK     *rk = (TS_RK*)ts->data;
558   RKTableau tab  = rk->tableau;
559 
560   PetscFunctionBegin;
561   if (ns) *ns = tab->s;
562   if (stagesensip) *stagesensip = rk->MatsFwdStageSensip;
563   PetscFunctionReturn(0);
564 }
565 
566 static PetscErrorCode TSForwardSetUp_RK(TS ts)
567 {
568   TS_RK          *rk = (TS_RK*)ts->data;
569   RKTableau      tab  = rk->tableau;
570   PetscInt       i;
571   PetscErrorCode ierr;
572 
573   PetscFunctionBegin;
574   /* backup sensitivity results for roll-backs */
575   ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatFwdSensip0);CHKERRQ(ierr);
576 
577   ierr = PetscMalloc1(tab->s,&rk->MatsFwdStageSensip);CHKERRQ(ierr);
578   ierr = PetscMalloc1(tab->s,&rk->MatsFwdSensipTemp);CHKERRQ(ierr);
579   for(i=0; i<tab->s; i++) {
580     ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatsFwdStageSensip[i]);CHKERRQ(ierr);
581     ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatsFwdSensipTemp[i]);CHKERRQ(ierr);
582   }
583   ierr = VecDuplicate(ts->vec_sol,&rk->VecDeltaFwdSensipCol);CHKERRQ(ierr);
584   PetscFunctionReturn(0);
585 }
586 
587 static PetscErrorCode TSForwardReset_RK(TS ts)
588 {
589   TS_RK          *rk = (TS_RK*)ts->data;
590   RKTableau      tab  = rk->tableau;
591   PetscInt       i;
592   PetscErrorCode ierr;
593 
594   PetscFunctionBegin;
595   ierr = MatDestroy(&rk->MatFwdSensip0);CHKERRQ(ierr);
596   if (rk->MatsFwdStageSensip) {
597     for (i=0; i<tab->s; i++) {
598       ierr = MatDestroy(&rk->MatsFwdStageSensip[i]);CHKERRQ(ierr);
599     }
600     ierr = PetscFree(rk->MatsFwdStageSensip);CHKERRQ(ierr);
601   }
602   if (rk->MatsFwdSensipTemp) {
603     for (i=0; i<tab->s; i++) {
604       ierr = MatDestroy(&rk->MatsFwdSensipTemp[i]);CHKERRQ(ierr);
605     }
606     ierr = PetscFree(rk->MatsFwdSensipTemp);CHKERRQ(ierr);
607   }
608   ierr = VecDestroy(&rk->VecDeltaFwdSensipCol);CHKERRQ(ierr);
609   PetscFunctionReturn(0);
610 }
611 
612 static PetscErrorCode TSStep_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,*c = tab->c;
618   PetscScalar     *w = rk->work;
619   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
620   PetscBool       FSAL = tab->FSAL;
621   TSAdapt         adapt;
622   PetscInt        i,j;
623   PetscInt        rejections = 0;
624   PetscBool       stageok,accept = PETSC_TRUE;
625   PetscReal       next_time_step = ts->time_step;
626   PetscErrorCode  ierr;
627 
628   PetscFunctionBegin;
629   if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE;
630   if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); }
631 
632   rk->status = TS_STEP_INCOMPLETE;
633   while (!ts->reason && rk->status != TS_STEP_COMPLETE) {
634     PetscReal t = ts->ptime;
635     PetscReal h = ts->time_step;
636     for (i=0; i<s; i++) {
637       rk->stage_time = t + h*c[i];
638       ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr);
639       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
640       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
641       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
642       ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
643       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
644       ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr);
645       if (!stageok) goto reject_step;
646       if (FSAL && !i) continue;
647       ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
648     }
649 
650     rk->status = TS_STEP_INCOMPLETE;
651     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
652     rk->status = TS_STEP_PENDING;
653     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
654     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
655     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr);
656     ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
657     rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
658     if (!accept) { /* Roll back the current step */
659       ierr = TSRollBack_RK(ts);CHKERRQ(ierr);
660       ts->time_step = next_time_step;
661       goto reject_step;
662     }
663 
664     if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
665       rk->ptime     = ts->ptime;
666       rk->time_step = ts->time_step;
667     }
668 
669     ts->ptime += ts->time_step;
670     ts->time_step = next_time_step;
671     break;
672 
673     reject_step:
674     ts->reject++; accept = PETSC_FALSE;
675     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
676       ts->reason = TS_DIVERGED_STEP_REJECTED;
677       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
678     }
679   }
680   PetscFunctionReturn(0);
681 }
682 
683 static PetscErrorCode TSAdjointSetUp_RK(TS ts)
684 {
685   TS_RK          *rk  = (TS_RK*)ts->data;
686   RKTableau      tab = rk->tableau;
687   PetscInt       s   = tab->s;
688   PetscErrorCode ierr;
689 
690   PetscFunctionBegin;
691   if (ts->adjointsetupcalled++) PetscFunctionReturn(0);
692   ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecsDeltaLam);CHKERRQ(ierr);
693   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&rk->VecsSensiTemp);CHKERRQ(ierr);
694   if(ts->vecs_sensip) {
695     ierr = VecDuplicate(ts->vecs_sensip[0],&rk->VecDeltaMu);CHKERRQ(ierr);
696   }
697   PetscFunctionReturn(0);
698 }
699 
700 static PetscErrorCode TSAdjointStep_RK(TS ts)
701 {
702   TS_RK            *rk = (TS_RK*)ts->data;
703   RKTableau        tab = rk->tableau;
704   Mat              J;
705   const PetscInt   s = tab->s;
706   const PetscReal  *A = tab->A,*b = tab->b,*c = tab->c;
707   PetscScalar      *w = rk->work;
708   Vec              *Y = rk->Y,*VecsDeltaLam = rk->VecsDeltaLam,VecDeltaMu = rk->VecDeltaMu,*VecsSensiTemp = rk->VecsSensiTemp;
709   PetscInt         i,j,nadj;
710   PetscReal        t = ts->ptime;
711   PetscReal        h = ts->time_step,stage_time;
712   PetscBool        zero;
713   PetscErrorCode   ierr;
714 
715   PetscFunctionBegin;
716   rk->status = TS_STEP_INCOMPLETE;
717 
718   for (i=s-1; i>=0; i--) {
719     if (tab->FSAL && i == s-1) {
720       /* VecsDeltaLam[nadj*s+s-1] are initialized with zeros and the values never change.*/
721       continue;
722     }
723     stage_time = t + h*(1.0-c[i]);
724     zero = PETSC_FALSE;
725     ierr = TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr);
726     ierr = TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);CHKERRQ(ierr);
727     if (ts->vec_costintegral) {
728       ierr = TSComputeDRDUFunction(ts,stage_time,Y[i],ts->vecs_drdu);CHKERRQ(ierr);
729     }
730     if (ts->vecs_sensip) {
731       ierr = TSComputeRHSJacobianP(ts,stage_time,Y[i],ts->Jacp);CHKERRQ(ierr); /* get f_p */
732       if (ts->vec_costintegral) {
733         ierr = TSComputeDRDPFunction(ts,stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr);
734       }
735     }
736 
737     if (b[i] == 0 && i == s-1) zero = PETSC_TRUE;
738 
739     for (j=i+1; j<s; j++) w[j-i-1]  = -h*A[j*s+i]; /* coefficients for computing VecsSensiTemp */
740 
741     for (nadj=0; nadj<ts->numcost; nadj++) {
742       /* Stage values of lambda */
743       if (!zero) {
744         ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); /* VecDeltaLam is an vec array of size s by numcost */
745         ierr = VecScale(VecsSensiTemp[nadj],-h*b[i]);CHKERRQ(ierr);
746         ierr = VecMAXPY(VecsSensiTemp[nadj],s-i-1,w,&VecsDeltaLam[nadj*s+i+1]);CHKERRQ(ierr);
747         ierr = MatMultTranspose(J,VecsSensiTemp[nadj],VecsDeltaLam[nadj*s+i]);CHKERRQ(ierr);
748       } else {
749         /* \sum_{j=i+1}^s a_{ji}*lambda_{s,j} */
750         ierr = VecSet(VecsSensiTemp[nadj],0);CHKERRQ(ierr);
751         ierr = VecMAXPY(VecsSensiTemp[nadj],s-i-1,w,&VecsDeltaLam[nadj*s+i+1]);CHKERRQ(ierr);
752         ierr = MatMultTranspose(J,VecsSensiTemp[nadj],VecsDeltaLam[nadj*s+i]);CHKERRQ(ierr);
753         ierr = VecScale(VecsDeltaLam[nadj*s+i],-h);CHKERRQ(ierr);
754       }
755       if (ts->vec_costintegral) {
756         ierr = VecAXPY(VecsDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdu[nadj]);CHKERRQ(ierr);
757       }
758 
759       /* Stage values of mu */
760       if (ts->vecs_sensip) {
761         if (!zero) {
762           ierr = MatMultTranspose(ts->Jacp,VecsSensiTemp[nadj],VecDeltaMu);CHKERRQ(ierr);
763         } else {
764           ierr = VecSet(VecDeltaMu,0);CHKERRQ(ierr);
765         }
766         if (ts->vec_costintegral) {
767           ierr = VecAXPY(VecDeltaMu,-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr);
768         }
769         ierr = VecAXPY(ts->vecs_sensip[nadj],1.,VecDeltaMu);CHKERRQ(ierr); /* update sensip for each stage */
770       }
771     }
772   }
773 
774   for (j=0; j<s; j++) w[j] = 1.0;
775   for (nadj=0; nadj<ts->numcost; nadj++) { /* no need to do this for mu's */
776     ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecsDeltaLam[nadj*s]);CHKERRQ(ierr);
777   }
778   rk->status = TS_STEP_COMPLETE;
779   PetscFunctionReturn(0);
780 }
781 
782 static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X)
783 {
784   TS_RK            *rk = (TS_RK*)ts->data;
785   PetscInt         s  = rk->tableau->s,p = rk->tableau->p,i,j;
786   PetscReal        h;
787   PetscReal        tt,t;
788   PetscScalar      *b;
789   const PetscReal  *B = rk->tableau->binterp;
790   PetscErrorCode   ierr;
791 
792   PetscFunctionBegin;
793   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
794 
795   switch (rk->status) {
796     case TS_STEP_INCOMPLETE:
797     case TS_STEP_PENDING:
798       h = ts->time_step;
799       t = (itime - ts->ptime)/h;
800       break;
801     case TS_STEP_COMPLETE:
802       h = ts->ptime - ts->ptime_prev;
803       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
804       break;
805     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
806   }
807   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
808   for (i=0; i<s; i++) b[i] = 0;
809   for (j=0,tt=t; j<p; j++,tt*=t) {
810     for (i=0; i<s; i++) {
811       b[i]  += h * B[i*p+j] * tt;
812     }
813   }
814   ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr);
815   ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr);
816   ierr = PetscFree(b);CHKERRQ(ierr);
817   PetscFunctionReturn(0);
818 }
819 
820 /*------------------------------------------------------------*/
821 
822 static PetscErrorCode TSRKTableauReset(TS ts)
823 {
824   TS_RK          *rk = (TS_RK*)ts->data;
825   RKTableau      tab = rk->tableau;
826   PetscErrorCode ierr;
827 
828   PetscFunctionBegin;
829   if (!tab) PetscFunctionReturn(0);
830   ierr = PetscFree(rk->work);CHKERRQ(ierr);
831   ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr);
832   ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr);
833   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam);CHKERRQ(ierr);
834   ierr = VecDestroyVecs(ts->numcost,&rk->VecsSensiTemp);CHKERRQ(ierr);
835   ierr = VecDestroy(&rk->VecDeltaMu);CHKERRQ(ierr);
836   PetscFunctionReturn(0);
837 }
838 
839 static PetscErrorCode TSReset_RK(TS ts)
840 {
841   TS_RK         *rk = (TS_RK*)ts->data;
842   PetscErrorCode ierr;
843 
844   PetscFunctionBegin;
845   ierr = TSRKTableauReset(ts);CHKERRQ(ierr);
846   ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr);
847   if (ts->use_splitrhsfunction) {
848     ierr = PetscTryMethod(ts,"TSReset_RK_MultirateSplit_C",(TS),(ts));CHKERRQ(ierr);
849   } else {
850     ierr = PetscTryMethod(ts,"TSReset_RK_MultirateNonsplit_C",(TS),(ts));CHKERRQ(ierr);
851   }
852   PetscFunctionReturn(0);
853 }
854 
855 static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx)
856 {
857   PetscFunctionBegin;
858   PetscFunctionReturn(0);
859 }
860 
861 static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
862 {
863   PetscFunctionBegin;
864   PetscFunctionReturn(0);
865 }
866 
867 
868 static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx)
869 {
870   PetscFunctionBegin;
871   PetscFunctionReturn(0);
872 }
873 
874 static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
875 {
876 
877   PetscFunctionBegin;
878   PetscFunctionReturn(0);
879 }
880 /*
881 static PetscErrorCode RKSetAdjCoe(RKTableau tab)
882 {
883   PetscReal *A,*b;
884   PetscInt        s,i,j;
885   PetscErrorCode  ierr;
886 
887   PetscFunctionBegin;
888   s    = tab->s;
889   ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr);
890 
891   for (i=0; i<s; i++)
892     for (j=0; j<s; j++) {
893       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];
894       ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr);
895     }
896 
897   for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i];
898 
899   ierr  = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
900   ierr  = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
901   ierr  = PetscFree2(A,b);CHKERRQ(ierr);
902   PetscFunctionReturn(0);
903 }
904 */
905 
906 static PetscErrorCode TSRKTableauSetUp(TS ts)
907 {
908   TS_RK          *rk  = (TS_RK*)ts->data;
909   RKTableau      tab = rk->tableau;
910   PetscErrorCode ierr;
911 
912   PetscFunctionBegin;
913   ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr);
914   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr);
915   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr);
916   PetscFunctionReturn(0);
917 }
918 
919 static PetscErrorCode TSSetUp_RK(TS ts)
920 {
921   TS_RK          *rk = (TS_RK*)ts->data;
922   PetscErrorCode ierr;
923   DM             dm;
924 
925   PetscFunctionBegin;
926   ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr);
927   ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);
928   if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */
929     ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr);
930   }
931   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
932   ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
933   ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
934   if (ts->use_splitrhsfunction) {
935     ierr = PetscTryMethod(ts,"TSSetUp_RK_MultirateSplit_C",(TS),(ts));CHKERRQ(ierr);
936   } else {
937     ierr = PetscTryMethod(ts,"TSSetUp_RK_MultirateNonsplit_C",(TS),(ts));CHKERRQ(ierr);
938   }
939   PetscFunctionReturn(0);
940 }
941 
942 static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts)
943 {
944   TS_RK          *rk = (TS_RK*)ts->data;
945   PetscErrorCode ierr;
946 
947   PetscFunctionBegin;
948   ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr);
949   {
950     RKTableauLink link;
951     PetscInt      count,choice;
952     PetscBool     flg,use_multirate = PETSC_FALSE;
953     const char    **namelist;
954 
955     for (link=RKTableauList,count=0; link; link=link->next,count++) ;
956     ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr);
957     for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
958     ierr = PetscOptionsBool("-ts_rk_multirate","Use interpolation-based multirate RK method","TSRKSetMultirate",rk->use_multirate,&use_multirate,&flg);CHKERRQ(ierr);
959     if (flg) {
960       ierr = TSRKSetMultirate(ts,use_multirate);CHKERRQ(ierr);
961     }
962     ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr);
963     if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);}
964     ierr = PetscFree(namelist);CHKERRQ(ierr);
965   }
966   ierr = PetscOptionsTail();CHKERRQ(ierr);
967   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Multirate methods options","");CHKERRQ(ierr);
968   ierr = PetscOptionsInt("-ts_rk_dtratio","time step ratio between slow and fast","",rk->dtratio,&rk->dtratio,NULL);CHKERRQ(ierr);
969   ierr = PetscOptionsEnd();CHKERRQ(ierr);
970   PetscFunctionReturn(0);
971 }
972 
973 static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
974 {
975   TS_RK          *rk = (TS_RK*)ts->data;
976   PetscBool      iascii;
977   PetscErrorCode ierr;
978 
979   PetscFunctionBegin;
980   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
981   if (iascii) {
982     RKTableau tab  = rk->tableau;
983     TSRKType  rktype;
984     char      buf[512];
985     ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr);
986     ierr = PetscViewerASCIIPrintf(viewer,"  RK type %s\n",rktype);CHKERRQ(ierr);
987     ierr = PetscViewerASCIIPrintf(viewer,"  Order: %D\n",tab->order);CHKERRQ(ierr);
988     ierr = PetscViewerASCIIPrintf(viewer,"  FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr);
989     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
990     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa c = %s\n",buf);CHKERRQ(ierr);
991   }
992   PetscFunctionReturn(0);
993 }
994 
995 static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer)
996 {
997   PetscErrorCode ierr;
998   TSAdapt        adapt;
999 
1000   PetscFunctionBegin;
1001   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
1002   ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr);
1003   PetscFunctionReturn(0);
1004 }
1005 
1006 /*@C
1007   TSRKSetType - Set the type of RK scheme
1008 
1009   Logically collective
1010 
1011   Input Parameter:
1012 +  ts - timestepping context
1013 -  rktype - type of RK-scheme
1014 
1015   Options Database:
1016 .   -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs>
1017 
1018   Level: intermediate
1019 
1020 .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS
1021 @*/
1022 PetscErrorCode TSRKSetType(TS ts,TSRKType rktype)
1023 {
1024   PetscErrorCode ierr;
1025 
1026   PetscFunctionBegin;
1027   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1028   PetscValidCharPointer(rktype,2);
1029   ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr);
1030   PetscFunctionReturn(0);
1031 }
1032 
1033 /*@C
1034   TSRKGetType - Get the type of RK scheme
1035 
1036   Logically collective
1037 
1038   Input Parameter:
1039 .  ts - timestepping context
1040 
1041   Output Parameter:
1042 .  rktype - type of RK-scheme
1043 
1044   Level: intermediate
1045 
1046 .seealso: TSRKGetType()
1047 @*/
1048 PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype)
1049 {
1050   PetscErrorCode ierr;
1051 
1052   PetscFunctionBegin;
1053   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1054   ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr);
1055   PetscFunctionReturn(0);
1056 }
1057 
1058 static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype)
1059 {
1060   TS_RK *rk = (TS_RK*)ts->data;
1061 
1062   PetscFunctionBegin;
1063   *rktype = rk->tableau->name;
1064   PetscFunctionReturn(0);
1065 }
1066 
1067 static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype)
1068 {
1069   TS_RK          *rk = (TS_RK*)ts->data;
1070   PetscErrorCode ierr;
1071   PetscBool      match;
1072   RKTableauLink  link;
1073 
1074   PetscFunctionBegin;
1075   if (rk->tableau) {
1076     ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr);
1077     if (match) PetscFunctionReturn(0);
1078   }
1079   for (link = RKTableauList; link; link=link->next) {
1080     ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr);
1081     if (match) {
1082       if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);}
1083       rk->tableau = &link->tab;
1084       if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);}
1085       ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
1086       PetscFunctionReturn(0);
1087     }
1088   }
1089   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype);
1090   PetscFunctionReturn(0);
1091 }
1092 
1093 static PetscErrorCode  TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y)
1094 {
1095   TS_RK *rk = (TS_RK*)ts->data;
1096 
1097   PetscFunctionBegin;
1098   if (ns) *ns = rk->tableau->s;
1099   if (Y)   *Y = rk->Y;
1100   PetscFunctionReturn(0);
1101 }
1102 
1103 static PetscErrorCode TSDestroy_RK(TS ts)
1104 {
1105   PetscErrorCode ierr;
1106 
1107   PetscFunctionBegin;
1108   ierr = TSReset_RK(ts);CHKERRQ(ierr);
1109   if (ts->dm) {
1110     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
1111     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
1112   }
1113   ierr = PetscFree(ts->data);CHKERRQ(ierr);
1114   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr);
1115   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr);
1116   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",NULL);CHKERRQ(ierr);
1117   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",NULL);CHKERRQ(ierr);
1118   PetscFunctionReturn(0);
1119 }
1120 
1121 /*@C
1122   TSRKSetMultirate - Use the interpolation-based multirate RK method
1123 
1124   Logically collective
1125 
1126   Input Parameter:
1127 +  ts - timestepping context
1128 -  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
1129 
1130   Options Database:
1131 .   -ts_rk_multirate - <true,false>
1132 
1133   Notes:
1134   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).
1135 
1136   Level: intermediate
1137 
1138 .seealso: TSRKGetMultirate()
1139 @*/
1140 PetscErrorCode TSRKSetMultirate(TS ts,PetscBool use_multirate)
1141 {
1142   PetscErrorCode ierr;
1143 
1144   PetscFunctionBegin;
1145   ierr = PetscTryMethod(ts,"TSRKSetMultirate_C",(TS,PetscBool),(ts,use_multirate));CHKERRQ(ierr);
1146   PetscFunctionReturn(0);
1147 }
1148 
1149 /*@C
1150   TSRKGetMultirate - Gets whether to Use the interpolation-based multirate RK method
1151 
1152   Not collective
1153 
1154   Input Parameter:
1155 .  ts - timestepping context
1156 
1157   Output Parameter:
1158 .  use_multirate - PETSC_TRUE if the multirate RK method is enabled, PETSC_FALSE otherwise
1159 
1160   Level: intermediate
1161 
1162 .seealso: TSRKSetMultirate()
1163 @*/
1164 PetscErrorCode TSRKGetMultirate(TS ts,PetscBool *use_multirate)
1165 {
1166   PetscErrorCode ierr;
1167 
1168   PetscFunctionBegin;
1169   ierr = PetscUseMethod(ts,"TSRKGetMultirate_C",(TS,PetscBool*),(ts,use_multirate));CHKERRQ(ierr);
1170   PetscFunctionReturn(0);
1171 }
1172 
1173 /*MC
1174       TSRK - ODE and DAE solver using Runge-Kutta schemes
1175 
1176   The user should provide the right hand side of the equation
1177   using TSSetRHSFunction().
1178 
1179   Notes:
1180   The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type
1181 
1182   Level: beginner
1183 
1184 .seealso:  TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3,
1185            TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister(), TSRKSetMultirate(), TSRKGetMultirate()
1186 
1187 M*/
1188 PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
1189 {
1190   TS_RK          *rk;
1191   PetscErrorCode ierr;
1192 
1193   PetscFunctionBegin;
1194   ierr = TSRKInitializePackage();CHKERRQ(ierr);
1195 
1196   ts->ops->reset          = TSReset_RK;
1197   ts->ops->destroy        = TSDestroy_RK;
1198   ts->ops->view           = TSView_RK;
1199   ts->ops->load           = TSLoad_RK;
1200   ts->ops->setup          = TSSetUp_RK;
1201   ts->ops->adjointsetup   = TSAdjointSetUp_RK;
1202   ts->ops->interpolate    = TSInterpolate_RK;
1203   ts->ops->step           = TSStep_RK;
1204   ts->ops->evaluatestep   = TSEvaluateStep_RK;
1205   ts->ops->rollback       = TSRollBack_RK;
1206   ts->ops->setfromoptions = TSSetFromOptions_RK;
1207   ts->ops->getstages      = TSGetStages_RK;
1208   ts->ops->adjointstep    = TSAdjointStep_RK;
1209 
1210   ts->ops->adjointintegral = TSAdjointCostIntegral_RK;
1211   ts->ops->forwardintegral = TSForwardCostIntegral_RK;
1212 
1213   ts->ops->forwardsetup    = TSForwardSetUp_RK;
1214   ts->ops->forwardreset    = TSForwardReset_RK;
1215   ts->ops->forwardstep     = TSForwardStep_RK;
1216   ts->ops->forwardgetstages= TSForwardGetStages_RK;
1217 
1218   ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr);
1219   ts->data = (void*)rk;
1220 
1221   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr);
1222   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr);
1223   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",TSRKSetMultirate_RK);CHKERRQ(ierr);
1224   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",TSRKGetMultirate_RK);CHKERRQ(ierr);
1225 
1226   ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr);
1227   rk->dtratio = 1;
1228   PetscFunctionReturn(0);
1229 }
1230