xref: /petsc/src/ts/impls/rosw/rosw.c (revision fe7e6d57163471cd46fdc637f0776ac5f202504b)
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 */
2861692a83SJed Brown   PetscReal *b;                 /* Step completion table */
29*fe7e6d57SJed Brown   PetscReal *bembed;            /* Step completion table for embedded method of order one less */
3061692a83SJed Brown   PetscReal *ASum;              /* Row sum of A */
3161692a83SJed Brown   PetscReal *GammaSum;          /* Row sum of Gamma, only needed for non-autonomous systems */
3261692a83SJed Brown   PetscReal *At;                /* Propagation table in transformed variables */
3361692a83SJed Brown   PetscReal *bt;                /* Step completion table in transformed variables */
34*fe7e6d57SJed Brown   PetscReal *bembedt;           /* Step completion table of order one less in transformed variables */
3561692a83SJed Brown   PetscReal *GammaInv;          /* Inverse of Gamma, used for transformed variables */
36e27a552bSJed Brown };
3761692a83SJed Brown typedef struct _RosWTableauLink *RosWTableauLink;
3861692a83SJed Brown struct _RosWTableauLink {
3961692a83SJed Brown   struct _RosWTableau tab;
4061692a83SJed Brown   RosWTableauLink next;
41e27a552bSJed Brown };
4261692a83SJed Brown static RosWTableauLink RosWTableauList;
43e27a552bSJed Brown 
44e27a552bSJed Brown typedef struct {
4561692a83SJed Brown   RosWTableau tableau;
4661692a83SJed Brown   Vec         *Y;               /* States computed during the step, used to complete the step */
47e27a552bSJed Brown   Vec         Ydot;             /* Work vector holding Ydot during residual evaluation */
4861692a83SJed Brown   Vec         Ystage;           /* Work vector for the state value at each stage */
4961692a83SJed Brown   Vec         Zdot;             /* Ydot = Zdot + shift*Y */
5061692a83SJed Brown   Vec         Zstage;           /* Y = Zstage + Y */
51e27a552bSJed Brown   PetscScalar *work;            /* Scalar work */
52e27a552bSJed Brown   PetscReal   shift;
53e27a552bSJed Brown   PetscReal   stage_time;
5461692a83SJed Brown   PetscBool   recompute_jacobian; /* Recompute the Jacobian at each stage, default is to freeze the Jacobian at the start of each step */
55e27a552bSJed Brown } TS_RosW;
56e27a552bSJed Brown 
57*fe7e6d57SJed Brown /*MC
58*fe7e6d57SJed Brown      TSROSW2M - Two stage second order L-stable Rosenbrock-W scheme.
59*fe7e6d57SJed Brown 
60*fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2P.
61*fe7e6d57SJed Brown 
62*fe7e6d57SJed Brown      Level: intermediate
63*fe7e6d57SJed Brown 
64*fe7e6d57SJed Brown .seealso: TSROSW
65*fe7e6d57SJed Brown M*/
66*fe7e6d57SJed Brown 
67*fe7e6d57SJed Brown /*MC
68*fe7e6d57SJed Brown      TSROSW2P - Two stage second order L-stable Rosenbrock-W scheme.
69*fe7e6d57SJed Brown 
70*fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2M.
71*fe7e6d57SJed Brown 
72*fe7e6d57SJed Brown      Level: intermediate
73*fe7e6d57SJed Brown 
74*fe7e6d57SJed Brown .seealso: TSROSW
75*fe7e6d57SJed Brown M*/
76*fe7e6d57SJed Brown 
77*fe7e6d57SJed Brown /*MC
78*fe7e6d57SJed Brown      TSROSWRA3PW - Three stage third order Rosenbrock-W scheme for PDAE of index 1.
79*fe7e6d57SJed Brown 
80*fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.
81*fe7e6d57SJed Brown 
82*fe7e6d57SJed 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.
83*fe7e6d57SJed Brown 
84*fe7e6d57SJed Brown      References:
85*fe7e6d57SJed Brown      Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005.
86*fe7e6d57SJed Brown 
87*fe7e6d57SJed Brown      Level: intermediate
88*fe7e6d57SJed Brown 
89*fe7e6d57SJed Brown .seealso: TSROSW
90*fe7e6d57SJed Brown M*/
91*fe7e6d57SJed Brown 
92*fe7e6d57SJed Brown /*MC
93*fe7e6d57SJed Brown      TSROSWRA34PW2 - Four stage third order L-stable Rosenbrock-W scheme for PDAE of index 1.
94*fe7e6d57SJed Brown 
95*fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.
96*fe7e6d57SJed Brown 
97*fe7e6d57SJed 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.
98*fe7e6d57SJed Brown 
99*fe7e6d57SJed Brown      References:
100*fe7e6d57SJed Brown      Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005.
101*fe7e6d57SJed Brown 
102*fe7e6d57SJed Brown      Level: intermediate
103*fe7e6d57SJed Brown 
104*fe7e6d57SJed Brown .seealso: TSROSW
105*fe7e6d57SJed Brown M*/
106*fe7e6d57SJed Brown 
107e27a552bSJed Brown #undef __FUNCT__
108e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterAll"
109e27a552bSJed Brown /*@C
110e27a552bSJed Brown   TSRosWRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSRosW
111e27a552bSJed Brown 
112e27a552bSJed Brown   Not Collective, but should be called by all processes which will need the schemes to be registered
113e27a552bSJed Brown 
114e27a552bSJed Brown   Level: advanced
115e27a552bSJed Brown 
116e27a552bSJed Brown .keywords: TS, TSRosW, register, all
117e27a552bSJed Brown 
118e27a552bSJed Brown .seealso:  TSRosWRegisterDestroy()
119e27a552bSJed Brown @*/
120e27a552bSJed Brown PetscErrorCode TSRosWRegisterAll(void)
121e27a552bSJed Brown {
122e27a552bSJed Brown   PetscErrorCode ierr;
123e27a552bSJed Brown 
124e27a552bSJed Brown   PetscFunctionBegin;
125e27a552bSJed Brown   if (TSRosWRegisterAllCalled) PetscFunctionReturn(0);
126e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_TRUE;
127e27a552bSJed Brown 
128e27a552bSJed Brown   {
12961692a83SJed Brown     const PetscReal g = 1. + 1./PetscSqrtReal(2.0);
130e27a552bSJed Brown     const PetscReal
13161692a83SJed Brown       A[2][2] = {{0,0}, {1.,0}},
13261692a83SJed Brown       Gamma[2][2] = {{g,0}, {-2.*g,g}},
13361692a83SJed Brown       b[2] = {0.5,0.5};
134*fe7e6d57SJed Brown     ierr = TSRosWRegister(TSROSW2P,2,2,&A[0][0],&Gamma[0][0],b,PETSC_NULL);CHKERRQ(ierr);
135e27a552bSJed Brown   }
136e27a552bSJed Brown   {
13761692a83SJed Brown     const PetscReal g = 1. - 1./PetscSqrtReal(2.0);
138e27a552bSJed Brown     const PetscReal
13961692a83SJed Brown       A[2][2] = {{0,0}, {1.,0}},
14061692a83SJed Brown       Gamma[2][2] = {{g,0}, {-2.*g,g}},
14161692a83SJed Brown       b[2] = {0.5,0.5};
142*fe7e6d57SJed Brown     ierr = TSRosWRegister(TSROSW2M,2,2,&A[0][0],&Gamma[0][0],b,PETSC_NULL);CHKERRQ(ierr);
143*fe7e6d57SJed Brown   }
144*fe7e6d57SJed Brown   {
145*fe7e6d57SJed Brown     const PetscReal g = 7.8867513459481287e-01;
146*fe7e6d57SJed Brown     const PetscReal
147*fe7e6d57SJed Brown       A[3][3] = {{0,0,0},
148*fe7e6d57SJed Brown                  {1.5773502691896257e+00,0,0},
149*fe7e6d57SJed Brown                  {0.5,0,0}},
150*fe7e6d57SJed Brown       Gamma[3][3] = {{g,0,0},
151*fe7e6d57SJed Brown                      {-1.5773502691896257e+00,g,0},
152*fe7e6d57SJed Brown                      {-6.7075317547305480e-01,1.7075317547305482e-01,g}},
153*fe7e6d57SJed Brown       b[3] = {1.0566243270259355e-01,4.9038105676657971e-02,8.4529946162074843e-01},
154*fe7e6d57SJed Brown       b2[3] = {-1.7863279495408180e-01,1./3.,8.4529946162074843e-01};
155*fe7e6d57SJed Brown     ierr = TSRosWRegister(TSROSWRA3PW,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
156*fe7e6d57SJed Brown   }
157*fe7e6d57SJed Brown   {
158*fe7e6d57SJed Brown     const PetscReal g = 4.3586652150845900e-01;
159*fe7e6d57SJed Brown     const PetscReal
160*fe7e6d57SJed Brown       A[4][4] = {{0,0,0,0},
161*fe7e6d57SJed Brown                  {8.7173304301691801e-01,0,0,0},
162*fe7e6d57SJed Brown                  {8.4457060015369423e-01,-1.1299064236484185e-01,0,0},
163*fe7e6d57SJed Brown                  {0,0,1.,0}},
164*fe7e6d57SJed Brown       Gamma[4][4] = {{g,0,0,0},
165*fe7e6d57SJed Brown                      {-8.7173304301691801e-01,g,0,0},
166*fe7e6d57SJed Brown                      {-9.0338057013044082e-01,5.4180672388095326e-02,g,0},
167*fe7e6d57SJed Brown                      {2.4212380706095346e-01,-1.2232505839045147e+00,5.4526025533510214e-01,g}},
168*fe7e6d57SJed Brown       b[4] = {2.4212380706095346e-01,-1.2232505839045147e+00,1.5452602553351020e+00,4.3586652150845900e-01},
169*fe7e6d57SJed Brown       b2[4] = {3.7810903145819369e-01,-9.6042292212423178e-02,5.0000000000000000e-01,2.1793326075422950e-01};
170*fe7e6d57SJed Brown     ierr = TSRosWRegister(TSROSWRA34PW2,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
171e27a552bSJed Brown   }
172e27a552bSJed Brown   PetscFunctionReturn(0);
173e27a552bSJed Brown }
174e27a552bSJed Brown 
175e27a552bSJed Brown #undef __FUNCT__
176e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterDestroy"
177e27a552bSJed Brown /*@C
178e27a552bSJed Brown    TSRosWRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister().
179e27a552bSJed Brown 
180e27a552bSJed Brown    Not Collective
181e27a552bSJed Brown 
182e27a552bSJed Brown    Level: advanced
183e27a552bSJed Brown 
184e27a552bSJed Brown .keywords: TSRosW, register, destroy
185e27a552bSJed Brown .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegisterDynamic()
186e27a552bSJed Brown @*/
187e27a552bSJed Brown PetscErrorCode TSRosWRegisterDestroy(void)
188e27a552bSJed Brown {
189e27a552bSJed Brown   PetscErrorCode ierr;
19061692a83SJed Brown   RosWTableauLink link;
191e27a552bSJed Brown 
192e27a552bSJed Brown   PetscFunctionBegin;
19361692a83SJed Brown   while ((link = RosWTableauList)) {
19461692a83SJed Brown     RosWTableau t = &link->tab;
19561692a83SJed Brown     RosWTableauList = link->next;
19661692a83SJed Brown     ierr = PetscFree5(t->A,t->Gamma,t->b,t->ASum,t->GammaSum);CHKERRQ(ierr);
19761692a83SJed Brown     ierr = PetscFree3(t->At,t->bt,t->GammaInv);CHKERRQ(ierr);
198*fe7e6d57SJed Brown     ierr = PetscFree2(t->bembed,t->bembedt);CHKERRQ(ierr);
199e27a552bSJed Brown     ierr = PetscFree(t->name);CHKERRQ(ierr);
200e27a552bSJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
201e27a552bSJed Brown   }
202e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_FALSE;
203e27a552bSJed Brown   PetscFunctionReturn(0);
204e27a552bSJed Brown }
205e27a552bSJed Brown 
206e27a552bSJed Brown #undef __FUNCT__
207e27a552bSJed Brown #define __FUNCT__ "TSRosWInitializePackage"
208e27a552bSJed Brown /*@C
209e27a552bSJed Brown   TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called
210e27a552bSJed Brown   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RosW()
211e27a552bSJed Brown   when using static libraries.
212e27a552bSJed Brown 
213e27a552bSJed Brown   Input Parameter:
214e27a552bSJed Brown   path - The dynamic library path, or PETSC_NULL
215e27a552bSJed Brown 
216e27a552bSJed Brown   Level: developer
217e27a552bSJed Brown 
218e27a552bSJed Brown .keywords: TS, TSRosW, initialize, package
219e27a552bSJed Brown .seealso: PetscInitialize()
220e27a552bSJed Brown @*/
221e27a552bSJed Brown PetscErrorCode TSRosWInitializePackage(const char path[])
222e27a552bSJed Brown {
223e27a552bSJed Brown   PetscErrorCode ierr;
224e27a552bSJed Brown 
225e27a552bSJed Brown   PetscFunctionBegin;
226e27a552bSJed Brown   if (TSRosWPackageInitialized) PetscFunctionReturn(0);
227e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_TRUE;
228e27a552bSJed Brown   ierr = TSRosWRegisterAll();CHKERRQ(ierr);
229e27a552bSJed Brown   ierr = PetscRegisterFinalize(TSRosWFinalizePackage);CHKERRQ(ierr);
230e27a552bSJed Brown   PetscFunctionReturn(0);
231e27a552bSJed Brown }
232e27a552bSJed Brown 
233e27a552bSJed Brown #undef __FUNCT__
234e27a552bSJed Brown #define __FUNCT__ "TSRosWFinalizePackage"
235e27a552bSJed Brown /*@C
236e27a552bSJed Brown   TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is
237e27a552bSJed Brown   called from PetscFinalize().
238e27a552bSJed Brown 
239e27a552bSJed Brown   Level: developer
240e27a552bSJed Brown 
241e27a552bSJed Brown .keywords: Petsc, destroy, package
242e27a552bSJed Brown .seealso: PetscFinalize()
243e27a552bSJed Brown @*/
244e27a552bSJed Brown PetscErrorCode TSRosWFinalizePackage(void)
245e27a552bSJed Brown {
246e27a552bSJed Brown   PetscErrorCode ierr;
247e27a552bSJed Brown 
248e27a552bSJed Brown   PetscFunctionBegin;
249e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_FALSE;
250e27a552bSJed Brown   ierr = TSRosWRegisterDestroy();CHKERRQ(ierr);
251e27a552bSJed Brown   PetscFunctionReturn(0);
252e27a552bSJed Brown }
253e27a552bSJed Brown 
254e27a552bSJed Brown #undef __FUNCT__
255e27a552bSJed Brown #define __FUNCT__ "TSRosWRegister"
256e27a552bSJed Brown /*@C
25761692a83SJed Brown    TSRosWRegister - register a Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
258e27a552bSJed Brown 
259e27a552bSJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
260e27a552bSJed Brown 
261e27a552bSJed Brown    Input Parameters:
262e27a552bSJed Brown +  name - identifier for method
263e27a552bSJed Brown .  order - approximation order of method
264e27a552bSJed Brown .  s - number of stages, this is the dimension of the matrices below
26561692a83SJed Brown .  A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular
26661692a83SJed Brown .  Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal
267*fe7e6d57SJed Brown .  b - Step completion table (dimension s)
268*fe7e6d57SJed Brown -  bembed - Step completion table for a scheme of order one less (dimension s, PETSC_NULL if no embedded scheme is available)
269e27a552bSJed Brown 
270e27a552bSJed Brown    Notes:
27161692a83SJed Brown    Several Rosenbrock W methods are provided, this function is only needed to create new methods.
272e27a552bSJed Brown 
273e27a552bSJed Brown    Level: advanced
274e27a552bSJed Brown 
275e27a552bSJed Brown .keywords: TS, register
276e27a552bSJed Brown 
277e27a552bSJed Brown .seealso: TSRosW
278e27a552bSJed Brown @*/
279e27a552bSJed Brown PetscErrorCode TSRosWRegister(const TSRosWType name,PetscInt order,PetscInt s,
280*fe7e6d57SJed Brown                               const PetscReal A[],const PetscReal Gamma[],const PetscReal b[],const PetscReal bembed[])
281e27a552bSJed Brown {
282e27a552bSJed Brown   PetscErrorCode ierr;
28361692a83SJed Brown   RosWTableauLink link;
28461692a83SJed Brown   RosWTableau t;
28561692a83SJed Brown   PetscInt i,j,k;
28661692a83SJed Brown   PetscScalar *GammaInv;
287e27a552bSJed Brown 
288e27a552bSJed Brown   PetscFunctionBegin;
289*fe7e6d57SJed Brown   PetscValidCharPointer(name,1);
290*fe7e6d57SJed Brown   PetscValidPointer(A,4);
291*fe7e6d57SJed Brown   PetscValidPointer(Gamma,5);
292*fe7e6d57SJed Brown   PetscValidPointer(b,6);
293*fe7e6d57SJed Brown   if (bembed) PetscValidPointer(bembed,7);
294*fe7e6d57SJed Brown 
295e27a552bSJed Brown   ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
296e27a552bSJed Brown   ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
297e27a552bSJed Brown   t = &link->tab;
298e27a552bSJed Brown   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
299e27a552bSJed Brown   t->order = order;
300e27a552bSJed Brown   t->s = s;
30161692a83SJed 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);
30261692a83SJed Brown   ierr = PetscMalloc3(s*s,PetscReal,&t->At,s,PetscReal,&t->bt,s*s,PetscReal,&t->GammaInv);CHKERRQ(ierr);
303e27a552bSJed Brown   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
30461692a83SJed Brown   ierr = PetscMemcpy(t->Gamma,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr);
30561692a83SJed Brown   ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
306*fe7e6d57SJed Brown   if (bembed) {
307*fe7e6d57SJed Brown     ierr = PetscMalloc2(s,PetscReal,&t->bembed,s,PetscReal,&t->bembedt);CHKERRQ(ierr);
308*fe7e6d57SJed Brown     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
309*fe7e6d57SJed Brown   }
31061692a83SJed Brown   for (i=0; i<s; i++) {
31161692a83SJed Brown     t->ASum[i] = 0;
31261692a83SJed Brown     t->GammaSum[i] = 0;
31361692a83SJed Brown     for (j=0; j<s; j++) {
31461692a83SJed Brown       t->ASum[i] += A[i*s+j];
315*fe7e6d57SJed Brown       t->GammaSum[i] += Gamma[i*s+j];
31661692a83SJed Brown     }
31761692a83SJed Brown   }
31861692a83SJed Brown   ierr = PetscMalloc(s*s*sizeof(PetscScalar),&GammaInv);CHKERRQ(ierr); /* Need to use Scalar for inverse, then convert back to Real */
31961692a83SJed Brown   for (i=0; i<s*s; i++) GammaInv[i] = Gamma[i];
32061692a83SJed Brown   switch (s) {
32161692a83SJed Brown   case 1: GammaInv[0] = 1./GammaInv[0]; break;
32261692a83SJed Brown   case 2: ierr = Kernel_A_gets_inverse_A_2(GammaInv,0);CHKERRQ(ierr); break;
32361692a83SJed Brown   case 3: ierr = Kernel_A_gets_inverse_A_3(GammaInv,0);CHKERRQ(ierr); break;
32461692a83SJed Brown   case 4: ierr = Kernel_A_gets_inverse_A_4(GammaInv,0);CHKERRQ(ierr); break;
32561692a83SJed Brown   case 5: {
32661692a83SJed Brown     PetscInt ipvt5[5];
32761692a83SJed Brown     MatScalar work5[5*5];
32861692a83SJed Brown     ierr = Kernel_A_gets_inverse_A_5(GammaInv,ipvt5,work5,0);CHKERRQ(ierr); break;
32961692a83SJed Brown   }
33061692a83SJed Brown   case 6: ierr = Kernel_A_gets_inverse_A_6(GammaInv,0);CHKERRQ(ierr); break;
33161692a83SJed Brown   case 7: ierr = Kernel_A_gets_inverse_A_7(GammaInv,0);CHKERRQ(ierr); break;
33261692a83SJed Brown   default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for %D stages",s);
33361692a83SJed Brown   }
33461692a83SJed Brown   for (i=0; i<s*s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]);
33561692a83SJed Brown   ierr = PetscFree(GammaInv);CHKERRQ(ierr);
33661692a83SJed Brown   for (i=0; i<s; i++) {
33761692a83SJed Brown     for (j=0; j<s; j++) {
33861692a83SJed Brown       t->At[i*s+j] = 0;
33961692a83SJed Brown       for (k=0; k<s; k++) {
34061692a83SJed Brown         t->At[i*s+j] += t->A[i*s+k] * t->GammaInv[k*s+j];
34161692a83SJed Brown       }
34261692a83SJed Brown     }
34361692a83SJed Brown     t->bt[i] = 0;
34461692a83SJed Brown     for (j=0; j<s; j++) {
34561692a83SJed Brown       t->bt[i] += t->b[j] * t->GammaInv[j*s+i];
34661692a83SJed Brown     }
347*fe7e6d57SJed Brown     if (bembed) {
348*fe7e6d57SJed Brown       t->bembedt[i] = 0;
349*fe7e6d57SJed Brown       for (j=0; j<s; j++) {
350*fe7e6d57SJed Brown         t->bembedt[i] += t->bembed[j] * t->GammaInv[j*s+i];
351*fe7e6d57SJed Brown       }
352*fe7e6d57SJed Brown     }
35361692a83SJed Brown   }
35461692a83SJed Brown   link->next = RosWTableauList;
35561692a83SJed Brown   RosWTableauList = link;
356e27a552bSJed Brown   PetscFunctionReturn(0);
357e27a552bSJed Brown }
358e27a552bSJed Brown 
359e27a552bSJed Brown #undef __FUNCT__
360e27a552bSJed Brown #define __FUNCT__ "TSStep_RosW"
361e27a552bSJed Brown static PetscErrorCode TSStep_RosW(TS ts)
362e27a552bSJed Brown {
36361692a83SJed Brown   TS_RosW         *ros = (TS_RosW*)ts->data;
36461692a83SJed Brown   RosWTableau     tab  = ros->tableau;
365e27a552bSJed Brown   const PetscInt  s    = tab->s;
36661692a83SJed Brown   const PetscReal *At  = tab->At,*Gamma = tab->Gamma,*bt = tab->bt,*ASum = tab->ASum,*GammaInv = tab->GammaInv;
36761692a83SJed Brown   PetscScalar     *w   = ros->work;
36861692a83SJed Brown   Vec             *Y   = ros->Y,Zdot = ros->Zdot,Zstage = ros->Zstage;
369e27a552bSJed Brown   SNES            snes;
370e27a552bSJed Brown   PetscInt        i,j,its,lits;
371cdbf8f93SLisandro Dalcin   PetscReal       next_time_step;
372e27a552bSJed Brown   PetscReal       h,t;
373e27a552bSJed Brown   PetscErrorCode  ierr;
374e27a552bSJed Brown 
375e27a552bSJed Brown   PetscFunctionBegin;
376e27a552bSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
377cdbf8f93SLisandro Dalcin   next_time_step = ts->time_step;
378cdbf8f93SLisandro Dalcin   h = ts->time_step;
379e27a552bSJed Brown   t = ts->ptime;
380e27a552bSJed Brown 
381e27a552bSJed Brown   for (i=0; i<s; i++) {
38261692a83SJed Brown     ros->stage_time = t + h*ASum[i];
38361692a83SJed Brown     ros->shift = 1./(h*Gamma[i*s+i]);
38461692a83SJed Brown 
38561692a83SJed Brown     ierr = VecCopy(ts->vec_sol,Zstage);CHKERRQ(ierr);
38661692a83SJed Brown     ierr = VecMAXPY(Zstage,i,&At[i*s+0],Y);CHKERRQ(ierr);
38761692a83SJed Brown 
38861692a83SJed Brown     for (j=0; j<i; j++) w[j] = 1./h * GammaInv[i*s+j];
38961692a83SJed Brown     ierr = VecZeroEntries(Zdot);CHKERRQ(ierr);
39061692a83SJed Brown     ierr = VecMAXPY(Zdot,i,w,Y);CHKERRQ(ierr);
39161692a83SJed Brown 
392e27a552bSJed Brown     /* Initial guess taken from last stage */
39361692a83SJed Brown     ierr = VecZeroEntries(Y[i]);CHKERRQ(ierr);
39461692a83SJed Brown 
39561692a83SJed Brown     if (!ros->recompute_jacobian && !i) {
39661692a83SJed Brown       ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */
39761692a83SJed Brown     }
39861692a83SJed Brown 
39961692a83SJed Brown     ierr = SNESSolve(snes,PETSC_NULL,Y[i]);CHKERRQ(ierr);
400e27a552bSJed Brown     ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
401e27a552bSJed Brown     ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
402e27a552bSJed Brown     ts->nonlinear_its += its; ts->linear_its += lits;
403e27a552bSJed Brown   }
40461692a83SJed Brown   ierr = VecMAXPY(ts->vec_sol,s,bt,Y);CHKERRQ(ierr);
405e27a552bSJed Brown 
406e27a552bSJed Brown   ts->ptime += ts->time_step;
407cdbf8f93SLisandro Dalcin   ts->time_step = next_time_step;
408e27a552bSJed Brown   ts->steps++;
409e27a552bSJed Brown   PetscFunctionReturn(0);
410e27a552bSJed Brown }
411e27a552bSJed Brown 
412e27a552bSJed Brown #undef __FUNCT__
413e27a552bSJed Brown #define __FUNCT__ "TSInterpolate_RosW"
414e27a552bSJed Brown static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec X)
415e27a552bSJed Brown {
41661692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
417e27a552bSJed Brown 
418e27a552bSJed Brown   PetscFunctionBegin;
41961692a83SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name);
420e27a552bSJed Brown   PetscFunctionReturn(0);
421e27a552bSJed Brown }
422e27a552bSJed Brown 
423e27a552bSJed Brown /*------------------------------------------------------------*/
424e27a552bSJed Brown #undef __FUNCT__
425e27a552bSJed Brown #define __FUNCT__ "TSReset_RosW"
426e27a552bSJed Brown static PetscErrorCode TSReset_RosW(TS ts)
427e27a552bSJed Brown {
42861692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
429e27a552bSJed Brown   PetscInt       s;
430e27a552bSJed Brown   PetscErrorCode ierr;
431e27a552bSJed Brown 
432e27a552bSJed Brown   PetscFunctionBegin;
43361692a83SJed Brown   if (!ros->tableau) PetscFunctionReturn(0);
43461692a83SJed Brown    s = ros->tableau->s;
43561692a83SJed Brown   ierr = VecDestroyVecs(s,&ros->Y);CHKERRQ(ierr);
43661692a83SJed Brown   ierr = VecDestroy(&ros->Ydot);CHKERRQ(ierr);
43761692a83SJed Brown   ierr = VecDestroy(&ros->Ystage);CHKERRQ(ierr);
43861692a83SJed Brown   ierr = VecDestroy(&ros->Zdot);CHKERRQ(ierr);
43961692a83SJed Brown   ierr = VecDestroy(&ros->Zstage);CHKERRQ(ierr);
44061692a83SJed Brown   ierr = PetscFree(ros->work);CHKERRQ(ierr);
441e27a552bSJed Brown   PetscFunctionReturn(0);
442e27a552bSJed Brown }
443e27a552bSJed Brown 
444e27a552bSJed Brown #undef __FUNCT__
445e27a552bSJed Brown #define __FUNCT__ "TSDestroy_RosW"
446e27a552bSJed Brown static PetscErrorCode TSDestroy_RosW(TS ts)
447e27a552bSJed Brown {
448e27a552bSJed Brown   PetscErrorCode  ierr;
449e27a552bSJed Brown 
450e27a552bSJed Brown   PetscFunctionBegin;
451e27a552bSJed Brown   ierr = TSReset_RosW(ts);CHKERRQ(ierr);
452e27a552bSJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
453e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","",PETSC_NULL);CHKERRQ(ierr);
454e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","",PETSC_NULL);CHKERRQ(ierr);
45561692a83SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","",PETSC_NULL);CHKERRQ(ierr);
456e27a552bSJed Brown   PetscFunctionReturn(0);
457e27a552bSJed Brown }
458e27a552bSJed Brown 
459e27a552bSJed Brown /*
460e27a552bSJed Brown   This defines the nonlinear equation that is to be solved with SNES
461e27a552bSJed Brown   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
462e27a552bSJed Brown */
463e27a552bSJed Brown #undef __FUNCT__
464e27a552bSJed Brown #define __FUNCT__ "SNESTSFormFunction_RosW"
465e27a552bSJed Brown static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec X,Vec F,TS ts)
466e27a552bSJed Brown {
46761692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
468e27a552bSJed Brown   PetscErrorCode ierr;
469e27a552bSJed Brown 
470e27a552bSJed Brown   PetscFunctionBegin;
47161692a83SJed Brown   ierr = VecWAXPY(ros->Ydot,ros->shift,X,ros->Zdot);CHKERRQ(ierr); /* Ydot = shift*X + Zdot */
47261692a83SJed Brown   ierr = VecWAXPY(ros->Ystage,1.0,X,ros->Zstage);CHKERRQ(ierr);    /* Ystage = X + Zstage */
47361692a83SJed Brown   ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,ros->Ydot,F,PETSC_FALSE);CHKERRQ(ierr);
474e27a552bSJed Brown   PetscFunctionReturn(0);
475e27a552bSJed Brown }
476e27a552bSJed Brown 
477e27a552bSJed Brown #undef __FUNCT__
478e27a552bSJed Brown #define __FUNCT__ "SNESTSFormJacobian_RosW"
479e27a552bSJed Brown static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
480e27a552bSJed Brown {
48161692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
482e27a552bSJed Brown   PetscErrorCode ierr;
483e27a552bSJed Brown 
484e27a552bSJed Brown   PetscFunctionBegin;
48561692a83SJed Brown   /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
48661692a83SJed Brown   ierr = TSComputeIJacobian(ts,ros->stage_time,ros->Ystage,ros->Ydot,ros->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr);
487e27a552bSJed Brown   PetscFunctionReturn(0);
488e27a552bSJed Brown }
489e27a552bSJed Brown 
490e27a552bSJed Brown #undef __FUNCT__
491e27a552bSJed Brown #define __FUNCT__ "TSSetUp_RosW"
492e27a552bSJed Brown static PetscErrorCode TSSetUp_RosW(TS ts)
493e27a552bSJed Brown {
49461692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
49561692a83SJed Brown   RosWTableau    tab  = ros->tableau;
496e27a552bSJed Brown   PetscInt       s    = tab->s;
497e27a552bSJed Brown   PetscErrorCode ierr;
498e27a552bSJed Brown 
499e27a552bSJed Brown   PetscFunctionBegin;
50061692a83SJed Brown   if (!ros->tableau) {
501e27a552bSJed Brown     ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);
502e27a552bSJed Brown   }
50361692a83SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ros->Y);CHKERRQ(ierr);
50461692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Ydot);CHKERRQ(ierr);
50561692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Ystage);CHKERRQ(ierr);
50661692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Zdot);CHKERRQ(ierr);
50761692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Zstage);CHKERRQ(ierr);
50861692a83SJed Brown   ierr = PetscMalloc(s*sizeof(ros->work[0]),&ros->work);CHKERRQ(ierr);
509e27a552bSJed Brown   PetscFunctionReturn(0);
510e27a552bSJed Brown }
511e27a552bSJed Brown /*------------------------------------------------------------*/
512e27a552bSJed Brown 
513e27a552bSJed Brown #undef __FUNCT__
514e27a552bSJed Brown #define __FUNCT__ "TSSetFromOptions_RosW"
515e27a552bSJed Brown static PetscErrorCode TSSetFromOptions_RosW(TS ts)
516e27a552bSJed Brown {
51761692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
518e27a552bSJed Brown   PetscErrorCode ierr;
51961692a83SJed Brown   char           rostype[256];
520e27a552bSJed Brown 
521e27a552bSJed Brown   PetscFunctionBegin;
522e27a552bSJed Brown   ierr = PetscOptionsHead("RosW ODE solver options");CHKERRQ(ierr);
523e27a552bSJed Brown   {
52461692a83SJed Brown     RosWTableauLink link;
525e27a552bSJed Brown     PetscInt count,choice;
526e27a552bSJed Brown     PetscBool flg;
527e27a552bSJed Brown     const char **namelist;
52861692a83SJed Brown     SNES snes;
52961692a83SJed Brown 
53061692a83SJed Brown     ierr = PetscStrncpy(rostype,TSRosWDefault,sizeof rostype);CHKERRQ(ierr);
53161692a83SJed Brown     for (link=RosWTableauList,count=0; link; link=link->next,count++) ;
532e27a552bSJed Brown     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
53361692a83SJed Brown     for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
53461692a83SJed Brown     ierr = PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,rostype,&choice,&flg);CHKERRQ(ierr);
53561692a83SJed Brown     ierr = TSRosWSetType(ts,flg ? namelist[choice] : rostype);CHKERRQ(ierr);
536e27a552bSJed Brown     ierr = PetscFree(namelist);CHKERRQ(ierr);
53761692a83SJed Brown 
53861692a83SJed Brown     ierr = PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,PETSC_NULL);CHKERRQ(ierr);
53961692a83SJed Brown 
54061692a83SJed Brown     /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
54161692a83SJed Brown     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
54261692a83SJed Brown     if (!((PetscObject)snes)->type_name) {
54361692a83SJed Brown       ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
54461692a83SJed Brown     }
54561692a83SJed Brown     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
546e27a552bSJed Brown   }
547e27a552bSJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
548e27a552bSJed Brown   PetscFunctionReturn(0);
549e27a552bSJed Brown }
550e27a552bSJed Brown 
551e27a552bSJed Brown #undef __FUNCT__
552e27a552bSJed Brown #define __FUNCT__ "PetscFormatRealArray"
553e27a552bSJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
554e27a552bSJed Brown {
555e27a552bSJed Brown   PetscErrorCode ierr;
556e408995aSJed Brown   PetscInt i;
557e408995aSJed Brown   size_t left,count;
558e27a552bSJed Brown   char *p;
559e27a552bSJed Brown 
560e27a552bSJed Brown   PetscFunctionBegin;
561e408995aSJed Brown   for (i=0,p=buf,left=len; i<n; i++) {
562e408995aSJed Brown     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
563e27a552bSJed Brown     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
564e27a552bSJed Brown     left -= count;
565e27a552bSJed Brown     p += count;
566e27a552bSJed Brown     *p++ = ' ';
567e27a552bSJed Brown   }
568e27a552bSJed Brown   p[i ? 0 : -1] = 0;
569e27a552bSJed Brown   PetscFunctionReturn(0);
570e27a552bSJed Brown }
571e27a552bSJed Brown 
572e27a552bSJed Brown #undef __FUNCT__
573e27a552bSJed Brown #define __FUNCT__ "TSView_RosW"
574e27a552bSJed Brown static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer)
575e27a552bSJed Brown {
57661692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
57761692a83SJed Brown   RosWTableau    tab  = ros->tableau;
578e27a552bSJed Brown   PetscBool      iascii;
579e27a552bSJed Brown   PetscErrorCode ierr;
580e27a552bSJed Brown 
581e27a552bSJed Brown   PetscFunctionBegin;
582e27a552bSJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
583e27a552bSJed Brown   if (iascii) {
58461692a83SJed Brown     const TSRosWType rostype;
585e408995aSJed Brown     PetscInt i;
586e408995aSJed Brown     PetscReal abscissa[512];
587e27a552bSJed Brown     char buf[512];
58861692a83SJed Brown     ierr = TSRosWGetType(ts,&rostype);CHKERRQ(ierr);
58961692a83SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Rosenbrock-W %s\n",rostype);CHKERRQ(ierr);
590e408995aSJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ASum);CHKERRQ(ierr);
59161692a83SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A       = %s\n",buf);CHKERRQ(ierr);
592e408995aSJed Brown     for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i];
593e408995aSJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,abscissa);CHKERRQ(ierr);
594e408995aSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A+Gamma = %s\n",buf);CHKERRQ(ierr);
595e27a552bSJed Brown   }
596e27a552bSJed Brown   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
597e27a552bSJed Brown   PetscFunctionReturn(0);
598e27a552bSJed Brown }
599e27a552bSJed Brown 
600e27a552bSJed Brown #undef __FUNCT__
601e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType"
602e27a552bSJed Brown /*@C
60361692a83SJed Brown   TSRosWSetType - Set the type of Rosenbrock-W scheme
604e27a552bSJed Brown 
605e27a552bSJed Brown   Logically collective
606e27a552bSJed Brown 
607e27a552bSJed Brown   Input Parameter:
608e27a552bSJed Brown +  ts - timestepping context
60961692a83SJed Brown -  rostype - type of Rosenbrock-W scheme
610e27a552bSJed Brown 
611e27a552bSJed Brown   Level: intermediate
612e27a552bSJed Brown 
613e27a552bSJed Brown .seealso: TSRosWGetType()
614e27a552bSJed Brown @*/
61561692a83SJed Brown PetscErrorCode TSRosWSetType(TS ts,const TSRosWType rostype)
616e27a552bSJed Brown {
617e27a552bSJed Brown   PetscErrorCode ierr;
618e27a552bSJed Brown 
619e27a552bSJed Brown   PetscFunctionBegin;
620e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
62161692a83SJed Brown   ierr = PetscTryMethod(ts,"TSRosWSetType_C",(TS,const TSRosWType),(ts,rostype));CHKERRQ(ierr);
622e27a552bSJed Brown   PetscFunctionReturn(0);
623e27a552bSJed Brown }
624e27a552bSJed Brown 
625e27a552bSJed Brown #undef __FUNCT__
626e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType"
627e27a552bSJed Brown /*@C
62861692a83SJed Brown   TSRosWGetType - Get the type of Rosenbrock-W scheme
629e27a552bSJed Brown 
630e27a552bSJed Brown   Logically collective
631e27a552bSJed Brown 
632e27a552bSJed Brown   Input Parameter:
633e27a552bSJed Brown .  ts - timestepping context
634e27a552bSJed Brown 
635e27a552bSJed Brown   Output Parameter:
63661692a83SJed Brown .  rostype - type of Rosenbrock-W scheme
637e27a552bSJed Brown 
638e27a552bSJed Brown   Level: intermediate
639e27a552bSJed Brown 
640e27a552bSJed Brown .seealso: TSRosWGetType()
641e27a552bSJed Brown @*/
64261692a83SJed Brown PetscErrorCode TSRosWGetType(TS ts,const TSRosWType *rostype)
643e27a552bSJed Brown {
644e27a552bSJed Brown   PetscErrorCode ierr;
645e27a552bSJed Brown 
646e27a552bSJed Brown   PetscFunctionBegin;
647e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
64861692a83SJed Brown   ierr = PetscUseMethod(ts,"TSRosWGetType_C",(TS,const TSRosWType*),(ts,rostype));CHKERRQ(ierr);
649e27a552bSJed Brown   PetscFunctionReturn(0);
650e27a552bSJed Brown }
651e27a552bSJed Brown 
652e27a552bSJed Brown #undef __FUNCT__
65361692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian"
654e27a552bSJed Brown /*@C
65561692a83SJed Brown   TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step.
656e27a552bSJed Brown 
657e27a552bSJed Brown   Logically collective
658e27a552bSJed Brown 
659e27a552bSJed Brown   Input Parameter:
660e27a552bSJed Brown +  ts - timestepping context
66161692a83SJed Brown -  flg - PETSC_TRUE to recompute the Jacobian at each stage
662e27a552bSJed Brown 
663e27a552bSJed Brown   Level: intermediate
664e27a552bSJed Brown 
665e27a552bSJed Brown .seealso: TSRosWGetType()
666e27a552bSJed Brown @*/
66761692a83SJed Brown PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg)
668e27a552bSJed Brown {
669e27a552bSJed Brown   PetscErrorCode ierr;
670e27a552bSJed Brown 
671e27a552bSJed Brown   PetscFunctionBegin;
672e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
67361692a83SJed Brown   ierr = PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
674e27a552bSJed Brown   PetscFunctionReturn(0);
675e27a552bSJed Brown }
676e27a552bSJed Brown 
677e27a552bSJed Brown EXTERN_C_BEGIN
678e27a552bSJed Brown #undef __FUNCT__
679e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType_RosW"
68061692a83SJed Brown PetscErrorCode  TSRosWGetType_RosW(TS ts,const TSRosWType *rostype)
681e27a552bSJed Brown {
68261692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
683e27a552bSJed Brown   PetscErrorCode ierr;
684e27a552bSJed Brown 
685e27a552bSJed Brown   PetscFunctionBegin;
68661692a83SJed Brown   if (!ros->tableau) {ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);}
68761692a83SJed Brown   *rostype = ros->tableau->name;
688e27a552bSJed Brown   PetscFunctionReturn(0);
689e27a552bSJed Brown }
690e27a552bSJed Brown #undef __FUNCT__
691e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType_RosW"
69261692a83SJed Brown PetscErrorCode  TSRosWSetType_RosW(TS ts,const TSRosWType rostype)
693e27a552bSJed Brown {
69461692a83SJed Brown   TS_RosW         *ros = (TS_RosW*)ts->data;
695e27a552bSJed Brown   PetscErrorCode  ierr;
696e27a552bSJed Brown   PetscBool       match;
69761692a83SJed Brown   RosWTableauLink link;
698e27a552bSJed Brown 
699e27a552bSJed Brown   PetscFunctionBegin;
70061692a83SJed Brown   if (ros->tableau) {
70161692a83SJed Brown     ierr = PetscStrcmp(ros->tableau->name,rostype,&match);CHKERRQ(ierr);
702e27a552bSJed Brown     if (match) PetscFunctionReturn(0);
703e27a552bSJed Brown   }
70461692a83SJed Brown   for (link = RosWTableauList; link; link=link->next) {
70561692a83SJed Brown     ierr = PetscStrcmp(link->tab.name,rostype,&match);CHKERRQ(ierr);
706e27a552bSJed Brown     if (match) {
707e27a552bSJed Brown       ierr = TSReset_RosW(ts);CHKERRQ(ierr);
70861692a83SJed Brown       ros->tableau = &link->tab;
709e27a552bSJed Brown       PetscFunctionReturn(0);
710e27a552bSJed Brown     }
711e27a552bSJed Brown   }
71261692a83SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype);
713e27a552bSJed Brown   PetscFunctionReturn(0);
714e27a552bSJed Brown }
71561692a83SJed Brown 
716e27a552bSJed Brown #undef __FUNCT__
71761692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian_RosW"
71861692a83SJed Brown PetscErrorCode  TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg)
719e27a552bSJed Brown {
72061692a83SJed Brown   TS_RosW *ros = (TS_RosW*)ts->data;
721e27a552bSJed Brown 
722e27a552bSJed Brown   PetscFunctionBegin;
72361692a83SJed Brown   ros->recompute_jacobian = flg;
724e27a552bSJed Brown   PetscFunctionReturn(0);
725e27a552bSJed Brown }
726e27a552bSJed Brown EXTERN_C_END
727e27a552bSJed Brown 
728e27a552bSJed Brown /* ------------------------------------------------------------ */
729e27a552bSJed Brown /*MC
730e27a552bSJed Brown       TSRosW - ODE solver using Rosenbrock-W schemes
731e27a552bSJed Brown 
732e27a552bSJed Brown   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
733e27a552bSJed Brown   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
734e27a552bSJed Brown   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
735e27a552bSJed Brown 
736e27a552bSJed Brown   Notes:
73761692a83SJed Brown   This method currently only works with autonomous ODE and DAE.
73861692a83SJed Brown 
73961692a83SJed Brown   Developer notes:
74061692a83SJed Brown   Rosenbrock-W methods are typically specified for autonomous ODE
74161692a83SJed Brown 
74261692a83SJed Brown $  xdot = f(x)
74361692a83SJed Brown 
74461692a83SJed Brown   by the stage equations
74561692a83SJed Brown 
74661692a83SJed Brown $  k_i = h f(x_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j
74761692a83SJed Brown 
74861692a83SJed Brown   and step completion formula
74961692a83SJed Brown 
75061692a83SJed Brown $  x_1 = x_0 + sum_j b_j k_j
75161692a83SJed Brown 
75261692a83SJed Brown   with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(x)
75361692a83SJed Brown   and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner,
75461692a83SJed Brown   we define new variables for the stage equations
75561692a83SJed Brown 
75661692a83SJed Brown $  y_i = gamma_ij k_j
75761692a83SJed Brown 
75861692a83SJed Brown   The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define
75961692a83SJed Brown 
76061692a83SJed Brown $  A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-i}
76161692a83SJed Brown 
76261692a83SJed Brown   to rewrite the method as
76361692a83SJed Brown 
76461692a83SJed 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
76561692a83SJed Brown $  x_1 = x_0 + sum_j bt_j y_j
76661692a83SJed Brown 
76761692a83SJed Brown    where we have introduced the mass matrix M. Continue by defining
76861692a83SJed Brown 
76961692a83SJed Brown $  ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j
77061692a83SJed Brown 
77161692a83SJed Brown    or, more compactly in tensor notation
77261692a83SJed Brown 
77361692a83SJed Brown $  Ydot = 1/h (Gamma^{-1} \otimes I) Y .
77461692a83SJed Brown 
77561692a83SJed Brown    Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current
77661692a83SJed Brown    stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the
77761692a83SJed Brown    equation
77861692a83SJed Brown 
77961692a83SJed Brown $  g(x_0 + sum_j a_ij y_j + y_i, ydot_i) = 0
78061692a83SJed Brown 
78161692a83SJed Brown    with initial guess y_i = 0.
782e27a552bSJed Brown 
783e27a552bSJed Brown   Level: beginner
784e27a552bSJed Brown 
785e27a552bSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSRosWRegister()
786e27a552bSJed Brown 
787e27a552bSJed Brown M*/
788e27a552bSJed Brown EXTERN_C_BEGIN
789e27a552bSJed Brown #undef __FUNCT__
790e27a552bSJed Brown #define __FUNCT__ "TSCreate_RosW"
791e27a552bSJed Brown PetscErrorCode  TSCreate_RosW(TS ts)
792e27a552bSJed Brown {
79361692a83SJed Brown   TS_RosW        *ros;
794e27a552bSJed Brown   PetscErrorCode ierr;
795e27a552bSJed Brown 
796e27a552bSJed Brown   PetscFunctionBegin;
797e27a552bSJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
798e27a552bSJed Brown   ierr = TSRosWInitializePackage(PETSC_NULL);CHKERRQ(ierr);
799e27a552bSJed Brown #endif
800e27a552bSJed Brown 
801e27a552bSJed Brown   ts->ops->reset          = TSReset_RosW;
802e27a552bSJed Brown   ts->ops->destroy        = TSDestroy_RosW;
803e27a552bSJed Brown   ts->ops->view           = TSView_RosW;
804e27a552bSJed Brown   ts->ops->setup          = TSSetUp_RosW;
805e27a552bSJed Brown   ts->ops->step           = TSStep_RosW;
806e27a552bSJed Brown   ts->ops->interpolate    = TSInterpolate_RosW;
807e27a552bSJed Brown   ts->ops->setfromoptions = TSSetFromOptions_RosW;
808e27a552bSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_RosW;
809e27a552bSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_RosW;
810e27a552bSJed Brown 
81161692a83SJed Brown   ierr = PetscNewLog(ts,TS_RosW,&ros);CHKERRQ(ierr);
81261692a83SJed Brown   ts->data = (void*)ros;
813e27a552bSJed Brown 
814e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","TSRosWGetType_RosW",TSRosWGetType_RosW);CHKERRQ(ierr);
815e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","TSRosWSetType_RosW",TSRosWSetType_RosW);CHKERRQ(ierr);
81661692a83SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","TSRosWSetRecomputeJacobian_RosW",TSRosWSetRecomputeJacobian_RosW);CHKERRQ(ierr);
817e27a552bSJed Brown   PetscFunctionReturn(0);
818e27a552bSJed Brown }
819e27a552bSJed Brown EXTERN_C_END
820