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