xref: /petsc/src/ts/impls/implicit/irk/irk.c (revision d2567f346889bc69f462cb815202b652ec2f297c)
1*d2567f34SHong Zhang /*
2*d2567f34SHong Zhang   Code for timestepping with implicit Runge-Kutta method
3*d2567f34SHong Zhang 
4*d2567f34SHong Zhang   Notes:
5*d2567f34SHong Zhang   The general system is written as
6*d2567f34SHong Zhang 
7*d2567f34SHong Zhang   F(t,U,Udot) = 0
8*d2567f34SHong Zhang 
9*d2567f34SHong Zhang */
10*d2567f34SHong Zhang #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
11*d2567f34SHong Zhang #include <petscdm.h>
12*d2567f34SHong Zhang #include <petscdt.h>
13*d2567f34SHong Zhang 
14*d2567f34SHong Zhang static TSIRKType         TSIRKDefault = TSIRKGAUSS;
15*d2567f34SHong Zhang static PetscBool         TSIRKRegisterAllCalled;
16*d2567f34SHong Zhang static PetscBool         TSIRKPackageInitialized;
17*d2567f34SHong Zhang static PetscFunctionList TSIRKList;
18*d2567f34SHong Zhang 
19*d2567f34SHong Zhang struct _IRKTableau{
20*d2567f34SHong Zhang   PetscReal   *A,*b,*c;
21*d2567f34SHong Zhang   PetscScalar *A_inv,*A_inv_rowsum,*I_s;
22*d2567f34SHong Zhang   PetscReal   *binterp;   /* Dense output formula */
23*d2567f34SHong Zhang };
24*d2567f34SHong Zhang 
25*d2567f34SHong Zhang typedef struct _IRKTableau *IRKTableau;
26*d2567f34SHong Zhang 
27*d2567f34SHong Zhang typedef struct {
28*d2567f34SHong Zhang   char         *method_name;
29*d2567f34SHong Zhang   PetscInt     order;            /* Classical approximation order of the method */
30*d2567f34SHong Zhang   PetscInt     nstages;          /* Number of stages */
31*d2567f34SHong Zhang   PetscBool    stiffly_accurate;
32*d2567f34SHong Zhang   PetscInt     pinterp;          /* Interpolation order */
33*d2567f34SHong Zhang   IRKTableau   tableau;
34*d2567f34SHong Zhang   Vec          U0;               /* Backup vector */
35*d2567f34SHong Zhang   Vec          Z;                /* Combined stage vector */
36*d2567f34SHong Zhang   Vec          *Y;               /* States computed during the step */
37*d2567f34SHong Zhang   Vec          Ydot;             /* Work vector holding time derivatives during residual evaluation */
38*d2567f34SHong Zhang   Vec          U;                /* U is used to compute Ydot = shift(Y-U) */
39*d2567f34SHong Zhang   Vec          *YdotI;           /* Work vectors to hold the residual evaluation */
40*d2567f34SHong Zhang   Mat          TJ;               /* KAIJ matrix for the Jacobian of the combined system */
41*d2567f34SHong Zhang   PetscScalar  *work;            /* Scalar work */
42*d2567f34SHong Zhang   TSStepStatus status;
43*d2567f34SHong Zhang   PetscBool    rebuild_completion;
44*d2567f34SHong Zhang   PetscReal    ccfl;
45*d2567f34SHong Zhang } TS_IRK;
46*d2567f34SHong Zhang 
47*d2567f34SHong Zhang /*@C
48*d2567f34SHong Zhang    TSIRKTableauCreate - create the tableau for TSIRK and provide the entries
49*d2567f34SHong Zhang 
50*d2567f34SHong Zhang    Not Collective
51*d2567f34SHong Zhang 
52*d2567f34SHong Zhang    Input Parameters:
53*d2567f34SHong Zhang +  ts - timestepping context
54*d2567f34SHong Zhang .  nstages - number of stages, this is the dimension of the matrices below
55*d2567f34SHong Zhang .  A - stage coefficients (dimension nstages*nstages, row-major)
56*d2567f34SHong Zhang .  b - step completion table (dimension nstages)
57*d2567f34SHong Zhang .  c - abscissa (dimension nstages)
58*d2567f34SHong Zhang -  binterp - coefficients of the interpolation formula (dimension nstages)
59*d2567f34SHong Zhang .  A_inv - inverse of A (dimension nstages*nstages, row-major)
60*d2567f34SHong Zhang .  A_inv_rowsum - row sum of the inverse of A (dimension nstages)
61*d2567f34SHong Zhang .  I_s - identity matrix (dimension nstages*nstages)
62*d2567f34SHong Zhang 
63*d2567f34SHong Zhang    Level: advanced
64*d2567f34SHong Zhang 
65*d2567f34SHong Zhang .seealso: TSIRK, TSIRKRegister()
66*d2567f34SHong Zhang @*/
67*d2567f34SHong Zhang PetscErrorCode TSIRKTableauCreate(TS ts,PetscInt nstages,const PetscReal *A,const PetscReal *b,const PetscReal *c,const PetscReal *binterp,const PetscScalar *A_inv,const PetscScalar *A_inv_rowsum,const PetscScalar *I_s)
68*d2567f34SHong Zhang {
69*d2567f34SHong Zhang   TS_IRK         *irk = (TS_IRK*)ts->data;
70*d2567f34SHong Zhang   IRKTableau     tab = irk->tableau;
71*d2567f34SHong Zhang   PetscErrorCode ierr;
72*d2567f34SHong Zhang 
73*d2567f34SHong Zhang   PetscFunctionBegin;
74*d2567f34SHong Zhang   irk->order = nstages;
75*d2567f34SHong Zhang   ierr = PetscMalloc3(PetscSqr(nstages),&tab->A,PetscSqr(nstages),&tab->A_inv,PetscSqr(nstages),&tab->I_s);CHKERRQ(ierr);
76*d2567f34SHong Zhang   ierr = PetscMalloc4(nstages,&tab->b,nstages,&tab->c,nstages,&tab->binterp,nstages,&tab->A_inv_rowsum);CHKERRQ(ierr);
77*d2567f34SHong Zhang   ierr = PetscArraycpy(tab->A,A,PetscSqr(nstages));CHKERRQ(ierr);
78*d2567f34SHong Zhang   ierr = PetscArraycpy(tab->b,b,nstages);CHKERRQ(ierr);
79*d2567f34SHong Zhang   ierr = PetscArraycpy(tab->c,c,nstages);CHKERRQ(ierr);
80*d2567f34SHong Zhang   /* optional coefficient arrays */
81*d2567f34SHong Zhang   if (binterp) {
82*d2567f34SHong Zhang     ierr = PetscArraycpy(tab->binterp,binterp,nstages);CHKERRQ(ierr);
83*d2567f34SHong Zhang   }
84*d2567f34SHong Zhang   if (A_inv) {
85*d2567f34SHong Zhang     ierr = PetscArraycpy(tab->A_inv,A_inv,PetscSqr(nstages));CHKERRQ(ierr);
86*d2567f34SHong Zhang   }
87*d2567f34SHong Zhang   if (A_inv_rowsum) {
88*d2567f34SHong Zhang     ierr = PetscArraycpy(tab->A_inv_rowsum,A_inv_rowsum,nstages);CHKERRQ(ierr);
89*d2567f34SHong Zhang   }
90*d2567f34SHong Zhang   if (I_s) {
91*d2567f34SHong Zhang     ierr = PetscArraycpy(tab->I_s,I_s,PetscSqr(nstages));CHKERRQ(ierr);
92*d2567f34SHong Zhang   }
93*d2567f34SHong Zhang   PetscFunctionReturn(0);
94*d2567f34SHong Zhang }
95*d2567f34SHong Zhang 
96*d2567f34SHong Zhang /* Arrays should be freed with PetscFree3(A,b,c) */
97*d2567f34SHong Zhang static PetscErrorCode TSIRKCreate_Gauss(TS ts)
98*d2567f34SHong Zhang {
99*d2567f34SHong Zhang   PetscInt       nstages;
100*d2567f34SHong Zhang   PetscReal      *gauss_A_real,*gauss_b,*b,*gauss_c;
101*d2567f34SHong Zhang   PetscScalar    *gauss_A,*gauss_A_inv,*gauss_A_inv_rowsum,*I_s;
102*d2567f34SHong Zhang   PetscScalar    *G0,*G1;
103*d2567f34SHong Zhang   PetscInt       i,j;
104*d2567f34SHong Zhang   Mat            G0mat,G1mat,Amat;
105*d2567f34SHong Zhang   PetscErrorCode ierr;
106*d2567f34SHong Zhang 
107*d2567f34SHong Zhang   PetscFunctionBegin;
108*d2567f34SHong Zhang   ierr = TSIRKGetNumStages(ts,&nstages);CHKERRQ(ierr);
109*d2567f34SHong Zhang   ierr = PetscMalloc3(PetscSqr(nstages),&gauss_A_real,nstages,&gauss_b,nstages,&gauss_c);CHKERRQ(ierr);
110*d2567f34SHong Zhang   ierr = PetscMalloc4(PetscSqr(nstages),&gauss_A,PetscSqr(nstages),&gauss_A_inv,nstages,&gauss_A_inv_rowsum,PetscSqr(nstages),&I_s);CHKERRQ(ierr);
111*d2567f34SHong Zhang   ierr = PetscMalloc3(nstages,&b,PetscSqr(nstages),&G0,PetscSqr(nstages),&G1);CHKERRQ(ierr);
112*d2567f34SHong Zhang   ierr = PetscDTGaussQuadrature(nstages,0.,1.,gauss_c,b);CHKERRQ(ierr);
113*d2567f34SHong Zhang   for (i=0; i<nstages; i++) gauss_b[i] = b[i]; /* copy to possibly-complex array */
114*d2567f34SHong Zhang 
115*d2567f34SHong Zhang   /* A^T = G0^{-1} G1 */
116*d2567f34SHong Zhang   for (i=0; i<nstages; i++) {
117*d2567f34SHong Zhang     for (j=0; j<nstages; j++) {
118*d2567f34SHong Zhang       G0[i*nstages+j] = PetscPowRealInt(gauss_c[i],j);
119*d2567f34SHong Zhang       G1[i*nstages+j] = PetscPowRealInt(gauss_c[i],j+1)/(j+1);
120*d2567f34SHong Zhang     }
121*d2567f34SHong Zhang   }
122*d2567f34SHong Zhang   /* The arrays above are row-aligned, but we create dense matrices as the transpose */
123*d2567f34SHong Zhang   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nstages,nstages,G0,&G0mat);CHKERRQ(ierr);
124*d2567f34SHong Zhang   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nstages,nstages,G1,&G1mat);CHKERRQ(ierr);
125*d2567f34SHong Zhang   ierr = MatCreateSeqDense(PETSC_COMM_SELF,nstages,nstages,gauss_A,&Amat);CHKERRQ(ierr);
126*d2567f34SHong Zhang   ierr = MatLUFactor(G0mat,NULL,NULL,NULL);CHKERRQ(ierr);
127*d2567f34SHong Zhang   ierr = MatMatSolve(G0mat,G1mat,Amat);CHKERRQ(ierr);
128*d2567f34SHong Zhang   ierr = MatTranspose(Amat,MAT_INPLACE_MATRIX,&Amat);CHKERRQ(ierr);
129*d2567f34SHong Zhang   for (i=0; i<nstages; i++)
130*d2567f34SHong Zhang     for (j=0; j<nstages; j++)
131*d2567f34SHong Zhang       gauss_A_real[i*nstages+j] = PetscRealPart(gauss_A[i*nstages+j]);
132*d2567f34SHong Zhang 
133*d2567f34SHong Zhang   ierr = MatDestroy(&G0mat);CHKERRQ(ierr);
134*d2567f34SHong Zhang   ierr = MatDestroy(&G1mat);CHKERRQ(ierr);
135*d2567f34SHong Zhang   ierr = MatDestroy(&Amat);CHKERRQ(ierr);
136*d2567f34SHong Zhang   ierr = PetscFree3(b,G0,G1);CHKERRQ(ierr);
137*d2567f34SHong Zhang 
138*d2567f34SHong Zhang   {/* Invert A */
139*d2567f34SHong Zhang     /* PETSc does not provide a routine to calculate the inverse of a general matrix.
140*d2567f34SHong Zhang      * To get the inverse of A, we form a sequential BAIJ matrix from it, consisting of a single block with block size
141*d2567f34SHong Zhang      * equal to the dimension of A, and then use MatInvertBlockDiagonal(). */
142*d2567f34SHong Zhang     Mat               A_baij;
143*d2567f34SHong Zhang     PetscInt          idxm[1]={0},idxn[1]={0};
144*d2567f34SHong Zhang     const PetscScalar *A_inv;
145*d2567f34SHong Zhang 
146*d2567f34SHong Zhang     ierr = MatCreateSeqBAIJ(PETSC_COMM_SELF,nstages,nstages,nstages,1,NULL,&A_baij);CHKERRQ(ierr);
147*d2567f34SHong Zhang     ierr = MatSetOption(A_baij,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
148*d2567f34SHong Zhang     ierr = MatSetValuesBlocked(A_baij,1,idxm,1,idxn,gauss_A,INSERT_VALUES);CHKERRQ(ierr);
149*d2567f34SHong Zhang     ierr = MatAssemblyBegin(A_baij,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
150*d2567f34SHong Zhang     ierr = MatAssemblyEnd(A_baij,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
151*d2567f34SHong Zhang     ierr = MatInvertBlockDiagonal(A_baij,&A_inv);CHKERRQ(ierr);
152*d2567f34SHong Zhang     ierr = PetscMemcpy(gauss_A_inv,A_inv,nstages*nstages*sizeof(PetscScalar));CHKERRQ(ierr);
153*d2567f34SHong Zhang     ierr = MatDestroy(&A_baij);CHKERRQ(ierr);
154*d2567f34SHong Zhang   }
155*d2567f34SHong Zhang 
156*d2567f34SHong Zhang   /* Compute row sums A_inv_rowsum and identity I_s */
157*d2567f34SHong Zhang   for (i=0; i<nstages; i++) {
158*d2567f34SHong Zhang     gauss_A_inv_rowsum[i] = 0;
159*d2567f34SHong Zhang     for (j=0; j<nstages; j++) {
160*d2567f34SHong Zhang       gauss_A_inv_rowsum[i] += gauss_A_inv[i+nstages*j];
161*d2567f34SHong Zhang       I_s[i+nstages*j] = 1.*(i == j);
162*d2567f34SHong Zhang     }
163*d2567f34SHong Zhang   }
164*d2567f34SHong Zhang   ierr = TSIRKTableauCreate(ts,nstages,gauss_A_real,gauss_b,gauss_c,NULL,gauss_A_inv,gauss_A_inv_rowsum,I_s);CHKERRQ(ierr);
165*d2567f34SHong Zhang   ierr = PetscFree3(gauss_A_real,gauss_b,gauss_c);CHKERRQ(ierr);
166*d2567f34SHong Zhang   ierr = PetscFree4(gauss_A,gauss_A_inv,gauss_A_inv_rowsum,I_s);CHKERRQ(ierr);
167*d2567f34SHong Zhang   PetscFunctionReturn(0);
168*d2567f34SHong Zhang }
169*d2567f34SHong Zhang 
170*d2567f34SHong Zhang /*@C
171*d2567f34SHong Zhang    TSIRKRegister -  adds a TSIRK implementation
172*d2567f34SHong Zhang 
173*d2567f34SHong Zhang    Not Collective
174*d2567f34SHong Zhang 
175*d2567f34SHong Zhang    Input Parameters:
176*d2567f34SHong Zhang +  sname - name of user-defined IRK scheme
177*d2567f34SHong Zhang -  function - function to create method context
178*d2567f34SHong Zhang 
179*d2567f34SHong Zhang    Notes:
180*d2567f34SHong Zhang    TSIRKRegister() may be called multiple times to add several user-defined families.
181*d2567f34SHong Zhang 
182*d2567f34SHong Zhang    Sample usage:
183*d2567f34SHong Zhang .vb
184*d2567f34SHong Zhang    TSIRKRegister("my_scheme",MySchemeCreate);
185*d2567f34SHong Zhang .ve
186*d2567f34SHong Zhang 
187*d2567f34SHong Zhang    Then, your scheme can be chosen with the procedural interface via
188*d2567f34SHong Zhang $     TSIRKSetType(ts,"my_scheme")
189*d2567f34SHong Zhang    or at runtime via the option
190*d2567f34SHong Zhang $     -ts_irk_type my_scheme
191*d2567f34SHong Zhang 
192*d2567f34SHong Zhang    Level: advanced
193*d2567f34SHong Zhang 
194*d2567f34SHong Zhang .seealso: TSIRKRegisterAll()
195*d2567f34SHong Zhang @*/
196*d2567f34SHong Zhang PetscErrorCode TSIRKRegister(const char sname[],PetscErrorCode (*function)(TS))
197*d2567f34SHong Zhang {
198*d2567f34SHong Zhang   PetscErrorCode ierr;
199*d2567f34SHong Zhang 
200*d2567f34SHong Zhang   PetscFunctionBegin;
201*d2567f34SHong Zhang   ierr = TSIRKInitializePackage();CHKERRQ(ierr);
202*d2567f34SHong Zhang   ierr = PetscFunctionListAdd(&TSIRKList,sname,function);CHKERRQ(ierr);
203*d2567f34SHong Zhang   PetscFunctionReturn(0);
204*d2567f34SHong Zhang }
205*d2567f34SHong Zhang 
206*d2567f34SHong Zhang /*@C
207*d2567f34SHong Zhang   TSIRKRegisterAll - Registers all of the implicit Runge-Kutta methods in TSIRK
208*d2567f34SHong Zhang 
209*d2567f34SHong Zhang   Not Collective, but should be called by all processes which will need the schemes to be registered
210*d2567f34SHong Zhang 
211*d2567f34SHong Zhang   Level: advanced
212*d2567f34SHong Zhang 
213*d2567f34SHong Zhang .seealso:  TSIRKRegisterDestroy()
214*d2567f34SHong Zhang @*/
215*d2567f34SHong Zhang PetscErrorCode TSIRKRegisterAll(void)
216*d2567f34SHong Zhang {
217*d2567f34SHong Zhang   PetscErrorCode ierr;
218*d2567f34SHong Zhang 
219*d2567f34SHong Zhang   PetscFunctionBegin;
220*d2567f34SHong Zhang   if (TSIRKRegisterAllCalled) PetscFunctionReturn(0);
221*d2567f34SHong Zhang   TSIRKRegisterAllCalled = PETSC_TRUE;
222*d2567f34SHong Zhang 
223*d2567f34SHong Zhang   ierr = TSIRKRegister(TSIRKGAUSS,TSIRKCreate_Gauss);CHKERRQ(ierr);
224*d2567f34SHong Zhang   PetscFunctionReturn(0);
225*d2567f34SHong Zhang }
226*d2567f34SHong Zhang 
227*d2567f34SHong Zhang /*@C
228*d2567f34SHong Zhang    TSIRKRegisterDestroy - Frees the list of schemes that were registered by TSIRKRegister().
229*d2567f34SHong Zhang 
230*d2567f34SHong Zhang    Not Collective
231*d2567f34SHong Zhang 
232*d2567f34SHong Zhang    Level: advanced
233*d2567f34SHong Zhang 
234*d2567f34SHong Zhang .seealso: TSIRKRegister(), TSIRKRegisterAll()
235*d2567f34SHong Zhang @*/
236*d2567f34SHong Zhang PetscErrorCode TSIRKRegisterDestroy(void)
237*d2567f34SHong Zhang {
238*d2567f34SHong Zhang   PetscFunctionBegin;
239*d2567f34SHong Zhang   TSIRKRegisterAllCalled = PETSC_FALSE;
240*d2567f34SHong Zhang   PetscFunctionReturn(0);
241*d2567f34SHong Zhang }
242*d2567f34SHong Zhang 
243*d2567f34SHong Zhang /*@C
244*d2567f34SHong Zhang   TSIRKInitializePackage - This function initializes everything in the TSIRK package. It is called
245*d2567f34SHong Zhang   from TSInitializePackage().
246*d2567f34SHong Zhang 
247*d2567f34SHong Zhang   Level: developer
248*d2567f34SHong Zhang 
249*d2567f34SHong Zhang .seealso: PetscInitialize()
250*d2567f34SHong Zhang @*/
251*d2567f34SHong Zhang PetscErrorCode TSIRKInitializePackage(void)
252*d2567f34SHong Zhang {
253*d2567f34SHong Zhang   PetscErrorCode ierr;
254*d2567f34SHong Zhang 
255*d2567f34SHong Zhang   PetscFunctionBegin;
256*d2567f34SHong Zhang   if (TSIRKPackageInitialized) PetscFunctionReturn(0);
257*d2567f34SHong Zhang   TSIRKPackageInitialized = PETSC_TRUE;
258*d2567f34SHong Zhang   ierr = TSIRKRegisterAll();CHKERRQ(ierr);
259*d2567f34SHong Zhang   ierr = PetscRegisterFinalize(TSIRKFinalizePackage);CHKERRQ(ierr);
260*d2567f34SHong Zhang   PetscFunctionReturn(0);
261*d2567f34SHong Zhang }
262*d2567f34SHong Zhang 
263*d2567f34SHong Zhang /*@C
264*d2567f34SHong Zhang   TSIRKFinalizePackage - This function destroys everything in the TSIRK package. It is
265*d2567f34SHong Zhang   called from PetscFinalize().
266*d2567f34SHong Zhang 
267*d2567f34SHong Zhang   Level: developer
268*d2567f34SHong Zhang 
269*d2567f34SHong Zhang .seealso: PetscFinalize()
270*d2567f34SHong Zhang @*/
271*d2567f34SHong Zhang PetscErrorCode TSIRKFinalizePackage(void)
272*d2567f34SHong Zhang {
273*d2567f34SHong Zhang   PetscErrorCode ierr;
274*d2567f34SHong Zhang 
275*d2567f34SHong Zhang   PetscFunctionBegin;
276*d2567f34SHong Zhang   ierr = PetscFunctionListDestroy(&TSIRKList);CHKERRQ(ierr);
277*d2567f34SHong Zhang   TSIRKPackageInitialized = PETSC_FALSE;
278*d2567f34SHong Zhang   PetscFunctionReturn(0);
279*d2567f34SHong Zhang }
280*d2567f34SHong Zhang 
281*d2567f34SHong Zhang /*
282*d2567f34SHong Zhang  This function can be called before or after ts->vec_sol has been updated.
283*d2567f34SHong Zhang */
284*d2567f34SHong Zhang static PetscErrorCode TSEvaluateStep_IRK(TS ts,PetscInt order,Vec U,PetscBool *done)
285*d2567f34SHong Zhang {
286*d2567f34SHong Zhang   TS_IRK         *irk = (TS_IRK*)ts->data;
287*d2567f34SHong Zhang   IRKTableau     tab = irk->tableau;
288*d2567f34SHong Zhang   Vec            *YdotI = irk->YdotI;
289*d2567f34SHong Zhang   PetscScalar    *w = irk->work;
290*d2567f34SHong Zhang   PetscReal      h;
291*d2567f34SHong Zhang   PetscInt       j;
292*d2567f34SHong Zhang   PetscErrorCode ierr;
293*d2567f34SHong Zhang 
294*d2567f34SHong Zhang   PetscFunctionBegin;
295*d2567f34SHong Zhang   switch (irk->status) {
296*d2567f34SHong Zhang   case TS_STEP_INCOMPLETE:
297*d2567f34SHong Zhang   case TS_STEP_PENDING:
298*d2567f34SHong Zhang     h = ts->time_step; break;
299*d2567f34SHong Zhang   case TS_STEP_COMPLETE:
300*d2567f34SHong Zhang     h = ts->ptime - ts->ptime_prev; break;
301*d2567f34SHong Zhang   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
302*d2567f34SHong Zhang   }
303*d2567f34SHong Zhang 
304*d2567f34SHong Zhang   ierr = VecCopy(ts->vec_sol,U);CHKERRQ(ierr);
305*d2567f34SHong Zhang   for (j=0; j<irk->nstages; j++) w[j] = h*tab->b[j];
306*d2567f34SHong Zhang   ierr = VecMAXPY(U,irk->nstages,w,YdotI);CHKERRQ(ierr);
307*d2567f34SHong Zhang   PetscFunctionReturn(0);
308*d2567f34SHong Zhang }
309*d2567f34SHong Zhang 
310*d2567f34SHong Zhang static PetscErrorCode TSRollBack_IRK(TS ts)
311*d2567f34SHong Zhang {
312*d2567f34SHong Zhang   TS_IRK         *irk = (TS_IRK*)ts->data;
313*d2567f34SHong Zhang   PetscErrorCode ierr;
314*d2567f34SHong Zhang 
315*d2567f34SHong Zhang   PetscFunctionBegin;
316*d2567f34SHong Zhang   ierr = VecCopy(irk->U0,ts->vec_sol);CHKERRQ(ierr);
317*d2567f34SHong Zhang   PetscFunctionReturn(0);
318*d2567f34SHong Zhang }
319*d2567f34SHong Zhang 
320*d2567f34SHong Zhang static PetscErrorCode TSStep_IRK(TS ts)
321*d2567f34SHong Zhang {
322*d2567f34SHong Zhang   TS_IRK          *irk = (TS_IRK*)ts->data;
323*d2567f34SHong Zhang   IRKTableau      tab = irk->tableau;
324*d2567f34SHong Zhang   PetscScalar     *A_inv = tab->A_inv,*A_inv_rowsum = tab->A_inv_rowsum;
325*d2567f34SHong Zhang   const PetscInt  nstages = irk->nstages;
326*d2567f34SHong Zhang   SNES            snes;
327*d2567f34SHong Zhang   PetscInt        i,j,its,lits,bs;
328*d2567f34SHong Zhang   TSAdapt         adapt;
329*d2567f34SHong Zhang   PetscInt        rejections = 0;
330*d2567f34SHong Zhang   PetscBool       accept = PETSC_TRUE;
331*d2567f34SHong Zhang   PetscReal       next_time_step = ts->time_step;
332*d2567f34SHong Zhang   PetscErrorCode  ierr;
333*d2567f34SHong Zhang 
334*d2567f34SHong Zhang   PetscFunctionBegin;
335*d2567f34SHong Zhang   if (!ts->steprollback) {
336*d2567f34SHong Zhang     ierr = VecCopy(ts->vec_sol,irk->U0);CHKERRQ(ierr);
337*d2567f34SHong Zhang   }
338*d2567f34SHong Zhang   ierr = VecGetBlockSize(ts->vec_sol,&bs);CHKERRQ(ierr);
339*d2567f34SHong Zhang   for (i=0; i<nstages; i++) {
340*d2567f34SHong Zhang     ierr = VecStrideScatter(ts->vec_sol,i*bs,irk->Z,INSERT_VALUES);CHKERRQ(ierr);
341*d2567f34SHong Zhang   }
342*d2567f34SHong Zhang 
343*d2567f34SHong Zhang   irk->status = TS_STEP_INCOMPLETE;
344*d2567f34SHong Zhang   while (!ts->reason && irk->status != TS_STEP_COMPLETE) {
345*d2567f34SHong Zhang     ierr = VecCopy(ts->vec_sol,irk->U);CHKERRQ(ierr);
346*d2567f34SHong Zhang     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
347*d2567f34SHong Zhang     ierr = SNESSolve(snes,NULL,irk->Z);CHKERRQ(ierr);
348*d2567f34SHong Zhang     ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
349*d2567f34SHong Zhang     ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
350*d2567f34SHong Zhang     ts->snes_its += its; ts->ksp_its += lits;
351*d2567f34SHong Zhang     ierr = VecStrideGatherAll(irk->Z,irk->Y,INSERT_VALUES);CHKERRQ(ierr);
352*d2567f34SHong Zhang     for (i=0; i<nstages; i++) {
353*d2567f34SHong Zhang       ierr = VecZeroEntries(irk->YdotI[i]);CHKERRQ(ierr);
354*d2567f34SHong Zhang       for (j=0; j<nstages; j++) {
355*d2567f34SHong Zhang         ierr = VecAXPY(irk->YdotI[i],A_inv[i+j*nstages]/ts->time_step,irk->Y[j]);CHKERRQ(ierr);
356*d2567f34SHong Zhang       }
357*d2567f34SHong Zhang       ierr = VecAXPY(irk->YdotI[i],-A_inv_rowsum[i]/ts->time_step,irk->U);CHKERRQ(ierr);
358*d2567f34SHong Zhang     }
359*d2567f34SHong Zhang     irk->status = TS_STEP_INCOMPLETE;
360*d2567f34SHong Zhang     ierr = TSEvaluateStep_IRK(ts,irk->order,ts->vec_sol,NULL);CHKERRQ(ierr);
361*d2567f34SHong Zhang     irk->status = TS_STEP_PENDING;
362*d2567f34SHong Zhang     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
363*d2567f34SHong Zhang     ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
364*d2567f34SHong Zhang     irk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
365*d2567f34SHong Zhang     if (!accept) {
366*d2567f34SHong Zhang       ierr = TSRollBack_IRK(ts);CHKERRQ(ierr);
367*d2567f34SHong Zhang       ts->time_step = next_time_step;
368*d2567f34SHong Zhang       goto reject_step;
369*d2567f34SHong Zhang     }
370*d2567f34SHong Zhang 
371*d2567f34SHong Zhang     ts->ptime += ts->time_step;
372*d2567f34SHong Zhang     ts->time_step = next_time_step;
373*d2567f34SHong Zhang     break;
374*d2567f34SHong Zhang   reject_step:
375*d2567f34SHong Zhang     ts->reject++; accept = PETSC_FALSE;
376*d2567f34SHong Zhang     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
377*d2567f34SHong Zhang       ts->reason = TS_DIVERGED_STEP_REJECTED;
378*d2567f34SHong Zhang       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
379*d2567f34SHong Zhang     }
380*d2567f34SHong Zhang   }
381*d2567f34SHong Zhang   PetscFunctionReturn(0);
382*d2567f34SHong Zhang }
383*d2567f34SHong Zhang 
384*d2567f34SHong Zhang static PetscErrorCode TSInterpolate_IRK(TS ts,PetscReal itime,Vec U)
385*d2567f34SHong Zhang {
386*d2567f34SHong Zhang   TS_IRK          *irk = (TS_IRK*)ts->data;
387*d2567f34SHong Zhang   PetscInt        nstages = irk->nstages,pinterp = irk->pinterp,i,j;
388*d2567f34SHong Zhang   PetscReal       h;
389*d2567f34SHong Zhang   PetscReal       tt,t;
390*d2567f34SHong Zhang   PetscScalar     *bt;
391*d2567f34SHong Zhang   const PetscReal *B = irk->tableau->binterp;
392*d2567f34SHong Zhang   PetscErrorCode  ierr;
393*d2567f34SHong Zhang 
394*d2567f34SHong Zhang   PetscFunctionBegin;
395*d2567f34SHong Zhang   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSIRK %s does not have an interpolation formula",irk->method_name);
396*d2567f34SHong Zhang   switch (irk->status) {
397*d2567f34SHong Zhang   case TS_STEP_INCOMPLETE:
398*d2567f34SHong Zhang   case TS_STEP_PENDING:
399*d2567f34SHong Zhang     h = ts->time_step;
400*d2567f34SHong Zhang     t = (itime - ts->ptime)/h;
401*d2567f34SHong Zhang     break;
402*d2567f34SHong Zhang   case TS_STEP_COMPLETE:
403*d2567f34SHong Zhang     h = ts->ptime - ts->ptime_prev;
404*d2567f34SHong Zhang     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
405*d2567f34SHong Zhang     break;
406*d2567f34SHong Zhang   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
407*d2567f34SHong Zhang   }
408*d2567f34SHong Zhang   ierr = PetscMalloc1(nstages,&bt);CHKERRQ(ierr);
409*d2567f34SHong Zhang   for (i=0; i<nstages; i++) bt[i] = 0;
410*d2567f34SHong Zhang   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
411*d2567f34SHong Zhang     for (i=0; i<nstages; i++) {
412*d2567f34SHong Zhang       bt[i] += h * B[i*pinterp+j] * tt;
413*d2567f34SHong Zhang     }
414*d2567f34SHong Zhang   }
415*d2567f34SHong Zhang   ierr = VecMAXPY(U,nstages,bt,irk->YdotI);CHKERRQ(ierr);
416*d2567f34SHong Zhang   PetscFunctionReturn(0);
417*d2567f34SHong Zhang }
418*d2567f34SHong Zhang 
419*d2567f34SHong Zhang static PetscErrorCode TSIRKTableauReset(TS ts)
420*d2567f34SHong Zhang {
421*d2567f34SHong Zhang   TS_IRK         *irk = (TS_IRK*)ts->data;
422*d2567f34SHong Zhang   IRKTableau     tab = irk->tableau;
423*d2567f34SHong Zhang   PetscErrorCode ierr;
424*d2567f34SHong Zhang 
425*d2567f34SHong Zhang   PetscFunctionBegin;
426*d2567f34SHong Zhang   if (!tab) PetscFunctionReturn(0);
427*d2567f34SHong Zhang   ierr = PetscFree3(tab->A,tab->A_inv,tab->I_s);CHKERRQ(ierr);
428*d2567f34SHong Zhang   ierr = PetscFree4(tab->b,tab->c,tab->binterp,tab->A_inv_rowsum);CHKERRQ(ierr);
429*d2567f34SHong Zhang   PetscFunctionReturn(0);
430*d2567f34SHong Zhang }
431*d2567f34SHong Zhang 
432*d2567f34SHong Zhang static PetscErrorCode TSReset_IRK(TS ts)
433*d2567f34SHong Zhang {
434*d2567f34SHong Zhang   TS_IRK         *irk = (TS_IRK*)ts->data;
435*d2567f34SHong Zhang   PetscErrorCode ierr;
436*d2567f34SHong Zhang 
437*d2567f34SHong Zhang   PetscFunctionBegin;
438*d2567f34SHong Zhang   ierr = TSIRKTableauReset(ts);CHKERRQ(ierr);
439*d2567f34SHong Zhang   if (irk->tableau) {
440*d2567f34SHong Zhang     ierr = PetscFree(irk->tableau);CHKERRQ(ierr);
441*d2567f34SHong Zhang   }
442*d2567f34SHong Zhang   if (irk->method_name) {
443*d2567f34SHong Zhang     ierr = PetscFree(irk->method_name);CHKERRQ(ierr);
444*d2567f34SHong Zhang   }
445*d2567f34SHong Zhang   if (irk->work) {
446*d2567f34SHong Zhang     ierr = PetscFree(irk->work);CHKERRQ(ierr);
447*d2567f34SHong Zhang   }
448*d2567f34SHong Zhang   ierr = VecDestroyVecs(irk->nstages,&irk->Y);CHKERRQ(ierr);
449*d2567f34SHong Zhang   ierr = VecDestroyVecs(irk->nstages,&irk->YdotI);CHKERRQ(ierr);
450*d2567f34SHong Zhang   ierr = VecDestroy(&irk->Ydot);CHKERRQ(ierr);
451*d2567f34SHong Zhang   ierr = VecDestroy(&irk->Z);CHKERRQ(ierr);
452*d2567f34SHong Zhang   ierr = VecDestroy(&irk->U);CHKERRQ(ierr);
453*d2567f34SHong Zhang   ierr = VecDestroy(&irk->U0);CHKERRQ(ierr);
454*d2567f34SHong Zhang   ierr = MatDestroy(&irk->TJ);CHKERRQ(ierr);
455*d2567f34SHong Zhang   PetscFunctionReturn(0);
456*d2567f34SHong Zhang }
457*d2567f34SHong Zhang 
458*d2567f34SHong Zhang static PetscErrorCode TSIRKGetVecs(TS ts,DM dm,Vec *U)
459*d2567f34SHong Zhang {
460*d2567f34SHong Zhang   TS_IRK         *irk = (TS_IRK*)ts->data;
461*d2567f34SHong Zhang   PetscErrorCode ierr;
462*d2567f34SHong Zhang 
463*d2567f34SHong Zhang   PetscFunctionBegin;
464*d2567f34SHong Zhang   if (U) {
465*d2567f34SHong Zhang     if (dm && dm != ts->dm) {
466*d2567f34SHong Zhang       ierr = DMGetNamedGlobalVector(dm,"TSIRK_U",U);CHKERRQ(ierr);
467*d2567f34SHong Zhang     } else *U = irk->U;
468*d2567f34SHong Zhang   }
469*d2567f34SHong Zhang   PetscFunctionReturn(0);
470*d2567f34SHong Zhang }
471*d2567f34SHong Zhang 
472*d2567f34SHong Zhang static PetscErrorCode TSIRKRestoreVecs(TS ts,DM dm,Vec *U)
473*d2567f34SHong Zhang {
474*d2567f34SHong Zhang   PetscErrorCode ierr;
475*d2567f34SHong Zhang 
476*d2567f34SHong Zhang   PetscFunctionBegin;
477*d2567f34SHong Zhang   if (U) {
478*d2567f34SHong Zhang     if (dm && dm != ts->dm) {
479*d2567f34SHong Zhang       ierr = DMRestoreNamedGlobalVector(dm,"TSIRK_U",U);CHKERRQ(ierr);
480*d2567f34SHong Zhang     }
481*d2567f34SHong Zhang   }
482*d2567f34SHong Zhang   PetscFunctionReturn(0);
483*d2567f34SHong Zhang }
484*d2567f34SHong Zhang 
485*d2567f34SHong Zhang /*
486*d2567f34SHong Zhang   This defines the nonlinear equations that is to be solved with SNES
487*d2567f34SHong Zhang     G[e\otimes t + C*dt, Z, Zdot] = 0
488*d2567f34SHong Zhang     Zdot = (In \otimes S)*Z - (In \otimes Se) U
489*d2567f34SHong Zhang   where S = 1/(dt*A)
490*d2567f34SHong Zhang */
491*d2567f34SHong Zhang static PetscErrorCode SNESTSFormFunction_IRK(SNES snes,Vec ZC,Vec FC,TS ts)
492*d2567f34SHong Zhang {
493*d2567f34SHong Zhang   TS_IRK            *irk = (TS_IRK*)ts->data;
494*d2567f34SHong Zhang   IRKTableau        tab  = irk->tableau;
495*d2567f34SHong Zhang   const PetscInt    nstages = irk->nstages;
496*d2567f34SHong Zhang   const PetscReal   *c = tab->c;
497*d2567f34SHong Zhang   const PetscScalar *A_inv = tab->A_inv,*A_inv_rowsum = tab->A_inv_rowsum;
498*d2567f34SHong Zhang   DM                dm,dmsave;
499*d2567f34SHong Zhang   Vec               U,*YdotI = irk->YdotI,Ydot = irk->Ydot,*Y = irk->Y;
500*d2567f34SHong Zhang   PetscReal         h = ts->time_step;
501*d2567f34SHong Zhang   PetscInt          i,j;
502*d2567f34SHong Zhang   PetscErrorCode    ierr;
503*d2567f34SHong Zhang 
504*d2567f34SHong Zhang   PetscFunctionBegin;
505*d2567f34SHong Zhang   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
506*d2567f34SHong Zhang   ierr = TSIRKGetVecs(ts,dm,&U);CHKERRQ(ierr);
507*d2567f34SHong Zhang   ierr = VecStrideGatherAll(ZC,Y,INSERT_VALUES);CHKERRQ(ierr);
508*d2567f34SHong Zhang   dmsave = ts->dm;
509*d2567f34SHong Zhang   ts->dm = dm;
510*d2567f34SHong Zhang   for (i=0; i<nstages; i++) {
511*d2567f34SHong Zhang     ierr = VecZeroEntries(Ydot);CHKERRQ(ierr);
512*d2567f34SHong Zhang     for (j=0; j<nstages; j++) {
513*d2567f34SHong Zhang       ierr = VecAXPY(Ydot,A_inv[j*nstages+i]/h,Y[j]);CHKERRQ(ierr);
514*d2567f34SHong Zhang     }
515*d2567f34SHong Zhang     ierr = VecAXPY(Ydot,-A_inv_rowsum[i]/h,U);CHKERRQ(ierr); /* Ydot = (S \otimes In)*Z - (Se \otimes In) U */
516*d2567f34SHong Zhang     ierr = TSComputeIFunction(ts,ts->ptime+ts->time_step*c[i],Y[i],Ydot,YdotI[i],PETSC_FALSE);CHKERRQ(ierr);
517*d2567f34SHong Zhang   }
518*d2567f34SHong Zhang   ierr = VecStrideScatterAll(YdotI,FC,INSERT_VALUES);CHKERRQ(ierr);
519*d2567f34SHong Zhang   ts->dm = dmsave;
520*d2567f34SHong Zhang   ierr   = TSIRKRestoreVecs(ts,dm,&U);CHKERRQ(ierr);
521*d2567f34SHong Zhang   PetscFunctionReturn(0);
522*d2567f34SHong Zhang }
523*d2567f34SHong Zhang 
524*d2567f34SHong Zhang /*
525*d2567f34SHong Zhang    For explicit ODE, the Jacobian is
526*d2567f34SHong Zhang      JC = I_n \otimes S - J \otimes I_s
527*d2567f34SHong Zhang    For DAE, the Jacobian is
528*d2567f34SHong Zhang      JC = M_n \otimes S - J \otimes I_s
529*d2567f34SHong Zhang */
530*d2567f34SHong Zhang static PetscErrorCode SNESTSFormJacobian_IRK(SNES snes,Vec ZC,Mat JC,Mat JCpre,TS ts)
531*d2567f34SHong Zhang {
532*d2567f34SHong Zhang   TS_IRK          *irk = (TS_IRK*)ts->data;
533*d2567f34SHong Zhang   IRKTableau      tab  = irk->tableau;
534*d2567f34SHong Zhang   const PetscInt  nstages = irk->nstages;
535*d2567f34SHong Zhang   const PetscReal *c = tab->c;
536*d2567f34SHong Zhang   DM              dm,dmsave;
537*d2567f34SHong Zhang   Vec             *Y = irk->Y,Ydot = irk->Ydot;
538*d2567f34SHong Zhang   Mat             J;
539*d2567f34SHong Zhang   PetscScalar     *S;
540*d2567f34SHong Zhang   PetscInt        i,j,bs;
541*d2567f34SHong Zhang   PetscErrorCode  ierr;
542*d2567f34SHong Zhang 
543*d2567f34SHong Zhang   PetscFunctionBegin;
544*d2567f34SHong Zhang   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
545*d2567f34SHong Zhang   /* irk->Ydot has already been computed in SNESTSFormFunction_IRK (SNES guarantees this) */
546*d2567f34SHong Zhang   dmsave = ts->dm;
547*d2567f34SHong Zhang   ts->dm = dm;
548*d2567f34SHong Zhang   ierr = VecGetBlockSize(Y[nstages-1],&bs);CHKERRQ(ierr);
549*d2567f34SHong Zhang   if (ts->equation_type <= TS_EQ_ODE_EXPLICIT) { /* Support explicit formulas only */
550*d2567f34SHong Zhang     ierr = VecStrideGather(ZC,(nstages-1)*bs,Y[nstages-1],INSERT_VALUES);CHKERRQ(ierr);
551*d2567f34SHong Zhang     ierr = MatKAIJGetAIJ(JC,&J);CHKERRQ(ierr);
552*d2567f34SHong Zhang     ierr = TSComputeIJacobian(ts,ts->ptime+ts->time_step*c[nstages-1],Y[nstages-1],Ydot,0,J,J,PETSC_FALSE);CHKERRQ(ierr);
553*d2567f34SHong Zhang     ierr = MatKAIJGetS(JC,NULL,NULL,&S);CHKERRQ(ierr);
554*d2567f34SHong Zhang     for (i=0; i<nstages; i++)
555*d2567f34SHong Zhang       for (j=0; j<nstages; j++)
556*d2567f34SHong Zhang         S[i+nstages*j] = tab->A_inv[i+nstages*j]/ts->time_step;
557*d2567f34SHong Zhang     ierr = MatKAIJRestoreS(JC,&S);CHKERRQ(ierr);
558*d2567f34SHong Zhang   } else
559*d2567f34SHong Zhang     SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSIRK %s does not support implicit formula",irk->method_name); /* ToDo: need the mass matrix for DAE  */
560*d2567f34SHong Zhang   ts->dm = dmsave;
561*d2567f34SHong Zhang   PetscFunctionReturn(0);
562*d2567f34SHong Zhang }
563*d2567f34SHong Zhang 
564*d2567f34SHong Zhang static PetscErrorCode DMCoarsenHook_TSIRK(DM fine,DM coarse,void *ctx)
565*d2567f34SHong Zhang {
566*d2567f34SHong Zhang   PetscFunctionBegin;
567*d2567f34SHong Zhang   PetscFunctionReturn(0);
568*d2567f34SHong Zhang }
569*d2567f34SHong Zhang 
570*d2567f34SHong Zhang static PetscErrorCode DMRestrictHook_TSIRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
571*d2567f34SHong Zhang {
572*d2567f34SHong Zhang   TS             ts = (TS)ctx;
573*d2567f34SHong Zhang   Vec            U,U_c;
574*d2567f34SHong Zhang   PetscErrorCode ierr;
575*d2567f34SHong Zhang 
576*d2567f34SHong Zhang   PetscFunctionBegin;
577*d2567f34SHong Zhang   ierr = TSIRKGetVecs(ts,fine,&U);CHKERRQ(ierr);
578*d2567f34SHong Zhang   ierr = TSIRKGetVecs(ts,coarse,&U_c);CHKERRQ(ierr);
579*d2567f34SHong Zhang   ierr = MatRestrict(restrct,U,U_c);CHKERRQ(ierr);
580*d2567f34SHong Zhang   ierr = VecPointwiseMult(U_c,rscale,U_c);CHKERRQ(ierr);
581*d2567f34SHong Zhang   ierr = TSIRKRestoreVecs(ts,fine,&U);CHKERRQ(ierr);
582*d2567f34SHong Zhang   ierr = TSIRKRestoreVecs(ts,coarse,&U_c);CHKERRQ(ierr);
583*d2567f34SHong Zhang   PetscFunctionReturn(0);
584*d2567f34SHong Zhang }
585*d2567f34SHong Zhang 
586*d2567f34SHong Zhang static PetscErrorCode DMSubDomainHook_TSIRK(DM dm,DM subdm,void *ctx)
587*d2567f34SHong Zhang {
588*d2567f34SHong Zhang   PetscFunctionBegin;
589*d2567f34SHong Zhang   PetscFunctionReturn(0);
590*d2567f34SHong Zhang }
591*d2567f34SHong Zhang 
592*d2567f34SHong Zhang static PetscErrorCode DMSubDomainRestrictHook_TSIRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
593*d2567f34SHong Zhang {
594*d2567f34SHong Zhang   TS             ts = (TS)ctx;
595*d2567f34SHong Zhang   Vec            U,U_c;
596*d2567f34SHong Zhang   PetscErrorCode ierr;
597*d2567f34SHong Zhang 
598*d2567f34SHong Zhang   PetscFunctionBegin;
599*d2567f34SHong Zhang   ierr = TSIRKGetVecs(ts,dm,&U);CHKERRQ(ierr);
600*d2567f34SHong Zhang   ierr = TSIRKGetVecs(ts,subdm,&U_c);CHKERRQ(ierr);
601*d2567f34SHong Zhang 
602*d2567f34SHong Zhang   ierr = VecScatterBegin(gscat,U,U_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
603*d2567f34SHong Zhang   ierr = VecScatterEnd(gscat,U,U_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
604*d2567f34SHong Zhang 
605*d2567f34SHong Zhang   ierr = TSIRKRestoreVecs(ts,dm,&U);CHKERRQ(ierr);
606*d2567f34SHong Zhang   ierr = TSIRKRestoreVecs(ts,subdm,&U_c);CHKERRQ(ierr);
607*d2567f34SHong Zhang   PetscFunctionReturn(0);
608*d2567f34SHong Zhang }
609*d2567f34SHong Zhang 
610*d2567f34SHong Zhang static PetscErrorCode TSSetUp_IRK(TS ts)
611*d2567f34SHong Zhang {
612*d2567f34SHong Zhang   TS_IRK         *irk = (TS_IRK*)ts->data;
613*d2567f34SHong Zhang   IRKTableau     tab = irk->tableau;
614*d2567f34SHong Zhang   DM             dm;
615*d2567f34SHong Zhang   Mat            J;
616*d2567f34SHong Zhang   Vec            R;
617*d2567f34SHong Zhang   const PetscInt nstages = irk->nstages;
618*d2567f34SHong Zhang   PetscInt       vsize,bs;
619*d2567f34SHong Zhang   PetscErrorCode ierr;
620*d2567f34SHong Zhang 
621*d2567f34SHong Zhang   PetscFunctionBegin;
622*d2567f34SHong Zhang   if (!irk->work) {
623*d2567f34SHong Zhang     ierr = PetscMalloc1(irk->nstages,&irk->work);CHKERRQ(ierr);
624*d2567f34SHong Zhang   }
625*d2567f34SHong Zhang   if (!irk->Y) {
626*d2567f34SHong Zhang     ierr = VecDuplicateVecs(ts->vec_sol,irk->nstages,&irk->Y);CHKERRQ(ierr);
627*d2567f34SHong Zhang   }
628*d2567f34SHong Zhang   if (!irk->YdotI) {
629*d2567f34SHong Zhang     ierr = VecDuplicateVecs(ts->vec_sol,irk->nstages,&irk->YdotI);CHKERRQ(ierr);
630*d2567f34SHong Zhang   }
631*d2567f34SHong Zhang   if (!irk->Ydot) {
632*d2567f34SHong Zhang     ierr = VecDuplicate(ts->vec_sol,&irk->Ydot);CHKERRQ(ierr);
633*d2567f34SHong Zhang   }
634*d2567f34SHong Zhang   if (!irk->U) {
635*d2567f34SHong Zhang     ierr = VecDuplicate(ts->vec_sol,&irk->U);CHKERRQ(ierr);
636*d2567f34SHong Zhang   }
637*d2567f34SHong Zhang   if (!irk->U0) {
638*d2567f34SHong Zhang     ierr = VecDuplicate(ts->vec_sol,&irk->U0);CHKERRQ(ierr);
639*d2567f34SHong Zhang   }
640*d2567f34SHong Zhang   if (!irk->Z) {
641*d2567f34SHong Zhang     ierr = VecCreate(PetscObjectComm((PetscObject)ts->vec_sol),&irk->Z);CHKERRQ(ierr);
642*d2567f34SHong Zhang     ierr = VecGetSize(ts->vec_sol,&vsize);CHKERRQ(ierr);
643*d2567f34SHong Zhang     ierr = VecSetSizes(irk->Z,PETSC_DECIDE,vsize*irk->nstages);CHKERRQ(ierr);
644*d2567f34SHong Zhang     ierr = VecGetBlockSize(ts->vec_sol,&bs);CHKERRQ(ierr);
645*d2567f34SHong Zhang     ierr = VecSetBlockSize(irk->Z,irk->nstages*bs);CHKERRQ(ierr);
646*d2567f34SHong Zhang     ierr = VecSetFromOptions(irk->Z);CHKERRQ(ierr);
647*d2567f34SHong Zhang   }
648*d2567f34SHong Zhang   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
649*d2567f34SHong Zhang   ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSIRK,DMRestrictHook_TSIRK,ts);CHKERRQ(ierr);
650*d2567f34SHong Zhang   ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSIRK,DMSubDomainRestrictHook_TSIRK,ts);CHKERRQ(ierr);
651*d2567f34SHong Zhang 
652*d2567f34SHong Zhang   ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr);
653*d2567f34SHong Zhang   ierr = VecDuplicate(irk->Z,&R);CHKERRQ(ierr);
654*d2567f34SHong Zhang   ierr = SNESSetFunction(ts->snes,R,SNESTSFormFunction,ts);CHKERRQ(ierr);
655*d2567f34SHong Zhang   ierr = SNESGetJacobian(ts->snes,&J,NULL,NULL,NULL);CHKERRQ(ierr);
656*d2567f34SHong Zhang   if (!irk->TJ) {
657*d2567f34SHong Zhang     /* Create the KAIJ matrix for solving the stages */
658*d2567f34SHong Zhang     ierr = MatCreateKAIJ(J,nstages,nstages,tab->A_inv,tab->I_s,&irk->TJ);CHKERRQ(ierr);
659*d2567f34SHong Zhang   }
660*d2567f34SHong Zhang   ierr = SNESSetJacobian(ts->snes,irk->TJ,irk->TJ,SNESTSFormJacobian,ts);CHKERRQ(ierr);
661*d2567f34SHong Zhang   ierr = VecDestroy(&R);CHKERRQ(ierr);
662*d2567f34SHong Zhang   PetscFunctionReturn(0);
663*d2567f34SHong Zhang }
664*d2567f34SHong Zhang 
665*d2567f34SHong Zhang static PetscErrorCode TSSetFromOptions_IRK(PetscOptionItems *PetscOptionsObject,TS ts)
666*d2567f34SHong Zhang {
667*d2567f34SHong Zhang   TS_IRK         *irk = (TS_IRK*)ts->data;
668*d2567f34SHong Zhang   char           tname[256] = TSIRKGAUSS;
669*d2567f34SHong Zhang   PetscErrorCode ierr;
670*d2567f34SHong Zhang 
671*d2567f34SHong Zhang   PetscFunctionBegin;
672*d2567f34SHong Zhang   ierr = PetscOptionsHead(PetscOptionsObject,"IRK ODE solver options");CHKERRQ(ierr);
673*d2567f34SHong Zhang   {
674*d2567f34SHong Zhang     PetscBool flg1,flg2;
675*d2567f34SHong Zhang     ierr = PetscOptionsInt("-ts_irk_nstages","Stages of the IRK method","TSIRKSetNumStages",irk->nstages,&irk->nstages,&flg1);CHKERRQ(ierr);
676*d2567f34SHong Zhang     ierr = PetscOptionsFList("-ts_irk_type","Type of IRK method","TSIRKSetType",TSIRKList,irk->method_name[0] ? irk->method_name : tname,tname,sizeof(tname),&flg2);CHKERRQ(ierr);
677*d2567f34SHong Zhang     if (flg1 ||flg2 || !irk->method_name[0]) { /* Create the method tableau after nstages or method is set */
678*d2567f34SHong Zhang       ierr = TSIRKSetType(ts,tname);CHKERRQ(ierr);
679*d2567f34SHong Zhang     }
680*d2567f34SHong Zhang   }
681*d2567f34SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
682*d2567f34SHong Zhang   PetscFunctionReturn(0);
683*d2567f34SHong Zhang }
684*d2567f34SHong Zhang 
685*d2567f34SHong Zhang static PetscErrorCode TSView_IRK(TS ts,PetscViewer viewer)
686*d2567f34SHong Zhang {
687*d2567f34SHong Zhang   TS_IRK         *irk = (TS_IRK*)ts->data;
688*d2567f34SHong Zhang   PetscBool      iascii;
689*d2567f34SHong Zhang   PetscErrorCode ierr;
690*d2567f34SHong Zhang 
691*d2567f34SHong Zhang   PetscFunctionBegin;
692*d2567f34SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
693*d2567f34SHong Zhang   if (iascii) {
694*d2567f34SHong Zhang     IRKTableau    tab = irk->tableau;
695*d2567f34SHong Zhang     TSIRKType irktype;
696*d2567f34SHong Zhang     char          buf[512];
697*d2567f34SHong Zhang 
698*d2567f34SHong Zhang     ierr = TSIRKGetType(ts,&irktype);CHKERRQ(ierr);
699*d2567f34SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  IRK type %s\n",irktype);CHKERRQ(ierr);
700*d2567f34SHong Zhang     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",irk->nstages,tab->c);CHKERRQ(ierr);
701*d2567f34SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa       c = %s\n",buf);CHKERRQ(ierr);
702*d2567f34SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s\n",irk->stiffly_accurate ? "yes" : "no");CHKERRQ(ierr);
703*d2567f34SHong Zhang     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",PetscSqr(irk->nstages),tab->A);CHKERRQ(ierr);
704*d2567f34SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  A coefficients       A = %s\n",buf);CHKERRQ(ierr);
705*d2567f34SHong Zhang   }
706*d2567f34SHong Zhang   PetscFunctionReturn(0);
707*d2567f34SHong Zhang }
708*d2567f34SHong Zhang 
709*d2567f34SHong Zhang static PetscErrorCode TSLoad_IRK(TS ts,PetscViewer viewer)
710*d2567f34SHong Zhang {
711*d2567f34SHong Zhang   PetscErrorCode ierr;
712*d2567f34SHong Zhang   SNES           snes;
713*d2567f34SHong Zhang   TSAdapt        adapt;
714*d2567f34SHong Zhang 
715*d2567f34SHong Zhang   PetscFunctionBegin;
716*d2567f34SHong Zhang   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
717*d2567f34SHong Zhang   ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr);
718*d2567f34SHong Zhang   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
719*d2567f34SHong Zhang   ierr = SNESLoad(snes,viewer);CHKERRQ(ierr);
720*d2567f34SHong Zhang   /* function and Jacobian context for SNES when used with TS is always ts object */
721*d2567f34SHong Zhang   ierr = SNESSetFunction(snes,NULL,NULL,ts);CHKERRQ(ierr);
722*d2567f34SHong Zhang   ierr = SNESSetJacobian(snes,NULL,NULL,NULL,ts);CHKERRQ(ierr);
723*d2567f34SHong Zhang   PetscFunctionReturn(0);
724*d2567f34SHong Zhang }
725*d2567f34SHong Zhang 
726*d2567f34SHong Zhang /*@C
727*d2567f34SHong Zhang   TSIRKSetType - Set the type of IRK scheme
728*d2567f34SHong Zhang 
729*d2567f34SHong Zhang   Logically collective
730*d2567f34SHong Zhang 
731*d2567f34SHong Zhang   Input Parameter:
732*d2567f34SHong Zhang +  ts - timestepping context
733*d2567f34SHong Zhang -  irktype - type of IRK scheme
734*d2567f34SHong Zhang 
735*d2567f34SHong Zhang   Options Database:
736*d2567f34SHong Zhang .  -ts_irk_type <gauss>
737*d2567f34SHong Zhang 
738*d2567f34SHong Zhang   Level: intermediate
739*d2567f34SHong Zhang 
740*d2567f34SHong Zhang .seealso: TSIRKGetType(), TSIRK, TSIRKType, TSIRKGAUSS
741*d2567f34SHong Zhang @*/
742*d2567f34SHong Zhang PetscErrorCode TSIRKSetType(TS ts,TSIRKType irktype)
743*d2567f34SHong Zhang {
744*d2567f34SHong Zhang   PetscErrorCode ierr;
745*d2567f34SHong Zhang 
746*d2567f34SHong Zhang   PetscFunctionBegin;
747*d2567f34SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
748*d2567f34SHong Zhang   PetscValidCharPointer(irktype,2);
749*d2567f34SHong Zhang   ierr = PetscTryMethod(ts,"TSIRKSetType_C",(TS,TSIRKType),(ts,irktype));CHKERRQ(ierr);
750*d2567f34SHong Zhang   PetscFunctionReturn(0);
751*d2567f34SHong Zhang }
752*d2567f34SHong Zhang 
753*d2567f34SHong Zhang /*@C
754*d2567f34SHong Zhang   TSIRKGetType - Get the type of IRK IMEX scheme
755*d2567f34SHong Zhang 
756*d2567f34SHong Zhang   Logically collective
757*d2567f34SHong Zhang 
758*d2567f34SHong Zhang   Input Parameter:
759*d2567f34SHong Zhang .  ts - timestepping context
760*d2567f34SHong Zhang 
761*d2567f34SHong Zhang   Output Parameter:
762*d2567f34SHong Zhang .  irktype - type of IRK-IMEX scheme
763*d2567f34SHong Zhang 
764*d2567f34SHong Zhang   Level: intermediate
765*d2567f34SHong Zhang 
766*d2567f34SHong Zhang .seealso: TSIRKGetType()
767*d2567f34SHong Zhang @*/
768*d2567f34SHong Zhang PetscErrorCode TSIRKGetType(TS ts,TSIRKType *irktype)
769*d2567f34SHong Zhang {
770*d2567f34SHong Zhang   PetscErrorCode ierr;
771*d2567f34SHong Zhang 
772*d2567f34SHong Zhang   PetscFunctionBegin;
773*d2567f34SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
774*d2567f34SHong Zhang   ierr = PetscUseMethod(ts,"TSIRKGetType_C",(TS,TSIRKType*),(ts,irktype));CHKERRQ(ierr);
775*d2567f34SHong Zhang   PetscFunctionReturn(0);
776*d2567f34SHong Zhang }
777*d2567f34SHong Zhang 
778*d2567f34SHong Zhang /*@C
779*d2567f34SHong Zhang   TSIRKSetNumStages - Set the number of stages of IRK scheme
780*d2567f34SHong Zhang 
781*d2567f34SHong Zhang   Logically collective
782*d2567f34SHong Zhang 
783*d2567f34SHong Zhang   Input Parameter:
784*d2567f34SHong Zhang +  ts - timestepping context
785*d2567f34SHong Zhang -  nstages - number of stages of IRK scheme
786*d2567f34SHong Zhang 
787*d2567f34SHong Zhang   Options Database:
788*d2567f34SHong Zhang .  -ts_irk_nstages <int>: Number of stages
789*d2567f34SHong Zhang 
790*d2567f34SHong Zhang   Level: intermediate
791*d2567f34SHong Zhang 
792*d2567f34SHong Zhang .seealso: TSIRKGetNumStages(), TSIRK
793*d2567f34SHong Zhang @*/
794*d2567f34SHong Zhang PetscErrorCode TSIRKSetNumStages(TS ts,PetscInt nstages)
795*d2567f34SHong Zhang {
796*d2567f34SHong Zhang   PetscErrorCode ierr;
797*d2567f34SHong Zhang 
798*d2567f34SHong Zhang   PetscFunctionBegin;
799*d2567f34SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
800*d2567f34SHong Zhang   ierr = PetscTryMethod(ts,"TSIRKSetNumStages_C",(TS,PetscInt),(ts,nstages));CHKERRQ(ierr);
801*d2567f34SHong Zhang   PetscFunctionReturn(0);
802*d2567f34SHong Zhang }
803*d2567f34SHong Zhang 
804*d2567f34SHong Zhang /*@C
805*d2567f34SHong Zhang   TSIRKGetNumStages - Get the number of stages of IRK scheme
806*d2567f34SHong Zhang 
807*d2567f34SHong Zhang   Logically collective
808*d2567f34SHong Zhang 
809*d2567f34SHong Zhang   Input Parameter:
810*d2567f34SHong Zhang +  ts - timestepping context
811*d2567f34SHong Zhang -  nstages - number of stages of IRK scheme
812*d2567f34SHong Zhang 
813*d2567f34SHong Zhang   Level: intermediate
814*d2567f34SHong Zhang 
815*d2567f34SHong Zhang .seealso: TSIRKSetNumStages(), TSIRK
816*d2567f34SHong Zhang @*/
817*d2567f34SHong Zhang PetscErrorCode TSIRKGetNumStages(TS ts,PetscInt *nstages)
818*d2567f34SHong Zhang {
819*d2567f34SHong Zhang   PetscErrorCode ierr;
820*d2567f34SHong Zhang 
821*d2567f34SHong Zhang   PetscFunctionBegin;
822*d2567f34SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
823*d2567f34SHong Zhang   PetscValidIntPointer(nstages,2);
824*d2567f34SHong Zhang   ierr = PetscTryMethod(ts,"TSIRKGetNumStages_C",(TS,PetscInt*),(ts,nstages));CHKERRQ(ierr);
825*d2567f34SHong Zhang   PetscFunctionReturn(0);
826*d2567f34SHong Zhang }
827*d2567f34SHong Zhang 
828*d2567f34SHong Zhang static PetscErrorCode TSIRKGetType_IRK(TS ts,TSIRKType *irktype)
829*d2567f34SHong Zhang {
830*d2567f34SHong Zhang   TS_IRK *irk = (TS_IRK*)ts->data;
831*d2567f34SHong Zhang 
832*d2567f34SHong Zhang   PetscFunctionBegin;
833*d2567f34SHong Zhang   *irktype = irk->method_name;
834*d2567f34SHong Zhang   PetscFunctionReturn(0);
835*d2567f34SHong Zhang }
836*d2567f34SHong Zhang 
837*d2567f34SHong Zhang static PetscErrorCode TSIRKSetType_IRK(TS ts,TSIRKType irktype)
838*d2567f34SHong Zhang {
839*d2567f34SHong Zhang   TS_IRK         *irk = (TS_IRK*)ts->data;
840*d2567f34SHong Zhang   PetscErrorCode ierr,(*irkcreate)(TS);
841*d2567f34SHong Zhang 
842*d2567f34SHong Zhang   PetscFunctionBegin;
843*d2567f34SHong Zhang   if (irk->method_name) {
844*d2567f34SHong Zhang     ierr = PetscFree(irk->method_name);CHKERRQ(ierr);
845*d2567f34SHong Zhang     ierr = TSIRKTableauReset(ts);CHKERRQ(ierr);
846*d2567f34SHong Zhang   }
847*d2567f34SHong Zhang   ierr = PetscFunctionListFind(TSIRKList,irktype,&irkcreate);CHKERRQ(ierr);
848*d2567f34SHong Zhang   if (!irkcreate) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSIRK type \"%s\" given",irktype);
849*d2567f34SHong Zhang   ierr = (*irkcreate)(ts);CHKERRQ(ierr);
850*d2567f34SHong Zhang   ierr = PetscStrallocpy(irktype,&irk->method_name);CHKERRQ(ierr);
851*d2567f34SHong Zhang   PetscFunctionReturn(0);
852*d2567f34SHong Zhang }
853*d2567f34SHong Zhang 
854*d2567f34SHong Zhang static PetscErrorCode TSIRKSetNumStages_IRK(TS ts,PetscInt nstages)
855*d2567f34SHong Zhang {
856*d2567f34SHong Zhang   TS_IRK *irk = (TS_IRK*)ts->data;
857*d2567f34SHong Zhang 
858*d2567f34SHong Zhang   PetscFunctionBegin;
859*d2567f34SHong Zhang   if (nstages<=0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"input argument, %d, out of range",nstages);
860*d2567f34SHong Zhang   irk->nstages = nstages;
861*d2567f34SHong Zhang   PetscFunctionReturn(0);
862*d2567f34SHong Zhang }
863*d2567f34SHong Zhang 
864*d2567f34SHong Zhang static PetscErrorCode TSIRKGetNumStages_IRK(TS ts,PetscInt *nstages)
865*d2567f34SHong Zhang {
866*d2567f34SHong Zhang   TS_IRK *irk = (TS_IRK*)ts->data;
867*d2567f34SHong Zhang 
868*d2567f34SHong Zhang   PetscFunctionBegin;
869*d2567f34SHong Zhang   PetscValidIntPointer(nstages,2);
870*d2567f34SHong Zhang   *nstages = irk->nstages;
871*d2567f34SHong Zhang   PetscFunctionReturn(0);
872*d2567f34SHong Zhang }
873*d2567f34SHong Zhang 
874*d2567f34SHong Zhang static PetscErrorCode TSDestroy_IRK(TS ts)
875*d2567f34SHong Zhang {
876*d2567f34SHong Zhang   PetscErrorCode ierr;
877*d2567f34SHong Zhang 
878*d2567f34SHong Zhang   PetscFunctionBegin;
879*d2567f34SHong Zhang   ierr = TSReset_IRK(ts);CHKERRQ(ierr);
880*d2567f34SHong Zhang   if (ts->dm) {
881*d2567f34SHong Zhang     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSIRK,DMRestrictHook_TSIRK,ts);CHKERRQ(ierr);
882*d2567f34SHong Zhang     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSIRK,DMSubDomainRestrictHook_TSIRK,ts);CHKERRQ(ierr);
883*d2567f34SHong Zhang   }
884*d2567f34SHong Zhang   ierr = PetscFree(ts->data);CHKERRQ(ierr);
885*d2567f34SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKSetType_C",NULL);CHKERRQ(ierr);
886*d2567f34SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKGetType_C",NULL);CHKERRQ(ierr);
887*d2567f34SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKSetNumStages_C",NULL);CHKERRQ(ierr);
888*d2567f34SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKGetNumStages_C",NULL);CHKERRQ(ierr);
889*d2567f34SHong Zhang   PetscFunctionReturn(0);
890*d2567f34SHong Zhang }
891*d2567f34SHong Zhang 
892*d2567f34SHong Zhang /*MC
893*d2567f34SHong Zhang       TSIRK - ODE and DAE solver using Implicit Runge-Kutta schemes
894*d2567f34SHong Zhang 
895*d2567f34SHong Zhang   Notes:
896*d2567f34SHong Zhang 
897*d2567f34SHong Zhang   TSIRK uses the sparse Kronecker product matrix implementation of MATKAIJ to achieve good arithmetic intensity.
898*d2567f34SHong Zhang 
899*d2567f34SHong Zhang   The default is TSIRK3, it can be changed with TSIRKSetType() or -ts_irk_type
900*d2567f34SHong Zhang 
901*d2567f34SHong Zhang   If the equation is implicit or a DAE, then TSSetEquationType() needs to be set accordingly. Refer to the manual for further information.
902*d2567f34SHong Zhang 
903*d2567f34SHong Zhang   Consider trying TSROSW if the stiff part is linear or weakly nonlinear.
904*d2567f34SHong Zhang 
905*d2567f34SHong Zhang   Level: beginner
906*d2567f34SHong Zhang 
907*d2567f34SHong Zhang .seealso:  TSCreate(), TS, TSSetType(), TSIRKSetType(), TSIRKGetType(),
908*d2567f34SHong Zhang            TSIRK1BEE, TSIRK2C, TSIRK2D, TSIRK2E, TSIRK3, TSIRKL2, TSIRKA2, TSIRKARS122,
909*d2567f34SHong Zhang            TSIRK4, TSIRK5, TSIRKPRSSP2, TSIRKARS443, TSIRKBPR3, TSIRKType, TSIRKRegister()
910*d2567f34SHong Zhang 
911*d2567f34SHong Zhang M*/
912*d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSCreate_IRK(TS ts)
913*d2567f34SHong Zhang {
914*d2567f34SHong Zhang   TS_IRK         *irk;
915*d2567f34SHong Zhang   PetscErrorCode ierr;
916*d2567f34SHong Zhang 
917*d2567f34SHong Zhang   PetscFunctionBegin;
918*d2567f34SHong Zhang   ierr = TSIRKInitializePackage();CHKERRQ(ierr);
919*d2567f34SHong Zhang 
920*d2567f34SHong Zhang   ts->ops->reset          = TSReset_IRK;
921*d2567f34SHong Zhang   ts->ops->destroy        = TSDestroy_IRK;
922*d2567f34SHong Zhang   ts->ops->view           = TSView_IRK;
923*d2567f34SHong Zhang   ts->ops->load           = TSLoad_IRK;
924*d2567f34SHong Zhang   ts->ops->setup          = TSSetUp_IRK;
925*d2567f34SHong Zhang   ts->ops->step           = TSStep_IRK;
926*d2567f34SHong Zhang   ts->ops->interpolate    = TSInterpolate_IRK;
927*d2567f34SHong Zhang   ts->ops->evaluatestep   = TSEvaluateStep_IRK;
928*d2567f34SHong Zhang   ts->ops->rollback       = TSRollBack_IRK;
929*d2567f34SHong Zhang   ts->ops->setfromoptions = TSSetFromOptions_IRK;
930*d2567f34SHong Zhang   ts->ops->snesfunction   = SNESTSFormFunction_IRK;
931*d2567f34SHong Zhang   ts->ops->snesjacobian   = SNESTSFormJacobian_IRK;
932*d2567f34SHong Zhang 
933*d2567f34SHong Zhang   ts->usessnes = PETSC_TRUE;
934*d2567f34SHong Zhang 
935*d2567f34SHong Zhang   ierr = PetscNewLog(ts,&irk);CHKERRQ(ierr);
936*d2567f34SHong Zhang   ts->data = (void*)irk;
937*d2567f34SHong Zhang 
938*d2567f34SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKSetType_C",TSIRKSetType_IRK);CHKERRQ(ierr);
939*d2567f34SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKGetType_C",TSIRKGetType_IRK);CHKERRQ(ierr);
940*d2567f34SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKSetNumStages_C",TSIRKSetNumStages_IRK);CHKERRQ(ierr);
941*d2567f34SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKGetNumStages_C",TSIRKGetNumStages_IRK);CHKERRQ(ierr);
942*d2567f34SHong Zhang   /* 3-stage IRK_Gauss is the default */
943*d2567f34SHong Zhang   ierr = PetscNew(&irk->tableau);CHKERRQ(ierr);
944*d2567f34SHong Zhang   irk->nstages = 3;
945*d2567f34SHong Zhang   ierr = TSIRKSetType(ts,TSIRKDefault);CHKERRQ(ierr);
946*d2567f34SHong Zhang   PetscFunctionReturn(0);
947*d2567f34SHong Zhang }
948