xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision 1cc06b555e92f8ec64db10330b8bbd830e5bc876)
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 
32*1cc06b55SBarry 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 
44*1cc06b55SBarry 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 
56*1cc06b55SBarry 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 
68*1cc06b55SBarry 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 
83*1cc06b55SBarry 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 
95*1cc06b55SBarry 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 
107*1cc06b55SBarry 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 
122*1cc06b55SBarry 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 
137*1cc06b55SBarry 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 
152*1cc06b55SBarry 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 
167*1cc06b55SBarry 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 
182*1cc06b55SBarry 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 
192*1cc06b55SBarry 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   {
2029371c9d4SSatish Balay     const PetscReal A[1][1] = {{0}}, b[1] = {RC(1.0)};
2039566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK1FE, 1, 1, &A[0][0], b, NULL, NULL, 0, NULL));
204e4dd521cSBarry Smith   }
205db046809SBarry Smith   {
2069371c9d4SSatish Balay     const PetscReal A[2][2] =
2079371c9d4SSatish Balay       {
2089371c9d4SSatish Balay         {0,       0},
2099371c9d4SSatish Balay         {RC(1.0), 0}
2109371c9d4SSatish Balay     },
2119371c9d4SSatish Balay                     b[2] = {RC(0.5), RC(0.5)}, bembed[2] = {RC(1.0), 0};
2129566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK2A, 2, 2, &A[0][0], b, NULL, bembed, 0, NULL));
213db046809SBarry Smith   }
214f68a32c8SEmil Constantinescu   {
2159371c9d4SSatish Balay     const PetscReal A[2][2] =
2169371c9d4SSatish Balay       {
2179371c9d4SSatish Balay         {0,       0},
2189371c9d4SSatish Balay         {RC(0.5), 0}
2199371c9d4SSatish Balay     },
2209958aef7SJacob Faibussowitsch                     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   {
2249371c9d4SSatish Balay     const PetscReal A[3][3] =
2259371c9d4SSatish Balay       {
2269371c9d4SSatish Balay         {0,                  0,       0},
2274e82626cSLisandro Dalcin         {RC(2.0) / RC(3.0),  0,       0},
2289371c9d4SSatish Balay         {RC(-1.0) / RC(3.0), RC(1.0), 0}
2299371c9d4SSatish Balay     },
2304e82626cSLisandro Dalcin                     b[3] = {RC(0.25), RC(0.5), RC(0.25)};
2319566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK3, 3, 3, &A[0][0], b, NULL, NULL, 0, NULL));
232fdefd5e5SDebojyoti Ghosh   }
233fdefd5e5SDebojyoti Ghosh   {
2349371c9d4SSatish Balay     const PetscReal A[4][4] =
2359371c9d4SSatish Balay       {
2369371c9d4SSatish Balay         {0,                 0,                 0,                 0},
2374e82626cSLisandro Dalcin         {RC(1.0) / RC(2.0), 0,                 0,                 0},
2384e82626cSLisandro Dalcin         {0,                 RC(3.0) / RC(4.0), 0,                 0},
2399371c9d4SSatish Balay         {RC(2.0) / RC(9.0), RC(1.0) / RC(3.0), RC(4.0) / RC(9.0), 0}
2409371c9d4SSatish Balay     },
2419371c9d4SSatish Balay                     b[4] = {RC(2.0) / RC(9.0), RC(1.0) / RC(3.0), RC(4.0) / RC(9.0), 0}, 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)};
2429566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK3BS, 3, 4, &A[0][0], b, NULL, bembed, 0, NULL));
243db046809SBarry Smith   }
244f68a32c8SEmil Constantinescu   {
2459371c9d4SSatish Balay     const PetscReal A[4][4] =
2469371c9d4SSatish Balay       {
2479371c9d4SSatish Balay         {0,       0,       0,       0},
2484e82626cSLisandro Dalcin         {RC(0.5), 0,       0,       0},
2494e82626cSLisandro Dalcin         {0,       RC(0.5), 0,       0},
2509371c9d4SSatish Balay         {0,       0,       RC(1.0), 0}
2519371c9d4SSatish Balay     },
2524e82626cSLisandro Dalcin                     b[4] = {RC(1.0) / RC(6.0), RC(1.0) / RC(3.0), RC(1.0) / RC(3.0), RC(1.0) / RC(6.0)};
2539566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK4, 4, 4, &A[0][0], b, NULL, NULL, 0, NULL));
254db046809SBarry Smith   }
255f68a32c8SEmil Constantinescu   {
2569371c9d4SSatish Balay     const PetscReal A[6][6] =
2579371c9d4SSatish Balay       {
2589371c9d4SSatish Balay         {0,                       0,                        0,                        0,                       0,                    0},
2594e82626cSLisandro Dalcin         {RC(0.25),                0,                        0,                        0,                       0,                    0},
2604e82626cSLisandro Dalcin         {RC(3.0) / RC(32.0),      RC(9.0) / RC(32.0),       0,                        0,                       0,                    0},
2614e82626cSLisandro Dalcin         {RC(1932.0) / RC(2197.0), RC(-7200.0) / RC(2197.0), RC(7296.0) / RC(2197.0),  0,                       0,                    0},
2624e82626cSLisandro Dalcin         {RC(439.0) / RC(216.0),   RC(-8.0),                 RC(3680.0) / RC(513.0),   RC(-845.0) / RC(4104.0), 0,                    0},
2639371c9d4SSatish 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}
2649371c9d4SSatish Balay     },
2654e82626cSLisandro Dalcin                     b[6]      = {RC(16.0) / RC(135.0), 0, RC(6656.0) / RC(12825.0), RC(28561.0) / RC(56430.0), RC(-9.0) / RC(50.0), RC(2.0) / RC(55.0)},
2664e82626cSLisandro Dalcin                     bembed[6] = {RC(25.0) / RC(216.0), 0, RC(1408.0) / RC(2565.0), RC(2197.0) / RC(4104.0), RC(-1.0) / RC(5.0), 0};
2679566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK5F, 5, 6, &A[0][0], b, NULL, bembed, 0, NULL));
268fdefd5e5SDebojyoti Ghosh   }
269fdefd5e5SDebojyoti Ghosh   {
2709371c9d4SSatish Balay     const PetscReal A[7][7] =
2719371c9d4SSatish Balay       {
2729371c9d4SSatish Balay         {0,                        0,                         0,                        0,                      0,                         0,                   0},
2734e82626cSLisandro Dalcin         {RC(1.0) / RC(5.0),        0,                         0,                        0,                      0,                         0,                   0},
2744e82626cSLisandro Dalcin         {RC(3.0) / RC(40.0),       RC(9.0) / RC(40.0),        0,                        0,                      0,                         0,                   0},
2754e82626cSLisandro Dalcin         {RC(44.0) / RC(45.0),      RC(-56.0) / RC(15.0),      RC(32.0) / RC(9.0),       0,                      0,                         0,                   0},
2764e82626cSLisandro 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},
2774e82626cSLisandro 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},
2789371c9d4SSatish 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}
2799371c9d4SSatish Balay     },
2804e82626cSLisandro Dalcin                     b[7] = {RC(35.0) / RC(384.0), 0, RC(500.0) / RC(1113.0), RC(125.0) / RC(192.0), RC(-2187.0) / RC(6784.0), RC(11.0) / RC(84.0), 0},
2819371c9d4SSatish Balay                     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)}, binterp[7][5] = {{RC(1.0), RC(-4034104133.0) / RC(1410260304.0), RC(105330401.0) / RC(33982176.0), RC(-13107642775.0) / RC(11282082432.0), RC(6542295.0) / RC(470086768.0)}, {0, 0, 0, 0, 0}, {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)}, {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)}, {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)}, {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)}, {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)}};
2829566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK5DP, 5, 7, &A[0][0], b, NULL, bembed, 5, binterp[0]));
283f68a32c8SEmil Constantinescu   }
28405e23783SLisandro Dalcin   {
2859371c9d4SSatish Balay     const PetscReal A[8][8] =
2869371c9d4SSatish Balay       {
2879371c9d4SSatish Balay         {0,                           0,                          0,                              0,                            0,                          0,                           0,                        0},
2884e82626cSLisandro Dalcin         {RC(1.0) / RC(6.0),           0,                          0,                              0,                            0,                          0,                           0,                        0},
2894e82626cSLisandro Dalcin         {RC(2.0) / RC(27.0),          RC(4.0) / RC(27.0),         0,                              0,                            0,                          0,                           0,                        0},
2904e82626cSLisandro Dalcin         {RC(183.0) / RC(1372.0),      RC(-162.0) / RC(343.0),     RC(1053.0) / RC(1372.0),        0,                            0,                          0,                           0,                        0},
2914e82626cSLisandro 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},
2924e82626cSLisandro 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},
2934e82626cSLisandro 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},
2949371c9d4SSatish 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}
2959371c9d4SSatish Balay     },
2964e82626cSLisandro Dalcin                     b[8]      = {RC(587.0) / RC(8064.0), 0, RC(4440339.0) / RC(15491840.0), RC(24353.0) / RC(124800.0), RC(387.0) / RC(44800.0), RC(2152.0) / RC(5985.0), RC(7267.0) / RC(94080.0), 0},
2974e82626cSLisandro Dalcin                     bembed[8] = {RC(2479.0) / RC(34992.0), 0, RC(123.0) / RC(416.0), RC(612941.0) / RC(3411720.0), RC(43.0) / RC(1440.0), RC(2272.0) / RC(6561.0), RC(79937.0) / RC(1113912.0), RC(3293.0) / RC(556956.0)};
2989566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK5BS, 5, 8, &A[0][0], b, NULL, bembed, 0, NULL));
29905e23783SLisandro Dalcin   }
30037acaa02SHendrik Ranocha   {
3019371c9d4SSatish Balay     const PetscReal A[9][9] =
3029371c9d4SSatish Balay       {
3039371c9d4SSatish Balay         {0,                                                   0,                                                  0,                                                   0,                                                   0,                                                   0,                                                   0,                                                  0,                                                  0},
30463fe90e8SHendrik Ranocha         {RC(1.8000000000000000000000000000000000000000e-01),  0,                                                  0,                                                   0,                                                   0,                                                   0,                                                   0,                                                  0,                                                  0},
30563fe90e8SHendrik Ranocha         {RC(8.9506172839506172839506172839506172839506e-02),  RC(7.7160493827160493827160493827160493827160e-02), 0,                                                   0,                                                   0,                                                   0,                                                   0,                                                  0,                                                  0},
30663fe90e8SHendrik Ranocha         {RC(6.2500000000000000000000000000000000000000e-02),  0,                                                  RC(1.8750000000000000000000000000000000000000e-01),  0,                                                   0,                                                   0,                                                   0,                                                  0,                                                  0},
30763fe90e8SHendrik Ranocha         {RC(3.1651600000000000000000000000000000000000e-01),  0,                                                  RC(-1.0449480000000000000000000000000000000000e+00), RC(1.2584320000000000000000000000000000000000e+00),  0,                                                   0,                                                   0,                                                  0,                                                  0},
30863fe90e8SHendrik Ranocha         {RC(2.7232612736485626257225065566674305502508e-01),  0,                                                  RC(-8.2513360323886639676113360323886639676113e-01), RC(1.0480917678812415654520917678812415654521e+00),  RC(1.0471570799276856873679117969088177628396e-01),  0,                                                   0,                                                  0,                                                  0},
30963fe90e8SHendrik 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},
31063fe90e8SHendrik 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},
3119371c9d4SSatish 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}
3129371c9d4SSatish Balay     },
3139371c9d4SSatish Balay                     b[9] = {RC(7.6388888888888888888888888888888888888889e-02), 0, 0, RC(3.6940836940836940836940836940836940836941e-01), 0, RC(2.4801587301587301587301587301587301587302e-01), RC(2.3674242424242424242424242424242424242424e-01),
3149371c9d4SSatish Balay                             RC(6.9444444444444444444444444444444444444444e-02), 0},
31563fe90e8SHendrik Ranocha                     bembed[9] = {RC(5.8700209643605870020964360587002096436059e-02), 0, 0, RC(4.8072562358276643990929705215419501133787e-01), RC(-8.5341242076919085578832094861228313083563e-01), RC(1.2046485260770975056689342403628117913832e+00), 0, RC(-5.9242373072160306202859394348756050883710e-02), RC(1.6858043453788134639198468985703028256220e-01)};
3169566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK6VR, 6, 9, &A[0][0], b, NULL, bembed, 0, NULL));
31737acaa02SHendrik Ranocha   }
31837acaa02SHendrik Ranocha   {
3199371c9d4SSatish Balay     const PetscReal A[10][10] =
3209371c9d4SSatish Balay       {
3219371c9d4SSatish Balay         {0,                                                   0,                                                  0,                                                   0,                                                   0,                                                   0,                                                   0,                                                  0,                                                  0, 0},
32263fe90e8SHendrik Ranocha         {RC(5.0000000000000000000000000000000000000000e-03),  0,                                                  0,                                                   0,                                                   0,                                                   0,                                                   0,                                                  0,                                                  0, 0},
32363fe90e8SHendrik Ranocha         {RC(-1.0767901234567901234567901234567901234568e+00), RC(1.1856790123456790123456790123456790123457e+00), 0,                                                   0,                                                   0,                                                   0,                                                   0,                                                  0,                                                  0, 0},
32463fe90e8SHendrik Ranocha         {RC(4.0833333333333333333333333333333333333333e-02),  0,                                                  RC(1.2250000000000000000000000000000000000000e-01),  0,                                                   0,                                                   0,                                                   0,                                                  0,                                                  0, 0},
32563fe90e8SHendrik Ranocha         {RC(6.3607142857142857142857142857142857142857e-01),  0,                                                  RC(-2.4444642857142857142857142857142857142857e+00), RC(2.2633928571428571428571428571428571428571e+00),  0,                                                   0,                                                   0,                                                  0,                                                  0, 0},
32663fe90e8SHendrik Ranocha         {RC(-2.5351211079349245229256383554660215487207e+00), 0,                                                  RC(1.0299374654449267920438514460756024913612e+01),  RC(-7.9513032885990579949493217458266876536482e+00), RC(7.9301148923100592201226014271115261823800e-01),  0,                                                   0,                                                  0,                                                  0, 0},
32763fe90e8SHendrik 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},
32863fe90e8SHendrik 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},
32963fe90e8SHendrik 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},
3309371c9d4SSatish 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}
3319371c9d4SSatish Balay     },
3329371c9d4SSatish Balay                     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}, bembed[10] = {RC(4.7485247699299631037531273805727961552268e-02), 0, 0, RC(2.5599412588690633297154918245905393870497e-01), RC(2.7058478081067688722530891099268135732387e-01), RC(1.2505618684425992913638822323746917920448e-01),
3339371c9d4SSatish Balay                                                                                                                                                                                                                                                                                                                                                                                                                                  RC(2.5204468723743860507184043820197442562182e-01), 0, 0, RC(4.8834971521418614557381971303093137592592e-02)};
3349566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK7VR, 7, 10, &A[0][0], b, NULL, bembed, 0, NULL));
33537acaa02SHendrik Ranocha   }
33637acaa02SHendrik Ranocha   {
3379371c9d4SSatish Balay     const PetscReal A[13][13] =
3389371c9d4SSatish Balay       {
3399371c9d4SSatish Balay         {0,                                                   0,                                                  0,                                                   0,                                                   0,                                                   0,                                                   0,                                                   0,                                                  0,                                                   0,                                                  0,                                                  0, 0},
34063fe90e8SHendrik Ranocha         {RC(2.5000000000000000000000000000000000000000e-01),  0,                                                  0,                                                   0,                                                   0,                                                   0,                                                   0,                                                   0,                                                  0,                                                   0,                                                  0,                                                  0, 0},
34163fe90e8SHendrik Ranocha         {RC(8.7400846504915232052686327594877411977046e-02),  RC(2.5487604938654321753087950620345685135815e-02), 0,                                                   0,                                                   0,                                                   0,                                                   0,                                                   0,                                                  0,                                                   0,                                                  0,                                                  0, 0},
34263fe90e8SHendrik Ranocha         {RC(4.2333169291338582677165354330708661417323e-02),  0,                                                  RC(1.2699950787401574803149606299212598425197e-01),  0,                                                   0,                                                   0,                                                   0,                                                   0,                                                  0,                                                   0,                                                  0,                                                  0, 0},
34363fe90e8SHendrik Ranocha         {RC(4.2609505888742261494881445237572274090942e-01),  0,                                                  RC(-1.5987952846591523265427733230657181117089e+00), RC(1.5967002257717297115939588706899953707994e+00),  0,                                                   0,                                                   0,                                                   0,                                                  0,                                                   0,                                                  0,                                                  0, 0},
34463fe90e8SHendrik Ranocha         {RC(5.0719337296713929515090618138513639239329e-02),  0,                                                  0,                                                   RC(2.5433377264600407582754714408877778031369e-01),  RC(2.0394689005728199465736223777270858044698e-01),  0,                                                   0,                                                   0,                                                  0,                                                   0,                                                  0,                                                  0, 0},
34563fe90e8SHendrik 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},
34663fe90e8SHendrik 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},
34763fe90e8SHendrik 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},
34863fe90e8SHendrik 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},
34963fe90e8SHendrik 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},
35063fe90e8SHendrik 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},
3519371c9d4SSatish 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}
3529371c9d4SSatish Balay     },
3539371c9d4SSatish Balay                     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}, 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)};
3549566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK8VR, 8, 13, &A[0][0], b, NULL, bembed, 0, NULL));
35537acaa02SHendrik Ranocha   }
3564e82626cSLisandro Dalcin #undef RC
3573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
358db046809SBarry Smith }
359db046809SBarry Smith 
360f68a32c8SEmil Constantinescu /*@C
361bcf0153eSBarry Smith    TSRKRegisterDestroy - Frees the list of schemes that were registered by `TSRKRegister()`.
362f68a32c8SEmil Constantinescu 
363f68a32c8SEmil Constantinescu    Not Collective
364f68a32c8SEmil Constantinescu 
365f68a32c8SEmil Constantinescu    Level: advanced
366f68a32c8SEmil Constantinescu 
367*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKRegister()`, `TSRKRegisterAll()`
368f68a32c8SEmil Constantinescu @*/
369d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKRegisterDestroy(void)
370d71ae5a4SJacob Faibussowitsch {
371f68a32c8SEmil Constantinescu   RKTableauLink link;
372f68a32c8SEmil Constantinescu 
373f68a32c8SEmil Constantinescu   PetscFunctionBegin;
374f68a32c8SEmil Constantinescu   while ((link = RKTableauList)) {
375f68a32c8SEmil Constantinescu     RKTableau t   = &link->tab;
376f68a32c8SEmil Constantinescu     RKTableauList = link->next;
3779566063dSJacob Faibussowitsch     PetscCall(PetscFree3(t->A, t->b, t->c));
3789566063dSJacob Faibussowitsch     PetscCall(PetscFree(t->bembed));
3799566063dSJacob Faibussowitsch     PetscCall(PetscFree(t->binterp));
3809566063dSJacob Faibussowitsch     PetscCall(PetscFree(t->name));
3819566063dSJacob Faibussowitsch     PetscCall(PetscFree(link));
382f68a32c8SEmil Constantinescu   }
383f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_FALSE;
3843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
385f68a32c8SEmil Constantinescu }
386f68a32c8SEmil Constantinescu 
387f68a32c8SEmil Constantinescu /*@C
388bcf0153eSBarry Smith   TSRKInitializePackage - This function initializes everything in the `TSRK` package. It is called
389bcf0153eSBarry Smith   from `TSInitializePackage()`.
390f68a32c8SEmil Constantinescu 
391f68a32c8SEmil Constantinescu   Level: developer
392f68a32c8SEmil Constantinescu 
393*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSInitializePackage()`, `PetscInitialize()`, `TSRKFinalizePackage()`
394f68a32c8SEmil Constantinescu @*/
395d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKInitializePackage(void)
396d71ae5a4SJacob Faibussowitsch {
397e4dd521cSBarry Smith   PetscFunctionBegin;
3983ba16761SJacob Faibussowitsch   if (TSRKPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
399f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_TRUE;
4009566063dSJacob Faibussowitsch   PetscCall(TSRKRegisterAll());
4019566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(TSRKFinalizePackage));
4023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
403f68a32c8SEmil Constantinescu }
404186e87acSLisandro Dalcin 
405f68a32c8SEmil Constantinescu /*@C
406bcf0153eSBarry Smith   TSRKFinalizePackage - This function destroys everything in the `TSRK` package. It is
407bcf0153eSBarry Smith   called from `PetscFinalize()`.
408186e87acSLisandro Dalcin 
409f68a32c8SEmil Constantinescu   Level: developer
410f68a32c8SEmil Constantinescu 
411*1cc06b55SBarry Smith .seealso: [](ch_ts), `PetscFinalize()`, `TSRKInitializePackage()`
412f68a32c8SEmil Constantinescu @*/
413d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKFinalizePackage(void)
414d71ae5a4SJacob Faibussowitsch {
415f68a32c8SEmil Constantinescu   PetscFunctionBegin;
416f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_FALSE;
4179566063dSJacob Faibussowitsch   PetscCall(TSRKRegisterDestroy());
4183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
419f68a32c8SEmil Constantinescu }
420f68a32c8SEmil Constantinescu 
421f68a32c8SEmil Constantinescu /*@C
422bcf0153eSBarry Smith    TSRKRegister - register an `TSRK` scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
423f68a32c8SEmil Constantinescu 
424f68a32c8SEmil Constantinescu    Not Collective, but the same schemes should be registered on all processes on which they will be used
425f68a32c8SEmil Constantinescu 
426f68a32c8SEmil Constantinescu    Input Parameters:
427f68a32c8SEmil Constantinescu +  name - identifier for method
428f68a32c8SEmil Constantinescu .  order - approximation order of method
429f68a32c8SEmil Constantinescu .  s - number of stages, this is the dimension of the matrices below
430f68a32c8SEmil Constantinescu .  A - stage coefficients (dimension s*s, row-major)
431f68a32c8SEmil Constantinescu .  b - step completion table (dimension s; NULL to use last row of A)
432f68a32c8SEmil Constantinescu .  c - abscissa (dimension s; NULL to use row sums of A)
433f68a32c8SEmil Constantinescu .  bembed - completion table for embedded method (dimension s; NULL if not available)
4343a8a9803SLisandro Dalcin .  p - Order of the interpolation scheme, equal to the number of columns of binterp
4353a8a9803SLisandro Dalcin -  binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1)
436f68a32c8SEmil Constantinescu 
437f68a32c8SEmil Constantinescu    Level: advanced
438f68a32c8SEmil Constantinescu 
439bcf0153eSBarry Smith    Note:
440bcf0153eSBarry Smith    Several `TSRK` methods are provided, this function is only needed to create new methods.
441bcf0153eSBarry Smith 
442*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`
443f68a32c8SEmil Constantinescu @*/
444d71ae5a4SJacob 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[])
445d71ae5a4SJacob Faibussowitsch {
446f68a32c8SEmil Constantinescu   RKTableauLink link;
447f68a32c8SEmil Constantinescu   RKTableau     t;
448f68a32c8SEmil Constantinescu   PetscInt      i, j;
449f68a32c8SEmil Constantinescu 
450f68a32c8SEmil Constantinescu   PetscFunctionBegin;
4513a8a9803SLisandro Dalcin   PetscValidCharPointer(name, 1);
4523a8a9803SLisandro Dalcin   PetscValidRealPointer(A, 4);
4533a8a9803SLisandro Dalcin   if (b) PetscValidRealPointer(b, 5);
4543a8a9803SLisandro Dalcin   if (c) PetscValidRealPointer(c, 6);
4553a8a9803SLisandro Dalcin   if (bembed) PetscValidRealPointer(bembed, 7);
4563a8a9803SLisandro Dalcin   if (binterp || p > 1) PetscValidRealPointer(binterp, 9);
4573a8a9803SLisandro Dalcin 
4589566063dSJacob Faibussowitsch   PetscCall(TSRKInitializePackage());
4599566063dSJacob Faibussowitsch   PetscCall(PetscNew(&link));
460f68a32c8SEmil Constantinescu   t = &link->tab;
4613a8a9803SLisandro Dalcin 
4629566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &t->name));
463f68a32c8SEmil Constantinescu   t->order = order;
464f68a32c8SEmil Constantinescu   t->s     = s;
4659566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(s * s, &t->A, s, &t->b, s, &t->c));
4669566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->A, A, s * s));
4679566063dSJacob Faibussowitsch   if (b) PetscCall(PetscArraycpy(t->b, b, s));
4689371c9d4SSatish Balay   else
4699371c9d4SSatish Balay     for (i = 0; i < s; i++) t->b[i] = A[(s - 1) * s + i];
4709566063dSJacob Faibussowitsch   if (c) PetscCall(PetscArraycpy(t->c, c, s));
4719371c9d4SSatish Balay   else
4729371c9d4SSatish Balay     for (i = 0; i < s; i++)
4739371c9d4SSatish Balay       for (j = 0, t->c[i] = 0; j < s; j++) t->c[i] += A[i * s + j];
474d760c35bSDebojyoti Ghosh   t->FSAL = PETSC_TRUE;
4759371c9d4SSatish Balay   for (i = 0; i < s; i++)
4769371c9d4SSatish Balay     if (t->A[(s - 1) * s + i] != t->b[i]) t->FSAL = PETSC_FALSE;
4773a8a9803SLisandro Dalcin 
478f68a32c8SEmil Constantinescu   if (bembed) {
4799566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(s, &t->bembed));
4809566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(t->bembed, bembed, s));
481f68a32c8SEmil Constantinescu   }
482f68a32c8SEmil Constantinescu 
4839371c9d4SSatish Balay   if (!binterp) {
4849371c9d4SSatish Balay     p       = 1;
4859371c9d4SSatish Balay     binterp = t->b;
4869371c9d4SSatish Balay   }
4873a8a9803SLisandro Dalcin   t->p = p;
4889566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s * p, &t->binterp));
4899566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->binterp, binterp, s * p));
4903a8a9803SLisandro Dalcin 
491f68a32c8SEmil Constantinescu   link->next    = RKTableauList;
492f68a32c8SEmil Constantinescu   RKTableauList = link;
4933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
494f68a32c8SEmil Constantinescu }
495f68a32c8SEmil Constantinescu 
496d71ae5a4SJacob 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)
497d71ae5a4SJacob Faibussowitsch {
4980f7a28cdSStefano Zampini   TS_RK    *rk  = (TS_RK *)ts->data;
4990f7a28cdSStefano Zampini   RKTableau tab = rk->tableau;
5000f7a28cdSStefano Zampini 
5010f7a28cdSStefano Zampini   PetscFunctionBegin;
5020f7a28cdSStefano Zampini   if (s) *s = tab->s;
5030f7a28cdSStefano Zampini   if (A) *A = tab->A;
5040f7a28cdSStefano Zampini   if (b) *b = tab->b;
5050f7a28cdSStefano Zampini   if (c) *c = tab->c;
5060f7a28cdSStefano Zampini   if (bembed) *bembed = tab->bembed;
5070f7a28cdSStefano Zampini   if (p) *p = tab->p;
5080f7a28cdSStefano Zampini   if (binterp) *binterp = tab->binterp;
5090f7a28cdSStefano Zampini   if (FSAL) *FSAL = tab->FSAL;
5103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5110f7a28cdSStefano Zampini }
5120f7a28cdSStefano Zampini 
5130f7a28cdSStefano Zampini /*@C
514bcf0153eSBarry Smith    TSRKGetTableau - Get info on the `TSRK` tableau
5150f7a28cdSStefano Zampini 
5160f7a28cdSStefano Zampini    Not Collective
5170f7a28cdSStefano Zampini 
518f899ff85SJose E. Roman    Input Parameter:
5190f7a28cdSStefano Zampini .  ts - timestepping context
5200f7a28cdSStefano Zampini 
5210f7a28cdSStefano Zampini    Output Parameters:
5220f7a28cdSStefano Zampini +  s - number of stages, this is the dimension of the matrices below
5230f7a28cdSStefano Zampini .  A - stage coefficients (dimension s*s, row-major)
5240f7a28cdSStefano Zampini .  b - step completion table (dimension s)
5250f7a28cdSStefano Zampini .  c - abscissa (dimension s)
5260f7a28cdSStefano Zampini .  bembed - completion table for embedded method (dimension s; NULL if not available)
5270f7a28cdSStefano Zampini .  p - Order of the interpolation scheme, equal to the number of columns of binterp
5280f7a28cdSStefano Zampini .  binterp - Coefficients of the interpolation formula (dimension s*p)
529aaa8cc7dSPierre Jolivet -  FSAL - whether or not the scheme has the First Same As Last property
5300f7a28cdSStefano Zampini 
5310f7a28cdSStefano Zampini    Level: developer
5320f7a28cdSStefano Zampini 
533*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKRegister()`, `TSRKSetType()`
5340f7a28cdSStefano Zampini @*/
535d71ae5a4SJacob 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)
536d71ae5a4SJacob Faibussowitsch {
5370f7a28cdSStefano Zampini   PetscFunctionBegin;
5380f7a28cdSStefano Zampini   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5399371c9d4SSatish 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));
5403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5410f7a28cdSStefano Zampini }
5420f7a28cdSStefano Zampini 
543e4dd521cSBarry Smith /*
5442c9cb062Sluzhanghpp  This is for single-step RK method
545f68a32c8SEmil Constantinescu  The step completion formula is
546e4dd521cSBarry Smith 
547f68a32c8SEmil Constantinescu  x1 = x0 + h b^T YdotRHS
548f68a32c8SEmil Constantinescu 
549f68a32c8SEmil Constantinescu  This function can be called before or after ts->vec_sol has been updated.
550f68a32c8SEmil Constantinescu  Suppose we have a completion formula (b) and an embedded formula (be) of different order.
551f68a32c8SEmil Constantinescu  We can write
552f68a32c8SEmil Constantinescu 
553f68a32c8SEmil Constantinescu  x1e = x0 + h be^T YdotRHS
554f68a32c8SEmil Constantinescu      = x1 - h b^T YdotRHS + h be^T YdotRHS
555f68a32c8SEmil Constantinescu      = x1 + h (be - b)^T YdotRHS
556f68a32c8SEmil Constantinescu 
557f68a32c8SEmil Constantinescu  so we can evaluate the method with different order even after the step has been optimistically completed.
558e4dd521cSBarry Smith */
559d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateStep_RK(TS ts, PetscInt order, Vec X, PetscBool *done)
560d71ae5a4SJacob Faibussowitsch {
561f68a32c8SEmil Constantinescu   TS_RK       *rk  = (TS_RK *)ts->data;
562f68a32c8SEmil Constantinescu   RKTableau    tab = rk->tableau;
563f68a32c8SEmil Constantinescu   PetscScalar *w   = rk->work;
564f68a32c8SEmil Constantinescu   PetscReal    h;
565f68a32c8SEmil Constantinescu   PetscInt     s = tab->s, j;
566f68a32c8SEmil Constantinescu 
567f68a32c8SEmil Constantinescu   PetscFunctionBegin;
568f68a32c8SEmil Constantinescu   switch (rk->status) {
569f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
570d71ae5a4SJacob Faibussowitsch   case TS_STEP_PENDING:
571d71ae5a4SJacob Faibussowitsch     h = ts->time_step;
572d71ae5a4SJacob Faibussowitsch     break;
573d71ae5a4SJacob Faibussowitsch   case TS_STEP_COMPLETE:
574d71ae5a4SJacob Faibussowitsch     h = ts->ptime - ts->ptime_prev;
575d71ae5a4SJacob Faibussowitsch     break;
576d71ae5a4SJacob Faibussowitsch   default:
577d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus");
578f68a32c8SEmil Constantinescu   }
579f68a32c8SEmil Constantinescu   if (order == tab->order) {
580f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) {
5819566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol, X));
58290b67129SHong Zhang       for (j = 0; j < s; j++) w[j] = h * tab->b[j] / rk->dtratio;
5839566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(X, s, w, rk->YdotRHS));
5849566063dSJacob Faibussowitsch     } else PetscCall(VecCopy(ts->vec_sol, X));
5853ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
586f68a32c8SEmil Constantinescu   } else if (order == tab->order - 1) {
587f68a32c8SEmil Constantinescu     if (!tab->bembed) goto unavailable;
588f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) { /*Complete with the embedded method (be)*/
5899566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol, X));
590f68a32c8SEmil Constantinescu       for (j = 0; j < s; j++) w[j] = h * tab->bembed[j];
5919566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(X, s, w, rk->YdotRHS));
592f68a32c8SEmil Constantinescu     } else { /*Rollback and re-complete using (be-b) */
5939566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol, X));
594f68a32c8SEmil Constantinescu       for (j = 0; j < s; j++) w[j] = h * (tab->bembed[j] - tab->b[j]);
5959566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(X, s, w, rk->YdotRHS));
596f68a32c8SEmil Constantinescu     }
597f68a32c8SEmil Constantinescu     if (done) *done = PETSC_TRUE;
5983ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
599f68a32c8SEmil Constantinescu   }
600f68a32c8SEmil Constantinescu unavailable:
601f68a32c8SEmil Constantinescu   if (done) *done = PETSC_FALSE;
6029371c9d4SSatish Balay   else
6039371c9d4SSatish 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);
6043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
605f68a32c8SEmil Constantinescu }
606f68a32c8SEmil Constantinescu 
607d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardCostIntegral_RK(TS ts)
608d71ae5a4SJacob Faibussowitsch {
6090f7a1166SHong Zhang   TS_RK           *rk     = (TS_RK *)ts->data;
610cd4cee2dSHong Zhang   TS               quadts = ts->quadraturets;
6110f7a1166SHong Zhang   RKTableau        tab    = rk->tableau;
6120f7a1166SHong Zhang   const PetscInt   s      = tab->s;
6130f7a1166SHong Zhang   const PetscReal *b = tab->b, *c = tab->c;
6140f7a1166SHong Zhang   Vec             *Y = rk->Y;
6150f7a1166SHong Zhang   PetscInt         i;
6160f7a1166SHong Zhang 
6170f7a1166SHong Zhang   PetscFunctionBegin;
618cd4cee2dSHong Zhang   /* No need to backup quadts->vec_sol since it can be reverted in TSRollBack_RK */
6190f7a1166SHong Zhang   for (i = s - 1; i >= 0; i--) {
620cd4cee2dSHong Zhang     /* Evolve quadrature TS solution to compute integrals */
6219566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(quadts, rk->ptime + rk->time_step * c[i], Y[i], ts->vec_costintegrand));
6229566063dSJacob Faibussowitsch     PetscCall(VecAXPY(quadts->vec_sol, rk->time_step * b[i], ts->vec_costintegrand));
6230f7a1166SHong Zhang   }
6243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6250f7a1166SHong Zhang }
6260f7a1166SHong Zhang 
627d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointCostIntegral_RK(TS ts)
628d71ae5a4SJacob Faibussowitsch {
6290f7a1166SHong Zhang   TS_RK           *rk     = (TS_RK *)ts->data;
6300f7a1166SHong Zhang   RKTableau        tab    = rk->tableau;
631cd4cee2dSHong Zhang   TS               quadts = ts->quadraturets;
6320f7a1166SHong Zhang   const PetscInt   s      = tab->s;
6330f7a1166SHong Zhang   const PetscReal *b = tab->b, *c = tab->c;
6340f7a1166SHong Zhang   Vec             *Y = rk->Y;
6350f7a1166SHong Zhang   PetscInt         i;
6360f7a1166SHong Zhang 
6370f7a1166SHong Zhang   PetscFunctionBegin;
6380f7a1166SHong Zhang   for (i = s - 1; i >= 0; i--) {
639cd4cee2dSHong Zhang     /* Evolve quadrature TS solution to compute integrals */
6409566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(quadts, ts->ptime + ts->time_step * (1.0 - c[i]), Y[i], ts->vec_costintegrand));
6419566063dSJacob Faibussowitsch     PetscCall(VecAXPY(quadts->vec_sol, -ts->time_step * b[i], ts->vec_costintegrand));
6420f7a1166SHong Zhang   }
6433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6440f7a1166SHong Zhang }
6450f7a1166SHong Zhang 
646d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRollBack_RK(TS ts)
647d71ae5a4SJacob Faibussowitsch {
648fecfb714SLisandro Dalcin   TS_RK           *rk     = (TS_RK *)ts->data;
649cd4cee2dSHong Zhang   TS               quadts = ts->quadraturets;
650fecfb714SLisandro Dalcin   RKTableau        tab    = rk->tableau;
651fecfb714SLisandro Dalcin   const PetscInt   s      = tab->s;
652cd4cee2dSHong Zhang   const PetscReal *b = tab->b, *c = tab->c;
653fecfb714SLisandro Dalcin   PetscScalar     *w = rk->work;
654cd4cee2dSHong Zhang   Vec             *Y = rk->Y, *YdotRHS = rk->YdotRHS;
655fecfb714SLisandro Dalcin   PetscInt         j;
656fecfb714SLisandro Dalcin   PetscReal        h;
657fecfb714SLisandro Dalcin 
658fecfb714SLisandro Dalcin   PetscFunctionBegin;
659fecfb714SLisandro Dalcin   switch (rk->status) {
660fecfb714SLisandro Dalcin   case TS_STEP_INCOMPLETE:
661d71ae5a4SJacob Faibussowitsch   case TS_STEP_PENDING:
662d71ae5a4SJacob Faibussowitsch     h = ts->time_step;
663d71ae5a4SJacob Faibussowitsch     break;
664d71ae5a4SJacob Faibussowitsch   case TS_STEP_COMPLETE:
665d71ae5a4SJacob Faibussowitsch     h = ts->ptime - ts->ptime_prev;
666d71ae5a4SJacob Faibussowitsch     break;
667d71ae5a4SJacob Faibussowitsch   default:
668d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus");
669fecfb714SLisandro Dalcin   }
670fecfb714SLisandro Dalcin   for (j = 0; j < s; j++) w[j] = -h * b[j];
6719566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(ts->vec_sol, s, w, YdotRHS));
672cd4cee2dSHong Zhang   if (quadts && ts->costintegralfwd) {
673cd4cee2dSHong Zhang     for (j = 0; j < s; j++) {
674cd4cee2dSHong Zhang       /* Revert the quadrature TS solution */
6759566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSFunction(quadts, rk->ptime + h * c[j], Y[j], ts->vec_costintegrand));
6769566063dSJacob Faibussowitsch       PetscCall(VecAXPY(quadts->vec_sol, -h * b[j], ts->vec_costintegrand));
677cd4cee2dSHong Zhang     }
678cd4cee2dSHong Zhang   }
6793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
680fecfb714SLisandro Dalcin }
681fecfb714SLisandro Dalcin 
682d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardStep_RK(TS ts)
683d71ae5a4SJacob Faibussowitsch {
684922a638cSHong Zhang   TS_RK           *rk  = (TS_RK *)ts->data;
685922a638cSHong Zhang   RKTableau        tab = rk->tableau;
686922a638cSHong Zhang   Mat              J, *MatsFwdSensipTemp = rk->MatsFwdSensipTemp;
687922a638cSHong Zhang   const PetscInt   s = tab->s;
688922a638cSHong Zhang   const PetscReal *A = tab->A, *c = tab->c, *b = tab->b;
689922a638cSHong Zhang   Vec             *Y = rk->Y;
690922a638cSHong Zhang   PetscInt         i, j;
691922a638cSHong Zhang   PetscReal        stage_time, h = ts->time_step;
692922a638cSHong Zhang   PetscBool        zero;
693922a638cSHong Zhang 
694922a638cSHong Zhang   PetscFunctionBegin;
6959566063dSJacob Faibussowitsch   PetscCall(MatCopy(ts->mat_sensip, rk->MatFwdSensip0, SAME_NONZERO_PATTERN));
6969566063dSJacob Faibussowitsch   PetscCall(TSGetRHSJacobian(ts, &J, NULL, NULL, NULL));
697922a638cSHong Zhang 
698922a638cSHong Zhang   for (i = 0; i < s; i++) {
699922a638cSHong Zhang     stage_time = ts->ptime + h * c[i];
700922a638cSHong Zhang     zero       = PETSC_FALSE;
701922a638cSHong Zhang     if (b[i] == 0 && i == s - 1) zero = PETSC_TRUE;
702922a638cSHong Zhang     /* TLM Stage values */
703922a638cSHong Zhang     if (!i) {
7049566063dSJacob Faibussowitsch       PetscCall(MatCopy(ts->mat_sensip, rk->MatsFwdStageSensip[i], SAME_NONZERO_PATTERN));
705922a638cSHong Zhang     } else if (!zero) {
7069566063dSJacob Faibussowitsch       PetscCall(MatZeroEntries(rk->MatsFwdStageSensip[i]));
70748a46eb9SPierre Jolivet       for (j = 0; j < i; j++) PetscCall(MatAXPY(rk->MatsFwdStageSensip[i], h * A[i * s + j], MatsFwdSensipTemp[j], SAME_NONZERO_PATTERN));
7089566063dSJacob Faibussowitsch       PetscCall(MatAXPY(rk->MatsFwdStageSensip[i], 1., ts->mat_sensip, SAME_NONZERO_PATTERN));
709922a638cSHong Zhang     } else {
7109566063dSJacob Faibussowitsch       PetscCall(MatZeroEntries(rk->MatsFwdStageSensip[i]));
711922a638cSHong Zhang     }
712922a638cSHong Zhang 
7139566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSJacobian(ts, stage_time, Y[i], J, J));
7149566063dSJacob Faibussowitsch     PetscCall(MatMatMult(J, rk->MatsFwdStageSensip[i], MAT_REUSE_MATRIX, PETSC_DEFAULT, &MatsFwdSensipTemp[i]));
71513af1a74SHong Zhang     if (ts->Jacprhs) {
7169566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(ts, stage_time, Y[i], ts->Jacprhs)); /* get f_p */
71713af1a74SHong Zhang       if (ts->vecs_sensi2p) {                                              /* TLM used for 2nd-order adjoint */
71813af1a74SHong Zhang         PetscScalar *xarr;
7199566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(MatsFwdSensipTemp[i], 0, &xarr));
7209566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(rk->VecDeltaFwdSensipCol, xarr));
7219566063dSJacob Faibussowitsch         PetscCall(MatMultAdd(ts->Jacprhs, ts->vec_dir, rk->VecDeltaFwdSensipCol, rk->VecDeltaFwdSensipCol));
7229566063dSJacob Faibussowitsch         PetscCall(VecResetArray(rk->VecDeltaFwdSensipCol));
7239566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(MatsFwdSensipTemp[i], &xarr));
72413af1a74SHong Zhang       } else {
7259566063dSJacob Faibussowitsch         PetscCall(MatAXPY(MatsFwdSensipTemp[i], 1., ts->Jacprhs, SUBSET_NONZERO_PATTERN));
72613af1a74SHong Zhang       }
727922a638cSHong Zhang     }
728922a638cSHong Zhang   }
729922a638cSHong Zhang 
73048a46eb9SPierre Jolivet   for (i = 0; i < s; i++) PetscCall(MatAXPY(ts->mat_sensip, h * b[i], rk->MatsFwdSensipTemp[i], SAME_NONZERO_PATTERN));
731922a638cSHong Zhang   rk->status = TS_STEP_COMPLETE;
7323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
733922a638cSHong Zhang }
734922a638cSHong Zhang 
735d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardGetStages_RK(TS ts, PetscInt *ns, Mat **stagesensip)
736d71ae5a4SJacob Faibussowitsch {
737922a638cSHong Zhang   TS_RK    *rk  = (TS_RK *)ts->data;
738922a638cSHong Zhang   RKTableau tab = rk->tableau;
739922a638cSHong Zhang 
740922a638cSHong Zhang   PetscFunctionBegin;
741922a638cSHong Zhang   if (ns) *ns = tab->s;
742922a638cSHong Zhang   if (stagesensip) *stagesensip = rk->MatsFwdStageSensip;
7433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
744922a638cSHong Zhang }
745922a638cSHong Zhang 
746d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardSetUp_RK(TS ts)
747d71ae5a4SJacob Faibussowitsch {
748922a638cSHong Zhang   TS_RK    *rk  = (TS_RK *)ts->data;
749922a638cSHong Zhang   RKTableau tab = rk->tableau;
750922a638cSHong Zhang   PetscInt  i;
751922a638cSHong Zhang 
752922a638cSHong Zhang   PetscFunctionBegin;
753922a638cSHong Zhang   /* backup sensitivity results for roll-backs */
7549566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &rk->MatFwdSensip0));
755922a638cSHong Zhang 
7569566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tab->s, &rk->MatsFwdStageSensip));
7579566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tab->s, &rk->MatsFwdSensipTemp));
758922a638cSHong Zhang   for (i = 0; i < tab->s; i++) {
7599566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &rk->MatsFwdStageSensip[i]));
7609566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &rk->MatsFwdSensipTemp[i]));
761922a638cSHong Zhang   }
7629566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &rk->VecDeltaFwdSensipCol));
7633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
764922a638cSHong Zhang }
765922a638cSHong Zhang 
766d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardReset_RK(TS ts)
767d71ae5a4SJacob Faibussowitsch {
768922a638cSHong Zhang   TS_RK    *rk  = (TS_RK *)ts->data;
769922a638cSHong Zhang   RKTableau tab = rk->tableau;
770922a638cSHong Zhang   PetscInt  i;
771922a638cSHong Zhang 
772922a638cSHong Zhang   PetscFunctionBegin;
7739566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&rk->MatFwdSensip0));
774922a638cSHong Zhang   if (rk->MatsFwdStageSensip) {
77548a46eb9SPierre Jolivet     for (i = 0; i < tab->s; i++) PetscCall(MatDestroy(&rk->MatsFwdStageSensip[i]));
7769566063dSJacob Faibussowitsch     PetscCall(PetscFree(rk->MatsFwdStageSensip));
777922a638cSHong Zhang   }
778922a638cSHong Zhang   if (rk->MatsFwdSensipTemp) {
77948a46eb9SPierre Jolivet     for (i = 0; i < tab->s; i++) PetscCall(MatDestroy(&rk->MatsFwdSensipTemp[i]));
7809566063dSJacob Faibussowitsch     PetscCall(PetscFree(rk->MatsFwdSensipTemp));
781922a638cSHong Zhang   }
7829566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&rk->VecDeltaFwdSensipCol));
7833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
784922a638cSHong Zhang }
785922a638cSHong Zhang 
786d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_RK(TS ts)
787d71ae5a4SJacob Faibussowitsch {
788f68a32c8SEmil Constantinescu   TS_RK           *rk  = (TS_RK *)ts->data;
789f68a32c8SEmil Constantinescu   RKTableau        tab = rk->tableau;
790f68a32c8SEmil Constantinescu   const PetscInt   s   = tab->s;
791fecfb714SLisandro Dalcin   const PetscReal *A = tab->A, *c = tab->c;
792f68a32c8SEmil Constantinescu   PetscScalar     *w = rk->work;
793f68a32c8SEmil Constantinescu   Vec             *Y = rk->Y, *YdotRHS = rk->YdotRHS;
794d1905564SLisandro Dalcin   PetscBool        FSAL = tab->FSAL;
795f68a32c8SEmil Constantinescu   TSAdapt          adapt;
796fecfb714SLisandro Dalcin   PetscInt         i, j;
797be5899b3SLisandro Dalcin   PetscInt         rejections = 0;
798be5899b3SLisandro Dalcin   PetscBool        stageok, accept = PETSC_TRUE;
799be5899b3SLisandro Dalcin   PetscReal        next_time_step = ts->time_step;
800f68a32c8SEmil Constantinescu 
801f68a32c8SEmil Constantinescu   PetscFunctionBegin;
802d1905564SLisandro Dalcin   if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE;
8039566063dSJacob Faibussowitsch   if (FSAL) PetscCall(VecCopy(YdotRHS[s - 1], YdotRHS[0]));
804d1905564SLisandro Dalcin 
805f68a32c8SEmil Constantinescu   rk->status = TS_STEP_INCOMPLETE;
806be5899b3SLisandro Dalcin   while (!ts->reason && rk->status != TS_STEP_COMPLETE) {
807be5899b3SLisandro Dalcin     PetscReal t = ts->ptime;
808f68a32c8SEmil Constantinescu     PetscReal h = ts->time_step;
809f68a32c8SEmil Constantinescu     for (i = 0; i < s; i++) {
8109be3e283SDebojyoti Ghosh       rk->stage_time = t + h * c[i];
8119566063dSJacob Faibussowitsch       PetscCall(TSPreStage(ts, rk->stage_time));
8129566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol, Y[i]));
813f68a32c8SEmil Constantinescu       for (j = 0; j < i; j++) w[j] = h * A[i * s + j];
8149566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(Y[i], i, w, YdotRHS));
8159566063dSJacob Faibussowitsch       PetscCall(TSPostStage(ts, rk->stage_time, i, Y));
8169566063dSJacob Faibussowitsch       PetscCall(TSGetAdapt(ts, &adapt));
8179566063dSJacob Faibussowitsch       PetscCall(TSAdaptCheckStage(adapt, ts, rk->stage_time, Y[i], &stageok));
818fecfb714SLisandro Dalcin       if (!stageok) goto reject_step;
8198f738a7cSLisandro Dalcin       if (FSAL && !i) continue;
8209566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSFunction(ts, t + h * c[i], Y[i], YdotRHS[i]));
821f68a32c8SEmil Constantinescu     }
822be5899b3SLisandro Dalcin 
823be5899b3SLisandro Dalcin     rk->status = TS_STEP_INCOMPLETE;
8249566063dSJacob Faibussowitsch     PetscCall(TSEvaluateStep(ts, tab->order, ts->vec_sol, NULL));
825f68a32c8SEmil Constantinescu     rk->status = TS_STEP_PENDING;
8269566063dSJacob Faibussowitsch     PetscCall(TSGetAdapt(ts, &adapt));
8279566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidatesClear(adapt));
8289566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidateAdd(adapt, tab->name, tab->order, 1, tab->ccfl, (PetscReal)tab->s, PETSC_TRUE));
8299566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
830be5899b3SLisandro Dalcin     rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
831be5899b3SLisandro Dalcin     if (!accept) { /* Roll back the current step */
8329566063dSJacob Faibussowitsch       PetscCall(TSRollBack_RK(ts));
833be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
834be5899b3SLisandro Dalcin       goto reject_step;
835be5899b3SLisandro Dalcin     }
836be5899b3SLisandro Dalcin 
8370f7a1166SHong Zhang     if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
8380f7a1166SHong Zhang       rk->ptime     = ts->ptime;
8390f7a1166SHong Zhang       rk->time_step = ts->time_step;
840493ed95fSHong Zhang     }
841be5899b3SLisandro Dalcin 
842f68a32c8SEmil Constantinescu     ts->ptime += ts->time_step;
843f68a32c8SEmil Constantinescu     ts->time_step = next_time_step;
844f68a32c8SEmil Constantinescu     break;
845be5899b3SLisandro Dalcin 
846be5899b3SLisandro Dalcin   reject_step:
8479371c9d4SSatish Balay     ts->reject++;
8489371c9d4SSatish Balay     accept = PETSC_FALSE;
849be5899b3SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
850be5899b3SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
85163a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
852f68a32c8SEmil Constantinescu     }
853f68a32c8SEmil Constantinescu   }
8543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
855e4dd521cSBarry Smith }
856e4dd521cSBarry Smith 
857d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointSetUp_RK(TS ts)
858d71ae5a4SJacob Faibussowitsch {
859f6a906c0SBarry Smith   TS_RK    *rk  = (TS_RK *)ts->data;
860be5899b3SLisandro Dalcin   RKTableau tab = rk->tableau;
861be5899b3SLisandro Dalcin   PetscInt  s   = tab->s;
862f6a906c0SBarry Smith 
863f6a906c0SBarry Smith   PetscFunctionBegin;
8643ba16761SJacob Faibussowitsch   if (ts->adjointsetupcalled++) PetscFunctionReturn(PETSC_SUCCESS);
8659566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], s * ts->numcost, &rk->VecsDeltaLam));
8669566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &rk->VecsSensiTemp));
86748a46eb9SPierre Jolivet   if (ts->vecs_sensip) PetscCall(VecDuplicate(ts->vecs_sensip[0], &rk->VecDeltaMu));
86813af1a74SHong Zhang   if (ts->vecs_sensi2) {
8699566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], s * ts->numcost, &rk->VecsDeltaLam2));
8709566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi2[0], ts->numcost, &rk->VecsSensi2Temp));
87113af1a74SHong Zhang   }
87248a46eb9SPierre Jolivet   if (ts->vecs_sensi2p) PetscCall(VecDuplicate(ts->vecs_sensi2p[0], &rk->VecDeltaMu2));
8733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
874f6a906c0SBarry Smith }
875f6a906c0SBarry Smith 
87635f5172aSHong Zhang /*
87735f5172aSHong Zhang   Assumptions:
87835f5172aSHong Zhang     - TSStep_RK() always evaluates the step with b, not bembed.
87935f5172aSHong Zhang */
880d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointStep_RK(TS ts)
881d71ae5a4SJacob Faibussowitsch {
882c235aa8dSHong Zhang   TS_RK           *rk     = (TS_RK *)ts->data;
883cd4cee2dSHong Zhang   TS               quadts = ts->quadraturets;
884c235aa8dSHong Zhang   RKTableau        tab    = rk->tableau;
8853ca0882eSHong Zhang   Mat              J, Jpre, Jquad;
886c235aa8dSHong Zhang   const PetscInt   s = tab->s;
887c235aa8dSHong Zhang   const PetscReal *A = tab->A, *b = tab->b, *c = tab->c;
88813af1a74SHong Zhang   PetscScalar     *w = rk->work, *xarr;
8892e7b7f96SHong Zhang   Vec             *Y = rk->Y, *VecsDeltaLam = rk->VecsDeltaLam, VecDeltaMu = rk->VecDeltaMu, *VecsSensiTemp = rk->VecsSensiTemp;
89013af1a74SHong Zhang   Vec             *VecsDeltaLam2 = rk->VecsDeltaLam2, VecDeltaMu2 = rk->VecDeltaMu2, *VecsSensi2Temp = rk->VecsSensi2Temp;
891cd4cee2dSHong Zhang   Vec              VecDRDUTransCol = ts->vec_drdu_col, VecDRDPTransCol = ts->vec_drdp_col;
892b18ea86cSHong Zhang   PetscInt         i, j, nadj;
8933d3eaba7SBarry Smith   PetscReal        t = ts->ptime;
8943ca0882eSHong Zhang   PetscReal        h = ts->time_step;
895c235aa8dSHong Zhang 
896d2daff3dSHong Zhang   PetscFunctionBegin;
897c235aa8dSHong Zhang   rk->status = TS_STEP_INCOMPLETE;
8983389c7d9SStefano Zampini 
8999566063dSJacob Faibussowitsch   PetscCall(TSGetRHSJacobian(ts, &J, &Jpre, NULL, NULL));
90048a46eb9SPierre Jolivet   if (quadts) PetscCall(TSGetRHSJacobian(quadts, &Jquad, NULL, NULL, NULL));
9012e7b7f96SHong Zhang   for (i = s - 1; i >= 0; i--) {
9026a1e4597SHong Zhang     if (tab->FSAL && i == s - 1) {
9036a1e4597SHong Zhang       /* VecsDeltaLam[nadj*s+s-1] are initialized with zeros and the values never change.*/
9046a1e4597SHong Zhang       continue;
9056a1e4597SHong Zhang     }
9063ca0882eSHong Zhang     rk->stage_time = t + h * (1.0 - c[i]);
9079566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts, Y[i], J, Jpre));
9089371c9d4SSatish Balay     if (quadts) { PetscCall(TSComputeRHSJacobian(quadts, rk->stage_time, Y[i], Jquad, Jquad)); /* get r_u^T */ }
9093389c7d9SStefano Zampini     if (ts->vecs_sensip) {
9109566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(ts, rk->stage_time, Y[i], ts->Jacprhs)); /* get f_p */
9119371c9d4SSatish Balay       if (quadts) { PetscCall(TSComputeRHSJacobianP(quadts, rk->stage_time, Y[i], quadts->Jacprhs)); /* get f_p for the quadrature */ }
91235f5172aSHong Zhang     }
9133389c7d9SStefano Zampini 
91435f5172aSHong Zhang     if (b[i]) {
91535f5172aSHong Zhang       for (j = i + 1; j < s; j++) w[j - i - 1] = A[j * s + i] / b[i]; /* coefficients for computing VecsSensiTemp */
91635f5172aSHong Zhang     } else {
91735f5172aSHong Zhang       for (j = i + 1; j < s; j++) w[j - i - 1] = A[j * s + i]; /* coefficients for computing VecsSensiTemp */
91835f5172aSHong Zhang     }
9192e7b7f96SHong Zhang 
9202e7b7f96SHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) {
9213389c7d9SStefano Zampini       /* Stage values of lambda */
92235f5172aSHong Zhang       if (b[i]) {
92335f5172aSHong Zhang         /* lambda_{n+1} + \sum_{j=i+1}^s a_{ji}/b[i]*lambda_{s,j} */
9249566063dSJacob Faibussowitsch         PetscCall(VecCopy(ts->vecs_sensi[nadj], VecsSensiTemp[nadj])); /* VecDeltaLam is an vec array of size s by numcost */
9259566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(VecsSensiTemp[nadj], s - i - 1, w, &VecsDeltaLam[nadj * s + i + 1]));
9269566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(J, VecsSensiTemp[nadj], VecsDeltaLam[nadj * s + i])); /* VecsSensiTemp will be reused by 2nd-order adjoint */
9279566063dSJacob Faibussowitsch         PetscCall(VecScale(VecsDeltaLam[nadj * s + i], -h * b[i]));
92835f5172aSHong Zhang         if (quadts) {
9299566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(Jquad, nadj, &xarr));
9309566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(VecDRDUTransCol, xarr));
9319566063dSJacob Faibussowitsch           PetscCall(VecAXPY(VecsDeltaLam[nadj * s + i], -h * b[i], VecDRDUTransCol));
9329566063dSJacob Faibussowitsch           PetscCall(VecResetArray(VecDRDUTransCol));
9339566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(Jquad, &xarr));
934cd4cee2dSHong Zhang         }
9353389c7d9SStefano Zampini       } else {
9366a1e4597SHong Zhang         /* \sum_{j=i+1}^s a_{ji}*lambda_{s,j} */
9379566063dSJacob Faibussowitsch         PetscCall(VecSet(VecsSensiTemp[nadj], 0));
9389566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(VecsSensiTemp[nadj], s - i - 1, w, &VecsDeltaLam[nadj * s + i + 1]));
9399566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(J, VecsSensiTemp[nadj], VecsDeltaLam[nadj * s + i]));
9409566063dSJacob Faibussowitsch         PetscCall(VecScale(VecsDeltaLam[nadj * s + i], -h));
9413389c7d9SStefano Zampini       }
9426497ce14SHong Zhang 
943ad8e2604SHong Zhang       /* Stage values of mu */
9446497ce14SHong Zhang       if (ts->vecs_sensip) {
94535f5172aSHong Zhang         if (b[i]) {
9469566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(ts->Jacprhs, VecsSensiTemp[nadj], VecDeltaMu));
9479566063dSJacob Faibussowitsch           PetscCall(VecScale(VecDeltaMu, -h * b[i]));
94835f5172aSHong Zhang           if (quadts) {
9499566063dSJacob Faibussowitsch             PetscCall(MatDenseGetColumn(quadts->Jacprhs, nadj, &xarr));
9509566063dSJacob Faibussowitsch             PetscCall(VecPlaceArray(VecDRDPTransCol, xarr));
9519566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecDeltaMu, -h * b[i], VecDRDPTransCol));
9529566063dSJacob Faibussowitsch             PetscCall(VecResetArray(VecDRDPTransCol));
9539566063dSJacob Faibussowitsch             PetscCall(MatDenseRestoreColumn(quadts->Jacprhs, &xarr));
954493ed95fSHong Zhang           }
95535f5172aSHong Zhang         } else {
9569566063dSJacob Faibussowitsch           PetscCall(VecScale(VecDeltaMu, -h));
95735f5172aSHong Zhang         }
9589566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj], 1., VecDeltaMu)); /* update sensip for each stage */
959ad8e2604SHong Zhang       }
960c235aa8dSHong Zhang     }
96113af1a74SHong Zhang 
96213af1a74SHong Zhang     if (ts->vecs_sensi2 && ts->forward_solve) { /* 2nd-order adjoint, TLM mode has to be turned on */
96313af1a74SHong Zhang       /* Get w1 at t_{n+1} from TLM matrix */
9649566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(rk->MatsFwdStageSensip[i], 0, &xarr));
9659566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
96613af1a74SHong Zhang       /* lambda_s^T F_UU w_1 */
9679566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSHessianProductFunctionUU(ts, rk->stage_time, Y[i], VecsSensiTemp, ts->vec_sensip_col, ts->vecs_guu));
96835f5172aSHong Zhang       if (quadts) {
96935f5172aSHong Zhang         /* R_UU w_1 */
9709566063dSJacob Faibussowitsch         PetscCall(TSComputeRHSHessianProductFunctionUU(quadts, rk->stage_time, Y[i], NULL, ts->vec_sensip_col, ts->vecs_guu));
97135f5172aSHong Zhang       }
97235f5172aSHong Zhang       if (ts->vecs_sensip) {
97313af1a74SHong Zhang         /* lambda_s^T F_UP w_2 */
9749566063dSJacob Faibussowitsch         PetscCall(TSComputeRHSHessianProductFunctionUP(ts, rk->stage_time, Y[i], VecsSensiTemp, ts->vec_dir, ts->vecs_gup));
97535f5172aSHong Zhang         if (quadts) {
97635f5172aSHong Zhang           /* R_UP w_2 */
9779566063dSJacob Faibussowitsch           PetscCall(TSComputeRHSHessianProductFunctionUP(quadts, rk->stage_time, Y[i], NULL, ts->vec_sensip_col, ts->vecs_gup));
97835f5172aSHong Zhang         }
97935f5172aSHong Zhang       }
98013af1a74SHong Zhang       if (ts->vecs_sensi2p) {
98113af1a74SHong Zhang         /* lambda_s^T F_PU w_1 */
9829566063dSJacob Faibussowitsch         PetscCall(TSComputeRHSHessianProductFunctionPU(ts, rk->stage_time, Y[i], VecsSensiTemp, ts->vec_sensip_col, ts->vecs_gpu));
98335f5172aSHong Zhang         /* lambda_s^T F_PP w_2 */
9849566063dSJacob Faibussowitsch         PetscCall(TSComputeRHSHessianProductFunctionPP(ts, rk->stage_time, Y[i], VecsSensiTemp, ts->vec_dir, ts->vecs_gpp));
98535f5172aSHong Zhang         if (b[i] && quadts) {
98635f5172aSHong Zhang           /* R_PU w_1 */
9879566063dSJacob Faibussowitsch           PetscCall(TSComputeRHSHessianProductFunctionPU(quadts, rk->stage_time, Y[i], NULL, ts->vec_sensip_col, ts->vecs_gpu));
98835f5172aSHong Zhang           /* R_PP w_2 */
9899566063dSJacob Faibussowitsch           PetscCall(TSComputeRHSHessianProductFunctionPP(quadts, rk->stage_time, Y[i], NULL, ts->vec_dir, ts->vecs_gpp));
99035f5172aSHong Zhang         }
99113af1a74SHong Zhang       }
9929566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_sensip_col));
9939566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(rk->MatsFwdStageSensip[i], &xarr));
99413af1a74SHong Zhang 
99513af1a74SHong Zhang       for (nadj = 0; nadj < ts->numcost; nadj++) {
99613af1a74SHong Zhang         /* Stage values of lambda */
99735f5172aSHong Zhang         if (b[i]) {
99835f5172aSHong Zhang           /* J_i^T*(Lambda_{n+1}+\sum_{j=i+1}^s a_{ji}/b_i*Lambda_{s,j} */
9999566063dSJacob Faibussowitsch           PetscCall(VecCopy(ts->vecs_sensi2[nadj], VecsSensi2Temp[nadj]));
10009566063dSJacob Faibussowitsch           PetscCall(VecMAXPY(VecsSensi2Temp[nadj], s - i - 1, w, &VecsDeltaLam2[nadj * s + i + 1]));
10019566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(J, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj * s + i]));
10029566063dSJacob Faibussowitsch           PetscCall(VecScale(VecsDeltaLam2[nadj * s + i], -h * b[i]));
10039566063dSJacob Faibussowitsch           PetscCall(VecAXPY(VecsDeltaLam2[nadj * s + i], -h * b[i], ts->vecs_guu[nadj]));
100448a46eb9SPierre Jolivet           if (ts->vecs_sensip) PetscCall(VecAXPY(VecsDeltaLam2[nadj * s + i], -h * b[i], ts->vecs_gup[nadj]));
100513af1a74SHong Zhang         } else {
100635f5172aSHong Zhang           /* \sum_{j=i+1}^s a_{ji}*Lambda_{s,j} */
10079566063dSJacob Faibussowitsch           PetscCall(VecSet(VecsDeltaLam2[nadj * s + i], 0));
10089566063dSJacob Faibussowitsch           PetscCall(VecMAXPY(VecsSensi2Temp[nadj], s - i - 1, w, &VecsDeltaLam2[nadj * s + i + 1]));
10099566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(J, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj * s + i]));
10109566063dSJacob Faibussowitsch           PetscCall(VecScale(VecsDeltaLam2[nadj * s + i], -h));
10119566063dSJacob Faibussowitsch           PetscCall(VecAXPY(VecsDeltaLam2[nadj * s + i], -h, ts->vecs_guu[nadj]));
101248a46eb9SPierre Jolivet           if (ts->vecs_sensip) PetscCall(VecAXPY(VecsDeltaLam2[nadj * s + i], -h, ts->vecs_gup[nadj]));
101335f5172aSHong Zhang         }
101435f5172aSHong Zhang         if (ts->vecs_sensi2p) { /* 2nd-order adjoint for parameters */
10159566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(ts->Jacprhs, VecsSensi2Temp[nadj], VecDeltaMu2));
101635f5172aSHong Zhang           if (b[i]) {
10179566063dSJacob Faibussowitsch             PetscCall(VecScale(VecDeltaMu2, -h * b[i]));
10189566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecDeltaMu2, -h * b[i], ts->vecs_gpu[nadj]));
10199566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecDeltaMu2, -h * b[i], ts->vecs_gpp[nadj]));
102013af1a74SHong Zhang           } else {
10219566063dSJacob Faibussowitsch             PetscCall(VecScale(VecDeltaMu2, -h));
10229566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecDeltaMu2, -h, ts->vecs_gpu[nadj]));
10239566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecDeltaMu2, -h, ts->vecs_gpp[nadj]));
102413af1a74SHong Zhang           }
10259566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], 1, VecDeltaMu2)); /* update sensi2p for each stage */
102613af1a74SHong Zhang         }
102713af1a74SHong Zhang       }
102813af1a74SHong Zhang     }
10296497ce14SHong Zhang   }
1030c235aa8dSHong Zhang 
1031c235aa8dSHong Zhang   for (j = 0; j < s; j++) w[j] = 1.0;
10322e7b7f96SHong Zhang   for (nadj = 0; nadj < ts->numcost; nadj++) { /* no need to do this for mu's */
10339566063dSJacob Faibussowitsch     PetscCall(VecMAXPY(ts->vecs_sensi[nadj], s, w, &VecsDeltaLam[nadj * s]));
103448a46eb9SPierre Jolivet     if (ts->vecs_sensi2) PetscCall(VecMAXPY(ts->vecs_sensi2[nadj], s, w, &VecsDeltaLam2[nadj * s]));
10356497ce14SHong Zhang   }
1036c235aa8dSHong Zhang   rk->status = TS_STEP_COMPLETE;
10373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1038d2daff3dSHong Zhang }
1039d2daff3dSHong Zhang 
1040d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointReset_RK(TS ts)
1041d71ae5a4SJacob Faibussowitsch {
104213af1a74SHong Zhang   TS_RK    *rk  = (TS_RK *)ts->data;
104313af1a74SHong Zhang   RKTableau tab = rk->tableau;
104413af1a74SHong Zhang 
104513af1a74SHong Zhang   PetscFunctionBegin;
10469566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(tab->s * ts->numcost, &rk->VecsDeltaLam));
10479566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &rk->VecsSensiTemp));
10489566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&rk->VecDeltaMu));
10499566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(tab->s * ts->numcost, &rk->VecsDeltaLam2));
10509566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&rk->VecDeltaMu2));
10519566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &rk->VecsSensi2Temp));
10523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
105313af1a74SHong Zhang }
105413af1a74SHong Zhang 
1055d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_RK(TS ts, PetscReal itime, Vec X)
1056d71ae5a4SJacob Faibussowitsch {
105755de54fcSHong Zhang   TS_RK           *rk = (TS_RK *)ts->data;
105855de54fcSHong Zhang   PetscInt         s = rk->tableau->s, p = rk->tableau->p, i, j;
105955de54fcSHong Zhang   PetscReal        h;
106055de54fcSHong Zhang   PetscReal        tt, t;
106155de54fcSHong Zhang   PetscScalar     *b;
106255de54fcSHong Zhang   const PetscReal *B = rk->tableau->binterp;
106355de54fcSHong Zhang 
106455de54fcSHong Zhang   PetscFunctionBegin;
10653c633725SBarry Smith   PetscCheck(B, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSRK %s does not have an interpolation formula", rk->tableau->name);
106655de54fcSHong Zhang 
106755de54fcSHong Zhang   switch (rk->status) {
106855de54fcSHong Zhang   case TS_STEP_INCOMPLETE:
106955de54fcSHong Zhang   case TS_STEP_PENDING:
107055de54fcSHong Zhang     h = ts->time_step;
107155de54fcSHong Zhang     t = (itime - ts->ptime) / h;
107255de54fcSHong Zhang     break;
107355de54fcSHong Zhang   case TS_STEP_COMPLETE:
107455de54fcSHong Zhang     h = ts->ptime - ts->ptime_prev;
107555de54fcSHong Zhang     t = (itime - ts->ptime) / h + 1; /* In the interval [0,1] */
107655de54fcSHong Zhang     break;
1077d71ae5a4SJacob Faibussowitsch   default:
1078d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus");
107955de54fcSHong Zhang   }
10809566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s, &b));
108155de54fcSHong Zhang   for (i = 0; i < s; i++) b[i] = 0;
108255de54fcSHong Zhang   for (j = 0, tt = t; j < p; j++, tt *= t) {
1083ad540459SPierre Jolivet     for (i = 0; i < s; i++) b[i] += h * B[i * p + j] * tt;
108455de54fcSHong Zhang   }
10859566063dSJacob Faibussowitsch   PetscCall(VecCopy(rk->Y[0], X));
10869566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(X, s, b, rk->YdotRHS));
10879566063dSJacob Faibussowitsch   PetscCall(PetscFree(b));
10883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
108955de54fcSHong Zhang }
109055de54fcSHong Zhang 
109155de54fcSHong Zhang /*------------------------------------------------------------*/
109255de54fcSHong Zhang 
1093d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRKTableauReset(TS ts)
1094d71ae5a4SJacob Faibussowitsch {
1095be5899b3SLisandro Dalcin   TS_RK    *rk  = (TS_RK *)ts->data;
1096be5899b3SLisandro Dalcin   RKTableau tab = rk->tableau;
1097be5899b3SLisandro Dalcin 
1098be5899b3SLisandro Dalcin   PetscFunctionBegin;
10993ba16761SJacob Faibussowitsch   if (!tab) PetscFunctionReturn(PETSC_SUCCESS);
11009566063dSJacob Faibussowitsch   PetscCall(PetscFree(rk->work));
11019566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(tab->s, &rk->Y));
11029566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(tab->s, &rk->YdotRHS));
11033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1104be5899b3SLisandro Dalcin }
1105be5899b3SLisandro Dalcin 
1106d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_RK(TS ts)
1107d71ae5a4SJacob Faibussowitsch {
1108e4dd521cSBarry Smith   PetscFunctionBegin;
11099566063dSJacob Faibussowitsch   PetscCall(TSRKTableauReset(ts));
11100fe4d17eSHong Zhang   if (ts->use_splitrhsfunction) {
1111cac4c232SBarry Smith     PetscTryMethod(ts, "TSReset_RK_MultirateSplit_C", (TS), (ts));
11120fe4d17eSHong Zhang   } else {
1113cac4c232SBarry Smith     PetscTryMethod(ts, "TSReset_RK_MultirateNonsplit_C", (TS), (ts));
11140fe4d17eSHong Zhang   }
11153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1116e4dd521cSBarry Smith }
1117277b19d0SLisandro Dalcin 
1118d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCoarsenHook_TSRK(DM fine, DM coarse, void *ctx)
1119d71ae5a4SJacob Faibussowitsch {
1120f68a32c8SEmil Constantinescu   PetscFunctionBegin;
11213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1122f68a32c8SEmil Constantinescu }
1123f68a32c8SEmil Constantinescu 
1124d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRestrictHook_TSRK(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx)
1125d71ae5a4SJacob Faibussowitsch {
1126f68a32c8SEmil Constantinescu   PetscFunctionBegin;
11273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1128f68a32c8SEmil Constantinescu }
1129f68a32c8SEmil Constantinescu 
1130d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainHook_TSRK(DM dm, DM subdm, void *ctx)
1131d71ae5a4SJacob Faibussowitsch {
1132f68a32c8SEmil Constantinescu   PetscFunctionBegin;
11333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1134f68a32c8SEmil Constantinescu }
1135f68a32c8SEmil Constantinescu 
1136d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx)
1137d71ae5a4SJacob Faibussowitsch {
1138f68a32c8SEmil Constantinescu   PetscFunctionBegin;
11393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1140f68a32c8SEmil Constantinescu }
1141be5899b3SLisandro Dalcin 
1142d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRKTableauSetUp(TS ts)
1143d71ae5a4SJacob Faibussowitsch {
1144be5899b3SLisandro Dalcin   TS_RK    *rk  = (TS_RK *)ts->data;
1145be5899b3SLisandro Dalcin   RKTableau tab = rk->tableau;
1146be5899b3SLisandro Dalcin 
1147be5899b3SLisandro Dalcin   PetscFunctionBegin;
11489566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tab->s, &rk->work));
11499566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &rk->Y));
11509566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &rk->YdotRHS));
11513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1152be5899b3SLisandro Dalcin }
1153be5899b3SLisandro Dalcin 
1154d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_RK(TS ts)
1155d71ae5a4SJacob Faibussowitsch {
1156cd4cee2dSHong Zhang   TS quadts = ts->quadraturets;
1157f68a32c8SEmil Constantinescu   DM dm;
1158f68a32c8SEmil Constantinescu 
1159f68a32c8SEmil Constantinescu   PetscFunctionBegin;
11609566063dSJacob Faibussowitsch   PetscCall(TSCheckImplicitTerm(ts));
11619566063dSJacob Faibussowitsch   PetscCall(TSRKTableauSetUp(ts));
1162cd4cee2dSHong Zhang   if (quadts && ts->costintegralfwd) {
1163cd4cee2dSHong Zhang     Mat Jquad;
11649566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobian(quadts, &Jquad, NULL, NULL, NULL));
11650f7a1166SHong Zhang   }
11669566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
11679566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSRK, DMRestrictHook_TSRK, ts));
11689566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSRK, DMSubDomainRestrictHook_TSRK, ts));
11690fe4d17eSHong Zhang   if (ts->use_splitrhsfunction) {
1170cac4c232SBarry Smith     PetscTryMethod(ts, "TSSetUp_RK_MultirateSplit_C", (TS), (ts));
11710fe4d17eSHong Zhang   } else {
1172cac4c232SBarry Smith     PetscTryMethod(ts, "TSSetUp_RK_MultirateNonsplit_C", (TS), (ts));
11730fe4d17eSHong Zhang   }
11743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1175e4dd521cSBarry Smith }
1176d2daff3dSHong Zhang 
1177d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_RK(TS ts, PetscOptionItems *PetscOptionsObject)
1178d71ae5a4SJacob Faibussowitsch {
1179be5899b3SLisandro Dalcin   TS_RK *rk = (TS_RK *)ts->data;
1180262ff7bbSBarry Smith 
1181e4dd521cSBarry Smith   PetscFunctionBegin;
1182d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "RK ODE solver options");
1183f68a32c8SEmil Constantinescu   {
1184f68a32c8SEmil Constantinescu     RKTableauLink link;
1185f68a32c8SEmil Constantinescu     PetscInt      count, choice;
11860fe4d17eSHong Zhang     PetscBool     flg, use_multirate = PETSC_FALSE;
1187f68a32c8SEmil Constantinescu     const char  **namelist;
11882c9cb062Sluzhanghpp 
11899371c9d4SSatish Balay     for (link = RKTableauList, count = 0; link; link = link->next, count++)
11909371c9d4SSatish Balay       ;
11919566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(count, (char ***)&namelist));
1192f68a32c8SEmil Constantinescu     for (link = RKTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name;
11939566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ts_rk_multirate", "Use interpolation-based multirate RK method", "TSRKSetMultirate", rk->use_multirate, &use_multirate, &flg));
11941baa6e33SBarry Smith     if (flg) PetscCall(TSRKSetMultirate(ts, use_multirate));
11959566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEList("-ts_rk_type", "Family of RK method", "TSRKSetType", (const char *const *)namelist, count, rk->tableau->name, &choice, &flg));
11969566063dSJacob Faibussowitsch     if (flg) PetscCall(TSRKSetType(ts, namelist[choice]));
11979566063dSJacob Faibussowitsch     PetscCall(PetscFree(namelist));
1198f68a32c8SEmil Constantinescu   }
1199d0609cedSBarry Smith   PetscOptionsHeadEnd();
1200d0609cedSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)ts), NULL, "Multirate methods options", "");
12019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-ts_rk_dtratio", "time step ratio between slow and fast", "", rk->dtratio, &rk->dtratio, NULL));
1202d0609cedSBarry Smith   PetscOptionsEnd();
12033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1204e4dd521cSBarry Smith }
1205e4dd521cSBarry Smith 
1206d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_RK(TS ts, PetscViewer viewer)
1207d71ae5a4SJacob Faibussowitsch {
12085f70b478SJed Brown   TS_RK    *rk = (TS_RK *)ts->data;
12098370ee3bSLisandro Dalcin   PetscBool iascii;
12102cdf8120Sasbjorn 
1211e4dd521cSBarry Smith   PetscFunctionBegin;
12129566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
12138370ee3bSLisandro Dalcin   if (iascii) {
12149c334d8fSLisandro Dalcin     RKTableau        tab = rk->tableau;
1215f68a32c8SEmil Constantinescu     TSRKType         rktype;
12160f7a28cdSStefano Zampini     const PetscReal *c;
12170f7a28cdSStefano Zampini     PetscInt         s;
1218f68a32c8SEmil Constantinescu     char             buf[512];
12190f7a28cdSStefano Zampini     PetscBool        FSAL;
12200f7a28cdSStefano Zampini 
12219566063dSJacob Faibussowitsch     PetscCall(TSRKGetType(ts, &rktype));
12229566063dSJacob Faibussowitsch     PetscCall(TSRKGetTableau(ts, &s, NULL, NULL, &c, NULL, NULL, NULL, &FSAL));
12239566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  RK type %s\n", rktype));
122463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Order: %" PetscInt_FMT "\n", tab->order));
12259566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  FSAL property: %s\n", FSAL ? "yes" : "no"));
12269566063dSJacob Faibussowitsch     PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", s, c));
12279566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Abscissa c = %s\n", buf));
12288370ee3bSLisandro Dalcin   }
12293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1230f68a32c8SEmil Constantinescu }
1231f68a32c8SEmil Constantinescu 
1232d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSLoad_RK(TS ts, PetscViewer viewer)
1233d71ae5a4SJacob Faibussowitsch {
12349c334d8fSLisandro Dalcin   TSAdapt adapt;
1235f68a32c8SEmil Constantinescu 
1236f68a32c8SEmil Constantinescu   PetscFunctionBegin;
12379566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &adapt));
12389566063dSJacob Faibussowitsch   PetscCall(TSAdaptLoad(adapt, viewer));
12393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1240f68a32c8SEmil Constantinescu }
1241f68a32c8SEmil Constantinescu 
12422ea87ba9SLisandro Dalcin /*@
1243bcf0153eSBarry Smith   TSRKGetOrder - Get the order of the `TSRK` scheme
12442ea87ba9SLisandro Dalcin 
124520f4b53cSBarry Smith   Not Collective
12462ea87ba9SLisandro Dalcin 
12472ea87ba9SLisandro Dalcin   Input Parameter:
12482ea87ba9SLisandro Dalcin .  ts - timestepping context
12492ea87ba9SLisandro Dalcin 
12502ea87ba9SLisandro Dalcin   Output Parameter:
1251bcf0153eSBarry Smith .  order - order of `TSRK` scheme
12522ea87ba9SLisandro Dalcin 
12532ea87ba9SLisandro Dalcin   Level: intermediate
12542ea87ba9SLisandro Dalcin 
1255*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKGetType()`
12562ea87ba9SLisandro Dalcin @*/
1257d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKGetOrder(TS ts, PetscInt *order)
1258d71ae5a4SJacob Faibussowitsch {
12592ea87ba9SLisandro Dalcin   PetscFunctionBegin;
12602ea87ba9SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
12612ea87ba9SLisandro Dalcin   PetscValidIntPointer(order, 2);
1262cac4c232SBarry Smith   PetscUseMethod(ts, "TSRKGetOrder_C", (TS, PetscInt *), (ts, order));
12633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12642ea87ba9SLisandro Dalcin }
12652ea87ba9SLisandro Dalcin 
1266f68a32c8SEmil Constantinescu /*@C
1267bcf0153eSBarry Smith   TSRKSetType - Set the type of the `TSRK` scheme
1268f68a32c8SEmil Constantinescu 
126920f4b53cSBarry Smith   Logically Collective
1270f68a32c8SEmil Constantinescu 
1271d8d19677SJose E. Roman   Input Parameters:
1272f68a32c8SEmil Constantinescu +  ts - timestepping context
1273bcf0153eSBarry Smith -  rktype - type of `TSRK` scheme
1274f68a32c8SEmil Constantinescu 
1275bcf0153eSBarry Smith   Options Database Key:
12769bd3e852SBarry Smith .   -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs>
1277287c2655SBarry Smith 
1278f68a32c8SEmil Constantinescu   Level: intermediate
1279f68a32c8SEmil Constantinescu 
1280*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSRKGetType()`, `TSRK`, `TSRKType`, `TSRK1FE`, `TSRK2A`, `TSRK2B`, `TSRK3`, `TSRK3BS`, `TSRK4`, `TSRK5F`, `TSRK5DP`, `TSRK5BS`, `TSRK6VR`, `TSRK7VR`, `TSRK8VR`
1281f68a32c8SEmil Constantinescu @*/
1282d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKSetType(TS ts, TSRKType rktype)
1283d71ae5a4SJacob Faibussowitsch {
1284f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1285f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1286b92453a8SLisandro Dalcin   PetscValidCharPointer(rktype, 2);
1287cac4c232SBarry Smith   PetscTryMethod(ts, "TSRKSetType_C", (TS, TSRKType), (ts, rktype));
12883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1289f68a32c8SEmil Constantinescu }
1290f68a32c8SEmil Constantinescu 
1291f68a32c8SEmil Constantinescu /*@C
1292bcf0153eSBarry Smith   TSRKGetType - Get the type of `TSRK` scheme
1293f68a32c8SEmil Constantinescu 
129420f4b53cSBarry Smith   Not Collective
1295f68a32c8SEmil Constantinescu 
1296f68a32c8SEmil Constantinescu   Input Parameter:
1297f68a32c8SEmil Constantinescu .  ts - timestepping context
1298f68a32c8SEmil Constantinescu 
1299f68a32c8SEmil Constantinescu   Output Parameter:
1300bcf0153eSBarry Smith .  rktype - type of `TSRK`-scheme
1301f68a32c8SEmil Constantinescu 
1302f68a32c8SEmil Constantinescu   Level: intermediate
1303f68a32c8SEmil Constantinescu 
1304*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSRKSetType()`
1305f68a32c8SEmil Constantinescu @*/
1306d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKGetType(TS ts, TSRKType *rktype)
1307d71ae5a4SJacob Faibussowitsch {
1308f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1309f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1310cac4c232SBarry Smith   PetscUseMethod(ts, "TSRKGetType_C", (TS, TSRKType *), (ts, rktype));
13113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1312f68a32c8SEmil Constantinescu }
1313f68a32c8SEmil Constantinescu 
1314d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRKGetOrder_RK(TS ts, PetscInt *order)
1315d71ae5a4SJacob Faibussowitsch {
13162ea87ba9SLisandro Dalcin   TS_RK *rk = (TS_RK *)ts->data;
13172ea87ba9SLisandro Dalcin 
13182ea87ba9SLisandro Dalcin   PetscFunctionBegin;
13192ea87ba9SLisandro Dalcin   *order = rk->tableau->order;
13203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13212ea87ba9SLisandro Dalcin }
13222ea87ba9SLisandro Dalcin 
1323d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRKGetType_RK(TS ts, TSRKType *rktype)
1324d71ae5a4SJacob Faibussowitsch {
1325f68a32c8SEmil Constantinescu   TS_RK *rk = (TS_RK *)ts->data;
1326f68a32c8SEmil Constantinescu 
1327f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1328f68a32c8SEmil Constantinescu   *rktype = rk->tableau->name;
13293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1330f68a32c8SEmil Constantinescu }
13312c9cb062Sluzhanghpp 
1332d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRKSetType_RK(TS ts, TSRKType rktype)
1333d71ae5a4SJacob Faibussowitsch {
1334f68a32c8SEmil Constantinescu   TS_RK        *rk = (TS_RK *)ts->data;
1335f68a32c8SEmil Constantinescu   PetscBool     match;
1336f68a32c8SEmil Constantinescu   RKTableauLink link;
1337f68a32c8SEmil Constantinescu 
1338f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1339f68a32c8SEmil Constantinescu   if (rk->tableau) {
13409566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(rk->tableau->name, rktype, &match));
13413ba16761SJacob Faibussowitsch     if (match) PetscFunctionReturn(PETSC_SUCCESS);
1342f68a32c8SEmil Constantinescu   }
1343f68a32c8SEmil Constantinescu   for (link = RKTableauList; link; link = link->next) {
13449566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(link->tab.name, rktype, &match));
1345f68a32c8SEmil Constantinescu     if (match) {
13469566063dSJacob Faibussowitsch       if (ts->setupcalled) PetscCall(TSRKTableauReset(ts));
1347f68a32c8SEmil Constantinescu       rk->tableau = &link->tab;
13489566063dSJacob Faibussowitsch       if (ts->setupcalled) PetscCall(TSRKTableauSetUp(ts));
13492ffb9264SLisandro Dalcin       ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
13503ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
1351f68a32c8SEmil Constantinescu     }
1352f68a32c8SEmil Constantinescu   }
135398921bdaSJacob Faibussowitsch   SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", rktype);
1354e4dd521cSBarry Smith }
1355e4dd521cSBarry Smith 
1356d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGetStages_RK(TS ts, PetscInt *ns, Vec **Y)
1357d71ae5a4SJacob Faibussowitsch {
1358ff22ae23SHong Zhang   TS_RK *rk = (TS_RK *)ts->data;
1359ff22ae23SHong Zhang 
1360ff22ae23SHong Zhang   PetscFunctionBegin;
13610429704eSStefano Zampini   if (ns) *ns = rk->tableau->s;
1362d2daff3dSHong Zhang   if (Y) *Y = rk->Y;
13633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1364ff22ae23SHong Zhang }
1365ff22ae23SHong Zhang 
1366d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_RK(TS ts)
1367d71ae5a4SJacob Faibussowitsch {
1368b3a6b972SJed Brown   PetscFunctionBegin;
13699566063dSJacob Faibussowitsch   PetscCall(TSReset_RK(ts));
1370b3a6b972SJed Brown   if (ts->dm) {
13719566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSRK, DMRestrictHook_TSRK, ts));
13729566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSRK, DMSubDomainRestrictHook_TSRK, ts));
1373b3a6b972SJed Brown   }
13749566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
13759566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetOrder_C", NULL));
13769566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetType_C", NULL));
13779566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKSetType_C", NULL));
13789566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetTableau_C", NULL));
13799566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKSetMultirate_C", NULL));
13809566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetMultirate_C", NULL));
13812e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSetUp_RK_MultirateSplit_C", NULL));
13822e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSReset_RK_MultirateSplit_C", NULL));
13832e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSetUp_RK_MultirateNonsplit_C", NULL));
13842e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSReset_RK_MultirateNonsplit_C", NULL));
13853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1386b3a6b972SJed Brown }
1387ff22ae23SHong Zhang 
13883ca0882eSHong Zhang /*
13893ca0882eSHong Zhang   This defines the nonlinear equation that is to be solved with SNES
13903ca0882eSHong Zhang   We do not need to solve the equation; we just use SNES to approximate the Jacobian
13913ca0882eSHong Zhang */
1392d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_RK(SNES snes, Vec x, Vec y, TS ts)
1393d71ae5a4SJacob Faibussowitsch {
13943ca0882eSHong Zhang   TS_RK *rk = (TS_RK *)ts->data;
13953ca0882eSHong Zhang   DM     dm, dmsave;
13963ca0882eSHong Zhang 
13973ca0882eSHong Zhang   PetscFunctionBegin;
13989566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
13993ca0882eSHong Zhang   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
14003ca0882eSHong Zhang   dmsave = ts->dm;
14013ca0882eSHong Zhang   ts->dm = dm;
14029566063dSJacob Faibussowitsch   PetscCall(TSComputeRHSFunction(ts, rk->stage_time, x, y));
14033ca0882eSHong Zhang   ts->dm = dmsave;
14043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14053ca0882eSHong Zhang }
14063ca0882eSHong Zhang 
1407d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_RK(SNES snes, Vec x, Mat A, Mat B, TS ts)
1408d71ae5a4SJacob Faibussowitsch {
14093ca0882eSHong Zhang   TS_RK *rk = (TS_RK *)ts->data;
14103ca0882eSHong Zhang   DM     dm, dmsave;
14113ca0882eSHong Zhang 
14123ca0882eSHong Zhang   PetscFunctionBegin;
14139566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
14143ca0882eSHong Zhang   dmsave = ts->dm;
14153ca0882eSHong Zhang   ts->dm = dm;
14169566063dSJacob Faibussowitsch   PetscCall(TSComputeRHSJacobian(ts, rk->stage_time, x, A, B));
14173ca0882eSHong Zhang   ts->dm = dmsave;
14183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14193ca0882eSHong Zhang }
14203ca0882eSHong Zhang 
142121052c0cSHong Zhang /*@C
1422bcf0153eSBarry Smith   TSRKSetMultirate - Use the interpolation-based multirate `TSRK` method
142321052c0cSHong Zhang 
142420f4b53cSBarry Smith   Logically Collective
142521052c0cSHong Zhang 
1426d8d19677SJose E. Roman   Input Parameters:
142721052c0cSHong Zhang +  ts - timestepping context
1428bcf0153eSBarry 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
142921052c0cSHong Zhang 
1430bcf0153eSBarry Smith   Options Database Key:
143121052c0cSHong Zhang .   -ts_rk_multirate - <true,false>
143221052c0cSHong Zhang 
143321052c0cSHong Zhang   Level: intermediate
143421052c0cSHong Zhang 
1435bcf0153eSBarry Smith   Note:
1436da81f932SPierre 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).
1437bcf0153eSBarry Smith 
1438*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKGetMultirate()`
143921052c0cSHong Zhang @*/
1440d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKSetMultirate(TS ts, PetscBool use_multirate)
1441d71ae5a4SJacob Faibussowitsch {
14428a4be4dbSHong Zhang   PetscFunctionBegin;
1443cac4c232SBarry Smith   PetscTryMethod(ts, "TSRKSetMultirate_C", (TS, PetscBool), (ts, use_multirate));
14443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
144521052c0cSHong Zhang }
144621052c0cSHong Zhang 
144721052c0cSHong Zhang /*@C
1448bcf0153eSBarry Smith   TSRKGetMultirate - Gets whether to use the interpolation-based multirate `TSRK` method
144921052c0cSHong Zhang 
145020f4b53cSBarry Smith   Not Collective
145121052c0cSHong Zhang 
145221052c0cSHong Zhang   Input Parameter:
145321052c0cSHong Zhang .  ts - timestepping context
145421052c0cSHong Zhang 
145521052c0cSHong Zhang   Output Parameter:
1456bcf0153eSBarry Smith .  use_multirate - `PETSC_TRUE` if the multirate RK method is enabled, `PETSC_FALSE` otherwise
145721052c0cSHong Zhang 
145821052c0cSHong Zhang   Level: intermediate
145921052c0cSHong Zhang 
1460*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKSetMultirate()`
146121052c0cSHong Zhang @*/
1462d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKGetMultirate(TS ts, PetscBool *use_multirate)
1463d71ae5a4SJacob Faibussowitsch {
14647dbacdf2SHong Zhang   PetscFunctionBegin;
1465cac4c232SBarry Smith   PetscUseMethod(ts, "TSRKGetMultirate_C", (TS, PetscBool *), (ts, use_multirate));
14663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
146721052c0cSHong Zhang }
146821052c0cSHong Zhang 
1469ebe3b25bSBarry Smith /*MC
1470f68a32c8SEmil Constantinescu       TSRK - ODE and DAE solver using Runge-Kutta schemes
1471ebe3b25bSBarry Smith 
1472f68a32c8SEmil Constantinescu   The user should provide the right hand side of the equation
1473bcf0153eSBarry Smith   using `TSSetRHSFunction()`.
1474ebe3b25bSBarry Smith 
1475d41222bbSBarry Smith   Level: beginner
1476d41222bbSBarry Smith 
1477bcf0153eSBarry Smith   Notes:
1478bcf0153eSBarry Smith   The default is `TSRK3BS`, it can be changed with `TSRKSetType()` or -ts_rk_type
1479ebe3b25bSBarry Smith 
1480*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSRK`, `TSSetType()`, `TSRKSetType()`, `TSRKGetType()`, `TSRK2D`, `TSRK2E`, `TSRK3`,
1481bcf0153eSBarry Smith           `TSRK4`, `TSRK5`, `TSRKPRSSP2`, `TSRKBPR3`, `TSRKType`, `TSRKRegister()`, `TSRKSetMultirate()`, `TSRKGetMultirate()`, `TSType`
1482ebe3b25bSBarry Smith M*/
1483d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
1484d71ae5a4SJacob Faibussowitsch {
14851566a47fSLisandro Dalcin   TS_RK *rk;
1486e4dd521cSBarry Smith 
1487e4dd521cSBarry Smith   PetscFunctionBegin;
14889566063dSJacob Faibussowitsch   PetscCall(TSRKInitializePackage());
1489f68a32c8SEmil Constantinescu 
1490f68a32c8SEmil Constantinescu   ts->ops->reset          = TSReset_RK;
14915f70b478SJed Brown   ts->ops->destroy        = TSDestroy_RK;
14925f70b478SJed Brown   ts->ops->view           = TSView_RK;
1493f68a32c8SEmil Constantinescu   ts->ops->load           = TSLoad_RK;
1494f68a32c8SEmil Constantinescu   ts->ops->setup          = TSSetUp_RK;
1495f68a32c8SEmil Constantinescu   ts->ops->interpolate    = TSInterpolate_RK;
14962c9cb062Sluzhanghpp   ts->ops->step           = TSStep_RK;
1497f68a32c8SEmil Constantinescu   ts->ops->evaluatestep   = TSEvaluateStep_RK;
1498fecfb714SLisandro Dalcin   ts->ops->rollback       = TSRollBack_RK;
1499f68a32c8SEmil Constantinescu   ts->ops->setfromoptions = TSSetFromOptions_RK;
1500ff22ae23SHong Zhang   ts->ops->getstages      = TSGetStages_RK;
1501e4dd521cSBarry Smith 
15023ca0882eSHong Zhang   ts->ops->snesfunction    = SNESTSFormFunction_RK;
15033ca0882eSHong Zhang   ts->ops->snesjacobian    = SNESTSFormJacobian_RK;
15040f7a1166SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_RK;
150513af1a74SHong Zhang   ts->ops->adjointsetup    = TSAdjointSetUp_RK;
150613af1a74SHong Zhang   ts->ops->adjointstep     = TSAdjointStep_RK;
150713af1a74SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_RK;
15080f7a1166SHong Zhang 
150913af1a74SHong Zhang   ts->ops->forwardintegral  = TSForwardCostIntegral_RK;
1510922a638cSHong Zhang   ts->ops->forwardsetup     = TSForwardSetUp_RK;
1511922a638cSHong Zhang   ts->ops->forwardreset     = TSForwardReset_RK;
1512922a638cSHong Zhang   ts->ops->forwardstep      = TSForwardStep_RK;
1513922a638cSHong Zhang   ts->ops->forwardgetstages = TSForwardGetStages_RK;
1514922a638cSHong Zhang 
15154dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&rk));
15161566a47fSLisandro Dalcin   ts->data = (void *)rk;
1517e4dd521cSBarry Smith 
15189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetOrder_C", TSRKGetOrder_RK));
15199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetType_C", TSRKGetType_RK));
15209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKSetType_C", TSRKSetType_RK));
15219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetTableau_C", TSRKGetTableau_RK));
15229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKSetMultirate_C", TSRKSetMultirate_RK));
15239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetMultirate_C", TSRKGetMultirate_RK));
1524be5899b3SLisandro Dalcin 
15259566063dSJacob Faibussowitsch   PetscCall(TSRKSetType(ts, TSRKDefault));
152690b67129SHong Zhang   rk->dtratio = 1;
15273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1528e4dd521cSBarry Smith }
1529