xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision 9566063d113dddea24716c546802770db7481bc0)
1e4dd521cSBarry Smith /*
2474dd773SHong Zhang   Code for time stepping with the Runge-Kutta method
3f68a32c8SEmil Constantinescu 
4f68a32c8SEmil Constantinescu   Notes:
5474dd773SHong Zhang   The general system is written as
6474dd773SHong Zhang 
7474dd773SHong Zhang   Udot = F(t,U)
8474dd773SHong Zhang 
9e4dd521cSBarry Smith */
102c9cb062Sluzhanghpp 
11af0996ceSBarry Smith #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
12f68a32c8SEmil Constantinescu #include <petscdm.h>
13474dd773SHong Zhang #include <../src/ts/impls/explicit/rk/rk.h>
1421052c0cSHong Zhang #include <../src/ts/impls/explicit/rk/mrk.h>
15f68a32c8SEmil Constantinescu 
16484bcdc7SDebojyoti Ghosh static TSRKType  TSRKDefault = TSRK3BS;
17f68a32c8SEmil Constantinescu static PetscBool TSRKRegisterAllCalled;
18f68a32c8SEmil Constantinescu static PetscBool TSRKPackageInitialized;
19f68a32c8SEmil Constantinescu 
20ab8985f5SHong Zhang static RKTableauLink RKTableauList;
21ab8985f5SHong Zhang 
22f68a32c8SEmil Constantinescu /*MC
23287c2655SBarry Smith      TSRK1FE - First order forward Euler scheme.
24262ff7bbSBarry Smith 
25f68a32c8SEmil Constantinescu      This method has one stage.
26f68a32c8SEmil Constantinescu 
27287c2655SBarry Smith      Options database:
2867b8a455SSatish Balay .     -ts_rk_type 1fe - use type 1fe
29287c2655SBarry Smith 
30f68a32c8SEmil Constantinescu      Level: advanced
31f68a32c8SEmil Constantinescu 
32287c2655SBarry Smith .seealso: 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 
39287c2655SBarry Smith      Options database:
4067b8a455SSatish Balay .     -ts_rk_type 2a - use type 2a
41287c2655SBarry Smith 
42f68a32c8SEmil Constantinescu      Level: advanced
43f68a32c8SEmil Constantinescu 
44287c2655SBarry Smith .seealso: 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 
51e7685601SHong Zhang      Options database:
5267b8a455SSatish Balay .     -ts_rk_type 2b - use type 2b
53e7685601SHong Zhang 
54e7685601SHong Zhang      Level: advanced
55e7685601SHong Zhang 
56e7685601SHong Zhang .seealso: 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 
63287c2655SBarry Smith      Options database:
6467b8a455SSatish Balay .     -ts_rk_type 3 - use type 3
65287c2655SBarry Smith 
66f68a32c8SEmil Constantinescu      Level: advanced
67f68a32c8SEmil Constantinescu 
68287c2655SBarry Smith .seealso: 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 
75287c2655SBarry Smith      Options database:
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 
83287c2655SBarry Smith .seealso: 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 
90287c2655SBarry Smith      Options database:
9167b8a455SSatish Balay .     -ts_rk_type 4 - use type 4
92287c2655SBarry Smith 
93f68a32c8SEmil Constantinescu      Level: advanced
94f68a32c8SEmil Constantinescu 
95287c2655SBarry Smith .seealso: 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 
102287c2655SBarry Smith      Options database:
10367b8a455SSatish Balay .     -ts_rk_type 5f - use type 5f
104287c2655SBarry Smith 
105f68a32c8SEmil Constantinescu      Level: advanced
106f68a32c8SEmil Constantinescu 
107287c2655SBarry Smith .seealso: 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 
114287c2655SBarry Smith      Options database:
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 
122287c2655SBarry Smith .seealso: 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 
129287c2655SBarry Smith      Options database:
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 
137287c2655SBarry Smith .seealso: 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 
14437acaa02SHendrik Ranocha      Options database:
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 
15237acaa02SHendrik Ranocha .seealso: 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 
15937acaa02SHendrik Ranocha      Options database:
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 
16737acaa02SHendrik Ranocha .seealso: TSRK, TSRKType, TSRKSetType()
16837acaa02SHendrik Ranocha M*/
16937acaa02SHendrik Ranocha /*MC
17037acaa02SHendrik Ranocha      TSRK8VR - Eigth order robust Verner RK scheme with seventh order embedded method.
17137acaa02SHendrik Ranocha 
1728f3d7ee2SCarsten Uphoff      This method has thirteen stages.
17337acaa02SHendrik Ranocha 
17437acaa02SHendrik Ranocha      Options database:
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 
18237acaa02SHendrik Ranocha .seealso: TSRK, TSRKType, TSRKSetType()
18337acaa02SHendrik Ranocha M*/
184262ff7bbSBarry Smith 
185f68a32c8SEmil Constantinescu /*@C
186f68a32c8SEmil Constantinescu   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 
192f68a32c8SEmil Constantinescu .seealso:  TSRKRegisterDestroy()
193262ff7bbSBarry Smith @*/
194f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void)
195262ff7bbSBarry Smith {
196262ff7bbSBarry Smith   PetscFunctionBegin;
197f68a32c8SEmil Constantinescu   if (TSRKRegisterAllCalled) PetscFunctionReturn(0);
198f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_TRUE;
199e4dd521cSBarry Smith 
2004e82626cSLisandro Dalcin #define RC PetscRealConstant
201e4dd521cSBarry Smith   {
202f68a32c8SEmil Constantinescu     const PetscReal
2034e82626cSLisandro Dalcin       A[1][1] = {{0}},
2044e82626cSLisandro Dalcin       b[1]    = {RC(1.0)};
205*9566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL));
206e4dd521cSBarry Smith   }
207db046809SBarry Smith   {
208f68a32c8SEmil Constantinescu     const PetscReal
2094e82626cSLisandro Dalcin       A[2][2]   = {{0,0},
2104e82626cSLisandro Dalcin                    {RC(1.0),0}},
2114e82626cSLisandro Dalcin       b[2]      =  {RC(0.5),RC(0.5)},
2124e82626cSLisandro Dalcin       bembed[2] =  {RC(1.0),0};
213*9566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL));
214db046809SBarry Smith   }
215f68a32c8SEmil Constantinescu   {
216f68a32c8SEmil Constantinescu     const PetscReal
217e7685601SHong Zhang       A[2][2]   = {{0,0},
218e7685601SHong Zhang                    {RC(0.5),0}},
2199958aef7SJacob Faibussowitsch       b[2]      =  {0,RC(1.0)};
220*9566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK2B,2,2,&A[0][0],b,NULL,NULL,0,NULL));
221e7685601SHong Zhang   }
222e7685601SHong Zhang   {
223e7685601SHong Zhang     const PetscReal
22417f6384fSBarry Smith       A[3][3] = {{0,0,0},
2254e82626cSLisandro Dalcin                  {RC(2.0)/RC(3.0),0,0},
2264e82626cSLisandro Dalcin                  {RC(-1.0)/RC(3.0),RC(1.0),0}},
2274e82626cSLisandro Dalcin       b[3]    =  {RC(0.25),RC(0.5),RC(0.25)};
228*9566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL));
229fdefd5e5SDebojyoti Ghosh   }
230fdefd5e5SDebojyoti Ghosh   {
231fdefd5e5SDebojyoti Ghosh     const PetscReal
23217f6384fSBarry Smith       A[4][4]   = {{0,0,0,0},
2334e82626cSLisandro Dalcin                    {RC(1.0)/RC(2.0),0,0,0},
2344e82626cSLisandro Dalcin                    {0,RC(3.0)/RC(4.0),0,0},
2354e82626cSLisandro Dalcin                    {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}},
2364e82626cSLisandro Dalcin       b[4]      =  {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0},
2374e82626cSLisandro Dalcin       bembed[4] =  {RC(7.0)/RC(24.0),RC(1.0)/RC(4.0),RC(1.0)/RC(3.0),RC(1.0)/RC(8.0)};
238*9566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL));
239db046809SBarry Smith   }
240f68a32c8SEmil Constantinescu   {
241f68a32c8SEmil Constantinescu     const PetscReal
242f68a32c8SEmil Constantinescu       A[4][4] = {{0,0,0,0},
2434e82626cSLisandro Dalcin                  {RC(0.5),0,0,0},
2444e82626cSLisandro Dalcin                  {0,RC(0.5),0,0},
2454e82626cSLisandro Dalcin                  {0,0,RC(1.0),0}},
2464e82626cSLisandro 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)};
247*9566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL));
248db046809SBarry Smith   }
249f68a32c8SEmil Constantinescu   {
250f68a32c8SEmil Constantinescu     const PetscReal
25117f6384fSBarry Smith       A[6][6]   = {{0,0,0,0,0,0},
2524e82626cSLisandro Dalcin                    {RC(0.25),0,0,0,0,0},
2534e82626cSLisandro Dalcin                    {RC(3.0)/RC(32.0),RC(9.0)/RC(32.0),0,0,0,0},
2544e82626cSLisandro Dalcin                    {RC(1932.0)/RC(2197.0),RC(-7200.0)/RC(2197.0),RC(7296.0)/RC(2197.0),0,0,0},
2554e82626cSLisandro Dalcin                    {RC(439.0)/RC(216.0),RC(-8.0),RC(3680.0)/RC(513.0),RC(-845.0)/RC(4104.0),0,0},
2564e82626cSLisandro Dalcin                    {RC(-8.0)/RC(27.0),RC(2.0),RC(-3544.0)/RC(2565.0),RC(1859.0)/RC(4104.0),RC(-11.0)/RC(40.0),0}},
2574e82626cSLisandro 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)},
2584e82626cSLisandro 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};
259*9566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL));
260fdefd5e5SDebojyoti Ghosh   }
261fdefd5e5SDebojyoti Ghosh   {
262fdefd5e5SDebojyoti Ghosh     const PetscReal
26317f6384fSBarry Smith       A[7][7]       = {{0,0,0,0,0,0,0},
2644e82626cSLisandro Dalcin                        {RC(1.0)/RC(5.0),0,0,0,0,0,0},
2654e82626cSLisandro Dalcin                        {RC(3.0)/RC(40.0),RC(9.0)/RC(40.0),0,0,0,0,0},
2664e82626cSLisandro Dalcin                        {RC(44.0)/RC(45.0),RC(-56.0)/RC(15.0),RC(32.0)/RC(9.0),0,0,0,0},
2674e82626cSLisandro 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},
2684e82626cSLisandro 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},
2694e82626cSLisandro Dalcin                        {RC(35.0)/RC(384.0),0,RC(500.0)/RC(1113.0),RC(125.0)/RC(192.0),RC(-2187.0)/RC(6784.0),RC(11.0)/RC(84.0),0}},
2704e82626cSLisandro 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},
271a685a763Sluzhanghpp       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)},
272a685a763Sluzhanghpp       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)},
273a685a763Sluzhanghpp                        {0,0,0,0,0},
274a685a763Sluzhanghpp                        {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)},
275a685a763Sluzhanghpp                        {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)},
276a685a763Sluzhanghpp                        {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)},
277a685a763Sluzhanghpp                        {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)},
278a685a763Sluzhanghpp                        {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)}};
279*9566063dSJacob Faibussowitsch       PetscCall(TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,5,binterp[0]));
280f68a32c8SEmil Constantinescu   }
28105e23783SLisandro Dalcin   {
28205e23783SLisandro Dalcin     const PetscReal
28317f6384fSBarry Smith       A[8][8]   = {{0,0,0,0,0,0,0,0},
2844e82626cSLisandro Dalcin                    {RC(1.0)/RC(6.0),0,0,0,0,0,0,0},
2854e82626cSLisandro Dalcin                    {RC(2.0)/RC(27.0),RC(4.0)/RC(27.0),0,0,0,0,0,0},
2864e82626cSLisandro Dalcin                    {RC(183.0)/RC(1372.0),RC(-162.0)/RC(343.0),RC(1053.0)/RC(1372.0),0,0,0,0,0},
2874e82626cSLisandro 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},
2884e82626cSLisandro 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},
2894e82626cSLisandro 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},
2904e82626cSLisandro Dalcin                    {RC(587.0)/RC(8064.0),0,RC(4440339.0)/RC(15491840.0),RC(24353.0)/RC(124800.0),RC(387.0)/RC(44800.0),RC(2152.0)/RC(5985.0),RC(7267.0)/RC(94080.0),0}},
2914e82626cSLisandro 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},
2924e82626cSLisandro 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)};
293*9566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL));
29405e23783SLisandro Dalcin   }
29537acaa02SHendrik Ranocha   {
29637acaa02SHendrik Ranocha     const PetscReal
29737acaa02SHendrik Ranocha       A[9][9]   = {{0,0,0,0,0,0,0,0,0},
29863fe90e8SHendrik Ranocha                    {RC(1.8000000000000000000000000000000000000000e-01),0,0,0,0,0,0,0,0},
29963fe90e8SHendrik Ranocha                    {RC(8.9506172839506172839506172839506172839506e-02),RC(7.7160493827160493827160493827160493827160e-02),0,0,0,0,0,0,0},
30063fe90e8SHendrik Ranocha                    {RC(6.2500000000000000000000000000000000000000e-02),0,RC(1.8750000000000000000000000000000000000000e-01),0,0,0,0,0,0},
30163fe90e8SHendrik Ranocha                    {RC(3.1651600000000000000000000000000000000000e-01),0,RC(-1.0449480000000000000000000000000000000000e+00),RC(1.2584320000000000000000000000000000000000e+00),0,0,0,0,0},
30263fe90e8SHendrik Ranocha                    {RC(2.7232612736485626257225065566674305502508e-01),0,RC(-8.2513360323886639676113360323886639676113e-01),RC(1.0480917678812415654520917678812415654521e+00),RC(1.0471570799276856873679117969088177628396e-01),0,0,0,0},
30363fe90e8SHendrik 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},
30463fe90e8SHendrik 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},
30563fe90e8SHendrik Ranocha                    {RC(7.6388888888888888888888888888888888888889e-02),0,0,RC(3.6940836940836940836940836940836940836941e-01),0,RC(2.4801587301587301587301587301587301587302e-01),RC(2.3674242424242424242424242424242424242424e-01),RC(6.9444444444444444444444444444444444444444e-02),0}},
30663fe90e8SHendrik Ranocha       b[9]      =  {RC(7.6388888888888888888888888888888888888889e-02),0,0,RC(3.6940836940836940836940836940836940836941e-01),0,RC(2.4801587301587301587301587301587301587302e-01),RC(2.3674242424242424242424242424242424242424e-01),RC(6.9444444444444444444444444444444444444444e-02),0},
30763fe90e8SHendrik 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)};
308*9566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK6VR,6,9,&A[0][0],b,NULL,bembed,0,NULL));
30937acaa02SHendrik Ranocha   }
31037acaa02SHendrik Ranocha   {
31137acaa02SHendrik Ranocha     const PetscReal
31237acaa02SHendrik Ranocha       A[10][10]  = {{0,0,0,0,0,0,0,0,0,0},
31363fe90e8SHendrik Ranocha                     {RC(5.0000000000000000000000000000000000000000e-03),0,0,0,0,0,0,0,0,0},
31463fe90e8SHendrik Ranocha                     {RC(-1.0767901234567901234567901234567901234568e+00),RC(1.1856790123456790123456790123456790123457e+00),0,0,0,0,0,0,0,0},
31563fe90e8SHendrik Ranocha                     {RC(4.0833333333333333333333333333333333333333e-02),0,RC(1.2250000000000000000000000000000000000000e-01),0,0,0,0,0,0,0},
31663fe90e8SHendrik Ranocha                     {RC(6.3607142857142857142857142857142857142857e-01),0,RC(-2.4444642857142857142857142857142857142857e+00),RC(2.2633928571428571428571428571428571428571e+00),0,0,0,0,0,0},
31763fe90e8SHendrik Ranocha                     {RC(-2.5351211079349245229256383554660215487207e+00),0,RC(1.0299374654449267920438514460756024913612e+01),RC(-7.9513032885990579949493217458266876536482e+00),RC(7.9301148923100592201226014271115261823800e-01),0,0,0,0,0},
31863fe90e8SHendrik 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},
31963fe90e8SHendrik 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},
32063fe90e8SHendrik 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},
32163fe90e8SHendrik Ranocha                     {RC(-1.4449518916777735137351003179355712360517e+00),0,RC(8.0318913859955919224117033223019560435041e+00),RC(-7.5831741663401346820798883023671588604984e+00),RC(3.5816169353190074211247685442452878696855e+00),RC(-2.4369722632199529411183809065693752383733e+00),RC(8.5158999992326179339689766032486142173390e-01),0,0,0}},
32263fe90e8SHendrik Ranocha       b[10]      =  {RC(4.7425837833706756083569172717574534698932e-02),0,0,RC(2.5622361659370562659961727458274623448160e-01),RC(2.6951376833074206619473817258075952886764e-01),RC(1.2686622409092782845989138364739173247882e-01),RC(2.4887225942060071622046449427647492767292e-01),RC(3.0744837408200631335304388479099184768645e-03),RC(4.8023809989496943308189063347143123323209e-02),0},
32363fe90e8SHendrik Ranocha       bembed[10] =  {RC(4.7485247699299631037531273805727961552268e-02),0,0,RC(2.5599412588690633297154918245905393870497e-01),RC(2.7058478081067688722530891099268135732387e-01),RC(1.2505618684425992913638822323746917920448e-01),RC(2.5204468723743860507184043820197442562182e-01),0,0,RC(4.8834971521418614557381971303093137592592e-02)};
324*9566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK7VR,7,10,&A[0][0],b,NULL,bembed,0,NULL));
32537acaa02SHendrik Ranocha   }
32637acaa02SHendrik Ranocha   {
32737acaa02SHendrik Ranocha     const PetscReal
32837acaa02SHendrik Ranocha       A[13][13]  = {{0,0,0,0,0,0,0,0,0,0,0,0,0},
32963fe90e8SHendrik Ranocha                     {RC(2.5000000000000000000000000000000000000000e-01),0,0,0,0,0,0,0,0,0,0,0,0},
33063fe90e8SHendrik Ranocha                     {RC(8.7400846504915232052686327594877411977046e-02),RC(2.5487604938654321753087950620345685135815e-02),0,0,0,0,0,0,0,0,0,0,0},
33163fe90e8SHendrik Ranocha                     {RC(4.2333169291338582677165354330708661417323e-02),0,RC(1.2699950787401574803149606299212598425197e-01),0,0,0,0,0,0,0,0,0,0},
33263fe90e8SHendrik Ranocha                     {RC(4.2609505888742261494881445237572274090942e-01),0,RC(-1.5987952846591523265427733230657181117089e+00),RC(1.5967002257717297115939588706899953707994e+00),0,0,0,0,0,0,0,0,0},
33363fe90e8SHendrik Ranocha                     {RC(5.0719337296713929515090618138513639239329e-02),0,0,RC(2.5433377264600407582754714408877778031369e-01),RC(2.0394689005728199465736223777270858044698e-01),0,0,0,0,0,0,0,0},
33463fe90e8SHendrik 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},
33563fe90e8SHendrik 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},
33663fe90e8SHendrik 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},
33763fe90e8SHendrik 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},
33863fe90e8SHendrik 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},
33963fe90e8SHendrik 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},
34063fe90e8SHendrik Ranocha                     {RC(8.8735762208534719663211694051981022704884e-01),0,0,RC(-2.9754597821085367558513632804709301581977e+00),RC(5.6007170094881630597990392548350098923829e+00),RC(-5.9156074505366744680014930189941657351840e+00),RC(2.2029689156134927016879142540807638331238e-01),RC(1.0155097824462216666143271340902996997549e-01),RC(1.1514345647386055909780397752125850553556e+00),RC(1.9297101665271239396134361900805843653065e+00),0,0,0}},
34163fe90e8SHendrik Ranocha       b[13]      =  {RC(4.4729564666695714203015840429049382466467e-02),0,0,0,0,RC(1.5691033527708199813368698010726645409175e-01),RC(1.8460973408151637740702451873526277892035e-01),RC(2.2516380602086991042479419400350721970920e-01),RC(1.4794615651970234687005179885449141753736e-01),RC(7.6055542444955825269798361910336491012732e-02),RC(1.2277290235018619610824346315921437388535e-01),RC(4.1811958638991631583384842800871882376786e-02),0},
34263fe90e8SHendrik Ranocha       bembed[13] =  {RC(4.5847111400495925878664730122010282095875e-02),0,0,0,0,RC(2.6231891404152387437443356584845803392392e-01),RC(1.9169372337852611904485738635688429008025e-01),RC(2.1709172327902618330978407422906448568196e-01),RC(1.2738189624833706796803169450656737867900e-01),RC(1.1510530385365326258240515750043192148894e-01),0,0,RC(4.0561327798437566841823391436583608050053e-02)};
343*9566063dSJacob Faibussowitsch     PetscCall(TSRKRegister(TSRK8VR,8,13,&A[0][0],b,NULL,bembed,0,NULL));
34437acaa02SHendrik Ranocha   }
3454e82626cSLisandro Dalcin #undef RC
346db046809SBarry Smith   PetscFunctionReturn(0);
347db046809SBarry Smith }
348db046809SBarry Smith 
349f68a32c8SEmil Constantinescu /*@C
350f68a32c8SEmil Constantinescu    TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister().
351f68a32c8SEmil Constantinescu 
352f68a32c8SEmil Constantinescu    Not Collective
353f68a32c8SEmil Constantinescu 
354f68a32c8SEmil Constantinescu    Level: advanced
355f68a32c8SEmil Constantinescu 
356f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll()
357f68a32c8SEmil Constantinescu @*/
358f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void)
359e4dd521cSBarry Smith {
360f68a32c8SEmil Constantinescu   RKTableauLink  link;
361f68a32c8SEmil Constantinescu 
362f68a32c8SEmil Constantinescu   PetscFunctionBegin;
363f68a32c8SEmil Constantinescu   while ((link = RKTableauList)) {
364f68a32c8SEmil Constantinescu     RKTableau t = &link->tab;
365f68a32c8SEmil Constantinescu     RKTableauList = link->next;
366*9566063dSJacob Faibussowitsch     PetscCall(PetscFree3(t->A,t->b,t->c));
367*9566063dSJacob Faibussowitsch     PetscCall(PetscFree(t->bembed));
368*9566063dSJacob Faibussowitsch     PetscCall(PetscFree(t->binterp));
369*9566063dSJacob Faibussowitsch     PetscCall(PetscFree(t->name));
370*9566063dSJacob Faibussowitsch     PetscCall(PetscFree(link));
371f68a32c8SEmil Constantinescu   }
372f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_FALSE;
373f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
374f68a32c8SEmil Constantinescu }
375f68a32c8SEmil Constantinescu 
376f68a32c8SEmil Constantinescu /*@C
377f68a32c8SEmil Constantinescu   TSRKInitializePackage - This function initializes everything in the TSRK package. It is called
3788a690491SBarry Smith   from TSInitializePackage().
379f68a32c8SEmil Constantinescu 
380f68a32c8SEmil Constantinescu   Level: developer
381f68a32c8SEmil Constantinescu 
382f68a32c8SEmil Constantinescu .seealso: PetscInitialize()
383f68a32c8SEmil Constantinescu @*/
384f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void)
385f68a32c8SEmil Constantinescu {
386e4dd521cSBarry Smith   PetscFunctionBegin;
387f68a32c8SEmil Constantinescu   if (TSRKPackageInitialized) PetscFunctionReturn(0);
388f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_TRUE;
389*9566063dSJacob Faibussowitsch   PetscCall(TSRKRegisterAll());
390*9566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(TSRKFinalizePackage));
391f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
392f68a32c8SEmil Constantinescu }
393186e87acSLisandro Dalcin 
394f68a32c8SEmil Constantinescu /*@C
395f68a32c8SEmil Constantinescu   TSRKFinalizePackage - This function destroys everything in the TSRK package. It is
396f68a32c8SEmil Constantinescu   called from PetscFinalize().
397186e87acSLisandro Dalcin 
398f68a32c8SEmil Constantinescu   Level: developer
399f68a32c8SEmil Constantinescu 
400f68a32c8SEmil Constantinescu .seealso: PetscFinalize()
401f68a32c8SEmil Constantinescu @*/
402f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void)
403f68a32c8SEmil Constantinescu {
404f68a32c8SEmil Constantinescu   PetscFunctionBegin;
405f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_FALSE;
406*9566063dSJacob Faibussowitsch   PetscCall(TSRKRegisterDestroy());
407f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
408f68a32c8SEmil Constantinescu }
409f68a32c8SEmil Constantinescu 
410f68a32c8SEmil Constantinescu /*@C
411f68a32c8SEmil Constantinescu    TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
412f68a32c8SEmil Constantinescu 
413f68a32c8SEmil Constantinescu    Not Collective, but the same schemes should be registered on all processes on which they will be used
414f68a32c8SEmil Constantinescu 
415f68a32c8SEmil Constantinescu    Input Parameters:
416f68a32c8SEmil Constantinescu +  name - identifier for method
417f68a32c8SEmil Constantinescu .  order - approximation order of method
418f68a32c8SEmil Constantinescu .  s - number of stages, this is the dimension of the matrices below
419f68a32c8SEmil Constantinescu .  A - stage coefficients (dimension s*s, row-major)
420f68a32c8SEmil Constantinescu .  b - step completion table (dimension s; NULL to use last row of A)
421f68a32c8SEmil Constantinescu .  c - abscissa (dimension s; NULL to use row sums of A)
422f68a32c8SEmil Constantinescu .  bembed - completion table for embedded method (dimension s; NULL if not available)
4233a8a9803SLisandro Dalcin .  p - Order of the interpolation scheme, equal to the number of columns of binterp
4243a8a9803SLisandro Dalcin -  binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1)
425f68a32c8SEmil Constantinescu 
426f68a32c8SEmil Constantinescu    Notes:
427f68a32c8SEmil Constantinescu    Several RK methods are provided, this function is only needed to create new methods.
428f68a32c8SEmil Constantinescu 
429f68a32c8SEmil Constantinescu    Level: advanced
430f68a32c8SEmil Constantinescu 
431f68a32c8SEmil Constantinescu .seealso: TSRK
432f68a32c8SEmil Constantinescu @*/
433f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s,
434f68a32c8SEmil Constantinescu                             const PetscReal A[],const PetscReal b[],const PetscReal c[],
4353a8a9803SLisandro Dalcin                             const PetscReal bembed[],PetscInt p,const PetscReal binterp[])
436f68a32c8SEmil Constantinescu {
437f68a32c8SEmil Constantinescu   RKTableauLink   link;
438f68a32c8SEmil Constantinescu   RKTableau       t;
439f68a32c8SEmil Constantinescu   PetscInt        i,j;
440f68a32c8SEmil Constantinescu 
441f68a32c8SEmil Constantinescu   PetscFunctionBegin;
4423a8a9803SLisandro Dalcin   PetscValidCharPointer(name,1);
4433a8a9803SLisandro Dalcin   PetscValidRealPointer(A,4);
4443a8a9803SLisandro Dalcin   if (b) PetscValidRealPointer(b,5);
4453a8a9803SLisandro Dalcin   if (c) PetscValidRealPointer(c,6);
4463a8a9803SLisandro Dalcin   if (bembed) PetscValidRealPointer(bembed,7);
4473a8a9803SLisandro Dalcin   if (binterp || p > 1) PetscValidRealPointer(binterp,9);
4483a8a9803SLisandro Dalcin 
449*9566063dSJacob Faibussowitsch   PetscCall(TSRKInitializePackage());
450*9566063dSJacob Faibussowitsch   PetscCall(PetscNew(&link));
451f68a32c8SEmil Constantinescu   t = &link->tab;
4523a8a9803SLisandro Dalcin 
453*9566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name,&t->name));
454f68a32c8SEmil Constantinescu   t->order = order;
455f68a32c8SEmil Constantinescu   t->s = s;
456*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c));
457*9566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->A,A,s*s));
458*9566063dSJacob Faibussowitsch   if (b)  PetscCall(PetscArraycpy(t->b,b,s));
459f68a32c8SEmil Constantinescu   else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i];
460*9566063dSJacob Faibussowitsch   if (c)  PetscCall(PetscArraycpy(t->c,c,s));
461f68a32c8SEmil Constantinescu   else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j];
462d760c35bSDebojyoti Ghosh   t->FSAL = PETSC_TRUE;
463d760c35bSDebojyoti Ghosh   for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE;
4643a8a9803SLisandro Dalcin 
465f68a32c8SEmil Constantinescu   if (bembed) {
466*9566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(s,&t->bembed));
467*9566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(t->bembed,bembed,s));
468f68a32c8SEmil Constantinescu   }
469f68a32c8SEmil Constantinescu 
4703a8a9803SLisandro Dalcin   if (!binterp) { p = 1; binterp = t->b; }
4713a8a9803SLisandro Dalcin   t->p = p;
472*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s*p,&t->binterp));
473*9566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->binterp,binterp,s*p));
4743a8a9803SLisandro Dalcin 
475f68a32c8SEmil Constantinescu   link->next = RKTableauList;
476f68a32c8SEmil Constantinescu   RKTableauList = link;
477f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
478f68a32c8SEmil Constantinescu }
479f68a32c8SEmil Constantinescu 
4800f7a28cdSStefano Zampini PetscErrorCode TSRKGetTableau_RK(TS ts, PetscInt *s, const PetscReal **A, const PetscReal **b, const PetscReal **c, const PetscReal **bembed,
4810f7a28cdSStefano Zampini                                         PetscInt *p, const PetscReal **binterp, PetscBool *FSAL)
4820f7a28cdSStefano Zampini {
4830f7a28cdSStefano Zampini   TS_RK     *rk   = (TS_RK*)ts->data;
4840f7a28cdSStefano Zampini   RKTableau tab  = rk->tableau;
4850f7a28cdSStefano Zampini 
4860f7a28cdSStefano Zampini   PetscFunctionBegin;
4870f7a28cdSStefano Zampini   if (s) *s = tab->s;
4880f7a28cdSStefano Zampini   if (A) *A = tab->A;
4890f7a28cdSStefano Zampini   if (b) *b = tab->b;
4900f7a28cdSStefano Zampini   if (c) *c = tab->c;
4910f7a28cdSStefano Zampini   if (bembed) *bembed = tab->bembed;
4920f7a28cdSStefano Zampini   if (p) *p = tab->p;
4930f7a28cdSStefano Zampini   if (binterp) *binterp = tab->binterp;
4940f7a28cdSStefano Zampini   if (FSAL) *FSAL = tab->FSAL;
4950f7a28cdSStefano Zampini   PetscFunctionReturn(0);
4960f7a28cdSStefano Zampini }
4970f7a28cdSStefano Zampini 
4980f7a28cdSStefano Zampini /*@C
4990f7a28cdSStefano Zampini    TSRKGetTableau - Get info on the RK tableau
5000f7a28cdSStefano Zampini 
5010f7a28cdSStefano Zampini    Not Collective
5020f7a28cdSStefano Zampini 
503f899ff85SJose E. Roman    Input Parameter:
5040f7a28cdSStefano Zampini .  ts - timestepping context
5050f7a28cdSStefano Zampini 
5060f7a28cdSStefano Zampini    Output Parameters:
5070f7a28cdSStefano Zampini +  s - number of stages, this is the dimension of the matrices below
5080f7a28cdSStefano Zampini .  A - stage coefficients (dimension s*s, row-major)
5090f7a28cdSStefano Zampini .  b - step completion table (dimension s)
5100f7a28cdSStefano Zampini .  c - abscissa (dimension s)
5110f7a28cdSStefano Zampini .  bembed - completion table for embedded method (dimension s; NULL if not available)
5120f7a28cdSStefano Zampini .  p - Order of the interpolation scheme, equal to the number of columns of binterp
5130f7a28cdSStefano Zampini .  binterp - Coefficients of the interpolation formula (dimension s*p)
5140f7a28cdSStefano Zampini -  FSAL - wheather or not the scheme has the First Same As Last property
5150f7a28cdSStefano Zampini 
5160f7a28cdSStefano Zampini    Level: developer
5170f7a28cdSStefano Zampini 
5180f7a28cdSStefano Zampini .seealso: TSRK
5190f7a28cdSStefano Zampini @*/
5200f7a28cdSStefano Zampini PetscErrorCode TSRKGetTableau(TS ts, PetscInt *s, const PetscReal **A, const PetscReal **b, const PetscReal **c, const PetscReal **bembed,
5210f7a28cdSStefano Zampini                                      PetscInt *p, const PetscReal **binterp, PetscBool *FSAL)
5220f7a28cdSStefano Zampini {
5230f7a28cdSStefano Zampini   PetscErrorCode ierr;
5240f7a28cdSStefano Zampini 
5250f7a28cdSStefano Zampini   PetscFunctionBegin;
5260f7a28cdSStefano Zampini   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5270f7a28cdSStefano Zampini   ierr = PetscUseMethod(ts,"TSRKGetTableau_C",(TS,PetscInt*,const PetscReal**,const PetscReal**,const PetscReal**,const PetscReal**,
528*9566063dSJacob Faibussowitsch                                                   PetscInt*,const PetscReal**,PetscBool*),(ts,s,A,b,c,bembed,p,binterp,FSAL));PetscCall(ierr);
5290f7a28cdSStefano Zampini   PetscFunctionReturn(0);
5300f7a28cdSStefano Zampini }
5310f7a28cdSStefano Zampini 
532e4dd521cSBarry Smith /*
5332c9cb062Sluzhanghpp  This is for single-step RK method
534f68a32c8SEmil Constantinescu  The step completion formula is
535e4dd521cSBarry Smith 
536f68a32c8SEmil Constantinescu  x1 = x0 + h b^T YdotRHS
537f68a32c8SEmil Constantinescu 
538f68a32c8SEmil Constantinescu  This function can be called before or after ts->vec_sol has been updated.
539f68a32c8SEmil Constantinescu  Suppose we have a completion formula (b) and an embedded formula (be) of different order.
540f68a32c8SEmil Constantinescu  We can write
541f68a32c8SEmil Constantinescu 
542f68a32c8SEmil Constantinescu  x1e = x0 + h be^T YdotRHS
543f68a32c8SEmil Constantinescu      = x1 - h b^T YdotRHS + h be^T YdotRHS
544f68a32c8SEmil Constantinescu      = x1 + h (be - b)^T YdotRHS
545f68a32c8SEmil Constantinescu 
546f68a32c8SEmil Constantinescu  so we can evaluate the method with different order even after the step has been optimistically completed.
547e4dd521cSBarry Smith */
548f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done)
549f68a32c8SEmil Constantinescu {
550f68a32c8SEmil Constantinescu   TS_RK          *rk   = (TS_RK*)ts->data;
551f68a32c8SEmil Constantinescu   RKTableau      tab  = rk->tableau;
552f68a32c8SEmil Constantinescu   PetscScalar    *w    = rk->work;
553f68a32c8SEmil Constantinescu   PetscReal      h;
554f68a32c8SEmil Constantinescu   PetscInt       s    = tab->s,j;
555f68a32c8SEmil Constantinescu 
556f68a32c8SEmil Constantinescu   PetscFunctionBegin;
557f68a32c8SEmil Constantinescu   switch (rk->status) {
558f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
559f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
560f68a32c8SEmil Constantinescu     h = ts->time_step; break;
561f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
562be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
563f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
564f68a32c8SEmil Constantinescu   }
565f68a32c8SEmil Constantinescu   if (order == tab->order) {
566f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) {
567*9566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol,X));
56890b67129SHong Zhang       for (j=0; j<s; j++) w[j] = h*tab->b[j]/rk->dtratio;
569*9566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(X,s,w,rk->YdotRHS));
570*9566063dSJacob Faibussowitsch     } else PetscCall(VecCopy(ts->vec_sol,X));
571f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
572f68a32c8SEmil Constantinescu   } else if (order == tab->order-1) {
573f68a32c8SEmil Constantinescu     if (!tab->bembed) goto unavailable;
574f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) { /*Complete with the embedded method (be)*/
575*9566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol,X));
576f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
577*9566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(X,s,w,rk->YdotRHS));
578f68a32c8SEmil Constantinescu     } else {  /*Rollback and re-complete using (be-b) */
579*9566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol,X));
580f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
581*9566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(X,s,w,rk->YdotRHS));
582f68a32c8SEmil Constantinescu     }
583f68a32c8SEmil Constantinescu     if (done) *done = PETSC_TRUE;
584f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
585f68a32c8SEmil Constantinescu   }
586f68a32c8SEmil Constantinescu unavailable:
587f68a32c8SEmil Constantinescu   if (done) *done = PETSC_FALSE;
58898921bdaSJacob Faibussowitsch   else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"RK '%s' of order %D cannot evaluate step at order %D. Consider using -ts_adapt_type none or a different method that has an embedded estimate.",tab->name,tab->order,order);
589f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
590f68a32c8SEmil Constantinescu }
591f68a32c8SEmil Constantinescu 
5920f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts)
5930f7a1166SHong Zhang {
5940f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
595cd4cee2dSHong Zhang   TS              quadts = ts->quadraturets;
5960f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
5970f7a1166SHong Zhang   const PetscInt  s = tab->s;
5980f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
5990f7a1166SHong Zhang   Vec             *Y = rk->Y;
6000f7a1166SHong Zhang   PetscInt        i;
6010f7a1166SHong Zhang 
6020f7a1166SHong Zhang   PetscFunctionBegin;
603cd4cee2dSHong Zhang   /* No need to backup quadts->vec_sol since it can be reverted in TSRollBack_RK */
6040f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
605cd4cee2dSHong Zhang     /* Evolve quadrature TS solution to compute integrals */
606*9566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(quadts,rk->ptime+rk->time_step*c[i],Y[i],ts->vec_costintegrand));
607*9566063dSJacob Faibussowitsch     PetscCall(VecAXPY(quadts->vec_sol,rk->time_step*b[i],ts->vec_costintegrand));
6080f7a1166SHong Zhang   }
6090f7a1166SHong Zhang   PetscFunctionReturn(0);
6100f7a1166SHong Zhang }
6110f7a1166SHong Zhang 
6120f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts)
6130f7a1166SHong Zhang {
6140f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
6150f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
616cd4cee2dSHong Zhang   TS              quadts = ts->quadraturets;
6170f7a1166SHong Zhang   const PetscInt  s = tab->s;
6180f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
6190f7a1166SHong Zhang   Vec             *Y = rk->Y;
6200f7a1166SHong Zhang   PetscInt        i;
6210f7a1166SHong Zhang 
6220f7a1166SHong Zhang   PetscFunctionBegin;
6230f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
624cd4cee2dSHong Zhang     /* Evolve quadrature TS solution to compute integrals */
625*9566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(quadts,ts->ptime+ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand));
626*9566063dSJacob Faibussowitsch     PetscCall(VecAXPY(quadts->vec_sol,-ts->time_step*b[i],ts->vec_costintegrand));
6270f7a1166SHong Zhang   }
6280f7a1166SHong Zhang   PetscFunctionReturn(0);
6290f7a1166SHong Zhang }
6300f7a1166SHong Zhang 
631fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts)
632fecfb714SLisandro Dalcin {
633fecfb714SLisandro Dalcin   TS_RK           *rk = (TS_RK*)ts->data;
634cd4cee2dSHong Zhang   TS              quadts = ts->quadraturets;
635fecfb714SLisandro Dalcin   RKTableau       tab = rk->tableau;
636fecfb714SLisandro Dalcin   const PetscInt  s  = tab->s;
637cd4cee2dSHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
638fecfb714SLisandro Dalcin   PetscScalar     *w = rk->work;
639cd4cee2dSHong Zhang   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
640fecfb714SLisandro Dalcin   PetscInt        j;
641fecfb714SLisandro Dalcin   PetscReal       h;
642fecfb714SLisandro Dalcin 
643fecfb714SLisandro Dalcin   PetscFunctionBegin;
644fecfb714SLisandro Dalcin   switch (rk->status) {
645fecfb714SLisandro Dalcin   case TS_STEP_INCOMPLETE:
646fecfb714SLisandro Dalcin   case TS_STEP_PENDING:
647fecfb714SLisandro Dalcin     h = ts->time_step; break;
648fecfb714SLisandro Dalcin   case TS_STEP_COMPLETE:
649fecfb714SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
650fecfb714SLisandro Dalcin   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
651fecfb714SLisandro Dalcin   }
652fecfb714SLisandro Dalcin   for (j=0; j<s; j++) w[j] = -h*b[j];
653*9566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(ts->vec_sol,s,w,YdotRHS));
654cd4cee2dSHong Zhang   if (quadts && ts->costintegralfwd) {
655cd4cee2dSHong Zhang     for (j=0; j<s; j++) {
656cd4cee2dSHong Zhang       /* Revert the quadrature TS solution */
657*9566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSFunction(quadts,rk->ptime+h*c[j],Y[j],ts->vec_costintegrand));
658*9566063dSJacob Faibussowitsch       PetscCall(VecAXPY(quadts->vec_sol,-h*b[j],ts->vec_costintegrand));
659cd4cee2dSHong Zhang     }
660cd4cee2dSHong Zhang   }
661fecfb714SLisandro Dalcin   PetscFunctionReturn(0);
662fecfb714SLisandro Dalcin }
663fecfb714SLisandro Dalcin 
664922a638cSHong Zhang static PetscErrorCode TSForwardStep_RK(TS ts)
665922a638cSHong Zhang {
666922a638cSHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
667922a638cSHong Zhang   RKTableau       tab = rk->tableau;
668922a638cSHong Zhang   Mat             J,*MatsFwdSensipTemp = rk->MatsFwdSensipTemp;
669922a638cSHong Zhang   const PetscInt  s = tab->s;
670922a638cSHong Zhang   const PetscReal *A = tab->A,*c = tab->c,*b = tab->b;
671922a638cSHong Zhang   Vec             *Y = rk->Y;
672922a638cSHong Zhang   PetscInt        i,j;
673922a638cSHong Zhang   PetscReal       stage_time,h = ts->time_step;
674922a638cSHong Zhang   PetscBool       zero;
675922a638cSHong Zhang 
676922a638cSHong Zhang   PetscFunctionBegin;
677*9566063dSJacob Faibussowitsch   PetscCall(MatCopy(ts->mat_sensip,rk->MatFwdSensip0,SAME_NONZERO_PATTERN));
678*9566063dSJacob Faibussowitsch   PetscCall(TSGetRHSJacobian(ts,&J,NULL,NULL,NULL));
679922a638cSHong Zhang 
680922a638cSHong Zhang   for (i=0; i<s; i++) {
681922a638cSHong Zhang     stage_time = ts->ptime + h*c[i];
682922a638cSHong Zhang     zero = PETSC_FALSE;
683922a638cSHong Zhang     if (b[i] == 0 && i == s-1) zero = PETSC_TRUE;
684922a638cSHong Zhang     /* TLM Stage values */
685922a638cSHong Zhang     if (!i) {
686*9566063dSJacob Faibussowitsch       PetscCall(MatCopy(ts->mat_sensip,rk->MatsFwdStageSensip[i],SAME_NONZERO_PATTERN));
687922a638cSHong Zhang     } else if (!zero) {
688*9566063dSJacob Faibussowitsch       PetscCall(MatZeroEntries(rk->MatsFwdStageSensip[i]));
689922a638cSHong Zhang       for (j=0; j<i; j++) {
690*9566063dSJacob Faibussowitsch         PetscCall(MatAXPY(rk->MatsFwdStageSensip[i],h*A[i*s+j],MatsFwdSensipTemp[j],SAME_NONZERO_PATTERN));
691922a638cSHong Zhang       }
692*9566063dSJacob Faibussowitsch       PetscCall(MatAXPY(rk->MatsFwdStageSensip[i],1.,ts->mat_sensip,SAME_NONZERO_PATTERN));
693922a638cSHong Zhang     } else {
694*9566063dSJacob Faibussowitsch       PetscCall(MatZeroEntries(rk->MatsFwdStageSensip[i]));
695922a638cSHong Zhang     }
696922a638cSHong Zhang 
697*9566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSJacobian(ts,stage_time,Y[i],J,J));
698*9566063dSJacob Faibussowitsch     PetscCall(MatMatMult(J,rk->MatsFwdStageSensip[i],MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatsFwdSensipTemp[i]));
69913af1a74SHong Zhang     if (ts->Jacprhs) {
700*9566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(ts,stage_time,Y[i],ts->Jacprhs)); /* get f_p */
70113af1a74SHong Zhang       if (ts->vecs_sensi2p) { /* TLM used for 2nd-order adjoint */
70213af1a74SHong Zhang         PetscScalar *xarr;
703*9566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(MatsFwdSensipTemp[i],0,&xarr));
704*9566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(rk->VecDeltaFwdSensipCol,xarr));
705*9566063dSJacob Faibussowitsch         PetscCall(MatMultAdd(ts->Jacprhs,ts->vec_dir,rk->VecDeltaFwdSensipCol,rk->VecDeltaFwdSensipCol));
706*9566063dSJacob Faibussowitsch         PetscCall(VecResetArray(rk->VecDeltaFwdSensipCol));
707*9566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(MatsFwdSensipTemp[i],&xarr));
70813af1a74SHong Zhang       } else {
709*9566063dSJacob Faibussowitsch         PetscCall(MatAXPY(MatsFwdSensipTemp[i],1.,ts->Jacprhs,SUBSET_NONZERO_PATTERN));
71013af1a74SHong Zhang       }
711922a638cSHong Zhang     }
712922a638cSHong Zhang   }
713922a638cSHong Zhang 
714922a638cSHong Zhang   for (i=0; i<s; i++) {
715*9566063dSJacob Faibussowitsch     PetscCall(MatAXPY(ts->mat_sensip,h*b[i],rk->MatsFwdSensipTemp[i],SAME_NONZERO_PATTERN));
716922a638cSHong Zhang   }
717922a638cSHong Zhang   rk->status = TS_STEP_COMPLETE;
718922a638cSHong Zhang   PetscFunctionReturn(0);
719922a638cSHong Zhang }
720922a638cSHong Zhang 
721922a638cSHong Zhang static PetscErrorCode TSForwardGetStages_RK(TS ts,PetscInt *ns,Mat **stagesensip)
722922a638cSHong Zhang {
723922a638cSHong Zhang   TS_RK     *rk = (TS_RK*)ts->data;
724922a638cSHong Zhang   RKTableau tab  = rk->tableau;
725922a638cSHong Zhang 
726922a638cSHong Zhang   PetscFunctionBegin;
727922a638cSHong Zhang   if (ns) *ns = tab->s;
728922a638cSHong Zhang   if (stagesensip) *stagesensip = rk->MatsFwdStageSensip;
729922a638cSHong Zhang   PetscFunctionReturn(0);
730922a638cSHong Zhang }
731922a638cSHong Zhang 
732922a638cSHong Zhang static PetscErrorCode TSForwardSetUp_RK(TS ts)
733922a638cSHong Zhang {
734922a638cSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
735922a638cSHong Zhang   RKTableau      tab  = rk->tableau;
736922a638cSHong Zhang   PetscInt       i;
737922a638cSHong Zhang 
738922a638cSHong Zhang   PetscFunctionBegin;
739922a638cSHong Zhang   /* backup sensitivity results for roll-backs */
740*9566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatFwdSensip0));
741922a638cSHong Zhang 
742*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tab->s,&rk->MatsFwdStageSensip));
743*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tab->s,&rk->MatsFwdSensipTemp));
744922a638cSHong Zhang   for (i=0; i<tab->s; i++) {
745*9566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatsFwdStageSensip[i]));
746*9566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatsFwdSensipTemp[i]));
747922a638cSHong Zhang   }
748*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol,&rk->VecDeltaFwdSensipCol));
749922a638cSHong Zhang   PetscFunctionReturn(0);
750922a638cSHong Zhang }
751922a638cSHong Zhang 
752922a638cSHong Zhang static PetscErrorCode TSForwardReset_RK(TS ts)
753922a638cSHong Zhang {
754922a638cSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
755922a638cSHong Zhang   RKTableau      tab  = rk->tableau;
756922a638cSHong Zhang   PetscInt       i;
757922a638cSHong Zhang 
758922a638cSHong Zhang   PetscFunctionBegin;
759*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&rk->MatFwdSensip0));
760922a638cSHong Zhang   if (rk->MatsFwdStageSensip) {
761922a638cSHong Zhang     for (i=0; i<tab->s; i++) {
762*9566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&rk->MatsFwdStageSensip[i]));
763922a638cSHong Zhang     }
764*9566063dSJacob Faibussowitsch     PetscCall(PetscFree(rk->MatsFwdStageSensip));
765922a638cSHong Zhang   }
766922a638cSHong Zhang   if (rk->MatsFwdSensipTemp) {
767922a638cSHong Zhang     for (i=0; i<tab->s; i++) {
768*9566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&rk->MatsFwdSensipTemp[i]));
769922a638cSHong Zhang     }
770*9566063dSJacob Faibussowitsch     PetscCall(PetscFree(rk->MatsFwdSensipTemp));
771922a638cSHong Zhang   }
772*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&rk->VecDeltaFwdSensipCol));
773922a638cSHong Zhang   PetscFunctionReturn(0);
774922a638cSHong Zhang }
775922a638cSHong Zhang 
776f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts)
777f68a32c8SEmil Constantinescu {
778f68a32c8SEmil Constantinescu   TS_RK           *rk  = (TS_RK*)ts->data;
779f68a32c8SEmil Constantinescu   RKTableau       tab  = rk->tableau;
780f68a32c8SEmil Constantinescu   const PetscInt  s = tab->s;
781fecfb714SLisandro Dalcin   const PetscReal *A = tab->A,*c = tab->c;
782f68a32c8SEmil Constantinescu   PetscScalar     *w = rk->work;
783f68a32c8SEmil Constantinescu   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
784d1905564SLisandro Dalcin   PetscBool       FSAL = tab->FSAL;
785f68a32c8SEmil Constantinescu   TSAdapt         adapt;
786fecfb714SLisandro Dalcin   PetscInt        i,j;
787be5899b3SLisandro Dalcin   PetscInt        rejections = 0;
788be5899b3SLisandro Dalcin   PetscBool       stageok,accept = PETSC_TRUE;
789be5899b3SLisandro Dalcin   PetscReal       next_time_step = ts->time_step;
790f68a32c8SEmil Constantinescu 
791f68a32c8SEmil Constantinescu   PetscFunctionBegin;
792d1905564SLisandro Dalcin   if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE;
793*9566063dSJacob Faibussowitsch   if (FSAL) PetscCall(VecCopy(YdotRHS[s-1],YdotRHS[0]));
794d1905564SLisandro Dalcin 
795f68a32c8SEmil Constantinescu   rk->status = TS_STEP_INCOMPLETE;
796be5899b3SLisandro Dalcin   while (!ts->reason && rk->status != TS_STEP_COMPLETE) {
797be5899b3SLisandro Dalcin     PetscReal t = ts->ptime;
798f68a32c8SEmil Constantinescu     PetscReal h = ts->time_step;
799f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
8009be3e283SDebojyoti Ghosh       rk->stage_time = t + h*c[i];
801*9566063dSJacob Faibussowitsch       PetscCall(TSPreStage(ts,rk->stage_time));
802*9566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol,Y[i]));
803f68a32c8SEmil Constantinescu       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
804*9566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(Y[i],i,w,YdotRHS));
805*9566063dSJacob Faibussowitsch       PetscCall(TSPostStage(ts,rk->stage_time,i,Y));
806*9566063dSJacob Faibussowitsch       PetscCall(TSGetAdapt(ts,&adapt));
807*9566063dSJacob Faibussowitsch       PetscCall(TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok));
808fecfb714SLisandro Dalcin       if (!stageok) goto reject_step;
8098f738a7cSLisandro Dalcin       if (FSAL && !i) continue;
810*9566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]));
811f68a32c8SEmil Constantinescu     }
812be5899b3SLisandro Dalcin 
813be5899b3SLisandro Dalcin     rk->status = TS_STEP_INCOMPLETE;
814*9566063dSJacob Faibussowitsch     PetscCall(TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL));
815f68a32c8SEmil Constantinescu     rk->status = TS_STEP_PENDING;
816*9566063dSJacob Faibussowitsch     PetscCall(TSGetAdapt(ts,&adapt));
817*9566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidatesClear(adapt));
818*9566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE));
819*9566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept));
820be5899b3SLisandro Dalcin     rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
821be5899b3SLisandro Dalcin     if (!accept) { /* Roll back the current step */
822*9566063dSJacob Faibussowitsch       PetscCall(TSRollBack_RK(ts));
823be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
824be5899b3SLisandro Dalcin       goto reject_step;
825be5899b3SLisandro Dalcin     }
826be5899b3SLisandro Dalcin 
8270f7a1166SHong Zhang     if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
8280f7a1166SHong Zhang       rk->ptime     = ts->ptime;
8290f7a1166SHong Zhang       rk->time_step = ts->time_step;
830493ed95fSHong Zhang     }
831be5899b3SLisandro Dalcin 
832f68a32c8SEmil Constantinescu     ts->ptime += ts->time_step;
833f68a32c8SEmil Constantinescu     ts->time_step = next_time_step;
834f68a32c8SEmil Constantinescu     break;
835be5899b3SLisandro Dalcin 
836be5899b3SLisandro Dalcin     reject_step:
837fecfb714SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
838be5899b3SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
839be5899b3SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
840*9566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections));
841f68a32c8SEmil Constantinescu     }
842f68a32c8SEmil Constantinescu   }
843f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
844e4dd521cSBarry Smith }
845e4dd521cSBarry Smith 
846f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts)
847f6a906c0SBarry Smith {
848f6a906c0SBarry Smith   TS_RK          *rk  = (TS_RK*)ts->data;
849be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
850be5899b3SLisandro Dalcin   PetscInt       s   = tab->s;
851f6a906c0SBarry Smith 
852f6a906c0SBarry Smith   PetscFunctionBegin;
853f6a906c0SBarry Smith   if (ts->adjointsetupcalled++) PetscFunctionReturn(0);
854*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecsDeltaLam));
855*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&rk->VecsSensiTemp));
856f6a906c0SBarry Smith   if (ts->vecs_sensip) {
857*9566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vecs_sensip[0],&rk->VecDeltaMu));
858f6a906c0SBarry Smith   }
85913af1a74SHong Zhang   if (ts->vecs_sensi2) {
860*9566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecsDeltaLam2));
861*9566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&rk->VecsSensi2Temp));
86213af1a74SHong Zhang   }
86313af1a74SHong Zhang   if (ts->vecs_sensi2p) {
864*9566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vecs_sensi2p[0],&rk->VecDeltaMu2));
86513af1a74SHong Zhang   }
866f6a906c0SBarry Smith   PetscFunctionReturn(0);
867f6a906c0SBarry Smith }
868f6a906c0SBarry Smith 
86935f5172aSHong Zhang /*
87035f5172aSHong Zhang   Assumptions:
87135f5172aSHong Zhang     - TSStep_RK() always evaluates the step with b, not bembed.
87235f5172aSHong Zhang */
87342f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts)
874d2daff3dSHong Zhang {
875c235aa8dSHong Zhang   TS_RK            *rk = (TS_RK*)ts->data;
876cd4cee2dSHong Zhang   TS               quadts = ts->quadraturets;
877c235aa8dSHong Zhang   RKTableau        tab = rk->tableau;
8783ca0882eSHong Zhang   Mat              J,Jpre,Jquad;
879c235aa8dSHong Zhang   const PetscInt   s = tab->s;
880c235aa8dSHong Zhang   const PetscReal  *A = tab->A,*b = tab->b,*c = tab->c;
88113af1a74SHong Zhang   PetscScalar      *w = rk->work,*xarr;
8822e7b7f96SHong Zhang   Vec              *Y = rk->Y,*VecsDeltaLam = rk->VecsDeltaLam,VecDeltaMu = rk->VecDeltaMu,*VecsSensiTemp = rk->VecsSensiTemp;
88313af1a74SHong Zhang   Vec              *VecsDeltaLam2 = rk->VecsDeltaLam2,VecDeltaMu2 = rk->VecDeltaMu2,*VecsSensi2Temp = rk->VecsSensi2Temp;
884cd4cee2dSHong Zhang   Vec              VecDRDUTransCol = ts->vec_drdu_col,VecDRDPTransCol = ts->vec_drdp_col;
885b18ea86cSHong Zhang   PetscInt         i,j,nadj;
8863d3eaba7SBarry Smith   PetscReal        t = ts->ptime;
8873ca0882eSHong Zhang   PetscReal        h = ts->time_step;
888c235aa8dSHong Zhang 
889d2daff3dSHong Zhang   PetscFunctionBegin;
890c235aa8dSHong Zhang   rk->status = TS_STEP_INCOMPLETE;
8913389c7d9SStefano Zampini 
892*9566063dSJacob Faibussowitsch   PetscCall(TSGetRHSJacobian(ts,&J,&Jpre,NULL,NULL));
89335f5172aSHong Zhang   if (quadts) {
894*9566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobian(quadts,&Jquad,NULL,NULL,NULL));
89535f5172aSHong Zhang   }
8962e7b7f96SHong Zhang   for (i=s-1; i>=0; i--) {
8976a1e4597SHong Zhang     if (tab->FSAL && i == s-1) {
8986a1e4597SHong Zhang       /* VecsDeltaLam[nadj*s+s-1] are initialized with zeros and the values never change.*/
8996a1e4597SHong Zhang       continue;
9006a1e4597SHong Zhang     }
9013ca0882eSHong Zhang     rk->stage_time = t + h*(1.0-c[i]);
902*9566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts,Y[i],J,Jpre));
90335f5172aSHong Zhang     if (quadts) {
904*9566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts,rk->stage_time,Y[i],Jquad,Jquad)); /* get r_u^T */
90535f5172aSHong Zhang     }
9063389c7d9SStefano Zampini     if (ts->vecs_sensip) {
907*9566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(ts,rk->stage_time,Y[i],ts->Jacprhs)); /* get f_p */
90835f5172aSHong Zhang       if (quadts) {
909*9566063dSJacob Faibussowitsch         PetscCall(TSComputeRHSJacobianP(quadts,rk->stage_time,Y[i],quadts->Jacprhs)); /* get f_p for the quadrature */
9103389c7d9SStefano Zampini       }
91135f5172aSHong Zhang     }
9123389c7d9SStefano Zampini 
91335f5172aSHong Zhang     if (b[i]) {
91435f5172aSHong Zhang       for (j=i+1; j<s; j++) w[j-i-1] = A[j*s+i]/b[i]; /* coefficients for computing VecsSensiTemp */
91535f5172aSHong Zhang     } else {
91635f5172aSHong Zhang       for (j=i+1; j<s; j++) w[j-i-1] = A[j*s+i]; /* coefficients for computing VecsSensiTemp */
91735f5172aSHong Zhang     }
9182e7b7f96SHong Zhang 
9192e7b7f96SHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
9203389c7d9SStefano Zampini       /* Stage values of lambda */
92135f5172aSHong Zhang       if (b[i]) {
92235f5172aSHong Zhang         /* lambda_{n+1} + \sum_{j=i+1}^s a_{ji}/b[i]*lambda_{s,j} */
923*9566063dSJacob Faibussowitsch         PetscCall(VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj])); /* VecDeltaLam is an vec array of size s by numcost */
924*9566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(VecsSensiTemp[nadj],s-i-1,w,&VecsDeltaLam[nadj*s+i+1]));
925*9566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(J,VecsSensiTemp[nadj],VecsDeltaLam[nadj*s+i])); /* VecsSensiTemp will be reused by 2nd-order adjoint */
926*9566063dSJacob Faibussowitsch         PetscCall(VecScale(VecsDeltaLam[nadj*s+i],-h*b[i]));
92735f5172aSHong Zhang         if (quadts) {
928*9566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(Jquad,nadj,&xarr));
929*9566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(VecDRDUTransCol,xarr));
930*9566063dSJacob Faibussowitsch           PetscCall(VecAXPY(VecsDeltaLam[nadj*s+i],-h*b[i],VecDRDUTransCol));
931*9566063dSJacob Faibussowitsch           PetscCall(VecResetArray(VecDRDUTransCol));
932*9566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(Jquad,&xarr));
933cd4cee2dSHong Zhang         }
9343389c7d9SStefano Zampini       } else {
9356a1e4597SHong Zhang         /* \sum_{j=i+1}^s a_{ji}*lambda_{s,j} */
936*9566063dSJacob Faibussowitsch         PetscCall(VecSet(VecsSensiTemp[nadj],0));
937*9566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(VecsSensiTemp[nadj],s-i-1,w,&VecsDeltaLam[nadj*s+i+1]));
938*9566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(J,VecsSensiTemp[nadj],VecsDeltaLam[nadj*s+i]));
939*9566063dSJacob Faibussowitsch         PetscCall(VecScale(VecsDeltaLam[nadj*s+i],-h));
9403389c7d9SStefano Zampini       }
9416497ce14SHong Zhang 
942ad8e2604SHong Zhang       /* Stage values of mu */
9436497ce14SHong Zhang       if (ts->vecs_sensip) {
94435f5172aSHong Zhang         if (b[i]) {
945*9566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(ts->Jacprhs,VecsSensiTemp[nadj],VecDeltaMu));
946*9566063dSJacob Faibussowitsch           PetscCall(VecScale(VecDeltaMu,-h*b[i]));
94735f5172aSHong Zhang           if (quadts) {
948*9566063dSJacob Faibussowitsch             PetscCall(MatDenseGetColumn(quadts->Jacprhs,nadj,&xarr));
949*9566063dSJacob Faibussowitsch             PetscCall(VecPlaceArray(VecDRDPTransCol,xarr));
950*9566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecDeltaMu,-h*b[i],VecDRDPTransCol));
951*9566063dSJacob Faibussowitsch             PetscCall(VecResetArray(VecDRDPTransCol));
952*9566063dSJacob Faibussowitsch             PetscCall(MatDenseRestoreColumn(quadts->Jacprhs,&xarr));
953493ed95fSHong Zhang           }
95435f5172aSHong Zhang         } else {
955*9566063dSJacob Faibussowitsch           PetscCall(VecScale(VecDeltaMu,-h));
95635f5172aSHong Zhang         }
957*9566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj],1.,VecDeltaMu)); /* update sensip for each stage */
958ad8e2604SHong Zhang       }
959c235aa8dSHong Zhang     }
96013af1a74SHong Zhang 
96113af1a74SHong Zhang     if (ts->vecs_sensi2 && ts->forward_solve) { /* 2nd-order adjoint, TLM mode has to be turned on */
96213af1a74SHong Zhang       /* Get w1 at t_{n+1} from TLM matrix */
963*9566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(rk->MatsFwdStageSensip[i],0,&xarr));
964*9566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_sensip_col,xarr));
96513af1a74SHong Zhang       /* lambda_s^T F_UU w_1 */
966*9566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSHessianProductFunctionUU(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_sensip_col,ts->vecs_guu));
96735f5172aSHong Zhang       if (quadts)  {
96835f5172aSHong Zhang         /* R_UU w_1 */
969*9566063dSJacob Faibussowitsch         PetscCall(TSComputeRHSHessianProductFunctionUU(quadts,rk->stage_time,Y[i],NULL,ts->vec_sensip_col,ts->vecs_guu));
97035f5172aSHong Zhang       }
97135f5172aSHong Zhang       if (ts->vecs_sensip) {
97213af1a74SHong Zhang         /* lambda_s^T F_UP w_2 */
973*9566063dSJacob Faibussowitsch         PetscCall(TSComputeRHSHessianProductFunctionUP(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_dir,ts->vecs_gup));
97435f5172aSHong Zhang         if (quadts)  {
97535f5172aSHong Zhang           /* R_UP w_2 */
976*9566063dSJacob Faibussowitsch           PetscCall(TSComputeRHSHessianProductFunctionUP(quadts,rk->stage_time,Y[i],NULL,ts->vec_sensip_col,ts->vecs_gup));
97735f5172aSHong Zhang         }
97835f5172aSHong Zhang       }
97913af1a74SHong Zhang       if (ts->vecs_sensi2p) {
98013af1a74SHong Zhang         /* lambda_s^T F_PU w_1 */
981*9566063dSJacob Faibussowitsch         PetscCall(TSComputeRHSHessianProductFunctionPU(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_sensip_col,ts->vecs_gpu));
98235f5172aSHong Zhang         /* lambda_s^T F_PP w_2 */
983*9566063dSJacob Faibussowitsch         PetscCall(TSComputeRHSHessianProductFunctionPP(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_dir,ts->vecs_gpp));
98435f5172aSHong Zhang         if (b[i] && quadts) {
98535f5172aSHong Zhang           /* R_PU w_1 */
986*9566063dSJacob Faibussowitsch           PetscCall(TSComputeRHSHessianProductFunctionPU(quadts,rk->stage_time,Y[i],NULL,ts->vec_sensip_col,ts->vecs_gpu));
98735f5172aSHong Zhang           /* R_PP w_2 */
988*9566063dSJacob Faibussowitsch           PetscCall(TSComputeRHSHessianProductFunctionPP(quadts,rk->stage_time,Y[i],NULL,ts->vec_dir,ts->vecs_gpp));
98935f5172aSHong Zhang         }
99013af1a74SHong Zhang       }
991*9566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_sensip_col));
992*9566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(rk->MatsFwdStageSensip[i],&xarr));
99313af1a74SHong Zhang 
99413af1a74SHong Zhang       for (nadj=0; nadj<ts->numcost; nadj++) {
99513af1a74SHong Zhang         /* Stage values of lambda */
99635f5172aSHong Zhang         if (b[i]) {
99735f5172aSHong Zhang           /* J_i^T*(Lambda_{n+1}+\sum_{j=i+1}^s a_{ji}/b_i*Lambda_{s,j} */
998*9566063dSJacob Faibussowitsch           PetscCall(VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]));
999*9566063dSJacob Faibussowitsch           PetscCall(VecMAXPY(VecsSensi2Temp[nadj],s-i-1,w,&VecsDeltaLam2[nadj*s+i+1]));
1000*9566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(J,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj*s+i]));
1001*9566063dSJacob Faibussowitsch           PetscCall(VecScale(VecsDeltaLam2[nadj*s+i],-h*b[i]));
1002*9566063dSJacob Faibussowitsch           PetscCall(VecAXPY(VecsDeltaLam2[nadj*s+i],-h*b[i],ts->vecs_guu[nadj]));
100335f5172aSHong Zhang           if (ts->vecs_sensip) {
1004*9566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecsDeltaLam2[nadj*s+i],-h*b[i],ts->vecs_gup[nadj]));
100513af1a74SHong Zhang           }
100613af1a74SHong Zhang         } else {
100735f5172aSHong Zhang           /* \sum_{j=i+1}^s a_{ji}*Lambda_{s,j} */
1008*9566063dSJacob Faibussowitsch           PetscCall(VecSet(VecsDeltaLam2[nadj*s+i],0));
1009*9566063dSJacob Faibussowitsch           PetscCall(VecMAXPY(VecsSensi2Temp[nadj],s-i-1,w,&VecsDeltaLam2[nadj*s+i+1]));
1010*9566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(J,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj*s+i]));
1011*9566063dSJacob Faibussowitsch           PetscCall(VecScale(VecsDeltaLam2[nadj*s+i],-h));
1012*9566063dSJacob Faibussowitsch           PetscCall(VecAXPY(VecsDeltaLam2[nadj*s+i],-h,ts->vecs_guu[nadj]));
101335f5172aSHong Zhang           if (ts->vecs_sensip) {
1014*9566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecsDeltaLam2[nadj*s+i],-h,ts->vecs_gup[nadj]));
101513af1a74SHong Zhang           }
101635f5172aSHong Zhang         }
101735f5172aSHong Zhang         if (ts->vecs_sensi2p) { /* 2nd-order adjoint for parameters */
1018*9566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(ts->Jacprhs,VecsSensi2Temp[nadj],VecDeltaMu2));
101935f5172aSHong Zhang           if (b[i]) {
1020*9566063dSJacob Faibussowitsch             PetscCall(VecScale(VecDeltaMu2,-h*b[i]));
1021*9566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecDeltaMu2,-h*b[i],ts->vecs_gpu[nadj]));
1022*9566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecDeltaMu2,-h*b[i],ts->vecs_gpp[nadj]));
102313af1a74SHong Zhang           } else {
1024*9566063dSJacob Faibussowitsch             PetscCall(VecScale(VecDeltaMu2,-h));
1025*9566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecDeltaMu2,-h,ts->vecs_gpu[nadj]));
1026*9566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecDeltaMu2,-h,ts->vecs_gpp[nadj]));
102713af1a74SHong Zhang           }
1028*9566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensi2p[nadj],1,VecDeltaMu2)); /* update sensi2p for each stage */
102913af1a74SHong Zhang         }
103013af1a74SHong Zhang       }
103113af1a74SHong Zhang     }
10326497ce14SHong Zhang   }
1033c235aa8dSHong Zhang 
1034c235aa8dSHong Zhang   for (j=0; j<s; j++) w[j] = 1.0;
10352e7b7f96SHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) { /* no need to do this for mu's */
1036*9566063dSJacob Faibussowitsch     PetscCall(VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecsDeltaLam[nadj*s]));
103713af1a74SHong Zhang     if (ts->vecs_sensi2) {
1038*9566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(ts->vecs_sensi2[nadj],s,w,&VecsDeltaLam2[nadj*s]));
103913af1a74SHong Zhang     }
10406497ce14SHong Zhang   }
1041c235aa8dSHong Zhang   rk->status = TS_STEP_COMPLETE;
1042d2daff3dSHong Zhang   PetscFunctionReturn(0);
1043d2daff3dSHong Zhang }
1044d2daff3dSHong Zhang 
104513af1a74SHong Zhang static PetscErrorCode TSAdjointReset_RK(TS ts)
104613af1a74SHong Zhang {
104713af1a74SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
104813af1a74SHong Zhang   RKTableau      tab = rk->tableau;
104913af1a74SHong Zhang 
105013af1a74SHong Zhang   PetscFunctionBegin;
1051*9566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam));
1052*9566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost,&rk->VecsSensiTemp));
1053*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&rk->VecDeltaMu));
1054*9566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam2));
1055*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&rk->VecDeltaMu2));
1056*9566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost,&rk->VecsSensi2Temp));
105713af1a74SHong Zhang   PetscFunctionReturn(0);
105813af1a74SHong Zhang }
105913af1a74SHong Zhang 
106055de54fcSHong Zhang static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X)
106155de54fcSHong Zhang {
106255de54fcSHong Zhang   TS_RK            *rk = (TS_RK*)ts->data;
106355de54fcSHong Zhang   PetscInt         s  = rk->tableau->s,p = rk->tableau->p,i,j;
106455de54fcSHong Zhang   PetscReal        h;
106555de54fcSHong Zhang   PetscReal        tt,t;
106655de54fcSHong Zhang   PetscScalar      *b;
106755de54fcSHong Zhang   const PetscReal  *B = rk->tableau->binterp;
106855de54fcSHong Zhang 
106955de54fcSHong Zhang   PetscFunctionBegin;
10703c633725SBarry Smith   PetscCheck(B,PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
107155de54fcSHong Zhang 
107255de54fcSHong Zhang   switch (rk->status) {
107355de54fcSHong Zhang     case TS_STEP_INCOMPLETE:
107455de54fcSHong Zhang     case TS_STEP_PENDING:
107555de54fcSHong Zhang       h = ts->time_step;
107655de54fcSHong Zhang       t = (itime - ts->ptime)/h;
107755de54fcSHong Zhang       break;
107855de54fcSHong Zhang     case TS_STEP_COMPLETE:
107955de54fcSHong Zhang       h = ts->ptime - ts->ptime_prev;
108055de54fcSHong Zhang       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
108155de54fcSHong Zhang       break;
108255de54fcSHong Zhang     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
108355de54fcSHong Zhang   }
1084*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s,&b));
108555de54fcSHong Zhang   for (i=0; i<s; i++) b[i] = 0;
108655de54fcSHong Zhang   for (j=0,tt=t; j<p; j++,tt*=t) {
108755de54fcSHong Zhang     for (i=0; i<s; i++) {
108855de54fcSHong Zhang       b[i]  += h * B[i*p+j] * tt;
108955de54fcSHong Zhang     }
109055de54fcSHong Zhang   }
1091*9566063dSJacob Faibussowitsch   PetscCall(VecCopy(rk->Y[0],X));
1092*9566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(X,s,b,rk->YdotRHS));
1093*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(b));
109455de54fcSHong Zhang   PetscFunctionReturn(0);
109555de54fcSHong Zhang }
109655de54fcSHong Zhang 
109755de54fcSHong Zhang /*------------------------------------------------------------*/
109855de54fcSHong Zhang 
1099be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts)
1100be5899b3SLisandro Dalcin {
1101be5899b3SLisandro Dalcin   TS_RK          *rk = (TS_RK*)ts->data;
1102be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
1103be5899b3SLisandro Dalcin 
1104be5899b3SLisandro Dalcin   PetscFunctionBegin;
1105be5899b3SLisandro Dalcin   if (!tab) PetscFunctionReturn(0);
1106*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(rk->work));
1107*9566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(tab->s,&rk->Y));
1108*9566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(tab->s,&rk->YdotRHS));
1109be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
1110be5899b3SLisandro Dalcin }
1111be5899b3SLisandro Dalcin 
1112277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts)
1113e4dd521cSBarry Smith {
1114e4dd521cSBarry Smith   PetscFunctionBegin;
1115*9566063dSJacob Faibussowitsch   PetscCall(TSRKTableauReset(ts));
11160fe4d17eSHong Zhang   if (ts->use_splitrhsfunction) {
1117*9566063dSJacob Faibussowitsch     PetscCall(PetscTryMethod(ts,"TSReset_RK_MultirateSplit_C",(TS),(ts)));
11180fe4d17eSHong Zhang   } else {
1119*9566063dSJacob Faibussowitsch     PetscCall(PetscTryMethod(ts,"TSReset_RK_MultirateNonsplit_C",(TS),(ts)));
11200fe4d17eSHong Zhang   }
1121277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
1122e4dd521cSBarry Smith }
1123277b19d0SLisandro Dalcin 
1124f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx)
1125f68a32c8SEmil Constantinescu {
1126f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1127f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1128f68a32c8SEmil Constantinescu }
1129f68a32c8SEmil Constantinescu 
1130f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
1131f68a32c8SEmil Constantinescu {
1132f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1133f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1134f68a32c8SEmil Constantinescu }
1135f68a32c8SEmil Constantinescu 
1136f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx)
1137f68a32c8SEmil Constantinescu {
1138f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1139f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1140f68a32c8SEmil Constantinescu }
1141f68a32c8SEmil Constantinescu 
1142f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
1143f68a32c8SEmil Constantinescu {
1144f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1145f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1146f68a32c8SEmil Constantinescu }
1147be5899b3SLisandro Dalcin 
1148be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts)
1149be5899b3SLisandro Dalcin {
1150be5899b3SLisandro Dalcin   TS_RK          *rk  = (TS_RK*)ts->data;
1151be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
1152be5899b3SLisandro Dalcin 
1153be5899b3SLisandro Dalcin   PetscFunctionBegin;
1154*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tab->s,&rk->work));
1155*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y));
1156*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS));
1157be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
1158be5899b3SLisandro Dalcin }
1159be5899b3SLisandro Dalcin 
1160f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts)
1161f68a32c8SEmil Constantinescu {
1162cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
1163f68a32c8SEmil Constantinescu   DM             dm;
1164f68a32c8SEmil Constantinescu 
1165f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1166*9566063dSJacob Faibussowitsch   PetscCall(TSCheckImplicitTerm(ts));
1167*9566063dSJacob Faibussowitsch   PetscCall(TSRKTableauSetUp(ts));
1168cd4cee2dSHong Zhang   if (quadts && ts->costintegralfwd) {
1169cd4cee2dSHong Zhang     Mat Jquad;
1170*9566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobian(quadts,&Jquad,NULL,NULL,NULL));
11710f7a1166SHong Zhang   }
1172*9566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&dm));
1173*9566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts));
1174*9566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts));
11750fe4d17eSHong Zhang   if (ts->use_splitrhsfunction) {
1176*9566063dSJacob Faibussowitsch     PetscCall(PetscTryMethod(ts,"TSSetUp_RK_MultirateSplit_C",(TS),(ts)));
11770fe4d17eSHong Zhang   } else {
1178*9566063dSJacob Faibussowitsch     PetscCall(PetscTryMethod(ts,"TSSetUp_RK_MultirateNonsplit_C",(TS),(ts)));
11790fe4d17eSHong Zhang   }
1180e4dd521cSBarry Smith   PetscFunctionReturn(0);
1181e4dd521cSBarry Smith }
1182d2daff3dSHong Zhang 
11834416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts)
1184e4dd521cSBarry Smith {
1185be5899b3SLisandro Dalcin   TS_RK          *rk = (TS_RK*)ts->data;
1186dfbe8321SBarry Smith   PetscErrorCode ierr;
1187262ff7bbSBarry Smith 
1188e4dd521cSBarry Smith   PetscFunctionBegin;
1189*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHead(PetscOptionsObject,"RK ODE solver options"));
1190f68a32c8SEmil Constantinescu   {
1191f68a32c8SEmil Constantinescu     RKTableauLink link;
1192f68a32c8SEmil Constantinescu     PetscInt      count,choice;
11930fe4d17eSHong Zhang     PetscBool     flg,use_multirate = PETSC_FALSE;
1194f68a32c8SEmil Constantinescu     const char    **namelist;
11952c9cb062Sluzhanghpp 
1196f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) ;
1197*9566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(count,(char***)&namelist));
1198f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
1199*9566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ts_rk_multirate","Use interpolation-based multirate RK method","TSRKSetMultirate",rk->use_multirate,&use_multirate,&flg));
12000fe4d17eSHong Zhang     if (flg) {
1201*9566063dSJacob Faibussowitsch       PetscCall(TSRKSetMultirate(ts,use_multirate));
12020fe4d17eSHong Zhang     }
1203*9566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg));
1204*9566063dSJacob Faibussowitsch     if (flg) PetscCall(TSRKSetType(ts,namelist[choice]));
1205*9566063dSJacob Faibussowitsch     PetscCall(PetscFree(namelist));
1206f68a32c8SEmil Constantinescu   }
1207*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsTail());
1208*9566063dSJacob Faibussowitsch   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)ts),NULL,"Multirate methods options","");PetscCall(ierr);
1209*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-ts_rk_dtratio","time step ratio between slow and fast","",rk->dtratio,&rk->dtratio,NULL));
1210*9566063dSJacob Faibussowitsch   ierr = PetscOptionsEnd();PetscCall(ierr);
1211e4dd521cSBarry Smith   PetscFunctionReturn(0);
1212e4dd521cSBarry Smith }
1213e4dd521cSBarry Smith 
12145f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
1215e4dd521cSBarry Smith {
12165f70b478SJed Brown   TS_RK          *rk = (TS_RK*)ts->data;
12178370ee3bSLisandro Dalcin   PetscBool      iascii;
12182cdf8120Sasbjorn 
1219e4dd521cSBarry Smith   PetscFunctionBegin;
1220*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
12218370ee3bSLisandro Dalcin   if (iascii) {
12229c334d8fSLisandro Dalcin     RKTableau       tab  = rk->tableau;
1223f68a32c8SEmil Constantinescu     TSRKType        rktype;
12240f7a28cdSStefano Zampini     const PetscReal *c;
12250f7a28cdSStefano Zampini     PetscInt        s;
1226f68a32c8SEmil Constantinescu     char            buf[512];
12270f7a28cdSStefano Zampini     PetscBool       FSAL;
12280f7a28cdSStefano Zampini 
1229*9566063dSJacob Faibussowitsch     PetscCall(TSRKGetType(ts,&rktype));
1230*9566063dSJacob Faibussowitsch     PetscCall(TSRKGetTableau(ts,&s,NULL,NULL,&c,NULL,NULL,NULL,&FSAL));
1231*9566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  RK type %s\n",rktype));
1232*9566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  Order: %D\n",tab->order));
1233*9566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  FSAL property: %s\n",FSAL ? "yes" : "no"));
1234*9566063dSJacob Faibussowitsch     PetscCall(PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",s,c));
1235*9566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  Abscissa c = %s\n",buf));
12368370ee3bSLisandro Dalcin   }
1237f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1238f68a32c8SEmil Constantinescu }
1239f68a32c8SEmil Constantinescu 
1240f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer)
1241f68a32c8SEmil Constantinescu {
12429c334d8fSLisandro Dalcin   TSAdapt        adapt;
1243f68a32c8SEmil Constantinescu 
1244f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1245*9566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts,&adapt));
1246*9566063dSJacob Faibussowitsch   PetscCall(TSAdaptLoad(adapt,viewer));
1247f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1248f68a32c8SEmil Constantinescu }
1249f68a32c8SEmil Constantinescu 
12502ea87ba9SLisandro Dalcin /*@
12512ea87ba9SLisandro Dalcin   TSRKGetOrder - Get the order of RK scheme
12522ea87ba9SLisandro Dalcin 
12532ea87ba9SLisandro Dalcin   Not collective
12542ea87ba9SLisandro Dalcin 
12552ea87ba9SLisandro Dalcin   Input Parameter:
12562ea87ba9SLisandro Dalcin .  ts - timestepping context
12572ea87ba9SLisandro Dalcin 
12582ea87ba9SLisandro Dalcin   Output Parameter:
12592ea87ba9SLisandro Dalcin .  order - order of RK-scheme
12602ea87ba9SLisandro Dalcin 
12612ea87ba9SLisandro Dalcin   Level: intermediate
12622ea87ba9SLisandro Dalcin 
12632ea87ba9SLisandro Dalcin .seealso: TSRKGetType()
12642ea87ba9SLisandro Dalcin @*/
12652ea87ba9SLisandro Dalcin PetscErrorCode TSRKGetOrder(TS ts,PetscInt *order)
12662ea87ba9SLisandro Dalcin {
12672ea87ba9SLisandro Dalcin   PetscFunctionBegin;
12682ea87ba9SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
12692ea87ba9SLisandro Dalcin   PetscValidIntPointer(order,2);
1270*9566063dSJacob Faibussowitsch   PetscCall(PetscUseMethod(ts,"TSRKGetOrder_C",(TS,PetscInt*),(ts,order)));
12712ea87ba9SLisandro Dalcin   PetscFunctionReturn(0);
12722ea87ba9SLisandro Dalcin }
12732ea87ba9SLisandro Dalcin 
1274f68a32c8SEmil Constantinescu /*@C
1275f68a32c8SEmil Constantinescu   TSRKSetType - Set the type of RK scheme
1276f68a32c8SEmil Constantinescu 
1277f68a32c8SEmil Constantinescu   Logically collective
1278f68a32c8SEmil Constantinescu 
1279d8d19677SJose E. Roman   Input Parameters:
1280f68a32c8SEmil Constantinescu +  ts - timestepping context
1281f68a32c8SEmil Constantinescu -  rktype - type of RK-scheme
1282f68a32c8SEmil Constantinescu 
1283287c2655SBarry Smith   Options Database:
12849bd3e852SBarry Smith .   -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs>
1285287c2655SBarry Smith 
1286f68a32c8SEmil Constantinescu   Level: intermediate
1287f68a32c8SEmil Constantinescu 
1288e7685601SHong Zhang .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK2B, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS, TSRK6VR, TSRK7VR, TSRK8VR
1289f68a32c8SEmil Constantinescu @*/
1290f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype)
1291f68a32c8SEmil Constantinescu {
1292f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1293f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1294b92453a8SLisandro Dalcin   PetscValidCharPointer(rktype,2);
1295*9566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype)));
1296f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1297f68a32c8SEmil Constantinescu }
1298f68a32c8SEmil Constantinescu 
1299f68a32c8SEmil Constantinescu /*@C
1300f68a32c8SEmil Constantinescu   TSRKGetType - Get the type of RK scheme
1301f68a32c8SEmil Constantinescu 
13022ea87ba9SLisandro Dalcin   Not collective
1303f68a32c8SEmil Constantinescu 
1304f68a32c8SEmil Constantinescu   Input Parameter:
1305f68a32c8SEmil Constantinescu .  ts - timestepping context
1306f68a32c8SEmil Constantinescu 
1307f68a32c8SEmil Constantinescu   Output Parameter:
1308f68a32c8SEmil Constantinescu .  rktype - type of RK-scheme
1309f68a32c8SEmil Constantinescu 
1310f68a32c8SEmil Constantinescu   Level: intermediate
1311f68a32c8SEmil Constantinescu 
13122ea87ba9SLisandro Dalcin .seealso: TSRKSetType()
1313f68a32c8SEmil Constantinescu @*/
1314f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype)
1315f68a32c8SEmil Constantinescu {
1316f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1317f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1318*9566063dSJacob Faibussowitsch   PetscCall(PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype)));
1319f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1320f68a32c8SEmil Constantinescu }
1321f68a32c8SEmil Constantinescu 
13222ea87ba9SLisandro Dalcin static PetscErrorCode TSRKGetOrder_RK(TS ts,PetscInt *order)
13232ea87ba9SLisandro Dalcin {
13242ea87ba9SLisandro Dalcin   TS_RK *rk = (TS_RK*)ts->data;
13252ea87ba9SLisandro Dalcin 
13262ea87ba9SLisandro Dalcin   PetscFunctionBegin;
13272ea87ba9SLisandro Dalcin   *order = rk->tableau->order;
13282ea87ba9SLisandro Dalcin   PetscFunctionReturn(0);
13292ea87ba9SLisandro Dalcin }
13302ea87ba9SLisandro Dalcin 
1331560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype)
1332f68a32c8SEmil Constantinescu {
1333f68a32c8SEmil Constantinescu   TS_RK *rk = (TS_RK*)ts->data;
1334f68a32c8SEmil Constantinescu 
1335f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1336f68a32c8SEmil Constantinescu   *rktype = rk->tableau->name;
1337f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1338f68a32c8SEmil Constantinescu }
13392c9cb062Sluzhanghpp 
1340560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype)
1341f68a32c8SEmil Constantinescu {
1342f68a32c8SEmil Constantinescu   TS_RK          *rk = (TS_RK*)ts->data;
1343f68a32c8SEmil Constantinescu   PetscBool      match;
1344f68a32c8SEmil Constantinescu   RKTableauLink  link;
1345f68a32c8SEmil Constantinescu 
1346f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1347f68a32c8SEmil Constantinescu   if (rk->tableau) {
1348*9566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(rk->tableau->name,rktype,&match));
1349f68a32c8SEmil Constantinescu     if (match) PetscFunctionReturn(0);
1350f68a32c8SEmil Constantinescu   }
1351f68a32c8SEmil Constantinescu   for (link = RKTableauList; link; link=link->next) {
1352*9566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(link->tab.name,rktype,&match));
1353f68a32c8SEmil Constantinescu     if (match) {
1354*9566063dSJacob Faibussowitsch       if (ts->setupcalled) PetscCall(TSRKTableauReset(ts));
1355f68a32c8SEmil Constantinescu       rk->tableau = &link->tab;
1356*9566063dSJacob Faibussowitsch       if (ts->setupcalled) PetscCall(TSRKTableauSetUp(ts));
13572ffb9264SLisandro Dalcin       ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
1358f68a32c8SEmil Constantinescu       PetscFunctionReturn(0);
1359f68a32c8SEmil Constantinescu     }
1360f68a32c8SEmil Constantinescu   }
136198921bdaSJacob Faibussowitsch   SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype);
1362e4dd521cSBarry Smith }
1363e4dd521cSBarry Smith 
1364ff22ae23SHong Zhang static PetscErrorCode  TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y)
1365ff22ae23SHong Zhang {
1366ff22ae23SHong Zhang   TS_RK *rk = (TS_RK*)ts->data;
1367ff22ae23SHong Zhang 
1368ff22ae23SHong Zhang   PetscFunctionBegin;
13690429704eSStefano Zampini   if (ns) *ns = rk->tableau->s;
1370d2daff3dSHong Zhang   if (Y)   *Y = rk->Y;
1371ff22ae23SHong Zhang   PetscFunctionReturn(0);
1372ff22ae23SHong Zhang }
1373ff22ae23SHong Zhang 
1374b3a6b972SJed Brown static PetscErrorCode TSDestroy_RK(TS ts)
1375b3a6b972SJed Brown {
1376b3a6b972SJed Brown   PetscFunctionBegin;
1377*9566063dSJacob Faibussowitsch   PetscCall(TSReset_RK(ts));
1378b3a6b972SJed Brown   if (ts->dm) {
1379*9566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts));
1380*9566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts));
1381b3a6b972SJed Brown   }
1382*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
1383*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKGetOrder_C",NULL));
1384*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL));
1385*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL));
1386*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKGetTableau_C",NULL));
1387*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",NULL));
1388*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",NULL));
1389b3a6b972SJed Brown   PetscFunctionReturn(0);
1390b3a6b972SJed Brown }
1391ff22ae23SHong Zhang 
13923ca0882eSHong Zhang /*
13933ca0882eSHong Zhang   This defines the nonlinear equation that is to be solved with SNES
13943ca0882eSHong Zhang   We do not need to solve the equation; we just use SNES to approximate the Jacobian
13953ca0882eSHong Zhang */
13963ca0882eSHong Zhang static PetscErrorCode SNESTSFormFunction_RK(SNES snes,Vec x,Vec y,TS ts)
13973ca0882eSHong Zhang {
13983ca0882eSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
13993ca0882eSHong Zhang   DM             dm,dmsave;
14003ca0882eSHong Zhang 
14013ca0882eSHong Zhang   PetscFunctionBegin;
1402*9566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
14033ca0882eSHong Zhang   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
14043ca0882eSHong Zhang   dmsave = ts->dm;
14053ca0882eSHong Zhang   ts->dm = dm;
1406*9566063dSJacob Faibussowitsch   PetscCall(TSComputeRHSFunction(ts,rk->stage_time,x,y));
14073ca0882eSHong Zhang   ts->dm = dmsave;
14083ca0882eSHong Zhang   PetscFunctionReturn(0);
14093ca0882eSHong Zhang }
14103ca0882eSHong Zhang 
14113ca0882eSHong Zhang static PetscErrorCode SNESTSFormJacobian_RK(SNES snes,Vec x,Mat A,Mat B,TS ts)
14123ca0882eSHong Zhang {
14133ca0882eSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
14143ca0882eSHong Zhang   DM             dm,dmsave;
14153ca0882eSHong Zhang 
14163ca0882eSHong Zhang   PetscFunctionBegin;
1417*9566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
14183ca0882eSHong Zhang   dmsave = ts->dm;
14193ca0882eSHong Zhang   ts->dm = dm;
1420*9566063dSJacob Faibussowitsch   PetscCall(TSComputeRHSJacobian(ts,rk->stage_time,x,A,B));
14213ca0882eSHong Zhang   ts->dm = dmsave;
14223ca0882eSHong Zhang   PetscFunctionReturn(0);
14233ca0882eSHong Zhang }
14243ca0882eSHong Zhang 
142521052c0cSHong Zhang /*@C
142621052c0cSHong Zhang   TSRKSetMultirate - Use the interpolation-based multirate RK method
142721052c0cSHong Zhang 
142821052c0cSHong Zhang   Logically collective
142921052c0cSHong Zhang 
1430d8d19677SJose E. Roman   Input Parameters:
143121052c0cSHong Zhang +  ts - timestepping context
143221052c0cSHong Zhang -  use_multirate - PETSC_TRUE enables the multirate RK method, sets the basic method to be RK2A and sets the ratio between slow stepsize and fast stepsize to be 2
143321052c0cSHong Zhang 
143421052c0cSHong Zhang   Options Database:
143521052c0cSHong Zhang .   -ts_rk_multirate - <true,false>
143621052c0cSHong Zhang 
143721052c0cSHong Zhang   Notes:
143821052c0cSHong Zhang   The multirate method requires interpolation. The default interpolation works for 1st- and 2nd- order RK, but not for high-order RKs except TSRK5DP which comes with the interpolation coeffcients (binterp).
143921052c0cSHong Zhang 
144021052c0cSHong Zhang   Level: intermediate
144121052c0cSHong Zhang 
144221052c0cSHong Zhang .seealso: TSRKGetMultirate()
144321052c0cSHong Zhang @*/
144421052c0cSHong Zhang PetscErrorCode TSRKSetMultirate(TS ts,PetscBool use_multirate)
144521052c0cSHong Zhang {
14468a4be4dbSHong Zhang   PetscFunctionBegin;
1447*9566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(ts,"TSRKSetMultirate_C",(TS,PetscBool),(ts,use_multirate)));
144821052c0cSHong Zhang   PetscFunctionReturn(0);
144921052c0cSHong Zhang }
145021052c0cSHong Zhang 
145121052c0cSHong Zhang /*@C
145221052c0cSHong Zhang   TSRKGetMultirate - Gets whether to Use the interpolation-based multirate RK method
145321052c0cSHong Zhang 
145421052c0cSHong Zhang   Not collective
145521052c0cSHong Zhang 
145621052c0cSHong Zhang   Input Parameter:
145721052c0cSHong Zhang .  ts - timestepping context
145821052c0cSHong Zhang 
145921052c0cSHong Zhang   Output Parameter:
146021052c0cSHong Zhang .  use_multirate - PETSC_TRUE if the multirate RK method is enabled, PETSC_FALSE otherwise
146121052c0cSHong Zhang 
146221052c0cSHong Zhang   Level: intermediate
146321052c0cSHong Zhang 
146421052c0cSHong Zhang .seealso: TSRKSetMultirate()
146521052c0cSHong Zhang @*/
146621052c0cSHong Zhang PetscErrorCode TSRKGetMultirate(TS ts,PetscBool *use_multirate)
146721052c0cSHong Zhang {
14687dbacdf2SHong Zhang   PetscFunctionBegin;
1469*9566063dSJacob Faibussowitsch   PetscCall(PetscUseMethod(ts,"TSRKGetMultirate_C",(TS,PetscBool*),(ts,use_multirate)));
147021052c0cSHong Zhang   PetscFunctionReturn(0);
147121052c0cSHong Zhang }
147221052c0cSHong Zhang 
1473ebe3b25bSBarry Smith /*MC
1474f68a32c8SEmil Constantinescu       TSRK - ODE and DAE solver using Runge-Kutta schemes
1475ebe3b25bSBarry Smith 
1476f68a32c8SEmil Constantinescu   The user should provide the right hand side of the equation
1477f68a32c8SEmil Constantinescu   using TSSetRHSFunction().
1478ebe3b25bSBarry Smith 
1479f68a32c8SEmil Constantinescu   Notes:
148098164e13SLisandro Dalcin   The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type
1481ebe3b25bSBarry Smith 
1482d41222bbSBarry Smith   Level: beginner
1483d41222bbSBarry Smith 
14845d1808f1SStefano Zampini .seealso:  TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRK2D, TTSRK2E, TSRK3,
14850fe4d17eSHong Zhang            TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister(), TSRKSetMultirate(), TSRKGetMultirate()
1486ebe3b25bSBarry Smith 
1487ebe3b25bSBarry Smith M*/
14888cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
1489e4dd521cSBarry Smith {
14901566a47fSLisandro Dalcin   TS_RK          *rk;
1491e4dd521cSBarry Smith 
1492e4dd521cSBarry Smith   PetscFunctionBegin;
1493*9566063dSJacob Faibussowitsch   PetscCall(TSRKInitializePackage());
1494f68a32c8SEmil Constantinescu 
1495f68a32c8SEmil Constantinescu   ts->ops->reset          = TSReset_RK;
14965f70b478SJed Brown   ts->ops->destroy        = TSDestroy_RK;
14975f70b478SJed Brown   ts->ops->view           = TSView_RK;
1498f68a32c8SEmil Constantinescu   ts->ops->load           = TSLoad_RK;
1499f68a32c8SEmil Constantinescu   ts->ops->setup          = TSSetUp_RK;
1500f68a32c8SEmil Constantinescu   ts->ops->interpolate    = TSInterpolate_RK;
15012c9cb062Sluzhanghpp   ts->ops->step           = TSStep_RK;
1502f68a32c8SEmil Constantinescu   ts->ops->evaluatestep   = TSEvaluateStep_RK;
1503fecfb714SLisandro Dalcin   ts->ops->rollback       = TSRollBack_RK;
1504f68a32c8SEmil Constantinescu   ts->ops->setfromoptions = TSSetFromOptions_RK;
1505ff22ae23SHong Zhang   ts->ops->getstages      = TSGetStages_RK;
1506e4dd521cSBarry Smith 
15073ca0882eSHong Zhang   ts->ops->snesfunction    = SNESTSFormFunction_RK;
15083ca0882eSHong Zhang   ts->ops->snesjacobian    = SNESTSFormJacobian_RK;
15090f7a1166SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_RK;
151013af1a74SHong Zhang   ts->ops->adjointsetup    = TSAdjointSetUp_RK;
151113af1a74SHong Zhang   ts->ops->adjointstep     = TSAdjointStep_RK;
151213af1a74SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_RK;
15130f7a1166SHong Zhang 
151413af1a74SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_RK;
1515922a638cSHong Zhang   ts->ops->forwardsetup    = TSForwardSetUp_RK;
1516922a638cSHong Zhang   ts->ops->forwardreset    = TSForwardReset_RK;
1517922a638cSHong Zhang   ts->ops->forwardstep     = TSForwardStep_RK;
1518922a638cSHong Zhang   ts->ops->forwardgetstages= TSForwardGetStages_RK;
1519922a638cSHong Zhang 
1520*9566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(ts,&rk));
15211566a47fSLisandro Dalcin   ts->data = (void*)rk;
1522e4dd521cSBarry Smith 
1523*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKGetOrder_C",TSRKGetOrder_RK));
1524*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK));
1525*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK));
1526*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKGetTableau_C",TSRKGetTableau_RK));
1527*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",TSRKSetMultirate_RK));
1528*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",TSRKGetMultirate_RK));
1529be5899b3SLisandro Dalcin 
1530*9566063dSJacob Faibussowitsch   PetscCall(TSRKSetType(ts,TSRKDefault));
153190b67129SHong Zhang   rk->dtratio = 1;
1532e4dd521cSBarry Smith   PetscFunctionReturn(0);
1533e4dd521cSBarry Smith }
1534