xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision b16ce868fadf09c63ae5be1e36b2479a88beb61d)
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 
27bcf0153eSBarry Smith      Options Database Key:
2867b8a455SSatish Balay .     -ts_rk_type 1fe - use type 1fe
29287c2655SBarry Smith 
30f68a32c8SEmil Constantinescu      Level: advanced
31f68a32c8SEmil Constantinescu 
321cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()`
33f68a32c8SEmil Constantinescu M*/
34f68a32c8SEmil Constantinescu /*MC
35e7685601SHong Zhang      TSRK2A - Second order RK scheme (Heun's method).
36f68a32c8SEmil Constantinescu 
37f68a32c8SEmil Constantinescu      This method has two stages.
38f68a32c8SEmil Constantinescu 
39bcf0153eSBarry Smith      Options Database Key:
4067b8a455SSatish Balay .     -ts_rk_type 2a - use type 2a
41287c2655SBarry Smith 
42f68a32c8SEmil Constantinescu      Level: advanced
43f68a32c8SEmil Constantinescu 
441cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()`
45f68a32c8SEmil Constantinescu M*/
46f68a32c8SEmil Constantinescu /*MC
47e7685601SHong Zhang      TSRK2B - Second order RK scheme (the midpoint method).
48e7685601SHong Zhang 
49e7685601SHong Zhang      This method has two stages.
50e7685601SHong Zhang 
51bcf0153eSBarry Smith      Options Database Key:
5267b8a455SSatish Balay .     -ts_rk_type 2b - use type 2b
53e7685601SHong Zhang 
54e7685601SHong Zhang      Level: advanced
55e7685601SHong Zhang 
561cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()`
57e7685601SHong Zhang M*/
58e7685601SHong Zhang /*MC
59f68a32c8SEmil Constantinescu      TSRK3 - Third order RK scheme.
60f68a32c8SEmil Constantinescu 
61f68a32c8SEmil Constantinescu      This method has three stages.
62f68a32c8SEmil Constantinescu 
63bcf0153eSBarry Smith      Options Database Key:
6467b8a455SSatish Balay .     -ts_rk_type 3 - use type 3
65287c2655SBarry Smith 
66f68a32c8SEmil Constantinescu      Level: advanced
67f68a32c8SEmil Constantinescu 
681cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()`
69f68a32c8SEmil Constantinescu M*/
70f68a32c8SEmil Constantinescu /*MC
712109b73fSDebojyoti Ghosh      TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method.
722109b73fSDebojyoti Ghosh 
73d1905564SLisandro Dalcin      This method has four stages with the First Same As Last (FSAL) property.
742109b73fSDebojyoti Ghosh 
75bcf0153eSBarry Smith      Options Database Key:
7667b8a455SSatish Balay .     -ts_rk_type 3bs - use type 3bs
77287c2655SBarry Smith 
782109b73fSDebojyoti Ghosh      Level: advanced
792109b73fSDebojyoti Ghosh 
80606c0280SSatish Balay      References:
81606c0280SSatish Balay . * - https://doi.org/10.1016/0893-9659(89)90079-7
82d1905564SLisandro Dalcin 
831cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()`
842109b73fSDebojyoti Ghosh M*/
852109b73fSDebojyoti Ghosh /*MC
86f68a32c8SEmil Constantinescu      TSRK4 - Fourth order RK scheme.
87f68a32c8SEmil Constantinescu 
882e698f8bSDebojyoti Ghosh      This is the classical Runge-Kutta method with four stages.
89f68a32c8SEmil Constantinescu 
90bcf0153eSBarry Smith      Options Database Key:
9167b8a455SSatish Balay .     -ts_rk_type 4 - use type 4
92287c2655SBarry Smith 
93f68a32c8SEmil Constantinescu      Level: advanced
94f68a32c8SEmil Constantinescu 
951cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()`
96f68a32c8SEmil Constantinescu M*/
97f68a32c8SEmil Constantinescu /*MC
982e698f8bSDebojyoti Ghosh      TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method.
99f68a32c8SEmil Constantinescu 
100f68a32c8SEmil Constantinescu      This method has six stages.
101f68a32c8SEmil Constantinescu 
102bcf0153eSBarry Smith      Options Database Key:
10367b8a455SSatish Balay .     -ts_rk_type 5f - use type 5f
104287c2655SBarry Smith 
105f68a32c8SEmil Constantinescu      Level: advanced
106f68a32c8SEmil Constantinescu 
1071cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()`
108f68a32c8SEmil Constantinescu M*/
1092109b73fSDebojyoti Ghosh /*MC
1102e698f8bSDebojyoti Ghosh      TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method.
1112109b73fSDebojyoti Ghosh 
112d1905564SLisandro Dalcin      This method has seven stages with the First Same As Last (FSAL) property.
1132109b73fSDebojyoti Ghosh 
114bcf0153eSBarry Smith      Options Database Key:
11567b8a455SSatish Balay .     -ts_rk_type 5dp - use type 5dp
116287c2655SBarry Smith 
1172109b73fSDebojyoti Ghosh      Level: advanced
1182109b73fSDebojyoti Ghosh 
119606c0280SSatish Balay      References:
120606c0280SSatish Balay . * - https://doi.org/10.1016/0771-050X(80)90013-3
121d1905564SLisandro Dalcin 
1221cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()`
1232109b73fSDebojyoti Ghosh M*/
12405e23783SLisandro Dalcin /*MC
12505e23783SLisandro Dalcin      TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method.
12605e23783SLisandro Dalcin 
12705e23783SLisandro Dalcin      This method has eight stages with the First Same As Last (FSAL) property.
12805e23783SLisandro Dalcin 
129bcf0153eSBarry Smith      Options Database Key:
13067b8a455SSatish Balay .     -ts_rk_type 5bs - use type 5bs
131287c2655SBarry Smith 
13205e23783SLisandro Dalcin      Level: advanced
13305e23783SLisandro Dalcin 
134606c0280SSatish Balay      References:
135606c0280SSatish Balay . * - https://doi.org/10.1016/0898-1221(96)00141-1
13605e23783SLisandro Dalcin 
1371cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()`
13805e23783SLisandro Dalcin M*/
13937acaa02SHendrik Ranocha /*MC
14037acaa02SHendrik Ranocha      TSRK6VR - Sixth order robust Verner RK scheme with fifth order embedded method.
14137acaa02SHendrik Ranocha 
14237acaa02SHendrik Ranocha      This method has nine stages with the First Same As Last (FSAL) property.
14337acaa02SHendrik Ranocha 
144bcf0153eSBarry Smith      Options Database Key:
14567b8a455SSatish Balay .     -ts_rk_type 6vr - use type 6vr
14637acaa02SHendrik Ranocha 
14737acaa02SHendrik Ranocha      Level: advanced
14837acaa02SHendrik Ranocha 
149606c0280SSatish Balay      References:
150606c0280SSatish Balay . * - http://people.math.sfu.ca/~jverner/RKV65.IIIXb.Robust.00010102836.081204.CoeffsOnlyRAT
15137acaa02SHendrik Ranocha 
1521cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()`
15337acaa02SHendrik Ranocha M*/
15437acaa02SHendrik Ranocha /*MC
15537acaa02SHendrik Ranocha      TSRK7VR - Seventh order robust Verner RK scheme with sixth order embedded method.
15637acaa02SHendrik Ranocha 
1578f3d7ee2SCarsten Uphoff      This method has ten stages.
15837acaa02SHendrik Ranocha 
159bcf0153eSBarry Smith      Options Database Key:
16067b8a455SSatish Balay .     -ts_rk_type 7vr - use type 7vr
16137acaa02SHendrik Ranocha 
16237acaa02SHendrik Ranocha      Level: advanced
16337acaa02SHendrik Ranocha 
164606c0280SSatish Balay      References:
165606c0280SSatish Balay . * - http://people.math.sfu.ca/~jverner/RKV76.IIa.Robust.000027015646.081206.CoeffsOnlyRAT
16637acaa02SHendrik Ranocha 
1671cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()`
16837acaa02SHendrik Ranocha M*/
16937acaa02SHendrik Ranocha /*MC
170d5b43468SJose E. Roman      TSRK8VR - Eighth order robust Verner RK scheme with seventh order embedded method.
17137acaa02SHendrik Ranocha 
1728f3d7ee2SCarsten Uphoff      This method has thirteen stages.
17337acaa02SHendrik Ranocha 
174bcf0153eSBarry Smith      Options Database Key:
17567b8a455SSatish Balay .     -ts_rk_type 8vr - use type 8vr
17637acaa02SHendrik Ranocha 
17737acaa02SHendrik Ranocha      Level: advanced
17837acaa02SHendrik Ranocha 
179606c0280SSatish Balay      References:
180606c0280SSatish Balay . * - http://people.math.sfu.ca/~jverner/RKV87.IIa.Robust.00000754677.081208.CoeffsOnlyRATandFLOAT
18137acaa02SHendrik Ranocha 
1821cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()`
18337acaa02SHendrik Ranocha M*/
184262ff7bbSBarry Smith 
185f68a32c8SEmil Constantinescu /*@C
186bcf0153eSBarry Smith   TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in `TSRK`
187262ff7bbSBarry Smith 
188f68a32c8SEmil Constantinescu   Not Collective, but should be called by all processes which will need the schemes to be registered
189262ff7bbSBarry Smith 
190f68a32c8SEmil Constantinescu   Level: advanced
191262ff7bbSBarry Smith 
1921cc06b55SBarry Smith .seealso: [](ch_ts), `TSRKRegisterDestroy()`, `TSRKRegister()`
193262ff7bbSBarry Smith @*/
194d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKRegisterAll(void)
195d71ae5a4SJacob Faibussowitsch {
196262ff7bbSBarry Smith   PetscFunctionBegin;
1973ba16761SJacob Faibussowitsch   if (TSRKRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
198f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_TRUE;
199e4dd521cSBarry Smith 
2004e82626cSLisandro Dalcin #define RC PetscRealConstant
201e4dd521cSBarry Smith   {
202*b16ce868SStefano Zampini     const PetscReal A[1][1] = {{0}};
203*b16ce868SStefano Zampini     const PetscReal b[1]    = {RC(1.0)};
2049566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK1FE, 1, 1, &A[0][0], b, NULL, NULL, 0, NULL));
205e4dd521cSBarry Smith   }
206db046809SBarry Smith   {
207*b16ce868SStefano Zampini     const PetscReal A[2][2] = {
2089371c9d4SSatish Balay       {0,       0},
2099371c9d4SSatish Balay       {RC(1.0), 0}
210*b16ce868SStefano Zampini     };
211*b16ce868SStefano Zampini     const PetscReal b[2]      = {RC(0.5), RC(0.5)};
212*b16ce868SStefano Zampini     const PetscReal bembed[2] = {RC(1.0), 0};
2139566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK2A, 2, 2, &A[0][0], b, NULL, bembed, 0, NULL));
214db046809SBarry Smith   }
215f68a32c8SEmil Constantinescu   {
216*b16ce868SStefano Zampini     const PetscReal A[2][2] = {
2179371c9d4SSatish Balay       {0,       0},
2189371c9d4SSatish Balay       {RC(0.5), 0}
219*b16ce868SStefano Zampini     };
220*b16ce868SStefano Zampini     const PetscReal b[2] = {0, RC(1.0)};
2219566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK2B, 2, 2, &A[0][0], b, NULL, NULL, 0, NULL));
222e7685601SHong Zhang   }
223e7685601SHong Zhang   {
224*b16ce868SStefano Zampini     const PetscReal A[3][3] = {
2259371c9d4SSatish Balay       {0,                  0,       0},
2264e82626cSLisandro Dalcin       {RC(2.0) / RC(3.0),  0,       0},
2279371c9d4SSatish Balay       {RC(-1.0) / RC(3.0), RC(1.0), 0}
228*b16ce868SStefano Zampini     };
229*b16ce868SStefano Zampini     const PetscReal b[3] = {RC(0.25), RC(0.5), RC(0.25)};
2309566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK3, 3, 3, &A[0][0], b, NULL, NULL, 0, NULL));
231fdefd5e5SDebojyoti Ghosh   }
232fdefd5e5SDebojyoti Ghosh   {
233*b16ce868SStefano Zampini     const PetscReal A[4][4] = {
2349371c9d4SSatish Balay       {0,                 0,                 0,                 0},
2354e82626cSLisandro Dalcin       {RC(1.0) / RC(2.0), 0,                 0,                 0},
2364e82626cSLisandro Dalcin       {0,                 RC(3.0) / RC(4.0), 0,                 0},
2379371c9d4SSatish Balay       {RC(2.0) / RC(9.0), RC(1.0) / RC(3.0), RC(4.0) / RC(9.0), 0}
238*b16ce868SStefano Zampini     };
239*b16ce868SStefano Zampini     const PetscReal b[4]      = {RC(2.0) / RC(9.0), RC(1.0) / RC(3.0), RC(4.0) / RC(9.0), 0};
240*b16ce868SStefano Zampini     const PetscReal 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)};
2419566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK3BS, 3, 4, &A[0][0], b, NULL, bembed, 0, NULL));
242db046809SBarry Smith   }
243f68a32c8SEmil Constantinescu   {
244*b16ce868SStefano Zampini     const PetscReal A[4][4] = {
2459371c9d4SSatish Balay       {0,       0,       0,       0},
2464e82626cSLisandro Dalcin       {RC(0.5), 0,       0,       0},
2474e82626cSLisandro Dalcin       {0,       RC(0.5), 0,       0},
2489371c9d4SSatish Balay       {0,       0,       RC(1.0), 0}
249*b16ce868SStefano Zampini     };
250*b16ce868SStefano Zampini     const PetscReal 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)};
2519566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK4, 4, 4, &A[0][0], b, NULL, NULL, 0, NULL));
252db046809SBarry Smith   }
253f68a32c8SEmil Constantinescu   {
254*b16ce868SStefano Zampini     const PetscReal A[6][6] = {
2559371c9d4SSatish Balay       {0,                       0,                        0,                        0,                       0,                    0},
2564e82626cSLisandro Dalcin       {RC(0.25),                0,                        0,                        0,                       0,                    0},
2574e82626cSLisandro Dalcin       {RC(3.0) / RC(32.0),      RC(9.0) / RC(32.0),       0,                        0,                       0,                    0},
2584e82626cSLisandro Dalcin       {RC(1932.0) / RC(2197.0), RC(-7200.0) / RC(2197.0), RC(7296.0) / RC(2197.0),  0,                       0,                    0},
2594e82626cSLisandro Dalcin       {RC(439.0) / RC(216.0),   RC(-8.0),                 RC(3680.0) / RC(513.0),   RC(-845.0) / RC(4104.0), 0,                    0},
2609371c9d4SSatish Balay       {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}
261*b16ce868SStefano Zampini     };
262*b16ce868SStefano Zampini     const PetscReal 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)};
263*b16ce868SStefano Zampini     const PetscReal 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};
2649566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK5F, 5, 6, &A[0][0], b, NULL, bembed, 0, NULL));
265fdefd5e5SDebojyoti Ghosh   }
266fdefd5e5SDebojyoti Ghosh   {
267*b16ce868SStefano Zampini     const PetscReal A[7][7] = {
2689371c9d4SSatish Balay       {0,                        0,                         0,                        0,                      0,                         0,                   0},
2694e82626cSLisandro Dalcin       {RC(1.0) / RC(5.0),        0,                         0,                        0,                      0,                         0,                   0},
2704e82626cSLisandro Dalcin       {RC(3.0) / RC(40.0),       RC(9.0) / RC(40.0),        0,                        0,                      0,                         0,                   0},
2714e82626cSLisandro Dalcin       {RC(44.0) / RC(45.0),      RC(-56.0) / RC(15.0),      RC(32.0) / RC(9.0),       0,                      0,                         0,                   0},
2724e82626cSLisandro 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},
2734e82626cSLisandro 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},
2749371c9d4SSatish Balay       {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}
275*b16ce868SStefano Zampini     };
276*b16ce868SStefano Zampini     const PetscReal 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};
277*b16ce868SStefano Zampini     const PetscReal 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)};
278*b16ce868SStefano Zampini     const PetscReal binterp[7][5] = {
279*b16ce868SStefano Zampini       {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)       },
280*b16ce868SStefano Zampini       {0,       0,                                      0,                                   0,                                       0                                     },
281*b16ce868SStefano Zampini       {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)  },
282*b16ce868SStefano Zampini       {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)      },
283*b16ce868SStefano Zampini       {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)},
284*b16ce868SStefano Zampini       {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)      },
285*b16ce868SStefano Zampini       {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)       }
286*b16ce868SStefano Zampini     };
2879566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK5DP, 5, 7, &A[0][0], b, NULL, bembed, 5, binterp[0]));
288f68a32c8SEmil Constantinescu   }
28905e23783SLisandro Dalcin   {
290*b16ce868SStefano Zampini     const PetscReal A[8][8] = {
2919371c9d4SSatish Balay       {0,                           0,                          0,                              0,                            0,                          0,                           0,                        0},
2924e82626cSLisandro Dalcin       {RC(1.0) / RC(6.0),           0,                          0,                              0,                            0,                          0,                           0,                        0},
2934e82626cSLisandro Dalcin       {RC(2.0) / RC(27.0),          RC(4.0) / RC(27.0),         0,                              0,                            0,                          0,                           0,                        0},
2944e82626cSLisandro Dalcin       {RC(183.0) / RC(1372.0),      RC(-162.0) / RC(343.0),     RC(1053.0) / RC(1372.0),        0,                            0,                          0,                           0,                        0},
2954e82626cSLisandro 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},
2964e82626cSLisandro 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},
2974e82626cSLisandro 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},
2989371c9d4SSatish Balay       {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}
299*b16ce868SStefano Zampini     };
300*b16ce868SStefano Zampini     const PetscReal 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};
301*b16ce868SStefano Zampini     const PetscReal 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)};
3029566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK5BS, 5, 8, &A[0][0], b, NULL, bembed, 0, NULL));
30305e23783SLisandro Dalcin   }
30437acaa02SHendrik Ranocha   {
305*b16ce868SStefano Zampini     const PetscReal A[9][9] = {
3069371c9d4SSatish Balay       {0,                                                   0,                                                  0,                                                   0,                                                   0,                                                   0,                                                   0,                                                  0,                                                  0},
30763fe90e8SHendrik Ranocha       {RC(1.8000000000000000000000000000000000000000e-01),  0,                                                  0,                                                   0,                                                   0,                                                   0,                                                   0,                                                  0,                                                  0},
30863fe90e8SHendrik Ranocha       {RC(8.9506172839506172839506172839506172839506e-02),  RC(7.7160493827160493827160493827160493827160e-02), 0,                                                   0,                                                   0,                                                   0,                                                   0,                                                  0,                                                  0},
30963fe90e8SHendrik Ranocha       {RC(6.2500000000000000000000000000000000000000e-02),  0,                                                  RC(1.8750000000000000000000000000000000000000e-01),  0,                                                   0,                                                   0,                                                   0,                                                  0,                                                  0},
31063fe90e8SHendrik Ranocha       {RC(3.1651600000000000000000000000000000000000e-01),  0,                                                  RC(-1.0449480000000000000000000000000000000000e+00), RC(1.2584320000000000000000000000000000000000e+00),  0,                                                   0,                                                   0,                                                  0,                                                  0},
31163fe90e8SHendrik Ranocha       {RC(2.7232612736485626257225065566674305502508e-01),  0,                                                  RC(-8.2513360323886639676113360323886639676113e-01), RC(1.0480917678812415654520917678812415654521e+00),  RC(1.0471570799276856873679117969088177628396e-01),  0,                                                   0,                                                  0,                                                  0},
31263fe90e8SHendrik 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},
31363fe90e8SHendrik 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},
3149371c9d4SSatish Balay       {RC(7.6388888888888888888888888888888888888889e-02),  0,                                                  0,                                                   RC(3.6940836940836940836940836940836940836941e-01),  0,                                                   RC(2.4801587301587301587301587301587301587302e-01),  RC(2.3674242424242424242424242424242424242424e-01), RC(6.9444444444444444444444444444444444444444e-02), 0}
315*b16ce868SStefano Zampini     };
316*b16ce868SStefano Zampini     const PetscReal b[9] = {RC(7.6388888888888888888888888888888888888889e-02), 0, 0, RC(3.6940836940836940836940836940836940836941e-01), 0, RC(2.4801587301587301587301587301587301587302e-01), RC(2.3674242424242424242424242424242424242424e-01),
317*b16ce868SStefano Zampini                             RC(6.9444444444444444444444444444444444444444e-02), 0};
318*b16ce868SStefano Zampini     const PetscReal 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)};
3199566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK6VR, 6, 9, &A[0][0], b, NULL, bembed, 0, NULL));
32037acaa02SHendrik Ranocha   }
32137acaa02SHendrik Ranocha   {
322*b16ce868SStefano Zampini     const PetscReal A[10][10] = {
3239371c9d4SSatish Balay       {0,                                                   0,                                                  0,                                                   0,                                                   0,                                                   0,                                                   0,                                                  0,                                                  0, 0},
32463fe90e8SHendrik Ranocha       {RC(5.0000000000000000000000000000000000000000e-03),  0,                                                  0,                                                   0,                                                   0,                                                   0,                                                   0,                                                  0,                                                  0, 0},
32563fe90e8SHendrik Ranocha       {RC(-1.0767901234567901234567901234567901234568e+00), RC(1.1856790123456790123456790123456790123457e+00), 0,                                                   0,                                                   0,                                                   0,                                                   0,                                                  0,                                                  0, 0},
32663fe90e8SHendrik Ranocha       {RC(4.0833333333333333333333333333333333333333e-02),  0,                                                  RC(1.2250000000000000000000000000000000000000e-01),  0,                                                   0,                                                   0,                                                   0,                                                  0,                                                  0, 0},
32763fe90e8SHendrik Ranocha       {RC(6.3607142857142857142857142857142857142857e-01),  0,                                                  RC(-2.4444642857142857142857142857142857142857e+00), RC(2.2633928571428571428571428571428571428571e+00),  0,                                                   0,                                                   0,                                                  0,                                                  0, 0},
32863fe90e8SHendrik Ranocha       {RC(-2.5351211079349245229256383554660215487207e+00), 0,                                                  RC(1.0299374654449267920438514460756024913612e+01),  RC(-7.9513032885990579949493217458266876536482e+00), RC(7.9301148923100592201226014271115261823800e-01),  0,                                                   0,                                                  0,                                                  0, 0},
32963fe90e8SHendrik 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},
33063fe90e8SHendrik 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},
33163fe90e8SHendrik 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},
3329371c9d4SSatish Balay       {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}
333*b16ce868SStefano Zampini     };
334*b16ce868SStefano Zampini     const PetscReal 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};
335*b16ce868SStefano Zampini     const PetscReal bembed[10] = {RC(4.7485247699299631037531273805727961552268e-02), 0, 0, RC(2.5599412588690633297154918245905393870497e-01), RC(2.7058478081067688722530891099268135732387e-01), RC(1.2505618684425992913638822323746917920448e-01),
3369371c9d4SSatish Balay                                   RC(2.5204468723743860507184043820197442562182e-01), 0, 0, RC(4.8834971521418614557381971303093137592592e-02)};
3379566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK7VR, 7, 10, &A[0][0], b, NULL, bembed, 0, NULL));
33837acaa02SHendrik Ranocha   }
33937acaa02SHendrik Ranocha   {
340*b16ce868SStefano Zampini     const PetscReal A[13][13] = {
3419371c9d4SSatish Balay       {0,                                                   0,                                                  0,                                                   0,                                                   0,                                                   0,                                                   0,                                                   0,                                                  0,                                                   0,                                                  0,                                                  0, 0},
34263fe90e8SHendrik Ranocha       {RC(2.5000000000000000000000000000000000000000e-01),  0,                                                  0,                                                   0,                                                   0,                                                   0,                                                   0,                                                   0,                                                  0,                                                   0,                                                  0,                                                  0, 0},
34363fe90e8SHendrik Ranocha       {RC(8.7400846504915232052686327594877411977046e-02),  RC(2.5487604938654321753087950620345685135815e-02), 0,                                                   0,                                                   0,                                                   0,                                                   0,                                                   0,                                                  0,                                                   0,                                                  0,                                                  0, 0},
34463fe90e8SHendrik Ranocha       {RC(4.2333169291338582677165354330708661417323e-02),  0,                                                  RC(1.2699950787401574803149606299212598425197e-01),  0,                                                   0,                                                   0,                                                   0,                                                   0,                                                  0,                                                   0,                                                  0,                                                  0, 0},
34563fe90e8SHendrik Ranocha       {RC(4.2609505888742261494881445237572274090942e-01),  0,                                                  RC(-1.5987952846591523265427733230657181117089e+00), RC(1.5967002257717297115939588706899953707994e+00),  0,                                                   0,                                                   0,                                                   0,                                                  0,                                                   0,                                                  0,                                                  0, 0},
34663fe90e8SHendrik Ranocha       {RC(5.0719337296713929515090618138513639239329e-02),  0,                                                  0,                                                   RC(2.5433377264600407582754714408877778031369e-01),  RC(2.0394689005728199465736223777270858044698e-01),  0,                                                   0,                                                   0,                                                  0,                                                   0,                                                  0,                                                  0, 0},
34763fe90e8SHendrik 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},
34863fe90e8SHendrik 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},
34963fe90e8SHendrik 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},
35063fe90e8SHendrik 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},
35163fe90e8SHendrik 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},
35263fe90e8SHendrik 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},
3539371c9d4SSatish Balay       {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}
354*b16ce868SStefano Zampini     };
355*b16ce868SStefano Zampini     const PetscReal 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};
356*b16ce868SStefano Zampini     const PetscReal 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)};
3579566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK8VR, 8, 13, &A[0][0], b, NULL, bembed, 0, NULL));
35837acaa02SHendrik Ranocha   }
3594e82626cSLisandro Dalcin #undef RC
3603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
361db046809SBarry Smith }
362db046809SBarry Smith 
363f68a32c8SEmil Constantinescu /*@C
364bcf0153eSBarry Smith   TSRKRegisterDestroy - Frees the list of schemes that were registered by `TSRKRegister()`.
365f68a32c8SEmil Constantinescu 
366f68a32c8SEmil Constantinescu   Not Collective
367f68a32c8SEmil Constantinescu 
368f68a32c8SEmil Constantinescu   Level: advanced
369f68a32c8SEmil Constantinescu 
3701cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKRegister()`, `TSRKRegisterAll()`
371f68a32c8SEmil Constantinescu @*/
372d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKRegisterDestroy(void)
373d71ae5a4SJacob Faibussowitsch {
374f68a32c8SEmil Constantinescu   RKTableauLink link;
375f68a32c8SEmil Constantinescu 
376f68a32c8SEmil Constantinescu   PetscFunctionBegin;
377f68a32c8SEmil Constantinescu   while ((link = RKTableauList)) {
378f68a32c8SEmil Constantinescu     RKTableau t   = &link->tab;
379f68a32c8SEmil Constantinescu     RKTableauList = link->next;
3809566063dSJacob Faibussowitsch     PetscCall(PetscFree3(t->A, t->b, t->c));
3819566063dSJacob Faibussowitsch     PetscCall(PetscFree(t->bembed));
3829566063dSJacob Faibussowitsch     PetscCall(PetscFree(t->binterp));
3839566063dSJacob Faibussowitsch     PetscCall(PetscFree(t->name));
3849566063dSJacob Faibussowitsch     PetscCall(PetscFree(link));
385f68a32c8SEmil Constantinescu   }
386f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_FALSE;
3873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
388f68a32c8SEmil Constantinescu }
389f68a32c8SEmil Constantinescu 
390f68a32c8SEmil Constantinescu /*@C
391bcf0153eSBarry Smith   TSRKInitializePackage - This function initializes everything in the `TSRK` package. It is called
392bcf0153eSBarry Smith   from `TSInitializePackage()`.
393f68a32c8SEmil Constantinescu 
394f68a32c8SEmil Constantinescu   Level: developer
395f68a32c8SEmil Constantinescu 
3961cc06b55SBarry Smith .seealso: [](ch_ts), `TSInitializePackage()`, `PetscInitialize()`, `TSRKFinalizePackage()`
397f68a32c8SEmil Constantinescu @*/
398d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKInitializePackage(void)
399d71ae5a4SJacob Faibussowitsch {
400e4dd521cSBarry Smith   PetscFunctionBegin;
4013ba16761SJacob Faibussowitsch   if (TSRKPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
402f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_TRUE;
4039566063dSJacob Faibussowitsch   PetscCall(TSRKRegisterAll());
4049566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(TSRKFinalizePackage));
4053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
406f68a32c8SEmil Constantinescu }
407186e87acSLisandro Dalcin 
408f68a32c8SEmil Constantinescu /*@C
409bcf0153eSBarry Smith   TSRKFinalizePackage - This function destroys everything in the `TSRK` package. It is
410bcf0153eSBarry Smith   called from `PetscFinalize()`.
411186e87acSLisandro Dalcin 
412f68a32c8SEmil Constantinescu   Level: developer
413f68a32c8SEmil Constantinescu 
4141cc06b55SBarry Smith .seealso: [](ch_ts), `PetscFinalize()`, `TSRKInitializePackage()`
415f68a32c8SEmil Constantinescu @*/
416d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKFinalizePackage(void)
417d71ae5a4SJacob Faibussowitsch {
418f68a32c8SEmil Constantinescu   PetscFunctionBegin;
419f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_FALSE;
4209566063dSJacob Faibussowitsch   PetscCall(TSRKRegisterDestroy());
4213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
422f68a32c8SEmil Constantinescu }
423f68a32c8SEmil Constantinescu 
424f68a32c8SEmil Constantinescu /*@C
425bcf0153eSBarry Smith   TSRKRegister - register an `TSRK` scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
426f68a32c8SEmil Constantinescu 
427f68a32c8SEmil Constantinescu   Not Collective, but the same schemes should be registered on all processes on which they will be used
428f68a32c8SEmil Constantinescu 
429f68a32c8SEmil Constantinescu   Input Parameters:
430f68a32c8SEmil Constantinescu + name    - identifier for method
431f68a32c8SEmil Constantinescu . order   - approximation order of method
432f68a32c8SEmil Constantinescu . s       - number of stages, this is the dimension of the matrices below
433f68a32c8SEmil Constantinescu . A       - stage coefficients (dimension s*s, row-major)
434f68a32c8SEmil Constantinescu . b       - step completion table (dimension s; NULL to use last row of A)
435f68a32c8SEmil Constantinescu . c       - abscissa (dimension s; NULL to use row sums of A)
436f68a32c8SEmil Constantinescu . bembed  - completion table for embedded method (dimension s; NULL if not available)
4373a8a9803SLisandro Dalcin . p       - Order of the interpolation scheme, equal to the number of columns of binterp
4383a8a9803SLisandro Dalcin - binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1)
439f68a32c8SEmil Constantinescu 
440f68a32c8SEmil Constantinescu   Level: advanced
441f68a32c8SEmil Constantinescu 
442bcf0153eSBarry Smith   Note:
443bcf0153eSBarry Smith   Several `TSRK` methods are provided, this function is only needed to create new methods.
444bcf0153eSBarry Smith 
4451cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`
446f68a32c8SEmil Constantinescu @*/
447d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKRegister(TSRKType name, PetscInt order, PetscInt s, const PetscReal A[], const PetscReal b[], const PetscReal c[], const PetscReal bembed[], PetscInt p, const PetscReal binterp[])
448d71ae5a4SJacob Faibussowitsch {
449f68a32c8SEmil Constantinescu   RKTableauLink link;
450f68a32c8SEmil Constantinescu   RKTableau     t;
451f68a32c8SEmil Constantinescu   PetscInt      i, j;
452f68a32c8SEmil Constantinescu 
453f68a32c8SEmil Constantinescu   PetscFunctionBegin;
4544f572ea9SToby Isaac   PetscAssertPointer(name, 1);
4554f572ea9SToby Isaac   PetscAssertPointer(A, 4);
4564f572ea9SToby Isaac   if (b) PetscAssertPointer(b, 5);
4574f572ea9SToby Isaac   if (c) PetscAssertPointer(c, 6);
4584f572ea9SToby Isaac   if (bembed) PetscAssertPointer(bembed, 7);
4594f572ea9SToby Isaac   if (binterp || p > 1) PetscAssertPointer(binterp, 9);
460095c51faSJed Brown   PetscCheck(s >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected number of stages s %" PetscInt_FMT " >= 0", s);
4613a8a9803SLisandro Dalcin 
4629566063dSJacob Faibussowitsch   PetscCall(TSRKInitializePackage());
4639566063dSJacob Faibussowitsch   PetscCall(PetscNew(&link));
464f68a32c8SEmil Constantinescu   t = &link->tab;
4653a8a9803SLisandro Dalcin 
4669566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &t->name));
467f68a32c8SEmil Constantinescu   t->order = order;
468f68a32c8SEmil Constantinescu   t->s     = s;
4699566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(s * s, &t->A, s, &t->b, s, &t->c));
4709566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->A, A, s * s));
4719566063dSJacob Faibussowitsch   if (b) PetscCall(PetscArraycpy(t->b, b, s));
4729371c9d4SSatish Balay   else
4739371c9d4SSatish Balay     for (i = 0; i < s; i++) t->b[i] = A[(s - 1) * s + i];
4749566063dSJacob Faibussowitsch   if (c) PetscCall(PetscArraycpy(t->c, c, s));
4759371c9d4SSatish Balay   else
4769371c9d4SSatish Balay     for (i = 0; i < s; i++)
4779371c9d4SSatish Balay       for (j = 0, t->c[i] = 0; j < s; j++) t->c[i] += A[i * s + j];
478d760c35bSDebojyoti Ghosh   t->FSAL = PETSC_TRUE;
4799371c9d4SSatish Balay   for (i = 0; i < s; i++)
4809371c9d4SSatish Balay     if (t->A[(s - 1) * s + i] != t->b[i]) t->FSAL = PETSC_FALSE;
4813a8a9803SLisandro Dalcin 
482f68a32c8SEmil Constantinescu   if (bembed) {
4839566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(s, &t->bembed));
4849566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(t->bembed, bembed, s));
485f68a32c8SEmil Constantinescu   }
486f68a32c8SEmil Constantinescu 
4879371c9d4SSatish Balay   if (!binterp) {
4889371c9d4SSatish Balay     p       = 1;
4899371c9d4SSatish Balay     binterp = t->b;
4909371c9d4SSatish Balay   }
4913a8a9803SLisandro Dalcin   t->p = p;
4929566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s * p, &t->binterp));
4939566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->binterp, binterp, s * p));
4943a8a9803SLisandro Dalcin 
495f68a32c8SEmil Constantinescu   link->next    = RKTableauList;
496f68a32c8SEmil Constantinescu   RKTableauList = link;
4973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
498f68a32c8SEmil Constantinescu }
499f68a32c8SEmil Constantinescu 
500d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKGetTableau_RK(TS ts, PetscInt *s, const PetscReal **A, const PetscReal **b, const PetscReal **c, const PetscReal **bembed, PetscInt *p, const PetscReal **binterp, PetscBool *FSAL)
501d71ae5a4SJacob Faibussowitsch {
5020f7a28cdSStefano Zampini   TS_RK    *rk  = (TS_RK *)ts->data;
5030f7a28cdSStefano Zampini   RKTableau tab = rk->tableau;
5040f7a28cdSStefano Zampini 
5050f7a28cdSStefano Zampini   PetscFunctionBegin;
5060f7a28cdSStefano Zampini   if (s) *s = tab->s;
5070f7a28cdSStefano Zampini   if (A) *A = tab->A;
5080f7a28cdSStefano Zampini   if (b) *b = tab->b;
5090f7a28cdSStefano Zampini   if (c) *c = tab->c;
5100f7a28cdSStefano Zampini   if (bembed) *bembed = tab->bembed;
5110f7a28cdSStefano Zampini   if (p) *p = tab->p;
5120f7a28cdSStefano Zampini   if (binterp) *binterp = tab->binterp;
5130f7a28cdSStefano Zampini   if (FSAL) *FSAL = tab->FSAL;
5143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5150f7a28cdSStefano Zampini }
5160f7a28cdSStefano Zampini 
5170f7a28cdSStefano Zampini /*@C
518bcf0153eSBarry Smith   TSRKGetTableau - Get info on the `TSRK` tableau
5190f7a28cdSStefano Zampini 
5200f7a28cdSStefano Zampini   Not Collective
5210f7a28cdSStefano Zampini 
522f899ff85SJose E. Roman   Input Parameter:
5230f7a28cdSStefano Zampini . ts - timestepping context
5240f7a28cdSStefano Zampini 
5250f7a28cdSStefano Zampini   Output Parameters:
5260f7a28cdSStefano Zampini + s       - number of stages, this is the dimension of the matrices below
5270f7a28cdSStefano Zampini . A       - stage coefficients (dimension s*s, row-major)
5280f7a28cdSStefano Zampini . b       - step completion table (dimension s)
5290f7a28cdSStefano Zampini . c       - abscissa (dimension s)
5300f7a28cdSStefano Zampini . bembed  - completion table for embedded method (dimension s; NULL if not available)
5310f7a28cdSStefano Zampini . p       - Order of the interpolation scheme, equal to the number of columns of binterp
5320f7a28cdSStefano Zampini . binterp - Coefficients of the interpolation formula (dimension s*p)
533aaa8cc7dSPierre Jolivet - FSAL    - whether or not the scheme has the First Same As Last property
5340f7a28cdSStefano Zampini 
5350f7a28cdSStefano Zampini   Level: developer
5360f7a28cdSStefano Zampini 
5371cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKRegister()`, `TSRKSetType()`
5380f7a28cdSStefano Zampini @*/
539d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKGetTableau(TS ts, PetscInt *s, const PetscReal **A, const PetscReal **b, const PetscReal **c, const PetscReal **bembed, PetscInt *p, const PetscReal **binterp, PetscBool *FSAL)
540d71ae5a4SJacob Faibussowitsch {
5410f7a28cdSStefano Zampini   PetscFunctionBegin;
5420f7a28cdSStefano Zampini   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5439371c9d4SSatish Balay   PetscUseMethod(ts, "TSRKGetTableau_C", (TS, PetscInt *, const PetscReal **, const PetscReal **, const PetscReal **, const PetscReal **, PetscInt *, const PetscReal **, PetscBool *), (ts, s, A, b, c, bembed, p, binterp, FSAL));
5443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5450f7a28cdSStefano Zampini }
5460f7a28cdSStefano Zampini 
547e4dd521cSBarry Smith /*
5482c9cb062Sluzhanghpp  This is for single-step RK method
549f68a32c8SEmil Constantinescu  The step completion formula is
550e4dd521cSBarry Smith 
551f68a32c8SEmil Constantinescu  x1 = x0 + h b^T YdotRHS
552f68a32c8SEmil Constantinescu 
553f68a32c8SEmil Constantinescu  This function can be called before or after ts->vec_sol has been updated.
554f68a32c8SEmil Constantinescu  Suppose we have a completion formula (b) and an embedded formula (be) of different order.
555f68a32c8SEmil Constantinescu  We can write
556f68a32c8SEmil Constantinescu 
557f68a32c8SEmil Constantinescu  x1e = x0 + h be^T YdotRHS
558f68a32c8SEmil Constantinescu      = x1 - h b^T YdotRHS + h be^T YdotRHS
559f68a32c8SEmil Constantinescu      = x1 + h (be - b)^T YdotRHS
560f68a32c8SEmil Constantinescu 
561f68a32c8SEmil Constantinescu  so we can evaluate the method with different order even after the step has been optimistically completed.
562e4dd521cSBarry Smith */
563d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateStep_RK(TS ts, PetscInt order, Vec X, PetscBool *done)
564d71ae5a4SJacob Faibussowitsch {
565f68a32c8SEmil Constantinescu   TS_RK       *rk  = (TS_RK *)ts->data;
566f68a32c8SEmil Constantinescu   RKTableau    tab = rk->tableau;
567f68a32c8SEmil Constantinescu   PetscScalar *w   = rk->work;
568f68a32c8SEmil Constantinescu   PetscReal    h;
569f68a32c8SEmil Constantinescu   PetscInt     s = tab->s, j;
570f68a32c8SEmil Constantinescu 
571f68a32c8SEmil Constantinescu   PetscFunctionBegin;
572f68a32c8SEmil Constantinescu   switch (rk->status) {
573f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
574d71ae5a4SJacob Faibussowitsch   case TS_STEP_PENDING:
575d71ae5a4SJacob Faibussowitsch     h = ts->time_step;
576d71ae5a4SJacob Faibussowitsch     break;
577d71ae5a4SJacob Faibussowitsch   case TS_STEP_COMPLETE:
578d71ae5a4SJacob Faibussowitsch     h = ts->ptime - ts->ptime_prev;
579d71ae5a4SJacob Faibussowitsch     break;
580d71ae5a4SJacob Faibussowitsch   default:
581d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus");
582f68a32c8SEmil Constantinescu   }
583f68a32c8SEmil Constantinescu   if (order == tab->order) {
584f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) {
5859566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol, X));
58690b67129SHong Zhang       for (j = 0; j < s; j++) w[j] = h * tab->b[j] / rk->dtratio;
5879566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(X, s, w, rk->YdotRHS));
5889566063dSJacob Faibussowitsch     } else PetscCall(VecCopy(ts->vec_sol, X));
5893ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
590f68a32c8SEmil Constantinescu   } else if (order == tab->order - 1) {
591f68a32c8SEmil Constantinescu     if (!tab->bembed) goto unavailable;
592f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) { /*Complete with the embedded method (be)*/
5939566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol, X));
594f68a32c8SEmil Constantinescu       for (j = 0; j < s; j++) w[j] = h * tab->bembed[j];
5959566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(X, s, w, rk->YdotRHS));
596f68a32c8SEmil Constantinescu     } else { /*Rollback and re-complete using (be-b) */
5979566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol, X));
598f68a32c8SEmil Constantinescu       for (j = 0; j < s; j++) w[j] = h * (tab->bembed[j] - tab->b[j]);
5999566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(X, s, w, rk->YdotRHS));
600f68a32c8SEmil Constantinescu     }
601f68a32c8SEmil Constantinescu     if (done) *done = PETSC_TRUE;
6023ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
603f68a32c8SEmil Constantinescu   }
604f68a32c8SEmil Constantinescu unavailable:
605f68a32c8SEmil Constantinescu   if (done) *done = PETSC_FALSE;
6069371c9d4SSatish Balay   else
6079371c9d4SSatish Balay     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "RK '%s' of order %" PetscInt_FMT " cannot evaluate step at order %" PetscInt_FMT ". Consider using -ts_adapt_type none or a different method that has an embedded estimate.", tab->name, tab->order, order);
6083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
609f68a32c8SEmil Constantinescu }
610f68a32c8SEmil Constantinescu 
611d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardCostIntegral_RK(TS ts)
612d71ae5a4SJacob Faibussowitsch {
6130f7a1166SHong Zhang   TS_RK           *rk     = (TS_RK *)ts->data;
614cd4cee2dSHong Zhang   TS               quadts = ts->quadraturets;
6150f7a1166SHong Zhang   RKTableau        tab    = rk->tableau;
6160f7a1166SHong Zhang   const PetscInt   s      = tab->s;
6170f7a1166SHong Zhang   const PetscReal *b = tab->b, *c = tab->c;
6180f7a1166SHong Zhang   Vec             *Y = rk->Y;
6190f7a1166SHong Zhang   PetscInt         i;
6200f7a1166SHong Zhang 
6210f7a1166SHong Zhang   PetscFunctionBegin;
622cd4cee2dSHong Zhang   /* No need to backup quadts->vec_sol since it can be reverted in TSRollBack_RK */
6230f7a1166SHong Zhang   for (i = s - 1; i >= 0; i--) {
624cd4cee2dSHong Zhang     /* Evolve quadrature TS solution to compute integrals */
6259566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(quadts, rk->ptime + rk->time_step * c[i], Y[i], ts->vec_costintegrand));
6269566063dSJacob Faibussowitsch     PetscCall(VecAXPY(quadts->vec_sol, rk->time_step * b[i], ts->vec_costintegrand));
6270f7a1166SHong Zhang   }
6283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6290f7a1166SHong Zhang }
6300f7a1166SHong Zhang 
631d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointCostIntegral_RK(TS ts)
632d71ae5a4SJacob Faibussowitsch {
6330f7a1166SHong Zhang   TS_RK           *rk     = (TS_RK *)ts->data;
6340f7a1166SHong Zhang   RKTableau        tab    = rk->tableau;
635cd4cee2dSHong Zhang   TS               quadts = ts->quadraturets;
6360f7a1166SHong Zhang   const PetscInt   s      = tab->s;
6370f7a1166SHong Zhang   const PetscReal *b = tab->b, *c = tab->c;
6380f7a1166SHong Zhang   Vec             *Y = rk->Y;
6390f7a1166SHong Zhang   PetscInt         i;
6400f7a1166SHong Zhang 
6410f7a1166SHong Zhang   PetscFunctionBegin;
6420f7a1166SHong Zhang   for (i = s - 1; i >= 0; i--) {
643cd4cee2dSHong Zhang     /* Evolve quadrature TS solution to compute integrals */
6449566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(quadts, ts->ptime + ts->time_step * (1.0 - c[i]), Y[i], ts->vec_costintegrand));
6459566063dSJacob Faibussowitsch     PetscCall(VecAXPY(quadts->vec_sol, -ts->time_step * b[i], ts->vec_costintegrand));
6460f7a1166SHong Zhang   }
6473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6480f7a1166SHong Zhang }
6490f7a1166SHong Zhang 
650d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRollBack_RK(TS ts)
651d71ae5a4SJacob Faibussowitsch {
652fecfb714SLisandro Dalcin   TS_RK           *rk     = (TS_RK *)ts->data;
653cd4cee2dSHong Zhang   TS               quadts = ts->quadraturets;
654fecfb714SLisandro Dalcin   RKTableau        tab    = rk->tableau;
655fecfb714SLisandro Dalcin   const PetscInt   s      = tab->s;
656cd4cee2dSHong Zhang   const PetscReal *b = tab->b, *c = tab->c;
657fecfb714SLisandro Dalcin   PetscScalar     *w = rk->work;
658cd4cee2dSHong Zhang   Vec             *Y = rk->Y, *YdotRHS = rk->YdotRHS;
659fecfb714SLisandro Dalcin   PetscInt         j;
660fecfb714SLisandro Dalcin   PetscReal        h;
661fecfb714SLisandro Dalcin 
662fecfb714SLisandro Dalcin   PetscFunctionBegin;
663fecfb714SLisandro Dalcin   switch (rk->status) {
664fecfb714SLisandro Dalcin   case TS_STEP_INCOMPLETE:
665d71ae5a4SJacob Faibussowitsch   case TS_STEP_PENDING:
666d71ae5a4SJacob Faibussowitsch     h = ts->time_step;
667d71ae5a4SJacob Faibussowitsch     break;
668d71ae5a4SJacob Faibussowitsch   case TS_STEP_COMPLETE:
669d71ae5a4SJacob Faibussowitsch     h = ts->ptime - ts->ptime_prev;
670d71ae5a4SJacob Faibussowitsch     break;
671d71ae5a4SJacob Faibussowitsch   default:
672d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus");
673fecfb714SLisandro Dalcin   }
674fecfb714SLisandro Dalcin   for (j = 0; j < s; j++) w[j] = -h * b[j];
6759566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(ts->vec_sol, s, w, YdotRHS));
676cd4cee2dSHong Zhang   if (quadts && ts->costintegralfwd) {
677cd4cee2dSHong Zhang     for (j = 0; j < s; j++) {
678cd4cee2dSHong Zhang       /* Revert the quadrature TS solution */
6799566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSFunction(quadts, rk->ptime + h * c[j], Y[j], ts->vec_costintegrand));
6809566063dSJacob Faibussowitsch       PetscCall(VecAXPY(quadts->vec_sol, -h * b[j], ts->vec_costintegrand));
681cd4cee2dSHong Zhang     }
682cd4cee2dSHong Zhang   }
6833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
684fecfb714SLisandro Dalcin }
685fecfb714SLisandro Dalcin 
686d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardStep_RK(TS ts)
687d71ae5a4SJacob Faibussowitsch {
688922a638cSHong Zhang   TS_RK           *rk  = (TS_RK *)ts->data;
689922a638cSHong Zhang   RKTableau        tab = rk->tableau;
690922a638cSHong Zhang   Mat              J, *MatsFwdSensipTemp = rk->MatsFwdSensipTemp;
691922a638cSHong Zhang   const PetscInt   s = tab->s;
692922a638cSHong Zhang   const PetscReal *A = tab->A, *c = tab->c, *b = tab->b;
693922a638cSHong Zhang   Vec             *Y = rk->Y;
694922a638cSHong Zhang   PetscInt         i, j;
695922a638cSHong Zhang   PetscReal        stage_time, h = ts->time_step;
696922a638cSHong Zhang   PetscBool        zero;
697922a638cSHong Zhang 
698922a638cSHong Zhang   PetscFunctionBegin;
6999566063dSJacob Faibussowitsch   PetscCall(MatCopy(ts->mat_sensip, rk->MatFwdSensip0, SAME_NONZERO_PATTERN));
7009566063dSJacob Faibussowitsch   PetscCall(TSGetRHSJacobian(ts, &J, NULL, NULL, NULL));
701922a638cSHong Zhang 
702922a638cSHong Zhang   for (i = 0; i < s; i++) {
703922a638cSHong Zhang     stage_time = ts->ptime + h * c[i];
704922a638cSHong Zhang     zero       = PETSC_FALSE;
705922a638cSHong Zhang     if (b[i] == 0 && i == s - 1) zero = PETSC_TRUE;
706922a638cSHong Zhang     /* TLM Stage values */
707922a638cSHong Zhang     if (!i) {
7089566063dSJacob Faibussowitsch       PetscCall(MatCopy(ts->mat_sensip, rk->MatsFwdStageSensip[i], SAME_NONZERO_PATTERN));
709922a638cSHong Zhang     } else if (!zero) {
7109566063dSJacob Faibussowitsch       PetscCall(MatZeroEntries(rk->MatsFwdStageSensip[i]));
71148a46eb9SPierre Jolivet       for (j = 0; j < i; j++) PetscCall(MatAXPY(rk->MatsFwdStageSensip[i], h * A[i * s + j], MatsFwdSensipTemp[j], SAME_NONZERO_PATTERN));
7129566063dSJacob Faibussowitsch       PetscCall(MatAXPY(rk->MatsFwdStageSensip[i], 1., ts->mat_sensip, SAME_NONZERO_PATTERN));
713922a638cSHong Zhang     } else {
7149566063dSJacob Faibussowitsch       PetscCall(MatZeroEntries(rk->MatsFwdStageSensip[i]));
715922a638cSHong Zhang     }
716922a638cSHong Zhang 
7179566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSJacobian(ts, stage_time, Y[i], J, J));
7189566063dSJacob Faibussowitsch     PetscCall(MatMatMult(J, rk->MatsFwdStageSensip[i], MAT_REUSE_MATRIX, PETSC_DEFAULT, &MatsFwdSensipTemp[i]));
71913af1a74SHong Zhang     if (ts->Jacprhs) {
7209566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(ts, stage_time, Y[i], ts->Jacprhs)); /* get f_p */
72113af1a74SHong Zhang       if (ts->vecs_sensi2p) {                                              /* TLM used for 2nd-order adjoint */
72213af1a74SHong Zhang         PetscScalar *xarr;
7239566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(MatsFwdSensipTemp[i], 0, &xarr));
7249566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(rk->VecDeltaFwdSensipCol, xarr));
7259566063dSJacob Faibussowitsch         PetscCall(MatMultAdd(ts->Jacprhs, ts->vec_dir, rk->VecDeltaFwdSensipCol, rk->VecDeltaFwdSensipCol));
7269566063dSJacob Faibussowitsch         PetscCall(VecResetArray(rk->VecDeltaFwdSensipCol));
7279566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(MatsFwdSensipTemp[i], &xarr));
72813af1a74SHong Zhang       } else {
7299566063dSJacob Faibussowitsch         PetscCall(MatAXPY(MatsFwdSensipTemp[i], 1., ts->Jacprhs, SUBSET_NONZERO_PATTERN));
73013af1a74SHong Zhang       }
731922a638cSHong Zhang     }
732922a638cSHong Zhang   }
733922a638cSHong Zhang 
73448a46eb9SPierre Jolivet   for (i = 0; i < s; i++) PetscCall(MatAXPY(ts->mat_sensip, h * b[i], rk->MatsFwdSensipTemp[i], SAME_NONZERO_PATTERN));
735922a638cSHong Zhang   rk->status = TS_STEP_COMPLETE;
7363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
737922a638cSHong Zhang }
738922a638cSHong Zhang 
739d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardGetStages_RK(TS ts, PetscInt *ns, Mat **stagesensip)
740d71ae5a4SJacob Faibussowitsch {
741922a638cSHong Zhang   TS_RK    *rk  = (TS_RK *)ts->data;
742922a638cSHong Zhang   RKTableau tab = rk->tableau;
743922a638cSHong Zhang 
744922a638cSHong Zhang   PetscFunctionBegin;
745922a638cSHong Zhang   if (ns) *ns = tab->s;
746922a638cSHong Zhang   if (stagesensip) *stagesensip = rk->MatsFwdStageSensip;
7473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
748922a638cSHong Zhang }
749922a638cSHong Zhang 
750d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardSetUp_RK(TS ts)
751d71ae5a4SJacob Faibussowitsch {
752922a638cSHong Zhang   TS_RK    *rk  = (TS_RK *)ts->data;
753922a638cSHong Zhang   RKTableau tab = rk->tableau;
754922a638cSHong Zhang   PetscInt  i;
755922a638cSHong Zhang 
756922a638cSHong Zhang   PetscFunctionBegin;
757922a638cSHong Zhang   /* backup sensitivity results for roll-backs */
7589566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &rk->MatFwdSensip0));
759922a638cSHong Zhang 
7609566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tab->s, &rk->MatsFwdStageSensip));
7619566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tab->s, &rk->MatsFwdSensipTemp));
762922a638cSHong Zhang   for (i = 0; i < tab->s; i++) {
7639566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &rk->MatsFwdStageSensip[i]));
7649566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &rk->MatsFwdSensipTemp[i]));
765922a638cSHong Zhang   }
7669566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &rk->VecDeltaFwdSensipCol));
7673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
768922a638cSHong Zhang }
769922a638cSHong Zhang 
770d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardReset_RK(TS ts)
771d71ae5a4SJacob Faibussowitsch {
772922a638cSHong Zhang   TS_RK    *rk  = (TS_RK *)ts->data;
773922a638cSHong Zhang   RKTableau tab = rk->tableau;
774922a638cSHong Zhang   PetscInt  i;
775922a638cSHong Zhang 
776922a638cSHong Zhang   PetscFunctionBegin;
7779566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&rk->MatFwdSensip0));
778922a638cSHong Zhang   if (rk->MatsFwdStageSensip) {
77948a46eb9SPierre Jolivet     for (i = 0; i < tab->s; i++) PetscCall(MatDestroy(&rk->MatsFwdStageSensip[i]));
7809566063dSJacob Faibussowitsch     PetscCall(PetscFree(rk->MatsFwdStageSensip));
781922a638cSHong Zhang   }
782922a638cSHong Zhang   if (rk->MatsFwdSensipTemp) {
78348a46eb9SPierre Jolivet     for (i = 0; i < tab->s; i++) PetscCall(MatDestroy(&rk->MatsFwdSensipTemp[i]));
7849566063dSJacob Faibussowitsch     PetscCall(PetscFree(rk->MatsFwdSensipTemp));
785922a638cSHong Zhang   }
7869566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&rk->VecDeltaFwdSensipCol));
7873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
788922a638cSHong Zhang }
789922a638cSHong Zhang 
790d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_RK(TS ts)
791d71ae5a4SJacob Faibussowitsch {
792f68a32c8SEmil Constantinescu   TS_RK           *rk  = (TS_RK *)ts->data;
793f68a32c8SEmil Constantinescu   RKTableau        tab = rk->tableau;
794f68a32c8SEmil Constantinescu   const PetscInt   s   = tab->s;
795fecfb714SLisandro Dalcin   const PetscReal *A = tab->A, *c = tab->c;
796f68a32c8SEmil Constantinescu   PetscScalar     *w = rk->work;
797f68a32c8SEmil Constantinescu   Vec             *Y = rk->Y, *YdotRHS = rk->YdotRHS;
798630f8c86SStefano Zampini   PetscBool        FSAL = (PetscBool)(tab->FSAL && !rk->newtableau);
799f68a32c8SEmil Constantinescu   TSAdapt          adapt;
800fecfb714SLisandro Dalcin   PetscInt         i, j;
801be5899b3SLisandro Dalcin   PetscInt         rejections = 0;
802be5899b3SLisandro Dalcin   PetscBool        stageok, accept = PETSC_TRUE;
803be5899b3SLisandro Dalcin   PetscReal        next_time_step = ts->time_step;
804f68a32c8SEmil Constantinescu 
805f68a32c8SEmil Constantinescu   PetscFunctionBegin;
806d1905564SLisandro Dalcin   if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE;
8079566063dSJacob Faibussowitsch   if (FSAL) PetscCall(VecCopy(YdotRHS[s - 1], YdotRHS[0]));
808630f8c86SStefano Zampini   rk->newtableau = PETSC_FALSE;
809d1905564SLisandro Dalcin 
810f68a32c8SEmil Constantinescu   rk->status = TS_STEP_INCOMPLETE;
811be5899b3SLisandro Dalcin   while (!ts->reason && rk->status != TS_STEP_COMPLETE) {
812be5899b3SLisandro Dalcin     PetscReal t = ts->ptime;
813f68a32c8SEmil Constantinescu     PetscReal h = ts->time_step;
814f68a32c8SEmil Constantinescu     for (i = 0; i < s; i++) {
8159be3e283SDebojyoti Ghosh       rk->stage_time = t + h * c[i];
8169566063dSJacob Faibussowitsch       PetscCall(TSPreStage(ts, rk->stage_time));
8179566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol, Y[i]));
818f68a32c8SEmil Constantinescu       for (j = 0; j < i; j++) w[j] = h * A[i * s + j];
8199566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(Y[i], i, w, YdotRHS));
8209566063dSJacob Faibussowitsch       PetscCall(TSPostStage(ts, rk->stage_time, i, Y));
8219566063dSJacob Faibussowitsch       PetscCall(TSGetAdapt(ts, &adapt));
8229566063dSJacob Faibussowitsch       PetscCall(TSAdaptCheckStage(adapt, ts, rk->stage_time, Y[i], &stageok));
823fecfb714SLisandro Dalcin       if (!stageok) goto reject_step;
8248f738a7cSLisandro Dalcin       if (FSAL && !i) continue;
8259566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSFunction(ts, t + h * c[i], Y[i], YdotRHS[i]));
826f68a32c8SEmil Constantinescu     }
827be5899b3SLisandro Dalcin 
828be5899b3SLisandro Dalcin     rk->status = TS_STEP_INCOMPLETE;
8299566063dSJacob Faibussowitsch     PetscCall(TSEvaluateStep(ts, tab->order, ts->vec_sol, NULL));
830f68a32c8SEmil Constantinescu     rk->status = TS_STEP_PENDING;
8319566063dSJacob Faibussowitsch     PetscCall(TSGetAdapt(ts, &adapt));
8329566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidatesClear(adapt));
8339566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidateAdd(adapt, tab->name, tab->order, 1, tab->ccfl, (PetscReal)tab->s, PETSC_TRUE));
8349566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
835be5899b3SLisandro Dalcin     rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
836be5899b3SLisandro Dalcin     if (!accept) { /* Roll back the current step */
8379566063dSJacob Faibussowitsch       PetscCall(TSRollBack_RK(ts));
838be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
839be5899b3SLisandro Dalcin       goto reject_step;
840be5899b3SLisandro Dalcin     }
841be5899b3SLisandro Dalcin 
8420f7a1166SHong Zhang     if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
8430f7a1166SHong Zhang       rk->ptime     = ts->ptime;
8440f7a1166SHong Zhang       rk->time_step = ts->time_step;
845493ed95fSHong Zhang     }
846be5899b3SLisandro Dalcin 
847f68a32c8SEmil Constantinescu     ts->ptime += ts->time_step;
848f68a32c8SEmil Constantinescu     ts->time_step = next_time_step;
849f68a32c8SEmil Constantinescu     break;
850be5899b3SLisandro Dalcin 
851be5899b3SLisandro Dalcin   reject_step:
8529371c9d4SSatish Balay     ts->reject++;
8539371c9d4SSatish Balay     accept = PETSC_FALSE;
854be5899b3SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
855be5899b3SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
85663a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
857f68a32c8SEmil Constantinescu     }
858f68a32c8SEmil Constantinescu   }
8593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
860e4dd521cSBarry Smith }
861e4dd521cSBarry Smith 
862d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointSetUp_RK(TS ts)
863d71ae5a4SJacob Faibussowitsch {
864f6a906c0SBarry Smith   TS_RK    *rk  = (TS_RK *)ts->data;
865be5899b3SLisandro Dalcin   RKTableau tab = rk->tableau;
866be5899b3SLisandro Dalcin   PetscInt  s   = tab->s;
867f6a906c0SBarry Smith 
868f6a906c0SBarry Smith   PetscFunctionBegin;
8693ba16761SJacob Faibussowitsch   if (ts->adjointsetupcalled++) PetscFunctionReturn(PETSC_SUCCESS);
8709566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], s * ts->numcost, &rk->VecsDeltaLam));
8719566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &rk->VecsSensiTemp));
87248a46eb9SPierre Jolivet   if (ts->vecs_sensip) PetscCall(VecDuplicate(ts->vecs_sensip[0], &rk->VecDeltaMu));
87313af1a74SHong Zhang   if (ts->vecs_sensi2) {
8749566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], s * ts->numcost, &rk->VecsDeltaLam2));
8759566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi2[0], ts->numcost, &rk->VecsSensi2Temp));
87613af1a74SHong Zhang   }
87748a46eb9SPierre Jolivet   if (ts->vecs_sensi2p) PetscCall(VecDuplicate(ts->vecs_sensi2p[0], &rk->VecDeltaMu2));
8783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
879f6a906c0SBarry Smith }
880f6a906c0SBarry Smith 
88135f5172aSHong Zhang /*
88235f5172aSHong Zhang   Assumptions:
88335f5172aSHong Zhang     - TSStep_RK() always evaluates the step with b, not bembed.
88435f5172aSHong Zhang */
885d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointStep_RK(TS ts)
886d71ae5a4SJacob Faibussowitsch {
887c235aa8dSHong Zhang   TS_RK           *rk     = (TS_RK *)ts->data;
888cd4cee2dSHong Zhang   TS               quadts = ts->quadraturets;
889c235aa8dSHong Zhang   RKTableau        tab    = rk->tableau;
8903ca0882eSHong Zhang   Mat              J, Jpre, Jquad;
891c235aa8dSHong Zhang   const PetscInt   s = tab->s;
892c235aa8dSHong Zhang   const PetscReal *A = tab->A, *b = tab->b, *c = tab->c;
89313af1a74SHong Zhang   PetscScalar     *w = rk->work, *xarr;
8942e7b7f96SHong Zhang   Vec             *Y = rk->Y, *VecsDeltaLam = rk->VecsDeltaLam, VecDeltaMu = rk->VecDeltaMu, *VecsSensiTemp = rk->VecsSensiTemp;
89513af1a74SHong Zhang   Vec             *VecsDeltaLam2 = rk->VecsDeltaLam2, VecDeltaMu2 = rk->VecDeltaMu2, *VecsSensi2Temp = rk->VecsSensi2Temp;
896cd4cee2dSHong Zhang   Vec              VecDRDUTransCol = ts->vec_drdu_col, VecDRDPTransCol = ts->vec_drdp_col;
897b18ea86cSHong Zhang   PetscInt         i, j, nadj;
8983d3eaba7SBarry Smith   PetscReal        t = ts->ptime;
8993ca0882eSHong Zhang   PetscReal        h = ts->time_step;
900c235aa8dSHong Zhang 
901d2daff3dSHong Zhang   PetscFunctionBegin;
902c235aa8dSHong Zhang   rk->status = TS_STEP_INCOMPLETE;
9033389c7d9SStefano Zampini 
9049566063dSJacob Faibussowitsch   PetscCall(TSGetRHSJacobian(ts, &J, &Jpre, NULL, NULL));
90548a46eb9SPierre Jolivet   if (quadts) PetscCall(TSGetRHSJacobian(quadts, &Jquad, NULL, NULL, NULL));
9062e7b7f96SHong Zhang   for (i = s - 1; i >= 0; i--) {
9076a1e4597SHong Zhang     if (tab->FSAL && i == s - 1) {
9086a1e4597SHong Zhang       /* VecsDeltaLam[nadj*s+s-1] are initialized with zeros and the values never change.*/
9096a1e4597SHong Zhang       continue;
9106a1e4597SHong Zhang     }
9113ca0882eSHong Zhang     rk->stage_time = t + h * (1.0 - c[i]);
9129566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts, Y[i], J, Jpre));
9139371c9d4SSatish Balay     if (quadts) { PetscCall(TSComputeRHSJacobian(quadts, rk->stage_time, Y[i], Jquad, Jquad)); /* get r_u^T */ }
9143389c7d9SStefano Zampini     if (ts->vecs_sensip) {
9159566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(ts, rk->stage_time, Y[i], ts->Jacprhs)); /* get f_p */
9169371c9d4SSatish Balay       if (quadts) { PetscCall(TSComputeRHSJacobianP(quadts, rk->stage_time, Y[i], quadts->Jacprhs)); /* get f_p for the quadrature */ }
91735f5172aSHong Zhang     }
9183389c7d9SStefano Zampini 
91935f5172aSHong Zhang     if (b[i]) {
92035f5172aSHong Zhang       for (j = i + 1; j < s; j++) w[j - i - 1] = A[j * s + i] / b[i]; /* coefficients for computing VecsSensiTemp */
92135f5172aSHong Zhang     } else {
92235f5172aSHong Zhang       for (j = i + 1; j < s; j++) w[j - i - 1] = A[j * s + i]; /* coefficients for computing VecsSensiTemp */
92335f5172aSHong Zhang     }
9242e7b7f96SHong Zhang 
9252e7b7f96SHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) {
9263389c7d9SStefano Zampini       /* Stage values of lambda */
92735f5172aSHong Zhang       if (b[i]) {
92835f5172aSHong Zhang         /* lambda_{n+1} + \sum_{j=i+1}^s a_{ji}/b[i]*lambda_{s,j} */
9299566063dSJacob Faibussowitsch         PetscCall(VecCopy(ts->vecs_sensi[nadj], VecsSensiTemp[nadj])); /* VecDeltaLam is an vec array of size s by numcost */
9309566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(VecsSensiTemp[nadj], s - i - 1, w, &VecsDeltaLam[nadj * s + i + 1]));
9319566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(J, VecsSensiTemp[nadj], VecsDeltaLam[nadj * s + i])); /* VecsSensiTemp will be reused by 2nd-order adjoint */
9329566063dSJacob Faibussowitsch         PetscCall(VecScale(VecsDeltaLam[nadj * s + i], -h * b[i]));
93335f5172aSHong Zhang         if (quadts) {
9349566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(Jquad, nadj, &xarr));
9359566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(VecDRDUTransCol, xarr));
9369566063dSJacob Faibussowitsch           PetscCall(VecAXPY(VecsDeltaLam[nadj * s + i], -h * b[i], VecDRDUTransCol));
9379566063dSJacob Faibussowitsch           PetscCall(VecResetArray(VecDRDUTransCol));
9389566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(Jquad, &xarr));
939cd4cee2dSHong Zhang         }
9403389c7d9SStefano Zampini       } else {
9416a1e4597SHong Zhang         /* \sum_{j=i+1}^s a_{ji}*lambda_{s,j} */
9429566063dSJacob Faibussowitsch         PetscCall(VecSet(VecsSensiTemp[nadj], 0));
9439566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(VecsSensiTemp[nadj], s - i - 1, w, &VecsDeltaLam[nadj * s + i + 1]));
9449566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(J, VecsSensiTemp[nadj], VecsDeltaLam[nadj * s + i]));
9459566063dSJacob Faibussowitsch         PetscCall(VecScale(VecsDeltaLam[nadj * s + i], -h));
9463389c7d9SStefano Zampini       }
9476497ce14SHong Zhang 
948ad8e2604SHong Zhang       /* Stage values of mu */
9496497ce14SHong Zhang       if (ts->vecs_sensip) {
95035f5172aSHong Zhang         if (b[i]) {
9519566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(ts->Jacprhs, VecsSensiTemp[nadj], VecDeltaMu));
9529566063dSJacob Faibussowitsch           PetscCall(VecScale(VecDeltaMu, -h * b[i]));
95335f5172aSHong Zhang           if (quadts) {
9549566063dSJacob Faibussowitsch             PetscCall(MatDenseGetColumn(quadts->Jacprhs, nadj, &xarr));
9559566063dSJacob Faibussowitsch             PetscCall(VecPlaceArray(VecDRDPTransCol, xarr));
9569566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecDeltaMu, -h * b[i], VecDRDPTransCol));
9579566063dSJacob Faibussowitsch             PetscCall(VecResetArray(VecDRDPTransCol));
9589566063dSJacob Faibussowitsch             PetscCall(MatDenseRestoreColumn(quadts->Jacprhs, &xarr));
959493ed95fSHong Zhang           }
96035f5172aSHong Zhang         } else {
9619566063dSJacob Faibussowitsch           PetscCall(VecScale(VecDeltaMu, -h));
96235f5172aSHong Zhang         }
9639566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj], 1., VecDeltaMu)); /* update sensip for each stage */
964ad8e2604SHong Zhang       }
965c235aa8dSHong Zhang     }
96613af1a74SHong Zhang 
96713af1a74SHong Zhang     if (ts->vecs_sensi2 && ts->forward_solve) { /* 2nd-order adjoint, TLM mode has to be turned on */
96813af1a74SHong Zhang       /* Get w1 at t_{n+1} from TLM matrix */
9699566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(rk->MatsFwdStageSensip[i], 0, &xarr));
9709566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
97113af1a74SHong Zhang       /* lambda_s^T F_UU w_1 */
9729566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSHessianProductFunctionUU(ts, rk->stage_time, Y[i], VecsSensiTemp, ts->vec_sensip_col, ts->vecs_guu));
97335f5172aSHong Zhang       if (quadts) {
97435f5172aSHong Zhang         /* R_UU w_1 */
9759566063dSJacob Faibussowitsch         PetscCall(TSComputeRHSHessianProductFunctionUU(quadts, rk->stage_time, Y[i], NULL, ts->vec_sensip_col, ts->vecs_guu));
97635f5172aSHong Zhang       }
97735f5172aSHong Zhang       if (ts->vecs_sensip) {
97813af1a74SHong Zhang         /* lambda_s^T F_UP w_2 */
9799566063dSJacob Faibussowitsch         PetscCall(TSComputeRHSHessianProductFunctionUP(ts, rk->stage_time, Y[i], VecsSensiTemp, ts->vec_dir, ts->vecs_gup));
98035f5172aSHong Zhang         if (quadts) {
98135f5172aSHong Zhang           /* R_UP w_2 */
9829566063dSJacob Faibussowitsch           PetscCall(TSComputeRHSHessianProductFunctionUP(quadts, rk->stage_time, Y[i], NULL, ts->vec_sensip_col, ts->vecs_gup));
98335f5172aSHong Zhang         }
98435f5172aSHong Zhang       }
98513af1a74SHong Zhang       if (ts->vecs_sensi2p) {
98613af1a74SHong Zhang         /* lambda_s^T F_PU w_1 */
9879566063dSJacob Faibussowitsch         PetscCall(TSComputeRHSHessianProductFunctionPU(ts, rk->stage_time, Y[i], VecsSensiTemp, ts->vec_sensip_col, ts->vecs_gpu));
98835f5172aSHong Zhang         /* lambda_s^T F_PP w_2 */
9899566063dSJacob Faibussowitsch         PetscCall(TSComputeRHSHessianProductFunctionPP(ts, rk->stage_time, Y[i], VecsSensiTemp, ts->vec_dir, ts->vecs_gpp));
99035f5172aSHong Zhang         if (b[i] && quadts) {
99135f5172aSHong Zhang           /* R_PU w_1 */
9929566063dSJacob Faibussowitsch           PetscCall(TSComputeRHSHessianProductFunctionPU(quadts, rk->stage_time, Y[i], NULL, ts->vec_sensip_col, ts->vecs_gpu));
99335f5172aSHong Zhang           /* R_PP w_2 */
9949566063dSJacob Faibussowitsch           PetscCall(TSComputeRHSHessianProductFunctionPP(quadts, rk->stage_time, Y[i], NULL, ts->vec_dir, ts->vecs_gpp));
99535f5172aSHong Zhang         }
99613af1a74SHong Zhang       }
9979566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_sensip_col));
9989566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(rk->MatsFwdStageSensip[i], &xarr));
99913af1a74SHong Zhang 
100013af1a74SHong Zhang       for (nadj = 0; nadj < ts->numcost; nadj++) {
100113af1a74SHong Zhang         /* Stage values of lambda */
100235f5172aSHong Zhang         if (b[i]) {
100335f5172aSHong Zhang           /* J_i^T*(Lambda_{n+1}+\sum_{j=i+1}^s a_{ji}/b_i*Lambda_{s,j} */
10049566063dSJacob Faibussowitsch           PetscCall(VecCopy(ts->vecs_sensi2[nadj], VecsSensi2Temp[nadj]));
10059566063dSJacob Faibussowitsch           PetscCall(VecMAXPY(VecsSensi2Temp[nadj], s - i - 1, w, &VecsDeltaLam2[nadj * s + i + 1]));
10069566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(J, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj * s + i]));
10079566063dSJacob Faibussowitsch           PetscCall(VecScale(VecsDeltaLam2[nadj * s + i], -h * b[i]));
10089566063dSJacob Faibussowitsch           PetscCall(VecAXPY(VecsDeltaLam2[nadj * s + i], -h * b[i], ts->vecs_guu[nadj]));
100948a46eb9SPierre Jolivet           if (ts->vecs_sensip) PetscCall(VecAXPY(VecsDeltaLam2[nadj * s + i], -h * b[i], ts->vecs_gup[nadj]));
101013af1a74SHong Zhang         } else {
101135f5172aSHong Zhang           /* \sum_{j=i+1}^s a_{ji}*Lambda_{s,j} */
10129566063dSJacob Faibussowitsch           PetscCall(VecSet(VecsDeltaLam2[nadj * s + i], 0));
10139566063dSJacob Faibussowitsch           PetscCall(VecMAXPY(VecsSensi2Temp[nadj], s - i - 1, w, &VecsDeltaLam2[nadj * s + i + 1]));
10149566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(J, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj * s + i]));
10159566063dSJacob Faibussowitsch           PetscCall(VecScale(VecsDeltaLam2[nadj * s + i], -h));
10169566063dSJacob Faibussowitsch           PetscCall(VecAXPY(VecsDeltaLam2[nadj * s + i], -h, ts->vecs_guu[nadj]));
101748a46eb9SPierre Jolivet           if (ts->vecs_sensip) PetscCall(VecAXPY(VecsDeltaLam2[nadj * s + i], -h, ts->vecs_gup[nadj]));
101835f5172aSHong Zhang         }
101935f5172aSHong Zhang         if (ts->vecs_sensi2p) { /* 2nd-order adjoint for parameters */
10209566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(ts->Jacprhs, VecsSensi2Temp[nadj], VecDeltaMu2));
102135f5172aSHong Zhang           if (b[i]) {
10229566063dSJacob Faibussowitsch             PetscCall(VecScale(VecDeltaMu2, -h * b[i]));
10239566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecDeltaMu2, -h * b[i], ts->vecs_gpu[nadj]));
10249566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecDeltaMu2, -h * b[i], ts->vecs_gpp[nadj]));
102513af1a74SHong Zhang           } else {
10269566063dSJacob Faibussowitsch             PetscCall(VecScale(VecDeltaMu2, -h));
10279566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecDeltaMu2, -h, ts->vecs_gpu[nadj]));
10289566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecDeltaMu2, -h, ts->vecs_gpp[nadj]));
102913af1a74SHong Zhang           }
10309566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], 1, VecDeltaMu2)); /* update sensi2p for each stage */
103113af1a74SHong Zhang         }
103213af1a74SHong Zhang       }
103313af1a74SHong Zhang     }
10346497ce14SHong Zhang   }
1035c235aa8dSHong Zhang 
1036c235aa8dSHong Zhang   for (j = 0; j < s; j++) w[j] = 1.0;
10372e7b7f96SHong Zhang   for (nadj = 0; nadj < ts->numcost; nadj++) { /* no need to do this for mu's */
10389566063dSJacob Faibussowitsch     PetscCall(VecMAXPY(ts->vecs_sensi[nadj], s, w, &VecsDeltaLam[nadj * s]));
103948a46eb9SPierre Jolivet     if (ts->vecs_sensi2) PetscCall(VecMAXPY(ts->vecs_sensi2[nadj], s, w, &VecsDeltaLam2[nadj * s]));
10406497ce14SHong Zhang   }
1041c235aa8dSHong Zhang   rk->status = TS_STEP_COMPLETE;
10423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1043d2daff3dSHong Zhang }
1044d2daff3dSHong Zhang 
1045d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointReset_RK(TS ts)
1046d71ae5a4SJacob Faibussowitsch {
104713af1a74SHong Zhang   TS_RK    *rk  = (TS_RK *)ts->data;
104813af1a74SHong Zhang   RKTableau tab = rk->tableau;
104913af1a74SHong Zhang 
105013af1a74SHong Zhang   PetscFunctionBegin;
10519566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(tab->s * ts->numcost, &rk->VecsDeltaLam));
10529566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &rk->VecsSensiTemp));
10539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&rk->VecDeltaMu));
10549566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(tab->s * ts->numcost, &rk->VecsDeltaLam2));
10559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&rk->VecDeltaMu2));
10569566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &rk->VecsSensi2Temp));
10573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
105813af1a74SHong Zhang }
105913af1a74SHong Zhang 
1060d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_RK(TS ts, PetscReal itime, Vec X)
1061d71ae5a4SJacob Faibussowitsch {
106255de54fcSHong Zhang   TS_RK           *rk = (TS_RK *)ts->data;
106355de54fcSHong Zhang   PetscInt         s = rk->tableau->s, p = rk->tableau->p, i, j;
106455de54fcSHong Zhang   PetscReal        h;
106555de54fcSHong Zhang   PetscReal        tt, t;
106655de54fcSHong Zhang   PetscScalar     *b;
106755de54fcSHong Zhang   const PetscReal *B = rk->tableau->binterp;
106855de54fcSHong Zhang 
106955de54fcSHong Zhang   PetscFunctionBegin;
10703c633725SBarry Smith   PetscCheck(B, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSRK %s does not have an interpolation formula", rk->tableau->name);
107155de54fcSHong Zhang 
107255de54fcSHong Zhang   switch (rk->status) {
107355de54fcSHong Zhang   case TS_STEP_INCOMPLETE:
107455de54fcSHong Zhang   case TS_STEP_PENDING:
107555de54fcSHong Zhang     h = ts->time_step;
107655de54fcSHong Zhang     t = (itime - ts->ptime) / h;
107755de54fcSHong Zhang     break;
107855de54fcSHong Zhang   case TS_STEP_COMPLETE:
107955de54fcSHong Zhang     h = ts->ptime - ts->ptime_prev;
108055de54fcSHong Zhang     t = (itime - ts->ptime) / h + 1; /* In the interval [0,1] */
108155de54fcSHong Zhang     break;
1082d71ae5a4SJacob Faibussowitsch   default:
1083d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus");
108455de54fcSHong Zhang   }
10859566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s, &b));
108655de54fcSHong Zhang   for (i = 0; i < s; i++) b[i] = 0;
108755de54fcSHong Zhang   for (j = 0, tt = t; j < p; j++, tt *= t) {
1088ad540459SPierre Jolivet     for (i = 0; i < s; i++) b[i] += h * B[i * p + j] * tt;
108955de54fcSHong Zhang   }
10909566063dSJacob Faibussowitsch   PetscCall(VecCopy(rk->Y[0], X));
10919566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(X, s, b, rk->YdotRHS));
10929566063dSJacob Faibussowitsch   PetscCall(PetscFree(b));
10933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
109455de54fcSHong Zhang }
109555de54fcSHong Zhang 
109655de54fcSHong Zhang /*------------------------------------------------------------*/
109755de54fcSHong Zhang 
1098d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRKTableauReset(TS ts)
1099d71ae5a4SJacob Faibussowitsch {
1100be5899b3SLisandro Dalcin   TS_RK    *rk  = (TS_RK *)ts->data;
1101be5899b3SLisandro Dalcin   RKTableau tab = rk->tableau;
1102be5899b3SLisandro Dalcin 
1103be5899b3SLisandro Dalcin   PetscFunctionBegin;
11043ba16761SJacob Faibussowitsch   if (!tab) PetscFunctionReturn(PETSC_SUCCESS);
11059566063dSJacob Faibussowitsch   PetscCall(PetscFree(rk->work));
11069566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(tab->s, &rk->Y));
11079566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(tab->s, &rk->YdotRHS));
11083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1109be5899b3SLisandro Dalcin }
1110be5899b3SLisandro Dalcin 
1111d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_RK(TS ts)
1112d71ae5a4SJacob Faibussowitsch {
1113e4dd521cSBarry Smith   PetscFunctionBegin;
11149566063dSJacob Faibussowitsch   PetscCall(TSRKTableauReset(ts));
11150fe4d17eSHong Zhang   if (ts->use_splitrhsfunction) {
1116cac4c232SBarry Smith     PetscTryMethod(ts, "TSReset_RK_MultirateSplit_C", (TS), (ts));
11170fe4d17eSHong Zhang   } else {
1118cac4c232SBarry Smith     PetscTryMethod(ts, "TSReset_RK_MultirateNonsplit_C", (TS), (ts));
11190fe4d17eSHong Zhang   }
11203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1121e4dd521cSBarry Smith }
1122277b19d0SLisandro Dalcin 
1123d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCoarsenHook_TSRK(DM fine, DM coarse, void *ctx)
1124d71ae5a4SJacob Faibussowitsch {
1125f68a32c8SEmil Constantinescu   PetscFunctionBegin;
11263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1127f68a32c8SEmil Constantinescu }
1128f68a32c8SEmil Constantinescu 
1129d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRestrictHook_TSRK(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx)
1130d71ae5a4SJacob Faibussowitsch {
1131f68a32c8SEmil Constantinescu   PetscFunctionBegin;
11323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1133f68a32c8SEmil Constantinescu }
1134f68a32c8SEmil Constantinescu 
1135d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainHook_TSRK(DM dm, DM subdm, void *ctx)
1136d71ae5a4SJacob Faibussowitsch {
1137f68a32c8SEmil Constantinescu   PetscFunctionBegin;
11383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1139f68a32c8SEmil Constantinescu }
1140f68a32c8SEmil Constantinescu 
1141d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx)
1142d71ae5a4SJacob Faibussowitsch {
1143f68a32c8SEmil Constantinescu   PetscFunctionBegin;
11443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1145f68a32c8SEmil Constantinescu }
1146be5899b3SLisandro Dalcin 
1147d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRKTableauSetUp(TS ts)
1148d71ae5a4SJacob Faibussowitsch {
1149be5899b3SLisandro Dalcin   TS_RK    *rk  = (TS_RK *)ts->data;
1150be5899b3SLisandro Dalcin   RKTableau tab = rk->tableau;
1151be5899b3SLisandro Dalcin 
1152be5899b3SLisandro Dalcin   PetscFunctionBegin;
11539566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tab->s, &rk->work));
11549566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &rk->Y));
11559566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &rk->YdotRHS));
1156630f8c86SStefano Zampini   rk->newtableau = PETSC_TRUE;
11573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1158be5899b3SLisandro Dalcin }
1159be5899b3SLisandro Dalcin 
1160d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_RK(TS ts)
1161d71ae5a4SJacob Faibussowitsch {
1162cd4cee2dSHong Zhang   TS quadts = ts->quadraturets;
1163f68a32c8SEmil Constantinescu   DM dm;
1164f68a32c8SEmil Constantinescu 
1165f68a32c8SEmil Constantinescu   PetscFunctionBegin;
11669566063dSJacob Faibussowitsch   PetscCall(TSCheckImplicitTerm(ts));
11679566063dSJacob Faibussowitsch   PetscCall(TSRKTableauSetUp(ts));
1168cd4cee2dSHong Zhang   if (quadts && ts->costintegralfwd) {
1169cd4cee2dSHong Zhang     Mat Jquad;
11709566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobian(quadts, &Jquad, NULL, NULL, NULL));
11710f7a1166SHong Zhang   }
11729566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
11739566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSRK, DMRestrictHook_TSRK, ts));
11749566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSRK, DMSubDomainRestrictHook_TSRK, ts));
11750fe4d17eSHong Zhang   if (ts->use_splitrhsfunction) {
1176cac4c232SBarry Smith     PetscTryMethod(ts, "TSSetUp_RK_MultirateSplit_C", (TS), (ts));
11770fe4d17eSHong Zhang   } else {
1178cac4c232SBarry Smith     PetscTryMethod(ts, "TSSetUp_RK_MultirateNonsplit_C", (TS), (ts));
11790fe4d17eSHong Zhang   }
11803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1181e4dd521cSBarry Smith }
1182d2daff3dSHong Zhang 
1183d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_RK(TS ts, PetscOptionItems *PetscOptionsObject)
1184d71ae5a4SJacob Faibussowitsch {
1185be5899b3SLisandro Dalcin   TS_RK *rk = (TS_RK *)ts->data;
1186262ff7bbSBarry Smith 
1187e4dd521cSBarry Smith   PetscFunctionBegin;
1188d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "RK ODE solver options");
1189f68a32c8SEmil Constantinescu   {
1190f68a32c8SEmil Constantinescu     RKTableauLink link;
1191f68a32c8SEmil Constantinescu     PetscInt      count, choice;
11920fe4d17eSHong Zhang     PetscBool     flg, use_multirate = PETSC_FALSE;
1193f68a32c8SEmil Constantinescu     const char  **namelist;
11942c9cb062Sluzhanghpp 
11959371c9d4SSatish Balay     for (link = RKTableauList, count = 0; link; link = link->next, count++)
11969371c9d4SSatish Balay       ;
11979566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(count, (char ***)&namelist));
1198f68a32c8SEmil Constantinescu     for (link = RKTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name;
11999566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ts_rk_multirate", "Use interpolation-based multirate RK method", "TSRKSetMultirate", rk->use_multirate, &use_multirate, &flg));
12001baa6e33SBarry Smith     if (flg) PetscCall(TSRKSetMultirate(ts, use_multirate));
12019566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEList("-ts_rk_type", "Family of RK method", "TSRKSetType", (const char *const *)namelist, count, rk->tableau->name, &choice, &flg));
12029566063dSJacob Faibussowitsch     if (flg) PetscCall(TSRKSetType(ts, namelist[choice]));
12039566063dSJacob Faibussowitsch     PetscCall(PetscFree(namelist));
1204f68a32c8SEmil Constantinescu   }
1205d0609cedSBarry Smith   PetscOptionsHeadEnd();
1206d0609cedSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)ts), NULL, "Multirate methods options", "");
12079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-ts_rk_dtratio", "time step ratio between slow and fast", "", rk->dtratio, &rk->dtratio, NULL));
1208d0609cedSBarry Smith   PetscOptionsEnd();
12093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1210e4dd521cSBarry Smith }
1211e4dd521cSBarry Smith 
1212d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_RK(TS ts, PetscViewer viewer)
1213d71ae5a4SJacob Faibussowitsch {
12145f70b478SJed Brown   TS_RK    *rk = (TS_RK *)ts->data;
12158370ee3bSLisandro Dalcin   PetscBool iascii;
12162cdf8120Sasbjorn 
1217e4dd521cSBarry Smith   PetscFunctionBegin;
12189566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
12198370ee3bSLisandro Dalcin   if (iascii) {
12209c334d8fSLisandro Dalcin     RKTableau        tab = rk->tableau;
1221f68a32c8SEmil Constantinescu     TSRKType         rktype;
12220f7a28cdSStefano Zampini     const PetscReal *c;
12230f7a28cdSStefano Zampini     PetscInt         s;
1224f68a32c8SEmil Constantinescu     char             buf[512];
12250f7a28cdSStefano Zampini     PetscBool        FSAL;
12260f7a28cdSStefano Zampini 
12279566063dSJacob Faibussowitsch     PetscCall(TSRKGetType(ts, &rktype));
12289566063dSJacob Faibussowitsch     PetscCall(TSRKGetTableau(ts, &s, NULL, NULL, &c, NULL, NULL, NULL, &FSAL));
12299566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  RK type %s\n", rktype));
123063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Order: %" PetscInt_FMT "\n", tab->order));
12319566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  FSAL property: %s\n", FSAL ? "yes" : "no"));
12329566063dSJacob Faibussowitsch     PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", s, c));
12339566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Abscissa c = %s\n", buf));
12348370ee3bSLisandro Dalcin   }
12353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1236f68a32c8SEmil Constantinescu }
1237f68a32c8SEmil Constantinescu 
1238d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSLoad_RK(TS ts, PetscViewer viewer)
1239d71ae5a4SJacob Faibussowitsch {
12409c334d8fSLisandro Dalcin   TSAdapt adapt;
1241f68a32c8SEmil Constantinescu 
1242f68a32c8SEmil Constantinescu   PetscFunctionBegin;
12439566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &adapt));
12449566063dSJacob Faibussowitsch   PetscCall(TSAdaptLoad(adapt, viewer));
12453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1246f68a32c8SEmil Constantinescu }
1247f68a32c8SEmil Constantinescu 
12482ea87ba9SLisandro Dalcin /*@
1249bcf0153eSBarry Smith   TSRKGetOrder - Get the order of the `TSRK` scheme
12502ea87ba9SLisandro Dalcin 
125120f4b53cSBarry Smith   Not Collective
12522ea87ba9SLisandro Dalcin 
12532ea87ba9SLisandro Dalcin   Input Parameter:
12542ea87ba9SLisandro Dalcin . ts - timestepping context
12552ea87ba9SLisandro Dalcin 
12562ea87ba9SLisandro Dalcin   Output Parameter:
1257bcf0153eSBarry Smith . order - order of `TSRK` scheme
12582ea87ba9SLisandro Dalcin 
12592ea87ba9SLisandro Dalcin   Level: intermediate
12602ea87ba9SLisandro Dalcin 
12611cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKGetType()`
12622ea87ba9SLisandro Dalcin @*/
1263d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKGetOrder(TS ts, PetscInt *order)
1264d71ae5a4SJacob Faibussowitsch {
12652ea87ba9SLisandro Dalcin   PetscFunctionBegin;
12662ea87ba9SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
12674f572ea9SToby Isaac   PetscAssertPointer(order, 2);
1268cac4c232SBarry Smith   PetscUseMethod(ts, "TSRKGetOrder_C", (TS, PetscInt *), (ts, order));
12693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12702ea87ba9SLisandro Dalcin }
12712ea87ba9SLisandro Dalcin 
1272f68a32c8SEmil Constantinescu /*@C
1273bcf0153eSBarry Smith   TSRKSetType - Set the type of the `TSRK` scheme
1274f68a32c8SEmil Constantinescu 
127520f4b53cSBarry Smith   Logically Collective
1276f68a32c8SEmil Constantinescu 
1277d8d19677SJose E. Roman   Input Parameters:
1278f68a32c8SEmil Constantinescu + ts     - timestepping context
1279bcf0153eSBarry Smith - rktype - type of `TSRK` scheme
1280f68a32c8SEmil Constantinescu 
1281bcf0153eSBarry Smith   Options Database Key:
12829bd3e852SBarry Smith . -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs>
1283287c2655SBarry Smith 
1284f68a32c8SEmil Constantinescu   Level: intermediate
1285f68a32c8SEmil Constantinescu 
12861cc06b55SBarry Smith .seealso: [](ch_ts), `TSRKGetType()`, `TSRK`, `TSRKType`, `TSRK1FE`, `TSRK2A`, `TSRK2B`, `TSRK3`, `TSRK3BS`, `TSRK4`, `TSRK5F`, `TSRK5DP`, `TSRK5BS`, `TSRK6VR`, `TSRK7VR`, `TSRK8VR`
1287f68a32c8SEmil Constantinescu @*/
1288d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKSetType(TS ts, TSRKType rktype)
1289d71ae5a4SJacob Faibussowitsch {
1290f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1291f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
12924f572ea9SToby Isaac   PetscAssertPointer(rktype, 2);
1293cac4c232SBarry Smith   PetscTryMethod(ts, "TSRKSetType_C", (TS, TSRKType), (ts, rktype));
12943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1295f68a32c8SEmil Constantinescu }
1296f68a32c8SEmil Constantinescu 
1297f68a32c8SEmil Constantinescu /*@C
1298bcf0153eSBarry Smith   TSRKGetType - Get the type of `TSRK` scheme
1299f68a32c8SEmil Constantinescu 
130020f4b53cSBarry Smith   Not Collective
1301f68a32c8SEmil Constantinescu 
1302f68a32c8SEmil Constantinescu   Input Parameter:
1303f68a32c8SEmil Constantinescu . ts - timestepping context
1304f68a32c8SEmil Constantinescu 
1305f68a32c8SEmil Constantinescu   Output Parameter:
1306bcf0153eSBarry Smith . rktype - type of `TSRK`-scheme
1307f68a32c8SEmil Constantinescu 
1308f68a32c8SEmil Constantinescu   Level: intermediate
1309f68a32c8SEmil Constantinescu 
13101cc06b55SBarry Smith .seealso: [](ch_ts), `TSRKSetType()`
1311f68a32c8SEmil Constantinescu @*/
1312d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKGetType(TS ts, TSRKType *rktype)
1313d71ae5a4SJacob Faibussowitsch {
1314f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1315f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1316cac4c232SBarry Smith   PetscUseMethod(ts, "TSRKGetType_C", (TS, TSRKType *), (ts, rktype));
13173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1318f68a32c8SEmil Constantinescu }
1319f68a32c8SEmil Constantinescu 
1320d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRKGetOrder_RK(TS ts, PetscInt *order)
1321d71ae5a4SJacob Faibussowitsch {
13222ea87ba9SLisandro Dalcin   TS_RK *rk = (TS_RK *)ts->data;
13232ea87ba9SLisandro Dalcin 
13242ea87ba9SLisandro Dalcin   PetscFunctionBegin;
13252ea87ba9SLisandro Dalcin   *order = rk->tableau->order;
13263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13272ea87ba9SLisandro Dalcin }
13282ea87ba9SLisandro Dalcin 
1329d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRKGetType_RK(TS ts, TSRKType *rktype)
1330d71ae5a4SJacob Faibussowitsch {
1331f68a32c8SEmil Constantinescu   TS_RK *rk = (TS_RK *)ts->data;
1332f68a32c8SEmil Constantinescu 
1333f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1334f68a32c8SEmil Constantinescu   *rktype = rk->tableau->name;
13353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1336f68a32c8SEmil Constantinescu }
13372c9cb062Sluzhanghpp 
1338d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRKSetType_RK(TS ts, TSRKType rktype)
1339d71ae5a4SJacob Faibussowitsch {
1340f68a32c8SEmil Constantinescu   TS_RK        *rk = (TS_RK *)ts->data;
1341f68a32c8SEmil Constantinescu   PetscBool     match;
1342f68a32c8SEmil Constantinescu   RKTableauLink link;
1343f68a32c8SEmil Constantinescu 
1344f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1345f68a32c8SEmil Constantinescu   if (rk->tableau) {
13469566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(rk->tableau->name, rktype, &match));
13473ba16761SJacob Faibussowitsch     if (match) PetscFunctionReturn(PETSC_SUCCESS);
1348f68a32c8SEmil Constantinescu   }
1349f68a32c8SEmil Constantinescu   for (link = RKTableauList; link; link = link->next) {
13509566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(link->tab.name, rktype, &match));
1351f68a32c8SEmil Constantinescu     if (match) {
13529566063dSJacob Faibussowitsch       if (ts->setupcalled) PetscCall(TSRKTableauReset(ts));
1353f68a32c8SEmil Constantinescu       rk->tableau = &link->tab;
13549566063dSJacob Faibussowitsch       if (ts->setupcalled) PetscCall(TSRKTableauSetUp(ts));
13552ffb9264SLisandro Dalcin       ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
13563ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
1357f68a32c8SEmil Constantinescu     }
1358f68a32c8SEmil Constantinescu   }
135998921bdaSJacob Faibussowitsch   SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", rktype);
1360e4dd521cSBarry Smith }
1361e4dd521cSBarry Smith 
1362d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGetStages_RK(TS ts, PetscInt *ns, Vec **Y)
1363d71ae5a4SJacob Faibussowitsch {
1364ff22ae23SHong Zhang   TS_RK *rk = (TS_RK *)ts->data;
1365ff22ae23SHong Zhang 
1366ff22ae23SHong Zhang   PetscFunctionBegin;
13670429704eSStefano Zampini   if (ns) *ns = rk->tableau->s;
1368d2daff3dSHong Zhang   if (Y) *Y = rk->Y;
13693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1370ff22ae23SHong Zhang }
1371ff22ae23SHong Zhang 
1372d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_RK(TS ts)
1373d71ae5a4SJacob Faibussowitsch {
1374b3a6b972SJed Brown   PetscFunctionBegin;
13759566063dSJacob Faibussowitsch   PetscCall(TSReset_RK(ts));
1376b3a6b972SJed Brown   if (ts->dm) {
13779566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSRK, DMRestrictHook_TSRK, ts));
13789566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSRK, DMSubDomainRestrictHook_TSRK, ts));
1379b3a6b972SJed Brown   }
13809566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
13819566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetOrder_C", NULL));
13829566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetType_C", NULL));
13839566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKSetType_C", NULL));
13849566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetTableau_C", NULL));
13859566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKSetMultirate_C", NULL));
13869566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetMultirate_C", NULL));
13872e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSetUp_RK_MultirateSplit_C", NULL));
13882e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSReset_RK_MultirateSplit_C", NULL));
13892e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSetUp_RK_MultirateNonsplit_C", NULL));
13902e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSReset_RK_MultirateNonsplit_C", NULL));
13913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1392b3a6b972SJed Brown }
1393ff22ae23SHong Zhang 
13943ca0882eSHong Zhang /*
13953ca0882eSHong Zhang   This defines the nonlinear equation that is to be solved with SNES
13963ca0882eSHong Zhang   We do not need to solve the equation; we just use SNES to approximate the Jacobian
13973ca0882eSHong Zhang */
1398d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_RK(SNES snes, Vec x, Vec y, TS ts)
1399d71ae5a4SJacob Faibussowitsch {
14003ca0882eSHong Zhang   TS_RK *rk = (TS_RK *)ts->data;
14013ca0882eSHong Zhang   DM     dm, dmsave;
14023ca0882eSHong Zhang 
14033ca0882eSHong Zhang   PetscFunctionBegin;
14049566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
14053ca0882eSHong Zhang   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
14063ca0882eSHong Zhang   dmsave = ts->dm;
14073ca0882eSHong Zhang   ts->dm = dm;
14089566063dSJacob Faibussowitsch   PetscCall(TSComputeRHSFunction(ts, rk->stage_time, x, y));
14093ca0882eSHong Zhang   ts->dm = dmsave;
14103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14113ca0882eSHong Zhang }
14123ca0882eSHong Zhang 
1413d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_RK(SNES snes, Vec x, Mat A, Mat B, TS ts)
1414d71ae5a4SJacob Faibussowitsch {
14153ca0882eSHong Zhang   TS_RK *rk = (TS_RK *)ts->data;
14163ca0882eSHong Zhang   DM     dm, dmsave;
14173ca0882eSHong Zhang 
14183ca0882eSHong Zhang   PetscFunctionBegin;
14199566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
14203ca0882eSHong Zhang   dmsave = ts->dm;
14213ca0882eSHong Zhang   ts->dm = dm;
14229566063dSJacob Faibussowitsch   PetscCall(TSComputeRHSJacobian(ts, rk->stage_time, x, A, B));
14233ca0882eSHong Zhang   ts->dm = dmsave;
14243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14253ca0882eSHong Zhang }
14263ca0882eSHong Zhang 
142721052c0cSHong Zhang /*@C
1428bcf0153eSBarry Smith   TSRKSetMultirate - Use the interpolation-based multirate `TSRK` method
142921052c0cSHong Zhang 
143020f4b53cSBarry Smith   Logically Collective
143121052c0cSHong Zhang 
1432d8d19677SJose E. Roman   Input Parameters:
143321052c0cSHong Zhang + ts            - timestepping context
1434bcf0153eSBarry Smith - use_multirate - `PETSC_TRUE` enables the multirate `TSRK` method, sets the basic method to be RK2A and sets the ratio between slow stepsize and fast stepsize to be 2
143521052c0cSHong Zhang 
1436bcf0153eSBarry Smith   Options Database Key:
143721052c0cSHong Zhang . -ts_rk_multirate - <true,false>
143821052c0cSHong Zhang 
143921052c0cSHong Zhang   Level: intermediate
144021052c0cSHong Zhang 
1441bcf0153eSBarry Smith   Note:
1442da81f932SPierre Jolivet   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 coefficients (binterp).
1443bcf0153eSBarry Smith 
14441cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKGetMultirate()`
144521052c0cSHong Zhang @*/
1446d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKSetMultirate(TS ts, PetscBool use_multirate)
1447d71ae5a4SJacob Faibussowitsch {
14488a4be4dbSHong Zhang   PetscFunctionBegin;
1449cac4c232SBarry Smith   PetscTryMethod(ts, "TSRKSetMultirate_C", (TS, PetscBool), (ts, use_multirate));
14503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
145121052c0cSHong Zhang }
145221052c0cSHong Zhang 
145321052c0cSHong Zhang /*@C
1454bcf0153eSBarry Smith   TSRKGetMultirate - Gets whether to use the interpolation-based multirate `TSRK` method
145521052c0cSHong Zhang 
145620f4b53cSBarry Smith   Not Collective
145721052c0cSHong Zhang 
145821052c0cSHong Zhang   Input Parameter:
145921052c0cSHong Zhang . ts - timestepping context
146021052c0cSHong Zhang 
146121052c0cSHong Zhang   Output Parameter:
1462bcf0153eSBarry Smith . use_multirate - `PETSC_TRUE` if the multirate RK method is enabled, `PETSC_FALSE` otherwise
146321052c0cSHong Zhang 
146421052c0cSHong Zhang   Level: intermediate
146521052c0cSHong Zhang 
14661cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKSetMultirate()`
146721052c0cSHong Zhang @*/
1468d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKGetMultirate(TS ts, PetscBool *use_multirate)
1469d71ae5a4SJacob Faibussowitsch {
14707dbacdf2SHong Zhang   PetscFunctionBegin;
1471cac4c232SBarry Smith   PetscUseMethod(ts, "TSRKGetMultirate_C", (TS, PetscBool *), (ts, use_multirate));
14723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
147321052c0cSHong Zhang }
147421052c0cSHong Zhang 
1475ebe3b25bSBarry Smith /*MC
1476f68a32c8SEmil Constantinescu       TSRK - ODE and DAE solver using Runge-Kutta schemes
1477ebe3b25bSBarry Smith 
1478f68a32c8SEmil Constantinescu   The user should provide the right hand side of the equation
1479bcf0153eSBarry Smith   using `TSSetRHSFunction()`.
1480ebe3b25bSBarry Smith 
1481d41222bbSBarry Smith   Level: beginner
1482d41222bbSBarry Smith 
1483bcf0153eSBarry Smith   Notes:
1484bcf0153eSBarry Smith   The default is `TSRK3BS`, it can be changed with `TSRKSetType()` or -ts_rk_type
1485ebe3b25bSBarry Smith 
14861cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSRK`, `TSSetType()`, `TSRKSetType()`, `TSRKGetType()`, `TSRK2D`, `TSRK2E`, `TSRK3`,
1487bcf0153eSBarry Smith           `TSRK4`, `TSRK5`, `TSRKPRSSP2`, `TSRKBPR3`, `TSRKType`, `TSRKRegister()`, `TSRKSetMultirate()`, `TSRKGetMultirate()`, `TSType`
1488ebe3b25bSBarry Smith M*/
1489d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
1490d71ae5a4SJacob Faibussowitsch {
14911566a47fSLisandro Dalcin   TS_RK *rk;
1492e4dd521cSBarry Smith 
1493e4dd521cSBarry Smith   PetscFunctionBegin;
14949566063dSJacob Faibussowitsch   PetscCall(TSRKInitializePackage());
1495f68a32c8SEmil Constantinescu 
1496f68a32c8SEmil Constantinescu   ts->ops->reset          = TSReset_RK;
14975f70b478SJed Brown   ts->ops->destroy        = TSDestroy_RK;
14985f70b478SJed Brown   ts->ops->view           = TSView_RK;
1499f68a32c8SEmil Constantinescu   ts->ops->load           = TSLoad_RK;
1500f68a32c8SEmil Constantinescu   ts->ops->setup          = TSSetUp_RK;
1501f68a32c8SEmil Constantinescu   ts->ops->interpolate    = TSInterpolate_RK;
15022c9cb062Sluzhanghpp   ts->ops->step           = TSStep_RK;
1503f68a32c8SEmil Constantinescu   ts->ops->evaluatestep   = TSEvaluateStep_RK;
1504fecfb714SLisandro Dalcin   ts->ops->rollback       = TSRollBack_RK;
1505f68a32c8SEmil Constantinescu   ts->ops->setfromoptions = TSSetFromOptions_RK;
1506ff22ae23SHong Zhang   ts->ops->getstages      = TSGetStages_RK;
1507e4dd521cSBarry Smith 
15083ca0882eSHong Zhang   ts->ops->snesfunction    = SNESTSFormFunction_RK;
15093ca0882eSHong Zhang   ts->ops->snesjacobian    = SNESTSFormJacobian_RK;
15100f7a1166SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_RK;
151113af1a74SHong Zhang   ts->ops->adjointsetup    = TSAdjointSetUp_RK;
151213af1a74SHong Zhang   ts->ops->adjointstep     = TSAdjointStep_RK;
151313af1a74SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_RK;
15140f7a1166SHong Zhang 
151513af1a74SHong Zhang   ts->ops->forwardintegral  = TSForwardCostIntegral_RK;
1516922a638cSHong Zhang   ts->ops->forwardsetup     = TSForwardSetUp_RK;
1517922a638cSHong Zhang   ts->ops->forwardreset     = TSForwardReset_RK;
1518922a638cSHong Zhang   ts->ops->forwardstep      = TSForwardStep_RK;
1519922a638cSHong Zhang   ts->ops->forwardgetstages = TSForwardGetStages_RK;
1520922a638cSHong Zhang 
15214dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&rk));
15221566a47fSLisandro Dalcin   ts->data = (void *)rk;
1523e4dd521cSBarry Smith 
15249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetOrder_C", TSRKGetOrder_RK));
15259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetType_C", TSRKGetType_RK));
15269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKSetType_C", TSRKSetType_RK));
15279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetTableau_C", TSRKGetTableau_RK));
15289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKSetMultirate_C", TSRKSetMultirate_RK));
15299566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetMultirate_C", TSRKGetMultirate_RK));
1530be5899b3SLisandro Dalcin 
15319566063dSJacob Faibussowitsch   PetscCall(TSRKSetType(ts, TSRKDefault));
153290b67129SHong Zhang   rk->dtratio = 1;
15333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1534e4dd521cSBarry Smith }
1535