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