1e4dd521cSBarry Smith /* 2474dd773SHong Zhang Code for time stepping with the Runge-Kutta method 3f68a32c8SEmil Constantinescu 4f68a32c8SEmil Constantinescu Notes: 5474dd773SHong Zhang The general system is written as 6474dd773SHong Zhang 7474dd773SHong Zhang Udot = F(t,U) 8474dd773SHong Zhang 9e4dd521cSBarry Smith */ 102c9cb062Sluzhanghpp 11af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 12f68a32c8SEmil Constantinescu #include <petscdm.h> 13474dd773SHong Zhang #include <../src/ts/impls/explicit/rk/rk.h> 1421052c0cSHong Zhang #include <../src/ts/impls/explicit/rk/mrk.h> 15f68a32c8SEmil Constantinescu 16484bcdc7SDebojyoti Ghosh static TSRKType TSRKDefault = TSRK3BS; 17f68a32c8SEmil Constantinescu static PetscBool TSRKRegisterAllCalled; 18f68a32c8SEmil Constantinescu static PetscBool TSRKPackageInitialized; 19f68a32c8SEmil Constantinescu 20ab8985f5SHong Zhang static RKTableauLink RKTableauList; 21ab8985f5SHong Zhang 22f68a32c8SEmil Constantinescu /*MC 23287c2655SBarry Smith TSRK1FE - First order forward Euler scheme. 24262ff7bbSBarry Smith 25f68a32c8SEmil Constantinescu This method has one stage. 26f68a32c8SEmil Constantinescu 27bcf0153eSBarry Smith Options Database Key: 2867b8a455SSatish Balay . -ts_rk_type 1fe - use type 1fe 29287c2655SBarry Smith 30f68a32c8SEmil Constantinescu Level: advanced 31f68a32c8SEmil Constantinescu 321cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 33f68a32c8SEmil Constantinescu M*/ 34f68a32c8SEmil Constantinescu /*MC 35e7685601SHong Zhang TSRK2A - Second order RK scheme (Heun's method). 36f68a32c8SEmil Constantinescu 37f68a32c8SEmil Constantinescu This method has two stages. 38f68a32c8SEmil Constantinescu 39bcf0153eSBarry Smith Options Database Key: 4067b8a455SSatish Balay . -ts_rk_type 2a - use type 2a 41287c2655SBarry Smith 42f68a32c8SEmil Constantinescu Level: advanced 43f68a32c8SEmil Constantinescu 441cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 45f68a32c8SEmil Constantinescu M*/ 46f68a32c8SEmil Constantinescu /*MC 47e7685601SHong Zhang TSRK2B - Second order RK scheme (the midpoint method). 48e7685601SHong Zhang 49e7685601SHong Zhang This method has two stages. 50e7685601SHong Zhang 51bcf0153eSBarry Smith Options Database Key: 5267b8a455SSatish Balay . -ts_rk_type 2b - use type 2b 53e7685601SHong Zhang 54e7685601SHong Zhang Level: advanced 55e7685601SHong Zhang 561cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 57e7685601SHong Zhang M*/ 58e7685601SHong Zhang /*MC 59f68a32c8SEmil Constantinescu TSRK3 - Third order RK scheme. 60f68a32c8SEmil Constantinescu 61f68a32c8SEmil Constantinescu This method has three stages. 62f68a32c8SEmil Constantinescu 63bcf0153eSBarry Smith Options Database Key: 6467b8a455SSatish Balay . -ts_rk_type 3 - use type 3 65287c2655SBarry Smith 66f68a32c8SEmil Constantinescu Level: advanced 67f68a32c8SEmil Constantinescu 681cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 69f68a32c8SEmil Constantinescu M*/ 70f68a32c8SEmil Constantinescu /*MC 712109b73fSDebojyoti Ghosh TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method. 722109b73fSDebojyoti Ghosh 73d1905564SLisandro Dalcin This method has four stages with the First Same As Last (FSAL) property. 742109b73fSDebojyoti Ghosh 75bcf0153eSBarry Smith Options Database Key: 7667b8a455SSatish Balay . -ts_rk_type 3bs - use type 3bs 77287c2655SBarry Smith 782109b73fSDebojyoti Ghosh Level: advanced 792109b73fSDebojyoti Ghosh 80606c0280SSatish Balay References: 81606c0280SSatish Balay . * - https://doi.org/10.1016/0893-9659(89)90079-7 82d1905564SLisandro Dalcin 831cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 842109b73fSDebojyoti Ghosh M*/ 852109b73fSDebojyoti Ghosh /*MC 86f68a32c8SEmil Constantinescu TSRK4 - Fourth order RK scheme. 87f68a32c8SEmil Constantinescu 882e698f8bSDebojyoti Ghosh This is the classical Runge-Kutta method with four stages. 89f68a32c8SEmil Constantinescu 90bcf0153eSBarry Smith Options Database Key: 9167b8a455SSatish Balay . -ts_rk_type 4 - use type 4 92287c2655SBarry Smith 93f68a32c8SEmil Constantinescu Level: advanced 94f68a32c8SEmil Constantinescu 951cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 96f68a32c8SEmil Constantinescu M*/ 97f68a32c8SEmil Constantinescu /*MC 982e698f8bSDebojyoti Ghosh TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method. 99f68a32c8SEmil Constantinescu 100f68a32c8SEmil Constantinescu This method has six stages. 101f68a32c8SEmil Constantinescu 102bcf0153eSBarry Smith Options Database Key: 10367b8a455SSatish Balay . -ts_rk_type 5f - use type 5f 104287c2655SBarry Smith 105f68a32c8SEmil Constantinescu Level: advanced 106f68a32c8SEmil Constantinescu 1071cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 108f68a32c8SEmil Constantinescu M*/ 1092109b73fSDebojyoti Ghosh /*MC 1102e698f8bSDebojyoti Ghosh TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method. 1112109b73fSDebojyoti Ghosh 112d1905564SLisandro Dalcin This method has seven stages with the First Same As Last (FSAL) property. 1132109b73fSDebojyoti Ghosh 114bcf0153eSBarry Smith Options Database Key: 11567b8a455SSatish Balay . -ts_rk_type 5dp - use type 5dp 116287c2655SBarry Smith 1172109b73fSDebojyoti Ghosh Level: advanced 1182109b73fSDebojyoti Ghosh 119606c0280SSatish Balay References: 120606c0280SSatish Balay . * - https://doi.org/10.1016/0771-050X(80)90013-3 121d1905564SLisandro Dalcin 1221cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 1232109b73fSDebojyoti Ghosh M*/ 12405e23783SLisandro Dalcin /*MC 12505e23783SLisandro Dalcin TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method. 12605e23783SLisandro Dalcin 12705e23783SLisandro Dalcin This method has eight stages with the First Same As Last (FSAL) property. 12805e23783SLisandro Dalcin 129bcf0153eSBarry Smith Options Database Key: 13067b8a455SSatish Balay . -ts_rk_type 5bs - use type 5bs 131287c2655SBarry Smith 13205e23783SLisandro Dalcin Level: advanced 13305e23783SLisandro Dalcin 134606c0280SSatish Balay References: 135606c0280SSatish Balay . * - https://doi.org/10.1016/0898-1221(96)00141-1 13605e23783SLisandro Dalcin 1371cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 13805e23783SLisandro Dalcin M*/ 13937acaa02SHendrik Ranocha /*MC 14037acaa02SHendrik Ranocha TSRK6VR - Sixth order robust Verner RK scheme with fifth order embedded method. 14137acaa02SHendrik Ranocha 14237acaa02SHendrik Ranocha This method has nine stages with the First Same As Last (FSAL) property. 14337acaa02SHendrik Ranocha 144bcf0153eSBarry Smith Options Database Key: 14567b8a455SSatish Balay . -ts_rk_type 6vr - use type 6vr 14637acaa02SHendrik Ranocha 14737acaa02SHendrik Ranocha Level: advanced 14837acaa02SHendrik Ranocha 149606c0280SSatish Balay References: 150606c0280SSatish Balay . * - http://people.math.sfu.ca/~jverner/RKV65.IIIXb.Robust.00010102836.081204.CoeffsOnlyRAT 15137acaa02SHendrik Ranocha 1521cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 15337acaa02SHendrik Ranocha M*/ 15437acaa02SHendrik Ranocha /*MC 15537acaa02SHendrik Ranocha TSRK7VR - Seventh order robust Verner RK scheme with sixth order embedded method. 15637acaa02SHendrik Ranocha 1578f3d7ee2SCarsten Uphoff This method has ten stages. 15837acaa02SHendrik Ranocha 159bcf0153eSBarry Smith Options Database Key: 16067b8a455SSatish Balay . -ts_rk_type 7vr - use type 7vr 16137acaa02SHendrik Ranocha 16237acaa02SHendrik Ranocha Level: advanced 16337acaa02SHendrik Ranocha 164606c0280SSatish Balay References: 165606c0280SSatish Balay . * - http://people.math.sfu.ca/~jverner/RKV76.IIa.Robust.000027015646.081206.CoeffsOnlyRAT 16637acaa02SHendrik Ranocha 1671cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 16837acaa02SHendrik Ranocha M*/ 16937acaa02SHendrik Ranocha /*MC 170d5b43468SJose E. Roman TSRK8VR - Eighth order robust Verner RK scheme with seventh order embedded method. 17137acaa02SHendrik Ranocha 1728f3d7ee2SCarsten Uphoff This method has thirteen stages. 17337acaa02SHendrik Ranocha 174bcf0153eSBarry Smith Options Database Key: 17567b8a455SSatish Balay . -ts_rk_type 8vr - use type 8vr 17637acaa02SHendrik Ranocha 17737acaa02SHendrik Ranocha Level: advanced 17837acaa02SHendrik Ranocha 179606c0280SSatish Balay References: 180606c0280SSatish Balay . * - http://people.math.sfu.ca/~jverner/RKV87.IIa.Robust.00000754677.081208.CoeffsOnlyRATandFLOAT 18137acaa02SHendrik Ranocha 1821cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKType`, `TSRKSetType()` 18337acaa02SHendrik Ranocha M*/ 184262ff7bbSBarry Smith 185f68a32c8SEmil Constantinescu /*@C 186bcf0153eSBarry Smith TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in `TSRK` 187262ff7bbSBarry Smith 188f68a32c8SEmil Constantinescu Not Collective, but should be called by all processes which will need the schemes to be registered 189262ff7bbSBarry Smith 190f68a32c8SEmil Constantinescu Level: advanced 191262ff7bbSBarry Smith 1921cc06b55SBarry Smith .seealso: [](ch_ts), `TSRKRegisterDestroy()`, `TSRKRegister()` 193262ff7bbSBarry Smith @*/ 194d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKRegisterAll(void) 195d71ae5a4SJacob Faibussowitsch { 196262ff7bbSBarry Smith PetscFunctionBegin; 1973ba16761SJacob Faibussowitsch if (TSRKRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS); 198f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_TRUE; 199e4dd521cSBarry Smith 2004e82626cSLisandro Dalcin #define RC PetscRealConstant 201e4dd521cSBarry Smith { 2029371c9d4SSatish Balay const PetscReal A[1][1] = {{0}}, b[1] = {RC(1.0)}; 2039566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK1FE, 1, 1, &A[0][0], b, NULL, NULL, 0, NULL)); 204e4dd521cSBarry Smith } 205db046809SBarry Smith { 2069371c9d4SSatish Balay const PetscReal A[2][2] = 2079371c9d4SSatish Balay { 2089371c9d4SSatish Balay {0, 0}, 2099371c9d4SSatish Balay {RC(1.0), 0} 2109371c9d4SSatish Balay }, 2119371c9d4SSatish Balay b[2] = {RC(0.5), RC(0.5)}, bembed[2] = {RC(1.0), 0}; 2129566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK2A, 2, 2, &A[0][0], b, NULL, bembed, 0, NULL)); 213db046809SBarry Smith } 214f68a32c8SEmil Constantinescu { 2159371c9d4SSatish Balay const PetscReal A[2][2] = 2169371c9d4SSatish Balay { 2179371c9d4SSatish Balay {0, 0}, 2189371c9d4SSatish Balay {RC(0.5), 0} 2199371c9d4SSatish Balay }, 2209958aef7SJacob Faibussowitsch b[2] = {0, RC(1.0)}; 2219566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK2B, 2, 2, &A[0][0], b, NULL, NULL, 0, NULL)); 222e7685601SHong Zhang } 223e7685601SHong Zhang { 2249371c9d4SSatish Balay const PetscReal A[3][3] = 2259371c9d4SSatish Balay { 2269371c9d4SSatish Balay {0, 0, 0}, 2274e82626cSLisandro Dalcin {RC(2.0) / RC(3.0), 0, 0}, 2289371c9d4SSatish Balay {RC(-1.0) / RC(3.0), RC(1.0), 0} 2299371c9d4SSatish Balay }, 2304e82626cSLisandro Dalcin b[3] = {RC(0.25), RC(0.5), RC(0.25)}; 2319566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK3, 3, 3, &A[0][0], b, NULL, NULL, 0, NULL)); 232fdefd5e5SDebojyoti Ghosh } 233fdefd5e5SDebojyoti Ghosh { 2349371c9d4SSatish Balay const PetscReal A[4][4] = 2359371c9d4SSatish Balay { 2369371c9d4SSatish Balay {0, 0, 0, 0}, 2374e82626cSLisandro Dalcin {RC(1.0) / RC(2.0), 0, 0, 0}, 2384e82626cSLisandro Dalcin {0, RC(3.0) / RC(4.0), 0, 0}, 2399371c9d4SSatish Balay {RC(2.0) / RC(9.0), RC(1.0) / RC(3.0), RC(4.0) / RC(9.0), 0} 2409371c9d4SSatish Balay }, 2419371c9d4SSatish Balay b[4] = {RC(2.0) / RC(9.0), RC(1.0) / RC(3.0), RC(4.0) / RC(9.0), 0}, bembed[4] = {RC(7.0) / RC(24.0), RC(1.0) / RC(4.0), RC(1.0) / RC(3.0), RC(1.0) / RC(8.0)}; 2429566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK3BS, 3, 4, &A[0][0], b, NULL, bembed, 0, NULL)); 243db046809SBarry Smith } 244f68a32c8SEmil Constantinescu { 2459371c9d4SSatish Balay const PetscReal A[4][4] = 2469371c9d4SSatish Balay { 2479371c9d4SSatish Balay {0, 0, 0, 0}, 2484e82626cSLisandro Dalcin {RC(0.5), 0, 0, 0}, 2494e82626cSLisandro Dalcin {0, RC(0.5), 0, 0}, 2509371c9d4SSatish Balay {0, 0, RC(1.0), 0} 2519371c9d4SSatish Balay }, 2524e82626cSLisandro Dalcin b[4] = {RC(1.0) / RC(6.0), RC(1.0) / RC(3.0), RC(1.0) / RC(3.0), RC(1.0) / RC(6.0)}; 2539566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK4, 4, 4, &A[0][0], b, NULL, NULL, 0, NULL)); 254db046809SBarry Smith } 255f68a32c8SEmil Constantinescu { 2569371c9d4SSatish Balay const PetscReal A[6][6] = 2579371c9d4SSatish Balay { 2589371c9d4SSatish Balay {0, 0, 0, 0, 0, 0}, 2594e82626cSLisandro Dalcin {RC(0.25), 0, 0, 0, 0, 0}, 2604e82626cSLisandro Dalcin {RC(3.0) / RC(32.0), RC(9.0) / RC(32.0), 0, 0, 0, 0}, 2614e82626cSLisandro Dalcin {RC(1932.0) / RC(2197.0), RC(-7200.0) / RC(2197.0), RC(7296.0) / RC(2197.0), 0, 0, 0}, 2624e82626cSLisandro Dalcin {RC(439.0) / RC(216.0), RC(-8.0), RC(3680.0) / RC(513.0), RC(-845.0) / RC(4104.0), 0, 0}, 2639371c9d4SSatish Balay {RC(-8.0) / RC(27.0), RC(2.0), RC(-3544.0) / RC(2565.0), RC(1859.0) / RC(4104.0), RC(-11.0) / RC(40.0), 0} 2649371c9d4SSatish Balay }, 2654e82626cSLisandro Dalcin b[6] = {RC(16.0) / RC(135.0), 0, RC(6656.0) / RC(12825.0), RC(28561.0) / RC(56430.0), RC(-9.0) / RC(50.0), RC(2.0) / RC(55.0)}, 2664e82626cSLisandro Dalcin bembed[6] = {RC(25.0) / RC(216.0), 0, RC(1408.0) / RC(2565.0), RC(2197.0) / RC(4104.0), RC(-1.0) / RC(5.0), 0}; 2679566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK5F, 5, 6, &A[0][0], b, NULL, bembed, 0, NULL)); 268fdefd5e5SDebojyoti Ghosh } 269fdefd5e5SDebojyoti Ghosh { 2709371c9d4SSatish Balay const PetscReal A[7][7] = 2719371c9d4SSatish Balay { 2729371c9d4SSatish Balay {0, 0, 0, 0, 0, 0, 0}, 2734e82626cSLisandro Dalcin {RC(1.0) / RC(5.0), 0, 0, 0, 0, 0, 0}, 2744e82626cSLisandro Dalcin {RC(3.0) / RC(40.0), RC(9.0) / RC(40.0), 0, 0, 0, 0, 0}, 2754e82626cSLisandro Dalcin {RC(44.0) / RC(45.0), RC(-56.0) / RC(15.0), RC(32.0) / RC(9.0), 0, 0, 0, 0}, 2764e82626cSLisandro Dalcin {RC(19372.0) / RC(6561.0), RC(-25360.0) / RC(2187.0), RC(64448.0) / RC(6561.0), RC(-212.0) / RC(729.0), 0, 0, 0}, 2774e82626cSLisandro Dalcin {RC(9017.0) / RC(3168.0), RC(-355.0) / RC(33.0), RC(46732.0) / RC(5247.0), RC(49.0) / RC(176.0), RC(-5103.0) / RC(18656.0), 0, 0}, 2789371c9d4SSatish Balay {RC(35.0) / RC(384.0), 0, RC(500.0) / RC(1113.0), RC(125.0) / RC(192.0), RC(-2187.0) / RC(6784.0), RC(11.0) / RC(84.0), 0} 2799371c9d4SSatish Balay }, 2804e82626cSLisandro Dalcin b[7] = {RC(35.0) / RC(384.0), 0, RC(500.0) / RC(1113.0), RC(125.0) / RC(192.0), RC(-2187.0) / RC(6784.0), RC(11.0) / RC(84.0), 0}, 2819371c9d4SSatish Balay bembed[7] = {RC(5179.0) / RC(57600.0), 0, RC(7571.0) / RC(16695.0), RC(393.0) / RC(640.0), RC(-92097.0) / RC(339200.0), RC(187.0) / RC(2100.0), RC(1.0) / RC(40.0)}, binterp[7][5] = {{RC(1.0), RC(-4034104133.0) / RC(1410260304.0), RC(105330401.0) / RC(33982176.0), RC(-13107642775.0) / RC(11282082432.0), RC(6542295.0) / RC(470086768.0)}, {0, 0, 0, 0, 0}, {0, RC(132343189600.0) / RC(32700410799.0), RC(-833316000.0) / RC(131326951.0), RC(91412856700.0) / RC(32700410799.0), RC(-523383600.0) / RC(10900136933.0)}, {0, RC(-115792950.0) / RC(29380423.0), RC(185270875.0) / RC(16991088.0), RC(-12653452475.0) / RC(1880347072.0), RC(98134425.0) / RC(235043384.0)}, {0, RC(70805911779.0) / RC(24914598704.0), RC(-4531260609.0) / RC(600351776.0), RC(988140236175.0) / RC(199316789632.0), RC(-14307999165.0) / RC(24914598704.0)}, {0, RC(-331320693.0) / RC(205662961.0), RC(31361737.0) / RC(7433601.0), RC(-2426908385.0) / RC(822651844.0), RC(97305120.0) / RC(205662961.0)}, {0, RC(44764047.0) / RC(29380423.0), RC(-1532549.0) / RC(353981.0), RC(90730570.0) / RC(29380423.0), RC(-8293050.0) / RC(29380423.0)}}; 2829566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK5DP, 5, 7, &A[0][0], b, NULL, bembed, 5, binterp[0])); 283f68a32c8SEmil Constantinescu } 28405e23783SLisandro Dalcin { 2859371c9d4SSatish Balay const PetscReal A[8][8] = 2869371c9d4SSatish Balay { 2879371c9d4SSatish Balay {0, 0, 0, 0, 0, 0, 0, 0}, 2884e82626cSLisandro Dalcin {RC(1.0) / RC(6.0), 0, 0, 0, 0, 0, 0, 0}, 2894e82626cSLisandro Dalcin {RC(2.0) / RC(27.0), RC(4.0) / RC(27.0), 0, 0, 0, 0, 0, 0}, 2904e82626cSLisandro Dalcin {RC(183.0) / RC(1372.0), RC(-162.0) / RC(343.0), RC(1053.0) / RC(1372.0), 0, 0, 0, 0, 0}, 2914e82626cSLisandro Dalcin {RC(68.0) / RC(297.0), RC(-4.0) / RC(11.0), RC(42.0) / RC(143.0), RC(1960.0) / RC(3861.0), 0, 0, 0, 0}, 2924e82626cSLisandro Dalcin {RC(597.0) / RC(22528.0), RC(81.0) / RC(352.0), RC(63099.0) / RC(585728.0), RC(58653.0) / RC(366080.0), RC(4617.0) / RC(20480.0), 0, 0, 0}, 2934e82626cSLisandro Dalcin {RC(174197.0) / RC(959244.0), RC(-30942.0) / RC(79937.0), RC(8152137.0) / RC(19744439.0), RC(666106.0) / RC(1039181.0), RC(-29421.0) / RC(29068.0), RC(482048.0) / RC(414219.0), 0, 0}, 2949371c9d4SSatish Balay {RC(587.0) / RC(8064.0), 0, RC(4440339.0) / RC(15491840.0), RC(24353.0) / RC(124800.0), RC(387.0) / RC(44800.0), RC(2152.0) / RC(5985.0), RC(7267.0) / RC(94080.0), 0} 2959371c9d4SSatish Balay }, 2964e82626cSLisandro Dalcin b[8] = {RC(587.0) / RC(8064.0), 0, RC(4440339.0) / RC(15491840.0), RC(24353.0) / RC(124800.0), RC(387.0) / RC(44800.0), RC(2152.0) / RC(5985.0), RC(7267.0) / RC(94080.0), 0}, 2974e82626cSLisandro Dalcin bembed[8] = {RC(2479.0) / RC(34992.0), 0, RC(123.0) / RC(416.0), RC(612941.0) / RC(3411720.0), RC(43.0) / RC(1440.0), RC(2272.0) / RC(6561.0), RC(79937.0) / RC(1113912.0), RC(3293.0) / RC(556956.0)}; 2989566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK5BS, 5, 8, &A[0][0], b, NULL, bembed, 0, NULL)); 29905e23783SLisandro Dalcin } 30037acaa02SHendrik Ranocha { 3019371c9d4SSatish Balay const PetscReal A[9][9] = 3029371c9d4SSatish Balay { 3039371c9d4SSatish Balay {0, 0, 0, 0, 0, 0, 0, 0, 0}, 30463fe90e8SHendrik Ranocha {RC(1.8000000000000000000000000000000000000000e-01), 0, 0, 0, 0, 0, 0, 0, 0}, 30563fe90e8SHendrik Ranocha {RC(8.9506172839506172839506172839506172839506e-02), RC(7.7160493827160493827160493827160493827160e-02), 0, 0, 0, 0, 0, 0, 0}, 30663fe90e8SHendrik Ranocha {RC(6.2500000000000000000000000000000000000000e-02), 0, RC(1.8750000000000000000000000000000000000000e-01), 0, 0, 0, 0, 0, 0}, 30763fe90e8SHendrik Ranocha {RC(3.1651600000000000000000000000000000000000e-01), 0, RC(-1.0449480000000000000000000000000000000000e+00), RC(1.2584320000000000000000000000000000000000e+00), 0, 0, 0, 0, 0}, 30863fe90e8SHendrik Ranocha {RC(2.7232612736485626257225065566674305502508e-01), 0, RC(-8.2513360323886639676113360323886639676113e-01), RC(1.0480917678812415654520917678812415654521e+00), RC(1.0471570799276856873679117969088177628396e-01), 0, 0, 0, 0}, 30963fe90e8SHendrik Ranocha {RC(-1.6699418599716514314329607278961797333198e-01), 0, RC(6.3170850202429149797570850202429149797571e-01), RC(1.7461044552773876082146758838488161796432e-01), RC(-1.0665356459086066122525194734018680677781e+00), RC(1.2272108843537414965986394557823129251701e+00), 0, 0, 0}, 31063fe90e8SHendrik Ranocha {RC(3.6423751686909581646423751686909581646424e-01), 0, RC(-2.0404858299595141700404858299595141700405e-01), RC(-3.4883737816068643136312309244640071707741e-01), RC(3.2619323032856867443333608747142581729048e+00), RC(-2.7551020408163265306122448979591836734694e+00), RC(6.8181818181818181818181818181818181818182e-01), 0, 0}, 3119371c9d4SSatish Balay {RC(7.6388888888888888888888888888888888888889e-02), 0, 0, RC(3.6940836940836940836940836940836940836941e-01), 0, RC(2.4801587301587301587301587301587301587302e-01), RC(2.3674242424242424242424242424242424242424e-01), RC(6.9444444444444444444444444444444444444444e-02), 0} 3129371c9d4SSatish Balay }, 3139371c9d4SSatish Balay b[9] = {RC(7.6388888888888888888888888888888888888889e-02), 0, 0, RC(3.6940836940836940836940836940836940836941e-01), 0, RC(2.4801587301587301587301587301587301587302e-01), RC(2.3674242424242424242424242424242424242424e-01), 3149371c9d4SSatish Balay RC(6.9444444444444444444444444444444444444444e-02), 0}, 31563fe90e8SHendrik Ranocha bembed[9] = {RC(5.8700209643605870020964360587002096436059e-02), 0, 0, RC(4.8072562358276643990929705215419501133787e-01), RC(-8.5341242076919085578832094861228313083563e-01), RC(1.2046485260770975056689342403628117913832e+00), 0, RC(-5.9242373072160306202859394348756050883710e-02), RC(1.6858043453788134639198468985703028256220e-01)}; 3169566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK6VR, 6, 9, &A[0][0], b, NULL, bembed, 0, NULL)); 31737acaa02SHendrik Ranocha } 31837acaa02SHendrik Ranocha { 3199371c9d4SSatish Balay const PetscReal A[10][10] = 3209371c9d4SSatish Balay { 3219371c9d4SSatish Balay {0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, 32263fe90e8SHendrik Ranocha {RC(5.0000000000000000000000000000000000000000e-03), 0, 0, 0, 0, 0, 0, 0, 0, 0}, 32363fe90e8SHendrik Ranocha {RC(-1.0767901234567901234567901234567901234568e+00), RC(1.1856790123456790123456790123456790123457e+00), 0, 0, 0, 0, 0, 0, 0, 0}, 32463fe90e8SHendrik Ranocha {RC(4.0833333333333333333333333333333333333333e-02), 0, RC(1.2250000000000000000000000000000000000000e-01), 0, 0, 0, 0, 0, 0, 0}, 32563fe90e8SHendrik Ranocha {RC(6.3607142857142857142857142857142857142857e-01), 0, RC(-2.4444642857142857142857142857142857142857e+00), RC(2.2633928571428571428571428571428571428571e+00), 0, 0, 0, 0, 0, 0}, 32663fe90e8SHendrik Ranocha {RC(-2.5351211079349245229256383554660215487207e+00), 0, RC(1.0299374654449267920438514460756024913612e+01), RC(-7.9513032885990579949493217458266876536482e+00), RC(7.9301148923100592201226014271115261823800e-01), 0, 0, 0, 0, 0}, 32763fe90e8SHendrik Ranocha {RC(1.0018765812524632961969196583094999808207e+00), 0, RC(-4.1665712824423798331313938005470971453189e+00), RC(3.8343432929128642412552665218251378665197e+00), RC(-5.0233333560710847547464330228611765612403e-01), RC(6.6768474388416077115385092269857695410259e-01), 0, 0, 0, 0}, 32863fe90e8SHendrik Ranocha {RC(2.7255018354630767130333963819175005717348e+01), 0, RC(-4.2004617278410638355318645443909295369611e+01), RC(-1.0535713126619489917921081600546526103722e+01), RC(8.0495536711411937147983652158926826634202e+01), RC(-6.7343882271790513468549075963212975640927e+01), RC(1.3048657610777937463471187029566964762710e+01), 0, 0, 0}, 32963fe90e8SHendrik Ranocha {RC(-3.0397378057114965146943658658755763226883e+00), 0, RC(1.0138161410329801111857946190709700150441e+01), RC(-6.4293056748647215721462825629555298064437e+00), RC(-1.5864371483408276587115312853798610579467e+00), RC(1.8921781841968424410864308909131353365021e+00), RC(1.9699335407608869061292360163336442838006e-02), RC(5.4416989827933235465102724247952572977903e-03), 0, 0}, 3309371c9d4SSatish Balay {RC(-1.4449518916777735137351003179355712360517e+00), 0, RC(8.0318913859955919224117033223019560435041e+00), RC(-7.5831741663401346820798883023671588604984e+00), RC(3.5816169353190074211247685442452878696855e+00), RC(-2.4369722632199529411183809065693752383733e+00), RC(8.5158999992326179339689766032486142173390e-01), 0, 0, 0} 3319371c9d4SSatish Balay }, 3329371c9d4SSatish Balay b[10] = {RC(4.7425837833706756083569172717574534698932e-02), 0, 0, RC(2.5622361659370562659961727458274623448160e-01), RC(2.6951376833074206619473817258075952886764e-01), RC(1.2686622409092782845989138364739173247882e-01), RC(2.4887225942060071622046449427647492767292e-01), RC(3.0744837408200631335304388479099184768645e-03), RC(4.8023809989496943308189063347143123323209e-02), 0}, bembed[10] = {RC(4.7485247699299631037531273805727961552268e-02), 0, 0, RC(2.5599412588690633297154918245905393870497e-01), RC(2.7058478081067688722530891099268135732387e-01), RC(1.2505618684425992913638822323746917920448e-01), 3339371c9d4SSatish Balay RC(2.5204468723743860507184043820197442562182e-01), 0, 0, RC(4.8834971521418614557381971303093137592592e-02)}; 3349566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK7VR, 7, 10, &A[0][0], b, NULL, bembed, 0, NULL)); 33537acaa02SHendrik Ranocha } 33637acaa02SHendrik Ranocha { 3379371c9d4SSatish Balay const PetscReal A[13][13] = 3389371c9d4SSatish Balay { 3399371c9d4SSatish Balay {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, 34063fe90e8SHendrik Ranocha {RC(2.5000000000000000000000000000000000000000e-01), 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, 34163fe90e8SHendrik Ranocha {RC(8.7400846504915232052686327594877411977046e-02), RC(2.5487604938654321753087950620345685135815e-02), 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, 34263fe90e8SHendrik Ranocha {RC(4.2333169291338582677165354330708661417323e-02), 0, RC(1.2699950787401574803149606299212598425197e-01), 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, 34363fe90e8SHendrik Ranocha {RC(4.2609505888742261494881445237572274090942e-01), 0, RC(-1.5987952846591523265427733230657181117089e+00), RC(1.5967002257717297115939588706899953707994e+00), 0, 0, 0, 0, 0, 0, 0, 0, 0}, 34463fe90e8SHendrik Ranocha {RC(5.0719337296713929515090618138513639239329e-02), 0, 0, RC(2.5433377264600407582754714408877778031369e-01), RC(2.0394689005728199465736223777270858044698e-01), 0, 0, 0, 0, 0, 0, 0, 0}, 34563fe90e8SHendrik Ranocha {RC(-2.9000374717523110970388379285425896124091e-01), 0, 0, RC(1.3441873910260789889438681109414337003184e+00), RC(-2.8647779433614427309611103827036562829470e+00), RC(2.6775942995105948517211260646164815438695e+00), 0, 0, 0, 0, 0, 0, 0}, 34663fe90e8SHendrik Ranocha {RC(9.8535011337993546469740402980727014284756e-02), 0, 0, 0, RC(2.2192680630751384842024036498197387903583e-01), RC(-1.8140622911806994312690338288073952457474e-01), RC(1.0944411472562548236922614918038631254153e-02), 0, 0, 0, 0, 0, 0}, 34763fe90e8SHendrik Ranocha {RC(3.8711052545731144679444618165166373405645e-01), 0, 0, RC(-1.4424454974855277571256745553077927767173e+00), RC(2.9053981890699509317691346449233848441744e+00), RC(-1.8537710696301059290843332675811978025183e+00), RC(1.4003648098728154269497325109771241479223e-01), RC(5.7273940811495816575746774624447706488753e-01), 0, 0, 0, 0, 0}, 34863fe90e8SHendrik Ranocha {RC(-1.6124403444439308100630016197913480595436e-01), 0, 0, RC(-1.7339602957358984083578404473962567894901e-01), RC(-1.3012892814065147406016812745172492529744e+00), RC(1.1379503751738617308558792131431003472124e+00), RC(-3.1747649663966880106923521138043024698980e-02), RC(9.3351293824933666439811064486056884856590e-01), RC(-8.3786318334733852703300855629616433201504e-02), 0, 0, 0, 0}, 34963fe90e8SHendrik Ranocha {RC(-1.9199444881589533281510804651483576073142e-02), 0, 0, RC(2.7330857265264284907942326254016124275617e-01), RC(-6.7534973206944372919691611210942380856240e-01), RC(3.4151849813846016071738489974728382711981e-01), RC(-6.7950064803375772478920516198524629391910e-02), RC(9.6591752247623878884265586491216376509746e-02), RC(1.3253082511182101180721038466545389951226e-01), RC(3.6854959360386113446906329951531666812946e-01), 0, 0, 0}, 35063fe90e8SHendrik Ranocha {RC(6.0918774036452898676888412111588817784584e-01), 0, 0, RC(-2.2725690858980016768999800931413088399719e+00), RC(4.7578983426940290068155255881914785497547e+00), RC(-5.5161067066927584824294689667844248244842e+00), RC(2.9005963696801192709095818565946174378180e-01), RC(5.6914239633590368229109858454801849145630e-01), RC(7.9267957603321670271339916205893327579951e-01), RC(1.5473720453288822894126190771849898232047e-01), RC(1.6149708956621816247083215106334544434974e+00), 0, 0}, 3519371c9d4SSatish Balay {RC(8.8735762208534719663211694051981022704884e-01), 0, 0, RC(-2.9754597821085367558513632804709301581977e+00), RC(5.6007170094881630597990392548350098923829e+00), RC(-5.9156074505366744680014930189941657351840e+00), RC(2.2029689156134927016879142540807638331238e-01), RC(1.0155097824462216666143271340902996997549e-01), RC(1.1514345647386055909780397752125850553556e+00), RC(1.9297101665271239396134361900805843653065e+00), 0, 0, 0} 3529371c9d4SSatish Balay }, 3539371c9d4SSatish Balay b[13] = {RC(4.4729564666695714203015840429049382466467e-02), 0, 0, 0, 0, RC(1.5691033527708199813368698010726645409175e-01), RC(1.8460973408151637740702451873526277892035e-01), RC(2.2516380602086991042479419400350721970920e-01), RC(1.4794615651970234687005179885449141753736e-01), RC(7.6055542444955825269798361910336491012732e-02), RC(1.2277290235018619610824346315921437388535e-01), RC(4.1811958638991631583384842800871882376786e-02), 0}, bembed[13] = {RC(4.5847111400495925878664730122010282095875e-02), 0, 0, 0, 0, RC(2.6231891404152387437443356584845803392392e-01), RC(1.9169372337852611904485738635688429008025e-01), RC(2.1709172327902618330978407422906448568196e-01), RC(1.2738189624833706796803169450656737867900e-01), RC(1.1510530385365326258240515750043192148894e-01), 0, 0, RC(4.0561327798437566841823391436583608050053e-02)}; 3549566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK8VR, 8, 13, &A[0][0], b, NULL, bembed, 0, NULL)); 35537acaa02SHendrik Ranocha } 3564e82626cSLisandro Dalcin #undef RC 3573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 358db046809SBarry Smith } 359db046809SBarry Smith 360f68a32c8SEmil Constantinescu /*@C 361bcf0153eSBarry Smith TSRKRegisterDestroy - Frees the list of schemes that were registered by `TSRKRegister()`. 362f68a32c8SEmil Constantinescu 363f68a32c8SEmil Constantinescu Not Collective 364f68a32c8SEmil Constantinescu 365f68a32c8SEmil Constantinescu Level: advanced 366f68a32c8SEmil Constantinescu 3671cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKRegister()`, `TSRKRegisterAll()` 368f68a32c8SEmil Constantinescu @*/ 369d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKRegisterDestroy(void) 370d71ae5a4SJacob Faibussowitsch { 371f68a32c8SEmil Constantinescu RKTableauLink link; 372f68a32c8SEmil Constantinescu 373f68a32c8SEmil Constantinescu PetscFunctionBegin; 374f68a32c8SEmil Constantinescu while ((link = RKTableauList)) { 375f68a32c8SEmil Constantinescu RKTableau t = &link->tab; 376f68a32c8SEmil Constantinescu RKTableauList = link->next; 3779566063dSJacob Faibussowitsch PetscCall(PetscFree3(t->A, t->b, t->c)); 3789566063dSJacob Faibussowitsch PetscCall(PetscFree(t->bembed)); 3799566063dSJacob Faibussowitsch PetscCall(PetscFree(t->binterp)); 3809566063dSJacob Faibussowitsch PetscCall(PetscFree(t->name)); 3819566063dSJacob Faibussowitsch PetscCall(PetscFree(link)); 382f68a32c8SEmil Constantinescu } 383f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_FALSE; 3843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 385f68a32c8SEmil Constantinescu } 386f68a32c8SEmil Constantinescu 387f68a32c8SEmil Constantinescu /*@C 388bcf0153eSBarry Smith TSRKInitializePackage - This function initializes everything in the `TSRK` package. It is called 389bcf0153eSBarry Smith from `TSInitializePackage()`. 390f68a32c8SEmil Constantinescu 391f68a32c8SEmil Constantinescu Level: developer 392f68a32c8SEmil Constantinescu 3931cc06b55SBarry Smith .seealso: [](ch_ts), `TSInitializePackage()`, `PetscInitialize()`, `TSRKFinalizePackage()` 394f68a32c8SEmil Constantinescu @*/ 395d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKInitializePackage(void) 396d71ae5a4SJacob Faibussowitsch { 397e4dd521cSBarry Smith PetscFunctionBegin; 3983ba16761SJacob Faibussowitsch if (TSRKPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS); 399f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_TRUE; 4009566063dSJacob Faibussowitsch PetscCall(TSRKRegisterAll()); 4019566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSRKFinalizePackage)); 4023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 403f68a32c8SEmil Constantinescu } 404186e87acSLisandro Dalcin 405f68a32c8SEmil Constantinescu /*@C 406bcf0153eSBarry Smith TSRKFinalizePackage - This function destroys everything in the `TSRK` package. It is 407bcf0153eSBarry Smith called from `PetscFinalize()`. 408186e87acSLisandro Dalcin 409f68a32c8SEmil Constantinescu Level: developer 410f68a32c8SEmil Constantinescu 4111cc06b55SBarry Smith .seealso: [](ch_ts), `PetscFinalize()`, `TSRKInitializePackage()` 412f68a32c8SEmil Constantinescu @*/ 413d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKFinalizePackage(void) 414d71ae5a4SJacob Faibussowitsch { 415f68a32c8SEmil Constantinescu PetscFunctionBegin; 416f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_FALSE; 4179566063dSJacob Faibussowitsch PetscCall(TSRKRegisterDestroy()); 4183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 419f68a32c8SEmil Constantinescu } 420f68a32c8SEmil Constantinescu 421f68a32c8SEmil Constantinescu /*@C 422bcf0153eSBarry Smith TSRKRegister - register an `TSRK` scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 423f68a32c8SEmil Constantinescu 424f68a32c8SEmil Constantinescu Not Collective, but the same schemes should be registered on all processes on which they will be used 425f68a32c8SEmil Constantinescu 426f68a32c8SEmil Constantinescu Input Parameters: 427f68a32c8SEmil Constantinescu + name - identifier for method 428f68a32c8SEmil Constantinescu . order - approximation order of method 429f68a32c8SEmil Constantinescu . s - number of stages, this is the dimension of the matrices below 430f68a32c8SEmil Constantinescu . A - stage coefficients (dimension s*s, row-major) 431f68a32c8SEmil Constantinescu . b - step completion table (dimension s; NULL to use last row of A) 432f68a32c8SEmil Constantinescu . c - abscissa (dimension s; NULL to use row sums of A) 433f68a32c8SEmil Constantinescu . bembed - completion table for embedded method (dimension s; NULL if not available) 4343a8a9803SLisandro Dalcin . p - Order of the interpolation scheme, equal to the number of columns of binterp 4353a8a9803SLisandro Dalcin - binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1) 436f68a32c8SEmil Constantinescu 437f68a32c8SEmil Constantinescu Level: advanced 438f68a32c8SEmil Constantinescu 439bcf0153eSBarry Smith Note: 440bcf0153eSBarry Smith Several `TSRK` methods are provided, this function is only needed to create new methods. 441bcf0153eSBarry Smith 4421cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK` 443f68a32c8SEmil Constantinescu @*/ 444d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKRegister(TSRKType name, PetscInt order, PetscInt s, const PetscReal A[], const PetscReal b[], const PetscReal c[], const PetscReal bembed[], PetscInt p, const PetscReal binterp[]) 445d71ae5a4SJacob Faibussowitsch { 446f68a32c8SEmil Constantinescu RKTableauLink link; 447f68a32c8SEmil Constantinescu RKTableau t; 448f68a32c8SEmil Constantinescu PetscInt i, j; 449f68a32c8SEmil Constantinescu 450f68a32c8SEmil Constantinescu PetscFunctionBegin; 4513a8a9803SLisandro Dalcin PetscValidCharPointer(name, 1); 4523a8a9803SLisandro Dalcin PetscValidRealPointer(A, 4); 4533a8a9803SLisandro Dalcin if (b) PetscValidRealPointer(b, 5); 4543a8a9803SLisandro Dalcin if (c) PetscValidRealPointer(c, 6); 4553a8a9803SLisandro Dalcin if (bembed) PetscValidRealPointer(bembed, 7); 4563a8a9803SLisandro Dalcin if (binterp || p > 1) PetscValidRealPointer(binterp, 9); 457095c51faSJed Brown PetscCheck(s >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected number of stages s %" PetscInt_FMT " >= 0", s); 4583a8a9803SLisandro Dalcin 4599566063dSJacob Faibussowitsch PetscCall(TSRKInitializePackage()); 4609566063dSJacob Faibussowitsch PetscCall(PetscNew(&link)); 461f68a32c8SEmil Constantinescu t = &link->tab; 4623a8a9803SLisandro Dalcin 4639566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &t->name)); 464f68a32c8SEmil Constantinescu t->order = order; 465f68a32c8SEmil Constantinescu t->s = s; 4669566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(s * s, &t->A, s, &t->b, s, &t->c)); 4679566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->A, A, s * s)); 4689566063dSJacob Faibussowitsch if (b) PetscCall(PetscArraycpy(t->b, b, s)); 4699371c9d4SSatish Balay else 4709371c9d4SSatish Balay for (i = 0; i < s; i++) t->b[i] = A[(s - 1) * s + i]; 4719566063dSJacob Faibussowitsch if (c) PetscCall(PetscArraycpy(t->c, c, s)); 4729371c9d4SSatish Balay else 4739371c9d4SSatish Balay for (i = 0; i < s; i++) 4749371c9d4SSatish Balay for (j = 0, t->c[i] = 0; j < s; j++) t->c[i] += A[i * s + j]; 475d760c35bSDebojyoti Ghosh t->FSAL = PETSC_TRUE; 4769371c9d4SSatish Balay for (i = 0; i < s; i++) 4779371c9d4SSatish Balay if (t->A[(s - 1) * s + i] != t->b[i]) t->FSAL = PETSC_FALSE; 4783a8a9803SLisandro Dalcin 479f68a32c8SEmil Constantinescu if (bembed) { 4809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s, &t->bembed)); 4819566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->bembed, bembed, s)); 482f68a32c8SEmil Constantinescu } 483f68a32c8SEmil Constantinescu 4849371c9d4SSatish Balay if (!binterp) { 4859371c9d4SSatish Balay p = 1; 4869371c9d4SSatish Balay binterp = t->b; 4879371c9d4SSatish Balay } 4883a8a9803SLisandro Dalcin t->p = p; 4899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s * p, &t->binterp)); 4909566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->binterp, binterp, s * p)); 4913a8a9803SLisandro Dalcin 492f68a32c8SEmil Constantinescu link->next = RKTableauList; 493f68a32c8SEmil Constantinescu RKTableauList = link; 4943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 495f68a32c8SEmil Constantinescu } 496f68a32c8SEmil Constantinescu 497d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKGetTableau_RK(TS ts, PetscInt *s, const PetscReal **A, const PetscReal **b, const PetscReal **c, const PetscReal **bembed, PetscInt *p, const PetscReal **binterp, PetscBool *FSAL) 498d71ae5a4SJacob Faibussowitsch { 4990f7a28cdSStefano Zampini TS_RK *rk = (TS_RK *)ts->data; 5000f7a28cdSStefano Zampini RKTableau tab = rk->tableau; 5010f7a28cdSStefano Zampini 5020f7a28cdSStefano Zampini PetscFunctionBegin; 5030f7a28cdSStefano Zampini if (s) *s = tab->s; 5040f7a28cdSStefano Zampini if (A) *A = tab->A; 5050f7a28cdSStefano Zampini if (b) *b = tab->b; 5060f7a28cdSStefano Zampini if (c) *c = tab->c; 5070f7a28cdSStefano Zampini if (bembed) *bembed = tab->bembed; 5080f7a28cdSStefano Zampini if (p) *p = tab->p; 5090f7a28cdSStefano Zampini if (binterp) *binterp = tab->binterp; 5100f7a28cdSStefano Zampini if (FSAL) *FSAL = tab->FSAL; 5113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5120f7a28cdSStefano Zampini } 5130f7a28cdSStefano Zampini 5140f7a28cdSStefano Zampini /*@C 515bcf0153eSBarry Smith TSRKGetTableau - Get info on the `TSRK` tableau 5160f7a28cdSStefano Zampini 5170f7a28cdSStefano Zampini Not Collective 5180f7a28cdSStefano Zampini 519f899ff85SJose E. Roman Input Parameter: 5200f7a28cdSStefano Zampini . ts - timestepping context 5210f7a28cdSStefano Zampini 5220f7a28cdSStefano Zampini Output Parameters: 5230f7a28cdSStefano Zampini + s - number of stages, this is the dimension of the matrices below 5240f7a28cdSStefano Zampini . A - stage coefficients (dimension s*s, row-major) 5250f7a28cdSStefano Zampini . b - step completion table (dimension s) 5260f7a28cdSStefano Zampini . c - abscissa (dimension s) 5270f7a28cdSStefano Zampini . bembed - completion table for embedded method (dimension s; NULL if not available) 5280f7a28cdSStefano Zampini . p - Order of the interpolation scheme, equal to the number of columns of binterp 5290f7a28cdSStefano Zampini . binterp - Coefficients of the interpolation formula (dimension s*p) 530aaa8cc7dSPierre Jolivet - FSAL - whether or not the scheme has the First Same As Last property 5310f7a28cdSStefano Zampini 5320f7a28cdSStefano Zampini Level: developer 5330f7a28cdSStefano Zampini 5341cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKRegister()`, `TSRKSetType()` 5350f7a28cdSStefano Zampini @*/ 536d71ae5a4SJacob 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) 537d71ae5a4SJacob Faibussowitsch { 5380f7a28cdSStefano Zampini PetscFunctionBegin; 5390f7a28cdSStefano Zampini PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 5409371c9d4SSatish 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)); 5413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5420f7a28cdSStefano Zampini } 5430f7a28cdSStefano Zampini 544e4dd521cSBarry Smith /* 5452c9cb062Sluzhanghpp This is for single-step RK method 546f68a32c8SEmil Constantinescu The step completion formula is 547e4dd521cSBarry Smith 548f68a32c8SEmil Constantinescu x1 = x0 + h b^T YdotRHS 549f68a32c8SEmil Constantinescu 550f68a32c8SEmil Constantinescu This function can be called before or after ts->vec_sol has been updated. 551f68a32c8SEmil Constantinescu Suppose we have a completion formula (b) and an embedded formula (be) of different order. 552f68a32c8SEmil Constantinescu We can write 553f68a32c8SEmil Constantinescu 554f68a32c8SEmil Constantinescu x1e = x0 + h be^T YdotRHS 555f68a32c8SEmil Constantinescu = x1 - h b^T YdotRHS + h be^T YdotRHS 556f68a32c8SEmil Constantinescu = x1 + h (be - b)^T YdotRHS 557f68a32c8SEmil Constantinescu 558f68a32c8SEmil Constantinescu so we can evaluate the method with different order even after the step has been optimistically completed. 559e4dd521cSBarry Smith */ 560d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateStep_RK(TS ts, PetscInt order, Vec X, PetscBool *done) 561d71ae5a4SJacob Faibussowitsch { 562f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK *)ts->data; 563f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 564f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 565f68a32c8SEmil Constantinescu PetscReal h; 566f68a32c8SEmil Constantinescu PetscInt s = tab->s, j; 567f68a32c8SEmil Constantinescu 568f68a32c8SEmil Constantinescu PetscFunctionBegin; 569f68a32c8SEmil Constantinescu switch (rk->status) { 570f68a32c8SEmil Constantinescu case TS_STEP_INCOMPLETE: 571d71ae5a4SJacob Faibussowitsch case TS_STEP_PENDING: 572d71ae5a4SJacob Faibussowitsch h = ts->time_step; 573d71ae5a4SJacob Faibussowitsch break; 574d71ae5a4SJacob Faibussowitsch case TS_STEP_COMPLETE: 575d71ae5a4SJacob Faibussowitsch h = ts->ptime - ts->ptime_prev; 576d71ae5a4SJacob Faibussowitsch break; 577d71ae5a4SJacob Faibussowitsch default: 578d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus"); 579f68a32c8SEmil Constantinescu } 580f68a32c8SEmil Constantinescu if (order == tab->order) { 581f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { 5829566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X)); 58390b67129SHong Zhang for (j = 0; j < s; j++) w[j] = h * tab->b[j] / rk->dtratio; 5849566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, rk->YdotRHS)); 5859566063dSJacob Faibussowitsch } else PetscCall(VecCopy(ts->vec_sol, X)); 5863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 587f68a32c8SEmil Constantinescu } else if (order == tab->order - 1) { 588f68a32c8SEmil Constantinescu if (!tab->bembed) goto unavailable; 589f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { /*Complete with the embedded method (be)*/ 5909566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X)); 591f68a32c8SEmil Constantinescu for (j = 0; j < s; j++) w[j] = h * tab->bembed[j]; 5929566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, rk->YdotRHS)); 593f68a32c8SEmil Constantinescu } else { /*Rollback and re-complete using (be-b) */ 5949566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X)); 595f68a32c8SEmil Constantinescu for (j = 0; j < s; j++) w[j] = h * (tab->bembed[j] - tab->b[j]); 5969566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, rk->YdotRHS)); 597f68a32c8SEmil Constantinescu } 598f68a32c8SEmil Constantinescu if (done) *done = PETSC_TRUE; 5993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 600f68a32c8SEmil Constantinescu } 601f68a32c8SEmil Constantinescu unavailable: 602f68a32c8SEmil Constantinescu if (done) *done = PETSC_FALSE; 6039371c9d4SSatish Balay else 6049371c9d4SSatish 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); 6053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 606f68a32c8SEmil Constantinescu } 607f68a32c8SEmil Constantinescu 608d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardCostIntegral_RK(TS ts) 609d71ae5a4SJacob Faibussowitsch { 6100f7a1166SHong Zhang TS_RK *rk = (TS_RK *)ts->data; 611cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 6120f7a1166SHong Zhang RKTableau tab = rk->tableau; 6130f7a1166SHong Zhang const PetscInt s = tab->s; 6140f7a1166SHong Zhang const PetscReal *b = tab->b, *c = tab->c; 6150f7a1166SHong Zhang Vec *Y = rk->Y; 6160f7a1166SHong Zhang PetscInt i; 6170f7a1166SHong Zhang 6180f7a1166SHong Zhang PetscFunctionBegin; 619cd4cee2dSHong Zhang /* No need to backup quadts->vec_sol since it can be reverted in TSRollBack_RK */ 6200f7a1166SHong Zhang for (i = s - 1; i >= 0; i--) { 621cd4cee2dSHong Zhang /* Evolve quadrature TS solution to compute integrals */ 6229566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(quadts, rk->ptime + rk->time_step * c[i], Y[i], ts->vec_costintegrand)); 6239566063dSJacob Faibussowitsch PetscCall(VecAXPY(quadts->vec_sol, rk->time_step * b[i], ts->vec_costintegrand)); 6240f7a1166SHong Zhang } 6253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6260f7a1166SHong Zhang } 6270f7a1166SHong Zhang 628d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointCostIntegral_RK(TS ts) 629d71ae5a4SJacob Faibussowitsch { 6300f7a1166SHong Zhang TS_RK *rk = (TS_RK *)ts->data; 6310f7a1166SHong Zhang RKTableau tab = rk->tableau; 632cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 6330f7a1166SHong Zhang const PetscInt s = tab->s; 6340f7a1166SHong Zhang const PetscReal *b = tab->b, *c = tab->c; 6350f7a1166SHong Zhang Vec *Y = rk->Y; 6360f7a1166SHong Zhang PetscInt i; 6370f7a1166SHong Zhang 6380f7a1166SHong Zhang PetscFunctionBegin; 6390f7a1166SHong Zhang for (i = s - 1; i >= 0; i--) { 640cd4cee2dSHong Zhang /* Evolve quadrature TS solution to compute integrals */ 6419566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(quadts, ts->ptime + ts->time_step * (1.0 - c[i]), Y[i], ts->vec_costintegrand)); 6429566063dSJacob Faibussowitsch PetscCall(VecAXPY(quadts->vec_sol, -ts->time_step * b[i], ts->vec_costintegrand)); 6430f7a1166SHong Zhang } 6443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6450f7a1166SHong Zhang } 6460f7a1166SHong Zhang 647d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRollBack_RK(TS ts) 648d71ae5a4SJacob Faibussowitsch { 649fecfb714SLisandro Dalcin TS_RK *rk = (TS_RK *)ts->data; 650cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 651fecfb714SLisandro Dalcin RKTableau tab = rk->tableau; 652fecfb714SLisandro Dalcin const PetscInt s = tab->s; 653cd4cee2dSHong Zhang const PetscReal *b = tab->b, *c = tab->c; 654fecfb714SLisandro Dalcin PetscScalar *w = rk->work; 655cd4cee2dSHong Zhang Vec *Y = rk->Y, *YdotRHS = rk->YdotRHS; 656fecfb714SLisandro Dalcin PetscInt j; 657fecfb714SLisandro Dalcin PetscReal h; 658fecfb714SLisandro Dalcin 659fecfb714SLisandro Dalcin PetscFunctionBegin; 660fecfb714SLisandro Dalcin switch (rk->status) { 661fecfb714SLisandro Dalcin case TS_STEP_INCOMPLETE: 662d71ae5a4SJacob Faibussowitsch case TS_STEP_PENDING: 663d71ae5a4SJacob Faibussowitsch h = ts->time_step; 664d71ae5a4SJacob Faibussowitsch break; 665d71ae5a4SJacob Faibussowitsch case TS_STEP_COMPLETE: 666d71ae5a4SJacob Faibussowitsch h = ts->ptime - ts->ptime_prev; 667d71ae5a4SJacob Faibussowitsch break; 668d71ae5a4SJacob Faibussowitsch default: 669d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus"); 670fecfb714SLisandro Dalcin } 671fecfb714SLisandro Dalcin for (j = 0; j < s; j++) w[j] = -h * b[j]; 6729566063dSJacob Faibussowitsch PetscCall(VecMAXPY(ts->vec_sol, s, w, YdotRHS)); 673cd4cee2dSHong Zhang if (quadts && ts->costintegralfwd) { 674cd4cee2dSHong Zhang for (j = 0; j < s; j++) { 675cd4cee2dSHong Zhang /* Revert the quadrature TS solution */ 6769566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(quadts, rk->ptime + h * c[j], Y[j], ts->vec_costintegrand)); 6779566063dSJacob Faibussowitsch PetscCall(VecAXPY(quadts->vec_sol, -h * b[j], ts->vec_costintegrand)); 678cd4cee2dSHong Zhang } 679cd4cee2dSHong Zhang } 6803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 681fecfb714SLisandro Dalcin } 682fecfb714SLisandro Dalcin 683d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardStep_RK(TS ts) 684d71ae5a4SJacob Faibussowitsch { 685922a638cSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 686922a638cSHong Zhang RKTableau tab = rk->tableau; 687922a638cSHong Zhang Mat J, *MatsFwdSensipTemp = rk->MatsFwdSensipTemp; 688922a638cSHong Zhang const PetscInt s = tab->s; 689922a638cSHong Zhang const PetscReal *A = tab->A, *c = tab->c, *b = tab->b; 690922a638cSHong Zhang Vec *Y = rk->Y; 691922a638cSHong Zhang PetscInt i, j; 692922a638cSHong Zhang PetscReal stage_time, h = ts->time_step; 693922a638cSHong Zhang PetscBool zero; 694922a638cSHong Zhang 695922a638cSHong Zhang PetscFunctionBegin; 6969566063dSJacob Faibussowitsch PetscCall(MatCopy(ts->mat_sensip, rk->MatFwdSensip0, SAME_NONZERO_PATTERN)); 6979566063dSJacob Faibussowitsch PetscCall(TSGetRHSJacobian(ts, &J, NULL, NULL, NULL)); 698922a638cSHong Zhang 699922a638cSHong Zhang for (i = 0; i < s; i++) { 700922a638cSHong Zhang stage_time = ts->ptime + h * c[i]; 701922a638cSHong Zhang zero = PETSC_FALSE; 702922a638cSHong Zhang if (b[i] == 0 && i == s - 1) zero = PETSC_TRUE; 703922a638cSHong Zhang /* TLM Stage values */ 704922a638cSHong Zhang if (!i) { 7059566063dSJacob Faibussowitsch PetscCall(MatCopy(ts->mat_sensip, rk->MatsFwdStageSensip[i], SAME_NONZERO_PATTERN)); 706922a638cSHong Zhang } else if (!zero) { 7079566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(rk->MatsFwdStageSensip[i])); 70848a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(MatAXPY(rk->MatsFwdStageSensip[i], h * A[i * s + j], MatsFwdSensipTemp[j], SAME_NONZERO_PATTERN)); 7099566063dSJacob Faibussowitsch PetscCall(MatAXPY(rk->MatsFwdStageSensip[i], 1., ts->mat_sensip, SAME_NONZERO_PATTERN)); 710922a638cSHong Zhang } else { 7119566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(rk->MatsFwdStageSensip[i])); 712922a638cSHong Zhang } 713922a638cSHong Zhang 7149566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobian(ts, stage_time, Y[i], J, J)); 7159566063dSJacob Faibussowitsch PetscCall(MatMatMult(J, rk->MatsFwdStageSensip[i], MAT_REUSE_MATRIX, PETSC_DEFAULT, &MatsFwdSensipTemp[i])); 71613af1a74SHong Zhang if (ts->Jacprhs) { 7179566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobianP(ts, stage_time, Y[i], ts->Jacprhs)); /* get f_p */ 71813af1a74SHong Zhang if (ts->vecs_sensi2p) { /* TLM used for 2nd-order adjoint */ 71913af1a74SHong Zhang PetscScalar *xarr; 7209566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(MatsFwdSensipTemp[i], 0, &xarr)); 7219566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(rk->VecDeltaFwdSensipCol, xarr)); 7229566063dSJacob Faibussowitsch PetscCall(MatMultAdd(ts->Jacprhs, ts->vec_dir, rk->VecDeltaFwdSensipCol, rk->VecDeltaFwdSensipCol)); 7239566063dSJacob Faibussowitsch PetscCall(VecResetArray(rk->VecDeltaFwdSensipCol)); 7249566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(MatsFwdSensipTemp[i], &xarr)); 72513af1a74SHong Zhang } else { 7269566063dSJacob Faibussowitsch PetscCall(MatAXPY(MatsFwdSensipTemp[i], 1., ts->Jacprhs, SUBSET_NONZERO_PATTERN)); 72713af1a74SHong Zhang } 728922a638cSHong Zhang } 729922a638cSHong Zhang } 730922a638cSHong Zhang 73148a46eb9SPierre Jolivet for (i = 0; i < s; i++) PetscCall(MatAXPY(ts->mat_sensip, h * b[i], rk->MatsFwdSensipTemp[i], SAME_NONZERO_PATTERN)); 732922a638cSHong Zhang rk->status = TS_STEP_COMPLETE; 7333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 734922a638cSHong Zhang } 735922a638cSHong Zhang 736d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardGetStages_RK(TS ts, PetscInt *ns, Mat **stagesensip) 737d71ae5a4SJacob Faibussowitsch { 738922a638cSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 739922a638cSHong Zhang RKTableau tab = rk->tableau; 740922a638cSHong Zhang 741922a638cSHong Zhang PetscFunctionBegin; 742922a638cSHong Zhang if (ns) *ns = tab->s; 743922a638cSHong Zhang if (stagesensip) *stagesensip = rk->MatsFwdStageSensip; 7443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 745922a638cSHong Zhang } 746922a638cSHong Zhang 747d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardSetUp_RK(TS ts) 748d71ae5a4SJacob Faibussowitsch { 749922a638cSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 750922a638cSHong Zhang RKTableau tab = rk->tableau; 751922a638cSHong Zhang PetscInt i; 752922a638cSHong Zhang 753922a638cSHong Zhang PetscFunctionBegin; 754922a638cSHong Zhang /* backup sensitivity results for roll-backs */ 7559566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &rk->MatFwdSensip0)); 756922a638cSHong Zhang 7579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &rk->MatsFwdStageSensip)); 7589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &rk->MatsFwdSensipTemp)); 759922a638cSHong Zhang for (i = 0; i < tab->s; i++) { 7609566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &rk->MatsFwdStageSensip[i])); 7619566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &rk->MatsFwdSensipTemp[i])); 762922a638cSHong Zhang } 7639566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &rk->VecDeltaFwdSensipCol)); 7643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 765922a638cSHong Zhang } 766922a638cSHong Zhang 767d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardReset_RK(TS ts) 768d71ae5a4SJacob Faibussowitsch { 769922a638cSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 770922a638cSHong Zhang RKTableau tab = rk->tableau; 771922a638cSHong Zhang PetscInt i; 772922a638cSHong Zhang 773922a638cSHong Zhang PetscFunctionBegin; 7749566063dSJacob Faibussowitsch PetscCall(MatDestroy(&rk->MatFwdSensip0)); 775922a638cSHong Zhang if (rk->MatsFwdStageSensip) { 77648a46eb9SPierre Jolivet for (i = 0; i < tab->s; i++) PetscCall(MatDestroy(&rk->MatsFwdStageSensip[i])); 7779566063dSJacob Faibussowitsch PetscCall(PetscFree(rk->MatsFwdStageSensip)); 778922a638cSHong Zhang } 779922a638cSHong Zhang if (rk->MatsFwdSensipTemp) { 78048a46eb9SPierre Jolivet for (i = 0; i < tab->s; i++) PetscCall(MatDestroy(&rk->MatsFwdSensipTemp[i])); 7819566063dSJacob Faibussowitsch PetscCall(PetscFree(rk->MatsFwdSensipTemp)); 782922a638cSHong Zhang } 7839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rk->VecDeltaFwdSensipCol)); 7843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 785922a638cSHong Zhang } 786922a638cSHong Zhang 787d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_RK(TS ts) 788d71ae5a4SJacob Faibussowitsch { 789f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK *)ts->data; 790f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 791f68a32c8SEmil Constantinescu const PetscInt s = tab->s; 792fecfb714SLisandro Dalcin const PetscReal *A = tab->A, *c = tab->c; 793f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 794f68a32c8SEmil Constantinescu Vec *Y = rk->Y, *YdotRHS = rk->YdotRHS; 795*630f8c86SStefano Zampini PetscBool FSAL = (PetscBool)(tab->FSAL && !rk->newtableau); 796f68a32c8SEmil Constantinescu TSAdapt adapt; 797fecfb714SLisandro Dalcin PetscInt i, j; 798be5899b3SLisandro Dalcin PetscInt rejections = 0; 799be5899b3SLisandro Dalcin PetscBool stageok, accept = PETSC_TRUE; 800be5899b3SLisandro Dalcin PetscReal next_time_step = ts->time_step; 801f68a32c8SEmil Constantinescu 802f68a32c8SEmil Constantinescu PetscFunctionBegin; 803d1905564SLisandro Dalcin if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE; 8049566063dSJacob Faibussowitsch if (FSAL) PetscCall(VecCopy(YdotRHS[s - 1], YdotRHS[0])); 805*630f8c86SStefano Zampini rk->newtableau = PETSC_FALSE; 806d1905564SLisandro Dalcin 807f68a32c8SEmil Constantinescu rk->status = TS_STEP_INCOMPLETE; 808be5899b3SLisandro Dalcin while (!ts->reason && rk->status != TS_STEP_COMPLETE) { 809be5899b3SLisandro Dalcin PetscReal t = ts->ptime; 810f68a32c8SEmil Constantinescu PetscReal h = ts->time_step; 811f68a32c8SEmil Constantinescu for (i = 0; i < s; i++) { 8129be3e283SDebojyoti Ghosh rk->stage_time = t + h * c[i]; 8139566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, rk->stage_time)); 8149566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, Y[i])); 815f68a32c8SEmil Constantinescu for (j = 0; j < i; j++) w[j] = h * A[i * s + j]; 8169566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Y[i], i, w, YdotRHS)); 8179566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, rk->stage_time, i, Y)); 8189566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 8199566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(adapt, ts, rk->stage_time, Y[i], &stageok)); 820fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 8218f738a7cSLisandro Dalcin if (FSAL && !i) continue; 8229566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, t + h * c[i], Y[i], YdotRHS[i])); 823f68a32c8SEmil Constantinescu } 824be5899b3SLisandro Dalcin 825be5899b3SLisandro Dalcin rk->status = TS_STEP_INCOMPLETE; 8269566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep(ts, tab->order, ts->vec_sol, NULL)); 827f68a32c8SEmil Constantinescu rk->status = TS_STEP_PENDING; 8289566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 8299566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(adapt)); 8309566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidateAdd(adapt, tab->name, tab->order, 1, tab->ccfl, (PetscReal)tab->s, PETSC_TRUE)); 8319566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 832be5899b3SLisandro Dalcin rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 833be5899b3SLisandro Dalcin if (!accept) { /* Roll back the current step */ 8349566063dSJacob Faibussowitsch PetscCall(TSRollBack_RK(ts)); 835be5899b3SLisandro Dalcin ts->time_step = next_time_step; 836be5899b3SLisandro Dalcin goto reject_step; 837be5899b3SLisandro Dalcin } 838be5899b3SLisandro Dalcin 8390f7a1166SHong Zhang if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 8400f7a1166SHong Zhang rk->ptime = ts->ptime; 8410f7a1166SHong Zhang rk->time_step = ts->time_step; 842493ed95fSHong Zhang } 843be5899b3SLisandro Dalcin 844f68a32c8SEmil Constantinescu ts->ptime += ts->time_step; 845f68a32c8SEmil Constantinescu ts->time_step = next_time_step; 846f68a32c8SEmil Constantinescu break; 847be5899b3SLisandro Dalcin 848be5899b3SLisandro Dalcin reject_step: 8499371c9d4SSatish Balay ts->reject++; 8509371c9d4SSatish Balay accept = PETSC_FALSE; 851be5899b3SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 852be5899b3SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 85363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections)); 854f68a32c8SEmil Constantinescu } 855f68a32c8SEmil Constantinescu } 8563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 857e4dd521cSBarry Smith } 858e4dd521cSBarry Smith 859d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointSetUp_RK(TS ts) 860d71ae5a4SJacob Faibussowitsch { 861f6a906c0SBarry Smith TS_RK *rk = (TS_RK *)ts->data; 862be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 863be5899b3SLisandro Dalcin PetscInt s = tab->s; 864f6a906c0SBarry Smith 865f6a906c0SBarry Smith PetscFunctionBegin; 8663ba16761SJacob Faibussowitsch if (ts->adjointsetupcalled++) PetscFunctionReturn(PETSC_SUCCESS); 8679566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], s * ts->numcost, &rk->VecsDeltaLam)); 8689566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &rk->VecsSensiTemp)); 86948a46eb9SPierre Jolivet if (ts->vecs_sensip) PetscCall(VecDuplicate(ts->vecs_sensip[0], &rk->VecDeltaMu)); 87013af1a74SHong Zhang if (ts->vecs_sensi2) { 8719566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], s * ts->numcost, &rk->VecsDeltaLam2)); 8729566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vecs_sensi2[0], ts->numcost, &rk->VecsSensi2Temp)); 87313af1a74SHong Zhang } 87448a46eb9SPierre Jolivet if (ts->vecs_sensi2p) PetscCall(VecDuplicate(ts->vecs_sensi2p[0], &rk->VecDeltaMu2)); 8753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 876f6a906c0SBarry Smith } 877f6a906c0SBarry Smith 87835f5172aSHong Zhang /* 87935f5172aSHong Zhang Assumptions: 88035f5172aSHong Zhang - TSStep_RK() always evaluates the step with b, not bembed. 88135f5172aSHong Zhang */ 882d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointStep_RK(TS ts) 883d71ae5a4SJacob Faibussowitsch { 884c235aa8dSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 885cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 886c235aa8dSHong Zhang RKTableau tab = rk->tableau; 8873ca0882eSHong Zhang Mat J, Jpre, Jquad; 888c235aa8dSHong Zhang const PetscInt s = tab->s; 889c235aa8dSHong Zhang const PetscReal *A = tab->A, *b = tab->b, *c = tab->c; 89013af1a74SHong Zhang PetscScalar *w = rk->work, *xarr; 8912e7b7f96SHong Zhang Vec *Y = rk->Y, *VecsDeltaLam = rk->VecsDeltaLam, VecDeltaMu = rk->VecDeltaMu, *VecsSensiTemp = rk->VecsSensiTemp; 89213af1a74SHong Zhang Vec *VecsDeltaLam2 = rk->VecsDeltaLam2, VecDeltaMu2 = rk->VecDeltaMu2, *VecsSensi2Temp = rk->VecsSensi2Temp; 893cd4cee2dSHong Zhang Vec VecDRDUTransCol = ts->vec_drdu_col, VecDRDPTransCol = ts->vec_drdp_col; 894b18ea86cSHong Zhang PetscInt i, j, nadj; 8953d3eaba7SBarry Smith PetscReal t = ts->ptime; 8963ca0882eSHong Zhang PetscReal h = ts->time_step; 897c235aa8dSHong Zhang 898d2daff3dSHong Zhang PetscFunctionBegin; 899c235aa8dSHong Zhang rk->status = TS_STEP_INCOMPLETE; 9003389c7d9SStefano Zampini 9019566063dSJacob Faibussowitsch PetscCall(TSGetRHSJacobian(ts, &J, &Jpre, NULL, NULL)); 90248a46eb9SPierre Jolivet if (quadts) PetscCall(TSGetRHSJacobian(quadts, &Jquad, NULL, NULL, NULL)); 9032e7b7f96SHong Zhang for (i = s - 1; i >= 0; i--) { 9046a1e4597SHong Zhang if (tab->FSAL && i == s - 1) { 9056a1e4597SHong Zhang /* VecsDeltaLam[nadj*s+s-1] are initialized with zeros and the values never change.*/ 9066a1e4597SHong Zhang continue; 9076a1e4597SHong Zhang } 9083ca0882eSHong Zhang rk->stage_time = t + h * (1.0 - c[i]); 9099566063dSJacob Faibussowitsch PetscCall(TSComputeSNESJacobian(ts, Y[i], J, Jpre)); 9109371c9d4SSatish Balay if (quadts) { PetscCall(TSComputeRHSJacobian(quadts, rk->stage_time, Y[i], Jquad, Jquad)); /* get r_u^T */ } 9113389c7d9SStefano Zampini if (ts->vecs_sensip) { 9129566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobianP(ts, rk->stage_time, Y[i], ts->Jacprhs)); /* get f_p */ 9139371c9d4SSatish Balay if (quadts) { PetscCall(TSComputeRHSJacobianP(quadts, rk->stage_time, Y[i], quadts->Jacprhs)); /* get f_p for the quadrature */ } 91435f5172aSHong Zhang } 9153389c7d9SStefano Zampini 91635f5172aSHong Zhang if (b[i]) { 91735f5172aSHong Zhang for (j = i + 1; j < s; j++) w[j - i - 1] = A[j * s + i] / b[i]; /* coefficients for computing VecsSensiTemp */ 91835f5172aSHong Zhang } else { 91935f5172aSHong Zhang for (j = i + 1; j < s; j++) w[j - i - 1] = A[j * s + i]; /* coefficients for computing VecsSensiTemp */ 92035f5172aSHong Zhang } 9212e7b7f96SHong Zhang 9222e7b7f96SHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 9233389c7d9SStefano Zampini /* Stage values of lambda */ 92435f5172aSHong Zhang if (b[i]) { 92535f5172aSHong Zhang /* lambda_{n+1} + \sum_{j=i+1}^s a_{ji}/b[i]*lambda_{s,j} */ 9269566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vecs_sensi[nadj], VecsSensiTemp[nadj])); /* VecDeltaLam is an vec array of size s by numcost */ 9279566063dSJacob Faibussowitsch PetscCall(VecMAXPY(VecsSensiTemp[nadj], s - i - 1, w, &VecsDeltaLam[nadj * s + i + 1])); 9289566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(J, VecsSensiTemp[nadj], VecsDeltaLam[nadj * s + i])); /* VecsSensiTemp will be reused by 2nd-order adjoint */ 9299566063dSJacob Faibussowitsch PetscCall(VecScale(VecsDeltaLam[nadj * s + i], -h * b[i])); 93035f5172aSHong Zhang if (quadts) { 9319566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(Jquad, nadj, &xarr)); 9329566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(VecDRDUTransCol, xarr)); 9339566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecsDeltaLam[nadj * s + i], -h * b[i], VecDRDUTransCol)); 9349566063dSJacob Faibussowitsch PetscCall(VecResetArray(VecDRDUTransCol)); 9359566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(Jquad, &xarr)); 936cd4cee2dSHong Zhang } 9373389c7d9SStefano Zampini } else { 9386a1e4597SHong Zhang /* \sum_{j=i+1}^s a_{ji}*lambda_{s,j} */ 9399566063dSJacob Faibussowitsch PetscCall(VecSet(VecsSensiTemp[nadj], 0)); 9409566063dSJacob Faibussowitsch PetscCall(VecMAXPY(VecsSensiTemp[nadj], s - i - 1, w, &VecsDeltaLam[nadj * s + i + 1])); 9419566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(J, VecsSensiTemp[nadj], VecsDeltaLam[nadj * s + i])); 9429566063dSJacob Faibussowitsch PetscCall(VecScale(VecsDeltaLam[nadj * s + i], -h)); 9433389c7d9SStefano Zampini } 9446497ce14SHong Zhang 945ad8e2604SHong Zhang /* Stage values of mu */ 9466497ce14SHong Zhang if (ts->vecs_sensip) { 94735f5172aSHong Zhang if (b[i]) { 9489566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ts->Jacprhs, VecsSensiTemp[nadj], VecDeltaMu)); 9499566063dSJacob Faibussowitsch PetscCall(VecScale(VecDeltaMu, -h * b[i])); 95035f5172aSHong Zhang if (quadts) { 9519566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(quadts->Jacprhs, nadj, &xarr)); 9529566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(VecDRDPTransCol, xarr)); 9539566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecDeltaMu, -h * b[i], VecDRDPTransCol)); 9549566063dSJacob Faibussowitsch PetscCall(VecResetArray(VecDRDPTransCol)); 9559566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(quadts->Jacprhs, &xarr)); 956493ed95fSHong Zhang } 95735f5172aSHong Zhang } else { 9589566063dSJacob Faibussowitsch PetscCall(VecScale(VecDeltaMu, -h)); 95935f5172aSHong Zhang } 9609566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensip[nadj], 1., VecDeltaMu)); /* update sensip for each stage */ 961ad8e2604SHong Zhang } 962c235aa8dSHong Zhang } 96313af1a74SHong Zhang 96413af1a74SHong Zhang if (ts->vecs_sensi2 && ts->forward_solve) { /* 2nd-order adjoint, TLM mode has to be turned on */ 96513af1a74SHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 9669566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(rk->MatsFwdStageSensip[i], 0, &xarr)); 9679566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr)); 96813af1a74SHong Zhang /* lambda_s^T F_UU w_1 */ 9699566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionUU(ts, rk->stage_time, Y[i], VecsSensiTemp, ts->vec_sensip_col, ts->vecs_guu)); 97035f5172aSHong Zhang if (quadts) { 97135f5172aSHong Zhang /* R_UU w_1 */ 9729566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionUU(quadts, rk->stage_time, Y[i], NULL, ts->vec_sensip_col, ts->vecs_guu)); 97335f5172aSHong Zhang } 97435f5172aSHong Zhang if (ts->vecs_sensip) { 97513af1a74SHong Zhang /* lambda_s^T F_UP w_2 */ 9769566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionUP(ts, rk->stage_time, Y[i], VecsSensiTemp, ts->vec_dir, ts->vecs_gup)); 97735f5172aSHong Zhang if (quadts) { 97835f5172aSHong Zhang /* R_UP w_2 */ 9799566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionUP(quadts, rk->stage_time, Y[i], NULL, ts->vec_sensip_col, ts->vecs_gup)); 98035f5172aSHong Zhang } 98135f5172aSHong Zhang } 98213af1a74SHong Zhang if (ts->vecs_sensi2p) { 98313af1a74SHong Zhang /* lambda_s^T F_PU w_1 */ 9849566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionPU(ts, rk->stage_time, Y[i], VecsSensiTemp, ts->vec_sensip_col, ts->vecs_gpu)); 98535f5172aSHong Zhang /* lambda_s^T F_PP w_2 */ 9869566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionPP(ts, rk->stage_time, Y[i], VecsSensiTemp, ts->vec_dir, ts->vecs_gpp)); 98735f5172aSHong Zhang if (b[i] && quadts) { 98835f5172aSHong Zhang /* R_PU w_1 */ 9899566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionPU(quadts, rk->stage_time, Y[i], NULL, ts->vec_sensip_col, ts->vecs_gpu)); 99035f5172aSHong Zhang /* R_PP w_2 */ 9919566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionPP(quadts, rk->stage_time, Y[i], NULL, ts->vec_dir, ts->vecs_gpp)); 99235f5172aSHong Zhang } 99313af1a74SHong Zhang } 9949566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_sensip_col)); 9959566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(rk->MatsFwdStageSensip[i], &xarr)); 99613af1a74SHong Zhang 99713af1a74SHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 99813af1a74SHong Zhang /* Stage values of lambda */ 99935f5172aSHong Zhang if (b[i]) { 100035f5172aSHong Zhang /* J_i^T*(Lambda_{n+1}+\sum_{j=i+1}^s a_{ji}/b_i*Lambda_{s,j} */ 10019566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vecs_sensi2[nadj], VecsSensi2Temp[nadj])); 10029566063dSJacob Faibussowitsch PetscCall(VecMAXPY(VecsSensi2Temp[nadj], s - i - 1, w, &VecsDeltaLam2[nadj * s + i + 1])); 10039566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(J, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj * s + i])); 10049566063dSJacob Faibussowitsch PetscCall(VecScale(VecsDeltaLam2[nadj * s + i], -h * b[i])); 10059566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecsDeltaLam2[nadj * s + i], -h * b[i], ts->vecs_guu[nadj])); 100648a46eb9SPierre Jolivet if (ts->vecs_sensip) PetscCall(VecAXPY(VecsDeltaLam2[nadj * s + i], -h * b[i], ts->vecs_gup[nadj])); 100713af1a74SHong Zhang } else { 100835f5172aSHong Zhang /* \sum_{j=i+1}^s a_{ji}*Lambda_{s,j} */ 10099566063dSJacob Faibussowitsch PetscCall(VecSet(VecsDeltaLam2[nadj * s + i], 0)); 10109566063dSJacob Faibussowitsch PetscCall(VecMAXPY(VecsSensi2Temp[nadj], s - i - 1, w, &VecsDeltaLam2[nadj * s + i + 1])); 10119566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(J, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj * s + i])); 10129566063dSJacob Faibussowitsch PetscCall(VecScale(VecsDeltaLam2[nadj * s + i], -h)); 10139566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecsDeltaLam2[nadj * s + i], -h, ts->vecs_guu[nadj])); 101448a46eb9SPierre Jolivet if (ts->vecs_sensip) PetscCall(VecAXPY(VecsDeltaLam2[nadj * s + i], -h, ts->vecs_gup[nadj])); 101535f5172aSHong Zhang } 101635f5172aSHong Zhang if (ts->vecs_sensi2p) { /* 2nd-order adjoint for parameters */ 10179566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ts->Jacprhs, VecsSensi2Temp[nadj], VecDeltaMu2)); 101835f5172aSHong Zhang if (b[i]) { 10199566063dSJacob Faibussowitsch PetscCall(VecScale(VecDeltaMu2, -h * b[i])); 10209566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecDeltaMu2, -h * b[i], ts->vecs_gpu[nadj])); 10219566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecDeltaMu2, -h * b[i], ts->vecs_gpp[nadj])); 102213af1a74SHong Zhang } else { 10239566063dSJacob Faibussowitsch PetscCall(VecScale(VecDeltaMu2, -h)); 10249566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecDeltaMu2, -h, ts->vecs_gpu[nadj])); 10259566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecDeltaMu2, -h, ts->vecs_gpp[nadj])); 102613af1a74SHong Zhang } 10279566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], 1, VecDeltaMu2)); /* update sensi2p for each stage */ 102813af1a74SHong Zhang } 102913af1a74SHong Zhang } 103013af1a74SHong Zhang } 10316497ce14SHong Zhang } 1032c235aa8dSHong Zhang 1033c235aa8dSHong Zhang for (j = 0; j < s; j++) w[j] = 1.0; 10342e7b7f96SHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { /* no need to do this for mu's */ 10359566063dSJacob Faibussowitsch PetscCall(VecMAXPY(ts->vecs_sensi[nadj], s, w, &VecsDeltaLam[nadj * s])); 103648a46eb9SPierre Jolivet if (ts->vecs_sensi2) PetscCall(VecMAXPY(ts->vecs_sensi2[nadj], s, w, &VecsDeltaLam2[nadj * s])); 10376497ce14SHong Zhang } 1038c235aa8dSHong Zhang rk->status = TS_STEP_COMPLETE; 10393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1040d2daff3dSHong Zhang } 1041d2daff3dSHong Zhang 1042d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointReset_RK(TS ts) 1043d71ae5a4SJacob Faibussowitsch { 104413af1a74SHong Zhang TS_RK *rk = (TS_RK *)ts->data; 104513af1a74SHong Zhang RKTableau tab = rk->tableau; 104613af1a74SHong Zhang 104713af1a74SHong Zhang PetscFunctionBegin; 10489566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s * ts->numcost, &rk->VecsDeltaLam)); 10499566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(ts->numcost, &rk->VecsSensiTemp)); 10509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rk->VecDeltaMu)); 10519566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s * ts->numcost, &rk->VecsDeltaLam2)); 10529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rk->VecDeltaMu2)); 10539566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(ts->numcost, &rk->VecsSensi2Temp)); 10543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 105513af1a74SHong Zhang } 105613af1a74SHong Zhang 1057d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_RK(TS ts, PetscReal itime, Vec X) 1058d71ae5a4SJacob Faibussowitsch { 105955de54fcSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 106055de54fcSHong Zhang PetscInt s = rk->tableau->s, p = rk->tableau->p, i, j; 106155de54fcSHong Zhang PetscReal h; 106255de54fcSHong Zhang PetscReal tt, t; 106355de54fcSHong Zhang PetscScalar *b; 106455de54fcSHong Zhang const PetscReal *B = rk->tableau->binterp; 106555de54fcSHong Zhang 106655de54fcSHong Zhang PetscFunctionBegin; 10673c633725SBarry Smith PetscCheck(B, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSRK %s does not have an interpolation formula", rk->tableau->name); 106855de54fcSHong Zhang 106955de54fcSHong Zhang switch (rk->status) { 107055de54fcSHong Zhang case TS_STEP_INCOMPLETE: 107155de54fcSHong Zhang case TS_STEP_PENDING: 107255de54fcSHong Zhang h = ts->time_step; 107355de54fcSHong Zhang t = (itime - ts->ptime) / h; 107455de54fcSHong Zhang break; 107555de54fcSHong Zhang case TS_STEP_COMPLETE: 107655de54fcSHong Zhang h = ts->ptime - ts->ptime_prev; 107755de54fcSHong Zhang t = (itime - ts->ptime) / h + 1; /* In the interval [0,1] */ 107855de54fcSHong Zhang break; 1079d71ae5a4SJacob Faibussowitsch default: 1080d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus"); 108155de54fcSHong Zhang } 10829566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s, &b)); 108355de54fcSHong Zhang for (i = 0; i < s; i++) b[i] = 0; 108455de54fcSHong Zhang for (j = 0, tt = t; j < p; j++, tt *= t) { 1085ad540459SPierre Jolivet for (i = 0; i < s; i++) b[i] += h * B[i * p + j] * tt; 108655de54fcSHong Zhang } 10879566063dSJacob Faibussowitsch PetscCall(VecCopy(rk->Y[0], X)); 10889566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, b, rk->YdotRHS)); 10899566063dSJacob Faibussowitsch PetscCall(PetscFree(b)); 10903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 109155de54fcSHong Zhang } 109255de54fcSHong Zhang 109355de54fcSHong Zhang /*------------------------------------------------------------*/ 109455de54fcSHong Zhang 1095d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRKTableauReset(TS ts) 1096d71ae5a4SJacob Faibussowitsch { 1097be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK *)ts->data; 1098be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 1099be5899b3SLisandro Dalcin 1100be5899b3SLisandro Dalcin PetscFunctionBegin; 11013ba16761SJacob Faibussowitsch if (!tab) PetscFunctionReturn(PETSC_SUCCESS); 11029566063dSJacob Faibussowitsch PetscCall(PetscFree(rk->work)); 11039566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &rk->Y)); 11049566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &rk->YdotRHS)); 11053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1106be5899b3SLisandro Dalcin } 1107be5899b3SLisandro Dalcin 1108d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_RK(TS ts) 1109d71ae5a4SJacob Faibussowitsch { 1110e4dd521cSBarry Smith PetscFunctionBegin; 11119566063dSJacob Faibussowitsch PetscCall(TSRKTableauReset(ts)); 11120fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 1113cac4c232SBarry Smith PetscTryMethod(ts, "TSReset_RK_MultirateSplit_C", (TS), (ts)); 11140fe4d17eSHong Zhang } else { 1115cac4c232SBarry Smith PetscTryMethod(ts, "TSReset_RK_MultirateNonsplit_C", (TS), (ts)); 11160fe4d17eSHong Zhang } 11173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1118e4dd521cSBarry Smith } 1119277b19d0SLisandro Dalcin 1120d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCoarsenHook_TSRK(DM fine, DM coarse, void *ctx) 1121d71ae5a4SJacob Faibussowitsch { 1122f68a32c8SEmil Constantinescu PetscFunctionBegin; 11233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1124f68a32c8SEmil Constantinescu } 1125f68a32c8SEmil Constantinescu 1126d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRestrictHook_TSRK(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) 1127d71ae5a4SJacob Faibussowitsch { 1128f68a32c8SEmil Constantinescu PetscFunctionBegin; 11293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1130f68a32c8SEmil Constantinescu } 1131f68a32c8SEmil Constantinescu 1132d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainHook_TSRK(DM dm, DM subdm, void *ctx) 1133d71ae5a4SJacob Faibussowitsch { 1134f68a32c8SEmil Constantinescu PetscFunctionBegin; 11353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1136f68a32c8SEmil Constantinescu } 1137f68a32c8SEmil Constantinescu 1138d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) 1139d71ae5a4SJacob Faibussowitsch { 1140f68a32c8SEmil Constantinescu PetscFunctionBegin; 11413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1142f68a32c8SEmil Constantinescu } 1143be5899b3SLisandro Dalcin 1144d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRKTableauSetUp(TS ts) 1145d71ae5a4SJacob Faibussowitsch { 1146be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK *)ts->data; 1147be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 1148be5899b3SLisandro Dalcin 1149be5899b3SLisandro Dalcin PetscFunctionBegin; 11509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &rk->work)); 11519566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &rk->Y)); 11529566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &rk->YdotRHS)); 1153*630f8c86SStefano Zampini rk->newtableau = PETSC_TRUE; 11543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1155be5899b3SLisandro Dalcin } 1156be5899b3SLisandro Dalcin 1157d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_RK(TS ts) 1158d71ae5a4SJacob Faibussowitsch { 1159cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 1160f68a32c8SEmil Constantinescu DM dm; 1161f68a32c8SEmil Constantinescu 1162f68a32c8SEmil Constantinescu PetscFunctionBegin; 11639566063dSJacob Faibussowitsch PetscCall(TSCheckImplicitTerm(ts)); 11649566063dSJacob Faibussowitsch PetscCall(TSRKTableauSetUp(ts)); 1165cd4cee2dSHong Zhang if (quadts && ts->costintegralfwd) { 1166cd4cee2dSHong Zhang Mat Jquad; 11679566063dSJacob Faibussowitsch PetscCall(TSGetRHSJacobian(quadts, &Jquad, NULL, NULL, NULL)); 11680f7a1166SHong Zhang } 11699566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 11709566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSRK, DMRestrictHook_TSRK, ts)); 11719566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSRK, DMSubDomainRestrictHook_TSRK, ts)); 11720fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 1173cac4c232SBarry Smith PetscTryMethod(ts, "TSSetUp_RK_MultirateSplit_C", (TS), (ts)); 11740fe4d17eSHong Zhang } else { 1175cac4c232SBarry Smith PetscTryMethod(ts, "TSSetUp_RK_MultirateNonsplit_C", (TS), (ts)); 11760fe4d17eSHong Zhang } 11773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1178e4dd521cSBarry Smith } 1179d2daff3dSHong Zhang 1180d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_RK(TS ts, PetscOptionItems *PetscOptionsObject) 1181d71ae5a4SJacob Faibussowitsch { 1182be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK *)ts->data; 1183262ff7bbSBarry Smith 1184e4dd521cSBarry Smith PetscFunctionBegin; 1185d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "RK ODE solver options"); 1186f68a32c8SEmil Constantinescu { 1187f68a32c8SEmil Constantinescu RKTableauLink link; 1188f68a32c8SEmil Constantinescu PetscInt count, choice; 11890fe4d17eSHong Zhang PetscBool flg, use_multirate = PETSC_FALSE; 1190f68a32c8SEmil Constantinescu const char **namelist; 11912c9cb062Sluzhanghpp 11929371c9d4SSatish Balay for (link = RKTableauList, count = 0; link; link = link->next, count++) 11939371c9d4SSatish Balay ; 11949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count, (char ***)&namelist)); 1195f68a32c8SEmil Constantinescu for (link = RKTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name; 11969566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_rk_multirate", "Use interpolation-based multirate RK method", "TSRKSetMultirate", rk->use_multirate, &use_multirate, &flg)); 11971baa6e33SBarry Smith if (flg) PetscCall(TSRKSetMultirate(ts, use_multirate)); 11989566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_rk_type", "Family of RK method", "TSRKSetType", (const char *const *)namelist, count, rk->tableau->name, &choice, &flg)); 11999566063dSJacob Faibussowitsch if (flg) PetscCall(TSRKSetType(ts, namelist[choice])); 12009566063dSJacob Faibussowitsch PetscCall(PetscFree(namelist)); 1201f68a32c8SEmil Constantinescu } 1202d0609cedSBarry Smith PetscOptionsHeadEnd(); 1203d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)ts), NULL, "Multirate methods options", ""); 12049566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_rk_dtratio", "time step ratio between slow and fast", "", rk->dtratio, &rk->dtratio, NULL)); 1205d0609cedSBarry Smith PetscOptionsEnd(); 12063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1207e4dd521cSBarry Smith } 1208e4dd521cSBarry Smith 1209d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_RK(TS ts, PetscViewer viewer) 1210d71ae5a4SJacob Faibussowitsch { 12115f70b478SJed Brown TS_RK *rk = (TS_RK *)ts->data; 12128370ee3bSLisandro Dalcin PetscBool iascii; 12132cdf8120Sasbjorn 1214e4dd521cSBarry Smith PetscFunctionBegin; 12159566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 12168370ee3bSLisandro Dalcin if (iascii) { 12179c334d8fSLisandro Dalcin RKTableau tab = rk->tableau; 1218f68a32c8SEmil Constantinescu TSRKType rktype; 12190f7a28cdSStefano Zampini const PetscReal *c; 12200f7a28cdSStefano Zampini PetscInt s; 1221f68a32c8SEmil Constantinescu char buf[512]; 12220f7a28cdSStefano Zampini PetscBool FSAL; 12230f7a28cdSStefano Zampini 12249566063dSJacob Faibussowitsch PetscCall(TSRKGetType(ts, &rktype)); 12259566063dSJacob Faibussowitsch PetscCall(TSRKGetTableau(ts, &s, NULL, NULL, &c, NULL, NULL, NULL, &FSAL)); 12269566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " RK type %s\n", rktype)); 122763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Order: %" PetscInt_FMT "\n", tab->order)); 12289566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " FSAL property: %s\n", FSAL ? "yes" : "no")); 12299566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", s, c)); 12309566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Abscissa c = %s\n", buf)); 12318370ee3bSLisandro Dalcin } 12323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1233f68a32c8SEmil Constantinescu } 1234f68a32c8SEmil Constantinescu 1235d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSLoad_RK(TS ts, PetscViewer viewer) 1236d71ae5a4SJacob Faibussowitsch { 12379c334d8fSLisandro Dalcin TSAdapt adapt; 1238f68a32c8SEmil Constantinescu 1239f68a32c8SEmil Constantinescu PetscFunctionBegin; 12409566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 12419566063dSJacob Faibussowitsch PetscCall(TSAdaptLoad(adapt, viewer)); 12423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1243f68a32c8SEmil Constantinescu } 1244f68a32c8SEmil Constantinescu 12452ea87ba9SLisandro Dalcin /*@ 1246bcf0153eSBarry Smith TSRKGetOrder - Get the order of the `TSRK` scheme 12472ea87ba9SLisandro Dalcin 124820f4b53cSBarry Smith Not Collective 12492ea87ba9SLisandro Dalcin 12502ea87ba9SLisandro Dalcin Input Parameter: 12512ea87ba9SLisandro Dalcin . ts - timestepping context 12522ea87ba9SLisandro Dalcin 12532ea87ba9SLisandro Dalcin Output Parameter: 1254bcf0153eSBarry Smith . order - order of `TSRK` scheme 12552ea87ba9SLisandro Dalcin 12562ea87ba9SLisandro Dalcin Level: intermediate 12572ea87ba9SLisandro Dalcin 12581cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKGetType()` 12592ea87ba9SLisandro Dalcin @*/ 1260d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKGetOrder(TS ts, PetscInt *order) 1261d71ae5a4SJacob Faibussowitsch { 12622ea87ba9SLisandro Dalcin PetscFunctionBegin; 12632ea87ba9SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 12642ea87ba9SLisandro Dalcin PetscValidIntPointer(order, 2); 1265cac4c232SBarry Smith PetscUseMethod(ts, "TSRKGetOrder_C", (TS, PetscInt *), (ts, order)); 12663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12672ea87ba9SLisandro Dalcin } 12682ea87ba9SLisandro Dalcin 1269f68a32c8SEmil Constantinescu /*@C 1270bcf0153eSBarry Smith TSRKSetType - Set the type of the `TSRK` scheme 1271f68a32c8SEmil Constantinescu 127220f4b53cSBarry Smith Logically Collective 1273f68a32c8SEmil Constantinescu 1274d8d19677SJose E. Roman Input Parameters: 1275f68a32c8SEmil Constantinescu + ts - timestepping context 1276bcf0153eSBarry Smith - rktype - type of `TSRK` scheme 1277f68a32c8SEmil Constantinescu 1278bcf0153eSBarry Smith Options Database Key: 12799bd3e852SBarry Smith . -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs> 1280287c2655SBarry Smith 1281f68a32c8SEmil Constantinescu Level: intermediate 1282f68a32c8SEmil Constantinescu 12831cc06b55SBarry Smith .seealso: [](ch_ts), `TSRKGetType()`, `TSRK`, `TSRKType`, `TSRK1FE`, `TSRK2A`, `TSRK2B`, `TSRK3`, `TSRK3BS`, `TSRK4`, `TSRK5F`, `TSRK5DP`, `TSRK5BS`, `TSRK6VR`, `TSRK7VR`, `TSRK8VR` 1284f68a32c8SEmil Constantinescu @*/ 1285d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKSetType(TS ts, TSRKType rktype) 1286d71ae5a4SJacob Faibussowitsch { 1287f68a32c8SEmil Constantinescu PetscFunctionBegin; 1288f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1289b92453a8SLisandro Dalcin PetscValidCharPointer(rktype, 2); 1290cac4c232SBarry Smith PetscTryMethod(ts, "TSRKSetType_C", (TS, TSRKType), (ts, rktype)); 12913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1292f68a32c8SEmil Constantinescu } 1293f68a32c8SEmil Constantinescu 1294f68a32c8SEmil Constantinescu /*@C 1295bcf0153eSBarry Smith TSRKGetType - Get the type of `TSRK` scheme 1296f68a32c8SEmil Constantinescu 129720f4b53cSBarry Smith Not Collective 1298f68a32c8SEmil Constantinescu 1299f68a32c8SEmil Constantinescu Input Parameter: 1300f68a32c8SEmil Constantinescu . ts - timestepping context 1301f68a32c8SEmil Constantinescu 1302f68a32c8SEmil Constantinescu Output Parameter: 1303bcf0153eSBarry Smith . rktype - type of `TSRK`-scheme 1304f68a32c8SEmil Constantinescu 1305f68a32c8SEmil Constantinescu Level: intermediate 1306f68a32c8SEmil Constantinescu 13071cc06b55SBarry Smith .seealso: [](ch_ts), `TSRKSetType()` 1308f68a32c8SEmil Constantinescu @*/ 1309d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKGetType(TS ts, TSRKType *rktype) 1310d71ae5a4SJacob Faibussowitsch { 1311f68a32c8SEmil Constantinescu PetscFunctionBegin; 1312f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1313cac4c232SBarry Smith PetscUseMethod(ts, "TSRKGetType_C", (TS, TSRKType *), (ts, rktype)); 13143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1315f68a32c8SEmil Constantinescu } 1316f68a32c8SEmil Constantinescu 1317d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRKGetOrder_RK(TS ts, PetscInt *order) 1318d71ae5a4SJacob Faibussowitsch { 13192ea87ba9SLisandro Dalcin TS_RK *rk = (TS_RK *)ts->data; 13202ea87ba9SLisandro Dalcin 13212ea87ba9SLisandro Dalcin PetscFunctionBegin; 13222ea87ba9SLisandro Dalcin *order = rk->tableau->order; 13233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13242ea87ba9SLisandro Dalcin } 13252ea87ba9SLisandro Dalcin 1326d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRKGetType_RK(TS ts, TSRKType *rktype) 1327d71ae5a4SJacob Faibussowitsch { 1328f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK *)ts->data; 1329f68a32c8SEmil Constantinescu 1330f68a32c8SEmil Constantinescu PetscFunctionBegin; 1331f68a32c8SEmil Constantinescu *rktype = rk->tableau->name; 13323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1333f68a32c8SEmil Constantinescu } 13342c9cb062Sluzhanghpp 1335d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRKSetType_RK(TS ts, TSRKType rktype) 1336d71ae5a4SJacob Faibussowitsch { 1337f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK *)ts->data; 1338f68a32c8SEmil Constantinescu PetscBool match; 1339f68a32c8SEmil Constantinescu RKTableauLink link; 1340f68a32c8SEmil Constantinescu 1341f68a32c8SEmil Constantinescu PetscFunctionBegin; 1342f68a32c8SEmil Constantinescu if (rk->tableau) { 13439566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(rk->tableau->name, rktype, &match)); 13443ba16761SJacob Faibussowitsch if (match) PetscFunctionReturn(PETSC_SUCCESS); 1345f68a32c8SEmil Constantinescu } 1346f68a32c8SEmil Constantinescu for (link = RKTableauList; link; link = link->next) { 13479566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(link->tab.name, rktype, &match)); 1348f68a32c8SEmil Constantinescu if (match) { 13499566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSRKTableauReset(ts)); 1350f68a32c8SEmil Constantinescu rk->tableau = &link->tab; 13519566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSRKTableauSetUp(ts)); 13522ffb9264SLisandro Dalcin ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 13533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1354f68a32c8SEmil Constantinescu } 1355f68a32c8SEmil Constantinescu } 135698921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", rktype); 1357e4dd521cSBarry Smith } 1358e4dd521cSBarry Smith 1359d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGetStages_RK(TS ts, PetscInt *ns, Vec **Y) 1360d71ae5a4SJacob Faibussowitsch { 1361ff22ae23SHong Zhang TS_RK *rk = (TS_RK *)ts->data; 1362ff22ae23SHong Zhang 1363ff22ae23SHong Zhang PetscFunctionBegin; 13640429704eSStefano Zampini if (ns) *ns = rk->tableau->s; 1365d2daff3dSHong Zhang if (Y) *Y = rk->Y; 13663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1367ff22ae23SHong Zhang } 1368ff22ae23SHong Zhang 1369d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_RK(TS ts) 1370d71ae5a4SJacob Faibussowitsch { 1371b3a6b972SJed Brown PetscFunctionBegin; 13729566063dSJacob Faibussowitsch PetscCall(TSReset_RK(ts)); 1373b3a6b972SJed Brown if (ts->dm) { 13749566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSRK, DMRestrictHook_TSRK, ts)); 13759566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSRK, DMSubDomainRestrictHook_TSRK, ts)); 1376b3a6b972SJed Brown } 13779566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 13789566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetOrder_C", NULL)); 13799566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetType_C", NULL)); 13809566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKSetType_C", NULL)); 13819566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetTableau_C", NULL)); 13829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKSetMultirate_C", NULL)); 13839566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetMultirate_C", NULL)); 13842e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSetUp_RK_MultirateSplit_C", NULL)); 13852e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSReset_RK_MultirateSplit_C", NULL)); 13862e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSetUp_RK_MultirateNonsplit_C", NULL)); 13872e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSReset_RK_MultirateNonsplit_C", NULL)); 13883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1389b3a6b972SJed Brown } 1390ff22ae23SHong Zhang 13913ca0882eSHong Zhang /* 13923ca0882eSHong Zhang This defines the nonlinear equation that is to be solved with SNES 13933ca0882eSHong Zhang We do not need to solve the equation; we just use SNES to approximate the Jacobian 13943ca0882eSHong Zhang */ 1395d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_RK(SNES snes, Vec x, Vec y, TS ts) 1396d71ae5a4SJacob Faibussowitsch { 13973ca0882eSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 13983ca0882eSHong Zhang DM dm, dmsave; 13993ca0882eSHong Zhang 14003ca0882eSHong Zhang PetscFunctionBegin; 14019566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 14023ca0882eSHong Zhang /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 14033ca0882eSHong Zhang dmsave = ts->dm; 14043ca0882eSHong Zhang ts->dm = dm; 14059566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, rk->stage_time, x, y)); 14063ca0882eSHong Zhang ts->dm = dmsave; 14073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14083ca0882eSHong Zhang } 14093ca0882eSHong Zhang 1410d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_RK(SNES snes, Vec x, Mat A, Mat B, TS ts) 1411d71ae5a4SJacob Faibussowitsch { 14123ca0882eSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 14133ca0882eSHong Zhang DM dm, dmsave; 14143ca0882eSHong Zhang 14153ca0882eSHong Zhang PetscFunctionBegin; 14169566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 14173ca0882eSHong Zhang dmsave = ts->dm; 14183ca0882eSHong Zhang ts->dm = dm; 14199566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobian(ts, rk->stage_time, x, A, B)); 14203ca0882eSHong Zhang ts->dm = dmsave; 14213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14223ca0882eSHong Zhang } 14233ca0882eSHong Zhang 142421052c0cSHong Zhang /*@C 1425bcf0153eSBarry Smith TSRKSetMultirate - Use the interpolation-based multirate `TSRK` method 142621052c0cSHong Zhang 142720f4b53cSBarry Smith Logically Collective 142821052c0cSHong Zhang 1429d8d19677SJose E. Roman Input Parameters: 143021052c0cSHong Zhang + ts - timestepping context 1431bcf0153eSBarry 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 143221052c0cSHong Zhang 1433bcf0153eSBarry Smith Options Database Key: 143421052c0cSHong Zhang . -ts_rk_multirate - <true,false> 143521052c0cSHong Zhang 143621052c0cSHong Zhang Level: intermediate 143721052c0cSHong Zhang 1438bcf0153eSBarry Smith Note: 1439da81f932SPierre 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). 1440bcf0153eSBarry Smith 14411cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKGetMultirate()` 144221052c0cSHong Zhang @*/ 1443d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKSetMultirate(TS ts, PetscBool use_multirate) 1444d71ae5a4SJacob Faibussowitsch { 14458a4be4dbSHong Zhang PetscFunctionBegin; 1446cac4c232SBarry Smith PetscTryMethod(ts, "TSRKSetMultirate_C", (TS, PetscBool), (ts, use_multirate)); 14473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 144821052c0cSHong Zhang } 144921052c0cSHong Zhang 145021052c0cSHong Zhang /*@C 1451bcf0153eSBarry Smith TSRKGetMultirate - Gets whether to use the interpolation-based multirate `TSRK` method 145221052c0cSHong Zhang 145320f4b53cSBarry Smith Not Collective 145421052c0cSHong Zhang 145521052c0cSHong Zhang Input Parameter: 145621052c0cSHong Zhang . ts - timestepping context 145721052c0cSHong Zhang 145821052c0cSHong Zhang Output Parameter: 1459bcf0153eSBarry Smith . use_multirate - `PETSC_TRUE` if the multirate RK method is enabled, `PETSC_FALSE` otherwise 146021052c0cSHong Zhang 146121052c0cSHong Zhang Level: intermediate 146221052c0cSHong Zhang 14631cc06b55SBarry Smith .seealso: [](ch_ts), `TSRK`, `TSRKSetMultirate()` 146421052c0cSHong Zhang @*/ 1465d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRKGetMultirate(TS ts, PetscBool *use_multirate) 1466d71ae5a4SJacob Faibussowitsch { 14677dbacdf2SHong Zhang PetscFunctionBegin; 1468cac4c232SBarry Smith PetscUseMethod(ts, "TSRKGetMultirate_C", (TS, PetscBool *), (ts, use_multirate)); 14693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 147021052c0cSHong Zhang } 147121052c0cSHong Zhang 1472ebe3b25bSBarry Smith /*MC 1473f68a32c8SEmil Constantinescu TSRK - ODE and DAE solver using Runge-Kutta schemes 1474ebe3b25bSBarry Smith 1475f68a32c8SEmil Constantinescu The user should provide the right hand side of the equation 1476bcf0153eSBarry Smith using `TSSetRHSFunction()`. 1477ebe3b25bSBarry Smith 1478d41222bbSBarry Smith Level: beginner 1479d41222bbSBarry Smith 1480bcf0153eSBarry Smith Notes: 1481bcf0153eSBarry Smith The default is `TSRK3BS`, it can be changed with `TSRKSetType()` or -ts_rk_type 1482ebe3b25bSBarry Smith 14831cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSRK`, `TSSetType()`, `TSRKSetType()`, `TSRKGetType()`, `TSRK2D`, `TSRK2E`, `TSRK3`, 1484bcf0153eSBarry Smith `TSRK4`, `TSRK5`, `TSRKPRSSP2`, `TSRKBPR3`, `TSRKType`, `TSRKRegister()`, `TSRKSetMultirate()`, `TSRKGetMultirate()`, `TSType` 1485ebe3b25bSBarry Smith M*/ 1486d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts) 1487d71ae5a4SJacob Faibussowitsch { 14881566a47fSLisandro Dalcin TS_RK *rk; 1489e4dd521cSBarry Smith 1490e4dd521cSBarry Smith PetscFunctionBegin; 14919566063dSJacob Faibussowitsch PetscCall(TSRKInitializePackage()); 1492f68a32c8SEmil Constantinescu 1493f68a32c8SEmil Constantinescu ts->ops->reset = TSReset_RK; 14945f70b478SJed Brown ts->ops->destroy = TSDestroy_RK; 14955f70b478SJed Brown ts->ops->view = TSView_RK; 1496f68a32c8SEmil Constantinescu ts->ops->load = TSLoad_RK; 1497f68a32c8SEmil Constantinescu ts->ops->setup = TSSetUp_RK; 1498f68a32c8SEmil Constantinescu ts->ops->interpolate = TSInterpolate_RK; 14992c9cb062Sluzhanghpp ts->ops->step = TSStep_RK; 1500f68a32c8SEmil Constantinescu ts->ops->evaluatestep = TSEvaluateStep_RK; 1501fecfb714SLisandro Dalcin ts->ops->rollback = TSRollBack_RK; 1502f68a32c8SEmil Constantinescu ts->ops->setfromoptions = TSSetFromOptions_RK; 1503ff22ae23SHong Zhang ts->ops->getstages = TSGetStages_RK; 1504e4dd521cSBarry Smith 15053ca0882eSHong Zhang ts->ops->snesfunction = SNESTSFormFunction_RK; 15063ca0882eSHong Zhang ts->ops->snesjacobian = SNESTSFormJacobian_RK; 15070f7a1166SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_RK; 150813af1a74SHong Zhang ts->ops->adjointsetup = TSAdjointSetUp_RK; 150913af1a74SHong Zhang ts->ops->adjointstep = TSAdjointStep_RK; 151013af1a74SHong Zhang ts->ops->adjointreset = TSAdjointReset_RK; 15110f7a1166SHong Zhang 151213af1a74SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_RK; 1513922a638cSHong Zhang ts->ops->forwardsetup = TSForwardSetUp_RK; 1514922a638cSHong Zhang ts->ops->forwardreset = TSForwardReset_RK; 1515922a638cSHong Zhang ts->ops->forwardstep = TSForwardStep_RK; 1516922a638cSHong Zhang ts->ops->forwardgetstages = TSForwardGetStages_RK; 1517922a638cSHong Zhang 15184dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&rk)); 15191566a47fSLisandro Dalcin ts->data = (void *)rk; 1520e4dd521cSBarry Smith 15219566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetOrder_C", TSRKGetOrder_RK)); 15229566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetType_C", TSRKGetType_RK)); 15239566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKSetType_C", TSRKSetType_RK)); 15249566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetTableau_C", TSRKGetTableau_RK)); 15259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKSetMultirate_C", TSRKSetMultirate_RK)); 15269566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetMultirate_C", TSRKGetMultirate_RK)); 1527be5899b3SLisandro Dalcin 15289566063dSJacob Faibussowitsch PetscCall(TSRKSetType(ts, TSRKDefault)); 152990b67129SHong Zhang rk->dtratio = 1; 15303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1531e4dd521cSBarry Smith } 1532