xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision 55de54fcc8bf4f46c3d8f3d7a760507b847cdb46)
1e4dd521cSBarry Smith /*
22c9cb062Sluzhanghpp   Code for time stepping with the Runge-Kutta method(single-rate RK and multi-rate RK)
3f68a32c8SEmil Constantinescu 
4f68a32c8SEmil Constantinescu   Notes:
52c9cb062Sluzhanghpp   1) The general system is written as
62c9cb062Sluzhanghpp      Udot = F(t,U) for single-rate RK;
72c9cb062Sluzhanghpp   2) The general system is written as
8*55de54fcSHong Zhang      Udot = F(t,U) for the nonsplit version of multi-rate RK,
92c9cb062Sluzhanghpp      user should give the indexes for both slow and fast components;
102c9cb062Sluzhanghpp   3) The general system is written as
112c9cb062Sluzhanghpp      Usdot = Fs(t,Us,Uf)
12*55de54fcSHong Zhang      Ufdot = Ff(t,Us,Uf) for multi-rate RK with RHS splits,
132c9cb062Sluzhanghpp      user should partioned RHS by themselves and also provide the indexes for both slow and fast components.
14e4dd521cSBarry Smith */
152c9cb062Sluzhanghpp 
16af0996ceSBarry Smith #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
17f68a32c8SEmil Constantinescu #include <petscdm.h>
18*55de54fcSHong Zhang #include <petscdmshell.h>
19f68a32c8SEmil Constantinescu 
20484bcdc7SDebojyoti Ghosh static TSRKType  TSRKDefault = TSRK3BS;
212c9cb062Sluzhanghpp static TSRKType  TSMRKDefault = TSRK2A;
22f68a32c8SEmil Constantinescu static PetscBool TSRKRegisterAllCalled;
23f68a32c8SEmil Constantinescu static PetscBool TSRKPackageInitialized;
24f68a32c8SEmil Constantinescu 
25f68a32c8SEmil Constantinescu typedef struct _RKTableau *RKTableau;
26f68a32c8SEmil Constantinescu struct _RKTableau {
27f68a32c8SEmil Constantinescu   char       *name;
28d760c35bSDebojyoti Ghosh   PetscInt   order;               /* Classical approximation order of the method i              */
29f68a32c8SEmil Constantinescu   PetscInt   s;                   /* Number of stages                                           */
303a8a9803SLisandro Dalcin   PetscInt   p;                   /* Interpolation order                                        */
31d760c35bSDebojyoti Ghosh   PetscBool  FSAL;                /* flag to indicate if tableau is FSAL                        */
32f68a32c8SEmil Constantinescu   PetscReal  *A,*b,*c;            /* Tableau                                                    */
33f68a32c8SEmil Constantinescu   PetscReal  *bembed;             /* Embedded formula of order one less (order-1)               */
34f68a32c8SEmil Constantinescu   PetscReal  *binterp;            /* Dense output formula                                       */
35f68a32c8SEmil Constantinescu   PetscReal  ccfl;                /* Placeholder for CFL coefficient relative to forward Euler  */
36f68a32c8SEmil Constantinescu };
37f68a32c8SEmil Constantinescu typedef struct _RKTableauLink *RKTableauLink;
38f68a32c8SEmil Constantinescu struct _RKTableauLink {
39f68a32c8SEmil Constantinescu   struct _RKTableau tab;
40f68a32c8SEmil Constantinescu   RKTableauLink     next;
41f68a32c8SEmil Constantinescu };
42f68a32c8SEmil Constantinescu static RKTableauLink RKTableauList;
43e4dd521cSBarry Smith 
44e4dd521cSBarry Smith typedef struct {
45f68a32c8SEmil Constantinescu   RKTableau    tableau;
46*55de54fcSHong Zhang   Vec          X0;
47f68a32c8SEmil Constantinescu   Vec          *Y;               /* States computed during the step                                              */
482c9cb062Sluzhanghpp   Vec          *YdotRHS;         /* Function evaluations for the non-stiff part and contains all components      */
492c9cb062Sluzhanghpp   Vec          *YdotRHS_fast;    /* Function evaluations for the non-stiff part and contains fast components     */
502c9cb062Sluzhanghpp   Vec          *YdotRHS_slow;    /* Function evaluations for the non-stiff part and contains slow components     */
51ad8e2604SHong Zhang   Vec          *VecDeltaLam;     /* Increment of the adjoint sensitivity w.r.t IC at stage                       */
52ad8e2604SHong Zhang   Vec          *VecDeltaMu;      /* Increment of the adjoint sensitivity w.r.t P at stage                        */
530f7a1166SHong Zhang   Vec          VecCostIntegral0; /* backup for roll-backs due to events                                          */
54f68a32c8SEmil Constantinescu   PetscScalar  *work;            /* Scalar work                                                                  */
552c9cb062Sluzhanghpp   PetscInt     slow;             /* flag indicates call slow components solver (0) or fast components solver (1) */
56f68a32c8SEmil Constantinescu   PetscReal    stage_time;
57f68a32c8SEmil Constantinescu   TSStepStatus status;
580f7a1166SHong Zhang   PetscReal    ptime;
590f7a1166SHong Zhang   PetscReal    time_step;
602c9cb062Sluzhanghpp   PetscInt     dtratio;          /* ratio between slow time step size and fast step size                         */
61*55de54fcSHong Zhang   IS           is_fast,is_slow;
62*55de54fcSHong Zhang   TS           subts_fast,subts_slow;
635f70b478SJed Brown } TS_RK;
64e4dd521cSBarry Smith 
652c9cb062Sluzhanghpp 
66f68a32c8SEmil Constantinescu /*MC
67287c2655SBarry Smith      TSRK1FE - First order forward Euler scheme.
68262ff7bbSBarry Smith 
69f68a32c8SEmil Constantinescu      This method has one stage.
70f68a32c8SEmil Constantinescu 
71287c2655SBarry Smith      Options database:
729bd3e852SBarry Smith .     -ts_rk_type 1fe
73287c2655SBarry Smith 
74f68a32c8SEmil Constantinescu      Level: advanced
75f68a32c8SEmil Constantinescu 
76287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
77f68a32c8SEmil Constantinescu M*/
78f68a32c8SEmil Constantinescu /*MC
792109b73fSDebojyoti Ghosh      TSRK2A - Second order RK scheme.
80f68a32c8SEmil Constantinescu 
81f68a32c8SEmil Constantinescu      This method has two stages.
82f68a32c8SEmil Constantinescu 
83287c2655SBarry Smith      Options database:
849bd3e852SBarry Smith .     -ts_rk_type 2a
85287c2655SBarry Smith 
86f68a32c8SEmil Constantinescu      Level: advanced
87f68a32c8SEmil Constantinescu 
88287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
89f68a32c8SEmil Constantinescu M*/
90f68a32c8SEmil Constantinescu /*MC
91f68a32c8SEmil Constantinescu      TSRK3 - Third order RK scheme.
92f68a32c8SEmil Constantinescu 
93f68a32c8SEmil Constantinescu      This method has three stages.
94f68a32c8SEmil Constantinescu 
95287c2655SBarry Smith      Options database:
969bd3e852SBarry Smith .     -ts_rk_type 3
97287c2655SBarry Smith 
98f68a32c8SEmil Constantinescu      Level: advanced
99f68a32c8SEmil Constantinescu 
100287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
101f68a32c8SEmil Constantinescu M*/
102f68a32c8SEmil Constantinescu /*MC
1032109b73fSDebojyoti Ghosh      TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method.
1042109b73fSDebojyoti Ghosh 
105d1905564SLisandro Dalcin      This method has four stages with the First Same As Last (FSAL) property.
1062109b73fSDebojyoti Ghosh 
107287c2655SBarry Smith      Options database:
1089bd3e852SBarry Smith .     -ts_rk_type 3bs
109287c2655SBarry Smith 
1102109b73fSDebojyoti Ghosh      Level: advanced
1112109b73fSDebojyoti Ghosh 
11298164e13SLisandro Dalcin      References: https://doi.org/10.1016/0893-9659(89)90079-7
113d1905564SLisandro Dalcin 
114287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
1152109b73fSDebojyoti Ghosh M*/
1162109b73fSDebojyoti Ghosh /*MC
117f68a32c8SEmil Constantinescu      TSRK4 - Fourth order RK scheme.
118f68a32c8SEmil Constantinescu 
1192e698f8bSDebojyoti Ghosh      This is the classical Runge-Kutta method with four stages.
120f68a32c8SEmil Constantinescu 
121287c2655SBarry Smith      Options database:
1229bd3e852SBarry Smith .     -ts_rk_type 4
123287c2655SBarry Smith 
124f68a32c8SEmil Constantinescu      Level: advanced
125f68a32c8SEmil Constantinescu 
126287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
127f68a32c8SEmil Constantinescu M*/
128f68a32c8SEmil Constantinescu /*MC
1292e698f8bSDebojyoti Ghosh      TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method.
130f68a32c8SEmil Constantinescu 
131f68a32c8SEmil Constantinescu      This method has six stages.
132f68a32c8SEmil Constantinescu 
133287c2655SBarry Smith      Options database:
1349bd3e852SBarry Smith .     -ts_rk_type 5f
135287c2655SBarry Smith 
136f68a32c8SEmil Constantinescu      Level: advanced
137f68a32c8SEmil Constantinescu 
138287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
139f68a32c8SEmil Constantinescu M*/
1402109b73fSDebojyoti Ghosh /*MC
1412e698f8bSDebojyoti Ghosh      TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method.
1422109b73fSDebojyoti Ghosh 
143d1905564SLisandro Dalcin      This method has seven stages with the First Same As Last (FSAL) property.
1442109b73fSDebojyoti Ghosh 
145287c2655SBarry Smith      Options database:
1469bd3e852SBarry Smith .     -ts_rk_type 5dp
147287c2655SBarry Smith 
1482109b73fSDebojyoti Ghosh      Level: advanced
1492109b73fSDebojyoti Ghosh 
15098164e13SLisandro Dalcin      References: https://doi.org/10.1016/0771-050X(80)90013-3
151d1905564SLisandro Dalcin 
152287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
1532109b73fSDebojyoti Ghosh M*/
15405e23783SLisandro Dalcin /*MC
15505e23783SLisandro Dalcin      TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method.
15605e23783SLisandro Dalcin 
15705e23783SLisandro Dalcin      This method has eight stages with the First Same As Last (FSAL) property.
15805e23783SLisandro Dalcin 
159287c2655SBarry Smith      Options database:
1609bd3e852SBarry Smith .     -ts_rk_type 5bs
161287c2655SBarry Smith 
16205e23783SLisandro Dalcin      Level: advanced
16305e23783SLisandro Dalcin 
16498164e13SLisandro Dalcin      References: https://doi.org/10.1016/0898-1221(96)00141-1
16505e23783SLisandro Dalcin 
166287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
16705e23783SLisandro Dalcin M*/
168262ff7bbSBarry Smith 
169f68a32c8SEmil Constantinescu /*@C
170f68a32c8SEmil Constantinescu   TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK
171262ff7bbSBarry Smith 
172f68a32c8SEmil Constantinescu   Not Collective, but should be called by all processes which will need the schemes to be registered
173262ff7bbSBarry Smith 
174f68a32c8SEmil Constantinescu   Level: advanced
175262ff7bbSBarry Smith 
176f68a32c8SEmil Constantinescu .keywords: TS, TSRK, register, all
177262ff7bbSBarry Smith 
178f68a32c8SEmil Constantinescu .seealso:  TSRKRegisterDestroy()
179262ff7bbSBarry Smith @*/
180f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void)
181262ff7bbSBarry Smith {
1824ac538c5SBarry Smith   PetscErrorCode ierr;
183262ff7bbSBarry Smith 
184262ff7bbSBarry Smith   PetscFunctionBegin;
185f68a32c8SEmil Constantinescu   if (TSRKRegisterAllCalled) PetscFunctionReturn(0);
186f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_TRUE;
187e4dd521cSBarry Smith 
1884e82626cSLisandro Dalcin #define RC PetscRealConstant
189e4dd521cSBarry Smith   {
190f68a32c8SEmil Constantinescu     const PetscReal
1914e82626cSLisandro Dalcin       A[1][1] = {{0}},
1924e82626cSLisandro Dalcin       b[1]    = {RC(1.0)};
1933a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
194e4dd521cSBarry Smith   }
195db046809SBarry Smith   {
196f68a32c8SEmil Constantinescu     const PetscReal
1974e82626cSLisandro Dalcin       A[2][2]   = {{0,0},
1984e82626cSLisandro Dalcin                    {RC(1.0),0}},
1994e82626cSLisandro Dalcin       b[2]      =  {RC(0.5),RC(0.5)},
2004e82626cSLisandro Dalcin       bembed[2] =  {RC(1.0),0};
2013a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
202db046809SBarry Smith   }
203f68a32c8SEmil Constantinescu   {
204f68a32c8SEmil Constantinescu     const PetscReal
20517f6384fSBarry Smith       A[3][3] = {{0,0,0},
2064e82626cSLisandro Dalcin                  {RC(2.0)/RC(3.0),0,0},
2074e82626cSLisandro Dalcin                  {RC(-1.0)/RC(3.0),RC(1.0),0}},
2084e82626cSLisandro Dalcin       b[3]    =  {RC(0.25),RC(0.5),RC(0.25)};
2093a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
210fdefd5e5SDebojyoti Ghosh   }
211fdefd5e5SDebojyoti Ghosh   {
212fdefd5e5SDebojyoti Ghosh     const PetscReal
21317f6384fSBarry Smith       A[4][4]   = {{0,0,0,0},
2144e82626cSLisandro Dalcin                    {RC(1.0)/RC(2.0),0,0,0},
2154e82626cSLisandro Dalcin                    {0,RC(3.0)/RC(4.0),0,0},
2164e82626cSLisandro Dalcin                    {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}},
2174e82626cSLisandro Dalcin       b[4]      =  {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0},
2184e82626cSLisandro Dalcin       bembed[4] =  {RC(7.0)/RC(24.0),RC(1.0)/RC(4.0),RC(1.0)/RC(3.0),RC(1.0)/RC(8.0)};
2193a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
220db046809SBarry Smith   }
221f68a32c8SEmil Constantinescu   {
222f68a32c8SEmil Constantinescu     const PetscReal
223f68a32c8SEmil Constantinescu       A[4][4] = {{0,0,0,0},
2244e82626cSLisandro Dalcin                  {RC(0.5),0,0,0},
2254e82626cSLisandro Dalcin                  {0,RC(0.5),0,0},
2264e82626cSLisandro Dalcin                  {0,0,RC(1.0),0}},
2274e82626cSLisandro Dalcin       b[4]    =  {RC(1.0)/RC(6.0),RC(1.0)/RC(3.0),RC(1.0)/RC(3.0),RC(1.0)/RC(6.0)};
2283a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
229db046809SBarry Smith   }
230f68a32c8SEmil Constantinescu   {
231f68a32c8SEmil Constantinescu     const PetscReal
23217f6384fSBarry Smith       A[6][6]   = {{0,0,0,0,0,0},
2334e82626cSLisandro Dalcin                    {RC(0.25),0,0,0,0,0},
2344e82626cSLisandro Dalcin                    {RC(3.0)/RC(32.0),RC(9.0)/RC(32.0),0,0,0,0},
2354e82626cSLisandro Dalcin                    {RC(1932.0)/RC(2197.0),RC(-7200.0)/RC(2197.0),RC(7296.0)/RC(2197.0),0,0,0},
2364e82626cSLisandro Dalcin                    {RC(439.0)/RC(216.0),RC(-8.0),RC(3680.0)/RC(513.0),RC(-845.0)/RC(4104.0),0,0},
2374e82626cSLisandro Dalcin                    {RC(-8.0)/RC(27.0),RC(2.0),RC(-3544.0)/RC(2565.0),RC(1859.0)/RC(4104.0),RC(-11.0)/RC(40.0),0}},
2384e82626cSLisandro Dalcin       b[6]      =  {RC(16.0)/RC(135.0),0,RC(6656.0)/RC(12825.0),RC(28561.0)/RC(56430.0),RC(-9.0)/RC(50.0),RC(2.0)/RC(55.0)},
2394e82626cSLisandro Dalcin       bembed[6] =  {RC(25.0)/RC(216.0),0,RC(1408.0)/RC(2565.0),RC(2197.0)/RC(4104.0),RC(-1.0)/RC(5.0),0};
2403a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
241fdefd5e5SDebojyoti Ghosh   }
242fdefd5e5SDebojyoti Ghosh   {
243fdefd5e5SDebojyoti Ghosh     const PetscReal
24417f6384fSBarry Smith       A[7][7]   = {{0,0,0,0,0,0,0},
2454e82626cSLisandro Dalcin                    {RC(1.0)/RC(5.0),0,0,0,0,0,0},
2464e82626cSLisandro Dalcin                    {RC(3.0)/RC(40.0),RC(9.0)/RC(40.0),0,0,0,0,0},
2474e82626cSLisandro Dalcin                    {RC(44.0)/RC(45.0),RC(-56.0)/RC(15.0),RC(32.0)/RC(9.0),0,0,0,0},
2484e82626cSLisandro Dalcin                    {RC(19372.0)/RC(6561.0),RC(-25360.0)/RC(2187.0),RC(64448.0)/RC(6561.0),RC(-212.0)/RC(729.0),0,0,0},
2494e82626cSLisandro Dalcin                    {RC(9017.0)/RC(3168.0),RC(-355.0)/RC(33.0),RC(46732.0)/RC(5247.0),RC(49.0)/RC(176.0),RC(-5103.0)/RC(18656.0),0,0},
2504e82626cSLisandro Dalcin                    {RC(35.0)/RC(384.0),0,RC(500.0)/RC(1113.0),RC(125.0)/RC(192.0),RC(-2187.0)/RC(6784.0),RC(11.0)/RC(84.0),0}},
2514e82626cSLisandro Dalcin       b[7]      =  {RC(35.0)/RC(384.0),0,RC(500.0)/RC(1113.0),RC(125.0)/RC(192.0),RC(-2187.0)/RC(6784.0),RC(11.0)/RC(84.0),0},
252a685a763Sluzhanghpp         bembed[7] =  {RC(5179.0)/RC(57600.0),0,RC(7571.0)/RC(16695.0),RC(393.0)/RC(640.0),RC(-92097.0)/RC(339200.0),RC(187.0)/RC(2100.0),RC(1.0)/RC(40.0)},
253a685a763Sluzhanghpp         binterp[7][5] =  {{RC(1.0),RC(-4034104133.0)/RC(1410260304.0),RC(105330401.0)/RC(33982176.0),RC(-13107642775.0)/RC(11282082432.0),RC(6542295.0)/RC(470086768.0)},
254a685a763Sluzhanghpp                     {0,0,0,0,0},
255a685a763Sluzhanghpp                     {0,RC(132343189600.0)/RC(32700410799.0),RC(-833316000.0)/RC(131326951.0),RC(91412856700.0)/RC(32700410799.0),RC(-523383600.0)/RC(10900136933.0)},
256a685a763Sluzhanghpp                     {0,RC(-115792950.0)/RC(29380423.0),RC(185270875.0)/RC(16991088.0),RC(-12653452475.0)/RC(1880347072.0),RC(98134425.0)/RC(235043384.0)},
257a685a763Sluzhanghpp                     {0,RC(70805911779.0)/RC(24914598704.0),RC(-4531260609.0)/RC(600351776.0),RC(988140236175.0)/RC(199316789632.0),RC(-14307999165.0)/RC(24914598704.0)},
258a685a763Sluzhanghpp                     {0,RC(-331320693.0)/RC(205662961.0),RC(31361737.0)/RC(7433601.0),RC(-2426908385.0)/RC(822651844.0),RC(97305120.0)/RC(205662961.0)},
259a685a763Sluzhanghpp                     {0,RC(44764047.0)/RC(29380423.0),RC(-1532549.0)/RC(353981.0),RC(90730570.0)/RC(29380423.0),RC(-8293050.0)/RC(29380423.0)}};
260a685a763Sluzhanghpp 
261a685a763Sluzhanghpp         ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,5,binterp[0]);CHKERRQ(ierr);
262f68a32c8SEmil Constantinescu   }
26305e23783SLisandro Dalcin   {
26405e23783SLisandro Dalcin     const PetscReal
26517f6384fSBarry Smith       A[8][8]   = {{0,0,0,0,0,0,0,0},
2664e82626cSLisandro Dalcin                    {RC(1.0)/RC(6.0),0,0,0,0,0,0,0},
2674e82626cSLisandro Dalcin                    {RC(2.0)/RC(27.0),RC(4.0)/RC(27.0),0,0,0,0,0,0},
2684e82626cSLisandro Dalcin                    {RC(183.0)/RC(1372.0),RC(-162.0)/RC(343.0),RC(1053.0)/RC(1372.0),0,0,0,0,0},
2694e82626cSLisandro Dalcin                    {RC(68.0)/RC(297.0),RC(-4.0)/RC(11.0),RC(42.0)/RC(143.0),RC(1960.0)/RC(3861.0),0,0,0,0},
2704e82626cSLisandro Dalcin                    {RC(597.0)/RC(22528.0),RC(81.0)/RC(352.0),RC(63099.0)/RC(585728.0),RC(58653.0)/RC(366080.0),RC(4617.0)/RC(20480.0),0,0,0},
2714e82626cSLisandro Dalcin                    {RC(174197.0)/RC(959244.0),RC(-30942.0)/RC(79937.0),RC(8152137.0)/RC(19744439.0),RC(666106.0)/RC(1039181.0),RC(-29421.0)/RC(29068.0),RC(482048.0)/RC(414219.0),0,0},
2724e82626cSLisandro Dalcin                    {RC(587.0)/RC(8064.0),0,RC(4440339.0)/RC(15491840.0),RC(24353.0)/RC(124800.0),RC(387.0)/RC(44800.0),RC(2152.0)/RC(5985.0),RC(7267.0)/RC(94080.0),0}},
2734e82626cSLisandro Dalcin       b[8]      =  {RC(587.0)/RC(8064.0),0,RC(4440339.0)/RC(15491840.0),RC(24353.0)/RC(124800.0),RC(387.0)/RC(44800.0),RC(2152.0)/RC(5985.0),RC(7267.0)/RC(94080.0),0},
2744e82626cSLisandro Dalcin       bembed[8] =  {RC(2479.0)/RC(34992.0),0,RC(123.0)/RC(416.0),RC(612941.0)/RC(3411720.0),RC(43.0)/RC(1440.0),RC(2272.0)/RC(6561.0),RC(79937.0)/RC(1113912.0),RC(3293.0)/RC(556956.0)};
27505e23783SLisandro Dalcin     ierr = TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
27605e23783SLisandro Dalcin   }
2774e82626cSLisandro Dalcin #undef RC
278db046809SBarry Smith   PetscFunctionReturn(0);
279db046809SBarry Smith }
280db046809SBarry Smith 
281f68a32c8SEmil Constantinescu /*@C
282f68a32c8SEmil Constantinescu    TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister().
283f68a32c8SEmil Constantinescu 
284f68a32c8SEmil Constantinescu    Not Collective
285f68a32c8SEmil Constantinescu 
286f68a32c8SEmil Constantinescu    Level: advanced
287f68a32c8SEmil Constantinescu 
288f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy
289f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll()
290f68a32c8SEmil Constantinescu @*/
291f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void)
292e4dd521cSBarry Smith {
293f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
294f68a32c8SEmil Constantinescu   RKTableauLink  link;
295f68a32c8SEmil Constantinescu 
296f68a32c8SEmil Constantinescu   PetscFunctionBegin;
297f68a32c8SEmil Constantinescu   while ((link = RKTableauList)) {
298f68a32c8SEmil Constantinescu     RKTableau t = &link->tab;
299f68a32c8SEmil Constantinescu     RKTableauList = link->next;
300f68a32c8SEmil Constantinescu     ierr = PetscFree3(t->A,t->b,t->c);  CHKERRQ(ierr);
301f68a32c8SEmil Constantinescu     ierr = PetscFree (t->bembed);       CHKERRQ(ierr);
302f68a32c8SEmil Constantinescu     ierr = PetscFree (t->binterp);      CHKERRQ(ierr);
303f68a32c8SEmil Constantinescu     ierr = PetscFree (t->name);         CHKERRQ(ierr);
304f68a32c8SEmil Constantinescu     ierr = PetscFree (link);            CHKERRQ(ierr);
305f68a32c8SEmil Constantinescu   }
306f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_FALSE;
307f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
308f68a32c8SEmil Constantinescu }
309f68a32c8SEmil Constantinescu 
310f68a32c8SEmil Constantinescu /*@C
311f68a32c8SEmil Constantinescu   TSRKInitializePackage - This function initializes everything in the TSRK package. It is called
3128a690491SBarry Smith   from TSInitializePackage().
313f68a32c8SEmil Constantinescu 
314f68a32c8SEmil Constantinescu   Level: developer
315f68a32c8SEmil Constantinescu 
316f68a32c8SEmil Constantinescu .keywords: TS, TSRK, initialize, package
317f68a32c8SEmil Constantinescu .seealso: PetscInitialize()
318f68a32c8SEmil Constantinescu @*/
319f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void)
320f68a32c8SEmil Constantinescu {
321186e87acSLisandro Dalcin   PetscErrorCode ierr;
322e4dd521cSBarry Smith 
323e4dd521cSBarry Smith   PetscFunctionBegin;
324f68a32c8SEmil Constantinescu   if (TSRKPackageInitialized) PetscFunctionReturn(0);
325f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_TRUE;
326f68a32c8SEmil Constantinescu   ierr = TSRKRegisterAll();CHKERRQ(ierr);
327f68a32c8SEmil Constantinescu   ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr);
328f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
329f68a32c8SEmil Constantinescu }
330186e87acSLisandro Dalcin 
331f68a32c8SEmil Constantinescu /*@C
332f68a32c8SEmil Constantinescu   TSRKFinalizePackage - This function destroys everything in the TSRK package. It is
333f68a32c8SEmil Constantinescu   called from PetscFinalize().
334186e87acSLisandro Dalcin 
335f68a32c8SEmil Constantinescu   Level: developer
336f68a32c8SEmil Constantinescu 
337f68a32c8SEmil Constantinescu .keywords: Petsc, destroy, package
338f68a32c8SEmil Constantinescu .seealso: PetscFinalize()
339f68a32c8SEmil Constantinescu @*/
340f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void)
341f68a32c8SEmil Constantinescu {
342f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
343f68a32c8SEmil Constantinescu 
344f68a32c8SEmil Constantinescu   PetscFunctionBegin;
345f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_FALSE;
346f68a32c8SEmil Constantinescu   ierr = TSRKRegisterDestroy();CHKERRQ(ierr);
347f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
348f68a32c8SEmil Constantinescu }
349f68a32c8SEmil Constantinescu 
350f68a32c8SEmil Constantinescu /*@C
351f68a32c8SEmil Constantinescu    TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
352f68a32c8SEmil Constantinescu 
353f68a32c8SEmil Constantinescu    Not Collective, but the same schemes should be registered on all processes on which they will be used
354f68a32c8SEmil Constantinescu 
355f68a32c8SEmil Constantinescu    Input Parameters:
356f68a32c8SEmil Constantinescu +  name - identifier for method
357f68a32c8SEmil Constantinescu .  order - approximation order of method
358f68a32c8SEmil Constantinescu .  s - number of stages, this is the dimension of the matrices below
359f68a32c8SEmil Constantinescu .  A - stage coefficients (dimension s*s, row-major)
360f68a32c8SEmil Constantinescu .  b - step completion table (dimension s; NULL to use last row of A)
361f68a32c8SEmil Constantinescu .  c - abscissa (dimension s; NULL to use row sums of A)
362f68a32c8SEmil Constantinescu .  bembed - completion table for embedded method (dimension s; NULL if not available)
3633a8a9803SLisandro Dalcin .  p - Order of the interpolation scheme, equal to the number of columns of binterp
3643a8a9803SLisandro Dalcin -  binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1)
365f68a32c8SEmil Constantinescu 
366f68a32c8SEmil Constantinescu    Notes:
367f68a32c8SEmil Constantinescu    Several RK methods are provided, this function is only needed to create new methods.
368f68a32c8SEmil Constantinescu 
369f68a32c8SEmil Constantinescu    Level: advanced
370f68a32c8SEmil Constantinescu 
371f68a32c8SEmil Constantinescu .keywords: TS, register
372f68a32c8SEmil Constantinescu 
373f68a32c8SEmil Constantinescu .seealso: TSRK
374f68a32c8SEmil Constantinescu @*/
375f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s,
376f68a32c8SEmil Constantinescu                             const PetscReal A[],const PetscReal b[],const PetscReal c[],
3773a8a9803SLisandro Dalcin                             const PetscReal bembed[],PetscInt p,const PetscReal binterp[])
378f68a32c8SEmil Constantinescu {
379f68a32c8SEmil Constantinescu   PetscErrorCode  ierr;
380f68a32c8SEmil Constantinescu   RKTableauLink   link;
381f68a32c8SEmil Constantinescu   RKTableau       t;
382f68a32c8SEmil Constantinescu   PetscInt        i,j;
383f68a32c8SEmil Constantinescu 
384f68a32c8SEmil Constantinescu   PetscFunctionBegin;
3853a8a9803SLisandro Dalcin   PetscValidCharPointer(name,1);
3863a8a9803SLisandro Dalcin   PetscValidRealPointer(A,4);
3873a8a9803SLisandro Dalcin   if (b) PetscValidRealPointer(b,5);
3883a8a9803SLisandro Dalcin   if (c) PetscValidRealPointer(c,6);
3893a8a9803SLisandro Dalcin   if (bembed) PetscValidRealPointer(bembed,7);
3903a8a9803SLisandro Dalcin   if (binterp || p > 1) PetscValidRealPointer(binterp,9);
3913a8a9803SLisandro Dalcin 
3921d36bdfdSBarry Smith   ierr = TSRKInitializePackage();CHKERRQ(ierr);
39395dccacaSBarry Smith   ierr = PetscNew(&link);CHKERRQ(ierr);
394f68a32c8SEmil Constantinescu   t = &link->tab;
3953a8a9803SLisandro Dalcin 
396f68a32c8SEmil Constantinescu   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
397f68a32c8SEmil Constantinescu   t->order = order;
398f68a32c8SEmil Constantinescu   t->s = s;
399dcca6d9dSJed Brown   ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr);
400f68a32c8SEmil Constantinescu   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
401f68a32c8SEmil Constantinescu   if (b)  { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); }
402f68a32c8SEmil Constantinescu   else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i];
403f68a32c8SEmil Constantinescu   if (c)  { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); }
404f68a32c8SEmil Constantinescu   else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j];
405d760c35bSDebojyoti Ghosh   t->FSAL = PETSC_TRUE;
406d760c35bSDebojyoti Ghosh   for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE;
4073a8a9803SLisandro Dalcin 
408f68a32c8SEmil Constantinescu   if (bembed) {
409785e854fSJed Brown     ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr);
410f68a32c8SEmil Constantinescu     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
411f68a32c8SEmil Constantinescu   }
412f68a32c8SEmil Constantinescu 
4133a8a9803SLisandro Dalcin   if (!binterp) { p = 1; binterp = t->b; }
4143a8a9803SLisandro Dalcin   t->p = p;
4153a8a9803SLisandro Dalcin   ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr);
4163a8a9803SLisandro Dalcin   ierr = PetscMemcpy(t->binterp,binterp,s*p*sizeof(binterp[0]));CHKERRQ(ierr);
4173a8a9803SLisandro Dalcin 
418f68a32c8SEmil Constantinescu   link->next = RKTableauList;
419f68a32c8SEmil Constantinescu   RKTableauList = link;
420f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
421f68a32c8SEmil Constantinescu }
422f68a32c8SEmil Constantinescu 
423e4dd521cSBarry Smith /*
4242c9cb062Sluzhanghpp  This is for single-step RK method
425f68a32c8SEmil Constantinescu  The step completion formula is
426e4dd521cSBarry Smith 
427f68a32c8SEmil Constantinescu  x1 = x0 + h b^T YdotRHS
428f68a32c8SEmil Constantinescu 
429f68a32c8SEmil Constantinescu  This function can be called before or after ts->vec_sol has been updated.
430f68a32c8SEmil Constantinescu  Suppose we have a completion formula (b) and an embedded formula (be) of different order.
431f68a32c8SEmil Constantinescu  We can write
432f68a32c8SEmil Constantinescu 
433f68a32c8SEmil Constantinescu  x1e = x0 + h be^T YdotRHS
434f68a32c8SEmil Constantinescu      = x1 - h b^T YdotRHS + h be^T YdotRHS
435f68a32c8SEmil Constantinescu      = x1 + h (be - b)^T YdotRHS
436f68a32c8SEmil Constantinescu 
437f68a32c8SEmil Constantinescu  so we can evaluate the method with different order even after the step has been optimistically completed.
438e4dd521cSBarry Smith */
439f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done)
440f68a32c8SEmil Constantinescu {
441f68a32c8SEmil Constantinescu   TS_RK          *rk   = (TS_RK*)ts->data;
442f68a32c8SEmil Constantinescu   RKTableau      tab  = rk->tableau;
443f68a32c8SEmil Constantinescu   PetscScalar    *w    = rk->work;
444f68a32c8SEmil Constantinescu   PetscReal      h;
445f68a32c8SEmil Constantinescu   PetscInt       s    = tab->s,j;
446f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
447f68a32c8SEmil Constantinescu 
448f68a32c8SEmil Constantinescu   PetscFunctionBegin;
449f68a32c8SEmil Constantinescu   switch (rk->status) {
450f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
451f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
452f68a32c8SEmil Constantinescu     h = ts->time_step; break;
453f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
454be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
455f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
456f68a32c8SEmil Constantinescu   }
457f68a32c8SEmil Constantinescu   if (order == tab->order) {
458f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) {
459f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
460f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->b[j];
461f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
462f68a32c8SEmil Constantinescu     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
463f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
464f68a32c8SEmil Constantinescu   } else if (order == tab->order-1) {
465f68a32c8SEmil Constantinescu     if (!tab->bembed) goto unavailable;
466f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) { /*Complete with the embedded method (be)*/
467f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
468f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
469f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
470f68a32c8SEmil Constantinescu     } else {  /*Rollback and re-complete using (be-b) */
471f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
472f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
473f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
474f68a32c8SEmil Constantinescu     }
475f68a32c8SEmil Constantinescu     if (done) *done = PETSC_TRUE;
476f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
477f68a32c8SEmil Constantinescu   }
478f68a32c8SEmil Constantinescu unavailable:
479f68a32c8SEmil Constantinescu   if (done) *done = PETSC_FALSE;
480a7fac7c2SEmil Constantinescu   else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"RK '%s' of order %D cannot evaluate step at order %D. Consider using -ts_adapt_type none or a different method that has an embedded estimate.",tab->name,tab->order,order);
481f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
482f68a32c8SEmil Constantinescu }
483f68a32c8SEmil Constantinescu 
4840f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts)
4850f7a1166SHong Zhang {
4860f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
4870f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
4880f7a1166SHong Zhang   const PetscInt  s = tab->s;
4890f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
4900f7a1166SHong Zhang   Vec             *Y = rk->Y;
4910f7a1166SHong Zhang   PetscInt        i;
4920f7a1166SHong Zhang   PetscErrorCode  ierr;
4930f7a1166SHong Zhang 
4940f7a1166SHong Zhang   PetscFunctionBegin;
4950f7a1166SHong Zhang   /* backup cost integral */
4960f7a1166SHong Zhang   ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr);
4970f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
4980f7a1166SHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
49951a7f651SHong Zhang     ierr = TSComputeCostIntegrand(ts,rk->ptime+rk->time_step*c[i],Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
5000f7a1166SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
5010f7a1166SHong Zhang   }
5020f7a1166SHong Zhang   PetscFunctionReturn(0);
5030f7a1166SHong Zhang }
5040f7a1166SHong Zhang 
5050f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts)
5060f7a1166SHong Zhang {
5070f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
5080f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
5090f7a1166SHong Zhang   const PetscInt  s = tab->s;
5100f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
5110f7a1166SHong Zhang   Vec             *Y = rk->Y;
5120f7a1166SHong Zhang   PetscInt        i;
5130f7a1166SHong Zhang   PetscErrorCode  ierr;
5140f7a1166SHong Zhang 
5150f7a1166SHong Zhang   PetscFunctionBegin;
5160f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
5170f7a1166SHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
51851a7f651SHong Zhang     ierr = TSComputeCostIntegrand(ts,ts->ptime+ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
5190f7a1166SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
5200f7a1166SHong Zhang   }
5210f7a1166SHong Zhang   PetscFunctionReturn(0);
5220f7a1166SHong Zhang }
5230f7a1166SHong Zhang 
524fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts)
525fecfb714SLisandro Dalcin {
526fecfb714SLisandro Dalcin   TS_RK           *rk = (TS_RK*)ts->data;
527fecfb714SLisandro Dalcin   RKTableau       tab = rk->tableau;
528fecfb714SLisandro Dalcin   const PetscInt  s  = tab->s;
529fecfb714SLisandro Dalcin   const PetscReal *b = tab->b;
530fecfb714SLisandro Dalcin   PetscScalar     *w = rk->work;
531fecfb714SLisandro Dalcin   Vec             *YdotRHS = rk->YdotRHS;
532fecfb714SLisandro Dalcin   PetscInt        j;
533fecfb714SLisandro Dalcin   PetscReal       h;
534fecfb714SLisandro Dalcin   PetscErrorCode  ierr;
535fecfb714SLisandro Dalcin 
536fecfb714SLisandro Dalcin   PetscFunctionBegin;
537fecfb714SLisandro Dalcin   switch (rk->status) {
538fecfb714SLisandro Dalcin   case TS_STEP_INCOMPLETE:
539fecfb714SLisandro Dalcin   case TS_STEP_PENDING:
540fecfb714SLisandro Dalcin     h = ts->time_step; break;
541fecfb714SLisandro Dalcin   case TS_STEP_COMPLETE:
542fecfb714SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
543fecfb714SLisandro Dalcin   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
544fecfb714SLisandro Dalcin   }
545fecfb714SLisandro Dalcin   for (j=0; j<s; j++) w[j] = -h*b[j];
546fecfb714SLisandro Dalcin   ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr);
547fecfb714SLisandro Dalcin   PetscFunctionReturn(0);
548fecfb714SLisandro Dalcin }
549fecfb714SLisandro Dalcin 
550f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts)
551f68a32c8SEmil Constantinescu {
552f68a32c8SEmil Constantinescu   TS_RK            *rk  = (TS_RK*)ts->data;
553f68a32c8SEmil Constantinescu   RKTableau        tab  = rk->tableau;
554f68a32c8SEmil Constantinescu   const PetscInt   s = tab->s;
555fecfb714SLisandro Dalcin   const PetscReal  *A = tab->A,*c = tab->c;
556f68a32c8SEmil Constantinescu   PetscScalar      *w = rk->work;
557f68a32c8SEmil Constantinescu   Vec              *Y = rk->Y,*YdotRHS = rk->YdotRHS;
558d1905564SLisandro Dalcin   PetscBool        FSAL = tab->FSAL;
559f68a32c8SEmil Constantinescu   TSAdapt          adapt;
560fecfb714SLisandro Dalcin   PetscInt         i,j;
561be5899b3SLisandro Dalcin   PetscInt         rejections = 0;
562be5899b3SLisandro Dalcin   PetscBool        stageok,accept = PETSC_TRUE;
563be5899b3SLisandro Dalcin   PetscReal        next_time_step = ts->time_step;
564f68a32c8SEmil Constantinescu   PetscErrorCode   ierr;
565f68a32c8SEmil Constantinescu 
566f68a32c8SEmil Constantinescu   PetscFunctionBegin;
567d1905564SLisandro Dalcin   if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE;
568d1905564SLisandro Dalcin   if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); }
569d1905564SLisandro Dalcin 
570f68a32c8SEmil Constantinescu   rk->status = TS_STEP_INCOMPLETE;
571be5899b3SLisandro Dalcin   while (!ts->reason && rk->status != TS_STEP_COMPLETE) {
572be5899b3SLisandro Dalcin     PetscReal t = ts->ptime;
573f68a32c8SEmil Constantinescu     PetscReal h = ts->time_step;
574f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
5759be3e283SDebojyoti Ghosh       rk->stage_time = t + h*c[i];
5769be3e283SDebojyoti Ghosh       ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr);
577f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
578f68a32c8SEmil Constantinescu       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
579f68a32c8SEmil Constantinescu       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
5809be3e283SDebojyoti Ghosh       ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
581f68a32c8SEmil Constantinescu       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
582be5899b3SLisandro Dalcin       ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr);
583fecfb714SLisandro Dalcin       if (!stageok) goto reject_step;
5848f738a7cSLisandro Dalcin       if (FSAL && !i) continue;
585f68a32c8SEmil Constantinescu       ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
586f68a32c8SEmil Constantinescu     }
587be5899b3SLisandro Dalcin 
588be5899b3SLisandro Dalcin     rk->status = TS_STEP_INCOMPLETE;
589f68a32c8SEmil Constantinescu     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
590f68a32c8SEmil Constantinescu     rk->status = TS_STEP_PENDING;
591f68a32c8SEmil Constantinescu     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
592f68a32c8SEmil Constantinescu     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
5931917a363SLisandro Dalcin     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr);
594fecfb714SLisandro Dalcin     ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
595be5899b3SLisandro Dalcin     rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
596be5899b3SLisandro Dalcin     if (!accept) {  /*Roll back the current step */
597fecfb714SLisandro Dalcin       ierr = TSRollBack_RK(ts);CHKERRQ(ierr);
598be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
599be5899b3SLisandro Dalcin       goto reject_step;
600be5899b3SLisandro Dalcin     }
601be5899b3SLisandro Dalcin 
6020f7a1166SHong Zhang     if (ts->costintegralfwd) {  /*Save the info for the later use in cost integral evaluation*/
6030f7a1166SHong Zhang       rk->ptime     = ts->ptime;
6040f7a1166SHong Zhang       rk->time_step = ts->time_step;
605493ed95fSHong Zhang     }
606be5899b3SLisandro Dalcin 
607f68a32c8SEmil Constantinescu     ts->ptime += ts->time_step;
608f68a32c8SEmil Constantinescu     ts->time_step = next_time_step;
609f68a32c8SEmil Constantinescu     break;
610be5899b3SLisandro Dalcin 
611be5899b3SLisandro Dalcin     reject_step:
612fecfb714SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
613be5899b3SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
614be5899b3SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
615be5899b3SLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
616f68a32c8SEmil Constantinescu     }
617f68a32c8SEmil Constantinescu   }
618f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
619e4dd521cSBarry Smith }
620e4dd521cSBarry Smith 
621f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts)
622f6a906c0SBarry Smith {
623f6a906c0SBarry Smith   TS_RK          *rk  = (TS_RK*)ts->data;
624be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
625be5899b3SLisandro Dalcin   PetscInt       s   = tab->s;
626f6a906c0SBarry Smith   PetscErrorCode ierr;
627f6a906c0SBarry Smith 
628f6a906c0SBarry Smith   PetscFunctionBegin;
629f6a906c0SBarry Smith   if (ts->adjointsetupcalled++) PetscFunctionReturn(0);
630abc2977eSBarry Smith   ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr);
631f6a906c0SBarry Smith   if(ts->vecs_sensip) {
632abc2977eSBarry Smith     ierr = VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr);
633f6a906c0SBarry Smith   }
634f6a906c0SBarry Smith   PetscFunctionReturn(0);
635f6a906c0SBarry Smith }
636f6a906c0SBarry Smith 
63742f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts)
638d2daff3dSHong Zhang {
639c235aa8dSHong Zhang   TS_RK            *rk  = (TS_RK*)ts->data;
640c235aa8dSHong Zhang   RKTableau        tab  = rk->tableau;
641c235aa8dSHong Zhang   const PetscInt   s    = tab->s;
642c235aa8dSHong Zhang   const PetscReal  *A = tab->A,*b = tab->b,*c = tab->c;
643c235aa8dSHong Zhang   PetscScalar      *w    = rk->work;
6443389c7d9SStefano Zampini   Vec              *Y    = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu;
645b18ea86cSHong Zhang   PetscInt         i,j,nadj;
6463d3eaba7SBarry Smith   PetscReal        t = ts->ptime;
647d2daff3dSHong Zhang   PetscErrorCode   ierr;
648cef76868SBarry Smith   PetscReal        h = ts->time_step;
649c235aa8dSHong Zhang 
650d2daff3dSHong Zhang   PetscFunctionBegin;
651c235aa8dSHong Zhang   rk->status = TS_STEP_INCOMPLETE;
652c235aa8dSHong Zhang   for (i=s-1; i>=0; i--) {
6533389c7d9SStefano Zampini     Mat       J;
6543389c7d9SStefano Zampini     PetscReal stage_time = t + h*(1.0-c[i]);
6553389c7d9SStefano Zampini     PetscBool zero = PETSC_FALSE;
6563389c7d9SStefano Zampini 
6573389c7d9SStefano Zampini     ierr = TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr);
6583389c7d9SStefano Zampini     ierr = TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);CHKERRQ(ierr);
659493ed95fSHong Zhang     if (ts->vec_costintegral) {
6603389c7d9SStefano Zampini       ierr = TSAdjointComputeDRDYFunction(ts,stage_time,Y[i],ts->vecs_drdy);CHKERRQ(ierr);
661493ed95fSHong Zhang     }
6623389c7d9SStefano Zampini     /* Stage values of mu */
6633389c7d9SStefano Zampini     if (ts->vecs_sensip) {
6643389c7d9SStefano Zampini       ierr = TSAdjointComputeRHSJacobian(ts,stage_time,Y[i],ts->Jacp);CHKERRQ(ierr);
6653389c7d9SStefano Zampini       if (ts->vec_costintegral) {
6663389c7d9SStefano Zampini         ierr = TSAdjointComputeDRDPFunction(ts,stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr);
6673389c7d9SStefano Zampini       }
6683389c7d9SStefano Zampini     }
6693389c7d9SStefano Zampini 
6703389c7d9SStefano Zampini     if (b[i] == 0 && i == s-1) zero = PETSC_TRUE;
671abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
6723389c7d9SStefano Zampini       DM  dm;
6733389c7d9SStefano Zampini       Vec VecSensiTemp;
6743389c7d9SStefano Zampini 
6753389c7d9SStefano Zampini       ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
6763389c7d9SStefano Zampini       ierr = DMGetGlobalVector(dm,&VecSensiTemp);CHKERRQ(ierr);
6773389c7d9SStefano Zampini       /* Stage values of lambda */
6783389c7d9SStefano Zampini       if (!zero) {
6793389c7d9SStefano Zampini         ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp);CHKERRQ(ierr);
6803389c7d9SStefano Zampini         ierr = VecScale(VecSensiTemp,-h*b[i]);CHKERRQ(ierr);
6813389c7d9SStefano Zampini         for (j=i+1; j<s; j++) w[j-i-1]  = -h*A[j*s+i];
6823389c7d9SStefano Zampini         ierr = VecMAXPY(VecSensiTemp,s-i-1,w,&VecDeltaLam[nadj*s+i+1]);CHKERRQ(ierr);
6833389c7d9SStefano Zampini         ierr = MatMultTranspose(J,VecSensiTemp,VecDeltaLam[nadj*s+i]);CHKERRQ(ierr);
6843389c7d9SStefano Zampini       } else {
6853389c7d9SStefano Zampini         ierr = VecSet(VecDeltaLam[nadj*s+i],0);CHKERRQ(ierr);
6863389c7d9SStefano Zampini       }
687493ed95fSHong Zhang       if (ts->vec_costintegral) {
688493ed95fSHong Zhang         ierr = VecAXPY(VecDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdy[nadj]);CHKERRQ(ierr);
689493ed95fSHong Zhang       }
6906497ce14SHong Zhang 
691ad8e2604SHong Zhang       /* Stage values of mu */
6926497ce14SHong Zhang       if (ts->vecs_sensip) {
6933389c7d9SStefano Zampini         if (!zero) {
6943389c7d9SStefano Zampini           ierr = MatMultTranspose(ts->Jacp,VecSensiTemp,VecDeltaMu[nadj*s+i]);CHKERRQ(ierr);
6953389c7d9SStefano Zampini         } else {
6963389c7d9SStefano Zampini           ierr = VecSet(VecDeltaMu[nadj*s+i],0);CHKERRQ(ierr);
697493ed95fSHong Zhang         }
698493ed95fSHong Zhang         if (ts->vec_costintegral) {
699493ed95fSHong Zhang           ierr = VecAXPY(VecDeltaMu[nadj*s+i],-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr);
700493ed95fSHong Zhang         }
701ad8e2604SHong Zhang       }
7023389c7d9SStefano Zampini       ierr = DMRestoreGlobalVector(dm,&VecSensiTemp);CHKERRQ(ierr);
703c235aa8dSHong Zhang     }
7046497ce14SHong Zhang   }
705c235aa8dSHong Zhang 
706c235aa8dSHong Zhang   for (j=0; j<s; j++) w[j] = 1.0;
707abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
708b18ea86cSHong Zhang     ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr);
7096497ce14SHong Zhang     if (ts->vecs_sensip) {
710ad8e2604SHong Zhang       ierr = VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);CHKERRQ(ierr);
711b18ea86cSHong Zhang     }
7126497ce14SHong Zhang   }
713c235aa8dSHong Zhang   rk->status = TS_STEP_COMPLETE;
714d2daff3dSHong Zhang   PetscFunctionReturn(0);
715d2daff3dSHong Zhang }
716d2daff3dSHong Zhang 
717*55de54fcSHong Zhang static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X)
718*55de54fcSHong Zhang {
719*55de54fcSHong Zhang   TS_RK            *rk = (TS_RK*)ts->data;
720*55de54fcSHong Zhang   PetscInt         s  = rk->tableau->s,p = rk->tableau->p,i,j;
721*55de54fcSHong Zhang   PetscReal        h;
722*55de54fcSHong Zhang   PetscReal        tt,t;
723*55de54fcSHong Zhang   PetscScalar      *b;
724*55de54fcSHong Zhang   const PetscReal  *B = rk->tableau->binterp;
725*55de54fcSHong Zhang   PetscErrorCode   ierr;
726*55de54fcSHong Zhang 
727*55de54fcSHong Zhang   PetscFunctionBegin;
728*55de54fcSHong Zhang   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
729*55de54fcSHong Zhang 
730*55de54fcSHong Zhang   switch (rk->status) {
731*55de54fcSHong Zhang     case TS_STEP_INCOMPLETE:
732*55de54fcSHong Zhang     case TS_STEP_PENDING:
733*55de54fcSHong Zhang       h = ts->time_step;
734*55de54fcSHong Zhang       t = (itime - ts->ptime)/h;
735*55de54fcSHong Zhang       break;
736*55de54fcSHong Zhang     case TS_STEP_COMPLETE:
737*55de54fcSHong Zhang       h = ts->ptime - ts->ptime_prev;
738*55de54fcSHong Zhang       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
739*55de54fcSHong Zhang       break;
740*55de54fcSHong Zhang     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
741*55de54fcSHong Zhang   }
742*55de54fcSHong Zhang   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
743*55de54fcSHong Zhang   for (i=0; i<s; i++) b[i] = 0;
744*55de54fcSHong Zhang   for (j=0,tt=t; j<p; j++,tt*=t) {
745*55de54fcSHong Zhang     for (i=0; i<s; i++) {
746*55de54fcSHong Zhang       b[i]  += h * B[i*p+j] * tt;
747*55de54fcSHong Zhang     }
748*55de54fcSHong Zhang   }
749*55de54fcSHong Zhang   ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr);
750*55de54fcSHong Zhang   ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr);
751*55de54fcSHong Zhang   ierr = PetscFree(b);CHKERRQ(ierr);
752*55de54fcSHong Zhang   PetscFunctionReturn(0);
753*55de54fcSHong Zhang }
754*55de54fcSHong Zhang 
755*55de54fcSHong Zhang /*------------------------------------------------------------*/
756*55de54fcSHong Zhang 
757be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts)
758be5899b3SLisandro Dalcin {
759be5899b3SLisandro Dalcin   TS_RK          *rk = (TS_RK*)ts->data;
760be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
761be5899b3SLisandro Dalcin   PetscErrorCode ierr;
762be5899b3SLisandro Dalcin 
763be5899b3SLisandro Dalcin   PetscFunctionBegin;
764be5899b3SLisandro Dalcin   if (!tab) PetscFunctionReturn(0);
765be5899b3SLisandro Dalcin   ierr = PetscFree(rk->work);CHKERRQ(ierr);
766be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr);
767be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr);
768be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr);
769be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr);
770*55de54fcSHong Zhang   ierr = PetscFree(rk->YdotRHS_fast);CHKERRQ(ierr);
771*55de54fcSHong Zhang   ierr = PetscFree(rk->YdotRHS_slow);CHKERRQ(ierr);
772be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
773be5899b3SLisandro Dalcin }
774be5899b3SLisandro Dalcin 
775277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts)
776e4dd521cSBarry Smith {
7775f70b478SJed Brown   TS_RK         *rk = (TS_RK*)ts->data;
7786849ba73SBarry Smith   PetscErrorCode ierr;
779e4dd521cSBarry Smith 
780e4dd521cSBarry Smith   PetscFunctionBegin;
781be5899b3SLisandro Dalcin   ierr = TSRKTableauReset(ts);CHKERRQ(ierr);
7820f7a1166SHong Zhang   ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr);
783*55de54fcSHong Zhang   ierr = PetscTryMethod(ts,"TSReset_MRKSPLIT_C",(TS),(ts));CHKERRQ(ierr);
784*55de54fcSHong Zhang   ierr = PetscTryMethod(ts,"TSReset_MRKNONSPLIT_C",(TS),(ts));CHKERRQ(ierr);
785277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
786e4dd521cSBarry Smith }
787277b19d0SLisandro Dalcin 
788f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx)
789f68a32c8SEmil Constantinescu {
790f68a32c8SEmil Constantinescu   PetscFunctionBegin;
791f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
792f68a32c8SEmil Constantinescu }
793f68a32c8SEmil Constantinescu 
794f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
795f68a32c8SEmil Constantinescu {
796f68a32c8SEmil Constantinescu   PetscFunctionBegin;
797f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
798f68a32c8SEmil Constantinescu }
799f68a32c8SEmil Constantinescu 
800f68a32c8SEmil Constantinescu 
801f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx)
802f68a32c8SEmil Constantinescu {
803f68a32c8SEmil Constantinescu   PetscFunctionBegin;
804f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
805f68a32c8SEmil Constantinescu }
806f68a32c8SEmil Constantinescu 
807f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
808f68a32c8SEmil Constantinescu {
809f68a32c8SEmil Constantinescu 
810f68a32c8SEmil Constantinescu   PetscFunctionBegin;
811f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
812f68a32c8SEmil Constantinescu }
813c235aa8dSHong Zhang /*
814d2daff3dSHong Zhang static PetscErrorCode RKSetAdjCoe(RKTableau tab)
815d2daff3dSHong Zhang {
816d2daff3dSHong Zhang   PetscReal *A,*b;
817d2daff3dSHong Zhang   PetscInt        s,i,j;
818d2daff3dSHong Zhang   PetscErrorCode  ierr;
819d2daff3dSHong Zhang 
820d2daff3dSHong Zhang   PetscFunctionBegin;
821d2daff3dSHong Zhang   s    = tab->s;
822d2daff3dSHong Zhang   ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr);
823d2daff3dSHong Zhang 
824d2daff3dSHong Zhang   for (i=0; i<s; i++)
825d2daff3dSHong Zhang     for (j=0; j<s; j++) {
826d2daff3dSHong Zhang       A[i*s+j] = (tab->b[s-1-i]==0) ? 0: -tab->A[s-1-i+(s-1-j)*s] * tab->b[s-1-j] / tab->b[s-1-i];
827d2daff3dSHong Zhang       ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr);
828d2daff3dSHong Zhang     }
829d2daff3dSHong Zhang 
830d2daff3dSHong Zhang   for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i];
831d2daff3dSHong Zhang 
832d2daff3dSHong Zhang   ierr  = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
833d2daff3dSHong Zhang   ierr  = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
834d2daff3dSHong Zhang   ierr  = PetscFree2(A,b);CHKERRQ(ierr);
835d2daff3dSHong Zhang   PetscFunctionReturn(0);
836d2daff3dSHong Zhang }
837c235aa8dSHong Zhang */
838be5899b3SLisandro Dalcin 
839be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts)
840be5899b3SLisandro Dalcin {
841be5899b3SLisandro Dalcin   TS_RK          *rk  = (TS_RK*)ts->data;
842be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
843be5899b3SLisandro Dalcin   PetscErrorCode ierr;
844be5899b3SLisandro Dalcin 
845be5899b3SLisandro Dalcin   PetscFunctionBegin;
846be5899b3SLisandro Dalcin   ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr);
847be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr);
848be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr);
849be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
850be5899b3SLisandro Dalcin }
851be5899b3SLisandro Dalcin 
852f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts)
853f68a32c8SEmil Constantinescu {
854f68a32c8SEmil Constantinescu   TS_RK          *rk = (TS_RK*)ts->data;
855f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
856f68a32c8SEmil Constantinescu   DM             dm;
857f68a32c8SEmil Constantinescu 
858f68a32c8SEmil Constantinescu   PetscFunctionBegin;
8593dd83b38SBarry Smith   ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr);
860be5899b3SLisandro Dalcin   ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);
8610f7a1166SHong Zhang   if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */
8620f7a1166SHong Zhang     ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr);
8630f7a1166SHong Zhang   }
864f68a32c8SEmil Constantinescu   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
865f68a32c8SEmil Constantinescu   ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
866f68a32c8SEmil Constantinescu   ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
867*55de54fcSHong Zhang   ierr = PetscTryMethod(ts,"TSSetUp_MRKSPLIT_C",(TS),(ts));CHKERRQ(ierr);
868*55de54fcSHong Zhang   ierr = PetscTryMethod(ts,"TSSetUp_MRKNONSPLIT_C",(TS),(ts));CHKERRQ(ierr);
869e4dd521cSBarry Smith   PetscFunctionReturn(0);
870e4dd521cSBarry Smith }
871d2daff3dSHong Zhang 
8722c9cb062Sluzhanghpp /*
873*55de54fcSHong Zhang   construct a database to choose single-step RK, or nonsplit version of mutirate RK or multirate RK with RHS splits
8742c9cb062Sluzhanghpp */
875*55de54fcSHong Zhang const char *const TSMRKTypes[] = {"NONE","NONSPLIT","SPLIT","TSMRKType","TSMRK",0};
876e4dd521cSBarry Smith 
8774416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts)
878e4dd521cSBarry Smith {
879be5899b3SLisandro Dalcin   TS_RK          *rk = (TS_RK*)ts->data;
880dfbe8321SBarry Smith   PetscErrorCode ierr;
881262ff7bbSBarry Smith 
882e4dd521cSBarry Smith   PetscFunctionBegin;
883e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr);
884f68a32c8SEmil Constantinescu   {
885f68a32c8SEmil Constantinescu     RKTableauLink link;
886f68a32c8SEmil Constantinescu     PetscInt      count,choice;
887f68a32c8SEmil Constantinescu     PetscBool     flg;
888f68a32c8SEmil Constantinescu     const char    **namelist;
889*55de54fcSHong Zhang     PetscInt      mrktype = 0;
8902c9cb062Sluzhanghpp 
891f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) ;
892f489ac74SBarry Smith     ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr);
893f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
894*55de54fcSHong Zhang     ierr = PetscOptionsEList("-ts_rk_multirate_type","Use Multirate or partioned RHS Multirate or single rate RK method","TSRKSetMultirateType",TSMRKTypes,3,TSMRKTypes[0],&mrktype,&flg);CHKERRQ(ierr);
895*55de54fcSHong Zhang     if (flg) {ierr = TSRKSetMultirateType(ts,mrktype);CHKERRQ(ierr);}
896be5899b3SLisandro Dalcin     ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr);
897be5899b3SLisandro Dalcin     if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);}
898f68a32c8SEmil Constantinescu     ierr = PetscFree(namelist);CHKERRQ(ierr);
899f68a32c8SEmil Constantinescu   }
900262ff7bbSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
9012c9cb062Sluzhanghpp   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Multirate methods options","");CHKERRQ(ierr);
9022c9cb062Sluzhanghpp   ierr = PetscOptionsInt("-ts_rk_dtratio","time step ratio between slow and fast","",rk->dtratio,&rk->dtratio,NULL);CHKERRQ(ierr);
9032c9cb062Sluzhanghpp   ierr = PetscOptionsEnd();CHKERRQ(ierr);
904e4dd521cSBarry Smith   PetscFunctionReturn(0);
905e4dd521cSBarry Smith }
906e4dd521cSBarry Smith 
9075f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
908e4dd521cSBarry Smith {
9095f70b478SJed Brown   TS_RK          *rk = (TS_RK*)ts->data;
9108370ee3bSLisandro Dalcin   PetscBool      iascii;
911dfbe8321SBarry Smith   PetscErrorCode ierr;
9122cdf8120Sasbjorn 
913e4dd521cSBarry Smith   PetscFunctionBegin;
914251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
9158370ee3bSLisandro Dalcin   if (iascii) {
9169c334d8fSLisandro Dalcin     RKTableau tab  = rk->tableau;
917f68a32c8SEmil Constantinescu     TSRKType  rktype;
918f68a32c8SEmil Constantinescu     char      buf[512];
919f68a32c8SEmil Constantinescu     ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr);
920efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  RK type %s\n",rktype);CHKERRQ(ierr);
9210c37eb8eSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  Order: %D\n",tab->order);CHKERRQ(ierr);
9220c37eb8eSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr);
923f68a32c8SEmil Constantinescu     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
924f68a32c8SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa c = %s\n",buf);CHKERRQ(ierr);
9258370ee3bSLisandro Dalcin   }
926f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
927f68a32c8SEmil Constantinescu }
928f68a32c8SEmil Constantinescu 
929f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer)
930f68a32c8SEmil Constantinescu {
931f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
9329c334d8fSLisandro Dalcin   TSAdapt        adapt;
933f68a32c8SEmil Constantinescu 
934f68a32c8SEmil Constantinescu   PetscFunctionBegin;
9359c334d8fSLisandro Dalcin   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
9369c334d8fSLisandro Dalcin   ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr);
937f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
938f68a32c8SEmil Constantinescu }
939f68a32c8SEmil Constantinescu 
940f68a32c8SEmil Constantinescu /*@C
941f68a32c8SEmil Constantinescu   TSRKSetType - Set the type of RK scheme
942f68a32c8SEmil Constantinescu 
943f68a32c8SEmil Constantinescu   Logically collective
944f68a32c8SEmil Constantinescu 
945f68a32c8SEmil Constantinescu   Input Parameter:
946f68a32c8SEmil Constantinescu +  ts - timestepping context
947f68a32c8SEmil Constantinescu -  rktype - type of RK-scheme
948f68a32c8SEmil Constantinescu 
949287c2655SBarry Smith   Options Database:
9509bd3e852SBarry Smith .   -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs>
951287c2655SBarry Smith 
952f68a32c8SEmil Constantinescu   Level: intermediate
953f68a32c8SEmil Constantinescu 
954287c2655SBarry Smith .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS
955f68a32c8SEmil Constantinescu @*/
956f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype)
957f68a32c8SEmil Constantinescu {
958f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
959f68a32c8SEmil Constantinescu 
960f68a32c8SEmil Constantinescu   PetscFunctionBegin;
961f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
962b92453a8SLisandro Dalcin   PetscValidCharPointer(rktype,2);
963f68a32c8SEmil Constantinescu   ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr);
964f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
965f68a32c8SEmil Constantinescu }
966f68a32c8SEmil Constantinescu 
967f68a32c8SEmil Constantinescu /*@C
968f68a32c8SEmil Constantinescu   TSRKGetType - Get the type of RK scheme
969f68a32c8SEmil Constantinescu 
970f68a32c8SEmil Constantinescu   Logically collective
971f68a32c8SEmil Constantinescu 
972f68a32c8SEmil Constantinescu   Input Parameter:
973f68a32c8SEmil Constantinescu .  ts - timestepping context
974f68a32c8SEmil Constantinescu 
975f68a32c8SEmil Constantinescu   Output Parameter:
976f68a32c8SEmil Constantinescu .  rktype - type of RK-scheme
977f68a32c8SEmil Constantinescu 
978f68a32c8SEmil Constantinescu   Level: intermediate
979f68a32c8SEmil Constantinescu 
980f68a32c8SEmil Constantinescu .seealso: TSRKGetType()
981f68a32c8SEmil Constantinescu @*/
982f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype)
983f68a32c8SEmil Constantinescu {
984f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
985f68a32c8SEmil Constantinescu 
986f68a32c8SEmil Constantinescu   PetscFunctionBegin;
987f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
988f68a32c8SEmil Constantinescu   ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr);
989f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
990f68a32c8SEmil Constantinescu }
991f68a32c8SEmil Constantinescu 
992560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype)
993f68a32c8SEmil Constantinescu {
994f68a32c8SEmil Constantinescu   TS_RK *rk = (TS_RK*)ts->data;
995f68a32c8SEmil Constantinescu 
996f68a32c8SEmil Constantinescu   PetscFunctionBegin;
997f68a32c8SEmil Constantinescu   *rktype = rk->tableau->name;
998f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
999f68a32c8SEmil Constantinescu }
10002c9cb062Sluzhanghpp 
1001560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype)
1002f68a32c8SEmil Constantinescu {
1003f68a32c8SEmil Constantinescu   TS_RK          *rk = (TS_RK*)ts->data;
1004f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
1005f68a32c8SEmil Constantinescu   PetscBool      match;
1006f68a32c8SEmil Constantinescu   RKTableauLink  link;
1007f68a32c8SEmil Constantinescu 
1008f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1009f68a32c8SEmil Constantinescu   if (rk->tableau) {
1010f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr);
1011f68a32c8SEmil Constantinescu     if (match) PetscFunctionReturn(0);
1012f68a32c8SEmil Constantinescu   }
1013f68a32c8SEmil Constantinescu   for (link = RKTableauList; link; link=link->next) {
1014f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr);
1015f68a32c8SEmil Constantinescu     if (match) {
1016be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);}
1017f68a32c8SEmil Constantinescu       rk->tableau = &link->tab;
1018be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);}
10192ffb9264SLisandro Dalcin       ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
1020f68a32c8SEmil Constantinescu       PetscFunctionReturn(0);
1021f68a32c8SEmil Constantinescu     }
1022f68a32c8SEmil Constantinescu   }
1023f68a32c8SEmil Constantinescu   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype);
1024e4dd521cSBarry Smith   PetscFunctionReturn(0);
1025e4dd521cSBarry Smith }
1026e4dd521cSBarry Smith 
1027ff22ae23SHong Zhang static PetscErrorCode  TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y)
1028ff22ae23SHong Zhang {
1029ff22ae23SHong Zhang   TS_RK *rk = (TS_RK*)ts->data;
1030ff22ae23SHong Zhang 
1031ff22ae23SHong Zhang   PetscFunctionBegin;
10320429704eSStefano Zampini   if (ns) *ns = rk->tableau->s;
1033d2daff3dSHong Zhang   if (Y)   *Y = rk->Y;
1034ff22ae23SHong Zhang   PetscFunctionReturn(0);
1035ff22ae23SHong Zhang }
1036ff22ae23SHong Zhang 
1037b3a6b972SJed Brown static PetscErrorCode TSDestroy_RK(TS ts)
1038b3a6b972SJed Brown {
1039b3a6b972SJed Brown   PetscErrorCode ierr;
1040b3a6b972SJed Brown 
1041b3a6b972SJed Brown   PetscFunctionBegin;
1042b3a6b972SJed Brown   ierr = TSReset_RK(ts);CHKERRQ(ierr);
1043b3a6b972SJed Brown   if (ts->dm) {
1044b3a6b972SJed Brown     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
1045b3a6b972SJed Brown     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
1046b3a6b972SJed Brown   }
1047b3a6b972SJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
1048b3a6b972SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr);
1049b3a6b972SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr);
10502c9cb062Sluzhanghpp   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirateType_C",NULL);CHKERRQ(ierr);
1051b3a6b972SJed Brown   PetscFunctionReturn(0);
1052b3a6b972SJed Brown }
1053ff22ae23SHong Zhang 
1054ebe3b25bSBarry Smith /*MC
1055f68a32c8SEmil Constantinescu       TSRK - ODE and DAE solver using Runge-Kutta schemes
1056ebe3b25bSBarry Smith 
1057f68a32c8SEmil Constantinescu   The user should provide the right hand side of the equation
1058f68a32c8SEmil Constantinescu   using TSSetRHSFunction().
1059ebe3b25bSBarry Smith 
1060f68a32c8SEmil Constantinescu   Notes:
106198164e13SLisandro Dalcin   The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type
1062ebe3b25bSBarry Smith 
1063d41222bbSBarry Smith   Level: beginner
1064d41222bbSBarry Smith 
1065f68a32c8SEmil Constantinescu .seealso:  TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3,
10662c9cb062Sluzhanghpp            TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister(), TSRKSetMultirateType()
1067ebe3b25bSBarry Smith 
1068ebe3b25bSBarry Smith M*/
10698cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
1070e4dd521cSBarry Smith {
10711566a47fSLisandro Dalcin   TS_RK          *rk;
1072dfbe8321SBarry Smith   PetscErrorCode ierr;
1073e4dd521cSBarry Smith 
1074e4dd521cSBarry Smith   PetscFunctionBegin;
1075f68a32c8SEmil Constantinescu   ierr = TSRKInitializePackage();CHKERRQ(ierr);
1076f68a32c8SEmil Constantinescu 
1077f68a32c8SEmil Constantinescu   ts->ops->reset          = TSReset_RK;
10785f70b478SJed Brown   ts->ops->destroy        = TSDestroy_RK;
10795f70b478SJed Brown   ts->ops->view           = TSView_RK;
1080f68a32c8SEmil Constantinescu   ts->ops->load           = TSLoad_RK;
1081f68a32c8SEmil Constantinescu   ts->ops->setup          = TSSetUp_RK;
108242f2b339SBarry Smith   ts->ops->adjointsetup   = TSAdjointSetUp_RK;
1083f68a32c8SEmil Constantinescu   ts->ops->interpolate    = TSInterpolate_RK;
10842c9cb062Sluzhanghpp   ts->ops->step           = TSStep_RK;
1085f68a32c8SEmil Constantinescu   ts->ops->evaluatestep   = TSEvaluateStep_RK;
1086fecfb714SLisandro Dalcin   ts->ops->rollback       = TSRollBack_RK;
1087f68a32c8SEmil Constantinescu   ts->ops->setfromoptions = TSSetFromOptions_RK;
1088ff22ae23SHong Zhang   ts->ops->getstages      = TSGetStages_RK;
108942f2b339SBarry Smith   ts->ops->adjointstep    = TSAdjointStep_RK;
1090e4dd521cSBarry Smith 
10910f7a1166SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_RK;
10920f7a1166SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_RK;
10930f7a1166SHong Zhang 
10941566a47fSLisandro Dalcin   ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr);
10951566a47fSLisandro Dalcin   ts->data = (void*)rk;
1096e4dd521cSBarry Smith 
1097f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr);
1098f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr);
10992c9cb062Sluzhanghpp   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirateType_C",TSRKSetMultirateType);CHKERRQ(ierr);
1100be5899b3SLisandro Dalcin 
1101be5899b3SLisandro Dalcin   ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr);
1102e4dd521cSBarry Smith   PetscFunctionReturn(0);
1103e4dd521cSBarry Smith }
1104*55de54fcSHong Zhang 
1105*55de54fcSHong Zhang 
1106*55de54fcSHong Zhang /*------------ Multirate support for partitioned systems --------------*/
1107*55de54fcSHong Zhang 
1108*55de54fcSHong Zhang static PetscErrorCode TSSetUp_MRKSPLIT(TS ts)
1109*55de54fcSHong Zhang {
1110*55de54fcSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
1111*55de54fcSHong Zhang   RKTableau      tab = rk->tableau;
1112*55de54fcSHong Zhang   DM             dm,subdm,newdm;
1113*55de54fcSHong Zhang   PetscErrorCode ierr;
1114*55de54fcSHong Zhang 
1115*55de54fcSHong Zhang   PetscFunctionBegin;
1116*55de54fcSHong Zhang   ierr = TSRHSSplitGetIS(ts,"slow",&rk->is_slow);CHKERRQ(ierr);
1117*55de54fcSHong Zhang   ierr = TSRHSSplitGetIS(ts,"fast",&rk->is_fast);CHKERRQ(ierr);
1118*55de54fcSHong Zhang   if (!rk->is_slow || !rk->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 bsi");
1119*55de54fcSHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"slow",&rk->subts_slow);CHKERRQ(ierr);
1120*55de54fcSHong Zhang   ierr = TSRHSSplitGetSubTS(ts,"fast",&rk->subts_fast);CHKERRQ(ierr);
1121*55de54fcSHong Zhang   if (!rk->subts_slow || !rk->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");
1122*55de54fcSHong Zhang 
1123*55de54fcSHong Zhang   /* Only copy */
1124*55de54fcSHong Zhang   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1125*55de54fcSHong Zhang   ierr = DMClone(dm,&newdm);CHKERRQ(ierr);
1126*55de54fcSHong Zhang   ierr = TSGetDM(rk->subts_fast,&subdm);CHKERRQ(ierr);
1127*55de54fcSHong Zhang   ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr);
1128*55de54fcSHong Zhang   ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr);
1129*55de54fcSHong Zhang   ierr = TSSetDM(rk->subts_fast,newdm);CHKERRQ(ierr);
1130*55de54fcSHong Zhang   ierr = DMDestroy(&newdm);CHKERRQ(ierr);
1131*55de54fcSHong Zhang   ierr = DMClone(dm,&newdm);CHKERRQ(ierr);
1132*55de54fcSHong Zhang   ierr = TSGetDM(rk->subts_slow,&subdm);CHKERRQ(ierr);
1133*55de54fcSHong Zhang   ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr);
1134*55de54fcSHong Zhang   ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr);
1135*55de54fcSHong Zhang   ierr = TSSetDM(rk->subts_slow,newdm);CHKERRQ(ierr);
1136*55de54fcSHong Zhang   ierr = DMDestroy(&newdm);CHKERRQ(ierr);
1137*55de54fcSHong Zhang 
1138*55de54fcSHong Zhang   ierr = PetscMalloc1(tab->s,&rk->YdotRHS_fast);CHKERRQ(ierr);
1139*55de54fcSHong Zhang   ierr = PetscMalloc1(tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr);
1140*55de54fcSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&rk->X0);CHKERRQ(ierr);
1141*55de54fcSHong Zhang   PetscFunctionReturn(0);
1142*55de54fcSHong Zhang }
1143*55de54fcSHong Zhang 
1144*55de54fcSHong Zhang static PetscErrorCode TSReset_MRKSPLIT(TS ts)
1145*55de54fcSHong Zhang {
1146*55de54fcSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
1147*55de54fcSHong Zhang   PetscErrorCode ierr;
1148*55de54fcSHong Zhang 
1149*55de54fcSHong Zhang   PetscFunctionBegin;
1150*55de54fcSHong Zhang   ierr = VecDestroy(&rk->X0);CHKERRQ(ierr);
1151*55de54fcSHong Zhang   PetscFunctionReturn(0);
1152*55de54fcSHong Zhang }
1153*55de54fcSHong Zhang 
1154*55de54fcSHong Zhang static PetscErrorCode TSSetUp_MRKNONSPLIT(TS ts)
1155*55de54fcSHong Zhang {
1156*55de54fcSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
1157*55de54fcSHong Zhang   RKTableau       tab  = rk->tableau;
1158*55de54fcSHong Zhang   PetscErrorCode ierr;
1159*55de54fcSHong Zhang 
1160*55de54fcSHong Zhang   PetscFunctionBegin;
1161*55de54fcSHong Zhang   ierr = TSRHSSplitGetIS(ts,"slow",&rk->is_slow);CHKERRQ(ierr);
1162*55de54fcSHong Zhang   ierr = TSRHSSplitGetIS(ts,"fast",&rk->is_fast);CHKERRQ(ierr);
1163*55de54fcSHong Zhang   if (!rk->is_slow || !rk->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 bsi");
1164*55de54fcSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&rk->X0);CHKERRQ(ierr);
1165*55de54fcSHong Zhang   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr);
1166*55de54fcSHong Zhang   PetscFunctionReturn(0);
1167*55de54fcSHong Zhang }
1168*55de54fcSHong Zhang 
1169*55de54fcSHong Zhang static PetscErrorCode TSReset_MRKNONSPLIT(TS ts)
1170*55de54fcSHong Zhang {
1171*55de54fcSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
1172*55de54fcSHong Zhang   RKTableau      tab  = rk->tableau;
1173*55de54fcSHong Zhang   PetscErrorCode ierr;
1174*55de54fcSHong Zhang 
1175*55de54fcSHong Zhang   PetscFunctionBegin;
1176*55de54fcSHong Zhang   ierr = VecDestroy(&rk->X0);CHKERRQ(ierr);
1177*55de54fcSHong Zhang   ierr = VecDestroyVecs(tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr);
1178*55de54fcSHong Zhang   PetscFunctionReturn(0);
1179*55de54fcSHong Zhang }
1180*55de54fcSHong Zhang 
1181*55de54fcSHong Zhang static PetscErrorCode TSInterpolate_MRKNONSPLIT(TS ts,PetscReal itime,Vec X)
1182*55de54fcSHong Zhang {
1183*55de54fcSHong Zhang   TS_RK            *rk = (TS_RK*)ts->data;
1184*55de54fcSHong Zhang   PetscInt         s  = rk->tableau->s,p = rk->tableau->p,i,j;
1185*55de54fcSHong Zhang   PetscReal        h;
1186*55de54fcSHong Zhang   PetscReal        tt,t;
1187*55de54fcSHong Zhang   PetscScalar      *b;
1188*55de54fcSHong Zhang   const PetscReal  *B = rk->tableau->binterp;
1189*55de54fcSHong Zhang   PetscErrorCode   ierr;
1190*55de54fcSHong Zhang 
1191*55de54fcSHong Zhang   PetscFunctionBegin;
1192*55de54fcSHong Zhang   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
1193*55de54fcSHong Zhang 
1194*55de54fcSHong Zhang   switch (rk->status) {
1195*55de54fcSHong Zhang     case TS_STEP_INCOMPLETE:
1196*55de54fcSHong Zhang     case TS_STEP_PENDING:
1197*55de54fcSHong Zhang       h = ts->time_step;
1198*55de54fcSHong Zhang       t = (itime - ts->ptime)/h;
1199*55de54fcSHong Zhang       break;
1200*55de54fcSHong Zhang     case TS_STEP_COMPLETE:
1201*55de54fcSHong Zhang       h = ts->ptime - ts->ptime_prev;
1202*55de54fcSHong Zhang       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
1203*55de54fcSHong Zhang       break;
1204*55de54fcSHong Zhang     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
1205*55de54fcSHong Zhang   }
1206*55de54fcSHong Zhang   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
1207*55de54fcSHong Zhang   for (i=0; i<s; i++) b[i] = 0;
1208*55de54fcSHong Zhang   for (j=0,tt=t; j<p; j++,tt*=t) {
1209*55de54fcSHong Zhang     for (i=0; i<s; i++) {
1210*55de54fcSHong Zhang       b[i]  += h * B[i*p+j] * tt;
1211*55de54fcSHong Zhang     }
1212*55de54fcSHong Zhang   }
1213*55de54fcSHong Zhang   ierr = VecCopy(rk->X0,X);CHKERRQ(ierr);
1214*55de54fcSHong Zhang   ierr = VecMAXPY(X,s,b,rk->YdotRHS_slow);CHKERRQ(ierr);
1215*55de54fcSHong Zhang   ierr = PetscFree(b);CHKERRQ(ierr);
1216*55de54fcSHong Zhang   PetscFunctionReturn(0);
1217*55de54fcSHong Zhang }
1218*55de54fcSHong Zhang 
1219*55de54fcSHong Zhang 
1220*55de54fcSHong Zhang static PetscErrorCode TSInterpolate_MRKSPLIT(TS ts,PetscReal itime,Vec X)
1221*55de54fcSHong Zhang {
1222*55de54fcSHong Zhang   TS_RK            *rk = (TS_RK*)ts->data;
1223*55de54fcSHong Zhang   PetscInt         s  = rk->tableau->s,p = rk->tableau->p,i,j;
1224*55de54fcSHong Zhang   Vec              Yslow;    /* vector holds the slow components in Y[0] */
1225*55de54fcSHong Zhang   PetscReal        h;
1226*55de54fcSHong Zhang   PetscReal        tt,t;
1227*55de54fcSHong Zhang   PetscScalar      *b;
1228*55de54fcSHong Zhang   const PetscReal  *B = rk->tableau->binterp;
1229*55de54fcSHong Zhang   PetscErrorCode   ierr;
1230*55de54fcSHong Zhang 
1231*55de54fcSHong Zhang   PetscFunctionBegin;
1232*55de54fcSHong Zhang   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
1233*55de54fcSHong Zhang 
1234*55de54fcSHong Zhang   switch (rk->status) {
1235*55de54fcSHong Zhang     case TS_STEP_INCOMPLETE:
1236*55de54fcSHong Zhang     case TS_STEP_PENDING:
1237*55de54fcSHong Zhang       h = ts->time_step;
1238*55de54fcSHong Zhang       t = (itime - ts->ptime)/h;
1239*55de54fcSHong Zhang       break;
1240*55de54fcSHong Zhang     case TS_STEP_COMPLETE:
1241*55de54fcSHong Zhang       h = ts->ptime - ts->ptime_prev;
1242*55de54fcSHong Zhang       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
1243*55de54fcSHong Zhang       break;
1244*55de54fcSHong Zhang     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
1245*55de54fcSHong Zhang   }
1246*55de54fcSHong Zhang   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
1247*55de54fcSHong Zhang   for (i=0; i<s; i++) b[i] = 0;
1248*55de54fcSHong Zhang   for (j=0,tt=t; j<p; j++,tt*=t) {
1249*55de54fcSHong Zhang     for (i=0; i<s; i++) {
1250*55de54fcSHong Zhang       b[i]  += h * B[i*p+j] * tt;
1251*55de54fcSHong Zhang     }
1252*55de54fcSHong Zhang   }
1253*55de54fcSHong Zhang   ierr = VecGetSubVector(rk->X0,rk->is_slow,&Yslow);CHKERRQ(ierr);
1254*55de54fcSHong Zhang   ierr = VecCopy(Yslow,X);CHKERRQ(ierr);
1255*55de54fcSHong Zhang   ierr = VecMAXPY(X,s,b,rk->YdotRHS_slow);CHKERRQ(ierr);
1256*55de54fcSHong Zhang   ierr = VecRestoreSubVector(rk->X0,rk->is_slow,&Yslow);CHKERRQ(ierr);
1257*55de54fcSHong Zhang   ierr = PetscFree(b);CHKERRQ(ierr);
1258*55de54fcSHong Zhang   PetscFunctionReturn(0);
1259*55de54fcSHong Zhang }
1260*55de54fcSHong Zhang 
1261*55de54fcSHong Zhang /*
1262*55de54fcSHong Zhang  This is for the nonsplit version of multirate RK method
1263*55de54fcSHong Zhang  The step completion formula is
1264*55de54fcSHong Zhang 
1265*55de54fcSHong Zhang  x1 = x0 + h b^T YdotRHS
1266*55de54fcSHong Zhang 
1267*55de54fcSHong Zhang */
1268*55de54fcSHong Zhang static PetscErrorCode TSEvaluateStep_MRKNONSPLIT(TS ts,PetscInt order,Vec X,PetscBool *done)
1269*55de54fcSHong Zhang {
1270*55de54fcSHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
1271*55de54fcSHong Zhang   RKTableau       tab = rk->tableau;
1272*55de54fcSHong Zhang   Vec             Xtemp;                          /* temporary work vector for X                                   */
1273*55de54fcSHong Zhang   PetscScalar     *w = rk->work;
1274*55de54fcSHong Zhang   PetscScalar     *x_ptr,*xtemp_ptr;              /* location to put the pointer to X and Xtemp respectively       */
1275*55de54fcSHong Zhang   PetscReal       h = ts->time_step;
1276*55de54fcSHong Zhang   PetscInt        s = tab->s,j;
1277*55de54fcSHong Zhang   PetscInt        len_slow,len_fast;              /* the number of slow components and fast components respectively */
1278*55de54fcSHong Zhang   const PetscInt  *is_slow,*is_fast;              /* the index for slow components and fast components respectively */
1279*55de54fcSHong Zhang   PetscErrorCode  ierr;
1280*55de54fcSHong Zhang 
1281*55de54fcSHong Zhang   PetscFunctionBegin;
1282*55de54fcSHong Zhang   ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
1283*55de54fcSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&Xtemp);CHKERRQ(ierr);
1284*55de54fcSHong Zhang   ierr = VecCopy(ts->vec_sol,Xtemp);CHKERRQ(ierr);
1285*55de54fcSHong Zhang   if (rk->slow) {
1286*55de54fcSHong Zhang     for (j=0; j<s; j++) w[j] = h*tab->b[j];
1287*55de54fcSHong Zhang     /* update the value of slow components, and discard the updating value of fast components */
1288*55de54fcSHong Zhang     ierr = VecMAXPY(Xtemp,s,w,rk->YdotRHS_slow);CHKERRQ(ierr);
1289*55de54fcSHong Zhang     ierr = VecGetArray(X,&x_ptr);CHKERRQ(ierr);
1290*55de54fcSHong Zhang     ierr = VecGetArray(Xtemp,&xtemp_ptr);CHKERRQ(ierr);
1291*55de54fcSHong Zhang     ierr = ISGetSize(rk->is_slow,&len_slow);CHKERRQ(ierr);
1292*55de54fcSHong Zhang     ierr = ISGetIndices(rk->is_slow,&is_slow);CHKERRQ(ierr);
1293*55de54fcSHong Zhang     ierr = ISRestoreIndices(rk->is_slow,&is_slow);CHKERRQ(ierr);
1294*55de54fcSHong Zhang     for (j=0; j<len_slow; j++) x_ptr[is_slow[j]] = xtemp_ptr[is_slow[j]];
1295*55de54fcSHong Zhang     ierr = VecRestoreArray(X,&x_ptr);CHKERRQ(ierr);
1296*55de54fcSHong Zhang     ierr = VecRestoreArray(Xtemp,&xtemp_ptr);CHKERRQ(ierr);
1297*55de54fcSHong Zhang   } else {
1298*55de54fcSHong Zhang     /* update the value of fast components, and discard the updating value of slow components */
1299*55de54fcSHong Zhang     for (j=0; j<s; j++) w[j] = h/rk->dtratio*tab->b[j];
1300*55de54fcSHong Zhang     ierr = VecMAXPY(Xtemp,s,w,rk->YdotRHS);CHKERRQ(ierr);
1301*55de54fcSHong Zhang     ierr = VecGetArray(X,&x_ptr);CHKERRQ(ierr);
1302*55de54fcSHong Zhang     ierr = VecGetArray(Xtemp,&xtemp_ptr);CHKERRQ(ierr);
1303*55de54fcSHong Zhang     ierr = ISGetSize(rk->is_fast,&len_fast);CHKERRQ(ierr);
1304*55de54fcSHong Zhang     ierr = ISGetIndices(rk->is_fast,&is_fast);CHKERRQ(ierr);
1305*55de54fcSHong Zhang     ierr = ISRestoreIndices(rk->is_fast,&is_fast);CHKERRQ(ierr);
1306*55de54fcSHong Zhang     for (j=0; j<len_fast; j++) x_ptr[is_fast[j]] = xtemp_ptr[is_fast[j]];
1307*55de54fcSHong Zhang     ierr = VecRestoreArray(X,&x_ptr);CHKERRQ(ierr);
1308*55de54fcSHong Zhang     ierr = VecRestoreArray(Xtemp,&xtemp_ptr);CHKERRQ(ierr);
1309*55de54fcSHong Zhang   }
1310*55de54fcSHong Zhang   /* free temporary vector space */
1311*55de54fcSHong Zhang   ierr = VecDestroy(&Xtemp);CHKERRQ(ierr);
1312*55de54fcSHong Zhang   PetscFunctionReturn(0);
1313*55de54fcSHong Zhang }
1314*55de54fcSHong Zhang 
1315*55de54fcSHong Zhang static PetscErrorCode TSStep_MRKNONSPLIT(TS ts)
1316*55de54fcSHong Zhang {
1317*55de54fcSHong Zhang   TS_RK             *rk = (TS_RK*)ts->data;
1318*55de54fcSHong Zhang   RKTableau         tab  = rk->tableau;
1319*55de54fcSHong Zhang   Vec               *Y = rk->Y,*YdotRHS = rk->YdotRHS,*YdotRHS_slow = rk->YdotRHS_slow;
1320*55de54fcSHong Zhang   Vec               stage_fast,stage_slow; /* vectors store the stage value of fast and slow components respectively */
1321*55de54fcSHong Zhang   Vec               Yslow,Sslow;           /* vectors store the previous and new time solution of slow components */
1322*55de54fcSHong Zhang   const PetscInt    s = tab->s;
1323*55de54fcSHong Zhang   const PetscReal   *A = tab->A,*c = tab->c;
1324*55de54fcSHong Zhang   PetscScalar       *w = rk->work;
1325*55de54fcSHong Zhang   PetscInt          i,j,k;
1326*55de54fcSHong Zhang   PetscReal         next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step;
1327*55de54fcSHong Zhang   PetscErrorCode    ierr;
1328*55de54fcSHong Zhang 
1329*55de54fcSHong Zhang   PetscFunctionBegin;
1330*55de54fcSHong Zhang   rk->status = TS_STEP_INCOMPLETE;
1331*55de54fcSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&stage_slow);CHKERRQ(ierr);
1332*55de54fcSHong Zhang   ierr = VecCopy(ts->vec_sol,rk->X0);CHKERRQ(ierr);
1333*55de54fcSHong Zhang   for (i=0; i<s; i++) {
1334*55de54fcSHong Zhang     rk->stage_time = t + h*c[i];
1335*55de54fcSHong Zhang     ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr);
1336*55de54fcSHong Zhang     ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
1337*55de54fcSHong Zhang     for (j=0; j<i; j++) w[j] = h*A[i*s+j];
1338*55de54fcSHong Zhang     ierr = VecMAXPY(Y[i],i,w,YdotRHS_slow);CHKERRQ(ierr);
1339*55de54fcSHong Zhang     ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
1340*55de54fcSHong Zhang     /* compute the stage RHS */
1341*55de54fcSHong Zhang     ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr);
1342*55de54fcSHong Zhang   }
1343*55de54fcSHong Zhang   rk->slow = PETSC_TRUE;
1344*55de54fcSHong Zhang   ierr = TSEvaluateStep_MRKNONSPLIT(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
1345*55de54fcSHong Zhang   for (k=0; k<rk->dtratio; k++) {
1346*55de54fcSHong Zhang     for (i=0; i<s; i++) {
1347*55de54fcSHong Zhang       rk->stage_time = t + h/rk->dtratio*c[i];
1348*55de54fcSHong Zhang       ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr);
1349*55de54fcSHong Zhang       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); /* Only need the values for fast components */
1350*55de54fcSHong Zhang       /* guarantee the stage RHS for slow components do not change while updating the stage RHS for fast components */
1351*55de54fcSHong Zhang       /* update the fast components stage value, slow components will be overwritten */
1352*55de54fcSHong Zhang       for (j=0; j<i; j++) w[j] = h/rk->dtratio*A[i*s+j];
1353*55de54fcSHong Zhang       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
1354*55de54fcSHong Zhang       /* update the slow components stage value */
1355*55de54fcSHong Zhang       ierr = TSInterpolate_MRKNONSPLIT(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],stage_slow);CHKERRQ(ierr);
1356*55de54fcSHong Zhang       /* combine the update fast components stage value and slow components stage value to Y[i] */
1357*55de54fcSHong Zhang       ierr = VecGetSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr);
1358*55de54fcSHong Zhang       ierr = VecGetSubVector(stage_slow,rk->is_slow,&Sslow);CHKERRQ(ierr);
1359*55de54fcSHong Zhang       ierr = VecCopy(Sslow,Yslow);CHKERRQ(ierr);
1360*55de54fcSHong Zhang       ierr = VecRestoreSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr);
1361*55de54fcSHong Zhang       ierr = VecRestoreSubVector(stage_slow,rk->is_slow,&Sslow);CHKERRQ(ierr);
1362*55de54fcSHong Zhang       ierr = TSPostStage(ts,rk->stage_time,i,Y);CHKERRQ(ierr);
1363*55de54fcSHong Zhang       /* compute the stage RHS */
1364*55de54fcSHong Zhang       ierr = TSComputeRHSFunction(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
1365*55de54fcSHong Zhang     }
1366*55de54fcSHong Zhang     rk->slow = PETSC_FALSE;
1367*55de54fcSHong Zhang     /* update the fast components */
1368*55de54fcSHong Zhang     ierr = TSEvaluateStep_MRKNONSPLIT(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
1369*55de54fcSHong Zhang   }
1370*55de54fcSHong Zhang   ts->ptime += ts->time_step;
1371*55de54fcSHong Zhang   ts->time_step = next_time_step;
1372*55de54fcSHong Zhang   rk->status = TS_STEP_COMPLETE;
1373*55de54fcSHong Zhang   /* free memory of work vectors */
1374*55de54fcSHong Zhang   ierr = VecDestroy(&stage_fast);CHKERRQ(ierr);
1375*55de54fcSHong Zhang   PetscFunctionReturn(0);
1376*55de54fcSHong Zhang }
1377*55de54fcSHong Zhang 
1378*55de54fcSHong Zhang /*
1379*55de54fcSHong Zhang  This is for partitioned RHS multirate RK method
1380*55de54fcSHong Zhang  The step completion formula is
1381*55de54fcSHong Zhang 
1382*55de54fcSHong Zhang  x1 = x0 + h b^T YdotRHS
1383*55de54fcSHong Zhang 
1384*55de54fcSHong Zhang */
1385*55de54fcSHong Zhang static PetscErrorCode TSEvaluateStep_MRKSPLIT(TS ts,PetscInt order,Vec X,PetscBool *done)
1386*55de54fcSHong Zhang {
1387*55de54fcSHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
1388*55de54fcSHong Zhang   RKTableau       tab  = rk->tableau;
1389*55de54fcSHong Zhang   Vec             Xslow,Xfast;                  /* subvectors of X which store slow components and fast components respectively */
1390*55de54fcSHong Zhang   PetscScalar     *w = rk->work;
1391*55de54fcSHong Zhang   PetscReal       h = ts->time_step;
1392*55de54fcSHong Zhang   PetscInt        s = tab->s,j;
1393*55de54fcSHong Zhang   PetscErrorCode  ierr;
1394*55de54fcSHong Zhang 
1395*55de54fcSHong Zhang   PetscFunctionBegin;
1396*55de54fcSHong Zhang   ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
1397*55de54fcSHong Zhang   if (rk->slow) {
1398*55de54fcSHong Zhang     for (j=0; j<s; j++) w[j] = h*tab->b[j];
1399*55de54fcSHong Zhang     ierr = VecGetSubVector(ts->vec_sol,rk->is_slow,&Xslow);CHKERRQ(ierr);
1400*55de54fcSHong Zhang     ierr = VecMAXPY(Xslow,s,w,rk->YdotRHS_slow);CHKERRQ(ierr);
1401*55de54fcSHong Zhang     ierr = VecRestoreSubVector(ts->vec_sol,rk->is_slow,&Xslow);CHKERRQ(ierr);;
1402*55de54fcSHong Zhang   } else {
1403*55de54fcSHong Zhang     for (j=0; j<s; j++) w[j] = h/rk->dtratio*tab->b[j];
1404*55de54fcSHong Zhang     ierr = VecGetSubVector(X,rk->is_fast,&Xfast);CHKERRQ(ierr);
1405*55de54fcSHong Zhang     ierr = VecMAXPY(Xfast,s,w,rk->YdotRHS_fast);CHKERRQ(ierr);
1406*55de54fcSHong Zhang     ierr = VecRestoreSubVector(X,rk->is_fast,&Xfast);CHKERRQ(ierr);
1407*55de54fcSHong Zhang   }
1408*55de54fcSHong Zhang   PetscFunctionReturn(0);
1409*55de54fcSHong Zhang }
1410*55de54fcSHong Zhang 
1411*55de54fcSHong Zhang static PetscErrorCode TSStep_MRKSPLIT(TS ts)
1412*55de54fcSHong Zhang {
1413*55de54fcSHong Zhang   TS_RK             *rk = (TS_RK*)ts->data;
1414*55de54fcSHong Zhang   RKTableau         tab = rk->tableau;
1415*55de54fcSHong Zhang   Vec               *Y = rk->Y,*YdotRHS = rk->YdotRHS;
1416*55de54fcSHong Zhang   Vec               *YdotRHS_fast = rk->YdotRHS_fast,*YdotRHS_slow = rk->YdotRHS_slow;
1417*55de54fcSHong Zhang   Vec               Yslow,Yfast;                         /* subvectors store the stges of slow components and fast components respectively                           */
1418*55de54fcSHong Zhang   const PetscInt    s = tab->s;
1419*55de54fcSHong Zhang   const PetscReal   *A = tab->A,*c = tab->c;
1420*55de54fcSHong Zhang   PetscScalar       *w = rk->work;
1421*55de54fcSHong Zhang   PetscInt          i,j,k;
1422*55de54fcSHong Zhang   PetscReal         next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step;
1423*55de54fcSHong Zhang   PetscErrorCode    ierr;
1424*55de54fcSHong Zhang 
1425*55de54fcSHong Zhang   PetscFunctionBegin;
1426*55de54fcSHong Zhang   rk->status = TS_STEP_INCOMPLETE;
1427*55de54fcSHong Zhang   for (i=0; i<s; i++) {
1428*55de54fcSHong Zhang     ierr = VecGetSubVector(YdotRHS[i],rk->is_slow,&YdotRHS_slow[i]);CHKERRQ(ierr);
1429*55de54fcSHong Zhang     ierr = VecGetSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr);
1430*55de54fcSHong Zhang   }
1431*55de54fcSHong Zhang   ierr = VecCopy(ts->vec_sol,rk->X0);CHKERRQ(ierr);
1432*55de54fcSHong Zhang   /* propogate both slow and fast components using large time steps */
1433*55de54fcSHong Zhang   for (i=0; i<s; i++) {
1434*55de54fcSHong Zhang     rk->stage_time = t + h*c[i];
1435*55de54fcSHong Zhang     ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr);
1436*55de54fcSHong Zhang     ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
1437*55de54fcSHong Zhang     /*ierr = VecCopy(ts->vec_sol,Yc[i]);CHKERRQ(ierr);*/
1438*55de54fcSHong Zhang     ierr = VecGetSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr);
1439*55de54fcSHong Zhang     ierr = VecGetSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
1440*55de54fcSHong Zhang     for (j=0; j<i; j++) w[j] = h*A[i*s+j];
1441*55de54fcSHong Zhang     ierr = VecMAXPY(Yslow,i,w,YdotRHS_slow);CHKERRQ(ierr);
1442*55de54fcSHong Zhang     ierr = VecMAXPY(Yfast,i,w,YdotRHS_fast);CHKERRQ(ierr);
1443*55de54fcSHong Zhang     ierr = VecRestoreSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr);
1444*55de54fcSHong Zhang     ierr = VecRestoreSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
1445*55de54fcSHong Zhang     ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
1446*55de54fcSHong Zhang     ierr = TSComputeRHSFunction(rk->subts_slow,t+h*c[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr);
1447*55de54fcSHong Zhang     ierr = TSComputeRHSFunction(rk->subts_fast,t+h*c[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr);
1448*55de54fcSHong Zhang   }
1449*55de54fcSHong Zhang   rk->slow = PETSC_TRUE;
1450*55de54fcSHong Zhang   ierr = TSEvaluateStep_MRKSPLIT(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
1451*55de54fcSHong Zhang   for (k=0; k<rk->dtratio; k++) {
1452*55de54fcSHong Zhang     /* propogate fast component using small time steps */
1453*55de54fcSHong Zhang     for (i=0; i<s; i++) {
1454*55de54fcSHong Zhang       rk->stage_time = t + h/rk->dtratio*c[i];
1455*55de54fcSHong Zhang       ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr);
1456*55de54fcSHong Zhang       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
1457*55de54fcSHong Zhang       /* stage value for fast components */
1458*55de54fcSHong Zhang       for (j=0; j<i; j++) w[j] = h/rk->dtratio*A[i*s+j];
1459*55de54fcSHong Zhang       ierr = VecGetSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
1460*55de54fcSHong Zhang       ierr = VecMAXPY(Yfast,i,w,YdotRHS_fast);CHKERRQ(ierr);
1461*55de54fcSHong Zhang       ierr = VecRestoreSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr);
1462*55de54fcSHong Zhang       /* stage value for slow components */
1463*55de54fcSHong Zhang       ierr = VecGetSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr);
1464*55de54fcSHong Zhang       ierr = TSInterpolate_MRKSPLIT(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Yslow);CHKERRQ(ierr);
1465*55de54fcSHong Zhang       ierr = VecRestoreSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr);
1466*55de54fcSHong Zhang       ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
1467*55de54fcSHong Zhang       /* compute the stage RHS for fast components */
1468*55de54fcSHong Zhang       ierr = TSComputeRHSFunction(rk->subts_fast,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr);
1469*55de54fcSHong Zhang     }
1470*55de54fcSHong Zhang     /* update the value of fast components when using fast time step */
1471*55de54fcSHong Zhang     rk->slow = PETSC_FALSE;
1472*55de54fcSHong Zhang     ierr = TSEvaluateStep_MRKSPLIT(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
1473*55de54fcSHong Zhang   }
1474*55de54fcSHong Zhang   for (i=0; i<s; i++) {
1475*55de54fcSHong Zhang     ierr = VecRestoreSubVector(YdotRHS[i],rk->is_slow,&YdotRHS_slow[i]);CHKERRQ(ierr);
1476*55de54fcSHong Zhang     ierr = VecRestoreSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr);
1477*55de54fcSHong Zhang   }
1478*55de54fcSHong Zhang   ts->ptime += ts->time_step;
1479*55de54fcSHong Zhang   ts->time_step = next_time_step;
1480*55de54fcSHong Zhang   rk->status = TS_STEP_COMPLETE;
1481*55de54fcSHong Zhang   PetscFunctionReturn(0);
1482*55de54fcSHong Zhang }
1483*55de54fcSHong Zhang 
1484*55de54fcSHong Zhang /*@C
1485*55de54fcSHong Zhang   TSRKSetMultirateType - Set the type of RK Multirate scheme
1486*55de54fcSHong Zhang 
1487*55de54fcSHong Zhang   Logically collective
1488*55de54fcSHong Zhang 
1489*55de54fcSHong Zhang   Input Parameter:
1490*55de54fcSHong Zhang +  ts - timestepping context
1491*55de54fcSHong Zhang -  mrktype - type of MRK-scheme
1492*55de54fcSHong Zhang 
1493*55de54fcSHong Zhang   Options Database:
1494*55de54fcSHong Zhang .   -ts_rk_multiarte_type - <none,nonsplit,split>
1495*55de54fcSHong Zhang 
1496*55de54fcSHong Zhang   Level: intermediate
1497*55de54fcSHong Zhang @*/
1498*55de54fcSHong Zhang PetscErrorCode TSRKSetMultirateType(TS ts, TSMRKType mrktype)
1499*55de54fcSHong Zhang {
1500*55de54fcSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
1501*55de54fcSHong Zhang   PetscErrorCode ierr;
1502*55de54fcSHong Zhang 
1503*55de54fcSHong Zhang   PetscFunctionBegin;
1504*55de54fcSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1505*55de54fcSHong Zhang   switch(mrktype){
1506*55de54fcSHong Zhang     case TSMRKNONE:
1507*55de54fcSHong Zhang       break;
1508*55de54fcSHong Zhang     case TSMRKNONSPLIT:
1509*55de54fcSHong Zhang       ts->ops->step           = TSStep_MRKNONSPLIT;
1510*55de54fcSHong Zhang       ts->ops->evaluatestep   = TSEvaluateStep_MRKNONSPLIT;
1511*55de54fcSHong Zhang       ts->ops->interpolate    = TSInterpolate_MRKNONSPLIT;
1512*55de54fcSHong Zhang       rk->dtratio             = 2;
1513*55de54fcSHong Zhang       ierr = TSRKSetType(ts,TSMRKDefault);CHKERRQ(ierr);
1514*55de54fcSHong Zhang       ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_MRKNONSPLIT_C",TSSetUp_MRKNONSPLIT);CHKERRQ(ierr);
1515*55de54fcSHong Zhang       ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_MRKNONSPLIT_C",TSReset_MRKNONSPLIT);CHKERRQ(ierr);
1516*55de54fcSHong Zhang       break;
1517*55de54fcSHong Zhang     case TSMRKSPLIT:
1518*55de54fcSHong Zhang       ts->ops->step           = TSStep_MRKSPLIT;
1519*55de54fcSHong Zhang       ts->ops->evaluatestep   = TSEvaluateStep_MRKSPLIT;
1520*55de54fcSHong Zhang       ts->ops->interpolate    = TSInterpolate_MRKSPLIT;
1521*55de54fcSHong Zhang       rk->dtratio             = 2;
1522*55de54fcSHong Zhang       ierr = TSRKSetType(ts,TSMRKDefault);CHKERRQ(ierr);
1523*55de54fcSHong Zhang       ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_MRKSPLIT_C",TSSetUp_MRKSPLIT);CHKERRQ(ierr);
1524*55de54fcSHong Zhang       ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_MRKSPLIT_C",TSReset_MRKSPLIT);CHKERRQ(ierr);
1525*55de54fcSHong Zhang       break;
1526*55de54fcSHong Zhang     default :
1527*55de54fcSHong Zhang       SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown type '%s'",mrktype);
1528*55de54fcSHong Zhang   }
1529*55de54fcSHong Zhang   PetscFunctionReturn(0);
1530*55de54fcSHong Zhang }
1531