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