1 /* 2 Code for timestepping with Rosenbrock W methods 3 4 Notes: 5 The general system is written as 6 7 F(t,U,Udot) = G(t,U) 8 9 where F represents the stiff part of the physics and G represents the non-stiff part. 10 This method is designed to be linearly implicit on F and can use an approximate and lagged Jacobian. 11 12 */ 13 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 14 #include <petscdm.h> 15 16 #include <petsc/private/kernels/blockinvert.h> 17 18 static TSRosWType TSRosWDefault = TSROSWRA34PW2; 19 static PetscBool TSRosWRegisterAllCalled; 20 static PetscBool TSRosWPackageInitialized; 21 22 typedef struct _RosWTableau *RosWTableau; 23 struct _RosWTableau { 24 char *name; 25 PetscInt order; /* Classical approximation order of the method */ 26 PetscInt s; /* Number of stages */ 27 PetscInt pinterp; /* Interpolation order */ 28 PetscReal *A; /* Propagation table, strictly lower triangular */ 29 PetscReal *Gamma; /* Stage table, lower triangular with nonzero diagonal */ 30 PetscBool *GammaZeroDiag; /* Diagonal entries that are zero in stage table Gamma, vector indicating explicit statages */ 31 PetscReal *GammaExplicitCorr; /* Coefficients for correction terms needed for explicit stages in transformed variables*/ 32 PetscReal *b; /* Step completion table */ 33 PetscReal *bembed; /* Step completion table for embedded method of order one less */ 34 PetscReal *ASum; /* Row sum of A */ 35 PetscReal *GammaSum; /* Row sum of Gamma, only needed for non-autonomous systems */ 36 PetscReal *At; /* Propagation table in transformed variables */ 37 PetscReal *bt; /* Step completion table in transformed variables */ 38 PetscReal *bembedt; /* Step completion table of order one less in transformed variables */ 39 PetscReal *GammaInv; /* Inverse of Gamma, used for transformed variables */ 40 PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */ 41 PetscReal *binterpt; /* Dense output formula */ 42 }; 43 typedef struct _RosWTableauLink *RosWTableauLink; 44 struct _RosWTableauLink { 45 struct _RosWTableau tab; 46 RosWTableauLink next; 47 }; 48 static RosWTableauLink RosWTableauList; 49 50 typedef struct { 51 RosWTableau tableau; 52 Vec *Y; /* States computed during the step, used to complete the step */ 53 Vec Ydot; /* Work vector holding Ydot during residual evaluation */ 54 Vec Ystage; /* Work vector for the state value at each stage */ 55 Vec Zdot; /* Ydot = Zdot + shift*Y */ 56 Vec Zstage; /* Y = Zstage + Y */ 57 Vec vec_sol_prev; /* Solution from the previous step (used for interpolation and rollback)*/ 58 PetscScalar *work; /* Scalar work space of length number of stages, used to prepare VecMAXPY() */ 59 PetscReal scoeff; /* shift = scoeff/dt */ 60 PetscReal stage_time; 61 PetscReal stage_explicit; /* Flag indicates that the current stage is explicit */ 62 PetscBool recompute_jacobian; /* Recompute the Jacobian at each stage, default is to freeze the Jacobian at the start of each step */ 63 TSStepStatus status; 64 } TS_RosW; 65 66 /*MC 67 TSROSWTHETA1 - One stage first order L-stable Rosenbrock-W scheme (aka theta method). 68 69 Only an approximate Jacobian is needed. 70 71 Level: intermediate 72 73 .seealso: [](ch_ts), `TSROSW` 74 M*/ 75 76 /*MC 77 TSROSWTHETA2 - One stage second order A-stable Rosenbrock-W scheme (aka theta method). 78 79 Only an approximate Jacobian is needed. 80 81 Level: intermediate 82 83 .seealso: [](ch_ts), `TSROSW` 84 M*/ 85 86 /*MC 87 TSROSW2M - Two stage second order L-stable Rosenbrock-W scheme. 88 89 Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2P. 90 91 Level: intermediate 92 93 .seealso: [](ch_ts), `TSROSW` 94 M*/ 95 96 /*MC 97 TSROSW2P - Two stage second order L-stable Rosenbrock-W scheme. 98 99 Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2M. 100 101 Level: intermediate 102 103 .seealso: [](ch_ts), `TSROSW` 104 M*/ 105 106 /*MC 107 TSROSWRA3PW - Three stage third order Rosenbrock-W scheme for PDAE of index 1 {cite}`rang_2005` 108 109 Only an approximate Jacobian is needed. By default, it is only recomputed once per step. 110 111 This is strongly A-stable with R(infty) = 0.73. The embedded method of order 2 is strongly A-stable with R(infty) = 0.73. 112 113 Level: intermediate 114 115 .seealso: [](ch_ts), `TSROSW` 116 M*/ 117 118 /*MC 119 TSROSWRA34PW2 - Four stage third order L-stable Rosenbrock-W scheme for PDAE of index 1 {cite}`rang_2005`. 120 121 Only an approximate Jacobian is needed. By default, it is only recomputed once per step. 122 123 This is strongly A-stable with R(infty) = 0. The embedded method of order 2 is strongly A-stable with R(infty) = 0.48. 124 125 Level: intermediate 126 127 .seealso: [](ch_ts), `TSROSW` 128 M*/ 129 130 /*MC 131 TSROSWRODAS3 - Four stage third order L-stable Rosenbrock scheme {cite}`sandu_1997` 132 133 By default, the Jacobian is only recomputed once per step. 134 135 Both the third order and embedded second order methods are stiffly accurate and L-stable. 136 137 Level: intermediate 138 139 .seealso: [](ch_ts), `TSROSW`, `TSROSWSANDU3` 140 M*/ 141 142 /*MC 143 TSROSWRODASPR - Six stage fourth order L-stable Rosenbrock scheme {cite}`rang2015improved` 144 145 By default, the Jacobian is only recomputed once per step. 146 147 Both the fourth order and embedded third order methods are stiffly accurate and L-stable. 148 The method is B_{PR} consistent of order 3, which ensures convergence order for non-stiff, medium stiff, and stiff problems. 149 150 Level: intermediate 151 152 .seealso: [](ch_ts), `TSROSW`, `TSROSWSANDU3` 153 M*/ 154 155 /*MC 156 TSROSWSANDU3 - Three stage third order L-stable Rosenbrock scheme {cite}`sandu_1997` 157 158 By default, the Jacobian is only recomputed once per step. 159 160 The third order method is L-stable, but not stiffly accurate. 161 The second order embedded method is strongly A-stable with R(infty) = 0.5. 162 The internal stages are L-stable. 163 This method is called ROS3 in {cite}`sandu_1997`. 164 165 Level: intermediate 166 167 .seealso: [](ch_ts), `TSROSW`, `TSROSWRODAS3` 168 M*/ 169 170 /*MC 171 TSROSWASSP3P3S1C - A-stable Rosenbrock-W method with SSP explicit part, third order, three stages 172 173 By default, the Jacobian is only recomputed once per step. 174 175 A-stable SPP explicit order 3, 3 stages, CFL 1 (eff = 1/3) 176 177 Level: intermediate 178 179 .seealso: [](ch_ts), `TSROSW`, `TSROSWLASSP3P4S2C`, `TSROSWLLSSP3P4S2C`, `SSP` 180 M*/ 181 182 /*MC 183 TSROSWLASSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, four stages 184 185 By default, the Jacobian is only recomputed once per step. 186 187 L-stable (A-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2) 188 189 Level: intermediate 190 191 .seealso: [](ch_ts), `TSROSW`, `TSROSWASSP3P3S1C`, `TSROSWLLSSP3P4S2C`, `TSSSP` 192 M*/ 193 194 /*MC 195 TSROSWLLSSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, four stages 196 197 By default, the Jacobian is only recomputed once per step. 198 199 L-stable (L-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2) 200 201 Level: intermediate 202 203 .seealso: [](ch_ts), `TSROSW`, `TSROSWASSP3P3S1C`, `TSROSWLASSP3P4S2C`, `TSSSP` 204 M*/ 205 206 /*MC 207 TSROSWGRK4T - four stage, fourth order Rosenbrock (not W) method from Kaps and Rentrop {cite}`kaps1979generalized` 208 209 By default, the Jacobian is only recomputed once per step. 210 211 A(89.3 degrees)-stable, |R(infty)| = 0.454. 212 213 This method does not provide a dense output formula. 214 215 Level: intermediate 216 217 Note: 218 See Section 4 Table 7.2 in {cite}`wanner1996solving` 219 220 .seealso: [](ch_ts), `TSROSW`, `TSROSWSHAMP4`, `TSROSWVELDD4`, `TSROSW4L` 221 M*/ 222 223 /*MC 224 TSROSWSHAMP4 - four stage, fourth order Rosenbrock (not W) method from Shampine {cite}`shampine1982implementation` 225 226 By default, the Jacobian is only recomputed once per step. 227 228 A-stable, |R(infty)| = 1/3. 229 230 This method does not provide a dense output formula. 231 232 Level: intermediate 233 234 Note: 235 See Section 4 Table 7.2 in {cite}`wanner1996solving` 236 237 .seealso: [](ch_ts), `TSROSW`, `TSROSWGRK4T`, `TSROSWVELDD4`, `TSROSW4L` 238 M*/ 239 240 /*MC 241 TSROSWVELDD4 - four stage, fourth order Rosenbrock (not W) method from van Veldhuizen {cite}`veldhuizen1984d` 242 243 By default, the Jacobian is only recomputed once per step. 244 245 A(89.5 degrees)-stable, |R(infty)| = 0.24. 246 247 This method does not provide a dense output formula. 248 249 Level: intermediate 250 251 Note: 252 See Section 4 Table 7.2 in {cite}`wanner1996solving` 253 254 .seealso: [](ch_ts), `TSROSW`, `TSROSWGRK4T`, `TSROSWSHAMP4`, `TSROSW4L` 255 M*/ 256 257 /*MC 258 TSROSW4L - four stage, fourth order Rosenbrock (not W) method 259 260 By default, the Jacobian is only recomputed once per step. 261 262 A-stable and L-stable 263 264 This method does not provide a dense output formula. 265 266 Level: intermediate 267 268 Note: 269 See Section 4 Table 7.2 in {cite}`wanner1996solving` 270 271 .seealso: [](ch_ts), `TSROSW`, `TSROSWGRK4T`, `TSROSWSHAMP4`, `TSROSW4L` 272 M*/ 273 274 /*@C 275 TSRosWRegisterAll - Registers all of the Rosenbrock-W methods in `TSROSW` 276 277 Not Collective, but should be called by all MPI processes which will need the schemes to be registered 278 279 Level: advanced 280 281 .seealso: [](ch_ts), `TSROSW`, `TSRosWRegisterDestroy()` 282 @*/ 283 PetscErrorCode TSRosWRegisterAll(void) 284 { 285 PetscFunctionBegin; 286 if (TSRosWRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS); 287 TSRosWRegisterAllCalled = PETSC_TRUE; 288 289 { 290 const PetscReal A = 0; 291 const PetscReal Gamma = 1; 292 const PetscReal b = 1; 293 const PetscReal binterpt = 1; 294 295 PetscCall(TSRosWRegister(TSROSWTHETA1, 1, 1, &A, &Gamma, &b, NULL, 1, &binterpt)); 296 } 297 298 { 299 const PetscReal A = 0; 300 const PetscReal Gamma = 0.5; 301 const PetscReal b = 1; 302 const PetscReal binterpt = 1; 303 304 PetscCall(TSRosWRegister(TSROSWTHETA2, 2, 1, &A, &Gamma, &b, NULL, 1, &binterpt)); 305 } 306 307 { 308 /*const PetscReal g = 1. + 1./PetscSqrtReal(2.0); Direct evaluation: 1.707106781186547524401. Used for setting up arrays of values known at compile time below. */ 309 const PetscReal A[2][2] = { 310 {0, 0}, 311 {1., 0} 312 }; 313 const PetscReal Gamma[2][2] = { 314 {1.707106781186547524401, 0 }, 315 {-2. * 1.707106781186547524401, 1.707106781186547524401} 316 }; 317 const PetscReal b[2] = {0.5, 0.5}; 318 const PetscReal b1[2] = {1.0, 0.0}; 319 PetscReal binterpt[2][2]; 320 binterpt[0][0] = 1.707106781186547524401 - 1.0; 321 binterpt[1][0] = 2.0 - 1.707106781186547524401; 322 binterpt[0][1] = 1.707106781186547524401 - 1.5; 323 binterpt[1][1] = 1.5 - 1.707106781186547524401; 324 325 PetscCall(TSRosWRegister(TSROSW2P, 2, 2, &A[0][0], &Gamma[0][0], b, b1, 2, &binterpt[0][0])); 326 } 327 { 328 /*const PetscReal g = 1. - 1./PetscSqrtReal(2.0); Direct evaluation: 0.2928932188134524755992. Used for setting up arrays of values known at compile time below. */ 329 const PetscReal A[2][2] = { 330 {0, 0}, 331 {1., 0} 332 }; 333 const PetscReal Gamma[2][2] = { 334 {0.2928932188134524755992, 0 }, 335 {-2. * 0.2928932188134524755992, 0.2928932188134524755992} 336 }; 337 const PetscReal b[2] = {0.5, 0.5}; 338 const PetscReal b1[2] = {1.0, 0.0}; 339 PetscReal binterpt[2][2]; 340 binterpt[0][0] = 0.2928932188134524755992 - 1.0; 341 binterpt[1][0] = 2.0 - 0.2928932188134524755992; 342 binterpt[0][1] = 0.2928932188134524755992 - 1.5; 343 binterpt[1][1] = 1.5 - 0.2928932188134524755992; 344 345 PetscCall(TSRosWRegister(TSROSW2M, 2, 2, &A[0][0], &Gamma[0][0], b, b1, 2, &binterpt[0][0])); 346 } 347 { 348 /*const PetscReal g = 7.8867513459481287e-01; Directly written in-place below */ 349 PetscReal binterpt[3][2]; 350 const PetscReal A[3][3] = { 351 {0, 0, 0}, 352 {1.5773502691896257e+00, 0, 0}, 353 {0.5, 0, 0} 354 }; 355 const PetscReal Gamma[3][3] = { 356 {7.8867513459481287e-01, 0, 0 }, 357 {-1.5773502691896257e+00, 7.8867513459481287e-01, 0 }, 358 {-6.7075317547305480e-01, -1.7075317547305482e-01, 7.8867513459481287e-01} 359 }; 360 const PetscReal b[3] = {1.0566243270259355e-01, 4.9038105676657971e-02, 8.4529946162074843e-01}; 361 const PetscReal b2[3] = {-1.7863279495408180e-01, 1. / 3., 8.4529946162074843e-01}; 362 363 binterpt[0][0] = -0.8094010767585034; 364 binterpt[1][0] = -0.5; 365 binterpt[2][0] = 2.3094010767585034; 366 binterpt[0][1] = 0.9641016151377548; 367 binterpt[1][1] = 0.5; 368 binterpt[2][1] = -1.4641016151377548; 369 370 PetscCall(TSRosWRegister(TSROSWRA3PW, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0])); 371 } 372 { 373 PetscReal binterpt[4][3]; 374 /*const PetscReal g = 4.3586652150845900e-01; Directly written in-place below */ 375 const PetscReal A[4][4] = { 376 {0, 0, 0, 0}, 377 {8.7173304301691801e-01, 0, 0, 0}, 378 {8.4457060015369423e-01, -1.1299064236484185e-01, 0, 0}, 379 {0, 0, 1., 0} 380 }; 381 const PetscReal Gamma[4][4] = { 382 {4.3586652150845900e-01, 0, 0, 0 }, 383 {-8.7173304301691801e-01, 4.3586652150845900e-01, 0, 0 }, 384 {-9.0338057013044082e-01, 5.4180672388095326e-02, 4.3586652150845900e-01, 0 }, 385 {2.4212380706095346e-01, -1.2232505839045147e+00, 5.4526025533510214e-01, 4.3586652150845900e-01} 386 }; 387 const PetscReal b[4] = {2.4212380706095346e-01, -1.2232505839045147e+00, 1.5452602553351020e+00, 4.3586652150845900e-01}; 388 const PetscReal b2[4] = {3.7810903145819369e-01, -9.6042292212423178e-02, 5.0000000000000000e-01, 2.1793326075422950e-01}; 389 390 binterpt[0][0] = 1.0564298455794094; 391 binterpt[1][0] = 2.296429974281067; 392 binterpt[2][0] = -1.307599564525376; 393 binterpt[3][0] = -1.045260255335102; 394 binterpt[0][1] = -1.3864882699759573; 395 binterpt[1][1] = -8.262611700275677; 396 binterpt[2][1] = 7.250979895056055; 397 binterpt[3][1] = 2.398120075195581; 398 binterpt[0][2] = 0.5721822314575016; 399 binterpt[1][2] = 4.742931142090097; 400 binterpt[2][2] = -4.398120075195578; 401 binterpt[3][2] = -0.9169932983520199; 402 403 PetscCall(TSRosWRegister(TSROSWRA34PW2, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0])); 404 } 405 { 406 /* const PetscReal g = 0.5; Directly written in-place below */ 407 const PetscReal A[4][4] = { 408 {0, 0, 0, 0}, 409 {0, 0, 0, 0}, 410 {1., 0, 0, 0}, 411 {0.75, -0.25, 0.5, 0} 412 }; 413 const PetscReal Gamma[4][4] = { 414 {0.5, 0, 0, 0 }, 415 {1., 0.5, 0, 0 }, 416 {-0.25, -0.25, 0.5, 0 }, 417 {1. / 12, 1. / 12, -2. / 3, 0.5} 418 }; 419 const PetscReal b[4] = {5. / 6, -1. / 6, -1. / 6, 0.5}; 420 const PetscReal b2[4] = {0.75, -0.25, 0.5, 0}; 421 422 PetscCall(TSRosWRegister(TSROSWRODAS3, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 0, NULL)); 423 } 424 { 425 /* const PetscReal g = 0.25; Directly written in-place below */ 426 const PetscReal A[6][6] = { 427 {0, 0, 0, 0, 0, 0}, 428 {0.75, 0, 0, 0, 0, 0}, 429 {7.5162877593868457e-02, 2.4837122406131545e-02, 0, 0, 0, 0}, 430 {1.6532708886396510e+00, 2.1545706385445562e-01, -1.3157488872766792e+00, 0, 0, 0}, 431 {1.9385003738039885e+01, 1.2007117225835324e+00, -1.9337924059522791e+01, -2.4779140110062559e-01, 0, 0}, 432 {-7.3844531665375115e+00, -3.0593419030174646e-01, 7.8622074209377981e+00, 5.7817993590145966e-01, 0.25, 0} 433 }; 434 const PetscReal Gamma[6][6] = { 435 {0.25, 0, 0, 0, 0, 0 }, 436 {-7.5000000000000000e-01, 0.25, 0, 0, 0, 0 }, 437 {-8.8644359075349941e-02, -2.8688974257983398e-02, 0.25, 0, 0, 0 }, 438 {-4.8470034585330284e+00, -3.1583244269672095e-01, 4.9536568360123221e+00, 0.25, 0, 0 }, 439 {-2.6769456904577400e+01, -1.5066459128852787e+00, 2.7200131480460591e+01, 8.2597133700208525e-01, 0.25, 0 }, 440 {6.5876206496361416e+00, 3.6807059172993878e-01, -6.7423520694658121e+00, -1.0619631475741095e-01, -3.5714285714285715e-01, 0.25} 441 }; 442 const PetscReal b[6] = {-7.9683251690137014e-01, 6.2136401428192344e-02, 1.1198553514719862e+00, 4.7198362114404874e-01, -1.0714285714285714e-01, 0.25}; 443 const PetscReal b2[6] = {-7.3844531665375115e+00, -3.0593419030174646e-01, 7.8622074209377981e+00, 5.7817993590145966e-01, 0.25, 0.0}; 444 445 PetscCall(TSRosWRegister(TSROSWRODASPR, 4, 6, &A[0][0], &Gamma[0][0], b, b2, 0, NULL)); 446 } 447 { 448 /*const PetscReal g = 0.43586652150845899941601945119356; Directly written in-place below */ 449 const PetscReal A[3][3] = { 450 {0, 0, 0}, 451 {0.43586652150845899941601945119356, 0, 0}, 452 {0.43586652150845899941601945119356, 0, 0} 453 }; 454 const PetscReal Gamma[3][3] = { 455 {0.43586652150845899941601945119356, 0, 0 }, 456 {-0.19294655696029095575009695436041, 0.43586652150845899941601945119356, 0 }, 457 {0, 1.74927148125794685173529749738960, 0.43586652150845899941601945119356} 458 }; 459 const PetscReal b[3] = {-0.75457412385404315829818998646589, 1.94100407061964420292840123379419, -0.18642994676560104463021124732829}; 460 const PetscReal b2[3] = {-1.53358745784149585370766523913002, 2.81745131148625772213931745457622, -0.28386385364476186843165221544619}; 461 462 PetscReal binterpt[3][2]; 463 binterpt[0][0] = 3.793692883777660870425141387941; 464 binterpt[1][0] = -2.918692883777660870425141387941; 465 binterpt[2][0] = 0.125; 466 binterpt[0][1] = -0.725741064379812106687651020584; 467 binterpt[1][1] = 0.559074397713145440020984353917; 468 binterpt[2][1] = 0.16666666666666666666666666666667; 469 470 PetscCall(TSRosWRegister(TSROSWSANDU3, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0])); 471 } 472 { 473 /*const PetscReal s3 = PetscSqrtReal(3.),g = (3.0+s3)/6.0; 474 * Direct evaluation: s3 = 1.732050807568877293527; 475 * g = 0.7886751345948128822546; 476 * Values are directly inserted below to ensure availability at compile time (compiler warnings otherwise...) */ 477 const PetscReal A[3][3] = { 478 {0, 0, 0}, 479 {1, 0, 0}, 480 {0.25, 0.25, 0} 481 }; 482 const PetscReal Gamma[3][3] = { 483 {0, 0, 0 }, 484 {(-3.0 - 1.732050807568877293527) / 6.0, 0.7886751345948128822546, 0 }, 485 {(-3.0 - 1.732050807568877293527) / 24.0, (-3.0 - 1.732050807568877293527) / 8.0, 0.7886751345948128822546} 486 }; 487 const PetscReal b[3] = {1. / 6., 1. / 6., 2. / 3.}; 488 const PetscReal b2[3] = {1. / 4., 1. / 4., 1. / 2.}; 489 PetscReal binterpt[3][2]; 490 491 binterpt[0][0] = 0.089316397477040902157517886164709; 492 binterpt[1][0] = -0.91068360252295909784248211383529; 493 binterpt[2][0] = 1.8213672050459181956849642276706; 494 binterpt[0][1] = 0.077350269189625764509148780501957; 495 binterpt[1][1] = 1.077350269189625764509148780502; 496 binterpt[2][1] = -1.1547005383792515290182975610039; 497 498 PetscCall(TSRosWRegister(TSROSWASSP3P3S1C, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0])); 499 } 500 501 { 502 const PetscReal A[4][4] = { 503 {0, 0, 0, 0}, 504 {1. / 2., 0, 0, 0}, 505 {1. / 2., 1. / 2., 0, 0}, 506 {1. / 6., 1. / 6., 1. / 6., 0} 507 }; 508 const PetscReal Gamma[4][4] = { 509 {1. / 2., 0, 0, 0}, 510 {0.0, 1. / 4., 0, 0}, 511 {-2., -2. / 3., 2. / 3., 0}, 512 {1. / 2., 5. / 36., -2. / 9, 0} 513 }; 514 const PetscReal b[4] = {1. / 6., 1. / 6., 1. / 6., 1. / 2.}; 515 const PetscReal b2[4] = {1. / 8., 3. / 4., 1. / 8., 0}; 516 PetscReal binterpt[4][3]; 517 518 binterpt[0][0] = 6.25; 519 binterpt[1][0] = -30.25; 520 binterpt[2][0] = 1.75; 521 binterpt[3][0] = 23.25; 522 binterpt[0][1] = -9.75; 523 binterpt[1][1] = 58.75; 524 binterpt[2][1] = -3.25; 525 binterpt[3][1] = -45.75; 526 binterpt[0][2] = 3.6666666666666666666666666666667; 527 binterpt[1][2] = -28.333333333333333333333333333333; 528 binterpt[2][2] = 1.6666666666666666666666666666667; 529 binterpt[3][2] = 23.; 530 531 PetscCall(TSRosWRegister(TSROSWLASSP3P4S2C, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0])); 532 } 533 534 { 535 const PetscReal A[4][4] = { 536 {0, 0, 0, 0}, 537 {1. / 2., 0, 0, 0}, 538 {1. / 2., 1. / 2., 0, 0}, 539 {1. / 6., 1. / 6., 1. / 6., 0} 540 }; 541 const PetscReal Gamma[4][4] = { 542 {1. / 2., 0, 0, 0}, 543 {0.0, 3. / 4., 0, 0}, 544 {-2. / 3., -23. / 9., 2. / 9., 0}, 545 {1. / 18., 65. / 108., -2. / 27, 0} 546 }; 547 const PetscReal b[4] = {1. / 6., 1. / 6., 1. / 6., 1. / 2.}; 548 const PetscReal b2[4] = {3. / 16., 10. / 16., 3. / 16., 0}; 549 PetscReal binterpt[4][3]; 550 551 binterpt[0][0] = 1.6911764705882352941176470588235; 552 binterpt[1][0] = 3.6813725490196078431372549019608; 553 binterpt[2][0] = 0.23039215686274509803921568627451; 554 binterpt[3][0] = -4.6029411764705882352941176470588; 555 binterpt[0][1] = -0.95588235294117647058823529411765; 556 binterpt[1][1] = -6.2401960784313725490196078431373; 557 binterpt[2][1] = -0.31862745098039215686274509803922; 558 binterpt[3][1] = 7.5147058823529411764705882352941; 559 binterpt[0][2] = -0.56862745098039215686274509803922; 560 binterpt[1][2] = 2.7254901960784313725490196078431; 561 binterpt[2][2] = 0.25490196078431372549019607843137; 562 binterpt[3][2] = -2.4117647058823529411764705882353; 563 564 PetscCall(TSRosWRegister(TSROSWLLSSP3P4S2C, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0])); 565 } 566 567 { 568 PetscReal A[4][4], Gamma[4][4], b[4], b2[4]; 569 PetscReal binterpt[4][3]; 570 571 Gamma[0][0] = 0.4358665215084589994160194475295062513822671686978816; 572 Gamma[0][1] = 0; 573 Gamma[0][2] = 0; 574 Gamma[0][3] = 0; 575 Gamma[1][0] = -1.997527830934941248426324674704153457289527280554476; 576 Gamma[1][1] = 0.4358665215084589994160194475295062513822671686978816; 577 Gamma[1][2] = 0; 578 Gamma[1][3] = 0; 579 Gamma[2][0] = -1.007948511795029620852002345345404191008352770119903; 580 Gamma[2][1] = -0.004648958462629345562774289390054679806993396798458131; 581 Gamma[2][2] = 0.4358665215084589994160194475295062513822671686978816; 582 Gamma[2][3] = 0; 583 Gamma[3][0] = -0.6685429734233467180451604600279552604364311322650783; 584 Gamma[3][1] = 0.6056625986449338476089525334450053439525178740492984; 585 Gamma[3][2] = -0.9717899277217721234705114616271378792182450260943198; 586 Gamma[3][3] = 0; 587 588 A[0][0] = 0; 589 A[0][1] = 0; 590 A[0][2] = 0; 591 A[0][3] = 0; 592 A[1][0] = 0.8717330430169179988320388950590125027645343373957631; 593 A[1][1] = 0; 594 A[1][2] = 0; 595 A[1][3] = 0; 596 A[2][0] = 0.5275890119763004115618079766722914408876108660811028; 597 A[2][1] = 0.07241098802369958843819203208518599088698057726988732; 598 A[2][2] = 0; 599 A[2][3] = 0; 600 A[3][0] = 0.3990960076760701320627260685975778145384666450351314; 601 A[3][1] = -0.4375576546135194437228463747348862825846903771419953; 602 A[3][2] = 1.038461646937449311660120300601880176655352737312713; 603 A[3][3] = 0; 604 605 b[0] = 0.1876410243467238251612921333138006734899663569186926; 606 b[1] = -0.5952974735769549480478230473706443582188442040780541; 607 b[2] = 0.9717899277217721234705114616271378792182450260943198; 608 b[3] = 0.4358665215084589994160194475295062513822671686978816; 609 610 b2[0] = 0.2147402862233891404862383521089097657790734483804460; 611 b2[1] = -0.4851622638849390928209050538171743017757490232519684; 612 b2[2] = 0.8687250025203875511662123688667549217531982787600080; 613 b2[3] = 0.4016969751411624011684543450940068201770721128357014; 614 615 binterpt[0][0] = 2.2565812720167954547104627844105; 616 binterpt[1][0] = 1.349166413351089573796243820819; 617 binterpt[2][0] = -2.4695174540533503758652847586647; 618 binterpt[3][0] = -0.13623023131453465264142184656474; 619 binterpt[0][1] = -3.0826699111559187902922463354557; 620 binterpt[1][1] = -2.4689115685996042534544925650515; 621 binterpt[2][1] = 5.7428279814696677152129332773553; 622 binterpt[3][1] = -0.19124650171414467146619437684812; 623 binterpt[0][2] = 1.0137296634858471607430756831148; 624 binterpt[1][2] = 0.52444768167155973161042570784064; 625 binterpt[2][2] = -2.3015205996945452158771370439586; 626 binterpt[3][2] = 0.76334325453713832352363565300308; 627 628 PetscCall(TSRosWRegister(TSROSWARK3, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0])); 629 } 630 PetscCall(TSRosWRegisterRos4(TSROSWGRK4T, 0.231, PETSC_DEFAULT, PETSC_DEFAULT, 0, -0.1282612945269037e+01)); 631 PetscCall(TSRosWRegisterRos4(TSROSWSHAMP4, 0.5, PETSC_DEFAULT, PETSC_DEFAULT, 0, 125. / 108.)); 632 PetscCall(TSRosWRegisterRos4(TSROSWVELDD4, 0.22570811482256823492, PETSC_DEFAULT, PETSC_DEFAULT, 0, -1.355958941201148)); 633 PetscCall(TSRosWRegisterRos4(TSROSW4L, 0.57282, PETSC_DEFAULT, PETSC_DEFAULT, 0, -1.093502252409163)); 634 PetscFunctionReturn(PETSC_SUCCESS); 635 } 636 637 /*@C 638 TSRosWRegisterDestroy - Frees the list of schemes that were registered by `TSRosWRegister()`. 639 640 Not Collective 641 642 Level: advanced 643 644 .seealso: [](ch_ts), `TSRosWRegister()`, `TSRosWRegisterAll()` 645 @*/ 646 PetscErrorCode TSRosWRegisterDestroy(void) 647 { 648 RosWTableauLink link; 649 650 PetscFunctionBegin; 651 while ((link = RosWTableauList)) { 652 RosWTableau t = &link->tab; 653 RosWTableauList = link->next; 654 PetscCall(PetscFree5(t->A, t->Gamma, t->b, t->ASum, t->GammaSum)); 655 PetscCall(PetscFree5(t->At, t->bt, t->GammaInv, t->GammaZeroDiag, t->GammaExplicitCorr)); 656 PetscCall(PetscFree2(t->bembed, t->bembedt)); 657 PetscCall(PetscFree(t->binterpt)); 658 PetscCall(PetscFree(t->name)); 659 PetscCall(PetscFree(link)); 660 } 661 TSRosWRegisterAllCalled = PETSC_FALSE; 662 PetscFunctionReturn(PETSC_SUCCESS); 663 } 664 665 /*@C 666 TSRosWInitializePackage - This function initializes everything in the `TSROSW` package. It is called 667 from `TSInitializePackage()`. 668 669 Level: developer 670 671 .seealso: [](ch_ts), `TSROSW`, `PetscInitialize()`, `TSRosWFinalizePackage()` 672 @*/ 673 PetscErrorCode TSRosWInitializePackage(void) 674 { 675 PetscFunctionBegin; 676 if (TSRosWPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS); 677 TSRosWPackageInitialized = PETSC_TRUE; 678 PetscCall(TSRosWRegisterAll()); 679 PetscCall(PetscRegisterFinalize(TSRosWFinalizePackage)); 680 PetscFunctionReturn(PETSC_SUCCESS); 681 } 682 683 /*@C 684 TSRosWFinalizePackage - This function destroys everything in the `TSROSW` package. It is 685 called from `PetscFinalize()`. 686 687 Level: developer 688 689 .seealso: [](ch_ts), `TSROSW`, `PetscFinalize()`, `TSRosWInitializePackage()` 690 @*/ 691 PetscErrorCode TSRosWFinalizePackage(void) 692 { 693 PetscFunctionBegin; 694 TSRosWPackageInitialized = PETSC_FALSE; 695 PetscCall(TSRosWRegisterDestroy()); 696 PetscFunctionReturn(PETSC_SUCCESS); 697 } 698 699 /*@C 700 TSRosWRegister - register a `TSROSW`, Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 701 702 Not Collective, but the same schemes should be registered on all processes on which they will be used 703 704 Input Parameters: 705 + name - identifier for method 706 . order - approximation order of method 707 . s - number of stages, this is the dimension of the matrices below 708 . A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular 709 . Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal 710 . b - Step completion table (dimension s) 711 . bembed - Step completion table for a scheme of order one less (dimension s, NULL if no embedded scheme is available) 712 . pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt 713 - binterpt - Coefficients of the interpolation formula (dimension s*pinterp) 714 715 Level: advanced 716 717 Note: 718 Several Rosenbrock W methods are provided, this function is only needed to create new methods. 719 720 .seealso: [](ch_ts), `TSROSW` 721 @*/ 722 PetscErrorCode TSRosWRegister(TSRosWType name, PetscInt order, PetscInt s, const PetscReal A[], const PetscReal Gamma[], const PetscReal b[], const PetscReal bembed[], PetscInt pinterp, const PetscReal binterpt[]) 723 { 724 RosWTableauLink link; 725 RosWTableau t; 726 PetscInt i, j, k; 727 PetscScalar *GammaInv; 728 729 PetscFunctionBegin; 730 PetscAssertPointer(name, 1); 731 PetscAssertPointer(A, 4); 732 PetscAssertPointer(Gamma, 5); 733 PetscAssertPointer(b, 6); 734 if (bembed) PetscAssertPointer(bembed, 7); 735 736 PetscCall(TSRosWInitializePackage()); 737 PetscCall(PetscNew(&link)); 738 t = &link->tab; 739 PetscCall(PetscStrallocpy(name, &t->name)); 740 t->order = order; 741 t->s = s; 742 PetscCall(PetscMalloc5(s * s, &t->A, s * s, &t->Gamma, s, &t->b, s, &t->ASum, s, &t->GammaSum)); 743 PetscCall(PetscMalloc5(s * s, &t->At, s, &t->bt, s * s, &t->GammaInv, s, &t->GammaZeroDiag, s * s, &t->GammaExplicitCorr)); 744 PetscCall(PetscArraycpy(t->A, A, s * s)); 745 PetscCall(PetscArraycpy(t->Gamma, Gamma, s * s)); 746 PetscCall(PetscArraycpy(t->GammaExplicitCorr, Gamma, s * s)); 747 PetscCall(PetscArraycpy(t->b, b, s)); 748 if (bembed) { 749 PetscCall(PetscMalloc2(s, &t->bembed, s, &t->bembedt)); 750 PetscCall(PetscArraycpy(t->bembed, bembed, s)); 751 } 752 for (i = 0; i < s; i++) { 753 t->ASum[i] = 0; 754 t->GammaSum[i] = 0; 755 for (j = 0; j < s; j++) { 756 t->ASum[i] += A[i * s + j]; 757 t->GammaSum[i] += Gamma[i * s + j]; 758 } 759 } 760 PetscCall(PetscMalloc1(s * s, &GammaInv)); /* Need to use Scalar for inverse, then convert back to Real */ 761 for (i = 0; i < s * s; i++) GammaInv[i] = Gamma[i]; 762 for (i = 0; i < s; i++) { 763 if (Gamma[i * s + i] == 0.0) { 764 GammaInv[i * s + i] = 1.0; 765 t->GammaZeroDiag[i] = PETSC_TRUE; 766 } else { 767 t->GammaZeroDiag[i] = PETSC_FALSE; 768 } 769 } 770 771 switch (s) { 772 case 1: 773 GammaInv[0] = 1. / GammaInv[0]; 774 break; 775 case 2: 776 PetscCall(PetscKernel_A_gets_inverse_A_2(GammaInv, 0, PETSC_FALSE, NULL)); 777 break; 778 case 3: 779 PetscCall(PetscKernel_A_gets_inverse_A_3(GammaInv, 0, PETSC_FALSE, NULL)); 780 break; 781 case 4: 782 PetscCall(PetscKernel_A_gets_inverse_A_4(GammaInv, 0, PETSC_FALSE, NULL)); 783 break; 784 case 5: { 785 PetscInt ipvt5[5]; 786 MatScalar work5[5 * 5]; 787 PetscCall(PetscKernel_A_gets_inverse_A_5(GammaInv, ipvt5, work5, 0, PETSC_FALSE, NULL)); 788 break; 789 } 790 case 6: 791 PetscCall(PetscKernel_A_gets_inverse_A_6(GammaInv, 0, PETSC_FALSE, NULL)); 792 break; 793 case 7: 794 PetscCall(PetscKernel_A_gets_inverse_A_7(GammaInv, 0, PETSC_FALSE, NULL)); 795 break; 796 default: 797 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not implemented for %" PetscInt_FMT " stages", s); 798 } 799 for (i = 0; i < s * s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]); 800 PetscCall(PetscFree(GammaInv)); 801 802 for (i = 0; i < s; i++) { 803 for (k = 0; k < i + 1; k++) { 804 t->GammaExplicitCorr[i * s + k] = (t->GammaExplicitCorr[i * s + k]) * (t->GammaInv[k * s + k]); 805 for (j = k + 1; j < i + 1; j++) t->GammaExplicitCorr[i * s + k] += (t->GammaExplicitCorr[i * s + j]) * (t->GammaInv[j * s + k]); 806 } 807 } 808 809 for (i = 0; i < s; i++) { 810 for (j = 0; j < s; j++) { 811 t->At[i * s + j] = 0; 812 for (k = 0; k < s; k++) t->At[i * s + j] += t->A[i * s + k] * t->GammaInv[k * s + j]; 813 } 814 t->bt[i] = 0; 815 for (j = 0; j < s; j++) t->bt[i] += t->b[j] * t->GammaInv[j * s + i]; 816 if (bembed) { 817 t->bembedt[i] = 0; 818 for (j = 0; j < s; j++) t->bembedt[i] += t->bembed[j] * t->GammaInv[j * s + i]; 819 } 820 } 821 t->ccfl = 1.0; /* Fix this */ 822 823 t->pinterp = pinterp; 824 PetscCall(PetscMalloc1(s * pinterp, &t->binterpt)); 825 PetscCall(PetscArraycpy(t->binterpt, binterpt, s * pinterp)); 826 link->next = RosWTableauList; 827 RosWTableauList = link; 828 PetscFunctionReturn(PETSC_SUCCESS); 829 } 830 831 /*@C 832 TSRosWRegisterRos4 - register a fourth order Rosenbrock scheme by providing parameter choices 833 834 Not Collective, but the same schemes should be registered on all processes on which they will be used 835 836 Input Parameters: 837 + name - identifier for method 838 . gamma - leading coefficient (diagonal entry) 839 . a2 - design parameter, see Table 7.2 of Hairer&Wanner 840 . a3 - design parameter or PETSC_DEFAULT to satisfy one of the order five conditions (Eq 7.22) 841 . b3 - design parameter, see Table 7.2 of Hairer&Wanner 842 - e4 - design parameter for embedded method, see coefficient E4 in ros4.f code from Hairer 843 844 Level: developer 845 846 Notes: 847 This routine encodes the design of fourth order Rosenbrock methods as described in Hairer and Wanner volume 2. 848 It is used here to implement several methods from the book and can be used to experiment with new methods. 849 It was written this way instead of by copying coefficients in order to provide better than double precision satisfaction of the order conditions. 850 851 .seealso: [](ch_ts), `TSRosW`, `TSRosWRegister()` 852 @*/ 853 PetscErrorCode TSRosWRegisterRos4(TSRosWType name, PetscReal gamma, PetscReal a2, PetscReal a3, PetscReal b3, PetscReal e4) 854 { 855 /* Declare numeric constants so they can be quad precision without being truncated at double */ 856 const PetscReal one = 1, two = 2, three = 3, four = 4, five = 5, six = 6, eight = 8, twelve = 12, twenty = 20, twentyfour = 24, p32 = one / six - gamma + gamma * gamma, p42 = one / eight - gamma / three, p43 = one / twelve - gamma / three, p44 = one / twentyfour - gamma / two + three / two * gamma * gamma - gamma * gamma * gamma, p56 = one / twenty - gamma / four; 857 PetscReal a4, a32, a42, a43, b1, b2, b4, beta2p, beta3p, beta4p, beta32, beta42, beta43, beta32beta2p, beta4jbetajp; 858 PetscReal A[4][4], Gamma[4][4], b[4], bm[4]; 859 PetscScalar M[3][3], rhs[3]; 860 861 PetscFunctionBegin; 862 /* Step 1: choose Gamma (input) */ 863 /* Step 2: choose a2,a3,a4; b1,b2,b3,b4 to satisfy order conditions */ 864 if (a3 == (PetscReal)PETSC_DEFAULT) a3 = (one / five - a2 / four) / (one / four - a2 / three); /* Eq 7.22 */ 865 a4 = a3; /* consequence of 7.20 */ 866 867 /* Solve order conditions 7.15a, 7.15c, 7.15e */ 868 M[0][0] = one; 869 M[0][1] = one; 870 M[0][2] = one; /* 7.15a */ 871 M[1][0] = 0.0; 872 M[1][1] = a2 * a2; 873 M[1][2] = a4 * a4; /* 7.15c */ 874 M[2][0] = 0.0; 875 M[2][1] = a2 * a2 * a2; 876 M[2][2] = a4 * a4 * a4; /* 7.15e */ 877 rhs[0] = one - b3; 878 rhs[1] = one / three - a3 * a3 * b3; 879 rhs[2] = one / four - a3 * a3 * a3 * b3; 880 PetscCall(PetscKernel_A_gets_inverse_A_3(&M[0][0], 0, PETSC_FALSE, NULL)); 881 b1 = PetscRealPart(M[0][0] * rhs[0] + M[0][1] * rhs[1] + M[0][2] * rhs[2]); 882 b2 = PetscRealPart(M[1][0] * rhs[0] + M[1][1] * rhs[1] + M[1][2] * rhs[2]); 883 b4 = PetscRealPart(M[2][0] * rhs[0] + M[2][1] * rhs[1] + M[2][2] * rhs[2]); 884 885 /* Step 3 */ 886 beta43 = (p56 - a2 * p43) / (b4 * a3 * a3 * (a3 - a2)); /* 7.21 */ 887 beta32beta2p = p44 / (b4 * beta43); /* 7.15h */ 888 beta4jbetajp = (p32 - b3 * beta32beta2p) / b4; 889 M[0][0] = b2; 890 M[0][1] = b3; 891 M[0][2] = b4; 892 M[1][0] = a4 * a4 * beta32beta2p - a3 * a3 * beta4jbetajp; 893 M[1][1] = a2 * a2 * beta4jbetajp; 894 M[1][2] = -a2 * a2 * beta32beta2p; 895 M[2][0] = b4 * beta43 * a3 * a3 - p43; 896 M[2][1] = -b4 * beta43 * a2 * a2; 897 M[2][2] = 0; 898 rhs[0] = one / two - gamma; 899 rhs[1] = 0; 900 rhs[2] = -a2 * a2 * p32; 901 PetscCall(PetscKernel_A_gets_inverse_A_3(&M[0][0], 0, PETSC_FALSE, NULL)); 902 beta2p = PetscRealPart(M[0][0] * rhs[0] + M[0][1] * rhs[1] + M[0][2] * rhs[2]); 903 beta3p = PetscRealPart(M[1][0] * rhs[0] + M[1][1] * rhs[1] + M[1][2] * rhs[2]); 904 beta4p = PetscRealPart(M[2][0] * rhs[0] + M[2][1] * rhs[1] + M[2][2] * rhs[2]); 905 906 /* Step 4: back-substitute */ 907 beta32 = beta32beta2p / beta2p; 908 beta42 = (beta4jbetajp - beta43 * beta3p) / beta2p; 909 910 /* Step 5: 7.15f and 7.20, then 7.16 */ 911 a43 = 0; 912 a32 = p42 / (b3 * a3 * beta2p + b4 * a4 * beta2p); 913 a42 = a32; 914 915 A[0][0] = 0; 916 A[0][1] = 0; 917 A[0][2] = 0; 918 A[0][3] = 0; 919 A[1][0] = a2; 920 A[1][1] = 0; 921 A[1][2] = 0; 922 A[1][3] = 0; 923 A[2][0] = a3 - a32; 924 A[2][1] = a32; 925 A[2][2] = 0; 926 A[2][3] = 0; 927 A[3][0] = a4 - a43 - a42; 928 A[3][1] = a42; 929 A[3][2] = a43; 930 A[3][3] = 0; 931 Gamma[0][0] = gamma; 932 Gamma[0][1] = 0; 933 Gamma[0][2] = 0; 934 Gamma[0][3] = 0; 935 Gamma[1][0] = beta2p - A[1][0]; 936 Gamma[1][1] = gamma; 937 Gamma[1][2] = 0; 938 Gamma[1][3] = 0; 939 Gamma[2][0] = beta3p - beta32 - A[2][0]; 940 Gamma[2][1] = beta32 - A[2][1]; 941 Gamma[2][2] = gamma; 942 Gamma[2][3] = 0; 943 Gamma[3][0] = beta4p - beta42 - beta43 - A[3][0]; 944 Gamma[3][1] = beta42 - A[3][1]; 945 Gamma[3][2] = beta43 - A[3][2]; 946 Gamma[3][3] = gamma; 947 b[0] = b1; 948 b[1] = b2; 949 b[2] = b3; 950 b[3] = b4; 951 952 /* Construct embedded formula using given e4. We are solving Equation 7.18. */ 953 bm[3] = b[3] - e4 * gamma; /* using definition of E4 */ 954 bm[2] = (p32 - beta4jbetajp * bm[3]) / (beta32 * beta2p); /* fourth row of 7.18 */ 955 bm[1] = (one / two - gamma - beta3p * bm[2] - beta4p * bm[3]) / beta2p; /* second row */ 956 bm[0] = one - bm[1] - bm[2] - bm[3]; /* first row */ 957 958 { 959 const PetscReal misfit = a2 * a2 * bm[1] + a3 * a3 * bm[2] + a4 * a4 * bm[3] - one / three; 960 PetscCheck(PetscAbs(misfit) <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_SUP, "Assumptions violated, could not construct a third order embedded method"); 961 } 962 PetscCall(TSRosWRegister(name, 4, 4, &A[0][0], &Gamma[0][0], b, bm, 0, NULL)); 963 PetscFunctionReturn(PETSC_SUCCESS); 964 } 965 966 /* 967 The step completion formula is 968 969 x1 = x0 + b^T Y 970 971 where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been 972 updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write 973 974 x1e = x0 + be^T Y 975 = x1 - b^T Y + be^T Y 976 = x1 + (be - b)^T Y 977 978 so we can evaluate the method of different order even after the step has been optimistically completed. 979 */ 980 static PetscErrorCode TSEvaluateStep_RosW(TS ts, PetscInt order, Vec U, PetscBool *done) 981 { 982 TS_RosW *ros = (TS_RosW *)ts->data; 983 RosWTableau tab = ros->tableau; 984 PetscScalar *w = ros->work; 985 PetscInt i; 986 987 PetscFunctionBegin; 988 if (order == tab->order) { 989 if (ros->status == TS_STEP_INCOMPLETE) { /* Use standard completion formula */ 990 PetscCall(VecCopy(ts->vec_sol, U)); 991 for (i = 0; i < tab->s; i++) w[i] = tab->bt[i]; 992 PetscCall(VecMAXPY(U, tab->s, w, ros->Y)); 993 } else PetscCall(VecCopy(ts->vec_sol, U)); 994 if (done) *done = PETSC_TRUE; 995 PetscFunctionReturn(PETSC_SUCCESS); 996 } else if (order == tab->order - 1) { 997 if (!tab->bembedt) goto unavailable; 998 if (ros->status == TS_STEP_INCOMPLETE) { /* Use embedded completion formula */ 999 PetscCall(VecCopy(ts->vec_sol, U)); 1000 for (i = 0; i < tab->s; i++) w[i] = tab->bembedt[i]; 1001 PetscCall(VecMAXPY(U, tab->s, w, ros->Y)); 1002 } else { /* Use rollback-and-recomplete formula (bembedt - bt) */ 1003 for (i = 0; i < tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i]; 1004 PetscCall(VecCopy(ts->vec_sol, U)); 1005 PetscCall(VecMAXPY(U, tab->s, w, ros->Y)); 1006 } 1007 if (done) *done = PETSC_TRUE; 1008 PetscFunctionReturn(PETSC_SUCCESS); 1009 } 1010 unavailable: 1011 if (done) *done = PETSC_FALSE; 1012 else 1013 SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Rosenbrock-W '%s' of order %" PetscInt_FMT " cannot evaluate step at order %" PetscInt_FMT ". Consider using -ts_adapt_type none or a different method that has an embedded estimate.", tab->name, 1014 tab->order, order); 1015 PetscFunctionReturn(PETSC_SUCCESS); 1016 } 1017 1018 static PetscErrorCode TSStep_RosW(TS ts) 1019 { 1020 TS_RosW *ros = (TS_RosW *)ts->data; 1021 RosWTableau tab = ros->tableau; 1022 const PetscInt s = tab->s; 1023 const PetscReal *At = tab->At, *Gamma = tab->Gamma, *ASum = tab->ASum, *GammaInv = tab->GammaInv; 1024 const PetscReal *GammaExplicitCorr = tab->GammaExplicitCorr; 1025 const PetscBool *GammaZeroDiag = tab->GammaZeroDiag; 1026 PetscScalar *w = ros->work; 1027 Vec *Y = ros->Y, Ydot = ros->Ydot, Zdot = ros->Zdot, Zstage = ros->Zstage; 1028 SNES snes; 1029 TSAdapt adapt; 1030 PetscInt i, j, its, lits; 1031 PetscInt rejections = 0; 1032 PetscBool stageok, accept = PETSC_TRUE; 1033 PetscReal next_time_step = ts->time_step; 1034 PetscInt lag; 1035 1036 PetscFunctionBegin; 1037 if (!ts->steprollback) PetscCall(VecCopy(ts->vec_sol, ros->vec_sol_prev)); 1038 1039 ros->status = TS_STEP_INCOMPLETE; 1040 while (!ts->reason && ros->status != TS_STEP_COMPLETE) { 1041 const PetscReal h = ts->time_step; 1042 for (i = 0; i < s; i++) { 1043 ros->stage_time = ts->ptime + h * ASum[i]; 1044 PetscCall(TSPreStage(ts, ros->stage_time)); 1045 if (GammaZeroDiag[i]) { 1046 ros->stage_explicit = PETSC_TRUE; 1047 ros->scoeff = 1.; 1048 } else { 1049 ros->stage_explicit = PETSC_FALSE; 1050 ros->scoeff = 1. / Gamma[i * s + i]; 1051 } 1052 1053 PetscCall(VecCopy(ts->vec_sol, Zstage)); 1054 for (j = 0; j < i; j++) w[j] = At[i * s + j]; 1055 PetscCall(VecMAXPY(Zstage, i, w, Y)); 1056 1057 for (j = 0; j < i; j++) w[j] = 1. / h * GammaInv[i * s + j]; 1058 PetscCall(VecZeroEntries(Zdot)); 1059 PetscCall(VecMAXPY(Zdot, i, w, Y)); 1060 1061 /* Initial guess taken from last stage */ 1062 PetscCall(VecZeroEntries(Y[i])); 1063 1064 if (!ros->stage_explicit) { 1065 PetscCall(TSGetSNES(ts, &snes)); 1066 if (!ros->recompute_jacobian && !i) { 1067 PetscCall(SNESGetLagJacobian(snes, &lag)); 1068 if (lag == 1) { /* use did not set a nontrivial lag, so lag over all stages */ 1069 PetscCall(SNESSetLagJacobian(snes, -2)); /* Recompute the Jacobian on this solve, but not again for the rest of the stages */ 1070 } 1071 } 1072 PetscCall(SNESSolve(snes, NULL, Y[i])); 1073 if (!ros->recompute_jacobian && i == s - 1 && lag == 1) { PetscCall(SNESSetLagJacobian(snes, lag)); /* Set lag back to 1 so we know user did not set it */ } 1074 PetscCall(SNESGetIterationNumber(snes, &its)); 1075 PetscCall(SNESGetLinearSolveIterations(snes, &lits)); 1076 ts->snes_its += its; 1077 ts->ksp_its += lits; 1078 } else { 1079 Mat J, Jp; 1080 PetscCall(VecZeroEntries(Ydot)); /* Evaluate Y[i]=G(t,Ydot=0,Zstage) */ 1081 PetscCall(TSComputeIFunction(ts, ros->stage_time, Zstage, Ydot, Y[i], PETSC_FALSE)); 1082 PetscCall(VecScale(Y[i], -1.0)); 1083 PetscCall(VecAXPY(Y[i], -1.0, Zdot)); /*Y[i] = F(Zstage)-Zdot[=GammaInv*Y]*/ 1084 1085 PetscCall(VecZeroEntries(Zstage)); /* Zstage = GammaExplicitCorr[i,j] * Y[j] */ 1086 for (j = 0; j < i; j++) w[j] = GammaExplicitCorr[i * s + j]; 1087 PetscCall(VecMAXPY(Zstage, i, w, Y)); 1088 1089 /* Y[i] = Y[i] + Jac*Zstage[=Jac*GammaExplicitCorr[i,j] * Y[j]] */ 1090 PetscCall(TSGetIJacobian(ts, &J, &Jp, NULL, NULL)); 1091 PetscCall(TSComputeIJacobian(ts, ros->stage_time, ts->vec_sol, Ydot, 0, J, Jp, PETSC_FALSE)); 1092 PetscCall(MatMult(J, Zstage, Zdot)); 1093 PetscCall(VecAXPY(Y[i], -1.0, Zdot)); 1094 ts->ksp_its += 1; 1095 1096 PetscCall(VecScale(Y[i], h)); 1097 } 1098 PetscCall(TSPostStage(ts, ros->stage_time, i, Y)); 1099 PetscCall(TSGetAdapt(ts, &adapt)); 1100 PetscCall(TSAdaptCheckStage(adapt, ts, ros->stage_time, Y[i], &stageok)); 1101 if (!stageok) goto reject_step; 1102 } 1103 1104 ros->status = TS_STEP_INCOMPLETE; 1105 PetscCall(TSEvaluateStep_RosW(ts, tab->order, ts->vec_sol, NULL)); 1106 ros->status = TS_STEP_PENDING; 1107 PetscCall(TSGetAdapt(ts, &adapt)); 1108 PetscCall(TSAdaptCandidatesClear(adapt)); 1109 PetscCall(TSAdaptCandidateAdd(adapt, tab->name, tab->order, 1, tab->ccfl, (PetscReal)tab->s, PETSC_TRUE)); 1110 PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 1111 ros->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 1112 if (!accept) { /* Roll back the current step */ 1113 PetscCall(VecCopy(ts->vec_sol0, ts->vec_sol)); 1114 ts->time_step = next_time_step; 1115 goto reject_step; 1116 } 1117 1118 ts->ptime += ts->time_step; 1119 ts->time_step = next_time_step; 1120 break; 1121 1122 reject_step: 1123 ts->reject++; 1124 accept = PETSC_FALSE; 1125 if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 1126 ts->reason = TS_DIVERGED_STEP_REJECTED; 1127 PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections)); 1128 } 1129 } 1130 PetscFunctionReturn(PETSC_SUCCESS); 1131 } 1132 1133 static PetscErrorCode TSInterpolate_RosW(TS ts, PetscReal itime, Vec U) 1134 { 1135 TS_RosW *ros = (TS_RosW *)ts->data; 1136 PetscInt s = ros->tableau->s, pinterp = ros->tableau->pinterp, i, j; 1137 PetscReal h; 1138 PetscReal tt, t; 1139 PetscScalar *bt; 1140 const PetscReal *Bt = ros->tableau->binterpt; 1141 const PetscReal *GammaInv = ros->tableau->GammaInv; 1142 PetscScalar *w = ros->work; 1143 Vec *Y = ros->Y; 1144 1145 PetscFunctionBegin; 1146 PetscCheck(Bt, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSRosW %s does not have an interpolation formula", ros->tableau->name); 1147 1148 switch (ros->status) { 1149 case TS_STEP_INCOMPLETE: 1150 case TS_STEP_PENDING: 1151 h = ts->time_step; 1152 t = (itime - ts->ptime) / h; 1153 break; 1154 case TS_STEP_COMPLETE: 1155 h = ts->ptime - ts->ptime_prev; 1156 t = (itime - ts->ptime) / h + 1; /* In the interval [0,1] */ 1157 break; 1158 default: 1159 SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus"); 1160 } 1161 PetscCall(PetscMalloc1(s, &bt)); 1162 for (i = 0; i < s; i++) bt[i] = 0; 1163 for (j = 0, tt = t; j < pinterp; j++, tt *= t) { 1164 for (i = 0; i < s; i++) bt[i] += Bt[i * pinterp + j] * tt; 1165 } 1166 1167 /* y(t+tt*h) = y(t) + Sum bt(tt) * GammaInv * Ydot */ 1168 /* U <- 0*/ 1169 PetscCall(VecZeroEntries(U)); 1170 /* U <- Sum bt_i * GammaInv(i,1:i) * Y(1:i) */ 1171 for (j = 0; j < s; j++) w[j] = 0; 1172 for (j = 0; j < s; j++) { 1173 for (i = j; i < s; i++) w[j] += bt[i] * GammaInv[i * s + j]; 1174 } 1175 PetscCall(VecMAXPY(U, i, w, Y)); 1176 /* U <- y(t) + U */ 1177 PetscCall(VecAXPY(U, 1, ros->vec_sol_prev)); 1178 1179 PetscCall(PetscFree(bt)); 1180 PetscFunctionReturn(PETSC_SUCCESS); 1181 } 1182 1183 /*------------------------------------------------------------*/ 1184 1185 static PetscErrorCode TSRosWTableauReset(TS ts) 1186 { 1187 TS_RosW *ros = (TS_RosW *)ts->data; 1188 RosWTableau tab = ros->tableau; 1189 1190 PetscFunctionBegin; 1191 if (!tab) PetscFunctionReturn(PETSC_SUCCESS); 1192 PetscCall(VecDestroyVecs(tab->s, &ros->Y)); 1193 PetscCall(PetscFree(ros->work)); 1194 PetscFunctionReturn(PETSC_SUCCESS); 1195 } 1196 1197 static PetscErrorCode TSReset_RosW(TS ts) 1198 { 1199 TS_RosW *ros = (TS_RosW *)ts->data; 1200 1201 PetscFunctionBegin; 1202 PetscCall(TSRosWTableauReset(ts)); 1203 PetscCall(VecDestroy(&ros->Ydot)); 1204 PetscCall(VecDestroy(&ros->Ystage)); 1205 PetscCall(VecDestroy(&ros->Zdot)); 1206 PetscCall(VecDestroy(&ros->Zstage)); 1207 PetscCall(VecDestroy(&ros->vec_sol_prev)); 1208 PetscFunctionReturn(PETSC_SUCCESS); 1209 } 1210 1211 static PetscErrorCode TSRosWGetVecs(TS ts, DM dm, Vec *Ydot, Vec *Zdot, Vec *Ystage, Vec *Zstage) 1212 { 1213 TS_RosW *rw = (TS_RosW *)ts->data; 1214 1215 PetscFunctionBegin; 1216 if (Ydot) { 1217 if (dm && dm != ts->dm) { 1218 PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Ydot", Ydot)); 1219 } else *Ydot = rw->Ydot; 1220 } 1221 if (Zdot) { 1222 if (dm && dm != ts->dm) { 1223 PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Zdot", Zdot)); 1224 } else *Zdot = rw->Zdot; 1225 } 1226 if (Ystage) { 1227 if (dm && dm != ts->dm) { 1228 PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Ystage", Ystage)); 1229 } else *Ystage = rw->Ystage; 1230 } 1231 if (Zstage) { 1232 if (dm && dm != ts->dm) { 1233 PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Zstage", Zstage)); 1234 } else *Zstage = rw->Zstage; 1235 } 1236 PetscFunctionReturn(PETSC_SUCCESS); 1237 } 1238 1239 static PetscErrorCode TSRosWRestoreVecs(TS ts, DM dm, Vec *Ydot, Vec *Zdot, Vec *Ystage, Vec *Zstage) 1240 { 1241 PetscFunctionBegin; 1242 if (Ydot) { 1243 if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Ydot", Ydot)); 1244 } 1245 if (Zdot) { 1246 if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Zdot", Zdot)); 1247 } 1248 if (Ystage) { 1249 if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Ystage", Ystage)); 1250 } 1251 if (Zstage) { 1252 if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Zstage", Zstage)); 1253 } 1254 PetscFunctionReturn(PETSC_SUCCESS); 1255 } 1256 1257 static PetscErrorCode DMCoarsenHook_TSRosW(DM fine, DM coarse, void *ctx) 1258 { 1259 PetscFunctionBegin; 1260 PetscFunctionReturn(PETSC_SUCCESS); 1261 } 1262 1263 static PetscErrorCode DMRestrictHook_TSRosW(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) 1264 { 1265 TS ts = (TS)ctx; 1266 Vec Ydot, Zdot, Ystage, Zstage; 1267 Vec Ydotc, Zdotc, Ystagec, Zstagec; 1268 1269 PetscFunctionBegin; 1270 PetscCall(TSRosWGetVecs(ts, fine, &Ydot, &Ystage, &Zdot, &Zstage)); 1271 PetscCall(TSRosWGetVecs(ts, coarse, &Ydotc, &Ystagec, &Zdotc, &Zstagec)); 1272 PetscCall(MatRestrict(restrct, Ydot, Ydotc)); 1273 PetscCall(VecPointwiseMult(Ydotc, rscale, Ydotc)); 1274 PetscCall(MatRestrict(restrct, Ystage, Ystagec)); 1275 PetscCall(VecPointwiseMult(Ystagec, rscale, Ystagec)); 1276 PetscCall(MatRestrict(restrct, Zdot, Zdotc)); 1277 PetscCall(VecPointwiseMult(Zdotc, rscale, Zdotc)); 1278 PetscCall(MatRestrict(restrct, Zstage, Zstagec)); 1279 PetscCall(VecPointwiseMult(Zstagec, rscale, Zstagec)); 1280 PetscCall(TSRosWRestoreVecs(ts, fine, &Ydot, &Ystage, &Zdot, &Zstage)); 1281 PetscCall(TSRosWRestoreVecs(ts, coarse, &Ydotc, &Ystagec, &Zdotc, &Zstagec)); 1282 PetscFunctionReturn(PETSC_SUCCESS); 1283 } 1284 1285 static PetscErrorCode DMSubDomainHook_TSRosW(DM fine, DM coarse, void *ctx) 1286 { 1287 PetscFunctionBegin; 1288 PetscFunctionReturn(PETSC_SUCCESS); 1289 } 1290 1291 static PetscErrorCode DMSubDomainRestrictHook_TSRosW(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) 1292 { 1293 TS ts = (TS)ctx; 1294 Vec Ydot, Zdot, Ystage, Zstage; 1295 Vec Ydots, Zdots, Ystages, Zstages; 1296 1297 PetscFunctionBegin; 1298 PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Ystage, &Zdot, &Zstage)); 1299 PetscCall(TSRosWGetVecs(ts, subdm, &Ydots, &Ystages, &Zdots, &Zstages)); 1300 1301 PetscCall(VecScatterBegin(gscat, Ydot, Ydots, INSERT_VALUES, SCATTER_FORWARD)); 1302 PetscCall(VecScatterEnd(gscat, Ydot, Ydots, INSERT_VALUES, SCATTER_FORWARD)); 1303 1304 PetscCall(VecScatterBegin(gscat, Ystage, Ystages, INSERT_VALUES, SCATTER_FORWARD)); 1305 PetscCall(VecScatterEnd(gscat, Ystage, Ystages, INSERT_VALUES, SCATTER_FORWARD)); 1306 1307 PetscCall(VecScatterBegin(gscat, Zdot, Zdots, INSERT_VALUES, SCATTER_FORWARD)); 1308 PetscCall(VecScatterEnd(gscat, Zdot, Zdots, INSERT_VALUES, SCATTER_FORWARD)); 1309 1310 PetscCall(VecScatterBegin(gscat, Zstage, Zstages, INSERT_VALUES, SCATTER_FORWARD)); 1311 PetscCall(VecScatterEnd(gscat, Zstage, Zstages, INSERT_VALUES, SCATTER_FORWARD)); 1312 1313 PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Ystage, &Zdot, &Zstage)); 1314 PetscCall(TSRosWRestoreVecs(ts, subdm, &Ydots, &Ystages, &Zdots, &Zstages)); 1315 PetscFunctionReturn(PETSC_SUCCESS); 1316 } 1317 1318 static PetscErrorCode SNESTSFormFunction_RosW(SNES snes, Vec U, Vec F, TS ts) 1319 { 1320 TS_RosW *ros = (TS_RosW *)ts->data; 1321 Vec Ydot, Zdot, Ystage, Zstage; 1322 PetscReal shift = ros->scoeff / ts->time_step; 1323 DM dm, dmsave; 1324 1325 PetscFunctionBegin; 1326 PetscCall(SNESGetDM(snes, &dm)); 1327 PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage)); 1328 PetscCall(VecWAXPY(Ydot, shift, U, Zdot)); /* Ydot = shift*U + Zdot */ 1329 PetscCall(VecWAXPY(Ystage, 1.0, U, Zstage)); /* Ystage = U + Zstage */ 1330 dmsave = ts->dm; 1331 ts->dm = dm; 1332 PetscCall(TSComputeIFunction(ts, ros->stage_time, Ystage, Ydot, F, PETSC_FALSE)); 1333 ts->dm = dmsave; 1334 PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage)); 1335 PetscFunctionReturn(PETSC_SUCCESS); 1336 } 1337 1338 static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes, Vec U, Mat A, Mat B, TS ts) 1339 { 1340 TS_RosW *ros = (TS_RosW *)ts->data; 1341 Vec Ydot, Zdot, Ystage, Zstage; 1342 PetscReal shift = ros->scoeff / ts->time_step; 1343 DM dm, dmsave; 1344 1345 PetscFunctionBegin; 1346 /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */ 1347 PetscCall(SNESGetDM(snes, &dm)); 1348 PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage)); 1349 dmsave = ts->dm; 1350 ts->dm = dm; 1351 PetscCall(TSComputeIJacobian(ts, ros->stage_time, Ystage, Ydot, shift, A, B, PETSC_TRUE)); 1352 ts->dm = dmsave; 1353 PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage)); 1354 PetscFunctionReturn(PETSC_SUCCESS); 1355 } 1356 1357 static PetscErrorCode TSRosWTableauSetUp(TS ts) 1358 { 1359 TS_RosW *ros = (TS_RosW *)ts->data; 1360 RosWTableau tab = ros->tableau; 1361 1362 PetscFunctionBegin; 1363 PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ros->Y)); 1364 PetscCall(PetscMalloc1(tab->s, &ros->work)); 1365 PetscFunctionReturn(PETSC_SUCCESS); 1366 } 1367 1368 static PetscErrorCode TSSetUp_RosW(TS ts) 1369 { 1370 TS_RosW *ros = (TS_RosW *)ts->data; 1371 DM dm; 1372 SNES snes; 1373 TSRHSJacobianFn *rhsjacobian; 1374 1375 PetscFunctionBegin; 1376 PetscCall(TSRosWTableauSetUp(ts)); 1377 PetscCall(VecDuplicate(ts->vec_sol, &ros->Ydot)); 1378 PetscCall(VecDuplicate(ts->vec_sol, &ros->Ystage)); 1379 PetscCall(VecDuplicate(ts->vec_sol, &ros->Zdot)); 1380 PetscCall(VecDuplicate(ts->vec_sol, &ros->Zstage)); 1381 PetscCall(VecDuplicate(ts->vec_sol, &ros->vec_sol_prev)); 1382 PetscCall(TSGetDM(ts, &dm)); 1383 PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSRosW, DMRestrictHook_TSRosW, ts)); 1384 PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSRosW, DMSubDomainRestrictHook_TSRosW, ts)); 1385 /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */ 1386 PetscCall(TSGetSNES(ts, &snes)); 1387 if (!((PetscObject)snes)->type_name) PetscCall(SNESSetType(snes, SNESKSPONLY)); 1388 PetscCall(DMTSGetRHSJacobian(dm, &rhsjacobian, NULL)); 1389 if (rhsjacobian == TSComputeRHSJacobianConstant) { 1390 Mat Amat, Pmat; 1391 1392 /* Set the SNES matrix to be different from the RHS matrix because there is no way to reconstruct shift*M-J */ 1393 PetscCall(SNESGetJacobian(snes, &Amat, &Pmat, NULL, NULL)); 1394 if (Amat && Amat == ts->Arhs) { 1395 if (Amat == Pmat) { 1396 PetscCall(MatDuplicate(ts->Arhs, MAT_COPY_VALUES, &Amat)); 1397 PetscCall(SNESSetJacobian(snes, Amat, Amat, NULL, NULL)); 1398 } else { 1399 PetscCall(MatDuplicate(ts->Arhs, MAT_COPY_VALUES, &Amat)); 1400 PetscCall(SNESSetJacobian(snes, Amat, NULL, NULL, NULL)); 1401 if (Pmat && Pmat == ts->Brhs) { 1402 PetscCall(MatDuplicate(ts->Brhs, MAT_COPY_VALUES, &Pmat)); 1403 PetscCall(SNESSetJacobian(snes, NULL, Pmat, NULL, NULL)); 1404 PetscCall(MatDestroy(&Pmat)); 1405 } 1406 } 1407 PetscCall(MatDestroy(&Amat)); 1408 } 1409 } 1410 PetscFunctionReturn(PETSC_SUCCESS); 1411 } 1412 /*------------------------------------------------------------*/ 1413 1414 static PetscErrorCode TSSetFromOptions_RosW(TS ts, PetscOptionItems *PetscOptionsObject) 1415 { 1416 TS_RosW *ros = (TS_RosW *)ts->data; 1417 SNES snes; 1418 1419 PetscFunctionBegin; 1420 PetscOptionsHeadBegin(PetscOptionsObject, "RosW ODE solver options"); 1421 { 1422 RosWTableauLink link; 1423 PetscInt count, choice; 1424 PetscBool flg; 1425 const char **namelist; 1426 1427 for (link = RosWTableauList, count = 0; link; link = link->next, count++); 1428 PetscCall(PetscMalloc1(count, (char ***)&namelist)); 1429 for (link = RosWTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name; 1430 PetscCall(PetscOptionsEList("-ts_rosw_type", "Family of Rosenbrock-W method", "TSRosWSetType", (const char *const *)namelist, count, ros->tableau->name, &choice, &flg)); 1431 if (flg) PetscCall(TSRosWSetType(ts, namelist[choice])); 1432 PetscCall(PetscFree(namelist)); 1433 1434 PetscCall(PetscOptionsBool("-ts_rosw_recompute_jacobian", "Recompute the Jacobian at each stage", "TSRosWSetRecomputeJacobian", ros->recompute_jacobian, &ros->recompute_jacobian, NULL)); 1435 } 1436 PetscOptionsHeadEnd(); 1437 /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */ 1438 PetscCall(TSGetSNES(ts, &snes)); 1439 if (!((PetscObject)snes)->type_name) PetscCall(SNESSetType(snes, SNESKSPONLY)); 1440 PetscFunctionReturn(PETSC_SUCCESS); 1441 } 1442 1443 static PetscErrorCode TSView_RosW(TS ts, PetscViewer viewer) 1444 { 1445 TS_RosW *ros = (TS_RosW *)ts->data; 1446 PetscBool iascii; 1447 1448 PetscFunctionBegin; 1449 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 1450 if (iascii) { 1451 RosWTableau tab = ros->tableau; 1452 TSRosWType rostype; 1453 char buf[512]; 1454 PetscInt i; 1455 PetscReal abscissa[512]; 1456 PetscCall(TSRosWGetType(ts, &rostype)); 1457 PetscCall(PetscViewerASCIIPrintf(viewer, " Rosenbrock-W %s\n", rostype)); 1458 PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, tab->ASum)); 1459 PetscCall(PetscViewerASCIIPrintf(viewer, " Abscissa of A = %s\n", buf)); 1460 for (i = 0; i < tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i]; 1461 PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, abscissa)); 1462 PetscCall(PetscViewerASCIIPrintf(viewer, " Abscissa of A+Gamma = %s\n", buf)); 1463 } 1464 PetscFunctionReturn(PETSC_SUCCESS); 1465 } 1466 1467 static PetscErrorCode TSLoad_RosW(TS ts, PetscViewer viewer) 1468 { 1469 SNES snes; 1470 TSAdapt adapt; 1471 1472 PetscFunctionBegin; 1473 PetscCall(TSGetAdapt(ts, &adapt)); 1474 PetscCall(TSAdaptLoad(adapt, viewer)); 1475 PetscCall(TSGetSNES(ts, &snes)); 1476 PetscCall(SNESLoad(snes, viewer)); 1477 /* function and Jacobian context for SNES when used with TS is always ts object */ 1478 PetscCall(SNESSetFunction(snes, NULL, NULL, ts)); 1479 PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, ts)); 1480 PetscFunctionReturn(PETSC_SUCCESS); 1481 } 1482 1483 /*@C 1484 TSRosWSetType - Set the type of Rosenbrock-W, `TSROSW`, scheme 1485 1486 Logically Collective 1487 1488 Input Parameters: 1489 + ts - timestepping context 1490 - roswtype - type of Rosenbrock-W scheme 1491 1492 Level: beginner 1493 1494 .seealso: [](ch_ts), `TSRosWGetType()`, `TSROSW`, `TSROSW2M`, `TSROSW2P`, `TSROSWRA3PW`, `TSROSWRA34PW2`, `TSROSWRODAS3`, `TSROSWSANDU3`, `TSROSWASSP3P3S1C`, `TSROSWLASSP3P4S2C`, `TSROSWLLSSP3P4S2C`, `TSROSWARK3` 1495 @*/ 1496 PetscErrorCode TSRosWSetType(TS ts, TSRosWType roswtype) 1497 { 1498 PetscFunctionBegin; 1499 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1500 PetscAssertPointer(roswtype, 2); 1501 PetscTryMethod(ts, "TSRosWSetType_C", (TS, TSRosWType), (ts, roswtype)); 1502 PetscFunctionReturn(PETSC_SUCCESS); 1503 } 1504 1505 /*@C 1506 TSRosWGetType - Get the type of Rosenbrock-W scheme 1507 1508 Logically Collective 1509 1510 Input Parameter: 1511 . ts - timestepping context 1512 1513 Output Parameter: 1514 . rostype - type of Rosenbrock-W scheme 1515 1516 Level: intermediate 1517 1518 .seealso: [](ch_ts), `TSRosWType`, `TSRosWSetType()` 1519 @*/ 1520 PetscErrorCode TSRosWGetType(TS ts, TSRosWType *rostype) 1521 { 1522 PetscFunctionBegin; 1523 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1524 PetscUseMethod(ts, "TSRosWGetType_C", (TS, TSRosWType *), (ts, rostype)); 1525 PetscFunctionReturn(PETSC_SUCCESS); 1526 } 1527 1528 /*@C 1529 TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step. 1530 1531 Logically Collective 1532 1533 Input Parameters: 1534 + ts - timestepping context 1535 - flg - `PETSC_TRUE` to recompute the Jacobian at each stage 1536 1537 Level: intermediate 1538 1539 .seealso: [](ch_ts), `TSRosWType`, `TSRosWGetType()` 1540 @*/ 1541 PetscErrorCode TSRosWSetRecomputeJacobian(TS ts, PetscBool flg) 1542 { 1543 PetscFunctionBegin; 1544 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1545 PetscTryMethod(ts, "TSRosWSetRecomputeJacobian_C", (TS, PetscBool), (ts, flg)); 1546 PetscFunctionReturn(PETSC_SUCCESS); 1547 } 1548 1549 static PetscErrorCode TSRosWGetType_RosW(TS ts, TSRosWType *rostype) 1550 { 1551 TS_RosW *ros = (TS_RosW *)ts->data; 1552 1553 PetscFunctionBegin; 1554 *rostype = ros->tableau->name; 1555 PetscFunctionReturn(PETSC_SUCCESS); 1556 } 1557 1558 static PetscErrorCode TSRosWSetType_RosW(TS ts, TSRosWType rostype) 1559 { 1560 TS_RosW *ros = (TS_RosW *)ts->data; 1561 PetscBool match; 1562 RosWTableauLink link; 1563 1564 PetscFunctionBegin; 1565 if (ros->tableau) { 1566 PetscCall(PetscStrcmp(ros->tableau->name, rostype, &match)); 1567 if (match) PetscFunctionReturn(PETSC_SUCCESS); 1568 } 1569 for (link = RosWTableauList; link; link = link->next) { 1570 PetscCall(PetscStrcmp(link->tab.name, rostype, &match)); 1571 if (match) { 1572 if (ts->setupcalled) PetscCall(TSRosWTableauReset(ts)); 1573 ros->tableau = &link->tab; 1574 if (ts->setupcalled) PetscCall(TSRosWTableauSetUp(ts)); 1575 ts->default_adapt_type = ros->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 1576 PetscFunctionReturn(PETSC_SUCCESS); 1577 } 1578 } 1579 SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", rostype); 1580 } 1581 1582 static PetscErrorCode TSRosWSetRecomputeJacobian_RosW(TS ts, PetscBool flg) 1583 { 1584 TS_RosW *ros = (TS_RosW *)ts->data; 1585 1586 PetscFunctionBegin; 1587 ros->recompute_jacobian = flg; 1588 PetscFunctionReturn(PETSC_SUCCESS); 1589 } 1590 1591 static PetscErrorCode TSDestroy_RosW(TS ts) 1592 { 1593 PetscFunctionBegin; 1594 PetscCall(TSReset_RosW(ts)); 1595 if (ts->dm) { 1596 PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSRosW, DMRestrictHook_TSRosW, ts)); 1597 PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSRosW, DMSubDomainRestrictHook_TSRosW, ts)); 1598 } 1599 PetscCall(PetscFree(ts->data)); 1600 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWGetType_C", NULL)); 1601 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetType_C", NULL)); 1602 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetRecomputeJacobian_C", NULL)); 1603 PetscFunctionReturn(PETSC_SUCCESS); 1604 } 1605 1606 /* ------------------------------------------------------------ */ 1607 /*MC 1608 TSROSW - ODE solver using Rosenbrock-W schemes 1609 1610 These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly 1611 nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part 1612 of the equation using `TSSetIFunction()` and the non-stiff part with `TSSetRHSFunction()`. 1613 1614 Level: beginner 1615 1616 Notes: 1617 This method currently only works with autonomous ODE and DAE. 1618 1619 Consider trying `TSARKIMEX` if the stiff part is strongly nonlinear. 1620 1621 Since this uses a single linear solve per time-step if you wish to lag the jacobian or preconditioner computation you must use also -snes_lag_jacobian_persists true or -snes_lag_jacobian_preconditioner true 1622 1623 Developer Notes: 1624 Rosenbrock-W methods are typically specified for autonomous ODE 1625 1626 $ udot = f(u) 1627 1628 by the stage equations 1629 1630 $ k_i = h f(u_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j 1631 1632 and step completion formula 1633 1634 $ u_1 = u_0 + sum_j b_j k_j 1635 1636 with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(u) 1637 and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner, 1638 we define new variables for the stage equations 1639 1640 $ y_i = gamma_ij k_j 1641 1642 The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define 1643 1644 $ A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-1} 1645 1646 to rewrite the method as 1647 1648 .vb 1649 [M/(h gamma_ii) - J] y_i = f(u_0 + sum_j a_ij y_j) + M sum_j (c_ij/h) y_j 1650 u_1 = u_0 + sum_j bt_j y_j 1651 .ve 1652 1653 where we have introduced the mass matrix M. Continue by defining 1654 1655 $ ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j 1656 1657 or, more compactly in tensor notation 1658 1659 $ Ydot = 1/h (Gamma^{-1} \otimes I) Y . 1660 1661 Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current 1662 stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the 1663 equation 1664 1665 $ g(u_0 + sum_j a_ij y_j + y_i, ydot_i) = 0 1666 1667 with initial guess y_i = 0. 1668 1669 .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSRosWSetType()`, `TSRosWRegister()`, `TSROSWTHETA1`, `TSROSWTHETA2`, `TSROSW2M`, `TSROSW2P`, `TSROSWRA3PW`, `TSROSWRA34PW2`, `TSROSWRODAS3`, 1670 `TSROSWSANDU3`, `TSROSWASSP3P3S1C`, `TSROSWLASSP3P4S2C`, `TSROSWLLSSP3P4S2C`, `TSROSWGRK4T`, `TSROSWSHAMP4`, `TSROSWVELDD4`, `TSROSW4L`, `TSType` 1671 M*/ 1672 PETSC_EXTERN PetscErrorCode TSCreate_RosW(TS ts) 1673 { 1674 TS_RosW *ros; 1675 1676 PetscFunctionBegin; 1677 PetscCall(TSRosWInitializePackage()); 1678 1679 ts->ops->reset = TSReset_RosW; 1680 ts->ops->destroy = TSDestroy_RosW; 1681 ts->ops->view = TSView_RosW; 1682 ts->ops->load = TSLoad_RosW; 1683 ts->ops->setup = TSSetUp_RosW; 1684 ts->ops->step = TSStep_RosW; 1685 ts->ops->interpolate = TSInterpolate_RosW; 1686 ts->ops->evaluatestep = TSEvaluateStep_RosW; 1687 ts->ops->setfromoptions = TSSetFromOptions_RosW; 1688 ts->ops->snesfunction = SNESTSFormFunction_RosW; 1689 ts->ops->snesjacobian = SNESTSFormJacobian_RosW; 1690 1691 ts->usessnes = PETSC_TRUE; 1692 1693 PetscCall(PetscNew(&ros)); 1694 ts->data = (void *)ros; 1695 1696 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWGetType_C", TSRosWGetType_RosW)); 1697 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetType_C", TSRosWSetType_RosW)); 1698 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetRecomputeJacobian_C", TSRosWSetRecomputeJacobian_RosW)); 1699 1700 PetscCall(TSRosWSetType(ts, TSRosWDefault)); 1701 PetscFunctionReturn(PETSC_SUCCESS); 1702 } 1703