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 73db781477SPatrick Sanan .seealso: `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 83db781477SPatrick Sanan .seealso: `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 93db781477SPatrick Sanan .seealso: `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 103db781477SPatrick Sanan .seealso: `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 113fe7e6d57SJed Brown References: 114606c0280SSatish Balay . * - Rang and Angermann, New Rosenbrock W methods of order 3 for partial differential algebraic equations of index 1, 2005. 115fe7e6d57SJed Brown 116fe7e6d57SJed Brown Level: intermediate 117fe7e6d57SJed Brown 118db781477SPatrick Sanan .seealso: `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 128fe7e6d57SJed Brown References: 129606c0280SSatish Balay . * - Rang and Angermann, New Rosenbrock W methods of order 3 for partial differential algebraic equations of index 1, 2005. 130fe7e6d57SJed Brown 131fe7e6d57SJed Brown Level: intermediate 132fe7e6d57SJed Brown 133db781477SPatrick Sanan .seealso: `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 143ef3c5b88SJed Brown References: 144606c0280SSatish Balay . * - Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997. 145ef3c5b88SJed Brown 146ef3c5b88SJed Brown Level: intermediate 147ef3c5b88SJed Brown 148db781477SPatrick Sanan .seealso: `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 161ef3c5b88SJed Brown References: 162606c0280SSatish Balay . * - Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997. 163ef3c5b88SJed Brown 164ef3c5b88SJed Brown Level: intermediate 165ef3c5b88SJed Brown 166db781477SPatrick Sanan .seealso: `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 176961f28d0SJed Brown References: 177606c0280SSatish Balay . * - Emil Constantinescu 178961f28d0SJed Brown 179961f28d0SJed Brown Level: intermediate 180961f28d0SJed Brown 181db781477SPatrick Sanan .seealso: `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 191961f28d0SJed Brown References: 192606c0280SSatish Balay . * - Emil Constantinescu 193961f28d0SJed Brown 194961f28d0SJed Brown Level: intermediate 195961f28d0SJed Brown 196db781477SPatrick Sanan .seealso: `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 206961f28d0SJed Brown References: 207606c0280SSatish Balay . * - Emil Constantinescu 208961f28d0SJed Brown 209961f28d0SJed Brown Level: intermediate 210961f28d0SJed Brown 211db781477SPatrick Sanan .seealso: `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 22342faf41dSJed Brown References: 224606c0280SSatish Balay + * - Kaps and Rentrop, Generalized Runge Kutta methods of order four with stepsize control for stiff ordinary differential equations, 1979. 225606c0280SSatish Balay - * - Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2. 22642faf41dSJed Brown 22742faf41dSJed Brown Hairer's code ros4.f 22842faf41dSJed Brown 22942faf41dSJed Brown Level: intermediate 23042faf41dSJed Brown 231db781477SPatrick Sanan .seealso: `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 24342faf41dSJed Brown References: 244606c0280SSatish Balay + * - Shampine, Implementation of Rosenbrock methods, 1982. 245606c0280SSatish Balay - * - Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2. 24642faf41dSJed Brown 24742faf41dSJed Brown Hairer's code ros4.f 24842faf41dSJed Brown 24942faf41dSJed Brown Level: intermediate 25042faf41dSJed Brown 251db781477SPatrick Sanan .seealso: `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 26342faf41dSJed Brown References: 264606c0280SSatish Balay + * - van Veldhuizen, D stability and Kaps Rentrop methods, 1984. 265606c0280SSatish Balay - * - Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2. 26642faf41dSJed Brown 26742faf41dSJed Brown Hairer's code ros4.f 26842faf41dSJed Brown 26942faf41dSJed Brown Level: intermediate 27042faf41dSJed Brown 271db781477SPatrick Sanan .seealso: `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 28342faf41dSJed Brown References: 284606c0280SSatish Balay . * - Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2. 28542faf41dSJed Brown 28642faf41dSJed Brown Hairer's code ros4.f 28742faf41dSJed Brown 28842faf41dSJed Brown Level: intermediate 28942faf41dSJed Brown 290db781477SPatrick Sanan .seealso: `TSROSW`, `TSROSWGRK4T`, `TSROSWSHAMP4`, `TSROSW4L` 29142faf41dSJed Brown M*/ 29242faf41dSJed Brown 293e27a552bSJed Brown /*@C 294be5899b3SLisandro Dalcin 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 300db781477SPatrick Sanan .seealso: `TSRosWRegisterDestroy()` 301e27a552bSJed Brown @*/ 3029371c9d4SSatish Balay PetscErrorCode TSRosWRegisterAll(void) { 303e27a552bSJed Brown PetscFunctionBegin; 304e27a552bSJed Brown if (TSRosWRegisterAllCalled) PetscFunctionReturn(0); 305e27a552bSJed Brown TSRosWRegisterAllCalled = PETSC_TRUE; 306e27a552bSJed Brown 307e27a552bSJed Brown { 308bbd56ea5SKarl Rupp const PetscReal A = 0; 309bbd56ea5SKarl Rupp const PetscReal Gamma = 1; 310bbd56ea5SKarl Rupp const PetscReal b = 1; 311bbd56ea5SKarl Rupp const PetscReal binterpt = 1; 3121f80e275SEmil Constantinescu 3139566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSWTHETA1, 1, 1, &A, &Gamma, &b, NULL, 1, &binterpt)); 3143606a31eSEmil Constantinescu } 3153606a31eSEmil Constantinescu 3163606a31eSEmil Constantinescu { 317bbd56ea5SKarl Rupp const PetscReal A = 0; 318bbd56ea5SKarl Rupp const PetscReal Gamma = 0.5; 319bbd56ea5SKarl Rupp const PetscReal b = 1; 320bbd56ea5SKarl Rupp const PetscReal binterpt = 1; 321bbd56ea5SKarl Rupp 3229566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSWTHETA2, 2, 1, &A, &Gamma, &b, NULL, 1, &binterpt)); 3233606a31eSEmil Constantinescu } 3243606a31eSEmil Constantinescu 3253606a31eSEmil Constantinescu { 326da80777bSKarl 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. */ 3279371c9d4SSatish Balay const PetscReal A[2][2] = 3289371c9d4SSatish Balay { 3299371c9d4SSatish Balay {0, 0}, 3309371c9d4SSatish Balay {1., 0} 3319371c9d4SSatish Balay }, 3329371c9d4SSatish Balay Gamma[2][2] = {{1.707106781186547524401, 0}, {-2. * 1.707106781186547524401, 1.707106781186547524401}}, b[2] = {0.5, 0.5}, b1[2] = {1.0, 0.0}; 3331f80e275SEmil Constantinescu PetscReal binterpt[2][2]; 334da80777bSKarl Rupp binterpt[0][0] = 1.707106781186547524401 - 1.0; 335da80777bSKarl Rupp binterpt[1][0] = 2.0 - 1.707106781186547524401; 336da80777bSKarl Rupp binterpt[0][1] = 1.707106781186547524401 - 1.5; 337da80777bSKarl Rupp binterpt[1][1] = 1.5 - 1.707106781186547524401; 338bbd56ea5SKarl Rupp 3399566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSW2P, 2, 2, &A[0][0], &Gamma[0][0], b, b1, 2, &binterpt[0][0])); 340e27a552bSJed Brown } 341e27a552bSJed Brown { 342da80777bSKarl 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. */ 3439371c9d4SSatish Balay const PetscReal A[2][2] = 3449371c9d4SSatish Balay { 3459371c9d4SSatish Balay {0, 0}, 3469371c9d4SSatish Balay {1., 0} 3479371c9d4SSatish Balay }, 3489371c9d4SSatish Balay Gamma[2][2] = {{0.2928932188134524755992, 0}, {-2. * 0.2928932188134524755992, 0.2928932188134524755992}}, b[2] = {0.5, 0.5}, b1[2] = {1.0, 0.0}; 3491f80e275SEmil Constantinescu PetscReal binterpt[2][2]; 350da80777bSKarl Rupp binterpt[0][0] = 0.2928932188134524755992 - 1.0; 351da80777bSKarl Rupp binterpt[1][0] = 2.0 - 0.2928932188134524755992; 352da80777bSKarl Rupp binterpt[0][1] = 0.2928932188134524755992 - 1.5; 353da80777bSKarl Rupp binterpt[1][1] = 1.5 - 0.2928932188134524755992; 354bbd56ea5SKarl Rupp 3559566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSW2M, 2, 2, &A[0][0], &Gamma[0][0], b, b1, 2, &binterpt[0][0])); 356fe7e6d57SJed Brown } 357fe7e6d57SJed Brown { 358da80777bSKarl Rupp /*const PetscReal g = 7.8867513459481287e-01; Directly written in-place below */ 3591f80e275SEmil Constantinescu PetscReal binterpt[3][2]; 3609371c9d4SSatish Balay const PetscReal A[3][3] = 3619371c9d4SSatish Balay { 3629371c9d4SSatish Balay {0, 0, 0}, 363fe7e6d57SJed Brown {1.5773502691896257e+00, 0, 0}, 3649371c9d4SSatish Balay {0.5, 0, 0} 3659371c9d4SSatish Balay }, 3669371c9d4SSatish 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}; 3671f80e275SEmil Constantinescu 3681f80e275SEmil Constantinescu binterpt[0][0] = -0.8094010767585034; 3691f80e275SEmil Constantinescu binterpt[1][0] = -0.5; 3701f80e275SEmil Constantinescu binterpt[2][0] = 2.3094010767585034; 3711f80e275SEmil Constantinescu binterpt[0][1] = 0.9641016151377548; 3721f80e275SEmil Constantinescu binterpt[1][1] = 0.5; 3731f80e275SEmil Constantinescu binterpt[2][1] = -1.4641016151377548; 374bbd56ea5SKarl Rupp 3759566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSWRA3PW, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0])); 376fe7e6d57SJed Brown } 377fe7e6d57SJed Brown { 3783ca35412SEmil Constantinescu PetscReal binterpt[4][3]; 379da80777bSKarl Rupp /*const PetscReal g = 4.3586652150845900e-01; Directly written in-place below */ 3809371c9d4SSatish Balay const PetscReal A[4][4] = 3819371c9d4SSatish Balay { 3829371c9d4SSatish Balay {0, 0, 0, 0}, 383fe7e6d57SJed Brown {8.7173304301691801e-01, 0, 0, 0}, 384fe7e6d57SJed Brown {8.4457060015369423e-01, -1.1299064236484185e-01, 0, 0}, 3859371c9d4SSatish Balay {0, 0, 1., 0} 3869371c9d4SSatish Balay }, 3879371c9d4SSatish 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}; 3883ca35412SEmil Constantinescu 3893ca35412SEmil Constantinescu binterpt[0][0] = 1.0564298455794094; 3903ca35412SEmil Constantinescu binterpt[1][0] = 2.296429974281067; 3913ca35412SEmil Constantinescu binterpt[2][0] = -1.307599564525376; 3923ca35412SEmil Constantinescu binterpt[3][0] = -1.045260255335102; 3933ca35412SEmil Constantinescu binterpt[0][1] = -1.3864882699759573; 3943ca35412SEmil Constantinescu binterpt[1][1] = -8.262611700275677; 3953ca35412SEmil Constantinescu binterpt[2][1] = 7.250979895056055; 3963ca35412SEmil Constantinescu binterpt[3][1] = 2.398120075195581; 3973ca35412SEmil Constantinescu binterpt[0][2] = 0.5721822314575016; 3983ca35412SEmil Constantinescu binterpt[1][2] = 4.742931142090097; 3993ca35412SEmil Constantinescu binterpt[2][2] = -4.398120075195578; 4003ca35412SEmil Constantinescu binterpt[3][2] = -0.9169932983520199; 4013ca35412SEmil Constantinescu 4029566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSWRA34PW2, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0])); 403e27a552bSJed Brown } 404ef3c5b88SJed Brown { 405da80777bSKarl Rupp /* const PetscReal g = 0.5; Directly written in-place below */ 4069371c9d4SSatish Balay const PetscReal A[4][4] = 4079371c9d4SSatish Balay { 4089371c9d4SSatish Balay {0, 0, 0, 0}, 409ef3c5b88SJed Brown {0, 0, 0, 0}, 410ef3c5b88SJed Brown {1., 0, 0, 0}, 4119371c9d4SSatish Balay {0.75, -0.25, 0.5, 0} 4129371c9d4SSatish Balay }, 4139371c9d4SSatish 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}; 414bbd56ea5SKarl Rupp 4159566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSWRODAS3, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 0, NULL)); 416ef3c5b88SJed Brown } 417ef3c5b88SJed Brown { 418da80777bSKarl Rupp /*const PetscReal g = 0.43586652150845899941601945119356; Directly written in-place below */ 4199371c9d4SSatish Balay const PetscReal A[3][3] = 4209371c9d4SSatish Balay { 4219371c9d4SSatish Balay {0, 0, 0}, 422da80777bSKarl Rupp {0.43586652150845899941601945119356, 0, 0}, 4239371c9d4SSatish Balay {0.43586652150845899941601945119356, 0, 0} 4249371c9d4SSatish Balay }, 4259371c9d4SSatish 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}; 4261f80e275SEmil Constantinescu 4271f80e275SEmil Constantinescu PetscReal binterpt[3][2]; 4281f80e275SEmil Constantinescu binterpt[0][0] = 3.793692883777660870425141387941; 4291f80e275SEmil Constantinescu binterpt[1][0] = -2.918692883777660870425141387941; 4301f80e275SEmil Constantinescu binterpt[2][0] = 0.125; 4311f80e275SEmil Constantinescu binterpt[0][1] = -0.725741064379812106687651020584; 4321f80e275SEmil Constantinescu binterpt[1][1] = 0.559074397713145440020984353917; 4331f80e275SEmil Constantinescu binterpt[2][1] = 0.16666666666666666666666666666667; 4341f80e275SEmil Constantinescu 4359566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSWSANDU3, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0])); 436ef3c5b88SJed Brown } 437b1c69cc3SEmil Constantinescu { 438da80777bSKarl Rupp /*const PetscReal s3 = PetscSqrtReal(3.),g = (3.0+s3)/6.0; 439da80777bSKarl Rupp * Direct evaluation: s3 = 1.732050807568877293527; 440da80777bSKarl Rupp * g = 0.7886751345948128822546; 441da80777bSKarl Rupp * Values are directly inserted below to ensure availability at compile time (compiler warnings otherwise...) */ 4429371c9d4SSatish Balay const PetscReal A[3][3] = 4439371c9d4SSatish Balay { 4449371c9d4SSatish Balay {0, 0, 0}, 445b1c69cc3SEmil Constantinescu {1, 0, 0}, 4469371c9d4SSatish Balay {0.25, 0.25, 0} 4479371c9d4SSatish Balay }, 4489371c9d4SSatish 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.}; 449c0cb691aSEmil Constantinescu PetscReal binterpt[3][2]; 450da80777bSKarl Rupp 451c0cb691aSEmil Constantinescu binterpt[0][0] = 0.089316397477040902157517886164709; 452c0cb691aSEmil Constantinescu binterpt[1][0] = -0.91068360252295909784248211383529; 453c0cb691aSEmil Constantinescu binterpt[2][0] = 1.8213672050459181956849642276706; 454c0cb691aSEmil Constantinescu binterpt[0][1] = 0.077350269189625764509148780501957; 455c0cb691aSEmil Constantinescu binterpt[1][1] = 1.077350269189625764509148780502; 456c0cb691aSEmil Constantinescu binterpt[2][1] = -1.1547005383792515290182975610039; 457bbd56ea5SKarl Rupp 4589566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSWASSP3P3S1C, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0])); 459b1c69cc3SEmil Constantinescu } 460b1c69cc3SEmil Constantinescu 461b1c69cc3SEmil Constantinescu { 4629371c9d4SSatish Balay const PetscReal A[4][4] = 4639371c9d4SSatish Balay { 4649371c9d4SSatish Balay {0, 0, 0, 0}, 465b1c69cc3SEmil Constantinescu {1. / 2., 0, 0, 0}, 466b1c69cc3SEmil Constantinescu {1. / 2., 1. / 2., 0, 0}, 4679371c9d4SSatish Balay {1. / 6., 1. / 6., 1. / 6., 0} 4689371c9d4SSatish Balay }, 4699371c9d4SSatish 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}; 470c0cb691aSEmil Constantinescu PetscReal binterpt[4][3]; 471da80777bSKarl Rupp 472c0cb691aSEmil Constantinescu binterpt[0][0] = 6.25; 473c0cb691aSEmil Constantinescu binterpt[1][0] = -30.25; 474c0cb691aSEmil Constantinescu binterpt[2][0] = 1.75; 475c0cb691aSEmil Constantinescu binterpt[3][0] = 23.25; 476c0cb691aSEmil Constantinescu binterpt[0][1] = -9.75; 477c0cb691aSEmil Constantinescu binterpt[1][1] = 58.75; 478c0cb691aSEmil Constantinescu binterpt[2][1] = -3.25; 479c0cb691aSEmil Constantinescu binterpt[3][1] = -45.75; 480c0cb691aSEmil Constantinescu binterpt[0][2] = 3.6666666666666666666666666666667; 481c0cb691aSEmil Constantinescu binterpt[1][2] = -28.333333333333333333333333333333; 482c0cb691aSEmil Constantinescu binterpt[2][2] = 1.6666666666666666666666666666667; 483c0cb691aSEmil Constantinescu binterpt[3][2] = 23.; 484bbd56ea5SKarl Rupp 4859566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSWLASSP3P4S2C, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0])); 486b1c69cc3SEmil Constantinescu } 487b1c69cc3SEmil Constantinescu 488b1c69cc3SEmil Constantinescu { 4899371c9d4SSatish Balay const PetscReal A[4][4] = 4909371c9d4SSatish Balay { 4919371c9d4SSatish Balay {0, 0, 0, 0}, 492b1c69cc3SEmil Constantinescu {1. / 2., 0, 0, 0}, 493b1c69cc3SEmil Constantinescu {1. / 2., 1. / 2., 0, 0}, 4949371c9d4SSatish Balay {1. / 6., 1. / 6., 1. / 6., 0} 4959371c9d4SSatish Balay }, 4969371c9d4SSatish 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}; 497c0cb691aSEmil Constantinescu PetscReal binterpt[4][3]; 498da80777bSKarl Rupp 499c0cb691aSEmil Constantinescu binterpt[0][0] = 1.6911764705882352941176470588235; 500c0cb691aSEmil Constantinescu binterpt[1][0] = 3.6813725490196078431372549019608; 501c0cb691aSEmil Constantinescu binterpt[2][0] = 0.23039215686274509803921568627451; 502c0cb691aSEmil Constantinescu binterpt[3][0] = -4.6029411764705882352941176470588; 503c0cb691aSEmil Constantinescu binterpt[0][1] = -0.95588235294117647058823529411765; 504c0cb691aSEmil Constantinescu binterpt[1][1] = -6.2401960784313725490196078431373; 505c0cb691aSEmil Constantinescu binterpt[2][1] = -0.31862745098039215686274509803922; 506c0cb691aSEmil Constantinescu binterpt[3][1] = 7.5147058823529411764705882352941; 507c0cb691aSEmil Constantinescu binterpt[0][2] = -0.56862745098039215686274509803922; 508c0cb691aSEmil Constantinescu binterpt[1][2] = 2.7254901960784313725490196078431; 509c0cb691aSEmil Constantinescu binterpt[2][2] = 0.25490196078431372549019607843137; 510c0cb691aSEmil Constantinescu binterpt[3][2] = -2.4117647058823529411764705882353; 511bbd56ea5SKarl Rupp 5129566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSWLLSSP3P4S2C, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0])); 513b1c69cc3SEmil Constantinescu } 514753f8adbSEmil Constantinescu 515753f8adbSEmil Constantinescu { 516753f8adbSEmil Constantinescu PetscReal A[4][4], Gamma[4][4], b[4], b2[4]; 5173ca35412SEmil Constantinescu PetscReal binterpt[4][3]; 518753f8adbSEmil Constantinescu 519753f8adbSEmil Constantinescu Gamma[0][0] = 0.4358665215084589994160194475295062513822671686978816; 5209371c9d4SSatish Balay Gamma[0][1] = 0; 5219371c9d4SSatish Balay Gamma[0][2] = 0; 5229371c9d4SSatish Balay Gamma[0][3] = 0; 523753f8adbSEmil Constantinescu Gamma[1][0] = -1.997527830934941248426324674704153457289527280554476; 524753f8adbSEmil Constantinescu Gamma[1][1] = 0.4358665215084589994160194475295062513822671686978816; 5259371c9d4SSatish Balay Gamma[1][2] = 0; 5269371c9d4SSatish Balay Gamma[1][3] = 0; 527753f8adbSEmil Constantinescu Gamma[2][0] = -1.007948511795029620852002345345404191008352770119903; 528753f8adbSEmil Constantinescu Gamma[2][1] = -0.004648958462629345562774289390054679806993396798458131; 529753f8adbSEmil Constantinescu Gamma[2][2] = 0.4358665215084589994160194475295062513822671686978816; 53005e8e825SJed Brown Gamma[2][3] = 0; 531753f8adbSEmil Constantinescu Gamma[3][0] = -0.6685429734233467180451604600279552604364311322650783; 532753f8adbSEmil Constantinescu Gamma[3][1] = 0.6056625986449338476089525334450053439525178740492984; 533753f8adbSEmil Constantinescu Gamma[3][2] = -0.9717899277217721234705114616271378792182450260943198; 534753f8adbSEmil Constantinescu Gamma[3][3] = 0; 535753f8adbSEmil Constantinescu 5369371c9d4SSatish Balay A[0][0] = 0; 5379371c9d4SSatish Balay A[0][1] = 0; 5389371c9d4SSatish Balay A[0][2] = 0; 5399371c9d4SSatish Balay A[0][3] = 0; 540753f8adbSEmil Constantinescu A[1][0] = 0.8717330430169179988320388950590125027645343373957631; 5419371c9d4SSatish Balay A[1][1] = 0; 5429371c9d4SSatish Balay A[1][2] = 0; 5439371c9d4SSatish Balay A[1][3] = 0; 544753f8adbSEmil Constantinescu A[2][0] = 0.5275890119763004115618079766722914408876108660811028; 545753f8adbSEmil Constantinescu A[2][1] = 0.07241098802369958843819203208518599088698057726988732; 5469371c9d4SSatish Balay A[2][2] = 0; 5479371c9d4SSatish Balay A[2][3] = 0; 548753f8adbSEmil Constantinescu A[3][0] = 0.3990960076760701320627260685975778145384666450351314; 549753f8adbSEmil Constantinescu A[3][1] = -0.4375576546135194437228463747348862825846903771419953; 550753f8adbSEmil Constantinescu A[3][2] = 1.038461646937449311660120300601880176655352737312713; 55105e8e825SJed Brown A[3][3] = 0; 552753f8adbSEmil Constantinescu 553753f8adbSEmil Constantinescu b[0] = 0.1876410243467238251612921333138006734899663569186926; 554753f8adbSEmil Constantinescu b[1] = -0.5952974735769549480478230473706443582188442040780541; 555753f8adbSEmil Constantinescu b[2] = 0.9717899277217721234705114616271378792182450260943198; 556753f8adbSEmil Constantinescu b[3] = 0.4358665215084589994160194475295062513822671686978816; 557753f8adbSEmil Constantinescu 558753f8adbSEmil Constantinescu b2[0] = 0.2147402862233891404862383521089097657790734483804460; 559753f8adbSEmil Constantinescu b2[1] = -0.4851622638849390928209050538171743017757490232519684; 560753f8adbSEmil Constantinescu b2[2] = 0.8687250025203875511662123688667549217531982787600080; 561753f8adbSEmil Constantinescu b2[3] = 0.4016969751411624011684543450940068201770721128357014; 562753f8adbSEmil Constantinescu 5633ca35412SEmil Constantinescu binterpt[0][0] = 2.2565812720167954547104627844105; 5643ca35412SEmil Constantinescu binterpt[1][0] = 1.349166413351089573796243820819; 5653ca35412SEmil Constantinescu binterpt[2][0] = -2.4695174540533503758652847586647; 5663ca35412SEmil Constantinescu binterpt[3][0] = -0.13623023131453465264142184656474; 5673ca35412SEmil Constantinescu binterpt[0][1] = -3.0826699111559187902922463354557; 5683ca35412SEmil Constantinescu binterpt[1][1] = -2.4689115685996042534544925650515; 5693ca35412SEmil Constantinescu binterpt[2][1] = 5.7428279814696677152129332773553; 5703ca35412SEmil Constantinescu binterpt[3][1] = -0.19124650171414467146619437684812; 5713ca35412SEmil Constantinescu binterpt[0][2] = 1.0137296634858471607430756831148; 5723ca35412SEmil Constantinescu binterpt[1][2] = 0.52444768167155973161042570784064; 5733ca35412SEmil Constantinescu binterpt[2][2] = -2.3015205996945452158771370439586; 5743ca35412SEmil Constantinescu binterpt[3][2] = 0.76334325453713832352363565300308; 575f4aed992SEmil Constantinescu 5769566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(TSROSWARK3, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0])); 577753f8adbSEmil Constantinescu } 5789566063dSJacob Faibussowitsch PetscCall(TSRosWRegisterRos4(TSROSWGRK4T, 0.231, PETSC_DEFAULT, PETSC_DEFAULT, 0, -0.1282612945269037e+01)); 5799566063dSJacob Faibussowitsch PetscCall(TSRosWRegisterRos4(TSROSWSHAMP4, 0.5, PETSC_DEFAULT, PETSC_DEFAULT, 0, 125. / 108.)); 5809566063dSJacob Faibussowitsch PetscCall(TSRosWRegisterRos4(TSROSWVELDD4, 0.22570811482256823492, PETSC_DEFAULT, PETSC_DEFAULT, 0, -1.355958941201148)); 5819566063dSJacob Faibussowitsch PetscCall(TSRosWRegisterRos4(TSROSW4L, 0.57282, PETSC_DEFAULT, PETSC_DEFAULT, 0, -1.093502252409163)); 582e27a552bSJed Brown PetscFunctionReturn(0); 583e27a552bSJed Brown } 584e27a552bSJed Brown 585e27a552bSJed Brown /*@C 586e27a552bSJed Brown TSRosWRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister(). 587e27a552bSJed Brown 588e27a552bSJed Brown Not Collective 589e27a552bSJed Brown 590e27a552bSJed Brown Level: advanced 591e27a552bSJed Brown 592db781477SPatrick Sanan .seealso: `TSRosWRegister()`, `TSRosWRegisterAll()` 593e27a552bSJed Brown @*/ 5949371c9d4SSatish Balay PetscErrorCode TSRosWRegisterDestroy(void) { 59561692a83SJed Brown RosWTableauLink link; 596e27a552bSJed Brown 597e27a552bSJed Brown PetscFunctionBegin; 59861692a83SJed Brown while ((link = RosWTableauList)) { 59961692a83SJed Brown RosWTableau t = &link->tab; 60061692a83SJed Brown RosWTableauList = link->next; 6019566063dSJacob Faibussowitsch PetscCall(PetscFree5(t->A, t->Gamma, t->b, t->ASum, t->GammaSum)); 6029566063dSJacob Faibussowitsch PetscCall(PetscFree5(t->At, t->bt, t->GammaInv, t->GammaZeroDiag, t->GammaExplicitCorr)); 6039566063dSJacob Faibussowitsch PetscCall(PetscFree2(t->bembed, t->bembedt)); 6049566063dSJacob Faibussowitsch PetscCall(PetscFree(t->binterpt)); 6059566063dSJacob Faibussowitsch PetscCall(PetscFree(t->name)); 6069566063dSJacob Faibussowitsch PetscCall(PetscFree(link)); 607e27a552bSJed Brown } 608e27a552bSJed Brown TSRosWRegisterAllCalled = PETSC_FALSE; 609e27a552bSJed Brown PetscFunctionReturn(0); 610e27a552bSJed Brown } 611e27a552bSJed Brown 612e27a552bSJed Brown /*@C 613e27a552bSJed Brown TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called 6148a690491SBarry Smith from TSInitializePackage(). 615e27a552bSJed Brown 616e27a552bSJed Brown Level: developer 617e27a552bSJed Brown 618db781477SPatrick Sanan .seealso: `PetscInitialize()` 619e27a552bSJed Brown @*/ 6209371c9d4SSatish Balay PetscErrorCode TSRosWInitializePackage(void) { 621e27a552bSJed Brown PetscFunctionBegin; 622e27a552bSJed Brown if (TSRosWPackageInitialized) PetscFunctionReturn(0); 623e27a552bSJed Brown TSRosWPackageInitialized = PETSC_TRUE; 6249566063dSJacob Faibussowitsch PetscCall(TSRosWRegisterAll()); 6259566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSRosWFinalizePackage)); 626e27a552bSJed Brown PetscFunctionReturn(0); 627e27a552bSJed Brown } 628e27a552bSJed Brown 629e27a552bSJed Brown /*@C 630e27a552bSJed Brown TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is 631e27a552bSJed Brown called from PetscFinalize(). 632e27a552bSJed Brown 633e27a552bSJed Brown Level: developer 634e27a552bSJed Brown 635db781477SPatrick Sanan .seealso: `PetscFinalize()` 636e27a552bSJed Brown @*/ 6379371c9d4SSatish Balay PetscErrorCode TSRosWFinalizePackage(void) { 638e27a552bSJed Brown PetscFunctionBegin; 639e27a552bSJed Brown TSRosWPackageInitialized = PETSC_FALSE; 6409566063dSJacob Faibussowitsch PetscCall(TSRosWRegisterDestroy()); 641e27a552bSJed Brown PetscFunctionReturn(0); 642e27a552bSJed Brown } 643e27a552bSJed Brown 644e27a552bSJed Brown /*@C 64561692a83SJed Brown TSRosWRegister - register a Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 646e27a552bSJed Brown 647e27a552bSJed Brown Not Collective, but the same schemes should be registered on all processes on which they will be used 648e27a552bSJed Brown 649e27a552bSJed Brown Input Parameters: 650e27a552bSJed Brown + name - identifier for method 651e27a552bSJed Brown . order - approximation order of method 652e27a552bSJed Brown . s - number of stages, this is the dimension of the matrices below 65361692a83SJed Brown . A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular 65461692a83SJed Brown . Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal 655fe7e6d57SJed Brown . b - Step completion table (dimension s) 6560298fd71SBarry Smith . bembed - Step completion table for a scheme of order one less (dimension s, NULL if no embedded scheme is available) 657f4aed992SEmil Constantinescu . pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt 65842faf41dSJed Brown - binterpt - Coefficients of the interpolation formula (dimension s*pinterp) 659e27a552bSJed Brown 660e27a552bSJed Brown Notes: 66161692a83SJed Brown Several Rosenbrock W methods are provided, this function is only needed to create new methods. 662e27a552bSJed Brown 663e27a552bSJed Brown Level: advanced 664e27a552bSJed Brown 665db781477SPatrick Sanan .seealso: `TSRosW` 666e27a552bSJed Brown @*/ 6679371c9d4SSatish Balay 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[]) { 66861692a83SJed Brown RosWTableauLink link; 66961692a83SJed Brown RosWTableau t; 67061692a83SJed Brown PetscInt i, j, k; 67161692a83SJed Brown PetscScalar *GammaInv; 672e27a552bSJed Brown 673e27a552bSJed Brown PetscFunctionBegin; 674fe7e6d57SJed Brown PetscValidCharPointer(name, 1); 675dadcf809SJacob Faibussowitsch PetscValidRealPointer(A, 4); 676dadcf809SJacob Faibussowitsch PetscValidRealPointer(Gamma, 5); 677dadcf809SJacob Faibussowitsch PetscValidRealPointer(b, 6); 678dadcf809SJacob Faibussowitsch if (bembed) PetscValidRealPointer(bembed, 7); 679fe7e6d57SJed Brown 6809566063dSJacob Faibussowitsch PetscCall(TSRosWInitializePackage()); 6819566063dSJacob Faibussowitsch PetscCall(PetscNew(&link)); 682e27a552bSJed Brown t = &link->tab; 6839566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &t->name)); 684e27a552bSJed Brown t->order = order; 685e27a552bSJed Brown t->s = s; 6869566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(s * s, &t->A, s * s, &t->Gamma, s, &t->b, s, &t->ASum, s, &t->GammaSum)); 6879566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(s * s, &t->At, s, &t->bt, s * s, &t->GammaInv, s, &t->GammaZeroDiag, s * s, &t->GammaExplicitCorr)); 6889566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->A, A, s * s)); 6899566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->Gamma, Gamma, s * s)); 6909566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->GammaExplicitCorr, Gamma, s * s)); 6919566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->b, b, s)); 692fe7e6d57SJed Brown if (bembed) { 6939566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(s, &t->bembed, s, &t->bembedt)); 6949566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->bembed, bembed, s)); 695fe7e6d57SJed Brown } 69661692a83SJed Brown for (i = 0; i < s; i++) { 69761692a83SJed Brown t->ASum[i] = 0; 69861692a83SJed Brown t->GammaSum[i] = 0; 69961692a83SJed Brown for (j = 0; j < s; j++) { 70061692a83SJed Brown t->ASum[i] += A[i * s + j]; 701fe7e6d57SJed Brown t->GammaSum[i] += Gamma[i * s + j]; 70261692a83SJed Brown } 70361692a83SJed Brown } 7049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s * s, &GammaInv)); /* Need to use Scalar for inverse, then convert back to Real */ 70561692a83SJed Brown for (i = 0; i < s * s; i++) GammaInv[i] = Gamma[i]; 706fd96d5b0SEmil Constantinescu for (i = 0; i < s; i++) { 707fd96d5b0SEmil Constantinescu if (Gamma[i * s + i] == 0.0) { 708fd96d5b0SEmil Constantinescu GammaInv[i * s + i] = 1.0; 709c17803e7SJed Brown t->GammaZeroDiag[i] = PETSC_TRUE; 710fd96d5b0SEmil Constantinescu } else { 711c17803e7SJed Brown t->GammaZeroDiag[i] = PETSC_FALSE; 712fd96d5b0SEmil Constantinescu } 713fd96d5b0SEmil Constantinescu } 714fd96d5b0SEmil Constantinescu 71561692a83SJed Brown switch (s) { 71661692a83SJed Brown case 1: GammaInv[0] = 1. / GammaInv[0]; break; 7179566063dSJacob Faibussowitsch case 2: PetscCall(PetscKernel_A_gets_inverse_A_2(GammaInv, 0, PETSC_FALSE, NULL)); break; 7189566063dSJacob Faibussowitsch case 3: PetscCall(PetscKernel_A_gets_inverse_A_3(GammaInv, 0, PETSC_FALSE, NULL)); break; 7199566063dSJacob Faibussowitsch case 4: PetscCall(PetscKernel_A_gets_inverse_A_4(GammaInv, 0, PETSC_FALSE, NULL)); break; 72061692a83SJed Brown case 5: { 72161692a83SJed Brown PetscInt ipvt5[5]; 72261692a83SJed Brown MatScalar work5[5 * 5]; 7239371c9d4SSatish Balay PetscCall(PetscKernel_A_gets_inverse_A_5(GammaInv, ipvt5, work5, 0, PETSC_FALSE, NULL)); 7249371c9d4SSatish Balay break; 72561692a83SJed Brown } 7269566063dSJacob Faibussowitsch case 6: PetscCall(PetscKernel_A_gets_inverse_A_6(GammaInv, 0, PETSC_FALSE, NULL)); break; 7279566063dSJacob Faibussowitsch case 7: PetscCall(PetscKernel_A_gets_inverse_A_7(GammaInv, 0, PETSC_FALSE, NULL)); break; 72863a3b9bcSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not implemented for %" PetscInt_FMT " stages", s); 72961692a83SJed Brown } 73061692a83SJed Brown for (i = 0; i < s * s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]); 7319566063dSJacob Faibussowitsch PetscCall(PetscFree(GammaInv)); 73243b21953SEmil Constantinescu 73343b21953SEmil Constantinescu for (i = 0; i < s; i++) { 73443b21953SEmil Constantinescu for (k = 0; k < i + 1; k++) { 73543b21953SEmil Constantinescu t->GammaExplicitCorr[i * s + k] = (t->GammaExplicitCorr[i * s + k]) * (t->GammaInv[k * s + k]); 7369371c9d4SSatish Balay for (j = k + 1; j < i + 1; j++) { t->GammaExplicitCorr[i * s + k] += (t->GammaExplicitCorr[i * s + j]) * (t->GammaInv[j * s + k]); } 73743b21953SEmil Constantinescu } 73843b21953SEmil Constantinescu } 73943b21953SEmil Constantinescu 74061692a83SJed Brown for (i = 0; i < s; i++) { 74161692a83SJed Brown for (j = 0; j < s; j++) { 74261692a83SJed Brown t->At[i * s + j] = 0; 7439371c9d4SSatish Balay for (k = 0; k < s; k++) { t->At[i * s + j] += t->A[i * s + k] * t->GammaInv[k * s + j]; } 74461692a83SJed Brown } 74561692a83SJed Brown t->bt[i] = 0; 7469371c9d4SSatish Balay for (j = 0; j < s; j++) { t->bt[i] += t->b[j] * t->GammaInv[j * s + i]; } 747fe7e6d57SJed Brown if (bembed) { 748fe7e6d57SJed Brown t->bembedt[i] = 0; 7499371c9d4SSatish Balay for (j = 0; j < s; j++) { t->bembedt[i] += t->bembed[j] * t->GammaInv[j * s + i]; } 750fe7e6d57SJed Brown } 75161692a83SJed Brown } 7528d59e960SJed Brown t->ccfl = 1.0; /* Fix this */ 7538d59e960SJed Brown 754f4aed992SEmil Constantinescu t->pinterp = pinterp; 7559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s * pinterp, &t->binterpt)); 7569566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->binterpt, binterpt, s * pinterp)); 75761692a83SJed Brown link->next = RosWTableauList; 75861692a83SJed Brown RosWTableauList = link; 759e27a552bSJed Brown PetscFunctionReturn(0); 760e27a552bSJed Brown } 761e27a552bSJed Brown 76242faf41dSJed Brown /*@C 763fd292e60Sprj- TSRosWRegisterRos4 - register a fourth order Rosenbrock scheme by providing parameter choices 76442faf41dSJed Brown 76542faf41dSJed Brown Not Collective, but the same schemes should be registered on all processes on which they will be used 76642faf41dSJed Brown 76742faf41dSJed Brown Input Parameters: 76842faf41dSJed Brown + name - identifier for method 76942faf41dSJed Brown . gamma - leading coefficient (diagonal entry) 77042faf41dSJed Brown . a2 - design parameter, see Table 7.2 of Hairer&Wanner 77142faf41dSJed Brown . a3 - design parameter or PETSC_DEFAULT to satisfy one of the order five conditions (Eq 7.22) 77242faf41dSJed Brown . b3 - design parameter, see Table 7.2 of Hairer&Wanner 77342faf41dSJed Brown . beta43 - design parameter or PETSC_DEFAULT to use Equation 7.21 of Hairer&Wanner 774a2b725a8SWilliam Gropp - e4 - design parameter for embedded method, see coefficient E4 in ros4.f code from Hairer 77542faf41dSJed Brown 77642faf41dSJed Brown Notes: 77742faf41dSJed Brown This routine encodes the design of fourth order Rosenbrock methods as described in Hairer and Wanner volume 2. 77842faf41dSJed Brown It is used here to implement several methods from the book and can be used to experiment with new methods. 77942faf41dSJed Brown It was written this way instead of by copying coefficients in order to provide better than double precision satisfaction of the order conditions. 78042faf41dSJed Brown 78142faf41dSJed Brown Level: developer 78242faf41dSJed Brown 783db781477SPatrick Sanan .seealso: `TSRosW`, `TSRosWRegister()` 78442faf41dSJed Brown @*/ 7859371c9d4SSatish Balay PetscErrorCode TSRosWRegisterRos4(TSRosWType name, PetscReal gamma, PetscReal a2, PetscReal a3, PetscReal b3, PetscReal e4) { 78642faf41dSJed Brown /* Declare numeric constants so they can be quad precision without being truncated at double */ 7879371c9d4SSatish 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; 78842faf41dSJed Brown PetscReal a4, a32, a42, a43, b1, b2, b4, beta2p, beta3p, beta4p, beta32, beta42, beta43, beta32beta2p, beta4jbetajp; 78942faf41dSJed Brown PetscReal A[4][4], Gamma[4][4], b[4], bm[4]; 79042faf41dSJed Brown PetscScalar M[3][3], rhs[3]; 79142faf41dSJed Brown 79242faf41dSJed Brown PetscFunctionBegin; 79342faf41dSJed Brown /* Step 1: choose Gamma (input) */ 79442faf41dSJed Brown /* Step 2: choose a2,a3,a4; b1,b2,b3,b4 to satisfy order conditions */ 79542faf41dSJed Brown if (a3 == PETSC_DEFAULT) a3 = (one / five - a2 / four) / (one / four - a2 / three); /* Eq 7.22 */ 79642faf41dSJed Brown a4 = a3; /* consequence of 7.20 */ 79742faf41dSJed Brown 79842faf41dSJed Brown /* Solve order conditions 7.15a, 7.15c, 7.15e */ 7999371c9d4SSatish Balay M[0][0] = one; 8009371c9d4SSatish Balay M[0][1] = one; 8019371c9d4SSatish Balay M[0][2] = one; /* 7.15a */ 8029371c9d4SSatish Balay M[1][0] = 0.0; 8039371c9d4SSatish Balay M[1][1] = a2 * a2; 8049371c9d4SSatish Balay M[1][2] = a4 * a4; /* 7.15c */ 8059371c9d4SSatish Balay M[2][0] = 0.0; 8069371c9d4SSatish Balay M[2][1] = a2 * a2 * a2; 8079371c9d4SSatish Balay M[2][2] = a4 * a4 * a4; /* 7.15e */ 80842faf41dSJed Brown rhs[0] = one - b3; 80942faf41dSJed Brown rhs[1] = one / three - a3 * a3 * b3; 81042faf41dSJed Brown rhs[2] = one / four - a3 * a3 * a3 * b3; 8119566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_3(&M[0][0], 0, PETSC_FALSE, NULL)); 81242faf41dSJed Brown b1 = PetscRealPart(M[0][0] * rhs[0] + M[0][1] * rhs[1] + M[0][2] * rhs[2]); 81342faf41dSJed Brown b2 = PetscRealPart(M[1][0] * rhs[0] + M[1][1] * rhs[1] + M[1][2] * rhs[2]); 81442faf41dSJed Brown b4 = PetscRealPart(M[2][0] * rhs[0] + M[2][1] * rhs[1] + M[2][2] * rhs[2]); 81542faf41dSJed Brown 81642faf41dSJed Brown /* Step 3 */ 81742faf41dSJed Brown beta43 = (p56 - a2 * p43) / (b4 * a3 * a3 * (a3 - a2)); /* 7.21 */ 81842faf41dSJed Brown beta32beta2p = p44 / (b4 * beta43); /* 7.15h */ 81942faf41dSJed Brown beta4jbetajp = (p32 - b3 * beta32beta2p) / b4; 8209371c9d4SSatish Balay M[0][0] = b2; 8219371c9d4SSatish Balay M[0][1] = b3; 8229371c9d4SSatish Balay M[0][2] = b4; 8239371c9d4SSatish Balay M[1][0] = a4 * a4 * beta32beta2p - a3 * a3 * beta4jbetajp; 8249371c9d4SSatish Balay M[1][1] = a2 * a2 * beta4jbetajp; 8259371c9d4SSatish Balay M[1][2] = -a2 * a2 * beta32beta2p; 8269371c9d4SSatish Balay M[2][0] = b4 * beta43 * a3 * a3 - p43; 8279371c9d4SSatish Balay M[2][1] = -b4 * beta43 * a2 * a2; 8289371c9d4SSatish Balay M[2][2] = 0; 8299371c9d4SSatish Balay rhs[0] = one / two - gamma; 8309371c9d4SSatish Balay rhs[1] = 0; 8319371c9d4SSatish Balay rhs[2] = -a2 * a2 * p32; 8329566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_3(&M[0][0], 0, PETSC_FALSE, NULL)); 83342faf41dSJed Brown beta2p = PetscRealPart(M[0][0] * rhs[0] + M[0][1] * rhs[1] + M[0][2] * rhs[2]); 83442faf41dSJed Brown beta3p = PetscRealPart(M[1][0] * rhs[0] + M[1][1] * rhs[1] + M[1][2] * rhs[2]); 83542faf41dSJed Brown beta4p = PetscRealPart(M[2][0] * rhs[0] + M[2][1] * rhs[1] + M[2][2] * rhs[2]); 83642faf41dSJed Brown 83742faf41dSJed Brown /* Step 4: back-substitute */ 83842faf41dSJed Brown beta32 = beta32beta2p / beta2p; 83942faf41dSJed Brown beta42 = (beta4jbetajp - beta43 * beta3p) / beta2p; 84042faf41dSJed Brown 84142faf41dSJed Brown /* Step 5: 7.15f and 7.20, then 7.16 */ 84242faf41dSJed Brown a43 = 0; 84342faf41dSJed Brown a32 = p42 / (b3 * a3 * beta2p + b4 * a4 * beta2p); 84442faf41dSJed Brown a42 = a32; 84542faf41dSJed Brown 8469371c9d4SSatish Balay A[0][0] = 0; 8479371c9d4SSatish Balay A[0][1] = 0; 8489371c9d4SSatish Balay A[0][2] = 0; 8499371c9d4SSatish Balay A[0][3] = 0; 8509371c9d4SSatish Balay A[1][0] = a2; 8519371c9d4SSatish Balay A[1][1] = 0; 8529371c9d4SSatish Balay A[1][2] = 0; 8539371c9d4SSatish Balay A[1][3] = 0; 8549371c9d4SSatish Balay A[2][0] = a3 - a32; 8559371c9d4SSatish Balay A[2][1] = a32; 8569371c9d4SSatish Balay A[2][2] = 0; 8579371c9d4SSatish Balay A[2][3] = 0; 8589371c9d4SSatish Balay A[3][0] = a4 - a43 - a42; 8599371c9d4SSatish Balay A[3][1] = a42; 8609371c9d4SSatish Balay A[3][2] = a43; 8619371c9d4SSatish Balay A[3][3] = 0; 8629371c9d4SSatish Balay Gamma[0][0] = gamma; 8639371c9d4SSatish Balay Gamma[0][1] = 0; 8649371c9d4SSatish Balay Gamma[0][2] = 0; 8659371c9d4SSatish Balay Gamma[0][3] = 0; 8669371c9d4SSatish Balay Gamma[1][0] = beta2p - A[1][0]; 8679371c9d4SSatish Balay Gamma[1][1] = gamma; 8689371c9d4SSatish Balay Gamma[1][2] = 0; 8699371c9d4SSatish Balay Gamma[1][3] = 0; 8709371c9d4SSatish Balay Gamma[2][0] = beta3p - beta32 - A[2][0]; 8719371c9d4SSatish Balay Gamma[2][1] = beta32 - A[2][1]; 8729371c9d4SSatish Balay Gamma[2][2] = gamma; 8739371c9d4SSatish Balay Gamma[2][3] = 0; 8749371c9d4SSatish Balay Gamma[3][0] = beta4p - beta42 - beta43 - A[3][0]; 8759371c9d4SSatish Balay Gamma[3][1] = beta42 - A[3][1]; 8769371c9d4SSatish Balay Gamma[3][2] = beta43 - A[3][2]; 8779371c9d4SSatish Balay Gamma[3][3] = gamma; 8789371c9d4SSatish Balay b[0] = b1; 8799371c9d4SSatish Balay b[1] = b2; 8809371c9d4SSatish Balay b[2] = b3; 8819371c9d4SSatish Balay b[3] = b4; 88242faf41dSJed Brown 88342faf41dSJed Brown /* Construct embedded formula using given e4. We are solving Equation 7.18. */ 88442faf41dSJed Brown bm[3] = b[3] - e4 * gamma; /* using definition of E4 */ 88542faf41dSJed Brown bm[2] = (p32 - beta4jbetajp * bm[3]) / (beta32 * beta2p); /* fourth row of 7.18 */ 88642faf41dSJed Brown bm[1] = (one / two - gamma - beta3p * bm[2] - beta4p * bm[3]) / beta2p; /* second row */ 88742faf41dSJed Brown bm[0] = one - bm[1] - bm[2] - bm[3]; /* first row */ 88842faf41dSJed Brown 88942faf41dSJed Brown { 89042faf41dSJed Brown const PetscReal misfit = a2 * a2 * bm[1] + a3 * a3 * bm[2] + a4 * a4 * bm[3] - one / three; 8913c633725SBarry Smith PetscCheck(PetscAbs(misfit) <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_SUP, "Assumptions violated, could not construct a third order embedded method"); 89242faf41dSJed Brown } 8939566063dSJacob Faibussowitsch PetscCall(TSRosWRegister(name, 4, 4, &A[0][0], &Gamma[0][0], b, bm, 0, NULL)); 89442faf41dSJed Brown PetscFunctionReturn(0); 89542faf41dSJed Brown } 89642faf41dSJed Brown 8971c3436cfSJed Brown /* 8981c3436cfSJed Brown The step completion formula is 8991c3436cfSJed Brown 9001c3436cfSJed Brown x1 = x0 + b^T Y 9011c3436cfSJed Brown 9021c3436cfSJed Brown where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been 9031c3436cfSJed Brown updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write 9041c3436cfSJed Brown 9051c3436cfSJed Brown x1e = x0 + be^T Y 9061c3436cfSJed Brown = x1 - b^T Y + be^T Y 9071c3436cfSJed Brown = x1 + (be - b)^T Y 9081c3436cfSJed Brown 9091c3436cfSJed Brown so we can evaluate the method of different order even after the step has been optimistically completed. 9101c3436cfSJed Brown */ 9119371c9d4SSatish Balay static PetscErrorCode TSEvaluateStep_RosW(TS ts, PetscInt order, Vec U, PetscBool *done) { 9121c3436cfSJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 9131c3436cfSJed Brown RosWTableau tab = ros->tableau; 9141c3436cfSJed Brown PetscScalar *w = ros->work; 9151c3436cfSJed Brown PetscInt i; 9161c3436cfSJed Brown 9171c3436cfSJed Brown PetscFunctionBegin; 9181c3436cfSJed Brown if (order == tab->order) { 919108c343cSJed Brown if (ros->status == TS_STEP_INCOMPLETE) { /* Use standard completion formula */ 9209566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, U)); 921de19f811SJed Brown for (i = 0; i < tab->s; i++) w[i] = tab->bt[i]; 9229566063dSJacob Faibussowitsch PetscCall(VecMAXPY(U, tab->s, w, ros->Y)); 9239566063dSJacob Faibussowitsch } else PetscCall(VecCopy(ts->vec_sol, U)); 9241c3436cfSJed Brown if (done) *done = PETSC_TRUE; 9251c3436cfSJed Brown PetscFunctionReturn(0); 9261c3436cfSJed Brown } else if (order == tab->order - 1) { 9271c3436cfSJed Brown if (!tab->bembedt) goto unavailable; 928108c343cSJed Brown if (ros->status == TS_STEP_INCOMPLETE) { /* Use embedded completion formula */ 9299566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, U)); 930de19f811SJed Brown for (i = 0; i < tab->s; i++) w[i] = tab->bembedt[i]; 9319566063dSJacob Faibussowitsch PetscCall(VecMAXPY(U, tab->s, w, ros->Y)); 932108c343cSJed Brown } else { /* Use rollback-and-recomplete formula (bembedt - bt) */ 933108c343cSJed Brown for (i = 0; i < tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i]; 9349566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, U)); 9359566063dSJacob Faibussowitsch PetscCall(VecMAXPY(U, tab->s, w, ros->Y)); 9361c3436cfSJed Brown } 9371c3436cfSJed Brown if (done) *done = PETSC_TRUE; 9381c3436cfSJed Brown PetscFunctionReturn(0); 9391c3436cfSJed Brown } 9401c3436cfSJed Brown unavailable: 9411c3436cfSJed Brown if (done) *done = PETSC_FALSE; 9429371c9d4SSatish Balay else 9439371c9d4SSatish 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, 9449371c9d4SSatish Balay tab->order, order); 9451c3436cfSJed Brown PetscFunctionReturn(0); 9461c3436cfSJed Brown } 9471c3436cfSJed Brown 9489371c9d4SSatish Balay static PetscErrorCode TSRollBack_RosW(TS ts) { 94924655328SShri TS_RosW *ros = (TS_RosW *)ts->data; 95024655328SShri 95124655328SShri PetscFunctionBegin; 9529566063dSJacob Faibussowitsch PetscCall(VecCopy(ros->vec_sol_prev, ts->vec_sol)); 95324655328SShri PetscFunctionReturn(0); 95424655328SShri } 95524655328SShri 9569371c9d4SSatish Balay static PetscErrorCode TSStep_RosW(TS ts) { 95761692a83SJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 95861692a83SJed Brown RosWTableau tab = ros->tableau; 959e27a552bSJed Brown const PetscInt s = tab->s; 9601c3436cfSJed Brown const PetscReal *At = tab->At, *Gamma = tab->Gamma, *ASum = tab->ASum, *GammaInv = tab->GammaInv; 9610feba352SEmil Constantinescu const PetscReal *GammaExplicitCorr = tab->GammaExplicitCorr; 962c17803e7SJed Brown const PetscBool *GammaZeroDiag = tab->GammaZeroDiag; 96361692a83SJed Brown PetscScalar *w = ros->work; 9647d4bf2deSEmil Constantinescu Vec *Y = ros->Y, Ydot = ros->Ydot, Zdot = ros->Zdot, Zstage = ros->Zstage; 965e27a552bSJed Brown SNES snes; 9661c3436cfSJed Brown TSAdapt adapt; 967fecfb714SLisandro Dalcin PetscInt i, j, its, lits; 968be5899b3SLisandro Dalcin PetscInt rejections = 0; 969b39943a6SLisandro Dalcin PetscBool stageok, accept = PETSC_TRUE; 970b39943a6SLisandro Dalcin PetscReal next_time_step = ts->time_step; 971f7f07198SBarry Smith PetscInt lag; 972e27a552bSJed Brown 973e27a552bSJed Brown PetscFunctionBegin; 974*48a46eb9SPierre Jolivet if (!ts->steprollback) PetscCall(VecCopy(ts->vec_sol, ros->vec_sol_prev)); 975e27a552bSJed Brown 976b39943a6SLisandro Dalcin ros->status = TS_STEP_INCOMPLETE; 977b39943a6SLisandro Dalcin while (!ts->reason && ros->status != TS_STEP_COMPLETE) { 9781c3436cfSJed Brown const PetscReal h = ts->time_step; 979e27a552bSJed Brown for (i = 0; i < s; i++) { 9801c3436cfSJed Brown ros->stage_time = ts->ptime + h * ASum[i]; 9819566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, ros->stage_time)); 982c17803e7SJed Brown if (GammaZeroDiag[i]) { 983c17803e7SJed Brown ros->stage_explicit = PETSC_TRUE; 984b296d7d5SJed Brown ros->scoeff = 1.; 985c17803e7SJed Brown } else { 986c17803e7SJed Brown ros->stage_explicit = PETSC_FALSE; 987b296d7d5SJed Brown ros->scoeff = 1. / Gamma[i * s + i]; 988fd96d5b0SEmil Constantinescu } 98961692a83SJed Brown 9909566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, Zstage)); 991de19f811SJed Brown for (j = 0; j < i; j++) w[j] = At[i * s + j]; 9929566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Zstage, i, w, Y)); 99361692a83SJed Brown 99461692a83SJed Brown for (j = 0; j < i; j++) w[j] = 1. / h * GammaInv[i * s + j]; 9959566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Zdot)); 9969566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Zdot, i, w, Y)); 99761692a83SJed Brown 998e27a552bSJed Brown /* Initial guess taken from last stage */ 9999566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Y[i])); 100061692a83SJed Brown 10017d4bf2deSEmil Constantinescu if (!ros->stage_explicit) { 10029566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 100361692a83SJed Brown if (!ros->recompute_jacobian && !i) { 10049566063dSJacob Faibussowitsch PetscCall(SNESGetLagJacobian(snes, &lag)); 10056aad120cSJose E. Roman if (lag == 1) { /* use did not set a nontrivial lag, so lag over all stages */ 10069566063dSJacob Faibussowitsch PetscCall(SNESSetLagJacobian(snes, -2)); /* Recompute the Jacobian on this solve, but not again for the rest of the stages */ 1007f7f07198SBarry Smith } 100861692a83SJed Brown } 10099566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, Y[i])); 10109371c9d4SSatish 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 */ } 10119566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 10129566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(snes, &lits)); 10139371c9d4SSatish Balay ts->snes_its += its; 10149371c9d4SSatish Balay ts->ksp_its += lits; 10157d4bf2deSEmil Constantinescu } else { 10161ce71dffSSatish Balay Mat J, Jp; 10179566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Ydot)); /* Evaluate Y[i]=G(t,Ydot=0,Zstage) */ 10189566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ros->stage_time, Zstage, Ydot, Y[i], PETSC_FALSE)); 10199566063dSJacob Faibussowitsch PetscCall(VecScale(Y[i], -1.0)); 10209566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y[i], -1.0, Zdot)); /*Y[i] = F(Zstage)-Zdot[=GammaInv*Y]*/ 10210feba352SEmil Constantinescu 10229566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Zstage)); /* Zstage = GammaExplicitCorr[i,j] * Y[j] */ 10230feba352SEmil Constantinescu for (j = 0; j < i; j++) w[j] = GammaExplicitCorr[i * s + j]; 10249566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Zstage, i, w, Y)); 1025fecfb714SLisandro Dalcin 1026fecfb714SLisandro Dalcin /* Y[i] = Y[i] + Jac*Zstage[=Jac*GammaExplicitCorr[i,j] * Y[j]] */ 10279566063dSJacob Faibussowitsch PetscCall(TSGetIJacobian(ts, &J, &Jp, NULL, NULL)); 10289566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, ros->stage_time, ts->vec_sol, Ydot, 0, J, Jp, PETSC_FALSE)); 10299566063dSJacob Faibussowitsch PetscCall(MatMult(J, Zstage, Zdot)); 10309566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y[i], -1.0, Zdot)); 10315ef26d82SJed Brown ts->ksp_its += 1; 1032fecfb714SLisandro Dalcin 10339566063dSJacob Faibussowitsch PetscCall(VecScale(Y[i], h)); 10347d4bf2deSEmil Constantinescu } 10359566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, ros->stage_time, i, Y)); 10369566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 10379566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(adapt, ts, ros->stage_time, Y[i], &stageok)); 1038fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 1039e27a552bSJed Brown } 1040e27a552bSJed Brown 1041b39943a6SLisandro Dalcin ros->status = TS_STEP_INCOMPLETE; 10429566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep_RosW(ts, tab->order, ts->vec_sol, NULL)); 1043b39943a6SLisandro Dalcin ros->status = TS_STEP_PENDING; 10449566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 10459566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(adapt)); 10469566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidateAdd(adapt, tab->name, tab->order, 1, tab->ccfl, (PetscReal)tab->s, PETSC_TRUE)); 10479566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 1048b39943a6SLisandro Dalcin ros->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 1049b39943a6SLisandro Dalcin if (!accept) { /* Roll back the current step */ 10509566063dSJacob Faibussowitsch PetscCall(TSRollBack_RosW(ts)); 1051be5899b3SLisandro Dalcin ts->time_step = next_time_step; 1052be5899b3SLisandro Dalcin goto reject_step; 1053b39943a6SLisandro Dalcin } 1054b39943a6SLisandro Dalcin 1055e27a552bSJed Brown ts->ptime += ts->time_step; 1056cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 10571c3436cfSJed Brown break; 1058b39943a6SLisandro Dalcin 1059b39943a6SLisandro Dalcin reject_step: 10609371c9d4SSatish Balay ts->reject++; 10619371c9d4SSatish Balay accept = PETSC_FALSE; 1062be5899b3SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 1063b39943a6SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 106463a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections)); 10651c3436cfSJed Brown } 10661c3436cfSJed Brown } 1067e27a552bSJed Brown PetscFunctionReturn(0); 1068e27a552bSJed Brown } 1069e27a552bSJed Brown 10709371c9d4SSatish Balay static PetscErrorCode TSInterpolate_RosW(TS ts, PetscReal itime, Vec U) { 107161692a83SJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 1072f4aed992SEmil Constantinescu PetscInt s = ros->tableau->s, pinterp = ros->tableau->pinterp, i, j; 1073f4aed992SEmil Constantinescu PetscReal h; 1074f4aed992SEmil Constantinescu PetscReal tt, t; 1075f4aed992SEmil Constantinescu PetscScalar *bt; 1076f4aed992SEmil Constantinescu const PetscReal *Bt = ros->tableau->binterpt; 1077f4aed992SEmil Constantinescu const PetscReal *GammaInv = ros->tableau->GammaInv; 1078f4aed992SEmil Constantinescu PetscScalar *w = ros->work; 1079f4aed992SEmil Constantinescu Vec *Y = ros->Y; 1080e27a552bSJed Brown 1081e27a552bSJed Brown PetscFunctionBegin; 10823c633725SBarry Smith PetscCheck(Bt, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSRosW %s does not have an interpolation formula", ros->tableau->name); 1083f4aed992SEmil Constantinescu 1084f4aed992SEmil Constantinescu switch (ros->status) { 1085f4aed992SEmil Constantinescu case TS_STEP_INCOMPLETE: 1086f4aed992SEmil Constantinescu case TS_STEP_PENDING: 1087f4aed992SEmil Constantinescu h = ts->time_step; 1088f4aed992SEmil Constantinescu t = (itime - ts->ptime) / h; 1089f4aed992SEmil Constantinescu break; 1090f4aed992SEmil Constantinescu case TS_STEP_COMPLETE: 1091be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; 1092f4aed992SEmil Constantinescu t = (itime - ts->ptime) / h + 1; /* In the interval [0,1] */ 1093f4aed992SEmil Constantinescu break; 1094ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus"); 1095f4aed992SEmil Constantinescu } 10969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s, &bt)); 1097f4aed992SEmil Constantinescu for (i = 0; i < s; i++) bt[i] = 0; 1098f4aed992SEmil Constantinescu for (j = 0, tt = t; j < pinterp; j++, tt *= t) { 10999371c9d4SSatish Balay for (i = 0; i < s; i++) { bt[i] += Bt[i * pinterp + j] * tt; } 1100f4aed992SEmil Constantinescu } 1101f4aed992SEmil Constantinescu 1102f4aed992SEmil Constantinescu /* y(t+tt*h) = y(t) + Sum bt(tt) * GammaInv * Ydot */ 1103f9c1d6abSBarry Smith /* U <- 0*/ 11049566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(U)); 1105f9c1d6abSBarry Smith /* U <- Sum bt_i * GammaInv(i,1:i) * Y(1:i) */ 11063ca35412SEmil Constantinescu for (j = 0; j < s; j++) w[j] = 0; 11073ca35412SEmil Constantinescu for (j = 0; j < s; j++) { 11089371c9d4SSatish Balay for (i = j; i < s; i++) { w[j] += bt[i] * GammaInv[i * s + j]; } 11093ca35412SEmil Constantinescu } 11109566063dSJacob Faibussowitsch PetscCall(VecMAXPY(U, i, w, Y)); 1111be5899b3SLisandro Dalcin /* U <- y(t) + U */ 11129566063dSJacob Faibussowitsch PetscCall(VecAXPY(U, 1, ros->vec_sol_prev)); 1113f4aed992SEmil Constantinescu 11149566063dSJacob Faibussowitsch PetscCall(PetscFree(bt)); 1115e27a552bSJed Brown PetscFunctionReturn(0); 1116e27a552bSJed Brown } 1117e27a552bSJed Brown 1118e27a552bSJed Brown /*------------------------------------------------------------*/ 1119b39943a6SLisandro Dalcin 11209371c9d4SSatish Balay static PetscErrorCode TSRosWTableauReset(TS ts) { 1121b39943a6SLisandro Dalcin TS_RosW *ros = (TS_RosW *)ts->data; 1122b39943a6SLisandro Dalcin RosWTableau tab = ros->tableau; 1123b39943a6SLisandro Dalcin 1124b39943a6SLisandro Dalcin PetscFunctionBegin; 1125b39943a6SLisandro Dalcin if (!tab) PetscFunctionReturn(0); 11269566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &ros->Y)); 11279566063dSJacob Faibussowitsch PetscCall(PetscFree(ros->work)); 1128b39943a6SLisandro Dalcin PetscFunctionReturn(0); 1129b39943a6SLisandro Dalcin } 1130b39943a6SLisandro Dalcin 11319371c9d4SSatish Balay static PetscErrorCode TSReset_RosW(TS ts) { 113261692a83SJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 1133e27a552bSJed Brown 1134e27a552bSJed Brown PetscFunctionBegin; 11359566063dSJacob Faibussowitsch PetscCall(TSRosWTableauReset(ts)); 11369566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ros->Ydot)); 11379566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ros->Ystage)); 11389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ros->Zdot)); 11399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ros->Zstage)); 11409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ros->vec_sol_prev)); 1141e27a552bSJed Brown PetscFunctionReturn(0); 1142e27a552bSJed Brown } 1143e27a552bSJed Brown 11449371c9d4SSatish Balay static PetscErrorCode TSRosWGetVecs(TS ts, DM dm, Vec *Ydot, Vec *Zdot, Vec *Ystage, Vec *Zstage) { 1145d5e6173cSPeter Brune TS_RosW *rw = (TS_RosW *)ts->data; 1146d5e6173cSPeter Brune 1147d5e6173cSPeter Brune PetscFunctionBegin; 1148d5e6173cSPeter Brune if (Ydot) { 1149d5e6173cSPeter Brune if (dm && dm != ts->dm) { 11509566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Ydot", Ydot)); 1151d5e6173cSPeter Brune } else *Ydot = rw->Ydot; 1152d5e6173cSPeter Brune } 1153d5e6173cSPeter Brune if (Zdot) { 1154d5e6173cSPeter Brune if (dm && dm != ts->dm) { 11559566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Zdot", Zdot)); 1156d5e6173cSPeter Brune } else *Zdot = rw->Zdot; 1157d5e6173cSPeter Brune } 1158d5e6173cSPeter Brune if (Ystage) { 1159d5e6173cSPeter Brune if (dm && dm != ts->dm) { 11609566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Ystage", Ystage)); 1161d5e6173cSPeter Brune } else *Ystage = rw->Ystage; 1162d5e6173cSPeter Brune } 1163d5e6173cSPeter Brune if (Zstage) { 1164d5e6173cSPeter Brune if (dm && dm != ts->dm) { 11659566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Zstage", Zstage)); 1166d5e6173cSPeter Brune } else *Zstage = rw->Zstage; 1167d5e6173cSPeter Brune } 1168d5e6173cSPeter Brune PetscFunctionReturn(0); 1169d5e6173cSPeter Brune } 1170d5e6173cSPeter Brune 11719371c9d4SSatish Balay static PetscErrorCode TSRosWRestoreVecs(TS ts, DM dm, Vec *Ydot, Vec *Zdot, Vec *Ystage, Vec *Zstage) { 1172d5e6173cSPeter Brune PetscFunctionBegin; 1173d5e6173cSPeter Brune if (Ydot) { 1174*48a46eb9SPierre Jolivet if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Ydot", Ydot)); 1175d5e6173cSPeter Brune } 1176d5e6173cSPeter Brune if (Zdot) { 1177*48a46eb9SPierre Jolivet if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Zdot", Zdot)); 1178d5e6173cSPeter Brune } 1179d5e6173cSPeter Brune if (Ystage) { 1180*48a46eb9SPierre Jolivet if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Ystage", Ystage)); 1181d5e6173cSPeter Brune } 1182d5e6173cSPeter Brune if (Zstage) { 1183*48a46eb9SPierre Jolivet if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Zstage", Zstage)); 1184d5e6173cSPeter Brune } 1185d5e6173cSPeter Brune PetscFunctionReturn(0); 1186d5e6173cSPeter Brune } 1187d5e6173cSPeter Brune 11889371c9d4SSatish Balay static PetscErrorCode DMCoarsenHook_TSRosW(DM fine, DM coarse, void *ctx) { 1189d5e6173cSPeter Brune PetscFunctionBegin; 1190d5e6173cSPeter Brune PetscFunctionReturn(0); 1191d5e6173cSPeter Brune } 1192d5e6173cSPeter Brune 11939371c9d4SSatish Balay static PetscErrorCode DMRestrictHook_TSRosW(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) { 1194d5e6173cSPeter Brune TS ts = (TS)ctx; 1195d5e6173cSPeter Brune Vec Ydot, Zdot, Ystage, Zstage; 1196d5e6173cSPeter Brune Vec Ydotc, Zdotc, Ystagec, Zstagec; 1197d5e6173cSPeter Brune 1198d5e6173cSPeter Brune PetscFunctionBegin; 11999566063dSJacob Faibussowitsch PetscCall(TSRosWGetVecs(ts, fine, &Ydot, &Ystage, &Zdot, &Zstage)); 12009566063dSJacob Faibussowitsch PetscCall(TSRosWGetVecs(ts, coarse, &Ydotc, &Ystagec, &Zdotc, &Zstagec)); 12019566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, Ydot, Ydotc)); 12029566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Ydotc, rscale, Ydotc)); 12039566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, Ystage, Ystagec)); 12049566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Ystagec, rscale, Ystagec)); 12059566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, Zdot, Zdotc)); 12069566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Zdotc, rscale, Zdotc)); 12079566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, Zstage, Zstagec)); 12089566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Zstagec, rscale, Zstagec)); 12099566063dSJacob Faibussowitsch PetscCall(TSRosWRestoreVecs(ts, fine, &Ydot, &Ystage, &Zdot, &Zstage)); 12109566063dSJacob Faibussowitsch PetscCall(TSRosWRestoreVecs(ts, coarse, &Ydotc, &Ystagec, &Zdotc, &Zstagec)); 1211d5e6173cSPeter Brune PetscFunctionReturn(0); 1212d5e6173cSPeter Brune } 1213d5e6173cSPeter Brune 12149371c9d4SSatish Balay static PetscErrorCode DMSubDomainHook_TSRosW(DM fine, DM coarse, void *ctx) { 1215258e1594SPeter Brune PetscFunctionBegin; 1216258e1594SPeter Brune PetscFunctionReturn(0); 1217258e1594SPeter Brune } 1218258e1594SPeter Brune 12199371c9d4SSatish Balay static PetscErrorCode DMSubDomainRestrictHook_TSRosW(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) { 1220258e1594SPeter Brune TS ts = (TS)ctx; 1221258e1594SPeter Brune Vec Ydot, Zdot, Ystage, Zstage; 1222258e1594SPeter Brune Vec Ydots, Zdots, Ystages, Zstages; 1223258e1594SPeter Brune 1224258e1594SPeter Brune PetscFunctionBegin; 12259566063dSJacob Faibussowitsch PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Ystage, &Zdot, &Zstage)); 12269566063dSJacob Faibussowitsch PetscCall(TSRosWGetVecs(ts, subdm, &Ydots, &Ystages, &Zdots, &Zstages)); 1227258e1594SPeter Brune 12289566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, Ydot, Ydots, INSERT_VALUES, SCATTER_FORWARD)); 12299566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, Ydot, Ydots, INSERT_VALUES, SCATTER_FORWARD)); 1230258e1594SPeter Brune 12319566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, Ystage, Ystages, INSERT_VALUES, SCATTER_FORWARD)); 12329566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, Ystage, Ystages, INSERT_VALUES, SCATTER_FORWARD)); 1233258e1594SPeter Brune 12349566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, Zdot, Zdots, INSERT_VALUES, SCATTER_FORWARD)); 12359566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, Zdot, Zdots, INSERT_VALUES, SCATTER_FORWARD)); 1236258e1594SPeter Brune 12379566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, Zstage, Zstages, INSERT_VALUES, SCATTER_FORWARD)); 12389566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, Zstage, Zstages, INSERT_VALUES, SCATTER_FORWARD)); 1239258e1594SPeter Brune 12409566063dSJacob Faibussowitsch PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Ystage, &Zdot, &Zstage)); 12419566063dSJacob Faibussowitsch PetscCall(TSRosWRestoreVecs(ts, subdm, &Ydots, &Ystages, &Zdots, &Zstages)); 1242258e1594SPeter Brune PetscFunctionReturn(0); 1243258e1594SPeter Brune } 1244258e1594SPeter Brune 1245e27a552bSJed Brown /* 1246e27a552bSJed Brown This defines the nonlinear equation that is to be solved with SNES 1247e27a552bSJed Brown G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 1248e27a552bSJed Brown */ 12499371c9d4SSatish Balay static PetscErrorCode SNESTSFormFunction_RosW(SNES snes, Vec U, Vec F, TS ts) { 125061692a83SJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 1251d5e6173cSPeter Brune Vec Ydot, Zdot, Ystage, Zstage; 1252b296d7d5SJed Brown PetscReal shift = ros->scoeff / ts->time_step; 1253d5e6173cSPeter Brune DM dm, dmsave; 1254e27a552bSJed Brown 1255e27a552bSJed Brown PetscFunctionBegin; 12569566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 12579566063dSJacob Faibussowitsch PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage)); 12589566063dSJacob Faibussowitsch PetscCall(VecWAXPY(Ydot, shift, U, Zdot)); /* Ydot = shift*U + Zdot */ 12599566063dSJacob Faibussowitsch PetscCall(VecWAXPY(Ystage, 1.0, U, Zstage)); /* Ystage = U + Zstage */ 1260d5e6173cSPeter Brune dmsave = ts->dm; 1261d5e6173cSPeter Brune ts->dm = dm; 12629566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ros->stage_time, Ystage, Ydot, F, PETSC_FALSE)); 1263d5e6173cSPeter Brune ts->dm = dmsave; 12649566063dSJacob Faibussowitsch PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage)); 1265e27a552bSJed Brown PetscFunctionReturn(0); 1266e27a552bSJed Brown } 1267e27a552bSJed Brown 12689371c9d4SSatish Balay static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes, Vec U, Mat A, Mat B, TS ts) { 126961692a83SJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 1270d5e6173cSPeter Brune Vec Ydot, Zdot, Ystage, Zstage; 1271b296d7d5SJed Brown PetscReal shift = ros->scoeff / ts->time_step; 1272d5e6173cSPeter Brune DM dm, dmsave; 1273e27a552bSJed Brown 1274e27a552bSJed Brown PetscFunctionBegin; 127561692a83SJed Brown /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */ 12769566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 12779566063dSJacob Faibussowitsch PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage)); 1278d5e6173cSPeter Brune dmsave = ts->dm; 1279d5e6173cSPeter Brune ts->dm = dm; 12809566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, ros->stage_time, Ystage, Ydot, shift, A, B, PETSC_TRUE)); 1281d5e6173cSPeter Brune ts->dm = dmsave; 12829566063dSJacob Faibussowitsch PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage)); 1283e27a552bSJed Brown PetscFunctionReturn(0); 1284e27a552bSJed Brown } 1285e27a552bSJed Brown 12869371c9d4SSatish Balay static PetscErrorCode TSRosWTableauSetUp(TS ts) { 1287b39943a6SLisandro Dalcin TS_RosW *ros = (TS_RosW *)ts->data; 1288b39943a6SLisandro Dalcin RosWTableau tab = ros->tableau; 1289b39943a6SLisandro Dalcin 1290b39943a6SLisandro Dalcin PetscFunctionBegin; 12919566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ros->Y)); 12929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &ros->work)); 1293b39943a6SLisandro Dalcin PetscFunctionReturn(0); 1294b39943a6SLisandro Dalcin } 1295b39943a6SLisandro Dalcin 12969371c9d4SSatish Balay static PetscErrorCode TSSetUp_RosW(TS ts) { 129761692a83SJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 1298d5e6173cSPeter Brune DM dm; 1299b39943a6SLisandro Dalcin SNES snes; 1300a3ab5968SHong Zhang TSRHSJacobian rhsjacobian; 1301e27a552bSJed Brown 1302e27a552bSJed Brown PetscFunctionBegin; 13039566063dSJacob Faibussowitsch PetscCall(TSRosWTableauSetUp(ts)); 13049566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &ros->Ydot)); 13059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &ros->Ystage)); 13069566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &ros->Zdot)); 13079566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &ros->Zstage)); 13089566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &ros->vec_sol_prev)); 13099566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 13109566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSRosW, DMRestrictHook_TSRosW, ts)); 13119566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSRosW, DMSubDomainRestrictHook_TSRosW, ts)); 1312b39943a6SLisandro Dalcin /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */ 13139566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 1314*48a46eb9SPierre Jolivet if (!((PetscObject)snes)->type_name) PetscCall(SNESSetType(snes, SNESKSPONLY)); 13159566063dSJacob Faibussowitsch PetscCall(DMTSGetRHSJacobian(dm, &rhsjacobian, NULL)); 1316a3ab5968SHong Zhang if (rhsjacobian == TSComputeRHSJacobianConstant) { 1317a3ab5968SHong Zhang Mat Amat, Pmat; 1318a3ab5968SHong Zhang 1319a3ab5968SHong Zhang /* Set the SNES matrix to be different from the RHS matrix because there is no way to reconstruct shift*M-J */ 13209566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(snes, &Amat, &Pmat, NULL, NULL)); 1321a3ab5968SHong Zhang if (Amat && Amat == ts->Arhs) { 1322a3ab5968SHong Zhang if (Amat == Pmat) { 13239566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ts->Arhs, MAT_COPY_VALUES, &Amat)); 13249566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, Amat, Amat, NULL, NULL)); 1325a3ab5968SHong Zhang } else { 13269566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ts->Arhs, MAT_COPY_VALUES, &Amat)); 13279566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, Amat, NULL, NULL, NULL)); 1328a3ab5968SHong Zhang if (Pmat && Pmat == ts->Brhs) { 13299566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ts->Brhs, MAT_COPY_VALUES, &Pmat)); 13309566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, NULL, Pmat, NULL, NULL)); 13319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Pmat)); 1332a3ab5968SHong Zhang } 1333a3ab5968SHong Zhang } 13349566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Amat)); 1335a3ab5968SHong Zhang } 1336a3ab5968SHong Zhang } 1337e27a552bSJed Brown PetscFunctionReturn(0); 1338e27a552bSJed Brown } 1339e27a552bSJed Brown /*------------------------------------------------------------*/ 1340e27a552bSJed Brown 13419371c9d4SSatish Balay static PetscErrorCode TSSetFromOptions_RosW(TS ts, PetscOptionItems *PetscOptionsObject) { 134261692a83SJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 1343b39943a6SLisandro Dalcin SNES snes; 1344e27a552bSJed Brown 1345e27a552bSJed Brown PetscFunctionBegin; 1346d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "RosW ODE solver options"); 1347e27a552bSJed Brown { 134861692a83SJed Brown RosWTableauLink link; 1349e27a552bSJed Brown PetscInt count, choice; 1350e27a552bSJed Brown PetscBool flg; 1351e27a552bSJed Brown const char **namelist; 135261692a83SJed Brown 13539371c9d4SSatish Balay for (link = RosWTableauList, count = 0; link; link = link->next, count++) 13549371c9d4SSatish Balay ; 13559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count, (char ***)&namelist)); 135661692a83SJed Brown for (link = RosWTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name; 13579566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_rosw_type", "Family of Rosenbrock-W method", "TSRosWSetType", (const char *const *)namelist, count, ros->tableau->name, &choice, &flg)); 13589566063dSJacob Faibussowitsch if (flg) PetscCall(TSRosWSetType(ts, namelist[choice])); 13599566063dSJacob Faibussowitsch PetscCall(PetscFree(namelist)); 136061692a83SJed Brown 13619566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_rosw_recompute_jacobian", "Recompute the Jacobian at each stage", "TSRosWSetRecomputeJacobian", ros->recompute_jacobian, &ros->recompute_jacobian, NULL)); 1362b39943a6SLisandro Dalcin } 1363d0609cedSBarry Smith PetscOptionsHeadEnd(); 136461692a83SJed Brown /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */ 13659566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 1366*48a46eb9SPierre Jolivet if (!((PetscObject)snes)->type_name) PetscCall(SNESSetType(snes, SNESKSPONLY)); 1367e27a552bSJed Brown PetscFunctionReturn(0); 1368e27a552bSJed Brown } 1369e27a552bSJed Brown 13709371c9d4SSatish Balay static PetscErrorCode TSView_RosW(TS ts, PetscViewer viewer) { 137161692a83SJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 1372e27a552bSJed Brown PetscBool iascii; 1373e27a552bSJed Brown 1374e27a552bSJed Brown PetscFunctionBegin; 13759566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 1376e27a552bSJed Brown if (iascii) { 13779c334d8fSLisandro Dalcin RosWTableau tab = ros->tableau; 137819fd82e9SBarry Smith TSRosWType rostype; 13799c334d8fSLisandro Dalcin char buf[512]; 1380e408995aSJed Brown PetscInt i; 1381e408995aSJed Brown PetscReal abscissa[512]; 13829566063dSJacob Faibussowitsch PetscCall(TSRosWGetType(ts, &rostype)); 13839566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Rosenbrock-W %s\n", rostype)); 13849566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, tab->ASum)); 13859566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Abscissa of A = %s\n", buf)); 1386e408995aSJed Brown for (i = 0; i < tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i]; 13879566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, abscissa)); 13889566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Abscissa of A+Gamma = %s\n", buf)); 1389e27a552bSJed Brown } 1390e27a552bSJed Brown PetscFunctionReturn(0); 1391e27a552bSJed Brown } 1392e27a552bSJed Brown 13939371c9d4SSatish Balay static PetscErrorCode TSLoad_RosW(TS ts, PetscViewer viewer) { 13949200755eSBarry Smith SNES snes; 13959c334d8fSLisandro Dalcin TSAdapt adapt; 13969200755eSBarry Smith 13979200755eSBarry Smith PetscFunctionBegin; 13989566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 13999566063dSJacob Faibussowitsch PetscCall(TSAdaptLoad(adapt, viewer)); 14009566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 14019566063dSJacob Faibussowitsch PetscCall(SNESLoad(snes, viewer)); 14029200755eSBarry Smith /* function and Jacobian context for SNES when used with TS is always ts object */ 14039566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, NULL, NULL, ts)); 14049566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, ts)); 14059200755eSBarry Smith PetscFunctionReturn(0); 14069200755eSBarry Smith } 14079200755eSBarry Smith 1408e27a552bSJed Brown /*@C 140961692a83SJed Brown TSRosWSetType - Set the type of Rosenbrock-W scheme 1410e27a552bSJed Brown 1411e27a552bSJed Brown Logically collective 1412e27a552bSJed Brown 1413d8d19677SJose E. Roman Input Parameters: 1414e27a552bSJed Brown + ts - timestepping context 1415b92453a8SLisandro Dalcin - roswtype - type of Rosenbrock-W scheme 1416e27a552bSJed Brown 1417020d8f30SJed Brown Level: beginner 1418e27a552bSJed Brown 1419db781477SPatrick Sanan .seealso: `TSRosWGetType()`, `TSROSW`, `TSROSW2M`, `TSROSW2P`, `TSROSWRA3PW`, `TSROSWRA34PW2`, `TSROSWRODAS3`, `TSROSWSANDU3`, `TSROSWASSP3P3S1C`, `TSROSWLASSP3P4S2C`, `TSROSWLLSSP3P4S2C`, `TSROSWARK3` 1420e27a552bSJed Brown @*/ 14219371c9d4SSatish Balay PetscErrorCode TSRosWSetType(TS ts, TSRosWType roswtype) { 1422e27a552bSJed Brown PetscFunctionBegin; 1423e27a552bSJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1424b92453a8SLisandro Dalcin PetscValidCharPointer(roswtype, 2); 1425cac4c232SBarry Smith PetscTryMethod(ts, "TSRosWSetType_C", (TS, TSRosWType), (ts, roswtype)); 1426e27a552bSJed Brown PetscFunctionReturn(0); 1427e27a552bSJed Brown } 1428e27a552bSJed Brown 1429e27a552bSJed Brown /*@C 143061692a83SJed Brown TSRosWGetType - Get the type of Rosenbrock-W scheme 1431e27a552bSJed Brown 1432e27a552bSJed Brown Logically collective 1433e27a552bSJed Brown 1434e27a552bSJed Brown Input Parameter: 1435e27a552bSJed Brown . ts - timestepping context 1436e27a552bSJed Brown 1437e27a552bSJed Brown Output Parameter: 143861692a83SJed Brown . rostype - type of Rosenbrock-W scheme 1439e27a552bSJed Brown 1440e27a552bSJed Brown Level: intermediate 1441e27a552bSJed Brown 1442db781477SPatrick Sanan .seealso: `TSRosWGetType()` 1443e27a552bSJed Brown @*/ 14449371c9d4SSatish Balay PetscErrorCode TSRosWGetType(TS ts, TSRosWType *rostype) { 1445e27a552bSJed Brown PetscFunctionBegin; 1446e27a552bSJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1447cac4c232SBarry Smith PetscUseMethod(ts, "TSRosWGetType_C", (TS, TSRosWType *), (ts, rostype)); 1448e27a552bSJed Brown PetscFunctionReturn(0); 1449e27a552bSJed Brown } 1450e27a552bSJed Brown 1451e27a552bSJed Brown /*@C 145261692a83SJed Brown TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step. 1453e27a552bSJed Brown 1454e27a552bSJed Brown Logically collective 1455e27a552bSJed Brown 1456d8d19677SJose E. Roman Input Parameters: 1457e27a552bSJed Brown + ts - timestepping context 145861692a83SJed Brown - flg - PETSC_TRUE to recompute the Jacobian at each stage 1459e27a552bSJed Brown 1460e27a552bSJed Brown Level: intermediate 1461e27a552bSJed Brown 1462db781477SPatrick Sanan .seealso: `TSRosWGetType()` 1463e27a552bSJed Brown @*/ 14649371c9d4SSatish Balay PetscErrorCode TSRosWSetRecomputeJacobian(TS ts, PetscBool flg) { 1465e27a552bSJed Brown PetscFunctionBegin; 1466e27a552bSJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1467cac4c232SBarry Smith PetscTryMethod(ts, "TSRosWSetRecomputeJacobian_C", (TS, PetscBool), (ts, flg)); 1468e27a552bSJed Brown PetscFunctionReturn(0); 1469e27a552bSJed Brown } 1470e27a552bSJed Brown 14719371c9d4SSatish Balay static PetscErrorCode TSRosWGetType_RosW(TS ts, TSRosWType *rostype) { 147261692a83SJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 1473e27a552bSJed Brown 1474e27a552bSJed Brown PetscFunctionBegin; 147561692a83SJed Brown *rostype = ros->tableau->name; 1476e27a552bSJed Brown PetscFunctionReturn(0); 1477e27a552bSJed Brown } 1478ef20d060SBarry Smith 14799371c9d4SSatish Balay static PetscErrorCode TSRosWSetType_RosW(TS ts, TSRosWType rostype) { 148061692a83SJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 1481e27a552bSJed Brown PetscBool match; 148261692a83SJed Brown RosWTableauLink link; 1483e27a552bSJed Brown 1484e27a552bSJed Brown PetscFunctionBegin; 148561692a83SJed Brown if (ros->tableau) { 14869566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(ros->tableau->name, rostype, &match)); 1487e27a552bSJed Brown if (match) PetscFunctionReturn(0); 1488e27a552bSJed Brown } 148961692a83SJed Brown for (link = RosWTableauList; link; link = link->next) { 14909566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(link->tab.name, rostype, &match)); 1491e27a552bSJed Brown if (match) { 14929566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSRosWTableauReset(ts)); 149361692a83SJed Brown ros->tableau = &link->tab; 14949566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSRosWTableauSetUp(ts)); 14952ffb9264SLisandro Dalcin ts->default_adapt_type = ros->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 1496e27a552bSJed Brown PetscFunctionReturn(0); 1497e27a552bSJed Brown } 1498e27a552bSJed Brown } 149998921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", rostype); 1500e27a552bSJed Brown } 150161692a83SJed Brown 15029371c9d4SSatish Balay static PetscErrorCode TSRosWSetRecomputeJacobian_RosW(TS ts, PetscBool flg) { 150361692a83SJed Brown TS_RosW *ros = (TS_RosW *)ts->data; 1504e27a552bSJed Brown 1505e27a552bSJed Brown PetscFunctionBegin; 150661692a83SJed Brown ros->recompute_jacobian = flg; 1507e27a552bSJed Brown PetscFunctionReturn(0); 1508e27a552bSJed Brown } 1509e27a552bSJed Brown 15109371c9d4SSatish Balay static PetscErrorCode TSDestroy_RosW(TS ts) { 1511b3a6b972SJed Brown PetscFunctionBegin; 15129566063dSJacob Faibussowitsch PetscCall(TSReset_RosW(ts)); 1513b3a6b972SJed Brown if (ts->dm) { 15149566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSRosW, DMRestrictHook_TSRosW, ts)); 15159566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSRosW, DMSubDomainRestrictHook_TSRosW, ts)); 1516b3a6b972SJed Brown } 15179566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 15189566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWGetType_C", NULL)); 15199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetType_C", NULL)); 15209566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetRecomputeJacobian_C", NULL)); 1521b3a6b972SJed Brown PetscFunctionReturn(0); 1522b3a6b972SJed Brown } 1523d5e6173cSPeter Brune 1524e27a552bSJed Brown /* ------------------------------------------------------------ */ 1525e27a552bSJed Brown /*MC 1526020d8f30SJed Brown TSROSW - ODE solver using Rosenbrock-W schemes 1527e27a552bSJed Brown 1528e27a552bSJed Brown These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly 1529e27a552bSJed Brown nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part 1530e27a552bSJed Brown of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction(). 1531e27a552bSJed Brown 1532e27a552bSJed Brown Notes: 153361692a83SJed Brown This method currently only works with autonomous ODE and DAE. 153461692a83SJed Brown 1535d0685a90SJed Brown Consider trying TSARKIMEX if the stiff part is strongly nonlinear. 1536d0685a90SJed Brown 15373d5a8a6aSBarry 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 15383d5a8a6aSBarry Smith 1539e94cfbe0SPatrick Sanan Developer Notes: 154061692a83SJed Brown Rosenbrock-W methods are typically specified for autonomous ODE 154161692a83SJed Brown 1542f9c1d6abSBarry Smith $ udot = f(u) 154361692a83SJed Brown 154461692a83SJed Brown by the stage equations 154561692a83SJed Brown 1546f9c1d6abSBarry Smith $ k_i = h f(u_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j 154761692a83SJed Brown 154861692a83SJed Brown and step completion formula 154961692a83SJed Brown 1550f9c1d6abSBarry Smith $ u_1 = u_0 + sum_j b_j k_j 155161692a83SJed Brown 1552f9c1d6abSBarry Smith with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(u) 155361692a83SJed Brown and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner, 155461692a83SJed Brown we define new variables for the stage equations 155561692a83SJed Brown 155661692a83SJed Brown $ y_i = gamma_ij k_j 155761692a83SJed Brown 155861692a83SJed Brown The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define 155961692a83SJed Brown 1560b70472e2SJed Brown $ A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-1} 156161692a83SJed Brown 156261692a83SJed Brown to rewrite the method as 156361692a83SJed Brown 1564f9c1d6abSBarry 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 1565f9c1d6abSBarry Smith $ u_1 = u_0 + sum_j bt_j y_j 156661692a83SJed Brown 156761692a83SJed Brown where we have introduced the mass matrix M. Continue by defining 156861692a83SJed Brown 156961692a83SJed Brown $ ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j 157061692a83SJed Brown 157161692a83SJed Brown or, more compactly in tensor notation 157261692a83SJed Brown 157361692a83SJed Brown $ Ydot = 1/h (Gamma^{-1} \otimes I) Y . 157461692a83SJed Brown 157561692a83SJed Brown Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current 157661692a83SJed Brown stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the 157761692a83SJed Brown equation 157861692a83SJed Brown 1579f9c1d6abSBarry Smith $ g(u_0 + sum_j a_ij y_j + y_i, ydot_i) = 0 158061692a83SJed Brown 158161692a83SJed Brown with initial guess y_i = 0. 1582e27a552bSJed Brown 1583e27a552bSJed Brown Level: beginner 1584e27a552bSJed Brown 1585db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSRosWSetType()`, `TSRosWRegister()`, `TSROSWTHETA1`, `TSROSWTHETA2`, `TSROSW2M`, `TSROSW2P`, `TSROSWRA3PW`, `TSROSWRA34PW2`, `TSROSWRODAS3`, 1586db781477SPatrick Sanan `TSROSWSANDU3`, `TSROSWASSP3P3S1C`, `TSROSWLASSP3P4S2C`, `TSROSWLLSSP3P4S2C`, `TSROSWGRK4T`, `TSROSWSHAMP4`, `TSROSWVELDD4`, `TSROSW4L` 1587e27a552bSJed Brown M*/ 15889371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TSCreate_RosW(TS ts) { 158961692a83SJed Brown TS_RosW *ros; 1590e27a552bSJed Brown 1591e27a552bSJed Brown PetscFunctionBegin; 15929566063dSJacob Faibussowitsch PetscCall(TSRosWInitializePackage()); 1593e27a552bSJed Brown 1594e27a552bSJed Brown ts->ops->reset = TSReset_RosW; 1595e27a552bSJed Brown ts->ops->destroy = TSDestroy_RosW; 1596e27a552bSJed Brown ts->ops->view = TSView_RosW; 15979200755eSBarry Smith ts->ops->load = TSLoad_RosW; 1598e27a552bSJed Brown ts->ops->setup = TSSetUp_RosW; 1599e27a552bSJed Brown ts->ops->step = TSStep_RosW; 1600e27a552bSJed Brown ts->ops->interpolate = TSInterpolate_RosW; 16011c3436cfSJed Brown ts->ops->evaluatestep = TSEvaluateStep_RosW; 160224655328SShri ts->ops->rollback = TSRollBack_RosW; 1603e27a552bSJed Brown ts->ops->setfromoptions = TSSetFromOptions_RosW; 1604e27a552bSJed Brown ts->ops->snesfunction = SNESTSFormFunction_RosW; 1605e27a552bSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_RosW; 1606e27a552bSJed Brown 1607efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 1608efd4aadfSBarry Smith 16099566063dSJacob Faibussowitsch PetscCall(PetscNewLog(ts, &ros)); 161061692a83SJed Brown ts->data = (void *)ros; 1611e27a552bSJed Brown 16129566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWGetType_C", TSRosWGetType_RosW)); 16139566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetType_C", TSRosWSetType_RosW)); 16149566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetRecomputeJacobian_C", TSRosWSetRecomputeJacobian_RosW)); 1615b39943a6SLisandro Dalcin 16169566063dSJacob Faibussowitsch PetscCall(TSRosWSetType(ts, TSRosWDefault)); 1617e27a552bSJed Brown PetscFunctionReturn(0); 1618e27a552bSJed Brown } 1619