1 /* 2 Code for time stepping with the Runge-Kutta method 3 4 Notes: 5 The general system is written as 6 7 Udot = F(t,U) 8 9 */ 10 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 11 #include <petscdm.h> 12 13 static TSRKType TSRKDefault = TSRK3BS; 14 static PetscBool TSRKRegisterAllCalled; 15 static PetscBool TSRKPackageInitialized; 16 17 typedef struct _RKTableau *RKTableau; 18 struct _RKTableau { 19 char *name; 20 PetscInt order; /* Classical approximation order of the method i */ 21 PetscInt s; /* Number of stages */ 22 PetscInt p; /* Interpolation order */ 23 PetscBool FSAL; /* flag to indicate if tableau is FSAL */ 24 PetscReal *A,*b,*c; /* Tableau */ 25 PetscReal *bembed; /* Embedded formula of order one less (order-1) */ 26 PetscReal *binterp; /* Dense output formula */ 27 PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */ 28 }; 29 typedef struct _RKTableauLink *RKTableauLink; 30 struct _RKTableauLink { 31 struct _RKTableau tab; 32 RKTableauLink next; 33 }; 34 static RKTableauLink RKTableauList; 35 36 typedef struct { 37 RKTableau tableau; 38 Vec *Y; /* States computed during the step */ 39 Vec *YdotRHS; /* Function evaluations for the non-stiff part */ 40 Vec *VecDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage */ 41 Vec *VecDeltaMu; /* Increment of the adjoint sensitivity w.r.t P at stage */ 42 Vec *VecSensiTemp; /* Vector to be timed with Jacobian transpose */ 43 Vec VecCostIntegral0; /* backup for roll-backs due to events */ 44 PetscScalar *work; /* Scalar work */ 45 PetscReal stage_time; 46 TSStepStatus status; 47 PetscReal ptime; 48 PetscReal time_step; 49 } TS_RK; 50 51 /*MC 52 TSRK1FE - First order forward Euler scheme. 53 54 This method has one stage. 55 56 Options database: 57 . -ts_tk_type 1fe 58 59 Level: advanced 60 61 .seealso: TSRK, TSRKType, TSRKSetType() 62 M*/ 63 /*MC 64 TSRK2A - Second order RK scheme. 65 66 This method has two stages. 67 68 Options database: 69 . -ts_tk_type 2a 70 71 Level: advanced 72 73 .seealso: TSRK, TSRKType, TSRKSetType() 74 M*/ 75 /*MC 76 TSRK3 - Third order RK scheme. 77 78 This method has three stages. 79 80 Options database: 81 . -ts_tk_type 3 82 83 Level: advanced 84 85 .seealso: TSRK, TSRKType, TSRKSetType() 86 M*/ 87 /*MC 88 TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method. 89 90 This method has four stages with the First Same As Last (FSAL) property. 91 92 Options database: 93 . -ts_tk_type 3bs 94 95 Level: advanced 96 97 References: https://doi.org/10.1016/0893-9659(89)90079-7 98 99 .seealso: TSRK, TSRKType, TSRKSetType() 100 M*/ 101 /*MC 102 TSRK4 - Fourth order RK scheme. 103 104 This is the classical Runge-Kutta method with four stages. 105 106 Options database: 107 . -ts_tk_type 4 108 109 Level: advanced 110 111 .seealso: TSRK, TSRKType, TSRKSetType() 112 M*/ 113 /*MC 114 TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method. 115 116 This method has six stages. 117 118 Options database: 119 . -ts_tk_type 5f 120 121 Level: advanced 122 123 .seealso: TSRK, TSRKType, TSRKSetType() 124 M*/ 125 /*MC 126 TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method. 127 128 This method has seven stages with the First Same As Last (FSAL) property. 129 130 Options database: 131 . -ts_tk_type 5dp 132 133 Level: advanced 134 135 References: https://doi.org/10.1016/0771-050X(80)90013-3 136 137 .seealso: TSRK, TSRKType, TSRKSetType() 138 M*/ 139 /*MC 140 TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method. 141 142 This method has eight stages with the First Same As Last (FSAL) property. 143 144 Options database: 145 . -ts_tk_type 5bs 146 147 Level: advanced 148 149 References: https://doi.org/10.1016/0898-1221(96)00141-1 150 151 .seealso: TSRK, TSRKType, TSRKSetType() 152 M*/ 153 154 /*@C 155 TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK 156 157 Not Collective, but should be called by all processes which will need the schemes to be registered 158 159 Level: advanced 160 161 .keywords: TS, TSRK, register, all 162 163 .seealso: TSRKRegisterDestroy() 164 @*/ 165 PetscErrorCode TSRKRegisterAll(void) 166 { 167 PetscErrorCode ierr; 168 169 PetscFunctionBegin; 170 if (TSRKRegisterAllCalled) PetscFunctionReturn(0); 171 TSRKRegisterAllCalled = PETSC_TRUE; 172 173 { 174 const PetscReal 175 A[1][1] = {{0.0}}, 176 b[1] = {1.0}; 177 ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 178 } 179 { 180 const PetscReal 181 A[2][2] = {{0.0,0.0}, 182 {1.0,0.0}}, 183 b[2] = {0.5,0.5}, 184 bembed[2] = {1.0,0}; 185 ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 186 } 187 { 188 const PetscReal 189 #if defined(PETSC_USE_REAL___FLOAT128) 190 A[3][3] = {{0,0,0}, 191 {2.0q/3.0q,0,0}, 192 {-1.0q/3.0q,1.0,0}}, 193 #else 194 A[3][3] = {{0,0,0}, 195 {2.0/3.0,0,0}, 196 {-1.0/3.0,1.0,0}}, 197 #endif 198 b[3] = {0.25,0.5,0.25}; 199 ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 200 } 201 { 202 const PetscReal 203 #if defined(PETSC_USE_REAL___FLOAT128) 204 A[4][4] = {{0,0,0,0}, 205 {1.0q/2.0q,0,0,0}, 206 {0,3.0q/4.0q,0,0}, 207 {2.0q/9.0q,1.0q/3.0q,4.0q/9.0q,0}}, 208 b[4] = {2.0q/9.0q,1.0q/3.0q,4.0q/9.0q,0}, 209 bembed[4] = {7.0q/24.0q,1.0q/4.0q,1.0q/3.0q,1.0q/8.0q}; 210 #else 211 A[4][4] = {{0,0,0,0}, 212 {1.0/2.0,0,0,0}, 213 {0,3.0/4.0,0,0}, 214 {2.0/9.0,1.0/3.0,4.0/9.0,0}}, 215 b[4] = {2.0/9.0,1.0/3.0,4.0/9.0,0}, 216 bembed[4] = {7.0/24.0,1.0/4.0,1.0/3.0,1.0/8.0}; 217 #endif 218 ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 219 } 220 { 221 const PetscReal 222 A[4][4] = {{0,0,0,0}, 223 {0.5,0,0,0}, 224 {0,0.5,0,0}, 225 {0,0,1.0,0}}, 226 #if defined(PETSC_USE_REAL___FLOAT128) 227 b[4] = {1.0q/6.0q,1.0q/3.0q,1.0q/3.0q,1.0q/6.0q}; 228 #else 229 b[4] = {1.0/6.0,1.0/3.0,1.0/3.0,1.0/6.0}; 230 #endif 231 ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 232 } 233 { 234 const PetscReal 235 #if defined(PETSC_USE_REAL___FLOAT128) 236 A[6][6] = {{0,0,0,0,0,0}, 237 {0.25q,0,0,0,0,0}, 238 {3.0q/32.0q,9.0q/32.0q,0,0,0,0}, 239 {1932.0q/2197.0q,-7200.0q/2197.0q,7296.0q/2197.0q,0,0,0}, 240 {439.0q/216.0q,-8.0q,3680.0q/513.0q,-845.0q/4104.0q,0,0}, 241 {-8.0q/27.0q,2.0q,-3544.0q/2565.0q,1859.0q/4104.0q,-11.0q/40.0q,0}}, 242 b[6] = {16.0q/135.0q,0,6656.0q/12825.0q,28561.0q/56430.0q,-9.0q/50.0q,2.0q/55.0q}, 243 bembed[6] = {25.0q/216.0q,0,1408.0q/2565.0q,2197.0q/4104.0q,-1.0q/5.0q,0}; 244 #else 245 A[6][6] = {{0,0,0,0,0,0}, 246 {0.25,0,0,0,0,0}, 247 {3.0/32.0,9.0/32.0,0,0,0,0}, 248 {1932.0/2197.0,-7200.0/2197.0,7296.0/2197.0,0,0,0}, 249 {439.0/216.0,-8.0,3680.0/513.0,-845.0/4104.0,0,0}, 250 {-8.0/27.0,2.0,-3544.0/2565.0,1859.0/4104.0,-11.0/40.0,0}}, 251 b[6] = {16.0/135.0,0,6656.0/12825.0,28561.0/56430.0,-9.0/50.0,2.0/55.0}, 252 bembed[6] = {25.0/216.0,0,1408.0/2565.0,2197.0/4104.0,-1.0/5.0,0}; 253 #endif 254 ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 255 } 256 { 257 const PetscReal 258 #if defined(PETSC_USE_REAL___FLOAT128) 259 A[7][7] = {{0,0,0,0,0,0,0}, 260 {1.0q/5.0q,0,0,0,0,0,0}, 261 {3.0q/40.0q,9.0q/40.0q,0,0,0,0,0}, 262 {44.0q/45.0q,-56.0q/15.0q,32.0q/9.0q,0,0,0,0}, 263 {19372.0q/6561.0q,-25360.0q/2187.0q,64448.0q/6561.0q,-212.0q/729.0q,0,0,0}, 264 {9017.0q/3168.0q,-355.0q/33.0q,46732.0q/5247.0q,49.0q/176.0q,-5103.0q/18656.0q,0,0}, 265 {35.0q/384.0q,0,500.0q/1113.0q,125.0q/192.0q,-2187.0q/6784.0q,11.0q/84.0q,0}}, 266 b[7] = {35.0q/384.0q,0,500.0q/1113.0q,125.0q/192.0q,-2187.0q/6784.0q,11.0q/84.0q,0}, 267 bembed[7] = {5179.0q/57600.0q,0,7571.0q/16695.0q,393.0q/640.0q,-92097.0q/339200.0q,187.0q/2100.0q,1.0q/40.0q}; 268 #else 269 A[7][7] = {{0,0,0,0,0,0,0}, 270 {1.0/5.0,0,0,0,0,0,0}, 271 {3.0/40.0,9.0/40.0,0,0,0,0,0}, 272 {44.0/45.0,-56.0/15.0,32.0/9.0,0,0,0,0}, 273 {19372.0/6561.0,-25360.0/2187.0,64448.0/6561.0,-212.0/729.0,0,0,0}, 274 {9017.0/3168.0,-355.0/33.0,46732.0/5247.0,49.0/176.0,-5103.0/18656.0,0,0}, 275 {35.0/384.0,0,500.0/1113.0,125.0/192.0,-2187.0/6784.0,11.0/84.0,0}}, 276 b[7] = {35.0/384.0,0,500.0/1113.0,125.0/192.0,-2187.0/6784.0,11.0/84.0,0}, 277 bembed[7] = {5179.0/57600.0,0,7571.0/16695.0,393.0/640.0,-92097.0/339200.0,187.0/2100.0,1.0/40.0}; 278 #endif 279 ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 280 } 281 { 282 const PetscReal 283 #if defined(PETSC_USE_REAL___FLOAT128) 284 A[8][8] = {{0,0,0,0,0,0,0,0}, 285 {1.0q/6.0q,0,0,0,0,0,0,0}, 286 {2.0q/27.0q,4.0q/27.0q,0,0,0,0,0,0}, 287 {183.0q/1372.0q,-162.0q/343.0q,1053.0q/1372.0q,0,0,0,0,0}, 288 {68.0q/297.0q,-4.0q/11.0q,42.0q/143.0q,1960.0q/3861.0q,0,0,0,0}, 289 {597.0q/22528.0q,81.0q/352.0q,63099.0q/585728.0q,58653.0q/366080.0q,4617.0q/20480.0q,0,0,0}, 290 {174197.0q/959244.0q,-30942.0q/79937.0q,8152137.0q/19744439.0q,666106.0q/1039181.0q,-29421.0q/29068.0q,482048.0q/414219.0q,0,0}, 291 {587.0q/8064.0q,0,4440339.0q/15491840.0q,24353.0q/124800.0q,387.0q/44800.0q,2152.0q/5985.0q,7267.0q/94080.0q,0}}, 292 b[8] = {587.0q/8064.0q,0,4440339.0q/15491840.0q,24353.0q/124800.0q,387.0q/44800.0q,2152.0q/5985.0q,7267.0q/94080.0q,0}, 293 bembed[8] = {2479.0q/34992.0q,0,123.0q/416.0q,612941.0q/3411720.0q,43.0q/1440.0q,2272.0q/6561.0q,79937.0q/1113912.0q,3293.0q/556956.0q}; 294 #else 295 A[8][8] = {{0,0,0,0,0,0,0,0}, 296 {1.0/6.0,0,0,0,0,0,0,0}, 297 {2.0/27.0,4.0/27.0,0,0,0,0,0,0}, 298 {183.0/1372.0,-162.0/343.0,1053.0/1372.0,0,0,0,0,0}, 299 {68.0/297.0,-4.0/11.0,42.0/143.0,1960.0/3861.0,0,0,0,0}, 300 {597.0/22528.0,81.0/352.0,63099.0/585728.0,58653.0/366080.0,4617.0/20480.0,0,0,0}, 301 {174197.0/959244.0,-30942.0/79937.0,8152137.0/19744439.0,666106.0/1039181.0,-29421.0/29068.0,482048.0/414219.0,0,0}, 302 {587.0/8064.0,0,4440339.0/15491840.0,24353.0/124800.0,387.0/44800.0,2152.0/5985.0,7267.0/94080.0,0}}, 303 b[8] = {587.0/8064.0,0,4440339.0/15491840.0,24353.0/124800.0,387.0/44800.0,2152.0/5985.0,7267.0/94080.0,0}, 304 bembed[8] = {2479.0/34992.0,0,123.0/416.0,612941.0/3411720.0,43.0/1440.0,2272.0/6561.0,79937.0/1113912.0,3293.0/556956.0}; 305 #endif 306 ierr = TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 307 } 308 PetscFunctionReturn(0); 309 } 310 311 /*@C 312 TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister(). 313 314 Not Collective 315 316 Level: advanced 317 318 .keywords: TSRK, register, destroy 319 .seealso: TSRKRegister(), TSRKRegisterAll() 320 @*/ 321 PetscErrorCode TSRKRegisterDestroy(void) 322 { 323 PetscErrorCode ierr; 324 RKTableauLink link; 325 326 PetscFunctionBegin; 327 while ((link = RKTableauList)) { 328 RKTableau t = &link->tab; 329 RKTableauList = link->next; 330 ierr = PetscFree3(t->A,t->b,t->c); CHKERRQ(ierr); 331 ierr = PetscFree (t->bembed); CHKERRQ(ierr); 332 ierr = PetscFree (t->binterp); CHKERRQ(ierr); 333 ierr = PetscFree (t->name); CHKERRQ(ierr); 334 ierr = PetscFree (link); CHKERRQ(ierr); 335 } 336 TSRKRegisterAllCalled = PETSC_FALSE; 337 PetscFunctionReturn(0); 338 } 339 340 /*@C 341 TSRKInitializePackage - This function initializes everything in the TSRK package. It is called 342 from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RK() 343 when using static libraries. 344 345 Level: developer 346 347 .keywords: TS, TSRK, initialize, package 348 .seealso: PetscInitialize() 349 @*/ 350 PetscErrorCode TSRKInitializePackage(void) 351 { 352 PetscErrorCode ierr; 353 354 PetscFunctionBegin; 355 if (TSRKPackageInitialized) PetscFunctionReturn(0); 356 TSRKPackageInitialized = PETSC_TRUE; 357 ierr = TSRKRegisterAll();CHKERRQ(ierr); 358 ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr); 359 PetscFunctionReturn(0); 360 } 361 362 /*@C 363 TSRKFinalizePackage - This function destroys everything in the TSRK package. It is 364 called from PetscFinalize(). 365 366 Level: developer 367 368 .keywords: Petsc, destroy, package 369 .seealso: PetscFinalize() 370 @*/ 371 PetscErrorCode TSRKFinalizePackage(void) 372 { 373 PetscErrorCode ierr; 374 375 PetscFunctionBegin; 376 TSRKPackageInitialized = PETSC_FALSE; 377 ierr = TSRKRegisterDestroy();CHKERRQ(ierr); 378 PetscFunctionReturn(0); 379 } 380 381 /*@C 382 TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 383 384 Not Collective, but the same schemes should be registered on all processes on which they will be used 385 386 Input Parameters: 387 + name - identifier for method 388 . order - approximation order of method 389 . s - number of stages, this is the dimension of the matrices below 390 . A - stage coefficients (dimension s*s, row-major) 391 . b - step completion table (dimension s; NULL to use last row of A) 392 . c - abscissa (dimension s; NULL to use row sums of A) 393 . bembed - completion table for embedded method (dimension s; NULL if not available) 394 . p - Order of the interpolation scheme, equal to the number of columns of binterp 395 - binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1) 396 397 Notes: 398 Several RK methods are provided, this function is only needed to create new methods. 399 400 Level: advanced 401 402 .keywords: TS, register 403 404 .seealso: TSRK 405 @*/ 406 PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s, 407 const PetscReal A[],const PetscReal b[],const PetscReal c[], 408 const PetscReal bembed[],PetscInt p,const PetscReal binterp[]) 409 { 410 PetscErrorCode ierr; 411 RKTableauLink link; 412 RKTableau t; 413 PetscInt i,j; 414 415 PetscFunctionBegin; 416 PetscValidCharPointer(name,1); 417 PetscValidRealPointer(A,4); 418 if (b) PetscValidRealPointer(b,5); 419 if (c) PetscValidRealPointer(c,6); 420 if (bembed) PetscValidRealPointer(bembed,7); 421 if (binterp || p > 1) PetscValidRealPointer(binterp,9); 422 423 ierr = PetscNew(&link);CHKERRQ(ierr); 424 t = &link->tab; 425 426 ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 427 t->order = order; 428 t->s = s; 429 ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr); 430 ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 431 if (b) { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); } 432 else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i]; 433 if (c) { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); } 434 else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j]; 435 t->FSAL = PETSC_TRUE; 436 for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE; 437 438 if (bembed) { 439 ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr); 440 ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr); 441 } 442 443 if (!binterp) { p = 1; binterp = t->b; } 444 t->p = p; 445 ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr); 446 ierr = PetscMemcpy(t->binterp,binterp,s*p*sizeof(binterp[0]));CHKERRQ(ierr); 447 448 link->next = RKTableauList; 449 RKTableauList = link; 450 PetscFunctionReturn(0); 451 } 452 453 /* 454 The step completion formula is 455 456 x1 = x0 + h b^T YdotRHS 457 458 This function can be called before or after ts->vec_sol has been updated. 459 Suppose we have a completion formula (b) and an embedded formula (be) of different order. 460 We can write 461 462 x1e = x0 + h be^T YdotRHS 463 = x1 - h b^T YdotRHS + h be^T YdotRHS 464 = x1 + h (be - b)^T YdotRHS 465 466 so we can evaluate the method with different order even after the step has been optimistically completed. 467 */ 468 static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done) 469 { 470 TS_RK *rk = (TS_RK*)ts->data; 471 RKTableau tab = rk->tableau; 472 PetscScalar *w = rk->work; 473 PetscReal h; 474 PetscInt s = tab->s,j; 475 PetscErrorCode ierr; 476 477 PetscFunctionBegin; 478 switch (rk->status) { 479 case TS_STEP_INCOMPLETE: 480 case TS_STEP_PENDING: 481 h = ts->time_step; break; 482 case TS_STEP_COMPLETE: 483 h = ts->ptime - ts->ptime_prev; break; 484 default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 485 } 486 if (order == tab->order) { 487 if (rk->status == TS_STEP_INCOMPLETE) { 488 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 489 for (j=0; j<s; j++) w[j] = h*tab->b[j]; 490 ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 491 } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 492 PetscFunctionReturn(0); 493 } else if (order == tab->order-1) { 494 if (!tab->bembed) goto unavailable; 495 if (rk->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (be) */ 496 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 497 for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 498 ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 499 } else { /* Rollback and re-complete using (be-b) */ 500 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 501 for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 502 ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 503 if (ts->vec_costintegral && ts->costintegralfwd) { 504 ierr = VecCopy(rk->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr); 505 } 506 } 507 if (done) *done = PETSC_TRUE; 508 PetscFunctionReturn(0); 509 } 510 unavailable: 511 if (done) *done = PETSC_FALSE; 512 else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"RK '%s' of order %D cannot evaluate step at order %D. Consider using -ts_adapt_type none or a different method that has an embedded estimate.",tab->name,tab->order,order); 513 PetscFunctionReturn(0); 514 } 515 516 static PetscErrorCode TSForwardCostIntegral_RK(TS ts) 517 { 518 TS_RK *rk = (TS_RK*)ts->data; 519 RKTableau tab = rk->tableau; 520 const PetscInt s = tab->s; 521 const PetscReal *b = tab->b,*c = tab->c; 522 Vec *Y = rk->Y; 523 PetscInt i; 524 PetscErrorCode ierr; 525 526 PetscFunctionBegin; 527 /* backup cost integral */ 528 ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr); 529 for (i=s-1; i>=0; i--) { 530 /* Evolve ts->vec_costintegral to compute integrals */ 531 ierr = TSAdjointComputeCostIntegrand(ts,rk->ptime+rk->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 532 ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 533 } 534 PetscFunctionReturn(0); 535 } 536 537 static PetscErrorCode TSAdjointCostIntegral_RK(TS ts) 538 { 539 TS_RK *rk = (TS_RK*)ts->data; 540 RKTableau tab = rk->tableau; 541 const PetscInt s = tab->s; 542 const PetscReal *b = tab->b,*c = tab->c; 543 Vec *Y = rk->Y; 544 PetscInt i; 545 PetscErrorCode ierr; 546 547 PetscFunctionBegin; 548 for (i=s-1; i>=0; i--) { 549 /* Evolve ts->vec_costintegral to compute integrals */ 550 ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime-ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 551 ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 552 } 553 PetscFunctionReturn(0); 554 } 555 556 static PetscErrorCode TSRollBack_RK(TS ts) 557 { 558 TS_RK *rk = (TS_RK*)ts->data; 559 RKTableau tab = rk->tableau; 560 const PetscInt s = tab->s; 561 const PetscReal *b = tab->b; 562 PetscScalar *w = rk->work; 563 Vec *YdotRHS = rk->YdotRHS; 564 PetscInt j; 565 PetscReal h; 566 PetscErrorCode ierr; 567 568 PetscFunctionBegin; 569 switch (rk->status) { 570 case TS_STEP_INCOMPLETE: 571 case TS_STEP_PENDING: 572 h = ts->time_step; break; 573 case TS_STEP_COMPLETE: 574 h = ts->ptime - ts->ptime_prev; break; 575 default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 576 } 577 for (j=0; j<s; j++) w[j] = -h*b[j]; 578 ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr); 579 PetscFunctionReturn(0); 580 } 581 582 static PetscErrorCode TSStep_RK(TS ts) 583 { 584 TS_RK *rk = (TS_RK*)ts->data; 585 RKTableau tab = rk->tableau; 586 const PetscInt s = tab->s; 587 const PetscReal *A = tab->A,*c = tab->c; 588 PetscScalar *w = rk->work; 589 Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 590 PetscBool FSAL = tab->FSAL; 591 TSAdapt adapt; 592 PetscInt i,j; 593 PetscInt rejections = 0; 594 PetscBool stageok,accept = PETSC_TRUE; 595 PetscReal next_time_step = ts->time_step; 596 PetscErrorCode ierr; 597 598 PetscFunctionBegin; 599 if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE; 600 if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); } 601 602 rk->status = TS_STEP_INCOMPLETE; 603 while (!ts->reason && rk->status != TS_STEP_COMPLETE) { 604 PetscReal t = ts->ptime; 605 PetscReal h = ts->time_step; 606 for (i=0; i<s; i++) { 607 rk->stage_time = t + h*c[i]; 608 ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr); 609 ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 610 for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 611 ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 612 ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 613 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 614 ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr); 615 if (!stageok) goto reject_step; 616 if (FSAL && !i) continue; 617 ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 618 } 619 620 rk->status = TS_STEP_INCOMPLETE; 621 ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 622 rk->status = TS_STEP_PENDING; 623 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 624 ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 625 ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr); 626 ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 627 rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 628 if (!accept) { /* Roll back the current step */ 629 ierr = TSRollBack_RK(ts);CHKERRQ(ierr); 630 ts->time_step = next_time_step; 631 goto reject_step; 632 } 633 634 if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation*/ 635 rk->ptime = ts->ptime; 636 rk->time_step = ts->time_step; 637 } 638 639 ts->ptime += ts->time_step; 640 ts->time_step = next_time_step; 641 break; 642 643 reject_step: 644 ts->reject++; accept = PETSC_FALSE; 645 if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 646 ts->reason = TS_DIVERGED_STEP_REJECTED; 647 ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 648 } 649 } 650 PetscFunctionReturn(0); 651 } 652 653 static PetscErrorCode TSAdjointSetUp_RK(TS ts) 654 { 655 TS_RK *rk = (TS_RK*)ts->data; 656 RKTableau tab = rk->tableau; 657 PetscInt s = tab->s; 658 PetscErrorCode ierr; 659 660 PetscFunctionBegin; 661 if (ts->adjointsetupcalled++) PetscFunctionReturn(0); 662 ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr); 663 if(ts->vecs_sensip) { 664 ierr = VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr); 665 } 666 ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr); 667 PetscFunctionReturn(0); 668 } 669 670 static PetscErrorCode TSAdjointStep_RK(TS ts) 671 { 672 TS_RK *rk = (TS_RK*)ts->data; 673 RKTableau tab = rk->tableau; 674 const PetscInt s = tab->s; 675 const PetscReal *A = tab->A,*b = tab->b,*c = tab->c; 676 PetscScalar *w = rk->work; 677 Vec *Y = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu,*VecSensiTemp = rk->VecSensiTemp; 678 PetscInt i,j,nadj; 679 PetscReal t = ts->ptime; 680 PetscErrorCode ierr; 681 PetscReal h = ts->time_step; 682 Mat J,Jp; 683 684 PetscFunctionBegin; 685 rk->status = TS_STEP_INCOMPLETE; 686 for (i=s-1; i>=0; i--) { 687 rk->stage_time = t + h*(1.0-c[i]); 688 for (nadj=0; nadj<ts->numcost; nadj++) { 689 ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr); 690 ierr = VecScale(VecSensiTemp[nadj],-h*b[i]);CHKERRQ(ierr); 691 for (j=i+1; j<s; j++) { 692 ierr = VecAXPY(VecSensiTemp[nadj],-h*A[j*s+i],VecDeltaLam[nadj*s+j]);CHKERRQ(ierr); 693 } 694 } 695 /* Stage values of lambda */ 696 ierr = TSGetRHSJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 697 ierr = TSComputeRHSJacobian(ts,rk->stage_time,Y[i],J,Jp);CHKERRQ(ierr); 698 if (ts->vec_costintegral) { 699 ierr = TSAdjointComputeDRDYFunction(ts,rk->stage_time,Y[i],ts->vecs_drdy);CHKERRQ(ierr); 700 } 701 for (nadj=0; nadj<ts->numcost; nadj++) { 702 ierr = MatMultTranspose(J,VecSensiTemp[nadj],VecDeltaLam[nadj*s+i]);CHKERRQ(ierr); 703 if (ts->vec_costintegral) { 704 ierr = VecAXPY(VecDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdy[nadj]);CHKERRQ(ierr); 705 } 706 } 707 708 /* Stage values of mu */ 709 if(ts->vecs_sensip) { 710 ierr = TSAdjointComputeRHSJacobian(ts,rk->stage_time,Y[i],ts->Jacp);CHKERRQ(ierr); 711 if (ts->vec_costintegral) { 712 ierr = TSAdjointComputeDRDPFunction(ts,rk->stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr); 713 } 714 715 for (nadj=0; nadj<ts->numcost; nadj++) { 716 ierr = MatMultTranspose(ts->Jacp,VecSensiTemp[nadj],VecDeltaMu[nadj*s+i]);CHKERRQ(ierr); 717 if (ts->vec_costintegral) { 718 ierr = VecAXPY(VecDeltaMu[nadj*s+i],-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr); 719 } 720 } 721 } 722 } 723 724 for (j=0; j<s; j++) w[j] = 1.0; 725 for (nadj=0; nadj<ts->numcost; nadj++) { 726 ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr); 727 if(ts->vecs_sensip) { 728 ierr = VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);CHKERRQ(ierr); 729 } 730 } 731 rk->status = TS_STEP_COMPLETE; 732 PetscFunctionReturn(0); 733 } 734 735 static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X) 736 { 737 TS_RK *rk = (TS_RK*)ts->data; 738 PetscInt s = rk->tableau->s,p = rk->tableau->p,i,j; 739 PetscReal h; 740 PetscReal tt,t; 741 PetscScalar *b; 742 const PetscReal *B = rk->tableau->binterp; 743 PetscErrorCode ierr; 744 745 PetscFunctionBegin; 746 if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 747 748 switch (rk->status) { 749 case TS_STEP_INCOMPLETE: 750 case TS_STEP_PENDING: 751 h = ts->time_step; 752 t = (itime - ts->ptime)/h; 753 break; 754 case TS_STEP_COMPLETE: 755 h = ts->ptime - ts->ptime_prev; 756 t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 757 break; 758 default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 759 } 760 ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 761 for (i=0; i<s; i++) b[i] = 0; 762 for (j=0,tt=t; j<p; j++,tt*=t) { 763 for (i=0; i<s; i++) { 764 b[i] += h * B[i*p+j] * tt; 765 } 766 } 767 768 ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr); 769 ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr); 770 771 ierr = PetscFree(b);CHKERRQ(ierr); 772 PetscFunctionReturn(0); 773 } 774 775 /*------------------------------------------------------------*/ 776 777 static PetscErrorCode TSRKTableauReset(TS ts) 778 { 779 TS_RK *rk = (TS_RK*)ts->data; 780 RKTableau tab = rk->tableau; 781 PetscErrorCode ierr; 782 783 PetscFunctionBegin; 784 if (!tab) PetscFunctionReturn(0); 785 ierr = PetscFree(rk->work);CHKERRQ(ierr); 786 ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr); 787 ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr); 788 ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr); 789 ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr); 790 PetscFunctionReturn(0); 791 } 792 793 static PetscErrorCode TSReset_RK(TS ts) 794 { 795 TS_RK *rk = (TS_RK*)ts->data; 796 PetscErrorCode ierr; 797 798 PetscFunctionBegin; 799 ierr = TSRKTableauReset(ts);CHKERRQ(ierr); 800 ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr); 801 ierr = VecDestroyVecs(ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr); 802 PetscFunctionReturn(0); 803 } 804 805 static PetscErrorCode TSDestroy_RK(TS ts) 806 { 807 PetscErrorCode ierr; 808 809 PetscFunctionBegin; 810 ierr = TSReset_RK(ts);CHKERRQ(ierr); 811 ierr = PetscFree(ts->data);CHKERRQ(ierr); 812 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr); 813 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr); 814 PetscFunctionReturn(0); 815 } 816 817 818 static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx) 819 { 820 PetscFunctionBegin; 821 PetscFunctionReturn(0); 822 } 823 824 static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 825 { 826 PetscFunctionBegin; 827 PetscFunctionReturn(0); 828 } 829 830 831 static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx) 832 { 833 PetscFunctionBegin; 834 PetscFunctionReturn(0); 835 } 836 837 static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 838 { 839 840 PetscFunctionBegin; 841 PetscFunctionReturn(0); 842 } 843 /* 844 static PetscErrorCode RKSetAdjCoe(RKTableau tab) 845 { 846 PetscReal *A,*b; 847 PetscInt s,i,j; 848 PetscErrorCode ierr; 849 850 PetscFunctionBegin; 851 s = tab->s; 852 ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr); 853 854 for (i=0; i<s; i++) 855 for (j=0; j<s; j++) { 856 A[i*s+j] = (tab->b[s-1-i]==0) ? 0: -tab->A[s-1-i+(s-1-j)*s] * tab->b[s-1-j] / tab->b[s-1-i]; 857 ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr); 858 } 859 860 for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i]; 861 862 ierr = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 863 ierr = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr); 864 ierr = PetscFree2(A,b);CHKERRQ(ierr); 865 PetscFunctionReturn(0); 866 } 867 */ 868 869 static PetscErrorCode TSRKTableauSetUp(TS ts) 870 { 871 TS_RK *rk = (TS_RK*)ts->data; 872 RKTableau tab = rk->tableau; 873 PetscErrorCode ierr; 874 875 PetscFunctionBegin; 876 ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr); 877 ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr); 878 ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr); 879 PetscFunctionReturn(0); 880 } 881 882 883 static PetscErrorCode TSSetUp_RK(TS ts) 884 { 885 TS_RK *rk = (TS_RK*)ts->data; 886 PetscErrorCode ierr; 887 DM dm; 888 889 PetscFunctionBegin; 890 ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr); 891 ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr); 892 if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 893 ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr); 894 } 895 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 896 if (dm) { 897 ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 898 ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 899 } 900 PetscFunctionReturn(0); 901 } 902 903 904 /*------------------------------------------------------------*/ 905 906 static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts) 907 { 908 TS_RK *rk = (TS_RK*)ts->data; 909 PetscErrorCode ierr; 910 911 PetscFunctionBegin; 912 ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr); 913 { 914 RKTableauLink link; 915 PetscInt count,choice; 916 PetscBool flg; 917 const char **namelist; 918 for (link=RKTableauList,count=0; link; link=link->next,count++) ; 919 ierr = PetscMalloc1(count,&namelist);CHKERRQ(ierr); 920 for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 921 ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr); 922 if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);} 923 ierr = PetscFree(namelist);CHKERRQ(ierr); 924 } 925 ierr = PetscOptionsTail();CHKERRQ(ierr); 926 PetscFunctionReturn(0); 927 } 928 929 static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[]) 930 { 931 PetscErrorCode ierr; 932 PetscInt i; 933 size_t left,count; 934 char *p; 935 936 PetscFunctionBegin; 937 for (i=0,p=buf,left=len; i<n; i++) { 938 ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr); 939 if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer"); 940 left -= count; 941 p += count; 942 *p++ = ' '; 943 } 944 p[i ? 0 : -1] = 0; 945 PetscFunctionReturn(0); 946 } 947 948 static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer) 949 { 950 TS_RK *rk = (TS_RK*)ts->data; 951 PetscBool iascii; 952 PetscErrorCode ierr; 953 954 PetscFunctionBegin; 955 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 956 if (iascii) { 957 RKTableau tab = rk->tableau; 958 TSRKType rktype; 959 char buf[512]; 960 ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr); 961 ierr = PetscViewerASCIIPrintf(viewer," RK: %s\n",rktype);CHKERRQ(ierr); 962 ierr = PetscViewerASCIIPrintf(viewer," Order: %D\n",tab->order);CHKERRQ(ierr); 963 ierr = PetscViewerASCIIPrintf(viewer," FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr); 964 ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 965 ierr = PetscViewerASCIIPrintf(viewer," Abscissa c = %s\n",buf);CHKERRQ(ierr); 966 } 967 if (ts->adapt) {ierr = TSAdaptView(ts->adapt,viewer);CHKERRQ(ierr);} 968 PetscFunctionReturn(0); 969 } 970 971 static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer) 972 { 973 PetscErrorCode ierr; 974 TSAdapt adapt; 975 976 PetscFunctionBegin; 977 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 978 ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr); 979 PetscFunctionReturn(0); 980 } 981 982 /*@C 983 TSRKSetType - Set the type of RK scheme 984 985 Logically collective 986 987 Input Parameter: 988 + ts - timestepping context 989 - rktype - type of RK-scheme 990 991 Options Database: 992 . -ts_tk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs> 993 994 Level: intermediate 995 996 .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS 997 @*/ 998 PetscErrorCode TSRKSetType(TS ts,TSRKType rktype) 999 { 1000 PetscErrorCode ierr; 1001 1002 PetscFunctionBegin; 1003 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1004 PetscValidCharPointer(rktype,2); 1005 ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr); 1006 PetscFunctionReturn(0); 1007 } 1008 1009 /*@C 1010 TSRKGetType - Get the type of RK scheme 1011 1012 Logically collective 1013 1014 Input Parameter: 1015 . ts - timestepping context 1016 1017 Output Parameter: 1018 . rktype - type of RK-scheme 1019 1020 Level: intermediate 1021 1022 .seealso: TSRKGetType() 1023 @*/ 1024 PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype) 1025 { 1026 PetscErrorCode ierr; 1027 1028 PetscFunctionBegin; 1029 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1030 ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr); 1031 PetscFunctionReturn(0); 1032 } 1033 1034 static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype) 1035 { 1036 TS_RK *rk = (TS_RK*)ts->data; 1037 1038 PetscFunctionBegin; 1039 *rktype = rk->tableau->name; 1040 PetscFunctionReturn(0); 1041 } 1042 static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype) 1043 { 1044 TS_RK *rk = (TS_RK*)ts->data; 1045 PetscErrorCode ierr; 1046 PetscBool match; 1047 RKTableauLink link; 1048 1049 PetscFunctionBegin; 1050 if (rk->tableau) { 1051 ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr); 1052 if (match) PetscFunctionReturn(0); 1053 } 1054 for (link = RKTableauList; link; link=link->next) { 1055 ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr); 1056 if (match) { 1057 if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);} 1058 rk->tableau = &link->tab; 1059 if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);} 1060 ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 1061 PetscFunctionReturn(0); 1062 } 1063 } 1064 SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype); 1065 PetscFunctionReturn(0); 1066 } 1067 1068 static PetscErrorCode TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y) 1069 { 1070 TS_RK *rk = (TS_RK*)ts->data; 1071 1072 PetscFunctionBegin; 1073 *ns = rk->tableau->s; 1074 if(Y) *Y = rk->Y; 1075 PetscFunctionReturn(0); 1076 } 1077 1078 1079 /* ------------------------------------------------------------ */ 1080 /*MC 1081 TSRK - ODE and DAE solver using Runge-Kutta schemes 1082 1083 The user should provide the right hand side of the equation 1084 using TSSetRHSFunction(). 1085 1086 Notes: 1087 The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type 1088 1089 Level: beginner 1090 1091 .seealso: TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3, 1092 TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister() 1093 1094 M*/ 1095 PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts) 1096 { 1097 TS_RK *rk; 1098 PetscErrorCode ierr; 1099 1100 PetscFunctionBegin; 1101 ierr = TSRKInitializePackage();CHKERRQ(ierr); 1102 1103 ts->ops->reset = TSReset_RK; 1104 ts->ops->destroy = TSDestroy_RK; 1105 ts->ops->view = TSView_RK; 1106 ts->ops->load = TSLoad_RK; 1107 ts->ops->setup = TSSetUp_RK; 1108 ts->ops->adjointsetup = TSAdjointSetUp_RK; 1109 ts->ops->step = TSStep_RK; 1110 ts->ops->interpolate = TSInterpolate_RK; 1111 ts->ops->evaluatestep = TSEvaluateStep_RK; 1112 ts->ops->rollback = TSRollBack_RK; 1113 ts->ops->setfromoptions = TSSetFromOptions_RK; 1114 ts->ops->getstages = TSGetStages_RK; 1115 ts->ops->adjointstep = TSAdjointStep_RK; 1116 1117 ts->ops->adjointintegral = TSAdjointCostIntegral_RK; 1118 ts->ops->forwardintegral = TSForwardCostIntegral_RK; 1119 1120 ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr); 1121 ts->data = (void*)rk; 1122 1123 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr); 1124 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr); 1125 1126 ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr); 1127 PetscFunctionReturn(0); 1128 } 1129