xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision cac4c232dc4f93991e342196e27ef7b0655dac7b)
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)};
2059566063dSJacob 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};
2139566063dSJacob 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)};
2209566063dSJacob 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)};
2289566063dSJacob 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)};
2389566063dSJacob 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)};
2479566063dSJacob 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};
2599566063dSJacob 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)}};
2799566063dSJacob 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)};
2939566063dSJacob 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)};
3089566063dSJacob 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)};
3249566063dSJacob 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)};
3439566063dSJacob 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;
3669566063dSJacob Faibussowitsch     PetscCall(PetscFree3(t->A,t->b,t->c));
3679566063dSJacob Faibussowitsch     PetscCall(PetscFree(t->bembed));
3689566063dSJacob Faibussowitsch     PetscCall(PetscFree(t->binterp));
3699566063dSJacob Faibussowitsch     PetscCall(PetscFree(t->name));
3709566063dSJacob 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;
3899566063dSJacob Faibussowitsch   PetscCall(TSRKRegisterAll());
3909566063dSJacob 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;
4069566063dSJacob 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 
4499566063dSJacob Faibussowitsch   PetscCall(TSRKInitializePackage());
4509566063dSJacob Faibussowitsch   PetscCall(PetscNew(&link));
451f68a32c8SEmil Constantinescu   t = &link->tab;
4523a8a9803SLisandro Dalcin 
4539566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name,&t->name));
454f68a32c8SEmil Constantinescu   t->order = order;
455f68a32c8SEmil Constantinescu   t->s = s;
4569566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c));
4579566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(t->A,A,s*s));
4589566063dSJacob 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];
4609566063dSJacob 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) {
4669566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(s,&t->bembed));
4679566063dSJacob 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;
4729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s*p,&t->binterp));
4739566063dSJacob 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   PetscFunctionBegin;
5240f7a28cdSStefano Zampini   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
525*cac4c232SBarry Smith   PetscUseMethod(ts,"TSRKGetTableau_C",(TS,PetscInt*,const PetscReal**,const PetscReal**,const PetscReal**,const PetscReal**,
526*cac4c232SBarry Smith                                         PetscInt*,const PetscReal**,PetscBool*),(ts,s,A,b,c,bembed,p,binterp,FSAL));
5270f7a28cdSStefano Zampini   PetscFunctionReturn(0);
5280f7a28cdSStefano Zampini }
5290f7a28cdSStefano Zampini 
530e4dd521cSBarry Smith /*
5312c9cb062Sluzhanghpp  This is for single-step RK method
532f68a32c8SEmil Constantinescu  The step completion formula is
533e4dd521cSBarry Smith 
534f68a32c8SEmil Constantinescu  x1 = x0 + h b^T YdotRHS
535f68a32c8SEmil Constantinescu 
536f68a32c8SEmil Constantinescu  This function can be called before or after ts->vec_sol has been updated.
537f68a32c8SEmil Constantinescu  Suppose we have a completion formula (b) and an embedded formula (be) of different order.
538f68a32c8SEmil Constantinescu  We can write
539f68a32c8SEmil Constantinescu 
540f68a32c8SEmil Constantinescu  x1e = x0 + h be^T YdotRHS
541f68a32c8SEmil Constantinescu      = x1 - h b^T YdotRHS + h be^T YdotRHS
542f68a32c8SEmil Constantinescu      = x1 + h (be - b)^T YdotRHS
543f68a32c8SEmil Constantinescu 
544f68a32c8SEmil Constantinescu  so we can evaluate the method with different order even after the step has been optimistically completed.
545e4dd521cSBarry Smith */
546f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done)
547f68a32c8SEmil Constantinescu {
548f68a32c8SEmil Constantinescu   TS_RK          *rk   = (TS_RK*)ts->data;
549f68a32c8SEmil Constantinescu   RKTableau      tab  = rk->tableau;
550f68a32c8SEmil Constantinescu   PetscScalar    *w    = rk->work;
551f68a32c8SEmil Constantinescu   PetscReal      h;
552f68a32c8SEmil Constantinescu   PetscInt       s    = tab->s,j;
553f68a32c8SEmil Constantinescu 
554f68a32c8SEmil Constantinescu   PetscFunctionBegin;
555f68a32c8SEmil Constantinescu   switch (rk->status) {
556f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
557f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
558f68a32c8SEmil Constantinescu     h = ts->time_step; break;
559f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
560be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
561f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
562f68a32c8SEmil Constantinescu   }
563f68a32c8SEmil Constantinescu   if (order == tab->order) {
564f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) {
5659566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol,X));
56690b67129SHong Zhang       for (j=0; j<s; j++) w[j] = h*tab->b[j]/rk->dtratio;
5679566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(X,s,w,rk->YdotRHS));
5689566063dSJacob Faibussowitsch     } else PetscCall(VecCopy(ts->vec_sol,X));
569f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
570f68a32c8SEmil Constantinescu   } else if (order == tab->order-1) {
571f68a32c8SEmil Constantinescu     if (!tab->bembed) goto unavailable;
572f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) { /*Complete with the embedded method (be)*/
5739566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol,X));
574f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
5759566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(X,s,w,rk->YdotRHS));
576f68a32c8SEmil Constantinescu     } else {  /*Rollback and re-complete using (be-b) */
5779566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol,X));
578f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
5799566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(X,s,w,rk->YdotRHS));
580f68a32c8SEmil Constantinescu     }
581f68a32c8SEmil Constantinescu     if (done) *done = PETSC_TRUE;
582f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
583f68a32c8SEmil Constantinescu   }
584f68a32c8SEmil Constantinescu unavailable:
585f68a32c8SEmil Constantinescu   if (done) *done = PETSC_FALSE;
58698921bdaSJacob 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);
587f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
588f68a32c8SEmil Constantinescu }
589f68a32c8SEmil Constantinescu 
5900f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts)
5910f7a1166SHong Zhang {
5920f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
593cd4cee2dSHong Zhang   TS              quadts = ts->quadraturets;
5940f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
5950f7a1166SHong Zhang   const PetscInt  s = tab->s;
5960f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
5970f7a1166SHong Zhang   Vec             *Y = rk->Y;
5980f7a1166SHong Zhang   PetscInt        i;
5990f7a1166SHong Zhang 
6000f7a1166SHong Zhang   PetscFunctionBegin;
601cd4cee2dSHong Zhang   /* No need to backup quadts->vec_sol since it can be reverted in TSRollBack_RK */
6020f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
603cd4cee2dSHong Zhang     /* Evolve quadrature TS solution to compute integrals */
6049566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(quadts,rk->ptime+rk->time_step*c[i],Y[i],ts->vec_costintegrand));
6059566063dSJacob Faibussowitsch     PetscCall(VecAXPY(quadts->vec_sol,rk->time_step*b[i],ts->vec_costintegrand));
6060f7a1166SHong Zhang   }
6070f7a1166SHong Zhang   PetscFunctionReturn(0);
6080f7a1166SHong Zhang }
6090f7a1166SHong Zhang 
6100f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts)
6110f7a1166SHong Zhang {
6120f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
6130f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
614cd4cee2dSHong Zhang   TS              quadts = ts->quadraturets;
6150f7a1166SHong Zhang   const PetscInt  s = tab->s;
6160f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
6170f7a1166SHong Zhang   Vec             *Y = rk->Y;
6180f7a1166SHong Zhang   PetscInt        i;
6190f7a1166SHong Zhang 
6200f7a1166SHong Zhang   PetscFunctionBegin;
6210f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
622cd4cee2dSHong Zhang     /* Evolve quadrature TS solution to compute integrals */
6239566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(quadts,ts->ptime+ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand));
6249566063dSJacob Faibussowitsch     PetscCall(VecAXPY(quadts->vec_sol,-ts->time_step*b[i],ts->vec_costintegrand));
6250f7a1166SHong Zhang   }
6260f7a1166SHong Zhang   PetscFunctionReturn(0);
6270f7a1166SHong Zhang }
6280f7a1166SHong Zhang 
629fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts)
630fecfb714SLisandro Dalcin {
631fecfb714SLisandro Dalcin   TS_RK           *rk = (TS_RK*)ts->data;
632cd4cee2dSHong Zhang   TS              quadts = ts->quadraturets;
633fecfb714SLisandro Dalcin   RKTableau       tab = rk->tableau;
634fecfb714SLisandro Dalcin   const PetscInt  s  = tab->s;
635cd4cee2dSHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
636fecfb714SLisandro Dalcin   PetscScalar     *w = rk->work;
637cd4cee2dSHong Zhang   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
638fecfb714SLisandro Dalcin   PetscInt        j;
639fecfb714SLisandro Dalcin   PetscReal       h;
640fecfb714SLisandro Dalcin 
641fecfb714SLisandro Dalcin   PetscFunctionBegin;
642fecfb714SLisandro Dalcin   switch (rk->status) {
643fecfb714SLisandro Dalcin   case TS_STEP_INCOMPLETE:
644fecfb714SLisandro Dalcin   case TS_STEP_PENDING:
645fecfb714SLisandro Dalcin     h = ts->time_step; break;
646fecfb714SLisandro Dalcin   case TS_STEP_COMPLETE:
647fecfb714SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
648fecfb714SLisandro Dalcin   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
649fecfb714SLisandro Dalcin   }
650fecfb714SLisandro Dalcin   for (j=0; j<s; j++) w[j] = -h*b[j];
6519566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(ts->vec_sol,s,w,YdotRHS));
652cd4cee2dSHong Zhang   if (quadts && ts->costintegralfwd) {
653cd4cee2dSHong Zhang     for (j=0; j<s; j++) {
654cd4cee2dSHong Zhang       /* Revert the quadrature TS solution */
6559566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSFunction(quadts,rk->ptime+h*c[j],Y[j],ts->vec_costintegrand));
6569566063dSJacob Faibussowitsch       PetscCall(VecAXPY(quadts->vec_sol,-h*b[j],ts->vec_costintegrand));
657cd4cee2dSHong Zhang     }
658cd4cee2dSHong Zhang   }
659fecfb714SLisandro Dalcin   PetscFunctionReturn(0);
660fecfb714SLisandro Dalcin }
661fecfb714SLisandro Dalcin 
662922a638cSHong Zhang static PetscErrorCode TSForwardStep_RK(TS ts)
663922a638cSHong Zhang {
664922a638cSHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
665922a638cSHong Zhang   RKTableau       tab = rk->tableau;
666922a638cSHong Zhang   Mat             J,*MatsFwdSensipTemp = rk->MatsFwdSensipTemp;
667922a638cSHong Zhang   const PetscInt  s = tab->s;
668922a638cSHong Zhang   const PetscReal *A = tab->A,*c = tab->c,*b = tab->b;
669922a638cSHong Zhang   Vec             *Y = rk->Y;
670922a638cSHong Zhang   PetscInt        i,j;
671922a638cSHong Zhang   PetscReal       stage_time,h = ts->time_step;
672922a638cSHong Zhang   PetscBool       zero;
673922a638cSHong Zhang 
674922a638cSHong Zhang   PetscFunctionBegin;
6759566063dSJacob Faibussowitsch   PetscCall(MatCopy(ts->mat_sensip,rk->MatFwdSensip0,SAME_NONZERO_PATTERN));
6769566063dSJacob Faibussowitsch   PetscCall(TSGetRHSJacobian(ts,&J,NULL,NULL,NULL));
677922a638cSHong Zhang 
678922a638cSHong Zhang   for (i=0; i<s; i++) {
679922a638cSHong Zhang     stage_time = ts->ptime + h*c[i];
680922a638cSHong Zhang     zero = PETSC_FALSE;
681922a638cSHong Zhang     if (b[i] == 0 && i == s-1) zero = PETSC_TRUE;
682922a638cSHong Zhang     /* TLM Stage values */
683922a638cSHong Zhang     if (!i) {
6849566063dSJacob Faibussowitsch       PetscCall(MatCopy(ts->mat_sensip,rk->MatsFwdStageSensip[i],SAME_NONZERO_PATTERN));
685922a638cSHong Zhang     } else if (!zero) {
6869566063dSJacob Faibussowitsch       PetscCall(MatZeroEntries(rk->MatsFwdStageSensip[i]));
687922a638cSHong Zhang       for (j=0; j<i; j++) {
6889566063dSJacob Faibussowitsch         PetscCall(MatAXPY(rk->MatsFwdStageSensip[i],h*A[i*s+j],MatsFwdSensipTemp[j],SAME_NONZERO_PATTERN));
689922a638cSHong Zhang       }
6909566063dSJacob Faibussowitsch       PetscCall(MatAXPY(rk->MatsFwdStageSensip[i],1.,ts->mat_sensip,SAME_NONZERO_PATTERN));
691922a638cSHong Zhang     } else {
6929566063dSJacob Faibussowitsch       PetscCall(MatZeroEntries(rk->MatsFwdStageSensip[i]));
693922a638cSHong Zhang     }
694922a638cSHong Zhang 
6959566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSJacobian(ts,stage_time,Y[i],J,J));
6969566063dSJacob Faibussowitsch     PetscCall(MatMatMult(J,rk->MatsFwdStageSensip[i],MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatsFwdSensipTemp[i]));
69713af1a74SHong Zhang     if (ts->Jacprhs) {
6989566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(ts,stage_time,Y[i],ts->Jacprhs)); /* get f_p */
69913af1a74SHong Zhang       if (ts->vecs_sensi2p) { /* TLM used for 2nd-order adjoint */
70013af1a74SHong Zhang         PetscScalar *xarr;
7019566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(MatsFwdSensipTemp[i],0,&xarr));
7029566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(rk->VecDeltaFwdSensipCol,xarr));
7039566063dSJacob Faibussowitsch         PetscCall(MatMultAdd(ts->Jacprhs,ts->vec_dir,rk->VecDeltaFwdSensipCol,rk->VecDeltaFwdSensipCol));
7049566063dSJacob Faibussowitsch         PetscCall(VecResetArray(rk->VecDeltaFwdSensipCol));
7059566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(MatsFwdSensipTemp[i],&xarr));
70613af1a74SHong Zhang       } else {
7079566063dSJacob Faibussowitsch         PetscCall(MatAXPY(MatsFwdSensipTemp[i],1.,ts->Jacprhs,SUBSET_NONZERO_PATTERN));
70813af1a74SHong Zhang       }
709922a638cSHong Zhang     }
710922a638cSHong Zhang   }
711922a638cSHong Zhang 
712922a638cSHong Zhang   for (i=0; i<s; i++) {
7139566063dSJacob Faibussowitsch     PetscCall(MatAXPY(ts->mat_sensip,h*b[i],rk->MatsFwdSensipTemp[i],SAME_NONZERO_PATTERN));
714922a638cSHong Zhang   }
715922a638cSHong Zhang   rk->status = TS_STEP_COMPLETE;
716922a638cSHong Zhang   PetscFunctionReturn(0);
717922a638cSHong Zhang }
718922a638cSHong Zhang 
719922a638cSHong Zhang static PetscErrorCode TSForwardGetStages_RK(TS ts,PetscInt *ns,Mat **stagesensip)
720922a638cSHong Zhang {
721922a638cSHong Zhang   TS_RK     *rk = (TS_RK*)ts->data;
722922a638cSHong Zhang   RKTableau tab  = rk->tableau;
723922a638cSHong Zhang 
724922a638cSHong Zhang   PetscFunctionBegin;
725922a638cSHong Zhang   if (ns) *ns = tab->s;
726922a638cSHong Zhang   if (stagesensip) *stagesensip = rk->MatsFwdStageSensip;
727922a638cSHong Zhang   PetscFunctionReturn(0);
728922a638cSHong Zhang }
729922a638cSHong Zhang 
730922a638cSHong Zhang static PetscErrorCode TSForwardSetUp_RK(TS ts)
731922a638cSHong Zhang {
732922a638cSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
733922a638cSHong Zhang   RKTableau      tab  = rk->tableau;
734922a638cSHong Zhang   PetscInt       i;
735922a638cSHong Zhang 
736922a638cSHong Zhang   PetscFunctionBegin;
737922a638cSHong Zhang   /* backup sensitivity results for roll-backs */
7389566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatFwdSensip0));
739922a638cSHong Zhang 
7409566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tab->s,&rk->MatsFwdStageSensip));
7419566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tab->s,&rk->MatsFwdSensipTemp));
742922a638cSHong Zhang   for (i=0; i<tab->s; i++) {
7439566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatsFwdStageSensip[i]));
7449566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatsFwdSensipTemp[i]));
745922a638cSHong Zhang   }
7469566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol,&rk->VecDeltaFwdSensipCol));
747922a638cSHong Zhang   PetscFunctionReturn(0);
748922a638cSHong Zhang }
749922a638cSHong Zhang 
750922a638cSHong Zhang static PetscErrorCode TSForwardReset_RK(TS ts)
751922a638cSHong Zhang {
752922a638cSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
753922a638cSHong Zhang   RKTableau      tab  = rk->tableau;
754922a638cSHong Zhang   PetscInt       i;
755922a638cSHong Zhang 
756922a638cSHong Zhang   PetscFunctionBegin;
7579566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&rk->MatFwdSensip0));
758922a638cSHong Zhang   if (rk->MatsFwdStageSensip) {
759922a638cSHong Zhang     for (i=0; i<tab->s; i++) {
7609566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&rk->MatsFwdStageSensip[i]));
761922a638cSHong Zhang     }
7629566063dSJacob Faibussowitsch     PetscCall(PetscFree(rk->MatsFwdStageSensip));
763922a638cSHong Zhang   }
764922a638cSHong Zhang   if (rk->MatsFwdSensipTemp) {
765922a638cSHong Zhang     for (i=0; i<tab->s; i++) {
7669566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&rk->MatsFwdSensipTemp[i]));
767922a638cSHong Zhang     }
7689566063dSJacob Faibussowitsch     PetscCall(PetscFree(rk->MatsFwdSensipTemp));
769922a638cSHong Zhang   }
7709566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&rk->VecDeltaFwdSensipCol));
771922a638cSHong Zhang   PetscFunctionReturn(0);
772922a638cSHong Zhang }
773922a638cSHong Zhang 
774f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts)
775f68a32c8SEmil Constantinescu {
776f68a32c8SEmil Constantinescu   TS_RK           *rk  = (TS_RK*)ts->data;
777f68a32c8SEmil Constantinescu   RKTableau       tab  = rk->tableau;
778f68a32c8SEmil Constantinescu   const PetscInt  s = tab->s;
779fecfb714SLisandro Dalcin   const PetscReal *A = tab->A,*c = tab->c;
780f68a32c8SEmil Constantinescu   PetscScalar     *w = rk->work;
781f68a32c8SEmil Constantinescu   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
782d1905564SLisandro Dalcin   PetscBool       FSAL = tab->FSAL;
783f68a32c8SEmil Constantinescu   TSAdapt         adapt;
784fecfb714SLisandro Dalcin   PetscInt        i,j;
785be5899b3SLisandro Dalcin   PetscInt        rejections = 0;
786be5899b3SLisandro Dalcin   PetscBool       stageok,accept = PETSC_TRUE;
787be5899b3SLisandro Dalcin   PetscReal       next_time_step = ts->time_step;
788f68a32c8SEmil Constantinescu 
789f68a32c8SEmil Constantinescu   PetscFunctionBegin;
790d1905564SLisandro Dalcin   if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE;
7919566063dSJacob Faibussowitsch   if (FSAL) PetscCall(VecCopy(YdotRHS[s-1],YdotRHS[0]));
792d1905564SLisandro Dalcin 
793f68a32c8SEmil Constantinescu   rk->status = TS_STEP_INCOMPLETE;
794be5899b3SLisandro Dalcin   while (!ts->reason && rk->status != TS_STEP_COMPLETE) {
795be5899b3SLisandro Dalcin     PetscReal t = ts->ptime;
796f68a32c8SEmil Constantinescu     PetscReal h = ts->time_step;
797f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
7989be3e283SDebojyoti Ghosh       rk->stage_time = t + h*c[i];
7999566063dSJacob Faibussowitsch       PetscCall(TSPreStage(ts,rk->stage_time));
8009566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vec_sol,Y[i]));
801f68a32c8SEmil Constantinescu       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
8029566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(Y[i],i,w,YdotRHS));
8039566063dSJacob Faibussowitsch       PetscCall(TSPostStage(ts,rk->stage_time,i,Y));
8049566063dSJacob Faibussowitsch       PetscCall(TSGetAdapt(ts,&adapt));
8059566063dSJacob Faibussowitsch       PetscCall(TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok));
806fecfb714SLisandro Dalcin       if (!stageok) goto reject_step;
8078f738a7cSLisandro Dalcin       if (FSAL && !i) continue;
8089566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]));
809f68a32c8SEmil Constantinescu     }
810be5899b3SLisandro Dalcin 
811be5899b3SLisandro Dalcin     rk->status = TS_STEP_INCOMPLETE;
8129566063dSJacob Faibussowitsch     PetscCall(TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL));
813f68a32c8SEmil Constantinescu     rk->status = TS_STEP_PENDING;
8149566063dSJacob Faibussowitsch     PetscCall(TSGetAdapt(ts,&adapt));
8159566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidatesClear(adapt));
8169566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE));
8179566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept));
818be5899b3SLisandro Dalcin     rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
819be5899b3SLisandro Dalcin     if (!accept) { /* Roll back the current step */
8209566063dSJacob Faibussowitsch       PetscCall(TSRollBack_RK(ts));
821be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
822be5899b3SLisandro Dalcin       goto reject_step;
823be5899b3SLisandro Dalcin     }
824be5899b3SLisandro Dalcin 
8250f7a1166SHong Zhang     if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
8260f7a1166SHong Zhang       rk->ptime     = ts->ptime;
8270f7a1166SHong Zhang       rk->time_step = ts->time_step;
828493ed95fSHong Zhang     }
829be5899b3SLisandro Dalcin 
830f68a32c8SEmil Constantinescu     ts->ptime += ts->time_step;
831f68a32c8SEmil Constantinescu     ts->time_step = next_time_step;
832f68a32c8SEmil Constantinescu     break;
833be5899b3SLisandro Dalcin 
834be5899b3SLisandro Dalcin     reject_step:
835fecfb714SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
836be5899b3SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
837be5899b3SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
8389566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections));
839f68a32c8SEmil Constantinescu     }
840f68a32c8SEmil Constantinescu   }
841f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
842e4dd521cSBarry Smith }
843e4dd521cSBarry Smith 
844f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts)
845f6a906c0SBarry Smith {
846f6a906c0SBarry Smith   TS_RK          *rk  = (TS_RK*)ts->data;
847be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
848be5899b3SLisandro Dalcin   PetscInt       s   = tab->s;
849f6a906c0SBarry Smith 
850f6a906c0SBarry Smith   PetscFunctionBegin;
851f6a906c0SBarry Smith   if (ts->adjointsetupcalled++) PetscFunctionReturn(0);
8529566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecsDeltaLam));
8539566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&rk->VecsSensiTemp));
854f6a906c0SBarry Smith   if (ts->vecs_sensip) {
8559566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vecs_sensip[0],&rk->VecDeltaMu));
856f6a906c0SBarry Smith   }
85713af1a74SHong Zhang   if (ts->vecs_sensi2) {
8589566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecsDeltaLam2));
8599566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&rk->VecsSensi2Temp));
86013af1a74SHong Zhang   }
86113af1a74SHong Zhang   if (ts->vecs_sensi2p) {
8629566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vecs_sensi2p[0],&rk->VecDeltaMu2));
86313af1a74SHong Zhang   }
864f6a906c0SBarry Smith   PetscFunctionReturn(0);
865f6a906c0SBarry Smith }
866f6a906c0SBarry Smith 
86735f5172aSHong Zhang /*
86835f5172aSHong Zhang   Assumptions:
86935f5172aSHong Zhang     - TSStep_RK() always evaluates the step with b, not bembed.
87035f5172aSHong Zhang */
87142f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts)
872d2daff3dSHong Zhang {
873c235aa8dSHong Zhang   TS_RK            *rk = (TS_RK*)ts->data;
874cd4cee2dSHong Zhang   TS               quadts = ts->quadraturets;
875c235aa8dSHong Zhang   RKTableau        tab = rk->tableau;
8763ca0882eSHong Zhang   Mat              J,Jpre,Jquad;
877c235aa8dSHong Zhang   const PetscInt   s = tab->s;
878c235aa8dSHong Zhang   const PetscReal  *A = tab->A,*b = tab->b,*c = tab->c;
87913af1a74SHong Zhang   PetscScalar      *w = rk->work,*xarr;
8802e7b7f96SHong Zhang   Vec              *Y = rk->Y,*VecsDeltaLam = rk->VecsDeltaLam,VecDeltaMu = rk->VecDeltaMu,*VecsSensiTemp = rk->VecsSensiTemp;
88113af1a74SHong Zhang   Vec              *VecsDeltaLam2 = rk->VecsDeltaLam2,VecDeltaMu2 = rk->VecDeltaMu2,*VecsSensi2Temp = rk->VecsSensi2Temp;
882cd4cee2dSHong Zhang   Vec              VecDRDUTransCol = ts->vec_drdu_col,VecDRDPTransCol = ts->vec_drdp_col;
883b18ea86cSHong Zhang   PetscInt         i,j,nadj;
8843d3eaba7SBarry Smith   PetscReal        t = ts->ptime;
8853ca0882eSHong Zhang   PetscReal        h = ts->time_step;
886c235aa8dSHong Zhang 
887d2daff3dSHong Zhang   PetscFunctionBegin;
888c235aa8dSHong Zhang   rk->status = TS_STEP_INCOMPLETE;
8893389c7d9SStefano Zampini 
8909566063dSJacob Faibussowitsch   PetscCall(TSGetRHSJacobian(ts,&J,&Jpre,NULL,NULL));
89135f5172aSHong Zhang   if (quadts) {
8929566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobian(quadts,&Jquad,NULL,NULL,NULL));
89335f5172aSHong Zhang   }
8942e7b7f96SHong Zhang   for (i=s-1; i>=0; i--) {
8956a1e4597SHong Zhang     if (tab->FSAL && i == s-1) {
8966a1e4597SHong Zhang       /* VecsDeltaLam[nadj*s+s-1] are initialized with zeros and the values never change.*/
8976a1e4597SHong Zhang       continue;
8986a1e4597SHong Zhang     }
8993ca0882eSHong Zhang     rk->stage_time = t + h*(1.0-c[i]);
9009566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts,Y[i],J,Jpre));
90135f5172aSHong Zhang     if (quadts) {
9029566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts,rk->stage_time,Y[i],Jquad,Jquad)); /* get r_u^T */
90335f5172aSHong Zhang     }
9043389c7d9SStefano Zampini     if (ts->vecs_sensip) {
9059566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(ts,rk->stage_time,Y[i],ts->Jacprhs)); /* get f_p */
90635f5172aSHong Zhang       if (quadts) {
9079566063dSJacob Faibussowitsch         PetscCall(TSComputeRHSJacobianP(quadts,rk->stage_time,Y[i],quadts->Jacprhs)); /* get f_p for the quadrature */
9083389c7d9SStefano Zampini       }
90935f5172aSHong Zhang     }
9103389c7d9SStefano Zampini 
91135f5172aSHong Zhang     if (b[i]) {
91235f5172aSHong Zhang       for (j=i+1; j<s; j++) w[j-i-1] = A[j*s+i]/b[i]; /* coefficients for computing VecsSensiTemp */
91335f5172aSHong Zhang     } else {
91435f5172aSHong Zhang       for (j=i+1; j<s; j++) w[j-i-1] = A[j*s+i]; /* coefficients for computing VecsSensiTemp */
91535f5172aSHong Zhang     }
9162e7b7f96SHong Zhang 
9172e7b7f96SHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
9183389c7d9SStefano Zampini       /* Stage values of lambda */
91935f5172aSHong Zhang       if (b[i]) {
92035f5172aSHong Zhang         /* lambda_{n+1} + \sum_{j=i+1}^s a_{ji}/b[i]*lambda_{s,j} */
9219566063dSJacob Faibussowitsch         PetscCall(VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj])); /* VecDeltaLam is an vec array of size s by numcost */
9229566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(VecsSensiTemp[nadj],s-i-1,w,&VecsDeltaLam[nadj*s+i+1]));
9239566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(J,VecsSensiTemp[nadj],VecsDeltaLam[nadj*s+i])); /* VecsSensiTemp will be reused by 2nd-order adjoint */
9249566063dSJacob Faibussowitsch         PetscCall(VecScale(VecsDeltaLam[nadj*s+i],-h*b[i]));
92535f5172aSHong Zhang         if (quadts) {
9269566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(Jquad,nadj,&xarr));
9279566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(VecDRDUTransCol,xarr));
9289566063dSJacob Faibussowitsch           PetscCall(VecAXPY(VecsDeltaLam[nadj*s+i],-h*b[i],VecDRDUTransCol));
9299566063dSJacob Faibussowitsch           PetscCall(VecResetArray(VecDRDUTransCol));
9309566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(Jquad,&xarr));
931cd4cee2dSHong Zhang         }
9323389c7d9SStefano Zampini       } else {
9336a1e4597SHong Zhang         /* \sum_{j=i+1}^s a_{ji}*lambda_{s,j} */
9349566063dSJacob Faibussowitsch         PetscCall(VecSet(VecsSensiTemp[nadj],0));
9359566063dSJacob Faibussowitsch         PetscCall(VecMAXPY(VecsSensiTemp[nadj],s-i-1,w,&VecsDeltaLam[nadj*s+i+1]));
9369566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(J,VecsSensiTemp[nadj],VecsDeltaLam[nadj*s+i]));
9379566063dSJacob Faibussowitsch         PetscCall(VecScale(VecsDeltaLam[nadj*s+i],-h));
9383389c7d9SStefano Zampini       }
9396497ce14SHong Zhang 
940ad8e2604SHong Zhang       /* Stage values of mu */
9416497ce14SHong Zhang       if (ts->vecs_sensip) {
94235f5172aSHong Zhang         if (b[i]) {
9439566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(ts->Jacprhs,VecsSensiTemp[nadj],VecDeltaMu));
9449566063dSJacob Faibussowitsch           PetscCall(VecScale(VecDeltaMu,-h*b[i]));
94535f5172aSHong Zhang           if (quadts) {
9469566063dSJacob Faibussowitsch             PetscCall(MatDenseGetColumn(quadts->Jacprhs,nadj,&xarr));
9479566063dSJacob Faibussowitsch             PetscCall(VecPlaceArray(VecDRDPTransCol,xarr));
9489566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecDeltaMu,-h*b[i],VecDRDPTransCol));
9499566063dSJacob Faibussowitsch             PetscCall(VecResetArray(VecDRDPTransCol));
9509566063dSJacob Faibussowitsch             PetscCall(MatDenseRestoreColumn(quadts->Jacprhs,&xarr));
951493ed95fSHong Zhang           }
95235f5172aSHong Zhang         } else {
9539566063dSJacob Faibussowitsch           PetscCall(VecScale(VecDeltaMu,-h));
95435f5172aSHong Zhang         }
9559566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj],1.,VecDeltaMu)); /* update sensip for each stage */
956ad8e2604SHong Zhang       }
957c235aa8dSHong Zhang     }
95813af1a74SHong Zhang 
95913af1a74SHong Zhang     if (ts->vecs_sensi2 && ts->forward_solve) { /* 2nd-order adjoint, TLM mode has to be turned on */
96013af1a74SHong Zhang       /* Get w1 at t_{n+1} from TLM matrix */
9619566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(rk->MatsFwdStageSensip[i],0,&xarr));
9629566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_sensip_col,xarr));
96313af1a74SHong Zhang       /* lambda_s^T F_UU w_1 */
9649566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSHessianProductFunctionUU(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_sensip_col,ts->vecs_guu));
96535f5172aSHong Zhang       if (quadts)  {
96635f5172aSHong Zhang         /* R_UU w_1 */
9679566063dSJacob Faibussowitsch         PetscCall(TSComputeRHSHessianProductFunctionUU(quadts,rk->stage_time,Y[i],NULL,ts->vec_sensip_col,ts->vecs_guu));
96835f5172aSHong Zhang       }
96935f5172aSHong Zhang       if (ts->vecs_sensip) {
97013af1a74SHong Zhang         /* lambda_s^T F_UP w_2 */
9719566063dSJacob Faibussowitsch         PetscCall(TSComputeRHSHessianProductFunctionUP(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_dir,ts->vecs_gup));
97235f5172aSHong Zhang         if (quadts)  {
97335f5172aSHong Zhang           /* R_UP w_2 */
9749566063dSJacob Faibussowitsch           PetscCall(TSComputeRHSHessianProductFunctionUP(quadts,rk->stage_time,Y[i],NULL,ts->vec_sensip_col,ts->vecs_gup));
97535f5172aSHong Zhang         }
97635f5172aSHong Zhang       }
97713af1a74SHong Zhang       if (ts->vecs_sensi2p) {
97813af1a74SHong Zhang         /* lambda_s^T F_PU w_1 */
9799566063dSJacob Faibussowitsch         PetscCall(TSComputeRHSHessianProductFunctionPU(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_sensip_col,ts->vecs_gpu));
98035f5172aSHong Zhang         /* lambda_s^T F_PP w_2 */
9819566063dSJacob Faibussowitsch         PetscCall(TSComputeRHSHessianProductFunctionPP(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_dir,ts->vecs_gpp));
98235f5172aSHong Zhang         if (b[i] && quadts) {
98335f5172aSHong Zhang           /* R_PU w_1 */
9849566063dSJacob Faibussowitsch           PetscCall(TSComputeRHSHessianProductFunctionPU(quadts,rk->stage_time,Y[i],NULL,ts->vec_sensip_col,ts->vecs_gpu));
98535f5172aSHong Zhang           /* R_PP w_2 */
9869566063dSJacob Faibussowitsch           PetscCall(TSComputeRHSHessianProductFunctionPP(quadts,rk->stage_time,Y[i],NULL,ts->vec_dir,ts->vecs_gpp));
98735f5172aSHong Zhang         }
98813af1a74SHong Zhang       }
9899566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_sensip_col));
9909566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(rk->MatsFwdStageSensip[i],&xarr));
99113af1a74SHong Zhang 
99213af1a74SHong Zhang       for (nadj=0; nadj<ts->numcost; nadj++) {
99313af1a74SHong Zhang         /* Stage values of lambda */
99435f5172aSHong Zhang         if (b[i]) {
99535f5172aSHong Zhang           /* J_i^T*(Lambda_{n+1}+\sum_{j=i+1}^s a_{ji}/b_i*Lambda_{s,j} */
9969566063dSJacob Faibussowitsch           PetscCall(VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]));
9979566063dSJacob Faibussowitsch           PetscCall(VecMAXPY(VecsSensi2Temp[nadj],s-i-1,w,&VecsDeltaLam2[nadj*s+i+1]));
9989566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(J,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj*s+i]));
9999566063dSJacob Faibussowitsch           PetscCall(VecScale(VecsDeltaLam2[nadj*s+i],-h*b[i]));
10009566063dSJacob Faibussowitsch           PetscCall(VecAXPY(VecsDeltaLam2[nadj*s+i],-h*b[i],ts->vecs_guu[nadj]));
100135f5172aSHong Zhang           if (ts->vecs_sensip) {
10029566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecsDeltaLam2[nadj*s+i],-h*b[i],ts->vecs_gup[nadj]));
100313af1a74SHong Zhang           }
100413af1a74SHong Zhang         } else {
100535f5172aSHong Zhang           /* \sum_{j=i+1}^s a_{ji}*Lambda_{s,j} */
10069566063dSJacob Faibussowitsch           PetscCall(VecSet(VecsDeltaLam2[nadj*s+i],0));
10079566063dSJacob Faibussowitsch           PetscCall(VecMAXPY(VecsSensi2Temp[nadj],s-i-1,w,&VecsDeltaLam2[nadj*s+i+1]));
10089566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(J,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj*s+i]));
10099566063dSJacob Faibussowitsch           PetscCall(VecScale(VecsDeltaLam2[nadj*s+i],-h));
10109566063dSJacob Faibussowitsch           PetscCall(VecAXPY(VecsDeltaLam2[nadj*s+i],-h,ts->vecs_guu[nadj]));
101135f5172aSHong Zhang           if (ts->vecs_sensip) {
10129566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecsDeltaLam2[nadj*s+i],-h,ts->vecs_gup[nadj]));
101313af1a74SHong Zhang           }
101435f5172aSHong Zhang         }
101535f5172aSHong Zhang         if (ts->vecs_sensi2p) { /* 2nd-order adjoint for parameters */
10169566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(ts->Jacprhs,VecsSensi2Temp[nadj],VecDeltaMu2));
101735f5172aSHong Zhang           if (b[i]) {
10189566063dSJacob Faibussowitsch             PetscCall(VecScale(VecDeltaMu2,-h*b[i]));
10199566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecDeltaMu2,-h*b[i],ts->vecs_gpu[nadj]));
10209566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecDeltaMu2,-h*b[i],ts->vecs_gpp[nadj]));
102113af1a74SHong Zhang           } else {
10229566063dSJacob Faibussowitsch             PetscCall(VecScale(VecDeltaMu2,-h));
10239566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecDeltaMu2,-h,ts->vecs_gpu[nadj]));
10249566063dSJacob Faibussowitsch             PetscCall(VecAXPY(VecDeltaMu2,-h,ts->vecs_gpp[nadj]));
102513af1a74SHong Zhang           }
10269566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensi2p[nadj],1,VecDeltaMu2)); /* update sensi2p for each stage */
102713af1a74SHong Zhang         }
102813af1a74SHong Zhang       }
102913af1a74SHong Zhang     }
10306497ce14SHong Zhang   }
1031c235aa8dSHong Zhang 
1032c235aa8dSHong Zhang   for (j=0; j<s; j++) w[j] = 1.0;
10332e7b7f96SHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) { /* no need to do this for mu's */
10349566063dSJacob Faibussowitsch     PetscCall(VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecsDeltaLam[nadj*s]));
103513af1a74SHong Zhang     if (ts->vecs_sensi2) {
10369566063dSJacob Faibussowitsch       PetscCall(VecMAXPY(ts->vecs_sensi2[nadj],s,w,&VecsDeltaLam2[nadj*s]));
103713af1a74SHong Zhang     }
10386497ce14SHong Zhang   }
1039c235aa8dSHong Zhang   rk->status = TS_STEP_COMPLETE;
1040d2daff3dSHong Zhang   PetscFunctionReturn(0);
1041d2daff3dSHong Zhang }
1042d2daff3dSHong Zhang 
104313af1a74SHong Zhang static PetscErrorCode TSAdjointReset_RK(TS ts)
104413af1a74SHong Zhang {
104513af1a74SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
104613af1a74SHong Zhang   RKTableau      tab = rk->tableau;
104713af1a74SHong Zhang 
104813af1a74SHong Zhang   PetscFunctionBegin;
10499566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam));
10509566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost,&rk->VecsSensiTemp));
10519566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&rk->VecDeltaMu));
10529566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam2));
10539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&rk->VecDeltaMu2));
10549566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost,&rk->VecsSensi2Temp));
105513af1a74SHong Zhang   PetscFunctionReturn(0);
105613af1a74SHong Zhang }
105713af1a74SHong Zhang 
105855de54fcSHong Zhang static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X)
105955de54fcSHong Zhang {
106055de54fcSHong Zhang   TS_RK            *rk = (TS_RK*)ts->data;
106155de54fcSHong Zhang   PetscInt         s  = rk->tableau->s,p = rk->tableau->p,i,j;
106255de54fcSHong Zhang   PetscReal        h;
106355de54fcSHong Zhang   PetscReal        tt,t;
106455de54fcSHong Zhang   PetscScalar      *b;
106555de54fcSHong Zhang   const PetscReal  *B = rk->tableau->binterp;
106655de54fcSHong Zhang 
106755de54fcSHong Zhang   PetscFunctionBegin;
10683c633725SBarry Smith   PetscCheck(B,PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
106955de54fcSHong Zhang 
107055de54fcSHong Zhang   switch (rk->status) {
107155de54fcSHong Zhang     case TS_STEP_INCOMPLETE:
107255de54fcSHong Zhang     case TS_STEP_PENDING:
107355de54fcSHong Zhang       h = ts->time_step;
107455de54fcSHong Zhang       t = (itime - ts->ptime)/h;
107555de54fcSHong Zhang       break;
107655de54fcSHong Zhang     case TS_STEP_COMPLETE:
107755de54fcSHong Zhang       h = ts->ptime - ts->ptime_prev;
107855de54fcSHong Zhang       t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
107955de54fcSHong Zhang       break;
108055de54fcSHong Zhang     default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
108155de54fcSHong Zhang   }
10829566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(s,&b));
108355de54fcSHong Zhang   for (i=0; i<s; i++) b[i] = 0;
108455de54fcSHong Zhang   for (j=0,tt=t; j<p; j++,tt*=t) {
108555de54fcSHong Zhang     for (i=0; i<s; i++) {
108655de54fcSHong Zhang       b[i]  += h * B[i*p+j] * tt;
108755de54fcSHong Zhang     }
108855de54fcSHong Zhang   }
10899566063dSJacob Faibussowitsch   PetscCall(VecCopy(rk->Y[0],X));
10909566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(X,s,b,rk->YdotRHS));
10919566063dSJacob Faibussowitsch   PetscCall(PetscFree(b));
109255de54fcSHong Zhang   PetscFunctionReturn(0);
109355de54fcSHong Zhang }
109455de54fcSHong Zhang 
109555de54fcSHong Zhang /*------------------------------------------------------------*/
109655de54fcSHong Zhang 
1097be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts)
1098be5899b3SLisandro Dalcin {
1099be5899b3SLisandro Dalcin   TS_RK          *rk = (TS_RK*)ts->data;
1100be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
1101be5899b3SLisandro Dalcin 
1102be5899b3SLisandro Dalcin   PetscFunctionBegin;
1103be5899b3SLisandro Dalcin   if (!tab) PetscFunctionReturn(0);
11049566063dSJacob Faibussowitsch   PetscCall(PetscFree(rk->work));
11059566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(tab->s,&rk->Y));
11069566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(tab->s,&rk->YdotRHS));
1107be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
1108be5899b3SLisandro Dalcin }
1109be5899b3SLisandro Dalcin 
1110277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts)
1111e4dd521cSBarry Smith {
1112e4dd521cSBarry Smith   PetscFunctionBegin;
11139566063dSJacob Faibussowitsch   PetscCall(TSRKTableauReset(ts));
11140fe4d17eSHong Zhang   if (ts->use_splitrhsfunction) {
1115*cac4c232SBarry Smith     PetscTryMethod(ts,"TSReset_RK_MultirateSplit_C",(TS),(ts));
11160fe4d17eSHong Zhang   } else {
1117*cac4c232SBarry Smith     PetscTryMethod(ts,"TSReset_RK_MultirateNonsplit_C",(TS),(ts));
11180fe4d17eSHong Zhang   }
1119277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
1120e4dd521cSBarry Smith }
1121277b19d0SLisandro Dalcin 
1122f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx)
1123f68a32c8SEmil Constantinescu {
1124f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1125f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1126f68a32c8SEmil Constantinescu }
1127f68a32c8SEmil Constantinescu 
1128f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
1129f68a32c8SEmil Constantinescu {
1130f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1131f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1132f68a32c8SEmil Constantinescu }
1133f68a32c8SEmil Constantinescu 
1134f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx)
1135f68a32c8SEmil Constantinescu {
1136f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1137f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1138f68a32c8SEmil Constantinescu }
1139f68a32c8SEmil Constantinescu 
1140f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
1141f68a32c8SEmil Constantinescu {
1142f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1143f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1144f68a32c8SEmil Constantinescu }
1145be5899b3SLisandro Dalcin 
1146be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts)
1147be5899b3SLisandro Dalcin {
1148be5899b3SLisandro Dalcin   TS_RK          *rk  = (TS_RK*)ts->data;
1149be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
1150be5899b3SLisandro Dalcin 
1151be5899b3SLisandro Dalcin   PetscFunctionBegin;
11529566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(tab->s,&rk->work));
11539566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y));
11549566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS));
1155be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
1156be5899b3SLisandro Dalcin }
1157be5899b3SLisandro Dalcin 
1158f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts)
1159f68a32c8SEmil Constantinescu {
1160cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
1161f68a32c8SEmil Constantinescu   DM             dm;
1162f68a32c8SEmil Constantinescu 
1163f68a32c8SEmil Constantinescu   PetscFunctionBegin;
11649566063dSJacob Faibussowitsch   PetscCall(TSCheckImplicitTerm(ts));
11659566063dSJacob Faibussowitsch   PetscCall(TSRKTableauSetUp(ts));
1166cd4cee2dSHong Zhang   if (quadts && ts->costintegralfwd) {
1167cd4cee2dSHong Zhang     Mat Jquad;
11689566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobian(quadts,&Jquad,NULL,NULL,NULL));
11690f7a1166SHong Zhang   }
11709566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&dm));
11719566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts));
11729566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts));
11730fe4d17eSHong Zhang   if (ts->use_splitrhsfunction) {
1174*cac4c232SBarry Smith     PetscTryMethod(ts,"TSSetUp_RK_MultirateSplit_C",(TS),(ts));
11750fe4d17eSHong Zhang   } else {
1176*cac4c232SBarry Smith     PetscTryMethod(ts,"TSSetUp_RK_MultirateNonsplit_C",(TS),(ts));
11770fe4d17eSHong Zhang   }
1178e4dd521cSBarry Smith   PetscFunctionReturn(0);
1179e4dd521cSBarry Smith }
1180d2daff3dSHong Zhang 
11814416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts)
1182e4dd521cSBarry Smith {
1183be5899b3SLisandro Dalcin   TS_RK          *rk = (TS_RK*)ts->data;
1184dfbe8321SBarry Smith   PetscErrorCode ierr;
1185262ff7bbSBarry Smith 
1186e4dd521cSBarry Smith   PetscFunctionBegin;
11879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHead(PetscOptionsObject,"RK ODE solver options"));
1188f68a32c8SEmil Constantinescu   {
1189f68a32c8SEmil Constantinescu     RKTableauLink link;
1190f68a32c8SEmil Constantinescu     PetscInt      count,choice;
11910fe4d17eSHong Zhang     PetscBool     flg,use_multirate = PETSC_FALSE;
1192f68a32c8SEmil Constantinescu     const char    **namelist;
11932c9cb062Sluzhanghpp 
1194f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) ;
11959566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(count,(char***)&namelist));
1196f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
11979566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ts_rk_multirate","Use interpolation-based multirate RK method","TSRKSetMultirate",rk->use_multirate,&use_multirate,&flg));
11980fe4d17eSHong Zhang     if (flg) {
11999566063dSJacob Faibussowitsch       PetscCall(TSRKSetMultirate(ts,use_multirate));
12000fe4d17eSHong Zhang     }
12019566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg));
12029566063dSJacob Faibussowitsch     if (flg) PetscCall(TSRKSetType(ts,namelist[choice]));
12039566063dSJacob Faibussowitsch     PetscCall(PetscFree(namelist));
1204f68a32c8SEmil Constantinescu   }
12059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsTail());
12069566063dSJacob Faibussowitsch   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)ts),NULL,"Multirate methods options","");PetscCall(ierr);
12079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-ts_rk_dtratio","time step ratio between slow and fast","",rk->dtratio,&rk->dtratio,NULL));
12089566063dSJacob Faibussowitsch   ierr = PetscOptionsEnd();PetscCall(ierr);
1209e4dd521cSBarry Smith   PetscFunctionReturn(0);
1210e4dd521cSBarry Smith }
1211e4dd521cSBarry Smith 
12125f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
1213e4dd521cSBarry Smith {
12145f70b478SJed Brown   TS_RK          *rk = (TS_RK*)ts->data;
12158370ee3bSLisandro Dalcin   PetscBool      iascii;
12162cdf8120Sasbjorn 
1217e4dd521cSBarry Smith   PetscFunctionBegin;
12189566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
12198370ee3bSLisandro Dalcin   if (iascii) {
12209c334d8fSLisandro Dalcin     RKTableau       tab  = rk->tableau;
1221f68a32c8SEmil Constantinescu     TSRKType        rktype;
12220f7a28cdSStefano Zampini     const PetscReal *c;
12230f7a28cdSStefano Zampini     PetscInt        s;
1224f68a32c8SEmil Constantinescu     char            buf[512];
12250f7a28cdSStefano Zampini     PetscBool       FSAL;
12260f7a28cdSStefano Zampini 
12279566063dSJacob Faibussowitsch     PetscCall(TSRKGetType(ts,&rktype));
12289566063dSJacob Faibussowitsch     PetscCall(TSRKGetTableau(ts,&s,NULL,NULL,&c,NULL,NULL,NULL,&FSAL));
12299566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  RK type %s\n",rktype));
12309566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  Order: %D\n",tab->order));
12319566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  FSAL property: %s\n",FSAL ? "yes" : "no"));
12329566063dSJacob Faibussowitsch     PetscCall(PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",s,c));
12339566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  Abscissa c = %s\n",buf));
12348370ee3bSLisandro Dalcin   }
1235f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1236f68a32c8SEmil Constantinescu }
1237f68a32c8SEmil Constantinescu 
1238f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer)
1239f68a32c8SEmil Constantinescu {
12409c334d8fSLisandro Dalcin   TSAdapt        adapt;
1241f68a32c8SEmil Constantinescu 
1242f68a32c8SEmil Constantinescu   PetscFunctionBegin;
12439566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts,&adapt));
12449566063dSJacob Faibussowitsch   PetscCall(TSAdaptLoad(adapt,viewer));
1245f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1246f68a32c8SEmil Constantinescu }
1247f68a32c8SEmil Constantinescu 
12482ea87ba9SLisandro Dalcin /*@
12492ea87ba9SLisandro Dalcin   TSRKGetOrder - Get the order of RK scheme
12502ea87ba9SLisandro Dalcin 
12512ea87ba9SLisandro Dalcin   Not collective
12522ea87ba9SLisandro Dalcin 
12532ea87ba9SLisandro Dalcin   Input Parameter:
12542ea87ba9SLisandro Dalcin .  ts - timestepping context
12552ea87ba9SLisandro Dalcin 
12562ea87ba9SLisandro Dalcin   Output Parameter:
12572ea87ba9SLisandro Dalcin .  order - order of RK-scheme
12582ea87ba9SLisandro Dalcin 
12592ea87ba9SLisandro Dalcin   Level: intermediate
12602ea87ba9SLisandro Dalcin 
12612ea87ba9SLisandro Dalcin .seealso: TSRKGetType()
12622ea87ba9SLisandro Dalcin @*/
12632ea87ba9SLisandro Dalcin PetscErrorCode TSRKGetOrder(TS ts,PetscInt *order)
12642ea87ba9SLisandro Dalcin {
12652ea87ba9SLisandro Dalcin   PetscFunctionBegin;
12662ea87ba9SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
12672ea87ba9SLisandro Dalcin   PetscValidIntPointer(order,2);
1268*cac4c232SBarry Smith   PetscUseMethod(ts,"TSRKGetOrder_C",(TS,PetscInt*),(ts,order));
12692ea87ba9SLisandro Dalcin   PetscFunctionReturn(0);
12702ea87ba9SLisandro Dalcin }
12712ea87ba9SLisandro Dalcin 
1272f68a32c8SEmil Constantinescu /*@C
1273f68a32c8SEmil Constantinescu   TSRKSetType - Set the type of RK scheme
1274f68a32c8SEmil Constantinescu 
1275f68a32c8SEmil Constantinescu   Logically collective
1276f68a32c8SEmil Constantinescu 
1277d8d19677SJose E. Roman   Input Parameters:
1278f68a32c8SEmil Constantinescu +  ts - timestepping context
1279f68a32c8SEmil Constantinescu -  rktype - type of RK-scheme
1280f68a32c8SEmil Constantinescu 
1281287c2655SBarry Smith   Options Database:
12829bd3e852SBarry Smith .   -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs>
1283287c2655SBarry Smith 
1284f68a32c8SEmil Constantinescu   Level: intermediate
1285f68a32c8SEmil Constantinescu 
1286e7685601SHong Zhang .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK2B, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS, TSRK6VR, TSRK7VR, TSRK8VR
1287f68a32c8SEmil Constantinescu @*/
1288f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype)
1289f68a32c8SEmil Constantinescu {
1290f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1291f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1292b92453a8SLisandro Dalcin   PetscValidCharPointer(rktype,2);
1293*cac4c232SBarry Smith   PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));
1294f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1295f68a32c8SEmil Constantinescu }
1296f68a32c8SEmil Constantinescu 
1297f68a32c8SEmil Constantinescu /*@C
1298f68a32c8SEmil Constantinescu   TSRKGetType - Get the type of RK scheme
1299f68a32c8SEmil Constantinescu 
13002ea87ba9SLisandro Dalcin   Not collective
1301f68a32c8SEmil Constantinescu 
1302f68a32c8SEmil Constantinescu   Input Parameter:
1303f68a32c8SEmil Constantinescu .  ts - timestepping context
1304f68a32c8SEmil Constantinescu 
1305f68a32c8SEmil Constantinescu   Output Parameter:
1306f68a32c8SEmil Constantinescu .  rktype - type of RK-scheme
1307f68a32c8SEmil Constantinescu 
1308f68a32c8SEmil Constantinescu   Level: intermediate
1309f68a32c8SEmil Constantinescu 
13102ea87ba9SLisandro Dalcin .seealso: TSRKSetType()
1311f68a32c8SEmil Constantinescu @*/
1312f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype)
1313f68a32c8SEmil Constantinescu {
1314f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1315f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1316*cac4c232SBarry Smith   PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));
1317f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1318f68a32c8SEmil Constantinescu }
1319f68a32c8SEmil Constantinescu 
13202ea87ba9SLisandro Dalcin static PetscErrorCode TSRKGetOrder_RK(TS ts,PetscInt *order)
13212ea87ba9SLisandro Dalcin {
13222ea87ba9SLisandro Dalcin   TS_RK *rk = (TS_RK*)ts->data;
13232ea87ba9SLisandro Dalcin 
13242ea87ba9SLisandro Dalcin   PetscFunctionBegin;
13252ea87ba9SLisandro Dalcin   *order = rk->tableau->order;
13262ea87ba9SLisandro Dalcin   PetscFunctionReturn(0);
13272ea87ba9SLisandro Dalcin }
13282ea87ba9SLisandro Dalcin 
1329560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype)
1330f68a32c8SEmil Constantinescu {
1331f68a32c8SEmil Constantinescu   TS_RK *rk = (TS_RK*)ts->data;
1332f68a32c8SEmil Constantinescu 
1333f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1334f68a32c8SEmil Constantinescu   *rktype = rk->tableau->name;
1335f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
1336f68a32c8SEmil Constantinescu }
13372c9cb062Sluzhanghpp 
1338560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype)
1339f68a32c8SEmil Constantinescu {
1340f68a32c8SEmil Constantinescu   TS_RK          *rk = (TS_RK*)ts->data;
1341f68a32c8SEmil Constantinescu   PetscBool      match;
1342f68a32c8SEmil Constantinescu   RKTableauLink  link;
1343f68a32c8SEmil Constantinescu 
1344f68a32c8SEmil Constantinescu   PetscFunctionBegin;
1345f68a32c8SEmil Constantinescu   if (rk->tableau) {
13469566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(rk->tableau->name,rktype,&match));
1347f68a32c8SEmil Constantinescu     if (match) PetscFunctionReturn(0);
1348f68a32c8SEmil Constantinescu   }
1349f68a32c8SEmil Constantinescu   for (link = RKTableauList; link; link=link->next) {
13509566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(link->tab.name,rktype,&match));
1351f68a32c8SEmil Constantinescu     if (match) {
13529566063dSJacob Faibussowitsch       if (ts->setupcalled) PetscCall(TSRKTableauReset(ts));
1353f68a32c8SEmil Constantinescu       rk->tableau = &link->tab;
13549566063dSJacob Faibussowitsch       if (ts->setupcalled) PetscCall(TSRKTableauSetUp(ts));
13552ffb9264SLisandro Dalcin       ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
1356f68a32c8SEmil Constantinescu       PetscFunctionReturn(0);
1357f68a32c8SEmil Constantinescu     }
1358f68a32c8SEmil Constantinescu   }
135998921bdaSJacob Faibussowitsch   SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype);
1360e4dd521cSBarry Smith }
1361e4dd521cSBarry Smith 
1362ff22ae23SHong Zhang static PetscErrorCode  TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y)
1363ff22ae23SHong Zhang {
1364ff22ae23SHong Zhang   TS_RK *rk = (TS_RK*)ts->data;
1365ff22ae23SHong Zhang 
1366ff22ae23SHong Zhang   PetscFunctionBegin;
13670429704eSStefano Zampini   if (ns) *ns = rk->tableau->s;
1368d2daff3dSHong Zhang   if (Y)   *Y = rk->Y;
1369ff22ae23SHong Zhang   PetscFunctionReturn(0);
1370ff22ae23SHong Zhang }
1371ff22ae23SHong Zhang 
1372b3a6b972SJed Brown static PetscErrorCode TSDestroy_RK(TS ts)
1373b3a6b972SJed Brown {
1374b3a6b972SJed Brown   PetscFunctionBegin;
13759566063dSJacob Faibussowitsch   PetscCall(TSReset_RK(ts));
1376b3a6b972SJed Brown   if (ts->dm) {
13779566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts));
13789566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts));
1379b3a6b972SJed Brown   }
13809566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
13819566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKGetOrder_C",NULL));
13829566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL));
13839566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL));
13849566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKGetTableau_C",NULL));
13859566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",NULL));
13869566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",NULL));
1387b3a6b972SJed Brown   PetscFunctionReturn(0);
1388b3a6b972SJed Brown }
1389ff22ae23SHong Zhang 
13903ca0882eSHong Zhang /*
13913ca0882eSHong Zhang   This defines the nonlinear equation that is to be solved with SNES
13923ca0882eSHong Zhang   We do not need to solve the equation; we just use SNES to approximate the Jacobian
13933ca0882eSHong Zhang */
13943ca0882eSHong Zhang static PetscErrorCode SNESTSFormFunction_RK(SNES snes,Vec x,Vec y,TS ts)
13953ca0882eSHong Zhang {
13963ca0882eSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
13973ca0882eSHong Zhang   DM             dm,dmsave;
13983ca0882eSHong Zhang 
13993ca0882eSHong Zhang   PetscFunctionBegin;
14009566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
14013ca0882eSHong Zhang   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
14023ca0882eSHong Zhang   dmsave = ts->dm;
14033ca0882eSHong Zhang   ts->dm = dm;
14049566063dSJacob Faibussowitsch   PetscCall(TSComputeRHSFunction(ts,rk->stage_time,x,y));
14053ca0882eSHong Zhang   ts->dm = dmsave;
14063ca0882eSHong Zhang   PetscFunctionReturn(0);
14073ca0882eSHong Zhang }
14083ca0882eSHong Zhang 
14093ca0882eSHong Zhang static PetscErrorCode SNESTSFormJacobian_RK(SNES snes,Vec x,Mat A,Mat B,TS ts)
14103ca0882eSHong Zhang {
14113ca0882eSHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
14123ca0882eSHong Zhang   DM             dm,dmsave;
14133ca0882eSHong Zhang 
14143ca0882eSHong Zhang   PetscFunctionBegin;
14159566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
14163ca0882eSHong Zhang   dmsave = ts->dm;
14173ca0882eSHong Zhang   ts->dm = dm;
14189566063dSJacob Faibussowitsch   PetscCall(TSComputeRHSJacobian(ts,rk->stage_time,x,A,B));
14193ca0882eSHong Zhang   ts->dm = dmsave;
14203ca0882eSHong Zhang   PetscFunctionReturn(0);
14213ca0882eSHong Zhang }
14223ca0882eSHong Zhang 
142321052c0cSHong Zhang /*@C
142421052c0cSHong Zhang   TSRKSetMultirate - Use the interpolation-based multirate RK method
142521052c0cSHong Zhang 
142621052c0cSHong Zhang   Logically collective
142721052c0cSHong Zhang 
1428d8d19677SJose E. Roman   Input Parameters:
142921052c0cSHong Zhang +  ts - timestepping context
143021052c0cSHong 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
143121052c0cSHong Zhang 
143221052c0cSHong Zhang   Options Database:
143321052c0cSHong Zhang .   -ts_rk_multirate - <true,false>
143421052c0cSHong Zhang 
143521052c0cSHong Zhang   Notes:
143621052c0cSHong 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).
143721052c0cSHong Zhang 
143821052c0cSHong Zhang   Level: intermediate
143921052c0cSHong Zhang 
144021052c0cSHong Zhang .seealso: TSRKGetMultirate()
144121052c0cSHong Zhang @*/
144221052c0cSHong Zhang PetscErrorCode TSRKSetMultirate(TS ts,PetscBool use_multirate)
144321052c0cSHong Zhang {
14448a4be4dbSHong Zhang   PetscFunctionBegin;
1445*cac4c232SBarry Smith   PetscTryMethod(ts,"TSRKSetMultirate_C",(TS,PetscBool),(ts,use_multirate));
144621052c0cSHong Zhang   PetscFunctionReturn(0);
144721052c0cSHong Zhang }
144821052c0cSHong Zhang 
144921052c0cSHong Zhang /*@C
145021052c0cSHong Zhang   TSRKGetMultirate - Gets whether to Use the interpolation-based multirate RK method
145121052c0cSHong Zhang 
145221052c0cSHong Zhang   Not collective
145321052c0cSHong Zhang 
145421052c0cSHong Zhang   Input Parameter:
145521052c0cSHong Zhang .  ts - timestepping context
145621052c0cSHong Zhang 
145721052c0cSHong Zhang   Output Parameter:
145821052c0cSHong Zhang .  use_multirate - PETSC_TRUE if the multirate RK method is enabled, PETSC_FALSE otherwise
145921052c0cSHong Zhang 
146021052c0cSHong Zhang   Level: intermediate
146121052c0cSHong Zhang 
146221052c0cSHong Zhang .seealso: TSRKSetMultirate()
146321052c0cSHong Zhang @*/
146421052c0cSHong Zhang PetscErrorCode TSRKGetMultirate(TS ts,PetscBool *use_multirate)
146521052c0cSHong Zhang {
14667dbacdf2SHong Zhang   PetscFunctionBegin;
1467*cac4c232SBarry Smith   PetscUseMethod(ts,"TSRKGetMultirate_C",(TS,PetscBool*),(ts,use_multirate));
146821052c0cSHong Zhang   PetscFunctionReturn(0);
146921052c0cSHong Zhang }
147021052c0cSHong Zhang 
1471ebe3b25bSBarry Smith /*MC
1472f68a32c8SEmil Constantinescu       TSRK - ODE and DAE solver using Runge-Kutta schemes
1473ebe3b25bSBarry Smith 
1474f68a32c8SEmil Constantinescu   The user should provide the right hand side of the equation
1475f68a32c8SEmil Constantinescu   using TSSetRHSFunction().
1476ebe3b25bSBarry Smith 
1477f68a32c8SEmil Constantinescu   Notes:
147898164e13SLisandro Dalcin   The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type
1479ebe3b25bSBarry Smith 
1480d41222bbSBarry Smith   Level: beginner
1481d41222bbSBarry Smith 
14825d1808f1SStefano Zampini .seealso:  TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRK2D, TTSRK2E, TSRK3,
14830fe4d17eSHong Zhang            TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister(), TSRKSetMultirate(), TSRKGetMultirate()
1484ebe3b25bSBarry Smith 
1485ebe3b25bSBarry Smith M*/
14868cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
1487e4dd521cSBarry Smith {
14881566a47fSLisandro Dalcin   TS_RK          *rk;
1489e4dd521cSBarry Smith 
1490e4dd521cSBarry Smith   PetscFunctionBegin;
14919566063dSJacob Faibussowitsch   PetscCall(TSRKInitializePackage());
1492f68a32c8SEmil Constantinescu 
1493f68a32c8SEmil Constantinescu   ts->ops->reset          = TSReset_RK;
14945f70b478SJed Brown   ts->ops->destroy        = TSDestroy_RK;
14955f70b478SJed Brown   ts->ops->view           = TSView_RK;
1496f68a32c8SEmil Constantinescu   ts->ops->load           = TSLoad_RK;
1497f68a32c8SEmil Constantinescu   ts->ops->setup          = TSSetUp_RK;
1498f68a32c8SEmil Constantinescu   ts->ops->interpolate    = TSInterpolate_RK;
14992c9cb062Sluzhanghpp   ts->ops->step           = TSStep_RK;
1500f68a32c8SEmil Constantinescu   ts->ops->evaluatestep   = TSEvaluateStep_RK;
1501fecfb714SLisandro Dalcin   ts->ops->rollback       = TSRollBack_RK;
1502f68a32c8SEmil Constantinescu   ts->ops->setfromoptions = TSSetFromOptions_RK;
1503ff22ae23SHong Zhang   ts->ops->getstages      = TSGetStages_RK;
1504e4dd521cSBarry Smith 
15053ca0882eSHong Zhang   ts->ops->snesfunction    = SNESTSFormFunction_RK;
15063ca0882eSHong Zhang   ts->ops->snesjacobian    = SNESTSFormJacobian_RK;
15070f7a1166SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_RK;
150813af1a74SHong Zhang   ts->ops->adjointsetup    = TSAdjointSetUp_RK;
150913af1a74SHong Zhang   ts->ops->adjointstep     = TSAdjointStep_RK;
151013af1a74SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_RK;
15110f7a1166SHong Zhang 
151213af1a74SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_RK;
1513922a638cSHong Zhang   ts->ops->forwardsetup    = TSForwardSetUp_RK;
1514922a638cSHong Zhang   ts->ops->forwardreset    = TSForwardReset_RK;
1515922a638cSHong Zhang   ts->ops->forwardstep     = TSForwardStep_RK;
1516922a638cSHong Zhang   ts->ops->forwardgetstages= TSForwardGetStages_RK;
1517922a638cSHong Zhang 
15189566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(ts,&rk));
15191566a47fSLisandro Dalcin   ts->data = (void*)rk;
1520e4dd521cSBarry Smith 
15219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKGetOrder_C",TSRKGetOrder_RK));
15229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK));
15239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK));
15249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKGetTableau_C",TSRKGetTableau_RK));
15259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",TSRKSetMultirate_RK));
15269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",TSRKGetMultirate_RK));
1527be5899b3SLisandro Dalcin 
15289566063dSJacob Faibussowitsch   PetscCall(TSRKSetType(ts,TSRKDefault));
152990b67129SHong Zhang   rk->dtratio = 1;
1530e4dd521cSBarry Smith   PetscFunctionReturn(0);
1531e4dd521cSBarry Smith }
1532