xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision cd4cee2d45328e4919add0f7279ebed5203e0566)
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*/
124262ff7bbSBarry Smith 
125f68a32c8SEmil Constantinescu /*@C
126f68a32c8SEmil Constantinescu   TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK
127262ff7bbSBarry Smith 
128f68a32c8SEmil Constantinescu   Not Collective, but should be called by all processes which will need the schemes to be registered
129262ff7bbSBarry Smith 
130f68a32c8SEmil Constantinescu   Level: advanced
131262ff7bbSBarry Smith 
132f68a32c8SEmil Constantinescu .keywords: TS, TSRK, register, all
133262ff7bbSBarry Smith 
134f68a32c8SEmil Constantinescu .seealso:  TSRKRegisterDestroy()
135262ff7bbSBarry Smith @*/
136f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void)
137262ff7bbSBarry Smith {
1384ac538c5SBarry Smith   PetscErrorCode ierr;
139262ff7bbSBarry Smith 
140262ff7bbSBarry Smith   PetscFunctionBegin;
141f68a32c8SEmil Constantinescu   if (TSRKRegisterAllCalled) PetscFunctionReturn(0);
142f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_TRUE;
143e4dd521cSBarry Smith 
1444e82626cSLisandro Dalcin #define RC PetscRealConstant
145e4dd521cSBarry Smith   {
146f68a32c8SEmil Constantinescu     const PetscReal
1474e82626cSLisandro Dalcin       A[1][1] = {{0}},
1484e82626cSLisandro Dalcin       b[1]    = {RC(1.0)};
1493a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
150e4dd521cSBarry Smith   }
151db046809SBarry Smith   {
152f68a32c8SEmil Constantinescu     const PetscReal
1534e82626cSLisandro Dalcin       A[2][2]   = {{0,0},
1544e82626cSLisandro Dalcin                    {RC(1.0),0}},
1554e82626cSLisandro Dalcin       b[2]      =  {RC(0.5),RC(0.5)},
1564e82626cSLisandro Dalcin       bembed[2] =  {RC(1.0),0};
1573a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
158db046809SBarry Smith   }
159f68a32c8SEmil Constantinescu   {
160f68a32c8SEmil Constantinescu     const PetscReal
16117f6384fSBarry Smith       A[3][3] = {{0,0,0},
1624e82626cSLisandro Dalcin                  {RC(2.0)/RC(3.0),0,0},
1634e82626cSLisandro Dalcin                  {RC(-1.0)/RC(3.0),RC(1.0),0}},
1644e82626cSLisandro Dalcin       b[3]    =  {RC(0.25),RC(0.5),RC(0.25)};
1653a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
166fdefd5e5SDebojyoti Ghosh   }
167fdefd5e5SDebojyoti Ghosh   {
168fdefd5e5SDebojyoti Ghosh     const PetscReal
16917f6384fSBarry Smith       A[4][4]   = {{0,0,0,0},
1704e82626cSLisandro Dalcin                    {RC(1.0)/RC(2.0),0,0,0},
1714e82626cSLisandro Dalcin                    {0,RC(3.0)/RC(4.0),0,0},
1724e82626cSLisandro Dalcin                    {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}},
1734e82626cSLisandro Dalcin       b[4]      =  {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0},
1744e82626cSLisandro 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)};
1753a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
176db046809SBarry Smith   }
177f68a32c8SEmil Constantinescu   {
178f68a32c8SEmil Constantinescu     const PetscReal
179f68a32c8SEmil Constantinescu       A[4][4] = {{0,0,0,0},
1804e82626cSLisandro Dalcin                  {RC(0.5),0,0,0},
1814e82626cSLisandro Dalcin                  {0,RC(0.5),0,0},
1824e82626cSLisandro Dalcin                  {0,0,RC(1.0),0}},
1834e82626cSLisandro 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)};
1843a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
185db046809SBarry Smith   }
186f68a32c8SEmil Constantinescu   {
187f68a32c8SEmil Constantinescu     const PetscReal
18817f6384fSBarry Smith       A[6][6]   = {{0,0,0,0,0,0},
1894e82626cSLisandro Dalcin                    {RC(0.25),0,0,0,0,0},
1904e82626cSLisandro Dalcin                    {RC(3.0)/RC(32.0),RC(9.0)/RC(32.0),0,0,0,0},
1914e82626cSLisandro Dalcin                    {RC(1932.0)/RC(2197.0),RC(-7200.0)/RC(2197.0),RC(7296.0)/RC(2197.0),0,0,0},
1924e82626cSLisandro Dalcin                    {RC(439.0)/RC(216.0),RC(-8.0),RC(3680.0)/RC(513.0),RC(-845.0)/RC(4104.0),0,0},
1934e82626cSLisandro 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}},
1944e82626cSLisandro 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)},
1954e82626cSLisandro 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};
1963a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
197fdefd5e5SDebojyoti Ghosh   }
198fdefd5e5SDebojyoti Ghosh   {
199fdefd5e5SDebojyoti Ghosh     const PetscReal
20017f6384fSBarry Smith       A[7][7]   = {{0,0,0,0,0,0,0},
2014e82626cSLisandro Dalcin                    {RC(1.0)/RC(5.0),0,0,0,0,0,0},
2024e82626cSLisandro Dalcin                    {RC(3.0)/RC(40.0),RC(9.0)/RC(40.0),0,0,0,0,0},
2034e82626cSLisandro Dalcin                    {RC(44.0)/RC(45.0),RC(-56.0)/RC(15.0),RC(32.0)/RC(9.0),0,0,0,0},
2044e82626cSLisandro 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},
2054e82626cSLisandro 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},
2064e82626cSLisandro 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}},
2074e82626cSLisandro 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},
208a685a763Sluzhanghpp         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)},
209a685a763Sluzhanghpp         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)},
210a685a763Sluzhanghpp                     {0,0,0,0,0},
211a685a763Sluzhanghpp                     {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)},
212a685a763Sluzhanghpp                     {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)},
213a685a763Sluzhanghpp                     {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)},
214a685a763Sluzhanghpp                     {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)},
215a685a763Sluzhanghpp                     {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)}};
216a685a763Sluzhanghpp 
217a685a763Sluzhanghpp         ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,5,binterp[0]);CHKERRQ(ierr);
218f68a32c8SEmil Constantinescu   }
21905e23783SLisandro Dalcin   {
22005e23783SLisandro Dalcin     const PetscReal
22117f6384fSBarry Smith       A[8][8]   = {{0,0,0,0,0,0,0,0},
2224e82626cSLisandro Dalcin                    {RC(1.0)/RC(6.0),0,0,0,0,0,0,0},
2234e82626cSLisandro Dalcin                    {RC(2.0)/RC(27.0),RC(4.0)/RC(27.0),0,0,0,0,0,0},
2244e82626cSLisandro Dalcin                    {RC(183.0)/RC(1372.0),RC(-162.0)/RC(343.0),RC(1053.0)/RC(1372.0),0,0,0,0,0},
2254e82626cSLisandro 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},
2264e82626cSLisandro 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},
2274e82626cSLisandro 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},
2284e82626cSLisandro 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}},
2294e82626cSLisandro 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},
2304e82626cSLisandro 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)};
23105e23783SLisandro Dalcin     ierr = TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
23205e23783SLisandro Dalcin   }
2334e82626cSLisandro Dalcin #undef RC
234db046809SBarry Smith   PetscFunctionReturn(0);
235db046809SBarry Smith }
236db046809SBarry Smith 
237f68a32c8SEmil Constantinescu /*@C
238f68a32c8SEmil Constantinescu    TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister().
239f68a32c8SEmil Constantinescu 
240f68a32c8SEmil Constantinescu    Not Collective
241f68a32c8SEmil Constantinescu 
242f68a32c8SEmil Constantinescu    Level: advanced
243f68a32c8SEmil Constantinescu 
244f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy
245f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll()
246f68a32c8SEmil Constantinescu @*/
247f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void)
248e4dd521cSBarry Smith {
249f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
250f68a32c8SEmil Constantinescu   RKTableauLink  link;
251f68a32c8SEmil Constantinescu 
252f68a32c8SEmil Constantinescu   PetscFunctionBegin;
253f68a32c8SEmil Constantinescu   while ((link = RKTableauList)) {
254f68a32c8SEmil Constantinescu     RKTableau t = &link->tab;
255f68a32c8SEmil Constantinescu     RKTableauList = link->next;
256f68a32c8SEmil Constantinescu     ierr = PetscFree3(t->A,t->b,t->c);  CHKERRQ(ierr);
257f68a32c8SEmil Constantinescu     ierr = PetscFree (t->bembed);       CHKERRQ(ierr);
258f68a32c8SEmil Constantinescu     ierr = PetscFree (t->binterp);      CHKERRQ(ierr);
259f68a32c8SEmil Constantinescu     ierr = PetscFree (t->name);         CHKERRQ(ierr);
260f68a32c8SEmil Constantinescu     ierr = PetscFree (link);            CHKERRQ(ierr);
261f68a32c8SEmil Constantinescu   }
262f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_FALSE;
263f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
264f68a32c8SEmil Constantinescu }
265f68a32c8SEmil Constantinescu 
266f68a32c8SEmil Constantinescu /*@C
267f68a32c8SEmil Constantinescu   TSRKInitializePackage - This function initializes everything in the TSRK package. It is called
2688a690491SBarry Smith   from TSInitializePackage().
269f68a32c8SEmil Constantinescu 
270f68a32c8SEmil Constantinescu   Level: developer
271f68a32c8SEmil Constantinescu 
272f68a32c8SEmil Constantinescu .keywords: TS, TSRK, initialize, package
273f68a32c8SEmil Constantinescu .seealso: PetscInitialize()
274f68a32c8SEmil Constantinescu @*/
275f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void)
276f68a32c8SEmil Constantinescu {
277186e87acSLisandro Dalcin   PetscErrorCode ierr;
278e4dd521cSBarry Smith 
279e4dd521cSBarry Smith   PetscFunctionBegin;
280f68a32c8SEmil Constantinescu   if (TSRKPackageInitialized) PetscFunctionReturn(0);
281f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_TRUE;
282f68a32c8SEmil Constantinescu   ierr = TSRKRegisterAll();CHKERRQ(ierr);
283f68a32c8SEmil Constantinescu   ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr);
284f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
285f68a32c8SEmil Constantinescu }
286186e87acSLisandro Dalcin 
287f68a32c8SEmil Constantinescu /*@C
288f68a32c8SEmil Constantinescu   TSRKFinalizePackage - This function destroys everything in the TSRK package. It is
289f68a32c8SEmil Constantinescu   called from PetscFinalize().
290186e87acSLisandro Dalcin 
291f68a32c8SEmil Constantinescu   Level: developer
292f68a32c8SEmil Constantinescu 
293f68a32c8SEmil Constantinescu .keywords: Petsc, destroy, package
294f68a32c8SEmil Constantinescu .seealso: PetscFinalize()
295f68a32c8SEmil Constantinescu @*/
296f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void)
297f68a32c8SEmil Constantinescu {
298f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
299f68a32c8SEmil Constantinescu 
300f68a32c8SEmil Constantinescu   PetscFunctionBegin;
301f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_FALSE;
302f68a32c8SEmil Constantinescu   ierr = TSRKRegisterDestroy();CHKERRQ(ierr);
303f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
304f68a32c8SEmil Constantinescu }
305f68a32c8SEmil Constantinescu 
306f68a32c8SEmil Constantinescu /*@C
307f68a32c8SEmil Constantinescu    TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
308f68a32c8SEmil Constantinescu 
309f68a32c8SEmil Constantinescu    Not Collective, but the same schemes should be registered on all processes on which they will be used
310f68a32c8SEmil Constantinescu 
311f68a32c8SEmil Constantinescu    Input Parameters:
312f68a32c8SEmil Constantinescu +  name - identifier for method
313f68a32c8SEmil Constantinescu .  order - approximation order of method
314f68a32c8SEmil Constantinescu .  s - number of stages, this is the dimension of the matrices below
315f68a32c8SEmil Constantinescu .  A - stage coefficients (dimension s*s, row-major)
316f68a32c8SEmil Constantinescu .  b - step completion table (dimension s; NULL to use last row of A)
317f68a32c8SEmil Constantinescu .  c - abscissa (dimension s; NULL to use row sums of A)
318f68a32c8SEmil Constantinescu .  bembed - completion table for embedded method (dimension s; NULL if not available)
3193a8a9803SLisandro Dalcin .  p - Order of the interpolation scheme, equal to the number of columns of binterp
3203a8a9803SLisandro Dalcin -  binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1)
321f68a32c8SEmil Constantinescu 
322f68a32c8SEmil Constantinescu    Notes:
323f68a32c8SEmil Constantinescu    Several RK methods are provided, this function is only needed to create new methods.
324f68a32c8SEmil Constantinescu 
325f68a32c8SEmil Constantinescu    Level: advanced
326f68a32c8SEmil Constantinescu 
327f68a32c8SEmil Constantinescu .keywords: TS, register
328f68a32c8SEmil Constantinescu 
329f68a32c8SEmil Constantinescu .seealso: TSRK
330f68a32c8SEmil Constantinescu @*/
331f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s,
332f68a32c8SEmil Constantinescu                             const PetscReal A[],const PetscReal b[],const PetscReal c[],
3333a8a9803SLisandro Dalcin                             const PetscReal bembed[],PetscInt p,const PetscReal binterp[])
334f68a32c8SEmil Constantinescu {
335f68a32c8SEmil Constantinescu   PetscErrorCode  ierr;
336f68a32c8SEmil Constantinescu   RKTableauLink   link;
337f68a32c8SEmil Constantinescu   RKTableau       t;
338f68a32c8SEmil Constantinescu   PetscInt        i,j;
339f68a32c8SEmil Constantinescu 
340f68a32c8SEmil Constantinescu   PetscFunctionBegin;
3413a8a9803SLisandro Dalcin   PetscValidCharPointer(name,1);
3423a8a9803SLisandro Dalcin   PetscValidRealPointer(A,4);
3433a8a9803SLisandro Dalcin   if (b) PetscValidRealPointer(b,5);
3443a8a9803SLisandro Dalcin   if (c) PetscValidRealPointer(c,6);
3453a8a9803SLisandro Dalcin   if (bembed) PetscValidRealPointer(bembed,7);
3463a8a9803SLisandro Dalcin   if (binterp || p > 1) PetscValidRealPointer(binterp,9);
3473a8a9803SLisandro Dalcin 
3481d36bdfdSBarry Smith   ierr = TSRKInitializePackage();CHKERRQ(ierr);
34995dccacaSBarry Smith   ierr = PetscNew(&link);CHKERRQ(ierr);
350f68a32c8SEmil Constantinescu   t = &link->tab;
3513a8a9803SLisandro Dalcin 
352f68a32c8SEmil Constantinescu   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
353f68a32c8SEmil Constantinescu   t->order = order;
354f68a32c8SEmil Constantinescu   t->s = s;
355dcca6d9dSJed Brown   ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr);
356f68a32c8SEmil Constantinescu   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
357f68a32c8SEmil Constantinescu   if (b)  { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); }
358f68a32c8SEmil Constantinescu   else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i];
359f68a32c8SEmil Constantinescu   if (c)  { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); }
360f68a32c8SEmil 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];
361d760c35bSDebojyoti Ghosh   t->FSAL = PETSC_TRUE;
362d760c35bSDebojyoti Ghosh   for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE;
3633a8a9803SLisandro Dalcin 
364f68a32c8SEmil Constantinescu   if (bembed) {
365785e854fSJed Brown     ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr);
366f68a32c8SEmil Constantinescu     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
367f68a32c8SEmil Constantinescu   }
368f68a32c8SEmil Constantinescu 
3693a8a9803SLisandro Dalcin   if (!binterp) { p = 1; binterp = t->b; }
3703a8a9803SLisandro Dalcin   t->p = p;
3713a8a9803SLisandro Dalcin   ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr);
3723a8a9803SLisandro Dalcin   ierr = PetscMemcpy(t->binterp,binterp,s*p*sizeof(binterp[0]));CHKERRQ(ierr);
3733a8a9803SLisandro Dalcin 
374f68a32c8SEmil Constantinescu   link->next = RKTableauList;
375f68a32c8SEmil Constantinescu   RKTableauList = link;
376f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
377f68a32c8SEmil Constantinescu }
378f68a32c8SEmil Constantinescu 
379e4dd521cSBarry Smith /*
3802c9cb062Sluzhanghpp  This is for single-step RK method
381f68a32c8SEmil Constantinescu  The step completion formula is
382e4dd521cSBarry Smith 
383f68a32c8SEmil Constantinescu  x1 = x0 + h b^T YdotRHS
384f68a32c8SEmil Constantinescu 
385f68a32c8SEmil Constantinescu  This function can be called before or after ts->vec_sol has been updated.
386f68a32c8SEmil Constantinescu  Suppose we have a completion formula (b) and an embedded formula (be) of different order.
387f68a32c8SEmil Constantinescu  We can write
388f68a32c8SEmil Constantinescu 
389f68a32c8SEmil Constantinescu  x1e = x0 + h be^T YdotRHS
390f68a32c8SEmil Constantinescu      = x1 - h b^T YdotRHS + h be^T YdotRHS
391f68a32c8SEmil Constantinescu      = x1 + h (be - b)^T YdotRHS
392f68a32c8SEmil Constantinescu 
393f68a32c8SEmil Constantinescu  so we can evaluate the method with different order even after the step has been optimistically completed.
394e4dd521cSBarry Smith */
395f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done)
396f68a32c8SEmil Constantinescu {
397f68a32c8SEmil Constantinescu   TS_RK          *rk   = (TS_RK*)ts->data;
398f68a32c8SEmil Constantinescu   RKTableau      tab  = rk->tableau;
399f68a32c8SEmil Constantinescu   PetscScalar    *w    = rk->work;
400f68a32c8SEmil Constantinescu   PetscReal      h;
401f68a32c8SEmil Constantinescu   PetscInt       s    = tab->s,j;
402f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
403f68a32c8SEmil Constantinescu 
404f68a32c8SEmil Constantinescu   PetscFunctionBegin;
405f68a32c8SEmil Constantinescu   switch (rk->status) {
406f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
407f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
408f68a32c8SEmil Constantinescu     h = ts->time_step; break;
409f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
410be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
411f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
412f68a32c8SEmil Constantinescu   }
413f68a32c8SEmil Constantinescu   if (order == tab->order) {
414f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) {
415f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
41690b67129SHong Zhang       for (j=0; j<s; j++) w[j] = h*tab->b[j]/rk->dtratio;
417f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
418f68a32c8SEmil Constantinescu     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
419f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
420f68a32c8SEmil Constantinescu   } else if (order == tab->order-1) {
421f68a32c8SEmil Constantinescu     if (!tab->bembed) goto unavailable;
422f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) { /*Complete with the embedded method (be)*/
423f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
424f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
425f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
426f68a32c8SEmil Constantinescu     } else {  /*Rollback and re-complete using (be-b) */
427f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
428f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
429f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
430f68a32c8SEmil Constantinescu     }
431f68a32c8SEmil Constantinescu     if (done) *done = PETSC_TRUE;
432f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
433f68a32c8SEmil Constantinescu   }
434f68a32c8SEmil Constantinescu unavailable:
435f68a32c8SEmil Constantinescu   if (done) *done = PETSC_FALSE;
436a7fac7c2SEmil 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);
437f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
438f68a32c8SEmil Constantinescu }
439f68a32c8SEmil Constantinescu 
4400f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts)
4410f7a1166SHong Zhang {
4420f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
443*cd4cee2dSHong Zhang   TS              quadts = ts->quadraturets;
4440f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
4450f7a1166SHong Zhang   const PetscInt  s = tab->s;
4460f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
4470f7a1166SHong Zhang   Vec             *Y = rk->Y;
4480f7a1166SHong Zhang   PetscInt        i;
4490f7a1166SHong Zhang   PetscErrorCode  ierr;
4500f7a1166SHong Zhang 
4510f7a1166SHong Zhang   PetscFunctionBegin;
452*cd4cee2dSHong Zhang   /* No need to backup quadts->vec_sol since it can be reverted in TSRollBack_RK */
4530f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
454*cd4cee2dSHong Zhang     /* Evolve quadrature TS solution to compute integrals */
455*cd4cee2dSHong Zhang     ierr = TSComputeRHSFunction(quadts,rk->ptime+rk->time_step*c[i],Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
456*cd4cee2dSHong Zhang     ierr = VecAXPY(quadts->vec_sol,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
4570f7a1166SHong Zhang   }
4580f7a1166SHong Zhang   PetscFunctionReturn(0);
4590f7a1166SHong Zhang }
4600f7a1166SHong Zhang 
4610f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts)
4620f7a1166SHong Zhang {
4630f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
4640f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
465*cd4cee2dSHong Zhang   TS              quadts = ts->quadraturets;
4660f7a1166SHong Zhang   const PetscInt  s = tab->s;
4670f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
4680f7a1166SHong Zhang   Vec             *Y = rk->Y;
4690f7a1166SHong Zhang   PetscInt        i;
4700f7a1166SHong Zhang   PetscErrorCode  ierr;
4710f7a1166SHong Zhang 
4720f7a1166SHong Zhang   PetscFunctionBegin;
4730f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
474*cd4cee2dSHong Zhang     /* Evolve quadrature TS solution to compute integrals */
475*cd4cee2dSHong Zhang     ierr = TSComputeRHSFunction(quadts,ts->ptime+ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
476*cd4cee2dSHong Zhang     ierr = VecAXPY(quadts->vec_sol,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
4770f7a1166SHong Zhang   }
4780f7a1166SHong Zhang   PetscFunctionReturn(0);
4790f7a1166SHong Zhang }
4800f7a1166SHong Zhang 
481fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts)
482fecfb714SLisandro Dalcin {
483fecfb714SLisandro Dalcin   TS_RK           *rk = (TS_RK*)ts->data;
484*cd4cee2dSHong Zhang   TS              quadts = ts->quadraturets;
485fecfb714SLisandro Dalcin   RKTableau       tab = rk->tableau;
486fecfb714SLisandro Dalcin   const PetscInt  s  = tab->s;
487*cd4cee2dSHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
488fecfb714SLisandro Dalcin   PetscScalar     *w = rk->work;
489*cd4cee2dSHong Zhang   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
490fecfb714SLisandro Dalcin   PetscInt        j;
491fecfb714SLisandro Dalcin   PetscReal       h;
492fecfb714SLisandro Dalcin   PetscErrorCode  ierr;
493fecfb714SLisandro Dalcin 
494fecfb714SLisandro Dalcin   PetscFunctionBegin;
495fecfb714SLisandro Dalcin   switch (rk->status) {
496fecfb714SLisandro Dalcin   case TS_STEP_INCOMPLETE:
497fecfb714SLisandro Dalcin   case TS_STEP_PENDING:
498fecfb714SLisandro Dalcin     h = ts->time_step; break;
499fecfb714SLisandro Dalcin   case TS_STEP_COMPLETE:
500fecfb714SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
501fecfb714SLisandro Dalcin   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
502fecfb714SLisandro Dalcin   }
503fecfb714SLisandro Dalcin   for (j=0; j<s; j++) w[j] = -h*b[j];
504fecfb714SLisandro Dalcin   ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr);
505*cd4cee2dSHong Zhang   if (quadts && ts->costintegralfwd) {
506*cd4cee2dSHong Zhang     for (j=0; j<s; j++) {
507*cd4cee2dSHong Zhang       /* Revert the quadrature TS solution */
508*cd4cee2dSHong Zhang       ierr = TSComputeRHSFunction(quadts,rk->ptime+h*c[j],Y[j],ts->vec_costintegrand);CHKERRQ(ierr);
509*cd4cee2dSHong Zhang       ierr = VecAXPY(quadts->vec_sol,-h*b[j],ts->vec_costintegrand);CHKERRQ(ierr);
510*cd4cee2dSHong Zhang     }
511*cd4cee2dSHong Zhang   }
512fecfb714SLisandro Dalcin   PetscFunctionReturn(0);
513fecfb714SLisandro Dalcin }
514fecfb714SLisandro Dalcin 
515922a638cSHong Zhang static PetscErrorCode TSForwardStep_RK(TS ts)
516922a638cSHong Zhang {
517922a638cSHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
518922a638cSHong Zhang   RKTableau       tab = rk->tableau;
519922a638cSHong Zhang   Mat             J,*MatsFwdSensipTemp = rk->MatsFwdSensipTemp;
520922a638cSHong Zhang   const PetscInt  s = tab->s;
521922a638cSHong Zhang   const PetscReal *A = tab->A,*c = tab->c,*b = tab->b;
522922a638cSHong Zhang   Vec             *Y = rk->Y;
523922a638cSHong Zhang   PetscInt        i,j;
524922a638cSHong Zhang   PetscReal       stage_time,h = ts->time_step;
525922a638cSHong Zhang   PetscBool       zero;
526922a638cSHong Zhang   PetscErrorCode  ierr;
527922a638cSHong Zhang 
528922a638cSHong Zhang   PetscFunctionBegin;
529922a638cSHong Zhang   ierr = MatCopy(ts->mat_sensip,rk->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
530922a638cSHong Zhang   ierr = TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr);
531922a638cSHong Zhang 
532922a638cSHong Zhang   for (i=0; i<s; i++) {
533922a638cSHong Zhang     stage_time = ts->ptime + h*c[i];
534922a638cSHong Zhang     zero = PETSC_FALSE;
535922a638cSHong Zhang     if (b[i] == 0 && i == s-1) zero = PETSC_TRUE;
536922a638cSHong Zhang     /* TLM Stage values */
537922a638cSHong Zhang     if(!i) {
538922a638cSHong Zhang       ierr = MatCopy(ts->mat_sensip,rk->MatsFwdStageSensip[i],SAME_NONZERO_PATTERN);CHKERRQ(ierr);
539922a638cSHong Zhang     } else if (!zero) {
540922a638cSHong Zhang       ierr = MatZeroEntries(rk->MatsFwdStageSensip[i]);CHKERRQ(ierr);
541922a638cSHong Zhang       for (j=0; j<i; j++) {
542922a638cSHong Zhang         ierr = MatAXPY(rk->MatsFwdStageSensip[i],h*A[i*s+j],MatsFwdSensipTemp[j],SAME_NONZERO_PATTERN);CHKERRQ(ierr);
543922a638cSHong Zhang       }
544922a638cSHong Zhang       ierr = MatAXPY(rk->MatsFwdStageSensip[i],1.,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
545922a638cSHong Zhang     } else {
546922a638cSHong Zhang       ierr = MatZeroEntries(rk->MatsFwdStageSensip[i]);CHKERRQ(ierr);
547922a638cSHong Zhang     }
548922a638cSHong Zhang 
549922a638cSHong Zhang     ierr = TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);CHKERRQ(ierr);
550922a638cSHong Zhang     ierr = MatMatMult(J,rk->MatsFwdStageSensip[i],MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatsFwdSensipTemp[i]);CHKERRQ(ierr);
55113af1a74SHong Zhang     if (ts->Jacprhs) {
55213af1a74SHong Zhang       ierr = TSComputeRHSJacobianP(ts,stage_time,Y[i],ts->Jacprhs);CHKERRQ(ierr); /* get f_p */
55313af1a74SHong Zhang       if (ts->vecs_sensi2p) { /* TLM used for 2nd-order adjoint */
55413af1a74SHong Zhang         PetscScalar *xarr;
55513af1a74SHong Zhang         ierr = MatDenseGetColumn(MatsFwdSensipTemp[i],0,&xarr);CHKERRQ(ierr);
55613af1a74SHong Zhang         ierr = VecPlaceArray(rk->VecDeltaFwdSensipCol,xarr);CHKERRQ(ierr);
55713af1a74SHong Zhang         ierr = MatMultAdd(ts->Jacprhs,ts->vec_dir,rk->VecDeltaFwdSensipCol,rk->VecDeltaFwdSensipCol);CHKERRQ(ierr);
55813af1a74SHong Zhang         ierr = VecResetArray(rk->VecDeltaFwdSensipCol);CHKERRQ(ierr);CHKERRQ(ierr);
55913af1a74SHong Zhang         ierr = MatDenseRestoreColumn(MatsFwdSensipTemp[i],&xarr);CHKERRQ(ierr);
56013af1a74SHong Zhang       } else {
56113af1a74SHong Zhang         ierr = MatAXPY(MatsFwdSensipTemp[i],1.,ts->Jacprhs,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
56213af1a74SHong Zhang       }
563922a638cSHong Zhang     }
564922a638cSHong Zhang   }
565922a638cSHong Zhang 
566922a638cSHong Zhang   for (i=0; i<s; i++) {
567922a638cSHong Zhang     ierr = MatAXPY(ts->mat_sensip,h*b[i],rk->MatsFwdSensipTemp[i],SAME_NONZERO_PATTERN);CHKERRQ(ierr);
568922a638cSHong Zhang   }
569922a638cSHong Zhang   rk->status = TS_STEP_COMPLETE;
570922a638cSHong Zhang   PetscFunctionReturn(0);
571922a638cSHong Zhang }
572922a638cSHong Zhang 
573922a638cSHong Zhang static PetscErrorCode TSForwardGetStages_RK(TS ts,PetscInt *ns,Mat **stagesensip)
574922a638cSHong Zhang {
575922a638cSHong Zhang   TS_RK     *rk = (TS_RK*)ts->data;
576922a638cSHong Zhang   RKTableau tab  = rk->tableau;
577922a638cSHong Zhang 
578922a638cSHong Zhang   PetscFunctionBegin;
579922a638cSHong Zhang   if (ns) *ns = tab->s;
580922a638cSHong Zhang   if (stagesensip) *stagesensip = rk->MatsFwdStageSensip;
581922a638cSHong Zhang   PetscFunctionReturn(0);
582922a638cSHong Zhang }
583922a638cSHong Zhang 
584922a638cSHong Zhang static PetscErrorCode TSForwardSetUp_RK(TS ts)
585922a638cSHong Zhang {
586922a638cSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
587922a638cSHong Zhang   RKTableau      tab  = rk->tableau;
588922a638cSHong Zhang   PetscInt       i;
589922a638cSHong Zhang   PetscErrorCode ierr;
590922a638cSHong Zhang 
591922a638cSHong Zhang   PetscFunctionBegin;
592922a638cSHong Zhang   /* backup sensitivity results for roll-backs */
593922a638cSHong Zhang   ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatFwdSensip0);CHKERRQ(ierr);
594922a638cSHong Zhang 
595922a638cSHong Zhang   ierr = PetscMalloc1(tab->s,&rk->MatsFwdStageSensip);CHKERRQ(ierr);
596922a638cSHong Zhang   ierr = PetscMalloc1(tab->s,&rk->MatsFwdSensipTemp);CHKERRQ(ierr);
597922a638cSHong Zhang   for(i=0; i<tab->s; i++) {
598922a638cSHong Zhang     ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatsFwdStageSensip[i]);CHKERRQ(ierr);
599922a638cSHong Zhang     ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatsFwdSensipTemp[i]);CHKERRQ(ierr);
600922a638cSHong Zhang   }
601922a638cSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&rk->VecDeltaFwdSensipCol);CHKERRQ(ierr);
602922a638cSHong Zhang   PetscFunctionReturn(0);
603922a638cSHong Zhang }
604922a638cSHong Zhang 
605922a638cSHong Zhang static PetscErrorCode TSForwardReset_RK(TS ts)
606922a638cSHong Zhang {
607922a638cSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
608922a638cSHong Zhang   RKTableau      tab  = rk->tableau;
609922a638cSHong Zhang   PetscInt       i;
610922a638cSHong Zhang   PetscErrorCode ierr;
611922a638cSHong Zhang 
612922a638cSHong Zhang   PetscFunctionBegin;
613922a638cSHong Zhang   ierr = MatDestroy(&rk->MatFwdSensip0);CHKERRQ(ierr);
614922a638cSHong Zhang   if (rk->MatsFwdStageSensip) {
615922a638cSHong Zhang     for (i=0; i<tab->s; i++) {
616922a638cSHong Zhang       ierr = MatDestroy(&rk->MatsFwdStageSensip[i]);CHKERRQ(ierr);
617922a638cSHong Zhang     }
618922a638cSHong Zhang     ierr = PetscFree(rk->MatsFwdStageSensip);CHKERRQ(ierr);
619922a638cSHong Zhang   }
620922a638cSHong Zhang   if (rk->MatsFwdSensipTemp) {
621922a638cSHong Zhang     for (i=0; i<tab->s; i++) {
622922a638cSHong Zhang       ierr = MatDestroy(&rk->MatsFwdSensipTemp[i]);CHKERRQ(ierr);
623922a638cSHong Zhang     }
624922a638cSHong Zhang     ierr = PetscFree(rk->MatsFwdSensipTemp);CHKERRQ(ierr);
625922a638cSHong Zhang   }
626922a638cSHong Zhang   ierr = VecDestroy(&rk->VecDeltaFwdSensipCol);CHKERRQ(ierr);
627922a638cSHong Zhang   PetscFunctionReturn(0);
628922a638cSHong Zhang }
629922a638cSHong Zhang 
630f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts)
631f68a32c8SEmil Constantinescu {
632f68a32c8SEmil Constantinescu   TS_RK           *rk  = (TS_RK*)ts->data;
633f68a32c8SEmil Constantinescu   RKTableau       tab  = rk->tableau;
634f68a32c8SEmil Constantinescu   const PetscInt  s = tab->s;
635fecfb714SLisandro Dalcin   const PetscReal *A = tab->A,*c = tab->c;
636f68a32c8SEmil Constantinescu   PetscScalar     *w = rk->work;
637f68a32c8SEmil Constantinescu   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
638d1905564SLisandro Dalcin   PetscBool       FSAL = tab->FSAL;
639f68a32c8SEmil Constantinescu   TSAdapt         adapt;
640fecfb714SLisandro Dalcin   PetscInt        i,j;
641be5899b3SLisandro Dalcin   PetscInt        rejections = 0;
642be5899b3SLisandro Dalcin   PetscBool       stageok,accept = PETSC_TRUE;
643be5899b3SLisandro Dalcin   PetscReal       next_time_step = ts->time_step;
644f68a32c8SEmil Constantinescu   PetscErrorCode  ierr;
645f68a32c8SEmil Constantinescu 
646f68a32c8SEmil Constantinescu   PetscFunctionBegin;
647d1905564SLisandro Dalcin   if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE;
648d1905564SLisandro Dalcin   if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); }
649d1905564SLisandro Dalcin 
650f68a32c8SEmil Constantinescu   rk->status = TS_STEP_INCOMPLETE;
651be5899b3SLisandro Dalcin   while (!ts->reason && rk->status != TS_STEP_COMPLETE) {
652be5899b3SLisandro Dalcin     PetscReal t = ts->ptime;
653f68a32c8SEmil Constantinescu     PetscReal h = ts->time_step;
654f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
6559be3e283SDebojyoti Ghosh       rk->stage_time = t + h*c[i];
6569be3e283SDebojyoti Ghosh       ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr);
657f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
658f68a32c8SEmil Constantinescu       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
659f68a32c8SEmil Constantinescu       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
6609be3e283SDebojyoti Ghosh       ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
661f68a32c8SEmil Constantinescu       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
662be5899b3SLisandro Dalcin       ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr);
663fecfb714SLisandro Dalcin       if (!stageok) goto reject_step;
6648f738a7cSLisandro Dalcin       if (FSAL && !i) continue;
665f68a32c8SEmil Constantinescu       ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
666f68a32c8SEmil Constantinescu     }
667be5899b3SLisandro Dalcin 
668be5899b3SLisandro Dalcin     rk->status = TS_STEP_INCOMPLETE;
669f68a32c8SEmil Constantinescu     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
670f68a32c8SEmil Constantinescu     rk->status = TS_STEP_PENDING;
671f68a32c8SEmil Constantinescu     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
672f68a32c8SEmil Constantinescu     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
6731917a363SLisandro Dalcin     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr);
674fecfb714SLisandro Dalcin     ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
675be5899b3SLisandro Dalcin     rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
676be5899b3SLisandro Dalcin     if (!accept) { /* Roll back the current step */
677fecfb714SLisandro Dalcin       ierr = TSRollBack_RK(ts);CHKERRQ(ierr);
678be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
679be5899b3SLisandro Dalcin       goto reject_step;
680be5899b3SLisandro Dalcin     }
681be5899b3SLisandro Dalcin 
6820f7a1166SHong Zhang     if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
6830f7a1166SHong Zhang       rk->ptime     = ts->ptime;
6840f7a1166SHong Zhang       rk->time_step = ts->time_step;
685493ed95fSHong Zhang     }
686be5899b3SLisandro Dalcin 
687f68a32c8SEmil Constantinescu     ts->ptime += ts->time_step;
688f68a32c8SEmil Constantinescu     ts->time_step = next_time_step;
689f68a32c8SEmil Constantinescu     break;
690be5899b3SLisandro Dalcin 
691be5899b3SLisandro Dalcin     reject_step:
692fecfb714SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
693be5899b3SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
694be5899b3SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
695be5899b3SLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
696f68a32c8SEmil Constantinescu     }
697f68a32c8SEmil Constantinescu   }
698f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
699e4dd521cSBarry Smith }
700e4dd521cSBarry Smith 
701f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts)
702f6a906c0SBarry Smith {
703f6a906c0SBarry Smith   TS_RK          *rk  = (TS_RK*)ts->data;
704be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
705be5899b3SLisandro Dalcin   PetscInt       s   = tab->s;
706f6a906c0SBarry Smith   PetscErrorCode ierr;
707f6a906c0SBarry Smith 
708f6a906c0SBarry Smith   PetscFunctionBegin;
709f6a906c0SBarry Smith   if (ts->adjointsetupcalled++) PetscFunctionReturn(0);
7102e7b7f96SHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecsDeltaLam);CHKERRQ(ierr);
7112e7b7f96SHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&rk->VecsSensiTemp);CHKERRQ(ierr);
712f6a906c0SBarry Smith   if(ts->vecs_sensip) {
7132e7b7f96SHong Zhang     ierr = VecDuplicate(ts->vecs_sensip[0],&rk->VecDeltaMu);CHKERRQ(ierr);
714f6a906c0SBarry Smith   }
71513af1a74SHong Zhang   if (ts->vecs_sensi2) {
71613af1a74SHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecsDeltaLam2);CHKERRQ(ierr);
71713af1a74SHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&rk->VecsSensi2Temp);CHKERRQ(ierr);
71813af1a74SHong Zhang   }
71913af1a74SHong Zhang   if (ts->vecs_sensi2p) {
72013af1a74SHong Zhang     ierr = VecDuplicate(ts->vecs_sensi2p[0],&rk->VecDeltaMu2);CHKERRQ(ierr);
72113af1a74SHong Zhang   }
722f6a906c0SBarry Smith   PetscFunctionReturn(0);
723f6a906c0SBarry Smith }
724f6a906c0SBarry Smith 
72542f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts)
726d2daff3dSHong Zhang {
727c235aa8dSHong Zhang   TS_RK            *rk = (TS_RK*)ts->data;
728*cd4cee2dSHong Zhang   TS               quadts = ts->quadraturets;
729c235aa8dSHong Zhang   RKTableau        tab = rk->tableau;
730*cd4cee2dSHong Zhang   Mat              J,Jquad;
731c235aa8dSHong Zhang   const PetscInt   s = tab->s;
732c235aa8dSHong Zhang   const PetscReal  *A = tab->A,*b = tab->b,*c = tab->c;
73313af1a74SHong Zhang   PetscScalar      *w = rk->work,*xarr;
7342e7b7f96SHong Zhang   Vec              *Y = rk->Y,*VecsDeltaLam = rk->VecsDeltaLam,VecDeltaMu = rk->VecDeltaMu,*VecsSensiTemp = rk->VecsSensiTemp;
73513af1a74SHong Zhang   Vec              *VecsDeltaLam2 = rk->VecsDeltaLam2,VecDeltaMu2 = rk->VecDeltaMu2,*VecsSensi2Temp = rk->VecsSensi2Temp;
736*cd4cee2dSHong Zhang   Vec              VecDRDUTransCol = ts->vec_drdu_col,VecDRDPTransCol = ts->vec_drdp_col;
737b18ea86cSHong Zhang   PetscInt         i,j,nadj;
7383d3eaba7SBarry Smith   PetscReal        t = ts->ptime;
7392e7b7f96SHong Zhang   PetscReal        h = ts->time_step,stage_time;
7402e7b7f96SHong Zhang   PetscBool        zero;
741d2daff3dSHong Zhang   PetscErrorCode   ierr;
742c235aa8dSHong Zhang 
743d2daff3dSHong Zhang   PetscFunctionBegin;
744c235aa8dSHong Zhang   rk->status = TS_STEP_INCOMPLETE;
7453389c7d9SStefano Zampini 
746*cd4cee2dSHong Zhang   ierr = TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr);
747*cd4cee2dSHong Zhang   ierr = TSGetRHSJacobian(quadts,&Jquad,NULL,NULL,NULL);CHKERRQ(ierr);
7482e7b7f96SHong Zhang   for (i=s-1; i>=0; i--) {
7496a1e4597SHong Zhang     if (tab->FSAL && i == s-1) {
7506a1e4597SHong Zhang       /* VecsDeltaLam[nadj*s+s-1] are initialized with zeros and the values never change.*/
7516a1e4597SHong Zhang       continue;
7526a1e4597SHong Zhang     }
7532e7b7f96SHong Zhang     stage_time = t + h*(1.0-c[i]);
7542e7b7f96SHong Zhang     zero = PETSC_FALSE;
7553389c7d9SStefano Zampini     ierr = TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);CHKERRQ(ierr);
756*cd4cee2dSHong Zhang     ierr = TSComputeRHSJacobian(quadts,stage_time,Y[i],Jquad,Jquad);CHKERRQ(ierr); /* get r_u^T */
7573389c7d9SStefano Zampini     if (ts->vecs_sensip) {
75813af1a74SHong Zhang       ierr = TSComputeRHSJacobianP(ts,stage_time,Y[i],ts->Jacprhs);CHKERRQ(ierr); /* get f_p */
759*cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobianP(quadts,stage_time,Y[i],quadts->Jacprhs);CHKERRQ(ierr); /* get f_p for the quadrature */
7603389c7d9SStefano Zampini     }
7613389c7d9SStefano Zampini 
7623389c7d9SStefano Zampini     if (b[i] == 0 && i == s-1) zero = PETSC_TRUE;
7633389c7d9SStefano Zampini 
7642e7b7f96SHong Zhang     for (j=i+1; j<s; j++) w[j-i-1]  = -h*A[j*s+i]; /* coefficients for computing VecsSensiTemp */
7652e7b7f96SHong Zhang 
7662e7b7f96SHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
7673389c7d9SStefano Zampini       /* Stage values of lambda */
7683389c7d9SStefano Zampini       if (!zero) {
76913af1a74SHong Zhang         /* b_i*lambda_{n+1} + \sum_{j=i+1}^s a_{ji}*lambda_{s,j} */
7702e7b7f96SHong Zhang         ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); /* VecDeltaLam is an vec array of size s by numcost */
7712e7b7f96SHong Zhang         ierr = VecScale(VecsSensiTemp[nadj],-h*b[i]);CHKERRQ(ierr);
7722e7b7f96SHong Zhang         ierr = VecMAXPY(VecsSensiTemp[nadj],s-i-1,w,&VecsDeltaLam[nadj*s+i+1]);CHKERRQ(ierr);
7732e7b7f96SHong Zhang         ierr = MatMultTranspose(J,VecsSensiTemp[nadj],VecsDeltaLam[nadj*s+i]);CHKERRQ(ierr);
774*cd4cee2dSHong Zhang         if (Jquad) {
775*cd4cee2dSHong Zhang           ierr = MatDenseGetColumn(Jquad,nadj,&xarr);CHKERRQ(ierr);
776*cd4cee2dSHong Zhang           ierr = VecPlaceArray(VecDRDUTransCol,xarr);CHKERRQ(ierr);
777*cd4cee2dSHong Zhang           ierr = VecAXPY(VecsDeltaLam[nadj*s+i],-h*b[i],VecDRDUTransCol);CHKERRQ(ierr);
778*cd4cee2dSHong Zhang           ierr = VecResetArray(VecDRDUTransCol);CHKERRQ(ierr);
779*cd4cee2dSHong Zhang           ierr = MatDenseRestoreColumn(Jquad,&xarr);CHKERRQ(ierr);
780*cd4cee2dSHong Zhang         }
7813389c7d9SStefano Zampini       } else {
7826a1e4597SHong Zhang         /* \sum_{j=i+1}^s a_{ji}*lambda_{s,j} */
7836a1e4597SHong Zhang         ierr = VecSet(VecsSensiTemp[nadj],0);CHKERRQ(ierr);
7846a1e4597SHong Zhang         ierr = VecMAXPY(VecsSensiTemp[nadj],s-i-1,w,&VecsDeltaLam[nadj*s+i+1]);CHKERRQ(ierr);
7856a1e4597SHong Zhang         ierr = MatMultTranspose(J,VecsSensiTemp[nadj],VecsDeltaLam[nadj*s+i]);CHKERRQ(ierr);
7866a1e4597SHong Zhang         ierr = VecScale(VecsDeltaLam[nadj*s+i],-h);CHKERRQ(ierr);
7873389c7d9SStefano Zampini       }
7886497ce14SHong Zhang 
789ad8e2604SHong Zhang       /* Stage values of mu */
7906497ce14SHong Zhang       if (ts->vecs_sensip) {
7913389c7d9SStefano Zampini         if (!zero) {
79213af1a74SHong Zhang           ierr = MatMultTranspose(ts->Jacprhs,VecsSensiTemp[nadj],VecDeltaMu);CHKERRQ(ierr);
7933389c7d9SStefano Zampini         } else {
7942e7b7f96SHong Zhang           ierr = VecSet(VecDeltaMu,0);CHKERRQ(ierr);
795493ed95fSHong Zhang         }
796*cd4cee2dSHong Zhang         if (quadts->Jacprhs) {
797*cd4cee2dSHong Zhang           ierr = MatDenseGetColumn(quadts->Jacprhs,nadj,&xarr);CHKERRQ(ierr);
798*cd4cee2dSHong Zhang           ierr = VecPlaceArray(VecDRDPTransCol,xarr);CHKERRQ(ierr);
799*cd4cee2dSHong Zhang           ierr = VecAXPY(VecDeltaMu,-h*b[i],VecDRDPTransCol);CHKERRQ(ierr);
800*cd4cee2dSHong Zhang           ierr = VecResetArray(VecDRDPTransCol);CHKERRQ(ierr);
801*cd4cee2dSHong Zhang           ierr = MatDenseRestoreColumn(quadts->Jacprhs,&xarr);CHKERRQ(ierr);
802493ed95fSHong Zhang         }
8032e7b7f96SHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],1.,VecDeltaMu);CHKERRQ(ierr); /* update sensip for each stage */
804ad8e2604SHong Zhang       }
805c235aa8dSHong Zhang     }
80613af1a74SHong Zhang 
80713af1a74SHong Zhang     if (ts->vecs_sensi2 && ts->forward_solve) { /* 2nd-order adjoint, TLM mode has to be turned on */
80813af1a74SHong Zhang       /* Get w1 at t_{n+1} from TLM matrix */
80913af1a74SHong Zhang       ierr = MatDenseGetColumn(rk->MatsFwdStageSensip[i],0,&xarr);CHKERRQ(ierr);
81013af1a74SHong Zhang       ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
81113af1a74SHong Zhang       /* lambda_s^T F_UU w_1 */
81213af1a74SHong Zhang       ierr = TSComputeRHSHessianProductFunction1(ts,stage_time,Y[i],VecsSensiTemp,ts->vec_sensip_col,ts->vecs_guu);CHKERRQ(ierr);
81313af1a74SHong Zhang       /* lambda_s^T F_UP w_2 */
81413af1a74SHong Zhang       ierr = TSComputeRHSHessianProductFunction2(ts,stage_time,Y[i],VecsSensiTemp,ts->vec_dir,ts->vecs_gup);CHKERRQ(ierr);
81513af1a74SHong Zhang       if (ts->vecs_sensi2p) {
81613af1a74SHong Zhang         /* lambda_s^T F_PU w_1 */
81713af1a74SHong Zhang         ierr = TSComputeRHSHessianProductFunction3(ts,stage_time,Y[i],VecsSensiTemp,ts->vec_sensip_col,ts->vecs_gpu);CHKERRQ(ierr);
81813af1a74SHong Zhang         /* lambda_s^T F_PU w_2 */
81913af1a74SHong Zhang         ierr = TSComputeRHSHessianProductFunction4(ts,stage_time,Y[i],VecsSensiTemp,ts->vec_dir,ts->vecs_gpp);CHKERRQ(ierr);
82013af1a74SHong Zhang       }
82113af1a74SHong Zhang       ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
82213af1a74SHong Zhang       ierr = MatDenseRestoreColumn(rk->MatsFwdStageSensip[i],&xarr);CHKERRQ(ierr);
82313af1a74SHong Zhang 
82413af1a74SHong Zhang       for (nadj=0; nadj<ts->numcost; nadj++) {
82513af1a74SHong Zhang         /* Stage values of lambda */
82613af1a74SHong Zhang         if (!zero) {
82713af1a74SHong Zhang           /* h*J_i^T*(b_i*Lambda_{n+1}+\sum_{j=i+1}^s a_{ji}*Lambda_{s,j} */
82813af1a74SHong Zhang           ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
82913af1a74SHong Zhang           ierr = VecScale(VecsSensi2Temp[nadj],-h*b[i]);CHKERRQ(ierr);
83013af1a74SHong Zhang           ierr = VecMAXPY(VecsSensi2Temp[nadj],s-i-1,w,&VecsDeltaLam2[nadj*s+i+1]);CHKERRQ(ierr);
83113af1a74SHong Zhang           ierr = MatMultTranspose(J,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj*s+i]);CHKERRQ(ierr);
83213af1a74SHong Zhang           ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],1.,ts->vecs_guu[nadj]);CHKERRQ(ierr);
83313af1a74SHong Zhang           if (ts->vecs_gup) {
83413af1a74SHong Zhang             ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],1.,ts->vecs_gup[nadj]);CHKERRQ(ierr);
83513af1a74SHong Zhang           }
83613af1a74SHong Zhang         } else {
83713af1a74SHong Zhang           ierr = VecSet(VecsDeltaLam2[nadj*s+i],0);CHKERRQ(ierr);
83813af1a74SHong Zhang         }
83913af1a74SHong Zhang         if (ts->vecs_sensi2p) { /* 2nd-order adjoint */
84013af1a74SHong Zhang           if (!zero) {
84113af1a74SHong Zhang             ierr = MatMultTranspose(ts->Jacprhs,VecsSensi2Temp[nadj],VecDeltaMu2);CHKERRQ(ierr);
84213af1a74SHong Zhang           } else {
84313af1a74SHong Zhang             ierr = VecSet(VecDeltaMu2,0);CHKERRQ(ierr);
84413af1a74SHong Zhang           }
84513af1a74SHong Zhang           if (ts->vecs_gpu) {
84613af1a74SHong Zhang             ierr = VecAXPY(VecDeltaMu2,1,ts->vecs_gpu[nadj]);CHKERRQ(ierr);
84713af1a74SHong Zhang           }
84813af1a74SHong Zhang           if (ts->vecs_gpp) {
84913af1a74SHong Zhang             ierr = VecAXPY(VecDeltaMu2,1,ts->vecs_gpp[nadj]);CHKERRQ(ierr);
85013af1a74SHong Zhang           }
85113af1a74SHong Zhang           ierr = VecAXPY(ts->vecs_sensi2p[nadj],1,VecDeltaMu2);CHKERRQ(ierr); /* update sensi2p for each stage */
85213af1a74SHong Zhang         }
85313af1a74SHong Zhang       }
85413af1a74SHong Zhang     }
8556497ce14SHong Zhang   }
856c235aa8dSHong Zhang 
857c235aa8dSHong Zhang   for (j=0; j<s; j++) w[j] = 1.0;
8582e7b7f96SHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) { /* no need to do this for mu's */
8592e7b7f96SHong Zhang     ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecsDeltaLam[nadj*s]);CHKERRQ(ierr);
86013af1a74SHong Zhang     if (ts->vecs_sensi2) {
86113af1a74SHong Zhang       ierr = VecMAXPY(ts->vecs_sensi2[nadj],s,w,&VecsDeltaLam2[nadj*s]);CHKERRQ(ierr);
86213af1a74SHong Zhang     }
8636497ce14SHong Zhang   }
864c235aa8dSHong Zhang   rk->status = TS_STEP_COMPLETE;
865d2daff3dSHong Zhang   PetscFunctionReturn(0);
866d2daff3dSHong Zhang }
867d2daff3dSHong Zhang 
86813af1a74SHong Zhang static PetscErrorCode TSAdjointReset_RK(TS ts)
86913af1a74SHong Zhang {
87013af1a74SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
87113af1a74SHong Zhang   RKTableau      tab = rk->tableau;
87213af1a74SHong Zhang   PetscErrorCode ierr;
87313af1a74SHong Zhang 
87413af1a74SHong Zhang   PetscFunctionBegin;
87513af1a74SHong Zhang   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam);CHKERRQ(ierr);
87613af1a74SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&rk->VecsSensiTemp);CHKERRQ(ierr);
87713af1a74SHong Zhang   ierr = VecDestroy(&rk->VecDeltaMu);CHKERRQ(ierr);
87813af1a74SHong Zhang   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam2);CHKERRQ(ierr);
87913af1a74SHong Zhang   ierr = VecDestroy(&rk->VecDeltaMu2);CHKERRQ(ierr);
88013af1a74SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&rk->VecsSensi2Temp);CHKERRQ(ierr);
88113af1a74SHong Zhang   PetscFunctionReturn(0);
88213af1a74SHong Zhang }
88313af1a74SHong Zhang 
88455de54fcSHong Zhang static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X)
88555de54fcSHong Zhang {
88655de54fcSHong Zhang   TS_RK            *rk = (TS_RK*)ts->data;
88755de54fcSHong Zhang   PetscInt         s  = rk->tableau->s,p = rk->tableau->p,i,j;
88855de54fcSHong Zhang   PetscReal        h;
88955de54fcSHong Zhang   PetscReal        tt,t;
89055de54fcSHong Zhang   PetscScalar      *b;
89155de54fcSHong Zhang   const PetscReal  *B = rk->tableau->binterp;
89255de54fcSHong Zhang   PetscErrorCode   ierr;
89355de54fcSHong Zhang 
89455de54fcSHong Zhang   PetscFunctionBegin;
89555de54fcSHong Zhang   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
89655de54fcSHong Zhang 
89755de54fcSHong Zhang   switch (rk->status) {
89855de54fcSHong Zhang     case TS_STEP_INCOMPLETE:
89955de54fcSHong Zhang     case TS_STEP_PENDING:
90055de54fcSHong Zhang       h = ts->time_step;
90155de54fcSHong Zhang       t = (itime - ts->ptime)/h;
90255de54fcSHong Zhang       break;
90355de54fcSHong Zhang     case TS_STEP_COMPLETE:
90455de54fcSHong Zhang       h = ts->ptime - ts->ptime_prev;
90555de54fcSHong Zhang       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
90655de54fcSHong Zhang       break;
90755de54fcSHong Zhang     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
90855de54fcSHong Zhang   }
90955de54fcSHong Zhang   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
91055de54fcSHong Zhang   for (i=0; i<s; i++) b[i] = 0;
91155de54fcSHong Zhang   for (j=0,tt=t; j<p; j++,tt*=t) {
91255de54fcSHong Zhang     for (i=0; i<s; i++) {
91355de54fcSHong Zhang       b[i]  += h * B[i*p+j] * tt;
91455de54fcSHong Zhang     }
91555de54fcSHong Zhang   }
91655de54fcSHong Zhang   ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr);
91755de54fcSHong Zhang   ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr);
91855de54fcSHong Zhang   ierr = PetscFree(b);CHKERRQ(ierr);
91955de54fcSHong Zhang   PetscFunctionReturn(0);
92055de54fcSHong Zhang }
92155de54fcSHong Zhang 
92255de54fcSHong Zhang /*------------------------------------------------------------*/
92355de54fcSHong Zhang 
924be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts)
925be5899b3SLisandro Dalcin {
926be5899b3SLisandro Dalcin   TS_RK          *rk = (TS_RK*)ts->data;
927be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
928be5899b3SLisandro Dalcin   PetscErrorCode ierr;
929be5899b3SLisandro Dalcin 
930be5899b3SLisandro Dalcin   PetscFunctionBegin;
931be5899b3SLisandro Dalcin   if (!tab) PetscFunctionReturn(0);
932be5899b3SLisandro Dalcin   ierr = PetscFree(rk->work);CHKERRQ(ierr);
933be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr);
934be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr);
9352e7b7f96SHong Zhang   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam);CHKERRQ(ierr);
9362e7b7f96SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&rk->VecsSensiTemp);CHKERRQ(ierr);
9372e7b7f96SHong Zhang   ierr = VecDestroy(&rk->VecDeltaMu);CHKERRQ(ierr);
938be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
939be5899b3SLisandro Dalcin }
940be5899b3SLisandro Dalcin 
941277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts)
942e4dd521cSBarry Smith {
9436849ba73SBarry Smith   PetscErrorCode ierr;
944e4dd521cSBarry Smith 
945e4dd521cSBarry Smith   PetscFunctionBegin;
946be5899b3SLisandro Dalcin   ierr = TSRKTableauReset(ts);CHKERRQ(ierr);
9470fe4d17eSHong Zhang   if (ts->use_splitrhsfunction) {
94802555c94SHong Zhang     ierr = PetscTryMethod(ts,"TSReset_RK_MultirateSplit_C",(TS),(ts));CHKERRQ(ierr);
9490fe4d17eSHong Zhang   } else {
95002555c94SHong Zhang     ierr = PetscTryMethod(ts,"TSReset_RK_MultirateNonsplit_C",(TS),(ts));CHKERRQ(ierr);
9510fe4d17eSHong Zhang   }
952277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
953e4dd521cSBarry Smith }
954277b19d0SLisandro Dalcin 
955f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx)
956f68a32c8SEmil Constantinescu {
957f68a32c8SEmil Constantinescu   PetscFunctionBegin;
958f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
959f68a32c8SEmil Constantinescu }
960f68a32c8SEmil Constantinescu 
961f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
962f68a32c8SEmil Constantinescu {
963f68a32c8SEmil Constantinescu   PetscFunctionBegin;
964f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
965f68a32c8SEmil Constantinescu }
966f68a32c8SEmil Constantinescu 
967f68a32c8SEmil Constantinescu 
968f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx)
969f68a32c8SEmil Constantinescu {
970f68a32c8SEmil Constantinescu   PetscFunctionBegin;
971f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
972f68a32c8SEmil Constantinescu }
973f68a32c8SEmil Constantinescu 
974f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
975f68a32c8SEmil Constantinescu {
976f68a32c8SEmil Constantinescu 
977f68a32c8SEmil Constantinescu   PetscFunctionBegin;
978f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
979f68a32c8SEmil Constantinescu }
980c235aa8dSHong Zhang /*
981d2daff3dSHong Zhang static PetscErrorCode RKSetAdjCoe(RKTableau tab)
982d2daff3dSHong Zhang {
983d2daff3dSHong Zhang   PetscReal *A,*b;
984d2daff3dSHong Zhang   PetscInt        s,i,j;
985d2daff3dSHong Zhang   PetscErrorCode  ierr;
986d2daff3dSHong Zhang 
987d2daff3dSHong Zhang   PetscFunctionBegin;
988d2daff3dSHong Zhang   s    = tab->s;
989d2daff3dSHong Zhang   ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr);
990d2daff3dSHong Zhang 
991d2daff3dSHong Zhang   for (i=0; i<s; i++)
992d2daff3dSHong Zhang     for (j=0; j<s; j++) {
993d2daff3dSHong Zhang       A[i*s+j] = (tab->b[s-1-i]==0) ? 0: -tab->A[s-1-i+(s-1-j)*s] * tab->b[s-1-j] / tab->b[s-1-i];
994d2daff3dSHong Zhang       ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr);
995d2daff3dSHong Zhang     }
996d2daff3dSHong Zhang 
997d2daff3dSHong Zhang   for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i];
998d2daff3dSHong Zhang 
999d2daff3dSHong Zhang   ierr  = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
1000d2daff3dSHong Zhang   ierr  = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
1001d2daff3dSHong Zhang   ierr  = PetscFree2(A,b);CHKERRQ(ierr);
1002d2daff3dSHong Zhang   PetscFunctionReturn(0);
1003d2daff3dSHong Zhang }
1004c235aa8dSHong Zhang */
1005be5899b3SLisandro Dalcin 
1006be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts)
1007be5899b3SLisandro Dalcin {
1008be5899b3SLisandro Dalcin   TS_RK          *rk  = (TS_RK*)ts->data;
1009be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
1010be5899b3SLisandro Dalcin   PetscErrorCode ierr;
1011be5899b3SLisandro Dalcin 
1012be5899b3SLisandro Dalcin   PetscFunctionBegin;
1013be5899b3SLisandro Dalcin   ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr);
1014be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr);
1015be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr);
1016be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
1017be5899b3SLisandro Dalcin }
1018be5899b3SLisandro Dalcin 
1019f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts)
1020f68a32c8SEmil Constantinescu {
1021*cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
1022f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
1023f68a32c8SEmil Constantinescu   DM             dm;
1024f68a32c8SEmil Constantinescu 
1025f68a32c8SEmil Constantinescu   PetscFunctionBegin;
10263dd83b38SBarry Smith   ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr);
1027be5899b3SLisandro Dalcin   ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);
1028*cd4cee2dSHong Zhang   if (quadts && ts->costintegralfwd) {
1029*cd4cee2dSHong Zhang     Mat Jquad;
1030*cd4cee2dSHong Zhang     ierr = TSGetRHSJacobian(quadts,&Jquad,NULL,NULL,NULL);CHKERRQ(ierr);
10310f7a1166SHong Zhang   }
1032f68a32c8SEmil Constantinescu   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1033f68a32c8SEmil Constantinescu   ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
1034f68a32c8SEmil Constantinescu   ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
10350fe4d17eSHong Zhang   if (ts->use_splitrhsfunction) {
103602555c94SHong Zhang     ierr = PetscTryMethod(ts,"TSSetUp_RK_MultirateSplit_C",(TS),(ts));CHKERRQ(ierr);
10370fe4d17eSHong Zhang   } else {
103802555c94SHong Zhang     ierr = PetscTryMethod(ts,"TSSetUp_RK_MultirateNonsplit_C",(TS),(ts));CHKERRQ(ierr);
10390fe4d17eSHong Zhang   }
1040e4dd521cSBarry Smith   PetscFunctionReturn(0);
1041e4dd521cSBarry Smith }
1042d2daff3dSHong Zhang 
10434416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts)
1044e4dd521cSBarry Smith {
1045be5899b3SLisandro Dalcin   TS_RK          *rk = (TS_RK*)ts->data;
1046dfbe8321SBarry Smith   PetscErrorCode ierr;
1047262ff7bbSBarry Smith 
1048e4dd521cSBarry Smith   PetscFunctionBegin;
1049e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr);
1050f68a32c8SEmil Constantinescu   {
1051f68a32c8SEmil Constantinescu     RKTableauLink link;
1052f68a32c8SEmil Constantinescu     PetscInt      count,choice;
10530fe4d17eSHong Zhang     PetscBool     flg,use_multirate = PETSC_FALSE;
1054f68a32c8SEmil Constantinescu     const char    **namelist;
10552c9cb062Sluzhanghpp 
1056f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) ;
1057f489ac74SBarry Smith     ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr);
1058f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
10590fe4d17eSHong Zhang     ierr = PetscOptionsBool("-ts_rk_multirate","Use interpolation-based multirate RK method","TSRKSetMultirate",rk->use_multirate,&use_multirate,&flg);CHKERRQ(ierr);
10600fe4d17eSHong Zhang     if (flg) {
10610fe4d17eSHong Zhang       ierr = TSRKSetMultirate(ts,use_multirate);CHKERRQ(ierr);
10620fe4d17eSHong Zhang     }
1063be5899b3SLisandro Dalcin     ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr);
1064be5899b3SLisandro Dalcin     if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);}
1065f68a32c8SEmil Constantinescu     ierr = PetscFree(namelist);CHKERRQ(ierr);
1066f68a32c8SEmil Constantinescu   }
1067262ff7bbSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
10682c9cb062Sluzhanghpp   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Multirate methods options","");CHKERRQ(ierr);
10692c9cb062Sluzhanghpp   ierr = PetscOptionsInt("-ts_rk_dtratio","time step ratio between slow and fast","",rk->dtratio,&rk->dtratio,NULL);CHKERRQ(ierr);
10702c9cb062Sluzhanghpp   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1071e4dd521cSBarry Smith   PetscFunctionReturn(0);
1072e4dd521cSBarry Smith }
1073e4dd521cSBarry Smith 
10745f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
1075e4dd521cSBarry Smith {
10765f70b478SJed Brown   TS_RK          *rk = (TS_RK*)ts->data;
10778370ee3bSLisandro Dalcin   PetscBool      iascii;
1078dfbe8321SBarry Smith   PetscErrorCode ierr;
10792cdf8120Sasbjorn 
1080e4dd521cSBarry Smith   PetscFunctionBegin;
1081251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
10828370ee3bSLisandro Dalcin   if (iascii) {
10839c334d8fSLisandro Dalcin     RKTableau tab  = rk->tableau;
1084f68a32c8SEmil Constantinescu     TSRKType  rktype;
1085f68a32c8SEmil Constantinescu     char      buf[512];
1086f68a32c8SEmil Constantinescu     ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr);
1087efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  RK type %s\n",rktype);CHKERRQ(ierr);
10880c37eb8eSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  Order: %D\n",tab->order);CHKERRQ(ierr);
10890c37eb8eSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr);
1090f68a32c8SEmil Constantinescu     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
1091f68a32c8SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa c = %s\n",buf);CHKERRQ(ierr);
10928370ee3bSLisandro Dalcin   }
1093f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1094f68a32c8SEmil Constantinescu }
1095f68a32c8SEmil Constantinescu 
1096f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer)
1097f68a32c8SEmil Constantinescu {
1098f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
10999c334d8fSLisandro Dalcin   TSAdapt        adapt;
1100f68a32c8SEmil Constantinescu 
1101f68a32c8SEmil Constantinescu   PetscFunctionBegin;
11029c334d8fSLisandro Dalcin   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
11039c334d8fSLisandro Dalcin   ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr);
1104f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1105f68a32c8SEmil Constantinescu }
1106f68a32c8SEmil Constantinescu 
1107f68a32c8SEmil Constantinescu /*@C
1108f68a32c8SEmil Constantinescu   TSRKSetType - Set the type of RK scheme
1109f68a32c8SEmil Constantinescu 
1110f68a32c8SEmil Constantinescu   Logically collective
1111f68a32c8SEmil Constantinescu 
1112f68a32c8SEmil Constantinescu   Input Parameter:
1113f68a32c8SEmil Constantinescu +  ts - timestepping context
1114f68a32c8SEmil Constantinescu -  rktype - type of RK-scheme
1115f68a32c8SEmil Constantinescu 
1116287c2655SBarry Smith   Options Database:
11179bd3e852SBarry Smith .   -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs>
1118287c2655SBarry Smith 
1119f68a32c8SEmil Constantinescu   Level: intermediate
1120f68a32c8SEmil Constantinescu 
1121287c2655SBarry Smith .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS
1122f68a32c8SEmil Constantinescu @*/
1123f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype)
1124f68a32c8SEmil Constantinescu {
1125f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
1126f68a32c8SEmil Constantinescu 
1127f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1128f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1129b92453a8SLisandro Dalcin   PetscValidCharPointer(rktype,2);
1130f68a32c8SEmil Constantinescu   ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr);
1131f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1132f68a32c8SEmil Constantinescu }
1133f68a32c8SEmil Constantinescu 
1134f68a32c8SEmil Constantinescu /*@C
1135f68a32c8SEmil Constantinescu   TSRKGetType - Get the type of RK scheme
1136f68a32c8SEmil Constantinescu 
1137f68a32c8SEmil Constantinescu   Logically collective
1138f68a32c8SEmil Constantinescu 
1139f68a32c8SEmil Constantinescu   Input Parameter:
1140f68a32c8SEmil Constantinescu .  ts - timestepping context
1141f68a32c8SEmil Constantinescu 
1142f68a32c8SEmil Constantinescu   Output Parameter:
1143f68a32c8SEmil Constantinescu .  rktype - type of RK-scheme
1144f68a32c8SEmil Constantinescu 
1145f68a32c8SEmil Constantinescu   Level: intermediate
1146f68a32c8SEmil Constantinescu 
1147f68a32c8SEmil Constantinescu .seealso: TSRKGetType()
1148f68a32c8SEmil Constantinescu @*/
1149f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype)
1150f68a32c8SEmil Constantinescu {
1151f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
1152f68a32c8SEmil Constantinescu 
1153f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1154f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1155f68a32c8SEmil Constantinescu   ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr);
1156f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1157f68a32c8SEmil Constantinescu }
1158f68a32c8SEmil Constantinescu 
1159560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype)
1160f68a32c8SEmil Constantinescu {
1161f68a32c8SEmil Constantinescu   TS_RK *rk = (TS_RK*)ts->data;
1162f68a32c8SEmil Constantinescu 
1163f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1164f68a32c8SEmil Constantinescu   *rktype = rk->tableau->name;
1165f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1166f68a32c8SEmil Constantinescu }
11672c9cb062Sluzhanghpp 
1168560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype)
1169f68a32c8SEmil Constantinescu {
1170f68a32c8SEmil Constantinescu   TS_RK          *rk = (TS_RK*)ts->data;
1171f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
1172f68a32c8SEmil Constantinescu   PetscBool      match;
1173f68a32c8SEmil Constantinescu   RKTableauLink  link;
1174f68a32c8SEmil Constantinescu 
1175f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1176f68a32c8SEmil Constantinescu   if (rk->tableau) {
1177f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr);
1178f68a32c8SEmil Constantinescu     if (match) PetscFunctionReturn(0);
1179f68a32c8SEmil Constantinescu   }
1180f68a32c8SEmil Constantinescu   for (link = RKTableauList; link; link=link->next) {
1181f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr);
1182f68a32c8SEmil Constantinescu     if (match) {
1183be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);}
1184f68a32c8SEmil Constantinescu       rk->tableau = &link->tab;
1185be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);}
11862ffb9264SLisandro Dalcin       ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
1187f68a32c8SEmil Constantinescu       PetscFunctionReturn(0);
1188f68a32c8SEmil Constantinescu     }
1189f68a32c8SEmil Constantinescu   }
1190f68a32c8SEmil Constantinescu   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype);
1191e4dd521cSBarry Smith   PetscFunctionReturn(0);
1192e4dd521cSBarry Smith }
1193e4dd521cSBarry Smith 
1194ff22ae23SHong Zhang static PetscErrorCode  TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y)
1195ff22ae23SHong Zhang {
1196ff22ae23SHong Zhang   TS_RK *rk = (TS_RK*)ts->data;
1197ff22ae23SHong Zhang 
1198ff22ae23SHong Zhang   PetscFunctionBegin;
11990429704eSStefano Zampini   if (ns) *ns = rk->tableau->s;
1200d2daff3dSHong Zhang   if (Y)   *Y = rk->Y;
1201ff22ae23SHong Zhang   PetscFunctionReturn(0);
1202ff22ae23SHong Zhang }
1203ff22ae23SHong Zhang 
1204b3a6b972SJed Brown static PetscErrorCode TSDestroy_RK(TS ts)
1205b3a6b972SJed Brown {
1206b3a6b972SJed Brown   PetscErrorCode ierr;
1207b3a6b972SJed Brown 
1208b3a6b972SJed Brown   PetscFunctionBegin;
1209b3a6b972SJed Brown   ierr = TSReset_RK(ts);CHKERRQ(ierr);
1210b3a6b972SJed Brown   if (ts->dm) {
1211b3a6b972SJed Brown     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
1212b3a6b972SJed Brown     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
1213b3a6b972SJed Brown   }
1214b3a6b972SJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
1215b3a6b972SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr);
1216b3a6b972SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr);
12170fe4d17eSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",NULL);CHKERRQ(ierr);
12180fe4d17eSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",NULL);CHKERRQ(ierr);
1219b3a6b972SJed Brown   PetscFunctionReturn(0);
1220b3a6b972SJed Brown }
1221ff22ae23SHong Zhang 
122221052c0cSHong Zhang /*@C
122321052c0cSHong Zhang   TSRKSetMultirate - Use the interpolation-based multirate RK method
122421052c0cSHong Zhang 
122521052c0cSHong Zhang   Logically collective
122621052c0cSHong Zhang 
122721052c0cSHong Zhang   Input Parameter:
122821052c0cSHong Zhang +  ts - timestepping context
122921052c0cSHong 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
123021052c0cSHong Zhang 
123121052c0cSHong Zhang   Options Database:
123221052c0cSHong Zhang .   -ts_rk_multirate - <true,false>
123321052c0cSHong Zhang 
123421052c0cSHong Zhang   Notes:
123521052c0cSHong 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).
123621052c0cSHong Zhang 
123721052c0cSHong Zhang   Level: intermediate
123821052c0cSHong Zhang 
123921052c0cSHong Zhang .seealso: TSRKGetMultirate()
124021052c0cSHong Zhang @*/
124121052c0cSHong Zhang PetscErrorCode TSRKSetMultirate(TS ts,PetscBool use_multirate)
124221052c0cSHong Zhang {
1243f06fb28fSHong Zhang   PetscErrorCode ierr;
12447dbacdf2SHong Zhang 
12458a4be4dbSHong Zhang   PetscFunctionBegin;
1246f06fb28fSHong Zhang   ierr = PetscTryMethod(ts,"TSRKSetMultirate_C",(TS,PetscBool),(ts,use_multirate));CHKERRQ(ierr);
124721052c0cSHong Zhang   PetscFunctionReturn(0);
124821052c0cSHong Zhang }
124921052c0cSHong Zhang 
125021052c0cSHong Zhang /*@C
125121052c0cSHong Zhang   TSRKGetMultirate - Gets whether to Use the interpolation-based multirate RK method
125221052c0cSHong Zhang 
125321052c0cSHong Zhang   Not collective
125421052c0cSHong Zhang 
125521052c0cSHong Zhang   Input Parameter:
125621052c0cSHong Zhang .  ts - timestepping context
125721052c0cSHong Zhang 
125821052c0cSHong Zhang   Output Parameter:
125921052c0cSHong Zhang .  use_multirate - PETSC_TRUE if the multirate RK method is enabled, PETSC_FALSE otherwise
126021052c0cSHong Zhang 
126121052c0cSHong Zhang   Level: intermediate
126221052c0cSHong Zhang 
126321052c0cSHong Zhang .seealso: TSRKSetMultirate()
126421052c0cSHong Zhang @*/
126521052c0cSHong Zhang PetscErrorCode TSRKGetMultirate(TS ts,PetscBool *use_multirate)
126621052c0cSHong Zhang {
1267f06fb28fSHong Zhang   PetscErrorCode ierr;
12687dbacdf2SHong Zhang 
12697dbacdf2SHong Zhang   PetscFunctionBegin;
1270f06fb28fSHong Zhang   ierr = PetscUseMethod(ts,"TSRKGetMultirate_C",(TS,PetscBool*),(ts,use_multirate));CHKERRQ(ierr);
127121052c0cSHong Zhang   PetscFunctionReturn(0);
127221052c0cSHong Zhang }
127321052c0cSHong Zhang 
1274ebe3b25bSBarry Smith /*MC
1275f68a32c8SEmil Constantinescu       TSRK - ODE and DAE solver using Runge-Kutta schemes
1276ebe3b25bSBarry Smith 
1277f68a32c8SEmil Constantinescu   The user should provide the right hand side of the equation
1278f68a32c8SEmil Constantinescu   using TSSetRHSFunction().
1279ebe3b25bSBarry Smith 
1280f68a32c8SEmil Constantinescu   Notes:
128198164e13SLisandro Dalcin   The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type
1282ebe3b25bSBarry Smith 
1283d41222bbSBarry Smith   Level: beginner
1284d41222bbSBarry Smith 
1285f68a32c8SEmil Constantinescu .seealso:  TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3,
12860fe4d17eSHong Zhang            TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister(), TSRKSetMultirate(), TSRKGetMultirate()
1287ebe3b25bSBarry Smith 
1288ebe3b25bSBarry Smith M*/
12898cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
1290e4dd521cSBarry Smith {
12911566a47fSLisandro Dalcin   TS_RK          *rk;
1292dfbe8321SBarry Smith   PetscErrorCode ierr;
1293e4dd521cSBarry Smith 
1294e4dd521cSBarry Smith   PetscFunctionBegin;
1295f68a32c8SEmil Constantinescu   ierr = TSRKInitializePackage();CHKERRQ(ierr);
1296f68a32c8SEmil Constantinescu 
1297f68a32c8SEmil Constantinescu   ts->ops->reset          = TSReset_RK;
12985f70b478SJed Brown   ts->ops->destroy        = TSDestroy_RK;
12995f70b478SJed Brown   ts->ops->view           = TSView_RK;
1300f68a32c8SEmil Constantinescu   ts->ops->load           = TSLoad_RK;
1301f68a32c8SEmil Constantinescu   ts->ops->setup          = TSSetUp_RK;
1302f68a32c8SEmil Constantinescu   ts->ops->interpolate    = TSInterpolate_RK;
13032c9cb062Sluzhanghpp   ts->ops->step           = TSStep_RK;
1304f68a32c8SEmil Constantinescu   ts->ops->evaluatestep   = TSEvaluateStep_RK;
1305fecfb714SLisandro Dalcin   ts->ops->rollback       = TSRollBack_RK;
1306f68a32c8SEmil Constantinescu   ts->ops->setfromoptions = TSSetFromOptions_RK;
1307ff22ae23SHong Zhang   ts->ops->getstages      = TSGetStages_RK;
1308e4dd521cSBarry Smith 
13090f7a1166SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_RK;
131013af1a74SHong Zhang   ts->ops->adjointsetup    = TSAdjointSetUp_RK;
131113af1a74SHong Zhang   ts->ops->adjointstep     = TSAdjointStep_RK;
131213af1a74SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_RK;
13130f7a1166SHong Zhang 
131413af1a74SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_RK;
1315922a638cSHong Zhang   ts->ops->forwardsetup    = TSForwardSetUp_RK;
1316922a638cSHong Zhang   ts->ops->forwardreset    = TSForwardReset_RK;
1317922a638cSHong Zhang   ts->ops->forwardstep     = TSForwardStep_RK;
1318922a638cSHong Zhang   ts->ops->forwardgetstages= TSForwardGetStages_RK;
1319922a638cSHong Zhang 
13201566a47fSLisandro Dalcin   ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr);
13211566a47fSLisandro Dalcin   ts->data = (void*)rk;
1322e4dd521cSBarry Smith 
1323f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr);
1324f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr);
132521052c0cSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",TSRKSetMultirate_RK);CHKERRQ(ierr);
132621052c0cSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",TSRKGetMultirate_RK);CHKERRQ(ierr);
1327be5899b3SLisandro Dalcin 
1328be5899b3SLisandro Dalcin   ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr);
132990b67129SHong Zhang   rk->dtratio = 1;
1330e4dd521cSBarry Smith   PetscFunctionReturn(0);
1331e4dd521cSBarry Smith }
1332