xref: /petsc/src/ts/impls/rosw/rosw.c (revision bcf0153e883cfed9568ef4557dcc209048fb58f7)
1e27a552bSJed Brown /*
261692a83SJed Brown   Code for timestepping with Rosenbrock W methods
3e27a552bSJed Brown 
4e27a552bSJed Brown   Notes:
5e27a552bSJed Brown   The general system is written as
6e27a552bSJed Brown 
7f9c1d6abSBarry Smith   F(t,U,Udot) = G(t,U)
8e27a552bSJed Brown 
9f9c1d6abSBarry Smith   where F represents the stiff part of the physics and G represents the non-stiff part.
10f9c1d6abSBarry Smith   This method is designed to be linearly implicit on F and can use an approximate and lagged Jacobian.
11e27a552bSJed Brown 
12e27a552bSJed Brown */
13af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I   "petscts.h"   I*/
141e25c274SJed Brown #include <petscdm.h>
15e27a552bSJed Brown 
16af0996ceSBarry Smith #include <petsc/private/kernels/blockinvert.h>
1761692a83SJed Brown 
1819fd82e9SBarry Smith static TSRosWType TSRosWDefault = TSROSWRA34PW2;
19e27a552bSJed Brown static PetscBool  TSRosWRegisterAllCalled;
20e27a552bSJed Brown static PetscBool  TSRosWPackageInitialized;
21e27a552bSJed Brown 
2261692a83SJed Brown typedef struct _RosWTableau *RosWTableau;
2361692a83SJed Brown struct _RosWTableau {
24e27a552bSJed Brown   char      *name;
25e27a552bSJed Brown   PetscInt   order;             /* Classical approximation order of the method */
26e27a552bSJed Brown   PetscInt   s;                 /* Number of stages */
27f4aed992SEmil Constantinescu   PetscInt   pinterp;           /* Interpolation order */
2861692a83SJed Brown   PetscReal *A;                 /* Propagation table, strictly lower triangular */
2961692a83SJed Brown   PetscReal *Gamma;             /* Stage table, lower triangular with nonzero diagonal */
30c17803e7SJed Brown   PetscBool *GammaZeroDiag;     /* Diagonal entries that are zero in stage table Gamma, vector indicating explicit statages */
3143b21953SEmil Constantinescu   PetscReal *GammaExplicitCorr; /* Coefficients for correction terms needed for explicit stages in transformed variables*/
3261692a83SJed Brown   PetscReal *b;                 /* Step completion table */
33fe7e6d57SJed Brown   PetscReal *bembed;            /* Step completion table for embedded method of order one less */
3461692a83SJed Brown   PetscReal *ASum;              /* Row sum of A */
3561692a83SJed Brown   PetscReal *GammaSum;          /* Row sum of Gamma, only needed for non-autonomous systems */
3661692a83SJed Brown   PetscReal *At;                /* Propagation table in transformed variables */
3761692a83SJed Brown   PetscReal *bt;                /* Step completion table in transformed variables */
38fe7e6d57SJed Brown   PetscReal *bembedt;           /* Step completion table of order one less in transformed variables */
3961692a83SJed Brown   PetscReal *GammaInv;          /* Inverse of Gamma, used for transformed variables */
408d59e960SJed Brown   PetscReal  ccfl;              /* Placeholder for CFL coefficient relative to forward Euler */
413ca35412SEmil Constantinescu   PetscReal *binterpt;          /* Dense output formula */
42e27a552bSJed Brown };
4361692a83SJed Brown typedef struct _RosWTableauLink *RosWTableauLink;
4461692a83SJed Brown struct _RosWTableauLink {
4561692a83SJed Brown   struct _RosWTableau tab;
4661692a83SJed Brown   RosWTableauLink     next;
47e27a552bSJed Brown };
4861692a83SJed Brown static RosWTableauLink RosWTableauList;
49e27a552bSJed Brown 
50e27a552bSJed Brown typedef struct {
5161692a83SJed Brown   RosWTableau  tableau;
5261692a83SJed Brown   Vec         *Y;            /* States computed during the step, used to complete the step */
53e27a552bSJed Brown   Vec          Ydot;         /* Work vector holding Ydot during residual evaluation */
5461692a83SJed Brown   Vec          Ystage;       /* Work vector for the state value at each stage */
5561692a83SJed Brown   Vec          Zdot;         /* Ydot = Zdot + shift*Y */
5661692a83SJed Brown   Vec          Zstage;       /* Y = Zstage + Y */
57be5899b3SLisandro Dalcin   Vec          vec_sol_prev; /* Solution from the previous step (used for interpolation and rollback)*/
581c3436cfSJed Brown   PetscScalar *work;         /* Scalar work space of length number of stages, used to prepare VecMAXPY() */
59b296d7d5SJed Brown   PetscReal    scoeff;       /* shift = scoeff/dt */
60e27a552bSJed Brown   PetscReal    stage_time;
61c17803e7SJed Brown   PetscReal    stage_explicit;     /* Flag indicates that the current stage is explicit */
6261692a83SJed Brown   PetscBool    recompute_jacobian; /* Recompute the Jacobian at each stage, default is to freeze the Jacobian at the start of each step */
63108c343cSJed Brown   TSStepStatus status;
64e27a552bSJed Brown } TS_RosW;
65e27a552bSJed Brown 
66fe7e6d57SJed Brown /*MC
673606a31eSEmil Constantinescu      TSROSWTHETA1 - One stage first order L-stable Rosenbrock-W scheme (aka theta method).
683606a31eSEmil Constantinescu 
693606a31eSEmil Constantinescu      Only an approximate Jacobian is needed.
703606a31eSEmil Constantinescu 
713606a31eSEmil Constantinescu      Level: intermediate
723606a31eSEmil Constantinescu 
73*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSROSW`
743606a31eSEmil Constantinescu M*/
753606a31eSEmil Constantinescu 
763606a31eSEmil Constantinescu /*MC
773606a31eSEmil Constantinescu      TSROSWTHETA2 - One stage second order A-stable Rosenbrock-W scheme (aka theta method).
783606a31eSEmil Constantinescu 
793606a31eSEmil Constantinescu      Only an approximate Jacobian is needed.
803606a31eSEmil Constantinescu 
813606a31eSEmil Constantinescu      Level: intermediate
823606a31eSEmil Constantinescu 
83*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSROSW`
843606a31eSEmil Constantinescu M*/
853606a31eSEmil Constantinescu 
863606a31eSEmil Constantinescu /*MC
87fe7e6d57SJed Brown      TSROSW2M - Two stage second order L-stable Rosenbrock-W scheme.
88fe7e6d57SJed Brown 
89fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2P.
90fe7e6d57SJed Brown 
91fe7e6d57SJed Brown      Level: intermediate
92fe7e6d57SJed Brown 
93*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSROSW`
94fe7e6d57SJed Brown M*/
95fe7e6d57SJed Brown 
96fe7e6d57SJed Brown /*MC
97fe7e6d57SJed Brown      TSROSW2P - Two stage second order L-stable Rosenbrock-W scheme.
98fe7e6d57SJed Brown 
99fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2M.
100fe7e6d57SJed Brown 
101fe7e6d57SJed Brown      Level: intermediate
102fe7e6d57SJed Brown 
103*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSROSW`
104fe7e6d57SJed Brown M*/
105fe7e6d57SJed Brown 
106fe7e6d57SJed Brown /*MC
107fe7e6d57SJed Brown      TSROSWRA3PW - Three stage third order Rosenbrock-W scheme for PDAE of index 1.
108fe7e6d57SJed Brown 
109fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.
110fe7e6d57SJed Brown 
111fe7e6d57SJed Brown      This is strongly A-stable with R(infty) = 0.73. The embedded method of order 2 is strongly A-stable with R(infty) = 0.73.
112fe7e6d57SJed Brown 
113*bcf0153eSBarry Smith      Level: intermediate
114*bcf0153eSBarry Smith 
115fe7e6d57SJed Brown      References:
116606c0280SSatish Balay .  * - Rang and Angermann, New Rosenbrock W methods of order 3 for partial differential algebraic equations of index 1, 2005.
117fe7e6d57SJed Brown 
118*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSROSW`
119fe7e6d57SJed Brown M*/
120fe7e6d57SJed Brown 
121fe7e6d57SJed Brown /*MC
122fe7e6d57SJed Brown      TSROSWRA34PW2 - Four stage third order L-stable Rosenbrock-W scheme for PDAE of index 1.
123fe7e6d57SJed Brown 
124fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.
125fe7e6d57SJed Brown 
126fe7e6d57SJed Brown      This is strongly A-stable with R(infty) = 0. The embedded method of order 2 is strongly A-stable with R(infty) = 0.48.
127fe7e6d57SJed Brown 
128*bcf0153eSBarry Smith      Level: intermediate
129*bcf0153eSBarry Smith 
130fe7e6d57SJed Brown      References:
131606c0280SSatish Balay .  * - Rang and Angermann, New Rosenbrock W methods of order 3 for partial differential algebraic equations of index 1, 2005.
132fe7e6d57SJed Brown 
133*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSROSW`
134fe7e6d57SJed Brown M*/
135fe7e6d57SJed Brown 
136ef3c5b88SJed Brown /*MC
137ef3c5b88SJed Brown      TSROSWRODAS3 - Four stage third order L-stable Rosenbrock scheme
138ef3c5b88SJed Brown 
139ef3c5b88SJed Brown      By default, the Jacobian is only recomputed once per step.
140ef3c5b88SJed Brown 
141ef3c5b88SJed Brown      Both the third order and embedded second order methods are stiffly accurate and L-stable.
142ef3c5b88SJed Brown 
143*bcf0153eSBarry Smith      Level: intermediate
144*bcf0153eSBarry Smith 
145ef3c5b88SJed Brown      References:
146606c0280SSatish Balay .  * - Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997.
147ef3c5b88SJed Brown 
148*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSROSW`, `TSROSWSANDU3`
149ef3c5b88SJed Brown M*/
150ef3c5b88SJed Brown 
151ef3c5b88SJed Brown /*MC
152ef3c5b88SJed Brown      TSROSWSANDU3 - Three stage third order L-stable Rosenbrock scheme
153ef3c5b88SJed Brown 
154ef3c5b88SJed Brown      By default, the Jacobian is only recomputed once per step.
155ef3c5b88SJed Brown 
156ef3c5b88SJed Brown      The third order method is L-stable, but not stiffly accurate.
157ef3c5b88SJed Brown      The second order embedded method is strongly A-stable with R(infty) = 0.5.
158ef3c5b88SJed Brown      The internal stages are L-stable.
159ef3c5b88SJed Brown      This method is called ROS3 in the paper.
160ef3c5b88SJed Brown 
161*bcf0153eSBarry Smith      Level: intermediate
162*bcf0153eSBarry Smith 
163ef3c5b88SJed Brown      References:
164606c0280SSatish Balay .  * - Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997.
165ef3c5b88SJed Brown 
166*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSROSW`, `TSROSWRODAS3`
167ef3c5b88SJed Brown M*/
168ef3c5b88SJed Brown 
169961f28d0SJed Brown /*MC
170961f28d0SJed Brown      TSROSWASSP3P3S1C - A-stable Rosenbrock-W method with SSP explicit part, third order, three stages
171961f28d0SJed Brown 
172961f28d0SJed Brown      By default, the Jacobian is only recomputed once per step.
173961f28d0SJed Brown 
174961f28d0SJed Brown      A-stable SPP explicit order 3, 3 stages, CFL 1 (eff = 1/3)
175961f28d0SJed Brown 
176*bcf0153eSBarry Smith      Level: intermediate
177*bcf0153eSBarry Smith 
178961f28d0SJed Brown      References:
179606c0280SSatish Balay . * - Emil Constantinescu
180961f28d0SJed Brown 
181*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSROSW`, `TSROSWLASSP3P4S2C`, `TSROSWLLSSP3P4S2C`, `SSP`
182961f28d0SJed Brown M*/
183961f28d0SJed Brown 
184961f28d0SJed Brown /*MC
185998eb97aSJed Brown      TSROSWLASSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, four stages
186961f28d0SJed Brown 
187961f28d0SJed Brown      By default, the Jacobian is only recomputed once per step.
188961f28d0SJed Brown 
189961f28d0SJed Brown      L-stable (A-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2)
190961f28d0SJed Brown 
191*bcf0153eSBarry Smith      Level: intermediate
192*bcf0153eSBarry Smith 
193961f28d0SJed Brown      References:
194606c0280SSatish Balay . * - Emil Constantinescu
195961f28d0SJed Brown 
196*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSROSW`, `TSROSWASSP3P3S1C`, `TSROSWLLSSP3P4S2C`, `TSSSP`
197961f28d0SJed Brown M*/
198961f28d0SJed Brown 
199961f28d0SJed Brown /*MC
200998eb97aSJed Brown      TSROSWLLSSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, four stages
201961f28d0SJed Brown 
202961f28d0SJed Brown      By default, the Jacobian is only recomputed once per step.
203961f28d0SJed Brown 
204961f28d0SJed Brown      L-stable (L-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2)
205961f28d0SJed Brown 
206*bcf0153eSBarry Smith      Level: intermediate
207*bcf0153eSBarry Smith 
208961f28d0SJed Brown      References:
209606c0280SSatish Balay . * - Emil Constantinescu
210961f28d0SJed Brown 
211*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSROSW`, `TSROSWASSP3P3S1C`, `TSROSWLASSP3P4S2C`, `TSSSP`
212961f28d0SJed Brown M*/
213961f28d0SJed Brown 
21442faf41dSJed Brown /*MC
21542faf41dSJed Brown      TSROSWGRK4T - four stage, fourth order Rosenbrock (not W) method from Kaps and Rentrop
21642faf41dSJed Brown 
21742faf41dSJed Brown      By default, the Jacobian is only recomputed once per step.
21842faf41dSJed Brown 
21942faf41dSJed Brown      A(89.3 degrees)-stable, |R(infty)| = 0.454.
22042faf41dSJed Brown 
22142faf41dSJed Brown      This method does not provide a dense output formula.
22242faf41dSJed Brown 
223*bcf0153eSBarry Smith      Level: intermediate
224*bcf0153eSBarry Smith 
22542faf41dSJed Brown      References:
226606c0280SSatish Balay +   * -  Kaps and Rentrop, Generalized Runge Kutta methods of order four with stepsize control for stiff ordinary differential equations, 1979.
227606c0280SSatish Balay -   * -  Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.
22842faf41dSJed Brown 
22942faf41dSJed Brown      Hairer's code ros4.f
23042faf41dSJed Brown 
231*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSROSW`, `TSROSWSHAMP4`, `TSROSWVELDD4`, `TSROSW4L`
23242faf41dSJed Brown M*/
23342faf41dSJed Brown 
23442faf41dSJed Brown /*MC
23542faf41dSJed Brown      TSROSWSHAMP4 - four stage, fourth order Rosenbrock (not W) method from Shampine
23642faf41dSJed Brown 
23742faf41dSJed Brown      By default, the Jacobian is only recomputed once per step.
23842faf41dSJed Brown 
23942faf41dSJed Brown      A-stable, |R(infty)| = 1/3.
24042faf41dSJed Brown 
24142faf41dSJed Brown      This method does not provide a dense output formula.
24242faf41dSJed Brown 
243*bcf0153eSBarry Smith      Level: intermediate
244*bcf0153eSBarry Smith 
24542faf41dSJed Brown      References:
246606c0280SSatish Balay +   * -  Shampine, Implementation of Rosenbrock methods, 1982.
247606c0280SSatish Balay -   * -  Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.
24842faf41dSJed Brown 
24942faf41dSJed Brown      Hairer's code ros4.f
25042faf41dSJed Brown 
251*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSROSW`, `TSROSWGRK4T`, `TSROSWVELDD4`, `TSROSW4L`
25242faf41dSJed Brown M*/
25342faf41dSJed Brown 
25442faf41dSJed Brown /*MC
25542faf41dSJed Brown      TSROSWVELDD4 - four stage, fourth order Rosenbrock (not W) method from van Veldhuizen
25642faf41dSJed Brown 
25742faf41dSJed Brown      By default, the Jacobian is only recomputed once per step.
25842faf41dSJed Brown 
25942faf41dSJed Brown      A(89.5 degrees)-stable, |R(infty)| = 0.24.
26042faf41dSJed Brown 
26142faf41dSJed Brown      This method does not provide a dense output formula.
26242faf41dSJed Brown 
263*bcf0153eSBarry Smith      Level: intermediate
264*bcf0153eSBarry Smith 
26542faf41dSJed Brown      References:
266606c0280SSatish Balay +   * -  van Veldhuizen, D stability and Kaps Rentrop methods, 1984.
267606c0280SSatish Balay -   * -  Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.
26842faf41dSJed Brown 
26942faf41dSJed Brown      Hairer's code ros4.f
27042faf41dSJed Brown 
271*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSROSW`, `TSROSWGRK4T`, `TSROSWSHAMP4`, `TSROSW4L`
27242faf41dSJed Brown M*/
27342faf41dSJed Brown 
27442faf41dSJed Brown /*MC
27542faf41dSJed Brown      TSROSW4L - four stage, fourth order Rosenbrock (not W) method
27642faf41dSJed Brown 
27742faf41dSJed Brown      By default, the Jacobian is only recomputed once per step.
27842faf41dSJed Brown 
27942faf41dSJed Brown      A-stable and L-stable
28042faf41dSJed Brown 
28142faf41dSJed Brown      This method does not provide a dense output formula.
28242faf41dSJed Brown 
283*bcf0153eSBarry Smith      Level: intermediate
284*bcf0153eSBarry Smith 
28542faf41dSJed Brown      References:
286606c0280SSatish Balay .  * -   Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.
28742faf41dSJed Brown 
28842faf41dSJed Brown      Hairer's code ros4.f
28942faf41dSJed Brown 
290*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSROSW`, `TSROSWGRK4T`, `TSROSWSHAMP4`, `TSROSW4L`
29142faf41dSJed Brown M*/
29242faf41dSJed Brown 
293e27a552bSJed Brown /*@C
294*bcf0153eSBarry Smith   TSRosWRegisterAll - Registers all of the Rosenbrock-W methods in `TSROSW`
295e27a552bSJed Brown 
296e27a552bSJed Brown   Not Collective, but should be called by all processes which will need the schemes to be registered
297e27a552bSJed Brown 
298e27a552bSJed Brown   Level: advanced
299e27a552bSJed Brown 
300*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSROSW`, `TSRosWRegisterDestroy()`
301e27a552bSJed Brown @*/
302d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWRegisterAll(void)
303d71ae5a4SJacob Faibussowitsch {
304e27a552bSJed Brown   PetscFunctionBegin;
305e27a552bSJed Brown   if (TSRosWRegisterAllCalled) PetscFunctionReturn(0);
306e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_TRUE;
307e27a552bSJed Brown 
308e27a552bSJed Brown   {
309bbd56ea5SKarl Rupp     const PetscReal A        = 0;
310bbd56ea5SKarl Rupp     const PetscReal Gamma    = 1;
311bbd56ea5SKarl Rupp     const PetscReal b        = 1;
312bbd56ea5SKarl Rupp     const PetscReal binterpt = 1;
3131f80e275SEmil Constantinescu 
3149566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWTHETA1, 1, 1, &A, &Gamma, &b, NULL, 1, &binterpt));
3153606a31eSEmil Constantinescu   }
3163606a31eSEmil Constantinescu 
3173606a31eSEmil Constantinescu   {
318bbd56ea5SKarl Rupp     const PetscReal A        = 0;
319bbd56ea5SKarl Rupp     const PetscReal Gamma    = 0.5;
320bbd56ea5SKarl Rupp     const PetscReal b        = 1;
321bbd56ea5SKarl Rupp     const PetscReal binterpt = 1;
322bbd56ea5SKarl Rupp 
3239566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWTHETA2, 2, 1, &A, &Gamma, &b, NULL, 1, &binterpt));
3243606a31eSEmil Constantinescu   }
3253606a31eSEmil Constantinescu 
3263606a31eSEmil Constantinescu   {
327da80777bSKarl Rupp     /*const PetscReal g = 1. + 1./PetscSqrtReal(2.0);   Direct evaluation: 1.707106781186547524401. Used for setting up arrays of values known at compile time below. */
3289371c9d4SSatish Balay     const PetscReal A[2][2] =
3299371c9d4SSatish Balay       {
3309371c9d4SSatish Balay         {0,  0},
3319371c9d4SSatish Balay         {1., 0}
3329371c9d4SSatish Balay     },
3339371c9d4SSatish Balay                     Gamma[2][2] = {{1.707106781186547524401, 0}, {-2. * 1.707106781186547524401, 1.707106781186547524401}}, b[2] = {0.5, 0.5}, b1[2] = {1.0, 0.0};
3341f80e275SEmil Constantinescu     PetscReal binterpt[2][2];
335da80777bSKarl Rupp     binterpt[0][0] = 1.707106781186547524401 - 1.0;
336da80777bSKarl Rupp     binterpt[1][0] = 2.0 - 1.707106781186547524401;
337da80777bSKarl Rupp     binterpt[0][1] = 1.707106781186547524401 - 1.5;
338da80777bSKarl Rupp     binterpt[1][1] = 1.5 - 1.707106781186547524401;
339bbd56ea5SKarl Rupp 
3409566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSW2P, 2, 2, &A[0][0], &Gamma[0][0], b, b1, 2, &binterpt[0][0]));
341e27a552bSJed Brown   }
342e27a552bSJed Brown   {
343da80777bSKarl Rupp     /*const PetscReal g = 1. - 1./PetscSqrtReal(2.0);   Direct evaluation: 0.2928932188134524755992. Used for setting up arrays of values known at compile time below. */
3449371c9d4SSatish Balay     const PetscReal A[2][2] =
3459371c9d4SSatish Balay       {
3469371c9d4SSatish Balay         {0,  0},
3479371c9d4SSatish Balay         {1., 0}
3489371c9d4SSatish Balay     },
3499371c9d4SSatish Balay                     Gamma[2][2] = {{0.2928932188134524755992, 0}, {-2. * 0.2928932188134524755992, 0.2928932188134524755992}}, b[2] = {0.5, 0.5}, b1[2] = {1.0, 0.0};
3501f80e275SEmil Constantinescu     PetscReal binterpt[2][2];
351da80777bSKarl Rupp     binterpt[0][0] = 0.2928932188134524755992 - 1.0;
352da80777bSKarl Rupp     binterpt[1][0] = 2.0 - 0.2928932188134524755992;
353da80777bSKarl Rupp     binterpt[0][1] = 0.2928932188134524755992 - 1.5;
354da80777bSKarl Rupp     binterpt[1][1] = 1.5 - 0.2928932188134524755992;
355bbd56ea5SKarl Rupp 
3569566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSW2M, 2, 2, &A[0][0], &Gamma[0][0], b, b1, 2, &binterpt[0][0]));
357fe7e6d57SJed Brown   }
358fe7e6d57SJed Brown   {
359da80777bSKarl Rupp     /*const PetscReal g = 7.8867513459481287e-01; Directly written in-place below */
3601f80e275SEmil Constantinescu     PetscReal       binterpt[3][2];
3619371c9d4SSatish Balay     const PetscReal A[3][3] =
3629371c9d4SSatish Balay       {
3639371c9d4SSatish Balay         {0,                      0, 0},
364fe7e6d57SJed Brown         {1.5773502691896257e+00, 0, 0},
3659371c9d4SSatish Balay         {0.5,                    0, 0}
3669371c9d4SSatish Balay     },
3679371c9d4SSatish Balay                     Gamma[3][3] = {{7.8867513459481287e-01, 0, 0}, {-1.5773502691896257e+00, 7.8867513459481287e-01, 0}, {-6.7075317547305480e-01, -1.7075317547305482e-01, 7.8867513459481287e-01}}, b[3] = {1.0566243270259355e-01, 4.9038105676657971e-02, 8.4529946162074843e-01}, b2[3] = {-1.7863279495408180e-01, 1. / 3., 8.4529946162074843e-01};
3681f80e275SEmil Constantinescu 
3691f80e275SEmil Constantinescu     binterpt[0][0] = -0.8094010767585034;
3701f80e275SEmil Constantinescu     binterpt[1][0] = -0.5;
3711f80e275SEmil Constantinescu     binterpt[2][0] = 2.3094010767585034;
3721f80e275SEmil Constantinescu     binterpt[0][1] = 0.9641016151377548;
3731f80e275SEmil Constantinescu     binterpt[1][1] = 0.5;
3741f80e275SEmil Constantinescu     binterpt[2][1] = -1.4641016151377548;
375bbd56ea5SKarl Rupp 
3769566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWRA3PW, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0]));
377fe7e6d57SJed Brown   }
378fe7e6d57SJed Brown   {
3793ca35412SEmil Constantinescu     PetscReal binterpt[4][3];
380da80777bSKarl Rupp     /*const PetscReal g = 4.3586652150845900e-01; Directly written in-place below */
3819371c9d4SSatish Balay     const PetscReal A[4][4] =
3829371c9d4SSatish Balay       {
3839371c9d4SSatish Balay         {0,                      0,                       0,  0},
384fe7e6d57SJed Brown         {8.7173304301691801e-01, 0,                       0,  0},
385fe7e6d57SJed Brown         {8.4457060015369423e-01, -1.1299064236484185e-01, 0,  0},
3869371c9d4SSatish Balay         {0,                      0,                       1., 0}
3879371c9d4SSatish Balay     },
3889371c9d4SSatish Balay                     Gamma[4][4] = {{4.3586652150845900e-01, 0, 0, 0}, {-8.7173304301691801e-01, 4.3586652150845900e-01, 0, 0}, {-9.0338057013044082e-01, 5.4180672388095326e-02, 4.3586652150845900e-01, 0}, {2.4212380706095346e-01, -1.2232505839045147e+00, 5.4526025533510214e-01, 4.3586652150845900e-01}}, b[4] = {2.4212380706095346e-01, -1.2232505839045147e+00, 1.5452602553351020e+00, 4.3586652150845900e-01}, b2[4] = {3.7810903145819369e-01, -9.6042292212423178e-02, 5.0000000000000000e-01, 2.1793326075422950e-01};
3893ca35412SEmil Constantinescu 
3903ca35412SEmil Constantinescu     binterpt[0][0] = 1.0564298455794094;
3913ca35412SEmil Constantinescu     binterpt[1][0] = 2.296429974281067;
3923ca35412SEmil Constantinescu     binterpt[2][0] = -1.307599564525376;
3933ca35412SEmil Constantinescu     binterpt[3][0] = -1.045260255335102;
3943ca35412SEmil Constantinescu     binterpt[0][1] = -1.3864882699759573;
3953ca35412SEmil Constantinescu     binterpt[1][1] = -8.262611700275677;
3963ca35412SEmil Constantinescu     binterpt[2][1] = 7.250979895056055;
3973ca35412SEmil Constantinescu     binterpt[3][1] = 2.398120075195581;
3983ca35412SEmil Constantinescu     binterpt[0][2] = 0.5721822314575016;
3993ca35412SEmil Constantinescu     binterpt[1][2] = 4.742931142090097;
4003ca35412SEmil Constantinescu     binterpt[2][2] = -4.398120075195578;
4013ca35412SEmil Constantinescu     binterpt[3][2] = -0.9169932983520199;
4023ca35412SEmil Constantinescu 
4039566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWRA34PW2, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0]));
404e27a552bSJed Brown   }
405ef3c5b88SJed Brown   {
406da80777bSKarl Rupp     /* const PetscReal g = 0.5;       Directly written in-place below */
4079371c9d4SSatish Balay     const PetscReal A[4][4] =
4089371c9d4SSatish Balay       {
4099371c9d4SSatish Balay         {0,    0,     0,   0},
410ef3c5b88SJed Brown         {0,    0,     0,   0},
411ef3c5b88SJed Brown         {1.,   0,     0,   0},
4129371c9d4SSatish Balay         {0.75, -0.25, 0.5, 0}
4139371c9d4SSatish Balay     },
4149371c9d4SSatish Balay                     Gamma[4][4] = {{0.5, 0, 0, 0}, {1., 0.5, 0, 0}, {-0.25, -0.25, 0.5, 0}, {1. / 12, 1. / 12, -2. / 3, 0.5}}, b[4] = {5. / 6, -1. / 6, -1. / 6, 0.5}, b2[4] = {0.75, -0.25, 0.5, 0};
415bbd56ea5SKarl Rupp 
4169566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWRODAS3, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 0, NULL));
417ef3c5b88SJed Brown   }
418ef3c5b88SJed Brown   {
419da80777bSKarl Rupp     /*const PetscReal g = 0.43586652150845899941601945119356;       Directly written in-place below */
4209371c9d4SSatish Balay     const PetscReal A[3][3] =
4219371c9d4SSatish Balay       {
4229371c9d4SSatish Balay         {0,                                  0, 0},
423da80777bSKarl Rupp         {0.43586652150845899941601945119356, 0, 0},
4249371c9d4SSatish Balay         {0.43586652150845899941601945119356, 0, 0}
4259371c9d4SSatish Balay     },
4269371c9d4SSatish Balay                     Gamma[3][3] = {{0.43586652150845899941601945119356, 0, 0}, {-0.19294655696029095575009695436041, 0.43586652150845899941601945119356, 0}, {0, 1.74927148125794685173529749738960, 0.43586652150845899941601945119356}}, b[3] = {-0.75457412385404315829818998646589, 1.94100407061964420292840123379419, -0.18642994676560104463021124732829}, b2[3] = {-1.53358745784149585370766523913002, 2.81745131148625772213931745457622, -0.28386385364476186843165221544619};
4271f80e275SEmil Constantinescu 
4281f80e275SEmil Constantinescu     PetscReal binterpt[3][2];
4291f80e275SEmil Constantinescu     binterpt[0][0] = 3.793692883777660870425141387941;
4301f80e275SEmil Constantinescu     binterpt[1][0] = -2.918692883777660870425141387941;
4311f80e275SEmil Constantinescu     binterpt[2][0] = 0.125;
4321f80e275SEmil Constantinescu     binterpt[0][1] = -0.725741064379812106687651020584;
4331f80e275SEmil Constantinescu     binterpt[1][1] = 0.559074397713145440020984353917;
4341f80e275SEmil Constantinescu     binterpt[2][1] = 0.16666666666666666666666666666667;
4351f80e275SEmil Constantinescu 
4369566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWSANDU3, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0]));
437ef3c5b88SJed Brown   }
438b1c69cc3SEmil Constantinescu   {
439da80777bSKarl Rupp     /*const PetscReal s3 = PetscSqrtReal(3.),g = (3.0+s3)/6.0;
440da80777bSKarl Rupp      * Direct evaluation: s3 = 1.732050807568877293527;
441da80777bSKarl Rupp      *                     g = 0.7886751345948128822546;
442da80777bSKarl Rupp      * Values are directly inserted below to ensure availability at compile time (compiler warnings otherwise...) */
4439371c9d4SSatish Balay     const PetscReal A[3][3] =
4449371c9d4SSatish Balay       {
4459371c9d4SSatish Balay         {0,    0,    0},
446b1c69cc3SEmil Constantinescu         {1,    0,    0},
4479371c9d4SSatish Balay         {0.25, 0.25, 0}
4489371c9d4SSatish Balay     },
4499371c9d4SSatish Balay                     Gamma[3][3] = {{0, 0, 0}, {(-3.0 - 1.732050807568877293527) / 6.0, 0.7886751345948128822546, 0}, {(-3.0 - 1.732050807568877293527) / 24.0, (-3.0 - 1.732050807568877293527) / 8.0, 0.7886751345948128822546}}, b[3] = {1. / 6., 1. / 6., 2. / 3.}, b2[3] = {1. / 4., 1. / 4., 1. / 2.};
450c0cb691aSEmil Constantinescu     PetscReal binterpt[3][2];
451da80777bSKarl Rupp 
452c0cb691aSEmil Constantinescu     binterpt[0][0] = 0.089316397477040902157517886164709;
453c0cb691aSEmil Constantinescu     binterpt[1][0] = -0.91068360252295909784248211383529;
454c0cb691aSEmil Constantinescu     binterpt[2][0] = 1.8213672050459181956849642276706;
455c0cb691aSEmil Constantinescu     binterpt[0][1] = 0.077350269189625764509148780501957;
456c0cb691aSEmil Constantinescu     binterpt[1][1] = 1.077350269189625764509148780502;
457c0cb691aSEmil Constantinescu     binterpt[2][1] = -1.1547005383792515290182975610039;
458bbd56ea5SKarl Rupp 
4599566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWASSP3P3S1C, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0]));
460b1c69cc3SEmil Constantinescu   }
461b1c69cc3SEmil Constantinescu 
462b1c69cc3SEmil Constantinescu   {
4639371c9d4SSatish Balay     const PetscReal A[4][4] =
4649371c9d4SSatish Balay       {
4659371c9d4SSatish Balay         {0,       0,       0,       0},
466b1c69cc3SEmil Constantinescu         {1. / 2., 0,       0,       0},
467b1c69cc3SEmil Constantinescu         {1. / 2., 1. / 2., 0,       0},
4689371c9d4SSatish Balay         {1. / 6., 1. / 6., 1. / 6., 0}
4699371c9d4SSatish Balay     },
4709371c9d4SSatish Balay                     Gamma[4][4] = {{1. / 2., 0, 0, 0}, {0.0, 1. / 4., 0, 0}, {-2., -2. / 3., 2. / 3., 0}, {1. / 2., 5. / 36., -2. / 9, 0}}, b[4] = {1. / 6., 1. / 6., 1. / 6., 1. / 2.}, b2[4] = {1. / 8., 3. / 4., 1. / 8., 0};
471c0cb691aSEmil Constantinescu     PetscReal binterpt[4][3];
472da80777bSKarl Rupp 
473c0cb691aSEmil Constantinescu     binterpt[0][0] = 6.25;
474c0cb691aSEmil Constantinescu     binterpt[1][0] = -30.25;
475c0cb691aSEmil Constantinescu     binterpt[2][0] = 1.75;
476c0cb691aSEmil Constantinescu     binterpt[3][0] = 23.25;
477c0cb691aSEmil Constantinescu     binterpt[0][1] = -9.75;
478c0cb691aSEmil Constantinescu     binterpt[1][1] = 58.75;
479c0cb691aSEmil Constantinescu     binterpt[2][1] = -3.25;
480c0cb691aSEmil Constantinescu     binterpt[3][1] = -45.75;
481c0cb691aSEmil Constantinescu     binterpt[0][2] = 3.6666666666666666666666666666667;
482c0cb691aSEmil Constantinescu     binterpt[1][2] = -28.333333333333333333333333333333;
483c0cb691aSEmil Constantinescu     binterpt[2][2] = 1.6666666666666666666666666666667;
484c0cb691aSEmil Constantinescu     binterpt[3][2] = 23.;
485bbd56ea5SKarl Rupp 
4869566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWLASSP3P4S2C, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0]));
487b1c69cc3SEmil Constantinescu   }
488b1c69cc3SEmil Constantinescu 
489b1c69cc3SEmil Constantinescu   {
4909371c9d4SSatish Balay     const PetscReal A[4][4] =
4919371c9d4SSatish Balay       {
4929371c9d4SSatish Balay         {0,       0,       0,       0},
493b1c69cc3SEmil Constantinescu         {1. / 2., 0,       0,       0},
494b1c69cc3SEmil Constantinescu         {1. / 2., 1. / 2., 0,       0},
4959371c9d4SSatish Balay         {1. / 6., 1. / 6., 1. / 6., 0}
4969371c9d4SSatish Balay     },
4979371c9d4SSatish Balay                     Gamma[4][4] = {{1. / 2., 0, 0, 0}, {0.0, 3. / 4., 0, 0}, {-2. / 3., -23. / 9., 2. / 9., 0}, {1. / 18., 65. / 108., -2. / 27, 0}}, b[4] = {1. / 6., 1. / 6., 1. / 6., 1. / 2.}, b2[4] = {3. / 16., 10. / 16., 3. / 16., 0};
498c0cb691aSEmil Constantinescu     PetscReal binterpt[4][3];
499da80777bSKarl Rupp 
500c0cb691aSEmil Constantinescu     binterpt[0][0] = 1.6911764705882352941176470588235;
501c0cb691aSEmil Constantinescu     binterpt[1][0] = 3.6813725490196078431372549019608;
502c0cb691aSEmil Constantinescu     binterpt[2][0] = 0.23039215686274509803921568627451;
503c0cb691aSEmil Constantinescu     binterpt[3][0] = -4.6029411764705882352941176470588;
504c0cb691aSEmil Constantinescu     binterpt[0][1] = -0.95588235294117647058823529411765;
505c0cb691aSEmil Constantinescu     binterpt[1][1] = -6.2401960784313725490196078431373;
506c0cb691aSEmil Constantinescu     binterpt[2][1] = -0.31862745098039215686274509803922;
507c0cb691aSEmil Constantinescu     binterpt[3][1] = 7.5147058823529411764705882352941;
508c0cb691aSEmil Constantinescu     binterpt[0][2] = -0.56862745098039215686274509803922;
509c0cb691aSEmil Constantinescu     binterpt[1][2] = 2.7254901960784313725490196078431;
510c0cb691aSEmil Constantinescu     binterpt[2][2] = 0.25490196078431372549019607843137;
511c0cb691aSEmil Constantinescu     binterpt[3][2] = -2.4117647058823529411764705882353;
512bbd56ea5SKarl Rupp 
5139566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWLLSSP3P4S2C, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0]));
514b1c69cc3SEmil Constantinescu   }
515753f8adbSEmil Constantinescu 
516753f8adbSEmil Constantinescu   {
517753f8adbSEmil Constantinescu     PetscReal A[4][4], Gamma[4][4], b[4], b2[4];
5183ca35412SEmil Constantinescu     PetscReal binterpt[4][3];
519753f8adbSEmil Constantinescu 
520753f8adbSEmil Constantinescu     Gamma[0][0] = 0.4358665215084589994160194475295062513822671686978816;
5219371c9d4SSatish Balay     Gamma[0][1] = 0;
5229371c9d4SSatish Balay     Gamma[0][2] = 0;
5239371c9d4SSatish Balay     Gamma[0][3] = 0;
524753f8adbSEmil Constantinescu     Gamma[1][0] = -1.997527830934941248426324674704153457289527280554476;
525753f8adbSEmil Constantinescu     Gamma[1][1] = 0.4358665215084589994160194475295062513822671686978816;
5269371c9d4SSatish Balay     Gamma[1][2] = 0;
5279371c9d4SSatish Balay     Gamma[1][3] = 0;
528753f8adbSEmil Constantinescu     Gamma[2][0] = -1.007948511795029620852002345345404191008352770119903;
529753f8adbSEmil Constantinescu     Gamma[2][1] = -0.004648958462629345562774289390054679806993396798458131;
530753f8adbSEmil Constantinescu     Gamma[2][2] = 0.4358665215084589994160194475295062513822671686978816;
53105e8e825SJed Brown     Gamma[2][3] = 0;
532753f8adbSEmil Constantinescu     Gamma[3][0] = -0.6685429734233467180451604600279552604364311322650783;
533753f8adbSEmil Constantinescu     Gamma[3][1] = 0.6056625986449338476089525334450053439525178740492984;
534753f8adbSEmil Constantinescu     Gamma[3][2] = -0.9717899277217721234705114616271378792182450260943198;
535753f8adbSEmil Constantinescu     Gamma[3][3] = 0;
536753f8adbSEmil Constantinescu 
5379371c9d4SSatish Balay     A[0][0] = 0;
5389371c9d4SSatish Balay     A[0][1] = 0;
5399371c9d4SSatish Balay     A[0][2] = 0;
5409371c9d4SSatish Balay     A[0][3] = 0;
541753f8adbSEmil Constantinescu     A[1][0] = 0.8717330430169179988320388950590125027645343373957631;
5429371c9d4SSatish Balay     A[1][1] = 0;
5439371c9d4SSatish Balay     A[1][2] = 0;
5449371c9d4SSatish Balay     A[1][3] = 0;
545753f8adbSEmil Constantinescu     A[2][0] = 0.5275890119763004115618079766722914408876108660811028;
546753f8adbSEmil Constantinescu     A[2][1] = 0.07241098802369958843819203208518599088698057726988732;
5479371c9d4SSatish Balay     A[2][2] = 0;
5489371c9d4SSatish Balay     A[2][3] = 0;
549753f8adbSEmil Constantinescu     A[3][0] = 0.3990960076760701320627260685975778145384666450351314;
550753f8adbSEmil Constantinescu     A[3][1] = -0.4375576546135194437228463747348862825846903771419953;
551753f8adbSEmil Constantinescu     A[3][2] = 1.038461646937449311660120300601880176655352737312713;
55205e8e825SJed Brown     A[3][3] = 0;
553753f8adbSEmil Constantinescu 
554753f8adbSEmil Constantinescu     b[0] = 0.1876410243467238251612921333138006734899663569186926;
555753f8adbSEmil Constantinescu     b[1] = -0.5952974735769549480478230473706443582188442040780541;
556753f8adbSEmil Constantinescu     b[2] = 0.9717899277217721234705114616271378792182450260943198;
557753f8adbSEmil Constantinescu     b[3] = 0.4358665215084589994160194475295062513822671686978816;
558753f8adbSEmil Constantinescu 
559753f8adbSEmil Constantinescu     b2[0] = 0.2147402862233891404862383521089097657790734483804460;
560753f8adbSEmil Constantinescu     b2[1] = -0.4851622638849390928209050538171743017757490232519684;
561753f8adbSEmil Constantinescu     b2[2] = 0.8687250025203875511662123688667549217531982787600080;
562753f8adbSEmil Constantinescu     b2[3] = 0.4016969751411624011684543450940068201770721128357014;
563753f8adbSEmil Constantinescu 
5643ca35412SEmil Constantinescu     binterpt[0][0] = 2.2565812720167954547104627844105;
5653ca35412SEmil Constantinescu     binterpt[1][0] = 1.349166413351089573796243820819;
5663ca35412SEmil Constantinescu     binterpt[2][0] = -2.4695174540533503758652847586647;
5673ca35412SEmil Constantinescu     binterpt[3][0] = -0.13623023131453465264142184656474;
5683ca35412SEmil Constantinescu     binterpt[0][1] = -3.0826699111559187902922463354557;
5693ca35412SEmil Constantinescu     binterpt[1][1] = -2.4689115685996042534544925650515;
5703ca35412SEmil Constantinescu     binterpt[2][1] = 5.7428279814696677152129332773553;
5713ca35412SEmil Constantinescu     binterpt[3][1] = -0.19124650171414467146619437684812;
5723ca35412SEmil Constantinescu     binterpt[0][2] = 1.0137296634858471607430756831148;
5733ca35412SEmil Constantinescu     binterpt[1][2] = 0.52444768167155973161042570784064;
5743ca35412SEmil Constantinescu     binterpt[2][2] = -2.3015205996945452158771370439586;
5753ca35412SEmil Constantinescu     binterpt[3][2] = 0.76334325453713832352363565300308;
576f4aed992SEmil Constantinescu 
5779566063dSJacob Faibussowitsch     PetscCall(TSRosWRegister(TSROSWARK3, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0]));
578753f8adbSEmil Constantinescu   }
5799566063dSJacob Faibussowitsch   PetscCall(TSRosWRegisterRos4(TSROSWGRK4T, 0.231, PETSC_DEFAULT, PETSC_DEFAULT, 0, -0.1282612945269037e+01));
5809566063dSJacob Faibussowitsch   PetscCall(TSRosWRegisterRos4(TSROSWSHAMP4, 0.5, PETSC_DEFAULT, PETSC_DEFAULT, 0, 125. / 108.));
5819566063dSJacob Faibussowitsch   PetscCall(TSRosWRegisterRos4(TSROSWVELDD4, 0.22570811482256823492, PETSC_DEFAULT, PETSC_DEFAULT, 0, -1.355958941201148));
5829566063dSJacob Faibussowitsch   PetscCall(TSRosWRegisterRos4(TSROSW4L, 0.57282, PETSC_DEFAULT, PETSC_DEFAULT, 0, -1.093502252409163));
583e27a552bSJed Brown   PetscFunctionReturn(0);
584e27a552bSJed Brown }
585e27a552bSJed Brown 
586e27a552bSJed Brown /*@C
587*bcf0153eSBarry Smith    TSRosWRegisterDestroy - Frees the list of schemes that were registered by `TSRosWRegister()`.
588e27a552bSJed Brown 
589e27a552bSJed Brown    Not Collective
590e27a552bSJed Brown 
591e27a552bSJed Brown    Level: advanced
592e27a552bSJed Brown 
593*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSRosWRegister()`, `TSRosWRegisterAll()`
594e27a552bSJed Brown @*/
595d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWRegisterDestroy(void)
596d71ae5a4SJacob Faibussowitsch {
59761692a83SJed Brown   RosWTableauLink link;
598e27a552bSJed Brown 
599e27a552bSJed Brown   PetscFunctionBegin;
60061692a83SJed Brown   while ((link = RosWTableauList)) {
60161692a83SJed Brown     RosWTableau t   = &link->tab;
60261692a83SJed Brown     RosWTableauList = link->next;
6039566063dSJacob Faibussowitsch     PetscCall(PetscFree5(t->A, t->Gamma, t->b, t->ASum, t->GammaSum));
6049566063dSJacob Faibussowitsch     PetscCall(PetscFree5(t->At, t->bt, t->GammaInv, t->GammaZeroDiag, t->GammaExplicitCorr));
6059566063dSJacob Faibussowitsch     PetscCall(PetscFree2(t->bembed, t->bembedt));
6069566063dSJacob Faibussowitsch     PetscCall(PetscFree(t->binterpt));
6079566063dSJacob Faibussowitsch     PetscCall(PetscFree(t->name));
6089566063dSJacob Faibussowitsch     PetscCall(PetscFree(link));
609e27a552bSJed Brown   }
610e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_FALSE;
611e27a552bSJed Brown   PetscFunctionReturn(0);
612e27a552bSJed Brown }
613e27a552bSJed Brown 
614e27a552bSJed Brown /*@C
615*bcf0153eSBarry Smith   TSRosWInitializePackage - This function initializes everything in the `TSROSW` package. It is called
616*bcf0153eSBarry Smith   from `TSInitializePackage()`.
617e27a552bSJed Brown 
618e27a552bSJed Brown   Level: developer
619e27a552bSJed Brown 
620*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSROSW`, `PetscInitialize()`, `TSRosWFinalizePackage()`
621e27a552bSJed Brown @*/
622d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWInitializePackage(void)
623d71ae5a4SJacob Faibussowitsch {
624e27a552bSJed Brown   PetscFunctionBegin;
625e27a552bSJed Brown   if (TSRosWPackageInitialized) PetscFunctionReturn(0);
626e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_TRUE;
6279566063dSJacob Faibussowitsch   PetscCall(TSRosWRegisterAll());
6289566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(TSRosWFinalizePackage));
629e27a552bSJed Brown   PetscFunctionReturn(0);
630e27a552bSJed Brown }
631e27a552bSJed Brown 
632e27a552bSJed Brown /*@C
633*bcf0153eSBarry Smith   TSRosWFinalizePackage - This function destroys everything in the `TSROSW` package. It is
634*bcf0153eSBarry Smith   called from `PetscFinalize()`.
635e27a552bSJed Brown 
636e27a552bSJed Brown   Level: developer
637e27a552bSJed Brown 
638*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSROSW`, `PetscFinalize()`, `TSRosWInitializePackage()`
639e27a552bSJed Brown @*/
640d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWFinalizePackage(void)
641d71ae5a4SJacob Faibussowitsch {
642e27a552bSJed Brown   PetscFunctionBegin;
643e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_FALSE;
6449566063dSJacob Faibussowitsch   PetscCall(TSRosWRegisterDestroy());
645e27a552bSJed Brown   PetscFunctionReturn(0);
646e27a552bSJed Brown }
647e27a552bSJed Brown 
648e27a552bSJed Brown /*@C
649*bcf0153eSBarry Smith    TSRosWRegister - register a `TSROSW`, Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
650e27a552bSJed Brown 
651e27a552bSJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
652e27a552bSJed Brown 
653e27a552bSJed Brown    Input Parameters:
654e27a552bSJed Brown +  name - identifier for method
655e27a552bSJed Brown .  order - approximation order of method
656e27a552bSJed Brown .  s - number of stages, this is the dimension of the matrices below
65761692a83SJed Brown .  A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular
65861692a83SJed Brown .  Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal
659fe7e6d57SJed Brown .  b - Step completion table (dimension s)
6600298fd71SBarry Smith .  bembed - Step completion table for a scheme of order one less (dimension s, NULL if no embedded scheme is available)
661f4aed992SEmil Constantinescu .  pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt
66242faf41dSJed Brown -  binterpt - Coefficients of the interpolation formula (dimension s*pinterp)
663e27a552bSJed Brown 
664e27a552bSJed Brown    Level: advanced
665e27a552bSJed Brown 
666*bcf0153eSBarry Smith    Note:
667*bcf0153eSBarry Smith    Several Rosenbrock W methods are provided, this function is only needed to create new methods.
668*bcf0153eSBarry Smith 
669*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSROSW`
670e27a552bSJed Brown @*/
671d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWRegister(TSRosWType name, PetscInt order, PetscInt s, const PetscReal A[], const PetscReal Gamma[], const PetscReal b[], const PetscReal bembed[], PetscInt pinterp, const PetscReal binterpt[])
672d71ae5a4SJacob Faibussowitsch {
67361692a83SJed Brown   RosWTableauLink link;
67461692a83SJed Brown   RosWTableau     t;
67561692a83SJed Brown   PetscInt        i, j, k;
67661692a83SJed Brown   PetscScalar    *GammaInv;
677e27a552bSJed Brown 
678e27a552bSJed Brown   PetscFunctionBegin;
679fe7e6d57SJed Brown   PetscValidCharPointer(name, 1);
680dadcf809SJacob Faibussowitsch   PetscValidRealPointer(A, 4);
681dadcf809SJacob Faibussowitsch   PetscValidRealPointer(Gamma, 5);
682dadcf809SJacob Faibussowitsch   PetscValidRealPointer(b, 6);
683dadcf809SJacob Faibussowitsch   if (bembed) PetscValidRealPointer(bembed, 7);
684fe7e6d57SJed Brown 
6859566063dSJacob Faibussowitsch   PetscCall(TSRosWInitializePackage());
6869566063dSJacob Faibussowitsch   PetscCall(PetscNew(&link));
687e27a552bSJed Brown   t = &link->tab;
6889566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &t->name));
689e27a552bSJed Brown   t->order = order;
690e27a552bSJed Brown   t->s     = s;
6919566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(s * s, &t->A, s * s, &t->Gamma, s, &t->b, s, &t->ASum, s, &t->GammaSum));
6929566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(s * s, &t->At, s, &t->bt, s * s, &t->GammaInv, s, &t->GammaZeroDiag, s * s, &t->GammaExplicitCorr));
6939566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->A, A, s * s));
6949566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->Gamma, Gamma, s * s));
6959566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->GammaExplicitCorr, Gamma, s * s));
6969566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->b, b, s));
697fe7e6d57SJed Brown   if (bembed) {
6989566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(s, &t->bembed, s, &t->bembedt));
6999566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(t->bembed, bembed, s));
700fe7e6d57SJed Brown   }
70161692a83SJed Brown   for (i = 0; i < s; i++) {
70261692a83SJed Brown     t->ASum[i]     = 0;
70361692a83SJed Brown     t->GammaSum[i] = 0;
70461692a83SJed Brown     for (j = 0; j < s; j++) {
70561692a83SJed Brown       t->ASum[i] += A[i * s + j];
706fe7e6d57SJed Brown       t->GammaSum[i] += Gamma[i * s + j];
70761692a83SJed Brown     }
70861692a83SJed Brown   }
7099566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s * s, &GammaInv)); /* Need to use Scalar for inverse, then convert back to Real */
71061692a83SJed Brown   for (i = 0; i < s * s; i++) GammaInv[i] = Gamma[i];
711fd96d5b0SEmil Constantinescu   for (i = 0; i < s; i++) {
712fd96d5b0SEmil Constantinescu     if (Gamma[i * s + i] == 0.0) {
713fd96d5b0SEmil Constantinescu       GammaInv[i * s + i] = 1.0;
714c17803e7SJed Brown       t->GammaZeroDiag[i] = PETSC_TRUE;
715fd96d5b0SEmil Constantinescu     } else {
716c17803e7SJed Brown       t->GammaZeroDiag[i] = PETSC_FALSE;
717fd96d5b0SEmil Constantinescu     }
718fd96d5b0SEmil Constantinescu   }
719fd96d5b0SEmil Constantinescu 
72061692a83SJed Brown   switch (s) {
721d71ae5a4SJacob Faibussowitsch   case 1:
722d71ae5a4SJacob Faibussowitsch     GammaInv[0] = 1. / GammaInv[0];
723d71ae5a4SJacob Faibussowitsch     break;
724d71ae5a4SJacob Faibussowitsch   case 2:
725d71ae5a4SJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_2(GammaInv, 0, PETSC_FALSE, NULL));
726d71ae5a4SJacob Faibussowitsch     break;
727d71ae5a4SJacob Faibussowitsch   case 3:
728d71ae5a4SJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_3(GammaInv, 0, PETSC_FALSE, NULL));
729d71ae5a4SJacob Faibussowitsch     break;
730d71ae5a4SJacob Faibussowitsch   case 4:
731d71ae5a4SJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_4(GammaInv, 0, PETSC_FALSE, NULL));
732d71ae5a4SJacob Faibussowitsch     break;
73361692a83SJed Brown   case 5: {
73461692a83SJed Brown     PetscInt  ipvt5[5];
73561692a83SJed Brown     MatScalar work5[5 * 5];
7369371c9d4SSatish Balay     PetscCall(PetscKernel_A_gets_inverse_A_5(GammaInv, ipvt5, work5, 0, PETSC_FALSE, NULL));
7379371c9d4SSatish Balay     break;
73861692a83SJed Brown   }
739d71ae5a4SJacob Faibussowitsch   case 6:
740d71ae5a4SJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_6(GammaInv, 0, PETSC_FALSE, NULL));
741d71ae5a4SJacob Faibussowitsch     break;
742d71ae5a4SJacob Faibussowitsch   case 7:
743d71ae5a4SJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_7(GammaInv, 0, PETSC_FALSE, NULL));
744d71ae5a4SJacob Faibussowitsch     break;
745d71ae5a4SJacob Faibussowitsch   default:
746d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not implemented for %" PetscInt_FMT " stages", s);
74761692a83SJed Brown   }
74861692a83SJed Brown   for (i = 0; i < s * s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]);
7499566063dSJacob Faibussowitsch   PetscCall(PetscFree(GammaInv));
75043b21953SEmil Constantinescu 
75143b21953SEmil Constantinescu   for (i = 0; i < s; i++) {
75243b21953SEmil Constantinescu     for (k = 0; k < i + 1; k++) {
75343b21953SEmil Constantinescu       t->GammaExplicitCorr[i * s + k] = (t->GammaExplicitCorr[i * s + k]) * (t->GammaInv[k * s + k]);
754ad540459SPierre Jolivet       for (j = k + 1; j < i + 1; j++) t->GammaExplicitCorr[i * s + k] += (t->GammaExplicitCorr[i * s + j]) * (t->GammaInv[j * s + k]);
75543b21953SEmil Constantinescu     }
75643b21953SEmil Constantinescu   }
75743b21953SEmil Constantinescu 
75861692a83SJed Brown   for (i = 0; i < s; i++) {
75961692a83SJed Brown     for (j = 0; j < s; j++) {
76061692a83SJed Brown       t->At[i * s + j] = 0;
761ad540459SPierre Jolivet       for (k = 0; k < s; k++) t->At[i * s + j] += t->A[i * s + k] * t->GammaInv[k * s + j];
76261692a83SJed Brown     }
76361692a83SJed Brown     t->bt[i] = 0;
764ad540459SPierre Jolivet     for (j = 0; j < s; j++) t->bt[i] += t->b[j] * t->GammaInv[j * s + i];
765fe7e6d57SJed Brown     if (bembed) {
766fe7e6d57SJed Brown       t->bembedt[i] = 0;
767ad540459SPierre Jolivet       for (j = 0; j < s; j++) t->bembedt[i] += t->bembed[j] * t->GammaInv[j * s + i];
768fe7e6d57SJed Brown     }
76961692a83SJed Brown   }
7708d59e960SJed Brown   t->ccfl = 1.0; /* Fix this */
7718d59e960SJed Brown 
772f4aed992SEmil Constantinescu   t->pinterp = pinterp;
7739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s * pinterp, &t->binterpt));
7749566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->binterpt, binterpt, s * pinterp));
77561692a83SJed Brown   link->next      = RosWTableauList;
77661692a83SJed Brown   RosWTableauList = link;
777e27a552bSJed Brown   PetscFunctionReturn(0);
778e27a552bSJed Brown }
779e27a552bSJed Brown 
78042faf41dSJed Brown /*@C
781fd292e60Sprj-    TSRosWRegisterRos4 - register a fourth order Rosenbrock scheme by providing parameter choices
78242faf41dSJed Brown 
78342faf41dSJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
78442faf41dSJed Brown 
78542faf41dSJed Brown    Input Parameters:
78642faf41dSJed Brown +  name - identifier for method
78742faf41dSJed Brown .  gamma - leading coefficient (diagonal entry)
78842faf41dSJed Brown .  a2 - design parameter, see Table 7.2 of Hairer&Wanner
78942faf41dSJed Brown .  a3 - design parameter or PETSC_DEFAULT to satisfy one of the order five conditions (Eq 7.22)
79042faf41dSJed Brown .  b3 - design parameter, see Table 7.2 of Hairer&Wanner
79142faf41dSJed Brown .  beta43 - design parameter or PETSC_DEFAULT to use Equation 7.21 of Hairer&Wanner
792a2b725a8SWilliam Gropp -  e4 - design parameter for embedded method, see coefficient E4 in ros4.f code from Hairer
79342faf41dSJed Brown 
794*bcf0153eSBarry Smith    Level: developer
795*bcf0153eSBarry Smith 
79642faf41dSJed Brown    Notes:
79742faf41dSJed Brown    This routine encodes the design of fourth order Rosenbrock methods as described in Hairer and Wanner volume 2.
79842faf41dSJed Brown    It is used here to implement several methods from the book and can be used to experiment with new methods.
79942faf41dSJed Brown    It was written this way instead of by copying coefficients in order to provide better than double precision satisfaction of the order conditions.
80042faf41dSJed Brown 
801*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSRosW`, `TSRosWRegister()`
80242faf41dSJed Brown @*/
803d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWRegisterRos4(TSRosWType name, PetscReal gamma, PetscReal a2, PetscReal a3, PetscReal b3, PetscReal e4)
804d71ae5a4SJacob Faibussowitsch {
80542faf41dSJed Brown   /* Declare numeric constants so they can be quad precision without being truncated at double */
8069371c9d4SSatish 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;
80742faf41dSJed Brown   PetscReal   a4, a32, a42, a43, b1, b2, b4, beta2p, beta3p, beta4p, beta32, beta42, beta43, beta32beta2p, beta4jbetajp;
80842faf41dSJed Brown   PetscReal   A[4][4], Gamma[4][4], b[4], bm[4];
80942faf41dSJed Brown   PetscScalar M[3][3], rhs[3];
81042faf41dSJed Brown 
81142faf41dSJed Brown   PetscFunctionBegin;
81242faf41dSJed Brown   /* Step 1: choose Gamma (input) */
81342faf41dSJed Brown   /* Step 2: choose a2,a3,a4; b1,b2,b3,b4 to satisfy order conditions */
81442faf41dSJed Brown   if (a3 == PETSC_DEFAULT) a3 = (one / five - a2 / four) / (one / four - a2 / three); /* Eq 7.22 */
81542faf41dSJed Brown   a4 = a3;                                                                            /* consequence of 7.20 */
81642faf41dSJed Brown 
81742faf41dSJed Brown   /* Solve order conditions 7.15a, 7.15c, 7.15e */
8189371c9d4SSatish Balay   M[0][0] = one;
8199371c9d4SSatish Balay   M[0][1] = one;
8209371c9d4SSatish Balay   M[0][2] = one; /* 7.15a */
8219371c9d4SSatish Balay   M[1][0] = 0.0;
8229371c9d4SSatish Balay   M[1][1] = a2 * a2;
8239371c9d4SSatish Balay   M[1][2] = a4 * a4; /* 7.15c */
8249371c9d4SSatish Balay   M[2][0] = 0.0;
8259371c9d4SSatish Balay   M[2][1] = a2 * a2 * a2;
8269371c9d4SSatish Balay   M[2][2] = a4 * a4 * a4; /* 7.15e */
82742faf41dSJed Brown   rhs[0]  = one - b3;
82842faf41dSJed Brown   rhs[1]  = one / three - a3 * a3 * b3;
82942faf41dSJed Brown   rhs[2]  = one / four - a3 * a3 * a3 * b3;
8309566063dSJacob Faibussowitsch   PetscCall(PetscKernel_A_gets_inverse_A_3(&M[0][0], 0, PETSC_FALSE, NULL));
83142faf41dSJed Brown   b1 = PetscRealPart(M[0][0] * rhs[0] + M[0][1] * rhs[1] + M[0][2] * rhs[2]);
83242faf41dSJed Brown   b2 = PetscRealPart(M[1][0] * rhs[0] + M[1][1] * rhs[1] + M[1][2] * rhs[2]);
83342faf41dSJed Brown   b4 = PetscRealPart(M[2][0] * rhs[0] + M[2][1] * rhs[1] + M[2][2] * rhs[2]);
83442faf41dSJed Brown 
83542faf41dSJed Brown   /* Step 3 */
83642faf41dSJed Brown   beta43       = (p56 - a2 * p43) / (b4 * a3 * a3 * (a3 - a2)); /* 7.21 */
83742faf41dSJed Brown   beta32beta2p = p44 / (b4 * beta43);                           /* 7.15h */
83842faf41dSJed Brown   beta4jbetajp = (p32 - b3 * beta32beta2p) / b4;
8399371c9d4SSatish Balay   M[0][0]      = b2;
8409371c9d4SSatish Balay   M[0][1]      = b3;
8419371c9d4SSatish Balay   M[0][2]      = b4;
8429371c9d4SSatish Balay   M[1][0]      = a4 * a4 * beta32beta2p - a3 * a3 * beta4jbetajp;
8439371c9d4SSatish Balay   M[1][1]      = a2 * a2 * beta4jbetajp;
8449371c9d4SSatish Balay   M[1][2]      = -a2 * a2 * beta32beta2p;
8459371c9d4SSatish Balay   M[2][0]      = b4 * beta43 * a3 * a3 - p43;
8469371c9d4SSatish Balay   M[2][1]      = -b4 * beta43 * a2 * a2;
8479371c9d4SSatish Balay   M[2][2]      = 0;
8489371c9d4SSatish Balay   rhs[0]       = one / two - gamma;
8499371c9d4SSatish Balay   rhs[1]       = 0;
8509371c9d4SSatish Balay   rhs[2]       = -a2 * a2 * p32;
8519566063dSJacob Faibussowitsch   PetscCall(PetscKernel_A_gets_inverse_A_3(&M[0][0], 0, PETSC_FALSE, NULL));
85242faf41dSJed Brown   beta2p = PetscRealPart(M[0][0] * rhs[0] + M[0][1] * rhs[1] + M[0][2] * rhs[2]);
85342faf41dSJed Brown   beta3p = PetscRealPart(M[1][0] * rhs[0] + M[1][1] * rhs[1] + M[1][2] * rhs[2]);
85442faf41dSJed Brown   beta4p = PetscRealPart(M[2][0] * rhs[0] + M[2][1] * rhs[1] + M[2][2] * rhs[2]);
85542faf41dSJed Brown 
85642faf41dSJed Brown   /* Step 4: back-substitute */
85742faf41dSJed Brown   beta32 = beta32beta2p / beta2p;
85842faf41dSJed Brown   beta42 = (beta4jbetajp - beta43 * beta3p) / beta2p;
85942faf41dSJed Brown 
86042faf41dSJed Brown   /* Step 5: 7.15f and 7.20, then 7.16 */
86142faf41dSJed Brown   a43 = 0;
86242faf41dSJed Brown   a32 = p42 / (b3 * a3 * beta2p + b4 * a4 * beta2p);
86342faf41dSJed Brown   a42 = a32;
86442faf41dSJed Brown 
8659371c9d4SSatish Balay   A[0][0]     = 0;
8669371c9d4SSatish Balay   A[0][1]     = 0;
8679371c9d4SSatish Balay   A[0][2]     = 0;
8689371c9d4SSatish Balay   A[0][3]     = 0;
8699371c9d4SSatish Balay   A[1][0]     = a2;
8709371c9d4SSatish Balay   A[1][1]     = 0;
8719371c9d4SSatish Balay   A[1][2]     = 0;
8729371c9d4SSatish Balay   A[1][3]     = 0;
8739371c9d4SSatish Balay   A[2][0]     = a3 - a32;
8749371c9d4SSatish Balay   A[2][1]     = a32;
8759371c9d4SSatish Balay   A[2][2]     = 0;
8769371c9d4SSatish Balay   A[2][3]     = 0;
8779371c9d4SSatish Balay   A[3][0]     = a4 - a43 - a42;
8789371c9d4SSatish Balay   A[3][1]     = a42;
8799371c9d4SSatish Balay   A[3][2]     = a43;
8809371c9d4SSatish Balay   A[3][3]     = 0;
8819371c9d4SSatish Balay   Gamma[0][0] = gamma;
8829371c9d4SSatish Balay   Gamma[0][1] = 0;
8839371c9d4SSatish Balay   Gamma[0][2] = 0;
8849371c9d4SSatish Balay   Gamma[0][3] = 0;
8859371c9d4SSatish Balay   Gamma[1][0] = beta2p - A[1][0];
8869371c9d4SSatish Balay   Gamma[1][1] = gamma;
8879371c9d4SSatish Balay   Gamma[1][2] = 0;
8889371c9d4SSatish Balay   Gamma[1][3] = 0;
8899371c9d4SSatish Balay   Gamma[2][0] = beta3p - beta32 - A[2][0];
8909371c9d4SSatish Balay   Gamma[2][1] = beta32 - A[2][1];
8919371c9d4SSatish Balay   Gamma[2][2] = gamma;
8929371c9d4SSatish Balay   Gamma[2][3] = 0;
8939371c9d4SSatish Balay   Gamma[3][0] = beta4p - beta42 - beta43 - A[3][0];
8949371c9d4SSatish Balay   Gamma[3][1] = beta42 - A[3][1];
8959371c9d4SSatish Balay   Gamma[3][2] = beta43 - A[3][2];
8969371c9d4SSatish Balay   Gamma[3][3] = gamma;
8979371c9d4SSatish Balay   b[0]        = b1;
8989371c9d4SSatish Balay   b[1]        = b2;
8999371c9d4SSatish Balay   b[2]        = b3;
9009371c9d4SSatish Balay   b[3]        = b4;
90142faf41dSJed Brown 
90242faf41dSJed Brown   /* Construct embedded formula using given e4. We are solving Equation 7.18. */
90342faf41dSJed Brown   bm[3] = b[3] - e4 * gamma;                                              /* using definition of E4 */
90442faf41dSJed Brown   bm[2] = (p32 - beta4jbetajp * bm[3]) / (beta32 * beta2p);               /* fourth row of 7.18 */
90542faf41dSJed Brown   bm[1] = (one / two - gamma - beta3p * bm[2] - beta4p * bm[3]) / beta2p; /* second row */
90642faf41dSJed Brown   bm[0] = one - bm[1] - bm[2] - bm[3];                                    /* first row */
90742faf41dSJed Brown 
90842faf41dSJed Brown   {
90942faf41dSJed Brown     const PetscReal misfit = a2 * a2 * bm[1] + a3 * a3 * bm[2] + a4 * a4 * bm[3] - one / three;
9103c633725SBarry Smith     PetscCheck(PetscAbs(misfit) <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_SUP, "Assumptions violated, could not construct a third order embedded method");
91142faf41dSJed Brown   }
9129566063dSJacob Faibussowitsch   PetscCall(TSRosWRegister(name, 4, 4, &A[0][0], &Gamma[0][0], b, bm, 0, NULL));
91342faf41dSJed Brown   PetscFunctionReturn(0);
91442faf41dSJed Brown }
91542faf41dSJed Brown 
9161c3436cfSJed Brown /*
9171c3436cfSJed Brown  The step completion formula is
9181c3436cfSJed Brown 
9191c3436cfSJed Brown  x1 = x0 + b^T Y
9201c3436cfSJed Brown 
9211c3436cfSJed Brown  where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been
9221c3436cfSJed Brown  updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write
9231c3436cfSJed Brown 
9241c3436cfSJed Brown  x1e = x0 + be^T Y
9251c3436cfSJed Brown      = x1 - b^T Y + be^T Y
9261c3436cfSJed Brown      = x1 + (be - b)^T Y
9271c3436cfSJed Brown 
9281c3436cfSJed Brown  so we can evaluate the method of different order even after the step has been optimistically completed.
9291c3436cfSJed Brown */
930d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateStep_RosW(TS ts, PetscInt order, Vec U, PetscBool *done)
931d71ae5a4SJacob Faibussowitsch {
9321c3436cfSJed Brown   TS_RosW     *ros = (TS_RosW *)ts->data;
9331c3436cfSJed Brown   RosWTableau  tab = ros->tableau;
9341c3436cfSJed Brown   PetscScalar *w   = ros->work;
9351c3436cfSJed Brown   PetscInt     i;
9361c3436cfSJed Brown 
9371c3436cfSJed Brown   PetscFunctionBegin;
9381c3436cfSJed Brown   if (order == tab->order) {
939108c343cSJed Brown     if (ros->status == TS_STEP_INCOMPLETE) { /* Use standard completion formula */
9409566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol, U));
941de19f811SJed Brown       for (i = 0; i < tab->s; i++) w[i] = tab->bt[i];
9429566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(U, tab->s, w, ros->Y));
9439566063dSJacob Faibussowitsch     } else PetscCall(VecCopy(ts->vec_sol, U));
9441c3436cfSJed Brown     if (done) *done = PETSC_TRUE;
9451c3436cfSJed Brown     PetscFunctionReturn(0);
9461c3436cfSJed Brown   } else if (order == tab->order - 1) {
9471c3436cfSJed Brown     if (!tab->bembedt) goto unavailable;
948108c343cSJed Brown     if (ros->status == TS_STEP_INCOMPLETE) { /* Use embedded completion formula */
9499566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol, U));
950de19f811SJed Brown       for (i = 0; i < tab->s; i++) w[i] = tab->bembedt[i];
9519566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(U, tab->s, w, ros->Y));
952108c343cSJed Brown     } else { /* Use rollback-and-recomplete formula (bembedt - bt) */
953108c343cSJed Brown       for (i = 0; i < tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i];
9549566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol, U));
9559566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(U, tab->s, w, ros->Y));
9561c3436cfSJed Brown     }
9571c3436cfSJed Brown     if (done) *done = PETSC_TRUE;
9581c3436cfSJed Brown     PetscFunctionReturn(0);
9591c3436cfSJed Brown   }
9601c3436cfSJed Brown unavailable:
9611c3436cfSJed Brown   if (done) *done = PETSC_FALSE;
9629371c9d4SSatish Balay   else
9639371c9d4SSatish 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,
9649371c9d4SSatish Balay             tab->order, order);
9651c3436cfSJed Brown   PetscFunctionReturn(0);
9661c3436cfSJed Brown }
9671c3436cfSJed Brown 
968d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRollBack_RosW(TS ts)
969d71ae5a4SJacob Faibussowitsch {
97024655328SShri   TS_RosW *ros = (TS_RosW *)ts->data;
97124655328SShri 
97224655328SShri   PetscFunctionBegin;
9739566063dSJacob Faibussowitsch   PetscCall(VecCopy(ros->vec_sol_prev, ts->vec_sol));
97424655328SShri   PetscFunctionReturn(0);
97524655328SShri }
97624655328SShri 
977d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_RosW(TS ts)
978d71ae5a4SJacob Faibussowitsch {
97961692a83SJed Brown   TS_RosW         *ros = (TS_RosW *)ts->data;
98061692a83SJed Brown   RosWTableau      tab = ros->tableau;
981e27a552bSJed Brown   const PetscInt   s   = tab->s;
9821c3436cfSJed Brown   const PetscReal *At = tab->At, *Gamma = tab->Gamma, *ASum = tab->ASum, *GammaInv = tab->GammaInv;
9830feba352SEmil Constantinescu   const PetscReal *GammaExplicitCorr = tab->GammaExplicitCorr;
984c17803e7SJed Brown   const PetscBool *GammaZeroDiag     = tab->GammaZeroDiag;
98561692a83SJed Brown   PetscScalar     *w                 = ros->work;
9867d4bf2deSEmil Constantinescu   Vec             *Y = ros->Y, Ydot = ros->Ydot, Zdot = ros->Zdot, Zstage = ros->Zstage;
987e27a552bSJed Brown   SNES             snes;
9881c3436cfSJed Brown   TSAdapt          adapt;
989fecfb714SLisandro Dalcin   PetscInt         i, j, its, lits;
990be5899b3SLisandro Dalcin   PetscInt         rejections = 0;
991b39943a6SLisandro Dalcin   PetscBool        stageok, accept = PETSC_TRUE;
992b39943a6SLisandro Dalcin   PetscReal        next_time_step = ts->time_step;
993f7f07198SBarry Smith   PetscInt         lag;
994e27a552bSJed Brown 
995e27a552bSJed Brown   PetscFunctionBegin;
99648a46eb9SPierre Jolivet   if (!ts->steprollback) PetscCall(VecCopy(ts->vec_sol, ros->vec_sol_prev));
997e27a552bSJed Brown 
998b39943a6SLisandro Dalcin   ros->status = TS_STEP_INCOMPLETE;
999b39943a6SLisandro Dalcin   while (!ts->reason && ros->status != TS_STEP_COMPLETE) {
10001c3436cfSJed Brown     const PetscReal h = ts->time_step;
1001e27a552bSJed Brown     for (i = 0; i < s; i++) {
10021c3436cfSJed Brown       ros->stage_time = ts->ptime + h * ASum[i];
10039566063dSJacob Faibussowitsch       PetscCall(TSPreStage(ts, ros->stage_time));
1004c17803e7SJed Brown       if (GammaZeroDiag[i]) {
1005c17803e7SJed Brown         ros->stage_explicit = PETSC_TRUE;
1006b296d7d5SJed Brown         ros->scoeff         = 1.;
1007c17803e7SJed Brown       } else {
1008c17803e7SJed Brown         ros->stage_explicit = PETSC_FALSE;
1009b296d7d5SJed Brown         ros->scoeff         = 1. / Gamma[i * s + i];
1010fd96d5b0SEmil Constantinescu       }
101161692a83SJed Brown 
10129566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol, Zstage));
1013de19f811SJed Brown       for (j = 0; j < i; j++) w[j] = At[i * s + j];
10149566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(Zstage, i, w, Y));
101561692a83SJed Brown 
101661692a83SJed Brown       for (j = 0; j < i; j++) w[j] = 1. / h * GammaInv[i * s + j];
10179566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(Zdot));
10189566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(Zdot, i, w, Y));
101961692a83SJed Brown 
1020e27a552bSJed Brown       /* Initial guess taken from last stage */
10219566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(Y[i]));
102261692a83SJed Brown 
10237d4bf2deSEmil Constantinescu       if (!ros->stage_explicit) {
10249566063dSJacob Faibussowitsch         PetscCall(TSGetSNES(ts, &snes));
102561692a83SJed Brown         if (!ros->recompute_jacobian && !i) {
10269566063dSJacob Faibussowitsch           PetscCall(SNESGetLagJacobian(snes, &lag));
10276aad120cSJose E. Roman           if (lag == 1) {                            /* use did not set a nontrivial lag, so lag over all stages */
10289566063dSJacob Faibussowitsch             PetscCall(SNESSetLagJacobian(snes, -2)); /* Recompute the Jacobian on this solve, but not again for the rest of the stages */
1029f7f07198SBarry Smith           }
103061692a83SJed Brown         }
10319566063dSJacob Faibussowitsch         PetscCall(SNESSolve(snes, NULL, Y[i]));
10329371c9d4SSatish 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 */ }
10339566063dSJacob Faibussowitsch         PetscCall(SNESGetIterationNumber(snes, &its));
10349566063dSJacob Faibussowitsch         PetscCall(SNESGetLinearSolveIterations(snes, &lits));
10359371c9d4SSatish Balay         ts->snes_its += its;
10369371c9d4SSatish Balay         ts->ksp_its += lits;
10377d4bf2deSEmil Constantinescu       } else {
10381ce71dffSSatish Balay         Mat J, Jp;
10399566063dSJacob Faibussowitsch         PetscCall(VecZeroEntries(Ydot)); /* Evaluate Y[i]=G(t,Ydot=0,Zstage) */
10409566063dSJacob Faibussowitsch         PetscCall(TSComputeIFunction(ts, ros->stage_time, Zstage, Ydot, Y[i], PETSC_FALSE));
10419566063dSJacob Faibussowitsch         PetscCall(VecScale(Y[i], -1.0));
10429566063dSJacob Faibussowitsch         PetscCall(VecAXPY(Y[i], -1.0, Zdot)); /*Y[i] = F(Zstage)-Zdot[=GammaInv*Y]*/
10430feba352SEmil Constantinescu 
10449566063dSJacob Faibussowitsch         PetscCall(VecZeroEntries(Zstage)); /* Zstage = GammaExplicitCorr[i,j] * Y[j] */
10450feba352SEmil Constantinescu         for (j = 0; j < i; j++) w[j] = GammaExplicitCorr[i * s + j];
10469566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(Zstage, i, w, Y));
1047fecfb714SLisandro Dalcin 
1048fecfb714SLisandro Dalcin         /* Y[i] = Y[i] + Jac*Zstage[=Jac*GammaExplicitCorr[i,j] * Y[j]] */
10499566063dSJacob Faibussowitsch         PetscCall(TSGetIJacobian(ts, &J, &Jp, NULL, NULL));
10509566063dSJacob Faibussowitsch         PetscCall(TSComputeIJacobian(ts, ros->stage_time, ts->vec_sol, Ydot, 0, J, Jp, PETSC_FALSE));
10519566063dSJacob Faibussowitsch         PetscCall(MatMult(J, Zstage, Zdot));
10529566063dSJacob Faibussowitsch         PetscCall(VecAXPY(Y[i], -1.0, Zdot));
10535ef26d82SJed Brown         ts->ksp_its += 1;
1054fecfb714SLisandro Dalcin 
10559566063dSJacob Faibussowitsch         PetscCall(VecScale(Y[i], h));
10567d4bf2deSEmil Constantinescu       }
10579566063dSJacob Faibussowitsch       PetscCall(TSPostStage(ts, ros->stage_time, i, Y));
10589566063dSJacob Faibussowitsch       PetscCall(TSGetAdapt(ts, &adapt));
10599566063dSJacob Faibussowitsch       PetscCall(TSAdaptCheckStage(adapt, ts, ros->stage_time, Y[i], &stageok));
1060fecfb714SLisandro Dalcin       if (!stageok) goto reject_step;
1061e27a552bSJed Brown     }
1062e27a552bSJed Brown 
1063b39943a6SLisandro Dalcin     ros->status = TS_STEP_INCOMPLETE;
10649566063dSJacob Faibussowitsch     PetscCall(TSEvaluateStep_RosW(ts, tab->order, ts->vec_sol, NULL));
1065b39943a6SLisandro Dalcin     ros->status = TS_STEP_PENDING;
10669566063dSJacob Faibussowitsch     PetscCall(TSGetAdapt(ts, &adapt));
10679566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidatesClear(adapt));
10689566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidateAdd(adapt, tab->name, tab->order, 1, tab->ccfl, (PetscReal)tab->s, PETSC_TRUE));
10699566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
1070b39943a6SLisandro Dalcin     ros->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
1071b39943a6SLisandro Dalcin     if (!accept) { /* Roll back the current step */
10729566063dSJacob Faibussowitsch       PetscCall(TSRollBack_RosW(ts));
1073be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
1074be5899b3SLisandro Dalcin       goto reject_step;
1075b39943a6SLisandro Dalcin     }
1076b39943a6SLisandro Dalcin 
1077e27a552bSJed Brown     ts->ptime += ts->time_step;
1078cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
10791c3436cfSJed Brown     break;
1080b39943a6SLisandro Dalcin 
1081b39943a6SLisandro Dalcin   reject_step:
10829371c9d4SSatish Balay     ts->reject++;
10839371c9d4SSatish Balay     accept = PETSC_FALSE;
1084be5899b3SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
1085b39943a6SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
108663a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
10871c3436cfSJed Brown     }
10881c3436cfSJed Brown   }
1089e27a552bSJed Brown   PetscFunctionReturn(0);
1090e27a552bSJed Brown }
1091e27a552bSJed Brown 
1092d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_RosW(TS ts, PetscReal itime, Vec U)
1093d71ae5a4SJacob Faibussowitsch {
109461692a83SJed Brown   TS_RosW         *ros = (TS_RosW *)ts->data;
1095f4aed992SEmil Constantinescu   PetscInt         s = ros->tableau->s, pinterp = ros->tableau->pinterp, i, j;
1096f4aed992SEmil Constantinescu   PetscReal        h;
1097f4aed992SEmil Constantinescu   PetscReal        tt, t;
1098f4aed992SEmil Constantinescu   PetscScalar     *bt;
1099f4aed992SEmil Constantinescu   const PetscReal *Bt       = ros->tableau->binterpt;
1100f4aed992SEmil Constantinescu   const PetscReal *GammaInv = ros->tableau->GammaInv;
1101f4aed992SEmil Constantinescu   PetscScalar     *w        = ros->work;
1102f4aed992SEmil Constantinescu   Vec             *Y        = ros->Y;
1103e27a552bSJed Brown 
1104e27a552bSJed Brown   PetscFunctionBegin;
11053c633725SBarry Smith   PetscCheck(Bt, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSRosW %s does not have an interpolation formula", ros->tableau->name);
1106f4aed992SEmil Constantinescu 
1107f4aed992SEmil Constantinescu   switch (ros->status) {
1108f4aed992SEmil Constantinescu   case TS_STEP_INCOMPLETE:
1109f4aed992SEmil Constantinescu   case TS_STEP_PENDING:
1110f4aed992SEmil Constantinescu     h = ts->time_step;
1111f4aed992SEmil Constantinescu     t = (itime - ts->ptime) / h;
1112f4aed992SEmil Constantinescu     break;
1113f4aed992SEmil Constantinescu   case TS_STEP_COMPLETE:
1114be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev;
1115f4aed992SEmil Constantinescu     t = (itime - ts->ptime) / h + 1; /* In the interval [0,1] */
1116f4aed992SEmil Constantinescu     break;
1117d71ae5a4SJacob Faibussowitsch   default:
1118d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus");
1119f4aed992SEmil Constantinescu   }
11209566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s, &bt));
1121f4aed992SEmil Constantinescu   for (i = 0; i < s; i++) bt[i] = 0;
1122f4aed992SEmil Constantinescu   for (j = 0, tt = t; j < pinterp; j++, tt *= t) {
1123ad540459SPierre Jolivet     for (i = 0; i < s; i++) bt[i] += Bt[i * pinterp + j] * tt;
1124f4aed992SEmil Constantinescu   }
1125f4aed992SEmil Constantinescu 
1126f4aed992SEmil Constantinescu   /* y(t+tt*h) = y(t) + Sum bt(tt) * GammaInv * Ydot */
1127f9c1d6abSBarry Smith   /* U <- 0*/
11289566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(U));
1129f9c1d6abSBarry Smith   /* U <- Sum bt_i * GammaInv(i,1:i) * Y(1:i) */
11303ca35412SEmil Constantinescu   for (j = 0; j < s; j++) w[j] = 0;
11313ca35412SEmil Constantinescu   for (j = 0; j < s; j++) {
1132ad540459SPierre Jolivet     for (i = j; i < s; i++) w[j] += bt[i] * GammaInv[i * s + j];
11333ca35412SEmil Constantinescu   }
11349566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(U, i, w, Y));
1135be5899b3SLisandro Dalcin   /* U <- y(t) + U */
11369566063dSJacob Faibussowitsch   PetscCall(VecAXPY(U, 1, ros->vec_sol_prev));
1137f4aed992SEmil Constantinescu 
11389566063dSJacob Faibussowitsch   PetscCall(PetscFree(bt));
1139e27a552bSJed Brown   PetscFunctionReturn(0);
1140e27a552bSJed Brown }
1141e27a552bSJed Brown 
1142e27a552bSJed Brown /*------------------------------------------------------------*/
1143b39943a6SLisandro Dalcin 
1144d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWTableauReset(TS ts)
1145d71ae5a4SJacob Faibussowitsch {
1146b39943a6SLisandro Dalcin   TS_RosW    *ros = (TS_RosW *)ts->data;
1147b39943a6SLisandro Dalcin   RosWTableau tab = ros->tableau;
1148b39943a6SLisandro Dalcin 
1149b39943a6SLisandro Dalcin   PetscFunctionBegin;
1150b39943a6SLisandro Dalcin   if (!tab) PetscFunctionReturn(0);
11519566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(tab->s, &ros->Y));
11529566063dSJacob Faibussowitsch   PetscCall(PetscFree(ros->work));
1153b39943a6SLisandro Dalcin   PetscFunctionReturn(0);
1154b39943a6SLisandro Dalcin }
1155b39943a6SLisandro Dalcin 
1156d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_RosW(TS ts)
1157d71ae5a4SJacob Faibussowitsch {
115861692a83SJed Brown   TS_RosW *ros = (TS_RosW *)ts->data;
1159e27a552bSJed Brown 
1160e27a552bSJed Brown   PetscFunctionBegin;
11619566063dSJacob Faibussowitsch   PetscCall(TSRosWTableauReset(ts));
11629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ros->Ydot));
11639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ros->Ystage));
11649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ros->Zdot));
11659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ros->Zstage));
11669566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ros->vec_sol_prev));
1167e27a552bSJed Brown   PetscFunctionReturn(0);
1168e27a552bSJed Brown }
1169e27a552bSJed Brown 
1170d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWGetVecs(TS ts, DM dm, Vec *Ydot, Vec *Zdot, Vec *Ystage, Vec *Zstage)
1171d71ae5a4SJacob Faibussowitsch {
1172d5e6173cSPeter Brune   TS_RosW *rw = (TS_RosW *)ts->data;
1173d5e6173cSPeter Brune 
1174d5e6173cSPeter Brune   PetscFunctionBegin;
1175d5e6173cSPeter Brune   if (Ydot) {
1176d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
11779566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Ydot", Ydot));
1178d5e6173cSPeter Brune     } else *Ydot = rw->Ydot;
1179d5e6173cSPeter Brune   }
1180d5e6173cSPeter Brune   if (Zdot) {
1181d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
11829566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Zdot", Zdot));
1183d5e6173cSPeter Brune     } else *Zdot = rw->Zdot;
1184d5e6173cSPeter Brune   }
1185d5e6173cSPeter Brune   if (Ystage) {
1186d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
11879566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Ystage", Ystage));
1188d5e6173cSPeter Brune     } else *Ystage = rw->Ystage;
1189d5e6173cSPeter Brune   }
1190d5e6173cSPeter Brune   if (Zstage) {
1191d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
11929566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Zstage", Zstage));
1193d5e6173cSPeter Brune     } else *Zstage = rw->Zstage;
1194d5e6173cSPeter Brune   }
1195d5e6173cSPeter Brune   PetscFunctionReturn(0);
1196d5e6173cSPeter Brune }
1197d5e6173cSPeter Brune 
1198d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWRestoreVecs(TS ts, DM dm, Vec *Ydot, Vec *Zdot, Vec *Ystage, Vec *Zstage)
1199d71ae5a4SJacob Faibussowitsch {
1200d5e6173cSPeter Brune   PetscFunctionBegin;
1201d5e6173cSPeter Brune   if (Ydot) {
120248a46eb9SPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Ydot", Ydot));
1203d5e6173cSPeter Brune   }
1204d5e6173cSPeter Brune   if (Zdot) {
120548a46eb9SPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Zdot", Zdot));
1206d5e6173cSPeter Brune   }
1207d5e6173cSPeter Brune   if (Ystage) {
120848a46eb9SPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Ystage", Ystage));
1209d5e6173cSPeter Brune   }
1210d5e6173cSPeter Brune   if (Zstage) {
121148a46eb9SPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Zstage", Zstage));
1212d5e6173cSPeter Brune   }
1213d5e6173cSPeter Brune   PetscFunctionReturn(0);
1214d5e6173cSPeter Brune }
1215d5e6173cSPeter Brune 
1216d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCoarsenHook_TSRosW(DM fine, DM coarse, void *ctx)
1217d71ae5a4SJacob Faibussowitsch {
1218d5e6173cSPeter Brune   PetscFunctionBegin;
1219d5e6173cSPeter Brune   PetscFunctionReturn(0);
1220d5e6173cSPeter Brune }
1221d5e6173cSPeter Brune 
1222d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRestrictHook_TSRosW(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx)
1223d71ae5a4SJacob Faibussowitsch {
1224d5e6173cSPeter Brune   TS  ts = (TS)ctx;
1225d5e6173cSPeter Brune   Vec Ydot, Zdot, Ystage, Zstage;
1226d5e6173cSPeter Brune   Vec Ydotc, Zdotc, Ystagec, Zstagec;
1227d5e6173cSPeter Brune 
1228d5e6173cSPeter Brune   PetscFunctionBegin;
12299566063dSJacob Faibussowitsch   PetscCall(TSRosWGetVecs(ts, fine, &Ydot, &Ystage, &Zdot, &Zstage));
12309566063dSJacob Faibussowitsch   PetscCall(TSRosWGetVecs(ts, coarse, &Ydotc, &Ystagec, &Zdotc, &Zstagec));
12319566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct, Ydot, Ydotc));
12329566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(Ydotc, rscale, Ydotc));
12339566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct, Ystage, Ystagec));
12349566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(Ystagec, rscale, Ystagec));
12359566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct, Zdot, Zdotc));
12369566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(Zdotc, rscale, Zdotc));
12379566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct, Zstage, Zstagec));
12389566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(Zstagec, rscale, Zstagec));
12399566063dSJacob Faibussowitsch   PetscCall(TSRosWRestoreVecs(ts, fine, &Ydot, &Ystage, &Zdot, &Zstage));
12409566063dSJacob Faibussowitsch   PetscCall(TSRosWRestoreVecs(ts, coarse, &Ydotc, &Ystagec, &Zdotc, &Zstagec));
1241d5e6173cSPeter Brune   PetscFunctionReturn(0);
1242d5e6173cSPeter Brune }
1243d5e6173cSPeter Brune 
1244d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainHook_TSRosW(DM fine, DM coarse, void *ctx)
1245d71ae5a4SJacob Faibussowitsch {
1246258e1594SPeter Brune   PetscFunctionBegin;
1247258e1594SPeter Brune   PetscFunctionReturn(0);
1248258e1594SPeter Brune }
1249258e1594SPeter Brune 
1250d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainRestrictHook_TSRosW(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx)
1251d71ae5a4SJacob Faibussowitsch {
1252258e1594SPeter Brune   TS  ts = (TS)ctx;
1253258e1594SPeter Brune   Vec Ydot, Zdot, Ystage, Zstage;
1254258e1594SPeter Brune   Vec Ydots, Zdots, Ystages, Zstages;
1255258e1594SPeter Brune 
1256258e1594SPeter Brune   PetscFunctionBegin;
12579566063dSJacob Faibussowitsch   PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Ystage, &Zdot, &Zstage));
12589566063dSJacob Faibussowitsch   PetscCall(TSRosWGetVecs(ts, subdm, &Ydots, &Ystages, &Zdots, &Zstages));
1259258e1594SPeter Brune 
12609566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(gscat, Ydot, Ydots, INSERT_VALUES, SCATTER_FORWARD));
12619566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(gscat, Ydot, Ydots, INSERT_VALUES, SCATTER_FORWARD));
1262258e1594SPeter Brune 
12639566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(gscat, Ystage, Ystages, INSERT_VALUES, SCATTER_FORWARD));
12649566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(gscat, Ystage, Ystages, INSERT_VALUES, SCATTER_FORWARD));
1265258e1594SPeter Brune 
12669566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(gscat, Zdot, Zdots, INSERT_VALUES, SCATTER_FORWARD));
12679566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(gscat, Zdot, Zdots, INSERT_VALUES, SCATTER_FORWARD));
1268258e1594SPeter Brune 
12699566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(gscat, Zstage, Zstages, INSERT_VALUES, SCATTER_FORWARD));
12709566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(gscat, Zstage, Zstages, INSERT_VALUES, SCATTER_FORWARD));
1271258e1594SPeter Brune 
12729566063dSJacob Faibussowitsch   PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Ystage, &Zdot, &Zstage));
12739566063dSJacob Faibussowitsch   PetscCall(TSRosWRestoreVecs(ts, subdm, &Ydots, &Ystages, &Zdots, &Zstages));
1274258e1594SPeter Brune   PetscFunctionReturn(0);
1275258e1594SPeter Brune }
1276258e1594SPeter Brune 
1277e27a552bSJed Brown /*
1278e27a552bSJed Brown   This defines the nonlinear equation that is to be solved with SNES
1279e27a552bSJed Brown   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
1280e27a552bSJed Brown */
1281d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_RosW(SNES snes, Vec U, Vec F, TS ts)
1282d71ae5a4SJacob Faibussowitsch {
128361692a83SJed Brown   TS_RosW  *ros = (TS_RosW *)ts->data;
1284d5e6173cSPeter Brune   Vec       Ydot, Zdot, Ystage, Zstage;
1285b296d7d5SJed Brown   PetscReal shift = ros->scoeff / ts->time_step;
1286d5e6173cSPeter Brune   DM        dm, dmsave;
1287e27a552bSJed Brown 
1288e27a552bSJed Brown   PetscFunctionBegin;
12899566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
12909566063dSJacob Faibussowitsch   PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage));
12919566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(Ydot, shift, U, Zdot));   /* Ydot = shift*U + Zdot */
12929566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(Ystage, 1.0, U, Zstage)); /* Ystage = U + Zstage */
1293d5e6173cSPeter Brune   dmsave = ts->dm;
1294d5e6173cSPeter Brune   ts->dm = dm;
12959566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts, ros->stage_time, Ystage, Ydot, F, PETSC_FALSE));
1296d5e6173cSPeter Brune   ts->dm = dmsave;
12979566063dSJacob Faibussowitsch   PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage));
1298e27a552bSJed Brown   PetscFunctionReturn(0);
1299e27a552bSJed Brown }
1300e27a552bSJed Brown 
1301d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes, Vec U, Mat A, Mat B, TS ts)
1302d71ae5a4SJacob Faibussowitsch {
130361692a83SJed Brown   TS_RosW  *ros = (TS_RosW *)ts->data;
1304d5e6173cSPeter Brune   Vec       Ydot, Zdot, Ystage, Zstage;
1305b296d7d5SJed Brown   PetscReal shift = ros->scoeff / ts->time_step;
1306d5e6173cSPeter Brune   DM        dm, dmsave;
1307e27a552bSJed Brown 
1308e27a552bSJed Brown   PetscFunctionBegin;
130961692a83SJed Brown   /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
13109566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
13119566063dSJacob Faibussowitsch   PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage));
1312d5e6173cSPeter Brune   dmsave = ts->dm;
1313d5e6173cSPeter Brune   ts->dm = dm;
13149566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts, ros->stage_time, Ystage, Ydot, shift, A, B, PETSC_TRUE));
1315d5e6173cSPeter Brune   ts->dm = dmsave;
13169566063dSJacob Faibussowitsch   PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage));
1317e27a552bSJed Brown   PetscFunctionReturn(0);
1318e27a552bSJed Brown }
1319e27a552bSJed Brown 
1320d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWTableauSetUp(TS ts)
1321d71ae5a4SJacob Faibussowitsch {
1322b39943a6SLisandro Dalcin   TS_RosW    *ros = (TS_RosW *)ts->data;
1323b39943a6SLisandro Dalcin   RosWTableau tab = ros->tableau;
1324b39943a6SLisandro Dalcin 
1325b39943a6SLisandro Dalcin   PetscFunctionBegin;
13269566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ros->Y));
13279566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tab->s, &ros->work));
1328b39943a6SLisandro Dalcin   PetscFunctionReturn(0);
1329b39943a6SLisandro Dalcin }
1330b39943a6SLisandro Dalcin 
1331d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_RosW(TS ts)
1332d71ae5a4SJacob Faibussowitsch {
133361692a83SJed Brown   TS_RosW      *ros = (TS_RosW *)ts->data;
1334d5e6173cSPeter Brune   DM            dm;
1335b39943a6SLisandro Dalcin   SNES          snes;
1336a3ab5968SHong Zhang   TSRHSJacobian rhsjacobian;
1337e27a552bSJed Brown 
1338e27a552bSJed Brown   PetscFunctionBegin;
13399566063dSJacob Faibussowitsch   PetscCall(TSRosWTableauSetUp(ts));
13409566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &ros->Ydot));
13419566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &ros->Ystage));
13429566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &ros->Zdot));
13439566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &ros->Zstage));
13449566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &ros->vec_sol_prev));
13459566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
13469566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSRosW, DMRestrictHook_TSRosW, ts));
13479566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSRosW, DMSubDomainRestrictHook_TSRosW, ts));
1348b39943a6SLisandro Dalcin   /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
13499566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
135048a46eb9SPierre Jolivet   if (!((PetscObject)snes)->type_name) PetscCall(SNESSetType(snes, SNESKSPONLY));
13519566063dSJacob Faibussowitsch   PetscCall(DMTSGetRHSJacobian(dm, &rhsjacobian, NULL));
1352a3ab5968SHong Zhang   if (rhsjacobian == TSComputeRHSJacobianConstant) {
1353a3ab5968SHong Zhang     Mat Amat, Pmat;
1354a3ab5968SHong Zhang 
1355a3ab5968SHong Zhang     /* Set the SNES matrix to be different from the RHS matrix because there is no way to reconstruct shift*M-J */
13569566063dSJacob Faibussowitsch     PetscCall(SNESGetJacobian(snes, &Amat, &Pmat, NULL, NULL));
1357a3ab5968SHong Zhang     if (Amat && Amat == ts->Arhs) {
1358a3ab5968SHong Zhang       if (Amat == Pmat) {
13599566063dSJacob Faibussowitsch         PetscCall(MatDuplicate(ts->Arhs, MAT_COPY_VALUES, &Amat));
13609566063dSJacob Faibussowitsch         PetscCall(SNESSetJacobian(snes, Amat, Amat, NULL, NULL));
1361a3ab5968SHong Zhang       } else {
13629566063dSJacob Faibussowitsch         PetscCall(MatDuplicate(ts->Arhs, MAT_COPY_VALUES, &Amat));
13639566063dSJacob Faibussowitsch         PetscCall(SNESSetJacobian(snes, Amat, NULL, NULL, NULL));
1364a3ab5968SHong Zhang         if (Pmat && Pmat == ts->Brhs) {
13659566063dSJacob Faibussowitsch           PetscCall(MatDuplicate(ts->Brhs, MAT_COPY_VALUES, &Pmat));
13669566063dSJacob Faibussowitsch           PetscCall(SNESSetJacobian(snes, NULL, Pmat, NULL, NULL));
13679566063dSJacob Faibussowitsch           PetscCall(MatDestroy(&Pmat));
1368a3ab5968SHong Zhang         }
1369a3ab5968SHong Zhang       }
13709566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Amat));
1371a3ab5968SHong Zhang     }
1372a3ab5968SHong Zhang   }
1373e27a552bSJed Brown   PetscFunctionReturn(0);
1374e27a552bSJed Brown }
1375e27a552bSJed Brown /*------------------------------------------------------------*/
1376e27a552bSJed Brown 
1377d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_RosW(TS ts, PetscOptionItems *PetscOptionsObject)
1378d71ae5a4SJacob Faibussowitsch {
137961692a83SJed Brown   TS_RosW *ros = (TS_RosW *)ts->data;
1380b39943a6SLisandro Dalcin   SNES     snes;
1381e27a552bSJed Brown 
1382e27a552bSJed Brown   PetscFunctionBegin;
1383d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "RosW ODE solver options");
1384e27a552bSJed Brown   {
138561692a83SJed Brown     RosWTableauLink link;
1386e27a552bSJed Brown     PetscInt        count, choice;
1387e27a552bSJed Brown     PetscBool       flg;
1388e27a552bSJed Brown     const char    **namelist;
138961692a83SJed Brown 
13909371c9d4SSatish Balay     for (link = RosWTableauList, count = 0; link; link = link->next, count++)
13919371c9d4SSatish Balay       ;
13929566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(count, (char ***)&namelist));
139361692a83SJed Brown     for (link = RosWTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name;
13949566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEList("-ts_rosw_type", "Family of Rosenbrock-W method", "TSRosWSetType", (const char *const *)namelist, count, ros->tableau->name, &choice, &flg));
13959566063dSJacob Faibussowitsch     if (flg) PetscCall(TSRosWSetType(ts, namelist[choice]));
13969566063dSJacob Faibussowitsch     PetscCall(PetscFree(namelist));
139761692a83SJed Brown 
13989566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ts_rosw_recompute_jacobian", "Recompute the Jacobian at each stage", "TSRosWSetRecomputeJacobian", ros->recompute_jacobian, &ros->recompute_jacobian, NULL));
1399b39943a6SLisandro Dalcin   }
1400d0609cedSBarry Smith   PetscOptionsHeadEnd();
140161692a83SJed Brown   /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
14029566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
140348a46eb9SPierre Jolivet   if (!((PetscObject)snes)->type_name) PetscCall(SNESSetType(snes, SNESKSPONLY));
1404e27a552bSJed Brown   PetscFunctionReturn(0);
1405e27a552bSJed Brown }
1406e27a552bSJed Brown 
1407d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_RosW(TS ts, PetscViewer viewer)
1408d71ae5a4SJacob Faibussowitsch {
140961692a83SJed Brown   TS_RosW  *ros = (TS_RosW *)ts->data;
1410e27a552bSJed Brown   PetscBool iascii;
1411e27a552bSJed Brown 
1412e27a552bSJed Brown   PetscFunctionBegin;
14139566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1414e27a552bSJed Brown   if (iascii) {
14159c334d8fSLisandro Dalcin     RosWTableau tab = ros->tableau;
141619fd82e9SBarry Smith     TSRosWType  rostype;
14179c334d8fSLisandro Dalcin     char        buf[512];
1418e408995aSJed Brown     PetscInt    i;
1419e408995aSJed Brown     PetscReal   abscissa[512];
14209566063dSJacob Faibussowitsch     PetscCall(TSRosWGetType(ts, &rostype));
14219566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Rosenbrock-W %s\n", rostype));
14229566063dSJacob Faibussowitsch     PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, tab->ASum));
14239566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Abscissa of A       = %s\n", buf));
1424e408995aSJed Brown     for (i = 0; i < tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i];
14259566063dSJacob Faibussowitsch     PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, abscissa));
14269566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Abscissa of A+Gamma = %s\n", buf));
1427e27a552bSJed Brown   }
1428e27a552bSJed Brown   PetscFunctionReturn(0);
1429e27a552bSJed Brown }
1430e27a552bSJed Brown 
1431d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSLoad_RosW(TS ts, PetscViewer viewer)
1432d71ae5a4SJacob Faibussowitsch {
14339200755eSBarry Smith   SNES    snes;
14349c334d8fSLisandro Dalcin   TSAdapt adapt;
14359200755eSBarry Smith 
14369200755eSBarry Smith   PetscFunctionBegin;
14379566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &adapt));
14389566063dSJacob Faibussowitsch   PetscCall(TSAdaptLoad(adapt, viewer));
14399566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
14409566063dSJacob Faibussowitsch   PetscCall(SNESLoad(snes, viewer));
14419200755eSBarry Smith   /* function and Jacobian context for SNES when used with TS is always ts object */
14429566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, NULL, NULL, ts));
14439566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, ts));
14449200755eSBarry Smith   PetscFunctionReturn(0);
14459200755eSBarry Smith }
14469200755eSBarry Smith 
1447e27a552bSJed Brown /*@C
1448*bcf0153eSBarry Smith   TSRosWSetType - Set the type of Rosenbrock-W, `TSROSW`, scheme
1449e27a552bSJed Brown 
1450e27a552bSJed Brown   Logically collective
1451e27a552bSJed Brown 
1452d8d19677SJose E. Roman   Input Parameters:
1453e27a552bSJed Brown +  ts - timestepping context
1454b92453a8SLisandro Dalcin -  roswtype - type of Rosenbrock-W scheme
1455e27a552bSJed Brown 
1456020d8f30SJed Brown   Level: beginner
1457e27a552bSJed Brown 
1458*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSRosWGetType()`, `TSROSW`, `TSROSW2M`, `TSROSW2P`, `TSROSWRA3PW`, `TSROSWRA34PW2`, `TSROSWRODAS3`, `TSROSWSANDU3`, `TSROSWASSP3P3S1C`, `TSROSWLASSP3P4S2C`, `TSROSWLLSSP3P4S2C`, `TSROSWARK3`
1459e27a552bSJed Brown @*/
1460d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWSetType(TS ts, TSRosWType roswtype)
1461d71ae5a4SJacob Faibussowitsch {
1462e27a552bSJed Brown   PetscFunctionBegin;
1463e27a552bSJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1464b92453a8SLisandro Dalcin   PetscValidCharPointer(roswtype, 2);
1465cac4c232SBarry Smith   PetscTryMethod(ts, "TSRosWSetType_C", (TS, TSRosWType), (ts, roswtype));
1466e27a552bSJed Brown   PetscFunctionReturn(0);
1467e27a552bSJed Brown }
1468e27a552bSJed Brown 
1469e27a552bSJed Brown /*@C
147061692a83SJed Brown   TSRosWGetType - Get the type of Rosenbrock-W scheme
1471e27a552bSJed Brown 
1472e27a552bSJed Brown   Logically collective
1473e27a552bSJed Brown 
1474e27a552bSJed Brown   Input Parameter:
1475e27a552bSJed Brown .  ts - timestepping context
1476e27a552bSJed Brown 
1477e27a552bSJed Brown   Output Parameter:
147861692a83SJed Brown .  rostype - type of Rosenbrock-W scheme
1479e27a552bSJed Brown 
1480e27a552bSJed Brown   Level: intermediate
1481e27a552bSJed Brown 
1482*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSRosWType`, `TSRosWSetType()`
1483e27a552bSJed Brown @*/
1484d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWGetType(TS ts, TSRosWType *rostype)
1485d71ae5a4SJacob Faibussowitsch {
1486e27a552bSJed Brown   PetscFunctionBegin;
1487e27a552bSJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1488cac4c232SBarry Smith   PetscUseMethod(ts, "TSRosWGetType_C", (TS, TSRosWType *), (ts, rostype));
1489e27a552bSJed Brown   PetscFunctionReturn(0);
1490e27a552bSJed Brown }
1491e27a552bSJed Brown 
1492e27a552bSJed Brown /*@C
149361692a83SJed Brown   TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step.
1494e27a552bSJed Brown 
1495e27a552bSJed Brown   Logically collective
1496e27a552bSJed Brown 
1497d8d19677SJose E. Roman   Input Parameters:
1498e27a552bSJed Brown +  ts - timestepping context
1499*bcf0153eSBarry Smith -  flg - `PETSC_TRUE` to recompute the Jacobian at each stage
1500e27a552bSJed Brown 
1501e27a552bSJed Brown   Level: intermediate
1502e27a552bSJed Brown 
1503*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSRosWType`, `TSRosWGetType()`
1504e27a552bSJed Brown @*/
1505d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRosWSetRecomputeJacobian(TS ts, PetscBool flg)
1506d71ae5a4SJacob Faibussowitsch {
1507e27a552bSJed Brown   PetscFunctionBegin;
1508e27a552bSJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1509cac4c232SBarry Smith   PetscTryMethod(ts, "TSRosWSetRecomputeJacobian_C", (TS, PetscBool), (ts, flg));
1510e27a552bSJed Brown   PetscFunctionReturn(0);
1511e27a552bSJed Brown }
1512e27a552bSJed Brown 
1513d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWGetType_RosW(TS ts, TSRosWType *rostype)
1514d71ae5a4SJacob Faibussowitsch {
151561692a83SJed Brown   TS_RosW *ros = (TS_RosW *)ts->data;
1516e27a552bSJed Brown 
1517e27a552bSJed Brown   PetscFunctionBegin;
151861692a83SJed Brown   *rostype = ros->tableau->name;
1519e27a552bSJed Brown   PetscFunctionReturn(0);
1520e27a552bSJed Brown }
1521ef20d060SBarry Smith 
1522d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWSetType_RosW(TS ts, TSRosWType rostype)
1523d71ae5a4SJacob Faibussowitsch {
152461692a83SJed Brown   TS_RosW        *ros = (TS_RosW *)ts->data;
1525e27a552bSJed Brown   PetscBool       match;
152661692a83SJed Brown   RosWTableauLink link;
1527e27a552bSJed Brown 
1528e27a552bSJed Brown   PetscFunctionBegin;
152961692a83SJed Brown   if (ros->tableau) {
15309566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(ros->tableau->name, rostype, &match));
1531e27a552bSJed Brown     if (match) PetscFunctionReturn(0);
1532e27a552bSJed Brown   }
153361692a83SJed Brown   for (link = RosWTableauList; link; link = link->next) {
15349566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(link->tab.name, rostype, &match));
1535e27a552bSJed Brown     if (match) {
15369566063dSJacob Faibussowitsch       if (ts->setupcalled) PetscCall(TSRosWTableauReset(ts));
153761692a83SJed Brown       ros->tableau = &link->tab;
15389566063dSJacob Faibussowitsch       if (ts->setupcalled) PetscCall(TSRosWTableauSetUp(ts));
15392ffb9264SLisandro Dalcin       ts->default_adapt_type = ros->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
1540e27a552bSJed Brown       PetscFunctionReturn(0);
1541e27a552bSJed Brown     }
1542e27a552bSJed Brown   }
154398921bdaSJacob Faibussowitsch   SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", rostype);
1544e27a552bSJed Brown }
154561692a83SJed Brown 
1546d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRosWSetRecomputeJacobian_RosW(TS ts, PetscBool flg)
1547d71ae5a4SJacob Faibussowitsch {
154861692a83SJed Brown   TS_RosW *ros = (TS_RosW *)ts->data;
1549e27a552bSJed Brown 
1550e27a552bSJed Brown   PetscFunctionBegin;
155161692a83SJed Brown   ros->recompute_jacobian = flg;
1552e27a552bSJed Brown   PetscFunctionReturn(0);
1553e27a552bSJed Brown }
1554e27a552bSJed Brown 
1555d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_RosW(TS ts)
1556d71ae5a4SJacob Faibussowitsch {
1557b3a6b972SJed Brown   PetscFunctionBegin;
15589566063dSJacob Faibussowitsch   PetscCall(TSReset_RosW(ts));
1559b3a6b972SJed Brown   if (ts->dm) {
15609566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSRosW, DMRestrictHook_TSRosW, ts));
15619566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSRosW, DMSubDomainRestrictHook_TSRosW, ts));
1562b3a6b972SJed Brown   }
15639566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
15649566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWGetType_C", NULL));
15659566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetType_C", NULL));
15669566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetRecomputeJacobian_C", NULL));
1567b3a6b972SJed Brown   PetscFunctionReturn(0);
1568b3a6b972SJed Brown }
1569d5e6173cSPeter Brune 
1570e27a552bSJed Brown /* ------------------------------------------------------------ */
1571e27a552bSJed Brown /*MC
1572020d8f30SJed Brown       TSROSW - ODE solver using Rosenbrock-W schemes
1573e27a552bSJed Brown 
1574e27a552bSJed Brown   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
1575e27a552bSJed Brown   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
1576*bcf0153eSBarry Smith   of the equation using `TSSetIFunction()` and the non-stiff part with `TSSetRHSFunction()`.
1577*bcf0153eSBarry Smith 
1578*bcf0153eSBarry Smith   Level: beginner
1579e27a552bSJed Brown 
1580e27a552bSJed Brown   Notes:
158161692a83SJed Brown   This method currently only works with autonomous ODE and DAE.
158261692a83SJed Brown 
1583*bcf0153eSBarry Smith   Consider trying `TSARKIMEX` if the stiff part is strongly nonlinear.
1584d0685a90SJed Brown 
15853d5a8a6aSBarry 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
15863d5a8a6aSBarry Smith 
1587e94cfbe0SPatrick Sanan   Developer Notes:
158861692a83SJed Brown   Rosenbrock-W methods are typically specified for autonomous ODE
158961692a83SJed Brown 
1590f9c1d6abSBarry Smith $  udot = f(u)
159161692a83SJed Brown 
159261692a83SJed Brown   by the stage equations
159361692a83SJed Brown 
1594f9c1d6abSBarry Smith $  k_i = h f(u_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j
159561692a83SJed Brown 
159661692a83SJed Brown   and step completion formula
159761692a83SJed Brown 
1598f9c1d6abSBarry Smith $  u_1 = u_0 + sum_j b_j k_j
159961692a83SJed Brown 
1600f9c1d6abSBarry Smith   with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(u)
160161692a83SJed Brown   and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner,
160261692a83SJed Brown   we define new variables for the stage equations
160361692a83SJed Brown 
160461692a83SJed Brown $  y_i = gamma_ij k_j
160561692a83SJed Brown 
160661692a83SJed Brown   The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define
160761692a83SJed Brown 
1608b70472e2SJed Brown $  A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-1}
160961692a83SJed Brown 
161061692a83SJed Brown   to rewrite the method as
161161692a83SJed Brown 
1612f9c1d6abSBarry 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
1613f9c1d6abSBarry Smith $  u_1 = u_0 + sum_j bt_j y_j
161461692a83SJed Brown 
161561692a83SJed Brown    where we have introduced the mass matrix M. Continue by defining
161661692a83SJed Brown 
161761692a83SJed Brown $  ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j
161861692a83SJed Brown 
161961692a83SJed Brown    or, more compactly in tensor notation
162061692a83SJed Brown 
162161692a83SJed Brown $  Ydot = 1/h (Gamma^{-1} \otimes I) Y .
162261692a83SJed Brown 
162361692a83SJed Brown    Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current
162461692a83SJed Brown    stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the
162561692a83SJed Brown    equation
162661692a83SJed Brown 
1627f9c1d6abSBarry Smith $  g(u_0 + sum_j a_ij y_j + y_i, ydot_i) = 0
162861692a83SJed Brown 
162961692a83SJed Brown    with initial guess y_i = 0.
1630e27a552bSJed Brown 
1631*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSRosWSetType()`, `TSRosWRegister()`, `TSROSWTHETA1`, `TSROSWTHETA2`, `TSROSW2M`, `TSROSW2P`, `TSROSWRA3PW`, `TSROSWRA34PW2`, `TSROSWRODAS3`,
1632*bcf0153eSBarry Smith           `TSROSWSANDU3`, `TSROSWASSP3P3S1C`, `TSROSWLASSP3P4S2C`, `TSROSWLLSSP3P4S2C`, `TSROSWGRK4T`, `TSROSWSHAMP4`, `TSROSWVELDD4`, `TSROSW4L`, `TSType`
1633e27a552bSJed Brown M*/
1634d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_RosW(TS ts)
1635d71ae5a4SJacob Faibussowitsch {
163661692a83SJed Brown   TS_RosW *ros;
1637e27a552bSJed Brown 
1638e27a552bSJed Brown   PetscFunctionBegin;
16399566063dSJacob Faibussowitsch   PetscCall(TSRosWInitializePackage());
1640e27a552bSJed Brown 
1641e27a552bSJed Brown   ts->ops->reset          = TSReset_RosW;
1642e27a552bSJed Brown   ts->ops->destroy        = TSDestroy_RosW;
1643e27a552bSJed Brown   ts->ops->view           = TSView_RosW;
16449200755eSBarry Smith   ts->ops->load           = TSLoad_RosW;
1645e27a552bSJed Brown   ts->ops->setup          = TSSetUp_RosW;
1646e27a552bSJed Brown   ts->ops->step           = TSStep_RosW;
1647e27a552bSJed Brown   ts->ops->interpolate    = TSInterpolate_RosW;
16481c3436cfSJed Brown   ts->ops->evaluatestep   = TSEvaluateStep_RosW;
164924655328SShri   ts->ops->rollback       = TSRollBack_RosW;
1650e27a552bSJed Brown   ts->ops->setfromoptions = TSSetFromOptions_RosW;
1651e27a552bSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_RosW;
1652e27a552bSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_RosW;
1653e27a552bSJed Brown 
1654efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
1655efd4aadfSBarry Smith 
16564dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&ros));
165761692a83SJed Brown   ts->data = (void *)ros;
1658e27a552bSJed Brown 
16599566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWGetType_C", TSRosWGetType_RosW));
16609566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetType_C", TSRosWSetType_RosW));
16619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetRecomputeJacobian_C", TSRosWSetRecomputeJacobian_RosW));
1662b39943a6SLisandro Dalcin 
16639566063dSJacob Faibussowitsch   PetscCall(TSRosWSetType(ts, TSRosWDefault));
1664e27a552bSJed Brown   PetscFunctionReturn(0);
1665e27a552bSJed Brown }
1666