163dd3a1aSKris Buschelman #define PETSCTS_DLL 263dd3a1aSKris Buschelman 3e4dd521cSBarry Smith /* 426edb414Sasbjorn * Code for Timestepping with Runge Kutta 526edb414Sasbjorn * 626edb414Sasbjorn * Written by 726edb414Sasbjorn * Asbjorn Hoiland Aarrestad 826edb414Sasbjorn * asbjorn@aarrestad.com 926edb414Sasbjorn * http://asbjorn.aarrestad.com/ 1026edb414Sasbjorn * 11e4dd521cSBarry Smith */ 12*7c4f633dSBarry Smith #include "private/tsimpl.h" /*I "petscts.h" I*/ 132cdf8120Sasbjorn #include "time.h" 14e4dd521cSBarry Smith 15e4dd521cSBarry Smith typedef struct { 16e4dd521cSBarry Smith Vec y1,y2; /* work wectors for the two rk permuations */ 17a7cc72afSBarry Smith PetscInt nok,nnok; /* counters for ok and not ok steps */ 18e4dd521cSBarry Smith PetscReal maxerror; /* variable to tell the maxerror allowed */ 19e4dd521cSBarry Smith PetscReal ferror; /* variable to tell (global maxerror)/(total time) */ 20262ff7bbSBarry Smith PetscReal tolerance; /* initial value set for maxerror by user */ 21e4dd521cSBarry Smith Vec tmp,tmp_y,*k; /* two temp vectors and the k vectors for rk */ 22e4dd521cSBarry Smith PetscScalar a[7][6]; /* rk scalars */ 23e4dd521cSBarry Smith PetscScalar b1[7],b2[7]; /* rk scalars */ 24e4dd521cSBarry Smith PetscReal c[7]; /* rk scalars */ 25a7cc72afSBarry Smith PetscInt p,s; /* variables to tell the size of the runge-kutta solver */ 26e4dd521cSBarry Smith } TS_Rk; 27e4dd521cSBarry Smith 28262ff7bbSBarry Smith EXTERN_C_BEGIN 29262ff7bbSBarry Smith #undef __FUNCT__ 30262ff7bbSBarry Smith #define __FUNCT__ "TSRKSetTolerance_RK" 3163dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSRKSetTolerance_RK(TS ts,PetscReal aabs) 32262ff7bbSBarry Smith { 33262ff7bbSBarry Smith TS_Rk *rk = (TS_Rk*)ts->data; 34262ff7bbSBarry Smith 35262ff7bbSBarry Smith PetscFunctionBegin; 36262ff7bbSBarry Smith rk->tolerance = aabs; 37262ff7bbSBarry Smith PetscFunctionReturn(0); 38262ff7bbSBarry Smith } 39262ff7bbSBarry Smith EXTERN_C_END 40262ff7bbSBarry Smith 41262ff7bbSBarry Smith #undef __FUNCT__ 42262ff7bbSBarry Smith #define __FUNCT__ "TSRKSetTolerance" 43262ff7bbSBarry Smith /*@ 44262ff7bbSBarry Smith TSRKSetTolerance - Sets the total error the RK explicit time integrators 45262ff7bbSBarry Smith will allow over the given time interval. 46262ff7bbSBarry Smith 47262ff7bbSBarry Smith Collective on TS 48262ff7bbSBarry Smith 49262ff7bbSBarry Smith Input parameters: 50262ff7bbSBarry Smith + ts - the time-step context 51262ff7bbSBarry Smith - aabs - the absolute tolerance 52262ff7bbSBarry Smith 53262ff7bbSBarry Smith Level: intermediate 54262ff7bbSBarry Smith 55262ff7bbSBarry Smith .keywords: RK, tolerance 56262ff7bbSBarry Smith 570f3b3ca1SHong Zhang .seealso: TSSundialsSetTolerance() 58262ff7bbSBarry Smith 59262ff7bbSBarry Smith @*/ 6063dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSRKSetTolerance(TS ts,PetscReal aabs) 61262ff7bbSBarry Smith { 62dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(TS,PetscReal); 63262ff7bbSBarry Smith 64262ff7bbSBarry Smith PetscFunctionBegin; 65262ff7bbSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)ts,"TSRKSetTolerance_C",(void (**)(void))&f);CHKERRQ(ierr); 66262ff7bbSBarry Smith if (f) { 67262ff7bbSBarry Smith ierr = (*f)(ts,aabs);CHKERRQ(ierr); 68262ff7bbSBarry Smith } 69262ff7bbSBarry Smith PetscFunctionReturn(0); 70262ff7bbSBarry Smith } 71e4dd521cSBarry Smith 72e4dd521cSBarry Smith 73e4dd521cSBarry Smith #undef __FUNCT__ 74e4dd521cSBarry Smith #define __FUNCT__ "TSSetUp_Rk" 756849ba73SBarry Smith static PetscErrorCode TSSetUp_Rk(TS ts) 76e4dd521cSBarry Smith { 77e4dd521cSBarry Smith TS_Rk *rk = (TS_Rk*)ts->data; 78a7cc72afSBarry Smith PetscErrorCode ierr; 79e4dd521cSBarry Smith 80e4dd521cSBarry Smith PetscFunctionBegin; 81e4dd521cSBarry Smith rk->nok = 0; 82e4dd521cSBarry Smith rk->nnok = 0; 83262ff7bbSBarry Smith rk->maxerror = rk->tolerance; 84e4dd521cSBarry Smith 85e4dd521cSBarry Smith /* fixing maxerror: global vs local */ 86e4dd521cSBarry Smith rk->ferror = rk->maxerror / (ts->max_time - ts->ptime); 87e4dd521cSBarry Smith 88e4dd521cSBarry Smith /* 34.0/45.0 gives double precision division */ 89e4dd521cSBarry Smith /* defining variables needed for Runge-Kutta computing */ 90e4dd521cSBarry Smith /* when changing below, please remember to change a, b1, b2 and c above! */ 91e4dd521cSBarry Smith /* Found in table on page 171: Dormand-Prince 5(4) */ 92e4dd521cSBarry Smith 93e4dd521cSBarry Smith /* are these right? */ 94e4dd521cSBarry Smith rk->p=6; 95e4dd521cSBarry Smith rk->s=7; 96e4dd521cSBarry Smith 97e4dd521cSBarry Smith rk->a[1][0]=1.0/5.0; 98e4dd521cSBarry Smith rk->a[2][0]=3.0/40.0; 99e4dd521cSBarry Smith rk->a[2][1]=9.0/40.0; 100e4dd521cSBarry Smith rk->a[3][0]=44.0/45.0; 101e4dd521cSBarry Smith rk->a[3][1]=-56.0/15.0; 102e4dd521cSBarry Smith rk->a[3][2]=32.0/9.0; 103e4dd521cSBarry Smith rk->a[4][0]=19372.0/6561.0; 104e4dd521cSBarry Smith rk->a[4][1]=-25360.0/2187.0; 105e4dd521cSBarry Smith rk->a[4][2]=64448.0/6561.0; 106e4dd521cSBarry Smith rk->a[4][3]=-212.0/729.0; 107e4dd521cSBarry Smith rk->a[5][0]=9017.0/3168.0; 108e4dd521cSBarry Smith rk->a[5][1]=-355.0/33.0; 109e4dd521cSBarry Smith rk->a[5][2]=46732.0/5247.0; 110e4dd521cSBarry Smith rk->a[5][3]=49.0/176.0; 111e4dd521cSBarry Smith rk->a[5][4]=-5103.0/18656.0; 112e4dd521cSBarry Smith rk->a[6][0]=35.0/384.0; 113e4dd521cSBarry Smith rk->a[6][1]=0.0; 114e4dd521cSBarry Smith rk->a[6][2]=500.0/1113.0; 115e4dd521cSBarry Smith rk->a[6][3]=125.0/192.0; 116e4dd521cSBarry Smith rk->a[6][4]=-2187.0/6784.0; 117e4dd521cSBarry Smith rk->a[6][5]=11.0/84.0; 118e4dd521cSBarry Smith 119e4dd521cSBarry Smith 120e4dd521cSBarry Smith rk->c[0]=0.0; 121e4dd521cSBarry Smith rk->c[1]=1.0/5.0; 122e4dd521cSBarry Smith rk->c[2]=3.0/10.0; 123e4dd521cSBarry Smith rk->c[3]=4.0/5.0; 124e4dd521cSBarry Smith rk->c[4]=8.0/9.0; 125e4dd521cSBarry Smith rk->c[5]=1.0; 126e4dd521cSBarry Smith rk->c[6]=1.0; 127e4dd521cSBarry Smith 128e4dd521cSBarry Smith rk->b1[0]=35.0/384.0; 129e4dd521cSBarry Smith rk->b1[1]=0.0; 130e4dd521cSBarry Smith rk->b1[2]=500.0/1113.0; 131e4dd521cSBarry Smith rk->b1[3]=125.0/192.0; 132e4dd521cSBarry Smith rk->b1[4]=-2187.0/6784.0; 133e4dd521cSBarry Smith rk->b1[5]=11.0/84.0; 134e4dd521cSBarry Smith rk->b1[6]=0.0; 135e4dd521cSBarry Smith 136e4dd521cSBarry Smith rk->b2[0]=5179.0/57600.0; 137e4dd521cSBarry Smith rk->b2[1]=0.0; 138e4dd521cSBarry Smith rk->b2[2]=7571.0/16695.0; 139e4dd521cSBarry Smith rk->b2[3]=393.0/640.0; 140e4dd521cSBarry Smith rk->b2[4]=-92097.0/339200.0; 141e4dd521cSBarry Smith rk->b2[5]=187.0/2100.0; 142e4dd521cSBarry Smith rk->b2[6]=1.0/40.0; 143e4dd521cSBarry Smith 144e4dd521cSBarry Smith 145e4dd521cSBarry Smith /* Found in table on page 170: Fehlberg 4(5) */ 146e4dd521cSBarry Smith /* 147e4dd521cSBarry Smith rk->p=5; 148e4dd521cSBarry Smith rk->s=6; 149e4dd521cSBarry Smith 150e4dd521cSBarry Smith rk->a[1][0]=1.0/4.0; 151e4dd521cSBarry Smith rk->a[2][0]=3.0/32.0; 152e4dd521cSBarry Smith rk->a[2][1]=9.0/32.0; 153e4dd521cSBarry Smith rk->a[3][0]=1932.0/2197.0; 154e4dd521cSBarry Smith rk->a[3][1]=-7200.0/2197.0; 155e4dd521cSBarry Smith rk->a[3][2]=7296.0/2197.0; 156e4dd521cSBarry Smith rk->a[4][0]=439.0/216.0; 157e4dd521cSBarry Smith rk->a[4][1]=-8.0; 158e4dd521cSBarry Smith rk->a[4][2]=3680.0/513.0; 159e4dd521cSBarry Smith rk->a[4][3]=-845.0/4104.0; 160e4dd521cSBarry Smith rk->a[5][0]=-8.0/27.0; 161e4dd521cSBarry Smith rk->a[5][1]=2.0; 162e4dd521cSBarry Smith rk->a[5][2]=-3544.0/2565.0; 163e4dd521cSBarry Smith rk->a[5][3]=1859.0/4104.0; 164e4dd521cSBarry Smith rk->a[5][4]=-11.0/40.0; 165e4dd521cSBarry Smith 166e4dd521cSBarry Smith rk->c[0]=0.0; 167e4dd521cSBarry Smith rk->c[1]=1.0/4.0; 168e4dd521cSBarry Smith rk->c[2]=3.0/8.0; 169e4dd521cSBarry Smith rk->c[3]=12.0/13.0; 170e4dd521cSBarry Smith rk->c[4]=1.0; 171e4dd521cSBarry Smith rk->c[5]=1.0/2.0; 172e4dd521cSBarry Smith 173e4dd521cSBarry Smith rk->b1[0]=25.0/216.0; 174e4dd521cSBarry Smith rk->b1[1]=0.0; 175e4dd521cSBarry Smith rk->b1[2]=1408.0/2565.0; 176e4dd521cSBarry Smith rk->b1[3]=2197.0/4104.0; 177e4dd521cSBarry Smith rk->b1[4]=-1.0/5.0; 178e4dd521cSBarry Smith rk->b1[5]=0.0; 179e4dd521cSBarry Smith 180e4dd521cSBarry Smith rk->b2[0]=16.0/135.0; 181e4dd521cSBarry Smith rk->b2[1]=0.0; 182e4dd521cSBarry Smith rk->b2[2]=6656.0/12825.0; 183e4dd521cSBarry Smith rk->b2[3]=28561.0/56430.0; 184e4dd521cSBarry Smith rk->b2[4]=-9.0/50.0; 185e4dd521cSBarry Smith rk->b2[5]=2.0/55.0; 186e4dd521cSBarry Smith */ 187e4dd521cSBarry Smith /* Found in table on page 169: Merson 4("5") */ 188e4dd521cSBarry Smith /* 189e4dd521cSBarry Smith rk->p=4; 190e4dd521cSBarry Smith rk->s=5; 191e4dd521cSBarry Smith rk->a[1][0] = 1.0/3.0; 192e4dd521cSBarry Smith rk->a[2][0] = 1.0/6.0; 193e4dd521cSBarry Smith rk->a[2][1] = 1.0/6.0; 194e4dd521cSBarry Smith rk->a[3][0] = 1.0/8.0; 195e4dd521cSBarry Smith rk->a[3][1] = 0.0; 196e4dd521cSBarry Smith rk->a[3][2] = 3.0/8.0; 197e4dd521cSBarry Smith rk->a[4][0] = 1.0/2.0; 198e4dd521cSBarry Smith rk->a[4][1] = 0.0; 199e4dd521cSBarry Smith rk->a[4][2] = -3.0/2.0; 200e4dd521cSBarry Smith rk->a[4][3] = 2.0; 201e4dd521cSBarry Smith 202e4dd521cSBarry Smith rk->c[0] = 0.0; 203e4dd521cSBarry Smith rk->c[1] = 1.0/3.0; 204e4dd521cSBarry Smith rk->c[2] = 1.0/3.0; 205e4dd521cSBarry Smith rk->c[3] = 0.5; 206e4dd521cSBarry Smith rk->c[4] = 1.0; 207e4dd521cSBarry Smith 208e4dd521cSBarry Smith rk->b1[0] = 1.0/2.0; 209e4dd521cSBarry Smith rk->b1[1] = 0.0; 210e4dd521cSBarry Smith rk->b1[2] = -3.0/2.0; 211e4dd521cSBarry Smith rk->b1[3] = 2.0; 212e4dd521cSBarry Smith rk->b1[4] = 0.0; 213e4dd521cSBarry Smith 214e4dd521cSBarry Smith rk->b2[0] = 1.0/6.0; 215e4dd521cSBarry Smith rk->b2[1] = 0.0; 216e4dd521cSBarry Smith rk->b2[2] = 0.0; 217e4dd521cSBarry Smith rk->b2[3] = 2.0/3.0; 218e4dd521cSBarry Smith rk->b2[4] = 1.0/6.0; 219e4dd521cSBarry Smith */ 220e4dd521cSBarry Smith 221e4dd521cSBarry Smith /* making b2 -> e=b1-b2 */ 222e4dd521cSBarry Smith /* 223e4dd521cSBarry Smith for(i=0;i<rk->s;i++){ 22426edb414Sasbjorn rk->b2[i] = (rk->b1[i]) - (rk->b2[i]); 225e4dd521cSBarry Smith } 226e4dd521cSBarry Smith */ 22726edb414Sasbjorn rk->b2[0]=71.0/57600.0; 22826edb414Sasbjorn rk->b2[1]=0.0; 22926edb414Sasbjorn rk->b2[2]=-71.0/16695.0; 23026edb414Sasbjorn rk->b2[3]=71.0/1920.0; 23126edb414Sasbjorn rk->b2[4]=-17253.0/339200.0; 23226edb414Sasbjorn rk->b2[5]=22.0/525.0; 23326edb414Sasbjorn rk->b2[6]=-1.0/40.0; 23426edb414Sasbjorn 235e4dd521cSBarry Smith /* initializing vectors */ 236e4dd521cSBarry Smith ierr = VecDuplicate(ts->vec_sol,&rk->y1);CHKERRQ(ierr); 237e4dd521cSBarry Smith ierr = VecDuplicate(ts->vec_sol,&rk->y2);CHKERRQ(ierr); 238e4dd521cSBarry Smith ierr = VecDuplicate(rk->y1,&rk->tmp);CHKERRQ(ierr); 239e4dd521cSBarry Smith ierr = VecDuplicate(rk->y1,&rk->tmp_y);CHKERRQ(ierr); 240e4dd521cSBarry Smith ierr = VecDuplicateVecs(rk->y1,rk->s,&rk->k);CHKERRQ(ierr); 241e4dd521cSBarry Smith 242e4dd521cSBarry Smith PetscFunctionReturn(0); 243e4dd521cSBarry Smith } 244e4dd521cSBarry Smith 245db046809SBarry Smith /*------------------------------------------------------------*/ 246db046809SBarry Smith #undef __FUNCT__ 247db046809SBarry Smith #define __FUNCT__ "TSRkqs" 248dfbe8321SBarry Smith PetscErrorCode TSRkqs(TS ts,PetscReal t,PetscReal h) 249db046809SBarry Smith { 250db046809SBarry Smith TS_Rk *rk = (TS_Rk*)ts->data; 2516849ba73SBarry Smith PetscErrorCode ierr; 252a7cc72afSBarry Smith PetscInt j,l; 253db046809SBarry Smith PetscReal tmp_t=t; 254b05257ddSBarry Smith PetscScalar hh=h; 255db046809SBarry Smith 256db046809SBarry Smith PetscFunctionBegin; 257db046809SBarry Smith /* k[0]=0 */ 258b05257ddSBarry Smith ierr = VecSet(rk->k[0],0.0);CHKERRQ(ierr); 259db046809SBarry Smith 260db046809SBarry Smith /* k[0] = derivs(t,y1) */ 261db046809SBarry Smith ierr = TSComputeRHSFunction(ts,t,rk->y1,rk->k[0]);CHKERRQ(ierr); 262db046809SBarry Smith /* looping over runge-kutta variables */ 263db046809SBarry Smith /* building the k - array of vectors */ 264db046809SBarry Smith for(j = 1 ; j < rk->s ; j++){ 265db046809SBarry Smith 266db046809SBarry Smith /* rk->tmp = 0 */ 267b05257ddSBarry Smith ierr = VecSet(rk->tmp,0.0);CHKERRQ(ierr); 268db046809SBarry Smith 269db046809SBarry Smith for(l=0;l<j;l++){ 270db046809SBarry Smith /* tmp += a(j,l)*k[l] */ 2712dcb1b2aSMatthew Knepley ierr = VecAXPY(rk->tmp,rk->a[j][l],rk->k[l]);CHKERRQ(ierr); 272db046809SBarry Smith } 273db046809SBarry Smith 274db046809SBarry Smith /* ierr = VecView(rk->tmp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 275db046809SBarry Smith 276db046809SBarry Smith /* k[j] = derivs(t+c(j)*h,y1+h*tmp,k(j)) */ 277db046809SBarry Smith /* I need the following helpers: 278db046809SBarry Smith PetscScalar tmp_t=t+c(j)*h 279db046809SBarry Smith Vec tmp_y=h*tmp+y1 280db046809SBarry Smith */ 281db046809SBarry Smith 282db046809SBarry Smith tmp_t = t + rk->c[j] * h; 283db046809SBarry Smith 284db046809SBarry Smith /* tmp_y = h * tmp + y1 */ 2852dcb1b2aSMatthew Knepley ierr = VecWAXPY(rk->tmp_y,hh,rk->tmp,rk->y1);CHKERRQ(ierr); 286db046809SBarry Smith 287db046809SBarry Smith /* rk->k[j]=0 */ 288b05257ddSBarry Smith ierr = VecSet(rk->k[j],0.0);CHKERRQ(ierr); 289db046809SBarry Smith ierr = TSComputeRHSFunction(ts,tmp_t,rk->tmp_y,rk->k[j]);CHKERRQ(ierr); 290db046809SBarry Smith } 291db046809SBarry Smith 292db046809SBarry Smith /* tmp=0 and tmp_y=0 */ 293b05257ddSBarry Smith ierr = VecSet(rk->tmp,0.0);CHKERRQ(ierr); 294b05257ddSBarry Smith ierr = VecSet(rk->tmp_y,0.0);CHKERRQ(ierr); 295db046809SBarry Smith 296db046809SBarry Smith for(j = 0 ; j < rk->s ; j++){ 297db046809SBarry Smith /* tmp=b1[j]*k[j]+tmp */ 2982dcb1b2aSMatthew Knepley ierr = VecAXPY(rk->tmp,rk->b1[j],rk->k[j]);CHKERRQ(ierr); 299db046809SBarry Smith /* tmp_y=b2[j]*k[j]+tmp_y */ 3002dcb1b2aSMatthew Knepley ierr = VecAXPY(rk->tmp_y,rk->b2[j],rk->k[j]);CHKERRQ(ierr); 301db046809SBarry Smith } 302db046809SBarry Smith 303db046809SBarry Smith /* y2 = hh * tmp_y */ 304b05257ddSBarry Smith ierr = VecSet(rk->y2,0.0);CHKERRQ(ierr); 3052dcb1b2aSMatthew Knepley ierr = VecAXPY(rk->y2,hh,rk->tmp_y);CHKERRQ(ierr); 306db046809SBarry Smith /* y1 = hh*tmp + y1 */ 3072dcb1b2aSMatthew Knepley ierr = VecAXPY(rk->y1,hh,rk->tmp);CHKERRQ(ierr); 308db046809SBarry Smith /* Finding difference between y1 and y2 */ 309db046809SBarry Smith PetscFunctionReturn(0); 310db046809SBarry Smith } 311db046809SBarry Smith 312e4dd521cSBarry Smith #undef __FUNCT__ 313e4dd521cSBarry Smith #define __FUNCT__ "TSStep_Rk" 314a7cc72afSBarry Smith static PetscErrorCode TSStep_Rk(TS ts,PetscInt *steps,PetscReal *ptime) 315e4dd521cSBarry Smith { 316e4dd521cSBarry Smith TS_Rk *rk = (TS_Rk*)ts->data; 317a7cc72afSBarry Smith PetscErrorCode ierr; 318b0dd9724SMatthew Knepley PetscReal norm=0.0,dt_fac=0.0,fac = 0.0/*,ttmp=0.0*/; 319e4dd521cSBarry Smith 320e4dd521cSBarry Smith PetscFunctionBegin; 321e4dd521cSBarry Smith ierr=VecCopy(ts->vec_sol,rk->y1);CHKERRQ(ierr); 322e4dd521cSBarry Smith *steps = -ts->steps; 323e4dd521cSBarry Smith /* trying to save the vector */ 324e4dd521cSBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 325e4dd521cSBarry Smith /* while loop to get from start to stop */ 326e4dd521cSBarry Smith while (ts->ptime < ts->max_time){ 327e4dd521cSBarry Smith /* calling rkqs */ 328e4dd521cSBarry Smith /* 329e4dd521cSBarry Smith -- input 330e4dd521cSBarry Smith ts - pointer to ts 331e4dd521cSBarry Smith ts->ptime - current time 332b05257ddSBarry Smith ts->time_step - try this timestep 333e4dd521cSBarry Smith y1 - solution for this step 334e4dd521cSBarry Smith 335e4dd521cSBarry Smith --output 336e4dd521cSBarry Smith y1 - suggested solution 337e4dd521cSBarry Smith y2 - check solution (runge - kutta second permutation) 338e4dd521cSBarry Smith */ 339b05257ddSBarry Smith ierr = TSRkqs(ts,ts->ptime,ts->time_step);CHKERRQ(ierr); 340e4dd521cSBarry Smith /* checking for maxerror */ 341e4dd521cSBarry Smith /* comparing difference to maxerror */ 342e4dd521cSBarry Smith ierr = VecNorm(rk->y2,NORM_2,&norm);CHKERRQ(ierr); 343e4dd521cSBarry Smith /* modifying maxerror to satisfy this timestep */ 344b05257ddSBarry Smith rk->maxerror = rk->ferror * ts->time_step; 345b05257ddSBarry Smith /* ierr = PetscPrintf(PETSC_COMM_WORLD,"norm err: %f maxerror: %f dt: %f",norm,rk->maxerror,ts->time_step);CHKERRQ(ierr); */ 346e4dd521cSBarry Smith 347e4dd521cSBarry Smith /* handling ok and not ok */ 348e4dd521cSBarry Smith if (norm < rk->maxerror){ 349e4dd521cSBarry Smith /* if ok: */ 350e4dd521cSBarry Smith ierr=VecCopy(rk->y1,ts->vec_sol);CHKERRQ(ierr); /* saves the suggested solution to current solution */ 351b05257ddSBarry Smith ts->ptime += ts->time_step; /* storing the new current time */ 352e4dd521cSBarry Smith rk->nok++; 353e4dd521cSBarry Smith fac=5.0; 354e4dd521cSBarry Smith /* trying to save the vector */ 355e4dd521cSBarry Smith /* calling monitor */ 356e4dd521cSBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 357e4dd521cSBarry Smith } else{ 358e4dd521cSBarry Smith /* if not OK */ 359e4dd521cSBarry Smith rk->nnok++; 360e4dd521cSBarry Smith fac=1.0; 361e4dd521cSBarry Smith ierr=VecCopy(ts->vec_sol,rk->y1);CHKERRQ(ierr); /* restores old solution */ 362e4dd521cSBarry Smith } 363e4dd521cSBarry Smith 364e4dd521cSBarry Smith /*Computing next stepsize. See page 167 in Solving ODE 1 365e4dd521cSBarry Smith * 366e4dd521cSBarry Smith * h_new = h * min( facmax , max( facmin , fac * (tol/err)^(1/(p+1)) ) ) 367e4dd521cSBarry Smith * facmax set above 368e4dd521cSBarry Smith * facmin 369e4dd521cSBarry Smith */ 370e4dd521cSBarry Smith dt_fac = exp(log((rk->maxerror) / norm) / ((rk->p) + 1) ) * 0.9 ; 371e4dd521cSBarry Smith 372e4dd521cSBarry Smith if (dt_fac > fac){ 3732cdf8120Sasbjorn /*ierr = PetscPrintf(PETSC_COMM_WORLD,"changing fac %f\n",fac);*/ 374e4dd521cSBarry Smith dt_fac = fac; 375e4dd521cSBarry Smith } 376e4dd521cSBarry Smith 377b05257ddSBarry Smith /* computing new ts->time_step */ 378b05257ddSBarry Smith ts->time_step = ts->time_step * dt_fac; 379e4dd521cSBarry Smith 380b05257ddSBarry Smith if (ts->ptime+ts->time_step > ts->max_time){ 381b05257ddSBarry Smith ts->time_step = ts->max_time - ts->ptime; 382e4dd521cSBarry Smith } 383e4dd521cSBarry Smith 384b05257ddSBarry Smith if (ts->time_step < 1e-14){ 385b05257ddSBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"Very small steps: %f\n",ts->time_step);CHKERRQ(ierr); 386b05257ddSBarry Smith ts->time_step = 1e-14; 387e4dd521cSBarry Smith } 388e4dd521cSBarry Smith 389e4dd521cSBarry Smith /* trying to purify h */ 390e4dd521cSBarry Smith /* (did not give any visible result) */ 391b05257ddSBarry Smith /* ttmp = ts->ptime + ts->time_step; 392b05257ddSBarry Smith ts->time_step = ttmp - ts->ptime; */ 393e4dd521cSBarry Smith 394e4dd521cSBarry Smith /* counting steps */ 395e4dd521cSBarry Smith ts->steps++; 396e4dd521cSBarry Smith } 397e4dd521cSBarry Smith 398e4dd521cSBarry Smith ierr=VecCopy(rk->y1,ts->vec_sol);CHKERRQ(ierr); 399e4dd521cSBarry Smith *steps += ts->steps; 400e4dd521cSBarry Smith *ptime = ts->ptime; 401e4dd521cSBarry Smith PetscFunctionReturn(0); 402e4dd521cSBarry Smith } 403e4dd521cSBarry Smith 404e4dd521cSBarry Smith /*------------------------------------------------------------*/ 405e4dd521cSBarry Smith #undef __FUNCT__ 406e4dd521cSBarry Smith #define __FUNCT__ "TSDestroy_Rk" 4076849ba73SBarry Smith static PetscErrorCode TSDestroy_Rk(TS ts) 408e4dd521cSBarry Smith { 409e4dd521cSBarry Smith TS_Rk *rk = (TS_Rk*)ts->data; 4106849ba73SBarry Smith PetscErrorCode ierr; 411a7cc72afSBarry Smith PetscInt i; 412e4dd521cSBarry Smith 413e4dd521cSBarry Smith /* REMEMBER TO DESTROY ALL */ 414e4dd521cSBarry Smith 415e4dd521cSBarry Smith PetscFunctionBegin; 416e4dd521cSBarry Smith if (rk->y1) {ierr = VecDestroy(rk->y1);CHKERRQ(ierr);} 417e4dd521cSBarry Smith if (rk->y2) {ierr = VecDestroy(rk->y2);CHKERRQ(ierr);} 418e4dd521cSBarry Smith if (rk->tmp) {ierr = VecDestroy(rk->tmp);CHKERRQ(ierr);} 419e4dd521cSBarry Smith if (rk->tmp_y) {ierr = VecDestroy(rk->tmp_y);CHKERRQ(ierr);} 420e4dd521cSBarry Smith for(i=0;i<rk->s;i++){ 421e4dd521cSBarry Smith if (rk->k[i]) {ierr = VecDestroy(rk->k[i]);CHKERRQ(ierr);} 422e4dd521cSBarry Smith } 423e4dd521cSBarry Smith ierr = PetscFree(rk);CHKERRQ(ierr); 424e4dd521cSBarry Smith PetscFunctionReturn(0); 425e4dd521cSBarry Smith } 426e4dd521cSBarry Smith /*------------------------------------------------------------*/ 427e4dd521cSBarry Smith 428e4dd521cSBarry Smith #undef __FUNCT__ 429e4dd521cSBarry Smith #define __FUNCT__ "TSSetFromOptions_Rk" 4306849ba73SBarry Smith static PetscErrorCode TSSetFromOptions_Rk(TS ts) 431e4dd521cSBarry Smith { 432262ff7bbSBarry Smith TS_Rk *rk = (TS_Rk*)ts->data; 433dfbe8321SBarry Smith PetscErrorCode ierr; 434262ff7bbSBarry Smith 435e4dd521cSBarry Smith PetscFunctionBegin; 436262ff7bbSBarry Smith ierr = PetscOptionsHead("RK ODE solver options");CHKERRQ(ierr); 437262ff7bbSBarry Smith ierr = PetscOptionsReal("-ts_rk_tol","Tolerance for convergence","TSRKSetTolerance",rk->tolerance,&rk->tolerance,PETSC_NULL);CHKERRQ(ierr); 438262ff7bbSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 439e4dd521cSBarry Smith PetscFunctionReturn(0); 440e4dd521cSBarry Smith } 441e4dd521cSBarry Smith 442e4dd521cSBarry Smith #undef __FUNCT__ 443e4dd521cSBarry Smith #define __FUNCT__ "TSView_Rk" 4446849ba73SBarry Smith static PetscErrorCode TSView_Rk(TS ts,PetscViewer viewer) 445e4dd521cSBarry Smith { 446e4dd521cSBarry Smith TS_Rk *rk = (TS_Rk*)ts->data; 447dfbe8321SBarry Smith PetscErrorCode ierr; 4482cdf8120Sasbjorn 449e4dd521cSBarry Smith PetscFunctionBegin; 45077431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD," number of ok steps: %D\n",rk->nok);CHKERRQ(ierr); 45177431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD," number of rejected steps: %D\n",rk->nnok);CHKERRQ(ierr); 452e4dd521cSBarry Smith PetscFunctionReturn(0); 453e4dd521cSBarry Smith } 454e4dd521cSBarry Smith 455e4dd521cSBarry Smith /* ------------------------------------------------------------ */ 456ebe3b25bSBarry Smith /*MC 457ebe3b25bSBarry Smith TS_RK - ODE solver using the explicit Runge-Kutta methods 458ebe3b25bSBarry Smith 459ebe3b25bSBarry Smith Options Database: 460ebe3b25bSBarry Smith . -ts_rk_tol <tol> Tolerance for convergence 461ebe3b25bSBarry Smith 462ebe3b25bSBarry Smith Contributed by: Asbjorn Hoiland Aarrestad, asbjorn@aarrestad.com, http://asbjorn.aarrestad.com/ 463ebe3b25bSBarry Smith 464d41222bbSBarry Smith Level: beginner 465d41222bbSBarry Smith 466ebe3b25bSBarry Smith .seealso: TSCreate(), TS, TSSetType(), TS_EULER, TSRKSetTolerance() 467ebe3b25bSBarry Smith 468ebe3b25bSBarry Smith M*/ 469ebe3b25bSBarry Smith 470e4dd521cSBarry Smith EXTERN_C_BEGIN 471e4dd521cSBarry Smith #undef __FUNCT__ 472e4dd521cSBarry Smith #define __FUNCT__ "TSCreate_Rk" 47363dd3a1aSKris Buschelman PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Rk(TS ts) 474e4dd521cSBarry Smith { 475e4dd521cSBarry Smith TS_Rk *rk; 476dfbe8321SBarry Smith PetscErrorCode ierr; 477e4dd521cSBarry Smith 478e4dd521cSBarry Smith PetscFunctionBegin; 479e4dd521cSBarry Smith ts->ops->setup = TSSetUp_Rk; 480e4dd521cSBarry Smith ts->ops->step = TSStep_Rk; 481e4dd521cSBarry Smith ts->ops->destroy = TSDestroy_Rk; 482e4dd521cSBarry Smith ts->ops->setfromoptions = TSSetFromOptions_Rk; 483e4dd521cSBarry Smith ts->ops->view = TSView_Rk; 484e4dd521cSBarry Smith 48538f2d2fdSLisandro Dalcin ierr = PetscNewLog(ts,TS_Rk,&rk);CHKERRQ(ierr); 486e4dd521cSBarry Smith ts->data = (void*)rk; 487e4dd521cSBarry Smith 488a7cc72afSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRKSetTolerance_C","TSRKSetTolerance_RK",TSRKSetTolerance_RK);CHKERRQ(ierr); 489262ff7bbSBarry Smith 490e4dd521cSBarry Smith PetscFunctionReturn(0); 491e4dd521cSBarry Smith } 492e4dd521cSBarry Smith EXTERN_C_END 493e4dd521cSBarry Smith 494e4dd521cSBarry Smith 495e4dd521cSBarry Smith 496e4dd521cSBarry Smith 497