xref: /petsc/src/ts/impls/multirate/mprk.c (revision 19c2959aedcbb17969a32880b9fa35aa707cd651)
14b84eec9SHong Zhang /*
24b84eec9SHong Zhang   Code for time stepping with the Partitioned Runge-Kutta method
34b84eec9SHong Zhang 
44b84eec9SHong Zhang   Notes:
54b84eec9SHong Zhang   1) The general system is written as
64b84eec9SHong Zhang      Udot = F(t,U) for nonsplit RHS multi-rate RK,
74b84eec9SHong Zhang      user should give the indexes for both slow and fast components;
84b84eec9SHong Zhang   2) The general system is written as
94b84eec9SHong Zhang      Usdot = Fs(t,Us,Uf)
104b84eec9SHong Zhang      Ufdot = Ff(t,Us,Uf) for split RHS multi-rate RK,
114b84eec9SHong Zhang      user should partioned RHS by themselves and also provide the indexes for both slow and fast components.
12*19c2959aSHong 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'.
13*19c2959aSHong 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.
14*19c2959aSHong Zhang 
15*19c2959aSHong Zhang   Reference:
16*19c2959aSHong Zhang   Emil M. Constantinescu, Adrian Sandu, Multirate Timestepping Methods for Hyperbolic Conservation Laws, Journal of Scientific Computing 2007
17*19c2959aSHong Zhang 
184b84eec9SHong Zhang */
194b84eec9SHong Zhang 
204b84eec9SHong Zhang #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
214b84eec9SHong Zhang #include <petscdm.h>
224b84eec9SHong Zhang 
23*19c2959aSHong 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 */
314b84eec9SHong Zhang   PetscInt  s;                              /* Number of stages */
32*19c2959aSHong Zhang   PetscInt  np;                             /* Number of partitions */
334b84eec9SHong Zhang   PetscReal *Af,*bf,*cf;                    /* Tableau for fast components */
34*19c2959aSHong Zhang   PetscReal *Amb,*bmb,*cmb;                 /* Tableau for medium components */
35*19c2959aSHong Zhang   PetscInt  *rmb;                           /* Array of flags for repeated stages in medium method */
36*19c2959aSHong Zhang   PetscReal *Asb,*bsb,*csb;                 /* Tableau for slow components */
37*19c2959aSHong Zhang   PetscInt  *rsb;                           /* Array of flags for repeated staged in slow method*/
384b84eec9SHong Zhang };
394b84eec9SHong Zhang typedef struct _MPRKTableauLink *MPRKTableauLink;
404b84eec9SHong Zhang struct _MPRKTableauLink {
414b84eec9SHong Zhang   struct _MPRKTableau tab;
424b84eec9SHong Zhang   MPRKTableauLink     next;
434b84eec9SHong Zhang };
444b84eec9SHong Zhang static MPRKTableauLink MPRKTableauList;
454b84eec9SHong Zhang 
464b84eec9SHong Zhang typedef struct {
474b84eec9SHong Zhang   MPRKTableau         tableau;
484b84eec9SHong Zhang   TSMPRKMultirateType mprkmtype;
494b84eec9SHong Zhang   Vec                 *Y;                          /* States computed during the step                           */
504b84eec9SHong Zhang   Vec                 Ytmp;
514b84eec9SHong Zhang   Vec                 *YdotRHS_slow;               /* Function evaluations by slow tableau for slow components  */
52*19c2959aSHong Zhang   Vec                 *YdotRHS_slowbuffer;         /* Function evaluations by slow tableau for slow components  */
53*19c2959aSHong Zhang   Vec                 *YdotRHS_medium;             /* Function evaluations by slow tableau for slow components  */
54*19c2959aSHong Zhang   Vec                 *YdotRHS_mediumbuffer;       /* Function evaluations by slow tableau for slow components  */
55*19c2959aSHong Zhang   Vec                 *YdotRHS_fast;               /* Function evaluations by fast tableau for fast components  */
564b84eec9SHong Zhang   PetscScalar         *work_slow;                  /* Scalar work_slow by slow tableau                          */
57*19c2959aSHong Zhang   PetscScalar         *work_slowbuffer;            /* Scalar work_slow by slow tableau                          */
58*19c2959aSHong Zhang   PetscScalar         *work_medium;                /* Scalar work_slow by medium tableau                        */
59*19c2959aSHong Zhang   PetscScalar         *work_mediumbuffer;          /* Scalar work_slow by medium tableau                        */
60*19c2959aSHong 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;
65*19c2959aSHong Zhang   IS                  is_slow,is_slowbuffer,is_medium,is_mediumbuffer,is_fast;
66*19c2959aSHong Zhang   TS                  subts_slow,subts_slowbuffer,subts_medium,subts_mediumbuffer,subts_fast;
674b84eec9SHong Zhang } TS_MPRK;
684b84eec9SHong Zhang 
69*19c2959aSHong Zhang static PetscErrorCode TSMPRKGenerateTableau2(PetscInt ratio,PetscInt s,const PetscReal Abase[],const PetscReal bbase[],PetscReal A1[],PetscReal b1[],PetscReal A2[],PetscReal b2[])
70*19c2959aSHong Zhang {
71*19c2959aSHong Zhang   PetscInt i,j,k,l;
72*19c2959aSHong Zhang 
73*19c2959aSHong Zhang   PetscFunctionBegin;
74*19c2959aSHong Zhang   for (k=0; k<ratio; k++) {
75*19c2959aSHong Zhang     /* diagonal blocks */
76*19c2959aSHong Zhang     for (i=0; i<s; i++)
77*19c2959aSHong Zhang       for (j=0; j<s; j++) {
78*19c2959aSHong Zhang         A1[(k*s+i)*ratio*s+k*s+j] = Abase[i*s+j];
79*19c2959aSHong Zhang         A2[(k*s+i)*ratio*s+k*s+j] = Abase[i*s+j]/ratio;
80*19c2959aSHong Zhang       }
81*19c2959aSHong Zhang     /* off diagonal blocks */
82*19c2959aSHong Zhang     for (l=0; l<k; l++)
83*19c2959aSHong Zhang       for (i=0; i<s; i++)
84*19c2959aSHong Zhang         for (j=0; j<s; j++)
85*19c2959aSHong Zhang           A2[(k*s+i)*ratio*s+l*s+j] = bbase[j]/ratio;
86*19c2959aSHong Zhang     for (j=0; j<s; j++) {
87*19c2959aSHong Zhang       b1[k*s+j] = bbase[j]/ratio;
88*19c2959aSHong Zhang       b2[k*s+j] = bbase[j]/ratio;
89*19c2959aSHong Zhang     }
90*19c2959aSHong Zhang   }
91*19c2959aSHong Zhang   PetscFunctionReturn(0);
92*19c2959aSHong Zhang }
93*19c2959aSHong Zhang 
94*19c2959aSHong Zhang 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[])
95*19c2959aSHong Zhang {
96*19c2959aSHong Zhang   PetscInt i,j,k,l,m,n;
97*19c2959aSHong Zhang 
98*19c2959aSHong Zhang   PetscFunctionBegin;
99*19c2959aSHong Zhang   for (k=0; k<ratio; k++) { /* diagonal blocks of size ratio*s by ratio*s */
100*19c2959aSHong Zhang     for (l=0; l<ratio; l++) /* diagonal sub-blocks of size s by s */
101*19c2959aSHong Zhang       for (i=0; i<s; i++)
102*19c2959aSHong Zhang         for (j=0; j<s; j++) {
103*19c2959aSHong Zhang           A1[((k*ratio+l)*s+i)*ratio*ratio*s+(k*ratio+l)*s+j] = Abase[i*s+j];
104*19c2959aSHong Zhang           A2[((k*ratio+l)*s+i)*ratio*ratio*s+(k*ratio+l)*s+j] = Abase[i*s+j]/ratio;
105*19c2959aSHong Zhang           A3[((k*ratio+l)*s+i)*ratio*ratio*s+(k*ratio+l)*s+j] = Abase[i*s+j]/ratio/ratio;
106*19c2959aSHong Zhang         }
107*19c2959aSHong Zhang     for (l=0; l<k; l++) /* off-diagonal blocks of size ratio*s by ratio*s */
108*19c2959aSHong Zhang       for (m=0; m<ratio; m++)
109*19c2959aSHong Zhang         for (n=0; n<ratio; n++)
110*19c2959aSHong Zhang           for (i=0; i<s; i++)
111*19c2959aSHong Zhang             for (j=0; j<s; j++) {
112*19c2959aSHong Zhang                A2[((k*ratio+m)*s+i)*ratio*ratio*s+(l*ratio+n)*s+j] = bbase[j]/ratio/ratio;
113*19c2959aSHong Zhang                A3[((k*ratio+m)*s+i)*ratio*ratio*s+(l*ratio+n)*s+j] = bbase[j]/ratio/ratio;
114*19c2959aSHong Zhang             }
115*19c2959aSHong Zhang     for (m=0; m<ratio; m++)
116*19c2959aSHong Zhang       for (n=0; n<m; n++) /* off-diagonal sub-blocks of size s by s in the diagonal blocks */
117*19c2959aSHong Zhang           for (i=0; i<s; i++)
118*19c2959aSHong Zhang             for (j=0; j<s; j++)
119*19c2959aSHong Zhang                A3[((k*ratio+m)*s+i)*ratio*ratio*s+(k*ratio+n)*s+j] = bbase[j]/ratio/ratio;
120*19c2959aSHong Zhang     for (n=0; n<ratio; n++)
121*19c2959aSHong Zhang       for (j=0; j<s; j++) {
122*19c2959aSHong Zhang         b1[(k*ratio+n)*s+j] = bbase[j]/ratio/ratio;
123*19c2959aSHong Zhang         b2[(k*ratio+n)*s+j] = bbase[j]/ratio/ratio;
124*19c2959aSHong Zhang         b3[(k*ratio+n)*s+j] = bbase[j]/ratio/ratio;
125*19c2959aSHong Zhang       }
126*19c2959aSHong Zhang   }
127*19c2959aSHong Zhang   PetscFunctionReturn(0);
128*19c2959aSHong Zhang }
129*19c2959aSHong Zhang 
1304b84eec9SHong Zhang /*MC
131*19c2959aSHong Zhang      TSMPRK2A22 - Second Order Partitioned Runge Kutta scheme based on RK2A.
1324b84eec9SHong Zhang 
133*19c2959aSHong Zhang      This method has four stages for slow and fast parts. The refinement factor of the stepsize is 2.
134*19c2959aSHong Zhang      r = 2, np = 2
1354b84eec9SHong Zhang      Options database:
136*19c2959aSHong Zhang .     -ts_mprk_type 2a22
1374b84eec9SHong Zhang 
1384b84eec9SHong Zhang      Level: advanced
1394b84eec9SHong Zhang 
1404b84eec9SHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType()
1414b84eec9SHong Zhang M*/
1424b84eec9SHong Zhang /*MC
143*19c2959aSHong Zhang      TSMPRK2A23 - Second Order Partitioned Runge Kutta scheme based on RK2A.
144*19c2959aSHong Zhang 
145*19c2959aSHong Zhang      This method has eight stages for slow and medium and fast parts. The refinement factor of the stepsize is 2.
146*19c2959aSHong Zhang      r = 2, np = 3
147*19c2959aSHong Zhang      Options database:
148*19c2959aSHong Zhang .     -ts_mprk_type 2a23
149*19c2959aSHong Zhang 
150*19c2959aSHong Zhang      Level: advanced
151*19c2959aSHong Zhang 
152*19c2959aSHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType()
153*19c2959aSHong Zhang M*/
154*19c2959aSHong Zhang /*MC
155*19c2959aSHong Zhang      TSMPRK2A32 - Second Order Partitioned Runge Kutta scheme based on RK2A.
156*19c2959aSHong Zhang 
157*19c2959aSHong Zhang      This method has four stages for slow and fast parts. The refinement factor of the stepsize is 3.
158*19c2959aSHong Zhang      r = 3, np = 2
159*19c2959aSHong Zhang      Options database:
160*19c2959aSHong Zhang .     -ts_mprk_type 2a32
161*19c2959aSHong Zhang 
162*19c2959aSHong Zhang      Level: advanced
163*19c2959aSHong Zhang 
164*19c2959aSHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType()
165*19c2959aSHong Zhang M*/
166*19c2959aSHong Zhang /*MC
167*19c2959aSHong Zhang      TSMPRK2A33 - Second Order Partitioned Runge Kutta scheme based on RK2A.
168*19c2959aSHong Zhang 
169*19c2959aSHong Zhang      This method has eight stages for slow and medium and fast parts. The refinement factor of the stepsize is 3.
170*19c2959aSHong Zhang      r = 3, np = 3
171*19c2959aSHong Zhang      Options database:
172*19c2959aSHong Zhang .     -ts_mprk_type 2a33
173*19c2959aSHong Zhang 
174*19c2959aSHong Zhang      Level: advanced
175*19c2959aSHong Zhang 
176*19c2959aSHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType()
177*19c2959aSHong Zhang M*/
178*19c2959aSHong Zhang /*MC
179*19c2959aSHong Zhang      TSMPRK3P2M - Third Order Partitioned Runge Kutta scheme.
1804b84eec9SHong Zhang 
1814b84eec9SHong Zhang      This method has eight stages for both slow and fast parts.
1824b84eec9SHong Zhang 
1834b84eec9SHong Zhang      Options database:
1844b84eec9SHong Zhang .     -ts_mprk_type pm3  (put here temporarily)
1854b84eec9SHong Zhang 
1864b84eec9SHong Zhang      Level: advanced
1874b84eec9SHong Zhang 
1884b84eec9SHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType()
1894b84eec9SHong Zhang M*/
1904b84eec9SHong Zhang /*MC
1914b84eec9SHong Zhang      TSMPRKP2 - Second Order Partitioned Runge Kutta scheme.
1924b84eec9SHong Zhang 
1934b84eec9SHong Zhang      This method has five stages for both slow and fast parts.
1944b84eec9SHong Zhang 
1954b84eec9SHong Zhang      Options database:
1964b84eec9SHong Zhang .     -ts_mprk_type p2
1974b84eec9SHong Zhang 
1984b84eec9SHong Zhang      Level: advanced
1994b84eec9SHong Zhang 
2004b84eec9SHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType()
2014b84eec9SHong Zhang M*/
2024b84eec9SHong Zhang /*MC
2034b84eec9SHong Zhang      TSMPRKP3 - Third Order Partitioned Runge Kutta scheme.
2044b84eec9SHong Zhang 
2054b84eec9SHong Zhang      This method has ten stages for both slow and fast parts.
2064b84eec9SHong Zhang 
2074b84eec9SHong Zhang      Options database:
2084b84eec9SHong Zhang .     -ts_mprk_type p3
2094b84eec9SHong Zhang 
2104b84eec9SHong Zhang      Level: advanced
2114b84eec9SHong Zhang 
2124b84eec9SHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType()
2134b84eec9SHong Zhang M*/
2144b84eec9SHong Zhang 
2154b84eec9SHong Zhang /*@C
2164b84eec9SHong Zhang   TSMPRKRegisterAll - Registers all of the Partirioned Runge-Kutta explicit methods in TSMPRK
2174b84eec9SHong Zhang 
2184b84eec9SHong Zhang   Not Collective, but should be called by all processes which will need the schemes to be registered
2194b84eec9SHong Zhang 
2204b84eec9SHong Zhang   Level: advanced
2214b84eec9SHong Zhang 
2224b84eec9SHong Zhang .keywords: TS, TSMPRK, register, all
2234b84eec9SHong Zhang 
2244b84eec9SHong Zhang .seealso:  TSMPRKRegisterDestroy()
2254b84eec9SHong Zhang @*/
2264b84eec9SHong Zhang PetscErrorCode TSMPRKRegisterAll(void)
2274b84eec9SHong Zhang {
2284b84eec9SHong Zhang   PetscErrorCode ierr;
2294b84eec9SHong Zhang 
2304b84eec9SHong Zhang   PetscFunctionBegin;
2314b84eec9SHong Zhang   if (TSMPRKRegisterAllCalled) PetscFunctionReturn(0);
2324b84eec9SHong Zhang   TSMPRKRegisterAllCalled = PETSC_TRUE;
2334b84eec9SHong Zhang 
2344b84eec9SHong Zhang #define RC PetscRealConstant
2354b84eec9SHong Zhang   {
2364b84eec9SHong Zhang     const PetscReal
237*19c2959aSHong Zhang       Abase[2][2] = {{0,0},
238*19c2959aSHong Zhang                      {RC(1.0),0}},
239*19c2959aSHong Zhang       bbase[2] = {RC(0.5),RC(0.5)};
240*19c2959aSHong Zhang     PetscReal
241*19c2959aSHong Zhang       Asb[4][4] = {{0}},Af[4][4] = {{0}},bsb[4] = {0},bf[4] = {0};
242*19c2959aSHong Zhang     PetscInt
243*19c2959aSHong Zhang       rsb[4] = {0,0,1,2};
244*19c2959aSHong Zhang     ierr = TSMPRKGenerateTableau2(2,2,&Abase[0][0],bbase,&Asb[0][0],bsb,&Af[0][0],bf);CHKERRQ(ierr);
245*19c2959aSHong Zhang     ierr = TSMPRKRegister(TSMPRK2A22,2,4,&Asb[0][0],bsb,NULL,rsb,NULL,NULL,NULL,NULL,&Af[0][0],bf,NULL);CHKERRQ(ierr);
2464b84eec9SHong Zhang   }
247*19c2959aSHong Zhang   {
248*19c2959aSHong Zhang     const PetscReal
249*19c2959aSHong Zhang       Abase[2][2] = {{0,0},
250*19c2959aSHong Zhang                      {RC(1.0),0}},
251*19c2959aSHong Zhang       bbase[2]    = {RC(0.5),RC(0.5)};
252*19c2959aSHong Zhang     PetscReal
253*19c2959aSHong Zhang       Asb[8][8] = {{0}},Amb[8][8] = {{0}},Af[8][8] = {{0}},bsb[8] ={0},bmb[8] = {0},bf[8] = {0};
254*19c2959aSHong Zhang     PetscInt
255*19c2959aSHong Zhang       rsb[8] = {0,0,1,2,1,2,1,2},rmb[8] = {0,0,1,2,0,0,5,6};
256*19c2959aSHong Zhang     ierr = TSMPRKGenerateTableau3(2,2,&Abase[0][0],bbase,&Asb[0][0],bsb,&Amb[0][0],bmb,&Af[0][0],bf);CHKERRQ(ierr);
257*19c2959aSHong Zhang     ierr = TSMPRKRegister(TSMPRK2A23,2,8,&Asb[0][0],bsb,NULL,rsb,&Amb[0][0],bmb,NULL,rmb,&Af[0][0],bf,NULL);CHKERRQ(ierr);
258*19c2959aSHong Zhang   }
259*19c2959aSHong Zhang   {
260*19c2959aSHong Zhang     const PetscReal
261*19c2959aSHong Zhang       Abase[2][2] = {{0,0},
262*19c2959aSHong Zhang                      {RC(1.0),0}},
263*19c2959aSHong Zhang       bbase[2]    = {RC(0.5),RC(0.5)};
264*19c2959aSHong Zhang     PetscReal
265*19c2959aSHong Zhang       Asb[6][6] = {{0}},Af[6][6] = {{0}},bsb[6] = {0},bf[6] = {0};
266*19c2959aSHong Zhang     PetscInt
267*19c2959aSHong Zhang       rsb[6] = {0,0,1,2,1,2};
268*19c2959aSHong Zhang     ierr = TSMPRKGenerateTableau2(3,2,&Abase[0][0],bbase,&Asb[0][0],bsb,&Af[0][0],bf);CHKERRQ(ierr);
269*19c2959aSHong Zhang     ierr = TSMPRKRegister(TSMPRK2A32,2,6,&Asb[0][0],bsb,NULL,rsb,NULL,NULL,NULL,NULL,&Af[0][0],bf,NULL);CHKERRQ(ierr);
270*19c2959aSHong Zhang   }
271*19c2959aSHong Zhang   {
272*19c2959aSHong Zhang     const PetscReal
273*19c2959aSHong Zhang       Abase[2][2] = {{0,0},
274*19c2959aSHong Zhang                      {RC(1.0),0}},
275*19c2959aSHong Zhang       bbase[2]    = {RC(0.5),RC(0.5)};
276*19c2959aSHong Zhang     PetscReal
277*19c2959aSHong Zhang       Asb[18][18] = {{0}},Amb[18][18] = {{0}},Af[18][18] = {{0}},bsb[18] ={0},bmb[18] = {0},bf[18] = {0};
278*19c2959aSHong Zhang     PetscInt
279*19c2959aSHong Zhang       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};
280*19c2959aSHong Zhang     ierr = TSMPRKGenerateTableau3(3,2,&Abase[0][0],bbase,&Asb[0][0],bsb,&Amb[0][0],bmb,&Af[0][0],bf);CHKERRQ(ierr);
281*19c2959aSHong Zhang     ierr = TSMPRKRegister(TSMPRK2A33,2,18,&Asb[0][0],bsb,NULL,rsb,&Amb[0][0],bmb,NULL,rmb,&Af[0][0],bf,NULL);CHKERRQ(ierr);
282*19c2959aSHong Zhang   }
283*19c2959aSHong Zhang /*
284*19c2959aSHong Zhang     PetscReal
285*19c2959aSHong Zhang       Asb[8][8] = {{Abase[0][0],Abase[0][1],0,0,0,0,0,0},
286*19c2959aSHong Zhang                    {Abase[1][0],Abase[1][1],0,0,0,0,0,0},
287*19c2959aSHong Zhang                    {0,0,Abase[0][0],Abase[0][1],0,0,0,0},
288*19c2959aSHong Zhang                    {0,0,Abase[1][0],Abase[1][1],0,0,0,0},
289*19c2959aSHong Zhang                    {0,0,0,0,Abase[0][0],Abase[0][1],0,0},
290*19c2959aSHong Zhang                    {0,0,0,0,Abase[1][0],Abase[1][1],0,0},
291*19c2959aSHong Zhang                    {0,0,0,0,0,0,Abase[0][0],Abase[0][1]},
292*19c2959aSHong Zhang                    {0,0,0,0,0,0,Abase[1][0],Abase[1][1]}},
293*19c2959aSHong Zhang       Amb[8][8] = {{Abase[0][0]/m,Abase[0][1]/m,0,0,0,0,0,0},
294*19c2959aSHong Zhang                    {Abase[1][0]/m,Abase[1][1]/m,0,0,0,0,0,0},
295*19c2959aSHong Zhang                    {0,0,Abase[0][0]/m,Abase[0][1]/m,0,0,0,0},
296*19c2959aSHong Zhang                    {0,0,Abase[1][0]/m,Abase[1][1]/m,0,0,0,0},
297*19c2959aSHong Zhang                    {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,Abase[0][0]/m,Abase[0][1]/m,0,0},
298*19c2959aSHong Zhang                    {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,Abase[1][0]/m,Abase[1][1]/m,0,0},
299*19c2959aSHong Zhang                    {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,0,0,Abase[0][0]/m,Abase[0][1]/m},
300*19c2959aSHong Zhang                    {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,0,0,Abase[1][0]/m,Abase[1][1]/m}},
301*19c2959aSHong Zhang       Af[8][8] = {{Abase[0][0]/m/m,Abase[0][1]/m/m,0,0,0,0,0,0},
302*19c2959aSHong Zhang                    {Abase[1][0]/m/m,Abase[1][1]/m/m,0,0,0,0,0,0},
303*19c2959aSHong Zhang                    {0,0,Abase[0][0]/m/m,Abase[0][1]/m/m,0,0,0,0},
304*19c2959aSHong Zhang                    {0,0,Abase[1][0]/m/m,Abase[1][1]/m/m,0,0,0,0},
305*19c2959aSHong 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},
306*19c2959aSHong 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},
307*19c2959aSHong 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},
308*19c2959aSHong 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}},
309*19c2959aSHong 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},
310*19c2959aSHong 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},
311*19c2959aSHong 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},
312*19c2959aSHong Zhang */
3134b84eec9SHong Zhang   /*{
3144b84eec9SHong Zhang       const PetscReal
3154b84eec9SHong Zhang         As[8][8] = {{0,0,0,0,0,0,0,0},
3164b84eec9SHong Zhang                     {RC(1.0)/RC(2.0),0,0,0,0,0,0,0},
3174b84eec9SHong Zhang                     {RC(-1.0)/RC(6.0),RC(2.0)/RC(3.0),0,0,0,0,0,0},
3184b84eec9SHong Zhang                     {RC(1.0)/RC(3.0),RC(-1.0)/RC(3.0),RC(1.0),0,0,0,0,0},
3194b84eec9SHong Zhang                     {0,0,0,0,0,0,0,0},
3204b84eec9SHong Zhang                     {0,0,0,0,RC(1.0)/RC(2.0),0,0,0},
3214b84eec9SHong Zhang                     {0,0,0,0,RC(-1.0)/RC(6.0),RC(2.0)/RC(3.0),0,0},
3224b84eec9SHong Zhang                     {0,0,0,0,RC(1.0)/RC(3.0),RC(-1.0)/RC(3.0),RC(1.0),0}},
3234b84eec9SHong Zhang          A[8][8] = {{0,0,0,0,0,0,0,0},
3244b84eec9SHong Zhang                     {RC(1.0)/RC(4.0),0,0,0,0,0,0,0},
3254b84eec9SHong Zhang                     {RC(-1.0)/RC(12.0),RC(1.0)/RC(3.0),0,0,0,0,0,0},
3264b84eec9SHong Zhang                     {RC(1.0)/RC(6.0),RC(-1.0)/RC(6.0),RC(1.0)/RC(2.0),0,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),0,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(4.0),0,0,0},
3294b84eec9SHong 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},
3304b84eec9SHong 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}},
3314b84eec9SHong 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)},
3324b84eec9SHong 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)};
3334b84eec9SHong Zhang            ierr = TSMPRKRegister(TSMPRKPM3,3,8,&As[0][0],bs,NULL,&A[0][0],b,NULL);CHKERRQ(ierr);
3344b84eec9SHong Zhang   }*/
3354b84eec9SHong Zhang 
3364b84eec9SHong Zhang   {
3374b84eec9SHong Zhang     const PetscReal
338*19c2959aSHong Zhang       Asb[5][5] = {{0,0,0,0,0},
3394b84eec9SHong Zhang                    {RC(1.0)/RC(2.0),0,0,0,0},
3404b84eec9SHong Zhang                    {RC(1.0)/RC(2.0),0,0,0,0},
3414b84eec9SHong Zhang                    {RC(1.0),0,0,0,0},
3424b84eec9SHong Zhang                    {RC(1.0),0,0,0,0}},
343*19c2959aSHong Zhang       Af[5][5]  = {{0,0,0,0,0},
3444b84eec9SHong Zhang                    {RC(1.0)/RC(2.0),0,0,0,0},
3454b84eec9SHong Zhang                    {RC(1.0)/RC(4.0),RC(1.0)/RC(4.0),0,0,0},
3464b84eec9SHong Zhang                    {RC(1.0)/RC(4.0),RC(1.0)/RC(4.0),RC(1.0)/RC(2.0),0,0},
3474b84eec9SHong Zhang                    {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}},
348*19c2959aSHong Zhang       bsb[5]    = {RC(1.0)/RC(2.0),0,0,0,RC(1.0)/RC(2.0)},
349*19c2959aSHong Zhang       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};
350*19c2959aSHong Zhang     const PetscInt
351*19c2959aSHong Zhang       rsb[5]    = {0,0,2,0,4};
352*19c2959aSHong Zhang     ierr = TSMPRKRegister(TSMPRKP2,2,5,&Asb[0][0],bsb,NULL,rsb,NULL,NULL,NULL,NULL,&Af[0][0],bf,NULL);CHKERRQ(ierr);
3534b84eec9SHong Zhang   }
3544b84eec9SHong Zhang 
3554b84eec9SHong Zhang   {
3564b84eec9SHong Zhang     const PetscReal
357*19c2959aSHong Zhang       Asb[10][10] = {{0,0,0,0,0,0,0,0,0,0},
3584b84eec9SHong Zhang                      {RC(1.0)/RC(4.0),0,0,0,0,0,0,0,0,0},
3594b84eec9SHong Zhang                      {RC(1.0)/RC(4.0),0,0,0,0,0,0,0,0,0},
3604b84eec9SHong Zhang                      {RC(1.0)/RC(2.0),0,0,0,0,0,0,0,0,0},
3614b84eec9SHong Zhang                      {RC(1.0)/RC(2.0),0,0,0,0,0,0,0,0,0},
3624b84eec9SHong Zhang                      {RC(-1.0)/RC(6.0),0,0,0,RC(2.0)/RC(3.0),0,0,0,0,0},
3634b84eec9SHong 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},
3644b84eec9SHong 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},
3654b84eec9SHong Zhang                      {RC(1.0)/RC(3.0),0,0,0,RC(-1.0)/RC(3.0),RC(1.0),0,0,0,0},
3664b84eec9SHong Zhang                      {RC(1.0)/RC(3.0),0,0,0,RC(-1.0)/RC(3.0),RC(1.0),0,0,0,0}},
367*19c2959aSHong Zhang       Af[10][10]  = {{0,0,0,0,0,0,0,0,0,0},
3684b84eec9SHong Zhang                      {RC(1.0)/RC(4.0),0,0,0,0,0,0,0,0,0},
3694b84eec9SHong Zhang                      {RC(-1.0)/RC(12.0),RC(1.0)/RC(3.0),0,0,0,0,0,0,0,0},
3704b84eec9SHong Zhang                      {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},
3714b84eec9SHong 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,0,0},
3724b84eec9SHong 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,0,0},
3734b84eec9SHong 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,RC(1.0)/RC(4.0),0,0,0,0},
3744b84eec9SHong 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,RC(-1.0)/RC(12.0),RC(1.0)/RC(3.0),0,0,0},
3754b84eec9SHong 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,RC(1.0)/RC(6.0),RC(-1.0)/RC(6.0),RC(1.0)/RC(2.0),0,0},
3764b84eec9SHong 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,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}},
377*19c2959aSHong Zhang       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)},
378*19c2959aSHong Zhang       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};
379*19c2959aSHong Zhang     const PetscInt
380*19c2959aSHong Zhang       rsb[10]     = {0,0,2,0,4,0,0,7,0,9};
381*19c2959aSHong Zhang     ierr = TSMPRKRegister(TSMPRKP3,3,10,&Asb[0][0],bsb,NULL,rsb,NULL,NULL,NULL,NULL,&Af[0][0],bf,NULL);CHKERRQ(ierr);
3824b84eec9SHong Zhang   }
3834b84eec9SHong Zhang #undef RC
3844b84eec9SHong Zhang   PetscFunctionReturn(0);
3854b84eec9SHong Zhang }
3864b84eec9SHong Zhang 
3874b84eec9SHong Zhang /*@C
3884b84eec9SHong Zhang    TSMPRKRegisterDestroy - Frees the list of schemes that were registered by TSMPRKRegister().
3894b84eec9SHong Zhang 
3904b84eec9SHong Zhang    Not Collective
3914b84eec9SHong Zhang 
3924b84eec9SHong Zhang    Level: advanced
3934b84eec9SHong Zhang 
3944b84eec9SHong Zhang .keywords: TSMPRK, register, destroy
3954b84eec9SHong Zhang .seealso: TSMPRKRegister(), TSMPRKRegisterAll()
3964b84eec9SHong Zhang @*/
3974b84eec9SHong Zhang PetscErrorCode TSMPRKRegisterDestroy(void)
3984b84eec9SHong Zhang {
3994b84eec9SHong Zhang   PetscErrorCode ierr;
4004b84eec9SHong Zhang   MPRKTableauLink link;
4014b84eec9SHong Zhang 
4024b84eec9SHong Zhang   PetscFunctionBegin;
4034b84eec9SHong Zhang   while ((link = MPRKTableauList)) {
4044b84eec9SHong Zhang     MPRKTableau t = &link->tab;
4054b84eec9SHong Zhang     MPRKTableauList = link->next;
406*19c2959aSHong Zhang     ierr = PetscFree3(t->Asb,t->bsb,t->csb);CHKERRQ(ierr);
407*19c2959aSHong Zhang     ierr = PetscFree3(t->Amb,t->bmb,t->cmb);CHKERRQ(ierr);
4084b84eec9SHong Zhang     ierr = PetscFree3(t->Af,t->bf,t->cf);CHKERRQ(ierr);
409*19c2959aSHong Zhang     ierr = PetscFree2(t->rsb,t->rmb);CHKERRQ(ierr);
4104b84eec9SHong Zhang     ierr = PetscFree (t->name);CHKERRQ(ierr);
4114b84eec9SHong Zhang     ierr = PetscFree (link);CHKERRQ(ierr);
4124b84eec9SHong Zhang   }
4134b84eec9SHong Zhang   TSMPRKRegisterAllCalled = PETSC_FALSE;
4144b84eec9SHong Zhang   PetscFunctionReturn(0);
4154b84eec9SHong Zhang }
4164b84eec9SHong Zhang 
4174b84eec9SHong Zhang /*@C
4184b84eec9SHong Zhang   TSMPRKInitializePackage - This function initializes everything in the TSMPRK package. It is called
4194b84eec9SHong Zhang   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_MPRK()
4204b84eec9SHong Zhang   when using static libraries.
4214b84eec9SHong Zhang 
4224b84eec9SHong Zhang   Level: developer
4234b84eec9SHong Zhang 
4244b84eec9SHong Zhang .keywords: TS, TSMPRK, initialize, package
4254b84eec9SHong Zhang .seealso: PetscInitialize()
4264b84eec9SHong Zhang @*/
4274b84eec9SHong Zhang PetscErrorCode TSMPRKInitializePackage(void)
4284b84eec9SHong Zhang {
4294b84eec9SHong Zhang   PetscErrorCode ierr;
4304b84eec9SHong Zhang 
4314b84eec9SHong Zhang   PetscFunctionBegin;
4324b84eec9SHong Zhang   if (TSMPRKPackageInitialized) PetscFunctionReturn(0);
4334b84eec9SHong Zhang   TSMPRKPackageInitialized = PETSC_TRUE;
4344b84eec9SHong Zhang   ierr = TSMPRKRegisterAll();CHKERRQ(ierr);
4354b84eec9SHong Zhang   ierr = PetscRegisterFinalize(TSMPRKFinalizePackage);CHKERRQ(ierr);
4364b84eec9SHong Zhang   PetscFunctionReturn(0);
4374b84eec9SHong Zhang }
4384b84eec9SHong Zhang 
4394b84eec9SHong Zhang /*@C
4404b84eec9SHong Zhang   TSRKFinalizePackage - This function destroys everything in the TSMPRK package. It is
4414b84eec9SHong Zhang   called from PetscFinalize().
4424b84eec9SHong Zhang 
4434b84eec9SHong Zhang   Level: developer
4444b84eec9SHong Zhang 
4454b84eec9SHong Zhang .keywords: Petsc, destroy, package
4464b84eec9SHong Zhang .seealso: PetscFinalize()
4474b84eec9SHong Zhang @*/
4484b84eec9SHong Zhang PetscErrorCode TSMPRKFinalizePackage(void)
4494b84eec9SHong Zhang {
4504b84eec9SHong Zhang   PetscErrorCode ierr;
4514b84eec9SHong Zhang 
4524b84eec9SHong Zhang   PetscFunctionBegin;
4534b84eec9SHong Zhang   TSMPRKPackageInitialized = PETSC_FALSE;
4544b84eec9SHong Zhang   ierr = TSMPRKRegisterDestroy();CHKERRQ(ierr);
4554b84eec9SHong Zhang   PetscFunctionReturn(0);
4564b84eec9SHong Zhang }
4574b84eec9SHong Zhang 
4584b84eec9SHong Zhang /*@C
4594b84eec9SHong Zhang    TSMPRKRegister - register a MPRK scheme by providing the entries in the Butcher tableau
4604b84eec9SHong Zhang 
4614b84eec9SHong Zhang    Not Collective, but the same schemes should be registered on all processes on which they will be used
4624b84eec9SHong Zhang 
4634b84eec9SHong Zhang    Input Parameters:
4644b84eec9SHong Zhang +  name - identifier for method
4654b84eec9SHong Zhang .  order - approximation order of method
4664b84eec9SHong Zhang .  s  - number of stages, this is the dimension of the matrices below
4674b84eec9SHong Zhang .  Af - stage coefficients for fast components(dimension s*s, row-major)
4684b84eec9SHong Zhang .  bf - step completion table for fast components(dimension s)
4694b84eec9SHong Zhang .  cf - abscissa for fast components(dimension s)
4704b84eec9SHong Zhang .  As - stage coefficients for slow components(dimension s*s, row-major)
4714b84eec9SHong Zhang .  bs - step completion table for slow components(dimension s)
4724b84eec9SHong Zhang -  cs - abscissa for slow components(dimension s)
4734b84eec9SHong Zhang 
4744b84eec9SHong Zhang    Notes:
4754b84eec9SHong Zhang    Several MPRK methods are provided, this function is only needed to create new methods.
4764b84eec9SHong Zhang 
4774b84eec9SHong Zhang    Level: advanced
4784b84eec9SHong Zhang 
4794b84eec9SHong Zhang .keywords: TS, register
4804b84eec9SHong Zhang 
4814b84eec9SHong Zhang .seealso: TSMPRK
4824b84eec9SHong Zhang @*/
4834b84eec9SHong Zhang PetscErrorCode TSMPRKRegister(TSMPRKType name,PetscInt order,PetscInt s,
484*19c2959aSHong Zhang                               const PetscReal Asb[],const PetscReal bsb[],const PetscReal csb[],const PetscInt rsb[],
485*19c2959aSHong Zhang                               const PetscReal Amb[],const PetscReal bmb[],const PetscReal cmb[],const PetscInt rmb[],
4864b84eec9SHong Zhang                               const PetscReal Af[],const PetscReal bf[],const PetscReal cf[])
4874b84eec9SHong Zhang {
4884b84eec9SHong Zhang   MPRKTableauLink link;
4894b84eec9SHong Zhang   MPRKTableau     t;
4904b84eec9SHong Zhang   PetscInt        i,j;
4914b84eec9SHong Zhang   PetscErrorCode  ierr;
4924b84eec9SHong Zhang 
4934b84eec9SHong Zhang   PetscFunctionBegin;
4944b84eec9SHong Zhang   PetscValidCharPointer(name,1);
495*19c2959aSHong Zhang   PetscValidRealPointer(Asb,4);
496*19c2959aSHong Zhang   if (bsb) PetscValidRealPointer(bsb,5);
497*19c2959aSHong Zhang   if (csb) PetscValidRealPointer(csb,6);
498*19c2959aSHong Zhang   if (rsb) PetscValidRealPointer(rsb,7);
499*19c2959aSHong Zhang   if (Amb) PetscValidRealPointer(Amb,8);
500*19c2959aSHong Zhang   if (bmb) PetscValidRealPointer(bmb,9);
501*19c2959aSHong Zhang   if (cmb) PetscValidRealPointer(cmb,10);
502*19c2959aSHong Zhang   if (rmb) PetscValidRealPointer(rmb,11);
503*19c2959aSHong Zhang   PetscValidRealPointer(Af,12);
5044b84eec9SHong Zhang   if (bf) PetscValidRealPointer(bf,8);
5054b84eec9SHong Zhang   if (cf) PetscValidRealPointer(cf,9);
5064b84eec9SHong Zhang 
5074b84eec9SHong Zhang   ierr = PetscNew(&link);CHKERRQ(ierr);
5084b84eec9SHong Zhang   t = &link->tab;
5094b84eec9SHong Zhang 
5104b84eec9SHong Zhang   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
5114b84eec9SHong Zhang   t->order = order;
5124b84eec9SHong Zhang   t->s  = s;
513*19c2959aSHong Zhang   t->np = 2;
5144b84eec9SHong Zhang   ierr = PetscMalloc3(s*s,&t->Af,s,&t->bf,s,&t->cf);CHKERRQ(ierr);
5154b84eec9SHong Zhang   ierr = PetscMemcpy(t->Af,Af,s*s*sizeof(Af[0]));CHKERRQ(ierr);
5164b84eec9SHong Zhang   if (bf) {
5174b84eec9SHong Zhang     ierr = PetscMemcpy(t->bf,bf,s*sizeof(bf[0]));CHKERRQ(ierr);
518*19c2959aSHong Zhang   } else
5194b84eec9SHong Zhang     for (i=0; i<s; i++) t->bf[i] = Af[(s-1)*s+i];
5204b84eec9SHong Zhang   if (cf) {
5214b84eec9SHong Zhang     ierr = PetscMemcpy(t->cf,cf,s*sizeof(cf[0]));CHKERRQ(ierr);
522*19c2959aSHong Zhang   } else {
5234b84eec9SHong Zhang     for (i=0; i<s; i++)
5244b84eec9SHong Zhang       for (j=0,t->cf[i]=0; j<s; j++)
5254b84eec9SHong Zhang         t->cf[i] += Af[i*s+j];
5264b84eec9SHong Zhang   }
527*19c2959aSHong Zhang 
528*19c2959aSHong Zhang   if (Amb) {
529*19c2959aSHong Zhang     t->np = 3;
530*19c2959aSHong Zhang     ierr = PetscMalloc3(s*s,&t->Amb,s,&t->bmb,s,&t->cmb);CHKERRQ(ierr);
531*19c2959aSHong Zhang     ierr = PetscCalloc1(s,&t->rmb);CHKERRQ(ierr);
532*19c2959aSHong Zhang     ierr = PetscMemcpy(t->Amb,Amb,s*s*sizeof(Amb[0]));CHKERRQ(ierr);
533*19c2959aSHong Zhang     if (bmb) {
534*19c2959aSHong Zhang       ierr = PetscMemcpy(t->bmb,bmb,s*sizeof(bmb[0]));CHKERRQ(ierr);
535*19c2959aSHong Zhang     } else {
536*19c2959aSHong Zhang       for (i=0; i<s; i++) t->bmb[i] = Amb[(s-1)*s+i];
5374b84eec9SHong Zhang     }
538*19c2959aSHong Zhang     if (cmb) {
539*19c2959aSHong Zhang       ierr = PetscMemcpy(t->cmb,cmb,s*sizeof(cmb[0]));CHKERRQ(ierr);
540*19c2959aSHong Zhang     } else {
5414b84eec9SHong Zhang       for (i=0; i<s; i++)
542*19c2959aSHong Zhang         for (j=0,t->cmb[i]=0; j<s; j++)
543*19c2959aSHong Zhang           t->cmb[i] += Amb[i*s+j];
544*19c2959aSHong Zhang     }
545*19c2959aSHong Zhang     if (rmb) {
546*19c2959aSHong Zhang       ierr = PetscMemcpy(t->rmb,rmb,s*sizeof(rmb[0]));CHKERRQ(ierr);
547*19c2959aSHong Zhang     }
548*19c2959aSHong Zhang   }
549*19c2959aSHong Zhang 
550*19c2959aSHong Zhang   ierr = PetscMalloc3(s*s,&t->Asb,s,&t->bsb,s,&t->csb);CHKERRQ(ierr);
551*19c2959aSHong Zhang   ierr = PetscMemcpy(t->Asb,Asb,s*s*sizeof(Asb[0]));CHKERRQ(ierr);
552*19c2959aSHong Zhang   ierr = PetscCalloc1(s,&t->rsb);CHKERRQ(ierr);
553*19c2959aSHong Zhang   if (bsb) {
554*19c2959aSHong Zhang     ierr = PetscMemcpy(t->bsb,bsb,s*sizeof(bsb[0]));CHKERRQ(ierr);
555*19c2959aSHong Zhang   } else
556*19c2959aSHong Zhang     for (i=0; i<s; i++) t->bsb[i] = Asb[(s-1)*s+i];
557*19c2959aSHong Zhang   if (csb) {
558*19c2959aSHong Zhang     ierr = PetscMemcpy(t->csb,csb,s*sizeof(csb[0]));CHKERRQ(ierr);
559*19c2959aSHong Zhang   } else {
560*19c2959aSHong Zhang     for (i=0; i<s; i++)
561*19c2959aSHong Zhang       for (j=0,t->csb[i]=0; j<s; j++)
562*19c2959aSHong Zhang         t->csb[i] += Asb[i*s+j];
563*19c2959aSHong Zhang   }
564*19c2959aSHong Zhang   if (rsb) {
565*19c2959aSHong Zhang     ierr = PetscMemcpy(t->rsb,rsb,s*sizeof(rsb[0]));CHKERRQ(ierr);
5664b84eec9SHong Zhang   }
5674b84eec9SHong Zhang   link->next = MPRKTableauList;
5684b84eec9SHong Zhang   MPRKTableauList = link;
5694b84eec9SHong Zhang   PetscFunctionReturn(0);
5704b84eec9SHong Zhang }
5714b84eec9SHong Zhang 
5724b84eec9SHong Zhang static PetscErrorCode TSMPRKSetSplits(TS ts)
5734b84eec9SHong Zhang {
5744b84eec9SHong Zhang   TS_MPRK        *mprk = (TS_MPRK*)ts->data;
575*19c2959aSHong Zhang   MPRKTableau    tab = mprk->tableau;
5764b84eec9SHong Zhang   DM             dm,subdm,newdm;
5774b84eec9SHong Zhang   PetscErrorCode ierr;
5784b84eec9SHong Zhang 
5794b84eec9SHong Zhang   PetscFunctionBegin;
5804b84eec9SHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"slow",&mprk->subts_slow);CHKERRQ(ierr);
5814b84eec9SHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"fast",&mprk->subts_fast);CHKERRQ(ierr);
5824b84eec9SHong Zhang   if (!mprk->subts_slow || !mprk->subts_fast) SETERRQ(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");
5834b84eec9SHong Zhang 
584*19c2959aSHong Zhang   /* Only copy the DM */
5854b84eec9SHong Zhang   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
586*19c2959aSHong Zhang 
587*19c2959aSHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"slowbuffer",&mprk->subts_slowbuffer);CHKERRQ(ierr);
588*19c2959aSHong Zhang   if (!mprk->subts_slowbuffer) {
589*19c2959aSHong Zhang     mprk->subts_slowbuffer = mprk->subts_slow;
590*19c2959aSHong Zhang     mprk->subts_slow       = NULL;
591*19c2959aSHong Zhang   }
592*19c2959aSHong Zhang   if (mprk->subts_slow) {
5934b84eec9SHong Zhang     ierr = DMClone(dm,&newdm);CHKERRQ(ierr);
5944b84eec9SHong Zhang     ierr = TSGetDM(mprk->subts_slow,&subdm);CHKERRQ(ierr);
5954b84eec9SHong Zhang     ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr);
5964b84eec9SHong Zhang     ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr);
5974b84eec9SHong Zhang     ierr = TSSetDM(mprk->subts_slow,newdm);CHKERRQ(ierr);
5984b84eec9SHong Zhang     ierr = DMDestroy(&newdm);CHKERRQ(ierr);
599*19c2959aSHong Zhang   }
600*19c2959aSHong Zhang   ierr = DMClone(dm,&newdm);CHKERRQ(ierr);
601*19c2959aSHong Zhang   ierr = TSGetDM(mprk->subts_slowbuffer,&subdm);CHKERRQ(ierr);
602*19c2959aSHong Zhang   ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr);
603*19c2959aSHong Zhang   ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr);
604*19c2959aSHong Zhang   ierr = TSSetDM(mprk->subts_slowbuffer,newdm);CHKERRQ(ierr);
605*19c2959aSHong Zhang   ierr = DMDestroy(&newdm);CHKERRQ(ierr);
606*19c2959aSHong Zhang 
607*19c2959aSHong Zhang   ierr = DMClone(dm,&newdm);CHKERRQ(ierr);
608*19c2959aSHong Zhang   ierr = TSGetDM(mprk->subts_fast,&subdm);CHKERRQ(ierr);
609*19c2959aSHong Zhang   ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr);
610*19c2959aSHong Zhang   ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr);
611*19c2959aSHong Zhang   ierr = TSSetDM(mprk->subts_fast,newdm);CHKERRQ(ierr);
612*19c2959aSHong Zhang   ierr = DMDestroy(&newdm);CHKERRQ(ierr);
613*19c2959aSHong Zhang 
614*19c2959aSHong Zhang   if (tab->np == 3) {
615*19c2959aSHong Zhang     ierr = TSRHSSplitGetSubTS(ts,"medium",&mprk->subts_medium);CHKERRQ(ierr);
616*19c2959aSHong Zhang     ierr = TSRHSSplitGetSubTS(ts,"mediumbuffer",&mprk->subts_mediumbuffer);CHKERRQ(ierr);
617*19c2959aSHong Zhang     if (mprk->subts_medium && !mprk->subts_mediumbuffer) {
618*19c2959aSHong Zhang       mprk->subts_mediumbuffer = mprk->subts_medium;
619*19c2959aSHong Zhang       mprk->subts_medium       = NULL;
620*19c2959aSHong Zhang     }
621*19c2959aSHong Zhang     if (mprk->subts_medium) {
622*19c2959aSHong Zhang       ierr = DMClone(dm,&newdm);CHKERRQ(ierr);
623*19c2959aSHong Zhang       ierr = TSGetDM(mprk->subts_medium,&subdm);CHKERRQ(ierr);
624*19c2959aSHong Zhang       ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr);
625*19c2959aSHong Zhang       ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr);
626*19c2959aSHong Zhang       ierr = TSSetDM(mprk->subts_medium,newdm);CHKERRQ(ierr);
627*19c2959aSHong Zhang       ierr = DMDestroy(&newdm);CHKERRQ(ierr);
628*19c2959aSHong Zhang     }
629*19c2959aSHong Zhang     ierr = DMClone(dm,&newdm);CHKERRQ(ierr);
630*19c2959aSHong Zhang     ierr = TSGetDM(mprk->subts_mediumbuffer,&subdm);CHKERRQ(ierr);
631*19c2959aSHong Zhang     ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr);
632*19c2959aSHong Zhang     ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr);
633*19c2959aSHong Zhang     ierr = TSSetDM(mprk->subts_mediumbuffer,newdm);CHKERRQ(ierr);
634*19c2959aSHong Zhang     ierr = DMDestroy(&newdm);CHKERRQ(ierr);
635*19c2959aSHong Zhang   }
6364b84eec9SHong Zhang   PetscFunctionReturn(0);
6374b84eec9SHong Zhang }
6384b84eec9SHong Zhang 
6394b84eec9SHong Zhang /*
6404b84eec9SHong Zhang  This if for nonsplit RHS MPRK
6414b84eec9SHong Zhang  The step completion formula is
6424b84eec9SHong Zhang 
6434b84eec9SHong Zhang  x1 = x0 + h b^T YdotRHS
6444b84eec9SHong Zhang 
6454b84eec9SHong Zhang */
6464b84eec9SHong Zhang static PetscErrorCode TSEvaluateStep_MPRK(TS ts,PetscInt order,Vec X,PetscBool *done)
6474b84eec9SHong Zhang {
6484b84eec9SHong Zhang   TS_MPRK        *mprk = (TS_MPRK*)ts->data;
6494b84eec9SHong Zhang   MPRKTableau    tab = mprk->tableau;
6504b84eec9SHong Zhang   Vec            Xtmp = mprk->Ytmp,Xslow,Xfast;
6514b84eec9SHong Zhang   PetscScalar    *wf = mprk->work_fast,*ws = mprk->work_slow;
6524b84eec9SHong Zhang   PetscReal      h = ts->time_step;
6534b84eec9SHong Zhang   PetscInt       s = tab->s,j;
6544b84eec9SHong Zhang   PetscErrorCode ierr;
6554b84eec9SHong Zhang 
6564b84eec9SHong Zhang   PetscFunctionBegin;
6574b84eec9SHong Zhang   ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
6584b84eec9SHong Zhang   for (j=0; j<s; j++) wf[j] = h*tab->bf[j];
659*19c2959aSHong Zhang   for (j=0; j<s; j++) ws[j] = h*tab->bsb[j];
6604b84eec9SHong Zhang   ierr = VecCopy(X,Xtmp);CHKERRQ(ierr);
6614b84eec9SHong Zhang   ierr = VecMAXPY(Xtmp,s,ws,mprk->YdotRHS_slow);CHKERRQ(ierr);
6624b84eec9SHong Zhang   ierr = VecGetSubVector(Xtmp,mprk->is_slow,&Xslow);CHKERRQ(ierr);
6634b84eec9SHong Zhang   ierr = VecISCopy(X,mprk->is_slow,SCATTER_FORWARD,Xslow);CHKERRQ(ierr);
6644b84eec9SHong Zhang   ierr = VecRestoreSubVector(Xtmp,mprk->is_slow,&Xslow);CHKERRQ(ierr);
6654b84eec9SHong Zhang 
6664b84eec9SHong Zhang   /* Update fast part of X, note that the slow part has been changed but is simply discarded here */
6674b84eec9SHong Zhang   ierr = VecCopy(X,Xtmp);CHKERRQ(ierr);
6684b84eec9SHong Zhang   ierr = VecMAXPY(Xtmp,s,wf,mprk->YdotRHS_fast);CHKERRQ(ierr);
6694b84eec9SHong Zhang   ierr = VecGetSubVector(Xtmp,mprk->is_fast,&Xfast);CHKERRQ(ierr);
6704b84eec9SHong Zhang   ierr = VecISCopy(X,mprk->is_fast,SCATTER_FORWARD,Xfast);CHKERRQ(ierr);
6714b84eec9SHong Zhang   ierr = VecRestoreSubVector(Xtmp,mprk->is_fast,&Xfast);CHKERRQ(ierr);
6724b84eec9SHong Zhang   PetscFunctionReturn(0);
6734b84eec9SHong Zhang }
6744b84eec9SHong Zhang 
6754b84eec9SHong Zhang static PetscErrorCode TSStep_MPRK(TS ts)
6764b84eec9SHong Zhang {
6774b84eec9SHong Zhang   TS_MPRK         *mprk = (TS_MPRK*)ts->data;
6784b84eec9SHong Zhang   Vec             *Y = mprk->Y,Ytmp = mprk->Ytmp,*YdotRHS_fast = mprk->YdotRHS_fast,*YdotRHS_slow = mprk->YdotRHS_slow;
6794b84eec9SHong Zhang   Vec             Yfast,Yslow;
6804b84eec9SHong Zhang   MPRKTableau     tab = mprk->tableau;
6814b84eec9SHong Zhang   const PetscInt  s   = tab->s;
682*19c2959aSHong Zhang   const PetscReal *Af = tab->Af,*cf = tab->cf,*Asb = tab->Asb,*csb = tab->csb;
683*19c2959aSHong Zhang   PetscScalar     *wf = mprk->work_fast,*wsb = mprk->work_slowbuffer;
6844b84eec9SHong Zhang   PetscInt        i,j;
6854b84eec9SHong Zhang   PetscReal       next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step;
6864b84eec9SHong Zhang   PetscErrorCode  ierr;
6874b84eec9SHong Zhang 
6884b84eec9SHong Zhang   PetscFunctionBegin;
6894b84eec9SHong Zhang   for (i=0; i<s; i++) {
6904b84eec9SHong Zhang     mprk->stage_time = t + h*cf[i];
6914b84eec9SHong Zhang     ierr = TSPreStage(ts,mprk->stage_time);CHKERRQ(ierr);
6924b84eec9SHong Zhang 
6934b84eec9SHong Zhang     /* update the satge value for all components by slow and fast tableau respectively */
694*19c2959aSHong Zhang     for (j=0; j<i; j++) wsb[j] = h*Asb[i*s+j];
6954b84eec9SHong Zhang     ierr = VecCopy(ts->vec_sol,Ytmp);CHKERRQ(ierr);
696*19c2959aSHong Zhang     ierr = VecMAXPY(Ytmp,i,wsb,YdotRHS_slow);CHKERRQ(ierr);
6974b84eec9SHong Zhang     ierr = VecGetSubVector(Ytmp,mprk->is_slow,&Yslow);CHKERRQ(ierr);
6984b84eec9SHong Zhang     ierr = VecISCopy(Y[i],mprk->is_slow,SCATTER_FORWARD,Yslow);CHKERRQ(ierr);
6994b84eec9SHong Zhang     ierr = VecRestoreSubVector(Ytmp,mprk->is_slow,&Yslow);CHKERRQ(ierr);
7004b84eec9SHong Zhang 
7014b84eec9SHong Zhang     for (j=0; j<i; j++) wf[j] = h*Af[i*s+j];
7024b84eec9SHong Zhang     ierr = VecCopy(ts->vec_sol,Ytmp);CHKERRQ(ierr);
7034b84eec9SHong Zhang     ierr = VecMAXPY(Ytmp,i,wf,YdotRHS_fast);CHKERRQ(ierr);
7044b84eec9SHong Zhang     ierr = VecGetSubVector(Ytmp,mprk->is_fast,&Yfast);CHKERRQ(ierr);
7054b84eec9SHong Zhang     ierr = VecISCopy(Y[i],mprk->is_fast,SCATTER_FORWARD,Yfast);CHKERRQ(ierr);
7064b84eec9SHong Zhang     ierr = VecRestoreSubVector(Ytmp,mprk->is_fast,&Yfast);CHKERRQ(ierr);
7074b84eec9SHong Zhang 
7084b84eec9SHong Zhang     ierr = TSPostStage(ts,mprk->stage_time,i,Y); CHKERRQ(ierr);
7094b84eec9SHong Zhang     /* compute the stage RHS by fast and slow tableau respectively */
710*19c2959aSHong Zhang     ierr = TSComputeRHSFunction(ts,t+h*csb[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr);
7114b84eec9SHong Zhang     ierr = TSComputeRHSFunction(ts,t+h*cf[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr);
7124b84eec9SHong Zhang   }
7134b84eec9SHong Zhang   ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
7144b84eec9SHong Zhang   ts->ptime += ts->time_step;
7154b84eec9SHong Zhang   ts->time_step = next_time_step;
7164b84eec9SHong Zhang   PetscFunctionReturn(0);
7174b84eec9SHong Zhang }
7184b84eec9SHong Zhang 
7194b84eec9SHong Zhang /*
7204b84eec9SHong Zhang  This if for partitioned RHS MPRK
7214b84eec9SHong Zhang  The step completion formula is
7224b84eec9SHong Zhang 
7234b84eec9SHong Zhang  x1 = x0 + h b^T YdotRHS
7244b84eec9SHong Zhang 
7254b84eec9SHong Zhang */
7264b84eec9SHong Zhang static PetscErrorCode TSEvaluateStep_MPRKSPLIT(TS ts,PetscInt order,Vec X,PetscBool *done)
7274b84eec9SHong Zhang {
7284b84eec9SHong Zhang   TS_MPRK        *mprk = (TS_MPRK*)ts->data;
7294b84eec9SHong Zhang   MPRKTableau    tab  = mprk->tableau;
730*19c2959aSHong Zhang   Vec            Xslow,Xfast,Xslowbuffer; /* subvectors for slow and fast componets in X respectively */
731*19c2959aSHong Zhang   PetscScalar    *wf = mprk->work_fast,*ws = mprk->work_slow,*wsb = mprk->work_slowbuffer;
7324b84eec9SHong Zhang   PetscReal      h = ts->time_step;
733*19c2959aSHong Zhang   PetscInt       s = tab->s,j,basestages;
7344b84eec9SHong Zhang   PetscErrorCode ierr;
7354b84eec9SHong Zhang 
7364b84eec9SHong Zhang   PetscFunctionBegin;
7374b84eec9SHong Zhang   ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
738*19c2959aSHong Zhang 
739*19c2959aSHong Zhang   /* slow region */
740*19c2959aSHong Zhang   if (mprk->is_slow) {
741*19c2959aSHong Zhang     basestages = 0;
742*19c2959aSHong Zhang     for (j=0; j<s; j++) {
743*19c2959aSHong Zhang       if (!tab->rsb[j]) ws[tab->rsb[j]-1] += h*tab->bsb[j];
744*19c2959aSHong Zhang       else ws[basestages++] = h*tab->bsb[j];
745*19c2959aSHong Zhang     }
7464b84eec9SHong Zhang     ierr = VecGetSubVector(X,mprk->is_slow,&Xslow);CHKERRQ(ierr);
747*19c2959aSHong Zhang     ierr = VecMAXPY(Xslow,basestages,ws,mprk->YdotRHS_slow);CHKERRQ(ierr);
748*19c2959aSHong Zhang     ierr = VecRestoreSubVector(X,mprk->is_slow,&Xslow);CHKERRQ(ierr);
749*19c2959aSHong Zhang   }
750*19c2959aSHong Zhang 
751*19c2959aSHong Zhang   /* slow buffer region */
752*19c2959aSHong Zhang   for (j=0; j<s; j++) wsb[j] = h*tab->bsb[j];
753*19c2959aSHong Zhang   ierr = VecGetSubVector(X,mprk->is_slowbuffer,&Xslowbuffer);CHKERRQ(ierr);
754*19c2959aSHong Zhang   ierr = VecMAXPY(Xslowbuffer,s,wsb,mprk->YdotRHS_slowbuffer);CHKERRQ(ierr);
755*19c2959aSHong Zhang   ierr = VecRestoreSubVector(X,mprk->is_slowbuffer,&Xslowbuffer);CHKERRQ(ierr);
756*19c2959aSHong Zhang 
757*19c2959aSHong Zhang   if (tab->np == 3) {
758*19c2959aSHong Zhang     Vec         Xmedium,Xmediumbuffer;
759*19c2959aSHong Zhang     PetscScalar *wm = mprk->work_medium,*wmb = mprk->work_mediumbuffer;
760*19c2959aSHong Zhang     /* medium region */
761*19c2959aSHong Zhang     if (mprk->is_medium) {
762*19c2959aSHong Zhang       basestages = 0;
763*19c2959aSHong Zhang       for (j=0; j<s; j++) {
764*19c2959aSHong Zhang         if (!tab->rmb[j]) wm[tab->rmb[j]-1] += h*tab->bmb[j];
765*19c2959aSHong Zhang         else wm[basestages++] = h*tab->bmb[j];
766*19c2959aSHong Zhang       }
767*19c2959aSHong Zhang       ierr = VecGetSubVector(X,mprk->is_medium,&Xmedium);CHKERRQ(ierr);
768*19c2959aSHong Zhang       ierr = VecMAXPY(Xmedium,basestages,wm,mprk->YdotRHS_medium);CHKERRQ(ierr);
769*19c2959aSHong Zhang       ierr = VecRestoreSubVector(X,mprk->is_medium,&Xmedium);CHKERRQ(ierr);
770*19c2959aSHong Zhang     }
771*19c2959aSHong Zhang     /* medium buffer region */
772*19c2959aSHong Zhang     for (j=0; j<s; j++) wmb[j] = h*tab->bmb[j];
773*19c2959aSHong Zhang     ierr = VecGetSubVector(X,mprk->is_mediumbuffer,&Xmediumbuffer);CHKERRQ(ierr);
774*19c2959aSHong Zhang     ierr = VecMAXPY(Xmediumbuffer,s,wmb,mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr);
775*19c2959aSHong Zhang     ierr = VecRestoreSubVector(X,mprk->is_mediumbuffer,&Xmediumbuffer);CHKERRQ(ierr);
776*19c2959aSHong Zhang   }
777*19c2959aSHong Zhang 
778*19c2959aSHong Zhang   /* fast region */
779*19c2959aSHong Zhang   for (j=0; j<s; j++) wf[j] = h*tab->bf[j];
7804b84eec9SHong Zhang   ierr = VecGetSubVector(X,mprk->is_fast,&Xfast);CHKERRQ(ierr);
7814b84eec9SHong Zhang   ierr = VecMAXPY(Xfast,s,wf,mprk->YdotRHS_fast);CHKERRQ(ierr);
782*19c2959aSHong Zhang   ierr = VecRestoreSubVector(X,mprk->is_fast,&Xfast);CHKERRQ(ierr);
7834b84eec9SHong Zhang   PetscFunctionReturn(0);
7844b84eec9SHong Zhang }
7854b84eec9SHong Zhang 
7864b84eec9SHong Zhang static PetscErrorCode TSStep_MPRKSPLIT(TS ts)
7874b84eec9SHong Zhang {
7884b84eec9SHong Zhang   TS_MPRK         *mprk = (TS_MPRK*)ts->data;
7894b84eec9SHong Zhang   MPRKTableau     tab = mprk->tableau;
790*19c2959aSHong Zhang   Vec             *Y = mprk->Y,*YdotRHS_fast = mprk->YdotRHS_fast,*YdotRHS_slow = mprk->YdotRHS_slow,*YdotRHS_slowbuffer = mprk->YdotRHS_slowbuffer;
791*19c2959aSHong Zhang   Vec             Yslow,Yslowbuffer,Yfast; /* subvectors for slow and fast components in Y[i] respectively */
792*19c2959aSHong Zhang   PetscInt        s = tab->s;
793*19c2959aSHong Zhang   const PetscReal *Af = tab->Af,*cf = tab->cf,*Asb = tab->Asb,*csb = tab->csb;
794*19c2959aSHong Zhang   PetscScalar     *wf = mprk->work_fast,*ws = mprk->work_slow,*wsb = mprk->work_slowbuffer;
795*19c2959aSHong Zhang   PetscInt        i,j,basestages;
7964b84eec9SHong Zhang   PetscReal       next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step;
7974b84eec9SHong Zhang   PetscErrorCode  ierr;
7984b84eec9SHong Zhang 
7994b84eec9SHong Zhang   PetscFunctionBegin;
8004b84eec9SHong Zhang   for (i=0; i<s; i++) {
8014b84eec9SHong Zhang     mprk->stage_time = t + h*cf[i];
8024b84eec9SHong Zhang     ierr = TSPreStage(ts,mprk->stage_time);CHKERRQ(ierr);
8034b84eec9SHong Zhang     /* calculate the stage value for fast and slow components respectively */
8044b84eec9SHong Zhang     ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
805*19c2959aSHong Zhang     for (j=0; j<i; j++) wsb[j] = h*Asb[i*s+j];
806*19c2959aSHong Zhang 
807*19c2959aSHong Zhang     /* slow buffer region */
808*19c2959aSHong Zhang     ierr = VecGetSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr);
809*19c2959aSHong Zhang     ierr = VecMAXPY(Yslowbuffer,i,wsb,YdotRHS_slowbuffer);CHKERRQ(ierr);
810*19c2959aSHong Zhang     ierr = VecRestoreSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr);
811*19c2959aSHong Zhang     /* slow region */
812*19c2959aSHong Zhang     if (!tab->rsb[i] && mprk->is_slow) { /* not a repeated stage */
813*19c2959aSHong Zhang       basestages = 0;
814*19c2959aSHong Zhang       for (j=0; j<i; j++) {
815*19c2959aSHong Zhang         if (!tab->rsb[j]) ws[tab->rsb[j]-1] += h*Asb[i*s+j];
816*19c2959aSHong Zhang         else ws[basestages++] = wsb[j];
817*19c2959aSHong Zhang       }
8184b84eec9SHong Zhang       ierr = VecGetSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr);
8194b84eec9SHong Zhang       ierr = VecMAXPY(Yslow,i,ws,YdotRHS_slow);CHKERRQ(ierr);
8204b84eec9SHong Zhang       ierr = VecRestoreSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr);
821*19c2959aSHong Zhang       /* only depends on the slow buffer region */
822*19c2959aSHong Zhang       ierr = TSComputeRHSFunction(mprk->subts_slow,t+h*csb[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr);
823*19c2959aSHong Zhang     }
824*19c2959aSHong Zhang 
825*19c2959aSHong Zhang     /* fast region */
826*19c2959aSHong Zhang     for (j=0; j<i; j++) wf[j] = h*Af[i*s+j];
827*19c2959aSHong Zhang     ierr = VecGetSubVector(Y[i],mprk->is_fast,&Yfast);CHKERRQ(ierr);
828*19c2959aSHong Zhang     ierr = VecMAXPY(Yfast,i,wf,YdotRHS_fast);CHKERRQ(ierr);
8294b84eec9SHong Zhang     ierr = VecRestoreSubVector(Y[i],mprk->is_fast,&Yfast);CHKERRQ(ierr);
830*19c2959aSHong Zhang 
831*19c2959aSHong Zhang     if (tab->np == 3) {
832*19c2959aSHong Zhang       Vec *YdotRHS_medium = mprk->YdotRHS_medium,*YdotRHS_mediumbuffer = mprk->YdotRHS_mediumbuffer;
833*19c2959aSHong Zhang       Vec Ymedium,Ymediumbuffer;
834*19c2959aSHong Zhang       const PetscReal *Amb = tab->Amb,*cmb = tab->cmb;
835*19c2959aSHong Zhang       PetscScalar *wm = mprk->work_medium,*wmb = mprk->work_mediumbuffer;
836*19c2959aSHong Zhang 
837*19c2959aSHong Zhang       for (j=0; j<i; j++) wmb[j] = h*Amb[i*s+j];
838*19c2959aSHong Zhang       /* medium buffer region */
839*19c2959aSHong Zhang       ierr = VecGetSubVector(Y[i],mprk->is_mediumbuffer,&Ymediumbuffer);CHKERRQ(ierr);
840*19c2959aSHong Zhang       ierr = VecMAXPY(Ymediumbuffer,i,wmb,YdotRHS_mediumbuffer);CHKERRQ(ierr);
841*19c2959aSHong Zhang       ierr = VecRestoreSubVector(Y[i],mprk->is_mediumbuffer,&Ymediumbuffer);CHKERRQ(ierr);
842*19c2959aSHong Zhang       /* medium region */
843*19c2959aSHong Zhang       if (!tab->rmb[i] && mprk->is_medium) { /* not a repeated stage */
844*19c2959aSHong Zhang         basestages = 0;
845*19c2959aSHong Zhang         for (j=0; j<i; j++) {
846*19c2959aSHong Zhang           if (tab->rmb[j]) wm[tab->rmb[j]-1] += h*Amb[i*s+j];
847*19c2959aSHong Zhang           else wm[basestages++] = wmb[j];
848*19c2959aSHong Zhang         }
849*19c2959aSHong Zhang         ierr = VecGetSubVector(Y[basestages],mprk->is_medium,&Ymedium);CHKERRQ(ierr);
850*19c2959aSHong Zhang         ierr = VecMAXPY(Ymedium,i,wm,YdotRHS_medium);CHKERRQ(ierr);
851*19c2959aSHong Zhang         ierr = VecRestoreSubVector(Y[i],mprk->is_medium,&Ymedium);CHKERRQ(ierr);
852*19c2959aSHong Zhang         /* only depends on the medium buffer region and the slow buffer region */
853*19c2959aSHong Zhang         ierr = TSComputeRHSFunction(mprk->subts_medium,t+h*cmb[i],Y[i],YdotRHS_medium[i]);CHKERRQ(ierr);
854*19c2959aSHong Zhang       }
855*19c2959aSHong Zhang       /* must be computed after fast region and slow region are updated in Y */
856*19c2959aSHong Zhang       ierr = TSComputeRHSFunction(mprk->subts_mediumbuffer,t+h*cmb[i],Y[i],YdotRHS_mediumbuffer[i]);CHKERRQ(ierr);
857*19c2959aSHong Zhang     }
858*19c2959aSHong Zhang     /* must be computed after all regions are updated in Y */
859*19c2959aSHong Zhang     ierr = TSComputeRHSFunction(mprk->subts_slowbuffer,t+h*csb[i],Y[i],YdotRHS_slowbuffer[i]);CHKERRQ(ierr);
8604b84eec9SHong Zhang     ierr = TSComputeRHSFunction(mprk->subts_fast,t+h*cf[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr);
8614b84eec9SHong Zhang   }
8624b84eec9SHong Zhang   ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
8634b84eec9SHong Zhang   ts->ptime += ts->time_step;
8644b84eec9SHong Zhang   ts->time_step = next_time_step;
8654b84eec9SHong Zhang   PetscFunctionReturn(0);
8664b84eec9SHong Zhang }
8674b84eec9SHong Zhang 
8684b84eec9SHong Zhang static PetscErrorCode TSMPRKTableauReset(TS ts)
8694b84eec9SHong Zhang {
8704b84eec9SHong Zhang   TS_MPRK        *mprk = (TS_MPRK*)ts->data;
8714b84eec9SHong Zhang   MPRKTableau    tab = mprk->tableau;
8724b84eec9SHong Zhang   PetscErrorCode ierr;
8734b84eec9SHong Zhang 
8744b84eec9SHong Zhang   PetscFunctionBegin;
8754b84eec9SHong Zhang   if (!tab) PetscFunctionReturn(0);
8764b84eec9SHong Zhang   ierr = PetscFree(mprk->work_fast);CHKERRQ(ierr);
8774b84eec9SHong Zhang   ierr = PetscFree(mprk->work_slow);CHKERRQ(ierr);
8784b84eec9SHong Zhang   ierr = VecDestroyVecs(tab->s,&mprk->Y);CHKERRQ(ierr);
8794b84eec9SHong Zhang   if (mprk->mprkmtype == TSMPRKNONSPLIT) {
8804b84eec9SHong Zhang     ierr = VecDestroy(&mprk->Ytmp);CHKERRQ(ierr);
8814b84eec9SHong Zhang   }
8824b84eec9SHong Zhang   ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_fast);CHKERRQ(ierr);
8834b84eec9SHong Zhang   ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_slow);CHKERRQ(ierr);
8844b84eec9SHong Zhang   PetscFunctionReturn(0);
8854b84eec9SHong Zhang }
8864b84eec9SHong Zhang 
8874b84eec9SHong Zhang static PetscErrorCode TSReset_MPRK(TS ts)
8884b84eec9SHong Zhang {
8894b84eec9SHong Zhang   PetscErrorCode ierr;
8904b84eec9SHong Zhang 
8914b84eec9SHong Zhang   PetscFunctionBegin;
8924b84eec9SHong Zhang   ierr = TSMPRKTableauReset(ts);CHKERRQ(ierr);
8934b84eec9SHong Zhang   PetscFunctionReturn(0);
8944b84eec9SHong Zhang }
8954b84eec9SHong Zhang 
8964b84eec9SHong Zhang static PetscErrorCode DMCoarsenHook_TSMPRK(DM fine,DM coarse,void *ctx)
8974b84eec9SHong Zhang {
8984b84eec9SHong Zhang   PetscFunctionBegin;
8994b84eec9SHong Zhang   PetscFunctionReturn(0);
9004b84eec9SHong Zhang }
9014b84eec9SHong Zhang 
9024b84eec9SHong Zhang static PetscErrorCode DMRestrictHook_TSMPRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
9034b84eec9SHong Zhang {
9044b84eec9SHong Zhang   PetscFunctionBegin;
9054b84eec9SHong Zhang   PetscFunctionReturn(0);
9064b84eec9SHong Zhang }
9074b84eec9SHong Zhang 
9084b84eec9SHong Zhang static PetscErrorCode DMSubDomainHook_TSMPRK(DM dm,DM subdm,void *ctx)
9094b84eec9SHong Zhang {
9104b84eec9SHong Zhang   PetscFunctionBegin;
9114b84eec9SHong Zhang   PetscFunctionReturn(0);
9124b84eec9SHong Zhang }
9134b84eec9SHong Zhang 
9144b84eec9SHong Zhang static PetscErrorCode DMSubDomainRestrictHook_TSMPRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
9154b84eec9SHong Zhang {
9164b84eec9SHong Zhang   PetscFunctionBegin;
9174b84eec9SHong Zhang   PetscFunctionReturn(0);
9184b84eec9SHong Zhang }
9194b84eec9SHong Zhang 
9204b84eec9SHong Zhang static PetscErrorCode TSMPRKTableauSetUp(TS ts)
9214b84eec9SHong Zhang {
9224b84eec9SHong Zhang   TS_MPRK        *mprk  = (TS_MPRK*)ts->data;
9234b84eec9SHong Zhang   MPRKTableau    tab = mprk->tableau;
924*19c2959aSHong Zhang   Vec            YdotRHS_slow,YdotRHS_slowbuffer,YdotRHS_medium,YdotRHS_mediumbuffer,YdotRHS_fast;
9254b84eec9SHong Zhang   PetscErrorCode ierr;
9264b84eec9SHong Zhang 
9274b84eec9SHong Zhang   PetscFunctionBegin;
9284b84eec9SHong Zhang   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&mprk->Y);CHKERRQ(ierr);
9294b84eec9SHong Zhang   if (mprk->mprkmtype == TSMPRKNONSPLIT) {
9304b84eec9SHong Zhang     ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&mprk->YdotRHS_slow);CHKERRQ(ierr);
9314b84eec9SHong Zhang     ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&mprk->YdotRHS_fast);CHKERRQ(ierr);
9324b84eec9SHong Zhang     ierr = VecDuplicate(ts->vec_sol,&mprk->Ytmp);CHKERRQ(ierr);
933*19c2959aSHong Zhang     ierr = PetscMalloc1(tab->s,&mprk->work_slow);CHKERRQ(ierr);
934*19c2959aSHong Zhang     ierr = PetscMalloc1(tab->s,&mprk->work_fast);CHKERRQ(ierr);
9354b84eec9SHong Zhang   }
9364b84eec9SHong Zhang   if (mprk->mprkmtype == TSMPRKSPLIT) {
937*19c2959aSHong Zhang     if (mprk->is_slow) {
9384b84eec9SHong Zhang       ierr = VecGetSubVector(ts->vec_sol,mprk->is_slow,&YdotRHS_slow);CHKERRQ(ierr);
9394b84eec9SHong Zhang       ierr = VecDuplicateVecs(YdotRHS_slow,tab->s,&mprk->YdotRHS_slow);CHKERRQ(ierr);
9404b84eec9SHong Zhang       ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_slow,&YdotRHS_slow);CHKERRQ(ierr);
941*19c2959aSHong Zhang       ierr = PetscMalloc1(tab->s,&mprk->work_slow);CHKERRQ(ierr);
942*19c2959aSHong Zhang     }
943*19c2959aSHong Zhang     ierr = VecGetSubVector(ts->vec_sol,mprk->is_slowbuffer,&YdotRHS_slowbuffer);CHKERRQ(ierr);
944*19c2959aSHong Zhang     ierr = VecDuplicateVecs(YdotRHS_slowbuffer,tab->s,&mprk->YdotRHS_slowbuffer);CHKERRQ(ierr);
945*19c2959aSHong Zhang     ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_slowbuffer,&YdotRHS_slowbuffer);CHKERRQ(ierr);
946*19c2959aSHong Zhang     ierr = PetscMalloc1(tab->s,&mprk->work_slowbuffer);CHKERRQ(ierr);
947*19c2959aSHong Zhang 
948*19c2959aSHong Zhang     if (tab->np == 3) {
949*19c2959aSHong Zhang       if (mprk->is_medium) {
950*19c2959aSHong Zhang         ierr = VecGetSubVector(ts->vec_sol,mprk->is_medium,&YdotRHS_medium);CHKERRQ(ierr);
951*19c2959aSHong Zhang         ierr = VecDuplicateVecs(YdotRHS_medium,tab->s,&mprk->YdotRHS_medium);CHKERRQ(ierr);
952*19c2959aSHong Zhang         ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_medium,&YdotRHS_medium);CHKERRQ(ierr);
953*19c2959aSHong Zhang         ierr = PetscMalloc1(tab->s,&mprk->work_medium);CHKERRQ(ierr);
954*19c2959aSHong Zhang       }
955*19c2959aSHong Zhang       ierr = VecGetSubVector(ts->vec_sol,mprk->is_mediumbuffer,&YdotRHS_mediumbuffer);CHKERRQ(ierr);
956*19c2959aSHong Zhang       ierr = VecDuplicateVecs(YdotRHS_mediumbuffer,tab->s,&mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr);
957*19c2959aSHong Zhang       ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_mediumbuffer,&YdotRHS_mediumbuffer);CHKERRQ(ierr);
958*19c2959aSHong Zhang       ierr = PetscMalloc1(tab->s,&mprk->work_mediumbuffer);CHKERRQ(ierr);
959*19c2959aSHong Zhang     }
960*19c2959aSHong Zhang     ierr = VecGetSubVector(ts->vec_sol,mprk->is_fast,&YdotRHS_fast);CHKERRQ(ierr);
961*19c2959aSHong Zhang     ierr = VecDuplicateVecs(YdotRHS_fast,tab->s,&mprk->YdotRHS_fast);CHKERRQ(ierr);
9624b84eec9SHong Zhang     ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_fast,&YdotRHS_fast);CHKERRQ(ierr);
963*19c2959aSHong Zhang     ierr = PetscMalloc1(tab->s,&mprk->work_fast);CHKERRQ(ierr);
9644b84eec9SHong Zhang   }
9654b84eec9SHong Zhang   PetscFunctionReturn(0);
9664b84eec9SHong Zhang }
9674b84eec9SHong Zhang 
9684b84eec9SHong Zhang static PetscErrorCode TSSetUp_MPRK(TS ts)
9694b84eec9SHong Zhang {
9704b84eec9SHong Zhang   TS_MPRK        *mprk = (TS_MPRK*)ts->data;
971*19c2959aSHong Zhang   MPRKTableau    tab = mprk->tableau;
9724b84eec9SHong Zhang   DM             dm;
9734b84eec9SHong Zhang   PetscErrorCode ierr;
9744b84eec9SHong Zhang 
9754b84eec9SHong Zhang   PetscFunctionBegin;
9764b84eec9SHong Zhang   ierr = TSRHSSplitGetIS(ts,"slow",&mprk->is_slow);CHKERRQ(ierr);
9774b84eec9SHong Zhang   ierr = TSRHSSplitGetIS(ts,"fast",&mprk->is_fast);CHKERRQ(ierr);
978*19c2959aSHong Zhang   if (!mprk->is_slow || !mprk->is_fast) SETERRQ1(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);
979*19c2959aSHong Zhang 
980*19c2959aSHong Zhang   if (tab->np == 3) {
981*19c2959aSHong Zhang     ierr = TSRHSSplitGetIS(ts,"medium",&mprk->is_medium);CHKERRQ(ierr);
982*19c2959aSHong Zhang     if (!mprk->is_medium) SETERRQ1(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);
983*19c2959aSHong Zhang     ierr = TSRHSSplitGetIS(ts,"mediumbuffer",&mprk->is_mediumbuffer);CHKERRQ(ierr);
984*19c2959aSHong Zhang     if (!mprk->is_mediumbuffer) { /* let medium buffer cover whole medium region */
985*19c2959aSHong Zhang       mprk->is_mediumbuffer = mprk->is_medium;
986*19c2959aSHong Zhang       mprk->is_medium = NULL;
987*19c2959aSHong Zhang     }
988*19c2959aSHong Zhang   }
989*19c2959aSHong Zhang 
990*19c2959aSHong Zhang   /* If users do not provide buffer region settings, the solver will do them automatically, but with a performance penalty */
991*19c2959aSHong Zhang   ierr = TSRHSSplitGetIS(ts,"slowbuffer",&mprk->is_slowbuffer);CHKERRQ(ierr);
992*19c2959aSHong Zhang   if (!mprk->is_slowbuffer) { /* let slow buffer cover whole slow region */
993*19c2959aSHong Zhang     mprk->is_slowbuffer = mprk->is_slow;
994*19c2959aSHong Zhang     mprk->is_slow = NULL;
995*19c2959aSHong Zhang   }
996*19c2959aSHong Zhang /*
997*19c2959aSHong Zhang   if (!mprk->is_medium) {
998*19c2959aSHong Zhang     mprk->is_medium = mprk->is_fast;
999*19c2959aSHong Zhang     mprk->is_fast = NULL;
1000*19c2959aSHong Zhang   } else {
1001*19c2959aSHong Zhang     ierr = TSRHSSplitGetIS(ts,"mediumbuffer",&mprk->is_mediumbuffer);CHKERRQ(ierr);
1002*19c2959aSHong Zhang   }
1003*19c2959aSHong Zhang */
10044b84eec9SHong Zhang   ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr);
10054b84eec9SHong Zhang   ierr = TSMPRKTableauSetUp(ts);CHKERRQ(ierr);
10064b84eec9SHong Zhang   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
10074b84eec9SHong Zhang   ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSMPRK,DMRestrictHook_TSMPRK,ts);CHKERRQ(ierr);
10084b84eec9SHong Zhang   ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSMPRK,DMSubDomainRestrictHook_TSMPRK,ts);CHKERRQ(ierr);
10094b84eec9SHong Zhang   ierr = PetscTryMethod(ts,"TSMPRKSetSplits_C",(TS),(ts));CHKERRQ(ierr);
10104b84eec9SHong Zhang   PetscFunctionReturn(0);
10114b84eec9SHong Zhang }
10124b84eec9SHong Zhang 
10134b84eec9SHong Zhang /* construct a database to chose nonsplit RHS mutirate mprk method or split RHS MPRK method */
10144b84eec9SHong Zhang const char *const TSMPRKMultirateTypes[] = {"NONSPLIT","SPLIT","TSMPRKMultirateType","TSMPRK",0};
10154b84eec9SHong Zhang 
10164b84eec9SHong Zhang static PetscErrorCode TSSetFromOptions_MPRK(PetscOptionItems *PetscOptionsObject,TS ts)
10174b84eec9SHong Zhang {
10184b84eec9SHong Zhang   TS_MPRK        *mprk = (TS_MPRK*)ts->data;
10194b84eec9SHong Zhang   PetscErrorCode ierr;
10204b84eec9SHong Zhang 
10214b84eec9SHong Zhang   PetscFunctionBegin;
10224b84eec9SHong Zhang   ierr = PetscOptionsHead(PetscOptionsObject,"PRK ODE solver options");CHKERRQ(ierr);
10234b84eec9SHong Zhang   {
10244b84eec9SHong Zhang     MPRKTableauLink link;
10254b84eec9SHong Zhang     PetscInt        count,choice;
10264b84eec9SHong Zhang     PetscBool       flg;
10274b84eec9SHong Zhang     const char      **namelist;
10284b84eec9SHong Zhang     PetscInt        mprkmtype = 0;
10294b84eec9SHong Zhang     for (link=MPRKTableauList,count=0; link; link=link->next,count++) ;
10304b84eec9SHong Zhang     ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr);
10314b84eec9SHong Zhang     for (link=MPRKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
10324b84eec9SHong Zhang     ierr = PetscOptionsEList("-ts_mprk_type","Family of MPRK method","TSMPRKSetType",(const char*const*)namelist,count,mprk->tableau->name,&choice,&flg);CHKERRQ(ierr);
10334b84eec9SHong Zhang     if (flg) {ierr = TSMPRKSetType(ts,namelist[choice]);CHKERRQ(ierr);}
10344b84eec9SHong Zhang     ierr = PetscFree(namelist);CHKERRQ(ierr);
10354b84eec9SHong Zhang     ierr = PetscOptionsEList("-ts_mprk_multirate_type","Use Combined RHS Multirate or Partioned RHS Multirat MPRK method","TSMPRKSetMultirateType",TSMPRKMultirateTypes,2,TSMPRKMultirateTypes[0],&mprkmtype,&flg);CHKERRQ(ierr);
10364b84eec9SHong Zhang      if (flg) {ierr = TSMPRKSetMultirateType(ts,mprkmtype);CHKERRQ(ierr);}
10374b84eec9SHong Zhang   }
10384b84eec9SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
10394b84eec9SHong Zhang   PetscFunctionReturn(0);
10404b84eec9SHong Zhang }
10414b84eec9SHong Zhang 
10424b84eec9SHong Zhang static PetscErrorCode TSView_MPRK(TS ts,PetscViewer viewer)
10434b84eec9SHong Zhang {
10444b84eec9SHong Zhang   TS_MPRK        *mprk = (TS_MPRK*)ts->data;
10454b84eec9SHong Zhang   PetscBool      iascii;
10464b84eec9SHong Zhang   PetscErrorCode ierr;
10474b84eec9SHong Zhang 
10484b84eec9SHong Zhang   PetscFunctionBegin;
10494b84eec9SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
10504b84eec9SHong Zhang   if (iascii) {
10514b84eec9SHong Zhang     MPRKTableau tab  = mprk->tableau;
10524b84eec9SHong Zhang     TSMPRKType  mprktype;
10534b84eec9SHong Zhang     char        fbuf[512];
10544b84eec9SHong Zhang     char        sbuf[512];
1055*19c2959aSHong Zhang     PetscInt    i;
10564b84eec9SHong Zhang     ierr = TSMPRKGetType(ts,&mprktype);CHKERRQ(ierr);
10574b84eec9SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  MPRK type %s\n",mprktype);CHKERRQ(ierr);
10584b84eec9SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Order: %D\n",tab->order);CHKERRQ(ierr);
1059*19c2959aSHong Zhang 
10604b84eec9SHong Zhang     ierr = PetscFormatRealArray(fbuf,sizeof(fbuf),"% 8.6f",tab->s,tab->cf);CHKERRQ(ierr);
10614b84eec9SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa cf = %s\n",fbuf);CHKERRQ(ierr);
1062*19c2959aSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Af = \n");CHKERRQ(ierr);
1063*19c2959aSHong Zhang     for (i=0; i<tab->s; i++) {
1064*19c2959aSHong Zhang       ierr = PetscFormatRealArray(fbuf,sizeof(fbuf),"% 8.6f",tab->s,&tab->Af[i*tab->s]);CHKERRQ(ierr);
1065*19c2959aSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"    %s\n",fbuf);CHKERRQ(ierr);
1066*19c2959aSHong Zhang     }
1067*19c2959aSHong Zhang     ierr = PetscFormatRealArray(fbuf,sizeof(fbuf),"% 8.6f",tab->s,tab->bf);CHKERRQ(ierr);
1068*19c2959aSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  bf = %s\n",fbuf);CHKERRQ(ierr);
1069*19c2959aSHong Zhang 
1070*19c2959aSHong Zhang     ierr = PetscFormatRealArray(sbuf,sizeof(sbuf),"% 8.6f",tab->s,tab->csb);CHKERRQ(ierr);
1071*19c2959aSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa csb = %s\n",sbuf);CHKERRQ(ierr);
1072*19c2959aSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Asb = \n");CHKERRQ(ierr);
1073*19c2959aSHong Zhang     for (i=0; i<tab->s; i++) {
1074*19c2959aSHong Zhang       ierr = PetscFormatRealArray(sbuf,sizeof(sbuf),"% 8.6f",tab->s,&tab->Asb[i*tab->s]);CHKERRQ(ierr);
1075*19c2959aSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"    %s\n",sbuf);CHKERRQ(ierr);
1076*19c2959aSHong Zhang     }
1077*19c2959aSHong Zhang     ierr = PetscFormatRealArray(sbuf,sizeof(sbuf),"% 8.6f",tab->s,tab->bsb);CHKERRQ(ierr);
1078*19c2959aSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  bsb = %s\n",sbuf);CHKERRQ(ierr);
1079*19c2959aSHong Zhang 
1080*19c2959aSHong Zhang     if (tab->np == 3) {
1081*19c2959aSHong Zhang       char mbuf[512];
1082*19c2959aSHong Zhang       ierr = PetscFormatRealArray(mbuf,sizeof(mbuf),"% 8.6f",tab->s,tab->cmb);CHKERRQ(ierr);
1083*19c2959aSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa cmb = %s\n",mbuf);CHKERRQ(ierr);
1084*19c2959aSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  Amb = \n");CHKERRQ(ierr);
1085*19c2959aSHong Zhang       for (i=0; i<tab->s; i++) {
1086*19c2959aSHong Zhang         ierr = PetscFormatRealArray(mbuf,sizeof(mbuf),"% 8.6f",tab->s,&tab->Amb[i*tab->s]);CHKERRQ(ierr);
1087*19c2959aSHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    %s\n",mbuf);CHKERRQ(ierr);
1088*19c2959aSHong Zhang       }
1089*19c2959aSHong Zhang       ierr = PetscFormatRealArray(mbuf,sizeof(mbuf),"% 8.6f",tab->s,tab->bmb);CHKERRQ(ierr);
1090*19c2959aSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  bmb = %s\n",mbuf);CHKERRQ(ierr);
1091*19c2959aSHong Zhang     }
10924b84eec9SHong Zhang   }
10934b84eec9SHong Zhang   PetscFunctionReturn(0);
10944b84eec9SHong Zhang }
10954b84eec9SHong Zhang 
10964b84eec9SHong Zhang static PetscErrorCode TSLoad_MPRK(TS ts,PetscViewer viewer)
10974b84eec9SHong Zhang {
10984b84eec9SHong Zhang   PetscErrorCode ierr;
10994b84eec9SHong Zhang   TSAdapt        adapt;
11004b84eec9SHong Zhang 
11014b84eec9SHong Zhang   PetscFunctionBegin;
11024b84eec9SHong Zhang   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
11034b84eec9SHong Zhang   ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr);
11044b84eec9SHong Zhang   PetscFunctionReturn(0);
11054b84eec9SHong Zhang }
11064b84eec9SHong Zhang 
11074b84eec9SHong Zhang /*@C
11084b84eec9SHong Zhang   TSMPRKSetType - Set the type of MPRK scheme
11094b84eec9SHong Zhang 
11104b84eec9SHong Zhang   Logically collective
11114b84eec9SHong Zhang 
11124b84eec9SHong Zhang   Input Parameter:
11134b84eec9SHong Zhang +  ts - timestepping context
11144b84eec9SHong Zhang -  mprktype - type of MPRK-scheme
11154b84eec9SHong Zhang 
11164b84eec9SHong Zhang   Options Database:
11174b84eec9SHong Zhang .   -ts_mprk_type - <pm2,p2,p3>
11184b84eec9SHong Zhang 
11194b84eec9SHong Zhang   Level: intermediate
11204b84eec9SHong Zhang 
11214b84eec9SHong Zhang .seealso: TSMPRKGetType(), TSMPRK, TSMPRKType
11224b84eec9SHong Zhang @*/
11234b84eec9SHong Zhang PetscErrorCode TSMPRKSetType(TS ts,TSMPRKType mprktype)
11244b84eec9SHong Zhang {
11254b84eec9SHong Zhang   PetscErrorCode ierr;
11264b84eec9SHong Zhang 
11274b84eec9SHong Zhang   PetscFunctionBegin;
11284b84eec9SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
11294b84eec9SHong Zhang   PetscValidCharPointer(mprktype,2);
11304b84eec9SHong Zhang   ierr = PetscTryMethod(ts,"TSMPRKSetType_C",(TS,TSMPRKType),(ts,mprktype));CHKERRQ(ierr);
11314b84eec9SHong Zhang   PetscFunctionReturn(0);
11324b84eec9SHong Zhang }
11334b84eec9SHong Zhang 
11344b84eec9SHong Zhang /*@C
11354b84eec9SHong Zhang   TSMPRKGetType - Get the type of MPRK scheme
11364b84eec9SHong Zhang 
11374b84eec9SHong Zhang   Logically collective
11384b84eec9SHong Zhang 
11394b84eec9SHong Zhang   Input Parameter:
11404b84eec9SHong Zhang .  ts - timestepping context
11414b84eec9SHong Zhang 
11424b84eec9SHong Zhang   Output Parameter:
11434b84eec9SHong Zhang .  mprktype - type of MPRK-scheme
11444b84eec9SHong Zhang 
11454b84eec9SHong Zhang   Level: intermediate
11464b84eec9SHong Zhang 
11474b84eec9SHong Zhang .seealso: TSMPRKGetType()
11484b84eec9SHong Zhang @*/
11494b84eec9SHong Zhang PetscErrorCode TSMPRKGetType(TS ts,TSMPRKType *mprktype)
11504b84eec9SHong Zhang {
11514b84eec9SHong Zhang   PetscErrorCode ierr;
11524b84eec9SHong Zhang 
11534b84eec9SHong Zhang   PetscFunctionBegin;
11544b84eec9SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
11554b84eec9SHong Zhang   ierr = PetscUseMethod(ts,"TSMPRKGetType_C",(TS,TSMPRKType*),(ts,mprktype));CHKERRQ(ierr);
11564b84eec9SHong Zhang   PetscFunctionReturn(0);
11574b84eec9SHong Zhang }
11584b84eec9SHong Zhang 
11594b84eec9SHong Zhang /*@C
11604b84eec9SHong Zhang   TSMPRKSetMultirateType - Set the type of MPRK multirate scheme
11614b84eec9SHong Zhang 
11624b84eec9SHong Zhang   Logically collective
11634b84eec9SHong Zhang 
11644b84eec9SHong Zhang   Input Parameter:
11654b84eec9SHong Zhang +  ts - timestepping context
11664b84eec9SHong Zhang -  mprkmtype - type of the multirate configuration
11674b84eec9SHong Zhang 
11684b84eec9SHong Zhang   Options Database:
11694b84eec9SHong Zhang .   -ts_mprk_multirate_type - <nonsplit,split>
11704b84eec9SHong Zhang 
11714b84eec9SHong Zhang   Level: intermediate
11724b84eec9SHong Zhang @*/
11734b84eec9SHong Zhang PetscErrorCode TSMPRKSetMultirateType(TS ts, TSMPRKMultirateType mprkmtype)
11744b84eec9SHong Zhang {
11754b84eec9SHong Zhang   TS_MPRK        *mprk = (TS_MPRK*)ts->data;
11764b84eec9SHong Zhang   PetscErrorCode ierr;
11774b84eec9SHong Zhang 
11784b84eec9SHong Zhang   PetscFunctionBegin;
11794b84eec9SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
11804b84eec9SHong Zhang   switch(mprkmtype){
11814b84eec9SHong Zhang     case TSMPRKNONSPLIT:
11824b84eec9SHong Zhang       ts->ops->step         = TSStep_MPRK;
11834b84eec9SHong Zhang       ts->ops->evaluatestep = TSEvaluateStep_MPRK;
11844b84eec9SHong Zhang       break;
11854b84eec9SHong Zhang     case TSMPRKSPLIT:
11864b84eec9SHong Zhang       ts->ops->step         = TSStep_MPRKSPLIT;
11874b84eec9SHong Zhang       ts->ops->evaluatestep = TSEvaluateStep_MPRKSPLIT;
11884b84eec9SHong Zhang       ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetSplits_C",TSMPRKSetSplits);CHKERRQ(ierr);
11894b84eec9SHong Zhang       break;
11904b84eec9SHong Zhang     default :
11914b84eec9SHong Zhang       SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown type '%s'",mprkmtype);
11924b84eec9SHong Zhang   }
11934b84eec9SHong Zhang   mprk->mprkmtype = mprkmtype;
11944b84eec9SHong Zhang   PetscFunctionReturn(0);
11954b84eec9SHong Zhang }
11964b84eec9SHong Zhang 
11974b84eec9SHong Zhang static PetscErrorCode TSMPRKGetType_MPRK(TS ts,TSMPRKType *mprktype)
11984b84eec9SHong Zhang {
11994b84eec9SHong Zhang   TS_MPRK *mprk = (TS_MPRK*)ts->data;
12004b84eec9SHong Zhang 
12014b84eec9SHong Zhang   PetscFunctionBegin;
12024b84eec9SHong Zhang   *mprktype = mprk->tableau->name;
12034b84eec9SHong Zhang   PetscFunctionReturn(0);
12044b84eec9SHong Zhang }
12054b84eec9SHong Zhang 
12064b84eec9SHong Zhang static PetscErrorCode TSMPRKSetType_MPRK(TS ts,TSMPRKType mprktype)
12074b84eec9SHong Zhang {
12084b84eec9SHong Zhang   TS_MPRK         *mprk = (TS_MPRK*)ts->data;
12094b84eec9SHong Zhang   PetscBool       match;
12104b84eec9SHong Zhang   MPRKTableauLink link;
12114b84eec9SHong Zhang   PetscErrorCode  ierr;
12124b84eec9SHong Zhang 
12134b84eec9SHong Zhang   PetscFunctionBegin;
12144b84eec9SHong Zhang   if (mprk->tableau) {
12154b84eec9SHong Zhang     ierr = PetscStrcmp(mprk->tableau->name,mprktype,&match);CHKERRQ(ierr);
12164b84eec9SHong Zhang     if (match) PetscFunctionReturn(0);
12174b84eec9SHong Zhang   }
12184b84eec9SHong Zhang   for (link = MPRKTableauList; link; link=link->next) {
12194b84eec9SHong Zhang     ierr = PetscStrcmp(link->tab.name,mprktype,&match);CHKERRQ(ierr);
12204b84eec9SHong Zhang     if (match) {
12214b84eec9SHong Zhang       if (ts->setupcalled) {ierr = TSMPRKTableauReset(ts);CHKERRQ(ierr);}
12224b84eec9SHong Zhang       mprk->tableau = &link->tab;
12234b84eec9SHong Zhang       if (ts->setupcalled) {ierr = TSMPRKTableauSetUp(ts);CHKERRQ(ierr);}
12244b84eec9SHong Zhang       PetscFunctionReturn(0);
12254b84eec9SHong Zhang     }
12264b84eec9SHong Zhang   }
12274b84eec9SHong Zhang   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",mprktype);
12284b84eec9SHong Zhang   PetscFunctionReturn(0);
12294b84eec9SHong Zhang }
12304b84eec9SHong Zhang 
12314b84eec9SHong Zhang static PetscErrorCode TSGetStages_MPRK(TS ts,PetscInt *ns,Vec **Y)
12324b84eec9SHong Zhang {
12334b84eec9SHong Zhang   TS_MPRK *mprk = (TS_MPRK*)ts->data;
12344b84eec9SHong Zhang 
12354b84eec9SHong Zhang   PetscFunctionBegin;
12364b84eec9SHong Zhang   *ns = mprk->tableau->s;
12374b84eec9SHong Zhang   if (Y) *Y = mprk->Y;
12384b84eec9SHong Zhang   PetscFunctionReturn(0);
12394b84eec9SHong Zhang }
12404b84eec9SHong Zhang 
12414b84eec9SHong Zhang static PetscErrorCode TSDestroy_MPRK(TS ts)
12424b84eec9SHong Zhang {
12434b84eec9SHong Zhang   PetscErrorCode ierr;
12444b84eec9SHong Zhang 
12454b84eec9SHong Zhang   PetscFunctionBegin;
12464b84eec9SHong Zhang   ierr = TSReset_MPRK(ts);CHKERRQ(ierr);
12474b84eec9SHong Zhang   if (ts->dm) {
12484b84eec9SHong Zhang     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSMPRK,DMRestrictHook_TSMPRK,ts);CHKERRQ(ierr);
12494b84eec9SHong Zhang     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSMPRK,DMSubDomainRestrictHook_TSMPRK,ts);CHKERRQ(ierr);
12504b84eec9SHong Zhang   }
12514b84eec9SHong Zhang   ierr = PetscFree(ts->data);CHKERRQ(ierr);
12524b84eec9SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKGetType_C",NULL);CHKERRQ(ierr);
12534b84eec9SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetType_C",NULL);CHKERRQ(ierr);
12544b84eec9SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetMultirateType_C",NULL);CHKERRQ(ierr);
12554b84eec9SHong Zhang   PetscFunctionReturn(0);
12564b84eec9SHong Zhang }
12574b84eec9SHong Zhang 
12584b84eec9SHong Zhang /*MC
12594b84eec9SHong Zhang       TSMPRK - ODE solver using Partitioned Runge-Kutta schemes
12604b84eec9SHong Zhang 
12614b84eec9SHong Zhang   The user should provide the right hand side of the equation
12624b84eec9SHong Zhang   using TSSetRHSFunction().
12634b84eec9SHong Zhang 
12644b84eec9SHong Zhang   Notes:
12654b84eec9SHong Zhang   The default is TSMPRKPM2, it can be changed with TSRKSetType() or -ts_mprk_type
12664b84eec9SHong Zhang 
12674b84eec9SHong Zhang   Level: beginner
12684b84eec9SHong Zhang 
12694b84eec9SHong Zhang .seealso:  TSCreate(), TS, TSSetType(), TSMPRKSetType(), TSMPRKGetType(), TSMPRKType, TSMPRKRegister(), TSMPRKSetMultirateType()
12704b84eec9SHong Zhang            TSMPRKM2, TSMPRKM3, TSMPRKRFSMR3, TSMPRKRFSMR2
12714b84eec9SHong Zhang 
12724b84eec9SHong Zhang M*/
12734b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSCreate_MPRK(TS ts)
12744b84eec9SHong Zhang {
12754b84eec9SHong Zhang   TS_MPRK        *mprk;
12764b84eec9SHong Zhang   PetscErrorCode ierr;
12774b84eec9SHong Zhang 
12784b84eec9SHong Zhang   PetscFunctionBegin;
12794b84eec9SHong Zhang   ierr = TSMPRKInitializePackage();CHKERRQ(ierr);
12804b84eec9SHong Zhang 
12814b84eec9SHong Zhang   ts->ops->reset          = TSReset_MPRK;
12824b84eec9SHong Zhang   ts->ops->destroy        = TSDestroy_MPRK;
12834b84eec9SHong Zhang   ts->ops->view           = TSView_MPRK;
12844b84eec9SHong Zhang   ts->ops->load           = TSLoad_MPRK;
12854b84eec9SHong Zhang   ts->ops->setup          = TSSetUp_MPRK;
12864b84eec9SHong Zhang   ts->ops->step           = TSStep_MPRK;
12874b84eec9SHong Zhang   ts->ops->evaluatestep   = TSEvaluateStep_MPRK;
12884b84eec9SHong Zhang   ts->ops->setfromoptions = TSSetFromOptions_MPRK;
12894b84eec9SHong Zhang   ts->ops->getstages      = TSGetStages_MPRK;
12904b84eec9SHong Zhang 
12914b84eec9SHong Zhang   ierr = PetscNewLog(ts,&mprk);CHKERRQ(ierr);
12924b84eec9SHong Zhang   ts->data = (void*)mprk;
12934b84eec9SHong Zhang 
12944b84eec9SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKGetType_C",TSMPRKGetType_MPRK);CHKERRQ(ierr);
12954b84eec9SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetType_C",TSMPRKSetType_MPRK);CHKERRQ(ierr);
12964b84eec9SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetMultirateType_C",TSMPRKSetMultirateType);CHKERRQ(ierr);
12974b84eec9SHong Zhang 
12984b84eec9SHong Zhang   ierr = TSMPRKSetType(ts,TSMPRKDefault);CHKERRQ(ierr);
12994b84eec9SHong Zhang   PetscFunctionReturn(0);
13004b84eec9SHong Zhang }
1301