xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision 3ca0882ebe6e6bc3c96f853b3244799ed8a662a6)
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 
14137acaa02SHendrik Ranocha      This method has ten stages with the First Same As Last (FSAL) property.
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 
15537acaa02SHendrik Ranocha      This method has thirteen stages with the First Same As Last (FSAL) property.
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 
257a685a763Sluzhanghpp         ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,5,binterp[0]);CHKERRQ(ierr);
258f68a32c8SEmil Constantinescu   }
25905e23783SLisandro Dalcin   {
26005e23783SLisandro Dalcin     const PetscReal
26117f6384fSBarry Smith       A[8][8]   = {{0,0,0,0,0,0,0,0},
2624e82626cSLisandro Dalcin                    {RC(1.0)/RC(6.0),0,0,0,0,0,0,0},
2634e82626cSLisandro Dalcin                    {RC(2.0)/RC(27.0),RC(4.0)/RC(27.0),0,0,0,0,0,0},
2644e82626cSLisandro Dalcin                    {RC(183.0)/RC(1372.0),RC(-162.0)/RC(343.0),RC(1053.0)/RC(1372.0),0,0,0,0,0},
2654e82626cSLisandro 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},
2664e82626cSLisandro 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},
2674e82626cSLisandro 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},
2684e82626cSLisandro 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}},
2694e82626cSLisandro 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},
2704e82626cSLisandro 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)};
27105e23783SLisandro Dalcin     ierr = TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
27205e23783SLisandro Dalcin   }
27337acaa02SHendrik Ranocha   {
27437acaa02SHendrik Ranocha     const PetscReal
27537acaa02SHendrik Ranocha       A[9][9]   = {{0,0,0,0,0,0,0,0,0},
27663fe90e8SHendrik Ranocha                    {RC(1.8000000000000000000000000000000000000000e-01),0,0,0,0,0,0,0,0},
27763fe90e8SHendrik Ranocha                    {RC(8.9506172839506172839506172839506172839506e-02),RC(7.7160493827160493827160493827160493827160e-02),0,0,0,0,0,0,0},
27863fe90e8SHendrik Ranocha                    {RC(6.2500000000000000000000000000000000000000e-02),0,RC(1.8750000000000000000000000000000000000000e-01),0,0,0,0,0,0},
27963fe90e8SHendrik Ranocha                    {RC(3.1651600000000000000000000000000000000000e-01),0,RC(-1.0449480000000000000000000000000000000000e+00),RC(1.2584320000000000000000000000000000000000e+00),0,0,0,0,0},
28063fe90e8SHendrik Ranocha                    {RC(2.7232612736485626257225065566674305502508e-01),0,RC(-8.2513360323886639676113360323886639676113e-01),RC(1.0480917678812415654520917678812415654521e+00),RC(1.0471570799276856873679117969088177628396e-01),0,0,0,0},
28163fe90e8SHendrik 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},
28263fe90e8SHendrik 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},
28363fe90e8SHendrik 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}},
28463fe90e8SHendrik 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},
28563fe90e8SHendrik 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)};
28637acaa02SHendrik Ranocha     ierr = TSRKRegister(TSRK6VR,6,9,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
28737acaa02SHendrik Ranocha   }
28837acaa02SHendrik Ranocha   {
28937acaa02SHendrik Ranocha     const PetscReal
29037acaa02SHendrik Ranocha       A[10][10]  = {{0,0,0,0,0,0,0,0,0,0},
29163fe90e8SHendrik Ranocha                     {RC(5.0000000000000000000000000000000000000000e-03),0,0,0,0,0,0,0,0,0},
29263fe90e8SHendrik Ranocha                     {RC(-1.0767901234567901234567901234567901234568e+00),RC(1.1856790123456790123456790123456790123457e+00),0,0,0,0,0,0,0,0},
29363fe90e8SHendrik Ranocha                     {RC(4.0833333333333333333333333333333333333333e-02),0,RC(1.2250000000000000000000000000000000000000e-01),0,0,0,0,0,0,0},
29463fe90e8SHendrik Ranocha                     {RC(6.3607142857142857142857142857142857142857e-01),0,RC(-2.4444642857142857142857142857142857142857e+00),RC(2.2633928571428571428571428571428571428571e+00),0,0,0,0,0,0},
29563fe90e8SHendrik Ranocha                     {RC(-2.5351211079349245229256383554660215487207e+00),0,RC(1.0299374654449267920438514460756024913612e+01),RC(-7.9513032885990579949493217458266876536482e+00),RC(7.9301148923100592201226014271115261823800e-01),0,0,0,0,0},
29663fe90e8SHendrik 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},
29763fe90e8SHendrik 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},
29863fe90e8SHendrik 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},
29963fe90e8SHendrik 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}},
30063fe90e8SHendrik 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},
30163fe90e8SHendrik 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)};
30237acaa02SHendrik Ranocha     ierr = TSRKRegister(TSRK7VR,7,10,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
30337acaa02SHendrik Ranocha   }
30437acaa02SHendrik Ranocha   {
30537acaa02SHendrik Ranocha     const PetscReal
30637acaa02SHendrik Ranocha       A[13][13]  = {{0,0,0,0,0,0,0,0,0,0,0,0,0},
30763fe90e8SHendrik Ranocha                     {RC(2.5000000000000000000000000000000000000000e-01),0,0,0,0,0,0,0,0,0,0,0,0},
30863fe90e8SHendrik Ranocha                     {RC(8.7400846504915232052686327594877411977046e-02),RC(2.5487604938654321753087950620345685135815e-02),0,0,0,0,0,0,0,0,0,0,0},
30963fe90e8SHendrik Ranocha                     {RC(4.2333169291338582677165354330708661417323e-02),0,RC(1.2699950787401574803149606299212598425197e-01),0,0,0,0,0,0,0,0,0,0},
31063fe90e8SHendrik Ranocha                     {RC(4.2609505888742261494881445237572274090942e-01),0,RC(-1.5987952846591523265427733230657181117089e+00),RC(1.5967002257717297115939588706899953707994e+00),0,0,0,0,0,0,0,0,0},
31163fe90e8SHendrik Ranocha                     {RC(5.0719337296713929515090618138513639239329e-02),0,0,RC(2.5433377264600407582754714408877778031369e-01),RC(2.0394689005728199465736223777270858044698e-01),0,0,0,0,0,0,0,0},
31263fe90e8SHendrik 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},
31363fe90e8SHendrik 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},
31463fe90e8SHendrik 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},
31563fe90e8SHendrik 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},
31663fe90e8SHendrik 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},
31763fe90e8SHendrik 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},
31863fe90e8SHendrik 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}},
31963fe90e8SHendrik 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},
32063fe90e8SHendrik 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)};
32137acaa02SHendrik Ranocha     ierr = TSRKRegister(TSRK8VR,8,13,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
32237acaa02SHendrik Ranocha   }
3234e82626cSLisandro Dalcin #undef RC
324db046809SBarry Smith   PetscFunctionReturn(0);
325db046809SBarry Smith }
326db046809SBarry Smith 
327f68a32c8SEmil Constantinescu /*@C
328f68a32c8SEmil Constantinescu    TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister().
329f68a32c8SEmil Constantinescu 
330f68a32c8SEmil Constantinescu    Not Collective
331f68a32c8SEmil Constantinescu 
332f68a32c8SEmil Constantinescu    Level: advanced
333f68a32c8SEmil Constantinescu 
334f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll()
335f68a32c8SEmil Constantinescu @*/
336f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void)
337e4dd521cSBarry Smith {
338f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
339f68a32c8SEmil Constantinescu   RKTableauLink  link;
340f68a32c8SEmil Constantinescu 
341f68a32c8SEmil Constantinescu   PetscFunctionBegin;
342f68a32c8SEmil Constantinescu   while ((link = RKTableauList)) {
343f68a32c8SEmil Constantinescu     RKTableau t = &link->tab;
344f68a32c8SEmil Constantinescu     RKTableauList = link->next;
345f68a32c8SEmil Constantinescu     ierr = PetscFree3(t->A,t->b,t->c);  CHKERRQ(ierr);
346f68a32c8SEmil Constantinescu     ierr = PetscFree (t->bembed);       CHKERRQ(ierr);
347f68a32c8SEmil Constantinescu     ierr = PetscFree (t->binterp);      CHKERRQ(ierr);
348f68a32c8SEmil Constantinescu     ierr = PetscFree (t->name);         CHKERRQ(ierr);
349f68a32c8SEmil Constantinescu     ierr = PetscFree (link);            CHKERRQ(ierr);
350f68a32c8SEmil Constantinescu   }
351f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_FALSE;
352f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
353f68a32c8SEmil Constantinescu }
354f68a32c8SEmil Constantinescu 
355f68a32c8SEmil Constantinescu /*@C
356f68a32c8SEmil Constantinescu   TSRKInitializePackage - This function initializes everything in the TSRK package. It is called
3578a690491SBarry Smith   from TSInitializePackage().
358f68a32c8SEmil Constantinescu 
359f68a32c8SEmil Constantinescu   Level: developer
360f68a32c8SEmil Constantinescu 
361f68a32c8SEmil Constantinescu .seealso: PetscInitialize()
362f68a32c8SEmil Constantinescu @*/
363f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void)
364f68a32c8SEmil Constantinescu {
365186e87acSLisandro Dalcin   PetscErrorCode ierr;
366e4dd521cSBarry Smith 
367e4dd521cSBarry Smith   PetscFunctionBegin;
368f68a32c8SEmil Constantinescu   if (TSRKPackageInitialized) PetscFunctionReturn(0);
369f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_TRUE;
370f68a32c8SEmil Constantinescu   ierr = TSRKRegisterAll();CHKERRQ(ierr);
371f68a32c8SEmil Constantinescu   ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr);
372f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
373f68a32c8SEmil Constantinescu }
374186e87acSLisandro Dalcin 
375f68a32c8SEmil Constantinescu /*@C
376f68a32c8SEmil Constantinescu   TSRKFinalizePackage - This function destroys everything in the TSRK package. It is
377f68a32c8SEmil Constantinescu   called from PetscFinalize().
378186e87acSLisandro Dalcin 
379f68a32c8SEmil Constantinescu   Level: developer
380f68a32c8SEmil Constantinescu 
381f68a32c8SEmil Constantinescu .seealso: PetscFinalize()
382f68a32c8SEmil Constantinescu @*/
383f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void)
384f68a32c8SEmil Constantinescu {
385f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
386f68a32c8SEmil Constantinescu 
387f68a32c8SEmil Constantinescu   PetscFunctionBegin;
388f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_FALSE;
389f68a32c8SEmil Constantinescu   ierr = TSRKRegisterDestroy();CHKERRQ(ierr);
390f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
391f68a32c8SEmil Constantinescu }
392f68a32c8SEmil Constantinescu 
393f68a32c8SEmil Constantinescu /*@C
394f68a32c8SEmil Constantinescu    TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
395f68a32c8SEmil Constantinescu 
396f68a32c8SEmil Constantinescu    Not Collective, but the same schemes should be registered on all processes on which they will be used
397f68a32c8SEmil Constantinescu 
398f68a32c8SEmil Constantinescu    Input Parameters:
399f68a32c8SEmil Constantinescu +  name - identifier for method
400f68a32c8SEmil Constantinescu .  order - approximation order of method
401f68a32c8SEmil Constantinescu .  s - number of stages, this is the dimension of the matrices below
402f68a32c8SEmil Constantinescu .  A - stage coefficients (dimension s*s, row-major)
403f68a32c8SEmil Constantinescu .  b - step completion table (dimension s; NULL to use last row of A)
404f68a32c8SEmil Constantinescu .  c - abscissa (dimension s; NULL to use row sums of A)
405f68a32c8SEmil Constantinescu .  bembed - completion table for embedded method (dimension s; NULL if not available)
4063a8a9803SLisandro Dalcin .  p - Order of the interpolation scheme, equal to the number of columns of binterp
4073a8a9803SLisandro Dalcin -  binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1)
408f68a32c8SEmil Constantinescu 
409f68a32c8SEmil Constantinescu    Notes:
410f68a32c8SEmil Constantinescu    Several RK methods are provided, this function is only needed to create new methods.
411f68a32c8SEmil Constantinescu 
412f68a32c8SEmil Constantinescu    Level: advanced
413f68a32c8SEmil Constantinescu 
414f68a32c8SEmil Constantinescu .seealso: TSRK
415f68a32c8SEmil Constantinescu @*/
416f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s,
417f68a32c8SEmil Constantinescu                             const PetscReal A[],const PetscReal b[],const PetscReal c[],
4183a8a9803SLisandro Dalcin                             const PetscReal bembed[],PetscInt p,const PetscReal binterp[])
419f68a32c8SEmil Constantinescu {
420f68a32c8SEmil Constantinescu   PetscErrorCode  ierr;
421f68a32c8SEmil Constantinescu   RKTableauLink   link;
422f68a32c8SEmil Constantinescu   RKTableau       t;
423f68a32c8SEmil Constantinescu   PetscInt        i,j;
424f68a32c8SEmil Constantinescu 
425f68a32c8SEmil Constantinescu   PetscFunctionBegin;
4263a8a9803SLisandro Dalcin   PetscValidCharPointer(name,1);
4273a8a9803SLisandro Dalcin   PetscValidRealPointer(A,4);
4283a8a9803SLisandro Dalcin   if (b) PetscValidRealPointer(b,5);
4293a8a9803SLisandro Dalcin   if (c) PetscValidRealPointer(c,6);
4303a8a9803SLisandro Dalcin   if (bembed) PetscValidRealPointer(bembed,7);
4313a8a9803SLisandro Dalcin   if (binterp || p > 1) PetscValidRealPointer(binterp,9);
4323a8a9803SLisandro Dalcin 
4331d36bdfdSBarry Smith   ierr = TSRKInitializePackage();CHKERRQ(ierr);
43495dccacaSBarry Smith   ierr = PetscNew(&link);CHKERRQ(ierr);
435f68a32c8SEmil Constantinescu   t = &link->tab;
4363a8a9803SLisandro Dalcin 
437f68a32c8SEmil Constantinescu   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
438f68a32c8SEmil Constantinescu   t->order = order;
439f68a32c8SEmil Constantinescu   t->s = s;
440dcca6d9dSJed Brown   ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr);
441580bdb30SBarry Smith   ierr = PetscArraycpy(t->A,A,s*s);CHKERRQ(ierr);
442580bdb30SBarry Smith   if (b)  { ierr = PetscArraycpy(t->b,b,s);CHKERRQ(ierr); }
443f68a32c8SEmil Constantinescu   else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i];
444580bdb30SBarry Smith   if (c)  { ierr = PetscArraycpy(t->c,c,s);CHKERRQ(ierr); }
445f68a32c8SEmil 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];
446d760c35bSDebojyoti Ghosh   t->FSAL = PETSC_TRUE;
447d760c35bSDebojyoti Ghosh   for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE;
4483a8a9803SLisandro Dalcin 
449f68a32c8SEmil Constantinescu   if (bembed) {
450785e854fSJed Brown     ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr);
451580bdb30SBarry Smith     ierr = PetscArraycpy(t->bembed,bembed,s);CHKERRQ(ierr);
452f68a32c8SEmil Constantinescu   }
453f68a32c8SEmil Constantinescu 
4543a8a9803SLisandro Dalcin   if (!binterp) { p = 1; binterp = t->b; }
4553a8a9803SLisandro Dalcin   t->p = p;
4563a8a9803SLisandro Dalcin   ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr);
457580bdb30SBarry Smith   ierr = PetscArraycpy(t->binterp,binterp,s*p);CHKERRQ(ierr);
4583a8a9803SLisandro Dalcin 
459f68a32c8SEmil Constantinescu   link->next = RKTableauList;
460f68a32c8SEmil Constantinescu   RKTableauList = link;
461f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
462f68a32c8SEmil Constantinescu }
463f68a32c8SEmil Constantinescu 
464e4dd521cSBarry Smith /*
4652c9cb062Sluzhanghpp  This is for single-step RK method
466f68a32c8SEmil Constantinescu  The step completion formula is
467e4dd521cSBarry Smith 
468f68a32c8SEmil Constantinescu  x1 = x0 + h b^T YdotRHS
469f68a32c8SEmil Constantinescu 
470f68a32c8SEmil Constantinescu  This function can be called before or after ts->vec_sol has been updated.
471f68a32c8SEmil Constantinescu  Suppose we have a completion formula (b) and an embedded formula (be) of different order.
472f68a32c8SEmil Constantinescu  We can write
473f68a32c8SEmil Constantinescu 
474f68a32c8SEmil Constantinescu  x1e = x0 + h be^T YdotRHS
475f68a32c8SEmil Constantinescu      = x1 - h b^T YdotRHS + h be^T YdotRHS
476f68a32c8SEmil Constantinescu      = x1 + h (be - b)^T YdotRHS
477f68a32c8SEmil Constantinescu 
478f68a32c8SEmil Constantinescu  so we can evaluate the method with different order even after the step has been optimistically completed.
479e4dd521cSBarry Smith */
480f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done)
481f68a32c8SEmil Constantinescu {
482f68a32c8SEmil Constantinescu   TS_RK          *rk   = (TS_RK*)ts->data;
483f68a32c8SEmil Constantinescu   RKTableau      tab  = rk->tableau;
484f68a32c8SEmil Constantinescu   PetscScalar    *w    = rk->work;
485f68a32c8SEmil Constantinescu   PetscReal      h;
486f68a32c8SEmil Constantinescu   PetscInt       s    = tab->s,j;
487f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
488f68a32c8SEmil Constantinescu 
489f68a32c8SEmil Constantinescu   PetscFunctionBegin;
490f68a32c8SEmil Constantinescu   switch (rk->status) {
491f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
492f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
493f68a32c8SEmil Constantinescu     h = ts->time_step; break;
494f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
495be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
496f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
497f68a32c8SEmil Constantinescu   }
498f68a32c8SEmil Constantinescu   if (order == tab->order) {
499f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) {
500f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
50190b67129SHong Zhang       for (j=0; j<s; j++) w[j] = h*tab->b[j]/rk->dtratio;
502f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
503f68a32c8SEmil Constantinescu     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
504f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
505f68a32c8SEmil Constantinescu   } else if (order == tab->order-1) {
506f68a32c8SEmil Constantinescu     if (!tab->bembed) goto unavailable;
507f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) { /*Complete with the embedded method (be)*/
508f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
509f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
510f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
511f68a32c8SEmil Constantinescu     } else {  /*Rollback and re-complete using (be-b) */
512f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
513f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
514f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
515f68a32c8SEmil Constantinescu     }
516f68a32c8SEmil Constantinescu     if (done) *done = PETSC_TRUE;
517f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
518f68a32c8SEmil Constantinescu   }
519f68a32c8SEmil Constantinescu unavailable:
520f68a32c8SEmil Constantinescu   if (done) *done = PETSC_FALSE;
521a7fac7c2SEmil 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);
522f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
523f68a32c8SEmil Constantinescu }
524f68a32c8SEmil Constantinescu 
5250f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts)
5260f7a1166SHong Zhang {
5270f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
528cd4cee2dSHong Zhang   TS              quadts = ts->quadraturets;
5290f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
5300f7a1166SHong Zhang   const PetscInt  s = tab->s;
5310f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
5320f7a1166SHong Zhang   Vec             *Y = rk->Y;
5330f7a1166SHong Zhang   PetscInt        i;
5340f7a1166SHong Zhang   PetscErrorCode  ierr;
5350f7a1166SHong Zhang 
5360f7a1166SHong Zhang   PetscFunctionBegin;
537cd4cee2dSHong Zhang   /* No need to backup quadts->vec_sol since it can be reverted in TSRollBack_RK */
5380f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
539cd4cee2dSHong Zhang     /* Evolve quadrature TS solution to compute integrals */
540cd4cee2dSHong Zhang     ierr = TSComputeRHSFunction(quadts,rk->ptime+rk->time_step*c[i],Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
541cd4cee2dSHong Zhang     ierr = VecAXPY(quadts->vec_sol,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
5420f7a1166SHong Zhang   }
5430f7a1166SHong Zhang   PetscFunctionReturn(0);
5440f7a1166SHong Zhang }
5450f7a1166SHong Zhang 
5460f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts)
5470f7a1166SHong Zhang {
5480f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
5490f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
550cd4cee2dSHong Zhang   TS              quadts = ts->quadraturets;
5510f7a1166SHong Zhang   const PetscInt  s = tab->s;
5520f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
5530f7a1166SHong Zhang   Vec             *Y = rk->Y;
5540f7a1166SHong Zhang   PetscInt        i;
5550f7a1166SHong Zhang   PetscErrorCode  ierr;
5560f7a1166SHong Zhang 
5570f7a1166SHong Zhang   PetscFunctionBegin;
5580f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
559cd4cee2dSHong Zhang     /* Evolve quadrature TS solution to compute integrals */
560cd4cee2dSHong Zhang     ierr = TSComputeRHSFunction(quadts,ts->ptime+ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
561cd4cee2dSHong Zhang     ierr = VecAXPY(quadts->vec_sol,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
5620f7a1166SHong Zhang   }
5630f7a1166SHong Zhang   PetscFunctionReturn(0);
5640f7a1166SHong Zhang }
5650f7a1166SHong Zhang 
566fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts)
567fecfb714SLisandro Dalcin {
568fecfb714SLisandro Dalcin   TS_RK           *rk = (TS_RK*)ts->data;
569cd4cee2dSHong Zhang   TS              quadts = ts->quadraturets;
570fecfb714SLisandro Dalcin   RKTableau       tab = rk->tableau;
571fecfb714SLisandro Dalcin   const PetscInt  s  = tab->s;
572cd4cee2dSHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
573fecfb714SLisandro Dalcin   PetscScalar     *w = rk->work;
574cd4cee2dSHong Zhang   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
575fecfb714SLisandro Dalcin   PetscInt        j;
576fecfb714SLisandro Dalcin   PetscReal       h;
577fecfb714SLisandro Dalcin   PetscErrorCode  ierr;
578fecfb714SLisandro Dalcin 
579fecfb714SLisandro Dalcin   PetscFunctionBegin;
580fecfb714SLisandro Dalcin   switch (rk->status) {
581fecfb714SLisandro Dalcin   case TS_STEP_INCOMPLETE:
582fecfb714SLisandro Dalcin   case TS_STEP_PENDING:
583fecfb714SLisandro Dalcin     h = ts->time_step; break;
584fecfb714SLisandro Dalcin   case TS_STEP_COMPLETE:
585fecfb714SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
586fecfb714SLisandro Dalcin   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
587fecfb714SLisandro Dalcin   }
588fecfb714SLisandro Dalcin   for (j=0; j<s; j++) w[j] = -h*b[j];
589fecfb714SLisandro Dalcin   ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr);
590cd4cee2dSHong Zhang   if (quadts && ts->costintegralfwd) {
591cd4cee2dSHong Zhang     for (j=0; j<s; j++) {
592cd4cee2dSHong Zhang       /* Revert the quadrature TS solution */
593cd4cee2dSHong Zhang       ierr = TSComputeRHSFunction(quadts,rk->ptime+h*c[j],Y[j],ts->vec_costintegrand);CHKERRQ(ierr);
594cd4cee2dSHong Zhang       ierr = VecAXPY(quadts->vec_sol,-h*b[j],ts->vec_costintegrand);CHKERRQ(ierr);
595cd4cee2dSHong Zhang     }
596cd4cee2dSHong Zhang   }
597fecfb714SLisandro Dalcin   PetscFunctionReturn(0);
598fecfb714SLisandro Dalcin }
599fecfb714SLisandro Dalcin 
600922a638cSHong Zhang static PetscErrorCode TSForwardStep_RK(TS ts)
601922a638cSHong Zhang {
602922a638cSHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
603922a638cSHong Zhang   RKTableau       tab = rk->tableau;
604922a638cSHong Zhang   Mat             J,*MatsFwdSensipTemp = rk->MatsFwdSensipTemp;
605922a638cSHong Zhang   const PetscInt  s = tab->s;
606922a638cSHong Zhang   const PetscReal *A = tab->A,*c = tab->c,*b = tab->b;
607922a638cSHong Zhang   Vec             *Y = rk->Y;
608922a638cSHong Zhang   PetscInt        i,j;
609922a638cSHong Zhang   PetscReal       stage_time,h = ts->time_step;
610922a638cSHong Zhang   PetscBool       zero;
611922a638cSHong Zhang   PetscErrorCode  ierr;
612922a638cSHong Zhang 
613922a638cSHong Zhang   PetscFunctionBegin;
614922a638cSHong Zhang   ierr = MatCopy(ts->mat_sensip,rk->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
615922a638cSHong Zhang   ierr = TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr);
616922a638cSHong Zhang 
617922a638cSHong Zhang   for (i=0; i<s; i++) {
618922a638cSHong Zhang     stage_time = ts->ptime + h*c[i];
619922a638cSHong Zhang     zero = PETSC_FALSE;
620922a638cSHong Zhang     if (b[i] == 0 && i == s-1) zero = PETSC_TRUE;
621922a638cSHong Zhang     /* TLM Stage values */
622922a638cSHong Zhang     if(!i) {
623922a638cSHong Zhang       ierr = MatCopy(ts->mat_sensip,rk->MatsFwdStageSensip[i],SAME_NONZERO_PATTERN);CHKERRQ(ierr);
624922a638cSHong Zhang     } else if (!zero) {
625922a638cSHong Zhang       ierr = MatZeroEntries(rk->MatsFwdStageSensip[i]);CHKERRQ(ierr);
626922a638cSHong Zhang       for (j=0; j<i; j++) {
627922a638cSHong Zhang         ierr = MatAXPY(rk->MatsFwdStageSensip[i],h*A[i*s+j],MatsFwdSensipTemp[j],SAME_NONZERO_PATTERN);CHKERRQ(ierr);
628922a638cSHong Zhang       }
629922a638cSHong Zhang       ierr = MatAXPY(rk->MatsFwdStageSensip[i],1.,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
630922a638cSHong Zhang     } else {
631922a638cSHong Zhang       ierr = MatZeroEntries(rk->MatsFwdStageSensip[i]);CHKERRQ(ierr);
632922a638cSHong Zhang     }
633922a638cSHong Zhang 
634922a638cSHong Zhang     ierr = TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);CHKERRQ(ierr);
635922a638cSHong Zhang     ierr = MatMatMult(J,rk->MatsFwdStageSensip[i],MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatsFwdSensipTemp[i]);CHKERRQ(ierr);
63613af1a74SHong Zhang     if (ts->Jacprhs) {
63713af1a74SHong Zhang       ierr = TSComputeRHSJacobianP(ts,stage_time,Y[i],ts->Jacprhs);CHKERRQ(ierr); /* get f_p */
63813af1a74SHong Zhang       if (ts->vecs_sensi2p) { /* TLM used for 2nd-order adjoint */
63913af1a74SHong Zhang         PetscScalar *xarr;
64013af1a74SHong Zhang         ierr = MatDenseGetColumn(MatsFwdSensipTemp[i],0,&xarr);CHKERRQ(ierr);
64113af1a74SHong Zhang         ierr = VecPlaceArray(rk->VecDeltaFwdSensipCol,xarr);CHKERRQ(ierr);
64213af1a74SHong Zhang         ierr = MatMultAdd(ts->Jacprhs,ts->vec_dir,rk->VecDeltaFwdSensipCol,rk->VecDeltaFwdSensipCol);CHKERRQ(ierr);
64313af1a74SHong Zhang         ierr = VecResetArray(rk->VecDeltaFwdSensipCol);CHKERRQ(ierr);CHKERRQ(ierr);
64413af1a74SHong Zhang         ierr = MatDenseRestoreColumn(MatsFwdSensipTemp[i],&xarr);CHKERRQ(ierr);
64513af1a74SHong Zhang       } else {
64613af1a74SHong Zhang         ierr = MatAXPY(MatsFwdSensipTemp[i],1.,ts->Jacprhs,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
64713af1a74SHong Zhang       }
648922a638cSHong Zhang     }
649922a638cSHong Zhang   }
650922a638cSHong Zhang 
651922a638cSHong Zhang   for (i=0; i<s; i++) {
652922a638cSHong Zhang     ierr = MatAXPY(ts->mat_sensip,h*b[i],rk->MatsFwdSensipTemp[i],SAME_NONZERO_PATTERN);CHKERRQ(ierr);
653922a638cSHong Zhang   }
654922a638cSHong Zhang   rk->status = TS_STEP_COMPLETE;
655922a638cSHong Zhang   PetscFunctionReturn(0);
656922a638cSHong Zhang }
657922a638cSHong Zhang 
658922a638cSHong Zhang static PetscErrorCode TSForwardGetStages_RK(TS ts,PetscInt *ns,Mat **stagesensip)
659922a638cSHong Zhang {
660922a638cSHong Zhang   TS_RK     *rk = (TS_RK*)ts->data;
661922a638cSHong Zhang   RKTableau tab  = rk->tableau;
662922a638cSHong Zhang 
663922a638cSHong Zhang   PetscFunctionBegin;
664922a638cSHong Zhang   if (ns) *ns = tab->s;
665922a638cSHong Zhang   if (stagesensip) *stagesensip = rk->MatsFwdStageSensip;
666922a638cSHong Zhang   PetscFunctionReturn(0);
667922a638cSHong Zhang }
668922a638cSHong Zhang 
669922a638cSHong Zhang static PetscErrorCode TSForwardSetUp_RK(TS ts)
670922a638cSHong Zhang {
671922a638cSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
672922a638cSHong Zhang   RKTableau      tab  = rk->tableau;
673922a638cSHong Zhang   PetscInt       i;
674922a638cSHong Zhang   PetscErrorCode ierr;
675922a638cSHong Zhang 
676922a638cSHong Zhang   PetscFunctionBegin;
677922a638cSHong Zhang   /* backup sensitivity results for roll-backs */
678922a638cSHong Zhang   ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatFwdSensip0);CHKERRQ(ierr);
679922a638cSHong Zhang 
680922a638cSHong Zhang   ierr = PetscMalloc1(tab->s,&rk->MatsFwdStageSensip);CHKERRQ(ierr);
681922a638cSHong Zhang   ierr = PetscMalloc1(tab->s,&rk->MatsFwdSensipTemp);CHKERRQ(ierr);
682922a638cSHong Zhang   for(i=0; i<tab->s; i++) {
683922a638cSHong Zhang     ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatsFwdStageSensip[i]);CHKERRQ(ierr);
684922a638cSHong Zhang     ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatsFwdSensipTemp[i]);CHKERRQ(ierr);
685922a638cSHong Zhang   }
686922a638cSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&rk->VecDeltaFwdSensipCol);CHKERRQ(ierr);
687922a638cSHong Zhang   PetscFunctionReturn(0);
688922a638cSHong Zhang }
689922a638cSHong Zhang 
690922a638cSHong Zhang static PetscErrorCode TSForwardReset_RK(TS ts)
691922a638cSHong Zhang {
692922a638cSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
693922a638cSHong Zhang   RKTableau      tab  = rk->tableau;
694922a638cSHong Zhang   PetscInt       i;
695922a638cSHong Zhang   PetscErrorCode ierr;
696922a638cSHong Zhang 
697922a638cSHong Zhang   PetscFunctionBegin;
698922a638cSHong Zhang   ierr = MatDestroy(&rk->MatFwdSensip0);CHKERRQ(ierr);
699922a638cSHong Zhang   if (rk->MatsFwdStageSensip) {
700922a638cSHong Zhang     for (i=0; i<tab->s; i++) {
701922a638cSHong Zhang       ierr = MatDestroy(&rk->MatsFwdStageSensip[i]);CHKERRQ(ierr);
702922a638cSHong Zhang     }
703922a638cSHong Zhang     ierr = PetscFree(rk->MatsFwdStageSensip);CHKERRQ(ierr);
704922a638cSHong Zhang   }
705922a638cSHong Zhang   if (rk->MatsFwdSensipTemp) {
706922a638cSHong Zhang     for (i=0; i<tab->s; i++) {
707922a638cSHong Zhang       ierr = MatDestroy(&rk->MatsFwdSensipTemp[i]);CHKERRQ(ierr);
708922a638cSHong Zhang     }
709922a638cSHong Zhang     ierr = PetscFree(rk->MatsFwdSensipTemp);CHKERRQ(ierr);
710922a638cSHong Zhang   }
711922a638cSHong Zhang   ierr = VecDestroy(&rk->VecDeltaFwdSensipCol);CHKERRQ(ierr);
712922a638cSHong Zhang   PetscFunctionReturn(0);
713922a638cSHong Zhang }
714922a638cSHong Zhang 
715f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts)
716f68a32c8SEmil Constantinescu {
717f68a32c8SEmil Constantinescu   TS_RK           *rk  = (TS_RK*)ts->data;
718f68a32c8SEmil Constantinescu   RKTableau       tab  = rk->tableau;
719f68a32c8SEmil Constantinescu   const PetscInt  s = tab->s;
720fecfb714SLisandro Dalcin   const PetscReal *A = tab->A,*c = tab->c;
721f68a32c8SEmil Constantinescu   PetscScalar     *w = rk->work;
722f68a32c8SEmil Constantinescu   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
723d1905564SLisandro Dalcin   PetscBool       FSAL = tab->FSAL;
724f68a32c8SEmil Constantinescu   TSAdapt         adapt;
725fecfb714SLisandro Dalcin   PetscInt        i,j;
726be5899b3SLisandro Dalcin   PetscInt        rejections = 0;
727be5899b3SLisandro Dalcin   PetscBool       stageok,accept = PETSC_TRUE;
728be5899b3SLisandro Dalcin   PetscReal       next_time_step = ts->time_step;
729f68a32c8SEmil Constantinescu   PetscErrorCode  ierr;
730f68a32c8SEmil Constantinescu 
731f68a32c8SEmil Constantinescu   PetscFunctionBegin;
732d1905564SLisandro Dalcin   if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE;
733d1905564SLisandro Dalcin   if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); }
734d1905564SLisandro Dalcin 
735f68a32c8SEmil Constantinescu   rk->status = TS_STEP_INCOMPLETE;
736be5899b3SLisandro Dalcin   while (!ts->reason && rk->status != TS_STEP_COMPLETE) {
737be5899b3SLisandro Dalcin     PetscReal t = ts->ptime;
738f68a32c8SEmil Constantinescu     PetscReal h = ts->time_step;
739f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
7409be3e283SDebojyoti Ghosh       rk->stage_time = t + h*c[i];
7419be3e283SDebojyoti Ghosh       ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr);
742f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
743f68a32c8SEmil Constantinescu       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
744f68a32c8SEmil Constantinescu       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
7459be3e283SDebojyoti Ghosh       ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
746f68a32c8SEmil Constantinescu       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
747be5899b3SLisandro Dalcin       ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr);
748fecfb714SLisandro Dalcin       if (!stageok) goto reject_step;
7498f738a7cSLisandro Dalcin       if (FSAL && !i) continue;
750f68a32c8SEmil Constantinescu       ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
751f68a32c8SEmil Constantinescu     }
752be5899b3SLisandro Dalcin 
753be5899b3SLisandro Dalcin     rk->status = TS_STEP_INCOMPLETE;
754f68a32c8SEmil Constantinescu     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
755f68a32c8SEmil Constantinescu     rk->status = TS_STEP_PENDING;
756f68a32c8SEmil Constantinescu     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
757f68a32c8SEmil Constantinescu     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
7581917a363SLisandro Dalcin     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr);
759fecfb714SLisandro Dalcin     ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
760be5899b3SLisandro Dalcin     rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
761be5899b3SLisandro Dalcin     if (!accept) { /* Roll back the current step */
762fecfb714SLisandro Dalcin       ierr = TSRollBack_RK(ts);CHKERRQ(ierr);
763be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
764be5899b3SLisandro Dalcin       goto reject_step;
765be5899b3SLisandro Dalcin     }
766be5899b3SLisandro Dalcin 
7670f7a1166SHong Zhang     if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
7680f7a1166SHong Zhang       rk->ptime     = ts->ptime;
7690f7a1166SHong Zhang       rk->time_step = ts->time_step;
770493ed95fSHong Zhang     }
771be5899b3SLisandro Dalcin 
772f68a32c8SEmil Constantinescu     ts->ptime += ts->time_step;
773f68a32c8SEmil Constantinescu     ts->time_step = next_time_step;
774f68a32c8SEmil Constantinescu     break;
775be5899b3SLisandro Dalcin 
776be5899b3SLisandro Dalcin     reject_step:
777fecfb714SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
778be5899b3SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
779be5899b3SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
780be5899b3SLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
781f68a32c8SEmil Constantinescu     }
782f68a32c8SEmil Constantinescu   }
783f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
784e4dd521cSBarry Smith }
785e4dd521cSBarry Smith 
786f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts)
787f6a906c0SBarry Smith {
788f6a906c0SBarry Smith   TS_RK          *rk  = (TS_RK*)ts->data;
789be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
790be5899b3SLisandro Dalcin   PetscInt       s   = tab->s;
791f6a906c0SBarry Smith   PetscErrorCode ierr;
792f6a906c0SBarry Smith 
793f6a906c0SBarry Smith   PetscFunctionBegin;
794f6a906c0SBarry Smith   if (ts->adjointsetupcalled++) PetscFunctionReturn(0);
7952e7b7f96SHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecsDeltaLam);CHKERRQ(ierr);
7962e7b7f96SHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&rk->VecsSensiTemp);CHKERRQ(ierr);
797f6a906c0SBarry Smith   if(ts->vecs_sensip) {
7982e7b7f96SHong Zhang     ierr = VecDuplicate(ts->vecs_sensip[0],&rk->VecDeltaMu);CHKERRQ(ierr);
799f6a906c0SBarry Smith   }
80013af1a74SHong Zhang   if (ts->vecs_sensi2) {
80113af1a74SHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecsDeltaLam2);CHKERRQ(ierr);
80213af1a74SHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&rk->VecsSensi2Temp);CHKERRQ(ierr);
80313af1a74SHong Zhang   }
80413af1a74SHong Zhang   if (ts->vecs_sensi2p) {
80513af1a74SHong Zhang     ierr = VecDuplicate(ts->vecs_sensi2p[0],&rk->VecDeltaMu2);CHKERRQ(ierr);
80613af1a74SHong Zhang   }
807f6a906c0SBarry Smith   PetscFunctionReturn(0);
808f6a906c0SBarry Smith }
809f6a906c0SBarry Smith 
81035f5172aSHong Zhang /*
811*3ca0882eSHong Zhang   Use SNES to compute the Jacobian so that finite differencing could be used when TS Jacobian is not available.
812*3ca0882eSHong Zhang */
813*3ca0882eSHong Zhang static PetscErrorCode TSFormRHSJacobian_Private(Vec x,Mat J,Mat Jpre,TS ts)
814*3ca0882eSHong Zhang {
815*3ca0882eSHong Zhang   SNES           snes = ts->snes;
816*3ca0882eSHong Zhang   PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*) = NULL;
817*3ca0882eSHong Zhang   PetscErrorCode ierr;
818*3ca0882eSHong Zhang 
819*3ca0882eSHong Zhang   PetscFunctionBegin;
820*3ca0882eSHong Zhang   /*
821*3ca0882eSHong Zhang     Unlike implicit methods, explicit methods do not have SNESMatFDColoring in the snes object
822*3ca0882eSHong Zhang     because SNESSolve() has not been called yet. Insteading of querying it, we check the
823*3ca0882eSHong Zhang     Jacobian compute function directly to determin if FD coloring is used.
824*3ca0882eSHong Zhang   */
825*3ca0882eSHong Zhang   ierr = SNESGetJacobian(snes,NULL,NULL,&jac,NULL);CHKERRQ(ierr);
826*3ca0882eSHong Zhang   if (jac == SNESComputeJacobianDefaultColor) {
827*3ca0882eSHong Zhang     Vec f;
828*3ca0882eSHong Zhang     ierr = SNESSetSolution(snes,x);CHKERRQ(ierr);
829*3ca0882eSHong Zhang     ierr = SNESGetFunction(snes,&f,NULL,NULL);CHKERRQ(ierr);
830*3ca0882eSHong Zhang     /* Force MatFDColoringApply to evaluate the SNES residual function for the base vector */
831*3ca0882eSHong Zhang     ierr = SNESComputeFunction(snes,x,f);CHKERRQ(ierr);
832*3ca0882eSHong Zhang   }
833*3ca0882eSHong Zhang   ierr = SNESComputeJacobian(snes,x,J,Jpre);CHKERRQ(ierr);
834*3ca0882eSHong Zhang   PetscFunctionReturn(0);
835*3ca0882eSHong Zhang }
836*3ca0882eSHong Zhang 
837*3ca0882eSHong Zhang /*
83835f5172aSHong Zhang   Assumptions:
83935f5172aSHong Zhang     - TSStep_RK() always evaluates the step with b, not bembed.
84035f5172aSHong Zhang */
84142f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts)
842d2daff3dSHong Zhang {
843c235aa8dSHong Zhang   TS_RK            *rk = (TS_RK*)ts->data;
844cd4cee2dSHong Zhang   TS               quadts = ts->quadraturets;
845c235aa8dSHong Zhang   RKTableau        tab = rk->tableau;
846*3ca0882eSHong Zhang   Mat              J,Jpre,Jquad;
847c235aa8dSHong Zhang   const PetscInt   s = tab->s;
848c235aa8dSHong Zhang   const PetscReal  *A = tab->A,*b = tab->b,*c = tab->c;
84913af1a74SHong Zhang   PetscScalar      *w = rk->work,*xarr;
8502e7b7f96SHong Zhang   Vec              *Y = rk->Y,*VecsDeltaLam = rk->VecsDeltaLam,VecDeltaMu = rk->VecDeltaMu,*VecsSensiTemp = rk->VecsSensiTemp;
85113af1a74SHong Zhang   Vec              *VecsDeltaLam2 = rk->VecsDeltaLam2,VecDeltaMu2 = rk->VecDeltaMu2,*VecsSensi2Temp = rk->VecsSensi2Temp;
852cd4cee2dSHong Zhang   Vec              VecDRDUTransCol = ts->vec_drdu_col,VecDRDPTransCol = ts->vec_drdp_col;
853b18ea86cSHong Zhang   PetscInt         i,j,nadj;
8543d3eaba7SBarry Smith   PetscReal        t = ts->ptime;
855*3ca0882eSHong Zhang   PetscReal        h = ts->time_step;
856d2daff3dSHong Zhang   PetscErrorCode   ierr;
857c235aa8dSHong Zhang 
858d2daff3dSHong Zhang   PetscFunctionBegin;
859c235aa8dSHong Zhang   rk->status = TS_STEP_INCOMPLETE;
8603389c7d9SStefano Zampini 
861*3ca0882eSHong Zhang   ierr = TSGetRHSJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr);
86235f5172aSHong Zhang   if (quadts) {
863cd4cee2dSHong Zhang     ierr = TSGetRHSJacobian(quadts,&Jquad,NULL,NULL,NULL);CHKERRQ(ierr);
86435f5172aSHong Zhang   }
8652e7b7f96SHong Zhang   for (i=s-1; i>=0; i--) {
8666a1e4597SHong Zhang     if (tab->FSAL && i == s-1) {
8676a1e4597SHong Zhang       /* VecsDeltaLam[nadj*s+s-1] are initialized with zeros and the values never change.*/
8686a1e4597SHong Zhang       continue;
8696a1e4597SHong Zhang     }
870*3ca0882eSHong Zhang     rk->stage_time = t + h*(1.0-c[i]);
871*3ca0882eSHong Zhang     ierr = TSFormRHSJacobian_Private(Y[i],J,Jpre,ts);CHKERRQ(ierr);
87235f5172aSHong Zhang     if (quadts) {
873*3ca0882eSHong Zhang       ierr = TSComputeRHSJacobian(quadts,rk->stage_time,Y[i],Jquad,Jquad);CHKERRQ(ierr); /* get r_u^T */
87435f5172aSHong Zhang     }
8753389c7d9SStefano Zampini     if (ts->vecs_sensip) {
876*3ca0882eSHong Zhang       ierr = TSComputeRHSJacobianP(ts,rk->stage_time,Y[i],ts->Jacprhs);CHKERRQ(ierr); /* get f_p */
87735f5172aSHong Zhang       if (quadts) {
878*3ca0882eSHong Zhang         ierr = TSComputeRHSJacobianP(quadts,rk->stage_time,Y[i],quadts->Jacprhs);CHKERRQ(ierr); /* get f_p for the quadrature */
8793389c7d9SStefano Zampini       }
88035f5172aSHong Zhang     }
8813389c7d9SStefano Zampini 
88235f5172aSHong Zhang     if (b[i]) {
88335f5172aSHong Zhang       for (j=i+1; j<s; j++) w[j-i-1] = A[j*s+i]/b[i]; /* coefficients for computing VecsSensiTemp */
88435f5172aSHong Zhang     } else {
88535f5172aSHong Zhang       for (j=i+1; j<s; j++) w[j-i-1] = A[j*s+i]; /* coefficients for computing VecsSensiTemp */
88635f5172aSHong Zhang     }
8872e7b7f96SHong Zhang 
8882e7b7f96SHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
8893389c7d9SStefano Zampini       /* Stage values of lambda */
89035f5172aSHong Zhang       if (b[i]) {
89135f5172aSHong Zhang         /* lambda_{n+1} + \sum_{j=i+1}^s a_{ji}/b[i]*lambda_{s,j} */
8922e7b7f96SHong Zhang         ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); /* VecDeltaLam is an vec array of size s by numcost */
8932e7b7f96SHong Zhang         ierr = VecMAXPY(VecsSensiTemp[nadj],s-i-1,w,&VecsDeltaLam[nadj*s+i+1]);CHKERRQ(ierr);
89435f5172aSHong Zhang         ierr = MatMultTranspose(J,VecsSensiTemp[nadj],VecsDeltaLam[nadj*s+i]);CHKERRQ(ierr); /* VecsSensiTemp will be reused by 2nd-order adjoint */
89535f5172aSHong Zhang         ierr = VecScale(VecsDeltaLam[nadj*s+i],-h*b[i]);CHKERRQ(ierr);
89635f5172aSHong Zhang         if (quadts) {
897cd4cee2dSHong Zhang           ierr = MatDenseGetColumn(Jquad,nadj,&xarr);CHKERRQ(ierr);
898cd4cee2dSHong Zhang           ierr = VecPlaceArray(VecDRDUTransCol,xarr);CHKERRQ(ierr);
899cd4cee2dSHong Zhang           ierr = VecAXPY(VecsDeltaLam[nadj*s+i],-h*b[i],VecDRDUTransCol);CHKERRQ(ierr);
900cd4cee2dSHong Zhang           ierr = VecResetArray(VecDRDUTransCol);CHKERRQ(ierr);
901cd4cee2dSHong Zhang           ierr = MatDenseRestoreColumn(Jquad,&xarr);CHKERRQ(ierr);
902cd4cee2dSHong Zhang         }
9033389c7d9SStefano Zampini       } else {
9046a1e4597SHong Zhang         /* \sum_{j=i+1}^s a_{ji}*lambda_{s,j} */
9056a1e4597SHong Zhang         ierr = VecSet(VecsSensiTemp[nadj],0);CHKERRQ(ierr);
9066a1e4597SHong Zhang         ierr = VecMAXPY(VecsSensiTemp[nadj],s-i-1,w,&VecsDeltaLam[nadj*s+i+1]);CHKERRQ(ierr);
9076a1e4597SHong Zhang         ierr = MatMultTranspose(J,VecsSensiTemp[nadj],VecsDeltaLam[nadj*s+i]);CHKERRQ(ierr);
9086a1e4597SHong Zhang         ierr = VecScale(VecsDeltaLam[nadj*s+i],-h);CHKERRQ(ierr);
9093389c7d9SStefano Zampini       }
9106497ce14SHong Zhang 
911ad8e2604SHong Zhang       /* Stage values of mu */
9126497ce14SHong Zhang       if (ts->vecs_sensip) {
91313af1a74SHong Zhang         ierr = MatMultTranspose(ts->Jacprhs,VecsSensiTemp[nadj],VecDeltaMu);CHKERRQ(ierr);
91435f5172aSHong Zhang         if (b[i]) {
91535f5172aSHong Zhang           ierr = VecScale(VecDeltaMu,-h*b[i]);CHKERRQ(ierr);
91635f5172aSHong Zhang           if (quadts) {
917cd4cee2dSHong Zhang             ierr = MatDenseGetColumn(quadts->Jacprhs,nadj,&xarr);CHKERRQ(ierr);
918cd4cee2dSHong Zhang             ierr = VecPlaceArray(VecDRDPTransCol,xarr);CHKERRQ(ierr);
919cd4cee2dSHong Zhang             ierr = VecAXPY(VecDeltaMu,-h*b[i],VecDRDPTransCol);CHKERRQ(ierr);
920cd4cee2dSHong Zhang             ierr = VecResetArray(VecDRDPTransCol);CHKERRQ(ierr);
921cd4cee2dSHong Zhang             ierr = MatDenseRestoreColumn(quadts->Jacprhs,&xarr);CHKERRQ(ierr);
922493ed95fSHong Zhang           }
92335f5172aSHong Zhang         } else {
92435f5172aSHong Zhang           ierr = VecScale(VecDeltaMu,-h);CHKERRQ(ierr);
92535f5172aSHong Zhang         }
9262e7b7f96SHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],1.,VecDeltaMu);CHKERRQ(ierr); /* update sensip for each stage */
927ad8e2604SHong Zhang       }
928c235aa8dSHong Zhang     }
92913af1a74SHong Zhang 
93013af1a74SHong Zhang     if (ts->vecs_sensi2 && ts->forward_solve) { /* 2nd-order adjoint, TLM mode has to be turned on */
93113af1a74SHong Zhang       /* Get w1 at t_{n+1} from TLM matrix */
93213af1a74SHong Zhang       ierr = MatDenseGetColumn(rk->MatsFwdStageSensip[i],0,&xarr);CHKERRQ(ierr);
93313af1a74SHong Zhang       ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
93413af1a74SHong Zhang       /* lambda_s^T F_UU w_1 */
935*3ca0882eSHong Zhang       ierr = TSComputeRHSHessianProductFunctionUU(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_sensip_col,ts->vecs_guu);CHKERRQ(ierr);
93635f5172aSHong Zhang       if (quadts)  {
93735f5172aSHong Zhang         /* R_UU w_1 */
938*3ca0882eSHong Zhang         ierr = TSComputeRHSHessianProductFunctionUU(quadts,rk->stage_time,Y[i],NULL,ts->vec_sensip_col,ts->vecs_guu);CHKERRQ(ierr);
93935f5172aSHong Zhang       }
94035f5172aSHong Zhang       if (ts->vecs_sensip) {
94113af1a74SHong Zhang         /* lambda_s^T F_UP w_2 */
942*3ca0882eSHong Zhang         ierr = TSComputeRHSHessianProductFunctionUP(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_dir,ts->vecs_gup);CHKERRQ(ierr);
94335f5172aSHong Zhang         if (quadts)  {
94435f5172aSHong Zhang           /* R_UP w_2 */
945*3ca0882eSHong Zhang           ierr = TSComputeRHSHessianProductFunctionUP(quadts,rk->stage_time,Y[i],NULL,ts->vec_sensip_col,ts->vecs_gup);CHKERRQ(ierr);
94635f5172aSHong Zhang         }
94735f5172aSHong Zhang       }
94813af1a74SHong Zhang       if (ts->vecs_sensi2p) {
94913af1a74SHong Zhang         /* lambda_s^T F_PU w_1 */
950*3ca0882eSHong Zhang         ierr = TSComputeRHSHessianProductFunctionPU(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_sensip_col,ts->vecs_gpu);CHKERRQ(ierr);
95135f5172aSHong Zhang         /* lambda_s^T F_PP w_2 */
952*3ca0882eSHong Zhang         ierr = TSComputeRHSHessianProductFunctionPP(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_dir,ts->vecs_gpp);CHKERRQ(ierr);
95335f5172aSHong Zhang         if (b[i] && quadts) {
95435f5172aSHong Zhang           /* R_PU w_1 */
955*3ca0882eSHong Zhang           ierr = TSComputeRHSHessianProductFunctionPU(quadts,rk->stage_time,Y[i],NULL,ts->vec_sensip_col,ts->vecs_gpu);CHKERRQ(ierr);
95635f5172aSHong Zhang           /* R_PP w_2 */
957*3ca0882eSHong Zhang           ierr = TSComputeRHSHessianProductFunctionPP(quadts,rk->stage_time,Y[i],NULL,ts->vec_dir,ts->vecs_gpp);CHKERRQ(ierr);
95835f5172aSHong Zhang         }
95913af1a74SHong Zhang       }
96013af1a74SHong Zhang       ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
96113af1a74SHong Zhang       ierr = MatDenseRestoreColumn(rk->MatsFwdStageSensip[i],&xarr);CHKERRQ(ierr);
96213af1a74SHong Zhang 
96313af1a74SHong Zhang       for (nadj=0; nadj<ts->numcost; nadj++) {
96413af1a74SHong Zhang         /* Stage values of lambda */
96535f5172aSHong Zhang         if (b[i]) {
96635f5172aSHong Zhang           /* J_i^T*(Lambda_{n+1}+\sum_{j=i+1}^s a_{ji}/b_i*Lambda_{s,j} */
96713af1a74SHong Zhang           ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
96813af1a74SHong Zhang           ierr = VecMAXPY(VecsSensi2Temp[nadj],s-i-1,w,&VecsDeltaLam2[nadj*s+i+1]);CHKERRQ(ierr);
96913af1a74SHong Zhang           ierr = MatMultTranspose(J,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj*s+i]);CHKERRQ(ierr);
97035f5172aSHong Zhang           ierr = VecScale(VecsDeltaLam2[nadj*s+i],-h*b[i]);CHKERRQ(ierr);
97135f5172aSHong Zhang           ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],-h*b[i],ts->vecs_guu[nadj]);CHKERRQ(ierr);
97235f5172aSHong Zhang           if (ts->vecs_sensip) {
97335f5172aSHong Zhang             ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],-h*b[i],ts->vecs_gup[nadj]);CHKERRQ(ierr);
97413af1a74SHong Zhang           }
97513af1a74SHong Zhang         } else {
97635f5172aSHong Zhang           /* \sum_{j=i+1}^s a_{ji}*Lambda_{s,j} */
97713af1a74SHong Zhang           ierr = VecSet(VecsDeltaLam2[nadj*s+i],0);CHKERRQ(ierr);
97835f5172aSHong Zhang           ierr = VecMAXPY(VecsSensi2Temp[nadj],s-i-1,w,&VecsDeltaLam2[nadj*s+i+1]);CHKERRQ(ierr);
97935f5172aSHong Zhang           ierr = MatMultTranspose(J,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj*s+i]);CHKERRQ(ierr);
98035f5172aSHong Zhang           ierr = VecScale(VecsDeltaLam2[nadj*s+i],-h);CHKERRQ(ierr);
98135f5172aSHong Zhang           ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],-h,ts->vecs_guu[nadj]);CHKERRQ(ierr);
98235f5172aSHong Zhang           if (ts->vecs_sensip) {
98335f5172aSHong Zhang             ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],-h,ts->vecs_gup[nadj]);CHKERRQ(ierr);
98413af1a74SHong Zhang           }
98535f5172aSHong Zhang         }
98635f5172aSHong Zhang         if (ts->vecs_sensi2p) { /* 2nd-order adjoint for parameters */
98713af1a74SHong Zhang           ierr = MatMultTranspose(ts->Jacprhs,VecsSensi2Temp[nadj],VecDeltaMu2);CHKERRQ(ierr);
98835f5172aSHong Zhang           if (b[i]) {
98935f5172aSHong Zhang             ierr = VecScale(VecDeltaMu2,-h*b[i]);CHKERRQ(ierr);
99035f5172aSHong Zhang             ierr = VecAXPY(VecDeltaMu2,-h*b[i],ts->vecs_gpu[nadj]);CHKERRQ(ierr);
99135f5172aSHong Zhang             ierr = VecAXPY(VecDeltaMu2,-h*b[i],ts->vecs_gpp[nadj]);CHKERRQ(ierr);
99213af1a74SHong Zhang           } else {
99335f5172aSHong Zhang             ierr = VecScale(VecDeltaMu2,-h);CHKERRQ(ierr);
99435f5172aSHong Zhang             ierr = VecAXPY(VecDeltaMu2,-h,ts->vecs_gpu[nadj]);CHKERRQ(ierr);
99535f5172aSHong Zhang             ierr = VecAXPY(VecDeltaMu2,-h,ts->vecs_gpp[nadj]);CHKERRQ(ierr);
99613af1a74SHong Zhang           }
99713af1a74SHong Zhang           ierr = VecAXPY(ts->vecs_sensi2p[nadj],1,VecDeltaMu2);CHKERRQ(ierr); /* update sensi2p for each stage */
99813af1a74SHong Zhang         }
99913af1a74SHong Zhang       }
100013af1a74SHong Zhang     }
10016497ce14SHong Zhang   }
1002c235aa8dSHong Zhang 
1003c235aa8dSHong Zhang   for (j=0; j<s; j++) w[j] = 1.0;
10042e7b7f96SHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) { /* no need to do this for mu's */
10052e7b7f96SHong Zhang     ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecsDeltaLam[nadj*s]);CHKERRQ(ierr);
100613af1a74SHong Zhang     if (ts->vecs_sensi2) {
100713af1a74SHong Zhang       ierr = VecMAXPY(ts->vecs_sensi2[nadj],s,w,&VecsDeltaLam2[nadj*s]);CHKERRQ(ierr);
100813af1a74SHong Zhang     }
10096497ce14SHong Zhang   }
1010c235aa8dSHong Zhang   rk->status = TS_STEP_COMPLETE;
1011d2daff3dSHong Zhang   PetscFunctionReturn(0);
1012d2daff3dSHong Zhang }
1013d2daff3dSHong Zhang 
101413af1a74SHong Zhang static PetscErrorCode TSAdjointReset_RK(TS ts)
101513af1a74SHong Zhang {
101613af1a74SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
101713af1a74SHong Zhang   RKTableau      tab = rk->tableau;
101813af1a74SHong Zhang   PetscErrorCode ierr;
101913af1a74SHong Zhang 
102013af1a74SHong Zhang   PetscFunctionBegin;
102113af1a74SHong Zhang   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam);CHKERRQ(ierr);
102213af1a74SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&rk->VecsSensiTemp);CHKERRQ(ierr);
102313af1a74SHong Zhang   ierr = VecDestroy(&rk->VecDeltaMu);CHKERRQ(ierr);
102413af1a74SHong Zhang   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam2);CHKERRQ(ierr);
102513af1a74SHong Zhang   ierr = VecDestroy(&rk->VecDeltaMu2);CHKERRQ(ierr);
102613af1a74SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&rk->VecsSensi2Temp);CHKERRQ(ierr);
102713af1a74SHong Zhang   PetscFunctionReturn(0);
102813af1a74SHong Zhang }
102913af1a74SHong Zhang 
103055de54fcSHong Zhang static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X)
103155de54fcSHong Zhang {
103255de54fcSHong Zhang   TS_RK            *rk = (TS_RK*)ts->data;
103355de54fcSHong Zhang   PetscInt         s  = rk->tableau->s,p = rk->tableau->p,i,j;
103455de54fcSHong Zhang   PetscReal        h;
103555de54fcSHong Zhang   PetscReal        tt,t;
103655de54fcSHong Zhang   PetscScalar      *b;
103755de54fcSHong Zhang   const PetscReal  *B = rk->tableau->binterp;
103855de54fcSHong Zhang   PetscErrorCode   ierr;
103955de54fcSHong Zhang 
104055de54fcSHong Zhang   PetscFunctionBegin;
104155de54fcSHong Zhang   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
104255de54fcSHong Zhang 
104355de54fcSHong Zhang   switch (rk->status) {
104455de54fcSHong Zhang     case TS_STEP_INCOMPLETE:
104555de54fcSHong Zhang     case TS_STEP_PENDING:
104655de54fcSHong Zhang       h = ts->time_step;
104755de54fcSHong Zhang       t = (itime - ts->ptime)/h;
104855de54fcSHong Zhang       break;
104955de54fcSHong Zhang     case TS_STEP_COMPLETE:
105055de54fcSHong Zhang       h = ts->ptime - ts->ptime_prev;
105155de54fcSHong Zhang       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
105255de54fcSHong Zhang       break;
105355de54fcSHong Zhang     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
105455de54fcSHong Zhang   }
105555de54fcSHong Zhang   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
105655de54fcSHong Zhang   for (i=0; i<s; i++) b[i] = 0;
105755de54fcSHong Zhang   for (j=0,tt=t; j<p; j++,tt*=t) {
105855de54fcSHong Zhang     for (i=0; i<s; i++) {
105955de54fcSHong Zhang       b[i]  += h * B[i*p+j] * tt;
106055de54fcSHong Zhang     }
106155de54fcSHong Zhang   }
106255de54fcSHong Zhang   ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr);
106355de54fcSHong Zhang   ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr);
106455de54fcSHong Zhang   ierr = PetscFree(b);CHKERRQ(ierr);
106555de54fcSHong Zhang   PetscFunctionReturn(0);
106655de54fcSHong Zhang }
106755de54fcSHong Zhang 
106855de54fcSHong Zhang /*------------------------------------------------------------*/
106955de54fcSHong Zhang 
1070be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts)
1071be5899b3SLisandro Dalcin {
1072be5899b3SLisandro Dalcin   TS_RK          *rk = (TS_RK*)ts->data;
1073be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
1074be5899b3SLisandro Dalcin   PetscErrorCode ierr;
1075be5899b3SLisandro Dalcin 
1076be5899b3SLisandro Dalcin   PetscFunctionBegin;
1077be5899b3SLisandro Dalcin   if (!tab) PetscFunctionReturn(0);
1078be5899b3SLisandro Dalcin   ierr = PetscFree(rk->work);CHKERRQ(ierr);
1079be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr);
1080be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr);
10812e7b7f96SHong Zhang   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam);CHKERRQ(ierr);
10822e7b7f96SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&rk->VecsSensiTemp);CHKERRQ(ierr);
10832e7b7f96SHong Zhang   ierr = VecDestroy(&rk->VecDeltaMu);CHKERRQ(ierr);
1084be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
1085be5899b3SLisandro Dalcin }
1086be5899b3SLisandro Dalcin 
1087277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts)
1088e4dd521cSBarry Smith {
10896849ba73SBarry Smith   PetscErrorCode ierr;
1090e4dd521cSBarry Smith 
1091e4dd521cSBarry Smith   PetscFunctionBegin;
1092be5899b3SLisandro Dalcin   ierr = TSRKTableauReset(ts);CHKERRQ(ierr);
10930fe4d17eSHong Zhang   if (ts->use_splitrhsfunction) {
109402555c94SHong Zhang     ierr = PetscTryMethod(ts,"TSReset_RK_MultirateSplit_C",(TS),(ts));CHKERRQ(ierr);
10950fe4d17eSHong Zhang   } else {
109602555c94SHong Zhang     ierr = PetscTryMethod(ts,"TSReset_RK_MultirateNonsplit_C",(TS),(ts));CHKERRQ(ierr);
10970fe4d17eSHong Zhang   }
1098277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
1099e4dd521cSBarry Smith }
1100277b19d0SLisandro Dalcin 
1101f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx)
1102f68a32c8SEmil Constantinescu {
1103f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1104f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1105f68a32c8SEmil Constantinescu }
1106f68a32c8SEmil Constantinescu 
1107f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
1108f68a32c8SEmil Constantinescu {
1109f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1110f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1111f68a32c8SEmil Constantinescu }
1112f68a32c8SEmil Constantinescu 
1113f68a32c8SEmil Constantinescu 
1114f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx)
1115f68a32c8SEmil Constantinescu {
1116f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1117f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1118f68a32c8SEmil Constantinescu }
1119f68a32c8SEmil Constantinescu 
1120f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
1121f68a32c8SEmil Constantinescu {
1122f68a32c8SEmil Constantinescu 
1123f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1124f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1125f68a32c8SEmil Constantinescu }
1126be5899b3SLisandro Dalcin 
1127be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts)
1128be5899b3SLisandro Dalcin {
1129be5899b3SLisandro Dalcin   TS_RK          *rk  = (TS_RK*)ts->data;
1130be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
1131be5899b3SLisandro Dalcin   PetscErrorCode ierr;
1132be5899b3SLisandro Dalcin 
1133be5899b3SLisandro Dalcin   PetscFunctionBegin;
1134be5899b3SLisandro Dalcin   ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr);
1135be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr);
1136be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr);
1137be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
1138be5899b3SLisandro Dalcin }
1139be5899b3SLisandro Dalcin 
1140f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts)
1141f68a32c8SEmil Constantinescu {
1142cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
1143f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
1144f68a32c8SEmil Constantinescu   DM             dm;
1145f68a32c8SEmil Constantinescu 
1146f68a32c8SEmil Constantinescu   PetscFunctionBegin;
11473dd83b38SBarry Smith   ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr);
1148be5899b3SLisandro Dalcin   ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);
1149cd4cee2dSHong Zhang   if (quadts && ts->costintegralfwd) {
1150cd4cee2dSHong Zhang     Mat Jquad;
1151cd4cee2dSHong Zhang     ierr = TSGetRHSJacobian(quadts,&Jquad,NULL,NULL,NULL);CHKERRQ(ierr);
11520f7a1166SHong Zhang   }
1153f68a32c8SEmil Constantinescu   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1154f68a32c8SEmil Constantinescu   ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
1155f68a32c8SEmil Constantinescu   ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
11560fe4d17eSHong Zhang   if (ts->use_splitrhsfunction) {
115702555c94SHong Zhang     ierr = PetscTryMethod(ts,"TSSetUp_RK_MultirateSplit_C",(TS),(ts));CHKERRQ(ierr);
11580fe4d17eSHong Zhang   } else {
115902555c94SHong Zhang     ierr = PetscTryMethod(ts,"TSSetUp_RK_MultirateNonsplit_C",(TS),(ts));CHKERRQ(ierr);
11600fe4d17eSHong Zhang   }
1161e4dd521cSBarry Smith   PetscFunctionReturn(0);
1162e4dd521cSBarry Smith }
1163d2daff3dSHong Zhang 
11644416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts)
1165e4dd521cSBarry Smith {
1166be5899b3SLisandro Dalcin   TS_RK          *rk = (TS_RK*)ts->data;
1167dfbe8321SBarry Smith   PetscErrorCode ierr;
1168262ff7bbSBarry Smith 
1169e4dd521cSBarry Smith   PetscFunctionBegin;
1170e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr);
1171f68a32c8SEmil Constantinescu   {
1172f68a32c8SEmil Constantinescu     RKTableauLink link;
1173f68a32c8SEmil Constantinescu     PetscInt      count,choice;
11740fe4d17eSHong Zhang     PetscBool     flg,use_multirate = PETSC_FALSE;
1175f68a32c8SEmil Constantinescu     const char    **namelist;
11762c9cb062Sluzhanghpp 
1177f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) ;
1178f489ac74SBarry Smith     ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr);
1179f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
11800fe4d17eSHong Zhang     ierr = PetscOptionsBool("-ts_rk_multirate","Use interpolation-based multirate RK method","TSRKSetMultirate",rk->use_multirate,&use_multirate,&flg);CHKERRQ(ierr);
11810fe4d17eSHong Zhang     if (flg) {
11820fe4d17eSHong Zhang       ierr = TSRKSetMultirate(ts,use_multirate);CHKERRQ(ierr);
11830fe4d17eSHong Zhang     }
1184be5899b3SLisandro Dalcin     ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr);
1185be5899b3SLisandro Dalcin     if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);}
1186f68a32c8SEmil Constantinescu     ierr = PetscFree(namelist);CHKERRQ(ierr);
1187f68a32c8SEmil Constantinescu   }
1188262ff7bbSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
11892c9cb062Sluzhanghpp   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Multirate methods options","");CHKERRQ(ierr);
11902c9cb062Sluzhanghpp   ierr = PetscOptionsInt("-ts_rk_dtratio","time step ratio between slow and fast","",rk->dtratio,&rk->dtratio,NULL);CHKERRQ(ierr);
11912c9cb062Sluzhanghpp   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1192e4dd521cSBarry Smith   PetscFunctionReturn(0);
1193e4dd521cSBarry Smith }
1194e4dd521cSBarry Smith 
11955f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
1196e4dd521cSBarry Smith {
11975f70b478SJed Brown   TS_RK          *rk = (TS_RK*)ts->data;
11988370ee3bSLisandro Dalcin   PetscBool      iascii;
1199dfbe8321SBarry Smith   PetscErrorCode ierr;
12002cdf8120Sasbjorn 
1201e4dd521cSBarry Smith   PetscFunctionBegin;
1202251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
12038370ee3bSLisandro Dalcin   if (iascii) {
12049c334d8fSLisandro Dalcin     RKTableau tab  = rk->tableau;
1205f68a32c8SEmil Constantinescu     TSRKType  rktype;
1206f68a32c8SEmil Constantinescu     char      buf[512];
1207f68a32c8SEmil Constantinescu     ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr);
1208efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  RK type %s\n",rktype);CHKERRQ(ierr);
12090c37eb8eSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  Order: %D\n",tab->order);CHKERRQ(ierr);
12100c37eb8eSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr);
1211f68a32c8SEmil Constantinescu     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
1212f68a32c8SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa c = %s\n",buf);CHKERRQ(ierr);
12138370ee3bSLisandro Dalcin   }
1214f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1215f68a32c8SEmil Constantinescu }
1216f68a32c8SEmil Constantinescu 
1217f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer)
1218f68a32c8SEmil Constantinescu {
1219f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
12209c334d8fSLisandro Dalcin   TSAdapt        adapt;
1221f68a32c8SEmil Constantinescu 
1222f68a32c8SEmil Constantinescu   PetscFunctionBegin;
12239c334d8fSLisandro Dalcin   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
12249c334d8fSLisandro Dalcin   ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr);
1225f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1226f68a32c8SEmil Constantinescu }
1227f68a32c8SEmil Constantinescu 
1228f68a32c8SEmil Constantinescu /*@C
1229f68a32c8SEmil Constantinescu   TSRKSetType - Set the type of RK scheme
1230f68a32c8SEmil Constantinescu 
1231f68a32c8SEmil Constantinescu   Logically collective
1232f68a32c8SEmil Constantinescu 
1233f68a32c8SEmil Constantinescu   Input Parameter:
1234f68a32c8SEmil Constantinescu +  ts - timestepping context
1235f68a32c8SEmil Constantinescu -  rktype - type of RK-scheme
1236f68a32c8SEmil Constantinescu 
1237287c2655SBarry Smith   Options Database:
12389bd3e852SBarry Smith .   -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs>
1239287c2655SBarry Smith 
1240f68a32c8SEmil Constantinescu   Level: intermediate
1241f68a32c8SEmil Constantinescu 
124237acaa02SHendrik Ranocha .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS, TSRK6VR, TSRK7VR, TSRK8VR
1243f68a32c8SEmil Constantinescu @*/
1244f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype)
1245f68a32c8SEmil Constantinescu {
1246f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
1247f68a32c8SEmil Constantinescu 
1248f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1249f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1250b92453a8SLisandro Dalcin   PetscValidCharPointer(rktype,2);
1251f68a32c8SEmil Constantinescu   ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr);
1252f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1253f68a32c8SEmil Constantinescu }
1254f68a32c8SEmil Constantinescu 
1255f68a32c8SEmil Constantinescu /*@C
1256f68a32c8SEmil Constantinescu   TSRKGetType - Get the type of RK scheme
1257f68a32c8SEmil Constantinescu 
1258f68a32c8SEmil Constantinescu   Logically collective
1259f68a32c8SEmil Constantinescu 
1260f68a32c8SEmil Constantinescu   Input Parameter:
1261f68a32c8SEmil Constantinescu .  ts - timestepping context
1262f68a32c8SEmil Constantinescu 
1263f68a32c8SEmil Constantinescu   Output Parameter:
1264f68a32c8SEmil Constantinescu .  rktype - type of RK-scheme
1265f68a32c8SEmil Constantinescu 
1266f68a32c8SEmil Constantinescu   Level: intermediate
1267f68a32c8SEmil Constantinescu 
1268f68a32c8SEmil Constantinescu .seealso: TSRKGetType()
1269f68a32c8SEmil Constantinescu @*/
1270f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype)
1271f68a32c8SEmil Constantinescu {
1272f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
1273f68a32c8SEmil Constantinescu 
1274f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1275f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1276f68a32c8SEmil Constantinescu   ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr);
1277f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1278f68a32c8SEmil Constantinescu }
1279f68a32c8SEmil Constantinescu 
1280560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype)
1281f68a32c8SEmil Constantinescu {
1282f68a32c8SEmil Constantinescu   TS_RK *rk = (TS_RK*)ts->data;
1283f68a32c8SEmil Constantinescu 
1284f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1285f68a32c8SEmil Constantinescu   *rktype = rk->tableau->name;
1286f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1287f68a32c8SEmil Constantinescu }
12882c9cb062Sluzhanghpp 
1289560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype)
1290f68a32c8SEmil Constantinescu {
1291f68a32c8SEmil Constantinescu   TS_RK          *rk = (TS_RK*)ts->data;
1292f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
1293f68a32c8SEmil Constantinescu   PetscBool      match;
1294f68a32c8SEmil Constantinescu   RKTableauLink  link;
1295f68a32c8SEmil Constantinescu 
1296f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1297f68a32c8SEmil Constantinescu   if (rk->tableau) {
1298f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr);
1299f68a32c8SEmil Constantinescu     if (match) PetscFunctionReturn(0);
1300f68a32c8SEmil Constantinescu   }
1301f68a32c8SEmil Constantinescu   for (link = RKTableauList; link; link=link->next) {
1302f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr);
1303f68a32c8SEmil Constantinescu     if (match) {
1304be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);}
1305f68a32c8SEmil Constantinescu       rk->tableau = &link->tab;
1306be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);}
13072ffb9264SLisandro Dalcin       ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
1308f68a32c8SEmil Constantinescu       PetscFunctionReturn(0);
1309f68a32c8SEmil Constantinescu     }
1310f68a32c8SEmil Constantinescu   }
1311f68a32c8SEmil Constantinescu   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype);
1312e4dd521cSBarry Smith   PetscFunctionReturn(0);
1313e4dd521cSBarry Smith }
1314e4dd521cSBarry Smith 
1315ff22ae23SHong Zhang static PetscErrorCode  TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y)
1316ff22ae23SHong Zhang {
1317ff22ae23SHong Zhang   TS_RK *rk = (TS_RK*)ts->data;
1318ff22ae23SHong Zhang 
1319ff22ae23SHong Zhang   PetscFunctionBegin;
13200429704eSStefano Zampini   if (ns) *ns = rk->tableau->s;
1321d2daff3dSHong Zhang   if (Y)   *Y = rk->Y;
1322ff22ae23SHong Zhang   PetscFunctionReturn(0);
1323ff22ae23SHong Zhang }
1324ff22ae23SHong Zhang 
1325b3a6b972SJed Brown static PetscErrorCode TSDestroy_RK(TS ts)
1326b3a6b972SJed Brown {
1327b3a6b972SJed Brown   PetscErrorCode ierr;
1328b3a6b972SJed Brown 
1329b3a6b972SJed Brown   PetscFunctionBegin;
1330b3a6b972SJed Brown   ierr = TSReset_RK(ts);CHKERRQ(ierr);
1331b3a6b972SJed Brown   if (ts->dm) {
1332b3a6b972SJed Brown     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
1333b3a6b972SJed Brown     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
1334b3a6b972SJed Brown   }
1335b3a6b972SJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
1336b3a6b972SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr);
1337b3a6b972SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr);
13380fe4d17eSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",NULL);CHKERRQ(ierr);
13390fe4d17eSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",NULL);CHKERRQ(ierr);
1340b3a6b972SJed Brown   PetscFunctionReturn(0);
1341b3a6b972SJed Brown }
1342ff22ae23SHong Zhang 
1343*3ca0882eSHong Zhang /*
1344*3ca0882eSHong Zhang   This defines the nonlinear equation that is to be solved with SNES
1345*3ca0882eSHong Zhang   We do not need to solve the equation; we just use SNES to approximate the Jacobian
1346*3ca0882eSHong Zhang */
1347*3ca0882eSHong Zhang static PetscErrorCode SNESTSFormFunction_RK(SNES snes,Vec x,Vec y,TS ts)
1348*3ca0882eSHong Zhang {
1349*3ca0882eSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
1350*3ca0882eSHong Zhang   DM             dm,dmsave;
1351*3ca0882eSHong Zhang   PetscErrorCode ierr;
1352*3ca0882eSHong Zhang 
1353*3ca0882eSHong Zhang   PetscFunctionBegin;
1354*3ca0882eSHong Zhang   ierr   = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1355*3ca0882eSHong Zhang   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
1356*3ca0882eSHong Zhang   dmsave = ts->dm;
1357*3ca0882eSHong Zhang   ts->dm = dm;
1358*3ca0882eSHong Zhang   ierr   = TSComputeRHSFunction(ts,rk->stage_time,x,y);CHKERRQ(ierr);
1359*3ca0882eSHong Zhang   ts->dm = dmsave;
1360*3ca0882eSHong Zhang   PetscFunctionReturn(0);
1361*3ca0882eSHong Zhang }
1362*3ca0882eSHong Zhang 
1363*3ca0882eSHong Zhang static PetscErrorCode SNESTSFormJacobian_RK(SNES snes,Vec x,Mat A,Mat B,TS ts)
1364*3ca0882eSHong Zhang {
1365*3ca0882eSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
1366*3ca0882eSHong Zhang   DM             dm,dmsave;
1367*3ca0882eSHong Zhang   PetscErrorCode ierr;
1368*3ca0882eSHong Zhang 
1369*3ca0882eSHong Zhang   PetscFunctionBegin;
1370*3ca0882eSHong Zhang   ierr   = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1371*3ca0882eSHong Zhang   dmsave = ts->dm;
1372*3ca0882eSHong Zhang   ts->dm = dm;
1373*3ca0882eSHong Zhang   ierr   = TSComputeRHSJacobian(ts,rk->stage_time,x,A,B);CHKERRQ(ierr);
1374*3ca0882eSHong Zhang   ts->dm = dmsave;
1375*3ca0882eSHong Zhang   PetscFunctionReturn(0);
1376*3ca0882eSHong Zhang }
1377*3ca0882eSHong Zhang 
137821052c0cSHong Zhang /*@C
137921052c0cSHong Zhang   TSRKSetMultirate - Use the interpolation-based multirate RK method
138021052c0cSHong Zhang 
138121052c0cSHong Zhang   Logically collective
138221052c0cSHong Zhang 
138321052c0cSHong Zhang   Input Parameter:
138421052c0cSHong Zhang +  ts - timestepping context
138521052c0cSHong 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
138621052c0cSHong Zhang 
138721052c0cSHong Zhang   Options Database:
138821052c0cSHong Zhang .   -ts_rk_multirate - <true,false>
138921052c0cSHong Zhang 
139021052c0cSHong Zhang   Notes:
139121052c0cSHong 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).
139221052c0cSHong Zhang 
139321052c0cSHong Zhang   Level: intermediate
139421052c0cSHong Zhang 
139521052c0cSHong Zhang .seealso: TSRKGetMultirate()
139621052c0cSHong Zhang @*/
139721052c0cSHong Zhang PetscErrorCode TSRKSetMultirate(TS ts,PetscBool use_multirate)
139821052c0cSHong Zhang {
1399f06fb28fSHong Zhang   PetscErrorCode ierr;
14007dbacdf2SHong Zhang 
14018a4be4dbSHong Zhang   PetscFunctionBegin;
1402f06fb28fSHong Zhang   ierr = PetscTryMethod(ts,"TSRKSetMultirate_C",(TS,PetscBool),(ts,use_multirate));CHKERRQ(ierr);
140321052c0cSHong Zhang   PetscFunctionReturn(0);
140421052c0cSHong Zhang }
140521052c0cSHong Zhang 
140621052c0cSHong Zhang /*@C
140721052c0cSHong Zhang   TSRKGetMultirate - Gets whether to Use the interpolation-based multirate RK method
140821052c0cSHong Zhang 
140921052c0cSHong Zhang   Not collective
141021052c0cSHong Zhang 
141121052c0cSHong Zhang   Input Parameter:
141221052c0cSHong Zhang .  ts - timestepping context
141321052c0cSHong Zhang 
141421052c0cSHong Zhang   Output Parameter:
141521052c0cSHong Zhang .  use_multirate - PETSC_TRUE if the multirate RK method is enabled, PETSC_FALSE otherwise
141621052c0cSHong Zhang 
141721052c0cSHong Zhang   Level: intermediate
141821052c0cSHong Zhang 
141921052c0cSHong Zhang .seealso: TSRKSetMultirate()
142021052c0cSHong Zhang @*/
142121052c0cSHong Zhang PetscErrorCode TSRKGetMultirate(TS ts,PetscBool *use_multirate)
142221052c0cSHong Zhang {
1423f06fb28fSHong Zhang   PetscErrorCode ierr;
14247dbacdf2SHong Zhang 
14257dbacdf2SHong Zhang   PetscFunctionBegin;
1426f06fb28fSHong Zhang   ierr = PetscUseMethod(ts,"TSRKGetMultirate_C",(TS,PetscBool*),(ts,use_multirate));CHKERRQ(ierr);
142721052c0cSHong Zhang   PetscFunctionReturn(0);
142821052c0cSHong Zhang }
142921052c0cSHong Zhang 
1430ebe3b25bSBarry Smith /*MC
1431f68a32c8SEmil Constantinescu       TSRK - ODE and DAE solver using Runge-Kutta schemes
1432ebe3b25bSBarry Smith 
1433f68a32c8SEmil Constantinescu   The user should provide the right hand side of the equation
1434f68a32c8SEmil Constantinescu   using TSSetRHSFunction().
1435ebe3b25bSBarry Smith 
1436f68a32c8SEmil Constantinescu   Notes:
143798164e13SLisandro Dalcin   The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type
1438ebe3b25bSBarry Smith 
1439d41222bbSBarry Smith   Level: beginner
1440d41222bbSBarry Smith 
1441f68a32c8SEmil Constantinescu .seealso:  TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3,
14420fe4d17eSHong Zhang            TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister(), TSRKSetMultirate(), TSRKGetMultirate()
1443ebe3b25bSBarry Smith 
1444ebe3b25bSBarry Smith M*/
14458cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
1446e4dd521cSBarry Smith {
14471566a47fSLisandro Dalcin   TS_RK          *rk;
1448dfbe8321SBarry Smith   PetscErrorCode ierr;
1449e4dd521cSBarry Smith 
1450e4dd521cSBarry Smith   PetscFunctionBegin;
1451f68a32c8SEmil Constantinescu   ierr = TSRKInitializePackage();CHKERRQ(ierr);
1452f68a32c8SEmil Constantinescu 
1453f68a32c8SEmil Constantinescu   ts->ops->reset          = TSReset_RK;
14545f70b478SJed Brown   ts->ops->destroy        = TSDestroy_RK;
14555f70b478SJed Brown   ts->ops->view           = TSView_RK;
1456f68a32c8SEmil Constantinescu   ts->ops->load           = TSLoad_RK;
1457f68a32c8SEmil Constantinescu   ts->ops->setup          = TSSetUp_RK;
1458f68a32c8SEmil Constantinescu   ts->ops->interpolate    = TSInterpolate_RK;
14592c9cb062Sluzhanghpp   ts->ops->step           = TSStep_RK;
1460f68a32c8SEmil Constantinescu   ts->ops->evaluatestep   = TSEvaluateStep_RK;
1461fecfb714SLisandro Dalcin   ts->ops->rollback       = TSRollBack_RK;
1462f68a32c8SEmil Constantinescu   ts->ops->setfromoptions = TSSetFromOptions_RK;
1463ff22ae23SHong Zhang   ts->ops->getstages      = TSGetStages_RK;
1464e4dd521cSBarry Smith 
1465*3ca0882eSHong Zhang   ts->ops->snesfunction    = SNESTSFormFunction_RK;
1466*3ca0882eSHong Zhang   ts->ops->snesjacobian    = SNESTSFormJacobian_RK;
14670f7a1166SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_RK;
146813af1a74SHong Zhang   ts->ops->adjointsetup    = TSAdjointSetUp_RK;
146913af1a74SHong Zhang   ts->ops->adjointstep     = TSAdjointStep_RK;
147013af1a74SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_RK;
14710f7a1166SHong Zhang 
147213af1a74SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_RK;
1473922a638cSHong Zhang   ts->ops->forwardsetup    = TSForwardSetUp_RK;
1474922a638cSHong Zhang   ts->ops->forwardreset    = TSForwardReset_RK;
1475922a638cSHong Zhang   ts->ops->forwardstep     = TSForwardStep_RK;
1476922a638cSHong Zhang   ts->ops->forwardgetstages= TSForwardGetStages_RK;
1477922a638cSHong Zhang 
14781566a47fSLisandro Dalcin   ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr);
14791566a47fSLisandro Dalcin   ts->data = (void*)rk;
1480e4dd521cSBarry Smith 
1481f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr);
1482f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr);
148321052c0cSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",TSRKSetMultirate_RK);CHKERRQ(ierr);
148421052c0cSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",TSRKGetMultirate_RK);CHKERRQ(ierr);
1485be5899b3SLisandro Dalcin 
1486be5899b3SLisandro Dalcin   ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr);
148790b67129SHong Zhang   rk->dtratio = 1;
1488e4dd521cSBarry Smith   PetscFunctionReturn(0);
1489e4dd521cSBarry Smith }
1490