xref: /petsc/src/ts/impls/rosw/rosw.c (revision f4aed992dd926c4e715d02cabab6ea04d956d070)
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 
761692a83SJed Brown   G(t,X,Xdot) = F(t,X)
8e27a552bSJed Brown 
961692a83SJed Brown   where G represents the stiff part of the physics and F represents the non-stiff part.
1061692a83SJed Brown   This method is designed to be linearly implicit on G and can use an approximate and lagged Jacobian.
11e27a552bSJed Brown 
12e27a552bSJed Brown */
13e27a552bSJed Brown #include <private/tsimpl.h>                /*I   "petscts.h"   I*/
14e27a552bSJed Brown 
1561692a83SJed Brown #include <../src/mat/blockinvert.h>
1661692a83SJed Brown 
17ae6f9490SJed Brown static const TSRosWType TSRosWDefault = TSROSWRA34PW2;
18e27a552bSJed Brown static PetscBool TSRosWRegisterAllCalled;
19e27a552bSJed Brown static PetscBool TSRosWPackageInitialized;
20e27a552bSJed Brown 
2161692a83SJed Brown typedef struct _RosWTableau *RosWTableau;
2261692a83SJed Brown struct _RosWTableau {
23e27a552bSJed Brown   char      *name;
24e27a552bSJed Brown   PetscInt  order;              /* Classical approximation order of the method */
25e27a552bSJed Brown   PetscInt  s;                  /* Number of stages */
26*f4aed992SEmil Constantinescu   PetscInt  pinterp;            /* Interpolation order */
2761692a83SJed Brown   PetscReal *A;                 /* Propagation table, strictly lower triangular */
2861692a83SJed Brown   PetscReal *Gamma;             /* Stage table, lower triangular with nonzero diagonal */
29c17803e7SJed Brown   PetscBool *GammaZeroDiag;     /* Diagonal entries that are zero in stage table Gamma, vector indicating explicit statages */
3043b21953SEmil Constantinescu   PetscReal *GammaExplicitCorr; /* Coefficients for correction terms needed for explicit stages in transformed variables*/
3161692a83SJed Brown   PetscReal *b;                 /* Step completion table */
32fe7e6d57SJed Brown   PetscReal *bembed;            /* Step completion table for embedded method of order one less */
3361692a83SJed Brown   PetscReal *ASum;              /* Row sum of A */
3461692a83SJed Brown   PetscReal *GammaSum;          /* Row sum of Gamma, only needed for non-autonomous systems */
3561692a83SJed Brown   PetscReal *At;                /* Propagation table in transformed variables */
3661692a83SJed Brown   PetscReal *bt;                /* Step completion table in transformed variables */
37fe7e6d57SJed Brown   PetscReal *bembedt;           /* Step completion table of order one less in transformed variables */
3861692a83SJed Brown   PetscReal *GammaInv;          /* Inverse of Gamma, used for transformed variables */
398d59e960SJed Brown   PetscReal ccfl;               /* Placeholder for CFL coefficient relative to forward Euler */
40*f4aed992SEmil Constantinescu   PetscReal *binterpt,*binterp; /* Dense output formula */
41e27a552bSJed Brown };
4261692a83SJed Brown typedef struct _RosWTableauLink *RosWTableauLink;
4361692a83SJed Brown struct _RosWTableauLink {
4461692a83SJed Brown   struct _RosWTableau tab;
4561692a83SJed Brown   RosWTableauLink next;
46e27a552bSJed Brown };
4761692a83SJed Brown static RosWTableauLink RosWTableauList;
48e27a552bSJed Brown 
49e27a552bSJed Brown typedef struct {
5061692a83SJed Brown   RosWTableau tableau;
5161692a83SJed Brown   Vec         *Y;               /* States computed during the step, used to complete the step */
52e27a552bSJed Brown   Vec         Ydot;             /* Work vector holding Ydot during residual evaluation */
5361692a83SJed Brown   Vec         Ystage;           /* Work vector for the state value at each stage */
5461692a83SJed Brown   Vec         Zdot;             /* Ydot = Zdot + shift*Y */
5561692a83SJed Brown   Vec         Zstage;           /* Y = Zstage + Y */
561c3436cfSJed Brown   PetscScalar *work;            /* Scalar work space of length number of stages, used to prepare VecMAXPY() */
57e27a552bSJed Brown   PetscReal   shift;
58e27a552bSJed Brown   PetscReal   stage_time;
59c17803e7SJed Brown   PetscReal   stage_explicit;     /* Flag indicates that the current stage is explicit */
6061692a83SJed Brown   PetscBool   recompute_jacobian; /* Recompute the Jacobian at each stage, default is to freeze the Jacobian at the start of each step */
61108c343cSJed Brown   TSStepStatus status;
62e27a552bSJed Brown } TS_RosW;
63e27a552bSJed Brown 
64fe7e6d57SJed Brown /*MC
653606a31eSEmil Constantinescu      TSROSWTHETA1 - One stage first order L-stable Rosenbrock-W scheme (aka theta method).
663606a31eSEmil Constantinescu 
673606a31eSEmil Constantinescu      Only an approximate Jacobian is needed.
683606a31eSEmil Constantinescu 
693606a31eSEmil Constantinescu      Level: intermediate
703606a31eSEmil Constantinescu 
713606a31eSEmil Constantinescu .seealso: TSROSW
723606a31eSEmil Constantinescu M*/
733606a31eSEmil Constantinescu 
743606a31eSEmil Constantinescu /*MC
753606a31eSEmil Constantinescu      TSROSWTHETA2 - One stage second order A-stable Rosenbrock-W scheme (aka theta method).
763606a31eSEmil Constantinescu 
773606a31eSEmil Constantinescu      Only an approximate Jacobian is needed.
783606a31eSEmil Constantinescu 
793606a31eSEmil Constantinescu      Level: intermediate
803606a31eSEmil Constantinescu 
813606a31eSEmil Constantinescu .seealso: TSROSW
823606a31eSEmil Constantinescu M*/
833606a31eSEmil Constantinescu 
843606a31eSEmil Constantinescu /*MC
85fe7e6d57SJed Brown      TSROSW2M - Two stage second order L-stable Rosenbrock-W scheme.
86fe7e6d57SJed Brown 
87fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2P.
88fe7e6d57SJed Brown 
89fe7e6d57SJed Brown      Level: intermediate
90fe7e6d57SJed Brown 
91fe7e6d57SJed Brown .seealso: TSROSW
92fe7e6d57SJed Brown M*/
93fe7e6d57SJed Brown 
94fe7e6d57SJed Brown /*MC
95fe7e6d57SJed Brown      TSROSW2P - Two stage second order L-stable Rosenbrock-W scheme.
96fe7e6d57SJed Brown 
97fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2M.
98fe7e6d57SJed Brown 
99fe7e6d57SJed Brown      Level: intermediate
100fe7e6d57SJed Brown 
101fe7e6d57SJed Brown .seealso: TSROSW
102fe7e6d57SJed Brown M*/
103fe7e6d57SJed Brown 
104fe7e6d57SJed Brown /*MC
105fe7e6d57SJed Brown      TSROSWRA3PW - Three stage third order Rosenbrock-W scheme for PDAE of index 1.
106fe7e6d57SJed Brown 
107fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.
108fe7e6d57SJed Brown 
109fe7e6d57SJed 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.
110fe7e6d57SJed Brown 
111fe7e6d57SJed Brown      References:
112fe7e6d57SJed Brown      Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005.
113fe7e6d57SJed Brown 
114fe7e6d57SJed Brown      Level: intermediate
115fe7e6d57SJed Brown 
116fe7e6d57SJed Brown .seealso: TSROSW
117fe7e6d57SJed Brown M*/
118fe7e6d57SJed Brown 
119fe7e6d57SJed Brown /*MC
120fe7e6d57SJed Brown      TSROSWRA34PW2 - Four stage third order L-stable Rosenbrock-W scheme for PDAE of index 1.
121fe7e6d57SJed Brown 
122fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.
123fe7e6d57SJed Brown 
124fe7e6d57SJed 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.
125fe7e6d57SJed Brown 
126fe7e6d57SJed Brown      References:
127fe7e6d57SJed Brown      Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005.
128fe7e6d57SJed Brown 
129fe7e6d57SJed Brown      Level: intermediate
130fe7e6d57SJed Brown 
131fe7e6d57SJed Brown .seealso: TSROSW
132fe7e6d57SJed Brown M*/
133fe7e6d57SJed Brown 
134ef3c5b88SJed Brown /*MC
135ef3c5b88SJed Brown      TSROSWRODAS3 - Four stage third order L-stable Rosenbrock scheme
136ef3c5b88SJed Brown 
137ef3c5b88SJed Brown      By default, the Jacobian is only recomputed once per step.
138ef3c5b88SJed Brown 
139ef3c5b88SJed Brown      Both the third order and embedded second order methods are stiffly accurate and L-stable.
140ef3c5b88SJed Brown 
141ef3c5b88SJed Brown      References:
142ef3c5b88SJed Brown      Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997.
143ef3c5b88SJed Brown 
144ef3c5b88SJed Brown      Level: intermediate
145ef3c5b88SJed Brown 
146ef3c5b88SJed Brown .seealso: TSROSW, TSROSWSANDU3
147ef3c5b88SJed Brown M*/
148ef3c5b88SJed Brown 
149ef3c5b88SJed Brown /*MC
150ef3c5b88SJed Brown      TSROSWSANDU3 - Three stage third order L-stable Rosenbrock scheme
151ef3c5b88SJed Brown 
152ef3c5b88SJed Brown      By default, the Jacobian is only recomputed once per step.
153ef3c5b88SJed Brown 
154ef3c5b88SJed Brown      The third order method is L-stable, but not stiffly accurate.
155ef3c5b88SJed Brown      The second order embedded method is strongly A-stable with R(infty) = 0.5.
156ef3c5b88SJed Brown      The internal stages are L-stable.
157ef3c5b88SJed Brown      This method is called ROS3 in the paper.
158ef3c5b88SJed Brown 
159ef3c5b88SJed Brown      References:
160ef3c5b88SJed Brown      Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997.
161ef3c5b88SJed Brown 
162ef3c5b88SJed Brown      Level: intermediate
163ef3c5b88SJed Brown 
164ef3c5b88SJed Brown .seealso: TSROSW, TSROSWRODAS3
165ef3c5b88SJed Brown M*/
166ef3c5b88SJed Brown 
167961f28d0SJed Brown /*MC
168961f28d0SJed Brown      TSROSWASSP3P3S1C - A-stable Rosenbrock-W method with SSP explicit part, third order, three stages
169961f28d0SJed Brown 
170961f28d0SJed Brown      By default, the Jacobian is only recomputed once per step.
171961f28d0SJed Brown 
172961f28d0SJed Brown      A-stable SPP explicit order 3, 3 stages, CFL 1 (eff = 1/3)
173961f28d0SJed Brown 
174961f28d0SJed Brown      References:
175961f28d0SJed Brown      Emil Constantinescu
176961f28d0SJed Brown 
177961f28d0SJed Brown      Level: intermediate
178961f28d0SJed Brown 
17943b21953SEmil Constantinescu .seealso: TSROSW, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, SSP
180961f28d0SJed Brown M*/
181961f28d0SJed Brown 
182961f28d0SJed Brown /*MC
183961f28d0SJed Brown      TSROSWLASSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, three stages
184961f28d0SJed Brown 
185961f28d0SJed Brown      By default, the Jacobian is only recomputed once per step.
186961f28d0SJed Brown 
187961f28d0SJed Brown      L-stable (A-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2)
188961f28d0SJed Brown 
189961f28d0SJed Brown      References:
190961f28d0SJed Brown      Emil Constantinescu
191961f28d0SJed Brown 
192961f28d0SJed Brown      Level: intermediate
193961f28d0SJed Brown 
19443b21953SEmil Constantinescu .seealso: TSROSW, TSROSWASSP3P3S1C, TSROSWLLSSP3P4S2C, TSSSP
195961f28d0SJed Brown M*/
196961f28d0SJed Brown 
197961f28d0SJed Brown /*MC
19843b21953SEmil Constantinescu      TSROSWLLSSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, three stages
199961f28d0SJed Brown 
200961f28d0SJed Brown      By default, the Jacobian is only recomputed once per step.
201961f28d0SJed Brown 
202961f28d0SJed Brown      L-stable (L-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2)
203961f28d0SJed Brown 
204961f28d0SJed Brown      References:
205961f28d0SJed Brown      Emil Constantinescu
206961f28d0SJed Brown 
207961f28d0SJed Brown      Level: intermediate
208961f28d0SJed Brown 
209961f28d0SJed Brown .seealso: TSROSW, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSSSP
210961f28d0SJed Brown M*/
211961f28d0SJed Brown 
212e27a552bSJed Brown #undef __FUNCT__
213e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterAll"
214e27a552bSJed Brown /*@C
215e27a552bSJed Brown   TSRosWRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSRosW
216e27a552bSJed Brown 
217e27a552bSJed Brown   Not Collective, but should be called by all processes which will need the schemes to be registered
218e27a552bSJed Brown 
219e27a552bSJed Brown   Level: advanced
220e27a552bSJed Brown 
221e27a552bSJed Brown .keywords: TS, TSRosW, register, all
222e27a552bSJed Brown 
223e27a552bSJed Brown .seealso:  TSRosWRegisterDestroy()
224e27a552bSJed Brown @*/
225e27a552bSJed Brown PetscErrorCode TSRosWRegisterAll(void)
226e27a552bSJed Brown {
227e27a552bSJed Brown   PetscErrorCode ierr;
228e27a552bSJed Brown 
229e27a552bSJed Brown   PetscFunctionBegin;
230e27a552bSJed Brown   if (TSRosWRegisterAllCalled) PetscFunctionReturn(0);
231e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_TRUE;
232e27a552bSJed Brown 
233e27a552bSJed Brown   {
2343606a31eSEmil Constantinescu     const PetscReal
2353606a31eSEmil Constantinescu       A = 0,
2363606a31eSEmil Constantinescu       Gamma = 1,
2373606a31eSEmil Constantinescu       b = 1;
238*f4aed992SEmil Constantinescu     ierr = TSRosWRegister(TSROSWTHETA1,1,1,&A,&Gamma,&b,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
2393606a31eSEmil Constantinescu   }
2403606a31eSEmil Constantinescu 
2413606a31eSEmil Constantinescu   {
2423606a31eSEmil Constantinescu     const PetscReal
2433606a31eSEmil Constantinescu       A= 0,
2443606a31eSEmil Constantinescu       Gamma = 0.5,
2453606a31eSEmil Constantinescu       b = 1;
246*f4aed992SEmil Constantinescu     ierr = TSRosWRegister(TSROSWTHETA2,2,1,&A,&Gamma,&b,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
2473606a31eSEmil Constantinescu   }
2483606a31eSEmil Constantinescu 
2493606a31eSEmil Constantinescu   {
25061692a83SJed Brown     const PetscReal g = 1. + 1./PetscSqrtReal(2.0);
251e27a552bSJed Brown     const PetscReal
25261692a83SJed Brown       A[2][2] = {{0,0}, {1.,0}},
25361692a83SJed Brown       Gamma[2][2] = {{g,0}, {-2.*g,g}},
2541c3436cfSJed Brown       b[2] = {0.5,0.5},
2551c3436cfSJed Brown       b1[2] = {1.0,0.0};
256*f4aed992SEmil Constantinescu     ierr = TSRosWRegister(TSROSW2P,2,2,&A[0][0],&Gamma[0][0],b,b1,0,PETSC_NULL);CHKERRQ(ierr);
257e27a552bSJed Brown   }
258e27a552bSJed Brown   {
25961692a83SJed Brown     const PetscReal g = 1. - 1./PetscSqrtReal(2.0);
260e27a552bSJed Brown     const PetscReal
26161692a83SJed Brown       A[2][2] = {{0,0}, {1.,0}},
26261692a83SJed Brown       Gamma[2][2] = {{g,0}, {-2.*g,g}},
2631c3436cfSJed Brown       b[2] = {0.5,0.5},
2641c3436cfSJed Brown       b1[2] = {1.0,0.0};
265*f4aed992SEmil Constantinescu     ierr = TSRosWRegister(TSROSW2M,2,2,&A[0][0],&Gamma[0][0],b,b1,0,PETSC_NULL);CHKERRQ(ierr);
266fe7e6d57SJed Brown   }
267fe7e6d57SJed Brown   {
268fe7e6d57SJed Brown     const PetscReal g = 7.8867513459481287e-01;
269fe7e6d57SJed Brown     const PetscReal
270fe7e6d57SJed Brown       A[3][3] = {{0,0,0},
271fe7e6d57SJed Brown                  {1.5773502691896257e+00,0,0},
272fe7e6d57SJed Brown                  {0.5,0,0}},
273fe7e6d57SJed Brown       Gamma[3][3] = {{g,0,0},
274fe7e6d57SJed Brown                      {-1.5773502691896257e+00,g,0},
27525833a80SEmil Constantinescu                      {-6.7075317547305480e-01,-1.7075317547305482e-01,g}},
276fe7e6d57SJed Brown       b[3] = {1.0566243270259355e-01,4.9038105676657971e-02,8.4529946162074843e-01},
277fe7e6d57SJed Brown       b2[3] = {-1.7863279495408180e-01,1./3.,8.4529946162074843e-01};
278*f4aed992SEmil Constantinescu     ierr = TSRosWRegister(TSROSWRA3PW,3,3,&A[0][0],&Gamma[0][0],b,b2,0,PETSC_NULL);CHKERRQ(ierr);
279fe7e6d57SJed Brown   }
280fe7e6d57SJed Brown   {
281fe7e6d57SJed Brown     const PetscReal g = 4.3586652150845900e-01;
282fe7e6d57SJed Brown     const PetscReal
283fe7e6d57SJed Brown       A[4][4] = {{0,0,0,0},
284fe7e6d57SJed Brown                  {8.7173304301691801e-01,0,0,0},
285fe7e6d57SJed Brown                  {8.4457060015369423e-01,-1.1299064236484185e-01,0,0},
286fe7e6d57SJed Brown                  {0,0,1.,0}},
287fe7e6d57SJed Brown       Gamma[4][4] = {{g,0,0,0},
288fe7e6d57SJed Brown                      {-8.7173304301691801e-01,g,0,0},
289fe7e6d57SJed Brown                      {-9.0338057013044082e-01,5.4180672388095326e-02,g,0},
290fe7e6d57SJed Brown                      {2.4212380706095346e-01,-1.2232505839045147e+00,5.4526025533510214e-01,g}},
291fe7e6d57SJed Brown         b[4] = {2.4212380706095346e-01,-1.2232505839045147e+00,1.5452602553351020e+00,4.3586652150845900e-01},
292*f4aed992SEmil Constantinescu           b2[4] = {3.7810903145819369e-01,-9.6042292212423178e-02,5.0000000000000000e-01,2.1793326075422950e-01},
293*f4aed992SEmil Constantinescu             binterpt[4][3]={{1.0564298455794094,-1.3864882699759573,0.5721822314575016},{2.296429974281067,-8.262611700275677,4.742931142090097},{-1.307599564525376,7.250979895056055,-4.398120075195578},{-1.045260255335102,2.398120075195581,-0.9169932983520199}};
294*f4aed992SEmil Constantinescu       ierr = TSRosWRegister(TSROSWRA34PW2,3,4,&A[0][0],&Gamma[0][0],b,b2,3,binterpt[0]);CHKERRQ(ierr);
295e27a552bSJed Brown   }
296ef3c5b88SJed Brown   {
297ef3c5b88SJed Brown     const PetscReal g = 0.5;
298ef3c5b88SJed Brown     const PetscReal
299ef3c5b88SJed Brown       A[4][4] = {{0,0,0,0},
300ef3c5b88SJed Brown                  {0,0,0,0},
301ef3c5b88SJed Brown                  {1.,0,0,0},
302ef3c5b88SJed Brown                  {0.75,-0.25,0.5,0}},
303ef3c5b88SJed Brown       Gamma[4][4] = {{g,0,0,0},
304ef3c5b88SJed Brown                      {1.,g,0,0},
305ef3c5b88SJed Brown                      {-0.25,-0.25,g,0},
306ef3c5b88SJed Brown                      {1./12,1./12,-2./3,g}},
307ef3c5b88SJed Brown       b[4] = {5./6,-1./6,-1./6,0.5},
308ef3c5b88SJed Brown       b2[4] = {0.75,-0.25,0.5,0};
309*f4aed992SEmil Constantinescu     ierr = TSRosWRegister(TSROSWRODAS3,3,4,&A[0][0],&Gamma[0][0],b,b2,0,PETSC_NULL);CHKERRQ(ierr);
310ef3c5b88SJed Brown   }
311ef3c5b88SJed Brown   {
312ef3c5b88SJed Brown     const PetscReal g = 0.43586652150845899941601945119356;
313ef3c5b88SJed Brown     const PetscReal
314ef3c5b88SJed Brown       A[3][3] = {{0,0,0},
315ef3c5b88SJed Brown                  {g,0,0},
316ef3c5b88SJed Brown                  {g,0,0}},
317ef3c5b88SJed Brown       Gamma[3][3] = {{g,0,0},
318ef3c5b88SJed Brown                      {-0.19294655696029095575009695436041,g,0},
319ef3c5b88SJed Brown                      {0,1.74927148125794685173529749738960,g}},
320ef3c5b88SJed Brown       b[3] = {-0.75457412385404315829818998646589,1.94100407061964420292840123379419,-0.18642994676560104463021124732829},
321ef3c5b88SJed Brown       b2[3] = {-1.53358745784149585370766523913002,2.81745131148625772213931745457622,-0.28386385364476186843165221544619};
322*f4aed992SEmil Constantinescu     ierr = TSRosWRegister(TSROSWSANDU3,3,3,&A[0][0],&Gamma[0][0],b,b2,0,PETSC_NULL);CHKERRQ(ierr);
323ef3c5b88SJed Brown   }
324b1c69cc3SEmil Constantinescu   {
325aaf9cf16SJed Brown     const PetscReal s3 = PetscSqrtReal(3.),g = (3.0+s3)/6.0;
326b1c69cc3SEmil Constantinescu     const PetscReal
327b1c69cc3SEmil Constantinescu       A[3][3] = {{0,0,0},
328b1c69cc3SEmil Constantinescu                  {1,0,0},
329b1c69cc3SEmil Constantinescu                  {0.25,0.25,0}},
330b1c69cc3SEmil Constantinescu       Gamma[3][3] = {{0,0,0},
331aaf9cf16SJed Brown                      {(-3.0-s3)/6.0,g,0},
332aaf9cf16SJed Brown                      {(-3.0-s3)/24.0,(-3.0-s3)/8.0,g}},
333b1c69cc3SEmil Constantinescu         b[3] = {1./6.,1./6.,2./3.},
334b1c69cc3SEmil Constantinescu           b2[3] = {1./4.,1./4.,1./2.};
335*f4aed992SEmil Constantinescu     ierr = TSRosWRegister(TSROSWASSP3P3S1C,3,3,&A[0][0],&Gamma[0][0],b,b2,0,PETSC_NULL);CHKERRQ(ierr);
336b1c69cc3SEmil Constantinescu   }
337b1c69cc3SEmil Constantinescu 
338b1c69cc3SEmil Constantinescu   {
339b1c69cc3SEmil Constantinescu     const PetscReal
340b1c69cc3SEmil Constantinescu       A[4][4] = {{0,0,0,0},
341b1c69cc3SEmil Constantinescu                  {1./2.,0,0,0},
342b1c69cc3SEmil Constantinescu                  {1./2.,1./2.,0,0},
343b1c69cc3SEmil Constantinescu                  {1./6.,1./6.,1./6.,0}},
344b1c69cc3SEmil Constantinescu       Gamma[4][4] = {{1./2.,0,0,0},
345b1c69cc3SEmil Constantinescu                      {0.0,1./4.,0,0},
346b1c69cc3SEmil Constantinescu                      {-2.,-2./3.,2./3.,0},
347b1c69cc3SEmil Constantinescu                      {1./2.,5./36.,-2./9,0}},
348b1c69cc3SEmil Constantinescu         b[4] = {1./6.,1./6.,1./6.,1./2.},
349b1c69cc3SEmil Constantinescu         b2[4] = {1./8.,3./4.,1./8.,0};
350*f4aed992SEmil Constantinescu      ierr = TSRosWRegister(TSROSWLASSP3P4S2C,3,4,&A[0][0],&Gamma[0][0],b,b2,0,PETSC_NULL);CHKERRQ(ierr);
351b1c69cc3SEmil Constantinescu   }
352b1c69cc3SEmil Constantinescu 
353b1c69cc3SEmil Constantinescu   {
354b1c69cc3SEmil Constantinescu     const PetscReal
355b1c69cc3SEmil Constantinescu       A[4][4] = {{0,0,0,0},
356b1c69cc3SEmil Constantinescu                  {1./2.,0,0,0},
357b1c69cc3SEmil Constantinescu                  {1./2.,1./2.,0,0},
358b1c69cc3SEmil Constantinescu                  {1./6.,1./6.,1./6.,0}},
359b1c69cc3SEmil Constantinescu       Gamma[4][4] = {{1./2.,0,0,0},
360b1c69cc3SEmil Constantinescu                      {0.0,3./4.,0,0},
361b1c69cc3SEmil Constantinescu                      {-2./3.,-23./9.,2./9.,0},
362b1c69cc3SEmil Constantinescu                      {1./18.,65./108.,-2./27,0}},
363b1c69cc3SEmil Constantinescu         b[4] = {1./6.,1./6.,1./6.,1./2.},
364b1c69cc3SEmil Constantinescu         b2[4] = {3./16.,10./16.,3./16.,0};
365*f4aed992SEmil Constantinescu      ierr = TSRosWRegister(TSROSWLLSSP3P4S2C,3,4,&A[0][0],&Gamma[0][0],b,b2,0,PETSC_NULL);CHKERRQ(ierr);
366b1c69cc3SEmil Constantinescu   }
367753f8adbSEmil Constantinescu 
368753f8adbSEmil Constantinescu  {
369753f8adbSEmil Constantinescu    PetscReal A[4][4],Gamma[4][4],b[4],b2[4];
370753f8adbSEmil Constantinescu 
371753f8adbSEmil Constantinescu    Gamma[0][0]=0.4358665215084589994160194475295062513822671686978816;
37205e8e825SJed Brown    Gamma[0][1]=0; Gamma[0][2]=0; Gamma[0][3]=0;
373753f8adbSEmil Constantinescu    Gamma[1][0]=-1.997527830934941248426324674704153457289527280554476;
374753f8adbSEmil Constantinescu    Gamma[1][1]=0.4358665215084589994160194475295062513822671686978816;
37505e8e825SJed Brown    Gamma[1][2]=0; Gamma[1][3]=0;
376753f8adbSEmil Constantinescu    Gamma[2][0]=-1.007948511795029620852002345345404191008352770119903;
377753f8adbSEmil Constantinescu    Gamma[2][1]=-0.004648958462629345562774289390054679806993396798458131;
378753f8adbSEmil Constantinescu    Gamma[2][2]=0.4358665215084589994160194475295062513822671686978816;
37905e8e825SJed Brown    Gamma[2][3]=0;
380753f8adbSEmil Constantinescu    Gamma[3][0]=-0.6685429734233467180451604600279552604364311322650783;
381753f8adbSEmil Constantinescu    Gamma[3][1]=0.6056625986449338476089525334450053439525178740492984;
382753f8adbSEmil Constantinescu    Gamma[3][2]=-0.9717899277217721234705114616271378792182450260943198;
383753f8adbSEmil Constantinescu    Gamma[3][3]=0;
384753f8adbSEmil Constantinescu 
38505e8e825SJed Brown    A[0][0]=0; A[0][1]=0; A[0][2]=0; A[0][3]=0;
386753f8adbSEmil Constantinescu    A[1][0]=0.8717330430169179988320388950590125027645343373957631;
38705e8e825SJed Brown    A[1][1]=0; A[1][2]=0; A[1][3]=0;
388753f8adbSEmil Constantinescu    A[2][0]=0.5275890119763004115618079766722914408876108660811028;
389753f8adbSEmil Constantinescu    A[2][1]=0.07241098802369958843819203208518599088698057726988732;
39005e8e825SJed Brown    A[2][2]=0; A[2][3]=0;
391753f8adbSEmil Constantinescu    A[3][0]=0.3990960076760701320627260685975778145384666450351314;
392753f8adbSEmil Constantinescu    A[3][1]=-0.4375576546135194437228463747348862825846903771419953;
393753f8adbSEmil Constantinescu    A[3][2]=1.038461646937449311660120300601880176655352737312713;
39405e8e825SJed Brown    A[3][3]=0;
395753f8adbSEmil Constantinescu 
396753f8adbSEmil Constantinescu    b[0]=0.1876410243467238251612921333138006734899663569186926;
397753f8adbSEmil Constantinescu    b[1]=-0.5952974735769549480478230473706443582188442040780541;
398753f8adbSEmil Constantinescu    b[2]=0.9717899277217721234705114616271378792182450260943198;
399753f8adbSEmil Constantinescu    b[3]=0.4358665215084589994160194475295062513822671686978816;
400753f8adbSEmil Constantinescu 
401753f8adbSEmil Constantinescu    b2[0]=0.2147402862233891404862383521089097657790734483804460;
402753f8adbSEmil Constantinescu    b2[1]=-0.4851622638849390928209050538171743017757490232519684;
403753f8adbSEmil Constantinescu    b2[2]=0.8687250025203875511662123688667549217531982787600080;
404753f8adbSEmil Constantinescu    b2[3]=0.4016969751411624011684543450940068201770721128357014;
405753f8adbSEmil Constantinescu 
406*f4aed992SEmil Constantinescu    PetscReal binterpt[4][3]={{2.2565812720167954547104627844105,-3.0826699111559187902922463354557,1.0137296634858471607430756831148},{1.349166413351089573796243820819,-2.4689115685996042534544925650515,0.52444768167155973161042570784064},{-2.4695174540533503758652847586647,5.7428279814696677152129332773553,-2.3015205996945452158771370439586},{-0.13623023131453465264142184656474,-0.19124650171414467146619437684812,0.76334325453713832352363565300308}};
407*f4aed992SEmil Constantinescu 
408*f4aed992SEmil Constantinescu    ierr = TSRosWRegister(TSROSWARK3,3,4,&A[0][0],&Gamma[0][0],b,b2,3,binterpt);CHKERRQ(ierr);
409753f8adbSEmil Constantinescu   }
410753f8adbSEmil Constantinescu 
411e27a552bSJed Brown   PetscFunctionReturn(0);
412e27a552bSJed Brown }
413e27a552bSJed Brown 
414e27a552bSJed Brown #undef __FUNCT__
415e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterDestroy"
416e27a552bSJed Brown /*@C
417e27a552bSJed Brown    TSRosWRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister().
418e27a552bSJed Brown 
419e27a552bSJed Brown    Not Collective
420e27a552bSJed Brown 
421e27a552bSJed Brown    Level: advanced
422e27a552bSJed Brown 
423e27a552bSJed Brown .keywords: TSRosW, register, destroy
424e27a552bSJed Brown .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegisterDynamic()
425e27a552bSJed Brown @*/
426e27a552bSJed Brown PetscErrorCode TSRosWRegisterDestroy(void)
427e27a552bSJed Brown {
428e27a552bSJed Brown   PetscErrorCode ierr;
42961692a83SJed Brown   RosWTableauLink link;
430e27a552bSJed Brown 
431e27a552bSJed Brown   PetscFunctionBegin;
43261692a83SJed Brown   while ((link = RosWTableauList)) {
43361692a83SJed Brown     RosWTableau t = &link->tab;
43461692a83SJed Brown     RosWTableauList = link->next;
43561692a83SJed Brown     ierr = PetscFree5(t->A,t->Gamma,t->b,t->ASum,t->GammaSum);CHKERRQ(ierr);
43643b21953SEmil Constantinescu     ierr = PetscFree5(t->At,t->bt,t->GammaInv,t->GammaZeroDiag,t->GammaExplicitCorr);CHKERRQ(ierr);
437fe7e6d57SJed Brown     ierr = PetscFree2(t->bembed,t->bembedt);CHKERRQ(ierr);
438*f4aed992SEmil Constantinescu     ierr = PetscFree(t->binterpt);CHKERRQ(ierr);
439e27a552bSJed Brown     ierr = PetscFree(t->name);CHKERRQ(ierr);
440e27a552bSJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
441e27a552bSJed Brown   }
442e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_FALSE;
443e27a552bSJed Brown   PetscFunctionReturn(0);
444e27a552bSJed Brown }
445e27a552bSJed Brown 
446e27a552bSJed Brown #undef __FUNCT__
447e27a552bSJed Brown #define __FUNCT__ "TSRosWInitializePackage"
448e27a552bSJed Brown /*@C
449e27a552bSJed Brown   TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called
450e27a552bSJed Brown   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RosW()
451e27a552bSJed Brown   when using static libraries.
452e27a552bSJed Brown 
453e27a552bSJed Brown   Input Parameter:
454e27a552bSJed Brown   path - The dynamic library path, or PETSC_NULL
455e27a552bSJed Brown 
456e27a552bSJed Brown   Level: developer
457e27a552bSJed Brown 
458e27a552bSJed Brown .keywords: TS, TSRosW, initialize, package
459e27a552bSJed Brown .seealso: PetscInitialize()
460e27a552bSJed Brown @*/
461e27a552bSJed Brown PetscErrorCode TSRosWInitializePackage(const char path[])
462e27a552bSJed Brown {
463e27a552bSJed Brown   PetscErrorCode ierr;
464e27a552bSJed Brown 
465e27a552bSJed Brown   PetscFunctionBegin;
466e27a552bSJed Brown   if (TSRosWPackageInitialized) PetscFunctionReturn(0);
467e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_TRUE;
468e27a552bSJed Brown   ierr = TSRosWRegisterAll();CHKERRQ(ierr);
469e27a552bSJed Brown   ierr = PetscRegisterFinalize(TSRosWFinalizePackage);CHKERRQ(ierr);
470e27a552bSJed Brown   PetscFunctionReturn(0);
471e27a552bSJed Brown }
472e27a552bSJed Brown 
473e27a552bSJed Brown #undef __FUNCT__
474e27a552bSJed Brown #define __FUNCT__ "TSRosWFinalizePackage"
475e27a552bSJed Brown /*@C
476e27a552bSJed Brown   TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is
477e27a552bSJed Brown   called from PetscFinalize().
478e27a552bSJed Brown 
479e27a552bSJed Brown   Level: developer
480e27a552bSJed Brown 
481e27a552bSJed Brown .keywords: Petsc, destroy, package
482e27a552bSJed Brown .seealso: PetscFinalize()
483e27a552bSJed Brown @*/
484e27a552bSJed Brown PetscErrorCode TSRosWFinalizePackage(void)
485e27a552bSJed Brown {
486e27a552bSJed Brown   PetscErrorCode ierr;
487e27a552bSJed Brown 
488e27a552bSJed Brown   PetscFunctionBegin;
489e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_FALSE;
490e27a552bSJed Brown   ierr = TSRosWRegisterDestroy();CHKERRQ(ierr);
491e27a552bSJed Brown   PetscFunctionReturn(0);
492e27a552bSJed Brown }
493e27a552bSJed Brown 
494e27a552bSJed Brown #undef __FUNCT__
495e27a552bSJed Brown #define __FUNCT__ "TSRosWRegister"
496e27a552bSJed Brown /*@C
49761692a83SJed Brown    TSRosWRegister - register a Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
498e27a552bSJed Brown 
499e27a552bSJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
500e27a552bSJed Brown 
501e27a552bSJed Brown    Input Parameters:
502e27a552bSJed Brown +  name - identifier for method
503e27a552bSJed Brown .  order - approximation order of method
504e27a552bSJed Brown .  s - number of stages, this is the dimension of the matrices below
50561692a83SJed Brown .  A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular
50661692a83SJed Brown .  Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal
507fe7e6d57SJed Brown .  b - Step completion table (dimension s)
508fe7e6d57SJed Brown -  bembed - Step completion table for a scheme of order one less (dimension s, PETSC_NULL if no embedded scheme is available)
509*f4aed992SEmil Constantinescu .  pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt
510*f4aed992SEmil Constantinescu .  binterpt - Coefficients of the interpolation formula (dimension s*pinterp)
511e27a552bSJed Brown 
512e27a552bSJed Brown    Notes:
51361692a83SJed Brown    Several Rosenbrock W methods are provided, this function is only needed to create new methods.
514e27a552bSJed Brown 
515e27a552bSJed Brown    Level: advanced
516e27a552bSJed Brown 
517e27a552bSJed Brown .keywords: TS, register
518e27a552bSJed Brown 
519e27a552bSJed Brown .seealso: TSRosW
520e27a552bSJed Brown @*/
521e27a552bSJed Brown PetscErrorCode TSRosWRegister(const TSRosWType name,PetscInt order,PetscInt s,
522*f4aed992SEmil Constantinescu                               const PetscReal A[],const PetscReal Gamma[],const PetscReal b[],const PetscReal bembed[],
523*f4aed992SEmil Constantinescu                                  PetscInt pinterp,const PetscReal binterpt[])
524e27a552bSJed Brown {
525e27a552bSJed Brown   PetscErrorCode ierr;
52661692a83SJed Brown   RosWTableauLink link;
52761692a83SJed Brown   RosWTableau t;
52861692a83SJed Brown   PetscInt i,j,k;
52961692a83SJed Brown   PetscScalar *GammaInv;
530e27a552bSJed Brown 
531e27a552bSJed Brown   PetscFunctionBegin;
532fe7e6d57SJed Brown   PetscValidCharPointer(name,1);
533fe7e6d57SJed Brown   PetscValidPointer(A,4);
534fe7e6d57SJed Brown   PetscValidPointer(Gamma,5);
535fe7e6d57SJed Brown   PetscValidPointer(b,6);
536fe7e6d57SJed Brown   if (bembed) PetscValidPointer(bembed,7);
537fe7e6d57SJed Brown 
538e27a552bSJed Brown   ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
539e27a552bSJed Brown   ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
540e27a552bSJed Brown   t = &link->tab;
541e27a552bSJed Brown   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
542e27a552bSJed Brown   t->order = order;
543e27a552bSJed Brown   t->s = s;
54461692a83SJed Brown   ierr = PetscMalloc5(s*s,PetscReal,&t->A,s*s,PetscReal,&t->Gamma,s,PetscReal,&t->b,s,PetscReal,&t->ASum,s,PetscReal,&t->GammaSum);CHKERRQ(ierr);
54543b21953SEmil Constantinescu   ierr = PetscMalloc5(s*s,PetscReal,&t->At,s,PetscReal,&t->bt,s*s,PetscReal,&t->GammaInv,s,PetscBool,&t->GammaZeroDiag,s*s,PetscReal,&t->GammaExplicitCorr);CHKERRQ(ierr);
546e27a552bSJed Brown   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
54761692a83SJed Brown   ierr = PetscMemcpy(t->Gamma,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr);
54843b21953SEmil Constantinescu   ierr = PetscMemcpy(t->GammaExplicitCorr,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr);
54961692a83SJed Brown   ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
550fe7e6d57SJed Brown   if (bembed) {
551fe7e6d57SJed Brown     ierr = PetscMalloc2(s,PetscReal,&t->bembed,s,PetscReal,&t->bembedt);CHKERRQ(ierr);
552fe7e6d57SJed Brown     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
553fe7e6d57SJed Brown   }
55461692a83SJed Brown   for (i=0; i<s; i++) {
55561692a83SJed Brown     t->ASum[i] = 0;
55661692a83SJed Brown     t->GammaSum[i] = 0;
55761692a83SJed Brown     for (j=0; j<s; j++) {
55861692a83SJed Brown       t->ASum[i] += A[i*s+j];
559fe7e6d57SJed Brown       t->GammaSum[i] += Gamma[i*s+j];
56061692a83SJed Brown     }
56161692a83SJed Brown   }
56261692a83SJed Brown   ierr = PetscMalloc(s*s*sizeof(PetscScalar),&GammaInv);CHKERRQ(ierr); /* Need to use Scalar for inverse, then convert back to Real */
56361692a83SJed Brown   for (i=0; i<s*s; i++) GammaInv[i] = Gamma[i];
564fd96d5b0SEmil Constantinescu   for (i=0; i<s; i++) {
565fd96d5b0SEmil Constantinescu     if (Gamma[i*s+i] == 0.0) {
566fd96d5b0SEmil Constantinescu       GammaInv[i*s+i] = 1.0;
567c17803e7SJed Brown       t->GammaZeroDiag[i] = PETSC_TRUE;
568fd96d5b0SEmil Constantinescu     } else {
569c17803e7SJed Brown       t->GammaZeroDiag[i] = PETSC_FALSE;
570fd96d5b0SEmil Constantinescu     }
571fd96d5b0SEmil Constantinescu   }
572fd96d5b0SEmil Constantinescu 
57361692a83SJed Brown   switch (s) {
57461692a83SJed Brown   case 1: GammaInv[0] = 1./GammaInv[0]; break;
57561692a83SJed Brown   case 2: ierr = Kernel_A_gets_inverse_A_2(GammaInv,0);CHKERRQ(ierr); break;
57661692a83SJed Brown   case 3: ierr = Kernel_A_gets_inverse_A_3(GammaInv,0);CHKERRQ(ierr); break;
57761692a83SJed Brown   case 4: ierr = Kernel_A_gets_inverse_A_4(GammaInv,0);CHKERRQ(ierr); break;
57861692a83SJed Brown   case 5: {
57961692a83SJed Brown     PetscInt ipvt5[5];
58061692a83SJed Brown     MatScalar work5[5*5];
58161692a83SJed Brown     ierr = Kernel_A_gets_inverse_A_5(GammaInv,ipvt5,work5,0);CHKERRQ(ierr); break;
58261692a83SJed Brown   }
58361692a83SJed Brown   case 6: ierr = Kernel_A_gets_inverse_A_6(GammaInv,0);CHKERRQ(ierr); break;
58461692a83SJed Brown   case 7: ierr = Kernel_A_gets_inverse_A_7(GammaInv,0);CHKERRQ(ierr); break;
58561692a83SJed Brown   default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for %D stages",s);
58661692a83SJed Brown   }
58761692a83SJed Brown   for (i=0; i<s*s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]);
58861692a83SJed Brown   ierr = PetscFree(GammaInv);CHKERRQ(ierr);
58943b21953SEmil Constantinescu 
59043b21953SEmil Constantinescu   for (i=0; i<s; i++) {
59143b21953SEmil Constantinescu     for (k=0; k<i+1; k++) {
59243b21953SEmil Constantinescu       t->GammaExplicitCorr[i*s+k]=(t->GammaExplicitCorr[i*s+k])*(t->GammaInv[k*s+k]);
59343b21953SEmil Constantinescu       for (j=k+1; j<i+1; j++) {
59443b21953SEmil Constantinescu         t->GammaExplicitCorr[i*s+k]+=(t->GammaExplicitCorr[i*s+j])*(t->GammaInv[j*s+k]);
59543b21953SEmil Constantinescu       }
59643b21953SEmil Constantinescu     }
59743b21953SEmil Constantinescu   }
59843b21953SEmil Constantinescu 
59961692a83SJed Brown   for (i=0; i<s; i++) {
60061692a83SJed Brown     for (j=0; j<s; j++) {
60161692a83SJed Brown       t->At[i*s+j] = 0;
60261692a83SJed Brown       for (k=0; k<s; k++) {
60361692a83SJed Brown         t->At[i*s+j] += t->A[i*s+k] * t->GammaInv[k*s+j];
60461692a83SJed Brown       }
60561692a83SJed Brown     }
60661692a83SJed Brown     t->bt[i] = 0;
60761692a83SJed Brown     for (j=0; j<s; j++) {
60861692a83SJed Brown       t->bt[i] += t->b[j] * t->GammaInv[j*s+i];
60961692a83SJed Brown     }
610fe7e6d57SJed Brown     if (bembed) {
611fe7e6d57SJed Brown       t->bembedt[i] = 0;
612fe7e6d57SJed Brown       for (j=0; j<s; j++) {
613fe7e6d57SJed Brown         t->bembedt[i] += t->bembed[j] * t->GammaInv[j*s+i];
614fe7e6d57SJed Brown       }
615fe7e6d57SJed Brown     }
61661692a83SJed Brown   }
6178d59e960SJed Brown   t->ccfl = 1.0;                /* Fix this */
6188d59e960SJed Brown 
619*f4aed992SEmil Constantinescu   t->pinterp = pinterp;
620*f4aed992SEmil Constantinescu   ierr = PetscMalloc(s*pinterp,&t->binterpt);CHKERRQ(ierr);
621*f4aed992SEmil Constantinescu 
62261692a83SJed Brown   link->next = RosWTableauList;
62361692a83SJed Brown   RosWTableauList = link;
624e27a552bSJed Brown   PetscFunctionReturn(0);
625e27a552bSJed Brown }
626e27a552bSJed Brown 
627e27a552bSJed Brown #undef __FUNCT__
6281c3436cfSJed Brown #define __FUNCT__ "TSEvaluateStep_RosW"
6291c3436cfSJed Brown /*
6301c3436cfSJed Brown  The step completion formula is
6311c3436cfSJed Brown 
6321c3436cfSJed Brown  x1 = x0 + b^T Y
6331c3436cfSJed Brown 
6341c3436cfSJed Brown  where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been
6351c3436cfSJed Brown  updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write
6361c3436cfSJed Brown 
6371c3436cfSJed Brown  x1e = x0 + be^T Y
6381c3436cfSJed Brown      = x1 - b^T Y + be^T Y
6391c3436cfSJed Brown      = x1 + (be - b)^T Y
6401c3436cfSJed Brown 
6411c3436cfSJed Brown  so we can evaluate the method of different order even after the step has been optimistically completed.
6421c3436cfSJed Brown */
6431c3436cfSJed Brown static PetscErrorCode TSEvaluateStep_RosW(TS ts,PetscInt order,Vec X,PetscBool *done)
6441c3436cfSJed Brown {
6451c3436cfSJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
6461c3436cfSJed Brown   RosWTableau    tab  = ros->tableau;
6471c3436cfSJed Brown   PetscScalar    *w = ros->work;
6481c3436cfSJed Brown   PetscInt       i;
6491c3436cfSJed Brown   PetscErrorCode ierr;
6501c3436cfSJed Brown 
6511c3436cfSJed Brown   PetscFunctionBegin;
6521c3436cfSJed Brown   if (order == tab->order) {
653108c343cSJed Brown     if (ros->status == TS_STEP_INCOMPLETE) { /* Use standard completion formula */
6541c3436cfSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
655de19f811SJed Brown       for (i=0; i<tab->s; i++) w[i] = tab->bt[i];
656de19f811SJed Brown       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
657108c343cSJed Brown     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
6581c3436cfSJed Brown     if (done) *done = PETSC_TRUE;
6591c3436cfSJed Brown     PetscFunctionReturn(0);
6601c3436cfSJed Brown   } else if (order == tab->order-1) {
6611c3436cfSJed Brown     if (!tab->bembedt) goto unavailable;
662108c343cSJed Brown     if (ros->status == TS_STEP_INCOMPLETE) { /* Use embedded completion formula */
6631c3436cfSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
664de19f811SJed Brown       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i];
665de19f811SJed Brown       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
666108c343cSJed Brown     } else {                    /* Use rollback-and-recomplete formula (bembedt - bt) */
667108c343cSJed Brown       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i];
668108c343cSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
669108c343cSJed Brown       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
6701c3436cfSJed Brown     }
6711c3436cfSJed Brown     if (done) *done = PETSC_TRUE;
6721c3436cfSJed Brown     PetscFunctionReturn(0);
6731c3436cfSJed Brown   }
6741c3436cfSJed Brown   unavailable:
6751c3436cfSJed Brown   if (done) *done = PETSC_FALSE;
6761c3436cfSJed Brown   else SETERRQ3(((PetscObject)ts)->comm,PETSC_ERR_SUP,"Rosenbrock-W '%s' of order %D cannot evaluate step at order %D",tab->name,tab->order,order);
6771c3436cfSJed Brown   PetscFunctionReturn(0);
6781c3436cfSJed Brown }
6791c3436cfSJed Brown 
6801c3436cfSJed Brown #undef __FUNCT__
681e27a552bSJed Brown #define __FUNCT__ "TSStep_RosW"
682e27a552bSJed Brown static PetscErrorCode TSStep_RosW(TS ts)
683e27a552bSJed Brown {
68461692a83SJed Brown   TS_RosW         *ros = (TS_RosW*)ts->data;
68561692a83SJed Brown   RosWTableau     tab  = ros->tableau;
686e27a552bSJed Brown   const PetscInt  s    = tab->s;
6871c3436cfSJed Brown   const PetscReal *At  = tab->At,*Gamma = tab->Gamma,*ASum = tab->ASum,*GammaInv = tab->GammaInv;
6880feba352SEmil Constantinescu   const PetscReal *GammaExplicitCorr = tab->GammaExplicitCorr;
689c17803e7SJed Brown   const PetscBool *GammaZeroDiag = tab->GammaZeroDiag;
69061692a83SJed Brown   PetscScalar     *w   = ros->work;
6917d4bf2deSEmil Constantinescu   Vec             *Y   = ros->Y,Ydot = ros->Ydot,Zdot = ros->Zdot,Zstage = ros->Zstage;
692e27a552bSJed Brown   SNES            snes;
6931c3436cfSJed Brown   TSAdapt         adapt;
6941c3436cfSJed Brown   PetscInt        i,j,its,lits,reject,next_scheme;
695cdbf8f93SLisandro Dalcin   PetscReal       next_time_step;
6961c3436cfSJed Brown   PetscBool       accept;
697e27a552bSJed Brown   PetscErrorCode  ierr;
6980feba352SEmil Constantinescu   MatStructure    str;
699e27a552bSJed Brown 
700e27a552bSJed Brown   PetscFunctionBegin;
701e27a552bSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
702cdbf8f93SLisandro Dalcin   next_time_step = ts->time_step;
7031c3436cfSJed Brown   accept = PETSC_TRUE;
704108c343cSJed Brown   ros->status = TS_STEP_INCOMPLETE;
705e27a552bSJed Brown 
70697335746SJed Brown   for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) {
7071c3436cfSJed Brown     const PetscReal h = ts->time_step;
708e27a552bSJed Brown     for (i=0; i<s; i++) {
7091c3436cfSJed Brown       ros->stage_time = ts->ptime + h*ASum[i];
710c17803e7SJed Brown       if (GammaZeroDiag[i]) {
711c17803e7SJed Brown         ros->stage_explicit = PETSC_TRUE;
712fd96d5b0SEmil Constantinescu         ros->shift = 1./h;
713c17803e7SJed Brown       } else {
714c17803e7SJed Brown         ros->stage_explicit = PETSC_FALSE;
71561692a83SJed Brown         ros->shift = 1./(h*Gamma[i*s+i]);
716fd96d5b0SEmil Constantinescu       }
71761692a83SJed Brown 
71861692a83SJed Brown       ierr = VecCopy(ts->vec_sol,Zstage);CHKERRQ(ierr);
719de19f811SJed Brown       for (j=0; j<i; j++) w[j] = At[i*s+j];
720de19f811SJed Brown       ierr = VecMAXPY(Zstage,i,w,Y);CHKERRQ(ierr);
72161692a83SJed Brown 
72261692a83SJed Brown       for (j=0; j<i; j++) w[j] = 1./h * GammaInv[i*s+j];
72361692a83SJed Brown       ierr = VecZeroEntries(Zdot);CHKERRQ(ierr);
72461692a83SJed Brown       ierr = VecMAXPY(Zdot,i,w,Y);CHKERRQ(ierr);
72561692a83SJed Brown 
726e27a552bSJed Brown       /* Initial guess taken from last stage */
72761692a83SJed Brown       ierr = VecZeroEntries(Y[i]);CHKERRQ(ierr);
72861692a83SJed Brown 
7297d4bf2deSEmil Constantinescu       if (!ros->stage_explicit) {
73061692a83SJed Brown         if (!ros->recompute_jacobian && !i) {
73161692a83SJed Brown           ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */
73261692a83SJed Brown         }
73361692a83SJed Brown         ierr = SNESSolve(snes,PETSC_NULL,Y[i]);CHKERRQ(ierr);
734e27a552bSJed Brown         ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
735e27a552bSJed Brown         ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
736e27a552bSJed Brown         ts->nonlinear_its += its; ts->linear_its += lits;
73797335746SJed Brown         ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
73897335746SJed Brown         ierr = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr);
73997335746SJed Brown         if (!accept) goto reject_step;
7407d4bf2deSEmil Constantinescu       } else {
7411ce71dffSSatish Balay         Mat J,Jp;
7420feba352SEmil Constantinescu         ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); /* Evaluate Y[i]=G(t,Ydot=0,Zstage) */
7430feba352SEmil Constantinescu         ierr = TSComputeIFunction(ts,ros->stage_time,Zstage,Ydot,Y[i],PETSC_FALSE);CHKERRQ(ierr);
7440feba352SEmil Constantinescu         ierr = VecScale(Y[i],-1.0);
7450feba352SEmil Constantinescu         ierr = VecAXPY(Y[i],-1.0,Zdot);CHKERRQ(ierr); /*Y[i]=F(Zstage)-Zdot[=GammaInv*Y]*/
7460feba352SEmil Constantinescu 
7470feba352SEmil Constantinescu         ierr = VecZeroEntries(Zstage);CHKERRQ(ierr); /* Zstage = GammaExplicitCorr[i,j] * Y[j] */
7480feba352SEmil Constantinescu         for (j=0; j<i; j++) w[j] = GammaExplicitCorr[i*s+j];
7490feba352SEmil Constantinescu         ierr = VecMAXPY(Zstage,i,w,Y);CHKERRQ(ierr);
7500feba352SEmil Constantinescu         /*Y[i] += Y[i] + Jac*Zstage[=Jac*GammaExplicitCorr[i,j] * Y[j]] */
7510feba352SEmil Constantinescu         str = SAME_NONZERO_PATTERN;
7520feba352SEmil Constantinescu         ierr = TSGetIJacobian(ts,&J,&Jp,PETSC_NULL,PETSC_NULL);
7530feba352SEmil Constantinescu         ierr = TSComputeIJacobian(ts,ros->stage_time,ts->vec_sol,Ydot,0,&J,&Jp,&str,PETSC_FALSE);CHKERRQ(ierr);
7540feba352SEmil Constantinescu         ierr = MatMult(J,Zstage,Zdot);
7550feba352SEmil Constantinescu 
7560feba352SEmil Constantinescu         ierr = VecAXPY(Y[i],-1.0,Zdot);CHKERRQ(ierr);
7570feba352SEmil Constantinescu         ierr = VecScale(Y[i],h);
7587d4bf2deSEmil Constantinescu         ts->linear_its += 1;
7597d4bf2deSEmil Constantinescu       }
760e27a552bSJed Brown     }
7611c3436cfSJed Brown     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
762108c343cSJed Brown     ros->status = TS_STEP_PENDING;
763e27a552bSJed Brown 
7641c3436cfSJed Brown     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
7651c3436cfSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
7661c3436cfSJed Brown     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
7678d59e960SJed Brown     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
7681c3436cfSJed Brown     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
7691c3436cfSJed Brown     if (accept) {
7701c3436cfSJed Brown       /* ignore next_scheme for now */
771e27a552bSJed Brown       ts->ptime += ts->time_step;
772cdbf8f93SLisandro Dalcin       ts->time_step = next_time_step;
773e27a552bSJed Brown       ts->steps++;
774108c343cSJed Brown       ros->status = TS_STEP_COMPLETE;
7751c3436cfSJed Brown       break;
7761c3436cfSJed Brown     } else {                    /* Roll back the current step */
7771c3436cfSJed Brown       for (i=0; i<s; i++) w[i] = -tab->bt[i];
7781c3436cfSJed Brown       ierr = VecMAXPY(ts->vec_sol,s,w,Y);CHKERRQ(ierr);
7791c3436cfSJed Brown       ts->time_step = next_time_step;
780108c343cSJed Brown       ros->status = TS_STEP_INCOMPLETE;
7811c3436cfSJed Brown     }
782476b6736SJed Brown     reject_step: continue;
7831c3436cfSJed Brown   }
784b2ce242eSJed Brown   if (ros->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
785e27a552bSJed Brown   PetscFunctionReturn(0);
786e27a552bSJed Brown }
787e27a552bSJed Brown 
788e27a552bSJed Brown #undef __FUNCT__
789e27a552bSJed Brown #define __FUNCT__ "TSInterpolate_RosW"
790e27a552bSJed Brown static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec X)
791e27a552bSJed Brown {
79261692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
793*f4aed992SEmil Constantinescu   PetscInt s = ros->tableau->s,pinterp = ros->tableau->pinterp,i,j;
794*f4aed992SEmil Constantinescu   PetscReal h;
795*f4aed992SEmil Constantinescu   PetscReal tt,t;
796*f4aed992SEmil Constantinescu   PetscScalar *bt;
797*f4aed992SEmil Constantinescu   const PetscReal *Bt = ros->tableau->binterpt;
798*f4aed992SEmil Constantinescu   PetscErrorCode ierr;
799*f4aed992SEmil Constantinescu   const PetscReal *GammaInv = ros->tableau->GammaInv;
800*f4aed992SEmil Constantinescu   PetscScalar     *w   = ros->work;
801*f4aed992SEmil Constantinescu   Vec             *Y   = ros->Y;
802e27a552bSJed Brown 
803e27a552bSJed Brown   PetscFunctionBegin;
804*f4aed992SEmil Constantinescu   //  SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name);
805*f4aed992SEmil Constantinescu   if (!Bt) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name);
806*f4aed992SEmil Constantinescu 
807*f4aed992SEmil Constantinescu   switch (ros->status) {
808*f4aed992SEmil Constantinescu   case TS_STEP_INCOMPLETE:
809*f4aed992SEmil Constantinescu   case TS_STEP_PENDING:
810*f4aed992SEmil Constantinescu     h = ts->time_step;
811*f4aed992SEmil Constantinescu     t = (itime - ts->ptime)/h;
812*f4aed992SEmil Constantinescu     break;
813*f4aed992SEmil Constantinescu   case TS_STEP_COMPLETE:
814*f4aed992SEmil Constantinescu     h = ts->time_step_prev;
815*f4aed992SEmil Constantinescu     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
816*f4aed992SEmil Constantinescu     break;
817*f4aed992SEmil Constantinescu   default: SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Invalid TSStepStatus");
818*f4aed992SEmil Constantinescu   }
819*f4aed992SEmil Constantinescu   ierr = PetscMalloc(s,&bt);CHKERRQ(ierr);
820*f4aed992SEmil Constantinescu   for (i=0; i<s; i++) bt[i] = 0;
821*f4aed992SEmil Constantinescu   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
822*f4aed992SEmil Constantinescu     for (i=0; i<s; i++) {
823*f4aed992SEmil Constantinescu       bt[i] += h * Bt[i*pinterp+j] * tt;
824*f4aed992SEmil Constantinescu     }
825*f4aed992SEmil Constantinescu   }
826*f4aed992SEmil Constantinescu 
827*f4aed992SEmil Constantinescu   /* y(t+tt*h) = y(t) + Sum bt(tt) * GammaInv * Ydot */
828*f4aed992SEmil Constantinescu   /*X<-0*/
829*f4aed992SEmil Constantinescu   ierr = VecZeroEntries(X);CHKERRQ(ierr);
830*f4aed992SEmil Constantinescu 
831*f4aed992SEmil Constantinescu   /*X<- Sum bt_i * GammaInv(i,1:i) * Y(1:i) */
832*f4aed992SEmil Constantinescu   for (i=0; i<s; i++) {
833*f4aed992SEmil Constantinescu     for (j=0; j<=i; j++) w[j] =  bt[i]*GammaInv[i*s+j];
834*f4aed992SEmil Constantinescu     ierr = VecMAXPY(X,i,w,Y);CHKERRQ(ierr);
835*f4aed992SEmil Constantinescu   }
836*f4aed992SEmil Constantinescu 
837*f4aed992SEmil Constantinescu   /*X<-y(t) + X*/
838*f4aed992SEmil Constantinescu   ierr = VecAXPY(X,1.0,ts->vec_sol);CHKERRQ(ierr);
839*f4aed992SEmil Constantinescu 
840*f4aed992SEmil Constantinescu   ierr = PetscFree(bt);CHKERRQ(ierr);
841*f4aed992SEmil Constantinescu 
842e27a552bSJed Brown   PetscFunctionReturn(0);
843e27a552bSJed Brown }
844e27a552bSJed Brown 
845e27a552bSJed Brown /*------------------------------------------------------------*/
846e27a552bSJed Brown #undef __FUNCT__
847e27a552bSJed Brown #define __FUNCT__ "TSReset_RosW"
848e27a552bSJed Brown static PetscErrorCode TSReset_RosW(TS ts)
849e27a552bSJed Brown {
85061692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
851e27a552bSJed Brown   PetscInt       s;
852e27a552bSJed Brown   PetscErrorCode ierr;
853e27a552bSJed Brown 
854e27a552bSJed Brown   PetscFunctionBegin;
85561692a83SJed Brown   if (!ros->tableau) PetscFunctionReturn(0);
85661692a83SJed Brown    s = ros->tableau->s;
85761692a83SJed Brown   ierr = VecDestroyVecs(s,&ros->Y);CHKERRQ(ierr);
85861692a83SJed Brown   ierr = VecDestroy(&ros->Ydot);CHKERRQ(ierr);
85961692a83SJed Brown   ierr = VecDestroy(&ros->Ystage);CHKERRQ(ierr);
86061692a83SJed Brown   ierr = VecDestroy(&ros->Zdot);CHKERRQ(ierr);
86161692a83SJed Brown   ierr = VecDestroy(&ros->Zstage);CHKERRQ(ierr);
86261692a83SJed Brown   ierr = PetscFree(ros->work);CHKERRQ(ierr);
863e27a552bSJed Brown   PetscFunctionReturn(0);
864e27a552bSJed Brown }
865e27a552bSJed Brown 
866e27a552bSJed Brown #undef __FUNCT__
867e27a552bSJed Brown #define __FUNCT__ "TSDestroy_RosW"
868e27a552bSJed Brown static PetscErrorCode TSDestroy_RosW(TS ts)
869e27a552bSJed Brown {
870e27a552bSJed Brown   PetscErrorCode  ierr;
871e27a552bSJed Brown 
872e27a552bSJed Brown   PetscFunctionBegin;
873e27a552bSJed Brown   ierr = TSReset_RosW(ts);CHKERRQ(ierr);
874e27a552bSJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
875e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","",PETSC_NULL);CHKERRQ(ierr);
876e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","",PETSC_NULL);CHKERRQ(ierr);
87761692a83SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","",PETSC_NULL);CHKERRQ(ierr);
878e27a552bSJed Brown   PetscFunctionReturn(0);
879e27a552bSJed Brown }
880e27a552bSJed Brown 
881e27a552bSJed Brown /*
882e27a552bSJed Brown   This defines the nonlinear equation that is to be solved with SNES
883e27a552bSJed Brown   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
884e27a552bSJed Brown */
885e27a552bSJed Brown #undef __FUNCT__
886e27a552bSJed Brown #define __FUNCT__ "SNESTSFormFunction_RosW"
887e27a552bSJed Brown static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec X,Vec F,TS ts)
888e27a552bSJed Brown {
88961692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
890e27a552bSJed Brown   PetscErrorCode ierr;
891e27a552bSJed Brown 
892e27a552bSJed Brown   PetscFunctionBegin;
89361692a83SJed Brown   ierr = VecWAXPY(ros->Ydot,ros->shift,X,ros->Zdot);CHKERRQ(ierr); /* Ydot = shift*X + Zdot */
89461692a83SJed Brown   ierr = VecWAXPY(ros->Ystage,1.0,X,ros->Zstage);CHKERRQ(ierr);    /* Ystage = X + Zstage */
89561692a83SJed Brown   ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,ros->Ydot,F,PETSC_FALSE);CHKERRQ(ierr);
896e27a552bSJed Brown   PetscFunctionReturn(0);
897e27a552bSJed Brown }
898e27a552bSJed Brown 
899e27a552bSJed Brown #undef __FUNCT__
900e27a552bSJed Brown #define __FUNCT__ "SNESTSFormJacobian_RosW"
901e27a552bSJed Brown static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
902e27a552bSJed Brown {
90361692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
904e27a552bSJed Brown   PetscErrorCode ierr;
905e27a552bSJed Brown 
906e27a552bSJed Brown   PetscFunctionBegin;
90761692a83SJed Brown   /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
90861692a83SJed Brown   ierr = TSComputeIJacobian(ts,ros->stage_time,ros->Ystage,ros->Ydot,ros->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr);
909e27a552bSJed Brown   PetscFunctionReturn(0);
910e27a552bSJed Brown }
911e27a552bSJed Brown 
912e27a552bSJed Brown #undef __FUNCT__
913e27a552bSJed Brown #define __FUNCT__ "TSSetUp_RosW"
914e27a552bSJed Brown static PetscErrorCode TSSetUp_RosW(TS ts)
915e27a552bSJed Brown {
91661692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
91761692a83SJed Brown   RosWTableau    tab  = ros->tableau;
918e27a552bSJed Brown   PetscInt       s    = tab->s;
919e27a552bSJed Brown   PetscErrorCode ierr;
920e27a552bSJed Brown 
921e27a552bSJed Brown   PetscFunctionBegin;
92261692a83SJed Brown   if (!ros->tableau) {
923e27a552bSJed Brown     ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);
924e27a552bSJed Brown   }
92561692a83SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ros->Y);CHKERRQ(ierr);
92661692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Ydot);CHKERRQ(ierr);
92761692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Ystage);CHKERRQ(ierr);
92861692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Zdot);CHKERRQ(ierr);
92961692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Zstage);CHKERRQ(ierr);
93061692a83SJed Brown   ierr = PetscMalloc(s*sizeof(ros->work[0]),&ros->work);CHKERRQ(ierr);
931e27a552bSJed Brown   PetscFunctionReturn(0);
932e27a552bSJed Brown }
933e27a552bSJed Brown /*------------------------------------------------------------*/
934e27a552bSJed Brown 
935e27a552bSJed Brown #undef __FUNCT__
936e27a552bSJed Brown #define __FUNCT__ "TSSetFromOptions_RosW"
937e27a552bSJed Brown static PetscErrorCode TSSetFromOptions_RosW(TS ts)
938e27a552bSJed Brown {
93961692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
940e27a552bSJed Brown   PetscErrorCode ierr;
94161692a83SJed Brown   char           rostype[256];
942e27a552bSJed Brown 
943e27a552bSJed Brown   PetscFunctionBegin;
944e27a552bSJed Brown   ierr = PetscOptionsHead("RosW ODE solver options");CHKERRQ(ierr);
945e27a552bSJed Brown   {
94661692a83SJed Brown     RosWTableauLink link;
947e27a552bSJed Brown     PetscInt count,choice;
948e27a552bSJed Brown     PetscBool flg;
949e27a552bSJed Brown     const char **namelist;
95061692a83SJed Brown     SNES snes;
95161692a83SJed Brown 
95261692a83SJed Brown     ierr = PetscStrncpy(rostype,TSRosWDefault,sizeof rostype);CHKERRQ(ierr);
95361692a83SJed Brown     for (link=RosWTableauList,count=0; link; link=link->next,count++) ;
954e27a552bSJed Brown     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
95561692a83SJed Brown     for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
95661692a83SJed Brown     ierr = PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,rostype,&choice,&flg);CHKERRQ(ierr);
95761692a83SJed Brown     ierr = TSRosWSetType(ts,flg ? namelist[choice] : rostype);CHKERRQ(ierr);
958e27a552bSJed Brown     ierr = PetscFree(namelist);CHKERRQ(ierr);
95961692a83SJed Brown 
96061692a83SJed Brown     ierr = PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,PETSC_NULL);CHKERRQ(ierr);
96161692a83SJed Brown 
96261692a83SJed Brown     /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
96361692a83SJed Brown     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
96461692a83SJed Brown     if (!((PetscObject)snes)->type_name) {
96561692a83SJed Brown       ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
96661692a83SJed Brown     }
96761692a83SJed Brown     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
968e27a552bSJed Brown   }
969e27a552bSJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
970e27a552bSJed Brown   PetscFunctionReturn(0);
971e27a552bSJed Brown }
972e27a552bSJed Brown 
973e27a552bSJed Brown #undef __FUNCT__
974e27a552bSJed Brown #define __FUNCT__ "PetscFormatRealArray"
975e27a552bSJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
976e27a552bSJed Brown {
977e27a552bSJed Brown   PetscErrorCode ierr;
978e408995aSJed Brown   PetscInt i;
979e408995aSJed Brown   size_t left,count;
980e27a552bSJed Brown   char *p;
981e27a552bSJed Brown 
982e27a552bSJed Brown   PetscFunctionBegin;
983e408995aSJed Brown   for (i=0,p=buf,left=len; i<n; i++) {
984e408995aSJed Brown     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
985e27a552bSJed Brown     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
986e27a552bSJed Brown     left -= count;
987e27a552bSJed Brown     p += count;
988e27a552bSJed Brown     *p++ = ' ';
989e27a552bSJed Brown   }
990e27a552bSJed Brown   p[i ? 0 : -1] = 0;
991e27a552bSJed Brown   PetscFunctionReturn(0);
992e27a552bSJed Brown }
993e27a552bSJed Brown 
994e27a552bSJed Brown #undef __FUNCT__
995e27a552bSJed Brown #define __FUNCT__ "TSView_RosW"
996e27a552bSJed Brown static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer)
997e27a552bSJed Brown {
99861692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
99961692a83SJed Brown   RosWTableau    tab  = ros->tableau;
1000e27a552bSJed Brown   PetscBool      iascii;
1001e27a552bSJed Brown   PetscErrorCode ierr;
1002e27a552bSJed Brown 
1003e27a552bSJed Brown   PetscFunctionBegin;
1004e27a552bSJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1005e27a552bSJed Brown   if (iascii) {
100661692a83SJed Brown     const TSRosWType rostype;
1007e408995aSJed Brown     PetscInt i;
1008e408995aSJed Brown     PetscReal abscissa[512];
1009e27a552bSJed Brown     char buf[512];
101061692a83SJed Brown     ierr = TSRosWGetType(ts,&rostype);CHKERRQ(ierr);
101161692a83SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Rosenbrock-W %s\n",rostype);CHKERRQ(ierr);
1012e408995aSJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ASum);CHKERRQ(ierr);
101361692a83SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A       = %s\n",buf);CHKERRQ(ierr);
1014e408995aSJed Brown     for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i];
1015e408995aSJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,abscissa);CHKERRQ(ierr);
1016e408995aSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A+Gamma = %s\n",buf);CHKERRQ(ierr);
1017e27a552bSJed Brown   }
1018e27a552bSJed Brown   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
1019e27a552bSJed Brown   PetscFunctionReturn(0);
1020e27a552bSJed Brown }
1021e27a552bSJed Brown 
1022e27a552bSJed Brown #undef __FUNCT__
1023e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType"
1024e27a552bSJed Brown /*@C
102561692a83SJed Brown   TSRosWSetType - Set the type of Rosenbrock-W scheme
1026e27a552bSJed Brown 
1027e27a552bSJed Brown   Logically collective
1028e27a552bSJed Brown 
1029e27a552bSJed Brown   Input Parameter:
1030e27a552bSJed Brown +  ts - timestepping context
103161692a83SJed Brown -  rostype - type of Rosenbrock-W scheme
1032e27a552bSJed Brown 
1033020d8f30SJed Brown   Level: beginner
1034e27a552bSJed Brown 
1035020d8f30SJed Brown .seealso: TSRosWGetType(), TSROSW, TSROSW2M, TSROSW2P, TSROSWRA3PW, TSROSWRA34PW2, TSROSWRODAS3, TSROSWSANDU3, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, TSROSWARK3
1036e27a552bSJed Brown @*/
103761692a83SJed Brown PetscErrorCode TSRosWSetType(TS ts,const TSRosWType rostype)
1038e27a552bSJed Brown {
1039e27a552bSJed Brown   PetscErrorCode ierr;
1040e27a552bSJed Brown 
1041e27a552bSJed Brown   PetscFunctionBegin;
1042e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
104361692a83SJed Brown   ierr = PetscTryMethod(ts,"TSRosWSetType_C",(TS,const TSRosWType),(ts,rostype));CHKERRQ(ierr);
1044e27a552bSJed Brown   PetscFunctionReturn(0);
1045e27a552bSJed Brown }
1046e27a552bSJed Brown 
1047e27a552bSJed Brown #undef __FUNCT__
1048e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType"
1049e27a552bSJed Brown /*@C
105061692a83SJed Brown   TSRosWGetType - Get the type of Rosenbrock-W scheme
1051e27a552bSJed Brown 
1052e27a552bSJed Brown   Logically collective
1053e27a552bSJed Brown 
1054e27a552bSJed Brown   Input Parameter:
1055e27a552bSJed Brown .  ts - timestepping context
1056e27a552bSJed Brown 
1057e27a552bSJed Brown   Output Parameter:
105861692a83SJed Brown .  rostype - type of Rosenbrock-W scheme
1059e27a552bSJed Brown 
1060e27a552bSJed Brown   Level: intermediate
1061e27a552bSJed Brown 
1062e27a552bSJed Brown .seealso: TSRosWGetType()
1063e27a552bSJed Brown @*/
106461692a83SJed Brown PetscErrorCode TSRosWGetType(TS ts,const TSRosWType *rostype)
1065e27a552bSJed Brown {
1066e27a552bSJed Brown   PetscErrorCode ierr;
1067e27a552bSJed Brown 
1068e27a552bSJed Brown   PetscFunctionBegin;
1069e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
107061692a83SJed Brown   ierr = PetscUseMethod(ts,"TSRosWGetType_C",(TS,const TSRosWType*),(ts,rostype));CHKERRQ(ierr);
1071e27a552bSJed Brown   PetscFunctionReturn(0);
1072e27a552bSJed Brown }
1073e27a552bSJed Brown 
1074e27a552bSJed Brown #undef __FUNCT__
107561692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian"
1076e27a552bSJed Brown /*@C
107761692a83SJed Brown   TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step.
1078e27a552bSJed Brown 
1079e27a552bSJed Brown   Logically collective
1080e27a552bSJed Brown 
1081e27a552bSJed Brown   Input Parameter:
1082e27a552bSJed Brown +  ts - timestepping context
108361692a83SJed Brown -  flg - PETSC_TRUE to recompute the Jacobian at each stage
1084e27a552bSJed Brown 
1085e27a552bSJed Brown   Level: intermediate
1086e27a552bSJed Brown 
1087e27a552bSJed Brown .seealso: TSRosWGetType()
1088e27a552bSJed Brown @*/
108961692a83SJed Brown PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg)
1090e27a552bSJed Brown {
1091e27a552bSJed Brown   PetscErrorCode ierr;
1092e27a552bSJed Brown 
1093e27a552bSJed Brown   PetscFunctionBegin;
1094e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
109561692a83SJed Brown   ierr = PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
1096e27a552bSJed Brown   PetscFunctionReturn(0);
1097e27a552bSJed Brown }
1098e27a552bSJed Brown 
1099e27a552bSJed Brown EXTERN_C_BEGIN
1100e27a552bSJed Brown #undef __FUNCT__
1101e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType_RosW"
110261692a83SJed Brown PetscErrorCode  TSRosWGetType_RosW(TS ts,const TSRosWType *rostype)
1103e27a552bSJed Brown {
110461692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
1105e27a552bSJed Brown   PetscErrorCode ierr;
1106e27a552bSJed Brown 
1107e27a552bSJed Brown   PetscFunctionBegin;
110861692a83SJed Brown   if (!ros->tableau) {ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);}
110961692a83SJed Brown   *rostype = ros->tableau->name;
1110e27a552bSJed Brown   PetscFunctionReturn(0);
1111e27a552bSJed Brown }
1112e27a552bSJed Brown #undef __FUNCT__
1113e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType_RosW"
111461692a83SJed Brown PetscErrorCode  TSRosWSetType_RosW(TS ts,const TSRosWType rostype)
1115e27a552bSJed Brown {
111661692a83SJed Brown   TS_RosW         *ros = (TS_RosW*)ts->data;
1117e27a552bSJed Brown   PetscErrorCode  ierr;
1118e27a552bSJed Brown   PetscBool       match;
111961692a83SJed Brown   RosWTableauLink link;
1120e27a552bSJed Brown 
1121e27a552bSJed Brown   PetscFunctionBegin;
112261692a83SJed Brown   if (ros->tableau) {
112361692a83SJed Brown     ierr = PetscStrcmp(ros->tableau->name,rostype,&match);CHKERRQ(ierr);
1124e27a552bSJed Brown     if (match) PetscFunctionReturn(0);
1125e27a552bSJed Brown   }
112661692a83SJed Brown   for (link = RosWTableauList; link; link=link->next) {
112761692a83SJed Brown     ierr = PetscStrcmp(link->tab.name,rostype,&match);CHKERRQ(ierr);
1128e27a552bSJed Brown     if (match) {
1129e27a552bSJed Brown       ierr = TSReset_RosW(ts);CHKERRQ(ierr);
113061692a83SJed Brown       ros->tableau = &link->tab;
1131e27a552bSJed Brown       PetscFunctionReturn(0);
1132e27a552bSJed Brown     }
1133e27a552bSJed Brown   }
113461692a83SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype);
1135e27a552bSJed Brown   PetscFunctionReturn(0);
1136e27a552bSJed Brown }
113761692a83SJed Brown 
1138e27a552bSJed Brown #undef __FUNCT__
113961692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian_RosW"
114061692a83SJed Brown PetscErrorCode  TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg)
1141e27a552bSJed Brown {
114261692a83SJed Brown   TS_RosW *ros = (TS_RosW*)ts->data;
1143e27a552bSJed Brown 
1144e27a552bSJed Brown   PetscFunctionBegin;
114561692a83SJed Brown   ros->recompute_jacobian = flg;
1146e27a552bSJed Brown   PetscFunctionReturn(0);
1147e27a552bSJed Brown }
1148e27a552bSJed Brown EXTERN_C_END
1149e27a552bSJed Brown 
1150e27a552bSJed Brown /* ------------------------------------------------------------ */
1151e27a552bSJed Brown /*MC
1152020d8f30SJed Brown       TSROSW - ODE solver using Rosenbrock-W schemes
1153e27a552bSJed Brown 
1154e27a552bSJed Brown   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
1155e27a552bSJed Brown   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
1156e27a552bSJed Brown   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
1157e27a552bSJed Brown 
1158e27a552bSJed Brown   Notes:
115961692a83SJed Brown   This method currently only works with autonomous ODE and DAE.
116061692a83SJed Brown 
116161692a83SJed Brown   Developer notes:
116261692a83SJed Brown   Rosenbrock-W methods are typically specified for autonomous ODE
116361692a83SJed Brown 
116461692a83SJed Brown $  xdot = f(x)
116561692a83SJed Brown 
116661692a83SJed Brown   by the stage equations
116761692a83SJed Brown 
116861692a83SJed Brown $  k_i = h f(x_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j
116961692a83SJed Brown 
117061692a83SJed Brown   and step completion formula
117161692a83SJed Brown 
117261692a83SJed Brown $  x_1 = x_0 + sum_j b_j k_j
117361692a83SJed Brown 
117461692a83SJed Brown   with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(x)
117561692a83SJed Brown   and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner,
117661692a83SJed Brown   we define new variables for the stage equations
117761692a83SJed Brown 
117861692a83SJed Brown $  y_i = gamma_ij k_j
117961692a83SJed Brown 
118061692a83SJed Brown   The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define
118161692a83SJed Brown 
118261692a83SJed Brown $  A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-i}
118361692a83SJed Brown 
118461692a83SJed Brown   to rewrite the method as
118561692a83SJed Brown 
118661692a83SJed Brown $  [M/(h gamma_ii) - J] y_i = f(x_0 + sum_j a_ij y_j) + M sum_j (c_ij/h) y_j
118761692a83SJed Brown $  x_1 = x_0 + sum_j bt_j y_j
118861692a83SJed Brown 
118961692a83SJed Brown    where we have introduced the mass matrix M. Continue by defining
119061692a83SJed Brown 
119161692a83SJed Brown $  ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j
119261692a83SJed Brown 
119361692a83SJed Brown    or, more compactly in tensor notation
119461692a83SJed Brown 
119561692a83SJed Brown $  Ydot = 1/h (Gamma^{-1} \otimes I) Y .
119661692a83SJed Brown 
119761692a83SJed Brown    Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current
119861692a83SJed Brown    stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the
119961692a83SJed Brown    equation
120061692a83SJed Brown 
120161692a83SJed Brown $  g(x_0 + sum_j a_ij y_j + y_i, ydot_i) = 0
120261692a83SJed Brown 
120361692a83SJed Brown    with initial guess y_i = 0.
1204e27a552bSJed Brown 
1205e27a552bSJed Brown   Level: beginner
1206e27a552bSJed Brown 
1207020d8f30SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSRosWSetType(), TSRosWRegister()
1208e27a552bSJed Brown 
1209e27a552bSJed Brown M*/
1210e27a552bSJed Brown EXTERN_C_BEGIN
1211e27a552bSJed Brown #undef __FUNCT__
1212e27a552bSJed Brown #define __FUNCT__ "TSCreate_RosW"
1213e27a552bSJed Brown PetscErrorCode  TSCreate_RosW(TS ts)
1214e27a552bSJed Brown {
121561692a83SJed Brown   TS_RosW        *ros;
1216e27a552bSJed Brown   PetscErrorCode ierr;
1217e27a552bSJed Brown 
1218e27a552bSJed Brown   PetscFunctionBegin;
1219e27a552bSJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
1220e27a552bSJed Brown   ierr = TSRosWInitializePackage(PETSC_NULL);CHKERRQ(ierr);
1221e27a552bSJed Brown #endif
1222e27a552bSJed Brown 
1223e27a552bSJed Brown   ts->ops->reset          = TSReset_RosW;
1224e27a552bSJed Brown   ts->ops->destroy        = TSDestroy_RosW;
1225e27a552bSJed Brown   ts->ops->view           = TSView_RosW;
1226e27a552bSJed Brown   ts->ops->setup          = TSSetUp_RosW;
1227e27a552bSJed Brown   ts->ops->step           = TSStep_RosW;
1228e27a552bSJed Brown   ts->ops->interpolate    = TSInterpolate_RosW;
12291c3436cfSJed Brown   ts->ops->evaluatestep   = TSEvaluateStep_RosW;
1230e27a552bSJed Brown   ts->ops->setfromoptions = TSSetFromOptions_RosW;
1231e27a552bSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_RosW;
1232e27a552bSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_RosW;
1233e27a552bSJed Brown 
123461692a83SJed Brown   ierr = PetscNewLog(ts,TS_RosW,&ros);CHKERRQ(ierr);
123561692a83SJed Brown   ts->data = (void*)ros;
1236e27a552bSJed Brown 
1237e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","TSRosWGetType_RosW",TSRosWGetType_RosW);CHKERRQ(ierr);
1238e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","TSRosWSetType_RosW",TSRosWSetType_RosW);CHKERRQ(ierr);
123961692a83SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","TSRosWSetRecomputeJacobian_RosW",TSRosWSetRecomputeJacobian_RosW);CHKERRQ(ierr);
1240e27a552bSJed Brown   PetscFunctionReturn(0);
1241e27a552bSJed Brown }
1242e27a552bSJed Brown EXTERN_C_END
1243