xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision 6b3eee04a67a44bd9a4e65b41b1dd96ffcd06ca3)
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     = PetscNew(&link);CHKERRQ(ierr);
320   t        = &link->tab;
321   ierr     = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
322   t->order = order;
323   t->s     = s;
324   ierr     = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr);
325   ierr     = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
326   if (b)  { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); }
327   else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i];
328   if (c)  { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); }
329   else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j];
330   t->FSAL = PETSC_TRUE;
331   for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE;
332   if (bembed) {
333     ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr);
334     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
335   }
336 
337   t->pinterp     = pinterp;
338   ierr           = PetscMalloc1(s*pinterp,&t->binterp);CHKERRQ(ierr);
339   ierr           = PetscMemcpy(t->binterp,binterp,s*pinterp*sizeof(binterp[0]));CHKERRQ(ierr);
340   link->next     = RKTableauList;
341   RKTableauList = link;
342   PetscFunctionReturn(0);
343 }
344 
345 #undef __FUNCT__
346 #define __FUNCT__ "TSEvaluateStep_RK"
347 /*
348  The step completion formula is
349 
350  x1 = x0 + h b^T YdotRHS
351 
352  This function can be called before or after ts->vec_sol has been updated.
353  Suppose we have a completion formula (b) and an embedded formula (be) of different order.
354  We can write
355 
356  x1e = x0 + h be^T YdotRHS
357      = x1 - h b^T YdotRHS + h be^T YdotRHS
358      = x1 + h (be - b)^T YdotRHS
359 
360  so we can evaluate the method with different order even after the step has been optimistically completed.
361 */
362 static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done)
363 {
364   TS_RK         *rk   = (TS_RK*)ts->data;
365   RKTableau      tab  = rk->tableau;
366   PetscScalar   *w    = rk->work;
367   PetscReal      h;
368   PetscInt       s    = tab->s,j;
369   PetscErrorCode ierr;
370 
371   PetscFunctionBegin;
372   switch (rk->status) {
373   case TS_STEP_INCOMPLETE:
374   case TS_STEP_PENDING:
375     h = ts->time_step; break;
376   case TS_STEP_COMPLETE:
377     h = ts->ptime - ts->ptime_prev; break;
378   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
379   }
380   if (order == tab->order) {
381     if (rk->status == TS_STEP_INCOMPLETE) {
382       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
383       for (j=0; j<s; j++) w[j] = h*tab->b[j];
384       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
385     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
386     PetscFunctionReturn(0);
387   } else if (order == tab->order-1) {
388     if (!tab->bembed) goto unavailable;
389     if (rk->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (be) */
390       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
391       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
392       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
393     } else { /* Rollback and re-complete using (be-b) */
394       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
395       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
396       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
397       if (ts->vec_costintegral && ts->costintegralfwd) {
398         ierr = VecCopy(rk->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr);
399       }
400     }
401     if (done) *done = PETSC_TRUE;
402     PetscFunctionReturn(0);
403   }
404 unavailable:
405   if (done) *done = PETSC_FALSE;
406   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);
407   PetscFunctionReturn(0);
408 }
409 
410 #undef __FUNCT__
411 #define __FUNCT__ "TSForwardCostIntegral_RK"
412 static PetscErrorCode TSForwardCostIntegral_RK(TS ts)
413 {
414   TS_RK           *rk = (TS_RK*)ts->data;
415   RKTableau       tab = rk->tableau;
416   const PetscInt  s = tab->s;
417   const PetscReal *b = tab->b,*c = tab->c;
418   Vec             *Y = rk->Y;
419   PetscInt        i;
420   PetscErrorCode  ierr;
421 
422   PetscFunctionBegin;
423   /* backup cost integral */
424   ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr);
425   for (i=s-1; i>=0; i--) {
426     /* Evolve ts->vec_costintegral to compute integrals */
427     ierr = TSAdjointComputeCostIntegrand(ts,rk->ptime+rk->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
428     ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
429   }
430   PetscFunctionReturn(0);
431 }
432 
433 #undef __FUNCT__
434 #define __FUNCT__ "TSAdjointCostIntegral_RK"
435 static PetscErrorCode TSAdjointCostIntegral_RK(TS ts)
436 {
437   TS_RK           *rk = (TS_RK*)ts->data;
438   RKTableau       tab = rk->tableau;
439   const PetscInt  s = tab->s;
440   const PetscReal *b = tab->b,*c = tab->c;
441   Vec             *Y = rk->Y;
442   PetscInt        i;
443   PetscErrorCode  ierr;
444 
445   PetscFunctionBegin;
446   for (i=s-1; i>=0; i--) {
447     /* Evolve ts->vec_costintegral to compute integrals */
448     ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime-ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
449     ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
450   }
451   PetscFunctionReturn(0);
452 }
453 
454 #undef __FUNCT__
455 #define __FUNCT__ "TSRollBack_RK"
456 static PetscErrorCode TSRollBack_RK(TS ts)
457 {
458   TS_RK           *rk = (TS_RK*)ts->data;
459   RKTableau       tab = rk->tableau;
460   const PetscInt  s  = tab->s;
461   const PetscReal *b = tab->b;
462   PetscScalar     *w = rk->work;
463   Vec             *YdotRHS = rk->YdotRHS;
464   PetscInt        j;
465   PetscReal       h;
466   PetscErrorCode  ierr;
467 
468   PetscFunctionBegin;
469   switch (rk->status) {
470   case TS_STEP_INCOMPLETE:
471   case TS_STEP_PENDING:
472     h = ts->time_step; break;
473   case TS_STEP_COMPLETE:
474     h = ts->ptime - ts->ptime_prev; break;
475   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
476   }
477   for (j=0; j<s; j++) w[j] = -h*b[j];
478   ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr);
479   PetscFunctionReturn(0);
480 }
481 
482 #undef __FUNCT__
483 #define __FUNCT__ "TSStep_RK"
484 static PetscErrorCode TSStep_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 *A = tab->A,*c = tab->c;
490   PetscScalar     *w = rk->work;
491   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
492   TSAdapt          adapt;
493   PetscInt         i,j;
494   PetscInt         rejections = 0;
495   PetscBool        stageok,accept = PETSC_TRUE;
496   PetscReal        next_time_step = ts->time_step;
497   PetscErrorCode   ierr;
498 
499   PetscFunctionBegin;
500 
501   rk->status = TS_STEP_INCOMPLETE;
502   while (!ts->reason && rk->status != TS_STEP_COMPLETE) {
503     PetscReal t = ts->ptime;
504     PetscReal h = ts->time_step;
505     for (i=0; i<s; i++) {
506       rk->stage_time = t + h*c[i];
507       ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr);
508       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
509       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
510       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
511       ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
512       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
513       ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr);
514       if (!stageok) goto reject_step;
515       ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
516     }
517 
518     rk->status = TS_STEP_INCOMPLETE;
519     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
520     rk->status = TS_STEP_PENDING;
521     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
522     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
523     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
524     ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
525     rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
526     if (!accept) { /* Roll back the current step */
527       ierr = TSRollBack_RK(ts);CHKERRQ(ierr);
528       ts->time_step = next_time_step;
529       goto reject_step;
530     }
531 
532     if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation*/
533       rk->ptime     = ts->ptime;
534       rk->time_step = ts->time_step;
535     }
536 
537     ts->ptime += ts->time_step;
538     ts->time_step = next_time_step;
539     break;
540 
541   reject_step:
542     ts->reject++; accept = PETSC_FALSE;
543     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
544       ts->reason = TS_DIVERGED_STEP_REJECTED;
545       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
546     }
547   }
548   PetscFunctionReturn(0);
549 }
550 
551 #undef __FUNCT__
552 #define __FUNCT__ "TSAdjointSetUp_RK"
553 static PetscErrorCode TSAdjointSetUp_RK(TS ts)
554 {
555   TS_RK         *rk  = (TS_RK*)ts->data;
556   RKTableau      tab = rk->tableau;
557   PetscInt       s   = tab->s;
558   PetscErrorCode ierr;
559 
560   PetscFunctionBegin;
561   if (ts->adjointsetupcalled++) PetscFunctionReturn(0);
562   ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr);
563   if(ts->vecs_sensip) {
564     ierr = VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr);
565   }
566   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr);
567   PetscFunctionReturn(0);
568 }
569 
570 #undef __FUNCT__
571 #define __FUNCT__ "TSAdjointStep_RK"
572 static PetscErrorCode TSAdjointStep_RK(TS ts)
573 {
574   TS_RK           *rk   = (TS_RK*)ts->data;
575   RKTableau        tab  = rk->tableau;
576   const PetscInt   s    = tab->s;
577   const PetscReal *A = tab->A,*b = tab->b,*c = tab->c;
578   PetscScalar     *w    = rk->work;
579   Vec             *Y    = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu,*VecSensiTemp = rk->VecSensiTemp;
580   PetscInt         i,j,nadj;
581   PetscReal        t = ts->ptime;
582   PetscErrorCode   ierr;
583   PetscReal        h = ts->time_step;
584   Mat              J,Jp;
585 
586   PetscFunctionBegin;
587   rk->status = TS_STEP_INCOMPLETE;
588   for (i=s-1; i>=0; i--) {
589     rk->stage_time = t + h*(1.0-c[i]);
590     for (nadj=0; nadj<ts->numcost; nadj++) {
591       ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr);
592       ierr = VecScale(VecSensiTemp[nadj],-h*b[i]);CHKERRQ(ierr);
593       for (j=i+1; j<s; j++) {
594         ierr = VecAXPY(VecSensiTemp[nadj],-h*A[j*s+i],VecDeltaLam[nadj*s+j]);CHKERRQ(ierr);
595       }
596     }
597     /* Stage values of lambda */
598     ierr = TSGetRHSJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
599     ierr = TSComputeRHSJacobian(ts,rk->stage_time,Y[i],J,Jp);CHKERRQ(ierr);
600     if (ts->vec_costintegral) {
601       ierr = TSAdjointComputeDRDYFunction(ts,rk->stage_time,Y[i],ts->vecs_drdy);CHKERRQ(ierr);
602     }
603     for (nadj=0; nadj<ts->numcost; nadj++) {
604       ierr = MatMultTranspose(J,VecSensiTemp[nadj],VecDeltaLam[nadj*s+i]);CHKERRQ(ierr);
605       if (ts->vec_costintegral) {
606         ierr = VecAXPY(VecDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdy[nadj]);CHKERRQ(ierr);
607       }
608     }
609 
610     /* Stage values of mu */
611     if(ts->vecs_sensip) {
612       ierr = TSAdjointComputeRHSJacobian(ts,rk->stage_time,Y[i],ts->Jacp);CHKERRQ(ierr);
613       if (ts->vec_costintegral) {
614         ierr = TSAdjointComputeDRDPFunction(ts,rk->stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr);
615       }
616 
617       for (nadj=0; nadj<ts->numcost; nadj++) {
618         ierr = MatMultTranspose(ts->Jacp,VecSensiTemp[nadj],VecDeltaMu[nadj*s+i]);CHKERRQ(ierr);
619         if (ts->vec_costintegral) {
620           ierr = VecAXPY(VecDeltaMu[nadj*s+i],-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr);
621         }
622       }
623     }
624   }
625 
626   for (j=0; j<s; j++) w[j] = 1.0;
627   for (nadj=0; nadj<ts->numcost; nadj++) {
628     ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr);
629     if(ts->vecs_sensip) {
630       ierr = VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);CHKERRQ(ierr);
631     }
632   }
633   rk->status = TS_STEP_COMPLETE;
634   PetscFunctionReturn(0);
635 }
636 
637 #undef __FUNCT__
638 #define __FUNCT__ "TSInterpolate_RK"
639 static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X)
640 {
641   TS_RK           *rk = (TS_RK*)ts->data;
642   PetscInt         s  = rk->tableau->s,pinterp = rk->tableau->pinterp,i,j;
643   PetscReal        h;
644   PetscReal        tt,t;
645   PetscScalar     *b;
646   const PetscReal *B = rk->tableau->binterp;
647   PetscErrorCode   ierr;
648 
649   PetscFunctionBegin;
650   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
651 
652   switch (rk->status) {
653   case TS_STEP_INCOMPLETE:
654   case TS_STEP_PENDING:
655     h = ts->time_step;
656     t = (itime - ts->ptime)/h;
657     break;
658   case TS_STEP_COMPLETE:
659     h = ts->ptime - ts->ptime_prev;
660     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
661     break;
662   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
663   }
664   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
665   for (i=0; i<s; i++) b[i] = 0;
666   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
667     for (i=0; i<s; i++) {
668       b[i]  += h * B[i*pinterp+j] * tt;
669     }
670   }
671 
672   ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr);
673   ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr);
674 
675   ierr = PetscFree(b);CHKERRQ(ierr);
676   PetscFunctionReturn(0);
677 }
678 
679 /*------------------------------------------------------------*/
680 
681 #undef __FUNCT__
682 #define __FUNCT__ "TSRKTableauReset"
683 static PetscErrorCode TSRKTableauReset(TS ts)
684 {
685   TS_RK          *rk = (TS_RK*)ts->data;
686   RKTableau      tab = rk->tableau;
687   PetscErrorCode ierr;
688 
689   PetscFunctionBegin;
690   if (!tab) PetscFunctionReturn(0);
691   ierr = PetscFree(rk->work);CHKERRQ(ierr);
692   ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr);
693   ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr);
694   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr);
695   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr);
696   PetscFunctionReturn(0);
697 }
698 
699 #undef __FUNCT__
700 #define __FUNCT__ "TSReset_RK"
701 static PetscErrorCode TSReset_RK(TS ts)
702 {
703   TS_RK         *rk = (TS_RK*)ts->data;
704   PetscErrorCode ierr;
705 
706   PetscFunctionBegin;
707   ierr = TSRKTableauReset(ts);CHKERRQ(ierr);
708   ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr);
709   ierr = VecDestroyVecs(ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr);
710   PetscFunctionReturn(0);
711 }
712 
713 #undef __FUNCT__
714 #define __FUNCT__ "TSDestroy_RK"
715 static PetscErrorCode TSDestroy_RK(TS ts)
716 {
717   PetscErrorCode ierr;
718 
719   PetscFunctionBegin;
720   ierr = TSReset_RK(ts);CHKERRQ(ierr);
721   ierr = PetscFree(ts->data);CHKERRQ(ierr);
722   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr);
723   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr);
724   PetscFunctionReturn(0);
725 }
726 
727 
728 #undef __FUNCT__
729 #define __FUNCT__ "DMCoarsenHook_TSRK"
730 static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx)
731 {
732   PetscFunctionBegin;
733   PetscFunctionReturn(0);
734 }
735 
736 #undef __FUNCT__
737 #define __FUNCT__ "DMRestrictHook_TSRK"
738 static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
739 {
740   PetscFunctionBegin;
741   PetscFunctionReturn(0);
742 }
743 
744 
745 #undef __FUNCT__
746 #define __FUNCT__ "DMSubDomainHook_TSRK"
747 static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx)
748 {
749   PetscFunctionBegin;
750   PetscFunctionReturn(0);
751 }
752 
753 #undef __FUNCT__
754 #define __FUNCT__ "DMSubDomainRestrictHook_TSRK"
755 static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
756 {
757 
758   PetscFunctionBegin;
759   PetscFunctionReturn(0);
760 }
761 /*
762 #undef __FUNCT__
763 #define __FUNCT__ "RKSetAdjCoe"
764 static PetscErrorCode RKSetAdjCoe(RKTableau tab)
765 {
766   PetscReal *A,*b;
767   PetscInt        s,i,j;
768   PetscErrorCode  ierr;
769 
770   PetscFunctionBegin;
771   s    = tab->s;
772   ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr);
773 
774   for (i=0; i<s; i++)
775     for (j=0; j<s; j++) {
776       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];
777       ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr);
778     }
779 
780   for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i];
781 
782   ierr  = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
783   ierr  = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
784   ierr  = PetscFree2(A,b);CHKERRQ(ierr);
785   PetscFunctionReturn(0);
786 }
787 */
788 
789 #undef __FUNCT__
790 #define __FUNCT__ "TSRKTableauSetUp"
791 static PetscErrorCode TSRKTableauSetUp(TS ts)
792 {
793   TS_RK         *rk  = (TS_RK*)ts->data;
794   RKTableau      tab = rk->tableau;
795   PetscErrorCode ierr;
796 
797   PetscFunctionBegin;
798   ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr);
799   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr);
800   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr);
801   PetscFunctionReturn(0);
802 }
803 
804 
805 #undef __FUNCT__
806 #define __FUNCT__ "TSSetUp_RK"
807 static PetscErrorCode TSSetUp_RK(TS ts)
808 {
809   TS_RK         *rk = (TS_RK*)ts->data;
810   PetscErrorCode ierr;
811   DM             dm;
812 
813   PetscFunctionBegin;
814   ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);
815   if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */
816     ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr);
817   }
818   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
819   if (dm) {
820     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
821     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
822   }
823   PetscFunctionReturn(0);
824 }
825 
826 
827 /*------------------------------------------------------------*/
828 
829 #undef __FUNCT__
830 #define __FUNCT__ "TSSetFromOptions_RK"
831 static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts)
832 {
833   TS_RK         *rk = (TS_RK*)ts->data;
834   PetscErrorCode ierr;
835 
836   PetscFunctionBegin;
837   ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr);
838   {
839     RKTableauLink  link;
840     PetscInt       count,choice;
841     PetscBool      flg;
842     const char   **namelist;
843     for (link=RKTableauList,count=0; link; link=link->next,count++) ;
844     ierr = PetscMalloc1(count,&namelist);CHKERRQ(ierr);
845     for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
846     ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr);
847     if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);}
848     ierr = PetscFree(namelist);CHKERRQ(ierr);
849   }
850   ierr = PetscOptionsTail();CHKERRQ(ierr);
851   PetscFunctionReturn(0);
852 }
853 
854 #undef __FUNCT__
855 #define __FUNCT__ "PetscFormatRealArray"
856 static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
857 {
858   PetscErrorCode ierr;
859   PetscInt       i;
860   size_t         left,count;
861   char           *p;
862 
863   PetscFunctionBegin;
864   for (i=0,p=buf,left=len; i<n; i++) {
865     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
866     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
867     left -= count;
868     p    += count;
869     *p++  = ' ';
870   }
871   p[i ? 0 : -1] = 0;
872   PetscFunctionReturn(0);
873 }
874 
875 #undef __FUNCT__
876 #define __FUNCT__ "TSView_RK"
877 static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
878 {
879   TS_RK          *rk = (TS_RK*)ts->data;
880   PetscBool      iascii;
881   PetscErrorCode ierr;
882 
883   PetscFunctionBegin;
884   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
885   if (iascii) {
886     RKTableau tab  = rk->tableau;
887     TSRKType  rktype;
888     char      buf[512];
889     ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr);
890     ierr = PetscViewerASCIIPrintf(viewer,"  RK %s\n",rktype);CHKERRQ(ierr);
891     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
892     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa     c = %s\n",buf);CHKERRQ(ierr);
893     ierr = PetscViewerASCIIPrintf(viewer,"FSAL: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr);
894   }
895   if (ts->adapt) {ierr = TSAdaptView(ts->adapt,viewer);CHKERRQ(ierr);}
896   PetscFunctionReturn(0);
897 }
898 
899 #undef __FUNCT__
900 #define __FUNCT__ "TSLoad_RK"
901 static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer)
902 {
903   PetscErrorCode ierr;
904   TSAdapt        adapt;
905 
906   PetscFunctionBegin;
907   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
908   ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr);
909   PetscFunctionReturn(0);
910 }
911 
912 #undef __FUNCT__
913 #define __FUNCT__ "TSRKSetType"
914 /*@C
915   TSRKSetType - Set the type of RK scheme
916 
917   Logically collective
918 
919   Input Parameter:
920 +  ts - timestepping context
921 -  rktype - type of RK-scheme
922 
923   Level: intermediate
924 
925 .seealso: TSRKGetType(), TSRK, TSRK2, TSRK3, TSRKPRSSP2, TSRK4, TSRK5
926 @*/
927 PetscErrorCode TSRKSetType(TS ts,TSRKType rktype)
928 {
929   PetscErrorCode ierr;
930 
931   PetscFunctionBegin;
932   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
933   ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr);
934   PetscFunctionReturn(0);
935 }
936 
937 #undef __FUNCT__
938 #define __FUNCT__ "TSRKGetType"
939 /*@C
940   TSRKGetType - Get the type of RK scheme
941 
942   Logically collective
943 
944   Input Parameter:
945 .  ts - timestepping context
946 
947   Output Parameter:
948 .  rktype - type of RK-scheme
949 
950   Level: intermediate
951 
952 .seealso: TSRKGetType()
953 @*/
954 PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype)
955 {
956   PetscErrorCode ierr;
957 
958   PetscFunctionBegin;
959   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
960   ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr);
961   PetscFunctionReturn(0);
962 }
963 
964 #undef __FUNCT__
965 #define __FUNCT__ "TSRKGetType_RK"
966 static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype)
967 {
968   TS_RK *rk = (TS_RK*)ts->data;
969 
970   PetscFunctionBegin;
971   *rktype = rk->tableau->name;
972   PetscFunctionReturn(0);
973 }
974 #undef __FUNCT__
975 #define __FUNCT__ "TSRKSetType_RK"
976 static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype)
977 {
978   TS_RK          *rk = (TS_RK*)ts->data;
979   PetscErrorCode ierr;
980   PetscBool      match;
981   RKTableauLink  link;
982 
983   PetscFunctionBegin;
984   if (rk->tableau) {
985     ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr);
986     if (match) PetscFunctionReturn(0);
987   }
988   for (link = RKTableauList; link; link=link->next) {
989     ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr);
990     if (match) {
991       if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);}
992       rk->tableau = &link->tab;
993       if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);}
994       PetscFunctionReturn(0);
995     }
996   }
997   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype);
998   PetscFunctionReturn(0);
999 }
1000 
1001 #undef __FUNCT__
1002 #define __FUNCT__ "TSGetStages_RK"
1003 static PetscErrorCode  TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y)
1004 {
1005   TS_RK          *rk = (TS_RK*)ts->data;
1006 
1007   PetscFunctionBegin;
1008   *ns = rk->tableau->s;
1009   if(Y) *Y  = rk->Y;
1010   PetscFunctionReturn(0);
1011 }
1012 
1013 
1014 /* ------------------------------------------------------------ */
1015 /*MC
1016       TSRK - ODE and DAE solver using Runge-Kutta schemes
1017 
1018   The user should provide the right hand side of the equation
1019   using TSSetRHSFunction().
1020 
1021   Notes:
1022   The default is TSRK3, it can be changed with TSRKSetType() or -ts_rk_type
1023 
1024   Level: beginner
1025 
1026 .seealso:  TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3,
1027            TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister()
1028 
1029 M*/
1030 #undef __FUNCT__
1031 #define __FUNCT__ "TSCreate_RK"
1032 PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
1033 {
1034   TS_RK          *rk;
1035   PetscErrorCode ierr;
1036 
1037   PetscFunctionBegin;
1038   ierr = TSRKInitializePackage();CHKERRQ(ierr);
1039 
1040   ts->ops->reset          = TSReset_RK;
1041   ts->ops->destroy        = TSDestroy_RK;
1042   ts->ops->view           = TSView_RK;
1043   ts->ops->load           = TSLoad_RK;
1044   ts->ops->setup          = TSSetUp_RK;
1045   ts->ops->adjointsetup   = TSAdjointSetUp_RK;
1046   ts->ops->step           = TSStep_RK;
1047   ts->ops->interpolate    = TSInterpolate_RK;
1048   ts->ops->evaluatestep   = TSEvaluateStep_RK;
1049   ts->ops->rollback       = TSRollBack_RK;
1050   ts->ops->setfromoptions = TSSetFromOptions_RK;
1051   ts->ops->getstages      = TSGetStages_RK;
1052   ts->ops->adjointstep    = TSAdjointStep_RK;
1053 
1054   ts->ops->adjointintegral = TSAdjointCostIntegral_RK;
1055   ts->ops->forwardintegral = TSForwardCostIntegral_RK;
1056 
1057   ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr);
1058   ts->data = (void*)rk;
1059 
1060   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr);
1061   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr);
1062 
1063   ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr);
1064   PetscFunctionReturn(0);
1065 }
1066