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