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