xref: /petsc/src/ts/impls/rosw/rosw.c (revision 1d27aa22b2f6148b2c4e3f06a75e0638d6493e09)
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
107*1d27aa22SBarry Smith      TSROSWRA3PW - Three stage third order Rosenbrock-W scheme for PDAE of index 1 {cite}`rang_2005`
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 
1151cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`
116fe7e6d57SJed Brown M*/
117fe7e6d57SJed Brown 
118fe7e6d57SJed Brown /*MC
119*1d27aa22SBarry Smith      TSROSWRA34PW2 - Four stage third order L-stable Rosenbrock-W scheme for PDAE of index 1 {cite}`rang_2005`.
120fe7e6d57SJed Brown 
121fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.
122fe7e6d57SJed Brown 
123fe7e6d57SJed 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.
124fe7e6d57SJed Brown 
125bcf0153eSBarry Smith      Level: intermediate
126bcf0153eSBarry Smith 
1271cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`
128fe7e6d57SJed Brown M*/
129fe7e6d57SJed Brown 
130ef3c5b88SJed Brown /*MC
131*1d27aa22SBarry Smith      TSROSWRODAS3 - Four stage third order L-stable Rosenbrock scheme {cite}`sandu_1997`
132ef3c5b88SJed Brown 
133ef3c5b88SJed Brown      By default, the Jacobian is only recomputed once per step.
134ef3c5b88SJed Brown 
135ef3c5b88SJed Brown      Both the third order and embedded second order methods are stiffly accurate and L-stable.
136ef3c5b88SJed Brown 
137bcf0153eSBarry Smith      Level: intermediate
138bcf0153eSBarry Smith 
1391cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWSANDU3`
140ef3c5b88SJed Brown M*/
141ef3c5b88SJed Brown 
142ef3c5b88SJed Brown /*MC
143*1d27aa22SBarry Smith      TSROSWSANDU3 - Three stage third order L-stable Rosenbrock scheme {cite}`sandu_1997`
144ef3c5b88SJed Brown 
145ef3c5b88SJed Brown      By default, the Jacobian is only recomputed once per step.
146ef3c5b88SJed Brown 
147ef3c5b88SJed Brown      The third order method is L-stable, but not stiffly accurate.
148ef3c5b88SJed Brown      The second order embedded method is strongly A-stable with R(infty) = 0.5.
149ef3c5b88SJed Brown      The internal stages are L-stable.
150*1d27aa22SBarry Smith      This method is called ROS3 in {cite}`sandu_1997`.
151ef3c5b88SJed Brown 
152bcf0153eSBarry Smith      Level: intermediate
153bcf0153eSBarry Smith 
1541cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWRODAS3`
155ef3c5b88SJed Brown M*/
156ef3c5b88SJed Brown 
157961f28d0SJed Brown /*MC
158961f28d0SJed Brown      TSROSWASSP3P3S1C - A-stable Rosenbrock-W method with SSP explicit part, third order, three stages
159961f28d0SJed Brown 
160961f28d0SJed Brown      By default, the Jacobian is only recomputed once per step.
161961f28d0SJed Brown 
162961f28d0SJed Brown      A-stable SPP explicit order 3, 3 stages, CFL 1 (eff = 1/3)
163961f28d0SJed Brown 
164bcf0153eSBarry Smith      Level: intermediate
165bcf0153eSBarry Smith 
1661cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWLASSP3P4S2C`, `TSROSWLLSSP3P4S2C`, `SSP`
167961f28d0SJed Brown M*/
168961f28d0SJed Brown 
169961f28d0SJed Brown /*MC
170998eb97aSJed Brown      TSROSWLASSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, four stages
171961f28d0SJed Brown 
172961f28d0SJed Brown      By default, the Jacobian is only recomputed once per step.
173961f28d0SJed Brown 
174961f28d0SJed Brown      L-stable (A-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2)
175961f28d0SJed Brown 
176bcf0153eSBarry Smith      Level: intermediate
177bcf0153eSBarry Smith 
1781cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWASSP3P3S1C`, `TSROSWLLSSP3P4S2C`, `TSSSP`
179961f28d0SJed Brown M*/
180961f28d0SJed Brown 
181961f28d0SJed Brown /*MC
182998eb97aSJed Brown      TSROSWLLSSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, four stages
183961f28d0SJed Brown 
184961f28d0SJed Brown      By default, the Jacobian is only recomputed once per step.
185961f28d0SJed Brown 
186961f28d0SJed Brown      L-stable (L-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2)
187961f28d0SJed Brown 
188bcf0153eSBarry Smith      Level: intermediate
189bcf0153eSBarry Smith 
1901cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWASSP3P3S1C`, `TSROSWLASSP3P4S2C`, `TSSSP`
191961f28d0SJed Brown M*/
192961f28d0SJed Brown 
19342faf41dSJed Brown /*MC
194*1d27aa22SBarry Smith      TSROSWGRK4T - four stage, fourth order Rosenbrock (not W) method from Kaps and Rentrop {cite}`kaps1979generalized`
19542faf41dSJed Brown 
19642faf41dSJed Brown      By default, the Jacobian is only recomputed once per step.
19742faf41dSJed Brown 
19842faf41dSJed Brown      A(89.3 degrees)-stable, |R(infty)| = 0.454.
19942faf41dSJed Brown 
20042faf41dSJed Brown      This method does not provide a dense output formula.
20142faf41dSJed Brown 
202bcf0153eSBarry Smith      Level: intermediate
203bcf0153eSBarry Smith 
204*1d27aa22SBarry Smith      Note:
205*1d27aa22SBarry Smith      See Section 4 Table 7.2 in {cite}`wanner1996solving`
20642faf41dSJed Brown 
2071cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWSHAMP4`, `TSROSWVELDD4`, `TSROSW4L`
20842faf41dSJed Brown M*/
20942faf41dSJed Brown 
21042faf41dSJed Brown /*MC
211*1d27aa22SBarry Smith      TSROSWSHAMP4 - four stage, fourth order Rosenbrock (not W) method from Shampine {cite}`shampine1982implementation`
21242faf41dSJed Brown 
21342faf41dSJed Brown      By default, the Jacobian is only recomputed once per step.
21442faf41dSJed Brown 
21542faf41dSJed Brown      A-stable, |R(infty)| = 1/3.
21642faf41dSJed Brown 
21742faf41dSJed Brown      This method does not provide a dense output formula.
21842faf41dSJed Brown 
219bcf0153eSBarry Smith      Level: intermediate
220bcf0153eSBarry Smith 
221*1d27aa22SBarry Smith      Note:
222*1d27aa22SBarry Smith      See Section 4 Table 7.2 in in {cite}`wanner1996solving`
22342faf41dSJed Brown 
2241cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWGRK4T`, `TSROSWVELDD4`, `TSROSW4L`
22542faf41dSJed Brown M*/
22642faf41dSJed Brown 
22742faf41dSJed Brown /*MC
228*1d27aa22SBarry Smith      TSROSWVELDD4 - four stage, fourth order Rosenbrock (not W) method from van Veldhuizen {cite}`veldhuizen1984d`
22942faf41dSJed Brown 
23042faf41dSJed Brown      By default, the Jacobian is only recomputed once per step.
23142faf41dSJed Brown 
23242faf41dSJed Brown      A(89.5 degrees)-stable, |R(infty)| = 0.24.
23342faf41dSJed Brown 
23442faf41dSJed Brown      This method does not provide a dense output formula.
23542faf41dSJed Brown 
236bcf0153eSBarry Smith      Level: intermediate
237bcf0153eSBarry Smith 
238*1d27aa22SBarry Smith      Note:
239*1d27aa22SBarry Smith      See Section 4 Table 7.2 in {cite}`wanner1996solving`
24042faf41dSJed Brown 
2411cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWGRK4T`, `TSROSWSHAMP4`, `TSROSW4L`
24242faf41dSJed Brown M*/
24342faf41dSJed Brown 
24442faf41dSJed Brown /*MC
24542faf41dSJed Brown      TSROSW4L - four stage, fourth order Rosenbrock (not W) method
24642faf41dSJed Brown 
24742faf41dSJed Brown      By default, the Jacobian is only recomputed once per step.
24842faf41dSJed Brown 
24942faf41dSJed Brown      A-stable and L-stable
25042faf41dSJed Brown 
25142faf41dSJed Brown      This method does not provide a dense output formula.
25242faf41dSJed Brown 
253bcf0153eSBarry Smith      Level: intermediate
254bcf0153eSBarry Smith 
255*1d27aa22SBarry Smith      Note:
256*1d27aa22SBarry Smith      See Section 4 Table 7.2 in in {cite}`wanner1996solving`
25742faf41dSJed Brown 
2581cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSROSWGRK4T`, `TSROSWSHAMP4`, `TSROSW4L`
25942faf41dSJed Brown M*/
26042faf41dSJed Brown 
261e27a552bSJed Brown /*@C
262bcf0153eSBarry Smith   TSRosWRegisterAll - Registers all of the Rosenbrock-W methods in `TSROSW`
263e27a552bSJed Brown 
264*1d27aa22SBarry Smith   Not Collective, but should be called by all MPI processes which will need the schemes to be registered
265e27a552bSJed Brown 
266e27a552bSJed Brown   Level: advanced
267e27a552bSJed Brown 
2681cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `TSRosWRegisterDestroy()`
269e27a552bSJed Brown @*/
270d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWRegisterAll(void)
271d71ae5a4SJacob Faibussowitsch {
272e27a552bSJed Brown   PetscFunctionBegin;
2733ba16761SJacob Faibussowitsch   if (TSRosWRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
274e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_TRUE;
275e27a552bSJed Brown 
276e27a552bSJed Brown   {
277bbd56ea5SKarl Rupp     const PetscReal A        = 0;
278bbd56ea5SKarl Rupp     const PetscReal Gamma    = 1;
279bbd56ea5SKarl Rupp     const PetscReal b        = 1;
280bbd56ea5SKarl Rupp     const PetscReal binterpt = 1;
2811f80e275SEmil Constantinescu 
2829566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWTHETA1, 1, 1, &A, &Gamma, &b, NULL, 1, &binterpt));
2833606a31eSEmil Constantinescu   }
2843606a31eSEmil Constantinescu 
2853606a31eSEmil Constantinescu   {
286bbd56ea5SKarl Rupp     const PetscReal A        = 0;
287bbd56ea5SKarl Rupp     const PetscReal Gamma    = 0.5;
288bbd56ea5SKarl Rupp     const PetscReal b        = 1;
289bbd56ea5SKarl Rupp     const PetscReal binterpt = 1;
290bbd56ea5SKarl Rupp 
2919566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWTHETA2, 2, 1, &A, &Gamma, &b, NULL, 1, &binterpt));
2923606a31eSEmil Constantinescu   }
2933606a31eSEmil Constantinescu 
2943606a31eSEmil Constantinescu   {
295da80777bSKarl 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. */
296b16ce868SStefano Zampini     const PetscReal A[2][2] = {
2979371c9d4SSatish Balay       {0,  0},
2989371c9d4SSatish Balay       {1., 0}
299b16ce868SStefano Zampini     };
300b16ce868SStefano Zampini     const PetscReal Gamma[2][2] = {
301b16ce868SStefano Zampini       {1.707106781186547524401,       0                      },
302b16ce868SStefano Zampini       {-2. * 1.707106781186547524401, 1.707106781186547524401}
303b16ce868SStefano Zampini     };
304b16ce868SStefano Zampini     const PetscReal b[2]  = {0.5, 0.5};
305b16ce868SStefano Zampini     const PetscReal b1[2] = {1.0, 0.0};
3061f80e275SEmil Constantinescu     PetscReal       binterpt[2][2];
307da80777bSKarl Rupp     binterpt[0][0] = 1.707106781186547524401 - 1.0;
308da80777bSKarl Rupp     binterpt[1][0] = 2.0 - 1.707106781186547524401;
309da80777bSKarl Rupp     binterpt[0][1] = 1.707106781186547524401 - 1.5;
310da80777bSKarl Rupp     binterpt[1][1] = 1.5 - 1.707106781186547524401;
311bbd56ea5SKarl Rupp 
3129566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSW2P, 2, 2, &A[0][0], &Gamma[0][0], b, b1, 2, &binterpt[0][0]));
313e27a552bSJed Brown   }
314e27a552bSJed Brown   {
315da80777bSKarl 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. */
316b16ce868SStefano Zampini     const PetscReal A[2][2] = {
3179371c9d4SSatish Balay       {0,  0},
3189371c9d4SSatish Balay       {1., 0}
319b16ce868SStefano Zampini     };
320b16ce868SStefano Zampini     const PetscReal Gamma[2][2] = {
321b16ce868SStefano Zampini       {0.2928932188134524755992,       0                       },
322b16ce868SStefano Zampini       {-2. * 0.2928932188134524755992, 0.2928932188134524755992}
323b16ce868SStefano Zampini     };
324b16ce868SStefano Zampini     const PetscReal b[2]  = {0.5, 0.5};
325b16ce868SStefano Zampini     const PetscReal b1[2] = {1.0, 0.0};
3261f80e275SEmil Constantinescu     PetscReal       binterpt[2][2];
327da80777bSKarl Rupp     binterpt[0][0] = 0.2928932188134524755992 - 1.0;
328da80777bSKarl Rupp     binterpt[1][0] = 2.0 - 0.2928932188134524755992;
329da80777bSKarl Rupp     binterpt[0][1] = 0.2928932188134524755992 - 1.5;
330da80777bSKarl Rupp     binterpt[1][1] = 1.5 - 0.2928932188134524755992;
331bbd56ea5SKarl Rupp 
3329566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSW2M, 2, 2, &A[0][0], &Gamma[0][0], b, b1, 2, &binterpt[0][0]));
333fe7e6d57SJed Brown   }
334fe7e6d57SJed Brown   {
335da80777bSKarl Rupp     /*const PetscReal g = 7.8867513459481287e-01; Directly written in-place below */
3361f80e275SEmil Constantinescu     PetscReal       binterpt[3][2];
337b16ce868SStefano Zampini     const PetscReal A[3][3] = {
3389371c9d4SSatish Balay       {0,                      0, 0},
339fe7e6d57SJed Brown       {1.5773502691896257e+00, 0, 0},
3409371c9d4SSatish Balay       {0.5,                    0, 0}
341b16ce868SStefano Zampini     };
342b16ce868SStefano Zampini     const PetscReal Gamma[3][3] = {
343b16ce868SStefano Zampini       {7.8867513459481287e-01,  0,                       0                     },
344b16ce868SStefano Zampini       {-1.5773502691896257e+00, 7.8867513459481287e-01,  0                     },
345b16ce868SStefano Zampini       {-6.7075317547305480e-01, -1.7075317547305482e-01, 7.8867513459481287e-01}
346b16ce868SStefano Zampini     };
347b16ce868SStefano Zampini     const PetscReal b[3]  = {1.0566243270259355e-01, 4.9038105676657971e-02, 8.4529946162074843e-01};
348b16ce868SStefano Zampini     const PetscReal b2[3] = {-1.7863279495408180e-01, 1. / 3., 8.4529946162074843e-01};
3491f80e275SEmil Constantinescu 
3501f80e275SEmil Constantinescu     binterpt[0][0] = -0.8094010767585034;
3511f80e275SEmil Constantinescu     binterpt[1][0] = -0.5;
3521f80e275SEmil Constantinescu     binterpt[2][0] = 2.3094010767585034;
3531f80e275SEmil Constantinescu     binterpt[0][1] = 0.9641016151377548;
3541f80e275SEmil Constantinescu     binterpt[1][1] = 0.5;
3551f80e275SEmil Constantinescu     binterpt[2][1] = -1.4641016151377548;
356bbd56ea5SKarl Rupp 
3579566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWRA3PW, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0]));
358fe7e6d57SJed Brown   }
359fe7e6d57SJed Brown   {
3603ca35412SEmil Constantinescu     PetscReal binterpt[4][3];
361da80777bSKarl Rupp     /*const PetscReal g = 4.3586652150845900e-01; Directly written in-place below */
362b16ce868SStefano Zampini     const PetscReal A[4][4] = {
3639371c9d4SSatish Balay       {0,                      0,                       0,  0},
364fe7e6d57SJed Brown       {8.7173304301691801e-01, 0,                       0,  0},
365fe7e6d57SJed Brown       {8.4457060015369423e-01, -1.1299064236484185e-01, 0,  0},
3669371c9d4SSatish Balay       {0,                      0,                       1., 0}
367b16ce868SStefano Zampini     };
368b16ce868SStefano Zampini     const PetscReal Gamma[4][4] = {
369b16ce868SStefano Zampini       {4.3586652150845900e-01,  0,                       0,                      0                     },
370b16ce868SStefano Zampini       {-8.7173304301691801e-01, 4.3586652150845900e-01,  0,                      0                     },
371b16ce868SStefano Zampini       {-9.0338057013044082e-01, 5.4180672388095326e-02,  4.3586652150845900e-01, 0                     },
372b16ce868SStefano Zampini       {2.4212380706095346e-01,  -1.2232505839045147e+00, 5.4526025533510214e-01, 4.3586652150845900e-01}
373b16ce868SStefano Zampini     };
374b16ce868SStefano Zampini     const PetscReal b[4]  = {2.4212380706095346e-01, -1.2232505839045147e+00, 1.5452602553351020e+00, 4.3586652150845900e-01};
375b16ce868SStefano Zampini     const PetscReal b2[4] = {3.7810903145819369e-01, -9.6042292212423178e-02, 5.0000000000000000e-01, 2.1793326075422950e-01};
3763ca35412SEmil Constantinescu 
3773ca35412SEmil Constantinescu     binterpt[0][0] = 1.0564298455794094;
3783ca35412SEmil Constantinescu     binterpt[1][0] = 2.296429974281067;
3793ca35412SEmil Constantinescu     binterpt[2][0] = -1.307599564525376;
3803ca35412SEmil Constantinescu     binterpt[3][0] = -1.045260255335102;
3813ca35412SEmil Constantinescu     binterpt[0][1] = -1.3864882699759573;
3823ca35412SEmil Constantinescu     binterpt[1][1] = -8.262611700275677;
3833ca35412SEmil Constantinescu     binterpt[2][1] = 7.250979895056055;
3843ca35412SEmil Constantinescu     binterpt[3][1] = 2.398120075195581;
3853ca35412SEmil Constantinescu     binterpt[0][2] = 0.5721822314575016;
3863ca35412SEmil Constantinescu     binterpt[1][2] = 4.742931142090097;
3873ca35412SEmil Constantinescu     binterpt[2][2] = -4.398120075195578;
3883ca35412SEmil Constantinescu     binterpt[3][2] = -0.9169932983520199;
3893ca35412SEmil Constantinescu 
3909566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWRA34PW2, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0]));
391e27a552bSJed Brown   }
392ef3c5b88SJed Brown   {
393da80777bSKarl Rupp     /* const PetscReal g = 0.5;       Directly written in-place below */
394b16ce868SStefano Zampini     const PetscReal A[4][4] = {
3959371c9d4SSatish Balay       {0,    0,     0,   0},
396ef3c5b88SJed Brown       {0,    0,     0,   0},
397ef3c5b88SJed Brown       {1.,   0,     0,   0},
3989371c9d4SSatish Balay       {0.75, -0.25, 0.5, 0}
399b16ce868SStefano Zampini     };
400b16ce868SStefano Zampini     const PetscReal Gamma[4][4] = {
401b16ce868SStefano Zampini       {0.5,     0,       0,       0  },
402b16ce868SStefano Zampini       {1.,      0.5,     0,       0  },
403b16ce868SStefano Zampini       {-0.25,   -0.25,   0.5,     0  },
404b16ce868SStefano Zampini       {1. / 12, 1. / 12, -2. / 3, 0.5}
405b16ce868SStefano Zampini     };
406b16ce868SStefano Zampini     const PetscReal b[4]  = {5. / 6, -1. / 6, -1. / 6, 0.5};
407b16ce868SStefano Zampini     const PetscReal b2[4] = {0.75, -0.25, 0.5, 0};
408bbd56ea5SKarl Rupp 
4099566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWRODAS3, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 0, NULL));
410ef3c5b88SJed Brown   }
411ef3c5b88SJed Brown   {
412da80777bSKarl Rupp     /*const PetscReal g = 0.43586652150845899941601945119356;       Directly written in-place below */
413b16ce868SStefano Zampini     const PetscReal A[3][3] = {
4149371c9d4SSatish Balay       {0,                                  0, 0},
415da80777bSKarl Rupp       {0.43586652150845899941601945119356, 0, 0},
4169371c9d4SSatish Balay       {0.43586652150845899941601945119356, 0, 0}
417b16ce868SStefano Zampini     };
418b16ce868SStefano Zampini     const PetscReal Gamma[3][3] = {
419b16ce868SStefano Zampini       {0.43586652150845899941601945119356,  0,                                  0                                 },
420b16ce868SStefano Zampini       {-0.19294655696029095575009695436041, 0.43586652150845899941601945119356, 0                                 },
421b16ce868SStefano Zampini       {0,                                   1.74927148125794685173529749738960, 0.43586652150845899941601945119356}
422b16ce868SStefano Zampini     };
423b16ce868SStefano Zampini     const PetscReal b[3]  = {-0.75457412385404315829818998646589, 1.94100407061964420292840123379419, -0.18642994676560104463021124732829};
424b16ce868SStefano Zampini     const PetscReal b2[3] = {-1.53358745784149585370766523913002, 2.81745131148625772213931745457622, -0.28386385364476186843165221544619};
4251f80e275SEmil Constantinescu 
4261f80e275SEmil Constantinescu     PetscReal binterpt[3][2];
4271f80e275SEmil Constantinescu     binterpt[0][0] = 3.793692883777660870425141387941;
4281f80e275SEmil Constantinescu     binterpt[1][0] = -2.918692883777660870425141387941;
4291f80e275SEmil Constantinescu     binterpt[2][0] = 0.125;
4301f80e275SEmil Constantinescu     binterpt[0][1] = -0.725741064379812106687651020584;
4311f80e275SEmil Constantinescu     binterpt[1][1] = 0.559074397713145440020984353917;
4321f80e275SEmil Constantinescu     binterpt[2][1] = 0.16666666666666666666666666666667;
4331f80e275SEmil Constantinescu 
4349566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWSANDU3, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0]));
435ef3c5b88SJed Brown   }
436b1c69cc3SEmil Constantinescu   {
437da80777bSKarl Rupp     /*const PetscReal s3 = PetscSqrtReal(3.),g = (3.0+s3)/6.0;
438da80777bSKarl Rupp      * Direct evaluation: s3 = 1.732050807568877293527;
439da80777bSKarl Rupp      *                     g = 0.7886751345948128822546;
440da80777bSKarl Rupp      * Values are directly inserted below to ensure availability at compile time (compiler warnings otherwise...) */
441b16ce868SStefano Zampini     const PetscReal A[3][3] = {
4429371c9d4SSatish Balay       {0,    0,    0},
443b1c69cc3SEmil Constantinescu       {1,    0,    0},
4449371c9d4SSatish Balay       {0.25, 0.25, 0}
445b16ce868SStefano Zampini     };
446b16ce868SStefano Zampini     const PetscReal Gamma[3][3] = {
447b16ce868SStefano Zampini       {0,                                       0,                                      0                       },
448b16ce868SStefano Zampini       {(-3.0 - 1.732050807568877293527) / 6.0,  0.7886751345948128822546,               0                       },
449b16ce868SStefano Zampini       {(-3.0 - 1.732050807568877293527) / 24.0, (-3.0 - 1.732050807568877293527) / 8.0, 0.7886751345948128822546}
450b16ce868SStefano Zampini     };
451b16ce868SStefano Zampini     const PetscReal b[3]  = {1. / 6., 1. / 6., 2. / 3.};
452b16ce868SStefano Zampini     const PetscReal b2[3] = {1. / 4., 1. / 4., 1. / 2.};
453c0cb691aSEmil Constantinescu     PetscReal       binterpt[3][2];
454da80777bSKarl Rupp 
455c0cb691aSEmil Constantinescu     binterpt[0][0] = 0.089316397477040902157517886164709;
456c0cb691aSEmil Constantinescu     binterpt[1][0] = -0.91068360252295909784248211383529;
457c0cb691aSEmil Constantinescu     binterpt[2][0] = 1.8213672050459181956849642276706;
458c0cb691aSEmil Constantinescu     binterpt[0][1] = 0.077350269189625764509148780501957;
459c0cb691aSEmil Constantinescu     binterpt[1][1] = 1.077350269189625764509148780502;
460c0cb691aSEmil Constantinescu     binterpt[2][1] = -1.1547005383792515290182975610039;
461bbd56ea5SKarl Rupp 
4629566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWASSP3P3S1C, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0]));
463b1c69cc3SEmil Constantinescu   }
464b1c69cc3SEmil Constantinescu 
465b1c69cc3SEmil Constantinescu   {
466b16ce868SStefano Zampini     const PetscReal A[4][4] = {
4679371c9d4SSatish Balay       {0,       0,       0,       0},
468b1c69cc3SEmil Constantinescu       {1. / 2., 0,       0,       0},
469b1c69cc3SEmil Constantinescu       {1. / 2., 1. / 2., 0,       0},
4709371c9d4SSatish Balay       {1. / 6., 1. / 6., 1. / 6., 0}
471b16ce868SStefano Zampini     };
472b16ce868SStefano Zampini     const PetscReal Gamma[4][4] = {
473b16ce868SStefano Zampini       {1. / 2., 0,        0,       0},
474b16ce868SStefano Zampini       {0.0,     1. / 4.,  0,       0},
475b16ce868SStefano Zampini       {-2.,     -2. / 3., 2. / 3., 0},
476b16ce868SStefano Zampini       {1. / 2., 5. / 36., -2. / 9, 0}
477b16ce868SStefano Zampini     };
478b16ce868SStefano Zampini     const PetscReal b[4]  = {1. / 6., 1. / 6., 1. / 6., 1. / 2.};
479b16ce868SStefano Zampini     const PetscReal b2[4] = {1. / 8., 3. / 4., 1. / 8., 0};
480c0cb691aSEmil Constantinescu     PetscReal       binterpt[4][3];
481da80777bSKarl Rupp 
482c0cb691aSEmil Constantinescu     binterpt[0][0] = 6.25;
483c0cb691aSEmil Constantinescu     binterpt[1][0] = -30.25;
484c0cb691aSEmil Constantinescu     binterpt[2][0] = 1.75;
485c0cb691aSEmil Constantinescu     binterpt[3][0] = 23.25;
486c0cb691aSEmil Constantinescu     binterpt[0][1] = -9.75;
487c0cb691aSEmil Constantinescu     binterpt[1][1] = 58.75;
488c0cb691aSEmil Constantinescu     binterpt[2][1] = -3.25;
489c0cb691aSEmil Constantinescu     binterpt[3][1] = -45.75;
490c0cb691aSEmil Constantinescu     binterpt[0][2] = 3.6666666666666666666666666666667;
491c0cb691aSEmil Constantinescu     binterpt[1][2] = -28.333333333333333333333333333333;
492c0cb691aSEmil Constantinescu     binterpt[2][2] = 1.6666666666666666666666666666667;
493c0cb691aSEmil Constantinescu     binterpt[3][2] = 23.;
494bbd56ea5SKarl Rupp 
4959566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWLASSP3P4S2C, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0]));
496b1c69cc3SEmil Constantinescu   }
497b1c69cc3SEmil Constantinescu 
498b1c69cc3SEmil Constantinescu   {
499b16ce868SStefano Zampini     const PetscReal A[4][4] = {
5009371c9d4SSatish Balay       {0,       0,       0,       0},
501b1c69cc3SEmil Constantinescu       {1. / 2., 0,       0,       0},
502b1c69cc3SEmil Constantinescu       {1. / 2., 1. / 2., 0,       0},
5039371c9d4SSatish Balay       {1. / 6., 1. / 6., 1. / 6., 0}
504b16ce868SStefano Zampini     };
505b16ce868SStefano Zampini     const PetscReal Gamma[4][4] = {
506b16ce868SStefano Zampini       {1. / 2.,  0,          0,        0},
507b16ce868SStefano Zampini       {0.0,      3. / 4.,    0,        0},
508b16ce868SStefano Zampini       {-2. / 3., -23. / 9.,  2. / 9.,  0},
509b16ce868SStefano Zampini       {1. / 18., 65. / 108., -2. / 27, 0}
510b16ce868SStefano Zampini     };
511b16ce868SStefano Zampini     const PetscReal b[4]  = {1. / 6., 1. / 6., 1. / 6., 1. / 2.};
512b16ce868SStefano Zampini     const PetscReal b2[4] = {3. / 16., 10. / 16., 3. / 16., 0};
513c0cb691aSEmil Constantinescu     PetscReal       binterpt[4][3];
514da80777bSKarl Rupp 
515c0cb691aSEmil Constantinescu     binterpt[0][0] = 1.6911764705882352941176470588235;
516c0cb691aSEmil Constantinescu     binterpt[1][0] = 3.6813725490196078431372549019608;
517c0cb691aSEmil Constantinescu     binterpt[2][0] = 0.23039215686274509803921568627451;
518c0cb691aSEmil Constantinescu     binterpt[3][0] = -4.6029411764705882352941176470588;
519c0cb691aSEmil Constantinescu     binterpt[0][1] = -0.95588235294117647058823529411765;
520c0cb691aSEmil Constantinescu     binterpt[1][1] = -6.2401960784313725490196078431373;
521c0cb691aSEmil Constantinescu     binterpt[2][1] = -0.31862745098039215686274509803922;
522c0cb691aSEmil Constantinescu     binterpt[3][1] = 7.5147058823529411764705882352941;
523c0cb691aSEmil Constantinescu     binterpt[0][2] = -0.56862745098039215686274509803922;
524c0cb691aSEmil Constantinescu     binterpt[1][2] = 2.7254901960784313725490196078431;
525c0cb691aSEmil Constantinescu     binterpt[2][2] = 0.25490196078431372549019607843137;
526c0cb691aSEmil Constantinescu     binterpt[3][2] = -2.4117647058823529411764705882353;
527bbd56ea5SKarl Rupp 
5289566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWLLSSP3P4S2C, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0]));
529b1c69cc3SEmil Constantinescu   }
530753f8adbSEmil Constantinescu 
531753f8adbSEmil Constantinescu   {
532753f8adbSEmil Constantinescu     PetscReal A[4][4], Gamma[4][4], b[4], b2[4];
5333ca35412SEmil Constantinescu     PetscReal binterpt[4][3];
534753f8adbSEmil Constantinescu 
535753f8adbSEmil Constantinescu     Gamma[0][0] = 0.4358665215084589994160194475295062513822671686978816;
5369371c9d4SSatish Balay     Gamma[0][1] = 0;
5379371c9d4SSatish Balay     Gamma[0][2] = 0;
5389371c9d4SSatish Balay     Gamma[0][3] = 0;
539753f8adbSEmil Constantinescu     Gamma[1][0] = -1.997527830934941248426324674704153457289527280554476;
540753f8adbSEmil Constantinescu     Gamma[1][1] = 0.4358665215084589994160194475295062513822671686978816;
5419371c9d4SSatish Balay     Gamma[1][2] = 0;
5429371c9d4SSatish Balay     Gamma[1][3] = 0;
543753f8adbSEmil Constantinescu     Gamma[2][0] = -1.007948511795029620852002345345404191008352770119903;
544753f8adbSEmil Constantinescu     Gamma[2][1] = -0.004648958462629345562774289390054679806993396798458131;
545753f8adbSEmil Constantinescu     Gamma[2][2] = 0.4358665215084589994160194475295062513822671686978816;
54605e8e825SJed Brown     Gamma[2][3] = 0;
547753f8adbSEmil Constantinescu     Gamma[3][0] = -0.6685429734233467180451604600279552604364311322650783;
548753f8adbSEmil Constantinescu     Gamma[3][1] = 0.6056625986449338476089525334450053439525178740492984;
549753f8adbSEmil Constantinescu     Gamma[3][2] = -0.9717899277217721234705114616271378792182450260943198;
550753f8adbSEmil Constantinescu     Gamma[3][3] = 0;
551753f8adbSEmil Constantinescu 
5529371c9d4SSatish Balay     A[0][0] = 0;
5539371c9d4SSatish Balay     A[0][1] = 0;
5549371c9d4SSatish Balay     A[0][2] = 0;
5559371c9d4SSatish Balay     A[0][3] = 0;
556753f8adbSEmil Constantinescu     A[1][0] = 0.8717330430169179988320388950590125027645343373957631;
5579371c9d4SSatish Balay     A[1][1] = 0;
5589371c9d4SSatish Balay     A[1][2] = 0;
5599371c9d4SSatish Balay     A[1][3] = 0;
560753f8adbSEmil Constantinescu     A[2][0] = 0.5275890119763004115618079766722914408876108660811028;
561753f8adbSEmil Constantinescu     A[2][1] = 0.07241098802369958843819203208518599088698057726988732;
5629371c9d4SSatish Balay     A[2][2] = 0;
5639371c9d4SSatish Balay     A[2][3] = 0;
564753f8adbSEmil Constantinescu     A[3][0] = 0.3990960076760701320627260685975778145384666450351314;
565753f8adbSEmil Constantinescu     A[3][1] = -0.4375576546135194437228463747348862825846903771419953;
566753f8adbSEmil Constantinescu     A[3][2] = 1.038461646937449311660120300601880176655352737312713;
56705e8e825SJed Brown     A[3][3] = 0;
568753f8adbSEmil Constantinescu 
569753f8adbSEmil Constantinescu     b[0] = 0.1876410243467238251612921333138006734899663569186926;
570753f8adbSEmil Constantinescu     b[1] = -0.5952974735769549480478230473706443582188442040780541;
571753f8adbSEmil Constantinescu     b[2] = 0.9717899277217721234705114616271378792182450260943198;
572753f8adbSEmil Constantinescu     b[3] = 0.4358665215084589994160194475295062513822671686978816;
573753f8adbSEmil Constantinescu 
574753f8adbSEmil Constantinescu     b2[0] = 0.2147402862233891404862383521089097657790734483804460;
575753f8adbSEmil Constantinescu     b2[1] = -0.4851622638849390928209050538171743017757490232519684;
576753f8adbSEmil Constantinescu     b2[2] = 0.8687250025203875511662123688667549217531982787600080;
577753f8adbSEmil Constantinescu     b2[3] = 0.4016969751411624011684543450940068201770721128357014;
578753f8adbSEmil Constantinescu 
5793ca35412SEmil Constantinescu     binterpt[0][0] = 2.2565812720167954547104627844105;
5803ca35412SEmil Constantinescu     binterpt[1][0] = 1.349166413351089573796243820819;
5813ca35412SEmil Constantinescu     binterpt[2][0] = -2.4695174540533503758652847586647;
5823ca35412SEmil Constantinescu     binterpt[3][0] = -0.13623023131453465264142184656474;
5833ca35412SEmil Constantinescu     binterpt[0][1] = -3.0826699111559187902922463354557;
5843ca35412SEmil Constantinescu     binterpt[1][1] = -2.4689115685996042534544925650515;
5853ca35412SEmil Constantinescu     binterpt[2][1] = 5.7428279814696677152129332773553;
5863ca35412SEmil Constantinescu     binterpt[3][1] = -0.19124650171414467146619437684812;
5873ca35412SEmil Constantinescu     binterpt[0][2] = 1.0137296634858471607430756831148;
5883ca35412SEmil Constantinescu     binterpt[1][2] = 0.52444768167155973161042570784064;
5893ca35412SEmil Constantinescu     binterpt[2][2] = -2.3015205996945452158771370439586;
5903ca35412SEmil Constantinescu     binterpt[3][2] = 0.76334325453713832352363565300308;
591f4aed992SEmil Constantinescu 
5929566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWARK3, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0]));
593753f8adbSEmil Constantinescu   }
5949566063dSJacob Faibussowitsch   PetscCall(TSRosWRegisterRos4(TSROSWGRK4T, 0.231, PETSC_DEFAULT, PETSC_DEFAULT, 0, -0.1282612945269037e+01));
5959566063dSJacob Faibussowitsch   PetscCall(TSRosWRegisterRos4(TSROSWSHAMP4, 0.5, PETSC_DEFAULT, PETSC_DEFAULT, 0, 125. / 108.));
5969566063dSJacob Faibussowitsch   PetscCall(TSRosWRegisterRos4(TSROSWVELDD4, 0.22570811482256823492, PETSC_DEFAULT, PETSC_DEFAULT, 0, -1.355958941201148));
5979566063dSJacob Faibussowitsch   PetscCall(TSRosWRegisterRos4(TSROSW4L, 0.57282, PETSC_DEFAULT, PETSC_DEFAULT, 0, -1.093502252409163));
5983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
599e27a552bSJed Brown }
600e27a552bSJed Brown 
601e27a552bSJed Brown /*@C
602bcf0153eSBarry Smith   TSRosWRegisterDestroy - Frees the list of schemes that were registered by `TSRosWRegister()`.
603e27a552bSJed Brown 
604e27a552bSJed Brown   Not Collective
605e27a552bSJed Brown 
606e27a552bSJed Brown   Level: advanced
607e27a552bSJed Brown 
6081cc06b55SBarry Smith .seealso: [](ch_ts), `TSRosWRegister()`, `TSRosWRegisterAll()`
609e27a552bSJed Brown @*/
610d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWRegisterDestroy(void)
611d71ae5a4SJacob Faibussowitsch {
61261692a83SJed Brown   RosWTableauLink link;
613e27a552bSJed Brown 
614e27a552bSJed Brown   PetscFunctionBegin;
61561692a83SJed Brown   while ((link = RosWTableauList)) {
61661692a83SJed Brown     RosWTableau t   = &link->tab;
61761692a83SJed Brown     RosWTableauList = link->next;
6189566063dSJacob Faibussowitsch     PetscCall(PetscFree5(t->A, t->Gamma, t->b, t->ASum, t->GammaSum));
6199566063dSJacob Faibussowitsch     PetscCall(PetscFree5(t->At, t->bt, t->GammaInv, t->GammaZeroDiag, t->GammaExplicitCorr));
6209566063dSJacob Faibussowitsch     PetscCall(PetscFree2(t->bembed, t->bembedt));
6219566063dSJacob Faibussowitsch     PetscCall(PetscFree(t->binterpt));
6229566063dSJacob Faibussowitsch     PetscCall(PetscFree(t->name));
6239566063dSJacob Faibussowitsch     PetscCall(PetscFree(link));
624e27a552bSJed Brown   }
625e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_FALSE;
6263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
627e27a552bSJed Brown }
628e27a552bSJed Brown 
629e27a552bSJed Brown /*@C
630bcf0153eSBarry Smith   TSRosWInitializePackage - This function initializes everything in the `TSROSW` package. It is called
631bcf0153eSBarry Smith   from `TSInitializePackage()`.
632e27a552bSJed Brown 
633e27a552bSJed Brown   Level: developer
634e27a552bSJed Brown 
6351cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `PetscInitialize()`, `TSRosWFinalizePackage()`
636e27a552bSJed Brown @*/
637d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWInitializePackage(void)
638d71ae5a4SJacob Faibussowitsch {
639e27a552bSJed Brown   PetscFunctionBegin;
6403ba16761SJacob Faibussowitsch   if (TSRosWPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
641e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_TRUE;
6429566063dSJacob Faibussowitsch   PetscCall(TSRosWRegisterAll());
6439566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(TSRosWFinalizePackage));
6443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
645e27a552bSJed Brown }
646e27a552bSJed Brown 
647e27a552bSJed Brown /*@C
648bcf0153eSBarry Smith   TSRosWFinalizePackage - This function destroys everything in the `TSROSW` package. It is
649bcf0153eSBarry Smith   called from `PetscFinalize()`.
650e27a552bSJed Brown 
651e27a552bSJed Brown   Level: developer
652e27a552bSJed Brown 
6531cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`, `PetscFinalize()`, `TSRosWInitializePackage()`
654e27a552bSJed Brown @*/
655d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWFinalizePackage(void)
656d71ae5a4SJacob Faibussowitsch {
657e27a552bSJed Brown   PetscFunctionBegin;
658e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_FALSE;
6599566063dSJacob Faibussowitsch   PetscCall(TSRosWRegisterDestroy());
6603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
661e27a552bSJed Brown }
662e27a552bSJed Brown 
663e27a552bSJed Brown /*@C
664bcf0153eSBarry Smith   TSRosWRegister - register a `TSROSW`, Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
665e27a552bSJed Brown 
666e27a552bSJed Brown   Not Collective, but the same schemes should be registered on all processes on which they will be used
667e27a552bSJed Brown 
668e27a552bSJed Brown   Input Parameters:
669e27a552bSJed Brown + name     - identifier for method
670e27a552bSJed Brown . order    - approximation order of method
671e27a552bSJed Brown . s        - number of stages, this is the dimension of the matrices below
67261692a83SJed Brown . A        - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular
67361692a83SJed Brown . Gamma    - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal
674fe7e6d57SJed Brown . b        - Step completion table (dimension s)
6750298fd71SBarry Smith . bembed   - Step completion table for a scheme of order one less (dimension s, NULL if no embedded scheme is available)
676f4aed992SEmil Constantinescu . pinterp  - Order of the interpolation scheme, equal to the number of columns of binterpt
67742faf41dSJed Brown - binterpt - Coefficients of the interpolation formula (dimension s*pinterp)
678e27a552bSJed Brown 
679e27a552bSJed Brown   Level: advanced
680e27a552bSJed Brown 
681bcf0153eSBarry Smith   Note:
682bcf0153eSBarry Smith   Several Rosenbrock W methods are provided, this function is only needed to create new methods.
683bcf0153eSBarry Smith 
6841cc06b55SBarry Smith .seealso: [](ch_ts), `TSROSW`
685e27a552bSJed Brown @*/
686d71ae5a4SJacob 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[])
687d71ae5a4SJacob Faibussowitsch {
68861692a83SJed Brown   RosWTableauLink link;
68961692a83SJed Brown   RosWTableau     t;
69061692a83SJed Brown   PetscInt        i, j, k;
69161692a83SJed Brown   PetscScalar    *GammaInv;
692e27a552bSJed Brown 
693e27a552bSJed Brown   PetscFunctionBegin;
6944f572ea9SToby Isaac   PetscAssertPointer(name, 1);
6954f572ea9SToby Isaac   PetscAssertPointer(A, 4);
6964f572ea9SToby Isaac   PetscAssertPointer(Gamma, 5);
6974f572ea9SToby Isaac   PetscAssertPointer(b, 6);
6984f572ea9SToby Isaac   if (bembed) PetscAssertPointer(bembed, 7);
699fe7e6d57SJed Brown 
7009566063dSJacob Faibussowitsch   PetscCall(TSRosWInitializePackage());
7019566063dSJacob Faibussowitsch   PetscCall(PetscNew(&link));
702e27a552bSJed Brown   t = &link->tab;
7039566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &t->name));
704e27a552bSJed Brown   t->order = order;
705e27a552bSJed Brown   t->s     = s;
7069566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(s * s, &t->A, s * s, &t->Gamma, s, &t->b, s, &t->ASum, s, &t->GammaSum));
7079566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(s * s, &t->At, s, &t->bt, s * s, &t->GammaInv, s, &t->GammaZeroDiag, s * s, &t->GammaExplicitCorr));
7089566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->A, A, s * s));
7099566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->Gamma, Gamma, s * s));
7109566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->GammaExplicitCorr, Gamma, s * s));
7119566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->b, b, s));
712fe7e6d57SJed Brown   if (bembed) {
7139566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(s, &t->bembed, s, &t->bembedt));
7149566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(t->bembed, bembed, s));
715fe7e6d57SJed Brown   }
71661692a83SJed Brown   for (i = 0; i < s; i++) {
71761692a83SJed Brown     t->ASum[i]     = 0;
71861692a83SJed Brown     t->GammaSum[i] = 0;
71961692a83SJed Brown     for (j = 0; j < s; j++) {
72061692a83SJed Brown       t->ASum[i] += A[i * s + j];
721fe7e6d57SJed Brown       t->GammaSum[i] += Gamma[i * s + j];
72261692a83SJed Brown     }
72361692a83SJed Brown   }
7249566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s * s, &GammaInv)); /* Need to use Scalar for inverse, then convert back to Real */
72561692a83SJed Brown   for (i = 0; i < s * s; i++) GammaInv[i] = Gamma[i];
726fd96d5b0SEmil Constantinescu   for (i = 0; i < s; i++) {
727fd96d5b0SEmil Constantinescu     if (Gamma[i * s + i] == 0.0) {
728fd96d5b0SEmil Constantinescu       GammaInv[i * s + i] = 1.0;
729c17803e7SJed Brown       t->GammaZeroDiag[i] = PETSC_TRUE;
730fd96d5b0SEmil Constantinescu     } else {
731c17803e7SJed Brown       t->GammaZeroDiag[i] = PETSC_FALSE;
732fd96d5b0SEmil Constantinescu     }
733fd96d5b0SEmil Constantinescu   }
734fd96d5b0SEmil Constantinescu 
73561692a83SJed Brown   switch (s) {
736d71ae5a4SJacob Faibussowitsch   case 1:
737d71ae5a4SJacob Faibussowitsch     GammaInv[0] = 1. / GammaInv[0];
738d71ae5a4SJacob Faibussowitsch     break;
739d71ae5a4SJacob Faibussowitsch   case 2:
740d71ae5a4SJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_2(GammaInv, 0, PETSC_FALSE, NULL));
741d71ae5a4SJacob Faibussowitsch     break;
742d71ae5a4SJacob Faibussowitsch   case 3:
743d71ae5a4SJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_3(GammaInv, 0, PETSC_FALSE, NULL));
744d71ae5a4SJacob Faibussowitsch     break;
745d71ae5a4SJacob Faibussowitsch   case 4:
746d71ae5a4SJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_4(GammaInv, 0, PETSC_FALSE, NULL));
747d71ae5a4SJacob Faibussowitsch     break;
74861692a83SJed Brown   case 5: {
74961692a83SJed Brown     PetscInt  ipvt5[5];
75061692a83SJed Brown     MatScalar work5[5 * 5];
7519371c9d4SSatish Balay     PetscCall(PetscKernel_A_gets_inverse_A_5(GammaInv, ipvt5, work5, 0, PETSC_FALSE, NULL));
7529371c9d4SSatish Balay     break;
75361692a83SJed Brown   }
754d71ae5a4SJacob Faibussowitsch   case 6:
755d71ae5a4SJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_6(GammaInv, 0, PETSC_FALSE, NULL));
756d71ae5a4SJacob Faibussowitsch     break;
757d71ae5a4SJacob Faibussowitsch   case 7:
758d71ae5a4SJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_7(GammaInv, 0, PETSC_FALSE, NULL));
759d71ae5a4SJacob Faibussowitsch     break;
760d71ae5a4SJacob Faibussowitsch   default:
761d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not implemented for %" PetscInt_FMT " stages", s);
76261692a83SJed Brown   }
76361692a83SJed Brown   for (i = 0; i < s * s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]);
7649566063dSJacob Faibussowitsch   PetscCall(PetscFree(GammaInv));
76543b21953SEmil Constantinescu 
76643b21953SEmil Constantinescu   for (i = 0; i < s; i++) {
76743b21953SEmil Constantinescu     for (k = 0; k < i + 1; k++) {
76843b21953SEmil Constantinescu       t->GammaExplicitCorr[i * s + k] = (t->GammaExplicitCorr[i * s + k]) * (t->GammaInv[k * s + k]);
769ad540459SPierre Jolivet       for (j = k + 1; j < i + 1; j++) t->GammaExplicitCorr[i * s + k] += (t->GammaExplicitCorr[i * s + j]) * (t->GammaInv[j * s + k]);
77043b21953SEmil Constantinescu     }
77143b21953SEmil Constantinescu   }
77243b21953SEmil Constantinescu 
77361692a83SJed Brown   for (i = 0; i < s; i++) {
77461692a83SJed Brown     for (j = 0; j < s; j++) {
77561692a83SJed Brown       t->At[i * s + j] = 0;
776ad540459SPierre Jolivet       for (k = 0; k < s; k++) t->At[i * s + j] += t->A[i * s + k] * t->GammaInv[k * s + j];
77761692a83SJed Brown     }
77861692a83SJed Brown     t->bt[i] = 0;
779ad540459SPierre Jolivet     for (j = 0; j < s; j++) t->bt[i] += t->b[j] * t->GammaInv[j * s + i];
780fe7e6d57SJed Brown     if (bembed) {
781fe7e6d57SJed Brown       t->bembedt[i] = 0;
782ad540459SPierre Jolivet       for (j = 0; j < s; j++) t->bembedt[i] += t->bembed[j] * t->GammaInv[j * s + i];
783fe7e6d57SJed Brown     }
78461692a83SJed Brown   }
7858d59e960SJed Brown   t->ccfl = 1.0; /* Fix this */
7868d59e960SJed Brown 
787f4aed992SEmil Constantinescu   t->pinterp = pinterp;
7889566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s * pinterp, &t->binterpt));
7899566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->binterpt, binterpt, s * pinterp));
79061692a83SJed Brown   link->next      = RosWTableauList;
79161692a83SJed Brown   RosWTableauList = link;
7923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
793e27a552bSJed Brown }
794e27a552bSJed Brown 
79542faf41dSJed Brown /*@C
796fd292e60Sprj-   TSRosWRegisterRos4 - register a fourth order Rosenbrock scheme by providing parameter choices
79742faf41dSJed Brown 
79842faf41dSJed Brown   Not Collective, but the same schemes should be registered on all processes on which they will be used
79942faf41dSJed Brown 
80042faf41dSJed Brown   Input Parameters:
80142faf41dSJed Brown + name  - identifier for method
80242faf41dSJed Brown . gamma - leading coefficient (diagonal entry)
80342faf41dSJed Brown . a2    - design parameter, see Table 7.2 of Hairer&Wanner
80442faf41dSJed Brown . a3    - design parameter or PETSC_DEFAULT to satisfy one of the order five conditions (Eq 7.22)
80542faf41dSJed Brown . b3    - design parameter, see Table 7.2 of Hairer&Wanner
806a2b725a8SWilliam Gropp - e4    - design parameter for embedded method, see coefficient E4 in ros4.f code from Hairer
80742faf41dSJed Brown 
808bcf0153eSBarry Smith   Level: developer
809bcf0153eSBarry Smith 
81042faf41dSJed Brown   Notes:
81142faf41dSJed Brown   This routine encodes the design of fourth order Rosenbrock methods as described in Hairer and Wanner volume 2.
81242faf41dSJed Brown   It is used here to implement several methods from the book and can be used to experiment with new methods.
81342faf41dSJed Brown   It was written this way instead of by copying coefficients in order to provide better than double precision satisfaction of the order conditions.
81442faf41dSJed Brown 
8151cc06b55SBarry Smith .seealso: [](ch_ts), `TSRosW`, `TSRosWRegister()`
81642faf41dSJed Brown @*/
817d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWRegisterRos4(TSRosWType name, PetscReal gamma, PetscReal a2, PetscReal a3, PetscReal b3, PetscReal e4)
818d71ae5a4SJacob Faibussowitsch {
81942faf41dSJed Brown   /* Declare numeric constants so they can be quad precision without being truncated at double */
8209371c9d4SSatish 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;
82142faf41dSJed Brown   PetscReal   a4, a32, a42, a43, b1, b2, b4, beta2p, beta3p, beta4p, beta32, beta42, beta43, beta32beta2p, beta4jbetajp;
82242faf41dSJed Brown   PetscReal   A[4][4], Gamma[4][4], b[4], bm[4];
82342faf41dSJed Brown   PetscScalar M[3][3], rhs[3];
82442faf41dSJed Brown 
82542faf41dSJed Brown   PetscFunctionBegin;
82642faf41dSJed Brown   /* Step 1: choose Gamma (input) */
82742faf41dSJed Brown   /* Step 2: choose a2,a3,a4; b1,b2,b3,b4 to satisfy order conditions */
82813bcc0bdSJacob Faibussowitsch   if (a3 == (PetscReal)PETSC_DEFAULT) a3 = (one / five - a2 / four) / (one / four - a2 / three); /* Eq 7.22 */
82942faf41dSJed Brown   a4 = a3;                                                                                       /* consequence of 7.20 */
83042faf41dSJed Brown 
83142faf41dSJed Brown   /* Solve order conditions 7.15a, 7.15c, 7.15e */
8329371c9d4SSatish Balay   M[0][0] = one;
8339371c9d4SSatish Balay   M[0][1] = one;
8349371c9d4SSatish Balay   M[0][2] = one; /* 7.15a */
8359371c9d4SSatish Balay   M[1][0] = 0.0;
8369371c9d4SSatish Balay   M[1][1] = a2 * a2;
8379371c9d4SSatish Balay   M[1][2] = a4 * a4; /* 7.15c */
8389371c9d4SSatish Balay   M[2][0] = 0.0;
8399371c9d4SSatish Balay   M[2][1] = a2 * a2 * a2;
8409371c9d4SSatish Balay   M[2][2] = a4 * a4 * a4; /* 7.15e */
84142faf41dSJed Brown   rhs[0]  = one - b3;
84242faf41dSJed Brown   rhs[1]  = one / three - a3 * a3 * b3;
84342faf41dSJed Brown   rhs[2]  = one / four - a3 * a3 * a3 * b3;
8449566063dSJacob Faibussowitsch   PetscCall(PetscKernel_A_gets_inverse_A_3(&M[0][0], 0, PETSC_FALSE, NULL));
84542faf41dSJed Brown   b1 = PetscRealPart(M[0][0] * rhs[0] + M[0][1] * rhs[1] + M[0][2] * rhs[2]);
84642faf41dSJed Brown   b2 = PetscRealPart(M[1][0] * rhs[0] + M[1][1] * rhs[1] + M[1][2] * rhs[2]);
84742faf41dSJed Brown   b4 = PetscRealPart(M[2][0] * rhs[0] + M[2][1] * rhs[1] + M[2][2] * rhs[2]);
84842faf41dSJed Brown 
84942faf41dSJed Brown   /* Step 3 */
85042faf41dSJed Brown   beta43       = (p56 - a2 * p43) / (b4 * a3 * a3 * (a3 - a2)); /* 7.21 */
85142faf41dSJed Brown   beta32beta2p = p44 / (b4 * beta43);                           /* 7.15h */
85242faf41dSJed Brown   beta4jbetajp = (p32 - b3 * beta32beta2p) / b4;
8539371c9d4SSatish Balay   M[0][0]      = b2;
8549371c9d4SSatish Balay   M[0][1]      = b3;
8559371c9d4SSatish Balay   M[0][2]      = b4;
8569371c9d4SSatish Balay   M[1][0]      = a4 * a4 * beta32beta2p - a3 * a3 * beta4jbetajp;
8579371c9d4SSatish Balay   M[1][1]      = a2 * a2 * beta4jbetajp;
8589371c9d4SSatish Balay   M[1][2]      = -a2 * a2 * beta32beta2p;
8599371c9d4SSatish Balay   M[2][0]      = b4 * beta43 * a3 * a3 - p43;
8609371c9d4SSatish Balay   M[2][1]      = -b4 * beta43 * a2 * a2;
8619371c9d4SSatish Balay   M[2][2]      = 0;
8629371c9d4SSatish Balay   rhs[0]       = one / two - gamma;
8639371c9d4SSatish Balay   rhs[1]       = 0;
8649371c9d4SSatish Balay   rhs[2]       = -a2 * a2 * p32;
8659566063dSJacob Faibussowitsch   PetscCall(PetscKernel_A_gets_inverse_A_3(&M[0][0], 0, PETSC_FALSE, NULL));
86642faf41dSJed Brown   beta2p = PetscRealPart(M[0][0] * rhs[0] + M[0][1] * rhs[1] + M[0][2] * rhs[2]);
86742faf41dSJed Brown   beta3p = PetscRealPart(M[1][0] * rhs[0] + M[1][1] * rhs[1] + M[1][2] * rhs[2]);
86842faf41dSJed Brown   beta4p = PetscRealPart(M[2][0] * rhs[0] + M[2][1] * rhs[1] + M[2][2] * rhs[2]);
86942faf41dSJed Brown 
87042faf41dSJed Brown   /* Step 4: back-substitute */
87142faf41dSJed Brown   beta32 = beta32beta2p / beta2p;
87242faf41dSJed Brown   beta42 = (beta4jbetajp - beta43 * beta3p) / beta2p;
87342faf41dSJed Brown 
87442faf41dSJed Brown   /* Step 5: 7.15f and 7.20, then 7.16 */
87542faf41dSJed Brown   a43 = 0;
87642faf41dSJed Brown   a32 = p42 / (b3 * a3 * beta2p + b4 * a4 * beta2p);
87742faf41dSJed Brown   a42 = a32;
87842faf41dSJed Brown 
8799371c9d4SSatish Balay   A[0][0]     = 0;
8809371c9d4SSatish Balay   A[0][1]     = 0;
8819371c9d4SSatish Balay   A[0][2]     = 0;
8829371c9d4SSatish Balay   A[0][3]     = 0;
8839371c9d4SSatish Balay   A[1][0]     = a2;
8849371c9d4SSatish Balay   A[1][1]     = 0;
8859371c9d4SSatish Balay   A[1][2]     = 0;
8869371c9d4SSatish Balay   A[1][3]     = 0;
8879371c9d4SSatish Balay   A[2][0]     = a3 - a32;
8889371c9d4SSatish Balay   A[2][1]     = a32;
8899371c9d4SSatish Balay   A[2][2]     = 0;
8909371c9d4SSatish Balay   A[2][3]     = 0;
8919371c9d4SSatish Balay   A[3][0]     = a4 - a43 - a42;
8929371c9d4SSatish Balay   A[3][1]     = a42;
8939371c9d4SSatish Balay   A[3][2]     = a43;
8949371c9d4SSatish Balay   A[3][3]     = 0;
8959371c9d4SSatish Balay   Gamma[0][0] = gamma;
8969371c9d4SSatish Balay   Gamma[0][1] = 0;
8979371c9d4SSatish Balay   Gamma[0][2] = 0;
8989371c9d4SSatish Balay   Gamma[0][3] = 0;
8999371c9d4SSatish Balay   Gamma[1][0] = beta2p - A[1][0];
9009371c9d4SSatish Balay   Gamma[1][1] = gamma;
9019371c9d4SSatish Balay   Gamma[1][2] = 0;
9029371c9d4SSatish Balay   Gamma[1][3] = 0;
9039371c9d4SSatish Balay   Gamma[2][0] = beta3p - beta32 - A[2][0];
9049371c9d4SSatish Balay   Gamma[2][1] = beta32 - A[2][1];
9059371c9d4SSatish Balay   Gamma[2][2] = gamma;
9069371c9d4SSatish Balay   Gamma[2][3] = 0;
9079371c9d4SSatish Balay   Gamma[3][0] = beta4p - beta42 - beta43 - A[3][0];
9089371c9d4SSatish Balay   Gamma[3][1] = beta42 - A[3][1];
9099371c9d4SSatish Balay   Gamma[3][2] = beta43 - A[3][2];
9109371c9d4SSatish Balay   Gamma[3][3] = gamma;
9119371c9d4SSatish Balay   b[0]        = b1;
9129371c9d4SSatish Balay   b[1]        = b2;
9139371c9d4SSatish Balay   b[2]        = b3;
9149371c9d4SSatish Balay   b[3]        = b4;
91542faf41dSJed Brown 
91642faf41dSJed Brown   /* Construct embedded formula using given e4. We are solving Equation 7.18. */
91742faf41dSJed Brown   bm[3] = b[3] - e4 * gamma;                                              /* using definition of E4 */
91842faf41dSJed Brown   bm[2] = (p32 - beta4jbetajp * bm[3]) / (beta32 * beta2p);               /* fourth row of 7.18 */
91942faf41dSJed Brown   bm[1] = (one / two - gamma - beta3p * bm[2] - beta4p * bm[3]) / beta2p; /* second row */
92042faf41dSJed Brown   bm[0] = one - bm[1] - bm[2] - bm[3];                                    /* first row */
92142faf41dSJed Brown 
92242faf41dSJed Brown   {
92342faf41dSJed Brown     const PetscReal misfit = a2 * a2 * bm[1] + a3 * a3 * bm[2] + a4 * a4 * bm[3] - one / three;
9243c633725SBarry Smith     PetscCheck(PetscAbs(misfit) <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_SUP, "Assumptions violated, could not construct a third order embedded method");
92542faf41dSJed Brown   }
9269566063dSJacob Faibussowitsch   PetscCall(TSRosWRegister(name, 4, 4, &A[0][0], &Gamma[0][0], b, bm, 0, NULL));
9273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
92842faf41dSJed Brown }
92942faf41dSJed Brown 
9301c3436cfSJed Brown /*
9311c3436cfSJed Brown  The step completion formula is
9321c3436cfSJed Brown 
9331c3436cfSJed Brown  x1 = x0 + b^T Y
9341c3436cfSJed Brown 
9351c3436cfSJed Brown  where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been
9361c3436cfSJed Brown  updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write
9371c3436cfSJed Brown 
9381c3436cfSJed Brown  x1e = x0 + be^T Y
9391c3436cfSJed Brown      = x1 - b^T Y + be^T Y
9401c3436cfSJed Brown      = x1 + (be - b)^T Y
9411c3436cfSJed Brown 
9421c3436cfSJed Brown  so we can evaluate the method of different order even after the step has been optimistically completed.
9431c3436cfSJed Brown */
944d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateStep_RosW(TS ts, PetscInt order, Vec U, PetscBool *done)
945d71ae5a4SJacob Faibussowitsch {
9461c3436cfSJed Brown   TS_RosW     *ros = (TS_RosW *)ts->data;
9471c3436cfSJed Brown   RosWTableau  tab = ros->tableau;
9481c3436cfSJed Brown   PetscScalar *w   = ros->work;
9491c3436cfSJed Brown   PetscInt     i;
9501c3436cfSJed Brown 
9511c3436cfSJed Brown   PetscFunctionBegin;
9521c3436cfSJed Brown   if (order == tab->order) {
953108c343cSJed Brown     if (ros->status == TS_STEP_INCOMPLETE) { /* Use standard completion formula */
9549566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol, U));
955de19f811SJed Brown       for (i = 0; i < tab->s; i++) w[i] = tab->bt[i];
9569566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(U, tab->s, w, ros->Y));
9579566063dSJacob Faibussowitsch     } else PetscCall(VecCopy(ts->vec_sol, U));
9581c3436cfSJed Brown     if (done) *done = PETSC_TRUE;
9593ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
9601c3436cfSJed Brown   } else if (order == tab->order - 1) {
9611c3436cfSJed Brown     if (!tab->bembedt) goto unavailable;
962108c343cSJed Brown     if (ros->status == TS_STEP_INCOMPLETE) { /* Use embedded completion formula */
9639566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol, U));
964de19f811SJed Brown       for (i = 0; i < tab->s; i++) w[i] = tab->bembedt[i];
9659566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(U, tab->s, w, ros->Y));
966108c343cSJed Brown     } else { /* Use rollback-and-recomplete formula (bembedt - bt) */
967108c343cSJed Brown       for (i = 0; i < tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i];
9689566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol, U));
9699566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(U, tab->s, w, ros->Y));
9701c3436cfSJed Brown     }
9711c3436cfSJed Brown     if (done) *done = PETSC_TRUE;
9723ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
9731c3436cfSJed Brown   }
9741c3436cfSJed Brown unavailable:
9751c3436cfSJed Brown   if (done) *done = PETSC_FALSE;
9769371c9d4SSatish Balay   else
9779371c9d4SSatish 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,
9789371c9d4SSatish Balay             tab->order, order);
9793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9801c3436cfSJed Brown }
9811c3436cfSJed Brown 
982d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRollBack_RosW(TS ts)
983d71ae5a4SJacob Faibussowitsch {
98424655328SShri   TS_RosW *ros = (TS_RosW *)ts->data;
98524655328SShri 
98624655328SShri   PetscFunctionBegin;
9879566063dSJacob Faibussowitsch   PetscCall(VecCopy(ros->vec_sol_prev, ts->vec_sol));
9883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
98924655328SShri }
99024655328SShri 
991d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_RosW(TS ts)
992d71ae5a4SJacob Faibussowitsch {
99361692a83SJed Brown   TS_RosW         *ros = (TS_RosW *)ts->data;
99461692a83SJed Brown   RosWTableau      tab = ros->tableau;
995e27a552bSJed Brown   const PetscInt   s   = tab->s;
9961c3436cfSJed Brown   const PetscReal *At = tab->At, *Gamma = tab->Gamma, *ASum = tab->ASum, *GammaInv = tab->GammaInv;
9970feba352SEmil Constantinescu   const PetscReal *GammaExplicitCorr = tab->GammaExplicitCorr;
998c17803e7SJed Brown   const PetscBool *GammaZeroDiag     = tab->GammaZeroDiag;
99961692a83SJed Brown   PetscScalar     *w                 = ros->work;
10007d4bf2deSEmil Constantinescu   Vec             *Y = ros->Y, Ydot = ros->Ydot, Zdot = ros->Zdot, Zstage = ros->Zstage;
1001e27a552bSJed Brown   SNES             snes;
10021c3436cfSJed Brown   TSAdapt          adapt;
1003fecfb714SLisandro Dalcin   PetscInt         i, j, its, lits;
1004be5899b3SLisandro Dalcin   PetscInt         rejections = 0;
1005b39943a6SLisandro Dalcin   PetscBool        stageok, accept = PETSC_TRUE;
1006b39943a6SLisandro Dalcin   PetscReal        next_time_step = ts->time_step;
1007f7f07198SBarry Smith   PetscInt         lag;
1008e27a552bSJed Brown 
1009e27a552bSJed Brown   PetscFunctionBegin;
101048a46eb9SPierre Jolivet   if (!ts->steprollback) PetscCall(VecCopy(ts->vec_sol, ros->vec_sol_prev));
1011e27a552bSJed Brown 
1012b39943a6SLisandro Dalcin   ros->status = TS_STEP_INCOMPLETE;
1013b39943a6SLisandro Dalcin   while (!ts->reason && ros->status != TS_STEP_COMPLETE) {
10141c3436cfSJed Brown     const PetscReal h = ts->time_step;
1015e27a552bSJed Brown     for (i = 0; i < s; i++) {
10161c3436cfSJed Brown       ros->stage_time = ts->ptime + h * ASum[i];
10179566063dSJacob Faibussowitsch       PetscCall(TSPreStage(ts, ros->stage_time));
1018c17803e7SJed Brown       if (GammaZeroDiag[i]) {
1019c17803e7SJed Brown         ros->stage_explicit = PETSC_TRUE;
1020b296d7d5SJed Brown         ros->scoeff         = 1.;
1021c17803e7SJed Brown       } else {
1022c17803e7SJed Brown         ros->stage_explicit = PETSC_FALSE;
1023b296d7d5SJed Brown         ros->scoeff         = 1. / Gamma[i * s + i];
1024fd96d5b0SEmil Constantinescu       }
102561692a83SJed Brown 
10269566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol, Zstage));
1027de19f811SJed Brown       for (j = 0; j < i; j++) w[j] = At[i * s + j];
10289566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(Zstage, i, w, Y));
102961692a83SJed Brown 
103061692a83SJed Brown       for (j = 0; j < i; j++) w[j] = 1. / h * GammaInv[i * s + j];
10319566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(Zdot));
10329566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(Zdot, i, w, Y));
103361692a83SJed Brown 
1034e27a552bSJed Brown       /* Initial guess taken from last stage */
10359566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(Y[i]));
103661692a83SJed Brown 
10377d4bf2deSEmil Constantinescu       if (!ros->stage_explicit) {
10389566063dSJacob Faibussowitsch         PetscCall(TSGetSNES(ts, &snes));
103961692a83SJed Brown         if (!ros->recompute_jacobian && !i) {
10409566063dSJacob Faibussowitsch           PetscCall(SNESGetLagJacobian(snes, &lag));
10416aad120cSJose E. Roman           if (lag == 1) {                            /* use did not set a nontrivial lag, so lag over all stages */
10429566063dSJacob Faibussowitsch             PetscCall(SNESSetLagJacobian(snes, -2)); /* Recompute the Jacobian on this solve, but not again for the rest of the stages */
1043f7f07198SBarry Smith           }
104461692a83SJed Brown         }
10459566063dSJacob Faibussowitsch         PetscCall(SNESSolve(snes, NULL, Y[i]));
10469371c9d4SSatish 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 */ }
10479566063dSJacob Faibussowitsch         PetscCall(SNESGetIterationNumber(snes, &its));
10489566063dSJacob Faibussowitsch         PetscCall(SNESGetLinearSolveIterations(snes, &lits));
10499371c9d4SSatish Balay         ts->snes_its += its;
10509371c9d4SSatish Balay         ts->ksp_its += lits;
10517d4bf2deSEmil Constantinescu       } else {
10521ce71dffSSatish Balay         Mat J, Jp;
10539566063dSJacob Faibussowitsch         PetscCall(VecZeroEntries(Ydot)); /* Evaluate Y[i]=G(t,Ydot=0,Zstage) */
10549566063dSJacob Faibussowitsch         PetscCall(TSComputeIFunction(ts, ros->stage_time, Zstage, Ydot, Y[i], PETSC_FALSE));
10559566063dSJacob Faibussowitsch         PetscCall(VecScale(Y[i], -1.0));
10569566063dSJacob Faibussowitsch         PetscCall(VecAXPY(Y[i], -1.0, Zdot)); /*Y[i] = F(Zstage)-Zdot[=GammaInv*Y]*/
10570feba352SEmil Constantinescu 
10589566063dSJacob Faibussowitsch         PetscCall(VecZeroEntries(Zstage)); /* Zstage = GammaExplicitCorr[i,j] * Y[j] */
10590feba352SEmil Constantinescu         for (j = 0; j < i; j++) w[j] = GammaExplicitCorr[i * s + j];
10609566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(Zstage, i, w, Y));
1061fecfb714SLisandro Dalcin 
1062fecfb714SLisandro Dalcin         /* Y[i] = Y[i] + Jac*Zstage[=Jac*GammaExplicitCorr[i,j] * Y[j]] */
10639566063dSJacob Faibussowitsch         PetscCall(TSGetIJacobian(ts, &J, &Jp, NULL, NULL));
10649566063dSJacob Faibussowitsch         PetscCall(TSComputeIJacobian(ts, ros->stage_time, ts->vec_sol, Ydot, 0, J, Jp, PETSC_FALSE));
10659566063dSJacob Faibussowitsch         PetscCall(MatMult(J, Zstage, Zdot));
10669566063dSJacob Faibussowitsch         PetscCall(VecAXPY(Y[i], -1.0, Zdot));
10675ef26d82SJed Brown         ts->ksp_its += 1;
1068fecfb714SLisandro Dalcin 
10699566063dSJacob Faibussowitsch         PetscCall(VecScale(Y[i], h));
10707d4bf2deSEmil Constantinescu       }
10719566063dSJacob Faibussowitsch       PetscCall(TSPostStage(ts, ros->stage_time, i, Y));
10729566063dSJacob Faibussowitsch       PetscCall(TSGetAdapt(ts, &adapt));
10739566063dSJacob Faibussowitsch       PetscCall(TSAdaptCheckStage(adapt, ts, ros->stage_time, Y[i], &stageok));
1074fecfb714SLisandro Dalcin       if (!stageok) goto reject_step;
1075e27a552bSJed Brown     }
1076e27a552bSJed Brown 
1077b39943a6SLisandro Dalcin     ros->status = TS_STEP_INCOMPLETE;
10789566063dSJacob Faibussowitsch     PetscCall(TSEvaluateStep_RosW(ts, tab->order, ts->vec_sol, NULL));
1079b39943a6SLisandro Dalcin     ros->status = TS_STEP_PENDING;
10809566063dSJacob Faibussowitsch     PetscCall(TSGetAdapt(ts, &adapt));
10819566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidatesClear(adapt));
10829566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidateAdd(adapt, tab->name, tab->order, 1, tab->ccfl, (PetscReal)tab->s, PETSC_TRUE));
10839566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
1084b39943a6SLisandro Dalcin     ros->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
1085b39943a6SLisandro Dalcin     if (!accept) { /* Roll back the current step */
10869566063dSJacob Faibussowitsch       PetscCall(TSRollBack_RosW(ts));
1087be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
1088be5899b3SLisandro Dalcin       goto reject_step;
1089b39943a6SLisandro Dalcin     }
1090b39943a6SLisandro Dalcin 
1091e27a552bSJed Brown     ts->ptime += ts->time_step;
1092cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
10931c3436cfSJed Brown     break;
1094b39943a6SLisandro Dalcin 
1095b39943a6SLisandro Dalcin   reject_step:
10969371c9d4SSatish Balay     ts->reject++;
10979371c9d4SSatish Balay     accept = PETSC_FALSE;
1098be5899b3SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
1099b39943a6SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
110063a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
11011c3436cfSJed Brown     }
11021c3436cfSJed Brown   }
11033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1104e27a552bSJed Brown }
1105e27a552bSJed Brown 
1106d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_RosW(TS ts, PetscReal itime, Vec U)
1107d71ae5a4SJacob Faibussowitsch {
110861692a83SJed Brown   TS_RosW         *ros = (TS_RosW *)ts->data;
1109f4aed992SEmil Constantinescu   PetscInt         s = ros->tableau->s, pinterp = ros->tableau->pinterp, i, j;
1110f4aed992SEmil Constantinescu   PetscReal        h;
1111f4aed992SEmil Constantinescu   PetscReal        tt, t;
1112f4aed992SEmil Constantinescu   PetscScalar     *bt;
1113f4aed992SEmil Constantinescu   const PetscReal *Bt       = ros->tableau->binterpt;
1114f4aed992SEmil Constantinescu   const PetscReal *GammaInv = ros->tableau->GammaInv;
1115f4aed992SEmil Constantinescu   PetscScalar     *w        = ros->work;
1116f4aed992SEmil Constantinescu   Vec             *Y        = ros->Y;
1117e27a552bSJed Brown 
1118e27a552bSJed Brown   PetscFunctionBegin;
11193c633725SBarry Smith   PetscCheck(Bt, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSRosW %s does not have an interpolation formula", ros->tableau->name);
1120f4aed992SEmil Constantinescu 
1121f4aed992SEmil Constantinescu   switch (ros->status) {
1122f4aed992SEmil Constantinescu   case TS_STEP_INCOMPLETE:
1123f4aed992SEmil Constantinescu   case TS_STEP_PENDING:
1124f4aed992SEmil Constantinescu     h = ts->time_step;
1125f4aed992SEmil Constantinescu     t = (itime - ts->ptime) / h;
1126f4aed992SEmil Constantinescu     break;
1127f4aed992SEmil Constantinescu   case TS_STEP_COMPLETE:
1128be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev;
1129f4aed992SEmil Constantinescu     t = (itime - ts->ptime) / h + 1; /* In the interval [0,1] */
1130f4aed992SEmil Constantinescu     break;
1131d71ae5a4SJacob Faibussowitsch   default:
1132d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus");
1133f4aed992SEmil Constantinescu   }
11349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s, &bt));
1135f4aed992SEmil Constantinescu   for (i = 0; i < s; i++) bt[i] = 0;
1136f4aed992SEmil Constantinescu   for (j = 0, tt = t; j < pinterp; j++, tt *= t) {
1137ad540459SPierre Jolivet     for (i = 0; i < s; i++) bt[i] += Bt[i * pinterp + j] * tt;
1138f4aed992SEmil Constantinescu   }
1139f4aed992SEmil Constantinescu 
1140f4aed992SEmil Constantinescu   /* y(t+tt*h) = y(t) + Sum bt(tt) * GammaInv * Ydot */
1141f9c1d6abSBarry Smith   /* U <- 0*/
11429566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(U));
1143f9c1d6abSBarry Smith   /* U <- Sum bt_i * GammaInv(i,1:i) * Y(1:i) */
11443ca35412SEmil Constantinescu   for (j = 0; j < s; j++) w[j] = 0;
11453ca35412SEmil Constantinescu   for (j = 0; j < s; j++) {
1146ad540459SPierre Jolivet     for (i = j; i < s; i++) w[j] += bt[i] * GammaInv[i * s + j];
11473ca35412SEmil Constantinescu   }
11489566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(U, i, w, Y));
1149be5899b3SLisandro Dalcin   /* U <- y(t) + U */
11509566063dSJacob Faibussowitsch   PetscCall(VecAXPY(U, 1, ros->vec_sol_prev));
1151f4aed992SEmil Constantinescu 
11529566063dSJacob Faibussowitsch   PetscCall(PetscFree(bt));
11533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1154e27a552bSJed Brown }
1155e27a552bSJed Brown 
1156e27a552bSJed Brown /*------------------------------------------------------------*/
1157b39943a6SLisandro Dalcin 
1158d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWTableauReset(TS ts)
1159d71ae5a4SJacob Faibussowitsch {
1160b39943a6SLisandro Dalcin   TS_RosW    *ros = (TS_RosW *)ts->data;
1161b39943a6SLisandro Dalcin   RosWTableau tab = ros->tableau;
1162b39943a6SLisandro Dalcin 
1163b39943a6SLisandro Dalcin   PetscFunctionBegin;
11643ba16761SJacob Faibussowitsch   if (!tab) PetscFunctionReturn(PETSC_SUCCESS);
11659566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(tab->s, &ros->Y));
11669566063dSJacob Faibussowitsch   PetscCall(PetscFree(ros->work));
11673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1168b39943a6SLisandro Dalcin }
1169b39943a6SLisandro Dalcin 
1170d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_RosW(TS ts)
1171d71ae5a4SJacob Faibussowitsch {
117261692a83SJed Brown   TS_RosW *ros = (TS_RosW *)ts->data;
1173e27a552bSJed Brown 
1174e27a552bSJed Brown   PetscFunctionBegin;
11759566063dSJacob Faibussowitsch   PetscCall(TSRosWTableauReset(ts));
11769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ros->Ydot));
11779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ros->Ystage));
11789566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ros->Zdot));
11799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ros->Zstage));
11809566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ros->vec_sol_prev));
11813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1182e27a552bSJed Brown }
1183e27a552bSJed Brown 
1184d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWGetVecs(TS ts, DM dm, Vec *Ydot, Vec *Zdot, Vec *Ystage, Vec *Zstage)
1185d71ae5a4SJacob Faibussowitsch {
1186d5e6173cSPeter Brune   TS_RosW *rw = (TS_RosW *)ts->data;
1187d5e6173cSPeter Brune 
1188d5e6173cSPeter Brune   PetscFunctionBegin;
1189d5e6173cSPeter Brune   if (Ydot) {
1190d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
11919566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Ydot", Ydot));
1192d5e6173cSPeter Brune     } else *Ydot = rw->Ydot;
1193d5e6173cSPeter Brune   }
1194d5e6173cSPeter Brune   if (Zdot) {
1195d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
11969566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Zdot", Zdot));
1197d5e6173cSPeter Brune     } else *Zdot = rw->Zdot;
1198d5e6173cSPeter Brune   }
1199d5e6173cSPeter Brune   if (Ystage) {
1200d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
12019566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Ystage", Ystage));
1202d5e6173cSPeter Brune     } else *Ystage = rw->Ystage;
1203d5e6173cSPeter Brune   }
1204d5e6173cSPeter Brune   if (Zstage) {
1205d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
12069566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Zstage", Zstage));
1207d5e6173cSPeter Brune     } else *Zstage = rw->Zstage;
1208d5e6173cSPeter Brune   }
12093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1210d5e6173cSPeter Brune }
1211d5e6173cSPeter Brune 
1212d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWRestoreVecs(TS ts, DM dm, Vec *Ydot, Vec *Zdot, Vec *Ystage, Vec *Zstage)
1213d71ae5a4SJacob Faibussowitsch {
1214d5e6173cSPeter Brune   PetscFunctionBegin;
1215d5e6173cSPeter Brune   if (Ydot) {
121648a46eb9SPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Ydot", Ydot));
1217d5e6173cSPeter Brune   }
1218d5e6173cSPeter Brune   if (Zdot) {
121948a46eb9SPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Zdot", Zdot));
1220d5e6173cSPeter Brune   }
1221d5e6173cSPeter Brune   if (Ystage) {
122248a46eb9SPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Ystage", Ystage));
1223d5e6173cSPeter Brune   }
1224d5e6173cSPeter Brune   if (Zstage) {
122548a46eb9SPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Zstage", Zstage));
1226d5e6173cSPeter Brune   }
12273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1228d5e6173cSPeter Brune }
1229d5e6173cSPeter Brune 
1230d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCoarsenHook_TSRosW(DM fine, DM coarse, void *ctx)
1231d71ae5a4SJacob Faibussowitsch {
1232d5e6173cSPeter Brune   PetscFunctionBegin;
12333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1234d5e6173cSPeter Brune }
1235d5e6173cSPeter Brune 
1236d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRestrictHook_TSRosW(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx)
1237d71ae5a4SJacob Faibussowitsch {
1238d5e6173cSPeter Brune   TS  ts = (TS)ctx;
1239d5e6173cSPeter Brune   Vec Ydot, Zdot, Ystage, Zstage;
1240d5e6173cSPeter Brune   Vec Ydotc, Zdotc, Ystagec, Zstagec;
1241d5e6173cSPeter Brune 
1242d5e6173cSPeter Brune   PetscFunctionBegin;
12439566063dSJacob Faibussowitsch   PetscCall(TSRosWGetVecs(ts, fine, &Ydot, &Ystage, &Zdot, &Zstage));
12449566063dSJacob Faibussowitsch   PetscCall(TSRosWGetVecs(ts, coarse, &Ydotc, &Ystagec, &Zdotc, &Zstagec));
12459566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct, Ydot, Ydotc));
12469566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(Ydotc, rscale, Ydotc));
12479566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct, Ystage, Ystagec));
12489566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(Ystagec, rscale, Ystagec));
12499566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct, Zdot, Zdotc));
12509566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(Zdotc, rscale, Zdotc));
12519566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct, Zstage, Zstagec));
12529566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(Zstagec, rscale, Zstagec));
12539566063dSJacob Faibussowitsch   PetscCall(TSRosWRestoreVecs(ts, fine, &Ydot, &Ystage, &Zdot, &Zstage));
12549566063dSJacob Faibussowitsch   PetscCall(TSRosWRestoreVecs(ts, coarse, &Ydotc, &Ystagec, &Zdotc, &Zstagec));
12553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1256d5e6173cSPeter Brune }
1257d5e6173cSPeter Brune 
1258d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainHook_TSRosW(DM fine, DM coarse, void *ctx)
1259d71ae5a4SJacob Faibussowitsch {
1260258e1594SPeter Brune   PetscFunctionBegin;
12613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1262258e1594SPeter Brune }
1263258e1594SPeter Brune 
1264d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainRestrictHook_TSRosW(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx)
1265d71ae5a4SJacob Faibussowitsch {
1266258e1594SPeter Brune   TS  ts = (TS)ctx;
1267258e1594SPeter Brune   Vec Ydot, Zdot, Ystage, Zstage;
1268258e1594SPeter Brune   Vec Ydots, Zdots, Ystages, Zstages;
1269258e1594SPeter Brune 
1270258e1594SPeter Brune   PetscFunctionBegin;
12719566063dSJacob Faibussowitsch   PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Ystage, &Zdot, &Zstage));
12729566063dSJacob Faibussowitsch   PetscCall(TSRosWGetVecs(ts, subdm, &Ydots, &Ystages, &Zdots, &Zstages));
1273258e1594SPeter Brune 
12749566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(gscat, Ydot, Ydots, INSERT_VALUES, SCATTER_FORWARD));
12759566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(gscat, Ydot, Ydots, INSERT_VALUES, SCATTER_FORWARD));
1276258e1594SPeter Brune 
12779566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(gscat, Ystage, Ystages, INSERT_VALUES, SCATTER_FORWARD));
12789566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(gscat, Ystage, Ystages, INSERT_VALUES, SCATTER_FORWARD));
1279258e1594SPeter Brune 
12809566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(gscat, Zdot, Zdots, INSERT_VALUES, SCATTER_FORWARD));
12819566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(gscat, Zdot, Zdots, INSERT_VALUES, SCATTER_FORWARD));
1282258e1594SPeter Brune 
12839566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(gscat, Zstage, Zstages, INSERT_VALUES, SCATTER_FORWARD));
12849566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(gscat, Zstage, Zstages, INSERT_VALUES, SCATTER_FORWARD));
1285258e1594SPeter Brune 
12869566063dSJacob Faibussowitsch   PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Ystage, &Zdot, &Zstage));
12879566063dSJacob Faibussowitsch   PetscCall(TSRosWRestoreVecs(ts, subdm, &Ydots, &Ystages, &Zdots, &Zstages));
12883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1289258e1594SPeter Brune }
1290258e1594SPeter Brune 
1291e27a552bSJed Brown /*
1292e27a552bSJed Brown   This defines the nonlinear equation that is to be solved with SNES
1293e27a552bSJed Brown   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
1294e27a552bSJed Brown */
1295d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_RosW(SNES snes, Vec U, Vec F, TS ts)
1296d71ae5a4SJacob Faibussowitsch {
129761692a83SJed Brown   TS_RosW  *ros = (TS_RosW *)ts->data;
1298d5e6173cSPeter Brune   Vec       Ydot, Zdot, Ystage, Zstage;
1299b296d7d5SJed Brown   PetscReal shift = ros->scoeff / ts->time_step;
1300d5e6173cSPeter Brune   DM        dm, dmsave;
1301e27a552bSJed Brown 
1302e27a552bSJed Brown   PetscFunctionBegin;
13039566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
13049566063dSJacob Faibussowitsch   PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage));
13059566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(Ydot, shift, U, Zdot));   /* Ydot = shift*U + Zdot */
13069566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(Ystage, 1.0, U, Zstage)); /* Ystage = U + Zstage */
1307d5e6173cSPeter Brune   dmsave = ts->dm;
1308d5e6173cSPeter Brune   ts->dm = dm;
13099566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts, ros->stage_time, Ystage, Ydot, F, PETSC_FALSE));
1310d5e6173cSPeter Brune   ts->dm = dmsave;
13119566063dSJacob Faibussowitsch   PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage));
13123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1313e27a552bSJed Brown }
1314e27a552bSJed Brown 
1315d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes, Vec U, Mat A, Mat B, TS ts)
1316d71ae5a4SJacob Faibussowitsch {
131761692a83SJed Brown   TS_RosW  *ros = (TS_RosW *)ts->data;
1318d5e6173cSPeter Brune   Vec       Ydot, Zdot, Ystage, Zstage;
1319b296d7d5SJed Brown   PetscReal shift = ros->scoeff / ts->time_step;
1320d5e6173cSPeter Brune   DM        dm, dmsave;
1321e27a552bSJed Brown 
1322e27a552bSJed Brown   PetscFunctionBegin;
132361692a83SJed Brown   /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
13249566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
13259566063dSJacob Faibussowitsch   PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage));
1326d5e6173cSPeter Brune   dmsave = ts->dm;
1327d5e6173cSPeter Brune   ts->dm = dm;
13289566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts, ros->stage_time, Ystage, Ydot, shift, A, B, PETSC_TRUE));
1329d5e6173cSPeter Brune   ts->dm = dmsave;
13309566063dSJacob Faibussowitsch   PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage));
13313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1332e27a552bSJed Brown }
1333e27a552bSJed Brown 
1334d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWTableauSetUp(TS ts)
1335d71ae5a4SJacob Faibussowitsch {
1336b39943a6SLisandro Dalcin   TS_RosW    *ros = (TS_RosW *)ts->data;
1337b39943a6SLisandro Dalcin   RosWTableau tab = ros->tableau;
1338b39943a6SLisandro Dalcin 
1339b39943a6SLisandro Dalcin   PetscFunctionBegin;
13409566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ros->Y));
13419566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tab->s, &ros->work));
13423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1343b39943a6SLisandro Dalcin }
1344b39943a6SLisandro Dalcin 
1345d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_RosW(TS ts)
1346d71ae5a4SJacob Faibussowitsch {
134761692a83SJed Brown   TS_RosW      *ros = (TS_RosW *)ts->data;
1348d5e6173cSPeter Brune   DM            dm;
1349b39943a6SLisandro Dalcin   SNES          snes;
1350a3ab5968SHong Zhang   TSRHSJacobian rhsjacobian;
1351e27a552bSJed Brown 
1352e27a552bSJed Brown   PetscFunctionBegin;
13539566063dSJacob Faibussowitsch   PetscCall(TSRosWTableauSetUp(ts));
13549566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &ros->Ydot));
13559566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &ros->Ystage));
13569566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &ros->Zdot));
13579566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &ros->Zstage));
13589566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &ros->vec_sol_prev));
13599566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
13609566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSRosW, DMRestrictHook_TSRosW, ts));
13619566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSRosW, DMSubDomainRestrictHook_TSRosW, ts));
1362b39943a6SLisandro Dalcin   /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
13639566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
136448a46eb9SPierre Jolivet   if (!((PetscObject)snes)->type_name) PetscCall(SNESSetType(snes, SNESKSPONLY));
13659566063dSJacob Faibussowitsch   PetscCall(DMTSGetRHSJacobian(dm, &rhsjacobian, NULL));
1366a3ab5968SHong Zhang   if (rhsjacobian == TSComputeRHSJacobianConstant) {
1367a3ab5968SHong Zhang     Mat Amat, Pmat;
1368a3ab5968SHong Zhang 
1369a3ab5968SHong Zhang     /* Set the SNES matrix to be different from the RHS matrix because there is no way to reconstruct shift*M-J */
13709566063dSJacob Faibussowitsch     PetscCall(SNESGetJacobian(snes, &Amat, &Pmat, NULL, NULL));
1371a3ab5968SHong Zhang     if (Amat && Amat == ts->Arhs) {
1372a3ab5968SHong Zhang       if (Amat == Pmat) {
13739566063dSJacob Faibussowitsch         PetscCall(MatDuplicate(ts->Arhs, MAT_COPY_VALUES, &Amat));
13749566063dSJacob Faibussowitsch         PetscCall(SNESSetJacobian(snes, Amat, Amat, NULL, NULL));
1375a3ab5968SHong Zhang       } else {
13769566063dSJacob Faibussowitsch         PetscCall(MatDuplicate(ts->Arhs, MAT_COPY_VALUES, &Amat));
13779566063dSJacob Faibussowitsch         PetscCall(SNESSetJacobian(snes, Amat, NULL, NULL, NULL));
1378a3ab5968SHong Zhang         if (Pmat && Pmat == ts->Brhs) {
13799566063dSJacob Faibussowitsch           PetscCall(MatDuplicate(ts->Brhs, MAT_COPY_VALUES, &Pmat));
13809566063dSJacob Faibussowitsch           PetscCall(SNESSetJacobian(snes, NULL, Pmat, NULL, NULL));
13819566063dSJacob Faibussowitsch           PetscCall(MatDestroy(&Pmat));
1382a3ab5968SHong Zhang         }
1383a3ab5968SHong Zhang       }
13849566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Amat));
1385a3ab5968SHong Zhang     }
1386a3ab5968SHong Zhang   }
13873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1388e27a552bSJed Brown }
1389e27a552bSJed Brown /*------------------------------------------------------------*/
1390e27a552bSJed Brown 
1391d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_RosW(TS ts, PetscOptionItems *PetscOptionsObject)
1392d71ae5a4SJacob Faibussowitsch {
139361692a83SJed Brown   TS_RosW *ros = (TS_RosW *)ts->data;
1394b39943a6SLisandro Dalcin   SNES     snes;
1395e27a552bSJed Brown 
1396e27a552bSJed Brown   PetscFunctionBegin;
1397d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "RosW ODE solver options");
1398e27a552bSJed Brown   {
139961692a83SJed Brown     RosWTableauLink link;
1400e27a552bSJed Brown     PetscInt        count, choice;
1401e27a552bSJed Brown     PetscBool       flg;
1402e27a552bSJed Brown     const char    **namelist;
140361692a83SJed Brown 
14049371c9d4SSatish Balay     for (link = RosWTableauList, count = 0; link; link = link->next, count++)
14059371c9d4SSatish Balay       ;
14069566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(count, (char ***)&namelist));
140761692a83SJed Brown     for (link = RosWTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name;
14089566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEList("-ts_rosw_type", "Family of Rosenbrock-W method", "TSRosWSetType", (const char *const *)namelist, count, ros->tableau->name, &choice, &flg));
14099566063dSJacob Faibussowitsch     if (flg) PetscCall(TSRosWSetType(ts, namelist[choice]));
14109566063dSJacob Faibussowitsch     PetscCall(PetscFree(namelist));
141161692a83SJed Brown 
14129566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ts_rosw_recompute_jacobian", "Recompute the Jacobian at each stage", "TSRosWSetRecomputeJacobian", ros->recompute_jacobian, &ros->recompute_jacobian, NULL));
1413b39943a6SLisandro Dalcin   }
1414d0609cedSBarry Smith   PetscOptionsHeadEnd();
141561692a83SJed Brown   /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
14169566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
141748a46eb9SPierre Jolivet   if (!((PetscObject)snes)->type_name) PetscCall(SNESSetType(snes, SNESKSPONLY));
14183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1419e27a552bSJed Brown }
1420e27a552bSJed Brown 
1421d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_RosW(TS ts, PetscViewer viewer)
1422d71ae5a4SJacob Faibussowitsch {
142361692a83SJed Brown   TS_RosW  *ros = (TS_RosW *)ts->data;
1424e27a552bSJed Brown   PetscBool iascii;
1425e27a552bSJed Brown 
1426e27a552bSJed Brown   PetscFunctionBegin;
14279566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1428e27a552bSJed Brown   if (iascii) {
14299c334d8fSLisandro Dalcin     RosWTableau tab = ros->tableau;
143019fd82e9SBarry Smith     TSRosWType  rostype;
14319c334d8fSLisandro Dalcin     char        buf[512];
1432e408995aSJed Brown     PetscInt    i;
1433e408995aSJed Brown     PetscReal   abscissa[512];
14349566063dSJacob Faibussowitsch     PetscCall(TSRosWGetType(ts, &rostype));
14359566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Rosenbrock-W %s\n", rostype));
14369566063dSJacob Faibussowitsch     PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, tab->ASum));
14379566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Abscissa of A       = %s\n", buf));
1438e408995aSJed Brown     for (i = 0; i < tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i];
14399566063dSJacob Faibussowitsch     PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, abscissa));
14409566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Abscissa of A+Gamma = %s\n", buf));
1441e27a552bSJed Brown   }
14423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1443e27a552bSJed Brown }
1444e27a552bSJed Brown 
1445d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSLoad_RosW(TS ts, PetscViewer viewer)
1446d71ae5a4SJacob Faibussowitsch {
14479200755eSBarry Smith   SNES    snes;
14489c334d8fSLisandro Dalcin   TSAdapt adapt;
14499200755eSBarry Smith 
14509200755eSBarry Smith   PetscFunctionBegin;
14519566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &adapt));
14529566063dSJacob Faibussowitsch   PetscCall(TSAdaptLoad(adapt, viewer));
14539566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
14549566063dSJacob Faibussowitsch   PetscCall(SNESLoad(snes, viewer));
14559200755eSBarry Smith   /* function and Jacobian context for SNES when used with TS is always ts object */
14569566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, NULL, NULL, ts));
14579566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, ts));
14583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14599200755eSBarry Smith }
14609200755eSBarry Smith 
1461e27a552bSJed Brown /*@C
1462bcf0153eSBarry Smith   TSRosWSetType - Set the type of Rosenbrock-W, `TSROSW`, scheme
1463e27a552bSJed Brown 
146420f4b53cSBarry Smith   Logically Collective
1465e27a552bSJed Brown 
1466d8d19677SJose E. Roman   Input Parameters:
1467e27a552bSJed Brown + ts       - timestepping context
1468b92453a8SLisandro Dalcin - roswtype - type of Rosenbrock-W scheme
1469e27a552bSJed Brown 
1470020d8f30SJed Brown   Level: beginner
1471e27a552bSJed Brown 
14721cc06b55SBarry Smith .seealso: [](ch_ts), `TSRosWGetType()`, `TSROSW`, `TSROSW2M`, `TSROSW2P`, `TSROSWRA3PW`, `TSROSWRA34PW2`, `TSROSWRODAS3`, `TSROSWSANDU3`, `TSROSWASSP3P3S1C`, `TSROSWLASSP3P4S2C`, `TSROSWLLSSP3P4S2C`, `TSROSWARK3`
1473e27a552bSJed Brown @*/
1474d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWSetType(TS ts, TSRosWType roswtype)
1475d71ae5a4SJacob Faibussowitsch {
1476e27a552bSJed Brown   PetscFunctionBegin;
1477e27a552bSJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
14784f572ea9SToby Isaac   PetscAssertPointer(roswtype, 2);
1479cac4c232SBarry Smith   PetscTryMethod(ts, "TSRosWSetType_C", (TS, TSRosWType), (ts, roswtype));
14803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1481e27a552bSJed Brown }
1482e27a552bSJed Brown 
1483e27a552bSJed Brown /*@C
148461692a83SJed Brown   TSRosWGetType - Get the type of Rosenbrock-W scheme
1485e27a552bSJed Brown 
148620f4b53cSBarry Smith   Logically Collective
1487e27a552bSJed Brown 
1488e27a552bSJed Brown   Input Parameter:
1489e27a552bSJed Brown . ts - timestepping context
1490e27a552bSJed Brown 
1491e27a552bSJed Brown   Output Parameter:
149261692a83SJed Brown . rostype - type of Rosenbrock-W scheme
1493e27a552bSJed Brown 
1494e27a552bSJed Brown   Level: intermediate
1495e27a552bSJed Brown 
14961cc06b55SBarry Smith .seealso: [](ch_ts), `TSRosWType`, `TSRosWSetType()`
1497e27a552bSJed Brown @*/
1498d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWGetType(TS ts, TSRosWType *rostype)
1499d71ae5a4SJacob Faibussowitsch {
1500e27a552bSJed Brown   PetscFunctionBegin;
1501e27a552bSJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1502cac4c232SBarry Smith   PetscUseMethod(ts, "TSRosWGetType_C", (TS, TSRosWType *), (ts, rostype));
15033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1504e27a552bSJed Brown }
1505e27a552bSJed Brown 
1506e27a552bSJed Brown /*@C
150761692a83SJed Brown   TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step.
1508e27a552bSJed Brown 
150920f4b53cSBarry Smith   Logically Collective
1510e27a552bSJed Brown 
1511d8d19677SJose E. Roman   Input Parameters:
1512e27a552bSJed Brown + ts  - timestepping context
1513bcf0153eSBarry Smith - flg - `PETSC_TRUE` to recompute the Jacobian at each stage
1514e27a552bSJed Brown 
1515e27a552bSJed Brown   Level: intermediate
1516e27a552bSJed Brown 
15171cc06b55SBarry Smith .seealso: [](ch_ts), `TSRosWType`, `TSRosWGetType()`
1518e27a552bSJed Brown @*/
1519d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWSetRecomputeJacobian(TS ts, PetscBool flg)
1520d71ae5a4SJacob Faibussowitsch {
1521e27a552bSJed Brown   PetscFunctionBegin;
1522e27a552bSJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1523cac4c232SBarry Smith   PetscTryMethod(ts, "TSRosWSetRecomputeJacobian_C", (TS, PetscBool), (ts, flg));
15243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1525e27a552bSJed Brown }
1526e27a552bSJed Brown 
1527d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWGetType_RosW(TS ts, TSRosWType *rostype)
1528d71ae5a4SJacob Faibussowitsch {
152961692a83SJed Brown   TS_RosW *ros = (TS_RosW *)ts->data;
1530e27a552bSJed Brown 
1531e27a552bSJed Brown   PetscFunctionBegin;
153261692a83SJed Brown   *rostype = ros->tableau->name;
15333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1534e27a552bSJed Brown }
1535ef20d060SBarry Smith 
1536d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWSetType_RosW(TS ts, TSRosWType rostype)
1537d71ae5a4SJacob Faibussowitsch {
153861692a83SJed Brown   TS_RosW        *ros = (TS_RosW *)ts->data;
1539e27a552bSJed Brown   PetscBool       match;
154061692a83SJed Brown   RosWTableauLink link;
1541e27a552bSJed Brown 
1542e27a552bSJed Brown   PetscFunctionBegin;
154361692a83SJed Brown   if (ros->tableau) {
15449566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(ros->tableau->name, rostype, &match));
15453ba16761SJacob Faibussowitsch     if (match) PetscFunctionReturn(PETSC_SUCCESS);
1546e27a552bSJed Brown   }
154761692a83SJed Brown   for (link = RosWTableauList; link; link = link->next) {
15489566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(link->tab.name, rostype, &match));
1549e27a552bSJed Brown     if (match) {
15509566063dSJacob Faibussowitsch       if (ts->setupcalled) PetscCall(TSRosWTableauReset(ts));
155161692a83SJed Brown       ros->tableau = &link->tab;
15529566063dSJacob Faibussowitsch       if (ts->setupcalled) PetscCall(TSRosWTableauSetUp(ts));
15532ffb9264SLisandro Dalcin       ts->default_adapt_type = ros->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
15543ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
1555e27a552bSJed Brown     }
1556e27a552bSJed Brown   }
155798921bdaSJacob Faibussowitsch   SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", rostype);
1558e27a552bSJed Brown }
155961692a83SJed Brown 
1560d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWSetRecomputeJacobian_RosW(TS ts, PetscBool flg)
1561d71ae5a4SJacob Faibussowitsch {
156261692a83SJed Brown   TS_RosW *ros = (TS_RosW *)ts->data;
1563e27a552bSJed Brown 
1564e27a552bSJed Brown   PetscFunctionBegin;
156561692a83SJed Brown   ros->recompute_jacobian = flg;
15663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1567e27a552bSJed Brown }
1568e27a552bSJed Brown 
1569d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_RosW(TS ts)
1570d71ae5a4SJacob Faibussowitsch {
1571b3a6b972SJed Brown   PetscFunctionBegin;
15729566063dSJacob Faibussowitsch   PetscCall(TSReset_RosW(ts));
1573b3a6b972SJed Brown   if (ts->dm) {
15749566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSRosW, DMRestrictHook_TSRosW, ts));
15759566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSRosW, DMSubDomainRestrictHook_TSRosW, ts));
1576b3a6b972SJed Brown   }
15779566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
15789566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWGetType_C", NULL));
15799566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetType_C", NULL));
15809566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetRecomputeJacobian_C", NULL));
15813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1582b3a6b972SJed Brown }
1583d5e6173cSPeter Brune 
1584e27a552bSJed Brown /* ------------------------------------------------------------ */
1585e27a552bSJed Brown /*MC
1586020d8f30SJed Brown       TSROSW - ODE solver using Rosenbrock-W schemes
1587e27a552bSJed Brown 
1588e27a552bSJed Brown   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
1589e27a552bSJed Brown   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
1590bcf0153eSBarry Smith   of the equation using `TSSetIFunction()` and the non-stiff part with `TSSetRHSFunction()`.
1591bcf0153eSBarry Smith 
1592bcf0153eSBarry Smith   Level: beginner
1593e27a552bSJed Brown 
1594e27a552bSJed Brown   Notes:
159561692a83SJed Brown   This method currently only works with autonomous ODE and DAE.
159661692a83SJed Brown 
1597bcf0153eSBarry Smith   Consider trying `TSARKIMEX` if the stiff part is strongly nonlinear.
1598d0685a90SJed Brown 
15993d5a8a6aSBarry 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
16003d5a8a6aSBarry Smith 
1601e94cfbe0SPatrick Sanan   Developer Notes:
160261692a83SJed Brown   Rosenbrock-W methods are typically specified for autonomous ODE
160361692a83SJed Brown 
1604f9c1d6abSBarry Smith $  udot = f(u)
160561692a83SJed Brown 
160661692a83SJed Brown   by the stage equations
160761692a83SJed Brown 
1608f9c1d6abSBarry Smith $  k_i = h f(u_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j
160961692a83SJed Brown 
161061692a83SJed Brown   and step completion formula
161161692a83SJed Brown 
1612f9c1d6abSBarry Smith $  u_1 = u_0 + sum_j b_j k_j
161361692a83SJed Brown 
1614f9c1d6abSBarry Smith   with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(u)
161561692a83SJed Brown   and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner,
161661692a83SJed Brown   we define new variables for the stage equations
161761692a83SJed Brown 
161861692a83SJed Brown $  y_i = gamma_ij k_j
161961692a83SJed Brown 
162061692a83SJed Brown   The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define
162161692a83SJed Brown 
1622b70472e2SJed Brown $  A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-1}
162361692a83SJed Brown 
162461692a83SJed Brown   to rewrite the method as
162561692a83SJed Brown 
162637fdd005SBarry Smith .vb
162737fdd005SBarry 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
162837fdd005SBarry Smith   u_1 = u_0 + sum_j bt_j y_j
162937fdd005SBarry Smith .ve
163061692a83SJed Brown 
163161692a83SJed Brown    where we have introduced the mass matrix M. Continue by defining
163261692a83SJed Brown 
163361692a83SJed Brown $  ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j
163461692a83SJed Brown 
163561692a83SJed Brown    or, more compactly in tensor notation
163661692a83SJed Brown 
163761692a83SJed Brown $  Ydot = 1/h (Gamma^{-1} \otimes I) Y .
163861692a83SJed Brown 
163961692a83SJed Brown    Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current
164061692a83SJed Brown    stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the
164161692a83SJed Brown    equation
164261692a83SJed Brown 
1643f9c1d6abSBarry Smith $  g(u_0 + sum_j a_ij y_j + y_i, ydot_i) = 0
164461692a83SJed Brown 
164561692a83SJed Brown    with initial guess y_i = 0.
1646e27a552bSJed Brown 
16471cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSRosWSetType()`, `TSRosWRegister()`, `TSROSWTHETA1`, `TSROSWTHETA2`, `TSROSW2M`, `TSROSW2P`, `TSROSWRA3PW`, `TSROSWRA34PW2`, `TSROSWRODAS3`,
1648bcf0153eSBarry Smith           `TSROSWSANDU3`, `TSROSWASSP3P3S1C`, `TSROSWLASSP3P4S2C`, `TSROSWLLSSP3P4S2C`, `TSROSWGRK4T`, `TSROSWSHAMP4`, `TSROSWVELDD4`, `TSROSW4L`, `TSType`
1649e27a552bSJed Brown M*/
1650d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_RosW(TS ts)
1651d71ae5a4SJacob Faibussowitsch {
165261692a83SJed Brown   TS_RosW *ros;
1653e27a552bSJed Brown 
1654e27a552bSJed Brown   PetscFunctionBegin;
16559566063dSJacob Faibussowitsch   PetscCall(TSRosWInitializePackage());
1656e27a552bSJed Brown 
1657e27a552bSJed Brown   ts->ops->reset          = TSReset_RosW;
1658e27a552bSJed Brown   ts->ops->destroy        = TSDestroy_RosW;
1659e27a552bSJed Brown   ts->ops->view           = TSView_RosW;
16609200755eSBarry Smith   ts->ops->load           = TSLoad_RosW;
1661e27a552bSJed Brown   ts->ops->setup          = TSSetUp_RosW;
1662e27a552bSJed Brown   ts->ops->step           = TSStep_RosW;
1663e27a552bSJed Brown   ts->ops->interpolate    = TSInterpolate_RosW;
16641c3436cfSJed Brown   ts->ops->evaluatestep   = TSEvaluateStep_RosW;
166524655328SShri   ts->ops->rollback       = TSRollBack_RosW;
1666e27a552bSJed Brown   ts->ops->setfromoptions = TSSetFromOptions_RosW;
1667e27a552bSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_RosW;
1668e27a552bSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_RosW;
1669e27a552bSJed Brown 
1670efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
1671efd4aadfSBarry Smith 
16724dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&ros));
167361692a83SJed Brown   ts->data = (void *)ros;
1674e27a552bSJed Brown 
16759566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWGetType_C", TSRosWGetType_RosW));
16769566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetType_C", TSRosWSetType_RosW));
16779566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetRecomputeJacobian_C", TSRosWSetRecomputeJacobian_RosW));
1678b39943a6SLisandro Dalcin 
16799566063dSJacob Faibussowitsch   PetscCall(TSRosWSetType(ts, TSRosWDefault));
16803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1681e27a552bSJed Brown }
1682