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