xref: /petsc/src/ts/impls/rosw/rosw.c (revision b16ce868fadf09c63ae5be1e36b2479a88beb61d)
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 
731cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`
743606a31eSEmil Constantinescu M*/
753606a31eSEmil Constantinescu 
763606a31eSEmil Constantinescu /*MC
773606a31eSEmil Constantinescu      TSROSWTHETA2 - One stage second order A-stable Rosenbrock-W scheme (aka theta method).
783606a31eSEmil Constantinescu 
793606a31eSEmil Constantinescu      Only an approximate Jacobian is needed.
803606a31eSEmil Constantinescu 
813606a31eSEmil Constantinescu      Level: intermediate
823606a31eSEmil Constantinescu 
831cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`
843606a31eSEmil Constantinescu M*/
853606a31eSEmil Constantinescu 
863606a31eSEmil Constantinescu /*MC
87fe7e6d57SJed Brown      TSROSW2M - Two stage second order L-stable Rosenbrock-W scheme.
88fe7e6d57SJed Brown 
89fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2P.
90fe7e6d57SJed Brown 
91fe7e6d57SJed Brown      Level: intermediate
92fe7e6d57SJed Brown 
931cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`
94fe7e6d57SJed Brown M*/
95fe7e6d57SJed Brown 
96fe7e6d57SJed Brown /*MC
97fe7e6d57SJed Brown      TSROSW2P - Two stage second order L-stable Rosenbrock-W scheme.
98fe7e6d57SJed Brown 
99fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2M.
100fe7e6d57SJed Brown 
101fe7e6d57SJed Brown      Level: intermediate
102fe7e6d57SJed Brown 
1031cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`
104fe7e6d57SJed Brown M*/
105fe7e6d57SJed Brown 
106fe7e6d57SJed Brown /*MC
107fe7e6d57SJed Brown      TSROSWRA3PW - Three stage third order Rosenbrock-W scheme for PDAE of index 1.
108fe7e6d57SJed Brown 
109fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.
110fe7e6d57SJed Brown 
111fe7e6d57SJed Brown      This is strongly A-stable with R(infty) = 0.73. The embedded method of order 2 is strongly A-stable with R(infty) = 0.73.
112fe7e6d57SJed Brown 
113bcf0153eSBarry Smith      Level: intermediate
114bcf0153eSBarry Smith 
115fe7e6d57SJed Brown      References:
116606c0280SSatish Balay .  * - Rang and Angermann, New Rosenbrock W methods of order 3 for partial differential algebraic equations of index 1, 2005.
117fe7e6d57SJed Brown 
1181cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`
119fe7e6d57SJed Brown M*/
120fe7e6d57SJed Brown 
121fe7e6d57SJed Brown /*MC
122fe7e6d57SJed Brown      TSROSWRA34PW2 - Four stage third order L-stable Rosenbrock-W scheme for PDAE of index 1.
123fe7e6d57SJed Brown 
124fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.
125fe7e6d57SJed Brown 
126fe7e6d57SJed Brown      This is strongly A-stable with R(infty) = 0. The embedded method of order 2 is strongly A-stable with R(infty) = 0.48.
127fe7e6d57SJed Brown 
128bcf0153eSBarry Smith      Level: intermediate
129bcf0153eSBarry Smith 
130fe7e6d57SJed Brown      References:
131606c0280SSatish Balay .  * - Rang and Angermann, New Rosenbrock W methods of order 3 for partial differential algebraic equations of index 1, 2005.
132fe7e6d57SJed Brown 
1331cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`
134fe7e6d57SJed Brown M*/
135fe7e6d57SJed Brown 
136ef3c5b88SJed Brown /*MC
137ef3c5b88SJed Brown      TSROSWRODAS3 - Four stage third order L-stable Rosenbrock scheme
138ef3c5b88SJed Brown 
139ef3c5b88SJed Brown      By default, the Jacobian is only recomputed once per step.
140ef3c5b88SJed Brown 
141ef3c5b88SJed Brown      Both the third order and embedded second order methods are stiffly accurate and L-stable.
142ef3c5b88SJed Brown 
143bcf0153eSBarry Smith      Level: intermediate
144bcf0153eSBarry Smith 
145ef3c5b88SJed Brown      References:
146606c0280SSatish Balay .  * - Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997.
147ef3c5b88SJed Brown 
1481cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWSANDU3`
149ef3c5b88SJed Brown M*/
150ef3c5b88SJed Brown 
151ef3c5b88SJed Brown /*MC
152ef3c5b88SJed Brown      TSROSWSANDU3 - Three stage third order L-stable Rosenbrock scheme
153ef3c5b88SJed Brown 
154ef3c5b88SJed Brown      By default, the Jacobian is only recomputed once per step.
155ef3c5b88SJed Brown 
156ef3c5b88SJed Brown      The third order method is L-stable, but not stiffly accurate.
157ef3c5b88SJed Brown      The second order embedded method is strongly A-stable with R(infty) = 0.5.
158ef3c5b88SJed Brown      The internal stages are L-stable.
159ef3c5b88SJed Brown      This method is called ROS3 in the paper.
160ef3c5b88SJed Brown 
161bcf0153eSBarry Smith      Level: intermediate
162bcf0153eSBarry Smith 
163ef3c5b88SJed Brown      References:
164606c0280SSatish Balay .  * - Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997.
165ef3c5b88SJed Brown 
1661cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWRODAS3`
167ef3c5b88SJed Brown M*/
168ef3c5b88SJed Brown 
169961f28d0SJed Brown /*MC
170961f28d0SJed Brown      TSROSWASSP3P3S1C - A-stable Rosenbrock-W method with SSP explicit part, third order, three stages
171961f28d0SJed Brown 
172961f28d0SJed Brown      By default, the Jacobian is only recomputed once per step.
173961f28d0SJed Brown 
174961f28d0SJed Brown      A-stable SPP explicit order 3, 3 stages, CFL 1 (eff = 1/3)
175961f28d0SJed Brown 
176bcf0153eSBarry Smith      Level: intermediate
177bcf0153eSBarry Smith 
178961f28d0SJed Brown      References:
179606c0280SSatish Balay . * - Emil Constantinescu
180961f28d0SJed Brown 
1811cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWLASSP3P4S2C`, `TSROSWLLSSP3P4S2C`, `SSP`
182961f28d0SJed Brown M*/
183961f28d0SJed Brown 
184961f28d0SJed Brown /*MC
185998eb97aSJed Brown      TSROSWLASSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, four stages
186961f28d0SJed Brown 
187961f28d0SJed Brown      By default, the Jacobian is only recomputed once per step.
188961f28d0SJed Brown 
189961f28d0SJed Brown      L-stable (A-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2)
190961f28d0SJed Brown 
191bcf0153eSBarry Smith      Level: intermediate
192bcf0153eSBarry Smith 
193961f28d0SJed Brown      References:
194606c0280SSatish Balay . * - Emil Constantinescu
195961f28d0SJed Brown 
1961cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWASSP3P3S1C`, `TSROSWLLSSP3P4S2C`, `TSSSP`
197961f28d0SJed Brown M*/
198961f28d0SJed Brown 
199961f28d0SJed Brown /*MC
200998eb97aSJed Brown      TSROSWLLSSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, four stages
201961f28d0SJed Brown 
202961f28d0SJed Brown      By default, the Jacobian is only recomputed once per step.
203961f28d0SJed Brown 
204961f28d0SJed Brown      L-stable (L-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2)
205961f28d0SJed Brown 
206bcf0153eSBarry Smith      Level: intermediate
207bcf0153eSBarry Smith 
208961f28d0SJed Brown      References:
209606c0280SSatish Balay . * - Emil Constantinescu
210961f28d0SJed Brown 
2111cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWASSP3P3S1C`, `TSROSWLASSP3P4S2C`, `TSSSP`
212961f28d0SJed Brown M*/
213961f28d0SJed Brown 
21442faf41dSJed Brown /*MC
21542faf41dSJed Brown      TSROSWGRK4T - four stage, fourth order Rosenbrock (not W) method from Kaps and Rentrop
21642faf41dSJed Brown 
21742faf41dSJed Brown      By default, the Jacobian is only recomputed once per step.
21842faf41dSJed Brown 
21942faf41dSJed Brown      A(89.3 degrees)-stable, |R(infty)| = 0.454.
22042faf41dSJed Brown 
22142faf41dSJed Brown      This method does not provide a dense output formula.
22242faf41dSJed Brown 
223bcf0153eSBarry Smith      Level: intermediate
224bcf0153eSBarry Smith 
22542faf41dSJed Brown      References:
226606c0280SSatish Balay +   * -  Kaps and Rentrop, Generalized Runge Kutta methods of order four with stepsize control for stiff ordinary differential equations, 1979.
227606c0280SSatish Balay -   * -  Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.
22842faf41dSJed Brown 
22942faf41dSJed Brown      Hairer's code ros4.f
23042faf41dSJed Brown 
2311cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWSHAMP4`, `TSROSWVELDD4`, `TSROSW4L`
23242faf41dSJed Brown M*/
23342faf41dSJed Brown 
23442faf41dSJed Brown /*MC
23542faf41dSJed Brown      TSROSWSHAMP4 - four stage, fourth order Rosenbrock (not W) method from Shampine
23642faf41dSJed Brown 
23742faf41dSJed Brown      By default, the Jacobian is only recomputed once per step.
23842faf41dSJed Brown 
23942faf41dSJed Brown      A-stable, |R(infty)| = 1/3.
24042faf41dSJed Brown 
24142faf41dSJed Brown      This method does not provide a dense output formula.
24242faf41dSJed Brown 
243bcf0153eSBarry Smith      Level: intermediate
244bcf0153eSBarry Smith 
24542faf41dSJed Brown      References:
246606c0280SSatish Balay +   * -  Shampine, Implementation of Rosenbrock methods, 1982.
247606c0280SSatish Balay -   * -  Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.
24842faf41dSJed Brown 
24942faf41dSJed Brown      Hairer's code ros4.f
25042faf41dSJed Brown 
2511cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWGRK4T`, `TSROSWVELDD4`, `TSROSW4L`
25242faf41dSJed Brown M*/
25342faf41dSJed Brown 
25442faf41dSJed Brown /*MC
25542faf41dSJed Brown      TSROSWVELDD4 - four stage, fourth order Rosenbrock (not W) method from van Veldhuizen
25642faf41dSJed Brown 
25742faf41dSJed Brown      By default, the Jacobian is only recomputed once per step.
25842faf41dSJed Brown 
25942faf41dSJed Brown      A(89.5 degrees)-stable, |R(infty)| = 0.24.
26042faf41dSJed Brown 
26142faf41dSJed Brown      This method does not provide a dense output formula.
26242faf41dSJed Brown 
263bcf0153eSBarry Smith      Level: intermediate
264bcf0153eSBarry Smith 
26542faf41dSJed Brown      References:
266606c0280SSatish Balay +   * -  van Veldhuizen, D stability and Kaps Rentrop methods, 1984.
267606c0280SSatish Balay -   * -  Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.
26842faf41dSJed Brown 
26942faf41dSJed Brown      Hairer's code ros4.f
27042faf41dSJed Brown 
2711cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWGRK4T`, `TSROSWSHAMP4`, `TSROSW4L`
27242faf41dSJed Brown M*/
27342faf41dSJed Brown 
27442faf41dSJed Brown /*MC
27542faf41dSJed Brown      TSROSW4L - four stage, fourth order Rosenbrock (not W) method
27642faf41dSJed Brown 
27742faf41dSJed Brown      By default, the Jacobian is only recomputed once per step.
27842faf41dSJed Brown 
27942faf41dSJed Brown      A-stable and L-stable
28042faf41dSJed Brown 
28142faf41dSJed Brown      This method does not provide a dense output formula.
28242faf41dSJed Brown 
283bcf0153eSBarry Smith      Level: intermediate
284bcf0153eSBarry Smith 
28542faf41dSJed Brown      References:
286606c0280SSatish Balay .  * -   Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.
28742faf41dSJed Brown 
28842faf41dSJed Brown      Hairer's code ros4.f
28942faf41dSJed Brown 
2901cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWGRK4T`, `TSROSWSHAMP4`, `TSROSW4L`
29142faf41dSJed Brown M*/
29242faf41dSJed Brown 
293e27a552bSJed Brown /*@C
294bcf0153eSBarry Smith   TSRosWRegisterAll - Registers all of the Rosenbrock-W methods in `TSROSW`
295e27a552bSJed Brown 
296e27a552bSJed Brown   Not Collective, but should be called by all processes which will need the schemes to be registered
297e27a552bSJed Brown 
298e27a552bSJed Brown   Level: advanced
299e27a552bSJed Brown 
3001cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSRosWRegisterDestroy()`
301e27a552bSJed Brown @*/
302d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWRegisterAll(void)
303d71ae5a4SJacob Faibussowitsch {
304e27a552bSJed Brown   PetscFunctionBegin;
3053ba16761SJacob Faibussowitsch   if (TSRosWRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
306e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_TRUE;
307e27a552bSJed Brown 
308e27a552bSJed Brown   {
309bbd56ea5SKarl Rupp     const PetscReal A        = 0;
310bbd56ea5SKarl Rupp     const PetscReal Gamma    = 1;
311bbd56ea5SKarl Rupp     const PetscReal b        = 1;
312bbd56ea5SKarl Rupp     const PetscReal binterpt = 1;
3131f80e275SEmil Constantinescu 
3149566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWTHETA1, 1, 1, &A, &Gamma, &b, NULL, 1, &binterpt));
3153606a31eSEmil Constantinescu   }
3163606a31eSEmil Constantinescu 
3173606a31eSEmil Constantinescu   {
318bbd56ea5SKarl Rupp     const PetscReal A        = 0;
319bbd56ea5SKarl Rupp     const PetscReal Gamma    = 0.5;
320bbd56ea5SKarl Rupp     const PetscReal b        = 1;
321bbd56ea5SKarl Rupp     const PetscReal binterpt = 1;
322bbd56ea5SKarl Rupp 
3239566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWTHETA2, 2, 1, &A, &Gamma, &b, NULL, 1, &binterpt));
3243606a31eSEmil Constantinescu   }
3253606a31eSEmil Constantinescu 
3263606a31eSEmil Constantinescu   {
327da80777bSKarl Rupp     /*const PetscReal g = 1. + 1./PetscSqrtReal(2.0);   Direct evaluation: 1.707106781186547524401. Used for setting up arrays of values known at compile time below. */
328*b16ce868SStefano Zampini     const PetscReal A[2][2] = {
3299371c9d4SSatish Balay       {0,  0},
3309371c9d4SSatish Balay       {1., 0}
331*b16ce868SStefano Zampini     };
332*b16ce868SStefano Zampini     const PetscReal Gamma[2][2] = {
333*b16ce868SStefano Zampini       {1.707106781186547524401,       0                      },
334*b16ce868SStefano Zampini       {-2. * 1.707106781186547524401, 1.707106781186547524401}
335*b16ce868SStefano Zampini     };
336*b16ce868SStefano Zampini     const PetscReal b[2]  = {0.5, 0.5};
337*b16ce868SStefano Zampini     const PetscReal b1[2] = {1.0, 0.0};
3381f80e275SEmil Constantinescu     PetscReal       binterpt[2][2];
339da80777bSKarl Rupp     binterpt[0][0] = 1.707106781186547524401 - 1.0;
340da80777bSKarl Rupp     binterpt[1][0] = 2.0 - 1.707106781186547524401;
341da80777bSKarl Rupp     binterpt[0][1] = 1.707106781186547524401 - 1.5;
342da80777bSKarl Rupp     binterpt[1][1] = 1.5 - 1.707106781186547524401;
343bbd56ea5SKarl Rupp 
3449566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSW2P, 2, 2, &A[0][0], &Gamma[0][0], b, b1, 2, &binterpt[0][0]));
345e27a552bSJed Brown   }
346e27a552bSJed Brown   {
347da80777bSKarl 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. */
348*b16ce868SStefano Zampini     const PetscReal A[2][2] = {
3499371c9d4SSatish Balay       {0,  0},
3509371c9d4SSatish Balay       {1., 0}
351*b16ce868SStefano Zampini     };
352*b16ce868SStefano Zampini     const PetscReal Gamma[2][2] = {
353*b16ce868SStefano Zampini       {0.2928932188134524755992,       0                       },
354*b16ce868SStefano Zampini       {-2. * 0.2928932188134524755992, 0.2928932188134524755992}
355*b16ce868SStefano Zampini     };
356*b16ce868SStefano Zampini     const PetscReal b[2]  = {0.5, 0.5};
357*b16ce868SStefano Zampini     const PetscReal b1[2] = {1.0, 0.0};
3581f80e275SEmil Constantinescu     PetscReal       binterpt[2][2];
359da80777bSKarl Rupp     binterpt[0][0] = 0.2928932188134524755992 - 1.0;
360da80777bSKarl Rupp     binterpt[1][0] = 2.0 - 0.2928932188134524755992;
361da80777bSKarl Rupp     binterpt[0][1] = 0.2928932188134524755992 - 1.5;
362da80777bSKarl Rupp     binterpt[1][1] = 1.5 - 0.2928932188134524755992;
363bbd56ea5SKarl Rupp 
3649566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSW2M, 2, 2, &A[0][0], &Gamma[0][0], b, b1, 2, &binterpt[0][0]));
365fe7e6d57SJed Brown   }
366fe7e6d57SJed Brown   {
367da80777bSKarl Rupp     /*const PetscReal g = 7.8867513459481287e-01; Directly written in-place below */
3681f80e275SEmil Constantinescu     PetscReal       binterpt[3][2];
369*b16ce868SStefano Zampini     const PetscReal A[3][3] = {
3709371c9d4SSatish Balay       {0,                      0, 0},
371fe7e6d57SJed Brown       {1.5773502691896257e+00, 0, 0},
3729371c9d4SSatish Balay       {0.5,                    0, 0}
373*b16ce868SStefano Zampini     };
374*b16ce868SStefano Zampini     const PetscReal Gamma[3][3] = {
375*b16ce868SStefano Zampini       {7.8867513459481287e-01,  0,                       0                     },
376*b16ce868SStefano Zampini       {-1.5773502691896257e+00, 7.8867513459481287e-01,  0                     },
377*b16ce868SStefano Zampini       {-6.7075317547305480e-01, -1.7075317547305482e-01, 7.8867513459481287e-01}
378*b16ce868SStefano Zampini     };
379*b16ce868SStefano Zampini     const PetscReal b[3]  = {1.0566243270259355e-01, 4.9038105676657971e-02, 8.4529946162074843e-01};
380*b16ce868SStefano Zampini     const PetscReal b2[3] = {-1.7863279495408180e-01, 1. / 3., 8.4529946162074843e-01};
3811f80e275SEmil Constantinescu 
3821f80e275SEmil Constantinescu     binterpt[0][0] = -0.8094010767585034;
3831f80e275SEmil Constantinescu     binterpt[1][0] = -0.5;
3841f80e275SEmil Constantinescu     binterpt[2][0] = 2.3094010767585034;
3851f80e275SEmil Constantinescu     binterpt[0][1] = 0.9641016151377548;
3861f80e275SEmil Constantinescu     binterpt[1][1] = 0.5;
3871f80e275SEmil Constantinescu     binterpt[2][1] = -1.4641016151377548;
388bbd56ea5SKarl Rupp 
3899566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWRA3PW, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0]));
390fe7e6d57SJed Brown   }
391fe7e6d57SJed Brown   {
3923ca35412SEmil Constantinescu     PetscReal binterpt[4][3];
393da80777bSKarl Rupp     /*const PetscReal g = 4.3586652150845900e-01; Directly written in-place below */
394*b16ce868SStefano Zampini     const PetscReal A[4][4] = {
3959371c9d4SSatish Balay       {0,                      0,                       0,  0},
396fe7e6d57SJed Brown       {8.7173304301691801e-01, 0,                       0,  0},
397fe7e6d57SJed Brown       {8.4457060015369423e-01, -1.1299064236484185e-01, 0,  0},
3989371c9d4SSatish Balay       {0,                      0,                       1., 0}
399*b16ce868SStefano Zampini     };
400*b16ce868SStefano Zampini     const PetscReal Gamma[4][4] = {
401*b16ce868SStefano Zampini       {4.3586652150845900e-01,  0,                       0,                      0                     },
402*b16ce868SStefano Zampini       {-8.7173304301691801e-01, 4.3586652150845900e-01,  0,                      0                     },
403*b16ce868SStefano Zampini       {-9.0338057013044082e-01, 5.4180672388095326e-02,  4.3586652150845900e-01, 0                     },
404*b16ce868SStefano Zampini       {2.4212380706095346e-01,  -1.2232505839045147e+00, 5.4526025533510214e-01, 4.3586652150845900e-01}
405*b16ce868SStefano Zampini     };
406*b16ce868SStefano Zampini     const PetscReal b[4]  = {2.4212380706095346e-01, -1.2232505839045147e+00, 1.5452602553351020e+00, 4.3586652150845900e-01};
407*b16ce868SStefano Zampini     const PetscReal b2[4] = {3.7810903145819369e-01, -9.6042292212423178e-02, 5.0000000000000000e-01, 2.1793326075422950e-01};
4083ca35412SEmil Constantinescu 
4093ca35412SEmil Constantinescu     binterpt[0][0] = 1.0564298455794094;
4103ca35412SEmil Constantinescu     binterpt[1][0] = 2.296429974281067;
4113ca35412SEmil Constantinescu     binterpt[2][0] = -1.307599564525376;
4123ca35412SEmil Constantinescu     binterpt[3][0] = -1.045260255335102;
4133ca35412SEmil Constantinescu     binterpt[0][1] = -1.3864882699759573;
4143ca35412SEmil Constantinescu     binterpt[1][1] = -8.262611700275677;
4153ca35412SEmil Constantinescu     binterpt[2][1] = 7.250979895056055;
4163ca35412SEmil Constantinescu     binterpt[3][1] = 2.398120075195581;
4173ca35412SEmil Constantinescu     binterpt[0][2] = 0.5721822314575016;
4183ca35412SEmil Constantinescu     binterpt[1][2] = 4.742931142090097;
4193ca35412SEmil Constantinescu     binterpt[2][2] = -4.398120075195578;
4203ca35412SEmil Constantinescu     binterpt[3][2] = -0.9169932983520199;
4213ca35412SEmil Constantinescu 
4229566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWRA34PW2, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0]));
423e27a552bSJed Brown   }
424ef3c5b88SJed Brown   {
425da80777bSKarl Rupp     /* const PetscReal g = 0.5;       Directly written in-place below */
426*b16ce868SStefano Zampini     const PetscReal A[4][4] = {
4279371c9d4SSatish Balay       {0,    0,     0,   0},
428ef3c5b88SJed Brown       {0,    0,     0,   0},
429ef3c5b88SJed Brown       {1.,   0,     0,   0},
4309371c9d4SSatish Balay       {0.75, -0.25, 0.5, 0}
431*b16ce868SStefano Zampini     };
432*b16ce868SStefano Zampini     const PetscReal Gamma[4][4] = {
433*b16ce868SStefano Zampini       {0.5,     0,       0,       0  },
434*b16ce868SStefano Zampini       {1.,      0.5,     0,       0  },
435*b16ce868SStefano Zampini       {-0.25,   -0.25,   0.5,     0  },
436*b16ce868SStefano Zampini       {1. / 12, 1. / 12, -2. / 3, 0.5}
437*b16ce868SStefano Zampini     };
438*b16ce868SStefano Zampini     const PetscReal b[4]  = {5. / 6, -1. / 6, -1. / 6, 0.5};
439*b16ce868SStefano Zampini     const PetscReal b2[4] = {0.75, -0.25, 0.5, 0};
440bbd56ea5SKarl Rupp 
4419566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWRODAS3, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 0, NULL));
442ef3c5b88SJed Brown   }
443ef3c5b88SJed Brown   {
444da80777bSKarl Rupp     /*const PetscReal g = 0.43586652150845899941601945119356;       Directly written in-place below */
445*b16ce868SStefano Zampini     const PetscReal A[3][3] = {
4469371c9d4SSatish Balay       {0,                                  0, 0},
447da80777bSKarl Rupp       {0.43586652150845899941601945119356, 0, 0},
4489371c9d4SSatish Balay       {0.43586652150845899941601945119356, 0, 0}
449*b16ce868SStefano Zampini     };
450*b16ce868SStefano Zampini     const PetscReal Gamma[3][3] = {
451*b16ce868SStefano Zampini       {0.43586652150845899941601945119356,  0,                                  0                                 },
452*b16ce868SStefano Zampini       {-0.19294655696029095575009695436041, 0.43586652150845899941601945119356, 0                                 },
453*b16ce868SStefano Zampini       {0,                                   1.74927148125794685173529749738960, 0.43586652150845899941601945119356}
454*b16ce868SStefano Zampini     };
455*b16ce868SStefano Zampini     const PetscReal b[3]  = {-0.75457412385404315829818998646589, 1.94100407061964420292840123379419, -0.18642994676560104463021124732829};
456*b16ce868SStefano Zampini     const PetscReal b2[3] = {-1.53358745784149585370766523913002, 2.81745131148625772213931745457622, -0.28386385364476186843165221544619};
4571f80e275SEmil Constantinescu 
4581f80e275SEmil Constantinescu     PetscReal binterpt[3][2];
4591f80e275SEmil Constantinescu     binterpt[0][0] = 3.793692883777660870425141387941;
4601f80e275SEmil Constantinescu     binterpt[1][0] = -2.918692883777660870425141387941;
4611f80e275SEmil Constantinescu     binterpt[2][0] = 0.125;
4621f80e275SEmil Constantinescu     binterpt[0][1] = -0.725741064379812106687651020584;
4631f80e275SEmil Constantinescu     binterpt[1][1] = 0.559074397713145440020984353917;
4641f80e275SEmil Constantinescu     binterpt[2][1] = 0.16666666666666666666666666666667;
4651f80e275SEmil Constantinescu 
4669566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWSANDU3, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0]));
467ef3c5b88SJed Brown   }
468b1c69cc3SEmil Constantinescu   {
469da80777bSKarl Rupp     /*const PetscReal s3 = PetscSqrtReal(3.),g = (3.0+s3)/6.0;
470da80777bSKarl Rupp      * Direct evaluation: s3 = 1.732050807568877293527;
471da80777bSKarl Rupp      *                     g = 0.7886751345948128822546;
472da80777bSKarl Rupp      * Values are directly inserted below to ensure availability at compile time (compiler warnings otherwise...) */
473*b16ce868SStefano Zampini     const PetscReal A[3][3] = {
4749371c9d4SSatish Balay       {0,    0,    0},
475b1c69cc3SEmil Constantinescu       {1,    0,    0},
4769371c9d4SSatish Balay       {0.25, 0.25, 0}
477*b16ce868SStefano Zampini     };
478*b16ce868SStefano Zampini     const PetscReal Gamma[3][3] = {
479*b16ce868SStefano Zampini       {0,                                       0,                                      0                       },
480*b16ce868SStefano Zampini       {(-3.0 - 1.732050807568877293527) / 6.0,  0.7886751345948128822546,               0                       },
481*b16ce868SStefano Zampini       {(-3.0 - 1.732050807568877293527) / 24.0, (-3.0 - 1.732050807568877293527) / 8.0, 0.7886751345948128822546}
482*b16ce868SStefano Zampini     };
483*b16ce868SStefano Zampini     const PetscReal b[3]  = {1. / 6., 1. / 6., 2. / 3.};
484*b16ce868SStefano Zampini     const PetscReal b2[3] = {1. / 4., 1. / 4., 1. / 2.};
485c0cb691aSEmil Constantinescu     PetscReal       binterpt[3][2];
486da80777bSKarl Rupp 
487c0cb691aSEmil Constantinescu     binterpt[0][0] = 0.089316397477040902157517886164709;
488c0cb691aSEmil Constantinescu     binterpt[1][0] = -0.91068360252295909784248211383529;
489c0cb691aSEmil Constantinescu     binterpt[2][0] = 1.8213672050459181956849642276706;
490c0cb691aSEmil Constantinescu     binterpt[0][1] = 0.077350269189625764509148780501957;
491c0cb691aSEmil Constantinescu     binterpt[1][1] = 1.077350269189625764509148780502;
492c0cb691aSEmil Constantinescu     binterpt[2][1] = -1.1547005383792515290182975610039;
493bbd56ea5SKarl Rupp 
4949566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWASSP3P3S1C, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0]));
495b1c69cc3SEmil Constantinescu   }
496b1c69cc3SEmil Constantinescu 
497b1c69cc3SEmil Constantinescu   {
498*b16ce868SStefano Zampini     const PetscReal A[4][4] = {
4999371c9d4SSatish Balay       {0,       0,       0,       0},
500b1c69cc3SEmil Constantinescu       {1. / 2., 0,       0,       0},
501b1c69cc3SEmil Constantinescu       {1. / 2., 1. / 2., 0,       0},
5029371c9d4SSatish Balay       {1. / 6., 1. / 6., 1. / 6., 0}
503*b16ce868SStefano Zampini     };
504*b16ce868SStefano Zampini     const PetscReal Gamma[4][4] = {
505*b16ce868SStefano Zampini       {1. / 2., 0,        0,       0},
506*b16ce868SStefano Zampini       {0.0,     1. / 4.,  0,       0},
507*b16ce868SStefano Zampini       {-2.,     -2. / 3., 2. / 3., 0},
508*b16ce868SStefano Zampini       {1. / 2., 5. / 36., -2. / 9, 0}
509*b16ce868SStefano Zampini     };
510*b16ce868SStefano Zampini     const PetscReal b[4]  = {1. / 6., 1. / 6., 1. / 6., 1. / 2.};
511*b16ce868SStefano Zampini     const PetscReal b2[4] = {1. / 8., 3. / 4., 1. / 8., 0};
512c0cb691aSEmil Constantinescu     PetscReal       binterpt[4][3];
513da80777bSKarl Rupp 
514c0cb691aSEmil Constantinescu     binterpt[0][0] = 6.25;
515c0cb691aSEmil Constantinescu     binterpt[1][0] = -30.25;
516c0cb691aSEmil Constantinescu     binterpt[2][0] = 1.75;
517c0cb691aSEmil Constantinescu     binterpt[3][0] = 23.25;
518c0cb691aSEmil Constantinescu     binterpt[0][1] = -9.75;
519c0cb691aSEmil Constantinescu     binterpt[1][1] = 58.75;
520c0cb691aSEmil Constantinescu     binterpt[2][1] = -3.25;
521c0cb691aSEmil Constantinescu     binterpt[3][1] = -45.75;
522c0cb691aSEmil Constantinescu     binterpt[0][2] = 3.6666666666666666666666666666667;
523c0cb691aSEmil Constantinescu     binterpt[1][2] = -28.333333333333333333333333333333;
524c0cb691aSEmil Constantinescu     binterpt[2][2] = 1.6666666666666666666666666666667;
525c0cb691aSEmil Constantinescu     binterpt[3][2] = 23.;
526bbd56ea5SKarl Rupp 
5279566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWLASSP3P4S2C, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0]));
528b1c69cc3SEmil Constantinescu   }
529b1c69cc3SEmil Constantinescu 
530b1c69cc3SEmil Constantinescu   {
531*b16ce868SStefano Zampini     const PetscReal A[4][4] = {
5329371c9d4SSatish Balay       {0,       0,       0,       0},
533b1c69cc3SEmil Constantinescu       {1. / 2., 0,       0,       0},
534b1c69cc3SEmil Constantinescu       {1. / 2., 1. / 2., 0,       0},
5359371c9d4SSatish Balay       {1. / 6., 1. / 6., 1. / 6., 0}
536*b16ce868SStefano Zampini     };
537*b16ce868SStefano Zampini     const PetscReal Gamma[4][4] = {
538*b16ce868SStefano Zampini       {1. / 2.,  0,          0,        0},
539*b16ce868SStefano Zampini       {0.0,      3. / 4.,    0,        0},
540*b16ce868SStefano Zampini       {-2. / 3., -23. / 9.,  2. / 9.,  0},
541*b16ce868SStefano Zampini       {1. / 18., 65. / 108., -2. / 27, 0}
542*b16ce868SStefano Zampini     };
543*b16ce868SStefano Zampini     const PetscReal b[4]  = {1. / 6., 1. / 6., 1. / 6., 1. / 2.};
544*b16ce868SStefano Zampini     const PetscReal b2[4] = {3. / 16., 10. / 16., 3. / 16., 0};
545c0cb691aSEmil Constantinescu     PetscReal       binterpt[4][3];
546da80777bSKarl Rupp 
547c0cb691aSEmil Constantinescu     binterpt[0][0] = 1.6911764705882352941176470588235;
548c0cb691aSEmil Constantinescu     binterpt[1][0] = 3.6813725490196078431372549019608;
549c0cb691aSEmil Constantinescu     binterpt[2][0] = 0.23039215686274509803921568627451;
550c0cb691aSEmil Constantinescu     binterpt[3][0] = -4.6029411764705882352941176470588;
551c0cb691aSEmil Constantinescu     binterpt[0][1] = -0.95588235294117647058823529411765;
552c0cb691aSEmil Constantinescu     binterpt[1][1] = -6.2401960784313725490196078431373;
553c0cb691aSEmil Constantinescu     binterpt[2][1] = -0.31862745098039215686274509803922;
554c0cb691aSEmil Constantinescu     binterpt[3][1] = 7.5147058823529411764705882352941;
555c0cb691aSEmil Constantinescu     binterpt[0][2] = -0.56862745098039215686274509803922;
556c0cb691aSEmil Constantinescu     binterpt[1][2] = 2.7254901960784313725490196078431;
557c0cb691aSEmil Constantinescu     binterpt[2][2] = 0.25490196078431372549019607843137;
558c0cb691aSEmil Constantinescu     binterpt[3][2] = -2.4117647058823529411764705882353;
559bbd56ea5SKarl Rupp 
5609566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWLLSSP3P4S2C, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0]));
561b1c69cc3SEmil Constantinescu   }
562753f8adbSEmil Constantinescu 
563753f8adbSEmil Constantinescu   {
564753f8adbSEmil Constantinescu     PetscReal A[4][4], Gamma[4][4], b[4], b2[4];
5653ca35412SEmil Constantinescu     PetscReal binterpt[4][3];
566753f8adbSEmil Constantinescu 
567753f8adbSEmil Constantinescu     Gamma[0][0] = 0.4358665215084589994160194475295062513822671686978816;
5689371c9d4SSatish Balay     Gamma[0][1] = 0;
5699371c9d4SSatish Balay     Gamma[0][2] = 0;
5709371c9d4SSatish Balay     Gamma[0][3] = 0;
571753f8adbSEmil Constantinescu     Gamma[1][0] = -1.997527830934941248426324674704153457289527280554476;
572753f8adbSEmil Constantinescu     Gamma[1][1] = 0.4358665215084589994160194475295062513822671686978816;
5739371c9d4SSatish Balay     Gamma[1][2] = 0;
5749371c9d4SSatish Balay     Gamma[1][3] = 0;
575753f8adbSEmil Constantinescu     Gamma[2][0] = -1.007948511795029620852002345345404191008352770119903;
576753f8adbSEmil Constantinescu     Gamma[2][1] = -0.004648958462629345562774289390054679806993396798458131;
577753f8adbSEmil Constantinescu     Gamma[2][2] = 0.4358665215084589994160194475295062513822671686978816;
57805e8e825SJed Brown     Gamma[2][3] = 0;
579753f8adbSEmil Constantinescu     Gamma[3][0] = -0.6685429734233467180451604600279552604364311322650783;
580753f8adbSEmil Constantinescu     Gamma[3][1] = 0.6056625986449338476089525334450053439525178740492984;
581753f8adbSEmil Constantinescu     Gamma[3][2] = -0.9717899277217721234705114616271378792182450260943198;
582753f8adbSEmil Constantinescu     Gamma[3][3] = 0;
583753f8adbSEmil Constantinescu 
5849371c9d4SSatish Balay     A[0][0] = 0;
5859371c9d4SSatish Balay     A[0][1] = 0;
5869371c9d4SSatish Balay     A[0][2] = 0;
5879371c9d4SSatish Balay     A[0][3] = 0;
588753f8adbSEmil Constantinescu     A[1][0] = 0.8717330430169179988320388950590125027645343373957631;
5899371c9d4SSatish Balay     A[1][1] = 0;
5909371c9d4SSatish Balay     A[1][2] = 0;
5919371c9d4SSatish Balay     A[1][3] = 0;
592753f8adbSEmil Constantinescu     A[2][0] = 0.5275890119763004115618079766722914408876108660811028;
593753f8adbSEmil Constantinescu     A[2][1] = 0.07241098802369958843819203208518599088698057726988732;
5949371c9d4SSatish Balay     A[2][2] = 0;
5959371c9d4SSatish Balay     A[2][3] = 0;
596753f8adbSEmil Constantinescu     A[3][0] = 0.3990960076760701320627260685975778145384666450351314;
597753f8adbSEmil Constantinescu     A[3][1] = -0.4375576546135194437228463747348862825846903771419953;
598753f8adbSEmil Constantinescu     A[3][2] = 1.038461646937449311660120300601880176655352737312713;
59905e8e825SJed Brown     A[3][3] = 0;
600753f8adbSEmil Constantinescu 
601753f8adbSEmil Constantinescu     b[0] = 0.1876410243467238251612921333138006734899663569186926;
602753f8adbSEmil Constantinescu     b[1] = -0.5952974735769549480478230473706443582188442040780541;
603753f8adbSEmil Constantinescu     b[2] = 0.9717899277217721234705114616271378792182450260943198;
604753f8adbSEmil Constantinescu     b[3] = 0.4358665215084589994160194475295062513822671686978816;
605753f8adbSEmil Constantinescu 
606753f8adbSEmil Constantinescu     b2[0] = 0.2147402862233891404862383521089097657790734483804460;
607753f8adbSEmil Constantinescu     b2[1] = -0.4851622638849390928209050538171743017757490232519684;
608753f8adbSEmil Constantinescu     b2[2] = 0.8687250025203875511662123688667549217531982787600080;
609753f8adbSEmil Constantinescu     b2[3] = 0.4016969751411624011684543450940068201770721128357014;
610753f8adbSEmil Constantinescu 
6113ca35412SEmil Constantinescu     binterpt[0][0] = 2.2565812720167954547104627844105;
6123ca35412SEmil Constantinescu     binterpt[1][0] = 1.349166413351089573796243820819;
6133ca35412SEmil Constantinescu     binterpt[2][0] = -2.4695174540533503758652847586647;
6143ca35412SEmil Constantinescu     binterpt[3][0] = -0.13623023131453465264142184656474;
6153ca35412SEmil Constantinescu     binterpt[0][1] = -3.0826699111559187902922463354557;
6163ca35412SEmil Constantinescu     binterpt[1][1] = -2.4689115685996042534544925650515;
6173ca35412SEmil Constantinescu     binterpt[2][1] = 5.7428279814696677152129332773553;
6183ca35412SEmil Constantinescu     binterpt[3][1] = -0.19124650171414467146619437684812;
6193ca35412SEmil Constantinescu     binterpt[0][2] = 1.0137296634858471607430756831148;
6203ca35412SEmil Constantinescu     binterpt[1][2] = 0.52444768167155973161042570784064;
6213ca35412SEmil Constantinescu     binterpt[2][2] = -2.3015205996945452158771370439586;
6223ca35412SEmil Constantinescu     binterpt[3][2] = 0.76334325453713832352363565300308;
623f4aed992SEmil Constantinescu 
6249566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWARK3, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0]));
625753f8adbSEmil Constantinescu   }
6269566063dSJacob Faibussowitsch   PetscCall(TSRosWRegisterRos4(TSROSWGRK4T, 0.231, PETSC_DEFAULT, PETSC_DEFAULT, 0, -0.1282612945269037e+01));
6279566063dSJacob Faibussowitsch   PetscCall(TSRosWRegisterRos4(TSROSWSHAMP4, 0.5, PETSC_DEFAULT, PETSC_DEFAULT, 0, 125. / 108.));
6289566063dSJacob Faibussowitsch   PetscCall(TSRosWRegisterRos4(TSROSWVELDD4, 0.22570811482256823492, PETSC_DEFAULT, PETSC_DEFAULT, 0, -1.355958941201148));
6299566063dSJacob Faibussowitsch   PetscCall(TSRosWRegisterRos4(TSROSW4L, 0.57282, PETSC_DEFAULT, PETSC_DEFAULT, 0, -1.093502252409163));
6303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
631e27a552bSJed Brown }
632e27a552bSJed Brown 
633e27a552bSJed Brown /*@C
634bcf0153eSBarry Smith   TSRosWRegisterDestroy - Frees the list of schemes that were registered by `TSRosWRegister()`.
635e27a552bSJed Brown 
636e27a552bSJed Brown   Not Collective
637e27a552bSJed Brown 
638e27a552bSJed Brown   Level: advanced
639e27a552bSJed Brown 
6401cc06b55SBarry Smith .seealso: [](ch_ts), `TSRosWRegister()`, `TSRosWRegisterAll()`
641e27a552bSJed Brown @*/
642d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWRegisterDestroy(void)
643d71ae5a4SJacob Faibussowitsch {
64461692a83SJed Brown   RosWTableauLink link;
645e27a552bSJed Brown 
646e27a552bSJed Brown   PetscFunctionBegin;
64761692a83SJed Brown   while ((link = RosWTableauList)) {
64861692a83SJed Brown     RosWTableau t   = &link->tab;
64961692a83SJed Brown     RosWTableauList = link->next;
6509566063dSJacob Faibussowitsch     PetscCall(PetscFree5(t->A, t->Gamma, t->b, t->ASum, t->GammaSum));
6519566063dSJacob Faibussowitsch     PetscCall(PetscFree5(t->At, t->bt, t->GammaInv, t->GammaZeroDiag, t->GammaExplicitCorr));
6529566063dSJacob Faibussowitsch     PetscCall(PetscFree2(t->bembed, t->bembedt));
6539566063dSJacob Faibussowitsch     PetscCall(PetscFree(t->binterpt));
6549566063dSJacob Faibussowitsch     PetscCall(PetscFree(t->name));
6559566063dSJacob Faibussowitsch     PetscCall(PetscFree(link));
656e27a552bSJed Brown   }
657e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_FALSE;
6583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
659e27a552bSJed Brown }
660e27a552bSJed Brown 
661e27a552bSJed Brown /*@C
662bcf0153eSBarry Smith   TSRosWInitializePackage - This function initializes everything in the `TSROSW` package. It is called
663bcf0153eSBarry Smith   from `TSInitializePackage()`.
664e27a552bSJed Brown 
665e27a552bSJed Brown   Level: developer
666e27a552bSJed Brown 
6671cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `PetscInitialize()`, `TSRosWFinalizePackage()`
668e27a552bSJed Brown @*/
669d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWInitializePackage(void)
670d71ae5a4SJacob Faibussowitsch {
671e27a552bSJed Brown   PetscFunctionBegin;
6723ba16761SJacob Faibussowitsch   if (TSRosWPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
673e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_TRUE;
6749566063dSJacob Faibussowitsch   PetscCall(TSRosWRegisterAll());
6759566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(TSRosWFinalizePackage));
6763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
677e27a552bSJed Brown }
678e27a552bSJed Brown 
679e27a552bSJed Brown /*@C
680bcf0153eSBarry Smith   TSRosWFinalizePackage - This function destroys everything in the `TSROSW` package. It is
681bcf0153eSBarry Smith   called from `PetscFinalize()`.
682e27a552bSJed Brown 
683e27a552bSJed Brown   Level: developer
684e27a552bSJed Brown 
6851cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `PetscFinalize()`, `TSRosWInitializePackage()`
686e27a552bSJed Brown @*/
687d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWFinalizePackage(void)
688d71ae5a4SJacob Faibussowitsch {
689e27a552bSJed Brown   PetscFunctionBegin;
690e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_FALSE;
6919566063dSJacob Faibussowitsch   PetscCall(TSRosWRegisterDestroy());
6923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
693e27a552bSJed Brown }
694e27a552bSJed Brown 
695e27a552bSJed Brown /*@C
696bcf0153eSBarry Smith   TSRosWRegister - register a `TSROSW`, Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
697e27a552bSJed Brown 
698e27a552bSJed Brown   Not Collective, but the same schemes should be registered on all processes on which they will be used
699e27a552bSJed Brown 
700e27a552bSJed Brown   Input Parameters:
701e27a552bSJed Brown + name     - identifier for method
702e27a552bSJed Brown . order    - approximation order of method
703e27a552bSJed Brown . s        - number of stages, this is the dimension of the matrices below
70461692a83SJed Brown . A        - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular
70561692a83SJed Brown . Gamma    - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal
706fe7e6d57SJed Brown . b        - Step completion table (dimension s)
7070298fd71SBarry Smith . bembed   - Step completion table for a scheme of order one less (dimension s, NULL if no embedded scheme is available)
708f4aed992SEmil Constantinescu . pinterp  - Order of the interpolation scheme, equal to the number of columns of binterpt
70942faf41dSJed Brown - binterpt - Coefficients of the interpolation formula (dimension s*pinterp)
710e27a552bSJed Brown 
711e27a552bSJed Brown   Level: advanced
712e27a552bSJed Brown 
713bcf0153eSBarry Smith   Note:
714bcf0153eSBarry Smith   Several Rosenbrock W methods are provided, this function is only needed to create new methods.
715bcf0153eSBarry Smith 
7161cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`
717e27a552bSJed Brown @*/
718d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWRegister(TSRosWType name, PetscInt order, PetscInt s, const PetscReal A[], const PetscReal Gamma[], const PetscReal b[], const PetscReal bembed[], PetscInt pinterp, const PetscReal binterpt[])
719d71ae5a4SJacob Faibussowitsch {
72061692a83SJed Brown   RosWTableauLink link;
72161692a83SJed Brown   RosWTableau     t;
72261692a83SJed Brown   PetscInt        i, j, k;
72361692a83SJed Brown   PetscScalar    *GammaInv;
724e27a552bSJed Brown 
725e27a552bSJed Brown   PetscFunctionBegin;
7264f572ea9SToby Isaac   PetscAssertPointer(name, 1);
7274f572ea9SToby Isaac   PetscAssertPointer(A, 4);
7284f572ea9SToby Isaac   PetscAssertPointer(Gamma, 5);
7294f572ea9SToby Isaac   PetscAssertPointer(b, 6);
7304f572ea9SToby Isaac   if (bembed) PetscAssertPointer(bembed, 7);
731fe7e6d57SJed Brown 
7329566063dSJacob Faibussowitsch   PetscCall(TSRosWInitializePackage());
7339566063dSJacob Faibussowitsch   PetscCall(PetscNew(&link));
734e27a552bSJed Brown   t = &link->tab;
7359566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &t->name));
736e27a552bSJed Brown   t->order = order;
737e27a552bSJed Brown   t->s     = s;
7389566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(s * s, &t->A, s * s, &t->Gamma, s, &t->b, s, &t->ASum, s, &t->GammaSum));
7399566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(s * s, &t->At, s, &t->bt, s * s, &t->GammaInv, s, &t->GammaZeroDiag, s * s, &t->GammaExplicitCorr));
7409566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->A, A, s * s));
7419566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->Gamma, Gamma, s * s));
7429566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->GammaExplicitCorr, Gamma, s * s));
7439566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->b, b, s));
744fe7e6d57SJed Brown   if (bembed) {
7459566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(s, &t->bembed, s, &t->bembedt));
7469566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(t->bembed, bembed, s));
747fe7e6d57SJed Brown   }
74861692a83SJed Brown   for (i = 0; i < s; i++) {
74961692a83SJed Brown     t->ASum[i]     = 0;
75061692a83SJed Brown     t->GammaSum[i] = 0;
75161692a83SJed Brown     for (j = 0; j < s; j++) {
75261692a83SJed Brown       t->ASum[i] += A[i * s + j];
753fe7e6d57SJed Brown       t->GammaSum[i] += Gamma[i * s + j];
75461692a83SJed Brown     }
75561692a83SJed Brown   }
7569566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s * s, &GammaInv)); /* Need to use Scalar for inverse, then convert back to Real */
75761692a83SJed Brown   for (i = 0; i < s * s; i++) GammaInv[i] = Gamma[i];
758fd96d5b0SEmil Constantinescu   for (i = 0; i < s; i++) {
759fd96d5b0SEmil Constantinescu     if (Gamma[i * s + i] == 0.0) {
760fd96d5b0SEmil Constantinescu       GammaInv[i * s + i] = 1.0;
761c17803e7SJed Brown       t->GammaZeroDiag[i] = PETSC_TRUE;
762fd96d5b0SEmil Constantinescu     } else {
763c17803e7SJed Brown       t->GammaZeroDiag[i] = PETSC_FALSE;
764fd96d5b0SEmil Constantinescu     }
765fd96d5b0SEmil Constantinescu   }
766fd96d5b0SEmil Constantinescu 
76761692a83SJed Brown   switch (s) {
768d71ae5a4SJacob Faibussowitsch   case 1:
769d71ae5a4SJacob Faibussowitsch     GammaInv[0] = 1. / GammaInv[0];
770d71ae5a4SJacob Faibussowitsch     break;
771d71ae5a4SJacob Faibussowitsch   case 2:
772d71ae5a4SJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_2(GammaInv, 0, PETSC_FALSE, NULL));
773d71ae5a4SJacob Faibussowitsch     break;
774d71ae5a4SJacob Faibussowitsch   case 3:
775d71ae5a4SJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_3(GammaInv, 0, PETSC_FALSE, NULL));
776d71ae5a4SJacob Faibussowitsch     break;
777d71ae5a4SJacob Faibussowitsch   case 4:
778d71ae5a4SJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_4(GammaInv, 0, PETSC_FALSE, NULL));
779d71ae5a4SJacob Faibussowitsch     break;
78061692a83SJed Brown   case 5: {
78161692a83SJed Brown     PetscInt  ipvt5[5];
78261692a83SJed Brown     MatScalar work5[5 * 5];
7839371c9d4SSatish Balay     PetscCall(PetscKernel_A_gets_inverse_A_5(GammaInv, ipvt5, work5, 0, PETSC_FALSE, NULL));
7849371c9d4SSatish Balay     break;
78561692a83SJed Brown   }
786d71ae5a4SJacob Faibussowitsch   case 6:
787d71ae5a4SJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_6(GammaInv, 0, PETSC_FALSE, NULL));
788d71ae5a4SJacob Faibussowitsch     break;
789d71ae5a4SJacob Faibussowitsch   case 7:
790d71ae5a4SJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_7(GammaInv, 0, PETSC_FALSE, NULL));
791d71ae5a4SJacob Faibussowitsch     break;
792d71ae5a4SJacob Faibussowitsch   default:
793d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not implemented for %" PetscInt_FMT " stages", s);
79461692a83SJed Brown   }
79561692a83SJed Brown   for (i = 0; i < s * s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]);
7969566063dSJacob Faibussowitsch   PetscCall(PetscFree(GammaInv));
79743b21953SEmil Constantinescu 
79843b21953SEmil Constantinescu   for (i = 0; i < s; i++) {
79943b21953SEmil Constantinescu     for (k = 0; k < i + 1; k++) {
80043b21953SEmil Constantinescu       t->GammaExplicitCorr[i * s + k] = (t->GammaExplicitCorr[i * s + k]) * (t->GammaInv[k * s + k]);
801ad540459SPierre Jolivet       for (j = k + 1; j < i + 1; j++) t->GammaExplicitCorr[i * s + k] += (t->GammaExplicitCorr[i * s + j]) * (t->GammaInv[j * s + k]);
80243b21953SEmil Constantinescu     }
80343b21953SEmil Constantinescu   }
80443b21953SEmil Constantinescu 
80561692a83SJed Brown   for (i = 0; i < s; i++) {
80661692a83SJed Brown     for (j = 0; j < s; j++) {
80761692a83SJed Brown       t->At[i * s + j] = 0;
808ad540459SPierre Jolivet       for (k = 0; k < s; k++) t->At[i * s + j] += t->A[i * s + k] * t->GammaInv[k * s + j];
80961692a83SJed Brown     }
81061692a83SJed Brown     t->bt[i] = 0;
811ad540459SPierre Jolivet     for (j = 0; j < s; j++) t->bt[i] += t->b[j] * t->GammaInv[j * s + i];
812fe7e6d57SJed Brown     if (bembed) {
813fe7e6d57SJed Brown       t->bembedt[i] = 0;
814ad540459SPierre Jolivet       for (j = 0; j < s; j++) t->bembedt[i] += t->bembed[j] * t->GammaInv[j * s + i];
815fe7e6d57SJed Brown     }
81661692a83SJed Brown   }
8178d59e960SJed Brown   t->ccfl = 1.0; /* Fix this */
8188d59e960SJed Brown 
819f4aed992SEmil Constantinescu   t->pinterp = pinterp;
8209566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s * pinterp, &t->binterpt));
8219566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->binterpt, binterpt, s * pinterp));
82261692a83SJed Brown   link->next      = RosWTableauList;
82361692a83SJed Brown   RosWTableauList = link;
8243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
825e27a552bSJed Brown }
826e27a552bSJed Brown 
82742faf41dSJed Brown /*@C
828fd292e60Sprj-   TSRosWRegisterRos4 - register a fourth order Rosenbrock scheme by providing parameter choices
82942faf41dSJed Brown 
83042faf41dSJed Brown   Not Collective, but the same schemes should be registered on all processes on which they will be used
83142faf41dSJed Brown 
83242faf41dSJed Brown   Input Parameters:
83342faf41dSJed Brown + name  - identifier for method
83442faf41dSJed Brown . gamma - leading coefficient (diagonal entry)
83542faf41dSJed Brown . a2    - design parameter, see Table 7.2 of Hairer&Wanner
83642faf41dSJed Brown . a3    - design parameter or PETSC_DEFAULT to satisfy one of the order five conditions (Eq 7.22)
83742faf41dSJed Brown . b3    - design parameter, see Table 7.2 of Hairer&Wanner
838a2b725a8SWilliam Gropp - e4    - design parameter for embedded method, see coefficient E4 in ros4.f code from Hairer
83942faf41dSJed Brown 
840bcf0153eSBarry Smith   Level: developer
841bcf0153eSBarry Smith 
84242faf41dSJed Brown   Notes:
84342faf41dSJed Brown   This routine encodes the design of fourth order Rosenbrock methods as described in Hairer and Wanner volume 2.
84442faf41dSJed Brown   It is used here to implement several methods from the book and can be used to experiment with new methods.
84542faf41dSJed Brown   It was written this way instead of by copying coefficients in order to provide better than double precision satisfaction of the order conditions.
84642faf41dSJed Brown 
8471cc06b55SBarry Smith .seealso: [](ch_ts), `TSRosW`, `TSRosWRegister()`
84842faf41dSJed Brown @*/
849d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWRegisterRos4(TSRosWType name, PetscReal gamma, PetscReal a2, PetscReal a3, PetscReal b3, PetscReal e4)
850d71ae5a4SJacob Faibussowitsch {
85142faf41dSJed Brown   /* Declare numeric constants so they can be quad precision without being truncated at double */
8529371c9d4SSatish 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;
85342faf41dSJed Brown   PetscReal   a4, a32, a42, a43, b1, b2, b4, beta2p, beta3p, beta4p, beta32, beta42, beta43, beta32beta2p, beta4jbetajp;
85442faf41dSJed Brown   PetscReal   A[4][4], Gamma[4][4], b[4], bm[4];
85542faf41dSJed Brown   PetscScalar M[3][3], rhs[3];
85642faf41dSJed Brown 
85742faf41dSJed Brown   PetscFunctionBegin;
85842faf41dSJed Brown   /* Step 1: choose Gamma (input) */
85942faf41dSJed Brown   /* Step 2: choose a2,a3,a4; b1,b2,b3,b4 to satisfy order conditions */
86013bcc0bdSJacob Faibussowitsch   if (a3 == (PetscReal)PETSC_DEFAULT) a3 = (one / five - a2 / four) / (one / four - a2 / three); /* Eq 7.22 */
86142faf41dSJed Brown   a4 = a3;                                                                                       /* consequence of 7.20 */
86242faf41dSJed Brown 
86342faf41dSJed Brown   /* Solve order conditions 7.15a, 7.15c, 7.15e */
8649371c9d4SSatish Balay   M[0][0] = one;
8659371c9d4SSatish Balay   M[0][1] = one;
8669371c9d4SSatish Balay   M[0][2] = one; /* 7.15a */
8679371c9d4SSatish Balay   M[1][0] = 0.0;
8689371c9d4SSatish Balay   M[1][1] = a2 * a2;
8699371c9d4SSatish Balay   M[1][2] = a4 * a4; /* 7.15c */
8709371c9d4SSatish Balay   M[2][0] = 0.0;
8719371c9d4SSatish Balay   M[2][1] = a2 * a2 * a2;
8729371c9d4SSatish Balay   M[2][2] = a4 * a4 * a4; /* 7.15e */
87342faf41dSJed Brown   rhs[0]  = one - b3;
87442faf41dSJed Brown   rhs[1]  = one / three - a3 * a3 * b3;
87542faf41dSJed Brown   rhs[2]  = one / four - a3 * a3 * a3 * b3;
8769566063dSJacob Faibussowitsch   PetscCall(PetscKernel_A_gets_inverse_A_3(&M[0][0], 0, PETSC_FALSE, NULL));
87742faf41dSJed Brown   b1 = PetscRealPart(M[0][0] * rhs[0] + M[0][1] * rhs[1] + M[0][2] * rhs[2]);
87842faf41dSJed Brown   b2 = PetscRealPart(M[1][0] * rhs[0] + M[1][1] * rhs[1] + M[1][2] * rhs[2]);
87942faf41dSJed Brown   b4 = PetscRealPart(M[2][0] * rhs[0] + M[2][1] * rhs[1] + M[2][2] * rhs[2]);
88042faf41dSJed Brown 
88142faf41dSJed Brown   /* Step 3 */
88242faf41dSJed Brown   beta43       = (p56 - a2 * p43) / (b4 * a3 * a3 * (a3 - a2)); /* 7.21 */
88342faf41dSJed Brown   beta32beta2p = p44 / (b4 * beta43);                           /* 7.15h */
88442faf41dSJed Brown   beta4jbetajp = (p32 - b3 * beta32beta2p) / b4;
8859371c9d4SSatish Balay   M[0][0]      = b2;
8869371c9d4SSatish Balay   M[0][1]      = b3;
8879371c9d4SSatish Balay   M[0][2]      = b4;
8889371c9d4SSatish Balay   M[1][0]      = a4 * a4 * beta32beta2p - a3 * a3 * beta4jbetajp;
8899371c9d4SSatish Balay   M[1][1]      = a2 * a2 * beta4jbetajp;
8909371c9d4SSatish Balay   M[1][2]      = -a2 * a2 * beta32beta2p;
8919371c9d4SSatish Balay   M[2][0]      = b4 * beta43 * a3 * a3 - p43;
8929371c9d4SSatish Balay   M[2][1]      = -b4 * beta43 * a2 * a2;
8939371c9d4SSatish Balay   M[2][2]      = 0;
8949371c9d4SSatish Balay   rhs[0]       = one / two - gamma;
8959371c9d4SSatish Balay   rhs[1]       = 0;
8969371c9d4SSatish Balay   rhs[2]       = -a2 * a2 * p32;
8979566063dSJacob Faibussowitsch   PetscCall(PetscKernel_A_gets_inverse_A_3(&M[0][0], 0, PETSC_FALSE, NULL));
89842faf41dSJed Brown   beta2p = PetscRealPart(M[0][0] * rhs[0] + M[0][1] * rhs[1] + M[0][2] * rhs[2]);
89942faf41dSJed Brown   beta3p = PetscRealPart(M[1][0] * rhs[0] + M[1][1] * rhs[1] + M[1][2] * rhs[2]);
90042faf41dSJed Brown   beta4p = PetscRealPart(M[2][0] * rhs[0] + M[2][1] * rhs[1] + M[2][2] * rhs[2]);
90142faf41dSJed Brown 
90242faf41dSJed Brown   /* Step 4: back-substitute */
90342faf41dSJed Brown   beta32 = beta32beta2p / beta2p;
90442faf41dSJed Brown   beta42 = (beta4jbetajp - beta43 * beta3p) / beta2p;
90542faf41dSJed Brown 
90642faf41dSJed Brown   /* Step 5: 7.15f and 7.20, then 7.16 */
90742faf41dSJed Brown   a43 = 0;
90842faf41dSJed Brown   a32 = p42 / (b3 * a3 * beta2p + b4 * a4 * beta2p);
90942faf41dSJed Brown   a42 = a32;
91042faf41dSJed Brown 
9119371c9d4SSatish Balay   A[0][0]     = 0;
9129371c9d4SSatish Balay   A[0][1]     = 0;
9139371c9d4SSatish Balay   A[0][2]     = 0;
9149371c9d4SSatish Balay   A[0][3]     = 0;
9159371c9d4SSatish Balay   A[1][0]     = a2;
9169371c9d4SSatish Balay   A[1][1]     = 0;
9179371c9d4SSatish Balay   A[1][2]     = 0;
9189371c9d4SSatish Balay   A[1][3]     = 0;
9199371c9d4SSatish Balay   A[2][0]     = a3 - a32;
9209371c9d4SSatish Balay   A[2][1]     = a32;
9219371c9d4SSatish Balay   A[2][2]     = 0;
9229371c9d4SSatish Balay   A[2][3]     = 0;
9239371c9d4SSatish Balay   A[3][0]     = a4 - a43 - a42;
9249371c9d4SSatish Balay   A[3][1]     = a42;
9259371c9d4SSatish Balay   A[3][2]     = a43;
9269371c9d4SSatish Balay   A[3][3]     = 0;
9279371c9d4SSatish Balay   Gamma[0][0] = gamma;
9289371c9d4SSatish Balay   Gamma[0][1] = 0;
9299371c9d4SSatish Balay   Gamma[0][2] = 0;
9309371c9d4SSatish Balay   Gamma[0][3] = 0;
9319371c9d4SSatish Balay   Gamma[1][0] = beta2p - A[1][0];
9329371c9d4SSatish Balay   Gamma[1][1] = gamma;
9339371c9d4SSatish Balay   Gamma[1][2] = 0;
9349371c9d4SSatish Balay   Gamma[1][3] = 0;
9359371c9d4SSatish Balay   Gamma[2][0] = beta3p - beta32 - A[2][0];
9369371c9d4SSatish Balay   Gamma[2][1] = beta32 - A[2][1];
9379371c9d4SSatish Balay   Gamma[2][2] = gamma;
9389371c9d4SSatish Balay   Gamma[2][3] = 0;
9399371c9d4SSatish Balay   Gamma[3][0] = beta4p - beta42 - beta43 - A[3][0];
9409371c9d4SSatish Balay   Gamma[3][1] = beta42 - A[3][1];
9419371c9d4SSatish Balay   Gamma[3][2] = beta43 - A[3][2];
9429371c9d4SSatish Balay   Gamma[3][3] = gamma;
9439371c9d4SSatish Balay   b[0]        = b1;
9449371c9d4SSatish Balay   b[1]        = b2;
9459371c9d4SSatish Balay   b[2]        = b3;
9469371c9d4SSatish Balay   b[3]        = b4;
94742faf41dSJed Brown 
94842faf41dSJed Brown   /* Construct embedded formula using given e4. We are solving Equation 7.18. */
94942faf41dSJed Brown   bm[3] = b[3] - e4 * gamma;                                              /* using definition of E4 */
95042faf41dSJed Brown   bm[2] = (p32 - beta4jbetajp * bm[3]) / (beta32 * beta2p);               /* fourth row of 7.18 */
95142faf41dSJed Brown   bm[1] = (one / two - gamma - beta3p * bm[2] - beta4p * bm[3]) / beta2p; /* second row */
95242faf41dSJed Brown   bm[0] = one - bm[1] - bm[2] - bm[3];                                    /* first row */
95342faf41dSJed Brown 
95442faf41dSJed Brown   {
95542faf41dSJed Brown     const PetscReal misfit = a2 * a2 * bm[1] + a3 * a3 * bm[2] + a4 * a4 * bm[3] - one / three;
9563c633725SBarry Smith     PetscCheck(PetscAbs(misfit) <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_SUP, "Assumptions violated, could not construct a third order embedded method");
95742faf41dSJed Brown   }
9589566063dSJacob Faibussowitsch   PetscCall(TSRosWRegister(name, 4, 4, &A[0][0], &Gamma[0][0], b, bm, 0, NULL));
9593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
96042faf41dSJed Brown }
96142faf41dSJed Brown 
9621c3436cfSJed Brown /*
9631c3436cfSJed Brown  The step completion formula is
9641c3436cfSJed Brown 
9651c3436cfSJed Brown  x1 = x0 + b^T Y
9661c3436cfSJed Brown 
9671c3436cfSJed Brown  where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been
9681c3436cfSJed Brown  updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write
9691c3436cfSJed Brown 
9701c3436cfSJed Brown  x1e = x0 + be^T Y
9711c3436cfSJed Brown      = x1 - b^T Y + be^T Y
9721c3436cfSJed Brown      = x1 + (be - b)^T Y
9731c3436cfSJed Brown 
9741c3436cfSJed Brown  so we can evaluate the method of different order even after the step has been optimistically completed.
9751c3436cfSJed Brown */
976d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateStep_RosW(TS ts, PetscInt order, Vec U, PetscBool *done)
977d71ae5a4SJacob Faibussowitsch {
9781c3436cfSJed Brown   TS_RosW     *ros = (TS_RosW *)ts->data;
9791c3436cfSJed Brown   RosWTableau  tab = ros->tableau;
9801c3436cfSJed Brown   PetscScalar *w   = ros->work;
9811c3436cfSJed Brown   PetscInt     i;
9821c3436cfSJed Brown 
9831c3436cfSJed Brown   PetscFunctionBegin;
9841c3436cfSJed Brown   if (order == tab->order) {
985108c343cSJed Brown     if (ros->status == TS_STEP_INCOMPLETE) { /* Use standard completion formula */
9869566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol, U));
987de19f811SJed Brown       for (i = 0; i < tab->s; i++) w[i] = tab->bt[i];
9889566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(U, tab->s, w, ros->Y));
9899566063dSJacob Faibussowitsch     } else PetscCall(VecCopy(ts->vec_sol, U));
9901c3436cfSJed Brown     if (done) *done = PETSC_TRUE;
9913ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
9921c3436cfSJed Brown   } else if (order == tab->order - 1) {
9931c3436cfSJed Brown     if (!tab->bembedt) goto unavailable;
994108c343cSJed Brown     if (ros->status == TS_STEP_INCOMPLETE) { /* Use embedded completion formula */
9959566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol, U));
996de19f811SJed Brown       for (i = 0; i < tab->s; i++) w[i] = tab->bembedt[i];
9979566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(U, tab->s, w, ros->Y));
998108c343cSJed Brown     } else { /* Use rollback-and-recomplete formula (bembedt - bt) */
999108c343cSJed Brown       for (i = 0; i < tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i];
10009566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol, U));
10019566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(U, tab->s, w, ros->Y));
10021c3436cfSJed Brown     }
10031c3436cfSJed Brown     if (done) *done = PETSC_TRUE;
10043ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
10051c3436cfSJed Brown   }
10061c3436cfSJed Brown unavailable:
10071c3436cfSJed Brown   if (done) *done = PETSC_FALSE;
10089371c9d4SSatish Balay   else
10099371c9d4SSatish 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,
10109371c9d4SSatish Balay             tab->order, order);
10113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10121c3436cfSJed Brown }
10131c3436cfSJed Brown 
1014d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRollBack_RosW(TS ts)
1015d71ae5a4SJacob Faibussowitsch {
101624655328SShri   TS_RosW *ros = (TS_RosW *)ts->data;
101724655328SShri 
101824655328SShri   PetscFunctionBegin;
10199566063dSJacob Faibussowitsch   PetscCall(VecCopy(ros->vec_sol_prev, ts->vec_sol));
10203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
102124655328SShri }
102224655328SShri 
1023d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_RosW(TS ts)
1024d71ae5a4SJacob Faibussowitsch {
102561692a83SJed Brown   TS_RosW         *ros = (TS_RosW *)ts->data;
102661692a83SJed Brown   RosWTableau      tab = ros->tableau;
1027e27a552bSJed Brown   const PetscInt   s   = tab->s;
10281c3436cfSJed Brown   const PetscReal *At = tab->At, *Gamma = tab->Gamma, *ASum = tab->ASum, *GammaInv = tab->GammaInv;
10290feba352SEmil Constantinescu   const PetscReal *GammaExplicitCorr = tab->GammaExplicitCorr;
1030c17803e7SJed Brown   const PetscBool *GammaZeroDiag     = tab->GammaZeroDiag;
103161692a83SJed Brown   PetscScalar     *w                 = ros->work;
10327d4bf2deSEmil Constantinescu   Vec             *Y = ros->Y, Ydot = ros->Ydot, Zdot = ros->Zdot, Zstage = ros->Zstage;
1033e27a552bSJed Brown   SNES             snes;
10341c3436cfSJed Brown   TSAdapt          adapt;
1035fecfb714SLisandro Dalcin   PetscInt         i, j, its, lits;
1036be5899b3SLisandro Dalcin   PetscInt         rejections = 0;
1037b39943a6SLisandro Dalcin   PetscBool        stageok, accept = PETSC_TRUE;
1038b39943a6SLisandro Dalcin   PetscReal        next_time_step = ts->time_step;
1039f7f07198SBarry Smith   PetscInt         lag;
1040e27a552bSJed Brown 
1041e27a552bSJed Brown   PetscFunctionBegin;
104248a46eb9SPierre Jolivet   if (!ts->steprollback) PetscCall(VecCopy(ts->vec_sol, ros->vec_sol_prev));
1043e27a552bSJed Brown 
1044b39943a6SLisandro Dalcin   ros->status = TS_STEP_INCOMPLETE;
1045b39943a6SLisandro Dalcin   while (!ts->reason && ros->status != TS_STEP_COMPLETE) {
10461c3436cfSJed Brown     const PetscReal h = ts->time_step;
1047e27a552bSJed Brown     for (i = 0; i < s; i++) {
10481c3436cfSJed Brown       ros->stage_time = ts->ptime + h * ASum[i];
10499566063dSJacob Faibussowitsch       PetscCall(TSPreStage(ts, ros->stage_time));
1050c17803e7SJed Brown       if (GammaZeroDiag[i]) {
1051c17803e7SJed Brown         ros->stage_explicit = PETSC_TRUE;
1052b296d7d5SJed Brown         ros->scoeff         = 1.;
1053c17803e7SJed Brown       } else {
1054c17803e7SJed Brown         ros->stage_explicit = PETSC_FALSE;
1055b296d7d5SJed Brown         ros->scoeff         = 1. / Gamma[i * s + i];
1056fd96d5b0SEmil Constantinescu       }
105761692a83SJed Brown 
10589566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol, Zstage));
1059de19f811SJed Brown       for (j = 0; j < i; j++) w[j] = At[i * s + j];
10609566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(Zstage, i, w, Y));
106161692a83SJed Brown 
106261692a83SJed Brown       for (j = 0; j < i; j++) w[j] = 1. / h * GammaInv[i * s + j];
10639566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(Zdot));
10649566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(Zdot, i, w, Y));
106561692a83SJed Brown 
1066e27a552bSJed Brown       /* Initial guess taken from last stage */
10679566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(Y[i]));
106861692a83SJed Brown 
10697d4bf2deSEmil Constantinescu       if (!ros->stage_explicit) {
10709566063dSJacob Faibussowitsch         PetscCall(TSGetSNES(ts, &snes));
107161692a83SJed Brown         if (!ros->recompute_jacobian && !i) {
10729566063dSJacob Faibussowitsch           PetscCall(SNESGetLagJacobian(snes, &lag));
10736aad120cSJose E. Roman           if (lag == 1) {                            /* use did not set a nontrivial lag, so lag over all stages */
10749566063dSJacob Faibussowitsch             PetscCall(SNESSetLagJacobian(snes, -2)); /* Recompute the Jacobian on this solve, but not again for the rest of the stages */
1075f7f07198SBarry Smith           }
107661692a83SJed Brown         }
10779566063dSJacob Faibussowitsch         PetscCall(SNESSolve(snes, NULL, Y[i]));
10789371c9d4SSatish 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 */ }
10799566063dSJacob Faibussowitsch         PetscCall(SNESGetIterationNumber(snes, &its));
10809566063dSJacob Faibussowitsch         PetscCall(SNESGetLinearSolveIterations(snes, &lits));
10819371c9d4SSatish Balay         ts->snes_its += its;
10829371c9d4SSatish Balay         ts->ksp_its += lits;
10837d4bf2deSEmil Constantinescu       } else {
10841ce71dffSSatish Balay         Mat J, Jp;
10859566063dSJacob Faibussowitsch         PetscCall(VecZeroEntries(Ydot)); /* Evaluate Y[i]=G(t,Ydot=0,Zstage) */
10869566063dSJacob Faibussowitsch         PetscCall(TSComputeIFunction(ts, ros->stage_time, Zstage, Ydot, Y[i], PETSC_FALSE));
10879566063dSJacob Faibussowitsch         PetscCall(VecScale(Y[i], -1.0));
10889566063dSJacob Faibussowitsch         PetscCall(VecAXPY(Y[i], -1.0, Zdot)); /*Y[i] = F(Zstage)-Zdot[=GammaInv*Y]*/
10890feba352SEmil Constantinescu 
10909566063dSJacob Faibussowitsch         PetscCall(VecZeroEntries(Zstage)); /* Zstage = GammaExplicitCorr[i,j] * Y[j] */
10910feba352SEmil Constantinescu         for (j = 0; j < i; j++) w[j] = GammaExplicitCorr[i * s + j];
10929566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(Zstage, i, w, Y));
1093fecfb714SLisandro Dalcin 
1094fecfb714SLisandro Dalcin         /* Y[i] = Y[i] + Jac*Zstage[=Jac*GammaExplicitCorr[i,j] * Y[j]] */
10959566063dSJacob Faibussowitsch         PetscCall(TSGetIJacobian(ts, &J, &Jp, NULL, NULL));
10969566063dSJacob Faibussowitsch         PetscCall(TSComputeIJacobian(ts, ros->stage_time, ts->vec_sol, Ydot, 0, J, Jp, PETSC_FALSE));
10979566063dSJacob Faibussowitsch         PetscCall(MatMult(J, Zstage, Zdot));
10989566063dSJacob Faibussowitsch         PetscCall(VecAXPY(Y[i], -1.0, Zdot));
10995ef26d82SJed Brown         ts->ksp_its += 1;
1100fecfb714SLisandro Dalcin 
11019566063dSJacob Faibussowitsch         PetscCall(VecScale(Y[i], h));
11027d4bf2deSEmil Constantinescu       }
11039566063dSJacob Faibussowitsch       PetscCall(TSPostStage(ts, ros->stage_time, i, Y));
11049566063dSJacob Faibussowitsch       PetscCall(TSGetAdapt(ts, &adapt));
11059566063dSJacob Faibussowitsch       PetscCall(TSAdaptCheckStage(adapt, ts, ros->stage_time, Y[i], &stageok));
1106fecfb714SLisandro Dalcin       if (!stageok) goto reject_step;
1107e27a552bSJed Brown     }
1108e27a552bSJed Brown 
1109b39943a6SLisandro Dalcin     ros->status = TS_STEP_INCOMPLETE;
11109566063dSJacob Faibussowitsch     PetscCall(TSEvaluateStep_RosW(ts, tab->order, ts->vec_sol, NULL));
1111b39943a6SLisandro Dalcin     ros->status = TS_STEP_PENDING;
11129566063dSJacob Faibussowitsch     PetscCall(TSGetAdapt(ts, &adapt));
11139566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidatesClear(adapt));
11149566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidateAdd(adapt, tab->name, tab->order, 1, tab->ccfl, (PetscReal)tab->s, PETSC_TRUE));
11159566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
1116b39943a6SLisandro Dalcin     ros->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
1117b39943a6SLisandro Dalcin     if (!accept) { /* Roll back the current step */
11189566063dSJacob Faibussowitsch       PetscCall(TSRollBack_RosW(ts));
1119be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
1120be5899b3SLisandro Dalcin       goto reject_step;
1121b39943a6SLisandro Dalcin     }
1122b39943a6SLisandro Dalcin 
1123e27a552bSJed Brown     ts->ptime += ts->time_step;
1124cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
11251c3436cfSJed Brown     break;
1126b39943a6SLisandro Dalcin 
1127b39943a6SLisandro Dalcin   reject_step:
11289371c9d4SSatish Balay     ts->reject++;
11299371c9d4SSatish Balay     accept = PETSC_FALSE;
1130be5899b3SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
1131b39943a6SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
113263a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
11331c3436cfSJed Brown     }
11341c3436cfSJed Brown   }
11353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1136e27a552bSJed Brown }
1137e27a552bSJed Brown 
1138d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_RosW(TS ts, PetscReal itime, Vec U)
1139d71ae5a4SJacob Faibussowitsch {
114061692a83SJed Brown   TS_RosW         *ros = (TS_RosW *)ts->data;
1141f4aed992SEmil Constantinescu   PetscInt         s = ros->tableau->s, pinterp = ros->tableau->pinterp, i, j;
1142f4aed992SEmil Constantinescu   PetscReal        h;
1143f4aed992SEmil Constantinescu   PetscReal        tt, t;
1144f4aed992SEmil Constantinescu   PetscScalar     *bt;
1145f4aed992SEmil Constantinescu   const PetscReal *Bt       = ros->tableau->binterpt;
1146f4aed992SEmil Constantinescu   const PetscReal *GammaInv = ros->tableau->GammaInv;
1147f4aed992SEmil Constantinescu   PetscScalar     *w        = ros->work;
1148f4aed992SEmil Constantinescu   Vec             *Y        = ros->Y;
1149e27a552bSJed Brown 
1150e27a552bSJed Brown   PetscFunctionBegin;
11513c633725SBarry Smith   PetscCheck(Bt, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSRosW %s does not have an interpolation formula", ros->tableau->name);
1152f4aed992SEmil Constantinescu 
1153f4aed992SEmil Constantinescu   switch (ros->status) {
1154f4aed992SEmil Constantinescu   case TS_STEP_INCOMPLETE:
1155f4aed992SEmil Constantinescu   case TS_STEP_PENDING:
1156f4aed992SEmil Constantinescu     h = ts->time_step;
1157f4aed992SEmil Constantinescu     t = (itime - ts->ptime) / h;
1158f4aed992SEmil Constantinescu     break;
1159f4aed992SEmil Constantinescu   case TS_STEP_COMPLETE:
1160be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev;
1161f4aed992SEmil Constantinescu     t = (itime - ts->ptime) / h + 1; /* In the interval [0,1] */
1162f4aed992SEmil Constantinescu     break;
1163d71ae5a4SJacob Faibussowitsch   default:
1164d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus");
1165f4aed992SEmil Constantinescu   }
11669566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s, &bt));
1167f4aed992SEmil Constantinescu   for (i = 0; i < s; i++) bt[i] = 0;
1168f4aed992SEmil Constantinescu   for (j = 0, tt = t; j < pinterp; j++, tt *= t) {
1169ad540459SPierre Jolivet     for (i = 0; i < s; i++) bt[i] += Bt[i * pinterp + j] * tt;
1170f4aed992SEmil Constantinescu   }
1171f4aed992SEmil Constantinescu 
1172f4aed992SEmil Constantinescu   /* y(t+tt*h) = y(t) + Sum bt(tt) * GammaInv * Ydot */
1173f9c1d6abSBarry Smith   /* U <- 0*/
11749566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(U));
1175f9c1d6abSBarry Smith   /* U <- Sum bt_i * GammaInv(i,1:i) * Y(1:i) */
11763ca35412SEmil Constantinescu   for (j = 0; j < s; j++) w[j] = 0;
11773ca35412SEmil Constantinescu   for (j = 0; j < s; j++) {
1178ad540459SPierre Jolivet     for (i = j; i < s; i++) w[j] += bt[i] * GammaInv[i * s + j];
11793ca35412SEmil Constantinescu   }
11809566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(U, i, w, Y));
1181be5899b3SLisandro Dalcin   /* U <- y(t) + U */
11829566063dSJacob Faibussowitsch   PetscCall(VecAXPY(U, 1, ros->vec_sol_prev));
1183f4aed992SEmil Constantinescu 
11849566063dSJacob Faibussowitsch   PetscCall(PetscFree(bt));
11853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1186e27a552bSJed Brown }
1187e27a552bSJed Brown 
1188e27a552bSJed Brown /*------------------------------------------------------------*/
1189b39943a6SLisandro Dalcin 
1190d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWTableauReset(TS ts)
1191d71ae5a4SJacob Faibussowitsch {
1192b39943a6SLisandro Dalcin   TS_RosW    *ros = (TS_RosW *)ts->data;
1193b39943a6SLisandro Dalcin   RosWTableau tab = ros->tableau;
1194b39943a6SLisandro Dalcin 
1195b39943a6SLisandro Dalcin   PetscFunctionBegin;
11963ba16761SJacob Faibussowitsch   if (!tab) PetscFunctionReturn(PETSC_SUCCESS);
11979566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(tab->s, &ros->Y));
11989566063dSJacob Faibussowitsch   PetscCall(PetscFree(ros->work));
11993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1200b39943a6SLisandro Dalcin }
1201b39943a6SLisandro Dalcin 
1202d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_RosW(TS ts)
1203d71ae5a4SJacob Faibussowitsch {
120461692a83SJed Brown   TS_RosW *ros = (TS_RosW *)ts->data;
1205e27a552bSJed Brown 
1206e27a552bSJed Brown   PetscFunctionBegin;
12079566063dSJacob Faibussowitsch   PetscCall(TSRosWTableauReset(ts));
12089566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ros->Ydot));
12099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ros->Ystage));
12109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ros->Zdot));
12119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ros->Zstage));
12129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ros->vec_sol_prev));
12133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1214e27a552bSJed Brown }
1215e27a552bSJed Brown 
1216d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWGetVecs(TS ts, DM dm, Vec *Ydot, Vec *Zdot, Vec *Ystage, Vec *Zstage)
1217d71ae5a4SJacob Faibussowitsch {
1218d5e6173cSPeter Brune   TS_RosW *rw = (TS_RosW *)ts->data;
1219d5e6173cSPeter Brune 
1220d5e6173cSPeter Brune   PetscFunctionBegin;
1221d5e6173cSPeter Brune   if (Ydot) {
1222d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
12239566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Ydot", Ydot));
1224d5e6173cSPeter Brune     } else *Ydot = rw->Ydot;
1225d5e6173cSPeter Brune   }
1226d5e6173cSPeter Brune   if (Zdot) {
1227d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
12289566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Zdot", Zdot));
1229d5e6173cSPeter Brune     } else *Zdot = rw->Zdot;
1230d5e6173cSPeter Brune   }
1231d5e6173cSPeter Brune   if (Ystage) {
1232d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
12339566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Ystage", Ystage));
1234d5e6173cSPeter Brune     } else *Ystage = rw->Ystage;
1235d5e6173cSPeter Brune   }
1236d5e6173cSPeter Brune   if (Zstage) {
1237d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
12389566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Zstage", Zstage));
1239d5e6173cSPeter Brune     } else *Zstage = rw->Zstage;
1240d5e6173cSPeter Brune   }
12413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1242d5e6173cSPeter Brune }
1243d5e6173cSPeter Brune 
1244d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWRestoreVecs(TS ts, DM dm, Vec *Ydot, Vec *Zdot, Vec *Ystage, Vec *Zstage)
1245d71ae5a4SJacob Faibussowitsch {
1246d5e6173cSPeter Brune   PetscFunctionBegin;
1247d5e6173cSPeter Brune   if (Ydot) {
124848a46eb9SPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Ydot", Ydot));
1249d5e6173cSPeter Brune   }
1250d5e6173cSPeter Brune   if (Zdot) {
125148a46eb9SPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Zdot", Zdot));
1252d5e6173cSPeter Brune   }
1253d5e6173cSPeter Brune   if (Ystage) {
125448a46eb9SPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Ystage", Ystage));
1255d5e6173cSPeter Brune   }
1256d5e6173cSPeter Brune   if (Zstage) {
125748a46eb9SPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Zstage", Zstage));
1258d5e6173cSPeter Brune   }
12593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1260d5e6173cSPeter Brune }
1261d5e6173cSPeter Brune 
1262d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCoarsenHook_TSRosW(DM fine, DM coarse, void *ctx)
1263d71ae5a4SJacob Faibussowitsch {
1264d5e6173cSPeter Brune   PetscFunctionBegin;
12653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1266d5e6173cSPeter Brune }
1267d5e6173cSPeter Brune 
1268d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRestrictHook_TSRosW(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx)
1269d71ae5a4SJacob Faibussowitsch {
1270d5e6173cSPeter Brune   TS  ts = (TS)ctx;
1271d5e6173cSPeter Brune   Vec Ydot, Zdot, Ystage, Zstage;
1272d5e6173cSPeter Brune   Vec Ydotc, Zdotc, Ystagec, Zstagec;
1273d5e6173cSPeter Brune 
1274d5e6173cSPeter Brune   PetscFunctionBegin;
12759566063dSJacob Faibussowitsch   PetscCall(TSRosWGetVecs(ts, fine, &Ydot, &Ystage, &Zdot, &Zstage));
12769566063dSJacob Faibussowitsch   PetscCall(TSRosWGetVecs(ts, coarse, &Ydotc, &Ystagec, &Zdotc, &Zstagec));
12779566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct, Ydot, Ydotc));
12789566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(Ydotc, rscale, Ydotc));
12799566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct, Ystage, Ystagec));
12809566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(Ystagec, rscale, Ystagec));
12819566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct, Zdot, Zdotc));
12829566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(Zdotc, rscale, Zdotc));
12839566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct, Zstage, Zstagec));
12849566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(Zstagec, rscale, Zstagec));
12859566063dSJacob Faibussowitsch   PetscCall(TSRosWRestoreVecs(ts, fine, &Ydot, &Ystage, &Zdot, &Zstage));
12869566063dSJacob Faibussowitsch   PetscCall(TSRosWRestoreVecs(ts, coarse, &Ydotc, &Ystagec, &Zdotc, &Zstagec));
12873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1288d5e6173cSPeter Brune }
1289d5e6173cSPeter Brune 
1290d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainHook_TSRosW(DM fine, DM coarse, void *ctx)
1291d71ae5a4SJacob Faibussowitsch {
1292258e1594SPeter Brune   PetscFunctionBegin;
12933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1294258e1594SPeter Brune }
1295258e1594SPeter Brune 
1296d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainRestrictHook_TSRosW(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx)
1297d71ae5a4SJacob Faibussowitsch {
1298258e1594SPeter Brune   TS  ts = (TS)ctx;
1299258e1594SPeter Brune   Vec Ydot, Zdot, Ystage, Zstage;
1300258e1594SPeter Brune   Vec Ydots, Zdots, Ystages, Zstages;
1301258e1594SPeter Brune 
1302258e1594SPeter Brune   PetscFunctionBegin;
13039566063dSJacob Faibussowitsch   PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Ystage, &Zdot, &Zstage));
13049566063dSJacob Faibussowitsch   PetscCall(TSRosWGetVecs(ts, subdm, &Ydots, &Ystages, &Zdots, &Zstages));
1305258e1594SPeter Brune 
13069566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(gscat, Ydot, Ydots, INSERT_VALUES, SCATTER_FORWARD));
13079566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(gscat, Ydot, Ydots, INSERT_VALUES, SCATTER_FORWARD));
1308258e1594SPeter Brune 
13099566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(gscat, Ystage, Ystages, INSERT_VALUES, SCATTER_FORWARD));
13109566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(gscat, Ystage, Ystages, INSERT_VALUES, SCATTER_FORWARD));
1311258e1594SPeter Brune 
13129566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(gscat, Zdot, Zdots, INSERT_VALUES, SCATTER_FORWARD));
13139566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(gscat, Zdot, Zdots, INSERT_VALUES, SCATTER_FORWARD));
1314258e1594SPeter Brune 
13159566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(gscat, Zstage, Zstages, INSERT_VALUES, SCATTER_FORWARD));
13169566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(gscat, Zstage, Zstages, INSERT_VALUES, SCATTER_FORWARD));
1317258e1594SPeter Brune 
13189566063dSJacob Faibussowitsch   PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Ystage, &Zdot, &Zstage));
13199566063dSJacob Faibussowitsch   PetscCall(TSRosWRestoreVecs(ts, subdm, &Ydots, &Ystages, &Zdots, &Zstages));
13203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1321258e1594SPeter Brune }
1322258e1594SPeter Brune 
1323e27a552bSJed Brown /*
1324e27a552bSJed Brown   This defines the nonlinear equation that is to be solved with SNES
1325e27a552bSJed Brown   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
1326e27a552bSJed Brown */
1327d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_RosW(SNES snes, Vec U, Vec F, TS ts)
1328d71ae5a4SJacob Faibussowitsch {
132961692a83SJed Brown   TS_RosW  *ros = (TS_RosW *)ts->data;
1330d5e6173cSPeter Brune   Vec       Ydot, Zdot, Ystage, Zstage;
1331b296d7d5SJed Brown   PetscReal shift = ros->scoeff / ts->time_step;
1332d5e6173cSPeter Brune   DM        dm, dmsave;
1333e27a552bSJed Brown 
1334e27a552bSJed Brown   PetscFunctionBegin;
13359566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
13369566063dSJacob Faibussowitsch   PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage));
13379566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(Ydot, shift, U, Zdot));   /* Ydot = shift*U + Zdot */
13389566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(Ystage, 1.0, U, Zstage)); /* Ystage = U + Zstage */
1339d5e6173cSPeter Brune   dmsave = ts->dm;
1340d5e6173cSPeter Brune   ts->dm = dm;
13419566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts, ros->stage_time, Ystage, Ydot, F, PETSC_FALSE));
1342d5e6173cSPeter Brune   ts->dm = dmsave;
13439566063dSJacob Faibussowitsch   PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage));
13443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1345e27a552bSJed Brown }
1346e27a552bSJed Brown 
1347d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes, Vec U, Mat A, Mat B, TS ts)
1348d71ae5a4SJacob Faibussowitsch {
134961692a83SJed Brown   TS_RosW  *ros = (TS_RosW *)ts->data;
1350d5e6173cSPeter Brune   Vec       Ydot, Zdot, Ystage, Zstage;
1351b296d7d5SJed Brown   PetscReal shift = ros->scoeff / ts->time_step;
1352d5e6173cSPeter Brune   DM        dm, dmsave;
1353e27a552bSJed Brown 
1354e27a552bSJed Brown   PetscFunctionBegin;
135561692a83SJed Brown   /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
13569566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
13579566063dSJacob Faibussowitsch   PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage));
1358d5e6173cSPeter Brune   dmsave = ts->dm;
1359d5e6173cSPeter Brune   ts->dm = dm;
13609566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts, ros->stage_time, Ystage, Ydot, shift, A, B, PETSC_TRUE));
1361d5e6173cSPeter Brune   ts->dm = dmsave;
13629566063dSJacob Faibussowitsch   PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage));
13633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1364e27a552bSJed Brown }
1365e27a552bSJed Brown 
1366d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWTableauSetUp(TS ts)
1367d71ae5a4SJacob Faibussowitsch {
1368b39943a6SLisandro Dalcin   TS_RosW    *ros = (TS_RosW *)ts->data;
1369b39943a6SLisandro Dalcin   RosWTableau tab = ros->tableau;
1370b39943a6SLisandro Dalcin 
1371b39943a6SLisandro Dalcin   PetscFunctionBegin;
13729566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ros->Y));
13739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tab->s, &ros->work));
13743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1375b39943a6SLisandro Dalcin }
1376b39943a6SLisandro Dalcin 
1377d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_RosW(TS ts)
1378d71ae5a4SJacob Faibussowitsch {
137961692a83SJed Brown   TS_RosW      *ros = (TS_RosW *)ts->data;
1380d5e6173cSPeter Brune   DM            dm;
1381b39943a6SLisandro Dalcin   SNES          snes;
1382a3ab5968SHong Zhang   TSRHSJacobian rhsjacobian;
1383e27a552bSJed Brown 
1384e27a552bSJed Brown   PetscFunctionBegin;
13859566063dSJacob Faibussowitsch   PetscCall(TSRosWTableauSetUp(ts));
13869566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &ros->Ydot));
13879566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &ros->Ystage));
13889566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &ros->Zdot));
13899566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &ros->Zstage));
13909566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &ros->vec_sol_prev));
13919566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
13929566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSRosW, DMRestrictHook_TSRosW, ts));
13939566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSRosW, DMSubDomainRestrictHook_TSRosW, ts));
1394b39943a6SLisandro Dalcin   /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
13959566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
139648a46eb9SPierre Jolivet   if (!((PetscObject)snes)->type_name) PetscCall(SNESSetType(snes, SNESKSPONLY));
13979566063dSJacob Faibussowitsch   PetscCall(DMTSGetRHSJacobian(dm, &rhsjacobian, NULL));
1398a3ab5968SHong Zhang   if (rhsjacobian == TSComputeRHSJacobianConstant) {
1399a3ab5968SHong Zhang     Mat Amat, Pmat;
1400a3ab5968SHong Zhang 
1401a3ab5968SHong Zhang     /* Set the SNES matrix to be different from the RHS matrix because there is no way to reconstruct shift*M-J */
14029566063dSJacob Faibussowitsch     PetscCall(SNESGetJacobian(snes, &Amat, &Pmat, NULL, NULL));
1403a3ab5968SHong Zhang     if (Amat && Amat == ts->Arhs) {
1404a3ab5968SHong Zhang       if (Amat == Pmat) {
14059566063dSJacob Faibussowitsch         PetscCall(MatDuplicate(ts->Arhs, MAT_COPY_VALUES, &Amat));
14069566063dSJacob Faibussowitsch         PetscCall(SNESSetJacobian(snes, Amat, Amat, NULL, NULL));
1407a3ab5968SHong Zhang       } else {
14089566063dSJacob Faibussowitsch         PetscCall(MatDuplicate(ts->Arhs, MAT_COPY_VALUES, &Amat));
14099566063dSJacob Faibussowitsch         PetscCall(SNESSetJacobian(snes, Amat, NULL, NULL, NULL));
1410a3ab5968SHong Zhang         if (Pmat && Pmat == ts->Brhs) {
14119566063dSJacob Faibussowitsch           PetscCall(MatDuplicate(ts->Brhs, MAT_COPY_VALUES, &Pmat));
14129566063dSJacob Faibussowitsch           PetscCall(SNESSetJacobian(snes, NULL, Pmat, NULL, NULL));
14139566063dSJacob Faibussowitsch           PetscCall(MatDestroy(&Pmat));
1414a3ab5968SHong Zhang         }
1415a3ab5968SHong Zhang       }
14169566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Amat));
1417a3ab5968SHong Zhang     }
1418a3ab5968SHong Zhang   }
14193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1420e27a552bSJed Brown }
1421e27a552bSJed Brown /*------------------------------------------------------------*/
1422e27a552bSJed Brown 
1423d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_RosW(TS ts, PetscOptionItems *PetscOptionsObject)
1424d71ae5a4SJacob Faibussowitsch {
142561692a83SJed Brown   TS_RosW *ros = (TS_RosW *)ts->data;
1426b39943a6SLisandro Dalcin   SNES     snes;
1427e27a552bSJed Brown 
1428e27a552bSJed Brown   PetscFunctionBegin;
1429d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "RosW ODE solver options");
1430e27a552bSJed Brown   {
143161692a83SJed Brown     RosWTableauLink link;
1432e27a552bSJed Brown     PetscInt        count, choice;
1433e27a552bSJed Brown     PetscBool       flg;
1434e27a552bSJed Brown     const char    **namelist;
143561692a83SJed Brown 
14369371c9d4SSatish Balay     for (link = RosWTableauList, count = 0; link; link = link->next, count++)
14379371c9d4SSatish Balay       ;
14389566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(count, (char ***)&namelist));
143961692a83SJed Brown     for (link = RosWTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name;
14409566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEList("-ts_rosw_type", "Family of Rosenbrock-W method", "TSRosWSetType", (const char *const *)namelist, count, ros->tableau->name, &choice, &flg));
14419566063dSJacob Faibussowitsch     if (flg) PetscCall(TSRosWSetType(ts, namelist[choice]));
14429566063dSJacob Faibussowitsch     PetscCall(PetscFree(namelist));
144361692a83SJed Brown 
14449566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ts_rosw_recompute_jacobian", "Recompute the Jacobian at each stage", "TSRosWSetRecomputeJacobian", ros->recompute_jacobian, &ros->recompute_jacobian, NULL));
1445b39943a6SLisandro Dalcin   }
1446d0609cedSBarry Smith   PetscOptionsHeadEnd();
144761692a83SJed Brown   /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
14489566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
144948a46eb9SPierre Jolivet   if (!((PetscObject)snes)->type_name) PetscCall(SNESSetType(snes, SNESKSPONLY));
14503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1451e27a552bSJed Brown }
1452e27a552bSJed Brown 
1453d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_RosW(TS ts, PetscViewer viewer)
1454d71ae5a4SJacob Faibussowitsch {
145561692a83SJed Brown   TS_RosW  *ros = (TS_RosW *)ts->data;
1456e27a552bSJed Brown   PetscBool iascii;
1457e27a552bSJed Brown 
1458e27a552bSJed Brown   PetscFunctionBegin;
14599566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1460e27a552bSJed Brown   if (iascii) {
14619c334d8fSLisandro Dalcin     RosWTableau tab = ros->tableau;
146219fd82e9SBarry Smith     TSRosWType  rostype;
14639c334d8fSLisandro Dalcin     char        buf[512];
1464e408995aSJed Brown     PetscInt    i;
1465e408995aSJed Brown     PetscReal   abscissa[512];
14669566063dSJacob Faibussowitsch     PetscCall(TSRosWGetType(ts, &rostype));
14679566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Rosenbrock-W %s\n", rostype));
14689566063dSJacob Faibussowitsch     PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, tab->ASum));
14699566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Abscissa of A       = %s\n", buf));
1470e408995aSJed Brown     for (i = 0; i < tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i];
14719566063dSJacob Faibussowitsch     PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, abscissa));
14729566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Abscissa of A+Gamma = %s\n", buf));
1473e27a552bSJed Brown   }
14743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1475e27a552bSJed Brown }
1476e27a552bSJed Brown 
1477d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSLoad_RosW(TS ts, PetscViewer viewer)
1478d71ae5a4SJacob Faibussowitsch {
14799200755eSBarry Smith   SNES    snes;
14809c334d8fSLisandro Dalcin   TSAdapt adapt;
14819200755eSBarry Smith 
14829200755eSBarry Smith   PetscFunctionBegin;
14839566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &adapt));
14849566063dSJacob Faibussowitsch   PetscCall(TSAdaptLoad(adapt, viewer));
14859566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
14869566063dSJacob Faibussowitsch   PetscCall(SNESLoad(snes, viewer));
14879200755eSBarry Smith   /* function and Jacobian context for SNES when used with TS is always ts object */
14889566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, NULL, NULL, ts));
14899566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, ts));
14903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14919200755eSBarry Smith }
14929200755eSBarry Smith 
1493e27a552bSJed Brown /*@C
1494bcf0153eSBarry Smith   TSRosWSetType - Set the type of Rosenbrock-W, `TSROSW`, scheme
1495e27a552bSJed Brown 
149620f4b53cSBarry Smith   Logically Collective
1497e27a552bSJed Brown 
1498d8d19677SJose E. Roman   Input Parameters:
1499e27a552bSJed Brown + ts       - timestepping context
1500b92453a8SLisandro Dalcin - roswtype - type of Rosenbrock-W scheme
1501e27a552bSJed Brown 
1502020d8f30SJed Brown   Level: beginner
1503e27a552bSJed Brown 
15041cc06b55SBarry Smith .seealso: [](ch_ts), `TSRosWGetType()`, `TSROSW`, `TSROSW2M`, `TSROSW2P`, `TSROSWRA3PW`, `TSROSWRA34PW2`, `TSROSWRODAS3`, `TSROSWSANDU3`, `TSROSWASSP3P3S1C`, `TSROSWLASSP3P4S2C`, `TSROSWLLSSP3P4S2C`, `TSROSWARK3`
1505e27a552bSJed Brown @*/
1506d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWSetType(TS ts, TSRosWType roswtype)
1507d71ae5a4SJacob Faibussowitsch {
1508e27a552bSJed Brown   PetscFunctionBegin;
1509e27a552bSJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
15104f572ea9SToby Isaac   PetscAssertPointer(roswtype, 2);
1511cac4c232SBarry Smith   PetscTryMethod(ts, "TSRosWSetType_C", (TS, TSRosWType), (ts, roswtype));
15123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1513e27a552bSJed Brown }
1514e27a552bSJed Brown 
1515e27a552bSJed Brown /*@C
151661692a83SJed Brown   TSRosWGetType - Get the type of Rosenbrock-W scheme
1517e27a552bSJed Brown 
151820f4b53cSBarry Smith   Logically Collective
1519e27a552bSJed Brown 
1520e27a552bSJed Brown   Input Parameter:
1521e27a552bSJed Brown . ts - timestepping context
1522e27a552bSJed Brown 
1523e27a552bSJed Brown   Output Parameter:
152461692a83SJed Brown . rostype - type of Rosenbrock-W scheme
1525e27a552bSJed Brown 
1526e27a552bSJed Brown   Level: intermediate
1527e27a552bSJed Brown 
15281cc06b55SBarry Smith .seealso: [](ch_ts), `TSRosWType`, `TSRosWSetType()`
1529e27a552bSJed Brown @*/
1530d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWGetType(TS ts, TSRosWType *rostype)
1531d71ae5a4SJacob Faibussowitsch {
1532e27a552bSJed Brown   PetscFunctionBegin;
1533e27a552bSJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1534cac4c232SBarry Smith   PetscUseMethod(ts, "TSRosWGetType_C", (TS, TSRosWType *), (ts, rostype));
15353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1536e27a552bSJed Brown }
1537e27a552bSJed Brown 
1538e27a552bSJed Brown /*@C
153961692a83SJed Brown   TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step.
1540e27a552bSJed Brown 
154120f4b53cSBarry Smith   Logically Collective
1542e27a552bSJed Brown 
1543d8d19677SJose E. Roman   Input Parameters:
1544e27a552bSJed Brown + ts  - timestepping context
1545bcf0153eSBarry Smith - flg - `PETSC_TRUE` to recompute the Jacobian at each stage
1546e27a552bSJed Brown 
1547e27a552bSJed Brown   Level: intermediate
1548e27a552bSJed Brown 
15491cc06b55SBarry Smith .seealso: [](ch_ts), `TSRosWType`, `TSRosWGetType()`
1550e27a552bSJed Brown @*/
1551d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWSetRecomputeJacobian(TS ts, PetscBool flg)
1552d71ae5a4SJacob Faibussowitsch {
1553e27a552bSJed Brown   PetscFunctionBegin;
1554e27a552bSJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1555cac4c232SBarry Smith   PetscTryMethod(ts, "TSRosWSetRecomputeJacobian_C", (TS, PetscBool), (ts, flg));
15563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1557e27a552bSJed Brown }
1558e27a552bSJed Brown 
1559d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWGetType_RosW(TS ts, TSRosWType *rostype)
1560d71ae5a4SJacob Faibussowitsch {
156161692a83SJed Brown   TS_RosW *ros = (TS_RosW *)ts->data;
1562e27a552bSJed Brown 
1563e27a552bSJed Brown   PetscFunctionBegin;
156461692a83SJed Brown   *rostype = ros->tableau->name;
15653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1566e27a552bSJed Brown }
1567ef20d060SBarry Smith 
1568d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWSetType_RosW(TS ts, TSRosWType rostype)
1569d71ae5a4SJacob Faibussowitsch {
157061692a83SJed Brown   TS_RosW        *ros = (TS_RosW *)ts->data;
1571e27a552bSJed Brown   PetscBool       match;
157261692a83SJed Brown   RosWTableauLink link;
1573e27a552bSJed Brown 
1574e27a552bSJed Brown   PetscFunctionBegin;
157561692a83SJed Brown   if (ros->tableau) {
15769566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(ros->tableau->name, rostype, &match));
15773ba16761SJacob Faibussowitsch     if (match) PetscFunctionReturn(PETSC_SUCCESS);
1578e27a552bSJed Brown   }
157961692a83SJed Brown   for (link = RosWTableauList; link; link = link->next) {
15809566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(link->tab.name, rostype, &match));
1581e27a552bSJed Brown     if (match) {
15829566063dSJacob Faibussowitsch       if (ts->setupcalled) PetscCall(TSRosWTableauReset(ts));
158361692a83SJed Brown       ros->tableau = &link->tab;
15849566063dSJacob Faibussowitsch       if (ts->setupcalled) PetscCall(TSRosWTableauSetUp(ts));
15852ffb9264SLisandro Dalcin       ts->default_adapt_type = ros->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
15863ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
1587e27a552bSJed Brown     }
1588e27a552bSJed Brown   }
158998921bdaSJacob Faibussowitsch   SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", rostype);
1590e27a552bSJed Brown }
159161692a83SJed Brown 
1592d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWSetRecomputeJacobian_RosW(TS ts, PetscBool flg)
1593d71ae5a4SJacob Faibussowitsch {
159461692a83SJed Brown   TS_RosW *ros = (TS_RosW *)ts->data;
1595e27a552bSJed Brown 
1596e27a552bSJed Brown   PetscFunctionBegin;
159761692a83SJed Brown   ros->recompute_jacobian = flg;
15983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1599e27a552bSJed Brown }
1600e27a552bSJed Brown 
1601d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_RosW(TS ts)
1602d71ae5a4SJacob Faibussowitsch {
1603b3a6b972SJed Brown   PetscFunctionBegin;
16049566063dSJacob Faibussowitsch   PetscCall(TSReset_RosW(ts));
1605b3a6b972SJed Brown   if (ts->dm) {
16069566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSRosW, DMRestrictHook_TSRosW, ts));
16079566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSRosW, DMSubDomainRestrictHook_TSRosW, ts));
1608b3a6b972SJed Brown   }
16099566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
16109566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWGetType_C", NULL));
16119566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetType_C", NULL));
16129566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetRecomputeJacobian_C", NULL));
16133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1614b3a6b972SJed Brown }
1615d5e6173cSPeter Brune 
1616e27a552bSJed Brown /* ------------------------------------------------------------ */
1617e27a552bSJed Brown /*MC
1618020d8f30SJed Brown       TSROSW - ODE solver using Rosenbrock-W schemes
1619e27a552bSJed Brown 
1620e27a552bSJed Brown   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
1621e27a552bSJed Brown   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
1622bcf0153eSBarry Smith   of the equation using `TSSetIFunction()` and the non-stiff part with `TSSetRHSFunction()`.
1623bcf0153eSBarry Smith 
1624bcf0153eSBarry Smith   Level: beginner
1625e27a552bSJed Brown 
1626e27a552bSJed Brown   Notes:
162761692a83SJed Brown   This method currently only works with autonomous ODE and DAE.
162861692a83SJed Brown 
1629bcf0153eSBarry Smith   Consider trying `TSARKIMEX` if the stiff part is strongly nonlinear.
1630d0685a90SJed Brown 
16313d5a8a6aSBarry 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
16323d5a8a6aSBarry Smith 
1633e94cfbe0SPatrick Sanan   Developer Notes:
163461692a83SJed Brown   Rosenbrock-W methods are typically specified for autonomous ODE
163561692a83SJed Brown 
1636f9c1d6abSBarry Smith $  udot = f(u)
163761692a83SJed Brown 
163861692a83SJed Brown   by the stage equations
163961692a83SJed Brown 
1640f9c1d6abSBarry Smith $  k_i = h f(u_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j
164161692a83SJed Brown 
164261692a83SJed Brown   and step completion formula
164361692a83SJed Brown 
1644f9c1d6abSBarry Smith $  u_1 = u_0 + sum_j b_j k_j
164561692a83SJed Brown 
1646f9c1d6abSBarry Smith   with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(u)
164761692a83SJed Brown   and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner,
164861692a83SJed Brown   we define new variables for the stage equations
164961692a83SJed Brown 
165061692a83SJed Brown $  y_i = gamma_ij k_j
165161692a83SJed Brown 
165261692a83SJed Brown   The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define
165361692a83SJed Brown 
1654b70472e2SJed Brown $  A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-1}
165561692a83SJed Brown 
165661692a83SJed Brown   to rewrite the method as
165761692a83SJed Brown 
165837fdd005SBarry Smith .vb
165937fdd005SBarry 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
166037fdd005SBarry Smith   u_1 = u_0 + sum_j bt_j y_j
166137fdd005SBarry Smith .ve
166261692a83SJed Brown 
166361692a83SJed Brown    where we have introduced the mass matrix M. Continue by defining
166461692a83SJed Brown 
166561692a83SJed Brown $  ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j
166661692a83SJed Brown 
166761692a83SJed Brown    or, more compactly in tensor notation
166861692a83SJed Brown 
166961692a83SJed Brown $  Ydot = 1/h (Gamma^{-1} \otimes I) Y .
167061692a83SJed Brown 
167161692a83SJed Brown    Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current
167261692a83SJed Brown    stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the
167361692a83SJed Brown    equation
167461692a83SJed Brown 
1675f9c1d6abSBarry Smith $  g(u_0 + sum_j a_ij y_j + y_i, ydot_i) = 0
167661692a83SJed Brown 
167761692a83SJed Brown    with initial guess y_i = 0.
1678e27a552bSJed Brown 
16791cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSRosWSetType()`, `TSRosWRegister()`, `TSROSWTHETA1`, `TSROSWTHETA2`, `TSROSW2M`, `TSROSW2P`, `TSROSWRA3PW`, `TSROSWRA34PW2`, `TSROSWRODAS3`,
1680bcf0153eSBarry Smith           `TSROSWSANDU3`, `TSROSWASSP3P3S1C`, `TSROSWLASSP3P4S2C`, `TSROSWLLSSP3P4S2C`, `TSROSWGRK4T`, `TSROSWSHAMP4`, `TSROSWVELDD4`, `TSROSW4L`, `TSType`
1681e27a552bSJed Brown M*/
1682d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_RosW(TS ts)
1683d71ae5a4SJacob Faibussowitsch {
168461692a83SJed Brown   TS_RosW *ros;
1685e27a552bSJed Brown 
1686e27a552bSJed Brown   PetscFunctionBegin;
16879566063dSJacob Faibussowitsch   PetscCall(TSRosWInitializePackage());
1688e27a552bSJed Brown 
1689e27a552bSJed Brown   ts->ops->reset          = TSReset_RosW;
1690e27a552bSJed Brown   ts->ops->destroy        = TSDestroy_RosW;
1691e27a552bSJed Brown   ts->ops->view           = TSView_RosW;
16929200755eSBarry Smith   ts->ops->load           = TSLoad_RosW;
1693e27a552bSJed Brown   ts->ops->setup          = TSSetUp_RosW;
1694e27a552bSJed Brown   ts->ops->step           = TSStep_RosW;
1695e27a552bSJed Brown   ts->ops->interpolate    = TSInterpolate_RosW;
16961c3436cfSJed Brown   ts->ops->evaluatestep   = TSEvaluateStep_RosW;
169724655328SShri   ts->ops->rollback       = TSRollBack_RosW;
1698e27a552bSJed Brown   ts->ops->setfromoptions = TSSetFromOptions_RosW;
1699e27a552bSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_RosW;
1700e27a552bSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_RosW;
1701e27a552bSJed Brown 
1702efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
1703efd4aadfSBarry Smith 
17044dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&ros));
170561692a83SJed Brown   ts->data = (void *)ros;
1706e27a552bSJed Brown 
17079566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWGetType_C", TSRosWGetType_RosW));
17089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetType_C", TSRosWSetType_RosW));
17099566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetRecomputeJacobian_C", TSRosWSetRecomputeJacobian_RosW));
1710b39943a6SLisandro Dalcin 
17119566063dSJacob Faibussowitsch   PetscCall(TSRosWSetType(ts, TSRosWDefault));
17123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1713e27a552bSJed Brown }
1714