163dd3a1aSKris Buschelman 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 */ 11*c6db04a5SJed Brown #include <private/tsimpl.h> /*I "petscts.h" I*/ 12*c6db04a5SJed Brown #include <time.h> 13e4dd521cSBarry Smith 14e4dd521cSBarry Smith typedef struct { 15e4dd521cSBarry Smith Vec y1,y2; /* work wectors for the two rk permuations */ 16a7cc72afSBarry Smith PetscInt 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) */ 19262ff7bbSBarry Smith PetscReal tolerance; /* initial value set for maxerror by user */ 20e4dd521cSBarry Smith Vec tmp,tmp_y,*k; /* two temp vectors and the k vectors for rk */ 21e4dd521cSBarry Smith PetscScalar a[7][6]; /* rk scalars */ 22e4dd521cSBarry Smith PetscScalar b1[7],b2[7]; /* rk scalars */ 23e4dd521cSBarry Smith PetscReal c[7]; /* rk scalars */ 24a7cc72afSBarry Smith PetscInt p,s; /* variables to tell the size of the runge-kutta solver */ 255f70b478SJed Brown } TS_RK; 26e4dd521cSBarry Smith 27262ff7bbSBarry Smith EXTERN_C_BEGIN 28262ff7bbSBarry Smith #undef __FUNCT__ 29262ff7bbSBarry Smith #define __FUNCT__ "TSRKSetTolerance_RK" 307087cfbeSBarry Smith PetscErrorCode TSRKSetTolerance_RK(TS ts,PetscReal aabs) 31262ff7bbSBarry Smith { 325f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 33262ff7bbSBarry Smith 34262ff7bbSBarry Smith PetscFunctionBegin; 35262ff7bbSBarry Smith rk->tolerance = aabs; 36262ff7bbSBarry Smith PetscFunctionReturn(0); 37262ff7bbSBarry Smith } 38262ff7bbSBarry Smith EXTERN_C_END 39262ff7bbSBarry Smith 40262ff7bbSBarry Smith #undef __FUNCT__ 41262ff7bbSBarry Smith #define __FUNCT__ "TSRKSetTolerance" 42262ff7bbSBarry Smith /*@ 43262ff7bbSBarry Smith TSRKSetTolerance - Sets the total error the RK explicit time integrators 44262ff7bbSBarry Smith will allow over the given time interval. 45262ff7bbSBarry Smith 463f9fe445SBarry Smith Logically Collective on TS 47262ff7bbSBarry Smith 48262ff7bbSBarry Smith Input parameters: 49262ff7bbSBarry Smith + ts - the time-step context 50262ff7bbSBarry Smith - aabs - the absolute tolerance 51262ff7bbSBarry Smith 52262ff7bbSBarry Smith Level: intermediate 53262ff7bbSBarry Smith 54262ff7bbSBarry Smith .keywords: RK, tolerance 55262ff7bbSBarry Smith 560f3b3ca1SHong Zhang .seealso: TSSundialsSetTolerance() 57262ff7bbSBarry Smith 58262ff7bbSBarry Smith @*/ 597087cfbeSBarry Smith PetscErrorCode TSRKSetTolerance(TS ts,PetscReal aabs) 60262ff7bbSBarry Smith { 614ac538c5SBarry Smith PetscErrorCode ierr; 62262ff7bbSBarry Smith 63262ff7bbSBarry Smith PetscFunctionBegin; 64c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,aabs,2); 654ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSRKSetTolerance_C",(TS,PetscReal),(ts,aabs));CHKERRQ(ierr); 66262ff7bbSBarry Smith PetscFunctionReturn(0); 67262ff7bbSBarry Smith } 68e4dd521cSBarry Smith 69e4dd521cSBarry Smith 70e4dd521cSBarry Smith #undef __FUNCT__ 715f70b478SJed Brown #define __FUNCT__ "TSSetUp_RK" 725f70b478SJed Brown static PetscErrorCode TSSetUp_RK(TS ts) 73e4dd521cSBarry Smith { 745f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 75a7cc72afSBarry Smith PetscErrorCode ierr; 76e4dd521cSBarry Smith 77e4dd521cSBarry Smith PetscFunctionBegin; 78e4dd521cSBarry Smith rk->nok = 0; 79e4dd521cSBarry Smith rk->nnok = 0; 80262ff7bbSBarry Smith rk->maxerror = rk->tolerance; 81e4dd521cSBarry Smith 82e4dd521cSBarry Smith /* fixing maxerror: global vs local */ 83e4dd521cSBarry Smith rk->ferror = rk->maxerror / (ts->max_time - ts->ptime); 84e4dd521cSBarry Smith 85e4dd521cSBarry Smith /* 34.0/45.0 gives double precision division */ 86e4dd521cSBarry Smith /* defining variables needed for Runge-Kutta computing */ 87e4dd521cSBarry Smith /* when changing below, please remember to change a, b1, b2 and c above! */ 88e4dd521cSBarry Smith /* Found in table on page 171: Dormand-Prince 5(4) */ 89e4dd521cSBarry Smith 90e4dd521cSBarry Smith /* are these right? */ 91e4dd521cSBarry Smith rk->p=6; 92e4dd521cSBarry Smith rk->s=7; 93e4dd521cSBarry Smith 94e4dd521cSBarry Smith rk->a[1][0]=1.0/5.0; 95e4dd521cSBarry Smith rk->a[2][0]=3.0/40.0; 96e4dd521cSBarry Smith rk->a[2][1]=9.0/40.0; 97e4dd521cSBarry Smith rk->a[3][0]=44.0/45.0; 98e4dd521cSBarry Smith rk->a[3][1]=-56.0/15.0; 99e4dd521cSBarry Smith rk->a[3][2]=32.0/9.0; 100e4dd521cSBarry Smith rk->a[4][0]=19372.0/6561.0; 101e4dd521cSBarry Smith rk->a[4][1]=-25360.0/2187.0; 102e4dd521cSBarry Smith rk->a[4][2]=64448.0/6561.0; 103e4dd521cSBarry Smith rk->a[4][3]=-212.0/729.0; 104e4dd521cSBarry Smith rk->a[5][0]=9017.0/3168.0; 105e4dd521cSBarry Smith rk->a[5][1]=-355.0/33.0; 106e4dd521cSBarry Smith rk->a[5][2]=46732.0/5247.0; 107e4dd521cSBarry Smith rk->a[5][3]=49.0/176.0; 108e4dd521cSBarry Smith rk->a[5][4]=-5103.0/18656.0; 109e4dd521cSBarry Smith rk->a[6][0]=35.0/384.0; 110e4dd521cSBarry Smith rk->a[6][1]=0.0; 111e4dd521cSBarry Smith rk->a[6][2]=500.0/1113.0; 112e4dd521cSBarry Smith rk->a[6][3]=125.0/192.0; 113e4dd521cSBarry Smith rk->a[6][4]=-2187.0/6784.0; 114e4dd521cSBarry Smith rk->a[6][5]=11.0/84.0; 115e4dd521cSBarry Smith 116e4dd521cSBarry Smith 117e4dd521cSBarry Smith rk->c[0]=0.0; 118e4dd521cSBarry Smith rk->c[1]=1.0/5.0; 119e4dd521cSBarry Smith rk->c[2]=3.0/10.0; 120e4dd521cSBarry Smith rk->c[3]=4.0/5.0; 121e4dd521cSBarry Smith rk->c[4]=8.0/9.0; 122e4dd521cSBarry Smith rk->c[5]=1.0; 123e4dd521cSBarry Smith rk->c[6]=1.0; 124e4dd521cSBarry Smith 125e4dd521cSBarry Smith rk->b1[0]=35.0/384.0; 126e4dd521cSBarry Smith rk->b1[1]=0.0; 127e4dd521cSBarry Smith rk->b1[2]=500.0/1113.0; 128e4dd521cSBarry Smith rk->b1[3]=125.0/192.0; 129e4dd521cSBarry Smith rk->b1[4]=-2187.0/6784.0; 130e4dd521cSBarry Smith rk->b1[5]=11.0/84.0; 131e4dd521cSBarry Smith rk->b1[6]=0.0; 132e4dd521cSBarry Smith 133e4dd521cSBarry Smith rk->b2[0]=5179.0/57600.0; 134e4dd521cSBarry Smith rk->b2[1]=0.0; 135e4dd521cSBarry Smith rk->b2[2]=7571.0/16695.0; 136e4dd521cSBarry Smith rk->b2[3]=393.0/640.0; 137e4dd521cSBarry Smith rk->b2[4]=-92097.0/339200.0; 138e4dd521cSBarry Smith rk->b2[5]=187.0/2100.0; 139e4dd521cSBarry Smith rk->b2[6]=1.0/40.0; 140e4dd521cSBarry Smith 141e4dd521cSBarry Smith 142e4dd521cSBarry Smith /* Found in table on page 170: Fehlberg 4(5) */ 143e4dd521cSBarry Smith /* 144e4dd521cSBarry Smith rk->p=5; 145e4dd521cSBarry Smith rk->s=6; 146e4dd521cSBarry Smith 147e4dd521cSBarry Smith rk->a[1][0]=1.0/4.0; 148e4dd521cSBarry Smith rk->a[2][0]=3.0/32.0; 149e4dd521cSBarry Smith rk->a[2][1]=9.0/32.0; 150e4dd521cSBarry Smith rk->a[3][0]=1932.0/2197.0; 151e4dd521cSBarry Smith rk->a[3][1]=-7200.0/2197.0; 152e4dd521cSBarry Smith rk->a[3][2]=7296.0/2197.0; 153e4dd521cSBarry Smith rk->a[4][0]=439.0/216.0; 154e4dd521cSBarry Smith rk->a[4][1]=-8.0; 155e4dd521cSBarry Smith rk->a[4][2]=3680.0/513.0; 156e4dd521cSBarry Smith rk->a[4][3]=-845.0/4104.0; 157e4dd521cSBarry Smith rk->a[5][0]=-8.0/27.0; 158e4dd521cSBarry Smith rk->a[5][1]=2.0; 159e4dd521cSBarry Smith rk->a[5][2]=-3544.0/2565.0; 160e4dd521cSBarry Smith rk->a[5][3]=1859.0/4104.0; 161e4dd521cSBarry Smith rk->a[5][4]=-11.0/40.0; 162e4dd521cSBarry Smith 163e4dd521cSBarry Smith rk->c[0]=0.0; 164e4dd521cSBarry Smith rk->c[1]=1.0/4.0; 165e4dd521cSBarry Smith rk->c[2]=3.0/8.0; 166e4dd521cSBarry Smith rk->c[3]=12.0/13.0; 167e4dd521cSBarry Smith rk->c[4]=1.0; 168e4dd521cSBarry Smith rk->c[5]=1.0/2.0; 169e4dd521cSBarry Smith 170e4dd521cSBarry Smith rk->b1[0]=25.0/216.0; 171e4dd521cSBarry Smith rk->b1[1]=0.0; 172e4dd521cSBarry Smith rk->b1[2]=1408.0/2565.0; 173e4dd521cSBarry Smith rk->b1[3]=2197.0/4104.0; 174e4dd521cSBarry Smith rk->b1[4]=-1.0/5.0; 175e4dd521cSBarry Smith rk->b1[5]=0.0; 176e4dd521cSBarry Smith 177e4dd521cSBarry Smith rk->b2[0]=16.0/135.0; 178e4dd521cSBarry Smith rk->b2[1]=0.0; 179e4dd521cSBarry Smith rk->b2[2]=6656.0/12825.0; 180e4dd521cSBarry Smith rk->b2[3]=28561.0/56430.0; 181e4dd521cSBarry Smith rk->b2[4]=-9.0/50.0; 182e4dd521cSBarry Smith rk->b2[5]=2.0/55.0; 183e4dd521cSBarry Smith */ 184e4dd521cSBarry Smith /* Found in table on page 169: Merson 4("5") */ 185e4dd521cSBarry Smith /* 186e4dd521cSBarry Smith rk->p=4; 187e4dd521cSBarry Smith rk->s=5; 188e4dd521cSBarry Smith rk->a[1][0] = 1.0/3.0; 189e4dd521cSBarry Smith rk->a[2][0] = 1.0/6.0; 190e4dd521cSBarry Smith rk->a[2][1] = 1.0/6.0; 191e4dd521cSBarry Smith rk->a[3][0] = 1.0/8.0; 192e4dd521cSBarry Smith rk->a[3][1] = 0.0; 193e4dd521cSBarry Smith rk->a[3][2] = 3.0/8.0; 194e4dd521cSBarry Smith rk->a[4][0] = 1.0/2.0; 195e4dd521cSBarry Smith rk->a[4][1] = 0.0; 196e4dd521cSBarry Smith rk->a[4][2] = -3.0/2.0; 197e4dd521cSBarry Smith rk->a[4][3] = 2.0; 198e4dd521cSBarry Smith 199e4dd521cSBarry Smith rk->c[0] = 0.0; 200e4dd521cSBarry Smith rk->c[1] = 1.0/3.0; 201e4dd521cSBarry Smith rk->c[2] = 1.0/3.0; 202e4dd521cSBarry Smith rk->c[3] = 0.5; 203e4dd521cSBarry Smith rk->c[4] = 1.0; 204e4dd521cSBarry Smith 205e4dd521cSBarry Smith rk->b1[0] = 1.0/2.0; 206e4dd521cSBarry Smith rk->b1[1] = 0.0; 207e4dd521cSBarry Smith rk->b1[2] = -3.0/2.0; 208e4dd521cSBarry Smith rk->b1[3] = 2.0; 209e4dd521cSBarry Smith rk->b1[4] = 0.0; 210e4dd521cSBarry Smith 211e4dd521cSBarry Smith rk->b2[0] = 1.0/6.0; 212e4dd521cSBarry Smith rk->b2[1] = 0.0; 213e4dd521cSBarry Smith rk->b2[2] = 0.0; 214e4dd521cSBarry Smith rk->b2[3] = 2.0/3.0; 215e4dd521cSBarry Smith rk->b2[4] = 1.0/6.0; 216e4dd521cSBarry Smith */ 217e4dd521cSBarry Smith 218e4dd521cSBarry Smith /* making b2 -> e=b1-b2 */ 219e4dd521cSBarry Smith /* 220e4dd521cSBarry Smith for(i=0;i<rk->s;i++){ 22126edb414Sasbjorn rk->b2[i] = (rk->b1[i]) - (rk->b2[i]); 222e4dd521cSBarry Smith } 223e4dd521cSBarry Smith */ 22426edb414Sasbjorn rk->b2[0]=71.0/57600.0; 22526edb414Sasbjorn rk->b2[1]=0.0; 22626edb414Sasbjorn rk->b2[2]=-71.0/16695.0; 22726edb414Sasbjorn rk->b2[3]=71.0/1920.0; 22826edb414Sasbjorn rk->b2[4]=-17253.0/339200.0; 22926edb414Sasbjorn rk->b2[5]=22.0/525.0; 23026edb414Sasbjorn rk->b2[6]=-1.0/40.0; 23126edb414Sasbjorn 232e4dd521cSBarry Smith /* initializing vectors */ 233e4dd521cSBarry Smith ierr = VecDuplicate(ts->vec_sol,&rk->y1);CHKERRQ(ierr); 234e4dd521cSBarry Smith ierr = VecDuplicate(ts->vec_sol,&rk->y2);CHKERRQ(ierr); 235e4dd521cSBarry Smith ierr = VecDuplicate(rk->y1,&rk->tmp);CHKERRQ(ierr); 236e4dd521cSBarry Smith ierr = VecDuplicate(rk->y1,&rk->tmp_y);CHKERRQ(ierr); 237e4dd521cSBarry Smith ierr = VecDuplicateVecs(rk->y1,rk->s,&rk->k);CHKERRQ(ierr); 238e4dd521cSBarry Smith 239e4dd521cSBarry Smith PetscFunctionReturn(0); 240e4dd521cSBarry Smith } 241e4dd521cSBarry Smith 242db046809SBarry Smith /*------------------------------------------------------------*/ 243db046809SBarry Smith #undef __FUNCT__ 2445f70b478SJed Brown #define __FUNCT__ "TSRKqs" 2455f70b478SJed Brown PetscErrorCode TSRKqs(TS ts,PetscReal t,PetscReal h) 246db046809SBarry Smith { 2475f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 2486849ba73SBarry Smith PetscErrorCode ierr; 249a7cc72afSBarry Smith PetscInt j,l; 250db046809SBarry Smith PetscReal tmp_t=t; 251b05257ddSBarry Smith PetscScalar hh=h; 252db046809SBarry Smith 253db046809SBarry Smith PetscFunctionBegin; 254db046809SBarry Smith /* k[0]=0 */ 255b05257ddSBarry Smith ierr = VecSet(rk->k[0],0.0);CHKERRQ(ierr); 256db046809SBarry Smith 257db046809SBarry Smith /* k[0] = derivs(t,y1) */ 258db046809SBarry Smith ierr = TSComputeRHSFunction(ts,t,rk->y1,rk->k[0]);CHKERRQ(ierr); 259db046809SBarry Smith /* looping over runge-kutta variables */ 260db046809SBarry Smith /* building the k - array of vectors */ 261db046809SBarry Smith for(j = 1 ; j < rk->s ; j++){ 262db046809SBarry Smith 263db046809SBarry Smith /* rk->tmp = 0 */ 264b05257ddSBarry Smith ierr = VecSet(rk->tmp,0.0);CHKERRQ(ierr); 265db046809SBarry Smith 266db046809SBarry Smith for(l=0;l<j;l++){ 267db046809SBarry Smith /* tmp += a(j,l)*k[l] */ 2682dcb1b2aSMatthew Knepley ierr = VecAXPY(rk->tmp,rk->a[j][l],rk->k[l]);CHKERRQ(ierr); 269db046809SBarry Smith } 270db046809SBarry Smith 271db046809SBarry Smith /* ierr = VecView(rk->tmp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 272db046809SBarry Smith 273db046809SBarry Smith /* k[j] = derivs(t+c(j)*h,y1+h*tmp,k(j)) */ 274db046809SBarry Smith /* I need the following helpers: 275db046809SBarry Smith PetscScalar tmp_t=t+c(j)*h 276db046809SBarry Smith Vec tmp_y=h*tmp+y1 277db046809SBarry Smith */ 278db046809SBarry Smith 279db046809SBarry Smith tmp_t = t + rk->c[j] * h; 280db046809SBarry Smith 281db046809SBarry Smith /* tmp_y = h * tmp + y1 */ 2822dcb1b2aSMatthew Knepley ierr = VecWAXPY(rk->tmp_y,hh,rk->tmp,rk->y1);CHKERRQ(ierr); 283db046809SBarry Smith 284db046809SBarry Smith /* rk->k[j]=0 */ 285b05257ddSBarry Smith ierr = VecSet(rk->k[j],0.0);CHKERRQ(ierr); 286db046809SBarry Smith ierr = TSComputeRHSFunction(ts,tmp_t,rk->tmp_y,rk->k[j]);CHKERRQ(ierr); 287db046809SBarry Smith } 288db046809SBarry Smith 289db046809SBarry Smith /* tmp=0 and tmp_y=0 */ 290b05257ddSBarry Smith ierr = VecSet(rk->tmp,0.0);CHKERRQ(ierr); 291b05257ddSBarry Smith ierr = VecSet(rk->tmp_y,0.0);CHKERRQ(ierr); 292db046809SBarry Smith 293db046809SBarry Smith for(j = 0 ; j < rk->s ; j++){ 294db046809SBarry Smith /* tmp=b1[j]*k[j]+tmp */ 2952dcb1b2aSMatthew Knepley ierr = VecAXPY(rk->tmp,rk->b1[j],rk->k[j]);CHKERRQ(ierr); 296db046809SBarry Smith /* tmp_y=b2[j]*k[j]+tmp_y */ 2972dcb1b2aSMatthew Knepley ierr = VecAXPY(rk->tmp_y,rk->b2[j],rk->k[j]);CHKERRQ(ierr); 298db046809SBarry Smith } 299db046809SBarry Smith 300db046809SBarry Smith /* y2 = hh * tmp_y */ 301b05257ddSBarry Smith ierr = VecSet(rk->y2,0.0);CHKERRQ(ierr); 3022dcb1b2aSMatthew Knepley ierr = VecAXPY(rk->y2,hh,rk->tmp_y);CHKERRQ(ierr); 303db046809SBarry Smith /* y1 = hh*tmp + y1 */ 3042dcb1b2aSMatthew Knepley ierr = VecAXPY(rk->y1,hh,rk->tmp);CHKERRQ(ierr); 305db046809SBarry Smith /* Finding difference between y1 and y2 */ 306db046809SBarry Smith PetscFunctionReturn(0); 307db046809SBarry Smith } 308db046809SBarry Smith 309e4dd521cSBarry Smith #undef __FUNCT__ 3105f70b478SJed Brown #define __FUNCT__ "TSStep_RK" 3115f70b478SJed Brown static PetscErrorCode TSStep_RK(TS ts,PetscInt *steps,PetscReal *ptime) 312e4dd521cSBarry Smith { 3135f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 314a7cc72afSBarry Smith PetscErrorCode ierr; 315b0dd9724SMatthew Knepley PetscReal norm=0.0,dt_fac=0.0,fac = 0.0/*,ttmp=0.0*/; 316e3caeda6SBarry Smith PetscInt i, max_steps = ts->max_steps; 317e4dd521cSBarry Smith 318e4dd521cSBarry Smith PetscFunctionBegin; 319e4dd521cSBarry Smith ierr=VecCopy(ts->vec_sol,rk->y1);CHKERRQ(ierr); 320e4dd521cSBarry Smith *steps = -ts->steps; 321e4dd521cSBarry Smith /* trying to save the vector */ 322e4dd521cSBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 323e4dd521cSBarry Smith /* while loop to get from start to stop */ 324e3caeda6SBarry Smith for (i = 0; i < max_steps; i++) { 3253f2090d5SJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); /* Note that this is called once per STEP, not once per STAGE. */ 326e4dd521cSBarry Smith /* calling rkqs */ 327e4dd521cSBarry Smith /* 328e4dd521cSBarry Smith -- input 329e4dd521cSBarry Smith ts - pointer to ts 330e4dd521cSBarry Smith ts->ptime - current time 331b05257ddSBarry Smith ts->time_step - try this timestep 332e4dd521cSBarry Smith y1 - solution for this step 333e4dd521cSBarry Smith 334e4dd521cSBarry Smith --output 335e4dd521cSBarry Smith y1 - suggested solution 336e4dd521cSBarry Smith y2 - check solution (runge - kutta second permutation) 337e4dd521cSBarry Smith */ 3385f70b478SJed Brown ierr = TSRKqs(ts,ts->ptime,ts->time_step);CHKERRQ(ierr); 339e3caeda6SBarry Smith /* counting steps */ 340e3caeda6SBarry Smith ts->steps++; 341e4dd521cSBarry Smith /* checking for maxerror */ 342e4dd521cSBarry Smith /* comparing difference to maxerror */ 343e4dd521cSBarry Smith ierr = VecNorm(rk->y2,NORM_2,&norm);CHKERRQ(ierr); 344e4dd521cSBarry Smith /* modifying maxerror to satisfy this timestep */ 345b05257ddSBarry Smith rk->maxerror = rk->ferror * ts->time_step; 346b05257ddSBarry Smith /* ierr = PetscPrintf(PETSC_COMM_WORLD,"norm err: %f maxerror: %f dt: %f",norm,rk->maxerror,ts->time_step);CHKERRQ(ierr); */ 347e4dd521cSBarry Smith 348e4dd521cSBarry Smith /* handling ok and not ok */ 349e4dd521cSBarry Smith if (norm < rk->maxerror){ 350e4dd521cSBarry Smith /* if ok: */ 351e4dd521cSBarry Smith ierr=VecCopy(rk->y1,ts->vec_sol);CHKERRQ(ierr); /* saves the suggested solution to current solution */ 352b05257ddSBarry Smith ts->ptime += ts->time_step; /* storing the new current time */ 353e4dd521cSBarry Smith rk->nok++; 354e4dd521cSBarry Smith fac=5.0; 355e4dd521cSBarry Smith /* trying to save the vector */ 3563f2090d5SJed Brown ierr = TSPostStep(ts);CHKERRQ(ierr); 357e4dd521cSBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 358e3caeda6SBarry Smith if (ts->ptime >= ts->max_time) break; 359e4dd521cSBarry Smith } else{ 360e4dd521cSBarry Smith /* if not OK */ 361e4dd521cSBarry Smith rk->nnok++; 362e4dd521cSBarry Smith fac=1.0; 363e4dd521cSBarry Smith ierr=VecCopy(ts->vec_sol,rk->y1);CHKERRQ(ierr); /* restores old solution */ 364e4dd521cSBarry Smith } 365e4dd521cSBarry Smith 366e4dd521cSBarry Smith /*Computing next stepsize. See page 167 in Solving ODE 1 367e4dd521cSBarry Smith * 368e4dd521cSBarry Smith * h_new = h * min( facmax , max( facmin , fac * (tol/err)^(1/(p+1)) ) ) 369e4dd521cSBarry Smith * facmax set above 370e4dd521cSBarry Smith * facmin 371e4dd521cSBarry Smith */ 372e4dd521cSBarry Smith dt_fac = exp(log((rk->maxerror) / norm) / ((rk->p) + 1) ) * 0.9 ; 373e4dd521cSBarry Smith 374e4dd521cSBarry Smith if (dt_fac > fac){ 3752cdf8120Sasbjorn /*ierr = PetscPrintf(PETSC_COMM_WORLD,"changing fac %f\n",fac);*/ 376e4dd521cSBarry Smith dt_fac = fac; 377e4dd521cSBarry Smith } 378e4dd521cSBarry Smith 379b05257ddSBarry Smith /* computing new ts->time_step */ 380b05257ddSBarry Smith ts->time_step = ts->time_step * dt_fac; 381e4dd521cSBarry Smith 382b05257ddSBarry Smith if (ts->ptime+ts->time_step > ts->max_time){ 383b05257ddSBarry Smith ts->time_step = ts->max_time - ts->ptime; 384e4dd521cSBarry Smith } 385e4dd521cSBarry Smith 386b05257ddSBarry Smith if (ts->time_step < 1e-14){ 387b05257ddSBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"Very small steps: %f\n",ts->time_step);CHKERRQ(ierr); 388b05257ddSBarry Smith ts->time_step = 1e-14; 389e4dd521cSBarry Smith } 390e4dd521cSBarry Smith 391e4dd521cSBarry Smith /* trying to purify h */ 392e4dd521cSBarry Smith /* (did not give any visible result) */ 393b05257ddSBarry Smith /* ttmp = ts->ptime + ts->time_step; 394b05257ddSBarry Smith ts->time_step = ttmp - ts->ptime; */ 395e4dd521cSBarry Smith 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__ 4065f70b478SJed Brown #define __FUNCT__ "TSDestroy_RK" 4075f70b478SJed Brown static PetscErrorCode TSDestroy_RK(TS ts) 408e4dd521cSBarry Smith { 4095f70b478SJed Brown 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__ 4295f70b478SJed Brown #define __FUNCT__ "TSSetFromOptions_RK" 4305f70b478SJed Brown static PetscErrorCode TSSetFromOptions_RK(TS ts) 431e4dd521cSBarry Smith { 4325f70b478SJed Brown 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__ 4435f70b478SJed Brown #define __FUNCT__ "TSView_RK" 4445f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer) 445e4dd521cSBarry Smith { 4465f70b478SJed Brown 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 45796f5712cSJed Brown TSRK - 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 4669596e0b4SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSRKSetTolerance() 467ebe3b25bSBarry Smith 468ebe3b25bSBarry Smith M*/ 469ebe3b25bSBarry Smith 470e4dd521cSBarry Smith EXTERN_C_BEGIN 471e4dd521cSBarry Smith #undef __FUNCT__ 4725f70b478SJed Brown #define __FUNCT__ "TSCreate_RK" 4737087cfbeSBarry Smith PetscErrorCode TSCreate_RK(TS ts) 474e4dd521cSBarry Smith { 4755f70b478SJed Brown TS_RK *rk; 476dfbe8321SBarry Smith PetscErrorCode ierr; 477e4dd521cSBarry Smith 478e4dd521cSBarry Smith PetscFunctionBegin; 4795f70b478SJed Brown ts->ops->setup = TSSetUp_RK; 4805f70b478SJed Brown ts->ops->step = TSStep_RK; 4815f70b478SJed Brown ts->ops->destroy = TSDestroy_RK; 4825f70b478SJed Brown ts->ops->setfromoptions = TSSetFromOptions_RK; 4835f70b478SJed Brown ts->ops->view = TSView_RK; 484e4dd521cSBarry Smith 4855f70b478SJed Brown 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