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