xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision 2c9cb0628d9ce4739c99b849a0530650236a11be)
1e4dd521cSBarry Smith /*
2*2c9cb062Sluzhanghpp   Code for time stepping with the Runge-Kutta method(single-rate RK and multi-rate RK)
3f68a32c8SEmil Constantinescu 
4f68a32c8SEmil Constantinescu   Notes:
5*2c9cb062Sluzhanghpp   1) The general system is written as
6*2c9cb062Sluzhanghpp      Udot = F(t,U) for single-rate RK;
7*2c9cb062Sluzhanghpp   2) The general system is written as
8*2c9cb062Sluzhanghpp      Udot = F(t,U) for combined RHS multi-rate RK,
9*2c9cb062Sluzhanghpp      user should give the indexes for both slow and fast components;
10*2c9cb062Sluzhanghpp   3) The general system is written as
11*2c9cb062Sluzhanghpp      Usdot = Fs(t,Us,Uf)
12*2c9cb062Sluzhanghpp      Ufdot = Ff(t,Us,Uf) for partitioned RHS multi-rate RK,
13*2c9cb062Sluzhanghpp      user should partioned RHS by themselves and also provide the indexes for both slow and fast components.
14e4dd521cSBarry Smith */
15*2c9cb062Sluzhanghpp 
16af0996ceSBarry Smith #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
17f68a32c8SEmil Constantinescu #include <petscdm.h>
18f68a32c8SEmil Constantinescu 
19484bcdc7SDebojyoti Ghosh static TSRKType  TSRKDefault = TSRK3BS;
20*2c9cb062Sluzhanghpp static TSRKType  TSMRKDefault = TSRK2A;
21f68a32c8SEmil Constantinescu static PetscBool TSRKRegisterAllCalled;
22f68a32c8SEmil Constantinescu static PetscBool TSRKPackageInitialized;
23f68a32c8SEmil Constantinescu 
24f68a32c8SEmil Constantinescu typedef struct _RKTableau *RKTableau;
25f68a32c8SEmil Constantinescu struct _RKTableau {
26f68a32c8SEmil Constantinescu   char       *name;
27d760c35bSDebojyoti Ghosh   PetscInt   order;               /* Classical approximation order of the method i              */
28f68a32c8SEmil Constantinescu   PetscInt   s;                   /* Number of stages                                           */
293a8a9803SLisandro Dalcin   PetscInt   p;                   /* Interpolation order                                        */
30d760c35bSDebojyoti Ghosh   PetscBool  FSAL;                /* flag to indicate if tableau is FSAL                        */
31f68a32c8SEmil Constantinescu   PetscReal  *A,*b,*c;            /* Tableau                                                    */
32f68a32c8SEmil Constantinescu   PetscReal  *bembed;             /* Embedded formula of order one less (order-1)               */
33f68a32c8SEmil Constantinescu   PetscReal  *binterp;            /* Dense output formula                                       */
34f68a32c8SEmil Constantinescu   PetscReal  ccfl;                /* Placeholder for CFL coefficient relative to forward Euler  */
35f68a32c8SEmil Constantinescu };
36f68a32c8SEmil Constantinescu typedef struct _RKTableauLink *RKTableauLink;
37f68a32c8SEmil Constantinescu struct _RKTableauLink {
38f68a32c8SEmil Constantinescu   struct _RKTableau tab;
39f68a32c8SEmil Constantinescu   RKTableauLink     next;
40f68a32c8SEmil Constantinescu };
41f68a32c8SEmil Constantinescu static RKTableauLink RKTableauList;
42e4dd521cSBarry Smith 
43e4dd521cSBarry Smith typedef struct {
44f68a32c8SEmil Constantinescu   RKTableau    tableau;
45f68a32c8SEmil Constantinescu   Vec          *Y;               /* States computed during the step                                              */
46*2c9cb062Sluzhanghpp   Vec          *YdotRHS;         /* Function evaluations for the non-stiff part and contains all components      */
47*2c9cb062Sluzhanghpp   Vec          *YdotRHS_fast;    /* Function evaluations for the non-stiff part and contains fast components     */
48*2c9cb062Sluzhanghpp   Vec          *YdotRHS_slow;    /* Function evaluations for the non-stiff part and contains slow components     */
49ad8e2604SHong Zhang   Vec          *VecDeltaLam;     /* Increment of the adjoint sensitivity w.r.t IC at stage                       */
50ad8e2604SHong Zhang   Vec          *VecDeltaMu;      /* Increment of the adjoint sensitivity w.r.t P at stage                        */
510f7a1166SHong Zhang   Vec          VecCostIntegral0; /* backup for roll-backs due to events                                          */
52f68a32c8SEmil Constantinescu   PetscScalar  *work;            /* Scalar work                                                                  */
53*2c9cb062Sluzhanghpp   PetscInt     slow;             /* flag indicates call slow components solver (0) or fast components solver (1) */
54f68a32c8SEmil Constantinescu   PetscReal    stage_time;
55f68a32c8SEmil Constantinescu   TSStepStatus status;
560f7a1166SHong Zhang   PetscReal    ptime;
570f7a1166SHong Zhang   PetscReal    time_step;
58*2c9cb062Sluzhanghpp   PetscInt     dtratio;          /* ratio between slow time step size and fast step size                         */
595f70b478SJed Brown } TS_RK;
60e4dd521cSBarry Smith 
61*2c9cb062Sluzhanghpp 
62f68a32c8SEmil Constantinescu /*MC
63287c2655SBarry Smith      TSRK1FE - First order forward Euler scheme.
64262ff7bbSBarry Smith 
65f68a32c8SEmil Constantinescu      This method has one stage.
66f68a32c8SEmil Constantinescu 
67287c2655SBarry Smith      Options database:
689bd3e852SBarry Smith .     -ts_rk_type 1fe
69287c2655SBarry Smith 
70f68a32c8SEmil Constantinescu      Level: advanced
71f68a32c8SEmil Constantinescu 
72287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
73f68a32c8SEmil Constantinescu M*/
74f68a32c8SEmil Constantinescu /*MC
752109b73fSDebojyoti Ghosh      TSRK2A - Second order RK scheme.
76f68a32c8SEmil Constantinescu 
77f68a32c8SEmil Constantinescu      This method has two stages.
78f68a32c8SEmil Constantinescu 
79287c2655SBarry Smith      Options database:
809bd3e852SBarry Smith .     -ts_rk_type 2a
81287c2655SBarry Smith 
82f68a32c8SEmil Constantinescu      Level: advanced
83f68a32c8SEmil Constantinescu 
84287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
85f68a32c8SEmil Constantinescu M*/
86f68a32c8SEmil Constantinescu /*MC
87f68a32c8SEmil Constantinescu      TSRK3 - Third order RK scheme.
88f68a32c8SEmil Constantinescu 
89f68a32c8SEmil Constantinescu      This method has three stages.
90f68a32c8SEmil Constantinescu 
91287c2655SBarry Smith      Options database:
929bd3e852SBarry Smith .     -ts_rk_type 3
93287c2655SBarry Smith 
94f68a32c8SEmil Constantinescu      Level: advanced
95f68a32c8SEmil Constantinescu 
96287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
97f68a32c8SEmil Constantinescu M*/
98f68a32c8SEmil Constantinescu /*MC
992109b73fSDebojyoti Ghosh      TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method.
1002109b73fSDebojyoti Ghosh 
101d1905564SLisandro Dalcin      This method has four stages with the First Same As Last (FSAL) property.
1022109b73fSDebojyoti Ghosh 
103287c2655SBarry Smith      Options database:
1049bd3e852SBarry Smith .     -ts_rk_type 3bs
105287c2655SBarry Smith 
1062109b73fSDebojyoti Ghosh      Level: advanced
1072109b73fSDebojyoti Ghosh 
10898164e13SLisandro Dalcin      References: https://doi.org/10.1016/0893-9659(89)90079-7
109d1905564SLisandro Dalcin 
110287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
1112109b73fSDebojyoti Ghosh M*/
1122109b73fSDebojyoti Ghosh /*MC
113f68a32c8SEmil Constantinescu      TSRK4 - Fourth order RK scheme.
114f68a32c8SEmil Constantinescu 
1152e698f8bSDebojyoti Ghosh      This is the classical Runge-Kutta method with four stages.
116f68a32c8SEmil Constantinescu 
117287c2655SBarry Smith      Options database:
1189bd3e852SBarry Smith .     -ts_rk_type 4
119287c2655SBarry Smith 
120f68a32c8SEmil Constantinescu      Level: advanced
121f68a32c8SEmil Constantinescu 
122287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
123f68a32c8SEmil Constantinescu M*/
124f68a32c8SEmil Constantinescu /*MC
1252e698f8bSDebojyoti Ghosh      TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method.
126f68a32c8SEmil Constantinescu 
127f68a32c8SEmil Constantinescu      This method has six stages.
128f68a32c8SEmil Constantinescu 
129287c2655SBarry Smith      Options database:
1309bd3e852SBarry Smith .     -ts_rk_type 5f
131287c2655SBarry Smith 
132f68a32c8SEmil Constantinescu      Level: advanced
133f68a32c8SEmil Constantinescu 
134287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
135f68a32c8SEmil Constantinescu M*/
1362109b73fSDebojyoti Ghosh /*MC
1372e698f8bSDebojyoti Ghosh      TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method.
1382109b73fSDebojyoti Ghosh 
139d1905564SLisandro Dalcin      This method has seven stages with the First Same As Last (FSAL) property.
1402109b73fSDebojyoti Ghosh 
141287c2655SBarry Smith      Options database:
1429bd3e852SBarry Smith .     -ts_rk_type 5dp
143287c2655SBarry Smith 
1442109b73fSDebojyoti Ghosh      Level: advanced
1452109b73fSDebojyoti Ghosh 
14698164e13SLisandro Dalcin      References: https://doi.org/10.1016/0771-050X(80)90013-3
147d1905564SLisandro Dalcin 
148287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
1492109b73fSDebojyoti Ghosh M*/
15005e23783SLisandro Dalcin /*MC
15105e23783SLisandro Dalcin      TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method.
15205e23783SLisandro Dalcin 
15305e23783SLisandro Dalcin      This method has eight stages with the First Same As Last (FSAL) property.
15405e23783SLisandro Dalcin 
155287c2655SBarry Smith      Options database:
1569bd3e852SBarry Smith .     -ts_rk_type 5bs
157287c2655SBarry Smith 
15805e23783SLisandro Dalcin      Level: advanced
15905e23783SLisandro Dalcin 
16098164e13SLisandro Dalcin      References: https://doi.org/10.1016/0898-1221(96)00141-1
16105e23783SLisandro Dalcin 
162287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
16305e23783SLisandro Dalcin M*/
164262ff7bbSBarry Smith 
165f68a32c8SEmil Constantinescu /*@C
166f68a32c8SEmil Constantinescu   TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK
167262ff7bbSBarry Smith 
168f68a32c8SEmil Constantinescu   Not Collective, but should be called by all processes which will need the schemes to be registered
169262ff7bbSBarry Smith 
170f68a32c8SEmil Constantinescu   Level: advanced
171262ff7bbSBarry Smith 
172f68a32c8SEmil Constantinescu .keywords: TS, TSRK, register, all
173262ff7bbSBarry Smith 
174f68a32c8SEmil Constantinescu .seealso:  TSRKRegisterDestroy()
175262ff7bbSBarry Smith @*/
176f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void)
177262ff7bbSBarry Smith {
1784ac538c5SBarry Smith   PetscErrorCode ierr;
179262ff7bbSBarry Smith 
180262ff7bbSBarry Smith   PetscFunctionBegin;
181f68a32c8SEmil Constantinescu   if (TSRKRegisterAllCalled) PetscFunctionReturn(0);
182f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_TRUE;
183e4dd521cSBarry Smith 
1844e82626cSLisandro Dalcin #define RC PetscRealConstant
185e4dd521cSBarry Smith   {
186f68a32c8SEmil Constantinescu     const PetscReal
1874e82626cSLisandro Dalcin       A[1][1] = {{0}},
1884e82626cSLisandro Dalcin       b[1]    = {RC(1.0)};
1893a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
190e4dd521cSBarry Smith   }
191db046809SBarry Smith   {
192f68a32c8SEmil Constantinescu     const PetscReal
1934e82626cSLisandro Dalcin       A[2][2]   = {{0,0},
1944e82626cSLisandro Dalcin                    {RC(1.0),0}},
1954e82626cSLisandro Dalcin       b[2]      =  {RC(0.5),RC(0.5)},
1964e82626cSLisandro Dalcin       bembed[2] =  {RC(1.0),0};
1973a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
198db046809SBarry Smith   }
199f68a32c8SEmil Constantinescu   {
200f68a32c8SEmil Constantinescu     const PetscReal
20117f6384fSBarry Smith       A[3][3] = {{0,0,0},
2024e82626cSLisandro Dalcin                  {RC(2.0)/RC(3.0),0,0},
2034e82626cSLisandro Dalcin                  {RC(-1.0)/RC(3.0),RC(1.0),0}},
2044e82626cSLisandro Dalcin       b[3]    =  {RC(0.25),RC(0.5),RC(0.25)};
2053a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
206fdefd5e5SDebojyoti Ghosh   }
207fdefd5e5SDebojyoti Ghosh   {
208fdefd5e5SDebojyoti Ghosh     const PetscReal
20917f6384fSBarry Smith       A[4][4]   = {{0,0,0,0},
2104e82626cSLisandro Dalcin                    {RC(1.0)/RC(2.0),0,0,0},
2114e82626cSLisandro Dalcin                    {0,RC(3.0)/RC(4.0),0,0},
2124e82626cSLisandro Dalcin                    {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}},
2134e82626cSLisandro Dalcin       b[4]      =  {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0},
2144e82626cSLisandro 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)};
2153a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
216db046809SBarry Smith   }
217f68a32c8SEmil Constantinescu   {
218f68a32c8SEmil Constantinescu     const PetscReal
219f68a32c8SEmil Constantinescu       A[4][4] = {{0,0,0,0},
2204e82626cSLisandro Dalcin                  {RC(0.5),0,0,0},
2214e82626cSLisandro Dalcin                  {0,RC(0.5),0,0},
2224e82626cSLisandro Dalcin                  {0,0,RC(1.0),0}},
2234e82626cSLisandro 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)};
2243a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
225db046809SBarry Smith   }
226f68a32c8SEmil Constantinescu   {
227f68a32c8SEmil Constantinescu     const PetscReal
22817f6384fSBarry Smith       A[6][6]   = {{0,0,0,0,0,0},
2294e82626cSLisandro Dalcin                    {RC(0.25),0,0,0,0,0},
2304e82626cSLisandro Dalcin                    {RC(3.0)/RC(32.0),RC(9.0)/RC(32.0),0,0,0,0},
2314e82626cSLisandro Dalcin                    {RC(1932.0)/RC(2197.0),RC(-7200.0)/RC(2197.0),RC(7296.0)/RC(2197.0),0,0,0},
2324e82626cSLisandro Dalcin                    {RC(439.0)/RC(216.0),RC(-8.0),RC(3680.0)/RC(513.0),RC(-845.0)/RC(4104.0),0,0},
2334e82626cSLisandro 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}},
2344e82626cSLisandro 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)},
2354e82626cSLisandro 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};
2363a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
237fdefd5e5SDebojyoti Ghosh   }
238fdefd5e5SDebojyoti Ghosh   {
239fdefd5e5SDebojyoti Ghosh     const PetscReal
24017f6384fSBarry Smith       A[7][7]   = {{0,0,0,0,0,0,0},
2414e82626cSLisandro Dalcin                    {RC(1.0)/RC(5.0),0,0,0,0,0,0},
2424e82626cSLisandro Dalcin                    {RC(3.0)/RC(40.0),RC(9.0)/RC(40.0),0,0,0,0,0},
2434e82626cSLisandro Dalcin                    {RC(44.0)/RC(45.0),RC(-56.0)/RC(15.0),RC(32.0)/RC(9.0),0,0,0,0},
2444e82626cSLisandro 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},
2454e82626cSLisandro 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},
2464e82626cSLisandro 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}},
2474e82626cSLisandro 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},
2484e82626cSLisandro Dalcin       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)};
2493a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
250f68a32c8SEmil Constantinescu   }
25105e23783SLisandro Dalcin   {
25205e23783SLisandro Dalcin     const PetscReal
25317f6384fSBarry Smith       A[8][8]   = {{0,0,0,0,0,0,0,0},
2544e82626cSLisandro Dalcin                    {RC(1.0)/RC(6.0),0,0,0,0,0,0,0},
2554e82626cSLisandro Dalcin                    {RC(2.0)/RC(27.0),RC(4.0)/RC(27.0),0,0,0,0,0,0},
2564e82626cSLisandro Dalcin                    {RC(183.0)/RC(1372.0),RC(-162.0)/RC(343.0),RC(1053.0)/RC(1372.0),0,0,0,0,0},
2574e82626cSLisandro 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},
2584e82626cSLisandro 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},
2594e82626cSLisandro 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},
2604e82626cSLisandro 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}},
2614e82626cSLisandro 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},
2624e82626cSLisandro 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)};
26305e23783SLisandro Dalcin     ierr = TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
26405e23783SLisandro Dalcin   }
2654e82626cSLisandro Dalcin #undef RC
266db046809SBarry Smith   PetscFunctionReturn(0);
267db046809SBarry Smith }
268db046809SBarry Smith 
269f68a32c8SEmil Constantinescu /*@C
270f68a32c8SEmil Constantinescu    TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister().
271f68a32c8SEmil Constantinescu 
272f68a32c8SEmil Constantinescu    Not Collective
273f68a32c8SEmil Constantinescu 
274f68a32c8SEmil Constantinescu    Level: advanced
275f68a32c8SEmil Constantinescu 
276f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy
277f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll()
278f68a32c8SEmil Constantinescu @*/
279f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void)
280e4dd521cSBarry Smith {
281f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
282f68a32c8SEmil Constantinescu   RKTableauLink  link;
283f68a32c8SEmil Constantinescu 
284f68a32c8SEmil Constantinescu   PetscFunctionBegin;
285f68a32c8SEmil Constantinescu   while ((link = RKTableauList)) {
286f68a32c8SEmil Constantinescu     RKTableau t = &link->tab;
287f68a32c8SEmil Constantinescu     RKTableauList = link->next;
288f68a32c8SEmil Constantinescu     ierr = PetscFree3(t->A,t->b,t->c);  CHKERRQ(ierr);
289f68a32c8SEmil Constantinescu     ierr = PetscFree (t->bembed);       CHKERRQ(ierr);
290f68a32c8SEmil Constantinescu     ierr = PetscFree (t->binterp);      CHKERRQ(ierr);
291f68a32c8SEmil Constantinescu     ierr = PetscFree (t->name);         CHKERRQ(ierr);
292f68a32c8SEmil Constantinescu     ierr = PetscFree (link);            CHKERRQ(ierr);
293f68a32c8SEmil Constantinescu   }
294f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_FALSE;
295f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
296f68a32c8SEmil Constantinescu }
297f68a32c8SEmil Constantinescu 
298f68a32c8SEmil Constantinescu /*@C
299f68a32c8SEmil Constantinescu   TSRKInitializePackage - This function initializes everything in the TSRK package. It is called
3008a690491SBarry Smith   from TSInitializePackage().
301f68a32c8SEmil Constantinescu 
302f68a32c8SEmil Constantinescu   Level: developer
303f68a32c8SEmil Constantinescu 
304f68a32c8SEmil Constantinescu .keywords: TS, TSRK, initialize, package
305f68a32c8SEmil Constantinescu .seealso: PetscInitialize()
306f68a32c8SEmil Constantinescu @*/
307f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void)
308f68a32c8SEmil Constantinescu {
309186e87acSLisandro Dalcin   PetscErrorCode ierr;
310e4dd521cSBarry Smith 
311e4dd521cSBarry Smith   PetscFunctionBegin;
312f68a32c8SEmil Constantinescu   if (TSRKPackageInitialized) PetscFunctionReturn(0);
313f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_TRUE;
314f68a32c8SEmil Constantinescu   ierr = TSRKRegisterAll();CHKERRQ(ierr);
315f68a32c8SEmil Constantinescu   ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr);
316f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
317f68a32c8SEmil Constantinescu }
318186e87acSLisandro Dalcin 
319f68a32c8SEmil Constantinescu /*@C
320f68a32c8SEmil Constantinescu   TSRKFinalizePackage - This function destroys everything in the TSRK package. It is
321f68a32c8SEmil Constantinescu   called from PetscFinalize().
322186e87acSLisandro Dalcin 
323f68a32c8SEmil Constantinescu   Level: developer
324f68a32c8SEmil Constantinescu 
325f68a32c8SEmil Constantinescu .keywords: Petsc, destroy, package
326f68a32c8SEmil Constantinescu .seealso: PetscFinalize()
327f68a32c8SEmil Constantinescu @*/
328f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void)
329f68a32c8SEmil Constantinescu {
330f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
331f68a32c8SEmil Constantinescu 
332f68a32c8SEmil Constantinescu   PetscFunctionBegin;
333f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_FALSE;
334f68a32c8SEmil Constantinescu   ierr = TSRKRegisterDestroy();CHKERRQ(ierr);
335f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
336f68a32c8SEmil Constantinescu }
337f68a32c8SEmil Constantinescu 
338f68a32c8SEmil Constantinescu /*@C
339f68a32c8SEmil Constantinescu    TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
340f68a32c8SEmil Constantinescu 
341f68a32c8SEmil Constantinescu    Not Collective, but the same schemes should be registered on all processes on which they will be used
342f68a32c8SEmil Constantinescu 
343f68a32c8SEmil Constantinescu    Input Parameters:
344f68a32c8SEmil Constantinescu +  name - identifier for method
345f68a32c8SEmil Constantinescu .  order - approximation order of method
346f68a32c8SEmil Constantinescu .  s - number of stages, this is the dimension of the matrices below
347f68a32c8SEmil Constantinescu .  A - stage coefficients (dimension s*s, row-major)
348f68a32c8SEmil Constantinescu .  b - step completion table (dimension s; NULL to use last row of A)
349f68a32c8SEmil Constantinescu .  c - abscissa (dimension s; NULL to use row sums of A)
350f68a32c8SEmil Constantinescu .  bembed - completion table for embedded method (dimension s; NULL if not available)
3513a8a9803SLisandro Dalcin .  p - Order of the interpolation scheme, equal to the number of columns of binterp
3523a8a9803SLisandro Dalcin -  binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1)
353f68a32c8SEmil Constantinescu 
354f68a32c8SEmil Constantinescu    Notes:
355f68a32c8SEmil Constantinescu    Several RK methods are provided, this function is only needed to create new methods.
356f68a32c8SEmil Constantinescu 
357f68a32c8SEmil Constantinescu    Level: advanced
358f68a32c8SEmil Constantinescu 
359f68a32c8SEmil Constantinescu .keywords: TS, register
360f68a32c8SEmil Constantinescu 
361f68a32c8SEmil Constantinescu .seealso: TSRK
362f68a32c8SEmil Constantinescu @*/
363f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s,
364f68a32c8SEmil Constantinescu                             const PetscReal A[],const PetscReal b[],const PetscReal c[],
3653a8a9803SLisandro Dalcin                             const PetscReal bembed[],PetscInt p,const PetscReal binterp[])
366f68a32c8SEmil Constantinescu {
367f68a32c8SEmil Constantinescu   PetscErrorCode  ierr;
368f68a32c8SEmil Constantinescu   RKTableauLink   link;
369f68a32c8SEmil Constantinescu   RKTableau       t;
370f68a32c8SEmil Constantinescu   PetscInt        i,j;
371f68a32c8SEmil Constantinescu 
372f68a32c8SEmil Constantinescu   PetscFunctionBegin;
3733a8a9803SLisandro Dalcin   PetscValidCharPointer(name,1);
3743a8a9803SLisandro Dalcin   PetscValidRealPointer(A,4);
3753a8a9803SLisandro Dalcin   if (b) PetscValidRealPointer(b,5);
3763a8a9803SLisandro Dalcin   if (c) PetscValidRealPointer(c,6);
3773a8a9803SLisandro Dalcin   if (bembed) PetscValidRealPointer(bembed,7);
3783a8a9803SLisandro Dalcin   if (binterp || p > 1) PetscValidRealPointer(binterp,9);
3793a8a9803SLisandro Dalcin 
3801d36bdfdSBarry Smith   ierr = TSRKInitializePackage();CHKERRQ(ierr);
38195dccacaSBarry Smith   ierr = PetscNew(&link);CHKERRQ(ierr);
382f68a32c8SEmil Constantinescu   t = &link->tab;
3833a8a9803SLisandro Dalcin 
384f68a32c8SEmil Constantinescu   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
385f68a32c8SEmil Constantinescu   t->order = order;
386f68a32c8SEmil Constantinescu   t->s = s;
387dcca6d9dSJed Brown   ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr);
388f68a32c8SEmil Constantinescu   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
389f68a32c8SEmil Constantinescu   if (b)  { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); }
390f68a32c8SEmil Constantinescu   else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i];
391f68a32c8SEmil Constantinescu   if (c)  { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); }
392f68a32c8SEmil 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];
393d760c35bSDebojyoti Ghosh   t->FSAL = PETSC_TRUE;
394d760c35bSDebojyoti Ghosh   for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE;
3953a8a9803SLisandro Dalcin 
396f68a32c8SEmil Constantinescu   if (bembed) {
397785e854fSJed Brown     ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr);
398f68a32c8SEmil Constantinescu     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
399f68a32c8SEmil Constantinescu   }
400f68a32c8SEmil Constantinescu 
4013a8a9803SLisandro Dalcin   if (!binterp) { p = 1; binterp = t->b; }
4023a8a9803SLisandro Dalcin   t->p = p;
4033a8a9803SLisandro Dalcin   ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr);
4043a8a9803SLisandro Dalcin   ierr = PetscMemcpy(t->binterp,binterp,s*p*sizeof(binterp[0]));CHKERRQ(ierr);
4053a8a9803SLisandro Dalcin 
406f68a32c8SEmil Constantinescu   link->next = RKTableauList;
407f68a32c8SEmil Constantinescu   RKTableauList = link;
408f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
409f68a32c8SEmil Constantinescu }
410f68a32c8SEmil Constantinescu 
411e4dd521cSBarry Smith /*
412*2c9cb062Sluzhanghpp  This is for single-step RK method
413f68a32c8SEmil Constantinescu  The step completion formula is
414e4dd521cSBarry Smith 
415f68a32c8SEmil Constantinescu  x1 = x0 + h b^T YdotRHS
416f68a32c8SEmil Constantinescu 
417f68a32c8SEmil Constantinescu  This function can be called before or after ts->vec_sol has been updated.
418f68a32c8SEmil Constantinescu  Suppose we have a completion formula (b) and an embedded formula (be) of different order.
419f68a32c8SEmil Constantinescu  We can write
420f68a32c8SEmil Constantinescu 
421f68a32c8SEmil Constantinescu  x1e = x0 + h be^T YdotRHS
422f68a32c8SEmil Constantinescu      = x1 - h b^T YdotRHS + h be^T YdotRHS
423f68a32c8SEmil Constantinescu      = x1 + h (be - b)^T YdotRHS
424f68a32c8SEmil Constantinescu 
425f68a32c8SEmil Constantinescu  so we can evaluate the method with different order even after the step has been optimistically completed.
426e4dd521cSBarry Smith */
427f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done)
428f68a32c8SEmil Constantinescu {
429f68a32c8SEmil Constantinescu   TS_RK          *rk   = (TS_RK*)ts->data;
430f68a32c8SEmil Constantinescu   RKTableau      tab  = rk->tableau;
431f68a32c8SEmil Constantinescu   PetscScalar    *w    = rk->work;
432f68a32c8SEmil Constantinescu   PetscReal      h;
433f68a32c8SEmil Constantinescu   PetscInt       s    = tab->s,j;
434f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
435f68a32c8SEmil Constantinescu 
436f68a32c8SEmil Constantinescu   PetscFunctionBegin;
437f68a32c8SEmil Constantinescu   switch (rk->status) {
438f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
439f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
440f68a32c8SEmil Constantinescu     h = ts->time_step; break;
441f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
442be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
443f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
444f68a32c8SEmil Constantinescu   }
445f68a32c8SEmil Constantinescu   if (order == tab->order) {
446f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) {
447f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
448f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->b[j];
449f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
450f68a32c8SEmil Constantinescu     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
451f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
452f68a32c8SEmil Constantinescu   } else if (order == tab->order-1) {
453f68a32c8SEmil Constantinescu     if (!tab->bembed) goto unavailable;
454f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) { /*Complete with the embedded method (be)*/
455f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
456f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
457f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
458f68a32c8SEmil Constantinescu     } else {  /*Rollback and re-complete using (be-b) */
459f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
460f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
461f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
462f68a32c8SEmil Constantinescu     }
463f68a32c8SEmil Constantinescu     if (done) *done = PETSC_TRUE;
464f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
465f68a32c8SEmil Constantinescu   }
466f68a32c8SEmil Constantinescu unavailable:
467f68a32c8SEmil Constantinescu   if (done) *done = PETSC_FALSE;
468a7fac7c2SEmil 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);
469f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
470f68a32c8SEmil Constantinescu }
471f68a32c8SEmil Constantinescu 
4720f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts)
4730f7a1166SHong Zhang {
4740f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
4750f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
4760f7a1166SHong Zhang   const PetscInt  s = tab->s;
4770f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
4780f7a1166SHong Zhang   Vec             *Y = rk->Y;
4790f7a1166SHong Zhang   PetscInt        i;
4800f7a1166SHong Zhang   PetscErrorCode  ierr;
4810f7a1166SHong Zhang 
4820f7a1166SHong Zhang   PetscFunctionBegin;
4830f7a1166SHong Zhang   /* backup cost integral */
4840f7a1166SHong Zhang   ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr);
4850f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
4860f7a1166SHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
48751a7f651SHong Zhang     ierr = TSComputeCostIntegrand(ts,rk->ptime+rk->time_step*c[i],Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
4880f7a1166SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
4890f7a1166SHong Zhang   }
4900f7a1166SHong Zhang   PetscFunctionReturn(0);
4910f7a1166SHong Zhang }
4920f7a1166SHong Zhang 
4930f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts)
4940f7a1166SHong Zhang {
4950f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
4960f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
4970f7a1166SHong Zhang   const PetscInt  s = tab->s;
4980f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
4990f7a1166SHong Zhang   Vec             *Y = rk->Y;
5000f7a1166SHong Zhang   PetscInt        i;
5010f7a1166SHong Zhang   PetscErrorCode  ierr;
5020f7a1166SHong Zhang 
5030f7a1166SHong Zhang   PetscFunctionBegin;
5040f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
5050f7a1166SHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
50651a7f651SHong Zhang     ierr = TSComputeCostIntegrand(ts,ts->ptime+ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
5070f7a1166SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
5080f7a1166SHong Zhang   }
5090f7a1166SHong Zhang   PetscFunctionReturn(0);
5100f7a1166SHong Zhang }
5110f7a1166SHong Zhang 
512fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts)
513fecfb714SLisandro Dalcin {
514fecfb714SLisandro Dalcin   TS_RK           *rk = (TS_RK*)ts->data;
515fecfb714SLisandro Dalcin   RKTableau       tab = rk->tableau;
516fecfb714SLisandro Dalcin   const PetscInt  s  = tab->s;
517fecfb714SLisandro Dalcin   const PetscReal *b = tab->b;
518fecfb714SLisandro Dalcin   PetscScalar     *w = rk->work;
519fecfb714SLisandro Dalcin   Vec             *YdotRHS = rk->YdotRHS;
520fecfb714SLisandro Dalcin   PetscInt        j;
521fecfb714SLisandro Dalcin   PetscReal       h;
522fecfb714SLisandro Dalcin   PetscErrorCode  ierr;
523fecfb714SLisandro Dalcin 
524fecfb714SLisandro Dalcin   PetscFunctionBegin;
525fecfb714SLisandro Dalcin   switch (rk->status) {
526fecfb714SLisandro Dalcin   case TS_STEP_INCOMPLETE:
527fecfb714SLisandro Dalcin   case TS_STEP_PENDING:
528fecfb714SLisandro Dalcin     h = ts->time_step; break;
529fecfb714SLisandro Dalcin   case TS_STEP_COMPLETE:
530fecfb714SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
531fecfb714SLisandro Dalcin   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
532fecfb714SLisandro Dalcin   }
533fecfb714SLisandro Dalcin   for (j=0; j<s; j++) w[j] = -h*b[j];
534fecfb714SLisandro Dalcin   ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr);
535fecfb714SLisandro Dalcin   PetscFunctionReturn(0);
536fecfb714SLisandro Dalcin }
537fecfb714SLisandro Dalcin 
538f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts)
539f68a32c8SEmil Constantinescu {
540f68a32c8SEmil Constantinescu   TS_RK            *rk  = (TS_RK*)ts->data;
541f68a32c8SEmil Constantinescu   RKTableau        tab  = rk->tableau;
542f68a32c8SEmil Constantinescu   const PetscInt   s = tab->s;
543fecfb714SLisandro Dalcin   const PetscReal  *A = tab->A,*c = tab->c;
544f68a32c8SEmil Constantinescu   PetscScalar      *w = rk->work;
545f68a32c8SEmil Constantinescu   Vec              *Y = rk->Y,*YdotRHS = rk->YdotRHS;
546d1905564SLisandro Dalcin   PetscBool        FSAL = tab->FSAL;
547f68a32c8SEmil Constantinescu   TSAdapt          adapt;
548fecfb714SLisandro Dalcin   PetscInt         i,j;
549be5899b3SLisandro Dalcin   PetscInt         rejections = 0;
550be5899b3SLisandro Dalcin   PetscBool        stageok,accept = PETSC_TRUE;
551be5899b3SLisandro Dalcin   PetscReal        next_time_step = ts->time_step;
552f68a32c8SEmil Constantinescu   PetscErrorCode   ierr;
553f68a32c8SEmil Constantinescu 
554f68a32c8SEmil Constantinescu   PetscFunctionBegin;
555d1905564SLisandro Dalcin   if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE;
556d1905564SLisandro Dalcin   if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); }
557d1905564SLisandro Dalcin 
558f68a32c8SEmil Constantinescu   rk->status = TS_STEP_INCOMPLETE;
559be5899b3SLisandro Dalcin   while (!ts->reason && rk->status != TS_STEP_COMPLETE) {
560be5899b3SLisandro Dalcin     PetscReal t = ts->ptime;
561f68a32c8SEmil Constantinescu     PetscReal h = ts->time_step;
562f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
5639be3e283SDebojyoti Ghosh       rk->stage_time = t + h*c[i];
5649be3e283SDebojyoti Ghosh       ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr);
565f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
566f68a32c8SEmil Constantinescu       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
567f68a32c8SEmil Constantinescu       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
5689be3e283SDebojyoti Ghosh       ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
569f68a32c8SEmil Constantinescu       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
570be5899b3SLisandro Dalcin       ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr);
571fecfb714SLisandro Dalcin       if (!stageok) goto reject_step;
5728f738a7cSLisandro Dalcin       if (FSAL && !i) continue;
573f68a32c8SEmil Constantinescu       ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
574f68a32c8SEmil Constantinescu     }
575be5899b3SLisandro Dalcin 
576be5899b3SLisandro Dalcin     rk->status = TS_STEP_INCOMPLETE;
577f68a32c8SEmil Constantinescu     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
578f68a32c8SEmil Constantinescu     rk->status = TS_STEP_PENDING;
579f68a32c8SEmil Constantinescu     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
580f68a32c8SEmil Constantinescu     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
5811917a363SLisandro Dalcin     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr);
582fecfb714SLisandro Dalcin     ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
583be5899b3SLisandro Dalcin     rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
584be5899b3SLisandro Dalcin     if (!accept) {  /*Roll back the current step */
585fecfb714SLisandro Dalcin       ierr = TSRollBack_RK(ts);CHKERRQ(ierr);
586be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
587be5899b3SLisandro Dalcin       goto reject_step;
588be5899b3SLisandro Dalcin     }
589be5899b3SLisandro Dalcin 
5900f7a1166SHong Zhang     if (ts->costintegralfwd) {  /*Save the info for the later use in cost integral evaluation*/
5910f7a1166SHong Zhang       rk->ptime     = ts->ptime;
5920f7a1166SHong Zhang       rk->time_step = ts->time_step;
593493ed95fSHong Zhang     }
594be5899b3SLisandro Dalcin 
595f68a32c8SEmil Constantinescu     ts->ptime += ts->time_step;
596f68a32c8SEmil Constantinescu     ts->time_step = next_time_step;
597f68a32c8SEmil Constantinescu     break;
598be5899b3SLisandro Dalcin 
599be5899b3SLisandro Dalcin     reject_step:
600fecfb714SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
601be5899b3SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
602be5899b3SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
603be5899b3SLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
604f68a32c8SEmil Constantinescu     }
605f68a32c8SEmil Constantinescu   }
606f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
607e4dd521cSBarry Smith }
608e4dd521cSBarry Smith 
609*2c9cb062Sluzhanghpp static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X)
610*2c9cb062Sluzhanghpp {
611*2c9cb062Sluzhanghpp   TS_RK            *rk = (TS_RK*)ts->data;
612*2c9cb062Sluzhanghpp   PetscInt         s  = rk->tableau->s,p = rk->tableau->p,i,j;
613*2c9cb062Sluzhanghpp   PetscReal        h;
614*2c9cb062Sluzhanghpp   PetscReal        tt,t;
615*2c9cb062Sluzhanghpp   PetscScalar      *b;
616*2c9cb062Sluzhanghpp   const PetscReal  *B = rk->tableau->binterp;
617*2c9cb062Sluzhanghpp   PetscErrorCode   ierr;
618*2c9cb062Sluzhanghpp 
619*2c9cb062Sluzhanghpp   PetscFunctionBegin;
620*2c9cb062Sluzhanghpp   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
621*2c9cb062Sluzhanghpp 
622*2c9cb062Sluzhanghpp   switch (rk->status) {
623*2c9cb062Sluzhanghpp     case TS_STEP_INCOMPLETE:
624*2c9cb062Sluzhanghpp     case TS_STEP_PENDING:
625*2c9cb062Sluzhanghpp       h = ts->time_step;
626*2c9cb062Sluzhanghpp       t = (itime - ts->ptime)/h;
627*2c9cb062Sluzhanghpp       break;
628*2c9cb062Sluzhanghpp     case TS_STEP_COMPLETE:
629*2c9cb062Sluzhanghpp       h = ts->ptime - ts->ptime_prev;
630*2c9cb062Sluzhanghpp       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
631*2c9cb062Sluzhanghpp       break;
632*2c9cb062Sluzhanghpp     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
633*2c9cb062Sluzhanghpp   }
634*2c9cb062Sluzhanghpp   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
635*2c9cb062Sluzhanghpp   for (i=0; i<s; i++) b[i] = 0;
636*2c9cb062Sluzhanghpp   for (j=0,tt=t; j<p; j++,tt*=t) {
637*2c9cb062Sluzhanghpp     for (i=0; i<s; i++) {
638*2c9cb062Sluzhanghpp       b[i]  += h * B[i*p+j] * tt;
639*2c9cb062Sluzhanghpp     }
640*2c9cb062Sluzhanghpp   }
641*2c9cb062Sluzhanghpp   ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr);
642*2c9cb062Sluzhanghpp   ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr);
643*2c9cb062Sluzhanghpp   ierr = PetscFree(b);CHKERRQ(ierr);
644*2c9cb062Sluzhanghpp   PetscFunctionReturn(0);
645*2c9cb062Sluzhanghpp }
646*2c9cb062Sluzhanghpp 
647*2c9cb062Sluzhanghpp /*
648*2c9cb062Sluzhanghpp   This is interpolate formula for partitioned RHS multirate RK method
649*2c9cb062Sluzhanghpp  */
650*2c9cb062Sluzhanghpp 
651*2c9cb062Sluzhanghpp static PetscErrorCode TSInterpolate_PARTITIONEDMRK(TS ts,PetscReal itime,Vec X)
652*2c9cb062Sluzhanghpp {
653*2c9cb062Sluzhanghpp   TS_RK            *rk = (TS_RK*)ts->data;
654*2c9cb062Sluzhanghpp   PetscInt         s  = rk->tableau->s,p = rk->tableau->p,i,j;
655*2c9cb062Sluzhanghpp   Vec              Yslow;    /* vector holds the slow components in Y[0] */
656*2c9cb062Sluzhanghpp   PetscReal        h;
657*2c9cb062Sluzhanghpp   PetscReal        tt,t;
658*2c9cb062Sluzhanghpp   PetscScalar      *b;
659*2c9cb062Sluzhanghpp   const PetscReal  *B = rk->tableau->binterp;
660*2c9cb062Sluzhanghpp   PetscErrorCode   ierr;
661*2c9cb062Sluzhanghpp 
662*2c9cb062Sluzhanghpp   PetscFunctionBegin;
663*2c9cb062Sluzhanghpp   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
664*2c9cb062Sluzhanghpp 
665*2c9cb062Sluzhanghpp   switch (rk->status) {
666*2c9cb062Sluzhanghpp     case TS_STEP_INCOMPLETE:
667*2c9cb062Sluzhanghpp     case TS_STEP_PENDING:
668*2c9cb062Sluzhanghpp       h = ts->time_step;
669*2c9cb062Sluzhanghpp       t = (itime - ts->ptime)/h;
670*2c9cb062Sluzhanghpp       break;
671*2c9cb062Sluzhanghpp     case TS_STEP_COMPLETE:
672*2c9cb062Sluzhanghpp       h = ts->ptime - ts->ptime_prev;
673*2c9cb062Sluzhanghpp       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
674*2c9cb062Sluzhanghpp       break;
675*2c9cb062Sluzhanghpp     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
676*2c9cb062Sluzhanghpp   }
677*2c9cb062Sluzhanghpp   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
678*2c9cb062Sluzhanghpp   for (i=0; i<s; i++) b[i] = 0;
679*2c9cb062Sluzhanghpp   for (j=0,tt=t; j<p; j++,tt*=t) {
680*2c9cb062Sluzhanghpp     for (i=0; i<s; i++) {
681*2c9cb062Sluzhanghpp       b[i]  += h * B[i*p+j] * tt;
682*2c9cb062Sluzhanghpp     }
683*2c9cb062Sluzhanghpp   }
684*2c9cb062Sluzhanghpp   ierr = VecGetSubVector(rk->Y[0],ts->iss,&Yslow);CHKERRQ(ierr);
685*2c9cb062Sluzhanghpp   ierr = VecCopy(Yslow,X);CHKERRQ(ierr);
686*2c9cb062Sluzhanghpp   ierr = VecMAXPY(X,s,b,rk->YdotRHS_slow);CHKERRQ(ierr);
687*2c9cb062Sluzhanghpp   ierr = VecRestoreSubVector(rk->Y[0],ts->iss,&Yslow);CHKERRQ(ierr);
688*2c9cb062Sluzhanghpp   ierr = PetscFree(b);CHKERRQ(ierr);
689*2c9cb062Sluzhanghpp   PetscFunctionReturn(0);
690*2c9cb062Sluzhanghpp }
691*2c9cb062Sluzhanghpp 
692*2c9cb062Sluzhanghpp /*
693*2c9cb062Sluzhanghpp  This is for combined RHS multirate RK method
694*2c9cb062Sluzhanghpp  The step completion formula is
695*2c9cb062Sluzhanghpp 
696*2c9cb062Sluzhanghpp  x1 = x0 + h b^T YdotRHS
697*2c9cb062Sluzhanghpp 
698*2c9cb062Sluzhanghpp */
699*2c9cb062Sluzhanghpp static PetscErrorCode TSEvaluateStep_RKMULTIRATE(TS ts,PetscInt order,Vec X,PetscBool *done)
700*2c9cb062Sluzhanghpp {
701*2c9cb062Sluzhanghpp   TS_RK           *rk = (TS_RK*)ts->data;
702*2c9cb062Sluzhanghpp   RKTableau       tab  = rk->tableau;
703*2c9cb062Sluzhanghpp   Vec             Xtemp;                          /* temporary work vector for X                                   */
704*2c9cb062Sluzhanghpp   PetscScalar     *w = rk->work;
705*2c9cb062Sluzhanghpp   PetscScalar     *x_ptr,*xtemp_ptr;              /* location to put the pointer to X and Xtemp respectively       */
706*2c9cb062Sluzhanghpp   PetscReal       h = ts->time_step;
707*2c9cb062Sluzhanghpp   PetscInt        s = tab->s,j;
708*2c9cb062Sluzhanghpp   PetscInt        len_slow,len_fast;              /* the number of slow components and fast components respectively */
709*2c9cb062Sluzhanghpp   const PetscInt  *is_slow,*is_fast;              /* the index for slow components and fast components respectively */
710*2c9cb062Sluzhanghpp   PetscErrorCode  ierr;
711*2c9cb062Sluzhanghpp 
712*2c9cb062Sluzhanghpp   PetscFunctionBegin;
713*2c9cb062Sluzhanghpp   ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
714*2c9cb062Sluzhanghpp   ierr = VecDuplicate(ts->vec_sol,&Xtemp);CHKERRQ(ierr);
715*2c9cb062Sluzhanghpp   ierr = VecCopy(ts->vec_sol,Xtemp);CHKERRQ(ierr);
716*2c9cb062Sluzhanghpp   if (!rk->slow){
717*2c9cb062Sluzhanghpp     for (j=0; j<s; j++) w[j] = h*tab->b[j];
718*2c9cb062Sluzhanghpp     /* update the value of slow components, and discard the updating value of fast components */
719*2c9cb062Sluzhanghpp     ierr = VecMAXPY(Xtemp,s,w,rk->YdotRHS);CHKERRQ(ierr);
720*2c9cb062Sluzhanghpp     ierr = VecGetArray(X,&x_ptr);CHKERRQ(ierr);
721*2c9cb062Sluzhanghpp     ierr = VecGetArray(Xtemp,&xtemp_ptr);CHKERRQ(ierr);
722*2c9cb062Sluzhanghpp     ierr = ISGetSize(ts->iss,&len_slow);CHKERRQ(ierr);
723*2c9cb062Sluzhanghpp     ierr = ISGetIndices(ts->iss,&is_slow);CHKERRQ(ierr);
724*2c9cb062Sluzhanghpp     ierr = ISRestoreIndices(ts->iss,&is_slow);CHKERRQ(ierr);
725*2c9cb062Sluzhanghpp     for (j=0; j<len_slow; j++) x_ptr[is_slow[j]] = xtemp_ptr[is_slow[j]];
726*2c9cb062Sluzhanghpp     ierr = VecRestoreArray(X,&x_ptr);CHKERRQ(ierr);
727*2c9cb062Sluzhanghpp     ierr = VecRestoreArray(Xtemp,&xtemp_ptr);CHKERRQ(ierr);
728*2c9cb062Sluzhanghpp   } else {
729*2c9cb062Sluzhanghpp     /* update the value of fast components, and discard the updating value of slow components */
730*2c9cb062Sluzhanghpp     for (j=0; j<s; j++) w[j] = h/rk->dtratio*tab->b[j];
731*2c9cb062Sluzhanghpp     ierr = VecMAXPY(Xtemp,s,w,rk->YdotRHS);CHKERRQ(ierr);
732*2c9cb062Sluzhanghpp     ierr = VecGetArray(X,&x_ptr);CHKERRQ(ierr);
733*2c9cb062Sluzhanghpp     ierr = VecGetArray(Xtemp,&xtemp_ptr);CHKERRQ(ierr);
734*2c9cb062Sluzhanghpp     ierr = ISGetSize(ts->isf,&len_fast);CHKERRQ(ierr);
735*2c9cb062Sluzhanghpp     ierr = ISGetIndices(ts->isf,&is_fast);CHKERRQ(ierr);
736*2c9cb062Sluzhanghpp     ierr = ISRestoreIndices(ts->isf,&is_fast);CHKERRQ(ierr);
737*2c9cb062Sluzhanghpp     for (j=0; j<len_fast; j++) x_ptr[is_fast[j]] = xtemp_ptr[is_fast[j]];
738*2c9cb062Sluzhanghpp     ierr = VecRestoreArray(X,&x_ptr);CHKERRQ(ierr);
739*2c9cb062Sluzhanghpp     ierr = VecRestoreArray(Xtemp,&xtemp_ptr);CHKERRQ(ierr);
740*2c9cb062Sluzhanghpp   }
741*2c9cb062Sluzhanghpp   /* free temporary vector space */
742*2c9cb062Sluzhanghpp   ierr = VecDestroy(&Xtemp);CHKERRQ(ierr);
743*2c9cb062Sluzhanghpp   PetscFunctionReturn(0);
744*2c9cb062Sluzhanghpp }
745*2c9cb062Sluzhanghpp 
746*2c9cb062Sluzhanghpp static PetscErrorCode TSStep_RKMULTIRATE(TS ts)
747*2c9cb062Sluzhanghpp {
748*2c9cb062Sluzhanghpp   TS_RK             *rk = (TS_RK*)ts->data;
749*2c9cb062Sluzhanghpp   RKTableau         tab  = rk->tableau;
750*2c9cb062Sluzhanghpp   Vec               *Y = rk->Y,*YdotRHS = rk->YdotRHS;
751*2c9cb062Sluzhanghpp   Vec               stage_fast,stage_slow;                              /* vectors store the stage value of fast and slow components respectively           */
752*2c9cb062Sluzhanghpp   Vec               presol_fast,newslo_fast;                            /* vectors store the previous and new time solution of fast components respectively */
753*2c9cb062Sluzhanghpp   Vec               YdotRHS_prev,prev_sol;                              /* vectors store the previous YdotRHS and solution respectively                     */
754*2c9cb062Sluzhanghpp   const PetscInt    s = tab->s,*is_slow,*is_fast;                       /* is_slow: index of slow components; is_fast: index of fast components             */
755*2c9cb062Sluzhanghpp   const PetscReal   *A = tab->A,*c = tab->c;
756*2c9cb062Sluzhanghpp   PetscScalar       *w = rk->work;
757*2c9cb062Sluzhanghpp   PetscScalar       *y_ptr,*faststage_ptr,*slowstage_ptr;               /* location to put the pointer to Y, stage_fast and stage_slow respectively         */
758*2c9cb062Sluzhanghpp   PetscScalar       *ydotrhsprev_ptr,*ydotrhs_ptr;                      /* location to put the pointer to YdotRHS_prev and YdotRHS respectively             */
759*2c9cb062Sluzhanghpp   PetscInt          i,j,k,len_slow,len_fast;                            /* len_slow: the number of slow comonents; len_fast: the number of fast components  */
760*2c9cb062Sluzhanghpp   PetscReal         next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step;
761*2c9cb062Sluzhanghpp   PetscErrorCode    ierr;
762*2c9cb062Sluzhanghpp 
763*2c9cb062Sluzhanghpp   PetscFunctionBegin;
764*2c9cb062Sluzhanghpp   rk->status = TS_STEP_INCOMPLETE;
765*2c9cb062Sluzhanghpp   ierr = VecDuplicate(ts->vec_sol,&stage_fast);CHKERRQ(ierr);
766*2c9cb062Sluzhanghpp   ierr = VecDuplicate(ts->vec_sol,&stage_slow);CHKERRQ(ierr);
767*2c9cb062Sluzhanghpp   ierr = VecDuplicate(ts->vec_sol,&YdotRHS_prev);CHKERRQ(ierr);
768*2c9cb062Sluzhanghpp   ierr = VecDuplicate(ts->vec_sol,&prev_sol);CHKERRQ(ierr);
769*2c9cb062Sluzhanghpp   ierr = ISGetSize(ts->iss,&len_slow);CHKERRQ(ierr);
770*2c9cb062Sluzhanghpp   ierr = ISGetSize(ts->isf,&len_fast);CHKERRQ(ierr);
771*2c9cb062Sluzhanghpp   ierr = ISGetIndices(ts->iss,&is_slow);CHKERRQ(ierr);
772*2c9cb062Sluzhanghpp   ierr = ISGetIndices(ts->isf,&is_fast);CHKERRQ(ierr);
773*2c9cb062Sluzhanghpp   ierr = ISRestoreIndices(ts->isf,&is_fast);CHKERRQ(ierr);
774*2c9cb062Sluzhanghpp   ierr = ISRestoreIndices(ts->iss,&is_slow);CHKERRQ(ierr);
775*2c9cb062Sluzhanghpp   for (i=0; i<s; i++) {
776*2c9cb062Sluzhanghpp     rk->stage_time = t + h*c[i];
777*2c9cb062Sluzhanghpp     ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr);
778*2c9cb062Sluzhanghpp     ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
779*2c9cb062Sluzhanghpp     for (j=0; j<i; j++) w[j] = h*A[i*s+j];
780*2c9cb062Sluzhanghpp     ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
781*2c9cb062Sluzhanghpp     ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
782*2c9cb062Sluzhanghpp     /* compute the stage RHS */
783*2c9cb062Sluzhanghpp     ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
784*2c9cb062Sluzhanghpp   }
785*2c9cb062Sluzhanghpp   ierr = VecCopy(ts->vec_sol,prev_sol);CHKERRQ(ierr);
786*2c9cb062Sluzhanghpp   rk->slow = 0;
787*2c9cb062Sluzhanghpp   ierr = TSEvaluateStep_RKMULTIRATE(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
788*2c9cb062Sluzhanghpp   for (k=0; k<rk->dtratio; k++){
789*2c9cb062Sluzhanghpp     for (i=0; i<s; i++) {
790*2c9cb062Sluzhanghpp       rk->stage_time = t + h/rk->dtratio*c[i];
791*2c9cb062Sluzhanghpp       ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr);
792*2c9cb062Sluzhanghpp       ierr = VecCopy(prev_sol,Y[i]);CHKERRQ(ierr);
793*2c9cb062Sluzhanghpp       ierr = VecCopy(prev_sol,stage_slow);CHKERRQ(ierr);
794*2c9cb062Sluzhanghpp       ierr = VecCopy(prev_sol,stage_fast);CHKERRQ(ierr);
795*2c9cb062Sluzhanghpp       /* guarantee the stage RHS for slow components do not change while updating the stage RHS for fast components */
796*2c9cb062Sluzhanghpp       ierr = VecCopy(YdotRHS[i],YdotRHS_prev);CHKERRQ(ierr);
797*2c9cb062Sluzhanghpp       /* update the fast components stage value */
798*2c9cb062Sluzhanghpp       for (j=0; j<i; j++) w[j] = h/rk->dtratio*A[i*s+j];
799*2c9cb062Sluzhanghpp       ierr = VecMAXPY(stage_fast,i,w,YdotRHS);CHKERRQ(ierr);
800*2c9cb062Sluzhanghpp       /* update the slow components stage value */
801*2c9cb062Sluzhanghpp       ierr = VecCopy(prev_sol,Y[0]);CHKERRQ(ierr);
802*2c9cb062Sluzhanghpp       ierr = TSInterpolate_RK(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],stage_slow);CHKERRQ(ierr);
803*2c9cb062Sluzhanghpp       /* combine the update fast components stage value and slow components stage value to Y[i] */
804*2c9cb062Sluzhanghpp       ierr = VecGetArray(Y[i],&y_ptr);CHKERRQ(ierr);
805*2c9cb062Sluzhanghpp       ierr = VecGetArray(stage_slow,&slowstage_ptr);CHKERRQ(ierr);
806*2c9cb062Sluzhanghpp       ierr = VecGetArray(stage_fast,&faststage_ptr);CHKERRQ(ierr);
807*2c9cb062Sluzhanghpp       for (j=0; j<len_slow; j++) y_ptr[is_slow[j]] = slowstage_ptr[is_slow[j]];
808*2c9cb062Sluzhanghpp       for (j=0; j<len_fast; j++) y_ptr[is_fast[j]] = faststage_ptr[is_fast[j]];
809*2c9cb062Sluzhanghpp       ierr = VecRestoreArray(Y[i],&y_ptr);CHKERRQ(ierr);
810*2c9cb062Sluzhanghpp       ierr = VecRestoreArray(stage_slow,&slowstage_ptr);CHKERRQ(ierr);
811*2c9cb062Sluzhanghpp       ierr = VecRestoreArray(stage_fast,&faststage_ptr);CHKERRQ(ierr);
812*2c9cb062Sluzhanghpp       ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
813*2c9cb062Sluzhanghpp       /* compute the stage RHS */
814*2c9cb062Sluzhanghpp       ierr = TSComputeRHSFunction(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
815*2c9cb062Sluzhanghpp       /* re-assign the value of slow components in YdotRHS[i] to the calculated value before working on smaller time step-size */
816*2c9cb062Sluzhanghpp       ierr = VecGetArray(YdotRHS[i],&ydotrhs_ptr);CHKERRQ(ierr);
817*2c9cb062Sluzhanghpp       ierr = VecGetArray(YdotRHS_prev,&ydotrhsprev_ptr);CHKERRQ(ierr);
818*2c9cb062Sluzhanghpp       for (j=0; j<len_slow; j++) ydotrhs_ptr[is_slow[j]] = ydotrhsprev_ptr[is_slow[j]];
819*2c9cb062Sluzhanghpp       ierr = VecRestoreArray(YdotRHS[i],&ydotrhs_ptr);CHKERRQ(ierr);
820*2c9cb062Sluzhanghpp       ierr = VecRestoreArray(YdotRHS_prev,&ydotrhsprev_ptr);CHKERRQ(ierr);
821*2c9cb062Sluzhanghpp     }
822*2c9cb062Sluzhanghpp     rk->slow = 1;
823*2c9cb062Sluzhanghpp     ierr = TSEvaluateStep_RKMULTIRATE(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
824*2c9cb062Sluzhanghpp     /* update the intial value for fast components */
825*2c9cb062Sluzhanghpp     ierr = VecGetSubVector(prev_sol,ts->isf,&presol_fast);CHKERRQ(ierr);
826*2c9cb062Sluzhanghpp     ierr = VecGetSubVector(ts->vec_sol,ts->isf,&newslo_fast);CHKERRQ(ierr);
827*2c9cb062Sluzhanghpp     ierr = VecCopy(newslo_fast,presol_fast);CHKERRQ(ierr);
828*2c9cb062Sluzhanghpp     ierr = VecRestoreSubVector(prev_sol,ts->isf,&presol_fast);CHKERRQ(ierr);
829*2c9cb062Sluzhanghpp     ierr = VecRestoreSubVector(ts->vec_sol,ts->isf,&newslo_fast);CHKERRQ(ierr);
830*2c9cb062Sluzhanghpp   }
831*2c9cb062Sluzhanghpp   ts->ptime += ts->time_step;
832*2c9cb062Sluzhanghpp   ts->time_step = next_time_step;
833*2c9cb062Sluzhanghpp   rk->slow = 0;
834*2c9cb062Sluzhanghpp   rk->status = TS_STEP_COMPLETE;
835*2c9cb062Sluzhanghpp   /* free memory of work vectors */
836*2c9cb062Sluzhanghpp   ierr = VecDestroy(&stage_fast);CHKERRQ(ierr);
837*2c9cb062Sluzhanghpp   ierr = VecDestroy(&stage_slow);CHKERRQ(ierr);
838*2c9cb062Sluzhanghpp   ierr = VecDestroy(&YdotRHS_prev);CHKERRQ(ierr);
839*2c9cb062Sluzhanghpp   ierr = VecDestroy(&prev_sol);CHKERRQ(ierr);
840*2c9cb062Sluzhanghpp   PetscFunctionReturn(0);
841*2c9cb062Sluzhanghpp }
842*2c9cb062Sluzhanghpp 
843*2c9cb062Sluzhanghpp /*
844*2c9cb062Sluzhanghpp  This is for partitioned RHS multirate RK method
845*2c9cb062Sluzhanghpp  The step completion formula is
846*2c9cb062Sluzhanghpp 
847*2c9cb062Sluzhanghpp  x1 = x0 + h b^T YdotRHS
848*2c9cb062Sluzhanghpp 
849*2c9cb062Sluzhanghpp */
850*2c9cb062Sluzhanghpp static PetscErrorCode TSEvaluateStep_RKPMULTIRATE(TS ts,PetscInt order,Vec X,PetscBool *done)
851*2c9cb062Sluzhanghpp {
852*2c9cb062Sluzhanghpp   TS_RK           *rk = (TS_RK*)ts->data;
853*2c9cb062Sluzhanghpp   RKTableau       tab  = rk->tableau;
854*2c9cb062Sluzhanghpp   Vec             Xslow,Xfast;                  /* subvectors of X which store slow components and fast components respectively */
855*2c9cb062Sluzhanghpp   PetscScalar     *w = rk->work;
856*2c9cb062Sluzhanghpp   PetscReal       h = ts->time_step;
857*2c9cb062Sluzhanghpp   PetscInt        s = tab->s,j;
858*2c9cb062Sluzhanghpp   PetscErrorCode  ierr;
859*2c9cb062Sluzhanghpp 
860*2c9cb062Sluzhanghpp   PetscFunctionBegin;
861*2c9cb062Sluzhanghpp   ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
862*2c9cb062Sluzhanghpp   if (!rk->slow){
863*2c9cb062Sluzhanghpp     for (j=0; j<s; j++) w[j] = h*tab->b[j];
864*2c9cb062Sluzhanghpp     ierr = VecGetSubVector(ts->vec_sol,ts->iss,&Xslow);CHKERRQ(ierr);
865*2c9cb062Sluzhanghpp     ierr = VecMAXPY(Xslow,s,w,rk->YdotRHS_slow);CHKERRQ(ierr);
866*2c9cb062Sluzhanghpp     ierr = VecRestoreSubVector(ts->vec_sol,ts->iss,&Xslow);CHKERRQ(ierr);;
867*2c9cb062Sluzhanghpp     } else {;
868*2c9cb062Sluzhanghpp     for (j=0; j<s; j++) w[j] = h/rk->dtratio*tab->b[j];
869*2c9cb062Sluzhanghpp     ierr = VecGetSubVector(X,ts->isf,&Xfast);CHKERRQ(ierr);
870*2c9cb062Sluzhanghpp     ierr = VecMAXPY(Xfast,s,w,rk->YdotRHS_fast);CHKERRQ(ierr);
871*2c9cb062Sluzhanghpp     ierr = VecRestoreSubVector(X,ts->isf,&Xfast);CHKERRQ(ierr);
872*2c9cb062Sluzhanghpp   }
873*2c9cb062Sluzhanghpp   PetscFunctionReturn(0);
874*2c9cb062Sluzhanghpp }
875*2c9cb062Sluzhanghpp 
876*2c9cb062Sluzhanghpp static PetscErrorCode TSStep_RKPMULTIRATE(TS ts)
877*2c9cb062Sluzhanghpp {
878*2c9cb062Sluzhanghpp   TS_RK             *rk = (TS_RK*)ts->data;
879*2c9cb062Sluzhanghpp   RKTableau         tab  = rk->tableau;
880*2c9cb062Sluzhanghpp   Vec               *Y = rk->Y;
881*2c9cb062Sluzhanghpp   Vec               *YdotRHS_fast = rk->YdotRHS_fast,*YdotRHS_slow = rk->YdotRHS_slow;
882*2c9cb062Sluzhanghpp   Vec               Yslow,Yfast;                         /* subvectors store the stges of slow components and fast components respectively                           */
883*2c9cb062Sluzhanghpp   Vec               Yfast_prev,Yfast_new,prev_sol;       /* subvectors store the previous and new solution of fast components and vector store the previous solution */
884*2c9cb062Sluzhanghpp   const PetscInt    s = tab->s;
885*2c9cb062Sluzhanghpp   const PetscReal   *A = tab->A,*c = tab->c;
886*2c9cb062Sluzhanghpp   PetscScalar       *w = rk->work;
887*2c9cb062Sluzhanghpp   PetscInt          i,j,k;
888*2c9cb062Sluzhanghpp   PetscReal         next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step;
889*2c9cb062Sluzhanghpp   PetscErrorCode    ierr;
890*2c9cb062Sluzhanghpp 
891*2c9cb062Sluzhanghpp   PetscFunctionBegin;
892*2c9cb062Sluzhanghpp   rk->status = TS_STEP_INCOMPLETE;
893*2c9cb062Sluzhanghpp   for (i=0; i<s; i++) {
894*2c9cb062Sluzhanghpp     rk->stage_time = t + h*c[i];
895*2c9cb062Sluzhanghpp     ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr);
896*2c9cb062Sluzhanghpp     ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
897*2c9cb062Sluzhanghpp     /*ierr = VecCopy(ts->vec_sol,Yc[i]);CHKERRQ(ierr);*/
898*2c9cb062Sluzhanghpp     ierr = VecGetSubVector(Y[i],ts->iss,&Yslow);CHKERRQ(ierr);
899*2c9cb062Sluzhanghpp     ierr = VecGetSubVector(Y[i],ts->isf,&Yfast);CHKERRQ(ierr);
900*2c9cb062Sluzhanghpp     for (j=0; j<i; j++) w[j] = h*A[i*s+j];
901*2c9cb062Sluzhanghpp     ierr = VecMAXPY(Yslow,i,w,YdotRHS_slow);CHKERRQ(ierr);
902*2c9cb062Sluzhanghpp     ierr = VecMAXPY(Yfast,i,w,YdotRHS_fast);CHKERRQ(ierr);
903*2c9cb062Sluzhanghpp     ierr = VecRestoreSubVector(Y[i],ts->iss,&Yslow);CHKERRQ(ierr);
904*2c9cb062Sluzhanghpp     ierr = VecRestoreSubVector(Y[i],ts->isf,&Yfast);CHKERRQ(ierr);
905*2c9cb062Sluzhanghpp     ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
906*2c9cb062Sluzhanghpp     /* compute the stage RHS for both slow and fast components */
907*2c9cb062Sluzhanghpp     ierr = TSComputeRHSFunctionslow(ts,t+h*c[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr);
908*2c9cb062Sluzhanghpp     ierr = TSComputeRHSFunctionfast(ts,t+h*c[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr);
909*2c9cb062Sluzhanghpp   }
910*2c9cb062Sluzhanghpp   /* store the value of slow components at previous time  */
911*2c9cb062Sluzhanghpp   ierr = VecDuplicate(ts->vec_sol,&prev_sol);CHKERRQ(ierr);
912*2c9cb062Sluzhanghpp   ierr = VecCopy(ts->vec_sol,prev_sol);CHKERRQ(ierr);
913*2c9cb062Sluzhanghpp   rk->slow = 0;
914*2c9cb062Sluzhanghpp   ierr = TSEvaluateStep_RKPMULTIRATE(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
915*2c9cb062Sluzhanghpp   for (k=0; k<rk->dtratio; k++){
916*2c9cb062Sluzhanghpp     for (i=0; i<s; i++) {
917*2c9cb062Sluzhanghpp     rk->stage_time = t + h/rk->dtratio*c[i];
918*2c9cb062Sluzhanghpp     ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr);
919*2c9cb062Sluzhanghpp     ierr = VecCopy(prev_sol,Y[i]);CHKERRQ(ierr);
920*2c9cb062Sluzhanghpp     /* stage value for fast components */
921*2c9cb062Sluzhanghpp     for (j=0; j<i; j++) w[j] = h/rk->dtratio*A[i*s+j];
922*2c9cb062Sluzhanghpp     ierr = VecGetSubVector(Y[i],ts->isf,&Yfast);CHKERRQ(ierr);
923*2c9cb062Sluzhanghpp     ierr = VecMAXPY(Yfast,i,w,YdotRHS_fast);CHKERRQ(ierr);
924*2c9cb062Sluzhanghpp     ierr = VecRestoreSubVector(Y[i],ts->isf,&Yfast);CHKERRQ(ierr);
925*2c9cb062Sluzhanghpp     /* stage value for slow components */
926*2c9cb062Sluzhanghpp     ierr = VecCopy(prev_sol,Y[0]);CHKERRQ(ierr);
927*2c9cb062Sluzhanghpp     ierr = VecGetSubVector(Y[i],ts->iss,&Yslow);CHKERRQ(ierr);
928*2c9cb062Sluzhanghpp     ierr = TSInterpolate_PARTITIONEDMRK(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Yslow);CHKERRQ(ierr);
929*2c9cb062Sluzhanghpp     ierr = VecRestoreSubVector(Y[i],ts->iss,&Yslow);CHKERRQ(ierr);
930*2c9cb062Sluzhanghpp     ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
931*2c9cb062Sluzhanghpp     /* compute the stage RHS for fast components */
932*2c9cb062Sluzhanghpp     ierr = TSComputeRHSFunctionfast(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr);
933*2c9cb062Sluzhanghpp    }
934*2c9cb062Sluzhanghpp     /* update the value of fast components whenusing fast time step */
935*2c9cb062Sluzhanghpp     rk->slow = 1;
936*2c9cb062Sluzhanghpp     ierr = TSEvaluateStep_RKPMULTIRATE(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
937*2c9cb062Sluzhanghpp     ierr = VecGetSubVector(prev_sol,ts->isf,&Yfast_prev);CHKERRQ(ierr);
938*2c9cb062Sluzhanghpp     ierr = VecGetSubVector(ts->vec_sol,ts->isf,&Yfast_new);CHKERRQ(ierr);
939*2c9cb062Sluzhanghpp     ierr = VecCopy(Yfast_new,Yfast_prev);CHKERRQ(ierr);
940*2c9cb062Sluzhanghpp     ierr = VecRestoreSubVector(prev_sol,ts->isf,&Yfast_prev);CHKERRQ(ierr);
941*2c9cb062Sluzhanghpp     ierr = VecRestoreSubVector(ts->vec_sol,ts->isf,&Yfast_new);CHKERRQ(ierr);
942*2c9cb062Sluzhanghpp   }
943*2c9cb062Sluzhanghpp   ts->ptime += ts->time_step;
944*2c9cb062Sluzhanghpp   ts->time_step = next_time_step;
945*2c9cb062Sluzhanghpp   rk->slow = 0;
946*2c9cb062Sluzhanghpp   rk->status = TS_STEP_COMPLETE;
947*2c9cb062Sluzhanghpp   ierr = VecDestroy(&prev_sol);CHKERRQ(ierr);
948*2c9cb062Sluzhanghpp   PetscFunctionReturn(0);
949*2c9cb062Sluzhanghpp }
950*2c9cb062Sluzhanghpp 
951f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts)
952f6a906c0SBarry Smith {
953f6a906c0SBarry Smith   TS_RK          *rk  = (TS_RK*)ts->data;
954be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
955be5899b3SLisandro Dalcin   PetscInt       s   = tab->s;
956f6a906c0SBarry Smith   PetscErrorCode ierr;
957f6a906c0SBarry Smith 
958f6a906c0SBarry Smith   PetscFunctionBegin;
959f6a906c0SBarry Smith   if (ts->adjointsetupcalled++) PetscFunctionReturn(0);
960abc2977eSBarry Smith   ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr);
961f6a906c0SBarry Smith   if(ts->vecs_sensip) {
962abc2977eSBarry Smith     ierr = VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr);
963f6a906c0SBarry Smith   }
964f6a906c0SBarry Smith   PetscFunctionReturn(0);
965f6a906c0SBarry Smith }
966f6a906c0SBarry Smith 
96742f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts)
968d2daff3dSHong Zhang {
969c235aa8dSHong Zhang   TS_RK            *rk  = (TS_RK*)ts->data;
970c235aa8dSHong Zhang   RKTableau        tab  = rk->tableau;
971c235aa8dSHong Zhang   const PetscInt   s    = tab->s;
972c235aa8dSHong Zhang   const PetscReal  *A = tab->A,*b = tab->b,*c = tab->c;
973c235aa8dSHong Zhang   PetscScalar      *w    = rk->work;
9743389c7d9SStefano Zampini   Vec              *Y    = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu;
975b18ea86cSHong Zhang   PetscInt         i,j,nadj;
9763d3eaba7SBarry Smith   PetscReal        t = ts->ptime;
977d2daff3dSHong Zhang   PetscErrorCode   ierr;
978cef76868SBarry Smith   PetscReal        h = ts->time_step;
979c235aa8dSHong Zhang 
980d2daff3dSHong Zhang   PetscFunctionBegin;
981c235aa8dSHong Zhang   rk->status = TS_STEP_INCOMPLETE;
982c235aa8dSHong Zhang   for (i=s-1; i>=0; i--) {
9833389c7d9SStefano Zampini     Mat       J;
9843389c7d9SStefano Zampini     PetscReal stage_time = t + h*(1.0-c[i]);
9853389c7d9SStefano Zampini     PetscBool zero = PETSC_FALSE;
9863389c7d9SStefano Zampini 
9873389c7d9SStefano Zampini     ierr = TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr);
9883389c7d9SStefano Zampini     ierr = TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);CHKERRQ(ierr);
989493ed95fSHong Zhang     if (ts->vec_costintegral) {
9903389c7d9SStefano Zampini       ierr = TSAdjointComputeDRDYFunction(ts,stage_time,Y[i],ts->vecs_drdy);CHKERRQ(ierr);
991493ed95fSHong Zhang     }
9923389c7d9SStefano Zampini     /* Stage values of mu */
9933389c7d9SStefano Zampini     if (ts->vecs_sensip) {
9943389c7d9SStefano Zampini       ierr = TSAdjointComputeRHSJacobian(ts,stage_time,Y[i],ts->Jacp);CHKERRQ(ierr);
9953389c7d9SStefano Zampini       if (ts->vec_costintegral) {
9963389c7d9SStefano Zampini         ierr = TSAdjointComputeDRDPFunction(ts,stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr);
9973389c7d9SStefano Zampini       }
9983389c7d9SStefano Zampini     }
9993389c7d9SStefano Zampini 
10003389c7d9SStefano Zampini     if (b[i] == 0 && i == s-1) zero = PETSC_TRUE;
1001abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
10023389c7d9SStefano Zampini       DM  dm;
10033389c7d9SStefano Zampini       Vec VecSensiTemp;
10043389c7d9SStefano Zampini 
10053389c7d9SStefano Zampini       ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
10063389c7d9SStefano Zampini       ierr = DMGetGlobalVector(dm,&VecSensiTemp);CHKERRQ(ierr);
10073389c7d9SStefano Zampini       /* Stage values of lambda */
10083389c7d9SStefano Zampini       if (!zero) {
10093389c7d9SStefano Zampini         ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp);CHKERRQ(ierr);
10103389c7d9SStefano Zampini         ierr = VecScale(VecSensiTemp,-h*b[i]);CHKERRQ(ierr);
10113389c7d9SStefano Zampini         for (j=i+1; j<s; j++) w[j-i-1]  = -h*A[j*s+i];
10123389c7d9SStefano Zampini         ierr = VecMAXPY(VecSensiTemp,s-i-1,w,&VecDeltaLam[nadj*s+i+1]);CHKERRQ(ierr);
10133389c7d9SStefano Zampini         ierr = MatMultTranspose(J,VecSensiTemp,VecDeltaLam[nadj*s+i]);CHKERRQ(ierr);
10143389c7d9SStefano Zampini       } else {
10153389c7d9SStefano Zampini         ierr = VecSet(VecDeltaLam[nadj*s+i],0);CHKERRQ(ierr);
10163389c7d9SStefano Zampini       }
1017493ed95fSHong Zhang       if (ts->vec_costintegral) {
1018493ed95fSHong Zhang         ierr = VecAXPY(VecDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdy[nadj]);CHKERRQ(ierr);
1019493ed95fSHong Zhang       }
10206497ce14SHong Zhang 
1021ad8e2604SHong Zhang       /* Stage values of mu */
10226497ce14SHong Zhang       if (ts->vecs_sensip) {
10233389c7d9SStefano Zampini         if (!zero) {
10243389c7d9SStefano Zampini           ierr = MatMultTranspose(ts->Jacp,VecSensiTemp,VecDeltaMu[nadj*s+i]);CHKERRQ(ierr);
10253389c7d9SStefano Zampini         } else {
10263389c7d9SStefano Zampini           ierr = VecSet(VecDeltaMu[nadj*s+i],0);CHKERRQ(ierr);
1027493ed95fSHong Zhang         }
1028493ed95fSHong Zhang         if (ts->vec_costintegral) {
1029493ed95fSHong Zhang           ierr = VecAXPY(VecDeltaMu[nadj*s+i],-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr);
1030493ed95fSHong Zhang         }
1031ad8e2604SHong Zhang       }
10323389c7d9SStefano Zampini       ierr = DMRestoreGlobalVector(dm,&VecSensiTemp);CHKERRQ(ierr);
1033c235aa8dSHong Zhang     }
10346497ce14SHong Zhang   }
1035c235aa8dSHong Zhang 
1036c235aa8dSHong Zhang   for (j=0; j<s; j++) w[j] = 1.0;
1037abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
1038b18ea86cSHong Zhang     ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr);
10396497ce14SHong Zhang     if (ts->vecs_sensip) {
1040ad8e2604SHong Zhang       ierr = VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);CHKERRQ(ierr);
1041b18ea86cSHong Zhang     }
10426497ce14SHong Zhang   }
1043c235aa8dSHong Zhang   rk->status = TS_STEP_COMPLETE;
1044d2daff3dSHong Zhang   PetscFunctionReturn(0);
1045d2daff3dSHong Zhang }
1046d2daff3dSHong Zhang 
1047be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts)
1048be5899b3SLisandro Dalcin {
1049be5899b3SLisandro Dalcin   TS_RK          *rk = (TS_RK*)ts->data;
1050be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
1051be5899b3SLisandro Dalcin   PetscErrorCode ierr;
1052be5899b3SLisandro Dalcin 
1053be5899b3SLisandro Dalcin   PetscFunctionBegin;
1054be5899b3SLisandro Dalcin   if (!tab) PetscFunctionReturn(0);
1055be5899b3SLisandro Dalcin   ierr = PetscFree(rk->work);CHKERRQ(ierr);
1056be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr);
1057be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr);
1058*2c9cb062Sluzhanghpp   ierr = VecDestroyVecs(tab->s,&rk->YdotRHS_fast);CHKERRQ(ierr);
1059*2c9cb062Sluzhanghpp   ierr = VecDestroyVecs(tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr);
1060be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr);
1061be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr);
1062be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
1063be5899b3SLisandro Dalcin }
1064be5899b3SLisandro Dalcin 
1065277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts)
1066e4dd521cSBarry Smith {
10675f70b478SJed Brown   TS_RK         *rk = (TS_RK*)ts->data;
10686849ba73SBarry Smith   PetscErrorCode ierr;
1069e4dd521cSBarry Smith 
1070e4dd521cSBarry Smith   PetscFunctionBegin;
1071be5899b3SLisandro Dalcin   ierr = TSRKTableauReset(ts);CHKERRQ(ierr);
10720f7a1166SHong Zhang   ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr);
1073277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
1074e4dd521cSBarry Smith }
1075277b19d0SLisandro Dalcin 
1076f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx)
1077f68a32c8SEmil Constantinescu {
1078f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1079f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1080f68a32c8SEmil Constantinescu }
1081f68a32c8SEmil Constantinescu 
1082f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
1083f68a32c8SEmil Constantinescu {
1084f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1085f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1086f68a32c8SEmil Constantinescu }
1087f68a32c8SEmil Constantinescu 
1088f68a32c8SEmil Constantinescu 
1089f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx)
1090f68a32c8SEmil Constantinescu {
1091f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1092f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1093f68a32c8SEmil Constantinescu }
1094f68a32c8SEmil Constantinescu 
1095f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
1096f68a32c8SEmil Constantinescu {
1097f68a32c8SEmil Constantinescu 
1098f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1099f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1100f68a32c8SEmil Constantinescu }
1101c235aa8dSHong Zhang /*
1102d2daff3dSHong Zhang static PetscErrorCode RKSetAdjCoe(RKTableau tab)
1103d2daff3dSHong Zhang {
1104d2daff3dSHong Zhang   PetscReal *A,*b;
1105d2daff3dSHong Zhang   PetscInt        s,i,j;
1106d2daff3dSHong Zhang   PetscErrorCode  ierr;
1107d2daff3dSHong Zhang 
1108d2daff3dSHong Zhang   PetscFunctionBegin;
1109d2daff3dSHong Zhang   s    = tab->s;
1110d2daff3dSHong Zhang   ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr);
1111d2daff3dSHong Zhang 
1112d2daff3dSHong Zhang   for (i=0; i<s; i++)
1113d2daff3dSHong Zhang     for (j=0; j<s; j++) {
1114d2daff3dSHong 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];
1115d2daff3dSHong Zhang       ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr);
1116d2daff3dSHong Zhang     }
1117d2daff3dSHong Zhang 
1118d2daff3dSHong Zhang   for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i];
1119d2daff3dSHong Zhang 
1120d2daff3dSHong Zhang   ierr  = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
1121d2daff3dSHong Zhang   ierr  = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
1122d2daff3dSHong Zhang   ierr  = PetscFree2(A,b);CHKERRQ(ierr);
1123d2daff3dSHong Zhang   PetscFunctionReturn(0);
1124d2daff3dSHong Zhang }
1125c235aa8dSHong Zhang */
1126be5899b3SLisandro Dalcin 
1127be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts)
1128be5899b3SLisandro Dalcin {
1129be5899b3SLisandro Dalcin   TS_RK          *rk  = (TS_RK*)ts->data;
1130be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
1131*2c9cb062Sluzhanghpp   Vec            YdotRHS_fast,YdotRHS_slow;
1132be5899b3SLisandro Dalcin   PetscErrorCode ierr;
1133be5899b3SLisandro Dalcin 
1134be5899b3SLisandro Dalcin   PetscFunctionBegin;
1135be5899b3SLisandro Dalcin   ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr);
1136be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr);
1137be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr);
1138*2c9cb062Sluzhanghpp   ierr = VecGetSubVector(ts->vec_sol,ts->isf,&YdotRHS_fast);CHKERRQ(ierr);
1139*2c9cb062Sluzhanghpp   ierr = VecDuplicateVecs(YdotRHS_fast,tab->s,&rk->YdotRHS_fast);CHKERRQ(ierr);
1140*2c9cb062Sluzhanghpp   ierr = VecRestoreSubVector(ts->vec_sol,ts->isf,&YdotRHS_fast);CHKERRQ(ierr);
1141*2c9cb062Sluzhanghpp   ierr = VecGetSubVector(ts->vec_sol,ts->iss,&YdotRHS_slow);CHKERRQ(ierr);
1142*2c9cb062Sluzhanghpp   ierr = VecDuplicateVecs(YdotRHS_slow,tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr);
1143*2c9cb062Sluzhanghpp   ierr = VecRestoreSubVector(ts->vec_sol,ts->iss,&YdotRHS_slow);CHKERRQ(ierr);
1144*2c9cb062Sluzhanghpp 
1145be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
1146be5899b3SLisandro Dalcin }
1147be5899b3SLisandro Dalcin 
1148be5899b3SLisandro Dalcin 
1149f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts)
1150f68a32c8SEmil Constantinescu {
1151f68a32c8SEmil Constantinescu   TS_RK          *rk = (TS_RK*)ts->data;
1152f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
1153f68a32c8SEmil Constantinescu   DM             dm;
1154f68a32c8SEmil Constantinescu 
1155f68a32c8SEmil Constantinescu   PetscFunctionBegin;
11563dd83b38SBarry Smith   ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr);
1157be5899b3SLisandro Dalcin   ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);
11580f7a1166SHong Zhang   if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */
11590f7a1166SHong Zhang     ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr);
11600f7a1166SHong Zhang   }
1161f68a32c8SEmil Constantinescu   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1162f68a32c8SEmil Constantinescu   ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
1163f68a32c8SEmil Constantinescu   ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
1164e4dd521cSBarry Smith   PetscFunctionReturn(0);
1165e4dd521cSBarry Smith }
1166d2daff3dSHong Zhang 
1167*2c9cb062Sluzhanghpp /*
1168*2c9cb062Sluzhanghpp   construct a database to choose single-step rk method, or combined rhs mutirate rk method or partitioned rhs rk method
1169*2c9cb062Sluzhanghpp */
1170*2c9cb062Sluzhanghpp const char *const TSRKMultirateTypes[] = {"NONE","COMBINED","PARTITIONED","TSRKMultirateType","RKM_",0};
1171e4dd521cSBarry Smith 
11724416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts)
1173e4dd521cSBarry Smith {
1174be5899b3SLisandro Dalcin   TS_RK          *rk = (TS_RK*)ts->data;
1175dfbe8321SBarry Smith   PetscErrorCode ierr;
1176262ff7bbSBarry Smith 
1177e4dd521cSBarry Smith   PetscFunctionBegin;
1178e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr);
1179f68a32c8SEmil Constantinescu   {
1180f68a32c8SEmil Constantinescu     RKTableauLink link;
1181f68a32c8SEmil Constantinescu     PetscInt      count,choice;
1182f68a32c8SEmil Constantinescu     PetscBool     flg;
1183f68a32c8SEmil Constantinescu     const char    **namelist;
1184*2c9cb062Sluzhanghpp     PetscInt      rkmtype = 0;
1185*2c9cb062Sluzhanghpp 
1186f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) ;
1187f489ac74SBarry Smith     ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr);
1188f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
1189*2c9cb062Sluzhanghpp     ierr = PetscOptionsEList("-ts_rk_multirate_type","Use Multirate or partioned RHS Multirate or single rate RK method","TSRKSetMultirateType",TSRKMultirateTypes,3,TSRKMultirateTypes[0],&rkmtype,&flg);CHKERRQ(ierr);
1190*2c9cb062Sluzhanghpp     if (flg) {ierr = TSRKSetMultirateType(ts,rkmtype);CHKERRQ(ierr);}
1191be5899b3SLisandro Dalcin     ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr);
1192be5899b3SLisandro Dalcin     if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);}
1193f68a32c8SEmil Constantinescu     ierr = PetscFree(namelist);CHKERRQ(ierr);
1194f68a32c8SEmil Constantinescu   }
1195262ff7bbSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
1196*2c9cb062Sluzhanghpp   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Multirate methods options","");CHKERRQ(ierr);
1197*2c9cb062Sluzhanghpp   ierr = PetscOptionsInt("-ts_rk_dtratio","time step ratio between slow and fast","",rk->dtratio,&rk->dtratio,NULL);CHKERRQ(ierr);
1198*2c9cb062Sluzhanghpp   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1199e4dd521cSBarry Smith   PetscFunctionReturn(0);
1200e4dd521cSBarry Smith }
1201e4dd521cSBarry Smith 
12025f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
1203e4dd521cSBarry Smith {
12045f70b478SJed Brown   TS_RK          *rk = (TS_RK*)ts->data;
12058370ee3bSLisandro Dalcin   PetscBool      iascii;
1206dfbe8321SBarry Smith   PetscErrorCode ierr;
12072cdf8120Sasbjorn 
1208e4dd521cSBarry Smith   PetscFunctionBegin;
1209251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
12108370ee3bSLisandro Dalcin   if (iascii) {
12119c334d8fSLisandro Dalcin     RKTableau tab  = rk->tableau;
1212f68a32c8SEmil Constantinescu     TSRKType  rktype;
1213f68a32c8SEmil Constantinescu     char      buf[512];
1214f68a32c8SEmil Constantinescu     ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr);
1215efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  RK type %s\n",rktype);CHKERRQ(ierr);
12160c37eb8eSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  Order: %D\n",tab->order);CHKERRQ(ierr);
12170c37eb8eSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr);
1218f68a32c8SEmil Constantinescu     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
1219f68a32c8SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa c = %s\n",buf);CHKERRQ(ierr);
12208370ee3bSLisandro Dalcin   }
1221f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1222f68a32c8SEmil Constantinescu }
1223f68a32c8SEmil Constantinescu 
1224f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer)
1225f68a32c8SEmil Constantinescu {
1226f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
12279c334d8fSLisandro Dalcin   TSAdapt        adapt;
1228f68a32c8SEmil Constantinescu 
1229f68a32c8SEmil Constantinescu   PetscFunctionBegin;
12309c334d8fSLisandro Dalcin   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
12319c334d8fSLisandro Dalcin   ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr);
1232f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1233f68a32c8SEmil Constantinescu }
1234f68a32c8SEmil Constantinescu 
1235f68a32c8SEmil Constantinescu /*@C
1236f68a32c8SEmil Constantinescu   TSRKSetType - Set the type of RK scheme
1237f68a32c8SEmil Constantinescu 
1238f68a32c8SEmil Constantinescu   Logically collective
1239f68a32c8SEmil Constantinescu 
1240f68a32c8SEmil Constantinescu   Input Parameter:
1241f68a32c8SEmil Constantinescu +  ts - timestepping context
1242f68a32c8SEmil Constantinescu -  rktype - type of RK-scheme
1243f68a32c8SEmil Constantinescu 
1244287c2655SBarry Smith   Options Database:
12459bd3e852SBarry Smith .   -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs>
1246287c2655SBarry Smith 
1247f68a32c8SEmil Constantinescu   Level: intermediate
1248f68a32c8SEmil Constantinescu 
1249287c2655SBarry Smith .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS
1250f68a32c8SEmil Constantinescu @*/
1251f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype)
1252f68a32c8SEmil Constantinescu {
1253f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
1254f68a32c8SEmil Constantinescu 
1255f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1256f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1257b92453a8SLisandro Dalcin   PetscValidCharPointer(rktype,2);
1258f68a32c8SEmil Constantinescu   ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr);
1259f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1260f68a32c8SEmil Constantinescu }
1261f68a32c8SEmil Constantinescu 
1262f68a32c8SEmil Constantinescu /*@C
1263*2c9cb062Sluzhanghpp   TSRKSetMultirateType - Set the type of RK Multirate scheme
1264*2c9cb062Sluzhanghpp 
1265*2c9cb062Sluzhanghpp   Logically collective
1266*2c9cb062Sluzhanghpp 
1267*2c9cb062Sluzhanghpp   Input Parameter:
1268*2c9cb062Sluzhanghpp +  ts - timestepping context
1269*2c9cb062Sluzhanghpp -  rkmtype - type of RKM-scheme
1270*2c9cb062Sluzhanghpp 
1271*2c9cb062Sluzhanghpp   Options Database:
1272*2c9cb062Sluzhanghpp .   -ts_rk_multiarte_type - <none,combined,partitioned>
1273*2c9cb062Sluzhanghpp 
1274*2c9cb062Sluzhanghpp   Level: intermediate
1275*2c9cb062Sluzhanghpp @*/
1276*2c9cb062Sluzhanghpp PetscErrorCode TSRKSetMultirateType(TS ts, TSRKMultirateType rkmtype)
1277*2c9cb062Sluzhanghpp {
1278*2c9cb062Sluzhanghpp   TS_RK *rk = (TS_RK*)ts->data;
1279*2c9cb062Sluzhanghpp   PetscErrorCode ierr;
1280*2c9cb062Sluzhanghpp 
1281*2c9cb062Sluzhanghpp   PetscFunctionBegin;
1282*2c9cb062Sluzhanghpp   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1283*2c9cb062Sluzhanghpp   switch(rkmtype){
1284*2c9cb062Sluzhanghpp     case RKM_NONE:
1285*2c9cb062Sluzhanghpp       break;
1286*2c9cb062Sluzhanghpp     case RKM_COMBINED:
1287*2c9cb062Sluzhanghpp       ts->ops->step           = TSStep_RKMULTIRATE;
1288*2c9cb062Sluzhanghpp       ts->ops->evaluatestep   = TSEvaluateStep_RKMULTIRATE;
1289*2c9cb062Sluzhanghpp       rk->dtratio             = 2;
1290*2c9cb062Sluzhanghpp       ierr = TSRKSetType(ts,TSMRKDefault);CHKERRQ(ierr);
1291*2c9cb062Sluzhanghpp       break;
1292*2c9cb062Sluzhanghpp     case RKM_PARTITIONED:
1293*2c9cb062Sluzhanghpp       ts->ops->step           = TSStep_RKPMULTIRATE;
1294*2c9cb062Sluzhanghpp       ts->ops->evaluatestep   = TSEvaluateStep_RKPMULTIRATE;
1295*2c9cb062Sluzhanghpp       ts->ops->interpolate    = TSInterpolate_PARTITIONEDMRK;
1296*2c9cb062Sluzhanghpp       rk->dtratio             = 2;
1297*2c9cb062Sluzhanghpp       ierr = TSRKSetType(ts,TSMRKDefault);CHKERRQ(ierr);
1298*2c9cb062Sluzhanghpp       break;
1299*2c9cb062Sluzhanghpp     default :
1300*2c9cb062Sluzhanghpp       SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown type '%s'",rkmtype);
1301*2c9cb062Sluzhanghpp   }
1302*2c9cb062Sluzhanghpp   PetscFunctionReturn(0);
1303*2c9cb062Sluzhanghpp }
1304*2c9cb062Sluzhanghpp 
1305*2c9cb062Sluzhanghpp /*@C
1306f68a32c8SEmil Constantinescu   TSRKGetType - Get the type of RK scheme
1307f68a32c8SEmil Constantinescu 
1308f68a32c8SEmil Constantinescu   Logically collective
1309f68a32c8SEmil Constantinescu 
1310f68a32c8SEmil Constantinescu   Input Parameter:
1311f68a32c8SEmil Constantinescu .  ts - timestepping context
1312f68a32c8SEmil Constantinescu 
1313f68a32c8SEmil Constantinescu   Output Parameter:
1314f68a32c8SEmil Constantinescu .  rktype - type of RK-scheme
1315f68a32c8SEmil Constantinescu 
1316f68a32c8SEmil Constantinescu   Level: intermediate
1317f68a32c8SEmil Constantinescu 
1318f68a32c8SEmil Constantinescu .seealso: TSRKGetType()
1319f68a32c8SEmil Constantinescu @*/
1320f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype)
1321f68a32c8SEmil Constantinescu {
1322f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
1323f68a32c8SEmil Constantinescu 
1324f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1325f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1326f68a32c8SEmil Constantinescu   ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr);
1327f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1328f68a32c8SEmil Constantinescu }
1329f68a32c8SEmil Constantinescu 
1330560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype)
1331f68a32c8SEmil Constantinescu {
1332f68a32c8SEmil Constantinescu   TS_RK *rk = (TS_RK*)ts->data;
1333f68a32c8SEmil Constantinescu 
1334f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1335f68a32c8SEmil Constantinescu   *rktype = rk->tableau->name;
1336f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1337f68a32c8SEmil Constantinescu }
1338*2c9cb062Sluzhanghpp 
1339560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype)
1340f68a32c8SEmil Constantinescu {
1341f68a32c8SEmil Constantinescu   TS_RK          *rk = (TS_RK*)ts->data;
1342f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
1343f68a32c8SEmil Constantinescu   PetscBool      match;
1344f68a32c8SEmil Constantinescu   RKTableauLink  link;
1345f68a32c8SEmil Constantinescu 
1346f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1347f68a32c8SEmil Constantinescu   if (rk->tableau) {
1348f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr);
1349f68a32c8SEmil Constantinescu     if (match) PetscFunctionReturn(0);
1350f68a32c8SEmil Constantinescu   }
1351f68a32c8SEmil Constantinescu   for (link = RKTableauList; link; link=link->next) {
1352f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr);
1353f68a32c8SEmil Constantinescu     if (match) {
1354be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);}
1355f68a32c8SEmil Constantinescu       rk->tableau = &link->tab;
1356be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);}
13572ffb9264SLisandro Dalcin       ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
1358f68a32c8SEmil Constantinescu       PetscFunctionReturn(0);
1359f68a32c8SEmil Constantinescu     }
1360f68a32c8SEmil Constantinescu   }
1361f68a32c8SEmil Constantinescu   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype);
1362e4dd521cSBarry Smith   PetscFunctionReturn(0);
1363e4dd521cSBarry Smith }
1364e4dd521cSBarry Smith 
1365ff22ae23SHong Zhang static PetscErrorCode  TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y)
1366ff22ae23SHong Zhang {
1367ff22ae23SHong Zhang   TS_RK *rk = (TS_RK*)ts->data;
1368ff22ae23SHong Zhang 
1369ff22ae23SHong Zhang   PetscFunctionBegin;
13700429704eSStefano Zampini   if (ns) *ns = rk->tableau->s;
1371d2daff3dSHong Zhang   if (Y)   *Y = rk->Y;
1372ff22ae23SHong Zhang   PetscFunctionReturn(0);
1373ff22ae23SHong Zhang }
1374ff22ae23SHong Zhang 
1375b3a6b972SJed Brown static PetscErrorCode TSDestroy_RK(TS ts)
1376b3a6b972SJed Brown {
1377b3a6b972SJed Brown   PetscErrorCode ierr;
1378b3a6b972SJed Brown 
1379b3a6b972SJed Brown   PetscFunctionBegin;
1380b3a6b972SJed Brown   ierr = TSReset_RK(ts);CHKERRQ(ierr);
1381b3a6b972SJed Brown   if (ts->dm) {
1382b3a6b972SJed Brown     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
1383b3a6b972SJed Brown     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
1384b3a6b972SJed Brown   }
1385b3a6b972SJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
1386b3a6b972SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr);
1387b3a6b972SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr);
1388*2c9cb062Sluzhanghpp   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirateType_C",NULL);CHKERRQ(ierr);
1389b3a6b972SJed Brown   PetscFunctionReturn(0);
1390b3a6b972SJed Brown }
1391ff22ae23SHong Zhang 
1392ebe3b25bSBarry Smith /*MC
1393f68a32c8SEmil Constantinescu       TSRK - ODE and DAE solver using Runge-Kutta schemes
1394ebe3b25bSBarry Smith 
1395f68a32c8SEmil Constantinescu   The user should provide the right hand side of the equation
1396f68a32c8SEmil Constantinescu   using TSSetRHSFunction().
1397ebe3b25bSBarry Smith 
1398f68a32c8SEmil Constantinescu   Notes:
139998164e13SLisandro Dalcin   The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type
1400ebe3b25bSBarry Smith 
1401d41222bbSBarry Smith   Level: beginner
1402d41222bbSBarry Smith 
1403f68a32c8SEmil Constantinescu .seealso:  TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3,
1404*2c9cb062Sluzhanghpp            TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister(), TSRKSetMultirateType()
1405ebe3b25bSBarry Smith 
1406ebe3b25bSBarry Smith M*/
14078cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
1408e4dd521cSBarry Smith {
14091566a47fSLisandro Dalcin   TS_RK          *rk;
1410dfbe8321SBarry Smith   PetscErrorCode ierr;
1411e4dd521cSBarry Smith 
1412e4dd521cSBarry Smith   PetscFunctionBegin;
1413f68a32c8SEmil Constantinescu   ierr = TSRKInitializePackage();CHKERRQ(ierr);
1414f68a32c8SEmil Constantinescu 
1415f68a32c8SEmil Constantinescu   ts->ops->reset          = TSReset_RK;
14165f70b478SJed Brown   ts->ops->destroy        = TSDestroy_RK;
14175f70b478SJed Brown   ts->ops->view           = TSView_RK;
1418f68a32c8SEmil Constantinescu   ts->ops->load           = TSLoad_RK;
1419f68a32c8SEmil Constantinescu   ts->ops->setup          = TSSetUp_RK;
142042f2b339SBarry Smith   ts->ops->adjointsetup   = TSAdjointSetUp_RK;
1421f68a32c8SEmil Constantinescu   ts->ops->interpolate    = TSInterpolate_RK;
1422*2c9cb062Sluzhanghpp   ts->ops->step           = TSStep_RK;
1423f68a32c8SEmil Constantinescu   ts->ops->evaluatestep   = TSEvaluateStep_RK;
1424fecfb714SLisandro Dalcin   ts->ops->rollback       = TSRollBack_RK;
1425f68a32c8SEmil Constantinescu   ts->ops->setfromoptions = TSSetFromOptions_RK;
1426ff22ae23SHong Zhang   ts->ops->getstages      = TSGetStages_RK;
142742f2b339SBarry Smith   ts->ops->adjointstep    = TSAdjointStep_RK;
1428e4dd521cSBarry Smith 
14290f7a1166SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_RK;
14300f7a1166SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_RK;
14310f7a1166SHong Zhang 
14321566a47fSLisandro Dalcin   ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr);
14331566a47fSLisandro Dalcin   ts->data = (void*)rk;
1434e4dd521cSBarry Smith 
1435f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr);
1436f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr);
1437*2c9cb062Sluzhanghpp   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirateType_C",TSRKSetMultirateType);CHKERRQ(ierr);
1438be5899b3SLisandro Dalcin 
1439be5899b3SLisandro Dalcin   ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr);
1440e4dd521cSBarry Smith   PetscFunctionReturn(0);
1441e4dd521cSBarry Smith }
1442