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