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