xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision d8d19677bbccf95218448bee62e6b87f4513e133)
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
352109b73fSDebojyoti Ghosh      TSRK2A - Second order RK scheme.
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
47f68a32c8SEmil Constantinescu      TSRK3 - Third order RK scheme.
48f68a32c8SEmil Constantinescu 
49f68a32c8SEmil Constantinescu      This method has three stages.
50f68a32c8SEmil Constantinescu 
51287c2655SBarry Smith      Options database:
529bd3e852SBarry Smith .     -ts_rk_type 3
53287c2655SBarry Smith 
54f68a32c8SEmil Constantinescu      Level: advanced
55f68a32c8SEmil Constantinescu 
56287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
57f68a32c8SEmil Constantinescu M*/
58f68a32c8SEmil Constantinescu /*MC
592109b73fSDebojyoti Ghosh      TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method.
602109b73fSDebojyoti Ghosh 
61d1905564SLisandro Dalcin      This method has four stages with the First Same As Last (FSAL) property.
622109b73fSDebojyoti Ghosh 
63287c2655SBarry Smith      Options database:
649bd3e852SBarry Smith .     -ts_rk_type 3bs
65287c2655SBarry Smith 
662109b73fSDebojyoti Ghosh      Level: advanced
672109b73fSDebojyoti Ghosh 
6898164e13SLisandro Dalcin      References: https://doi.org/10.1016/0893-9659(89)90079-7
69d1905564SLisandro Dalcin 
70287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
712109b73fSDebojyoti Ghosh M*/
722109b73fSDebojyoti Ghosh /*MC
73f68a32c8SEmil Constantinescu      TSRK4 - Fourth order RK scheme.
74f68a32c8SEmil Constantinescu 
752e698f8bSDebojyoti Ghosh      This is the classical Runge-Kutta method with four stages.
76f68a32c8SEmil Constantinescu 
77287c2655SBarry Smith      Options database:
789bd3e852SBarry Smith .     -ts_rk_type 4
79287c2655SBarry Smith 
80f68a32c8SEmil Constantinescu      Level: advanced
81f68a32c8SEmil Constantinescu 
82287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
83f68a32c8SEmil Constantinescu M*/
84f68a32c8SEmil Constantinescu /*MC
852e698f8bSDebojyoti Ghosh      TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method.
86f68a32c8SEmil Constantinescu 
87f68a32c8SEmil Constantinescu      This method has six stages.
88f68a32c8SEmil Constantinescu 
89287c2655SBarry Smith      Options database:
909bd3e852SBarry Smith .     -ts_rk_type 5f
91287c2655SBarry Smith 
92f68a32c8SEmil Constantinescu      Level: advanced
93f68a32c8SEmil Constantinescu 
94287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
95f68a32c8SEmil Constantinescu M*/
962109b73fSDebojyoti Ghosh /*MC
972e698f8bSDebojyoti Ghosh      TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method.
982109b73fSDebojyoti Ghosh 
99d1905564SLisandro Dalcin      This method has seven stages with the First Same As Last (FSAL) property.
1002109b73fSDebojyoti Ghosh 
101287c2655SBarry Smith      Options database:
1029bd3e852SBarry Smith .     -ts_rk_type 5dp
103287c2655SBarry Smith 
1042109b73fSDebojyoti Ghosh      Level: advanced
1052109b73fSDebojyoti Ghosh 
10698164e13SLisandro Dalcin      References: https://doi.org/10.1016/0771-050X(80)90013-3
107d1905564SLisandro Dalcin 
108287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
1092109b73fSDebojyoti Ghosh M*/
11005e23783SLisandro Dalcin /*MC
11105e23783SLisandro Dalcin      TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method.
11205e23783SLisandro Dalcin 
11305e23783SLisandro Dalcin      This method has eight stages with the First Same As Last (FSAL) property.
11405e23783SLisandro Dalcin 
115287c2655SBarry Smith      Options database:
1169bd3e852SBarry Smith .     -ts_rk_type 5bs
117287c2655SBarry Smith 
11805e23783SLisandro Dalcin      Level: advanced
11905e23783SLisandro Dalcin 
12098164e13SLisandro Dalcin      References: https://doi.org/10.1016/0898-1221(96)00141-1
12105e23783SLisandro Dalcin 
122287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
12305e23783SLisandro Dalcin M*/
12437acaa02SHendrik Ranocha /*MC
12537acaa02SHendrik Ranocha      TSRK6VR - Sixth order robust Verner RK scheme with fifth order embedded method.
12637acaa02SHendrik Ranocha 
12737acaa02SHendrik Ranocha      This method has nine stages with the First Same As Last (FSAL) property.
12837acaa02SHendrik Ranocha 
12937acaa02SHendrik Ranocha      Options database:
13037acaa02SHendrik Ranocha .     -ts_rk_type 6vr
13137acaa02SHendrik Ranocha 
13237acaa02SHendrik Ranocha      Level: advanced
13337acaa02SHendrik Ranocha 
13437acaa02SHendrik Ranocha      References: http://people.math.sfu.ca/~jverner/RKV65.IIIXb.Robust.00010102836.081204.CoeffsOnlyRAT
13537acaa02SHendrik Ranocha 
13637acaa02SHendrik Ranocha .seealso: TSRK, TSRKType, TSRKSetType()
13737acaa02SHendrik Ranocha M*/
13837acaa02SHendrik Ranocha /*MC
13937acaa02SHendrik Ranocha      TSRK7VR - Seventh order robust Verner RK scheme with sixth order embedded method.
14037acaa02SHendrik Ranocha 
1418f3d7ee2SCarsten Uphoff      This method has ten stages.
14237acaa02SHendrik Ranocha 
14337acaa02SHendrik Ranocha      Options database:
14437acaa02SHendrik Ranocha .     -ts_rk_type 7vr
14537acaa02SHendrik Ranocha 
14637acaa02SHendrik Ranocha      Level: advanced
14737acaa02SHendrik Ranocha 
14837acaa02SHendrik Ranocha      References: http://people.math.sfu.ca/~jverner/RKV76.IIa.Robust.000027015646.081206.CoeffsOnlyRAT
14937acaa02SHendrik Ranocha 
15037acaa02SHendrik Ranocha .seealso: TSRK, TSRKType, TSRKSetType()
15137acaa02SHendrik Ranocha M*/
15237acaa02SHendrik Ranocha /*MC
15337acaa02SHendrik Ranocha      TSRK8VR - Eigth order robust Verner RK scheme with seventh order embedded method.
15437acaa02SHendrik Ranocha 
1558f3d7ee2SCarsten Uphoff      This method has thirteen stages.
15637acaa02SHendrik Ranocha 
15737acaa02SHendrik Ranocha      Options database:
15837acaa02SHendrik Ranocha .     -ts_rk_type 8vr
15937acaa02SHendrik Ranocha 
16037acaa02SHendrik Ranocha      Level: advanced
16137acaa02SHendrik Ranocha 
16237acaa02SHendrik Ranocha      References: http://people.math.sfu.ca/~jverner/RKV87.IIa.Robust.00000754677.081208.CoeffsOnlyRATandFLOAT
16337acaa02SHendrik Ranocha 
16437acaa02SHendrik Ranocha .seealso: TSRK, TSRKType, TSRKSetType()
16537acaa02SHendrik Ranocha M*/
166262ff7bbSBarry Smith 
167f68a32c8SEmil Constantinescu /*@C
168f68a32c8SEmil Constantinescu   TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK
169262ff7bbSBarry Smith 
170f68a32c8SEmil Constantinescu   Not Collective, but should be called by all processes which will need the schemes to be registered
171262ff7bbSBarry Smith 
172f68a32c8SEmil Constantinescu   Level: advanced
173262ff7bbSBarry Smith 
174f68a32c8SEmil Constantinescu .seealso:  TSRKRegisterDestroy()
175262ff7bbSBarry Smith @*/
176f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void)
177262ff7bbSBarry Smith {
1784ac538c5SBarry Smith   PetscErrorCode ierr;
179262ff7bbSBarry Smith 
180262ff7bbSBarry Smith   PetscFunctionBegin;
181f68a32c8SEmil Constantinescu   if (TSRKRegisterAllCalled) PetscFunctionReturn(0);
182f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_TRUE;
183e4dd521cSBarry Smith 
1844e82626cSLisandro Dalcin #define RC PetscRealConstant
185e4dd521cSBarry Smith   {
186f68a32c8SEmil Constantinescu     const PetscReal
1874e82626cSLisandro Dalcin       A[1][1] = {{0}},
1884e82626cSLisandro Dalcin       b[1]    = {RC(1.0)};
1893a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
190e4dd521cSBarry Smith   }
191db046809SBarry Smith   {
192f68a32c8SEmil Constantinescu     const PetscReal
1934e82626cSLisandro Dalcin       A[2][2]   = {{0,0},
1944e82626cSLisandro Dalcin                    {RC(1.0),0}},
1954e82626cSLisandro Dalcin       b[2]      =  {RC(0.5),RC(0.5)},
1964e82626cSLisandro Dalcin       bembed[2] =  {RC(1.0),0};
1973a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
198db046809SBarry Smith   }
199f68a32c8SEmil Constantinescu   {
200f68a32c8SEmil Constantinescu     const PetscReal
20117f6384fSBarry Smith       A[3][3] = {{0,0,0},
2024e82626cSLisandro Dalcin                  {RC(2.0)/RC(3.0),0,0},
2034e82626cSLisandro Dalcin                  {RC(-1.0)/RC(3.0),RC(1.0),0}},
2044e82626cSLisandro Dalcin       b[3]    =  {RC(0.25),RC(0.5),RC(0.25)};
2053a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
206fdefd5e5SDebojyoti Ghosh   }
207fdefd5e5SDebojyoti Ghosh   {
208fdefd5e5SDebojyoti Ghosh     const PetscReal
20917f6384fSBarry Smith       A[4][4]   = {{0,0,0,0},
2104e82626cSLisandro Dalcin                    {RC(1.0)/RC(2.0),0,0,0},
2114e82626cSLisandro Dalcin                    {0,RC(3.0)/RC(4.0),0,0},
2124e82626cSLisandro Dalcin                    {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}},
2134e82626cSLisandro Dalcin       b[4]      =  {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0},
2144e82626cSLisandro Dalcin       bembed[4] =  {RC(7.0)/RC(24.0),RC(1.0)/RC(4.0),RC(1.0)/RC(3.0),RC(1.0)/RC(8.0)};
2153a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
216db046809SBarry Smith   }
217f68a32c8SEmil Constantinescu   {
218f68a32c8SEmil Constantinescu     const PetscReal
219f68a32c8SEmil Constantinescu       A[4][4] = {{0,0,0,0},
2204e82626cSLisandro Dalcin                  {RC(0.5),0,0,0},
2214e82626cSLisandro Dalcin                  {0,RC(0.5),0,0},
2224e82626cSLisandro Dalcin                  {0,0,RC(1.0),0}},
2234e82626cSLisandro Dalcin       b[4]    =  {RC(1.0)/RC(6.0),RC(1.0)/RC(3.0),RC(1.0)/RC(3.0),RC(1.0)/RC(6.0)};
2243a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
225db046809SBarry Smith   }
226f68a32c8SEmil Constantinescu   {
227f68a32c8SEmil Constantinescu     const PetscReal
22817f6384fSBarry Smith       A[6][6]   = {{0,0,0,0,0,0},
2294e82626cSLisandro Dalcin                    {RC(0.25),0,0,0,0,0},
2304e82626cSLisandro Dalcin                    {RC(3.0)/RC(32.0),RC(9.0)/RC(32.0),0,0,0,0},
2314e82626cSLisandro Dalcin                    {RC(1932.0)/RC(2197.0),RC(-7200.0)/RC(2197.0),RC(7296.0)/RC(2197.0),0,0,0},
2324e82626cSLisandro Dalcin                    {RC(439.0)/RC(216.0),RC(-8.0),RC(3680.0)/RC(513.0),RC(-845.0)/RC(4104.0),0,0},
2334e82626cSLisandro Dalcin                    {RC(-8.0)/RC(27.0),RC(2.0),RC(-3544.0)/RC(2565.0),RC(1859.0)/RC(4104.0),RC(-11.0)/RC(40.0),0}},
2344e82626cSLisandro Dalcin       b[6]      =  {RC(16.0)/RC(135.0),0,RC(6656.0)/RC(12825.0),RC(28561.0)/RC(56430.0),RC(-9.0)/RC(50.0),RC(2.0)/RC(55.0)},
2354e82626cSLisandro Dalcin       bembed[6] =  {RC(25.0)/RC(216.0),0,RC(1408.0)/RC(2565.0),RC(2197.0)/RC(4104.0),RC(-1.0)/RC(5.0),0};
2363a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
237fdefd5e5SDebojyoti Ghosh   }
238fdefd5e5SDebojyoti Ghosh   {
239fdefd5e5SDebojyoti Ghosh     const PetscReal
24017f6384fSBarry Smith       A[7][7]       = {{0,0,0,0,0,0,0},
2414e82626cSLisandro Dalcin                        {RC(1.0)/RC(5.0),0,0,0,0,0,0},
2424e82626cSLisandro Dalcin                        {RC(3.0)/RC(40.0),RC(9.0)/RC(40.0),0,0,0,0,0},
2434e82626cSLisandro Dalcin                        {RC(44.0)/RC(45.0),RC(-56.0)/RC(15.0),RC(32.0)/RC(9.0),0,0,0,0},
2444e82626cSLisandro Dalcin                        {RC(19372.0)/RC(6561.0),RC(-25360.0)/RC(2187.0),RC(64448.0)/RC(6561.0),RC(-212.0)/RC(729.0),0,0,0},
2454e82626cSLisandro Dalcin                        {RC(9017.0)/RC(3168.0),RC(-355.0)/RC(33.0),RC(46732.0)/RC(5247.0),RC(49.0)/RC(176.0),RC(-5103.0)/RC(18656.0),0,0},
2464e82626cSLisandro Dalcin                        {RC(35.0)/RC(384.0),0,RC(500.0)/RC(1113.0),RC(125.0)/RC(192.0),RC(-2187.0)/RC(6784.0),RC(11.0)/RC(84.0),0}},
2474e82626cSLisandro Dalcin       b[7]          =  {RC(35.0)/RC(384.0),0,RC(500.0)/RC(1113.0),RC(125.0)/RC(192.0),RC(-2187.0)/RC(6784.0),RC(11.0)/RC(84.0),0},
248a685a763Sluzhanghpp       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)},
249a685a763Sluzhanghpp       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)},
250a685a763Sluzhanghpp                        {0,0,0,0,0},
251a685a763Sluzhanghpp                        {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)},
252a685a763Sluzhanghpp                        {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)},
253a685a763Sluzhanghpp                        {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)},
254a685a763Sluzhanghpp                        {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)},
255a685a763Sluzhanghpp                        {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)}};
256a685a763Sluzhanghpp       ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,5,binterp[0]);CHKERRQ(ierr);
257f68a32c8SEmil Constantinescu   }
25805e23783SLisandro Dalcin   {
25905e23783SLisandro Dalcin     const PetscReal
26017f6384fSBarry Smith       A[8][8]   = {{0,0,0,0,0,0,0,0},
2614e82626cSLisandro Dalcin                    {RC(1.0)/RC(6.0),0,0,0,0,0,0,0},
2624e82626cSLisandro Dalcin                    {RC(2.0)/RC(27.0),RC(4.0)/RC(27.0),0,0,0,0,0,0},
2634e82626cSLisandro Dalcin                    {RC(183.0)/RC(1372.0),RC(-162.0)/RC(343.0),RC(1053.0)/RC(1372.0),0,0,0,0,0},
2644e82626cSLisandro 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},
2654e82626cSLisandro 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},
2664e82626cSLisandro 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},
2674e82626cSLisandro 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}},
2684e82626cSLisandro 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},
2694e82626cSLisandro 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)};
27005e23783SLisandro Dalcin     ierr = TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
27105e23783SLisandro Dalcin   }
27237acaa02SHendrik Ranocha   {
27337acaa02SHendrik Ranocha     const PetscReal
27437acaa02SHendrik Ranocha       A[9][9]   = {{0,0,0,0,0,0,0,0,0},
27563fe90e8SHendrik Ranocha                    {RC(1.8000000000000000000000000000000000000000e-01),0,0,0,0,0,0,0,0},
27663fe90e8SHendrik Ranocha                    {RC(8.9506172839506172839506172839506172839506e-02),RC(7.7160493827160493827160493827160493827160e-02),0,0,0,0,0,0,0},
27763fe90e8SHendrik Ranocha                    {RC(6.2500000000000000000000000000000000000000e-02),0,RC(1.8750000000000000000000000000000000000000e-01),0,0,0,0,0,0},
27863fe90e8SHendrik Ranocha                    {RC(3.1651600000000000000000000000000000000000e-01),0,RC(-1.0449480000000000000000000000000000000000e+00),RC(1.2584320000000000000000000000000000000000e+00),0,0,0,0,0},
27963fe90e8SHendrik Ranocha                    {RC(2.7232612736485626257225065566674305502508e-01),0,RC(-8.2513360323886639676113360323886639676113e-01),RC(1.0480917678812415654520917678812415654521e+00),RC(1.0471570799276856873679117969088177628396e-01),0,0,0,0},
28063fe90e8SHendrik 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},
28163fe90e8SHendrik 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},
28263fe90e8SHendrik 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}},
28363fe90e8SHendrik 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},
28463fe90e8SHendrik 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)};
28537acaa02SHendrik Ranocha     ierr = TSRKRegister(TSRK6VR,6,9,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
28637acaa02SHendrik Ranocha   }
28737acaa02SHendrik Ranocha   {
28837acaa02SHendrik Ranocha     const PetscReal
28937acaa02SHendrik Ranocha       A[10][10]  = {{0,0,0,0,0,0,0,0,0,0},
29063fe90e8SHendrik Ranocha                     {RC(5.0000000000000000000000000000000000000000e-03),0,0,0,0,0,0,0,0,0},
29163fe90e8SHendrik Ranocha                     {RC(-1.0767901234567901234567901234567901234568e+00),RC(1.1856790123456790123456790123456790123457e+00),0,0,0,0,0,0,0,0},
29263fe90e8SHendrik Ranocha                     {RC(4.0833333333333333333333333333333333333333e-02),0,RC(1.2250000000000000000000000000000000000000e-01),0,0,0,0,0,0,0},
29363fe90e8SHendrik Ranocha                     {RC(6.3607142857142857142857142857142857142857e-01),0,RC(-2.4444642857142857142857142857142857142857e+00),RC(2.2633928571428571428571428571428571428571e+00),0,0,0,0,0,0},
29463fe90e8SHendrik Ranocha                     {RC(-2.5351211079349245229256383554660215487207e+00),0,RC(1.0299374654449267920438514460756024913612e+01),RC(-7.9513032885990579949493217458266876536482e+00),RC(7.9301148923100592201226014271115261823800e-01),0,0,0,0,0},
29563fe90e8SHendrik 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},
29663fe90e8SHendrik 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},
29763fe90e8SHendrik 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},
29863fe90e8SHendrik 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}},
29963fe90e8SHendrik 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},
30063fe90e8SHendrik 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)};
30137acaa02SHendrik Ranocha     ierr = TSRKRegister(TSRK7VR,7,10,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
30237acaa02SHendrik Ranocha   }
30337acaa02SHendrik Ranocha   {
30437acaa02SHendrik Ranocha     const PetscReal
30537acaa02SHendrik Ranocha       A[13][13]  = {{0,0,0,0,0,0,0,0,0,0,0,0,0},
30663fe90e8SHendrik Ranocha                     {RC(2.5000000000000000000000000000000000000000e-01),0,0,0,0,0,0,0,0,0,0,0,0},
30763fe90e8SHendrik Ranocha                     {RC(8.7400846504915232052686327594877411977046e-02),RC(2.5487604938654321753087950620345685135815e-02),0,0,0,0,0,0,0,0,0,0,0},
30863fe90e8SHendrik Ranocha                     {RC(4.2333169291338582677165354330708661417323e-02),0,RC(1.2699950787401574803149606299212598425197e-01),0,0,0,0,0,0,0,0,0,0},
30963fe90e8SHendrik Ranocha                     {RC(4.2609505888742261494881445237572274090942e-01),0,RC(-1.5987952846591523265427733230657181117089e+00),RC(1.5967002257717297115939588706899953707994e+00),0,0,0,0,0,0,0,0,0},
31063fe90e8SHendrik Ranocha                     {RC(5.0719337296713929515090618138513639239329e-02),0,0,RC(2.5433377264600407582754714408877778031369e-01),RC(2.0394689005728199465736223777270858044698e-01),0,0,0,0,0,0,0,0},
31163fe90e8SHendrik 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},
31263fe90e8SHendrik 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},
31363fe90e8SHendrik 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},
31463fe90e8SHendrik 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},
31563fe90e8SHendrik 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},
31663fe90e8SHendrik 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},
31763fe90e8SHendrik 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}},
31863fe90e8SHendrik 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},
31963fe90e8SHendrik 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)};
32037acaa02SHendrik Ranocha     ierr = TSRKRegister(TSRK8VR,8,13,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
32137acaa02SHendrik Ranocha   }
3224e82626cSLisandro Dalcin #undef RC
323db046809SBarry Smith   PetscFunctionReturn(0);
324db046809SBarry Smith }
325db046809SBarry Smith 
326f68a32c8SEmil Constantinescu /*@C
327f68a32c8SEmil Constantinescu    TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister().
328f68a32c8SEmil Constantinescu 
329f68a32c8SEmil Constantinescu    Not Collective
330f68a32c8SEmil Constantinescu 
331f68a32c8SEmil Constantinescu    Level: advanced
332f68a32c8SEmil Constantinescu 
333f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll()
334f68a32c8SEmil Constantinescu @*/
335f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void)
336e4dd521cSBarry Smith {
337f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
338f68a32c8SEmil Constantinescu   RKTableauLink  link;
339f68a32c8SEmil Constantinescu 
340f68a32c8SEmil Constantinescu   PetscFunctionBegin;
341f68a32c8SEmil Constantinescu   while ((link = RKTableauList)) {
342f68a32c8SEmil Constantinescu     RKTableau t = &link->tab;
343f68a32c8SEmil Constantinescu     RKTableauList = link->next;
344f68a32c8SEmil Constantinescu     ierr = PetscFree3(t->A,t->b,t->c);CHKERRQ(ierr);
345f68a32c8SEmil Constantinescu     ierr = PetscFree(t->bembed);CHKERRQ(ierr);
346f68a32c8SEmil Constantinescu     ierr = PetscFree(t->binterp);CHKERRQ(ierr);
347f68a32c8SEmil Constantinescu     ierr = PetscFree(t->name);CHKERRQ(ierr);
348f68a32c8SEmil Constantinescu     ierr = PetscFree(link);CHKERRQ(ierr);
349f68a32c8SEmil Constantinescu   }
350f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_FALSE;
351f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
352f68a32c8SEmil Constantinescu }
353f68a32c8SEmil Constantinescu 
354f68a32c8SEmil Constantinescu /*@C
355f68a32c8SEmil Constantinescu   TSRKInitializePackage - This function initializes everything in the TSRK package. It is called
3568a690491SBarry Smith   from TSInitializePackage().
357f68a32c8SEmil Constantinescu 
358f68a32c8SEmil Constantinescu   Level: developer
359f68a32c8SEmil Constantinescu 
360f68a32c8SEmil Constantinescu .seealso: PetscInitialize()
361f68a32c8SEmil Constantinescu @*/
362f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void)
363f68a32c8SEmil Constantinescu {
364186e87acSLisandro Dalcin   PetscErrorCode ierr;
365e4dd521cSBarry Smith 
366e4dd521cSBarry Smith   PetscFunctionBegin;
367f68a32c8SEmil Constantinescu   if (TSRKPackageInitialized) PetscFunctionReturn(0);
368f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_TRUE;
369f68a32c8SEmil Constantinescu   ierr = TSRKRegisterAll();CHKERRQ(ierr);
370f68a32c8SEmil Constantinescu   ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr);
371f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
372f68a32c8SEmil Constantinescu }
373186e87acSLisandro Dalcin 
374f68a32c8SEmil Constantinescu /*@C
375f68a32c8SEmil Constantinescu   TSRKFinalizePackage - This function destroys everything in the TSRK package. It is
376f68a32c8SEmil Constantinescu   called from PetscFinalize().
377186e87acSLisandro Dalcin 
378f68a32c8SEmil Constantinescu   Level: developer
379f68a32c8SEmil Constantinescu 
380f68a32c8SEmil Constantinescu .seealso: PetscFinalize()
381f68a32c8SEmil Constantinescu @*/
382f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void)
383f68a32c8SEmil Constantinescu {
384f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
385f68a32c8SEmil Constantinescu 
386f68a32c8SEmil Constantinescu   PetscFunctionBegin;
387f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_FALSE;
388f68a32c8SEmil Constantinescu   ierr = TSRKRegisterDestroy();CHKERRQ(ierr);
389f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
390f68a32c8SEmil Constantinescu }
391f68a32c8SEmil Constantinescu 
392f68a32c8SEmil Constantinescu /*@C
393f68a32c8SEmil Constantinescu    TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
394f68a32c8SEmil Constantinescu 
395f68a32c8SEmil Constantinescu    Not Collective, but the same schemes should be registered on all processes on which they will be used
396f68a32c8SEmil Constantinescu 
397f68a32c8SEmil Constantinescu    Input Parameters:
398f68a32c8SEmil Constantinescu +  name - identifier for method
399f68a32c8SEmil Constantinescu .  order - approximation order of method
400f68a32c8SEmil Constantinescu .  s - number of stages, this is the dimension of the matrices below
401f68a32c8SEmil Constantinescu .  A - stage coefficients (dimension s*s, row-major)
402f68a32c8SEmil Constantinescu .  b - step completion table (dimension s; NULL to use last row of A)
403f68a32c8SEmil Constantinescu .  c - abscissa (dimension s; NULL to use row sums of A)
404f68a32c8SEmil Constantinescu .  bembed - completion table for embedded method (dimension s; NULL if not available)
4053a8a9803SLisandro Dalcin .  p - Order of the interpolation scheme, equal to the number of columns of binterp
4063a8a9803SLisandro Dalcin -  binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1)
407f68a32c8SEmil Constantinescu 
408f68a32c8SEmil Constantinescu    Notes:
409f68a32c8SEmil Constantinescu    Several RK methods are provided, this function is only needed to create new methods.
410f68a32c8SEmil Constantinescu 
411f68a32c8SEmil Constantinescu    Level: advanced
412f68a32c8SEmil Constantinescu 
413f68a32c8SEmil Constantinescu .seealso: TSRK
414f68a32c8SEmil Constantinescu @*/
415f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s,
416f68a32c8SEmil Constantinescu                             const PetscReal A[],const PetscReal b[],const PetscReal c[],
4173a8a9803SLisandro Dalcin                             const PetscReal bembed[],PetscInt p,const PetscReal binterp[])
418f68a32c8SEmil Constantinescu {
419f68a32c8SEmil Constantinescu   PetscErrorCode  ierr;
420f68a32c8SEmil Constantinescu   RKTableauLink   link;
421f68a32c8SEmil Constantinescu   RKTableau       t;
422f68a32c8SEmil Constantinescu   PetscInt        i,j;
423f68a32c8SEmil Constantinescu 
424f68a32c8SEmil Constantinescu   PetscFunctionBegin;
4253a8a9803SLisandro Dalcin   PetscValidCharPointer(name,1);
4263a8a9803SLisandro Dalcin   PetscValidRealPointer(A,4);
4273a8a9803SLisandro Dalcin   if (b) PetscValidRealPointer(b,5);
4283a8a9803SLisandro Dalcin   if (c) PetscValidRealPointer(c,6);
4293a8a9803SLisandro Dalcin   if (bembed) PetscValidRealPointer(bembed,7);
4303a8a9803SLisandro Dalcin   if (binterp || p > 1) PetscValidRealPointer(binterp,9);
4313a8a9803SLisandro Dalcin 
4321d36bdfdSBarry Smith   ierr = TSRKInitializePackage();CHKERRQ(ierr);
43395dccacaSBarry Smith   ierr = PetscNew(&link);CHKERRQ(ierr);
434f68a32c8SEmil Constantinescu   t = &link->tab;
4353a8a9803SLisandro Dalcin 
436f68a32c8SEmil Constantinescu   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
437f68a32c8SEmil Constantinescu   t->order = order;
438f68a32c8SEmil Constantinescu   t->s = s;
439dcca6d9dSJed Brown   ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr);
440580bdb30SBarry Smith   ierr = PetscArraycpy(t->A,A,s*s);CHKERRQ(ierr);
441580bdb30SBarry Smith   if (b)  { ierr = PetscArraycpy(t->b,b,s);CHKERRQ(ierr); }
442f68a32c8SEmil Constantinescu   else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i];
443580bdb30SBarry Smith   if (c)  { ierr = PetscArraycpy(t->c,c,s);CHKERRQ(ierr); }
444f68a32c8SEmil 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];
445d760c35bSDebojyoti Ghosh   t->FSAL = PETSC_TRUE;
446d760c35bSDebojyoti Ghosh   for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE;
4473a8a9803SLisandro Dalcin 
448f68a32c8SEmil Constantinescu   if (bembed) {
449785e854fSJed Brown     ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr);
450580bdb30SBarry Smith     ierr = PetscArraycpy(t->bembed,bembed,s);CHKERRQ(ierr);
451f68a32c8SEmil Constantinescu   }
452f68a32c8SEmil Constantinescu 
4533a8a9803SLisandro Dalcin   if (!binterp) { p = 1; binterp = t->b; }
4543a8a9803SLisandro Dalcin   t->p = p;
4553a8a9803SLisandro Dalcin   ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr);
456580bdb30SBarry Smith   ierr = PetscArraycpy(t->binterp,binterp,s*p);CHKERRQ(ierr);
4573a8a9803SLisandro Dalcin 
458f68a32c8SEmil Constantinescu   link->next = RKTableauList;
459f68a32c8SEmil Constantinescu   RKTableauList = link;
460f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
461f68a32c8SEmil Constantinescu }
462f68a32c8SEmil Constantinescu 
4630f7a28cdSStefano Zampini PetscErrorCode TSRKGetTableau_RK(TS ts, PetscInt *s, const PetscReal **A, const PetscReal **b, const PetscReal **c, const PetscReal **bembed,
4640f7a28cdSStefano Zampini                                         PetscInt *p, const PetscReal **binterp, PetscBool *FSAL)
4650f7a28cdSStefano Zampini {
4660f7a28cdSStefano Zampini   TS_RK     *rk   = (TS_RK*)ts->data;
4670f7a28cdSStefano Zampini   RKTableau tab  = rk->tableau;
4680f7a28cdSStefano Zampini 
4690f7a28cdSStefano Zampini   PetscFunctionBegin;
4700f7a28cdSStefano Zampini   if (s) *s = tab->s;
4710f7a28cdSStefano Zampini   if (A) *A = tab->A;
4720f7a28cdSStefano Zampini   if (b) *b = tab->b;
4730f7a28cdSStefano Zampini   if (c) *c = tab->c;
4740f7a28cdSStefano Zampini   if (bembed) *bembed = tab->bembed;
4750f7a28cdSStefano Zampini   if (p) *p = tab->p;
4760f7a28cdSStefano Zampini   if (binterp) *binterp = tab->binterp;
4770f7a28cdSStefano Zampini   if (FSAL) *FSAL = tab->FSAL;
4780f7a28cdSStefano Zampini   PetscFunctionReturn(0);
4790f7a28cdSStefano Zampini }
4800f7a28cdSStefano Zampini 
4810f7a28cdSStefano Zampini /*@C
4820f7a28cdSStefano Zampini    TSRKGetTableau - Get info on the RK tableau
4830f7a28cdSStefano Zampini 
4840f7a28cdSStefano Zampini    Not Collective
4850f7a28cdSStefano Zampini 
4860f7a28cdSStefano Zampini    Input Parameters:
4870f7a28cdSStefano Zampini .  ts - timestepping context
4880f7a28cdSStefano Zampini 
4890f7a28cdSStefano Zampini    Output Parameters:
4900f7a28cdSStefano Zampini +  s - number of stages, this is the dimension of the matrices below
4910f7a28cdSStefano Zampini .  A - stage coefficients (dimension s*s, row-major)
4920f7a28cdSStefano Zampini .  b - step completion table (dimension s)
4930f7a28cdSStefano Zampini .  c - abscissa (dimension s)
4940f7a28cdSStefano Zampini .  bembed - completion table for embedded method (dimension s; NULL if not available)
4950f7a28cdSStefano Zampini .  p - Order of the interpolation scheme, equal to the number of columns of binterp
4960f7a28cdSStefano Zampini .  binterp - Coefficients of the interpolation formula (dimension s*p)
4970f7a28cdSStefano Zampini -  FSAL - wheather or not the scheme has the First Same As Last property
4980f7a28cdSStefano Zampini 
4990f7a28cdSStefano Zampini    Level: developer
5000f7a28cdSStefano Zampini 
5010f7a28cdSStefano Zampini .seealso: TSRK
5020f7a28cdSStefano Zampini @*/
5030f7a28cdSStefano Zampini PetscErrorCode TSRKGetTableau(TS ts, PetscInt *s, const PetscReal **A, const PetscReal **b, const PetscReal **c, const PetscReal **bembed,
5040f7a28cdSStefano Zampini                                      PetscInt *p, const PetscReal **binterp, PetscBool *FSAL)
5050f7a28cdSStefano Zampini {
5060f7a28cdSStefano Zampini   PetscErrorCode ierr;
5070f7a28cdSStefano Zampini 
5080f7a28cdSStefano Zampini   PetscFunctionBegin;
5090f7a28cdSStefano Zampini   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5100f7a28cdSStefano Zampini   ierr = PetscUseMethod(ts,"TSRKGetTableau_C",(TS,PetscInt*,const PetscReal**,const PetscReal**,const PetscReal**,const PetscReal**,
5110f7a28cdSStefano Zampini                                                   PetscInt*,const PetscReal**,PetscBool*),(ts,s,A,b,c,bembed,p,binterp,FSAL));CHKERRQ(ierr);
5120f7a28cdSStefano Zampini   PetscFunctionReturn(0);
5130f7a28cdSStefano Zampini }
5140f7a28cdSStefano Zampini 
515e4dd521cSBarry Smith /*
5162c9cb062Sluzhanghpp  This is for single-step RK method
517f68a32c8SEmil Constantinescu  The step completion formula is
518e4dd521cSBarry Smith 
519f68a32c8SEmil Constantinescu  x1 = x0 + h b^T YdotRHS
520f68a32c8SEmil Constantinescu 
521f68a32c8SEmil Constantinescu  This function can be called before or after ts->vec_sol has been updated.
522f68a32c8SEmil Constantinescu  Suppose we have a completion formula (b) and an embedded formula (be) of different order.
523f68a32c8SEmil Constantinescu  We can write
524f68a32c8SEmil Constantinescu 
525f68a32c8SEmil Constantinescu  x1e = x0 + h be^T YdotRHS
526f68a32c8SEmil Constantinescu      = x1 - h b^T YdotRHS + h be^T YdotRHS
527f68a32c8SEmil Constantinescu      = x1 + h (be - b)^T YdotRHS
528f68a32c8SEmil Constantinescu 
529f68a32c8SEmil Constantinescu  so we can evaluate the method with different order even after the step has been optimistically completed.
530e4dd521cSBarry Smith */
531f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done)
532f68a32c8SEmil Constantinescu {
533f68a32c8SEmil Constantinescu   TS_RK          *rk   = (TS_RK*)ts->data;
534f68a32c8SEmil Constantinescu   RKTableau      tab  = rk->tableau;
535f68a32c8SEmil Constantinescu   PetscScalar    *w    = rk->work;
536f68a32c8SEmil Constantinescu   PetscReal      h;
537f68a32c8SEmil Constantinescu   PetscInt       s    = tab->s,j;
538f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
539f68a32c8SEmil Constantinescu 
540f68a32c8SEmil Constantinescu   PetscFunctionBegin;
541f68a32c8SEmil Constantinescu   switch (rk->status) {
542f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
543f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
544f68a32c8SEmil Constantinescu     h = ts->time_step; break;
545f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
546be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
547f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
548f68a32c8SEmil Constantinescu   }
549f68a32c8SEmil Constantinescu   if (order == tab->order) {
550f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) {
551f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
55290b67129SHong Zhang       for (j=0; j<s; j++) w[j] = h*tab->b[j]/rk->dtratio;
553f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
554f68a32c8SEmil Constantinescu     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
555f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
556f68a32c8SEmil Constantinescu   } else if (order == tab->order-1) {
557f68a32c8SEmil Constantinescu     if (!tab->bembed) goto unavailable;
558f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) { /*Complete with the embedded method (be)*/
559f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
560f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
561f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
562f68a32c8SEmil Constantinescu     } else {  /*Rollback and re-complete using (be-b) */
563f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
564f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
565f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
566f68a32c8SEmil Constantinescu     }
567f68a32c8SEmil Constantinescu     if (done) *done = PETSC_TRUE;
568f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
569f68a32c8SEmil Constantinescu   }
570f68a32c8SEmil Constantinescu unavailable:
571f68a32c8SEmil Constantinescu   if (done) *done = PETSC_FALSE;
572a7fac7c2SEmil Constantinescu   else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"RK '%s' of order %D cannot evaluate step at order %D. Consider using -ts_adapt_type none or a different method that has an embedded estimate.",tab->name,tab->order,order);
573f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
574f68a32c8SEmil Constantinescu }
575f68a32c8SEmil Constantinescu 
5760f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts)
5770f7a1166SHong Zhang {
5780f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
579cd4cee2dSHong Zhang   TS              quadts = ts->quadraturets;
5800f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
5810f7a1166SHong Zhang   const PetscInt  s = tab->s;
5820f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
5830f7a1166SHong Zhang   Vec             *Y = rk->Y;
5840f7a1166SHong Zhang   PetscInt        i;
5850f7a1166SHong Zhang   PetscErrorCode  ierr;
5860f7a1166SHong Zhang 
5870f7a1166SHong Zhang   PetscFunctionBegin;
588cd4cee2dSHong Zhang   /* No need to backup quadts->vec_sol since it can be reverted in TSRollBack_RK */
5890f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
590cd4cee2dSHong Zhang     /* Evolve quadrature TS solution to compute integrals */
591cd4cee2dSHong Zhang     ierr = TSComputeRHSFunction(quadts,rk->ptime+rk->time_step*c[i],Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
592cd4cee2dSHong Zhang     ierr = VecAXPY(quadts->vec_sol,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
5930f7a1166SHong Zhang   }
5940f7a1166SHong Zhang   PetscFunctionReturn(0);
5950f7a1166SHong Zhang }
5960f7a1166SHong Zhang 
5970f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts)
5980f7a1166SHong Zhang {
5990f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
6000f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
601cd4cee2dSHong Zhang   TS              quadts = ts->quadraturets;
6020f7a1166SHong Zhang   const PetscInt  s = tab->s;
6030f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
6040f7a1166SHong Zhang   Vec             *Y = rk->Y;
6050f7a1166SHong Zhang   PetscInt        i;
6060f7a1166SHong Zhang   PetscErrorCode  ierr;
6070f7a1166SHong Zhang 
6080f7a1166SHong Zhang   PetscFunctionBegin;
6090f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
610cd4cee2dSHong Zhang     /* Evolve quadrature TS solution to compute integrals */
611cd4cee2dSHong Zhang     ierr = TSComputeRHSFunction(quadts,ts->ptime+ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
612cd4cee2dSHong Zhang     ierr = VecAXPY(quadts->vec_sol,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
6130f7a1166SHong Zhang   }
6140f7a1166SHong Zhang   PetscFunctionReturn(0);
6150f7a1166SHong Zhang }
6160f7a1166SHong Zhang 
617fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts)
618fecfb714SLisandro Dalcin {
619fecfb714SLisandro Dalcin   TS_RK           *rk = (TS_RK*)ts->data;
620cd4cee2dSHong Zhang   TS              quadts = ts->quadraturets;
621fecfb714SLisandro Dalcin   RKTableau       tab = rk->tableau;
622fecfb714SLisandro Dalcin   const PetscInt  s  = tab->s;
623cd4cee2dSHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
624fecfb714SLisandro Dalcin   PetscScalar     *w = rk->work;
625cd4cee2dSHong Zhang   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
626fecfb714SLisandro Dalcin   PetscInt        j;
627fecfb714SLisandro Dalcin   PetscReal       h;
628fecfb714SLisandro Dalcin   PetscErrorCode  ierr;
629fecfb714SLisandro Dalcin 
630fecfb714SLisandro Dalcin   PetscFunctionBegin;
631fecfb714SLisandro Dalcin   switch (rk->status) {
632fecfb714SLisandro Dalcin   case TS_STEP_INCOMPLETE:
633fecfb714SLisandro Dalcin   case TS_STEP_PENDING:
634fecfb714SLisandro Dalcin     h = ts->time_step; break;
635fecfb714SLisandro Dalcin   case TS_STEP_COMPLETE:
636fecfb714SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
637fecfb714SLisandro Dalcin   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
638fecfb714SLisandro Dalcin   }
639fecfb714SLisandro Dalcin   for (j=0; j<s; j++) w[j] = -h*b[j];
640fecfb714SLisandro Dalcin   ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr);
641cd4cee2dSHong Zhang   if (quadts && ts->costintegralfwd) {
642cd4cee2dSHong Zhang     for (j=0; j<s; j++) {
643cd4cee2dSHong Zhang       /* Revert the quadrature TS solution */
644cd4cee2dSHong Zhang       ierr = TSComputeRHSFunction(quadts,rk->ptime+h*c[j],Y[j],ts->vec_costintegrand);CHKERRQ(ierr);
645cd4cee2dSHong Zhang       ierr = VecAXPY(quadts->vec_sol,-h*b[j],ts->vec_costintegrand);CHKERRQ(ierr);
646cd4cee2dSHong Zhang     }
647cd4cee2dSHong Zhang   }
648fecfb714SLisandro Dalcin   PetscFunctionReturn(0);
649fecfb714SLisandro Dalcin }
650fecfb714SLisandro Dalcin 
651922a638cSHong Zhang static PetscErrorCode TSForwardStep_RK(TS ts)
652922a638cSHong Zhang {
653922a638cSHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
654922a638cSHong Zhang   RKTableau       tab = rk->tableau;
655922a638cSHong Zhang   Mat             J,*MatsFwdSensipTemp = rk->MatsFwdSensipTemp;
656922a638cSHong Zhang   const PetscInt  s = tab->s;
657922a638cSHong Zhang   const PetscReal *A = tab->A,*c = tab->c,*b = tab->b;
658922a638cSHong Zhang   Vec             *Y = rk->Y;
659922a638cSHong Zhang   PetscInt        i,j;
660922a638cSHong Zhang   PetscReal       stage_time,h = ts->time_step;
661922a638cSHong Zhang   PetscBool       zero;
662922a638cSHong Zhang   PetscErrorCode  ierr;
663922a638cSHong Zhang 
664922a638cSHong Zhang   PetscFunctionBegin;
665922a638cSHong Zhang   ierr = MatCopy(ts->mat_sensip,rk->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
666922a638cSHong Zhang   ierr = TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr);
667922a638cSHong Zhang 
668922a638cSHong Zhang   for (i=0; i<s; i++) {
669922a638cSHong Zhang     stage_time = ts->ptime + h*c[i];
670922a638cSHong Zhang     zero = PETSC_FALSE;
671922a638cSHong Zhang     if (b[i] == 0 && i == s-1) zero = PETSC_TRUE;
672922a638cSHong Zhang     /* TLM Stage values */
673922a638cSHong Zhang     if (!i) {
674922a638cSHong Zhang       ierr = MatCopy(ts->mat_sensip,rk->MatsFwdStageSensip[i],SAME_NONZERO_PATTERN);CHKERRQ(ierr);
675922a638cSHong Zhang     } else if (!zero) {
676922a638cSHong Zhang       ierr = MatZeroEntries(rk->MatsFwdStageSensip[i]);CHKERRQ(ierr);
677922a638cSHong Zhang       for (j=0; j<i; j++) {
678922a638cSHong Zhang         ierr = MatAXPY(rk->MatsFwdStageSensip[i],h*A[i*s+j],MatsFwdSensipTemp[j],SAME_NONZERO_PATTERN);CHKERRQ(ierr);
679922a638cSHong Zhang       }
680922a638cSHong Zhang       ierr = MatAXPY(rk->MatsFwdStageSensip[i],1.,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
681922a638cSHong Zhang     } else {
682922a638cSHong Zhang       ierr = MatZeroEntries(rk->MatsFwdStageSensip[i]);CHKERRQ(ierr);
683922a638cSHong Zhang     }
684922a638cSHong Zhang 
685922a638cSHong Zhang     ierr = TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);CHKERRQ(ierr);
686922a638cSHong Zhang     ierr = MatMatMult(J,rk->MatsFwdStageSensip[i],MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatsFwdSensipTemp[i]);CHKERRQ(ierr);
68713af1a74SHong Zhang     if (ts->Jacprhs) {
68813af1a74SHong Zhang       ierr = TSComputeRHSJacobianP(ts,stage_time,Y[i],ts->Jacprhs);CHKERRQ(ierr); /* get f_p */
68913af1a74SHong Zhang       if (ts->vecs_sensi2p) { /* TLM used for 2nd-order adjoint */
69013af1a74SHong Zhang         PetscScalar *xarr;
69113af1a74SHong Zhang         ierr = MatDenseGetColumn(MatsFwdSensipTemp[i],0,&xarr);CHKERRQ(ierr);
69213af1a74SHong Zhang         ierr = VecPlaceArray(rk->VecDeltaFwdSensipCol,xarr);CHKERRQ(ierr);
69313af1a74SHong Zhang         ierr = MatMultAdd(ts->Jacprhs,ts->vec_dir,rk->VecDeltaFwdSensipCol,rk->VecDeltaFwdSensipCol);CHKERRQ(ierr);
694c3b366b1Sprj-         ierr = VecResetArray(rk->VecDeltaFwdSensipCol);CHKERRQ(ierr);
69513af1a74SHong Zhang         ierr = MatDenseRestoreColumn(MatsFwdSensipTemp[i],&xarr);CHKERRQ(ierr);
69613af1a74SHong Zhang       } else {
69713af1a74SHong Zhang         ierr = MatAXPY(MatsFwdSensipTemp[i],1.,ts->Jacprhs,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
69813af1a74SHong Zhang       }
699922a638cSHong Zhang     }
700922a638cSHong Zhang   }
701922a638cSHong Zhang 
702922a638cSHong Zhang   for (i=0; i<s; i++) {
703922a638cSHong Zhang     ierr = MatAXPY(ts->mat_sensip,h*b[i],rk->MatsFwdSensipTemp[i],SAME_NONZERO_PATTERN);CHKERRQ(ierr);
704922a638cSHong Zhang   }
705922a638cSHong Zhang   rk->status = TS_STEP_COMPLETE;
706922a638cSHong Zhang   PetscFunctionReturn(0);
707922a638cSHong Zhang }
708922a638cSHong Zhang 
709922a638cSHong Zhang static PetscErrorCode TSForwardGetStages_RK(TS ts,PetscInt *ns,Mat **stagesensip)
710922a638cSHong Zhang {
711922a638cSHong Zhang   TS_RK     *rk = (TS_RK*)ts->data;
712922a638cSHong Zhang   RKTableau tab  = rk->tableau;
713922a638cSHong Zhang 
714922a638cSHong Zhang   PetscFunctionBegin;
715922a638cSHong Zhang   if (ns) *ns = tab->s;
716922a638cSHong Zhang   if (stagesensip) *stagesensip = rk->MatsFwdStageSensip;
717922a638cSHong Zhang   PetscFunctionReturn(0);
718922a638cSHong Zhang }
719922a638cSHong Zhang 
720922a638cSHong Zhang static PetscErrorCode TSForwardSetUp_RK(TS ts)
721922a638cSHong Zhang {
722922a638cSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
723922a638cSHong Zhang   RKTableau      tab  = rk->tableau;
724922a638cSHong Zhang   PetscInt       i;
725922a638cSHong Zhang   PetscErrorCode ierr;
726922a638cSHong Zhang 
727922a638cSHong Zhang   PetscFunctionBegin;
728922a638cSHong Zhang   /* backup sensitivity results for roll-backs */
729922a638cSHong Zhang   ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatFwdSensip0);CHKERRQ(ierr);
730922a638cSHong Zhang 
731922a638cSHong Zhang   ierr = PetscMalloc1(tab->s,&rk->MatsFwdStageSensip);CHKERRQ(ierr);
732922a638cSHong Zhang   ierr = PetscMalloc1(tab->s,&rk->MatsFwdSensipTemp);CHKERRQ(ierr);
733922a638cSHong Zhang   for (i=0; i<tab->s; i++) {
734922a638cSHong Zhang     ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatsFwdStageSensip[i]);CHKERRQ(ierr);
735922a638cSHong Zhang     ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatsFwdSensipTemp[i]);CHKERRQ(ierr);
736922a638cSHong Zhang   }
737922a638cSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&rk->VecDeltaFwdSensipCol);CHKERRQ(ierr);
738922a638cSHong Zhang   PetscFunctionReturn(0);
739922a638cSHong Zhang }
740922a638cSHong Zhang 
741922a638cSHong Zhang static PetscErrorCode TSForwardReset_RK(TS ts)
742922a638cSHong Zhang {
743922a638cSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
744922a638cSHong Zhang   RKTableau      tab  = rk->tableau;
745922a638cSHong Zhang   PetscInt       i;
746922a638cSHong Zhang   PetscErrorCode ierr;
747922a638cSHong Zhang 
748922a638cSHong Zhang   PetscFunctionBegin;
749922a638cSHong Zhang   ierr = MatDestroy(&rk->MatFwdSensip0);CHKERRQ(ierr);
750922a638cSHong Zhang   if (rk->MatsFwdStageSensip) {
751922a638cSHong Zhang     for (i=0; i<tab->s; i++) {
752922a638cSHong Zhang       ierr = MatDestroy(&rk->MatsFwdStageSensip[i]);CHKERRQ(ierr);
753922a638cSHong Zhang     }
754922a638cSHong Zhang     ierr = PetscFree(rk->MatsFwdStageSensip);CHKERRQ(ierr);
755922a638cSHong Zhang   }
756922a638cSHong Zhang   if (rk->MatsFwdSensipTemp) {
757922a638cSHong Zhang     for (i=0; i<tab->s; i++) {
758922a638cSHong Zhang       ierr = MatDestroy(&rk->MatsFwdSensipTemp[i]);CHKERRQ(ierr);
759922a638cSHong Zhang     }
760922a638cSHong Zhang     ierr = PetscFree(rk->MatsFwdSensipTemp);CHKERRQ(ierr);
761922a638cSHong Zhang   }
762922a638cSHong Zhang   ierr = VecDestroy(&rk->VecDeltaFwdSensipCol);CHKERRQ(ierr);
763922a638cSHong Zhang   PetscFunctionReturn(0);
764922a638cSHong Zhang }
765922a638cSHong Zhang 
766f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts)
767f68a32c8SEmil Constantinescu {
768f68a32c8SEmil Constantinescu   TS_RK           *rk  = (TS_RK*)ts->data;
769f68a32c8SEmil Constantinescu   RKTableau       tab  = rk->tableau;
770f68a32c8SEmil Constantinescu   const PetscInt  s = tab->s;
771fecfb714SLisandro Dalcin   const PetscReal *A = tab->A,*c = tab->c;
772f68a32c8SEmil Constantinescu   PetscScalar     *w = rk->work;
773f68a32c8SEmil Constantinescu   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
774d1905564SLisandro Dalcin   PetscBool       FSAL = tab->FSAL;
775f68a32c8SEmil Constantinescu   TSAdapt         adapt;
776fecfb714SLisandro Dalcin   PetscInt        i,j;
777be5899b3SLisandro Dalcin   PetscInt        rejections = 0;
778be5899b3SLisandro Dalcin   PetscBool       stageok,accept = PETSC_TRUE;
779be5899b3SLisandro Dalcin   PetscReal       next_time_step = ts->time_step;
780f68a32c8SEmil Constantinescu   PetscErrorCode  ierr;
781f68a32c8SEmil Constantinescu 
782f68a32c8SEmil Constantinescu   PetscFunctionBegin;
783d1905564SLisandro Dalcin   if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE;
784d1905564SLisandro Dalcin   if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); }
785d1905564SLisandro Dalcin 
786f68a32c8SEmil Constantinescu   rk->status = TS_STEP_INCOMPLETE;
787be5899b3SLisandro Dalcin   while (!ts->reason && rk->status != TS_STEP_COMPLETE) {
788be5899b3SLisandro Dalcin     PetscReal t = ts->ptime;
789f68a32c8SEmil Constantinescu     PetscReal h = ts->time_step;
790f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
7919be3e283SDebojyoti Ghosh       rk->stage_time = t + h*c[i];
7929be3e283SDebojyoti Ghosh       ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr);
793f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
794f68a32c8SEmil Constantinescu       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
795f68a32c8SEmil Constantinescu       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
7969be3e283SDebojyoti Ghosh       ierr = TSPostStage(ts,rk->stage_time,i,Y);CHKERRQ(ierr);
797f68a32c8SEmil Constantinescu       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
798be5899b3SLisandro Dalcin       ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr);
799fecfb714SLisandro Dalcin       if (!stageok) goto reject_step;
8008f738a7cSLisandro Dalcin       if (FSAL && !i) continue;
801f68a32c8SEmil Constantinescu       ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
802f68a32c8SEmil Constantinescu     }
803be5899b3SLisandro Dalcin 
804be5899b3SLisandro Dalcin     rk->status = TS_STEP_INCOMPLETE;
805f68a32c8SEmil Constantinescu     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
806f68a32c8SEmil Constantinescu     rk->status = TS_STEP_PENDING;
807f68a32c8SEmil Constantinescu     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
808f68a32c8SEmil Constantinescu     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
8091917a363SLisandro Dalcin     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr);
810fecfb714SLisandro Dalcin     ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
811be5899b3SLisandro Dalcin     rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
812be5899b3SLisandro Dalcin     if (!accept) { /* Roll back the current step */
813fecfb714SLisandro Dalcin       ierr = TSRollBack_RK(ts);CHKERRQ(ierr);
814be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
815be5899b3SLisandro Dalcin       goto reject_step;
816be5899b3SLisandro Dalcin     }
817be5899b3SLisandro Dalcin 
8180f7a1166SHong Zhang     if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
8190f7a1166SHong Zhang       rk->ptime     = ts->ptime;
8200f7a1166SHong Zhang       rk->time_step = ts->time_step;
821493ed95fSHong Zhang     }
822be5899b3SLisandro Dalcin 
823f68a32c8SEmil Constantinescu     ts->ptime += ts->time_step;
824f68a32c8SEmil Constantinescu     ts->time_step = next_time_step;
825f68a32c8SEmil Constantinescu     break;
826be5899b3SLisandro Dalcin 
827be5899b3SLisandro Dalcin     reject_step:
828fecfb714SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
829be5899b3SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
830be5899b3SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
831be5899b3SLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
832f68a32c8SEmil Constantinescu     }
833f68a32c8SEmil Constantinescu   }
834f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
835e4dd521cSBarry Smith }
836e4dd521cSBarry Smith 
837f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts)
838f6a906c0SBarry Smith {
839f6a906c0SBarry Smith   TS_RK          *rk  = (TS_RK*)ts->data;
840be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
841be5899b3SLisandro Dalcin   PetscInt       s   = tab->s;
842f6a906c0SBarry Smith   PetscErrorCode ierr;
843f6a906c0SBarry Smith 
844f6a906c0SBarry Smith   PetscFunctionBegin;
845f6a906c0SBarry Smith   if (ts->adjointsetupcalled++) PetscFunctionReturn(0);
8462e7b7f96SHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecsDeltaLam);CHKERRQ(ierr);
8472e7b7f96SHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&rk->VecsSensiTemp);CHKERRQ(ierr);
848f6a906c0SBarry Smith   if (ts->vecs_sensip) {
8492e7b7f96SHong Zhang     ierr = VecDuplicate(ts->vecs_sensip[0],&rk->VecDeltaMu);CHKERRQ(ierr);
850f6a906c0SBarry Smith   }
85113af1a74SHong Zhang   if (ts->vecs_sensi2) {
85213af1a74SHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecsDeltaLam2);CHKERRQ(ierr);
85313af1a74SHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&rk->VecsSensi2Temp);CHKERRQ(ierr);
85413af1a74SHong Zhang   }
85513af1a74SHong Zhang   if (ts->vecs_sensi2p) {
85613af1a74SHong Zhang     ierr = VecDuplicate(ts->vecs_sensi2p[0],&rk->VecDeltaMu2);CHKERRQ(ierr);
85713af1a74SHong Zhang   }
858f6a906c0SBarry Smith   PetscFunctionReturn(0);
859f6a906c0SBarry Smith }
860f6a906c0SBarry Smith 
86135f5172aSHong Zhang /*
86235f5172aSHong Zhang   Assumptions:
86335f5172aSHong Zhang     - TSStep_RK() always evaluates the step with b, not bembed.
86435f5172aSHong Zhang */
86542f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts)
866d2daff3dSHong Zhang {
867c235aa8dSHong Zhang   TS_RK            *rk = (TS_RK*)ts->data;
868cd4cee2dSHong Zhang   TS               quadts = ts->quadraturets;
869c235aa8dSHong Zhang   RKTableau        tab = rk->tableau;
8703ca0882eSHong Zhang   Mat              J,Jpre,Jquad;
871c235aa8dSHong Zhang   const PetscInt   s = tab->s;
872c235aa8dSHong Zhang   const PetscReal  *A = tab->A,*b = tab->b,*c = tab->c;
87313af1a74SHong Zhang   PetscScalar      *w = rk->work,*xarr;
8742e7b7f96SHong Zhang   Vec              *Y = rk->Y,*VecsDeltaLam = rk->VecsDeltaLam,VecDeltaMu = rk->VecDeltaMu,*VecsSensiTemp = rk->VecsSensiTemp;
87513af1a74SHong Zhang   Vec              *VecsDeltaLam2 = rk->VecsDeltaLam2,VecDeltaMu2 = rk->VecDeltaMu2,*VecsSensi2Temp = rk->VecsSensi2Temp;
876cd4cee2dSHong Zhang   Vec              VecDRDUTransCol = ts->vec_drdu_col,VecDRDPTransCol = ts->vec_drdp_col;
877b18ea86cSHong Zhang   PetscInt         i,j,nadj;
8783d3eaba7SBarry Smith   PetscReal        t = ts->ptime;
8793ca0882eSHong Zhang   PetscReal        h = ts->time_step;
880d2daff3dSHong Zhang   PetscErrorCode   ierr;
881c235aa8dSHong Zhang 
882d2daff3dSHong Zhang   PetscFunctionBegin;
883c235aa8dSHong Zhang   rk->status = TS_STEP_INCOMPLETE;
8843389c7d9SStefano Zampini 
8853ca0882eSHong Zhang   ierr = TSGetRHSJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr);
88635f5172aSHong Zhang   if (quadts) {
887cd4cee2dSHong Zhang     ierr = TSGetRHSJacobian(quadts,&Jquad,NULL,NULL,NULL);CHKERRQ(ierr);
88835f5172aSHong Zhang   }
8892e7b7f96SHong Zhang   for (i=s-1; i>=0; i--) {
8906a1e4597SHong Zhang     if (tab->FSAL && i == s-1) {
8916a1e4597SHong Zhang       /* VecsDeltaLam[nadj*s+s-1] are initialized with zeros and the values never change.*/
8926a1e4597SHong Zhang       continue;
8936a1e4597SHong Zhang     }
8943ca0882eSHong Zhang     rk->stage_time = t + h*(1.0-c[i]);
895fd850964SHong Zhang     ierr = TSComputeSNESJacobian(ts,Y[i],J,Jpre);CHKERRQ(ierr);
89635f5172aSHong Zhang     if (quadts) {
8973ca0882eSHong Zhang       ierr = TSComputeRHSJacobian(quadts,rk->stage_time,Y[i],Jquad,Jquad);CHKERRQ(ierr); /* get r_u^T */
89835f5172aSHong Zhang     }
8993389c7d9SStefano Zampini     if (ts->vecs_sensip) {
9003ca0882eSHong Zhang       ierr = TSComputeRHSJacobianP(ts,rk->stage_time,Y[i],ts->Jacprhs);CHKERRQ(ierr); /* get f_p */
90135f5172aSHong Zhang       if (quadts) {
9023ca0882eSHong Zhang         ierr = TSComputeRHSJacobianP(quadts,rk->stage_time,Y[i],quadts->Jacprhs);CHKERRQ(ierr); /* get f_p for the quadrature */
9033389c7d9SStefano Zampini       }
90435f5172aSHong Zhang     }
9053389c7d9SStefano Zampini 
90635f5172aSHong Zhang     if (b[i]) {
90735f5172aSHong Zhang       for (j=i+1; j<s; j++) w[j-i-1] = A[j*s+i]/b[i]; /* coefficients for computing VecsSensiTemp */
90835f5172aSHong Zhang     } else {
90935f5172aSHong Zhang       for (j=i+1; j<s; j++) w[j-i-1] = A[j*s+i]; /* coefficients for computing VecsSensiTemp */
91035f5172aSHong Zhang     }
9112e7b7f96SHong Zhang 
9122e7b7f96SHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
9133389c7d9SStefano Zampini       /* Stage values of lambda */
91435f5172aSHong Zhang       if (b[i]) {
91535f5172aSHong Zhang         /* lambda_{n+1} + \sum_{j=i+1}^s a_{ji}/b[i]*lambda_{s,j} */
9162e7b7f96SHong Zhang         ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); /* VecDeltaLam is an vec array of size s by numcost */
9172e7b7f96SHong Zhang         ierr = VecMAXPY(VecsSensiTemp[nadj],s-i-1,w,&VecsDeltaLam[nadj*s+i+1]);CHKERRQ(ierr);
91835f5172aSHong Zhang         ierr = MatMultTranspose(J,VecsSensiTemp[nadj],VecsDeltaLam[nadj*s+i]);CHKERRQ(ierr); /* VecsSensiTemp will be reused by 2nd-order adjoint */
91935f5172aSHong Zhang         ierr = VecScale(VecsDeltaLam[nadj*s+i],-h*b[i]);CHKERRQ(ierr);
92035f5172aSHong Zhang         if (quadts) {
921cd4cee2dSHong Zhang           ierr = MatDenseGetColumn(Jquad,nadj,&xarr);CHKERRQ(ierr);
922cd4cee2dSHong Zhang           ierr = VecPlaceArray(VecDRDUTransCol,xarr);CHKERRQ(ierr);
923cd4cee2dSHong Zhang           ierr = VecAXPY(VecsDeltaLam[nadj*s+i],-h*b[i],VecDRDUTransCol);CHKERRQ(ierr);
924cd4cee2dSHong Zhang           ierr = VecResetArray(VecDRDUTransCol);CHKERRQ(ierr);
925cd4cee2dSHong Zhang           ierr = MatDenseRestoreColumn(Jquad,&xarr);CHKERRQ(ierr);
926cd4cee2dSHong Zhang         }
9273389c7d9SStefano Zampini       } else {
9286a1e4597SHong Zhang         /* \sum_{j=i+1}^s a_{ji}*lambda_{s,j} */
9296a1e4597SHong Zhang         ierr = VecSet(VecsSensiTemp[nadj],0);CHKERRQ(ierr);
9306a1e4597SHong Zhang         ierr = VecMAXPY(VecsSensiTemp[nadj],s-i-1,w,&VecsDeltaLam[nadj*s+i+1]);CHKERRQ(ierr);
9316a1e4597SHong Zhang         ierr = MatMultTranspose(J,VecsSensiTemp[nadj],VecsDeltaLam[nadj*s+i]);CHKERRQ(ierr);
9326a1e4597SHong Zhang         ierr = VecScale(VecsDeltaLam[nadj*s+i],-h);CHKERRQ(ierr);
9333389c7d9SStefano Zampini       }
9346497ce14SHong Zhang 
935ad8e2604SHong Zhang       /* Stage values of mu */
9366497ce14SHong Zhang       if (ts->vecs_sensip) {
93713af1a74SHong Zhang         ierr = MatMultTranspose(ts->Jacprhs,VecsSensiTemp[nadj],VecDeltaMu);CHKERRQ(ierr);
93835f5172aSHong Zhang         if (b[i]) {
93935f5172aSHong Zhang           ierr = VecScale(VecDeltaMu,-h*b[i]);CHKERRQ(ierr);
94035f5172aSHong Zhang           if (quadts) {
941cd4cee2dSHong Zhang             ierr = MatDenseGetColumn(quadts->Jacprhs,nadj,&xarr);CHKERRQ(ierr);
942cd4cee2dSHong Zhang             ierr = VecPlaceArray(VecDRDPTransCol,xarr);CHKERRQ(ierr);
943cd4cee2dSHong Zhang             ierr = VecAXPY(VecDeltaMu,-h*b[i],VecDRDPTransCol);CHKERRQ(ierr);
944cd4cee2dSHong Zhang             ierr = VecResetArray(VecDRDPTransCol);CHKERRQ(ierr);
945cd4cee2dSHong Zhang             ierr = MatDenseRestoreColumn(quadts->Jacprhs,&xarr);CHKERRQ(ierr);
946493ed95fSHong Zhang           }
94735f5172aSHong Zhang         } else {
94835f5172aSHong Zhang           ierr = VecScale(VecDeltaMu,-h);CHKERRQ(ierr);
94935f5172aSHong Zhang         }
9502e7b7f96SHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],1.,VecDeltaMu);CHKERRQ(ierr); /* update sensip for each stage */
951ad8e2604SHong Zhang       }
952c235aa8dSHong Zhang     }
95313af1a74SHong Zhang 
95413af1a74SHong Zhang     if (ts->vecs_sensi2 && ts->forward_solve) { /* 2nd-order adjoint, TLM mode has to be turned on */
95513af1a74SHong Zhang       /* Get w1 at t_{n+1} from TLM matrix */
95613af1a74SHong Zhang       ierr = MatDenseGetColumn(rk->MatsFwdStageSensip[i],0,&xarr);CHKERRQ(ierr);
95713af1a74SHong Zhang       ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
95813af1a74SHong Zhang       /* lambda_s^T F_UU w_1 */
9593ca0882eSHong Zhang       ierr = TSComputeRHSHessianProductFunctionUU(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_sensip_col,ts->vecs_guu);CHKERRQ(ierr);
96035f5172aSHong Zhang       if (quadts)  {
96135f5172aSHong Zhang         /* R_UU w_1 */
9623ca0882eSHong Zhang         ierr = TSComputeRHSHessianProductFunctionUU(quadts,rk->stage_time,Y[i],NULL,ts->vec_sensip_col,ts->vecs_guu);CHKERRQ(ierr);
96335f5172aSHong Zhang       }
96435f5172aSHong Zhang       if (ts->vecs_sensip) {
96513af1a74SHong Zhang         /* lambda_s^T F_UP w_2 */
9663ca0882eSHong Zhang         ierr = TSComputeRHSHessianProductFunctionUP(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_dir,ts->vecs_gup);CHKERRQ(ierr);
96735f5172aSHong Zhang         if (quadts)  {
96835f5172aSHong Zhang           /* R_UP w_2 */
9693ca0882eSHong Zhang           ierr = TSComputeRHSHessianProductFunctionUP(quadts,rk->stage_time,Y[i],NULL,ts->vec_sensip_col,ts->vecs_gup);CHKERRQ(ierr);
97035f5172aSHong Zhang         }
97135f5172aSHong Zhang       }
97213af1a74SHong Zhang       if (ts->vecs_sensi2p) {
97313af1a74SHong Zhang         /* lambda_s^T F_PU w_1 */
9743ca0882eSHong Zhang         ierr = TSComputeRHSHessianProductFunctionPU(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_sensip_col,ts->vecs_gpu);CHKERRQ(ierr);
97535f5172aSHong Zhang         /* lambda_s^T F_PP w_2 */
9763ca0882eSHong Zhang         ierr = TSComputeRHSHessianProductFunctionPP(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_dir,ts->vecs_gpp);CHKERRQ(ierr);
97735f5172aSHong Zhang         if (b[i] && quadts) {
97835f5172aSHong Zhang           /* R_PU w_1 */
9793ca0882eSHong Zhang           ierr = TSComputeRHSHessianProductFunctionPU(quadts,rk->stage_time,Y[i],NULL,ts->vec_sensip_col,ts->vecs_gpu);CHKERRQ(ierr);
98035f5172aSHong Zhang           /* R_PP w_2 */
9813ca0882eSHong Zhang           ierr = TSComputeRHSHessianProductFunctionPP(quadts,rk->stage_time,Y[i],NULL,ts->vec_dir,ts->vecs_gpp);CHKERRQ(ierr);
98235f5172aSHong Zhang         }
98313af1a74SHong Zhang       }
98413af1a74SHong Zhang       ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
98513af1a74SHong Zhang       ierr = MatDenseRestoreColumn(rk->MatsFwdStageSensip[i],&xarr);CHKERRQ(ierr);
98613af1a74SHong Zhang 
98713af1a74SHong Zhang       for (nadj=0; nadj<ts->numcost; nadj++) {
98813af1a74SHong Zhang         /* Stage values of lambda */
98935f5172aSHong Zhang         if (b[i]) {
99035f5172aSHong Zhang           /* J_i^T*(Lambda_{n+1}+\sum_{j=i+1}^s a_{ji}/b_i*Lambda_{s,j} */
99113af1a74SHong Zhang           ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
99213af1a74SHong Zhang           ierr = VecMAXPY(VecsSensi2Temp[nadj],s-i-1,w,&VecsDeltaLam2[nadj*s+i+1]);CHKERRQ(ierr);
99313af1a74SHong Zhang           ierr = MatMultTranspose(J,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj*s+i]);CHKERRQ(ierr);
99435f5172aSHong Zhang           ierr = VecScale(VecsDeltaLam2[nadj*s+i],-h*b[i]);CHKERRQ(ierr);
99535f5172aSHong Zhang           ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],-h*b[i],ts->vecs_guu[nadj]);CHKERRQ(ierr);
99635f5172aSHong Zhang           if (ts->vecs_sensip) {
99735f5172aSHong Zhang             ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],-h*b[i],ts->vecs_gup[nadj]);CHKERRQ(ierr);
99813af1a74SHong Zhang           }
99913af1a74SHong Zhang         } else {
100035f5172aSHong Zhang           /* \sum_{j=i+1}^s a_{ji}*Lambda_{s,j} */
100113af1a74SHong Zhang           ierr = VecSet(VecsDeltaLam2[nadj*s+i],0);CHKERRQ(ierr);
100235f5172aSHong Zhang           ierr = VecMAXPY(VecsSensi2Temp[nadj],s-i-1,w,&VecsDeltaLam2[nadj*s+i+1]);CHKERRQ(ierr);
100335f5172aSHong Zhang           ierr = MatMultTranspose(J,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj*s+i]);CHKERRQ(ierr);
100435f5172aSHong Zhang           ierr = VecScale(VecsDeltaLam2[nadj*s+i],-h);CHKERRQ(ierr);
100535f5172aSHong Zhang           ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],-h,ts->vecs_guu[nadj]);CHKERRQ(ierr);
100635f5172aSHong Zhang           if (ts->vecs_sensip) {
100735f5172aSHong Zhang             ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],-h,ts->vecs_gup[nadj]);CHKERRQ(ierr);
100813af1a74SHong Zhang           }
100935f5172aSHong Zhang         }
101035f5172aSHong Zhang         if (ts->vecs_sensi2p) { /* 2nd-order adjoint for parameters */
101113af1a74SHong Zhang           ierr = MatMultTranspose(ts->Jacprhs,VecsSensi2Temp[nadj],VecDeltaMu2);CHKERRQ(ierr);
101235f5172aSHong Zhang           if (b[i]) {
101335f5172aSHong Zhang             ierr = VecScale(VecDeltaMu2,-h*b[i]);CHKERRQ(ierr);
101435f5172aSHong Zhang             ierr = VecAXPY(VecDeltaMu2,-h*b[i],ts->vecs_gpu[nadj]);CHKERRQ(ierr);
101535f5172aSHong Zhang             ierr = VecAXPY(VecDeltaMu2,-h*b[i],ts->vecs_gpp[nadj]);CHKERRQ(ierr);
101613af1a74SHong Zhang           } else {
101735f5172aSHong Zhang             ierr = VecScale(VecDeltaMu2,-h);CHKERRQ(ierr);
101835f5172aSHong Zhang             ierr = VecAXPY(VecDeltaMu2,-h,ts->vecs_gpu[nadj]);CHKERRQ(ierr);
101935f5172aSHong Zhang             ierr = VecAXPY(VecDeltaMu2,-h,ts->vecs_gpp[nadj]);CHKERRQ(ierr);
102013af1a74SHong Zhang           }
102113af1a74SHong Zhang           ierr = VecAXPY(ts->vecs_sensi2p[nadj],1,VecDeltaMu2);CHKERRQ(ierr); /* update sensi2p for each stage */
102213af1a74SHong Zhang         }
102313af1a74SHong Zhang       }
102413af1a74SHong Zhang     }
10256497ce14SHong Zhang   }
1026c235aa8dSHong Zhang 
1027c235aa8dSHong Zhang   for (j=0; j<s; j++) w[j] = 1.0;
10282e7b7f96SHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) { /* no need to do this for mu's */
10292e7b7f96SHong Zhang     ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecsDeltaLam[nadj*s]);CHKERRQ(ierr);
103013af1a74SHong Zhang     if (ts->vecs_sensi2) {
103113af1a74SHong Zhang       ierr = VecMAXPY(ts->vecs_sensi2[nadj],s,w,&VecsDeltaLam2[nadj*s]);CHKERRQ(ierr);
103213af1a74SHong Zhang     }
10336497ce14SHong Zhang   }
1034c235aa8dSHong Zhang   rk->status = TS_STEP_COMPLETE;
1035d2daff3dSHong Zhang   PetscFunctionReturn(0);
1036d2daff3dSHong Zhang }
1037d2daff3dSHong Zhang 
103813af1a74SHong Zhang static PetscErrorCode TSAdjointReset_RK(TS ts)
103913af1a74SHong Zhang {
104013af1a74SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
104113af1a74SHong Zhang   RKTableau      tab = rk->tableau;
104213af1a74SHong Zhang   PetscErrorCode ierr;
104313af1a74SHong Zhang 
104413af1a74SHong Zhang   PetscFunctionBegin;
104513af1a74SHong Zhang   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam);CHKERRQ(ierr);
104613af1a74SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&rk->VecsSensiTemp);CHKERRQ(ierr);
104713af1a74SHong Zhang   ierr = VecDestroy(&rk->VecDeltaMu);CHKERRQ(ierr);
104813af1a74SHong Zhang   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam2);CHKERRQ(ierr);
104913af1a74SHong Zhang   ierr = VecDestroy(&rk->VecDeltaMu2);CHKERRQ(ierr);
105013af1a74SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&rk->VecsSensi2Temp);CHKERRQ(ierr);
105113af1a74SHong Zhang   PetscFunctionReturn(0);
105213af1a74SHong Zhang }
105313af1a74SHong Zhang 
105455de54fcSHong Zhang static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X)
105555de54fcSHong Zhang {
105655de54fcSHong Zhang   TS_RK            *rk = (TS_RK*)ts->data;
105755de54fcSHong Zhang   PetscInt         s  = rk->tableau->s,p = rk->tableau->p,i,j;
105855de54fcSHong Zhang   PetscReal        h;
105955de54fcSHong Zhang   PetscReal        tt,t;
106055de54fcSHong Zhang   PetscScalar      *b;
106155de54fcSHong Zhang   const PetscReal  *B = rk->tableau->binterp;
106255de54fcSHong Zhang   PetscErrorCode   ierr;
106355de54fcSHong Zhang 
106455de54fcSHong Zhang   PetscFunctionBegin;
106555de54fcSHong Zhang   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
106655de54fcSHong Zhang 
106755de54fcSHong Zhang   switch (rk->status) {
106855de54fcSHong Zhang     case TS_STEP_INCOMPLETE:
106955de54fcSHong Zhang     case TS_STEP_PENDING:
107055de54fcSHong Zhang       h = ts->time_step;
107155de54fcSHong Zhang       t = (itime - ts->ptime)/h;
107255de54fcSHong Zhang       break;
107355de54fcSHong Zhang     case TS_STEP_COMPLETE:
107455de54fcSHong Zhang       h = ts->ptime - ts->ptime_prev;
107555de54fcSHong Zhang       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
107655de54fcSHong Zhang       break;
107755de54fcSHong Zhang     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
107855de54fcSHong Zhang   }
107955de54fcSHong Zhang   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
108055de54fcSHong Zhang   for (i=0; i<s; i++) b[i] = 0;
108155de54fcSHong Zhang   for (j=0,tt=t; j<p; j++,tt*=t) {
108255de54fcSHong Zhang     for (i=0; i<s; i++) {
108355de54fcSHong Zhang       b[i]  += h * B[i*p+j] * tt;
108455de54fcSHong Zhang     }
108555de54fcSHong Zhang   }
108655de54fcSHong Zhang   ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr);
108755de54fcSHong Zhang   ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr);
108855de54fcSHong Zhang   ierr = PetscFree(b);CHKERRQ(ierr);
108955de54fcSHong Zhang   PetscFunctionReturn(0);
109055de54fcSHong Zhang }
109155de54fcSHong Zhang 
109255de54fcSHong Zhang /*------------------------------------------------------------*/
109355de54fcSHong Zhang 
1094be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts)
1095be5899b3SLisandro Dalcin {
1096be5899b3SLisandro Dalcin   TS_RK          *rk = (TS_RK*)ts->data;
1097be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
1098be5899b3SLisandro Dalcin   PetscErrorCode ierr;
1099be5899b3SLisandro Dalcin 
1100be5899b3SLisandro Dalcin   PetscFunctionBegin;
1101be5899b3SLisandro Dalcin   if (!tab) PetscFunctionReturn(0);
1102be5899b3SLisandro Dalcin   ierr = PetscFree(rk->work);CHKERRQ(ierr);
1103be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr);
1104be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr);
11052e7b7f96SHong Zhang   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam);CHKERRQ(ierr);
11062e7b7f96SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&rk->VecsSensiTemp);CHKERRQ(ierr);
11072e7b7f96SHong Zhang   ierr = VecDestroy(&rk->VecDeltaMu);CHKERRQ(ierr);
1108be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
1109be5899b3SLisandro Dalcin }
1110be5899b3SLisandro Dalcin 
1111277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts)
1112e4dd521cSBarry Smith {
11136849ba73SBarry Smith   PetscErrorCode ierr;
1114e4dd521cSBarry Smith 
1115e4dd521cSBarry Smith   PetscFunctionBegin;
1116be5899b3SLisandro Dalcin   ierr = TSRKTableauReset(ts);CHKERRQ(ierr);
11170fe4d17eSHong Zhang   if (ts->use_splitrhsfunction) {
111802555c94SHong Zhang     ierr = PetscTryMethod(ts,"TSReset_RK_MultirateSplit_C",(TS),(ts));CHKERRQ(ierr);
11190fe4d17eSHong Zhang   } else {
112002555c94SHong Zhang     ierr = PetscTryMethod(ts,"TSReset_RK_MultirateNonsplit_C",(TS),(ts));CHKERRQ(ierr);
11210fe4d17eSHong Zhang   }
1122277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
1123e4dd521cSBarry Smith }
1124277b19d0SLisandro Dalcin 
1125f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx)
1126f68a32c8SEmil Constantinescu {
1127f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1128f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1129f68a32c8SEmil Constantinescu }
1130f68a32c8SEmil Constantinescu 
1131f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
1132f68a32c8SEmil Constantinescu {
1133f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1134f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1135f68a32c8SEmil Constantinescu }
1136f68a32c8SEmil Constantinescu 
1137f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx)
1138f68a32c8SEmil Constantinescu {
1139f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1140f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1141f68a32c8SEmil Constantinescu }
1142f68a32c8SEmil Constantinescu 
1143f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
1144f68a32c8SEmil Constantinescu {
1145f68a32c8SEmil Constantinescu 
1146f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1147f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1148f68a32c8SEmil Constantinescu }
1149be5899b3SLisandro Dalcin 
1150be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts)
1151be5899b3SLisandro Dalcin {
1152be5899b3SLisandro Dalcin   TS_RK          *rk  = (TS_RK*)ts->data;
1153be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
1154be5899b3SLisandro Dalcin   PetscErrorCode ierr;
1155be5899b3SLisandro Dalcin 
1156be5899b3SLisandro Dalcin   PetscFunctionBegin;
1157be5899b3SLisandro Dalcin   ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr);
1158be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr);
1159be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr);
1160be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
1161be5899b3SLisandro Dalcin }
1162be5899b3SLisandro Dalcin 
1163f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts)
1164f68a32c8SEmil Constantinescu {
1165cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
1166f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
1167f68a32c8SEmil Constantinescu   DM             dm;
1168f68a32c8SEmil Constantinescu 
1169f68a32c8SEmil Constantinescu   PetscFunctionBegin;
11703dd83b38SBarry Smith   ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr);
1171be5899b3SLisandro Dalcin   ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);
1172cd4cee2dSHong Zhang   if (quadts && ts->costintegralfwd) {
1173cd4cee2dSHong Zhang     Mat Jquad;
1174cd4cee2dSHong Zhang     ierr = TSGetRHSJacobian(quadts,&Jquad,NULL,NULL,NULL);CHKERRQ(ierr);
11750f7a1166SHong Zhang   }
1176f68a32c8SEmil Constantinescu   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1177f68a32c8SEmil Constantinescu   ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
1178f68a32c8SEmil Constantinescu   ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
11790fe4d17eSHong Zhang   if (ts->use_splitrhsfunction) {
118002555c94SHong Zhang     ierr = PetscTryMethod(ts,"TSSetUp_RK_MultirateSplit_C",(TS),(ts));CHKERRQ(ierr);
11810fe4d17eSHong Zhang   } else {
118202555c94SHong Zhang     ierr = PetscTryMethod(ts,"TSSetUp_RK_MultirateNonsplit_C",(TS),(ts));CHKERRQ(ierr);
11830fe4d17eSHong Zhang   }
1184e4dd521cSBarry Smith   PetscFunctionReturn(0);
1185e4dd521cSBarry Smith }
1186d2daff3dSHong Zhang 
11874416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts)
1188e4dd521cSBarry Smith {
1189be5899b3SLisandro Dalcin   TS_RK          *rk = (TS_RK*)ts->data;
1190dfbe8321SBarry Smith   PetscErrorCode ierr;
1191262ff7bbSBarry Smith 
1192e4dd521cSBarry Smith   PetscFunctionBegin;
1193e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr);
1194f68a32c8SEmil Constantinescu   {
1195f68a32c8SEmil Constantinescu     RKTableauLink link;
1196f68a32c8SEmil Constantinescu     PetscInt      count,choice;
11970fe4d17eSHong Zhang     PetscBool     flg,use_multirate = PETSC_FALSE;
1198f68a32c8SEmil Constantinescu     const char    **namelist;
11992c9cb062Sluzhanghpp 
1200f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) ;
1201f489ac74SBarry Smith     ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr);
1202f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
12030fe4d17eSHong Zhang     ierr = PetscOptionsBool("-ts_rk_multirate","Use interpolation-based multirate RK method","TSRKSetMultirate",rk->use_multirate,&use_multirate,&flg);CHKERRQ(ierr);
12040fe4d17eSHong Zhang     if (flg) {
12050fe4d17eSHong Zhang       ierr = TSRKSetMultirate(ts,use_multirate);CHKERRQ(ierr);
12060fe4d17eSHong Zhang     }
1207be5899b3SLisandro Dalcin     ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr);
1208be5899b3SLisandro Dalcin     if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);}
1209f68a32c8SEmil Constantinescu     ierr = PetscFree(namelist);CHKERRQ(ierr);
1210f68a32c8SEmil Constantinescu   }
1211262ff7bbSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
1212ea13f565SStefano Zampini   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)ts),NULL,"Multirate methods options","");CHKERRQ(ierr);
12132c9cb062Sluzhanghpp   ierr = PetscOptionsInt("-ts_rk_dtratio","time step ratio between slow and fast","",rk->dtratio,&rk->dtratio,NULL);CHKERRQ(ierr);
12142c9cb062Sluzhanghpp   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1215e4dd521cSBarry Smith   PetscFunctionReturn(0);
1216e4dd521cSBarry Smith }
1217e4dd521cSBarry Smith 
12185f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
1219e4dd521cSBarry Smith {
12205f70b478SJed Brown   TS_RK          *rk = (TS_RK*)ts->data;
12218370ee3bSLisandro Dalcin   PetscBool      iascii;
1222dfbe8321SBarry Smith   PetscErrorCode ierr;
12232cdf8120Sasbjorn 
1224e4dd521cSBarry Smith   PetscFunctionBegin;
1225251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
12268370ee3bSLisandro Dalcin   if (iascii) {
12279c334d8fSLisandro Dalcin     RKTableau       tab  = rk->tableau;
1228f68a32c8SEmil Constantinescu     TSRKType        rktype;
12290f7a28cdSStefano Zampini     const PetscReal *c;
12300f7a28cdSStefano Zampini     PetscInt        s;
1231f68a32c8SEmil Constantinescu     char            buf[512];
12320f7a28cdSStefano Zampini     PetscBool       FSAL;
12330f7a28cdSStefano Zampini 
1234f68a32c8SEmil Constantinescu     ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr);
12350f7a28cdSStefano Zampini     ierr = TSRKGetTableau(ts,&s,NULL,NULL,&c,NULL,NULL,NULL,&FSAL);CHKERRQ(ierr);
1236efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  RK type %s\n",rktype);CHKERRQ(ierr);
12370c37eb8eSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  Order: %D\n",tab->order);CHKERRQ(ierr);
12380f7a28cdSStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"  FSAL property: %s\n",FSAL ? "yes" : "no");CHKERRQ(ierr);
12390f7a28cdSStefano Zampini     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",s,c);CHKERRQ(ierr);
1240f68a32c8SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa c = %s\n",buf);CHKERRQ(ierr);
12418370ee3bSLisandro Dalcin   }
1242f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1243f68a32c8SEmil Constantinescu }
1244f68a32c8SEmil Constantinescu 
1245f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer)
1246f68a32c8SEmil Constantinescu {
1247f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
12489c334d8fSLisandro Dalcin   TSAdapt        adapt;
1249f68a32c8SEmil Constantinescu 
1250f68a32c8SEmil Constantinescu   PetscFunctionBegin;
12519c334d8fSLisandro Dalcin   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
12529c334d8fSLisandro Dalcin   ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr);
1253f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1254f68a32c8SEmil Constantinescu }
1255f68a32c8SEmil Constantinescu 
12562ea87ba9SLisandro Dalcin /*@
12572ea87ba9SLisandro Dalcin   TSRKGetOrder - Get the order of RK scheme
12582ea87ba9SLisandro Dalcin 
12592ea87ba9SLisandro Dalcin   Not collective
12602ea87ba9SLisandro Dalcin 
12612ea87ba9SLisandro Dalcin   Input Parameter:
12622ea87ba9SLisandro Dalcin .  ts - timestepping context
12632ea87ba9SLisandro Dalcin 
12642ea87ba9SLisandro Dalcin   Output Parameter:
12652ea87ba9SLisandro Dalcin .  order - order of RK-scheme
12662ea87ba9SLisandro Dalcin 
12672ea87ba9SLisandro Dalcin   Level: intermediate
12682ea87ba9SLisandro Dalcin 
12692ea87ba9SLisandro Dalcin .seealso: TSRKGetType()
12702ea87ba9SLisandro Dalcin @*/
12712ea87ba9SLisandro Dalcin PetscErrorCode TSRKGetOrder(TS ts,PetscInt *order)
12722ea87ba9SLisandro Dalcin {
12732ea87ba9SLisandro Dalcin   PetscErrorCode ierr;
12742ea87ba9SLisandro Dalcin 
12752ea87ba9SLisandro Dalcin   PetscFunctionBegin;
12762ea87ba9SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
12772ea87ba9SLisandro Dalcin   PetscValidIntPointer(order,2);
12782ea87ba9SLisandro Dalcin   ierr = PetscUseMethod(ts,"TSRKGetOrder_C",(TS,PetscInt*),(ts,order));CHKERRQ(ierr);
12792ea87ba9SLisandro Dalcin   PetscFunctionReturn(0);
12802ea87ba9SLisandro Dalcin }
12812ea87ba9SLisandro Dalcin 
1282f68a32c8SEmil Constantinescu /*@C
1283f68a32c8SEmil Constantinescu   TSRKSetType - Set the type of RK scheme
1284f68a32c8SEmil Constantinescu 
1285f68a32c8SEmil Constantinescu   Logically collective
1286f68a32c8SEmil Constantinescu 
1287*d8d19677SJose E. Roman   Input Parameters:
1288f68a32c8SEmil Constantinescu +  ts - timestepping context
1289f68a32c8SEmil Constantinescu -  rktype - type of RK-scheme
1290f68a32c8SEmil Constantinescu 
1291287c2655SBarry Smith   Options Database:
12929bd3e852SBarry Smith .   -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs>
1293287c2655SBarry Smith 
1294f68a32c8SEmil Constantinescu   Level: intermediate
1295f68a32c8SEmil Constantinescu 
129637acaa02SHendrik Ranocha .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS, TSRK6VR, TSRK7VR, TSRK8VR
1297f68a32c8SEmil Constantinescu @*/
1298f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype)
1299f68a32c8SEmil Constantinescu {
1300f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
1301f68a32c8SEmil Constantinescu 
1302f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1303f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1304b92453a8SLisandro Dalcin   PetscValidCharPointer(rktype,2);
1305f68a32c8SEmil Constantinescu   ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr);
1306f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1307f68a32c8SEmil Constantinescu }
1308f68a32c8SEmil Constantinescu 
1309f68a32c8SEmil Constantinescu /*@C
1310f68a32c8SEmil Constantinescu   TSRKGetType - Get the type of RK scheme
1311f68a32c8SEmil Constantinescu 
13122ea87ba9SLisandro Dalcin   Not collective
1313f68a32c8SEmil Constantinescu 
1314f68a32c8SEmil Constantinescu   Input Parameter:
1315f68a32c8SEmil Constantinescu .  ts - timestepping context
1316f68a32c8SEmil Constantinescu 
1317f68a32c8SEmil Constantinescu   Output Parameter:
1318f68a32c8SEmil Constantinescu .  rktype - type of RK-scheme
1319f68a32c8SEmil Constantinescu 
1320f68a32c8SEmil Constantinescu   Level: intermediate
1321f68a32c8SEmil Constantinescu 
13222ea87ba9SLisandro Dalcin .seealso: TSRKSetType()
1323f68a32c8SEmil Constantinescu @*/
1324f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype)
1325f68a32c8SEmil Constantinescu {
1326f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
1327f68a32c8SEmil Constantinescu 
1328f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1329f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1330f68a32c8SEmil Constantinescu   ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr);
1331f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1332f68a32c8SEmil Constantinescu }
1333f68a32c8SEmil Constantinescu 
13342ea87ba9SLisandro Dalcin static PetscErrorCode TSRKGetOrder_RK(TS ts,PetscInt *order)
13352ea87ba9SLisandro Dalcin {
13362ea87ba9SLisandro Dalcin   TS_RK *rk = (TS_RK*)ts->data;
13372ea87ba9SLisandro Dalcin 
13382ea87ba9SLisandro Dalcin   PetscFunctionBegin;
13392ea87ba9SLisandro Dalcin   *order = rk->tableau->order;
13402ea87ba9SLisandro Dalcin   PetscFunctionReturn(0);
13412ea87ba9SLisandro Dalcin }
13422ea87ba9SLisandro Dalcin 
1343560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype)
1344f68a32c8SEmil Constantinescu {
1345f68a32c8SEmil Constantinescu   TS_RK *rk = (TS_RK*)ts->data;
1346f68a32c8SEmil Constantinescu 
1347f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1348f68a32c8SEmil Constantinescu   *rktype = rk->tableau->name;
1349f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1350f68a32c8SEmil Constantinescu }
13512c9cb062Sluzhanghpp 
1352560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype)
1353f68a32c8SEmil Constantinescu {
1354f68a32c8SEmil Constantinescu   TS_RK          *rk = (TS_RK*)ts->data;
1355f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
1356f68a32c8SEmil Constantinescu   PetscBool      match;
1357f68a32c8SEmil Constantinescu   RKTableauLink  link;
1358f68a32c8SEmil Constantinescu 
1359f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1360f68a32c8SEmil Constantinescu   if (rk->tableau) {
1361f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr);
1362f68a32c8SEmil Constantinescu     if (match) PetscFunctionReturn(0);
1363f68a32c8SEmil Constantinescu   }
1364f68a32c8SEmil Constantinescu   for (link = RKTableauList; link; link=link->next) {
1365f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr);
1366f68a32c8SEmil Constantinescu     if (match) {
1367be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);}
1368f68a32c8SEmil Constantinescu       rk->tableau = &link->tab;
1369be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);}
13702ffb9264SLisandro Dalcin       ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
1371f68a32c8SEmil Constantinescu       PetscFunctionReturn(0);
1372f68a32c8SEmil Constantinescu     }
1373f68a32c8SEmil Constantinescu   }
1374f68a32c8SEmil Constantinescu   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype);
1375e4dd521cSBarry Smith }
1376e4dd521cSBarry Smith 
1377ff22ae23SHong Zhang static PetscErrorCode  TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y)
1378ff22ae23SHong Zhang {
1379ff22ae23SHong Zhang   TS_RK *rk = (TS_RK*)ts->data;
1380ff22ae23SHong Zhang 
1381ff22ae23SHong Zhang   PetscFunctionBegin;
13820429704eSStefano Zampini   if (ns) *ns = rk->tableau->s;
1383d2daff3dSHong Zhang   if (Y)   *Y = rk->Y;
1384ff22ae23SHong Zhang   PetscFunctionReturn(0);
1385ff22ae23SHong Zhang }
1386ff22ae23SHong Zhang 
1387b3a6b972SJed Brown static PetscErrorCode TSDestroy_RK(TS ts)
1388b3a6b972SJed Brown {
1389b3a6b972SJed Brown   PetscErrorCode ierr;
1390b3a6b972SJed Brown 
1391b3a6b972SJed Brown   PetscFunctionBegin;
1392b3a6b972SJed Brown   ierr = TSReset_RK(ts);CHKERRQ(ierr);
1393b3a6b972SJed Brown   if (ts->dm) {
1394b3a6b972SJed Brown     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
1395b3a6b972SJed Brown     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
1396b3a6b972SJed Brown   }
1397b3a6b972SJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
13982ea87ba9SLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetOrder_C",NULL);CHKERRQ(ierr);
1399b3a6b972SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr);
1400b3a6b972SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr);
14010f7a28cdSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetTableau_C",NULL);CHKERRQ(ierr);
14020fe4d17eSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",NULL);CHKERRQ(ierr);
14030fe4d17eSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",NULL);CHKERRQ(ierr);
1404b3a6b972SJed Brown   PetscFunctionReturn(0);
1405b3a6b972SJed Brown }
1406ff22ae23SHong Zhang 
14073ca0882eSHong Zhang /*
14083ca0882eSHong Zhang   This defines the nonlinear equation that is to be solved with SNES
14093ca0882eSHong Zhang   We do not need to solve the equation; we just use SNES to approximate the Jacobian
14103ca0882eSHong Zhang */
14113ca0882eSHong Zhang static PetscErrorCode SNESTSFormFunction_RK(SNES snes,Vec x,Vec y,TS ts)
14123ca0882eSHong Zhang {
14133ca0882eSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
14143ca0882eSHong Zhang   DM             dm,dmsave;
14153ca0882eSHong Zhang   PetscErrorCode ierr;
14163ca0882eSHong Zhang 
14173ca0882eSHong Zhang   PetscFunctionBegin;
14183ca0882eSHong Zhang   ierr   = SNESGetDM(snes,&dm);CHKERRQ(ierr);
14193ca0882eSHong Zhang   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
14203ca0882eSHong Zhang   dmsave = ts->dm;
14213ca0882eSHong Zhang   ts->dm = dm;
14223ca0882eSHong Zhang   ierr   = TSComputeRHSFunction(ts,rk->stage_time,x,y);CHKERRQ(ierr);
14233ca0882eSHong Zhang   ts->dm = dmsave;
14243ca0882eSHong Zhang   PetscFunctionReturn(0);
14253ca0882eSHong Zhang }
14263ca0882eSHong Zhang 
14273ca0882eSHong Zhang static PetscErrorCode SNESTSFormJacobian_RK(SNES snes,Vec x,Mat A,Mat B,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   dmsave = ts->dm;
14363ca0882eSHong Zhang   ts->dm = dm;
14373ca0882eSHong Zhang   ierr   = TSComputeRHSJacobian(ts,rk->stage_time,x,A,B);CHKERRQ(ierr);
14383ca0882eSHong Zhang   ts->dm = dmsave;
14393ca0882eSHong Zhang   PetscFunctionReturn(0);
14403ca0882eSHong Zhang }
14413ca0882eSHong Zhang 
144221052c0cSHong Zhang /*@C
144321052c0cSHong Zhang   TSRKSetMultirate - Use the interpolation-based multirate RK method
144421052c0cSHong Zhang 
144521052c0cSHong Zhang   Logically collective
144621052c0cSHong Zhang 
1447*d8d19677SJose E. Roman   Input Parameters:
144821052c0cSHong Zhang +  ts - timestepping context
144921052c0cSHong 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
145021052c0cSHong Zhang 
145121052c0cSHong Zhang   Options Database:
145221052c0cSHong Zhang .   -ts_rk_multirate - <true,false>
145321052c0cSHong Zhang 
145421052c0cSHong Zhang   Notes:
145521052c0cSHong 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).
145621052c0cSHong Zhang 
145721052c0cSHong Zhang   Level: intermediate
145821052c0cSHong Zhang 
145921052c0cSHong Zhang .seealso: TSRKGetMultirate()
146021052c0cSHong Zhang @*/
146121052c0cSHong Zhang PetscErrorCode TSRKSetMultirate(TS ts,PetscBool use_multirate)
146221052c0cSHong Zhang {
1463f06fb28fSHong Zhang   PetscErrorCode ierr;
14647dbacdf2SHong Zhang 
14658a4be4dbSHong Zhang   PetscFunctionBegin;
1466f06fb28fSHong Zhang   ierr = PetscTryMethod(ts,"TSRKSetMultirate_C",(TS,PetscBool),(ts,use_multirate));CHKERRQ(ierr);
146721052c0cSHong Zhang   PetscFunctionReturn(0);
146821052c0cSHong Zhang }
146921052c0cSHong Zhang 
147021052c0cSHong Zhang /*@C
147121052c0cSHong Zhang   TSRKGetMultirate - Gets whether to Use the interpolation-based multirate RK method
147221052c0cSHong Zhang 
147321052c0cSHong Zhang   Not collective
147421052c0cSHong Zhang 
147521052c0cSHong Zhang   Input Parameter:
147621052c0cSHong Zhang .  ts - timestepping context
147721052c0cSHong Zhang 
147821052c0cSHong Zhang   Output Parameter:
147921052c0cSHong Zhang .  use_multirate - PETSC_TRUE if the multirate RK method is enabled, PETSC_FALSE otherwise
148021052c0cSHong Zhang 
148121052c0cSHong Zhang   Level: intermediate
148221052c0cSHong Zhang 
148321052c0cSHong Zhang .seealso: TSRKSetMultirate()
148421052c0cSHong Zhang @*/
148521052c0cSHong Zhang PetscErrorCode TSRKGetMultirate(TS ts,PetscBool *use_multirate)
148621052c0cSHong Zhang {
1487f06fb28fSHong Zhang   PetscErrorCode ierr;
14887dbacdf2SHong Zhang 
14897dbacdf2SHong Zhang   PetscFunctionBegin;
1490f06fb28fSHong Zhang   ierr = PetscUseMethod(ts,"TSRKGetMultirate_C",(TS,PetscBool*),(ts,use_multirate));CHKERRQ(ierr);
149121052c0cSHong Zhang   PetscFunctionReturn(0);
149221052c0cSHong Zhang }
149321052c0cSHong Zhang 
1494ebe3b25bSBarry Smith /*MC
1495f68a32c8SEmil Constantinescu       TSRK - ODE and DAE solver using Runge-Kutta schemes
1496ebe3b25bSBarry Smith 
1497f68a32c8SEmil Constantinescu   The user should provide the right hand side of the equation
1498f68a32c8SEmil Constantinescu   using TSSetRHSFunction().
1499ebe3b25bSBarry Smith 
1500f68a32c8SEmil Constantinescu   Notes:
150198164e13SLisandro Dalcin   The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type
1502ebe3b25bSBarry Smith 
1503d41222bbSBarry Smith   Level: beginner
1504d41222bbSBarry Smith 
15055d1808f1SStefano Zampini .seealso:  TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRK2D, TTSRK2E, TSRK3,
15060fe4d17eSHong Zhang            TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister(), TSRKSetMultirate(), TSRKGetMultirate()
1507ebe3b25bSBarry Smith 
1508ebe3b25bSBarry Smith M*/
15098cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
1510e4dd521cSBarry Smith {
15111566a47fSLisandro Dalcin   TS_RK          *rk;
1512dfbe8321SBarry Smith   PetscErrorCode ierr;
1513e4dd521cSBarry Smith 
1514e4dd521cSBarry Smith   PetscFunctionBegin;
1515f68a32c8SEmil Constantinescu   ierr = TSRKInitializePackage();CHKERRQ(ierr);
1516f68a32c8SEmil Constantinescu 
1517f68a32c8SEmil Constantinescu   ts->ops->reset          = TSReset_RK;
15185f70b478SJed Brown   ts->ops->destroy        = TSDestroy_RK;
15195f70b478SJed Brown   ts->ops->view           = TSView_RK;
1520f68a32c8SEmil Constantinescu   ts->ops->load           = TSLoad_RK;
1521f68a32c8SEmil Constantinescu   ts->ops->setup          = TSSetUp_RK;
1522f68a32c8SEmil Constantinescu   ts->ops->interpolate    = TSInterpolate_RK;
15232c9cb062Sluzhanghpp   ts->ops->step           = TSStep_RK;
1524f68a32c8SEmil Constantinescu   ts->ops->evaluatestep   = TSEvaluateStep_RK;
1525fecfb714SLisandro Dalcin   ts->ops->rollback       = TSRollBack_RK;
1526f68a32c8SEmil Constantinescu   ts->ops->setfromoptions = TSSetFromOptions_RK;
1527ff22ae23SHong Zhang   ts->ops->getstages      = TSGetStages_RK;
1528e4dd521cSBarry Smith 
15293ca0882eSHong Zhang   ts->ops->snesfunction    = SNESTSFormFunction_RK;
15303ca0882eSHong Zhang   ts->ops->snesjacobian    = SNESTSFormJacobian_RK;
15310f7a1166SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_RK;
153213af1a74SHong Zhang   ts->ops->adjointsetup    = TSAdjointSetUp_RK;
153313af1a74SHong Zhang   ts->ops->adjointstep     = TSAdjointStep_RK;
153413af1a74SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_RK;
15350f7a1166SHong Zhang 
153613af1a74SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_RK;
1537922a638cSHong Zhang   ts->ops->forwardsetup    = TSForwardSetUp_RK;
1538922a638cSHong Zhang   ts->ops->forwardreset    = TSForwardReset_RK;
1539922a638cSHong Zhang   ts->ops->forwardstep     = TSForwardStep_RK;
1540922a638cSHong Zhang   ts->ops->forwardgetstages= TSForwardGetStages_RK;
1541922a638cSHong Zhang 
15421566a47fSLisandro Dalcin   ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr);
15431566a47fSLisandro Dalcin   ts->data = (void*)rk;
1544e4dd521cSBarry Smith 
15452ea87ba9SLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetOrder_C",TSRKGetOrder_RK);CHKERRQ(ierr);
1546f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr);
1547f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr);
15480f7a28cdSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetTableau_C",TSRKGetTableau_RK);CHKERRQ(ierr);
154921052c0cSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",TSRKSetMultirate_RK);CHKERRQ(ierr);
155021052c0cSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",TSRKGetMultirate_RK);CHKERRQ(ierr);
1551be5899b3SLisandro Dalcin 
1552be5899b3SLisandro Dalcin   ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr);
155390b67129SHong Zhang   rk->dtratio = 1;
1554e4dd521cSBarry Smith   PetscFunctionReturn(0);
1555e4dd521cSBarry Smith }
1556