xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision d61e7eb36d836a8e2bf1ad3f41344afcd0d5cf20)
1 #define PETSCTS_DLL
2 
3 /*
4  * Code for Timestepping with Runge Kutta
5  *
6  * Written by
7  * Asbjorn Hoiland Aarrestad
8  * asbjorn@aarrestad.com
9  * http://asbjorn.aarrestad.com/
10  *
11  */
12 #include "private/tsimpl.h"                /*I   "petscts.h"   I*/
13 #include "time.h"
14 
15 typedef struct {
16    Vec          y1,y2;  /* work wectors for the two rk permuations */
17    PetscInt     nok,nnok; /* counters for ok and not ok steps */
18    PetscReal    maxerror; /* variable to tell the maxerror allowed */
19    PetscReal    ferror; /* variable to tell (global maxerror)/(total time) */
20    PetscReal    tolerance; /* initial value set for maxerror by user */
21    Vec          tmp,tmp_y,*k; /* two temp vectors and the k vectors for rk */
22    PetscScalar  a[7][6]; /* rk scalars */
23    PetscScalar  b1[7],b2[7]; /* rk scalars */
24    PetscReal    c[7]; /* rk scalars */
25    PetscInt     p,s; /* variables to tell the size of the runge-kutta solver */
26 } TS_RK;
27 
28 EXTERN_C_BEGIN
29 #undef __FUNCT__
30 #define __FUNCT__ "TSRKSetTolerance_RK"
31 PetscErrorCode PETSCTS_DLLEXPORT TSRKSetTolerance_RK(TS ts,PetscReal aabs)
32 {
33   TS_RK *rk = (TS_RK*)ts->data;
34 
35   PetscFunctionBegin;
36   rk->tolerance = aabs;
37   PetscFunctionReturn(0);
38 }
39 EXTERN_C_END
40 
41 #undef __FUNCT__
42 #define __FUNCT__ "TSRKSetTolerance"
43 /*@
44    TSRKSetTolerance - Sets the total error the RK explicit time integrators
45                       will allow over the given time interval.
46 
47    Logically Collective on TS
48 
49    Input parameters:
50 +    ts  - the time-step context
51 -    aabs - the absolute tolerance
52 
53    Level: intermediate
54 
55 .keywords: RK, tolerance
56 
57 .seealso: TSSundialsSetTolerance()
58 
59 @*/
60 PetscErrorCode PETSCTS_DLLEXPORT TSRKSetTolerance(TS ts,PetscReal aabs)
61 {
62   PetscErrorCode ierr;
63 
64   PetscFunctionBegin;
65   PetscValidLogicalCollectiveReal(ts,aabs,2);
66   ierr = PetscTryMethod(ts,"TSRKSetTolerance_C",(TS,PetscReal),(ts,aabs));CHKERRQ(ierr);
67   PetscFunctionReturn(0);
68 }
69 
70 
71 #undef __FUNCT__
72 #define __FUNCT__ "TSSetUp_RK"
73 static PetscErrorCode TSSetUp_RK(TS ts)
74 {
75   TS_RK          *rk = (TS_RK*)ts->data;
76   PetscErrorCode ierr;
77 
78   PetscFunctionBegin;
79   rk->nok      = 0;
80   rk->nnok     = 0;
81   rk->maxerror = rk->tolerance;
82 
83   /* fixing maxerror: global vs local */
84   rk->ferror = rk->maxerror / (ts->max_time - ts->ptime);
85 
86   /* 34.0/45.0 gives double precision division */
87   /* defining variables needed for Runge-Kutta computing */
88   /* when changing below, please remember to change a, b1, b2 and c above! */
89   /* Found in table on page 171: Dormand-Prince 5(4) */
90 
91   /* are these right? */
92   rk->p=6;
93   rk->s=7;
94 
95   rk->a[1][0]=1.0/5.0;
96   rk->a[2][0]=3.0/40.0;
97   rk->a[2][1]=9.0/40.0;
98   rk->a[3][0]=44.0/45.0;
99   rk->a[3][1]=-56.0/15.0;
100   rk->a[3][2]=32.0/9.0;
101   rk->a[4][0]=19372.0/6561.0;
102   rk->a[4][1]=-25360.0/2187.0;
103   rk->a[4][2]=64448.0/6561.0;
104   rk->a[4][3]=-212.0/729.0;
105   rk->a[5][0]=9017.0/3168.0;
106   rk->a[5][1]=-355.0/33.0;
107   rk->a[5][2]=46732.0/5247.0;
108   rk->a[5][3]=49.0/176.0;
109   rk->a[5][4]=-5103.0/18656.0;
110   rk->a[6][0]=35.0/384.0;
111   rk->a[6][1]=0.0;
112   rk->a[6][2]=500.0/1113.0;
113   rk->a[6][3]=125.0/192.0;
114   rk->a[6][4]=-2187.0/6784.0;
115   rk->a[6][5]=11.0/84.0;
116 
117 
118   rk->c[0]=0.0;
119   rk->c[1]=1.0/5.0;
120   rk->c[2]=3.0/10.0;
121   rk->c[3]=4.0/5.0;
122   rk->c[4]=8.0/9.0;
123   rk->c[5]=1.0;
124   rk->c[6]=1.0;
125 
126   rk->b1[0]=35.0/384.0;
127   rk->b1[1]=0.0;
128   rk->b1[2]=500.0/1113.0;
129   rk->b1[3]=125.0/192.0;
130   rk->b1[4]=-2187.0/6784.0;
131   rk->b1[5]=11.0/84.0;
132   rk->b1[6]=0.0;
133 
134   rk->b2[0]=5179.0/57600.0;
135   rk->b2[1]=0.0;
136   rk->b2[2]=7571.0/16695.0;
137   rk->b2[3]=393.0/640.0;
138   rk->b2[4]=-92097.0/339200.0;
139   rk->b2[5]=187.0/2100.0;
140   rk->b2[6]=1.0/40.0;
141 
142 
143   /* Found in table on page 170: Fehlberg 4(5) */
144   /*
145   rk->p=5;
146   rk->s=6;
147 
148   rk->a[1][0]=1.0/4.0;
149   rk->a[2][0]=3.0/32.0;
150   rk->a[2][1]=9.0/32.0;
151   rk->a[3][0]=1932.0/2197.0;
152   rk->a[3][1]=-7200.0/2197.0;
153   rk->a[3][2]=7296.0/2197.0;
154   rk->a[4][0]=439.0/216.0;
155   rk->a[4][1]=-8.0;
156   rk->a[4][2]=3680.0/513.0;
157   rk->a[4][3]=-845.0/4104.0;
158   rk->a[5][0]=-8.0/27.0;
159   rk->a[5][1]=2.0;
160   rk->a[5][2]=-3544.0/2565.0;
161   rk->a[5][3]=1859.0/4104.0;
162   rk->a[5][4]=-11.0/40.0;
163 
164   rk->c[0]=0.0;
165   rk->c[1]=1.0/4.0;
166   rk->c[2]=3.0/8.0;
167   rk->c[3]=12.0/13.0;
168   rk->c[4]=1.0;
169   rk->c[5]=1.0/2.0;
170 
171   rk->b1[0]=25.0/216.0;
172   rk->b1[1]=0.0;
173   rk->b1[2]=1408.0/2565.0;
174   rk->b1[3]=2197.0/4104.0;
175   rk->b1[4]=-1.0/5.0;
176   rk->b1[5]=0.0;
177 
178   rk->b2[0]=16.0/135.0;
179   rk->b2[1]=0.0;
180   rk->b2[2]=6656.0/12825.0;
181   rk->b2[3]=28561.0/56430.0;
182   rk->b2[4]=-9.0/50.0;
183   rk->b2[5]=2.0/55.0;
184   */
185   /* Found in table on page 169: Merson 4("5") */
186   /*
187   rk->p=4;
188   rk->s=5;
189   rk->a[1][0] = 1.0/3.0;
190   rk->a[2][0] = 1.0/6.0;
191   rk->a[2][1] = 1.0/6.0;
192   rk->a[3][0] = 1.0/8.0;
193   rk->a[3][1] = 0.0;
194   rk->a[3][2] = 3.0/8.0;
195   rk->a[4][0] = 1.0/2.0;
196   rk->a[4][1] = 0.0;
197   rk->a[4][2] = -3.0/2.0;
198   rk->a[4][3] = 2.0;
199 
200   rk->c[0] = 0.0;
201   rk->c[1] = 1.0/3.0;
202   rk->c[2] = 1.0/3.0;
203   rk->c[3] = 0.5;
204   rk->c[4] = 1.0;
205 
206   rk->b1[0] = 1.0/2.0;
207   rk->b1[1] = 0.0;
208   rk->b1[2] = -3.0/2.0;
209   rk->b1[3] = 2.0;
210   rk->b1[4] = 0.0;
211 
212   rk->b2[0] = 1.0/6.0;
213   rk->b2[1] = 0.0;
214   rk->b2[2] = 0.0;
215   rk->b2[3] = 2.0/3.0;
216   rk->b2[4] = 1.0/6.0;
217   */
218 
219   /* making b2 -> e=b1-b2 */
220   /*
221     for(i=0;i<rk->s;i++){
222      rk->b2[i] = (rk->b1[i]) - (rk->b2[i]);
223   }
224   */
225   rk->b2[0]=71.0/57600.0;
226   rk->b2[1]=0.0;
227   rk->b2[2]=-71.0/16695.0;
228   rk->b2[3]=71.0/1920.0;
229   rk->b2[4]=-17253.0/339200.0;
230   rk->b2[5]=22.0/525.0;
231   rk->b2[6]=-1.0/40.0;
232 
233   /* initializing vectors */
234   ierr = VecDuplicate(ts->vec_sol,&rk->y1);CHKERRQ(ierr);
235   ierr = VecDuplicate(ts->vec_sol,&rk->y2);CHKERRQ(ierr);
236   ierr = VecDuplicate(rk->y1,&rk->tmp);CHKERRQ(ierr);
237   ierr = VecDuplicate(rk->y1,&rk->tmp_y);CHKERRQ(ierr);
238   ierr = VecDuplicateVecs(rk->y1,rk->s,&rk->k);CHKERRQ(ierr);
239 
240   PetscFunctionReturn(0);
241 }
242 
243 /*------------------------------------------------------------*/
244 #undef __FUNCT__
245 #define __FUNCT__ "TSRKqs"
246 PetscErrorCode TSRKqs(TS ts,PetscReal t,PetscReal h)
247 {
248   TS_RK          *rk = (TS_RK*)ts->data;
249   PetscErrorCode ierr;
250   PetscInt       j,l;
251   PetscReal      tmp_t=t;
252   PetscScalar    hh=h;
253 
254   PetscFunctionBegin;
255   /* k[0]=0  */
256   ierr = VecSet(rk->k[0],0.0);CHKERRQ(ierr);
257 
258   /* k[0] = derivs(t,y1) */
259   ierr = TSComputeRHSFunction(ts,t,rk->y1,rk->k[0]);CHKERRQ(ierr);
260   /* looping over runge-kutta variables */
261   /* building the k - array of vectors */
262   for(j = 1 ; j < rk->s ; j++){
263 
264      /* rk->tmp = 0 */
265      ierr = VecSet(rk->tmp,0.0);CHKERRQ(ierr);
266 
267      for(l=0;l<j;l++){
268         /* tmp += a(j,l)*k[l] */
269        ierr = VecAXPY(rk->tmp,rk->a[j][l],rk->k[l]);CHKERRQ(ierr);
270      }
271 
272      /* ierr = VecView(rk->tmp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
273 
274      /* k[j] = derivs(t+c(j)*h,y1+h*tmp,k(j)) */
275      /* I need the following helpers:
276         PetscScalar  tmp_t=t+c(j)*h
277         Vec          tmp_y=h*tmp+y1
278      */
279 
280      tmp_t = t + rk->c[j] * h;
281 
282      /* tmp_y = h * tmp + y1 */
283      ierr = VecWAXPY(rk->tmp_y,hh,rk->tmp,rk->y1);CHKERRQ(ierr);
284 
285      /* rk->k[j]=0 */
286      ierr = VecSet(rk->k[j],0.0);CHKERRQ(ierr);
287      ierr = TSComputeRHSFunction(ts,tmp_t,rk->tmp_y,rk->k[j]);CHKERRQ(ierr);
288   }
289 
290   /* tmp=0 and tmp_y=0 */
291   ierr = VecSet(rk->tmp,0.0);CHKERRQ(ierr);
292   ierr = VecSet(rk->tmp_y,0.0);CHKERRQ(ierr);
293 
294   for(j = 0 ; j < rk->s ; j++){
295      /* tmp=b1[j]*k[j]+tmp  */
296     ierr = VecAXPY(rk->tmp,rk->b1[j],rk->k[j]);CHKERRQ(ierr);
297      /* tmp_y=b2[j]*k[j]+tmp_y */
298     ierr = VecAXPY(rk->tmp_y,rk->b2[j],rk->k[j]);CHKERRQ(ierr);
299   }
300 
301   /* y2 = hh * tmp_y */
302   ierr = VecSet(rk->y2,0.0);CHKERRQ(ierr);
303   ierr = VecAXPY(rk->y2,hh,rk->tmp_y);CHKERRQ(ierr);
304   /* y1 = hh*tmp + y1 */
305   ierr = VecAXPY(rk->y1,hh,rk->tmp);CHKERRQ(ierr);
306   /* Finding difference between y1 and y2 */
307   PetscFunctionReturn(0);
308 }
309 
310 #undef __FUNCT__
311 #define __FUNCT__ "TSStep_RK"
312 static PetscErrorCode TSStep_RK(TS ts,PetscInt *steps,PetscReal *ptime)
313 {
314   TS_RK          *rk = (TS_RK*)ts->data;
315   PetscErrorCode ierr;
316   PetscReal      norm=0.0,dt_fac=0.0,fac = 0.0/*,ttmp=0.0*/;
317   PetscInt       i, max_steps = ts->max_steps;
318 
319   PetscFunctionBegin;
320   ierr=VecCopy(ts->vec_sol,rk->y1);CHKERRQ(ierr);
321   *steps = -ts->steps;
322   /* trying to save the vector */
323   ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
324   /* while loop to get from start to stop */
325   for (i = 0; i < max_steps; i++) {
326     ierr = TSPreStep(ts);CHKERRQ(ierr); /* Note that this is called once per STEP, not once per STAGE. */
327    /* calling rkqs */
328      /*
329        -- input
330        ts        - pointer to ts
331        ts->ptime - current time
332        ts->time_step        - try this timestep
333        y1        - solution for this step
334 
335        --output
336        y1        - suggested solution
337        y2        - check solution (runge - kutta second permutation)
338      */
339      ierr = TSRKqs(ts,ts->ptime,ts->time_step);CHKERRQ(ierr);
340      /* counting steps */
341      ts->steps++;
342    /* checking for maxerror */
343      /* comparing difference to maxerror */
344      ierr = VecNorm(rk->y2,NORM_2,&norm);CHKERRQ(ierr);
345      /* modifying maxerror to satisfy this timestep */
346      rk->maxerror = rk->ferror * ts->time_step;
347      /* ierr = PetscPrintf(PETSC_COMM_WORLD,"norm err: %f maxerror: %f dt: %f",norm,rk->maxerror,ts->time_step);CHKERRQ(ierr); */
348 
349    /* handling ok and not ok */
350      if (norm < rk->maxerror){
351         /* if ok: */
352         ierr=VecCopy(rk->y1,ts->vec_sol);CHKERRQ(ierr); /* saves the suggested solution to current solution */
353         ts->ptime += ts->time_step; /* storing the new current time */
354         rk->nok++;
355         fac=5.0;
356         /* trying to save the vector */
357         ierr = TSPostStep(ts);CHKERRQ(ierr);
358         ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
359         if (ts->ptime >= ts->max_time) break;
360      } else{
361         /* if not OK */
362         rk->nnok++;
363         fac=1.0;
364         ierr=VecCopy(ts->vec_sol,rk->y1);CHKERRQ(ierr);  /* restores old solution */
365      }
366 
367      /*Computing next stepsize. See page 167 in Solving ODE 1
368       *
369       * h_new = h * min( facmax , max( facmin , fac * (tol/err)^(1/(p+1)) ) )
370       * facmax set above
371       * facmin
372       */
373      dt_fac = exp(log((rk->maxerror) / norm) / ((rk->p) + 1) ) * 0.9 ;
374 
375      if (dt_fac > fac){
376         /*ierr = PetscPrintf(PETSC_COMM_WORLD,"changing fac %f\n",fac);*/
377         dt_fac = fac;
378      }
379 
380      /* computing new ts->time_step */
381      ts->time_step = ts->time_step * dt_fac;
382 
383      if (ts->ptime+ts->time_step > ts->max_time){
384         ts->time_step = ts->max_time - ts->ptime;
385      }
386 
387      if (ts->time_step < 1e-14){
388         ierr = PetscPrintf(PETSC_COMM_WORLD,"Very small steps: %f\n",ts->time_step);CHKERRQ(ierr);
389         ts->time_step = 1e-14;
390      }
391 
392      /* trying to purify h */
393      /* (did not give any visible result) */
394      /* ttmp = ts->ptime + ts->time_step;
395         ts->time_step = ttmp - ts->ptime; */
396 
397   }
398 
399   ierr=VecCopy(rk->y1,ts->vec_sol);CHKERRQ(ierr);
400   *steps += ts->steps;
401   *ptime  = ts->ptime;
402   PetscFunctionReturn(0);
403 }
404 
405 /*------------------------------------------------------------*/
406 #undef __FUNCT__
407 #define __FUNCT__ "TSDestroy_RK"
408 static PetscErrorCode TSDestroy_RK(TS ts)
409 {
410   TS_RK          *rk = (TS_RK*)ts->data;
411   PetscErrorCode ierr;
412   PetscInt       i;
413 
414   /* REMEMBER TO DESTROY ALL */
415 
416   PetscFunctionBegin;
417   if (rk->y1) {ierr = VecDestroy(rk->y1);CHKERRQ(ierr);}
418   if (rk->y2) {ierr = VecDestroy(rk->y2);CHKERRQ(ierr);}
419   if (rk->tmp) {ierr = VecDestroy(rk->tmp);CHKERRQ(ierr);}
420   if (rk->tmp_y) {ierr = VecDestroy(rk->tmp_y);CHKERRQ(ierr);}
421   for(i=0;i<rk->s;i++){
422      if (rk->k[i]) {ierr = VecDestroy(rk->k[i]);CHKERRQ(ierr);}
423   }
424   ierr = PetscFree(rk);CHKERRQ(ierr);
425   PetscFunctionReturn(0);
426 }
427 /*------------------------------------------------------------*/
428 
429 #undef __FUNCT__
430 #define __FUNCT__ "TSSetFromOptions_RK"
431 static PetscErrorCode TSSetFromOptions_RK(TS ts)
432 {
433   TS_RK          *rk = (TS_RK*)ts->data;
434   PetscErrorCode ierr;
435 
436   PetscFunctionBegin;
437   ierr = PetscOptionsHead("RK ODE solver options");CHKERRQ(ierr);
438     ierr = PetscOptionsReal("-ts_rk_tol","Tolerance for convergence","TSRKSetTolerance",rk->tolerance,&rk->tolerance,PETSC_NULL);CHKERRQ(ierr);
439   ierr = PetscOptionsTail();CHKERRQ(ierr);
440   PetscFunctionReturn(0);
441 }
442 
443 #undef __FUNCT__
444 #define __FUNCT__ "TSView_RK"
445 static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
446 {
447    TS_RK          *rk = (TS_RK*)ts->data;
448    PetscErrorCode ierr;
449 
450    PetscFunctionBegin;
451    ierr = PetscPrintf(PETSC_COMM_WORLD,"  number of ok steps: %D\n",rk->nok);CHKERRQ(ierr);
452    ierr = PetscPrintf(PETSC_COMM_WORLD,"  number of rejected steps: %D\n",rk->nnok);CHKERRQ(ierr);
453    PetscFunctionReturn(0);
454 }
455 
456 /* ------------------------------------------------------------ */
457 /*MC
458       TSRK - ODE solver using the explicit Runge-Kutta methods
459 
460    Options Database:
461 .  -ts_rk_tol <tol> Tolerance for convergence
462 
463   Contributed by: Asbjorn Hoiland Aarrestad, asbjorn@aarrestad.com, http://asbjorn.aarrestad.com/
464 
465   Level: beginner
466 
467 .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSRKSetTolerance()
468 
469 M*/
470 
471 EXTERN_C_BEGIN
472 #undef __FUNCT__
473 #define __FUNCT__ "TSCreate_RK"
474 PetscErrorCode PETSCTS_DLLEXPORT TSCreate_RK(TS ts)
475 {
476   TS_RK          *rk;
477   PetscErrorCode ierr;
478 
479   PetscFunctionBegin;
480   ts->ops->setup           = TSSetUp_RK;
481   ts->ops->step            = TSStep_RK;
482   ts->ops->destroy         = TSDestroy_RK;
483   ts->ops->setfromoptions  = TSSetFromOptions_RK;
484   ts->ops->view            = TSView_RK;
485 
486   ierr = PetscNewLog(ts,TS_RK,&rk);CHKERRQ(ierr);
487   ts->data = (void*)rk;
488 
489   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRKSetTolerance_C","TSRKSetTolerance_RK",TSRKSetTolerance_RK);CHKERRQ(ierr);
490 
491   PetscFunctionReturn(0);
492 }
493 EXTERN_C_END
494 
495 
496 
497 
498