xref: /petsc/src/ts/impls/rosw/rosw.c (revision 61692a839f98c39afd201bc4752f07afc59b0300)
1e27a552bSJed Brown /*
2*61692a83SJed Brown   Code for timestepping with Rosenbrock W methods
3e27a552bSJed Brown 
4e27a552bSJed Brown   Notes:
5e27a552bSJed Brown   The general system is written as
6e27a552bSJed Brown 
7*61692a83SJed Brown   G(t,X,Xdot) = F(t,X)
8e27a552bSJed Brown 
9*61692a83SJed Brown   where G represents the stiff part of the physics and F represents the non-stiff part.
10*61692a83SJed 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 
15*61692a83SJed Brown #include <../src/mat/blockinvert.h>
16*61692a83SJed Brown 
17*61692a83SJed Brown static const TSRosWType TSRosWDefault = TSROSW2P;
18e27a552bSJed Brown static PetscBool TSRosWRegisterAllCalled;
19e27a552bSJed Brown static PetscBool TSRosWPackageInitialized;
20e27a552bSJed Brown 
21*61692a83SJed Brown typedef struct _RosWTableau *RosWTableau;
22*61692a83SJed Brown struct _RosWTableau {
23e27a552bSJed Brown   char *name;
24e27a552bSJed Brown   PetscInt order;               /* Classical approximation order of the method */
25e27a552bSJed Brown   PetscInt s;                   /* Number of stages */
26e27a552bSJed Brown   PetscInt pinterp;             /* Interpolation order */
27*61692a83SJed Brown   PetscReal *A;                 /* Propagation table, strictly lower triangular */
28*61692a83SJed Brown   PetscReal *Gamma;             /* Stage table, lower triangular with nonzero diagonal */
29*61692a83SJed Brown   PetscReal *b;                 /* Step completion table */
30*61692a83SJed Brown   PetscReal *ASum;              /* Row sum of A */
31*61692a83SJed Brown   PetscReal *GammaSum;          /* Row sum of Gamma, only needed for non-autonomous systems */
32*61692a83SJed Brown   PetscReal *At;                /* Propagation table in transformed variables */
33*61692a83SJed Brown   PetscReal *bt;                /* Step completion table in transformed variables  */
34*61692a83SJed Brown   PetscReal *GammaInv;          /* Inverse of Gamma, used for transformed variables */
35e27a552bSJed Brown };
36*61692a83SJed Brown typedef struct _RosWTableauLink *RosWTableauLink;
37*61692a83SJed Brown struct _RosWTableauLink {
38*61692a83SJed Brown   struct _RosWTableau tab;
39*61692a83SJed Brown   RosWTableauLink next;
40e27a552bSJed Brown };
41*61692a83SJed Brown static RosWTableauLink RosWTableauList;
42e27a552bSJed Brown 
43e27a552bSJed Brown typedef struct {
44*61692a83SJed Brown   RosWTableau tableau;
45*61692a83SJed Brown   Vec         *Y;               /* States computed during the step, used to complete the step */
46e27a552bSJed Brown   Vec         Ydot;             /* Work vector holding Ydot during residual evaluation */
47*61692a83SJed Brown   Vec         Ystage;           /* Work vector for the state value at each stage */
48*61692a83SJed Brown   Vec         Zdot;             /* Ydot = Zdot + shift*Y */
49*61692a83SJed Brown   Vec         Zstage;           /* Y = Zstage + Y */
50e27a552bSJed Brown   PetscScalar *work;            /* Scalar work */
51e27a552bSJed Brown   PetscReal   shift;
52e27a552bSJed Brown   PetscReal   stage_time;
53*61692a83SJed Brown   PetscBool   recompute_jacobian; /* Recompute the Jacobian at each stage, default is to freeze the Jacobian at the start of each step */
54e27a552bSJed Brown } TS_RosW;
55e27a552bSJed Brown 
56e27a552bSJed Brown #undef __FUNCT__
57e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterAll"
58e27a552bSJed Brown /*@C
59e27a552bSJed Brown   TSRosWRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSRosW
60e27a552bSJed Brown 
61e27a552bSJed Brown   Not Collective, but should be called by all processes which will need the schemes to be registered
62e27a552bSJed Brown 
63e27a552bSJed Brown   Level: advanced
64e27a552bSJed Brown 
65e27a552bSJed Brown .keywords: TS, TSRosW, register, all
66e27a552bSJed Brown 
67e27a552bSJed Brown .seealso:  TSRosWRegisterDestroy()
68e27a552bSJed Brown @*/
69e27a552bSJed Brown PetscErrorCode TSRosWRegisterAll(void)
70e27a552bSJed Brown {
71e27a552bSJed Brown   PetscErrorCode ierr;
72e27a552bSJed Brown 
73e27a552bSJed Brown   PetscFunctionBegin;
74e27a552bSJed Brown   if (TSRosWRegisterAllCalled) PetscFunctionReturn(0);
75e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_TRUE;
76e27a552bSJed Brown 
77e27a552bSJed Brown   {
78*61692a83SJed Brown     const PetscReal g = 1. + 1./PetscSqrtReal(2.0);
79e27a552bSJed Brown     const PetscReal
80*61692a83SJed Brown       A[2][2] = {{0,0}, {1.,0}},
81*61692a83SJed Brown       Gamma[2][2] = {{g,0}, {-2.*g,g}},
82*61692a83SJed Brown       b[2] = {0.5,0.5};
83*61692a83SJed Brown     ierr = TSRosWRegister(TSROSW2P,2,2,&A[0][0],&Gamma[0][0],b);CHKERRQ(ierr);
84e27a552bSJed Brown   }
85e27a552bSJed Brown   {
86*61692a83SJed Brown     const PetscReal g = 1. - 1./PetscSqrtReal(2.0);
87e27a552bSJed Brown     const PetscReal
88*61692a83SJed Brown       A[2][2] = {{0,0}, {1.,0}},
89*61692a83SJed Brown       Gamma[2][2] = {{g,0}, {-2.*g,g}},
90*61692a83SJed Brown       b[2] = {0.5,0.5};
91*61692a83SJed Brown     ierr = TSRosWRegister(TSROSW2M,2,2,&A[0][0],&Gamma[0][0],b);CHKERRQ(ierr);
92e27a552bSJed Brown   }
93e27a552bSJed Brown   PetscFunctionReturn(0);
94e27a552bSJed Brown }
95e27a552bSJed Brown 
96e27a552bSJed Brown #undef __FUNCT__
97e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterDestroy"
98e27a552bSJed Brown /*@C
99e27a552bSJed Brown    TSRosWRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister().
100e27a552bSJed Brown 
101e27a552bSJed Brown    Not Collective
102e27a552bSJed Brown 
103e27a552bSJed Brown    Level: advanced
104e27a552bSJed Brown 
105e27a552bSJed Brown .keywords: TSRosW, register, destroy
106e27a552bSJed Brown .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegisterDynamic()
107e27a552bSJed Brown @*/
108e27a552bSJed Brown PetscErrorCode TSRosWRegisterDestroy(void)
109e27a552bSJed Brown {
110e27a552bSJed Brown   PetscErrorCode ierr;
111*61692a83SJed Brown   RosWTableauLink link;
112e27a552bSJed Brown 
113e27a552bSJed Brown   PetscFunctionBegin;
114*61692a83SJed Brown   while ((link = RosWTableauList)) {
115*61692a83SJed Brown     RosWTableau t = &link->tab;
116*61692a83SJed Brown     RosWTableauList = link->next;
117*61692a83SJed Brown     ierr = PetscFree5(t->A,t->Gamma,t->b,t->ASum,t->GammaSum);CHKERRQ(ierr);
118*61692a83SJed Brown     ierr = PetscFree3(t->At,t->bt,t->GammaInv);CHKERRQ(ierr);
119e27a552bSJed Brown     ierr = PetscFree(t->name);CHKERRQ(ierr);
120e27a552bSJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
121e27a552bSJed Brown   }
122e27a552bSJed Brown   TSRosWRegisterAllCalled = PETSC_FALSE;
123e27a552bSJed Brown   PetscFunctionReturn(0);
124e27a552bSJed Brown }
125e27a552bSJed Brown 
126e27a552bSJed Brown #undef __FUNCT__
127e27a552bSJed Brown #define __FUNCT__ "TSRosWInitializePackage"
128e27a552bSJed Brown /*@C
129e27a552bSJed Brown   TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called
130e27a552bSJed Brown   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RosW()
131e27a552bSJed Brown   when using static libraries.
132e27a552bSJed Brown 
133e27a552bSJed Brown   Input Parameter:
134e27a552bSJed Brown   path - The dynamic library path, or PETSC_NULL
135e27a552bSJed Brown 
136e27a552bSJed Brown   Level: developer
137e27a552bSJed Brown 
138e27a552bSJed Brown .keywords: TS, TSRosW, initialize, package
139e27a552bSJed Brown .seealso: PetscInitialize()
140e27a552bSJed Brown @*/
141e27a552bSJed Brown PetscErrorCode TSRosWInitializePackage(const char path[])
142e27a552bSJed Brown {
143e27a552bSJed Brown   PetscErrorCode ierr;
144e27a552bSJed Brown 
145e27a552bSJed Brown   PetscFunctionBegin;
146e27a552bSJed Brown   if (TSRosWPackageInitialized) PetscFunctionReturn(0);
147e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_TRUE;
148e27a552bSJed Brown   ierr = TSRosWRegisterAll();CHKERRQ(ierr);
149e27a552bSJed Brown   ierr = PetscRegisterFinalize(TSRosWFinalizePackage);CHKERRQ(ierr);
150e27a552bSJed Brown   PetscFunctionReturn(0);
151e27a552bSJed Brown }
152e27a552bSJed Brown 
153e27a552bSJed Brown #undef __FUNCT__
154e27a552bSJed Brown #define __FUNCT__ "TSRosWFinalizePackage"
155e27a552bSJed Brown /*@C
156e27a552bSJed Brown   TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is
157e27a552bSJed Brown   called from PetscFinalize().
158e27a552bSJed Brown 
159e27a552bSJed Brown   Level: developer
160e27a552bSJed Brown 
161e27a552bSJed Brown .keywords: Petsc, destroy, package
162e27a552bSJed Brown .seealso: PetscFinalize()
163e27a552bSJed Brown @*/
164e27a552bSJed Brown PetscErrorCode TSRosWFinalizePackage(void)
165e27a552bSJed Brown {
166e27a552bSJed Brown   PetscErrorCode ierr;
167e27a552bSJed Brown 
168e27a552bSJed Brown   PetscFunctionBegin;
169e27a552bSJed Brown   TSRosWPackageInitialized = PETSC_FALSE;
170e27a552bSJed Brown   ierr = TSRosWRegisterDestroy();CHKERRQ(ierr);
171e27a552bSJed Brown   PetscFunctionReturn(0);
172e27a552bSJed Brown }
173e27a552bSJed Brown 
174e27a552bSJed Brown #undef __FUNCT__
175e27a552bSJed Brown #define __FUNCT__ "TSRosWRegister"
176e27a552bSJed Brown /*@C
177*61692a83SJed Brown    TSRosWRegister - register a Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
178e27a552bSJed Brown 
179e27a552bSJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
180e27a552bSJed Brown 
181e27a552bSJed Brown    Input Parameters:
182e27a552bSJed Brown +  name - identifier for method
183e27a552bSJed Brown .  order - approximation order of method
184e27a552bSJed Brown .  s - number of stages, this is the dimension of the matrices below
185*61692a83SJed Brown .  A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular
186*61692a83SJed Brown .  Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal
187*61692a83SJed Brown -  b - Step completion table (dimension s)
188e27a552bSJed Brown 
189e27a552bSJed Brown    Notes:
190*61692a83SJed Brown    Several Rosenbrock W methods are provided, this function is only needed to create new methods.
191e27a552bSJed Brown 
192e27a552bSJed Brown    Level: advanced
193e27a552bSJed Brown 
194e27a552bSJed Brown .keywords: TS, register
195e27a552bSJed Brown 
196e27a552bSJed Brown .seealso: TSRosW
197e27a552bSJed Brown @*/
198e27a552bSJed Brown PetscErrorCode TSRosWRegister(const TSRosWType name,PetscInt order,PetscInt s,
199*61692a83SJed Brown                               const PetscReal A[],const PetscReal Gamma[],const PetscReal b[])
200e27a552bSJed Brown {
201e27a552bSJed Brown   PetscErrorCode ierr;
202*61692a83SJed Brown   RosWTableauLink link;
203*61692a83SJed Brown   RosWTableau t;
204*61692a83SJed Brown   PetscInt i,j,k;
205*61692a83SJed Brown   PetscScalar *GammaInv;
206e27a552bSJed Brown 
207e27a552bSJed Brown   PetscFunctionBegin;
208e27a552bSJed Brown   ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
209e27a552bSJed Brown   ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
210e27a552bSJed Brown   t = &link->tab;
211e27a552bSJed Brown   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
212e27a552bSJed Brown   t->order = order;
213e27a552bSJed Brown   t->s = s;
214*61692a83SJed 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);
215*61692a83SJed Brown   ierr = PetscMalloc3(s*s,PetscReal,&t->At,s,PetscReal,&t->bt,s*s,PetscReal,&t->GammaInv);CHKERRQ(ierr);
216e27a552bSJed Brown   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
217*61692a83SJed Brown   ierr = PetscMemcpy(t->Gamma,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr);
218*61692a83SJed Brown   ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
219*61692a83SJed Brown   for (i=0; i<s; i++) {
220*61692a83SJed Brown     t->ASum[i] = 0;
221*61692a83SJed Brown     t->GammaSum[i] = 0;
222*61692a83SJed Brown     for (j=0; j<s; j++) {
223*61692a83SJed Brown       t->ASum[i] += A[i*s+j];
224*61692a83SJed Brown       t->GammaSum[i] = Gamma[i*s+j];
225*61692a83SJed Brown     }
226*61692a83SJed Brown   }
227*61692a83SJed Brown   ierr = PetscMalloc(s*s*sizeof(PetscScalar),&GammaInv);CHKERRQ(ierr); /* Need to use Scalar for inverse, then convert back to Real */
228*61692a83SJed Brown   for (i=0; i<s*s; i++) GammaInv[i] = Gamma[i];
229*61692a83SJed Brown   switch (s) {
230*61692a83SJed Brown   case 1: GammaInv[0] = 1./GammaInv[0]; break;
231*61692a83SJed Brown   case 2: ierr = Kernel_A_gets_inverse_A_2(GammaInv,0);CHKERRQ(ierr); break;
232*61692a83SJed Brown   case 3: ierr = Kernel_A_gets_inverse_A_3(GammaInv,0);CHKERRQ(ierr); break;
233*61692a83SJed Brown   case 4: ierr = Kernel_A_gets_inverse_A_4(GammaInv,0);CHKERRQ(ierr); break;
234*61692a83SJed Brown   case 5: {
235*61692a83SJed Brown     PetscInt ipvt5[5];
236*61692a83SJed Brown     MatScalar work5[5*5];
237*61692a83SJed Brown     ierr = Kernel_A_gets_inverse_A_5(GammaInv,ipvt5,work5,0);CHKERRQ(ierr); break;
238*61692a83SJed Brown   }
239*61692a83SJed Brown   case 6: ierr = Kernel_A_gets_inverse_A_6(GammaInv,0);CHKERRQ(ierr); break;
240*61692a83SJed Brown   case 7: ierr = Kernel_A_gets_inverse_A_7(GammaInv,0);CHKERRQ(ierr); break;
241*61692a83SJed Brown   default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for %D stages",s);
242*61692a83SJed Brown   }
243*61692a83SJed Brown   for (i=0; i<s*s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]);
244*61692a83SJed Brown   ierr = PetscFree(GammaInv);CHKERRQ(ierr);
245*61692a83SJed Brown   for (i=0; i<s; i++) {
246*61692a83SJed Brown     for (j=0; j<s; j++) {
247*61692a83SJed Brown       t->At[i*s+j] = 0;
248*61692a83SJed Brown       for (k=0; k<s; k++) {
249*61692a83SJed Brown         t->At[i*s+j] += t->A[i*s+k] * t->GammaInv[k*s+j];
250*61692a83SJed Brown       }
251*61692a83SJed Brown     }
252*61692a83SJed Brown     t->bt[i] = 0;
253*61692a83SJed Brown     for (j=0; j<s; j++) {
254*61692a83SJed Brown       t->bt[i] += t->b[j] * t->GammaInv[j*s+i];
255*61692a83SJed Brown     }
256*61692a83SJed Brown   }
257*61692a83SJed Brown   link->next = RosWTableauList;
258*61692a83SJed Brown   RosWTableauList = link;
259e27a552bSJed Brown   PetscFunctionReturn(0);
260e27a552bSJed Brown }
261e27a552bSJed Brown 
262e27a552bSJed Brown #undef __FUNCT__
263e27a552bSJed Brown #define __FUNCT__ "TSStep_RosW"
264e27a552bSJed Brown static PetscErrorCode TSStep_RosW(TS ts)
265e27a552bSJed Brown {
266*61692a83SJed Brown   TS_RosW         *ros = (TS_RosW*)ts->data;
267*61692a83SJed Brown   RosWTableau     tab  = ros->tableau;
268e27a552bSJed Brown   const PetscInt  s    = tab->s;
269*61692a83SJed Brown   const PetscReal *At  = tab->At,*Gamma = tab->Gamma,*bt = tab->bt,*ASum = tab->ASum,*GammaInv = tab->GammaInv;
270*61692a83SJed Brown   PetscScalar     *w   = ros->work;
271*61692a83SJed Brown   Vec             *Y   = ros->Y,Zdot = ros->Zdot,Zstage = ros->Zstage;
272e27a552bSJed Brown   SNES            snes;
273e27a552bSJed Brown   PetscInt        i,j,its,lits;
274cdbf8f93SLisandro Dalcin   PetscReal       next_time_step;
275e27a552bSJed Brown   PetscReal       h,t;
276e27a552bSJed Brown   PetscErrorCode  ierr;
277e27a552bSJed Brown 
278e27a552bSJed Brown   PetscFunctionBegin;
279e27a552bSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
280cdbf8f93SLisandro Dalcin   next_time_step = ts->time_step;
281cdbf8f93SLisandro Dalcin   h = ts->time_step;
282e27a552bSJed Brown   t = ts->ptime;
283e27a552bSJed Brown 
284e27a552bSJed Brown   for (i=0; i<s; i++) {
285*61692a83SJed Brown     ros->stage_time = t + h*ASum[i];
286*61692a83SJed Brown     ros->shift = 1./(h*Gamma[i*s+i]);
287*61692a83SJed Brown 
288*61692a83SJed Brown     ierr = VecCopy(ts->vec_sol,Zstage);CHKERRQ(ierr);
289*61692a83SJed Brown     ierr = VecMAXPY(Zstage,i,&At[i*s+0],Y);CHKERRQ(ierr);
290*61692a83SJed Brown 
291*61692a83SJed Brown     for (j=0; j<i; j++) w[j] = 1./h * GammaInv[i*s+j];
292*61692a83SJed Brown     ierr = VecZeroEntries(Zdot);CHKERRQ(ierr);
293*61692a83SJed Brown     ierr = VecMAXPY(Zdot,i,w,Y);CHKERRQ(ierr);
294*61692a83SJed Brown 
295e27a552bSJed Brown     /* Initial guess taken from last stage */
296*61692a83SJed Brown     ierr = VecZeroEntries(Y[i]);CHKERRQ(ierr);
297*61692a83SJed Brown 
298*61692a83SJed Brown     if (!ros->recompute_jacobian && !i) {
299*61692a83SJed Brown       ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */
300*61692a83SJed Brown     }
301*61692a83SJed Brown 
302*61692a83SJed Brown     ierr = SNESSolve(snes,PETSC_NULL,Y[i]);CHKERRQ(ierr);
303e27a552bSJed Brown     ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
304e27a552bSJed Brown     ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
305e27a552bSJed Brown     ts->nonlinear_its += its; ts->linear_its += lits;
306e27a552bSJed Brown   }
307*61692a83SJed Brown   ierr = VecMAXPY(ts->vec_sol,s,bt,Y);CHKERRQ(ierr);
308e27a552bSJed Brown 
309e27a552bSJed Brown   ts->ptime += ts->time_step;
310cdbf8f93SLisandro Dalcin   ts->time_step = next_time_step;
311e27a552bSJed Brown   ts->steps++;
312e27a552bSJed Brown   PetscFunctionReturn(0);
313e27a552bSJed Brown }
314e27a552bSJed Brown 
315e27a552bSJed Brown #undef __FUNCT__
316e27a552bSJed Brown #define __FUNCT__ "TSInterpolate_RosW"
317e27a552bSJed Brown static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec X)
318e27a552bSJed Brown {
319*61692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
320e27a552bSJed Brown 
321e27a552bSJed Brown   PetscFunctionBegin;
322*61692a83SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name);
323e27a552bSJed Brown   PetscFunctionReturn(0);
324e27a552bSJed Brown }
325e27a552bSJed Brown 
326e27a552bSJed Brown /*------------------------------------------------------------*/
327e27a552bSJed Brown #undef __FUNCT__
328e27a552bSJed Brown #define __FUNCT__ "TSReset_RosW"
329e27a552bSJed Brown static PetscErrorCode TSReset_RosW(TS ts)
330e27a552bSJed Brown {
331*61692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
332e27a552bSJed Brown   PetscInt       s;
333e27a552bSJed Brown   PetscErrorCode ierr;
334e27a552bSJed Brown 
335e27a552bSJed Brown   PetscFunctionBegin;
336*61692a83SJed Brown   if (!ros->tableau) PetscFunctionReturn(0);
337*61692a83SJed Brown    s = ros->tableau->s;
338*61692a83SJed Brown   ierr = VecDestroyVecs(s,&ros->Y);CHKERRQ(ierr);
339*61692a83SJed Brown   ierr = VecDestroy(&ros->Ydot);CHKERRQ(ierr);
340*61692a83SJed Brown   ierr = VecDestroy(&ros->Ystage);CHKERRQ(ierr);
341*61692a83SJed Brown   ierr = VecDestroy(&ros->Zdot);CHKERRQ(ierr);
342*61692a83SJed Brown   ierr = VecDestroy(&ros->Zstage);CHKERRQ(ierr);
343*61692a83SJed Brown   ierr = PetscFree(ros->work);CHKERRQ(ierr);
344e27a552bSJed Brown   PetscFunctionReturn(0);
345e27a552bSJed Brown }
346e27a552bSJed Brown 
347e27a552bSJed Brown #undef __FUNCT__
348e27a552bSJed Brown #define __FUNCT__ "TSDestroy_RosW"
349e27a552bSJed Brown static PetscErrorCode TSDestroy_RosW(TS ts)
350e27a552bSJed Brown {
351e27a552bSJed Brown   PetscErrorCode  ierr;
352e27a552bSJed Brown 
353e27a552bSJed Brown   PetscFunctionBegin;
354e27a552bSJed Brown   ierr = TSReset_RosW(ts);CHKERRQ(ierr);
355e27a552bSJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
356e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","",PETSC_NULL);CHKERRQ(ierr);
357e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","",PETSC_NULL);CHKERRQ(ierr);
358*61692a83SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","",PETSC_NULL);CHKERRQ(ierr);
359e27a552bSJed Brown   PetscFunctionReturn(0);
360e27a552bSJed Brown }
361e27a552bSJed Brown 
362e27a552bSJed Brown /*
363e27a552bSJed Brown   This defines the nonlinear equation that is to be solved with SNES
364e27a552bSJed Brown   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
365e27a552bSJed Brown */
366e27a552bSJed Brown #undef __FUNCT__
367e27a552bSJed Brown #define __FUNCT__ "SNESTSFormFunction_RosW"
368e27a552bSJed Brown static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec X,Vec F,TS ts)
369e27a552bSJed Brown {
370*61692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
371e27a552bSJed Brown   PetscErrorCode ierr;
372e27a552bSJed Brown 
373e27a552bSJed Brown   PetscFunctionBegin;
374*61692a83SJed Brown   ierr = VecWAXPY(ros->Ydot,ros->shift,X,ros->Zdot);CHKERRQ(ierr); /* Ydot = shift*X + Zdot */
375*61692a83SJed Brown   ierr = VecWAXPY(ros->Ystage,1.0,X,ros->Zstage);CHKERRQ(ierr);    /* Ystage = X + Zstage */
376*61692a83SJed Brown   ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,ros->Ydot,F,PETSC_FALSE);CHKERRQ(ierr);
377e27a552bSJed Brown   PetscFunctionReturn(0);
378e27a552bSJed Brown }
379e27a552bSJed Brown 
380e27a552bSJed Brown #undef __FUNCT__
381e27a552bSJed Brown #define __FUNCT__ "SNESTSFormJacobian_RosW"
382e27a552bSJed Brown static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
383e27a552bSJed Brown {
384*61692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
385e27a552bSJed Brown   PetscErrorCode ierr;
386e27a552bSJed Brown 
387e27a552bSJed Brown   PetscFunctionBegin;
388*61692a83SJed Brown   /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
389*61692a83SJed Brown   ierr = TSComputeIJacobian(ts,ros->stage_time,ros->Ystage,ros->Ydot,ros->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr);
390e27a552bSJed Brown   PetscFunctionReturn(0);
391e27a552bSJed Brown }
392e27a552bSJed Brown 
393e27a552bSJed Brown #undef __FUNCT__
394e27a552bSJed Brown #define __FUNCT__ "TSSetUp_RosW"
395e27a552bSJed Brown static PetscErrorCode TSSetUp_RosW(TS ts)
396e27a552bSJed Brown {
397*61692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
398*61692a83SJed Brown   RosWTableau    tab  = ros->tableau;
399e27a552bSJed Brown   PetscInt       s    = tab->s;
400e27a552bSJed Brown   PetscErrorCode ierr;
401e27a552bSJed Brown 
402e27a552bSJed Brown   PetscFunctionBegin;
403*61692a83SJed Brown   if (!ros->tableau) {
404e27a552bSJed Brown     ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);
405e27a552bSJed Brown   }
406*61692a83SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ros->Y);CHKERRQ(ierr);
407*61692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Ydot);CHKERRQ(ierr);
408*61692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Ystage);CHKERRQ(ierr);
409*61692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Zdot);CHKERRQ(ierr);
410*61692a83SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ros->Zstage);CHKERRQ(ierr);
411*61692a83SJed Brown   ierr = PetscMalloc(s*sizeof(ros->work[0]),&ros->work);CHKERRQ(ierr);
412e27a552bSJed Brown   PetscFunctionReturn(0);
413e27a552bSJed Brown }
414e27a552bSJed Brown /*------------------------------------------------------------*/
415e27a552bSJed Brown 
416e27a552bSJed Brown #undef __FUNCT__
417e27a552bSJed Brown #define __FUNCT__ "TSSetFromOptions_RosW"
418e27a552bSJed Brown static PetscErrorCode TSSetFromOptions_RosW(TS ts)
419e27a552bSJed Brown {
420*61692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
421e27a552bSJed Brown   PetscErrorCode ierr;
422*61692a83SJed Brown   char           rostype[256];
423e27a552bSJed Brown 
424e27a552bSJed Brown   PetscFunctionBegin;
425e27a552bSJed Brown   ierr = PetscOptionsHead("RosW ODE solver options");CHKERRQ(ierr);
426e27a552bSJed Brown   {
427*61692a83SJed Brown     RosWTableauLink link;
428e27a552bSJed Brown     PetscInt count,choice;
429e27a552bSJed Brown     PetscBool flg;
430e27a552bSJed Brown     const char **namelist;
431*61692a83SJed Brown     SNES snes;
432*61692a83SJed Brown 
433*61692a83SJed Brown     ierr = PetscStrncpy(rostype,TSRosWDefault,sizeof rostype);CHKERRQ(ierr);
434*61692a83SJed Brown     for (link=RosWTableauList,count=0; link; link=link->next,count++) ;
435e27a552bSJed Brown     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
436*61692a83SJed Brown     for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
437*61692a83SJed Brown     ierr = PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,rostype,&choice,&flg);CHKERRQ(ierr);
438*61692a83SJed Brown     ierr = TSRosWSetType(ts,flg ? namelist[choice] : rostype);CHKERRQ(ierr);
439e27a552bSJed Brown     ierr = PetscFree(namelist);CHKERRQ(ierr);
440*61692a83SJed Brown 
441*61692a83SJed Brown     ierr = PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,PETSC_NULL);CHKERRQ(ierr);
442*61692a83SJed Brown 
443*61692a83SJed Brown     /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
444*61692a83SJed Brown     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
445*61692a83SJed Brown     if (!((PetscObject)snes)->type_name) {
446*61692a83SJed Brown       ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);
447*61692a83SJed Brown     }
448*61692a83SJed Brown     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
449e27a552bSJed Brown   }
450e27a552bSJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
451e27a552bSJed Brown   PetscFunctionReturn(0);
452e27a552bSJed Brown }
453e27a552bSJed Brown 
454e27a552bSJed Brown #undef __FUNCT__
455e27a552bSJed Brown #define __FUNCT__ "PetscFormatRealArray"
456e27a552bSJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
457e27a552bSJed Brown {
458e27a552bSJed Brown   PetscErrorCode ierr;
459e27a552bSJed Brown   int i,left,count;
460e27a552bSJed Brown   char *p;
461e27a552bSJed Brown 
462e27a552bSJed Brown   PetscFunctionBegin;
463e27a552bSJed Brown   for (i=0,p=buf,left=(int)len; i<n; i++) {
464e27a552bSJed Brown     ierr = PetscSNPrintf(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
465e27a552bSJed Brown     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
466e27a552bSJed Brown     left -= count;
467e27a552bSJed Brown     p += count;
468e27a552bSJed Brown     *p++ = ' ';
469e27a552bSJed Brown   }
470e27a552bSJed Brown   p[i ? 0 : -1] = 0;
471e27a552bSJed Brown   PetscFunctionReturn(0);
472e27a552bSJed Brown }
473e27a552bSJed Brown 
474e27a552bSJed Brown #undef __FUNCT__
475e27a552bSJed Brown #define __FUNCT__ "TSView_RosW"
476e27a552bSJed Brown static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer)
477e27a552bSJed Brown {
478*61692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
479*61692a83SJed Brown   RosWTableau    tab  = ros->tableau;
480e27a552bSJed Brown   PetscBool      iascii;
481e27a552bSJed Brown   PetscErrorCode ierr;
482e27a552bSJed Brown 
483e27a552bSJed Brown   PetscFunctionBegin;
484e27a552bSJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
485e27a552bSJed Brown   if (iascii) {
486*61692a83SJed Brown     const TSRosWType rostype;
487e27a552bSJed Brown     char buf[512];
488*61692a83SJed Brown     ierr = TSRosWGetType(ts,&rostype);CHKERRQ(ierr);
489*61692a83SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Rosenbrock-W %s\n",rostype);CHKERRQ(ierr);
490*61692a83SJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6F",tab->s,tab->ASum);CHKERRQ(ierr);
491*61692a83SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of A     = %s\n",buf);CHKERRQ(ierr);
492*61692a83SJed Brown     ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6F",tab->s,tab->GammaSum);CHKERRQ(ierr);
493*61692a83SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa of Gamma = %s\n",buf);CHKERRQ(ierr);
494e27a552bSJed Brown   }
495e27a552bSJed Brown   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
496e27a552bSJed Brown   PetscFunctionReturn(0);
497e27a552bSJed Brown }
498e27a552bSJed Brown 
499e27a552bSJed Brown #undef __FUNCT__
500e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType"
501e27a552bSJed Brown /*@C
502*61692a83SJed Brown   TSRosWSetType - Set the type of Rosenbrock-W scheme
503e27a552bSJed Brown 
504e27a552bSJed Brown   Logically collective
505e27a552bSJed Brown 
506e27a552bSJed Brown   Input Parameter:
507e27a552bSJed Brown +  ts - timestepping context
508*61692a83SJed Brown -  rostype - type of Rosenbrock-W scheme
509e27a552bSJed Brown 
510e27a552bSJed Brown   Level: intermediate
511e27a552bSJed Brown 
512e27a552bSJed Brown .seealso: TSRosWGetType()
513e27a552bSJed Brown @*/
514*61692a83SJed Brown PetscErrorCode TSRosWSetType(TS ts,const TSRosWType rostype)
515e27a552bSJed Brown {
516e27a552bSJed Brown   PetscErrorCode ierr;
517e27a552bSJed Brown 
518e27a552bSJed Brown   PetscFunctionBegin;
519e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
520*61692a83SJed Brown   ierr = PetscTryMethod(ts,"TSRosWSetType_C",(TS,const TSRosWType),(ts,rostype));CHKERRQ(ierr);
521e27a552bSJed Brown   PetscFunctionReturn(0);
522e27a552bSJed Brown }
523e27a552bSJed Brown 
524e27a552bSJed Brown #undef __FUNCT__
525e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType"
526e27a552bSJed Brown /*@C
527*61692a83SJed Brown   TSRosWGetType - Get the type of Rosenbrock-W scheme
528e27a552bSJed Brown 
529e27a552bSJed Brown   Logically collective
530e27a552bSJed Brown 
531e27a552bSJed Brown   Input Parameter:
532e27a552bSJed Brown .  ts - timestepping context
533e27a552bSJed Brown 
534e27a552bSJed Brown   Output Parameter:
535*61692a83SJed Brown .  rostype - type of Rosenbrock-W scheme
536e27a552bSJed Brown 
537e27a552bSJed Brown   Level: intermediate
538e27a552bSJed Brown 
539e27a552bSJed Brown .seealso: TSRosWGetType()
540e27a552bSJed Brown @*/
541*61692a83SJed Brown PetscErrorCode TSRosWGetType(TS ts,const TSRosWType *rostype)
542e27a552bSJed Brown {
543e27a552bSJed Brown   PetscErrorCode ierr;
544e27a552bSJed Brown 
545e27a552bSJed Brown   PetscFunctionBegin;
546e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
547*61692a83SJed Brown   ierr = PetscUseMethod(ts,"TSRosWGetType_C",(TS,const TSRosWType*),(ts,rostype));CHKERRQ(ierr);
548e27a552bSJed Brown   PetscFunctionReturn(0);
549e27a552bSJed Brown }
550e27a552bSJed Brown 
551e27a552bSJed Brown #undef __FUNCT__
552*61692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian"
553e27a552bSJed Brown /*@C
554*61692a83SJed Brown   TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step.
555e27a552bSJed Brown 
556e27a552bSJed Brown   Logically collective
557e27a552bSJed Brown 
558e27a552bSJed Brown   Input Parameter:
559e27a552bSJed Brown +  ts - timestepping context
560*61692a83SJed Brown -  flg - PETSC_TRUE to recompute the Jacobian at each stage
561e27a552bSJed Brown 
562e27a552bSJed Brown   Level: intermediate
563e27a552bSJed Brown 
564e27a552bSJed Brown .seealso: TSRosWGetType()
565e27a552bSJed Brown @*/
566*61692a83SJed Brown PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg)
567e27a552bSJed Brown {
568e27a552bSJed Brown   PetscErrorCode ierr;
569e27a552bSJed Brown 
570e27a552bSJed Brown   PetscFunctionBegin;
571e27a552bSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
572*61692a83SJed Brown   ierr = PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
573e27a552bSJed Brown   PetscFunctionReturn(0);
574e27a552bSJed Brown }
575e27a552bSJed Brown 
576e27a552bSJed Brown EXTERN_C_BEGIN
577e27a552bSJed Brown #undef __FUNCT__
578e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType_RosW"
579*61692a83SJed Brown PetscErrorCode  TSRosWGetType_RosW(TS ts,const TSRosWType *rostype)
580e27a552bSJed Brown {
581*61692a83SJed Brown   TS_RosW        *ros = (TS_RosW*)ts->data;
582e27a552bSJed Brown   PetscErrorCode ierr;
583e27a552bSJed Brown 
584e27a552bSJed Brown   PetscFunctionBegin;
585*61692a83SJed Brown   if (!ros->tableau) {ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);}
586*61692a83SJed Brown   *rostype = ros->tableau->name;
587e27a552bSJed Brown   PetscFunctionReturn(0);
588e27a552bSJed Brown }
589e27a552bSJed Brown #undef __FUNCT__
590e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType_RosW"
591*61692a83SJed Brown PetscErrorCode  TSRosWSetType_RosW(TS ts,const TSRosWType rostype)
592e27a552bSJed Brown {
593*61692a83SJed Brown   TS_RosW         *ros = (TS_RosW*)ts->data;
594e27a552bSJed Brown   PetscErrorCode  ierr;
595e27a552bSJed Brown   PetscBool       match;
596*61692a83SJed Brown   RosWTableauLink link;
597e27a552bSJed Brown 
598e27a552bSJed Brown   PetscFunctionBegin;
599*61692a83SJed Brown   if (ros->tableau) {
600*61692a83SJed Brown     ierr = PetscStrcmp(ros->tableau->name,rostype,&match);CHKERRQ(ierr);
601e27a552bSJed Brown     if (match) PetscFunctionReturn(0);
602e27a552bSJed Brown   }
603*61692a83SJed Brown   for (link = RosWTableauList; link; link=link->next) {
604*61692a83SJed Brown     ierr = PetscStrcmp(link->tab.name,rostype,&match);CHKERRQ(ierr);
605e27a552bSJed Brown     if (match) {
606e27a552bSJed Brown       ierr = TSReset_RosW(ts);CHKERRQ(ierr);
607*61692a83SJed Brown       ros->tableau = &link->tab;
608e27a552bSJed Brown       PetscFunctionReturn(0);
609e27a552bSJed Brown     }
610e27a552bSJed Brown   }
611*61692a83SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype);
612e27a552bSJed Brown   PetscFunctionReturn(0);
613e27a552bSJed Brown }
614*61692a83SJed Brown 
615e27a552bSJed Brown #undef __FUNCT__
616*61692a83SJed Brown #define __FUNCT__ "TSRosWSetRecomputeJacobian_RosW"
617*61692a83SJed Brown PetscErrorCode  TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg)
618e27a552bSJed Brown {
619*61692a83SJed Brown   TS_RosW *ros = (TS_RosW*)ts->data;
620e27a552bSJed Brown 
621e27a552bSJed Brown   PetscFunctionBegin;
622*61692a83SJed Brown   ros->recompute_jacobian = flg;
623e27a552bSJed Brown   PetscFunctionReturn(0);
624e27a552bSJed Brown }
625e27a552bSJed Brown EXTERN_C_END
626e27a552bSJed Brown 
627e27a552bSJed Brown /* ------------------------------------------------------------ */
628e27a552bSJed Brown /*MC
629e27a552bSJed Brown       TSRosW - ODE solver using Rosenbrock-W schemes
630e27a552bSJed Brown 
631e27a552bSJed Brown   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
632e27a552bSJed Brown   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
633e27a552bSJed Brown   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
634e27a552bSJed Brown 
635e27a552bSJed Brown   Notes:
636*61692a83SJed Brown   This method currently only works with autonomous ODE and DAE.
637*61692a83SJed Brown 
638*61692a83SJed Brown   Developer notes:
639*61692a83SJed Brown   Rosenbrock-W methods are typically specified for autonomous ODE
640*61692a83SJed Brown 
641*61692a83SJed Brown $  xdot = f(x)
642*61692a83SJed Brown 
643*61692a83SJed Brown   by the stage equations
644*61692a83SJed Brown 
645*61692a83SJed Brown $  k_i = h f(x_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j
646*61692a83SJed Brown 
647*61692a83SJed Brown   and step completion formula
648*61692a83SJed Brown 
649*61692a83SJed Brown $  x_1 = x_0 + sum_j b_j k_j
650*61692a83SJed Brown 
651*61692a83SJed Brown   with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(x)
652*61692a83SJed Brown   and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner,
653*61692a83SJed Brown   we define new variables for the stage equations
654*61692a83SJed Brown 
655*61692a83SJed Brown $  y_i = gamma_ij k_j
656*61692a83SJed Brown 
657*61692a83SJed Brown   The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define
658*61692a83SJed Brown 
659*61692a83SJed Brown $  A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-i}
660*61692a83SJed Brown 
661*61692a83SJed Brown   to rewrite the method as
662*61692a83SJed Brown 
663*61692a83SJed 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
664*61692a83SJed Brown $  x_1 = x_0 + sum_j bt_j y_j
665*61692a83SJed Brown 
666*61692a83SJed Brown    where we have introduced the mass matrix M. Continue by defining
667*61692a83SJed Brown 
668*61692a83SJed Brown $  ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j
669*61692a83SJed Brown 
670*61692a83SJed Brown    or, more compactly in tensor notation
671*61692a83SJed Brown 
672*61692a83SJed Brown $  Ydot = 1/h (Gamma^{-1} \otimes I) Y .
673*61692a83SJed Brown 
674*61692a83SJed Brown    Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current
675*61692a83SJed Brown    stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the
676*61692a83SJed Brown    equation
677*61692a83SJed Brown 
678*61692a83SJed Brown $  g(x_0 + sum_j a_ij y_j + y_i, ydot_i) = 0
679*61692a83SJed Brown 
680*61692a83SJed Brown    with initial guess y_i = 0.
681e27a552bSJed Brown 
682e27a552bSJed Brown   Level: beginner
683e27a552bSJed Brown 
684e27a552bSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSRosWRegister()
685e27a552bSJed Brown 
686e27a552bSJed Brown M*/
687e27a552bSJed Brown EXTERN_C_BEGIN
688e27a552bSJed Brown #undef __FUNCT__
689e27a552bSJed Brown #define __FUNCT__ "TSCreate_RosW"
690e27a552bSJed Brown PetscErrorCode  TSCreate_RosW(TS ts)
691e27a552bSJed Brown {
692*61692a83SJed Brown   TS_RosW        *ros;
693e27a552bSJed Brown   PetscErrorCode ierr;
694e27a552bSJed Brown 
695e27a552bSJed Brown   PetscFunctionBegin;
696e27a552bSJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
697e27a552bSJed Brown   ierr = TSRosWInitializePackage(PETSC_NULL);CHKERRQ(ierr);
698e27a552bSJed Brown #endif
699e27a552bSJed Brown 
700e27a552bSJed Brown   ts->ops->reset          = TSReset_RosW;
701e27a552bSJed Brown   ts->ops->destroy        = TSDestroy_RosW;
702e27a552bSJed Brown   ts->ops->view           = TSView_RosW;
703e27a552bSJed Brown   ts->ops->setup          = TSSetUp_RosW;
704e27a552bSJed Brown   ts->ops->step           = TSStep_RosW;
705e27a552bSJed Brown   ts->ops->interpolate    = TSInterpolate_RosW;
706e27a552bSJed Brown   ts->ops->setfromoptions = TSSetFromOptions_RosW;
707e27a552bSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_RosW;
708e27a552bSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_RosW;
709e27a552bSJed Brown 
710*61692a83SJed Brown   ierr = PetscNewLog(ts,TS_RosW,&ros);CHKERRQ(ierr);
711*61692a83SJed Brown   ts->data = (void*)ros;
712e27a552bSJed Brown 
713e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","TSRosWGetType_RosW",TSRosWGetType_RosW);CHKERRQ(ierr);
714e27a552bSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","TSRosWSetType_RosW",TSRosWSetType_RosW);CHKERRQ(ierr);
715*61692a83SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","TSRosWSetRecomputeJacobian_RosW",TSRosWSetRecomputeJacobian_RosW);CHKERRQ(ierr);
716e27a552bSJed Brown   PetscFunctionReturn(0);
717e27a552bSJed Brown }
718e27a552bSJed Brown EXTERN_C_END
719