xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision 98921bda46e76d7aaed9e0138c5ff9d0ce93f355)
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>
1421052c0cSHong 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 
20ab8985f5SHong Zhang static RKTableauLink RKTableauList;
21ab8985f5SHong Zhang 
22f68a32c8SEmil Constantinescu /*MC
23287c2655SBarry Smith      TSRK1FE - First order forward Euler scheme.
24262ff7bbSBarry Smith 
25f68a32c8SEmil Constantinescu      This method has one stage.
26f68a32c8SEmil Constantinescu 
27287c2655SBarry Smith      Options database:
289bd3e852SBarry Smith .     -ts_rk_type 1fe
29287c2655SBarry Smith 
30f68a32c8SEmil Constantinescu      Level: advanced
31f68a32c8SEmil Constantinescu 
32287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
33f68a32c8SEmil Constantinescu M*/
34f68a32c8SEmil Constantinescu /*MC
35e7685601SHong Zhang      TSRK2A - Second order RK scheme (Heun's method).
36f68a32c8SEmil Constantinescu 
37f68a32c8SEmil Constantinescu      This method has two stages.
38f68a32c8SEmil Constantinescu 
39287c2655SBarry Smith      Options database:
409bd3e852SBarry Smith .     -ts_rk_type 2a
41287c2655SBarry Smith 
42f68a32c8SEmil Constantinescu      Level: advanced
43f68a32c8SEmil Constantinescu 
44287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
45f68a32c8SEmil Constantinescu M*/
46f68a32c8SEmil Constantinescu /*MC
47e7685601SHong Zhang      TSRK2B - Second order RK scheme (the midpoint method).
48e7685601SHong Zhang 
49e7685601SHong Zhang      This method has two stages.
50e7685601SHong Zhang 
51e7685601SHong Zhang      Options database:
52e7685601SHong Zhang .     -ts_rk_type 2b
53e7685601SHong Zhang 
54e7685601SHong Zhang      Level: advanced
55e7685601SHong Zhang 
56e7685601SHong Zhang .seealso: TSRK, TSRKType, TSRKSetType()
57e7685601SHong Zhang M*/
58e7685601SHong Zhang /*MC
59f68a32c8SEmil Constantinescu      TSRK3 - Third order RK scheme.
60f68a32c8SEmil Constantinescu 
61f68a32c8SEmil Constantinescu      This method has three stages.
62f68a32c8SEmil Constantinescu 
63287c2655SBarry Smith      Options database:
649bd3e852SBarry Smith .     -ts_rk_type 3
65287c2655SBarry Smith 
66f68a32c8SEmil Constantinescu      Level: advanced
67f68a32c8SEmil Constantinescu 
68287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
69f68a32c8SEmil Constantinescu M*/
70f68a32c8SEmil Constantinescu /*MC
712109b73fSDebojyoti Ghosh      TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method.
722109b73fSDebojyoti Ghosh 
73d1905564SLisandro Dalcin      This method has four stages with the First Same As Last (FSAL) property.
742109b73fSDebojyoti Ghosh 
75287c2655SBarry Smith      Options database:
769bd3e852SBarry Smith .     -ts_rk_type 3bs
77287c2655SBarry Smith 
782109b73fSDebojyoti Ghosh      Level: advanced
792109b73fSDebojyoti Ghosh 
8098164e13SLisandro Dalcin      References: https://doi.org/10.1016/0893-9659(89)90079-7
81d1905564SLisandro Dalcin 
82287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
832109b73fSDebojyoti Ghosh M*/
842109b73fSDebojyoti Ghosh /*MC
85f68a32c8SEmil Constantinescu      TSRK4 - Fourth order RK scheme.
86f68a32c8SEmil Constantinescu 
872e698f8bSDebojyoti Ghosh      This is the classical Runge-Kutta method with four stages.
88f68a32c8SEmil Constantinescu 
89287c2655SBarry Smith      Options database:
909bd3e852SBarry Smith .     -ts_rk_type 4
91287c2655SBarry Smith 
92f68a32c8SEmil Constantinescu      Level: advanced
93f68a32c8SEmil Constantinescu 
94287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
95f68a32c8SEmil Constantinescu M*/
96f68a32c8SEmil Constantinescu /*MC
972e698f8bSDebojyoti Ghosh      TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method.
98f68a32c8SEmil Constantinescu 
99f68a32c8SEmil Constantinescu      This method has six stages.
100f68a32c8SEmil Constantinescu 
101287c2655SBarry Smith      Options database:
1029bd3e852SBarry Smith .     -ts_rk_type 5f
103287c2655SBarry Smith 
104f68a32c8SEmil Constantinescu      Level: advanced
105f68a32c8SEmil Constantinescu 
106287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
107f68a32c8SEmil Constantinescu M*/
1082109b73fSDebojyoti Ghosh /*MC
1092e698f8bSDebojyoti Ghosh      TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method.
1102109b73fSDebojyoti Ghosh 
111d1905564SLisandro Dalcin      This method has seven stages with the First Same As Last (FSAL) property.
1122109b73fSDebojyoti Ghosh 
113287c2655SBarry Smith      Options database:
1149bd3e852SBarry Smith .     -ts_rk_type 5dp
115287c2655SBarry Smith 
1162109b73fSDebojyoti Ghosh      Level: advanced
1172109b73fSDebojyoti Ghosh 
11898164e13SLisandro Dalcin      References: https://doi.org/10.1016/0771-050X(80)90013-3
119d1905564SLisandro Dalcin 
120287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
1212109b73fSDebojyoti Ghosh M*/
12205e23783SLisandro Dalcin /*MC
12305e23783SLisandro Dalcin      TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method.
12405e23783SLisandro Dalcin 
12505e23783SLisandro Dalcin      This method has eight stages with the First Same As Last (FSAL) property.
12605e23783SLisandro Dalcin 
127287c2655SBarry Smith      Options database:
1289bd3e852SBarry Smith .     -ts_rk_type 5bs
129287c2655SBarry Smith 
13005e23783SLisandro Dalcin      Level: advanced
13105e23783SLisandro Dalcin 
13298164e13SLisandro Dalcin      References: https://doi.org/10.1016/0898-1221(96)00141-1
13305e23783SLisandro Dalcin 
134287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
13505e23783SLisandro Dalcin M*/
13637acaa02SHendrik Ranocha /*MC
13737acaa02SHendrik Ranocha      TSRK6VR - Sixth order robust Verner RK scheme with fifth order embedded method.
13837acaa02SHendrik Ranocha 
13937acaa02SHendrik Ranocha      This method has nine stages with the First Same As Last (FSAL) property.
14037acaa02SHendrik Ranocha 
14137acaa02SHendrik Ranocha      Options database:
14237acaa02SHendrik Ranocha .     -ts_rk_type 6vr
14337acaa02SHendrik Ranocha 
14437acaa02SHendrik Ranocha      Level: advanced
14537acaa02SHendrik Ranocha 
14637acaa02SHendrik Ranocha      References: http://people.math.sfu.ca/~jverner/RKV65.IIIXb.Robust.00010102836.081204.CoeffsOnlyRAT
14737acaa02SHendrik Ranocha 
14837acaa02SHendrik Ranocha .seealso: TSRK, TSRKType, TSRKSetType()
14937acaa02SHendrik Ranocha M*/
15037acaa02SHendrik Ranocha /*MC
15137acaa02SHendrik Ranocha      TSRK7VR - Seventh order robust Verner RK scheme with sixth order embedded method.
15237acaa02SHendrik Ranocha 
1538f3d7ee2SCarsten Uphoff      This method has ten stages.
15437acaa02SHendrik Ranocha 
15537acaa02SHendrik Ranocha      Options database:
15637acaa02SHendrik Ranocha .     -ts_rk_type 7vr
15737acaa02SHendrik Ranocha 
15837acaa02SHendrik Ranocha      Level: advanced
15937acaa02SHendrik Ranocha 
16037acaa02SHendrik Ranocha      References: http://people.math.sfu.ca/~jverner/RKV76.IIa.Robust.000027015646.081206.CoeffsOnlyRAT
16137acaa02SHendrik Ranocha 
16237acaa02SHendrik Ranocha .seealso: TSRK, TSRKType, TSRKSetType()
16337acaa02SHendrik Ranocha M*/
16437acaa02SHendrik Ranocha /*MC
16537acaa02SHendrik Ranocha      TSRK8VR - Eigth order robust Verner RK scheme with seventh order embedded method.
16637acaa02SHendrik Ranocha 
1678f3d7ee2SCarsten Uphoff      This method has thirteen stages.
16837acaa02SHendrik Ranocha 
16937acaa02SHendrik Ranocha      Options database:
17037acaa02SHendrik Ranocha .     -ts_rk_type 8vr
17137acaa02SHendrik Ranocha 
17237acaa02SHendrik Ranocha      Level: advanced
17337acaa02SHendrik Ranocha 
17437acaa02SHendrik Ranocha      References: http://people.math.sfu.ca/~jverner/RKV87.IIa.Robust.00000754677.081208.CoeffsOnlyRATandFLOAT
17537acaa02SHendrik Ranocha 
17637acaa02SHendrik Ranocha .seealso: TSRK, TSRKType, TSRKSetType()
17737acaa02SHendrik Ranocha M*/
178262ff7bbSBarry Smith 
179f68a32c8SEmil Constantinescu /*@C
180f68a32c8SEmil Constantinescu   TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK
181262ff7bbSBarry Smith 
182f68a32c8SEmil Constantinescu   Not Collective, but should be called by all processes which will need the schemes to be registered
183262ff7bbSBarry Smith 
184f68a32c8SEmil Constantinescu   Level: advanced
185262ff7bbSBarry Smith 
186f68a32c8SEmil Constantinescu .seealso:  TSRKRegisterDestroy()
187262ff7bbSBarry Smith @*/
188f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void)
189262ff7bbSBarry Smith {
1904ac538c5SBarry Smith   PetscErrorCode ierr;
191262ff7bbSBarry Smith 
192262ff7bbSBarry Smith   PetscFunctionBegin;
193f68a32c8SEmil Constantinescu   if (TSRKRegisterAllCalled) PetscFunctionReturn(0);
194f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_TRUE;
195e4dd521cSBarry Smith 
1964e82626cSLisandro Dalcin #define RC PetscRealConstant
197e4dd521cSBarry Smith   {
198f68a32c8SEmil Constantinescu     const PetscReal
1994e82626cSLisandro Dalcin       A[1][1] = {{0}},
2004e82626cSLisandro Dalcin       b[1]    = {RC(1.0)};
2013a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
202e4dd521cSBarry Smith   }
203db046809SBarry Smith   {
204f68a32c8SEmil Constantinescu     const PetscReal
2054e82626cSLisandro Dalcin       A[2][2]   = {{0,0},
2064e82626cSLisandro Dalcin                    {RC(1.0),0}},
2074e82626cSLisandro Dalcin       b[2]      =  {RC(0.5),RC(0.5)},
2084e82626cSLisandro Dalcin       bembed[2] =  {RC(1.0),0};
2093a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
210db046809SBarry Smith   }
211f68a32c8SEmil Constantinescu   {
212f68a32c8SEmil Constantinescu     const PetscReal
213e7685601SHong Zhang       A[2][2]   = {{0,0},
214e7685601SHong Zhang                    {RC(0.5),0}},
215e7685601SHong Zhang       b[2]      =  {0,RC(1.0)},
216e7685601SHong Zhang     ierr = TSRKRegister(TSRK2B,2,2,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
217e7685601SHong Zhang   }
218e7685601SHong Zhang   {
219e7685601SHong Zhang     const PetscReal
22017f6384fSBarry Smith       A[3][3] = {{0,0,0},
2214e82626cSLisandro Dalcin                  {RC(2.0)/RC(3.0),0,0},
2224e82626cSLisandro Dalcin                  {RC(-1.0)/RC(3.0),RC(1.0),0}},
2234e82626cSLisandro Dalcin       b[3]    =  {RC(0.25),RC(0.5),RC(0.25)};
2243a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
225fdefd5e5SDebojyoti Ghosh   }
226fdefd5e5SDebojyoti Ghosh   {
227fdefd5e5SDebojyoti Ghosh     const PetscReal
22817f6384fSBarry Smith       A[4][4]   = {{0,0,0,0},
2294e82626cSLisandro Dalcin                    {RC(1.0)/RC(2.0),0,0,0},
2304e82626cSLisandro Dalcin                    {0,RC(3.0)/RC(4.0),0,0},
2314e82626cSLisandro Dalcin                    {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}},
2324e82626cSLisandro Dalcin       b[4]      =  {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0},
2334e82626cSLisandro 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)};
2343a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
235db046809SBarry Smith   }
236f68a32c8SEmil Constantinescu   {
237f68a32c8SEmil Constantinescu     const PetscReal
238f68a32c8SEmil Constantinescu       A[4][4] = {{0,0,0,0},
2394e82626cSLisandro Dalcin                  {RC(0.5),0,0,0},
2404e82626cSLisandro Dalcin                  {0,RC(0.5),0,0},
2414e82626cSLisandro Dalcin                  {0,0,RC(1.0),0}},
2424e82626cSLisandro 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)};
2433a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
244db046809SBarry Smith   }
245f68a32c8SEmil Constantinescu   {
246f68a32c8SEmil Constantinescu     const PetscReal
24717f6384fSBarry Smith       A[6][6]   = {{0,0,0,0,0,0},
2484e82626cSLisandro Dalcin                    {RC(0.25),0,0,0,0,0},
2494e82626cSLisandro Dalcin                    {RC(3.0)/RC(32.0),RC(9.0)/RC(32.0),0,0,0,0},
2504e82626cSLisandro Dalcin                    {RC(1932.0)/RC(2197.0),RC(-7200.0)/RC(2197.0),RC(7296.0)/RC(2197.0),0,0,0},
2514e82626cSLisandro Dalcin                    {RC(439.0)/RC(216.0),RC(-8.0),RC(3680.0)/RC(513.0),RC(-845.0)/RC(4104.0),0,0},
2524e82626cSLisandro 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}},
2534e82626cSLisandro 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)},
2544e82626cSLisandro 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};
2553a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
256fdefd5e5SDebojyoti Ghosh   }
257fdefd5e5SDebojyoti Ghosh   {
258fdefd5e5SDebojyoti Ghosh     const PetscReal
25917f6384fSBarry Smith       A[7][7]       = {{0,0,0,0,0,0,0},
2604e82626cSLisandro Dalcin                        {RC(1.0)/RC(5.0),0,0,0,0,0,0},
2614e82626cSLisandro Dalcin                        {RC(3.0)/RC(40.0),RC(9.0)/RC(40.0),0,0,0,0,0},
2624e82626cSLisandro Dalcin                        {RC(44.0)/RC(45.0),RC(-56.0)/RC(15.0),RC(32.0)/RC(9.0),0,0,0,0},
2634e82626cSLisandro 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},
2644e82626cSLisandro 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},
2654e82626cSLisandro 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}},
2664e82626cSLisandro 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},
267a685a763Sluzhanghpp       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)},
268a685a763Sluzhanghpp       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)},
269a685a763Sluzhanghpp                        {0,0,0,0,0},
270a685a763Sluzhanghpp                        {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)},
271a685a763Sluzhanghpp                        {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)},
272a685a763Sluzhanghpp                        {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)},
273a685a763Sluzhanghpp                        {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)},
274a685a763Sluzhanghpp                        {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)}};
275a685a763Sluzhanghpp       ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,5,binterp[0]);CHKERRQ(ierr);
276f68a32c8SEmil Constantinescu   }
27705e23783SLisandro Dalcin   {
27805e23783SLisandro Dalcin     const PetscReal
27917f6384fSBarry Smith       A[8][8]   = {{0,0,0,0,0,0,0,0},
2804e82626cSLisandro Dalcin                    {RC(1.0)/RC(6.0),0,0,0,0,0,0,0},
2814e82626cSLisandro Dalcin                    {RC(2.0)/RC(27.0),RC(4.0)/RC(27.0),0,0,0,0,0,0},
2824e82626cSLisandro Dalcin                    {RC(183.0)/RC(1372.0),RC(-162.0)/RC(343.0),RC(1053.0)/RC(1372.0),0,0,0,0,0},
2834e82626cSLisandro 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},
2844e82626cSLisandro 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},
2854e82626cSLisandro 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},
2864e82626cSLisandro 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}},
2874e82626cSLisandro 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},
2884e82626cSLisandro 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)};
28905e23783SLisandro Dalcin     ierr = TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
29005e23783SLisandro Dalcin   }
29137acaa02SHendrik Ranocha   {
29237acaa02SHendrik Ranocha     const PetscReal
29337acaa02SHendrik Ranocha       A[9][9]   = {{0,0,0,0,0,0,0,0,0},
29463fe90e8SHendrik Ranocha                    {RC(1.8000000000000000000000000000000000000000e-01),0,0,0,0,0,0,0,0},
29563fe90e8SHendrik Ranocha                    {RC(8.9506172839506172839506172839506172839506e-02),RC(7.7160493827160493827160493827160493827160e-02),0,0,0,0,0,0,0},
29663fe90e8SHendrik Ranocha                    {RC(6.2500000000000000000000000000000000000000e-02),0,RC(1.8750000000000000000000000000000000000000e-01),0,0,0,0,0,0},
29763fe90e8SHendrik Ranocha                    {RC(3.1651600000000000000000000000000000000000e-01),0,RC(-1.0449480000000000000000000000000000000000e+00),RC(1.2584320000000000000000000000000000000000e+00),0,0,0,0,0},
29863fe90e8SHendrik Ranocha                    {RC(2.7232612736485626257225065566674305502508e-01),0,RC(-8.2513360323886639676113360323886639676113e-01),RC(1.0480917678812415654520917678812415654521e+00),RC(1.0471570799276856873679117969088177628396e-01),0,0,0,0},
29963fe90e8SHendrik Ranocha                    {RC(-1.6699418599716514314329607278961797333198e-01),0,RC(6.3170850202429149797570850202429149797571e-01),RC(1.7461044552773876082146758838488161796432e-01),RC(-1.0665356459086066122525194734018680677781e+00),RC(1.2272108843537414965986394557823129251701e+00),0,0,0},
30063fe90e8SHendrik Ranocha                    {RC(3.6423751686909581646423751686909581646424e-01),0,RC(-2.0404858299595141700404858299595141700405e-01),RC(-3.4883737816068643136312309244640071707741e-01),RC(3.2619323032856867443333608747142581729048e+00),RC(-2.7551020408163265306122448979591836734694e+00),RC(6.8181818181818181818181818181818181818182e-01),0,0},
30163fe90e8SHendrik Ranocha                    {RC(7.6388888888888888888888888888888888888889e-02),0,0,RC(3.6940836940836940836940836940836940836941e-01),0,RC(2.4801587301587301587301587301587301587302e-01),RC(2.3674242424242424242424242424242424242424e-01),RC(6.9444444444444444444444444444444444444444e-02),0}},
30263fe90e8SHendrik Ranocha       b[9]      =  {RC(7.6388888888888888888888888888888888888889e-02),0,0,RC(3.6940836940836940836940836940836940836941e-01),0,RC(2.4801587301587301587301587301587301587302e-01),RC(2.3674242424242424242424242424242424242424e-01),RC(6.9444444444444444444444444444444444444444e-02),0},
30363fe90e8SHendrik Ranocha       bembed[9] =  {RC(5.8700209643605870020964360587002096436059e-02),0,0,RC(4.8072562358276643990929705215419501133787e-01),RC(-8.5341242076919085578832094861228313083563e-01),RC(1.2046485260770975056689342403628117913832e+00),0,RC(-5.9242373072160306202859394348756050883710e-02),RC(1.6858043453788134639198468985703028256220e-01)};
30437acaa02SHendrik Ranocha     ierr = TSRKRegister(TSRK6VR,6,9,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
30537acaa02SHendrik Ranocha   }
30637acaa02SHendrik Ranocha   {
30737acaa02SHendrik Ranocha     const PetscReal
30837acaa02SHendrik Ranocha       A[10][10]  = {{0,0,0,0,0,0,0,0,0,0},
30963fe90e8SHendrik Ranocha                     {RC(5.0000000000000000000000000000000000000000e-03),0,0,0,0,0,0,0,0,0},
31063fe90e8SHendrik Ranocha                     {RC(-1.0767901234567901234567901234567901234568e+00),RC(1.1856790123456790123456790123456790123457e+00),0,0,0,0,0,0,0,0},
31163fe90e8SHendrik Ranocha                     {RC(4.0833333333333333333333333333333333333333e-02),0,RC(1.2250000000000000000000000000000000000000e-01),0,0,0,0,0,0,0},
31263fe90e8SHendrik Ranocha                     {RC(6.3607142857142857142857142857142857142857e-01),0,RC(-2.4444642857142857142857142857142857142857e+00),RC(2.2633928571428571428571428571428571428571e+00),0,0,0,0,0,0},
31363fe90e8SHendrik Ranocha                     {RC(-2.5351211079349245229256383554660215487207e+00),0,RC(1.0299374654449267920438514460756024913612e+01),RC(-7.9513032885990579949493217458266876536482e+00),RC(7.9301148923100592201226014271115261823800e-01),0,0,0,0,0},
31463fe90e8SHendrik Ranocha                     {RC(1.0018765812524632961969196583094999808207e+00),0,RC(-4.1665712824423798331313938005470971453189e+00),RC(3.8343432929128642412552665218251378665197e+00),RC(-5.0233333560710847547464330228611765612403e-01),RC(6.6768474388416077115385092269857695410259e-01),0,0,0,0},
31563fe90e8SHendrik Ranocha                     {RC(2.7255018354630767130333963819175005717348e+01),0,RC(-4.2004617278410638355318645443909295369611e+01),RC(-1.0535713126619489917921081600546526103722e+01),RC(8.0495536711411937147983652158926826634202e+01),RC(-6.7343882271790513468549075963212975640927e+01),RC(1.3048657610777937463471187029566964762710e+01),0,0,0},
31663fe90e8SHendrik Ranocha                     {RC(-3.0397378057114965146943658658755763226883e+00),0,RC(1.0138161410329801111857946190709700150441e+01),RC(-6.4293056748647215721462825629555298064437e+00),RC(-1.5864371483408276587115312853798610579467e+00),RC(1.8921781841968424410864308909131353365021e+00),RC(1.9699335407608869061292360163336442838006e-02),RC(5.4416989827933235465102724247952572977903e-03),0,0},
31763fe90e8SHendrik Ranocha                     {RC(-1.4449518916777735137351003179355712360517e+00),0,RC(8.0318913859955919224117033223019560435041e+00),RC(-7.5831741663401346820798883023671588604984e+00),RC(3.5816169353190074211247685442452878696855e+00),RC(-2.4369722632199529411183809065693752383733e+00),RC(8.5158999992326179339689766032486142173390e-01),0,0,0}},
31863fe90e8SHendrik Ranocha       b[10]      =  {RC(4.7425837833706756083569172717574534698932e-02),0,0,RC(2.5622361659370562659961727458274623448160e-01),RC(2.6951376833074206619473817258075952886764e-01),RC(1.2686622409092782845989138364739173247882e-01),RC(2.4887225942060071622046449427647492767292e-01),RC(3.0744837408200631335304388479099184768645e-03),RC(4.8023809989496943308189063347143123323209e-02),0},
31963fe90e8SHendrik Ranocha       bembed[10] =  {RC(4.7485247699299631037531273805727961552268e-02),0,0,RC(2.5599412588690633297154918245905393870497e-01),RC(2.7058478081067688722530891099268135732387e-01),RC(1.2505618684425992913638822323746917920448e-01),RC(2.5204468723743860507184043820197442562182e-01),0,0,RC(4.8834971521418614557381971303093137592592e-02)};
32037acaa02SHendrik Ranocha     ierr = TSRKRegister(TSRK7VR,7,10,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
32137acaa02SHendrik Ranocha   }
32237acaa02SHendrik Ranocha   {
32337acaa02SHendrik Ranocha     const PetscReal
32437acaa02SHendrik Ranocha       A[13][13]  = {{0,0,0,0,0,0,0,0,0,0,0,0,0},
32563fe90e8SHendrik Ranocha                     {RC(2.5000000000000000000000000000000000000000e-01),0,0,0,0,0,0,0,0,0,0,0,0},
32663fe90e8SHendrik Ranocha                     {RC(8.7400846504915232052686327594877411977046e-02),RC(2.5487604938654321753087950620345685135815e-02),0,0,0,0,0,0,0,0,0,0,0},
32763fe90e8SHendrik Ranocha                     {RC(4.2333169291338582677165354330708661417323e-02),0,RC(1.2699950787401574803149606299212598425197e-01),0,0,0,0,0,0,0,0,0,0},
32863fe90e8SHendrik Ranocha                     {RC(4.2609505888742261494881445237572274090942e-01),0,RC(-1.5987952846591523265427733230657181117089e+00),RC(1.5967002257717297115939588706899953707994e+00),0,0,0,0,0,0,0,0,0},
32963fe90e8SHendrik Ranocha                     {RC(5.0719337296713929515090618138513639239329e-02),0,0,RC(2.5433377264600407582754714408877778031369e-01),RC(2.0394689005728199465736223777270858044698e-01),0,0,0,0,0,0,0,0},
33063fe90e8SHendrik Ranocha                     {RC(-2.9000374717523110970388379285425896124091e-01),0,0,RC(1.3441873910260789889438681109414337003184e+00),RC(-2.8647779433614427309611103827036562829470e+00),RC(2.6775942995105948517211260646164815438695e+00),0,0,0,0,0,0,0},
33163fe90e8SHendrik Ranocha                     {RC(9.8535011337993546469740402980727014284756e-02),0,0,0,RC(2.2192680630751384842024036498197387903583e-01),RC(-1.8140622911806994312690338288073952457474e-01),RC(1.0944411472562548236922614918038631254153e-02),0,0,0,0,0,0},
33263fe90e8SHendrik Ranocha                     {RC(3.8711052545731144679444618165166373405645e-01),0,0,RC(-1.4424454974855277571256745553077927767173e+00),RC(2.9053981890699509317691346449233848441744e+00),RC(-1.8537710696301059290843332675811978025183e+00),RC(1.4003648098728154269497325109771241479223e-01),RC(5.7273940811495816575746774624447706488753e-01),0,0,0,0,0},
33363fe90e8SHendrik Ranocha                     {RC(-1.6124403444439308100630016197913480595436e-01),0,0,RC(-1.7339602957358984083578404473962567894901e-01),RC(-1.3012892814065147406016812745172492529744e+00),RC(1.1379503751738617308558792131431003472124e+00),RC(-3.1747649663966880106923521138043024698980e-02),RC(9.3351293824933666439811064486056884856590e-01),RC(-8.3786318334733852703300855629616433201504e-02),0,0,0,0},
33463fe90e8SHendrik Ranocha                     {RC(-1.9199444881589533281510804651483576073142e-02),0,0,RC(2.7330857265264284907942326254016124275617e-01),RC(-6.7534973206944372919691611210942380856240e-01),RC(3.4151849813846016071738489974728382711981e-01),RC(-6.7950064803375772478920516198524629391910e-02),RC(9.6591752247623878884265586491216376509746e-02),RC(1.3253082511182101180721038466545389951226e-01),RC(3.6854959360386113446906329951531666812946e-01),0,0,0},
33563fe90e8SHendrik Ranocha                     {RC(6.0918774036452898676888412111588817784584e-01),0,0,RC(-2.2725690858980016768999800931413088399719e+00),RC(4.7578983426940290068155255881914785497547e+00),RC(-5.5161067066927584824294689667844248244842e+00),RC(2.9005963696801192709095818565946174378180e-01),RC(5.6914239633590368229109858454801849145630e-01),RC(7.9267957603321670271339916205893327579951e-01),RC(1.5473720453288822894126190771849898232047e-01),RC(1.6149708956621816247083215106334544434974e+00),0,0},
33663fe90e8SHendrik Ranocha                     {RC(8.8735762208534719663211694051981022704884e-01),0,0,RC(-2.9754597821085367558513632804709301581977e+00),RC(5.6007170094881630597990392548350098923829e+00),RC(-5.9156074505366744680014930189941657351840e+00),RC(2.2029689156134927016879142540807638331238e-01),RC(1.0155097824462216666143271340902996997549e-01),RC(1.1514345647386055909780397752125850553556e+00),RC(1.9297101665271239396134361900805843653065e+00),0,0,0}},
33763fe90e8SHendrik Ranocha       b[13]      =  {RC(4.4729564666695714203015840429049382466467e-02),0,0,0,0,RC(1.5691033527708199813368698010726645409175e-01),RC(1.8460973408151637740702451873526277892035e-01),RC(2.2516380602086991042479419400350721970920e-01),RC(1.4794615651970234687005179885449141753736e-01),RC(7.6055542444955825269798361910336491012732e-02),RC(1.2277290235018619610824346315921437388535e-01),RC(4.1811958638991631583384842800871882376786e-02),0},
33863fe90e8SHendrik Ranocha       bembed[13] =  {RC(4.5847111400495925878664730122010282095875e-02),0,0,0,0,RC(2.6231891404152387437443356584845803392392e-01),RC(1.9169372337852611904485738635688429008025e-01),RC(2.1709172327902618330978407422906448568196e-01),RC(1.2738189624833706796803169450656737867900e-01),RC(1.1510530385365326258240515750043192148894e-01),0,0,RC(4.0561327798437566841823391436583608050053e-02)};
33937acaa02SHendrik Ranocha     ierr = TSRKRegister(TSRK8VR,8,13,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
34037acaa02SHendrik Ranocha   }
3414e82626cSLisandro Dalcin #undef RC
342db046809SBarry Smith   PetscFunctionReturn(0);
343db046809SBarry Smith }
344db046809SBarry Smith 
345f68a32c8SEmil Constantinescu /*@C
346f68a32c8SEmil Constantinescu    TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister().
347f68a32c8SEmil Constantinescu 
348f68a32c8SEmil Constantinescu    Not Collective
349f68a32c8SEmil Constantinescu 
350f68a32c8SEmil Constantinescu    Level: advanced
351f68a32c8SEmil Constantinescu 
352f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll()
353f68a32c8SEmil Constantinescu @*/
354f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void)
355e4dd521cSBarry Smith {
356f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
357f68a32c8SEmil Constantinescu   RKTableauLink  link;
358f68a32c8SEmil Constantinescu 
359f68a32c8SEmil Constantinescu   PetscFunctionBegin;
360f68a32c8SEmil Constantinescu   while ((link = RKTableauList)) {
361f68a32c8SEmil Constantinescu     RKTableau t = &link->tab;
362f68a32c8SEmil Constantinescu     RKTableauList = link->next;
363f68a32c8SEmil Constantinescu     ierr = PetscFree3(t->A,t->b,t->c);CHKERRQ(ierr);
364f68a32c8SEmil Constantinescu     ierr = PetscFree(t->bembed);CHKERRQ(ierr);
365f68a32c8SEmil Constantinescu     ierr = PetscFree(t->binterp);CHKERRQ(ierr);
366f68a32c8SEmil Constantinescu     ierr = PetscFree(t->name);CHKERRQ(ierr);
367f68a32c8SEmil Constantinescu     ierr = PetscFree(link);CHKERRQ(ierr);
368f68a32c8SEmil Constantinescu   }
369f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_FALSE;
370f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
371f68a32c8SEmil Constantinescu }
372f68a32c8SEmil Constantinescu 
373f68a32c8SEmil Constantinescu /*@C
374f68a32c8SEmil Constantinescu   TSRKInitializePackage - This function initializes everything in the TSRK package. It is called
3758a690491SBarry Smith   from TSInitializePackage().
376f68a32c8SEmil Constantinescu 
377f68a32c8SEmil Constantinescu   Level: developer
378f68a32c8SEmil Constantinescu 
379f68a32c8SEmil Constantinescu .seealso: PetscInitialize()
380f68a32c8SEmil Constantinescu @*/
381f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void)
382f68a32c8SEmil Constantinescu {
383186e87acSLisandro Dalcin   PetscErrorCode ierr;
384e4dd521cSBarry Smith 
385e4dd521cSBarry Smith   PetscFunctionBegin;
386f68a32c8SEmil Constantinescu   if (TSRKPackageInitialized) PetscFunctionReturn(0);
387f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_TRUE;
388f68a32c8SEmil Constantinescu   ierr = TSRKRegisterAll();CHKERRQ(ierr);
389f68a32c8SEmil Constantinescu   ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr);
390f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
391f68a32c8SEmil Constantinescu }
392186e87acSLisandro Dalcin 
393f68a32c8SEmil Constantinescu /*@C
394f68a32c8SEmil Constantinescu   TSRKFinalizePackage - This function destroys everything in the TSRK package. It is
395f68a32c8SEmil Constantinescu   called from PetscFinalize().
396186e87acSLisandro Dalcin 
397f68a32c8SEmil Constantinescu   Level: developer
398f68a32c8SEmil Constantinescu 
399f68a32c8SEmil Constantinescu .seealso: PetscFinalize()
400f68a32c8SEmil Constantinescu @*/
401f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void)
402f68a32c8SEmil Constantinescu {
403f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
404f68a32c8SEmil Constantinescu 
405f68a32c8SEmil Constantinescu   PetscFunctionBegin;
406f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_FALSE;
407f68a32c8SEmil Constantinescu   ierr = TSRKRegisterDestroy();CHKERRQ(ierr);
408f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
409f68a32c8SEmil Constantinescu }
410f68a32c8SEmil Constantinescu 
411f68a32c8SEmil Constantinescu /*@C
412f68a32c8SEmil Constantinescu    TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
413f68a32c8SEmil Constantinescu 
414f68a32c8SEmil Constantinescu    Not Collective, but the same schemes should be registered on all processes on which they will be used
415f68a32c8SEmil Constantinescu 
416f68a32c8SEmil Constantinescu    Input Parameters:
417f68a32c8SEmil Constantinescu +  name - identifier for method
418f68a32c8SEmil Constantinescu .  order - approximation order of method
419f68a32c8SEmil Constantinescu .  s - number of stages, this is the dimension of the matrices below
420f68a32c8SEmil Constantinescu .  A - stage coefficients (dimension s*s, row-major)
421f68a32c8SEmil Constantinescu .  b - step completion table (dimension s; NULL to use last row of A)
422f68a32c8SEmil Constantinescu .  c - abscissa (dimension s; NULL to use row sums of A)
423f68a32c8SEmil Constantinescu .  bembed - completion table for embedded method (dimension s; NULL if not available)
4243a8a9803SLisandro Dalcin .  p - Order of the interpolation scheme, equal to the number of columns of binterp
4253a8a9803SLisandro Dalcin -  binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1)
426f68a32c8SEmil Constantinescu 
427f68a32c8SEmil Constantinescu    Notes:
428f68a32c8SEmil Constantinescu    Several RK methods are provided, this function is only needed to create new methods.
429f68a32c8SEmil Constantinescu 
430f68a32c8SEmil Constantinescu    Level: advanced
431f68a32c8SEmil Constantinescu 
432f68a32c8SEmil Constantinescu .seealso: TSRK
433f68a32c8SEmil Constantinescu @*/
434f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s,
435f68a32c8SEmil Constantinescu                             const PetscReal A[],const PetscReal b[],const PetscReal c[],
4363a8a9803SLisandro Dalcin                             const PetscReal bembed[],PetscInt p,const PetscReal binterp[])
437f68a32c8SEmil Constantinescu {
438f68a32c8SEmil Constantinescu   PetscErrorCode  ierr;
439f68a32c8SEmil Constantinescu   RKTableauLink   link;
440f68a32c8SEmil Constantinescu   RKTableau       t;
441f68a32c8SEmil Constantinescu   PetscInt        i,j;
442f68a32c8SEmil Constantinescu 
443f68a32c8SEmil Constantinescu   PetscFunctionBegin;
4443a8a9803SLisandro Dalcin   PetscValidCharPointer(name,1);
4453a8a9803SLisandro Dalcin   PetscValidRealPointer(A,4);
4463a8a9803SLisandro Dalcin   if (b) PetscValidRealPointer(b,5);
4473a8a9803SLisandro Dalcin   if (c) PetscValidRealPointer(c,6);
4483a8a9803SLisandro Dalcin   if (bembed) PetscValidRealPointer(bembed,7);
4493a8a9803SLisandro Dalcin   if (binterp || p > 1) PetscValidRealPointer(binterp,9);
4503a8a9803SLisandro Dalcin 
4511d36bdfdSBarry Smith   ierr = TSRKInitializePackage();CHKERRQ(ierr);
45295dccacaSBarry Smith   ierr = PetscNew(&link);CHKERRQ(ierr);
453f68a32c8SEmil Constantinescu   t = &link->tab;
4543a8a9803SLisandro Dalcin 
455f68a32c8SEmil Constantinescu   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
456f68a32c8SEmil Constantinescu   t->order = order;
457f68a32c8SEmil Constantinescu   t->s = s;
458dcca6d9dSJed Brown   ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr);
459580bdb30SBarry Smith   ierr = PetscArraycpy(t->A,A,s*s);CHKERRQ(ierr);
460580bdb30SBarry Smith   if (b)  { ierr = PetscArraycpy(t->b,b,s);CHKERRQ(ierr); }
461f68a32c8SEmil Constantinescu   else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i];
462580bdb30SBarry Smith   if (c)  { ierr = PetscArraycpy(t->c,c,s);CHKERRQ(ierr); }
463f68a32c8SEmil 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];
464d760c35bSDebojyoti Ghosh   t->FSAL = PETSC_TRUE;
465d760c35bSDebojyoti Ghosh   for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE;
4663a8a9803SLisandro Dalcin 
467f68a32c8SEmil Constantinescu   if (bembed) {
468785e854fSJed Brown     ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr);
469580bdb30SBarry Smith     ierr = PetscArraycpy(t->bembed,bembed,s);CHKERRQ(ierr);
470f68a32c8SEmil Constantinescu   }
471f68a32c8SEmil Constantinescu 
4723a8a9803SLisandro Dalcin   if (!binterp) { p = 1; binterp = t->b; }
4733a8a9803SLisandro Dalcin   t->p = p;
4743a8a9803SLisandro Dalcin   ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr);
475580bdb30SBarry Smith   ierr = PetscArraycpy(t->binterp,binterp,s*p);CHKERRQ(ierr);
4763a8a9803SLisandro Dalcin 
477f68a32c8SEmil Constantinescu   link->next = RKTableauList;
478f68a32c8SEmil Constantinescu   RKTableauList = link;
479f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
480f68a32c8SEmil Constantinescu }
481f68a32c8SEmil Constantinescu 
4820f7a28cdSStefano Zampini PetscErrorCode TSRKGetTableau_RK(TS ts, PetscInt *s, const PetscReal **A, const PetscReal **b, const PetscReal **c, const PetscReal **bembed,
4830f7a28cdSStefano Zampini                                         PetscInt *p, const PetscReal **binterp, PetscBool *FSAL)
4840f7a28cdSStefano Zampini {
4850f7a28cdSStefano Zampini   TS_RK     *rk   = (TS_RK*)ts->data;
4860f7a28cdSStefano Zampini   RKTableau tab  = rk->tableau;
4870f7a28cdSStefano Zampini 
4880f7a28cdSStefano Zampini   PetscFunctionBegin;
4890f7a28cdSStefano Zampini   if (s) *s = tab->s;
4900f7a28cdSStefano Zampini   if (A) *A = tab->A;
4910f7a28cdSStefano Zampini   if (b) *b = tab->b;
4920f7a28cdSStefano Zampini   if (c) *c = tab->c;
4930f7a28cdSStefano Zampini   if (bembed) *bembed = tab->bembed;
4940f7a28cdSStefano Zampini   if (p) *p = tab->p;
4950f7a28cdSStefano Zampini   if (binterp) *binterp = tab->binterp;
4960f7a28cdSStefano Zampini   if (FSAL) *FSAL = tab->FSAL;
4970f7a28cdSStefano Zampini   PetscFunctionReturn(0);
4980f7a28cdSStefano Zampini }
4990f7a28cdSStefano Zampini 
5000f7a28cdSStefano Zampini /*@C
5010f7a28cdSStefano Zampini    TSRKGetTableau - Get info on the RK tableau
5020f7a28cdSStefano Zampini 
5030f7a28cdSStefano Zampini    Not Collective
5040f7a28cdSStefano Zampini 
505f899ff85SJose E. Roman    Input Parameter:
5060f7a28cdSStefano Zampini .  ts - timestepping context
5070f7a28cdSStefano Zampini 
5080f7a28cdSStefano Zampini    Output Parameters:
5090f7a28cdSStefano Zampini +  s - number of stages, this is the dimension of the matrices below
5100f7a28cdSStefano Zampini .  A - stage coefficients (dimension s*s, row-major)
5110f7a28cdSStefano Zampini .  b - step completion table (dimension s)
5120f7a28cdSStefano Zampini .  c - abscissa (dimension s)
5130f7a28cdSStefano Zampini .  bembed - completion table for embedded method (dimension s; NULL if not available)
5140f7a28cdSStefano Zampini .  p - Order of the interpolation scheme, equal to the number of columns of binterp
5150f7a28cdSStefano Zampini .  binterp - Coefficients of the interpolation formula (dimension s*p)
5160f7a28cdSStefano Zampini -  FSAL - wheather or not the scheme has the First Same As Last property
5170f7a28cdSStefano Zampini 
5180f7a28cdSStefano Zampini    Level: developer
5190f7a28cdSStefano Zampini 
5200f7a28cdSStefano Zampini .seealso: TSRK
5210f7a28cdSStefano Zampini @*/
5220f7a28cdSStefano Zampini PetscErrorCode TSRKGetTableau(TS ts, PetscInt *s, const PetscReal **A, const PetscReal **b, const PetscReal **c, const PetscReal **bembed,
5230f7a28cdSStefano Zampini                                      PetscInt *p, const PetscReal **binterp, PetscBool *FSAL)
5240f7a28cdSStefano Zampini {
5250f7a28cdSStefano Zampini   PetscErrorCode ierr;
5260f7a28cdSStefano Zampini 
5270f7a28cdSStefano Zampini   PetscFunctionBegin;
5280f7a28cdSStefano Zampini   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5290f7a28cdSStefano Zampini   ierr = PetscUseMethod(ts,"TSRKGetTableau_C",(TS,PetscInt*,const PetscReal**,const PetscReal**,const PetscReal**,const PetscReal**,
5300f7a28cdSStefano Zampini                                                   PetscInt*,const PetscReal**,PetscBool*),(ts,s,A,b,c,bembed,p,binterp,FSAL));CHKERRQ(ierr);
5310f7a28cdSStefano Zampini   PetscFunctionReturn(0);
5320f7a28cdSStefano Zampini }
5330f7a28cdSStefano Zampini 
534e4dd521cSBarry Smith /*
5352c9cb062Sluzhanghpp  This is for single-step RK method
536f68a32c8SEmil Constantinescu  The step completion formula is
537e4dd521cSBarry Smith 
538f68a32c8SEmil Constantinescu  x1 = x0 + h b^T YdotRHS
539f68a32c8SEmil Constantinescu 
540f68a32c8SEmil Constantinescu  This function can be called before or after ts->vec_sol has been updated.
541f68a32c8SEmil Constantinescu  Suppose we have a completion formula (b) and an embedded formula (be) of different order.
542f68a32c8SEmil Constantinescu  We can write
543f68a32c8SEmil Constantinescu 
544f68a32c8SEmil Constantinescu  x1e = x0 + h be^T YdotRHS
545f68a32c8SEmil Constantinescu      = x1 - h b^T YdotRHS + h be^T YdotRHS
546f68a32c8SEmil Constantinescu      = x1 + h (be - b)^T YdotRHS
547f68a32c8SEmil Constantinescu 
548f68a32c8SEmil Constantinescu  so we can evaluate the method with different order even after the step has been optimistically completed.
549e4dd521cSBarry Smith */
550f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done)
551f68a32c8SEmil Constantinescu {
552f68a32c8SEmil Constantinescu   TS_RK          *rk   = (TS_RK*)ts->data;
553f68a32c8SEmil Constantinescu   RKTableau      tab  = rk->tableau;
554f68a32c8SEmil Constantinescu   PetscScalar    *w    = rk->work;
555f68a32c8SEmil Constantinescu   PetscReal      h;
556f68a32c8SEmil Constantinescu   PetscInt       s    = tab->s,j;
557f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
558f68a32c8SEmil Constantinescu 
559f68a32c8SEmil Constantinescu   PetscFunctionBegin;
560f68a32c8SEmil Constantinescu   switch (rk->status) {
561f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
562f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
563f68a32c8SEmil Constantinescu     h = ts->time_step; break;
564f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
565be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
566f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
567f68a32c8SEmil Constantinescu   }
568f68a32c8SEmil Constantinescu   if (order == tab->order) {
569f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) {
570f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
57190b67129SHong Zhang       for (j=0; j<s; j++) w[j] = h*tab->b[j]/rk->dtratio;
572f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
573f68a32c8SEmil Constantinescu     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
574f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
575f68a32c8SEmil Constantinescu   } else if (order == tab->order-1) {
576f68a32c8SEmil Constantinescu     if (!tab->bembed) goto unavailable;
577f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) { /*Complete with the embedded method (be)*/
578f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
579f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
580f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
581f68a32c8SEmil Constantinescu     } else {  /*Rollback and re-complete using (be-b) */
582f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
583f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
584f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
585f68a32c8SEmil Constantinescu     }
586f68a32c8SEmil Constantinescu     if (done) *done = PETSC_TRUE;
587f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
588f68a32c8SEmil Constantinescu   }
589f68a32c8SEmil Constantinescu unavailable:
590f68a32c8SEmil Constantinescu   if (done) *done = PETSC_FALSE;
591*98921bdaSJacob Faibussowitsch   else SETERRQ(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);
592f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
593f68a32c8SEmil Constantinescu }
594f68a32c8SEmil Constantinescu 
5950f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts)
5960f7a1166SHong Zhang {
5970f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
598cd4cee2dSHong Zhang   TS              quadts = ts->quadraturets;
5990f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
6000f7a1166SHong Zhang   const PetscInt  s = tab->s;
6010f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
6020f7a1166SHong Zhang   Vec             *Y = rk->Y;
6030f7a1166SHong Zhang   PetscInt        i;
6040f7a1166SHong Zhang   PetscErrorCode  ierr;
6050f7a1166SHong Zhang 
6060f7a1166SHong Zhang   PetscFunctionBegin;
607cd4cee2dSHong Zhang   /* No need to backup quadts->vec_sol since it can be reverted in TSRollBack_RK */
6080f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
609cd4cee2dSHong Zhang     /* Evolve quadrature TS solution to compute integrals */
610cd4cee2dSHong Zhang     ierr = TSComputeRHSFunction(quadts,rk->ptime+rk->time_step*c[i],Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
611cd4cee2dSHong Zhang     ierr = VecAXPY(quadts->vec_sol,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
6120f7a1166SHong Zhang   }
6130f7a1166SHong Zhang   PetscFunctionReturn(0);
6140f7a1166SHong Zhang }
6150f7a1166SHong Zhang 
6160f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts)
6170f7a1166SHong Zhang {
6180f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
6190f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
620cd4cee2dSHong Zhang   TS              quadts = ts->quadraturets;
6210f7a1166SHong Zhang   const PetscInt  s = tab->s;
6220f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
6230f7a1166SHong Zhang   Vec             *Y = rk->Y;
6240f7a1166SHong Zhang   PetscInt        i;
6250f7a1166SHong Zhang   PetscErrorCode  ierr;
6260f7a1166SHong Zhang 
6270f7a1166SHong Zhang   PetscFunctionBegin;
6280f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
629cd4cee2dSHong Zhang     /* Evolve quadrature TS solution to compute integrals */
630cd4cee2dSHong Zhang     ierr = TSComputeRHSFunction(quadts,ts->ptime+ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
631cd4cee2dSHong Zhang     ierr = VecAXPY(quadts->vec_sol,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
6320f7a1166SHong Zhang   }
6330f7a1166SHong Zhang   PetscFunctionReturn(0);
6340f7a1166SHong Zhang }
6350f7a1166SHong Zhang 
636fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts)
637fecfb714SLisandro Dalcin {
638fecfb714SLisandro Dalcin   TS_RK           *rk = (TS_RK*)ts->data;
639cd4cee2dSHong Zhang   TS              quadts = ts->quadraturets;
640fecfb714SLisandro Dalcin   RKTableau       tab = rk->tableau;
641fecfb714SLisandro Dalcin   const PetscInt  s  = tab->s;
642cd4cee2dSHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
643fecfb714SLisandro Dalcin   PetscScalar     *w = rk->work;
644cd4cee2dSHong Zhang   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
645fecfb714SLisandro Dalcin   PetscInt        j;
646fecfb714SLisandro Dalcin   PetscReal       h;
647fecfb714SLisandro Dalcin   PetscErrorCode  ierr;
648fecfb714SLisandro Dalcin 
649fecfb714SLisandro Dalcin   PetscFunctionBegin;
650fecfb714SLisandro Dalcin   switch (rk->status) {
651fecfb714SLisandro Dalcin   case TS_STEP_INCOMPLETE:
652fecfb714SLisandro Dalcin   case TS_STEP_PENDING:
653fecfb714SLisandro Dalcin     h = ts->time_step; break;
654fecfb714SLisandro Dalcin   case TS_STEP_COMPLETE:
655fecfb714SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
656fecfb714SLisandro Dalcin   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
657fecfb714SLisandro Dalcin   }
658fecfb714SLisandro Dalcin   for (j=0; j<s; j++) w[j] = -h*b[j];
659fecfb714SLisandro Dalcin   ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr);
660cd4cee2dSHong Zhang   if (quadts && ts->costintegralfwd) {
661cd4cee2dSHong Zhang     for (j=0; j<s; j++) {
662cd4cee2dSHong Zhang       /* Revert the quadrature TS solution */
663cd4cee2dSHong Zhang       ierr = TSComputeRHSFunction(quadts,rk->ptime+h*c[j],Y[j],ts->vec_costintegrand);CHKERRQ(ierr);
664cd4cee2dSHong Zhang       ierr = VecAXPY(quadts->vec_sol,-h*b[j],ts->vec_costintegrand);CHKERRQ(ierr);
665cd4cee2dSHong Zhang     }
666cd4cee2dSHong Zhang   }
667fecfb714SLisandro Dalcin   PetscFunctionReturn(0);
668fecfb714SLisandro Dalcin }
669fecfb714SLisandro Dalcin 
670922a638cSHong Zhang static PetscErrorCode TSForwardStep_RK(TS ts)
671922a638cSHong Zhang {
672922a638cSHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
673922a638cSHong Zhang   RKTableau       tab = rk->tableau;
674922a638cSHong Zhang   Mat             J,*MatsFwdSensipTemp = rk->MatsFwdSensipTemp;
675922a638cSHong Zhang   const PetscInt  s = tab->s;
676922a638cSHong Zhang   const PetscReal *A = tab->A,*c = tab->c,*b = tab->b;
677922a638cSHong Zhang   Vec             *Y = rk->Y;
678922a638cSHong Zhang   PetscInt        i,j;
679922a638cSHong Zhang   PetscReal       stage_time,h = ts->time_step;
680922a638cSHong Zhang   PetscBool       zero;
681922a638cSHong Zhang   PetscErrorCode  ierr;
682922a638cSHong Zhang 
683922a638cSHong Zhang   PetscFunctionBegin;
684922a638cSHong Zhang   ierr = MatCopy(ts->mat_sensip,rk->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
685922a638cSHong Zhang   ierr = TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr);
686922a638cSHong Zhang 
687922a638cSHong Zhang   for (i=0; i<s; i++) {
688922a638cSHong Zhang     stage_time = ts->ptime + h*c[i];
689922a638cSHong Zhang     zero = PETSC_FALSE;
690922a638cSHong Zhang     if (b[i] == 0 && i == s-1) zero = PETSC_TRUE;
691922a638cSHong Zhang     /* TLM Stage values */
692922a638cSHong Zhang     if (!i) {
693922a638cSHong Zhang       ierr = MatCopy(ts->mat_sensip,rk->MatsFwdStageSensip[i],SAME_NONZERO_PATTERN);CHKERRQ(ierr);
694922a638cSHong Zhang     } else if (!zero) {
695922a638cSHong Zhang       ierr = MatZeroEntries(rk->MatsFwdStageSensip[i]);CHKERRQ(ierr);
696922a638cSHong Zhang       for (j=0; j<i; j++) {
697922a638cSHong Zhang         ierr = MatAXPY(rk->MatsFwdStageSensip[i],h*A[i*s+j],MatsFwdSensipTemp[j],SAME_NONZERO_PATTERN);CHKERRQ(ierr);
698922a638cSHong Zhang       }
699922a638cSHong Zhang       ierr = MatAXPY(rk->MatsFwdStageSensip[i],1.,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
700922a638cSHong Zhang     } else {
701922a638cSHong Zhang       ierr = MatZeroEntries(rk->MatsFwdStageSensip[i]);CHKERRQ(ierr);
702922a638cSHong Zhang     }
703922a638cSHong Zhang 
704922a638cSHong Zhang     ierr = TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);CHKERRQ(ierr);
705922a638cSHong Zhang     ierr = MatMatMult(J,rk->MatsFwdStageSensip[i],MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatsFwdSensipTemp[i]);CHKERRQ(ierr);
70613af1a74SHong Zhang     if (ts->Jacprhs) {
70713af1a74SHong Zhang       ierr = TSComputeRHSJacobianP(ts,stage_time,Y[i],ts->Jacprhs);CHKERRQ(ierr); /* get f_p */
70813af1a74SHong Zhang       if (ts->vecs_sensi2p) { /* TLM used for 2nd-order adjoint */
70913af1a74SHong Zhang         PetscScalar *xarr;
71013af1a74SHong Zhang         ierr = MatDenseGetColumn(MatsFwdSensipTemp[i],0,&xarr);CHKERRQ(ierr);
71113af1a74SHong Zhang         ierr = VecPlaceArray(rk->VecDeltaFwdSensipCol,xarr);CHKERRQ(ierr);
71213af1a74SHong Zhang         ierr = MatMultAdd(ts->Jacprhs,ts->vec_dir,rk->VecDeltaFwdSensipCol,rk->VecDeltaFwdSensipCol);CHKERRQ(ierr);
713c3b366b1Sprj-         ierr = VecResetArray(rk->VecDeltaFwdSensipCol);CHKERRQ(ierr);
71413af1a74SHong Zhang         ierr = MatDenseRestoreColumn(MatsFwdSensipTemp[i],&xarr);CHKERRQ(ierr);
71513af1a74SHong Zhang       } else {
71613af1a74SHong Zhang         ierr = MatAXPY(MatsFwdSensipTemp[i],1.,ts->Jacprhs,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
71713af1a74SHong Zhang       }
718922a638cSHong Zhang     }
719922a638cSHong Zhang   }
720922a638cSHong Zhang 
721922a638cSHong Zhang   for (i=0; i<s; i++) {
722922a638cSHong Zhang     ierr = MatAXPY(ts->mat_sensip,h*b[i],rk->MatsFwdSensipTemp[i],SAME_NONZERO_PATTERN);CHKERRQ(ierr);
723922a638cSHong Zhang   }
724922a638cSHong Zhang   rk->status = TS_STEP_COMPLETE;
725922a638cSHong Zhang   PetscFunctionReturn(0);
726922a638cSHong Zhang }
727922a638cSHong Zhang 
728922a638cSHong Zhang static PetscErrorCode TSForwardGetStages_RK(TS ts,PetscInt *ns,Mat **stagesensip)
729922a638cSHong Zhang {
730922a638cSHong Zhang   TS_RK     *rk = (TS_RK*)ts->data;
731922a638cSHong Zhang   RKTableau tab  = rk->tableau;
732922a638cSHong Zhang 
733922a638cSHong Zhang   PetscFunctionBegin;
734922a638cSHong Zhang   if (ns) *ns = tab->s;
735922a638cSHong Zhang   if (stagesensip) *stagesensip = rk->MatsFwdStageSensip;
736922a638cSHong Zhang   PetscFunctionReturn(0);
737922a638cSHong Zhang }
738922a638cSHong Zhang 
739922a638cSHong Zhang static PetscErrorCode TSForwardSetUp_RK(TS ts)
740922a638cSHong Zhang {
741922a638cSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
742922a638cSHong Zhang   RKTableau      tab  = rk->tableau;
743922a638cSHong Zhang   PetscInt       i;
744922a638cSHong Zhang   PetscErrorCode ierr;
745922a638cSHong Zhang 
746922a638cSHong Zhang   PetscFunctionBegin;
747922a638cSHong Zhang   /* backup sensitivity results for roll-backs */
748922a638cSHong Zhang   ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatFwdSensip0);CHKERRQ(ierr);
749922a638cSHong Zhang 
750922a638cSHong Zhang   ierr = PetscMalloc1(tab->s,&rk->MatsFwdStageSensip);CHKERRQ(ierr);
751922a638cSHong Zhang   ierr = PetscMalloc1(tab->s,&rk->MatsFwdSensipTemp);CHKERRQ(ierr);
752922a638cSHong Zhang   for (i=0; i<tab->s; i++) {
753922a638cSHong Zhang     ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatsFwdStageSensip[i]);CHKERRQ(ierr);
754922a638cSHong Zhang     ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatsFwdSensipTemp[i]);CHKERRQ(ierr);
755922a638cSHong Zhang   }
756922a638cSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&rk->VecDeltaFwdSensipCol);CHKERRQ(ierr);
757922a638cSHong Zhang   PetscFunctionReturn(0);
758922a638cSHong Zhang }
759922a638cSHong Zhang 
760922a638cSHong Zhang static PetscErrorCode TSForwardReset_RK(TS ts)
761922a638cSHong Zhang {
762922a638cSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
763922a638cSHong Zhang   RKTableau      tab  = rk->tableau;
764922a638cSHong Zhang   PetscInt       i;
765922a638cSHong Zhang   PetscErrorCode ierr;
766922a638cSHong Zhang 
767922a638cSHong Zhang   PetscFunctionBegin;
768922a638cSHong Zhang   ierr = MatDestroy(&rk->MatFwdSensip0);CHKERRQ(ierr);
769922a638cSHong Zhang   if (rk->MatsFwdStageSensip) {
770922a638cSHong Zhang     for (i=0; i<tab->s; i++) {
771922a638cSHong Zhang       ierr = MatDestroy(&rk->MatsFwdStageSensip[i]);CHKERRQ(ierr);
772922a638cSHong Zhang     }
773922a638cSHong Zhang     ierr = PetscFree(rk->MatsFwdStageSensip);CHKERRQ(ierr);
774922a638cSHong Zhang   }
775922a638cSHong Zhang   if (rk->MatsFwdSensipTemp) {
776922a638cSHong Zhang     for (i=0; i<tab->s; i++) {
777922a638cSHong Zhang       ierr = MatDestroy(&rk->MatsFwdSensipTemp[i]);CHKERRQ(ierr);
778922a638cSHong Zhang     }
779922a638cSHong Zhang     ierr = PetscFree(rk->MatsFwdSensipTemp);CHKERRQ(ierr);
780922a638cSHong Zhang   }
781922a638cSHong Zhang   ierr = VecDestroy(&rk->VecDeltaFwdSensipCol);CHKERRQ(ierr);
782922a638cSHong Zhang   PetscFunctionReturn(0);
783922a638cSHong Zhang }
784922a638cSHong Zhang 
785f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts)
786f68a32c8SEmil Constantinescu {
787f68a32c8SEmil Constantinescu   TS_RK           *rk  = (TS_RK*)ts->data;
788f68a32c8SEmil Constantinescu   RKTableau       tab  = rk->tableau;
789f68a32c8SEmil Constantinescu   const PetscInt  s = tab->s;
790fecfb714SLisandro Dalcin   const PetscReal *A = tab->A,*c = tab->c;
791f68a32c8SEmil Constantinescu   PetscScalar     *w = rk->work;
792f68a32c8SEmil Constantinescu   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
793d1905564SLisandro Dalcin   PetscBool       FSAL = tab->FSAL;
794f68a32c8SEmil Constantinescu   TSAdapt         adapt;
795fecfb714SLisandro Dalcin   PetscInt        i,j;
796be5899b3SLisandro Dalcin   PetscInt        rejections = 0;
797be5899b3SLisandro Dalcin   PetscBool       stageok,accept = PETSC_TRUE;
798be5899b3SLisandro Dalcin   PetscReal       next_time_step = ts->time_step;
799f68a32c8SEmil Constantinescu   PetscErrorCode  ierr;
800f68a32c8SEmil Constantinescu 
801f68a32c8SEmil Constantinescu   PetscFunctionBegin;
802d1905564SLisandro Dalcin   if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE;
803d1905564SLisandro Dalcin   if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); }
804d1905564SLisandro Dalcin 
805f68a32c8SEmil Constantinescu   rk->status = TS_STEP_INCOMPLETE;
806be5899b3SLisandro Dalcin   while (!ts->reason && rk->status != TS_STEP_COMPLETE) {
807be5899b3SLisandro Dalcin     PetscReal t = ts->ptime;
808f68a32c8SEmil Constantinescu     PetscReal h = ts->time_step;
809f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
8109be3e283SDebojyoti Ghosh       rk->stage_time = t + h*c[i];
8119be3e283SDebojyoti Ghosh       ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr);
812f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
813f68a32c8SEmil Constantinescu       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
814f68a32c8SEmil Constantinescu       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
8159be3e283SDebojyoti Ghosh       ierr = TSPostStage(ts,rk->stage_time,i,Y);CHKERRQ(ierr);
816f68a32c8SEmil Constantinescu       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
817be5899b3SLisandro Dalcin       ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr);
818fecfb714SLisandro Dalcin       if (!stageok) goto reject_step;
8198f738a7cSLisandro Dalcin       if (FSAL && !i) continue;
820f68a32c8SEmil Constantinescu       ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
821f68a32c8SEmil Constantinescu     }
822be5899b3SLisandro Dalcin 
823be5899b3SLisandro Dalcin     rk->status = TS_STEP_INCOMPLETE;
824f68a32c8SEmil Constantinescu     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
825f68a32c8SEmil Constantinescu     rk->status = TS_STEP_PENDING;
826f68a32c8SEmil Constantinescu     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
827f68a32c8SEmil Constantinescu     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
8281917a363SLisandro Dalcin     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr);
829fecfb714SLisandro Dalcin     ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
830be5899b3SLisandro Dalcin     rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
831be5899b3SLisandro Dalcin     if (!accept) { /* Roll back the current step */
832fecfb714SLisandro Dalcin       ierr = TSRollBack_RK(ts);CHKERRQ(ierr);
833be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
834be5899b3SLisandro Dalcin       goto reject_step;
835be5899b3SLisandro Dalcin     }
836be5899b3SLisandro Dalcin 
8370f7a1166SHong Zhang     if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
8380f7a1166SHong Zhang       rk->ptime     = ts->ptime;
8390f7a1166SHong Zhang       rk->time_step = ts->time_step;
840493ed95fSHong Zhang     }
841be5899b3SLisandro Dalcin 
842f68a32c8SEmil Constantinescu     ts->ptime += ts->time_step;
843f68a32c8SEmil Constantinescu     ts->time_step = next_time_step;
844f68a32c8SEmil Constantinescu     break;
845be5899b3SLisandro Dalcin 
846be5899b3SLisandro Dalcin     reject_step:
847fecfb714SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
848be5899b3SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
849be5899b3SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
850be5899b3SLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
851f68a32c8SEmil Constantinescu     }
852f68a32c8SEmil Constantinescu   }
853f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
854e4dd521cSBarry Smith }
855e4dd521cSBarry Smith 
856f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts)
857f6a906c0SBarry Smith {
858f6a906c0SBarry Smith   TS_RK          *rk  = (TS_RK*)ts->data;
859be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
860be5899b3SLisandro Dalcin   PetscInt       s   = tab->s;
861f6a906c0SBarry Smith   PetscErrorCode ierr;
862f6a906c0SBarry Smith 
863f6a906c0SBarry Smith   PetscFunctionBegin;
864f6a906c0SBarry Smith   if (ts->adjointsetupcalled++) PetscFunctionReturn(0);
8652e7b7f96SHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecsDeltaLam);CHKERRQ(ierr);
8662e7b7f96SHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&rk->VecsSensiTemp);CHKERRQ(ierr);
867f6a906c0SBarry Smith   if (ts->vecs_sensip) {
8682e7b7f96SHong Zhang     ierr = VecDuplicate(ts->vecs_sensip[0],&rk->VecDeltaMu);CHKERRQ(ierr);
869f6a906c0SBarry Smith   }
87013af1a74SHong Zhang   if (ts->vecs_sensi2) {
87113af1a74SHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecsDeltaLam2);CHKERRQ(ierr);
87213af1a74SHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&rk->VecsSensi2Temp);CHKERRQ(ierr);
87313af1a74SHong Zhang   }
87413af1a74SHong Zhang   if (ts->vecs_sensi2p) {
87513af1a74SHong Zhang     ierr = VecDuplicate(ts->vecs_sensi2p[0],&rk->VecDeltaMu2);CHKERRQ(ierr);
87613af1a74SHong Zhang   }
877f6a906c0SBarry Smith   PetscFunctionReturn(0);
878f6a906c0SBarry Smith }
879f6a906c0SBarry Smith 
88035f5172aSHong Zhang /*
88135f5172aSHong Zhang   Assumptions:
88235f5172aSHong Zhang     - TSStep_RK() always evaluates the step with b, not bembed.
88335f5172aSHong Zhang */
88442f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts)
885d2daff3dSHong Zhang {
886c235aa8dSHong Zhang   TS_RK            *rk = (TS_RK*)ts->data;
887cd4cee2dSHong Zhang   TS               quadts = ts->quadraturets;
888c235aa8dSHong Zhang   RKTableau        tab = rk->tableau;
8893ca0882eSHong Zhang   Mat              J,Jpre,Jquad;
890c235aa8dSHong Zhang   const PetscInt   s = tab->s;
891c235aa8dSHong Zhang   const PetscReal  *A = tab->A,*b = tab->b,*c = tab->c;
89213af1a74SHong Zhang   PetscScalar      *w = rk->work,*xarr;
8932e7b7f96SHong Zhang   Vec              *Y = rk->Y,*VecsDeltaLam = rk->VecsDeltaLam,VecDeltaMu = rk->VecDeltaMu,*VecsSensiTemp = rk->VecsSensiTemp;
89413af1a74SHong Zhang   Vec              *VecsDeltaLam2 = rk->VecsDeltaLam2,VecDeltaMu2 = rk->VecDeltaMu2,*VecsSensi2Temp = rk->VecsSensi2Temp;
895cd4cee2dSHong Zhang   Vec              VecDRDUTransCol = ts->vec_drdu_col,VecDRDPTransCol = ts->vec_drdp_col;
896b18ea86cSHong Zhang   PetscInt         i,j,nadj;
8973d3eaba7SBarry Smith   PetscReal        t = ts->ptime;
8983ca0882eSHong Zhang   PetscReal        h = ts->time_step;
899d2daff3dSHong Zhang   PetscErrorCode   ierr;
900c235aa8dSHong Zhang 
901d2daff3dSHong Zhang   PetscFunctionBegin;
902c235aa8dSHong Zhang   rk->status = TS_STEP_INCOMPLETE;
9033389c7d9SStefano Zampini 
9043ca0882eSHong Zhang   ierr = TSGetRHSJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr);
90535f5172aSHong Zhang   if (quadts) {
906cd4cee2dSHong Zhang     ierr = TSGetRHSJacobian(quadts,&Jquad,NULL,NULL,NULL);CHKERRQ(ierr);
90735f5172aSHong Zhang   }
9082e7b7f96SHong Zhang   for (i=s-1; i>=0; i--) {
9096a1e4597SHong Zhang     if (tab->FSAL && i == s-1) {
9106a1e4597SHong Zhang       /* VecsDeltaLam[nadj*s+s-1] are initialized with zeros and the values never change.*/
9116a1e4597SHong Zhang       continue;
9126a1e4597SHong Zhang     }
9133ca0882eSHong Zhang     rk->stage_time = t + h*(1.0-c[i]);
914fd850964SHong Zhang     ierr = TSComputeSNESJacobian(ts,Y[i],J,Jpre);CHKERRQ(ierr);
91535f5172aSHong Zhang     if (quadts) {
9163ca0882eSHong Zhang       ierr = TSComputeRHSJacobian(quadts,rk->stage_time,Y[i],Jquad,Jquad);CHKERRQ(ierr); /* get r_u^T */
91735f5172aSHong Zhang     }
9183389c7d9SStefano Zampini     if (ts->vecs_sensip) {
9193ca0882eSHong Zhang       ierr = TSComputeRHSJacobianP(ts,rk->stage_time,Y[i],ts->Jacprhs);CHKERRQ(ierr); /* get f_p */
92035f5172aSHong Zhang       if (quadts) {
9213ca0882eSHong Zhang         ierr = TSComputeRHSJacobianP(quadts,rk->stage_time,Y[i],quadts->Jacprhs);CHKERRQ(ierr); /* get f_p for the quadrature */
9223389c7d9SStefano Zampini       }
92335f5172aSHong Zhang     }
9243389c7d9SStefano Zampini 
92535f5172aSHong Zhang     if (b[i]) {
92635f5172aSHong Zhang       for (j=i+1; j<s; j++) w[j-i-1] = A[j*s+i]/b[i]; /* coefficients for computing VecsSensiTemp */
92735f5172aSHong Zhang     } else {
92835f5172aSHong Zhang       for (j=i+1; j<s; j++) w[j-i-1] = A[j*s+i]; /* coefficients for computing VecsSensiTemp */
92935f5172aSHong Zhang     }
9302e7b7f96SHong Zhang 
9312e7b7f96SHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
9323389c7d9SStefano Zampini       /* Stage values of lambda */
93335f5172aSHong Zhang       if (b[i]) {
93435f5172aSHong Zhang         /* lambda_{n+1} + \sum_{j=i+1}^s a_{ji}/b[i]*lambda_{s,j} */
9352e7b7f96SHong Zhang         ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); /* VecDeltaLam is an vec array of size s by numcost */
9362e7b7f96SHong Zhang         ierr = VecMAXPY(VecsSensiTemp[nadj],s-i-1,w,&VecsDeltaLam[nadj*s+i+1]);CHKERRQ(ierr);
93735f5172aSHong Zhang         ierr = MatMultTranspose(J,VecsSensiTemp[nadj],VecsDeltaLam[nadj*s+i]);CHKERRQ(ierr); /* VecsSensiTemp will be reused by 2nd-order adjoint */
93835f5172aSHong Zhang         ierr = VecScale(VecsDeltaLam[nadj*s+i],-h*b[i]);CHKERRQ(ierr);
93935f5172aSHong Zhang         if (quadts) {
940cd4cee2dSHong Zhang           ierr = MatDenseGetColumn(Jquad,nadj,&xarr);CHKERRQ(ierr);
941cd4cee2dSHong Zhang           ierr = VecPlaceArray(VecDRDUTransCol,xarr);CHKERRQ(ierr);
942cd4cee2dSHong Zhang           ierr = VecAXPY(VecsDeltaLam[nadj*s+i],-h*b[i],VecDRDUTransCol);CHKERRQ(ierr);
943cd4cee2dSHong Zhang           ierr = VecResetArray(VecDRDUTransCol);CHKERRQ(ierr);
944cd4cee2dSHong Zhang           ierr = MatDenseRestoreColumn(Jquad,&xarr);CHKERRQ(ierr);
945cd4cee2dSHong Zhang         }
9463389c7d9SStefano Zampini       } else {
9476a1e4597SHong Zhang         /* \sum_{j=i+1}^s a_{ji}*lambda_{s,j} */
9486a1e4597SHong Zhang         ierr = VecSet(VecsSensiTemp[nadj],0);CHKERRQ(ierr);
9496a1e4597SHong Zhang         ierr = VecMAXPY(VecsSensiTemp[nadj],s-i-1,w,&VecsDeltaLam[nadj*s+i+1]);CHKERRQ(ierr);
9506a1e4597SHong Zhang         ierr = MatMultTranspose(J,VecsSensiTemp[nadj],VecsDeltaLam[nadj*s+i]);CHKERRQ(ierr);
9516a1e4597SHong Zhang         ierr = VecScale(VecsDeltaLam[nadj*s+i],-h);CHKERRQ(ierr);
9523389c7d9SStefano Zampini       }
9536497ce14SHong Zhang 
954ad8e2604SHong Zhang       /* Stage values of mu */
9556497ce14SHong Zhang       if (ts->vecs_sensip) {
95613af1a74SHong Zhang         ierr = MatMultTranspose(ts->Jacprhs,VecsSensiTemp[nadj],VecDeltaMu);CHKERRQ(ierr);
95735f5172aSHong Zhang         if (b[i]) {
95835f5172aSHong Zhang           ierr = VecScale(VecDeltaMu,-h*b[i]);CHKERRQ(ierr);
95935f5172aSHong Zhang           if (quadts) {
960cd4cee2dSHong Zhang             ierr = MatDenseGetColumn(quadts->Jacprhs,nadj,&xarr);CHKERRQ(ierr);
961cd4cee2dSHong Zhang             ierr = VecPlaceArray(VecDRDPTransCol,xarr);CHKERRQ(ierr);
962cd4cee2dSHong Zhang             ierr = VecAXPY(VecDeltaMu,-h*b[i],VecDRDPTransCol);CHKERRQ(ierr);
963cd4cee2dSHong Zhang             ierr = VecResetArray(VecDRDPTransCol);CHKERRQ(ierr);
964cd4cee2dSHong Zhang             ierr = MatDenseRestoreColumn(quadts->Jacprhs,&xarr);CHKERRQ(ierr);
965493ed95fSHong Zhang           }
96635f5172aSHong Zhang         } else {
96735f5172aSHong Zhang           ierr = VecScale(VecDeltaMu,-h);CHKERRQ(ierr);
96835f5172aSHong Zhang         }
9692e7b7f96SHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],1.,VecDeltaMu);CHKERRQ(ierr); /* update sensip for each stage */
970ad8e2604SHong Zhang       }
971c235aa8dSHong Zhang     }
97213af1a74SHong Zhang 
97313af1a74SHong Zhang     if (ts->vecs_sensi2 && ts->forward_solve) { /* 2nd-order adjoint, TLM mode has to be turned on */
97413af1a74SHong Zhang       /* Get w1 at t_{n+1} from TLM matrix */
97513af1a74SHong Zhang       ierr = MatDenseGetColumn(rk->MatsFwdStageSensip[i],0,&xarr);CHKERRQ(ierr);
97613af1a74SHong Zhang       ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
97713af1a74SHong Zhang       /* lambda_s^T F_UU w_1 */
9783ca0882eSHong Zhang       ierr = TSComputeRHSHessianProductFunctionUU(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_sensip_col,ts->vecs_guu);CHKERRQ(ierr);
97935f5172aSHong Zhang       if (quadts)  {
98035f5172aSHong Zhang         /* R_UU w_1 */
9813ca0882eSHong Zhang         ierr = TSComputeRHSHessianProductFunctionUU(quadts,rk->stage_time,Y[i],NULL,ts->vec_sensip_col,ts->vecs_guu);CHKERRQ(ierr);
98235f5172aSHong Zhang       }
98335f5172aSHong Zhang       if (ts->vecs_sensip) {
98413af1a74SHong Zhang         /* lambda_s^T F_UP w_2 */
9853ca0882eSHong Zhang         ierr = TSComputeRHSHessianProductFunctionUP(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_dir,ts->vecs_gup);CHKERRQ(ierr);
98635f5172aSHong Zhang         if (quadts)  {
98735f5172aSHong Zhang           /* R_UP w_2 */
9883ca0882eSHong Zhang           ierr = TSComputeRHSHessianProductFunctionUP(quadts,rk->stage_time,Y[i],NULL,ts->vec_sensip_col,ts->vecs_gup);CHKERRQ(ierr);
98935f5172aSHong Zhang         }
99035f5172aSHong Zhang       }
99113af1a74SHong Zhang       if (ts->vecs_sensi2p) {
99213af1a74SHong Zhang         /* lambda_s^T F_PU w_1 */
9933ca0882eSHong Zhang         ierr = TSComputeRHSHessianProductFunctionPU(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_sensip_col,ts->vecs_gpu);CHKERRQ(ierr);
99435f5172aSHong Zhang         /* lambda_s^T F_PP w_2 */
9953ca0882eSHong Zhang         ierr = TSComputeRHSHessianProductFunctionPP(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_dir,ts->vecs_gpp);CHKERRQ(ierr);
99635f5172aSHong Zhang         if (b[i] && quadts) {
99735f5172aSHong Zhang           /* R_PU w_1 */
9983ca0882eSHong Zhang           ierr = TSComputeRHSHessianProductFunctionPU(quadts,rk->stage_time,Y[i],NULL,ts->vec_sensip_col,ts->vecs_gpu);CHKERRQ(ierr);
99935f5172aSHong Zhang           /* R_PP w_2 */
10003ca0882eSHong Zhang           ierr = TSComputeRHSHessianProductFunctionPP(quadts,rk->stage_time,Y[i],NULL,ts->vec_dir,ts->vecs_gpp);CHKERRQ(ierr);
100135f5172aSHong Zhang         }
100213af1a74SHong Zhang       }
100313af1a74SHong Zhang       ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
100413af1a74SHong Zhang       ierr = MatDenseRestoreColumn(rk->MatsFwdStageSensip[i],&xarr);CHKERRQ(ierr);
100513af1a74SHong Zhang 
100613af1a74SHong Zhang       for (nadj=0; nadj<ts->numcost; nadj++) {
100713af1a74SHong Zhang         /* Stage values of lambda */
100835f5172aSHong Zhang         if (b[i]) {
100935f5172aSHong Zhang           /* J_i^T*(Lambda_{n+1}+\sum_{j=i+1}^s a_{ji}/b_i*Lambda_{s,j} */
101013af1a74SHong Zhang           ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
101113af1a74SHong Zhang           ierr = VecMAXPY(VecsSensi2Temp[nadj],s-i-1,w,&VecsDeltaLam2[nadj*s+i+1]);CHKERRQ(ierr);
101213af1a74SHong Zhang           ierr = MatMultTranspose(J,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj*s+i]);CHKERRQ(ierr);
101335f5172aSHong Zhang           ierr = VecScale(VecsDeltaLam2[nadj*s+i],-h*b[i]);CHKERRQ(ierr);
101435f5172aSHong Zhang           ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],-h*b[i],ts->vecs_guu[nadj]);CHKERRQ(ierr);
101535f5172aSHong Zhang           if (ts->vecs_sensip) {
101635f5172aSHong Zhang             ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],-h*b[i],ts->vecs_gup[nadj]);CHKERRQ(ierr);
101713af1a74SHong Zhang           }
101813af1a74SHong Zhang         } else {
101935f5172aSHong Zhang           /* \sum_{j=i+1}^s a_{ji}*Lambda_{s,j} */
102013af1a74SHong Zhang           ierr = VecSet(VecsDeltaLam2[nadj*s+i],0);CHKERRQ(ierr);
102135f5172aSHong Zhang           ierr = VecMAXPY(VecsSensi2Temp[nadj],s-i-1,w,&VecsDeltaLam2[nadj*s+i+1]);CHKERRQ(ierr);
102235f5172aSHong Zhang           ierr = MatMultTranspose(J,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj*s+i]);CHKERRQ(ierr);
102335f5172aSHong Zhang           ierr = VecScale(VecsDeltaLam2[nadj*s+i],-h);CHKERRQ(ierr);
102435f5172aSHong Zhang           ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],-h,ts->vecs_guu[nadj]);CHKERRQ(ierr);
102535f5172aSHong Zhang           if (ts->vecs_sensip) {
102635f5172aSHong Zhang             ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],-h,ts->vecs_gup[nadj]);CHKERRQ(ierr);
102713af1a74SHong Zhang           }
102835f5172aSHong Zhang         }
102935f5172aSHong Zhang         if (ts->vecs_sensi2p) { /* 2nd-order adjoint for parameters */
103013af1a74SHong Zhang           ierr = MatMultTranspose(ts->Jacprhs,VecsSensi2Temp[nadj],VecDeltaMu2);CHKERRQ(ierr);
103135f5172aSHong Zhang           if (b[i]) {
103235f5172aSHong Zhang             ierr = VecScale(VecDeltaMu2,-h*b[i]);CHKERRQ(ierr);
103335f5172aSHong Zhang             ierr = VecAXPY(VecDeltaMu2,-h*b[i],ts->vecs_gpu[nadj]);CHKERRQ(ierr);
103435f5172aSHong Zhang             ierr = VecAXPY(VecDeltaMu2,-h*b[i],ts->vecs_gpp[nadj]);CHKERRQ(ierr);
103513af1a74SHong Zhang           } else {
103635f5172aSHong Zhang             ierr = VecScale(VecDeltaMu2,-h);CHKERRQ(ierr);
103735f5172aSHong Zhang             ierr = VecAXPY(VecDeltaMu2,-h,ts->vecs_gpu[nadj]);CHKERRQ(ierr);
103835f5172aSHong Zhang             ierr = VecAXPY(VecDeltaMu2,-h,ts->vecs_gpp[nadj]);CHKERRQ(ierr);
103913af1a74SHong Zhang           }
104013af1a74SHong Zhang           ierr = VecAXPY(ts->vecs_sensi2p[nadj],1,VecDeltaMu2);CHKERRQ(ierr); /* update sensi2p for each stage */
104113af1a74SHong Zhang         }
104213af1a74SHong Zhang       }
104313af1a74SHong Zhang     }
10446497ce14SHong Zhang   }
1045c235aa8dSHong Zhang 
1046c235aa8dSHong Zhang   for (j=0; j<s; j++) w[j] = 1.0;
10472e7b7f96SHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) { /* no need to do this for mu's */
10482e7b7f96SHong Zhang     ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecsDeltaLam[nadj*s]);CHKERRQ(ierr);
104913af1a74SHong Zhang     if (ts->vecs_sensi2) {
105013af1a74SHong Zhang       ierr = VecMAXPY(ts->vecs_sensi2[nadj],s,w,&VecsDeltaLam2[nadj*s]);CHKERRQ(ierr);
105113af1a74SHong Zhang     }
10526497ce14SHong Zhang   }
1053c235aa8dSHong Zhang   rk->status = TS_STEP_COMPLETE;
1054d2daff3dSHong Zhang   PetscFunctionReturn(0);
1055d2daff3dSHong Zhang }
1056d2daff3dSHong Zhang 
105713af1a74SHong Zhang static PetscErrorCode TSAdjointReset_RK(TS ts)
105813af1a74SHong Zhang {
105913af1a74SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
106013af1a74SHong Zhang   RKTableau      tab = rk->tableau;
106113af1a74SHong Zhang   PetscErrorCode ierr;
106213af1a74SHong Zhang 
106313af1a74SHong Zhang   PetscFunctionBegin;
106413af1a74SHong Zhang   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam);CHKERRQ(ierr);
106513af1a74SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&rk->VecsSensiTemp);CHKERRQ(ierr);
106613af1a74SHong Zhang   ierr = VecDestroy(&rk->VecDeltaMu);CHKERRQ(ierr);
106713af1a74SHong Zhang   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam2);CHKERRQ(ierr);
106813af1a74SHong Zhang   ierr = VecDestroy(&rk->VecDeltaMu2);CHKERRQ(ierr);
106913af1a74SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&rk->VecsSensi2Temp);CHKERRQ(ierr);
107013af1a74SHong Zhang   PetscFunctionReturn(0);
107113af1a74SHong Zhang }
107213af1a74SHong Zhang 
107355de54fcSHong Zhang static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X)
107455de54fcSHong Zhang {
107555de54fcSHong Zhang   TS_RK            *rk = (TS_RK*)ts->data;
107655de54fcSHong Zhang   PetscInt         s  = rk->tableau->s,p = rk->tableau->p,i,j;
107755de54fcSHong Zhang   PetscReal        h;
107855de54fcSHong Zhang   PetscReal        tt,t;
107955de54fcSHong Zhang   PetscScalar      *b;
108055de54fcSHong Zhang   const PetscReal  *B = rk->tableau->binterp;
108155de54fcSHong Zhang   PetscErrorCode   ierr;
108255de54fcSHong Zhang 
108355de54fcSHong Zhang   PetscFunctionBegin;
1084*98921bdaSJacob Faibussowitsch   if (!B) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
108555de54fcSHong Zhang 
108655de54fcSHong Zhang   switch (rk->status) {
108755de54fcSHong Zhang     case TS_STEP_INCOMPLETE:
108855de54fcSHong Zhang     case TS_STEP_PENDING:
108955de54fcSHong Zhang       h = ts->time_step;
109055de54fcSHong Zhang       t = (itime - ts->ptime)/h;
109155de54fcSHong Zhang       break;
109255de54fcSHong Zhang     case TS_STEP_COMPLETE:
109355de54fcSHong Zhang       h = ts->ptime - ts->ptime_prev;
109455de54fcSHong Zhang       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
109555de54fcSHong Zhang       break;
109655de54fcSHong Zhang     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
109755de54fcSHong Zhang   }
109855de54fcSHong Zhang   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
109955de54fcSHong Zhang   for (i=0; i<s; i++) b[i] = 0;
110055de54fcSHong Zhang   for (j=0,tt=t; j<p; j++,tt*=t) {
110155de54fcSHong Zhang     for (i=0; i<s; i++) {
110255de54fcSHong Zhang       b[i]  += h * B[i*p+j] * tt;
110355de54fcSHong Zhang     }
110455de54fcSHong Zhang   }
110555de54fcSHong Zhang   ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr);
110655de54fcSHong Zhang   ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr);
110755de54fcSHong Zhang   ierr = PetscFree(b);CHKERRQ(ierr);
110855de54fcSHong Zhang   PetscFunctionReturn(0);
110955de54fcSHong Zhang }
111055de54fcSHong Zhang 
111155de54fcSHong Zhang /*------------------------------------------------------------*/
111255de54fcSHong Zhang 
1113be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts)
1114be5899b3SLisandro Dalcin {
1115be5899b3SLisandro Dalcin   TS_RK          *rk = (TS_RK*)ts->data;
1116be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
1117be5899b3SLisandro Dalcin   PetscErrorCode ierr;
1118be5899b3SLisandro Dalcin 
1119be5899b3SLisandro Dalcin   PetscFunctionBegin;
1120be5899b3SLisandro Dalcin   if (!tab) PetscFunctionReturn(0);
1121be5899b3SLisandro Dalcin   ierr = PetscFree(rk->work);CHKERRQ(ierr);
1122be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr);
1123be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr);
1124be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
1125be5899b3SLisandro Dalcin }
1126be5899b3SLisandro Dalcin 
1127277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts)
1128e4dd521cSBarry Smith {
11296849ba73SBarry Smith   PetscErrorCode ierr;
1130e4dd521cSBarry Smith 
1131e4dd521cSBarry Smith   PetscFunctionBegin;
1132be5899b3SLisandro Dalcin   ierr = TSRKTableauReset(ts);CHKERRQ(ierr);
11330fe4d17eSHong Zhang   if (ts->use_splitrhsfunction) {
113402555c94SHong Zhang     ierr = PetscTryMethod(ts,"TSReset_RK_MultirateSplit_C",(TS),(ts));CHKERRQ(ierr);
11350fe4d17eSHong Zhang   } else {
113602555c94SHong Zhang     ierr = PetscTryMethod(ts,"TSReset_RK_MultirateNonsplit_C",(TS),(ts));CHKERRQ(ierr);
11370fe4d17eSHong Zhang   }
1138277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
1139e4dd521cSBarry Smith }
1140277b19d0SLisandro Dalcin 
1141f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx)
1142f68a32c8SEmil Constantinescu {
1143f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1144f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1145f68a32c8SEmil Constantinescu }
1146f68a32c8SEmil Constantinescu 
1147f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
1148f68a32c8SEmil Constantinescu {
1149f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1150f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1151f68a32c8SEmil Constantinescu }
1152f68a32c8SEmil Constantinescu 
1153f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx)
1154f68a32c8SEmil Constantinescu {
1155f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1156f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1157f68a32c8SEmil Constantinescu }
1158f68a32c8SEmil Constantinescu 
1159f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
1160f68a32c8SEmil Constantinescu {
1161f68a32c8SEmil Constantinescu 
1162f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1163f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1164f68a32c8SEmil Constantinescu }
1165be5899b3SLisandro Dalcin 
1166be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts)
1167be5899b3SLisandro Dalcin {
1168be5899b3SLisandro Dalcin   TS_RK          *rk  = (TS_RK*)ts->data;
1169be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
1170be5899b3SLisandro Dalcin   PetscErrorCode ierr;
1171be5899b3SLisandro Dalcin 
1172be5899b3SLisandro Dalcin   PetscFunctionBegin;
1173be5899b3SLisandro Dalcin   ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr);
1174be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr);
1175be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr);
1176be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
1177be5899b3SLisandro Dalcin }
1178be5899b3SLisandro Dalcin 
1179f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts)
1180f68a32c8SEmil Constantinescu {
1181cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
1182f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
1183f68a32c8SEmil Constantinescu   DM             dm;
1184f68a32c8SEmil Constantinescu 
1185f68a32c8SEmil Constantinescu   PetscFunctionBegin;
11863dd83b38SBarry Smith   ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr);
1187be5899b3SLisandro Dalcin   ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);
1188cd4cee2dSHong Zhang   if (quadts && ts->costintegralfwd) {
1189cd4cee2dSHong Zhang     Mat Jquad;
1190cd4cee2dSHong Zhang     ierr = TSGetRHSJacobian(quadts,&Jquad,NULL,NULL,NULL);CHKERRQ(ierr);
11910f7a1166SHong Zhang   }
1192f68a32c8SEmil Constantinescu   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1193f68a32c8SEmil Constantinescu   ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
1194f68a32c8SEmil Constantinescu   ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
11950fe4d17eSHong Zhang   if (ts->use_splitrhsfunction) {
119602555c94SHong Zhang     ierr = PetscTryMethod(ts,"TSSetUp_RK_MultirateSplit_C",(TS),(ts));CHKERRQ(ierr);
11970fe4d17eSHong Zhang   } else {
119802555c94SHong Zhang     ierr = PetscTryMethod(ts,"TSSetUp_RK_MultirateNonsplit_C",(TS),(ts));CHKERRQ(ierr);
11990fe4d17eSHong Zhang   }
1200e4dd521cSBarry Smith   PetscFunctionReturn(0);
1201e4dd521cSBarry Smith }
1202d2daff3dSHong Zhang 
12034416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts)
1204e4dd521cSBarry Smith {
1205be5899b3SLisandro Dalcin   TS_RK          *rk = (TS_RK*)ts->data;
1206dfbe8321SBarry Smith   PetscErrorCode ierr;
1207262ff7bbSBarry Smith 
1208e4dd521cSBarry Smith   PetscFunctionBegin;
1209e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr);
1210f68a32c8SEmil Constantinescu   {
1211f68a32c8SEmil Constantinescu     RKTableauLink link;
1212f68a32c8SEmil Constantinescu     PetscInt      count,choice;
12130fe4d17eSHong Zhang     PetscBool     flg,use_multirate = PETSC_FALSE;
1214f68a32c8SEmil Constantinescu     const char    **namelist;
12152c9cb062Sluzhanghpp 
1216f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) ;
1217f489ac74SBarry Smith     ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr);
1218f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
12190fe4d17eSHong Zhang     ierr = PetscOptionsBool("-ts_rk_multirate","Use interpolation-based multirate RK method","TSRKSetMultirate",rk->use_multirate,&use_multirate,&flg);CHKERRQ(ierr);
12200fe4d17eSHong Zhang     if (flg) {
12210fe4d17eSHong Zhang       ierr = TSRKSetMultirate(ts,use_multirate);CHKERRQ(ierr);
12220fe4d17eSHong Zhang     }
1223be5899b3SLisandro Dalcin     ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr);
1224be5899b3SLisandro Dalcin     if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);}
1225f68a32c8SEmil Constantinescu     ierr = PetscFree(namelist);CHKERRQ(ierr);
1226f68a32c8SEmil Constantinescu   }
1227262ff7bbSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
1228ea13f565SStefano Zampini   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)ts),NULL,"Multirate methods options","");CHKERRQ(ierr);
12292c9cb062Sluzhanghpp   ierr = PetscOptionsInt("-ts_rk_dtratio","time step ratio between slow and fast","",rk->dtratio,&rk->dtratio,NULL);CHKERRQ(ierr);
12302c9cb062Sluzhanghpp   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1231e4dd521cSBarry Smith   PetscFunctionReturn(0);
1232e4dd521cSBarry Smith }
1233e4dd521cSBarry Smith 
12345f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
1235e4dd521cSBarry Smith {
12365f70b478SJed Brown   TS_RK          *rk = (TS_RK*)ts->data;
12378370ee3bSLisandro Dalcin   PetscBool      iascii;
1238dfbe8321SBarry Smith   PetscErrorCode ierr;
12392cdf8120Sasbjorn 
1240e4dd521cSBarry Smith   PetscFunctionBegin;
1241251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
12428370ee3bSLisandro Dalcin   if (iascii) {
12439c334d8fSLisandro Dalcin     RKTableau       tab  = rk->tableau;
1244f68a32c8SEmil Constantinescu     TSRKType        rktype;
12450f7a28cdSStefano Zampini     const PetscReal *c;
12460f7a28cdSStefano Zampini     PetscInt        s;
1247f68a32c8SEmil Constantinescu     char            buf[512];
12480f7a28cdSStefano Zampini     PetscBool       FSAL;
12490f7a28cdSStefano Zampini 
1250f68a32c8SEmil Constantinescu     ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr);
12510f7a28cdSStefano Zampini     ierr = TSRKGetTableau(ts,&s,NULL,NULL,&c,NULL,NULL,NULL,&FSAL);CHKERRQ(ierr);
1252efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  RK type %s\n",rktype);CHKERRQ(ierr);
12530c37eb8eSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  Order: %D\n",tab->order);CHKERRQ(ierr);
12540f7a28cdSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"  FSAL property: %s\n",FSAL ? "yes" : "no");CHKERRQ(ierr);
12550f7a28cdSStefano Zampini     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",s,c);CHKERRQ(ierr);
1256f68a32c8SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa c = %s\n",buf);CHKERRQ(ierr);
12578370ee3bSLisandro Dalcin   }
1258f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1259f68a32c8SEmil Constantinescu }
1260f68a32c8SEmil Constantinescu 
1261f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer)
1262f68a32c8SEmil Constantinescu {
1263f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
12649c334d8fSLisandro Dalcin   TSAdapt        adapt;
1265f68a32c8SEmil Constantinescu 
1266f68a32c8SEmil Constantinescu   PetscFunctionBegin;
12679c334d8fSLisandro Dalcin   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
12689c334d8fSLisandro Dalcin   ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr);
1269f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1270f68a32c8SEmil Constantinescu }
1271f68a32c8SEmil Constantinescu 
12722ea87ba9SLisandro Dalcin /*@
12732ea87ba9SLisandro Dalcin   TSRKGetOrder - Get the order of RK scheme
12742ea87ba9SLisandro Dalcin 
12752ea87ba9SLisandro Dalcin   Not collective
12762ea87ba9SLisandro Dalcin 
12772ea87ba9SLisandro Dalcin   Input Parameter:
12782ea87ba9SLisandro Dalcin .  ts - timestepping context
12792ea87ba9SLisandro Dalcin 
12802ea87ba9SLisandro Dalcin   Output Parameter:
12812ea87ba9SLisandro Dalcin .  order - order of RK-scheme
12822ea87ba9SLisandro Dalcin 
12832ea87ba9SLisandro Dalcin   Level: intermediate
12842ea87ba9SLisandro Dalcin 
12852ea87ba9SLisandro Dalcin .seealso: TSRKGetType()
12862ea87ba9SLisandro Dalcin @*/
12872ea87ba9SLisandro Dalcin PetscErrorCode TSRKGetOrder(TS ts,PetscInt *order)
12882ea87ba9SLisandro Dalcin {
12892ea87ba9SLisandro Dalcin   PetscErrorCode ierr;
12902ea87ba9SLisandro Dalcin 
12912ea87ba9SLisandro Dalcin   PetscFunctionBegin;
12922ea87ba9SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
12932ea87ba9SLisandro Dalcin   PetscValidIntPointer(order,2);
12942ea87ba9SLisandro Dalcin   ierr = PetscUseMethod(ts,"TSRKGetOrder_C",(TS,PetscInt*),(ts,order));CHKERRQ(ierr);
12952ea87ba9SLisandro Dalcin   PetscFunctionReturn(0);
12962ea87ba9SLisandro Dalcin }
12972ea87ba9SLisandro Dalcin 
1298f68a32c8SEmil Constantinescu /*@C
1299f68a32c8SEmil Constantinescu   TSRKSetType - Set the type of RK scheme
1300f68a32c8SEmil Constantinescu 
1301f68a32c8SEmil Constantinescu   Logically collective
1302f68a32c8SEmil Constantinescu 
1303d8d19677SJose E. Roman   Input Parameters:
1304f68a32c8SEmil Constantinescu +  ts - timestepping context
1305f68a32c8SEmil Constantinescu -  rktype - type of RK-scheme
1306f68a32c8SEmil Constantinescu 
1307287c2655SBarry Smith   Options Database:
13089bd3e852SBarry Smith .   -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs>
1309287c2655SBarry Smith 
1310f68a32c8SEmil Constantinescu   Level: intermediate
1311f68a32c8SEmil Constantinescu 
1312e7685601SHong Zhang .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK2B, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS, TSRK6VR, TSRK7VR, TSRK8VR
1313f68a32c8SEmil Constantinescu @*/
1314f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype)
1315f68a32c8SEmil Constantinescu {
1316f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
1317f68a32c8SEmil Constantinescu 
1318f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1319f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1320b92453a8SLisandro Dalcin   PetscValidCharPointer(rktype,2);
1321f68a32c8SEmil Constantinescu   ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr);
1322f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1323f68a32c8SEmil Constantinescu }
1324f68a32c8SEmil Constantinescu 
1325f68a32c8SEmil Constantinescu /*@C
1326f68a32c8SEmil Constantinescu   TSRKGetType - Get the type of RK scheme
1327f68a32c8SEmil Constantinescu 
13282ea87ba9SLisandro Dalcin   Not collective
1329f68a32c8SEmil Constantinescu 
1330f68a32c8SEmil Constantinescu   Input Parameter:
1331f68a32c8SEmil Constantinescu .  ts - timestepping context
1332f68a32c8SEmil Constantinescu 
1333f68a32c8SEmil Constantinescu   Output Parameter:
1334f68a32c8SEmil Constantinescu .  rktype - type of RK-scheme
1335f68a32c8SEmil Constantinescu 
1336f68a32c8SEmil Constantinescu   Level: intermediate
1337f68a32c8SEmil Constantinescu 
13382ea87ba9SLisandro Dalcin .seealso: TSRKSetType()
1339f68a32c8SEmil Constantinescu @*/
1340f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype)
1341f68a32c8SEmil Constantinescu {
1342f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
1343f68a32c8SEmil Constantinescu 
1344f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1345f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1346f68a32c8SEmil Constantinescu   ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr);
1347f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1348f68a32c8SEmil Constantinescu }
1349f68a32c8SEmil Constantinescu 
13502ea87ba9SLisandro Dalcin static PetscErrorCode TSRKGetOrder_RK(TS ts,PetscInt *order)
13512ea87ba9SLisandro Dalcin {
13522ea87ba9SLisandro Dalcin   TS_RK *rk = (TS_RK*)ts->data;
13532ea87ba9SLisandro Dalcin 
13542ea87ba9SLisandro Dalcin   PetscFunctionBegin;
13552ea87ba9SLisandro Dalcin   *order = rk->tableau->order;
13562ea87ba9SLisandro Dalcin   PetscFunctionReturn(0);
13572ea87ba9SLisandro Dalcin }
13582ea87ba9SLisandro Dalcin 
1359560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype)
1360f68a32c8SEmil Constantinescu {
1361f68a32c8SEmil Constantinescu   TS_RK *rk = (TS_RK*)ts->data;
1362f68a32c8SEmil Constantinescu 
1363f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1364f68a32c8SEmil Constantinescu   *rktype = rk->tableau->name;
1365f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1366f68a32c8SEmil Constantinescu }
13672c9cb062Sluzhanghpp 
1368560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype)
1369f68a32c8SEmil Constantinescu {
1370f68a32c8SEmil Constantinescu   TS_RK          *rk = (TS_RK*)ts->data;
1371f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
1372f68a32c8SEmil Constantinescu   PetscBool      match;
1373f68a32c8SEmil Constantinescu   RKTableauLink  link;
1374f68a32c8SEmil Constantinescu 
1375f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1376f68a32c8SEmil Constantinescu   if (rk->tableau) {
1377f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr);
1378f68a32c8SEmil Constantinescu     if (match) PetscFunctionReturn(0);
1379f68a32c8SEmil Constantinescu   }
1380f68a32c8SEmil Constantinescu   for (link = RKTableauList; link; link=link->next) {
1381f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr);
1382f68a32c8SEmil Constantinescu     if (match) {
1383be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);}
1384f68a32c8SEmil Constantinescu       rk->tableau = &link->tab;
1385be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);}
13862ffb9264SLisandro Dalcin       ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
1387f68a32c8SEmil Constantinescu       PetscFunctionReturn(0);
1388f68a32c8SEmil Constantinescu     }
1389f68a32c8SEmil Constantinescu   }
1390*98921bdaSJacob Faibussowitsch   SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype);
1391e4dd521cSBarry Smith }
1392e4dd521cSBarry Smith 
1393ff22ae23SHong Zhang static PetscErrorCode  TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y)
1394ff22ae23SHong Zhang {
1395ff22ae23SHong Zhang   TS_RK *rk = (TS_RK*)ts->data;
1396ff22ae23SHong Zhang 
1397ff22ae23SHong Zhang   PetscFunctionBegin;
13980429704eSStefano Zampini   if (ns) *ns = rk->tableau->s;
1399d2daff3dSHong Zhang   if (Y)   *Y = rk->Y;
1400ff22ae23SHong Zhang   PetscFunctionReturn(0);
1401ff22ae23SHong Zhang }
1402ff22ae23SHong Zhang 
1403b3a6b972SJed Brown static PetscErrorCode TSDestroy_RK(TS ts)
1404b3a6b972SJed Brown {
1405b3a6b972SJed Brown   PetscErrorCode ierr;
1406b3a6b972SJed Brown 
1407b3a6b972SJed Brown   PetscFunctionBegin;
1408b3a6b972SJed Brown   ierr = TSReset_RK(ts);CHKERRQ(ierr);
1409b3a6b972SJed Brown   if (ts->dm) {
1410b3a6b972SJed Brown     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
1411b3a6b972SJed Brown     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
1412b3a6b972SJed Brown   }
1413b3a6b972SJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
14142ea87ba9SLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetOrder_C",NULL);CHKERRQ(ierr);
1415b3a6b972SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr);
1416b3a6b972SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr);
14170f7a28cdSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetTableau_C",NULL);CHKERRQ(ierr);
14180fe4d17eSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",NULL);CHKERRQ(ierr);
14190fe4d17eSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",NULL);CHKERRQ(ierr);
1420b3a6b972SJed Brown   PetscFunctionReturn(0);
1421b3a6b972SJed Brown }
1422ff22ae23SHong Zhang 
14233ca0882eSHong Zhang /*
14243ca0882eSHong Zhang   This defines the nonlinear equation that is to be solved with SNES
14253ca0882eSHong Zhang   We do not need to solve the equation; we just use SNES to approximate the Jacobian
14263ca0882eSHong Zhang */
14273ca0882eSHong Zhang static PetscErrorCode SNESTSFormFunction_RK(SNES snes,Vec x,Vec y,TS ts)
14283ca0882eSHong Zhang {
14293ca0882eSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
14303ca0882eSHong Zhang   DM             dm,dmsave;
14313ca0882eSHong Zhang   PetscErrorCode ierr;
14323ca0882eSHong Zhang 
14333ca0882eSHong Zhang   PetscFunctionBegin;
14343ca0882eSHong Zhang   ierr   = SNESGetDM(snes,&dm);CHKERRQ(ierr);
14353ca0882eSHong Zhang   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
14363ca0882eSHong Zhang   dmsave = ts->dm;
14373ca0882eSHong Zhang   ts->dm = dm;
14383ca0882eSHong Zhang   ierr   = TSComputeRHSFunction(ts,rk->stage_time,x,y);CHKERRQ(ierr);
14393ca0882eSHong Zhang   ts->dm = dmsave;
14403ca0882eSHong Zhang   PetscFunctionReturn(0);
14413ca0882eSHong Zhang }
14423ca0882eSHong Zhang 
14433ca0882eSHong Zhang static PetscErrorCode SNESTSFormJacobian_RK(SNES snes,Vec x,Mat A,Mat B,TS ts)
14443ca0882eSHong Zhang {
14453ca0882eSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
14463ca0882eSHong Zhang   DM             dm,dmsave;
14473ca0882eSHong Zhang   PetscErrorCode ierr;
14483ca0882eSHong Zhang 
14493ca0882eSHong Zhang   PetscFunctionBegin;
14503ca0882eSHong Zhang   ierr   = SNESGetDM(snes,&dm);CHKERRQ(ierr);
14513ca0882eSHong Zhang   dmsave = ts->dm;
14523ca0882eSHong Zhang   ts->dm = dm;
14533ca0882eSHong Zhang   ierr   = TSComputeRHSJacobian(ts,rk->stage_time,x,A,B);CHKERRQ(ierr);
14543ca0882eSHong Zhang   ts->dm = dmsave;
14553ca0882eSHong Zhang   PetscFunctionReturn(0);
14563ca0882eSHong Zhang }
14573ca0882eSHong Zhang 
145821052c0cSHong Zhang /*@C
145921052c0cSHong Zhang   TSRKSetMultirate - Use the interpolation-based multirate RK method
146021052c0cSHong Zhang 
146121052c0cSHong Zhang   Logically collective
146221052c0cSHong Zhang 
1463d8d19677SJose E. Roman   Input Parameters:
146421052c0cSHong Zhang +  ts - timestepping context
146521052c0cSHong 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
146621052c0cSHong Zhang 
146721052c0cSHong Zhang   Options Database:
146821052c0cSHong Zhang .   -ts_rk_multirate - <true,false>
146921052c0cSHong Zhang 
147021052c0cSHong Zhang   Notes:
147121052c0cSHong 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).
147221052c0cSHong Zhang 
147321052c0cSHong Zhang   Level: intermediate
147421052c0cSHong Zhang 
147521052c0cSHong Zhang .seealso: TSRKGetMultirate()
147621052c0cSHong Zhang @*/
147721052c0cSHong Zhang PetscErrorCode TSRKSetMultirate(TS ts,PetscBool use_multirate)
147821052c0cSHong Zhang {
1479f06fb28fSHong Zhang   PetscErrorCode ierr;
14807dbacdf2SHong Zhang 
14818a4be4dbSHong Zhang   PetscFunctionBegin;
1482f06fb28fSHong Zhang   ierr = PetscTryMethod(ts,"TSRKSetMultirate_C",(TS,PetscBool),(ts,use_multirate));CHKERRQ(ierr);
148321052c0cSHong Zhang   PetscFunctionReturn(0);
148421052c0cSHong Zhang }
148521052c0cSHong Zhang 
148621052c0cSHong Zhang /*@C
148721052c0cSHong Zhang   TSRKGetMultirate - Gets whether to Use the interpolation-based multirate RK method
148821052c0cSHong Zhang 
148921052c0cSHong Zhang   Not collective
149021052c0cSHong Zhang 
149121052c0cSHong Zhang   Input Parameter:
149221052c0cSHong Zhang .  ts - timestepping context
149321052c0cSHong Zhang 
149421052c0cSHong Zhang   Output Parameter:
149521052c0cSHong Zhang .  use_multirate - PETSC_TRUE if the multirate RK method is enabled, PETSC_FALSE otherwise
149621052c0cSHong Zhang 
149721052c0cSHong Zhang   Level: intermediate
149821052c0cSHong Zhang 
149921052c0cSHong Zhang .seealso: TSRKSetMultirate()
150021052c0cSHong Zhang @*/
150121052c0cSHong Zhang PetscErrorCode TSRKGetMultirate(TS ts,PetscBool *use_multirate)
150221052c0cSHong Zhang {
1503f06fb28fSHong Zhang   PetscErrorCode ierr;
15047dbacdf2SHong Zhang 
15057dbacdf2SHong Zhang   PetscFunctionBegin;
1506f06fb28fSHong Zhang   ierr = PetscUseMethod(ts,"TSRKGetMultirate_C",(TS,PetscBool*),(ts,use_multirate));CHKERRQ(ierr);
150721052c0cSHong Zhang   PetscFunctionReturn(0);
150821052c0cSHong Zhang }
150921052c0cSHong Zhang 
1510ebe3b25bSBarry Smith /*MC
1511f68a32c8SEmil Constantinescu       TSRK - ODE and DAE solver using Runge-Kutta schemes
1512ebe3b25bSBarry Smith 
1513f68a32c8SEmil Constantinescu   The user should provide the right hand side of the equation
1514f68a32c8SEmil Constantinescu   using TSSetRHSFunction().
1515ebe3b25bSBarry Smith 
1516f68a32c8SEmil Constantinescu   Notes:
151798164e13SLisandro Dalcin   The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type
1518ebe3b25bSBarry Smith 
1519d41222bbSBarry Smith   Level: beginner
1520d41222bbSBarry Smith 
15215d1808f1SStefano Zampini .seealso:  TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRK2D, TTSRK2E, TSRK3,
15220fe4d17eSHong Zhang            TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister(), TSRKSetMultirate(), TSRKGetMultirate()
1523ebe3b25bSBarry Smith 
1524ebe3b25bSBarry Smith M*/
15258cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
1526e4dd521cSBarry Smith {
15271566a47fSLisandro Dalcin   TS_RK          *rk;
1528dfbe8321SBarry Smith   PetscErrorCode ierr;
1529e4dd521cSBarry Smith 
1530e4dd521cSBarry Smith   PetscFunctionBegin;
1531f68a32c8SEmil Constantinescu   ierr = TSRKInitializePackage();CHKERRQ(ierr);
1532f68a32c8SEmil Constantinescu 
1533f68a32c8SEmil Constantinescu   ts->ops->reset          = TSReset_RK;
15345f70b478SJed Brown   ts->ops->destroy        = TSDestroy_RK;
15355f70b478SJed Brown   ts->ops->view           = TSView_RK;
1536f68a32c8SEmil Constantinescu   ts->ops->load           = TSLoad_RK;
1537f68a32c8SEmil Constantinescu   ts->ops->setup          = TSSetUp_RK;
1538f68a32c8SEmil Constantinescu   ts->ops->interpolate    = TSInterpolate_RK;
15392c9cb062Sluzhanghpp   ts->ops->step           = TSStep_RK;
1540f68a32c8SEmil Constantinescu   ts->ops->evaluatestep   = TSEvaluateStep_RK;
1541fecfb714SLisandro Dalcin   ts->ops->rollback       = TSRollBack_RK;
1542f68a32c8SEmil Constantinescu   ts->ops->setfromoptions = TSSetFromOptions_RK;
1543ff22ae23SHong Zhang   ts->ops->getstages      = TSGetStages_RK;
1544e4dd521cSBarry Smith 
15453ca0882eSHong Zhang   ts->ops->snesfunction    = SNESTSFormFunction_RK;
15463ca0882eSHong Zhang   ts->ops->snesjacobian    = SNESTSFormJacobian_RK;
15470f7a1166SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_RK;
154813af1a74SHong Zhang   ts->ops->adjointsetup    = TSAdjointSetUp_RK;
154913af1a74SHong Zhang   ts->ops->adjointstep     = TSAdjointStep_RK;
155013af1a74SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_RK;
15510f7a1166SHong Zhang 
155213af1a74SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_RK;
1553922a638cSHong Zhang   ts->ops->forwardsetup    = TSForwardSetUp_RK;
1554922a638cSHong Zhang   ts->ops->forwardreset    = TSForwardReset_RK;
1555922a638cSHong Zhang   ts->ops->forwardstep     = TSForwardStep_RK;
1556922a638cSHong Zhang   ts->ops->forwardgetstages= TSForwardGetStages_RK;
1557922a638cSHong Zhang 
15581566a47fSLisandro Dalcin   ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr);
15591566a47fSLisandro Dalcin   ts->data = (void*)rk;
1560e4dd521cSBarry Smith 
15612ea87ba9SLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetOrder_C",TSRKGetOrder_RK);CHKERRQ(ierr);
1562f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr);
1563f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr);
15640f7a28cdSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetTableau_C",TSRKGetTableau_RK);CHKERRQ(ierr);
156521052c0cSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",TSRKSetMultirate_RK);CHKERRQ(ierr);
156621052c0cSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",TSRKGetMultirate_RK);CHKERRQ(ierr);
1567be5899b3SLisandro Dalcin 
1568be5899b3SLisandro Dalcin   ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr);
156990b67129SHong Zhang   rk->dtratio = 1;
1570e4dd521cSBarry Smith   PetscFunctionReturn(0);
1571e4dd521cSBarry Smith }
1572