xref: /petsc/src/ts/impls/rosw/rosw.c (revision 8d59e96071807dd8d22c4ae77f4b4dc641017dcf)
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 */
37*8d59e960SJed 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 
111e27a552bSJed Brown #undef __FUNCT__
112e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterAll"
113e27a552bSJed Brown /*@C
114e27a552bSJed Brown   TSRosWRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSRosW
115e27a552bSJed Brown 
116e27a552bSJed Brown   Not Collective, but should be called by all processes which will need the schemes to be registered
117e27a552bSJed Brown 
118e27a552bSJed Brown   Level: advanced
119e27a552bSJed Brown 
120e27a552bSJed Brown .keywords: TS, TSRosW, register, all
121e27a552bSJed Brown 
122e27a552bSJed Brown .seealso:  TSRosWRegisterDestroy()
123e27a552bSJed Brown @*/
124e27a552bSJed Brown PetscErrorCode TSRosWRegisterAll(void)
125e27a552bSJed Brown {
126e27a552bSJed Brown   PetscErrorCode ierr;
127e27a552bSJed Brown 
128e27a552bSJed Brown   PetscFunctionBegin;
129e27a552bSJed Brown   if (TSRosWRegisterAllCalled) PetscFunctionReturn(0);
130e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_TRUE;
131e27a552bSJed Brown 
132e27a552bSJed Brown   {
13361692a83SJed Brown     const PetscReal g = 1. + 1./PetscSqrtReal(2.0);
134e27a552bSJed Brown     const PetscReal
13561692a83SJed Brown       A[2][2] = {{0,0}, {1.,0}},
13661692a83SJed Brown       Gamma[2][2] = {{g,0}, {-2.*g,g}},
1371c3436cfSJed Brown       b[2] = {0.5,0.5},
1381c3436cfSJed Brown       b1[2] = {1.0,0.0};
1391c3436cfSJed Brown     ierr = TSRosWRegister(TSROSW2P,2,2,&A[0][0],&Gamma[0][0],b,b1);CHKERRQ(ierr);
140e27a552bSJed Brown   }
141e27a552bSJed Brown   {
14261692a83SJed Brown     const PetscReal g = 1. - 1./PetscSqrtReal(2.0);
143e27a552bSJed Brown     const PetscReal
14461692a83SJed Brown       A[2][2] = {{0,0}, {1.,0}},
14561692a83SJed Brown       Gamma[2][2] = {{g,0}, {-2.*g,g}},
1461c3436cfSJed Brown       b[2] = {0.5,0.5},
1471c3436cfSJed Brown       b1[2] = {1.0,0.0};
1481c3436cfSJed Brown     ierr = TSRosWRegister(TSROSW2M,2,2,&A[0][0],&Gamma[0][0],b,b1);CHKERRQ(ierr);
149fe7e6d57SJed Brown   }
150fe7e6d57SJed Brown   {
151fe7e6d57SJed Brown     const PetscReal g = 7.8867513459481287e-01;
152fe7e6d57SJed Brown     const PetscReal
153fe7e6d57SJed Brown       A[3][3] = {{0,0,0},
154fe7e6d57SJed Brown                  {1.5773502691896257e+00,0,0},
155fe7e6d57SJed Brown                  {0.5,0,0}},
156fe7e6d57SJed Brown       Gamma[3][3] = {{g,0,0},
157fe7e6d57SJed Brown                      {-1.5773502691896257e+00,g,0},
158fe7e6d57SJed Brown                      {-6.7075317547305480e-01,1.7075317547305482e-01,g}},
159fe7e6d57SJed Brown       b[3] = {1.0566243270259355e-01,4.9038105676657971e-02,8.4529946162074843e-01},
160fe7e6d57SJed Brown       b2[3] = {-1.7863279495408180e-01,1./3.,8.4529946162074843e-01};
161fe7e6d57SJed Brown     ierr = TSRosWRegister(TSROSWRA3PW,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
162fe7e6d57SJed Brown   }
163fe7e6d57SJed Brown   {
164fe7e6d57SJed Brown     const PetscReal g = 4.3586652150845900e-01;
165fe7e6d57SJed Brown     const PetscReal
166fe7e6d57SJed Brown       A[4][4] = {{0,0,0,0},
167fe7e6d57SJed Brown                  {8.7173304301691801e-01,0,0,0},
168fe7e6d57SJed Brown                  {8.4457060015369423e-01,-1.1299064236484185e-01,0,0},
169fe7e6d57SJed Brown                  {0,0,1.,0}},
170fe7e6d57SJed Brown       Gamma[4][4] = {{g,0,0,0},
171fe7e6d57SJed Brown                      {-8.7173304301691801e-01,g,0,0},
172fe7e6d57SJed Brown                      {-9.0338057013044082e-01,5.4180672388095326e-02,g,0},
173fe7e6d57SJed Brown                      {2.4212380706095346e-01,-1.2232505839045147e+00,5.4526025533510214e-01,g}},
174fe7e6d57SJed Brown       b[4] = {2.4212380706095346e-01,-1.2232505839045147e+00,1.5452602553351020e+00,4.3586652150845900e-01},
175fe7e6d57SJed Brown       b2[4] = {3.7810903145819369e-01,-9.6042292212423178e-02,5.0000000000000000e-01,2.1793326075422950e-01};
176fe7e6d57SJed Brown     ierr = TSRosWRegister(TSROSWRA34PW2,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
177e27a552bSJed Brown   }
178e27a552bSJed Brown   PetscFunctionReturn(0);
179e27a552bSJed Brown }
180e27a552bSJed Brown 
181e27a552bSJed Brown #undef __FUNCT__
182e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterDestroy"
183e27a552bSJed Brown /*@C
184e27a552bSJed Brown    TSRosWRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister().
185e27a552bSJed Brown 
186e27a552bSJed Brown    Not Collective
187e27a552bSJed Brown 
188e27a552bSJed Brown    Level: advanced
189e27a552bSJed Brown 
190e27a552bSJed Brown .keywords: TSRosW, register, destroy
191e27a552bSJed Brown .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegisterDynamic()
192e27a552bSJed Brown @*/
193e27a552bSJed Brown PetscErrorCode TSRosWRegisterDestroy(void)
194e27a552bSJed Brown {
195e27a552bSJed Brown   PetscErrorCode ierr;
19661692a83SJed Brown   RosWTableauLink link;
197e27a552bSJed Brown 
198e27a552bSJed Brown   PetscFunctionBegin;
19961692a83SJed Brown   while ((link = RosWTableauList)) {
20061692a83SJed Brown     RosWTableau t = &link->tab;
20161692a83SJed Brown     RosWTableauList = link->next;
20261692a83SJed Brown     ierr = PetscFree5(t->A,t->Gamma,t->b,t->ASum,t->GammaSum);CHKERRQ(ierr);
203c17803e7SJed Brown     ierr = PetscFree4(t->At,t->bt,t->GammaInv,t->GammaZeroDiag);CHKERRQ(ierr);
204fe7e6d57SJed Brown     ierr = PetscFree2(t->bembed,t->bembedt);CHKERRQ(ierr);
205e27a552bSJed Brown     ierr = PetscFree(t->name);CHKERRQ(ierr);
206e27a552bSJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
207e27a552bSJed Brown   }
208e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_FALSE;
209e27a552bSJed Brown   PetscFunctionReturn(0);
210e27a552bSJed Brown }
211e27a552bSJed Brown 
212e27a552bSJed Brown #undef __FUNCT__
213e27a552bSJed Brown #define __FUNCT__ "TSRosWInitializePackage"
214e27a552bSJed Brown /*@C
215e27a552bSJed Brown   TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called
216e27a552bSJed Brown   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RosW()
217e27a552bSJed Brown   when using static libraries.
218e27a552bSJed Brown 
219e27a552bSJed Brown   Input Parameter:
220e27a552bSJed Brown   path - The dynamic library path, or PETSC_NULL
221e27a552bSJed Brown 
222e27a552bSJed Brown   Level: developer
223e27a552bSJed Brown 
224e27a552bSJed Brown .keywords: TS, TSRosW, initialize, package
225e27a552bSJed Brown .seealso: PetscInitialize()
226e27a552bSJed Brown @*/
227e27a552bSJed Brown PetscErrorCode TSRosWInitializePackage(const char path[])
228e27a552bSJed Brown {
229e27a552bSJed Brown   PetscErrorCode ierr;
230e27a552bSJed Brown 
231e27a552bSJed Brown   PetscFunctionBegin;
232e27a552bSJed Brown   if (TSRosWPackageInitialized) PetscFunctionReturn(0);
233e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_TRUE;
234e27a552bSJed Brown   ierr = TSRosWRegisterAll();CHKERRQ(ierr);
235e27a552bSJed Brown   ierr = PetscRegisterFinalize(TSRosWFinalizePackage);CHKERRQ(ierr);
236e27a552bSJed Brown   PetscFunctionReturn(0);
237e27a552bSJed Brown }
238e27a552bSJed Brown 
239e27a552bSJed Brown #undef __FUNCT__
240e27a552bSJed Brown #define __FUNCT__ "TSRosWFinalizePackage"
241e27a552bSJed Brown /*@C
242e27a552bSJed Brown   TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is
243e27a552bSJed Brown   called from PetscFinalize().
244e27a552bSJed Brown 
245e27a552bSJed Brown   Level: developer
246e27a552bSJed Brown 
247e27a552bSJed Brown .keywords: Petsc, destroy, package
248e27a552bSJed Brown .seealso: PetscFinalize()
249e27a552bSJed Brown @*/
250e27a552bSJed Brown PetscErrorCode TSRosWFinalizePackage(void)
251e27a552bSJed Brown {
252e27a552bSJed Brown   PetscErrorCode ierr;
253e27a552bSJed Brown 
254e27a552bSJed Brown   PetscFunctionBegin;
255e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_FALSE;
256e27a552bSJed Brown   ierr = TSRosWRegisterDestroy();CHKERRQ(ierr);
257e27a552bSJed Brown   PetscFunctionReturn(0);
258e27a552bSJed Brown }
259e27a552bSJed Brown 
260e27a552bSJed Brown #undef __FUNCT__
261e27a552bSJed Brown #define __FUNCT__ "TSRosWRegister"
262e27a552bSJed Brown /*@C
26361692a83SJed Brown    TSRosWRegister - register a Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
264e27a552bSJed Brown 
265e27a552bSJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
266e27a552bSJed Brown 
267e27a552bSJed Brown    Input Parameters:
268e27a552bSJed Brown +  name - identifier for method
269e27a552bSJed Brown .  order - approximation order of method
270e27a552bSJed Brown .  s - number of stages, this is the dimension of the matrices below
27161692a83SJed Brown .  A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular
27261692a83SJed Brown .  Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal
273fe7e6d57SJed Brown .  b - Step completion table (dimension s)
274fe7e6d57SJed Brown -  bembed - Step completion table for a scheme of order one less (dimension s, PETSC_NULL if no embedded scheme is available)
275e27a552bSJed Brown 
276e27a552bSJed Brown    Notes:
27761692a83SJed Brown    Several Rosenbrock W methods are provided, this function is only needed to create new methods.
278e27a552bSJed Brown 
279e27a552bSJed Brown    Level: advanced
280e27a552bSJed Brown 
281e27a552bSJed Brown .keywords: TS, register
282e27a552bSJed Brown 
283e27a552bSJed Brown .seealso: TSRosW
284e27a552bSJed Brown @*/
285e27a552bSJed Brown PetscErrorCode TSRosWRegister(const TSRosWType name,PetscInt order,PetscInt s,
286fe7e6d57SJed Brown                               const PetscReal A[],const PetscReal Gamma[],const PetscReal b[],const PetscReal bembed[])
287e27a552bSJed Brown {
288e27a552bSJed Brown   PetscErrorCode ierr;
28961692a83SJed Brown   RosWTableauLink link;
29061692a83SJed Brown   RosWTableau t;
29161692a83SJed Brown   PetscInt i,j,k;
29261692a83SJed Brown   PetscScalar *GammaInv;
293e27a552bSJed Brown 
294e27a552bSJed Brown   PetscFunctionBegin;
295fe7e6d57SJed Brown   PetscValidCharPointer(name,1);
296fe7e6d57SJed Brown   PetscValidPointer(A,4);
297fe7e6d57SJed Brown   PetscValidPointer(Gamma,5);
298fe7e6d57SJed Brown   PetscValidPointer(b,6);
299fe7e6d57SJed Brown   if (bembed) PetscValidPointer(bembed,7);
300fe7e6d57SJed Brown 
301e27a552bSJed Brown   ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
302e27a552bSJed Brown   ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
303e27a552bSJed Brown   t = &link->tab;
304e27a552bSJed Brown   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
305e27a552bSJed Brown   t->order = order;
306e27a552bSJed Brown   t->s = s;
30761692a83SJed 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);
308c17803e7SJed Brown   ierr = PetscMalloc4(s*s,PetscReal,&t->At,s,PetscReal,&t->bt,s*s,PetscReal,&t->GammaInv,s,PetscBool,&t->GammaZeroDiag);CHKERRQ(ierr);
309e27a552bSJed Brown   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
31061692a83SJed Brown   ierr = PetscMemcpy(t->Gamma,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr);
31161692a83SJed Brown   ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
312fe7e6d57SJed Brown   if (bembed) {
313fe7e6d57SJed Brown     ierr = PetscMalloc2(s,PetscReal,&t->bembed,s,PetscReal,&t->bembedt);CHKERRQ(ierr);
314fe7e6d57SJed Brown     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
315fe7e6d57SJed Brown   }
31661692a83SJed Brown   for (i=0; i<s; i++) {
31761692a83SJed Brown     t->ASum[i] = 0;
31861692a83SJed Brown     t->GammaSum[i] = 0;
31961692a83SJed Brown     for (j=0; j<s; j++) {
32061692a83SJed Brown       t->ASum[i] += A[i*s+j];
321fe7e6d57SJed Brown       t->GammaSum[i] += Gamma[i*s+j];
32261692a83SJed Brown     }
32361692a83SJed Brown   }
32461692a83SJed Brown   ierr = PetscMalloc(s*s*sizeof(PetscScalar),&GammaInv);CHKERRQ(ierr); /* Need to use Scalar for inverse, then convert back to Real */
32561692a83SJed Brown   for (i=0; i<s*s; i++) GammaInv[i] = Gamma[i];
326fd96d5b0SEmil Constantinescu   for (i=0; i<s; i++) {
327fd96d5b0SEmil Constantinescu     if (Gamma[i*s+i] == 0.0) {
328fd96d5b0SEmil Constantinescu       GammaInv[i*s+i] = 1.0;
329c17803e7SJed Brown       t->GammaZeroDiag[i] = PETSC_TRUE;
330fd96d5b0SEmil Constantinescu     } else {
331c17803e7SJed Brown       t->GammaZeroDiag[i] = PETSC_FALSE;
332fd96d5b0SEmil Constantinescu     }
333fd96d5b0SEmil Constantinescu   }
334fd96d5b0SEmil Constantinescu 
33561692a83SJed Brown   switch (s) {
33661692a83SJed Brown   case 1: GammaInv[0] = 1./GammaInv[0]; break;
33761692a83SJed Brown   case 2: ierr = Kernel_A_gets_inverse_A_2(GammaInv,0);CHKERRQ(ierr); break;
33861692a83SJed Brown   case 3: ierr = Kernel_A_gets_inverse_A_3(GammaInv,0);CHKERRQ(ierr); break;
33961692a83SJed Brown   case 4: ierr = Kernel_A_gets_inverse_A_4(GammaInv,0);CHKERRQ(ierr); break;
34061692a83SJed Brown   case 5: {
34161692a83SJed Brown     PetscInt ipvt5[5];
34261692a83SJed Brown     MatScalar work5[5*5];
34361692a83SJed Brown     ierr = Kernel_A_gets_inverse_A_5(GammaInv,ipvt5,work5,0);CHKERRQ(ierr); break;
34461692a83SJed Brown   }
34561692a83SJed Brown   case 6: ierr = Kernel_A_gets_inverse_A_6(GammaInv,0);CHKERRQ(ierr); break;
34661692a83SJed Brown   case 7: ierr = Kernel_A_gets_inverse_A_7(GammaInv,0);CHKERRQ(ierr); break;
34761692a83SJed Brown   default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for %D stages",s);
34861692a83SJed Brown   }
34961692a83SJed Brown   for (i=0; i<s*s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]);
35061692a83SJed Brown   ierr = PetscFree(GammaInv);CHKERRQ(ierr);
35161692a83SJed Brown   for (i=0; i<s; i++) {
35261692a83SJed Brown     for (j=0; j<s; j++) {
35361692a83SJed Brown       t->At[i*s+j] = 0;
35461692a83SJed Brown       for (k=0; k<s; k++) {
35561692a83SJed Brown         t->At[i*s+j] += t->A[i*s+k] * t->GammaInv[k*s+j];
35661692a83SJed Brown       }
35761692a83SJed Brown     }
35861692a83SJed Brown     t->bt[i] = 0;
35961692a83SJed Brown     for (j=0; j<s; j++) {
36061692a83SJed Brown       t->bt[i] += t->b[j] * t->GammaInv[j*s+i];
36161692a83SJed Brown     }
362fe7e6d57SJed Brown     if (bembed) {
363fe7e6d57SJed Brown       t->bembedt[i] = 0;
364fe7e6d57SJed Brown       for (j=0; j<s; j++) {
365fe7e6d57SJed Brown         t->bembedt[i] += t->bembed[j] * t->GammaInv[j*s+i];
366fe7e6d57SJed Brown       }
367fe7e6d57SJed Brown     }
36861692a83SJed Brown   }
369*8d59e960SJed Brown   t->ccfl = 1.0;                /* Fix this */
370*8d59e960SJed Brown 
37161692a83SJed Brown   link->next = RosWTableauList;
37261692a83SJed Brown   RosWTableauList = link;
373e27a552bSJed Brown   PetscFunctionReturn(0);
374e27a552bSJed Brown }
375e27a552bSJed Brown 
376e27a552bSJed Brown #undef __FUNCT__
3771c3436cfSJed Brown #define __FUNCT__ "TSEvaluateStep_RosW"
3781c3436cfSJed Brown /*
3791c3436cfSJed Brown  The step completion formula is
3801c3436cfSJed Brown 
3811c3436cfSJed Brown  x1 = x0 + b^T Y
3821c3436cfSJed Brown 
3831c3436cfSJed Brown  where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been
3841c3436cfSJed Brown  updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write
3851c3436cfSJed Brown 
3861c3436cfSJed Brown  x1e = x0 + be^T Y
3871c3436cfSJed Brown      = x1 - b^T Y + be^T Y
3881c3436cfSJed Brown      = x1 + (be - b)^T Y
3891c3436cfSJed Brown 
3901c3436cfSJed Brown  so we can evaluate the method of different order even after the step has been optimistically completed.
3911c3436cfSJed Brown */
3921c3436cfSJed Brown static PetscErrorCode TSEvaluateStep_RosW(TS ts,PetscInt order,Vec X,PetscBool *done)
3931c3436cfSJed Brown {
3941c3436cfSJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
3951c3436cfSJed Brown   RosWTableau    tab  = ros->tableau;
3961c3436cfSJed Brown   PetscScalar    *w = ros->work;
3971c3436cfSJed Brown   PetscInt       i;
3981c3436cfSJed Brown   PetscErrorCode ierr;
3991c3436cfSJed Brown 
4001c3436cfSJed Brown   PetscFunctionBegin;
4011c3436cfSJed Brown   if (order == tab->order) {
4021c3436cfSJed Brown     if (ros->step_taken) {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
4031c3436cfSJed Brown     else {
4041c3436cfSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
4051c3436cfSJed Brown       ierr = VecMAXPY(X,tab->s,tab->bt,ros->Y);CHKERRQ(ierr);
4061c3436cfSJed Brown     }
4071c3436cfSJed Brown     if (done) *done = PETSC_TRUE;
4081c3436cfSJed Brown     PetscFunctionReturn(0);
4091c3436cfSJed Brown   } else if (order == tab->order-1) {
4101c3436cfSJed Brown     if (!tab->bembedt) goto unavailable;
4111c3436cfSJed Brown     if (ros->step_taken) {
4121c3436cfSJed Brown       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i];
4131c3436cfSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
4141c3436cfSJed Brown       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
4151c3436cfSJed Brown     } else {
4161c3436cfSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
4171c3436cfSJed Brown       ierr = VecMAXPY(X,tab->s,tab->bembedt,ros->Y);CHKERRQ(ierr);
4181c3436cfSJed Brown     }
4191c3436cfSJed Brown     if (done) *done = PETSC_TRUE;
4201c3436cfSJed Brown     PetscFunctionReturn(0);
4211c3436cfSJed Brown   }
4221c3436cfSJed Brown   unavailable:
4231c3436cfSJed Brown   if (done) *done = PETSC_FALSE;
4241c3436cfSJed 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);
4251c3436cfSJed Brown   PetscFunctionReturn(0);
4261c3436cfSJed Brown }
4271c3436cfSJed Brown 
4281c3436cfSJed Brown #undef __FUNCT__
429e27a552bSJed Brown #define __FUNCT__ "TSStep_RosW"
430e27a552bSJed Brown static PetscErrorCode TSStep_RosW(TS ts)
431e27a552bSJed Brown {
43261692a83SJed Brown   TS_RosW         *ros = (TS_RosW*)ts->data;
43361692a83SJed Brown   RosWTableau     tab  = ros->tableau;
434e27a552bSJed Brown   const PetscInt  s    = tab->s;
4351c3436cfSJed Brown   const PetscReal *At  = tab->At,*Gamma = tab->Gamma,*ASum = tab->ASum,*GammaInv = tab->GammaInv;
436c17803e7SJed Brown   const PetscBool *GammaZeroDiag = tab->GammaZeroDiag;
43761692a83SJed Brown   PetscScalar     *w   = ros->work;
43861692a83SJed Brown   Vec             *Y   = ros->Y,Zdot = ros->Zdot,Zstage = ros->Zstage;
439e27a552bSJed Brown   SNES            snes;
4401c3436cfSJed Brown   TSAdapt         adapt;
4411c3436cfSJed Brown   PetscInt        i,j,its,lits,reject,next_scheme;
442cdbf8f93SLisandro Dalcin   PetscReal       next_time_step;
4431c3436cfSJed Brown   PetscBool       accept;
444e27a552bSJed Brown   PetscErrorCode  ierr;
445e27a552bSJed Brown 
446e27a552bSJed Brown   PetscFunctionBegin;
447e27a552bSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
448cdbf8f93SLisandro Dalcin   next_time_step = ts->time_step;
4491c3436cfSJed Brown   accept = PETSC_TRUE;
4501c3436cfSJed Brown   ros->step_taken = PETSC_FALSE;
451e27a552bSJed Brown 
4521c3436cfSJed Brown   for (reject=0; reject<ts->max_reject; reject++,ts->reject++) {
4531c3436cfSJed Brown     const PetscReal h = ts->time_step;
454e27a552bSJed Brown     for (i=0; i<s; i++) {
4551c3436cfSJed Brown       ros->stage_time = ts->ptime + h*ASum[i];
456c17803e7SJed Brown       if (GammaZeroDiag[i]) {
457c17803e7SJed Brown         ros->stage_explicit = PETSC_TRUE;
458fd96d5b0SEmil Constantinescu         ros->shift = 1./h;
459c17803e7SJed Brown       } else {
460c17803e7SJed Brown         ros->stage_explicit = PETSC_FALSE;
46161692a83SJed Brown         ros->shift = 1./(h*Gamma[i*s+i]);
462fd96d5b0SEmil Constantinescu       }
46361692a83SJed Brown 
46461692a83SJed Brown       ierr = VecCopy(ts->vec_sol,Zstage);CHKERRQ(ierr);
46561692a83SJed Brown       ierr = VecMAXPY(Zstage,i,&At[i*s+0],Y);CHKERRQ(ierr);
46661692a83SJed Brown 
46761692a83SJed Brown       for (j=0; j<i; j++) w[j] = 1./h * GammaInv[i*s+j];
46861692a83SJed Brown       ierr = VecZeroEntries(Zdot);CHKERRQ(ierr);
46961692a83SJed Brown       ierr = VecMAXPY(Zdot,i,w,Y);CHKERRQ(ierr);
47061692a83SJed Brown 
471e27a552bSJed Brown       /* Initial guess taken from last stage */
47261692a83SJed Brown       ierr = VecZeroEntries(Y[i]);CHKERRQ(ierr);
47361692a83SJed Brown 
47461692a83SJed Brown       if (!ros->recompute_jacobian && !i) {
47561692a83SJed Brown         ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */
47661692a83SJed Brown       }
47761692a83SJed Brown 
47861692a83SJed Brown       ierr = SNESSolve(snes,PETSC_NULL,Y[i]);CHKERRQ(ierr);
479e27a552bSJed Brown       ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
480e27a552bSJed Brown       ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
481e27a552bSJed Brown       ts->nonlinear_its += its; ts->linear_its += lits;
482e27a552bSJed Brown     }
4831c3436cfSJed Brown     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
4841c3436cfSJed Brown     ros->step_taken = PETSC_TRUE;
485e27a552bSJed Brown 
4861c3436cfSJed Brown     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
4871c3436cfSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
4881c3436cfSJed Brown     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
489*8d59e960SJed Brown     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
4901c3436cfSJed Brown     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
4911c3436cfSJed Brown     if (accept) {
4921c3436cfSJed Brown       /* ignore next_scheme for now */
493e27a552bSJed Brown       ts->ptime += ts->time_step;
494cdbf8f93SLisandro Dalcin       ts->time_step = next_time_step;
495e27a552bSJed Brown       ts->steps++;
4961c3436cfSJed Brown       break;
4971c3436cfSJed Brown     } else {                    /* Roll back the current step */
4981c3436cfSJed Brown       for (i=0; i<s; i++) w[i] = -tab->bt[i];
4991c3436cfSJed Brown       ierr = VecMAXPY(ts->vec_sol,s,w,Y);CHKERRQ(ierr);
5001c3436cfSJed Brown       ts->time_step = next_time_step;
5011c3436cfSJed Brown       ros->step_taken = PETSC_FALSE;
5021c3436cfSJed Brown     }
5031c3436cfSJed Brown   }
5041c3436cfSJed Brown 
505e27a552bSJed Brown   PetscFunctionReturn(0);
506e27a552bSJed Brown }
507e27a552bSJed Brown 
508e27a552bSJed Brown #undef __FUNCT__
509e27a552bSJed Brown #define __FUNCT__ "TSInterpolate_RosW"
510e27a552bSJed Brown static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec X)
511e27a552bSJed Brown {
51261692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
513e27a552bSJed Brown 
514e27a552bSJed Brown   PetscFunctionBegin;
51561692a83SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name);
516e27a552bSJed Brown   PetscFunctionReturn(0);
517e27a552bSJed Brown }
518e27a552bSJed Brown 
519e27a552bSJed Brown /*------------------------------------------------------------*/
520e27a552bSJed Brown #undef __FUNCT__
521e27a552bSJed Brown #define __FUNCT__ "TSReset_RosW"
522e27a552bSJed Brown static PetscErrorCode TSReset_RosW(TS ts)
523e27a552bSJed Brown {
52461692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
525e27a552bSJed Brown   PetscInt       s;
526e27a552bSJed Brown   PetscErrorCode ierr;
527e27a552bSJed Brown 
528e27a552bSJed Brown   PetscFunctionBegin;
52961692a83SJed Brown   if (!ros->tableau) PetscFunctionReturn(0);
53061692a83SJed Brown    s = ros->tableau->s;
53161692a83SJed Brown   ierr = VecDestroyVecs(s,&ros->Y);CHKERRQ(ierr);
53261692a83SJed Brown   ierr = VecDestroy(&ros->Ydot);CHKERRQ(ierr);
53361692a83SJed Brown   ierr = VecDestroy(&ros->Ystage);CHKERRQ(ierr);
53461692a83SJed Brown   ierr = VecDestroy(&ros->Zdot);CHKERRQ(ierr);
53561692a83SJed Brown   ierr = VecDestroy(&ros->Zstage);CHKERRQ(ierr);
53661692a83SJed Brown   ierr = PetscFree(ros->work);CHKERRQ(ierr);
537e27a552bSJed Brown   PetscFunctionReturn(0);
538e27a552bSJed Brown }
539e27a552bSJed Brown 
540e27a552bSJed Brown #undef __FUNCT__
541e27a552bSJed Brown #define __FUNCT__ "TSDestroy_RosW"
542e27a552bSJed Brown static PetscErrorCode TSDestroy_RosW(TS ts)
543e27a552bSJed Brown {
544e27a552bSJed Brown   PetscErrorCode  ierr;
545e27a552bSJed Brown 
546e27a552bSJed Brown   PetscFunctionBegin;
547e27a552bSJed Brown   ierr = TSReset_RosW(ts);CHKERRQ(ierr);
548e27a552bSJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
549e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","",PETSC_NULL);CHKERRQ(ierr);
550e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","",PETSC_NULL);CHKERRQ(ierr);
55161692a83SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","",PETSC_NULL);CHKERRQ(ierr);
552e27a552bSJed Brown   PetscFunctionReturn(0);
553e27a552bSJed Brown }
554e27a552bSJed Brown 
555e27a552bSJed Brown /*
556e27a552bSJed Brown   This defines the nonlinear equation that is to be solved with SNES
557e27a552bSJed Brown   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
558e27a552bSJed Brown */
559e27a552bSJed Brown #undef __FUNCT__
560e27a552bSJed Brown #define __FUNCT__ "SNESTSFormFunction_RosW"
561e27a552bSJed Brown static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec X,Vec F,TS ts)
562e27a552bSJed Brown {
56361692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
564e27a552bSJed Brown   PetscErrorCode ierr;
565e27a552bSJed Brown 
566e27a552bSJed Brown   PetscFunctionBegin;
567c17803e7SJed Brown   if (ros->stage_explicit) {
568c17803e7SJed Brown     ierr = VecAXPBY(ros->Ydot,ros->shift,0.0,X);CHKERRQ(ierr);       /* Ydot = shift*X*/
569c17803e7SJed Brown   } else {
57061692a83SJed Brown     ierr = VecWAXPY(ros->Ydot,ros->shift,X,ros->Zdot);CHKERRQ(ierr); /* Ydot = shift*X + Zdot */
571c17803e7SJed Brown   }
57261692a83SJed Brown   ierr = VecWAXPY(ros->Ystage,1.0,X,ros->Zstage);CHKERRQ(ierr);    /* Ystage = X + Zstage */
57361692a83SJed Brown   ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,ros->Ydot,F,PETSC_FALSE);CHKERRQ(ierr);
574e27a552bSJed Brown   PetscFunctionReturn(0);
575e27a552bSJed Brown }
576e27a552bSJed Brown 
577e27a552bSJed Brown #undef __FUNCT__
578e27a552bSJed Brown #define __FUNCT__ "SNESTSFormJacobian_RosW"
579e27a552bSJed Brown static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
580e27a552bSJed Brown {
58161692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
582e27a552bSJed Brown   PetscErrorCode ierr;
583e27a552bSJed Brown 
584e27a552bSJed Brown   PetscFunctionBegin;
58561692a83SJed Brown   /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
58661692a83SJed Brown   ierr = TSComputeIJacobian(ts,ros->stage_time,ros->Ystage,ros->Ydot,ros->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr);
587e27a552bSJed Brown   PetscFunctionReturn(0);
588e27a552bSJed Brown }
589e27a552bSJed Brown 
590e27a552bSJed Brown #undef __FUNCT__
591e27a552bSJed Brown #define __FUNCT__ "TSSetUp_RosW"
592e27a552bSJed Brown static PetscErrorCode TSSetUp_RosW(TS ts)
593e27a552bSJed Brown {
59461692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
59561692a83SJed Brown   RosWTableau    tab  = ros->tableau;
596e27a552bSJed Brown   PetscInt       s    = tab->s;
597e27a552bSJed Brown   PetscErrorCode ierr;
598e27a552bSJed Brown 
599e27a552bSJed Brown   PetscFunctionBegin;
60061692a83SJed Brown   if (!ros->tableau) {
601e27a552bSJed Brown     ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);
602e27a552bSJed Brown   }
60361692a83SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ros->Y);CHKERRQ(ierr);
60461692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Ydot);CHKERRQ(ierr);
60561692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Ystage);CHKERRQ(ierr);
60661692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Zdot);CHKERRQ(ierr);
60761692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Zstage);CHKERRQ(ierr);
60861692a83SJed Brown   ierr = PetscMalloc(s*sizeof(ros->work[0]),&ros->work);CHKERRQ(ierr);
609e27a552bSJed Brown   PetscFunctionReturn(0);
610e27a552bSJed Brown }
611e27a552bSJed Brown /*------------------------------------------------------------*/
612e27a552bSJed Brown 
613e27a552bSJed Brown #undef __FUNCT__
614e27a552bSJed Brown #define __FUNCT__ "TSSetFromOptions_RosW"
615e27a552bSJed Brown static PetscErrorCode TSSetFromOptions_RosW(TS ts)
616e27a552bSJed Brown {
61761692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
618e27a552bSJed Brown   PetscErrorCode ierr;
61961692a83SJed Brown   char           rostype[256];
620e27a552bSJed Brown 
621e27a552bSJed Brown   PetscFunctionBegin;
622e27a552bSJed Brown   ierr = PetscOptionsHead("RosW ODE solver options");CHKERRQ(ierr);
623e27a552bSJed Brown   {
62461692a83SJed Brown     RosWTableauLink link;
625e27a552bSJed Brown     PetscInt count,choice;
626e27a552bSJed Brown     PetscBool flg;
627e27a552bSJed Brown     const char **namelist;
62861692a83SJed Brown     SNES snes;
62961692a83SJed Brown 
63061692a83SJed Brown     ierr = PetscStrncpy(rostype,TSRosWDefault,sizeof rostype);CHKERRQ(ierr);
63161692a83SJed Brown     for (link=RosWTableauList,count=0; link; link=link->next,count++) ;
632e27a552bSJed Brown     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
63361692a83SJed Brown     for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
63461692a83SJed Brown     ierr = PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,rostype,&choice,&flg);CHKERRQ(ierr);
63561692a83SJed Brown     ierr = TSRosWSetType(ts,flg ? namelist[choice] : rostype);CHKERRQ(ierr);
636e27a552bSJed Brown     ierr = PetscFree(namelist);CHKERRQ(ierr);
63761692a83SJed Brown 
63861692a83SJed Brown     ierr = PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,PETSC_NULL);CHKERRQ(ierr);
63961692a83SJed Brown 
64061692a83SJed Brown     /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
64161692a83SJed Brown     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
64261692a83SJed Brown     if (!((PetscObject)snes)->type_name) {
64361692a83SJed Brown       ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
64461692a83SJed Brown     }
64561692a83SJed Brown     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
646e27a552bSJed Brown   }
647e27a552bSJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
648e27a552bSJed Brown   PetscFunctionReturn(0);
649e27a552bSJed Brown }
650e27a552bSJed Brown 
651e27a552bSJed Brown #undef __FUNCT__
652e27a552bSJed Brown #define __FUNCT__ "PetscFormatRealArray"
653e27a552bSJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
654e27a552bSJed Brown {
655e27a552bSJed Brown   PetscErrorCode ierr;
656e408995aSJed Brown   PetscInt i;
657e408995aSJed Brown   size_t left,count;
658e27a552bSJed Brown   char *p;
659e27a552bSJed Brown 
660e27a552bSJed Brown   PetscFunctionBegin;
661e408995aSJed Brown   for (i=0,p=buf,left=len; i<n; i++) {
662e408995aSJed Brown     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
663e27a552bSJed Brown     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
664e27a552bSJed Brown     left -= count;
665e27a552bSJed Brown     p += count;
666e27a552bSJed Brown     *p++ = ' ';
667e27a552bSJed Brown   }
668e27a552bSJed Brown   p[i ? 0 : -1] = 0;
669e27a552bSJed Brown   PetscFunctionReturn(0);
670e27a552bSJed Brown }
671e27a552bSJed Brown 
672e27a552bSJed Brown #undef __FUNCT__
673e27a552bSJed Brown #define __FUNCT__ "TSView_RosW"
674e27a552bSJed Brown static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer)
675e27a552bSJed Brown {
67661692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
67761692a83SJed Brown   RosWTableau    tab  = ros->tableau;
678e27a552bSJed Brown   PetscBool      iascii;
679e27a552bSJed Brown   PetscErrorCode ierr;
680e27a552bSJed Brown 
681e27a552bSJed Brown   PetscFunctionBegin;
682e27a552bSJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
683e27a552bSJed Brown   if (iascii) {
68461692a83SJed Brown     const TSRosWType rostype;
685e408995aSJed Brown     PetscInt i;
686e408995aSJed Brown     PetscReal abscissa[512];
687e27a552bSJed Brown     char buf[512];
68861692a83SJed Brown     ierr = TSRosWGetType(ts,&rostype);CHKERRQ(ierr);
68961692a83SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Rosenbrock-W %s\n",rostype);CHKERRQ(ierr);
690e408995aSJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ASum);CHKERRQ(ierr);
69161692a83SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A       = %s\n",buf);CHKERRQ(ierr);
692e408995aSJed Brown     for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i];
693e408995aSJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,abscissa);CHKERRQ(ierr);
694e408995aSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A+Gamma = %s\n",buf);CHKERRQ(ierr);
695e27a552bSJed Brown   }
696e27a552bSJed Brown   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
697e27a552bSJed Brown   PetscFunctionReturn(0);
698e27a552bSJed Brown }
699e27a552bSJed Brown 
700e27a552bSJed Brown #undef __FUNCT__
701e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType"
702e27a552bSJed Brown /*@C
70361692a83SJed Brown   TSRosWSetType - Set the type of Rosenbrock-W scheme
704e27a552bSJed Brown 
705e27a552bSJed Brown   Logically collective
706e27a552bSJed Brown 
707e27a552bSJed Brown   Input Parameter:
708e27a552bSJed Brown +  ts - timestepping context
70961692a83SJed Brown -  rostype - type of Rosenbrock-W scheme
710e27a552bSJed Brown 
711e27a552bSJed Brown   Level: intermediate
712e27a552bSJed Brown 
713e27a552bSJed Brown .seealso: TSRosWGetType()
714e27a552bSJed Brown @*/
71561692a83SJed Brown PetscErrorCode TSRosWSetType(TS ts,const TSRosWType rostype)
716e27a552bSJed Brown {
717e27a552bSJed Brown   PetscErrorCode ierr;
718e27a552bSJed Brown 
719e27a552bSJed Brown   PetscFunctionBegin;
720e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
72161692a83SJed Brown   ierr = PetscTryMethod(ts,"TSRosWSetType_C",(TS,const TSRosWType),(ts,rostype));CHKERRQ(ierr);
722e27a552bSJed Brown   PetscFunctionReturn(0);
723e27a552bSJed Brown }
724e27a552bSJed Brown 
725e27a552bSJed Brown #undef __FUNCT__
726e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType"
727e27a552bSJed Brown /*@C
72861692a83SJed Brown   TSRosWGetType - Get the type of Rosenbrock-W scheme
729e27a552bSJed Brown 
730e27a552bSJed Brown   Logically collective
731e27a552bSJed Brown 
732e27a552bSJed Brown   Input Parameter:
733e27a552bSJed Brown .  ts - timestepping context
734e27a552bSJed Brown 
735e27a552bSJed Brown   Output Parameter:
73661692a83SJed Brown .  rostype - type of Rosenbrock-W scheme
737e27a552bSJed Brown 
738e27a552bSJed Brown   Level: intermediate
739e27a552bSJed Brown 
740e27a552bSJed Brown .seealso: TSRosWGetType()
741e27a552bSJed Brown @*/
74261692a83SJed Brown PetscErrorCode TSRosWGetType(TS ts,const TSRosWType *rostype)
743e27a552bSJed Brown {
744e27a552bSJed Brown   PetscErrorCode ierr;
745e27a552bSJed Brown 
746e27a552bSJed Brown   PetscFunctionBegin;
747e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
74861692a83SJed Brown   ierr = PetscUseMethod(ts,"TSRosWGetType_C",(TS,const TSRosWType*),(ts,rostype));CHKERRQ(ierr);
749e27a552bSJed Brown   PetscFunctionReturn(0);
750e27a552bSJed Brown }
751e27a552bSJed Brown 
752e27a552bSJed Brown #undef __FUNCT__
75361692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian"
754e27a552bSJed Brown /*@C
75561692a83SJed Brown   TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step.
756e27a552bSJed Brown 
757e27a552bSJed Brown   Logically collective
758e27a552bSJed Brown 
759e27a552bSJed Brown   Input Parameter:
760e27a552bSJed Brown +  ts - timestepping context
76161692a83SJed Brown -  flg - PETSC_TRUE to recompute the Jacobian at each stage
762e27a552bSJed Brown 
763e27a552bSJed Brown   Level: intermediate
764e27a552bSJed Brown 
765e27a552bSJed Brown .seealso: TSRosWGetType()
766e27a552bSJed Brown @*/
76761692a83SJed Brown PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg)
768e27a552bSJed Brown {
769e27a552bSJed Brown   PetscErrorCode ierr;
770e27a552bSJed Brown 
771e27a552bSJed Brown   PetscFunctionBegin;
772e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
77361692a83SJed Brown   ierr = PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
774e27a552bSJed Brown   PetscFunctionReturn(0);
775e27a552bSJed Brown }
776e27a552bSJed Brown 
777e27a552bSJed Brown EXTERN_C_BEGIN
778e27a552bSJed Brown #undef __FUNCT__
779e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType_RosW"
78061692a83SJed Brown PetscErrorCode  TSRosWGetType_RosW(TS ts,const TSRosWType *rostype)
781e27a552bSJed Brown {
78261692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
783e27a552bSJed Brown   PetscErrorCode ierr;
784e27a552bSJed Brown 
785e27a552bSJed Brown   PetscFunctionBegin;
78661692a83SJed Brown   if (!ros->tableau) {ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);}
78761692a83SJed Brown   *rostype = ros->tableau->name;
788e27a552bSJed Brown   PetscFunctionReturn(0);
789e27a552bSJed Brown }
790e27a552bSJed Brown #undef __FUNCT__
791e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType_RosW"
79261692a83SJed Brown PetscErrorCode  TSRosWSetType_RosW(TS ts,const TSRosWType rostype)
793e27a552bSJed Brown {
79461692a83SJed Brown   TS_RosW         *ros = (TS_RosW*)ts->data;
795e27a552bSJed Brown   PetscErrorCode  ierr;
796e27a552bSJed Brown   PetscBool       match;
79761692a83SJed Brown   RosWTableauLink link;
798e27a552bSJed Brown 
799e27a552bSJed Brown   PetscFunctionBegin;
80061692a83SJed Brown   if (ros->tableau) {
80161692a83SJed Brown     ierr = PetscStrcmp(ros->tableau->name,rostype,&match);CHKERRQ(ierr);
802e27a552bSJed Brown     if (match) PetscFunctionReturn(0);
803e27a552bSJed Brown   }
80461692a83SJed Brown   for (link = RosWTableauList; link; link=link->next) {
80561692a83SJed Brown     ierr = PetscStrcmp(link->tab.name,rostype,&match);CHKERRQ(ierr);
806e27a552bSJed Brown     if (match) {
807e27a552bSJed Brown       ierr = TSReset_RosW(ts);CHKERRQ(ierr);
80861692a83SJed Brown       ros->tableau = &link->tab;
809e27a552bSJed Brown       PetscFunctionReturn(0);
810e27a552bSJed Brown     }
811e27a552bSJed Brown   }
81261692a83SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype);
813e27a552bSJed Brown   PetscFunctionReturn(0);
814e27a552bSJed Brown }
81561692a83SJed Brown 
816e27a552bSJed Brown #undef __FUNCT__
81761692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian_RosW"
81861692a83SJed Brown PetscErrorCode  TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg)
819e27a552bSJed Brown {
82061692a83SJed Brown   TS_RosW *ros = (TS_RosW*)ts->data;
821e27a552bSJed Brown 
822e27a552bSJed Brown   PetscFunctionBegin;
82361692a83SJed Brown   ros->recompute_jacobian = flg;
824e27a552bSJed Brown   PetscFunctionReturn(0);
825e27a552bSJed Brown }
826e27a552bSJed Brown EXTERN_C_END
827e27a552bSJed Brown 
828e27a552bSJed Brown /* ------------------------------------------------------------ */
829e27a552bSJed Brown /*MC
830e27a552bSJed Brown       TSRosW - ODE solver using Rosenbrock-W schemes
831e27a552bSJed Brown 
832e27a552bSJed Brown   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
833e27a552bSJed Brown   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
834e27a552bSJed Brown   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
835e27a552bSJed Brown 
836e27a552bSJed Brown   Notes:
83761692a83SJed Brown   This method currently only works with autonomous ODE and DAE.
83861692a83SJed Brown 
83961692a83SJed Brown   Developer notes:
84061692a83SJed Brown   Rosenbrock-W methods are typically specified for autonomous ODE
84161692a83SJed Brown 
84261692a83SJed Brown $  xdot = f(x)
84361692a83SJed Brown 
84461692a83SJed Brown   by the stage equations
84561692a83SJed Brown 
84661692a83SJed Brown $  k_i = h f(x_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j
84761692a83SJed Brown 
84861692a83SJed Brown   and step completion formula
84961692a83SJed Brown 
85061692a83SJed Brown $  x_1 = x_0 + sum_j b_j k_j
85161692a83SJed Brown 
85261692a83SJed Brown   with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(x)
85361692a83SJed Brown   and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner,
85461692a83SJed Brown   we define new variables for the stage equations
85561692a83SJed Brown 
85661692a83SJed Brown $  y_i = gamma_ij k_j
85761692a83SJed Brown 
85861692a83SJed Brown   The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define
85961692a83SJed Brown 
86061692a83SJed Brown $  A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-i}
86161692a83SJed Brown 
86261692a83SJed Brown   to rewrite the method as
86361692a83SJed Brown 
86461692a83SJed 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
86561692a83SJed Brown $  x_1 = x_0 + sum_j bt_j y_j
86661692a83SJed Brown 
86761692a83SJed Brown    where we have introduced the mass matrix M. Continue by defining
86861692a83SJed Brown 
86961692a83SJed Brown $  ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j
87061692a83SJed Brown 
87161692a83SJed Brown    or, more compactly in tensor notation
87261692a83SJed Brown 
87361692a83SJed Brown $  Ydot = 1/h (Gamma^{-1} \otimes I) Y .
87461692a83SJed Brown 
87561692a83SJed Brown    Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current
87661692a83SJed Brown    stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the
87761692a83SJed Brown    equation
87861692a83SJed Brown 
87961692a83SJed Brown $  g(x_0 + sum_j a_ij y_j + y_i, ydot_i) = 0
88061692a83SJed Brown 
88161692a83SJed Brown    with initial guess y_i = 0.
882e27a552bSJed Brown 
883e27a552bSJed Brown   Level: beginner
884e27a552bSJed Brown 
885e27a552bSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSRosWRegister()
886e27a552bSJed Brown 
887e27a552bSJed Brown M*/
888e27a552bSJed Brown EXTERN_C_BEGIN
889e27a552bSJed Brown #undef __FUNCT__
890e27a552bSJed Brown #define __FUNCT__ "TSCreate_RosW"
891e27a552bSJed Brown PetscErrorCode  TSCreate_RosW(TS ts)
892e27a552bSJed Brown {
89361692a83SJed Brown   TS_RosW        *ros;
894e27a552bSJed Brown   PetscErrorCode ierr;
895e27a552bSJed Brown 
896e27a552bSJed Brown   PetscFunctionBegin;
897e27a552bSJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
898e27a552bSJed Brown   ierr = TSRosWInitializePackage(PETSC_NULL);CHKERRQ(ierr);
899e27a552bSJed Brown #endif
900e27a552bSJed Brown 
901e27a552bSJed Brown   ts->ops->reset          = TSReset_RosW;
902e27a552bSJed Brown   ts->ops->destroy        = TSDestroy_RosW;
903e27a552bSJed Brown   ts->ops->view           = TSView_RosW;
904e27a552bSJed Brown   ts->ops->setup          = TSSetUp_RosW;
905e27a552bSJed Brown   ts->ops->step           = TSStep_RosW;
906e27a552bSJed Brown   ts->ops->interpolate    = TSInterpolate_RosW;
9071c3436cfSJed Brown   ts->ops->evaluatestep   = TSEvaluateStep_RosW;
908e27a552bSJed Brown   ts->ops->setfromoptions = TSSetFromOptions_RosW;
909e27a552bSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_RosW;
910e27a552bSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_RosW;
911e27a552bSJed Brown 
91261692a83SJed Brown   ierr = PetscNewLog(ts,TS_RosW,&ros);CHKERRQ(ierr);
91361692a83SJed Brown   ts->data = (void*)ros;
914e27a552bSJed Brown 
915e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","TSRosWGetType_RosW",TSRosWGetType_RosW);CHKERRQ(ierr);
916e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","TSRosWSetType_RosW",TSRosWSetType_RosW);CHKERRQ(ierr);
91761692a83SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","TSRosWSetRecomputeJacobian_RosW",TSRosWSetRecomputeJacobian_RosW);CHKERRQ(ierr);
918e27a552bSJed Brown   PetscFunctionReturn(0);
919e27a552bSJed Brown }
920e27a552bSJed Brown EXTERN_C_END
921