1 2 /* 3 * Code for Timestepping with Runge Kutta 4 * 5 * Written by 6 * Asbjorn Hoiland Aarrestad 7 * asbjorn@aarrestad.com 8 * http://asbjorn.aarrestad.com/ 9 * 10 */ 11 #include "private/tsimpl.h" /*I "petscts.h" I*/ 12 #include "time.h" 13 14 typedef struct { 15 Vec y1,y2; /* work wectors for the two rk permuations */ 16 PetscInt nok,nnok; /* counters for ok and not ok steps */ 17 PetscReal maxerror; /* variable to tell the maxerror allowed */ 18 PetscReal ferror; /* variable to tell (global maxerror)/(total time) */ 19 PetscReal tolerance; /* initial value set for maxerror by user */ 20 Vec tmp,tmp_y,*k; /* two temp vectors and the k vectors for rk */ 21 PetscScalar a[7][6]; /* rk scalars */ 22 PetscScalar b1[7],b2[7]; /* rk scalars */ 23 PetscReal c[7]; /* rk scalars */ 24 PetscInt p,s; /* variables to tell the size of the runge-kutta solver */ 25 } TS_RK; 26 27 EXTERN_C_BEGIN 28 #undef __FUNCT__ 29 #define __FUNCT__ "TSRKSetTolerance_RK" 30 PetscErrorCode TSRKSetTolerance_RK(TS ts,PetscReal aabs) 31 { 32 TS_RK *rk = (TS_RK*)ts->data; 33 34 PetscFunctionBegin; 35 rk->tolerance = aabs; 36 PetscFunctionReturn(0); 37 } 38 EXTERN_C_END 39 40 #undef __FUNCT__ 41 #define __FUNCT__ "TSRKSetTolerance" 42 /*@ 43 TSRKSetTolerance - Sets the total error the RK explicit time integrators 44 will allow over the given time interval. 45 46 Logically Collective on TS 47 48 Input parameters: 49 + ts - the time-step context 50 - aabs - the absolute tolerance 51 52 Level: intermediate 53 54 .keywords: RK, tolerance 55 56 .seealso: TSSundialsSetTolerance() 57 58 @*/ 59 PetscErrorCode TSRKSetTolerance(TS ts,PetscReal aabs) 60 { 61 PetscErrorCode ierr; 62 63 PetscFunctionBegin; 64 PetscValidLogicalCollectiveReal(ts,aabs,2); 65 ierr = PetscTryMethod(ts,"TSRKSetTolerance_C",(TS,PetscReal),(ts,aabs));CHKERRQ(ierr); 66 PetscFunctionReturn(0); 67 } 68 69 70 #undef __FUNCT__ 71 #define __FUNCT__ "TSSetUp_RK" 72 static PetscErrorCode TSSetUp_RK(TS ts) 73 { 74 TS_RK *rk = (TS_RK*)ts->data; 75 PetscErrorCode ierr; 76 77 PetscFunctionBegin; 78 rk->nok = 0; 79 rk->nnok = 0; 80 rk->maxerror = rk->tolerance; 81 82 /* fixing maxerror: global vs local */ 83 rk->ferror = rk->maxerror / (ts->max_time - ts->ptime); 84 85 /* 34.0/45.0 gives double precision division */ 86 /* defining variables needed for Runge-Kutta computing */ 87 /* when changing below, please remember to change a, b1, b2 and c above! */ 88 /* Found in table on page 171: Dormand-Prince 5(4) */ 89 90 /* are these right? */ 91 rk->p=6; 92 rk->s=7; 93 94 rk->a[1][0]=1.0/5.0; 95 rk->a[2][0]=3.0/40.0; 96 rk->a[2][1]=9.0/40.0; 97 rk->a[3][0]=44.0/45.0; 98 rk->a[3][1]=-56.0/15.0; 99 rk->a[3][2]=32.0/9.0; 100 rk->a[4][0]=19372.0/6561.0; 101 rk->a[4][1]=-25360.0/2187.0; 102 rk->a[4][2]=64448.0/6561.0; 103 rk->a[4][3]=-212.0/729.0; 104 rk->a[5][0]=9017.0/3168.0; 105 rk->a[5][1]=-355.0/33.0; 106 rk->a[5][2]=46732.0/5247.0; 107 rk->a[5][3]=49.0/176.0; 108 rk->a[5][4]=-5103.0/18656.0; 109 rk->a[6][0]=35.0/384.0; 110 rk->a[6][1]=0.0; 111 rk->a[6][2]=500.0/1113.0; 112 rk->a[6][3]=125.0/192.0; 113 rk->a[6][4]=-2187.0/6784.0; 114 rk->a[6][5]=11.0/84.0; 115 116 117 rk->c[0]=0.0; 118 rk->c[1]=1.0/5.0; 119 rk->c[2]=3.0/10.0; 120 rk->c[3]=4.0/5.0; 121 rk->c[4]=8.0/9.0; 122 rk->c[5]=1.0; 123 rk->c[6]=1.0; 124 125 rk->b1[0]=35.0/384.0; 126 rk->b1[1]=0.0; 127 rk->b1[2]=500.0/1113.0; 128 rk->b1[3]=125.0/192.0; 129 rk->b1[4]=-2187.0/6784.0; 130 rk->b1[5]=11.0/84.0; 131 rk->b1[6]=0.0; 132 133 rk->b2[0]=5179.0/57600.0; 134 rk->b2[1]=0.0; 135 rk->b2[2]=7571.0/16695.0; 136 rk->b2[3]=393.0/640.0; 137 rk->b2[4]=-92097.0/339200.0; 138 rk->b2[5]=187.0/2100.0; 139 rk->b2[6]=1.0/40.0; 140 141 142 /* Found in table on page 170: Fehlberg 4(5) */ 143 /* 144 rk->p=5; 145 rk->s=6; 146 147 rk->a[1][0]=1.0/4.0; 148 rk->a[2][0]=3.0/32.0; 149 rk->a[2][1]=9.0/32.0; 150 rk->a[3][0]=1932.0/2197.0; 151 rk->a[3][1]=-7200.0/2197.0; 152 rk->a[3][2]=7296.0/2197.0; 153 rk->a[4][0]=439.0/216.0; 154 rk->a[4][1]=-8.0; 155 rk->a[4][2]=3680.0/513.0; 156 rk->a[4][3]=-845.0/4104.0; 157 rk->a[5][0]=-8.0/27.0; 158 rk->a[5][1]=2.0; 159 rk->a[5][2]=-3544.0/2565.0; 160 rk->a[5][3]=1859.0/4104.0; 161 rk->a[5][4]=-11.0/40.0; 162 163 rk->c[0]=0.0; 164 rk->c[1]=1.0/4.0; 165 rk->c[2]=3.0/8.0; 166 rk->c[3]=12.0/13.0; 167 rk->c[4]=1.0; 168 rk->c[5]=1.0/2.0; 169 170 rk->b1[0]=25.0/216.0; 171 rk->b1[1]=0.0; 172 rk->b1[2]=1408.0/2565.0; 173 rk->b1[3]=2197.0/4104.0; 174 rk->b1[4]=-1.0/5.0; 175 rk->b1[5]=0.0; 176 177 rk->b2[0]=16.0/135.0; 178 rk->b2[1]=0.0; 179 rk->b2[2]=6656.0/12825.0; 180 rk->b2[3]=28561.0/56430.0; 181 rk->b2[4]=-9.0/50.0; 182 rk->b2[5]=2.0/55.0; 183 */ 184 /* Found in table on page 169: Merson 4("5") */ 185 /* 186 rk->p=4; 187 rk->s=5; 188 rk->a[1][0] = 1.0/3.0; 189 rk->a[2][0] = 1.0/6.0; 190 rk->a[2][1] = 1.0/6.0; 191 rk->a[3][0] = 1.0/8.0; 192 rk->a[3][1] = 0.0; 193 rk->a[3][2] = 3.0/8.0; 194 rk->a[4][0] = 1.0/2.0; 195 rk->a[4][1] = 0.0; 196 rk->a[4][2] = -3.0/2.0; 197 rk->a[4][3] = 2.0; 198 199 rk->c[0] = 0.0; 200 rk->c[1] = 1.0/3.0; 201 rk->c[2] = 1.0/3.0; 202 rk->c[3] = 0.5; 203 rk->c[4] = 1.0; 204 205 rk->b1[0] = 1.0/2.0; 206 rk->b1[1] = 0.0; 207 rk->b1[2] = -3.0/2.0; 208 rk->b1[3] = 2.0; 209 rk->b1[4] = 0.0; 210 211 rk->b2[0] = 1.0/6.0; 212 rk->b2[1] = 0.0; 213 rk->b2[2] = 0.0; 214 rk->b2[3] = 2.0/3.0; 215 rk->b2[4] = 1.0/6.0; 216 */ 217 218 /* making b2 -> e=b1-b2 */ 219 /* 220 for(i=0;i<rk->s;i++){ 221 rk->b2[i] = (rk->b1[i]) - (rk->b2[i]); 222 } 223 */ 224 rk->b2[0]=71.0/57600.0; 225 rk->b2[1]=0.0; 226 rk->b2[2]=-71.0/16695.0; 227 rk->b2[3]=71.0/1920.0; 228 rk->b2[4]=-17253.0/339200.0; 229 rk->b2[5]=22.0/525.0; 230 rk->b2[6]=-1.0/40.0; 231 232 /* initializing vectors */ 233 ierr = VecDuplicate(ts->vec_sol,&rk->y1);CHKERRQ(ierr); 234 ierr = VecDuplicate(ts->vec_sol,&rk->y2);CHKERRQ(ierr); 235 ierr = VecDuplicate(rk->y1,&rk->tmp);CHKERRQ(ierr); 236 ierr = VecDuplicate(rk->y1,&rk->tmp_y);CHKERRQ(ierr); 237 ierr = VecDuplicateVecs(rk->y1,rk->s,&rk->k);CHKERRQ(ierr); 238 239 PetscFunctionReturn(0); 240 } 241 242 /*------------------------------------------------------------*/ 243 #undef __FUNCT__ 244 #define __FUNCT__ "TSRKqs" 245 PetscErrorCode TSRKqs(TS ts,PetscReal t,PetscReal h) 246 { 247 TS_RK *rk = (TS_RK*)ts->data; 248 PetscErrorCode ierr; 249 PetscInt j,l; 250 PetscReal tmp_t=t; 251 PetscScalar hh=h; 252 253 PetscFunctionBegin; 254 /* k[0]=0 */ 255 ierr = VecSet(rk->k[0],0.0);CHKERRQ(ierr); 256 257 /* k[0] = derivs(t,y1) */ 258 ierr = TSComputeRHSFunction(ts,t,rk->y1,rk->k[0]);CHKERRQ(ierr); 259 /* looping over runge-kutta variables */ 260 /* building the k - array of vectors */ 261 for(j = 1 ; j < rk->s ; j++){ 262 263 /* rk->tmp = 0 */ 264 ierr = VecSet(rk->tmp,0.0);CHKERRQ(ierr); 265 266 for(l=0;l<j;l++){ 267 /* tmp += a(j,l)*k[l] */ 268 ierr = VecAXPY(rk->tmp,rk->a[j][l],rk->k[l]);CHKERRQ(ierr); 269 } 270 271 /* ierr = VecView(rk->tmp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 272 273 /* k[j] = derivs(t+c(j)*h,y1+h*tmp,k(j)) */ 274 /* I need the following helpers: 275 PetscScalar tmp_t=t+c(j)*h 276 Vec tmp_y=h*tmp+y1 277 */ 278 279 tmp_t = t + rk->c[j] * h; 280 281 /* tmp_y = h * tmp + y1 */ 282 ierr = VecWAXPY(rk->tmp_y,hh,rk->tmp,rk->y1);CHKERRQ(ierr); 283 284 /* rk->k[j]=0 */ 285 ierr = VecSet(rk->k[j],0.0);CHKERRQ(ierr); 286 ierr = TSComputeRHSFunction(ts,tmp_t,rk->tmp_y,rk->k[j]);CHKERRQ(ierr); 287 } 288 289 /* tmp=0 and tmp_y=0 */ 290 ierr = VecSet(rk->tmp,0.0);CHKERRQ(ierr); 291 ierr = VecSet(rk->tmp_y,0.0);CHKERRQ(ierr); 292 293 for(j = 0 ; j < rk->s ; j++){ 294 /* tmp=b1[j]*k[j]+tmp */ 295 ierr = VecAXPY(rk->tmp,rk->b1[j],rk->k[j]);CHKERRQ(ierr); 296 /* tmp_y=b2[j]*k[j]+tmp_y */ 297 ierr = VecAXPY(rk->tmp_y,rk->b2[j],rk->k[j]);CHKERRQ(ierr); 298 } 299 300 /* y2 = hh * tmp_y */ 301 ierr = VecSet(rk->y2,0.0);CHKERRQ(ierr); 302 ierr = VecAXPY(rk->y2,hh,rk->tmp_y);CHKERRQ(ierr); 303 /* y1 = hh*tmp + y1 */ 304 ierr = VecAXPY(rk->y1,hh,rk->tmp);CHKERRQ(ierr); 305 /* Finding difference between y1 and y2 */ 306 PetscFunctionReturn(0); 307 } 308 309 #undef __FUNCT__ 310 #define __FUNCT__ "TSStep_RK" 311 static PetscErrorCode TSStep_RK(TS ts,PetscInt *steps,PetscReal *ptime) 312 { 313 TS_RK *rk = (TS_RK*)ts->data; 314 PetscErrorCode ierr; 315 PetscReal norm=0.0,dt_fac=0.0,fac = 0.0/*,ttmp=0.0*/; 316 PetscInt i, max_steps = ts->max_steps; 317 318 PetscFunctionBegin; 319 ierr=VecCopy(ts->vec_sol,rk->y1);CHKERRQ(ierr); 320 *steps = -ts->steps; 321 /* trying to save the vector */ 322 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 323 /* while loop to get from start to stop */ 324 for (i = 0; i < max_steps; i++) { 325 ierr = TSPreStep(ts);CHKERRQ(ierr); /* Note that this is called once per STEP, not once per STAGE. */ 326 /* calling rkqs */ 327 /* 328 -- input 329 ts - pointer to ts 330 ts->ptime - current time 331 ts->time_step - try this timestep 332 y1 - solution for this step 333 334 --output 335 y1 - suggested solution 336 y2 - check solution (runge - kutta second permutation) 337 */ 338 ierr = TSRKqs(ts,ts->ptime,ts->time_step);CHKERRQ(ierr); 339 /* counting steps */ 340 ts->steps++; 341 /* checking for maxerror */ 342 /* comparing difference to maxerror */ 343 ierr = VecNorm(rk->y2,NORM_2,&norm);CHKERRQ(ierr); 344 /* modifying maxerror to satisfy this timestep */ 345 rk->maxerror = rk->ferror * ts->time_step; 346 /* ierr = PetscPrintf(PETSC_COMM_WORLD,"norm err: %f maxerror: %f dt: %f",norm,rk->maxerror,ts->time_step);CHKERRQ(ierr); */ 347 348 /* handling ok and not ok */ 349 if (norm < rk->maxerror){ 350 /* if ok: */ 351 ierr=VecCopy(rk->y1,ts->vec_sol);CHKERRQ(ierr); /* saves the suggested solution to current solution */ 352 ts->ptime += ts->time_step; /* storing the new current time */ 353 rk->nok++; 354 fac=5.0; 355 /* trying to save the vector */ 356 ierr = TSPostStep(ts);CHKERRQ(ierr); 357 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 358 if (ts->ptime >= ts->max_time) break; 359 } else{ 360 /* if not OK */ 361 rk->nnok++; 362 fac=1.0; 363 ierr=VecCopy(ts->vec_sol,rk->y1);CHKERRQ(ierr); /* restores old solution */ 364 } 365 366 /*Computing next stepsize. See page 167 in Solving ODE 1 367 * 368 * h_new = h * min( facmax , max( facmin , fac * (tol/err)^(1/(p+1)) ) ) 369 * facmax set above 370 * facmin 371 */ 372 dt_fac = exp(log((rk->maxerror) / norm) / ((rk->p) + 1) ) * 0.9 ; 373 374 if (dt_fac > fac){ 375 /*ierr = PetscPrintf(PETSC_COMM_WORLD,"changing fac %f\n",fac);*/ 376 dt_fac = fac; 377 } 378 379 /* computing new ts->time_step */ 380 ts->time_step = ts->time_step * dt_fac; 381 382 if (ts->ptime+ts->time_step > ts->max_time){ 383 ts->time_step = ts->max_time - ts->ptime; 384 } 385 386 if (ts->time_step < 1e-14){ 387 ierr = PetscPrintf(PETSC_COMM_WORLD,"Very small steps: %f\n",ts->time_step);CHKERRQ(ierr); 388 ts->time_step = 1e-14; 389 } 390 391 /* trying to purify h */ 392 /* (did not give any visible result) */ 393 /* ttmp = ts->ptime + ts->time_step; 394 ts->time_step = ttmp - ts->ptime; */ 395 396 } 397 398 ierr=VecCopy(rk->y1,ts->vec_sol);CHKERRQ(ierr); 399 *steps += ts->steps; 400 *ptime = ts->ptime; 401 PetscFunctionReturn(0); 402 } 403 404 /*------------------------------------------------------------*/ 405 #undef __FUNCT__ 406 #define __FUNCT__ "TSDestroy_RK" 407 static PetscErrorCode TSDestroy_RK(TS ts) 408 { 409 TS_RK *rk = (TS_RK*)ts->data; 410 PetscErrorCode ierr; 411 PetscInt i; 412 413 /* REMEMBER TO DESTROY ALL */ 414 415 PetscFunctionBegin; 416 if (rk->y1) {ierr = VecDestroy(rk->y1);CHKERRQ(ierr);} 417 if (rk->y2) {ierr = VecDestroy(rk->y2);CHKERRQ(ierr);} 418 if (rk->tmp) {ierr = VecDestroy(rk->tmp);CHKERRQ(ierr);} 419 if (rk->tmp_y) {ierr = VecDestroy(rk->tmp_y);CHKERRQ(ierr);} 420 for(i=0;i<rk->s;i++){ 421 if (rk->k[i]) {ierr = VecDestroy(rk->k[i]);CHKERRQ(ierr);} 422 } 423 ierr = PetscFree(rk);CHKERRQ(ierr); 424 PetscFunctionReturn(0); 425 } 426 /*------------------------------------------------------------*/ 427 428 #undef __FUNCT__ 429 #define __FUNCT__ "TSSetFromOptions_RK" 430 static PetscErrorCode TSSetFromOptions_RK(TS ts) 431 { 432 TS_RK *rk = (TS_RK*)ts->data; 433 PetscErrorCode ierr; 434 435 PetscFunctionBegin; 436 ierr = PetscOptionsHead("RK ODE solver options");CHKERRQ(ierr); 437 ierr = PetscOptionsReal("-ts_rk_tol","Tolerance for convergence","TSRKSetTolerance",rk->tolerance,&rk->tolerance,PETSC_NULL);CHKERRQ(ierr); 438 ierr = PetscOptionsTail();CHKERRQ(ierr); 439 PetscFunctionReturn(0); 440 } 441 442 #undef __FUNCT__ 443 #define __FUNCT__ "TSView_RK" 444 static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer) 445 { 446 TS_RK *rk = (TS_RK*)ts->data; 447 PetscErrorCode ierr; 448 449 PetscFunctionBegin; 450 ierr = PetscPrintf(PETSC_COMM_WORLD," number of ok steps: %D\n",rk->nok);CHKERRQ(ierr); 451 ierr = PetscPrintf(PETSC_COMM_WORLD," number of rejected steps: %D\n",rk->nnok);CHKERRQ(ierr); 452 PetscFunctionReturn(0); 453 } 454 455 /* ------------------------------------------------------------ */ 456 /*MC 457 TSRK - ODE solver using the explicit Runge-Kutta methods 458 459 Options Database: 460 . -ts_rk_tol <tol> Tolerance for convergence 461 462 Contributed by: Asbjorn Hoiland Aarrestad, asbjorn@aarrestad.com, http://asbjorn.aarrestad.com/ 463 464 Level: beginner 465 466 .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSRKSetTolerance() 467 468 M*/ 469 470 EXTERN_C_BEGIN 471 #undef __FUNCT__ 472 #define __FUNCT__ "TSCreate_RK" 473 PetscErrorCode TSCreate_RK(TS ts) 474 { 475 TS_RK *rk; 476 PetscErrorCode ierr; 477 478 PetscFunctionBegin; 479 ts->ops->setup = TSSetUp_RK; 480 ts->ops->step = TSStep_RK; 481 ts->ops->destroy = TSDestroy_RK; 482 ts->ops->setfromoptions = TSSetFromOptions_RK; 483 ts->ops->view = TSView_RK; 484 485 ierr = PetscNewLog(ts,TS_RK,&rk);CHKERRQ(ierr); 486 ts->data = (void*)rk; 487 488 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRKSetTolerance_C","TSRKSetTolerance_RK",TSRKSetTolerance_RK);CHKERRQ(ierr); 489 490 PetscFunctionReturn(0); 491 } 492 EXTERN_C_END 493 494 495 496 497