xref: /petsc/src/ts/impls/implicit/irk/irk.c (revision 1baa6e3354bfe224b33a0c158f482508792a8a6e)
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 
65db781477SPatrick Sanan .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;
749566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(PetscSqr(nstages),&tab->A,PetscSqr(nstages),&tab->A_inv,PetscSqr(nstages),&tab->I_s));
759566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(nstages,&tab->b,nstages,&tab->c,nstages,&tab->binterp,nstages,&tab->A_inv_rowsum));
769566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(tab->A,A,PetscSqr(nstages)));
779566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(tab->b,b,nstages));
789566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(tab->c,c,nstages));
79d2567f34SHong Zhang   /* optional coefficient arrays */
80*1baa6e33SBarry Smith   if (binterp) PetscCall(PetscArraycpy(tab->binterp,binterp,nstages));
81*1baa6e33SBarry Smith   if (A_inv) PetscCall(PetscArraycpy(tab->A_inv,A_inv,PetscSqr(nstages)));
82*1baa6e33SBarry Smith   if (A_inv_rowsum) PetscCall(PetscArraycpy(tab->A_inv_rowsum,A_inv_rowsum,nstages));
83*1baa6e33SBarry Smith   if (I_s) PetscCall(PetscArraycpy(tab->I_s,I_s,PetscSqr(nstages)));
84d2567f34SHong Zhang   PetscFunctionReturn(0);
85d2567f34SHong Zhang }
86d2567f34SHong Zhang 
87d2567f34SHong Zhang /* Arrays should be freed with PetscFree3(A,b,c) */
88d2567f34SHong Zhang static PetscErrorCode TSIRKCreate_Gauss(TS ts)
89d2567f34SHong Zhang {
90d2567f34SHong Zhang   PetscInt       nstages;
91d2567f34SHong Zhang   PetscReal      *gauss_A_real,*gauss_b,*b,*gauss_c;
92d2567f34SHong Zhang   PetscScalar    *gauss_A,*gauss_A_inv,*gauss_A_inv_rowsum,*I_s;
93d2567f34SHong Zhang   PetscScalar    *G0,*G1;
94d2567f34SHong Zhang   PetscInt       i,j;
95d2567f34SHong Zhang   Mat            G0mat,G1mat,Amat;
96d2567f34SHong Zhang 
97d2567f34SHong Zhang   PetscFunctionBegin;
989566063dSJacob Faibussowitsch   PetscCall(TSIRKGetNumStages(ts,&nstages));
999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(PetscSqr(nstages),&gauss_A_real,nstages,&gauss_b,nstages,&gauss_c));
1009566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(PetscSqr(nstages),&gauss_A,PetscSqr(nstages),&gauss_A_inv,nstages,&gauss_A_inv_rowsum,PetscSqr(nstages),&I_s));
1019566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(nstages,&b,PetscSqr(nstages),&G0,PetscSqr(nstages),&G1));
1029566063dSJacob Faibussowitsch   PetscCall(PetscDTGaussQuadrature(nstages,0.,1.,gauss_c,b));
103d2567f34SHong Zhang   for (i=0; i<nstages; i++) gauss_b[i] = b[i]; /* copy to possibly-complex array */
104d2567f34SHong Zhang 
105d2567f34SHong Zhang   /* A^T = G0^{-1} G1 */
106d2567f34SHong Zhang   for (i=0; i<nstages; i++) {
107d2567f34SHong Zhang     for (j=0; j<nstages; j++) {
108d2567f34SHong Zhang       G0[i*nstages+j] = PetscPowRealInt(gauss_c[i],j);
109d2567f34SHong Zhang       G1[i*nstages+j] = PetscPowRealInt(gauss_c[i],j+1)/(j+1);
110d2567f34SHong Zhang     }
111d2567f34SHong Zhang   }
112d2567f34SHong Zhang   /* The arrays above are row-aligned, but we create dense matrices as the transpose */
1139566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,nstages,nstages,G0,&G0mat));
1149566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,nstages,nstages,G1,&G1mat));
1159566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,nstages,nstages,gauss_A,&Amat));
1169566063dSJacob Faibussowitsch   PetscCall(MatLUFactor(G0mat,NULL,NULL,NULL));
1179566063dSJacob Faibussowitsch   PetscCall(MatMatSolve(G0mat,G1mat,Amat));
1189566063dSJacob Faibussowitsch   PetscCall(MatTranspose(Amat,MAT_INPLACE_MATRIX,&Amat));
119d2567f34SHong Zhang   for (i=0; i<nstages; i++)
120d2567f34SHong Zhang     for (j=0; j<nstages; j++)
121d2567f34SHong Zhang       gauss_A_real[i*nstages+j] = PetscRealPart(gauss_A[i*nstages+j]);
122d2567f34SHong Zhang 
1239566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&G0mat));
1249566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&G1mat));
1259566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Amat));
1269566063dSJacob Faibussowitsch   PetscCall(PetscFree3(b,G0,G1));
127d2567f34SHong Zhang 
128d2567f34SHong Zhang   {/* Invert A */
129d2567f34SHong Zhang     /* PETSc does not provide a routine to calculate the inverse of a general matrix.
130d2567f34SHong Zhang      * To get the inverse of A, we form a sequential BAIJ matrix from it, consisting of a single block with block size
131d2567f34SHong Zhang      * equal to the dimension of A, and then use MatInvertBlockDiagonal(). */
132d2567f34SHong Zhang     Mat               A_baij;
133d2567f34SHong Zhang     PetscInt          idxm[1]={0},idxn[1]={0};
134d2567f34SHong Zhang     const PetscScalar *A_inv;
135d2567f34SHong Zhang 
1369566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF,nstages,nstages,nstages,1,NULL,&A_baij));
1379566063dSJacob Faibussowitsch     PetscCall(MatSetOption(A_baij,MAT_ROW_ORIENTED,PETSC_FALSE));
1389566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(A_baij,1,idxm,1,idxn,gauss_A,INSERT_VALUES));
1399566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A_baij,MAT_FINAL_ASSEMBLY));
1409566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A_baij,MAT_FINAL_ASSEMBLY));
1419566063dSJacob Faibussowitsch     PetscCall(MatInvertBlockDiagonal(A_baij,&A_inv));
1429566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(gauss_A_inv,A_inv,nstages*nstages*sizeof(PetscScalar)));
1439566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A_baij));
144d2567f34SHong Zhang   }
145d2567f34SHong Zhang 
146d2567f34SHong Zhang   /* Compute row sums A_inv_rowsum and identity I_s */
147d2567f34SHong Zhang   for (i=0; i<nstages; i++) {
148d2567f34SHong Zhang     gauss_A_inv_rowsum[i] = 0;
149d2567f34SHong Zhang     for (j=0; j<nstages; j++) {
150d2567f34SHong Zhang       gauss_A_inv_rowsum[i] += gauss_A_inv[i+nstages*j];
151d2567f34SHong Zhang       I_s[i+nstages*j] = 1.*(i == j);
152d2567f34SHong Zhang     }
153d2567f34SHong Zhang   }
1549566063dSJacob Faibussowitsch   PetscCall(TSIRKTableauCreate(ts,nstages,gauss_A_real,gauss_b,gauss_c,NULL,gauss_A_inv,gauss_A_inv_rowsum,I_s));
1559566063dSJacob Faibussowitsch   PetscCall(PetscFree3(gauss_A_real,gauss_b,gauss_c));
1569566063dSJacob Faibussowitsch   PetscCall(PetscFree4(gauss_A,gauss_A_inv,gauss_A_inv_rowsum,I_s));
157d2567f34SHong Zhang   PetscFunctionReturn(0);
158d2567f34SHong Zhang }
159d2567f34SHong Zhang 
160d2567f34SHong Zhang /*@C
161d2567f34SHong Zhang    TSIRKRegister -  adds a TSIRK implementation
162d2567f34SHong Zhang 
163d2567f34SHong Zhang    Not Collective
164d2567f34SHong Zhang 
165d2567f34SHong Zhang    Input Parameters:
166d2567f34SHong Zhang +  sname - name of user-defined IRK scheme
167d2567f34SHong Zhang -  function - function to create method context
168d2567f34SHong Zhang 
169d2567f34SHong Zhang    Notes:
170d2567f34SHong Zhang    TSIRKRegister() may be called multiple times to add several user-defined families.
171d2567f34SHong Zhang 
172d2567f34SHong Zhang    Sample usage:
173d2567f34SHong Zhang .vb
174d2567f34SHong Zhang    TSIRKRegister("my_scheme",MySchemeCreate);
175d2567f34SHong Zhang .ve
176d2567f34SHong Zhang 
177d2567f34SHong Zhang    Then, your scheme can be chosen with the procedural interface via
178d2567f34SHong Zhang $     TSIRKSetType(ts,"my_scheme")
179d2567f34SHong Zhang    or at runtime via the option
180d2567f34SHong Zhang $     -ts_irk_type my_scheme
181d2567f34SHong Zhang 
182d2567f34SHong Zhang    Level: advanced
183d2567f34SHong Zhang 
184db781477SPatrick Sanan .seealso: `TSIRKRegisterAll()`
185d2567f34SHong Zhang @*/
186d2567f34SHong Zhang PetscErrorCode TSIRKRegister(const char sname[],PetscErrorCode (*function)(TS))
187d2567f34SHong Zhang {
188d2567f34SHong Zhang   PetscFunctionBegin;
1899566063dSJacob Faibussowitsch   PetscCall(TSIRKInitializePackage());
1909566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&TSIRKList,sname,function));
191d2567f34SHong Zhang   PetscFunctionReturn(0);
192d2567f34SHong Zhang }
193d2567f34SHong Zhang 
194d2567f34SHong Zhang /*@C
195d2567f34SHong Zhang   TSIRKRegisterAll - Registers all of the implicit Runge-Kutta methods in TSIRK
196d2567f34SHong Zhang 
197d2567f34SHong Zhang   Not Collective, but should be called by all processes which will need the schemes to be registered
198d2567f34SHong Zhang 
199d2567f34SHong Zhang   Level: advanced
200d2567f34SHong Zhang 
201db781477SPatrick Sanan .seealso: `TSIRKRegisterDestroy()`
202d2567f34SHong Zhang @*/
203d2567f34SHong Zhang PetscErrorCode TSIRKRegisterAll(void)
204d2567f34SHong Zhang {
205d2567f34SHong Zhang   PetscFunctionBegin;
206d2567f34SHong Zhang   if (TSIRKRegisterAllCalled) PetscFunctionReturn(0);
207d2567f34SHong Zhang   TSIRKRegisterAllCalled = PETSC_TRUE;
208d2567f34SHong Zhang 
2099566063dSJacob Faibussowitsch   PetscCall(TSIRKRegister(TSIRKGAUSS,TSIRKCreate_Gauss));
210d2567f34SHong Zhang   PetscFunctionReturn(0);
211d2567f34SHong Zhang }
212d2567f34SHong Zhang 
213d2567f34SHong Zhang /*@C
214d2567f34SHong Zhang    TSIRKRegisterDestroy - Frees the list of schemes that were registered by TSIRKRegister().
215d2567f34SHong Zhang 
216d2567f34SHong Zhang    Not Collective
217d2567f34SHong Zhang 
218d2567f34SHong Zhang    Level: advanced
219d2567f34SHong Zhang 
220db781477SPatrick Sanan .seealso: `TSIRKRegister()`, `TSIRKRegisterAll()`
221d2567f34SHong Zhang @*/
222d2567f34SHong Zhang PetscErrorCode TSIRKRegisterDestroy(void)
223d2567f34SHong Zhang {
224d2567f34SHong Zhang   PetscFunctionBegin;
225d2567f34SHong Zhang   TSIRKRegisterAllCalled = PETSC_FALSE;
226d2567f34SHong Zhang   PetscFunctionReturn(0);
227d2567f34SHong Zhang }
228d2567f34SHong Zhang 
229d2567f34SHong Zhang /*@C
230d2567f34SHong Zhang   TSIRKInitializePackage - This function initializes everything in the TSIRK package. It is called
231d2567f34SHong Zhang   from TSInitializePackage().
232d2567f34SHong Zhang 
233d2567f34SHong Zhang   Level: developer
234d2567f34SHong Zhang 
235db781477SPatrick Sanan .seealso: `PetscInitialize()`
236d2567f34SHong Zhang @*/
237d2567f34SHong Zhang PetscErrorCode TSIRKInitializePackage(void)
238d2567f34SHong Zhang {
239d2567f34SHong Zhang   PetscFunctionBegin;
240d2567f34SHong Zhang   if (TSIRKPackageInitialized) PetscFunctionReturn(0);
241d2567f34SHong Zhang   TSIRKPackageInitialized = PETSC_TRUE;
2429566063dSJacob Faibussowitsch   PetscCall(TSIRKRegisterAll());
2439566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(TSIRKFinalizePackage));
244d2567f34SHong Zhang   PetscFunctionReturn(0);
245d2567f34SHong Zhang }
246d2567f34SHong Zhang 
247d2567f34SHong Zhang /*@C
248d2567f34SHong Zhang   TSIRKFinalizePackage - This function destroys everything in the TSIRK package. It is
249d2567f34SHong Zhang   called from PetscFinalize().
250d2567f34SHong Zhang 
251d2567f34SHong Zhang   Level: developer
252d2567f34SHong Zhang 
253db781477SPatrick Sanan .seealso: `PetscFinalize()`
254d2567f34SHong Zhang @*/
255d2567f34SHong Zhang PetscErrorCode TSIRKFinalizePackage(void)
256d2567f34SHong Zhang {
257d2567f34SHong Zhang   PetscFunctionBegin;
2589566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&TSIRKList));
259d2567f34SHong Zhang   TSIRKPackageInitialized = PETSC_FALSE;
260d2567f34SHong Zhang   PetscFunctionReturn(0);
261d2567f34SHong Zhang }
262d2567f34SHong Zhang 
263d2567f34SHong Zhang /*
264d2567f34SHong Zhang  This function can be called before or after ts->vec_sol has been updated.
265d2567f34SHong Zhang */
266d2567f34SHong Zhang static PetscErrorCode TSEvaluateStep_IRK(TS ts,PetscInt order,Vec U,PetscBool *done)
267d2567f34SHong Zhang {
268d2567f34SHong Zhang   TS_IRK         *irk = (TS_IRK*)ts->data;
269d2567f34SHong Zhang   IRKTableau     tab = irk->tableau;
270d2567f34SHong Zhang   Vec            *YdotI = irk->YdotI;
271d2567f34SHong Zhang   PetscScalar    *w = irk->work;
272d2567f34SHong Zhang   PetscReal      h;
273d2567f34SHong Zhang   PetscInt       j;
274d2567f34SHong Zhang 
275d2567f34SHong Zhang   PetscFunctionBegin;
276d2567f34SHong Zhang   switch (irk->status) {
277d2567f34SHong Zhang   case TS_STEP_INCOMPLETE:
278d2567f34SHong Zhang   case TS_STEP_PENDING:
279d2567f34SHong Zhang     h = ts->time_step; break;
280d2567f34SHong Zhang   case TS_STEP_COMPLETE:
281d2567f34SHong Zhang     h = ts->ptime - ts->ptime_prev; break;
282d2567f34SHong Zhang   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
283d2567f34SHong Zhang   }
284d2567f34SHong Zhang 
2859566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol,U));
286d2567f34SHong Zhang   for (j=0; j<irk->nstages; j++) w[j] = h*tab->b[j];
2879566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(U,irk->nstages,w,YdotI));
288d2567f34SHong Zhang   PetscFunctionReturn(0);
289d2567f34SHong Zhang }
290d2567f34SHong Zhang 
291d2567f34SHong Zhang static PetscErrorCode TSRollBack_IRK(TS ts)
292d2567f34SHong Zhang {
293d2567f34SHong Zhang   TS_IRK         *irk = (TS_IRK*)ts->data;
294d2567f34SHong Zhang 
295d2567f34SHong Zhang   PetscFunctionBegin;
2969566063dSJacob Faibussowitsch   PetscCall(VecCopy(irk->U0,ts->vec_sol));
297d2567f34SHong Zhang   PetscFunctionReturn(0);
298d2567f34SHong Zhang }
299d2567f34SHong Zhang 
300d2567f34SHong Zhang static PetscErrorCode TSStep_IRK(TS ts)
301d2567f34SHong Zhang {
302d2567f34SHong Zhang   TS_IRK          *irk = (TS_IRK*)ts->data;
303d2567f34SHong Zhang   IRKTableau      tab = irk->tableau;
304d2567f34SHong Zhang   PetscScalar     *A_inv = tab->A_inv,*A_inv_rowsum = tab->A_inv_rowsum;
305d2567f34SHong Zhang   const PetscInt  nstages = irk->nstages;
306d2567f34SHong Zhang   SNES            snes;
307d2567f34SHong Zhang   PetscInt        i,j,its,lits,bs;
308d2567f34SHong Zhang   TSAdapt         adapt;
309d2567f34SHong Zhang   PetscInt        rejections = 0;
310d2567f34SHong Zhang   PetscBool       accept = PETSC_TRUE;
311d2567f34SHong Zhang   PetscReal       next_time_step = ts->time_step;
312d2567f34SHong Zhang 
313d2567f34SHong Zhang   PetscFunctionBegin;
314d2567f34SHong Zhang   if (!ts->steprollback) {
3159566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vec_sol,irk->U0));
316d2567f34SHong Zhang   }
3179566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(ts->vec_sol,&bs));
318d2567f34SHong Zhang   for (i=0; i<nstages; i++) {
3199566063dSJacob Faibussowitsch     PetscCall(VecStrideScatter(ts->vec_sol,i*bs,irk->Z,INSERT_VALUES));
320d2567f34SHong Zhang   }
321d2567f34SHong Zhang 
322d2567f34SHong Zhang   irk->status = TS_STEP_INCOMPLETE;
323d2567f34SHong Zhang   while (!ts->reason && irk->status != TS_STEP_COMPLETE) {
3249566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vec_sol,irk->U));
3259566063dSJacob Faibussowitsch     PetscCall(TSGetSNES(ts,&snes));
3269566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes,NULL,irk->Z));
3279566063dSJacob Faibussowitsch     PetscCall(SNESGetIterationNumber(snes,&its));
3289566063dSJacob Faibussowitsch     PetscCall(SNESGetLinearSolveIterations(snes,&lits));
329d2567f34SHong Zhang     ts->snes_its += its; ts->ksp_its += lits;
3309566063dSJacob Faibussowitsch     PetscCall(VecStrideGatherAll(irk->Z,irk->Y,INSERT_VALUES));
331d2567f34SHong Zhang     for (i=0; i<nstages; i++) {
3329566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(irk->YdotI[i]));
333d2567f34SHong Zhang       for (j=0; j<nstages; j++) {
3349566063dSJacob Faibussowitsch         PetscCall(VecAXPY(irk->YdotI[i],A_inv[i+j*nstages]/ts->time_step,irk->Y[j]));
335d2567f34SHong Zhang       }
3369566063dSJacob Faibussowitsch       PetscCall(VecAXPY(irk->YdotI[i],-A_inv_rowsum[i]/ts->time_step,irk->U));
337d2567f34SHong Zhang     }
338d2567f34SHong Zhang     irk->status = TS_STEP_INCOMPLETE;
3399566063dSJacob Faibussowitsch     PetscCall(TSEvaluateStep_IRK(ts,irk->order,ts->vec_sol,NULL));
340d2567f34SHong Zhang     irk->status = TS_STEP_PENDING;
3419566063dSJacob Faibussowitsch     PetscCall(TSGetAdapt(ts,&adapt));
3429566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept));
343d2567f34SHong Zhang     irk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
344d2567f34SHong Zhang     if (!accept) {
3459566063dSJacob Faibussowitsch       PetscCall(TSRollBack_IRK(ts));
346d2567f34SHong Zhang       ts->time_step = next_time_step;
347d2567f34SHong Zhang       goto reject_step;
348d2567f34SHong Zhang     }
349d2567f34SHong Zhang 
350d2567f34SHong Zhang     ts->ptime += ts->time_step;
351d2567f34SHong Zhang     ts->time_step = next_time_step;
352d2567f34SHong Zhang     break;
353d2567f34SHong Zhang   reject_step:
354d2567f34SHong Zhang     ts->reject++; accept = PETSC_FALSE;
355d2567f34SHong Zhang     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
356d2567f34SHong Zhang       ts->reason = TS_DIVERGED_STEP_REJECTED;
35763a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts,"Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n",ts->steps,rejections));
358d2567f34SHong Zhang     }
359d2567f34SHong Zhang   }
360d2567f34SHong Zhang   PetscFunctionReturn(0);
361d2567f34SHong Zhang }
362d2567f34SHong Zhang 
363d2567f34SHong Zhang static PetscErrorCode TSInterpolate_IRK(TS ts,PetscReal itime,Vec U)
364d2567f34SHong Zhang {
365d2567f34SHong Zhang   TS_IRK          *irk = (TS_IRK*)ts->data;
366d2567f34SHong Zhang   PetscInt        nstages = irk->nstages,pinterp = irk->pinterp,i,j;
367d2567f34SHong Zhang   PetscReal       h;
368d2567f34SHong Zhang   PetscReal       tt,t;
369d2567f34SHong Zhang   PetscScalar     *bt;
370d2567f34SHong Zhang   const PetscReal *B = irk->tableau->binterp;
371d2567f34SHong Zhang 
372d2567f34SHong Zhang   PetscFunctionBegin;
3733c633725SBarry Smith   PetscCheck(B,PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSIRK %s does not have an interpolation formula",irk->method_name);
374d2567f34SHong Zhang   switch (irk->status) {
375d2567f34SHong Zhang   case TS_STEP_INCOMPLETE:
376d2567f34SHong Zhang   case TS_STEP_PENDING:
377d2567f34SHong Zhang     h = ts->time_step;
378d2567f34SHong Zhang     t = (itime - ts->ptime)/h;
379d2567f34SHong Zhang     break;
380d2567f34SHong Zhang   case TS_STEP_COMPLETE:
381d2567f34SHong Zhang     h = ts->ptime - ts->ptime_prev;
382d2567f34SHong Zhang     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
383d2567f34SHong Zhang     break;
384d2567f34SHong Zhang   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
385d2567f34SHong Zhang   }
3869566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nstages,&bt));
387d2567f34SHong Zhang   for (i=0; i<nstages; i++) bt[i] = 0;
388d2567f34SHong Zhang   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
389d2567f34SHong Zhang     for (i=0; i<nstages; i++) {
390d2567f34SHong Zhang       bt[i] += h * B[i*pinterp+j] * tt;
391d2567f34SHong Zhang     }
392d2567f34SHong Zhang   }
3939566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(U,nstages,bt,irk->YdotI));
394d2567f34SHong Zhang   PetscFunctionReturn(0);
395d2567f34SHong Zhang }
396d2567f34SHong Zhang 
397d2567f34SHong Zhang static PetscErrorCode TSIRKTableauReset(TS ts)
398d2567f34SHong Zhang {
399d2567f34SHong Zhang   TS_IRK         *irk = (TS_IRK*)ts->data;
400d2567f34SHong Zhang   IRKTableau     tab = irk->tableau;
401d2567f34SHong Zhang 
402d2567f34SHong Zhang   PetscFunctionBegin;
403d2567f34SHong Zhang   if (!tab) PetscFunctionReturn(0);
4049566063dSJacob Faibussowitsch   PetscCall(PetscFree3(tab->A,tab->A_inv,tab->I_s));
4059566063dSJacob Faibussowitsch   PetscCall(PetscFree4(tab->b,tab->c,tab->binterp,tab->A_inv_rowsum));
406d2567f34SHong Zhang   PetscFunctionReturn(0);
407d2567f34SHong Zhang }
408d2567f34SHong Zhang 
409d2567f34SHong Zhang static PetscErrorCode TSReset_IRK(TS ts)
410d2567f34SHong Zhang {
411d2567f34SHong Zhang   TS_IRK         *irk = (TS_IRK*)ts->data;
412d2567f34SHong Zhang 
413d2567f34SHong Zhang   PetscFunctionBegin;
4149566063dSJacob Faibussowitsch   PetscCall(TSIRKTableauReset(ts));
415*1baa6e33SBarry Smith   if (irk->tableau) PetscCall(PetscFree(irk->tableau));
416*1baa6e33SBarry Smith   if (irk->method_name) PetscCall(PetscFree(irk->method_name));
417*1baa6e33SBarry Smith   if (irk->work) PetscCall(PetscFree(irk->work));
4189566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(irk->nstages,&irk->Y));
4199566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(irk->nstages,&irk->YdotI));
4209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&irk->Ydot));
4219566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&irk->Z));
4229566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&irk->U));
4239566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&irk->U0));
4249566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&irk->TJ));
425d2567f34SHong Zhang   PetscFunctionReturn(0);
426d2567f34SHong Zhang }
427d2567f34SHong Zhang 
428d2567f34SHong Zhang static PetscErrorCode TSIRKGetVecs(TS ts,DM dm,Vec *U)
429d2567f34SHong Zhang {
430d2567f34SHong Zhang   TS_IRK         *irk = (TS_IRK*)ts->data;
431d2567f34SHong Zhang 
432d2567f34SHong Zhang   PetscFunctionBegin;
433d2567f34SHong Zhang   if (U) {
434d2567f34SHong Zhang     if (dm && dm != ts->dm) {
4359566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm,"TSIRK_U",U));
436d2567f34SHong Zhang     } else *U = irk->U;
437d2567f34SHong Zhang   }
438d2567f34SHong Zhang   PetscFunctionReturn(0);
439d2567f34SHong Zhang }
440d2567f34SHong Zhang 
441d2567f34SHong Zhang static PetscErrorCode TSIRKRestoreVecs(TS ts,DM dm,Vec *U)
442d2567f34SHong Zhang {
443d2567f34SHong Zhang   PetscFunctionBegin;
444d2567f34SHong Zhang   if (U) {
445d2567f34SHong Zhang     if (dm && dm != ts->dm) {
4469566063dSJacob Faibussowitsch       PetscCall(DMRestoreNamedGlobalVector(dm,"TSIRK_U",U));
447d2567f34SHong Zhang     }
448d2567f34SHong Zhang   }
449d2567f34SHong Zhang   PetscFunctionReturn(0);
450d2567f34SHong Zhang }
451d2567f34SHong Zhang 
452d2567f34SHong Zhang /*
453d2567f34SHong Zhang   This defines the nonlinear equations that is to be solved with SNES
454d2567f34SHong Zhang     G[e\otimes t + C*dt, Z, Zdot] = 0
455d2567f34SHong Zhang     Zdot = (In \otimes S)*Z - (In \otimes Se) U
456d2567f34SHong Zhang   where S = 1/(dt*A)
457d2567f34SHong Zhang */
458d2567f34SHong Zhang static PetscErrorCode SNESTSFormFunction_IRK(SNES snes,Vec ZC,Vec FC,TS ts)
459d2567f34SHong Zhang {
460d2567f34SHong Zhang   TS_IRK            *irk = (TS_IRK*)ts->data;
461d2567f34SHong Zhang   IRKTableau        tab  = irk->tableau;
462d2567f34SHong Zhang   const PetscInt    nstages = irk->nstages;
463d2567f34SHong Zhang   const PetscReal   *c = tab->c;
464d2567f34SHong Zhang   const PetscScalar *A_inv = tab->A_inv,*A_inv_rowsum = tab->A_inv_rowsum;
465d2567f34SHong Zhang   DM                dm,dmsave;
466d2567f34SHong Zhang   Vec               U,*YdotI = irk->YdotI,Ydot = irk->Ydot,*Y = irk->Y;
467d2567f34SHong Zhang   PetscReal         h = ts->time_step;
468d2567f34SHong Zhang   PetscInt          i,j;
469d2567f34SHong Zhang 
470d2567f34SHong Zhang   PetscFunctionBegin;
4719566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
4729566063dSJacob Faibussowitsch   PetscCall(TSIRKGetVecs(ts,dm,&U));
4739566063dSJacob Faibussowitsch   PetscCall(VecStrideGatherAll(ZC,Y,INSERT_VALUES));
474d2567f34SHong Zhang   dmsave = ts->dm;
475d2567f34SHong Zhang   ts->dm = dm;
476d2567f34SHong Zhang   for (i=0; i<nstages; i++) {
4779566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Ydot));
478d2567f34SHong Zhang     for (j=0; j<nstages; j++) {
4799566063dSJacob Faibussowitsch       PetscCall(VecAXPY(Ydot,A_inv[j*nstages+i]/h,Y[j]));
480d2567f34SHong Zhang     }
4819566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Ydot,-A_inv_rowsum[i]/h,U)); /* Ydot = (S \otimes In)*Z - (Se \otimes In) U */
4829566063dSJacob Faibussowitsch     PetscCall(TSComputeIFunction(ts,ts->ptime+ts->time_step*c[i],Y[i],Ydot,YdotI[i],PETSC_FALSE));
483d2567f34SHong Zhang   }
4849566063dSJacob Faibussowitsch   PetscCall(VecStrideScatterAll(YdotI,FC,INSERT_VALUES));
485d2567f34SHong Zhang   ts->dm = dmsave;
4869566063dSJacob Faibussowitsch   PetscCall(TSIRKRestoreVecs(ts,dm,&U));
487d2567f34SHong Zhang   PetscFunctionReturn(0);
488d2567f34SHong Zhang }
489d2567f34SHong Zhang 
490d2567f34SHong Zhang /*
491d2567f34SHong Zhang    For explicit ODE, the Jacobian is
492d2567f34SHong Zhang      JC = I_n \otimes S - J \otimes I_s
493d2567f34SHong Zhang    For DAE, the Jacobian is
494d2567f34SHong Zhang      JC = M_n \otimes S - J \otimes I_s
495d2567f34SHong Zhang */
496d2567f34SHong Zhang static PetscErrorCode SNESTSFormJacobian_IRK(SNES snes,Vec ZC,Mat JC,Mat JCpre,TS ts)
497d2567f34SHong Zhang {
498d2567f34SHong Zhang   TS_IRK          *irk = (TS_IRK*)ts->data;
499d2567f34SHong Zhang   IRKTableau      tab  = irk->tableau;
500d2567f34SHong Zhang   const PetscInt  nstages = irk->nstages;
501d2567f34SHong Zhang   const PetscReal *c = tab->c;
502d2567f34SHong Zhang   DM              dm,dmsave;
503d2567f34SHong Zhang   Vec             *Y = irk->Y,Ydot = irk->Ydot;
504d2567f34SHong Zhang   Mat             J;
505d2567f34SHong Zhang   PetscScalar     *S;
506d2567f34SHong Zhang   PetscInt        i,j,bs;
507d2567f34SHong Zhang 
508d2567f34SHong Zhang   PetscFunctionBegin;
5099566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
510d2567f34SHong Zhang   /* irk->Ydot has already been computed in SNESTSFormFunction_IRK (SNES guarantees this) */
511d2567f34SHong Zhang   dmsave = ts->dm;
512d2567f34SHong Zhang   ts->dm = dm;
5139566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(Y[nstages-1],&bs));
514d2567f34SHong Zhang   if (ts->equation_type <= TS_EQ_ODE_EXPLICIT) { /* Support explicit formulas only */
5159566063dSJacob Faibussowitsch     PetscCall(VecStrideGather(ZC,(nstages-1)*bs,Y[nstages-1],INSERT_VALUES));
5169566063dSJacob Faibussowitsch     PetscCall(MatKAIJGetAIJ(JC,&J));
5179566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts,ts->ptime+ts->time_step*c[nstages-1],Y[nstages-1],Ydot,0,J,J,PETSC_FALSE));
5189566063dSJacob Faibussowitsch     PetscCall(MatKAIJGetS(JC,NULL,NULL,&S));
519d2567f34SHong Zhang     for (i=0; i<nstages; i++)
520d2567f34SHong Zhang       for (j=0; j<nstages; j++)
521d2567f34SHong Zhang         S[i+nstages*j] = tab->A_inv[i+nstages*j]/ts->time_step;
5229566063dSJacob Faibussowitsch     PetscCall(MatKAIJRestoreS(JC,&S));
52398921bdaSJacob 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  */
524d2567f34SHong Zhang   ts->dm = dmsave;
525d2567f34SHong Zhang   PetscFunctionReturn(0);
526d2567f34SHong Zhang }
527d2567f34SHong Zhang 
528d2567f34SHong Zhang static PetscErrorCode DMCoarsenHook_TSIRK(DM fine,DM coarse,void *ctx)
529d2567f34SHong Zhang {
530d2567f34SHong Zhang   PetscFunctionBegin;
531d2567f34SHong Zhang   PetscFunctionReturn(0);
532d2567f34SHong Zhang }
533d2567f34SHong Zhang 
534d2567f34SHong Zhang static PetscErrorCode DMRestrictHook_TSIRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
535d2567f34SHong Zhang {
536d2567f34SHong Zhang   TS             ts = (TS)ctx;
537d2567f34SHong Zhang   Vec            U,U_c;
538d2567f34SHong Zhang 
539d2567f34SHong Zhang   PetscFunctionBegin;
5409566063dSJacob Faibussowitsch   PetscCall(TSIRKGetVecs(ts,fine,&U));
5419566063dSJacob Faibussowitsch   PetscCall(TSIRKGetVecs(ts,coarse,&U_c));
5429566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct,U,U_c));
5439566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(U_c,rscale,U_c));
5449566063dSJacob Faibussowitsch   PetscCall(TSIRKRestoreVecs(ts,fine,&U));
5459566063dSJacob Faibussowitsch   PetscCall(TSIRKRestoreVecs(ts,coarse,&U_c));
546d2567f34SHong Zhang   PetscFunctionReturn(0);
547d2567f34SHong Zhang }
548d2567f34SHong Zhang 
549d2567f34SHong Zhang static PetscErrorCode DMSubDomainHook_TSIRK(DM dm,DM subdm,void *ctx)
550d2567f34SHong Zhang {
551d2567f34SHong Zhang   PetscFunctionBegin;
552d2567f34SHong Zhang   PetscFunctionReturn(0);
553d2567f34SHong Zhang }
554d2567f34SHong Zhang 
555d2567f34SHong Zhang static PetscErrorCode DMSubDomainRestrictHook_TSIRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
556d2567f34SHong Zhang {
557d2567f34SHong Zhang   TS             ts = (TS)ctx;
558d2567f34SHong Zhang   Vec            U,U_c;
559d2567f34SHong Zhang 
560d2567f34SHong Zhang   PetscFunctionBegin;
5619566063dSJacob Faibussowitsch   PetscCall(TSIRKGetVecs(ts,dm,&U));
5629566063dSJacob Faibussowitsch   PetscCall(TSIRKGetVecs(ts,subdm,&U_c));
563d2567f34SHong Zhang 
5649566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(gscat,U,U_c,INSERT_VALUES,SCATTER_FORWARD));
5659566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(gscat,U,U_c,INSERT_VALUES,SCATTER_FORWARD));
566d2567f34SHong Zhang 
5679566063dSJacob Faibussowitsch   PetscCall(TSIRKRestoreVecs(ts,dm,&U));
5689566063dSJacob Faibussowitsch   PetscCall(TSIRKRestoreVecs(ts,subdm,&U_c));
569d2567f34SHong Zhang   PetscFunctionReturn(0);
570d2567f34SHong Zhang }
571d2567f34SHong Zhang 
572d2567f34SHong Zhang static PetscErrorCode TSSetUp_IRK(TS ts)
573d2567f34SHong Zhang {
574d2567f34SHong Zhang   TS_IRK         *irk = (TS_IRK*)ts->data;
575d2567f34SHong Zhang   IRKTableau     tab = irk->tableau;
576d2567f34SHong Zhang   DM             dm;
577d2567f34SHong Zhang   Mat            J;
578d2567f34SHong Zhang   Vec            R;
579d2567f34SHong Zhang   const PetscInt nstages = irk->nstages;
580d2567f34SHong Zhang   PetscInt       vsize,bs;
581d2567f34SHong Zhang 
582d2567f34SHong Zhang   PetscFunctionBegin;
583d2567f34SHong Zhang   if (!irk->work) {
5849566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(irk->nstages,&irk->work));
585d2567f34SHong Zhang   }
586d2567f34SHong Zhang   if (!irk->Y) {
5879566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vec_sol,irk->nstages,&irk->Y));
588d2567f34SHong Zhang   }
589d2567f34SHong Zhang   if (!irk->YdotI) {
5909566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vec_sol,irk->nstages,&irk->YdotI));
591d2567f34SHong Zhang   }
592d2567f34SHong Zhang   if (!irk->Ydot) {
5939566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol,&irk->Ydot));
594d2567f34SHong Zhang   }
595d2567f34SHong Zhang   if (!irk->U) {
5969566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol,&irk->U));
597d2567f34SHong Zhang   }
598d2567f34SHong Zhang   if (!irk->U0) {
5999566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol,&irk->U0));
600d2567f34SHong Zhang   }
601d2567f34SHong Zhang   if (!irk->Z) {
6029566063dSJacob Faibussowitsch     PetscCall(VecCreate(PetscObjectComm((PetscObject)ts->vec_sol),&irk->Z));
6039566063dSJacob Faibussowitsch     PetscCall(VecGetSize(ts->vec_sol,&vsize));
6049566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(irk->Z,PETSC_DECIDE,vsize*irk->nstages));
6059566063dSJacob Faibussowitsch     PetscCall(VecGetBlockSize(ts->vec_sol,&bs));
6069566063dSJacob Faibussowitsch     PetscCall(VecSetBlockSize(irk->Z,irk->nstages*bs));
6079566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(irk->Z));
608d2567f34SHong Zhang   }
6099566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&dm));
6109566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(dm,DMCoarsenHook_TSIRK,DMRestrictHook_TSIRK,ts));
6119566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(dm,DMSubDomainHook_TSIRK,DMSubDomainRestrictHook_TSIRK,ts));
612d2567f34SHong Zhang 
6139566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts,&ts->snes));
6149566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(irk->Z,&R));
6159566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(ts->snes,R,SNESTSFormFunction,ts));
6169566063dSJacob Faibussowitsch   PetscCall(TSGetIJacobian(ts,&J,NULL,NULL,NULL));
617d2567f34SHong Zhang   if (!irk->TJ) {
618d2567f34SHong Zhang     /* Create the KAIJ matrix for solving the stages */
6199566063dSJacob Faibussowitsch     PetscCall(MatCreateKAIJ(J,nstages,nstages,tab->A_inv,tab->I_s,&irk->TJ));
620d2567f34SHong Zhang   }
6219566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(ts->snes,irk->TJ,irk->TJ,SNESTSFormJacobian,ts));
6229566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&R));
623d2567f34SHong Zhang   PetscFunctionReturn(0);
624d2567f34SHong Zhang }
625d2567f34SHong Zhang 
626d2567f34SHong Zhang static PetscErrorCode TSSetFromOptions_IRK(PetscOptionItems *PetscOptionsObject,TS ts)
627d2567f34SHong Zhang {
628d2567f34SHong Zhang   TS_IRK         *irk = (TS_IRK*)ts->data;
629d2567f34SHong Zhang   char           tname[256] = TSIRKGAUSS;
630d2567f34SHong Zhang 
631d2567f34SHong Zhang   PetscFunctionBegin;
632d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"IRK ODE solver options");
633d2567f34SHong Zhang   {
634d2567f34SHong Zhang     PetscBool flg1,flg2;
6359566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-ts_irk_nstages","Stages of the IRK method","TSIRKSetNumStages",irk->nstages,&irk->nstages,&flg1));
6369566063dSJacob Faibussowitsch     PetscCall(PetscOptionsFList("-ts_irk_type","Type of IRK method","TSIRKSetType",TSIRKList,irk->method_name[0] ? irk->method_name : tname,tname,sizeof(tname),&flg2));
637d2567f34SHong Zhang     if (flg1 ||flg2 || !irk->method_name[0]) { /* Create the method tableau after nstages or method is set */
6389566063dSJacob Faibussowitsch       PetscCall(TSIRKSetType(ts,tname));
639d2567f34SHong Zhang     }
640d2567f34SHong Zhang   }
641d0609cedSBarry Smith   PetscOptionsHeadEnd();
642d2567f34SHong Zhang   PetscFunctionReturn(0);
643d2567f34SHong Zhang }
644d2567f34SHong Zhang 
645d2567f34SHong Zhang static PetscErrorCode TSView_IRK(TS ts,PetscViewer viewer)
646d2567f34SHong Zhang {
647d2567f34SHong Zhang   TS_IRK         *irk = (TS_IRK*)ts->data;
648d2567f34SHong Zhang   PetscBool      iascii;
649d2567f34SHong Zhang 
650d2567f34SHong Zhang   PetscFunctionBegin;
6519566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
652d2567f34SHong Zhang   if (iascii) {
653d2567f34SHong Zhang     IRKTableau    tab = irk->tableau;
654d2567f34SHong Zhang     TSIRKType irktype;
655d2567f34SHong Zhang     char          buf[512];
656d2567f34SHong Zhang 
6579566063dSJacob Faibussowitsch     PetscCall(TSIRKGetType(ts,&irktype));
6589566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  IRK type %s\n",irktype));
6599566063dSJacob Faibussowitsch     PetscCall(PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",irk->nstages,tab->c));
6609566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  Abscissa       c = %s\n",buf));
6619566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s\n",irk->stiffly_accurate ? "yes" : "no"));
6629566063dSJacob Faibussowitsch     PetscCall(PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",PetscSqr(irk->nstages),tab->A));
6639566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  A coefficients       A = %s\n",buf));
664d2567f34SHong Zhang   }
665d2567f34SHong Zhang   PetscFunctionReturn(0);
666d2567f34SHong Zhang }
667d2567f34SHong Zhang 
668d2567f34SHong Zhang static PetscErrorCode TSLoad_IRK(TS ts,PetscViewer viewer)
669d2567f34SHong Zhang {
670d2567f34SHong Zhang   SNES           snes;
671d2567f34SHong Zhang   TSAdapt        adapt;
672d2567f34SHong Zhang 
673d2567f34SHong Zhang   PetscFunctionBegin;
6749566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts,&adapt));
6759566063dSJacob Faibussowitsch   PetscCall(TSAdaptLoad(adapt,viewer));
6769566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts,&snes));
6779566063dSJacob Faibussowitsch   PetscCall(SNESLoad(snes,viewer));
678d2567f34SHong Zhang   /* function and Jacobian context for SNES when used with TS is always ts object */
6799566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes,NULL,NULL,ts));
6809566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes,NULL,NULL,NULL,ts));
681d2567f34SHong Zhang   PetscFunctionReturn(0);
682d2567f34SHong Zhang }
683d2567f34SHong Zhang 
684d2567f34SHong Zhang /*@C
685d2567f34SHong Zhang   TSIRKSetType - Set the type of IRK scheme
686d2567f34SHong Zhang 
687d2567f34SHong Zhang   Logically collective
688d2567f34SHong Zhang 
68902024617SSatish Balay   Input Parameters:
690d2567f34SHong Zhang +  ts - timestepping context
691d2567f34SHong Zhang -  irktype - type of IRK scheme
692d2567f34SHong Zhang 
693d2567f34SHong Zhang   Options Database:
69467b8a455SSatish Balay .  -ts_irk_type <gauss> - set irk type
695d2567f34SHong Zhang 
696d2567f34SHong Zhang   Level: intermediate
697d2567f34SHong Zhang 
698db781477SPatrick Sanan .seealso: `TSIRKGetType()`, `TSIRK`, `TSIRKType`, `TSIRKGAUSS`
699d2567f34SHong Zhang @*/
700d2567f34SHong Zhang PetscErrorCode TSIRKSetType(TS ts,TSIRKType irktype)
701d2567f34SHong Zhang {
702d2567f34SHong Zhang   PetscFunctionBegin;
703d2567f34SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
704d2567f34SHong Zhang   PetscValidCharPointer(irktype,2);
705cac4c232SBarry Smith   PetscTryMethod(ts,"TSIRKSetType_C",(TS,TSIRKType),(ts,irktype));
706d2567f34SHong Zhang   PetscFunctionReturn(0);
707d2567f34SHong Zhang }
708d2567f34SHong Zhang 
709d2567f34SHong Zhang /*@C
710d2567f34SHong Zhang   TSIRKGetType - Get the type of IRK IMEX scheme
711d2567f34SHong Zhang 
712d2567f34SHong Zhang   Logically collective
713d2567f34SHong Zhang 
714d2567f34SHong Zhang   Input Parameter:
715d2567f34SHong Zhang .  ts - timestepping context
716d2567f34SHong Zhang 
717d2567f34SHong Zhang   Output Parameter:
718d2567f34SHong Zhang .  irktype - type of IRK-IMEX scheme
719d2567f34SHong Zhang 
720d2567f34SHong Zhang   Level: intermediate
721d2567f34SHong Zhang 
722db781477SPatrick Sanan .seealso: `TSIRKGetType()`
723d2567f34SHong Zhang @*/
724d2567f34SHong Zhang PetscErrorCode TSIRKGetType(TS ts,TSIRKType *irktype)
725d2567f34SHong Zhang {
726d2567f34SHong Zhang   PetscFunctionBegin;
727d2567f34SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
728cac4c232SBarry Smith   PetscUseMethod(ts,"TSIRKGetType_C",(TS,TSIRKType*),(ts,irktype));
729d2567f34SHong Zhang   PetscFunctionReturn(0);
730d2567f34SHong Zhang }
731d2567f34SHong Zhang 
732d2567f34SHong Zhang /*@C
733d2567f34SHong Zhang   TSIRKSetNumStages - Set the number of stages of IRK scheme
734d2567f34SHong Zhang 
735d2567f34SHong Zhang   Logically collective
736d2567f34SHong Zhang 
73702024617SSatish Balay   Input Parameters:
738d2567f34SHong Zhang +  ts - timestepping context
739d2567f34SHong Zhang -  nstages - number of stages of IRK scheme
740d2567f34SHong Zhang 
741d2567f34SHong Zhang   Options Database:
74267b8a455SSatish Balay .  -ts_irk_nstages <int> - set number of stages
743d2567f34SHong Zhang 
744d2567f34SHong Zhang   Level: intermediate
745d2567f34SHong Zhang 
746db781477SPatrick Sanan .seealso: `TSIRKGetNumStages()`, `TSIRK`
747d2567f34SHong Zhang @*/
748d2567f34SHong Zhang PetscErrorCode TSIRKSetNumStages(TS ts,PetscInt nstages)
749d2567f34SHong Zhang {
750d2567f34SHong Zhang   PetscFunctionBegin;
751d2567f34SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
752cac4c232SBarry Smith   PetscTryMethod(ts,"TSIRKSetNumStages_C",(TS,PetscInt),(ts,nstages));
753d2567f34SHong Zhang   PetscFunctionReturn(0);
754d2567f34SHong Zhang }
755d2567f34SHong Zhang 
756d2567f34SHong Zhang /*@C
757d2567f34SHong Zhang   TSIRKGetNumStages - Get the number of stages of IRK scheme
758d2567f34SHong Zhang 
759d2567f34SHong Zhang   Logically collective
760d2567f34SHong Zhang 
76102024617SSatish Balay   Input Parameters:
762d2567f34SHong Zhang +  ts - timestepping context
763d2567f34SHong Zhang -  nstages - number of stages of IRK scheme
764d2567f34SHong Zhang 
765d2567f34SHong Zhang   Level: intermediate
766d2567f34SHong Zhang 
767db781477SPatrick Sanan .seealso: `TSIRKSetNumStages()`, `TSIRK`
768d2567f34SHong Zhang @*/
769d2567f34SHong Zhang PetscErrorCode TSIRKGetNumStages(TS ts,PetscInt *nstages)
770d2567f34SHong Zhang {
771d2567f34SHong Zhang   PetscFunctionBegin;
772d2567f34SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
773d2567f34SHong Zhang   PetscValidIntPointer(nstages,2);
774cac4c232SBarry Smith   PetscTryMethod(ts,"TSIRKGetNumStages_C",(TS,PetscInt*),(ts,nstages));
775d2567f34SHong Zhang   PetscFunctionReturn(0);
776d2567f34SHong Zhang }
777d2567f34SHong Zhang 
778d2567f34SHong Zhang static PetscErrorCode TSIRKGetType_IRK(TS ts,TSIRKType *irktype)
779d2567f34SHong Zhang {
780d2567f34SHong Zhang   TS_IRK *irk = (TS_IRK*)ts->data;
781d2567f34SHong Zhang 
782d2567f34SHong Zhang   PetscFunctionBegin;
783d2567f34SHong Zhang   *irktype = irk->method_name;
784d2567f34SHong Zhang   PetscFunctionReturn(0);
785d2567f34SHong Zhang }
786d2567f34SHong Zhang 
787d2567f34SHong Zhang static PetscErrorCode TSIRKSetType_IRK(TS ts,TSIRKType irktype)
788d2567f34SHong Zhang {
789d2567f34SHong Zhang   TS_IRK         *irk = (TS_IRK*)ts->data;
7905f80ce2aSJacob Faibussowitsch   PetscErrorCode (*irkcreate)(TS);
791d2567f34SHong Zhang 
792d2567f34SHong Zhang   PetscFunctionBegin;
793d2567f34SHong Zhang   if (irk->method_name) {
7949566063dSJacob Faibussowitsch     PetscCall(PetscFree(irk->method_name));
7959566063dSJacob Faibussowitsch     PetscCall(TSIRKTableauReset(ts));
796d2567f34SHong Zhang   }
7979566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(TSIRKList,irktype,&irkcreate));
7983c633725SBarry Smith   PetscCheck(irkcreate,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSIRK type \"%s\" given",irktype);
7999566063dSJacob Faibussowitsch   PetscCall((*irkcreate)(ts));
8009566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(irktype,&irk->method_name));
801d2567f34SHong Zhang   PetscFunctionReturn(0);
802d2567f34SHong Zhang }
803d2567f34SHong Zhang 
804d2567f34SHong Zhang static PetscErrorCode TSIRKSetNumStages_IRK(TS ts,PetscInt nstages)
805d2567f34SHong Zhang {
806d2567f34SHong Zhang   TS_IRK *irk = (TS_IRK*)ts->data;
807d2567f34SHong Zhang 
808d2567f34SHong Zhang   PetscFunctionBegin;
80963a3b9bcSJacob Faibussowitsch   PetscCheck(nstages>0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"input argument, %" PetscInt_FMT ", out of range",nstages);
810d2567f34SHong Zhang   irk->nstages = nstages;
811d2567f34SHong Zhang   PetscFunctionReturn(0);
812d2567f34SHong Zhang }
813d2567f34SHong Zhang 
814d2567f34SHong Zhang static PetscErrorCode TSIRKGetNumStages_IRK(TS ts,PetscInt *nstages)
815d2567f34SHong Zhang {
816d2567f34SHong Zhang   TS_IRK *irk = (TS_IRK*)ts->data;
817d2567f34SHong Zhang 
818d2567f34SHong Zhang   PetscFunctionBegin;
819d2567f34SHong Zhang   PetscValidIntPointer(nstages,2);
820d2567f34SHong Zhang   *nstages = irk->nstages;
821d2567f34SHong Zhang   PetscFunctionReturn(0);
822d2567f34SHong Zhang }
823d2567f34SHong Zhang 
824d2567f34SHong Zhang static PetscErrorCode TSDestroy_IRK(TS ts)
825d2567f34SHong Zhang {
826d2567f34SHong Zhang   PetscFunctionBegin;
8279566063dSJacob Faibussowitsch   PetscCall(TSReset_IRK(ts));
828d2567f34SHong Zhang   if (ts->dm) {
8299566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSIRK,DMRestrictHook_TSIRK,ts));
8309566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSIRK,DMSubDomainRestrictHook_TSIRK,ts));
831d2567f34SHong Zhang   }
8329566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
8339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSIRKSetType_C",NULL));
8349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSIRKGetType_C",NULL));
8359566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSIRKSetNumStages_C",NULL));
8369566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSIRKGetNumStages_C",NULL));
837d2567f34SHong Zhang   PetscFunctionReturn(0);
838d2567f34SHong Zhang }
839d2567f34SHong Zhang 
840d2567f34SHong Zhang /*MC
841d2567f34SHong Zhang       TSIRK - ODE and DAE solver using Implicit Runge-Kutta schemes
842d2567f34SHong Zhang 
843d2567f34SHong Zhang   Notes:
844d2567f34SHong Zhang 
845d2567f34SHong Zhang   TSIRK uses the sparse Kronecker product matrix implementation of MATKAIJ to achieve good arithmetic intensity.
846d2567f34SHong Zhang 
847517f343bSHong 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().
848d2567f34SHong Zhang 
849d2567f34SHong Zhang   Level: beginner
850d2567f34SHong Zhang 
851db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSIRKSetType()`, `TSIRKGetType()`, `TSIRKGAUSS`, `TSIRKRegister()`, `TSIRKSetNumStages()`
852d2567f34SHong Zhang 
853d2567f34SHong Zhang M*/
854d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSCreate_IRK(TS ts)
855d2567f34SHong Zhang {
856d2567f34SHong Zhang   TS_IRK         *irk;
857d2567f34SHong Zhang 
858d2567f34SHong Zhang   PetscFunctionBegin;
8599566063dSJacob Faibussowitsch   PetscCall(TSIRKInitializePackage());
860d2567f34SHong Zhang 
861d2567f34SHong Zhang   ts->ops->reset          = TSReset_IRK;
862d2567f34SHong Zhang   ts->ops->destroy        = TSDestroy_IRK;
863d2567f34SHong Zhang   ts->ops->view           = TSView_IRK;
864d2567f34SHong Zhang   ts->ops->load           = TSLoad_IRK;
865d2567f34SHong Zhang   ts->ops->setup          = TSSetUp_IRK;
866d2567f34SHong Zhang   ts->ops->step           = TSStep_IRK;
867d2567f34SHong Zhang   ts->ops->interpolate    = TSInterpolate_IRK;
868d2567f34SHong Zhang   ts->ops->evaluatestep   = TSEvaluateStep_IRK;
869d2567f34SHong Zhang   ts->ops->rollback       = TSRollBack_IRK;
870d2567f34SHong Zhang   ts->ops->setfromoptions = TSSetFromOptions_IRK;
871d2567f34SHong Zhang   ts->ops->snesfunction   = SNESTSFormFunction_IRK;
872d2567f34SHong Zhang   ts->ops->snesjacobian   = SNESTSFormJacobian_IRK;
873d2567f34SHong Zhang 
874d2567f34SHong Zhang   ts->usessnes = PETSC_TRUE;
875d2567f34SHong Zhang 
8769566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(ts,&irk));
877d2567f34SHong Zhang   ts->data = (void*)irk;
878d2567f34SHong Zhang 
8799566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSIRKSetType_C",TSIRKSetType_IRK));
8809566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSIRKGetType_C",TSIRKGetType_IRK));
8819566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSIRKSetNumStages_C",TSIRKSetNumStages_IRK));
8829566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSIRKGetNumStages_C",TSIRKGetNumStages_IRK));
883d2567f34SHong Zhang   /* 3-stage IRK_Gauss is the default */
8849566063dSJacob Faibussowitsch   PetscCall(PetscNew(&irk->tableau));
885d2567f34SHong Zhang   irk->nstages = 3;
8869566063dSJacob Faibussowitsch   PetscCall(TSIRKSetType(ts,TSIRKDefault));
887d2567f34SHong Zhang   PetscFunctionReturn(0);
888d2567f34SHong Zhang }
889