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 */ 10b45d2f2cSJed Brown #include <petsc-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 #undef __FUNCT__ 27262ff7bbSBarry Smith #define __FUNCT__ "TSRKSetTolerance_RK" 287087cfbeSBarry Smith PetscErrorCode TSRKSetTolerance_RK(TS ts,PetscReal aabs) 29262ff7bbSBarry Smith { 305f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 31262ff7bbSBarry Smith 32262ff7bbSBarry Smith PetscFunctionBegin; 33262ff7bbSBarry Smith rk->tolerance = aabs; 34262ff7bbSBarry Smith PetscFunctionReturn(0); 35262ff7bbSBarry Smith } 36262ff7bbSBarry Smith 37262ff7bbSBarry Smith #undef __FUNCT__ 38262ff7bbSBarry Smith #define __FUNCT__ "TSRKSetTolerance" 39262ff7bbSBarry Smith /*@ 40262ff7bbSBarry Smith TSRKSetTolerance - Sets the total error the RK explicit time integrators 41262ff7bbSBarry Smith will allow over the given time interval. 42262ff7bbSBarry Smith 433f9fe445SBarry Smith Logically Collective on TS 44262ff7bbSBarry Smith 45262ff7bbSBarry Smith Input parameters: 46262ff7bbSBarry Smith + ts - the time-step context 47262ff7bbSBarry Smith - aabs - the absolute tolerance 48262ff7bbSBarry Smith 49262ff7bbSBarry Smith Level: intermediate 50262ff7bbSBarry Smith 51262ff7bbSBarry Smith .keywords: RK, tolerance 52262ff7bbSBarry Smith 530f3b3ca1SHong Zhang .seealso: TSSundialsSetTolerance() 54262ff7bbSBarry Smith 55262ff7bbSBarry Smith @*/ 567087cfbeSBarry Smith PetscErrorCode TSRKSetTolerance(TS ts,PetscReal aabs) 57262ff7bbSBarry Smith { 584ac538c5SBarry Smith PetscErrorCode ierr; 59262ff7bbSBarry Smith 60262ff7bbSBarry Smith PetscFunctionBegin; 61c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(ts,aabs,2); 624ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSRKSetTolerance_C",(TS,PetscReal),(ts,aabs));CHKERRQ(ierr); 63262ff7bbSBarry Smith PetscFunctionReturn(0); 64262ff7bbSBarry Smith } 65e4dd521cSBarry Smith 66e4dd521cSBarry Smith 67e4dd521cSBarry Smith #undef __FUNCT__ 685f70b478SJed Brown #define __FUNCT__ "TSSetUp_RK" 695f70b478SJed Brown static PetscErrorCode TSSetUp_RK(TS ts) 70e4dd521cSBarry Smith { 715f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 72a7cc72afSBarry Smith PetscErrorCode ierr; 73e4dd521cSBarry Smith 74e4dd521cSBarry Smith PetscFunctionBegin; 75e4dd521cSBarry Smith rk->nok = 0; 76e4dd521cSBarry Smith rk->nnok = 0; 77262ff7bbSBarry Smith rk->maxerror = rk->tolerance; 78e4dd521cSBarry Smith 79e4dd521cSBarry Smith /* fixing maxerror: global vs local */ 80e4dd521cSBarry Smith rk->ferror = rk->maxerror / (ts->max_time - ts->ptime); 81e4dd521cSBarry Smith 82e4dd521cSBarry Smith /* 34.0/45.0 gives double precision division */ 83e4dd521cSBarry Smith /* defining variables needed for Runge-Kutta computing */ 84e4dd521cSBarry Smith /* when changing below, please remember to change a, b1, b2 and c above! */ 85e4dd521cSBarry Smith /* Found in table on page 171: Dormand-Prince 5(4) */ 86e4dd521cSBarry Smith 87e4dd521cSBarry Smith /* are these right? */ 88e4dd521cSBarry Smith rk->p=6; 89e4dd521cSBarry Smith rk->s=7; 90e4dd521cSBarry Smith 91e4dd521cSBarry Smith rk->a[1][0]=1.0/5.0; 92e4dd521cSBarry Smith rk->a[2][0]=3.0/40.0; 93e4dd521cSBarry Smith rk->a[2][1]=9.0/40.0; 94e4dd521cSBarry Smith rk->a[3][0]=44.0/45.0; 95e4dd521cSBarry Smith rk->a[3][1]=-56.0/15.0; 96e4dd521cSBarry Smith rk->a[3][2]=32.0/9.0; 97e4dd521cSBarry Smith rk->a[4][0]=19372.0/6561.0; 98e4dd521cSBarry Smith rk->a[4][1]=-25360.0/2187.0; 99e4dd521cSBarry Smith rk->a[4][2]=64448.0/6561.0; 100e4dd521cSBarry Smith rk->a[4][3]=-212.0/729.0; 101e4dd521cSBarry Smith rk->a[5][0]=9017.0/3168.0; 102e4dd521cSBarry Smith rk->a[5][1]=-355.0/33.0; 103e4dd521cSBarry Smith rk->a[5][2]=46732.0/5247.0; 104e4dd521cSBarry Smith rk->a[5][3]=49.0/176.0; 105e4dd521cSBarry Smith rk->a[5][4]=-5103.0/18656.0; 106e4dd521cSBarry Smith rk->a[6][0]=35.0/384.0; 107e4dd521cSBarry Smith rk->a[6][1]=0.0; 108e4dd521cSBarry Smith rk->a[6][2]=500.0/1113.0; 109e4dd521cSBarry Smith rk->a[6][3]=125.0/192.0; 110e4dd521cSBarry Smith rk->a[6][4]=-2187.0/6784.0; 111e4dd521cSBarry Smith rk->a[6][5]=11.0/84.0; 112e4dd521cSBarry Smith 113e4dd521cSBarry Smith 114e4dd521cSBarry Smith rk->c[0]=0.0; 115e4dd521cSBarry Smith rk->c[1]=1.0/5.0; 116e4dd521cSBarry Smith rk->c[2]=3.0/10.0; 117e4dd521cSBarry Smith rk->c[3]=4.0/5.0; 118e4dd521cSBarry Smith rk->c[4]=8.0/9.0; 119e4dd521cSBarry Smith rk->c[5]=1.0; 120e4dd521cSBarry Smith rk->c[6]=1.0; 121e4dd521cSBarry Smith 122e4dd521cSBarry Smith rk->b1[0]=35.0/384.0; 123e4dd521cSBarry Smith rk->b1[1]=0.0; 124e4dd521cSBarry Smith rk->b1[2]=500.0/1113.0; 125e4dd521cSBarry Smith rk->b1[3]=125.0/192.0; 126e4dd521cSBarry Smith rk->b1[4]=-2187.0/6784.0; 127e4dd521cSBarry Smith rk->b1[5]=11.0/84.0; 128e4dd521cSBarry Smith rk->b1[6]=0.0; 129e4dd521cSBarry Smith 130e4dd521cSBarry Smith rk->b2[0]=5179.0/57600.0; 131e4dd521cSBarry Smith rk->b2[1]=0.0; 132e4dd521cSBarry Smith rk->b2[2]=7571.0/16695.0; 133e4dd521cSBarry Smith rk->b2[3]=393.0/640.0; 134e4dd521cSBarry Smith rk->b2[4]=-92097.0/339200.0; 135e4dd521cSBarry Smith rk->b2[5]=187.0/2100.0; 136e4dd521cSBarry Smith rk->b2[6]=1.0/40.0; 137e4dd521cSBarry Smith 138e4dd521cSBarry Smith 139e4dd521cSBarry Smith /* Found in table on page 170: Fehlberg 4(5) */ 140e4dd521cSBarry Smith /* 141e4dd521cSBarry Smith rk->p=5; 142e4dd521cSBarry Smith rk->s=6; 143e4dd521cSBarry Smith 144e4dd521cSBarry Smith rk->a[1][0]=1.0/4.0; 145e4dd521cSBarry Smith rk->a[2][0]=3.0/32.0; 146e4dd521cSBarry Smith rk->a[2][1]=9.0/32.0; 147e4dd521cSBarry Smith rk->a[3][0]=1932.0/2197.0; 148e4dd521cSBarry Smith rk->a[3][1]=-7200.0/2197.0; 149e4dd521cSBarry Smith rk->a[3][2]=7296.0/2197.0; 150e4dd521cSBarry Smith rk->a[4][0]=439.0/216.0; 151e4dd521cSBarry Smith rk->a[4][1]=-8.0; 152e4dd521cSBarry Smith rk->a[4][2]=3680.0/513.0; 153e4dd521cSBarry Smith rk->a[4][3]=-845.0/4104.0; 154e4dd521cSBarry Smith rk->a[5][0]=-8.0/27.0; 155e4dd521cSBarry Smith rk->a[5][1]=2.0; 156e4dd521cSBarry Smith rk->a[5][2]=-3544.0/2565.0; 157e4dd521cSBarry Smith rk->a[5][3]=1859.0/4104.0; 158e4dd521cSBarry Smith rk->a[5][4]=-11.0/40.0; 159e4dd521cSBarry Smith 160e4dd521cSBarry Smith rk->c[0]=0.0; 161e4dd521cSBarry Smith rk->c[1]=1.0/4.0; 162e4dd521cSBarry Smith rk->c[2]=3.0/8.0; 163e4dd521cSBarry Smith rk->c[3]=12.0/13.0; 164e4dd521cSBarry Smith rk->c[4]=1.0; 165e4dd521cSBarry Smith rk->c[5]=1.0/2.0; 166e4dd521cSBarry Smith 167e4dd521cSBarry Smith rk->b1[0]=25.0/216.0; 168e4dd521cSBarry Smith rk->b1[1]=0.0; 169e4dd521cSBarry Smith rk->b1[2]=1408.0/2565.0; 170e4dd521cSBarry Smith rk->b1[3]=2197.0/4104.0; 171e4dd521cSBarry Smith rk->b1[4]=-1.0/5.0; 172e4dd521cSBarry Smith rk->b1[5]=0.0; 173e4dd521cSBarry Smith 174e4dd521cSBarry Smith rk->b2[0]=16.0/135.0; 175e4dd521cSBarry Smith rk->b2[1]=0.0; 176e4dd521cSBarry Smith rk->b2[2]=6656.0/12825.0; 177e4dd521cSBarry Smith rk->b2[3]=28561.0/56430.0; 178e4dd521cSBarry Smith rk->b2[4]=-9.0/50.0; 179e4dd521cSBarry Smith rk->b2[5]=2.0/55.0; 180e4dd521cSBarry Smith */ 181e4dd521cSBarry Smith /* Found in table on page 169: Merson 4("5") */ 182e4dd521cSBarry Smith /* 183e4dd521cSBarry Smith rk->p=4; 184e4dd521cSBarry Smith rk->s=5; 185e4dd521cSBarry Smith rk->a[1][0] = 1.0/3.0; 186e4dd521cSBarry Smith rk->a[2][0] = 1.0/6.0; 187e4dd521cSBarry Smith rk->a[2][1] = 1.0/6.0; 188e4dd521cSBarry Smith rk->a[3][0] = 1.0/8.0; 189e4dd521cSBarry Smith rk->a[3][1] = 0.0; 190e4dd521cSBarry Smith rk->a[3][2] = 3.0/8.0; 191e4dd521cSBarry Smith rk->a[4][0] = 1.0/2.0; 192e4dd521cSBarry Smith rk->a[4][1] = 0.0; 193e4dd521cSBarry Smith rk->a[4][2] = -3.0/2.0; 194e4dd521cSBarry Smith rk->a[4][3] = 2.0; 195e4dd521cSBarry Smith 196e4dd521cSBarry Smith rk->c[0] = 0.0; 197e4dd521cSBarry Smith rk->c[1] = 1.0/3.0; 198e4dd521cSBarry Smith rk->c[2] = 1.0/3.0; 199e4dd521cSBarry Smith rk->c[3] = 0.5; 200e4dd521cSBarry Smith rk->c[4] = 1.0; 201e4dd521cSBarry Smith 202e4dd521cSBarry Smith rk->b1[0] = 1.0/2.0; 203e4dd521cSBarry Smith rk->b1[1] = 0.0; 204e4dd521cSBarry Smith rk->b1[2] = -3.0/2.0; 205e4dd521cSBarry Smith rk->b1[3] = 2.0; 206e4dd521cSBarry Smith rk->b1[4] = 0.0; 207e4dd521cSBarry Smith 208e4dd521cSBarry Smith rk->b2[0] = 1.0/6.0; 209e4dd521cSBarry Smith rk->b2[1] = 0.0; 210e4dd521cSBarry Smith rk->b2[2] = 0.0; 211e4dd521cSBarry Smith rk->b2[3] = 2.0/3.0; 212e4dd521cSBarry Smith rk->b2[4] = 1.0/6.0; 213e4dd521cSBarry Smith */ 214e4dd521cSBarry Smith 215e4dd521cSBarry Smith /* making b2 -> e=b1-b2 */ 216e4dd521cSBarry Smith /* 217e4dd521cSBarry Smith for (i=0;i<rk->s;i++) { 21826edb414Sasbjorn rk->b2[i] = (rk->b1[i]) - (rk->b2[i]); 219e4dd521cSBarry Smith } 220e4dd521cSBarry Smith */ 22126edb414Sasbjorn rk->b2[0]=71.0/57600.0; 22226edb414Sasbjorn rk->b2[1]=0.0; 22326edb414Sasbjorn rk->b2[2]=-71.0/16695.0; 22426edb414Sasbjorn rk->b2[3]=71.0/1920.0; 22526edb414Sasbjorn rk->b2[4]=-17253.0/339200.0; 22626edb414Sasbjorn rk->b2[5]=22.0/525.0; 22726edb414Sasbjorn rk->b2[6]=-1.0/40.0; 22826edb414Sasbjorn 229e4dd521cSBarry Smith /* initializing vectors */ 230e4dd521cSBarry Smith ierr = VecDuplicate(ts->vec_sol,&rk->y1);CHKERRQ(ierr); 231e4dd521cSBarry Smith ierr = VecDuplicate(ts->vec_sol,&rk->y2);CHKERRQ(ierr); 232e4dd521cSBarry Smith ierr = VecDuplicate(rk->y1,&rk->tmp);CHKERRQ(ierr); 233e4dd521cSBarry Smith ierr = VecDuplicate(rk->y1,&rk->tmp_y);CHKERRQ(ierr); 234e4dd521cSBarry Smith ierr = VecDuplicateVecs(rk->y1,rk->s,&rk->k);CHKERRQ(ierr); 235e4dd521cSBarry Smith PetscFunctionReturn(0); 236e4dd521cSBarry Smith } 237e4dd521cSBarry Smith 238db046809SBarry Smith /*------------------------------------------------------------*/ 239db046809SBarry Smith #undef __FUNCT__ 2405f70b478SJed Brown #define __FUNCT__ "TSRKqs" 2415f70b478SJed Brown PetscErrorCode TSRKqs(TS ts,PetscReal t,PetscReal h) 242db046809SBarry Smith { 2435f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 2446849ba73SBarry Smith PetscErrorCode ierr; 245a7cc72afSBarry Smith PetscInt j,l; 246db046809SBarry Smith PetscReal tmp_t = t; 247b05257ddSBarry Smith PetscScalar hh = h; 248db046809SBarry Smith 249db046809SBarry Smith PetscFunctionBegin; 250db046809SBarry Smith /* k[0]=0 */ 251b05257ddSBarry Smith ierr = VecSet(rk->k[0],0.0);CHKERRQ(ierr); 252db046809SBarry Smith 253db046809SBarry Smith /* k[0] = derivs(t,y1) */ 254db046809SBarry Smith ierr = TSComputeRHSFunction(ts,t,rk->y1,rk->k[0]);CHKERRQ(ierr); 255db046809SBarry Smith /* looping over runge-kutta variables */ 256db046809SBarry Smith /* building the k - array of vectors */ 257db046809SBarry Smith for (j = 1; j < rk->s; j++) { 258db046809SBarry Smith 259db046809SBarry Smith /* rk->tmp = 0 */ 260b05257ddSBarry Smith ierr = VecSet(rk->tmp,0.0);CHKERRQ(ierr); 261db046809SBarry Smith 262db046809SBarry Smith for (l=0; l<j; l++) { 263db046809SBarry Smith /* tmp += a(j,l)*k[l] */ 2642dcb1b2aSMatthew Knepley ierr = VecAXPY(rk->tmp,rk->a[j][l],rk->k[l]);CHKERRQ(ierr); 265db046809SBarry Smith } 266db046809SBarry Smith 267db046809SBarry Smith /* ierr = VecView(rk->tmp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 268db046809SBarry Smith 269db046809SBarry Smith /* k[j] = derivs(t+c(j)*h,y1+h*tmp,k(j)) */ 270db046809SBarry Smith /* I need the following helpers: 271db046809SBarry Smith PetscScalar tmp_t=t+c(j)*h 272db046809SBarry Smith Vec tmp_y=h*tmp+y1 273db046809SBarry Smith */ 274db046809SBarry Smith 275db046809SBarry Smith tmp_t = t + rk->c[j] * h; 276db046809SBarry Smith 277db046809SBarry Smith /* tmp_y = h * tmp + y1 */ 2782dcb1b2aSMatthew Knepley ierr = VecWAXPY(rk->tmp_y,hh,rk->tmp,rk->y1);CHKERRQ(ierr); 279db046809SBarry Smith 280db046809SBarry Smith /* rk->k[j]=0 */ 281b05257ddSBarry Smith ierr = VecSet(rk->k[j],0.0);CHKERRQ(ierr); 282db046809SBarry Smith ierr = TSComputeRHSFunction(ts,tmp_t,rk->tmp_y,rk->k[j]);CHKERRQ(ierr); 283db046809SBarry Smith } 284db046809SBarry Smith 285db046809SBarry Smith /* tmp=0 and tmp_y=0 */ 286b05257ddSBarry Smith ierr = VecSet(rk->tmp,0.0);CHKERRQ(ierr); 287b05257ddSBarry Smith ierr = VecSet(rk->tmp_y,0.0);CHKERRQ(ierr); 288db046809SBarry Smith 289db046809SBarry Smith for (j = 0; j < rk->s; j++) { 290db046809SBarry Smith /* tmp=b1[j]*k[j]+tmp */ 2912dcb1b2aSMatthew Knepley ierr = VecAXPY(rk->tmp,rk->b1[j],rk->k[j]);CHKERRQ(ierr); 292db046809SBarry Smith /* tmp_y=b2[j]*k[j]+tmp_y */ 2932dcb1b2aSMatthew Knepley ierr = VecAXPY(rk->tmp_y,rk->b2[j],rk->k[j]);CHKERRQ(ierr); 294db046809SBarry Smith } 295db046809SBarry Smith 296db046809SBarry Smith /* y2 = hh * tmp_y */ 297b05257ddSBarry Smith ierr = VecSet(rk->y2,0.0);CHKERRQ(ierr); 2982dcb1b2aSMatthew Knepley ierr = VecAXPY(rk->y2,hh,rk->tmp_y);CHKERRQ(ierr); 299db046809SBarry Smith /* y1 = hh*tmp + y1 */ 3002dcb1b2aSMatthew Knepley ierr = VecAXPY(rk->y1,hh,rk->tmp);CHKERRQ(ierr); 301db046809SBarry Smith /* Finding difference between y1 and y2 */ 302db046809SBarry Smith PetscFunctionReturn(0); 303db046809SBarry Smith } 304db046809SBarry Smith 305e4dd521cSBarry Smith #undef __FUNCT__ 306193ac0bcSJed Brown #define __FUNCT__ "TSSolve_RK" 307193ac0bcSJed Brown static PetscErrorCode TSSolve_RK(TS ts) 308e4dd521cSBarry Smith { 3095f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 310b0dd9724SMatthew Knepley PetscReal norm=0.0,dt_fac=0.0,fac = 0.0 /*,ttmp=0.0*/; 311186e87acSLisandro Dalcin PetscInt i; 312186e87acSLisandro Dalcin PetscErrorCode ierr; 313e4dd521cSBarry Smith 314e4dd521cSBarry Smith PetscFunctionBegin; 315186e87acSLisandro Dalcin ierr = VecCopy(ts->vec_sol,rk->y1);CHKERRQ(ierr); 316186e87acSLisandro Dalcin 317e4dd521cSBarry Smith /* while loop to get from start to stop */ 318186e87acSLisandro Dalcin for (i = 0; i < ts->max_steps; i++) { 3193f2090d5SJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); /* Note that this is called once per STEP, not once per STAGE. */ 320186e87acSLisandro Dalcin 321e4dd521cSBarry Smith /* calling rkqs */ 322e4dd521cSBarry Smith /* 323e4dd521cSBarry Smith -- input 324e4dd521cSBarry Smith ts - pointer to ts 325e4dd521cSBarry Smith ts->ptime - current time 326b05257ddSBarry Smith ts->time_step - try this timestep 327e4dd521cSBarry Smith y1 - solution for this step 328e4dd521cSBarry Smith 329e4dd521cSBarry Smith --output 330e4dd521cSBarry Smith y1 - suggested solution 331e4dd521cSBarry Smith y2 - check solution (runge - kutta second permutation) 332e4dd521cSBarry Smith */ 3335f70b478SJed Brown ierr = TSRKqs(ts,ts->ptime,ts->time_step);CHKERRQ(ierr); 334e3caeda6SBarry Smith /* counting steps */ 335e3caeda6SBarry Smith ts->steps++; 336e4dd521cSBarry Smith /* checking for maxerror */ 337e4dd521cSBarry Smith /* comparing difference to maxerror */ 338e4dd521cSBarry Smith ierr = VecNorm(rk->y2,NORM_2,&norm);CHKERRQ(ierr); 339e4dd521cSBarry Smith /* modifying maxerror to satisfy this timestep */ 340b05257ddSBarry Smith rk->maxerror = rk->ferror * ts->time_step; 341b05257ddSBarry Smith /* ierr = PetscPrintf(PETSC_COMM_WORLD,"norm err: %f maxerror: %f dt: %f",norm,rk->maxerror,ts->time_step);CHKERRQ(ierr); */ 342e4dd521cSBarry Smith 343e4dd521cSBarry Smith /* handling ok and not ok */ 344e4dd521cSBarry Smith if (norm < rk->maxerror) { 345e4dd521cSBarry Smith /* if ok: */ 346e4dd521cSBarry Smith ierr = VecCopy(rk->y1,ts->vec_sol);CHKERRQ(ierr); /* saves the suggested solution to current solution */ 347b05257ddSBarry Smith ts->ptime += ts->time_step; /* storing the new current time */ 348e4dd521cSBarry Smith rk->nok++; 349e4dd521cSBarry Smith fac=5.0; 350e4dd521cSBarry Smith /* trying to save the vector */ 3513f2090d5SJed Brown ierr = TSPostStep(ts);CHKERRQ(ierr); 352e4dd521cSBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 353e3caeda6SBarry Smith if (ts->ptime >= ts->max_time) break; 354e4dd521cSBarry Smith } else { 355e4dd521cSBarry Smith /* if not OK */ 356e4dd521cSBarry Smith rk->nnok++; 357e4dd521cSBarry Smith fac =1.0; 358e4dd521cSBarry Smith ierr=VecCopy(ts->vec_sol,rk->y1);CHKERRQ(ierr); /* restores old solution */ 359e4dd521cSBarry Smith } 360e4dd521cSBarry Smith 361e4dd521cSBarry Smith /*Computing next stepsize. See page 167 in Solving ODE 1 362e4dd521cSBarry Smith * 363e4dd521cSBarry Smith * h_new = h * min(facmax , max(facmin , fac * (tol/err)^(1/(p+1)))) 364e4dd521cSBarry Smith * facmax set above 365e4dd521cSBarry Smith * facmin 366e4dd521cSBarry Smith */ 367e4dd521cSBarry Smith dt_fac = exp(log((rk->maxerror) / norm) / ((rk->p) + 1)) * 0.9; 368e4dd521cSBarry Smith 369bbd56ea5SKarl Rupp if (dt_fac > fac) dt_fac = fac; 370bbd56ea5SKarl Rupp 371e4dd521cSBarry Smith 372b05257ddSBarry Smith /* computing new ts->time_step */ 373b05257ddSBarry Smith ts->time_step = ts->time_step * dt_fac; 374e4dd521cSBarry Smith 375bbd56ea5SKarl Rupp if (ts->ptime+ts->time_step > ts->max_time) ts->time_step = ts->max_time - ts->ptime; 376e4dd521cSBarry Smith 377b05257ddSBarry Smith if (ts->time_step < 1e-14) { 378b05257ddSBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"Very small steps: %f\n",ts->time_step);CHKERRQ(ierr); 379b05257ddSBarry Smith ts->time_step = 1e-14; 380e4dd521cSBarry Smith } 381e4dd521cSBarry Smith 382e4dd521cSBarry Smith /* trying to purify h */ 383e4dd521cSBarry Smith /* (did not give any visible result) */ 384b05257ddSBarry Smith /* ttmp = ts->ptime + ts->time_step; 385b05257ddSBarry Smith ts->time_step = ttmp - ts->ptime; */ 386e4dd521cSBarry Smith 387e4dd521cSBarry Smith } 388e4dd521cSBarry Smith 389e4dd521cSBarry Smith ierr=VecCopy(rk->y1,ts->vec_sol);CHKERRQ(ierr); 390e4dd521cSBarry Smith PetscFunctionReturn(0); 391e4dd521cSBarry Smith } 392e4dd521cSBarry Smith 393e4dd521cSBarry Smith /*------------------------------------------------------------*/ 394e4dd521cSBarry Smith #undef __FUNCT__ 395277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_RK" 396277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts) 397e4dd521cSBarry Smith { 3985f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 3996849ba73SBarry Smith PetscErrorCode ierr; 400e4dd521cSBarry Smith 401e4dd521cSBarry Smith PetscFunctionBegin; 4026bf464f9SBarry Smith ierr = VecDestroy(&rk->y1);CHKERRQ(ierr); 4036bf464f9SBarry Smith ierr = VecDestroy(&rk->y2);CHKERRQ(ierr); 4046bf464f9SBarry Smith ierr = VecDestroy(&rk->tmp);CHKERRQ(ierr); 4056bf464f9SBarry Smith ierr = VecDestroy(&rk->tmp_y);CHKERRQ(ierr); 406277b19d0SLisandro Dalcin if (rk->k) {ierr = VecDestroyVecs(rk->s,&rk->k);CHKERRQ(ierr);} 407277b19d0SLisandro Dalcin PetscFunctionReturn(0); 408e4dd521cSBarry Smith } 409277b19d0SLisandro Dalcin 410277b19d0SLisandro Dalcin #undef __FUNCT__ 411277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_RK" 412277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_RK(TS ts) 413277b19d0SLisandro Dalcin { 414277b19d0SLisandro Dalcin PetscErrorCode ierr; 415277b19d0SLisandro Dalcin 416277b19d0SLisandro Dalcin PetscFunctionBegin; 417277b19d0SLisandro Dalcin ierr = TSReset_RK(ts);CHKERRQ(ierr); 418277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 41900de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetTolerance_C","",NULL);CHKERRQ(ierr); 420e4dd521cSBarry Smith PetscFunctionReturn(0); 421e4dd521cSBarry Smith } 422e4dd521cSBarry Smith /*------------------------------------------------------------*/ 423e4dd521cSBarry Smith 424e4dd521cSBarry Smith #undef __FUNCT__ 4255f70b478SJed Brown #define __FUNCT__ "TSSetFromOptions_RK" 4265f70b478SJed Brown static PetscErrorCode TSSetFromOptions_RK(TS ts) 427e4dd521cSBarry Smith { 4285f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 429dfbe8321SBarry Smith PetscErrorCode ierr; 430262ff7bbSBarry Smith 431e4dd521cSBarry Smith PetscFunctionBegin; 432262ff7bbSBarry Smith ierr = PetscOptionsHead("RK ODE solver options");CHKERRQ(ierr); 4330298fd71SBarry Smith ierr = PetscOptionsReal("-ts_rk_tol","Tolerance for convergence","TSRKSetTolerance",rk->tolerance,&rk->tolerance,NULL);CHKERRQ(ierr); 434262ff7bbSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 435e4dd521cSBarry Smith PetscFunctionReturn(0); 436e4dd521cSBarry Smith } 437e4dd521cSBarry Smith 438e4dd521cSBarry Smith #undef __FUNCT__ 4395f70b478SJed Brown #define __FUNCT__ "TSView_RK" 4405f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer) 441e4dd521cSBarry Smith { 4425f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 4438370ee3bSLisandro Dalcin PetscBool iascii; 444dfbe8321SBarry Smith PetscErrorCode ierr; 4452cdf8120Sasbjorn 446e4dd521cSBarry Smith PetscFunctionBegin; 447251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 4488370ee3bSLisandro Dalcin if (iascii) { 4498370ee3bSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer,"number of ok steps: %D\n",rk->nok);CHKERRQ(ierr); 4508370ee3bSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer,"number of rejected steps: %D\n",rk->nnok);CHKERRQ(ierr); 4518370ee3bSLisandro Dalcin } 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 #undef __FUNCT__ 4715f70b478SJed Brown #define __FUNCT__ "TSCreate_RK" 472*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts) 473e4dd521cSBarry Smith { 4745f70b478SJed Brown TS_RK *rk; 475dfbe8321SBarry Smith PetscErrorCode ierr; 476e4dd521cSBarry Smith 477e4dd521cSBarry Smith PetscFunctionBegin; 4785f70b478SJed Brown ts->ops->setup = TSSetUp_RK; 479193ac0bcSJed Brown ts->ops->solve = TSSolve_RK; 4805f70b478SJed Brown ts->ops->destroy = TSDestroy_RK; 4815f70b478SJed Brown ts->ops->setfromoptions = TSSetFromOptions_RK; 4825f70b478SJed Brown ts->ops->view = TSView_RK; 483e4dd521cSBarry Smith 4845f70b478SJed Brown ierr = PetscNewLog(ts,TS_RK,&rk);CHKERRQ(ierr); 485e4dd521cSBarry Smith ts->data = (void*)rk; 486e4dd521cSBarry Smith 48700de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetTolerance_C","TSRKSetTolerance_RK",TSRKSetTolerance_RK);CHKERRQ(ierr); 488e4dd521cSBarry Smith PetscFunctionReturn(0); 489e4dd521cSBarry Smith } 490