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