1e4dd521cSBarry Smith /* 2474dd773SHong Zhang Code for time stepping with the Runge-Kutta method 3f68a32c8SEmil Constantinescu 4f68a32c8SEmil Constantinescu Notes: 5474dd773SHong Zhang The general system is written as 6474dd773SHong Zhang 7474dd773SHong Zhang Udot = F(t,U) 8474dd773SHong Zhang 9e4dd521cSBarry Smith */ 102c9cb062Sluzhanghpp 11af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 12f68a32c8SEmil Constantinescu #include <petscdm.h> 13474dd773SHong Zhang #include <../src/ts/impls/explicit/rk/rk.h> 1421052c0cSHong Zhang #include <../src/ts/impls/explicit/rk/mrk.h> 15f68a32c8SEmil Constantinescu 16484bcdc7SDebojyoti Ghosh static TSRKType TSRKDefault = TSRK3BS; 17f68a32c8SEmil Constantinescu static PetscBool TSRKRegisterAllCalled; 18f68a32c8SEmil Constantinescu static PetscBool TSRKPackageInitialized; 19f68a32c8SEmil Constantinescu 20ab8985f5SHong Zhang static RKTableauLink RKTableauList; 21ab8985f5SHong Zhang 22f68a32c8SEmil Constantinescu /*MC 23287c2655SBarry Smith TSRK1FE - First order forward Euler scheme. 24262ff7bbSBarry Smith 25f68a32c8SEmil Constantinescu This method has one stage. 26f68a32c8SEmil Constantinescu 27bcf0153eSBarry Smith Options Database Key: 2867b8a455SSatish Balay . -ts_rk_type 1fe - use type 1fe 29287c2655SBarry Smith 30f68a32c8SEmil Constantinescu Level: advanced 31f68a32c8SEmil Constantinescu 321cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 33f68a32c8SEmil Constantinescu M*/ 34f68a32c8SEmil Constantinescu /*MC 35e7685601SHong Zhang TSRK2A - Second order RK scheme (Heun's method). 36f68a32c8SEmil Constantinescu 37f68a32c8SEmil Constantinescu This method has two stages. 38f68a32c8SEmil Constantinescu 39bcf0153eSBarry Smith Options Database Key: 4067b8a455SSatish Balay . -ts_rk_type 2a - use type 2a 41287c2655SBarry Smith 42f68a32c8SEmil Constantinescu Level: advanced 43f68a32c8SEmil Constantinescu 441cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 45f68a32c8SEmil Constantinescu M*/ 46f68a32c8SEmil Constantinescu /*MC 47e7685601SHong Zhang TSRK2B - Second order RK scheme (the midpoint method). 48e7685601SHong Zhang 49e7685601SHong Zhang This method has two stages. 50e7685601SHong Zhang 51bcf0153eSBarry Smith Options Database Key: 5267b8a455SSatish Balay . -ts_rk_type 2b - use type 2b 53e7685601SHong Zhang 54e7685601SHong Zhang Level: advanced 55e7685601SHong Zhang 561cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 57e7685601SHong Zhang M*/ 58e7685601SHong Zhang /*MC 59f68a32c8SEmil Constantinescu TSRK3 - Third order RK scheme. 60f68a32c8SEmil Constantinescu 61f68a32c8SEmil Constantinescu This method has three stages. 62f68a32c8SEmil Constantinescu 63bcf0153eSBarry Smith Options Database Key: 6467b8a455SSatish Balay . -ts_rk_type 3 - use type 3 65287c2655SBarry Smith 66f68a32c8SEmil Constantinescu Level: advanced 67f68a32c8SEmil Constantinescu 681cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 69f68a32c8SEmil Constantinescu M*/ 70f68a32c8SEmil Constantinescu /*MC 711d27aa22SBarry Smith TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method <https://doi.org/10.1016/0893-9659(89)90079-7> 722109b73fSDebojyoti Ghosh 73d1905564SLisandro Dalcin This method has four stages with the First Same As Last (FSAL) property. 742109b73fSDebojyoti Ghosh 75bcf0153eSBarry Smith Options Database Key: 7667b8a455SSatish Balay . -ts_rk_type 3bs - use type 3bs 77287c2655SBarry Smith 782109b73fSDebojyoti Ghosh Level: advanced 792109b73fSDebojyoti Ghosh 801cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 812109b73fSDebojyoti Ghosh M*/ 822109b73fSDebojyoti Ghosh /*MC 83f68a32c8SEmil Constantinescu TSRK4 - Fourth order RK scheme. 84f68a32c8SEmil Constantinescu 852e698f8bSDebojyoti Ghosh This is the classical Runge-Kutta method with four stages. 86f68a32c8SEmil Constantinescu 87bcf0153eSBarry Smith Options Database Key: 8867b8a455SSatish Balay . -ts_rk_type 4 - use type 4 89287c2655SBarry Smith 90f68a32c8SEmil Constantinescu Level: advanced 91f68a32c8SEmil Constantinescu 921cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 93f68a32c8SEmil Constantinescu M*/ 94f68a32c8SEmil Constantinescu /*MC 952e698f8bSDebojyoti Ghosh TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method. 96f68a32c8SEmil Constantinescu 97f68a32c8SEmil Constantinescu This method has six stages. 98f68a32c8SEmil Constantinescu 99bcf0153eSBarry Smith Options Database Key: 10067b8a455SSatish Balay . -ts_rk_type 5f - use type 5f 101287c2655SBarry Smith 102f68a32c8SEmil Constantinescu Level: advanced 103f68a32c8SEmil Constantinescu 1041cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 105f68a32c8SEmil Constantinescu M*/ 1062109b73fSDebojyoti Ghosh /*MC 1071d27aa22SBarry Smith TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method <https://doi.org/10.1016/0771-050X(80)90013-3> 1082109b73fSDebojyoti Ghosh 109d1905564SLisandro Dalcin This method has seven stages with the First Same As Last (FSAL) property. 1102109b73fSDebojyoti Ghosh 111bcf0153eSBarry Smith Options Database Key: 11267b8a455SSatish Balay . -ts_rk_type 5dp - use type 5dp 113287c2655SBarry Smith 1142109b73fSDebojyoti Ghosh Level: advanced 1152109b73fSDebojyoti Ghosh 1161cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 1172109b73fSDebojyoti Ghosh M*/ 11805e23783SLisandro Dalcin /*MC 1191d27aa22SBarry Smith TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method <https://doi.org/10.1016/0898-1221(96)00141-1> 12005e23783SLisandro Dalcin 12105e23783SLisandro Dalcin This method has eight stages with the First Same As Last (FSAL) property. 12205e23783SLisandro Dalcin 123bcf0153eSBarry Smith Options Database Key: 12467b8a455SSatish Balay . -ts_rk_type 5bs - use type 5bs 125287c2655SBarry Smith 12605e23783SLisandro Dalcin Level: advanced 12705e23783SLisandro Dalcin 1281cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 12905e23783SLisandro Dalcin M*/ 13037acaa02SHendrik Ranocha /*MC 13137acaa02SHendrik Ranocha TSRK6VR - Sixth order robust Verner RK scheme with fifth order embedded method. 1321d27aa22SBarry Smith <http://people.math.sfu.ca/~jverner/RKV65.IIIXb.Robust.00010102836.081204.CoeffsOnlyRAT> 13337acaa02SHendrik Ranocha 13437acaa02SHendrik Ranocha This method has nine stages with the First Same As Last (FSAL) property. 13537acaa02SHendrik Ranocha 136bcf0153eSBarry Smith Options Database Key: 13767b8a455SSatish Balay . -ts_rk_type 6vr - use type 6vr 13837acaa02SHendrik Ranocha 13937acaa02SHendrik Ranocha Level: advanced 14037acaa02SHendrik Ranocha 1411cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 14237acaa02SHendrik Ranocha M*/ 14337acaa02SHendrik Ranocha /*MC 14437acaa02SHendrik Ranocha TSRK7VR - Seventh order robust Verner RK scheme with sixth order embedded method. 1451d27aa22SBarry Smith <http://people.math.sfu.ca/~jverner/RKV65.IIIXb.Robust.00010102836.081204.CoeffsOnlyRAT> 14637acaa02SHendrik Ranocha 1478f3d7ee2SCarsten Uphoff This method has ten stages. 14837acaa02SHendrik Ranocha 149bcf0153eSBarry Smith Options Database Key: 15067b8a455SSatish Balay . -ts_rk_type 7vr - use type 7vr 15137acaa02SHendrik Ranocha 15237acaa02SHendrik Ranocha Level: advanced 15337acaa02SHendrik Ranocha 1541cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 15537acaa02SHendrik Ranocha M*/ 15637acaa02SHendrik Ranocha /*MC 157d5b43468SJose E. Roman TSRK8VR - Eighth order robust Verner RK scheme with seventh order embedded method. 1581d27aa22SBarry Smith <http://people.math.sfu.ca/~jverner/RKV87.IIa.Robust.00000754677.081208.CoeffsOnlyRATandFLOAT> 15937acaa02SHendrik Ranocha 1608f3d7ee2SCarsten Uphoff This method has thirteen stages. 16137acaa02SHendrik Ranocha 162bcf0153eSBarry Smith Options Database Key: 16367b8a455SSatish Balay . -ts_rk_type 8vr - use type 8vr 16437acaa02SHendrik Ranocha 16537acaa02SHendrik Ranocha Level: advanced 16637acaa02SHendrik Ranocha 1671cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 16837acaa02SHendrik Ranocha M*/ 169262ff7bbSBarry Smith 170f68a32c8SEmil Constantinescu /*@C 171bcf0153eSBarry Smith TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in `TSRK` 172262ff7bbSBarry Smith 173f68a32c8SEmil Constantinescu Not Collective, but should be called by all processes which will need the schemes to be registered 174262ff7bbSBarry Smith 175f68a32c8SEmil Constantinescu Level: advanced 176262ff7bbSBarry Smith 1771cc06b55SBarry Smith .seealso: [](ch_ts), `TSRKRegisterDestroy()`, `TSRKRegister()` 178262ff7bbSBarry Smith @*/ 179d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKRegisterAll(void) 180d71ae5a4SJacob Faibussowitsch { 181262ff7bbSBarry Smith PetscFunctionBegin; 1823ba16761SJacob Faibussowitsch if (TSRKRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS); 183f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_TRUE; 184e4dd521cSBarry Smith 1854e82626cSLisandro Dalcin #define RC PetscRealConstant 186e4dd521cSBarry Smith { 187b16ce868SStefano Zampini const PetscReal A[1][1] = {{0}}; 188b16ce868SStefano Zampini const PetscReal b[1] = {RC(1.0)}; 1899566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK1FE, 1, 1, &A[0][0], b, NULL, NULL, 0, NULL)); 190e4dd521cSBarry Smith } 191db046809SBarry Smith { 192b16ce868SStefano Zampini const PetscReal A[2][2] = { 1939371c9d4SSatish Balay {0, 0}, 1949371c9d4SSatish Balay {RC(1.0), 0} 195b16ce868SStefano Zampini }; 196b16ce868SStefano Zampini const PetscReal b[2] = {RC(0.5), RC(0.5)}; 197b16ce868SStefano Zampini const PetscReal bembed[2] = {RC(1.0), 0}; 1989566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK2A, 2, 2, &A[0][0], b, NULL, bembed, 0, NULL)); 199db046809SBarry Smith } 200f68a32c8SEmil Constantinescu { 201b16ce868SStefano Zampini const PetscReal A[2][2] = { 2029371c9d4SSatish Balay {0, 0}, 2039371c9d4SSatish Balay {RC(0.5), 0} 204b16ce868SStefano Zampini }; 205b16ce868SStefano Zampini const PetscReal b[2] = {0, RC(1.0)}; 2069566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK2B, 2, 2, &A[0][0], b, NULL, NULL, 0, NULL)); 207e7685601SHong Zhang } 208e7685601SHong Zhang { 209b16ce868SStefano Zampini const PetscReal A[3][3] = { 2109371c9d4SSatish Balay {0, 0, 0}, 2114e82626cSLisandro Dalcin {RC(2.0) / RC(3.0), 0, 0}, 2129371c9d4SSatish Balay {RC(-1.0) / RC(3.0), RC(1.0), 0} 213b16ce868SStefano Zampini }; 214b16ce868SStefano Zampini const PetscReal b[3] = {RC(0.25), RC(0.5), RC(0.25)}; 2159566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK3, 3, 3, &A[0][0], b, NULL, NULL, 0, NULL)); 216fdefd5e5SDebojyoti Ghosh } 217fdefd5e5SDebojyoti Ghosh { 218b16ce868SStefano Zampini const PetscReal A[4][4] = { 2199371c9d4SSatish Balay {0, 0, 0, 0}, 2204e82626cSLisandro Dalcin {RC(1.0) / RC(2.0), 0, 0, 0}, 2214e82626cSLisandro Dalcin {0, RC(3.0) / RC(4.0), 0, 0}, 2229371c9d4SSatish Balay {RC(2.0) / RC(9.0), RC(1.0) / RC(3.0), RC(4.0) / RC(9.0), 0} 223b16ce868SStefano Zampini }; 224b16ce868SStefano Zampini const PetscReal b[4] = {RC(2.0) / RC(9.0), RC(1.0) / RC(3.0), RC(4.0) / RC(9.0), 0}; 225b16ce868SStefano Zampini const PetscReal bembed[4] = {RC(7.0) / RC(24.0), RC(1.0) / RC(4.0), RC(1.0) / RC(3.0), RC(1.0) / RC(8.0)}; 2269566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK3BS, 3, 4, &A[0][0], b, NULL, bembed, 0, NULL)); 227db046809SBarry Smith } 228f68a32c8SEmil Constantinescu { 229b16ce868SStefano Zampini const PetscReal A[4][4] = { 2309371c9d4SSatish Balay {0, 0, 0, 0}, 2314e82626cSLisandro Dalcin {RC(0.5), 0, 0, 0}, 2324e82626cSLisandro Dalcin {0, RC(0.5), 0, 0}, 2339371c9d4SSatish Balay {0, 0, RC(1.0), 0} 234b16ce868SStefano Zampini }; 235b16ce868SStefano Zampini const PetscReal b[4] = {RC(1.0) / RC(6.0), RC(1.0) / RC(3.0), RC(1.0) / RC(3.0), RC(1.0) / RC(6.0)}; 2369566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK4, 4, 4, &A[0][0], b, NULL, NULL, 0, NULL)); 237db046809SBarry Smith } 238f68a32c8SEmil Constantinescu { 239b16ce868SStefano Zampini const PetscReal A[6][6] = { 2409371c9d4SSatish Balay {0, 0, 0, 0, 0, 0}, 2414e82626cSLisandro Dalcin {RC(0.25), 0, 0, 0, 0, 0}, 2424e82626cSLisandro Dalcin {RC(3.0) / RC(32.0), RC(9.0) / RC(32.0), 0, 0, 0, 0}, 2434e82626cSLisandro Dalcin {RC(1932.0) / RC(2197.0), RC(-7200.0) / RC(2197.0), RC(7296.0) / RC(2197.0), 0, 0, 0}, 2444e82626cSLisandro Dalcin {RC(439.0) / RC(216.0), RC(-8.0), RC(3680.0) / RC(513.0), RC(-845.0) / RC(4104.0), 0, 0}, 2459371c9d4SSatish Balay {RC(-8.0) / RC(27.0), RC(2.0), RC(-3544.0) / RC(2565.0), RC(1859.0) / RC(4104.0), RC(-11.0) / RC(40.0), 0} 246b16ce868SStefano Zampini }; 247b16ce868SStefano Zampini const PetscReal b[6] = {RC(16.0) / RC(135.0), 0, RC(6656.0) / RC(12825.0), RC(28561.0) / RC(56430.0), RC(-9.0) / RC(50.0), RC(2.0) / RC(55.0)}; 248b16ce868SStefano Zampini const PetscReal bembed[6] = {RC(25.0) / RC(216.0), 0, RC(1408.0) / RC(2565.0), RC(2197.0) / RC(4104.0), RC(-1.0) / RC(5.0), 0}; 2499566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK5F, 5, 6, &A[0][0], b, NULL, bembed, 0, NULL)); 250fdefd5e5SDebojyoti Ghosh } 251fdefd5e5SDebojyoti Ghosh { 252b16ce868SStefano Zampini const PetscReal A[7][7] = { 2539371c9d4SSatish Balay {0, 0, 0, 0, 0, 0, 0}, 2544e82626cSLisandro Dalcin {RC(1.0) / RC(5.0), 0, 0, 0, 0, 0, 0}, 2554e82626cSLisandro Dalcin {RC(3.0) / RC(40.0), RC(9.0) / RC(40.0), 0, 0, 0, 0, 0}, 2564e82626cSLisandro Dalcin {RC(44.0) / RC(45.0), RC(-56.0) / RC(15.0), RC(32.0) / RC(9.0), 0, 0, 0, 0}, 2574e82626cSLisandro 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}, 2584e82626cSLisandro 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}, 2599371c9d4SSatish Balay {RC(35.0) / RC(384.0), 0, RC(500.0) / RC(1113.0), RC(125.0) / RC(192.0), RC(-2187.0) / RC(6784.0), RC(11.0) / RC(84.0), 0} 260b16ce868SStefano Zampini }; 261b16ce868SStefano Zampini const PetscReal b[7] = {RC(35.0) / RC(384.0), 0, RC(500.0) / RC(1113.0), RC(125.0) / RC(192.0), RC(-2187.0) / RC(6784.0), RC(11.0) / RC(84.0), 0}; 262b16ce868SStefano Zampini const PetscReal bembed[7] = {RC(5179.0) / RC(57600.0), 0, RC(7571.0) / RC(16695.0), RC(393.0) / RC(640.0), RC(-92097.0) / RC(339200.0), RC(187.0) / RC(2100.0), RC(1.0) / RC(40.0)}; 263b16ce868SStefano Zampini const PetscReal binterp[7][5] = { 264b16ce868SStefano Zampini {RC(1.0), RC(-4034104133.0) / RC(1410260304.0), RC(105330401.0) / RC(33982176.0), RC(-13107642775.0) / RC(11282082432.0), RC(6542295.0) / RC(470086768.0) }, 265b16ce868SStefano Zampini {0, 0, 0, 0, 0 }, 266b16ce868SStefano Zampini {0, RC(132343189600.0) / RC(32700410799.0), RC(-833316000.0) / RC(131326951.0), RC(91412856700.0) / RC(32700410799.0), RC(-523383600.0) / RC(10900136933.0) }, 267b16ce868SStefano Zampini {0, RC(-115792950.0) / RC(29380423.0), RC(185270875.0) / RC(16991088.0), RC(-12653452475.0) / RC(1880347072.0), RC(98134425.0) / RC(235043384.0) }, 268b16ce868SStefano Zampini {0, RC(70805911779.0) / RC(24914598704.0), RC(-4531260609.0) / RC(600351776.0), RC(988140236175.0) / RC(199316789632.0), RC(-14307999165.0) / RC(24914598704.0)}, 269b16ce868SStefano Zampini {0, RC(-331320693.0) / RC(205662961.0), RC(31361737.0) / RC(7433601.0), RC(-2426908385.0) / RC(822651844.0), RC(97305120.0) / RC(205662961.0) }, 270b16ce868SStefano Zampini {0, RC(44764047.0) / RC(29380423.0), RC(-1532549.0) / RC(353981.0), RC(90730570.0) / RC(29380423.0), RC(-8293050.0) / RC(29380423.0) } 271b16ce868SStefano Zampini }; 2729566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK5DP, 5, 7, &A[0][0], b, NULL, bembed, 5, binterp[0])); 273f68a32c8SEmil Constantinescu } 27405e23783SLisandro Dalcin { 275b16ce868SStefano Zampini const PetscReal A[8][8] = { 2769371c9d4SSatish Balay {0, 0, 0, 0, 0, 0, 0, 0}, 2774e82626cSLisandro Dalcin {RC(1.0) / RC(6.0), 0, 0, 0, 0, 0, 0, 0}, 2784e82626cSLisandro Dalcin {RC(2.0) / RC(27.0), RC(4.0) / RC(27.0), 0, 0, 0, 0, 0, 0}, 2794e82626cSLisandro Dalcin {RC(183.0) / RC(1372.0), RC(-162.0) / RC(343.0), RC(1053.0) / RC(1372.0), 0, 0, 0, 0, 0}, 2804e82626cSLisandro 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}, 2814e82626cSLisandro 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}, 2824e82626cSLisandro 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}, 2839371c9d4SSatish Balay {RC(587.0) / RC(8064.0), 0, RC(4440339.0) / RC(15491840.0), RC(24353.0) / RC(124800.0), RC(387.0) / RC(44800.0), RC(2152.0) / RC(5985.0), RC(7267.0) / RC(94080.0), 0} 284b16ce868SStefano Zampini }; 285b16ce868SStefano Zampini const PetscReal b[8] = {RC(587.0) / RC(8064.0), 0, RC(4440339.0) / RC(15491840.0), RC(24353.0) / RC(124800.0), RC(387.0) / RC(44800.0), RC(2152.0) / RC(5985.0), RC(7267.0) / RC(94080.0), 0}; 286b16ce868SStefano Zampini const PetscReal bembed[8] = {RC(2479.0) / RC(34992.0), 0, RC(123.0) / RC(416.0), RC(612941.0) / RC(3411720.0), RC(43.0) / RC(1440.0), RC(2272.0) / RC(6561.0), RC(79937.0) / RC(1113912.0), RC(3293.0) / RC(556956.0)}; 2879566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK5BS, 5, 8, &A[0][0], b, NULL, bembed, 0, NULL)); 28805e23783SLisandro Dalcin } 28937acaa02SHendrik Ranocha { 290b16ce868SStefano Zampini const PetscReal A[9][9] = { 2919371c9d4SSatish Balay {0, 0, 0, 0, 0, 0, 0, 0, 0}, 29263fe90e8SHendrik Ranocha {RC(1.8000000000000000000000000000000000000000e-01), 0, 0, 0, 0, 0, 0, 0, 0}, 29363fe90e8SHendrik Ranocha {RC(8.9506172839506172839506172839506172839506e-02), RC(7.7160493827160493827160493827160493827160e-02), 0, 0, 0, 0, 0, 0, 0}, 29463fe90e8SHendrik Ranocha {RC(6.2500000000000000000000000000000000000000e-02), 0, RC(1.8750000000000000000000000000000000000000e-01), 0, 0, 0, 0, 0, 0}, 29563fe90e8SHendrik Ranocha {RC(3.1651600000000000000000000000000000000000e-01), 0, RC(-1.0449480000000000000000000000000000000000e+00), RC(1.2584320000000000000000000000000000000000e+00), 0, 0, 0, 0, 0}, 29663fe90e8SHendrik Ranocha {RC(2.7232612736485626257225065566674305502508e-01), 0, RC(-8.2513360323886639676113360323886639676113e-01), RC(1.0480917678812415654520917678812415654521e+00), RC(1.0471570799276856873679117969088177628396e-01), 0, 0, 0, 0}, 29763fe90e8SHendrik 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}, 29863fe90e8SHendrik 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}, 2999371c9d4SSatish Balay {RC(7.6388888888888888888888888888888888888889e-02), 0, 0, RC(3.6940836940836940836940836940836940836941e-01), 0, RC(2.4801587301587301587301587301587301587302e-01), RC(2.3674242424242424242424242424242424242424e-01), RC(6.9444444444444444444444444444444444444444e-02), 0} 300b16ce868SStefano Zampini }; 301b16ce868SStefano Zampini const PetscReal b[9] = {RC(7.6388888888888888888888888888888888888889e-02), 0, 0, RC(3.6940836940836940836940836940836940836941e-01), 0, RC(2.4801587301587301587301587301587301587302e-01), RC(2.3674242424242424242424242424242424242424e-01), 302b16ce868SStefano Zampini RC(6.9444444444444444444444444444444444444444e-02), 0}; 303b16ce868SStefano Zampini const PetscReal bembed[9] = {RC(5.8700209643605870020964360587002096436059e-02), 0, 0, RC(4.8072562358276643990929705215419501133787e-01), RC(-8.5341242076919085578832094861228313083563e-01), RC(1.2046485260770975056689342403628117913832e+00), 0, RC(-5.9242373072160306202859394348756050883710e-02), RC(1.6858043453788134639198468985703028256220e-01)}; 3049566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK6VR, 6, 9, &A[0][0], b, NULL, bembed, 0, NULL)); 30537acaa02SHendrik Ranocha } 30637acaa02SHendrik Ranocha { 307b16ce868SStefano Zampini const PetscReal A[10][10] = { 3089371c9d4SSatish Balay {0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, 30963fe90e8SHendrik Ranocha {RC(5.0000000000000000000000000000000000000000e-03), 0, 0, 0, 0, 0, 0, 0, 0, 0}, 31063fe90e8SHendrik Ranocha {RC(-1.0767901234567901234567901234567901234568e+00), RC(1.1856790123456790123456790123456790123457e+00), 0, 0, 0, 0, 0, 0, 0, 0}, 31163fe90e8SHendrik Ranocha {RC(4.0833333333333333333333333333333333333333e-02), 0, RC(1.2250000000000000000000000000000000000000e-01), 0, 0, 0, 0, 0, 0, 0}, 31263fe90e8SHendrik Ranocha {RC(6.3607142857142857142857142857142857142857e-01), 0, RC(-2.4444642857142857142857142857142857142857e+00), RC(2.2633928571428571428571428571428571428571e+00), 0, 0, 0, 0, 0, 0}, 31363fe90e8SHendrik Ranocha {RC(-2.5351211079349245229256383554660215487207e+00), 0, RC(1.0299374654449267920438514460756024913612e+01), RC(-7.9513032885990579949493217458266876536482e+00), RC(7.9301148923100592201226014271115261823800e-01), 0, 0, 0, 0, 0}, 31463fe90e8SHendrik 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}, 31563fe90e8SHendrik 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}, 31663fe90e8SHendrik 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}, 3179371c9d4SSatish Balay {RC(-1.4449518916777735137351003179355712360517e+00), 0, RC(8.0318913859955919224117033223019560435041e+00), RC(-7.5831741663401346820798883023671588604984e+00), RC(3.5816169353190074211247685442452878696855e+00), RC(-2.4369722632199529411183809065693752383733e+00), RC(8.5158999992326179339689766032486142173390e-01), 0, 0, 0} 318b16ce868SStefano Zampini }; 319b16ce868SStefano Zampini const PetscReal b[10] = {RC(4.7425837833706756083569172717574534698932e-02), 0, 0, RC(2.5622361659370562659961727458274623448160e-01), RC(2.6951376833074206619473817258075952886764e-01), RC(1.2686622409092782845989138364739173247882e-01), RC(2.4887225942060071622046449427647492767292e-01), RC(3.0744837408200631335304388479099184768645e-03), RC(4.8023809989496943308189063347143123323209e-02), 0}; 320b16ce868SStefano Zampini const PetscReal bembed[10] = {RC(4.7485247699299631037531273805727961552268e-02), 0, 0, RC(2.5599412588690633297154918245905393870497e-01), RC(2.7058478081067688722530891099268135732387e-01), RC(1.2505618684425992913638822323746917920448e-01), 3219371c9d4SSatish Balay RC(2.5204468723743860507184043820197442562182e-01), 0, 0, RC(4.8834971521418614557381971303093137592592e-02)}; 3229566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK7VR, 7, 10, &A[0][0], b, NULL, bembed, 0, NULL)); 32337acaa02SHendrik Ranocha } 32437acaa02SHendrik Ranocha { 325b16ce868SStefano Zampini const PetscReal A[13][13] = { 3269371c9d4SSatish Balay {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, 32763fe90e8SHendrik Ranocha {RC(2.5000000000000000000000000000000000000000e-01), 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, 32863fe90e8SHendrik Ranocha {RC(8.7400846504915232052686327594877411977046e-02), RC(2.5487604938654321753087950620345685135815e-02), 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, 32963fe90e8SHendrik Ranocha {RC(4.2333169291338582677165354330708661417323e-02), 0, RC(1.2699950787401574803149606299212598425197e-01), 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, 33063fe90e8SHendrik Ranocha {RC(4.2609505888742261494881445237572274090942e-01), 0, RC(-1.5987952846591523265427733230657181117089e+00), RC(1.5967002257717297115939588706899953707994e+00), 0, 0, 0, 0, 0, 0, 0, 0, 0}, 33163fe90e8SHendrik Ranocha {RC(5.0719337296713929515090618138513639239329e-02), 0, 0, RC(2.5433377264600407582754714408877778031369e-01), RC(2.0394689005728199465736223777270858044698e-01), 0, 0, 0, 0, 0, 0, 0, 0}, 33263fe90e8SHendrik 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}, 33363fe90e8SHendrik 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}, 33463fe90e8SHendrik 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}, 33563fe90e8SHendrik 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}, 33663fe90e8SHendrik 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}, 33763fe90e8SHendrik 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}, 3389371c9d4SSatish Balay {RC(8.8735762208534719663211694051981022704884e-01), 0, 0, RC(-2.9754597821085367558513632804709301581977e+00), RC(5.6007170094881630597990392548350098923829e+00), RC(-5.9156074505366744680014930189941657351840e+00), RC(2.2029689156134927016879142540807638331238e-01), RC(1.0155097824462216666143271340902996997549e-01), RC(1.1514345647386055909780397752125850553556e+00), RC(1.9297101665271239396134361900805843653065e+00), 0, 0, 0} 339b16ce868SStefano Zampini }; 340b16ce868SStefano Zampini const PetscReal b[13] = {RC(4.4729564666695714203015840429049382466467e-02), 0, 0, 0, 0, RC(1.5691033527708199813368698010726645409175e-01), RC(1.8460973408151637740702451873526277892035e-01), RC(2.2516380602086991042479419400350721970920e-01), RC(1.4794615651970234687005179885449141753736e-01), RC(7.6055542444955825269798361910336491012732e-02), RC(1.2277290235018619610824346315921437388535e-01), RC(4.1811958638991631583384842800871882376786e-02), 0}; 341b16ce868SStefano Zampini const PetscReal bembed[13] = {RC(4.5847111400495925878664730122010282095875e-02), 0, 0, 0, 0, RC(2.6231891404152387437443356584845803392392e-01), RC(1.9169372337852611904485738635688429008025e-01), RC(2.1709172327902618330978407422906448568196e-01), RC(1.2738189624833706796803169450656737867900e-01), RC(1.1510530385365326258240515750043192148894e-01), 0, 0, RC(4.0561327798437566841823391436583608050053e-02)}; 3429566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK8VR, 8, 13, &A[0][0], b, NULL, bembed, 0, NULL)); 34337acaa02SHendrik Ranocha } 3444e82626cSLisandro Dalcin #undef RC 3453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 346db046809SBarry Smith } 347db046809SBarry Smith 348f68a32c8SEmil Constantinescu /*@C 349bcf0153eSBarry Smith TSRKRegisterDestroy - Frees the list of schemes that were registered by `TSRKRegister()`. 350f68a32c8SEmil Constantinescu 351f68a32c8SEmil Constantinescu Not Collective 352f68a32c8SEmil Constantinescu 353f68a32c8SEmil Constantinescu Level: advanced 354f68a32c8SEmil Constantinescu 3551cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKRegister()`, `TSRKRegisterAll()` 356f68a32c8SEmil Constantinescu @*/ 357d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKRegisterDestroy(void) 358d71ae5a4SJacob Faibussowitsch { 359f68a32c8SEmil Constantinescu RKTableauLink link; 360f68a32c8SEmil Constantinescu 361f68a32c8SEmil Constantinescu PetscFunctionBegin; 362f68a32c8SEmil Constantinescu while ((link = RKTableauList)) { 363f68a32c8SEmil Constantinescu RKTableau t = &link->tab; 364f68a32c8SEmil Constantinescu RKTableauList = link->next; 3659566063dSJacob Faibussowitsch PetscCall(PetscFree3(t->A, t->b, t->c)); 3669566063dSJacob Faibussowitsch PetscCall(PetscFree(t->bembed)); 3679566063dSJacob Faibussowitsch PetscCall(PetscFree(t->binterp)); 3689566063dSJacob Faibussowitsch PetscCall(PetscFree(t->name)); 3699566063dSJacob Faibussowitsch PetscCall(PetscFree(link)); 370f68a32c8SEmil Constantinescu } 371f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_FALSE; 3723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 373f68a32c8SEmil Constantinescu } 374f68a32c8SEmil Constantinescu 375f68a32c8SEmil Constantinescu /*@C 376bcf0153eSBarry Smith TSRKInitializePackage - This function initializes everything in the `TSRK` package. It is called 377bcf0153eSBarry Smith from `TSInitializePackage()`. 378f68a32c8SEmil Constantinescu 379f68a32c8SEmil Constantinescu Level: developer 380f68a32c8SEmil Constantinescu 3811cc06b55SBarry Smith .seealso: [](ch_ts), `TSInitializePackage()`, `PetscInitialize()`, `TSRKFinalizePackage()` 382f68a32c8SEmil Constantinescu @*/ 383d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKInitializePackage(void) 384d71ae5a4SJacob Faibussowitsch { 385e4dd521cSBarry Smith PetscFunctionBegin; 3863ba16761SJacob Faibussowitsch if (TSRKPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS); 387f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_TRUE; 3889566063dSJacob Faibussowitsch PetscCall(TSRKRegisterAll()); 3899566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSRKFinalizePackage)); 3903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 391f68a32c8SEmil Constantinescu } 392186e87acSLisandro Dalcin 393f68a32c8SEmil Constantinescu /*@C 394bcf0153eSBarry Smith TSRKFinalizePackage - This function destroys everything in the `TSRK` package. It is 395bcf0153eSBarry Smith called from `PetscFinalize()`. 396186e87acSLisandro Dalcin 397f68a32c8SEmil Constantinescu Level: developer 398f68a32c8SEmil Constantinescu 3991cc06b55SBarry Smith .seealso: [](ch_ts), `PetscFinalize()`, `TSRKInitializePackage()` 400f68a32c8SEmil Constantinescu @*/ 401d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKFinalizePackage(void) 402d71ae5a4SJacob Faibussowitsch { 403f68a32c8SEmil Constantinescu PetscFunctionBegin; 404f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_FALSE; 4059566063dSJacob Faibussowitsch PetscCall(TSRKRegisterDestroy()); 4063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 407f68a32c8SEmil Constantinescu } 408f68a32c8SEmil Constantinescu 409f68a32c8SEmil Constantinescu /*@C 410bcf0153eSBarry Smith TSRKRegister - register an `TSRK` scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 411f68a32c8SEmil Constantinescu 412*cc4c1da9SBarry Smith Not Collective, but the same schemes should be registered on all processes on which they will be used, No Fortran Support 413f68a32c8SEmil Constantinescu 414f68a32c8SEmil Constantinescu Input Parameters: 415f68a32c8SEmil Constantinescu + name - identifier for method 416f68a32c8SEmil Constantinescu . order - approximation order of method 417f68a32c8SEmil Constantinescu . s - number of stages, this is the dimension of the matrices below 418f68a32c8SEmil Constantinescu . A - stage coefficients (dimension s*s, row-major) 419f68a32c8SEmil Constantinescu . b - step completion table (dimension s; NULL to use last row of A) 420f68a32c8SEmil Constantinescu . c - abscissa (dimension s; NULL to use row sums of A) 421f68a32c8SEmil Constantinescu . bembed - completion table for embedded method (dimension s; NULL if not available) 4223a8a9803SLisandro Dalcin . p - Order of the interpolation scheme, equal to the number of columns of binterp 4233a8a9803SLisandro Dalcin - binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1) 424f68a32c8SEmil Constantinescu 425f68a32c8SEmil Constantinescu Level: advanced 426f68a32c8SEmil Constantinescu 427bcf0153eSBarry Smith Note: 428bcf0153eSBarry Smith Several `TSRK` methods are provided, this function is only needed to create new methods. 429bcf0153eSBarry Smith 4301cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK` 431f68a32c8SEmil Constantinescu @*/ 432d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKRegister(TSRKType name, PetscInt order, PetscInt s, const PetscReal A[], const PetscReal b[], const PetscReal c[], const PetscReal bembed[], PetscInt p, const PetscReal binterp[]) 433d71ae5a4SJacob Faibussowitsch { 434f68a32c8SEmil Constantinescu RKTableauLink link; 435f68a32c8SEmil Constantinescu RKTableau t; 436f68a32c8SEmil Constantinescu PetscInt i, j; 437f68a32c8SEmil Constantinescu 438f68a32c8SEmil Constantinescu PetscFunctionBegin; 4394f572ea9SToby Isaac PetscAssertPointer(name, 1); 4404f572ea9SToby Isaac PetscAssertPointer(A, 4); 4414f572ea9SToby Isaac if (b) PetscAssertPointer(b, 5); 4424f572ea9SToby Isaac if (c) PetscAssertPointer(c, 6); 4434f572ea9SToby Isaac if (bembed) PetscAssertPointer(bembed, 7); 4444f572ea9SToby Isaac if (binterp || p > 1) PetscAssertPointer(binterp, 9); 445095c51faSJed Brown PetscCheck(s >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected number of stages s %" PetscInt_FMT " >= 0", s); 4463a8a9803SLisandro Dalcin 4479566063dSJacob Faibussowitsch PetscCall(TSRKInitializePackage()); 4489566063dSJacob Faibussowitsch PetscCall(PetscNew(&link)); 449f68a32c8SEmil Constantinescu t = &link->tab; 4503a8a9803SLisandro Dalcin 4519566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &t->name)); 452f68a32c8SEmil Constantinescu t->order = order; 453f68a32c8SEmil Constantinescu t->s = s; 4549566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(s * s, &t->A, s, &t->b, s, &t->c)); 4559566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->A, A, s * s)); 4569566063dSJacob Faibussowitsch if (b) PetscCall(PetscArraycpy(t->b, b, s)); 4579371c9d4SSatish Balay else 4589371c9d4SSatish Balay for (i = 0; i < s; i++) t->b[i] = A[(s - 1) * s + i]; 4599566063dSJacob Faibussowitsch if (c) PetscCall(PetscArraycpy(t->c, c, s)); 4609371c9d4SSatish Balay else 4619371c9d4SSatish Balay for (i = 0; i < s; i++) 4629371c9d4SSatish Balay for (j = 0, t->c[i] = 0; j < s; j++) t->c[i] += A[i * s + j]; 463d760c35bSDebojyoti Ghosh t->FSAL = PETSC_TRUE; 4649371c9d4SSatish Balay for (i = 0; i < s; i++) 4659371c9d4SSatish Balay if (t->A[(s - 1) * s + i] != t->b[i]) t->FSAL = PETSC_FALSE; 4663a8a9803SLisandro Dalcin 467f68a32c8SEmil Constantinescu if (bembed) { 4689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s, &t->bembed)); 4699566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->bembed, bembed, s)); 470f68a32c8SEmil Constantinescu } 471f68a32c8SEmil Constantinescu 4729371c9d4SSatish Balay if (!binterp) { 4739371c9d4SSatish Balay p = 1; 4749371c9d4SSatish Balay binterp = t->b; 4759371c9d4SSatish Balay } 4763a8a9803SLisandro Dalcin t->p = p; 4779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s * p, &t->binterp)); 4789566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->binterp, binterp, s * p)); 4793a8a9803SLisandro Dalcin 480f68a32c8SEmil Constantinescu link->next = RKTableauList; 481f68a32c8SEmil Constantinescu RKTableauList = link; 4823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 483f68a32c8SEmil Constantinescu } 484f68a32c8SEmil Constantinescu 48566976f2fSJacob Faibussowitsch static PetscErrorCode TSRKGetTableau_RK(TS ts, PetscInt *s, const PetscReal **A, const PetscReal **b, const PetscReal **c, const PetscReal **bembed, PetscInt *p, const PetscReal **binterp, PetscBool *FSAL) 486d71ae5a4SJacob Faibussowitsch { 4870f7a28cdSStefano Zampini TS_RK *rk = (TS_RK *)ts->data; 4880f7a28cdSStefano Zampini RKTableau tab = rk->tableau; 4890f7a28cdSStefano Zampini 4900f7a28cdSStefano Zampini PetscFunctionBegin; 4910f7a28cdSStefano Zampini if (s) *s = tab->s; 4920f7a28cdSStefano Zampini if (A) *A = tab->A; 4930f7a28cdSStefano Zampini if (b) *b = tab->b; 4940f7a28cdSStefano Zampini if (c) *c = tab->c; 4950f7a28cdSStefano Zampini if (bembed) *bembed = tab->bembed; 4960f7a28cdSStefano Zampini if (p) *p = tab->p; 4970f7a28cdSStefano Zampini if (binterp) *binterp = tab->binterp; 4980f7a28cdSStefano Zampini if (FSAL) *FSAL = tab->FSAL; 4993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5000f7a28cdSStefano Zampini } 5010f7a28cdSStefano Zampini 5020f7a28cdSStefano Zampini /*@C 503bcf0153eSBarry Smith TSRKGetTableau - Get info on the `TSRK` tableau 5040f7a28cdSStefano Zampini 5050f7a28cdSStefano Zampini Not Collective 5060f7a28cdSStefano Zampini 507f899ff85SJose E. Roman Input Parameter: 5080f7a28cdSStefano Zampini . ts - timestepping context 5090f7a28cdSStefano Zampini 5100f7a28cdSStefano Zampini Output Parameters: 5110f7a28cdSStefano Zampini + s - number of stages, this is the dimension of the matrices below 5120f7a28cdSStefano Zampini . A - stage coefficients (dimension s*s, row-major) 5130f7a28cdSStefano Zampini . b - step completion table (dimension s) 5140f7a28cdSStefano Zampini . c - abscissa (dimension s) 5150f7a28cdSStefano Zampini . bembed - completion table for embedded method (dimension s; NULL if not available) 5160f7a28cdSStefano Zampini . p - Order of the interpolation scheme, equal to the number of columns of binterp 5170f7a28cdSStefano Zampini . binterp - Coefficients of the interpolation formula (dimension s*p) 518aaa8cc7dSPierre Jolivet - FSAL - whether or not the scheme has the First Same As Last property 5190f7a28cdSStefano Zampini 5200f7a28cdSStefano Zampini Level: developer 5210f7a28cdSStefano Zampini 5221cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKRegister()`, `TSRKSetType()` 5230f7a28cdSStefano Zampini @*/ 524d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKGetTableau(TS ts, PetscInt *s, const PetscReal **A, const PetscReal **b, const PetscReal **c, const PetscReal **bembed, PetscInt *p, const PetscReal **binterp, PetscBool *FSAL) 525d71ae5a4SJacob Faibussowitsch { 5260f7a28cdSStefano Zampini PetscFunctionBegin; 5270f7a28cdSStefano Zampini PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 5289371c9d4SSatish Balay PetscUseMethod(ts, "TSRKGetTableau_C", (TS, PetscInt *, const PetscReal **, const PetscReal **, const PetscReal **, const PetscReal **, PetscInt *, const PetscReal **, PetscBool *), (ts, s, A, b, c, bembed, p, binterp, FSAL)); 5293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5300f7a28cdSStefano Zampini } 5310f7a28cdSStefano Zampini 532e4dd521cSBarry Smith /* 5332c9cb062Sluzhanghpp This is for single-step RK method 534f68a32c8SEmil Constantinescu The step completion formula is 535e4dd521cSBarry Smith 536f68a32c8SEmil Constantinescu x1 = x0 + h b^T YdotRHS 537f68a32c8SEmil Constantinescu 538f68a32c8SEmil Constantinescu This function can be called before or after ts->vec_sol has been updated. 539f68a32c8SEmil Constantinescu Suppose we have a completion formula (b) and an embedded formula (be) of different order. 540f68a32c8SEmil Constantinescu We can write 541f68a32c8SEmil Constantinescu 542f68a32c8SEmil Constantinescu x1e = x0 + h be^T YdotRHS 543f68a32c8SEmil Constantinescu = x1 - h b^T YdotRHS + h be^T YdotRHS 544f68a32c8SEmil Constantinescu = x1 + h (be - b)^T YdotRHS 545f68a32c8SEmil Constantinescu 546f68a32c8SEmil Constantinescu so we can evaluate the method with different order even after the step has been optimistically completed. 547e4dd521cSBarry Smith */ 548d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateStep_RK(TS ts, PetscInt order, Vec X, PetscBool *done) 549d71ae5a4SJacob Faibussowitsch { 550f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK *)ts->data; 551f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 552f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 553f68a32c8SEmil Constantinescu PetscReal h; 554f68a32c8SEmil Constantinescu PetscInt s = tab->s, j; 555f68a32c8SEmil Constantinescu 556f68a32c8SEmil Constantinescu PetscFunctionBegin; 557f68a32c8SEmil Constantinescu switch (rk->status) { 558f68a32c8SEmil Constantinescu case TS_STEP_INCOMPLETE: 559d71ae5a4SJacob Faibussowitsch case TS_STEP_PENDING: 560d71ae5a4SJacob Faibussowitsch h = ts->time_step; 561d71ae5a4SJacob Faibussowitsch break; 562d71ae5a4SJacob Faibussowitsch case TS_STEP_COMPLETE: 563d71ae5a4SJacob Faibussowitsch h = ts->ptime - ts->ptime_prev; 564d71ae5a4SJacob Faibussowitsch break; 565d71ae5a4SJacob Faibussowitsch default: 566d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus"); 567f68a32c8SEmil Constantinescu } 568f68a32c8SEmil Constantinescu if (order == tab->order) { 569f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { 5709566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X)); 57190b67129SHong Zhang for (j = 0; j < s; j++) w[j] = h * tab->b[j] / rk->dtratio; 5729566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, rk->YdotRHS)); 5739566063dSJacob Faibussowitsch } else PetscCall(VecCopy(ts->vec_sol, X)); 5743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 575f68a32c8SEmil Constantinescu } else if (order == tab->order - 1) { 576f68a32c8SEmil Constantinescu if (!tab->bembed) goto unavailable; 577f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { /*Complete with the embedded method (be)*/ 5789566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X)); 579f68a32c8SEmil Constantinescu for (j = 0; j < s; j++) w[j] = h * tab->bembed[j]; 5809566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, rk->YdotRHS)); 581f68a32c8SEmil Constantinescu } else { /*Rollback and re-complete using (be-b) */ 5829566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X)); 583f68a32c8SEmil Constantinescu for (j = 0; j < s; j++) w[j] = h * (tab->bembed[j] - tab->b[j]); 5849566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, rk->YdotRHS)); 585f68a32c8SEmil Constantinescu } 586f68a32c8SEmil Constantinescu if (done) *done = PETSC_TRUE; 5873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 588f68a32c8SEmil Constantinescu } 589f68a32c8SEmil Constantinescu unavailable: 590f68a32c8SEmil Constantinescu if (done) *done = PETSC_FALSE; 5919371c9d4SSatish Balay else 5929371c9d4SSatish Balay SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "RK '%s' of order %" PetscInt_FMT " cannot evaluate step at order %" PetscInt_FMT ". Consider using -ts_adapt_type none or a different method that has an embedded estimate.", tab->name, tab->order, order); 5933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 594f68a32c8SEmil Constantinescu } 595f68a32c8SEmil Constantinescu 596d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardCostIntegral_RK(TS ts) 597d71ae5a4SJacob Faibussowitsch { 5980f7a1166SHong Zhang TS_RK *rk = (TS_RK *)ts->data; 599cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 6000f7a1166SHong Zhang RKTableau tab = rk->tableau; 6010f7a1166SHong Zhang const PetscInt s = tab->s; 6020f7a1166SHong Zhang const PetscReal *b = tab->b, *c = tab->c; 6030f7a1166SHong Zhang Vec *Y = rk->Y; 6040f7a1166SHong Zhang PetscInt i; 6050f7a1166SHong Zhang 6060f7a1166SHong Zhang PetscFunctionBegin; 607cd4cee2dSHong Zhang /* No need to backup quadts->vec_sol since it can be reverted in TSRollBack_RK */ 6080f7a1166SHong Zhang for (i = s - 1; i >= 0; i--) { 609cd4cee2dSHong Zhang /* Evolve quadrature TS solution to compute integrals */ 6109566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(quadts, rk->ptime + rk->time_step * c[i], Y[i], ts->vec_costintegrand)); 6119566063dSJacob Faibussowitsch PetscCall(VecAXPY(quadts->vec_sol, rk->time_step * b[i], ts->vec_costintegrand)); 6120f7a1166SHong Zhang } 6133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6140f7a1166SHong Zhang } 6150f7a1166SHong Zhang 616d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointCostIntegral_RK(TS ts) 617d71ae5a4SJacob Faibussowitsch { 6180f7a1166SHong Zhang TS_RK *rk = (TS_RK *)ts->data; 6190f7a1166SHong Zhang RKTableau tab = rk->tableau; 620cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 6210f7a1166SHong Zhang const PetscInt s = tab->s; 6220f7a1166SHong Zhang const PetscReal *b = tab->b, *c = tab->c; 6230f7a1166SHong Zhang Vec *Y = rk->Y; 6240f7a1166SHong Zhang PetscInt i; 6250f7a1166SHong Zhang 6260f7a1166SHong Zhang PetscFunctionBegin; 6270f7a1166SHong Zhang for (i = s - 1; i >= 0; i--) { 628cd4cee2dSHong Zhang /* Evolve quadrature TS solution to compute integrals */ 6299566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(quadts, ts->ptime + ts->time_step * (1.0 - c[i]), Y[i], ts->vec_costintegrand)); 6309566063dSJacob Faibussowitsch PetscCall(VecAXPY(quadts->vec_sol, -ts->time_step * b[i], ts->vec_costintegrand)); 6310f7a1166SHong Zhang } 6323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6330f7a1166SHong Zhang } 6340f7a1166SHong Zhang 635d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRollBack_RK(TS ts) 636d71ae5a4SJacob Faibussowitsch { 637fecfb714SLisandro Dalcin TS_RK *rk = (TS_RK *)ts->data; 638cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 639fecfb714SLisandro Dalcin RKTableau tab = rk->tableau; 640fecfb714SLisandro Dalcin const PetscInt s = tab->s; 641cd4cee2dSHong Zhang const PetscReal *b = tab->b, *c = tab->c; 642fecfb714SLisandro Dalcin PetscScalar *w = rk->work; 643cd4cee2dSHong Zhang Vec *Y = rk->Y, *YdotRHS = rk->YdotRHS; 644fecfb714SLisandro Dalcin PetscInt j; 645fecfb714SLisandro Dalcin PetscReal h; 646fecfb714SLisandro Dalcin 647fecfb714SLisandro Dalcin PetscFunctionBegin; 648fecfb714SLisandro Dalcin switch (rk->status) { 649fecfb714SLisandro Dalcin case TS_STEP_INCOMPLETE: 650d71ae5a4SJacob Faibussowitsch case TS_STEP_PENDING: 651d71ae5a4SJacob Faibussowitsch h = ts->time_step; 652d71ae5a4SJacob Faibussowitsch break; 653d71ae5a4SJacob Faibussowitsch case TS_STEP_COMPLETE: 654d71ae5a4SJacob Faibussowitsch h = ts->ptime - ts->ptime_prev; 655d71ae5a4SJacob Faibussowitsch break; 656d71ae5a4SJacob Faibussowitsch default: 657d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus"); 658fecfb714SLisandro Dalcin } 659fecfb714SLisandro Dalcin for (j = 0; j < s; j++) w[j] = -h * b[j]; 6609566063dSJacob Faibussowitsch PetscCall(VecMAXPY(ts->vec_sol, s, w, YdotRHS)); 661cd4cee2dSHong Zhang if (quadts && ts->costintegralfwd) { 662cd4cee2dSHong Zhang for (j = 0; j < s; j++) { 663cd4cee2dSHong Zhang /* Revert the quadrature TS solution */ 6649566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(quadts, rk->ptime + h * c[j], Y[j], ts->vec_costintegrand)); 6659566063dSJacob Faibussowitsch PetscCall(VecAXPY(quadts->vec_sol, -h * b[j], ts->vec_costintegrand)); 666cd4cee2dSHong Zhang } 667cd4cee2dSHong Zhang } 6683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 669fecfb714SLisandro Dalcin } 670fecfb714SLisandro Dalcin 671d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardStep_RK(TS ts) 672d71ae5a4SJacob Faibussowitsch { 673922a638cSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 674922a638cSHong Zhang RKTableau tab = rk->tableau; 675922a638cSHong Zhang Mat J, *MatsFwdSensipTemp = rk->MatsFwdSensipTemp; 676922a638cSHong Zhang const PetscInt s = tab->s; 677922a638cSHong Zhang const PetscReal *A = tab->A, *c = tab->c, *b = tab->b; 678922a638cSHong Zhang Vec *Y = rk->Y; 679922a638cSHong Zhang PetscInt i, j; 680922a638cSHong Zhang PetscReal stage_time, h = ts->time_step; 681922a638cSHong Zhang PetscBool zero; 682922a638cSHong Zhang 683922a638cSHong Zhang PetscFunctionBegin; 6849566063dSJacob Faibussowitsch PetscCall(MatCopy(ts->mat_sensip, rk->MatFwdSensip0, SAME_NONZERO_PATTERN)); 6859566063dSJacob Faibussowitsch PetscCall(TSGetRHSJacobian(ts, &J, NULL, NULL, NULL)); 686922a638cSHong Zhang 687922a638cSHong Zhang for (i = 0; i < s; i++) { 688922a638cSHong Zhang stage_time = ts->ptime + h * c[i]; 689922a638cSHong Zhang zero = PETSC_FALSE; 690922a638cSHong Zhang if (b[i] == 0 && i == s - 1) zero = PETSC_TRUE; 691922a638cSHong Zhang /* TLM Stage values */ 692922a638cSHong Zhang if (!i) { 6939566063dSJacob Faibussowitsch PetscCall(MatCopy(ts->mat_sensip, rk->MatsFwdStageSensip[i], SAME_NONZERO_PATTERN)); 694922a638cSHong Zhang } else if (!zero) { 6959566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(rk->MatsFwdStageSensip[i])); 69648a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(MatAXPY(rk->MatsFwdStageSensip[i], h * A[i * s + j], MatsFwdSensipTemp[j], SAME_NONZERO_PATTERN)); 6979566063dSJacob Faibussowitsch PetscCall(MatAXPY(rk->MatsFwdStageSensip[i], 1., ts->mat_sensip, SAME_NONZERO_PATTERN)); 698922a638cSHong Zhang } else { 6999566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(rk->MatsFwdStageSensip[i])); 700922a638cSHong Zhang } 701922a638cSHong Zhang 7029566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobian(ts, stage_time, Y[i], J, J)); 7039566063dSJacob Faibussowitsch PetscCall(MatMatMult(J, rk->MatsFwdStageSensip[i], MAT_REUSE_MATRIX, PETSC_DEFAULT, &MatsFwdSensipTemp[i])); 70413af1a74SHong Zhang if (ts->Jacprhs) { 7059566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobianP(ts, stage_time, Y[i], ts->Jacprhs)); /* get f_p */ 70613af1a74SHong Zhang if (ts->vecs_sensi2p) { /* TLM used for 2nd-order adjoint */ 70713af1a74SHong Zhang PetscScalar *xarr; 7089566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(MatsFwdSensipTemp[i], 0, &xarr)); 7099566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(rk->VecDeltaFwdSensipCol, xarr)); 7109566063dSJacob Faibussowitsch PetscCall(MatMultAdd(ts->Jacprhs, ts->vec_dir, rk->VecDeltaFwdSensipCol, rk->VecDeltaFwdSensipCol)); 7119566063dSJacob Faibussowitsch PetscCall(VecResetArray(rk->VecDeltaFwdSensipCol)); 7129566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(MatsFwdSensipTemp[i], &xarr)); 71313af1a74SHong Zhang } else { 7149566063dSJacob Faibussowitsch PetscCall(MatAXPY(MatsFwdSensipTemp[i], 1., ts->Jacprhs, SUBSET_NONZERO_PATTERN)); 71513af1a74SHong Zhang } 716922a638cSHong Zhang } 717922a638cSHong Zhang } 718922a638cSHong Zhang 71948a46eb9SPierre Jolivet for (i = 0; i < s; i++) PetscCall(MatAXPY(ts->mat_sensip, h * b[i], rk->MatsFwdSensipTemp[i], SAME_NONZERO_PATTERN)); 720922a638cSHong Zhang rk->status = TS_STEP_COMPLETE; 7213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 722922a638cSHong Zhang } 723922a638cSHong Zhang 724d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardGetStages_RK(TS ts, PetscInt *ns, Mat **stagesensip) 725d71ae5a4SJacob Faibussowitsch { 726922a638cSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 727922a638cSHong Zhang RKTableau tab = rk->tableau; 728922a638cSHong Zhang 729922a638cSHong Zhang PetscFunctionBegin; 730922a638cSHong Zhang if (ns) *ns = tab->s; 731922a638cSHong Zhang if (stagesensip) *stagesensip = rk->MatsFwdStageSensip; 7323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 733922a638cSHong Zhang } 734922a638cSHong Zhang 735d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardSetUp_RK(TS ts) 736d71ae5a4SJacob Faibussowitsch { 737922a638cSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 738922a638cSHong Zhang RKTableau tab = rk->tableau; 739922a638cSHong Zhang PetscInt i; 740922a638cSHong Zhang 741922a638cSHong Zhang PetscFunctionBegin; 742922a638cSHong Zhang /* backup sensitivity results for roll-backs */ 7439566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &rk->MatFwdSensip0)); 744922a638cSHong Zhang 7459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &rk->MatsFwdStageSensip)); 7469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &rk->MatsFwdSensipTemp)); 747922a638cSHong Zhang for (i = 0; i < tab->s; i++) { 7489566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &rk->MatsFwdStageSensip[i])); 7499566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &rk->MatsFwdSensipTemp[i])); 750922a638cSHong Zhang } 7519566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &rk->VecDeltaFwdSensipCol)); 7523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 753922a638cSHong Zhang } 754922a638cSHong Zhang 755d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardReset_RK(TS ts) 756d71ae5a4SJacob Faibussowitsch { 757922a638cSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 758922a638cSHong Zhang RKTableau tab = rk->tableau; 759922a638cSHong Zhang PetscInt i; 760922a638cSHong Zhang 761922a638cSHong Zhang PetscFunctionBegin; 7629566063dSJacob Faibussowitsch PetscCall(MatDestroy(&rk->MatFwdSensip0)); 763922a638cSHong Zhang if (rk->MatsFwdStageSensip) { 76448a46eb9SPierre Jolivet for (i = 0; i < tab->s; i++) PetscCall(MatDestroy(&rk->MatsFwdStageSensip[i])); 7659566063dSJacob Faibussowitsch PetscCall(PetscFree(rk->MatsFwdStageSensip)); 766922a638cSHong Zhang } 767922a638cSHong Zhang if (rk->MatsFwdSensipTemp) { 76848a46eb9SPierre Jolivet for (i = 0; i < tab->s; i++) PetscCall(MatDestroy(&rk->MatsFwdSensipTemp[i])); 7699566063dSJacob Faibussowitsch PetscCall(PetscFree(rk->MatsFwdSensipTemp)); 770922a638cSHong Zhang } 7719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rk->VecDeltaFwdSensipCol)); 7723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 773922a638cSHong Zhang } 774922a638cSHong Zhang 775d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_RK(TS ts) 776d71ae5a4SJacob Faibussowitsch { 777f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK *)ts->data; 778f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 779f68a32c8SEmil Constantinescu const PetscInt s = tab->s; 780fecfb714SLisandro Dalcin const PetscReal *A = tab->A, *c = tab->c; 781f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 782f68a32c8SEmil Constantinescu Vec *Y = rk->Y, *YdotRHS = rk->YdotRHS; 783630f8c86SStefano Zampini PetscBool FSAL = (PetscBool)(tab->FSAL && !rk->newtableau); 784f68a32c8SEmil Constantinescu TSAdapt adapt; 785fecfb714SLisandro Dalcin PetscInt i, j; 786be5899b3SLisandro Dalcin PetscInt rejections = 0; 787be5899b3SLisandro Dalcin PetscBool stageok, accept = PETSC_TRUE; 788be5899b3SLisandro Dalcin PetscReal next_time_step = ts->time_step; 789f68a32c8SEmil Constantinescu 790f68a32c8SEmil Constantinescu PetscFunctionBegin; 791d1905564SLisandro Dalcin if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE; 7929566063dSJacob Faibussowitsch if (FSAL) PetscCall(VecCopy(YdotRHS[s - 1], YdotRHS[0])); 793630f8c86SStefano Zampini rk->newtableau = PETSC_FALSE; 794d1905564SLisandro Dalcin 795f68a32c8SEmil Constantinescu rk->status = TS_STEP_INCOMPLETE; 796be5899b3SLisandro Dalcin while (!ts->reason && rk->status != TS_STEP_COMPLETE) { 797be5899b3SLisandro Dalcin PetscReal t = ts->ptime; 798f68a32c8SEmil Constantinescu PetscReal h = ts->time_step; 799f68a32c8SEmil Constantinescu for (i = 0; i < s; i++) { 8009be3e283SDebojyoti Ghosh rk->stage_time = t + h * c[i]; 8019566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, rk->stage_time)); 8029566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, Y[i])); 803f68a32c8SEmil Constantinescu for (j = 0; j < i; j++) w[j] = h * A[i * s + j]; 8049566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Y[i], i, w, YdotRHS)); 8059566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, rk->stage_time, i, Y)); 8069566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 8079566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(adapt, ts, rk->stage_time, Y[i], &stageok)); 808fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 8098f738a7cSLisandro Dalcin if (FSAL && !i) continue; 8109566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, t + h * c[i], Y[i], YdotRHS[i])); 811f68a32c8SEmil Constantinescu } 812be5899b3SLisandro Dalcin 813be5899b3SLisandro Dalcin rk->status = TS_STEP_INCOMPLETE; 8149566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep(ts, tab->order, ts->vec_sol, NULL)); 815f68a32c8SEmil Constantinescu rk->status = TS_STEP_PENDING; 8169566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 8179566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(adapt)); 8189566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidateAdd(adapt, tab->name, tab->order, 1, tab->ccfl, (PetscReal)tab->s, PETSC_TRUE)); 8199566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 820be5899b3SLisandro Dalcin rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 821be5899b3SLisandro Dalcin if (!accept) { /* Roll back the current step */ 8229566063dSJacob Faibussowitsch PetscCall(TSRollBack_RK(ts)); 823be5899b3SLisandro Dalcin ts->time_step = next_time_step; 824be5899b3SLisandro Dalcin goto reject_step; 825be5899b3SLisandro Dalcin } 826be5899b3SLisandro Dalcin 8270f7a1166SHong Zhang if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 8280f7a1166SHong Zhang rk->ptime = ts->ptime; 8290f7a1166SHong Zhang rk->time_step = ts->time_step; 830493ed95fSHong Zhang } 831be5899b3SLisandro Dalcin 832f68a32c8SEmil Constantinescu ts->ptime += ts->time_step; 833f68a32c8SEmil Constantinescu ts->time_step = next_time_step; 834f68a32c8SEmil Constantinescu break; 835be5899b3SLisandro Dalcin 836be5899b3SLisandro Dalcin reject_step: 8379371c9d4SSatish Balay ts->reject++; 8389371c9d4SSatish Balay accept = PETSC_FALSE; 839be5899b3SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 840be5899b3SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 84163a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections)); 842f68a32c8SEmil Constantinescu } 843f68a32c8SEmil Constantinescu } 8443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 845e4dd521cSBarry Smith } 846e4dd521cSBarry Smith 847d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointSetUp_RK(TS ts) 848d71ae5a4SJacob Faibussowitsch { 849f6a906c0SBarry Smith TS_RK *rk = (TS_RK *)ts->data; 850be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 851be5899b3SLisandro Dalcin PetscInt s = tab->s; 852f6a906c0SBarry Smith 853f6a906c0SBarry Smith PetscFunctionBegin; 8543ba16761SJacob Faibussowitsch if (ts->adjointsetupcalled++) PetscFunctionReturn(PETSC_SUCCESS); 8559566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], s * ts->numcost, &rk->VecsDeltaLam)); 8569566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &rk->VecsSensiTemp)); 85748a46eb9SPierre Jolivet if (ts->vecs_sensip) PetscCall(VecDuplicate(ts->vecs_sensip[0], &rk->VecDeltaMu)); 85813af1a74SHong Zhang if (ts->vecs_sensi2) { 8599566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], s * ts->numcost, &rk->VecsDeltaLam2)); 8609566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vecs_sensi2[0], ts->numcost, &rk->VecsSensi2Temp)); 86113af1a74SHong Zhang } 86248a46eb9SPierre Jolivet if (ts->vecs_sensi2p) PetscCall(VecDuplicate(ts->vecs_sensi2p[0], &rk->VecDeltaMu2)); 8633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 864f6a906c0SBarry Smith } 865f6a906c0SBarry Smith 86635f5172aSHong Zhang /* 86735f5172aSHong Zhang Assumptions: 86835f5172aSHong Zhang - TSStep_RK() always evaluates the step with b, not bembed. 86935f5172aSHong Zhang */ 870d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointStep_RK(TS ts) 871d71ae5a4SJacob Faibussowitsch { 872c235aa8dSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 873cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 874c235aa8dSHong Zhang RKTableau tab = rk->tableau; 8753ca0882eSHong Zhang Mat J, Jpre, Jquad; 876c235aa8dSHong Zhang const PetscInt s = tab->s; 877c235aa8dSHong Zhang const PetscReal *A = tab->A, *b = tab->b, *c = tab->c; 87813af1a74SHong Zhang PetscScalar *w = rk->work, *xarr; 8792e7b7f96SHong Zhang Vec *Y = rk->Y, *VecsDeltaLam = rk->VecsDeltaLam, VecDeltaMu = rk->VecDeltaMu, *VecsSensiTemp = rk->VecsSensiTemp; 88013af1a74SHong Zhang Vec *VecsDeltaLam2 = rk->VecsDeltaLam2, VecDeltaMu2 = rk->VecDeltaMu2, *VecsSensi2Temp = rk->VecsSensi2Temp; 881cd4cee2dSHong Zhang Vec VecDRDUTransCol = ts->vec_drdu_col, VecDRDPTransCol = ts->vec_drdp_col; 882b18ea86cSHong Zhang PetscInt i, j, nadj; 8833d3eaba7SBarry Smith PetscReal t = ts->ptime; 8843ca0882eSHong Zhang PetscReal h = ts->time_step; 885c235aa8dSHong Zhang 886d2daff3dSHong Zhang PetscFunctionBegin; 887c235aa8dSHong Zhang rk->status = TS_STEP_INCOMPLETE; 8883389c7d9SStefano Zampini 8899566063dSJacob Faibussowitsch PetscCall(TSGetRHSJacobian(ts, &J, &Jpre, NULL, NULL)); 89048a46eb9SPierre Jolivet if (quadts) PetscCall(TSGetRHSJacobian(quadts, &Jquad, NULL, NULL, NULL)); 8912e7b7f96SHong Zhang for (i = s - 1; i >= 0; i--) { 8926a1e4597SHong Zhang if (tab->FSAL && i == s - 1) { 8936a1e4597SHong Zhang /* VecsDeltaLam[nadj*s+s-1] are initialized with zeros and the values never change.*/ 8946a1e4597SHong Zhang continue; 8956a1e4597SHong Zhang } 8963ca0882eSHong Zhang rk->stage_time = t + h * (1.0 - c[i]); 8979566063dSJacob Faibussowitsch PetscCall(TSComputeSNESJacobian(ts, Y[i], J, Jpre)); 8989371c9d4SSatish Balay if (quadts) { PetscCall(TSComputeRHSJacobian(quadts, rk->stage_time, Y[i], Jquad, Jquad)); /* get r_u^T */ } 8993389c7d9SStefano Zampini if (ts->vecs_sensip) { 9009566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobianP(ts, rk->stage_time, Y[i], ts->Jacprhs)); /* get f_p */ 9019371c9d4SSatish Balay if (quadts) { PetscCall(TSComputeRHSJacobianP(quadts, rk->stage_time, Y[i], quadts->Jacprhs)); /* get f_p for the quadrature */ } 90235f5172aSHong Zhang } 9033389c7d9SStefano Zampini 90435f5172aSHong Zhang if (b[i]) { 90535f5172aSHong Zhang for (j = i + 1; j < s; j++) w[j - i - 1] = A[j * s + i] / b[i]; /* coefficients for computing VecsSensiTemp */ 90635f5172aSHong Zhang } else { 90735f5172aSHong Zhang for (j = i + 1; j < s; j++) w[j - i - 1] = A[j * s + i]; /* coefficients for computing VecsSensiTemp */ 90835f5172aSHong Zhang } 9092e7b7f96SHong Zhang 9102e7b7f96SHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 9113389c7d9SStefano Zampini /* Stage values of lambda */ 91235f5172aSHong Zhang if (b[i]) { 91335f5172aSHong Zhang /* lambda_{n+1} + \sum_{j=i+1}^s a_{ji}/b[i]*lambda_{s,j} */ 9149566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vecs_sensi[nadj], VecsSensiTemp[nadj])); /* VecDeltaLam is an vec array of size s by numcost */ 9159566063dSJacob Faibussowitsch PetscCall(VecMAXPY(VecsSensiTemp[nadj], s - i - 1, w, &VecsDeltaLam[nadj * s + i + 1])); 9169566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(J, VecsSensiTemp[nadj], VecsDeltaLam[nadj * s + i])); /* VecsSensiTemp will be reused by 2nd-order adjoint */ 9179566063dSJacob Faibussowitsch PetscCall(VecScale(VecsDeltaLam[nadj * s + i], -h * b[i])); 91835f5172aSHong Zhang if (quadts) { 9199566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(Jquad, nadj, &xarr)); 9209566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(VecDRDUTransCol, xarr)); 9219566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecsDeltaLam[nadj * s + i], -h * b[i], VecDRDUTransCol)); 9229566063dSJacob Faibussowitsch PetscCall(VecResetArray(VecDRDUTransCol)); 9239566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(Jquad, &xarr)); 924cd4cee2dSHong Zhang } 9253389c7d9SStefano Zampini } else { 9266a1e4597SHong Zhang /* \sum_{j=i+1}^s a_{ji}*lambda_{s,j} */ 9279566063dSJacob Faibussowitsch PetscCall(VecSet(VecsSensiTemp[nadj], 0)); 9289566063dSJacob Faibussowitsch PetscCall(VecMAXPY(VecsSensiTemp[nadj], s - i - 1, w, &VecsDeltaLam[nadj * s + i + 1])); 9299566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(J, VecsSensiTemp[nadj], VecsDeltaLam[nadj * s + i])); 9309566063dSJacob Faibussowitsch PetscCall(VecScale(VecsDeltaLam[nadj * s + i], -h)); 9313389c7d9SStefano Zampini } 9326497ce14SHong Zhang 933ad8e2604SHong Zhang /* Stage values of mu */ 9346497ce14SHong Zhang if (ts->vecs_sensip) { 93535f5172aSHong Zhang if (b[i]) { 9369566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ts->Jacprhs, VecsSensiTemp[nadj], VecDeltaMu)); 9379566063dSJacob Faibussowitsch PetscCall(VecScale(VecDeltaMu, -h * b[i])); 93835f5172aSHong Zhang if (quadts) { 9399566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(quadts->Jacprhs, nadj, &xarr)); 9409566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(VecDRDPTransCol, xarr)); 9419566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecDeltaMu, -h * b[i], VecDRDPTransCol)); 9429566063dSJacob Faibussowitsch PetscCall(VecResetArray(VecDRDPTransCol)); 9439566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(quadts->Jacprhs, &xarr)); 944493ed95fSHong Zhang } 94535f5172aSHong Zhang } else { 9469566063dSJacob Faibussowitsch PetscCall(VecScale(VecDeltaMu, -h)); 94735f5172aSHong Zhang } 9489566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensip[nadj], 1., VecDeltaMu)); /* update sensip for each stage */ 949ad8e2604SHong Zhang } 950c235aa8dSHong Zhang } 95113af1a74SHong Zhang 95213af1a74SHong Zhang if (ts->vecs_sensi2 && ts->forward_solve) { /* 2nd-order adjoint, TLM mode has to be turned on */ 95313af1a74SHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 9549566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(rk->MatsFwdStageSensip[i], 0, &xarr)); 9559566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr)); 95613af1a74SHong Zhang /* lambda_s^T F_UU w_1 */ 9579566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionUU(ts, rk->stage_time, Y[i], VecsSensiTemp, ts->vec_sensip_col, ts->vecs_guu)); 95835f5172aSHong Zhang if (quadts) { 95935f5172aSHong Zhang /* R_UU w_1 */ 9609566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionUU(quadts, rk->stage_time, Y[i], NULL, ts->vec_sensip_col, ts->vecs_guu)); 96135f5172aSHong Zhang } 96235f5172aSHong Zhang if (ts->vecs_sensip) { 96313af1a74SHong Zhang /* lambda_s^T F_UP w_2 */ 9649566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionUP(ts, rk->stage_time, Y[i], VecsSensiTemp, ts->vec_dir, ts->vecs_gup)); 96535f5172aSHong Zhang if (quadts) { 96635f5172aSHong Zhang /* R_UP w_2 */ 9679566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionUP(quadts, rk->stage_time, Y[i], NULL, ts->vec_sensip_col, ts->vecs_gup)); 96835f5172aSHong Zhang } 96935f5172aSHong Zhang } 97013af1a74SHong Zhang if (ts->vecs_sensi2p) { 97113af1a74SHong Zhang /* lambda_s^T F_PU w_1 */ 9729566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionPU(ts, rk->stage_time, Y[i], VecsSensiTemp, ts->vec_sensip_col, ts->vecs_gpu)); 97335f5172aSHong Zhang /* lambda_s^T F_PP w_2 */ 9749566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionPP(ts, rk->stage_time, Y[i], VecsSensiTemp, ts->vec_dir, ts->vecs_gpp)); 97535f5172aSHong Zhang if (b[i] && quadts) { 97635f5172aSHong Zhang /* R_PU w_1 */ 9779566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionPU(quadts, rk->stage_time, Y[i], NULL, ts->vec_sensip_col, ts->vecs_gpu)); 97835f5172aSHong Zhang /* R_PP w_2 */ 9799566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionPP(quadts, rk->stage_time, Y[i], NULL, ts->vec_dir, ts->vecs_gpp)); 98035f5172aSHong Zhang } 98113af1a74SHong Zhang } 9829566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_sensip_col)); 9839566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(rk->MatsFwdStageSensip[i], &xarr)); 98413af1a74SHong Zhang 98513af1a74SHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 98613af1a74SHong Zhang /* Stage values of lambda */ 98735f5172aSHong Zhang if (b[i]) { 98835f5172aSHong Zhang /* J_i^T*(Lambda_{n+1}+\sum_{j=i+1}^s a_{ji}/b_i*Lambda_{s,j} */ 9899566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vecs_sensi2[nadj], VecsSensi2Temp[nadj])); 9909566063dSJacob Faibussowitsch PetscCall(VecMAXPY(VecsSensi2Temp[nadj], s - i - 1, w, &VecsDeltaLam2[nadj * s + i + 1])); 9919566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(J, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj * s + i])); 9929566063dSJacob Faibussowitsch PetscCall(VecScale(VecsDeltaLam2[nadj * s + i], -h * b[i])); 9939566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecsDeltaLam2[nadj * s + i], -h * b[i], ts->vecs_guu[nadj])); 99448a46eb9SPierre Jolivet if (ts->vecs_sensip) PetscCall(VecAXPY(VecsDeltaLam2[nadj * s + i], -h * b[i], ts->vecs_gup[nadj])); 99513af1a74SHong Zhang } else { 99635f5172aSHong Zhang /* \sum_{j=i+1}^s a_{ji}*Lambda_{s,j} */ 9979566063dSJacob Faibussowitsch PetscCall(VecSet(VecsDeltaLam2[nadj * s + i], 0)); 9989566063dSJacob Faibussowitsch PetscCall(VecMAXPY(VecsSensi2Temp[nadj], s - i - 1, w, &VecsDeltaLam2[nadj * s + i + 1])); 9999566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(J, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj * s + i])); 10009566063dSJacob Faibussowitsch PetscCall(VecScale(VecsDeltaLam2[nadj * s + i], -h)); 10019566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecsDeltaLam2[nadj * s + i], -h, ts->vecs_guu[nadj])); 100248a46eb9SPierre Jolivet if (ts->vecs_sensip) PetscCall(VecAXPY(VecsDeltaLam2[nadj * s + i], -h, ts->vecs_gup[nadj])); 100335f5172aSHong Zhang } 100435f5172aSHong Zhang if (ts->vecs_sensi2p) { /* 2nd-order adjoint for parameters */ 10059566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ts->Jacprhs, VecsSensi2Temp[nadj], VecDeltaMu2)); 100635f5172aSHong Zhang if (b[i]) { 10079566063dSJacob Faibussowitsch PetscCall(VecScale(VecDeltaMu2, -h * b[i])); 10089566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecDeltaMu2, -h * b[i], ts->vecs_gpu[nadj])); 10099566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecDeltaMu2, -h * b[i], ts->vecs_gpp[nadj])); 101013af1a74SHong Zhang } else { 10119566063dSJacob Faibussowitsch PetscCall(VecScale(VecDeltaMu2, -h)); 10129566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecDeltaMu2, -h, ts->vecs_gpu[nadj])); 10139566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecDeltaMu2, -h, ts->vecs_gpp[nadj])); 101413af1a74SHong Zhang } 10159566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], 1, VecDeltaMu2)); /* update sensi2p for each stage */ 101613af1a74SHong Zhang } 101713af1a74SHong Zhang } 101813af1a74SHong Zhang } 10196497ce14SHong Zhang } 1020c235aa8dSHong Zhang 1021c235aa8dSHong Zhang for (j = 0; j < s; j++) w[j] = 1.0; 10222e7b7f96SHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { /* no need to do this for mu's */ 10239566063dSJacob Faibussowitsch PetscCall(VecMAXPY(ts->vecs_sensi[nadj], s, w, &VecsDeltaLam[nadj * s])); 102448a46eb9SPierre Jolivet if (ts->vecs_sensi2) PetscCall(VecMAXPY(ts->vecs_sensi2[nadj], s, w, &VecsDeltaLam2[nadj * s])); 10256497ce14SHong Zhang } 1026c235aa8dSHong Zhang rk->status = TS_STEP_COMPLETE; 10273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1028d2daff3dSHong Zhang } 1029d2daff3dSHong Zhang 1030d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointReset_RK(TS ts) 1031d71ae5a4SJacob Faibussowitsch { 103213af1a74SHong Zhang TS_RK *rk = (TS_RK *)ts->data; 103313af1a74SHong Zhang RKTableau tab = rk->tableau; 103413af1a74SHong Zhang 103513af1a74SHong Zhang PetscFunctionBegin; 10369566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s * ts->numcost, &rk->VecsDeltaLam)); 10379566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(ts->numcost, &rk->VecsSensiTemp)); 10389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rk->VecDeltaMu)); 10399566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s * ts->numcost, &rk->VecsDeltaLam2)); 10409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rk->VecDeltaMu2)); 10419566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(ts->numcost, &rk->VecsSensi2Temp)); 10423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 104313af1a74SHong Zhang } 104413af1a74SHong Zhang 1045d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_RK(TS ts, PetscReal itime, Vec X) 1046d71ae5a4SJacob Faibussowitsch { 104755de54fcSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 104855de54fcSHong Zhang PetscInt s = rk->tableau->s, p = rk->tableau->p, i, j; 104955de54fcSHong Zhang PetscReal h; 105055de54fcSHong Zhang PetscReal tt, t; 105155de54fcSHong Zhang PetscScalar *b; 105255de54fcSHong Zhang const PetscReal *B = rk->tableau->binterp; 105355de54fcSHong Zhang 105455de54fcSHong Zhang PetscFunctionBegin; 10553c633725SBarry Smith PetscCheck(B, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSRK %s does not have an interpolation formula", rk->tableau->name); 105655de54fcSHong Zhang 105755de54fcSHong Zhang switch (rk->status) { 105855de54fcSHong Zhang case TS_STEP_INCOMPLETE: 105955de54fcSHong Zhang case TS_STEP_PENDING: 106055de54fcSHong Zhang h = ts->time_step; 106155de54fcSHong Zhang t = (itime - ts->ptime) / h; 106255de54fcSHong Zhang break; 106355de54fcSHong Zhang case TS_STEP_COMPLETE: 106455de54fcSHong Zhang h = ts->ptime - ts->ptime_prev; 106555de54fcSHong Zhang t = (itime - ts->ptime) / h + 1; /* In the interval [0,1] */ 106655de54fcSHong Zhang break; 1067d71ae5a4SJacob Faibussowitsch default: 1068d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus"); 106955de54fcSHong Zhang } 10709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s, &b)); 107155de54fcSHong Zhang for (i = 0; i < s; i++) b[i] = 0; 107255de54fcSHong Zhang for (j = 0, tt = t; j < p; j++, tt *= t) { 1073ad540459SPierre Jolivet for (i = 0; i < s; i++) b[i] += h * B[i * p + j] * tt; 107455de54fcSHong Zhang } 10759566063dSJacob Faibussowitsch PetscCall(VecCopy(rk->Y[0], X)); 10769566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, b, rk->YdotRHS)); 10779566063dSJacob Faibussowitsch PetscCall(PetscFree(b)); 10783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 107955de54fcSHong Zhang } 108055de54fcSHong Zhang 108155de54fcSHong Zhang /*------------------------------------------------------------*/ 108255de54fcSHong Zhang 1083d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRKTableauReset(TS ts) 1084d71ae5a4SJacob Faibussowitsch { 1085be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK *)ts->data; 1086be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 1087be5899b3SLisandro Dalcin 1088be5899b3SLisandro Dalcin PetscFunctionBegin; 10893ba16761SJacob Faibussowitsch if (!tab) PetscFunctionReturn(PETSC_SUCCESS); 10909566063dSJacob Faibussowitsch PetscCall(PetscFree(rk->work)); 10919566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &rk->Y)); 10929566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &rk->YdotRHS)); 10933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1094be5899b3SLisandro Dalcin } 1095be5899b3SLisandro Dalcin 1096d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_RK(TS ts) 1097d71ae5a4SJacob Faibussowitsch { 1098e4dd521cSBarry Smith PetscFunctionBegin; 10999566063dSJacob Faibussowitsch PetscCall(TSRKTableauReset(ts)); 11000fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 1101cac4c232SBarry Smith PetscTryMethod(ts, "TSReset_RK_MultirateSplit_C", (TS), (ts)); 11020fe4d17eSHong Zhang } else { 1103cac4c232SBarry Smith PetscTryMethod(ts, "TSReset_RK_MultirateNonsplit_C", (TS), (ts)); 11040fe4d17eSHong Zhang } 11053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1106e4dd521cSBarry Smith } 1107277b19d0SLisandro Dalcin 1108d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCoarsenHook_TSRK(DM fine, DM coarse, void *ctx) 1109d71ae5a4SJacob Faibussowitsch { 1110f68a32c8SEmil Constantinescu PetscFunctionBegin; 11113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1112f68a32c8SEmil Constantinescu } 1113f68a32c8SEmil Constantinescu 1114d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRestrictHook_TSRK(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) 1115d71ae5a4SJacob Faibussowitsch { 1116f68a32c8SEmil Constantinescu PetscFunctionBegin; 11173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1118f68a32c8SEmil Constantinescu } 1119f68a32c8SEmil Constantinescu 1120d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainHook_TSRK(DM dm, DM subdm, void *ctx) 1121d71ae5a4SJacob Faibussowitsch { 1122f68a32c8SEmil Constantinescu PetscFunctionBegin; 11233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1124f68a32c8SEmil Constantinescu } 1125f68a32c8SEmil Constantinescu 1126d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) 1127d71ae5a4SJacob Faibussowitsch { 1128f68a32c8SEmil Constantinescu PetscFunctionBegin; 11293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1130f68a32c8SEmil Constantinescu } 1131be5899b3SLisandro Dalcin 1132d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRKTableauSetUp(TS ts) 1133d71ae5a4SJacob Faibussowitsch { 1134be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK *)ts->data; 1135be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 1136be5899b3SLisandro Dalcin 1137be5899b3SLisandro Dalcin PetscFunctionBegin; 11389566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &rk->work)); 11399566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &rk->Y)); 11409566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &rk->YdotRHS)); 1141630f8c86SStefano Zampini rk->newtableau = PETSC_TRUE; 11423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1143be5899b3SLisandro Dalcin } 1144be5899b3SLisandro Dalcin 1145d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_RK(TS ts) 1146d71ae5a4SJacob Faibussowitsch { 1147cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 1148f68a32c8SEmil Constantinescu DM dm; 1149f68a32c8SEmil Constantinescu 1150f68a32c8SEmil Constantinescu PetscFunctionBegin; 11519566063dSJacob Faibussowitsch PetscCall(TSCheckImplicitTerm(ts)); 11529566063dSJacob Faibussowitsch PetscCall(TSRKTableauSetUp(ts)); 1153cd4cee2dSHong Zhang if (quadts && ts->costintegralfwd) { 1154cd4cee2dSHong Zhang Mat Jquad; 11559566063dSJacob Faibussowitsch PetscCall(TSGetRHSJacobian(quadts, &Jquad, NULL, NULL, NULL)); 11560f7a1166SHong Zhang } 11579566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 11589566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSRK, DMRestrictHook_TSRK, ts)); 11599566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSRK, DMSubDomainRestrictHook_TSRK, ts)); 11600fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 1161cac4c232SBarry Smith PetscTryMethod(ts, "TSSetUp_RK_MultirateSplit_C", (TS), (ts)); 11620fe4d17eSHong Zhang } else { 1163cac4c232SBarry Smith PetscTryMethod(ts, "TSSetUp_RK_MultirateNonsplit_C", (TS), (ts)); 11640fe4d17eSHong Zhang } 11653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1166e4dd521cSBarry Smith } 1167d2daff3dSHong Zhang 1168d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_RK(TS ts, PetscOptionItems *PetscOptionsObject) 1169d71ae5a4SJacob Faibussowitsch { 1170be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK *)ts->data; 1171262ff7bbSBarry Smith 1172e4dd521cSBarry Smith PetscFunctionBegin; 1173d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "RK ODE solver options"); 1174f68a32c8SEmil Constantinescu { 1175f68a32c8SEmil Constantinescu RKTableauLink link; 1176f68a32c8SEmil Constantinescu PetscInt count, choice; 11770fe4d17eSHong Zhang PetscBool flg, use_multirate = PETSC_FALSE; 1178f68a32c8SEmil Constantinescu const char **namelist; 11792c9cb062Sluzhanghpp 1180fbccb6d4SPierre Jolivet for (link = RKTableauList, count = 0; link; link = link->next, count++); 11819566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count, (char ***)&namelist)); 1182f68a32c8SEmil Constantinescu for (link = RKTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name; 11839566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_rk_multirate", "Use interpolation-based multirate RK method", "TSRKSetMultirate", rk->use_multirate, &use_multirate, &flg)); 11841baa6e33SBarry Smith if (flg) PetscCall(TSRKSetMultirate(ts, use_multirate)); 11859566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_rk_type", "Family of RK method", "TSRKSetType", (const char *const *)namelist, count, rk->tableau->name, &choice, &flg)); 11869566063dSJacob Faibussowitsch if (flg) PetscCall(TSRKSetType(ts, namelist[choice])); 11879566063dSJacob Faibussowitsch PetscCall(PetscFree(namelist)); 1188f68a32c8SEmil Constantinescu } 1189d0609cedSBarry Smith PetscOptionsHeadEnd(); 1190d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)ts), NULL, "Multirate methods options", ""); 11919566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_rk_dtratio", "time step ratio between slow and fast", "", rk->dtratio, &rk->dtratio, NULL)); 1192d0609cedSBarry Smith PetscOptionsEnd(); 11933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1194e4dd521cSBarry Smith } 1195e4dd521cSBarry Smith 1196d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_RK(TS ts, PetscViewer viewer) 1197d71ae5a4SJacob Faibussowitsch { 11985f70b478SJed Brown TS_RK *rk = (TS_RK *)ts->data; 11998370ee3bSLisandro Dalcin PetscBool iascii; 12002cdf8120Sasbjorn 1201e4dd521cSBarry Smith PetscFunctionBegin; 12029566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 12038370ee3bSLisandro Dalcin if (iascii) { 12049c334d8fSLisandro Dalcin RKTableau tab = rk->tableau; 1205f68a32c8SEmil Constantinescu TSRKType rktype; 12060f7a28cdSStefano Zampini const PetscReal *c; 12070f7a28cdSStefano Zampini PetscInt s; 1208f68a32c8SEmil Constantinescu char buf[512]; 12090f7a28cdSStefano Zampini PetscBool FSAL; 12100f7a28cdSStefano Zampini 12119566063dSJacob Faibussowitsch PetscCall(TSRKGetType(ts, &rktype)); 12129566063dSJacob Faibussowitsch PetscCall(TSRKGetTableau(ts, &s, NULL, NULL, &c, NULL, NULL, NULL, &FSAL)); 12139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " RK type %s\n", rktype)); 121463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Order: %" PetscInt_FMT "\n", tab->order)); 12159566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " FSAL property: %s\n", FSAL ? "yes" : "no")); 12169566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", s, c)); 12179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Abscissa c = %s\n", buf)); 12188370ee3bSLisandro Dalcin } 12193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1220f68a32c8SEmil Constantinescu } 1221f68a32c8SEmil Constantinescu 1222d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSLoad_RK(TS ts, PetscViewer viewer) 1223d71ae5a4SJacob Faibussowitsch { 12249c334d8fSLisandro Dalcin TSAdapt adapt; 1225f68a32c8SEmil Constantinescu 1226f68a32c8SEmil Constantinescu PetscFunctionBegin; 12279566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 12289566063dSJacob Faibussowitsch PetscCall(TSAdaptLoad(adapt, viewer)); 12293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1230f68a32c8SEmil Constantinescu } 1231f68a32c8SEmil Constantinescu 12322ea87ba9SLisandro Dalcin /*@ 1233bcf0153eSBarry Smith TSRKGetOrder - Get the order of the `TSRK` scheme 12342ea87ba9SLisandro Dalcin 123520f4b53cSBarry Smith Not Collective 12362ea87ba9SLisandro Dalcin 12372ea87ba9SLisandro Dalcin Input Parameter: 12382ea87ba9SLisandro Dalcin . ts - timestepping context 12392ea87ba9SLisandro Dalcin 12402ea87ba9SLisandro Dalcin Output Parameter: 1241bcf0153eSBarry Smith . order - order of `TSRK` scheme 12422ea87ba9SLisandro Dalcin 12432ea87ba9SLisandro Dalcin Level: intermediate 12442ea87ba9SLisandro Dalcin 12451cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKGetType()` 12462ea87ba9SLisandro Dalcin @*/ 1247d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKGetOrder(TS ts, PetscInt *order) 1248d71ae5a4SJacob Faibussowitsch { 12492ea87ba9SLisandro Dalcin PetscFunctionBegin; 12502ea87ba9SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 12514f572ea9SToby Isaac PetscAssertPointer(order, 2); 1252cac4c232SBarry Smith PetscUseMethod(ts, "TSRKGetOrder_C", (TS, PetscInt *), (ts, order)); 12533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12542ea87ba9SLisandro Dalcin } 12552ea87ba9SLisandro Dalcin 1256*cc4c1da9SBarry Smith /*@ 1257bcf0153eSBarry Smith TSRKSetType - Set the type of the `TSRK` scheme 1258f68a32c8SEmil Constantinescu 125920f4b53cSBarry Smith Logically Collective 1260f68a32c8SEmil Constantinescu 1261d8d19677SJose E. Roman Input Parameters: 1262f68a32c8SEmil Constantinescu + ts - timestepping context 1263bcf0153eSBarry Smith - rktype - type of `TSRK` scheme 1264f68a32c8SEmil Constantinescu 1265bcf0153eSBarry Smith Options Database Key: 12669bd3e852SBarry Smith . -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs> 1267287c2655SBarry Smith 1268f68a32c8SEmil Constantinescu Level: intermediate 1269f68a32c8SEmil Constantinescu 12701cc06b55SBarry Smith .seealso: [](ch_ts), `TSRKGetType()`, `TSRK`, `TSRKType`, `TSRK1FE`, `TSRK2A`, `TSRK2B`, `TSRK3`, `TSRK3BS`, `TSRK4`, `TSRK5F`, `TSRK5DP`, `TSRK5BS`, `TSRK6VR`, `TSRK7VR`, `TSRK8VR` 1271f68a32c8SEmil Constantinescu @*/ 1272d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKSetType(TS ts, TSRKType rktype) 1273d71ae5a4SJacob Faibussowitsch { 1274f68a32c8SEmil Constantinescu PetscFunctionBegin; 1275f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 12764f572ea9SToby Isaac PetscAssertPointer(rktype, 2); 1277cac4c232SBarry Smith PetscTryMethod(ts, "TSRKSetType_C", (TS, TSRKType), (ts, rktype)); 12783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1279f68a32c8SEmil Constantinescu } 1280f68a32c8SEmil Constantinescu 1281*cc4c1da9SBarry Smith /*@ 1282bcf0153eSBarry Smith TSRKGetType - Get the type of `TSRK` scheme 1283f68a32c8SEmil Constantinescu 128420f4b53cSBarry Smith Not Collective 1285f68a32c8SEmil Constantinescu 1286f68a32c8SEmil Constantinescu Input Parameter: 1287f68a32c8SEmil Constantinescu . ts - timestepping context 1288f68a32c8SEmil Constantinescu 1289f68a32c8SEmil Constantinescu Output Parameter: 1290bcf0153eSBarry Smith . rktype - type of `TSRK`-scheme 1291f68a32c8SEmil Constantinescu 1292f68a32c8SEmil Constantinescu Level: intermediate 1293f68a32c8SEmil Constantinescu 12941cc06b55SBarry Smith .seealso: [](ch_ts), `TSRKSetType()` 1295f68a32c8SEmil Constantinescu @*/ 1296d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKGetType(TS ts, TSRKType *rktype) 1297d71ae5a4SJacob Faibussowitsch { 1298f68a32c8SEmil Constantinescu PetscFunctionBegin; 1299f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1300cac4c232SBarry Smith PetscUseMethod(ts, "TSRKGetType_C", (TS, TSRKType *), (ts, rktype)); 13013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1302f68a32c8SEmil Constantinescu } 1303f68a32c8SEmil Constantinescu 1304d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRKGetOrder_RK(TS ts, PetscInt *order) 1305d71ae5a4SJacob Faibussowitsch { 13062ea87ba9SLisandro Dalcin TS_RK *rk = (TS_RK *)ts->data; 13072ea87ba9SLisandro Dalcin 13082ea87ba9SLisandro Dalcin PetscFunctionBegin; 13092ea87ba9SLisandro Dalcin *order = rk->tableau->order; 13103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13112ea87ba9SLisandro Dalcin } 13122ea87ba9SLisandro Dalcin 1313d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRKGetType_RK(TS ts, TSRKType *rktype) 1314d71ae5a4SJacob Faibussowitsch { 1315f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK *)ts->data; 1316f68a32c8SEmil Constantinescu 1317f68a32c8SEmil Constantinescu PetscFunctionBegin; 1318f68a32c8SEmil Constantinescu *rktype = rk->tableau->name; 13193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1320f68a32c8SEmil Constantinescu } 13212c9cb062Sluzhanghpp 1322d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRKSetType_RK(TS ts, TSRKType rktype) 1323d71ae5a4SJacob Faibussowitsch { 1324f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK *)ts->data; 1325f68a32c8SEmil Constantinescu PetscBool match; 1326f68a32c8SEmil Constantinescu RKTableauLink link; 1327f68a32c8SEmil Constantinescu 1328f68a32c8SEmil Constantinescu PetscFunctionBegin; 1329f68a32c8SEmil Constantinescu if (rk->tableau) { 13309566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(rk->tableau->name, rktype, &match)); 13313ba16761SJacob Faibussowitsch if (match) PetscFunctionReturn(PETSC_SUCCESS); 1332f68a32c8SEmil Constantinescu } 1333f68a32c8SEmil Constantinescu for (link = RKTableauList; link; link = link->next) { 13349566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(link->tab.name, rktype, &match)); 1335f68a32c8SEmil Constantinescu if (match) { 13369566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSRKTableauReset(ts)); 1337f68a32c8SEmil Constantinescu rk->tableau = &link->tab; 13389566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSRKTableauSetUp(ts)); 13392ffb9264SLisandro Dalcin ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 13403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1341f68a32c8SEmil Constantinescu } 1342f68a32c8SEmil Constantinescu } 134398921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", rktype); 1344e4dd521cSBarry Smith } 1345e4dd521cSBarry Smith 1346d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGetStages_RK(TS ts, PetscInt *ns, Vec **Y) 1347d71ae5a4SJacob Faibussowitsch { 1348ff22ae23SHong Zhang TS_RK *rk = (TS_RK *)ts->data; 1349ff22ae23SHong Zhang 1350ff22ae23SHong Zhang PetscFunctionBegin; 13510429704eSStefano Zampini if (ns) *ns = rk->tableau->s; 1352d2daff3dSHong Zhang if (Y) *Y = rk->Y; 13533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1354ff22ae23SHong Zhang } 1355ff22ae23SHong Zhang 1356d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_RK(TS ts) 1357d71ae5a4SJacob Faibussowitsch { 1358b3a6b972SJed Brown PetscFunctionBegin; 13599566063dSJacob Faibussowitsch PetscCall(TSReset_RK(ts)); 1360b3a6b972SJed Brown if (ts->dm) { 13619566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSRK, DMRestrictHook_TSRK, ts)); 13629566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSRK, DMSubDomainRestrictHook_TSRK, ts)); 1363b3a6b972SJed Brown } 13649566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 13659566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetOrder_C", NULL)); 13669566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetType_C", NULL)); 13679566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKSetType_C", NULL)); 13689566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetTableau_C", NULL)); 13699566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKSetMultirate_C", NULL)); 13709566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetMultirate_C", NULL)); 13712e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSetUp_RK_MultirateSplit_C", NULL)); 13722e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSReset_RK_MultirateSplit_C", NULL)); 13732e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSetUp_RK_MultirateNonsplit_C", NULL)); 13742e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSReset_RK_MultirateNonsplit_C", NULL)); 13753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1376b3a6b972SJed Brown } 1377ff22ae23SHong Zhang 13783ca0882eSHong Zhang /* 13793ca0882eSHong Zhang This defines the nonlinear equation that is to be solved with SNES 13803ca0882eSHong Zhang We do not need to solve the equation; we just use SNES to approximate the Jacobian 13813ca0882eSHong Zhang */ 1382d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_RK(SNES snes, Vec x, Vec y, TS ts) 1383d71ae5a4SJacob Faibussowitsch { 13843ca0882eSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 13853ca0882eSHong Zhang DM dm, dmsave; 13863ca0882eSHong Zhang 13873ca0882eSHong Zhang PetscFunctionBegin; 13889566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 13893ca0882eSHong Zhang /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 13903ca0882eSHong Zhang dmsave = ts->dm; 13913ca0882eSHong Zhang ts->dm = dm; 13929566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, rk->stage_time, x, y)); 13933ca0882eSHong Zhang ts->dm = dmsave; 13943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13953ca0882eSHong Zhang } 13963ca0882eSHong Zhang 1397d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_RK(SNES snes, Vec x, Mat A, Mat B, TS ts) 1398d71ae5a4SJacob Faibussowitsch { 13993ca0882eSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 14003ca0882eSHong Zhang DM dm, dmsave; 14013ca0882eSHong Zhang 14023ca0882eSHong Zhang PetscFunctionBegin; 14039566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 14043ca0882eSHong Zhang dmsave = ts->dm; 14053ca0882eSHong Zhang ts->dm = dm; 14069566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobian(ts, rk->stage_time, x, A, B)); 14073ca0882eSHong Zhang ts->dm = dmsave; 14083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14093ca0882eSHong Zhang } 14103ca0882eSHong Zhang 1411*cc4c1da9SBarry Smith /*@ 1412bcf0153eSBarry Smith TSRKSetMultirate - Use the interpolation-based multirate `TSRK` method 141321052c0cSHong Zhang 141420f4b53cSBarry Smith Logically Collective 141521052c0cSHong Zhang 1416d8d19677SJose E. Roman Input Parameters: 141721052c0cSHong Zhang + ts - timestepping context 1418bcf0153eSBarry Smith - use_multirate - `PETSC_TRUE` enables the multirate `TSRK` method, sets the basic method to be RK2A and sets the ratio between slow stepsize and fast stepsize to be 2 141921052c0cSHong Zhang 1420bcf0153eSBarry Smith Options Database Key: 142121052c0cSHong Zhang . -ts_rk_multirate - <true,false> 142221052c0cSHong Zhang 142321052c0cSHong Zhang Level: intermediate 142421052c0cSHong Zhang 1425bcf0153eSBarry Smith Note: 1426da81f932SPierre Jolivet The multirate method requires interpolation. The default interpolation works for 1st- and 2nd- order RK, but not for high-order RKs except `TSRK5DP` which comes with the interpolation coefficients (binterp). 1427bcf0153eSBarry Smith 14281cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKGetMultirate()` 142921052c0cSHong Zhang @*/ 1430d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKSetMultirate(TS ts, PetscBool use_multirate) 1431d71ae5a4SJacob Faibussowitsch { 14328a4be4dbSHong Zhang PetscFunctionBegin; 1433cac4c232SBarry Smith PetscTryMethod(ts, "TSRKSetMultirate_C", (TS, PetscBool), (ts, use_multirate)); 14343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 143521052c0cSHong Zhang } 143621052c0cSHong Zhang 1437*cc4c1da9SBarry Smith /*@ 1438bcf0153eSBarry Smith TSRKGetMultirate - Gets whether to use the interpolation-based multirate `TSRK` method 143921052c0cSHong Zhang 144020f4b53cSBarry Smith Not Collective 144121052c0cSHong Zhang 144221052c0cSHong Zhang Input Parameter: 144321052c0cSHong Zhang . ts - timestepping context 144421052c0cSHong Zhang 144521052c0cSHong Zhang Output Parameter: 1446bcf0153eSBarry Smith . use_multirate - `PETSC_TRUE` if the multirate RK method is enabled, `PETSC_FALSE` otherwise 144721052c0cSHong Zhang 144821052c0cSHong Zhang Level: intermediate 144921052c0cSHong Zhang 14501cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKSetMultirate()` 145121052c0cSHong Zhang @*/ 1452d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKGetMultirate(TS ts, PetscBool *use_multirate) 1453d71ae5a4SJacob Faibussowitsch { 14547dbacdf2SHong Zhang PetscFunctionBegin; 1455cac4c232SBarry Smith PetscUseMethod(ts, "TSRKGetMultirate_C", (TS, PetscBool *), (ts, use_multirate)); 14563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 145721052c0cSHong Zhang } 145821052c0cSHong Zhang 1459ebe3b25bSBarry Smith /*MC 1460f68a32c8SEmil Constantinescu TSRK - ODE and DAE solver using Runge-Kutta schemes 1461ebe3b25bSBarry Smith 1462dd8e379bSPierre Jolivet The user should provide the right-hand side of the equation 1463bcf0153eSBarry Smith using `TSSetRHSFunction()`. 1464ebe3b25bSBarry Smith 1465d41222bbSBarry Smith Level: beginner 1466d41222bbSBarry Smith 1467bcf0153eSBarry Smith Notes: 1468bcf0153eSBarry Smith The default is `TSRK3BS`, it can be changed with `TSRKSetType()` or -ts_rk_type 1469ebe3b25bSBarry Smith 14701cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSRK`, `TSSetType()`, `TSRKSetType()`, `TSRKGetType()`, `TSRK2D`, `TSRK2E`, `TSRK3`, 1471bcf0153eSBarry Smith `TSRK4`, `TSRK5`, `TSRKPRSSP2`, `TSRKBPR3`, `TSRKType`, `TSRKRegister()`, `TSRKSetMultirate()`, `TSRKGetMultirate()`, `TSType` 1472ebe3b25bSBarry Smith M*/ 1473d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts) 1474d71ae5a4SJacob Faibussowitsch { 14751566a47fSLisandro Dalcin TS_RK *rk; 1476e4dd521cSBarry Smith 1477e4dd521cSBarry Smith PetscFunctionBegin; 14789566063dSJacob Faibussowitsch PetscCall(TSRKInitializePackage()); 1479f68a32c8SEmil Constantinescu 1480f68a32c8SEmil Constantinescu ts->ops->reset = TSReset_RK; 14815f70b478SJed Brown ts->ops->destroy = TSDestroy_RK; 14825f70b478SJed Brown ts->ops->view = TSView_RK; 1483f68a32c8SEmil Constantinescu ts->ops->load = TSLoad_RK; 1484f68a32c8SEmil Constantinescu ts->ops->setup = TSSetUp_RK; 1485f68a32c8SEmil Constantinescu ts->ops->interpolate = TSInterpolate_RK; 14862c9cb062Sluzhanghpp ts->ops->step = TSStep_RK; 1487f68a32c8SEmil Constantinescu ts->ops->evaluatestep = TSEvaluateStep_RK; 1488fecfb714SLisandro Dalcin ts->ops->rollback = TSRollBack_RK; 1489f68a32c8SEmil Constantinescu ts->ops->setfromoptions = TSSetFromOptions_RK; 1490ff22ae23SHong Zhang ts->ops->getstages = TSGetStages_RK; 1491e4dd521cSBarry Smith 14923ca0882eSHong Zhang ts->ops->snesfunction = SNESTSFormFunction_RK; 14933ca0882eSHong Zhang ts->ops->snesjacobian = SNESTSFormJacobian_RK; 14940f7a1166SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_RK; 149513af1a74SHong Zhang ts->ops->adjointsetup = TSAdjointSetUp_RK; 149613af1a74SHong Zhang ts->ops->adjointstep = TSAdjointStep_RK; 149713af1a74SHong Zhang ts->ops->adjointreset = TSAdjointReset_RK; 14980f7a1166SHong Zhang 149913af1a74SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_RK; 1500922a638cSHong Zhang ts->ops->forwardsetup = TSForwardSetUp_RK; 1501922a638cSHong Zhang ts->ops->forwardreset = TSForwardReset_RK; 1502922a638cSHong Zhang ts->ops->forwardstep = TSForwardStep_RK; 1503922a638cSHong Zhang ts->ops->forwardgetstages = TSForwardGetStages_RK; 1504922a638cSHong Zhang 15054dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&rk)); 15061566a47fSLisandro Dalcin ts->data = (void *)rk; 1507e4dd521cSBarry Smith 15089566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetOrder_C", TSRKGetOrder_RK)); 15099566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetType_C", TSRKGetType_RK)); 15109566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKSetType_C", TSRKSetType_RK)); 15119566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetTableau_C", TSRKGetTableau_RK)); 15129566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKSetMultirate_C", TSRKSetMultirate_RK)); 15139566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetMultirate_C", TSRKGetMultirate_RK)); 1514be5899b3SLisandro Dalcin 15159566063dSJacob Faibussowitsch PetscCall(TSRKSetType(ts, TSRKDefault)); 151690b67129SHong Zhang rk->dtratio = 1; 15173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1518e4dd521cSBarry Smith } 1519