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