xref: /petsc/src/ts/impls/rosw/rosw.c (revision b2ce242ef4d0d7b1f624d096cdff676ee62d037d)
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 */
2661692a83SJed Brown   PetscReal *A;                 /* Propagation table, strictly lower triangular */
2761692a83SJed Brown   PetscReal *Gamma;             /* Stage table, lower triangular with nonzero diagonal */
28c17803e7SJed Brown   PetscBool *GammaZeroDiag;     /* Diagonal entries that are zero in stage table Gamma, vector indicating explicit statages */
2943b21953SEmil Constantinescu   PetscReal *GammaExplicitCorr; /* Coefficients for correction terms needed for explicit stages in transformed variables*/
3061692a83SJed Brown   PetscReal *b;                 /* Step completion table */
31fe7e6d57SJed Brown   PetscReal *bembed;            /* Step completion table for embedded method of order one less */
3261692a83SJed Brown   PetscReal *ASum;              /* Row sum of A */
3361692a83SJed Brown   PetscReal *GammaSum;          /* Row sum of Gamma, only needed for non-autonomous systems */
3461692a83SJed Brown   PetscReal *At;                /* Propagation table in transformed variables */
3561692a83SJed Brown   PetscReal *bt;                /* Step completion table in transformed variables */
36fe7e6d57SJed Brown   PetscReal *bembedt;           /* Step completion table of order one less in transformed variables */
3761692a83SJed Brown   PetscReal *GammaInv;          /* Inverse of Gamma, used for transformed variables */
388d59e960SJed Brown   PetscReal ccfl;               /* Placeholder for CFL coefficient relative to forward Euler */
39e27a552bSJed Brown };
4061692a83SJed Brown typedef struct _RosWTableauLink *RosWTableauLink;
4161692a83SJed Brown struct _RosWTableauLink {
4261692a83SJed Brown   struct _RosWTableau tab;
4361692a83SJed Brown   RosWTableauLink next;
44e27a552bSJed Brown };
4561692a83SJed Brown static RosWTableauLink RosWTableauList;
46e27a552bSJed Brown 
47e27a552bSJed Brown typedef struct {
4861692a83SJed Brown   RosWTableau tableau;
4961692a83SJed Brown   Vec         *Y;               /* States computed during the step, used to complete the step */
50e27a552bSJed Brown   Vec         Ydot;             /* Work vector holding Ydot during residual evaluation */
5161692a83SJed Brown   Vec         Ystage;           /* Work vector for the state value at each stage */
5261692a83SJed Brown   Vec         Zdot;             /* Ydot = Zdot + shift*Y */
5361692a83SJed Brown   Vec         Zstage;           /* Y = Zstage + Y */
541c3436cfSJed Brown   PetscScalar *work;            /* Scalar work space of length number of stages, used to prepare VecMAXPY() */
55e27a552bSJed Brown   PetscReal   shift;
56e27a552bSJed Brown   PetscReal   stage_time;
57c17803e7SJed Brown   PetscReal   stage_explicit;     /* Flag indicates that the current stage is explicit */
5861692a83SJed Brown   PetscBool   recompute_jacobian; /* Recompute the Jacobian at each stage, default is to freeze the Jacobian at the start of each step */
59108c343cSJed Brown   TSStepStatus status;
60e27a552bSJed Brown } TS_RosW;
61e27a552bSJed Brown 
62fe7e6d57SJed Brown /*MC
633606a31eSEmil Constantinescu      TSROSWTHETA1 - One stage first order L-stable Rosenbrock-W scheme (aka theta method).
643606a31eSEmil Constantinescu 
653606a31eSEmil Constantinescu      Only an approximate Jacobian is needed.
663606a31eSEmil Constantinescu 
673606a31eSEmil Constantinescu      Level: intermediate
683606a31eSEmil Constantinescu 
693606a31eSEmil Constantinescu .seealso: TSROSW
703606a31eSEmil Constantinescu M*/
713606a31eSEmil Constantinescu 
723606a31eSEmil Constantinescu /*MC
733606a31eSEmil Constantinescu      TSROSWTHETA2 - One stage second order A-stable Rosenbrock-W scheme (aka theta method).
743606a31eSEmil Constantinescu 
753606a31eSEmil Constantinescu      Only an approximate Jacobian is needed.
763606a31eSEmil Constantinescu 
773606a31eSEmil Constantinescu      Level: intermediate
783606a31eSEmil Constantinescu 
793606a31eSEmil Constantinescu .seealso: TSROSW
803606a31eSEmil Constantinescu M*/
813606a31eSEmil Constantinescu 
823606a31eSEmil Constantinescu /*MC
83fe7e6d57SJed Brown      TSROSW2M - Two stage second order L-stable Rosenbrock-W scheme.
84fe7e6d57SJed Brown 
85fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2P.
86fe7e6d57SJed Brown 
87fe7e6d57SJed Brown      Level: intermediate
88fe7e6d57SJed Brown 
89fe7e6d57SJed Brown .seealso: TSROSW
90fe7e6d57SJed Brown M*/
91fe7e6d57SJed Brown 
92fe7e6d57SJed Brown /*MC
93fe7e6d57SJed Brown      TSROSW2P - Two stage second order L-stable Rosenbrock-W scheme.
94fe7e6d57SJed Brown 
95fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2M.
96fe7e6d57SJed Brown 
97fe7e6d57SJed Brown      Level: intermediate
98fe7e6d57SJed Brown 
99fe7e6d57SJed Brown .seealso: TSROSW
100fe7e6d57SJed Brown M*/
101fe7e6d57SJed Brown 
102fe7e6d57SJed Brown /*MC
103fe7e6d57SJed Brown      TSROSWRA3PW - Three stage third order Rosenbrock-W scheme for PDAE of index 1.
104fe7e6d57SJed Brown 
105fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.
106fe7e6d57SJed Brown 
107fe7e6d57SJed 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.
108fe7e6d57SJed Brown 
109fe7e6d57SJed Brown      References:
110fe7e6d57SJed Brown      Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005.
111fe7e6d57SJed Brown 
112fe7e6d57SJed Brown      Level: intermediate
113fe7e6d57SJed Brown 
114fe7e6d57SJed Brown .seealso: TSROSW
115fe7e6d57SJed Brown M*/
116fe7e6d57SJed Brown 
117fe7e6d57SJed Brown /*MC
118fe7e6d57SJed Brown      TSROSWRA34PW2 - Four stage third order L-stable Rosenbrock-W scheme for PDAE of index 1.
119fe7e6d57SJed Brown 
120fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.
121fe7e6d57SJed Brown 
122fe7e6d57SJed 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.
123fe7e6d57SJed Brown 
124fe7e6d57SJed Brown      References:
125fe7e6d57SJed Brown      Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005.
126fe7e6d57SJed Brown 
127fe7e6d57SJed Brown      Level: intermediate
128fe7e6d57SJed Brown 
129fe7e6d57SJed Brown .seealso: TSROSW
130fe7e6d57SJed Brown M*/
131fe7e6d57SJed Brown 
132ef3c5b88SJed Brown /*MC
133ef3c5b88SJed Brown      TSROSWRODAS3 - Four stage third order L-stable Rosenbrock scheme
134ef3c5b88SJed Brown 
135ef3c5b88SJed Brown      By default, the Jacobian is only recomputed once per step.
136ef3c5b88SJed Brown 
137ef3c5b88SJed Brown      Both the third order and embedded second order methods are stiffly accurate and L-stable.
138ef3c5b88SJed Brown 
139ef3c5b88SJed Brown      References:
140ef3c5b88SJed Brown      Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997.
141ef3c5b88SJed Brown 
142ef3c5b88SJed Brown      Level: intermediate
143ef3c5b88SJed Brown 
144ef3c5b88SJed Brown .seealso: TSROSW, TSROSWSANDU3
145ef3c5b88SJed Brown M*/
146ef3c5b88SJed Brown 
147ef3c5b88SJed Brown /*MC
148ef3c5b88SJed Brown      TSROSWSANDU3 - Three stage third order L-stable Rosenbrock scheme
149ef3c5b88SJed Brown 
150ef3c5b88SJed Brown      By default, the Jacobian is only recomputed once per step.
151ef3c5b88SJed Brown 
152ef3c5b88SJed Brown      The third order method is L-stable, but not stiffly accurate.
153ef3c5b88SJed Brown      The second order embedded method is strongly A-stable with R(infty) = 0.5.
154ef3c5b88SJed Brown      The internal stages are L-stable.
155ef3c5b88SJed Brown      This method is called ROS3 in the paper.
156ef3c5b88SJed Brown 
157ef3c5b88SJed Brown      References:
158ef3c5b88SJed Brown      Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997.
159ef3c5b88SJed Brown 
160ef3c5b88SJed Brown      Level: intermediate
161ef3c5b88SJed Brown 
162ef3c5b88SJed Brown .seealso: TSROSW, TSROSWRODAS3
163ef3c5b88SJed Brown M*/
164ef3c5b88SJed Brown 
165961f28d0SJed Brown /*MC
166961f28d0SJed Brown      TSROSWASSP3P3S1C - A-stable Rosenbrock-W method with SSP explicit part, third order, three stages
167961f28d0SJed Brown 
168961f28d0SJed Brown      By default, the Jacobian is only recomputed once per step.
169961f28d0SJed Brown 
170961f28d0SJed Brown      A-stable SPP explicit order 3, 3 stages, CFL 1 (eff = 1/3)
171961f28d0SJed Brown 
172961f28d0SJed Brown      References:
173961f28d0SJed Brown      Emil Constantinescu
174961f28d0SJed Brown 
175961f28d0SJed Brown      Level: intermediate
176961f28d0SJed Brown 
17743b21953SEmil Constantinescu .seealso: TSROSW, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, SSP
178961f28d0SJed Brown M*/
179961f28d0SJed Brown 
180961f28d0SJed Brown /*MC
181961f28d0SJed Brown      TSROSWLASSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, three stages
182961f28d0SJed Brown 
183961f28d0SJed Brown      By default, the Jacobian is only recomputed once per step.
184961f28d0SJed Brown 
185961f28d0SJed Brown      L-stable (A-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2)
186961f28d0SJed Brown 
187961f28d0SJed Brown      References:
188961f28d0SJed Brown      Emil Constantinescu
189961f28d0SJed Brown 
190961f28d0SJed Brown      Level: intermediate
191961f28d0SJed Brown 
19243b21953SEmil Constantinescu .seealso: TSROSW, TSROSWASSP3P3S1C, TSROSWLLSSP3P4S2C, TSSSP
193961f28d0SJed Brown M*/
194961f28d0SJed Brown 
195961f28d0SJed Brown /*MC
19643b21953SEmil Constantinescu      TSROSWLLSSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, three stages
197961f28d0SJed Brown 
198961f28d0SJed Brown      By default, the Jacobian is only recomputed once per step.
199961f28d0SJed Brown 
200961f28d0SJed Brown      L-stable (L-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2)
201961f28d0SJed Brown 
202961f28d0SJed Brown      References:
203961f28d0SJed Brown      Emil Constantinescu
204961f28d0SJed Brown 
205961f28d0SJed Brown      Level: intermediate
206961f28d0SJed Brown 
207961f28d0SJed Brown .seealso: TSROSW, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSSSP
208961f28d0SJed Brown M*/
209961f28d0SJed Brown 
210e27a552bSJed Brown #undef __FUNCT__
211e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterAll"
212e27a552bSJed Brown /*@C
213e27a552bSJed Brown   TSRosWRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSRosW
214e27a552bSJed Brown 
215e27a552bSJed Brown   Not Collective, but should be called by all processes which will need the schemes to be registered
216e27a552bSJed Brown 
217e27a552bSJed Brown   Level: advanced
218e27a552bSJed Brown 
219e27a552bSJed Brown .keywords: TS, TSRosW, register, all
220e27a552bSJed Brown 
221e27a552bSJed Brown .seealso:  TSRosWRegisterDestroy()
222e27a552bSJed Brown @*/
223e27a552bSJed Brown PetscErrorCode TSRosWRegisterAll(void)
224e27a552bSJed Brown {
225e27a552bSJed Brown   PetscErrorCode ierr;
226e27a552bSJed Brown 
227e27a552bSJed Brown   PetscFunctionBegin;
228e27a552bSJed Brown   if (TSRosWRegisterAllCalled) PetscFunctionReturn(0);
229e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_TRUE;
230e27a552bSJed Brown 
231e27a552bSJed Brown   {
2323606a31eSEmil Constantinescu     const PetscReal
2333606a31eSEmil Constantinescu       A = 0,
2343606a31eSEmil Constantinescu       Gamma = 1,
2353606a31eSEmil Constantinescu       b = 1;
2363606a31eSEmil Constantinescu     ierr = TSRosWRegister(TSROSWTHETA1,1,1,&A,&Gamma,&b,PETSC_NULL);CHKERRQ(ierr);
2373606a31eSEmil Constantinescu   }
2383606a31eSEmil Constantinescu 
2393606a31eSEmil Constantinescu   {
2403606a31eSEmil Constantinescu     const PetscReal
2413606a31eSEmil Constantinescu       A= 0,
2423606a31eSEmil Constantinescu       Gamma = 0.5,
2433606a31eSEmil Constantinescu       b = 1;
2443606a31eSEmil Constantinescu     ierr = TSRosWRegister(TSROSWTHETA2,2,1,&A,&Gamma,&b,PETSC_NULL);CHKERRQ(ierr);
2453606a31eSEmil Constantinescu   }
2463606a31eSEmil Constantinescu 
2473606a31eSEmil Constantinescu   {
24861692a83SJed Brown     const PetscReal g = 1. + 1./PetscSqrtReal(2.0);
249e27a552bSJed Brown     const PetscReal
25061692a83SJed Brown       A[2][2] = {{0,0}, {1.,0}},
25161692a83SJed Brown       Gamma[2][2] = {{g,0}, {-2.*g,g}},
2521c3436cfSJed Brown       b[2] = {0.5,0.5},
2531c3436cfSJed Brown       b1[2] = {1.0,0.0};
2541c3436cfSJed Brown     ierr = TSRosWRegister(TSROSW2P,2,2,&A[0][0],&Gamma[0][0],b,b1);CHKERRQ(ierr);
255e27a552bSJed Brown   }
256e27a552bSJed Brown   {
25761692a83SJed Brown     const PetscReal g = 1. - 1./PetscSqrtReal(2.0);
258e27a552bSJed Brown     const PetscReal
25961692a83SJed Brown       A[2][2] = {{0,0}, {1.,0}},
26061692a83SJed Brown       Gamma[2][2] = {{g,0}, {-2.*g,g}},
2611c3436cfSJed Brown       b[2] = {0.5,0.5},
2621c3436cfSJed Brown       b1[2] = {1.0,0.0};
2631c3436cfSJed Brown     ierr = TSRosWRegister(TSROSW2M,2,2,&A[0][0],&Gamma[0][0],b,b1);CHKERRQ(ierr);
264fe7e6d57SJed Brown   }
265fe7e6d57SJed Brown   {
266fe7e6d57SJed Brown     const PetscReal g = 7.8867513459481287e-01;
267fe7e6d57SJed Brown     const PetscReal
268fe7e6d57SJed Brown       A[3][3] = {{0,0,0},
269fe7e6d57SJed Brown                  {1.5773502691896257e+00,0,0},
270fe7e6d57SJed Brown                  {0.5,0,0}},
271fe7e6d57SJed Brown       Gamma[3][3] = {{g,0,0},
272fe7e6d57SJed Brown                      {-1.5773502691896257e+00,g,0},
273fe7e6d57SJed Brown                      {-6.7075317547305480e-01,1.7075317547305482e-01,g}},
274fe7e6d57SJed Brown       b[3] = {1.0566243270259355e-01,4.9038105676657971e-02,8.4529946162074843e-01},
275fe7e6d57SJed Brown       b2[3] = {-1.7863279495408180e-01,1./3.,8.4529946162074843e-01};
276fe7e6d57SJed Brown     ierr = TSRosWRegister(TSROSWRA3PW,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
277fe7e6d57SJed Brown   }
278fe7e6d57SJed Brown   {
279fe7e6d57SJed Brown     const PetscReal g = 4.3586652150845900e-01;
280fe7e6d57SJed Brown     const PetscReal
281fe7e6d57SJed Brown       A[4][4] = {{0,0,0,0},
282fe7e6d57SJed Brown                  {8.7173304301691801e-01,0,0,0},
283fe7e6d57SJed Brown                  {8.4457060015369423e-01,-1.1299064236484185e-01,0,0},
284fe7e6d57SJed Brown                  {0,0,1.,0}},
285fe7e6d57SJed Brown       Gamma[4][4] = {{g,0,0,0},
286fe7e6d57SJed Brown                      {-8.7173304301691801e-01,g,0,0},
287fe7e6d57SJed Brown                      {-9.0338057013044082e-01,5.4180672388095326e-02,g,0},
288fe7e6d57SJed Brown                      {2.4212380706095346e-01,-1.2232505839045147e+00,5.4526025533510214e-01,g}},
289fe7e6d57SJed Brown       b[4] = {2.4212380706095346e-01,-1.2232505839045147e+00,1.5452602553351020e+00,4.3586652150845900e-01},
290fe7e6d57SJed Brown       b2[4] = {3.7810903145819369e-01,-9.6042292212423178e-02,5.0000000000000000e-01,2.1793326075422950e-01};
291fe7e6d57SJed Brown     ierr = TSRosWRegister(TSROSWRA34PW2,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
292e27a552bSJed Brown   }
293ef3c5b88SJed Brown   {
294ef3c5b88SJed Brown     const PetscReal g = 0.5;
295ef3c5b88SJed Brown     const PetscReal
296ef3c5b88SJed Brown       A[4][4] = {{0,0,0,0},
297ef3c5b88SJed Brown                  {0,0,0,0},
298ef3c5b88SJed Brown                  {1.,0,0,0},
299ef3c5b88SJed Brown                  {0.75,-0.25,0.5,0}},
300ef3c5b88SJed Brown       Gamma[4][4] = {{g,0,0,0},
301ef3c5b88SJed Brown                      {1.,g,0,0},
302ef3c5b88SJed Brown                      {-0.25,-0.25,g,0},
303ef3c5b88SJed Brown                      {1./12,1./12,-2./3,g}},
304ef3c5b88SJed Brown       b[4] = {5./6,-1./6,-1./6,0.5},
305ef3c5b88SJed Brown       b2[4] = {0.75,-0.25,0.5,0};
306ef3c5b88SJed Brown     ierr = TSRosWRegister(TSROSWRODAS3,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
307ef3c5b88SJed Brown   }
308ef3c5b88SJed Brown   {
309ef3c5b88SJed Brown     const PetscReal g = 0.43586652150845899941601945119356;
310ef3c5b88SJed Brown     const PetscReal
311ef3c5b88SJed Brown       A[3][3] = {{0,0,0},
312ef3c5b88SJed Brown                  {g,0,0},
313ef3c5b88SJed Brown                  {g,0,0}},
314ef3c5b88SJed Brown       Gamma[3][3] = {{g,0,0},
315ef3c5b88SJed Brown                      {-0.19294655696029095575009695436041,g,0},
316ef3c5b88SJed Brown                      {0,1.74927148125794685173529749738960,g}},
317ef3c5b88SJed Brown       b[3] = {-0.75457412385404315829818998646589,1.94100407061964420292840123379419,-0.18642994676560104463021124732829},
318ef3c5b88SJed Brown       b2[3] = {-1.53358745784149585370766523913002,2.81745131148625772213931745457622,-0.28386385364476186843165221544619};
319ef3c5b88SJed Brown     ierr = TSRosWRegister(TSROSWSANDU3,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
320ef3c5b88SJed Brown   }
321b1c69cc3SEmil Constantinescu   {
322*aaf9cf16SJed Brown     const PetscReal s3 = PetscSqrtReal(3.),g = (3.0+s3)/6.0;
323b1c69cc3SEmil Constantinescu     const PetscReal
324b1c69cc3SEmil Constantinescu       A[3][3] = {{0,0,0},
325b1c69cc3SEmil Constantinescu                  {1,0,0},
326b1c69cc3SEmil Constantinescu                  {0.25,0.25,0}},
327b1c69cc3SEmil Constantinescu       Gamma[3][3] = {{0,0,0},
328*aaf9cf16SJed Brown                      {(-3.0-s3)/6.0,g,0},
329*aaf9cf16SJed Brown                      {(-3.0-s3)/24.0,(-3.0-s3)/8.0,g}},
330b1c69cc3SEmil Constantinescu         b[3] = {1./6.,1./6.,2./3.},
331b1c69cc3SEmil Constantinescu           b2[3] = {1./4.,1./4.,1./2.};
332b1c69cc3SEmil Constantinescu     ierr = TSRosWRegister(TSROSWASSP3P3S1C,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
333b1c69cc3SEmil Constantinescu   }
334b1c69cc3SEmil Constantinescu 
335b1c69cc3SEmil Constantinescu   {
336b1c69cc3SEmil Constantinescu     const PetscReal
337b1c69cc3SEmil Constantinescu       A[4][4] = {{0,0,0,0},
338b1c69cc3SEmil Constantinescu                  {1./2.,0,0,0},
339b1c69cc3SEmil Constantinescu                  {1./2.,1./2.,0,0},
340b1c69cc3SEmil Constantinescu                  {1./6.,1./6.,1./6.,0}},
341b1c69cc3SEmil Constantinescu       Gamma[4][4] = {{1./2.,0,0,0},
342b1c69cc3SEmil Constantinescu                      {0.0,1./4.,0,0},
343b1c69cc3SEmil Constantinescu                      {-2.,-2./3.,2./3.,0},
344b1c69cc3SEmil Constantinescu                      {1./2.,5./36.,-2./9,0}},
345b1c69cc3SEmil Constantinescu         b[4] = {1./6.,1./6.,1./6.,1./2.},
346b1c69cc3SEmil Constantinescu         b2[4] = {1./8.,3./4.,1./8.,0};
347b1c69cc3SEmil Constantinescu      ierr = TSRosWRegister(TSROSWLASSP3P4S2C,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
348b1c69cc3SEmil Constantinescu   }
349b1c69cc3SEmil Constantinescu 
350b1c69cc3SEmil Constantinescu   {
351b1c69cc3SEmil Constantinescu     const PetscReal
352b1c69cc3SEmil Constantinescu       A[4][4] = {{0,0,0,0},
353b1c69cc3SEmil Constantinescu                  {1./2.,0,0,0},
354b1c69cc3SEmil Constantinescu                  {1./2.,1./2.,0,0},
355b1c69cc3SEmil Constantinescu                  {1./6.,1./6.,1./6.,0}},
356b1c69cc3SEmil Constantinescu       Gamma[4][4] = {{1./2.,0,0,0},
357b1c69cc3SEmil Constantinescu                      {0.0,3./4.,0,0},
358b1c69cc3SEmil Constantinescu                      {-2./3.,-23./9.,2./9.,0},
359b1c69cc3SEmil Constantinescu                      {1./18.,65./108.,-2./27,0}},
360b1c69cc3SEmil Constantinescu         b[4] = {1./6.,1./6.,1./6.,1./2.},
361b1c69cc3SEmil Constantinescu         b2[4] = {3./16.,10./16.,3./16.,0};
36243b21953SEmil Constantinescu      ierr = TSRosWRegister(TSROSWLLSSP3P4S2C,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
363b1c69cc3SEmil Constantinescu   }
364753f8adbSEmil Constantinescu 
365753f8adbSEmil Constantinescu  {
366753f8adbSEmil Constantinescu    PetscReal A[4][4],Gamma[4][4],b[4],b2[4];
367753f8adbSEmil Constantinescu 
368753f8adbSEmil Constantinescu    Gamma[0][0]=0.4358665215084589994160194475295062513822671686978816;
36905e8e825SJed Brown    Gamma[0][1]=0; Gamma[0][2]=0; Gamma[0][3]=0;
370753f8adbSEmil Constantinescu    Gamma[1][0]=-1.997527830934941248426324674704153457289527280554476;
371753f8adbSEmil Constantinescu    Gamma[1][1]=0.4358665215084589994160194475295062513822671686978816;
37205e8e825SJed Brown    Gamma[1][2]=0; Gamma[1][3]=0;
373753f8adbSEmil Constantinescu    Gamma[2][0]=-1.007948511795029620852002345345404191008352770119903;
374753f8adbSEmil Constantinescu    Gamma[2][1]=-0.004648958462629345562774289390054679806993396798458131;
375753f8adbSEmil Constantinescu    Gamma[2][2]=0.4358665215084589994160194475295062513822671686978816;
37605e8e825SJed Brown    Gamma[2][3]=0;
377753f8adbSEmil Constantinescu    Gamma[3][0]=-0.6685429734233467180451604600279552604364311322650783;
378753f8adbSEmil Constantinescu    Gamma[3][1]=0.6056625986449338476089525334450053439525178740492984;
379753f8adbSEmil Constantinescu    Gamma[3][2]=-0.9717899277217721234705114616271378792182450260943198;
380753f8adbSEmil Constantinescu    Gamma[3][3]=0;
381753f8adbSEmil Constantinescu 
38205e8e825SJed Brown    A[0][0]=0; A[0][1]=0; A[0][2]=0; A[0][3]=0;
383753f8adbSEmil Constantinescu    A[1][0]=0.8717330430169179988320388950590125027645343373957631;
38405e8e825SJed Brown    A[1][1]=0; A[1][2]=0; A[1][3]=0;
385753f8adbSEmil Constantinescu    A[2][0]=0.5275890119763004115618079766722914408876108660811028;
386753f8adbSEmil Constantinescu    A[2][1]=0.07241098802369958843819203208518599088698057726988732;
38705e8e825SJed Brown    A[2][2]=0; A[2][3]=0;
388753f8adbSEmil Constantinescu    A[3][0]=0.3990960076760701320627260685975778145384666450351314;
389753f8adbSEmil Constantinescu    A[3][1]=-0.4375576546135194437228463747348862825846903771419953;
390753f8adbSEmil Constantinescu    A[3][2]=1.038461646937449311660120300601880176655352737312713;
39105e8e825SJed Brown    A[3][3]=0;
392753f8adbSEmil Constantinescu 
393753f8adbSEmil Constantinescu    b[0]=0.1876410243467238251612921333138006734899663569186926;
394753f8adbSEmil Constantinescu    b[1]=-0.5952974735769549480478230473706443582188442040780541;
395753f8adbSEmil Constantinescu    b[2]=0.9717899277217721234705114616271378792182450260943198;
396753f8adbSEmil Constantinescu    b[3]=0.4358665215084589994160194475295062513822671686978816;
397753f8adbSEmil Constantinescu 
398753f8adbSEmil Constantinescu    b2[0]=0.2147402862233891404862383521089097657790734483804460;
399753f8adbSEmil Constantinescu    b2[1]=-0.4851622638849390928209050538171743017757490232519684;
400753f8adbSEmil Constantinescu    b2[2]=0.8687250025203875511662123688667549217531982787600080;
401753f8adbSEmil Constantinescu    b2[3]=0.4016969751411624011684543450940068201770721128357014;
402753f8adbSEmil Constantinescu 
403753f8adbSEmil Constantinescu    ierr = TSRosWRegister(TSROSWARK3,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
404753f8adbSEmil Constantinescu   }
405753f8adbSEmil Constantinescu 
406e27a552bSJed Brown   PetscFunctionReturn(0);
407e27a552bSJed Brown }
408e27a552bSJed Brown 
409e27a552bSJed Brown #undef __FUNCT__
410e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterDestroy"
411e27a552bSJed Brown /*@C
412e27a552bSJed Brown    TSRosWRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister().
413e27a552bSJed Brown 
414e27a552bSJed Brown    Not Collective
415e27a552bSJed Brown 
416e27a552bSJed Brown    Level: advanced
417e27a552bSJed Brown 
418e27a552bSJed Brown .keywords: TSRosW, register, destroy
419e27a552bSJed Brown .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegisterDynamic()
420e27a552bSJed Brown @*/
421e27a552bSJed Brown PetscErrorCode TSRosWRegisterDestroy(void)
422e27a552bSJed Brown {
423e27a552bSJed Brown   PetscErrorCode ierr;
42461692a83SJed Brown   RosWTableauLink link;
425e27a552bSJed Brown 
426e27a552bSJed Brown   PetscFunctionBegin;
42761692a83SJed Brown   while ((link = RosWTableauList)) {
42861692a83SJed Brown     RosWTableau t = &link->tab;
42961692a83SJed Brown     RosWTableauList = link->next;
43061692a83SJed Brown     ierr = PetscFree5(t->A,t->Gamma,t->b,t->ASum,t->GammaSum);CHKERRQ(ierr);
43143b21953SEmil Constantinescu     ierr = PetscFree5(t->At,t->bt,t->GammaInv,t->GammaZeroDiag,t->GammaExplicitCorr);CHKERRQ(ierr);
432fe7e6d57SJed Brown     ierr = PetscFree2(t->bembed,t->bembedt);CHKERRQ(ierr);
433e27a552bSJed Brown     ierr = PetscFree(t->name);CHKERRQ(ierr);
434e27a552bSJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
435e27a552bSJed Brown   }
436e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_FALSE;
437e27a552bSJed Brown   PetscFunctionReturn(0);
438e27a552bSJed Brown }
439e27a552bSJed Brown 
440e27a552bSJed Brown #undef __FUNCT__
441e27a552bSJed Brown #define __FUNCT__ "TSRosWInitializePackage"
442e27a552bSJed Brown /*@C
443e27a552bSJed Brown   TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called
444e27a552bSJed Brown   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RosW()
445e27a552bSJed Brown   when using static libraries.
446e27a552bSJed Brown 
447e27a552bSJed Brown   Input Parameter:
448e27a552bSJed Brown   path - The dynamic library path, or PETSC_NULL
449e27a552bSJed Brown 
450e27a552bSJed Brown   Level: developer
451e27a552bSJed Brown 
452e27a552bSJed Brown .keywords: TS, TSRosW, initialize, package
453e27a552bSJed Brown .seealso: PetscInitialize()
454e27a552bSJed Brown @*/
455e27a552bSJed Brown PetscErrorCode TSRosWInitializePackage(const char path[])
456e27a552bSJed Brown {
457e27a552bSJed Brown   PetscErrorCode ierr;
458e27a552bSJed Brown 
459e27a552bSJed Brown   PetscFunctionBegin;
460e27a552bSJed Brown   if (TSRosWPackageInitialized) PetscFunctionReturn(0);
461e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_TRUE;
462e27a552bSJed Brown   ierr = TSRosWRegisterAll();CHKERRQ(ierr);
463e27a552bSJed Brown   ierr = PetscRegisterFinalize(TSRosWFinalizePackage);CHKERRQ(ierr);
464e27a552bSJed Brown   PetscFunctionReturn(0);
465e27a552bSJed Brown }
466e27a552bSJed Brown 
467e27a552bSJed Brown #undef __FUNCT__
468e27a552bSJed Brown #define __FUNCT__ "TSRosWFinalizePackage"
469e27a552bSJed Brown /*@C
470e27a552bSJed Brown   TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is
471e27a552bSJed Brown   called from PetscFinalize().
472e27a552bSJed Brown 
473e27a552bSJed Brown   Level: developer
474e27a552bSJed Brown 
475e27a552bSJed Brown .keywords: Petsc, destroy, package
476e27a552bSJed Brown .seealso: PetscFinalize()
477e27a552bSJed Brown @*/
478e27a552bSJed Brown PetscErrorCode TSRosWFinalizePackage(void)
479e27a552bSJed Brown {
480e27a552bSJed Brown   PetscErrorCode ierr;
481e27a552bSJed Brown 
482e27a552bSJed Brown   PetscFunctionBegin;
483e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_FALSE;
484e27a552bSJed Brown   ierr = TSRosWRegisterDestroy();CHKERRQ(ierr);
485e27a552bSJed Brown   PetscFunctionReturn(0);
486e27a552bSJed Brown }
487e27a552bSJed Brown 
488e27a552bSJed Brown #undef __FUNCT__
489e27a552bSJed Brown #define __FUNCT__ "TSRosWRegister"
490e27a552bSJed Brown /*@C
49161692a83SJed Brown    TSRosWRegister - register a Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
492e27a552bSJed Brown 
493e27a552bSJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
494e27a552bSJed Brown 
495e27a552bSJed Brown    Input Parameters:
496e27a552bSJed Brown +  name - identifier for method
497e27a552bSJed Brown .  order - approximation order of method
498e27a552bSJed Brown .  s - number of stages, this is the dimension of the matrices below
49961692a83SJed Brown .  A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular
50061692a83SJed Brown .  Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal
501fe7e6d57SJed Brown .  b - Step completion table (dimension s)
502fe7e6d57SJed Brown -  bembed - Step completion table for a scheme of order one less (dimension s, PETSC_NULL if no embedded scheme is available)
503e27a552bSJed Brown 
504e27a552bSJed Brown    Notes:
50561692a83SJed Brown    Several Rosenbrock W methods are provided, this function is only needed to create new methods.
506e27a552bSJed Brown 
507e27a552bSJed Brown    Level: advanced
508e27a552bSJed Brown 
509e27a552bSJed Brown .keywords: TS, register
510e27a552bSJed Brown 
511e27a552bSJed Brown .seealso: TSRosW
512e27a552bSJed Brown @*/
513e27a552bSJed Brown PetscErrorCode TSRosWRegister(const TSRosWType name,PetscInt order,PetscInt s,
514fe7e6d57SJed Brown                               const PetscReal A[],const PetscReal Gamma[],const PetscReal b[],const PetscReal bembed[])
515e27a552bSJed Brown {
516e27a552bSJed Brown   PetscErrorCode ierr;
51761692a83SJed Brown   RosWTableauLink link;
51861692a83SJed Brown   RosWTableau t;
51961692a83SJed Brown   PetscInt i,j,k;
52061692a83SJed Brown   PetscScalar *GammaInv;
521e27a552bSJed Brown 
522e27a552bSJed Brown   PetscFunctionBegin;
523fe7e6d57SJed Brown   PetscValidCharPointer(name,1);
524fe7e6d57SJed Brown   PetscValidPointer(A,4);
525fe7e6d57SJed Brown   PetscValidPointer(Gamma,5);
526fe7e6d57SJed Brown   PetscValidPointer(b,6);
527fe7e6d57SJed Brown   if (bembed) PetscValidPointer(bembed,7);
528fe7e6d57SJed Brown 
529e27a552bSJed Brown   ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
530e27a552bSJed Brown   ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
531e27a552bSJed Brown   t = &link->tab;
532e27a552bSJed Brown   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
533e27a552bSJed Brown   t->order = order;
534e27a552bSJed Brown   t->s = s;
53561692a83SJed 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);
53643b21953SEmil 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);
537e27a552bSJed Brown   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
53861692a83SJed Brown   ierr = PetscMemcpy(t->Gamma,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr);
53943b21953SEmil Constantinescu   ierr = PetscMemcpy(t->GammaExplicitCorr,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr);
54061692a83SJed Brown   ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
541fe7e6d57SJed Brown   if (bembed) {
542fe7e6d57SJed Brown     ierr = PetscMalloc2(s,PetscReal,&t->bembed,s,PetscReal,&t->bembedt);CHKERRQ(ierr);
543fe7e6d57SJed Brown     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
544fe7e6d57SJed Brown   }
54561692a83SJed Brown   for (i=0; i<s; i++) {
54661692a83SJed Brown     t->ASum[i] = 0;
54761692a83SJed Brown     t->GammaSum[i] = 0;
54861692a83SJed Brown     for (j=0; j<s; j++) {
54961692a83SJed Brown       t->ASum[i] += A[i*s+j];
550fe7e6d57SJed Brown       t->GammaSum[i] += Gamma[i*s+j];
55161692a83SJed Brown     }
55261692a83SJed Brown   }
55361692a83SJed Brown   ierr = PetscMalloc(s*s*sizeof(PetscScalar),&GammaInv);CHKERRQ(ierr); /* Need to use Scalar for inverse, then convert back to Real */
55461692a83SJed Brown   for (i=0; i<s*s; i++) GammaInv[i] = Gamma[i];
555fd96d5b0SEmil Constantinescu   for (i=0; i<s; i++) {
556fd96d5b0SEmil Constantinescu     if (Gamma[i*s+i] == 0.0) {
557fd96d5b0SEmil Constantinescu       GammaInv[i*s+i] = 1.0;
558c17803e7SJed Brown       t->GammaZeroDiag[i] = PETSC_TRUE;
559fd96d5b0SEmil Constantinescu     } else {
560c17803e7SJed Brown       t->GammaZeroDiag[i] = PETSC_FALSE;
561fd96d5b0SEmil Constantinescu     }
562fd96d5b0SEmil Constantinescu   }
563fd96d5b0SEmil Constantinescu 
56461692a83SJed Brown   switch (s) {
56561692a83SJed Brown   case 1: GammaInv[0] = 1./GammaInv[0]; break;
56661692a83SJed Brown   case 2: ierr = Kernel_A_gets_inverse_A_2(GammaInv,0);CHKERRQ(ierr); break;
56761692a83SJed Brown   case 3: ierr = Kernel_A_gets_inverse_A_3(GammaInv,0);CHKERRQ(ierr); break;
56861692a83SJed Brown   case 4: ierr = Kernel_A_gets_inverse_A_4(GammaInv,0);CHKERRQ(ierr); break;
56961692a83SJed Brown   case 5: {
57061692a83SJed Brown     PetscInt ipvt5[5];
57161692a83SJed Brown     MatScalar work5[5*5];
57261692a83SJed Brown     ierr = Kernel_A_gets_inverse_A_5(GammaInv,ipvt5,work5,0);CHKERRQ(ierr); break;
57361692a83SJed Brown   }
57461692a83SJed Brown   case 6: ierr = Kernel_A_gets_inverse_A_6(GammaInv,0);CHKERRQ(ierr); break;
57561692a83SJed Brown   case 7: ierr = Kernel_A_gets_inverse_A_7(GammaInv,0);CHKERRQ(ierr); break;
57661692a83SJed Brown   default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for %D stages",s);
57761692a83SJed Brown   }
57861692a83SJed Brown   for (i=0; i<s*s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]);
57961692a83SJed Brown   ierr = PetscFree(GammaInv);CHKERRQ(ierr);
58043b21953SEmil Constantinescu 
58143b21953SEmil Constantinescu   for (i=0; i<s; i++) {
58243b21953SEmil Constantinescu     for (k=0; k<i+1; k++) {
58343b21953SEmil Constantinescu       t->GammaExplicitCorr[i*s+k]=(t->GammaExplicitCorr[i*s+k])*(t->GammaInv[k*s+k]);
58443b21953SEmil Constantinescu       for (j=k+1; j<i+1; j++) {
58543b21953SEmil Constantinescu         t->GammaExplicitCorr[i*s+k]+=(t->GammaExplicitCorr[i*s+j])*(t->GammaInv[j*s+k]);
58643b21953SEmil Constantinescu       }
58743b21953SEmil Constantinescu     }
58843b21953SEmil Constantinescu   }
58943b21953SEmil Constantinescu 
59061692a83SJed Brown   for (i=0; i<s; i++) {
59161692a83SJed Brown     for (j=0; j<s; j++) {
59261692a83SJed Brown       t->At[i*s+j] = 0;
59361692a83SJed Brown       for (k=0; k<s; k++) {
59461692a83SJed Brown         t->At[i*s+j] += t->A[i*s+k] * t->GammaInv[k*s+j];
59561692a83SJed Brown       }
59661692a83SJed Brown     }
59761692a83SJed Brown     t->bt[i] = 0;
59861692a83SJed Brown     for (j=0; j<s; j++) {
59961692a83SJed Brown       t->bt[i] += t->b[j] * t->GammaInv[j*s+i];
60061692a83SJed Brown     }
601fe7e6d57SJed Brown     if (bembed) {
602fe7e6d57SJed Brown       t->bembedt[i] = 0;
603fe7e6d57SJed Brown       for (j=0; j<s; j++) {
604fe7e6d57SJed Brown         t->bembedt[i] += t->bembed[j] * t->GammaInv[j*s+i];
605fe7e6d57SJed Brown       }
606fe7e6d57SJed Brown     }
60761692a83SJed Brown   }
6088d59e960SJed Brown   t->ccfl = 1.0;                /* Fix this */
6098d59e960SJed Brown 
61061692a83SJed Brown   link->next = RosWTableauList;
61161692a83SJed Brown   RosWTableauList = link;
612e27a552bSJed Brown   PetscFunctionReturn(0);
613e27a552bSJed Brown }
614e27a552bSJed Brown 
615e27a552bSJed Brown #undef __FUNCT__
6161c3436cfSJed Brown #define __FUNCT__ "TSEvaluateStep_RosW"
6171c3436cfSJed Brown /*
6181c3436cfSJed Brown  The step completion formula is
6191c3436cfSJed Brown 
6201c3436cfSJed Brown  x1 = x0 + b^T Y
6211c3436cfSJed Brown 
6221c3436cfSJed Brown  where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been
6231c3436cfSJed Brown  updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write
6241c3436cfSJed Brown 
6251c3436cfSJed Brown  x1e = x0 + be^T Y
6261c3436cfSJed Brown      = x1 - b^T Y + be^T Y
6271c3436cfSJed Brown      = x1 + (be - b)^T Y
6281c3436cfSJed Brown 
6291c3436cfSJed Brown  so we can evaluate the method of different order even after the step has been optimistically completed.
6301c3436cfSJed Brown */
6311c3436cfSJed Brown static PetscErrorCode TSEvaluateStep_RosW(TS ts,PetscInt order,Vec X,PetscBool *done)
6321c3436cfSJed Brown {
6331c3436cfSJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
6341c3436cfSJed Brown   RosWTableau    tab  = ros->tableau;
6351c3436cfSJed Brown   PetscScalar    *w = ros->work;
6361c3436cfSJed Brown   PetscInt       i;
6371c3436cfSJed Brown   PetscErrorCode ierr;
6381c3436cfSJed Brown 
6391c3436cfSJed Brown   PetscFunctionBegin;
6401c3436cfSJed Brown   if (order == tab->order) {
641108c343cSJed Brown     if (ros->status == TS_STEP_INCOMPLETE) { /* Use standard completion formula */
6421c3436cfSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
643de19f811SJed Brown       for (i=0; i<tab->s; i++) w[i] = tab->bt[i];
644de19f811SJed Brown       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
645108c343cSJed Brown     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
6461c3436cfSJed Brown     if (done) *done = PETSC_TRUE;
6471c3436cfSJed Brown     PetscFunctionReturn(0);
6481c3436cfSJed Brown   } else if (order == tab->order-1) {
6491c3436cfSJed Brown     if (!tab->bembedt) goto unavailable;
650108c343cSJed Brown     if (ros->status == TS_STEP_INCOMPLETE) { /* Use embedded completion formula */
6511c3436cfSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
652de19f811SJed Brown       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i];
653de19f811SJed Brown       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
654108c343cSJed Brown     } else {                    /* Use rollback-and-recomplete formula (bembedt - bt) */
655108c343cSJed Brown       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i];
656108c343cSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
657108c343cSJed Brown       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
6581c3436cfSJed Brown     }
6591c3436cfSJed Brown     if (done) *done = PETSC_TRUE;
6601c3436cfSJed Brown     PetscFunctionReturn(0);
6611c3436cfSJed Brown   }
6621c3436cfSJed Brown   unavailable:
6631c3436cfSJed Brown   if (done) *done = PETSC_FALSE;
6641c3436cfSJed 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);
6651c3436cfSJed Brown   PetscFunctionReturn(0);
6661c3436cfSJed Brown }
6671c3436cfSJed Brown 
6681c3436cfSJed Brown #undef __FUNCT__
669e27a552bSJed Brown #define __FUNCT__ "TSStep_RosW"
670e27a552bSJed Brown static PetscErrorCode TSStep_RosW(TS ts)
671e27a552bSJed Brown {
67261692a83SJed Brown   TS_RosW         *ros = (TS_RosW*)ts->data;
67361692a83SJed Brown   RosWTableau     tab  = ros->tableau;
674e27a552bSJed Brown   const PetscInt  s    = tab->s;
6751c3436cfSJed Brown   const PetscReal *At  = tab->At,*Gamma = tab->Gamma,*ASum = tab->ASum,*GammaInv = tab->GammaInv;
6760feba352SEmil Constantinescu   const PetscReal *GammaExplicitCorr = tab->GammaExplicitCorr;
677c17803e7SJed Brown   const PetscBool *GammaZeroDiag = tab->GammaZeroDiag;
67861692a83SJed Brown   PetscScalar     *w   = ros->work;
6797d4bf2deSEmil Constantinescu   Vec             *Y   = ros->Y,Ydot = ros->Ydot,Zdot = ros->Zdot,Zstage = ros->Zstage;
680e27a552bSJed Brown   SNES            snes;
6811c3436cfSJed Brown   TSAdapt         adapt;
6821c3436cfSJed Brown   PetscInt        i,j,its,lits,reject,next_scheme;
683cdbf8f93SLisandro Dalcin   PetscReal       next_time_step;
6841c3436cfSJed Brown   PetscBool       accept;
685e27a552bSJed Brown   PetscErrorCode  ierr;
6860feba352SEmil Constantinescu   MatStructure    str;
687e27a552bSJed Brown 
688e27a552bSJed Brown   PetscFunctionBegin;
689e27a552bSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
690cdbf8f93SLisandro Dalcin   next_time_step = ts->time_step;
6911c3436cfSJed Brown   accept = PETSC_TRUE;
692108c343cSJed Brown   ros->status = TS_STEP_INCOMPLETE;
693e27a552bSJed Brown 
69497335746SJed Brown   for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) {
6951c3436cfSJed Brown     const PetscReal h = ts->time_step;
696e27a552bSJed Brown     for (i=0; i<s; i++) {
6971c3436cfSJed Brown       ros->stage_time = ts->ptime + h*ASum[i];
698c17803e7SJed Brown       if (GammaZeroDiag[i]) {
699c17803e7SJed Brown         ros->stage_explicit = PETSC_TRUE;
700fd96d5b0SEmil Constantinescu         ros->shift = 1./h;
701c17803e7SJed Brown       } else {
702c17803e7SJed Brown         ros->stage_explicit = PETSC_FALSE;
70361692a83SJed Brown         ros->shift = 1./(h*Gamma[i*s+i]);
704fd96d5b0SEmil Constantinescu       }
70561692a83SJed Brown 
70661692a83SJed Brown       ierr = VecCopy(ts->vec_sol,Zstage);CHKERRQ(ierr);
707de19f811SJed Brown       for (j=0; j<i; j++) w[j] = At[i*s+j];
708de19f811SJed Brown       ierr = VecMAXPY(Zstage,i,w,Y);CHKERRQ(ierr);
70961692a83SJed Brown 
71061692a83SJed Brown       for (j=0; j<i; j++) w[j] = 1./h * GammaInv[i*s+j];
71161692a83SJed Brown       ierr = VecZeroEntries(Zdot);CHKERRQ(ierr);
71261692a83SJed Brown       ierr = VecMAXPY(Zdot,i,w,Y);CHKERRQ(ierr);
71361692a83SJed Brown 
714e27a552bSJed Brown       /* Initial guess taken from last stage */
71561692a83SJed Brown       ierr = VecZeroEntries(Y[i]);CHKERRQ(ierr);
71661692a83SJed Brown 
7177d4bf2deSEmil Constantinescu       if (!ros->stage_explicit) {
71861692a83SJed Brown         if (!ros->recompute_jacobian && !i) {
71961692a83SJed Brown           ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */
72061692a83SJed Brown         }
72161692a83SJed Brown         ierr = SNESSolve(snes,PETSC_NULL,Y[i]);CHKERRQ(ierr);
722e27a552bSJed Brown         ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
723e27a552bSJed Brown         ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
724e27a552bSJed Brown         ts->nonlinear_its += its; ts->linear_its += lits;
72597335746SJed Brown         ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
72697335746SJed Brown         ierr = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr);
72797335746SJed Brown         if (!accept) goto reject_step;
7287d4bf2deSEmil Constantinescu       } else {
7291ce71dffSSatish Balay         Mat J,Jp;
7300feba352SEmil Constantinescu         ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); /* Evaluate Y[i]=G(t,Ydot=0,Zstage) */
7310feba352SEmil Constantinescu         ierr = TSComputeIFunction(ts,ros->stage_time,Zstage,Ydot,Y[i],PETSC_FALSE);CHKERRQ(ierr);
7320feba352SEmil Constantinescu         ierr = VecScale(Y[i],-1.0);
7330feba352SEmil Constantinescu         ierr = VecAXPY(Y[i],-1.0,Zdot);CHKERRQ(ierr); /*Y[i]=F(Zstage)-Zdot[=GammaInv*Y]*/
7340feba352SEmil Constantinescu 
7350feba352SEmil Constantinescu         ierr = VecZeroEntries(Zstage);CHKERRQ(ierr); /* Zstage = GammaExplicitCorr[i,j] * Y[j] */
7360feba352SEmil Constantinescu         for (j=0; j<i; j++) w[j] = GammaExplicitCorr[i*s+j];
7370feba352SEmil Constantinescu         ierr = VecMAXPY(Zstage,i,w,Y);CHKERRQ(ierr);
7380feba352SEmil Constantinescu         /*Y[i] += Y[i] + Jac*Zstage[=Jac*GammaExplicitCorr[i,j] * Y[j]] */
7390feba352SEmil Constantinescu         str = SAME_NONZERO_PATTERN;
7400feba352SEmil Constantinescu         ierr = TSGetIJacobian(ts,&J,&Jp,PETSC_NULL,PETSC_NULL);
7410feba352SEmil Constantinescu         ierr = TSComputeIJacobian(ts,ros->stage_time,ts->vec_sol,Ydot,0,&J,&Jp,&str,PETSC_FALSE);CHKERRQ(ierr);
7420feba352SEmil Constantinescu         ierr = MatMult(J,Zstage,Zdot);
7430feba352SEmil Constantinescu 
7440feba352SEmil Constantinescu         ierr = VecAXPY(Y[i],-1.0,Zdot);CHKERRQ(ierr);
7450feba352SEmil Constantinescu         ierr = VecScale(Y[i],h);
7467d4bf2deSEmil Constantinescu         ts->linear_its += 1;
7477d4bf2deSEmil Constantinescu       }
748e27a552bSJed Brown     }
7491c3436cfSJed Brown     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
750108c343cSJed Brown     ros->status = TS_STEP_PENDING;
751e27a552bSJed Brown 
7521c3436cfSJed Brown     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
7531c3436cfSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
7541c3436cfSJed Brown     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
7558d59e960SJed Brown     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
7561c3436cfSJed Brown     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
7571c3436cfSJed Brown     if (accept) {
7581c3436cfSJed Brown       /* ignore next_scheme for now */
759e27a552bSJed Brown       ts->ptime += ts->time_step;
760cdbf8f93SLisandro Dalcin       ts->time_step = next_time_step;
761e27a552bSJed Brown       ts->steps++;
762108c343cSJed Brown       ros->status = TS_STEP_COMPLETE;
7631c3436cfSJed Brown       break;
7641c3436cfSJed Brown     } else {                    /* Roll back the current step */
7651c3436cfSJed Brown       for (i=0; i<s; i++) w[i] = -tab->bt[i];
7661c3436cfSJed Brown       ierr = VecMAXPY(ts->vec_sol,s,w,Y);CHKERRQ(ierr);
7671c3436cfSJed Brown       ts->time_step = next_time_step;
768108c343cSJed Brown       ros->status = TS_STEP_INCOMPLETE;
7691c3436cfSJed Brown     }
770476b6736SJed Brown     reject_step: continue;
7711c3436cfSJed Brown   }
772b2ce242eSJed Brown   if (ros->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
773e27a552bSJed Brown   PetscFunctionReturn(0);
774e27a552bSJed Brown }
775e27a552bSJed Brown 
776e27a552bSJed Brown #undef __FUNCT__
777e27a552bSJed Brown #define __FUNCT__ "TSInterpolate_RosW"
778e27a552bSJed Brown static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec X)
779e27a552bSJed Brown {
78061692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
781e27a552bSJed Brown 
782e27a552bSJed Brown   PetscFunctionBegin;
78361692a83SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name);
784e27a552bSJed Brown   PetscFunctionReturn(0);
785e27a552bSJed Brown }
786e27a552bSJed Brown 
787e27a552bSJed Brown /*------------------------------------------------------------*/
788e27a552bSJed Brown #undef __FUNCT__
789e27a552bSJed Brown #define __FUNCT__ "TSReset_RosW"
790e27a552bSJed Brown static PetscErrorCode TSReset_RosW(TS ts)
791e27a552bSJed Brown {
79261692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
793e27a552bSJed Brown   PetscInt       s;
794e27a552bSJed Brown   PetscErrorCode ierr;
795e27a552bSJed Brown 
796e27a552bSJed Brown   PetscFunctionBegin;
79761692a83SJed Brown   if (!ros->tableau) PetscFunctionReturn(0);
79861692a83SJed Brown    s = ros->tableau->s;
79961692a83SJed Brown   ierr = VecDestroyVecs(s,&ros->Y);CHKERRQ(ierr);
80061692a83SJed Brown   ierr = VecDestroy(&ros->Ydot);CHKERRQ(ierr);
80161692a83SJed Brown   ierr = VecDestroy(&ros->Ystage);CHKERRQ(ierr);
80261692a83SJed Brown   ierr = VecDestroy(&ros->Zdot);CHKERRQ(ierr);
80361692a83SJed Brown   ierr = VecDestroy(&ros->Zstage);CHKERRQ(ierr);
80461692a83SJed Brown   ierr = PetscFree(ros->work);CHKERRQ(ierr);
805e27a552bSJed Brown   PetscFunctionReturn(0);
806e27a552bSJed Brown }
807e27a552bSJed Brown 
808e27a552bSJed Brown #undef __FUNCT__
809e27a552bSJed Brown #define __FUNCT__ "TSDestroy_RosW"
810e27a552bSJed Brown static PetscErrorCode TSDestroy_RosW(TS ts)
811e27a552bSJed Brown {
812e27a552bSJed Brown   PetscErrorCode  ierr;
813e27a552bSJed Brown 
814e27a552bSJed Brown   PetscFunctionBegin;
815e27a552bSJed Brown   ierr = TSReset_RosW(ts);CHKERRQ(ierr);
816e27a552bSJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
817e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","",PETSC_NULL);CHKERRQ(ierr);
818e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","",PETSC_NULL);CHKERRQ(ierr);
81961692a83SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","",PETSC_NULL);CHKERRQ(ierr);
820e27a552bSJed Brown   PetscFunctionReturn(0);
821e27a552bSJed Brown }
822e27a552bSJed Brown 
823e27a552bSJed Brown /*
824e27a552bSJed Brown   This defines the nonlinear equation that is to be solved with SNES
825e27a552bSJed Brown   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
826e27a552bSJed Brown */
827e27a552bSJed Brown #undef __FUNCT__
828e27a552bSJed Brown #define __FUNCT__ "SNESTSFormFunction_RosW"
829e27a552bSJed Brown static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec X,Vec F,TS ts)
830e27a552bSJed Brown {
83161692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
832e27a552bSJed Brown   PetscErrorCode ierr;
833e27a552bSJed Brown 
834e27a552bSJed Brown   PetscFunctionBegin;
83561692a83SJed Brown   ierr = VecWAXPY(ros->Ydot,ros->shift,X,ros->Zdot);CHKERRQ(ierr); /* Ydot = shift*X + Zdot */
83661692a83SJed Brown   ierr = VecWAXPY(ros->Ystage,1.0,X,ros->Zstage);CHKERRQ(ierr);    /* Ystage = X + Zstage */
83761692a83SJed Brown   ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,ros->Ydot,F,PETSC_FALSE);CHKERRQ(ierr);
838e27a552bSJed Brown   PetscFunctionReturn(0);
839e27a552bSJed Brown }
840e27a552bSJed Brown 
841e27a552bSJed Brown #undef __FUNCT__
842e27a552bSJed Brown #define __FUNCT__ "SNESTSFormJacobian_RosW"
843e27a552bSJed Brown static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
844e27a552bSJed Brown {
84561692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
846e27a552bSJed Brown   PetscErrorCode ierr;
847e27a552bSJed Brown 
848e27a552bSJed Brown   PetscFunctionBegin;
84961692a83SJed Brown   /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
85061692a83SJed Brown   ierr = TSComputeIJacobian(ts,ros->stage_time,ros->Ystage,ros->Ydot,ros->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr);
851e27a552bSJed Brown   PetscFunctionReturn(0);
852e27a552bSJed Brown }
853e27a552bSJed Brown 
854e27a552bSJed Brown #undef __FUNCT__
855e27a552bSJed Brown #define __FUNCT__ "TSSetUp_RosW"
856e27a552bSJed Brown static PetscErrorCode TSSetUp_RosW(TS ts)
857e27a552bSJed Brown {
85861692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
85961692a83SJed Brown   RosWTableau    tab  = ros->tableau;
860e27a552bSJed Brown   PetscInt       s    = tab->s;
861e27a552bSJed Brown   PetscErrorCode ierr;
862e27a552bSJed Brown 
863e27a552bSJed Brown   PetscFunctionBegin;
86461692a83SJed Brown   if (!ros->tableau) {
865e27a552bSJed Brown     ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);
866e27a552bSJed Brown   }
86761692a83SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ros->Y);CHKERRQ(ierr);
86861692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Ydot);CHKERRQ(ierr);
86961692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Ystage);CHKERRQ(ierr);
87061692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Zdot);CHKERRQ(ierr);
87161692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Zstage);CHKERRQ(ierr);
87261692a83SJed Brown   ierr = PetscMalloc(s*sizeof(ros->work[0]),&ros->work);CHKERRQ(ierr);
873e27a552bSJed Brown   PetscFunctionReturn(0);
874e27a552bSJed Brown }
875e27a552bSJed Brown /*------------------------------------------------------------*/
876e27a552bSJed Brown 
877e27a552bSJed Brown #undef __FUNCT__
878e27a552bSJed Brown #define __FUNCT__ "TSSetFromOptions_RosW"
879e27a552bSJed Brown static PetscErrorCode TSSetFromOptions_RosW(TS ts)
880e27a552bSJed Brown {
88161692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
882e27a552bSJed Brown   PetscErrorCode ierr;
88361692a83SJed Brown   char           rostype[256];
884e27a552bSJed Brown 
885e27a552bSJed Brown   PetscFunctionBegin;
886e27a552bSJed Brown   ierr = PetscOptionsHead("RosW ODE solver options");CHKERRQ(ierr);
887e27a552bSJed Brown   {
88861692a83SJed Brown     RosWTableauLink link;
889e27a552bSJed Brown     PetscInt count,choice;
890e27a552bSJed Brown     PetscBool flg;
891e27a552bSJed Brown     const char **namelist;
89261692a83SJed Brown     SNES snes;
89361692a83SJed Brown 
89461692a83SJed Brown     ierr = PetscStrncpy(rostype,TSRosWDefault,sizeof rostype);CHKERRQ(ierr);
89561692a83SJed Brown     for (link=RosWTableauList,count=0; link; link=link->next,count++) ;
896e27a552bSJed Brown     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
89761692a83SJed Brown     for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
89861692a83SJed Brown     ierr = PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,rostype,&choice,&flg);CHKERRQ(ierr);
89961692a83SJed Brown     ierr = TSRosWSetType(ts,flg ? namelist[choice] : rostype);CHKERRQ(ierr);
900e27a552bSJed Brown     ierr = PetscFree(namelist);CHKERRQ(ierr);
90161692a83SJed Brown 
90261692a83SJed Brown     ierr = PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,PETSC_NULL);CHKERRQ(ierr);
90361692a83SJed Brown 
90461692a83SJed Brown     /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
90561692a83SJed Brown     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
90661692a83SJed Brown     if (!((PetscObject)snes)->type_name) {
90761692a83SJed Brown       ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
90861692a83SJed Brown     }
90961692a83SJed Brown     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
910e27a552bSJed Brown   }
911e27a552bSJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
912e27a552bSJed Brown   PetscFunctionReturn(0);
913e27a552bSJed Brown }
914e27a552bSJed Brown 
915e27a552bSJed Brown #undef __FUNCT__
916e27a552bSJed Brown #define __FUNCT__ "PetscFormatRealArray"
917e27a552bSJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
918e27a552bSJed Brown {
919e27a552bSJed Brown   PetscErrorCode ierr;
920e408995aSJed Brown   PetscInt i;
921e408995aSJed Brown   size_t left,count;
922e27a552bSJed Brown   char *p;
923e27a552bSJed Brown 
924e27a552bSJed Brown   PetscFunctionBegin;
925e408995aSJed Brown   for (i=0,p=buf,left=len; i<n; i++) {
926e408995aSJed Brown     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
927e27a552bSJed Brown     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
928e27a552bSJed Brown     left -= count;
929e27a552bSJed Brown     p += count;
930e27a552bSJed Brown     *p++ = ' ';
931e27a552bSJed Brown   }
932e27a552bSJed Brown   p[i ? 0 : -1] = 0;
933e27a552bSJed Brown   PetscFunctionReturn(0);
934e27a552bSJed Brown }
935e27a552bSJed Brown 
936e27a552bSJed Brown #undef __FUNCT__
937e27a552bSJed Brown #define __FUNCT__ "TSView_RosW"
938e27a552bSJed Brown static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer)
939e27a552bSJed Brown {
94061692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
94161692a83SJed Brown   RosWTableau    tab  = ros->tableau;
942e27a552bSJed Brown   PetscBool      iascii;
943e27a552bSJed Brown   PetscErrorCode ierr;
944e27a552bSJed Brown 
945e27a552bSJed Brown   PetscFunctionBegin;
946e27a552bSJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
947e27a552bSJed Brown   if (iascii) {
94861692a83SJed Brown     const TSRosWType rostype;
949e408995aSJed Brown     PetscInt i;
950e408995aSJed Brown     PetscReal abscissa[512];
951e27a552bSJed Brown     char buf[512];
95261692a83SJed Brown     ierr = TSRosWGetType(ts,&rostype);CHKERRQ(ierr);
95361692a83SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Rosenbrock-W %s\n",rostype);CHKERRQ(ierr);
954e408995aSJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ASum);CHKERRQ(ierr);
95561692a83SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A       = %s\n",buf);CHKERRQ(ierr);
956e408995aSJed Brown     for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i];
957e408995aSJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,abscissa);CHKERRQ(ierr);
958e408995aSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A+Gamma = %s\n",buf);CHKERRQ(ierr);
959e27a552bSJed Brown   }
960e27a552bSJed Brown   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
961e27a552bSJed Brown   PetscFunctionReturn(0);
962e27a552bSJed Brown }
963e27a552bSJed Brown 
964e27a552bSJed Brown #undef __FUNCT__
965e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType"
966e27a552bSJed Brown /*@C
96761692a83SJed Brown   TSRosWSetType - Set the type of Rosenbrock-W scheme
968e27a552bSJed Brown 
969e27a552bSJed Brown   Logically collective
970e27a552bSJed Brown 
971e27a552bSJed Brown   Input Parameter:
972e27a552bSJed Brown +  ts - timestepping context
97361692a83SJed Brown -  rostype - type of Rosenbrock-W scheme
974e27a552bSJed Brown 
975020d8f30SJed Brown   Level: beginner
976e27a552bSJed Brown 
977020d8f30SJed Brown .seealso: TSRosWGetType(), TSROSW, TSROSW2M, TSROSW2P, TSROSWRA3PW, TSROSWRA34PW2, TSROSWRODAS3, TSROSWSANDU3, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, TSROSWARK3
978e27a552bSJed Brown @*/
97961692a83SJed Brown PetscErrorCode TSRosWSetType(TS ts,const TSRosWType rostype)
980e27a552bSJed Brown {
981e27a552bSJed Brown   PetscErrorCode ierr;
982e27a552bSJed Brown 
983e27a552bSJed Brown   PetscFunctionBegin;
984e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
98561692a83SJed Brown   ierr = PetscTryMethod(ts,"TSRosWSetType_C",(TS,const TSRosWType),(ts,rostype));CHKERRQ(ierr);
986e27a552bSJed Brown   PetscFunctionReturn(0);
987e27a552bSJed Brown }
988e27a552bSJed Brown 
989e27a552bSJed Brown #undef __FUNCT__
990e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType"
991e27a552bSJed Brown /*@C
99261692a83SJed Brown   TSRosWGetType - Get the type of Rosenbrock-W scheme
993e27a552bSJed Brown 
994e27a552bSJed Brown   Logically collective
995e27a552bSJed Brown 
996e27a552bSJed Brown   Input Parameter:
997e27a552bSJed Brown .  ts - timestepping context
998e27a552bSJed Brown 
999e27a552bSJed Brown   Output Parameter:
100061692a83SJed Brown .  rostype - type of Rosenbrock-W scheme
1001e27a552bSJed Brown 
1002e27a552bSJed Brown   Level: intermediate
1003e27a552bSJed Brown 
1004e27a552bSJed Brown .seealso: TSRosWGetType()
1005e27a552bSJed Brown @*/
100661692a83SJed Brown PetscErrorCode TSRosWGetType(TS ts,const TSRosWType *rostype)
1007e27a552bSJed Brown {
1008e27a552bSJed Brown   PetscErrorCode ierr;
1009e27a552bSJed Brown 
1010e27a552bSJed Brown   PetscFunctionBegin;
1011e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
101261692a83SJed Brown   ierr = PetscUseMethod(ts,"TSRosWGetType_C",(TS,const TSRosWType*),(ts,rostype));CHKERRQ(ierr);
1013e27a552bSJed Brown   PetscFunctionReturn(0);
1014e27a552bSJed Brown }
1015e27a552bSJed Brown 
1016e27a552bSJed Brown #undef __FUNCT__
101761692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian"
1018e27a552bSJed Brown /*@C
101961692a83SJed Brown   TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step.
1020e27a552bSJed Brown 
1021e27a552bSJed Brown   Logically collective
1022e27a552bSJed Brown 
1023e27a552bSJed Brown   Input Parameter:
1024e27a552bSJed Brown +  ts - timestepping context
102561692a83SJed Brown -  flg - PETSC_TRUE to recompute the Jacobian at each stage
1026e27a552bSJed Brown 
1027e27a552bSJed Brown   Level: intermediate
1028e27a552bSJed Brown 
1029e27a552bSJed Brown .seealso: TSRosWGetType()
1030e27a552bSJed Brown @*/
103161692a83SJed Brown PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg)
1032e27a552bSJed Brown {
1033e27a552bSJed Brown   PetscErrorCode ierr;
1034e27a552bSJed Brown 
1035e27a552bSJed Brown   PetscFunctionBegin;
1036e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
103761692a83SJed Brown   ierr = PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
1038e27a552bSJed Brown   PetscFunctionReturn(0);
1039e27a552bSJed Brown }
1040e27a552bSJed Brown 
1041e27a552bSJed Brown EXTERN_C_BEGIN
1042e27a552bSJed Brown #undef __FUNCT__
1043e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType_RosW"
104461692a83SJed Brown PetscErrorCode  TSRosWGetType_RosW(TS ts,const TSRosWType *rostype)
1045e27a552bSJed Brown {
104661692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
1047e27a552bSJed Brown   PetscErrorCode ierr;
1048e27a552bSJed Brown 
1049e27a552bSJed Brown   PetscFunctionBegin;
105061692a83SJed Brown   if (!ros->tableau) {ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);}
105161692a83SJed Brown   *rostype = ros->tableau->name;
1052e27a552bSJed Brown   PetscFunctionReturn(0);
1053e27a552bSJed Brown }
1054e27a552bSJed Brown #undef __FUNCT__
1055e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType_RosW"
105661692a83SJed Brown PetscErrorCode  TSRosWSetType_RosW(TS ts,const TSRosWType rostype)
1057e27a552bSJed Brown {
105861692a83SJed Brown   TS_RosW         *ros = (TS_RosW*)ts->data;
1059e27a552bSJed Brown   PetscErrorCode  ierr;
1060e27a552bSJed Brown   PetscBool       match;
106161692a83SJed Brown   RosWTableauLink link;
1062e27a552bSJed Brown 
1063e27a552bSJed Brown   PetscFunctionBegin;
106461692a83SJed Brown   if (ros->tableau) {
106561692a83SJed Brown     ierr = PetscStrcmp(ros->tableau->name,rostype,&match);CHKERRQ(ierr);
1066e27a552bSJed Brown     if (match) PetscFunctionReturn(0);
1067e27a552bSJed Brown   }
106861692a83SJed Brown   for (link = RosWTableauList; link; link=link->next) {
106961692a83SJed Brown     ierr = PetscStrcmp(link->tab.name,rostype,&match);CHKERRQ(ierr);
1070e27a552bSJed Brown     if (match) {
1071e27a552bSJed Brown       ierr = TSReset_RosW(ts);CHKERRQ(ierr);
107261692a83SJed Brown       ros->tableau = &link->tab;
1073e27a552bSJed Brown       PetscFunctionReturn(0);
1074e27a552bSJed Brown     }
1075e27a552bSJed Brown   }
107661692a83SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype);
1077e27a552bSJed Brown   PetscFunctionReturn(0);
1078e27a552bSJed Brown }
107961692a83SJed Brown 
1080e27a552bSJed Brown #undef __FUNCT__
108161692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian_RosW"
108261692a83SJed Brown PetscErrorCode  TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg)
1083e27a552bSJed Brown {
108461692a83SJed Brown   TS_RosW *ros = (TS_RosW*)ts->data;
1085e27a552bSJed Brown 
1086e27a552bSJed Brown   PetscFunctionBegin;
108761692a83SJed Brown   ros->recompute_jacobian = flg;
1088e27a552bSJed Brown   PetscFunctionReturn(0);
1089e27a552bSJed Brown }
1090e27a552bSJed Brown EXTERN_C_END
1091e27a552bSJed Brown 
1092e27a552bSJed Brown /* ------------------------------------------------------------ */
1093e27a552bSJed Brown /*MC
1094020d8f30SJed Brown       TSROSW - ODE solver using Rosenbrock-W schemes
1095e27a552bSJed Brown 
1096e27a552bSJed Brown   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
1097e27a552bSJed Brown   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
1098e27a552bSJed Brown   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
1099e27a552bSJed Brown 
1100e27a552bSJed Brown   Notes:
110161692a83SJed Brown   This method currently only works with autonomous ODE and DAE.
110261692a83SJed Brown 
110361692a83SJed Brown   Developer notes:
110461692a83SJed Brown   Rosenbrock-W methods are typically specified for autonomous ODE
110561692a83SJed Brown 
110661692a83SJed Brown $  xdot = f(x)
110761692a83SJed Brown 
110861692a83SJed Brown   by the stage equations
110961692a83SJed Brown 
111061692a83SJed Brown $  k_i = h f(x_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j
111161692a83SJed Brown 
111261692a83SJed Brown   and step completion formula
111361692a83SJed Brown 
111461692a83SJed Brown $  x_1 = x_0 + sum_j b_j k_j
111561692a83SJed Brown 
111661692a83SJed Brown   with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(x)
111761692a83SJed Brown   and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner,
111861692a83SJed Brown   we define new variables for the stage equations
111961692a83SJed Brown 
112061692a83SJed Brown $  y_i = gamma_ij k_j
112161692a83SJed Brown 
112261692a83SJed Brown   The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define
112361692a83SJed Brown 
112461692a83SJed Brown $  A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-i}
112561692a83SJed Brown 
112661692a83SJed Brown   to rewrite the method as
112761692a83SJed Brown 
112861692a83SJed 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
112961692a83SJed Brown $  x_1 = x_0 + sum_j bt_j y_j
113061692a83SJed Brown 
113161692a83SJed Brown    where we have introduced the mass matrix M. Continue by defining
113261692a83SJed Brown 
113361692a83SJed Brown $  ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j
113461692a83SJed Brown 
113561692a83SJed Brown    or, more compactly in tensor notation
113661692a83SJed Brown 
113761692a83SJed Brown $  Ydot = 1/h (Gamma^{-1} \otimes I) Y .
113861692a83SJed Brown 
113961692a83SJed Brown    Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current
114061692a83SJed Brown    stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the
114161692a83SJed Brown    equation
114261692a83SJed Brown 
114361692a83SJed Brown $  g(x_0 + sum_j a_ij y_j + y_i, ydot_i) = 0
114461692a83SJed Brown 
114561692a83SJed Brown    with initial guess y_i = 0.
1146e27a552bSJed Brown 
1147e27a552bSJed Brown   Level: beginner
1148e27a552bSJed Brown 
1149020d8f30SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSRosWSetType(), TSRosWRegister()
1150e27a552bSJed Brown 
1151e27a552bSJed Brown M*/
1152e27a552bSJed Brown EXTERN_C_BEGIN
1153e27a552bSJed Brown #undef __FUNCT__
1154e27a552bSJed Brown #define __FUNCT__ "TSCreate_RosW"
1155e27a552bSJed Brown PetscErrorCode  TSCreate_RosW(TS ts)
1156e27a552bSJed Brown {
115761692a83SJed Brown   TS_RosW        *ros;
1158e27a552bSJed Brown   PetscErrorCode ierr;
1159e27a552bSJed Brown 
1160e27a552bSJed Brown   PetscFunctionBegin;
1161e27a552bSJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
1162e27a552bSJed Brown   ierr = TSRosWInitializePackage(PETSC_NULL);CHKERRQ(ierr);
1163e27a552bSJed Brown #endif
1164e27a552bSJed Brown 
1165e27a552bSJed Brown   ts->ops->reset          = TSReset_RosW;
1166e27a552bSJed Brown   ts->ops->destroy        = TSDestroy_RosW;
1167e27a552bSJed Brown   ts->ops->view           = TSView_RosW;
1168e27a552bSJed Brown   ts->ops->setup          = TSSetUp_RosW;
1169e27a552bSJed Brown   ts->ops->step           = TSStep_RosW;
1170e27a552bSJed Brown   ts->ops->interpolate    = TSInterpolate_RosW;
11711c3436cfSJed Brown   ts->ops->evaluatestep   = TSEvaluateStep_RosW;
1172e27a552bSJed Brown   ts->ops->setfromoptions = TSSetFromOptions_RosW;
1173e27a552bSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_RosW;
1174e27a552bSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_RosW;
1175e27a552bSJed Brown 
117661692a83SJed Brown   ierr = PetscNewLog(ts,TS_RosW,&ros);CHKERRQ(ierr);
117761692a83SJed Brown   ts->data = (void*)ros;
1178e27a552bSJed Brown 
1179e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","TSRosWGetType_RosW",TSRosWGetType_RosW);CHKERRQ(ierr);
1180e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","TSRosWSetType_RosW",TSRosWSetType_RosW);CHKERRQ(ierr);
118161692a83SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","TSRosWSetRecomputeJacobian_RosW",TSRosWSetRecomputeJacobian_RosW);CHKERRQ(ierr);
1182e27a552bSJed Brown   PetscFunctionReturn(0);
1183e27a552bSJed Brown }
1184e27a552bSJed Brown EXTERN_C_END
1185