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 "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 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 ierr = PetscObjectQueryFunction((PetscObject)ts,"TSRKSetTolerance_C",(void (**)(void))&f);CHKERRQ(ierr); 66 if (f) { 67 ierr = (*f)(ts,aabs);CHKERRQ(ierr); 68 } 69 PetscFunctionReturn(0); 70 } 71 72 73 #undef __FUNCT__ 74 #define __FUNCT__ "TSSetUp_Rk" 75 static PetscErrorCode TSSetUp_Rk(TS ts) 76 { 77 TS_Rk *rk = (TS_Rk*)ts->data; 78 PetscErrorCode ierr; 79 80 PetscFunctionBegin; 81 rk->nok = 0; 82 rk->nnok = 0; 83 rk->maxerror = rk->tolerance; 84 85 /* fixing maxerror: global vs local */ 86 rk->ferror = rk->maxerror / (ts->max_time - ts->ptime); 87 88 /* 34.0/45.0 gives double precision division */ 89 /* defining variables needed for Runge-Kutta computing */ 90 /* when changing below, please remember to change a, b1, b2 and c above! */ 91 /* Found in table on page 171: Dormand-Prince 5(4) */ 92 93 /* are these right? */ 94 rk->p=6; 95 rk->s=7; 96 97 rk->a[1][0]=1.0/5.0; 98 rk->a[2][0]=3.0/40.0; 99 rk->a[2][1]=9.0/40.0; 100 rk->a[3][0]=44.0/45.0; 101 rk->a[3][1]=-56.0/15.0; 102 rk->a[3][2]=32.0/9.0; 103 rk->a[4][0]=19372.0/6561.0; 104 rk->a[4][1]=-25360.0/2187.0; 105 rk->a[4][2]=64448.0/6561.0; 106 rk->a[4][3]=-212.0/729.0; 107 rk->a[5][0]=9017.0/3168.0; 108 rk->a[5][1]=-355.0/33.0; 109 rk->a[5][2]=46732.0/5247.0; 110 rk->a[5][3]=49.0/176.0; 111 rk->a[5][4]=-5103.0/18656.0; 112 rk->a[6][0]=35.0/384.0; 113 rk->a[6][1]=0.0; 114 rk->a[6][2]=500.0/1113.0; 115 rk->a[6][3]=125.0/192.0; 116 rk->a[6][4]=-2187.0/6784.0; 117 rk->a[6][5]=11.0/84.0; 118 119 120 rk->c[0]=0.0; 121 rk->c[1]=1.0/5.0; 122 rk->c[2]=3.0/10.0; 123 rk->c[3]=4.0/5.0; 124 rk->c[4]=8.0/9.0; 125 rk->c[5]=1.0; 126 rk->c[6]=1.0; 127 128 rk->b1[0]=35.0/384.0; 129 rk->b1[1]=0.0; 130 rk->b1[2]=500.0/1113.0; 131 rk->b1[3]=125.0/192.0; 132 rk->b1[4]=-2187.0/6784.0; 133 rk->b1[5]=11.0/84.0; 134 rk->b1[6]=0.0; 135 136 rk->b2[0]=5179.0/57600.0; 137 rk->b2[1]=0.0; 138 rk->b2[2]=7571.0/16695.0; 139 rk->b2[3]=393.0/640.0; 140 rk->b2[4]=-92097.0/339200.0; 141 rk->b2[5]=187.0/2100.0; 142 rk->b2[6]=1.0/40.0; 143 144 145 /* Found in table on page 170: Fehlberg 4(5) */ 146 /* 147 rk->p=5; 148 rk->s=6; 149 150 rk->a[1][0]=1.0/4.0; 151 rk->a[2][0]=3.0/32.0; 152 rk->a[2][1]=9.0/32.0; 153 rk->a[3][0]=1932.0/2197.0; 154 rk->a[3][1]=-7200.0/2197.0; 155 rk->a[3][2]=7296.0/2197.0; 156 rk->a[4][0]=439.0/216.0; 157 rk->a[4][1]=-8.0; 158 rk->a[4][2]=3680.0/513.0; 159 rk->a[4][3]=-845.0/4104.0; 160 rk->a[5][0]=-8.0/27.0; 161 rk->a[5][1]=2.0; 162 rk->a[5][2]=-3544.0/2565.0; 163 rk->a[5][3]=1859.0/4104.0; 164 rk->a[5][4]=-11.0/40.0; 165 166 rk->c[0]=0.0; 167 rk->c[1]=1.0/4.0; 168 rk->c[2]=3.0/8.0; 169 rk->c[3]=12.0/13.0; 170 rk->c[4]=1.0; 171 rk->c[5]=1.0/2.0; 172 173 rk->b1[0]=25.0/216.0; 174 rk->b1[1]=0.0; 175 rk->b1[2]=1408.0/2565.0; 176 rk->b1[3]=2197.0/4104.0; 177 rk->b1[4]=-1.0/5.0; 178 rk->b1[5]=0.0; 179 180 rk->b2[0]=16.0/135.0; 181 rk->b2[1]=0.0; 182 rk->b2[2]=6656.0/12825.0; 183 rk->b2[3]=28561.0/56430.0; 184 rk->b2[4]=-9.0/50.0; 185 rk->b2[5]=2.0/55.0; 186 */ 187 /* Found in table on page 169: Merson 4("5") */ 188 /* 189 rk->p=4; 190 rk->s=5; 191 rk->a[1][0] = 1.0/3.0; 192 rk->a[2][0] = 1.0/6.0; 193 rk->a[2][1] = 1.0/6.0; 194 rk->a[3][0] = 1.0/8.0; 195 rk->a[3][1] = 0.0; 196 rk->a[3][2] = 3.0/8.0; 197 rk->a[4][0] = 1.0/2.0; 198 rk->a[4][1] = 0.0; 199 rk->a[4][2] = -3.0/2.0; 200 rk->a[4][3] = 2.0; 201 202 rk->c[0] = 0.0; 203 rk->c[1] = 1.0/3.0; 204 rk->c[2] = 1.0/3.0; 205 rk->c[3] = 0.5; 206 rk->c[4] = 1.0; 207 208 rk->b1[0] = 1.0/2.0; 209 rk->b1[1] = 0.0; 210 rk->b1[2] = -3.0/2.0; 211 rk->b1[3] = 2.0; 212 rk->b1[4] = 0.0; 213 214 rk->b2[0] = 1.0/6.0; 215 rk->b2[1] = 0.0; 216 rk->b2[2] = 0.0; 217 rk->b2[3] = 2.0/3.0; 218 rk->b2[4] = 1.0/6.0; 219 */ 220 221 /* making b2 -> e=b1-b2 */ 222 /* 223 for(i=0;i<rk->s;i++){ 224 rk->b2[i] = (rk->b1[i]) - (rk->b2[i]); 225 } 226 */ 227 rk->b2[0]=71.0/57600.0; 228 rk->b2[1]=0.0; 229 rk->b2[2]=-71.0/16695.0; 230 rk->b2[3]=71.0/1920.0; 231 rk->b2[4]=-17253.0/339200.0; 232 rk->b2[5]=22.0/525.0; 233 rk->b2[6]=-1.0/40.0; 234 235 /* initializing vectors */ 236 ierr = VecDuplicate(ts->vec_sol,&rk->y1);CHKERRQ(ierr); 237 ierr = VecDuplicate(ts->vec_sol,&rk->y2);CHKERRQ(ierr); 238 ierr = VecDuplicate(rk->y1,&rk->tmp);CHKERRQ(ierr); 239 ierr = VecDuplicate(rk->y1,&rk->tmp_y);CHKERRQ(ierr); 240 ierr = VecDuplicateVecs(rk->y1,rk->s,&rk->k);CHKERRQ(ierr); 241 242 PetscFunctionReturn(0); 243 } 244 245 /*------------------------------------------------------------*/ 246 #undef __FUNCT__ 247 #define __FUNCT__ "TSRkqs" 248 PetscErrorCode TSRkqs(TS ts,PetscReal t,PetscReal h) 249 { 250 TS_Rk *rk = (TS_Rk*)ts->data; 251 PetscErrorCode ierr; 252 PetscInt j,l; 253 PetscReal tmp_t=t; 254 PetscScalar hh=h; 255 256 PetscFunctionBegin; 257 /* k[0]=0 */ 258 ierr = VecSet(rk->k[0],0.0);CHKERRQ(ierr); 259 260 /* k[0] = derivs(t,y1) */ 261 ierr = TSComputeRHSFunction(ts,t,rk->y1,rk->k[0]);CHKERRQ(ierr); 262 /* looping over runge-kutta variables */ 263 /* building the k - array of vectors */ 264 for(j = 1 ; j < rk->s ; j++){ 265 266 /* rk->tmp = 0 */ 267 ierr = VecSet(rk->tmp,0.0);CHKERRQ(ierr); 268 269 for(l=0;l<j;l++){ 270 /* tmp += a(j,l)*k[l] */ 271 ierr = VecAXPY(rk->tmp,rk->a[j][l],rk->k[l]);CHKERRQ(ierr); 272 } 273 274 /* ierr = VecView(rk->tmp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 275 276 /* k[j] = derivs(t+c(j)*h,y1+h*tmp,k(j)) */ 277 /* I need the following helpers: 278 PetscScalar tmp_t=t+c(j)*h 279 Vec tmp_y=h*tmp+y1 280 */ 281 282 tmp_t = t + rk->c[j] * h; 283 284 /* tmp_y = h * tmp + y1 */ 285 ierr = VecWAXPY(rk->tmp_y,hh,rk->tmp,rk->y1);CHKERRQ(ierr); 286 287 /* rk->k[j]=0 */ 288 ierr = VecSet(rk->k[j],0.0);CHKERRQ(ierr); 289 ierr = TSComputeRHSFunction(ts,tmp_t,rk->tmp_y,rk->k[j]);CHKERRQ(ierr); 290 } 291 292 /* tmp=0 and tmp_y=0 */ 293 ierr = VecSet(rk->tmp,0.0);CHKERRQ(ierr); 294 ierr = VecSet(rk->tmp_y,0.0);CHKERRQ(ierr); 295 296 for(j = 0 ; j < rk->s ; j++){ 297 /* tmp=b1[j]*k[j]+tmp */ 298 ierr = VecAXPY(rk->tmp,rk->b1[j],rk->k[j]);CHKERRQ(ierr); 299 /* tmp_y=b2[j]*k[j]+tmp_y */ 300 ierr = VecAXPY(rk->tmp_y,rk->b2[j],rk->k[j]);CHKERRQ(ierr); 301 } 302 303 /* y2 = hh * tmp_y */ 304 ierr = VecSet(rk->y2,0.0);CHKERRQ(ierr); 305 ierr = VecAXPY(rk->y2,hh,rk->tmp_y);CHKERRQ(ierr); 306 /* y1 = hh*tmp + y1 */ 307 ierr = VecAXPY(rk->y1,hh,rk->tmp);CHKERRQ(ierr); 308 /* Finding difference between y1 and y2 */ 309 PetscFunctionReturn(0); 310 } 311 312 #undef __FUNCT__ 313 #define __FUNCT__ "TSStep_Rk" 314 static PetscErrorCode TSStep_Rk(TS ts,PetscInt *steps,PetscReal *ptime) 315 { 316 TS_Rk *rk = (TS_Rk*)ts->data; 317 PetscErrorCode ierr; 318 PetscReal norm=0.0,dt_fac=0.0,fac = 0.0/*,ttmp=0.0*/; 319 320 PetscFunctionBegin; 321 ierr=VecCopy(ts->vec_sol,rk->y1);CHKERRQ(ierr); 322 *steps = -ts->steps; 323 /* trying to save the vector */ 324 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 325 /* while loop to get from start to stop */ 326 while (ts->ptime < ts->max_time){ 327 /* calling rkqs */ 328 /* 329 -- input 330 ts - pointer to ts 331 ts->ptime - current time 332 ts->time_step - try this timestep 333 y1 - solution for this step 334 335 --output 336 y1 - suggested solution 337 y2 - check solution (runge - kutta second permutation) 338 */ 339 ierr = TSRkqs(ts,ts->ptime,ts->time_step);CHKERRQ(ierr); 340 /* checking for maxerror */ 341 /* comparing difference to maxerror */ 342 ierr = VecNorm(rk->y2,NORM_2,&norm);CHKERRQ(ierr); 343 /* modifying maxerror to satisfy this timestep */ 344 rk->maxerror = rk->ferror * ts->time_step; 345 /* ierr = PetscPrintf(PETSC_COMM_WORLD,"norm err: %f maxerror: %f dt: %f",norm,rk->maxerror,ts->time_step);CHKERRQ(ierr); */ 346 347 /* handling ok and not ok */ 348 if (norm < rk->maxerror){ 349 /* if ok: */ 350 ierr=VecCopy(rk->y1,ts->vec_sol);CHKERRQ(ierr); /* saves the suggested solution to current solution */ 351 ts->ptime += ts->time_step; /* storing the new current time */ 352 rk->nok++; 353 fac=5.0; 354 /* trying to save the vector */ 355 /* calling monitor */ 356 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 357 } else{ 358 /* if not OK */ 359 rk->nnok++; 360 fac=1.0; 361 ierr=VecCopy(ts->vec_sol,rk->y1);CHKERRQ(ierr); /* restores old solution */ 362 } 363 364 /*Computing next stepsize. See page 167 in Solving ODE 1 365 * 366 * h_new = h * min( facmax , max( facmin , fac * (tol/err)^(1/(p+1)) ) ) 367 * facmax set above 368 * facmin 369 */ 370 dt_fac = exp(log((rk->maxerror) / norm) / ((rk->p) + 1) ) * 0.9 ; 371 372 if (dt_fac > fac){ 373 /*ierr = PetscPrintf(PETSC_COMM_WORLD,"changing fac %f\n",fac);*/ 374 dt_fac = fac; 375 } 376 377 /* computing new ts->time_step */ 378 ts->time_step = ts->time_step * dt_fac; 379 380 if (ts->ptime+ts->time_step > ts->max_time){ 381 ts->time_step = ts->max_time - ts->ptime; 382 } 383 384 if (ts->time_step < 1e-14){ 385 ierr = PetscPrintf(PETSC_COMM_WORLD,"Very small steps: %f\n",ts->time_step);CHKERRQ(ierr); 386 ts->time_step = 1e-14; 387 } 388 389 /* trying to purify h */ 390 /* (did not give any visible result) */ 391 /* ttmp = ts->ptime + ts->time_step; 392 ts->time_step = ttmp - ts->ptime; */ 393 394 /* counting steps */ 395 ts->steps++; 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 TS_RK - 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(), TS_EULER, TSRKSetTolerance() 467 468 M*/ 469 470 EXTERN_C_BEGIN 471 #undef __FUNCT__ 472 #define __FUNCT__ "TSCreate_Rk" 473 PetscErrorCode PETSCTS_DLLEXPORT 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