xref: /petsc/src/ts/impls/multirate/mprk.c (revision 79bc8a7798545193596939ed5feadbaa18d20db4)
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.
1219c2959aSHong 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'.
1319c2959aSHong 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.
1419c2959aSHong Zhang 
1519c2959aSHong Zhang   Reference:
1619c2959aSHong Zhang   Emil M. Constantinescu, Adrian Sandu, Multirate Timestepping Methods for Hyperbolic Conservation Laws, Journal of Scientific Computing 2007
1719c2959aSHong Zhang 
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 */
31*79bc8a77SHong 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   TSMPRKMultirateType mprkmtype;
504b84eec9SHong Zhang   Vec                 *Y;                          /* States computed during the step                           */
517c0df07dSHong Zhang   Vec                 *YdotRHS;
524b84eec9SHong Zhang   Vec                 *YdotRHS_slow;               /* Function evaluations by slow tableau for slow components  */
5319c2959aSHong Zhang   Vec                 *YdotRHS_slowbuffer;         /* Function evaluations by slow tableau for slow components  */
5419c2959aSHong Zhang   Vec                 *YdotRHS_medium;             /* Function evaluations by slow tableau for slow components  */
5519c2959aSHong Zhang   Vec                 *YdotRHS_mediumbuffer;       /* Function evaluations by slow tableau for slow components  */
5619c2959aSHong Zhang   Vec                 *YdotRHS_fast;               /* Function evaluations by fast tableau for fast components  */
574b84eec9SHong Zhang   PetscScalar         *work_slow;                  /* Scalar work_slow by slow tableau                          */
5819c2959aSHong Zhang   PetscScalar         *work_slowbuffer;            /* Scalar work_slow by slow tableau                          */
5919c2959aSHong Zhang   PetscScalar         *work_medium;                /* Scalar work_slow by medium tableau                        */
6019c2959aSHong Zhang   PetscScalar         *work_mediumbuffer;          /* Scalar work_slow by medium tableau                        */
6119c2959aSHong Zhang   PetscScalar         *work_fast;                  /* Scalar work_fast by fast tableau                          */
624b84eec9SHong Zhang   PetscReal           stage_time;
634b84eec9SHong Zhang   TSStepStatus        status;
644b84eec9SHong Zhang   PetscReal           ptime;
654b84eec9SHong Zhang   PetscReal           time_step;
6619c2959aSHong Zhang   IS                  is_slow,is_slowbuffer,is_medium,is_mediumbuffer,is_fast;
6719c2959aSHong Zhang   TS                  subts_slow,subts_slowbuffer,subts_medium,subts_mediumbuffer,subts_fast;
684b84eec9SHong Zhang } TS_MPRK;
694b84eec9SHong Zhang 
7019c2959aSHong Zhang static PetscErrorCode TSMPRKGenerateTableau2(PetscInt ratio,PetscInt s,const PetscReal Abase[],const PetscReal bbase[],PetscReal A1[],PetscReal b1[],PetscReal A2[],PetscReal b2[])
7119c2959aSHong Zhang {
7219c2959aSHong Zhang   PetscInt i,j,k,l;
7319c2959aSHong Zhang 
7419c2959aSHong Zhang   PetscFunctionBegin;
7519c2959aSHong Zhang   for (k=0; k<ratio; k++) {
7619c2959aSHong Zhang     /* diagonal blocks */
7719c2959aSHong Zhang     for (i=0; i<s; i++)
7819c2959aSHong Zhang       for (j=0; j<s; j++) {
7919c2959aSHong Zhang         A1[(k*s+i)*ratio*s+k*s+j] = Abase[i*s+j];
8019c2959aSHong Zhang         A2[(k*s+i)*ratio*s+k*s+j] = Abase[i*s+j]/ratio;
8119c2959aSHong Zhang       }
8219c2959aSHong Zhang     /* off diagonal blocks */
8319c2959aSHong Zhang     for (l=0; l<k; l++)
8419c2959aSHong Zhang       for (i=0; i<s; i++)
8519c2959aSHong Zhang         for (j=0; j<s; j++)
8619c2959aSHong Zhang           A2[(k*s+i)*ratio*s+l*s+j] = bbase[j]/ratio;
8719c2959aSHong Zhang     for (j=0; j<s; j++) {
8819c2959aSHong Zhang       b1[k*s+j] = bbase[j]/ratio;
8919c2959aSHong Zhang       b2[k*s+j] = bbase[j]/ratio;
9019c2959aSHong Zhang     }
9119c2959aSHong Zhang   }
9219c2959aSHong Zhang   PetscFunctionReturn(0);
9319c2959aSHong Zhang }
9419c2959aSHong Zhang 
9519c2959aSHong 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[])
9619c2959aSHong Zhang {
9719c2959aSHong Zhang   PetscInt i,j,k,l,m,n;
9819c2959aSHong Zhang 
9919c2959aSHong Zhang   PetscFunctionBegin;
10019c2959aSHong Zhang   for (k=0; k<ratio; k++) { /* diagonal blocks of size ratio*s by ratio*s */
10119c2959aSHong Zhang     for (l=0; l<ratio; l++) /* diagonal sub-blocks of size s by s */
10219c2959aSHong Zhang       for (i=0; i<s; i++)
10319c2959aSHong Zhang         for (j=0; j<s; j++) {
10419c2959aSHong Zhang           A1[((k*ratio+l)*s+i)*ratio*ratio*s+(k*ratio+l)*s+j] = Abase[i*s+j];
10519c2959aSHong Zhang           A2[((k*ratio+l)*s+i)*ratio*ratio*s+(k*ratio+l)*s+j] = Abase[i*s+j]/ratio;
10619c2959aSHong Zhang           A3[((k*ratio+l)*s+i)*ratio*ratio*s+(k*ratio+l)*s+j] = Abase[i*s+j]/ratio/ratio;
10719c2959aSHong Zhang         }
10819c2959aSHong Zhang     for (l=0; l<k; l++) /* off-diagonal blocks of size ratio*s by ratio*s */
10919c2959aSHong Zhang       for (m=0; m<ratio; m++)
11019c2959aSHong Zhang         for (n=0; n<ratio; n++)
11119c2959aSHong Zhang           for (i=0; i<s; i++)
11219c2959aSHong Zhang             for (j=0; j<s; j++) {
11319c2959aSHong Zhang                A2[((k*ratio+m)*s+i)*ratio*ratio*s+(l*ratio+n)*s+j] = bbase[j]/ratio/ratio;
11419c2959aSHong Zhang                A3[((k*ratio+m)*s+i)*ratio*ratio*s+(l*ratio+n)*s+j] = bbase[j]/ratio/ratio;
11519c2959aSHong Zhang             }
11619c2959aSHong Zhang     for (m=0; m<ratio; m++)
11719c2959aSHong Zhang       for (n=0; n<m; n++) /* off-diagonal sub-blocks of size s by s in the diagonal blocks */
11819c2959aSHong Zhang           for (i=0; i<s; i++)
11919c2959aSHong Zhang             for (j=0; j<s; j++)
12019c2959aSHong Zhang                A3[((k*ratio+m)*s+i)*ratio*ratio*s+(k*ratio+n)*s+j] = bbase[j]/ratio/ratio;
12119c2959aSHong Zhang     for (n=0; n<ratio; n++)
12219c2959aSHong Zhang       for (j=0; j<s; j++) {
12319c2959aSHong Zhang         b1[(k*ratio+n)*s+j] = bbase[j]/ratio/ratio;
12419c2959aSHong Zhang         b2[(k*ratio+n)*s+j] = bbase[j]/ratio/ratio;
12519c2959aSHong Zhang         b3[(k*ratio+n)*s+j] = bbase[j]/ratio/ratio;
12619c2959aSHong Zhang       }
12719c2959aSHong Zhang   }
12819c2959aSHong Zhang   PetscFunctionReturn(0);
12919c2959aSHong Zhang }
13019c2959aSHong Zhang 
1314b84eec9SHong Zhang /*MC
13219c2959aSHong Zhang      TSMPRK2A22 - Second Order Partitioned Runge Kutta scheme based on RK2A.
1334b84eec9SHong Zhang 
13419c2959aSHong Zhang      This method has four stages for slow and fast parts. The refinement factor of the stepsize is 2.
13519c2959aSHong Zhang      r = 2, np = 2
1364b84eec9SHong Zhang      Options database:
13719c2959aSHong Zhang .     -ts_mprk_type 2a22
1384b84eec9SHong Zhang 
1394b84eec9SHong Zhang      Level: advanced
1404b84eec9SHong Zhang 
1414b84eec9SHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType()
1424b84eec9SHong Zhang M*/
1434b84eec9SHong Zhang /*MC
14419c2959aSHong Zhang      TSMPRK2A23 - Second Order Partitioned Runge Kutta scheme based on RK2A.
14519c2959aSHong Zhang 
14619c2959aSHong Zhang      This method has eight stages for slow and medium and fast parts. The refinement factor of the stepsize is 2.
14719c2959aSHong Zhang      r = 2, np = 3
14819c2959aSHong Zhang      Options database:
14919c2959aSHong Zhang .     -ts_mprk_type 2a23
15019c2959aSHong Zhang 
15119c2959aSHong Zhang      Level: advanced
15219c2959aSHong Zhang 
15319c2959aSHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType()
15419c2959aSHong Zhang M*/
15519c2959aSHong Zhang /*MC
15619c2959aSHong Zhang      TSMPRK2A32 - Second Order Partitioned Runge Kutta scheme based on RK2A.
15719c2959aSHong Zhang 
15819c2959aSHong Zhang      This method has four stages for slow and fast parts. The refinement factor of the stepsize is 3.
15919c2959aSHong Zhang      r = 3, np = 2
16019c2959aSHong Zhang      Options database:
16119c2959aSHong Zhang .     -ts_mprk_type 2a32
16219c2959aSHong Zhang 
16319c2959aSHong Zhang      Level: advanced
16419c2959aSHong Zhang 
16519c2959aSHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType()
16619c2959aSHong Zhang M*/
16719c2959aSHong Zhang /*MC
16819c2959aSHong Zhang      TSMPRK2A33 - Second Order 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
17219c2959aSHong Zhang      Options database:
17319c2959aSHong Zhang .     -ts_mprk_type 2a33
17419c2959aSHong Zhang 
17519c2959aSHong Zhang      Level: advanced
17619c2959aSHong Zhang 
17719c2959aSHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType()
17819c2959aSHong Zhang M*/
17919c2959aSHong Zhang /*MC
18019c2959aSHong Zhang      TSMPRK3P2M - Third Order Partitioned Runge Kutta scheme.
1814b84eec9SHong Zhang 
1824b84eec9SHong Zhang      This method has eight stages for both slow and fast parts.
1834b84eec9SHong Zhang 
1844b84eec9SHong Zhang      Options database:
1854b84eec9SHong Zhang .     -ts_mprk_type pm3  (put here temporarily)
1864b84eec9SHong Zhang 
1874b84eec9SHong Zhang      Level: advanced
1884b84eec9SHong Zhang 
1894b84eec9SHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType()
1904b84eec9SHong Zhang M*/
1914b84eec9SHong Zhang /*MC
1924b84eec9SHong Zhang      TSMPRKP2 - Second Order Partitioned Runge Kutta scheme.
1934b84eec9SHong Zhang 
1944b84eec9SHong Zhang      This method has five stages for both slow and fast parts.
1954b84eec9SHong Zhang 
1964b84eec9SHong Zhang      Options database:
1974b84eec9SHong Zhang .     -ts_mprk_type p2
1984b84eec9SHong Zhang 
1994b84eec9SHong Zhang      Level: advanced
2004b84eec9SHong Zhang 
2014b84eec9SHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType()
2024b84eec9SHong Zhang M*/
2034b84eec9SHong Zhang /*MC
2044b84eec9SHong Zhang      TSMPRKP3 - Third Order Partitioned Runge Kutta scheme.
2054b84eec9SHong Zhang 
2064b84eec9SHong Zhang      This method has ten stages for both slow and fast parts.
2074b84eec9SHong Zhang 
2084b84eec9SHong Zhang      Options database:
2094b84eec9SHong Zhang .     -ts_mprk_type p3
2104b84eec9SHong Zhang 
2114b84eec9SHong Zhang      Level: advanced
2124b84eec9SHong Zhang 
2134b84eec9SHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType()
2144b84eec9SHong Zhang M*/
2154b84eec9SHong Zhang 
2164b84eec9SHong Zhang /*@C
2174b84eec9SHong Zhang   TSMPRKRegisterAll - Registers all of the Partirioned Runge-Kutta explicit methods in TSMPRK
2184b84eec9SHong Zhang 
2194b84eec9SHong Zhang   Not Collective, but should be called by all processes which will need the schemes to be registered
2204b84eec9SHong Zhang 
2214b84eec9SHong Zhang   Level: advanced
2224b84eec9SHong Zhang 
2234b84eec9SHong Zhang .keywords: TS, TSMPRK, register, all
2244b84eec9SHong Zhang 
2254b84eec9SHong Zhang .seealso:  TSMPRKRegisterDestroy()
2264b84eec9SHong Zhang @*/
2274b84eec9SHong Zhang PetscErrorCode TSMPRKRegisterAll(void)
2284b84eec9SHong Zhang {
2294b84eec9SHong Zhang   PetscErrorCode ierr;
2304b84eec9SHong Zhang 
2314b84eec9SHong Zhang   PetscFunctionBegin;
2324b84eec9SHong Zhang   if (TSMPRKRegisterAllCalled) PetscFunctionReturn(0);
2334b84eec9SHong Zhang   TSMPRKRegisterAllCalled = PETSC_TRUE;
2344b84eec9SHong Zhang 
2354b84eec9SHong Zhang #define RC PetscRealConstant
2364b84eec9SHong Zhang   {
2374b84eec9SHong Zhang     const PetscReal
23819c2959aSHong Zhang       Abase[2][2] = {{0,0},
23919c2959aSHong Zhang                      {RC(1.0),0}},
24019c2959aSHong Zhang       bbase[2] = {RC(0.5),RC(0.5)};
24119c2959aSHong Zhang     PetscReal
24219c2959aSHong Zhang       Asb[4][4] = {{0}},Af[4][4] = {{0}},bsb[4] = {0},bf[4] = {0};
24319c2959aSHong Zhang     PetscInt
24419c2959aSHong Zhang       rsb[4] = {0,0,1,2};
24519c2959aSHong Zhang     ierr = TSMPRKGenerateTableau2(2,2,&Abase[0][0],bbase,&Asb[0][0],bsb,&Af[0][0],bf);CHKERRQ(ierr);
246*79bc8a77SHong Zhang     ierr = TSMPRKRegister(TSMPRK2A22,2,2,2,1,&Asb[0][0],bsb,NULL,rsb,NULL,NULL,NULL,NULL,&Af[0][0],bf,NULL);CHKERRQ(ierr);
2474b84eec9SHong Zhang   }
24819c2959aSHong Zhang   {
24919c2959aSHong Zhang     const PetscReal
25019c2959aSHong Zhang       Abase[2][2] = {{0,0},
25119c2959aSHong Zhang                      {RC(1.0),0}},
25219c2959aSHong Zhang       bbase[2]    = {RC(0.5),RC(0.5)};
25319c2959aSHong Zhang     PetscReal
25419c2959aSHong Zhang       Asb[8][8] = {{0}},Amb[8][8] = {{0}},Af[8][8] = {{0}},bsb[8] ={0},bmb[8] = {0},bf[8] = {0};
25519c2959aSHong Zhang     PetscInt
25619c2959aSHong Zhang       rsb[8] = {0,0,1,2,1,2,1,2},rmb[8] = {0,0,1,2,0,0,5,6};
25719c2959aSHong Zhang     ierr = TSMPRKGenerateTableau3(2,2,&Abase[0][0],bbase,&Asb[0][0],bsb,&Amb[0][0],bmb,&Af[0][0],bf);CHKERRQ(ierr);
258*79bc8a77SHong Zhang     ierr = TSMPRKRegister(TSMPRK2A23,2,2,2,2,&Asb[0][0],bsb,NULL,rsb,&Amb[0][0],bmb,NULL,rmb,&Af[0][0],bf,NULL);CHKERRQ(ierr);
25919c2959aSHong Zhang   }
26019c2959aSHong Zhang   {
26119c2959aSHong Zhang     const PetscReal
26219c2959aSHong Zhang       Abase[2][2] = {{0,0},
26319c2959aSHong Zhang                      {RC(1.0),0}},
26419c2959aSHong Zhang       bbase[2]    = {RC(0.5),RC(0.5)};
26519c2959aSHong Zhang     PetscReal
26619c2959aSHong Zhang       Asb[6][6] = {{0}},Af[6][6] = {{0}},bsb[6] = {0},bf[6] = {0};
26719c2959aSHong Zhang     PetscInt
26819c2959aSHong Zhang       rsb[6] = {0,0,1,2,1,2};
26919c2959aSHong Zhang     ierr = TSMPRKGenerateTableau2(3,2,&Abase[0][0],bbase,&Asb[0][0],bsb,&Af[0][0],bf);CHKERRQ(ierr);
270*79bc8a77SHong Zhang     ierr = TSMPRKRegister(TSMPRK2A32,2,2,3,1,&Asb[0][0],bsb,NULL,rsb,NULL,NULL,NULL,NULL,&Af[0][0],bf,NULL);CHKERRQ(ierr);
27119c2959aSHong Zhang   }
27219c2959aSHong Zhang   {
27319c2959aSHong Zhang     const PetscReal
27419c2959aSHong Zhang       Abase[2][2] = {{0,0},
27519c2959aSHong Zhang                      {RC(1.0),0}},
27619c2959aSHong Zhang       bbase[2]    = {RC(0.5),RC(0.5)};
27719c2959aSHong Zhang     PetscReal
27819c2959aSHong Zhang       Asb[18][18] = {{0}},Amb[18][18] = {{0}},Af[18][18] = {{0}},bsb[18] ={0},bmb[18] = {0},bf[18] = {0};
27919c2959aSHong Zhang     PetscInt
28019c2959aSHong 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};
28119c2959aSHong Zhang     ierr = TSMPRKGenerateTableau3(3,2,&Abase[0][0],bbase,&Asb[0][0],bsb,&Amb[0][0],bmb,&Af[0][0],bf);CHKERRQ(ierr);
282*79bc8a77SHong Zhang     ierr = TSMPRKRegister(TSMPRK2A33,2,2,3,3,&Asb[0][0],bsb,NULL,rsb,&Amb[0][0],bmb,NULL,rmb,&Af[0][0],bf,NULL);CHKERRQ(ierr);
28319c2959aSHong Zhang   }
28419c2959aSHong Zhang /*
28519c2959aSHong Zhang     PetscReal
28619c2959aSHong Zhang       Asb[8][8] = {{Abase[0][0],Abase[0][1],0,0,0,0,0,0},
28719c2959aSHong Zhang                    {Abase[1][0],Abase[1][1],0,0,0,0,0,0},
28819c2959aSHong Zhang                    {0,0,Abase[0][0],Abase[0][1],0,0,0,0},
28919c2959aSHong Zhang                    {0,0,Abase[1][0],Abase[1][1],0,0,0,0},
29019c2959aSHong Zhang                    {0,0,0,0,Abase[0][0],Abase[0][1],0,0},
29119c2959aSHong Zhang                    {0,0,0,0,Abase[1][0],Abase[1][1],0,0},
29219c2959aSHong Zhang                    {0,0,0,0,0,0,Abase[0][0],Abase[0][1]},
29319c2959aSHong Zhang                    {0,0,0,0,0,0,Abase[1][0],Abase[1][1]}},
29419c2959aSHong Zhang       Amb[8][8] = {{Abase[0][0]/m,Abase[0][1]/m,0,0,0,0,0,0},
29519c2959aSHong Zhang                    {Abase[1][0]/m,Abase[1][1]/m,0,0,0,0,0,0},
29619c2959aSHong Zhang                    {0,0,Abase[0][0]/m,Abase[0][1]/m,0,0,0,0},
29719c2959aSHong Zhang                    {0,0,Abase[1][0]/m,Abase[1][1]/m,0,0,0,0},
29819c2959aSHong Zhang                    {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,Abase[0][0]/m,Abase[0][1]/m,0,0},
29919c2959aSHong Zhang                    {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,Abase[1][0]/m,Abase[1][1]/m,0,0},
30019c2959aSHong Zhang                    {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,0,0,Abase[0][0]/m,Abase[0][1]/m},
30119c2959aSHong Zhang                    {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,0,0,Abase[1][0]/m,Abase[1][1]/m}},
30219c2959aSHong Zhang       Af[8][8] = {{Abase[0][0]/m/m,Abase[0][1]/m/m,0,0,0,0,0,0},
30319c2959aSHong Zhang                    {Abase[1][0]/m/m,Abase[1][1]/m/m,0,0,0,0,0,0},
30419c2959aSHong Zhang                    {0,0,Abase[0][0]/m/m,Abase[0][1]/m/m,0,0,0,0},
30519c2959aSHong Zhang                    {0,0,Abase[1][0]/m/m,Abase[1][1]/m/m,0,0,0,0},
30619c2959aSHong 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},
30719c2959aSHong 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},
30819c2959aSHong 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},
30919c2959aSHong 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}},
31019c2959aSHong 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},
31119c2959aSHong 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},
31219c2959aSHong 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},
31319c2959aSHong Zhang */
3144b84eec9SHong Zhang   /*{
3154b84eec9SHong Zhang       const PetscReal
3164b84eec9SHong Zhang         As[8][8] = {{0,0,0,0,0,0,0,0},
3174b84eec9SHong Zhang                     {RC(1.0)/RC(2.0),0,0,0,0,0,0,0},
3184b84eec9SHong Zhang                     {RC(-1.0)/RC(6.0),RC(2.0)/RC(3.0),0,0,0,0,0,0},
3194b84eec9SHong Zhang                     {RC(1.0)/RC(3.0),RC(-1.0)/RC(3.0),RC(1.0),0,0,0,0,0},
3204b84eec9SHong Zhang                     {0,0,0,0,0,0,0,0},
3214b84eec9SHong Zhang                     {0,0,0,0,RC(1.0)/RC(2.0),0,0,0},
3224b84eec9SHong Zhang                     {0,0,0,0,RC(-1.0)/RC(6.0),RC(2.0)/RC(3.0),0,0},
3234b84eec9SHong Zhang                     {0,0,0,0,RC(1.0)/RC(3.0),RC(-1.0)/RC(3.0),RC(1.0),0}},
3244b84eec9SHong Zhang          A[8][8] = {{0,0,0,0,0,0,0,0},
3254b84eec9SHong Zhang                     {RC(1.0)/RC(4.0),0,0,0,0,0,0,0},
3264b84eec9SHong Zhang                     {RC(-1.0)/RC(12.0),RC(1.0)/RC(3.0),0,0,0,0,0,0},
3274b84eec9SHong Zhang                     {RC(1.0)/RC(6.0),RC(-1.0)/RC(6.0),RC(1.0)/RC(2.0),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),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(4.0),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(12.0),RC(1.0)/RC(3.0),0,0},
3314b84eec9SHong 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}},
3324b84eec9SHong 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)},
3334b84eec9SHong 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)};
3344b84eec9SHong Zhang            ierr = TSMPRKRegister(TSMPRKPM3,3,8,&As[0][0],bs,NULL,&A[0][0],b,NULL);CHKERRQ(ierr);
3354b84eec9SHong Zhang   }*/
3364b84eec9SHong Zhang 
3374b84eec9SHong Zhang   {
3384b84eec9SHong Zhang     const PetscReal
33919c2959aSHong Zhang       Asb[5][5] = {{0,0,0,0,0},
3404b84eec9SHong Zhang                    {RC(1.0)/RC(2.0),0,0,0,0},
3414b84eec9SHong Zhang                    {RC(1.0)/RC(2.0),0,0,0,0},
3424b84eec9SHong Zhang                    {RC(1.0),0,0,0,0},
3434b84eec9SHong Zhang                    {RC(1.0),0,0,0,0}},
34419c2959aSHong Zhang       Af[5][5]  = {{0,0,0,0,0},
3454b84eec9SHong Zhang                    {RC(1.0)/RC(2.0),0,0,0,0},
3464b84eec9SHong Zhang                    {RC(1.0)/RC(4.0),RC(1.0)/RC(4.0),0,0,0},
3474b84eec9SHong Zhang                    {RC(1.0)/RC(4.0),RC(1.0)/RC(4.0),RC(1.0)/RC(2.0),0,0},
3484b84eec9SHong 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}},
34919c2959aSHong Zhang       bsb[5]    = {RC(1.0)/RC(2.0),0,0,0,RC(1.0)/RC(2.0)},
35019c2959aSHong 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};
35119c2959aSHong Zhang     const PetscInt
35219c2959aSHong Zhang       rsb[5]    = {0,0,2,0,4};
353*79bc8a77SHong Zhang     ierr = TSMPRKRegister(TSMPRKP2,2,5,1,1,&Asb[0][0],bsb,NULL,rsb,NULL,NULL,NULL,NULL,&Af[0][0],bf,NULL);CHKERRQ(ierr);
3544b84eec9SHong Zhang   }
3554b84eec9SHong Zhang 
3564b84eec9SHong Zhang   {
3574b84eec9SHong Zhang     const PetscReal
35819c2959aSHong Zhang       Asb[10][10] = {{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(4.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(2.0),0,0,0,0,0,0,0,0,0},
3634b84eec9SHong Zhang                      {RC(-1.0)/RC(6.0),0,0,0,RC(2.0)/RC(3.0),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(12.0),0,0,0,RC(1.0)/RC(6.0),RC(1.0)/RC(2.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},
3674b84eec9SHong Zhang                      {RC(1.0)/RC(3.0),0,0,0,RC(-1.0)/RC(3.0),RC(1.0),0,0,0,0}},
36819c2959aSHong Zhang       Af[10][10]  = {{0,0,0,0,0,0,0,0,0,0},
3694b84eec9SHong Zhang                      {RC(1.0)/RC(4.0),0,0,0,0,0,0,0,0,0},
3704b84eec9SHong Zhang                      {RC(-1.0)/RC(12.0),RC(1.0)/RC(3.0),0,0,0,0,0,0,0,0},
3714b84eec9SHong 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},
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,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(4.0),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(12.0),RC(1.0)/RC(3.0),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(6.0),RC(-1.0)/RC(6.0),RC(1.0)/RC(2.0),0,0},
3774b84eec9SHong 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}},
37819c2959aSHong 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)},
37919c2959aSHong 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};
38019c2959aSHong Zhang     const PetscInt
38119c2959aSHong Zhang       rsb[10]     = {0,0,2,0,4,0,0,7,0,9};
382*79bc8a77SHong Zhang     ierr = TSMPRKRegister(TSMPRKP3,3,5,2,1,&Asb[0][0],bsb,NULL,rsb,NULL,NULL,NULL,NULL,&Af[0][0],bf,NULL);CHKERRQ(ierr);
3834b84eec9SHong Zhang   }
3844b84eec9SHong Zhang #undef RC
3854b84eec9SHong Zhang   PetscFunctionReturn(0);
3864b84eec9SHong Zhang }
3874b84eec9SHong Zhang 
3884b84eec9SHong Zhang /*@C
3894b84eec9SHong Zhang    TSMPRKRegisterDestroy - Frees the list of schemes that were registered by TSMPRKRegister().
3904b84eec9SHong Zhang 
3914b84eec9SHong Zhang    Not Collective
3924b84eec9SHong Zhang 
3934b84eec9SHong Zhang    Level: advanced
3944b84eec9SHong Zhang 
3954b84eec9SHong Zhang .keywords: TSMPRK, register, destroy
3964b84eec9SHong Zhang .seealso: TSMPRKRegister(), TSMPRKRegisterAll()
3974b84eec9SHong Zhang @*/
3984b84eec9SHong Zhang PetscErrorCode TSMPRKRegisterDestroy(void)
3994b84eec9SHong Zhang {
4004b84eec9SHong Zhang   PetscErrorCode ierr;
4014b84eec9SHong Zhang   MPRKTableauLink link;
4024b84eec9SHong Zhang 
4034b84eec9SHong Zhang   PetscFunctionBegin;
4044b84eec9SHong Zhang   while ((link = MPRKTableauList)) {
4054b84eec9SHong Zhang     MPRKTableau t = &link->tab;
4064b84eec9SHong Zhang     MPRKTableauList = link->next;
40719c2959aSHong Zhang     ierr = PetscFree3(t->Asb,t->bsb,t->csb);CHKERRQ(ierr);
40819c2959aSHong Zhang     ierr = PetscFree3(t->Amb,t->bmb,t->cmb);CHKERRQ(ierr);
4094b84eec9SHong Zhang     ierr = PetscFree3(t->Af,t->bf,t->cf);CHKERRQ(ierr);
41019c2959aSHong Zhang     ierr = PetscFree2(t->rsb,t->rmb);CHKERRQ(ierr);
4114b84eec9SHong Zhang     ierr = PetscFree (t->name);CHKERRQ(ierr);
4124b84eec9SHong Zhang     ierr = PetscFree (link);CHKERRQ(ierr);
4134b84eec9SHong Zhang   }
4144b84eec9SHong Zhang   TSMPRKRegisterAllCalled = PETSC_FALSE;
4154b84eec9SHong Zhang   PetscFunctionReturn(0);
4164b84eec9SHong Zhang }
4174b84eec9SHong Zhang 
4184b84eec9SHong Zhang /*@C
4194b84eec9SHong Zhang   TSMPRKInitializePackage - This function initializes everything in the TSMPRK package. It is called
4204b84eec9SHong Zhang   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_MPRK()
4214b84eec9SHong Zhang   when using static libraries.
4224b84eec9SHong Zhang 
4234b84eec9SHong Zhang   Level: developer
4244b84eec9SHong Zhang 
4254b84eec9SHong Zhang .keywords: TS, TSMPRK, initialize, package
4264b84eec9SHong Zhang .seealso: PetscInitialize()
4274b84eec9SHong Zhang @*/
4284b84eec9SHong Zhang PetscErrorCode TSMPRKInitializePackage(void)
4294b84eec9SHong Zhang {
4304b84eec9SHong Zhang   PetscErrorCode ierr;
4314b84eec9SHong Zhang 
4324b84eec9SHong Zhang   PetscFunctionBegin;
4334b84eec9SHong Zhang   if (TSMPRKPackageInitialized) PetscFunctionReturn(0);
4344b84eec9SHong Zhang   TSMPRKPackageInitialized = PETSC_TRUE;
4354b84eec9SHong Zhang   ierr = TSMPRKRegisterAll();CHKERRQ(ierr);
4364b84eec9SHong Zhang   ierr = PetscRegisterFinalize(TSMPRKFinalizePackage);CHKERRQ(ierr);
4374b84eec9SHong Zhang   PetscFunctionReturn(0);
4384b84eec9SHong Zhang }
4394b84eec9SHong Zhang 
4404b84eec9SHong Zhang /*@C
4414b84eec9SHong Zhang   TSRKFinalizePackage - This function destroys everything in the TSMPRK package. It is
4424b84eec9SHong Zhang   called from PetscFinalize().
4434b84eec9SHong Zhang 
4444b84eec9SHong Zhang   Level: developer
4454b84eec9SHong Zhang 
4464b84eec9SHong Zhang .keywords: Petsc, destroy, package
4474b84eec9SHong Zhang .seealso: PetscFinalize()
4484b84eec9SHong Zhang @*/
4494b84eec9SHong Zhang PetscErrorCode TSMPRKFinalizePackage(void)
4504b84eec9SHong Zhang {
4514b84eec9SHong Zhang   PetscErrorCode ierr;
4524b84eec9SHong Zhang 
4534b84eec9SHong Zhang   PetscFunctionBegin;
4544b84eec9SHong Zhang   TSMPRKPackageInitialized = PETSC_FALSE;
4554b84eec9SHong Zhang   ierr = TSMPRKRegisterDestroy();CHKERRQ(ierr);
4564b84eec9SHong Zhang   PetscFunctionReturn(0);
4574b84eec9SHong Zhang }
4584b84eec9SHong Zhang 
4594b84eec9SHong Zhang /*@C
4604b84eec9SHong Zhang    TSMPRKRegister - register a MPRK scheme by providing the entries in the Butcher tableau
4614b84eec9SHong Zhang 
4624b84eec9SHong Zhang    Not Collective, but the same schemes should be registered on all processes on which they will be used
4634b84eec9SHong Zhang 
4644b84eec9SHong Zhang    Input Parameters:
4654b84eec9SHong Zhang +  name - identifier for method
4664b84eec9SHong Zhang .  order - approximation order of method
467*79bc8a77SHong Zhang .  s  - number of stages in the base methods
468*79bc8a77SHong Zhang .  ratio1 - stepsize ratio at 1st level (e.g. slow/medium)
469*79bc8a77SHong Zhang .  ratio2 - stepsize ratio at 2nd level (e.g. medium/fast)
4704b84eec9SHong Zhang .  Af - stage coefficients for fast components(dimension s*s, row-major)
4714b84eec9SHong Zhang .  bf - step completion table for fast components(dimension s)
4724b84eec9SHong Zhang .  cf - abscissa for fast components(dimension s)
4734b84eec9SHong Zhang .  As - stage coefficients for slow components(dimension s*s, row-major)
4744b84eec9SHong Zhang .  bs - step completion table for slow components(dimension s)
4754b84eec9SHong Zhang -  cs - abscissa for slow components(dimension s)
4764b84eec9SHong Zhang 
4774b84eec9SHong Zhang    Notes:
4784b84eec9SHong Zhang    Several MPRK methods are provided, this function is only needed to create new methods.
4794b84eec9SHong Zhang 
4804b84eec9SHong Zhang    Level: advanced
4814b84eec9SHong Zhang 
4824b84eec9SHong Zhang .keywords: TS, register
4834b84eec9SHong Zhang 
4844b84eec9SHong Zhang .seealso: TSMPRK
4854b84eec9SHong Zhang @*/
486*79bc8a77SHong Zhang PetscErrorCode TSMPRKRegister(TSMPRKType name,PetscInt order,
487*79bc8a77SHong Zhang                               PetscInt sbase,PetscInt ratio1,PetscInt ratio2,
48819c2959aSHong Zhang                               const PetscReal Asb[],const PetscReal bsb[],const PetscReal csb[],const PetscInt rsb[],
48919c2959aSHong Zhang                               const PetscReal Amb[],const PetscReal bmb[],const PetscReal cmb[],const PetscInt rmb[],
4904b84eec9SHong Zhang                               const PetscReal Af[],const PetscReal bf[],const PetscReal cf[])
4914b84eec9SHong Zhang {
4924b84eec9SHong Zhang   MPRKTableauLink link;
4934b84eec9SHong Zhang   MPRKTableau     t;
494*79bc8a77SHong Zhang   PetscInt        s,i,j;
4954b84eec9SHong Zhang   PetscErrorCode  ierr;
4964b84eec9SHong Zhang 
4974b84eec9SHong Zhang   PetscFunctionBegin;
4984b84eec9SHong Zhang   PetscValidCharPointer(name,1);
49919c2959aSHong Zhang   PetscValidRealPointer(Asb,4);
500*79bc8a77SHong Zhang   if (bsb) PetscValidRealPointer(bsb,7);
501*79bc8a77SHong Zhang   if (csb) PetscValidRealPointer(csb,8);
502*79bc8a77SHong Zhang   if (rsb) PetscValidRealPointer(rsb,9);
503*79bc8a77SHong Zhang   if (Amb) PetscValidRealPointer(Amb,10);
504*79bc8a77SHong Zhang   if (bmb) PetscValidRealPointer(bmb,11);
505*79bc8a77SHong Zhang   if (cmb) PetscValidRealPointer(cmb,12);
506*79bc8a77SHong Zhang   if (rmb) PetscValidRealPointer(rmb,13);
507*79bc8a77SHong Zhang   PetscValidRealPointer(Af,14);
508*79bc8a77SHong Zhang   if (bf) PetscValidRealPointer(bf,15);
509*79bc8a77SHong Zhang   if (cf) PetscValidRealPointer(cf,16);
5104b84eec9SHong Zhang 
5114b84eec9SHong Zhang   ierr = PetscNew(&link);CHKERRQ(ierr);
5124b84eec9SHong Zhang   t = &link->tab;
5134b84eec9SHong Zhang 
5144b84eec9SHong Zhang   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
515*79bc8a77SHong Zhang   s = sbase*ratio1*ratio2; /*  this is the dimension of the matrices below */
5164b84eec9SHong Zhang   t->order = order;
517*79bc8a77SHong Zhang   t->sbase = sbase;
5184b84eec9SHong Zhang   t->s  = s;
51919c2959aSHong Zhang   t->np = 2;
520*79bc8a77SHong Zhang 
5214b84eec9SHong Zhang   ierr = PetscMalloc3(s*s,&t->Af,s,&t->bf,s,&t->cf);CHKERRQ(ierr);
5224b84eec9SHong Zhang   ierr = PetscMemcpy(t->Af,Af,s*s*sizeof(Af[0]));CHKERRQ(ierr);
5234b84eec9SHong Zhang   if (bf) {
5244b84eec9SHong Zhang     ierr = PetscMemcpy(t->bf,bf,s*sizeof(bf[0]));CHKERRQ(ierr);
52519c2959aSHong Zhang   } else
5264b84eec9SHong Zhang     for (i=0; i<s; i++) t->bf[i] = Af[(s-1)*s+i];
5274b84eec9SHong Zhang   if (cf) {
5284b84eec9SHong Zhang     ierr = PetscMemcpy(t->cf,cf,s*sizeof(cf[0]));CHKERRQ(ierr);
52919c2959aSHong Zhang   } else {
5304b84eec9SHong Zhang     for (i=0; i<s; i++)
5314b84eec9SHong Zhang       for (j=0,t->cf[i]=0; j<s; j++)
5324b84eec9SHong Zhang         t->cf[i] += Af[i*s+j];
5334b84eec9SHong Zhang   }
53419c2959aSHong Zhang 
53519c2959aSHong Zhang   if (Amb) {
53619c2959aSHong Zhang     t->np = 3;
53719c2959aSHong Zhang     ierr = PetscMalloc3(s*s,&t->Amb,s,&t->bmb,s,&t->cmb);CHKERRQ(ierr);
53819c2959aSHong Zhang     ierr = PetscCalloc1(s,&t->rmb);CHKERRQ(ierr);
53919c2959aSHong Zhang     ierr = PetscMemcpy(t->Amb,Amb,s*s*sizeof(Amb[0]));CHKERRQ(ierr);
54019c2959aSHong Zhang     if (bmb) {
54119c2959aSHong Zhang       ierr = PetscMemcpy(t->bmb,bmb,s*sizeof(bmb[0]));CHKERRQ(ierr);
54219c2959aSHong Zhang     } else {
54319c2959aSHong Zhang       for (i=0; i<s; i++) t->bmb[i] = Amb[(s-1)*s+i];
5444b84eec9SHong Zhang     }
54519c2959aSHong Zhang     if (cmb) {
54619c2959aSHong Zhang       ierr = PetscMemcpy(t->cmb,cmb,s*sizeof(cmb[0]));CHKERRQ(ierr);
54719c2959aSHong Zhang     } else {
5484b84eec9SHong Zhang       for (i=0; i<s; i++)
54919c2959aSHong Zhang         for (j=0,t->cmb[i]=0; j<s; j++)
55019c2959aSHong Zhang           t->cmb[i] += Amb[i*s+j];
55119c2959aSHong Zhang     }
55219c2959aSHong Zhang     if (rmb) {
55319c2959aSHong Zhang       ierr = PetscMemcpy(t->rmb,rmb,s*sizeof(rmb[0]));CHKERRQ(ierr);
55419c2959aSHong Zhang     }
55519c2959aSHong Zhang   }
55619c2959aSHong Zhang 
55719c2959aSHong Zhang   ierr = PetscMalloc3(s*s,&t->Asb,s,&t->bsb,s,&t->csb);CHKERRQ(ierr);
55819c2959aSHong Zhang   ierr = PetscMemcpy(t->Asb,Asb,s*s*sizeof(Asb[0]));CHKERRQ(ierr);
55919c2959aSHong Zhang   ierr = PetscCalloc1(s,&t->rsb);CHKERRQ(ierr);
56019c2959aSHong Zhang   if (bsb) {
56119c2959aSHong Zhang     ierr = PetscMemcpy(t->bsb,bsb,s*sizeof(bsb[0]));CHKERRQ(ierr);
56219c2959aSHong Zhang   } else
56319c2959aSHong Zhang     for (i=0; i<s; i++) t->bsb[i] = Asb[(s-1)*s+i];
56419c2959aSHong Zhang   if (csb) {
56519c2959aSHong Zhang     ierr = PetscMemcpy(t->csb,csb,s*sizeof(csb[0]));CHKERRQ(ierr);
56619c2959aSHong Zhang   } else {
56719c2959aSHong Zhang     for (i=0; i<s; i++)
56819c2959aSHong Zhang       for (j=0,t->csb[i]=0; j<s; j++)
56919c2959aSHong Zhang         t->csb[i] += Asb[i*s+j];
57019c2959aSHong Zhang   }
57119c2959aSHong Zhang   if (rsb) {
57219c2959aSHong Zhang     ierr = PetscMemcpy(t->rsb,rsb,s*sizeof(rsb[0]));CHKERRQ(ierr);
5734b84eec9SHong Zhang   }
5744b84eec9SHong Zhang   link->next = MPRKTableauList;
5754b84eec9SHong Zhang   MPRKTableauList = link;
5764b84eec9SHong Zhang   PetscFunctionReturn(0);
5774b84eec9SHong Zhang }
5784b84eec9SHong Zhang 
5794b84eec9SHong Zhang static PetscErrorCode TSMPRKSetSplits(TS ts)
5804b84eec9SHong Zhang {
5814b84eec9SHong Zhang   TS_MPRK        *mprk = (TS_MPRK*)ts->data;
58219c2959aSHong Zhang   MPRKTableau    tab = mprk->tableau;
5834b84eec9SHong Zhang   DM             dm,subdm,newdm;
5844b84eec9SHong Zhang   PetscErrorCode ierr;
5854b84eec9SHong Zhang 
5864b84eec9SHong Zhang   PetscFunctionBegin;
5874b84eec9SHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"slow",&mprk->subts_slow);CHKERRQ(ierr);
5884b84eec9SHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"fast",&mprk->subts_fast);CHKERRQ(ierr);
5894b84eec9SHong 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");
5904b84eec9SHong Zhang 
59119c2959aSHong Zhang   /* Only copy the DM */
5924b84eec9SHong Zhang   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
59319c2959aSHong Zhang 
59419c2959aSHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"slowbuffer",&mprk->subts_slowbuffer);CHKERRQ(ierr);
59519c2959aSHong Zhang   if (!mprk->subts_slowbuffer) {
59619c2959aSHong Zhang     mprk->subts_slowbuffer = mprk->subts_slow;
59719c2959aSHong Zhang     mprk->subts_slow       = NULL;
59819c2959aSHong Zhang   }
59919c2959aSHong Zhang   if (mprk->subts_slow) {
6004b84eec9SHong Zhang     ierr = DMClone(dm,&newdm);CHKERRQ(ierr);
6014b84eec9SHong Zhang     ierr = TSGetDM(mprk->subts_slow,&subdm);CHKERRQ(ierr);
6024b84eec9SHong Zhang     ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr);
6034b84eec9SHong Zhang     ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr);
6044b84eec9SHong Zhang     ierr = TSSetDM(mprk->subts_slow,newdm);CHKERRQ(ierr);
6054b84eec9SHong Zhang     ierr = DMDestroy(&newdm);CHKERRQ(ierr);
60619c2959aSHong Zhang   }
60719c2959aSHong Zhang   ierr = DMClone(dm,&newdm);CHKERRQ(ierr);
60819c2959aSHong Zhang   ierr = TSGetDM(mprk->subts_slowbuffer,&subdm);CHKERRQ(ierr);
60919c2959aSHong Zhang   ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr);
61019c2959aSHong Zhang   ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr);
61119c2959aSHong Zhang   ierr = TSSetDM(mprk->subts_slowbuffer,newdm);CHKERRQ(ierr);
61219c2959aSHong Zhang   ierr = DMDestroy(&newdm);CHKERRQ(ierr);
61319c2959aSHong Zhang 
61419c2959aSHong Zhang   ierr = DMClone(dm,&newdm);CHKERRQ(ierr);
61519c2959aSHong Zhang   ierr = TSGetDM(mprk->subts_fast,&subdm);CHKERRQ(ierr);
61619c2959aSHong Zhang   ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr);
61719c2959aSHong Zhang   ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr);
61819c2959aSHong Zhang   ierr = TSSetDM(mprk->subts_fast,newdm);CHKERRQ(ierr);
61919c2959aSHong Zhang   ierr = DMDestroy(&newdm);CHKERRQ(ierr);
62019c2959aSHong Zhang 
62119c2959aSHong Zhang   if (tab->np == 3) {
62219c2959aSHong Zhang     ierr = TSRHSSplitGetSubTS(ts,"medium",&mprk->subts_medium);CHKERRQ(ierr);
62319c2959aSHong Zhang     ierr = TSRHSSplitGetSubTS(ts,"mediumbuffer",&mprk->subts_mediumbuffer);CHKERRQ(ierr);
62419c2959aSHong Zhang     if (mprk->subts_medium && !mprk->subts_mediumbuffer) {
62519c2959aSHong Zhang       mprk->subts_mediumbuffer = mprk->subts_medium;
62619c2959aSHong Zhang       mprk->subts_medium       = NULL;
62719c2959aSHong Zhang     }
62819c2959aSHong Zhang     if (mprk->subts_medium) {
62919c2959aSHong Zhang       ierr = DMClone(dm,&newdm);CHKERRQ(ierr);
63019c2959aSHong Zhang       ierr = TSGetDM(mprk->subts_medium,&subdm);CHKERRQ(ierr);
63119c2959aSHong Zhang       ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr);
63219c2959aSHong Zhang       ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr);
63319c2959aSHong Zhang       ierr = TSSetDM(mprk->subts_medium,newdm);CHKERRQ(ierr);
63419c2959aSHong Zhang       ierr = DMDestroy(&newdm);CHKERRQ(ierr);
63519c2959aSHong Zhang     }
63619c2959aSHong Zhang     ierr = DMClone(dm,&newdm);CHKERRQ(ierr);
63719c2959aSHong Zhang     ierr = TSGetDM(mprk->subts_mediumbuffer,&subdm);CHKERRQ(ierr);
63819c2959aSHong Zhang     ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr);
63919c2959aSHong Zhang     ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr);
64019c2959aSHong Zhang     ierr = TSSetDM(mprk->subts_mediumbuffer,newdm);CHKERRQ(ierr);
64119c2959aSHong Zhang     ierr = DMDestroy(&newdm);CHKERRQ(ierr);
64219c2959aSHong Zhang   }
6434b84eec9SHong Zhang   PetscFunctionReturn(0);
6444b84eec9SHong Zhang }
6454b84eec9SHong Zhang 
6464b84eec9SHong Zhang /*
6474b84eec9SHong Zhang  This if for nonsplit RHS MPRK
6484b84eec9SHong Zhang  The step completion formula is
6494b84eec9SHong Zhang 
6504b84eec9SHong Zhang  x1 = x0 + h b^T YdotRHS
6514b84eec9SHong Zhang 
6524b84eec9SHong Zhang */
6534b84eec9SHong Zhang static PetscErrorCode TSEvaluateStep_MPRK(TS ts,PetscInt order,Vec X,PetscBool *done)
6544b84eec9SHong Zhang {
6554b84eec9SHong Zhang   TS_MPRK        *mprk = (TS_MPRK*)ts->data;
6564b84eec9SHong Zhang   MPRKTableau    tab = mprk->tableau;
6577c0df07dSHong Zhang   PetscScalar    *wf = mprk->work_fast;
6584b84eec9SHong Zhang   PetscReal      h = ts->time_step;
6594b84eec9SHong Zhang   PetscInt       s = tab->s,j;
6604b84eec9SHong Zhang   PetscErrorCode ierr;
6614b84eec9SHong Zhang 
6624b84eec9SHong Zhang   PetscFunctionBegin;
6634b84eec9SHong Zhang   for (j=0; j<s; j++) wf[j] = h*tab->bf[j];
6647c0df07dSHong Zhang   ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
6657c0df07dSHong Zhang   ierr = VecMAXPY(X,s,wf,mprk->YdotRHS);CHKERRQ(ierr);
6664b84eec9SHong Zhang   PetscFunctionReturn(0);
6674b84eec9SHong Zhang }
6684b84eec9SHong Zhang 
6694b84eec9SHong Zhang static PetscErrorCode TSStep_MPRK(TS ts)
6704b84eec9SHong Zhang {
6714b84eec9SHong Zhang   TS_MPRK         *mprk = (TS_MPRK*)ts->data;
6727c0df07dSHong Zhang   Vec             *Y = mprk->Y,*YdotRHS = mprk->YdotRHS,*YdotRHS_fast = mprk->YdotRHS_fast,*YdotRHS_slow = mprk->YdotRHS_slow,*YdotRHS_slowbuffer = mprk->YdotRHS_slowbuffer;
6737c0df07dSHong Zhang   Vec             Yslow,Yslowbuffer,Yfast;
6744b84eec9SHong Zhang   MPRKTableau     tab = mprk->tableau;
6754b84eec9SHong Zhang   const PetscInt  s = tab->s;
67619c2959aSHong Zhang   const PetscReal *Af = tab->Af,*cf = tab->cf,*Asb = tab->Asb,*csb = tab->csb;
6777c0df07dSHong Zhang   PetscScalar     *wf = mprk->work_fast,*ws = mprk->work_slow,*wsb = mprk->work_slowbuffer;
678*79bc8a77SHong Zhang   PetscInt        i,j,computedstages;
6794b84eec9SHong Zhang   PetscReal       next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step;
6804b84eec9SHong Zhang   PetscErrorCode  ierr;
6814b84eec9SHong Zhang 
6824b84eec9SHong Zhang   PetscFunctionBegin;
6834b84eec9SHong Zhang   for (i=0; i<s; i++) {
6844b84eec9SHong Zhang     mprk->stage_time = t + h*cf[i];
6854b84eec9SHong Zhang     ierr = TSPreStage(ts,mprk->stage_time);CHKERRQ(ierr);
6867c0df07dSHong Zhang     ierr = VecCopy(ts->vec_sol, Y[i]);CHKERRQ(ierr);
6874b84eec9SHong Zhang 
6887c0df07dSHong Zhang     /* slow buffer region */
68919c2959aSHong Zhang     for (j=0; j<i; j++) wsb[j] = h*Asb[i*s+j];
6909849be05SHong Zhang     for (j=0; j<i; j++) {
6917c0df07dSHong Zhang       ierr = VecGetSubVector(YdotRHS[j],mprk->is_slowbuffer,&YdotRHS_slowbuffer[j]);CHKERRQ(ierr);
6927c0df07dSHong Zhang     }
6937c0df07dSHong Zhang     ierr = VecGetSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr);
6947c0df07dSHong Zhang     ierr = VecMAXPY(Yslowbuffer,i,wsb,mprk->YdotRHS_slowbuffer);CHKERRQ(ierr);
6957c0df07dSHong Zhang     ierr = VecRestoreSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr);
6969849be05SHong Zhang     for (j=0; j<i; j++) {
6977c0df07dSHong Zhang       ierr = VecRestoreSubVector(YdotRHS[j],mprk->is_slowbuffer,&YdotRHS_slowbuffer[j]);CHKERRQ(ierr);
6987c0df07dSHong Zhang     }
6997c0df07dSHong Zhang     /* slow region */
7007c0df07dSHong Zhang     if (mprk->is_slow) {
701*79bc8a77SHong Zhang       computedstages = 0;
7029849be05SHong Zhang       for (j=0; j<i; j++) {
7039849be05SHong Zhang         if (tab->rsb[j]) ws[tab->rsb[j]-1] += h*tab->Asb[i*s+j];
704*79bc8a77SHong Zhang         else ws[computedstages++] = h*tab->Asb[i*s+j];
7057c0df07dSHong Zhang       }
706*79bc8a77SHong Zhang       for (j=0; j<computedstages; j++) {
7077c0df07dSHong Zhang         ierr = VecGetSubVector(YdotRHS[j],mprk->is_slow,&YdotRHS_slow[j]);CHKERRQ(ierr);
7087c0df07dSHong Zhang       }
7097c0df07dSHong Zhang       ierr = VecGetSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr);
710*79bc8a77SHong Zhang       ierr = VecMAXPY(Yslow,computedstages,ws,mprk->YdotRHS_slow);CHKERRQ(ierr);
7117c0df07dSHong Zhang       ierr = VecRestoreSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr);
712*79bc8a77SHong Zhang       for (j=0; j<computedstages; j++) {
7137c0df07dSHong Zhang         ierr = VecRestoreSubVector(YdotRHS[j],mprk->is_slow,&YdotRHS_slow[j]);CHKERRQ(ierr);
7147c0df07dSHong Zhang       }
7157c0df07dSHong Zhang     }
7164b84eec9SHong Zhang 
7177c0df07dSHong Zhang     /* fast region */
7184b84eec9SHong Zhang     for (j=0; j<i; j++) wf[j] = h*Af[i*s+j];
7199849be05SHong Zhang     for (j=0; j<i; j++) {
7207c0df07dSHong Zhang       ierr = VecGetSubVector(YdotRHS[j],mprk->is_fast,&YdotRHS_fast[j]);CHKERRQ(ierr);
7217c0df07dSHong Zhang     }
7227c0df07dSHong Zhang     ierr = VecGetSubVector(Y[i],mprk->is_fast,&Yfast);CHKERRQ(ierr);
7237c0df07dSHong Zhang     ierr = VecMAXPY(Yfast,i,wf,mprk->YdotRHS_fast);CHKERRQ(ierr);
7247c0df07dSHong Zhang     ierr = VecRestoreSubVector(Y[i],mprk->is_fast,&Yfast);CHKERRQ(ierr);
7259849be05SHong Zhang     for (j=0; j<i; j++) {
7267c0df07dSHong Zhang       ierr = VecRestoreSubVector(YdotRHS[j],mprk->is_fast,&YdotRHS_fast[j]);CHKERRQ(ierr);
7277c0df07dSHong Zhang     }
7287c0df07dSHong Zhang     if (tab->np == 3) {
7297c0df07dSHong Zhang       Vec         *YdotRHS_medium = mprk->YdotRHS_medium,*YdotRHS_mediumbuffer = mprk->YdotRHS_mediumbuffer;
7307c0df07dSHong Zhang       Vec         Ymedium,Ymediumbuffer;
7317c0df07dSHong Zhang       PetscScalar *wm = mprk->work_medium,*wmb = mprk->work_mediumbuffer;
7327c0df07dSHong Zhang       /* medium region */
7337c0df07dSHong Zhang       if (mprk->is_medium) {
734*79bc8a77SHong Zhang         computedstages = 0;
7357c0df07dSHong Zhang         for (j=0; j<i; j++) {
736*79bc8a77SHong Zhang           if (tab->rmb[j]) wm[computedstages-tab->sbase+(tab->rmb[j]-1)%tab->sbase] += h*tab->Amb[i*s+j];
737*79bc8a77SHong Zhang           else wm[computedstages++] = h*tab->Amb[i*s+j];
7387c0df07dSHong Zhang         }
739*79bc8a77SHong Zhang         for (j=0; j<computedstages; j++) {
7407c0df07dSHong Zhang           ierr = VecGetSubVector(YdotRHS[j],mprk->is_medium,&YdotRHS_medium[j]);CHKERRQ(ierr);
7417c0df07dSHong Zhang         }
7427c0df07dSHong Zhang         ierr = VecGetSubVector(Y[i],mprk->is_medium,&Ymedium);CHKERRQ(ierr);
743*79bc8a77SHong Zhang         ierr = VecMAXPY(Ymedium,computedstages,wm,mprk->YdotRHS_medium);CHKERRQ(ierr);
7447c0df07dSHong Zhang         ierr = VecRestoreSubVector(Y[i],mprk->is_medium,&Ymedium);CHKERRQ(ierr);
745*79bc8a77SHong Zhang         for (j=0; j<computedstages; j++) {
7467c0df07dSHong Zhang           ierr = VecRestoreSubVector(YdotRHS[j],mprk->is_medium,&YdotRHS_medium[j]);CHKERRQ(ierr);
7477c0df07dSHong Zhang         }
7487c0df07dSHong Zhang       }
7497c0df07dSHong Zhang       /* medium buffer region */
7507c0df07dSHong Zhang       for (j=0; j<i; j++) wmb[j] = h*tab->Amb[i*s+j];
7519849be05SHong Zhang       for (j=0; j<i; j++) {
7527c0df07dSHong Zhang         ierr = VecGetSubVector(YdotRHS[j],mprk->is_mediumbuffer,&YdotRHS_mediumbuffer[j]);CHKERRQ(ierr);
7537c0df07dSHong Zhang       }
7547c0df07dSHong Zhang       ierr = VecGetSubVector(Y[i],mprk->is_mediumbuffer,&Ymediumbuffer);CHKERRQ(ierr);
7557c0df07dSHong Zhang       ierr = VecMAXPY(Ymediumbuffer,i,wmb,mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr);
7567c0df07dSHong Zhang       ierr = VecRestoreSubVector(Y[i],mprk->is_mediumbuffer,&Ymediumbuffer);CHKERRQ(ierr);
7579849be05SHong Zhang       for (j=0; j<i; j++) {
7587c0df07dSHong Zhang         ierr = VecRestoreSubVector(YdotRHS[j],mprk->is_mediumbuffer,&YdotRHS_mediumbuffer[j]);CHKERRQ(ierr);
7597c0df07dSHong Zhang       }
7607c0df07dSHong Zhang     }
7614b84eec9SHong Zhang     ierr = TSPostStage(ts,mprk->stage_time,i,Y); CHKERRQ(ierr);
7624b84eec9SHong Zhang     /* compute the stage RHS by fast and slow tableau respectively */
7637c0df07dSHong Zhang     ierr = TSComputeRHSFunction(ts,t+h*csb[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
7644b84eec9SHong Zhang   }
7654b84eec9SHong Zhang   ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
7664b84eec9SHong Zhang   ts->ptime += ts->time_step;
7674b84eec9SHong Zhang   ts->time_step = next_time_step;
7684b84eec9SHong Zhang   PetscFunctionReturn(0);
7694b84eec9SHong Zhang }
7704b84eec9SHong Zhang 
7714b84eec9SHong Zhang /*
7724b84eec9SHong Zhang  This if for partitioned RHS MPRK
7734b84eec9SHong Zhang  The step completion formula is
7744b84eec9SHong Zhang 
7754b84eec9SHong Zhang  x1 = x0 + h b^T YdotRHS
7764b84eec9SHong Zhang 
7774b84eec9SHong Zhang */
7784b84eec9SHong Zhang static PetscErrorCode TSEvaluateStep_MPRKSPLIT(TS ts,PetscInt order,Vec X,PetscBool *done)
7794b84eec9SHong Zhang {
7804b84eec9SHong Zhang   TS_MPRK        *mprk = (TS_MPRK*)ts->data;
7814b84eec9SHong Zhang   MPRKTableau    tab  = mprk->tableau;
78219c2959aSHong Zhang   Vec            Xslow,Xfast,Xslowbuffer; /* subvectors for slow and fast componets in X respectively */
78319c2959aSHong Zhang   PetscScalar    *wf = mprk->work_fast,*ws = mprk->work_slow,*wsb = mprk->work_slowbuffer;
7844b84eec9SHong Zhang   PetscReal      h = ts->time_step;
785*79bc8a77SHong Zhang   PetscInt       s = tab->s,j,computedstages;
7864b84eec9SHong Zhang   PetscErrorCode ierr;
7874b84eec9SHong Zhang 
7884b84eec9SHong Zhang   PetscFunctionBegin;
7894b84eec9SHong Zhang   ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
79019c2959aSHong Zhang 
79119c2959aSHong Zhang   /* slow region */
79219c2959aSHong Zhang   if (mprk->is_slow) {
793*79bc8a77SHong Zhang     computedstages = 0;
79419c2959aSHong Zhang     for (j=0; j<s; j++) {
7959849be05SHong Zhang       if (tab->rsb[j]) ws[tab->rsb[j]-1] += h*tab->bsb[j];
796*79bc8a77SHong Zhang       else ws[computedstages++] = h*tab->bsb[j];
79719c2959aSHong Zhang     }
7984b84eec9SHong Zhang     ierr = VecGetSubVector(X,mprk->is_slow,&Xslow);CHKERRQ(ierr);
799*79bc8a77SHong Zhang     ierr = VecMAXPY(Xslow,computedstages,ws,mprk->YdotRHS_slow);CHKERRQ(ierr);
80019c2959aSHong Zhang     ierr = VecRestoreSubVector(X,mprk->is_slow,&Xslow);CHKERRQ(ierr);
80119c2959aSHong Zhang   }
80219c2959aSHong Zhang 
80319c2959aSHong Zhang   /* slow buffer region */
80419c2959aSHong Zhang   for (j=0; j<s; j++) wsb[j] = h*tab->bsb[j];
80519c2959aSHong Zhang   ierr = VecGetSubVector(X,mprk->is_slowbuffer,&Xslowbuffer);CHKERRQ(ierr);
80619c2959aSHong Zhang   ierr = VecMAXPY(Xslowbuffer,s,wsb,mprk->YdotRHS_slowbuffer);CHKERRQ(ierr);
80719c2959aSHong Zhang   ierr = VecRestoreSubVector(X,mprk->is_slowbuffer,&Xslowbuffer);CHKERRQ(ierr);
80819c2959aSHong Zhang 
80919c2959aSHong Zhang   if (tab->np == 3) {
81019c2959aSHong Zhang     Vec         Xmedium,Xmediumbuffer;
81119c2959aSHong Zhang     PetscScalar *wm = mprk->work_medium,*wmb = mprk->work_mediumbuffer;
81219c2959aSHong Zhang     /* medium region */
81319c2959aSHong Zhang     if (mprk->is_medium) {
814*79bc8a77SHong Zhang       computedstages = 0;
81519c2959aSHong Zhang       for (j=0; j<s; j++) {
816*79bc8a77SHong Zhang         if (tab->rmb[j]) wm[computedstages-tab->sbase+(tab->rmb[j]-1)%tab->sbase] += h*tab->bmb[j];
817*79bc8a77SHong Zhang         else wm[computedstages++] = h*tab->bmb[j];
81819c2959aSHong Zhang       }
81919c2959aSHong Zhang       ierr = VecGetSubVector(X,mprk->is_medium,&Xmedium);CHKERRQ(ierr);
820*79bc8a77SHong Zhang       ierr = VecMAXPY(Xmedium,computedstages,wm,mprk->YdotRHS_medium);CHKERRQ(ierr);
82119c2959aSHong Zhang       ierr = VecRestoreSubVector(X,mprk->is_medium,&Xmedium);CHKERRQ(ierr);
82219c2959aSHong Zhang     }
82319c2959aSHong Zhang     /* medium buffer region */
82419c2959aSHong Zhang     for (j=0; j<s; j++) wmb[j] = h*tab->bmb[j];
82519c2959aSHong Zhang     ierr = VecGetSubVector(X,mprk->is_mediumbuffer,&Xmediumbuffer);CHKERRQ(ierr);
82619c2959aSHong Zhang     ierr = VecMAXPY(Xmediumbuffer,s,wmb,mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr);
82719c2959aSHong Zhang     ierr = VecRestoreSubVector(X,mprk->is_mediumbuffer,&Xmediumbuffer);CHKERRQ(ierr);
82819c2959aSHong Zhang   }
82919c2959aSHong Zhang   /* fast region */
83019c2959aSHong Zhang   for (j=0; j<s; j++) wf[j] = h*tab->bf[j];
8314b84eec9SHong Zhang   ierr = VecGetSubVector(X,mprk->is_fast,&Xfast);CHKERRQ(ierr);
8324b84eec9SHong Zhang   ierr = VecMAXPY(Xfast,s,wf,mprk->YdotRHS_fast);CHKERRQ(ierr);
83319c2959aSHong Zhang   ierr = VecRestoreSubVector(X,mprk->is_fast,&Xfast);CHKERRQ(ierr);
8344b84eec9SHong Zhang   PetscFunctionReturn(0);
8354b84eec9SHong Zhang }
8364b84eec9SHong Zhang 
8374b84eec9SHong Zhang static PetscErrorCode TSStep_MPRKSPLIT(TS ts)
8384b84eec9SHong Zhang {
8394b84eec9SHong Zhang   TS_MPRK         *mprk = (TS_MPRK*)ts->data;
8404b84eec9SHong Zhang   MPRKTableau     tab = mprk->tableau;
84119c2959aSHong Zhang   Vec             *Y = mprk->Y,*YdotRHS_fast = mprk->YdotRHS_fast,*YdotRHS_slow = mprk->YdotRHS_slow,*YdotRHS_slowbuffer = mprk->YdotRHS_slowbuffer;
84219c2959aSHong Zhang   Vec             Yslow,Yslowbuffer,Yfast; /* subvectors for slow and fast components in Y[i] respectively */
84319c2959aSHong Zhang   PetscInt        s = tab->s;
84419c2959aSHong Zhang   const PetscReal *Af = tab->Af,*cf = tab->cf,*Asb = tab->Asb,*csb = tab->csb;
84519c2959aSHong Zhang   PetscScalar     *wf = mprk->work_fast,*ws = mprk->work_slow,*wsb = mprk->work_slowbuffer;
846*79bc8a77SHong Zhang   PetscInt        i,j,computedstages;
8474b84eec9SHong Zhang   PetscReal       next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step;
8484b84eec9SHong Zhang   PetscErrorCode  ierr;
8494b84eec9SHong Zhang 
8504b84eec9SHong Zhang   PetscFunctionBegin;
8514b84eec9SHong Zhang   for (i=0; i<s; i++) {
8524b84eec9SHong Zhang     mprk->stage_time = t + h*cf[i];
8534b84eec9SHong Zhang     ierr = TSPreStage(ts,mprk->stage_time);CHKERRQ(ierr);
8544b84eec9SHong Zhang     /* calculate the stage value for fast and slow components respectively */
8554b84eec9SHong Zhang     ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
85619c2959aSHong Zhang     for (j=0; j<i; j++) wsb[j] = h*Asb[i*s+j];
85719c2959aSHong Zhang 
85819c2959aSHong Zhang     /* slow buffer region */
85919c2959aSHong Zhang     ierr = VecGetSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr);
86019c2959aSHong Zhang     ierr = VecMAXPY(Yslowbuffer,i,wsb,YdotRHS_slowbuffer);CHKERRQ(ierr);
86119c2959aSHong Zhang     ierr = VecRestoreSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr);
8629849be05SHong Zhang 
86319c2959aSHong Zhang     /* slow region */
8649849be05SHong Zhang     if (mprk->is_slow) {
8659849be05SHong Zhang       if (tab->rsb[i]) { /* repeat previous stage */
8669849be05SHong Zhang         ierr = VecGetSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr);
8679849be05SHong Zhang         ierr = VecISCopy(Y[tab->rsb[i]-1],mprk->is_slow,SCATTER_REVERSE,Yslow);CHKERRQ(ierr);
8689849be05SHong Zhang         ierr = VecRestoreSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr);
8699849be05SHong Zhang       } else {
870*79bc8a77SHong Zhang         computedstages = 0;
87119c2959aSHong Zhang         for (j=0; j<i; j++) {
872*79bc8a77SHong Zhang           if (tab->rsb[j]) ws[tab->rsb[j]-1] += wsb[j];
873*79bc8a77SHong Zhang           else ws[computedstages++] = wsb[j];
87419c2959aSHong Zhang         }
8754b84eec9SHong Zhang         ierr = VecGetSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr);
876*79bc8a77SHong Zhang         ierr = VecMAXPY(Yslow,computedstages,ws,YdotRHS_slow);CHKERRQ(ierr);
8774b84eec9SHong Zhang         ierr = VecRestoreSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr);
87819c2959aSHong Zhang         /* only depends on the slow buffer region */
879*79bc8a77SHong Zhang         ierr = TSComputeRHSFunction(mprk->subts_slow,t+h*csb[i],Y[i],YdotRHS_slow[computedstages]);CHKERRQ(ierr);
88019c2959aSHong Zhang       }
8819849be05SHong Zhang     }
88219c2959aSHong Zhang 
88319c2959aSHong Zhang     /* fast region */
88419c2959aSHong Zhang     for (j=0; j<i; j++) wf[j] = h*Af[i*s+j];
88519c2959aSHong Zhang     ierr = VecGetSubVector(Y[i],mprk->is_fast,&Yfast);CHKERRQ(ierr);
88619c2959aSHong Zhang     ierr = VecMAXPY(Yfast,i,wf,YdotRHS_fast);CHKERRQ(ierr);
8874b84eec9SHong Zhang     ierr = VecRestoreSubVector(Y[i],mprk->is_fast,&Yfast);CHKERRQ(ierr);
88819c2959aSHong Zhang 
88919c2959aSHong Zhang     if (tab->np == 3) {
89019c2959aSHong Zhang       Vec *YdotRHS_medium = mprk->YdotRHS_medium,*YdotRHS_mediumbuffer = mprk->YdotRHS_mediumbuffer;
89119c2959aSHong Zhang       Vec Ymedium,Ymediumbuffer;
89219c2959aSHong Zhang       const PetscReal *Amb = tab->Amb,*cmb = tab->cmb;
89319c2959aSHong Zhang       PetscScalar *wm = mprk->work_medium,*wmb = mprk->work_mediumbuffer;
89419c2959aSHong Zhang 
89519c2959aSHong Zhang       for (j=0; j<i; j++) wmb[j] = h*Amb[i*s+j];
89619c2959aSHong Zhang       /* medium buffer region */
89719c2959aSHong Zhang       ierr = VecGetSubVector(Y[i],mprk->is_mediumbuffer,&Ymediumbuffer);CHKERRQ(ierr);
89819c2959aSHong Zhang       ierr = VecMAXPY(Ymediumbuffer,i,wmb,YdotRHS_mediumbuffer);CHKERRQ(ierr);
89919c2959aSHong Zhang       ierr = VecRestoreSubVector(Y[i],mprk->is_mediumbuffer,&Ymediumbuffer);CHKERRQ(ierr);
90019c2959aSHong Zhang       /* medium region */
901*79bc8a77SHong Zhang       if (mprk->is_medium) {
902*79bc8a77SHong Zhang         if (tab->rmb[i]) { /* repeat previous stage */
903*79bc8a77SHong Zhang           ierr = VecGetSubVector(Y[i],mprk->is_medium,&Ymedium);CHKERRQ(ierr);
904*79bc8a77SHong Zhang           ierr = VecISCopy(Y[tab->rmb[i]-1],mprk->is_medium,SCATTER_REVERSE,Ymedium);CHKERRQ(ierr);
90519c2959aSHong Zhang           ierr = VecRestoreSubVector(Y[i],mprk->is_medium,&Ymedium);CHKERRQ(ierr);
906*79bc8a77SHong Zhang         } else {
907*79bc8a77SHong Zhang           computedstages = 0;
908*79bc8a77SHong Zhang           for (j=0; j<i; j++) {
909*79bc8a77SHong Zhang             if (tab->rmb[j]) wm[computedstages-tab->sbase+(tab->rmb[j]-1)%tab->sbase] += wmb[j];
910*79bc8a77SHong Zhang             else wm[computedstages++] = wmb[j];
911*79bc8a77SHong Zhang 
912*79bc8a77SHong Zhang           }
913*79bc8a77SHong Zhang           ierr = VecGetSubVector(Y[i],mprk->is_medium,&Ymedium);CHKERRQ(ierr);
914*79bc8a77SHong Zhang           ierr = VecMAXPY(Ymedium,computedstages,wm,YdotRHS_medium);CHKERRQ(ierr);
915*79bc8a77SHong Zhang           ierr = VecRestoreSubVector(Y[i],mprk->is_medium,&Ymedium);CHKERRQ(ierr);
916*79bc8a77SHong Zhang           /* only depends on the medium buffer region */
917*79bc8a77SHong Zhang           ierr = TSComputeRHSFunction(mprk->subts_medium,t+h*cmb[i],Y[i],YdotRHS_medium[computedstages]);CHKERRQ(ierr);
918*79bc8a77SHong Zhang         }
91919c2959aSHong Zhang       }
92019c2959aSHong Zhang       /* must be computed after fast region and slow region are updated in Y */
92119c2959aSHong Zhang       ierr = TSComputeRHSFunction(mprk->subts_mediumbuffer,t+h*cmb[i],Y[i],YdotRHS_mediumbuffer[i]);CHKERRQ(ierr);
92219c2959aSHong Zhang     }
92319c2959aSHong Zhang     /* must be computed after all regions are updated in Y */
92419c2959aSHong Zhang     ierr = TSComputeRHSFunction(mprk->subts_slowbuffer,t+h*csb[i],Y[i],YdotRHS_slowbuffer[i]);CHKERRQ(ierr);
9257c0df07dSHong Zhang     ierr = TSComputeRHSFunction(mprk->subts_fast,t+h*cf[i],Y[i],YdotRHS_fast[i]);
9264b84eec9SHong Zhang   }
927*79bc8a77SHong Zhang 
9284b84eec9SHong Zhang   ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
9294b84eec9SHong Zhang   ts->ptime += ts->time_step;
9304b84eec9SHong Zhang   ts->time_step = next_time_step;
9314b84eec9SHong Zhang   PetscFunctionReturn(0);
9324b84eec9SHong Zhang }
9334b84eec9SHong Zhang 
9344b84eec9SHong Zhang static PetscErrorCode TSMPRKTableauReset(TS ts)
9354b84eec9SHong Zhang {
9364b84eec9SHong Zhang   TS_MPRK        *mprk = (TS_MPRK*)ts->data;
9374b84eec9SHong Zhang   MPRKTableau    tab = mprk->tableau;
9384b84eec9SHong Zhang   PetscErrorCode ierr;
9394b84eec9SHong Zhang 
9404b84eec9SHong Zhang   PetscFunctionBegin;
9414b84eec9SHong Zhang   if (!tab) PetscFunctionReturn(0);
9424b84eec9SHong Zhang   ierr = PetscFree(mprk->work_fast);CHKERRQ(ierr);
9434b84eec9SHong Zhang   ierr = PetscFree(mprk->work_slow);CHKERRQ(ierr);
9447c0df07dSHong Zhang   ierr = PetscFree(mprk->work_slowbuffer);CHKERRQ(ierr);
9457c0df07dSHong Zhang   ierr = PetscFree(mprk->work_medium);CHKERRQ(ierr);
9467c0df07dSHong Zhang   ierr = PetscFree(mprk->work_mediumbuffer);CHKERRQ(ierr);
9474b84eec9SHong Zhang   ierr = VecDestroyVecs(tab->s,&mprk->Y);CHKERRQ(ierr);
9484b84eec9SHong Zhang   if (mprk->mprkmtype == TSMPRKNONSPLIT) {
9497c0df07dSHong Zhang     ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS);CHKERRQ(ierr);
9507c0df07dSHong Zhang     if (mprk->is_slow) {
9517c0df07dSHong Zhang       ierr = PetscFree(mprk->YdotRHS_slow);CHKERRQ(ierr);
9524b84eec9SHong Zhang     }
9537c0df07dSHong Zhang     ierr = PetscFree(mprk->YdotRHS_slowbuffer);CHKERRQ(ierr);
9547c0df07dSHong Zhang     if (tab->np == 3) {
9557c0df07dSHong Zhang       if (mprk->is_medium) {
9567c0df07dSHong Zhang         ierr = PetscFree(mprk->YdotRHS_medium);CHKERRQ(ierr);
9577c0df07dSHong Zhang       }
9587c0df07dSHong Zhang       ierr = PetscFree(mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr);
9597c0df07dSHong Zhang     }
9607c0df07dSHong Zhang     ierr = PetscFree(mprk->YdotRHS_fast);CHKERRQ(ierr);
9617c0df07dSHong Zhang 
9627c0df07dSHong Zhang   } else {
9634b84eec9SHong Zhang     ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_fast);CHKERRQ(ierr);
9644b84eec9SHong Zhang     ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_slow);CHKERRQ(ierr);
9657c0df07dSHong Zhang     ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_slowbuffer);CHKERRQ(ierr);
9667c0df07dSHong Zhang     ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_medium);CHKERRQ(ierr);
9677c0df07dSHong Zhang     ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr);
9687c0df07dSHong Zhang   }
9694b84eec9SHong Zhang   PetscFunctionReturn(0);
9704b84eec9SHong Zhang }
9714b84eec9SHong Zhang 
9724b84eec9SHong Zhang static PetscErrorCode TSReset_MPRK(TS ts)
9734b84eec9SHong Zhang {
9744b84eec9SHong Zhang   PetscErrorCode ierr;
9754b84eec9SHong Zhang 
9764b84eec9SHong Zhang   PetscFunctionBegin;
9774b84eec9SHong Zhang   ierr = TSMPRKTableauReset(ts);CHKERRQ(ierr);
9784b84eec9SHong Zhang   PetscFunctionReturn(0);
9794b84eec9SHong Zhang }
9804b84eec9SHong Zhang 
9814b84eec9SHong Zhang static PetscErrorCode DMCoarsenHook_TSMPRK(DM fine,DM coarse,void *ctx)
9824b84eec9SHong Zhang {
9834b84eec9SHong Zhang   PetscFunctionBegin;
9844b84eec9SHong Zhang   PetscFunctionReturn(0);
9854b84eec9SHong Zhang }
9864b84eec9SHong Zhang 
9874b84eec9SHong Zhang static PetscErrorCode DMRestrictHook_TSMPRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
9884b84eec9SHong Zhang {
9894b84eec9SHong Zhang   PetscFunctionBegin;
9904b84eec9SHong Zhang   PetscFunctionReturn(0);
9914b84eec9SHong Zhang }
9924b84eec9SHong Zhang 
9934b84eec9SHong Zhang static PetscErrorCode DMSubDomainHook_TSMPRK(DM dm,DM subdm,void *ctx)
9944b84eec9SHong Zhang {
9954b84eec9SHong Zhang   PetscFunctionBegin;
9964b84eec9SHong Zhang   PetscFunctionReturn(0);
9974b84eec9SHong Zhang }
9984b84eec9SHong Zhang 
9994b84eec9SHong Zhang static PetscErrorCode DMSubDomainRestrictHook_TSMPRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
10004b84eec9SHong Zhang {
10014b84eec9SHong Zhang   PetscFunctionBegin;
10024b84eec9SHong Zhang   PetscFunctionReturn(0);
10034b84eec9SHong Zhang }
10044b84eec9SHong Zhang 
10054b84eec9SHong Zhang static PetscErrorCode TSMPRKTableauSetUp(TS ts)
10064b84eec9SHong Zhang {
10074b84eec9SHong Zhang   TS_MPRK        *mprk  = (TS_MPRK*)ts->data;
10084b84eec9SHong Zhang   MPRKTableau    tab = mprk->tableau;
100919c2959aSHong Zhang   Vec            YdotRHS_slow,YdotRHS_slowbuffer,YdotRHS_medium,YdotRHS_mediumbuffer,YdotRHS_fast;
10104b84eec9SHong Zhang   PetscErrorCode ierr;
10114b84eec9SHong Zhang 
10124b84eec9SHong Zhang   PetscFunctionBegin;
10134b84eec9SHong Zhang   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&mprk->Y);CHKERRQ(ierr);
10147c0df07dSHong Zhang   if (mprk->is_slow) {
101519c2959aSHong Zhang     ierr = PetscMalloc1(tab->s,&mprk->work_slow);CHKERRQ(ierr);
10164b84eec9SHong Zhang   }
10177c0df07dSHong Zhang   ierr = PetscMalloc1(tab->s,&mprk->work_slowbuffer);CHKERRQ(ierr);
10187c0df07dSHong Zhang   if (tab->np == 3) {
10197c0df07dSHong Zhang     if (mprk->is_medium) {
10207c0df07dSHong Zhang       ierr = PetscMalloc1(tab->s,&mprk->work_medium);CHKERRQ(ierr);
10217c0df07dSHong Zhang     }
10227c0df07dSHong Zhang     ierr = PetscMalloc1(tab->s,&mprk->work_mediumbuffer);CHKERRQ(ierr);
10237c0df07dSHong Zhang   }
10247c0df07dSHong Zhang   ierr = PetscMalloc1(tab->s,&mprk->work_fast);CHKERRQ(ierr);
10257c0df07dSHong Zhang 
10267c0df07dSHong Zhang   if (mprk->mprkmtype == TSMPRKNONSPLIT) {
10277c0df07dSHong Zhang     ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&mprk->YdotRHS);CHKERRQ(ierr);
10287c0df07dSHong Zhang     if (mprk->is_slow) {
10297c0df07dSHong Zhang       ierr = PetscMalloc1(tab->s,&mprk->YdotRHS_slow);CHKERRQ(ierr);
10307c0df07dSHong Zhang     }
10317c0df07dSHong Zhang     ierr = PetscMalloc1(tab->s,&mprk->YdotRHS_slowbuffer);CHKERRQ(ierr);
10327c0df07dSHong Zhang     if (tab->np == 3) {
10337c0df07dSHong Zhang       if (mprk->is_medium) {
10347c0df07dSHong Zhang         ierr = PetscMalloc1(tab->s,&mprk->YdotRHS_medium);CHKERRQ(ierr);
10357c0df07dSHong Zhang       }
10367c0df07dSHong Zhang       ierr = PetscMalloc1(tab->s,&mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr);
10377c0df07dSHong Zhang     }
10387c0df07dSHong Zhang     ierr = PetscMalloc1(tab->s,&mprk->YdotRHS_fast);CHKERRQ(ierr);
10397c0df07dSHong Zhang   } else {
104019c2959aSHong Zhang     if (mprk->is_slow) {
10414b84eec9SHong Zhang       ierr = VecGetSubVector(ts->vec_sol,mprk->is_slow,&YdotRHS_slow);CHKERRQ(ierr);
10424b84eec9SHong Zhang       ierr = VecDuplicateVecs(YdotRHS_slow,tab->s,&mprk->YdotRHS_slow);CHKERRQ(ierr);
10434b84eec9SHong Zhang       ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_slow,&YdotRHS_slow);CHKERRQ(ierr);
104419c2959aSHong Zhang     }
104519c2959aSHong Zhang     ierr = VecGetSubVector(ts->vec_sol,mprk->is_slowbuffer,&YdotRHS_slowbuffer);CHKERRQ(ierr);
104619c2959aSHong Zhang     ierr = VecDuplicateVecs(YdotRHS_slowbuffer,tab->s,&mprk->YdotRHS_slowbuffer);CHKERRQ(ierr);
104719c2959aSHong Zhang     ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_slowbuffer,&YdotRHS_slowbuffer);CHKERRQ(ierr);
104819c2959aSHong Zhang     if (tab->np == 3) {
104919c2959aSHong Zhang       if (mprk->is_medium) {
105019c2959aSHong Zhang         ierr = VecGetSubVector(ts->vec_sol,mprk->is_medium,&YdotRHS_medium);CHKERRQ(ierr);
105119c2959aSHong Zhang         ierr = VecDuplicateVecs(YdotRHS_medium,tab->s,&mprk->YdotRHS_medium);CHKERRQ(ierr);
105219c2959aSHong Zhang         ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_medium,&YdotRHS_medium);CHKERRQ(ierr);
105319c2959aSHong Zhang       }
105419c2959aSHong Zhang       ierr = VecGetSubVector(ts->vec_sol,mprk->is_mediumbuffer,&YdotRHS_mediumbuffer);CHKERRQ(ierr);
105519c2959aSHong Zhang       ierr = VecDuplicateVecs(YdotRHS_mediumbuffer,tab->s,&mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr);
105619c2959aSHong Zhang       ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_mediumbuffer,&YdotRHS_mediumbuffer);CHKERRQ(ierr);
105719c2959aSHong Zhang     }
105819c2959aSHong Zhang     ierr = VecGetSubVector(ts->vec_sol,mprk->is_fast,&YdotRHS_fast);CHKERRQ(ierr);
105919c2959aSHong Zhang     ierr = VecDuplicateVecs(YdotRHS_fast,tab->s,&mprk->YdotRHS_fast);CHKERRQ(ierr);
10604b84eec9SHong Zhang     ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_fast,&YdotRHS_fast);CHKERRQ(ierr);
10614b84eec9SHong Zhang   }
10624b84eec9SHong Zhang   PetscFunctionReturn(0);
10634b84eec9SHong Zhang }
10644b84eec9SHong Zhang 
10654b84eec9SHong Zhang static PetscErrorCode TSSetUp_MPRK(TS ts)
10664b84eec9SHong Zhang {
10674b84eec9SHong Zhang   TS_MPRK        *mprk = (TS_MPRK*)ts->data;
106819c2959aSHong Zhang   MPRKTableau    tab = mprk->tableau;
10694b84eec9SHong Zhang   DM             dm;
10704b84eec9SHong Zhang   PetscErrorCode ierr;
10714b84eec9SHong Zhang 
10724b84eec9SHong Zhang   PetscFunctionBegin;
10734b84eec9SHong Zhang   ierr = TSRHSSplitGetIS(ts,"slow",&mprk->is_slow);CHKERRQ(ierr);
10744b84eec9SHong Zhang   ierr = TSRHSSplitGetIS(ts,"fast",&mprk->is_fast);CHKERRQ(ierr);
107519c2959aSHong 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);
107619c2959aSHong Zhang 
107719c2959aSHong Zhang   if (tab->np == 3) {
107819c2959aSHong Zhang     ierr = TSRHSSplitGetIS(ts,"medium",&mprk->is_medium);CHKERRQ(ierr);
107919c2959aSHong 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);
108019c2959aSHong Zhang     ierr = TSRHSSplitGetIS(ts,"mediumbuffer",&mprk->is_mediumbuffer);CHKERRQ(ierr);
108119c2959aSHong Zhang     if (!mprk->is_mediumbuffer) { /* let medium buffer cover whole medium region */
108219c2959aSHong Zhang       mprk->is_mediumbuffer = mprk->is_medium;
108319c2959aSHong Zhang       mprk->is_medium = NULL;
108419c2959aSHong Zhang     }
108519c2959aSHong Zhang   }
108619c2959aSHong Zhang 
108719c2959aSHong Zhang   /* If users do not provide buffer region settings, the solver will do them automatically, but with a performance penalty */
108819c2959aSHong Zhang   ierr = TSRHSSplitGetIS(ts,"slowbuffer",&mprk->is_slowbuffer);CHKERRQ(ierr);
108919c2959aSHong Zhang   if (!mprk->is_slowbuffer) { /* let slow buffer cover whole slow region */
109019c2959aSHong Zhang     mprk->is_slowbuffer = mprk->is_slow;
109119c2959aSHong Zhang     mprk->is_slow = NULL;
109219c2959aSHong Zhang   }
109319c2959aSHong Zhang /*
109419c2959aSHong Zhang   if (!mprk->is_medium) {
109519c2959aSHong Zhang     mprk->is_medium = mprk->is_fast;
109619c2959aSHong Zhang     mprk->is_fast = NULL;
109719c2959aSHong Zhang   } else {
109819c2959aSHong Zhang     ierr = TSRHSSplitGetIS(ts,"mediumbuffer",&mprk->is_mediumbuffer);CHKERRQ(ierr);
109919c2959aSHong Zhang   }
110019c2959aSHong Zhang */
11014b84eec9SHong Zhang   ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr);
11024b84eec9SHong Zhang   ierr = TSMPRKTableauSetUp(ts);CHKERRQ(ierr);
11034b84eec9SHong Zhang   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
11044b84eec9SHong Zhang   ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSMPRK,DMRestrictHook_TSMPRK,ts);CHKERRQ(ierr);
11054b84eec9SHong Zhang   ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSMPRK,DMSubDomainRestrictHook_TSMPRK,ts);CHKERRQ(ierr);
11064b84eec9SHong Zhang   ierr = PetscTryMethod(ts,"TSMPRKSetSplits_C",(TS),(ts));CHKERRQ(ierr);
11074b84eec9SHong Zhang   PetscFunctionReturn(0);
11084b84eec9SHong Zhang }
11094b84eec9SHong Zhang 
11104b84eec9SHong Zhang /* construct a database to chose nonsplit RHS mutirate mprk method or split RHS MPRK method */
11114b84eec9SHong Zhang const char *const TSMPRKMultirateTypes[] = {"NONSPLIT","SPLIT","TSMPRKMultirateType","TSMPRK",0};
11124b84eec9SHong Zhang 
11134b84eec9SHong Zhang static PetscErrorCode TSSetFromOptions_MPRK(PetscOptionItems *PetscOptionsObject,TS ts)
11144b84eec9SHong Zhang {
11154b84eec9SHong Zhang   TS_MPRK        *mprk = (TS_MPRK*)ts->data;
11164b84eec9SHong Zhang   PetscErrorCode ierr;
11174b84eec9SHong Zhang 
11184b84eec9SHong Zhang   PetscFunctionBegin;
11194b84eec9SHong Zhang   ierr = PetscOptionsHead(PetscOptionsObject,"PRK ODE solver options");CHKERRQ(ierr);
11204b84eec9SHong Zhang   {
11214b84eec9SHong Zhang     MPRKTableauLink link;
11224b84eec9SHong Zhang     PetscInt        count,choice;
11234b84eec9SHong Zhang     PetscBool       flg;
11244b84eec9SHong Zhang     const char      **namelist;
11254b84eec9SHong Zhang     PetscInt        mprkmtype = 0;
11264b84eec9SHong Zhang     for (link=MPRKTableauList,count=0; link; link=link->next,count++) ;
11274b84eec9SHong Zhang     ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr);
11284b84eec9SHong Zhang     for (link=MPRKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
11294b84eec9SHong Zhang     ierr = PetscOptionsEList("-ts_mprk_type","Family of MPRK method","TSMPRKSetType",(const char*const*)namelist,count,mprk->tableau->name,&choice,&flg);CHKERRQ(ierr);
11304b84eec9SHong Zhang     if (flg) {ierr = TSMPRKSetType(ts,namelist[choice]);CHKERRQ(ierr);}
11314b84eec9SHong Zhang     ierr = PetscFree(namelist);CHKERRQ(ierr);
11324b84eec9SHong 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);
11334b84eec9SHong Zhang      if (flg) {ierr = TSMPRKSetMultirateType(ts,mprkmtype);CHKERRQ(ierr);}
11344b84eec9SHong Zhang   }
11354b84eec9SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
11364b84eec9SHong Zhang   PetscFunctionReturn(0);
11374b84eec9SHong Zhang }
11384b84eec9SHong Zhang 
11394b84eec9SHong Zhang static PetscErrorCode TSView_MPRK(TS ts,PetscViewer viewer)
11404b84eec9SHong Zhang {
11414b84eec9SHong Zhang   TS_MPRK        *mprk = (TS_MPRK*)ts->data;
11424b84eec9SHong Zhang   PetscBool      iascii;
11434b84eec9SHong Zhang   PetscErrorCode ierr;
11444b84eec9SHong Zhang 
11454b84eec9SHong Zhang   PetscFunctionBegin;
11464b84eec9SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
11474b84eec9SHong Zhang   if (iascii) {
11484b84eec9SHong Zhang     MPRKTableau tab  = mprk->tableau;
11494b84eec9SHong Zhang     TSMPRKType  mprktype;
11504b84eec9SHong Zhang     char        fbuf[512];
11514b84eec9SHong Zhang     char        sbuf[512];
115219c2959aSHong Zhang     PetscInt    i;
11534b84eec9SHong Zhang     ierr = TSMPRKGetType(ts,&mprktype);CHKERRQ(ierr);
11544b84eec9SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  MPRK type %s\n",mprktype);CHKERRQ(ierr);
11554b84eec9SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Order: %D\n",tab->order);CHKERRQ(ierr);
115619c2959aSHong Zhang 
11574b84eec9SHong Zhang     ierr = PetscFormatRealArray(fbuf,sizeof(fbuf),"% 8.6f",tab->s,tab->cf);CHKERRQ(ierr);
11584b84eec9SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa cf = %s\n",fbuf);CHKERRQ(ierr);
115919c2959aSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Af = \n");CHKERRQ(ierr);
116019c2959aSHong Zhang     for (i=0; i<tab->s; i++) {
116119c2959aSHong Zhang       ierr = PetscFormatRealArray(fbuf,sizeof(fbuf),"% 8.6f",tab->s,&tab->Af[i*tab->s]);CHKERRQ(ierr);
116219c2959aSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"    %s\n",fbuf);CHKERRQ(ierr);
116319c2959aSHong Zhang     }
116419c2959aSHong Zhang     ierr = PetscFormatRealArray(fbuf,sizeof(fbuf),"% 8.6f",tab->s,tab->bf);CHKERRQ(ierr);
116519c2959aSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  bf = %s\n",fbuf);CHKERRQ(ierr);
116619c2959aSHong Zhang 
116719c2959aSHong Zhang     ierr = PetscFormatRealArray(sbuf,sizeof(sbuf),"% 8.6f",tab->s,tab->csb);CHKERRQ(ierr);
116819c2959aSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa csb = %s\n",sbuf);CHKERRQ(ierr);
116919c2959aSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Asb = \n");CHKERRQ(ierr);
117019c2959aSHong Zhang     for (i=0; i<tab->s; i++) {
117119c2959aSHong Zhang       ierr = PetscFormatRealArray(sbuf,sizeof(sbuf),"% 8.6f",tab->s,&tab->Asb[i*tab->s]);CHKERRQ(ierr);
117219c2959aSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"    %s\n",sbuf);CHKERRQ(ierr);
117319c2959aSHong Zhang     }
117419c2959aSHong Zhang     ierr = PetscFormatRealArray(sbuf,sizeof(sbuf),"% 8.6f",tab->s,tab->bsb);CHKERRQ(ierr);
117519c2959aSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  bsb = %s\n",sbuf);CHKERRQ(ierr);
117619c2959aSHong Zhang 
117719c2959aSHong Zhang     if (tab->np == 3) {
117819c2959aSHong Zhang       char mbuf[512];
117919c2959aSHong Zhang       ierr = PetscFormatRealArray(mbuf,sizeof(mbuf),"% 8.6f",tab->s,tab->cmb);CHKERRQ(ierr);
118019c2959aSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa cmb = %s\n",mbuf);CHKERRQ(ierr);
118119c2959aSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  Amb = \n");CHKERRQ(ierr);
118219c2959aSHong Zhang       for (i=0; i<tab->s; i++) {
118319c2959aSHong Zhang         ierr = PetscFormatRealArray(mbuf,sizeof(mbuf),"% 8.6f",tab->s,&tab->Amb[i*tab->s]);CHKERRQ(ierr);
118419c2959aSHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    %s\n",mbuf);CHKERRQ(ierr);
118519c2959aSHong Zhang       }
118619c2959aSHong Zhang       ierr = PetscFormatRealArray(mbuf,sizeof(mbuf),"% 8.6f",tab->s,tab->bmb);CHKERRQ(ierr);
118719c2959aSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  bmb = %s\n",mbuf);CHKERRQ(ierr);
118819c2959aSHong Zhang     }
11894b84eec9SHong Zhang   }
11904b84eec9SHong Zhang   PetscFunctionReturn(0);
11914b84eec9SHong Zhang }
11924b84eec9SHong Zhang 
11934b84eec9SHong Zhang static PetscErrorCode TSLoad_MPRK(TS ts,PetscViewer viewer)
11944b84eec9SHong Zhang {
11954b84eec9SHong Zhang   PetscErrorCode ierr;
11964b84eec9SHong Zhang   TSAdapt        adapt;
11974b84eec9SHong Zhang 
11984b84eec9SHong Zhang   PetscFunctionBegin;
11994b84eec9SHong Zhang   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
12004b84eec9SHong Zhang   ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr);
12014b84eec9SHong Zhang   PetscFunctionReturn(0);
12024b84eec9SHong Zhang }
12034b84eec9SHong Zhang 
12044b84eec9SHong Zhang /*@C
12054b84eec9SHong Zhang   TSMPRKSetType - Set the type of MPRK scheme
12064b84eec9SHong Zhang 
12074b84eec9SHong Zhang   Logically collective
12084b84eec9SHong Zhang 
12094b84eec9SHong Zhang   Input Parameter:
12104b84eec9SHong Zhang +  ts - timestepping context
12114b84eec9SHong Zhang -  mprktype - type of MPRK-scheme
12124b84eec9SHong Zhang 
12134b84eec9SHong Zhang   Options Database:
12144b84eec9SHong Zhang .   -ts_mprk_type - <pm2,p2,p3>
12154b84eec9SHong Zhang 
12164b84eec9SHong Zhang   Level: intermediate
12174b84eec9SHong Zhang 
12184b84eec9SHong Zhang .seealso: TSMPRKGetType(), TSMPRK, TSMPRKType
12194b84eec9SHong Zhang @*/
12204b84eec9SHong Zhang PetscErrorCode TSMPRKSetType(TS ts,TSMPRKType mprktype)
12214b84eec9SHong Zhang {
12224b84eec9SHong Zhang   PetscErrorCode ierr;
12234b84eec9SHong Zhang 
12244b84eec9SHong Zhang   PetscFunctionBegin;
12254b84eec9SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
12264b84eec9SHong Zhang   PetscValidCharPointer(mprktype,2);
12274b84eec9SHong Zhang   ierr = PetscTryMethod(ts,"TSMPRKSetType_C",(TS,TSMPRKType),(ts,mprktype));CHKERRQ(ierr);
12284b84eec9SHong Zhang   PetscFunctionReturn(0);
12294b84eec9SHong Zhang }
12304b84eec9SHong Zhang 
12314b84eec9SHong Zhang /*@C
12324b84eec9SHong Zhang   TSMPRKGetType - Get the type of MPRK scheme
12334b84eec9SHong Zhang 
12344b84eec9SHong Zhang   Logically collective
12354b84eec9SHong Zhang 
12364b84eec9SHong Zhang   Input Parameter:
12374b84eec9SHong Zhang .  ts - timestepping context
12384b84eec9SHong Zhang 
12394b84eec9SHong Zhang   Output Parameter:
12404b84eec9SHong Zhang .  mprktype - type of MPRK-scheme
12414b84eec9SHong Zhang 
12424b84eec9SHong Zhang   Level: intermediate
12434b84eec9SHong Zhang 
12444b84eec9SHong Zhang .seealso: TSMPRKGetType()
12454b84eec9SHong Zhang @*/
12464b84eec9SHong Zhang PetscErrorCode TSMPRKGetType(TS ts,TSMPRKType *mprktype)
12474b84eec9SHong Zhang {
12484b84eec9SHong Zhang   PetscErrorCode ierr;
12494b84eec9SHong Zhang 
12504b84eec9SHong Zhang   PetscFunctionBegin;
12514b84eec9SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
12524b84eec9SHong Zhang   ierr = PetscUseMethod(ts,"TSMPRKGetType_C",(TS,TSMPRKType*),(ts,mprktype));CHKERRQ(ierr);
12534b84eec9SHong Zhang   PetscFunctionReturn(0);
12544b84eec9SHong Zhang }
12554b84eec9SHong Zhang 
12564b84eec9SHong Zhang /*@C
12574b84eec9SHong Zhang   TSMPRKSetMultirateType - Set the type of MPRK multirate scheme
12584b84eec9SHong Zhang 
12594b84eec9SHong Zhang   Logically collective
12604b84eec9SHong Zhang 
12614b84eec9SHong Zhang   Input Parameter:
12624b84eec9SHong Zhang +  ts - timestepping context
12634b84eec9SHong Zhang -  mprkmtype - type of the multirate configuration
12644b84eec9SHong Zhang 
12654b84eec9SHong Zhang   Options Database:
12664b84eec9SHong Zhang .   -ts_mprk_multirate_type - <nonsplit,split>
12674b84eec9SHong Zhang 
12684b84eec9SHong Zhang   Level: intermediate
12694b84eec9SHong Zhang @*/
12704b84eec9SHong Zhang PetscErrorCode TSMPRKSetMultirateType(TS ts, TSMPRKMultirateType mprkmtype)
12714b84eec9SHong Zhang {
12724b84eec9SHong Zhang   TS_MPRK        *mprk = (TS_MPRK*)ts->data;
12734b84eec9SHong Zhang   PetscErrorCode ierr;
12744b84eec9SHong Zhang 
12754b84eec9SHong Zhang   PetscFunctionBegin;
12764b84eec9SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
12774b84eec9SHong Zhang   switch(mprkmtype){
12784b84eec9SHong Zhang     case TSMPRKNONSPLIT:
12794b84eec9SHong Zhang       ts->ops->step         = TSStep_MPRK;
12804b84eec9SHong Zhang       ts->ops->evaluatestep = TSEvaluateStep_MPRK;
12814b84eec9SHong Zhang       break;
12824b84eec9SHong Zhang     case TSMPRKSPLIT:
12834b84eec9SHong Zhang       ts->ops->step         = TSStep_MPRKSPLIT;
12844b84eec9SHong Zhang       ts->ops->evaluatestep = TSEvaluateStep_MPRKSPLIT;
12854b84eec9SHong Zhang       ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetSplits_C",TSMPRKSetSplits);CHKERRQ(ierr);
12864b84eec9SHong Zhang       break;
12874b84eec9SHong Zhang     default :
12884b84eec9SHong Zhang       SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown type '%s'",mprkmtype);
12894b84eec9SHong Zhang   }
12904b84eec9SHong Zhang   mprk->mprkmtype = mprkmtype;
12914b84eec9SHong Zhang   PetscFunctionReturn(0);
12924b84eec9SHong Zhang }
12934b84eec9SHong Zhang 
12944b84eec9SHong Zhang static PetscErrorCode TSMPRKGetType_MPRK(TS ts,TSMPRKType *mprktype)
12954b84eec9SHong Zhang {
12964b84eec9SHong Zhang   TS_MPRK *mprk = (TS_MPRK*)ts->data;
12974b84eec9SHong Zhang 
12984b84eec9SHong Zhang   PetscFunctionBegin;
12994b84eec9SHong Zhang   *mprktype = mprk->tableau->name;
13004b84eec9SHong Zhang   PetscFunctionReturn(0);
13014b84eec9SHong Zhang }
13024b84eec9SHong Zhang 
13034b84eec9SHong Zhang static PetscErrorCode TSMPRKSetType_MPRK(TS ts,TSMPRKType mprktype)
13044b84eec9SHong Zhang {
13054b84eec9SHong Zhang   TS_MPRK         *mprk = (TS_MPRK*)ts->data;
13064b84eec9SHong Zhang   PetscBool       match;
13074b84eec9SHong Zhang   MPRKTableauLink link;
13084b84eec9SHong Zhang   PetscErrorCode  ierr;
13094b84eec9SHong Zhang 
13104b84eec9SHong Zhang   PetscFunctionBegin;
13114b84eec9SHong Zhang   if (mprk->tableau) {
13124b84eec9SHong Zhang     ierr = PetscStrcmp(mprk->tableau->name,mprktype,&match);CHKERRQ(ierr);
13134b84eec9SHong Zhang     if (match) PetscFunctionReturn(0);
13144b84eec9SHong Zhang   }
13154b84eec9SHong Zhang   for (link = MPRKTableauList; link; link=link->next) {
13164b84eec9SHong Zhang     ierr = PetscStrcmp(link->tab.name,mprktype,&match);CHKERRQ(ierr);
13174b84eec9SHong Zhang     if (match) {
13184b84eec9SHong Zhang       if (ts->setupcalled) {ierr = TSMPRKTableauReset(ts);CHKERRQ(ierr);}
13194b84eec9SHong Zhang       mprk->tableau = &link->tab;
13204b84eec9SHong Zhang       if (ts->setupcalled) {ierr = TSMPRKTableauSetUp(ts);CHKERRQ(ierr);}
13214b84eec9SHong Zhang       PetscFunctionReturn(0);
13224b84eec9SHong Zhang     }
13234b84eec9SHong Zhang   }
13244b84eec9SHong Zhang   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",mprktype);
13254b84eec9SHong Zhang   PetscFunctionReturn(0);
13264b84eec9SHong Zhang }
13274b84eec9SHong Zhang 
13284b84eec9SHong Zhang static PetscErrorCode TSGetStages_MPRK(TS ts,PetscInt *ns,Vec **Y)
13294b84eec9SHong Zhang {
13304b84eec9SHong Zhang   TS_MPRK *mprk = (TS_MPRK*)ts->data;
13314b84eec9SHong Zhang 
13324b84eec9SHong Zhang   PetscFunctionBegin;
13334b84eec9SHong Zhang   *ns = mprk->tableau->s;
13344b84eec9SHong Zhang   if (Y) *Y = mprk->Y;
13354b84eec9SHong Zhang   PetscFunctionReturn(0);
13364b84eec9SHong Zhang }
13374b84eec9SHong Zhang 
13384b84eec9SHong Zhang static PetscErrorCode TSDestroy_MPRK(TS ts)
13394b84eec9SHong Zhang {
13404b84eec9SHong Zhang   PetscErrorCode ierr;
13414b84eec9SHong Zhang 
13424b84eec9SHong Zhang   PetscFunctionBegin;
13434b84eec9SHong Zhang   ierr = TSReset_MPRK(ts);CHKERRQ(ierr);
13444b84eec9SHong Zhang   if (ts->dm) {
13454b84eec9SHong Zhang     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSMPRK,DMRestrictHook_TSMPRK,ts);CHKERRQ(ierr);
13464b84eec9SHong Zhang     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSMPRK,DMSubDomainRestrictHook_TSMPRK,ts);CHKERRQ(ierr);
13474b84eec9SHong Zhang   }
13484b84eec9SHong Zhang   ierr = PetscFree(ts->data);CHKERRQ(ierr);
13494b84eec9SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKGetType_C",NULL);CHKERRQ(ierr);
13504b84eec9SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetType_C",NULL);CHKERRQ(ierr);
13514b84eec9SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetMultirateType_C",NULL);CHKERRQ(ierr);
13524b84eec9SHong Zhang   PetscFunctionReturn(0);
13534b84eec9SHong Zhang }
13544b84eec9SHong Zhang 
13554b84eec9SHong Zhang /*MC
13564b84eec9SHong Zhang       TSMPRK - ODE solver using Partitioned Runge-Kutta schemes
13574b84eec9SHong Zhang 
13584b84eec9SHong Zhang   The user should provide the right hand side of the equation
13594b84eec9SHong Zhang   using TSSetRHSFunction().
13604b84eec9SHong Zhang 
13614b84eec9SHong Zhang   Notes:
13624b84eec9SHong Zhang   The default is TSMPRKPM2, it can be changed with TSRKSetType() or -ts_mprk_type
13634b84eec9SHong Zhang 
13644b84eec9SHong Zhang   Level: beginner
13654b84eec9SHong Zhang 
13664b84eec9SHong Zhang .seealso:  TSCreate(), TS, TSSetType(), TSMPRKSetType(), TSMPRKGetType(), TSMPRKType, TSMPRKRegister(), TSMPRKSetMultirateType()
13674b84eec9SHong Zhang            TSMPRKM2, TSMPRKM3, TSMPRKRFSMR3, TSMPRKRFSMR2
13684b84eec9SHong Zhang 
13694b84eec9SHong Zhang M*/
13704b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSCreate_MPRK(TS ts)
13714b84eec9SHong Zhang {
13724b84eec9SHong Zhang   TS_MPRK        *mprk;
13734b84eec9SHong Zhang   PetscErrorCode ierr;
13744b84eec9SHong Zhang 
13754b84eec9SHong Zhang   PetscFunctionBegin;
13764b84eec9SHong Zhang   ierr = TSMPRKInitializePackage();CHKERRQ(ierr);
13774b84eec9SHong Zhang 
13784b84eec9SHong Zhang   ts->ops->reset          = TSReset_MPRK;
13794b84eec9SHong Zhang   ts->ops->destroy        = TSDestroy_MPRK;
13804b84eec9SHong Zhang   ts->ops->view           = TSView_MPRK;
13814b84eec9SHong Zhang   ts->ops->load           = TSLoad_MPRK;
13824b84eec9SHong Zhang   ts->ops->setup          = TSSetUp_MPRK;
13834b84eec9SHong Zhang   ts->ops->step           = TSStep_MPRK;
13844b84eec9SHong Zhang   ts->ops->evaluatestep   = TSEvaluateStep_MPRK;
13854b84eec9SHong Zhang   ts->ops->setfromoptions = TSSetFromOptions_MPRK;
13864b84eec9SHong Zhang   ts->ops->getstages      = TSGetStages_MPRK;
13874b84eec9SHong Zhang 
13884b84eec9SHong Zhang   ierr = PetscNewLog(ts,&mprk);CHKERRQ(ierr);
13894b84eec9SHong Zhang   ts->data = (void*)mprk;
13904b84eec9SHong Zhang 
13914b84eec9SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKGetType_C",TSMPRKGetType_MPRK);CHKERRQ(ierr);
13924b84eec9SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetType_C",TSMPRKSetType_MPRK);CHKERRQ(ierr);
13934b84eec9SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetMultirateType_C",TSMPRKSetMultirateType);CHKERRQ(ierr);
13944b84eec9SHong Zhang 
13954b84eec9SHong Zhang   ierr = TSMPRKSetType(ts,TSMPRKDefault);CHKERRQ(ierr);
13964b84eec9SHong Zhang   PetscFunctionReturn(0);
13974b84eec9SHong Zhang }
1398