14b84eec9SHong Zhang /*
2f944a40eSHong Zhang Code for time stepping with the Multirate Partitioned Runge-Kutta method
34b84eec9SHong Zhang
44b84eec9SHong Zhang Notes:
54b84eec9SHong Zhang 1) The general system is written as
6f944a40eSHong Zhang Udot = F(t,U)
7f944a40eSHong Zhang if one does not split the RHS function, but gives the indexes for both slow and fast components;
84b84eec9SHong Zhang 2) The general system is written as
94b84eec9SHong Zhang Usdot = Fs(t,Us,Uf)
10f944a40eSHong Zhang Ufdot = Ff(t,Us,Uf)
11f944a40eSHong Zhang for component-wise partitioned system,
12f944a40eSHong Zhang users should split the RHS function themselves and also provide the indexes for both slow and fast components.
1319c2959aSHong Zhang 3) To correct The confusing terminology in the paper, we use 'slow method', 'slow buffer method' and 'fast method' to denote the methods applied to 'slow region', 'slow buffer region' and 'fast region' respectively. The 'slow method' in the original paper actually means the 'slow buffer method'.
1419c2959aSHong Zhang 4) Why does the buffer region have to be inside the slow region? The buffer region is treated with a slow method essentially. Applying the slow method to a region with a fast characteristic time scale is apparently not a good choice.
1519c2959aSHong Zhang
1619c2959aSHong Zhang Reference:
1719c2959aSHong Zhang Emil M. Constantinescu, Adrian Sandu, Multirate Timestepping Methods for Hyperbolic Conservation Laws, Journal of Scientific Computing 2007
184b84eec9SHong Zhang */
194b84eec9SHong Zhang
204b84eec9SHong Zhang #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
214b84eec9SHong Zhang #include <petscdm.h>
224b84eec9SHong Zhang
2319c2959aSHong Zhang static TSMPRKType TSMPRKDefault = TSMPRK2A22;
244b84eec9SHong Zhang static PetscBool TSMPRKRegisterAllCalled;
254b84eec9SHong Zhang static PetscBool TSMPRKPackageInitialized;
264b84eec9SHong Zhang
274b84eec9SHong Zhang typedef struct _MPRKTableau *MPRKTableau;
284b84eec9SHong Zhang struct _MPRKTableau {
294b84eec9SHong Zhang char *name;
304b84eec9SHong Zhang PetscInt order; /* Classical approximation order of the method i */
3179bc8a77SHong Zhang PetscInt sbase; /* Number of stages in the base method*/
324b84eec9SHong Zhang PetscInt s; /* Number of stages */
3319c2959aSHong Zhang PetscInt np; /* Number of partitions */
344b84eec9SHong Zhang PetscReal *Af, *bf, *cf; /* Tableau for fast components */
3519c2959aSHong Zhang PetscReal *Amb, *bmb, *cmb; /* Tableau for medium components */
3619c2959aSHong Zhang PetscInt *rmb; /* Array of flags for repeated stages in medium method */
3719c2959aSHong Zhang PetscReal *Asb, *bsb, *csb; /* Tableau for slow components */
3819c2959aSHong Zhang PetscInt *rsb; /* Array of flags for repeated staged in slow method*/
394b84eec9SHong Zhang };
404b84eec9SHong Zhang typedef struct _MPRKTableauLink *MPRKTableauLink;
414b84eec9SHong Zhang struct _MPRKTableauLink {
424b84eec9SHong Zhang struct _MPRKTableau tab;
434b84eec9SHong Zhang MPRKTableauLink next;
444b84eec9SHong Zhang };
454b84eec9SHong Zhang static MPRKTableauLink MPRKTableauList;
464b84eec9SHong Zhang
474b84eec9SHong Zhang typedef struct {
484b84eec9SHong Zhang MPRKTableau tableau;
494b84eec9SHong Zhang Vec *Y; /* States computed during the step */
507c0df07dSHong Zhang Vec *YdotRHS;
514b84eec9SHong Zhang Vec *YdotRHS_slow; /* Function evaluations by slow tableau for slow components */
5214d0ab18SJacob Faibussowitsch Vec *YdotRHS_slowbuffer; /* Function evaluations by medium tableau for slow components */
5314d0ab18SJacob Faibussowitsch Vec *YdotRHS_medium; /* Function evaluations by medium tableau for medium components*/
5414d0ab18SJacob Faibussowitsch Vec *YdotRHS_mediumbuffer; /* Function evaluations by fast tableau for medium components */
5519c2959aSHong Zhang Vec *YdotRHS_fast; /* Function evaluations by fast tableau for fast components */
564b84eec9SHong Zhang PetscScalar *work_slow; /* Scalar work_slow by slow tableau */
5719c2959aSHong Zhang PetscScalar *work_slowbuffer; /* Scalar work_slow by slow tableau */
5819c2959aSHong Zhang PetscScalar *work_medium; /* Scalar work_slow by medium tableau */
5914d0ab18SJacob Faibussowitsch PetscScalar *work_mediumbuffer; /* Scalar work_mediumbuffer by fast tableau */
6019c2959aSHong Zhang PetscScalar *work_fast; /* Scalar work_fast by fast tableau */
614b84eec9SHong Zhang PetscReal stage_time;
624b84eec9SHong Zhang TSStepStatus status;
634b84eec9SHong Zhang PetscReal ptime;
644b84eec9SHong Zhang PetscReal time_step;
6519c2959aSHong Zhang IS is_slow, is_slowbuffer, is_medium, is_mediumbuffer, is_fast;
6619c2959aSHong Zhang TS subts_slow, subts_slowbuffer, subts_medium, subts_mediumbuffer, subts_fast;
674b84eec9SHong Zhang } TS_MPRK;
684b84eec9SHong Zhang
TSMPRKGenerateTableau2(PetscInt ratio,PetscInt s,const PetscReal Abase[],const PetscReal bbase[],PetscReal A1[],PetscReal b1[],PetscReal A2[],PetscReal b2[])69d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSMPRKGenerateTableau2(PetscInt ratio, PetscInt s, const PetscReal Abase[], const PetscReal bbase[], PetscReal A1[], PetscReal b1[], PetscReal A2[], PetscReal b2[])
70d71ae5a4SJacob Faibussowitsch {
7119c2959aSHong Zhang PetscInt i, j, k, l;
7219c2959aSHong Zhang
7319c2959aSHong Zhang PetscFunctionBegin;
7419c2959aSHong Zhang for (k = 0; k < ratio; k++) {
7519c2959aSHong Zhang /* diagonal blocks */
7619c2959aSHong Zhang for (i = 0; i < s; i++)
7719c2959aSHong Zhang for (j = 0; j < s; j++) {
7819c2959aSHong Zhang A1[(k * s + i) * ratio * s + k * s + j] = Abase[i * s + j];
7919c2959aSHong Zhang A2[(k * s + i) * ratio * s + k * s + j] = Abase[i * s + j] / ratio;
8019c2959aSHong Zhang }
814cf0e950SBarry Smith /* off-diagonal blocks */
8219c2959aSHong Zhang for (l = 0; l < k; l++)
8319c2959aSHong Zhang for (i = 0; i < s; i++)
849371c9d4SSatish Balay for (j = 0; j < s; j++) A2[(k * s + i) * ratio * s + l * s + j] = bbase[j] / ratio;
8519c2959aSHong Zhang for (j = 0; j < s; j++) {
8619c2959aSHong Zhang b1[k * s + j] = bbase[j] / ratio;
8719c2959aSHong Zhang b2[k * s + j] = bbase[j] / ratio;
8819c2959aSHong Zhang }
8919c2959aSHong Zhang }
903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9119c2959aSHong Zhang }
9219c2959aSHong Zhang
TSMPRKGenerateTableau3(PetscInt ratio,PetscInt s,const PetscReal Abase[],const PetscReal bbase[],PetscReal A1[],PetscReal b1[],PetscReal A2[],PetscReal b2[],PetscReal A3[],PetscReal b3[])93d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSMPRKGenerateTableau3(PetscInt ratio, PetscInt s, const PetscReal Abase[], const PetscReal bbase[], PetscReal A1[], PetscReal b1[], PetscReal A2[], PetscReal b2[], PetscReal A3[], PetscReal b3[])
94d71ae5a4SJacob Faibussowitsch {
9519c2959aSHong Zhang PetscInt i, j, k, l, m, n;
9619c2959aSHong Zhang
9719c2959aSHong Zhang PetscFunctionBegin;
9819c2959aSHong Zhang for (k = 0; k < ratio; k++) { /* diagonal blocks of size ratio*s by ratio*s */
9919c2959aSHong Zhang for (l = 0; l < ratio; l++) /* diagonal sub-blocks of size s by s */
10019c2959aSHong Zhang for (i = 0; i < s; i++)
10119c2959aSHong Zhang for (j = 0; j < s; j++) {
10219c2959aSHong Zhang A1[((k * ratio + l) * s + i) * ratio * ratio * s + (k * ratio + l) * s + j] = Abase[i * s + j];
10319c2959aSHong Zhang A2[((k * ratio + l) * s + i) * ratio * ratio * s + (k * ratio + l) * s + j] = Abase[i * s + j] / ratio;
10419c2959aSHong Zhang A3[((k * ratio + l) * s + i) * ratio * ratio * s + (k * ratio + l) * s + j] = Abase[i * s + j] / ratio / ratio;
10519c2959aSHong Zhang }
10619c2959aSHong Zhang for (l = 0; l < k; l++) /* off-diagonal blocks of size ratio*s by ratio*s */
10719c2959aSHong Zhang for (m = 0; m < ratio; m++)
10819c2959aSHong Zhang for (n = 0; n < ratio; n++)
10919c2959aSHong Zhang for (i = 0; i < s; i++)
11019c2959aSHong Zhang for (j = 0; j < s; j++) {
11119c2959aSHong Zhang A2[((k * ratio + m) * s + i) * ratio * ratio * s + (l * ratio + n) * s + j] = bbase[j] / ratio / ratio;
11219c2959aSHong Zhang A3[((k * ratio + m) * s + i) * ratio * ratio * s + (l * ratio + n) * s + j] = bbase[j] / ratio / ratio;
11319c2959aSHong Zhang }
11419c2959aSHong Zhang for (m = 0; m < ratio; m++)
11519c2959aSHong Zhang for (n = 0; n < m; n++) /* off-diagonal sub-blocks of size s by s in the diagonal blocks */
11619c2959aSHong Zhang for (i = 0; i < s; i++)
1179371c9d4SSatish Balay for (j = 0; j < s; j++) A3[((k * ratio + m) * s + i) * ratio * ratio * s + (k * ratio + n) * s + j] = bbase[j] / ratio / ratio;
11819c2959aSHong Zhang for (n = 0; n < ratio; n++)
11919c2959aSHong Zhang for (j = 0; j < s; j++) {
12019c2959aSHong Zhang b1[(k * ratio + n) * s + j] = bbase[j] / ratio / ratio;
12119c2959aSHong Zhang b2[(k * ratio + n) * s + j] = bbase[j] / ratio / ratio;
12219c2959aSHong Zhang b3[(k * ratio + n) * s + j] = bbase[j] / ratio / ratio;
12319c2959aSHong Zhang }
12419c2959aSHong Zhang }
1253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
12619c2959aSHong Zhang }
12719c2959aSHong Zhang
1284b84eec9SHong Zhang /*MC
129f944a40eSHong Zhang TSMPRK2A22 - Second Order Multirate Partitioned Runge Kutta scheme based on RK2A.
1304b84eec9SHong Zhang
13119c2959aSHong Zhang This method has four stages for slow and fast parts. The refinement factor of the stepsize is 2.
13219c2959aSHong Zhang r = 2, np = 2
133147403d9SBarry Smith
134bcf0153eSBarry Smith Options Database Key:
135147403d9SBarry Smith . -ts_mprk_type 2a22 - select this scheme
1364b84eec9SHong Zhang
1374b84eec9SHong Zhang Level: advanced
1384b84eec9SHong Zhang
1391cc06b55SBarry Smith .seealso: [](ch_ts), `TSMPRK`, `TSMPRKType`, `TSMPRKSetType()`
1404b84eec9SHong Zhang M*/
1414b84eec9SHong Zhang /*MC
142f944a40eSHong Zhang TSMPRK2A23 - Second Order Multirate Partitioned Runge-Kutta scheme based on RK2A.
14319c2959aSHong Zhang
14419c2959aSHong Zhang This method has eight stages for slow and medium and fast parts. The refinement factor of the stepsize is 2.
14519c2959aSHong Zhang r = 2, np = 3
146147403d9SBarry Smith
147bcf0153eSBarry Smith Options Database Key:
148147403d9SBarry Smith . -ts_mprk_type 2a23 - select this scheme
14919c2959aSHong Zhang
15019c2959aSHong Zhang Level: advanced
15119c2959aSHong Zhang
1521cc06b55SBarry Smith .seealso: [](ch_ts), `TSMPRK`, `TSMPRKType`, `TSMPRKSetType()`
15319c2959aSHong Zhang M*/
15419c2959aSHong Zhang /*MC
155f944a40eSHong Zhang TSMPRK2A32 - Second Order Multirate Partitioned Runge-Kutta scheme based on RK2A.
15619c2959aSHong Zhang
15719c2959aSHong Zhang This method has four stages for slow and fast parts. The refinement factor of the stepsize is 3.
15819c2959aSHong Zhang r = 3, np = 2
159147403d9SBarry Smith
160bcf0153eSBarry Smith Options Database Key:
161147403d9SBarry Smith . -ts_mprk_type 2a32 - select this scheme
16219c2959aSHong Zhang
16319c2959aSHong Zhang Level: advanced
16419c2959aSHong Zhang
1651cc06b55SBarry Smith .seealso: [](ch_ts), `TSMPRK`, `TSMPRKType`, `TSMPRKSetType()`
16619c2959aSHong Zhang M*/
16719c2959aSHong Zhang /*MC
168f944a40eSHong Zhang TSMPRK2A33 - Second Order Multirate Partitioned Runge-Kutta scheme based on RK2A.
16919c2959aSHong Zhang
17019c2959aSHong Zhang This method has eight stages for slow and medium and fast parts. The refinement factor of the stepsize is 3.
17119c2959aSHong Zhang r = 3, np = 3
172147403d9SBarry Smith
173bcf0153eSBarry Smith Options Database Key:
174147403d9SBarry Smith . -ts_mprk_type 2a33- select this scheme
17519c2959aSHong Zhang
17619c2959aSHong Zhang Level: advanced
17719c2959aSHong Zhang
1781cc06b55SBarry Smith .seealso: [](ch_ts), `TSMPRK`, `TSMPRKType`, `TSMPRKSetType()`
17919c2959aSHong Zhang M*/
18019c2959aSHong Zhang /*MC
181f944a40eSHong Zhang TSMPRK3P2M - Third Order Multirate Partitioned Runge-Kutta scheme.
1824b84eec9SHong Zhang
1834b84eec9SHong Zhang This method has eight stages for both slow and fast parts.
1844b84eec9SHong Zhang
185bcf0153eSBarry Smith Options Database Key:
186147403d9SBarry Smith . -ts_mprk_type pm3 - select this scheme
1874b84eec9SHong Zhang
1884b84eec9SHong Zhang Level: advanced
1894b84eec9SHong Zhang
1901cc06b55SBarry Smith .seealso: [](ch_ts), `TSMPRK`, `TSMPRKType`, `TSMPRKSetType()`
1914b84eec9SHong Zhang M*/
1924b84eec9SHong Zhang /*MC
193f944a40eSHong Zhang TSMPRKP2 - Second Order Multirate Partitioned Runge-Kutta scheme.
1944b84eec9SHong Zhang
1954b84eec9SHong Zhang This method has five stages for both slow and fast parts.
1964b84eec9SHong Zhang
197bcf0153eSBarry Smith Options Database Key:
198147403d9SBarry Smith . -ts_mprk_type p2 - select this scheme
1994b84eec9SHong Zhang
2004b84eec9SHong Zhang Level: advanced
2014b84eec9SHong Zhang
2021cc06b55SBarry Smith .seealso: [](ch_ts), `TSMPRK`, `TSMPRKType`, `TSMPRKSetType()`
2034b84eec9SHong Zhang M*/
2044b84eec9SHong Zhang /*MC
205f944a40eSHong Zhang TSMPRKP3 - Third Order Multirate Partitioned Runge-Kutta scheme.
2064b84eec9SHong Zhang
2074b84eec9SHong Zhang This method has ten stages for both slow and fast parts.
2084b84eec9SHong Zhang
209bcf0153eSBarry Smith Options Database Key:
210147403d9SBarry Smith . -ts_mprk_type p3 - select this scheme
2114b84eec9SHong Zhang
2124b84eec9SHong Zhang Level: advanced
2134b84eec9SHong Zhang
2141cc06b55SBarry Smith .seealso: [](ch_ts), `TSMPRK`, `TSMPRKType`, `TSMPRKSetType()`
2154b84eec9SHong Zhang M*/
2164b84eec9SHong Zhang
2174b84eec9SHong Zhang /*@C
218d5b43468SJose E. Roman TSMPRKRegisterAll - Registers all of the Partitioned Runge-Kutta explicit methods in `TSMPRK`
2194b84eec9SHong Zhang
2204b84eec9SHong Zhang Not Collective, but should be called by all processes which will need the schemes to be registered
2214b84eec9SHong Zhang
2224b84eec9SHong Zhang Level: advanced
2234b84eec9SHong Zhang
2241cc06b55SBarry Smith .seealso: [](ch_ts), `TSMPRK`, `TSMPRKRegisterDestroy()`
2254b84eec9SHong Zhang @*/
TSMPRKRegisterAll(void)226d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMPRKRegisterAll(void)
227d71ae5a4SJacob Faibussowitsch {
2284b84eec9SHong Zhang PetscFunctionBegin;
2293ba16761SJacob Faibussowitsch if (TSMPRKRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
2304b84eec9SHong Zhang TSMPRKRegisterAllCalled = PETSC_TRUE;
2314b84eec9SHong Zhang
2324b84eec9SHong Zhang #define RC PetscRealConstant
2334b84eec9SHong Zhang {
2349371c9d4SSatish Balay const PetscReal Abase[2][2] =
2359371c9d4SSatish Balay {
2369371c9d4SSatish Balay {0, 0},
2379371c9d4SSatish Balay {RC(1.0), 0}
2389371c9d4SSatish Balay },
23919c2959aSHong Zhang bbase[2] = {RC(0.5), RC(0.5)};
2409371c9d4SSatish Balay PetscReal Asb[4][4] = {{0}}, Af[4][4] = {{0}}, bsb[4] = {0}, bf[4] = {0};
2419371c9d4SSatish Balay PetscInt rsb[4] = {0, 0, 1, 2};
2429566063dSJacob Faibussowitsch PetscCall(TSMPRKGenerateTableau2(2, 2, &Abase[0][0], bbase, &Asb[0][0], bsb, &Af[0][0], bf));
2439566063dSJacob Faibussowitsch PetscCall(TSMPRKRegister(TSMPRK2A22, 2, 2, 2, 1, &Asb[0][0], bsb, NULL, rsb, NULL, NULL, NULL, NULL, &Af[0][0], bf, NULL));
2444b84eec9SHong Zhang }
24519c2959aSHong Zhang {
2469371c9d4SSatish Balay const PetscReal Abase[2][2] =
2479371c9d4SSatish Balay {
2489371c9d4SSatish Balay {0, 0},
2499371c9d4SSatish Balay {RC(1.0), 0}
2509371c9d4SSatish Balay },
25119c2959aSHong Zhang bbase[2] = {RC(0.5), RC(0.5)};
2529371c9d4SSatish Balay PetscReal Asb[8][8] = {{0}}, Amb[8][8] = {{0}}, Af[8][8] = {{0}}, bsb[8] = {0}, bmb[8] = {0}, bf[8] = {0};
2539371c9d4SSatish Balay PetscInt rsb[8] = {0, 0, 1, 2, 1, 2, 1, 2}, rmb[8] = {0, 0, 1, 2, 0, 0, 5, 6};
2549566063dSJacob Faibussowitsch PetscCall(TSMPRKGenerateTableau3(2, 2, &Abase[0][0], bbase, &Asb[0][0], bsb, &Amb[0][0], bmb, &Af[0][0], bf));
2559566063dSJacob Faibussowitsch PetscCall(TSMPRKRegister(TSMPRK2A23, 2, 2, 2, 2, &Asb[0][0], bsb, NULL, rsb, &Amb[0][0], bmb, NULL, rmb, &Af[0][0], bf, NULL));
25619c2959aSHong Zhang }
25719c2959aSHong Zhang {
2589371c9d4SSatish Balay const PetscReal Abase[2][2] =
2599371c9d4SSatish Balay {
2609371c9d4SSatish Balay {0, 0},
2619371c9d4SSatish Balay {RC(1.0), 0}
2629371c9d4SSatish Balay },
26319c2959aSHong Zhang bbase[2] = {RC(0.5), RC(0.5)};
2649371c9d4SSatish Balay PetscReal Asb[6][6] = {{0}}, Af[6][6] = {{0}}, bsb[6] = {0}, bf[6] = {0};
2659371c9d4SSatish Balay PetscInt rsb[6] = {0, 0, 1, 2, 1, 2};
2669566063dSJacob Faibussowitsch PetscCall(TSMPRKGenerateTableau2(3, 2, &Abase[0][0], bbase, &Asb[0][0], bsb, &Af[0][0], bf));
2679566063dSJacob Faibussowitsch PetscCall(TSMPRKRegister(TSMPRK2A32, 2, 2, 3, 1, &Asb[0][0], bsb, NULL, rsb, NULL, NULL, NULL, NULL, &Af[0][0], bf, NULL));
26819c2959aSHong Zhang }
26919c2959aSHong Zhang {
2709371c9d4SSatish Balay const PetscReal Abase[2][2] =
2719371c9d4SSatish Balay {
2729371c9d4SSatish Balay {0, 0},
2739371c9d4SSatish Balay {RC(1.0), 0}
2749371c9d4SSatish Balay },
27519c2959aSHong Zhang bbase[2] = {RC(0.5), RC(0.5)};
2769371c9d4SSatish Balay PetscReal Asb[18][18] = {{0}}, Amb[18][18] = {{0}}, Af[18][18] = {{0}}, bsb[18] = {0}, bmb[18] = {0}, bf[18] = {0};
2779371c9d4SSatish Balay PetscInt rsb[18] = {0, 0, 1, 2, 1, 2, 1, 2, 1, 2, 1, 2, 1, 2, 1, 2, 1, 2}, rmb[18] = {0, 0, 1, 2, 1, 2, 0, 0, 7, 8, 7, 8, 0, 0, 13, 14, 13, 14};
2789566063dSJacob Faibussowitsch PetscCall(TSMPRKGenerateTableau3(3, 2, &Abase[0][0], bbase, &Asb[0][0], bsb, &Amb[0][0], bmb, &Af[0][0], bf));
2799566063dSJacob Faibussowitsch PetscCall(TSMPRKRegister(TSMPRK2A33, 2, 2, 3, 3, &Asb[0][0], bsb, NULL, rsb, &Amb[0][0], bmb, NULL, rmb, &Af[0][0], bf, NULL));
28019c2959aSHong Zhang }
28119c2959aSHong Zhang /*
28219c2959aSHong Zhang PetscReal
28319c2959aSHong Zhang Asb[8][8] = {{Abase[0][0],Abase[0][1],0,0,0,0,0,0},
28419c2959aSHong Zhang {Abase[1][0],Abase[1][1],0,0,0,0,0,0},
28519c2959aSHong Zhang {0,0,Abase[0][0],Abase[0][1],0,0,0,0},
28619c2959aSHong Zhang {0,0,Abase[1][0],Abase[1][1],0,0,0,0},
28719c2959aSHong Zhang {0,0,0,0,Abase[0][0],Abase[0][1],0,0},
28819c2959aSHong Zhang {0,0,0,0,Abase[1][0],Abase[1][1],0,0},
28919c2959aSHong Zhang {0,0,0,0,0,0,Abase[0][0],Abase[0][1]},
29019c2959aSHong Zhang {0,0,0,0,0,0,Abase[1][0],Abase[1][1]}},
29119c2959aSHong Zhang Amb[8][8] = {{Abase[0][0]/m,Abase[0][1]/m,0,0,0,0,0,0},
29219c2959aSHong Zhang {Abase[1][0]/m,Abase[1][1]/m,0,0,0,0,0,0},
29319c2959aSHong Zhang {0,0,Abase[0][0]/m,Abase[0][1]/m,0,0,0,0},
29419c2959aSHong Zhang {0,0,Abase[1][0]/m,Abase[1][1]/m,0,0,0,0},
29519c2959aSHong Zhang {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,Abase[0][0]/m,Abase[0][1]/m,0,0},
29619c2959aSHong Zhang {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,Abase[1][0]/m,Abase[1][1]/m,0,0},
29719c2959aSHong Zhang {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,0,0,Abase[0][0]/m,Abase[0][1]/m},
29819c2959aSHong Zhang {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,0,0,Abase[1][0]/m,Abase[1][1]/m}},
29919c2959aSHong Zhang Af[8][8] = {{Abase[0][0]/m/m,Abase[0][1]/m/m,0,0,0,0,0,0},
30019c2959aSHong Zhang {Abase[1][0]/m/m,Abase[1][1]/m/m,0,0,0,0,0,0},
30119c2959aSHong Zhang {0,0,Abase[0][0]/m/m,Abase[0][1]/m/m,0,0,0,0},
30219c2959aSHong Zhang {0,0,Abase[1][0]/m/m,Abase[1][1]/m/m,0,0,0,0},
30319c2959aSHong Zhang {bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m,Abase[0][0]/m/m,Abase[0][1]/m/m,0,0},
30419c2959aSHong Zhang {bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m,Abase[1][0]/m/m,Abase[1][1]/m/m,0,0},
30519c2959aSHong Zhang {bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m,Abase[0][0]/m,Abase[0][1]/m},
30619c2959aSHong Zhang {bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m,Abase[1][0]/m,Abase[1][1]/m}},
30719c2959aSHong Zhang bsb[8] = {bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m},
30819c2959aSHong Zhang bmb[8] = {bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m,bbase[1]/m/m},
30919c2959aSHong Zhang bf[8] = {bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m,bbase[0]/m/m,bbase[1]/m/m},
31019c2959aSHong Zhang */
3114b84eec9SHong Zhang /*{
3124b84eec9SHong Zhang const PetscReal
3134b84eec9SHong Zhang As[8][8] = {{0,0,0,0,0,0,0,0},
3144b84eec9SHong Zhang {RC(1.0)/RC(2.0),0,0,0,0,0,0,0},
3154b84eec9SHong Zhang {RC(-1.0)/RC(6.0),RC(2.0)/RC(3.0),0,0,0,0,0,0},
3164b84eec9SHong Zhang {RC(1.0)/RC(3.0),RC(-1.0)/RC(3.0),RC(1.0),0,0,0,0,0},
3174b84eec9SHong Zhang {0,0,0,0,0,0,0,0},
3184b84eec9SHong Zhang {0,0,0,0,RC(1.0)/RC(2.0),0,0,0},
3194b84eec9SHong Zhang {0,0,0,0,RC(-1.0)/RC(6.0),RC(2.0)/RC(3.0),0,0},
3204b84eec9SHong Zhang {0,0,0,0,RC(1.0)/RC(3.0),RC(-1.0)/RC(3.0),RC(1.0),0}},
3214b84eec9SHong Zhang A[8][8] = {{0,0,0,0,0,0,0,0},
3224b84eec9SHong Zhang {RC(1.0)/RC(4.0),0,0,0,0,0,0,0},
3234b84eec9SHong Zhang {RC(-1.0)/RC(12.0),RC(1.0)/RC(3.0),0,0,0,0,0,0},
3244b84eec9SHong Zhang {RC(1.0)/RC(6.0),RC(-1.0)/RC(6.0),RC(1.0)/RC(2.0),0,0,0,0,0},
3254b84eec9SHong Zhang {RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0),0,0,0,0},
3264b84eec9SHong Zhang {RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0),RC(1.0)/RC(4.0),0,0,0},
3274b84eec9SHong Zhang {RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0),RC(-1.0)/RC(12.0),RC(1.0)/RC(3.0),0,0},
3284b84eec9SHong Zhang {RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(-1.0)/RC(6.0),RC(1.0)/RC(2.0),0}},
3294b84eec9SHong Zhang bs[8] = {RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0),RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0)},
3304b84eec9SHong Zhang b[8] = {RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0),RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0)};
3319566063dSJacob Faibussowitsch PetscCall(TSMPRKRegister(TSMPRKPM3,3,8,&As[0][0],bs,NULL,&A[0][0],b,NULL));
3324b84eec9SHong Zhang }*/
3334b84eec9SHong Zhang
3344b84eec9SHong Zhang {
3359371c9d4SSatish Balay const PetscReal Asb[5][5] =
3369371c9d4SSatish Balay {
3379371c9d4SSatish Balay {0, 0, 0, 0, 0},
3384b84eec9SHong Zhang {RC(1.0) / RC(2.0), 0, 0, 0, 0},
3394b84eec9SHong Zhang {RC(1.0) / RC(2.0), 0, 0, 0, 0},
3404b84eec9SHong Zhang {RC(1.0), 0, 0, 0, 0},
3419371c9d4SSatish Balay {RC(1.0), 0, 0, 0, 0}
3429371c9d4SSatish Balay },
3439371c9d4SSatish Balay Af[5][5] = {{0, 0, 0, 0, 0}, {RC(1.0) / RC(2.0), 0, 0, 0, 0}, {RC(1.0) / RC(4.0), RC(1.0) / RC(4.0), 0, 0, 0}, {RC(1.0) / RC(4.0), RC(1.0) / RC(4.0), RC(1.0) / RC(2.0), 0, 0}, {RC(1.0) / RC(4.0), RC(1.0) / RC(4.0), RC(1.0) / RC(4.0), RC(1.0) / RC(4.0), 0}}, bsb[5] = {RC(1.0) / RC(2.0), 0, 0, 0, RC(1.0) / RC(2.0)}, bf[5] = {RC(1.0) / RC(4.0), RC(1.0) / RC(4.0), RC(1.0) / RC(4.0), RC(1.0) / RC(4.0), 0};
3449371c9d4SSatish Balay const PetscInt rsb[5] = {0, 0, 2, 0, 4};
3459566063dSJacob Faibussowitsch PetscCall(TSMPRKRegister(TSMPRKP2, 2, 5, 1, 1, &Asb[0][0], bsb, NULL, rsb, NULL, NULL, NULL, NULL, &Af[0][0], bf, NULL));
3464b84eec9SHong Zhang }
3474b84eec9SHong Zhang
3484b84eec9SHong Zhang {
3499371c9d4SSatish Balay const PetscReal Asb[10][10] =
3509371c9d4SSatish Balay {
3519371c9d4SSatish Balay {0, 0, 0, 0, 0, 0, 0, 0, 0, 0},
3524b84eec9SHong Zhang {RC(1.0) / RC(4.0), 0, 0, 0, 0, 0, 0, 0, 0, 0},
3534b84eec9SHong Zhang {RC(1.0) / RC(4.0), 0, 0, 0, 0, 0, 0, 0, 0, 0},
3544b84eec9SHong Zhang {RC(1.0) / RC(2.0), 0, 0, 0, 0, 0, 0, 0, 0, 0},
3554b84eec9SHong Zhang {RC(1.0) / RC(2.0), 0, 0, 0, 0, 0, 0, 0, 0, 0},
3564b84eec9SHong Zhang {RC(-1.0) / RC(6.0), 0, 0, 0, RC(2.0) / RC(3.0), 0, 0, 0, 0, 0},
3574b84eec9SHong Zhang {RC(1.0) / RC(12.0), 0, 0, 0, RC(1.0) / RC(6.0), RC(1.0) / RC(2.0), 0, 0, 0, 0},
3584b84eec9SHong Zhang {RC(1.0) / RC(12.0), 0, 0, 0, RC(1.0) / RC(6.0), RC(1.0) / RC(2.0), 0, 0, 0, 0},
3594b84eec9SHong Zhang {RC(1.0) / RC(3.0), 0, 0, 0, RC(-1.0) / RC(3.0), RC(1.0), 0, 0, 0, 0},
3609371c9d4SSatish Balay {RC(1.0) / RC(3.0), 0, 0, 0, RC(-1.0) / RC(3.0), RC(1.0), 0, 0, 0, 0}
3619371c9d4SSatish Balay },
3629371c9d4SSatish Balay Af[10][10] = {{0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, {RC(1.0) / RC(4.0), 0, 0, 0, 0, 0, 0, 0, 0, 0}, {RC(-1.0) / RC(12.0), RC(1.0) / RC(3.0), 0, 0, 0, 0, 0, 0, 0, 0}, {RC(1.0) / RC(6.0), RC(-1.0) / RC(6.0), RC(1.0) / RC(2.0), 0, 0, 0, 0, 0, 0, 0}, {RC(1.0) / RC(12.0), RC(1.0) / RC(6.0), RC(1.0) / RC(6.0), RC(1.0) / RC(12.0), 0, 0, 0, 0, 0, 0}, {RC(1.0) / RC(12.0), RC(1.0) / RC(6.0), RC(1.0) / RC(6.0), RC(1.0) / RC(12.0), 0, 0, 0, 0, 0, 0}, {RC(1.0) / RC(12.0), RC(1.0) / RC(6.0), RC(1.0) / RC(6.0), RC(1.0) / RC(12.0), 0, RC(1.0) / RC(4.0), 0, 0, 0, 0}, {RC(1.0) / RC(12.0), RC(1.0) / RC(6.0), RC(1.0) / RC(6.0), RC(1.0) / RC(12.0), 0, RC(-1.0) / RC(12.0), RC(1.0) / RC(3.0), 0, 0, 0}, {RC(1.0) / RC(12.0), RC(1.0) / RC(6.0), RC(1.0) / RC(6.0), RC(1.0) / RC(12.0), 0, RC(1.0) / RC(6.0), RC(-1.0) / RC(6.0), RC(1.0) / RC(2.0), 0, 0}, {RC(1.0) / RC(12.0), RC(1.0) / RC(6.0), RC(1.0) / RC(6.0), RC(1.0) / RC(12.0), 0, RC(1.0) / RC(12.0), RC(1.0) / RC(6.0), RC(1.0) / RC(6.0), RC(1.0) / RC(12.0), 0}}, bsb[10] = {RC(1.0) / RC(6.0), 0, 0, 0, RC(1.0) / RC(3.0), RC(1.0) / RC(3.0), 0, 0, 0, RC(1.0) / RC(6.0)}, bf[10] = {RC(1.0) / RC(12.0), RC(1.0) / RC(6.0), RC(1.0) / RC(6.0), RC(1.0) / RC(12.0), 0, RC(1.0) / RC(12.0), RC(1.0) / RC(6.0), RC(1.0) / RC(6.0), RC(1.0) / RC(12.0), 0};
3639371c9d4SSatish Balay const PetscInt rsb[10] = {0, 0, 2, 0, 4, 0, 0, 7, 0, 9};
3649566063dSJacob Faibussowitsch PetscCall(TSMPRKRegister(TSMPRKP3, 3, 5, 2, 1, &Asb[0][0], bsb, NULL, rsb, NULL, NULL, NULL, NULL, &Af[0][0], bf, NULL));
3654b84eec9SHong Zhang }
3664b84eec9SHong Zhang #undef RC
3673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3684b84eec9SHong Zhang }
3694b84eec9SHong Zhang
3704b84eec9SHong Zhang /*@C
371bcf0153eSBarry Smith TSMPRKRegisterDestroy - Frees the list of schemes that were registered by `TSMPRKRegister()`.
3724b84eec9SHong Zhang
3734b84eec9SHong Zhang Not Collective
3744b84eec9SHong Zhang
3754b84eec9SHong Zhang Level: advanced
3764b84eec9SHong Zhang
3771cc06b55SBarry Smith .seealso: [](ch_ts), `TSMPRK`, `TSMPRKRegister()`, `TSMPRKRegisterAll()`
3784b84eec9SHong Zhang @*/
TSMPRKRegisterDestroy(void)379d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMPRKRegisterDestroy(void)
380d71ae5a4SJacob Faibussowitsch {
3814b84eec9SHong Zhang MPRKTableauLink link;
3824b84eec9SHong Zhang
3834b84eec9SHong Zhang PetscFunctionBegin;
3844b84eec9SHong Zhang while ((link = MPRKTableauList)) {
3854b84eec9SHong Zhang MPRKTableau t = &link->tab;
3864b84eec9SHong Zhang MPRKTableauList = link->next;
3879566063dSJacob Faibussowitsch PetscCall(PetscFree3(t->Asb, t->bsb, t->csb));
3889566063dSJacob Faibussowitsch PetscCall(PetscFree3(t->Amb, t->bmb, t->cmb));
3899566063dSJacob Faibussowitsch PetscCall(PetscFree3(t->Af, t->bf, t->cf));
3909566063dSJacob Faibussowitsch PetscCall(PetscFree(t->rsb));
3919566063dSJacob Faibussowitsch PetscCall(PetscFree(t->rmb));
3929566063dSJacob Faibussowitsch PetscCall(PetscFree(t->name));
3939566063dSJacob Faibussowitsch PetscCall(PetscFree(link));
3944b84eec9SHong Zhang }
3954b84eec9SHong Zhang TSMPRKRegisterAllCalled = PETSC_FALSE;
3963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3974b84eec9SHong Zhang }
3984b84eec9SHong Zhang
3994b84eec9SHong Zhang /*@C
400bcf0153eSBarry Smith TSMPRKInitializePackage - This function initializes everything in the `TSMPRK` package. It is called
401bcf0153eSBarry Smith from `PetscDLLibraryRegister()` when using dynamic libraries, and on the first call to `TSCreate_MPRK()`
4024b84eec9SHong Zhang when using static libraries.
4034b84eec9SHong Zhang
4044b84eec9SHong Zhang Level: developer
4054b84eec9SHong Zhang
4061cc06b55SBarry Smith .seealso: [](ch_ts), `TSMPRK`, `PetscInitialize()`
4074b84eec9SHong Zhang @*/
TSMPRKInitializePackage(void)408d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMPRKInitializePackage(void)
409d71ae5a4SJacob Faibussowitsch {
4104b84eec9SHong Zhang PetscFunctionBegin;
4113ba16761SJacob Faibussowitsch if (TSMPRKPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
4124b84eec9SHong Zhang TSMPRKPackageInitialized = PETSC_TRUE;
4139566063dSJacob Faibussowitsch PetscCall(TSMPRKRegisterAll());
4149566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSMPRKFinalizePackage));
4153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4164b84eec9SHong Zhang }
4174b84eec9SHong Zhang
4184b84eec9SHong Zhang /*@C
419bcf0153eSBarry Smith TSMPRKFinalizePackage - This function destroys everything in the `TSMPRK` package. It is
420bcf0153eSBarry Smith called from `PetscFinalize()`.
4214b84eec9SHong Zhang
4224b84eec9SHong Zhang Level: developer
4234b84eec9SHong Zhang
4241cc06b55SBarry Smith .seealso: [](ch_ts), `TSMPRK`, `PetscFinalize()`
4254b84eec9SHong Zhang @*/
TSMPRKFinalizePackage(void)426d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMPRKFinalizePackage(void)
427d71ae5a4SJacob Faibussowitsch {
4284b84eec9SHong Zhang PetscFunctionBegin;
4294b84eec9SHong Zhang TSMPRKPackageInitialized = PETSC_FALSE;
4309566063dSJacob Faibussowitsch PetscCall(TSMPRKRegisterDestroy());
4313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4324b84eec9SHong Zhang }
4334b84eec9SHong Zhang
4344b84eec9SHong Zhang /*@C
435bcf0153eSBarry Smith TSMPRKRegister - register a `TSMPRK` scheme by providing the entries in the Butcher tableau
4364b84eec9SHong Zhang
437cc4c1da9SBarry Smith Not Collective, but the same schemes should be registered on all processes on which they will be used, No Fortran Support
4384b84eec9SHong Zhang
4394b84eec9SHong Zhang Input Parameters:
4404b84eec9SHong Zhang + name - identifier for method
4414b84eec9SHong Zhang . order - approximation order of method
4422fe279fdSBarry Smith . sbase - number of stages in the base methods
44379bc8a77SHong Zhang . ratio1 - stepsize ratio at 1st level (e.g. slow/medium)
44479bc8a77SHong Zhang . ratio2 - stepsize ratio at 2nd level (e.g. medium/fast)
44514d0ab18SJacob Faibussowitsch . Asb - stage coefficients for slow components(dimension s*s, row-major)
44614d0ab18SJacob Faibussowitsch . bsb - step completion table for slow components(dimension s)
44714d0ab18SJacob Faibussowitsch . csb - abscissa for slow components(dimension s)
44814d0ab18SJacob Faibussowitsch . rsb - array of flags for repeated stages for slow components (dimension s)
44914d0ab18SJacob Faibussowitsch . Amb - stage coefficients for medium components(dimension s*s, row-major)
45014d0ab18SJacob Faibussowitsch . bmb - step completion table for medium components(dimension s)
45114d0ab18SJacob Faibussowitsch . cmb - abscissa for medium components(dimension s)
45214d0ab18SJacob Faibussowitsch . rmb - array of flags for repeated stages for medium components (dimension s)
4534b84eec9SHong Zhang . Af - stage coefficients for fast components(dimension s*s, row-major)
4544b84eec9SHong Zhang . bf - step completion table for fast components(dimension s)
45514d0ab18SJacob Faibussowitsch - cf - abscissa for fast components(dimension s)
4564b84eec9SHong Zhang
4574b84eec9SHong Zhang Level: advanced
4584b84eec9SHong Zhang
459bcf0153eSBarry Smith Note:
460bcf0153eSBarry Smith Several `TSMPRK` methods are provided, this function is only needed to create new methods.
461bcf0153eSBarry Smith
4621cc06b55SBarry Smith .seealso: [](ch_ts), `TSMPRK`
4634b84eec9SHong Zhang @*/
TSMPRKRegister(TSMPRKType name,PetscInt order,PetscInt sbase,PetscInt ratio1,PetscInt ratio2,const PetscReal Asb[],const PetscReal bsb[],const PetscReal csb[],const PetscInt rsb[],const PetscReal Amb[],const PetscReal bmb[],const PetscReal cmb[],const PetscInt rmb[],const PetscReal Af[],const PetscReal bf[],const PetscReal cf[])464d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMPRKRegister(TSMPRKType name, PetscInt order, PetscInt sbase, PetscInt ratio1, PetscInt ratio2, const PetscReal Asb[], const PetscReal bsb[], const PetscReal csb[], const PetscInt rsb[], const PetscReal Amb[], const PetscReal bmb[], const PetscReal cmb[], const PetscInt rmb[], const PetscReal Af[], const PetscReal bf[], const PetscReal cf[])
465d71ae5a4SJacob Faibussowitsch {
4664b84eec9SHong Zhang MPRKTableauLink link;
4674b84eec9SHong Zhang MPRKTableau t;
46879bc8a77SHong Zhang PetscInt s, i, j;
4694b84eec9SHong Zhang
4704b84eec9SHong Zhang PetscFunctionBegin;
4714f572ea9SToby Isaac PetscAssertPointer(name, 1);
4724f572ea9SToby Isaac PetscAssertPointer(Asb, 6);
4734f572ea9SToby Isaac if (bsb) PetscAssertPointer(bsb, 7);
4744f572ea9SToby Isaac if (csb) PetscAssertPointer(csb, 8);
4754f572ea9SToby Isaac if (rsb) PetscAssertPointer(rsb, 9);
4764f572ea9SToby Isaac if (Amb) PetscAssertPointer(Amb, 10);
4774f572ea9SToby Isaac if (bmb) PetscAssertPointer(bmb, 11);
4784f572ea9SToby Isaac if (cmb) PetscAssertPointer(cmb, 12);
4794f572ea9SToby Isaac if (rmb) PetscAssertPointer(rmb, 13);
4804f572ea9SToby Isaac PetscAssertPointer(Af, 14);
4814f572ea9SToby Isaac if (bf) PetscAssertPointer(bf, 15);
4824f572ea9SToby Isaac if (cf) PetscAssertPointer(cf, 16);
4834b84eec9SHong Zhang
4849566063dSJacob Faibussowitsch PetscCall(PetscNew(&link));
4854b84eec9SHong Zhang t = &link->tab;
4864b84eec9SHong Zhang
4879566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &t->name));
48879bc8a77SHong Zhang s = sbase * ratio1 * ratio2; /* this is the dimension of the matrices below */
4894b84eec9SHong Zhang t->order = order;
49079bc8a77SHong Zhang t->sbase = sbase;
4914b84eec9SHong Zhang t->s = s;
49219c2959aSHong Zhang t->np = 2;
49379bc8a77SHong Zhang
4949566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(s * s, &t->Af, s, &t->bf, s, &t->cf));
4959566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->Af, Af, s * s));
4964b84eec9SHong Zhang if (bf) {
4979566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->bf, bf, s));
49819c2959aSHong Zhang } else
4994b84eec9SHong Zhang for (i = 0; i < s; i++) t->bf[i] = Af[(s - 1) * s + i];
5004b84eec9SHong Zhang if (cf) {
5019566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->cf, cf, s));
50219c2959aSHong Zhang } else {
5034b84eec9SHong Zhang for (i = 0; i < s; i++)
5049371c9d4SSatish Balay for (j = 0, t->cf[i] = 0; j < s; j++) t->cf[i] += Af[i * s + j];
5054b84eec9SHong Zhang }
50619c2959aSHong Zhang
50719c2959aSHong Zhang if (Amb) {
50819c2959aSHong Zhang t->np = 3;
5099566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(s * s, &t->Amb, s, &t->bmb, s, &t->cmb));
5109566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->Amb, Amb, s * s));
51119c2959aSHong Zhang if (bmb) {
5129566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->bmb, bmb, s));
51319c2959aSHong Zhang } else {
51419c2959aSHong Zhang for (i = 0; i < s; i++) t->bmb[i] = Amb[(s - 1) * s + i];
5154b84eec9SHong Zhang }
51619c2959aSHong Zhang if (cmb) {
5179566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->cmb, cmb, s));
51819c2959aSHong Zhang } else {
5194b84eec9SHong Zhang for (i = 0; i < s; i++)
5209371c9d4SSatish Balay for (j = 0, t->cmb[i] = 0; j < s; j++) t->cmb[i] += Amb[i * s + j];
52119c2959aSHong Zhang }
52219c2959aSHong Zhang if (rmb) {
5239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s, &t->rmb));
5249566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->rmb, rmb, s));
525071fcb05SBarry Smith } else {
5269566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(s, &t->rmb));
52719c2959aSHong Zhang }
52819c2959aSHong Zhang }
52919c2959aSHong Zhang
5309566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(s * s, &t->Asb, s, &t->bsb, s, &t->csb));
5319566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->Asb, Asb, s * s));
53219c2959aSHong Zhang if (bsb) {
5339566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->bsb, bsb, s));
53419c2959aSHong Zhang } else
53519c2959aSHong Zhang for (i = 0; i < s; i++) t->bsb[i] = Asb[(s - 1) * s + i];
53619c2959aSHong Zhang if (csb) {
5379566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->csb, csb, s));
53819c2959aSHong Zhang } else {
53919c2959aSHong Zhang for (i = 0; i < s; i++)
5409371c9d4SSatish Balay for (j = 0, t->csb[i] = 0; j < s; j++) t->csb[i] += Asb[i * s + j];
54119c2959aSHong Zhang }
54219c2959aSHong Zhang if (rsb) {
5439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s, &t->rsb));
5449566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->rsb, rsb, s));
545071fcb05SBarry Smith } else {
5469566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(s, &t->rsb));
5474b84eec9SHong Zhang }
5484b84eec9SHong Zhang link->next = MPRKTableauList;
5494b84eec9SHong Zhang MPRKTableauList = link;
5503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5514b84eec9SHong Zhang }
5524b84eec9SHong Zhang
TSMPRKSetSplits(TS ts)553d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSMPRKSetSplits(TS ts)
554d71ae5a4SJacob Faibussowitsch {
5554b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data;
55619c2959aSHong Zhang MPRKTableau tab = mprk->tableau;
5574b84eec9SHong Zhang DM dm, subdm, newdm;
5584b84eec9SHong Zhang
5594b84eec9SHong Zhang PetscFunctionBegin;
5609566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts, "slow", &mprk->subts_slow));
5619566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts, "fast", &mprk->subts_fast));
5623c633725SBarry Smith PetscCheck(mprk->subts_slow && mprk->subts_fast, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must set up the RHSFunctions for 'slow' and 'fast' components using TSRHSSplitSetRHSFunction() or calling TSSetRHSFunction() for each sub-TS");
5634b84eec9SHong Zhang
56419c2959aSHong Zhang /* Only copy the DM */
5659566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm));
56619c2959aSHong Zhang
5679566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts, "slowbuffer", &mprk->subts_slowbuffer));
56819c2959aSHong Zhang if (!mprk->subts_slowbuffer) {
56919c2959aSHong Zhang mprk->subts_slowbuffer = mprk->subts_slow;
57019c2959aSHong Zhang mprk->subts_slow = NULL;
57119c2959aSHong Zhang }
57219c2959aSHong Zhang if (mprk->subts_slow) {
5739566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &newdm));
5749566063dSJacob Faibussowitsch PetscCall(TSGetDM(mprk->subts_slow, &subdm));
5759566063dSJacob Faibussowitsch PetscCall(DMCopyDMTS(subdm, newdm));
5769566063dSJacob Faibussowitsch PetscCall(DMCopyDMSNES(subdm, newdm));
5779566063dSJacob Faibussowitsch PetscCall(TSSetDM(mprk->subts_slow, newdm));
5789566063dSJacob Faibussowitsch PetscCall(DMDestroy(&newdm));
57919c2959aSHong Zhang }
5809566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &newdm));
5819566063dSJacob Faibussowitsch PetscCall(TSGetDM(mprk->subts_slowbuffer, &subdm));
5829566063dSJacob Faibussowitsch PetscCall(DMCopyDMTS(subdm, newdm));
5839566063dSJacob Faibussowitsch PetscCall(DMCopyDMSNES(subdm, newdm));
5849566063dSJacob Faibussowitsch PetscCall(TSSetDM(mprk->subts_slowbuffer, newdm));
5859566063dSJacob Faibussowitsch PetscCall(DMDestroy(&newdm));
58619c2959aSHong Zhang
5879566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &newdm));
5889566063dSJacob Faibussowitsch PetscCall(TSGetDM(mprk->subts_fast, &subdm));
5899566063dSJacob Faibussowitsch PetscCall(DMCopyDMTS(subdm, newdm));
5909566063dSJacob Faibussowitsch PetscCall(DMCopyDMSNES(subdm, newdm));
5919566063dSJacob Faibussowitsch PetscCall(TSSetDM(mprk->subts_fast, newdm));
5929566063dSJacob Faibussowitsch PetscCall(DMDestroy(&newdm));
59319c2959aSHong Zhang
59419c2959aSHong Zhang if (tab->np == 3) {
5959566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts, "medium", &mprk->subts_medium));
5969566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts, "mediumbuffer", &mprk->subts_mediumbuffer));
59719c2959aSHong Zhang if (mprk->subts_medium && !mprk->subts_mediumbuffer) {
59819c2959aSHong Zhang mprk->subts_mediumbuffer = mprk->subts_medium;
59919c2959aSHong Zhang mprk->subts_medium = NULL;
60019c2959aSHong Zhang }
60119c2959aSHong Zhang if (mprk->subts_medium) {
6029566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &newdm));
6039566063dSJacob Faibussowitsch PetscCall(TSGetDM(mprk->subts_medium, &subdm));
6049566063dSJacob Faibussowitsch PetscCall(DMCopyDMTS(subdm, newdm));
6059566063dSJacob Faibussowitsch PetscCall(DMCopyDMSNES(subdm, newdm));
6069566063dSJacob Faibussowitsch PetscCall(TSSetDM(mprk->subts_medium, newdm));
6079566063dSJacob Faibussowitsch PetscCall(DMDestroy(&newdm));
60819c2959aSHong Zhang }
6099566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &newdm));
6109566063dSJacob Faibussowitsch PetscCall(TSGetDM(mprk->subts_mediumbuffer, &subdm));
6119566063dSJacob Faibussowitsch PetscCall(DMCopyDMTS(subdm, newdm));
6129566063dSJacob Faibussowitsch PetscCall(DMCopyDMSNES(subdm, newdm));
6139566063dSJacob Faibussowitsch PetscCall(TSSetDM(mprk->subts_mediumbuffer, newdm));
6149566063dSJacob Faibussowitsch PetscCall(DMDestroy(&newdm));
61519c2959aSHong Zhang }
6163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6174b84eec9SHong Zhang }
6184b84eec9SHong Zhang
6194b84eec9SHong Zhang /*
6204b84eec9SHong Zhang This if for nonsplit RHS MPRK
6214b84eec9SHong Zhang The step completion formula is
6224b84eec9SHong Zhang
6234b84eec9SHong Zhang x1 = x0 + h b^T YdotRHS
6244b84eec9SHong Zhang
6254b84eec9SHong Zhang */
TSEvaluateStep_MPRK(TS ts,PetscInt order,Vec X,PetscBool * done)626d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateStep_MPRK(TS ts, PetscInt order, Vec X, PetscBool *done)
627d71ae5a4SJacob Faibussowitsch {
6284b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data;
6294b84eec9SHong Zhang MPRKTableau tab = mprk->tableau;
6307c0df07dSHong Zhang PetscScalar *wf = mprk->work_fast;
6314b84eec9SHong Zhang PetscReal h = ts->time_step;
6324b84eec9SHong Zhang PetscInt s = tab->s, j;
6334b84eec9SHong Zhang
6344b84eec9SHong Zhang PetscFunctionBegin;
6354b84eec9SHong Zhang for (j = 0; j < s; j++) wf[j] = h * tab->bf[j];
6369566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X));
6379566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, wf, mprk->YdotRHS));
6383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6394b84eec9SHong Zhang }
6404b84eec9SHong Zhang
TSStep_MPRK(TS ts)641d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_MPRK(TS ts)
642d71ae5a4SJacob Faibussowitsch {
6434b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data;
6447c0df07dSHong Zhang Vec *Y = mprk->Y, *YdotRHS = mprk->YdotRHS, *YdotRHS_fast = mprk->YdotRHS_fast, *YdotRHS_slow = mprk->YdotRHS_slow, *YdotRHS_slowbuffer = mprk->YdotRHS_slowbuffer;
6457c0df07dSHong Zhang Vec Yslow, Yslowbuffer, Yfast;
6464b84eec9SHong Zhang MPRKTableau tab = mprk->tableau;
6474b84eec9SHong Zhang const PetscInt s = tab->s;
64819c2959aSHong Zhang const PetscReal *Af = tab->Af, *cf = tab->cf, *Asb = tab->Asb, *csb = tab->csb;
649ebd5ed4eSHong Zhang PetscScalar *wf = mprk->work_fast, *wsb = mprk->work_slowbuffer;
650ebd5ed4eSHong Zhang PetscInt i, j;
6514b84eec9SHong Zhang PetscReal next_time_step = ts->time_step, t = ts->ptime, h = ts->time_step;
6524b84eec9SHong Zhang
6534b84eec9SHong Zhang PetscFunctionBegin;
6544b84eec9SHong Zhang for (i = 0; i < s; i++) {
6554b84eec9SHong Zhang mprk->stage_time = t + h * cf[i];
6569566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, mprk->stage_time));
6579566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, Y[i]));
6584b84eec9SHong Zhang
6597c0df07dSHong Zhang /* slow buffer region */
66019c2959aSHong Zhang for (j = 0; j < i; j++) wsb[j] = h * Asb[i * s + j];
66148a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(VecGetSubVector(YdotRHS[j], mprk->is_slowbuffer, &YdotRHS_slowbuffer[j]));
6629566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_slowbuffer, &Yslowbuffer));
6639566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yslowbuffer, i, wsb, mprk->YdotRHS_slowbuffer));
6649566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_slowbuffer, &Yslowbuffer));
66548a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(VecRestoreSubVector(YdotRHS[j], mprk->is_slowbuffer, &YdotRHS_slowbuffer[j]));
6667c0df07dSHong Zhang /* slow region */
6677c0df07dSHong Zhang if (mprk->is_slow) {
66848a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(VecGetSubVector(YdotRHS[j], mprk->is_slow, &YdotRHS_slow[j]));
6699566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_slow, &Yslow));
6709566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yslow, i, wsb, mprk->YdotRHS_slow));
6719566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_slow, &Yslow));
67248a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(VecRestoreSubVector(YdotRHS[j], mprk->is_slow, &YdotRHS_slow[j]));
6737c0df07dSHong Zhang }
6744b84eec9SHong Zhang
6757c0df07dSHong Zhang /* fast region */
6764b84eec9SHong Zhang for (j = 0; j < i; j++) wf[j] = h * Af[i * s + j];
67748a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(VecGetSubVector(YdotRHS[j], mprk->is_fast, &YdotRHS_fast[j]));
6789566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_fast, &Yfast));
6799566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yfast, i, wf, mprk->YdotRHS_fast));
6809566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_fast, &Yfast));
68148a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(VecRestoreSubVector(YdotRHS[j], mprk->is_fast, &YdotRHS_fast[j]));
6827c0df07dSHong Zhang if (tab->np == 3) {
6837c0df07dSHong Zhang Vec *YdotRHS_medium = mprk->YdotRHS_medium, *YdotRHS_mediumbuffer = mprk->YdotRHS_mediumbuffer;
6847c0df07dSHong Zhang Vec Ymedium, Ymediumbuffer;
685ebd5ed4eSHong Zhang PetscScalar *wmb = mprk->work_mediumbuffer;
686ebd5ed4eSHong Zhang
687ebd5ed4eSHong Zhang for (j = 0; j < i; j++) wmb[j] = h * tab->Amb[i * s + j];
6887c0df07dSHong Zhang /* medium region */
6897c0df07dSHong Zhang if (mprk->is_medium) {
69048a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(VecGetSubVector(YdotRHS[j], mprk->is_medium, &YdotRHS_medium[j]));
6919566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_medium, &Ymedium));
6929566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Ymedium, i, wmb, mprk->YdotRHS_medium));
6939566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_medium, &Ymedium));
69448a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(VecRestoreSubVector(YdotRHS[j], mprk->is_medium, &YdotRHS_medium[j]));
6957c0df07dSHong Zhang }
6967c0df07dSHong Zhang /* medium buffer region */
69748a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(VecGetSubVector(YdotRHS[j], mprk->is_mediumbuffer, &YdotRHS_mediumbuffer[j]));
6989566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_mediumbuffer, &Ymediumbuffer));
6999566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Ymediumbuffer, i, wmb, mprk->YdotRHS_mediumbuffer));
7009566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_mediumbuffer, &Ymediumbuffer));
70148a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(VecRestoreSubVector(YdotRHS[j], mprk->is_mediumbuffer, &YdotRHS_mediumbuffer[j]));
7027c0df07dSHong Zhang }
7039566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, mprk->stage_time, i, Y));
7044b84eec9SHong Zhang /* compute the stage RHS by fast and slow tableau respectively */
7059566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, t + h * csb[i], Y[i], YdotRHS[i]));
7064b84eec9SHong Zhang }
7079566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep(ts, tab->order, ts->vec_sol, NULL));
7084b84eec9SHong Zhang ts->ptime += ts->time_step;
7094b84eec9SHong Zhang ts->time_step = next_time_step;
7103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7114b84eec9SHong Zhang }
7124b84eec9SHong Zhang
7134b84eec9SHong Zhang /*
714f944a40eSHong Zhang This if for the case when split RHS is used
7154b84eec9SHong Zhang The step completion formula is
7164b84eec9SHong Zhang x1 = x0 + h b^T YdotRHS
7174b84eec9SHong Zhang */
TSEvaluateStep_MPRKSPLIT(TS ts,PetscInt order,Vec X,PetscBool * done)718d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateStep_MPRKSPLIT(TS ts, PetscInt order, Vec X, PetscBool *done)
719d71ae5a4SJacob Faibussowitsch {
7204b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data;
7214b84eec9SHong Zhang MPRKTableau tab = mprk->tableau;
7226aad120cSJose E. Roman Vec Xslow, Xfast, Xslowbuffer; /* subvectors for slow and fast components in X respectively */
72319c2959aSHong Zhang PetscScalar *wf = mprk->work_fast, *ws = mprk->work_slow, *wsb = mprk->work_slowbuffer;
7244b84eec9SHong Zhang PetscReal h = ts->time_step;
72579bc8a77SHong Zhang PetscInt s = tab->s, j, computedstages;
7264b84eec9SHong Zhang
7274b84eec9SHong Zhang PetscFunctionBegin;
7289566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X));
72919c2959aSHong Zhang
73019c2959aSHong Zhang /* slow region */
73119c2959aSHong Zhang if (mprk->is_slow) {
73279bc8a77SHong Zhang computedstages = 0;
73319c2959aSHong Zhang for (j = 0; j < s; j++) {
7349849be05SHong Zhang if (tab->rsb[j]) ws[tab->rsb[j] - 1] += h * tab->bsb[j];
73579bc8a77SHong Zhang else ws[computedstages++] = h * tab->bsb[j];
73619c2959aSHong Zhang }
7379566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X, mprk->is_slow, &Xslow));
7389566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Xslow, computedstages, ws, mprk->YdotRHS_slow));
7399566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X, mprk->is_slow, &Xslow));
74019c2959aSHong Zhang }
74119c2959aSHong Zhang
7429d6e09e9SHong Zhang if (tab->np == 3 && mprk->is_medium) {
7439d6e09e9SHong Zhang computedstages = 0;
7449d6e09e9SHong Zhang for (j = 0; j < s; j++) {
7459d6e09e9SHong Zhang if (tab->rmb[j]) wsb[computedstages - tab->sbase + (tab->rmb[j] - 1) % tab->sbase] += h * tab->bsb[j];
7469d6e09e9SHong Zhang else wsb[computedstages++] = h * tab->bsb[j];
7479d6e09e9SHong Zhang }
7489566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X, mprk->is_slowbuffer, &Xslowbuffer));
7499566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Xslowbuffer, computedstages, wsb, mprk->YdotRHS_slowbuffer));
7509566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X, mprk->is_slowbuffer, &Xslowbuffer));
7519d6e09e9SHong Zhang } else {
75219c2959aSHong Zhang /* slow buffer region */
75319c2959aSHong Zhang for (j = 0; j < s; j++) wsb[j] = h * tab->bsb[j];
7549566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X, mprk->is_slowbuffer, &Xslowbuffer));
7559566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Xslowbuffer, s, wsb, mprk->YdotRHS_slowbuffer));
7569566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X, mprk->is_slowbuffer, &Xslowbuffer));
7579d6e09e9SHong Zhang }
75819c2959aSHong Zhang if (tab->np == 3) {
75919c2959aSHong Zhang Vec Xmedium, Xmediumbuffer;
76019c2959aSHong Zhang PetscScalar *wm = mprk->work_medium, *wmb = mprk->work_mediumbuffer;
7619d6e09e9SHong Zhang /* medium region and slow buffer region */
76219c2959aSHong Zhang if (mprk->is_medium) {
76379bc8a77SHong Zhang computedstages = 0;
76419c2959aSHong Zhang for (j = 0; j < s; j++) {
76579bc8a77SHong Zhang if (tab->rmb[j]) wm[computedstages - tab->sbase + (tab->rmb[j] - 1) % tab->sbase] += h * tab->bmb[j];
76679bc8a77SHong Zhang else wm[computedstages++] = h * tab->bmb[j];
76719c2959aSHong Zhang }
7689566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X, mprk->is_medium, &Xmedium));
7699566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Xmedium, computedstages, wm, mprk->YdotRHS_medium));
7709566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X, mprk->is_medium, &Xmedium));
77119c2959aSHong Zhang }
77219c2959aSHong Zhang /* medium buffer region */
77319c2959aSHong Zhang for (j = 0; j < s; j++) wmb[j] = h * tab->bmb[j];
7749566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X, mprk->is_mediumbuffer, &Xmediumbuffer));
7759566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Xmediumbuffer, s, wmb, mprk->YdotRHS_mediumbuffer));
7769566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X, mprk->is_mediumbuffer, &Xmediumbuffer));
77719c2959aSHong Zhang }
77819c2959aSHong Zhang /* fast region */
77919c2959aSHong Zhang for (j = 0; j < s; j++) wf[j] = h * tab->bf[j];
7809566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X, mprk->is_fast, &Xfast));
7819566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Xfast, s, wf, mprk->YdotRHS_fast));
7829566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X, mprk->is_fast, &Xfast));
7833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7844b84eec9SHong Zhang }
7854b84eec9SHong Zhang
TSStep_MPRKSPLIT(TS ts)786d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_MPRKSPLIT(TS ts)
787d71ae5a4SJacob Faibussowitsch {
7884b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data;
7894b84eec9SHong Zhang MPRKTableau tab = mprk->tableau;
79019c2959aSHong Zhang Vec *Y = mprk->Y, *YdotRHS_fast = mprk->YdotRHS_fast, *YdotRHS_slow = mprk->YdotRHS_slow, *YdotRHS_slowbuffer = mprk->YdotRHS_slowbuffer;
79119c2959aSHong Zhang Vec Yslow, Yslowbuffer, Yfast; /* subvectors for slow and fast components in Y[i] respectively */
79219c2959aSHong Zhang PetscInt s = tab->s;
79319c2959aSHong Zhang const PetscReal *Af = tab->Af, *cf = tab->cf, *Asb = tab->Asb, *csb = tab->csb;
79419c2959aSHong Zhang PetscScalar *wf = mprk->work_fast, *ws = mprk->work_slow, *wsb = mprk->work_slowbuffer;
79579bc8a77SHong Zhang PetscInt i, j, computedstages;
7964b84eec9SHong Zhang PetscReal next_time_step = ts->time_step, t = ts->ptime, h = ts->time_step;
7974b84eec9SHong Zhang
7984b84eec9SHong Zhang PetscFunctionBegin;
7994b84eec9SHong Zhang for (i = 0; i < s; i++) {
8004b84eec9SHong Zhang mprk->stage_time = t + h * cf[i];
8019566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, mprk->stage_time));
8024b84eec9SHong Zhang /* calculate the stage value for fast and slow components respectively */
8039566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, Y[i]));
80419c2959aSHong Zhang for (j = 0; j < i; j++) wsb[j] = h * Asb[i * s + j];
80519c2959aSHong Zhang
80619c2959aSHong Zhang /* slow buffer region */
8079d6e09e9SHong Zhang if (tab->np == 3 && mprk->is_medium) {
8089d6e09e9SHong Zhang if (tab->rmb[i]) {
8099566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_slowbuffer, &Yslowbuffer));
8109566063dSJacob Faibussowitsch PetscCall(VecISCopy(Y[tab->rmb[i] - 1], mprk->is_slowbuffer, SCATTER_REVERSE, Yslowbuffer));
8119566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_slowbuffer, &Yslowbuffer));
8129d6e09e9SHong Zhang } else {
8139d6e09e9SHong Zhang PetscScalar *wm = mprk->work_medium;
8149d6e09e9SHong Zhang computedstages = 0;
8159d6e09e9SHong Zhang for (j = 0; j < i; j++) {
8169d6e09e9SHong Zhang if (tab->rmb[j]) wm[computedstages - tab->sbase + (tab->rmb[j] - 1) % tab->sbase] += wsb[j];
8179d6e09e9SHong Zhang else wm[computedstages++] = wsb[j];
8189d6e09e9SHong Zhang }
8199566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_slowbuffer, &Yslowbuffer));
8209566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yslowbuffer, computedstages, wm, YdotRHS_slowbuffer));
8219566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_slowbuffer, &Yslowbuffer));
8229d6e09e9SHong Zhang }
8239d6e09e9SHong Zhang } else {
8249566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_slowbuffer, &Yslowbuffer));
8259566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yslowbuffer, i, wsb, YdotRHS_slowbuffer));
8269566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_slowbuffer, &Yslowbuffer));
8279d6e09e9SHong Zhang }
8289849be05SHong Zhang
82919c2959aSHong Zhang /* slow region */
8309849be05SHong Zhang if (mprk->is_slow) {
8319849be05SHong Zhang if (tab->rsb[i]) { /* repeat previous stage */
8329566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_slow, &Yslow));
8339566063dSJacob Faibussowitsch PetscCall(VecISCopy(Y[tab->rsb[i] - 1], mprk->is_slow, SCATTER_REVERSE, Yslow));
8349566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_slow, &Yslow));
8359849be05SHong Zhang } else {
83679bc8a77SHong Zhang computedstages = 0;
83719c2959aSHong Zhang for (j = 0; j < i; j++) {
83879bc8a77SHong Zhang if (tab->rsb[j]) ws[tab->rsb[j] - 1] += wsb[j];
83979bc8a77SHong Zhang else ws[computedstages++] = wsb[j];
84019c2959aSHong Zhang }
8419566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_slow, &Yslow));
8429566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yslow, computedstages, ws, YdotRHS_slow));
8439566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_slow, &Yslow));
84419c2959aSHong Zhang /* only depends on the slow buffer region */
8459566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(mprk->subts_slow, t + h * csb[i], Y[i], YdotRHS_slow[computedstages]));
84619c2959aSHong Zhang }
8479849be05SHong Zhang }
84819c2959aSHong Zhang
84919c2959aSHong Zhang /* fast region */
85019c2959aSHong Zhang for (j = 0; j < i; j++) wf[j] = h * Af[i * s + j];
8519566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_fast, &Yfast));
8529566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yfast, i, wf, YdotRHS_fast));
8539566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_fast, &Yfast));
85419c2959aSHong Zhang
85519c2959aSHong Zhang if (tab->np == 3) {
85619c2959aSHong Zhang Vec *YdotRHS_medium = mprk->YdotRHS_medium, *YdotRHS_mediumbuffer = mprk->YdotRHS_mediumbuffer;
85719c2959aSHong Zhang Vec Ymedium, Ymediumbuffer;
85819c2959aSHong Zhang const PetscReal *Amb = tab->Amb, *cmb = tab->cmb;
85919c2959aSHong Zhang PetscScalar *wm = mprk->work_medium, *wmb = mprk->work_mediumbuffer;
86019c2959aSHong Zhang
86119c2959aSHong Zhang for (j = 0; j < i; j++) wmb[j] = h * Amb[i * s + j];
86219c2959aSHong Zhang /* medium buffer region */
8639566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_mediumbuffer, &Ymediumbuffer));
8649566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Ymediumbuffer, i, wmb, YdotRHS_mediumbuffer));
8659566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_mediumbuffer, &Ymediumbuffer));
86619c2959aSHong Zhang /* medium region */
86779bc8a77SHong Zhang if (mprk->is_medium) {
86879bc8a77SHong Zhang if (tab->rmb[i]) { /* repeat previous stage */
8699566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_medium, &Ymedium));
8709566063dSJacob Faibussowitsch PetscCall(VecISCopy(Y[tab->rmb[i] - 1], mprk->is_medium, SCATTER_REVERSE, Ymedium));
8719566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_medium, &Ymedium));
87279bc8a77SHong Zhang } else {
87379bc8a77SHong Zhang computedstages = 0;
87479bc8a77SHong Zhang for (j = 0; j < i; j++) {
87579bc8a77SHong Zhang if (tab->rmb[j]) wm[computedstages - tab->sbase + (tab->rmb[j] - 1) % tab->sbase] += wmb[j];
87679bc8a77SHong Zhang else wm[computedstages++] = wmb[j];
87779bc8a77SHong Zhang }
8789566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_medium, &Ymedium));
8799566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Ymedium, computedstages, wm, YdotRHS_medium));
8809566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_medium, &Ymedium));
88179bc8a77SHong Zhang /* only depends on the medium buffer region */
8829566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(mprk->subts_medium, t + h * cmb[i], Y[i], YdotRHS_medium[computedstages]));
8839d6e09e9SHong Zhang /* must be computed after all regions are updated in Y */
8849566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(mprk->subts_slowbuffer, t + h * csb[i], Y[i], YdotRHS_slowbuffer[computedstages]));
88579bc8a77SHong Zhang }
88619c2959aSHong Zhang }
88719c2959aSHong Zhang /* must be computed after fast region and slow region are updated in Y */
8889566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(mprk->subts_mediumbuffer, t + h * cmb[i], Y[i], YdotRHS_mediumbuffer[i]));
88919c2959aSHong Zhang }
89048a46eb9SPierre Jolivet if (!(tab->np == 3 && mprk->is_medium)) PetscCall(TSComputeRHSFunction(mprk->subts_slowbuffer, t + h * csb[i], Y[i], YdotRHS_slowbuffer[i]));
8919566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(mprk->subts_fast, t + h * cf[i], Y[i], YdotRHS_fast[i]));
8924b84eec9SHong Zhang }
89379bc8a77SHong Zhang
8949566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep(ts, tab->order, ts->vec_sol, NULL));
8954b84eec9SHong Zhang ts->ptime += ts->time_step;
8964b84eec9SHong Zhang ts->time_step = next_time_step;
8973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
8984b84eec9SHong Zhang }
8994b84eec9SHong Zhang
TSMPRKTableauReset(TS ts)900d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSMPRKTableauReset(TS ts)
901d71ae5a4SJacob Faibussowitsch {
9024b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data;
9034b84eec9SHong Zhang MPRKTableau tab = mprk->tableau;
9044b84eec9SHong Zhang
9054b84eec9SHong Zhang PetscFunctionBegin;
9063ba16761SJacob Faibussowitsch if (!tab) PetscFunctionReturn(PETSC_SUCCESS);
9079566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->work_fast));
9089566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->work_slow));
9099566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->work_slowbuffer));
9109566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->work_medium));
9119566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->work_mediumbuffer));
9129566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &mprk->Y));
9130fe4d17eSHong Zhang if (ts->use_splitrhsfunction) {
9149566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &mprk->YdotRHS_fast));
9159566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &mprk->YdotRHS_slow));
9169566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &mprk->YdotRHS_slowbuffer));
9179566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &mprk->YdotRHS_medium));
9189566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &mprk->YdotRHS_mediumbuffer));
9190fe4d17eSHong Zhang } else {
9209566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &mprk->YdotRHS));
9211baa6e33SBarry Smith if (mprk->is_slow) PetscCall(PetscFree(mprk->YdotRHS_slow));
9229566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->YdotRHS_slowbuffer));
9237c0df07dSHong Zhang if (tab->np == 3) {
9241baa6e33SBarry Smith if (mprk->is_medium) PetscCall(PetscFree(mprk->YdotRHS_medium));
9259566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->YdotRHS_mediumbuffer));
9267c0df07dSHong Zhang }
9279566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->YdotRHS_fast));
9287c0df07dSHong Zhang }
9293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9304b84eec9SHong Zhang }
9314b84eec9SHong Zhang
TSReset_MPRK(TS ts)932d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_MPRK(TS ts)
933d71ae5a4SJacob Faibussowitsch {
9344b84eec9SHong Zhang PetscFunctionBegin;
9359566063dSJacob Faibussowitsch PetscCall(TSMPRKTableauReset(ts));
9363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9374b84eec9SHong Zhang }
9384b84eec9SHong Zhang
DMCoarsenHook_TSMPRK(DM fine,DM coarse,PetscCtx ctx)939*2a8381b2SBarry Smith static PetscErrorCode DMCoarsenHook_TSMPRK(DM fine, DM coarse, PetscCtx ctx)
940d71ae5a4SJacob Faibussowitsch {
9414b84eec9SHong Zhang PetscFunctionBegin;
9423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9434b84eec9SHong Zhang }
9444b84eec9SHong Zhang
DMRestrictHook_TSMPRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,PetscCtx ctx)945*2a8381b2SBarry Smith static PetscErrorCode DMRestrictHook_TSMPRK(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, PetscCtx ctx)
946d71ae5a4SJacob Faibussowitsch {
9474b84eec9SHong Zhang PetscFunctionBegin;
9483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9494b84eec9SHong Zhang }
9504b84eec9SHong Zhang
DMSubDomainHook_TSMPRK(DM dm,DM subdm,PetscCtx ctx)951*2a8381b2SBarry Smith static PetscErrorCode DMSubDomainHook_TSMPRK(DM dm, DM subdm, PetscCtx ctx)
952d71ae5a4SJacob Faibussowitsch {
9534b84eec9SHong Zhang PetscFunctionBegin;
9543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9554b84eec9SHong Zhang }
9564b84eec9SHong Zhang
DMSubDomainRestrictHook_TSMPRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,PetscCtx ctx)957*2a8381b2SBarry Smith static PetscErrorCode DMSubDomainRestrictHook_TSMPRK(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, PetscCtx ctx)
958d71ae5a4SJacob Faibussowitsch {
9594b84eec9SHong Zhang PetscFunctionBegin;
9603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9614b84eec9SHong Zhang }
9624b84eec9SHong Zhang
TSMPRKTableauSetUp(TS ts)963d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSMPRKTableauSetUp(TS ts)
964d71ae5a4SJacob Faibussowitsch {
9654b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data;
9664b84eec9SHong Zhang MPRKTableau tab = mprk->tableau;
96719c2959aSHong Zhang Vec YdotRHS_slow, YdotRHS_slowbuffer, YdotRHS_medium, YdotRHS_mediumbuffer, YdotRHS_fast;
9684b84eec9SHong Zhang
9694b84eec9SHong Zhang PetscFunctionBegin;
9709566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &mprk->Y));
97148a46eb9SPierre Jolivet if (mprk->is_slow) PetscCall(PetscMalloc1(tab->s, &mprk->work_slow));
9729566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &mprk->work_slowbuffer));
9737c0df07dSHong Zhang if (tab->np == 3) {
97448a46eb9SPierre Jolivet if (mprk->is_medium) PetscCall(PetscMalloc1(tab->s, &mprk->work_medium));
9759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &mprk->work_mediumbuffer));
9767c0df07dSHong Zhang }
9779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &mprk->work_fast));
9787c0df07dSHong Zhang
9790fe4d17eSHong Zhang if (ts->use_splitrhsfunction) {
98019c2959aSHong Zhang if (mprk->is_slow) {
9819566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ts->vec_sol, mprk->is_slow, &YdotRHS_slow));
9829566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(YdotRHS_slow, tab->s, &mprk->YdotRHS_slow));
9839566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ts->vec_sol, mprk->is_slow, &YdotRHS_slow));
98419c2959aSHong Zhang }
9859566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ts->vec_sol, mprk->is_slowbuffer, &YdotRHS_slowbuffer));
9869566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(YdotRHS_slowbuffer, tab->s, &mprk->YdotRHS_slowbuffer));
9879566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ts->vec_sol, mprk->is_slowbuffer, &YdotRHS_slowbuffer));
98819c2959aSHong Zhang if (tab->np == 3) {
98919c2959aSHong Zhang if (mprk->is_medium) {
9909566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ts->vec_sol, mprk->is_medium, &YdotRHS_medium));
9919566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(YdotRHS_medium, tab->s, &mprk->YdotRHS_medium));
9929566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ts->vec_sol, mprk->is_medium, &YdotRHS_medium));
99319c2959aSHong Zhang }
9949566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ts->vec_sol, mprk->is_mediumbuffer, &YdotRHS_mediumbuffer));
9959566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(YdotRHS_mediumbuffer, tab->s, &mprk->YdotRHS_mediumbuffer));
9969566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ts->vec_sol, mprk->is_mediumbuffer, &YdotRHS_mediumbuffer));
99719c2959aSHong Zhang }
9989566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ts->vec_sol, mprk->is_fast, &YdotRHS_fast));
9999566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(YdotRHS_fast, tab->s, &mprk->YdotRHS_fast));
10009566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ts->vec_sol, mprk->is_fast, &YdotRHS_fast));
10010fe4d17eSHong Zhang } else {
10029566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &mprk->YdotRHS));
100348a46eb9SPierre Jolivet if (mprk->is_slow) PetscCall(PetscMalloc1(tab->s, &mprk->YdotRHS_slow));
10049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &mprk->YdotRHS_slowbuffer));
10050fe4d17eSHong Zhang if (tab->np == 3) {
100648a46eb9SPierre Jolivet if (mprk->is_medium) PetscCall(PetscMalloc1(tab->s, &mprk->YdotRHS_medium));
10079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &mprk->YdotRHS_mediumbuffer));
10080fe4d17eSHong Zhang }
10099566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &mprk->YdotRHS_fast));
10104b84eec9SHong Zhang }
10113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
10124b84eec9SHong Zhang }
10134b84eec9SHong Zhang
TSSetUp_MPRK(TS ts)1014d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_MPRK(TS ts)
1015d71ae5a4SJacob Faibussowitsch {
10164b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data;
101719c2959aSHong Zhang MPRKTableau tab = mprk->tableau;
10184b84eec9SHong Zhang DM dm;
10194b84eec9SHong Zhang
10204b84eec9SHong Zhang PetscFunctionBegin;
10219566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts, "slow", &mprk->is_slow));
10229566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts, "fast", &mprk->is_fast));
10233c633725SBarry Smith PetscCheck(mprk->is_slow && mprk->is_fast, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must set up RHSSplits with TSRHSSplitSetIS() using split names 'slow' and 'fast' respectively in order to use the method '%s'", tab->name);
102419c2959aSHong Zhang
102519c2959aSHong Zhang if (tab->np == 3) {
10269566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts, "medium", &mprk->is_medium));
10273c633725SBarry Smith PetscCheck(mprk->is_medium, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must set up RHSSplits with TSRHSSplitSetIS() using split names 'slow' and 'medium' and 'fast' respectively in order to use the method '%s'", tab->name);
10289566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts, "mediumbuffer", &mprk->is_mediumbuffer));
102919c2959aSHong Zhang if (!mprk->is_mediumbuffer) { /* let medium buffer cover whole medium region */
103019c2959aSHong Zhang mprk->is_mediumbuffer = mprk->is_medium;
103119c2959aSHong Zhang mprk->is_medium = NULL;
103219c2959aSHong Zhang }
103319c2959aSHong Zhang }
103419c2959aSHong Zhang
103519c2959aSHong Zhang /* If users do not provide buffer region settings, the solver will do them automatically, but with a performance penalty */
10369566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts, "slowbuffer", &mprk->is_slowbuffer));
103719c2959aSHong Zhang if (!mprk->is_slowbuffer) { /* let slow buffer cover whole slow region */
103819c2959aSHong Zhang mprk->is_slowbuffer = mprk->is_slow;
103919c2959aSHong Zhang mprk->is_slow = NULL;
104019c2959aSHong Zhang }
10419566063dSJacob Faibussowitsch PetscCall(TSCheckImplicitTerm(ts));
10429566063dSJacob Faibussowitsch PetscCall(TSMPRKTableauSetUp(ts));
10439566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm));
10449566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSMPRK, DMRestrictHook_TSMPRK, ts));
10459566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSMPRK, DMSubDomainRestrictHook_TSMPRK, ts));
10460fe4d17eSHong Zhang if (ts->use_splitrhsfunction) {
10470fe4d17eSHong Zhang ts->ops->step = TSStep_MPRKSPLIT;
10480fe4d17eSHong Zhang ts->ops->evaluatestep = TSEvaluateStep_MPRKSPLIT;
10499566063dSJacob Faibussowitsch PetscCall(TSMPRKSetSplits(ts));
10500fe4d17eSHong Zhang } else {
10510fe4d17eSHong Zhang ts->ops->step = TSStep_MPRK;
10520fe4d17eSHong Zhang ts->ops->evaluatestep = TSEvaluateStep_MPRK;
10530fe4d17eSHong Zhang }
10543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
10554b84eec9SHong Zhang }
10564b84eec9SHong Zhang
TSSetFromOptions_MPRK(TS ts,PetscOptionItems PetscOptionsObject)1057ce78bad3SBarry Smith static PetscErrorCode TSSetFromOptions_MPRK(TS ts, PetscOptionItems PetscOptionsObject)
1058d71ae5a4SJacob Faibussowitsch {
10594b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data;
10604b84eec9SHong Zhang
10614b84eec9SHong Zhang PetscFunctionBegin;
1062d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "PRK ODE solver options");
10634b84eec9SHong Zhang {
10644b84eec9SHong Zhang MPRKTableauLink link;
10654b84eec9SHong Zhang PetscInt count, choice;
10664b84eec9SHong Zhang PetscBool flg;
10674b84eec9SHong Zhang const char **namelist;
1068fbccb6d4SPierre Jolivet for (link = MPRKTableauList, count = 0; link; link = link->next, count++);
10699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count, (char ***)&namelist));
10704b84eec9SHong Zhang for (link = MPRKTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name;
10719566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_mprk_type", "Family of MPRK method", "TSMPRKSetType", (const char *const *)namelist, count, mprk->tableau->name, &choice, &flg));
10729566063dSJacob Faibussowitsch if (flg) PetscCall(TSMPRKSetType(ts, namelist[choice]));
10739566063dSJacob Faibussowitsch PetscCall(PetscFree(namelist));
10744b84eec9SHong Zhang }
1075d0609cedSBarry Smith PetscOptionsHeadEnd();
10763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
10774b84eec9SHong Zhang }
10784b84eec9SHong Zhang
TSView_MPRK(TS ts,PetscViewer viewer)1079d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_MPRK(TS ts, PetscViewer viewer)
1080d71ae5a4SJacob Faibussowitsch {
10814b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data;
10829f196a02SMartin Diehl PetscBool isascii;
10834b84eec9SHong Zhang
10844b84eec9SHong Zhang PetscFunctionBegin;
10859f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
10869f196a02SMartin Diehl if (isascii) {
10874b84eec9SHong Zhang MPRKTableau tab = mprk->tableau;
10884b84eec9SHong Zhang TSMPRKType mprktype;
10894b84eec9SHong Zhang char fbuf[512];
10904b84eec9SHong Zhang char sbuf[512];
109119c2959aSHong Zhang PetscInt i;
10929566063dSJacob Faibussowitsch PetscCall(TSMPRKGetType(ts, &mprktype));
10939566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " MPRK type %s\n", mprktype));
109463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Order: %" PetscInt_FMT "\n", tab->order));
109519c2959aSHong Zhang
10969566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(fbuf, sizeof(fbuf), "% 8.6f", tab->s, tab->cf));
10979566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Abscissa cf = %s\n", fbuf));
10989566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Af = \n"));
109919c2959aSHong Zhang for (i = 0; i < tab->s; i++) {
11009566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(fbuf, sizeof(fbuf), "% 8.6f", tab->s, &tab->Af[i * tab->s]));
11019566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " %s\n", fbuf));
110219c2959aSHong Zhang }
11039566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(fbuf, sizeof(fbuf), "% 8.6f", tab->s, tab->bf));
11049566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " bf = %s\n", fbuf));
110519c2959aSHong Zhang
11069566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(sbuf, sizeof(sbuf), "% 8.6f", tab->s, tab->csb));
11079566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Abscissa csb = %s\n", sbuf));
11089566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Asb = \n"));
110919c2959aSHong Zhang for (i = 0; i < tab->s; i++) {
11109566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(sbuf, sizeof(sbuf), "% 8.6f", tab->s, &tab->Asb[i * tab->s]));
11119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " %s\n", sbuf));
111219c2959aSHong Zhang }
11139566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(sbuf, sizeof(sbuf), "% 8.6f", tab->s, tab->bsb));
11149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " bsb = %s\n", sbuf));
111519c2959aSHong Zhang
111619c2959aSHong Zhang if (tab->np == 3) {
111719c2959aSHong Zhang char mbuf[512];
11189566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(mbuf, sizeof(mbuf), "% 8.6f", tab->s, tab->cmb));
11199566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Abscissa cmb = %s\n", mbuf));
11209566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Amb = \n"));
112119c2959aSHong Zhang for (i = 0; i < tab->s; i++) {
11229566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(mbuf, sizeof(mbuf), "% 8.6f", tab->s, &tab->Amb[i * tab->s]));
11239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " %s\n", mbuf));
112419c2959aSHong Zhang }
11259566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(mbuf, sizeof(mbuf), "% 8.6f", tab->s, tab->bmb));
11269566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " bmb = %s\n", mbuf));
112719c2959aSHong Zhang }
11284b84eec9SHong Zhang }
11293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11304b84eec9SHong Zhang }
11314b84eec9SHong Zhang
TSLoad_MPRK(TS ts,PetscViewer viewer)1132d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSLoad_MPRK(TS ts, PetscViewer viewer)
1133d71ae5a4SJacob Faibussowitsch {
11344b84eec9SHong Zhang TSAdapt adapt;
11354b84eec9SHong Zhang
11364b84eec9SHong Zhang PetscFunctionBegin;
11379566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt));
11389566063dSJacob Faibussowitsch PetscCall(TSAdaptLoad(adapt, viewer));
11393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11404b84eec9SHong Zhang }
11414b84eec9SHong Zhang
1142cc4c1da9SBarry Smith /*@
1143bcf0153eSBarry Smith TSMPRKSetType - Set the type of `TSMPRK` scheme
11444b84eec9SHong Zhang
114520f4b53cSBarry Smith Not Collective
11464b84eec9SHong Zhang
1147d8d19677SJose E. Roman Input Parameters:
11484b84eec9SHong Zhang + ts - timestepping context
1149bcf0153eSBarry Smith - mprktype - type of `TSMPRK` scheme
11504b84eec9SHong Zhang
1151bcf0153eSBarry Smith Options Database Key:
1152147403d9SBarry Smith . -ts_mprk_type - <pm2,p2,p3> - select the specific scheme
11534b84eec9SHong Zhang
11544b84eec9SHong Zhang Level: intermediate
11554b84eec9SHong Zhang
11561cc06b55SBarry Smith .seealso: [](ch_ts), `TSMPRKGetType()`, `TSMPRK`, `TSMPRKType`
11574b84eec9SHong Zhang @*/
TSMPRKSetType(TS ts,TSMPRKType mprktype)1158d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMPRKSetType(TS ts, TSMPRKType mprktype)
1159d71ae5a4SJacob Faibussowitsch {
11604b84eec9SHong Zhang PetscFunctionBegin;
11614b84eec9SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
11624f572ea9SToby Isaac PetscAssertPointer(mprktype, 2);
1163cac4c232SBarry Smith PetscTryMethod(ts, "TSMPRKSetType_C", (TS, TSMPRKType), (ts, mprktype));
11643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11654b84eec9SHong Zhang }
11664b84eec9SHong Zhang
1167cc4c1da9SBarry Smith /*@
1168bcf0153eSBarry Smith TSMPRKGetType - Get the type of `TSMPRK` scheme
11694b84eec9SHong Zhang
117020f4b53cSBarry Smith Not Collective
11714b84eec9SHong Zhang
11724b84eec9SHong Zhang Input Parameter:
11734b84eec9SHong Zhang . ts - timestepping context
11744b84eec9SHong Zhang
11754b84eec9SHong Zhang Output Parameter:
1176bcf0153eSBarry Smith . mprktype - type of `TSMPRK` scheme
11774b84eec9SHong Zhang
11784b84eec9SHong Zhang Level: intermediate
11794b84eec9SHong Zhang
118042747ad1SJacob Faibussowitsch .seealso: [](ch_ts), `TSMPRK`
11814b84eec9SHong Zhang @*/
TSMPRKGetType(TS ts,TSMPRKType * mprktype)1182d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMPRKGetType(TS ts, TSMPRKType *mprktype)
1183d71ae5a4SJacob Faibussowitsch {
11844b84eec9SHong Zhang PetscFunctionBegin;
11854b84eec9SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1186cac4c232SBarry Smith PetscUseMethod(ts, "TSMPRKGetType_C", (TS, TSMPRKType *), (ts, mprktype));
11873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11884b84eec9SHong Zhang }
11894b84eec9SHong Zhang
TSMPRKGetType_MPRK(TS ts,TSMPRKType * mprktype)1190d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSMPRKGetType_MPRK(TS ts, TSMPRKType *mprktype)
1191d71ae5a4SJacob Faibussowitsch {
11924b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data;
11934b84eec9SHong Zhang
11944b84eec9SHong Zhang PetscFunctionBegin;
11954b84eec9SHong Zhang *mprktype = mprk->tableau->name;
11963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11974b84eec9SHong Zhang }
11984b84eec9SHong Zhang
TSMPRKSetType_MPRK(TS ts,TSMPRKType mprktype)1199d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSMPRKSetType_MPRK(TS ts, TSMPRKType mprktype)
1200d71ae5a4SJacob Faibussowitsch {
12014b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data;
12024b84eec9SHong Zhang PetscBool match;
12034b84eec9SHong Zhang MPRKTableauLink link;
12044b84eec9SHong Zhang
12054b84eec9SHong Zhang PetscFunctionBegin;
12064b84eec9SHong Zhang if (mprk->tableau) {
12079566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(mprk->tableau->name, mprktype, &match));
12083ba16761SJacob Faibussowitsch if (match) PetscFunctionReturn(PETSC_SUCCESS);
12094b84eec9SHong Zhang }
12104b84eec9SHong Zhang for (link = MPRKTableauList; link; link = link->next) {
12119566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(link->tab.name, mprktype, &match));
12124b84eec9SHong Zhang if (match) {
12139566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSMPRKTableauReset(ts));
12144b84eec9SHong Zhang mprk->tableau = &link->tab;
12159566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSMPRKTableauSetUp(ts));
12163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
12174b84eec9SHong Zhang }
12184b84eec9SHong Zhang }
121998921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", mprktype);
12204b84eec9SHong Zhang }
12214b84eec9SHong Zhang
TSGetStages_MPRK(TS ts,PetscInt * ns,Vec ** Y)1222d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGetStages_MPRK(TS ts, PetscInt *ns, Vec **Y)
1223d71ae5a4SJacob Faibussowitsch {
12244b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data;
12254b84eec9SHong Zhang
12264b84eec9SHong Zhang PetscFunctionBegin;
12274b84eec9SHong Zhang *ns = mprk->tableau->s;
12284b84eec9SHong Zhang if (Y) *Y = mprk->Y;
12293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
12304b84eec9SHong Zhang }
12314b84eec9SHong Zhang
TSDestroy_MPRK(TS ts)1232d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_MPRK(TS ts)
1233d71ae5a4SJacob Faibussowitsch {
12344b84eec9SHong Zhang PetscFunctionBegin;
12359566063dSJacob Faibussowitsch PetscCall(TSReset_MPRK(ts));
12364b84eec9SHong Zhang if (ts->dm) {
12379566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSMPRK, DMRestrictHook_TSMPRK, ts));
12389566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSMPRK, DMSubDomainRestrictHook_TSMPRK, ts));
12394b84eec9SHong Zhang }
12409566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data));
12419566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSMPRKGetType_C", NULL));
12429566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSMPRKSetType_C", NULL));
12433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
12444b84eec9SHong Zhang }
12454b84eec9SHong Zhang
12464b84eec9SHong Zhang /*MC
1247f944a40eSHong Zhang TSMPRK - ODE solver using Multirate Partitioned Runge-Kutta schemes
12484b84eec9SHong Zhang
1249dd8e379bSPierre Jolivet The user should provide the right-hand side of the equation using `TSSetRHSFunction()`.
12504b84eec9SHong Zhang
12514b84eec9SHong Zhang Level: beginner
12524b84eec9SHong Zhang
1253bcf0153eSBarry Smith Note:
1254bcf0153eSBarry Smith The default is `TSMPRKPM2`, it can be changed with `TSMPRKSetType()` or -ts_mprk_type
12554b84eec9SHong Zhang
12561cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSMPRKSetType()`, `TSMPRKGetType()`, `TSMPRKType`, `TSMPRKRegister()`, `TSMPRKSetMultirateType()`
1257bcf0153eSBarry Smith `TSMPRKM2`, `TSMPRKM3`, `TSMPRKRFSMR3`, `TSMPRKRFSMR2`, `TSType`
12584b84eec9SHong Zhang M*/
TSCreate_MPRK(TS ts)1259d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_MPRK(TS ts)
1260d71ae5a4SJacob Faibussowitsch {
12614b84eec9SHong Zhang TS_MPRK *mprk;
12624b84eec9SHong Zhang
12634b84eec9SHong Zhang PetscFunctionBegin;
12649566063dSJacob Faibussowitsch PetscCall(TSMPRKInitializePackage());
12654b84eec9SHong Zhang
12664b84eec9SHong Zhang ts->ops->reset = TSReset_MPRK;
12674b84eec9SHong Zhang ts->ops->destroy = TSDestroy_MPRK;
12684b84eec9SHong Zhang ts->ops->view = TSView_MPRK;
12694b84eec9SHong Zhang ts->ops->load = TSLoad_MPRK;
12704b84eec9SHong Zhang ts->ops->setup = TSSetUp_MPRK;
12714b84eec9SHong Zhang ts->ops->step = TSStep_MPRK;
12724b84eec9SHong Zhang ts->ops->evaluatestep = TSEvaluateStep_MPRK;
12734b84eec9SHong Zhang ts->ops->setfromoptions = TSSetFromOptions_MPRK;
12744b84eec9SHong Zhang ts->ops->getstages = TSGetStages_MPRK;
12754b84eec9SHong Zhang
12764dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&mprk));
12774b84eec9SHong Zhang ts->data = (void *)mprk;
12784b84eec9SHong Zhang
12799566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSMPRKGetType_C", TSMPRKGetType_MPRK));
12809566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSMPRKSetType_C", TSMPRKSetType_MPRK));
12814b84eec9SHong Zhang
12829566063dSJacob Faibussowitsch PetscCall(TSMPRKSetType(ts, TSMPRKDefault));
12833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
12844b84eec9SHong Zhang }
1285