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