1e27a552bSJed Brown /* 261692a83SJed Brown Code for timestepping with Rosenbrock W methods 3e27a552bSJed Brown 4e27a552bSJed Brown Notes: 5e27a552bSJed Brown The general system is written as 6e27a552bSJed Brown 7f9c1d6abSBarry Smith F(t,U,Udot) = G(t,U) 8e27a552bSJed Brown 9f9c1d6abSBarry Smith where F represents the stiff part of the physics and G represents the non-stiff part. 10f9c1d6abSBarry Smith This method is designed to be linearly implicit on F and can use an approximate and lagged Jacobian. 11e27a552bSJed Brown 12e27a552bSJed Brown */ 13af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 141e25c274SJed Brown #include <petscdm.h> 15e27a552bSJed Brown 16af0996ceSBarry Smith #include <petsc/private/kernels/blockinvert.h> 1761692a83SJed Brown 1819fd82e9SBarry Smith static TSRosWType TSRosWDefault = TSROSWRA34PW2; 19e27a552bSJed Brown static PetscBool TSRosWRegisterAllCalled; 20e27a552bSJed Brown static PetscBool TSRosWPackageInitialized; 21e27a552bSJed Brown 2261692a83SJed Brown typedef struct _RosWTableau *RosWTableau; 2361692a83SJed Brown struct _RosWTableau { 24e27a552bSJed Brown char *name; 25e27a552bSJed Brown PetscInt order; /* Classical approximation order of the method */ 26e27a552bSJed Brown PetscInt s; /* Number of stages */ 27f4aed992SEmil Constantinescu PetscInt pinterp; /* Interpolation order */ 2861692a83SJed Brown PetscReal *A; /* Propagation table, strictly lower triangular */ 2961692a83SJed Brown PetscReal *Gamma; /* Stage table, lower triangular with nonzero diagonal */ 30c17803e7SJed Brown PetscBool *GammaZeroDiag; /* Diagonal entries that are zero in stage table Gamma, vector indicating explicit statages */ 3143b21953SEmil Constantinescu PetscReal *GammaExplicitCorr; /* Coefficients for correction terms needed for explicit stages in transformed variables*/ 3261692a83SJed Brown PetscReal *b; /* Step completion table */ 33fe7e6d57SJed Brown PetscReal *bembed; /* Step completion table for embedded method of order one less */ 3461692a83SJed Brown PetscReal *ASum; /* Row sum of A */ 3561692a83SJed Brown PetscReal *GammaSum; /* Row sum of Gamma, only needed for non-autonomous systems */ 3661692a83SJed Brown PetscReal *At; /* Propagation table in transformed variables */ 3761692a83SJed Brown PetscReal *bt; /* Step completion table in transformed variables */ 38fe7e6d57SJed Brown PetscReal *bembedt; /* Step completion table of order one less in transformed variables */ 3961692a83SJed Brown PetscReal *GammaInv; /* Inverse of Gamma, used for transformed variables */ 408d59e960SJed Brown PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */ 413ca35412SEmil Constantinescu PetscReal *binterpt; /* Dense output formula */ 42e27a552bSJed Brown }; 4361692a83SJed Brown typedef struct _RosWTableauLink *RosWTableauLink; 4461692a83SJed Brown struct _RosWTableauLink { 4561692a83SJed Brown struct _RosWTableau tab; 4661692a83SJed Brown RosWTableauLink next; 47e27a552bSJed Brown }; 4861692a83SJed Brown static RosWTableauLink RosWTableauList; 49e27a552bSJed Brown 50e27a552bSJed Brown typedef struct { 5161692a83SJed Brown RosWTableau tableau; 5261692a83SJed Brown Vec *Y; /* States computed during the step, used to complete the step */ 53e27a552bSJed Brown Vec Ydot; /* Work vector holding Ydot during residual evaluation */ 5461692a83SJed Brown Vec Ystage; /* Work vector for the state value at each stage */ 5561692a83SJed Brown Vec Zdot; /* Ydot = Zdot + shift*Y */ 5661692a83SJed Brown Vec Zstage; /* Y = Zstage + Y */ 57be5899b3SLisandro Dalcin Vec vec_sol_prev; /* Solution from the previous step (used for interpolation and rollback)*/ 581c3436cfSJed Brown PetscScalar *work; /* Scalar work space of length number of stages, used to prepare VecMAXPY() */ 59b296d7d5SJed Brown PetscReal scoeff; /* shift = scoeff/dt */ 60e27a552bSJed Brown PetscReal stage_time; 61c17803e7SJed Brown PetscReal stage_explicit; /* Flag indicates that the current stage is explicit */ 6261692a83SJed Brown PetscBool recompute_jacobian; /* Recompute the Jacobian at each stage, default is to freeze the Jacobian at the start of each step */ 63108c343cSJed Brown TSStepStatus status; 64e27a552bSJed Brown } TS_RosW; 65e27a552bSJed Brown 66fe7e6d57SJed Brown /*MC 673606a31eSEmil Constantinescu TSROSWTHETA1 - One stage first order L-stable Rosenbrock-W scheme (aka theta method). 683606a31eSEmil Constantinescu 693606a31eSEmil Constantinescu Only an approximate Jacobian is needed. 703606a31eSEmil Constantinescu 713606a31eSEmil Constantinescu Level: intermediate 723606a31eSEmil Constantinescu 73*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW` 743606a31eSEmil Constantinescu M*/ 753606a31eSEmil Constantinescu 763606a31eSEmil Constantinescu /*MC 773606a31eSEmil Constantinescu TSROSWTHETA2 - One stage second order A-stable Rosenbrock-W scheme (aka theta method). 783606a31eSEmil Constantinescu 793606a31eSEmil Constantinescu Only an approximate Jacobian is needed. 803606a31eSEmil Constantinescu 813606a31eSEmil Constantinescu Level: intermediate 823606a31eSEmil Constantinescu 83*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW` 843606a31eSEmil Constantinescu M*/ 853606a31eSEmil Constantinescu 863606a31eSEmil Constantinescu /*MC 87fe7e6d57SJed Brown TSROSW2M - Two stage second order L-stable Rosenbrock-W scheme. 88fe7e6d57SJed Brown 89fe7e6d57SJed Brown Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2P. 90fe7e6d57SJed Brown 91fe7e6d57SJed Brown Level: intermediate 92fe7e6d57SJed Brown 93*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW` 94fe7e6d57SJed Brown M*/ 95fe7e6d57SJed Brown 96fe7e6d57SJed Brown /*MC 97fe7e6d57SJed Brown TSROSW2P - Two stage second order L-stable Rosenbrock-W scheme. 98fe7e6d57SJed Brown 99fe7e6d57SJed Brown Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2M. 100fe7e6d57SJed Brown 101fe7e6d57SJed Brown Level: intermediate 102fe7e6d57SJed Brown 103*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW` 104fe7e6d57SJed Brown M*/ 105fe7e6d57SJed Brown 106fe7e6d57SJed Brown /*MC 107fe7e6d57SJed Brown TSROSWRA3PW - Three stage third order Rosenbrock-W scheme for PDAE of index 1. 108fe7e6d57SJed Brown 109fe7e6d57SJed Brown Only an approximate Jacobian is needed. By default, it is only recomputed once per step. 110fe7e6d57SJed Brown 111fe7e6d57SJed Brown 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. 112fe7e6d57SJed Brown 113bcf0153eSBarry Smith Level: intermediate 114bcf0153eSBarry Smith 115fe7e6d57SJed Brown References: 116606c0280SSatish Balay . * - Rang and Angermann, New Rosenbrock W methods of order 3 for partial differential algebraic equations of index 1, 2005. 117fe7e6d57SJed Brown 118*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW` 119fe7e6d57SJed Brown M*/ 120fe7e6d57SJed Brown 121fe7e6d57SJed Brown /*MC 122fe7e6d57SJed Brown TSROSWRA34PW2 - Four stage third order L-stable Rosenbrock-W scheme for PDAE of index 1. 123fe7e6d57SJed Brown 124fe7e6d57SJed Brown Only an approximate Jacobian is needed. By default, it is only recomputed once per step. 125fe7e6d57SJed Brown 126fe7e6d57SJed Brown This is strongly A-stable with R(infty) = 0. The embedded method of order 2 is strongly A-stable with R(infty) = 0.48. 127fe7e6d57SJed Brown 128bcf0153eSBarry Smith Level: intermediate 129bcf0153eSBarry Smith 130fe7e6d57SJed Brown References: 131606c0280SSatish Balay . * - Rang and Angermann, New Rosenbrock W methods of order 3 for partial differential algebraic equations of index 1, 2005. 132fe7e6d57SJed Brown 133*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW` 134fe7e6d57SJed Brown M*/ 135fe7e6d57SJed Brown 136ef3c5b88SJed Brown /*MC 137ef3c5b88SJed Brown TSROSWRODAS3 - Four stage third order L-stable Rosenbrock scheme 138ef3c5b88SJed Brown 139ef3c5b88SJed Brown By default, the Jacobian is only recomputed once per step. 140ef3c5b88SJed Brown 141ef3c5b88SJed Brown Both the third order and embedded second order methods are stiffly accurate and L-stable. 142ef3c5b88SJed Brown 143bcf0153eSBarry Smith Level: intermediate 144bcf0153eSBarry Smith 145ef3c5b88SJed Brown References: 146606c0280SSatish Balay . * - Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997. 147ef3c5b88SJed Brown 148*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWSANDU3` 149ef3c5b88SJed Brown M*/ 150ef3c5b88SJed Brown 151ef3c5b88SJed Brown /*MC 152ef3c5b88SJed Brown TSROSWSANDU3 - Three stage third order L-stable Rosenbrock scheme 153ef3c5b88SJed Brown 154ef3c5b88SJed Brown By default, the Jacobian is only recomputed once per step. 155ef3c5b88SJed Brown 156ef3c5b88SJed Brown The third order method is L-stable, but not stiffly accurate. 157ef3c5b88SJed Brown The second order embedded method is strongly A-stable with R(infty) = 0.5. 158ef3c5b88SJed Brown The internal stages are L-stable. 159ef3c5b88SJed Brown This method is called ROS3 in the paper. 160ef3c5b88SJed Brown 161bcf0153eSBarry Smith Level: intermediate 162bcf0153eSBarry Smith 163ef3c5b88SJed Brown References: 164606c0280SSatish Balay . * - Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997. 165ef3c5b88SJed Brown 166*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWRODAS3` 167ef3c5b88SJed Brown M*/ 168ef3c5b88SJed Brown 169961f28d0SJed Brown /*MC 170961f28d0SJed Brown TSROSWASSP3P3S1C - A-stable Rosenbrock-W method with SSP explicit part, third order, three stages 171961f28d0SJed Brown 172961f28d0SJed Brown By default, the Jacobian is only recomputed once per step. 173961f28d0SJed Brown 174961f28d0SJed Brown A-stable SPP explicit order 3, 3 stages, CFL 1 (eff = 1/3) 175961f28d0SJed Brown 176bcf0153eSBarry Smith Level: intermediate 177bcf0153eSBarry Smith 178961f28d0SJed Brown References: 179606c0280SSatish Balay . * - Emil Constantinescu 180961f28d0SJed Brown 181*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWLASSP3P4S2C`, `TSROSWLLSSP3P4S2C`, `SSP` 182961f28d0SJed Brown M*/ 183961f28d0SJed Brown 184961f28d0SJed Brown /*MC 185998eb97aSJed Brown TSROSWLASSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, four stages 186961f28d0SJed Brown 187961f28d0SJed Brown By default, the Jacobian is only recomputed once per step. 188961f28d0SJed Brown 189961f28d0SJed Brown L-stable (A-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2) 190961f28d0SJed Brown 191bcf0153eSBarry Smith Level: intermediate 192bcf0153eSBarry Smith 193961f28d0SJed Brown References: 194606c0280SSatish Balay . * - Emil Constantinescu 195961f28d0SJed Brown 196*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWASSP3P3S1C`, `TSROSWLLSSP3P4S2C`, `TSSSP` 197961f28d0SJed Brown M*/ 198961f28d0SJed Brown 199961f28d0SJed Brown /*MC 200998eb97aSJed Brown TSROSWLLSSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, four stages 201961f28d0SJed Brown 202961f28d0SJed Brown By default, the Jacobian is only recomputed once per step. 203961f28d0SJed Brown 204961f28d0SJed Brown L-stable (L-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2) 205961f28d0SJed Brown 206bcf0153eSBarry Smith Level: intermediate 207bcf0153eSBarry Smith 208961f28d0SJed Brown References: 209606c0280SSatish Balay . * - Emil Constantinescu 210961f28d0SJed Brown 211*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWASSP3P3S1C`, `TSROSWLASSP3P4S2C`, `TSSSP` 212961f28d0SJed Brown M*/ 213961f28d0SJed Brown 21442faf41dSJed Brown /*MC 21542faf41dSJed Brown TSROSWGRK4T - four stage, fourth order Rosenbrock (not W) method from Kaps and Rentrop 21642faf41dSJed Brown 21742faf41dSJed Brown By default, the Jacobian is only recomputed once per step. 21842faf41dSJed Brown 21942faf41dSJed Brown A(89.3 degrees)-stable, |R(infty)| = 0.454. 22042faf41dSJed Brown 22142faf41dSJed Brown This method does not provide a dense output formula. 22242faf41dSJed Brown 223bcf0153eSBarry Smith Level: intermediate 224bcf0153eSBarry Smith 22542faf41dSJed Brown References: 226606c0280SSatish Balay + * - Kaps and Rentrop, Generalized Runge Kutta methods of order four with stepsize control for stiff ordinary differential equations, 1979. 227606c0280SSatish Balay - * - Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2. 22842faf41dSJed Brown 22942faf41dSJed Brown Hairer's code ros4.f 23042faf41dSJed Brown 231*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWSHAMP4`, `TSROSWVELDD4`, `TSROSW4L` 23242faf41dSJed Brown M*/ 23342faf41dSJed Brown 23442faf41dSJed Brown /*MC 23542faf41dSJed Brown TSROSWSHAMP4 - four stage, fourth order Rosenbrock (not W) method from Shampine 23642faf41dSJed Brown 23742faf41dSJed Brown By default, the Jacobian is only recomputed once per step. 23842faf41dSJed Brown 23942faf41dSJed Brown A-stable, |R(infty)| = 1/3. 24042faf41dSJed Brown 24142faf41dSJed Brown This method does not provide a dense output formula. 24242faf41dSJed Brown 243bcf0153eSBarry Smith Level: intermediate 244bcf0153eSBarry Smith 24542faf41dSJed Brown References: 246606c0280SSatish Balay + * - Shampine, Implementation of Rosenbrock methods, 1982. 247606c0280SSatish Balay - * - Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2. 24842faf41dSJed Brown 24942faf41dSJed Brown Hairer's code ros4.f 25042faf41dSJed Brown 251*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWGRK4T`, `TSROSWVELDD4`, `TSROSW4L` 25242faf41dSJed Brown M*/ 25342faf41dSJed Brown 25442faf41dSJed Brown /*MC 25542faf41dSJed Brown TSROSWVELDD4 - four stage, fourth order Rosenbrock (not W) method from van Veldhuizen 25642faf41dSJed Brown 25742faf41dSJed Brown By default, the Jacobian is only recomputed once per step. 25842faf41dSJed Brown 25942faf41dSJed Brown A(89.5 degrees)-stable, |R(infty)| = 0.24. 26042faf41dSJed Brown 26142faf41dSJed Brown This method does not provide a dense output formula. 26242faf41dSJed Brown 263bcf0153eSBarry Smith Level: intermediate 264bcf0153eSBarry Smith 26542faf41dSJed Brown References: 266606c0280SSatish Balay + * - van Veldhuizen, D stability and Kaps Rentrop methods, 1984. 267606c0280SSatish Balay - * - Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2. 26842faf41dSJed Brown 26942faf41dSJed Brown Hairer's code ros4.f 27042faf41dSJed Brown 271*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWGRK4T`, `TSROSWSHAMP4`, `TSROSW4L` 27242faf41dSJed Brown M*/ 27342faf41dSJed Brown 27442faf41dSJed Brown /*MC 27542faf41dSJed Brown TSROSW4L - four stage, fourth order Rosenbrock (not W) method 27642faf41dSJed Brown 27742faf41dSJed Brown By default, the Jacobian is only recomputed once per step. 27842faf41dSJed Brown 27942faf41dSJed Brown A-stable and L-stable 28042faf41dSJed Brown 28142faf41dSJed Brown This method does not provide a dense output formula. 28242faf41dSJed Brown 283bcf0153eSBarry Smith Level: intermediate 284bcf0153eSBarry Smith 28542faf41dSJed Brown References: 286606c0280SSatish Balay . * - Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2. 28742faf41dSJed Brown 28842faf41dSJed Brown Hairer's code ros4.f 28942faf41dSJed Brown 290*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWGRK4T`, `TSROSWSHAMP4`, `TSROSW4L` 29142faf41dSJed Brown M*/ 29242faf41dSJed Brown 293e27a552bSJed Brown /*@C 294bcf0153eSBarry Smith TSRosWRegisterAll - Registers all of the Rosenbrock-W methods in `TSROSW` 295e27a552bSJed Brown 296e27a552bSJed Brown Not Collective, but should be called by all processes which will need the schemes to be registered 297e27a552bSJed Brown 298e27a552bSJed Brown Level: advanced 299e27a552bSJed Brown 300*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSRosWRegisterDestroy()` 301e27a552bSJed Brown @*/ 302d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWRegisterAll(void) 303d71ae5a4SJacob Faibussowitsch { 304e27a552bSJed Brown PetscFunctionBegin; 3053ba16761SJacob Faibussowitsch if (TSRosWRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS); 306e27a552bSJed Brown TSRosWRegisterAllCalled = PETSC_TRUE; 307e27a552bSJed Brown 308e27a552bSJed Brown { 309bbd56ea5SKarl Rupp const PetscReal A = 0; 310bbd56ea5SKarl Rupp const PetscReal Gamma = 1; 311bbd56ea5SKarl Rupp const PetscReal b = 1; 312bbd56ea5SKarl Rupp const PetscReal binterpt = 1; 3131f80e275SEmil Constantinescu 3149566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSWTHETA1, 1, 1, &A, &Gamma, &b, NULL, 1, &binterpt)); 3153606a31eSEmil Constantinescu } 3163606a31eSEmil Constantinescu 3173606a31eSEmil Constantinescu { 318bbd56ea5SKarl Rupp const PetscReal A = 0; 319bbd56ea5SKarl Rupp const PetscReal Gamma = 0.5; 320bbd56ea5SKarl Rupp const PetscReal b = 1; 321bbd56ea5SKarl Rupp const PetscReal binterpt = 1; 322bbd56ea5SKarl Rupp 3239566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSWTHETA2, 2, 1, &A, &Gamma, &b, NULL, 1, &binterpt)); 3243606a31eSEmil Constantinescu } 3253606a31eSEmil Constantinescu 3263606a31eSEmil Constantinescu { 327da80777bSKarl Rupp /*const PetscReal g = 1. + 1./PetscSqrtReal(2.0); Direct evaluation: 1.707106781186547524401. Used for setting up arrays of values known at compile time below. */ 3289371c9d4SSatish Balay const PetscReal A[2][2] = 3299371c9d4SSatish Balay { 3309371c9d4SSatish Balay {0, 0}, 3319371c9d4SSatish Balay {1., 0} 3329371c9d4SSatish Balay }, 3339371c9d4SSatish Balay Gamma[2][2] = {{1.707106781186547524401, 0}, {-2. * 1.707106781186547524401, 1.707106781186547524401}}, b[2] = {0.5, 0.5}, b1[2] = {1.0, 0.0}; 3341f80e275SEmil Constantinescu PetscReal binterpt[2][2]; 335da80777bSKarl Rupp binterpt[0][0] = 1.707106781186547524401 - 1.0; 336da80777bSKarl Rupp binterpt[1][0] = 2.0 - 1.707106781186547524401; 337da80777bSKarl Rupp binterpt[0][1] = 1.707106781186547524401 - 1.5; 338da80777bSKarl Rupp binterpt[1][1] = 1.5 - 1.707106781186547524401; 339bbd56ea5SKarl Rupp 3409566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSW2P, 2, 2, &A[0][0], &Gamma[0][0], b, b1, 2, &binterpt[0][0])); 341e27a552bSJed Brown } 342e27a552bSJed Brown { 343da80777bSKarl Rupp /*const PetscReal g = 1. - 1./PetscSqrtReal(2.0); Direct evaluation: 0.2928932188134524755992. Used for setting up arrays of values known at compile time below. */ 3449371c9d4SSatish Balay const PetscReal A[2][2] = 3459371c9d4SSatish Balay { 3469371c9d4SSatish Balay {0, 0}, 3479371c9d4SSatish Balay {1., 0} 3489371c9d4SSatish Balay }, 3499371c9d4SSatish Balay Gamma[2][2] = {{0.2928932188134524755992, 0}, {-2. * 0.2928932188134524755992, 0.2928932188134524755992}}, b[2] = {0.5, 0.5}, b1[2] = {1.0, 0.0}; 3501f80e275SEmil Constantinescu PetscReal binterpt[2][2]; 351da80777bSKarl Rupp binterpt[0][0] = 0.2928932188134524755992 - 1.0; 352da80777bSKarl Rupp binterpt[1][0] = 2.0 - 0.2928932188134524755992; 353da80777bSKarl Rupp binterpt[0][1] = 0.2928932188134524755992 - 1.5; 354da80777bSKarl Rupp binterpt[1][1] = 1.5 - 0.2928932188134524755992; 355bbd56ea5SKarl Rupp 3569566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSW2M, 2, 2, &A[0][0], &Gamma[0][0], b, b1, 2, &binterpt[0][0])); 357fe7e6d57SJed Brown } 358fe7e6d57SJed Brown { 359da80777bSKarl Rupp /*const PetscReal g = 7.8867513459481287e-01; Directly written in-place below */ 3601f80e275SEmil Constantinescu PetscReal binterpt[3][2]; 3619371c9d4SSatish Balay const PetscReal A[3][3] = 3629371c9d4SSatish Balay { 3639371c9d4SSatish Balay {0, 0, 0}, 364fe7e6d57SJed Brown {1.5773502691896257e+00, 0, 0}, 3659371c9d4SSatish Balay {0.5, 0, 0} 3669371c9d4SSatish Balay }, 3679371c9d4SSatish Balay Gamma[3][3] = {{7.8867513459481287e-01, 0, 0}, {-1.5773502691896257e+00, 7.8867513459481287e-01, 0}, {-6.7075317547305480e-01, -1.7075317547305482e-01, 7.8867513459481287e-01}}, b[3] = {1.0566243270259355e-01, 4.9038105676657971e-02, 8.4529946162074843e-01}, b2[3] = {-1.7863279495408180e-01, 1. / 3., 8.4529946162074843e-01}; 3681f80e275SEmil Constantinescu 3691f80e275SEmil Constantinescu binterpt[0][0] = -0.8094010767585034; 3701f80e275SEmil Constantinescu binterpt[1][0] = -0.5; 3711f80e275SEmil Constantinescu binterpt[2][0] = 2.3094010767585034; 3721f80e275SEmil Constantinescu binterpt[0][1] = 0.9641016151377548; 3731f80e275SEmil Constantinescu binterpt[1][1] = 0.5; 3741f80e275SEmil Constantinescu binterpt[2][1] = -1.4641016151377548; 375bbd56ea5SKarl Rupp 3769566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSWRA3PW, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0])); 377fe7e6d57SJed Brown } 378fe7e6d57SJed Brown { 3793ca35412SEmil Constantinescu PetscReal binterpt[4][3]; 380da80777bSKarl Rupp /*const PetscReal g = 4.3586652150845900e-01; Directly written in-place below */ 3819371c9d4SSatish Balay const PetscReal A[4][4] = 3829371c9d4SSatish Balay { 3839371c9d4SSatish Balay {0, 0, 0, 0}, 384fe7e6d57SJed Brown {8.7173304301691801e-01, 0, 0, 0}, 385fe7e6d57SJed Brown {8.4457060015369423e-01, -1.1299064236484185e-01, 0, 0}, 3869371c9d4SSatish Balay {0, 0, 1., 0} 3879371c9d4SSatish Balay }, 3889371c9d4SSatish Balay Gamma[4][4] = {{4.3586652150845900e-01, 0, 0, 0}, {-8.7173304301691801e-01, 4.3586652150845900e-01, 0, 0}, {-9.0338057013044082e-01, 5.4180672388095326e-02, 4.3586652150845900e-01, 0}, {2.4212380706095346e-01, -1.2232505839045147e+00, 5.4526025533510214e-01, 4.3586652150845900e-01}}, b[4] = {2.4212380706095346e-01, -1.2232505839045147e+00, 1.5452602553351020e+00, 4.3586652150845900e-01}, b2[4] = {3.7810903145819369e-01, -9.6042292212423178e-02, 5.0000000000000000e-01, 2.1793326075422950e-01}; 3893ca35412SEmil Constantinescu 3903ca35412SEmil Constantinescu binterpt[0][0] = 1.0564298455794094; 3913ca35412SEmil Constantinescu binterpt[1][0] = 2.296429974281067; 3923ca35412SEmil Constantinescu binterpt[2][0] = -1.307599564525376; 3933ca35412SEmil Constantinescu binterpt[3][0] = -1.045260255335102; 3943ca35412SEmil Constantinescu binterpt[0][1] = -1.3864882699759573; 3953ca35412SEmil Constantinescu binterpt[1][1] = -8.262611700275677; 3963ca35412SEmil Constantinescu binterpt[2][1] = 7.250979895056055; 3973ca35412SEmil Constantinescu binterpt[3][1] = 2.398120075195581; 3983ca35412SEmil Constantinescu binterpt[0][2] = 0.5721822314575016; 3993ca35412SEmil Constantinescu binterpt[1][2] = 4.742931142090097; 4003ca35412SEmil Constantinescu binterpt[2][2] = -4.398120075195578; 4013ca35412SEmil Constantinescu binterpt[3][2] = -0.9169932983520199; 4023ca35412SEmil Constantinescu 4039566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSWRA34PW2, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0])); 404e27a552bSJed Brown } 405ef3c5b88SJed Brown { 406da80777bSKarl Rupp /* const PetscReal g = 0.5; Directly written in-place below */ 4079371c9d4SSatish Balay const PetscReal A[4][4] = 4089371c9d4SSatish Balay { 4099371c9d4SSatish Balay {0, 0, 0, 0}, 410ef3c5b88SJed Brown {0, 0, 0, 0}, 411ef3c5b88SJed Brown {1., 0, 0, 0}, 4129371c9d4SSatish Balay {0.75, -0.25, 0.5, 0} 4139371c9d4SSatish Balay }, 4149371c9d4SSatish Balay Gamma[4][4] = {{0.5, 0, 0, 0}, {1., 0.5, 0, 0}, {-0.25, -0.25, 0.5, 0}, {1. / 12, 1. / 12, -2. / 3, 0.5}}, b[4] = {5. / 6, -1. / 6, -1. / 6, 0.5}, b2[4] = {0.75, -0.25, 0.5, 0}; 415bbd56ea5SKarl Rupp 4169566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSWRODAS3, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 0, NULL)); 417ef3c5b88SJed Brown } 418ef3c5b88SJed Brown { 419da80777bSKarl Rupp /*const PetscReal g = 0.43586652150845899941601945119356; Directly written in-place below */ 4209371c9d4SSatish Balay const PetscReal A[3][3] = 4219371c9d4SSatish Balay { 4229371c9d4SSatish Balay {0, 0, 0}, 423da80777bSKarl Rupp {0.43586652150845899941601945119356, 0, 0}, 4249371c9d4SSatish Balay {0.43586652150845899941601945119356, 0, 0} 4259371c9d4SSatish Balay }, 4269371c9d4SSatish Balay Gamma[3][3] = {{0.43586652150845899941601945119356, 0, 0}, {-0.19294655696029095575009695436041, 0.43586652150845899941601945119356, 0}, {0, 1.74927148125794685173529749738960, 0.43586652150845899941601945119356}}, b[3] = {-0.75457412385404315829818998646589, 1.94100407061964420292840123379419, -0.18642994676560104463021124732829}, b2[3] = {-1.53358745784149585370766523913002, 2.81745131148625772213931745457622, -0.28386385364476186843165221544619}; 4271f80e275SEmil Constantinescu 4281f80e275SEmil Constantinescu PetscReal binterpt[3][2]; 4291f80e275SEmil Constantinescu binterpt[0][0] = 3.793692883777660870425141387941; 4301f80e275SEmil Constantinescu binterpt[1][0] = -2.918692883777660870425141387941; 4311f80e275SEmil Constantinescu binterpt[2][0] = 0.125; 4321f80e275SEmil Constantinescu binterpt[0][1] = -0.725741064379812106687651020584; 4331f80e275SEmil Constantinescu binterpt[1][1] = 0.559074397713145440020984353917; 4341f80e275SEmil Constantinescu binterpt[2][1] = 0.16666666666666666666666666666667; 4351f80e275SEmil Constantinescu 4369566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSWSANDU3, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0])); 437ef3c5b88SJed Brown } 438b1c69cc3SEmil Constantinescu { 439da80777bSKarl Rupp /*const PetscReal s3 = PetscSqrtReal(3.),g = (3.0+s3)/6.0; 440da80777bSKarl Rupp * Direct evaluation: s3 = 1.732050807568877293527; 441da80777bSKarl Rupp * g = 0.7886751345948128822546; 442da80777bSKarl Rupp * Values are directly inserted below to ensure availability at compile time (compiler warnings otherwise...) */ 4439371c9d4SSatish Balay const PetscReal A[3][3] = 4449371c9d4SSatish Balay { 4459371c9d4SSatish Balay {0, 0, 0}, 446b1c69cc3SEmil Constantinescu {1, 0, 0}, 4479371c9d4SSatish Balay {0.25, 0.25, 0} 4489371c9d4SSatish Balay }, 4499371c9d4SSatish Balay Gamma[3][3] = {{0, 0, 0}, {(-3.0 - 1.732050807568877293527) / 6.0, 0.7886751345948128822546, 0}, {(-3.0 - 1.732050807568877293527) / 24.0, (-3.0 - 1.732050807568877293527) / 8.0, 0.7886751345948128822546}}, b[3] = {1. / 6., 1. / 6., 2. / 3.}, b2[3] = {1. / 4., 1. / 4., 1. / 2.}; 450c0cb691aSEmil Constantinescu PetscReal binterpt[3][2]; 451da80777bSKarl Rupp 452c0cb691aSEmil Constantinescu binterpt[0][0] = 0.089316397477040902157517886164709; 453c0cb691aSEmil Constantinescu binterpt[1][0] = -0.91068360252295909784248211383529; 454c0cb691aSEmil Constantinescu binterpt[2][0] = 1.8213672050459181956849642276706; 455c0cb691aSEmil Constantinescu binterpt[0][1] = 0.077350269189625764509148780501957; 456c0cb691aSEmil Constantinescu binterpt[1][1] = 1.077350269189625764509148780502; 457c0cb691aSEmil Constantinescu binterpt[2][1] = -1.1547005383792515290182975610039; 458bbd56ea5SKarl Rupp 4599566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSWASSP3P3S1C, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0])); 460b1c69cc3SEmil Constantinescu } 461b1c69cc3SEmil Constantinescu 462b1c69cc3SEmil Constantinescu { 4639371c9d4SSatish Balay const PetscReal A[4][4] = 4649371c9d4SSatish Balay { 4659371c9d4SSatish Balay {0, 0, 0, 0}, 466b1c69cc3SEmil Constantinescu {1. / 2., 0, 0, 0}, 467b1c69cc3SEmil Constantinescu {1. / 2., 1. / 2., 0, 0}, 4689371c9d4SSatish Balay {1. / 6., 1. / 6., 1. / 6., 0} 4699371c9d4SSatish Balay }, 4709371c9d4SSatish Balay Gamma[4][4] = {{1. / 2., 0, 0, 0}, {0.0, 1. / 4., 0, 0}, {-2., -2. / 3., 2. / 3., 0}, {1. / 2., 5. / 36., -2. / 9, 0}}, b[4] = {1. / 6., 1. / 6., 1. / 6., 1. / 2.}, b2[4] = {1. / 8., 3. / 4., 1. / 8., 0}; 471c0cb691aSEmil Constantinescu PetscReal binterpt[4][3]; 472da80777bSKarl Rupp 473c0cb691aSEmil Constantinescu binterpt[0][0] = 6.25; 474c0cb691aSEmil Constantinescu binterpt[1][0] = -30.25; 475c0cb691aSEmil Constantinescu binterpt[2][0] = 1.75; 476c0cb691aSEmil Constantinescu binterpt[3][0] = 23.25; 477c0cb691aSEmil Constantinescu binterpt[0][1] = -9.75; 478c0cb691aSEmil Constantinescu binterpt[1][1] = 58.75; 479c0cb691aSEmil Constantinescu binterpt[2][1] = -3.25; 480c0cb691aSEmil Constantinescu binterpt[3][1] = -45.75; 481c0cb691aSEmil Constantinescu binterpt[0][2] = 3.6666666666666666666666666666667; 482c0cb691aSEmil Constantinescu binterpt[1][2] = -28.333333333333333333333333333333; 483c0cb691aSEmil Constantinescu binterpt[2][2] = 1.6666666666666666666666666666667; 484c0cb691aSEmil Constantinescu binterpt[3][2] = 23.; 485bbd56ea5SKarl Rupp 4869566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSWLASSP3P4S2C, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0])); 487b1c69cc3SEmil Constantinescu } 488b1c69cc3SEmil Constantinescu 489b1c69cc3SEmil Constantinescu { 4909371c9d4SSatish Balay const PetscReal A[4][4] = 4919371c9d4SSatish Balay { 4929371c9d4SSatish Balay {0, 0, 0, 0}, 493b1c69cc3SEmil Constantinescu {1. / 2., 0, 0, 0}, 494b1c69cc3SEmil Constantinescu {1. / 2., 1. / 2., 0, 0}, 4959371c9d4SSatish Balay {1. / 6., 1. / 6., 1. / 6., 0} 4969371c9d4SSatish Balay }, 4979371c9d4SSatish Balay Gamma[4][4] = {{1. / 2., 0, 0, 0}, {0.0, 3. / 4., 0, 0}, {-2. / 3., -23. / 9., 2. / 9., 0}, {1. / 18., 65. / 108., -2. / 27, 0}}, b[4] = {1. / 6., 1. / 6., 1. / 6., 1. / 2.}, b2[4] = {3. / 16., 10. / 16., 3. / 16., 0}; 498c0cb691aSEmil Constantinescu PetscReal binterpt[4][3]; 499da80777bSKarl Rupp 500c0cb691aSEmil Constantinescu binterpt[0][0] = 1.6911764705882352941176470588235; 501c0cb691aSEmil Constantinescu binterpt[1][0] = 3.6813725490196078431372549019608; 502c0cb691aSEmil Constantinescu binterpt[2][0] = 0.23039215686274509803921568627451; 503c0cb691aSEmil Constantinescu binterpt[3][0] = -4.6029411764705882352941176470588; 504c0cb691aSEmil Constantinescu binterpt[0][1] = -0.95588235294117647058823529411765; 505c0cb691aSEmil Constantinescu binterpt[1][1] = -6.2401960784313725490196078431373; 506c0cb691aSEmil Constantinescu binterpt[2][1] = -0.31862745098039215686274509803922; 507c0cb691aSEmil Constantinescu binterpt[3][1] = 7.5147058823529411764705882352941; 508c0cb691aSEmil Constantinescu binterpt[0][2] = -0.56862745098039215686274509803922; 509c0cb691aSEmil Constantinescu binterpt[1][2] = 2.7254901960784313725490196078431; 510c0cb691aSEmil Constantinescu binterpt[2][2] = 0.25490196078431372549019607843137; 511c0cb691aSEmil Constantinescu binterpt[3][2] = -2.4117647058823529411764705882353; 512bbd56ea5SKarl Rupp 5139566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSWLLSSP3P4S2C, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0])); 514b1c69cc3SEmil Constantinescu } 515753f8adbSEmil Constantinescu 516753f8adbSEmil Constantinescu { 517753f8adbSEmil Constantinescu PetscReal A[4][4], Gamma[4][4], b[4], b2[4]; 5183ca35412SEmil Constantinescu PetscReal binterpt[4][3]; 519753f8adbSEmil Constantinescu 520753f8adbSEmil Constantinescu Gamma[0][0] = 0.4358665215084589994160194475295062513822671686978816; 5219371c9d4SSatish Balay Gamma[0][1] = 0; 5229371c9d4SSatish Balay Gamma[0][2] = 0; 5239371c9d4SSatish Balay Gamma[0][3] = 0; 524753f8adbSEmil Constantinescu Gamma[1][0] = -1.997527830934941248426324674704153457289527280554476; 525753f8adbSEmil Constantinescu Gamma[1][1] = 0.4358665215084589994160194475295062513822671686978816; 5269371c9d4SSatish Balay Gamma[1][2] = 0; 5279371c9d4SSatish Balay Gamma[1][3] = 0; 528753f8adbSEmil Constantinescu Gamma[2][0] = -1.007948511795029620852002345345404191008352770119903; 529753f8adbSEmil Constantinescu Gamma[2][1] = -0.004648958462629345562774289390054679806993396798458131; 530753f8adbSEmil Constantinescu Gamma[2][2] = 0.4358665215084589994160194475295062513822671686978816; 53105e8e825SJed Brown Gamma[2][3] = 0; 532753f8adbSEmil Constantinescu Gamma[3][0] = -0.6685429734233467180451604600279552604364311322650783; 533753f8adbSEmil Constantinescu Gamma[3][1] = 0.6056625986449338476089525334450053439525178740492984; 534753f8adbSEmil Constantinescu Gamma[3][2] = -0.9717899277217721234705114616271378792182450260943198; 535753f8adbSEmil Constantinescu Gamma[3][3] = 0; 536753f8adbSEmil Constantinescu 5379371c9d4SSatish Balay A[0][0] = 0; 5389371c9d4SSatish Balay A[0][1] = 0; 5399371c9d4SSatish Balay A[0][2] = 0; 5409371c9d4SSatish Balay A[0][3] = 0; 541753f8adbSEmil Constantinescu A[1][0] = 0.8717330430169179988320388950590125027645343373957631; 5429371c9d4SSatish Balay A[1][1] = 0; 5439371c9d4SSatish Balay A[1][2] = 0; 5449371c9d4SSatish Balay A[1][3] = 0; 545753f8adbSEmil Constantinescu A[2][0] = 0.5275890119763004115618079766722914408876108660811028; 546753f8adbSEmil Constantinescu A[2][1] = 0.07241098802369958843819203208518599088698057726988732; 5479371c9d4SSatish Balay A[2][2] = 0; 5489371c9d4SSatish Balay A[2][3] = 0; 549753f8adbSEmil Constantinescu A[3][0] = 0.3990960076760701320627260685975778145384666450351314; 550753f8adbSEmil Constantinescu A[3][1] = -0.4375576546135194437228463747348862825846903771419953; 551753f8adbSEmil Constantinescu A[3][2] = 1.038461646937449311660120300601880176655352737312713; 55205e8e825SJed Brown A[3][3] = 0; 553753f8adbSEmil Constantinescu 554753f8adbSEmil Constantinescu b[0] = 0.1876410243467238251612921333138006734899663569186926; 555753f8adbSEmil Constantinescu b[1] = -0.5952974735769549480478230473706443582188442040780541; 556753f8adbSEmil Constantinescu b[2] = 0.9717899277217721234705114616271378792182450260943198; 557753f8adbSEmil Constantinescu b[3] = 0.4358665215084589994160194475295062513822671686978816; 558753f8adbSEmil Constantinescu 559753f8adbSEmil Constantinescu b2[0] = 0.2147402862233891404862383521089097657790734483804460; 560753f8adbSEmil Constantinescu b2[1] = -0.4851622638849390928209050538171743017757490232519684; 561753f8adbSEmil Constantinescu b2[2] = 0.8687250025203875511662123688667549217531982787600080; 562753f8adbSEmil Constantinescu b2[3] = 0.4016969751411624011684543450940068201770721128357014; 563753f8adbSEmil Constantinescu 5643ca35412SEmil Constantinescu binterpt[0][0] = 2.2565812720167954547104627844105; 5653ca35412SEmil Constantinescu binterpt[1][0] = 1.349166413351089573796243820819; 5663ca35412SEmil Constantinescu binterpt[2][0] = -2.4695174540533503758652847586647; 5673ca35412SEmil Constantinescu binterpt[3][0] = -0.13623023131453465264142184656474; 5683ca35412SEmil Constantinescu binterpt[0][1] = -3.0826699111559187902922463354557; 5693ca35412SEmil Constantinescu binterpt[1][1] = -2.4689115685996042534544925650515; 5703ca35412SEmil Constantinescu binterpt[2][1] = 5.7428279814696677152129332773553; 5713ca35412SEmil Constantinescu binterpt[3][1] = -0.19124650171414467146619437684812; 5723ca35412SEmil Constantinescu binterpt[0][2] = 1.0137296634858471607430756831148; 5733ca35412SEmil Constantinescu binterpt[1][2] = 0.52444768167155973161042570784064; 5743ca35412SEmil Constantinescu binterpt[2][2] = -2.3015205996945452158771370439586; 5753ca35412SEmil Constantinescu binterpt[3][2] = 0.76334325453713832352363565300308; 576f4aed992SEmil Constantinescu 5779566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSWARK3, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0])); 578753f8adbSEmil Constantinescu } 5799566063dSJacob Faibussowitsch PetscCall(TSRosWRegisterRos4(TSROSWGRK4T, 0.231, PETSC_DEFAULT, PETSC_DEFAULT, 0, -0.1282612945269037e+01)); 5809566063dSJacob Faibussowitsch PetscCall(TSRosWRegisterRos4(TSROSWSHAMP4, 0.5, PETSC_DEFAULT, PETSC_DEFAULT, 0, 125. / 108.)); 5819566063dSJacob Faibussowitsch PetscCall(TSRosWRegisterRos4(TSROSWVELDD4, 0.22570811482256823492, PETSC_DEFAULT, PETSC_DEFAULT, 0, -1.355958941201148)); 5829566063dSJacob Faibussowitsch PetscCall(TSRosWRegisterRos4(TSROSW4L, 0.57282, PETSC_DEFAULT, PETSC_DEFAULT, 0, -1.093502252409163)); 5833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 584e27a552bSJed Brown } 585e27a552bSJed Brown 586e27a552bSJed Brown /*@C 587bcf0153eSBarry Smith TSRosWRegisterDestroy - Frees the list of schemes that were registered by `TSRosWRegister()`. 588e27a552bSJed Brown 589e27a552bSJed Brown Not Collective 590e27a552bSJed Brown 591e27a552bSJed Brown Level: advanced 592e27a552bSJed Brown 593*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSRosWRegister()`, `TSRosWRegisterAll()` 594e27a552bSJed Brown @*/ 595d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWRegisterDestroy(void) 596d71ae5a4SJacob Faibussowitsch { 59761692a83SJed Brown RosWTableauLink link; 598e27a552bSJed Brown 599e27a552bSJed Brown PetscFunctionBegin; 60061692a83SJed Brown while ((link = RosWTableauList)) { 60161692a83SJed Brown RosWTableau t = &link->tab; 60261692a83SJed Brown RosWTableauList = link->next; 6039566063dSJacob Faibussowitsch PetscCall(PetscFree5(t->A, t->Gamma, t->b, t->ASum, t->GammaSum)); 6049566063dSJacob Faibussowitsch PetscCall(PetscFree5(t->At, t->bt, t->GammaInv, t->GammaZeroDiag, t->GammaExplicitCorr)); 6059566063dSJacob Faibussowitsch PetscCall(PetscFree2(t->bembed, t->bembedt)); 6069566063dSJacob Faibussowitsch PetscCall(PetscFree(t->binterpt)); 6079566063dSJacob Faibussowitsch PetscCall(PetscFree(t->name)); 6089566063dSJacob Faibussowitsch PetscCall(PetscFree(link)); 609e27a552bSJed Brown } 610e27a552bSJed Brown TSRosWRegisterAllCalled = PETSC_FALSE; 6113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 612e27a552bSJed Brown } 613e27a552bSJed Brown 614e27a552bSJed Brown /*@C 615bcf0153eSBarry Smith TSRosWInitializePackage - This function initializes everything in the `TSROSW` package. It is called 616bcf0153eSBarry Smith from `TSInitializePackage()`. 617e27a552bSJed Brown 618e27a552bSJed Brown Level: developer 619e27a552bSJed Brown 620*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `PetscInitialize()`, `TSRosWFinalizePackage()` 621e27a552bSJed Brown @*/ 622d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWInitializePackage(void) 623d71ae5a4SJacob Faibussowitsch { 624e27a552bSJed Brown PetscFunctionBegin; 6253ba16761SJacob Faibussowitsch if (TSRosWPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS); 626e27a552bSJed Brown TSRosWPackageInitialized = PETSC_TRUE; 6279566063dSJacob Faibussowitsch PetscCall(TSRosWRegisterAll()); 6289566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSRosWFinalizePackage)); 6293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 630e27a552bSJed Brown } 631e27a552bSJed Brown 632e27a552bSJed Brown /*@C 633bcf0153eSBarry Smith TSRosWFinalizePackage - This function destroys everything in the `TSROSW` package. It is 634bcf0153eSBarry Smith called from `PetscFinalize()`. 635e27a552bSJed Brown 636e27a552bSJed Brown Level: developer 637e27a552bSJed Brown 638*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `PetscFinalize()`, `TSRosWInitializePackage()` 639e27a552bSJed Brown @*/ 640d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWFinalizePackage(void) 641d71ae5a4SJacob Faibussowitsch { 642e27a552bSJed Brown PetscFunctionBegin; 643e27a552bSJed Brown TSRosWPackageInitialized = PETSC_FALSE; 6449566063dSJacob Faibussowitsch PetscCall(TSRosWRegisterDestroy()); 6453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 646e27a552bSJed Brown } 647e27a552bSJed Brown 648e27a552bSJed Brown /*@C 649bcf0153eSBarry Smith TSRosWRegister - register a `TSROSW`, Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 650e27a552bSJed Brown 651e27a552bSJed Brown Not Collective, but the same schemes should be registered on all processes on which they will be used 652e27a552bSJed Brown 653e27a552bSJed Brown Input Parameters: 654e27a552bSJed Brown + name - identifier for method 655e27a552bSJed Brown . order - approximation order of method 656e27a552bSJed Brown . s - number of stages, this is the dimension of the matrices below 65761692a83SJed Brown . A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular 65861692a83SJed Brown . Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal 659fe7e6d57SJed Brown . b - Step completion table (dimension s) 6600298fd71SBarry Smith . bembed - Step completion table for a scheme of order one less (dimension s, NULL if no embedded scheme is available) 661f4aed992SEmil Constantinescu . pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt 66242faf41dSJed Brown - binterpt - Coefficients of the interpolation formula (dimension s*pinterp) 663e27a552bSJed Brown 664e27a552bSJed Brown Level: advanced 665e27a552bSJed Brown 666bcf0153eSBarry Smith Note: 667bcf0153eSBarry Smith Several Rosenbrock W methods are provided, this function is only needed to create new methods. 668bcf0153eSBarry Smith 669*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW` 670e27a552bSJed Brown @*/ 671d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWRegister(TSRosWType name, PetscInt order, PetscInt s, const PetscReal A[], const PetscReal Gamma[], const PetscReal b[], const PetscReal bembed[], PetscInt pinterp, const PetscReal binterpt[]) 672d71ae5a4SJacob Faibussowitsch { 67361692a83SJed Brown RosWTableauLink link; 67461692a83SJed Brown RosWTableau t; 67561692a83SJed Brown PetscInt i, j, k; 67661692a83SJed Brown PetscScalar *GammaInv; 677e27a552bSJed Brown 678e27a552bSJed Brown PetscFunctionBegin; 679fe7e6d57SJed Brown PetscValidCharPointer(name, 1); 680dadcf809SJacob Faibussowitsch PetscValidRealPointer(A, 4); 681dadcf809SJacob Faibussowitsch PetscValidRealPointer(Gamma, 5); 682dadcf809SJacob Faibussowitsch PetscValidRealPointer(b, 6); 683dadcf809SJacob Faibussowitsch if (bembed) PetscValidRealPointer(bembed, 7); 684fe7e6d57SJed Brown 6859566063dSJacob Faibussowitsch PetscCall(TSRosWInitializePackage()); 6869566063dSJacob Faibussowitsch PetscCall(PetscNew(&link)); 687e27a552bSJed Brown t = &link->tab; 6889566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &t->name)); 689e27a552bSJed Brown t->order = order; 690e27a552bSJed Brown t->s = s; 6919566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(s * s, &t->A, s * s, &t->Gamma, s, &t->b, s, &t->ASum, s, &t->GammaSum)); 6929566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(s * s, &t->At, s, &t->bt, s * s, &t->GammaInv, s, &t->GammaZeroDiag, s * s, &t->GammaExplicitCorr)); 6939566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->A, A, s * s)); 6949566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->Gamma, Gamma, s * s)); 6959566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->GammaExplicitCorr, Gamma, s * s)); 6969566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->b, b, s)); 697fe7e6d57SJed Brown if (bembed) { 6989566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(s, &t->bembed, s, &t->bembedt)); 6999566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->bembed, bembed, s)); 700fe7e6d57SJed Brown } 70161692a83SJed Brown for (i = 0; i < s; i++) { 70261692a83SJed Brown t->ASum[i] = 0; 70361692a83SJed Brown t->GammaSum[i] = 0; 70461692a83SJed Brown for (j = 0; j < s; j++) { 70561692a83SJed Brown t->ASum[i] += A[i * s + j]; 706fe7e6d57SJed Brown t->GammaSum[i] += Gamma[i * s + j]; 70761692a83SJed Brown } 70861692a83SJed Brown } 7099566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s * s, &GammaInv)); /* Need to use Scalar for inverse, then convert back to Real */ 71061692a83SJed Brown for (i = 0; i < s * s; i++) GammaInv[i] = Gamma[i]; 711fd96d5b0SEmil Constantinescu for (i = 0; i < s; i++) { 712fd96d5b0SEmil Constantinescu if (Gamma[i * s + i] == 0.0) { 713fd96d5b0SEmil Constantinescu GammaInv[i * s + i] = 1.0; 714c17803e7SJed Brown t->GammaZeroDiag[i] = PETSC_TRUE; 715fd96d5b0SEmil Constantinescu } else { 716c17803e7SJed Brown t->GammaZeroDiag[i] = PETSC_FALSE; 717fd96d5b0SEmil Constantinescu } 718fd96d5b0SEmil Constantinescu } 719fd96d5b0SEmil Constantinescu 72061692a83SJed Brown switch (s) { 721d71ae5a4SJacob Faibussowitsch case 1: 722d71ae5a4SJacob Faibussowitsch GammaInv[0] = 1. / GammaInv[0]; 723d71ae5a4SJacob Faibussowitsch break; 724d71ae5a4SJacob Faibussowitsch case 2: 725d71ae5a4SJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_2(GammaInv, 0, PETSC_FALSE, NULL)); 726d71ae5a4SJacob Faibussowitsch break; 727d71ae5a4SJacob Faibussowitsch case 3: 728d71ae5a4SJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_3(GammaInv, 0, PETSC_FALSE, NULL)); 729d71ae5a4SJacob Faibussowitsch break; 730d71ae5a4SJacob Faibussowitsch case 4: 731d71ae5a4SJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_4(GammaInv, 0, PETSC_FALSE, NULL)); 732d71ae5a4SJacob Faibussowitsch break; 73361692a83SJed Brown case 5: { 73461692a83SJed Brown PetscInt ipvt5[5]; 73561692a83SJed Brown MatScalar work5[5 * 5]; 7369371c9d4SSatish Balay PetscCall(PetscKernel_A_gets_inverse_A_5(GammaInv, ipvt5, work5, 0, PETSC_FALSE, NULL)); 7379371c9d4SSatish Balay break; 73861692a83SJed Brown } 739d71ae5a4SJacob Faibussowitsch case 6: 740d71ae5a4SJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_6(GammaInv, 0, PETSC_FALSE, NULL)); 741d71ae5a4SJacob Faibussowitsch break; 742d71ae5a4SJacob Faibussowitsch case 7: 743d71ae5a4SJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_7(GammaInv, 0, PETSC_FALSE, NULL)); 744d71ae5a4SJacob Faibussowitsch break; 745d71ae5a4SJacob Faibussowitsch default: 746d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not implemented for %" PetscInt_FMT " stages", s); 74761692a83SJed Brown } 74861692a83SJed Brown for (i = 0; i < s * s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]); 7499566063dSJacob Faibussowitsch PetscCall(PetscFree(GammaInv)); 75043b21953SEmil Constantinescu 75143b21953SEmil Constantinescu for (i = 0; i < s; i++) { 75243b21953SEmil Constantinescu for (k = 0; k < i + 1; k++) { 75343b21953SEmil Constantinescu t->GammaExplicitCorr[i * s + k] = (t->GammaExplicitCorr[i * s + k]) * (t->GammaInv[k * s + k]); 754ad540459SPierre Jolivet for (j = k + 1; j < i + 1; j++) t->GammaExplicitCorr[i * s + k] += (t->GammaExplicitCorr[i * s + j]) * (t->GammaInv[j * s + k]); 75543b21953SEmil Constantinescu } 75643b21953SEmil Constantinescu } 75743b21953SEmil Constantinescu 75861692a83SJed Brown for (i = 0; i < s; i++) { 75961692a83SJed Brown for (j = 0; j < s; j++) { 76061692a83SJed Brown t->At[i * s + j] = 0; 761ad540459SPierre Jolivet for (k = 0; k < s; k++) t->At[i * s + j] += t->A[i * s + k] * t->GammaInv[k * s + j]; 76261692a83SJed Brown } 76361692a83SJed Brown t->bt[i] = 0; 764ad540459SPierre Jolivet for (j = 0; j < s; j++) t->bt[i] += t->b[j] * t->GammaInv[j * s + i]; 765fe7e6d57SJed Brown if (bembed) { 766fe7e6d57SJed Brown t->bembedt[i] = 0; 767ad540459SPierre Jolivet for (j = 0; j < s; j++) t->bembedt[i] += t->bembed[j] * t->GammaInv[j * s + i]; 768fe7e6d57SJed Brown } 76961692a83SJed Brown } 7708d59e960SJed Brown t->ccfl = 1.0; /* Fix this */ 7718d59e960SJed Brown 772f4aed992SEmil Constantinescu t->pinterp = pinterp; 7739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s * pinterp, &t->binterpt)); 7749566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->binterpt, binterpt, s * pinterp)); 77561692a83SJed Brown link->next = RosWTableauList; 77661692a83SJed Brown RosWTableauList = link; 7773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 778e27a552bSJed Brown } 779e27a552bSJed Brown 78042faf41dSJed Brown /*@C 781fd292e60Sprj- TSRosWRegisterRos4 - register a fourth order Rosenbrock scheme by providing parameter choices 78242faf41dSJed Brown 78342faf41dSJed Brown Not Collective, but the same schemes should be registered on all processes on which they will be used 78442faf41dSJed Brown 78542faf41dSJed Brown Input Parameters: 78642faf41dSJed Brown + name - identifier for method 78742faf41dSJed Brown . gamma - leading coefficient (diagonal entry) 78842faf41dSJed Brown . a2 - design parameter, see Table 7.2 of Hairer&Wanner 78942faf41dSJed Brown . a3 - design parameter or PETSC_DEFAULT to satisfy one of the order five conditions (Eq 7.22) 79042faf41dSJed Brown . b3 - design parameter, see Table 7.2 of Hairer&Wanner 791a2b725a8SWilliam Gropp - e4 - design parameter for embedded method, see coefficient E4 in ros4.f code from Hairer 79242faf41dSJed Brown 793bcf0153eSBarry Smith Level: developer 794bcf0153eSBarry Smith 79542faf41dSJed Brown Notes: 79642faf41dSJed Brown This routine encodes the design of fourth order Rosenbrock methods as described in Hairer and Wanner volume 2. 79742faf41dSJed Brown It is used here to implement several methods from the book and can be used to experiment with new methods. 79842faf41dSJed Brown It was written this way instead of by copying coefficients in order to provide better than double precision satisfaction of the order conditions. 79942faf41dSJed Brown 800*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSRosW`, `TSRosWRegister()` 80142faf41dSJed Brown @*/ 802d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWRegisterRos4(TSRosWType name, PetscReal gamma, PetscReal a2, PetscReal a3, PetscReal b3, PetscReal e4) 803d71ae5a4SJacob Faibussowitsch { 80442faf41dSJed Brown /* Declare numeric constants so they can be quad precision without being truncated at double */ 8059371c9d4SSatish Balay const PetscReal one = 1, two = 2, three = 3, four = 4, five = 5, six = 6, eight = 8, twelve = 12, twenty = 20, twentyfour = 24, p32 = one / six - gamma + gamma * gamma, p42 = one / eight - gamma / three, p43 = one / twelve - gamma / three, p44 = one / twentyfour - gamma / two + three / two * gamma * gamma - gamma * gamma * gamma, p56 = one / twenty - gamma / four; 80642faf41dSJed Brown PetscReal a4, a32, a42, a43, b1, b2, b4, beta2p, beta3p, beta4p, beta32, beta42, beta43, beta32beta2p, beta4jbetajp; 80742faf41dSJed Brown PetscReal A[4][4], Gamma[4][4], b[4], bm[4]; 80842faf41dSJed Brown PetscScalar M[3][3], rhs[3]; 80942faf41dSJed Brown 81042faf41dSJed Brown PetscFunctionBegin; 81142faf41dSJed Brown /* Step 1: choose Gamma (input) */ 81242faf41dSJed Brown /* Step 2: choose a2,a3,a4; b1,b2,b3,b4 to satisfy order conditions */ 81313bcc0bdSJacob Faibussowitsch if (a3 == (PetscReal)PETSC_DEFAULT) a3 = (one / five - a2 / four) / (one / four - a2 / three); /* Eq 7.22 */ 81442faf41dSJed Brown a4 = a3; /* consequence of 7.20 */ 81542faf41dSJed Brown 81642faf41dSJed Brown /* Solve order conditions 7.15a, 7.15c, 7.15e */ 8179371c9d4SSatish Balay M[0][0] = one; 8189371c9d4SSatish Balay M[0][1] = one; 8199371c9d4SSatish Balay M[0][2] = one; /* 7.15a */ 8209371c9d4SSatish Balay M[1][0] = 0.0; 8219371c9d4SSatish Balay M[1][1] = a2 * a2; 8229371c9d4SSatish Balay M[1][2] = a4 * a4; /* 7.15c */ 8239371c9d4SSatish Balay M[2][0] = 0.0; 8249371c9d4SSatish Balay M[2][1] = a2 * a2 * a2; 8259371c9d4SSatish Balay M[2][2] = a4 * a4 * a4; /* 7.15e */ 82642faf41dSJed Brown rhs[0] = one - b3; 82742faf41dSJed Brown rhs[1] = one / three - a3 * a3 * b3; 82842faf41dSJed Brown rhs[2] = one / four - a3 * a3 * a3 * b3; 8299566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_3(&M[0][0], 0, PETSC_FALSE, NULL)); 83042faf41dSJed Brown b1 = PetscRealPart(M[0][0] * rhs[0] + M[0][1] * rhs[1] + M[0][2] * rhs[2]); 83142faf41dSJed Brown b2 = PetscRealPart(M[1][0] * rhs[0] + M[1][1] * rhs[1] + M[1][2] * rhs[2]); 83242faf41dSJed Brown b4 = PetscRealPart(M[2][0] * rhs[0] + M[2][1] * rhs[1] + M[2][2] * rhs[2]); 83342faf41dSJed Brown 83442faf41dSJed Brown /* Step 3 */ 83542faf41dSJed Brown beta43 = (p56 - a2 * p43) / (b4 * a3 * a3 * (a3 - a2)); /* 7.21 */ 83642faf41dSJed Brown beta32beta2p = p44 / (b4 * beta43); /* 7.15h */ 83742faf41dSJed Brown beta4jbetajp = (p32 - b3 * beta32beta2p) / b4; 8389371c9d4SSatish Balay M[0][0] = b2; 8399371c9d4SSatish Balay M[0][1] = b3; 8409371c9d4SSatish Balay M[0][2] = b4; 8419371c9d4SSatish Balay M[1][0] = a4 * a4 * beta32beta2p - a3 * a3 * beta4jbetajp; 8429371c9d4SSatish Balay M[1][1] = a2 * a2 * beta4jbetajp; 8439371c9d4SSatish Balay M[1][2] = -a2 * a2 * beta32beta2p; 8449371c9d4SSatish Balay M[2][0] = b4 * beta43 * a3 * a3 - p43; 8459371c9d4SSatish Balay M[2][1] = -b4 * beta43 * a2 * a2; 8469371c9d4SSatish Balay M[2][2] = 0; 8479371c9d4SSatish Balay rhs[0] = one / two - gamma; 8489371c9d4SSatish Balay rhs[1] = 0; 8499371c9d4SSatish Balay rhs[2] = -a2 * a2 * p32; 8509566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_3(&M[0][0], 0, PETSC_FALSE, NULL)); 85142faf41dSJed Brown beta2p = PetscRealPart(M[0][0] * rhs[0] + M[0][1] * rhs[1] + M[0][2] * rhs[2]); 85242faf41dSJed Brown beta3p = PetscRealPart(M[1][0] * rhs[0] + M[1][1] * rhs[1] + M[1][2] * rhs[2]); 85342faf41dSJed Brown beta4p = PetscRealPart(M[2][0] * rhs[0] + M[2][1] * rhs[1] + M[2][2] * rhs[2]); 85442faf41dSJed Brown 85542faf41dSJed Brown /* Step 4: back-substitute */ 85642faf41dSJed Brown beta32 = beta32beta2p / beta2p; 85742faf41dSJed Brown beta42 = (beta4jbetajp - beta43 * beta3p) / beta2p; 85842faf41dSJed Brown 85942faf41dSJed Brown /* Step 5: 7.15f and 7.20, then 7.16 */ 86042faf41dSJed Brown a43 = 0; 86142faf41dSJed Brown a32 = p42 / (b3 * a3 * beta2p + b4 * a4 * beta2p); 86242faf41dSJed Brown a42 = a32; 86342faf41dSJed Brown 8649371c9d4SSatish Balay A[0][0] = 0; 8659371c9d4SSatish Balay A[0][1] = 0; 8669371c9d4SSatish Balay A[0][2] = 0; 8679371c9d4SSatish Balay A[0][3] = 0; 8689371c9d4SSatish Balay A[1][0] = a2; 8699371c9d4SSatish Balay A[1][1] = 0; 8709371c9d4SSatish Balay A[1][2] = 0; 8719371c9d4SSatish Balay A[1][3] = 0; 8729371c9d4SSatish Balay A[2][0] = a3 - a32; 8739371c9d4SSatish Balay A[2][1] = a32; 8749371c9d4SSatish Balay A[2][2] = 0; 8759371c9d4SSatish Balay A[2][3] = 0; 8769371c9d4SSatish Balay A[3][0] = a4 - a43 - a42; 8779371c9d4SSatish Balay A[3][1] = a42; 8789371c9d4SSatish Balay A[3][2] = a43; 8799371c9d4SSatish Balay A[3][3] = 0; 8809371c9d4SSatish Balay Gamma[0][0] = gamma; 8819371c9d4SSatish Balay Gamma[0][1] = 0; 8829371c9d4SSatish Balay Gamma[0][2] = 0; 8839371c9d4SSatish Balay Gamma[0][3] = 0; 8849371c9d4SSatish Balay Gamma[1][0] = beta2p - A[1][0]; 8859371c9d4SSatish Balay Gamma[1][1] = gamma; 8869371c9d4SSatish Balay Gamma[1][2] = 0; 8879371c9d4SSatish Balay Gamma[1][3] = 0; 8889371c9d4SSatish Balay Gamma[2][0] = beta3p - beta32 - A[2][0]; 8899371c9d4SSatish Balay Gamma[2][1] = beta32 - A[2][1]; 8909371c9d4SSatish Balay Gamma[2][2] = gamma; 8919371c9d4SSatish Balay Gamma[2][3] = 0; 8929371c9d4SSatish Balay Gamma[3][0] = beta4p - beta42 - beta43 - A[3][0]; 8939371c9d4SSatish Balay Gamma[3][1] = beta42 - A[3][1]; 8949371c9d4SSatish Balay Gamma[3][2] = beta43 - A[3][2]; 8959371c9d4SSatish Balay Gamma[3][3] = gamma; 8969371c9d4SSatish Balay b[0] = b1; 8979371c9d4SSatish Balay b[1] = b2; 8989371c9d4SSatish Balay b[2] = b3; 8999371c9d4SSatish Balay b[3] = b4; 90042faf41dSJed Brown 90142faf41dSJed Brown /* Construct embedded formula using given e4. We are solving Equation 7.18. */ 90242faf41dSJed Brown bm[3] = b[3] - e4 * gamma; /* using definition of E4 */ 90342faf41dSJed Brown bm[2] = (p32 - beta4jbetajp * bm[3]) / (beta32 * beta2p); /* fourth row of 7.18 */ 90442faf41dSJed Brown bm[1] = (one / two - gamma - beta3p * bm[2] - beta4p * bm[3]) / beta2p; /* second row */ 90542faf41dSJed Brown bm[0] = one - bm[1] - bm[2] - bm[3]; /* first row */ 90642faf41dSJed Brown 90742faf41dSJed Brown { 90842faf41dSJed Brown const PetscReal misfit = a2 * a2 * bm[1] + a3 * a3 * bm[2] + a4 * a4 * bm[3] - one / three; 9093c633725SBarry Smith PetscCheck(PetscAbs(misfit) <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_SUP, "Assumptions violated, could not construct a third order embedded method"); 91042faf41dSJed Brown } 9119566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(name, 4, 4, &A[0][0], &Gamma[0][0], b, bm, 0, NULL)); 9123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 91342faf41dSJed Brown } 91442faf41dSJed Brown 9151c3436cfSJed Brown /* 9161c3436cfSJed Brown The step completion formula is 9171c3436cfSJed Brown 9181c3436cfSJed Brown x1 = x0 + b^T Y 9191c3436cfSJed Brown 9201c3436cfSJed Brown where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been 9211c3436cfSJed Brown updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write 9221c3436cfSJed Brown 9231c3436cfSJed Brown x1e = x0 + be^T Y 9241c3436cfSJed Brown = x1 - b^T Y + be^T Y 9251c3436cfSJed Brown = x1 + (be - b)^T Y 9261c3436cfSJed Brown 9271c3436cfSJed Brown so we can evaluate the method of different order even after the step has been optimistically completed. 9281c3436cfSJed Brown */ 929d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateStep_RosW(TS ts, PetscInt order, Vec U, PetscBool *done) 930d71ae5a4SJacob Faibussowitsch { 9311c3436cfSJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 9321c3436cfSJed Brown RosWTableau tab = ros->tableau; 9331c3436cfSJed Brown PetscScalar *w = ros->work; 9341c3436cfSJed Brown PetscInt i; 9351c3436cfSJed Brown 9361c3436cfSJed Brown PetscFunctionBegin; 9371c3436cfSJed Brown if (order == tab->order) { 938108c343cSJed Brown if (ros->status == TS_STEP_INCOMPLETE) { /* Use standard completion formula */ 9399566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, U)); 940de19f811SJed Brown for (i = 0; i < tab->s; i++) w[i] = tab->bt[i]; 9419566063dSJacob Faibussowitsch PetscCall(VecMAXPY(U, tab->s, w, ros->Y)); 9429566063dSJacob Faibussowitsch } else PetscCall(VecCopy(ts->vec_sol, U)); 9431c3436cfSJed Brown if (done) *done = PETSC_TRUE; 9443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9451c3436cfSJed Brown } else if (order == tab->order - 1) { 9461c3436cfSJed Brown if (!tab->bembedt) goto unavailable; 947108c343cSJed Brown if (ros->status == TS_STEP_INCOMPLETE) { /* Use embedded completion formula */ 9489566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, U)); 949de19f811SJed Brown for (i = 0; i < tab->s; i++) w[i] = tab->bembedt[i]; 9509566063dSJacob Faibussowitsch PetscCall(VecMAXPY(U, tab->s, w, ros->Y)); 951108c343cSJed Brown } else { /* Use rollback-and-recomplete formula (bembedt - bt) */ 952108c343cSJed Brown for (i = 0; i < tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i]; 9539566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, U)); 9549566063dSJacob Faibussowitsch PetscCall(VecMAXPY(U, tab->s, w, ros->Y)); 9551c3436cfSJed Brown } 9561c3436cfSJed Brown if (done) *done = PETSC_TRUE; 9573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9581c3436cfSJed Brown } 9591c3436cfSJed Brown unavailable: 9601c3436cfSJed Brown if (done) *done = PETSC_FALSE; 9619371c9d4SSatish Balay else 9629371c9d4SSatish Balay SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Rosenbrock-W '%s' of order %" PetscInt_FMT " cannot evaluate step at order %" PetscInt_FMT ". Consider using -ts_adapt_type none or a different method that has an embedded estimate.", tab->name, 9639371c9d4SSatish Balay tab->order, order); 9643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9651c3436cfSJed Brown } 9661c3436cfSJed Brown 967d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRollBack_RosW(TS ts) 968d71ae5a4SJacob Faibussowitsch { 96924655328SShri TS_RosW *ros = (TS_RosW *)ts->data; 97024655328SShri 97124655328SShri PetscFunctionBegin; 9729566063dSJacob Faibussowitsch PetscCall(VecCopy(ros->vec_sol_prev, ts->vec_sol)); 9733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 97424655328SShri } 97524655328SShri 976d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_RosW(TS ts) 977d71ae5a4SJacob Faibussowitsch { 97861692a83SJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 97961692a83SJed Brown RosWTableau tab = ros->tableau; 980e27a552bSJed Brown const PetscInt s = tab->s; 9811c3436cfSJed Brown const PetscReal *At = tab->At, *Gamma = tab->Gamma, *ASum = tab->ASum, *GammaInv = tab->GammaInv; 9820feba352SEmil Constantinescu const PetscReal *GammaExplicitCorr = tab->GammaExplicitCorr; 983c17803e7SJed Brown const PetscBool *GammaZeroDiag = tab->GammaZeroDiag; 98461692a83SJed Brown PetscScalar *w = ros->work; 9857d4bf2deSEmil Constantinescu Vec *Y = ros->Y, Ydot = ros->Ydot, Zdot = ros->Zdot, Zstage = ros->Zstage; 986e27a552bSJed Brown SNES snes; 9871c3436cfSJed Brown TSAdapt adapt; 988fecfb714SLisandro Dalcin PetscInt i, j, its, lits; 989be5899b3SLisandro Dalcin PetscInt rejections = 0; 990b39943a6SLisandro Dalcin PetscBool stageok, accept = PETSC_TRUE; 991b39943a6SLisandro Dalcin PetscReal next_time_step = ts->time_step; 992f7f07198SBarry Smith PetscInt lag; 993e27a552bSJed Brown 994e27a552bSJed Brown PetscFunctionBegin; 99548a46eb9SPierre Jolivet if (!ts->steprollback) PetscCall(VecCopy(ts->vec_sol, ros->vec_sol_prev)); 996e27a552bSJed Brown 997b39943a6SLisandro Dalcin ros->status = TS_STEP_INCOMPLETE; 998b39943a6SLisandro Dalcin while (!ts->reason && ros->status != TS_STEP_COMPLETE) { 9991c3436cfSJed Brown const PetscReal h = ts->time_step; 1000e27a552bSJed Brown for (i = 0; i < s; i++) { 10011c3436cfSJed Brown ros->stage_time = ts->ptime + h * ASum[i]; 10029566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, ros->stage_time)); 1003c17803e7SJed Brown if (GammaZeroDiag[i]) { 1004c17803e7SJed Brown ros->stage_explicit = PETSC_TRUE; 1005b296d7d5SJed Brown ros->scoeff = 1.; 1006c17803e7SJed Brown } else { 1007c17803e7SJed Brown ros->stage_explicit = PETSC_FALSE; 1008b296d7d5SJed Brown ros->scoeff = 1. / Gamma[i * s + i]; 1009fd96d5b0SEmil Constantinescu } 101061692a83SJed Brown 10119566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, Zstage)); 1012de19f811SJed Brown for (j = 0; j < i; j++) w[j] = At[i * s + j]; 10139566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Zstage, i, w, Y)); 101461692a83SJed Brown 101561692a83SJed Brown for (j = 0; j < i; j++) w[j] = 1. / h * GammaInv[i * s + j]; 10169566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Zdot)); 10179566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Zdot, i, w, Y)); 101861692a83SJed Brown 1019e27a552bSJed Brown /* Initial guess taken from last stage */ 10209566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Y[i])); 102161692a83SJed Brown 10227d4bf2deSEmil Constantinescu if (!ros->stage_explicit) { 10239566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 102461692a83SJed Brown if (!ros->recompute_jacobian && !i) { 10259566063dSJacob Faibussowitsch PetscCall(SNESGetLagJacobian(snes, &lag)); 10266aad120cSJose E. Roman if (lag == 1) { /* use did not set a nontrivial lag, so lag over all stages */ 10279566063dSJacob Faibussowitsch PetscCall(SNESSetLagJacobian(snes, -2)); /* Recompute the Jacobian on this solve, but not again for the rest of the stages */ 1028f7f07198SBarry Smith } 102961692a83SJed Brown } 10309566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, Y[i])); 10319371c9d4SSatish Balay if (!ros->recompute_jacobian && i == s - 1 && lag == 1) { PetscCall(SNESSetLagJacobian(snes, lag)); /* Set lag back to 1 so we know user did not set it */ } 10329566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 10339566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(snes, &lits)); 10349371c9d4SSatish Balay ts->snes_its += its; 10359371c9d4SSatish Balay ts->ksp_its += lits; 10367d4bf2deSEmil Constantinescu } else { 10371ce71dffSSatish Balay Mat J, Jp; 10389566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Ydot)); /* Evaluate Y[i]=G(t,Ydot=0,Zstage) */ 10399566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ros->stage_time, Zstage, Ydot, Y[i], PETSC_FALSE)); 10409566063dSJacob Faibussowitsch PetscCall(VecScale(Y[i], -1.0)); 10419566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y[i], -1.0, Zdot)); /*Y[i] = F(Zstage)-Zdot[=GammaInv*Y]*/ 10420feba352SEmil Constantinescu 10439566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Zstage)); /* Zstage = GammaExplicitCorr[i,j] * Y[j] */ 10440feba352SEmil Constantinescu for (j = 0; j < i; j++) w[j] = GammaExplicitCorr[i * s + j]; 10459566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Zstage, i, w, Y)); 1046fecfb714SLisandro Dalcin 1047fecfb714SLisandro Dalcin /* Y[i] = Y[i] + Jac*Zstage[=Jac*GammaExplicitCorr[i,j] * Y[j]] */ 10489566063dSJacob Faibussowitsch PetscCall(TSGetIJacobian(ts, &J, &Jp, NULL, NULL)); 10499566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, ros->stage_time, ts->vec_sol, Ydot, 0, J, Jp, PETSC_FALSE)); 10509566063dSJacob Faibussowitsch PetscCall(MatMult(J, Zstage, Zdot)); 10519566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y[i], -1.0, Zdot)); 10525ef26d82SJed Brown ts->ksp_its += 1; 1053fecfb714SLisandro Dalcin 10549566063dSJacob Faibussowitsch PetscCall(VecScale(Y[i], h)); 10557d4bf2deSEmil Constantinescu } 10569566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, ros->stage_time, i, Y)); 10579566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 10589566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(adapt, ts, ros->stage_time, Y[i], &stageok)); 1059fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 1060e27a552bSJed Brown } 1061e27a552bSJed Brown 1062b39943a6SLisandro Dalcin ros->status = TS_STEP_INCOMPLETE; 10639566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep_RosW(ts, tab->order, ts->vec_sol, NULL)); 1064b39943a6SLisandro Dalcin ros->status = TS_STEP_PENDING; 10659566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 10669566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(adapt)); 10679566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidateAdd(adapt, tab->name, tab->order, 1, tab->ccfl, (PetscReal)tab->s, PETSC_TRUE)); 10689566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 1069b39943a6SLisandro Dalcin ros->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 1070b39943a6SLisandro Dalcin if (!accept) { /* Roll back the current step */ 10719566063dSJacob Faibussowitsch PetscCall(TSRollBack_RosW(ts)); 1072be5899b3SLisandro Dalcin ts->time_step = next_time_step; 1073be5899b3SLisandro Dalcin goto reject_step; 1074b39943a6SLisandro Dalcin } 1075b39943a6SLisandro Dalcin 1076e27a552bSJed Brown ts->ptime += ts->time_step; 1077cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 10781c3436cfSJed Brown break; 1079b39943a6SLisandro Dalcin 1080b39943a6SLisandro Dalcin reject_step: 10819371c9d4SSatish Balay ts->reject++; 10829371c9d4SSatish Balay accept = PETSC_FALSE; 1083be5899b3SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 1084b39943a6SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 108563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections)); 10861c3436cfSJed Brown } 10871c3436cfSJed Brown } 10883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1089e27a552bSJed Brown } 1090e27a552bSJed Brown 1091d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_RosW(TS ts, PetscReal itime, Vec U) 1092d71ae5a4SJacob Faibussowitsch { 109361692a83SJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 1094f4aed992SEmil Constantinescu PetscInt s = ros->tableau->s, pinterp = ros->tableau->pinterp, i, j; 1095f4aed992SEmil Constantinescu PetscReal h; 1096f4aed992SEmil Constantinescu PetscReal tt, t; 1097f4aed992SEmil Constantinescu PetscScalar *bt; 1098f4aed992SEmil Constantinescu const PetscReal *Bt = ros->tableau->binterpt; 1099f4aed992SEmil Constantinescu const PetscReal *GammaInv = ros->tableau->GammaInv; 1100f4aed992SEmil Constantinescu PetscScalar *w = ros->work; 1101f4aed992SEmil Constantinescu Vec *Y = ros->Y; 1102e27a552bSJed Brown 1103e27a552bSJed Brown PetscFunctionBegin; 11043c633725SBarry Smith PetscCheck(Bt, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSRosW %s does not have an interpolation formula", ros->tableau->name); 1105f4aed992SEmil Constantinescu 1106f4aed992SEmil Constantinescu switch (ros->status) { 1107f4aed992SEmil Constantinescu case TS_STEP_INCOMPLETE: 1108f4aed992SEmil Constantinescu case TS_STEP_PENDING: 1109f4aed992SEmil Constantinescu h = ts->time_step; 1110f4aed992SEmil Constantinescu t = (itime - ts->ptime) / h; 1111f4aed992SEmil Constantinescu break; 1112f4aed992SEmil Constantinescu case TS_STEP_COMPLETE: 1113be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; 1114f4aed992SEmil Constantinescu t = (itime - ts->ptime) / h + 1; /* In the interval [0,1] */ 1115f4aed992SEmil Constantinescu break; 1116d71ae5a4SJacob Faibussowitsch default: 1117d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus"); 1118f4aed992SEmil Constantinescu } 11199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s, &bt)); 1120f4aed992SEmil Constantinescu for (i = 0; i < s; i++) bt[i] = 0; 1121f4aed992SEmil Constantinescu for (j = 0, tt = t; j < pinterp; j++, tt *= t) { 1122ad540459SPierre Jolivet for (i = 0; i < s; i++) bt[i] += Bt[i * pinterp + j] * tt; 1123f4aed992SEmil Constantinescu } 1124f4aed992SEmil Constantinescu 1125f4aed992SEmil Constantinescu /* y(t+tt*h) = y(t) + Sum bt(tt) * GammaInv * Ydot */ 1126f9c1d6abSBarry Smith /* U <- 0*/ 11279566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(U)); 1128f9c1d6abSBarry Smith /* U <- Sum bt_i * GammaInv(i,1:i) * Y(1:i) */ 11293ca35412SEmil Constantinescu for (j = 0; j < s; j++) w[j] = 0; 11303ca35412SEmil Constantinescu for (j = 0; j < s; j++) { 1131ad540459SPierre Jolivet for (i = j; i < s; i++) w[j] += bt[i] * GammaInv[i * s + j]; 11323ca35412SEmil Constantinescu } 11339566063dSJacob Faibussowitsch PetscCall(VecMAXPY(U, i, w, Y)); 1134be5899b3SLisandro Dalcin /* U <- y(t) + U */ 11359566063dSJacob Faibussowitsch PetscCall(VecAXPY(U, 1, ros->vec_sol_prev)); 1136f4aed992SEmil Constantinescu 11379566063dSJacob Faibussowitsch PetscCall(PetscFree(bt)); 11383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1139e27a552bSJed Brown } 1140e27a552bSJed Brown 1141e27a552bSJed Brown /*------------------------------------------------------------*/ 1142b39943a6SLisandro Dalcin 1143d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWTableauReset(TS ts) 1144d71ae5a4SJacob Faibussowitsch { 1145b39943a6SLisandro Dalcin TS_RosW *ros = (TS_RosW *)ts->data; 1146b39943a6SLisandro Dalcin RosWTableau tab = ros->tableau; 1147b39943a6SLisandro Dalcin 1148b39943a6SLisandro Dalcin PetscFunctionBegin; 11493ba16761SJacob Faibussowitsch if (!tab) PetscFunctionReturn(PETSC_SUCCESS); 11509566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &ros->Y)); 11519566063dSJacob Faibussowitsch PetscCall(PetscFree(ros->work)); 11523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1153b39943a6SLisandro Dalcin } 1154b39943a6SLisandro Dalcin 1155d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_RosW(TS ts) 1156d71ae5a4SJacob Faibussowitsch { 115761692a83SJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 1158e27a552bSJed Brown 1159e27a552bSJed Brown PetscFunctionBegin; 11609566063dSJacob Faibussowitsch PetscCall(TSRosWTableauReset(ts)); 11619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ros->Ydot)); 11629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ros->Ystage)); 11639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ros->Zdot)); 11649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ros->Zstage)); 11659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ros->vec_sol_prev)); 11663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1167e27a552bSJed Brown } 1168e27a552bSJed Brown 1169d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWGetVecs(TS ts, DM dm, Vec *Ydot, Vec *Zdot, Vec *Ystage, Vec *Zstage) 1170d71ae5a4SJacob Faibussowitsch { 1171d5e6173cSPeter Brune TS_RosW *rw = (TS_RosW *)ts->data; 1172d5e6173cSPeter Brune 1173d5e6173cSPeter Brune PetscFunctionBegin; 1174d5e6173cSPeter Brune if (Ydot) { 1175d5e6173cSPeter Brune if (dm && dm != ts->dm) { 11769566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Ydot", Ydot)); 1177d5e6173cSPeter Brune } else *Ydot = rw->Ydot; 1178d5e6173cSPeter Brune } 1179d5e6173cSPeter Brune if (Zdot) { 1180d5e6173cSPeter Brune if (dm && dm != ts->dm) { 11819566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Zdot", Zdot)); 1182d5e6173cSPeter Brune } else *Zdot = rw->Zdot; 1183d5e6173cSPeter Brune } 1184d5e6173cSPeter Brune if (Ystage) { 1185d5e6173cSPeter Brune if (dm && dm != ts->dm) { 11869566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Ystage", Ystage)); 1187d5e6173cSPeter Brune } else *Ystage = rw->Ystage; 1188d5e6173cSPeter Brune } 1189d5e6173cSPeter Brune if (Zstage) { 1190d5e6173cSPeter Brune if (dm && dm != ts->dm) { 11919566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Zstage", Zstage)); 1192d5e6173cSPeter Brune } else *Zstage = rw->Zstage; 1193d5e6173cSPeter Brune } 11943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1195d5e6173cSPeter Brune } 1196d5e6173cSPeter Brune 1197d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWRestoreVecs(TS ts, DM dm, Vec *Ydot, Vec *Zdot, Vec *Ystage, Vec *Zstage) 1198d71ae5a4SJacob Faibussowitsch { 1199d5e6173cSPeter Brune PetscFunctionBegin; 1200d5e6173cSPeter Brune if (Ydot) { 120148a46eb9SPierre Jolivet if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Ydot", Ydot)); 1202d5e6173cSPeter Brune } 1203d5e6173cSPeter Brune if (Zdot) { 120448a46eb9SPierre Jolivet if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Zdot", Zdot)); 1205d5e6173cSPeter Brune } 1206d5e6173cSPeter Brune if (Ystage) { 120748a46eb9SPierre Jolivet if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Ystage", Ystage)); 1208d5e6173cSPeter Brune } 1209d5e6173cSPeter Brune if (Zstage) { 121048a46eb9SPierre Jolivet if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Zstage", Zstage)); 1211d5e6173cSPeter Brune } 12123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1213d5e6173cSPeter Brune } 1214d5e6173cSPeter Brune 1215d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCoarsenHook_TSRosW(DM fine, DM coarse, void *ctx) 1216d71ae5a4SJacob Faibussowitsch { 1217d5e6173cSPeter Brune PetscFunctionBegin; 12183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1219d5e6173cSPeter Brune } 1220d5e6173cSPeter Brune 1221d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRestrictHook_TSRosW(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) 1222d71ae5a4SJacob Faibussowitsch { 1223d5e6173cSPeter Brune TS ts = (TS)ctx; 1224d5e6173cSPeter Brune Vec Ydot, Zdot, Ystage, Zstage; 1225d5e6173cSPeter Brune Vec Ydotc, Zdotc, Ystagec, Zstagec; 1226d5e6173cSPeter Brune 1227d5e6173cSPeter Brune PetscFunctionBegin; 12289566063dSJacob Faibussowitsch PetscCall(TSRosWGetVecs(ts, fine, &Ydot, &Ystage, &Zdot, &Zstage)); 12299566063dSJacob Faibussowitsch PetscCall(TSRosWGetVecs(ts, coarse, &Ydotc, &Ystagec, &Zdotc, &Zstagec)); 12309566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, Ydot, Ydotc)); 12319566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Ydotc, rscale, Ydotc)); 12329566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, Ystage, Ystagec)); 12339566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Ystagec, rscale, Ystagec)); 12349566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, Zdot, Zdotc)); 12359566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Zdotc, rscale, Zdotc)); 12369566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, Zstage, Zstagec)); 12379566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Zstagec, rscale, Zstagec)); 12389566063dSJacob Faibussowitsch PetscCall(TSRosWRestoreVecs(ts, fine, &Ydot, &Ystage, &Zdot, &Zstage)); 12399566063dSJacob Faibussowitsch PetscCall(TSRosWRestoreVecs(ts, coarse, &Ydotc, &Ystagec, &Zdotc, &Zstagec)); 12403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1241d5e6173cSPeter Brune } 1242d5e6173cSPeter Brune 1243d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainHook_TSRosW(DM fine, DM coarse, void *ctx) 1244d71ae5a4SJacob Faibussowitsch { 1245258e1594SPeter Brune PetscFunctionBegin; 12463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1247258e1594SPeter Brune } 1248258e1594SPeter Brune 1249d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainRestrictHook_TSRosW(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) 1250d71ae5a4SJacob Faibussowitsch { 1251258e1594SPeter Brune TS ts = (TS)ctx; 1252258e1594SPeter Brune Vec Ydot, Zdot, Ystage, Zstage; 1253258e1594SPeter Brune Vec Ydots, Zdots, Ystages, Zstages; 1254258e1594SPeter Brune 1255258e1594SPeter Brune PetscFunctionBegin; 12569566063dSJacob Faibussowitsch PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Ystage, &Zdot, &Zstage)); 12579566063dSJacob Faibussowitsch PetscCall(TSRosWGetVecs(ts, subdm, &Ydots, &Ystages, &Zdots, &Zstages)); 1258258e1594SPeter Brune 12599566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, Ydot, Ydots, INSERT_VALUES, SCATTER_FORWARD)); 12609566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, Ydot, Ydots, INSERT_VALUES, SCATTER_FORWARD)); 1261258e1594SPeter Brune 12629566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, Ystage, Ystages, INSERT_VALUES, SCATTER_FORWARD)); 12639566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, Ystage, Ystages, INSERT_VALUES, SCATTER_FORWARD)); 1264258e1594SPeter Brune 12659566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, Zdot, Zdots, INSERT_VALUES, SCATTER_FORWARD)); 12669566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, Zdot, Zdots, INSERT_VALUES, SCATTER_FORWARD)); 1267258e1594SPeter Brune 12689566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, Zstage, Zstages, INSERT_VALUES, SCATTER_FORWARD)); 12699566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, Zstage, Zstages, INSERT_VALUES, SCATTER_FORWARD)); 1270258e1594SPeter Brune 12719566063dSJacob Faibussowitsch PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Ystage, &Zdot, &Zstage)); 12729566063dSJacob Faibussowitsch PetscCall(TSRosWRestoreVecs(ts, subdm, &Ydots, &Ystages, &Zdots, &Zstages)); 12733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1274258e1594SPeter Brune } 1275258e1594SPeter Brune 1276e27a552bSJed Brown /* 1277e27a552bSJed Brown This defines the nonlinear equation that is to be solved with SNES 1278e27a552bSJed Brown G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 1279e27a552bSJed Brown */ 1280d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_RosW(SNES snes, Vec U, Vec F, TS ts) 1281d71ae5a4SJacob Faibussowitsch { 128261692a83SJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 1283d5e6173cSPeter Brune Vec Ydot, Zdot, Ystage, Zstage; 1284b296d7d5SJed Brown PetscReal shift = ros->scoeff / ts->time_step; 1285d5e6173cSPeter Brune DM dm, dmsave; 1286e27a552bSJed Brown 1287e27a552bSJed Brown PetscFunctionBegin; 12889566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 12899566063dSJacob Faibussowitsch PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage)); 12909566063dSJacob Faibussowitsch PetscCall(VecWAXPY(Ydot, shift, U, Zdot)); /* Ydot = shift*U + Zdot */ 12919566063dSJacob Faibussowitsch PetscCall(VecWAXPY(Ystage, 1.0, U, Zstage)); /* Ystage = U + Zstage */ 1292d5e6173cSPeter Brune dmsave = ts->dm; 1293d5e6173cSPeter Brune ts->dm = dm; 12949566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ros->stage_time, Ystage, Ydot, F, PETSC_FALSE)); 1295d5e6173cSPeter Brune ts->dm = dmsave; 12969566063dSJacob Faibussowitsch PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage)); 12973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1298e27a552bSJed Brown } 1299e27a552bSJed Brown 1300d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes, Vec U, Mat A, Mat B, TS ts) 1301d71ae5a4SJacob Faibussowitsch { 130261692a83SJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 1303d5e6173cSPeter Brune Vec Ydot, Zdot, Ystage, Zstage; 1304b296d7d5SJed Brown PetscReal shift = ros->scoeff / ts->time_step; 1305d5e6173cSPeter Brune DM dm, dmsave; 1306e27a552bSJed Brown 1307e27a552bSJed Brown PetscFunctionBegin; 130861692a83SJed Brown /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */ 13099566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 13109566063dSJacob Faibussowitsch PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage)); 1311d5e6173cSPeter Brune dmsave = ts->dm; 1312d5e6173cSPeter Brune ts->dm = dm; 13139566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, ros->stage_time, Ystage, Ydot, shift, A, B, PETSC_TRUE)); 1314d5e6173cSPeter Brune ts->dm = dmsave; 13159566063dSJacob Faibussowitsch PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage)); 13163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1317e27a552bSJed Brown } 1318e27a552bSJed Brown 1319d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWTableauSetUp(TS ts) 1320d71ae5a4SJacob Faibussowitsch { 1321b39943a6SLisandro Dalcin TS_RosW *ros = (TS_RosW *)ts->data; 1322b39943a6SLisandro Dalcin RosWTableau tab = ros->tableau; 1323b39943a6SLisandro Dalcin 1324b39943a6SLisandro Dalcin PetscFunctionBegin; 13259566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ros->Y)); 13269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &ros->work)); 13273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1328b39943a6SLisandro Dalcin } 1329b39943a6SLisandro Dalcin 1330d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_RosW(TS ts) 1331d71ae5a4SJacob Faibussowitsch { 133261692a83SJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 1333d5e6173cSPeter Brune DM dm; 1334b39943a6SLisandro Dalcin SNES snes; 1335a3ab5968SHong Zhang TSRHSJacobian rhsjacobian; 1336e27a552bSJed Brown 1337e27a552bSJed Brown PetscFunctionBegin; 13389566063dSJacob Faibussowitsch PetscCall(TSRosWTableauSetUp(ts)); 13399566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &ros->Ydot)); 13409566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &ros->Ystage)); 13419566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &ros->Zdot)); 13429566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &ros->Zstage)); 13439566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &ros->vec_sol_prev)); 13449566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 13459566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSRosW, DMRestrictHook_TSRosW, ts)); 13469566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSRosW, DMSubDomainRestrictHook_TSRosW, ts)); 1347b39943a6SLisandro Dalcin /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */ 13489566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 134948a46eb9SPierre Jolivet if (!((PetscObject)snes)->type_name) PetscCall(SNESSetType(snes, SNESKSPONLY)); 13509566063dSJacob Faibussowitsch PetscCall(DMTSGetRHSJacobian(dm, &rhsjacobian, NULL)); 1351a3ab5968SHong Zhang if (rhsjacobian == TSComputeRHSJacobianConstant) { 1352a3ab5968SHong Zhang Mat Amat, Pmat; 1353a3ab5968SHong Zhang 1354a3ab5968SHong Zhang /* Set the SNES matrix to be different from the RHS matrix because there is no way to reconstruct shift*M-J */ 13559566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(snes, &Amat, &Pmat, NULL, NULL)); 1356a3ab5968SHong Zhang if (Amat && Amat == ts->Arhs) { 1357a3ab5968SHong Zhang if (Amat == Pmat) { 13589566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ts->Arhs, MAT_COPY_VALUES, &Amat)); 13599566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, Amat, Amat, NULL, NULL)); 1360a3ab5968SHong Zhang } else { 13619566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ts->Arhs, MAT_COPY_VALUES, &Amat)); 13629566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, Amat, NULL, NULL, NULL)); 1363a3ab5968SHong Zhang if (Pmat && Pmat == ts->Brhs) { 13649566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ts->Brhs, MAT_COPY_VALUES, &Pmat)); 13659566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, NULL, Pmat, NULL, NULL)); 13669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Pmat)); 1367a3ab5968SHong Zhang } 1368a3ab5968SHong Zhang } 13699566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Amat)); 1370a3ab5968SHong Zhang } 1371a3ab5968SHong Zhang } 13723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1373e27a552bSJed Brown } 1374e27a552bSJed Brown /*------------------------------------------------------------*/ 1375e27a552bSJed Brown 1376d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_RosW(TS ts, PetscOptionItems *PetscOptionsObject) 1377d71ae5a4SJacob Faibussowitsch { 137861692a83SJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 1379b39943a6SLisandro Dalcin SNES snes; 1380e27a552bSJed Brown 1381e27a552bSJed Brown PetscFunctionBegin; 1382d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "RosW ODE solver options"); 1383e27a552bSJed Brown { 138461692a83SJed Brown RosWTableauLink link; 1385e27a552bSJed Brown PetscInt count, choice; 1386e27a552bSJed Brown PetscBool flg; 1387e27a552bSJed Brown const char **namelist; 138861692a83SJed Brown 13899371c9d4SSatish Balay for (link = RosWTableauList, count = 0; link; link = link->next, count++) 13909371c9d4SSatish Balay ; 13919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count, (char ***)&namelist)); 139261692a83SJed Brown for (link = RosWTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name; 13939566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_rosw_type", "Family of Rosenbrock-W method", "TSRosWSetType", (const char *const *)namelist, count, ros->tableau->name, &choice, &flg)); 13949566063dSJacob Faibussowitsch if (flg) PetscCall(TSRosWSetType(ts, namelist[choice])); 13959566063dSJacob Faibussowitsch PetscCall(PetscFree(namelist)); 139661692a83SJed Brown 13979566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_rosw_recompute_jacobian", "Recompute the Jacobian at each stage", "TSRosWSetRecomputeJacobian", ros->recompute_jacobian, &ros->recompute_jacobian, NULL)); 1398b39943a6SLisandro Dalcin } 1399d0609cedSBarry Smith PetscOptionsHeadEnd(); 140061692a83SJed Brown /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */ 14019566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 140248a46eb9SPierre Jolivet if (!((PetscObject)snes)->type_name) PetscCall(SNESSetType(snes, SNESKSPONLY)); 14033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1404e27a552bSJed Brown } 1405e27a552bSJed Brown 1406d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_RosW(TS ts, PetscViewer viewer) 1407d71ae5a4SJacob Faibussowitsch { 140861692a83SJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 1409e27a552bSJed Brown PetscBool iascii; 1410e27a552bSJed Brown 1411e27a552bSJed Brown PetscFunctionBegin; 14129566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 1413e27a552bSJed Brown if (iascii) { 14149c334d8fSLisandro Dalcin RosWTableau tab = ros->tableau; 141519fd82e9SBarry Smith TSRosWType rostype; 14169c334d8fSLisandro Dalcin char buf[512]; 1417e408995aSJed Brown PetscInt i; 1418e408995aSJed Brown PetscReal abscissa[512]; 14199566063dSJacob Faibussowitsch PetscCall(TSRosWGetType(ts, &rostype)); 14209566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Rosenbrock-W %s\n", rostype)); 14219566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, tab->ASum)); 14229566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Abscissa of A = %s\n", buf)); 1423e408995aSJed Brown for (i = 0; i < tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i]; 14249566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, abscissa)); 14259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Abscissa of A+Gamma = %s\n", buf)); 1426e27a552bSJed Brown } 14273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1428e27a552bSJed Brown } 1429e27a552bSJed Brown 1430d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSLoad_RosW(TS ts, PetscViewer viewer) 1431d71ae5a4SJacob Faibussowitsch { 14329200755eSBarry Smith SNES snes; 14339c334d8fSLisandro Dalcin TSAdapt adapt; 14349200755eSBarry Smith 14359200755eSBarry Smith PetscFunctionBegin; 14369566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 14379566063dSJacob Faibussowitsch PetscCall(TSAdaptLoad(adapt, viewer)); 14389566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 14399566063dSJacob Faibussowitsch PetscCall(SNESLoad(snes, viewer)); 14409200755eSBarry Smith /* function and Jacobian context for SNES when used with TS is always ts object */ 14419566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, NULL, NULL, ts)); 14429566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, ts)); 14433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14449200755eSBarry Smith } 14459200755eSBarry Smith 1446e27a552bSJed Brown /*@C 1447bcf0153eSBarry Smith TSRosWSetType - Set the type of Rosenbrock-W, `TSROSW`, scheme 1448e27a552bSJed Brown 144920f4b53cSBarry Smith Logically Collective 1450e27a552bSJed Brown 1451d8d19677SJose E. Roman Input Parameters: 1452e27a552bSJed Brown + ts - timestepping context 1453b92453a8SLisandro Dalcin - roswtype - type of Rosenbrock-W scheme 1454e27a552bSJed Brown 1455020d8f30SJed Brown Level: beginner 1456e27a552bSJed Brown 1457*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSRosWGetType()`, `TSROSW`, `TSROSW2M`, `TSROSW2P`, `TSROSWRA3PW`, `TSROSWRA34PW2`, `TSROSWRODAS3`, `TSROSWSANDU3`, `TSROSWASSP3P3S1C`, `TSROSWLASSP3P4S2C`, `TSROSWLLSSP3P4S2C`, `TSROSWARK3` 1458e27a552bSJed Brown @*/ 1459d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWSetType(TS ts, TSRosWType roswtype) 1460d71ae5a4SJacob Faibussowitsch { 1461e27a552bSJed Brown PetscFunctionBegin; 1462e27a552bSJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1463b92453a8SLisandro Dalcin PetscValidCharPointer(roswtype, 2); 1464cac4c232SBarry Smith PetscTryMethod(ts, "TSRosWSetType_C", (TS, TSRosWType), (ts, roswtype)); 14653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1466e27a552bSJed Brown } 1467e27a552bSJed Brown 1468e27a552bSJed Brown /*@C 146961692a83SJed Brown TSRosWGetType - Get the type of Rosenbrock-W scheme 1470e27a552bSJed Brown 147120f4b53cSBarry Smith Logically Collective 1472e27a552bSJed Brown 1473e27a552bSJed Brown Input Parameter: 1474e27a552bSJed Brown . ts - timestepping context 1475e27a552bSJed Brown 1476e27a552bSJed Brown Output Parameter: 147761692a83SJed Brown . rostype - type of Rosenbrock-W scheme 1478e27a552bSJed Brown 1479e27a552bSJed Brown Level: intermediate 1480e27a552bSJed Brown 1481*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSRosWType`, `TSRosWSetType()` 1482e27a552bSJed Brown @*/ 1483d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWGetType(TS ts, TSRosWType *rostype) 1484d71ae5a4SJacob Faibussowitsch { 1485e27a552bSJed Brown PetscFunctionBegin; 1486e27a552bSJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1487cac4c232SBarry Smith PetscUseMethod(ts, "TSRosWGetType_C", (TS, TSRosWType *), (ts, rostype)); 14883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1489e27a552bSJed Brown } 1490e27a552bSJed Brown 1491e27a552bSJed Brown /*@C 149261692a83SJed Brown TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step. 1493e27a552bSJed Brown 149420f4b53cSBarry Smith Logically Collective 1495e27a552bSJed Brown 1496d8d19677SJose E. Roman Input Parameters: 1497e27a552bSJed Brown + ts - timestepping context 1498bcf0153eSBarry Smith - flg - `PETSC_TRUE` to recompute the Jacobian at each stage 1499e27a552bSJed Brown 1500e27a552bSJed Brown Level: intermediate 1501e27a552bSJed Brown 1502*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSRosWType`, `TSRosWGetType()` 1503e27a552bSJed Brown @*/ 1504d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWSetRecomputeJacobian(TS ts, PetscBool flg) 1505d71ae5a4SJacob Faibussowitsch { 1506e27a552bSJed Brown PetscFunctionBegin; 1507e27a552bSJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1508cac4c232SBarry Smith PetscTryMethod(ts, "TSRosWSetRecomputeJacobian_C", (TS, PetscBool), (ts, flg)); 15093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1510e27a552bSJed Brown } 1511e27a552bSJed Brown 1512d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWGetType_RosW(TS ts, TSRosWType *rostype) 1513d71ae5a4SJacob Faibussowitsch { 151461692a83SJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 1515e27a552bSJed Brown 1516e27a552bSJed Brown PetscFunctionBegin; 151761692a83SJed Brown *rostype = ros->tableau->name; 15183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1519e27a552bSJed Brown } 1520ef20d060SBarry Smith 1521d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWSetType_RosW(TS ts, TSRosWType rostype) 1522d71ae5a4SJacob Faibussowitsch { 152361692a83SJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 1524e27a552bSJed Brown PetscBool match; 152561692a83SJed Brown RosWTableauLink link; 1526e27a552bSJed Brown 1527e27a552bSJed Brown PetscFunctionBegin; 152861692a83SJed Brown if (ros->tableau) { 15299566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(ros->tableau->name, rostype, &match)); 15303ba16761SJacob Faibussowitsch if (match) PetscFunctionReturn(PETSC_SUCCESS); 1531e27a552bSJed Brown } 153261692a83SJed Brown for (link = RosWTableauList; link; link = link->next) { 15339566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(link->tab.name, rostype, &match)); 1534e27a552bSJed Brown if (match) { 15359566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSRosWTableauReset(ts)); 153661692a83SJed Brown ros->tableau = &link->tab; 15379566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSRosWTableauSetUp(ts)); 15382ffb9264SLisandro Dalcin ts->default_adapt_type = ros->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 15393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1540e27a552bSJed Brown } 1541e27a552bSJed Brown } 154298921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", rostype); 1543e27a552bSJed Brown } 154461692a83SJed Brown 1545d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWSetRecomputeJacobian_RosW(TS ts, PetscBool flg) 1546d71ae5a4SJacob Faibussowitsch { 154761692a83SJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 1548e27a552bSJed Brown 1549e27a552bSJed Brown PetscFunctionBegin; 155061692a83SJed Brown ros->recompute_jacobian = flg; 15513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1552e27a552bSJed Brown } 1553e27a552bSJed Brown 1554d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_RosW(TS ts) 1555d71ae5a4SJacob Faibussowitsch { 1556b3a6b972SJed Brown PetscFunctionBegin; 15579566063dSJacob Faibussowitsch PetscCall(TSReset_RosW(ts)); 1558b3a6b972SJed Brown if (ts->dm) { 15599566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSRosW, DMRestrictHook_TSRosW, ts)); 15609566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSRosW, DMSubDomainRestrictHook_TSRosW, ts)); 1561b3a6b972SJed Brown } 15629566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 15639566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWGetType_C", NULL)); 15649566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetType_C", NULL)); 15659566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetRecomputeJacobian_C", NULL)); 15663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1567b3a6b972SJed Brown } 1568d5e6173cSPeter Brune 1569e27a552bSJed Brown /* ------------------------------------------------------------ */ 1570e27a552bSJed Brown /*MC 1571020d8f30SJed Brown TSROSW - ODE solver using Rosenbrock-W schemes 1572e27a552bSJed Brown 1573e27a552bSJed Brown These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly 1574e27a552bSJed Brown nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part 1575bcf0153eSBarry Smith of the equation using `TSSetIFunction()` and the non-stiff part with `TSSetRHSFunction()`. 1576bcf0153eSBarry Smith 1577bcf0153eSBarry Smith Level: beginner 1578e27a552bSJed Brown 1579e27a552bSJed Brown Notes: 158061692a83SJed Brown This method currently only works with autonomous ODE and DAE. 158161692a83SJed Brown 1582bcf0153eSBarry Smith Consider trying `TSARKIMEX` if the stiff part is strongly nonlinear. 1583d0685a90SJed Brown 15843d5a8a6aSBarry Smith Since this uses a single linear solve per time-step if you wish to lag the jacobian or preconditioner computation you must use also -snes_lag_jacobian_persists true or -snes_lag_jacobian_preconditioner true 15853d5a8a6aSBarry Smith 1586e94cfbe0SPatrick Sanan Developer Notes: 158761692a83SJed Brown Rosenbrock-W methods are typically specified for autonomous ODE 158861692a83SJed Brown 1589f9c1d6abSBarry Smith $ udot = f(u) 159061692a83SJed Brown 159161692a83SJed Brown by the stage equations 159261692a83SJed Brown 1593f9c1d6abSBarry Smith $ k_i = h f(u_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j 159461692a83SJed Brown 159561692a83SJed Brown and step completion formula 159661692a83SJed Brown 1597f9c1d6abSBarry Smith $ u_1 = u_0 + sum_j b_j k_j 159861692a83SJed Brown 1599f9c1d6abSBarry Smith with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(u) 160061692a83SJed Brown and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner, 160161692a83SJed Brown we define new variables for the stage equations 160261692a83SJed Brown 160361692a83SJed Brown $ y_i = gamma_ij k_j 160461692a83SJed Brown 160561692a83SJed Brown The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define 160661692a83SJed Brown 1607b70472e2SJed Brown $ A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-1} 160861692a83SJed Brown 160961692a83SJed Brown to rewrite the method as 161061692a83SJed Brown 161137fdd005SBarry Smith .vb 161237fdd005SBarry Smith [M/(h gamma_ii) - J] y_i = f(u_0 + sum_j a_ij y_j) + M sum_j (c_ij/h) y_j 161337fdd005SBarry Smith u_1 = u_0 + sum_j bt_j y_j 161437fdd005SBarry Smith .ve 161561692a83SJed Brown 161661692a83SJed Brown where we have introduced the mass matrix M. Continue by defining 161761692a83SJed Brown 161861692a83SJed Brown $ ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j 161961692a83SJed Brown 162061692a83SJed Brown or, more compactly in tensor notation 162161692a83SJed Brown 162261692a83SJed Brown $ Ydot = 1/h (Gamma^{-1} \otimes I) Y . 162361692a83SJed Brown 162461692a83SJed Brown Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current 162561692a83SJed Brown stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the 162661692a83SJed Brown equation 162761692a83SJed Brown 1628f9c1d6abSBarry Smith $ g(u_0 + sum_j a_ij y_j + y_i, ydot_i) = 0 162961692a83SJed Brown 163061692a83SJed Brown with initial guess y_i = 0. 1631e27a552bSJed Brown 1632*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSRosWSetType()`, `TSRosWRegister()`, `TSROSWTHETA1`, `TSROSWTHETA2`, `TSROSW2M`, `TSROSW2P`, `TSROSWRA3PW`, `TSROSWRA34PW2`, `TSROSWRODAS3`, 1633bcf0153eSBarry Smith `TSROSWSANDU3`, `TSROSWASSP3P3S1C`, `TSROSWLASSP3P4S2C`, `TSROSWLLSSP3P4S2C`, `TSROSWGRK4T`, `TSROSWSHAMP4`, `TSROSWVELDD4`, `TSROSW4L`, `TSType` 1634e27a552bSJed Brown M*/ 1635d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_RosW(TS ts) 1636d71ae5a4SJacob Faibussowitsch { 163761692a83SJed Brown TS_RosW *ros; 1638e27a552bSJed Brown 1639e27a552bSJed Brown PetscFunctionBegin; 16409566063dSJacob Faibussowitsch PetscCall(TSRosWInitializePackage()); 1641e27a552bSJed Brown 1642e27a552bSJed Brown ts->ops->reset = TSReset_RosW; 1643e27a552bSJed Brown ts->ops->destroy = TSDestroy_RosW; 1644e27a552bSJed Brown ts->ops->view = TSView_RosW; 16459200755eSBarry Smith ts->ops->load = TSLoad_RosW; 1646e27a552bSJed Brown ts->ops->setup = TSSetUp_RosW; 1647e27a552bSJed Brown ts->ops->step = TSStep_RosW; 1648e27a552bSJed Brown ts->ops->interpolate = TSInterpolate_RosW; 16491c3436cfSJed Brown ts->ops->evaluatestep = TSEvaluateStep_RosW; 165024655328SShri ts->ops->rollback = TSRollBack_RosW; 1651e27a552bSJed Brown ts->ops->setfromoptions = TSSetFromOptions_RosW; 1652e27a552bSJed Brown ts->ops->snesfunction = SNESTSFormFunction_RosW; 1653e27a552bSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_RosW; 1654e27a552bSJed Brown 1655efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 1656efd4aadfSBarry Smith 16574dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&ros)); 165861692a83SJed Brown ts->data = (void *)ros; 1659e27a552bSJed Brown 16609566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWGetType_C", TSRosWGetType_RosW)); 16619566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetType_C", TSRosWSetType_RosW)); 16629566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetRecomputeJacobian_C", TSRosWSetRecomputeJacobian_RosW)); 1663b39943a6SLisandro Dalcin 16649566063dSJacob Faibussowitsch PetscCall(TSRosWSetType(ts, TSRosWDefault)); 16653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1666e27a552bSJed Brown } 1667