xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision fdefd5e5da6933eb9789981da0487a06cd30939a)
1 /*
2   Code for timestepping with Runge-Kutta method
3 
4   Notes:
5   The general system is written as
6 
7   F(t,U,Udot) = G(t,U)
8 
9   where F represents the stiff part of the physics and G represents the non-stiff part.
10 
11 */
12 #include <petsc-private/tsimpl.h>                /*I   "petscts.h"   I*/
13 #include <petscdm.h>
14 
15 static TSRKType TSRKDefault = TSRK3;
16 static PetscBool     TSRKRegisterAllCalled;
17 static PetscBool     TSRKPackageInitialized;
18 static PetscInt      explicit_stage_time_id;
19 
20 typedef struct _RKTableau *RKTableau;
21 struct _RKTableau {
22   char      *name;
23   PetscInt   order;               /* Classical approximation order of the method */
24   PetscInt   s;                   /* Number of stages */
25   PetscInt   pinterp;             /* Interpolation order */
26   PetscReal *A,*b,*c;             /* Tableau */
27   PetscReal *bembed;              /* Embedded formula of order one less (order-1) */
28   PetscReal *binterp;             /* Dense output formula */
29   PetscReal  ccfl;                /* Placeholder for CFL coefficient relative to forward Euler */
30 };
31 typedef struct _RKTableauLink *RKTableauLink;
32 struct _RKTableauLink {
33   struct _RKTableau tab;
34   RKTableauLink     next;
35 };
36 static RKTableauLink RKTableauList;
37 
38 typedef struct {
39   RKTableau   tableau;
40   Vec          *Y;               /* States computed during the step */
41   Vec          *YdotRHS;         /* Function evaluations for the non-stiff part */
42   PetscScalar  *work;            /* Scalar work */
43   PetscReal    stage_time;
44   TSStepStatus status;
45 } TS_RK;
46 
47 /*MC
48      TSRK1 - First order forward Euler scheme.
49 
50      This method has one stage.
51 
52      Level: advanced
53 
54 .seealso: TSRK
55 M*/
56 /*MC
57      TSRK2 - Second order RK scheme.
58 
59      This method has two stages.
60 
61      Level: advanced
62 
63 .seealso: TSRK
64 M*/
65 /*MC
66      TSRK3 - Third order RK scheme.
67 
68      This method has three stages.
69 
70      Level: advanced
71 
72 .seealso: TSRK
73 M*/
74 /*MC
75      TSRK4 - Fourth order RK scheme.
76 
77      This method has four stages.
78 
79      Level: advanced
80 
81 .seealso: TSRK
82 M*/
83 /*MC
84      TSRK5F - Fifth order Fehlberg RK scheme with 4th order embedded method.
85 
86      This method has six stages.
87 
88      Level: advanced
89 
90 .seealso: TSRK
91 M*/
92 
93 #undef __FUNCT__
94 #define __FUNCT__ "TSRKRegisterAll"
95 /*@C
96   TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK
97 
98   Not Collective, but should be called by all processes which will need the schemes to be registered
99 
100   Level: advanced
101 
102 .keywords: TS, TSRK, register, all
103 
104 .seealso:  TSRKRegisterDestroy()
105 @*/
106 PetscErrorCode TSRKRegisterAll(void)
107 {
108   PetscErrorCode ierr;
109 
110   PetscFunctionBegin;
111   if (TSRKRegisterAllCalled) PetscFunctionReturn(0);
112   TSRKRegisterAllCalled = PETSC_TRUE;
113 
114   {
115     const PetscReal
116       A[1][1] = {{0.0}},
117       b[1]    = {1.0};
118     ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,1,b);CHKERRQ(ierr);
119   }
120   {
121     const PetscReal
122       A[2][2]     = {{0.0,0.0},
123                     {1.0,0.0}},
124       b[2]        = {0.5,0.5},
125       bembed[2]   = {1.0,0};
126     ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,2,b);CHKERRQ(ierr);
127   }
128   {
129     const PetscReal
130       A[3][3] = {{0,0,0},
131                  {2.0/3.0,0,0},
132                  {-1.0/3.0,1.0,0}},
133       b[3]    = {0.25,0.5,0.25};
134     ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,3,b);CHKERRQ(ierr);
135   }
136   {
137     const PetscReal
138       A[4][4] = {{0,0,0,0},
139                  {1.0/2.0,0,0,0},
140                  {0,3.0/4.0,0,0},
141                  {2.0/9.0,1.0/3.0,4.0/9.0,0}},
142       b[4]    = {2.0/9.0,1.0/3.0,4.0/9.0,0},
143       bembed[4] = {7.0/24.0,1.0/4.0,1.0/3.0,1.0/8.0};
144     ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,3,b);CHKERRQ(ierr);
145   }
146   {
147     const PetscReal
148       A[4][4] = {{0,0,0,0},
149                  {0.5,0,0,0},
150                  {0,0.5,0,0},
151                  {0,0,1.0,0}},
152       b[4]    = {1.0/6.0,1.0/3.0,1.0/3.0,1.0/6.0};
153     ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,4,b);CHKERRQ(ierr);
154   }
155   {
156     const PetscReal
157       A[6][6]   = {{0,0,0,0,0,0},
158                    {0.25,0,0,0,0,0},
159                    {3.0/32.0,9.0/32.0,0,0,0,0},
160                    {1932.0/2197.0,-7200.0/2197.0,7296.0/2197.0,0,0,0},
161                    {439.0/216.0,-8.0,3680.0/513.0,-845.0/4104.0,0,0},
162                    {-8.0/27.0,2.0,-3544.0/2565.0,1859.0/4104.0,-11.0/40.0,0}},
163       b[6]      = {16.0/135.0,0,6656.0/12825.0,28561.0/56430.0,-9.0/50.0,2.0/55.0},
164       bembed[6] = {25.0/216.0,0,1408.0/2565.0,2197.0/4104.0,-1.0/5.0,0};
165     ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,5,b);CHKERRQ(ierr);
166   }
167   {
168     const PetscReal
169       A[7][7]   = {{0,0,0,0,0,0,0},
170                    {1.0/5.0,0,0,0,0,0,0},
171                    {3.0/40.0,9.0/40.0,0,0,0,0,0},
172                    {44.0/45.0,-56.0/15.0,32.0/9.0,0,0,0,0},
173                    {19372.0/6561.0,-25360.0/2187.0,64448.0/6561.0,-212.0/729.0,0,0,0},
174                    {9017.0/3168.0,-355.0/33.0,46732.0/5247.0,49.0/176.0,-5103.0/18656.0,0,0},
175                    {35.0/384.0,0,500.0/1113.0,125.0/192.0,-2187.0/6784.0,11.0/84.0,0}},
176       b[7]      = {35.0/384.0,0,500.0/1113.0,125.0/192.0,-2187.0/6784.0,11.0/84.0,0},
177       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};
178     ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,5,b);CHKERRQ(ierr);
179   }
180   PetscFunctionReturn(0);
181 }
182 
183 #undef __FUNCT__
184 #define __FUNCT__ "TSRKRegisterDestroy"
185 /*@C
186    TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister().
187 
188    Not Collective
189 
190    Level: advanced
191 
192 .keywords: TSRK, register, destroy
193 .seealso: TSRKRegister(), TSRKRegisterAll()
194 @*/
195 PetscErrorCode TSRKRegisterDestroy(void)
196 {
197   PetscErrorCode ierr;
198   RKTableauLink link;
199 
200   PetscFunctionBegin;
201   while ((link = RKTableauList)) {
202     RKTableau t = &link->tab;
203     RKTableauList = link->next;
204     ierr = PetscFree3(t->A,t->b,t->c);  CHKERRQ(ierr);
205     ierr = PetscFree (t->bembed);       CHKERRQ(ierr);
206     ierr = PetscFree (t->binterp);      CHKERRQ(ierr);
207     ierr = PetscFree (t->name);         CHKERRQ(ierr);
208     ierr = PetscFree (link);            CHKERRQ(ierr);
209   }
210   TSRKRegisterAllCalled = PETSC_FALSE;
211   PetscFunctionReturn(0);
212 }
213 
214 #undef __FUNCT__
215 #define __FUNCT__ "TSRKInitializePackage"
216 /*@C
217   TSRKInitializePackage - This function initializes everything in the TSRK package. It is called
218   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RK()
219   when using static libraries.
220 
221   Level: developer
222 
223 .keywords: TS, TSRK, initialize, package
224 .seealso: PetscInitialize()
225 @*/
226 PetscErrorCode TSRKInitializePackage(void)
227 {
228   PetscErrorCode ierr;
229 
230   PetscFunctionBegin;
231   if (TSRKPackageInitialized) PetscFunctionReturn(0);
232   TSRKPackageInitialized = PETSC_TRUE;
233   ierr = TSRKRegisterAll();CHKERRQ(ierr);
234   ierr = PetscObjectComposedDataRegister(&explicit_stage_time_id);CHKERRQ(ierr);
235   ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr);
236   PetscFunctionReturn(0);
237 }
238 
239 #undef __FUNCT__
240 #define __FUNCT__ "TSRKFinalizePackage"
241 /*@C
242   TSRKFinalizePackage - This function destroys everything in the TSRK package. It is
243   called from PetscFinalize().
244 
245   Level: developer
246 
247 .keywords: Petsc, destroy, package
248 .seealso: PetscFinalize()
249 @*/
250 PetscErrorCode TSRKFinalizePackage(void)
251 {
252   PetscErrorCode ierr;
253 
254   PetscFunctionBegin;
255   TSRKPackageInitialized = PETSC_FALSE;
256   ierr = TSRKRegisterDestroy();CHKERRQ(ierr);
257   PetscFunctionReturn(0);
258 }
259 
260 #undef __FUNCT__
261 #define __FUNCT__ "TSRKRegister"
262 /*@C
263    TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
264 
265    Not Collective, but the same schemes should be registered on all processes on which they will be used
266 
267    Input Parameters:
268 +  name - identifier for method
269 .  order - approximation order of method
270 .  s - number of stages, this is the dimension of the matrices below
271 .  A - stage coefficients (dimension s*s, row-major)
272 .  b - step completion table (dimension s; NULL to use last row of A)
273 .  c - abscissa (dimension s; NULL to use row sums of A)
274 .  bembed - completion table for embedded method (dimension s; NULL if not available)
275 .  pinterp - Order of the interpolation scheme, equal to the number of columns of binterp
276 -  binterp - Coefficients of the interpolation formula (dimension s*pinterp; NULL to reuse binterpt)
277 
278    Notes:
279    Several RK methods are provided, this function is only needed to create new methods.
280 
281    Level: advanced
282 
283 .keywords: TS, register
284 
285 .seealso: TSRK
286 @*/
287 PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s,
288                                  const PetscReal A[],const PetscReal b[],const PetscReal c[],
289                                  const PetscReal bembed[],
290                                  PetscInt pinterp,const PetscReal binterp[])
291 {
292   PetscErrorCode  ierr;
293   RKTableauLink   link;
294   RKTableau       t;
295   PetscInt        i,j;
296 
297   PetscFunctionBegin;
298   ierr     = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
299   ierr     = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
300   t        = &link->tab;
301   ierr     = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
302   t->order = order;
303   t->s     = s;
304   ierr     = PetscMalloc3(s*s,PetscReal,&t->A,s,PetscReal,&t->b,s,PetscReal,&t->c);CHKERRQ(ierr);
305   ierr     = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
306   if (b)  { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); }
307   else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i];
308   if (c)  { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); }
309   else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j];
310 
311   if (bembed) {
312     ierr = PetscMalloc(s*sizeof(PetscReal),&t->bembed);CHKERRQ(ierr);
313     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
314   }
315 
316   t->pinterp     = pinterp;
317   ierr           = PetscMalloc(s*pinterp*sizeof(PetscReal),&t->binterp);CHKERRQ(ierr);
318   ierr           = PetscMemcpy(t->binterp,binterp,s*pinterp*sizeof(binterp[0]));CHKERRQ(ierr);
319   link->next     = RKTableauList;
320   RKTableauList = link;
321   PetscFunctionReturn(0);
322 }
323 
324 #undef __FUNCT__
325 #define __FUNCT__ "TSEvaluateStep_RK"
326 /*
327  The step completion formula is
328 
329  x1 = x0 + h b^T YdotRHS
330 
331  This function can be called before or after ts->vec_sol has been updated.
332  Suppose we have a completion formula (b) and an embedded formula (be) of different order.
333  We can write
334 
335  x1e = x0 + h be^T YdotRHS
336      = x1 - h b^T YdotRHS + h be^T YdotRHS
337      = x1 + h (be - b)^T YdotRHS
338 
339  so we can evaluate the method with different order even after the step has been optimistically completed.
340 */
341 static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done)
342 {
343   TS_RK         *rk   = (TS_RK*)ts->data;
344   RKTableau      tab  = rk->tableau;
345   PetscScalar   *w    = rk->work;
346   PetscReal      h;
347   PetscInt       s    = tab->s,j;
348   PetscErrorCode ierr;
349 
350   PetscFunctionBegin;
351   switch (rk->status) {
352   case TS_STEP_INCOMPLETE:
353   case TS_STEP_PENDING:
354     h = ts->time_step; break;
355   case TS_STEP_COMPLETE:
356     h = ts->time_step_prev; break;
357   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
358   }
359   if (order == tab->order) {
360     if (rk->status == TS_STEP_INCOMPLETE) {
361       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
362       for (j=0; j<s; j++) w[j] = h*tab->b[j];
363       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
364     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
365     if (done) *done = PETSC_TRUE;
366     PetscFunctionReturn(0);
367   } else if (order == tab->order-1) {
368     if (!tab->bembed) goto unavailable;
369     if (rk->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (be) */
370       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
371       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
372       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
373     } else {                    /* Rollback and re-complete using (be-b) */
374       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
375       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
376       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
377     }
378     if (done) *done = PETSC_TRUE;
379     PetscFunctionReturn(0);
380   }
381 unavailable:
382   if (done) *done = PETSC_FALSE;
383   else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"RK '%s' of order %D cannot evaluate step at order %D",tab->name,tab->order,order);
384   PetscFunctionReturn(0);
385 }
386 
387 #undef __FUNCT__
388 #define __FUNCT__ "TSStep_RK"
389 static PetscErrorCode TSStep_RK(TS ts)
390 {
391   TS_RK           *rk   = (TS_RK*)ts->data;
392   RKTableau        tab  = rk->tableau;
393   const PetscInt   s    = tab->s;
394   const PetscReal *A = tab->A,*b = tab->b,*c = tab->c;
395   PetscScalar     *w    = rk->work;
396   Vec             *Y    = rk->Y,*YdotRHS = rk->YdotRHS;
397   TSAdapt          adapt;
398   PetscInt         i,j,reject,next_scheme;
399   PetscReal        next_time_step;
400   PetscReal        t;
401   PetscBool        accept;
402   PetscErrorCode   ierr;
403 
404   PetscFunctionBegin;
405 
406   next_time_step = ts->time_step;
407   t              = ts->ptime;
408   accept         = PETSC_TRUE;
409   rk->status     = TS_STEP_INCOMPLETE;
410 
411 
412   for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) {
413     PetscReal h = ts->time_step;
414     ierr = TSPreStep(ts);CHKERRQ(ierr);
415     for (i=0; i<s; i++) {
416       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
417       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
418       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
419       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
420       ierr = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr);
421       if (!accept) goto reject_step;
422       ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
423     }
424     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
425     rk->status = TS_STEP_PENDING;
426 
427     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
428     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
429     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
430     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
431     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
432     if (accept) {
433       /* ignore next_scheme for now */
434       ts->ptime    += ts->time_step;
435       ts->time_step = next_time_step;
436       ts->steps++;
437       rk->status = TS_STEP_COMPLETE;
438       ierr = PetscObjectComposedDataSetReal((PetscObject)ts->vec_sol,explicit_stage_time_id,ts->ptime);CHKERRQ(ierr);
439 
440       break;
441     } else {                    /* Roll back the current step */
442       for (j=0; j<s; j++) w[j] = -h*b[j];
443       ierr = VecMAXPY(ts->vec_sol,s,w,rk->YdotRHS);CHKERRQ(ierr);
444       ts->time_step = next_time_step;
445       rk->status   = TS_STEP_INCOMPLETE;
446     }
447 reject_step: continue;
448   }
449   if (rk->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
450   PetscFunctionReturn(0);
451 }
452 
453 #undef __FUNCT__
454 #define __FUNCT__ "TSInterpolate_RK"
455 static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X)
456 {
457   TS_RK           *rk = (TS_RK*)ts->data;
458   PetscInt         s  = rk->tableau->s,pinterp = rk->tableau->pinterp,i,j;
459   PetscReal        h;
460   PetscReal        tt,t;
461   PetscScalar     *b;
462   const PetscReal *B = rk->tableau->binterp;
463   PetscErrorCode   ierr;
464 
465   PetscFunctionBegin;
466   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
467   switch (rk->status) {
468   case TS_STEP_INCOMPLETE:
469   case TS_STEP_PENDING:
470     h = ts->time_step;
471     t = (itime - ts->ptime)/h;
472     break;
473   case TS_STEP_COMPLETE:
474     h = ts->time_step_prev;
475     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
476     break;
477   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
478   }
479   ierr = PetscMalloc(s*sizeof(PetscScalar),&b);CHKERRQ(ierr);
480   for (i=0; i<s; i++) b[i] = 0;
481   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
482     for (i=0; i<s; i++) {
483       b[i]  += h * B[i*pinterp+j] * tt;
484     }
485   }
486   ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr);
487   ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr);
488   ierr = PetscFree(b);CHKERRQ(ierr);
489   PetscFunctionReturn(0);
490 }
491 
492 /*------------------------------------------------------------*/
493 #undef __FUNCT__
494 #define __FUNCT__ "TSReset_RK"
495 static PetscErrorCode TSReset_RK(TS ts)
496 {
497   TS_RK         *rk = (TS_RK*)ts->data;
498   PetscInt       s;
499   PetscErrorCode ierr;
500 
501   PetscFunctionBegin;
502   if (!rk->tableau) PetscFunctionReturn(0);
503   s    = rk->tableau->s;
504   ierr = VecDestroyVecs(s,&rk->Y);CHKERRQ(ierr);
505   ierr = VecDestroyVecs(s,&rk->YdotRHS);CHKERRQ(ierr);
506   ierr = PetscFree(rk->work);CHKERRQ(ierr);
507   PetscFunctionReturn(0);
508 }
509 
510 #undef __FUNCT__
511 #define __FUNCT__ "TSDestroy_RK"
512 static PetscErrorCode TSDestroy_RK(TS ts)
513 {
514   PetscErrorCode ierr;
515 
516   PetscFunctionBegin;
517   ierr = TSReset_RK(ts);CHKERRQ(ierr);
518   ierr = PetscFree(ts->data);CHKERRQ(ierr);
519   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr);
520   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr);
521   PetscFunctionReturn(0);
522 }
523 
524 
525 #undef __FUNCT__
526 #define __FUNCT__ "DMCoarsenHook_TSRK"
527 static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx)
528 {
529   PetscFunctionBegin;
530   PetscFunctionReturn(0);
531 }
532 
533 #undef __FUNCT__
534 #define __FUNCT__ "DMRestrictHook_TSRK"
535 static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
536 {
537   PetscFunctionBegin;
538   PetscFunctionReturn(0);
539 }
540 
541 
542 #undef __FUNCT__
543 #define __FUNCT__ "DMSubDomainHook_TSRK"
544 static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx)
545 {
546   PetscFunctionBegin;
547   PetscFunctionReturn(0);
548 }
549 
550 #undef __FUNCT__
551 #define __FUNCT__ "DMSubDomainRestrictHook_TSRK"
552 static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
553 {
554 
555   PetscFunctionBegin;
556   PetscFunctionReturn(0);
557 }
558 
559 #undef __FUNCT__
560 #define __FUNCT__ "TSSetUp_RK"
561 static PetscErrorCode TSSetUp_RK(TS ts)
562 {
563   TS_RK         *rk = (TS_RK*)ts->data;
564   RKTableau      tab;
565   PetscInt       s;
566   PetscErrorCode ierr;
567   DM             dm;
568 
569   PetscFunctionBegin;
570   if (!rk->tableau) {
571     ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr);
572   }
573   tab  = rk->tableau;
574   s    = tab->s;
575   ierr = VecDuplicateVecs(ts->vec_sol,s,&rk->Y);CHKERRQ(ierr);
576   ierr = VecDuplicateVecs(ts->vec_sol,s,&rk->YdotRHS);CHKERRQ(ierr);
577   ierr = PetscMalloc(s*sizeof(rk->work[0]),&rk->work);CHKERRQ(ierr);
578   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
579   if (dm) {
580     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
581     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
582   }
583   PetscFunctionReturn(0);
584 }
585 /*------------------------------------------------------------*/
586 
587 #undef __FUNCT__
588 #define __FUNCT__ "TSSetFromOptions_RK"
589 static PetscErrorCode TSSetFromOptions_RK(TS ts)
590 {
591   PetscErrorCode ierr;
592   char           rktype[256];
593 
594   PetscFunctionBegin;
595   ierr = PetscOptionsHead("RK ODE solver options");CHKERRQ(ierr);
596   {
597     RKTableauLink  link;
598     PetscInt       count,choice;
599     PetscBool      flg;
600     const char   **namelist;
601     ierr = PetscStrncpy(rktype,TSRKDefault,sizeof(rktype));CHKERRQ(ierr);
602     for (link=RKTableauList,count=0; link; link=link->next,count++) ;
603     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
604     for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
605     ierr      = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rktype,&choice,&flg);CHKERRQ(ierr);
606     ierr      = TSRKSetType(ts,flg ? namelist[choice] : rktype);CHKERRQ(ierr);
607     ierr      = PetscFree(namelist);CHKERRQ(ierr);
608   }
609   ierr = PetscOptionsTail();CHKERRQ(ierr);
610   PetscFunctionReturn(0);
611 }
612 
613 #undef __FUNCT__
614 #define __FUNCT__ "PetscFormatRealArray"
615 static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
616 {
617   PetscErrorCode ierr;
618   PetscInt       i;
619   size_t         left,count;
620   char           *p;
621 
622   PetscFunctionBegin;
623   for (i=0,p=buf,left=len; i<n; i++) {
624     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
625     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
626     left -= count;
627     p    += count;
628     *p++  = ' ';
629   }
630   p[i ? 0 : -1] = 0;
631   PetscFunctionReturn(0);
632 }
633 
634 #undef __FUNCT__
635 #define __FUNCT__ "TSView_RK"
636 static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
637 {
638   TS_RK         *rk   = (TS_RK*)ts->data;
639   RKTableau      tab  = rk->tableau;
640   PetscBool      iascii;
641   PetscErrorCode ierr;
642   TSAdapt        adapt;
643 
644   PetscFunctionBegin;
645   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
646   if (iascii) {
647     TSRKType rktype;
648     char     buf[512];
649     ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr);
650     ierr = PetscViewerASCIIPrintf(viewer,"  RK %s\n",rktype);CHKERRQ(ierr);
651     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
652     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa     c = %s\n",buf);CHKERRQ(ierr);
653   }
654   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
655   ierr = TSAdaptView(adapt,viewer);CHKERRQ(ierr);
656   PetscFunctionReturn(0);
657 }
658 
659 #undef __FUNCT__
660 #define __FUNCT__ "TSLoad_RK"
661 static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer)
662 {
663   PetscErrorCode ierr;
664   TSAdapt        tsadapt;
665 
666   PetscFunctionBegin;
667   ierr = TSGetAdapt(ts,&tsadapt);CHKERRQ(ierr);
668   ierr = TSAdaptLoad(tsadapt,viewer);CHKERRQ(ierr);
669   PetscFunctionReturn(0);
670 }
671 
672 #undef __FUNCT__
673 #define __FUNCT__ "TSRKSetType"
674 /*@C
675   TSRKSetType - Set the type of RK scheme
676 
677   Logically collective
678 
679   Input Parameter:
680 +  ts - timestepping context
681 -  rktype - type of RK-scheme
682 
683   Level: intermediate
684 
685 .seealso: TSRKGetType(), TSRK, TSRK2, TSRK3, TSRKPRSSP2, TSRK4, TSRK5
686 @*/
687 PetscErrorCode TSRKSetType(TS ts,TSRKType rktype)
688 {
689   PetscErrorCode ierr;
690 
691   PetscFunctionBegin;
692   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
693   ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr);
694   PetscFunctionReturn(0);
695 }
696 
697 #undef __FUNCT__
698 #define __FUNCT__ "TSRKGetType"
699 /*@C
700   TSRKGetType - Get the type of RK scheme
701 
702   Logically collective
703 
704   Input Parameter:
705 .  ts - timestepping context
706 
707   Output Parameter:
708 .  rktype - type of RK-scheme
709 
710   Level: intermediate
711 
712 .seealso: TSRKGetType()
713 @*/
714 PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype)
715 {
716   PetscErrorCode ierr;
717 
718   PetscFunctionBegin;
719   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
720   ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr);
721   PetscFunctionReturn(0);
722 }
723 
724 #undef __FUNCT__
725 #define __FUNCT__ "TSRKGetType_RK"
726 PetscErrorCode  TSRKGetType_RK(TS ts,TSRKType *rktype)
727 {
728   TS_RK     *rk = (TS_RK*)ts->data;
729   PetscErrorCode ierr;
730 
731   PetscFunctionBegin;
732   if (!rk->tableau) {
733     ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr);
734   }
735   *rktype = rk->tableau->name;
736   PetscFunctionReturn(0);
737 }
738 #undef __FUNCT__
739 #define __FUNCT__ "TSRKSetType_RK"
740 PetscErrorCode  TSRKSetType_RK(TS ts,TSRKType rktype)
741 {
742   TS_RK     *rk = (TS_RK*)ts->data;
743   PetscErrorCode ierr;
744   PetscBool      match;
745   RKTableauLink link;
746 
747   PetscFunctionBegin;
748   if (rk->tableau) {
749     ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr);
750     if (match) PetscFunctionReturn(0);
751   }
752   for (link = RKTableauList; link; link=link->next) {
753     ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr);
754     if (match) {
755       ierr = TSReset_RK(ts);CHKERRQ(ierr);
756       rk->tableau = &link->tab;
757       PetscFunctionReturn(0);
758     }
759   }
760   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype);
761   PetscFunctionReturn(0);
762 }
763 
764 /* ------------------------------------------------------------ */
765 /*MC
766       TSRK - ODE and DAE solver using Runge-Kutta schemes
767 
768   The user should provide the right hand side of the equation
769   using TSSetRHSFunction().
770 
771   Notes:
772   The default is TSRK3, it can be changed with TSRKSetType() or -ts_rk_type
773 
774   Level: beginner
775 
776 .seealso:  TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3,
777            TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister()
778 
779 M*/
780 #undef __FUNCT__
781 #define __FUNCT__ "TSCreate_RK"
782 PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
783 {
784   TS_RK     *th;
785   PetscErrorCode ierr;
786 
787   PetscFunctionBegin;
788 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
789   ierr = TSRKInitializePackage();CHKERRQ(ierr);
790 #endif
791 
792   ts->ops->reset          = TSReset_RK;
793   ts->ops->destroy        = TSDestroy_RK;
794   ts->ops->view           = TSView_RK;
795   ts->ops->load           = TSLoad_RK;
796   ts->ops->setup          = TSSetUp_RK;
797   ts->ops->step           = TSStep_RK;
798   ts->ops->interpolate    = TSInterpolate_RK;
799   ts->ops->evaluatestep   = TSEvaluateStep_RK;
800   ts->ops->setfromoptions = TSSetFromOptions_RK;
801 
802   ierr = PetscNewLog(ts,TS_RK,&th);CHKERRQ(ierr);
803   ts->data = (void*)th;
804 
805   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr);
806   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr);
807   PetscFunctionReturn(0);
808 }
809