xref: /petsc/src/ts/impls/rosw/rosw.c (revision b1c69cc3815f73ced67073ff25ec5313da0a61d1)
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 
144e27a552bSJed Brown #undef __FUNCT__
145e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterAll"
146e27a552bSJed Brown /*@C
147e27a552bSJed Brown   TSRosWRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSRosW
148e27a552bSJed Brown 
149e27a552bSJed Brown   Not Collective, but should be called by all processes which will need the schemes to be registered
150e27a552bSJed Brown 
151e27a552bSJed Brown   Level: advanced
152e27a552bSJed Brown 
153e27a552bSJed Brown .keywords: TS, TSRosW, register, all
154e27a552bSJed Brown 
155e27a552bSJed Brown .seealso:  TSRosWRegisterDestroy()
156e27a552bSJed Brown @*/
157e27a552bSJed Brown PetscErrorCode TSRosWRegisterAll(void)
158e27a552bSJed Brown {
159e27a552bSJed Brown   PetscErrorCode ierr;
160e27a552bSJed Brown 
161e27a552bSJed Brown   PetscFunctionBegin;
162e27a552bSJed Brown   if (TSRosWRegisterAllCalled) PetscFunctionReturn(0);
163e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_TRUE;
164e27a552bSJed Brown 
165e27a552bSJed Brown   {
16661692a83SJed Brown     const PetscReal g = 1. + 1./PetscSqrtReal(2.0);
167e27a552bSJed Brown     const PetscReal
16861692a83SJed Brown       A[2][2] = {{0,0}, {1.,0}},
16961692a83SJed Brown       Gamma[2][2] = {{g,0}, {-2.*g,g}},
1701c3436cfSJed Brown       b[2] = {0.5,0.5},
1711c3436cfSJed Brown       b1[2] = {1.0,0.0};
1721c3436cfSJed Brown     ierr = TSRosWRegister(TSROSW2P,2,2,&A[0][0],&Gamma[0][0],b,b1);CHKERRQ(ierr);
173e27a552bSJed Brown   }
174e27a552bSJed Brown   {
17561692a83SJed Brown     const PetscReal g = 1. - 1./PetscSqrtReal(2.0);
176e27a552bSJed Brown     const PetscReal
17761692a83SJed Brown       A[2][2] = {{0,0}, {1.,0}},
17861692a83SJed Brown       Gamma[2][2] = {{g,0}, {-2.*g,g}},
1791c3436cfSJed Brown       b[2] = {0.5,0.5},
1801c3436cfSJed Brown       b1[2] = {1.0,0.0};
1811c3436cfSJed Brown     ierr = TSRosWRegister(TSROSW2M,2,2,&A[0][0],&Gamma[0][0],b,b1);CHKERRQ(ierr);
182fe7e6d57SJed Brown   }
183fe7e6d57SJed Brown   {
184fe7e6d57SJed Brown     const PetscReal g = 7.8867513459481287e-01;
185fe7e6d57SJed Brown     const PetscReal
186fe7e6d57SJed Brown       A[3][3] = {{0,0,0},
187fe7e6d57SJed Brown                  {1.5773502691896257e+00,0,0},
188fe7e6d57SJed Brown                  {0.5,0,0}},
189fe7e6d57SJed Brown       Gamma[3][3] = {{g,0,0},
190fe7e6d57SJed Brown                      {-1.5773502691896257e+00,g,0},
191fe7e6d57SJed Brown                      {-6.7075317547305480e-01,1.7075317547305482e-01,g}},
192fe7e6d57SJed Brown       b[3] = {1.0566243270259355e-01,4.9038105676657971e-02,8.4529946162074843e-01},
193fe7e6d57SJed Brown       b2[3] = {-1.7863279495408180e-01,1./3.,8.4529946162074843e-01};
194fe7e6d57SJed Brown     ierr = TSRosWRegister(TSROSWRA3PW,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
195fe7e6d57SJed Brown   }
196fe7e6d57SJed Brown   {
197fe7e6d57SJed Brown     const PetscReal g = 4.3586652150845900e-01;
198fe7e6d57SJed Brown     const PetscReal
199fe7e6d57SJed Brown       A[4][4] = {{0,0,0,0},
200fe7e6d57SJed Brown                  {8.7173304301691801e-01,0,0,0},
201fe7e6d57SJed Brown                  {8.4457060015369423e-01,-1.1299064236484185e-01,0,0},
202fe7e6d57SJed Brown                  {0,0,1.,0}},
203fe7e6d57SJed Brown       Gamma[4][4] = {{g,0,0,0},
204fe7e6d57SJed Brown                      {-8.7173304301691801e-01,g,0,0},
205fe7e6d57SJed Brown                      {-9.0338057013044082e-01,5.4180672388095326e-02,g,0},
206fe7e6d57SJed Brown                      {2.4212380706095346e-01,-1.2232505839045147e+00,5.4526025533510214e-01,g}},
207fe7e6d57SJed Brown       b[4] = {2.4212380706095346e-01,-1.2232505839045147e+00,1.5452602553351020e+00,4.3586652150845900e-01},
208fe7e6d57SJed Brown       b2[4] = {3.7810903145819369e-01,-9.6042292212423178e-02,5.0000000000000000e-01,2.1793326075422950e-01};
209fe7e6d57SJed Brown     ierr = TSRosWRegister(TSROSWRA34PW2,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
210e27a552bSJed Brown   }
211ef3c5b88SJed Brown   {
212ef3c5b88SJed Brown     const PetscReal g = 0.5;
213ef3c5b88SJed Brown     const PetscReal
214ef3c5b88SJed Brown       A[4][4] = {{0,0,0,0},
215ef3c5b88SJed Brown                  {0,0,0,0},
216ef3c5b88SJed Brown                  {1.,0,0,0},
217ef3c5b88SJed Brown                  {0.75,-0.25,0.5,0}},
218ef3c5b88SJed Brown       Gamma[4][4] = {{g,0,0,0},
219ef3c5b88SJed Brown                      {1.,g,0,0},
220ef3c5b88SJed Brown                      {-0.25,-0.25,g,0},
221ef3c5b88SJed Brown                      {1./12,1./12,-2./3,g}},
222ef3c5b88SJed Brown       b[4] = {5./6,-1./6,-1./6,0.5},
223ef3c5b88SJed Brown       b2[4] = {0.75,-0.25,0.5,0};
224ef3c5b88SJed Brown     ierr = TSRosWRegister(TSROSWRODAS3,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
225ef3c5b88SJed Brown   }
226ef3c5b88SJed Brown   {
227ef3c5b88SJed Brown     const PetscReal g = 0.43586652150845899941601945119356;
228ef3c5b88SJed Brown     const PetscReal
229ef3c5b88SJed Brown       A[3][3] = {{0,0,0},
230ef3c5b88SJed Brown                  {g,0,0},
231ef3c5b88SJed Brown                  {g,0,0}},
232ef3c5b88SJed Brown       Gamma[3][3] = {{g,0,0},
233ef3c5b88SJed Brown                      {-0.19294655696029095575009695436041,g,0},
234ef3c5b88SJed Brown                      {0,1.74927148125794685173529749738960,g}},
235ef3c5b88SJed Brown       b[3] = {-0.75457412385404315829818998646589,1.94100407061964420292840123379419,-0.18642994676560104463021124732829},
236ef3c5b88SJed Brown       b2[3] = {-1.53358745784149585370766523913002,2.81745131148625772213931745457622,-0.28386385364476186843165221544619};
237ef3c5b88SJed Brown     ierr = TSRosWRegister(TSROSWSANDU3,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
238ef3c5b88SJed Brown   }
239*b1c69cc3SEmil Constantinescu   {
240*b1c69cc3SEmil Constantinescu     const PetscReal g = (3.0+sqrt(3.0))/6.0;
241*b1c69cc3SEmil Constantinescu     const PetscReal
242*b1c69cc3SEmil Constantinescu       A[3][3] = {{0,0,0},
243*b1c69cc3SEmil Constantinescu                  {1,0,0},
244*b1c69cc3SEmil Constantinescu                  {0.25,0.25,0}},
245*b1c69cc3SEmil Constantinescu       Gamma[3][3] = {{0,0,0},
246*b1c69cc3SEmil Constantinescu                      {(-3.0-sqrt(3.0))/6.0,g,0},
247*b1c69cc3SEmil Constantinescu                      {(-3.0-sqrt(3.0))/24.0,(-3.0-sqrt(3.0))/8.0,g}},
248*b1c69cc3SEmil Constantinescu         b[3] = {1./6.,1./6.,2./3.},
249*b1c69cc3SEmil Constantinescu           b2[3] = {1./4.,1./4.,1./2.};
250*b1c69cc3SEmil Constantinescu     ierr = TSRosWRegister(TSROSWASSP3P3S1C,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
251*b1c69cc3SEmil Constantinescu   }
252*b1c69cc3SEmil Constantinescu 
253*b1c69cc3SEmil Constantinescu   {
254*b1c69cc3SEmil Constantinescu     const PetscReal
255*b1c69cc3SEmil Constantinescu       A[4][4] = {{0,0,0,0},
256*b1c69cc3SEmil Constantinescu                  {1./2.,0,0,0},
257*b1c69cc3SEmil Constantinescu                  {1./2.,1./2.,0,0},
258*b1c69cc3SEmil Constantinescu                  {1./6.,1./6.,1./6.,0}},
259*b1c69cc3SEmil Constantinescu       Gamma[4][4] = {{1./2.,0,0,0},
260*b1c69cc3SEmil Constantinescu                      {0.0,1./4.,0,0},
261*b1c69cc3SEmil Constantinescu                      {-2.,-2./3.,2./3.,0},
262*b1c69cc3SEmil Constantinescu                      {1./2.,5./36.,-2./9,0}},
263*b1c69cc3SEmil Constantinescu         b[4] = {1./6.,1./6.,1./6.,1./2.},
264*b1c69cc3SEmil Constantinescu         b2[4] = {1./8.,3./4.,1./8.,0};
265*b1c69cc3SEmil Constantinescu      ierr = TSRosWRegister(TSROSWLASSP3P4S2C,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
266*b1c69cc3SEmil Constantinescu   }
267*b1c69cc3SEmil Constantinescu 
268*b1c69cc3SEmil Constantinescu   {
269*b1c69cc3SEmil Constantinescu     const PetscReal
270*b1c69cc3SEmil Constantinescu       A[4][4] = {{0,0,0,0},
271*b1c69cc3SEmil Constantinescu                  {1./2.,0,0,0},
272*b1c69cc3SEmil Constantinescu                  {1./2.,1./2.,0,0},
273*b1c69cc3SEmil Constantinescu                  {1./6.,1./6.,1./6.,0}},
274*b1c69cc3SEmil Constantinescu       Gamma[4][4] = {{1./2.,0,0,0},
275*b1c69cc3SEmil Constantinescu                      {0.0,3./4.,0,0},
276*b1c69cc3SEmil Constantinescu                      {-2./3.,-23./9.,2./9.,0},
277*b1c69cc3SEmil Constantinescu                      {1./18.,65./108.,-2./27,0}},
278*b1c69cc3SEmil Constantinescu         b[4] = {1./6.,1./6.,1./6.,1./2.},
279*b1c69cc3SEmil Constantinescu         b2[4] = {3./16.,10./16.,3./16.,0};
280*b1c69cc3SEmil Constantinescu      ierr = TSRosWRegister(TSROSWLLSSP3P3S2C,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
281*b1c69cc3SEmil Constantinescu   }
282e27a552bSJed Brown   PetscFunctionReturn(0);
283e27a552bSJed Brown }
284e27a552bSJed Brown 
285e27a552bSJed Brown #undef __FUNCT__
286e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterDestroy"
287e27a552bSJed Brown /*@C
288e27a552bSJed Brown    TSRosWRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister().
289e27a552bSJed Brown 
290e27a552bSJed Brown    Not Collective
291e27a552bSJed Brown 
292e27a552bSJed Brown    Level: advanced
293e27a552bSJed Brown 
294e27a552bSJed Brown .keywords: TSRosW, register, destroy
295e27a552bSJed Brown .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegisterDynamic()
296e27a552bSJed Brown @*/
297e27a552bSJed Brown PetscErrorCode TSRosWRegisterDestroy(void)
298e27a552bSJed Brown {
299e27a552bSJed Brown   PetscErrorCode ierr;
30061692a83SJed Brown   RosWTableauLink link;
301e27a552bSJed Brown 
302e27a552bSJed Brown   PetscFunctionBegin;
30361692a83SJed Brown   while ((link = RosWTableauList)) {
30461692a83SJed Brown     RosWTableau t = &link->tab;
30561692a83SJed Brown     RosWTableauList = link->next;
30661692a83SJed Brown     ierr = PetscFree5(t->A,t->Gamma,t->b,t->ASum,t->GammaSum);CHKERRQ(ierr);
307c17803e7SJed Brown     ierr = PetscFree4(t->At,t->bt,t->GammaInv,t->GammaZeroDiag);CHKERRQ(ierr);
308fe7e6d57SJed Brown     ierr = PetscFree2(t->bembed,t->bembedt);CHKERRQ(ierr);
309e27a552bSJed Brown     ierr = PetscFree(t->name);CHKERRQ(ierr);
310e27a552bSJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
311e27a552bSJed Brown   }
312e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_FALSE;
313e27a552bSJed Brown   PetscFunctionReturn(0);
314e27a552bSJed Brown }
315e27a552bSJed Brown 
316e27a552bSJed Brown #undef __FUNCT__
317e27a552bSJed Brown #define __FUNCT__ "TSRosWInitializePackage"
318e27a552bSJed Brown /*@C
319e27a552bSJed Brown   TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called
320e27a552bSJed Brown   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RosW()
321e27a552bSJed Brown   when using static libraries.
322e27a552bSJed Brown 
323e27a552bSJed Brown   Input Parameter:
324e27a552bSJed Brown   path - The dynamic library path, or PETSC_NULL
325e27a552bSJed Brown 
326e27a552bSJed Brown   Level: developer
327e27a552bSJed Brown 
328e27a552bSJed Brown .keywords: TS, TSRosW, initialize, package
329e27a552bSJed Brown .seealso: PetscInitialize()
330e27a552bSJed Brown @*/
331e27a552bSJed Brown PetscErrorCode TSRosWInitializePackage(const char path[])
332e27a552bSJed Brown {
333e27a552bSJed Brown   PetscErrorCode ierr;
334e27a552bSJed Brown 
335e27a552bSJed Brown   PetscFunctionBegin;
336e27a552bSJed Brown   if (TSRosWPackageInitialized) PetscFunctionReturn(0);
337e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_TRUE;
338e27a552bSJed Brown   ierr = TSRosWRegisterAll();CHKERRQ(ierr);
339e27a552bSJed Brown   ierr = PetscRegisterFinalize(TSRosWFinalizePackage);CHKERRQ(ierr);
340e27a552bSJed Brown   PetscFunctionReturn(0);
341e27a552bSJed Brown }
342e27a552bSJed Brown 
343e27a552bSJed Brown #undef __FUNCT__
344e27a552bSJed Brown #define __FUNCT__ "TSRosWFinalizePackage"
345e27a552bSJed Brown /*@C
346e27a552bSJed Brown   TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is
347e27a552bSJed Brown   called from PetscFinalize().
348e27a552bSJed Brown 
349e27a552bSJed Brown   Level: developer
350e27a552bSJed Brown 
351e27a552bSJed Brown .keywords: Petsc, destroy, package
352e27a552bSJed Brown .seealso: PetscFinalize()
353e27a552bSJed Brown @*/
354e27a552bSJed Brown PetscErrorCode TSRosWFinalizePackage(void)
355e27a552bSJed Brown {
356e27a552bSJed Brown   PetscErrorCode ierr;
357e27a552bSJed Brown 
358e27a552bSJed Brown   PetscFunctionBegin;
359e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_FALSE;
360e27a552bSJed Brown   ierr = TSRosWRegisterDestroy();CHKERRQ(ierr);
361e27a552bSJed Brown   PetscFunctionReturn(0);
362e27a552bSJed Brown }
363e27a552bSJed Brown 
364e27a552bSJed Brown #undef __FUNCT__
365e27a552bSJed Brown #define __FUNCT__ "TSRosWRegister"
366e27a552bSJed Brown /*@C
36761692a83SJed Brown    TSRosWRegister - register a Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
368e27a552bSJed Brown 
369e27a552bSJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
370e27a552bSJed Brown 
371e27a552bSJed Brown    Input Parameters:
372e27a552bSJed Brown +  name - identifier for method
373e27a552bSJed Brown .  order - approximation order of method
374e27a552bSJed Brown .  s - number of stages, this is the dimension of the matrices below
37561692a83SJed Brown .  A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular
37661692a83SJed Brown .  Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal
377fe7e6d57SJed Brown .  b - Step completion table (dimension s)
378fe7e6d57SJed Brown -  bembed - Step completion table for a scheme of order one less (dimension s, PETSC_NULL if no embedded scheme is available)
379e27a552bSJed Brown 
380e27a552bSJed Brown    Notes:
38161692a83SJed Brown    Several Rosenbrock W methods are provided, this function is only needed to create new methods.
382e27a552bSJed Brown 
383e27a552bSJed Brown    Level: advanced
384e27a552bSJed Brown 
385e27a552bSJed Brown .keywords: TS, register
386e27a552bSJed Brown 
387e27a552bSJed Brown .seealso: TSRosW
388e27a552bSJed Brown @*/
389e27a552bSJed Brown PetscErrorCode TSRosWRegister(const TSRosWType name,PetscInt order,PetscInt s,
390fe7e6d57SJed Brown                               const PetscReal A[],const PetscReal Gamma[],const PetscReal b[],const PetscReal bembed[])
391e27a552bSJed Brown {
392e27a552bSJed Brown   PetscErrorCode ierr;
39361692a83SJed Brown   RosWTableauLink link;
39461692a83SJed Brown   RosWTableau t;
39561692a83SJed Brown   PetscInt i,j,k;
39661692a83SJed Brown   PetscScalar *GammaInv;
397e27a552bSJed Brown 
398e27a552bSJed Brown   PetscFunctionBegin;
399fe7e6d57SJed Brown   PetscValidCharPointer(name,1);
400fe7e6d57SJed Brown   PetscValidPointer(A,4);
401fe7e6d57SJed Brown   PetscValidPointer(Gamma,5);
402fe7e6d57SJed Brown   PetscValidPointer(b,6);
403fe7e6d57SJed Brown   if (bembed) PetscValidPointer(bembed,7);
404fe7e6d57SJed Brown 
405e27a552bSJed Brown   ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
406e27a552bSJed Brown   ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
407e27a552bSJed Brown   t = &link->tab;
408e27a552bSJed Brown   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
409e27a552bSJed Brown   t->order = order;
410e27a552bSJed Brown   t->s = s;
41161692a83SJed 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);
412c17803e7SJed Brown   ierr = PetscMalloc4(s*s,PetscReal,&t->At,s,PetscReal,&t->bt,s*s,PetscReal,&t->GammaInv,s,PetscBool,&t->GammaZeroDiag);CHKERRQ(ierr);
413e27a552bSJed Brown   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
41461692a83SJed Brown   ierr = PetscMemcpy(t->Gamma,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr);
41561692a83SJed Brown   ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
416fe7e6d57SJed Brown   if (bembed) {
417fe7e6d57SJed Brown     ierr = PetscMalloc2(s,PetscReal,&t->bembed,s,PetscReal,&t->bembedt);CHKERRQ(ierr);
418fe7e6d57SJed Brown     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
419fe7e6d57SJed Brown   }
42061692a83SJed Brown   for (i=0; i<s; i++) {
42161692a83SJed Brown     t->ASum[i] = 0;
42261692a83SJed Brown     t->GammaSum[i] = 0;
42361692a83SJed Brown     for (j=0; j<s; j++) {
42461692a83SJed Brown       t->ASum[i] += A[i*s+j];
425fe7e6d57SJed Brown       t->GammaSum[i] += Gamma[i*s+j];
42661692a83SJed Brown     }
42761692a83SJed Brown   }
42861692a83SJed Brown   ierr = PetscMalloc(s*s*sizeof(PetscScalar),&GammaInv);CHKERRQ(ierr); /* Need to use Scalar for inverse, then convert back to Real */
42961692a83SJed Brown   for (i=0; i<s*s; i++) GammaInv[i] = Gamma[i];
430fd96d5b0SEmil Constantinescu   for (i=0; i<s; i++) {
431fd96d5b0SEmil Constantinescu     if (Gamma[i*s+i] == 0.0) {
432fd96d5b0SEmil Constantinescu       GammaInv[i*s+i] = 1.0;
433c17803e7SJed Brown       t->GammaZeroDiag[i] = PETSC_TRUE;
434fd96d5b0SEmil Constantinescu     } else {
435c17803e7SJed Brown       t->GammaZeroDiag[i] = PETSC_FALSE;
436fd96d5b0SEmil Constantinescu     }
437fd96d5b0SEmil Constantinescu   }
438fd96d5b0SEmil Constantinescu 
43961692a83SJed Brown   switch (s) {
44061692a83SJed Brown   case 1: GammaInv[0] = 1./GammaInv[0]; break;
44161692a83SJed Brown   case 2: ierr = Kernel_A_gets_inverse_A_2(GammaInv,0);CHKERRQ(ierr); break;
44261692a83SJed Brown   case 3: ierr = Kernel_A_gets_inverse_A_3(GammaInv,0);CHKERRQ(ierr); break;
44361692a83SJed Brown   case 4: ierr = Kernel_A_gets_inverse_A_4(GammaInv,0);CHKERRQ(ierr); break;
44461692a83SJed Brown   case 5: {
44561692a83SJed Brown     PetscInt ipvt5[5];
44661692a83SJed Brown     MatScalar work5[5*5];
44761692a83SJed Brown     ierr = Kernel_A_gets_inverse_A_5(GammaInv,ipvt5,work5,0);CHKERRQ(ierr); break;
44861692a83SJed Brown   }
44961692a83SJed Brown   case 6: ierr = Kernel_A_gets_inverse_A_6(GammaInv,0);CHKERRQ(ierr); break;
45061692a83SJed Brown   case 7: ierr = Kernel_A_gets_inverse_A_7(GammaInv,0);CHKERRQ(ierr); break;
45161692a83SJed Brown   default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for %D stages",s);
45261692a83SJed Brown   }
45361692a83SJed Brown   for (i=0; i<s*s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]);
45461692a83SJed Brown   ierr = PetscFree(GammaInv);CHKERRQ(ierr);
45561692a83SJed Brown   for (i=0; i<s; i++) {
45661692a83SJed Brown     for (j=0; j<s; j++) {
45761692a83SJed Brown       t->At[i*s+j] = 0;
45861692a83SJed Brown       for (k=0; k<s; k++) {
45961692a83SJed Brown         t->At[i*s+j] += t->A[i*s+k] * t->GammaInv[k*s+j];
46061692a83SJed Brown       }
46161692a83SJed Brown     }
46261692a83SJed Brown     t->bt[i] = 0;
46361692a83SJed Brown     for (j=0; j<s; j++) {
46461692a83SJed Brown       t->bt[i] += t->b[j] * t->GammaInv[j*s+i];
46561692a83SJed Brown     }
466fe7e6d57SJed Brown     if (bembed) {
467fe7e6d57SJed Brown       t->bembedt[i] = 0;
468fe7e6d57SJed Brown       for (j=0; j<s; j++) {
469fe7e6d57SJed Brown         t->bembedt[i] += t->bembed[j] * t->GammaInv[j*s+i];
470fe7e6d57SJed Brown       }
471fe7e6d57SJed Brown     }
47261692a83SJed Brown   }
4738d59e960SJed Brown   t->ccfl = 1.0;                /* Fix this */
4748d59e960SJed Brown 
47561692a83SJed Brown   link->next = RosWTableauList;
47661692a83SJed Brown   RosWTableauList = link;
477e27a552bSJed Brown   PetscFunctionReturn(0);
478e27a552bSJed Brown }
479e27a552bSJed Brown 
480e27a552bSJed Brown #undef __FUNCT__
4811c3436cfSJed Brown #define __FUNCT__ "TSEvaluateStep_RosW"
4821c3436cfSJed Brown /*
4831c3436cfSJed Brown  The step completion formula is
4841c3436cfSJed Brown 
4851c3436cfSJed Brown  x1 = x0 + b^T Y
4861c3436cfSJed Brown 
4871c3436cfSJed Brown  where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been
4881c3436cfSJed Brown  updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write
4891c3436cfSJed Brown 
4901c3436cfSJed Brown  x1e = x0 + be^T Y
4911c3436cfSJed Brown      = x1 - b^T Y + be^T Y
4921c3436cfSJed Brown      = x1 + (be - b)^T Y
4931c3436cfSJed Brown 
4941c3436cfSJed Brown  so we can evaluate the method of different order even after the step has been optimistically completed.
4951c3436cfSJed Brown */
4961c3436cfSJed Brown static PetscErrorCode TSEvaluateStep_RosW(TS ts,PetscInt order,Vec X,PetscBool *done)
4971c3436cfSJed Brown {
4981c3436cfSJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
4991c3436cfSJed Brown   RosWTableau    tab  = ros->tableau;
5001c3436cfSJed Brown   PetscScalar    *w = ros->work;
5011c3436cfSJed Brown   PetscInt       i;
5021c3436cfSJed Brown   PetscErrorCode ierr;
5031c3436cfSJed Brown 
5041c3436cfSJed Brown   PetscFunctionBegin;
5051c3436cfSJed Brown   if (order == tab->order) {
5061c3436cfSJed Brown     if (ros->step_taken) {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
5071c3436cfSJed Brown     else {
5081c3436cfSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
5091c3436cfSJed Brown       ierr = VecMAXPY(X,tab->s,tab->bt,ros->Y);CHKERRQ(ierr);
5101c3436cfSJed Brown     }
5111c3436cfSJed Brown     if (done) *done = PETSC_TRUE;
5121c3436cfSJed Brown     PetscFunctionReturn(0);
5131c3436cfSJed Brown   } else if (order == tab->order-1) {
5141c3436cfSJed Brown     if (!tab->bembedt) goto unavailable;
5151c3436cfSJed Brown     if (ros->step_taken) {
5161c3436cfSJed Brown       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i];
5171c3436cfSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
5181c3436cfSJed Brown       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
5191c3436cfSJed Brown     } else {
5201c3436cfSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
5211c3436cfSJed Brown       ierr = VecMAXPY(X,tab->s,tab->bembedt,ros->Y);CHKERRQ(ierr);
5221c3436cfSJed Brown     }
5231c3436cfSJed Brown     if (done) *done = PETSC_TRUE;
5241c3436cfSJed Brown     PetscFunctionReturn(0);
5251c3436cfSJed Brown   }
5261c3436cfSJed Brown   unavailable:
5271c3436cfSJed Brown   if (done) *done = PETSC_FALSE;
5281c3436cfSJed 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);
5291c3436cfSJed Brown   PetscFunctionReturn(0);
5301c3436cfSJed Brown }
5311c3436cfSJed Brown 
5321c3436cfSJed Brown #undef __FUNCT__
533e27a552bSJed Brown #define __FUNCT__ "TSStep_RosW"
534e27a552bSJed Brown static PetscErrorCode TSStep_RosW(TS ts)
535e27a552bSJed Brown {
53661692a83SJed Brown   TS_RosW         *ros = (TS_RosW*)ts->data;
53761692a83SJed Brown   RosWTableau     tab  = ros->tableau;
538e27a552bSJed Brown   const PetscInt  s    = tab->s;
5391c3436cfSJed Brown   const PetscReal *At  = tab->At,*Gamma = tab->Gamma,*ASum = tab->ASum,*GammaInv = tab->GammaInv;
540c17803e7SJed Brown   const PetscBool *GammaZeroDiag = tab->GammaZeroDiag;
54161692a83SJed Brown   PetscScalar     *w   = ros->work;
54261692a83SJed Brown   Vec             *Y   = ros->Y,Zdot = ros->Zdot,Zstage = ros->Zstage;
543e27a552bSJed Brown   SNES            snes;
5441c3436cfSJed Brown   TSAdapt         adapt;
5451c3436cfSJed Brown   PetscInt        i,j,its,lits,reject,next_scheme;
546cdbf8f93SLisandro Dalcin   PetscReal       next_time_step;
5471c3436cfSJed Brown   PetscBool       accept;
548e27a552bSJed Brown   PetscErrorCode  ierr;
549e27a552bSJed Brown 
550e27a552bSJed Brown   PetscFunctionBegin;
551e27a552bSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
552cdbf8f93SLisandro Dalcin   next_time_step = ts->time_step;
5531c3436cfSJed Brown   accept = PETSC_TRUE;
5541c3436cfSJed Brown   ros->step_taken = PETSC_FALSE;
555e27a552bSJed Brown 
5561c3436cfSJed Brown   for (reject=0; reject<ts->max_reject; reject++,ts->reject++) {
5571c3436cfSJed Brown     const PetscReal h = ts->time_step;
558e27a552bSJed Brown     for (i=0; i<s; i++) {
5591c3436cfSJed Brown       ros->stage_time = ts->ptime + h*ASum[i];
560c17803e7SJed Brown       if (GammaZeroDiag[i]) {
561c17803e7SJed Brown         ros->stage_explicit = PETSC_TRUE;
562fd96d5b0SEmil Constantinescu         ros->shift = 1./h;
563c17803e7SJed Brown       } else {
564c17803e7SJed Brown         ros->stage_explicit = PETSC_FALSE;
56561692a83SJed Brown         ros->shift = 1./(h*Gamma[i*s+i]);
566fd96d5b0SEmil Constantinescu       }
56761692a83SJed Brown 
56861692a83SJed Brown       ierr = VecCopy(ts->vec_sol,Zstage);CHKERRQ(ierr);
56961692a83SJed Brown       ierr = VecMAXPY(Zstage,i,&At[i*s+0],Y);CHKERRQ(ierr);
57061692a83SJed Brown 
57161692a83SJed Brown       for (j=0; j<i; j++) w[j] = 1./h * GammaInv[i*s+j];
57261692a83SJed Brown       ierr = VecZeroEntries(Zdot);CHKERRQ(ierr);
57361692a83SJed Brown       ierr = VecMAXPY(Zdot,i,w,Y);CHKERRQ(ierr);
57461692a83SJed Brown 
575e27a552bSJed Brown       /* Initial guess taken from last stage */
57661692a83SJed Brown       ierr = VecZeroEntries(Y[i]);CHKERRQ(ierr);
57761692a83SJed Brown 
57861692a83SJed Brown       if (!ros->recompute_jacobian && !i) {
57961692a83SJed Brown         ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */
58061692a83SJed Brown       }
58161692a83SJed Brown 
58261692a83SJed Brown       ierr = SNESSolve(snes,PETSC_NULL,Y[i]);CHKERRQ(ierr);
583e27a552bSJed Brown       ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
584e27a552bSJed Brown       ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
585e27a552bSJed Brown       ts->nonlinear_its += its; ts->linear_its += lits;
586e27a552bSJed Brown     }
5871c3436cfSJed Brown     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
5881c3436cfSJed Brown     ros->step_taken = PETSC_TRUE;
589e27a552bSJed Brown 
5901c3436cfSJed Brown     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
5911c3436cfSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
5921c3436cfSJed Brown     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
5938d59e960SJed Brown     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
5941c3436cfSJed Brown     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
5951c3436cfSJed Brown     if (accept) {
5961c3436cfSJed Brown       /* ignore next_scheme for now */
597e27a552bSJed Brown       ts->ptime += ts->time_step;
598cdbf8f93SLisandro Dalcin       ts->time_step = next_time_step;
599e27a552bSJed Brown       ts->steps++;
6001c3436cfSJed Brown       break;
6011c3436cfSJed Brown     } else {                    /* Roll back the current step */
6021c3436cfSJed Brown       for (i=0; i<s; i++) w[i] = -tab->bt[i];
6031c3436cfSJed Brown       ierr = VecMAXPY(ts->vec_sol,s,w,Y);CHKERRQ(ierr);
6041c3436cfSJed Brown       ts->time_step = next_time_step;
6051c3436cfSJed Brown       ros->step_taken = PETSC_FALSE;
6061c3436cfSJed Brown     }
6071c3436cfSJed Brown   }
6081c3436cfSJed Brown 
609e27a552bSJed Brown   PetscFunctionReturn(0);
610e27a552bSJed Brown }
611e27a552bSJed Brown 
612e27a552bSJed Brown #undef __FUNCT__
613e27a552bSJed Brown #define __FUNCT__ "TSInterpolate_RosW"
614e27a552bSJed Brown static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec X)
615e27a552bSJed Brown {
61661692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
617e27a552bSJed Brown 
618e27a552bSJed Brown   PetscFunctionBegin;
61961692a83SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name);
620e27a552bSJed Brown   PetscFunctionReturn(0);
621e27a552bSJed Brown }
622e27a552bSJed Brown 
623e27a552bSJed Brown /*------------------------------------------------------------*/
624e27a552bSJed Brown #undef __FUNCT__
625e27a552bSJed Brown #define __FUNCT__ "TSReset_RosW"
626e27a552bSJed Brown static PetscErrorCode TSReset_RosW(TS ts)
627e27a552bSJed Brown {
62861692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
629e27a552bSJed Brown   PetscInt       s;
630e27a552bSJed Brown   PetscErrorCode ierr;
631e27a552bSJed Brown 
632e27a552bSJed Brown   PetscFunctionBegin;
63361692a83SJed Brown   if (!ros->tableau) PetscFunctionReturn(0);
63461692a83SJed Brown    s = ros->tableau->s;
63561692a83SJed Brown   ierr = VecDestroyVecs(s,&ros->Y);CHKERRQ(ierr);
63661692a83SJed Brown   ierr = VecDestroy(&ros->Ydot);CHKERRQ(ierr);
63761692a83SJed Brown   ierr = VecDestroy(&ros->Ystage);CHKERRQ(ierr);
63861692a83SJed Brown   ierr = VecDestroy(&ros->Zdot);CHKERRQ(ierr);
63961692a83SJed Brown   ierr = VecDestroy(&ros->Zstage);CHKERRQ(ierr);
64061692a83SJed Brown   ierr = PetscFree(ros->work);CHKERRQ(ierr);
641e27a552bSJed Brown   PetscFunctionReturn(0);
642e27a552bSJed Brown }
643e27a552bSJed Brown 
644e27a552bSJed Brown #undef __FUNCT__
645e27a552bSJed Brown #define __FUNCT__ "TSDestroy_RosW"
646e27a552bSJed Brown static PetscErrorCode TSDestroy_RosW(TS ts)
647e27a552bSJed Brown {
648e27a552bSJed Brown   PetscErrorCode  ierr;
649e27a552bSJed Brown 
650e27a552bSJed Brown   PetscFunctionBegin;
651e27a552bSJed Brown   ierr = TSReset_RosW(ts);CHKERRQ(ierr);
652e27a552bSJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
653e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","",PETSC_NULL);CHKERRQ(ierr);
654e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","",PETSC_NULL);CHKERRQ(ierr);
65561692a83SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","",PETSC_NULL);CHKERRQ(ierr);
656e27a552bSJed Brown   PetscFunctionReturn(0);
657e27a552bSJed Brown }
658e27a552bSJed Brown 
659e27a552bSJed Brown /*
660e27a552bSJed Brown   This defines the nonlinear equation that is to be solved with SNES
661e27a552bSJed Brown   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
662e27a552bSJed Brown */
663e27a552bSJed Brown #undef __FUNCT__
664e27a552bSJed Brown #define __FUNCT__ "SNESTSFormFunction_RosW"
665e27a552bSJed Brown static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec X,Vec F,TS ts)
666e27a552bSJed Brown {
66761692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
668e27a552bSJed Brown   PetscErrorCode ierr;
669e27a552bSJed Brown 
670e27a552bSJed Brown   PetscFunctionBegin;
671c17803e7SJed Brown   if (ros->stage_explicit) {
672c17803e7SJed Brown     ierr = VecAXPBY(ros->Ydot,ros->shift,0.0,X);CHKERRQ(ierr);       /* Ydot = shift*X*/
673c17803e7SJed Brown   } else {
67461692a83SJed Brown     ierr = VecWAXPY(ros->Ydot,ros->shift,X,ros->Zdot);CHKERRQ(ierr); /* Ydot = shift*X + Zdot */
675c17803e7SJed Brown   }
67661692a83SJed Brown   ierr = VecWAXPY(ros->Ystage,1.0,X,ros->Zstage);CHKERRQ(ierr);    /* Ystage = X + Zstage */
67761692a83SJed Brown   ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,ros->Ydot,F,PETSC_FALSE);CHKERRQ(ierr);
678e27a552bSJed Brown   PetscFunctionReturn(0);
679e27a552bSJed Brown }
680e27a552bSJed Brown 
681e27a552bSJed Brown #undef __FUNCT__
682e27a552bSJed Brown #define __FUNCT__ "SNESTSFormJacobian_RosW"
683e27a552bSJed Brown static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
684e27a552bSJed Brown {
68561692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
686e27a552bSJed Brown   PetscErrorCode ierr;
687e27a552bSJed Brown 
688e27a552bSJed Brown   PetscFunctionBegin;
68961692a83SJed Brown   /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
69061692a83SJed Brown   ierr = TSComputeIJacobian(ts,ros->stage_time,ros->Ystage,ros->Ydot,ros->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr);
691e27a552bSJed Brown   PetscFunctionReturn(0);
692e27a552bSJed Brown }
693e27a552bSJed Brown 
694e27a552bSJed Brown #undef __FUNCT__
695e27a552bSJed Brown #define __FUNCT__ "TSSetUp_RosW"
696e27a552bSJed Brown static PetscErrorCode TSSetUp_RosW(TS ts)
697e27a552bSJed Brown {
69861692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
69961692a83SJed Brown   RosWTableau    tab  = ros->tableau;
700e27a552bSJed Brown   PetscInt       s    = tab->s;
701e27a552bSJed Brown   PetscErrorCode ierr;
702e27a552bSJed Brown 
703e27a552bSJed Brown   PetscFunctionBegin;
70461692a83SJed Brown   if (!ros->tableau) {
705e27a552bSJed Brown     ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);
706e27a552bSJed Brown   }
70761692a83SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ros->Y);CHKERRQ(ierr);
70861692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Ydot);CHKERRQ(ierr);
70961692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Ystage);CHKERRQ(ierr);
71061692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Zdot);CHKERRQ(ierr);
71161692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Zstage);CHKERRQ(ierr);
71261692a83SJed Brown   ierr = PetscMalloc(s*sizeof(ros->work[0]),&ros->work);CHKERRQ(ierr);
713e27a552bSJed Brown   PetscFunctionReturn(0);
714e27a552bSJed Brown }
715e27a552bSJed Brown /*------------------------------------------------------------*/
716e27a552bSJed Brown 
717e27a552bSJed Brown #undef __FUNCT__
718e27a552bSJed Brown #define __FUNCT__ "TSSetFromOptions_RosW"
719e27a552bSJed Brown static PetscErrorCode TSSetFromOptions_RosW(TS ts)
720e27a552bSJed Brown {
72161692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
722e27a552bSJed Brown   PetscErrorCode ierr;
72361692a83SJed Brown   char           rostype[256];
724e27a552bSJed Brown 
725e27a552bSJed Brown   PetscFunctionBegin;
726e27a552bSJed Brown   ierr = PetscOptionsHead("RosW ODE solver options");CHKERRQ(ierr);
727e27a552bSJed Brown   {
72861692a83SJed Brown     RosWTableauLink link;
729e27a552bSJed Brown     PetscInt count,choice;
730e27a552bSJed Brown     PetscBool flg;
731e27a552bSJed Brown     const char **namelist;
73261692a83SJed Brown     SNES snes;
73361692a83SJed Brown 
73461692a83SJed Brown     ierr = PetscStrncpy(rostype,TSRosWDefault,sizeof rostype);CHKERRQ(ierr);
73561692a83SJed Brown     for (link=RosWTableauList,count=0; link; link=link->next,count++) ;
736e27a552bSJed Brown     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
73761692a83SJed Brown     for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
73861692a83SJed Brown     ierr = PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,rostype,&choice,&flg);CHKERRQ(ierr);
73961692a83SJed Brown     ierr = TSRosWSetType(ts,flg ? namelist[choice] : rostype);CHKERRQ(ierr);
740e27a552bSJed Brown     ierr = PetscFree(namelist);CHKERRQ(ierr);
74161692a83SJed Brown 
74261692a83SJed Brown     ierr = PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,PETSC_NULL);CHKERRQ(ierr);
74361692a83SJed Brown 
74461692a83SJed Brown     /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
74561692a83SJed Brown     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
74661692a83SJed Brown     if (!((PetscObject)snes)->type_name) {
74761692a83SJed Brown       ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
74861692a83SJed Brown     }
74961692a83SJed Brown     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
750e27a552bSJed Brown   }
751e27a552bSJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
752e27a552bSJed Brown   PetscFunctionReturn(0);
753e27a552bSJed Brown }
754e27a552bSJed Brown 
755e27a552bSJed Brown #undef __FUNCT__
756e27a552bSJed Brown #define __FUNCT__ "PetscFormatRealArray"
757e27a552bSJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
758e27a552bSJed Brown {
759e27a552bSJed Brown   PetscErrorCode ierr;
760e408995aSJed Brown   PetscInt i;
761e408995aSJed Brown   size_t left,count;
762e27a552bSJed Brown   char *p;
763e27a552bSJed Brown 
764e27a552bSJed Brown   PetscFunctionBegin;
765e408995aSJed Brown   for (i=0,p=buf,left=len; i<n; i++) {
766e408995aSJed Brown     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
767e27a552bSJed Brown     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
768e27a552bSJed Brown     left -= count;
769e27a552bSJed Brown     p += count;
770e27a552bSJed Brown     *p++ = ' ';
771e27a552bSJed Brown   }
772e27a552bSJed Brown   p[i ? 0 : -1] = 0;
773e27a552bSJed Brown   PetscFunctionReturn(0);
774e27a552bSJed Brown }
775e27a552bSJed Brown 
776e27a552bSJed Brown #undef __FUNCT__
777e27a552bSJed Brown #define __FUNCT__ "TSView_RosW"
778e27a552bSJed Brown static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer)
779e27a552bSJed Brown {
78061692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
78161692a83SJed Brown   RosWTableau    tab  = ros->tableau;
782e27a552bSJed Brown   PetscBool      iascii;
783e27a552bSJed Brown   PetscErrorCode ierr;
784e27a552bSJed Brown 
785e27a552bSJed Brown   PetscFunctionBegin;
786e27a552bSJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
787e27a552bSJed Brown   if (iascii) {
78861692a83SJed Brown     const TSRosWType rostype;
789e408995aSJed Brown     PetscInt i;
790e408995aSJed Brown     PetscReal abscissa[512];
791e27a552bSJed Brown     char buf[512];
79261692a83SJed Brown     ierr = TSRosWGetType(ts,&rostype);CHKERRQ(ierr);
79361692a83SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Rosenbrock-W %s\n",rostype);CHKERRQ(ierr);
794e408995aSJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ASum);CHKERRQ(ierr);
79561692a83SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A       = %s\n",buf);CHKERRQ(ierr);
796e408995aSJed Brown     for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i];
797e408995aSJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,abscissa);CHKERRQ(ierr);
798e408995aSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A+Gamma = %s\n",buf);CHKERRQ(ierr);
799e27a552bSJed Brown   }
800e27a552bSJed Brown   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
801e27a552bSJed Brown   PetscFunctionReturn(0);
802e27a552bSJed Brown }
803e27a552bSJed Brown 
804e27a552bSJed Brown #undef __FUNCT__
805e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType"
806e27a552bSJed Brown /*@C
80761692a83SJed Brown   TSRosWSetType - Set the type of Rosenbrock-W scheme
808e27a552bSJed Brown 
809e27a552bSJed Brown   Logically collective
810e27a552bSJed Brown 
811e27a552bSJed Brown   Input Parameter:
812e27a552bSJed Brown +  ts - timestepping context
81361692a83SJed Brown -  rostype - type of Rosenbrock-W scheme
814e27a552bSJed Brown 
815e27a552bSJed Brown   Level: intermediate
816e27a552bSJed Brown 
817e27a552bSJed Brown .seealso: TSRosWGetType()
818e27a552bSJed Brown @*/
81961692a83SJed Brown PetscErrorCode TSRosWSetType(TS ts,const TSRosWType rostype)
820e27a552bSJed Brown {
821e27a552bSJed Brown   PetscErrorCode ierr;
822e27a552bSJed Brown 
823e27a552bSJed Brown   PetscFunctionBegin;
824e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
82561692a83SJed Brown   ierr = PetscTryMethod(ts,"TSRosWSetType_C",(TS,const TSRosWType),(ts,rostype));CHKERRQ(ierr);
826e27a552bSJed Brown   PetscFunctionReturn(0);
827e27a552bSJed Brown }
828e27a552bSJed Brown 
829e27a552bSJed Brown #undef __FUNCT__
830e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType"
831e27a552bSJed Brown /*@C
83261692a83SJed Brown   TSRosWGetType - Get the type of Rosenbrock-W scheme
833e27a552bSJed Brown 
834e27a552bSJed Brown   Logically collective
835e27a552bSJed Brown 
836e27a552bSJed Brown   Input Parameter:
837e27a552bSJed Brown .  ts - timestepping context
838e27a552bSJed Brown 
839e27a552bSJed Brown   Output Parameter:
84061692a83SJed Brown .  rostype - type of Rosenbrock-W scheme
841e27a552bSJed Brown 
842e27a552bSJed Brown   Level: intermediate
843e27a552bSJed Brown 
844e27a552bSJed Brown .seealso: TSRosWGetType()
845e27a552bSJed Brown @*/
84661692a83SJed Brown PetscErrorCode TSRosWGetType(TS ts,const TSRosWType *rostype)
847e27a552bSJed Brown {
848e27a552bSJed Brown   PetscErrorCode ierr;
849e27a552bSJed Brown 
850e27a552bSJed Brown   PetscFunctionBegin;
851e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
85261692a83SJed Brown   ierr = PetscUseMethod(ts,"TSRosWGetType_C",(TS,const TSRosWType*),(ts,rostype));CHKERRQ(ierr);
853e27a552bSJed Brown   PetscFunctionReturn(0);
854e27a552bSJed Brown }
855e27a552bSJed Brown 
856e27a552bSJed Brown #undef __FUNCT__
85761692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian"
858e27a552bSJed Brown /*@C
85961692a83SJed Brown   TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step.
860e27a552bSJed Brown 
861e27a552bSJed Brown   Logically collective
862e27a552bSJed Brown 
863e27a552bSJed Brown   Input Parameter:
864e27a552bSJed Brown +  ts - timestepping context
86561692a83SJed Brown -  flg - PETSC_TRUE to recompute the Jacobian at each stage
866e27a552bSJed Brown 
867e27a552bSJed Brown   Level: intermediate
868e27a552bSJed Brown 
869e27a552bSJed Brown .seealso: TSRosWGetType()
870e27a552bSJed Brown @*/
87161692a83SJed Brown PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg)
872e27a552bSJed Brown {
873e27a552bSJed Brown   PetscErrorCode ierr;
874e27a552bSJed Brown 
875e27a552bSJed Brown   PetscFunctionBegin;
876e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
87761692a83SJed Brown   ierr = PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
878e27a552bSJed Brown   PetscFunctionReturn(0);
879e27a552bSJed Brown }
880e27a552bSJed Brown 
881e27a552bSJed Brown EXTERN_C_BEGIN
882e27a552bSJed Brown #undef __FUNCT__
883e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType_RosW"
88461692a83SJed Brown PetscErrorCode  TSRosWGetType_RosW(TS ts,const TSRosWType *rostype)
885e27a552bSJed Brown {
88661692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
887e27a552bSJed Brown   PetscErrorCode ierr;
888e27a552bSJed Brown 
889e27a552bSJed Brown   PetscFunctionBegin;
89061692a83SJed Brown   if (!ros->tableau) {ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);}
89161692a83SJed Brown   *rostype = ros->tableau->name;
892e27a552bSJed Brown   PetscFunctionReturn(0);
893e27a552bSJed Brown }
894e27a552bSJed Brown #undef __FUNCT__
895e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType_RosW"
89661692a83SJed Brown PetscErrorCode  TSRosWSetType_RosW(TS ts,const TSRosWType rostype)
897e27a552bSJed Brown {
89861692a83SJed Brown   TS_RosW         *ros = (TS_RosW*)ts->data;
899e27a552bSJed Brown   PetscErrorCode  ierr;
900e27a552bSJed Brown   PetscBool       match;
90161692a83SJed Brown   RosWTableauLink link;
902e27a552bSJed Brown 
903e27a552bSJed Brown   PetscFunctionBegin;
90461692a83SJed Brown   if (ros->tableau) {
90561692a83SJed Brown     ierr = PetscStrcmp(ros->tableau->name,rostype,&match);CHKERRQ(ierr);
906e27a552bSJed Brown     if (match) PetscFunctionReturn(0);
907e27a552bSJed Brown   }
90861692a83SJed Brown   for (link = RosWTableauList; link; link=link->next) {
90961692a83SJed Brown     ierr = PetscStrcmp(link->tab.name,rostype,&match);CHKERRQ(ierr);
910e27a552bSJed Brown     if (match) {
911e27a552bSJed Brown       ierr = TSReset_RosW(ts);CHKERRQ(ierr);
91261692a83SJed Brown       ros->tableau = &link->tab;
913e27a552bSJed Brown       PetscFunctionReturn(0);
914e27a552bSJed Brown     }
915e27a552bSJed Brown   }
91661692a83SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype);
917e27a552bSJed Brown   PetscFunctionReturn(0);
918e27a552bSJed Brown }
91961692a83SJed Brown 
920e27a552bSJed Brown #undef __FUNCT__
92161692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian_RosW"
92261692a83SJed Brown PetscErrorCode  TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg)
923e27a552bSJed Brown {
92461692a83SJed Brown   TS_RosW *ros = (TS_RosW*)ts->data;
925e27a552bSJed Brown 
926e27a552bSJed Brown   PetscFunctionBegin;
92761692a83SJed Brown   ros->recompute_jacobian = flg;
928e27a552bSJed Brown   PetscFunctionReturn(0);
929e27a552bSJed Brown }
930e27a552bSJed Brown EXTERN_C_END
931e27a552bSJed Brown 
932e27a552bSJed Brown /* ------------------------------------------------------------ */
933e27a552bSJed Brown /*MC
934e27a552bSJed Brown       TSRosW - ODE solver using Rosenbrock-W schemes
935e27a552bSJed Brown 
936e27a552bSJed Brown   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
937e27a552bSJed Brown   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
938e27a552bSJed Brown   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
939e27a552bSJed Brown 
940e27a552bSJed Brown   Notes:
94161692a83SJed Brown   This method currently only works with autonomous ODE and DAE.
94261692a83SJed Brown 
94361692a83SJed Brown   Developer notes:
94461692a83SJed Brown   Rosenbrock-W methods are typically specified for autonomous ODE
94561692a83SJed Brown 
94661692a83SJed Brown $  xdot = f(x)
94761692a83SJed Brown 
94861692a83SJed Brown   by the stage equations
94961692a83SJed Brown 
95061692a83SJed Brown $  k_i = h f(x_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j
95161692a83SJed Brown 
95261692a83SJed Brown   and step completion formula
95361692a83SJed Brown 
95461692a83SJed Brown $  x_1 = x_0 + sum_j b_j k_j
95561692a83SJed Brown 
95661692a83SJed Brown   with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(x)
95761692a83SJed Brown   and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner,
95861692a83SJed Brown   we define new variables for the stage equations
95961692a83SJed Brown 
96061692a83SJed Brown $  y_i = gamma_ij k_j
96161692a83SJed Brown 
96261692a83SJed Brown   The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define
96361692a83SJed Brown 
96461692a83SJed Brown $  A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-i}
96561692a83SJed Brown 
96661692a83SJed Brown   to rewrite the method as
96761692a83SJed Brown 
96861692a83SJed 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
96961692a83SJed Brown $  x_1 = x_0 + sum_j bt_j y_j
97061692a83SJed Brown 
97161692a83SJed Brown    where we have introduced the mass matrix M. Continue by defining
97261692a83SJed Brown 
97361692a83SJed Brown $  ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j
97461692a83SJed Brown 
97561692a83SJed Brown    or, more compactly in tensor notation
97661692a83SJed Brown 
97761692a83SJed Brown $  Ydot = 1/h (Gamma^{-1} \otimes I) Y .
97861692a83SJed Brown 
97961692a83SJed Brown    Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current
98061692a83SJed Brown    stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the
98161692a83SJed Brown    equation
98261692a83SJed Brown 
98361692a83SJed Brown $  g(x_0 + sum_j a_ij y_j + y_i, ydot_i) = 0
98461692a83SJed Brown 
98561692a83SJed Brown    with initial guess y_i = 0.
986e27a552bSJed Brown 
987e27a552bSJed Brown   Level: beginner
988e27a552bSJed Brown 
989e27a552bSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSRosWRegister()
990e27a552bSJed Brown 
991e27a552bSJed Brown M*/
992e27a552bSJed Brown EXTERN_C_BEGIN
993e27a552bSJed Brown #undef __FUNCT__
994e27a552bSJed Brown #define __FUNCT__ "TSCreate_RosW"
995e27a552bSJed Brown PetscErrorCode  TSCreate_RosW(TS ts)
996e27a552bSJed Brown {
99761692a83SJed Brown   TS_RosW        *ros;
998e27a552bSJed Brown   PetscErrorCode ierr;
999e27a552bSJed Brown 
1000e27a552bSJed Brown   PetscFunctionBegin;
1001e27a552bSJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
1002e27a552bSJed Brown   ierr = TSRosWInitializePackage(PETSC_NULL);CHKERRQ(ierr);
1003e27a552bSJed Brown #endif
1004e27a552bSJed Brown 
1005e27a552bSJed Brown   ts->ops->reset          = TSReset_RosW;
1006e27a552bSJed Brown   ts->ops->destroy        = TSDestroy_RosW;
1007e27a552bSJed Brown   ts->ops->view           = TSView_RosW;
1008e27a552bSJed Brown   ts->ops->setup          = TSSetUp_RosW;
1009e27a552bSJed Brown   ts->ops->step           = TSStep_RosW;
1010e27a552bSJed Brown   ts->ops->interpolate    = TSInterpolate_RosW;
10111c3436cfSJed Brown   ts->ops->evaluatestep   = TSEvaluateStep_RosW;
1012e27a552bSJed Brown   ts->ops->setfromoptions = TSSetFromOptions_RosW;
1013e27a552bSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_RosW;
1014e27a552bSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_RosW;
1015e27a552bSJed Brown 
101661692a83SJed Brown   ierr = PetscNewLog(ts,TS_RosW,&ros);CHKERRQ(ierr);
101761692a83SJed Brown   ts->data = (void*)ros;
1018e27a552bSJed Brown 
1019e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","TSRosWGetType_RosW",TSRosWGetType_RosW);CHKERRQ(ierr);
1020e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","TSRosWSetType_RosW",TSRosWSetType_RosW);CHKERRQ(ierr);
102161692a83SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","TSRosWSetRecomputeJacobian_RosW",TSRosWSetRecomputeJacobian_RosW);CHKERRQ(ierr);
1022e27a552bSJed Brown   PetscFunctionReturn(0);
1023e27a552bSJed Brown }
1024e27a552bSJed Brown EXTERN_C_END
1025