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