xref: /petsc/src/ts/impls/implicit/irk/irk.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 @*/
679371c9d4SSatish Balay 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   TS_IRK    *irk = (TS_IRK *)ts->data;
69d2567f34SHong Zhang   IRKTableau tab = irk->tableau;
70d2567f34SHong Zhang 
71d2567f34SHong Zhang   PetscFunctionBegin;
72d2567f34SHong Zhang   irk->order = nstages;
739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(PetscSqr(nstages), &tab->A, PetscSqr(nstages), &tab->A_inv, PetscSqr(nstages), &tab->I_s));
749566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(nstages, &tab->b, nstages, &tab->c, nstages, &tab->binterp, nstages, &tab->A_inv_rowsum));
759566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(tab->A, A, PetscSqr(nstages)));
769566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(tab->b, b, nstages));
779566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(tab->c, c, nstages));
78d2567f34SHong Zhang   /* optional coefficient arrays */
791baa6e33SBarry Smith   if (binterp) PetscCall(PetscArraycpy(tab->binterp, binterp, nstages));
801baa6e33SBarry Smith   if (A_inv) PetscCall(PetscArraycpy(tab->A_inv, A_inv, PetscSqr(nstages)));
811baa6e33SBarry Smith   if (A_inv_rowsum) PetscCall(PetscArraycpy(tab->A_inv_rowsum, A_inv_rowsum, nstages));
821baa6e33SBarry Smith   if (I_s) PetscCall(PetscArraycpy(tab->I_s, I_s, PetscSqr(nstages)));
83d2567f34SHong Zhang   PetscFunctionReturn(0);
84d2567f34SHong Zhang }
85d2567f34SHong Zhang 
86d2567f34SHong Zhang /* Arrays should be freed with PetscFree3(A,b,c) */
879371c9d4SSatish Balay static PetscErrorCode TSIRKCreate_Gauss(TS ts) {
88d2567f34SHong Zhang   PetscInt     nstages;
89d2567f34SHong Zhang   PetscReal   *gauss_A_real, *gauss_b, *b, *gauss_c;
90d2567f34SHong Zhang   PetscScalar *gauss_A, *gauss_A_inv, *gauss_A_inv_rowsum, *I_s;
91d2567f34SHong Zhang   PetscScalar *G0, *G1;
92d2567f34SHong Zhang   PetscInt     i, j;
93d2567f34SHong Zhang   Mat          G0mat, G1mat, Amat;
94d2567f34SHong Zhang 
95d2567f34SHong Zhang   PetscFunctionBegin;
969566063dSJacob Faibussowitsch   PetscCall(TSIRKGetNumStages(ts, &nstages));
979566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(PetscSqr(nstages), &gauss_A_real, nstages, &gauss_b, nstages, &gauss_c));
989566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(PetscSqr(nstages), &gauss_A, PetscSqr(nstages), &gauss_A_inv, nstages, &gauss_A_inv_rowsum, PetscSqr(nstages), &I_s));
999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(nstages, &b, PetscSqr(nstages), &G0, PetscSqr(nstages), &G1));
1009566063dSJacob Faibussowitsch   PetscCall(PetscDTGaussQuadrature(nstages, 0., 1., gauss_c, b));
101d2567f34SHong Zhang   for (i = 0; i < nstages; i++) gauss_b[i] = b[i]; /* copy to possibly-complex array */
102d2567f34SHong Zhang 
103d2567f34SHong Zhang   /* A^T = G0^{-1} G1 */
104d2567f34SHong Zhang   for (i = 0; i < nstages; i++) {
105d2567f34SHong Zhang     for (j = 0; j < nstages; j++) {
106d2567f34SHong Zhang       G0[i * nstages + j] = PetscPowRealInt(gauss_c[i], j);
107d2567f34SHong Zhang       G1[i * nstages + j] = PetscPowRealInt(gauss_c[i], j + 1) / (j + 1);
108d2567f34SHong Zhang     }
109d2567f34SHong Zhang   }
110d2567f34SHong Zhang   /* The arrays above are row-aligned, but we create dense matrices as the transpose */
1119566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, nstages, nstages, G0, &G0mat));
1129566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, nstages, nstages, G1, &G1mat));
1139566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, nstages, nstages, gauss_A, &Amat));
1149566063dSJacob Faibussowitsch   PetscCall(MatLUFactor(G0mat, NULL, NULL, NULL));
1159566063dSJacob Faibussowitsch   PetscCall(MatMatSolve(G0mat, G1mat, Amat));
1169566063dSJacob Faibussowitsch   PetscCall(MatTranspose(Amat, MAT_INPLACE_MATRIX, &Amat));
117d2567f34SHong Zhang   for (i = 0; i < nstages; i++)
1189371c9d4SSatish Balay     for (j = 0; j < nstages; j++) gauss_A_real[i * nstages + j] = PetscRealPart(gauss_A[i * nstages + j]);
119d2567f34SHong Zhang 
1209566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&G0mat));
1219566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&G1mat));
1229566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Amat));
1239566063dSJacob Faibussowitsch   PetscCall(PetscFree3(b, G0, G1));
124d2567f34SHong Zhang 
125d2567f34SHong Zhang   { /* Invert A */
126d2567f34SHong Zhang     /* PETSc does not provide a routine to calculate the inverse of a general matrix.
127d2567f34SHong Zhang      * To get the inverse of A, we form a sequential BAIJ matrix from it, consisting of a single block with block size
128d2567f34SHong Zhang      * equal to the dimension of A, and then use MatInvertBlockDiagonal(). */
129d2567f34SHong Zhang     Mat                A_baij;
130d2567f34SHong Zhang     PetscInt           idxm[1] = {0}, idxn[1] = {0};
131d2567f34SHong Zhang     const PetscScalar *A_inv;
132d2567f34SHong Zhang 
1339566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF, nstages, nstages, nstages, 1, NULL, &A_baij));
1349566063dSJacob Faibussowitsch     PetscCall(MatSetOption(A_baij, MAT_ROW_ORIENTED, PETSC_FALSE));
1359566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(A_baij, 1, idxm, 1, idxn, gauss_A, INSERT_VALUES));
1369566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A_baij, MAT_FINAL_ASSEMBLY));
1379566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A_baij, MAT_FINAL_ASSEMBLY));
1389566063dSJacob Faibussowitsch     PetscCall(MatInvertBlockDiagonal(A_baij, &A_inv));
1399566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(gauss_A_inv, A_inv, nstages * nstages * sizeof(PetscScalar)));
1409566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A_baij));
141d2567f34SHong Zhang   }
142d2567f34SHong Zhang 
143d2567f34SHong Zhang   /* Compute row sums A_inv_rowsum and identity I_s */
144d2567f34SHong Zhang   for (i = 0; i < nstages; i++) {
145d2567f34SHong Zhang     gauss_A_inv_rowsum[i] = 0;
146d2567f34SHong Zhang     for (j = 0; j < nstages; j++) {
147d2567f34SHong Zhang       gauss_A_inv_rowsum[i] += gauss_A_inv[i + nstages * j];
148d2567f34SHong Zhang       I_s[i + nstages * j] = 1. * (i == j);
149d2567f34SHong Zhang     }
150d2567f34SHong Zhang   }
1519566063dSJacob Faibussowitsch   PetscCall(TSIRKTableauCreate(ts, nstages, gauss_A_real, gauss_b, gauss_c, NULL, gauss_A_inv, gauss_A_inv_rowsum, I_s));
1529566063dSJacob Faibussowitsch   PetscCall(PetscFree3(gauss_A_real, gauss_b, gauss_c));
1539566063dSJacob Faibussowitsch   PetscCall(PetscFree4(gauss_A, gauss_A_inv, gauss_A_inv_rowsum, I_s));
154d2567f34SHong Zhang   PetscFunctionReturn(0);
155d2567f34SHong Zhang }
156d2567f34SHong Zhang 
157d2567f34SHong Zhang /*@C
158d2567f34SHong Zhang    TSIRKRegister -  adds a TSIRK implementation
159d2567f34SHong Zhang 
160d2567f34SHong Zhang    Not Collective
161d2567f34SHong Zhang 
162d2567f34SHong Zhang    Input Parameters:
163d2567f34SHong Zhang +  sname - name of user-defined IRK scheme
164d2567f34SHong Zhang -  function - function to create method context
165d2567f34SHong Zhang 
166d2567f34SHong Zhang    Notes:
167d2567f34SHong Zhang    TSIRKRegister() may be called multiple times to add several user-defined families.
168d2567f34SHong Zhang 
169d2567f34SHong Zhang    Sample usage:
170d2567f34SHong Zhang .vb
171d2567f34SHong Zhang    TSIRKRegister("my_scheme",MySchemeCreate);
172d2567f34SHong Zhang .ve
173d2567f34SHong Zhang 
174d2567f34SHong Zhang    Then, your scheme can be chosen with the procedural interface via
175d2567f34SHong Zhang $     TSIRKSetType(ts,"my_scheme")
176d2567f34SHong Zhang    or at runtime via the option
177d2567f34SHong Zhang $     -ts_irk_type my_scheme
178d2567f34SHong Zhang 
179d2567f34SHong Zhang    Level: advanced
180d2567f34SHong Zhang 
181db781477SPatrick Sanan .seealso: `TSIRKRegisterAll()`
182d2567f34SHong Zhang @*/
1839371c9d4SSatish Balay PetscErrorCode TSIRKRegister(const char sname[], PetscErrorCode (*function)(TS)) {
184d2567f34SHong Zhang   PetscFunctionBegin;
1859566063dSJacob Faibussowitsch   PetscCall(TSIRKInitializePackage());
1869566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&TSIRKList, sname, function));
187d2567f34SHong Zhang   PetscFunctionReturn(0);
188d2567f34SHong Zhang }
189d2567f34SHong Zhang 
190d2567f34SHong Zhang /*@C
191d2567f34SHong Zhang   TSIRKRegisterAll - Registers all of the implicit Runge-Kutta methods in TSIRK
192d2567f34SHong Zhang 
193d2567f34SHong Zhang   Not Collective, but should be called by all processes which will need the schemes to be registered
194d2567f34SHong Zhang 
195d2567f34SHong Zhang   Level: advanced
196d2567f34SHong Zhang 
197db781477SPatrick Sanan .seealso: `TSIRKRegisterDestroy()`
198d2567f34SHong Zhang @*/
1999371c9d4SSatish Balay PetscErrorCode TSIRKRegisterAll(void) {
200d2567f34SHong Zhang   PetscFunctionBegin;
201d2567f34SHong Zhang   if (TSIRKRegisterAllCalled) PetscFunctionReturn(0);
202d2567f34SHong Zhang   TSIRKRegisterAllCalled = PETSC_TRUE;
203d2567f34SHong Zhang 
2049566063dSJacob Faibussowitsch   PetscCall(TSIRKRegister(TSIRKGAUSS, TSIRKCreate_Gauss));
205d2567f34SHong Zhang   PetscFunctionReturn(0);
206d2567f34SHong Zhang }
207d2567f34SHong Zhang 
208d2567f34SHong Zhang /*@C
209d2567f34SHong Zhang    TSIRKRegisterDestroy - Frees the list of schemes that were registered by TSIRKRegister().
210d2567f34SHong Zhang 
211d2567f34SHong Zhang    Not Collective
212d2567f34SHong Zhang 
213d2567f34SHong Zhang    Level: advanced
214d2567f34SHong Zhang 
215db781477SPatrick Sanan .seealso: `TSIRKRegister()`, `TSIRKRegisterAll()`
216d2567f34SHong Zhang @*/
2179371c9d4SSatish Balay PetscErrorCode TSIRKRegisterDestroy(void) {
218d2567f34SHong Zhang   PetscFunctionBegin;
219d2567f34SHong Zhang   TSIRKRegisterAllCalled = PETSC_FALSE;
220d2567f34SHong Zhang   PetscFunctionReturn(0);
221d2567f34SHong Zhang }
222d2567f34SHong Zhang 
223d2567f34SHong Zhang /*@C
224d2567f34SHong Zhang   TSIRKInitializePackage - This function initializes everything in the TSIRK package. It is called
225d2567f34SHong Zhang   from TSInitializePackage().
226d2567f34SHong Zhang 
227d2567f34SHong Zhang   Level: developer
228d2567f34SHong Zhang 
229db781477SPatrick Sanan .seealso: `PetscInitialize()`
230d2567f34SHong Zhang @*/
2319371c9d4SSatish Balay PetscErrorCode TSIRKInitializePackage(void) {
232d2567f34SHong Zhang   PetscFunctionBegin;
233d2567f34SHong Zhang   if (TSIRKPackageInitialized) PetscFunctionReturn(0);
234d2567f34SHong Zhang   TSIRKPackageInitialized = PETSC_TRUE;
2359566063dSJacob Faibussowitsch   PetscCall(TSIRKRegisterAll());
2369566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(TSIRKFinalizePackage));
237d2567f34SHong Zhang   PetscFunctionReturn(0);
238d2567f34SHong Zhang }
239d2567f34SHong Zhang 
240d2567f34SHong Zhang /*@C
241d2567f34SHong Zhang   TSIRKFinalizePackage - This function destroys everything in the TSIRK package. It is
242d2567f34SHong Zhang   called from PetscFinalize().
243d2567f34SHong Zhang 
244d2567f34SHong Zhang   Level: developer
245d2567f34SHong Zhang 
246db781477SPatrick Sanan .seealso: `PetscFinalize()`
247d2567f34SHong Zhang @*/
2489371c9d4SSatish Balay PetscErrorCode TSIRKFinalizePackage(void) {
249d2567f34SHong Zhang   PetscFunctionBegin;
2509566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&TSIRKList));
251d2567f34SHong Zhang   TSIRKPackageInitialized = PETSC_FALSE;
252d2567f34SHong Zhang   PetscFunctionReturn(0);
253d2567f34SHong Zhang }
254d2567f34SHong Zhang 
255d2567f34SHong Zhang /*
256d2567f34SHong Zhang  This function can be called before or after ts->vec_sol has been updated.
257d2567f34SHong Zhang */
2589371c9d4SSatish Balay static PetscErrorCode TSEvaluateStep_IRK(TS ts, PetscInt order, Vec U, PetscBool *done) {
259d2567f34SHong Zhang   TS_IRK      *irk   = (TS_IRK *)ts->data;
260d2567f34SHong Zhang   IRKTableau   tab   = irk->tableau;
261d2567f34SHong Zhang   Vec         *YdotI = irk->YdotI;
262d2567f34SHong Zhang   PetscScalar *w     = irk->work;
263d2567f34SHong Zhang   PetscReal    h;
264d2567f34SHong Zhang   PetscInt     j;
265d2567f34SHong Zhang 
266d2567f34SHong Zhang   PetscFunctionBegin;
267d2567f34SHong Zhang   switch (irk->status) {
268d2567f34SHong Zhang   case TS_STEP_INCOMPLETE:
2699371c9d4SSatish Balay   case TS_STEP_PENDING: h = ts->time_step; break;
2709371c9d4SSatish Balay   case TS_STEP_COMPLETE: h = ts->ptime - ts->ptime_prev; break;
271d2567f34SHong Zhang   default: SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus");
272d2567f34SHong Zhang   }
273d2567f34SHong Zhang 
2749566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol, U));
275d2567f34SHong Zhang   for (j = 0; j < irk->nstages; j++) w[j] = h * tab->b[j];
2769566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(U, irk->nstages, w, YdotI));
277d2567f34SHong Zhang   PetscFunctionReturn(0);
278d2567f34SHong Zhang }
279d2567f34SHong Zhang 
2809371c9d4SSatish Balay static PetscErrorCode TSRollBack_IRK(TS ts) {
281d2567f34SHong Zhang   TS_IRK *irk = (TS_IRK *)ts->data;
282d2567f34SHong Zhang 
283d2567f34SHong Zhang   PetscFunctionBegin;
2849566063dSJacob Faibussowitsch   PetscCall(VecCopy(irk->U0, ts->vec_sol));
285d2567f34SHong Zhang   PetscFunctionReturn(0);
286d2567f34SHong Zhang }
287d2567f34SHong Zhang 
2889371c9d4SSatish Balay static PetscErrorCode TSStep_IRK(TS ts) {
289d2567f34SHong Zhang   TS_IRK        *irk   = (TS_IRK *)ts->data;
290d2567f34SHong Zhang   IRKTableau     tab   = irk->tableau;
291d2567f34SHong Zhang   PetscScalar   *A_inv = tab->A_inv, *A_inv_rowsum = tab->A_inv_rowsum;
292d2567f34SHong Zhang   const PetscInt nstages = irk->nstages;
293d2567f34SHong Zhang   SNES           snes;
294d2567f34SHong Zhang   PetscInt       i, j, its, lits, bs;
295d2567f34SHong Zhang   TSAdapt        adapt;
296d2567f34SHong Zhang   PetscInt       rejections     = 0;
297d2567f34SHong Zhang   PetscBool      accept         = PETSC_TRUE;
298d2567f34SHong Zhang   PetscReal      next_time_step = ts->time_step;
299d2567f34SHong Zhang 
300d2567f34SHong Zhang   PetscFunctionBegin;
301*48a46eb9SPierre Jolivet   if (!ts->steprollback) PetscCall(VecCopy(ts->vec_sol, irk->U0));
3029566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(ts->vec_sol, &bs));
303*48a46eb9SPierre Jolivet   for (i = 0; i < nstages; i++) PetscCall(VecStrideScatter(ts->vec_sol, i * bs, irk->Z, INSERT_VALUES));
304d2567f34SHong Zhang 
305d2567f34SHong Zhang   irk->status = TS_STEP_INCOMPLETE;
306d2567f34SHong Zhang   while (!ts->reason && irk->status != TS_STEP_COMPLETE) {
3079566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vec_sol, irk->U));
3089566063dSJacob Faibussowitsch     PetscCall(TSGetSNES(ts, &snes));
3099566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes, NULL, irk->Z));
3109566063dSJacob Faibussowitsch     PetscCall(SNESGetIterationNumber(snes, &its));
3119566063dSJacob Faibussowitsch     PetscCall(SNESGetLinearSolveIterations(snes, &lits));
3129371c9d4SSatish Balay     ts->snes_its += its;
3139371c9d4SSatish Balay     ts->ksp_its += lits;
3149566063dSJacob Faibussowitsch     PetscCall(VecStrideGatherAll(irk->Z, irk->Y, INSERT_VALUES));
315d2567f34SHong Zhang     for (i = 0; i < nstages; i++) {
3169566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(irk->YdotI[i]));
317*48a46eb9SPierre Jolivet       for (j = 0; j < nstages; j++) PetscCall(VecAXPY(irk->YdotI[i], A_inv[i + j * nstages] / ts->time_step, irk->Y[j]));
3189566063dSJacob Faibussowitsch       PetscCall(VecAXPY(irk->YdotI[i], -A_inv_rowsum[i] / ts->time_step, irk->U));
319d2567f34SHong Zhang     }
320d2567f34SHong Zhang     irk->status = TS_STEP_INCOMPLETE;
3219566063dSJacob Faibussowitsch     PetscCall(TSEvaluateStep_IRK(ts, irk->order, ts->vec_sol, NULL));
322d2567f34SHong Zhang     irk->status = TS_STEP_PENDING;
3239566063dSJacob Faibussowitsch     PetscCall(TSGetAdapt(ts, &adapt));
3249566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
325d2567f34SHong Zhang     irk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
326d2567f34SHong Zhang     if (!accept) {
3279566063dSJacob Faibussowitsch       PetscCall(TSRollBack_IRK(ts));
328d2567f34SHong Zhang       ts->time_step = next_time_step;
329d2567f34SHong Zhang       goto reject_step;
330d2567f34SHong Zhang     }
331d2567f34SHong Zhang 
332d2567f34SHong Zhang     ts->ptime += ts->time_step;
333d2567f34SHong Zhang     ts->time_step = next_time_step;
334d2567f34SHong Zhang     break;
335d2567f34SHong Zhang   reject_step:
3369371c9d4SSatish Balay     ts->reject++;
3379371c9d4SSatish Balay     accept = PETSC_FALSE;
338d2567f34SHong Zhang     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
339d2567f34SHong Zhang       ts->reason = TS_DIVERGED_STEP_REJECTED;
34063a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
341d2567f34SHong Zhang     }
342d2567f34SHong Zhang   }
343d2567f34SHong Zhang   PetscFunctionReturn(0);
344d2567f34SHong Zhang }
345d2567f34SHong Zhang 
3469371c9d4SSatish Balay static PetscErrorCode TSInterpolate_IRK(TS ts, PetscReal itime, Vec U) {
347d2567f34SHong Zhang   TS_IRK          *irk     = (TS_IRK *)ts->data;
348d2567f34SHong Zhang   PetscInt         nstages = irk->nstages, pinterp = irk->pinterp, i, j;
349d2567f34SHong Zhang   PetscReal        h;
350d2567f34SHong Zhang   PetscReal        tt, t;
351d2567f34SHong Zhang   PetscScalar     *bt;
352d2567f34SHong Zhang   const PetscReal *B = irk->tableau->binterp;
353d2567f34SHong Zhang 
354d2567f34SHong Zhang   PetscFunctionBegin;
3553c633725SBarry Smith   PetscCheck(B, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSIRK %s does not have an interpolation formula", irk->method_name);
356d2567f34SHong Zhang   switch (irk->status) {
357d2567f34SHong Zhang   case TS_STEP_INCOMPLETE:
358d2567f34SHong Zhang   case TS_STEP_PENDING:
359d2567f34SHong Zhang     h = ts->time_step;
360d2567f34SHong Zhang     t = (itime - ts->ptime) / h;
361d2567f34SHong Zhang     break;
362d2567f34SHong Zhang   case TS_STEP_COMPLETE:
363d2567f34SHong Zhang     h = ts->ptime - ts->ptime_prev;
364d2567f34SHong Zhang     t = (itime - ts->ptime) / h + 1; /* In the interval [0,1] */
365d2567f34SHong Zhang     break;
366d2567f34SHong Zhang   default: SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus");
367d2567f34SHong Zhang   }
3689566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nstages, &bt));
369d2567f34SHong Zhang   for (i = 0; i < nstages; i++) bt[i] = 0;
370d2567f34SHong Zhang   for (j = 0, tt = t; j < pinterp; j++, tt *= t) {
3719371c9d4SSatish Balay     for (i = 0; i < nstages; i++) { bt[i] += h * B[i * pinterp + j] * tt; }
372d2567f34SHong Zhang   }
3739566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(U, nstages, bt, irk->YdotI));
374d2567f34SHong Zhang   PetscFunctionReturn(0);
375d2567f34SHong Zhang }
376d2567f34SHong Zhang 
3779371c9d4SSatish Balay static PetscErrorCode TSIRKTableauReset(TS ts) {
378d2567f34SHong Zhang   TS_IRK    *irk = (TS_IRK *)ts->data;
379d2567f34SHong Zhang   IRKTableau tab = irk->tableau;
380d2567f34SHong Zhang 
381d2567f34SHong Zhang   PetscFunctionBegin;
382d2567f34SHong Zhang   if (!tab) PetscFunctionReturn(0);
3839566063dSJacob Faibussowitsch   PetscCall(PetscFree3(tab->A, tab->A_inv, tab->I_s));
3849566063dSJacob Faibussowitsch   PetscCall(PetscFree4(tab->b, tab->c, tab->binterp, tab->A_inv_rowsum));
385d2567f34SHong Zhang   PetscFunctionReturn(0);
386d2567f34SHong Zhang }
387d2567f34SHong Zhang 
3889371c9d4SSatish Balay static PetscErrorCode TSReset_IRK(TS ts) {
389d2567f34SHong Zhang   TS_IRK *irk = (TS_IRK *)ts->data;
390d2567f34SHong Zhang 
391d2567f34SHong Zhang   PetscFunctionBegin;
3929566063dSJacob Faibussowitsch   PetscCall(TSIRKTableauReset(ts));
3931baa6e33SBarry Smith   if (irk->tableau) PetscCall(PetscFree(irk->tableau));
3941baa6e33SBarry Smith   if (irk->method_name) PetscCall(PetscFree(irk->method_name));
3951baa6e33SBarry Smith   if (irk->work) PetscCall(PetscFree(irk->work));
3969566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(irk->nstages, &irk->Y));
3979566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(irk->nstages, &irk->YdotI));
3989566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&irk->Ydot));
3999566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&irk->Z));
4009566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&irk->U));
4019566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&irk->U0));
4029566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&irk->TJ));
403d2567f34SHong Zhang   PetscFunctionReturn(0);
404d2567f34SHong Zhang }
405d2567f34SHong Zhang 
4069371c9d4SSatish Balay static PetscErrorCode TSIRKGetVecs(TS ts, DM dm, Vec *U) {
407d2567f34SHong Zhang   TS_IRK *irk = (TS_IRK *)ts->data;
408d2567f34SHong Zhang 
409d2567f34SHong Zhang   PetscFunctionBegin;
410d2567f34SHong Zhang   if (U) {
411d2567f34SHong Zhang     if (dm && dm != ts->dm) {
4129566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm, "TSIRK_U", U));
413d2567f34SHong Zhang     } else *U = irk->U;
414d2567f34SHong Zhang   }
415d2567f34SHong Zhang   PetscFunctionReturn(0);
416d2567f34SHong Zhang }
417d2567f34SHong Zhang 
4189371c9d4SSatish Balay static PetscErrorCode TSIRKRestoreVecs(TS ts, DM dm, Vec *U) {
419d2567f34SHong Zhang   PetscFunctionBegin;
420d2567f34SHong Zhang   if (U) {
421*48a46eb9SPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSIRK_U", U));
422d2567f34SHong Zhang   }
423d2567f34SHong Zhang   PetscFunctionReturn(0);
424d2567f34SHong Zhang }
425d2567f34SHong Zhang 
426d2567f34SHong Zhang /*
427d2567f34SHong Zhang   This defines the nonlinear equations that is to be solved with SNES
428d2567f34SHong Zhang     G[e\otimes t + C*dt, Z, Zdot] = 0
429d2567f34SHong Zhang     Zdot = (In \otimes S)*Z - (In \otimes Se) U
430d2567f34SHong Zhang   where S = 1/(dt*A)
431d2567f34SHong Zhang */
4329371c9d4SSatish Balay static PetscErrorCode SNESTSFormFunction_IRK(SNES snes, Vec ZC, Vec FC, TS ts) {
433d2567f34SHong Zhang   TS_IRK            *irk     = (TS_IRK *)ts->data;
434d2567f34SHong Zhang   IRKTableau         tab     = irk->tableau;
435d2567f34SHong Zhang   const PetscInt     nstages = irk->nstages;
436d2567f34SHong Zhang   const PetscReal   *c       = tab->c;
437d2567f34SHong Zhang   const PetscScalar *A_inv = tab->A_inv, *A_inv_rowsum = tab->A_inv_rowsum;
438d2567f34SHong Zhang   DM                 dm, dmsave;
439d2567f34SHong Zhang   Vec                U, *YdotI = irk->YdotI, Ydot = irk->Ydot, *Y = irk->Y;
440d2567f34SHong Zhang   PetscReal          h = ts->time_step;
441d2567f34SHong Zhang   PetscInt           i, j;
442d2567f34SHong Zhang 
443d2567f34SHong Zhang   PetscFunctionBegin;
4449566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
4459566063dSJacob Faibussowitsch   PetscCall(TSIRKGetVecs(ts, dm, &U));
4469566063dSJacob Faibussowitsch   PetscCall(VecStrideGatherAll(ZC, Y, INSERT_VALUES));
447d2567f34SHong Zhang   dmsave = ts->dm;
448d2567f34SHong Zhang   ts->dm = dm;
449d2567f34SHong Zhang   for (i = 0; i < nstages; i++) {
4509566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Ydot));
451*48a46eb9SPierre Jolivet     for (j = 0; j < nstages; j++) PetscCall(VecAXPY(Ydot, A_inv[j * nstages + i] / h, Y[j]));
4529566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Ydot, -A_inv_rowsum[i] / h, U)); /* Ydot = (S \otimes In)*Z - (Se \otimes In) U */
4539566063dSJacob Faibussowitsch     PetscCall(TSComputeIFunction(ts, ts->ptime + ts->time_step * c[i], Y[i], Ydot, YdotI[i], PETSC_FALSE));
454d2567f34SHong Zhang   }
4559566063dSJacob Faibussowitsch   PetscCall(VecStrideScatterAll(YdotI, FC, INSERT_VALUES));
456d2567f34SHong Zhang   ts->dm = dmsave;
4579566063dSJacob Faibussowitsch   PetscCall(TSIRKRestoreVecs(ts, dm, &U));
458d2567f34SHong Zhang   PetscFunctionReturn(0);
459d2567f34SHong Zhang }
460d2567f34SHong Zhang 
461d2567f34SHong Zhang /*
462d2567f34SHong Zhang    For explicit ODE, the Jacobian is
463d2567f34SHong Zhang      JC = I_n \otimes S - J \otimes I_s
464d2567f34SHong Zhang    For DAE, the Jacobian is
465d2567f34SHong Zhang      JC = M_n \otimes S - J \otimes I_s
466d2567f34SHong Zhang */
4679371c9d4SSatish Balay static PetscErrorCode SNESTSFormJacobian_IRK(SNES snes, Vec ZC, Mat JC, Mat JCpre, TS ts) {
468d2567f34SHong Zhang   TS_IRK          *irk     = (TS_IRK *)ts->data;
469d2567f34SHong Zhang   IRKTableau       tab     = irk->tableau;
470d2567f34SHong Zhang   const PetscInt   nstages = irk->nstages;
471d2567f34SHong Zhang   const PetscReal *c       = tab->c;
472d2567f34SHong Zhang   DM               dm, dmsave;
473d2567f34SHong Zhang   Vec             *Y = irk->Y, Ydot = irk->Ydot;
474d2567f34SHong Zhang   Mat              J;
475d2567f34SHong Zhang   PetscScalar     *S;
476d2567f34SHong Zhang   PetscInt         i, j, bs;
477d2567f34SHong Zhang 
478d2567f34SHong Zhang   PetscFunctionBegin;
4799566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
480d2567f34SHong Zhang   /* irk->Ydot has already been computed in SNESTSFormFunction_IRK (SNES guarantees this) */
481d2567f34SHong Zhang   dmsave = ts->dm;
482d2567f34SHong Zhang   ts->dm = dm;
4839566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(Y[nstages - 1], &bs));
484d2567f34SHong Zhang   if (ts->equation_type <= TS_EQ_ODE_EXPLICIT) { /* Support explicit formulas only */
4859566063dSJacob Faibussowitsch     PetscCall(VecStrideGather(ZC, (nstages - 1) * bs, Y[nstages - 1], INSERT_VALUES));
4869566063dSJacob Faibussowitsch     PetscCall(MatKAIJGetAIJ(JC, &J));
4879566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, ts->ptime + ts->time_step * c[nstages - 1], Y[nstages - 1], Ydot, 0, J, J, PETSC_FALSE));
4889566063dSJacob Faibussowitsch     PetscCall(MatKAIJGetS(JC, NULL, NULL, &S));
489d2567f34SHong Zhang     for (i = 0; i < nstages; i++)
4909371c9d4SSatish Balay       for (j = 0; j < nstages; j++) S[i + nstages * j] = tab->A_inv[i + nstages * j] / ts->time_step;
4919566063dSJacob Faibussowitsch     PetscCall(MatKAIJRestoreS(JC, &S));
49298921bdaSJacob 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  */
493d2567f34SHong Zhang   ts->dm = dmsave;
494d2567f34SHong Zhang   PetscFunctionReturn(0);
495d2567f34SHong Zhang }
496d2567f34SHong Zhang 
4979371c9d4SSatish Balay static PetscErrorCode DMCoarsenHook_TSIRK(DM fine, DM coarse, void *ctx) {
498d2567f34SHong Zhang   PetscFunctionBegin;
499d2567f34SHong Zhang   PetscFunctionReturn(0);
500d2567f34SHong Zhang }
501d2567f34SHong Zhang 
5029371c9d4SSatish Balay static PetscErrorCode DMRestrictHook_TSIRK(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) {
503d2567f34SHong Zhang   TS  ts = (TS)ctx;
504d2567f34SHong Zhang   Vec U, U_c;
505d2567f34SHong Zhang 
506d2567f34SHong Zhang   PetscFunctionBegin;
5079566063dSJacob Faibussowitsch   PetscCall(TSIRKGetVecs(ts, fine, &U));
5089566063dSJacob Faibussowitsch   PetscCall(TSIRKGetVecs(ts, coarse, &U_c));
5099566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct, U, U_c));
5109566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(U_c, rscale, U_c));
5119566063dSJacob Faibussowitsch   PetscCall(TSIRKRestoreVecs(ts, fine, &U));
5129566063dSJacob Faibussowitsch   PetscCall(TSIRKRestoreVecs(ts, coarse, &U_c));
513d2567f34SHong Zhang   PetscFunctionReturn(0);
514d2567f34SHong Zhang }
515d2567f34SHong Zhang 
5169371c9d4SSatish Balay static PetscErrorCode DMSubDomainHook_TSIRK(DM dm, DM subdm, void *ctx) {
517d2567f34SHong Zhang   PetscFunctionBegin;
518d2567f34SHong Zhang   PetscFunctionReturn(0);
519d2567f34SHong Zhang }
520d2567f34SHong Zhang 
5219371c9d4SSatish Balay static PetscErrorCode DMSubDomainRestrictHook_TSIRK(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) {
522d2567f34SHong Zhang   TS  ts = (TS)ctx;
523d2567f34SHong Zhang   Vec U, U_c;
524d2567f34SHong Zhang 
525d2567f34SHong Zhang   PetscFunctionBegin;
5269566063dSJacob Faibussowitsch   PetscCall(TSIRKGetVecs(ts, dm, &U));
5279566063dSJacob Faibussowitsch   PetscCall(TSIRKGetVecs(ts, subdm, &U_c));
528d2567f34SHong Zhang 
5299566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(gscat, U, U_c, INSERT_VALUES, SCATTER_FORWARD));
5309566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(gscat, U, U_c, INSERT_VALUES, SCATTER_FORWARD));
531d2567f34SHong Zhang 
5329566063dSJacob Faibussowitsch   PetscCall(TSIRKRestoreVecs(ts, dm, &U));
5339566063dSJacob Faibussowitsch   PetscCall(TSIRKRestoreVecs(ts, subdm, &U_c));
534d2567f34SHong Zhang   PetscFunctionReturn(0);
535d2567f34SHong Zhang }
536d2567f34SHong Zhang 
5379371c9d4SSatish Balay static PetscErrorCode TSSetUp_IRK(TS ts) {
538d2567f34SHong Zhang   TS_IRK        *irk = (TS_IRK *)ts->data;
539d2567f34SHong Zhang   IRKTableau     tab = irk->tableau;
540d2567f34SHong Zhang   DM             dm;
541d2567f34SHong Zhang   Mat            J;
542d2567f34SHong Zhang   Vec            R;
543d2567f34SHong Zhang   const PetscInt nstages = irk->nstages;
544d2567f34SHong Zhang   PetscInt       vsize, bs;
545d2567f34SHong Zhang 
546d2567f34SHong Zhang   PetscFunctionBegin;
547*48a46eb9SPierre Jolivet   if (!irk->work) PetscCall(PetscMalloc1(irk->nstages, &irk->work));
548*48a46eb9SPierre Jolivet   if (!irk->Y) PetscCall(VecDuplicateVecs(ts->vec_sol, irk->nstages, &irk->Y));
549*48a46eb9SPierre Jolivet   if (!irk->YdotI) PetscCall(VecDuplicateVecs(ts->vec_sol, irk->nstages, &irk->YdotI));
550*48a46eb9SPierre Jolivet   if (!irk->Ydot) PetscCall(VecDuplicate(ts->vec_sol, &irk->Ydot));
551*48a46eb9SPierre Jolivet   if (!irk->U) PetscCall(VecDuplicate(ts->vec_sol, &irk->U));
552*48a46eb9SPierre Jolivet   if (!irk->U0) PetscCall(VecDuplicate(ts->vec_sol, &irk->U0));
553d2567f34SHong Zhang   if (!irk->Z) {
5549566063dSJacob Faibussowitsch     PetscCall(VecCreate(PetscObjectComm((PetscObject)ts->vec_sol), &irk->Z));
5559566063dSJacob Faibussowitsch     PetscCall(VecGetSize(ts->vec_sol, &vsize));
5569566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(irk->Z, PETSC_DECIDE, vsize * irk->nstages));
5579566063dSJacob Faibussowitsch     PetscCall(VecGetBlockSize(ts->vec_sol, &bs));
5589566063dSJacob Faibussowitsch     PetscCall(VecSetBlockSize(irk->Z, irk->nstages * bs));
5599566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(irk->Z));
560d2567f34SHong Zhang   }
5619566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
5629566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSIRK, DMRestrictHook_TSIRK, ts));
5639566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSIRK, DMSubDomainRestrictHook_TSIRK, ts));
564d2567f34SHong Zhang 
5659566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &ts->snes));
5669566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(irk->Z, &R));
5679566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(ts->snes, R, SNESTSFormFunction, ts));
5689566063dSJacob Faibussowitsch   PetscCall(TSGetIJacobian(ts, &J, NULL, NULL, NULL));
569d2567f34SHong Zhang   if (!irk->TJ) {
570d2567f34SHong Zhang     /* Create the KAIJ matrix for solving the stages */
5719566063dSJacob Faibussowitsch     PetscCall(MatCreateKAIJ(J, nstages, nstages, tab->A_inv, tab->I_s, &irk->TJ));
572d2567f34SHong Zhang   }
5739566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(ts->snes, irk->TJ, irk->TJ, SNESTSFormJacobian, ts));
5749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&R));
575d2567f34SHong Zhang   PetscFunctionReturn(0);
576d2567f34SHong Zhang }
577d2567f34SHong Zhang 
5789371c9d4SSatish Balay static PetscErrorCode TSSetFromOptions_IRK(TS ts, PetscOptionItems *PetscOptionsObject) {
579d2567f34SHong Zhang   TS_IRK *irk        = (TS_IRK *)ts->data;
580d2567f34SHong Zhang   char    tname[256] = TSIRKGAUSS;
581d2567f34SHong Zhang 
582d2567f34SHong Zhang   PetscFunctionBegin;
583d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "IRK ODE solver options");
584d2567f34SHong Zhang   {
585d2567f34SHong Zhang     PetscBool flg1, flg2;
5869566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-ts_irk_nstages", "Stages of the IRK method", "TSIRKSetNumStages", irk->nstages, &irk->nstages, &flg1));
5879566063dSJacob Faibussowitsch     PetscCall(PetscOptionsFList("-ts_irk_type", "Type of IRK method", "TSIRKSetType", TSIRKList, irk->method_name[0] ? irk->method_name : tname, tname, sizeof(tname), &flg2));
588d2567f34SHong Zhang     if (flg1 || flg2 || !irk->method_name[0]) { /* Create the method tableau after nstages or method is set */
5899566063dSJacob Faibussowitsch       PetscCall(TSIRKSetType(ts, tname));
590d2567f34SHong Zhang     }
591d2567f34SHong Zhang   }
592d0609cedSBarry Smith   PetscOptionsHeadEnd();
593d2567f34SHong Zhang   PetscFunctionReturn(0);
594d2567f34SHong Zhang }
595d2567f34SHong Zhang 
5969371c9d4SSatish Balay static PetscErrorCode TSView_IRK(TS ts, PetscViewer viewer) {
597d2567f34SHong Zhang   TS_IRK   *irk = (TS_IRK *)ts->data;
598d2567f34SHong Zhang   PetscBool iascii;
599d2567f34SHong Zhang 
600d2567f34SHong Zhang   PetscFunctionBegin;
6019566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
602d2567f34SHong Zhang   if (iascii) {
603d2567f34SHong Zhang     IRKTableau tab = irk->tableau;
604d2567f34SHong Zhang     TSIRKType  irktype;
605d2567f34SHong Zhang     char       buf[512];
606d2567f34SHong Zhang 
6079566063dSJacob Faibussowitsch     PetscCall(TSIRKGetType(ts, &irktype));
6089566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  IRK type %s\n", irktype));
6099566063dSJacob Faibussowitsch     PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", irk->nstages, tab->c));
6109566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Abscissa       c = %s\n", buf));
6119566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Stiffly accurate: %s\n", irk->stiffly_accurate ? "yes" : "no"));
6129566063dSJacob Faibussowitsch     PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", PetscSqr(irk->nstages), tab->A));
6139566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  A coefficients       A = %s\n", buf));
614d2567f34SHong Zhang   }
615d2567f34SHong Zhang   PetscFunctionReturn(0);
616d2567f34SHong Zhang }
617d2567f34SHong Zhang 
6189371c9d4SSatish Balay static PetscErrorCode TSLoad_IRK(TS ts, PetscViewer viewer) {
619d2567f34SHong Zhang   SNES    snes;
620d2567f34SHong Zhang   TSAdapt adapt;
621d2567f34SHong Zhang 
622d2567f34SHong Zhang   PetscFunctionBegin;
6239566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &adapt));
6249566063dSJacob Faibussowitsch   PetscCall(TSAdaptLoad(adapt, viewer));
6259566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
6269566063dSJacob Faibussowitsch   PetscCall(SNESLoad(snes, viewer));
627d2567f34SHong Zhang   /* function and Jacobian context for SNES when used with TS is always ts object */
6289566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, NULL, NULL, ts));
6299566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, ts));
630d2567f34SHong Zhang   PetscFunctionReturn(0);
631d2567f34SHong Zhang }
632d2567f34SHong Zhang 
633d2567f34SHong Zhang /*@C
634d2567f34SHong Zhang   TSIRKSetType - Set the type of IRK scheme
635d2567f34SHong Zhang 
636d2567f34SHong Zhang   Logically collective
637d2567f34SHong Zhang 
63802024617SSatish Balay   Input Parameters:
639d2567f34SHong Zhang +  ts - timestepping context
640d2567f34SHong Zhang -  irktype - type of IRK scheme
641d2567f34SHong Zhang 
642d2567f34SHong Zhang   Options Database:
64367b8a455SSatish Balay .  -ts_irk_type <gauss> - set irk type
644d2567f34SHong Zhang 
645d2567f34SHong Zhang   Level: intermediate
646d2567f34SHong Zhang 
647db781477SPatrick Sanan .seealso: `TSIRKGetType()`, `TSIRK`, `TSIRKType`, `TSIRKGAUSS`
648d2567f34SHong Zhang @*/
6499371c9d4SSatish Balay PetscErrorCode TSIRKSetType(TS ts, TSIRKType irktype) {
650d2567f34SHong Zhang   PetscFunctionBegin;
651d2567f34SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
652d2567f34SHong Zhang   PetscValidCharPointer(irktype, 2);
653cac4c232SBarry Smith   PetscTryMethod(ts, "TSIRKSetType_C", (TS, TSIRKType), (ts, irktype));
654d2567f34SHong Zhang   PetscFunctionReturn(0);
655d2567f34SHong Zhang }
656d2567f34SHong Zhang 
657d2567f34SHong Zhang /*@C
658d2567f34SHong Zhang   TSIRKGetType - Get the type of IRK IMEX scheme
659d2567f34SHong Zhang 
660d2567f34SHong Zhang   Logically collective
661d2567f34SHong Zhang 
662d2567f34SHong Zhang   Input Parameter:
663d2567f34SHong Zhang .  ts - timestepping context
664d2567f34SHong Zhang 
665d2567f34SHong Zhang   Output Parameter:
666d2567f34SHong Zhang .  irktype - type of IRK-IMEX scheme
667d2567f34SHong Zhang 
668d2567f34SHong Zhang   Level: intermediate
669d2567f34SHong Zhang 
670db781477SPatrick Sanan .seealso: `TSIRKGetType()`
671d2567f34SHong Zhang @*/
6729371c9d4SSatish Balay PetscErrorCode TSIRKGetType(TS ts, TSIRKType *irktype) {
673d2567f34SHong Zhang   PetscFunctionBegin;
674d2567f34SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
675cac4c232SBarry Smith   PetscUseMethod(ts, "TSIRKGetType_C", (TS, TSIRKType *), (ts, irktype));
676d2567f34SHong Zhang   PetscFunctionReturn(0);
677d2567f34SHong Zhang }
678d2567f34SHong Zhang 
679d2567f34SHong Zhang /*@C
680d2567f34SHong Zhang   TSIRKSetNumStages - Set the number of stages of IRK scheme
681d2567f34SHong Zhang 
682d2567f34SHong Zhang   Logically collective
683d2567f34SHong Zhang 
68402024617SSatish Balay   Input Parameters:
685d2567f34SHong Zhang +  ts - timestepping context
686d2567f34SHong Zhang -  nstages - number of stages of IRK scheme
687d2567f34SHong Zhang 
688d2567f34SHong Zhang   Options Database:
68967b8a455SSatish Balay .  -ts_irk_nstages <int> - set number of stages
690d2567f34SHong Zhang 
691d2567f34SHong Zhang   Level: intermediate
692d2567f34SHong Zhang 
693db781477SPatrick Sanan .seealso: `TSIRKGetNumStages()`, `TSIRK`
694d2567f34SHong Zhang @*/
6959371c9d4SSatish Balay PetscErrorCode TSIRKSetNumStages(TS ts, PetscInt nstages) {
696d2567f34SHong Zhang   PetscFunctionBegin;
697d2567f34SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
698cac4c232SBarry Smith   PetscTryMethod(ts, "TSIRKSetNumStages_C", (TS, PetscInt), (ts, nstages));
699d2567f34SHong Zhang   PetscFunctionReturn(0);
700d2567f34SHong Zhang }
701d2567f34SHong Zhang 
702d2567f34SHong Zhang /*@C
703d2567f34SHong Zhang   TSIRKGetNumStages - Get the number of stages of IRK scheme
704d2567f34SHong Zhang 
705d2567f34SHong Zhang   Logically collective
706d2567f34SHong Zhang 
70702024617SSatish Balay   Input Parameters:
708d2567f34SHong Zhang +  ts - timestepping context
709d2567f34SHong Zhang -  nstages - number of stages of IRK scheme
710d2567f34SHong Zhang 
711d2567f34SHong Zhang   Level: intermediate
712d2567f34SHong Zhang 
713db781477SPatrick Sanan .seealso: `TSIRKSetNumStages()`, `TSIRK`
714d2567f34SHong Zhang @*/
7159371c9d4SSatish Balay PetscErrorCode TSIRKGetNumStages(TS ts, PetscInt *nstages) {
716d2567f34SHong Zhang   PetscFunctionBegin;
717d2567f34SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
718d2567f34SHong Zhang   PetscValidIntPointer(nstages, 2);
719cac4c232SBarry Smith   PetscTryMethod(ts, "TSIRKGetNumStages_C", (TS, PetscInt *), (ts, nstages));
720d2567f34SHong Zhang   PetscFunctionReturn(0);
721d2567f34SHong Zhang }
722d2567f34SHong Zhang 
7239371c9d4SSatish Balay static PetscErrorCode TSIRKGetType_IRK(TS ts, TSIRKType *irktype) {
724d2567f34SHong Zhang   TS_IRK *irk = (TS_IRK *)ts->data;
725d2567f34SHong Zhang 
726d2567f34SHong Zhang   PetscFunctionBegin;
727d2567f34SHong Zhang   *irktype = irk->method_name;
728d2567f34SHong Zhang   PetscFunctionReturn(0);
729d2567f34SHong Zhang }
730d2567f34SHong Zhang 
7319371c9d4SSatish Balay static PetscErrorCode TSIRKSetType_IRK(TS ts, TSIRKType irktype) {
732d2567f34SHong Zhang   TS_IRK *irk = (TS_IRK *)ts->data;
7335f80ce2aSJacob Faibussowitsch   PetscErrorCode (*irkcreate)(TS);
734d2567f34SHong Zhang 
735d2567f34SHong Zhang   PetscFunctionBegin;
736d2567f34SHong Zhang   if (irk->method_name) {
7379566063dSJacob Faibussowitsch     PetscCall(PetscFree(irk->method_name));
7389566063dSJacob Faibussowitsch     PetscCall(TSIRKTableauReset(ts));
739d2567f34SHong Zhang   }
7409566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(TSIRKList, irktype, &irkcreate));
7413c633725SBarry Smith   PetscCheck(irkcreate, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TSIRK type \"%s\" given", irktype);
7429566063dSJacob Faibussowitsch   PetscCall((*irkcreate)(ts));
7439566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(irktype, &irk->method_name));
744d2567f34SHong Zhang   PetscFunctionReturn(0);
745d2567f34SHong Zhang }
746d2567f34SHong Zhang 
7479371c9d4SSatish Balay static PetscErrorCode TSIRKSetNumStages_IRK(TS ts, PetscInt nstages) {
748d2567f34SHong Zhang   TS_IRK *irk = (TS_IRK *)ts->data;
749d2567f34SHong Zhang 
750d2567f34SHong Zhang   PetscFunctionBegin;
75163a3b9bcSJacob Faibussowitsch   PetscCheck(nstages > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "input argument, %" PetscInt_FMT ", out of range", nstages);
752d2567f34SHong Zhang   irk->nstages = nstages;
753d2567f34SHong Zhang   PetscFunctionReturn(0);
754d2567f34SHong Zhang }
755d2567f34SHong Zhang 
7569371c9d4SSatish Balay static PetscErrorCode TSIRKGetNumStages_IRK(TS ts, PetscInt *nstages) {
757d2567f34SHong Zhang   TS_IRK *irk = (TS_IRK *)ts->data;
758d2567f34SHong Zhang 
759d2567f34SHong Zhang   PetscFunctionBegin;
760d2567f34SHong Zhang   PetscValidIntPointer(nstages, 2);
761d2567f34SHong Zhang   *nstages = irk->nstages;
762d2567f34SHong Zhang   PetscFunctionReturn(0);
763d2567f34SHong Zhang }
764d2567f34SHong Zhang 
7659371c9d4SSatish Balay static PetscErrorCode TSDestroy_IRK(TS ts) {
766d2567f34SHong Zhang   PetscFunctionBegin;
7679566063dSJacob Faibussowitsch   PetscCall(TSReset_IRK(ts));
768d2567f34SHong Zhang   if (ts->dm) {
7699566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSIRK, DMRestrictHook_TSIRK, ts));
7709566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSIRK, DMSubDomainRestrictHook_TSIRK, ts));
771d2567f34SHong Zhang   }
7729566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
7739566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSIRKSetType_C", NULL));
7749566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSIRKGetType_C", NULL));
7759566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSIRKSetNumStages_C", NULL));
7769566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSIRKGetNumStages_C", NULL));
777d2567f34SHong Zhang   PetscFunctionReturn(0);
778d2567f34SHong Zhang }
779d2567f34SHong Zhang 
780d2567f34SHong Zhang /*MC
781d2567f34SHong Zhang       TSIRK - ODE and DAE solver using Implicit Runge-Kutta schemes
782d2567f34SHong Zhang 
783d2567f34SHong Zhang   Notes:
784d2567f34SHong Zhang 
785d2567f34SHong Zhang   TSIRK uses the sparse Kronecker product matrix implementation of MATKAIJ to achieve good arithmetic intensity.
786d2567f34SHong Zhang 
787517f343bSHong 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().
788d2567f34SHong Zhang 
789d2567f34SHong Zhang   Level: beginner
790d2567f34SHong Zhang 
791db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSIRKSetType()`, `TSIRKGetType()`, `TSIRKGAUSS`, `TSIRKRegister()`, `TSIRKSetNumStages()`
792d2567f34SHong Zhang 
793d2567f34SHong Zhang M*/
7949371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TSCreate_IRK(TS ts) {
795d2567f34SHong Zhang   TS_IRK *irk;
796d2567f34SHong Zhang 
797d2567f34SHong Zhang   PetscFunctionBegin;
7989566063dSJacob Faibussowitsch   PetscCall(TSIRKInitializePackage());
799d2567f34SHong Zhang 
800d2567f34SHong Zhang   ts->ops->reset          = TSReset_IRK;
801d2567f34SHong Zhang   ts->ops->destroy        = TSDestroy_IRK;
802d2567f34SHong Zhang   ts->ops->view           = TSView_IRK;
803d2567f34SHong Zhang   ts->ops->load           = TSLoad_IRK;
804d2567f34SHong Zhang   ts->ops->setup          = TSSetUp_IRK;
805d2567f34SHong Zhang   ts->ops->step           = TSStep_IRK;
806d2567f34SHong Zhang   ts->ops->interpolate    = TSInterpolate_IRK;
807d2567f34SHong Zhang   ts->ops->evaluatestep   = TSEvaluateStep_IRK;
808d2567f34SHong Zhang   ts->ops->rollback       = TSRollBack_IRK;
809d2567f34SHong Zhang   ts->ops->setfromoptions = TSSetFromOptions_IRK;
810d2567f34SHong Zhang   ts->ops->snesfunction   = SNESTSFormFunction_IRK;
811d2567f34SHong Zhang   ts->ops->snesjacobian   = SNESTSFormJacobian_IRK;
812d2567f34SHong Zhang 
813d2567f34SHong Zhang   ts->usessnes = PETSC_TRUE;
814d2567f34SHong Zhang 
8159566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(ts, &irk));
816d2567f34SHong Zhang   ts->data = (void *)irk;
817d2567f34SHong Zhang 
8189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSIRKSetType_C", TSIRKSetType_IRK));
8199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSIRKGetType_C", TSIRKGetType_IRK));
8209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSIRKSetNumStages_C", TSIRKSetNumStages_IRK));
8219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSIRKGetNumStages_C", TSIRKGetNumStages_IRK));
822d2567f34SHong Zhang   /* 3-stage IRK_Gauss is the default */
8239566063dSJacob Faibussowitsch   PetscCall(PetscNew(&irk->tableau));
824d2567f34SHong Zhang   irk->nstages = 3;
8259566063dSJacob Faibussowitsch   PetscCall(TSIRKSetType(ts, TSIRKDefault));
826d2567f34SHong Zhang   PetscFunctionReturn(0);
827d2567f34SHong Zhang }
828