xref: /petsc/src/ts/impls/rosw/rosw.c (revision c17803e75cc0ba2a59653a0c06a54a91a25fd724)
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 */
28*c17803e7SJed 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 */
37e27a552bSJed Brown };
3861692a83SJed Brown typedef struct _RosWTableauLink *RosWTableauLink;
3961692a83SJed Brown struct _RosWTableauLink {
4061692a83SJed Brown   struct _RosWTableau tab;
4161692a83SJed Brown   RosWTableauLink next;
42e27a552bSJed Brown };
4361692a83SJed Brown static RosWTableauLink RosWTableauList;
44e27a552bSJed Brown 
45e27a552bSJed Brown typedef struct {
4661692a83SJed Brown   RosWTableau tableau;
4761692a83SJed Brown   Vec         *Y;               /* States computed during the step, used to complete the step */
48e27a552bSJed Brown   Vec         Ydot;             /* Work vector holding Ydot during residual evaluation */
4961692a83SJed Brown   Vec         Ystage;           /* Work vector for the state value at each stage */
5061692a83SJed Brown   Vec         Zdot;             /* Ydot = Zdot + shift*Y */
5161692a83SJed Brown   Vec         Zstage;           /* Y = Zstage + Y */
521c3436cfSJed Brown   PetscScalar *work;            /* Scalar work space of length number of stages, used to prepare VecMAXPY() */
53e27a552bSJed Brown   PetscReal   shift;
54e27a552bSJed Brown   PetscReal   stage_time;
55*c17803e7SJed Brown   PetscReal   stage_explicit;     /* Flag indicates that the current stage is explicit */
5661692a83SJed Brown   PetscBool   recompute_jacobian; /* Recompute the Jacobian at each stage, default is to freeze the Jacobian at the start of each step */
571c3436cfSJed Brown   PetscBool   step_taken;         /* ts->vec_sol has been advanced to the end of the current time step */
58e27a552bSJed Brown } TS_RosW;
59e27a552bSJed Brown 
60fe7e6d57SJed Brown /*MC
61fe7e6d57SJed Brown      TSROSW2M - Two stage second order L-stable Rosenbrock-W scheme.
62fe7e6d57SJed Brown 
63fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2P.
64fe7e6d57SJed Brown 
65fe7e6d57SJed Brown      Level: intermediate
66fe7e6d57SJed Brown 
67fe7e6d57SJed Brown .seealso: TSROSW
68fe7e6d57SJed Brown M*/
69fe7e6d57SJed Brown 
70fe7e6d57SJed Brown /*MC
71fe7e6d57SJed Brown      TSROSW2P - Two stage second order L-stable Rosenbrock-W scheme.
72fe7e6d57SJed Brown 
73fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2M.
74fe7e6d57SJed Brown 
75fe7e6d57SJed Brown      Level: intermediate
76fe7e6d57SJed Brown 
77fe7e6d57SJed Brown .seealso: TSROSW
78fe7e6d57SJed Brown M*/
79fe7e6d57SJed Brown 
80fe7e6d57SJed Brown /*MC
81fe7e6d57SJed Brown      TSROSWRA3PW - Three stage third order Rosenbrock-W scheme for PDAE of index 1.
82fe7e6d57SJed Brown 
83fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.
84fe7e6d57SJed Brown 
85fe7e6d57SJed 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.
86fe7e6d57SJed Brown 
87fe7e6d57SJed Brown      References:
88fe7e6d57SJed Brown      Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005.
89fe7e6d57SJed Brown 
90fe7e6d57SJed Brown      Level: intermediate
91fe7e6d57SJed Brown 
92fe7e6d57SJed Brown .seealso: TSROSW
93fe7e6d57SJed Brown M*/
94fe7e6d57SJed Brown 
95fe7e6d57SJed Brown /*MC
96fe7e6d57SJed Brown      TSROSWRA34PW2 - Four stage third order L-stable Rosenbrock-W scheme for PDAE of index 1.
97fe7e6d57SJed Brown 
98fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.
99fe7e6d57SJed Brown 
100fe7e6d57SJed 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.
101fe7e6d57SJed Brown 
102fe7e6d57SJed Brown      References:
103fe7e6d57SJed Brown      Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005.
104fe7e6d57SJed Brown 
105fe7e6d57SJed Brown      Level: intermediate
106fe7e6d57SJed Brown 
107fe7e6d57SJed Brown .seealso: TSROSW
108fe7e6d57SJed Brown M*/
109fe7e6d57SJed Brown 
110e27a552bSJed Brown #undef __FUNCT__
111e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterAll"
112e27a552bSJed Brown /*@C
113e27a552bSJed Brown   TSRosWRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSRosW
114e27a552bSJed Brown 
115e27a552bSJed Brown   Not Collective, but should be called by all processes which will need the schemes to be registered
116e27a552bSJed Brown 
117e27a552bSJed Brown   Level: advanced
118e27a552bSJed Brown 
119e27a552bSJed Brown .keywords: TS, TSRosW, register, all
120e27a552bSJed Brown 
121e27a552bSJed Brown .seealso:  TSRosWRegisterDestroy()
122e27a552bSJed Brown @*/
123e27a552bSJed Brown PetscErrorCode TSRosWRegisterAll(void)
124e27a552bSJed Brown {
125e27a552bSJed Brown   PetscErrorCode ierr;
126e27a552bSJed Brown 
127e27a552bSJed Brown   PetscFunctionBegin;
128e27a552bSJed Brown   if (TSRosWRegisterAllCalled) PetscFunctionReturn(0);
129e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_TRUE;
130e27a552bSJed Brown 
131e27a552bSJed Brown   {
13261692a83SJed Brown     const PetscReal g = 1. + 1./PetscSqrtReal(2.0);
133e27a552bSJed Brown     const PetscReal
13461692a83SJed Brown       A[2][2] = {{0,0}, {1.,0}},
13561692a83SJed Brown       Gamma[2][2] = {{g,0}, {-2.*g,g}},
1361c3436cfSJed Brown       b[2] = {0.5,0.5},
1371c3436cfSJed Brown       b1[2] = {1.0,0.0};
1381c3436cfSJed Brown     ierr = TSRosWRegister(TSROSW2P,2,2,&A[0][0],&Gamma[0][0],b,b1);CHKERRQ(ierr);
139e27a552bSJed Brown   }
140e27a552bSJed Brown   {
14161692a83SJed Brown     const PetscReal g = 1. - 1./PetscSqrtReal(2.0);
142e27a552bSJed Brown     const PetscReal
14361692a83SJed Brown       A[2][2] = {{0,0}, {1.,0}},
14461692a83SJed Brown       Gamma[2][2] = {{g,0}, {-2.*g,g}},
1451c3436cfSJed Brown       b[2] = {0.5,0.5},
1461c3436cfSJed Brown       b1[2] = {1.0,0.0};
1471c3436cfSJed Brown     ierr = TSRosWRegister(TSROSW2M,2,2,&A[0][0],&Gamma[0][0],b,b1);CHKERRQ(ierr);
148fe7e6d57SJed Brown   }
149fe7e6d57SJed Brown   {
150fe7e6d57SJed Brown     const PetscReal g = 7.8867513459481287e-01;
151fe7e6d57SJed Brown     const PetscReal
152fe7e6d57SJed Brown       A[3][3] = {{0,0,0},
153fe7e6d57SJed Brown                  {1.5773502691896257e+00,0,0},
154fe7e6d57SJed Brown                  {0.5,0,0}},
155fe7e6d57SJed Brown       Gamma[3][3] = {{g,0,0},
156fe7e6d57SJed Brown                      {-1.5773502691896257e+00,g,0},
157fe7e6d57SJed Brown                      {-6.7075317547305480e-01,1.7075317547305482e-01,g}},
158fe7e6d57SJed Brown       b[3] = {1.0566243270259355e-01,4.9038105676657971e-02,8.4529946162074843e-01},
159fe7e6d57SJed Brown       b2[3] = {-1.7863279495408180e-01,1./3.,8.4529946162074843e-01};
160fe7e6d57SJed Brown     ierr = TSRosWRegister(TSROSWRA3PW,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
161fe7e6d57SJed Brown   }
162fe7e6d57SJed Brown   {
163fe7e6d57SJed Brown     const PetscReal g = 4.3586652150845900e-01;
164fe7e6d57SJed Brown     const PetscReal
165fe7e6d57SJed Brown       A[4][4] = {{0,0,0,0},
166fe7e6d57SJed Brown                  {8.7173304301691801e-01,0,0,0},
167fe7e6d57SJed Brown                  {8.4457060015369423e-01,-1.1299064236484185e-01,0,0},
168fe7e6d57SJed Brown                  {0,0,1.,0}},
169fe7e6d57SJed Brown       Gamma[4][4] = {{g,0,0,0},
170fe7e6d57SJed Brown                      {-8.7173304301691801e-01,g,0,0},
171fe7e6d57SJed Brown                      {-9.0338057013044082e-01,5.4180672388095326e-02,g,0},
172fe7e6d57SJed Brown                      {2.4212380706095346e-01,-1.2232505839045147e+00,5.4526025533510214e-01,g}},
173fe7e6d57SJed Brown       b[4] = {2.4212380706095346e-01,-1.2232505839045147e+00,1.5452602553351020e+00,4.3586652150845900e-01},
174fe7e6d57SJed Brown       b2[4] = {3.7810903145819369e-01,-9.6042292212423178e-02,5.0000000000000000e-01,2.1793326075422950e-01};
175fe7e6d57SJed Brown     ierr = TSRosWRegister(TSROSWRA34PW2,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
176e27a552bSJed Brown   }
177e27a552bSJed Brown   PetscFunctionReturn(0);
178e27a552bSJed Brown }
179e27a552bSJed Brown 
180e27a552bSJed Brown #undef __FUNCT__
181e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterDestroy"
182e27a552bSJed Brown /*@C
183e27a552bSJed Brown    TSRosWRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister().
184e27a552bSJed Brown 
185e27a552bSJed Brown    Not Collective
186e27a552bSJed Brown 
187e27a552bSJed Brown    Level: advanced
188e27a552bSJed Brown 
189e27a552bSJed Brown .keywords: TSRosW, register, destroy
190e27a552bSJed Brown .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegisterDynamic()
191e27a552bSJed Brown @*/
192e27a552bSJed Brown PetscErrorCode TSRosWRegisterDestroy(void)
193e27a552bSJed Brown {
194e27a552bSJed Brown   PetscErrorCode ierr;
19561692a83SJed Brown   RosWTableauLink link;
196e27a552bSJed Brown 
197e27a552bSJed Brown   PetscFunctionBegin;
19861692a83SJed Brown   while ((link = RosWTableauList)) {
19961692a83SJed Brown     RosWTableau t = &link->tab;
20061692a83SJed Brown     RosWTableauList = link->next;
20161692a83SJed Brown     ierr = PetscFree5(t->A,t->Gamma,t->b,t->ASum,t->GammaSum);CHKERRQ(ierr);
202*c17803e7SJed Brown     ierr = PetscFree4(t->At,t->bt,t->GammaInv,t->GammaZeroDiag);CHKERRQ(ierr);
203fe7e6d57SJed Brown     ierr = PetscFree2(t->bembed,t->bembedt);CHKERRQ(ierr);
204e27a552bSJed Brown     ierr = PetscFree(t->name);CHKERRQ(ierr);
205e27a552bSJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
206e27a552bSJed Brown   }
207e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_FALSE;
208e27a552bSJed Brown   PetscFunctionReturn(0);
209e27a552bSJed Brown }
210e27a552bSJed Brown 
211e27a552bSJed Brown #undef __FUNCT__
212e27a552bSJed Brown #define __FUNCT__ "TSRosWInitializePackage"
213e27a552bSJed Brown /*@C
214e27a552bSJed Brown   TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called
215e27a552bSJed Brown   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RosW()
216e27a552bSJed Brown   when using static libraries.
217e27a552bSJed Brown 
218e27a552bSJed Brown   Input Parameter:
219e27a552bSJed Brown   path - The dynamic library path, or PETSC_NULL
220e27a552bSJed Brown 
221e27a552bSJed Brown   Level: developer
222e27a552bSJed Brown 
223e27a552bSJed Brown .keywords: TS, TSRosW, initialize, package
224e27a552bSJed Brown .seealso: PetscInitialize()
225e27a552bSJed Brown @*/
226e27a552bSJed Brown PetscErrorCode TSRosWInitializePackage(const char path[])
227e27a552bSJed Brown {
228e27a552bSJed Brown   PetscErrorCode ierr;
229e27a552bSJed Brown 
230e27a552bSJed Brown   PetscFunctionBegin;
231e27a552bSJed Brown   if (TSRosWPackageInitialized) PetscFunctionReturn(0);
232e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_TRUE;
233e27a552bSJed Brown   ierr = TSRosWRegisterAll();CHKERRQ(ierr);
234e27a552bSJed Brown   ierr = PetscRegisterFinalize(TSRosWFinalizePackage);CHKERRQ(ierr);
235e27a552bSJed Brown   PetscFunctionReturn(0);
236e27a552bSJed Brown }
237e27a552bSJed Brown 
238e27a552bSJed Brown #undef __FUNCT__
239e27a552bSJed Brown #define __FUNCT__ "TSRosWFinalizePackage"
240e27a552bSJed Brown /*@C
241e27a552bSJed Brown   TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is
242e27a552bSJed Brown   called from PetscFinalize().
243e27a552bSJed Brown 
244e27a552bSJed Brown   Level: developer
245e27a552bSJed Brown 
246e27a552bSJed Brown .keywords: Petsc, destroy, package
247e27a552bSJed Brown .seealso: PetscFinalize()
248e27a552bSJed Brown @*/
249e27a552bSJed Brown PetscErrorCode TSRosWFinalizePackage(void)
250e27a552bSJed Brown {
251e27a552bSJed Brown   PetscErrorCode ierr;
252e27a552bSJed Brown 
253e27a552bSJed Brown   PetscFunctionBegin;
254e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_FALSE;
255e27a552bSJed Brown   ierr = TSRosWRegisterDestroy();CHKERRQ(ierr);
256e27a552bSJed Brown   PetscFunctionReturn(0);
257e27a552bSJed Brown }
258e27a552bSJed Brown 
259e27a552bSJed Brown #undef __FUNCT__
260e27a552bSJed Brown #define __FUNCT__ "TSRosWRegister"
261e27a552bSJed Brown /*@C
26261692a83SJed Brown    TSRosWRegister - register a Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
263e27a552bSJed Brown 
264e27a552bSJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
265e27a552bSJed Brown 
266e27a552bSJed Brown    Input Parameters:
267e27a552bSJed Brown +  name - identifier for method
268e27a552bSJed Brown .  order - approximation order of method
269e27a552bSJed Brown .  s - number of stages, this is the dimension of the matrices below
27061692a83SJed Brown .  A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular
27161692a83SJed Brown .  Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal
272fe7e6d57SJed Brown .  b - Step completion table (dimension s)
273fe7e6d57SJed Brown -  bembed - Step completion table for a scheme of order one less (dimension s, PETSC_NULL if no embedded scheme is available)
274e27a552bSJed Brown 
275e27a552bSJed Brown    Notes:
27661692a83SJed Brown    Several Rosenbrock W methods are provided, this function is only needed to create new methods.
277e27a552bSJed Brown 
278e27a552bSJed Brown    Level: advanced
279e27a552bSJed Brown 
280e27a552bSJed Brown .keywords: TS, register
281e27a552bSJed Brown 
282e27a552bSJed Brown .seealso: TSRosW
283e27a552bSJed Brown @*/
284e27a552bSJed Brown PetscErrorCode TSRosWRegister(const TSRosWType name,PetscInt order,PetscInt s,
285fe7e6d57SJed Brown                               const PetscReal A[],const PetscReal Gamma[],const PetscReal b[],const PetscReal bembed[])
286e27a552bSJed Brown {
287e27a552bSJed Brown   PetscErrorCode ierr;
28861692a83SJed Brown   RosWTableauLink link;
28961692a83SJed Brown   RosWTableau t;
29061692a83SJed Brown   PetscInt i,j,k;
29161692a83SJed Brown   PetscScalar *GammaInv;
292e27a552bSJed Brown 
293e27a552bSJed Brown   PetscFunctionBegin;
294fe7e6d57SJed Brown   PetscValidCharPointer(name,1);
295fe7e6d57SJed Brown   PetscValidPointer(A,4);
296fe7e6d57SJed Brown   PetscValidPointer(Gamma,5);
297fe7e6d57SJed Brown   PetscValidPointer(b,6);
298fe7e6d57SJed Brown   if (bembed) PetscValidPointer(bembed,7);
299fe7e6d57SJed Brown 
300e27a552bSJed Brown   ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
301e27a552bSJed Brown   ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
302e27a552bSJed Brown   t = &link->tab;
303e27a552bSJed Brown   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
304e27a552bSJed Brown   t->order = order;
305e27a552bSJed Brown   t->s = s;
30661692a83SJed 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);
307*c17803e7SJed Brown   ierr = PetscMalloc4(s*s,PetscReal,&t->At,s,PetscReal,&t->bt,s*s,PetscReal,&t->GammaInv,s,PetscBool,&t->GammaZeroDiag);CHKERRQ(ierr);
308e27a552bSJed Brown   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
30961692a83SJed Brown   ierr = PetscMemcpy(t->Gamma,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr);
31061692a83SJed Brown   ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
311fe7e6d57SJed Brown   if (bembed) {
312fe7e6d57SJed Brown     ierr = PetscMalloc2(s,PetscReal,&t->bembed,s,PetscReal,&t->bembedt);CHKERRQ(ierr);
313fe7e6d57SJed Brown     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
314fe7e6d57SJed Brown   }
31561692a83SJed Brown   for (i=0; i<s; i++) {
31661692a83SJed Brown     t->ASum[i] = 0;
31761692a83SJed Brown     t->GammaSum[i] = 0;
31861692a83SJed Brown     for (j=0; j<s; j++) {
31961692a83SJed Brown       t->ASum[i] += A[i*s+j];
320fe7e6d57SJed Brown       t->GammaSum[i] += Gamma[i*s+j];
32161692a83SJed Brown     }
32261692a83SJed Brown   }
32361692a83SJed Brown   ierr = PetscMalloc(s*s*sizeof(PetscScalar),&GammaInv);CHKERRQ(ierr); /* Need to use Scalar for inverse, then convert back to Real */
32461692a83SJed Brown   for (i=0; i<s*s; i++) GammaInv[i] = Gamma[i];
325fd96d5b0SEmil Constantinescu   for (i=0; i<s; i++) {
326fd96d5b0SEmil Constantinescu     if (Gamma[i*s+i] == 0.0) {
327fd96d5b0SEmil Constantinescu       GammaInv[i*s+i] = 1.0;
328*c17803e7SJed Brown       t->GammaZeroDiag[i] = PETSC_TRUE;
329fd96d5b0SEmil Constantinescu     } else {
330*c17803e7SJed Brown       t->GammaZeroDiag[i] = PETSC_FALSE;
331fd96d5b0SEmil Constantinescu     }
332fd96d5b0SEmil Constantinescu   }
333fd96d5b0SEmil Constantinescu 
33461692a83SJed Brown   switch (s) {
33561692a83SJed Brown   case 1: GammaInv[0] = 1./GammaInv[0]; break;
33661692a83SJed Brown   case 2: ierr = Kernel_A_gets_inverse_A_2(GammaInv,0);CHKERRQ(ierr); break;
33761692a83SJed Brown   case 3: ierr = Kernel_A_gets_inverse_A_3(GammaInv,0);CHKERRQ(ierr); break;
33861692a83SJed Brown   case 4: ierr = Kernel_A_gets_inverse_A_4(GammaInv,0);CHKERRQ(ierr); break;
33961692a83SJed Brown   case 5: {
34061692a83SJed Brown     PetscInt ipvt5[5];
34161692a83SJed Brown     MatScalar work5[5*5];
34261692a83SJed Brown     ierr = Kernel_A_gets_inverse_A_5(GammaInv,ipvt5,work5,0);CHKERRQ(ierr); break;
34361692a83SJed Brown   }
34461692a83SJed Brown   case 6: ierr = Kernel_A_gets_inverse_A_6(GammaInv,0);CHKERRQ(ierr); break;
34561692a83SJed Brown   case 7: ierr = Kernel_A_gets_inverse_A_7(GammaInv,0);CHKERRQ(ierr); break;
34661692a83SJed Brown   default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for %D stages",s);
34761692a83SJed Brown   }
34861692a83SJed Brown   for (i=0; i<s*s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]);
34961692a83SJed Brown   ierr = PetscFree(GammaInv);CHKERRQ(ierr);
35061692a83SJed Brown   for (i=0; i<s; i++) {
35161692a83SJed Brown     for (j=0; j<s; j++) {
35261692a83SJed Brown       t->At[i*s+j] = 0;
35361692a83SJed Brown       for (k=0; k<s; k++) {
35461692a83SJed Brown         t->At[i*s+j] += t->A[i*s+k] * t->GammaInv[k*s+j];
35561692a83SJed Brown       }
35661692a83SJed Brown     }
35761692a83SJed Brown     t->bt[i] = 0;
35861692a83SJed Brown     for (j=0; j<s; j++) {
35961692a83SJed Brown       t->bt[i] += t->b[j] * t->GammaInv[j*s+i];
36061692a83SJed Brown     }
361fe7e6d57SJed Brown     if (bembed) {
362fe7e6d57SJed Brown       t->bembedt[i] = 0;
363fe7e6d57SJed Brown       for (j=0; j<s; j++) {
364fe7e6d57SJed Brown         t->bembedt[i] += t->bembed[j] * t->GammaInv[j*s+i];
365fe7e6d57SJed Brown       }
366fe7e6d57SJed Brown     }
36761692a83SJed Brown   }
36861692a83SJed Brown   link->next = RosWTableauList;
36961692a83SJed Brown   RosWTableauList = link;
370e27a552bSJed Brown   PetscFunctionReturn(0);
371e27a552bSJed Brown }
372e27a552bSJed Brown 
373e27a552bSJed Brown #undef __FUNCT__
3741c3436cfSJed Brown #define __FUNCT__ "TSEvaluateStep_RosW"
3751c3436cfSJed Brown /*
3761c3436cfSJed Brown  The step completion formula is
3771c3436cfSJed Brown 
3781c3436cfSJed Brown  x1 = x0 + b^T Y
3791c3436cfSJed Brown 
3801c3436cfSJed Brown  where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been
3811c3436cfSJed Brown  updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write
3821c3436cfSJed Brown 
3831c3436cfSJed Brown  x1e = x0 + be^T Y
3841c3436cfSJed Brown      = x1 - b^T Y + be^T Y
3851c3436cfSJed Brown      = x1 + (be - b)^T Y
3861c3436cfSJed Brown 
3871c3436cfSJed Brown  so we can evaluate the method of different order even after the step has been optimistically completed.
3881c3436cfSJed Brown */
3891c3436cfSJed Brown static PetscErrorCode TSEvaluateStep_RosW(TS ts,PetscInt order,Vec X,PetscBool *done)
3901c3436cfSJed Brown {
3911c3436cfSJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
3921c3436cfSJed Brown   RosWTableau    tab  = ros->tableau;
3931c3436cfSJed Brown   PetscScalar    *w = ros->work;
3941c3436cfSJed Brown   PetscInt       i;
3951c3436cfSJed Brown   PetscErrorCode ierr;
3961c3436cfSJed Brown 
3971c3436cfSJed Brown   PetscFunctionBegin;
3981c3436cfSJed Brown   if (order == tab->order) {
3991c3436cfSJed Brown     if (ros->step_taken) {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
4001c3436cfSJed Brown     else {
4011c3436cfSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
4021c3436cfSJed Brown       ierr = VecMAXPY(X,tab->s,tab->bt,ros->Y);CHKERRQ(ierr);
4031c3436cfSJed Brown     }
4041c3436cfSJed Brown     if (done) *done = PETSC_TRUE;
4051c3436cfSJed Brown     PetscFunctionReturn(0);
4061c3436cfSJed Brown   } else if (order == tab->order-1) {
4071c3436cfSJed Brown     if (!tab->bembedt) goto unavailable;
4081c3436cfSJed Brown     if (ros->step_taken) {
4091c3436cfSJed Brown       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i];
4101c3436cfSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
4111c3436cfSJed Brown       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
4121c3436cfSJed Brown     } else {
4131c3436cfSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
4141c3436cfSJed Brown       ierr = VecMAXPY(X,tab->s,tab->bembedt,ros->Y);CHKERRQ(ierr);
4151c3436cfSJed Brown     }
4161c3436cfSJed Brown     if (done) *done = PETSC_TRUE;
4171c3436cfSJed Brown     PetscFunctionReturn(0);
4181c3436cfSJed Brown   }
4191c3436cfSJed Brown   unavailable:
4201c3436cfSJed Brown   if (done) *done = PETSC_FALSE;
4211c3436cfSJed 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);
4221c3436cfSJed Brown   PetscFunctionReturn(0);
4231c3436cfSJed Brown }
4241c3436cfSJed Brown 
4251c3436cfSJed Brown #undef __FUNCT__
426e27a552bSJed Brown #define __FUNCT__ "TSStep_RosW"
427e27a552bSJed Brown static PetscErrorCode TSStep_RosW(TS ts)
428e27a552bSJed Brown {
42961692a83SJed Brown   TS_RosW         *ros = (TS_RosW*)ts->data;
43061692a83SJed Brown   RosWTableau     tab  = ros->tableau;
431e27a552bSJed Brown   const PetscInt  s    = tab->s;
4321c3436cfSJed Brown   const PetscReal *At  = tab->At,*Gamma = tab->Gamma,*ASum = tab->ASum,*GammaInv = tab->GammaInv;
433*c17803e7SJed Brown   const PetscBool *GammaZeroDiag = tab->GammaZeroDiag;
43461692a83SJed Brown   PetscScalar     *w   = ros->work;
43561692a83SJed Brown   Vec             *Y   = ros->Y,Zdot = ros->Zdot,Zstage = ros->Zstage;
436e27a552bSJed Brown   SNES            snes;
4371c3436cfSJed Brown   TSAdapt         adapt;
4381c3436cfSJed Brown   PetscInt        i,j,its,lits,reject,next_scheme;
439cdbf8f93SLisandro Dalcin   PetscReal       next_time_step;
4401c3436cfSJed Brown   PetscBool       accept;
441e27a552bSJed Brown   PetscErrorCode  ierr;
442e27a552bSJed Brown 
443e27a552bSJed Brown   PetscFunctionBegin;
444e27a552bSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
445cdbf8f93SLisandro Dalcin   next_time_step = ts->time_step;
4461c3436cfSJed Brown   accept = PETSC_TRUE;
4471c3436cfSJed Brown   ros->step_taken = PETSC_FALSE;
448e27a552bSJed Brown 
4491c3436cfSJed Brown   for (reject=0; reject<ts->max_reject; reject++,ts->reject++) {
4501c3436cfSJed Brown     const PetscReal h = ts->time_step;
451e27a552bSJed Brown     for (i=0; i<s; i++) {
4521c3436cfSJed Brown       ros->stage_time = ts->ptime + h*ASum[i];
453*c17803e7SJed Brown       if (GammaZeroDiag[i]) {
454*c17803e7SJed Brown         ros->stage_explicit = PETSC_TRUE;
455fd96d5b0SEmil Constantinescu         ros->shift = 1./h;
456*c17803e7SJed Brown       } else {
457*c17803e7SJed Brown         ros->stage_explicit = PETSC_FALSE;
45861692a83SJed Brown         ros->shift = 1./(h*Gamma[i*s+i]);
459fd96d5b0SEmil Constantinescu       }
46061692a83SJed Brown 
46161692a83SJed Brown       ierr = VecCopy(ts->vec_sol,Zstage);CHKERRQ(ierr);
46261692a83SJed Brown       ierr = VecMAXPY(Zstage,i,&At[i*s+0],Y);CHKERRQ(ierr);
46361692a83SJed Brown 
46461692a83SJed Brown       for (j=0; j<i; j++) w[j] = 1./h * GammaInv[i*s+j];
46561692a83SJed Brown       ierr = VecZeroEntries(Zdot);CHKERRQ(ierr);
46661692a83SJed Brown       ierr = VecMAXPY(Zdot,i,w,Y);CHKERRQ(ierr);
46761692a83SJed Brown 
468e27a552bSJed Brown       /* Initial guess taken from last stage */
46961692a83SJed Brown       ierr = VecZeroEntries(Y[i]);CHKERRQ(ierr);
47061692a83SJed Brown 
47161692a83SJed Brown       if (!ros->recompute_jacobian && !i) {
47261692a83SJed Brown         ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */
47361692a83SJed Brown       }
47461692a83SJed Brown 
47561692a83SJed Brown       ierr = SNESSolve(snes,PETSC_NULL,Y[i]);CHKERRQ(ierr);
476e27a552bSJed Brown       ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
477e27a552bSJed Brown       ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
478e27a552bSJed Brown       ts->nonlinear_its += its; ts->linear_its += lits;
479e27a552bSJed Brown     }
4801c3436cfSJed Brown     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
4811c3436cfSJed Brown     ros->step_taken = PETSC_TRUE;
482e27a552bSJed Brown 
4831c3436cfSJed Brown     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
4841c3436cfSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
4851c3436cfSJed Brown     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
4861c3436cfSJed Brown     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,0.0,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
4871c3436cfSJed Brown     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
4881c3436cfSJed Brown     if (accept) {
4891c3436cfSJed Brown       /* ignore next_scheme for now */
490e27a552bSJed Brown       ts->ptime += ts->time_step;
491cdbf8f93SLisandro Dalcin       ts->time_step = next_time_step;
492e27a552bSJed Brown       ts->steps++;
4931c3436cfSJed Brown       break;
4941c3436cfSJed Brown     } else {                    /* Roll back the current step */
4951c3436cfSJed Brown       for (i=0; i<s; i++) w[i] = -tab->bt[i];
4961c3436cfSJed Brown       ierr = VecMAXPY(ts->vec_sol,s,w,Y);CHKERRQ(ierr);
4971c3436cfSJed Brown       ts->time_step = next_time_step;
4981c3436cfSJed Brown       ros->step_taken = PETSC_FALSE;
4991c3436cfSJed Brown     }
5001c3436cfSJed Brown   }
5011c3436cfSJed Brown 
502e27a552bSJed Brown   PetscFunctionReturn(0);
503e27a552bSJed Brown }
504e27a552bSJed Brown 
505e27a552bSJed Brown #undef __FUNCT__
506e27a552bSJed Brown #define __FUNCT__ "TSInterpolate_RosW"
507e27a552bSJed Brown static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec X)
508e27a552bSJed Brown {
50961692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
510e27a552bSJed Brown 
511e27a552bSJed Brown   PetscFunctionBegin;
51261692a83SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name);
513e27a552bSJed Brown   PetscFunctionReturn(0);
514e27a552bSJed Brown }
515e27a552bSJed Brown 
516e27a552bSJed Brown /*------------------------------------------------------------*/
517e27a552bSJed Brown #undef __FUNCT__
518e27a552bSJed Brown #define __FUNCT__ "TSReset_RosW"
519e27a552bSJed Brown static PetscErrorCode TSReset_RosW(TS ts)
520e27a552bSJed Brown {
52161692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
522e27a552bSJed Brown   PetscInt       s;
523e27a552bSJed Brown   PetscErrorCode ierr;
524e27a552bSJed Brown 
525e27a552bSJed Brown   PetscFunctionBegin;
52661692a83SJed Brown   if (!ros->tableau) PetscFunctionReturn(0);
52761692a83SJed Brown    s = ros->tableau->s;
52861692a83SJed Brown   ierr = VecDestroyVecs(s,&ros->Y);CHKERRQ(ierr);
52961692a83SJed Brown   ierr = VecDestroy(&ros->Ydot);CHKERRQ(ierr);
53061692a83SJed Brown   ierr = VecDestroy(&ros->Ystage);CHKERRQ(ierr);
53161692a83SJed Brown   ierr = VecDestroy(&ros->Zdot);CHKERRQ(ierr);
53261692a83SJed Brown   ierr = VecDestroy(&ros->Zstage);CHKERRQ(ierr);
53361692a83SJed Brown   ierr = PetscFree(ros->work);CHKERRQ(ierr);
534e27a552bSJed Brown   PetscFunctionReturn(0);
535e27a552bSJed Brown }
536e27a552bSJed Brown 
537e27a552bSJed Brown #undef __FUNCT__
538e27a552bSJed Brown #define __FUNCT__ "TSDestroy_RosW"
539e27a552bSJed Brown static PetscErrorCode TSDestroy_RosW(TS ts)
540e27a552bSJed Brown {
541e27a552bSJed Brown   PetscErrorCode  ierr;
542e27a552bSJed Brown 
543e27a552bSJed Brown   PetscFunctionBegin;
544e27a552bSJed Brown   ierr = TSReset_RosW(ts);CHKERRQ(ierr);
545e27a552bSJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
546e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","",PETSC_NULL);CHKERRQ(ierr);
547e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","",PETSC_NULL);CHKERRQ(ierr);
54861692a83SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","",PETSC_NULL);CHKERRQ(ierr);
549e27a552bSJed Brown   PetscFunctionReturn(0);
550e27a552bSJed Brown }
551e27a552bSJed Brown 
552e27a552bSJed Brown /*
553e27a552bSJed Brown   This defines the nonlinear equation that is to be solved with SNES
554e27a552bSJed Brown   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
555e27a552bSJed Brown */
556e27a552bSJed Brown #undef __FUNCT__
557e27a552bSJed Brown #define __FUNCT__ "SNESTSFormFunction_RosW"
558e27a552bSJed Brown static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec X,Vec F,TS ts)
559e27a552bSJed Brown {
56061692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
561e27a552bSJed Brown   PetscErrorCode ierr;
562e27a552bSJed Brown 
563e27a552bSJed Brown   PetscFunctionBegin;
564*c17803e7SJed Brown   if (ros->stage_explicit) {
565*c17803e7SJed Brown     ierr = VecAXPBY(ros->Ydot,ros->shift,0.0,X);CHKERRQ(ierr);       /* Ydot = shift*X*/
566*c17803e7SJed Brown   } else {
56761692a83SJed Brown     ierr = VecWAXPY(ros->Ydot,ros->shift,X,ros->Zdot);CHKERRQ(ierr); /* Ydot = shift*X + Zdot */
568*c17803e7SJed Brown   }
56961692a83SJed Brown   ierr = VecWAXPY(ros->Ystage,1.0,X,ros->Zstage);CHKERRQ(ierr);    /* Ystage = X + Zstage */
57061692a83SJed Brown   ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,ros->Ydot,F,PETSC_FALSE);CHKERRQ(ierr);
571e27a552bSJed Brown   PetscFunctionReturn(0);
572e27a552bSJed Brown }
573e27a552bSJed Brown 
574e27a552bSJed Brown #undef __FUNCT__
575e27a552bSJed Brown #define __FUNCT__ "SNESTSFormJacobian_RosW"
576e27a552bSJed Brown static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
577e27a552bSJed Brown {
57861692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
579e27a552bSJed Brown   PetscErrorCode ierr;
580e27a552bSJed Brown 
581e27a552bSJed Brown   PetscFunctionBegin;
58261692a83SJed Brown   /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
58361692a83SJed Brown   ierr = TSComputeIJacobian(ts,ros->stage_time,ros->Ystage,ros->Ydot,ros->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr);
584e27a552bSJed Brown   PetscFunctionReturn(0);
585e27a552bSJed Brown }
586e27a552bSJed Brown 
587e27a552bSJed Brown #undef __FUNCT__
588e27a552bSJed Brown #define __FUNCT__ "TSSetUp_RosW"
589e27a552bSJed Brown static PetscErrorCode TSSetUp_RosW(TS ts)
590e27a552bSJed Brown {
59161692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
59261692a83SJed Brown   RosWTableau    tab  = ros->tableau;
593e27a552bSJed Brown   PetscInt       s    = tab->s;
594e27a552bSJed Brown   PetscErrorCode ierr;
595e27a552bSJed Brown 
596e27a552bSJed Brown   PetscFunctionBegin;
59761692a83SJed Brown   if (!ros->tableau) {
598e27a552bSJed Brown     ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);
599e27a552bSJed Brown   }
60061692a83SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ros->Y);CHKERRQ(ierr);
60161692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Ydot);CHKERRQ(ierr);
60261692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Ystage);CHKERRQ(ierr);
60361692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Zdot);CHKERRQ(ierr);
60461692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Zstage);CHKERRQ(ierr);
60561692a83SJed Brown   ierr = PetscMalloc(s*sizeof(ros->work[0]),&ros->work);CHKERRQ(ierr);
606e27a552bSJed Brown   PetscFunctionReturn(0);
607e27a552bSJed Brown }
608e27a552bSJed Brown /*------------------------------------------------------------*/
609e27a552bSJed Brown 
610e27a552bSJed Brown #undef __FUNCT__
611e27a552bSJed Brown #define __FUNCT__ "TSSetFromOptions_RosW"
612e27a552bSJed Brown static PetscErrorCode TSSetFromOptions_RosW(TS ts)
613e27a552bSJed Brown {
61461692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
615e27a552bSJed Brown   PetscErrorCode ierr;
61661692a83SJed Brown   char           rostype[256];
617e27a552bSJed Brown 
618e27a552bSJed Brown   PetscFunctionBegin;
619e27a552bSJed Brown   ierr = PetscOptionsHead("RosW ODE solver options");CHKERRQ(ierr);
620e27a552bSJed Brown   {
62161692a83SJed Brown     RosWTableauLink link;
622e27a552bSJed Brown     PetscInt count,choice;
623e27a552bSJed Brown     PetscBool flg;
624e27a552bSJed Brown     const char **namelist;
62561692a83SJed Brown     SNES snes;
62661692a83SJed Brown 
62761692a83SJed Brown     ierr = PetscStrncpy(rostype,TSRosWDefault,sizeof rostype);CHKERRQ(ierr);
62861692a83SJed Brown     for (link=RosWTableauList,count=0; link; link=link->next,count++) ;
629e27a552bSJed Brown     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
63061692a83SJed Brown     for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
63161692a83SJed Brown     ierr = PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,rostype,&choice,&flg);CHKERRQ(ierr);
63261692a83SJed Brown     ierr = TSRosWSetType(ts,flg ? namelist[choice] : rostype);CHKERRQ(ierr);
633e27a552bSJed Brown     ierr = PetscFree(namelist);CHKERRQ(ierr);
63461692a83SJed Brown 
63561692a83SJed Brown     ierr = PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,PETSC_NULL);CHKERRQ(ierr);
63661692a83SJed Brown 
63761692a83SJed Brown     /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
63861692a83SJed Brown     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
63961692a83SJed Brown     if (!((PetscObject)snes)->type_name) {
64061692a83SJed Brown       ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
64161692a83SJed Brown     }
64261692a83SJed Brown     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
643e27a552bSJed Brown   }
644e27a552bSJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
645e27a552bSJed Brown   PetscFunctionReturn(0);
646e27a552bSJed Brown }
647e27a552bSJed Brown 
648e27a552bSJed Brown #undef __FUNCT__
649e27a552bSJed Brown #define __FUNCT__ "PetscFormatRealArray"
650e27a552bSJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
651e27a552bSJed Brown {
652e27a552bSJed Brown   PetscErrorCode ierr;
653e408995aSJed Brown   PetscInt i;
654e408995aSJed Brown   size_t left,count;
655e27a552bSJed Brown   char *p;
656e27a552bSJed Brown 
657e27a552bSJed Brown   PetscFunctionBegin;
658e408995aSJed Brown   for (i=0,p=buf,left=len; i<n; i++) {
659e408995aSJed Brown     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
660e27a552bSJed Brown     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
661e27a552bSJed Brown     left -= count;
662e27a552bSJed Brown     p += count;
663e27a552bSJed Brown     *p++ = ' ';
664e27a552bSJed Brown   }
665e27a552bSJed Brown   p[i ? 0 : -1] = 0;
666e27a552bSJed Brown   PetscFunctionReturn(0);
667e27a552bSJed Brown }
668e27a552bSJed Brown 
669e27a552bSJed Brown #undef __FUNCT__
670e27a552bSJed Brown #define __FUNCT__ "TSView_RosW"
671e27a552bSJed Brown static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer)
672e27a552bSJed Brown {
67361692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
67461692a83SJed Brown   RosWTableau    tab  = ros->tableau;
675e27a552bSJed Brown   PetscBool      iascii;
676e27a552bSJed Brown   PetscErrorCode ierr;
677e27a552bSJed Brown 
678e27a552bSJed Brown   PetscFunctionBegin;
679e27a552bSJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
680e27a552bSJed Brown   if (iascii) {
68161692a83SJed Brown     const TSRosWType rostype;
682e408995aSJed Brown     PetscInt i;
683e408995aSJed Brown     PetscReal abscissa[512];
684e27a552bSJed Brown     char buf[512];
68561692a83SJed Brown     ierr = TSRosWGetType(ts,&rostype);CHKERRQ(ierr);
68661692a83SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Rosenbrock-W %s\n",rostype);CHKERRQ(ierr);
687e408995aSJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ASum);CHKERRQ(ierr);
68861692a83SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A       = %s\n",buf);CHKERRQ(ierr);
689e408995aSJed Brown     for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i];
690e408995aSJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,abscissa);CHKERRQ(ierr);
691e408995aSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A+Gamma = %s\n",buf);CHKERRQ(ierr);
692e27a552bSJed Brown   }
693e27a552bSJed Brown   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
694e27a552bSJed Brown   PetscFunctionReturn(0);
695e27a552bSJed Brown }
696e27a552bSJed Brown 
697e27a552bSJed Brown #undef __FUNCT__
698e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType"
699e27a552bSJed Brown /*@C
70061692a83SJed Brown   TSRosWSetType - Set the type of Rosenbrock-W scheme
701e27a552bSJed Brown 
702e27a552bSJed Brown   Logically collective
703e27a552bSJed Brown 
704e27a552bSJed Brown   Input Parameter:
705e27a552bSJed Brown +  ts - timestepping context
70661692a83SJed Brown -  rostype - type of Rosenbrock-W scheme
707e27a552bSJed Brown 
708e27a552bSJed Brown   Level: intermediate
709e27a552bSJed Brown 
710e27a552bSJed Brown .seealso: TSRosWGetType()
711e27a552bSJed Brown @*/
71261692a83SJed Brown PetscErrorCode TSRosWSetType(TS ts,const TSRosWType rostype)
713e27a552bSJed Brown {
714e27a552bSJed Brown   PetscErrorCode ierr;
715e27a552bSJed Brown 
716e27a552bSJed Brown   PetscFunctionBegin;
717e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
71861692a83SJed Brown   ierr = PetscTryMethod(ts,"TSRosWSetType_C",(TS,const TSRosWType),(ts,rostype));CHKERRQ(ierr);
719e27a552bSJed Brown   PetscFunctionReturn(0);
720e27a552bSJed Brown }
721e27a552bSJed Brown 
722e27a552bSJed Brown #undef __FUNCT__
723e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType"
724e27a552bSJed Brown /*@C
72561692a83SJed Brown   TSRosWGetType - Get the type of Rosenbrock-W scheme
726e27a552bSJed Brown 
727e27a552bSJed Brown   Logically collective
728e27a552bSJed Brown 
729e27a552bSJed Brown   Input Parameter:
730e27a552bSJed Brown .  ts - timestepping context
731e27a552bSJed Brown 
732e27a552bSJed Brown   Output Parameter:
73361692a83SJed Brown .  rostype - type of Rosenbrock-W scheme
734e27a552bSJed Brown 
735e27a552bSJed Brown   Level: intermediate
736e27a552bSJed Brown 
737e27a552bSJed Brown .seealso: TSRosWGetType()
738e27a552bSJed Brown @*/
73961692a83SJed Brown PetscErrorCode TSRosWGetType(TS ts,const TSRosWType *rostype)
740e27a552bSJed Brown {
741e27a552bSJed Brown   PetscErrorCode ierr;
742e27a552bSJed Brown 
743e27a552bSJed Brown   PetscFunctionBegin;
744e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
74561692a83SJed Brown   ierr = PetscUseMethod(ts,"TSRosWGetType_C",(TS,const TSRosWType*),(ts,rostype));CHKERRQ(ierr);
746e27a552bSJed Brown   PetscFunctionReturn(0);
747e27a552bSJed Brown }
748e27a552bSJed Brown 
749e27a552bSJed Brown #undef __FUNCT__
75061692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian"
751e27a552bSJed Brown /*@C
75261692a83SJed Brown   TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step.
753e27a552bSJed Brown 
754e27a552bSJed Brown   Logically collective
755e27a552bSJed Brown 
756e27a552bSJed Brown   Input Parameter:
757e27a552bSJed Brown +  ts - timestepping context
75861692a83SJed Brown -  flg - PETSC_TRUE to recompute the Jacobian at each stage
759e27a552bSJed Brown 
760e27a552bSJed Brown   Level: intermediate
761e27a552bSJed Brown 
762e27a552bSJed Brown .seealso: TSRosWGetType()
763e27a552bSJed Brown @*/
76461692a83SJed Brown PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg)
765e27a552bSJed Brown {
766e27a552bSJed Brown   PetscErrorCode ierr;
767e27a552bSJed Brown 
768e27a552bSJed Brown   PetscFunctionBegin;
769e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
77061692a83SJed Brown   ierr = PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
771e27a552bSJed Brown   PetscFunctionReturn(0);
772e27a552bSJed Brown }
773e27a552bSJed Brown 
774e27a552bSJed Brown EXTERN_C_BEGIN
775e27a552bSJed Brown #undef __FUNCT__
776e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType_RosW"
77761692a83SJed Brown PetscErrorCode  TSRosWGetType_RosW(TS ts,const TSRosWType *rostype)
778e27a552bSJed Brown {
77961692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
780e27a552bSJed Brown   PetscErrorCode ierr;
781e27a552bSJed Brown 
782e27a552bSJed Brown   PetscFunctionBegin;
78361692a83SJed Brown   if (!ros->tableau) {ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);}
78461692a83SJed Brown   *rostype = ros->tableau->name;
785e27a552bSJed Brown   PetscFunctionReturn(0);
786e27a552bSJed Brown }
787e27a552bSJed Brown #undef __FUNCT__
788e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType_RosW"
78961692a83SJed Brown PetscErrorCode  TSRosWSetType_RosW(TS ts,const TSRosWType rostype)
790e27a552bSJed Brown {
79161692a83SJed Brown   TS_RosW         *ros = (TS_RosW*)ts->data;
792e27a552bSJed Brown   PetscErrorCode  ierr;
793e27a552bSJed Brown   PetscBool       match;
79461692a83SJed Brown   RosWTableauLink link;
795e27a552bSJed Brown 
796e27a552bSJed Brown   PetscFunctionBegin;
79761692a83SJed Brown   if (ros->tableau) {
79861692a83SJed Brown     ierr = PetscStrcmp(ros->tableau->name,rostype,&match);CHKERRQ(ierr);
799e27a552bSJed Brown     if (match) PetscFunctionReturn(0);
800e27a552bSJed Brown   }
80161692a83SJed Brown   for (link = RosWTableauList; link; link=link->next) {
80261692a83SJed Brown     ierr = PetscStrcmp(link->tab.name,rostype,&match);CHKERRQ(ierr);
803e27a552bSJed Brown     if (match) {
804e27a552bSJed Brown       ierr = TSReset_RosW(ts);CHKERRQ(ierr);
80561692a83SJed Brown       ros->tableau = &link->tab;
806e27a552bSJed Brown       PetscFunctionReturn(0);
807e27a552bSJed Brown     }
808e27a552bSJed Brown   }
80961692a83SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype);
810e27a552bSJed Brown   PetscFunctionReturn(0);
811e27a552bSJed Brown }
81261692a83SJed Brown 
813e27a552bSJed Brown #undef __FUNCT__
81461692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian_RosW"
81561692a83SJed Brown PetscErrorCode  TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg)
816e27a552bSJed Brown {
81761692a83SJed Brown   TS_RosW *ros = (TS_RosW*)ts->data;
818e27a552bSJed Brown 
819e27a552bSJed Brown   PetscFunctionBegin;
82061692a83SJed Brown   ros->recompute_jacobian = flg;
821e27a552bSJed Brown   PetscFunctionReturn(0);
822e27a552bSJed Brown }
823e27a552bSJed Brown EXTERN_C_END
824e27a552bSJed Brown 
825e27a552bSJed Brown /* ------------------------------------------------------------ */
826e27a552bSJed Brown /*MC
827e27a552bSJed Brown       TSRosW - ODE solver using Rosenbrock-W schemes
828e27a552bSJed Brown 
829e27a552bSJed Brown   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
830e27a552bSJed Brown   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
831e27a552bSJed Brown   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
832e27a552bSJed Brown 
833e27a552bSJed Brown   Notes:
83461692a83SJed Brown   This method currently only works with autonomous ODE and DAE.
83561692a83SJed Brown 
83661692a83SJed Brown   Developer notes:
83761692a83SJed Brown   Rosenbrock-W methods are typically specified for autonomous ODE
83861692a83SJed Brown 
83961692a83SJed Brown $  xdot = f(x)
84061692a83SJed Brown 
84161692a83SJed Brown   by the stage equations
84261692a83SJed Brown 
84361692a83SJed Brown $  k_i = h f(x_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j
84461692a83SJed Brown 
84561692a83SJed Brown   and step completion formula
84661692a83SJed Brown 
84761692a83SJed Brown $  x_1 = x_0 + sum_j b_j k_j
84861692a83SJed Brown 
84961692a83SJed Brown   with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(x)
85061692a83SJed Brown   and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner,
85161692a83SJed Brown   we define new variables for the stage equations
85261692a83SJed Brown 
85361692a83SJed Brown $  y_i = gamma_ij k_j
85461692a83SJed Brown 
85561692a83SJed Brown   The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define
85661692a83SJed Brown 
85761692a83SJed Brown $  A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-i}
85861692a83SJed Brown 
85961692a83SJed Brown   to rewrite the method as
86061692a83SJed Brown 
86161692a83SJed 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
86261692a83SJed Brown $  x_1 = x_0 + sum_j bt_j y_j
86361692a83SJed Brown 
86461692a83SJed Brown    where we have introduced the mass matrix M. Continue by defining
86561692a83SJed Brown 
86661692a83SJed Brown $  ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j
86761692a83SJed Brown 
86861692a83SJed Brown    or, more compactly in tensor notation
86961692a83SJed Brown 
87061692a83SJed Brown $  Ydot = 1/h (Gamma^{-1} \otimes I) Y .
87161692a83SJed Brown 
87261692a83SJed Brown    Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current
87361692a83SJed Brown    stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the
87461692a83SJed Brown    equation
87561692a83SJed Brown 
87661692a83SJed Brown $  g(x_0 + sum_j a_ij y_j + y_i, ydot_i) = 0
87761692a83SJed Brown 
87861692a83SJed Brown    with initial guess y_i = 0.
879e27a552bSJed Brown 
880e27a552bSJed Brown   Level: beginner
881e27a552bSJed Brown 
882e27a552bSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSRosWRegister()
883e27a552bSJed Brown 
884e27a552bSJed Brown M*/
885e27a552bSJed Brown EXTERN_C_BEGIN
886e27a552bSJed Brown #undef __FUNCT__
887e27a552bSJed Brown #define __FUNCT__ "TSCreate_RosW"
888e27a552bSJed Brown PetscErrorCode  TSCreate_RosW(TS ts)
889e27a552bSJed Brown {
89061692a83SJed Brown   TS_RosW        *ros;
891e27a552bSJed Brown   PetscErrorCode ierr;
892e27a552bSJed Brown 
893e27a552bSJed Brown   PetscFunctionBegin;
894e27a552bSJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
895e27a552bSJed Brown   ierr = TSRosWInitializePackage(PETSC_NULL);CHKERRQ(ierr);
896e27a552bSJed Brown #endif
897e27a552bSJed Brown 
898e27a552bSJed Brown   ts->ops->reset          = TSReset_RosW;
899e27a552bSJed Brown   ts->ops->destroy        = TSDestroy_RosW;
900e27a552bSJed Brown   ts->ops->view           = TSView_RosW;
901e27a552bSJed Brown   ts->ops->setup          = TSSetUp_RosW;
902e27a552bSJed Brown   ts->ops->step           = TSStep_RosW;
903e27a552bSJed Brown   ts->ops->interpolate    = TSInterpolate_RosW;
9041c3436cfSJed Brown   ts->ops->evaluatestep   = TSEvaluateStep_RosW;
905e27a552bSJed Brown   ts->ops->setfromoptions = TSSetFromOptions_RosW;
906e27a552bSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_RosW;
907e27a552bSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_RosW;
908e27a552bSJed Brown 
90961692a83SJed Brown   ierr = PetscNewLog(ts,TS_RosW,&ros);CHKERRQ(ierr);
91061692a83SJed Brown   ts->data = (void*)ros;
911e27a552bSJed Brown 
912e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","TSRosWGetType_RosW",TSRosWGetType_RosW);CHKERRQ(ierr);
913e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","TSRosWSetType_RosW",TSRosWSetType_RosW);CHKERRQ(ierr);
91461692a83SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","TSRosWSetRecomputeJacobian_RosW",TSRosWSetRecomputeJacobian_RosW);CHKERRQ(ierr);
915e27a552bSJed Brown   PetscFunctionReturn(0);
916e27a552bSJed Brown }
917e27a552bSJed Brown EXTERN_C_END
918