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