xref: /petsc/src/ts/impls/rosw/rosw.c (revision 1c3436cf04c11f6e65a053bdffa658cacfbd163a)
1e27a552bSJed Brown /*
261692a83SJed Brown   Code for timestepping with Rosenbrock W methods
3e27a552bSJed Brown 
4e27a552bSJed Brown   Notes:
5e27a552bSJed Brown   The general system is written as
6e27a552bSJed Brown 
761692a83SJed Brown   G(t,X,Xdot) = F(t,X)
8e27a552bSJed Brown 
961692a83SJed Brown   where G represents the stiff part of the physics and F represents the non-stiff part.
1061692a83SJed Brown   This method is designed to be linearly implicit on G and can use an approximate and lagged Jacobian.
11e27a552bSJed Brown 
12e27a552bSJed Brown */
13e27a552bSJed Brown #include <private/tsimpl.h>                /*I   "petscts.h"   I*/
14e27a552bSJed Brown 
1561692a83SJed Brown #include <../src/mat/blockinvert.h>
1661692a83SJed Brown 
1761692a83SJed Brown static const TSRosWType TSRosWDefault = TSROSW2P;
18e27a552bSJed Brown static PetscBool TSRosWRegisterAllCalled;
19e27a552bSJed Brown static PetscBool TSRosWPackageInitialized;
20e27a552bSJed Brown 
2161692a83SJed Brown typedef struct _RosWTableau *RosWTableau;
2261692a83SJed Brown struct _RosWTableau {
23e27a552bSJed Brown   char *name;
24e27a552bSJed Brown   PetscInt order;               /* Classical approximation order of the method */
25e27a552bSJed Brown   PetscInt s;                   /* Number of stages */
2661692a83SJed Brown   PetscReal *A;                 /* Propagation table, strictly lower triangular */
2761692a83SJed Brown   PetscReal *Gamma;             /* Stage table, lower triangular with nonzero diagonal */
2861692a83SJed Brown   PetscReal *b;                 /* Step completion table */
29fe7e6d57SJed Brown   PetscReal *bembed;            /* Step completion table for embedded method of order one less */
3061692a83SJed Brown   PetscReal *ASum;              /* Row sum of A */
3161692a83SJed Brown   PetscReal *GammaSum;          /* Row sum of Gamma, only needed for non-autonomous systems */
3261692a83SJed Brown   PetscReal *At;                /* Propagation table in transformed variables */
3361692a83SJed Brown   PetscReal *bt;                /* Step completion table in transformed variables */
34fe7e6d57SJed Brown   PetscReal *bembedt;           /* Step completion table of order one less in transformed variables */
3561692a83SJed Brown   PetscReal *GammaInv;          /* Inverse of Gamma, used for transformed variables */
36e27a552bSJed Brown };
3761692a83SJed Brown typedef struct _RosWTableauLink *RosWTableauLink;
3861692a83SJed Brown struct _RosWTableauLink {
3961692a83SJed Brown   struct _RosWTableau tab;
4061692a83SJed Brown   RosWTableauLink next;
41e27a552bSJed Brown };
4261692a83SJed Brown static RosWTableauLink RosWTableauList;
43e27a552bSJed Brown 
44e27a552bSJed Brown typedef struct {
4561692a83SJed Brown   RosWTableau tableau;
4661692a83SJed Brown   Vec         *Y;               /* States computed during the step, used to complete the step */
47e27a552bSJed Brown   Vec         Ydot;             /* Work vector holding Ydot during residual evaluation */
4861692a83SJed Brown   Vec         Ystage;           /* Work vector for the state value at each stage */
4961692a83SJed Brown   Vec         Zdot;             /* Ydot = Zdot + shift*Y */
5061692a83SJed Brown   Vec         Zstage;           /* Y = Zstage + Y */
51*1c3436cfSJed Brown   PetscScalar *work;            /* Scalar work space of length number of stages, used to prepare VecMAXPY() */
52e27a552bSJed Brown   PetscReal   shift;
53e27a552bSJed Brown   PetscReal   stage_time;
5461692a83SJed Brown   PetscBool   recompute_jacobian; /* Recompute the Jacobian at each stage, default is to freeze the Jacobian at the start of each step */
55*1c3436cfSJed Brown   PetscBool   step_taken;         /* ts->vec_sol has been advanced to the end of the current time step */
56e27a552bSJed Brown } TS_RosW;
57e27a552bSJed Brown 
58fe7e6d57SJed Brown /*MC
59fe7e6d57SJed Brown      TSROSW2M - Two stage second order L-stable Rosenbrock-W scheme.
60fe7e6d57SJed Brown 
61fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2P.
62fe7e6d57SJed Brown 
63fe7e6d57SJed Brown      Level: intermediate
64fe7e6d57SJed Brown 
65fe7e6d57SJed Brown .seealso: TSROSW
66fe7e6d57SJed Brown M*/
67fe7e6d57SJed Brown 
68fe7e6d57SJed Brown /*MC
69fe7e6d57SJed Brown      TSROSW2P - Two stage second order L-stable Rosenbrock-W scheme.
70fe7e6d57SJed Brown 
71fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2M.
72fe7e6d57SJed Brown 
73fe7e6d57SJed Brown      Level: intermediate
74fe7e6d57SJed Brown 
75fe7e6d57SJed Brown .seealso: TSROSW
76fe7e6d57SJed Brown M*/
77fe7e6d57SJed Brown 
78fe7e6d57SJed Brown /*MC
79fe7e6d57SJed Brown      TSROSWRA3PW - Three stage third order Rosenbrock-W scheme for PDAE of index 1.
80fe7e6d57SJed Brown 
81fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.
82fe7e6d57SJed Brown 
83fe7e6d57SJed 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.
84fe7e6d57SJed Brown 
85fe7e6d57SJed Brown      References:
86fe7e6d57SJed Brown      Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005.
87fe7e6d57SJed Brown 
88fe7e6d57SJed Brown      Level: intermediate
89fe7e6d57SJed Brown 
90fe7e6d57SJed Brown .seealso: TSROSW
91fe7e6d57SJed Brown M*/
92fe7e6d57SJed Brown 
93fe7e6d57SJed Brown /*MC
94fe7e6d57SJed Brown      TSROSWRA34PW2 - Four stage third order L-stable Rosenbrock-W scheme for PDAE of index 1.
95fe7e6d57SJed Brown 
96fe7e6d57SJed Brown      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.
97fe7e6d57SJed Brown 
98fe7e6d57SJed 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.
99fe7e6d57SJed Brown 
100fe7e6d57SJed Brown      References:
101fe7e6d57SJed Brown      Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005.
102fe7e6d57SJed Brown 
103fe7e6d57SJed Brown      Level: intermediate
104fe7e6d57SJed Brown 
105fe7e6d57SJed Brown .seealso: TSROSW
106fe7e6d57SJed Brown M*/
107fe7e6d57SJed Brown 
108e27a552bSJed Brown #undef __FUNCT__
109e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterAll"
110e27a552bSJed Brown /*@C
111e27a552bSJed Brown   TSRosWRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSRosW
112e27a552bSJed Brown 
113e27a552bSJed Brown   Not Collective, but should be called by all processes which will need the schemes to be registered
114e27a552bSJed Brown 
115e27a552bSJed Brown   Level: advanced
116e27a552bSJed Brown 
117e27a552bSJed Brown .keywords: TS, TSRosW, register, all
118e27a552bSJed Brown 
119e27a552bSJed Brown .seealso:  TSRosWRegisterDestroy()
120e27a552bSJed Brown @*/
121e27a552bSJed Brown PetscErrorCode TSRosWRegisterAll(void)
122e27a552bSJed Brown {
123e27a552bSJed Brown   PetscErrorCode ierr;
124e27a552bSJed Brown 
125e27a552bSJed Brown   PetscFunctionBegin;
126e27a552bSJed Brown   if (TSRosWRegisterAllCalled) PetscFunctionReturn(0);
127e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_TRUE;
128e27a552bSJed Brown 
129e27a552bSJed Brown   {
13061692a83SJed Brown     const PetscReal g = 1. + 1./PetscSqrtReal(2.0);
131e27a552bSJed Brown     const PetscReal
13261692a83SJed Brown       A[2][2] = {{0,0}, {1.,0}},
13361692a83SJed Brown       Gamma[2][2] = {{g,0}, {-2.*g,g}},
134*1c3436cfSJed Brown       b[2] = {0.5,0.5},
135*1c3436cfSJed Brown       b1[2] = {1.0,0.0};
136*1c3436cfSJed Brown     ierr = TSRosWRegister(TSROSW2P,2,2,&A[0][0],&Gamma[0][0],b,b1);CHKERRQ(ierr);
137e27a552bSJed Brown   }
138e27a552bSJed Brown   {
13961692a83SJed Brown     const PetscReal g = 1. - 1./PetscSqrtReal(2.0);
140e27a552bSJed Brown     const PetscReal
14161692a83SJed Brown       A[2][2] = {{0,0}, {1.,0}},
14261692a83SJed Brown       Gamma[2][2] = {{g,0}, {-2.*g,g}},
143*1c3436cfSJed Brown       b[2] = {0.5,0.5},
144*1c3436cfSJed Brown       b1[2] = {1.0,0.0};
145*1c3436cfSJed Brown     ierr = TSRosWRegister(TSROSW2M,2,2,&A[0][0],&Gamma[0][0],b,b1);CHKERRQ(ierr);
146fe7e6d57SJed Brown   }
147fe7e6d57SJed Brown   {
148fe7e6d57SJed Brown     const PetscReal g = 7.8867513459481287e-01;
149fe7e6d57SJed Brown     const PetscReal
150fe7e6d57SJed Brown       A[3][3] = {{0,0,0},
151fe7e6d57SJed Brown                  {1.5773502691896257e+00,0,0},
152fe7e6d57SJed Brown                  {0.5,0,0}},
153fe7e6d57SJed Brown       Gamma[3][3] = {{g,0,0},
154fe7e6d57SJed Brown                      {-1.5773502691896257e+00,g,0},
155fe7e6d57SJed Brown                      {-6.7075317547305480e-01,1.7075317547305482e-01,g}},
156fe7e6d57SJed Brown       b[3] = {1.0566243270259355e-01,4.9038105676657971e-02,8.4529946162074843e-01},
157fe7e6d57SJed Brown       b2[3] = {-1.7863279495408180e-01,1./3.,8.4529946162074843e-01};
158fe7e6d57SJed Brown     ierr = TSRosWRegister(TSROSWRA3PW,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
159fe7e6d57SJed Brown   }
160fe7e6d57SJed Brown   {
161fe7e6d57SJed Brown     const PetscReal g = 4.3586652150845900e-01;
162fe7e6d57SJed Brown     const PetscReal
163fe7e6d57SJed Brown       A[4][4] = {{0,0,0,0},
164fe7e6d57SJed Brown                  {8.7173304301691801e-01,0,0,0},
165fe7e6d57SJed Brown                  {8.4457060015369423e-01,-1.1299064236484185e-01,0,0},
166fe7e6d57SJed Brown                  {0,0,1.,0}},
167fe7e6d57SJed Brown       Gamma[4][4] = {{g,0,0,0},
168fe7e6d57SJed Brown                      {-8.7173304301691801e-01,g,0,0},
169fe7e6d57SJed Brown                      {-9.0338057013044082e-01,5.4180672388095326e-02,g,0},
170fe7e6d57SJed Brown                      {2.4212380706095346e-01,-1.2232505839045147e+00,5.4526025533510214e-01,g}},
171fe7e6d57SJed Brown       b[4] = {2.4212380706095346e-01,-1.2232505839045147e+00,1.5452602553351020e+00,4.3586652150845900e-01},
172fe7e6d57SJed Brown       b2[4] = {3.7810903145819369e-01,-9.6042292212423178e-02,5.0000000000000000e-01,2.1793326075422950e-01};
173fe7e6d57SJed Brown     ierr = TSRosWRegister(TSROSWRA34PW2,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
174e27a552bSJed Brown   }
175e27a552bSJed Brown   PetscFunctionReturn(0);
176e27a552bSJed Brown }
177e27a552bSJed Brown 
178e27a552bSJed Brown #undef __FUNCT__
179e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterDestroy"
180e27a552bSJed Brown /*@C
181e27a552bSJed Brown    TSRosWRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister().
182e27a552bSJed Brown 
183e27a552bSJed Brown    Not Collective
184e27a552bSJed Brown 
185e27a552bSJed Brown    Level: advanced
186e27a552bSJed Brown 
187e27a552bSJed Brown .keywords: TSRosW, register, destroy
188e27a552bSJed Brown .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegisterDynamic()
189e27a552bSJed Brown @*/
190e27a552bSJed Brown PetscErrorCode TSRosWRegisterDestroy(void)
191e27a552bSJed Brown {
192e27a552bSJed Brown   PetscErrorCode ierr;
19361692a83SJed Brown   RosWTableauLink link;
194e27a552bSJed Brown 
195e27a552bSJed Brown   PetscFunctionBegin;
19661692a83SJed Brown   while ((link = RosWTableauList)) {
19761692a83SJed Brown     RosWTableau t = &link->tab;
19861692a83SJed Brown     RosWTableauList = link->next;
19961692a83SJed Brown     ierr = PetscFree5(t->A,t->Gamma,t->b,t->ASum,t->GammaSum);CHKERRQ(ierr);
20061692a83SJed Brown     ierr = PetscFree3(t->At,t->bt,t->GammaInv);CHKERRQ(ierr);
201fe7e6d57SJed Brown     ierr = PetscFree2(t->bembed,t->bembedt);CHKERRQ(ierr);
202e27a552bSJed Brown     ierr = PetscFree(t->name);CHKERRQ(ierr);
203e27a552bSJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
204e27a552bSJed Brown   }
205e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_FALSE;
206e27a552bSJed Brown   PetscFunctionReturn(0);
207e27a552bSJed Brown }
208e27a552bSJed Brown 
209e27a552bSJed Brown #undef __FUNCT__
210e27a552bSJed Brown #define __FUNCT__ "TSRosWInitializePackage"
211e27a552bSJed Brown /*@C
212e27a552bSJed Brown   TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called
213e27a552bSJed Brown   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RosW()
214e27a552bSJed Brown   when using static libraries.
215e27a552bSJed Brown 
216e27a552bSJed Brown   Input Parameter:
217e27a552bSJed Brown   path - The dynamic library path, or PETSC_NULL
218e27a552bSJed Brown 
219e27a552bSJed Brown   Level: developer
220e27a552bSJed Brown 
221e27a552bSJed Brown .keywords: TS, TSRosW, initialize, package
222e27a552bSJed Brown .seealso: PetscInitialize()
223e27a552bSJed Brown @*/
224e27a552bSJed Brown PetscErrorCode TSRosWInitializePackage(const char path[])
225e27a552bSJed Brown {
226e27a552bSJed Brown   PetscErrorCode ierr;
227e27a552bSJed Brown 
228e27a552bSJed Brown   PetscFunctionBegin;
229e27a552bSJed Brown   if (TSRosWPackageInitialized) PetscFunctionReturn(0);
230e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_TRUE;
231e27a552bSJed Brown   ierr = TSRosWRegisterAll();CHKERRQ(ierr);
232e27a552bSJed Brown   ierr = PetscRegisterFinalize(TSRosWFinalizePackage);CHKERRQ(ierr);
233e27a552bSJed Brown   PetscFunctionReturn(0);
234e27a552bSJed Brown }
235e27a552bSJed Brown 
236e27a552bSJed Brown #undef __FUNCT__
237e27a552bSJed Brown #define __FUNCT__ "TSRosWFinalizePackage"
238e27a552bSJed Brown /*@C
239e27a552bSJed Brown   TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is
240e27a552bSJed Brown   called from PetscFinalize().
241e27a552bSJed Brown 
242e27a552bSJed Brown   Level: developer
243e27a552bSJed Brown 
244e27a552bSJed Brown .keywords: Petsc, destroy, package
245e27a552bSJed Brown .seealso: PetscFinalize()
246e27a552bSJed Brown @*/
247e27a552bSJed Brown PetscErrorCode TSRosWFinalizePackage(void)
248e27a552bSJed Brown {
249e27a552bSJed Brown   PetscErrorCode ierr;
250e27a552bSJed Brown 
251e27a552bSJed Brown   PetscFunctionBegin;
252e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_FALSE;
253e27a552bSJed Brown   ierr = TSRosWRegisterDestroy();CHKERRQ(ierr);
254e27a552bSJed Brown   PetscFunctionReturn(0);
255e27a552bSJed Brown }
256e27a552bSJed Brown 
257e27a552bSJed Brown #undef __FUNCT__
258e27a552bSJed Brown #define __FUNCT__ "TSRosWRegister"
259e27a552bSJed Brown /*@C
26061692a83SJed Brown    TSRosWRegister - register a Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
261e27a552bSJed Brown 
262e27a552bSJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
263e27a552bSJed Brown 
264e27a552bSJed Brown    Input Parameters:
265e27a552bSJed Brown +  name - identifier for method
266e27a552bSJed Brown .  order - approximation order of method
267e27a552bSJed Brown .  s - number of stages, this is the dimension of the matrices below
26861692a83SJed Brown .  A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular
26961692a83SJed Brown .  Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal
270fe7e6d57SJed Brown .  b - Step completion table (dimension s)
271fe7e6d57SJed Brown -  bembed - Step completion table for a scheme of order one less (dimension s, PETSC_NULL if no embedded scheme is available)
272e27a552bSJed Brown 
273e27a552bSJed Brown    Notes:
27461692a83SJed Brown    Several Rosenbrock W methods are provided, this function is only needed to create new methods.
275e27a552bSJed Brown 
276e27a552bSJed Brown    Level: advanced
277e27a552bSJed Brown 
278e27a552bSJed Brown .keywords: TS, register
279e27a552bSJed Brown 
280e27a552bSJed Brown .seealso: TSRosW
281e27a552bSJed Brown @*/
282e27a552bSJed Brown PetscErrorCode TSRosWRegister(const TSRosWType name,PetscInt order,PetscInt s,
283fe7e6d57SJed Brown                               const PetscReal A[],const PetscReal Gamma[],const PetscReal b[],const PetscReal bembed[])
284e27a552bSJed Brown {
285e27a552bSJed Brown   PetscErrorCode ierr;
28661692a83SJed Brown   RosWTableauLink link;
28761692a83SJed Brown   RosWTableau t;
28861692a83SJed Brown   PetscInt i,j,k;
28961692a83SJed Brown   PetscScalar *GammaInv;
290e27a552bSJed Brown 
291e27a552bSJed Brown   PetscFunctionBegin;
292fe7e6d57SJed Brown   PetscValidCharPointer(name,1);
293fe7e6d57SJed Brown   PetscValidPointer(A,4);
294fe7e6d57SJed Brown   PetscValidPointer(Gamma,5);
295fe7e6d57SJed Brown   PetscValidPointer(b,6);
296fe7e6d57SJed Brown   if (bembed) PetscValidPointer(bembed,7);
297fe7e6d57SJed Brown 
298e27a552bSJed Brown   ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
299e27a552bSJed Brown   ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
300e27a552bSJed Brown   t = &link->tab;
301e27a552bSJed Brown   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
302e27a552bSJed Brown   t->order = order;
303e27a552bSJed Brown   t->s = s;
30461692a83SJed 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);
30561692a83SJed Brown   ierr = PetscMalloc3(s*s,PetscReal,&t->At,s,PetscReal,&t->bt,s*s,PetscReal,&t->GammaInv);CHKERRQ(ierr);
306e27a552bSJed Brown   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
30761692a83SJed Brown   ierr = PetscMemcpy(t->Gamma,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr);
30861692a83SJed Brown   ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
309fe7e6d57SJed Brown   if (bembed) {
310fe7e6d57SJed Brown     ierr = PetscMalloc2(s,PetscReal,&t->bembed,s,PetscReal,&t->bembedt);CHKERRQ(ierr);
311fe7e6d57SJed Brown     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
312fe7e6d57SJed Brown   }
31361692a83SJed Brown   for (i=0; i<s; i++) {
31461692a83SJed Brown     t->ASum[i] = 0;
31561692a83SJed Brown     t->GammaSum[i] = 0;
31661692a83SJed Brown     for (j=0; j<s; j++) {
31761692a83SJed Brown       t->ASum[i] += A[i*s+j];
318fe7e6d57SJed Brown       t->GammaSum[i] += Gamma[i*s+j];
31961692a83SJed Brown     }
32061692a83SJed Brown   }
32161692a83SJed Brown   ierr = PetscMalloc(s*s*sizeof(PetscScalar),&GammaInv);CHKERRQ(ierr); /* Need to use Scalar for inverse, then convert back to Real */
32261692a83SJed Brown   for (i=0; i<s*s; i++) GammaInv[i] = Gamma[i];
32361692a83SJed Brown   switch (s) {
32461692a83SJed Brown   case 1: GammaInv[0] = 1./GammaInv[0]; break;
32561692a83SJed Brown   case 2: ierr = Kernel_A_gets_inverse_A_2(GammaInv,0);CHKERRQ(ierr); break;
32661692a83SJed Brown   case 3: ierr = Kernel_A_gets_inverse_A_3(GammaInv,0);CHKERRQ(ierr); break;
32761692a83SJed Brown   case 4: ierr = Kernel_A_gets_inverse_A_4(GammaInv,0);CHKERRQ(ierr); break;
32861692a83SJed Brown   case 5: {
32961692a83SJed Brown     PetscInt ipvt5[5];
33061692a83SJed Brown     MatScalar work5[5*5];
33161692a83SJed Brown     ierr = Kernel_A_gets_inverse_A_5(GammaInv,ipvt5,work5,0);CHKERRQ(ierr); break;
33261692a83SJed Brown   }
33361692a83SJed Brown   case 6: ierr = Kernel_A_gets_inverse_A_6(GammaInv,0);CHKERRQ(ierr); break;
33461692a83SJed Brown   case 7: ierr = Kernel_A_gets_inverse_A_7(GammaInv,0);CHKERRQ(ierr); break;
33561692a83SJed Brown   default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for %D stages",s);
33661692a83SJed Brown   }
33761692a83SJed Brown   for (i=0; i<s*s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]);
33861692a83SJed Brown   ierr = PetscFree(GammaInv);CHKERRQ(ierr);
33961692a83SJed Brown   for (i=0; i<s; i++) {
34061692a83SJed Brown     for (j=0; j<s; j++) {
34161692a83SJed Brown       t->At[i*s+j] = 0;
34261692a83SJed Brown       for (k=0; k<s; k++) {
34361692a83SJed Brown         t->At[i*s+j] += t->A[i*s+k] * t->GammaInv[k*s+j];
34461692a83SJed Brown       }
34561692a83SJed Brown     }
34661692a83SJed Brown     t->bt[i] = 0;
34761692a83SJed Brown     for (j=0; j<s; j++) {
34861692a83SJed Brown       t->bt[i] += t->b[j] * t->GammaInv[j*s+i];
34961692a83SJed Brown     }
350fe7e6d57SJed Brown     if (bembed) {
351fe7e6d57SJed Brown       t->bembedt[i] = 0;
352fe7e6d57SJed Brown       for (j=0; j<s; j++) {
353fe7e6d57SJed Brown         t->bembedt[i] += t->bembed[j] * t->GammaInv[j*s+i];
354fe7e6d57SJed Brown       }
355fe7e6d57SJed Brown     }
35661692a83SJed Brown   }
35761692a83SJed Brown   link->next = RosWTableauList;
35861692a83SJed Brown   RosWTableauList = link;
359e27a552bSJed Brown   PetscFunctionReturn(0);
360e27a552bSJed Brown }
361e27a552bSJed Brown 
362e27a552bSJed Brown #undef __FUNCT__
363*1c3436cfSJed Brown #define __FUNCT__ "TSEvaluateStep_RosW"
364*1c3436cfSJed Brown /*
365*1c3436cfSJed Brown  The step completion formula is
366*1c3436cfSJed Brown 
367*1c3436cfSJed Brown  x1 = x0 + b^T Y
368*1c3436cfSJed Brown 
369*1c3436cfSJed Brown  where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been
370*1c3436cfSJed Brown  updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write
371*1c3436cfSJed Brown 
372*1c3436cfSJed Brown  x1e = x0 + be^T Y
373*1c3436cfSJed Brown      = x1 - b^T Y + be^T Y
374*1c3436cfSJed Brown      = x1 + (be - b)^T Y
375*1c3436cfSJed Brown 
376*1c3436cfSJed Brown  so we can evaluate the method of different order even after the step has been optimistically completed.
377*1c3436cfSJed Brown */
378*1c3436cfSJed Brown static PetscErrorCode TSEvaluateStep_RosW(TS ts,PetscInt order,Vec X,PetscBool *done)
379*1c3436cfSJed Brown {
380*1c3436cfSJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
381*1c3436cfSJed Brown   RosWTableau    tab  = ros->tableau;
382*1c3436cfSJed Brown   PetscScalar    *w = ros->work;
383*1c3436cfSJed Brown   PetscInt       i;
384*1c3436cfSJed Brown   PetscErrorCode ierr;
385*1c3436cfSJed Brown 
386*1c3436cfSJed Brown   PetscFunctionBegin;
387*1c3436cfSJed Brown   if (order == tab->order) {
388*1c3436cfSJed Brown     if (ros->step_taken) {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
389*1c3436cfSJed Brown     else {
390*1c3436cfSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
391*1c3436cfSJed Brown       ierr = VecMAXPY(X,tab->s,tab->bt,ros->Y);CHKERRQ(ierr);
392*1c3436cfSJed Brown     }
393*1c3436cfSJed Brown     if (done) *done = PETSC_TRUE;
394*1c3436cfSJed Brown     PetscFunctionReturn(0);
395*1c3436cfSJed Brown   } else if (order == tab->order-1) {
396*1c3436cfSJed Brown     if (!tab->bembedt) goto unavailable;
397*1c3436cfSJed Brown     if (ros->step_taken) {
398*1c3436cfSJed Brown       for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i];
399*1c3436cfSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
400*1c3436cfSJed Brown       ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr);
401*1c3436cfSJed Brown     } else {
402*1c3436cfSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
403*1c3436cfSJed Brown       ierr = VecMAXPY(X,tab->s,tab->bembedt,ros->Y);CHKERRQ(ierr);
404*1c3436cfSJed Brown     }
405*1c3436cfSJed Brown     if (done) *done = PETSC_TRUE;
406*1c3436cfSJed Brown     PetscFunctionReturn(0);
407*1c3436cfSJed Brown   }
408*1c3436cfSJed Brown   unavailable:
409*1c3436cfSJed Brown   if (done) *done = PETSC_FALSE;
410*1c3436cfSJed 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);
411*1c3436cfSJed Brown   PetscFunctionReturn(0);
412*1c3436cfSJed Brown }
413*1c3436cfSJed Brown 
414*1c3436cfSJed Brown #undef __FUNCT__
415e27a552bSJed Brown #define __FUNCT__ "TSStep_RosW"
416e27a552bSJed Brown static PetscErrorCode TSStep_RosW(TS ts)
417e27a552bSJed Brown {
41861692a83SJed Brown   TS_RosW         *ros = (TS_RosW*)ts->data;
41961692a83SJed Brown   RosWTableau     tab  = ros->tableau;
420e27a552bSJed Brown   const PetscInt  s    = tab->s;
421*1c3436cfSJed Brown   const PetscReal *At  = tab->At,*Gamma = tab->Gamma,*ASum = tab->ASum,*GammaInv = tab->GammaInv;
42261692a83SJed Brown   PetscScalar     *w   = ros->work;
42361692a83SJed Brown   Vec             *Y   = ros->Y,Zdot = ros->Zdot,Zstage = ros->Zstage;
424e27a552bSJed Brown   SNES            snes;
425*1c3436cfSJed Brown   TSAdapt         adapt;
426*1c3436cfSJed Brown   PetscInt        i,j,its,lits,reject,next_scheme;
427cdbf8f93SLisandro Dalcin   PetscReal       next_time_step;
428*1c3436cfSJed Brown   PetscBool       accept;
429e27a552bSJed Brown   PetscErrorCode  ierr;
430e27a552bSJed Brown 
431e27a552bSJed Brown   PetscFunctionBegin;
432e27a552bSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
433cdbf8f93SLisandro Dalcin   next_time_step = ts->time_step;
434*1c3436cfSJed Brown   accept = PETSC_TRUE;
435*1c3436cfSJed Brown   ros->step_taken = PETSC_FALSE;
436e27a552bSJed Brown 
437*1c3436cfSJed Brown   for (reject=0; reject<ts->max_reject; reject++,ts->reject++) {
438*1c3436cfSJed Brown     const PetscReal h = ts->time_step;
439e27a552bSJed Brown     for (i=0; i<s; i++) {
440*1c3436cfSJed Brown       ros->stage_time = ts->ptime + h*ASum[i];
44161692a83SJed Brown       ros->shift = 1./(h*Gamma[i*s+i]);
44261692a83SJed Brown 
44361692a83SJed Brown       ierr = VecCopy(ts->vec_sol,Zstage);CHKERRQ(ierr);
44461692a83SJed Brown       ierr = VecMAXPY(Zstage,i,&At[i*s+0],Y);CHKERRQ(ierr);
44561692a83SJed Brown 
44661692a83SJed Brown       for (j=0; j<i; j++) w[j] = 1./h * GammaInv[i*s+j];
44761692a83SJed Brown       ierr = VecZeroEntries(Zdot);CHKERRQ(ierr);
44861692a83SJed Brown       ierr = VecMAXPY(Zdot,i,w,Y);CHKERRQ(ierr);
44961692a83SJed Brown 
450e27a552bSJed Brown       /* Initial guess taken from last stage */
45161692a83SJed Brown       ierr = VecZeroEntries(Y[i]);CHKERRQ(ierr);
45261692a83SJed Brown 
45361692a83SJed Brown       if (!ros->recompute_jacobian && !i) {
45461692a83SJed Brown         ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */
45561692a83SJed Brown       }
45661692a83SJed Brown 
45761692a83SJed Brown       ierr = SNESSolve(snes,PETSC_NULL,Y[i]);CHKERRQ(ierr);
458e27a552bSJed Brown       ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
459e27a552bSJed Brown       ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
460e27a552bSJed Brown       ts->nonlinear_its += its; ts->linear_its += lits;
461e27a552bSJed Brown     }
462*1c3436cfSJed Brown     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
463*1c3436cfSJed Brown     ros->step_taken = PETSC_TRUE;
464e27a552bSJed Brown 
465*1c3436cfSJed Brown     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
466*1c3436cfSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
467*1c3436cfSJed Brown     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
468*1c3436cfSJed Brown     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,0.0,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
469*1c3436cfSJed Brown     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
470*1c3436cfSJed Brown     if (accept) {
471*1c3436cfSJed Brown       /* ignore next_scheme for now */
472e27a552bSJed Brown       ts->ptime += ts->time_step;
473cdbf8f93SLisandro Dalcin       ts->time_step = next_time_step;
474e27a552bSJed Brown       ts->steps++;
475*1c3436cfSJed Brown       break;
476*1c3436cfSJed Brown     } else {                    /* Roll back the current step */
477*1c3436cfSJed Brown       for (i=0; i<s; i++) w[i] = -tab->bt[i];
478*1c3436cfSJed Brown       ierr = VecMAXPY(ts->vec_sol,s,w,Y);CHKERRQ(ierr);
479*1c3436cfSJed Brown       ts->time_step = next_time_step;
480*1c3436cfSJed Brown       ros->step_taken = PETSC_FALSE;
481*1c3436cfSJed Brown     }
482*1c3436cfSJed Brown   }
483*1c3436cfSJed Brown 
484e27a552bSJed Brown   PetscFunctionReturn(0);
485e27a552bSJed Brown }
486e27a552bSJed Brown 
487e27a552bSJed Brown #undef __FUNCT__
488e27a552bSJed Brown #define __FUNCT__ "TSInterpolate_RosW"
489e27a552bSJed Brown static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec X)
490e27a552bSJed Brown {
49161692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
492e27a552bSJed Brown 
493e27a552bSJed Brown   PetscFunctionBegin;
49461692a83SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name);
495e27a552bSJed Brown   PetscFunctionReturn(0);
496e27a552bSJed Brown }
497e27a552bSJed Brown 
498e27a552bSJed Brown /*------------------------------------------------------------*/
499e27a552bSJed Brown #undef __FUNCT__
500e27a552bSJed Brown #define __FUNCT__ "TSReset_RosW"
501e27a552bSJed Brown static PetscErrorCode TSReset_RosW(TS ts)
502e27a552bSJed Brown {
50361692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
504e27a552bSJed Brown   PetscInt       s;
505e27a552bSJed Brown   PetscErrorCode ierr;
506e27a552bSJed Brown 
507e27a552bSJed Brown   PetscFunctionBegin;
50861692a83SJed Brown   if (!ros->tableau) PetscFunctionReturn(0);
50961692a83SJed Brown    s = ros->tableau->s;
51061692a83SJed Brown   ierr = VecDestroyVecs(s,&ros->Y);CHKERRQ(ierr);
51161692a83SJed Brown   ierr = VecDestroy(&ros->Ydot);CHKERRQ(ierr);
51261692a83SJed Brown   ierr = VecDestroy(&ros->Ystage);CHKERRQ(ierr);
51361692a83SJed Brown   ierr = VecDestroy(&ros->Zdot);CHKERRQ(ierr);
51461692a83SJed Brown   ierr = VecDestroy(&ros->Zstage);CHKERRQ(ierr);
51561692a83SJed Brown   ierr = PetscFree(ros->work);CHKERRQ(ierr);
516e27a552bSJed Brown   PetscFunctionReturn(0);
517e27a552bSJed Brown }
518e27a552bSJed Brown 
519e27a552bSJed Brown #undef __FUNCT__
520e27a552bSJed Brown #define __FUNCT__ "TSDestroy_RosW"
521e27a552bSJed Brown static PetscErrorCode TSDestroy_RosW(TS ts)
522e27a552bSJed Brown {
523e27a552bSJed Brown   PetscErrorCode  ierr;
524e27a552bSJed Brown 
525e27a552bSJed Brown   PetscFunctionBegin;
526e27a552bSJed Brown   ierr = TSReset_RosW(ts);CHKERRQ(ierr);
527e27a552bSJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
528e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","",PETSC_NULL);CHKERRQ(ierr);
529e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","",PETSC_NULL);CHKERRQ(ierr);
53061692a83SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","",PETSC_NULL);CHKERRQ(ierr);
531e27a552bSJed Brown   PetscFunctionReturn(0);
532e27a552bSJed Brown }
533e27a552bSJed Brown 
534e27a552bSJed Brown /*
535e27a552bSJed Brown   This defines the nonlinear equation that is to be solved with SNES
536e27a552bSJed Brown   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
537e27a552bSJed Brown */
538e27a552bSJed Brown #undef __FUNCT__
539e27a552bSJed Brown #define __FUNCT__ "SNESTSFormFunction_RosW"
540e27a552bSJed Brown static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec X,Vec F,TS ts)
541e27a552bSJed Brown {
54261692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
543e27a552bSJed Brown   PetscErrorCode ierr;
544e27a552bSJed Brown 
545e27a552bSJed Brown   PetscFunctionBegin;
54661692a83SJed Brown   ierr = VecWAXPY(ros->Ydot,ros->shift,X,ros->Zdot);CHKERRQ(ierr); /* Ydot = shift*X + Zdot */
54761692a83SJed Brown   ierr = VecWAXPY(ros->Ystage,1.0,X,ros->Zstage);CHKERRQ(ierr);    /* Ystage = X + Zstage */
54861692a83SJed Brown   ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,ros->Ydot,F,PETSC_FALSE);CHKERRQ(ierr);
549e27a552bSJed Brown   PetscFunctionReturn(0);
550e27a552bSJed Brown }
551e27a552bSJed Brown 
552e27a552bSJed Brown #undef __FUNCT__
553e27a552bSJed Brown #define __FUNCT__ "SNESTSFormJacobian_RosW"
554e27a552bSJed Brown static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
555e27a552bSJed Brown {
55661692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
557e27a552bSJed Brown   PetscErrorCode ierr;
558e27a552bSJed Brown 
559e27a552bSJed Brown   PetscFunctionBegin;
56061692a83SJed Brown   /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
56161692a83SJed Brown   ierr = TSComputeIJacobian(ts,ros->stage_time,ros->Ystage,ros->Ydot,ros->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr);
562e27a552bSJed Brown   PetscFunctionReturn(0);
563e27a552bSJed Brown }
564e27a552bSJed Brown 
565e27a552bSJed Brown #undef __FUNCT__
566e27a552bSJed Brown #define __FUNCT__ "TSSetUp_RosW"
567e27a552bSJed Brown static PetscErrorCode TSSetUp_RosW(TS ts)
568e27a552bSJed Brown {
56961692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
57061692a83SJed Brown   RosWTableau    tab  = ros->tableau;
571e27a552bSJed Brown   PetscInt       s    = tab->s;
572e27a552bSJed Brown   PetscErrorCode ierr;
573e27a552bSJed Brown 
574e27a552bSJed Brown   PetscFunctionBegin;
57561692a83SJed Brown   if (!ros->tableau) {
576e27a552bSJed Brown     ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);
577e27a552bSJed Brown   }
57861692a83SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ros->Y);CHKERRQ(ierr);
57961692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Ydot);CHKERRQ(ierr);
58061692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Ystage);CHKERRQ(ierr);
58161692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Zdot);CHKERRQ(ierr);
58261692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Zstage);CHKERRQ(ierr);
58361692a83SJed Brown   ierr = PetscMalloc(s*sizeof(ros->work[0]),&ros->work);CHKERRQ(ierr);
584e27a552bSJed Brown   PetscFunctionReturn(0);
585e27a552bSJed Brown }
586e27a552bSJed Brown /*------------------------------------------------------------*/
587e27a552bSJed Brown 
588e27a552bSJed Brown #undef __FUNCT__
589e27a552bSJed Brown #define __FUNCT__ "TSSetFromOptions_RosW"
590e27a552bSJed Brown static PetscErrorCode TSSetFromOptions_RosW(TS ts)
591e27a552bSJed Brown {
59261692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
593e27a552bSJed Brown   PetscErrorCode ierr;
59461692a83SJed Brown   char           rostype[256];
595e27a552bSJed Brown 
596e27a552bSJed Brown   PetscFunctionBegin;
597e27a552bSJed Brown   ierr = PetscOptionsHead("RosW ODE solver options");CHKERRQ(ierr);
598e27a552bSJed Brown   {
59961692a83SJed Brown     RosWTableauLink link;
600e27a552bSJed Brown     PetscInt count,choice;
601e27a552bSJed Brown     PetscBool flg;
602e27a552bSJed Brown     const char **namelist;
60361692a83SJed Brown     SNES snes;
60461692a83SJed Brown 
60561692a83SJed Brown     ierr = PetscStrncpy(rostype,TSRosWDefault,sizeof rostype);CHKERRQ(ierr);
60661692a83SJed Brown     for (link=RosWTableauList,count=0; link; link=link->next,count++) ;
607e27a552bSJed Brown     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
60861692a83SJed Brown     for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
60961692a83SJed Brown     ierr = PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,rostype,&choice,&flg);CHKERRQ(ierr);
61061692a83SJed Brown     ierr = TSRosWSetType(ts,flg ? namelist[choice] : rostype);CHKERRQ(ierr);
611e27a552bSJed Brown     ierr = PetscFree(namelist);CHKERRQ(ierr);
61261692a83SJed Brown 
61361692a83SJed Brown     ierr = PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,PETSC_NULL);CHKERRQ(ierr);
61461692a83SJed Brown 
61561692a83SJed Brown     /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
61661692a83SJed Brown     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
61761692a83SJed Brown     if (!((PetscObject)snes)->type_name) {
61861692a83SJed Brown       ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
61961692a83SJed Brown     }
62061692a83SJed Brown     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
621e27a552bSJed Brown   }
622e27a552bSJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
623e27a552bSJed Brown   PetscFunctionReturn(0);
624e27a552bSJed Brown }
625e27a552bSJed Brown 
626e27a552bSJed Brown #undef __FUNCT__
627e27a552bSJed Brown #define __FUNCT__ "PetscFormatRealArray"
628e27a552bSJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
629e27a552bSJed Brown {
630e27a552bSJed Brown   PetscErrorCode ierr;
631e408995aSJed Brown   PetscInt i;
632e408995aSJed Brown   size_t left,count;
633e27a552bSJed Brown   char *p;
634e27a552bSJed Brown 
635e27a552bSJed Brown   PetscFunctionBegin;
636e408995aSJed Brown   for (i=0,p=buf,left=len; i<n; i++) {
637e408995aSJed Brown     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
638e27a552bSJed Brown     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
639e27a552bSJed Brown     left -= count;
640e27a552bSJed Brown     p += count;
641e27a552bSJed Brown     *p++ = ' ';
642e27a552bSJed Brown   }
643e27a552bSJed Brown   p[i ? 0 : -1] = 0;
644e27a552bSJed Brown   PetscFunctionReturn(0);
645e27a552bSJed Brown }
646e27a552bSJed Brown 
647e27a552bSJed Brown #undef __FUNCT__
648e27a552bSJed Brown #define __FUNCT__ "TSView_RosW"
649e27a552bSJed Brown static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer)
650e27a552bSJed Brown {
65161692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
65261692a83SJed Brown   RosWTableau    tab  = ros->tableau;
653e27a552bSJed Brown   PetscBool      iascii;
654e27a552bSJed Brown   PetscErrorCode ierr;
655e27a552bSJed Brown 
656e27a552bSJed Brown   PetscFunctionBegin;
657e27a552bSJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
658e27a552bSJed Brown   if (iascii) {
65961692a83SJed Brown     const TSRosWType rostype;
660e408995aSJed Brown     PetscInt i;
661e408995aSJed Brown     PetscReal abscissa[512];
662e27a552bSJed Brown     char buf[512];
66361692a83SJed Brown     ierr = TSRosWGetType(ts,&rostype);CHKERRQ(ierr);
66461692a83SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Rosenbrock-W %s\n",rostype);CHKERRQ(ierr);
665e408995aSJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ASum);CHKERRQ(ierr);
66661692a83SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A       = %s\n",buf);CHKERRQ(ierr);
667e408995aSJed Brown     for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i];
668e408995aSJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,abscissa);CHKERRQ(ierr);
669e408995aSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A+Gamma = %s\n",buf);CHKERRQ(ierr);
670e27a552bSJed Brown   }
671e27a552bSJed Brown   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
672e27a552bSJed Brown   PetscFunctionReturn(0);
673e27a552bSJed Brown }
674e27a552bSJed Brown 
675e27a552bSJed Brown #undef __FUNCT__
676e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType"
677e27a552bSJed Brown /*@C
67861692a83SJed Brown   TSRosWSetType - Set the type of Rosenbrock-W scheme
679e27a552bSJed Brown 
680e27a552bSJed Brown   Logically collective
681e27a552bSJed Brown 
682e27a552bSJed Brown   Input Parameter:
683e27a552bSJed Brown +  ts - timestepping context
68461692a83SJed Brown -  rostype - type of Rosenbrock-W scheme
685e27a552bSJed Brown 
686e27a552bSJed Brown   Level: intermediate
687e27a552bSJed Brown 
688e27a552bSJed Brown .seealso: TSRosWGetType()
689e27a552bSJed Brown @*/
69061692a83SJed Brown PetscErrorCode TSRosWSetType(TS ts,const TSRosWType rostype)
691e27a552bSJed Brown {
692e27a552bSJed Brown   PetscErrorCode ierr;
693e27a552bSJed Brown 
694e27a552bSJed Brown   PetscFunctionBegin;
695e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
69661692a83SJed Brown   ierr = PetscTryMethod(ts,"TSRosWSetType_C",(TS,const TSRosWType),(ts,rostype));CHKERRQ(ierr);
697e27a552bSJed Brown   PetscFunctionReturn(0);
698e27a552bSJed Brown }
699e27a552bSJed Brown 
700e27a552bSJed Brown #undef __FUNCT__
701e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType"
702e27a552bSJed Brown /*@C
70361692a83SJed Brown   TSRosWGetType - Get the type of Rosenbrock-W scheme
704e27a552bSJed Brown 
705e27a552bSJed Brown   Logically collective
706e27a552bSJed Brown 
707e27a552bSJed Brown   Input Parameter:
708e27a552bSJed Brown .  ts - timestepping context
709e27a552bSJed Brown 
710e27a552bSJed Brown   Output Parameter:
71161692a83SJed Brown .  rostype - type of Rosenbrock-W scheme
712e27a552bSJed Brown 
713e27a552bSJed Brown   Level: intermediate
714e27a552bSJed Brown 
715e27a552bSJed Brown .seealso: TSRosWGetType()
716e27a552bSJed Brown @*/
71761692a83SJed Brown PetscErrorCode TSRosWGetType(TS ts,const TSRosWType *rostype)
718e27a552bSJed Brown {
719e27a552bSJed Brown   PetscErrorCode ierr;
720e27a552bSJed Brown 
721e27a552bSJed Brown   PetscFunctionBegin;
722e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
72361692a83SJed Brown   ierr = PetscUseMethod(ts,"TSRosWGetType_C",(TS,const TSRosWType*),(ts,rostype));CHKERRQ(ierr);
724e27a552bSJed Brown   PetscFunctionReturn(0);
725e27a552bSJed Brown }
726e27a552bSJed Brown 
727e27a552bSJed Brown #undef __FUNCT__
72861692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian"
729e27a552bSJed Brown /*@C
73061692a83SJed Brown   TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step.
731e27a552bSJed Brown 
732e27a552bSJed Brown   Logically collective
733e27a552bSJed Brown 
734e27a552bSJed Brown   Input Parameter:
735e27a552bSJed Brown +  ts - timestepping context
73661692a83SJed Brown -  flg - PETSC_TRUE to recompute the Jacobian at each stage
737e27a552bSJed Brown 
738e27a552bSJed Brown   Level: intermediate
739e27a552bSJed Brown 
740e27a552bSJed Brown .seealso: TSRosWGetType()
741e27a552bSJed Brown @*/
74261692a83SJed Brown PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg)
743e27a552bSJed Brown {
744e27a552bSJed Brown   PetscErrorCode ierr;
745e27a552bSJed Brown 
746e27a552bSJed Brown   PetscFunctionBegin;
747e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
74861692a83SJed Brown   ierr = PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
749e27a552bSJed Brown   PetscFunctionReturn(0);
750e27a552bSJed Brown }
751e27a552bSJed Brown 
752e27a552bSJed Brown EXTERN_C_BEGIN
753e27a552bSJed Brown #undef __FUNCT__
754e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType_RosW"
75561692a83SJed Brown PetscErrorCode  TSRosWGetType_RosW(TS ts,const TSRosWType *rostype)
756e27a552bSJed Brown {
75761692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
758e27a552bSJed Brown   PetscErrorCode ierr;
759e27a552bSJed Brown 
760e27a552bSJed Brown   PetscFunctionBegin;
76161692a83SJed Brown   if (!ros->tableau) {ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);}
76261692a83SJed Brown   *rostype = ros->tableau->name;
763e27a552bSJed Brown   PetscFunctionReturn(0);
764e27a552bSJed Brown }
765e27a552bSJed Brown #undef __FUNCT__
766e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType_RosW"
76761692a83SJed Brown PetscErrorCode  TSRosWSetType_RosW(TS ts,const TSRosWType rostype)
768e27a552bSJed Brown {
76961692a83SJed Brown   TS_RosW         *ros = (TS_RosW*)ts->data;
770e27a552bSJed Brown   PetscErrorCode  ierr;
771e27a552bSJed Brown   PetscBool       match;
77261692a83SJed Brown   RosWTableauLink link;
773e27a552bSJed Brown 
774e27a552bSJed Brown   PetscFunctionBegin;
77561692a83SJed Brown   if (ros->tableau) {
77661692a83SJed Brown     ierr = PetscStrcmp(ros->tableau->name,rostype,&match);CHKERRQ(ierr);
777e27a552bSJed Brown     if (match) PetscFunctionReturn(0);
778e27a552bSJed Brown   }
77961692a83SJed Brown   for (link = RosWTableauList; link; link=link->next) {
78061692a83SJed Brown     ierr = PetscStrcmp(link->tab.name,rostype,&match);CHKERRQ(ierr);
781e27a552bSJed Brown     if (match) {
782e27a552bSJed Brown       ierr = TSReset_RosW(ts);CHKERRQ(ierr);
78361692a83SJed Brown       ros->tableau = &link->tab;
784e27a552bSJed Brown       PetscFunctionReturn(0);
785e27a552bSJed Brown     }
786e27a552bSJed Brown   }
78761692a83SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype);
788e27a552bSJed Brown   PetscFunctionReturn(0);
789e27a552bSJed Brown }
79061692a83SJed Brown 
791e27a552bSJed Brown #undef __FUNCT__
79261692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian_RosW"
79361692a83SJed Brown PetscErrorCode  TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg)
794e27a552bSJed Brown {
79561692a83SJed Brown   TS_RosW *ros = (TS_RosW*)ts->data;
796e27a552bSJed Brown 
797e27a552bSJed Brown   PetscFunctionBegin;
79861692a83SJed Brown   ros->recompute_jacobian = flg;
799e27a552bSJed Brown   PetscFunctionReturn(0);
800e27a552bSJed Brown }
801e27a552bSJed Brown EXTERN_C_END
802e27a552bSJed Brown 
803e27a552bSJed Brown /* ------------------------------------------------------------ */
804e27a552bSJed Brown /*MC
805e27a552bSJed Brown       TSRosW - ODE solver using Rosenbrock-W schemes
806e27a552bSJed Brown 
807e27a552bSJed Brown   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
808e27a552bSJed Brown   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
809e27a552bSJed Brown   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
810e27a552bSJed Brown 
811e27a552bSJed Brown   Notes:
81261692a83SJed Brown   This method currently only works with autonomous ODE and DAE.
81361692a83SJed Brown 
81461692a83SJed Brown   Developer notes:
81561692a83SJed Brown   Rosenbrock-W methods are typically specified for autonomous ODE
81661692a83SJed Brown 
81761692a83SJed Brown $  xdot = f(x)
81861692a83SJed Brown 
81961692a83SJed Brown   by the stage equations
82061692a83SJed Brown 
82161692a83SJed Brown $  k_i = h f(x_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j
82261692a83SJed Brown 
82361692a83SJed Brown   and step completion formula
82461692a83SJed Brown 
82561692a83SJed Brown $  x_1 = x_0 + sum_j b_j k_j
82661692a83SJed Brown 
82761692a83SJed Brown   with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(x)
82861692a83SJed Brown   and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner,
82961692a83SJed Brown   we define new variables for the stage equations
83061692a83SJed Brown 
83161692a83SJed Brown $  y_i = gamma_ij k_j
83261692a83SJed Brown 
83361692a83SJed Brown   The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define
83461692a83SJed Brown 
83561692a83SJed Brown $  A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-i}
83661692a83SJed Brown 
83761692a83SJed Brown   to rewrite the method as
83861692a83SJed Brown 
83961692a83SJed 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
84061692a83SJed Brown $  x_1 = x_0 + sum_j bt_j y_j
84161692a83SJed Brown 
84261692a83SJed Brown    where we have introduced the mass matrix M. Continue by defining
84361692a83SJed Brown 
84461692a83SJed Brown $  ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j
84561692a83SJed Brown 
84661692a83SJed Brown    or, more compactly in tensor notation
84761692a83SJed Brown 
84861692a83SJed Brown $  Ydot = 1/h (Gamma^{-1} \otimes I) Y .
84961692a83SJed Brown 
85061692a83SJed Brown    Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current
85161692a83SJed Brown    stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the
85261692a83SJed Brown    equation
85361692a83SJed Brown 
85461692a83SJed Brown $  g(x_0 + sum_j a_ij y_j + y_i, ydot_i) = 0
85561692a83SJed Brown 
85661692a83SJed Brown    with initial guess y_i = 0.
857e27a552bSJed Brown 
858e27a552bSJed Brown   Level: beginner
859e27a552bSJed Brown 
860e27a552bSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSRosWRegister()
861e27a552bSJed Brown 
862e27a552bSJed Brown M*/
863e27a552bSJed Brown EXTERN_C_BEGIN
864e27a552bSJed Brown #undef __FUNCT__
865e27a552bSJed Brown #define __FUNCT__ "TSCreate_RosW"
866e27a552bSJed Brown PetscErrorCode  TSCreate_RosW(TS ts)
867e27a552bSJed Brown {
86861692a83SJed Brown   TS_RosW        *ros;
869e27a552bSJed Brown   PetscErrorCode ierr;
870e27a552bSJed Brown 
871e27a552bSJed Brown   PetscFunctionBegin;
872e27a552bSJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
873e27a552bSJed Brown   ierr = TSRosWInitializePackage(PETSC_NULL);CHKERRQ(ierr);
874e27a552bSJed Brown #endif
875e27a552bSJed Brown 
876e27a552bSJed Brown   ts->ops->reset          = TSReset_RosW;
877e27a552bSJed Brown   ts->ops->destroy        = TSDestroy_RosW;
878e27a552bSJed Brown   ts->ops->view           = TSView_RosW;
879e27a552bSJed Brown   ts->ops->setup          = TSSetUp_RosW;
880e27a552bSJed Brown   ts->ops->step           = TSStep_RosW;
881e27a552bSJed Brown   ts->ops->interpolate    = TSInterpolate_RosW;
882*1c3436cfSJed Brown   ts->ops->evaluatestep   = TSEvaluateStep_RosW;
883e27a552bSJed Brown   ts->ops->setfromoptions = TSSetFromOptions_RosW;
884e27a552bSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_RosW;
885e27a552bSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_RosW;
886e27a552bSJed Brown 
88761692a83SJed Brown   ierr = PetscNewLog(ts,TS_RosW,&ros);CHKERRQ(ierr);
88861692a83SJed Brown   ts->data = (void*)ros;
889e27a552bSJed Brown 
890e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","TSRosWGetType_RosW",TSRosWGetType_RosW);CHKERRQ(ierr);
891e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","TSRosWSetType_RosW",TSRosWSetType_RosW);CHKERRQ(ierr);
89261692a83SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","TSRosWSetRecomputeJacobian_RosW",TSRosWSetRecomputeJacobian_RosW);CHKERRQ(ierr);
893e27a552bSJed Brown   PetscFunctionReturn(0);
894e27a552bSJed Brown }
895e27a552bSJed Brown EXTERN_C_END
896