xref: /petsc/src/ts/impls/rosw/rosw.c (revision 961f28d0a49ee45ff1f844d187d293ed8e84b6a1)
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 
1761692a83SJed Brown static const TSRosWType TSRosWDefault = TSROSW2P;
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 */
2961692a83SJed Brown   PetscReal *b;                 /* Step completion table */
30fe7e6d57SJed Brown   PetscReal *bembed;            /* Step completion table for embedded method of order one less */
3161692a83SJed Brown   PetscReal *ASum;              /* Row sum of A */
3261692a83SJed Brown   PetscReal *GammaSum;          /* Row sum of Gamma, only needed for non-autonomous systems */
3361692a83SJed Brown   PetscReal *At;                /* Propagation table in transformed variables */
3461692a83SJed Brown   PetscReal *bt;                /* Step completion table in transformed variables */
35fe7e6d57SJed Brown   PetscReal *bembedt;           /* Step completion table of order one less in transformed variables */
3661692a83SJed Brown   PetscReal *GammaInv;          /* Inverse of Gamma, used for transformed variables */
378d59e960SJed Brown   PetscReal  ccfl;              /* Placeholder for CFL coefficient relative to forward Euler */
38e27a552bSJed Brown };
3961692a83SJed Brown typedef struct _RosWTableauLink *RosWTableauLink;
4061692a83SJed Brown struct _RosWTableauLink {
4161692a83SJed Brown   struct _RosWTableau tab;
4261692a83SJed Brown   RosWTableauLink next;
43e27a552bSJed Brown };
4461692a83SJed Brown static RosWTableauLink RosWTableauList;
45e27a552bSJed Brown 
46e27a552bSJed Brown typedef struct {
4761692a83SJed Brown   RosWTableau tableau;
4861692a83SJed Brown   Vec         *Y;               /* States computed during the step, used to complete the step */
49e27a552bSJed Brown   Vec         Ydot;             /* Work vector holding Ydot during residual evaluation */
5061692a83SJed Brown   Vec         Ystage;           /* Work vector for the state value at each stage */
5161692a83SJed Brown   Vec         Zdot;             /* Ydot = Zdot + shift*Y */
5261692a83SJed Brown   Vec         Zstage;           /* Y = Zstage + Y */
531c3436cfSJed Brown   PetscScalar *work;            /* Scalar work space of length number of stages, used to prepare VecMAXPY() */
54e27a552bSJed Brown   PetscReal   shift;
55e27a552bSJed Brown   PetscReal   stage_time;
56c17803e7SJed Brown   PetscReal   stage_explicit;     /* Flag indicates that the current stage is explicit */
5761692a83SJed Brown   PetscBool   recompute_jacobian; /* Recompute the Jacobian at each stage, default is to freeze the Jacobian at the start of each step */
581c3436cfSJed Brown   PetscBool   step_taken;         /* ts->vec_sol has been advanced to the end of the current time step */
59e27a552bSJed Brown } TS_RosW;
60e27a552bSJed Brown 
61fe7e6d57SJed Brown /*MC
62fe7e6d57SJed Brown      TSROSW2M - Two stage second order L-stable Rosenbrock-W scheme.
63fe7e6d57SJed Brown 
64fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2P.
65fe7e6d57SJed Brown 
66fe7e6d57SJed Brown      Level: intermediate
67fe7e6d57SJed Brown 
68fe7e6d57SJed Brown .seealso: TSROSW
69fe7e6d57SJed Brown M*/
70fe7e6d57SJed Brown 
71fe7e6d57SJed Brown /*MC
72fe7e6d57SJed Brown      TSROSW2P - Two stage second order L-stable Rosenbrock-W scheme.
73fe7e6d57SJed Brown 
74fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2M.
75fe7e6d57SJed Brown 
76fe7e6d57SJed Brown      Level: intermediate
77fe7e6d57SJed Brown 
78fe7e6d57SJed Brown .seealso: TSROSW
79fe7e6d57SJed Brown M*/
80fe7e6d57SJed Brown 
81fe7e6d57SJed Brown /*MC
82fe7e6d57SJed Brown      TSROSWRA3PW - Three stage third order Rosenbrock-W scheme for PDAE of index 1.
83fe7e6d57SJed Brown 
84fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.
85fe7e6d57SJed Brown 
86fe7e6d57SJed 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.
87fe7e6d57SJed Brown 
88fe7e6d57SJed Brown      References:
89fe7e6d57SJed Brown      Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005.
90fe7e6d57SJed Brown 
91fe7e6d57SJed Brown      Level: intermediate
92fe7e6d57SJed Brown 
93fe7e6d57SJed Brown .seealso: TSROSW
94fe7e6d57SJed Brown M*/
95fe7e6d57SJed Brown 
96fe7e6d57SJed Brown /*MC
97fe7e6d57SJed Brown      TSROSWRA34PW2 - Four stage third order L-stable Rosenbrock-W scheme for PDAE of index 1.
98fe7e6d57SJed Brown 
99fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.
100fe7e6d57SJed Brown 
101fe7e6d57SJed 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.
102fe7e6d57SJed Brown 
103fe7e6d57SJed Brown      References:
104fe7e6d57SJed Brown      Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005.
105fe7e6d57SJed Brown 
106fe7e6d57SJed Brown      Level: intermediate
107fe7e6d57SJed Brown 
108fe7e6d57SJed Brown .seealso: TSROSW
109fe7e6d57SJed Brown M*/
110fe7e6d57SJed Brown 
111ef3c5b88SJed Brown /*MC
112ef3c5b88SJed Brown      TSROSWRODAS3 - Four stage third order L-stable Rosenbrock scheme
113ef3c5b88SJed Brown 
114ef3c5b88SJed Brown      By default, the Jacobian is only recomputed once per step.
115ef3c5b88SJed Brown 
116ef3c5b88SJed Brown      Both the third order and embedded second order methods are stiffly accurate and L-stable.
117ef3c5b88SJed Brown 
118ef3c5b88SJed Brown      References:
119ef3c5b88SJed Brown      Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997.
120ef3c5b88SJed Brown 
121ef3c5b88SJed Brown      Level: intermediate
122ef3c5b88SJed Brown 
123ef3c5b88SJed Brown .seealso: TSROSW, TSROSWSANDU3
124ef3c5b88SJed Brown M*/
125ef3c5b88SJed Brown 
126ef3c5b88SJed Brown /*MC
127ef3c5b88SJed Brown      TSROSWSANDU3 - Three stage third order L-stable Rosenbrock scheme
128ef3c5b88SJed Brown 
129ef3c5b88SJed Brown      By default, the Jacobian is only recomputed once per step.
130ef3c5b88SJed Brown 
131ef3c5b88SJed Brown      The third order method is L-stable, but not stiffly accurate.
132ef3c5b88SJed Brown      The second order embedded method is strongly A-stable with R(infty) = 0.5.
133ef3c5b88SJed Brown      The internal stages are L-stable.
134ef3c5b88SJed Brown      This method is called ROS3 in the paper.
135ef3c5b88SJed Brown 
136ef3c5b88SJed Brown      References:
137ef3c5b88SJed Brown      Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997.
138ef3c5b88SJed Brown 
139ef3c5b88SJed Brown      Level: intermediate
140ef3c5b88SJed Brown 
141ef3c5b88SJed Brown .seealso: TSROSW, TSROSWRODAS3
142ef3c5b88SJed Brown M*/
143ef3c5b88SJed Brown 
144*961f28d0SJed Brown /*MC
145*961f28d0SJed Brown      TSROSWASSP3P3S1C - A-stable Rosenbrock-W method with SSP explicit part, third order, three stages
146*961f28d0SJed Brown 
147*961f28d0SJed Brown      By default, the Jacobian is only recomputed once per step.
148*961f28d0SJed Brown 
149*961f28d0SJed Brown      A-stable SPP explicit order 3, 3 stages, CFL 1 (eff = 1/3)
150*961f28d0SJed Brown 
151*961f28d0SJed Brown      References:
152*961f28d0SJed Brown      Emil Constantinescu
153*961f28d0SJed Brown 
154*961f28d0SJed Brown      Level: intermediate
155*961f28d0SJed Brown 
156*961f28d0SJed Brown .seealso: TSROSW, TSROSWLASSP3P4S2C, TSROSWLLSSP3P3S2C, SSP
157*961f28d0SJed Brown M*/
158*961f28d0SJed Brown 
159*961f28d0SJed Brown /*MC
160*961f28d0SJed Brown      TSROSWLASSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, three stages
161*961f28d0SJed Brown 
162*961f28d0SJed Brown      By default, the Jacobian is only recomputed once per step.
163*961f28d0SJed Brown 
164*961f28d0SJed Brown      L-stable (A-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2)
165*961f28d0SJed Brown 
166*961f28d0SJed Brown      References:
167*961f28d0SJed Brown      Emil Constantinescu
168*961f28d0SJed Brown 
169*961f28d0SJed Brown      Level: intermediate
170*961f28d0SJed Brown 
171*961f28d0SJed Brown .seealso: TSROSW, TSROSWASSP3P3S1C, TSROSWLLSSP3P3S2C, TSSSP
172*961f28d0SJed Brown M*/
173*961f28d0SJed Brown 
174*961f28d0SJed Brown /*MC
175*961f28d0SJed Brown      TSROSWLLSSP3P3S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, three stages
176*961f28d0SJed Brown 
177*961f28d0SJed Brown      By default, the Jacobian is only recomputed once per step.
178*961f28d0SJed Brown 
179*961f28d0SJed Brown      L-stable (L-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2)
180*961f28d0SJed Brown 
181*961f28d0SJed Brown      References:
182*961f28d0SJed Brown      Emil Constantinescu
183*961f28d0SJed Brown 
184*961f28d0SJed Brown      Level: intermediate
185*961f28d0SJed Brown 
186*961f28d0SJed Brown .seealso: TSROSW, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSSSP
187*961f28d0SJed Brown M*/
188*961f28d0SJed Brown 
189e27a552bSJed Brown #undef __FUNCT__
190e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterAll"
191e27a552bSJed Brown /*@C
192e27a552bSJed Brown   TSRosWRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSRosW
193e27a552bSJed Brown 
194e27a552bSJed Brown   Not Collective, but should be called by all processes which will need the schemes to be registered
195e27a552bSJed Brown 
196e27a552bSJed Brown   Level: advanced
197e27a552bSJed Brown 
198e27a552bSJed Brown .keywords: TS, TSRosW, register, all
199e27a552bSJed Brown 
200e27a552bSJed Brown .seealso:  TSRosWRegisterDestroy()
201e27a552bSJed Brown @*/
202e27a552bSJed Brown PetscErrorCode TSRosWRegisterAll(void)
203e27a552bSJed Brown {
204e27a552bSJed Brown   PetscErrorCode ierr;
205e27a552bSJed Brown 
206e27a552bSJed Brown   PetscFunctionBegin;
207e27a552bSJed Brown   if (TSRosWRegisterAllCalled) PetscFunctionReturn(0);
208e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_TRUE;
209e27a552bSJed Brown 
210e27a552bSJed Brown   {
21161692a83SJed Brown     const PetscReal g = 1. + 1./PetscSqrtReal(2.0);
212e27a552bSJed Brown     const PetscReal
21361692a83SJed Brown       A[2][2] = {{0,0}, {1.,0}},
21461692a83SJed Brown       Gamma[2][2] = {{g,0}, {-2.*g,g}},
2151c3436cfSJed Brown       b[2] = {0.5,0.5},
2161c3436cfSJed Brown       b1[2] = {1.0,0.0};
2171c3436cfSJed Brown     ierr = TSRosWRegister(TSROSW2P,2,2,&A[0][0],&Gamma[0][0],b,b1);CHKERRQ(ierr);
218e27a552bSJed Brown   }
219e27a552bSJed Brown   {
22061692a83SJed Brown     const PetscReal g = 1. - 1./PetscSqrtReal(2.0);
221e27a552bSJed Brown     const PetscReal
22261692a83SJed Brown       A[2][2] = {{0,0}, {1.,0}},
22361692a83SJed Brown       Gamma[2][2] = {{g,0}, {-2.*g,g}},
2241c3436cfSJed Brown       b[2] = {0.5,0.5},
2251c3436cfSJed Brown       b1[2] = {1.0,0.0};
2261c3436cfSJed Brown     ierr = TSRosWRegister(TSROSW2M,2,2,&A[0][0],&Gamma[0][0],b,b1);CHKERRQ(ierr);
227fe7e6d57SJed Brown   }
228fe7e6d57SJed Brown   {
229fe7e6d57SJed Brown     const PetscReal g = 7.8867513459481287e-01;
230fe7e6d57SJed Brown     const PetscReal
231fe7e6d57SJed Brown       A[3][3] = {{0,0,0},
232fe7e6d57SJed Brown                  {1.5773502691896257e+00,0,0},
233fe7e6d57SJed Brown                  {0.5,0,0}},
234fe7e6d57SJed Brown       Gamma[3][3] = {{g,0,0},
235fe7e6d57SJed Brown                      {-1.5773502691896257e+00,g,0},
236fe7e6d57SJed Brown                      {-6.7075317547305480e-01,1.7075317547305482e-01,g}},
237fe7e6d57SJed Brown       b[3] = {1.0566243270259355e-01,4.9038105676657971e-02,8.4529946162074843e-01},
238fe7e6d57SJed Brown       b2[3] = {-1.7863279495408180e-01,1./3.,8.4529946162074843e-01};
239fe7e6d57SJed Brown     ierr = TSRosWRegister(TSROSWRA3PW,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
240fe7e6d57SJed Brown   }
241fe7e6d57SJed Brown   {
242fe7e6d57SJed Brown     const PetscReal g = 4.3586652150845900e-01;
243fe7e6d57SJed Brown     const PetscReal
244fe7e6d57SJed Brown       A[4][4] = {{0,0,0,0},
245fe7e6d57SJed Brown                  {8.7173304301691801e-01,0,0,0},
246fe7e6d57SJed Brown                  {8.4457060015369423e-01,-1.1299064236484185e-01,0,0},
247fe7e6d57SJed Brown                  {0,0,1.,0}},
248fe7e6d57SJed Brown       Gamma[4][4] = {{g,0,0,0},
249fe7e6d57SJed Brown                      {-8.7173304301691801e-01,g,0,0},
250fe7e6d57SJed Brown                      {-9.0338057013044082e-01,5.4180672388095326e-02,g,0},
251fe7e6d57SJed Brown                      {2.4212380706095346e-01,-1.2232505839045147e+00,5.4526025533510214e-01,g}},
252fe7e6d57SJed Brown       b[4] = {2.4212380706095346e-01,-1.2232505839045147e+00,1.5452602553351020e+00,4.3586652150845900e-01},
253fe7e6d57SJed Brown       b2[4] = {3.7810903145819369e-01,-9.6042292212423178e-02,5.0000000000000000e-01,2.1793326075422950e-01};
254fe7e6d57SJed Brown     ierr = TSRosWRegister(TSROSWRA34PW2,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
255e27a552bSJed Brown   }
256ef3c5b88SJed Brown   {
257ef3c5b88SJed Brown     const PetscReal g = 0.5;
258ef3c5b88SJed Brown     const PetscReal
259ef3c5b88SJed Brown       A[4][4] = {{0,0,0,0},
260ef3c5b88SJed Brown                  {0,0,0,0},
261ef3c5b88SJed Brown                  {1.,0,0,0},
262ef3c5b88SJed Brown                  {0.75,-0.25,0.5,0}},
263ef3c5b88SJed Brown       Gamma[4][4] = {{g,0,0,0},
264ef3c5b88SJed Brown                      {1.,g,0,0},
265ef3c5b88SJed Brown                      {-0.25,-0.25,g,0},
266ef3c5b88SJed Brown                      {1./12,1./12,-2./3,g}},
267ef3c5b88SJed Brown       b[4] = {5./6,-1./6,-1./6,0.5},
268ef3c5b88SJed Brown       b2[4] = {0.75,-0.25,0.5,0};
269ef3c5b88SJed Brown     ierr = TSRosWRegister(TSROSWRODAS3,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
270ef3c5b88SJed Brown   }
271ef3c5b88SJed Brown   {
272ef3c5b88SJed Brown     const PetscReal g = 0.43586652150845899941601945119356;
273ef3c5b88SJed Brown     const PetscReal
274ef3c5b88SJed Brown       A[3][3] = {{0,0,0},
275ef3c5b88SJed Brown                  {g,0,0},
276ef3c5b88SJed Brown                  {g,0,0}},
277ef3c5b88SJed Brown       Gamma[3][3] = {{g,0,0},
278ef3c5b88SJed Brown                      {-0.19294655696029095575009695436041,g,0},
279ef3c5b88SJed Brown                      {0,1.74927148125794685173529749738960,g}},
280ef3c5b88SJed Brown       b[3] = {-0.75457412385404315829818998646589,1.94100407061964420292840123379419,-0.18642994676560104463021124732829},
281ef3c5b88SJed Brown       b2[3] = {-1.53358745784149585370766523913002,2.81745131148625772213931745457622,-0.28386385364476186843165221544619};
282ef3c5b88SJed Brown     ierr = TSRosWRegister(TSROSWSANDU3,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
283ef3c5b88SJed Brown   }
284b1c69cc3SEmil Constantinescu   {
285b1c69cc3SEmil Constantinescu     const PetscReal g = (3.0+sqrt(3.0))/6.0;
286b1c69cc3SEmil Constantinescu     const PetscReal
287b1c69cc3SEmil Constantinescu       A[3][3] = {{0,0,0},
288b1c69cc3SEmil Constantinescu                  {1,0,0},
289b1c69cc3SEmil Constantinescu                  {0.25,0.25,0}},
290b1c69cc3SEmil Constantinescu       Gamma[3][3] = {{0,0,0},
291b1c69cc3SEmil Constantinescu                      {(-3.0-sqrt(3.0))/6.0,g,0},
292b1c69cc3SEmil Constantinescu                      {(-3.0-sqrt(3.0))/24.0,(-3.0-sqrt(3.0))/8.0,g}},
293b1c69cc3SEmil Constantinescu         b[3] = {1./6.,1./6.,2./3.},
294b1c69cc3SEmil Constantinescu           b2[3] = {1./4.,1./4.,1./2.};
295b1c69cc3SEmil Constantinescu     ierr = TSRosWRegister(TSROSWASSP3P3S1C,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
296b1c69cc3SEmil Constantinescu   }
297b1c69cc3SEmil Constantinescu 
298b1c69cc3SEmil Constantinescu   {
299b1c69cc3SEmil Constantinescu     const PetscReal
300b1c69cc3SEmil Constantinescu       A[4][4] = {{0,0,0,0},
301b1c69cc3SEmil Constantinescu                  {1./2.,0,0,0},
302b1c69cc3SEmil Constantinescu                  {1./2.,1./2.,0,0},
303b1c69cc3SEmil Constantinescu                  {1./6.,1./6.,1./6.,0}},
304b1c69cc3SEmil Constantinescu       Gamma[4][4] = {{1./2.,0,0,0},
305b1c69cc3SEmil Constantinescu                      {0.0,1./4.,0,0},
306b1c69cc3SEmil Constantinescu                      {-2.,-2./3.,2./3.,0},
307b1c69cc3SEmil Constantinescu                      {1./2.,5./36.,-2./9,0}},
308b1c69cc3SEmil Constantinescu         b[4] = {1./6.,1./6.,1./6.,1./2.},
309b1c69cc3SEmil Constantinescu         b2[4] = {1./8.,3./4.,1./8.,0};
310b1c69cc3SEmil Constantinescu      ierr = TSRosWRegister(TSROSWLASSP3P4S2C,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
311b1c69cc3SEmil Constantinescu   }
312b1c69cc3SEmil Constantinescu 
313b1c69cc3SEmil Constantinescu   {
314b1c69cc3SEmil Constantinescu     const PetscReal
315b1c69cc3SEmil Constantinescu       A[4][4] = {{0,0,0,0},
316b1c69cc3SEmil Constantinescu                  {1./2.,0,0,0},
317b1c69cc3SEmil Constantinescu                  {1./2.,1./2.,0,0},
318b1c69cc3SEmil Constantinescu                  {1./6.,1./6.,1./6.,0}},
319b1c69cc3SEmil Constantinescu       Gamma[4][4] = {{1./2.,0,0,0},
320b1c69cc3SEmil Constantinescu                      {0.0,3./4.,0,0},
321b1c69cc3SEmil Constantinescu                      {-2./3.,-23./9.,2./9.,0},
322b1c69cc3SEmil Constantinescu                      {1./18.,65./108.,-2./27,0}},
323b1c69cc3SEmil Constantinescu         b[4] = {1./6.,1./6.,1./6.,1./2.},
324b1c69cc3SEmil Constantinescu         b2[4] = {3./16.,10./16.,3./16.,0};
325b1c69cc3SEmil Constantinescu      ierr = TSRosWRegister(TSROSWLLSSP3P3S2C,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
326b1c69cc3SEmil Constantinescu   }
327e27a552bSJed Brown   PetscFunctionReturn(0);
328e27a552bSJed Brown }
329e27a552bSJed Brown 
330e27a552bSJed Brown #undef __FUNCT__
331e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterDestroy"
332e27a552bSJed Brown /*@C
333e27a552bSJed Brown    TSRosWRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister().
334e27a552bSJed Brown 
335e27a552bSJed Brown    Not Collective
336e27a552bSJed Brown 
337e27a552bSJed Brown    Level: advanced
338e27a552bSJed Brown 
339e27a552bSJed Brown .keywords: TSRosW, register, destroy
340e27a552bSJed Brown .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegisterDynamic()
341e27a552bSJed Brown @*/
342e27a552bSJed Brown PetscErrorCode TSRosWRegisterDestroy(void)
343e27a552bSJed Brown {
344e27a552bSJed Brown   PetscErrorCode ierr;
34561692a83SJed Brown   RosWTableauLink link;
346e27a552bSJed Brown 
347e27a552bSJed Brown   PetscFunctionBegin;
34861692a83SJed Brown   while ((link = RosWTableauList)) {
34961692a83SJed Brown     RosWTableau t = &link->tab;
35061692a83SJed Brown     RosWTableauList = link->next;
35161692a83SJed Brown     ierr = PetscFree5(t->A,t->Gamma,t->b,t->ASum,t->GammaSum);CHKERRQ(ierr);
352c17803e7SJed Brown     ierr = PetscFree4(t->At,t->bt,t->GammaInv,t->GammaZeroDiag);CHKERRQ(ierr);
353fe7e6d57SJed Brown     ierr = PetscFree2(t->bembed,t->bembedt);CHKERRQ(ierr);
354e27a552bSJed Brown     ierr = PetscFree(t->name);CHKERRQ(ierr);
355e27a552bSJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
356e27a552bSJed Brown   }
357e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_FALSE;
358e27a552bSJed Brown   PetscFunctionReturn(0);
359e27a552bSJed Brown }
360e27a552bSJed Brown 
361e27a552bSJed Brown #undef __FUNCT__
362e27a552bSJed Brown #define __FUNCT__ "TSRosWInitializePackage"
363e27a552bSJed Brown /*@C
364e27a552bSJed Brown   TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called
365e27a552bSJed Brown   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RosW()
366e27a552bSJed Brown   when using static libraries.
367e27a552bSJed Brown 
368e27a552bSJed Brown   Input Parameter:
369e27a552bSJed Brown   path - The dynamic library path, or PETSC_NULL
370e27a552bSJed Brown 
371e27a552bSJed Brown   Level: developer
372e27a552bSJed Brown 
373e27a552bSJed Brown .keywords: TS, TSRosW, initialize, package
374e27a552bSJed Brown .seealso: PetscInitialize()
375e27a552bSJed Brown @*/
376e27a552bSJed Brown PetscErrorCode TSRosWInitializePackage(const char path[])
377e27a552bSJed Brown {
378e27a552bSJed Brown   PetscErrorCode ierr;
379e27a552bSJed Brown 
380e27a552bSJed Brown   PetscFunctionBegin;
381e27a552bSJed Brown   if (TSRosWPackageInitialized) PetscFunctionReturn(0);
382e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_TRUE;
383e27a552bSJed Brown   ierr = TSRosWRegisterAll();CHKERRQ(ierr);
384e27a552bSJed Brown   ierr = PetscRegisterFinalize(TSRosWFinalizePackage);CHKERRQ(ierr);
385e27a552bSJed Brown   PetscFunctionReturn(0);
386e27a552bSJed Brown }
387e27a552bSJed Brown 
388e27a552bSJed Brown #undef __FUNCT__
389e27a552bSJed Brown #define __FUNCT__ "TSRosWFinalizePackage"
390e27a552bSJed Brown /*@C
391e27a552bSJed Brown   TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is
392e27a552bSJed Brown   called from PetscFinalize().
393e27a552bSJed Brown 
394e27a552bSJed Brown   Level: developer
395e27a552bSJed Brown 
396e27a552bSJed Brown .keywords: Petsc, destroy, package
397e27a552bSJed Brown .seealso: PetscFinalize()
398e27a552bSJed Brown @*/
399e27a552bSJed Brown PetscErrorCode TSRosWFinalizePackage(void)
400e27a552bSJed Brown {
401e27a552bSJed Brown   PetscErrorCode ierr;
402e27a552bSJed Brown 
403e27a552bSJed Brown   PetscFunctionBegin;
404e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_FALSE;
405e27a552bSJed Brown   ierr = TSRosWRegisterDestroy();CHKERRQ(ierr);
406e27a552bSJed Brown   PetscFunctionReturn(0);
407e27a552bSJed Brown }
408e27a552bSJed Brown 
409e27a552bSJed Brown #undef __FUNCT__
410e27a552bSJed Brown #define __FUNCT__ "TSRosWRegister"
411e27a552bSJed Brown /*@C
41261692a83SJed Brown    TSRosWRegister - register a Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
413e27a552bSJed Brown 
414e27a552bSJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
415e27a552bSJed Brown 
416e27a552bSJed Brown    Input Parameters:
417e27a552bSJed Brown +  name - identifier for method
418e27a552bSJed Brown .  order - approximation order of method
419e27a552bSJed Brown .  s - number of stages, this is the dimension of the matrices below
42061692a83SJed Brown .  A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular
42161692a83SJed Brown .  Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal
422fe7e6d57SJed Brown .  b - Step completion table (dimension s)
423fe7e6d57SJed Brown -  bembed - Step completion table for a scheme of order one less (dimension s, PETSC_NULL if no embedded scheme is available)
424e27a552bSJed Brown 
425e27a552bSJed Brown    Notes:
42661692a83SJed Brown    Several Rosenbrock W methods are provided, this function is only needed to create new methods.
427e27a552bSJed Brown 
428e27a552bSJed Brown    Level: advanced
429e27a552bSJed Brown 
430e27a552bSJed Brown .keywords: TS, register
431e27a552bSJed Brown 
432e27a552bSJed Brown .seealso: TSRosW
433e27a552bSJed Brown @*/
434e27a552bSJed Brown PetscErrorCode TSRosWRegister(const TSRosWType name,PetscInt order,PetscInt s,
435fe7e6d57SJed Brown                               const PetscReal A[],const PetscReal Gamma[],const PetscReal b[],const PetscReal bembed[])
436e27a552bSJed Brown {
437e27a552bSJed Brown   PetscErrorCode ierr;
43861692a83SJed Brown   RosWTableauLink link;
43961692a83SJed Brown   RosWTableau t;
44061692a83SJed Brown   PetscInt i,j,k;
44161692a83SJed Brown   PetscScalar *GammaInv;
442e27a552bSJed Brown 
443e27a552bSJed Brown   PetscFunctionBegin;
444fe7e6d57SJed Brown   PetscValidCharPointer(name,1);
445fe7e6d57SJed Brown   PetscValidPointer(A,4);
446fe7e6d57SJed Brown   PetscValidPointer(Gamma,5);
447fe7e6d57SJed Brown   PetscValidPointer(b,6);
448fe7e6d57SJed Brown   if (bembed) PetscValidPointer(bembed,7);
449fe7e6d57SJed Brown 
450e27a552bSJed Brown   ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
451e27a552bSJed Brown   ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
452e27a552bSJed Brown   t = &link->tab;
453e27a552bSJed Brown   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
454e27a552bSJed Brown   t->order = order;
455e27a552bSJed Brown   t->s = s;
45661692a83SJed 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);
457c17803e7SJed Brown   ierr = PetscMalloc4(s*s,PetscReal,&t->At,s,PetscReal,&t->bt,s*s,PetscReal,&t->GammaInv,s,PetscBool,&t->GammaZeroDiag);CHKERRQ(ierr);
458e27a552bSJed Brown   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
45961692a83SJed Brown   ierr = PetscMemcpy(t->Gamma,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr);
46061692a83SJed Brown   ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
461fe7e6d57SJed Brown   if (bembed) {
462fe7e6d57SJed Brown     ierr = PetscMalloc2(s,PetscReal,&t->bembed,s,PetscReal,&t->bembedt);CHKERRQ(ierr);
463fe7e6d57SJed Brown     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
464fe7e6d57SJed Brown   }
46561692a83SJed Brown   for (i=0; i<s; i++) {
46661692a83SJed Brown     t->ASum[i] = 0;
46761692a83SJed Brown     t->GammaSum[i] = 0;
46861692a83SJed Brown     for (j=0; j<s; j++) {
46961692a83SJed Brown       t->ASum[i] += A[i*s+j];
470fe7e6d57SJed Brown       t->GammaSum[i] += Gamma[i*s+j];
47161692a83SJed Brown     }
47261692a83SJed Brown   }
47361692a83SJed Brown   ierr = PetscMalloc(s*s*sizeof(PetscScalar),&GammaInv);CHKERRQ(ierr); /* Need to use Scalar for inverse, then convert back to Real */
47461692a83SJed Brown   for (i=0; i<s*s; i++) GammaInv[i] = Gamma[i];
475fd96d5b0SEmil Constantinescu   for (i=0; i<s; i++) {
476fd96d5b0SEmil Constantinescu     if (Gamma[i*s+i] == 0.0) {
477fd96d5b0SEmil Constantinescu       GammaInv[i*s+i] = 1.0;
478c17803e7SJed Brown       t->GammaZeroDiag[i] = PETSC_TRUE;
479fd96d5b0SEmil Constantinescu     } else {
480c17803e7SJed Brown       t->GammaZeroDiag[i] = PETSC_FALSE;
481fd96d5b0SEmil Constantinescu     }
482fd96d5b0SEmil Constantinescu   }
483fd96d5b0SEmil Constantinescu 
48461692a83SJed Brown   switch (s) {
48561692a83SJed Brown   case 1: GammaInv[0] = 1./GammaInv[0]; break;
48661692a83SJed Brown   case 2: ierr = Kernel_A_gets_inverse_A_2(GammaInv,0);CHKERRQ(ierr); break;
48761692a83SJed Brown   case 3: ierr = Kernel_A_gets_inverse_A_3(GammaInv,0);CHKERRQ(ierr); break;
48861692a83SJed Brown   case 4: ierr = Kernel_A_gets_inverse_A_4(GammaInv,0);CHKERRQ(ierr); break;
48961692a83SJed Brown   case 5: {
49061692a83SJed Brown     PetscInt ipvt5[5];
49161692a83SJed Brown     MatScalar work5[5*5];
49261692a83SJed Brown     ierr = Kernel_A_gets_inverse_A_5(GammaInv,ipvt5,work5,0);CHKERRQ(ierr); break;
49361692a83SJed Brown   }
49461692a83SJed Brown   case 6: ierr = Kernel_A_gets_inverse_A_6(GammaInv,0);CHKERRQ(ierr); break;
49561692a83SJed Brown   case 7: ierr = Kernel_A_gets_inverse_A_7(GammaInv,0);CHKERRQ(ierr); break;
49661692a83SJed Brown   default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for %D stages",s);
49761692a83SJed Brown   }
49861692a83SJed Brown   for (i=0; i<s*s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]);
49961692a83SJed Brown   ierr = PetscFree(GammaInv);CHKERRQ(ierr);
50061692a83SJed Brown   for (i=0; i<s; i++) {
50161692a83SJed Brown     for (j=0; j<s; j++) {
50261692a83SJed Brown       t->At[i*s+j] = 0;
50361692a83SJed Brown       for (k=0; k<s; k++) {
50461692a83SJed Brown         t->At[i*s+j] += t->A[i*s+k] * t->GammaInv[k*s+j];
50561692a83SJed Brown       }
50661692a83SJed Brown     }
50761692a83SJed Brown     t->bt[i] = 0;
50861692a83SJed Brown     for (j=0; j<s; j++) {
50961692a83SJed Brown       t->bt[i] += t->b[j] * t->GammaInv[j*s+i];
51061692a83SJed Brown     }
511fe7e6d57SJed Brown     if (bembed) {
512fe7e6d57SJed Brown       t->bembedt[i] = 0;
513fe7e6d57SJed Brown       for (j=0; j<s; j++) {
514fe7e6d57SJed Brown         t->bembedt[i] += t->bembed[j] * t->GammaInv[j*s+i];
515fe7e6d57SJed Brown       }
516fe7e6d57SJed Brown     }
51761692a83SJed Brown   }
5188d59e960SJed Brown   t->ccfl = 1.0;                /* Fix this */
5198d59e960SJed Brown 
52061692a83SJed Brown   link->next = RosWTableauList;
52161692a83SJed Brown   RosWTableauList = link;
522e27a552bSJed Brown   PetscFunctionReturn(0);
523e27a552bSJed Brown }
524e27a552bSJed Brown 
525e27a552bSJed Brown #undef __FUNCT__
5261c3436cfSJed Brown #define __FUNCT__ "TSEvaluateStep_RosW"
5271c3436cfSJed Brown /*
5281c3436cfSJed Brown  The step completion formula is
5291c3436cfSJed Brown 
5301c3436cfSJed Brown  x1 = x0 + b^T Y
5311c3436cfSJed Brown 
5321c3436cfSJed Brown  where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been
5331c3436cfSJed Brown  updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write
5341c3436cfSJed Brown 
5351c3436cfSJed Brown  x1e = x0 + be^T Y
5361c3436cfSJed Brown      = x1 - b^T Y + be^T Y
5371c3436cfSJed Brown      = x1 + (be - b)^T Y
5381c3436cfSJed Brown 
5391c3436cfSJed Brown  so we can evaluate the method of different order even after the step has been optimistically completed.
5401c3436cfSJed Brown */
5411c3436cfSJed Brown static PetscErrorCode TSEvaluateStep_RosW(TS ts,PetscInt order,Vec X,PetscBool *done)
5421c3436cfSJed Brown {
5431c3436cfSJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
5441c3436cfSJed Brown   RosWTableau    tab  = ros->tableau;
5451c3436cfSJed Brown   PetscScalar    *w = ros->work;
5461c3436cfSJed Brown   PetscInt       i;
5471c3436cfSJed Brown   PetscErrorCode ierr;
5481c3436cfSJed Brown 
5491c3436cfSJed Brown   PetscFunctionBegin;
5501c3436cfSJed Brown   if (order == tab->order) {
5511c3436cfSJed Brown     if (ros->step_taken) {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
5521c3436cfSJed Brown     else {
5531c3436cfSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
5541c3436cfSJed Brown       ierr = VecMAXPY(X,tab->s,tab->bt,ros->Y);CHKERRQ(ierr);
5551c3436cfSJed Brown     }
5561c3436cfSJed Brown     if (done) *done = PETSC_TRUE;
5571c3436cfSJed Brown     PetscFunctionReturn(0);
5581c3436cfSJed Brown   } else if (order == tab->order-1) {
5591c3436cfSJed Brown     if (!tab->bembedt) goto unavailable;
5601c3436cfSJed Brown     if (ros->step_taken) {
5611c3436cfSJed Brown       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i];
5621c3436cfSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
5631c3436cfSJed Brown       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
5641c3436cfSJed Brown     } else {
5651c3436cfSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
5661c3436cfSJed Brown       ierr = VecMAXPY(X,tab->s,tab->bembedt,ros->Y);CHKERRQ(ierr);
5671c3436cfSJed Brown     }
5681c3436cfSJed Brown     if (done) *done = PETSC_TRUE;
5691c3436cfSJed Brown     PetscFunctionReturn(0);
5701c3436cfSJed Brown   }
5711c3436cfSJed Brown   unavailable:
5721c3436cfSJed Brown   if (done) *done = PETSC_FALSE;
5731c3436cfSJed 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);
5741c3436cfSJed Brown   PetscFunctionReturn(0);
5751c3436cfSJed Brown }
5761c3436cfSJed Brown 
5771c3436cfSJed Brown #undef __FUNCT__
578e27a552bSJed Brown #define __FUNCT__ "TSStep_RosW"
579e27a552bSJed Brown static PetscErrorCode TSStep_RosW(TS ts)
580e27a552bSJed Brown {
58161692a83SJed Brown   TS_RosW         *ros = (TS_RosW*)ts->data;
58261692a83SJed Brown   RosWTableau     tab  = ros->tableau;
583e27a552bSJed Brown   const PetscInt  s    = tab->s;
5841c3436cfSJed Brown   const PetscReal *At  = tab->At,*Gamma = tab->Gamma,*ASum = tab->ASum,*GammaInv = tab->GammaInv;
585c17803e7SJed Brown   const PetscBool *GammaZeroDiag = tab->GammaZeroDiag;
58661692a83SJed Brown   PetscScalar     *w   = ros->work;
58761692a83SJed Brown   Vec             *Y   = ros->Y,Zdot = ros->Zdot,Zstage = ros->Zstage;
588e27a552bSJed Brown   SNES            snes;
5891c3436cfSJed Brown   TSAdapt         adapt;
5901c3436cfSJed Brown   PetscInt        i,j,its,lits,reject,next_scheme;
591cdbf8f93SLisandro Dalcin   PetscReal       next_time_step;
5921c3436cfSJed Brown   PetscBool       accept;
593e27a552bSJed Brown   PetscErrorCode  ierr;
594e27a552bSJed Brown 
595e27a552bSJed Brown   PetscFunctionBegin;
596e27a552bSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
597cdbf8f93SLisandro Dalcin   next_time_step = ts->time_step;
5981c3436cfSJed Brown   accept = PETSC_TRUE;
5991c3436cfSJed Brown   ros->step_taken = PETSC_FALSE;
600e27a552bSJed Brown 
6011c3436cfSJed Brown   for (reject=0; reject<ts->max_reject; reject++,ts->reject++) {
6021c3436cfSJed Brown     const PetscReal h = ts->time_step;
603e27a552bSJed Brown     for (i=0; i<s; i++) {
6041c3436cfSJed Brown       ros->stage_time = ts->ptime + h*ASum[i];
605c17803e7SJed Brown       if (GammaZeroDiag[i]) {
606c17803e7SJed Brown         ros->stage_explicit = PETSC_TRUE;
607fd96d5b0SEmil Constantinescu         ros->shift = 1./h;
608c17803e7SJed Brown       } else {
609c17803e7SJed Brown         ros->stage_explicit = PETSC_FALSE;
61061692a83SJed Brown         ros->shift = 1./(h*Gamma[i*s+i]);
611fd96d5b0SEmil Constantinescu       }
61261692a83SJed Brown 
61361692a83SJed Brown       ierr = VecCopy(ts->vec_sol,Zstage);CHKERRQ(ierr);
61461692a83SJed Brown       ierr = VecMAXPY(Zstage,i,&At[i*s+0],Y);CHKERRQ(ierr);
61561692a83SJed Brown 
61661692a83SJed Brown       for (j=0; j<i; j++) w[j] = 1./h * GammaInv[i*s+j];
61761692a83SJed Brown       ierr = VecZeroEntries(Zdot);CHKERRQ(ierr);
61861692a83SJed Brown       ierr = VecMAXPY(Zdot,i,w,Y);CHKERRQ(ierr);
61961692a83SJed Brown 
620e27a552bSJed Brown       /* Initial guess taken from last stage */
62161692a83SJed Brown       ierr = VecZeroEntries(Y[i]);CHKERRQ(ierr);
62261692a83SJed Brown 
62361692a83SJed Brown       if (!ros->recompute_jacobian && !i) {
62461692a83SJed Brown         ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */
62561692a83SJed Brown       }
62661692a83SJed Brown 
62761692a83SJed Brown       ierr = SNESSolve(snes,PETSC_NULL,Y[i]);CHKERRQ(ierr);
628e27a552bSJed Brown       ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
629e27a552bSJed Brown       ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
630e27a552bSJed Brown       ts->nonlinear_its += its; ts->linear_its += lits;
631e27a552bSJed Brown     }
6321c3436cfSJed Brown     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
6331c3436cfSJed Brown     ros->step_taken = PETSC_TRUE;
634e27a552bSJed Brown 
6351c3436cfSJed Brown     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
6361c3436cfSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
6371c3436cfSJed Brown     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
6388d59e960SJed Brown     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
6391c3436cfSJed Brown     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
6401c3436cfSJed Brown     if (accept) {
6411c3436cfSJed Brown       /* ignore next_scheme for now */
642e27a552bSJed Brown       ts->ptime += ts->time_step;
643cdbf8f93SLisandro Dalcin       ts->time_step = next_time_step;
644e27a552bSJed Brown       ts->steps++;
6451c3436cfSJed Brown       break;
6461c3436cfSJed Brown     } else {                    /* Roll back the current step */
6471c3436cfSJed Brown       for (i=0; i<s; i++) w[i] = -tab->bt[i];
6481c3436cfSJed Brown       ierr = VecMAXPY(ts->vec_sol,s,w,Y);CHKERRQ(ierr);
6491c3436cfSJed Brown       ts->time_step = next_time_step;
6501c3436cfSJed Brown       ros->step_taken = PETSC_FALSE;
6511c3436cfSJed Brown     }
6521c3436cfSJed Brown   }
6531c3436cfSJed Brown 
654e27a552bSJed Brown   PetscFunctionReturn(0);
655e27a552bSJed Brown }
656e27a552bSJed Brown 
657e27a552bSJed Brown #undef __FUNCT__
658e27a552bSJed Brown #define __FUNCT__ "TSInterpolate_RosW"
659e27a552bSJed Brown static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec X)
660e27a552bSJed Brown {
66161692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
662e27a552bSJed Brown 
663e27a552bSJed Brown   PetscFunctionBegin;
66461692a83SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name);
665e27a552bSJed Brown   PetscFunctionReturn(0);
666e27a552bSJed Brown }
667e27a552bSJed Brown 
668e27a552bSJed Brown /*------------------------------------------------------------*/
669e27a552bSJed Brown #undef __FUNCT__
670e27a552bSJed Brown #define __FUNCT__ "TSReset_RosW"
671e27a552bSJed Brown static PetscErrorCode TSReset_RosW(TS ts)
672e27a552bSJed Brown {
67361692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
674e27a552bSJed Brown   PetscInt       s;
675e27a552bSJed Brown   PetscErrorCode ierr;
676e27a552bSJed Brown 
677e27a552bSJed Brown   PetscFunctionBegin;
67861692a83SJed Brown   if (!ros->tableau) PetscFunctionReturn(0);
67961692a83SJed Brown    s = ros->tableau->s;
68061692a83SJed Brown   ierr = VecDestroyVecs(s,&ros->Y);CHKERRQ(ierr);
68161692a83SJed Brown   ierr = VecDestroy(&ros->Ydot);CHKERRQ(ierr);
68261692a83SJed Brown   ierr = VecDestroy(&ros->Ystage);CHKERRQ(ierr);
68361692a83SJed Brown   ierr = VecDestroy(&ros->Zdot);CHKERRQ(ierr);
68461692a83SJed Brown   ierr = VecDestroy(&ros->Zstage);CHKERRQ(ierr);
68561692a83SJed Brown   ierr = PetscFree(ros->work);CHKERRQ(ierr);
686e27a552bSJed Brown   PetscFunctionReturn(0);
687e27a552bSJed Brown }
688e27a552bSJed Brown 
689e27a552bSJed Brown #undef __FUNCT__
690e27a552bSJed Brown #define __FUNCT__ "TSDestroy_RosW"
691e27a552bSJed Brown static PetscErrorCode TSDestroy_RosW(TS ts)
692e27a552bSJed Brown {
693e27a552bSJed Brown   PetscErrorCode  ierr;
694e27a552bSJed Brown 
695e27a552bSJed Brown   PetscFunctionBegin;
696e27a552bSJed Brown   ierr = TSReset_RosW(ts);CHKERRQ(ierr);
697e27a552bSJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
698e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","",PETSC_NULL);CHKERRQ(ierr);
699e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","",PETSC_NULL);CHKERRQ(ierr);
70061692a83SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","",PETSC_NULL);CHKERRQ(ierr);
701e27a552bSJed Brown   PetscFunctionReturn(0);
702e27a552bSJed Brown }
703e27a552bSJed Brown 
704e27a552bSJed Brown /*
705e27a552bSJed Brown   This defines the nonlinear equation that is to be solved with SNES
706e27a552bSJed Brown   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
707e27a552bSJed Brown */
708e27a552bSJed Brown #undef __FUNCT__
709e27a552bSJed Brown #define __FUNCT__ "SNESTSFormFunction_RosW"
710e27a552bSJed Brown static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec X,Vec F,TS ts)
711e27a552bSJed Brown {
71261692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
713e27a552bSJed Brown   PetscErrorCode ierr;
714e27a552bSJed Brown 
715e27a552bSJed Brown   PetscFunctionBegin;
716c17803e7SJed Brown   if (ros->stage_explicit) {
717c17803e7SJed Brown     ierr = VecAXPBY(ros->Ydot,ros->shift,0.0,X);CHKERRQ(ierr);       /* Ydot = shift*X*/
718c17803e7SJed Brown   } else {
71961692a83SJed Brown     ierr = VecWAXPY(ros->Ydot,ros->shift,X,ros->Zdot);CHKERRQ(ierr); /* Ydot = shift*X + Zdot */
720c17803e7SJed Brown   }
72161692a83SJed Brown   ierr = VecWAXPY(ros->Ystage,1.0,X,ros->Zstage);CHKERRQ(ierr);    /* Ystage = X + Zstage */
72261692a83SJed Brown   ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,ros->Ydot,F,PETSC_FALSE);CHKERRQ(ierr);
723e27a552bSJed Brown   PetscFunctionReturn(0);
724e27a552bSJed Brown }
725e27a552bSJed Brown 
726e27a552bSJed Brown #undef __FUNCT__
727e27a552bSJed Brown #define __FUNCT__ "SNESTSFormJacobian_RosW"
728e27a552bSJed Brown static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
729e27a552bSJed Brown {
73061692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
731e27a552bSJed Brown   PetscErrorCode ierr;
732e27a552bSJed Brown 
733e27a552bSJed Brown   PetscFunctionBegin;
73461692a83SJed Brown   /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
73561692a83SJed Brown   ierr = TSComputeIJacobian(ts,ros->stage_time,ros->Ystage,ros->Ydot,ros->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr);
736e27a552bSJed Brown   PetscFunctionReturn(0);
737e27a552bSJed Brown }
738e27a552bSJed Brown 
739e27a552bSJed Brown #undef __FUNCT__
740e27a552bSJed Brown #define __FUNCT__ "TSSetUp_RosW"
741e27a552bSJed Brown static PetscErrorCode TSSetUp_RosW(TS ts)
742e27a552bSJed Brown {
74361692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
74461692a83SJed Brown   RosWTableau    tab  = ros->tableau;
745e27a552bSJed Brown   PetscInt       s    = tab->s;
746e27a552bSJed Brown   PetscErrorCode ierr;
747e27a552bSJed Brown 
748e27a552bSJed Brown   PetscFunctionBegin;
74961692a83SJed Brown   if (!ros->tableau) {
750e27a552bSJed Brown     ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);
751e27a552bSJed Brown   }
75261692a83SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ros->Y);CHKERRQ(ierr);
75361692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Ydot);CHKERRQ(ierr);
75461692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Ystage);CHKERRQ(ierr);
75561692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Zdot);CHKERRQ(ierr);
75661692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Zstage);CHKERRQ(ierr);
75761692a83SJed Brown   ierr = PetscMalloc(s*sizeof(ros->work[0]),&ros->work);CHKERRQ(ierr);
758e27a552bSJed Brown   PetscFunctionReturn(0);
759e27a552bSJed Brown }
760e27a552bSJed Brown /*------------------------------------------------------------*/
761e27a552bSJed Brown 
762e27a552bSJed Brown #undef __FUNCT__
763e27a552bSJed Brown #define __FUNCT__ "TSSetFromOptions_RosW"
764e27a552bSJed Brown static PetscErrorCode TSSetFromOptions_RosW(TS ts)
765e27a552bSJed Brown {
76661692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
767e27a552bSJed Brown   PetscErrorCode ierr;
76861692a83SJed Brown   char           rostype[256];
769e27a552bSJed Brown 
770e27a552bSJed Brown   PetscFunctionBegin;
771e27a552bSJed Brown   ierr = PetscOptionsHead("RosW ODE solver options");CHKERRQ(ierr);
772e27a552bSJed Brown   {
77361692a83SJed Brown     RosWTableauLink link;
774e27a552bSJed Brown     PetscInt count,choice;
775e27a552bSJed Brown     PetscBool flg;
776e27a552bSJed Brown     const char **namelist;
77761692a83SJed Brown     SNES snes;
77861692a83SJed Brown 
77961692a83SJed Brown     ierr = PetscStrncpy(rostype,TSRosWDefault,sizeof rostype);CHKERRQ(ierr);
78061692a83SJed Brown     for (link=RosWTableauList,count=0; link; link=link->next,count++) ;
781e27a552bSJed Brown     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
78261692a83SJed Brown     for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
78361692a83SJed Brown     ierr = PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,rostype,&choice,&flg);CHKERRQ(ierr);
78461692a83SJed Brown     ierr = TSRosWSetType(ts,flg ? namelist[choice] : rostype);CHKERRQ(ierr);
785e27a552bSJed Brown     ierr = PetscFree(namelist);CHKERRQ(ierr);
78661692a83SJed Brown 
78761692a83SJed Brown     ierr = PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,PETSC_NULL);CHKERRQ(ierr);
78861692a83SJed Brown 
78961692a83SJed Brown     /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
79061692a83SJed Brown     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
79161692a83SJed Brown     if (!((PetscObject)snes)->type_name) {
79261692a83SJed Brown       ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
79361692a83SJed Brown     }
79461692a83SJed Brown     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
795e27a552bSJed Brown   }
796e27a552bSJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
797e27a552bSJed Brown   PetscFunctionReturn(0);
798e27a552bSJed Brown }
799e27a552bSJed Brown 
800e27a552bSJed Brown #undef __FUNCT__
801e27a552bSJed Brown #define __FUNCT__ "PetscFormatRealArray"
802e27a552bSJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
803e27a552bSJed Brown {
804e27a552bSJed Brown   PetscErrorCode ierr;
805e408995aSJed Brown   PetscInt i;
806e408995aSJed Brown   size_t left,count;
807e27a552bSJed Brown   char *p;
808e27a552bSJed Brown 
809e27a552bSJed Brown   PetscFunctionBegin;
810e408995aSJed Brown   for (i=0,p=buf,left=len; i<n; i++) {
811e408995aSJed Brown     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
812e27a552bSJed Brown     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
813e27a552bSJed Brown     left -= count;
814e27a552bSJed Brown     p += count;
815e27a552bSJed Brown     *p++ = ' ';
816e27a552bSJed Brown   }
817e27a552bSJed Brown   p[i ? 0 : -1] = 0;
818e27a552bSJed Brown   PetscFunctionReturn(0);
819e27a552bSJed Brown }
820e27a552bSJed Brown 
821e27a552bSJed Brown #undef __FUNCT__
822e27a552bSJed Brown #define __FUNCT__ "TSView_RosW"
823e27a552bSJed Brown static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer)
824e27a552bSJed Brown {
82561692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
82661692a83SJed Brown   RosWTableau    tab  = ros->tableau;
827e27a552bSJed Brown   PetscBool      iascii;
828e27a552bSJed Brown   PetscErrorCode ierr;
829e27a552bSJed Brown 
830e27a552bSJed Brown   PetscFunctionBegin;
831e27a552bSJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
832e27a552bSJed Brown   if (iascii) {
83361692a83SJed Brown     const TSRosWType rostype;
834e408995aSJed Brown     PetscInt i;
835e408995aSJed Brown     PetscReal abscissa[512];
836e27a552bSJed Brown     char buf[512];
83761692a83SJed Brown     ierr = TSRosWGetType(ts,&rostype);CHKERRQ(ierr);
83861692a83SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Rosenbrock-W %s\n",rostype);CHKERRQ(ierr);
839e408995aSJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ASum);CHKERRQ(ierr);
84061692a83SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A       = %s\n",buf);CHKERRQ(ierr);
841e408995aSJed Brown     for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i];
842e408995aSJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,abscissa);CHKERRQ(ierr);
843e408995aSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A+Gamma = %s\n",buf);CHKERRQ(ierr);
844e27a552bSJed Brown   }
845e27a552bSJed Brown   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
846e27a552bSJed Brown   PetscFunctionReturn(0);
847e27a552bSJed Brown }
848e27a552bSJed Brown 
849e27a552bSJed Brown #undef __FUNCT__
850e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType"
851e27a552bSJed Brown /*@C
85261692a83SJed Brown   TSRosWSetType - Set the type of Rosenbrock-W scheme
853e27a552bSJed Brown 
854e27a552bSJed Brown   Logically collective
855e27a552bSJed Brown 
856e27a552bSJed Brown   Input Parameter:
857e27a552bSJed Brown +  ts - timestepping context
85861692a83SJed Brown -  rostype - type of Rosenbrock-W scheme
859e27a552bSJed Brown 
860e27a552bSJed Brown   Level: intermediate
861e27a552bSJed Brown 
862e27a552bSJed Brown .seealso: TSRosWGetType()
863e27a552bSJed Brown @*/
86461692a83SJed Brown PetscErrorCode TSRosWSetType(TS ts,const TSRosWType rostype)
865e27a552bSJed Brown {
866e27a552bSJed Brown   PetscErrorCode ierr;
867e27a552bSJed Brown 
868e27a552bSJed Brown   PetscFunctionBegin;
869e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
87061692a83SJed Brown   ierr = PetscTryMethod(ts,"TSRosWSetType_C",(TS,const TSRosWType),(ts,rostype));CHKERRQ(ierr);
871e27a552bSJed Brown   PetscFunctionReturn(0);
872e27a552bSJed Brown }
873e27a552bSJed Brown 
874e27a552bSJed Brown #undef __FUNCT__
875e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType"
876e27a552bSJed Brown /*@C
87761692a83SJed Brown   TSRosWGetType - Get the type of Rosenbrock-W scheme
878e27a552bSJed Brown 
879e27a552bSJed Brown   Logically collective
880e27a552bSJed Brown 
881e27a552bSJed Brown   Input Parameter:
882e27a552bSJed Brown .  ts - timestepping context
883e27a552bSJed Brown 
884e27a552bSJed Brown   Output Parameter:
88561692a83SJed Brown .  rostype - type of Rosenbrock-W scheme
886e27a552bSJed Brown 
887e27a552bSJed Brown   Level: intermediate
888e27a552bSJed Brown 
889e27a552bSJed Brown .seealso: TSRosWGetType()
890e27a552bSJed Brown @*/
89161692a83SJed Brown PetscErrorCode TSRosWGetType(TS ts,const TSRosWType *rostype)
892e27a552bSJed Brown {
893e27a552bSJed Brown   PetscErrorCode ierr;
894e27a552bSJed Brown 
895e27a552bSJed Brown   PetscFunctionBegin;
896e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
89761692a83SJed Brown   ierr = PetscUseMethod(ts,"TSRosWGetType_C",(TS,const TSRosWType*),(ts,rostype));CHKERRQ(ierr);
898e27a552bSJed Brown   PetscFunctionReturn(0);
899e27a552bSJed Brown }
900e27a552bSJed Brown 
901e27a552bSJed Brown #undef __FUNCT__
90261692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian"
903e27a552bSJed Brown /*@C
90461692a83SJed Brown   TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step.
905e27a552bSJed Brown 
906e27a552bSJed Brown   Logically collective
907e27a552bSJed Brown 
908e27a552bSJed Brown   Input Parameter:
909e27a552bSJed Brown +  ts - timestepping context
91061692a83SJed Brown -  flg - PETSC_TRUE to recompute the Jacobian at each stage
911e27a552bSJed Brown 
912e27a552bSJed Brown   Level: intermediate
913e27a552bSJed Brown 
914e27a552bSJed Brown .seealso: TSRosWGetType()
915e27a552bSJed Brown @*/
91661692a83SJed Brown PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg)
917e27a552bSJed Brown {
918e27a552bSJed Brown   PetscErrorCode ierr;
919e27a552bSJed Brown 
920e27a552bSJed Brown   PetscFunctionBegin;
921e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
92261692a83SJed Brown   ierr = PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
923e27a552bSJed Brown   PetscFunctionReturn(0);
924e27a552bSJed Brown }
925e27a552bSJed Brown 
926e27a552bSJed Brown EXTERN_C_BEGIN
927e27a552bSJed Brown #undef __FUNCT__
928e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType_RosW"
92961692a83SJed Brown PetscErrorCode  TSRosWGetType_RosW(TS ts,const TSRosWType *rostype)
930e27a552bSJed Brown {
93161692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
932e27a552bSJed Brown   PetscErrorCode ierr;
933e27a552bSJed Brown 
934e27a552bSJed Brown   PetscFunctionBegin;
93561692a83SJed Brown   if (!ros->tableau) {ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);}
93661692a83SJed Brown   *rostype = ros->tableau->name;
937e27a552bSJed Brown   PetscFunctionReturn(0);
938e27a552bSJed Brown }
939e27a552bSJed Brown #undef __FUNCT__
940e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType_RosW"
94161692a83SJed Brown PetscErrorCode  TSRosWSetType_RosW(TS ts,const TSRosWType rostype)
942e27a552bSJed Brown {
94361692a83SJed Brown   TS_RosW         *ros = (TS_RosW*)ts->data;
944e27a552bSJed Brown   PetscErrorCode  ierr;
945e27a552bSJed Brown   PetscBool       match;
94661692a83SJed Brown   RosWTableauLink link;
947e27a552bSJed Brown 
948e27a552bSJed Brown   PetscFunctionBegin;
94961692a83SJed Brown   if (ros->tableau) {
95061692a83SJed Brown     ierr = PetscStrcmp(ros->tableau->name,rostype,&match);CHKERRQ(ierr);
951e27a552bSJed Brown     if (match) PetscFunctionReturn(0);
952e27a552bSJed Brown   }
95361692a83SJed Brown   for (link = RosWTableauList; link; link=link->next) {
95461692a83SJed Brown     ierr = PetscStrcmp(link->tab.name,rostype,&match);CHKERRQ(ierr);
955e27a552bSJed Brown     if (match) {
956e27a552bSJed Brown       ierr = TSReset_RosW(ts);CHKERRQ(ierr);
95761692a83SJed Brown       ros->tableau = &link->tab;
958e27a552bSJed Brown       PetscFunctionReturn(0);
959e27a552bSJed Brown     }
960e27a552bSJed Brown   }
96161692a83SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype);
962e27a552bSJed Brown   PetscFunctionReturn(0);
963e27a552bSJed Brown }
96461692a83SJed Brown 
965e27a552bSJed Brown #undef __FUNCT__
96661692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian_RosW"
96761692a83SJed Brown PetscErrorCode  TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg)
968e27a552bSJed Brown {
96961692a83SJed Brown   TS_RosW *ros = (TS_RosW*)ts->data;
970e27a552bSJed Brown 
971e27a552bSJed Brown   PetscFunctionBegin;
97261692a83SJed Brown   ros->recompute_jacobian = flg;
973e27a552bSJed Brown   PetscFunctionReturn(0);
974e27a552bSJed Brown }
975e27a552bSJed Brown EXTERN_C_END
976e27a552bSJed Brown 
977e27a552bSJed Brown /* ------------------------------------------------------------ */
978e27a552bSJed Brown /*MC
979e27a552bSJed Brown       TSRosW - ODE solver using Rosenbrock-W schemes
980e27a552bSJed Brown 
981e27a552bSJed Brown   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
982e27a552bSJed Brown   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
983e27a552bSJed Brown   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
984e27a552bSJed Brown 
985e27a552bSJed Brown   Notes:
98661692a83SJed Brown   This method currently only works with autonomous ODE and DAE.
98761692a83SJed Brown 
98861692a83SJed Brown   Developer notes:
98961692a83SJed Brown   Rosenbrock-W methods are typically specified for autonomous ODE
99061692a83SJed Brown 
99161692a83SJed Brown $  xdot = f(x)
99261692a83SJed Brown 
99361692a83SJed Brown   by the stage equations
99461692a83SJed Brown 
99561692a83SJed Brown $  k_i = h f(x_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j
99661692a83SJed Brown 
99761692a83SJed Brown   and step completion formula
99861692a83SJed Brown 
99961692a83SJed Brown $  x_1 = x_0 + sum_j b_j k_j
100061692a83SJed Brown 
100161692a83SJed Brown   with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(x)
100261692a83SJed Brown   and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner,
100361692a83SJed Brown   we define new variables for the stage equations
100461692a83SJed Brown 
100561692a83SJed Brown $  y_i = gamma_ij k_j
100661692a83SJed Brown 
100761692a83SJed Brown   The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define
100861692a83SJed Brown 
100961692a83SJed Brown $  A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-i}
101061692a83SJed Brown 
101161692a83SJed Brown   to rewrite the method as
101261692a83SJed Brown 
101361692a83SJed 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
101461692a83SJed Brown $  x_1 = x_0 + sum_j bt_j y_j
101561692a83SJed Brown 
101661692a83SJed Brown    where we have introduced the mass matrix M. Continue by defining
101761692a83SJed Brown 
101861692a83SJed Brown $  ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j
101961692a83SJed Brown 
102061692a83SJed Brown    or, more compactly in tensor notation
102161692a83SJed Brown 
102261692a83SJed Brown $  Ydot = 1/h (Gamma^{-1} \otimes I) Y .
102361692a83SJed Brown 
102461692a83SJed Brown    Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current
102561692a83SJed Brown    stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the
102661692a83SJed Brown    equation
102761692a83SJed Brown 
102861692a83SJed Brown $  g(x_0 + sum_j a_ij y_j + y_i, ydot_i) = 0
102961692a83SJed Brown 
103061692a83SJed Brown    with initial guess y_i = 0.
1031e27a552bSJed Brown 
1032e27a552bSJed Brown   Level: beginner
1033e27a552bSJed Brown 
1034e27a552bSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSRosWRegister()
1035e27a552bSJed Brown 
1036e27a552bSJed Brown M*/
1037e27a552bSJed Brown EXTERN_C_BEGIN
1038e27a552bSJed Brown #undef __FUNCT__
1039e27a552bSJed Brown #define __FUNCT__ "TSCreate_RosW"
1040e27a552bSJed Brown PetscErrorCode  TSCreate_RosW(TS ts)
1041e27a552bSJed Brown {
104261692a83SJed Brown   TS_RosW        *ros;
1043e27a552bSJed Brown   PetscErrorCode ierr;
1044e27a552bSJed Brown 
1045e27a552bSJed Brown   PetscFunctionBegin;
1046e27a552bSJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
1047e27a552bSJed Brown   ierr = TSRosWInitializePackage(PETSC_NULL);CHKERRQ(ierr);
1048e27a552bSJed Brown #endif
1049e27a552bSJed Brown 
1050e27a552bSJed Brown   ts->ops->reset          = TSReset_RosW;
1051e27a552bSJed Brown   ts->ops->destroy        = TSDestroy_RosW;
1052e27a552bSJed Brown   ts->ops->view           = TSView_RosW;
1053e27a552bSJed Brown   ts->ops->setup          = TSSetUp_RosW;
1054e27a552bSJed Brown   ts->ops->step           = TSStep_RosW;
1055e27a552bSJed Brown   ts->ops->interpolate    = TSInterpolate_RosW;
10561c3436cfSJed Brown   ts->ops->evaluatestep   = TSEvaluateStep_RosW;
1057e27a552bSJed Brown   ts->ops->setfromoptions = TSSetFromOptions_RosW;
1058e27a552bSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_RosW;
1059e27a552bSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_RosW;
1060e27a552bSJed Brown 
106161692a83SJed Brown   ierr = PetscNewLog(ts,TS_RosW,&ros);CHKERRQ(ierr);
106261692a83SJed Brown   ts->data = (void*)ros;
1063e27a552bSJed Brown 
1064e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","TSRosWGetType_RosW",TSRosWGetType_RosW);CHKERRQ(ierr);
1065e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","TSRosWSetType_RosW",TSRosWSetType_RosW);CHKERRQ(ierr);
106661692a83SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","TSRosWSetRecomputeJacobian_RosW",TSRosWSetRecomputeJacobian_RosW);CHKERRQ(ierr);
1067e27a552bSJed Brown   PetscFunctionReturn(0);
1068e27a552bSJed Brown }
1069e27a552bSJed Brown EXTERN_C_END
1070