xref: /petsc/src/ts/impls/rosw/rosw.c (revision ef3c5b884be35d09441e6283d7315bbe6c271bd2)
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 
111*ef3c5b88SJed Brown /*MC
112*ef3c5b88SJed Brown      TSROSWRODAS3 - Four stage third order L-stable Rosenbrock scheme
113*ef3c5b88SJed Brown 
114*ef3c5b88SJed Brown      By default, the Jacobian is only recomputed once per step.
115*ef3c5b88SJed Brown 
116*ef3c5b88SJed Brown      Both the third order and embedded second order methods are stiffly accurate and L-stable.
117*ef3c5b88SJed Brown 
118*ef3c5b88SJed Brown      References:
119*ef3c5b88SJed Brown      Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997.
120*ef3c5b88SJed Brown 
121*ef3c5b88SJed Brown      Level: intermediate
122*ef3c5b88SJed Brown 
123*ef3c5b88SJed Brown .seealso: TSROSW, TSROSWSANDU3
124*ef3c5b88SJed Brown M*/
125*ef3c5b88SJed Brown 
126*ef3c5b88SJed Brown /*MC
127*ef3c5b88SJed Brown      TSROSWSANDU3 - Three stage third order L-stable Rosenbrock scheme
128*ef3c5b88SJed Brown 
129*ef3c5b88SJed Brown      By default, the Jacobian is only recomputed once per step.
130*ef3c5b88SJed Brown 
131*ef3c5b88SJed Brown      The third order method is L-stable, but not stiffly accurate.
132*ef3c5b88SJed Brown      The second order embedded method is strongly A-stable with R(infty) = 0.5.
133*ef3c5b88SJed Brown      The internal stages are L-stable.
134*ef3c5b88SJed Brown      This method is called ROS3 in the paper.
135*ef3c5b88SJed Brown 
136*ef3c5b88SJed Brown      References:
137*ef3c5b88SJed Brown      Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997.
138*ef3c5b88SJed Brown 
139*ef3c5b88SJed Brown      Level: intermediate
140*ef3c5b88SJed Brown 
141*ef3c5b88SJed Brown .seealso: TSROSW, TSROSWRODAS3
142*ef3c5b88SJed Brown M*/
143*ef3c5b88SJed 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   }
211*ef3c5b88SJed Brown   {
212*ef3c5b88SJed Brown     const PetscReal g = 0.5;
213*ef3c5b88SJed Brown     const PetscReal
214*ef3c5b88SJed Brown       A[4][4] = {{0,0,0,0},
215*ef3c5b88SJed Brown                  {0,0,0,0},
216*ef3c5b88SJed Brown                  {1.,0,0,0},
217*ef3c5b88SJed Brown                  {0.75,-0.25,0.5,0}},
218*ef3c5b88SJed Brown       Gamma[4][4] = {{g,0,0,0},
219*ef3c5b88SJed Brown                      {1.,g,0,0},
220*ef3c5b88SJed Brown                      {-0.25,-0.25,g,0},
221*ef3c5b88SJed Brown                      {1./12,1./12,-2./3,g}},
222*ef3c5b88SJed Brown       b[4] = {5./6,-1./6,-1./6,0.5},
223*ef3c5b88SJed Brown       b2[4] = {0.75,-0.25,0.5,0};
224*ef3c5b88SJed Brown     ierr = TSRosWRegister(TSROSWRODAS3,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
225*ef3c5b88SJed Brown   }
226*ef3c5b88SJed Brown   {
227*ef3c5b88SJed Brown     const PetscReal g = 0.43586652150845899941601945119356;
228*ef3c5b88SJed Brown     const PetscReal
229*ef3c5b88SJed Brown       A[3][3] = {{0,0,0},
230*ef3c5b88SJed Brown                  {g,0,0},
231*ef3c5b88SJed Brown                  {g,0,0}},
232*ef3c5b88SJed Brown       Gamma[3][3] = {{g,0,0},
233*ef3c5b88SJed Brown                      {-0.19294655696029095575009695436041,g,0},
234*ef3c5b88SJed Brown                      {0,1.74927148125794685173529749738960,g}},
235*ef3c5b88SJed Brown       b[3] = {-0.75457412385404315829818998646589,1.94100407061964420292840123379419,-0.18642994676560104463021124732829},
236*ef3c5b88SJed Brown       b2[3] = {-1.53358745784149585370766523913002,2.81745131148625772213931745457622,-0.28386385364476186843165221544619};
237*ef3c5b88SJed Brown     ierr = TSRosWRegister(TSROSWSANDU3,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
238*ef3c5b88SJed Brown   }
239e27a552bSJed Brown   PetscFunctionReturn(0);
240e27a552bSJed Brown }
241e27a552bSJed Brown 
242e27a552bSJed Brown #undef __FUNCT__
243e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterDestroy"
244e27a552bSJed Brown /*@C
245e27a552bSJed Brown    TSRosWRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister().
246e27a552bSJed Brown 
247e27a552bSJed Brown    Not Collective
248e27a552bSJed Brown 
249e27a552bSJed Brown    Level: advanced
250e27a552bSJed Brown 
251e27a552bSJed Brown .keywords: TSRosW, register, destroy
252e27a552bSJed Brown .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegisterDynamic()
253e27a552bSJed Brown @*/
254e27a552bSJed Brown PetscErrorCode TSRosWRegisterDestroy(void)
255e27a552bSJed Brown {
256e27a552bSJed Brown   PetscErrorCode ierr;
25761692a83SJed Brown   RosWTableauLink link;
258e27a552bSJed Brown 
259e27a552bSJed Brown   PetscFunctionBegin;
26061692a83SJed Brown   while ((link = RosWTableauList)) {
26161692a83SJed Brown     RosWTableau t = &link->tab;
26261692a83SJed Brown     RosWTableauList = link->next;
26361692a83SJed Brown     ierr = PetscFree5(t->A,t->Gamma,t->b,t->ASum,t->GammaSum);CHKERRQ(ierr);
264c17803e7SJed Brown     ierr = PetscFree4(t->At,t->bt,t->GammaInv,t->GammaZeroDiag);CHKERRQ(ierr);
265fe7e6d57SJed Brown     ierr = PetscFree2(t->bembed,t->bembedt);CHKERRQ(ierr);
266e27a552bSJed Brown     ierr = PetscFree(t->name);CHKERRQ(ierr);
267e27a552bSJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
268e27a552bSJed Brown   }
269e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_FALSE;
270e27a552bSJed Brown   PetscFunctionReturn(0);
271e27a552bSJed Brown }
272e27a552bSJed Brown 
273e27a552bSJed Brown #undef __FUNCT__
274e27a552bSJed Brown #define __FUNCT__ "TSRosWInitializePackage"
275e27a552bSJed Brown /*@C
276e27a552bSJed Brown   TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called
277e27a552bSJed Brown   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RosW()
278e27a552bSJed Brown   when using static libraries.
279e27a552bSJed Brown 
280e27a552bSJed Brown   Input Parameter:
281e27a552bSJed Brown   path - The dynamic library path, or PETSC_NULL
282e27a552bSJed Brown 
283e27a552bSJed Brown   Level: developer
284e27a552bSJed Brown 
285e27a552bSJed Brown .keywords: TS, TSRosW, initialize, package
286e27a552bSJed Brown .seealso: PetscInitialize()
287e27a552bSJed Brown @*/
288e27a552bSJed Brown PetscErrorCode TSRosWInitializePackage(const char path[])
289e27a552bSJed Brown {
290e27a552bSJed Brown   PetscErrorCode ierr;
291e27a552bSJed Brown 
292e27a552bSJed Brown   PetscFunctionBegin;
293e27a552bSJed Brown   if (TSRosWPackageInitialized) PetscFunctionReturn(0);
294e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_TRUE;
295e27a552bSJed Brown   ierr = TSRosWRegisterAll();CHKERRQ(ierr);
296e27a552bSJed Brown   ierr = PetscRegisterFinalize(TSRosWFinalizePackage);CHKERRQ(ierr);
297e27a552bSJed Brown   PetscFunctionReturn(0);
298e27a552bSJed Brown }
299e27a552bSJed Brown 
300e27a552bSJed Brown #undef __FUNCT__
301e27a552bSJed Brown #define __FUNCT__ "TSRosWFinalizePackage"
302e27a552bSJed Brown /*@C
303e27a552bSJed Brown   TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is
304e27a552bSJed Brown   called from PetscFinalize().
305e27a552bSJed Brown 
306e27a552bSJed Brown   Level: developer
307e27a552bSJed Brown 
308e27a552bSJed Brown .keywords: Petsc, destroy, package
309e27a552bSJed Brown .seealso: PetscFinalize()
310e27a552bSJed Brown @*/
311e27a552bSJed Brown PetscErrorCode TSRosWFinalizePackage(void)
312e27a552bSJed Brown {
313e27a552bSJed Brown   PetscErrorCode ierr;
314e27a552bSJed Brown 
315e27a552bSJed Brown   PetscFunctionBegin;
316e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_FALSE;
317e27a552bSJed Brown   ierr = TSRosWRegisterDestroy();CHKERRQ(ierr);
318e27a552bSJed Brown   PetscFunctionReturn(0);
319e27a552bSJed Brown }
320e27a552bSJed Brown 
321e27a552bSJed Brown #undef __FUNCT__
322e27a552bSJed Brown #define __FUNCT__ "TSRosWRegister"
323e27a552bSJed Brown /*@C
32461692a83SJed Brown    TSRosWRegister - register a Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
325e27a552bSJed Brown 
326e27a552bSJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
327e27a552bSJed Brown 
328e27a552bSJed Brown    Input Parameters:
329e27a552bSJed Brown +  name - identifier for method
330e27a552bSJed Brown .  order - approximation order of method
331e27a552bSJed Brown .  s - number of stages, this is the dimension of the matrices below
33261692a83SJed Brown .  A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular
33361692a83SJed Brown .  Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal
334fe7e6d57SJed Brown .  b - Step completion table (dimension s)
335fe7e6d57SJed Brown -  bembed - Step completion table for a scheme of order one less (dimension s, PETSC_NULL if no embedded scheme is available)
336e27a552bSJed Brown 
337e27a552bSJed Brown    Notes:
33861692a83SJed Brown    Several Rosenbrock W methods are provided, this function is only needed to create new methods.
339e27a552bSJed Brown 
340e27a552bSJed Brown    Level: advanced
341e27a552bSJed Brown 
342e27a552bSJed Brown .keywords: TS, register
343e27a552bSJed Brown 
344e27a552bSJed Brown .seealso: TSRosW
345e27a552bSJed Brown @*/
346e27a552bSJed Brown PetscErrorCode TSRosWRegister(const TSRosWType name,PetscInt order,PetscInt s,
347fe7e6d57SJed Brown                               const PetscReal A[],const PetscReal Gamma[],const PetscReal b[],const PetscReal bembed[])
348e27a552bSJed Brown {
349e27a552bSJed Brown   PetscErrorCode ierr;
35061692a83SJed Brown   RosWTableauLink link;
35161692a83SJed Brown   RosWTableau t;
35261692a83SJed Brown   PetscInt i,j,k;
35361692a83SJed Brown   PetscScalar *GammaInv;
354e27a552bSJed Brown 
355e27a552bSJed Brown   PetscFunctionBegin;
356fe7e6d57SJed Brown   PetscValidCharPointer(name,1);
357fe7e6d57SJed Brown   PetscValidPointer(A,4);
358fe7e6d57SJed Brown   PetscValidPointer(Gamma,5);
359fe7e6d57SJed Brown   PetscValidPointer(b,6);
360fe7e6d57SJed Brown   if (bembed) PetscValidPointer(bembed,7);
361fe7e6d57SJed Brown 
362e27a552bSJed Brown   ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
363e27a552bSJed Brown   ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
364e27a552bSJed Brown   t = &link->tab;
365e27a552bSJed Brown   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
366e27a552bSJed Brown   t->order = order;
367e27a552bSJed Brown   t->s = s;
36861692a83SJed 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);
369c17803e7SJed Brown   ierr = PetscMalloc4(s*s,PetscReal,&t->At,s,PetscReal,&t->bt,s*s,PetscReal,&t->GammaInv,s,PetscBool,&t->GammaZeroDiag);CHKERRQ(ierr);
370e27a552bSJed Brown   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
37161692a83SJed Brown   ierr = PetscMemcpy(t->Gamma,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr);
37261692a83SJed Brown   ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
373fe7e6d57SJed Brown   if (bembed) {
374fe7e6d57SJed Brown     ierr = PetscMalloc2(s,PetscReal,&t->bembed,s,PetscReal,&t->bembedt);CHKERRQ(ierr);
375fe7e6d57SJed Brown     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
376fe7e6d57SJed Brown   }
37761692a83SJed Brown   for (i=0; i<s; i++) {
37861692a83SJed Brown     t->ASum[i] = 0;
37961692a83SJed Brown     t->GammaSum[i] = 0;
38061692a83SJed Brown     for (j=0; j<s; j++) {
38161692a83SJed Brown       t->ASum[i] += A[i*s+j];
382fe7e6d57SJed Brown       t->GammaSum[i] += Gamma[i*s+j];
38361692a83SJed Brown     }
38461692a83SJed Brown   }
38561692a83SJed Brown   ierr = PetscMalloc(s*s*sizeof(PetscScalar),&GammaInv);CHKERRQ(ierr); /* Need to use Scalar for inverse, then convert back to Real */
38661692a83SJed Brown   for (i=0; i<s*s; i++) GammaInv[i] = Gamma[i];
387fd96d5b0SEmil Constantinescu   for (i=0; i<s; i++) {
388fd96d5b0SEmil Constantinescu     if (Gamma[i*s+i] == 0.0) {
389fd96d5b0SEmil Constantinescu       GammaInv[i*s+i] = 1.0;
390c17803e7SJed Brown       t->GammaZeroDiag[i] = PETSC_TRUE;
391fd96d5b0SEmil Constantinescu     } else {
392c17803e7SJed Brown       t->GammaZeroDiag[i] = PETSC_FALSE;
393fd96d5b0SEmil Constantinescu     }
394fd96d5b0SEmil Constantinescu   }
395fd96d5b0SEmil Constantinescu 
39661692a83SJed Brown   switch (s) {
39761692a83SJed Brown   case 1: GammaInv[0] = 1./GammaInv[0]; break;
39861692a83SJed Brown   case 2: ierr = Kernel_A_gets_inverse_A_2(GammaInv,0);CHKERRQ(ierr); break;
39961692a83SJed Brown   case 3: ierr = Kernel_A_gets_inverse_A_3(GammaInv,0);CHKERRQ(ierr); break;
40061692a83SJed Brown   case 4: ierr = Kernel_A_gets_inverse_A_4(GammaInv,0);CHKERRQ(ierr); break;
40161692a83SJed Brown   case 5: {
40261692a83SJed Brown     PetscInt ipvt5[5];
40361692a83SJed Brown     MatScalar work5[5*5];
40461692a83SJed Brown     ierr = Kernel_A_gets_inverse_A_5(GammaInv,ipvt5,work5,0);CHKERRQ(ierr); break;
40561692a83SJed Brown   }
40661692a83SJed Brown   case 6: ierr = Kernel_A_gets_inverse_A_6(GammaInv,0);CHKERRQ(ierr); break;
40761692a83SJed Brown   case 7: ierr = Kernel_A_gets_inverse_A_7(GammaInv,0);CHKERRQ(ierr); break;
40861692a83SJed Brown   default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for %D stages",s);
40961692a83SJed Brown   }
41061692a83SJed Brown   for (i=0; i<s*s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]);
41161692a83SJed Brown   ierr = PetscFree(GammaInv);CHKERRQ(ierr);
41261692a83SJed Brown   for (i=0; i<s; i++) {
41361692a83SJed Brown     for (j=0; j<s; j++) {
41461692a83SJed Brown       t->At[i*s+j] = 0;
41561692a83SJed Brown       for (k=0; k<s; k++) {
41661692a83SJed Brown         t->At[i*s+j] += t->A[i*s+k] * t->GammaInv[k*s+j];
41761692a83SJed Brown       }
41861692a83SJed Brown     }
41961692a83SJed Brown     t->bt[i] = 0;
42061692a83SJed Brown     for (j=0; j<s; j++) {
42161692a83SJed Brown       t->bt[i] += t->b[j] * t->GammaInv[j*s+i];
42261692a83SJed Brown     }
423fe7e6d57SJed Brown     if (bembed) {
424fe7e6d57SJed Brown       t->bembedt[i] = 0;
425fe7e6d57SJed Brown       for (j=0; j<s; j++) {
426fe7e6d57SJed Brown         t->bembedt[i] += t->bembed[j] * t->GammaInv[j*s+i];
427fe7e6d57SJed Brown       }
428fe7e6d57SJed Brown     }
42961692a83SJed Brown   }
4308d59e960SJed Brown   t->ccfl = 1.0;                /* Fix this */
4318d59e960SJed Brown 
43261692a83SJed Brown   link->next = RosWTableauList;
43361692a83SJed Brown   RosWTableauList = link;
434e27a552bSJed Brown   PetscFunctionReturn(0);
435e27a552bSJed Brown }
436e27a552bSJed Brown 
437e27a552bSJed Brown #undef __FUNCT__
4381c3436cfSJed Brown #define __FUNCT__ "TSEvaluateStep_RosW"
4391c3436cfSJed Brown /*
4401c3436cfSJed Brown  The step completion formula is
4411c3436cfSJed Brown 
4421c3436cfSJed Brown  x1 = x0 + b^T Y
4431c3436cfSJed Brown 
4441c3436cfSJed Brown  where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been
4451c3436cfSJed Brown  updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write
4461c3436cfSJed Brown 
4471c3436cfSJed Brown  x1e = x0 + be^T Y
4481c3436cfSJed Brown      = x1 - b^T Y + be^T Y
4491c3436cfSJed Brown      = x1 + (be - b)^T Y
4501c3436cfSJed Brown 
4511c3436cfSJed Brown  so we can evaluate the method of different order even after the step has been optimistically completed.
4521c3436cfSJed Brown */
4531c3436cfSJed Brown static PetscErrorCode TSEvaluateStep_RosW(TS ts,PetscInt order,Vec X,PetscBool *done)
4541c3436cfSJed Brown {
4551c3436cfSJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
4561c3436cfSJed Brown   RosWTableau    tab  = ros->tableau;
4571c3436cfSJed Brown   PetscScalar    *w = ros->work;
4581c3436cfSJed Brown   PetscInt       i;
4591c3436cfSJed Brown   PetscErrorCode ierr;
4601c3436cfSJed Brown 
4611c3436cfSJed Brown   PetscFunctionBegin;
4621c3436cfSJed Brown   if (order == tab->order) {
4631c3436cfSJed Brown     if (ros->step_taken) {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
4641c3436cfSJed Brown     else {
4651c3436cfSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
4661c3436cfSJed Brown       ierr = VecMAXPY(X,tab->s,tab->bt,ros->Y);CHKERRQ(ierr);
4671c3436cfSJed Brown     }
4681c3436cfSJed Brown     if (done) *done = PETSC_TRUE;
4691c3436cfSJed Brown     PetscFunctionReturn(0);
4701c3436cfSJed Brown   } else if (order == tab->order-1) {
4711c3436cfSJed Brown     if (!tab->bembedt) goto unavailable;
4721c3436cfSJed Brown     if (ros->step_taken) {
4731c3436cfSJed Brown       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i];
4741c3436cfSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
4751c3436cfSJed Brown       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
4761c3436cfSJed Brown     } else {
4771c3436cfSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
4781c3436cfSJed Brown       ierr = VecMAXPY(X,tab->s,tab->bembedt,ros->Y);CHKERRQ(ierr);
4791c3436cfSJed Brown     }
4801c3436cfSJed Brown     if (done) *done = PETSC_TRUE;
4811c3436cfSJed Brown     PetscFunctionReturn(0);
4821c3436cfSJed Brown   }
4831c3436cfSJed Brown   unavailable:
4841c3436cfSJed Brown   if (done) *done = PETSC_FALSE;
4851c3436cfSJed 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);
4861c3436cfSJed Brown   PetscFunctionReturn(0);
4871c3436cfSJed Brown }
4881c3436cfSJed Brown 
4891c3436cfSJed Brown #undef __FUNCT__
490e27a552bSJed Brown #define __FUNCT__ "TSStep_RosW"
491e27a552bSJed Brown static PetscErrorCode TSStep_RosW(TS ts)
492e27a552bSJed Brown {
49361692a83SJed Brown   TS_RosW         *ros = (TS_RosW*)ts->data;
49461692a83SJed Brown   RosWTableau     tab  = ros->tableau;
495e27a552bSJed Brown   const PetscInt  s    = tab->s;
4961c3436cfSJed Brown   const PetscReal *At  = tab->At,*Gamma = tab->Gamma,*ASum = tab->ASum,*GammaInv = tab->GammaInv;
497c17803e7SJed Brown   const PetscBool *GammaZeroDiag = tab->GammaZeroDiag;
49861692a83SJed Brown   PetscScalar     *w   = ros->work;
49961692a83SJed Brown   Vec             *Y   = ros->Y,Zdot = ros->Zdot,Zstage = ros->Zstage;
500e27a552bSJed Brown   SNES            snes;
5011c3436cfSJed Brown   TSAdapt         adapt;
5021c3436cfSJed Brown   PetscInt        i,j,its,lits,reject,next_scheme;
503cdbf8f93SLisandro Dalcin   PetscReal       next_time_step;
5041c3436cfSJed Brown   PetscBool       accept;
505e27a552bSJed Brown   PetscErrorCode  ierr;
506e27a552bSJed Brown 
507e27a552bSJed Brown   PetscFunctionBegin;
508e27a552bSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
509cdbf8f93SLisandro Dalcin   next_time_step = ts->time_step;
5101c3436cfSJed Brown   accept = PETSC_TRUE;
5111c3436cfSJed Brown   ros->step_taken = PETSC_FALSE;
512e27a552bSJed Brown 
5131c3436cfSJed Brown   for (reject=0; reject<ts->max_reject; reject++,ts->reject++) {
5141c3436cfSJed Brown     const PetscReal h = ts->time_step;
515e27a552bSJed Brown     for (i=0; i<s; i++) {
5161c3436cfSJed Brown       ros->stage_time = ts->ptime + h*ASum[i];
517c17803e7SJed Brown       if (GammaZeroDiag[i]) {
518c17803e7SJed Brown         ros->stage_explicit = PETSC_TRUE;
519fd96d5b0SEmil Constantinescu         ros->shift = 1./h;
520c17803e7SJed Brown       } else {
521c17803e7SJed Brown         ros->stage_explicit = PETSC_FALSE;
52261692a83SJed Brown         ros->shift = 1./(h*Gamma[i*s+i]);
523fd96d5b0SEmil Constantinescu       }
52461692a83SJed Brown 
52561692a83SJed Brown       ierr = VecCopy(ts->vec_sol,Zstage);CHKERRQ(ierr);
52661692a83SJed Brown       ierr = VecMAXPY(Zstage,i,&At[i*s+0],Y);CHKERRQ(ierr);
52761692a83SJed Brown 
52861692a83SJed Brown       for (j=0; j<i; j++) w[j] = 1./h * GammaInv[i*s+j];
52961692a83SJed Brown       ierr = VecZeroEntries(Zdot);CHKERRQ(ierr);
53061692a83SJed Brown       ierr = VecMAXPY(Zdot,i,w,Y);CHKERRQ(ierr);
53161692a83SJed Brown 
532e27a552bSJed Brown       /* Initial guess taken from last stage */
53361692a83SJed Brown       ierr = VecZeroEntries(Y[i]);CHKERRQ(ierr);
53461692a83SJed Brown 
53561692a83SJed Brown       if (!ros->recompute_jacobian && !i) {
53661692a83SJed Brown         ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */
53761692a83SJed Brown       }
53861692a83SJed Brown 
53961692a83SJed Brown       ierr = SNESSolve(snes,PETSC_NULL,Y[i]);CHKERRQ(ierr);
540e27a552bSJed Brown       ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
541e27a552bSJed Brown       ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
542e27a552bSJed Brown       ts->nonlinear_its += its; ts->linear_its += lits;
543e27a552bSJed Brown     }
5441c3436cfSJed Brown     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
5451c3436cfSJed Brown     ros->step_taken = PETSC_TRUE;
546e27a552bSJed Brown 
5471c3436cfSJed Brown     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
5481c3436cfSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
5491c3436cfSJed Brown     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
5508d59e960SJed Brown     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
5511c3436cfSJed Brown     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
5521c3436cfSJed Brown     if (accept) {
5531c3436cfSJed Brown       /* ignore next_scheme for now */
554e27a552bSJed Brown       ts->ptime += ts->time_step;
555cdbf8f93SLisandro Dalcin       ts->time_step = next_time_step;
556e27a552bSJed Brown       ts->steps++;
5571c3436cfSJed Brown       break;
5581c3436cfSJed Brown     } else {                    /* Roll back the current step */
5591c3436cfSJed Brown       for (i=0; i<s; i++) w[i] = -tab->bt[i];
5601c3436cfSJed Brown       ierr = VecMAXPY(ts->vec_sol,s,w,Y);CHKERRQ(ierr);
5611c3436cfSJed Brown       ts->time_step = next_time_step;
5621c3436cfSJed Brown       ros->step_taken = PETSC_FALSE;
5631c3436cfSJed Brown     }
5641c3436cfSJed Brown   }
5651c3436cfSJed Brown 
566e27a552bSJed Brown   PetscFunctionReturn(0);
567e27a552bSJed Brown }
568e27a552bSJed Brown 
569e27a552bSJed Brown #undef __FUNCT__
570e27a552bSJed Brown #define __FUNCT__ "TSInterpolate_RosW"
571e27a552bSJed Brown static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec X)
572e27a552bSJed Brown {
57361692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
574e27a552bSJed Brown 
575e27a552bSJed Brown   PetscFunctionBegin;
57661692a83SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name);
577e27a552bSJed Brown   PetscFunctionReturn(0);
578e27a552bSJed Brown }
579e27a552bSJed Brown 
580e27a552bSJed Brown /*------------------------------------------------------------*/
581e27a552bSJed Brown #undef __FUNCT__
582e27a552bSJed Brown #define __FUNCT__ "TSReset_RosW"
583e27a552bSJed Brown static PetscErrorCode TSReset_RosW(TS ts)
584e27a552bSJed Brown {
58561692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
586e27a552bSJed Brown   PetscInt       s;
587e27a552bSJed Brown   PetscErrorCode ierr;
588e27a552bSJed Brown 
589e27a552bSJed Brown   PetscFunctionBegin;
59061692a83SJed Brown   if (!ros->tableau) PetscFunctionReturn(0);
59161692a83SJed Brown    s = ros->tableau->s;
59261692a83SJed Brown   ierr = VecDestroyVecs(s,&ros->Y);CHKERRQ(ierr);
59361692a83SJed Brown   ierr = VecDestroy(&ros->Ydot);CHKERRQ(ierr);
59461692a83SJed Brown   ierr = VecDestroy(&ros->Ystage);CHKERRQ(ierr);
59561692a83SJed Brown   ierr = VecDestroy(&ros->Zdot);CHKERRQ(ierr);
59661692a83SJed Brown   ierr = VecDestroy(&ros->Zstage);CHKERRQ(ierr);
59761692a83SJed Brown   ierr = PetscFree(ros->work);CHKERRQ(ierr);
598e27a552bSJed Brown   PetscFunctionReturn(0);
599e27a552bSJed Brown }
600e27a552bSJed Brown 
601e27a552bSJed Brown #undef __FUNCT__
602e27a552bSJed Brown #define __FUNCT__ "TSDestroy_RosW"
603e27a552bSJed Brown static PetscErrorCode TSDestroy_RosW(TS ts)
604e27a552bSJed Brown {
605e27a552bSJed Brown   PetscErrorCode  ierr;
606e27a552bSJed Brown 
607e27a552bSJed Brown   PetscFunctionBegin;
608e27a552bSJed Brown   ierr = TSReset_RosW(ts);CHKERRQ(ierr);
609e27a552bSJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
610e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","",PETSC_NULL);CHKERRQ(ierr);
611e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","",PETSC_NULL);CHKERRQ(ierr);
61261692a83SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","",PETSC_NULL);CHKERRQ(ierr);
613e27a552bSJed Brown   PetscFunctionReturn(0);
614e27a552bSJed Brown }
615e27a552bSJed Brown 
616e27a552bSJed Brown /*
617e27a552bSJed Brown   This defines the nonlinear equation that is to be solved with SNES
618e27a552bSJed Brown   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
619e27a552bSJed Brown */
620e27a552bSJed Brown #undef __FUNCT__
621e27a552bSJed Brown #define __FUNCT__ "SNESTSFormFunction_RosW"
622e27a552bSJed Brown static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec X,Vec F,TS ts)
623e27a552bSJed Brown {
62461692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
625e27a552bSJed Brown   PetscErrorCode ierr;
626e27a552bSJed Brown 
627e27a552bSJed Brown   PetscFunctionBegin;
628c17803e7SJed Brown   if (ros->stage_explicit) {
629c17803e7SJed Brown     ierr = VecAXPBY(ros->Ydot,ros->shift,0.0,X);CHKERRQ(ierr);       /* Ydot = shift*X*/
630c17803e7SJed Brown   } else {
63161692a83SJed Brown     ierr = VecWAXPY(ros->Ydot,ros->shift,X,ros->Zdot);CHKERRQ(ierr); /* Ydot = shift*X + Zdot */
632c17803e7SJed Brown   }
63361692a83SJed Brown   ierr = VecWAXPY(ros->Ystage,1.0,X,ros->Zstage);CHKERRQ(ierr);    /* Ystage = X + Zstage */
63461692a83SJed Brown   ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,ros->Ydot,F,PETSC_FALSE);CHKERRQ(ierr);
635e27a552bSJed Brown   PetscFunctionReturn(0);
636e27a552bSJed Brown }
637e27a552bSJed Brown 
638e27a552bSJed Brown #undef __FUNCT__
639e27a552bSJed Brown #define __FUNCT__ "SNESTSFormJacobian_RosW"
640e27a552bSJed Brown static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
641e27a552bSJed Brown {
64261692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
643e27a552bSJed Brown   PetscErrorCode ierr;
644e27a552bSJed Brown 
645e27a552bSJed Brown   PetscFunctionBegin;
64661692a83SJed Brown   /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
64761692a83SJed Brown   ierr = TSComputeIJacobian(ts,ros->stage_time,ros->Ystage,ros->Ydot,ros->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr);
648e27a552bSJed Brown   PetscFunctionReturn(0);
649e27a552bSJed Brown }
650e27a552bSJed Brown 
651e27a552bSJed Brown #undef __FUNCT__
652e27a552bSJed Brown #define __FUNCT__ "TSSetUp_RosW"
653e27a552bSJed Brown static PetscErrorCode TSSetUp_RosW(TS ts)
654e27a552bSJed Brown {
65561692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
65661692a83SJed Brown   RosWTableau    tab  = ros->tableau;
657e27a552bSJed Brown   PetscInt       s    = tab->s;
658e27a552bSJed Brown   PetscErrorCode ierr;
659e27a552bSJed Brown 
660e27a552bSJed Brown   PetscFunctionBegin;
66161692a83SJed Brown   if (!ros->tableau) {
662e27a552bSJed Brown     ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);
663e27a552bSJed Brown   }
66461692a83SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ros->Y);CHKERRQ(ierr);
66561692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Ydot);CHKERRQ(ierr);
66661692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Ystage);CHKERRQ(ierr);
66761692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Zdot);CHKERRQ(ierr);
66861692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Zstage);CHKERRQ(ierr);
66961692a83SJed Brown   ierr = PetscMalloc(s*sizeof(ros->work[0]),&ros->work);CHKERRQ(ierr);
670e27a552bSJed Brown   PetscFunctionReturn(0);
671e27a552bSJed Brown }
672e27a552bSJed Brown /*------------------------------------------------------------*/
673e27a552bSJed Brown 
674e27a552bSJed Brown #undef __FUNCT__
675e27a552bSJed Brown #define __FUNCT__ "TSSetFromOptions_RosW"
676e27a552bSJed Brown static PetscErrorCode TSSetFromOptions_RosW(TS ts)
677e27a552bSJed Brown {
67861692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
679e27a552bSJed Brown   PetscErrorCode ierr;
68061692a83SJed Brown   char           rostype[256];
681e27a552bSJed Brown 
682e27a552bSJed Brown   PetscFunctionBegin;
683e27a552bSJed Brown   ierr = PetscOptionsHead("RosW ODE solver options");CHKERRQ(ierr);
684e27a552bSJed Brown   {
68561692a83SJed Brown     RosWTableauLink link;
686e27a552bSJed Brown     PetscInt count,choice;
687e27a552bSJed Brown     PetscBool flg;
688e27a552bSJed Brown     const char **namelist;
68961692a83SJed Brown     SNES snes;
69061692a83SJed Brown 
69161692a83SJed Brown     ierr = PetscStrncpy(rostype,TSRosWDefault,sizeof rostype);CHKERRQ(ierr);
69261692a83SJed Brown     for (link=RosWTableauList,count=0; link; link=link->next,count++) ;
693e27a552bSJed Brown     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
69461692a83SJed Brown     for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
69561692a83SJed Brown     ierr = PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,rostype,&choice,&flg);CHKERRQ(ierr);
69661692a83SJed Brown     ierr = TSRosWSetType(ts,flg ? namelist[choice] : rostype);CHKERRQ(ierr);
697e27a552bSJed Brown     ierr = PetscFree(namelist);CHKERRQ(ierr);
69861692a83SJed Brown 
69961692a83SJed Brown     ierr = PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,PETSC_NULL);CHKERRQ(ierr);
70061692a83SJed Brown 
70161692a83SJed Brown     /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
70261692a83SJed Brown     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
70361692a83SJed Brown     if (!((PetscObject)snes)->type_name) {
70461692a83SJed Brown       ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
70561692a83SJed Brown     }
70661692a83SJed Brown     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
707e27a552bSJed Brown   }
708e27a552bSJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
709e27a552bSJed Brown   PetscFunctionReturn(0);
710e27a552bSJed Brown }
711e27a552bSJed Brown 
712e27a552bSJed Brown #undef __FUNCT__
713e27a552bSJed Brown #define __FUNCT__ "PetscFormatRealArray"
714e27a552bSJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
715e27a552bSJed Brown {
716e27a552bSJed Brown   PetscErrorCode ierr;
717e408995aSJed Brown   PetscInt i;
718e408995aSJed Brown   size_t left,count;
719e27a552bSJed Brown   char *p;
720e27a552bSJed Brown 
721e27a552bSJed Brown   PetscFunctionBegin;
722e408995aSJed Brown   for (i=0,p=buf,left=len; i<n; i++) {
723e408995aSJed Brown     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
724e27a552bSJed Brown     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
725e27a552bSJed Brown     left -= count;
726e27a552bSJed Brown     p += count;
727e27a552bSJed Brown     *p++ = ' ';
728e27a552bSJed Brown   }
729e27a552bSJed Brown   p[i ? 0 : -1] = 0;
730e27a552bSJed Brown   PetscFunctionReturn(0);
731e27a552bSJed Brown }
732e27a552bSJed Brown 
733e27a552bSJed Brown #undef __FUNCT__
734e27a552bSJed Brown #define __FUNCT__ "TSView_RosW"
735e27a552bSJed Brown static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer)
736e27a552bSJed Brown {
73761692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
73861692a83SJed Brown   RosWTableau    tab  = ros->tableau;
739e27a552bSJed Brown   PetscBool      iascii;
740e27a552bSJed Brown   PetscErrorCode ierr;
741e27a552bSJed Brown 
742e27a552bSJed Brown   PetscFunctionBegin;
743e27a552bSJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
744e27a552bSJed Brown   if (iascii) {
74561692a83SJed Brown     const TSRosWType rostype;
746e408995aSJed Brown     PetscInt i;
747e408995aSJed Brown     PetscReal abscissa[512];
748e27a552bSJed Brown     char buf[512];
74961692a83SJed Brown     ierr = TSRosWGetType(ts,&rostype);CHKERRQ(ierr);
75061692a83SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Rosenbrock-W %s\n",rostype);CHKERRQ(ierr);
751e408995aSJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ASum);CHKERRQ(ierr);
75261692a83SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A       = %s\n",buf);CHKERRQ(ierr);
753e408995aSJed Brown     for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i];
754e408995aSJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,abscissa);CHKERRQ(ierr);
755e408995aSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A+Gamma = %s\n",buf);CHKERRQ(ierr);
756e27a552bSJed Brown   }
757e27a552bSJed Brown   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
758e27a552bSJed Brown   PetscFunctionReturn(0);
759e27a552bSJed Brown }
760e27a552bSJed Brown 
761e27a552bSJed Brown #undef __FUNCT__
762e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType"
763e27a552bSJed Brown /*@C
76461692a83SJed Brown   TSRosWSetType - Set the type of Rosenbrock-W scheme
765e27a552bSJed Brown 
766e27a552bSJed Brown   Logically collective
767e27a552bSJed Brown 
768e27a552bSJed Brown   Input Parameter:
769e27a552bSJed Brown +  ts - timestepping context
77061692a83SJed Brown -  rostype - type of Rosenbrock-W scheme
771e27a552bSJed Brown 
772e27a552bSJed Brown   Level: intermediate
773e27a552bSJed Brown 
774e27a552bSJed Brown .seealso: TSRosWGetType()
775e27a552bSJed Brown @*/
77661692a83SJed Brown PetscErrorCode TSRosWSetType(TS ts,const TSRosWType rostype)
777e27a552bSJed Brown {
778e27a552bSJed Brown   PetscErrorCode ierr;
779e27a552bSJed Brown 
780e27a552bSJed Brown   PetscFunctionBegin;
781e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
78261692a83SJed Brown   ierr = PetscTryMethod(ts,"TSRosWSetType_C",(TS,const TSRosWType),(ts,rostype));CHKERRQ(ierr);
783e27a552bSJed Brown   PetscFunctionReturn(0);
784e27a552bSJed Brown }
785e27a552bSJed Brown 
786e27a552bSJed Brown #undef __FUNCT__
787e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType"
788e27a552bSJed Brown /*@C
78961692a83SJed Brown   TSRosWGetType - Get the type of Rosenbrock-W scheme
790e27a552bSJed Brown 
791e27a552bSJed Brown   Logically collective
792e27a552bSJed Brown 
793e27a552bSJed Brown   Input Parameter:
794e27a552bSJed Brown .  ts - timestepping context
795e27a552bSJed Brown 
796e27a552bSJed Brown   Output Parameter:
79761692a83SJed Brown .  rostype - type of Rosenbrock-W scheme
798e27a552bSJed Brown 
799e27a552bSJed Brown   Level: intermediate
800e27a552bSJed Brown 
801e27a552bSJed Brown .seealso: TSRosWGetType()
802e27a552bSJed Brown @*/
80361692a83SJed Brown PetscErrorCode TSRosWGetType(TS ts,const TSRosWType *rostype)
804e27a552bSJed Brown {
805e27a552bSJed Brown   PetscErrorCode ierr;
806e27a552bSJed Brown 
807e27a552bSJed Brown   PetscFunctionBegin;
808e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
80961692a83SJed Brown   ierr = PetscUseMethod(ts,"TSRosWGetType_C",(TS,const TSRosWType*),(ts,rostype));CHKERRQ(ierr);
810e27a552bSJed Brown   PetscFunctionReturn(0);
811e27a552bSJed Brown }
812e27a552bSJed Brown 
813e27a552bSJed Brown #undef __FUNCT__
81461692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian"
815e27a552bSJed Brown /*@C
81661692a83SJed Brown   TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step.
817e27a552bSJed Brown 
818e27a552bSJed Brown   Logically collective
819e27a552bSJed Brown 
820e27a552bSJed Brown   Input Parameter:
821e27a552bSJed Brown +  ts - timestepping context
82261692a83SJed Brown -  flg - PETSC_TRUE to recompute the Jacobian at each stage
823e27a552bSJed Brown 
824e27a552bSJed Brown   Level: intermediate
825e27a552bSJed Brown 
826e27a552bSJed Brown .seealso: TSRosWGetType()
827e27a552bSJed Brown @*/
82861692a83SJed Brown PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg)
829e27a552bSJed Brown {
830e27a552bSJed Brown   PetscErrorCode ierr;
831e27a552bSJed Brown 
832e27a552bSJed Brown   PetscFunctionBegin;
833e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
83461692a83SJed Brown   ierr = PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
835e27a552bSJed Brown   PetscFunctionReturn(0);
836e27a552bSJed Brown }
837e27a552bSJed Brown 
838e27a552bSJed Brown EXTERN_C_BEGIN
839e27a552bSJed Brown #undef __FUNCT__
840e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType_RosW"
84161692a83SJed Brown PetscErrorCode  TSRosWGetType_RosW(TS ts,const TSRosWType *rostype)
842e27a552bSJed Brown {
84361692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
844e27a552bSJed Brown   PetscErrorCode ierr;
845e27a552bSJed Brown 
846e27a552bSJed Brown   PetscFunctionBegin;
84761692a83SJed Brown   if (!ros->tableau) {ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);}
84861692a83SJed Brown   *rostype = ros->tableau->name;
849e27a552bSJed Brown   PetscFunctionReturn(0);
850e27a552bSJed Brown }
851e27a552bSJed Brown #undef __FUNCT__
852e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType_RosW"
85361692a83SJed Brown PetscErrorCode  TSRosWSetType_RosW(TS ts,const TSRosWType rostype)
854e27a552bSJed Brown {
85561692a83SJed Brown   TS_RosW         *ros = (TS_RosW*)ts->data;
856e27a552bSJed Brown   PetscErrorCode  ierr;
857e27a552bSJed Brown   PetscBool       match;
85861692a83SJed Brown   RosWTableauLink link;
859e27a552bSJed Brown 
860e27a552bSJed Brown   PetscFunctionBegin;
86161692a83SJed Brown   if (ros->tableau) {
86261692a83SJed Brown     ierr = PetscStrcmp(ros->tableau->name,rostype,&match);CHKERRQ(ierr);
863e27a552bSJed Brown     if (match) PetscFunctionReturn(0);
864e27a552bSJed Brown   }
86561692a83SJed Brown   for (link = RosWTableauList; link; link=link->next) {
86661692a83SJed Brown     ierr = PetscStrcmp(link->tab.name,rostype,&match);CHKERRQ(ierr);
867e27a552bSJed Brown     if (match) {
868e27a552bSJed Brown       ierr = TSReset_RosW(ts);CHKERRQ(ierr);
86961692a83SJed Brown       ros->tableau = &link->tab;
870e27a552bSJed Brown       PetscFunctionReturn(0);
871e27a552bSJed Brown     }
872e27a552bSJed Brown   }
87361692a83SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype);
874e27a552bSJed Brown   PetscFunctionReturn(0);
875e27a552bSJed Brown }
87661692a83SJed Brown 
877e27a552bSJed Brown #undef __FUNCT__
87861692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian_RosW"
87961692a83SJed Brown PetscErrorCode  TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg)
880e27a552bSJed Brown {
88161692a83SJed Brown   TS_RosW *ros = (TS_RosW*)ts->data;
882e27a552bSJed Brown 
883e27a552bSJed Brown   PetscFunctionBegin;
88461692a83SJed Brown   ros->recompute_jacobian = flg;
885e27a552bSJed Brown   PetscFunctionReturn(0);
886e27a552bSJed Brown }
887e27a552bSJed Brown EXTERN_C_END
888e27a552bSJed Brown 
889e27a552bSJed Brown /* ------------------------------------------------------------ */
890e27a552bSJed Brown /*MC
891e27a552bSJed Brown       TSRosW - ODE solver using Rosenbrock-W schemes
892e27a552bSJed Brown 
893e27a552bSJed Brown   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
894e27a552bSJed Brown   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
895e27a552bSJed Brown   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
896e27a552bSJed Brown 
897e27a552bSJed Brown   Notes:
89861692a83SJed Brown   This method currently only works with autonomous ODE and DAE.
89961692a83SJed Brown 
90061692a83SJed Brown   Developer notes:
90161692a83SJed Brown   Rosenbrock-W methods are typically specified for autonomous ODE
90261692a83SJed Brown 
90361692a83SJed Brown $  xdot = f(x)
90461692a83SJed Brown 
90561692a83SJed Brown   by the stage equations
90661692a83SJed Brown 
90761692a83SJed Brown $  k_i = h f(x_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j
90861692a83SJed Brown 
90961692a83SJed Brown   and step completion formula
91061692a83SJed Brown 
91161692a83SJed Brown $  x_1 = x_0 + sum_j b_j k_j
91261692a83SJed Brown 
91361692a83SJed Brown   with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(x)
91461692a83SJed Brown   and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner,
91561692a83SJed Brown   we define new variables for the stage equations
91661692a83SJed Brown 
91761692a83SJed Brown $  y_i = gamma_ij k_j
91861692a83SJed Brown 
91961692a83SJed Brown   The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define
92061692a83SJed Brown 
92161692a83SJed Brown $  A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-i}
92261692a83SJed Brown 
92361692a83SJed Brown   to rewrite the method as
92461692a83SJed Brown 
92561692a83SJed 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
92661692a83SJed Brown $  x_1 = x_0 + sum_j bt_j y_j
92761692a83SJed Brown 
92861692a83SJed Brown    where we have introduced the mass matrix M. Continue by defining
92961692a83SJed Brown 
93061692a83SJed Brown $  ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j
93161692a83SJed Brown 
93261692a83SJed Brown    or, more compactly in tensor notation
93361692a83SJed Brown 
93461692a83SJed Brown $  Ydot = 1/h (Gamma^{-1} \otimes I) Y .
93561692a83SJed Brown 
93661692a83SJed Brown    Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current
93761692a83SJed Brown    stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the
93861692a83SJed Brown    equation
93961692a83SJed Brown 
94061692a83SJed Brown $  g(x_0 + sum_j a_ij y_j + y_i, ydot_i) = 0
94161692a83SJed Brown 
94261692a83SJed Brown    with initial guess y_i = 0.
943e27a552bSJed Brown 
944e27a552bSJed Brown   Level: beginner
945e27a552bSJed Brown 
946e27a552bSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSRosWRegister()
947e27a552bSJed Brown 
948e27a552bSJed Brown M*/
949e27a552bSJed Brown EXTERN_C_BEGIN
950e27a552bSJed Brown #undef __FUNCT__
951e27a552bSJed Brown #define __FUNCT__ "TSCreate_RosW"
952e27a552bSJed Brown PetscErrorCode  TSCreate_RosW(TS ts)
953e27a552bSJed Brown {
95461692a83SJed Brown   TS_RosW        *ros;
955e27a552bSJed Brown   PetscErrorCode ierr;
956e27a552bSJed Brown 
957e27a552bSJed Brown   PetscFunctionBegin;
958e27a552bSJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
959e27a552bSJed Brown   ierr = TSRosWInitializePackage(PETSC_NULL);CHKERRQ(ierr);
960e27a552bSJed Brown #endif
961e27a552bSJed Brown 
962e27a552bSJed Brown   ts->ops->reset          = TSReset_RosW;
963e27a552bSJed Brown   ts->ops->destroy        = TSDestroy_RosW;
964e27a552bSJed Brown   ts->ops->view           = TSView_RosW;
965e27a552bSJed Brown   ts->ops->setup          = TSSetUp_RosW;
966e27a552bSJed Brown   ts->ops->step           = TSStep_RosW;
967e27a552bSJed Brown   ts->ops->interpolate    = TSInterpolate_RosW;
9681c3436cfSJed Brown   ts->ops->evaluatestep   = TSEvaluateStep_RosW;
969e27a552bSJed Brown   ts->ops->setfromoptions = TSSetFromOptions_RosW;
970e27a552bSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_RosW;
971e27a552bSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_RosW;
972e27a552bSJed Brown 
97361692a83SJed Brown   ierr = PetscNewLog(ts,TS_RosW,&ros);CHKERRQ(ierr);
97461692a83SJed Brown   ts->data = (void*)ros;
975e27a552bSJed Brown 
976e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","TSRosWGetType_RosW",TSRosWGetType_RosW);CHKERRQ(ierr);
977e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","TSRosWSetType_RosW",TSRosWSetType_RosW);CHKERRQ(ierr);
97861692a83SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","TSRosWSetRecomputeJacobian_RosW",TSRosWSetRecomputeJacobian_RosW);CHKERRQ(ierr);
979e27a552bSJed Brown   PetscFunctionReturn(0);
980e27a552bSJed Brown }
981e27a552bSJed Brown EXTERN_C_END
982