xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision 186e87aca15291e38b957d29c5c6a67aa326bcab)
1e4dd521cSBarry Smith /*
226edb414Sasbjorn  * Code for Timestepping with Runge Kutta
326edb414Sasbjorn  *
426edb414Sasbjorn  * Written by
526edb414Sasbjorn  * Asbjorn Hoiland Aarrestad
626edb414Sasbjorn  * asbjorn@aarrestad.com
726edb414Sasbjorn  * http://asbjorn.aarrestad.com/
826edb414Sasbjorn  *
9e4dd521cSBarry Smith  */
10c6db04a5SJed Brown #include <private/tsimpl.h>                /*I   "petscts.h"   I*/
11c6db04a5SJed Brown #include <time.h>
12e4dd521cSBarry Smith 
13e4dd521cSBarry Smith typedef struct {
14e4dd521cSBarry Smith    Vec          y1,y2;  /* work wectors for the two rk permuations */
15a7cc72afSBarry Smith    PetscInt     nok,nnok; /* counters for ok and not ok steps */
16e4dd521cSBarry Smith    PetscReal    maxerror; /* variable to tell the maxerror allowed */
17e4dd521cSBarry Smith    PetscReal    ferror; /* variable to tell (global maxerror)/(total time) */
18262ff7bbSBarry Smith    PetscReal    tolerance; /* initial value set for maxerror by user */
19e4dd521cSBarry Smith    Vec          tmp,tmp_y,*k; /* two temp vectors and the k vectors for rk */
20e4dd521cSBarry Smith    PetscScalar  a[7][6]; /* rk scalars */
21e4dd521cSBarry Smith    PetscScalar  b1[7],b2[7]; /* rk scalars */
22e4dd521cSBarry Smith    PetscReal    c[7]; /* rk scalars */
23a7cc72afSBarry Smith    PetscInt     p,s; /* variables to tell the size of the runge-kutta solver */
245f70b478SJed Brown } TS_RK;
25e4dd521cSBarry Smith 
26262ff7bbSBarry Smith EXTERN_C_BEGIN
27262ff7bbSBarry Smith #undef __FUNCT__
28262ff7bbSBarry Smith #define __FUNCT__ "TSRKSetTolerance_RK"
297087cfbeSBarry Smith PetscErrorCode  TSRKSetTolerance_RK(TS ts,PetscReal aabs)
30262ff7bbSBarry Smith {
315f70b478SJed Brown   TS_RK *rk = (TS_RK*)ts->data;
32262ff7bbSBarry Smith 
33262ff7bbSBarry Smith   PetscFunctionBegin;
34262ff7bbSBarry Smith   rk->tolerance = aabs;
35262ff7bbSBarry Smith   PetscFunctionReturn(0);
36262ff7bbSBarry Smith }
37262ff7bbSBarry Smith EXTERN_C_END
38262ff7bbSBarry Smith 
39262ff7bbSBarry Smith #undef __FUNCT__
40262ff7bbSBarry Smith #define __FUNCT__ "TSRKSetTolerance"
41262ff7bbSBarry Smith /*@
42262ff7bbSBarry Smith    TSRKSetTolerance - Sets the total error the RK explicit time integrators
43262ff7bbSBarry Smith                       will allow over the given time interval.
44262ff7bbSBarry Smith 
453f9fe445SBarry Smith    Logically Collective on TS
46262ff7bbSBarry Smith 
47262ff7bbSBarry Smith    Input parameters:
48262ff7bbSBarry Smith +    ts  - the time-step context
49262ff7bbSBarry Smith -    aabs - the absolute tolerance
50262ff7bbSBarry Smith 
51262ff7bbSBarry Smith    Level: intermediate
52262ff7bbSBarry Smith 
53262ff7bbSBarry Smith .keywords: RK, tolerance
54262ff7bbSBarry Smith 
550f3b3ca1SHong Zhang .seealso: TSSundialsSetTolerance()
56262ff7bbSBarry Smith 
57262ff7bbSBarry Smith @*/
587087cfbeSBarry Smith PetscErrorCode  TSRKSetTolerance(TS ts,PetscReal aabs)
59262ff7bbSBarry Smith {
604ac538c5SBarry Smith   PetscErrorCode ierr;
61262ff7bbSBarry Smith 
62262ff7bbSBarry Smith   PetscFunctionBegin;
63c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,aabs,2);
644ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSRKSetTolerance_C",(TS,PetscReal),(ts,aabs));CHKERRQ(ierr);
65262ff7bbSBarry Smith   PetscFunctionReturn(0);
66262ff7bbSBarry Smith }
67e4dd521cSBarry Smith 
68e4dd521cSBarry Smith 
69e4dd521cSBarry Smith #undef __FUNCT__
705f70b478SJed Brown #define __FUNCT__ "TSSetUp_RK"
715f70b478SJed Brown static PetscErrorCode TSSetUp_RK(TS ts)
72e4dd521cSBarry Smith {
735f70b478SJed Brown   TS_RK          *rk = (TS_RK*)ts->data;
74a7cc72afSBarry Smith   PetscErrorCode ierr;
75e4dd521cSBarry Smith 
76e4dd521cSBarry Smith   PetscFunctionBegin;
77e4dd521cSBarry Smith   rk->nok      = 0;
78e4dd521cSBarry Smith   rk->nnok     = 0;
79262ff7bbSBarry Smith   rk->maxerror = rk->tolerance;
80e4dd521cSBarry Smith 
81e4dd521cSBarry Smith   /* fixing maxerror: global vs local */
82e4dd521cSBarry Smith   rk->ferror = rk->maxerror / (ts->max_time - ts->ptime);
83e4dd521cSBarry Smith 
84e4dd521cSBarry Smith   /* 34.0/45.0 gives double precision division */
85e4dd521cSBarry Smith   /* defining variables needed for Runge-Kutta computing */
86e4dd521cSBarry Smith   /* when changing below, please remember to change a, b1, b2 and c above! */
87e4dd521cSBarry Smith   /* Found in table on page 171: Dormand-Prince 5(4) */
88e4dd521cSBarry Smith 
89e4dd521cSBarry Smith   /* are these right? */
90e4dd521cSBarry Smith   rk->p=6;
91e4dd521cSBarry Smith   rk->s=7;
92e4dd521cSBarry Smith 
93e4dd521cSBarry Smith   rk->a[1][0]=1.0/5.0;
94e4dd521cSBarry Smith   rk->a[2][0]=3.0/40.0;
95e4dd521cSBarry Smith   rk->a[2][1]=9.0/40.0;
96e4dd521cSBarry Smith   rk->a[3][0]=44.0/45.0;
97e4dd521cSBarry Smith   rk->a[3][1]=-56.0/15.0;
98e4dd521cSBarry Smith   rk->a[3][2]=32.0/9.0;
99e4dd521cSBarry Smith   rk->a[4][0]=19372.0/6561.0;
100e4dd521cSBarry Smith   rk->a[4][1]=-25360.0/2187.0;
101e4dd521cSBarry Smith   rk->a[4][2]=64448.0/6561.0;
102e4dd521cSBarry Smith   rk->a[4][3]=-212.0/729.0;
103e4dd521cSBarry Smith   rk->a[5][0]=9017.0/3168.0;
104e4dd521cSBarry Smith   rk->a[5][1]=-355.0/33.0;
105e4dd521cSBarry Smith   rk->a[5][2]=46732.0/5247.0;
106e4dd521cSBarry Smith   rk->a[5][3]=49.0/176.0;
107e4dd521cSBarry Smith   rk->a[5][4]=-5103.0/18656.0;
108e4dd521cSBarry Smith   rk->a[6][0]=35.0/384.0;
109e4dd521cSBarry Smith   rk->a[6][1]=0.0;
110e4dd521cSBarry Smith   rk->a[6][2]=500.0/1113.0;
111e4dd521cSBarry Smith   rk->a[6][3]=125.0/192.0;
112e4dd521cSBarry Smith   rk->a[6][4]=-2187.0/6784.0;
113e4dd521cSBarry Smith   rk->a[6][5]=11.0/84.0;
114e4dd521cSBarry Smith 
115e4dd521cSBarry Smith 
116e4dd521cSBarry Smith   rk->c[0]=0.0;
117e4dd521cSBarry Smith   rk->c[1]=1.0/5.0;
118e4dd521cSBarry Smith   rk->c[2]=3.0/10.0;
119e4dd521cSBarry Smith   rk->c[3]=4.0/5.0;
120e4dd521cSBarry Smith   rk->c[4]=8.0/9.0;
121e4dd521cSBarry Smith   rk->c[5]=1.0;
122e4dd521cSBarry Smith   rk->c[6]=1.0;
123e4dd521cSBarry Smith 
124e4dd521cSBarry Smith   rk->b1[0]=35.0/384.0;
125e4dd521cSBarry Smith   rk->b1[1]=0.0;
126e4dd521cSBarry Smith   rk->b1[2]=500.0/1113.0;
127e4dd521cSBarry Smith   rk->b1[3]=125.0/192.0;
128e4dd521cSBarry Smith   rk->b1[4]=-2187.0/6784.0;
129e4dd521cSBarry Smith   rk->b1[5]=11.0/84.0;
130e4dd521cSBarry Smith   rk->b1[6]=0.0;
131e4dd521cSBarry Smith 
132e4dd521cSBarry Smith   rk->b2[0]=5179.0/57600.0;
133e4dd521cSBarry Smith   rk->b2[1]=0.0;
134e4dd521cSBarry Smith   rk->b2[2]=7571.0/16695.0;
135e4dd521cSBarry Smith   rk->b2[3]=393.0/640.0;
136e4dd521cSBarry Smith   rk->b2[4]=-92097.0/339200.0;
137e4dd521cSBarry Smith   rk->b2[5]=187.0/2100.0;
138e4dd521cSBarry Smith   rk->b2[6]=1.0/40.0;
139e4dd521cSBarry Smith 
140e4dd521cSBarry Smith 
141e4dd521cSBarry Smith   /* Found in table on page 170: Fehlberg 4(5) */
142e4dd521cSBarry Smith   /*
143e4dd521cSBarry Smith   rk->p=5;
144e4dd521cSBarry Smith   rk->s=6;
145e4dd521cSBarry Smith 
146e4dd521cSBarry Smith   rk->a[1][0]=1.0/4.0;
147e4dd521cSBarry Smith   rk->a[2][0]=3.0/32.0;
148e4dd521cSBarry Smith   rk->a[2][1]=9.0/32.0;
149e4dd521cSBarry Smith   rk->a[3][0]=1932.0/2197.0;
150e4dd521cSBarry Smith   rk->a[3][1]=-7200.0/2197.0;
151e4dd521cSBarry Smith   rk->a[3][2]=7296.0/2197.0;
152e4dd521cSBarry Smith   rk->a[4][0]=439.0/216.0;
153e4dd521cSBarry Smith   rk->a[4][1]=-8.0;
154e4dd521cSBarry Smith   rk->a[4][2]=3680.0/513.0;
155e4dd521cSBarry Smith   rk->a[4][3]=-845.0/4104.0;
156e4dd521cSBarry Smith   rk->a[5][0]=-8.0/27.0;
157e4dd521cSBarry Smith   rk->a[5][1]=2.0;
158e4dd521cSBarry Smith   rk->a[5][2]=-3544.0/2565.0;
159e4dd521cSBarry Smith   rk->a[5][3]=1859.0/4104.0;
160e4dd521cSBarry Smith   rk->a[5][4]=-11.0/40.0;
161e4dd521cSBarry Smith 
162e4dd521cSBarry Smith   rk->c[0]=0.0;
163e4dd521cSBarry Smith   rk->c[1]=1.0/4.0;
164e4dd521cSBarry Smith   rk->c[2]=3.0/8.0;
165e4dd521cSBarry Smith   rk->c[3]=12.0/13.0;
166e4dd521cSBarry Smith   rk->c[4]=1.0;
167e4dd521cSBarry Smith   rk->c[5]=1.0/2.0;
168e4dd521cSBarry Smith 
169e4dd521cSBarry Smith   rk->b1[0]=25.0/216.0;
170e4dd521cSBarry Smith   rk->b1[1]=0.0;
171e4dd521cSBarry Smith   rk->b1[2]=1408.0/2565.0;
172e4dd521cSBarry Smith   rk->b1[3]=2197.0/4104.0;
173e4dd521cSBarry Smith   rk->b1[4]=-1.0/5.0;
174e4dd521cSBarry Smith   rk->b1[5]=0.0;
175e4dd521cSBarry Smith 
176e4dd521cSBarry Smith   rk->b2[0]=16.0/135.0;
177e4dd521cSBarry Smith   rk->b2[1]=0.0;
178e4dd521cSBarry Smith   rk->b2[2]=6656.0/12825.0;
179e4dd521cSBarry Smith   rk->b2[3]=28561.0/56430.0;
180e4dd521cSBarry Smith   rk->b2[4]=-9.0/50.0;
181e4dd521cSBarry Smith   rk->b2[5]=2.0/55.0;
182e4dd521cSBarry Smith   */
183e4dd521cSBarry Smith   /* Found in table on page 169: Merson 4("5") */
184e4dd521cSBarry Smith   /*
185e4dd521cSBarry Smith   rk->p=4;
186e4dd521cSBarry Smith   rk->s=5;
187e4dd521cSBarry Smith   rk->a[1][0] = 1.0/3.0;
188e4dd521cSBarry Smith   rk->a[2][0] = 1.0/6.0;
189e4dd521cSBarry Smith   rk->a[2][1] = 1.0/6.0;
190e4dd521cSBarry Smith   rk->a[3][0] = 1.0/8.0;
191e4dd521cSBarry Smith   rk->a[3][1] = 0.0;
192e4dd521cSBarry Smith   rk->a[3][2] = 3.0/8.0;
193e4dd521cSBarry Smith   rk->a[4][0] = 1.0/2.0;
194e4dd521cSBarry Smith   rk->a[4][1] = 0.0;
195e4dd521cSBarry Smith   rk->a[4][2] = -3.0/2.0;
196e4dd521cSBarry Smith   rk->a[4][3] = 2.0;
197e4dd521cSBarry Smith 
198e4dd521cSBarry Smith   rk->c[0] = 0.0;
199e4dd521cSBarry Smith   rk->c[1] = 1.0/3.0;
200e4dd521cSBarry Smith   rk->c[2] = 1.0/3.0;
201e4dd521cSBarry Smith   rk->c[3] = 0.5;
202e4dd521cSBarry Smith   rk->c[4] = 1.0;
203e4dd521cSBarry Smith 
204e4dd521cSBarry Smith   rk->b1[0] = 1.0/2.0;
205e4dd521cSBarry Smith   rk->b1[1] = 0.0;
206e4dd521cSBarry Smith   rk->b1[2] = -3.0/2.0;
207e4dd521cSBarry Smith   rk->b1[3] = 2.0;
208e4dd521cSBarry Smith   rk->b1[4] = 0.0;
209e4dd521cSBarry Smith 
210e4dd521cSBarry Smith   rk->b2[0] = 1.0/6.0;
211e4dd521cSBarry Smith   rk->b2[1] = 0.0;
212e4dd521cSBarry Smith   rk->b2[2] = 0.0;
213e4dd521cSBarry Smith   rk->b2[3] = 2.0/3.0;
214e4dd521cSBarry Smith   rk->b2[4] = 1.0/6.0;
215e4dd521cSBarry Smith   */
216e4dd521cSBarry Smith 
217e4dd521cSBarry Smith   /* making b2 -> e=b1-b2 */
218e4dd521cSBarry Smith   /*
219e4dd521cSBarry Smith     for(i=0;i<rk->s;i++){
22026edb414Sasbjorn      rk->b2[i] = (rk->b1[i]) - (rk->b2[i]);
221e4dd521cSBarry Smith   }
222e4dd521cSBarry Smith   */
22326edb414Sasbjorn   rk->b2[0]=71.0/57600.0;
22426edb414Sasbjorn   rk->b2[1]=0.0;
22526edb414Sasbjorn   rk->b2[2]=-71.0/16695.0;
22626edb414Sasbjorn   rk->b2[3]=71.0/1920.0;
22726edb414Sasbjorn   rk->b2[4]=-17253.0/339200.0;
22826edb414Sasbjorn   rk->b2[5]=22.0/525.0;
22926edb414Sasbjorn   rk->b2[6]=-1.0/40.0;
23026edb414Sasbjorn 
231e4dd521cSBarry Smith   /* initializing vectors */
232e4dd521cSBarry Smith   ierr = VecDuplicate(ts->vec_sol,&rk->y1);CHKERRQ(ierr);
233e4dd521cSBarry Smith   ierr = VecDuplicate(ts->vec_sol,&rk->y2);CHKERRQ(ierr);
234e4dd521cSBarry Smith   ierr = VecDuplicate(rk->y1,&rk->tmp);CHKERRQ(ierr);
235e4dd521cSBarry Smith   ierr = VecDuplicate(rk->y1,&rk->tmp_y);CHKERRQ(ierr);
236e4dd521cSBarry Smith   ierr = VecDuplicateVecs(rk->y1,rk->s,&rk->k);CHKERRQ(ierr);
237e4dd521cSBarry Smith 
238e4dd521cSBarry Smith   PetscFunctionReturn(0);
239e4dd521cSBarry Smith }
240e4dd521cSBarry Smith 
241db046809SBarry Smith /*------------------------------------------------------------*/
242db046809SBarry Smith #undef __FUNCT__
2435f70b478SJed Brown #define __FUNCT__ "TSRKqs"
2445f70b478SJed Brown PetscErrorCode TSRKqs(TS ts,PetscReal t,PetscReal h)
245db046809SBarry Smith {
2465f70b478SJed Brown   TS_RK          *rk = (TS_RK*)ts->data;
2476849ba73SBarry Smith   PetscErrorCode ierr;
248a7cc72afSBarry Smith   PetscInt       j,l;
249db046809SBarry Smith   PetscReal      tmp_t=t;
250b05257ddSBarry Smith   PetscScalar    hh=h;
251db046809SBarry Smith 
252db046809SBarry Smith   PetscFunctionBegin;
253db046809SBarry Smith   /* k[0]=0  */
254b05257ddSBarry Smith   ierr = VecSet(rk->k[0],0.0);CHKERRQ(ierr);
255db046809SBarry Smith 
256db046809SBarry Smith   /* k[0] = derivs(t,y1) */
257db046809SBarry Smith   ierr = TSComputeRHSFunction(ts,t,rk->y1,rk->k[0]);CHKERRQ(ierr);
258db046809SBarry Smith   /* looping over runge-kutta variables */
259db046809SBarry Smith   /* building the k - array of vectors */
260db046809SBarry Smith   for(j = 1 ; j < rk->s ; j++){
261db046809SBarry Smith 
262db046809SBarry Smith      /* rk->tmp = 0 */
263b05257ddSBarry Smith      ierr = VecSet(rk->tmp,0.0);CHKERRQ(ierr);
264db046809SBarry Smith 
265db046809SBarry Smith      for(l=0;l<j;l++){
266db046809SBarry Smith         /* tmp += a(j,l)*k[l] */
2672dcb1b2aSMatthew Knepley        ierr = VecAXPY(rk->tmp,rk->a[j][l],rk->k[l]);CHKERRQ(ierr);
268db046809SBarry Smith      }
269db046809SBarry Smith 
270db046809SBarry Smith      /* ierr = VecView(rk->tmp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
271db046809SBarry Smith 
272db046809SBarry Smith      /* k[j] = derivs(t+c(j)*h,y1+h*tmp,k(j)) */
273db046809SBarry Smith      /* I need the following helpers:
274db046809SBarry Smith         PetscScalar  tmp_t=t+c(j)*h
275db046809SBarry Smith         Vec          tmp_y=h*tmp+y1
276db046809SBarry Smith      */
277db046809SBarry Smith 
278db046809SBarry Smith      tmp_t = t + rk->c[j] * h;
279db046809SBarry Smith 
280db046809SBarry Smith      /* tmp_y = h * tmp + y1 */
2812dcb1b2aSMatthew Knepley      ierr = VecWAXPY(rk->tmp_y,hh,rk->tmp,rk->y1);CHKERRQ(ierr);
282db046809SBarry Smith 
283db046809SBarry Smith      /* rk->k[j]=0 */
284b05257ddSBarry Smith      ierr = VecSet(rk->k[j],0.0);CHKERRQ(ierr);
285db046809SBarry Smith      ierr = TSComputeRHSFunction(ts,tmp_t,rk->tmp_y,rk->k[j]);CHKERRQ(ierr);
286db046809SBarry Smith   }
287db046809SBarry Smith 
288db046809SBarry Smith   /* tmp=0 and tmp_y=0 */
289b05257ddSBarry Smith   ierr = VecSet(rk->tmp,0.0);CHKERRQ(ierr);
290b05257ddSBarry Smith   ierr = VecSet(rk->tmp_y,0.0);CHKERRQ(ierr);
291db046809SBarry Smith 
292db046809SBarry Smith   for(j = 0 ; j < rk->s ; j++){
293db046809SBarry Smith      /* tmp=b1[j]*k[j]+tmp  */
2942dcb1b2aSMatthew Knepley     ierr = VecAXPY(rk->tmp,rk->b1[j],rk->k[j]);CHKERRQ(ierr);
295db046809SBarry Smith      /* tmp_y=b2[j]*k[j]+tmp_y */
2962dcb1b2aSMatthew Knepley     ierr = VecAXPY(rk->tmp_y,rk->b2[j],rk->k[j]);CHKERRQ(ierr);
297db046809SBarry Smith   }
298db046809SBarry Smith 
299db046809SBarry Smith   /* y2 = hh * tmp_y */
300b05257ddSBarry Smith   ierr = VecSet(rk->y2,0.0);CHKERRQ(ierr);
3012dcb1b2aSMatthew Knepley   ierr = VecAXPY(rk->y2,hh,rk->tmp_y);CHKERRQ(ierr);
302db046809SBarry Smith   /* y1 = hh*tmp + y1 */
3032dcb1b2aSMatthew Knepley   ierr = VecAXPY(rk->y1,hh,rk->tmp);CHKERRQ(ierr);
304db046809SBarry Smith   /* Finding difference between y1 and y2 */
305db046809SBarry Smith   PetscFunctionReturn(0);
306db046809SBarry Smith }
307db046809SBarry Smith 
308e4dd521cSBarry Smith #undef __FUNCT__
3095f70b478SJed Brown #define __FUNCT__ "TSStep_RK"
3105f70b478SJed Brown static PetscErrorCode TSStep_RK(TS ts,PetscInt *steps,PetscReal *ptime)
311e4dd521cSBarry Smith {
3125f70b478SJed Brown   TS_RK          *rk = (TS_RK*)ts->data;
313b0dd9724SMatthew Knepley   PetscReal      norm=0.0,dt_fac=0.0,fac = 0.0/*,ttmp=0.0*/;
314*186e87acSLisandro Dalcin   PetscInt       i;
315*186e87acSLisandro Dalcin   PetscErrorCode ierr;
316e4dd521cSBarry Smith 
317e4dd521cSBarry Smith   PetscFunctionBegin;
318e4dd521cSBarry Smith   *steps = -ts->steps;
319*186e87acSLisandro Dalcin   *ptime  = ts->ptime;
320*186e87acSLisandro Dalcin 
321e4dd521cSBarry Smith   ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
322*186e87acSLisandro Dalcin 
323*186e87acSLisandro Dalcin   ierr = VecCopy(ts->vec_sol,rk->y1);CHKERRQ(ierr);
324*186e87acSLisandro Dalcin 
325e4dd521cSBarry Smith   /* while loop to get from start to stop */
326*186e87acSLisandro Dalcin   for (i = 0; i < ts->max_steps; i++) {
3273f2090d5SJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr); /* Note that this is called once per STEP, not once per STAGE. */
328*186e87acSLisandro Dalcin 
329e4dd521cSBarry Smith    /* calling rkqs */
330e4dd521cSBarry Smith      /*
331e4dd521cSBarry Smith        -- input
332e4dd521cSBarry Smith        ts        - pointer to ts
333e4dd521cSBarry Smith        ts->ptime - current time
334b05257ddSBarry Smith        ts->time_step        - try this timestep
335e4dd521cSBarry Smith        y1        - solution for this step
336e4dd521cSBarry Smith 
337e4dd521cSBarry Smith        --output
338e4dd521cSBarry Smith        y1        - suggested solution
339e4dd521cSBarry Smith        y2        - check solution (runge - kutta second permutation)
340e4dd521cSBarry Smith      */
3415f70b478SJed Brown      ierr = TSRKqs(ts,ts->ptime,ts->time_step);CHKERRQ(ierr);
342e3caeda6SBarry Smith      /* counting steps */
343e3caeda6SBarry Smith      ts->steps++;
344e4dd521cSBarry Smith    /* checking for maxerror */
345e4dd521cSBarry Smith      /* comparing difference to maxerror */
346e4dd521cSBarry Smith      ierr = VecNorm(rk->y2,NORM_2,&norm);CHKERRQ(ierr);
347e4dd521cSBarry Smith      /* modifying maxerror to satisfy this timestep */
348b05257ddSBarry Smith      rk->maxerror = rk->ferror * ts->time_step;
349b05257ddSBarry Smith      /* ierr = PetscPrintf(PETSC_COMM_WORLD,"norm err: %f maxerror: %f dt: %f",norm,rk->maxerror,ts->time_step);CHKERRQ(ierr); */
350e4dd521cSBarry Smith 
351e4dd521cSBarry Smith    /* handling ok and not ok */
352e4dd521cSBarry Smith      if (norm < rk->maxerror){
353e4dd521cSBarry Smith         /* if ok: */
354e4dd521cSBarry Smith         ierr=VecCopy(rk->y1,ts->vec_sol);CHKERRQ(ierr); /* saves the suggested solution to current solution */
355b05257ddSBarry Smith         ts->ptime += ts->time_step; /* storing the new current time */
356e4dd521cSBarry Smith         rk->nok++;
357e4dd521cSBarry Smith         fac=5.0;
358e4dd521cSBarry Smith         /* trying to save the vector */
3593f2090d5SJed Brown         ierr = TSPostStep(ts);CHKERRQ(ierr);
360e4dd521cSBarry Smith         ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
361e3caeda6SBarry Smith         if (ts->ptime >= ts->max_time) break;
362e4dd521cSBarry Smith      } else{
363e4dd521cSBarry Smith         /* if not OK */
364e4dd521cSBarry Smith         rk->nnok++;
365e4dd521cSBarry Smith         fac=1.0;
366e4dd521cSBarry Smith         ierr=VecCopy(ts->vec_sol,rk->y1);CHKERRQ(ierr);  /* restores old solution */
367e4dd521cSBarry Smith      }
368e4dd521cSBarry Smith 
369e4dd521cSBarry Smith      /*Computing next stepsize. See page 167 in Solving ODE 1
370e4dd521cSBarry Smith       *
371e4dd521cSBarry Smith       * h_new = h * min( facmax , max( facmin , fac * (tol/err)^(1/(p+1)) ) )
372e4dd521cSBarry Smith       * facmax set above
373e4dd521cSBarry Smith       * facmin
374e4dd521cSBarry Smith       */
375e4dd521cSBarry Smith      dt_fac = exp(log((rk->maxerror) / norm) / ((rk->p) + 1) ) * 0.9 ;
376e4dd521cSBarry Smith 
377e4dd521cSBarry Smith      if (dt_fac > fac){
3782cdf8120Sasbjorn         /*ierr = PetscPrintf(PETSC_COMM_WORLD,"changing fac %f\n",fac);*/
379e4dd521cSBarry Smith         dt_fac = fac;
380e4dd521cSBarry Smith      }
381e4dd521cSBarry Smith 
382b05257ddSBarry Smith      /* computing new ts->time_step */
383b05257ddSBarry Smith      ts->time_step = ts->time_step * dt_fac;
384e4dd521cSBarry Smith 
385b05257ddSBarry Smith      if (ts->ptime+ts->time_step > ts->max_time){
386b05257ddSBarry Smith         ts->time_step = ts->max_time - ts->ptime;
387e4dd521cSBarry Smith      }
388e4dd521cSBarry Smith 
389b05257ddSBarry Smith      if (ts->time_step < 1e-14){
390b05257ddSBarry Smith         ierr = PetscPrintf(PETSC_COMM_WORLD,"Very small steps: %f\n",ts->time_step);CHKERRQ(ierr);
391b05257ddSBarry Smith         ts->time_step = 1e-14;
392e4dd521cSBarry Smith      }
393e4dd521cSBarry Smith 
394e4dd521cSBarry Smith      /* trying to purify h */
395e4dd521cSBarry Smith      /* (did not give any visible result) */
396b05257ddSBarry Smith      /* ttmp = ts->ptime + ts->time_step;
397b05257ddSBarry Smith         ts->time_step = ttmp - ts->ptime; */
398e4dd521cSBarry Smith 
399e4dd521cSBarry Smith   }
400e4dd521cSBarry Smith 
401e4dd521cSBarry Smith   ierr=VecCopy(rk->y1,ts->vec_sol);CHKERRQ(ierr);
402*186e87acSLisandro Dalcin 
403e4dd521cSBarry Smith   *steps += ts->steps;
404e4dd521cSBarry Smith   *ptime  = ts->ptime;
405e4dd521cSBarry Smith   PetscFunctionReturn(0);
406e4dd521cSBarry Smith }
407e4dd521cSBarry Smith 
408e4dd521cSBarry Smith /*------------------------------------------------------------*/
409e4dd521cSBarry Smith #undef __FUNCT__
410277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_RK"
411277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts)
412e4dd521cSBarry Smith {
4135f70b478SJed Brown   TS_RK          *rk = (TS_RK*)ts->data;
4146849ba73SBarry Smith   PetscErrorCode ierr;
415e4dd521cSBarry Smith 
416e4dd521cSBarry Smith   PetscFunctionBegin;
417e4dd521cSBarry Smith   if (rk->y1)    {ierr = VecDestroy(rk->y1);CHKERRQ(ierr);}
418e4dd521cSBarry Smith   if (rk->y2)    {ierr = VecDestroy(rk->y2);CHKERRQ(ierr);}
419e4dd521cSBarry Smith   if (rk->tmp)   {ierr = VecDestroy(rk->tmp);CHKERRQ(ierr);}
420e4dd521cSBarry Smith   if (rk->tmp_y) {ierr = VecDestroy(rk->tmp_y);CHKERRQ(ierr);}
421277b19d0SLisandro Dalcin   if (rk->k)     {ierr = VecDestroyVecs(rk->s,&rk->k);CHKERRQ(ierr);}
422277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
423e4dd521cSBarry Smith }
424277b19d0SLisandro Dalcin 
425277b19d0SLisandro Dalcin #undef __FUNCT__
426277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_RK"
427277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_RK(TS ts)
428277b19d0SLisandro Dalcin {
429277b19d0SLisandro Dalcin   PetscErrorCode ierr;
430277b19d0SLisandro Dalcin 
431277b19d0SLisandro Dalcin   PetscFunctionBegin;
432277b19d0SLisandro Dalcin   ierr = TSReset_RK(ts);CHKERRQ(ierr);
433277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
434e4dd521cSBarry Smith   PetscFunctionReturn(0);
435e4dd521cSBarry Smith }
436e4dd521cSBarry Smith /*------------------------------------------------------------*/
437e4dd521cSBarry Smith 
438e4dd521cSBarry Smith #undef __FUNCT__
4395f70b478SJed Brown #define __FUNCT__ "TSSetFromOptions_RK"
4405f70b478SJed Brown static PetscErrorCode TSSetFromOptions_RK(TS ts)
441e4dd521cSBarry Smith {
4425f70b478SJed Brown   TS_RK          *rk = (TS_RK*)ts->data;
443dfbe8321SBarry Smith   PetscErrorCode ierr;
444262ff7bbSBarry Smith 
445e4dd521cSBarry Smith   PetscFunctionBegin;
446262ff7bbSBarry Smith   ierr = PetscOptionsHead("RK ODE solver options");CHKERRQ(ierr);
447262ff7bbSBarry Smith     ierr = PetscOptionsReal("-ts_rk_tol","Tolerance for convergence","TSRKSetTolerance",rk->tolerance,&rk->tolerance,PETSC_NULL);CHKERRQ(ierr);
448262ff7bbSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
449e4dd521cSBarry Smith   PetscFunctionReturn(0);
450e4dd521cSBarry Smith }
451e4dd521cSBarry Smith 
452e4dd521cSBarry Smith #undef __FUNCT__
4535f70b478SJed Brown #define __FUNCT__ "TSView_RK"
4545f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
455e4dd521cSBarry Smith {
4565f70b478SJed Brown    TS_RK          *rk = (TS_RK*)ts->data;
457dfbe8321SBarry Smith    PetscErrorCode ierr;
4582cdf8120Sasbjorn 
459e4dd521cSBarry Smith    PetscFunctionBegin;
46077431f27SBarry Smith    ierr = PetscPrintf(PETSC_COMM_WORLD,"  number of ok steps: %D\n",rk->nok);CHKERRQ(ierr);
46177431f27SBarry Smith    ierr = PetscPrintf(PETSC_COMM_WORLD,"  number of rejected steps: %D\n",rk->nnok);CHKERRQ(ierr);
462e4dd521cSBarry Smith    PetscFunctionReturn(0);
463e4dd521cSBarry Smith }
464e4dd521cSBarry Smith 
465e4dd521cSBarry Smith /* ------------------------------------------------------------ */
466ebe3b25bSBarry Smith /*MC
46796f5712cSJed Brown       TSRK - ODE solver using the explicit Runge-Kutta methods
468ebe3b25bSBarry Smith 
469ebe3b25bSBarry Smith    Options Database:
470ebe3b25bSBarry Smith .  -ts_rk_tol <tol> Tolerance for convergence
471ebe3b25bSBarry Smith 
472ebe3b25bSBarry Smith   Contributed by: Asbjorn Hoiland Aarrestad, asbjorn@aarrestad.com, http://asbjorn.aarrestad.com/
473ebe3b25bSBarry Smith 
474d41222bbSBarry Smith   Level: beginner
475d41222bbSBarry Smith 
4769596e0b4SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSRKSetTolerance()
477ebe3b25bSBarry Smith 
478ebe3b25bSBarry Smith M*/
479ebe3b25bSBarry Smith 
480e4dd521cSBarry Smith EXTERN_C_BEGIN
481e4dd521cSBarry Smith #undef __FUNCT__
4825f70b478SJed Brown #define __FUNCT__ "TSCreate_RK"
4837087cfbeSBarry Smith PetscErrorCode  TSCreate_RK(TS ts)
484e4dd521cSBarry Smith {
4855f70b478SJed Brown   TS_RK          *rk;
486dfbe8321SBarry Smith   PetscErrorCode ierr;
487e4dd521cSBarry Smith 
488e4dd521cSBarry Smith   PetscFunctionBegin;
4895f70b478SJed Brown   ts->ops->setup           = TSSetUp_RK;
4905f70b478SJed Brown   ts->ops->step            = TSStep_RK;
4915f70b478SJed Brown   ts->ops->destroy         = TSDestroy_RK;
4925f70b478SJed Brown   ts->ops->setfromoptions  = TSSetFromOptions_RK;
4935f70b478SJed Brown   ts->ops->view            = TSView_RK;
494e4dd521cSBarry Smith 
4955f70b478SJed Brown   ierr = PetscNewLog(ts,TS_RK,&rk);CHKERRQ(ierr);
496e4dd521cSBarry Smith   ts->data = (void*)rk;
497e4dd521cSBarry Smith 
498a7cc72afSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRKSetTolerance_C","TSRKSetTolerance_RK",TSRKSetTolerance_RK);CHKERRQ(ierr);
499262ff7bbSBarry Smith 
500e4dd521cSBarry Smith   PetscFunctionReturn(0);
501e4dd521cSBarry Smith }
502e4dd521cSBarry Smith EXTERN_C_END
503