xref: /petsc/src/ts/impls/multirate/mprk.c (revision 4b84eec9ecbe260473731d36f19df2cab502356b)
1*4b84eec9SHong Zhang /*
2*4b84eec9SHong Zhang   Code for time stepping with the Partitioned Runge-Kutta method
3*4b84eec9SHong Zhang 
4*4b84eec9SHong Zhang   Notes:
5*4b84eec9SHong Zhang   1) The general system is written as
6*4b84eec9SHong Zhang      Udot = F(t,U) for nonsplit RHS multi-rate RK,
7*4b84eec9SHong Zhang      user should give the indexes for both slow and fast components;
8*4b84eec9SHong Zhang   2) The general system is written as
9*4b84eec9SHong Zhang      Usdot = Fs(t,Us,Uf)
10*4b84eec9SHong Zhang      Ufdot = Ff(t,Us,Uf) for split  RHS multi-rate RK,
11*4b84eec9SHong Zhang      user should partioned RHS by themselves and also provide the indexes for both slow and fast components.
12*4b84eec9SHong Zhang */
13*4b84eec9SHong Zhang 
14*4b84eec9SHong Zhang #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
15*4b84eec9SHong Zhang #include <petscdm.h>
16*4b84eec9SHong Zhang 
17*4b84eec9SHong Zhang static TSMPRKType TSMPRKDefault = TSMPRKPM2;
18*4b84eec9SHong Zhang static PetscBool TSMPRKRegisterAllCalled;
19*4b84eec9SHong Zhang static PetscBool TSMPRKPackageInitialized;
20*4b84eec9SHong Zhang 
21*4b84eec9SHong Zhang typedef struct _MPRKTableau *MPRKTableau;
22*4b84eec9SHong Zhang struct _MPRKTableau {
23*4b84eec9SHong Zhang   char       *name;
24*4b84eec9SHong Zhang   PetscInt   order;                          /* Classical approximation order of the method i  */
25*4b84eec9SHong Zhang   PetscInt   s;                              /* Number of stages                               */
26*4b84eec9SHong Zhang   PetscReal  *Af,*bf,*cf;                    /* Tableau for fast components                    */
27*4b84eec9SHong Zhang   PetscReal  *As,*bs,*cs;                    /* Tableau for slow components                    */
28*4b84eec9SHong Zhang };
29*4b84eec9SHong Zhang typedef struct _MPRKTableauLink *MPRKTableauLink;
30*4b84eec9SHong Zhang struct _MPRKTableauLink {
31*4b84eec9SHong Zhang   struct _MPRKTableau tab;
32*4b84eec9SHong Zhang   MPRKTableauLink     next;
33*4b84eec9SHong Zhang };
34*4b84eec9SHong Zhang static MPRKTableauLink MPRKTableauList;
35*4b84eec9SHong Zhang 
36*4b84eec9SHong Zhang typedef struct {
37*4b84eec9SHong Zhang   MPRKTableau          tableau;
38*4b84eec9SHong Zhang   TSMPRKMultirateType mprkmtype;
39*4b84eec9SHong Zhang   Vec                 *Y;                          /* States computed during the step                           */
40*4b84eec9SHong Zhang   Vec                 Ytmp;
41*4b84eec9SHong Zhang   Vec                 *YdotRHS_fast;               /* Function evaluations by fast tableau for fast components  */
42*4b84eec9SHong Zhang   Vec                 *YdotRHS_slow;               /* Function evaluations by slow tableau for slow components  */
43*4b84eec9SHong Zhang   PetscScalar         *work_fast;                  /* Scalar work_fast by fast tableau                          */
44*4b84eec9SHong Zhang   PetscScalar         *work_slow;                  /* Scalar work_slow by slow tableau                          */
45*4b84eec9SHong Zhang   PetscReal           stage_time;
46*4b84eec9SHong Zhang   TSStepStatus        status;
47*4b84eec9SHong Zhang   PetscReal           ptime;
48*4b84eec9SHong Zhang   PetscReal           time_step;
49*4b84eec9SHong Zhang   IS                  is_slow,is_fast;
50*4b84eec9SHong Zhang   TS                  subts_slow,subts_fast;
51*4b84eec9SHong Zhang } TS_MPRK;
52*4b84eec9SHong Zhang 
53*4b84eec9SHong Zhang /*MC
54*4b84eec9SHong Zhang      TSMPRKPM2 - Second Order Partitioned Runge Kutta scheme.
55*4b84eec9SHong Zhang 
56*4b84eec9SHong Zhang      This method has four stages for both slow and fast parts.
57*4b84eec9SHong Zhang 
58*4b84eec9SHong Zhang      Options database:
59*4b84eec9SHong Zhang .     -ts_mprk_type pm2
60*4b84eec9SHong Zhang 
61*4b84eec9SHong Zhang      Level: advanced
62*4b84eec9SHong Zhang 
63*4b84eec9SHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType()
64*4b84eec9SHong Zhang M*/
65*4b84eec9SHong Zhang /*MC
66*4b84eec9SHong Zhang      TSMPRKPM3 - Third Order Partitioned Runge Kutta scheme.
67*4b84eec9SHong Zhang 
68*4b84eec9SHong Zhang      This method has eight stages for both slow and fast parts.
69*4b84eec9SHong Zhang 
70*4b84eec9SHong Zhang      Options database:
71*4b84eec9SHong Zhang .     -ts_mprk_type pm3  (put here temporarily)
72*4b84eec9SHong Zhang 
73*4b84eec9SHong Zhang      Level: advanced
74*4b84eec9SHong Zhang 
75*4b84eec9SHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType()
76*4b84eec9SHong Zhang M*/
77*4b84eec9SHong Zhang /*MC
78*4b84eec9SHong Zhang      TSMPRKP2 - Second Order Partitioned Runge Kutta scheme.
79*4b84eec9SHong Zhang 
80*4b84eec9SHong Zhang      This method has five stages for both slow and fast parts.
81*4b84eec9SHong Zhang 
82*4b84eec9SHong Zhang      Options database:
83*4b84eec9SHong Zhang .     -ts_mprk_type p2
84*4b84eec9SHong Zhang 
85*4b84eec9SHong Zhang      Level: advanced
86*4b84eec9SHong Zhang 
87*4b84eec9SHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType()
88*4b84eec9SHong Zhang M*/
89*4b84eec9SHong Zhang /*MC
90*4b84eec9SHong Zhang      TSMPRKP3 - Third Order Partitioned Runge Kutta scheme.
91*4b84eec9SHong Zhang 
92*4b84eec9SHong Zhang      This method has ten stages for both slow and fast parts.
93*4b84eec9SHong Zhang 
94*4b84eec9SHong Zhang      Options database:
95*4b84eec9SHong Zhang .     -ts_mprk_type p3
96*4b84eec9SHong Zhang 
97*4b84eec9SHong Zhang      Level: advanced
98*4b84eec9SHong Zhang 
99*4b84eec9SHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType()
100*4b84eec9SHong Zhang M*/
101*4b84eec9SHong Zhang 
102*4b84eec9SHong Zhang /*@C
103*4b84eec9SHong Zhang   TSMPRKRegisterAll - Registers all of the Partirioned Runge-Kutta explicit methods in TSMPRK
104*4b84eec9SHong Zhang 
105*4b84eec9SHong Zhang   Not Collective, but should be called by all processes which will need the schemes to be registered
106*4b84eec9SHong Zhang 
107*4b84eec9SHong Zhang   Level: advanced
108*4b84eec9SHong Zhang 
109*4b84eec9SHong Zhang .keywords: TS, TSMPRK, register, all
110*4b84eec9SHong Zhang 
111*4b84eec9SHong Zhang .seealso:  TSMPRKRegisterDestroy()
112*4b84eec9SHong Zhang @*/
113*4b84eec9SHong Zhang PetscErrorCode TSMPRKRegisterAll(void)
114*4b84eec9SHong Zhang {
115*4b84eec9SHong Zhang   PetscErrorCode ierr;
116*4b84eec9SHong Zhang 
117*4b84eec9SHong Zhang   PetscFunctionBegin;
118*4b84eec9SHong Zhang   if (TSMPRKRegisterAllCalled) PetscFunctionReturn(0);
119*4b84eec9SHong Zhang   TSMPRKRegisterAllCalled = PETSC_TRUE;
120*4b84eec9SHong Zhang 
121*4b84eec9SHong Zhang #define RC PetscRealConstant
122*4b84eec9SHong Zhang   {
123*4b84eec9SHong Zhang     const PetscReal
124*4b84eec9SHong Zhang       As[4][4] = {{0,0,0,0},
125*4b84eec9SHong Zhang                   {RC(1.0),0,0,0},
126*4b84eec9SHong Zhang                   {0,0,0,0},
127*4b84eec9SHong Zhang                   {0,0,RC(1.0),0}},
128*4b84eec9SHong Zhang       A[4][4]  = {{0,0,0,0},
129*4b84eec9SHong Zhang                   {RC(0.5),0,0,0},
130*4b84eec9SHong Zhang                   {RC(0.25),RC(0.25),0,0},
131*4b84eec9SHong Zhang                   {RC(0.25),RC(0.25),RC(0.5),0}},
132*4b84eec9SHong Zhang       bs[4]    = {RC(0.25),RC(0.25),RC(0.25),RC(0.25)},
133*4b84eec9SHong Zhang       b[4]     = {RC(0.25),RC(0.25),RC(0.25),RC(0.25)};
134*4b84eec9SHong Zhang     ierr = TSMPRKRegister(TSMPRKPM2,2,4,&As[0][0],bs,NULL,&A[0][0],b,NULL);CHKERRQ(ierr);
135*4b84eec9SHong Zhang   }
136*4b84eec9SHong Zhang 
137*4b84eec9SHong Zhang   /*{
138*4b84eec9SHong Zhang       const PetscReal
139*4b84eec9SHong Zhang         As[8][8] = {{0,0,0,0,0,0,0,0},
140*4b84eec9SHong Zhang                     {RC(1.0)/RC(2.0),0,0,0,0,0,0,0},
141*4b84eec9SHong Zhang                     {RC(-1.0)/RC(6.0),RC(2.0)/RC(3.0),0,0,0,0,0,0},
142*4b84eec9SHong Zhang                     {RC(1.0)/RC(3.0),RC(-1.0)/RC(3.0),RC(1.0),0,0,0,0,0},
143*4b84eec9SHong Zhang                     {0,0,0,0,0,0,0,0},
144*4b84eec9SHong Zhang                     {0,0,0,0,RC(1.0)/RC(2.0),0,0,0},
145*4b84eec9SHong Zhang                     {0,0,0,0,RC(-1.0)/RC(6.0),RC(2.0)/RC(3.0),0,0},
146*4b84eec9SHong Zhang                     {0,0,0,0,RC(1.0)/RC(3.0),RC(-1.0)/RC(3.0),RC(1.0),0}},
147*4b84eec9SHong Zhang          A[8][8] = {{0,0,0,0,0,0,0,0},
148*4b84eec9SHong Zhang                     {RC(1.0)/RC(4.0),0,0,0,0,0,0,0},
149*4b84eec9SHong Zhang                     {RC(-1.0)/RC(12.0),RC(1.0)/RC(3.0),0,0,0,0,0,0},
150*4b84eec9SHong Zhang                     {RC(1.0)/RC(6.0),RC(-1.0)/RC(6.0),RC(1.0)/RC(2.0),0,0,0,0,0},
151*4b84eec9SHong 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},
152*4b84eec9SHong 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},
153*4b84eec9SHong 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},
154*4b84eec9SHong 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}},
155*4b84eec9SHong 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)},
156*4b84eec9SHong 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)};
157*4b84eec9SHong Zhang            ierr = TSMPRKRegister(TSMPRKPM3,3,8,&As[0][0],bs,NULL,&A[0][0],b,NULL);CHKERRQ(ierr);
158*4b84eec9SHong Zhang   }*/
159*4b84eec9SHong Zhang 
160*4b84eec9SHong Zhang   {
161*4b84eec9SHong Zhang     const PetscReal
162*4b84eec9SHong Zhang       As[5][5] = {{0,0,0,0,0},
163*4b84eec9SHong Zhang                   {RC(1.0)/RC(2.0),0,0,0,0},
164*4b84eec9SHong Zhang                   {RC(1.0)/RC(2.0),0,0,0,0},
165*4b84eec9SHong Zhang                   {RC(1.0),0,0,0,0},
166*4b84eec9SHong Zhang                   {RC(1.0),0,0,0,0}},
167*4b84eec9SHong Zhang       A[5][5]  = {{0,0,0,0,0},
168*4b84eec9SHong Zhang                   {RC(1.0)/RC(2.0),0,0,0,0},
169*4b84eec9SHong Zhang                   {RC(1.0)/RC(4.0),RC(1.0)/RC(4.0),0,0,0},
170*4b84eec9SHong Zhang                   {RC(1.0)/RC(4.0),RC(1.0)/RC(4.0),RC(1.0)/RC(2.0),0,0},
171*4b84eec9SHong 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}},
172*4b84eec9SHong Zhang       bs[5]     = {RC(1.0)/RC(2.0),0,0,0,RC(1.0)/RC(2.0)},
173*4b84eec9SHong Zhang       b[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};
174*4b84eec9SHong Zhang     ierr = TSMPRKRegister(TSMPRKP2,2,5,&As[0][0],bs,NULL,&A[0][0],b,NULL);CHKERRQ(ierr);
175*4b84eec9SHong Zhang   }
176*4b84eec9SHong Zhang 
177*4b84eec9SHong Zhang   {
178*4b84eec9SHong Zhang     const PetscReal
179*4b84eec9SHong Zhang       As[10][10] = {{0,0,0,0,0,0,0,0,0,0},
180*4b84eec9SHong Zhang                     {RC(1.0)/RC(4.0),0,0,0,0,0,0,0,0,0},
181*4b84eec9SHong Zhang                     {RC(1.0)/RC(4.0),0,0,0,0,0,0,0,0,0},
182*4b84eec9SHong Zhang                     {RC(1.0)/RC(2.0),0,0,0,0,0,0,0,0,0},
183*4b84eec9SHong Zhang                     {RC(1.0)/RC(2.0),0,0,0,0,0,0,0,0,0},
184*4b84eec9SHong Zhang                     {RC(-1.0)/RC(6.0),0,0,0,RC(2.0)/RC(3.0),0,0,0,0,0},
185*4b84eec9SHong 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},
186*4b84eec9SHong 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},
187*4b84eec9SHong Zhang                     {RC(1.0)/RC(3.0),0,0,0,RC(-1.0)/RC(3.0),RC(1.0),0,0,0,0},
188*4b84eec9SHong Zhang                     {RC(1.0)/RC(3.0),0,0,0,RC(-1.0)/RC(3.0),RC(1.0),0,0,0,0}},
189*4b84eec9SHong Zhang       A[10][10]  = {{0,0,0,0,0,0,0,0,0,0},
190*4b84eec9SHong Zhang                     {RC(1.0)/RC(4.0),0,0,0,0,0,0,0,0,0},
191*4b84eec9SHong Zhang                     {RC(-1.0)/RC(12.0),RC(1.0)/RC(3.0),0,0,0,0,0,0,0,0},
192*4b84eec9SHong 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},
193*4b84eec9SHong 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},
194*4b84eec9SHong 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},
195*4b84eec9SHong 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},
196*4b84eec9SHong 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},
197*4b84eec9SHong 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},
198*4b84eec9SHong 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}},
199*4b84eec9SHong Zhang       bs[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)},
200*4b84eec9SHong Zhang       b[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};
201*4b84eec9SHong Zhang     ierr = TSMPRKRegister(TSMPRKP3,3,10,&As[0][0],bs,NULL,&A[0][0],b,NULL);CHKERRQ(ierr);
202*4b84eec9SHong Zhang   }
203*4b84eec9SHong Zhang #undef RC
204*4b84eec9SHong Zhang   PetscFunctionReturn(0);
205*4b84eec9SHong Zhang }
206*4b84eec9SHong Zhang 
207*4b84eec9SHong Zhang /*@C
208*4b84eec9SHong Zhang    TSMPRKRegisterDestroy - Frees the list of schemes that were registered by TSMPRKRegister().
209*4b84eec9SHong Zhang 
210*4b84eec9SHong Zhang    Not Collective
211*4b84eec9SHong Zhang 
212*4b84eec9SHong Zhang    Level: advanced
213*4b84eec9SHong Zhang 
214*4b84eec9SHong Zhang .keywords: TSMPRK, register, destroy
215*4b84eec9SHong Zhang .seealso: TSMPRKRegister(), TSMPRKRegisterAll()
216*4b84eec9SHong Zhang @*/
217*4b84eec9SHong Zhang PetscErrorCode TSMPRKRegisterDestroy(void)
218*4b84eec9SHong Zhang {
219*4b84eec9SHong Zhang   PetscErrorCode ierr;
220*4b84eec9SHong Zhang   MPRKTableauLink link;
221*4b84eec9SHong Zhang 
222*4b84eec9SHong Zhang   PetscFunctionBegin;
223*4b84eec9SHong Zhang   while ((link = MPRKTableauList)) {
224*4b84eec9SHong Zhang     MPRKTableau t = &link->tab;
225*4b84eec9SHong Zhang     MPRKTableauList = link->next;
226*4b84eec9SHong Zhang     ierr = PetscFree3(t->Af,t->bf,t->cf);CHKERRQ(ierr);
227*4b84eec9SHong Zhang     ierr = PetscFree3(t->As,t->bs,t->cs);CHKERRQ(ierr);
228*4b84eec9SHong Zhang     ierr = PetscFree (t->name);CHKERRQ(ierr);
229*4b84eec9SHong Zhang     ierr = PetscFree (link);CHKERRQ(ierr);
230*4b84eec9SHong Zhang   }
231*4b84eec9SHong Zhang   TSMPRKRegisterAllCalled = PETSC_FALSE;
232*4b84eec9SHong Zhang   PetscFunctionReturn(0);
233*4b84eec9SHong Zhang }
234*4b84eec9SHong Zhang 
235*4b84eec9SHong Zhang /*@C
236*4b84eec9SHong Zhang   TSMPRKInitializePackage - This function initializes everything in the TSMPRK package. It is called
237*4b84eec9SHong Zhang   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_MPRK()
238*4b84eec9SHong Zhang   when using static libraries.
239*4b84eec9SHong Zhang 
240*4b84eec9SHong Zhang   Level: developer
241*4b84eec9SHong Zhang 
242*4b84eec9SHong Zhang .keywords: TS, TSMPRK, initialize, package
243*4b84eec9SHong Zhang .seealso: PetscInitialize()
244*4b84eec9SHong Zhang @*/
245*4b84eec9SHong Zhang PetscErrorCode TSMPRKInitializePackage(void)
246*4b84eec9SHong Zhang {
247*4b84eec9SHong Zhang   PetscErrorCode ierr;
248*4b84eec9SHong Zhang 
249*4b84eec9SHong Zhang   PetscFunctionBegin;
250*4b84eec9SHong Zhang   if (TSMPRKPackageInitialized) PetscFunctionReturn(0);
251*4b84eec9SHong Zhang   TSMPRKPackageInitialized = PETSC_TRUE;
252*4b84eec9SHong Zhang   ierr = TSMPRKRegisterAll();CHKERRQ(ierr);
253*4b84eec9SHong Zhang   ierr = PetscRegisterFinalize(TSMPRKFinalizePackage);CHKERRQ(ierr);
254*4b84eec9SHong Zhang   PetscFunctionReturn(0);
255*4b84eec9SHong Zhang }
256*4b84eec9SHong Zhang 
257*4b84eec9SHong Zhang /*@C
258*4b84eec9SHong Zhang   TSRKFinalizePackage - This function destroys everything in the TSMPRK package. It is
259*4b84eec9SHong Zhang   called from PetscFinalize().
260*4b84eec9SHong Zhang 
261*4b84eec9SHong Zhang   Level: developer
262*4b84eec9SHong Zhang 
263*4b84eec9SHong Zhang .keywords: Petsc, destroy, package
264*4b84eec9SHong Zhang .seealso: PetscFinalize()
265*4b84eec9SHong Zhang @*/
266*4b84eec9SHong Zhang PetscErrorCode TSMPRKFinalizePackage(void)
267*4b84eec9SHong Zhang {
268*4b84eec9SHong Zhang   PetscErrorCode ierr;
269*4b84eec9SHong Zhang 
270*4b84eec9SHong Zhang   PetscFunctionBegin;
271*4b84eec9SHong Zhang   TSMPRKPackageInitialized = PETSC_FALSE;
272*4b84eec9SHong Zhang   ierr = TSMPRKRegisterDestroy();CHKERRQ(ierr);
273*4b84eec9SHong Zhang   PetscFunctionReturn(0);
274*4b84eec9SHong Zhang }
275*4b84eec9SHong Zhang 
276*4b84eec9SHong Zhang /*@C
277*4b84eec9SHong Zhang    TSMPRKRegister - register a MPRK scheme by providing the entries in the Butcher tableau
278*4b84eec9SHong Zhang 
279*4b84eec9SHong Zhang    Not Collective, but the same schemes should be registered on all processes on which they will be used
280*4b84eec9SHong Zhang 
281*4b84eec9SHong Zhang    Input Parameters:
282*4b84eec9SHong Zhang +  name - identifier for method
283*4b84eec9SHong Zhang .  order - approximation order of method
284*4b84eec9SHong Zhang .  s  - number of stages, this is the dimension of the matrices below
285*4b84eec9SHong Zhang .  Af - stage coefficients for fast components(dimension s*s, row-major)
286*4b84eec9SHong Zhang .  bf - step completion table for fast components(dimension s)
287*4b84eec9SHong Zhang .  cf - abscissa for fast components(dimension s)
288*4b84eec9SHong Zhang .  As - stage coefficients for slow components(dimension s*s, row-major)
289*4b84eec9SHong Zhang .  bs - step completion table for slow components(dimension s)
290*4b84eec9SHong Zhang -  cs - abscissa for slow components(dimension s)
291*4b84eec9SHong Zhang 
292*4b84eec9SHong Zhang    Notes:
293*4b84eec9SHong Zhang    Several MPRK methods are provided, this function is only needed to create new methods.
294*4b84eec9SHong Zhang 
295*4b84eec9SHong Zhang    Level: advanced
296*4b84eec9SHong Zhang 
297*4b84eec9SHong Zhang .keywords: TS, register
298*4b84eec9SHong Zhang 
299*4b84eec9SHong Zhang .seealso: TSMPRK
300*4b84eec9SHong Zhang @*/
301*4b84eec9SHong Zhang PetscErrorCode TSMPRKRegister(TSMPRKType name,PetscInt order,PetscInt s,
302*4b84eec9SHong Zhang                               const PetscReal As[],const PetscReal bs[],const PetscReal cs[],
303*4b84eec9SHong Zhang                               const PetscReal Af[],const PetscReal bf[],const PetscReal cf[])
304*4b84eec9SHong Zhang {
305*4b84eec9SHong Zhang   MPRKTableauLink link;
306*4b84eec9SHong Zhang   MPRKTableau     t;
307*4b84eec9SHong Zhang   PetscInt        i,j;
308*4b84eec9SHong Zhang   PetscErrorCode  ierr;
309*4b84eec9SHong Zhang 
310*4b84eec9SHong Zhang   PetscFunctionBegin;
311*4b84eec9SHong Zhang   PetscValidCharPointer(name,1);
312*4b84eec9SHong Zhang   PetscValidRealPointer(Af,7);
313*4b84eec9SHong Zhang   if (bf) PetscValidRealPointer(bf,8);
314*4b84eec9SHong Zhang   if (cf) PetscValidRealPointer(cf,9);
315*4b84eec9SHong Zhang   PetscValidRealPointer(As,4);
316*4b84eec9SHong Zhang   if (bs) PetscValidRealPointer(bs,5);
317*4b84eec9SHong Zhang   if (cs) PetscValidRealPointer(cs,6);
318*4b84eec9SHong Zhang 
319*4b84eec9SHong Zhang   ierr = PetscNew(&link);CHKERRQ(ierr);
320*4b84eec9SHong Zhang   t = &link->tab;
321*4b84eec9SHong Zhang 
322*4b84eec9SHong Zhang   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
323*4b84eec9SHong Zhang   t->order = order;
324*4b84eec9SHong Zhang   t->s = s;
325*4b84eec9SHong Zhang   ierr = PetscMalloc3(s*s,&t->Af,s,&t->bf,s,&t->cf);CHKERRQ(ierr);
326*4b84eec9SHong Zhang   ierr = PetscMemcpy(t->Af,Af,s*s*sizeof(Af[0]));CHKERRQ(ierr);
327*4b84eec9SHong Zhang   if (bf) {
328*4b84eec9SHong Zhang     ierr = PetscMemcpy(t->bf,bf,s*sizeof(bf[0]));CHKERRQ(ierr);
329*4b84eec9SHong Zhang   }
330*4b84eec9SHong Zhang   else
331*4b84eec9SHong Zhang     for (i=0; i<s; i++) t->bf[i] = Af[(s-1)*s+i];
332*4b84eec9SHong Zhang   if (cf) {
333*4b84eec9SHong Zhang     ierr = PetscMemcpy(t->cf,cf,s*sizeof(cf[0]));CHKERRQ(ierr);
334*4b84eec9SHong Zhang   }
335*4b84eec9SHong Zhang   else {
336*4b84eec9SHong Zhang     for (i=0; i<s; i++)
337*4b84eec9SHong Zhang       for (j=0,t->cf[i]=0; j<s; j++)
338*4b84eec9SHong Zhang         t->cf[i] += Af[i*s+j];
339*4b84eec9SHong Zhang   }
340*4b84eec9SHong Zhang   ierr = PetscMalloc3(s*s,&t->As,s,&t->bs,s,&t->cs);CHKERRQ(ierr);
341*4b84eec9SHong Zhang   ierr = PetscMemcpy(t->As,As,s*s*sizeof(As[0]));CHKERRQ(ierr);
342*4b84eec9SHong Zhang   if (bs) {
343*4b84eec9SHong Zhang     ierr = PetscMemcpy(t->bs,bs,s*sizeof(bs[0]));CHKERRQ(ierr);
344*4b84eec9SHong Zhang   }
345*4b84eec9SHong Zhang   else
346*4b84eec9SHong Zhang     for (i=0; i<s; i++) t->bs[i] = As[(s-1)*s+i];
347*4b84eec9SHong Zhang   if (cs) {
348*4b84eec9SHong Zhang     ierr = PetscMemcpy(t->cs,cs,s*sizeof(cs[0]));CHKERRQ(ierr);
349*4b84eec9SHong Zhang   }
350*4b84eec9SHong Zhang   else {
351*4b84eec9SHong Zhang     for (i=0; i<s; i++)
352*4b84eec9SHong Zhang       for (j=0,t->cs[i]=0; j<s; j++)
353*4b84eec9SHong Zhang         t->cs[i] += As[i*s+j];
354*4b84eec9SHong Zhang   }
355*4b84eec9SHong Zhang   link->next = MPRKTableauList;
356*4b84eec9SHong Zhang   MPRKTableauList = link;
357*4b84eec9SHong Zhang   PetscFunctionReturn(0);
358*4b84eec9SHong Zhang }
359*4b84eec9SHong Zhang 
360*4b84eec9SHong Zhang static PetscErrorCode TSMPRKSetSplits(TS ts)
361*4b84eec9SHong Zhang {
362*4b84eec9SHong Zhang   TS_MPRK        *mprk = (TS_MPRK*)ts->data;
363*4b84eec9SHong Zhang   DM             dm,subdm,newdm;
364*4b84eec9SHong Zhang   PetscErrorCode ierr;
365*4b84eec9SHong Zhang 
366*4b84eec9SHong Zhang   PetscFunctionBegin;
367*4b84eec9SHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"slow",&mprk->subts_slow);CHKERRQ(ierr);
368*4b84eec9SHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"fast",&mprk->subts_fast);CHKERRQ(ierr);
369*4b84eec9SHong 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");
370*4b84eec9SHong Zhang 
371*4b84eec9SHong Zhang   /* Only copy */
372*4b84eec9SHong Zhang   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
373*4b84eec9SHong Zhang   ierr = DMClone(dm,&newdm);CHKERRQ(ierr);
374*4b84eec9SHong Zhang   ierr = TSGetDM(mprk->subts_fast,&subdm);CHKERRQ(ierr);
375*4b84eec9SHong Zhang   ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr);
376*4b84eec9SHong Zhang   ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr);
377*4b84eec9SHong Zhang   ierr = TSSetDM(mprk->subts_fast,newdm);CHKERRQ(ierr);
378*4b84eec9SHong Zhang   ierr = DMDestroy(&newdm);CHKERRQ(ierr);
379*4b84eec9SHong Zhang   ierr = DMClone(dm,&newdm);CHKERRQ(ierr);
380*4b84eec9SHong Zhang   ierr = TSGetDM(mprk->subts_slow,&subdm);CHKERRQ(ierr);
381*4b84eec9SHong Zhang   ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr);
382*4b84eec9SHong Zhang   ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr);
383*4b84eec9SHong Zhang   ierr = TSSetDM(mprk->subts_slow,newdm);CHKERRQ(ierr);
384*4b84eec9SHong Zhang   ierr = DMDestroy(&newdm);CHKERRQ(ierr);
385*4b84eec9SHong Zhang   PetscFunctionReturn(0);
386*4b84eec9SHong Zhang }
387*4b84eec9SHong Zhang 
388*4b84eec9SHong Zhang /*
389*4b84eec9SHong Zhang  This if for nonsplit RHS MPRK
390*4b84eec9SHong Zhang  The step completion formula is
391*4b84eec9SHong Zhang 
392*4b84eec9SHong Zhang  x1 = x0 + h b^T YdotRHS
393*4b84eec9SHong Zhang 
394*4b84eec9SHong Zhang */
395*4b84eec9SHong Zhang static PetscErrorCode TSEvaluateStep_MPRK(TS ts,PetscInt order,Vec X,PetscBool *done)
396*4b84eec9SHong Zhang {
397*4b84eec9SHong Zhang   TS_MPRK        *mprk = (TS_MPRK*)ts->data;
398*4b84eec9SHong Zhang   MPRKTableau    tab = mprk->tableau;
399*4b84eec9SHong Zhang   Vec            Xtmp = mprk->Ytmp,Xslow,Xfast;
400*4b84eec9SHong Zhang   PetscScalar    *wf = mprk->work_fast,*ws = mprk->work_slow;
401*4b84eec9SHong Zhang   PetscReal      h = ts->time_step;
402*4b84eec9SHong Zhang   PetscInt       s = tab->s,j;
403*4b84eec9SHong Zhang   PetscErrorCode ierr;
404*4b84eec9SHong Zhang 
405*4b84eec9SHong Zhang   PetscFunctionBegin;
406*4b84eec9SHong Zhang   ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
407*4b84eec9SHong Zhang   for (j=0; j<s; j++) wf[j] = h*tab->bf[j];
408*4b84eec9SHong Zhang   for (j=0; j<s; j++) ws[j] = h*tab->bs[j];
409*4b84eec9SHong Zhang   ierr = VecCopy(X,Xtmp);CHKERRQ(ierr);
410*4b84eec9SHong Zhang   ierr = VecMAXPY(Xtmp,s,ws,mprk->YdotRHS_slow);CHKERRQ(ierr);
411*4b84eec9SHong Zhang   ierr = VecGetSubVector(Xtmp,mprk->is_slow,&Xslow);CHKERRQ(ierr);
412*4b84eec9SHong Zhang   ierr = VecISCopy(X,mprk->is_slow,SCATTER_FORWARD,Xslow);CHKERRQ(ierr);
413*4b84eec9SHong Zhang   ierr = VecRestoreSubVector(Xtmp,mprk->is_slow,&Xslow);CHKERRQ(ierr);
414*4b84eec9SHong Zhang 
415*4b84eec9SHong Zhang   /* Update fast part of X, note that the slow part has been changed but is simply discarded here */
416*4b84eec9SHong Zhang   ierr = VecCopy(X,Xtmp);CHKERRQ(ierr);
417*4b84eec9SHong Zhang   ierr = VecMAXPY(Xtmp,s,wf,mprk->YdotRHS_fast);CHKERRQ(ierr);
418*4b84eec9SHong Zhang   ierr = VecGetSubVector(Xtmp,mprk->is_fast,&Xfast);CHKERRQ(ierr);
419*4b84eec9SHong Zhang   ierr = VecISCopy(X,mprk->is_fast,SCATTER_FORWARD,Xfast);CHKERRQ(ierr);
420*4b84eec9SHong Zhang   ierr = VecRestoreSubVector(Xtmp,mprk->is_fast,&Xfast);CHKERRQ(ierr);
421*4b84eec9SHong Zhang   PetscFunctionReturn(0);
422*4b84eec9SHong Zhang }
423*4b84eec9SHong Zhang 
424*4b84eec9SHong Zhang static PetscErrorCode TSStep_MPRK(TS ts)
425*4b84eec9SHong Zhang {
426*4b84eec9SHong Zhang   TS_MPRK         *mprk = (TS_MPRK*)ts->data;
427*4b84eec9SHong Zhang   Vec             *Y = mprk->Y,Ytmp = mprk->Ytmp,*YdotRHS_fast = mprk->YdotRHS_fast,*YdotRHS_slow = mprk->YdotRHS_slow;
428*4b84eec9SHong Zhang   Vec             Yfast,Yslow;
429*4b84eec9SHong Zhang   MPRKTableau     tab = mprk->tableau;
430*4b84eec9SHong Zhang   const PetscInt  s   = tab->s;
431*4b84eec9SHong Zhang   const PetscReal *Af = tab->Af,*cf = tab->cf,*As = tab->As,*cs = tab->cs;
432*4b84eec9SHong Zhang   PetscScalar     *wf = mprk->work_fast, *ws = mprk->work_slow;
433*4b84eec9SHong Zhang   PetscInt        i,j;
434*4b84eec9SHong Zhang   PetscReal       next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step;
435*4b84eec9SHong Zhang   PetscErrorCode  ierr;
436*4b84eec9SHong Zhang 
437*4b84eec9SHong Zhang   PetscFunctionBegin;
438*4b84eec9SHong Zhang   for (i=0; i<s; i++) {
439*4b84eec9SHong Zhang     mprk->stage_time = t + h*cf[i];
440*4b84eec9SHong Zhang     ierr = TSPreStage(ts,mprk->stage_time);CHKERRQ(ierr);
441*4b84eec9SHong Zhang 
442*4b84eec9SHong Zhang     /* update the satge value for all components by slow and fast tableau respectively */
443*4b84eec9SHong Zhang     for (j=0; j<i; j++) ws[j] = h*As[i*s+j];
444*4b84eec9SHong Zhang     ierr = VecCopy(ts->vec_sol,Ytmp);CHKERRQ(ierr);
445*4b84eec9SHong Zhang     ierr = VecMAXPY(Ytmp,i,ws,YdotRHS_slow);CHKERRQ(ierr);
446*4b84eec9SHong Zhang     ierr = VecGetSubVector(Ytmp,mprk->is_slow,&Yslow);CHKERRQ(ierr);
447*4b84eec9SHong Zhang     ierr = VecISCopy(Y[i],mprk->is_slow,SCATTER_FORWARD,Yslow);CHKERRQ(ierr);
448*4b84eec9SHong Zhang     ierr = VecRestoreSubVector(Ytmp,mprk->is_slow,&Yslow);CHKERRQ(ierr);
449*4b84eec9SHong Zhang 
450*4b84eec9SHong Zhang     for (j=0; j<i; j++) wf[j] = h*Af[i*s+j];
451*4b84eec9SHong Zhang     ierr = VecCopy(ts->vec_sol,Ytmp);CHKERRQ(ierr);
452*4b84eec9SHong Zhang     ierr = VecMAXPY(Ytmp,i,wf,YdotRHS_fast);CHKERRQ(ierr);
453*4b84eec9SHong Zhang     ierr = VecGetSubVector(Ytmp,mprk->is_fast,&Yfast);CHKERRQ(ierr);
454*4b84eec9SHong Zhang     ierr = VecISCopy(Y[i],mprk->is_fast,SCATTER_FORWARD,Yfast);CHKERRQ(ierr);
455*4b84eec9SHong Zhang     ierr = VecRestoreSubVector(Ytmp,mprk->is_fast,&Yfast);CHKERRQ(ierr);
456*4b84eec9SHong Zhang 
457*4b84eec9SHong Zhang     ierr = TSPostStage(ts,mprk->stage_time,i,Y); CHKERRQ(ierr);
458*4b84eec9SHong Zhang     /* compute the stage RHS by fast and slow tableau respectively */
459*4b84eec9SHong Zhang     ierr = TSComputeRHSFunction(ts,t+h*cs[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr);
460*4b84eec9SHong Zhang     ierr = TSComputeRHSFunction(ts,t+h*cf[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr);
461*4b84eec9SHong Zhang   }
462*4b84eec9SHong Zhang   ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
463*4b84eec9SHong Zhang   ts->ptime += ts->time_step;
464*4b84eec9SHong Zhang   ts->time_step = next_time_step;
465*4b84eec9SHong Zhang   PetscFunctionReturn(0);
466*4b84eec9SHong Zhang }
467*4b84eec9SHong Zhang 
468*4b84eec9SHong Zhang /*
469*4b84eec9SHong Zhang  This if for partitioned RHS MPRK
470*4b84eec9SHong Zhang  The step completion formula is
471*4b84eec9SHong Zhang 
472*4b84eec9SHong Zhang  x1 = x0 + h b^T YdotRHS
473*4b84eec9SHong Zhang 
474*4b84eec9SHong Zhang */
475*4b84eec9SHong Zhang static PetscErrorCode TSEvaluateStep_MPRKSPLIT(TS ts,PetscInt order,Vec X,PetscBool *done)
476*4b84eec9SHong Zhang {
477*4b84eec9SHong Zhang   TS_MPRK        *mprk = (TS_MPRK*)ts->data;
478*4b84eec9SHong Zhang   MPRKTableau    tab  = mprk->tableau;
479*4b84eec9SHong Zhang   Vec            Xslow,Xfast; /* subvectors for slow and fast componets in X respectively */
480*4b84eec9SHong Zhang   PetscScalar    *wf = mprk->work_fast,*ws = mprk->work_slow;
481*4b84eec9SHong Zhang   PetscReal      h = ts->time_step;
482*4b84eec9SHong Zhang   PetscInt       s = tab->s,j;
483*4b84eec9SHong Zhang   PetscErrorCode ierr;
484*4b84eec9SHong Zhang 
485*4b84eec9SHong Zhang   PetscFunctionBegin;
486*4b84eec9SHong Zhang   ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
487*4b84eec9SHong Zhang   for (j=0; j<s; j++) wf[j] = h*tab->bf[j];
488*4b84eec9SHong Zhang   for (j=0; j<s; j++) ws[j] = h*tab->bs[j];
489*4b84eec9SHong Zhang   ierr = VecGetSubVector(X,mprk->is_slow,&Xslow);CHKERRQ(ierr);
490*4b84eec9SHong Zhang   ierr = VecGetSubVector(X,mprk->is_fast,&Xfast);CHKERRQ(ierr);
491*4b84eec9SHong Zhang   ierr = VecMAXPY(Xslow,s,ws,mprk->YdotRHS_slow);CHKERRQ(ierr);
492*4b84eec9SHong Zhang   ierr = VecMAXPY(Xfast,s,wf,mprk->YdotRHS_fast);CHKERRQ(ierr);
493*4b84eec9SHong Zhang   ierr = VecRestoreSubVector(X,mprk->is_slow,&Xfast);CHKERRQ(ierr);
494*4b84eec9SHong Zhang   ierr = VecRestoreSubVector(X,mprk->is_fast,&Xslow);CHKERRQ(ierr);
495*4b84eec9SHong Zhang   PetscFunctionReturn(0);
496*4b84eec9SHong Zhang }
497*4b84eec9SHong Zhang 
498*4b84eec9SHong Zhang static PetscErrorCode TSStep_MPRKSPLIT(TS ts)
499*4b84eec9SHong Zhang {
500*4b84eec9SHong Zhang   TS_MPRK         *mprk = (TS_MPRK*)ts->data;
501*4b84eec9SHong Zhang   MPRKTableau     tab = mprk->tableau;
502*4b84eec9SHong Zhang   Vec             *Y = mprk->Y,*YdotRHS_fast = mprk->YdotRHS_fast, *YdotRHS_slow = mprk->YdotRHS_slow;
503*4b84eec9SHong Zhang   Vec             Yslow,Yfast; /* subvectors for slow and fast components in Y[i] respectively */
504*4b84eec9SHong Zhang   const PetscInt  s = tab->s;
505*4b84eec9SHong Zhang   const PetscReal *Af = tab->Af,*cf = tab->cf,*As = tab->As,*cs = tab->cs;
506*4b84eec9SHong Zhang   PetscScalar     *wf = mprk->work_fast, *ws = mprk->work_slow;
507*4b84eec9SHong Zhang   PetscInt        i,j;
508*4b84eec9SHong Zhang   PetscReal       next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step;
509*4b84eec9SHong Zhang   PetscErrorCode  ierr;
510*4b84eec9SHong Zhang 
511*4b84eec9SHong Zhang   PetscFunctionBegin;
512*4b84eec9SHong Zhang   for (i=0; i<s; i++) {
513*4b84eec9SHong Zhang     mprk->stage_time = t + h*cf[i];
514*4b84eec9SHong Zhang     ierr = TSPreStage(ts,mprk->stage_time);CHKERRQ(ierr);
515*4b84eec9SHong Zhang     /* calculate the stage value for fast and slow components respectively */
516*4b84eec9SHong Zhang     ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
517*4b84eec9SHong Zhang     for (j=0; j<i; j++) wf[j] = h*Af[i*s+j];
518*4b84eec9SHong Zhang     for (j=0; j<i; j++) ws[j] = h*As[i*s+j];
519*4b84eec9SHong Zhang     ierr = VecGetSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr);
520*4b84eec9SHong Zhang     ierr = VecGetSubVector(Y[i],mprk->is_fast,&Yfast);CHKERRQ(ierr);
521*4b84eec9SHong Zhang     ierr = VecMAXPY(Yslow,i,ws,YdotRHS_slow);CHKERRQ(ierr);
522*4b84eec9SHong Zhang     ierr = VecMAXPY(Yfast,i,wf,YdotRHS_fast);CHKERRQ(ierr);
523*4b84eec9SHong Zhang     ierr = VecRestoreSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr);
524*4b84eec9SHong Zhang     ierr = VecRestoreSubVector(Y[i],mprk->is_fast,&Yfast);CHKERRQ(ierr);
525*4b84eec9SHong Zhang     ierr = TSPostStage(ts,mprk->stage_time,i,Y); CHKERRQ(ierr);
526*4b84eec9SHong Zhang     /* calculate the stage RHS for slow and fast components respectively */
527*4b84eec9SHong Zhang     ierr = TSComputeRHSFunction(mprk->subts_slow,t+h*cs[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr);
528*4b84eec9SHong Zhang     ierr = TSComputeRHSFunction(mprk->subts_fast,t+h*cf[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr);
529*4b84eec9SHong Zhang   }
530*4b84eec9SHong Zhang   ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
531*4b84eec9SHong Zhang   ts->ptime += ts->time_step;
532*4b84eec9SHong Zhang   ts->time_step = next_time_step;
533*4b84eec9SHong Zhang   PetscFunctionReturn(0);
534*4b84eec9SHong Zhang }
535*4b84eec9SHong Zhang 
536*4b84eec9SHong Zhang static PetscErrorCode TSMPRKTableauReset(TS ts)
537*4b84eec9SHong Zhang {
538*4b84eec9SHong Zhang   TS_MPRK        *mprk = (TS_MPRK*)ts->data;
539*4b84eec9SHong Zhang   MPRKTableau    tab = mprk->tableau;
540*4b84eec9SHong Zhang   PetscErrorCode ierr;
541*4b84eec9SHong Zhang 
542*4b84eec9SHong Zhang   PetscFunctionBegin;
543*4b84eec9SHong Zhang   if (!tab) PetscFunctionReturn(0);
544*4b84eec9SHong Zhang   ierr = PetscFree(mprk->work_fast);CHKERRQ(ierr);
545*4b84eec9SHong Zhang   ierr = PetscFree(mprk->work_slow);CHKERRQ(ierr);
546*4b84eec9SHong Zhang   ierr = VecDestroyVecs(tab->s,&mprk->Y);CHKERRQ(ierr);
547*4b84eec9SHong Zhang   if (mprk->mprkmtype == TSMPRKNONSPLIT) {
548*4b84eec9SHong Zhang     ierr = VecDestroy(&mprk->Ytmp);CHKERRQ(ierr);
549*4b84eec9SHong Zhang   }
550*4b84eec9SHong Zhang   ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_fast);CHKERRQ(ierr);
551*4b84eec9SHong Zhang   ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_slow);CHKERRQ(ierr);
552*4b84eec9SHong Zhang   PetscFunctionReturn(0);
553*4b84eec9SHong Zhang }
554*4b84eec9SHong Zhang 
555*4b84eec9SHong Zhang static PetscErrorCode TSReset_MPRK(TS ts)
556*4b84eec9SHong Zhang {
557*4b84eec9SHong Zhang   PetscErrorCode ierr;
558*4b84eec9SHong Zhang 
559*4b84eec9SHong Zhang   PetscFunctionBegin;
560*4b84eec9SHong Zhang   ierr = TSMPRKTableauReset(ts);CHKERRQ(ierr);
561*4b84eec9SHong Zhang   PetscFunctionReturn(0);
562*4b84eec9SHong Zhang }
563*4b84eec9SHong Zhang 
564*4b84eec9SHong Zhang static PetscErrorCode DMCoarsenHook_TSMPRK(DM fine,DM coarse,void *ctx)
565*4b84eec9SHong Zhang {
566*4b84eec9SHong Zhang   PetscFunctionBegin;
567*4b84eec9SHong Zhang   PetscFunctionReturn(0);
568*4b84eec9SHong Zhang }
569*4b84eec9SHong Zhang 
570*4b84eec9SHong Zhang static PetscErrorCode DMRestrictHook_TSMPRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
571*4b84eec9SHong Zhang {
572*4b84eec9SHong Zhang   PetscFunctionBegin;
573*4b84eec9SHong Zhang   PetscFunctionReturn(0);
574*4b84eec9SHong Zhang }
575*4b84eec9SHong Zhang 
576*4b84eec9SHong Zhang static PetscErrorCode DMSubDomainHook_TSMPRK(DM dm,DM subdm,void *ctx)
577*4b84eec9SHong Zhang {
578*4b84eec9SHong Zhang   PetscFunctionBegin;
579*4b84eec9SHong Zhang   PetscFunctionReturn(0);
580*4b84eec9SHong Zhang }
581*4b84eec9SHong Zhang 
582*4b84eec9SHong Zhang static PetscErrorCode DMSubDomainRestrictHook_TSMPRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
583*4b84eec9SHong Zhang {
584*4b84eec9SHong Zhang   PetscFunctionBegin;
585*4b84eec9SHong Zhang   PetscFunctionReturn(0);
586*4b84eec9SHong Zhang }
587*4b84eec9SHong Zhang 
588*4b84eec9SHong Zhang static PetscErrorCode TSMPRKTableauSetUp(TS ts)
589*4b84eec9SHong Zhang {
590*4b84eec9SHong Zhang   TS_MPRK        *mprk  = (TS_MPRK*)ts->data;
591*4b84eec9SHong Zhang   MPRKTableau    tab = mprk->tableau;
592*4b84eec9SHong Zhang   Vec            YdotRHS_fast,YdotRHS_slow;
593*4b84eec9SHong Zhang   PetscErrorCode ierr;
594*4b84eec9SHong Zhang 
595*4b84eec9SHong Zhang   PetscFunctionBegin;
596*4b84eec9SHong Zhang   ierr = PetscMalloc1(tab->s,&mprk->work_fast);CHKERRQ(ierr);
597*4b84eec9SHong Zhang   ierr = PetscMalloc1(tab->s,&mprk->work_slow);CHKERRQ(ierr);
598*4b84eec9SHong Zhang   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&mprk->Y);CHKERRQ(ierr);
599*4b84eec9SHong Zhang   if (mprk->mprkmtype == TSMPRKNONSPLIT) {
600*4b84eec9SHong Zhang     ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&mprk->YdotRHS_slow);CHKERRQ(ierr);
601*4b84eec9SHong Zhang     ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&mprk->YdotRHS_fast);CHKERRQ(ierr);
602*4b84eec9SHong Zhang     ierr = VecDuplicate(ts->vec_sol,&mprk->Ytmp);CHKERRQ(ierr);
603*4b84eec9SHong Zhang   }
604*4b84eec9SHong Zhang   if (mprk->mprkmtype == TSMPRKSPLIT) {
605*4b84eec9SHong Zhang     ierr = VecGetSubVector(ts->vec_sol,mprk->is_slow,&YdotRHS_slow);CHKERRQ(ierr);
606*4b84eec9SHong Zhang     ierr = VecGetSubVector(ts->vec_sol,mprk->is_fast,&YdotRHS_fast);CHKERRQ(ierr);
607*4b84eec9SHong Zhang     ierr = VecDuplicateVecs(YdotRHS_slow,tab->s,&mprk->YdotRHS_slow);CHKERRQ(ierr);
608*4b84eec9SHong Zhang     ierr = VecDuplicateVecs(YdotRHS_fast,tab->s,&mprk->YdotRHS_fast);CHKERRQ(ierr);
609*4b84eec9SHong Zhang     ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_slow,&YdotRHS_slow);CHKERRQ(ierr);
610*4b84eec9SHong Zhang     ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_fast,&YdotRHS_fast);CHKERRQ(ierr);
611*4b84eec9SHong Zhang   }
612*4b84eec9SHong Zhang   PetscFunctionReturn(0);
613*4b84eec9SHong Zhang }
614*4b84eec9SHong Zhang 
615*4b84eec9SHong Zhang static PetscErrorCode TSSetUp_MPRK(TS ts)
616*4b84eec9SHong Zhang {
617*4b84eec9SHong Zhang   TS_MPRK         *mprk = (TS_MPRK*)ts->data;
618*4b84eec9SHong Zhang   DM             dm;
619*4b84eec9SHong Zhang   PetscErrorCode ierr;
620*4b84eec9SHong Zhang 
621*4b84eec9SHong Zhang   PetscFunctionBegin;
622*4b84eec9SHong Zhang   ierr = TSRHSSplitGetIS(ts,"slow",&mprk->is_slow);CHKERRQ(ierr);
623*4b84eec9SHong Zhang   ierr = TSRHSSplitGetIS(ts,"fast",&mprk->is_fast);CHKERRQ(ierr);
624*4b84eec9SHong Zhang   if (!mprk->is_slow || !mprk->is_fast) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set up RHSSplits with TSRHSSplitSetIS() using split names 'slow' and 'fast' respectively in order to use -ts_type mprk");
625*4b84eec9SHong Zhang   ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr);
626*4b84eec9SHong Zhang   ierr = TSMPRKTableauSetUp(ts);CHKERRQ(ierr);
627*4b84eec9SHong Zhang   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
628*4b84eec9SHong Zhang   ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSMPRK,DMRestrictHook_TSMPRK,ts);CHKERRQ(ierr);
629*4b84eec9SHong Zhang   ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSMPRK,DMSubDomainRestrictHook_TSMPRK,ts);CHKERRQ(ierr);
630*4b84eec9SHong Zhang   ierr = PetscTryMethod(ts,"TSMPRKSetSplits_C",(TS),(ts));CHKERRQ(ierr);
631*4b84eec9SHong Zhang   PetscFunctionReturn(0);
632*4b84eec9SHong Zhang }
633*4b84eec9SHong Zhang 
634*4b84eec9SHong Zhang /* construct a database to chose nonsplit RHS mutirate mprk method or split RHS MPRK method */
635*4b84eec9SHong Zhang const char *const TSMPRKMultirateTypes[] = {"NONSPLIT","SPLIT","TSMPRKMultirateType","TSMPRK",0};
636*4b84eec9SHong Zhang 
637*4b84eec9SHong Zhang static PetscErrorCode TSSetFromOptions_MPRK(PetscOptionItems *PetscOptionsObject,TS ts)
638*4b84eec9SHong Zhang {
639*4b84eec9SHong Zhang   TS_MPRK        *mprk = (TS_MPRK*)ts->data;
640*4b84eec9SHong Zhang   PetscErrorCode ierr;
641*4b84eec9SHong Zhang 
642*4b84eec9SHong Zhang   PetscFunctionBegin;
643*4b84eec9SHong Zhang   ierr = PetscOptionsHead(PetscOptionsObject,"PRK ODE solver options");CHKERRQ(ierr);
644*4b84eec9SHong Zhang   {
645*4b84eec9SHong Zhang     MPRKTableauLink link;
646*4b84eec9SHong Zhang     PetscInt        count,choice;
647*4b84eec9SHong Zhang     PetscBool       flg;
648*4b84eec9SHong Zhang     const char      **namelist;
649*4b84eec9SHong Zhang     PetscInt        mprkmtype = 0;
650*4b84eec9SHong Zhang     for (link=MPRKTableauList,count=0; link; link=link->next,count++) ;
651*4b84eec9SHong Zhang     ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr);
652*4b84eec9SHong Zhang     for (link=MPRKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
653*4b84eec9SHong Zhang     ierr = PetscOptionsEList("-ts_mprk_type","Family of MPRK method","TSMPRKSetType",(const char*const*)namelist,count,mprk->tableau->name,&choice,&flg);CHKERRQ(ierr);
654*4b84eec9SHong Zhang     if (flg) {ierr = TSMPRKSetType(ts,namelist[choice]);CHKERRQ(ierr);}
655*4b84eec9SHong Zhang     ierr = PetscFree(namelist);CHKERRQ(ierr);
656*4b84eec9SHong 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);
657*4b84eec9SHong Zhang      if (flg) {ierr = TSMPRKSetMultirateType(ts,mprkmtype);CHKERRQ(ierr);}
658*4b84eec9SHong Zhang   }
659*4b84eec9SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
660*4b84eec9SHong Zhang   PetscFunctionReturn(0);
661*4b84eec9SHong Zhang }
662*4b84eec9SHong Zhang 
663*4b84eec9SHong Zhang static PetscErrorCode TSView_MPRK(TS ts,PetscViewer viewer)
664*4b84eec9SHong Zhang {
665*4b84eec9SHong Zhang   TS_MPRK        *mprk = (TS_MPRK*)ts->data;
666*4b84eec9SHong Zhang   PetscBool      iascii;
667*4b84eec9SHong Zhang   PetscErrorCode ierr;
668*4b84eec9SHong Zhang 
669*4b84eec9SHong Zhang   PetscFunctionBegin;
670*4b84eec9SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
671*4b84eec9SHong Zhang   if (iascii) {
672*4b84eec9SHong Zhang     MPRKTableau tab  = mprk->tableau;
673*4b84eec9SHong Zhang     TSMPRKType  mprktype;
674*4b84eec9SHong Zhang     char        fbuf[512];
675*4b84eec9SHong Zhang     char        sbuf[512];
676*4b84eec9SHong Zhang     ierr = TSMPRKGetType(ts,&mprktype);CHKERRQ(ierr);
677*4b84eec9SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  MPRK type %s\n",mprktype);CHKERRQ(ierr);
678*4b84eec9SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Order: %D\n",tab->order);CHKERRQ(ierr);
679*4b84eec9SHong Zhang     ierr = PetscFormatRealArray(fbuf,sizeof(fbuf),"% 8.6f",tab->s,tab->cf);CHKERRQ(ierr);
680*4b84eec9SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa cf = %s\n",fbuf);CHKERRQ(ierr);
681*4b84eec9SHong Zhang     ierr = PetscFormatRealArray(sbuf,sizeof(sbuf),"% 8.6f",tab->s,tab->cs);CHKERRQ(ierr);
682*4b84eec9SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa cs = %s\n",sbuf);CHKERRQ(ierr);
683*4b84eec9SHong Zhang   }
684*4b84eec9SHong Zhang   PetscFunctionReturn(0);
685*4b84eec9SHong Zhang }
686*4b84eec9SHong Zhang 
687*4b84eec9SHong Zhang static PetscErrorCode TSLoad_MPRK(TS ts,PetscViewer viewer)
688*4b84eec9SHong Zhang {
689*4b84eec9SHong Zhang   PetscErrorCode ierr;
690*4b84eec9SHong Zhang   TSAdapt        adapt;
691*4b84eec9SHong Zhang 
692*4b84eec9SHong Zhang   PetscFunctionBegin;
693*4b84eec9SHong Zhang   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
694*4b84eec9SHong Zhang   ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr);
695*4b84eec9SHong Zhang   PetscFunctionReturn(0);
696*4b84eec9SHong Zhang }
697*4b84eec9SHong Zhang 
698*4b84eec9SHong Zhang /*@C
699*4b84eec9SHong Zhang   TSMPRKSetType - Set the type of MPRK scheme
700*4b84eec9SHong Zhang 
701*4b84eec9SHong Zhang   Logically collective
702*4b84eec9SHong Zhang 
703*4b84eec9SHong Zhang   Input Parameter:
704*4b84eec9SHong Zhang +  ts - timestepping context
705*4b84eec9SHong Zhang -  mprktype - type of MPRK-scheme
706*4b84eec9SHong Zhang 
707*4b84eec9SHong Zhang   Options Database:
708*4b84eec9SHong Zhang .   -ts_mprk_type - <pm2,p2,p3>
709*4b84eec9SHong Zhang 
710*4b84eec9SHong Zhang   Level: intermediate
711*4b84eec9SHong Zhang 
712*4b84eec9SHong Zhang .seealso: TSMPRKGetType(), TSMPRK, TSMPRKType
713*4b84eec9SHong Zhang @*/
714*4b84eec9SHong Zhang PetscErrorCode TSMPRKSetType(TS ts,TSMPRKType mprktype)
715*4b84eec9SHong Zhang {
716*4b84eec9SHong Zhang   PetscErrorCode ierr;
717*4b84eec9SHong Zhang 
718*4b84eec9SHong Zhang   PetscFunctionBegin;
719*4b84eec9SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
720*4b84eec9SHong Zhang   PetscValidCharPointer(mprktype,2);
721*4b84eec9SHong Zhang   ierr = PetscTryMethod(ts,"TSMPRKSetType_C",(TS,TSMPRKType),(ts,mprktype));CHKERRQ(ierr);
722*4b84eec9SHong Zhang   PetscFunctionReturn(0);
723*4b84eec9SHong Zhang }
724*4b84eec9SHong Zhang 
725*4b84eec9SHong Zhang /*@C
726*4b84eec9SHong Zhang   TSMPRKGetType - Get the type of MPRK scheme
727*4b84eec9SHong Zhang 
728*4b84eec9SHong Zhang   Logically collective
729*4b84eec9SHong Zhang 
730*4b84eec9SHong Zhang   Input Parameter:
731*4b84eec9SHong Zhang .  ts - timestepping context
732*4b84eec9SHong Zhang 
733*4b84eec9SHong Zhang   Output Parameter:
734*4b84eec9SHong Zhang .  mprktype - type of MPRK-scheme
735*4b84eec9SHong Zhang 
736*4b84eec9SHong Zhang   Level: intermediate
737*4b84eec9SHong Zhang 
738*4b84eec9SHong Zhang .seealso: TSMPRKGetType()
739*4b84eec9SHong Zhang @*/
740*4b84eec9SHong Zhang PetscErrorCode TSMPRKGetType(TS ts,TSMPRKType *mprktype)
741*4b84eec9SHong Zhang {
742*4b84eec9SHong Zhang   PetscErrorCode ierr;
743*4b84eec9SHong Zhang 
744*4b84eec9SHong Zhang   PetscFunctionBegin;
745*4b84eec9SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
746*4b84eec9SHong Zhang   ierr = PetscUseMethod(ts,"TSMPRKGetType_C",(TS,TSMPRKType*),(ts,mprktype));CHKERRQ(ierr);
747*4b84eec9SHong Zhang   PetscFunctionReturn(0);
748*4b84eec9SHong Zhang }
749*4b84eec9SHong Zhang 
750*4b84eec9SHong Zhang /*@C
751*4b84eec9SHong Zhang   TSMPRKSetMultirateType - Set the type of MPRK multirate scheme
752*4b84eec9SHong Zhang 
753*4b84eec9SHong Zhang   Logically collective
754*4b84eec9SHong Zhang 
755*4b84eec9SHong Zhang   Input Parameter:
756*4b84eec9SHong Zhang +  ts - timestepping context
757*4b84eec9SHong Zhang -  mprkmtype - type of the multirate configuration
758*4b84eec9SHong Zhang 
759*4b84eec9SHong Zhang   Options Database:
760*4b84eec9SHong Zhang .   -ts_mprk_multirate_type - <nonsplit,split>
761*4b84eec9SHong Zhang 
762*4b84eec9SHong Zhang   Level: intermediate
763*4b84eec9SHong Zhang @*/
764*4b84eec9SHong Zhang PetscErrorCode TSMPRKSetMultirateType(TS ts, TSMPRKMultirateType mprkmtype)
765*4b84eec9SHong Zhang {
766*4b84eec9SHong Zhang   TS_MPRK        *mprk = (TS_MPRK*)ts->data;
767*4b84eec9SHong Zhang   PetscErrorCode ierr;
768*4b84eec9SHong Zhang 
769*4b84eec9SHong Zhang   PetscFunctionBegin;
770*4b84eec9SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
771*4b84eec9SHong Zhang   switch(mprkmtype){
772*4b84eec9SHong Zhang     case TSMPRKNONSPLIT:
773*4b84eec9SHong Zhang       ts->ops->step         = TSStep_MPRK;
774*4b84eec9SHong Zhang       ts->ops->evaluatestep = TSEvaluateStep_MPRK;
775*4b84eec9SHong Zhang       break;
776*4b84eec9SHong Zhang     case TSMPRKSPLIT:
777*4b84eec9SHong Zhang       ts->ops->step         = TSStep_MPRKSPLIT;
778*4b84eec9SHong Zhang       ts->ops->evaluatestep = TSEvaluateStep_MPRKSPLIT;
779*4b84eec9SHong Zhang       ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetSplits_C",TSMPRKSetSplits);CHKERRQ(ierr);
780*4b84eec9SHong Zhang       break;
781*4b84eec9SHong Zhang     default :
782*4b84eec9SHong Zhang       SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown type '%s'",mprkmtype);
783*4b84eec9SHong Zhang   }
784*4b84eec9SHong Zhang   mprk->mprkmtype = mprkmtype;
785*4b84eec9SHong Zhang   PetscFunctionReturn(0);
786*4b84eec9SHong Zhang }
787*4b84eec9SHong Zhang 
788*4b84eec9SHong Zhang static PetscErrorCode TSMPRKGetType_MPRK(TS ts,TSMPRKType *mprktype)
789*4b84eec9SHong Zhang {
790*4b84eec9SHong Zhang   TS_MPRK *mprk = (TS_MPRK*)ts->data;
791*4b84eec9SHong Zhang 
792*4b84eec9SHong Zhang   PetscFunctionBegin;
793*4b84eec9SHong Zhang   *mprktype = mprk->tableau->name;
794*4b84eec9SHong Zhang   PetscFunctionReturn(0);
795*4b84eec9SHong Zhang }
796*4b84eec9SHong Zhang 
797*4b84eec9SHong Zhang static PetscErrorCode TSMPRKSetType_MPRK(TS ts,TSMPRKType mprktype)
798*4b84eec9SHong Zhang {
799*4b84eec9SHong Zhang   TS_MPRK         *mprk = (TS_MPRK*)ts->data;
800*4b84eec9SHong Zhang   PetscBool       match;
801*4b84eec9SHong Zhang   MPRKTableauLink link;
802*4b84eec9SHong Zhang   PetscErrorCode  ierr;
803*4b84eec9SHong Zhang 
804*4b84eec9SHong Zhang   PetscFunctionBegin;
805*4b84eec9SHong Zhang   if (mprk->tableau) {
806*4b84eec9SHong Zhang     ierr = PetscStrcmp(mprk->tableau->name,mprktype,&match);CHKERRQ(ierr);
807*4b84eec9SHong Zhang     if (match) PetscFunctionReturn(0);
808*4b84eec9SHong Zhang   }
809*4b84eec9SHong Zhang   for (link = MPRKTableauList; link; link=link->next) {
810*4b84eec9SHong Zhang     ierr = PetscStrcmp(link->tab.name,mprktype,&match);CHKERRQ(ierr);
811*4b84eec9SHong Zhang     if (match) {
812*4b84eec9SHong Zhang       if (ts->setupcalled) {ierr = TSMPRKTableauReset(ts);CHKERRQ(ierr);}
813*4b84eec9SHong Zhang       mprk->tableau = &link->tab;
814*4b84eec9SHong Zhang       if (ts->setupcalled) {ierr = TSMPRKTableauSetUp(ts);CHKERRQ(ierr);}
815*4b84eec9SHong Zhang       PetscFunctionReturn(0);
816*4b84eec9SHong Zhang     }
817*4b84eec9SHong Zhang   }
818*4b84eec9SHong Zhang   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",mprktype);
819*4b84eec9SHong Zhang   PetscFunctionReturn(0);
820*4b84eec9SHong Zhang }
821*4b84eec9SHong Zhang 
822*4b84eec9SHong Zhang static PetscErrorCode TSGetStages_MPRK(TS ts,PetscInt *ns,Vec **Y)
823*4b84eec9SHong Zhang {
824*4b84eec9SHong Zhang   TS_MPRK *mprk = (TS_MPRK*)ts->data;
825*4b84eec9SHong Zhang 
826*4b84eec9SHong Zhang   PetscFunctionBegin;
827*4b84eec9SHong Zhang   *ns = mprk->tableau->s;
828*4b84eec9SHong Zhang   if (Y) *Y = mprk->Y;
829*4b84eec9SHong Zhang   PetscFunctionReturn(0);
830*4b84eec9SHong Zhang }
831*4b84eec9SHong Zhang 
832*4b84eec9SHong Zhang static PetscErrorCode TSDestroy_MPRK(TS ts)
833*4b84eec9SHong Zhang {
834*4b84eec9SHong Zhang   PetscErrorCode ierr;
835*4b84eec9SHong Zhang 
836*4b84eec9SHong Zhang   PetscFunctionBegin;
837*4b84eec9SHong Zhang   ierr = TSReset_MPRK(ts);CHKERRQ(ierr);
838*4b84eec9SHong Zhang   if (ts->dm) {
839*4b84eec9SHong Zhang     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSMPRK,DMRestrictHook_TSMPRK,ts);CHKERRQ(ierr);
840*4b84eec9SHong Zhang     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSMPRK,DMSubDomainRestrictHook_TSMPRK,ts);CHKERRQ(ierr);
841*4b84eec9SHong Zhang   }
842*4b84eec9SHong Zhang   ierr = PetscFree(ts->data);CHKERRQ(ierr);
843*4b84eec9SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKGetType_C",NULL);CHKERRQ(ierr);
844*4b84eec9SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetType_C",NULL);CHKERRQ(ierr);
845*4b84eec9SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetMultirateType_C",NULL);CHKERRQ(ierr);
846*4b84eec9SHong Zhang   PetscFunctionReturn(0);
847*4b84eec9SHong Zhang }
848*4b84eec9SHong Zhang 
849*4b84eec9SHong Zhang /*MC
850*4b84eec9SHong Zhang       TSMPRK - ODE solver using Partitioned Runge-Kutta schemes
851*4b84eec9SHong Zhang 
852*4b84eec9SHong Zhang   The user should provide the right hand side of the equation
853*4b84eec9SHong Zhang   using TSSetRHSFunction().
854*4b84eec9SHong Zhang 
855*4b84eec9SHong Zhang   Notes:
856*4b84eec9SHong Zhang   The default is TSMPRKPM2, it can be changed with TSRKSetType() or -ts_mprk_type
857*4b84eec9SHong Zhang 
858*4b84eec9SHong Zhang   Level: beginner
859*4b84eec9SHong Zhang 
860*4b84eec9SHong Zhang .seealso:  TSCreate(), TS, TSSetType(), TSMPRKSetType(), TSMPRKGetType(), TSMPRKType, TSMPRKRegister(), TSMPRKSetMultirateType()
861*4b84eec9SHong Zhang            TSMPRKM2, TSMPRKM3, TSMPRKRFSMR3, TSMPRKRFSMR2
862*4b84eec9SHong Zhang 
863*4b84eec9SHong Zhang M*/
864*4b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSCreate_MPRK(TS ts)
865*4b84eec9SHong Zhang {
866*4b84eec9SHong Zhang   TS_MPRK        *mprk;
867*4b84eec9SHong Zhang   PetscErrorCode ierr;
868*4b84eec9SHong Zhang 
869*4b84eec9SHong Zhang   PetscFunctionBegin;
870*4b84eec9SHong Zhang   ierr = TSMPRKInitializePackage();CHKERRQ(ierr);
871*4b84eec9SHong Zhang 
872*4b84eec9SHong Zhang   ts->ops->reset          = TSReset_MPRK;
873*4b84eec9SHong Zhang   ts->ops->destroy        = TSDestroy_MPRK;
874*4b84eec9SHong Zhang   ts->ops->view           = TSView_MPRK;
875*4b84eec9SHong Zhang   ts->ops->load           = TSLoad_MPRK;
876*4b84eec9SHong Zhang   ts->ops->setup          = TSSetUp_MPRK;
877*4b84eec9SHong Zhang   ts->ops->step           = TSStep_MPRK;
878*4b84eec9SHong Zhang   ts->ops->evaluatestep   = TSEvaluateStep_MPRK;
879*4b84eec9SHong Zhang   ts->ops->setfromoptions = TSSetFromOptions_MPRK;
880*4b84eec9SHong Zhang   ts->ops->getstages      = TSGetStages_MPRK;
881*4b84eec9SHong Zhang 
882*4b84eec9SHong Zhang   ierr = PetscNewLog(ts,&mprk);CHKERRQ(ierr);
883*4b84eec9SHong Zhang   ts->data = (void*)mprk;
884*4b84eec9SHong Zhang 
885*4b84eec9SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKGetType_C",TSMPRKGetType_MPRK);CHKERRQ(ierr);
886*4b84eec9SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetType_C",TSMPRKSetType_MPRK);CHKERRQ(ierr);
887*4b84eec9SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetMultirateType_C",TSMPRKSetMultirateType);CHKERRQ(ierr);
888*4b84eec9SHong Zhang 
889*4b84eec9SHong Zhang   ierr = TSMPRKSetType(ts,TSMPRKDefault);CHKERRQ(ierr);
890*4b84eec9SHong Zhang   PetscFunctionReturn(0);
891*4b84eec9SHong Zhang }
892