xref: /petsc/src/ts/impls/implicit/irk/irk.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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
48bcf0153eSBarry Smith   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
55efa39862SBarry Smith . A            - stage coefficients (dimension `nstages` * `nstages`, row-major)
56efa39862SBarry Smith . b            - step completion table (dimension `nstages`)
57efa39862SBarry Smith . c            - abscissa (dimension `nstages`)
58efa39862SBarry Smith . binterp      - coefficients of the interpolation formula (dimension `nstages`), optional (use `NULL` to skip)
59efa39862SBarry Smith . A_inv        - inverse of `A` (dimension `nstages` * `nstages`, row-major), optional (use `NULL` to skip)
60efa39862SBarry Smith . A_inv_rowsum - row sum of the inverse of `A` (dimension `nstages`), optional (use `NULL` to skip)
61efa39862SBarry Smith - I_s          - identity matrix (dimension `nstages` * `nstages`), optional (use `NULL` to skip)
62d2567f34SHong Zhang 
63d2567f34SHong Zhang   Level: advanced
64d2567f34SHong Zhang 
651cc06b55SBarry Smith .seealso: [](ch_ts), `TSIRK`, `TSIRKRegister()`
66d2567f34SHong Zhang @*/
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)67d71ae5a4SJacob Faibussowitsch 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)
68d71ae5a4SJacob Faibussowitsch {
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 */
801baa6e33SBarry Smith   if (binterp) PetscCall(PetscArraycpy(tab->binterp, binterp, nstages));
811baa6e33SBarry Smith   if (A_inv) PetscCall(PetscArraycpy(tab->A_inv, A_inv, PetscSqr(nstages)));
821baa6e33SBarry Smith   if (A_inv_rowsum) PetscCall(PetscArraycpy(tab->A_inv_rowsum, A_inv_rowsum, nstages));
831baa6e33SBarry Smith   if (I_s) PetscCall(PetscArraycpy(tab->I_s, I_s, PetscSqr(nstages)));
843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
85d2567f34SHong Zhang }
86d2567f34SHong Zhang 
87d2567f34SHong Zhang /* Arrays should be freed with PetscFree3(A,b,c) */
TSIRKCreate_Gauss(TS ts)88d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSIRKCreate_Gauss(TS ts)
89d71ae5a4SJacob Faibussowitsch {
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++)
1209371c9d4SSatish Balay     for (j = 0; j < nstages; j++) gauss_A_real[i * nstages + j] = PetscRealPart(gauss_A[i * nstages + j]);
121d2567f34SHong Zhang 
1229566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&G0mat));
1239566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&G1mat));
1249566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Amat));
1259566063dSJacob Faibussowitsch   PetscCall(PetscFree3(b, G0, G1));
126d2567f34SHong Zhang 
127d2567f34SHong Zhang   { /* Invert A */
128d2567f34SHong Zhang     /* PETSc does not provide a routine to calculate the inverse of a general matrix.
129d2567f34SHong Zhang      * To get the inverse of A, we form a sequential BAIJ matrix from it, consisting of a single block with block size
130d2567f34SHong Zhang      * equal to the dimension of A, and then use MatInvertBlockDiagonal(). */
131d2567f34SHong Zhang     Mat                A_baij;
132d2567f34SHong Zhang     PetscInt           idxm[1] = {0}, idxn[1] = {0};
133d2567f34SHong Zhang     const PetscScalar *A_inv;
134d2567f34SHong Zhang 
1359566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF, nstages, nstages, nstages, 1, NULL, &A_baij));
1369566063dSJacob Faibussowitsch     PetscCall(MatSetOption(A_baij, MAT_ROW_ORIENTED, PETSC_FALSE));
1379566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(A_baij, 1, idxm, 1, idxn, gauss_A, INSERT_VALUES));
1389566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A_baij, MAT_FINAL_ASSEMBLY));
1399566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A_baij, MAT_FINAL_ASSEMBLY));
1409566063dSJacob Faibussowitsch     PetscCall(MatInvertBlockDiagonal(A_baij, &A_inv));
141418fb43bSPierre Jolivet     PetscCall(PetscArraycpy(gauss_A_inv, A_inv, nstages * nstages));
1429566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A_baij));
143d2567f34SHong Zhang   }
144d2567f34SHong Zhang 
145d2567f34SHong Zhang   /* Compute row sums A_inv_rowsum and identity I_s */
146d2567f34SHong Zhang   for (i = 0; i < nstages; i++) {
147d2567f34SHong Zhang     gauss_A_inv_rowsum[i] = 0;
148d2567f34SHong Zhang     for (j = 0; j < nstages; j++) {
149d2567f34SHong Zhang       gauss_A_inv_rowsum[i] += gauss_A_inv[i + nstages * j];
150d2567f34SHong Zhang       I_s[i + nstages * j] = 1. * (i == j);
151d2567f34SHong Zhang     }
152d2567f34SHong Zhang   }
1539566063dSJacob Faibussowitsch   PetscCall(TSIRKTableauCreate(ts, nstages, gauss_A_real, gauss_b, gauss_c, NULL, gauss_A_inv, gauss_A_inv_rowsum, I_s));
1549566063dSJacob Faibussowitsch   PetscCall(PetscFree3(gauss_A_real, gauss_b, gauss_c));
1559566063dSJacob Faibussowitsch   PetscCall(PetscFree4(gauss_A, gauss_A_inv, gauss_A_inv_rowsum, I_s));
1563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
157d2567f34SHong Zhang }
158d2567f34SHong Zhang 
159d2567f34SHong Zhang /*@C
160bcf0153eSBarry Smith   TSIRKRegister -  adds a `TSIRK` implementation
161d2567f34SHong Zhang 
162cc4c1da9SBarry Smith   Not Collective, No Fortran Support
163d2567f34SHong Zhang 
164d2567f34SHong Zhang   Input Parameters:
165d2567f34SHong Zhang + sname    - name of user-defined IRK scheme
166d2567f34SHong Zhang - function - function to create method context
167d2567f34SHong Zhang 
168bcf0153eSBarry Smith   Level: advanced
169bcf0153eSBarry Smith 
170bcf0153eSBarry Smith   Note:
171bcf0153eSBarry Smith   `TSIRKRegister()` may be called multiple times to add several user-defined families.
172d2567f34SHong Zhang 
173b43aa488SJacob Faibussowitsch   Example Usage:
174d2567f34SHong Zhang .vb
175d2567f34SHong Zhang    TSIRKRegister("my_scheme", MySchemeCreate);
176d2567f34SHong Zhang .ve
177d2567f34SHong Zhang 
178d2567f34SHong Zhang   Then, your scheme can be chosen with the procedural interface via
179b44f4de4SBarry Smith .vb
180b44f4de4SBarry Smith   TSIRKSetType(ts, "my_scheme")
181b44f4de4SBarry Smith .ve
182d2567f34SHong Zhang   or at runtime via the option
183b44f4de4SBarry Smith .vb
184b44f4de4SBarry Smith   -ts_irk_type my_scheme
185b44f4de4SBarry Smith .ve
186d2567f34SHong Zhang 
1871cc06b55SBarry Smith .seealso: [](ch_ts), `TSIRK`, `TSIRKRegisterAll()`
188d2567f34SHong Zhang @*/
TSIRKRegister(const char sname[],PetscErrorCode (* function)(TS))189d71ae5a4SJacob Faibussowitsch PetscErrorCode TSIRKRegister(const char sname[], PetscErrorCode (*function)(TS))
190d71ae5a4SJacob Faibussowitsch {
191d2567f34SHong Zhang   PetscFunctionBegin;
1929566063dSJacob Faibussowitsch   PetscCall(TSIRKInitializePackage());
1939566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&TSIRKList, sname, function));
1943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
195d2567f34SHong Zhang }
196d2567f34SHong Zhang 
197d2567f34SHong Zhang /*@C
198bcf0153eSBarry Smith   TSIRKRegisterAll - Registers all of the implicit Runge-Kutta methods in `TSIRK`
199d2567f34SHong Zhang 
200d2567f34SHong Zhang   Not Collective, but should be called by all processes which will need the schemes to be registered
201d2567f34SHong Zhang 
202d2567f34SHong Zhang   Level: advanced
203d2567f34SHong Zhang 
2041cc06b55SBarry Smith .seealso: [](ch_ts), `TSIRK`, `TSIRKRegisterDestroy()`
205d2567f34SHong Zhang @*/
TSIRKRegisterAll(void)206d71ae5a4SJacob Faibussowitsch PetscErrorCode TSIRKRegisterAll(void)
207d71ae5a4SJacob Faibussowitsch {
208d2567f34SHong Zhang   PetscFunctionBegin;
2093ba16761SJacob Faibussowitsch   if (TSIRKRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
210d2567f34SHong Zhang   TSIRKRegisterAllCalled = PETSC_TRUE;
211d2567f34SHong Zhang 
2129566063dSJacob Faibussowitsch   PetscCall(TSIRKRegister(TSIRKGAUSS, TSIRKCreate_Gauss));
2133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
214d2567f34SHong Zhang }
215d2567f34SHong Zhang 
216d2567f34SHong Zhang /*@C
217bcf0153eSBarry Smith   TSIRKRegisterDestroy - Frees the list of schemes that were registered by `TSIRKRegister()`.
218d2567f34SHong Zhang 
219d2567f34SHong Zhang   Not Collective
220d2567f34SHong Zhang 
221d2567f34SHong Zhang   Level: advanced
222d2567f34SHong Zhang 
2231cc06b55SBarry Smith .seealso: [](ch_ts), `TSIRK`, `TSIRKRegister()`, `TSIRKRegisterAll()`
224d2567f34SHong Zhang @*/
TSIRKRegisterDestroy(void)225d71ae5a4SJacob Faibussowitsch PetscErrorCode TSIRKRegisterDestroy(void)
226d71ae5a4SJacob Faibussowitsch {
227d2567f34SHong Zhang   PetscFunctionBegin;
228d2567f34SHong Zhang   TSIRKRegisterAllCalled = PETSC_FALSE;
2293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
230d2567f34SHong Zhang }
231d2567f34SHong Zhang 
232d2567f34SHong Zhang /*@C
233bcf0153eSBarry Smith   TSIRKInitializePackage - This function initializes everything in the `TSIRK` package. It is called
234bcf0153eSBarry Smith   from `TSInitializePackage()`.
235d2567f34SHong Zhang 
236d2567f34SHong Zhang   Level: developer
237d2567f34SHong Zhang 
2381cc06b55SBarry Smith .seealso: [](ch_ts), `TSIRK`, `PetscInitialize()`, `TSIRKFinalizePackage()`, `TSInitializePackage()`
239d2567f34SHong Zhang @*/
TSIRKInitializePackage(void)240d71ae5a4SJacob Faibussowitsch PetscErrorCode TSIRKInitializePackage(void)
241d71ae5a4SJacob Faibussowitsch {
242d2567f34SHong Zhang   PetscFunctionBegin;
2433ba16761SJacob Faibussowitsch   if (TSIRKPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
244d2567f34SHong Zhang   TSIRKPackageInitialized = PETSC_TRUE;
2459566063dSJacob Faibussowitsch   PetscCall(TSIRKRegisterAll());
2469566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(TSIRKFinalizePackage));
2473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
248d2567f34SHong Zhang }
249d2567f34SHong Zhang 
250d2567f34SHong Zhang /*@C
251bcf0153eSBarry Smith   TSIRKFinalizePackage - This function destroys everything in the `TSIRK` package. It is
252bcf0153eSBarry Smith   called from `PetscFinalize()`.
253d2567f34SHong Zhang 
254d2567f34SHong Zhang   Level: developer
255d2567f34SHong Zhang 
25642747ad1SJacob Faibussowitsch .seealso: [](ch_ts), `TSIRK`, `PetscFinalize()`, `TSInitializePackage()`
257d2567f34SHong Zhang @*/
TSIRKFinalizePackage(void)258d71ae5a4SJacob Faibussowitsch PetscErrorCode TSIRKFinalizePackage(void)
259d71ae5a4SJacob Faibussowitsch {
260d2567f34SHong Zhang   PetscFunctionBegin;
2619566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&TSIRKList));
262d2567f34SHong Zhang   TSIRKPackageInitialized = PETSC_FALSE;
2633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
264d2567f34SHong Zhang }
265d2567f34SHong Zhang 
266d2567f34SHong Zhang /*
267d2567f34SHong Zhang  This function can be called before or after ts->vec_sol has been updated.
268d2567f34SHong Zhang */
TSEvaluateStep_IRK(TS ts,PetscInt order,Vec U,PetscBool * done)269d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateStep_IRK(TS ts, PetscInt order, Vec U, PetscBool *done)
270d71ae5a4SJacob Faibussowitsch {
271d2567f34SHong Zhang   TS_IRK      *irk   = (TS_IRK *)ts->data;
272d2567f34SHong Zhang   IRKTableau   tab   = irk->tableau;
273d2567f34SHong Zhang   Vec         *YdotI = irk->YdotI;
274d2567f34SHong Zhang   PetscScalar *w     = irk->work;
275d2567f34SHong Zhang   PetscReal    h;
276d2567f34SHong Zhang   PetscInt     j;
277d2567f34SHong Zhang 
278d2567f34SHong Zhang   PetscFunctionBegin;
279d2567f34SHong Zhang   switch (irk->status) {
280d2567f34SHong Zhang   case TS_STEP_INCOMPLETE:
281d71ae5a4SJacob Faibussowitsch   case TS_STEP_PENDING:
282d71ae5a4SJacob Faibussowitsch     h = ts->time_step;
283d71ae5a4SJacob Faibussowitsch     break;
284d71ae5a4SJacob Faibussowitsch   case TS_STEP_COMPLETE:
285d71ae5a4SJacob Faibussowitsch     h = ts->ptime - ts->ptime_prev;
286d71ae5a4SJacob Faibussowitsch     break;
287d71ae5a4SJacob Faibussowitsch   default:
288d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus");
289d2567f34SHong Zhang   }
290d2567f34SHong Zhang 
2919566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol, U));
292d2567f34SHong Zhang   for (j = 0; j < irk->nstages; j++) w[j] = h * tab->b[j];
2939566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(U, irk->nstages, w, YdotI));
2943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
295d2567f34SHong Zhang }
296d2567f34SHong Zhang 
TSRollBack_IRK(TS ts)297d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRollBack_IRK(TS ts)
298d71ae5a4SJacob Faibussowitsch {
299d2567f34SHong Zhang   TS_IRK *irk = (TS_IRK *)ts->data;
300d2567f34SHong Zhang 
301d2567f34SHong Zhang   PetscFunctionBegin;
3029566063dSJacob Faibussowitsch   PetscCall(VecCopy(irk->U0, ts->vec_sol));
3033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
304d2567f34SHong Zhang }
305d2567f34SHong Zhang 
TSStep_IRK(TS ts)306d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_IRK(TS ts)
307d71ae5a4SJacob Faibussowitsch {
308d2567f34SHong Zhang   TS_IRK        *irk   = (TS_IRK *)ts->data;
309d2567f34SHong Zhang   IRKTableau     tab   = irk->tableau;
310d2567f34SHong Zhang   PetscScalar   *A_inv = tab->A_inv, *A_inv_rowsum = tab->A_inv_rowsum;
311d2567f34SHong Zhang   const PetscInt nstages = irk->nstages;
312d2567f34SHong Zhang   SNES           snes;
313d2567f34SHong Zhang   PetscInt       i, j, its, lits, bs;
314d2567f34SHong Zhang   TSAdapt        adapt;
315d2567f34SHong Zhang   PetscInt       rejections     = 0;
316d2567f34SHong Zhang   PetscBool      accept         = PETSC_TRUE;
317d2567f34SHong Zhang   PetscReal      next_time_step = ts->time_step;
318d2567f34SHong Zhang 
319d2567f34SHong Zhang   PetscFunctionBegin;
32048a46eb9SPierre Jolivet   if (!ts->steprollback) PetscCall(VecCopy(ts->vec_sol, irk->U0));
3219566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(ts->vec_sol, &bs));
32248a46eb9SPierre Jolivet   for (i = 0; i < nstages; i++) PetscCall(VecStrideScatter(ts->vec_sol, i * bs, irk->Z, INSERT_VALUES));
323d2567f34SHong Zhang 
324d2567f34SHong Zhang   irk->status = TS_STEP_INCOMPLETE;
325d2567f34SHong Zhang   while (!ts->reason && irk->status != TS_STEP_COMPLETE) {
3269566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vec_sol, irk->U));
3279566063dSJacob Faibussowitsch     PetscCall(TSGetSNES(ts, &snes));
3289566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes, NULL, irk->Z));
3299566063dSJacob Faibussowitsch     PetscCall(SNESGetIterationNumber(snes, &its));
3309566063dSJacob Faibussowitsch     PetscCall(SNESGetLinearSolveIterations(snes, &lits));
3319371c9d4SSatish Balay     ts->snes_its += its;
3329371c9d4SSatish Balay     ts->ksp_its += lits;
3339566063dSJacob Faibussowitsch     PetscCall(VecStrideGatherAll(irk->Z, irk->Y, INSERT_VALUES));
334d2567f34SHong Zhang     for (i = 0; i < nstages; i++) {
3359566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(irk->YdotI[i]));
33648a46eb9SPierre Jolivet       for (j = 0; j < nstages; j++) PetscCall(VecAXPY(irk->YdotI[i], A_inv[i + j * nstages] / ts->time_step, irk->Y[j]));
3379566063dSJacob Faibussowitsch       PetscCall(VecAXPY(irk->YdotI[i], -A_inv_rowsum[i] / ts->time_step, irk->U));
338d2567f34SHong Zhang     }
339d2567f34SHong Zhang     irk->status = TS_STEP_INCOMPLETE;
3409566063dSJacob Faibussowitsch     PetscCall(TSEvaluateStep_IRK(ts, irk->order, ts->vec_sol, NULL));
341d2567f34SHong Zhang     irk->status = TS_STEP_PENDING;
3429566063dSJacob Faibussowitsch     PetscCall(TSGetAdapt(ts, &adapt));
3439566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
344d2567f34SHong Zhang     irk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
345d2567f34SHong Zhang     if (!accept) {
3469566063dSJacob Faibussowitsch       PetscCall(TSRollBack_IRK(ts));
347d2567f34SHong Zhang       ts->time_step = next_time_step;
348d2567f34SHong Zhang       goto reject_step;
349d2567f34SHong Zhang     }
350d2567f34SHong Zhang 
351d2567f34SHong Zhang     ts->ptime += ts->time_step;
352d2567f34SHong Zhang     ts->time_step = next_time_step;
353d2567f34SHong Zhang     break;
354d2567f34SHong Zhang   reject_step:
3559371c9d4SSatish Balay     ts->reject++;
3569371c9d4SSatish Balay     accept = PETSC_FALSE;
357d2567f34SHong Zhang     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
358d2567f34SHong Zhang       ts->reason = TS_DIVERGED_STEP_REJECTED;
35963a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
360d2567f34SHong Zhang     }
361d2567f34SHong Zhang   }
3623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
363d2567f34SHong Zhang }
364d2567f34SHong Zhang 
TSInterpolate_IRK(TS ts,PetscReal itime,Vec U)365d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_IRK(TS ts, PetscReal itime, Vec U)
366d71ae5a4SJacob Faibussowitsch {
367d2567f34SHong Zhang   TS_IRK          *irk     = (TS_IRK *)ts->data;
368d2567f34SHong Zhang   PetscInt         nstages = irk->nstages, pinterp = irk->pinterp, i, j;
369d2567f34SHong Zhang   PetscReal        h;
370d2567f34SHong Zhang   PetscReal        tt, t;
371d2567f34SHong Zhang   PetscScalar     *bt;
372d2567f34SHong Zhang   const PetscReal *B = irk->tableau->binterp;
373d2567f34SHong Zhang 
374d2567f34SHong Zhang   PetscFunctionBegin;
3753c633725SBarry Smith   PetscCheck(B, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSIRK %s does not have an interpolation formula", irk->method_name);
376d2567f34SHong Zhang   switch (irk->status) {
377d2567f34SHong Zhang   case TS_STEP_INCOMPLETE:
378d2567f34SHong Zhang   case TS_STEP_PENDING:
379d2567f34SHong Zhang     h = ts->time_step;
380d2567f34SHong Zhang     t = (itime - ts->ptime) / h;
381d2567f34SHong Zhang     break;
382d2567f34SHong Zhang   case TS_STEP_COMPLETE:
383d2567f34SHong Zhang     h = ts->ptime - ts->ptime_prev;
384d2567f34SHong Zhang     t = (itime - ts->ptime) / h + 1; /* In the interval [0,1] */
385d2567f34SHong Zhang     break;
386d71ae5a4SJacob Faibussowitsch   default:
387d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus");
388d2567f34SHong Zhang   }
3899566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nstages, &bt));
390d2567f34SHong Zhang   for (i = 0; i < nstages; i++) bt[i] = 0;
391d2567f34SHong Zhang   for (j = 0, tt = t; j < pinterp; j++, tt *= t) {
392ad540459SPierre Jolivet     for (i = 0; i < nstages; i++) bt[i] += h * B[i * pinterp + j] * tt;
393d2567f34SHong Zhang   }
3949566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(U, nstages, bt, irk->YdotI));
3953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
396d2567f34SHong Zhang }
397d2567f34SHong Zhang 
TSIRKTableauReset(TS ts)398d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSIRKTableauReset(TS ts)
399d71ae5a4SJacob Faibussowitsch {
400d2567f34SHong Zhang   TS_IRK    *irk = (TS_IRK *)ts->data;
401d2567f34SHong Zhang   IRKTableau tab = irk->tableau;
402d2567f34SHong Zhang 
403d2567f34SHong Zhang   PetscFunctionBegin;
4043ba16761SJacob Faibussowitsch   if (!tab) PetscFunctionReturn(PETSC_SUCCESS);
4059566063dSJacob Faibussowitsch   PetscCall(PetscFree3(tab->A, tab->A_inv, tab->I_s));
4069566063dSJacob Faibussowitsch   PetscCall(PetscFree4(tab->b, tab->c, tab->binterp, tab->A_inv_rowsum));
4073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
408d2567f34SHong Zhang }
409d2567f34SHong Zhang 
TSReset_IRK(TS ts)410d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_IRK(TS ts)
411d71ae5a4SJacob Faibussowitsch {
412d2567f34SHong Zhang   TS_IRK *irk = (TS_IRK *)ts->data;
413d2567f34SHong Zhang 
414d2567f34SHong Zhang   PetscFunctionBegin;
4159566063dSJacob Faibussowitsch   PetscCall(TSIRKTableauReset(ts));
4161baa6e33SBarry Smith   if (irk->tableau) PetscCall(PetscFree(irk->tableau));
4171baa6e33SBarry Smith   if (irk->method_name) PetscCall(PetscFree(irk->method_name));
4181baa6e33SBarry Smith   if (irk->work) PetscCall(PetscFree(irk->work));
4199566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(irk->nstages, &irk->Y));
4209566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(irk->nstages, &irk->YdotI));
4219566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&irk->Ydot));
4229566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&irk->Z));
4239566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&irk->U));
4249566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&irk->U0));
4259566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&irk->TJ));
4263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
427d2567f34SHong Zhang }
428d2567f34SHong Zhang 
TSIRKGetVecs(TS ts,DM dm,Vec * U)429d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSIRKGetVecs(TS ts, DM dm, Vec *U)
430d71ae5a4SJacob Faibussowitsch {
431d2567f34SHong Zhang   TS_IRK *irk = (TS_IRK *)ts->data;
432d2567f34SHong Zhang 
433d2567f34SHong Zhang   PetscFunctionBegin;
434d2567f34SHong Zhang   if (U) {
435ac530a7eSPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMGetNamedGlobalVector(dm, "TSIRK_U", U));
436ac530a7eSPierre Jolivet     else *U = irk->U;
437d2567f34SHong Zhang   }
4383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
439d2567f34SHong Zhang }
440d2567f34SHong Zhang 
TSIRKRestoreVecs(TS ts,DM dm,Vec * U)441d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSIRKRestoreVecs(TS ts, DM dm, Vec *U)
442d71ae5a4SJacob Faibussowitsch {
443d2567f34SHong Zhang   PetscFunctionBegin;
444d2567f34SHong Zhang   if (U) {
44548a46eb9SPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSIRK_U", U));
446d2567f34SHong Zhang   }
4473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
448d2567f34SHong Zhang }
449d2567f34SHong Zhang 
450d2567f34SHong Zhang /*
451d2567f34SHong Zhang   This defines the nonlinear equations that is to be solved with SNES
452d2567f34SHong Zhang     G[e\otimes t + C*dt, Z, Zdot] = 0
453d2567f34SHong Zhang     Zdot = (In \otimes S)*Z - (In \otimes Se) U
454d2567f34SHong Zhang   where S = 1/(dt*A)
455d2567f34SHong Zhang */
SNESTSFormFunction_IRK(SNES snes,Vec ZC,Vec FC,TS ts)456d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_IRK(SNES snes, Vec ZC, Vec FC, TS ts)
457d71ae5a4SJacob Faibussowitsch {
458d2567f34SHong Zhang   TS_IRK            *irk     = (TS_IRK *)ts->data;
459d2567f34SHong Zhang   IRKTableau         tab     = irk->tableau;
460d2567f34SHong Zhang   const PetscInt     nstages = irk->nstages;
461d2567f34SHong Zhang   const PetscReal   *c       = tab->c;
462d2567f34SHong Zhang   const PetscScalar *A_inv = tab->A_inv, *A_inv_rowsum = tab->A_inv_rowsum;
463d2567f34SHong Zhang   DM                 dm, dmsave;
464d2567f34SHong Zhang   Vec                U, *YdotI = irk->YdotI, Ydot = irk->Ydot, *Y = irk->Y;
465d2567f34SHong Zhang   PetscReal          h = ts->time_step;
466d2567f34SHong Zhang   PetscInt           i, j;
467d2567f34SHong Zhang 
468d2567f34SHong Zhang   PetscFunctionBegin;
4699566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
4709566063dSJacob Faibussowitsch   PetscCall(TSIRKGetVecs(ts, dm, &U));
4719566063dSJacob Faibussowitsch   PetscCall(VecStrideGatherAll(ZC, Y, INSERT_VALUES));
472d2567f34SHong Zhang   dmsave = ts->dm;
473d2567f34SHong Zhang   ts->dm = dm;
474d2567f34SHong Zhang   for (i = 0; i < nstages; i++) {
4759566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Ydot));
47648a46eb9SPierre Jolivet     for (j = 0; j < nstages; j++) PetscCall(VecAXPY(Ydot, A_inv[j * nstages + i] / h, Y[j]));
4779566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Ydot, -A_inv_rowsum[i] / h, U)); /* Ydot = (S \otimes In)*Z - (Se \otimes In) U */
4789566063dSJacob Faibussowitsch     PetscCall(TSComputeIFunction(ts, ts->ptime + ts->time_step * c[i], Y[i], Ydot, YdotI[i], PETSC_FALSE));
479d2567f34SHong Zhang   }
4809566063dSJacob Faibussowitsch   PetscCall(VecStrideScatterAll(YdotI, FC, INSERT_VALUES));
481d2567f34SHong Zhang   ts->dm = dmsave;
4829566063dSJacob Faibussowitsch   PetscCall(TSIRKRestoreVecs(ts, dm, &U));
4833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
484d2567f34SHong Zhang }
485d2567f34SHong Zhang 
486d2567f34SHong Zhang /*
487d2567f34SHong Zhang    For explicit ODE, the Jacobian is
488d2567f34SHong Zhang      JC = I_n \otimes S - J \otimes I_s
489d2567f34SHong Zhang    For DAE, the Jacobian is
490d2567f34SHong Zhang      JC = M_n \otimes S - J \otimes I_s
491d2567f34SHong Zhang */
SNESTSFormJacobian_IRK(SNES snes,Vec ZC,Mat JC,Mat JCpre,TS ts)492d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_IRK(SNES snes, Vec ZC, Mat JC, Mat JCpre, TS ts)
493d71ae5a4SJacob Faibussowitsch {
494d2567f34SHong Zhang   TS_IRK          *irk     = (TS_IRK *)ts->data;
495d2567f34SHong Zhang   IRKTableau       tab     = irk->tableau;
496d2567f34SHong Zhang   const PetscInt   nstages = irk->nstages;
497d2567f34SHong Zhang   const PetscReal *c       = tab->c;
498d2567f34SHong Zhang   DM               dm, dmsave;
499d2567f34SHong Zhang   Vec             *Y = irk->Y, Ydot = irk->Ydot;
500d2567f34SHong Zhang   Mat              J;
501d2567f34SHong Zhang   PetscScalar     *S;
502d2567f34SHong Zhang   PetscInt         i, j, bs;
503d2567f34SHong Zhang 
504d2567f34SHong Zhang   PetscFunctionBegin;
5059566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
506d2567f34SHong Zhang   /* irk->Ydot has already been computed in SNESTSFormFunction_IRK (SNES guarantees this) */
507d2567f34SHong Zhang   dmsave = ts->dm;
508d2567f34SHong Zhang   ts->dm = dm;
5099566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(Y[nstages - 1], &bs));
510966bd95aSPierre Jolivet   PetscCheck(ts->equation_type <= TS_EQ_ODE_EXPLICIT, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSIRK %s does not support implicit formula", irk->method_name); /* TODO: need the mass matrix for DAE  */
511966bd95aSPierre Jolivet   /* Support explicit formulas only */
5129566063dSJacob Faibussowitsch   PetscCall(VecStrideGather(ZC, (nstages - 1) * bs, Y[nstages - 1], INSERT_VALUES));
5139566063dSJacob Faibussowitsch   PetscCall(MatKAIJGetAIJ(JC, &J));
5149566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts, ts->ptime + ts->time_step * c[nstages - 1], Y[nstages - 1], Ydot, 0, J, J, PETSC_FALSE));
5159566063dSJacob Faibussowitsch   PetscCall(MatKAIJGetS(JC, NULL, NULL, &S));
516d2567f34SHong Zhang   for (i = 0; i < nstages; i++)
5179371c9d4SSatish Balay     for (j = 0; j < nstages; j++) S[i + nstages * j] = tab->A_inv[i + nstages * j] / ts->time_step;
5189566063dSJacob Faibussowitsch   PetscCall(MatKAIJRestoreS(JC, &S));
519d2567f34SHong Zhang   ts->dm = dmsave;
5203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
521d2567f34SHong Zhang }
522d2567f34SHong Zhang 
DMCoarsenHook_TSIRK(DM fine,DM coarse,PetscCtx ctx)523*2a8381b2SBarry Smith static PetscErrorCode DMCoarsenHook_TSIRK(DM fine, DM coarse, PetscCtx ctx)
524d71ae5a4SJacob Faibussowitsch {
525d2567f34SHong Zhang   PetscFunctionBegin;
5263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
527d2567f34SHong Zhang }
528d2567f34SHong Zhang 
DMRestrictHook_TSIRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,PetscCtx ctx)529*2a8381b2SBarry Smith static PetscErrorCode DMRestrictHook_TSIRK(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, PetscCtx ctx)
530d71ae5a4SJacob Faibussowitsch {
531d2567f34SHong Zhang   TS  ts = (TS)ctx;
532d2567f34SHong Zhang   Vec U, U_c;
533d2567f34SHong Zhang 
534d2567f34SHong Zhang   PetscFunctionBegin;
5359566063dSJacob Faibussowitsch   PetscCall(TSIRKGetVecs(ts, fine, &U));
5369566063dSJacob Faibussowitsch   PetscCall(TSIRKGetVecs(ts, coarse, &U_c));
5379566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct, U, U_c));
5389566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(U_c, rscale, U_c));
5399566063dSJacob Faibussowitsch   PetscCall(TSIRKRestoreVecs(ts, fine, &U));
5409566063dSJacob Faibussowitsch   PetscCall(TSIRKRestoreVecs(ts, coarse, &U_c));
5413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
542d2567f34SHong Zhang }
543d2567f34SHong Zhang 
DMSubDomainHook_TSIRK(DM dm,DM subdm,PetscCtx ctx)544*2a8381b2SBarry Smith static PetscErrorCode DMSubDomainHook_TSIRK(DM dm, DM subdm, PetscCtx ctx)
545d71ae5a4SJacob Faibussowitsch {
546d2567f34SHong Zhang   PetscFunctionBegin;
5473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
548d2567f34SHong Zhang }
549d2567f34SHong Zhang 
DMSubDomainRestrictHook_TSIRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,PetscCtx ctx)550*2a8381b2SBarry Smith static PetscErrorCode DMSubDomainRestrictHook_TSIRK(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, PetscCtx ctx)
551d71ae5a4SJacob Faibussowitsch {
552d2567f34SHong Zhang   TS  ts = (TS)ctx;
553d2567f34SHong Zhang   Vec U, U_c;
554d2567f34SHong Zhang 
555d2567f34SHong Zhang   PetscFunctionBegin;
5569566063dSJacob Faibussowitsch   PetscCall(TSIRKGetVecs(ts, dm, &U));
5579566063dSJacob Faibussowitsch   PetscCall(TSIRKGetVecs(ts, subdm, &U_c));
558d2567f34SHong Zhang 
5599566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(gscat, U, U_c, INSERT_VALUES, SCATTER_FORWARD));
5609566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(gscat, U, U_c, INSERT_VALUES, SCATTER_FORWARD));
561d2567f34SHong Zhang 
5629566063dSJacob Faibussowitsch   PetscCall(TSIRKRestoreVecs(ts, dm, &U));
5639566063dSJacob Faibussowitsch   PetscCall(TSIRKRestoreVecs(ts, subdm, &U_c));
5643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
565d2567f34SHong Zhang }
566d2567f34SHong Zhang 
TSSetUp_IRK(TS ts)567d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_IRK(TS ts)
568d71ae5a4SJacob Faibussowitsch {
569d2567f34SHong Zhang   TS_IRK        *irk = (TS_IRK *)ts->data;
570d2567f34SHong Zhang   IRKTableau     tab = irk->tableau;
571d2567f34SHong Zhang   DM             dm;
572d2567f34SHong Zhang   Mat            J;
573d2567f34SHong Zhang   Vec            R;
574d2567f34SHong Zhang   const PetscInt nstages = irk->nstages;
575d2567f34SHong Zhang   PetscInt       vsize, bs;
576d2567f34SHong Zhang 
577d2567f34SHong Zhang   PetscFunctionBegin;
57848a46eb9SPierre Jolivet   if (!irk->work) PetscCall(PetscMalloc1(irk->nstages, &irk->work));
57948a46eb9SPierre Jolivet   if (!irk->Y) PetscCall(VecDuplicateVecs(ts->vec_sol, irk->nstages, &irk->Y));
58048a46eb9SPierre Jolivet   if (!irk->YdotI) PetscCall(VecDuplicateVecs(ts->vec_sol, irk->nstages, &irk->YdotI));
58148a46eb9SPierre Jolivet   if (!irk->Ydot) PetscCall(VecDuplicate(ts->vec_sol, &irk->Ydot));
58248a46eb9SPierre Jolivet   if (!irk->U) PetscCall(VecDuplicate(ts->vec_sol, &irk->U));
58348a46eb9SPierre Jolivet   if (!irk->U0) PetscCall(VecDuplicate(ts->vec_sol, &irk->U0));
584d2567f34SHong Zhang   if (!irk->Z) {
5859566063dSJacob Faibussowitsch     PetscCall(VecCreate(PetscObjectComm((PetscObject)ts->vec_sol), &irk->Z));
5869566063dSJacob Faibussowitsch     PetscCall(VecGetSize(ts->vec_sol, &vsize));
5879566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(irk->Z, PETSC_DECIDE, vsize * irk->nstages));
5889566063dSJacob Faibussowitsch     PetscCall(VecGetBlockSize(ts->vec_sol, &bs));
5899566063dSJacob Faibussowitsch     PetscCall(VecSetBlockSize(irk->Z, irk->nstages * bs));
5909566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(irk->Z));
591d2567f34SHong Zhang   }
5929566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
5939566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSIRK, DMRestrictHook_TSIRK, ts));
5949566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSIRK, DMSubDomainRestrictHook_TSIRK, ts));
595d2567f34SHong Zhang 
5969566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &ts->snes));
5979566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(irk->Z, &R));
5989566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(ts->snes, R, SNESTSFormFunction, ts));
5999566063dSJacob Faibussowitsch   PetscCall(TSGetIJacobian(ts, &J, NULL, NULL, NULL));
600d2567f34SHong Zhang   if (!irk->TJ) {
601d2567f34SHong Zhang     /* Create the KAIJ matrix for solving the stages */
6029566063dSJacob Faibussowitsch     PetscCall(MatCreateKAIJ(J, nstages, nstages, tab->A_inv, tab->I_s, &irk->TJ));
603d2567f34SHong Zhang   }
6049566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(ts->snes, irk->TJ, irk->TJ, SNESTSFormJacobian, ts));
6059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&R));
6063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
607d2567f34SHong Zhang }
608d2567f34SHong Zhang 
TSSetFromOptions_IRK(TS ts,PetscOptionItems PetscOptionsObject)609ce78bad3SBarry Smith static PetscErrorCode TSSetFromOptions_IRK(TS ts, PetscOptionItems PetscOptionsObject)
610d71ae5a4SJacob Faibussowitsch {
611d2567f34SHong Zhang   TS_IRK *irk        = (TS_IRK *)ts->data;
612d2567f34SHong Zhang   char    tname[256] = TSIRKGAUSS;
613d2567f34SHong Zhang 
614d2567f34SHong Zhang   PetscFunctionBegin;
615d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "IRK ODE solver options");
616d2567f34SHong Zhang   {
617d2567f34SHong Zhang     PetscBool flg1, flg2;
6189566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-ts_irk_nstages", "Stages of the IRK method", "TSIRKSetNumStages", irk->nstages, &irk->nstages, &flg1));
6199566063dSJacob Faibussowitsch     PetscCall(PetscOptionsFList("-ts_irk_type", "Type of IRK method", "TSIRKSetType", TSIRKList, irk->method_name[0] ? irk->method_name : tname, tname, sizeof(tname), &flg2));
620d2567f34SHong Zhang     if (flg1 || flg2 || !irk->method_name[0]) { /* Create the method tableau after nstages or method is set */
6219566063dSJacob Faibussowitsch       PetscCall(TSIRKSetType(ts, tname));
622d2567f34SHong Zhang     }
623d2567f34SHong Zhang   }
624d0609cedSBarry Smith   PetscOptionsHeadEnd();
6253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
626d2567f34SHong Zhang }
627d2567f34SHong Zhang 
TSView_IRK(TS ts,PetscViewer viewer)628d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_IRK(TS ts, PetscViewer viewer)
629d71ae5a4SJacob Faibussowitsch {
630d2567f34SHong Zhang   TS_IRK   *irk = (TS_IRK *)ts->data;
6319f196a02SMartin Diehl   PetscBool isascii;
632d2567f34SHong Zhang 
633d2567f34SHong Zhang   PetscFunctionBegin;
6349f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
6359f196a02SMartin Diehl   if (isascii) {
636d2567f34SHong Zhang     IRKTableau tab = irk->tableau;
637d2567f34SHong Zhang     TSIRKType  irktype;
638d2567f34SHong Zhang     char       buf[512];
639d2567f34SHong Zhang 
6409566063dSJacob Faibussowitsch     PetscCall(TSIRKGetType(ts, &irktype));
6419566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  IRK type %s\n", irktype));
6429566063dSJacob Faibussowitsch     PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", irk->nstages, tab->c));
6439566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Abscissa       c = %s\n", buf));
6449566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Stiffly accurate: %s\n", irk->stiffly_accurate ? "yes" : "no"));
6459566063dSJacob Faibussowitsch     PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", PetscSqr(irk->nstages), tab->A));
6469566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  A coefficients       A = %s\n", buf));
647d2567f34SHong Zhang   }
6483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
649d2567f34SHong Zhang }
650d2567f34SHong Zhang 
TSLoad_IRK(TS ts,PetscViewer viewer)651d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSLoad_IRK(TS ts, PetscViewer viewer)
652d71ae5a4SJacob Faibussowitsch {
653d2567f34SHong Zhang   SNES    snes;
654d2567f34SHong Zhang   TSAdapt adapt;
655d2567f34SHong Zhang 
656d2567f34SHong Zhang   PetscFunctionBegin;
6579566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &adapt));
6589566063dSJacob Faibussowitsch   PetscCall(TSAdaptLoad(adapt, viewer));
6599566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
6609566063dSJacob Faibussowitsch   PetscCall(SNESLoad(snes, viewer));
661d2567f34SHong Zhang   /* function and Jacobian context for SNES when used with TS is always ts object */
6629566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, NULL, NULL, ts));
6639566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, ts));
6643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
665d2567f34SHong Zhang }
666d2567f34SHong Zhang 
667cc4c1da9SBarry Smith /*@
668bcf0153eSBarry Smith   TSIRKSetType - Set the type of `TSIRK` scheme to use
669d2567f34SHong Zhang 
67020f4b53cSBarry Smith   Logically Collective
671d2567f34SHong Zhang 
67202024617SSatish Balay   Input Parameters:
673d2567f34SHong Zhang + ts      - timestepping context
674bcf0153eSBarry Smith - irktype - type of `TSIRK` scheme
675d2567f34SHong Zhang 
676bcf0153eSBarry Smith   Options Database Key:
67767b8a455SSatish Balay . -ts_irk_type <gauss> - set irk type
678d2567f34SHong Zhang 
679d2567f34SHong Zhang   Level: intermediate
680d2567f34SHong Zhang 
6811cc06b55SBarry Smith .seealso: [](ch_ts), `TSIRKGetType()`, `TSIRK`, `TSIRKType`, `TSIRKGAUSS`
682d2567f34SHong Zhang @*/
TSIRKSetType(TS ts,TSIRKType irktype)683d71ae5a4SJacob Faibussowitsch PetscErrorCode TSIRKSetType(TS ts, TSIRKType irktype)
684d71ae5a4SJacob Faibussowitsch {
685d2567f34SHong Zhang   PetscFunctionBegin;
686d2567f34SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
6874f572ea9SToby Isaac   PetscAssertPointer(irktype, 2);
688cac4c232SBarry Smith   PetscTryMethod(ts, "TSIRKSetType_C", (TS, TSIRKType), (ts, irktype));
6893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
690d2567f34SHong Zhang }
691d2567f34SHong Zhang 
692cc4c1da9SBarry Smith /*@
693bcf0153eSBarry Smith   TSIRKGetType - Get the type of `TSIRK` IMEX scheme being used
694d2567f34SHong Zhang 
69520f4b53cSBarry Smith   Logically Collective
696d2567f34SHong Zhang 
697d2567f34SHong Zhang   Input Parameter:
698d2567f34SHong Zhang . ts - timestepping context
699d2567f34SHong Zhang 
700d2567f34SHong Zhang   Output Parameter:
701bcf0153eSBarry Smith . irktype - type of `TSIRK` IMEX scheme
702d2567f34SHong Zhang 
703d2567f34SHong Zhang   Level: intermediate
704d2567f34SHong Zhang 
70542747ad1SJacob Faibussowitsch .seealso: [](ch_ts), `TSIRK`, `TSIRKType`, `TSIRKGAUSS`
706d2567f34SHong Zhang @*/
TSIRKGetType(TS ts,TSIRKType * irktype)707d71ae5a4SJacob Faibussowitsch PetscErrorCode TSIRKGetType(TS ts, TSIRKType *irktype)
708d71ae5a4SJacob Faibussowitsch {
709d2567f34SHong Zhang   PetscFunctionBegin;
710d2567f34SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
711cac4c232SBarry Smith   PetscUseMethod(ts, "TSIRKGetType_C", (TS, TSIRKType *), (ts, irktype));
7123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
713d2567f34SHong Zhang }
714d2567f34SHong Zhang 
715cc4c1da9SBarry Smith /*@
716bcf0153eSBarry Smith   TSIRKSetNumStages - Set the number of stages of `TSIRK` scheme to use
717d2567f34SHong Zhang 
71820f4b53cSBarry Smith   Logically Collective
719d2567f34SHong Zhang 
72002024617SSatish Balay   Input Parameters:
721d2567f34SHong Zhang + ts      - timestepping context
722bcf0153eSBarry Smith - nstages - number of stages of `TSIRK` scheme
723d2567f34SHong Zhang 
724bcf0153eSBarry Smith   Options Database Key:
72567b8a455SSatish Balay . -ts_irk_nstages <int> - set number of stages
726d2567f34SHong Zhang 
727d2567f34SHong Zhang   Level: intermediate
728d2567f34SHong Zhang 
7291cc06b55SBarry Smith .seealso: [](ch_ts), `TSIRKGetNumStages()`, `TSIRK`
730d2567f34SHong Zhang @*/
TSIRKSetNumStages(TS ts,PetscInt nstages)731d71ae5a4SJacob Faibussowitsch PetscErrorCode TSIRKSetNumStages(TS ts, PetscInt nstages)
732d71ae5a4SJacob Faibussowitsch {
733d2567f34SHong Zhang   PetscFunctionBegin;
734d2567f34SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
735cac4c232SBarry Smith   PetscTryMethod(ts, "TSIRKSetNumStages_C", (TS, PetscInt), (ts, nstages));
7363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
737d2567f34SHong Zhang }
738d2567f34SHong Zhang 
739cc4c1da9SBarry Smith /*@
740bcf0153eSBarry Smith   TSIRKGetNumStages - Get the number of stages of `TSIRK` scheme
741d2567f34SHong Zhang 
74220f4b53cSBarry Smith   Logically Collective
743d2567f34SHong Zhang 
74402024617SSatish Balay   Input Parameters:
745d2567f34SHong Zhang + ts      - timestepping context
746bcf0153eSBarry Smith - nstages - number of stages of `TSIRK` scheme
747d2567f34SHong Zhang 
748d2567f34SHong Zhang   Level: intermediate
749d2567f34SHong Zhang 
7501cc06b55SBarry Smith .seealso: [](ch_ts), `TSIRKSetNumStages()`, `TSIRK`
751d2567f34SHong Zhang @*/
TSIRKGetNumStages(TS ts,PetscInt * nstages)752d71ae5a4SJacob Faibussowitsch PetscErrorCode TSIRKGetNumStages(TS ts, PetscInt *nstages)
753d71ae5a4SJacob Faibussowitsch {
754d2567f34SHong Zhang   PetscFunctionBegin;
755d2567f34SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
7564f572ea9SToby Isaac   PetscAssertPointer(nstages, 2);
757cac4c232SBarry Smith   PetscTryMethod(ts, "TSIRKGetNumStages_C", (TS, PetscInt *), (ts, nstages));
7583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
759d2567f34SHong Zhang }
760d2567f34SHong Zhang 
TSIRKGetType_IRK(TS ts,TSIRKType * irktype)761d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSIRKGetType_IRK(TS ts, TSIRKType *irktype)
762d71ae5a4SJacob Faibussowitsch {
763d2567f34SHong Zhang   TS_IRK *irk = (TS_IRK *)ts->data;
764d2567f34SHong Zhang 
765d2567f34SHong Zhang   PetscFunctionBegin;
766d2567f34SHong Zhang   *irktype = irk->method_name;
7673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
768d2567f34SHong Zhang }
769d2567f34SHong Zhang 
TSIRKSetType_IRK(TS ts,TSIRKType irktype)770d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSIRKSetType_IRK(TS ts, TSIRKType irktype)
771d71ae5a4SJacob Faibussowitsch {
772d2567f34SHong Zhang   TS_IRK *irk = (TS_IRK *)ts->data;
7735f80ce2aSJacob Faibussowitsch   PetscErrorCode (*irkcreate)(TS);
774d2567f34SHong Zhang 
775d2567f34SHong Zhang   PetscFunctionBegin;
776d2567f34SHong Zhang   if (irk->method_name) {
7779566063dSJacob Faibussowitsch     PetscCall(PetscFree(irk->method_name));
7789566063dSJacob Faibussowitsch     PetscCall(TSIRKTableauReset(ts));
779d2567f34SHong Zhang   }
7809566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(TSIRKList, irktype, &irkcreate));
7816adde796SStefano Zampini   PetscCheck(irkcreate, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TSIRK type \"%s\" given", irktype);
7829566063dSJacob Faibussowitsch   PetscCall((*irkcreate)(ts));
7839566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(irktype, &irk->method_name));
7843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
785d2567f34SHong Zhang }
786d2567f34SHong Zhang 
TSIRKSetNumStages_IRK(TS ts,PetscInt nstages)787d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSIRKSetNumStages_IRK(TS ts, PetscInt nstages)
788d71ae5a4SJacob Faibussowitsch {
789d2567f34SHong Zhang   TS_IRK *irk = (TS_IRK *)ts->data;
790d2567f34SHong Zhang 
791d2567f34SHong Zhang   PetscFunctionBegin;
79263a3b9bcSJacob Faibussowitsch   PetscCheck(nstages > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "input argument, %" PetscInt_FMT ", out of range", nstages);
793d2567f34SHong Zhang   irk->nstages = nstages;
7943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
795d2567f34SHong Zhang }
796d2567f34SHong Zhang 
TSIRKGetNumStages_IRK(TS ts,PetscInt * nstages)797d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSIRKGetNumStages_IRK(TS ts, PetscInt *nstages)
798d71ae5a4SJacob Faibussowitsch {
799d2567f34SHong Zhang   TS_IRK *irk = (TS_IRK *)ts->data;
800d2567f34SHong Zhang 
801d2567f34SHong Zhang   PetscFunctionBegin;
8024f572ea9SToby Isaac   PetscAssertPointer(nstages, 2);
803d2567f34SHong Zhang   *nstages = irk->nstages;
8043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
805d2567f34SHong Zhang }
806d2567f34SHong Zhang 
TSDestroy_IRK(TS ts)807d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_IRK(TS ts)
808d71ae5a4SJacob Faibussowitsch {
809d2567f34SHong Zhang   PetscFunctionBegin;
8109566063dSJacob Faibussowitsch   PetscCall(TSReset_IRK(ts));
811d2567f34SHong Zhang   if (ts->dm) {
8129566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSIRK, DMRestrictHook_TSIRK, ts));
8139566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSIRK, DMSubDomainRestrictHook_TSIRK, ts));
814d2567f34SHong Zhang   }
8159566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
8169566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSIRKSetType_C", NULL));
8179566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSIRKGetType_C", NULL));
8189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSIRKSetNumStages_C", NULL));
8199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSIRKGetNumStages_C", NULL));
8203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
821d2567f34SHong Zhang }
822d2567f34SHong Zhang 
823d2567f34SHong Zhang /*MC
824d2567f34SHong Zhang   TSIRK - ODE and DAE solver using Implicit Runge-Kutta schemes
825d2567f34SHong Zhang 
826d2567f34SHong Zhang   Level: beginner
827d2567f34SHong Zhang 
828bcf0153eSBarry Smith   Notes:
829bcf0153eSBarry Smith   `TSIRK` uses the sparse Kronecker product matrix implementation of `MATKAIJ` to achieve good arithmetic intensity.
830d2567f34SHong Zhang 
831bcf0153eSBarry Smith   Gauss-Legrendre methods are currently supported. These are A-stable symplectic methods with an arbitrary number of stages. The order of accuracy is 2s
832bcf0153eSBarry Smith   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
833efa39862SBarry Smith   `-ts_irk_nstages` or `TSIRKSetNumStages()`.
834bcf0153eSBarry Smith 
8351cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSIRKSetType()`, `TSIRKGetType()`, `TSIRKGAUSS`, `TSIRKRegister()`, `TSIRKSetNumStages()`, `TSType`
836d2567f34SHong Zhang M*/
TSCreate_IRK(TS ts)837d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_IRK(TS ts)
838d71ae5a4SJacob Faibussowitsch {
839d2567f34SHong Zhang   TS_IRK *irk;
840d2567f34SHong Zhang 
841d2567f34SHong Zhang   PetscFunctionBegin;
8429566063dSJacob Faibussowitsch   PetscCall(TSIRKInitializePackage());
843d2567f34SHong Zhang 
844d2567f34SHong Zhang   ts->ops->reset          = TSReset_IRK;
845d2567f34SHong Zhang   ts->ops->destroy        = TSDestroy_IRK;
846d2567f34SHong Zhang   ts->ops->view           = TSView_IRK;
847d2567f34SHong Zhang   ts->ops->load           = TSLoad_IRK;
848d2567f34SHong Zhang   ts->ops->setup          = TSSetUp_IRK;
849d2567f34SHong Zhang   ts->ops->step           = TSStep_IRK;
850d2567f34SHong Zhang   ts->ops->interpolate    = TSInterpolate_IRK;
851d2567f34SHong Zhang   ts->ops->evaluatestep   = TSEvaluateStep_IRK;
852d2567f34SHong Zhang   ts->ops->rollback       = TSRollBack_IRK;
853d2567f34SHong Zhang   ts->ops->setfromoptions = TSSetFromOptions_IRK;
854d2567f34SHong Zhang   ts->ops->snesfunction   = SNESTSFormFunction_IRK;
855d2567f34SHong Zhang   ts->ops->snesjacobian   = SNESTSFormJacobian_IRK;
856d2567f34SHong Zhang 
857d2567f34SHong Zhang   ts->usessnes = PETSC_TRUE;
858d2567f34SHong Zhang 
8594dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&irk));
860d2567f34SHong Zhang   ts->data = (void *)irk;
861d2567f34SHong Zhang 
8629566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSIRKSetType_C", TSIRKSetType_IRK));
8639566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSIRKGetType_C", TSIRKGetType_IRK));
8649566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSIRKSetNumStages_C", TSIRKSetNumStages_IRK));
8659566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSIRKGetNumStages_C", TSIRKGetNumStages_IRK));
866d2567f34SHong Zhang   /* 3-stage IRK_Gauss is the default */
8679566063dSJacob Faibussowitsch   PetscCall(PetscNew(&irk->tableau));
868d2567f34SHong Zhang   irk->nstages = 3;
8699566063dSJacob Faibussowitsch   PetscCall(TSIRKSetType(ts, TSIRKDefault));
8703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
871d2567f34SHong Zhang }
872