xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision 21052c0c12488caedf8314cf02959ecfca93a77d)
1e4dd521cSBarry Smith /*
2474dd773SHong Zhang   Code for time stepping with the Runge-Kutta method
3f68a32c8SEmil Constantinescu 
4f68a32c8SEmil Constantinescu   Notes:
5474dd773SHong Zhang   The general system is written as
6474dd773SHong Zhang 
7474dd773SHong Zhang   Udot = F(t,U)
8474dd773SHong Zhang 
9e4dd521cSBarry Smith */
102c9cb062Sluzhanghpp 
11af0996ceSBarry Smith #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
12f68a32c8SEmil Constantinescu #include <petscdm.h>
13474dd773SHong Zhang #include <../src/ts/impls/explicit/rk/rk.h>
14*21052c0cSHong Zhang #include <../src/ts/impls/explicit/rk/mrk.h>
15f68a32c8SEmil Constantinescu 
16484bcdc7SDebojyoti Ghosh static TSRKType  TSRKDefault = TSRK3BS;
17f68a32c8SEmil Constantinescu static PetscBool TSRKRegisterAllCalled;
18f68a32c8SEmil Constantinescu static PetscBool TSRKPackageInitialized;
19f68a32c8SEmil Constantinescu 
20f68a32c8SEmil Constantinescu /*MC
21287c2655SBarry Smith      TSRK1FE - First order forward Euler scheme.
22262ff7bbSBarry Smith 
23f68a32c8SEmil Constantinescu      This method has one stage.
24f68a32c8SEmil Constantinescu 
25287c2655SBarry Smith      Options database:
269bd3e852SBarry Smith .     -ts_rk_type 1fe
27287c2655SBarry Smith 
28f68a32c8SEmil Constantinescu      Level: advanced
29f68a32c8SEmil Constantinescu 
30287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
31f68a32c8SEmil Constantinescu M*/
32f68a32c8SEmil Constantinescu /*MC
332109b73fSDebojyoti Ghosh      TSRK2A - Second order RK scheme.
34f68a32c8SEmil Constantinescu 
35f68a32c8SEmil Constantinescu      This method has two stages.
36f68a32c8SEmil Constantinescu 
37287c2655SBarry Smith      Options database:
389bd3e852SBarry Smith .     -ts_rk_type 2a
39287c2655SBarry Smith 
40f68a32c8SEmil Constantinescu      Level: advanced
41f68a32c8SEmil Constantinescu 
42287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
43f68a32c8SEmil Constantinescu M*/
44f68a32c8SEmil Constantinescu /*MC
45f68a32c8SEmil Constantinescu      TSRK3 - Third order RK scheme.
46f68a32c8SEmil Constantinescu 
47f68a32c8SEmil Constantinescu      This method has three stages.
48f68a32c8SEmil Constantinescu 
49287c2655SBarry Smith      Options database:
509bd3e852SBarry Smith .     -ts_rk_type 3
51287c2655SBarry Smith 
52f68a32c8SEmil Constantinescu      Level: advanced
53f68a32c8SEmil Constantinescu 
54287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
55f68a32c8SEmil Constantinescu M*/
56f68a32c8SEmil Constantinescu /*MC
572109b73fSDebojyoti Ghosh      TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method.
582109b73fSDebojyoti Ghosh 
59d1905564SLisandro Dalcin      This method has four stages with the First Same As Last (FSAL) property.
602109b73fSDebojyoti Ghosh 
61287c2655SBarry Smith      Options database:
629bd3e852SBarry Smith .     -ts_rk_type 3bs
63287c2655SBarry Smith 
642109b73fSDebojyoti Ghosh      Level: advanced
652109b73fSDebojyoti Ghosh 
6698164e13SLisandro Dalcin      References: https://doi.org/10.1016/0893-9659(89)90079-7
67d1905564SLisandro Dalcin 
68287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
692109b73fSDebojyoti Ghosh M*/
702109b73fSDebojyoti Ghosh /*MC
71f68a32c8SEmil Constantinescu      TSRK4 - Fourth order RK scheme.
72f68a32c8SEmil Constantinescu 
732e698f8bSDebojyoti Ghosh      This is the classical Runge-Kutta method with four stages.
74f68a32c8SEmil Constantinescu 
75287c2655SBarry Smith      Options database:
769bd3e852SBarry Smith .     -ts_rk_type 4
77287c2655SBarry Smith 
78f68a32c8SEmil Constantinescu      Level: advanced
79f68a32c8SEmil Constantinescu 
80287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
81f68a32c8SEmil Constantinescu M*/
82f68a32c8SEmil Constantinescu /*MC
832e698f8bSDebojyoti Ghosh      TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method.
84f68a32c8SEmil Constantinescu 
85f68a32c8SEmil Constantinescu      This method has six stages.
86f68a32c8SEmil Constantinescu 
87287c2655SBarry Smith      Options database:
889bd3e852SBarry Smith .     -ts_rk_type 5f
89287c2655SBarry Smith 
90f68a32c8SEmil Constantinescu      Level: advanced
91f68a32c8SEmil Constantinescu 
92287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
93f68a32c8SEmil Constantinescu M*/
942109b73fSDebojyoti Ghosh /*MC
952e698f8bSDebojyoti Ghosh      TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method.
962109b73fSDebojyoti Ghosh 
97d1905564SLisandro Dalcin      This method has seven stages with the First Same As Last (FSAL) property.
982109b73fSDebojyoti Ghosh 
99287c2655SBarry Smith      Options database:
1009bd3e852SBarry Smith .     -ts_rk_type 5dp
101287c2655SBarry Smith 
1022109b73fSDebojyoti Ghosh      Level: advanced
1032109b73fSDebojyoti Ghosh 
10498164e13SLisandro Dalcin      References: https://doi.org/10.1016/0771-050X(80)90013-3
105d1905564SLisandro Dalcin 
106287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
1072109b73fSDebojyoti Ghosh M*/
10805e23783SLisandro Dalcin /*MC
10905e23783SLisandro Dalcin      TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method.
11005e23783SLisandro Dalcin 
11105e23783SLisandro Dalcin      This method has eight stages with the First Same As Last (FSAL) property.
11205e23783SLisandro Dalcin 
113287c2655SBarry Smith      Options database:
1149bd3e852SBarry Smith .     -ts_rk_type 5bs
115287c2655SBarry Smith 
11605e23783SLisandro Dalcin      Level: advanced
11705e23783SLisandro Dalcin 
11898164e13SLisandro Dalcin      References: https://doi.org/10.1016/0898-1221(96)00141-1
11905e23783SLisandro Dalcin 
120287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
12105e23783SLisandro Dalcin M*/
122262ff7bbSBarry Smith 
123f68a32c8SEmil Constantinescu /*@C
124f68a32c8SEmil Constantinescu   TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK
125262ff7bbSBarry Smith 
126f68a32c8SEmil Constantinescu   Not Collective, but should be called by all processes which will need the schemes to be registered
127262ff7bbSBarry Smith 
128f68a32c8SEmil Constantinescu   Level: advanced
129262ff7bbSBarry Smith 
130f68a32c8SEmil Constantinescu .keywords: TS, TSRK, register, all
131262ff7bbSBarry Smith 
132f68a32c8SEmil Constantinescu .seealso:  TSRKRegisterDestroy()
133262ff7bbSBarry Smith @*/
134f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void)
135262ff7bbSBarry Smith {
1364ac538c5SBarry Smith   PetscErrorCode ierr;
137262ff7bbSBarry Smith 
138262ff7bbSBarry Smith   PetscFunctionBegin;
139f68a32c8SEmil Constantinescu   if (TSRKRegisterAllCalled) PetscFunctionReturn(0);
140f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_TRUE;
141e4dd521cSBarry Smith 
1424e82626cSLisandro Dalcin #define RC PetscRealConstant
143e4dd521cSBarry Smith   {
144f68a32c8SEmil Constantinescu     const PetscReal
1454e82626cSLisandro Dalcin       A[1][1] = {{0}},
1464e82626cSLisandro Dalcin       b[1]    = {RC(1.0)};
1473a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
148e4dd521cSBarry Smith   }
149db046809SBarry Smith   {
150f68a32c8SEmil Constantinescu     const PetscReal
1514e82626cSLisandro Dalcin       A[2][2]   = {{0,0},
1524e82626cSLisandro Dalcin                    {RC(1.0),0}},
1534e82626cSLisandro Dalcin       b[2]      =  {RC(0.5),RC(0.5)},
1544e82626cSLisandro Dalcin       bembed[2] =  {RC(1.0),0};
1553a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
156db046809SBarry Smith   }
157f68a32c8SEmil Constantinescu   {
158f68a32c8SEmil Constantinescu     const PetscReal
15917f6384fSBarry Smith       A[3][3] = {{0,0,0},
1604e82626cSLisandro Dalcin                  {RC(2.0)/RC(3.0),0,0},
1614e82626cSLisandro Dalcin                  {RC(-1.0)/RC(3.0),RC(1.0),0}},
1624e82626cSLisandro Dalcin       b[3]    =  {RC(0.25),RC(0.5),RC(0.25)};
1633a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
164fdefd5e5SDebojyoti Ghosh   }
165fdefd5e5SDebojyoti Ghosh   {
166fdefd5e5SDebojyoti Ghosh     const PetscReal
16717f6384fSBarry Smith       A[4][4]   = {{0,0,0,0},
1684e82626cSLisandro Dalcin                    {RC(1.0)/RC(2.0),0,0,0},
1694e82626cSLisandro Dalcin                    {0,RC(3.0)/RC(4.0),0,0},
1704e82626cSLisandro Dalcin                    {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}},
1714e82626cSLisandro Dalcin       b[4]      =  {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0},
1724e82626cSLisandro 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)};
1733a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
174db046809SBarry Smith   }
175f68a32c8SEmil Constantinescu   {
176f68a32c8SEmil Constantinescu     const PetscReal
177f68a32c8SEmil Constantinescu       A[4][4] = {{0,0,0,0},
1784e82626cSLisandro Dalcin                  {RC(0.5),0,0,0},
1794e82626cSLisandro Dalcin                  {0,RC(0.5),0,0},
1804e82626cSLisandro Dalcin                  {0,0,RC(1.0),0}},
1814e82626cSLisandro 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)};
1823a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
183db046809SBarry Smith   }
184f68a32c8SEmil Constantinescu   {
185f68a32c8SEmil Constantinescu     const PetscReal
18617f6384fSBarry Smith       A[6][6]   = {{0,0,0,0,0,0},
1874e82626cSLisandro Dalcin                    {RC(0.25),0,0,0,0,0},
1884e82626cSLisandro Dalcin                    {RC(3.0)/RC(32.0),RC(9.0)/RC(32.0),0,0,0,0},
1894e82626cSLisandro Dalcin                    {RC(1932.0)/RC(2197.0),RC(-7200.0)/RC(2197.0),RC(7296.0)/RC(2197.0),0,0,0},
1904e82626cSLisandro Dalcin                    {RC(439.0)/RC(216.0),RC(-8.0),RC(3680.0)/RC(513.0),RC(-845.0)/RC(4104.0),0,0},
1914e82626cSLisandro 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}},
1924e82626cSLisandro 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)},
1934e82626cSLisandro 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};
1943a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
195fdefd5e5SDebojyoti Ghosh   }
196fdefd5e5SDebojyoti Ghosh   {
197fdefd5e5SDebojyoti Ghosh     const PetscReal
19817f6384fSBarry Smith       A[7][7]   = {{0,0,0,0,0,0,0},
1994e82626cSLisandro Dalcin                    {RC(1.0)/RC(5.0),0,0,0,0,0,0},
2004e82626cSLisandro Dalcin                    {RC(3.0)/RC(40.0),RC(9.0)/RC(40.0),0,0,0,0,0},
2014e82626cSLisandro Dalcin                    {RC(44.0)/RC(45.0),RC(-56.0)/RC(15.0),RC(32.0)/RC(9.0),0,0,0,0},
2024e82626cSLisandro 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},
2034e82626cSLisandro 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},
2044e82626cSLisandro 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}},
2054e82626cSLisandro 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},
206a685a763Sluzhanghpp         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)},
207a685a763Sluzhanghpp         binterp[7][5] =  {{RC(1.0),RC(-4034104133.0)/RC(1410260304.0),RC(105330401.0)/RC(33982176.0),RC(-13107642775.0)/RC(11282082432.0),RC(6542295.0)/RC(470086768.0)},
208a685a763Sluzhanghpp                     {0,0,0,0,0},
209a685a763Sluzhanghpp                     {0,RC(132343189600.0)/RC(32700410799.0),RC(-833316000.0)/RC(131326951.0),RC(91412856700.0)/RC(32700410799.0),RC(-523383600.0)/RC(10900136933.0)},
210a685a763Sluzhanghpp                     {0,RC(-115792950.0)/RC(29380423.0),RC(185270875.0)/RC(16991088.0),RC(-12653452475.0)/RC(1880347072.0),RC(98134425.0)/RC(235043384.0)},
211a685a763Sluzhanghpp                     {0,RC(70805911779.0)/RC(24914598704.0),RC(-4531260609.0)/RC(600351776.0),RC(988140236175.0)/RC(199316789632.0),RC(-14307999165.0)/RC(24914598704.0)},
212a685a763Sluzhanghpp                     {0,RC(-331320693.0)/RC(205662961.0),RC(31361737.0)/RC(7433601.0),RC(-2426908385.0)/RC(822651844.0),RC(97305120.0)/RC(205662961.0)},
213a685a763Sluzhanghpp                     {0,RC(44764047.0)/RC(29380423.0),RC(-1532549.0)/RC(353981.0),RC(90730570.0)/RC(29380423.0),RC(-8293050.0)/RC(29380423.0)}};
214a685a763Sluzhanghpp 
215a685a763Sluzhanghpp         ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,5,binterp[0]);CHKERRQ(ierr);
216f68a32c8SEmil Constantinescu   }
21705e23783SLisandro Dalcin   {
21805e23783SLisandro Dalcin     const PetscReal
21917f6384fSBarry Smith       A[8][8]   = {{0,0,0,0,0,0,0,0},
2204e82626cSLisandro Dalcin                    {RC(1.0)/RC(6.0),0,0,0,0,0,0,0},
2214e82626cSLisandro Dalcin                    {RC(2.0)/RC(27.0),RC(4.0)/RC(27.0),0,0,0,0,0,0},
2224e82626cSLisandro Dalcin                    {RC(183.0)/RC(1372.0),RC(-162.0)/RC(343.0),RC(1053.0)/RC(1372.0),0,0,0,0,0},
2234e82626cSLisandro 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},
2244e82626cSLisandro 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},
2254e82626cSLisandro 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},
2264e82626cSLisandro 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}},
2274e82626cSLisandro 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},
2284e82626cSLisandro 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)};
22905e23783SLisandro Dalcin     ierr = TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
23005e23783SLisandro Dalcin   }
2314e82626cSLisandro Dalcin #undef RC
232db046809SBarry Smith   PetscFunctionReturn(0);
233db046809SBarry Smith }
234db046809SBarry Smith 
235f68a32c8SEmil Constantinescu /*@C
236f68a32c8SEmil Constantinescu    TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister().
237f68a32c8SEmil Constantinescu 
238f68a32c8SEmil Constantinescu    Not Collective
239f68a32c8SEmil Constantinescu 
240f68a32c8SEmil Constantinescu    Level: advanced
241f68a32c8SEmil Constantinescu 
242f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy
243f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll()
244f68a32c8SEmil Constantinescu @*/
245f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void)
246e4dd521cSBarry Smith {
247f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
248f68a32c8SEmil Constantinescu   RKTableauLink  link;
249f68a32c8SEmil Constantinescu 
250f68a32c8SEmil Constantinescu   PetscFunctionBegin;
251f68a32c8SEmil Constantinescu   while ((link = RKTableauList)) {
252f68a32c8SEmil Constantinescu     RKTableau t = &link->tab;
253f68a32c8SEmil Constantinescu     RKTableauList = link->next;
254f68a32c8SEmil Constantinescu     ierr = PetscFree3(t->A,t->b,t->c);  CHKERRQ(ierr);
255f68a32c8SEmil Constantinescu     ierr = PetscFree (t->bembed);       CHKERRQ(ierr);
256f68a32c8SEmil Constantinescu     ierr = PetscFree (t->binterp);      CHKERRQ(ierr);
257f68a32c8SEmil Constantinescu     ierr = PetscFree (t->name);         CHKERRQ(ierr);
258f68a32c8SEmil Constantinescu     ierr = PetscFree (link);            CHKERRQ(ierr);
259f68a32c8SEmil Constantinescu   }
260f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_FALSE;
261f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
262f68a32c8SEmil Constantinescu }
263f68a32c8SEmil Constantinescu 
264f68a32c8SEmil Constantinescu /*@C
265f68a32c8SEmil Constantinescu   TSRKInitializePackage - This function initializes everything in the TSRK package. It is called
2668a690491SBarry Smith   from TSInitializePackage().
267f68a32c8SEmil Constantinescu 
268f68a32c8SEmil Constantinescu   Level: developer
269f68a32c8SEmil Constantinescu 
270f68a32c8SEmil Constantinescu .keywords: TS, TSRK, initialize, package
271f68a32c8SEmil Constantinescu .seealso: PetscInitialize()
272f68a32c8SEmil Constantinescu @*/
273f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void)
274f68a32c8SEmil Constantinescu {
275186e87acSLisandro Dalcin   PetscErrorCode ierr;
276e4dd521cSBarry Smith 
277e4dd521cSBarry Smith   PetscFunctionBegin;
278f68a32c8SEmil Constantinescu   if (TSRKPackageInitialized) PetscFunctionReturn(0);
279f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_TRUE;
280f68a32c8SEmil Constantinescu   ierr = TSRKRegisterAll();CHKERRQ(ierr);
281f68a32c8SEmil Constantinescu   ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr);
282f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
283f68a32c8SEmil Constantinescu }
284186e87acSLisandro Dalcin 
285f68a32c8SEmil Constantinescu /*@C
286f68a32c8SEmil Constantinescu   TSRKFinalizePackage - This function destroys everything in the TSRK package. It is
287f68a32c8SEmil Constantinescu   called from PetscFinalize().
288186e87acSLisandro Dalcin 
289f68a32c8SEmil Constantinescu   Level: developer
290f68a32c8SEmil Constantinescu 
291f68a32c8SEmil Constantinescu .keywords: Petsc, destroy, package
292f68a32c8SEmil Constantinescu .seealso: PetscFinalize()
293f68a32c8SEmil Constantinescu @*/
294f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void)
295f68a32c8SEmil Constantinescu {
296f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
297f68a32c8SEmil Constantinescu 
298f68a32c8SEmil Constantinescu   PetscFunctionBegin;
299f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_FALSE;
300f68a32c8SEmil Constantinescu   ierr = TSRKRegisterDestroy();CHKERRQ(ierr);
301f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
302f68a32c8SEmil Constantinescu }
303f68a32c8SEmil Constantinescu 
304f68a32c8SEmil Constantinescu /*@C
305f68a32c8SEmil Constantinescu    TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
306f68a32c8SEmil Constantinescu 
307f68a32c8SEmil Constantinescu    Not Collective, but the same schemes should be registered on all processes on which they will be used
308f68a32c8SEmil Constantinescu 
309f68a32c8SEmil Constantinescu    Input Parameters:
310f68a32c8SEmil Constantinescu +  name - identifier for method
311f68a32c8SEmil Constantinescu .  order - approximation order of method
312f68a32c8SEmil Constantinescu .  s - number of stages, this is the dimension of the matrices below
313f68a32c8SEmil Constantinescu .  A - stage coefficients (dimension s*s, row-major)
314f68a32c8SEmil Constantinescu .  b - step completion table (dimension s; NULL to use last row of A)
315f68a32c8SEmil Constantinescu .  c - abscissa (dimension s; NULL to use row sums of A)
316f68a32c8SEmil Constantinescu .  bembed - completion table for embedded method (dimension s; NULL if not available)
3173a8a9803SLisandro Dalcin .  p - Order of the interpolation scheme, equal to the number of columns of binterp
3183a8a9803SLisandro Dalcin -  binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1)
319f68a32c8SEmil Constantinescu 
320f68a32c8SEmil Constantinescu    Notes:
321f68a32c8SEmil Constantinescu    Several RK methods are provided, this function is only needed to create new methods.
322f68a32c8SEmil Constantinescu 
323f68a32c8SEmil Constantinescu    Level: advanced
324f68a32c8SEmil Constantinescu 
325f68a32c8SEmil Constantinescu .keywords: TS, register
326f68a32c8SEmil Constantinescu 
327f68a32c8SEmil Constantinescu .seealso: TSRK
328f68a32c8SEmil Constantinescu @*/
329f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s,
330f68a32c8SEmil Constantinescu                             const PetscReal A[],const PetscReal b[],const PetscReal c[],
3313a8a9803SLisandro Dalcin                             const PetscReal bembed[],PetscInt p,const PetscReal binterp[])
332f68a32c8SEmil Constantinescu {
333f68a32c8SEmil Constantinescu   PetscErrorCode  ierr;
334f68a32c8SEmil Constantinescu   RKTableauLink   link;
335f68a32c8SEmil Constantinescu   RKTableau       t;
336f68a32c8SEmil Constantinescu   PetscInt        i,j;
337f68a32c8SEmil Constantinescu 
338f68a32c8SEmil Constantinescu   PetscFunctionBegin;
3393a8a9803SLisandro Dalcin   PetscValidCharPointer(name,1);
3403a8a9803SLisandro Dalcin   PetscValidRealPointer(A,4);
3413a8a9803SLisandro Dalcin   if (b) PetscValidRealPointer(b,5);
3423a8a9803SLisandro Dalcin   if (c) PetscValidRealPointer(c,6);
3433a8a9803SLisandro Dalcin   if (bembed) PetscValidRealPointer(bembed,7);
3443a8a9803SLisandro Dalcin   if (binterp || p > 1) PetscValidRealPointer(binterp,9);
3453a8a9803SLisandro Dalcin 
3461d36bdfdSBarry Smith   ierr = TSRKInitializePackage();CHKERRQ(ierr);
34795dccacaSBarry Smith   ierr = PetscNew(&link);CHKERRQ(ierr);
348f68a32c8SEmil Constantinescu   t = &link->tab;
3493a8a9803SLisandro Dalcin 
350f68a32c8SEmil Constantinescu   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
351f68a32c8SEmil Constantinescu   t->order = order;
352f68a32c8SEmil Constantinescu   t->s = s;
353dcca6d9dSJed Brown   ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr);
354f68a32c8SEmil Constantinescu   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
355f68a32c8SEmil Constantinescu   if (b)  { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); }
356f68a32c8SEmil Constantinescu   else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i];
357f68a32c8SEmil Constantinescu   if (c)  { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); }
358f68a32c8SEmil 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];
359d760c35bSDebojyoti Ghosh   t->FSAL = PETSC_TRUE;
360d760c35bSDebojyoti Ghosh   for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE;
3613a8a9803SLisandro Dalcin 
362f68a32c8SEmil Constantinescu   if (bembed) {
363785e854fSJed Brown     ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr);
364f68a32c8SEmil Constantinescu     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
365f68a32c8SEmil Constantinescu   }
366f68a32c8SEmil Constantinescu 
3673a8a9803SLisandro Dalcin   if (!binterp) { p = 1; binterp = t->b; }
3683a8a9803SLisandro Dalcin   t->p = p;
3693a8a9803SLisandro Dalcin   ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr);
3703a8a9803SLisandro Dalcin   ierr = PetscMemcpy(t->binterp,binterp,s*p*sizeof(binterp[0]));CHKERRQ(ierr);
3713a8a9803SLisandro Dalcin 
372f68a32c8SEmil Constantinescu   link->next = RKTableauList;
373f68a32c8SEmil Constantinescu   RKTableauList = link;
374f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
375f68a32c8SEmil Constantinescu }
376f68a32c8SEmil Constantinescu 
377e4dd521cSBarry Smith /*
3782c9cb062Sluzhanghpp  This is for single-step RK method
379f68a32c8SEmil Constantinescu  The step completion formula is
380e4dd521cSBarry Smith 
381f68a32c8SEmil Constantinescu  x1 = x0 + h b^T YdotRHS
382f68a32c8SEmil Constantinescu 
383f68a32c8SEmil Constantinescu  This function can be called before or after ts->vec_sol has been updated.
384f68a32c8SEmil Constantinescu  Suppose we have a completion formula (b) and an embedded formula (be) of different order.
385f68a32c8SEmil Constantinescu  We can write
386f68a32c8SEmil Constantinescu 
387f68a32c8SEmil Constantinescu  x1e = x0 + h be^T YdotRHS
388f68a32c8SEmil Constantinescu      = x1 - h b^T YdotRHS + h be^T YdotRHS
389f68a32c8SEmil Constantinescu      = x1 + h (be - b)^T YdotRHS
390f68a32c8SEmil Constantinescu 
391f68a32c8SEmil Constantinescu  so we can evaluate the method with different order even after the step has been optimistically completed.
392e4dd521cSBarry Smith */
393f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done)
394f68a32c8SEmil Constantinescu {
395f68a32c8SEmil Constantinescu   TS_RK          *rk   = (TS_RK*)ts->data;
396f68a32c8SEmil Constantinescu   RKTableau      tab  = rk->tableau;
397f68a32c8SEmil Constantinescu   PetscScalar    *w    = rk->work;
398f68a32c8SEmil Constantinescu   PetscReal      h;
399f68a32c8SEmil Constantinescu   PetscInt       s    = tab->s,j;
400f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
401f68a32c8SEmil Constantinescu 
402f68a32c8SEmil Constantinescu   PetscFunctionBegin;
403f68a32c8SEmil Constantinescu   switch (rk->status) {
404f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
405f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
406f68a32c8SEmil Constantinescu     h = ts->time_step; break;
407f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
408be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
409f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
410f68a32c8SEmil Constantinescu   }
411f68a32c8SEmil Constantinescu   if (order == tab->order) {
412f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) {
413f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
41490b67129SHong Zhang       for (j=0; j<s; j++) w[j] = h*tab->b[j]/rk->dtratio;
415f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
416f68a32c8SEmil Constantinescu     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
417f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
418f68a32c8SEmil Constantinescu   } else if (order == tab->order-1) {
419f68a32c8SEmil Constantinescu     if (!tab->bembed) goto unavailable;
420f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) { /*Complete with the embedded method (be)*/
421f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
422f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
423f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
424f68a32c8SEmil Constantinescu     } else {  /*Rollback and re-complete using (be-b) */
425f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
426f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
427f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
428f68a32c8SEmil Constantinescu     }
429f68a32c8SEmil Constantinescu     if (done) *done = PETSC_TRUE;
430f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
431f68a32c8SEmil Constantinescu   }
432f68a32c8SEmil Constantinescu unavailable:
433f68a32c8SEmil Constantinescu   if (done) *done = PETSC_FALSE;
434a7fac7c2SEmil 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);
435f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
436f68a32c8SEmil Constantinescu }
437f68a32c8SEmil Constantinescu 
4380f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts)
4390f7a1166SHong Zhang {
4400f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
4410f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
4420f7a1166SHong Zhang   const PetscInt  s = tab->s;
4430f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
4440f7a1166SHong Zhang   Vec             *Y = rk->Y;
4450f7a1166SHong Zhang   PetscInt        i;
4460f7a1166SHong Zhang   PetscErrorCode  ierr;
4470f7a1166SHong Zhang 
4480f7a1166SHong Zhang   PetscFunctionBegin;
4490f7a1166SHong Zhang   /* backup cost integral */
4500f7a1166SHong Zhang   ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr);
4510f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
4520f7a1166SHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
45351a7f651SHong Zhang     ierr = TSComputeCostIntegrand(ts,rk->ptime+rk->time_step*c[i],Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
4540f7a1166SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
4550f7a1166SHong Zhang   }
4560f7a1166SHong Zhang   PetscFunctionReturn(0);
4570f7a1166SHong Zhang }
4580f7a1166SHong Zhang 
4590f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts)
4600f7a1166SHong Zhang {
4610f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
4620f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
4630f7a1166SHong Zhang   const PetscInt  s = tab->s;
4640f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
4650f7a1166SHong Zhang   Vec             *Y = rk->Y;
4660f7a1166SHong Zhang   PetscInt        i;
4670f7a1166SHong Zhang   PetscErrorCode  ierr;
4680f7a1166SHong Zhang 
4690f7a1166SHong Zhang   PetscFunctionBegin;
4700f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
4710f7a1166SHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
47251a7f651SHong Zhang     ierr = TSComputeCostIntegrand(ts,ts->ptime+ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
4730f7a1166SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
4740f7a1166SHong Zhang   }
4750f7a1166SHong Zhang   PetscFunctionReturn(0);
4760f7a1166SHong Zhang }
4770f7a1166SHong Zhang 
478fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts)
479fecfb714SLisandro Dalcin {
480fecfb714SLisandro Dalcin   TS_RK           *rk = (TS_RK*)ts->data;
481fecfb714SLisandro Dalcin   RKTableau       tab = rk->tableau;
482fecfb714SLisandro Dalcin   const PetscInt  s  = tab->s;
483fecfb714SLisandro Dalcin   const PetscReal *b = tab->b;
484fecfb714SLisandro Dalcin   PetscScalar     *w = rk->work;
485fecfb714SLisandro Dalcin   Vec             *YdotRHS = rk->YdotRHS;
486fecfb714SLisandro Dalcin   PetscInt        j;
487fecfb714SLisandro Dalcin   PetscReal       h;
488fecfb714SLisandro Dalcin   PetscErrorCode  ierr;
489fecfb714SLisandro Dalcin 
490fecfb714SLisandro Dalcin   PetscFunctionBegin;
491fecfb714SLisandro Dalcin   switch (rk->status) {
492fecfb714SLisandro Dalcin   case TS_STEP_INCOMPLETE:
493fecfb714SLisandro Dalcin   case TS_STEP_PENDING:
494fecfb714SLisandro Dalcin     h = ts->time_step; break;
495fecfb714SLisandro Dalcin   case TS_STEP_COMPLETE:
496fecfb714SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
497fecfb714SLisandro Dalcin   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
498fecfb714SLisandro Dalcin   }
499fecfb714SLisandro Dalcin   for (j=0; j<s; j++) w[j] = -h*b[j];
500fecfb714SLisandro Dalcin   ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr);
501fecfb714SLisandro Dalcin   PetscFunctionReturn(0);
502fecfb714SLisandro Dalcin }
503fecfb714SLisandro Dalcin 
504f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts)
505f68a32c8SEmil Constantinescu {
506f68a32c8SEmil Constantinescu   TS_RK           *rk  = (TS_RK*)ts->data;
507f68a32c8SEmil Constantinescu   RKTableau       tab  = rk->tableau;
508f68a32c8SEmil Constantinescu   const PetscInt  s = tab->s;
509fecfb714SLisandro Dalcin   const PetscReal *A = tab->A,*c = tab->c;
510f68a32c8SEmil Constantinescu   PetscScalar     *w = rk->work;
511f68a32c8SEmil Constantinescu   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
512d1905564SLisandro Dalcin   PetscBool       FSAL = tab->FSAL;
513f68a32c8SEmil Constantinescu   TSAdapt         adapt;
514fecfb714SLisandro Dalcin   PetscInt        i,j;
515be5899b3SLisandro Dalcin   PetscInt        rejections = 0;
516be5899b3SLisandro Dalcin   PetscBool       stageok,accept = PETSC_TRUE;
517be5899b3SLisandro Dalcin   PetscReal       next_time_step = ts->time_step;
518f68a32c8SEmil Constantinescu   PetscErrorCode  ierr;
519f68a32c8SEmil Constantinescu 
520f68a32c8SEmil Constantinescu   PetscFunctionBegin;
521d1905564SLisandro Dalcin   if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE;
522d1905564SLisandro Dalcin   if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); }
523d1905564SLisandro Dalcin 
524f68a32c8SEmil Constantinescu   rk->status = TS_STEP_INCOMPLETE;
525be5899b3SLisandro Dalcin   while (!ts->reason && rk->status != TS_STEP_COMPLETE) {
526be5899b3SLisandro Dalcin     PetscReal t = ts->ptime;
527f68a32c8SEmil Constantinescu     PetscReal h = ts->time_step;
528f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
5299be3e283SDebojyoti Ghosh       rk->stage_time = t + h*c[i];
5309be3e283SDebojyoti Ghosh       ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr);
531f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
532f68a32c8SEmil Constantinescu       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
533f68a32c8SEmil Constantinescu       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
5349be3e283SDebojyoti Ghosh       ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
535f68a32c8SEmil Constantinescu       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
536be5899b3SLisandro Dalcin       ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr);
537fecfb714SLisandro Dalcin       if (!stageok) goto reject_step;
5388f738a7cSLisandro Dalcin       if (FSAL && !i) continue;
539f68a32c8SEmil Constantinescu       ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
540f68a32c8SEmil Constantinescu     }
541be5899b3SLisandro Dalcin 
542be5899b3SLisandro Dalcin     rk->status = TS_STEP_INCOMPLETE;
543f68a32c8SEmil Constantinescu     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
544f68a32c8SEmil Constantinescu     rk->status = TS_STEP_PENDING;
545f68a32c8SEmil Constantinescu     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
546f68a32c8SEmil Constantinescu     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
5471917a363SLisandro Dalcin     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr);
548fecfb714SLisandro Dalcin     ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
549be5899b3SLisandro Dalcin     rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
550be5899b3SLisandro Dalcin     if (!accept) { /* Roll back the current step */
551fecfb714SLisandro Dalcin       ierr = TSRollBack_RK(ts);CHKERRQ(ierr);
552be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
553be5899b3SLisandro Dalcin       goto reject_step;
554be5899b3SLisandro Dalcin     }
555be5899b3SLisandro Dalcin 
5560f7a1166SHong Zhang     if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
5570f7a1166SHong Zhang       rk->ptime     = ts->ptime;
5580f7a1166SHong Zhang       rk->time_step = ts->time_step;
559493ed95fSHong Zhang     }
560be5899b3SLisandro Dalcin 
561f68a32c8SEmil Constantinescu     ts->ptime += ts->time_step;
562f68a32c8SEmil Constantinescu     ts->time_step = next_time_step;
563f68a32c8SEmil Constantinescu     break;
564be5899b3SLisandro Dalcin 
565be5899b3SLisandro Dalcin     reject_step:
566fecfb714SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
567be5899b3SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
568be5899b3SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
569be5899b3SLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
570f68a32c8SEmil Constantinescu     }
571f68a32c8SEmil Constantinescu   }
572f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
573e4dd521cSBarry Smith }
574e4dd521cSBarry Smith 
575f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts)
576f6a906c0SBarry Smith {
577f6a906c0SBarry Smith   TS_RK          *rk  = (TS_RK*)ts->data;
578be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
579be5899b3SLisandro Dalcin   PetscInt       s   = tab->s;
580f6a906c0SBarry Smith   PetscErrorCode ierr;
581f6a906c0SBarry Smith 
582f6a906c0SBarry Smith   PetscFunctionBegin;
583f6a906c0SBarry Smith   if (ts->adjointsetupcalled++) PetscFunctionReturn(0);
584abc2977eSBarry Smith   ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr);
585f6a906c0SBarry Smith   if(ts->vecs_sensip) {
586abc2977eSBarry Smith     ierr = VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr);
587f6a906c0SBarry Smith   }
588f6a906c0SBarry Smith   PetscFunctionReturn(0);
589f6a906c0SBarry Smith }
590f6a906c0SBarry Smith 
59142f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts)
592d2daff3dSHong Zhang {
593c235aa8dSHong Zhang   TS_RK            *rk  = (TS_RK*)ts->data;
594c235aa8dSHong Zhang   RKTableau        tab  = rk->tableau;
595c235aa8dSHong Zhang   const PetscInt   s    = tab->s;
596c235aa8dSHong Zhang   const PetscReal  *A = tab->A,*b = tab->b,*c = tab->c;
597c235aa8dSHong Zhang   PetscScalar      *w    = rk->work;
5983389c7d9SStefano Zampini   Vec              *Y    = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu;
599b18ea86cSHong Zhang   PetscInt         i,j,nadj;
6003d3eaba7SBarry Smith   PetscReal        t = ts->ptime;
601d2daff3dSHong Zhang   PetscErrorCode   ierr;
602cef76868SBarry Smith   PetscReal        h = ts->time_step;
603c235aa8dSHong Zhang 
604d2daff3dSHong Zhang   PetscFunctionBegin;
605c235aa8dSHong Zhang   rk->status = TS_STEP_INCOMPLETE;
606c235aa8dSHong Zhang   for (i=s-1; i>=0; i--) {
6073389c7d9SStefano Zampini     Mat       J;
6083389c7d9SStefano Zampini     PetscReal stage_time = t + h*(1.0-c[i]);
6093389c7d9SStefano Zampini     PetscBool zero = PETSC_FALSE;
6103389c7d9SStefano Zampini 
6113389c7d9SStefano Zampini     ierr = TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr);
6123389c7d9SStefano Zampini     ierr = TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);CHKERRQ(ierr);
613493ed95fSHong Zhang     if (ts->vec_costintegral) {
6143389c7d9SStefano Zampini       ierr = TSAdjointComputeDRDYFunction(ts,stage_time,Y[i],ts->vecs_drdy);CHKERRQ(ierr);
615493ed95fSHong Zhang     }
6163389c7d9SStefano Zampini     /* Stage values of mu */
6173389c7d9SStefano Zampini     if (ts->vecs_sensip) {
6183389c7d9SStefano Zampini       ierr = TSAdjointComputeRHSJacobian(ts,stage_time,Y[i],ts->Jacp);CHKERRQ(ierr);
6193389c7d9SStefano Zampini       if (ts->vec_costintegral) {
6203389c7d9SStefano Zampini         ierr = TSAdjointComputeDRDPFunction(ts,stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr);
6213389c7d9SStefano Zampini       }
6223389c7d9SStefano Zampini     }
6233389c7d9SStefano Zampini 
6243389c7d9SStefano Zampini     if (b[i] == 0 && i == s-1) zero = PETSC_TRUE;
625abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
6263389c7d9SStefano Zampini       DM  dm;
6273389c7d9SStefano Zampini       Vec VecSensiTemp;
6283389c7d9SStefano Zampini 
6293389c7d9SStefano Zampini       ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
6303389c7d9SStefano Zampini       ierr = DMGetGlobalVector(dm,&VecSensiTemp);CHKERRQ(ierr);
6313389c7d9SStefano Zampini       /* Stage values of lambda */
6323389c7d9SStefano Zampini       if (!zero) {
6333389c7d9SStefano Zampini         ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp);CHKERRQ(ierr);
6343389c7d9SStefano Zampini         ierr = VecScale(VecSensiTemp,-h*b[i]);CHKERRQ(ierr);
6353389c7d9SStefano Zampini         for (j=i+1; j<s; j++) w[j-i-1]  = -h*A[j*s+i];
6363389c7d9SStefano Zampini         ierr = VecMAXPY(VecSensiTemp,s-i-1,w,&VecDeltaLam[nadj*s+i+1]);CHKERRQ(ierr);
6373389c7d9SStefano Zampini         ierr = MatMultTranspose(J,VecSensiTemp,VecDeltaLam[nadj*s+i]);CHKERRQ(ierr);
6383389c7d9SStefano Zampini       } else {
6393389c7d9SStefano Zampini         ierr = VecSet(VecDeltaLam[nadj*s+i],0);CHKERRQ(ierr);
6403389c7d9SStefano Zampini       }
641493ed95fSHong Zhang       if (ts->vec_costintegral) {
642493ed95fSHong Zhang         ierr = VecAXPY(VecDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdy[nadj]);CHKERRQ(ierr);
643493ed95fSHong Zhang       }
6446497ce14SHong Zhang 
645ad8e2604SHong Zhang       /* Stage values of mu */
6466497ce14SHong Zhang       if (ts->vecs_sensip) {
6473389c7d9SStefano Zampini         if (!zero) {
6483389c7d9SStefano Zampini           ierr = MatMultTranspose(ts->Jacp,VecSensiTemp,VecDeltaMu[nadj*s+i]);CHKERRQ(ierr);
6493389c7d9SStefano Zampini         } else {
6503389c7d9SStefano Zampini           ierr = VecSet(VecDeltaMu[nadj*s+i],0);CHKERRQ(ierr);
651493ed95fSHong Zhang         }
652493ed95fSHong Zhang         if (ts->vec_costintegral) {
653493ed95fSHong Zhang           ierr = VecAXPY(VecDeltaMu[nadj*s+i],-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr);
654493ed95fSHong Zhang         }
655ad8e2604SHong Zhang       }
6563389c7d9SStefano Zampini       ierr = DMRestoreGlobalVector(dm,&VecSensiTemp);CHKERRQ(ierr);
657c235aa8dSHong Zhang     }
6586497ce14SHong Zhang   }
659c235aa8dSHong Zhang 
660c235aa8dSHong Zhang   for (j=0; j<s; j++) w[j] = 1.0;
661abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
662b18ea86cSHong Zhang     ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr);
6636497ce14SHong Zhang     if (ts->vecs_sensip) {
664ad8e2604SHong Zhang       ierr = VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);CHKERRQ(ierr);
665b18ea86cSHong Zhang     }
6666497ce14SHong Zhang   }
667c235aa8dSHong Zhang   rk->status = TS_STEP_COMPLETE;
668d2daff3dSHong Zhang   PetscFunctionReturn(0);
669d2daff3dSHong Zhang }
670d2daff3dSHong Zhang 
67155de54fcSHong Zhang static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X)
67255de54fcSHong Zhang {
67355de54fcSHong Zhang   TS_RK            *rk = (TS_RK*)ts->data;
67455de54fcSHong Zhang   PetscInt         s  = rk->tableau->s,p = rk->tableau->p,i,j;
67555de54fcSHong Zhang   PetscReal        h;
67655de54fcSHong Zhang   PetscReal        tt,t;
67755de54fcSHong Zhang   PetscScalar      *b;
67855de54fcSHong Zhang   const PetscReal  *B = rk->tableau->binterp;
67955de54fcSHong Zhang   PetscErrorCode   ierr;
68055de54fcSHong Zhang 
68155de54fcSHong Zhang   PetscFunctionBegin;
68255de54fcSHong Zhang   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
68355de54fcSHong Zhang 
68455de54fcSHong Zhang   switch (rk->status) {
68555de54fcSHong Zhang     case TS_STEP_INCOMPLETE:
68655de54fcSHong Zhang     case TS_STEP_PENDING:
68755de54fcSHong Zhang       h = ts->time_step;
68855de54fcSHong Zhang       t = (itime - ts->ptime)/h;
68955de54fcSHong Zhang       break;
69055de54fcSHong Zhang     case TS_STEP_COMPLETE:
69155de54fcSHong Zhang       h = ts->ptime - ts->ptime_prev;
69255de54fcSHong Zhang       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
69355de54fcSHong Zhang       break;
69455de54fcSHong Zhang     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
69555de54fcSHong Zhang   }
69655de54fcSHong Zhang   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
69755de54fcSHong Zhang   for (i=0; i<s; i++) b[i] = 0;
69855de54fcSHong Zhang   for (j=0,tt=t; j<p; j++,tt*=t) {
69955de54fcSHong Zhang     for (i=0; i<s; i++) {
70055de54fcSHong Zhang       b[i]  += h * B[i*p+j] * tt;
70155de54fcSHong Zhang     }
70255de54fcSHong Zhang   }
70355de54fcSHong Zhang   ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr);
70455de54fcSHong Zhang   ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr);
70555de54fcSHong Zhang   ierr = PetscFree(b);CHKERRQ(ierr);
70655de54fcSHong Zhang   PetscFunctionReturn(0);
70755de54fcSHong Zhang }
70855de54fcSHong Zhang 
70955de54fcSHong Zhang /*------------------------------------------------------------*/
71055de54fcSHong Zhang 
711be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts)
712be5899b3SLisandro Dalcin {
713be5899b3SLisandro Dalcin   TS_RK          *rk = (TS_RK*)ts->data;
714be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
715be5899b3SLisandro Dalcin   PetscErrorCode ierr;
716be5899b3SLisandro Dalcin 
717be5899b3SLisandro Dalcin   PetscFunctionBegin;
718be5899b3SLisandro Dalcin   if (!tab) PetscFunctionReturn(0);
719be5899b3SLisandro Dalcin   ierr = PetscFree(rk->work);CHKERRQ(ierr);
720be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr);
721be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr);
722be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr);
723be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr);
724be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
725be5899b3SLisandro Dalcin }
726be5899b3SLisandro Dalcin 
727277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts)
728e4dd521cSBarry Smith {
7295f70b478SJed Brown   TS_RK         *rk = (TS_RK*)ts->data;
7306849ba73SBarry Smith   PetscErrorCode ierr;
731e4dd521cSBarry Smith 
732e4dd521cSBarry Smith   PetscFunctionBegin;
733be5899b3SLisandro Dalcin   ierr = TSRKTableauReset(ts);CHKERRQ(ierr);
7340f7a1166SHong Zhang   ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr);
7350fe4d17eSHong Zhang   if (ts->use_splitrhsfunction) {
73655de54fcSHong Zhang     ierr = PetscTryMethod(ts,"TSReset_MRKSPLIT_C",(TS),(ts));CHKERRQ(ierr);
7370fe4d17eSHong Zhang   } else {
73855de54fcSHong Zhang     ierr = PetscTryMethod(ts,"TSReset_MRKNONSPLIT_C",(TS),(ts));CHKERRQ(ierr);
7390fe4d17eSHong Zhang   }
740277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
741e4dd521cSBarry Smith }
742277b19d0SLisandro Dalcin 
743f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx)
744f68a32c8SEmil Constantinescu {
745f68a32c8SEmil Constantinescu   PetscFunctionBegin;
746f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
747f68a32c8SEmil Constantinescu }
748f68a32c8SEmil Constantinescu 
749f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
750f68a32c8SEmil Constantinescu {
751f68a32c8SEmil Constantinescu   PetscFunctionBegin;
752f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
753f68a32c8SEmil Constantinescu }
754f68a32c8SEmil Constantinescu 
755f68a32c8SEmil Constantinescu 
756f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx)
757f68a32c8SEmil Constantinescu {
758f68a32c8SEmil Constantinescu   PetscFunctionBegin;
759f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
760f68a32c8SEmil Constantinescu }
761f68a32c8SEmil Constantinescu 
762f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
763f68a32c8SEmil Constantinescu {
764f68a32c8SEmil Constantinescu 
765f68a32c8SEmil Constantinescu   PetscFunctionBegin;
766f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
767f68a32c8SEmil Constantinescu }
768c235aa8dSHong Zhang /*
769d2daff3dSHong Zhang static PetscErrorCode RKSetAdjCoe(RKTableau tab)
770d2daff3dSHong Zhang {
771d2daff3dSHong Zhang   PetscReal *A,*b;
772d2daff3dSHong Zhang   PetscInt        s,i,j;
773d2daff3dSHong Zhang   PetscErrorCode  ierr;
774d2daff3dSHong Zhang 
775d2daff3dSHong Zhang   PetscFunctionBegin;
776d2daff3dSHong Zhang   s    = tab->s;
777d2daff3dSHong Zhang   ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr);
778d2daff3dSHong Zhang 
779d2daff3dSHong Zhang   for (i=0; i<s; i++)
780d2daff3dSHong Zhang     for (j=0; j<s; j++) {
781d2daff3dSHong 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];
782d2daff3dSHong Zhang       ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr);
783d2daff3dSHong Zhang     }
784d2daff3dSHong Zhang 
785d2daff3dSHong Zhang   for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i];
786d2daff3dSHong Zhang 
787d2daff3dSHong Zhang   ierr  = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
788d2daff3dSHong Zhang   ierr  = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
789d2daff3dSHong Zhang   ierr  = PetscFree2(A,b);CHKERRQ(ierr);
790d2daff3dSHong Zhang   PetscFunctionReturn(0);
791d2daff3dSHong Zhang }
792c235aa8dSHong Zhang */
793be5899b3SLisandro Dalcin 
794be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts)
795be5899b3SLisandro Dalcin {
796be5899b3SLisandro Dalcin   TS_RK          *rk  = (TS_RK*)ts->data;
797be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
798be5899b3SLisandro Dalcin   PetscErrorCode ierr;
799be5899b3SLisandro Dalcin 
800be5899b3SLisandro Dalcin   PetscFunctionBegin;
801be5899b3SLisandro Dalcin   ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr);
802be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr);
803be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr);
804be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
805be5899b3SLisandro Dalcin }
806be5899b3SLisandro Dalcin 
807f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts)
808f68a32c8SEmil Constantinescu {
809f68a32c8SEmil Constantinescu   TS_RK          *rk = (TS_RK*)ts->data;
810f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
811f68a32c8SEmil Constantinescu   DM             dm;
812f68a32c8SEmil Constantinescu 
813f68a32c8SEmil Constantinescu   PetscFunctionBegin;
8143dd83b38SBarry Smith   ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr);
815be5899b3SLisandro Dalcin   ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);
8160f7a1166SHong Zhang   if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */
8170f7a1166SHong Zhang     ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr);
8180f7a1166SHong Zhang   }
819f68a32c8SEmil Constantinescu   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
820f68a32c8SEmil Constantinescu   ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
821f68a32c8SEmil Constantinescu   ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
8220fe4d17eSHong Zhang   if (ts->use_splitrhsfunction) {
82355de54fcSHong Zhang     ierr = PetscTryMethod(ts,"TSSetUp_MRKSPLIT_C",(TS),(ts));CHKERRQ(ierr);
8240fe4d17eSHong Zhang   } else {
82555de54fcSHong Zhang     ierr = PetscTryMethod(ts,"TSSetUp_MRKNONSPLIT_C",(TS),(ts));CHKERRQ(ierr);
8260fe4d17eSHong Zhang   }
827e4dd521cSBarry Smith   PetscFunctionReturn(0);
828e4dd521cSBarry Smith }
829d2daff3dSHong Zhang 
8304416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts)
831e4dd521cSBarry Smith {
832be5899b3SLisandro Dalcin   TS_RK          *rk = (TS_RK*)ts->data;
833dfbe8321SBarry Smith   PetscErrorCode ierr;
834262ff7bbSBarry Smith 
835e4dd521cSBarry Smith   PetscFunctionBegin;
836e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr);
837f68a32c8SEmil Constantinescu   {
838f68a32c8SEmil Constantinescu     RKTableauLink link;
839f68a32c8SEmil Constantinescu     PetscInt      count,choice;
8400fe4d17eSHong Zhang     PetscBool     flg,use_multirate = PETSC_FALSE;
841f68a32c8SEmil Constantinescu     const char    **namelist;
8422c9cb062Sluzhanghpp 
843f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) ;
844f489ac74SBarry Smith     ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr);
845f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
8460fe4d17eSHong Zhang     ierr = PetscOptionsBool("-ts_rk_multirate","Use interpolation-based multirate RK method","TSRKSetMultirate",rk->use_multirate,&use_multirate,&flg);CHKERRQ(ierr);
8470fe4d17eSHong Zhang     if (flg) {
8480fe4d17eSHong Zhang       ierr = TSRKSetMultirate(ts,use_multirate);CHKERRQ(ierr);
8490fe4d17eSHong Zhang     }
850be5899b3SLisandro Dalcin     ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr);
851be5899b3SLisandro Dalcin     if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);}
852f68a32c8SEmil Constantinescu     ierr = PetscFree(namelist);CHKERRQ(ierr);
853f68a32c8SEmil Constantinescu   }
854262ff7bbSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
8552c9cb062Sluzhanghpp   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Multirate methods options","");CHKERRQ(ierr);
8562c9cb062Sluzhanghpp   ierr = PetscOptionsInt("-ts_rk_dtratio","time step ratio between slow and fast","",rk->dtratio,&rk->dtratio,NULL);CHKERRQ(ierr);
8572c9cb062Sluzhanghpp   ierr = PetscOptionsEnd();CHKERRQ(ierr);
858e4dd521cSBarry Smith   PetscFunctionReturn(0);
859e4dd521cSBarry Smith }
860e4dd521cSBarry Smith 
8615f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
862e4dd521cSBarry Smith {
8635f70b478SJed Brown   TS_RK          *rk = (TS_RK*)ts->data;
8648370ee3bSLisandro Dalcin   PetscBool      iascii;
865dfbe8321SBarry Smith   PetscErrorCode ierr;
8662cdf8120Sasbjorn 
867e4dd521cSBarry Smith   PetscFunctionBegin;
868251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
8698370ee3bSLisandro Dalcin   if (iascii) {
8709c334d8fSLisandro Dalcin     RKTableau tab  = rk->tableau;
871f68a32c8SEmil Constantinescu     TSRKType  rktype;
872f68a32c8SEmil Constantinescu     char      buf[512];
873f68a32c8SEmil Constantinescu     ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr);
874efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  RK type %s\n",rktype);CHKERRQ(ierr);
8750c37eb8eSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  Order: %D\n",tab->order);CHKERRQ(ierr);
8760c37eb8eSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr);
877f68a32c8SEmil Constantinescu     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
878f68a32c8SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa c = %s\n",buf);CHKERRQ(ierr);
8798370ee3bSLisandro Dalcin   }
880f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
881f68a32c8SEmil Constantinescu }
882f68a32c8SEmil Constantinescu 
883f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer)
884f68a32c8SEmil Constantinescu {
885f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
8869c334d8fSLisandro Dalcin   TSAdapt        adapt;
887f68a32c8SEmil Constantinescu 
888f68a32c8SEmil Constantinescu   PetscFunctionBegin;
8899c334d8fSLisandro Dalcin   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
8909c334d8fSLisandro Dalcin   ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr);
891f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
892f68a32c8SEmil Constantinescu }
893f68a32c8SEmil Constantinescu 
894f68a32c8SEmil Constantinescu /*@C
895f68a32c8SEmil Constantinescu   TSRKSetType - Set the type of RK scheme
896f68a32c8SEmil Constantinescu 
897f68a32c8SEmil Constantinescu   Logically collective
898f68a32c8SEmil Constantinescu 
899f68a32c8SEmil Constantinescu   Input Parameter:
900f68a32c8SEmil Constantinescu +  ts - timestepping context
901f68a32c8SEmil Constantinescu -  rktype - type of RK-scheme
902f68a32c8SEmil Constantinescu 
903287c2655SBarry Smith   Options Database:
9049bd3e852SBarry Smith .   -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs>
905287c2655SBarry Smith 
906f68a32c8SEmil Constantinescu   Level: intermediate
907f68a32c8SEmil Constantinescu 
908287c2655SBarry Smith .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS
909f68a32c8SEmil Constantinescu @*/
910f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype)
911f68a32c8SEmil Constantinescu {
912f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
913f68a32c8SEmil Constantinescu 
914f68a32c8SEmil Constantinescu   PetscFunctionBegin;
915f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
916b92453a8SLisandro Dalcin   PetscValidCharPointer(rktype,2);
917f68a32c8SEmil Constantinescu   ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr);
918f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
919f68a32c8SEmil Constantinescu }
920f68a32c8SEmil Constantinescu 
921f68a32c8SEmil Constantinescu /*@C
922f68a32c8SEmil Constantinescu   TSRKGetType - Get the type of RK scheme
923f68a32c8SEmil Constantinescu 
924f68a32c8SEmil Constantinescu   Logically collective
925f68a32c8SEmil Constantinescu 
926f68a32c8SEmil Constantinescu   Input Parameter:
927f68a32c8SEmil Constantinescu .  ts - timestepping context
928f68a32c8SEmil Constantinescu 
929f68a32c8SEmil Constantinescu   Output Parameter:
930f68a32c8SEmil Constantinescu .  rktype - type of RK-scheme
931f68a32c8SEmil Constantinescu 
932f68a32c8SEmil Constantinescu   Level: intermediate
933f68a32c8SEmil Constantinescu 
934f68a32c8SEmil Constantinescu .seealso: TSRKGetType()
935f68a32c8SEmil Constantinescu @*/
936f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype)
937f68a32c8SEmil Constantinescu {
938f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
939f68a32c8SEmil Constantinescu 
940f68a32c8SEmil Constantinescu   PetscFunctionBegin;
941f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
942f68a32c8SEmil Constantinescu   ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr);
943f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
944f68a32c8SEmil Constantinescu }
945f68a32c8SEmil Constantinescu 
946560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype)
947f68a32c8SEmil Constantinescu {
948f68a32c8SEmil Constantinescu   TS_RK *rk = (TS_RK*)ts->data;
949f68a32c8SEmil Constantinescu 
950f68a32c8SEmil Constantinescu   PetscFunctionBegin;
951f68a32c8SEmil Constantinescu   *rktype = rk->tableau->name;
952f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
953f68a32c8SEmil Constantinescu }
9542c9cb062Sluzhanghpp 
955560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype)
956f68a32c8SEmil Constantinescu {
957f68a32c8SEmil Constantinescu   TS_RK          *rk = (TS_RK*)ts->data;
958f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
959f68a32c8SEmil Constantinescu   PetscBool      match;
960f68a32c8SEmil Constantinescu   RKTableauLink  link;
961f68a32c8SEmil Constantinescu 
962f68a32c8SEmil Constantinescu   PetscFunctionBegin;
963f68a32c8SEmil Constantinescu   if (rk->tableau) {
964f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr);
965f68a32c8SEmil Constantinescu     if (match) PetscFunctionReturn(0);
966f68a32c8SEmil Constantinescu   }
967f68a32c8SEmil Constantinescu   for (link = RKTableauList; link; link=link->next) {
968f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr);
969f68a32c8SEmil Constantinescu     if (match) {
970be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);}
971f68a32c8SEmil Constantinescu       rk->tableau = &link->tab;
972be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);}
9732ffb9264SLisandro Dalcin       ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
974f68a32c8SEmil Constantinescu       PetscFunctionReturn(0);
975f68a32c8SEmil Constantinescu     }
976f68a32c8SEmil Constantinescu   }
977f68a32c8SEmil Constantinescu   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype);
978e4dd521cSBarry Smith   PetscFunctionReturn(0);
979e4dd521cSBarry Smith }
980e4dd521cSBarry Smith 
981ff22ae23SHong Zhang static PetscErrorCode  TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y)
982ff22ae23SHong Zhang {
983ff22ae23SHong Zhang   TS_RK *rk = (TS_RK*)ts->data;
984ff22ae23SHong Zhang 
985ff22ae23SHong Zhang   PetscFunctionBegin;
9860429704eSStefano Zampini   if (ns) *ns = rk->tableau->s;
987d2daff3dSHong Zhang   if (Y)   *Y = rk->Y;
988ff22ae23SHong Zhang   PetscFunctionReturn(0);
989ff22ae23SHong Zhang }
990ff22ae23SHong Zhang 
991b3a6b972SJed Brown static PetscErrorCode TSDestroy_RK(TS ts)
992b3a6b972SJed Brown {
993b3a6b972SJed Brown   PetscErrorCode ierr;
994b3a6b972SJed Brown 
995b3a6b972SJed Brown   PetscFunctionBegin;
996b3a6b972SJed Brown   ierr = TSReset_RK(ts);CHKERRQ(ierr);
997b3a6b972SJed Brown   if (ts->dm) {
998b3a6b972SJed Brown     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
999b3a6b972SJed Brown     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
1000b3a6b972SJed Brown   }
1001b3a6b972SJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
1002b3a6b972SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr);
1003b3a6b972SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr);
10040fe4d17eSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",NULL);CHKERRQ(ierr);
10050fe4d17eSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",NULL);CHKERRQ(ierr);
1006b3a6b972SJed Brown   PetscFunctionReturn(0);
1007b3a6b972SJed Brown }
1008ff22ae23SHong Zhang 
1009*21052c0cSHong Zhang /*@C
1010*21052c0cSHong Zhang   TSRKSetMultirate - Use the interpolation-based multirate RK method
1011*21052c0cSHong Zhang 
1012*21052c0cSHong Zhang   Logically collective
1013*21052c0cSHong Zhang 
1014*21052c0cSHong Zhang   Input Parameter:
1015*21052c0cSHong Zhang +  ts - timestepping context
1016*21052c0cSHong Zhang -  use_multirate - PETSC_TRUE enables the multirate RK method, sets the basic method to be RK2A and sets the ratio between slow stepsize and fast stepsize to be 2
1017*21052c0cSHong Zhang 
1018*21052c0cSHong Zhang   Options Database:
1019*21052c0cSHong Zhang .   -ts_rk_multirate - <true,false>
1020*21052c0cSHong Zhang 
1021*21052c0cSHong Zhang   Notes:
1022*21052c0cSHong Zhang   The multirate method requires interpolation. The default interpolation works for 1st- and 2nd- order RK, but not for high-order RKs except TSRK5DP which comes with the interpolation coeffcients (binterp).
1023*21052c0cSHong Zhang 
1024*21052c0cSHong Zhang   Level: intermediate
1025*21052c0cSHong Zhang 
1026*21052c0cSHong Zhang .seealso: TSRKGetMultirate()
1027*21052c0cSHong Zhang @*/
1028*21052c0cSHong Zhang PetscErrorCode TSRKSetMultirate(TS ts,PetscBool use_multirate)
1029*21052c0cSHong Zhang {
1030*21052c0cSHong Zhang   PetscTryMethod(ts,"TSRKSetMultirate_C",(TS,PetscBool),(ts,use_multirate));
1031*21052c0cSHong Zhang   PetscFunctionReturn(0);
1032*21052c0cSHong Zhang }
1033*21052c0cSHong Zhang 
1034*21052c0cSHong Zhang /*@C
1035*21052c0cSHong Zhang   TSRKGetMultirate - Gets whether to Use the interpolation-based multirate RK method
1036*21052c0cSHong Zhang 
1037*21052c0cSHong Zhang   Not collective
1038*21052c0cSHong Zhang 
1039*21052c0cSHong Zhang   Input Parameter:
1040*21052c0cSHong Zhang .  ts - timestepping context
1041*21052c0cSHong Zhang 
1042*21052c0cSHong Zhang   Output Parameter:
1043*21052c0cSHong Zhang .  use_multirate - PETSC_TRUE if the multirate RK method is enabled, PETSC_FALSE otherwise
1044*21052c0cSHong Zhang 
1045*21052c0cSHong Zhang   Level: intermediate
1046*21052c0cSHong Zhang 
1047*21052c0cSHong Zhang .seealso: TSRKSetMultirate()
1048*21052c0cSHong Zhang @*/
1049*21052c0cSHong Zhang PetscErrorCode TSRKGetMultirate(TS ts,PetscBool *use_multirate)
1050*21052c0cSHong Zhang {
1051*21052c0cSHong Zhang   PetscUseMethod(ts,"TSRKGetMultirate_C",(TS,PetscBool*),(ts,use_multirate));
1052*21052c0cSHong Zhang   PetscFunctionReturn(0);
1053*21052c0cSHong Zhang }
1054*21052c0cSHong Zhang 
1055ebe3b25bSBarry Smith /*MC
1056f68a32c8SEmil Constantinescu       TSRK - ODE and DAE solver using Runge-Kutta schemes
1057ebe3b25bSBarry Smith 
1058f68a32c8SEmil Constantinescu   The user should provide the right hand side of the equation
1059f68a32c8SEmil Constantinescu   using TSSetRHSFunction().
1060ebe3b25bSBarry Smith 
1061f68a32c8SEmil Constantinescu   Notes:
106298164e13SLisandro Dalcin   The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type
1063ebe3b25bSBarry Smith 
1064d41222bbSBarry Smith   Level: beginner
1065d41222bbSBarry Smith 
1066f68a32c8SEmil Constantinescu .seealso:  TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3,
10670fe4d17eSHong Zhang            TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister(), TSRKSetMultirate(), TSRKGetMultirate()
1068ebe3b25bSBarry Smith 
1069ebe3b25bSBarry Smith M*/
10708cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
1071e4dd521cSBarry Smith {
10721566a47fSLisandro Dalcin   TS_RK          *rk;
1073dfbe8321SBarry Smith   PetscErrorCode ierr;
1074e4dd521cSBarry Smith 
1075e4dd521cSBarry Smith   PetscFunctionBegin;
1076f68a32c8SEmil Constantinescu   ierr = TSRKInitializePackage();CHKERRQ(ierr);
1077f68a32c8SEmil Constantinescu 
1078f68a32c8SEmil Constantinescu   ts->ops->reset          = TSReset_RK;
10795f70b478SJed Brown   ts->ops->destroy        = TSDestroy_RK;
10805f70b478SJed Brown   ts->ops->view           = TSView_RK;
1081f68a32c8SEmil Constantinescu   ts->ops->load           = TSLoad_RK;
1082f68a32c8SEmil Constantinescu   ts->ops->setup          = TSSetUp_RK;
108342f2b339SBarry Smith   ts->ops->adjointsetup   = TSAdjointSetUp_RK;
1084f68a32c8SEmil Constantinescu   ts->ops->interpolate    = TSInterpolate_RK;
10852c9cb062Sluzhanghpp   ts->ops->step           = TSStep_RK;
1086f68a32c8SEmil Constantinescu   ts->ops->evaluatestep   = TSEvaluateStep_RK;
1087fecfb714SLisandro Dalcin   ts->ops->rollback       = TSRollBack_RK;
1088f68a32c8SEmil Constantinescu   ts->ops->setfromoptions = TSSetFromOptions_RK;
1089ff22ae23SHong Zhang   ts->ops->getstages      = TSGetStages_RK;
109042f2b339SBarry Smith   ts->ops->adjointstep    = TSAdjointStep_RK;
1091e4dd521cSBarry Smith 
10920f7a1166SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_RK;
10930f7a1166SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_RK;
10940f7a1166SHong Zhang 
10951566a47fSLisandro Dalcin   ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr);
10961566a47fSLisandro Dalcin   ts->data = (void*)rk;
1097e4dd521cSBarry Smith 
1098f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr);
1099f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr);
1100*21052c0cSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",TSRKSetMultirate_RK);CHKERRQ(ierr);
1101*21052c0cSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",TSRKGetMultirate_RK);CHKERRQ(ierr);
1102be5899b3SLisandro Dalcin 
1103be5899b3SLisandro Dalcin   ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr);
110490b67129SHong Zhang   rk->dtratio = 1;
1105e4dd521cSBarry Smith   PetscFunctionReturn(0);
1106e4dd521cSBarry Smith }
1107