xref: /petsc/src/ts/impls/rosw/rosw.c (revision fd96d5b03496a934fb926c5b5e2d2bda53726738)
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 */
27*fd96d5b0SEmil Constantinescu   PetscReal *Gamma;             /* Stage table, lower triangular*/
28*fd96d5b0SEmil Constantinescu   PetscInt *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 */
52e27a552bSJed Brown   PetscScalar *work;            /* Scalar work */
53e27a552bSJed Brown   PetscReal   shift;
54e27a552bSJed Brown   PetscReal   stage_time;
5561692a83SJed Brown   PetscBool   recompute_jacobian; /* Recompute the Jacobian at each stage, default is to freeze the Jacobian at the start of each 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}},
13461692a83SJed Brown       b[2] = {0.5,0.5};
135fe7e6d57SJed Brown     ierr = TSRosWRegister(TSROSW2P,2,2,&A[0][0],&Gamma[0][0],b,PETSC_NULL);CHKERRQ(ierr);
136e27a552bSJed Brown   }
137e27a552bSJed Brown   {
13861692a83SJed Brown     const PetscReal g = 1. - 1./PetscSqrtReal(2.0);
139e27a552bSJed Brown     const PetscReal
14061692a83SJed Brown       A[2][2] = {{0,0}, {1.,0}},
14161692a83SJed Brown       Gamma[2][2] = {{g,0}, {-2.*g,g}},
14261692a83SJed Brown       b[2] = {0.5,0.5};
143fe7e6d57SJed Brown     ierr = TSRosWRegister(TSROSW2M,2,2,&A[0][0],&Gamma[0][0],b,PETSC_NULL);CHKERRQ(ierr);
144fe7e6d57SJed Brown   }
145fe7e6d57SJed Brown   {
146fe7e6d57SJed Brown     const PetscReal g = 7.8867513459481287e-01;
147fe7e6d57SJed Brown     const PetscReal
148fe7e6d57SJed Brown       A[3][3] = {{0,0,0},
149fe7e6d57SJed Brown                  {1.5773502691896257e+00,0,0},
150fe7e6d57SJed Brown                  {0.5,0,0}},
151fe7e6d57SJed Brown       Gamma[3][3] = {{g,0,0},
152fe7e6d57SJed Brown                      {-1.5773502691896257e+00,g,0},
153fe7e6d57SJed Brown                      {-6.7075317547305480e-01,1.7075317547305482e-01,g}},
154fe7e6d57SJed Brown       b[3] = {1.0566243270259355e-01,4.9038105676657971e-02,8.4529946162074843e-01},
155fe7e6d57SJed Brown       b2[3] = {-1.7863279495408180e-01,1./3.,8.4529946162074843e-01};
156fe7e6d57SJed Brown     ierr = TSRosWRegister(TSROSWRA3PW,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
157fe7e6d57SJed Brown   }
158fe7e6d57SJed Brown   {
159fe7e6d57SJed Brown     const PetscReal g = 4.3586652150845900e-01;
160fe7e6d57SJed Brown     const PetscReal
161fe7e6d57SJed Brown       A[4][4] = {{0,0,0,0},
162fe7e6d57SJed Brown                  {8.7173304301691801e-01,0,0,0},
163fe7e6d57SJed Brown                  {8.4457060015369423e-01,-1.1299064236484185e-01,0,0},
164fe7e6d57SJed Brown                  {0,0,1.,0}},
165fe7e6d57SJed Brown       Gamma[4][4] = {{g,0,0,0},
166fe7e6d57SJed Brown                      {-8.7173304301691801e-01,g,0,0},
167fe7e6d57SJed Brown                      {-9.0338057013044082e-01,5.4180672388095326e-02,g,0},
168fe7e6d57SJed Brown                      {2.4212380706095346e-01,-1.2232505839045147e+00,5.4526025533510214e-01,g}},
169fe7e6d57SJed Brown       b[4] = {2.4212380706095346e-01,-1.2232505839045147e+00,1.5452602553351020e+00,4.3586652150845900e-01},
170fe7e6d57SJed Brown       b2[4] = {3.7810903145819369e-01,-9.6042292212423178e-02,5.0000000000000000e-01,2.1793326075422950e-01};
171fe7e6d57SJed Brown     ierr = TSRosWRegister(TSROSWRA34PW2,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr);
172e27a552bSJed Brown   }
173e27a552bSJed Brown   PetscFunctionReturn(0);
174e27a552bSJed Brown }
175e27a552bSJed Brown 
176e27a552bSJed Brown #undef __FUNCT__
177e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterDestroy"
178e27a552bSJed Brown /*@C
179e27a552bSJed Brown    TSRosWRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister().
180e27a552bSJed Brown 
181e27a552bSJed Brown    Not Collective
182e27a552bSJed Brown 
183e27a552bSJed Brown    Level: advanced
184e27a552bSJed Brown 
185e27a552bSJed Brown .keywords: TSRosW, register, destroy
186e27a552bSJed Brown .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegisterDynamic()
187e27a552bSJed Brown @*/
188e27a552bSJed Brown PetscErrorCode TSRosWRegisterDestroy(void)
189e27a552bSJed Brown {
190e27a552bSJed Brown   PetscErrorCode ierr;
19161692a83SJed Brown   RosWTableauLink link;
192e27a552bSJed Brown 
193e27a552bSJed Brown   PetscFunctionBegin;
19461692a83SJed Brown   while ((link = RosWTableauList)) {
19561692a83SJed Brown     RosWTableau t = &link->tab;
19661692a83SJed Brown     RosWTableauList = link->next;
19761692a83SJed Brown     ierr = PetscFree5(t->A,t->Gamma,t->b,t->ASum,t->GammaSum);CHKERRQ(ierr);
19861692a83SJed Brown     ierr = PetscFree3(t->At,t->bt,t->GammaInv);CHKERRQ(ierr);
199fe7e6d57SJed Brown     ierr = PetscFree2(t->bembed,t->bembedt);CHKERRQ(ierr);
200*fd96d5b0SEmil Constantinescu     ierr = PetscFree(t->GammaZeroDiag);CHKERRQ(ierr);
201e27a552bSJed Brown     ierr = PetscFree(t->name);CHKERRQ(ierr);
202e27a552bSJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
203e27a552bSJed Brown   }
204e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_FALSE;
205e27a552bSJed Brown   PetscFunctionReturn(0);
206e27a552bSJed Brown }
207e27a552bSJed Brown 
208e27a552bSJed Brown #undef __FUNCT__
209e27a552bSJed Brown #define __FUNCT__ "TSRosWInitializePackage"
210e27a552bSJed Brown /*@C
211e27a552bSJed Brown   TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called
212e27a552bSJed Brown   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RosW()
213e27a552bSJed Brown   when using static libraries.
214e27a552bSJed Brown 
215e27a552bSJed Brown   Input Parameter:
216e27a552bSJed Brown   path - The dynamic library path, or PETSC_NULL
217e27a552bSJed Brown 
218e27a552bSJed Brown   Level: developer
219e27a552bSJed Brown 
220e27a552bSJed Brown .keywords: TS, TSRosW, initialize, package
221e27a552bSJed Brown .seealso: PetscInitialize()
222e27a552bSJed Brown @*/
223e27a552bSJed Brown PetscErrorCode TSRosWInitializePackage(const char path[])
224e27a552bSJed Brown {
225e27a552bSJed Brown   PetscErrorCode ierr;
226e27a552bSJed Brown 
227e27a552bSJed Brown   PetscFunctionBegin;
228e27a552bSJed Brown   if (TSRosWPackageInitialized) PetscFunctionReturn(0);
229e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_TRUE;
230e27a552bSJed Brown   ierr = TSRosWRegisterAll();CHKERRQ(ierr);
231e27a552bSJed Brown   ierr = PetscRegisterFinalize(TSRosWFinalizePackage);CHKERRQ(ierr);
232e27a552bSJed Brown   PetscFunctionReturn(0);
233e27a552bSJed Brown }
234e27a552bSJed Brown 
235e27a552bSJed Brown #undef __FUNCT__
236e27a552bSJed Brown #define __FUNCT__ "TSRosWFinalizePackage"
237e27a552bSJed Brown /*@C
238e27a552bSJed Brown   TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is
239e27a552bSJed Brown   called from PetscFinalize().
240e27a552bSJed Brown 
241e27a552bSJed Brown   Level: developer
242e27a552bSJed Brown 
243e27a552bSJed Brown .keywords: Petsc, destroy, package
244e27a552bSJed Brown .seealso: PetscFinalize()
245e27a552bSJed Brown @*/
246e27a552bSJed Brown PetscErrorCode TSRosWFinalizePackage(void)
247e27a552bSJed Brown {
248e27a552bSJed Brown   PetscErrorCode ierr;
249e27a552bSJed Brown 
250e27a552bSJed Brown   PetscFunctionBegin;
251e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_FALSE;
252e27a552bSJed Brown   ierr = TSRosWRegisterDestroy();CHKERRQ(ierr);
253e27a552bSJed Brown   PetscFunctionReturn(0);
254e27a552bSJed Brown }
255e27a552bSJed Brown 
256e27a552bSJed Brown #undef __FUNCT__
257e27a552bSJed Brown #define __FUNCT__ "TSRosWRegister"
258e27a552bSJed Brown /*@C
25961692a83SJed Brown    TSRosWRegister - register a Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
260e27a552bSJed Brown 
261e27a552bSJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
262e27a552bSJed Brown 
263e27a552bSJed Brown    Input Parameters:
264e27a552bSJed Brown +  name - identifier for method
265e27a552bSJed Brown .  order - approximation order of method
266e27a552bSJed Brown .  s - number of stages, this is the dimension of the matrices below
26761692a83SJed Brown .  A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular
26861692a83SJed Brown .  Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal
269fe7e6d57SJed Brown .  b - Step completion table (dimension s)
270fe7e6d57SJed Brown -  bembed - Step completion table for a scheme of order one less (dimension s, PETSC_NULL if no embedded scheme is available)
271e27a552bSJed Brown 
272e27a552bSJed Brown    Notes:
27361692a83SJed Brown    Several Rosenbrock W methods are provided, this function is only needed to create new methods.
274e27a552bSJed Brown 
275e27a552bSJed Brown    Level: advanced
276e27a552bSJed Brown 
277e27a552bSJed Brown .keywords: TS, register
278e27a552bSJed Brown 
279e27a552bSJed Brown .seealso: TSRosW
280e27a552bSJed Brown @*/
281e27a552bSJed Brown PetscErrorCode TSRosWRegister(const TSRosWType name,PetscInt order,PetscInt s,
282fe7e6d57SJed Brown                               const PetscReal A[],const PetscReal Gamma[],const PetscReal b[],const PetscReal bembed[])
283e27a552bSJed Brown {
284e27a552bSJed Brown   PetscErrorCode ierr;
28561692a83SJed Brown   RosWTableauLink link;
28661692a83SJed Brown   RosWTableau t;
28761692a83SJed Brown   PetscInt i,j,k;
28861692a83SJed Brown   PetscScalar *GammaInv;
289e27a552bSJed Brown 
290e27a552bSJed Brown   PetscFunctionBegin;
291fe7e6d57SJed Brown   PetscValidCharPointer(name,1);
292fe7e6d57SJed Brown   PetscValidPointer(A,4);
293fe7e6d57SJed Brown   PetscValidPointer(Gamma,5);
294fe7e6d57SJed Brown   PetscValidPointer(b,6);
295fe7e6d57SJed Brown   if (bembed) PetscValidPointer(bembed,7);
296fe7e6d57SJed Brown 
297e27a552bSJed Brown   ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
298e27a552bSJed Brown   ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
299e27a552bSJed Brown   t = &link->tab;
300e27a552bSJed Brown   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
301e27a552bSJed Brown   t->order = order;
302e27a552bSJed Brown   t->s = s;
30361692a83SJed 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);
304*fd96d5b0SEmil Constantinescu   ierr = PetscMalloc4(s*s,PetscReal,&t->At,s,PetscReal,&t->bt,s*s,PetscReal,&t->GammaInv,s,PetscInt,&t->GammaZeroDiag);CHKERRQ(ierr);
305e27a552bSJed Brown   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
30661692a83SJed Brown   ierr = PetscMemcpy(t->Gamma,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr);
30761692a83SJed Brown   ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
308fe7e6d57SJed Brown   if (bembed) {
309fe7e6d57SJed Brown     ierr = PetscMalloc2(s,PetscReal,&t->bembed,s,PetscReal,&t->bembedt);CHKERRQ(ierr);
310fe7e6d57SJed Brown     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
311fe7e6d57SJed Brown   }
31261692a83SJed Brown   for (i=0; i<s; i++) {
31361692a83SJed Brown     t->ASum[i] = 0;
31461692a83SJed Brown     t->GammaSum[i] = 0;
31561692a83SJed Brown     for (j=0; j<s; j++) {
31661692a83SJed Brown       t->ASum[i] += A[i*s+j];
317fe7e6d57SJed Brown       t->GammaSum[i] += Gamma[i*s+j];
31861692a83SJed Brown     }
31961692a83SJed Brown   }
320*fd96d5b0SEmil Constantinescu 
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];
323*fd96d5b0SEmil Constantinescu   for (i=0; i<s; i++) {
324*fd96d5b0SEmil Constantinescu     if(Gamma[i*s+i]==0.0){
325*fd96d5b0SEmil Constantinescu       GammaInv[i*s+i]=1.0;
326*fd96d5b0SEmil Constantinescu       t->GammaZeroDiag[i]=1;
327*fd96d5b0SEmil Constantinescu     } else {
328*fd96d5b0SEmil Constantinescu       t->GammaZeroDiag[i]=0;
329*fd96d5b0SEmil Constantinescu     }
330*fd96d5b0SEmil Constantinescu   }
331*fd96d5b0SEmil Constantinescu 
33261692a83SJed Brown   switch (s) {
33361692a83SJed Brown   case 1: GammaInv[0] = 1./GammaInv[0]; break;
33461692a83SJed Brown   case 2: ierr = Kernel_A_gets_inverse_A_2(GammaInv,0);CHKERRQ(ierr); break;
33561692a83SJed Brown   case 3: ierr = Kernel_A_gets_inverse_A_3(GammaInv,0);CHKERRQ(ierr); break;
33661692a83SJed Brown   case 4: ierr = Kernel_A_gets_inverse_A_4(GammaInv,0);CHKERRQ(ierr); break;
33761692a83SJed Brown   case 5: {
33861692a83SJed Brown     PetscInt ipvt5[5];
33961692a83SJed Brown     MatScalar work5[5*5];
34061692a83SJed Brown     ierr = Kernel_A_gets_inverse_A_5(GammaInv,ipvt5,work5,0);CHKERRQ(ierr); break;
34161692a83SJed Brown   }
34261692a83SJed Brown   case 6: ierr = Kernel_A_gets_inverse_A_6(GammaInv,0);CHKERRQ(ierr); break;
34361692a83SJed Brown   case 7: ierr = Kernel_A_gets_inverse_A_7(GammaInv,0);CHKERRQ(ierr); break;
34461692a83SJed Brown   default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for %D stages",s);
34561692a83SJed Brown   }
34661692a83SJed Brown   for (i=0; i<s*s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]);
34761692a83SJed Brown   ierr = PetscFree(GammaInv);CHKERRQ(ierr);
34861692a83SJed Brown   for (i=0; i<s; i++) {
34961692a83SJed Brown     for (j=0; j<s; j++) {
35061692a83SJed Brown       t->At[i*s+j] = 0;
35161692a83SJed Brown       for (k=0; k<s; k++) {
35261692a83SJed Brown         t->At[i*s+j] += t->A[i*s+k] * t->GammaInv[k*s+j];
35361692a83SJed Brown       }
35461692a83SJed Brown     }
35561692a83SJed Brown     t->bt[i] = 0;
35661692a83SJed Brown     for (j=0; j<s; j++) {
35761692a83SJed Brown       t->bt[i] += t->b[j] * t->GammaInv[j*s+i];
35861692a83SJed Brown     }
359fe7e6d57SJed Brown     if (bembed) {
360fe7e6d57SJed Brown       t->bembedt[i] = 0;
361fe7e6d57SJed Brown       for (j=0; j<s; j++) {
362fe7e6d57SJed Brown         t->bembedt[i] += t->bembed[j] * t->GammaInv[j*s+i];
363fe7e6d57SJed Brown       }
364fe7e6d57SJed Brown     }
36561692a83SJed Brown   }
36661692a83SJed Brown   link->next = RosWTableauList;
36761692a83SJed Brown   RosWTableauList = link;
368e27a552bSJed Brown   PetscFunctionReturn(0);
369e27a552bSJed Brown }
370e27a552bSJed Brown 
371*fd96d5b0SEmil Constantinescu /*
372*fd96d5b0SEmil Constantinescu   This defines the nonlinear equation that is to be solved with SNES
373*fd96d5b0SEmil Constantinescu   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
374*fd96d5b0SEmil Constantinescu */
375*fd96d5b0SEmil Constantinescu #undef __FUNCT__
376*fd96d5b0SEmil Constantinescu #define __FUNCT__ "SNESTSFormFunction_RosW"
377*fd96d5b0SEmil Constantinescu static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec X,Vec F,TS ts)
378*fd96d5b0SEmil Constantinescu {
379*fd96d5b0SEmil Constantinescu   TS_RosW        *ros = (TS_RosW*)ts->data;
380*fd96d5b0SEmil Constantinescu   PetscErrorCode ierr;
381*fd96d5b0SEmil Constantinescu 
382*fd96d5b0SEmil Constantinescu   PetscFunctionBegin;
383*fd96d5b0SEmil Constantinescu   ierr = VecWAXPY(ros->Ydot,ros->shift,X,ros->Zdot);CHKERRQ(ierr); /* Ydot = shift*X + Zdot */
384*fd96d5b0SEmil Constantinescu   ierr = VecWAXPY(ros->Ystage,1.0,X,ros->Zstage);CHKERRQ(ierr);    /* Ystage = X + Zstage */
385*fd96d5b0SEmil Constantinescu   ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,ros->Ydot,F,PETSC_FALSE);CHKERRQ(ierr);
386*fd96d5b0SEmil Constantinescu   PetscFunctionReturn(0);
387*fd96d5b0SEmil Constantinescu }
388*fd96d5b0SEmil Constantinescu 
389*fd96d5b0SEmil Constantinescu /*
390*fd96d5b0SEmil Constantinescu   This defines the nonlinear equation that is to be solved with SNES for explicit stages
391*fd96d5b0SEmil Constantinescu   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
392*fd96d5b0SEmil Constantinescu */
393*fd96d5b0SEmil Constantinescu #undef __FUNCT__
394*fd96d5b0SEmil Constantinescu #define __FUNCT__ "SNESTSFormFunction_RosWExplicit"
395*fd96d5b0SEmil Constantinescu static PetscErrorCode SNESTSFormFunction_RosWExplicit(SNES snes,Vec X,Vec F,TS ts)
396*fd96d5b0SEmil Constantinescu {
397*fd96d5b0SEmil Constantinescu   TS_RosW        *ros = (TS_RosW*)ts->data;
398*fd96d5b0SEmil Constantinescu   PetscErrorCode ierr;
399*fd96d5b0SEmil Constantinescu 
400*fd96d5b0SEmil Constantinescu   PetscFunctionBegin;
401*fd96d5b0SEmil Constantinescu   ierr = VecCopy(X,ros->Ydot);CHKERRQ(ierr);
402*fd96d5b0SEmil Constantinescu   ierr = VecScale(ros->Ydot,ros->shift);CHKERRQ(ierr); /* Ydot = shift*X*/
403*fd96d5b0SEmil Constantinescu   ierr = VecWAXPY(ros->Ystage,1.0,X,ros->Zstage);CHKERRQ(ierr);    /* Ystage = X + Zstage */
404*fd96d5b0SEmil Constantinescu   ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,ros->Ydot,F,PETSC_FALSE);CHKERRQ(ierr);
405*fd96d5b0SEmil Constantinescu   PetscFunctionReturn(0);
406*fd96d5b0SEmil Constantinescu }
407*fd96d5b0SEmil Constantinescu 
408*fd96d5b0SEmil Constantinescu 
409e27a552bSJed Brown #undef __FUNCT__
410e27a552bSJed Brown #define __FUNCT__ "TSStep_RosW"
411e27a552bSJed Brown static PetscErrorCode TSStep_RosW(TS ts)
412e27a552bSJed Brown {
41361692a83SJed Brown   TS_RosW         *ros = (TS_RosW*)ts->data;
41461692a83SJed Brown   RosWTableau     tab  = ros->tableau;
415e27a552bSJed Brown   const PetscInt  s    = tab->s;
41661692a83SJed Brown   const PetscReal *At  = tab->At,*Gamma = tab->Gamma,*bt = tab->bt,*ASum = tab->ASum,*GammaInv = tab->GammaInv;
417*fd96d5b0SEmil Constantinescu   const PetscInt  *GammaZeroDiag = tab->GammaZeroDiag;
41861692a83SJed Brown   PetscScalar     *w   = ros->work;
41961692a83SJed Brown   Vec             *Y   = ros->Y,Zdot = ros->Zdot,Zstage = ros->Zstage;
420e27a552bSJed Brown   SNES            snes;
421e27a552bSJed Brown   PetscInt        i,j,its,lits;
422cdbf8f93SLisandro Dalcin   PetscReal       next_time_step;
423e27a552bSJed Brown   PetscReal       h,t;
424e27a552bSJed Brown   PetscErrorCode  ierr;
425e27a552bSJed Brown 
426e27a552bSJed Brown   PetscFunctionBegin;
427e27a552bSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
428cdbf8f93SLisandro Dalcin   next_time_step = ts->time_step;
429cdbf8f93SLisandro Dalcin   h = ts->time_step;
430e27a552bSJed Brown   t = ts->ptime;
431e27a552bSJed Brown 
432e27a552bSJed Brown   for (i=0; i<s; i++) {
43361692a83SJed Brown     ros->stage_time = t + h*ASum[i];
434*fd96d5b0SEmil Constantinescu     if(GammaZeroDiag[i]==0){
43561692a83SJed Brown       ros->shift = 1./(h*Gamma[i*s+i]);
436*fd96d5b0SEmil Constantinescu     } else {
437*fd96d5b0SEmil Constantinescu       ros->shift = 1./h;
438*fd96d5b0SEmil Constantinescu     }
43961692a83SJed Brown     ierr = VecCopy(ts->vec_sol,Zstage);CHKERRQ(ierr);
44061692a83SJed Brown     ierr = VecMAXPY(Zstage,i,&At[i*s+0],Y);CHKERRQ(ierr);
44161692a83SJed Brown 
44261692a83SJed Brown     for (j=0; j<i; j++) w[j] = 1./h * GammaInv[i*s+j];
44361692a83SJed Brown     ierr = VecZeroEntries(Zdot);CHKERRQ(ierr);
44461692a83SJed Brown     ierr = VecMAXPY(Zdot,i,w,Y);CHKERRQ(ierr);
44561692a83SJed Brown 
446*fd96d5b0SEmil Constantinescu     if(GammaZeroDiag[i]==1){ /* computing an explicit stage */
447*fd96d5b0SEmil Constantinescu       ts->ops->snesfunction   = SNESTSFormFunction_RosWExplicit;
448*fd96d5b0SEmil Constantinescu     }
449e27a552bSJed Brown     /* Initial guess taken from last stage */
45061692a83SJed Brown     ierr = VecZeroEntries(Y[i]);CHKERRQ(ierr);
45161692a83SJed Brown 
45261692a83SJed Brown     if (!ros->recompute_jacobian && !i) {
45361692a83SJed Brown       ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */
45461692a83SJed Brown     }
45561692a83SJed Brown 
45661692a83SJed Brown     ierr = SNESSolve(snes,PETSC_NULL,Y[i]);CHKERRQ(ierr);
457e27a552bSJed Brown     ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
458e27a552bSJed Brown     ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
459*fd96d5b0SEmil Constantinescu 
460*fd96d5b0SEmil Constantinescu     if(GammaZeroDiag[i]==1){ /* switching back to the implicit stage function */
461*fd96d5b0SEmil Constantinescu       ts->ops->snesfunction   = SNESTSFormFunction_RosW;
462*fd96d5b0SEmil Constantinescu     }
463e27a552bSJed Brown     ts->nonlinear_its += its; ts->linear_its += lits;
464e27a552bSJed Brown   }
46561692a83SJed Brown   ierr = VecMAXPY(ts->vec_sol,s,bt,Y);CHKERRQ(ierr);
466e27a552bSJed Brown 
467e27a552bSJed Brown   ts->ptime += ts->time_step;
468cdbf8f93SLisandro Dalcin   ts->time_step = next_time_step;
469e27a552bSJed Brown   ts->steps++;
470e27a552bSJed Brown   PetscFunctionReturn(0);
471e27a552bSJed Brown }
472e27a552bSJed Brown 
473e27a552bSJed Brown #undef __FUNCT__
474e27a552bSJed Brown #define __FUNCT__ "TSInterpolate_RosW"
475e27a552bSJed Brown static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec X)
476e27a552bSJed Brown {
47761692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
478e27a552bSJed Brown 
479e27a552bSJed Brown   PetscFunctionBegin;
48061692a83SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name);
481e27a552bSJed Brown   PetscFunctionReturn(0);
482e27a552bSJed Brown }
483e27a552bSJed Brown 
484e27a552bSJed Brown /*------------------------------------------------------------*/
485e27a552bSJed Brown #undef __FUNCT__
486e27a552bSJed Brown #define __FUNCT__ "TSReset_RosW"
487e27a552bSJed Brown static PetscErrorCode TSReset_RosW(TS ts)
488e27a552bSJed Brown {
48961692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
490e27a552bSJed Brown   PetscInt       s;
491e27a552bSJed Brown   PetscErrorCode ierr;
492e27a552bSJed Brown 
493e27a552bSJed Brown   PetscFunctionBegin;
49461692a83SJed Brown   if (!ros->tableau) PetscFunctionReturn(0);
49561692a83SJed Brown    s = ros->tableau->s;
49661692a83SJed Brown   ierr = VecDestroyVecs(s,&ros->Y);CHKERRQ(ierr);
49761692a83SJed Brown   ierr = VecDestroy(&ros->Ydot);CHKERRQ(ierr);
49861692a83SJed Brown   ierr = VecDestroy(&ros->Ystage);CHKERRQ(ierr);
49961692a83SJed Brown   ierr = VecDestroy(&ros->Zdot);CHKERRQ(ierr);
50061692a83SJed Brown   ierr = VecDestroy(&ros->Zstage);CHKERRQ(ierr);
50161692a83SJed Brown   ierr = PetscFree(ros->work);CHKERRQ(ierr);
502e27a552bSJed Brown   PetscFunctionReturn(0);
503e27a552bSJed Brown }
504e27a552bSJed Brown 
505e27a552bSJed Brown #undef __FUNCT__
506e27a552bSJed Brown #define __FUNCT__ "TSDestroy_RosW"
507e27a552bSJed Brown static PetscErrorCode TSDestroy_RosW(TS ts)
508e27a552bSJed Brown {
509e27a552bSJed Brown   PetscErrorCode  ierr;
510e27a552bSJed Brown 
511e27a552bSJed Brown   PetscFunctionBegin;
512e27a552bSJed Brown   ierr = TSReset_RosW(ts);CHKERRQ(ierr);
513e27a552bSJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
514e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","",PETSC_NULL);CHKERRQ(ierr);
515e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","",PETSC_NULL);CHKERRQ(ierr);
51661692a83SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","",PETSC_NULL);CHKERRQ(ierr);
517e27a552bSJed Brown   PetscFunctionReturn(0);
518e27a552bSJed Brown }
519e27a552bSJed Brown 
520e27a552bSJed Brown #undef __FUNCT__
521e27a552bSJed Brown #define __FUNCT__ "SNESTSFormJacobian_RosW"
522e27a552bSJed Brown static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
523e27a552bSJed Brown {
52461692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
525e27a552bSJed Brown   PetscErrorCode ierr;
526e27a552bSJed Brown 
527e27a552bSJed Brown   PetscFunctionBegin;
52861692a83SJed Brown   /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
52961692a83SJed Brown   ierr = TSComputeIJacobian(ts,ros->stage_time,ros->Ystage,ros->Ydot,ros->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr);
530e27a552bSJed Brown   PetscFunctionReturn(0);
531e27a552bSJed Brown }
532e27a552bSJed Brown 
533e27a552bSJed Brown #undef __FUNCT__
534e27a552bSJed Brown #define __FUNCT__ "TSSetUp_RosW"
535e27a552bSJed Brown static PetscErrorCode TSSetUp_RosW(TS ts)
536e27a552bSJed Brown {
53761692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
53861692a83SJed Brown   RosWTableau    tab  = ros->tableau;
539e27a552bSJed Brown   PetscInt       s    = tab->s;
540e27a552bSJed Brown   PetscErrorCode ierr;
541e27a552bSJed Brown 
542e27a552bSJed Brown   PetscFunctionBegin;
54361692a83SJed Brown   if (!ros->tableau) {
544e27a552bSJed Brown     ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);
545e27a552bSJed Brown   }
54661692a83SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ros->Y);CHKERRQ(ierr);
54761692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Ydot);CHKERRQ(ierr);
54861692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Ystage);CHKERRQ(ierr);
54961692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Zdot);CHKERRQ(ierr);
55061692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Zstage);CHKERRQ(ierr);
55161692a83SJed Brown   ierr = PetscMalloc(s*sizeof(ros->work[0]),&ros->work);CHKERRQ(ierr);
552e27a552bSJed Brown   PetscFunctionReturn(0);
553e27a552bSJed Brown }
554e27a552bSJed Brown /*------------------------------------------------------------*/
555e27a552bSJed Brown 
556e27a552bSJed Brown #undef __FUNCT__
557e27a552bSJed Brown #define __FUNCT__ "TSSetFromOptions_RosW"
558e27a552bSJed Brown static PetscErrorCode TSSetFromOptions_RosW(TS ts)
559e27a552bSJed Brown {
56061692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
561e27a552bSJed Brown   PetscErrorCode ierr;
56261692a83SJed Brown   char           rostype[256];
563e27a552bSJed Brown 
564e27a552bSJed Brown   PetscFunctionBegin;
565e27a552bSJed Brown   ierr = PetscOptionsHead("RosW ODE solver options");CHKERRQ(ierr);
566e27a552bSJed Brown   {
56761692a83SJed Brown     RosWTableauLink link;
568e27a552bSJed Brown     PetscInt count,choice;
569e27a552bSJed Brown     PetscBool flg;
570e27a552bSJed Brown     const char **namelist;
57161692a83SJed Brown     SNES snes;
57261692a83SJed Brown 
57361692a83SJed Brown     ierr = PetscStrncpy(rostype,TSRosWDefault,sizeof rostype);CHKERRQ(ierr);
57461692a83SJed Brown     for (link=RosWTableauList,count=0; link; link=link->next,count++) ;
575e27a552bSJed Brown     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
57661692a83SJed Brown     for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
57761692a83SJed Brown     ierr = PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,rostype,&choice,&flg);CHKERRQ(ierr);
57861692a83SJed Brown     ierr = TSRosWSetType(ts,flg ? namelist[choice] : rostype);CHKERRQ(ierr);
579e27a552bSJed Brown     ierr = PetscFree(namelist);CHKERRQ(ierr);
58061692a83SJed Brown 
58161692a83SJed Brown     ierr = PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,PETSC_NULL);CHKERRQ(ierr);
58261692a83SJed Brown 
58361692a83SJed Brown     /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
58461692a83SJed Brown     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
58561692a83SJed Brown     if (!((PetscObject)snes)->type_name) {
58661692a83SJed Brown       ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
58761692a83SJed Brown     }
58861692a83SJed Brown     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
589e27a552bSJed Brown   }
590e27a552bSJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
591e27a552bSJed Brown   PetscFunctionReturn(0);
592e27a552bSJed Brown }
593e27a552bSJed Brown 
594e27a552bSJed Brown #undef __FUNCT__
595e27a552bSJed Brown #define __FUNCT__ "PetscFormatRealArray"
596e27a552bSJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
597e27a552bSJed Brown {
598e27a552bSJed Brown   PetscErrorCode ierr;
599e408995aSJed Brown   PetscInt i;
600e408995aSJed Brown   size_t left,count;
601e27a552bSJed Brown   char *p;
602e27a552bSJed Brown 
603e27a552bSJed Brown   PetscFunctionBegin;
604e408995aSJed Brown   for (i=0,p=buf,left=len; i<n; i++) {
605e408995aSJed Brown     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
606e27a552bSJed Brown     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
607e27a552bSJed Brown     left -= count;
608e27a552bSJed Brown     p += count;
609e27a552bSJed Brown     *p++ = ' ';
610e27a552bSJed Brown   }
611e27a552bSJed Brown   p[i ? 0 : -1] = 0;
612e27a552bSJed Brown   PetscFunctionReturn(0);
613e27a552bSJed Brown }
614e27a552bSJed Brown 
615e27a552bSJed Brown #undef __FUNCT__
616e27a552bSJed Brown #define __FUNCT__ "TSView_RosW"
617e27a552bSJed Brown static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer)
618e27a552bSJed Brown {
61961692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
62061692a83SJed Brown   RosWTableau    tab  = ros->tableau;
621e27a552bSJed Brown   PetscBool      iascii;
622e27a552bSJed Brown   PetscErrorCode ierr;
623e27a552bSJed Brown 
624e27a552bSJed Brown   PetscFunctionBegin;
625e27a552bSJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
626e27a552bSJed Brown   if (iascii) {
62761692a83SJed Brown     const TSRosWType rostype;
628e408995aSJed Brown     PetscInt i;
629e408995aSJed Brown     PetscReal abscissa[512];
630e27a552bSJed Brown     char buf[512];
63161692a83SJed Brown     ierr = TSRosWGetType(ts,&rostype);CHKERRQ(ierr);
63261692a83SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Rosenbrock-W %s\n",rostype);CHKERRQ(ierr);
633e408995aSJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ASum);CHKERRQ(ierr);
63461692a83SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A       = %s\n",buf);CHKERRQ(ierr);
635e408995aSJed Brown     for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i];
636e408995aSJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,abscissa);CHKERRQ(ierr);
637e408995aSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A+Gamma = %s\n",buf);CHKERRQ(ierr);
638e27a552bSJed Brown   }
639e27a552bSJed Brown   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
640e27a552bSJed Brown   PetscFunctionReturn(0);
641e27a552bSJed Brown }
642e27a552bSJed Brown 
643e27a552bSJed Brown #undef __FUNCT__
644e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType"
645e27a552bSJed Brown /*@C
64661692a83SJed Brown   TSRosWSetType - Set the type of Rosenbrock-W scheme
647e27a552bSJed Brown 
648e27a552bSJed Brown   Logically collective
649e27a552bSJed Brown 
650e27a552bSJed Brown   Input Parameter:
651e27a552bSJed Brown +  ts - timestepping context
65261692a83SJed Brown -  rostype - type of Rosenbrock-W scheme
653e27a552bSJed Brown 
654e27a552bSJed Brown   Level: intermediate
655e27a552bSJed Brown 
656e27a552bSJed Brown .seealso: TSRosWGetType()
657e27a552bSJed Brown @*/
65861692a83SJed Brown PetscErrorCode TSRosWSetType(TS ts,const TSRosWType rostype)
659e27a552bSJed Brown {
660e27a552bSJed Brown   PetscErrorCode ierr;
661e27a552bSJed Brown 
662e27a552bSJed Brown   PetscFunctionBegin;
663e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
66461692a83SJed Brown   ierr = PetscTryMethod(ts,"TSRosWSetType_C",(TS,const TSRosWType),(ts,rostype));CHKERRQ(ierr);
665e27a552bSJed Brown   PetscFunctionReturn(0);
666e27a552bSJed Brown }
667e27a552bSJed Brown 
668e27a552bSJed Brown #undef __FUNCT__
669e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType"
670e27a552bSJed Brown /*@C
67161692a83SJed Brown   TSRosWGetType - Get the type of Rosenbrock-W scheme
672e27a552bSJed Brown 
673e27a552bSJed Brown   Logically collective
674e27a552bSJed Brown 
675e27a552bSJed Brown   Input Parameter:
676e27a552bSJed Brown .  ts - timestepping context
677e27a552bSJed Brown 
678e27a552bSJed Brown   Output Parameter:
67961692a83SJed Brown .  rostype - type of Rosenbrock-W scheme
680e27a552bSJed Brown 
681e27a552bSJed Brown   Level: intermediate
682e27a552bSJed Brown 
683e27a552bSJed Brown .seealso: TSRosWGetType()
684e27a552bSJed Brown @*/
68561692a83SJed Brown PetscErrorCode TSRosWGetType(TS ts,const TSRosWType *rostype)
686e27a552bSJed Brown {
687e27a552bSJed Brown   PetscErrorCode ierr;
688e27a552bSJed Brown 
689e27a552bSJed Brown   PetscFunctionBegin;
690e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
69161692a83SJed Brown   ierr = PetscUseMethod(ts,"TSRosWGetType_C",(TS,const TSRosWType*),(ts,rostype));CHKERRQ(ierr);
692e27a552bSJed Brown   PetscFunctionReturn(0);
693e27a552bSJed Brown }
694e27a552bSJed Brown 
695e27a552bSJed Brown #undef __FUNCT__
69661692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian"
697e27a552bSJed Brown /*@C
69861692a83SJed Brown   TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step.
699e27a552bSJed Brown 
700e27a552bSJed Brown   Logically collective
701e27a552bSJed Brown 
702e27a552bSJed Brown   Input Parameter:
703e27a552bSJed Brown +  ts - timestepping context
70461692a83SJed Brown -  flg - PETSC_TRUE to recompute the Jacobian at each stage
705e27a552bSJed Brown 
706e27a552bSJed Brown   Level: intermediate
707e27a552bSJed Brown 
708e27a552bSJed Brown .seealso: TSRosWGetType()
709e27a552bSJed Brown @*/
71061692a83SJed Brown PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg)
711e27a552bSJed Brown {
712e27a552bSJed Brown   PetscErrorCode ierr;
713e27a552bSJed Brown 
714e27a552bSJed Brown   PetscFunctionBegin;
715e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
71661692a83SJed Brown   ierr = PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
717e27a552bSJed Brown   PetscFunctionReturn(0);
718e27a552bSJed Brown }
719e27a552bSJed Brown 
720e27a552bSJed Brown EXTERN_C_BEGIN
721e27a552bSJed Brown #undef __FUNCT__
722e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType_RosW"
72361692a83SJed Brown PetscErrorCode  TSRosWGetType_RosW(TS ts,const TSRosWType *rostype)
724e27a552bSJed Brown {
72561692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
726e27a552bSJed Brown   PetscErrorCode ierr;
727e27a552bSJed Brown 
728e27a552bSJed Brown   PetscFunctionBegin;
72961692a83SJed Brown   if (!ros->tableau) {ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);}
73061692a83SJed Brown   *rostype = ros->tableau->name;
731e27a552bSJed Brown   PetscFunctionReturn(0);
732e27a552bSJed Brown }
733e27a552bSJed Brown #undef __FUNCT__
734e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType_RosW"
73561692a83SJed Brown PetscErrorCode  TSRosWSetType_RosW(TS ts,const TSRosWType rostype)
736e27a552bSJed Brown {
73761692a83SJed Brown   TS_RosW         *ros = (TS_RosW*)ts->data;
738e27a552bSJed Brown   PetscErrorCode  ierr;
739e27a552bSJed Brown   PetscBool       match;
74061692a83SJed Brown   RosWTableauLink link;
741e27a552bSJed Brown 
742e27a552bSJed Brown   PetscFunctionBegin;
74361692a83SJed Brown   if (ros->tableau) {
74461692a83SJed Brown     ierr = PetscStrcmp(ros->tableau->name,rostype,&match);CHKERRQ(ierr);
745e27a552bSJed Brown     if (match) PetscFunctionReturn(0);
746e27a552bSJed Brown   }
74761692a83SJed Brown   for (link = RosWTableauList; link; link=link->next) {
74861692a83SJed Brown     ierr = PetscStrcmp(link->tab.name,rostype,&match);CHKERRQ(ierr);
749e27a552bSJed Brown     if (match) {
750e27a552bSJed Brown       ierr = TSReset_RosW(ts);CHKERRQ(ierr);
75161692a83SJed Brown       ros->tableau = &link->tab;
752e27a552bSJed Brown       PetscFunctionReturn(0);
753e27a552bSJed Brown     }
754e27a552bSJed Brown   }
75561692a83SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype);
756e27a552bSJed Brown   PetscFunctionReturn(0);
757e27a552bSJed Brown }
75861692a83SJed Brown 
759e27a552bSJed Brown #undef __FUNCT__
76061692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian_RosW"
76161692a83SJed Brown PetscErrorCode  TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg)
762e27a552bSJed Brown {
76361692a83SJed Brown   TS_RosW *ros = (TS_RosW*)ts->data;
764e27a552bSJed Brown 
765e27a552bSJed Brown   PetscFunctionBegin;
76661692a83SJed Brown   ros->recompute_jacobian = flg;
767e27a552bSJed Brown   PetscFunctionReturn(0);
768e27a552bSJed Brown }
769e27a552bSJed Brown EXTERN_C_END
770e27a552bSJed Brown 
771e27a552bSJed Brown /* ------------------------------------------------------------ */
772e27a552bSJed Brown /*MC
773e27a552bSJed Brown       TSRosW - ODE solver using Rosenbrock-W schemes
774e27a552bSJed Brown 
775e27a552bSJed Brown   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
776e27a552bSJed Brown   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
777e27a552bSJed Brown   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
778e27a552bSJed Brown 
779e27a552bSJed Brown   Notes:
78061692a83SJed Brown   This method currently only works with autonomous ODE and DAE.
78161692a83SJed Brown 
78261692a83SJed Brown   Developer notes:
78361692a83SJed Brown   Rosenbrock-W methods are typically specified for autonomous ODE
78461692a83SJed Brown 
78561692a83SJed Brown $  xdot = f(x)
78661692a83SJed Brown 
78761692a83SJed Brown   by the stage equations
78861692a83SJed Brown 
78961692a83SJed Brown $  k_i = h f(x_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j
79061692a83SJed Brown 
79161692a83SJed Brown   and step completion formula
79261692a83SJed Brown 
79361692a83SJed Brown $  x_1 = x_0 + sum_j b_j k_j
79461692a83SJed Brown 
79561692a83SJed Brown   with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(x)
79661692a83SJed Brown   and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner,
79761692a83SJed Brown   we define new variables for the stage equations
79861692a83SJed Brown 
79961692a83SJed Brown $  y_i = gamma_ij k_j
80061692a83SJed Brown 
80161692a83SJed Brown   The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define
80261692a83SJed Brown 
80361692a83SJed Brown $  A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-i}
80461692a83SJed Brown 
80561692a83SJed Brown   to rewrite the method as
80661692a83SJed Brown 
80761692a83SJed 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
80861692a83SJed Brown $  x_1 = x_0 + sum_j bt_j y_j
80961692a83SJed Brown 
81061692a83SJed Brown    where we have introduced the mass matrix M. Continue by defining
81161692a83SJed Brown 
81261692a83SJed Brown $  ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j
81361692a83SJed Brown 
81461692a83SJed Brown    or, more compactly in tensor notation
81561692a83SJed Brown 
81661692a83SJed Brown $  Ydot = 1/h (Gamma^{-1} \otimes I) Y .
81761692a83SJed Brown 
81861692a83SJed Brown    Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current
81961692a83SJed Brown    stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the
82061692a83SJed Brown    equation
82161692a83SJed Brown 
82261692a83SJed Brown $  g(x_0 + sum_j a_ij y_j + y_i, ydot_i) = 0
82361692a83SJed Brown 
82461692a83SJed Brown    with initial guess y_i = 0.
825e27a552bSJed Brown 
826e27a552bSJed Brown   Level: beginner
827e27a552bSJed Brown 
828e27a552bSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSRosWRegister()
829e27a552bSJed Brown 
830e27a552bSJed Brown M*/
831e27a552bSJed Brown EXTERN_C_BEGIN
832e27a552bSJed Brown #undef __FUNCT__
833e27a552bSJed Brown #define __FUNCT__ "TSCreate_RosW"
834e27a552bSJed Brown PetscErrorCode  TSCreate_RosW(TS ts)
835e27a552bSJed Brown {
83661692a83SJed Brown   TS_RosW        *ros;
837e27a552bSJed Brown   PetscErrorCode ierr;
838e27a552bSJed Brown 
839e27a552bSJed Brown   PetscFunctionBegin;
840e27a552bSJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
841e27a552bSJed Brown   ierr = TSRosWInitializePackage(PETSC_NULL);CHKERRQ(ierr);
842e27a552bSJed Brown #endif
843e27a552bSJed Brown 
844e27a552bSJed Brown   ts->ops->reset          = TSReset_RosW;
845e27a552bSJed Brown   ts->ops->destroy        = TSDestroy_RosW;
846e27a552bSJed Brown   ts->ops->view           = TSView_RosW;
847e27a552bSJed Brown   ts->ops->setup          = TSSetUp_RosW;
848e27a552bSJed Brown   ts->ops->step           = TSStep_RosW;
849e27a552bSJed Brown   ts->ops->interpolate    = TSInterpolate_RosW;
850e27a552bSJed Brown   ts->ops->setfromoptions = TSSetFromOptions_RosW;
851e27a552bSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_RosW;
852e27a552bSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_RosW;
853e27a552bSJed Brown 
85461692a83SJed Brown   ierr = PetscNewLog(ts,TS_RosW,&ros);CHKERRQ(ierr);
85561692a83SJed Brown   ts->data = (void*)ros;
856e27a552bSJed Brown 
857e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","TSRosWGetType_RosW",TSRosWGetType_RosW);CHKERRQ(ierr);
858e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","TSRosWSetType_RosW",TSRosWSetType_RosW);CHKERRQ(ierr);
85961692a83SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","TSRosWSetRecomputeJacobian_RosW",TSRosWSetRecomputeJacobian_RosW);CHKERRQ(ierr);
860e27a552bSJed Brown   PetscFunctionReturn(0);
861e27a552bSJed Brown }
862e27a552bSJed Brown EXTERN_C_END
863