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