1e4dd521cSBarry Smith /* 2474dd773SHong Zhang Code for time stepping with the Runge-Kutta method 3f68a32c8SEmil Constantinescu 4f68a32c8SEmil Constantinescu Notes: 5474dd773SHong Zhang The general system is written as 6474dd773SHong Zhang 7474dd773SHong Zhang Udot = F(t,U) 8474dd773SHong Zhang 9e4dd521cSBarry Smith */ 102c9cb062Sluzhanghpp 11af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 12f68a32c8SEmil Constantinescu #include <petscdm.h> 13474dd773SHong Zhang #include <../src/ts/impls/explicit/rk/rk.h> 1421052c0cSHong Zhang #include <../src/ts/impls/explicit/rk/mrk.h> 15f68a32c8SEmil Constantinescu 16484bcdc7SDebojyoti Ghosh static TSRKType TSRKDefault = TSRK3BS; 17f68a32c8SEmil Constantinescu static PetscBool TSRKRegisterAllCalled; 18f68a32c8SEmil Constantinescu static PetscBool TSRKPackageInitialized; 19f68a32c8SEmil Constantinescu 20ab8985f5SHong Zhang static RKTableauLink RKTableauList; 21ab8985f5SHong Zhang 22f68a32c8SEmil Constantinescu /*MC 23287c2655SBarry Smith TSRK1FE - First order forward Euler scheme. 24262ff7bbSBarry Smith 25f68a32c8SEmil Constantinescu This method has one stage. 26f68a32c8SEmil Constantinescu 27287c2655SBarry Smith Options database: 2867b8a455SSatish Balay . -ts_rk_type 1fe - use type 1fe 29287c2655SBarry Smith 30f68a32c8SEmil Constantinescu Level: advanced 31f68a32c8SEmil Constantinescu 32db781477SPatrick Sanan .seealso: `TSRK`, `TSRKType`, `TSRKSetType()` 33f68a32c8SEmil Constantinescu M*/ 34f68a32c8SEmil Constantinescu /*MC 35e7685601SHong Zhang TSRK2A - Second order RK scheme (Heun's method). 36f68a32c8SEmil Constantinescu 37f68a32c8SEmil Constantinescu This method has two stages. 38f68a32c8SEmil Constantinescu 39287c2655SBarry Smith Options database: 4067b8a455SSatish Balay . -ts_rk_type 2a - use type 2a 41287c2655SBarry Smith 42f68a32c8SEmil Constantinescu Level: advanced 43f68a32c8SEmil Constantinescu 44db781477SPatrick Sanan .seealso: `TSRK`, `TSRKType`, `TSRKSetType()` 45f68a32c8SEmil Constantinescu M*/ 46f68a32c8SEmil Constantinescu /*MC 47e7685601SHong Zhang TSRK2B - Second order RK scheme (the midpoint method). 48e7685601SHong Zhang 49e7685601SHong Zhang This method has two stages. 50e7685601SHong Zhang 51e7685601SHong Zhang Options database: 5267b8a455SSatish Balay . -ts_rk_type 2b - use type 2b 53e7685601SHong Zhang 54e7685601SHong Zhang Level: advanced 55e7685601SHong Zhang 56db781477SPatrick Sanan .seealso: `TSRK`, `TSRKType`, `TSRKSetType()` 57e7685601SHong Zhang M*/ 58e7685601SHong Zhang /*MC 59f68a32c8SEmil Constantinescu TSRK3 - Third order RK scheme. 60f68a32c8SEmil Constantinescu 61f68a32c8SEmil Constantinescu This method has three stages. 62f68a32c8SEmil Constantinescu 63287c2655SBarry Smith Options database: 6467b8a455SSatish Balay . -ts_rk_type 3 - use type 3 65287c2655SBarry Smith 66f68a32c8SEmil Constantinescu Level: advanced 67f68a32c8SEmil Constantinescu 68db781477SPatrick Sanan .seealso: `TSRK`, `TSRKType`, `TSRKSetType()` 69f68a32c8SEmil Constantinescu M*/ 70f68a32c8SEmil Constantinescu /*MC 712109b73fSDebojyoti Ghosh TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method. 722109b73fSDebojyoti Ghosh 73d1905564SLisandro Dalcin This method has four stages with the First Same As Last (FSAL) property. 742109b73fSDebojyoti Ghosh 75287c2655SBarry Smith Options database: 7667b8a455SSatish Balay . -ts_rk_type 3bs - use type 3bs 77287c2655SBarry Smith 782109b73fSDebojyoti Ghosh Level: advanced 792109b73fSDebojyoti Ghosh 80606c0280SSatish Balay References: 81606c0280SSatish Balay . * - https://doi.org/10.1016/0893-9659(89)90079-7 82d1905564SLisandro Dalcin 83db781477SPatrick Sanan .seealso: `TSRK`, `TSRKType`, `TSRKSetType()` 842109b73fSDebojyoti Ghosh M*/ 852109b73fSDebojyoti Ghosh /*MC 86f68a32c8SEmil Constantinescu TSRK4 - Fourth order RK scheme. 87f68a32c8SEmil Constantinescu 882e698f8bSDebojyoti Ghosh This is the classical Runge-Kutta method with four stages. 89f68a32c8SEmil Constantinescu 90287c2655SBarry Smith Options database: 9167b8a455SSatish Balay . -ts_rk_type 4 - use type 4 92287c2655SBarry Smith 93f68a32c8SEmil Constantinescu Level: advanced 94f68a32c8SEmil Constantinescu 95db781477SPatrick Sanan .seealso: `TSRK`, `TSRKType`, `TSRKSetType()` 96f68a32c8SEmil Constantinescu M*/ 97f68a32c8SEmil Constantinescu /*MC 982e698f8bSDebojyoti Ghosh TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method. 99f68a32c8SEmil Constantinescu 100f68a32c8SEmil Constantinescu This method has six stages. 101f68a32c8SEmil Constantinescu 102287c2655SBarry Smith Options database: 10367b8a455SSatish Balay . -ts_rk_type 5f - use type 5f 104287c2655SBarry Smith 105f68a32c8SEmil Constantinescu Level: advanced 106f68a32c8SEmil Constantinescu 107db781477SPatrick Sanan .seealso: `TSRK`, `TSRKType`, `TSRKSetType()` 108f68a32c8SEmil Constantinescu M*/ 1092109b73fSDebojyoti Ghosh /*MC 1102e698f8bSDebojyoti Ghosh TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method. 1112109b73fSDebojyoti Ghosh 112d1905564SLisandro Dalcin This method has seven stages with the First Same As Last (FSAL) property. 1132109b73fSDebojyoti Ghosh 114287c2655SBarry Smith Options database: 11567b8a455SSatish Balay . -ts_rk_type 5dp - use type 5dp 116287c2655SBarry Smith 1172109b73fSDebojyoti Ghosh Level: advanced 1182109b73fSDebojyoti Ghosh 119606c0280SSatish Balay References: 120606c0280SSatish Balay . * - https://doi.org/10.1016/0771-050X(80)90013-3 121d1905564SLisandro Dalcin 122db781477SPatrick Sanan .seealso: `TSRK`, `TSRKType`, `TSRKSetType()` 1232109b73fSDebojyoti Ghosh M*/ 12405e23783SLisandro Dalcin /*MC 12505e23783SLisandro Dalcin TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method. 12605e23783SLisandro Dalcin 12705e23783SLisandro Dalcin This method has eight stages with the First Same As Last (FSAL) property. 12805e23783SLisandro Dalcin 129287c2655SBarry Smith Options database: 13067b8a455SSatish Balay . -ts_rk_type 5bs - use type 5bs 131287c2655SBarry Smith 13205e23783SLisandro Dalcin Level: advanced 13305e23783SLisandro Dalcin 134606c0280SSatish Balay References: 135606c0280SSatish Balay . * - https://doi.org/10.1016/0898-1221(96)00141-1 13605e23783SLisandro Dalcin 137db781477SPatrick Sanan .seealso: `TSRK`, `TSRKType`, `TSRKSetType()` 13805e23783SLisandro Dalcin M*/ 13937acaa02SHendrik Ranocha /*MC 14037acaa02SHendrik Ranocha TSRK6VR - Sixth order robust Verner RK scheme with fifth order embedded method. 14137acaa02SHendrik Ranocha 14237acaa02SHendrik Ranocha This method has nine stages with the First Same As Last (FSAL) property. 14337acaa02SHendrik Ranocha 14437acaa02SHendrik Ranocha Options database: 14567b8a455SSatish Balay . -ts_rk_type 6vr - use type 6vr 14637acaa02SHendrik Ranocha 14737acaa02SHendrik Ranocha Level: advanced 14837acaa02SHendrik Ranocha 149606c0280SSatish Balay References: 150606c0280SSatish Balay . * - http://people.math.sfu.ca/~jverner/RKV65.IIIXb.Robust.00010102836.081204.CoeffsOnlyRAT 15137acaa02SHendrik Ranocha 152db781477SPatrick Sanan .seealso: `TSRK`, `TSRKType`, `TSRKSetType()` 15337acaa02SHendrik Ranocha M*/ 15437acaa02SHendrik Ranocha /*MC 15537acaa02SHendrik Ranocha TSRK7VR - Seventh order robust Verner RK scheme with sixth order embedded method. 15637acaa02SHendrik Ranocha 1578f3d7ee2SCarsten Uphoff This method has ten stages. 15837acaa02SHendrik Ranocha 15937acaa02SHendrik Ranocha Options database: 16067b8a455SSatish Balay . -ts_rk_type 7vr - use type 7vr 16137acaa02SHendrik Ranocha 16237acaa02SHendrik Ranocha Level: advanced 16337acaa02SHendrik Ranocha 164606c0280SSatish Balay References: 165606c0280SSatish Balay . * - http://people.math.sfu.ca/~jverner/RKV76.IIa.Robust.000027015646.081206.CoeffsOnlyRAT 16637acaa02SHendrik Ranocha 167db781477SPatrick Sanan .seealso: `TSRK`, `TSRKType`, `TSRKSetType()` 16837acaa02SHendrik Ranocha M*/ 16937acaa02SHendrik Ranocha /*MC 17037acaa02SHendrik Ranocha TSRK8VR - Eigth order robust Verner RK scheme with seventh order embedded method. 17137acaa02SHendrik Ranocha 1728f3d7ee2SCarsten Uphoff This method has thirteen stages. 17337acaa02SHendrik Ranocha 17437acaa02SHendrik Ranocha Options database: 17567b8a455SSatish Balay . -ts_rk_type 8vr - use type 8vr 17637acaa02SHendrik Ranocha 17737acaa02SHendrik Ranocha Level: advanced 17837acaa02SHendrik Ranocha 179606c0280SSatish Balay References: 180606c0280SSatish Balay . * - http://people.math.sfu.ca/~jverner/RKV87.IIa.Robust.00000754677.081208.CoeffsOnlyRATandFLOAT 18137acaa02SHendrik Ranocha 182db781477SPatrick Sanan .seealso: `TSRK`, `TSRKType`, `TSRKSetType()` 18337acaa02SHendrik Ranocha M*/ 184262ff7bbSBarry Smith 185f68a32c8SEmil Constantinescu /*@C 186f68a32c8SEmil Constantinescu TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK 187262ff7bbSBarry Smith 188f68a32c8SEmil Constantinescu Not Collective, but should be called by all processes which will need the schemes to be registered 189262ff7bbSBarry Smith 190f68a32c8SEmil Constantinescu Level: advanced 191262ff7bbSBarry Smith 192db781477SPatrick Sanan .seealso: `TSRKRegisterDestroy()` 193262ff7bbSBarry Smith @*/ 1949371c9d4SSatish Balay PetscErrorCode TSRKRegisterAll(void) { 195262ff7bbSBarry Smith PetscFunctionBegin; 196f68a32c8SEmil Constantinescu if (TSRKRegisterAllCalled) PetscFunctionReturn(0); 197f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_TRUE; 198e4dd521cSBarry Smith 1994e82626cSLisandro Dalcin #define RC PetscRealConstant 200e4dd521cSBarry Smith { 2019371c9d4SSatish Balay const PetscReal A[1][1] = {{0}}, b[1] = {RC(1.0)}; 2029566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK1FE, 1, 1, &A[0][0], b, NULL, NULL, 0, NULL)); 203e4dd521cSBarry Smith } 204db046809SBarry Smith { 2059371c9d4SSatish Balay const PetscReal A[2][2] = 2069371c9d4SSatish Balay { 2079371c9d4SSatish Balay {0, 0}, 2089371c9d4SSatish Balay {RC(1.0), 0} 2099371c9d4SSatish Balay }, 2109371c9d4SSatish Balay b[2] = {RC(0.5), RC(0.5)}, bembed[2] = {RC(1.0), 0}; 2119566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK2A, 2, 2, &A[0][0], b, NULL, bembed, 0, NULL)); 212db046809SBarry Smith } 213f68a32c8SEmil Constantinescu { 2149371c9d4SSatish Balay const PetscReal A[2][2] = 2159371c9d4SSatish Balay { 2169371c9d4SSatish Balay {0, 0}, 2179371c9d4SSatish Balay {RC(0.5), 0} 2189371c9d4SSatish Balay }, 2199958aef7SJacob Faibussowitsch b[2] = {0, RC(1.0)}; 2209566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK2B, 2, 2, &A[0][0], b, NULL, NULL, 0, NULL)); 221e7685601SHong Zhang } 222e7685601SHong Zhang { 2239371c9d4SSatish Balay const PetscReal A[3][3] = 2249371c9d4SSatish Balay { 2259371c9d4SSatish Balay {0, 0, 0}, 2264e82626cSLisandro Dalcin {RC(2.0) / RC(3.0), 0, 0}, 2279371c9d4SSatish Balay {RC(-1.0) / RC(3.0), RC(1.0), 0} 2289371c9d4SSatish Balay }, 2294e82626cSLisandro Dalcin b[3] = {RC(0.25), RC(0.5), RC(0.25)}; 2309566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK3, 3, 3, &A[0][0], b, NULL, NULL, 0, NULL)); 231fdefd5e5SDebojyoti Ghosh } 232fdefd5e5SDebojyoti Ghosh { 2339371c9d4SSatish Balay const PetscReal A[4][4] = 2349371c9d4SSatish Balay { 2359371c9d4SSatish Balay {0, 0, 0, 0}, 2364e82626cSLisandro Dalcin {RC(1.0) / RC(2.0), 0, 0, 0}, 2374e82626cSLisandro Dalcin {0, RC(3.0) / RC(4.0), 0, 0}, 2389371c9d4SSatish Balay {RC(2.0) / RC(9.0), RC(1.0) / RC(3.0), RC(4.0) / RC(9.0), 0} 2399371c9d4SSatish Balay }, 2409371c9d4SSatish 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)}; 2419566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK3BS, 3, 4, &A[0][0], b, NULL, bembed, 0, NULL)); 242db046809SBarry Smith } 243f68a32c8SEmil Constantinescu { 2449371c9d4SSatish Balay const PetscReal A[4][4] = 2459371c9d4SSatish Balay { 2469371c9d4SSatish Balay {0, 0, 0, 0}, 2474e82626cSLisandro Dalcin {RC(0.5), 0, 0, 0}, 2484e82626cSLisandro Dalcin {0, RC(0.5), 0, 0}, 2499371c9d4SSatish Balay {0, 0, RC(1.0), 0} 2509371c9d4SSatish Balay }, 2514e82626cSLisandro 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)}; 2529566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK4, 4, 4, &A[0][0], b, NULL, NULL, 0, NULL)); 253db046809SBarry Smith } 254f68a32c8SEmil Constantinescu { 2559371c9d4SSatish Balay const PetscReal A[6][6] = 2569371c9d4SSatish Balay { 2579371c9d4SSatish Balay {0, 0, 0, 0, 0, 0}, 2584e82626cSLisandro Dalcin {RC(0.25), 0, 0, 0, 0, 0}, 2594e82626cSLisandro Dalcin {RC(3.0) / RC(32.0), RC(9.0) / RC(32.0), 0, 0, 0, 0}, 2604e82626cSLisandro Dalcin {RC(1932.0) / RC(2197.0), RC(-7200.0) / RC(2197.0), RC(7296.0) / RC(2197.0), 0, 0, 0}, 2614e82626cSLisandro Dalcin {RC(439.0) / RC(216.0), RC(-8.0), RC(3680.0) / RC(513.0), RC(-845.0) / RC(4104.0), 0, 0}, 2629371c9d4SSatish 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} 2639371c9d4SSatish Balay }, 2644e82626cSLisandro 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)}, 2654e82626cSLisandro 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}; 2669566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK5F, 5, 6, &A[0][0], b, NULL, bembed, 0, NULL)); 267fdefd5e5SDebojyoti Ghosh } 268fdefd5e5SDebojyoti Ghosh { 2699371c9d4SSatish Balay const PetscReal A[7][7] = 2709371c9d4SSatish Balay { 2719371c9d4SSatish Balay {0, 0, 0, 0, 0, 0, 0}, 2724e82626cSLisandro Dalcin {RC(1.0) / RC(5.0), 0, 0, 0, 0, 0, 0}, 2734e82626cSLisandro Dalcin {RC(3.0) / RC(40.0), RC(9.0) / RC(40.0), 0, 0, 0, 0, 0}, 2744e82626cSLisandro Dalcin {RC(44.0) / RC(45.0), RC(-56.0) / RC(15.0), RC(32.0) / RC(9.0), 0, 0, 0, 0}, 2754e82626cSLisandro 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}, 2764e82626cSLisandro 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}, 2779371c9d4SSatish 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} 2789371c9d4SSatish Balay }, 2794e82626cSLisandro 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}, 2809371c9d4SSatish 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)}}; 2819566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK5DP, 5, 7, &A[0][0], b, NULL, bembed, 5, binterp[0])); 282f68a32c8SEmil Constantinescu } 28305e23783SLisandro Dalcin { 2849371c9d4SSatish Balay const PetscReal A[8][8] = 2859371c9d4SSatish Balay { 2869371c9d4SSatish Balay {0, 0, 0, 0, 0, 0, 0, 0}, 2874e82626cSLisandro Dalcin {RC(1.0) / RC(6.0), 0, 0, 0, 0, 0, 0, 0}, 2884e82626cSLisandro Dalcin {RC(2.0) / RC(27.0), RC(4.0) / RC(27.0), 0, 0, 0, 0, 0, 0}, 2894e82626cSLisandro Dalcin {RC(183.0) / RC(1372.0), RC(-162.0) / RC(343.0), RC(1053.0) / RC(1372.0), 0, 0, 0, 0, 0}, 2904e82626cSLisandro 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}, 2914e82626cSLisandro 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}, 2924e82626cSLisandro 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}, 2939371c9d4SSatish 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} 2949371c9d4SSatish Balay }, 2954e82626cSLisandro 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}, 2964e82626cSLisandro 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)}; 2979566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK5BS, 5, 8, &A[0][0], b, NULL, bembed, 0, NULL)); 29805e23783SLisandro Dalcin } 29937acaa02SHendrik Ranocha { 3009371c9d4SSatish Balay const PetscReal A[9][9] = 3019371c9d4SSatish Balay { 3029371c9d4SSatish Balay {0, 0, 0, 0, 0, 0, 0, 0, 0}, 30363fe90e8SHendrik Ranocha {RC(1.8000000000000000000000000000000000000000e-01), 0, 0, 0, 0, 0, 0, 0, 0}, 30463fe90e8SHendrik Ranocha {RC(8.9506172839506172839506172839506172839506e-02), RC(7.7160493827160493827160493827160493827160e-02), 0, 0, 0, 0, 0, 0, 0}, 30563fe90e8SHendrik Ranocha {RC(6.2500000000000000000000000000000000000000e-02), 0, RC(1.8750000000000000000000000000000000000000e-01), 0, 0, 0, 0, 0, 0}, 30663fe90e8SHendrik Ranocha {RC(3.1651600000000000000000000000000000000000e-01), 0, RC(-1.0449480000000000000000000000000000000000e+00), RC(1.2584320000000000000000000000000000000000e+00), 0, 0, 0, 0, 0}, 30763fe90e8SHendrik Ranocha {RC(2.7232612736485626257225065566674305502508e-01), 0, RC(-8.2513360323886639676113360323886639676113e-01), RC(1.0480917678812415654520917678812415654521e+00), RC(1.0471570799276856873679117969088177628396e-01), 0, 0, 0, 0}, 30863fe90e8SHendrik 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}, 30963fe90e8SHendrik 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}, 3109371c9d4SSatish 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} 3119371c9d4SSatish Balay }, 3129371c9d4SSatish Balay b[9] = {RC(7.6388888888888888888888888888888888888889e-02), 0, 0, RC(3.6940836940836940836940836940836940836941e-01), 0, RC(2.4801587301587301587301587301587301587302e-01), RC(2.3674242424242424242424242424242424242424e-01), 3139371c9d4SSatish Balay RC(6.9444444444444444444444444444444444444444e-02), 0}, 31463fe90e8SHendrik 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)}; 3159566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK6VR, 6, 9, &A[0][0], b, NULL, bembed, 0, NULL)); 31637acaa02SHendrik Ranocha } 31737acaa02SHendrik Ranocha { 3189371c9d4SSatish Balay const PetscReal A[10][10] = 3199371c9d4SSatish Balay { 3209371c9d4SSatish Balay {0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, 32163fe90e8SHendrik Ranocha {RC(5.0000000000000000000000000000000000000000e-03), 0, 0, 0, 0, 0, 0, 0, 0, 0}, 32263fe90e8SHendrik Ranocha {RC(-1.0767901234567901234567901234567901234568e+00), RC(1.1856790123456790123456790123456790123457e+00), 0, 0, 0, 0, 0, 0, 0, 0}, 32363fe90e8SHendrik Ranocha {RC(4.0833333333333333333333333333333333333333e-02), 0, RC(1.2250000000000000000000000000000000000000e-01), 0, 0, 0, 0, 0, 0, 0}, 32463fe90e8SHendrik Ranocha {RC(6.3607142857142857142857142857142857142857e-01), 0, RC(-2.4444642857142857142857142857142857142857e+00), RC(2.2633928571428571428571428571428571428571e+00), 0, 0, 0, 0, 0, 0}, 32563fe90e8SHendrik Ranocha {RC(-2.5351211079349245229256383554660215487207e+00), 0, RC(1.0299374654449267920438514460756024913612e+01), RC(-7.9513032885990579949493217458266876536482e+00), RC(7.9301148923100592201226014271115261823800e-01), 0, 0, 0, 0, 0}, 32663fe90e8SHendrik 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}, 32763fe90e8SHendrik 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}, 32863fe90e8SHendrik 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}, 3299371c9d4SSatish 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} 3309371c9d4SSatish Balay }, 3319371c9d4SSatish 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), 3329371c9d4SSatish Balay RC(2.5204468723743860507184043820197442562182e-01), 0, 0, RC(4.8834971521418614557381971303093137592592e-02)}; 3339566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK7VR, 7, 10, &A[0][0], b, NULL, bembed, 0, NULL)); 33437acaa02SHendrik Ranocha } 33537acaa02SHendrik Ranocha { 3369371c9d4SSatish Balay const PetscReal A[13][13] = 3379371c9d4SSatish Balay { 3389371c9d4SSatish Balay {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, 33963fe90e8SHendrik Ranocha {RC(2.5000000000000000000000000000000000000000e-01), 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, 34063fe90e8SHendrik Ranocha {RC(8.7400846504915232052686327594877411977046e-02), RC(2.5487604938654321753087950620345685135815e-02), 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, 34163fe90e8SHendrik Ranocha {RC(4.2333169291338582677165354330708661417323e-02), 0, RC(1.2699950787401574803149606299212598425197e-01), 0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, 34263fe90e8SHendrik Ranocha {RC(4.2609505888742261494881445237572274090942e-01), 0, RC(-1.5987952846591523265427733230657181117089e+00), RC(1.5967002257717297115939588706899953707994e+00), 0, 0, 0, 0, 0, 0, 0, 0, 0}, 34363fe90e8SHendrik Ranocha {RC(5.0719337296713929515090618138513639239329e-02), 0, 0, RC(2.5433377264600407582754714408877778031369e-01), RC(2.0394689005728199465736223777270858044698e-01), 0, 0, 0, 0, 0, 0, 0, 0}, 34463fe90e8SHendrik 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}, 34563fe90e8SHendrik 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}, 34663fe90e8SHendrik 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}, 34763fe90e8SHendrik 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}, 34863fe90e8SHendrik 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}, 34963fe90e8SHendrik 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}, 3509371c9d4SSatish 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} 3519371c9d4SSatish Balay }, 3529371c9d4SSatish 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)}; 3539566063dSJacob Faibussowitsch PetscCall(TSRKRegister(TSRK8VR, 8, 13, &A[0][0], b, NULL, bembed, 0, NULL)); 35437acaa02SHendrik Ranocha } 3554e82626cSLisandro Dalcin #undef RC 356db046809SBarry Smith PetscFunctionReturn(0); 357db046809SBarry Smith } 358db046809SBarry Smith 359f68a32c8SEmil Constantinescu /*@C 360f68a32c8SEmil Constantinescu TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister(). 361f68a32c8SEmil Constantinescu 362f68a32c8SEmil Constantinescu Not Collective 363f68a32c8SEmil Constantinescu 364f68a32c8SEmil Constantinescu Level: advanced 365f68a32c8SEmil Constantinescu 366db781477SPatrick Sanan .seealso: `TSRKRegister()`, `TSRKRegisterAll()` 367f68a32c8SEmil Constantinescu @*/ 3689371c9d4SSatish Balay PetscErrorCode TSRKRegisterDestroy(void) { 369f68a32c8SEmil Constantinescu RKTableauLink link; 370f68a32c8SEmil Constantinescu 371f68a32c8SEmil Constantinescu PetscFunctionBegin; 372f68a32c8SEmil Constantinescu while ((link = RKTableauList)) { 373f68a32c8SEmil Constantinescu RKTableau t = &link->tab; 374f68a32c8SEmil Constantinescu RKTableauList = link->next; 3759566063dSJacob Faibussowitsch PetscCall(PetscFree3(t->A, t->b, t->c)); 3769566063dSJacob Faibussowitsch PetscCall(PetscFree(t->bembed)); 3779566063dSJacob Faibussowitsch PetscCall(PetscFree(t->binterp)); 3789566063dSJacob Faibussowitsch PetscCall(PetscFree(t->name)); 3799566063dSJacob Faibussowitsch PetscCall(PetscFree(link)); 380f68a32c8SEmil Constantinescu } 381f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_FALSE; 382f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 383f68a32c8SEmil Constantinescu } 384f68a32c8SEmil Constantinescu 385f68a32c8SEmil Constantinescu /*@C 386f68a32c8SEmil Constantinescu TSRKInitializePackage - This function initializes everything in the TSRK package. It is called 3878a690491SBarry Smith from TSInitializePackage(). 388f68a32c8SEmil Constantinescu 389f68a32c8SEmil Constantinescu Level: developer 390f68a32c8SEmil Constantinescu 391db781477SPatrick Sanan .seealso: `PetscInitialize()` 392f68a32c8SEmil Constantinescu @*/ 3939371c9d4SSatish Balay PetscErrorCode TSRKInitializePackage(void) { 394e4dd521cSBarry Smith PetscFunctionBegin; 395f68a32c8SEmil Constantinescu if (TSRKPackageInitialized) PetscFunctionReturn(0); 396f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_TRUE; 3979566063dSJacob Faibussowitsch PetscCall(TSRKRegisterAll()); 3989566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSRKFinalizePackage)); 399f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 400f68a32c8SEmil Constantinescu } 401186e87acSLisandro Dalcin 402f68a32c8SEmil Constantinescu /*@C 403f68a32c8SEmil Constantinescu TSRKFinalizePackage - This function destroys everything in the TSRK package. It is 404f68a32c8SEmil Constantinescu called from PetscFinalize(). 405186e87acSLisandro Dalcin 406f68a32c8SEmil Constantinescu Level: developer 407f68a32c8SEmil Constantinescu 408db781477SPatrick Sanan .seealso: `PetscFinalize()` 409f68a32c8SEmil Constantinescu @*/ 4109371c9d4SSatish Balay PetscErrorCode TSRKFinalizePackage(void) { 411f68a32c8SEmil Constantinescu PetscFunctionBegin; 412f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_FALSE; 4139566063dSJacob Faibussowitsch PetscCall(TSRKRegisterDestroy()); 414f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 415f68a32c8SEmil Constantinescu } 416f68a32c8SEmil Constantinescu 417f68a32c8SEmil Constantinescu /*@C 418f68a32c8SEmil Constantinescu TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 419f68a32c8SEmil Constantinescu 420f68a32c8SEmil Constantinescu Not Collective, but the same schemes should be registered on all processes on which they will be used 421f68a32c8SEmil Constantinescu 422f68a32c8SEmil Constantinescu Input Parameters: 423f68a32c8SEmil Constantinescu + name - identifier for method 424f68a32c8SEmil Constantinescu . order - approximation order of method 425f68a32c8SEmil Constantinescu . s - number of stages, this is the dimension of the matrices below 426f68a32c8SEmil Constantinescu . A - stage coefficients (dimension s*s, row-major) 427f68a32c8SEmil Constantinescu . b - step completion table (dimension s; NULL to use last row of A) 428f68a32c8SEmil Constantinescu . c - abscissa (dimension s; NULL to use row sums of A) 429f68a32c8SEmil Constantinescu . bembed - completion table for embedded method (dimension s; NULL if not available) 4303a8a9803SLisandro Dalcin . p - Order of the interpolation scheme, equal to the number of columns of binterp 4313a8a9803SLisandro Dalcin - binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1) 432f68a32c8SEmil Constantinescu 433f68a32c8SEmil Constantinescu Notes: 434f68a32c8SEmil Constantinescu Several RK methods are provided, this function is only needed to create new methods. 435f68a32c8SEmil Constantinescu 436f68a32c8SEmil Constantinescu Level: advanced 437f68a32c8SEmil Constantinescu 438db781477SPatrick Sanan .seealso: `TSRK` 439f68a32c8SEmil Constantinescu @*/ 4409371c9d4SSatish Balay 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[]) { 441f68a32c8SEmil Constantinescu RKTableauLink link; 442f68a32c8SEmil Constantinescu RKTableau t; 443f68a32c8SEmil Constantinescu PetscInt i, j; 444f68a32c8SEmil Constantinescu 445f68a32c8SEmil Constantinescu PetscFunctionBegin; 4463a8a9803SLisandro Dalcin PetscValidCharPointer(name, 1); 4473a8a9803SLisandro Dalcin PetscValidRealPointer(A, 4); 4483a8a9803SLisandro Dalcin if (b) PetscValidRealPointer(b, 5); 4493a8a9803SLisandro Dalcin if (c) PetscValidRealPointer(c, 6); 4503a8a9803SLisandro Dalcin if (bembed) PetscValidRealPointer(bembed, 7); 4513a8a9803SLisandro Dalcin if (binterp || p > 1) PetscValidRealPointer(binterp, 9); 4523a8a9803SLisandro Dalcin 4539566063dSJacob Faibussowitsch PetscCall(TSRKInitializePackage()); 4549566063dSJacob Faibussowitsch PetscCall(PetscNew(&link)); 455f68a32c8SEmil Constantinescu t = &link->tab; 4563a8a9803SLisandro Dalcin 4579566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &t->name)); 458f68a32c8SEmil Constantinescu t->order = order; 459f68a32c8SEmil Constantinescu t->s = s; 4609566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(s * s, &t->A, s, &t->b, s, &t->c)); 4619566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->A, A, s * s)); 4629566063dSJacob Faibussowitsch if (b) PetscCall(PetscArraycpy(t->b, b, s)); 4639371c9d4SSatish Balay else 4649371c9d4SSatish Balay for (i = 0; i < s; i++) t->b[i] = A[(s - 1) * s + i]; 4659566063dSJacob Faibussowitsch if (c) PetscCall(PetscArraycpy(t->c, c, s)); 4669371c9d4SSatish Balay else 4679371c9d4SSatish Balay for (i = 0; i < s; i++) 4689371c9d4SSatish Balay for (j = 0, t->c[i] = 0; j < s; j++) t->c[i] += A[i * s + j]; 469d760c35bSDebojyoti Ghosh t->FSAL = PETSC_TRUE; 4709371c9d4SSatish Balay for (i = 0; i < s; i++) 4719371c9d4SSatish Balay if (t->A[(s - 1) * s + i] != t->b[i]) t->FSAL = PETSC_FALSE; 4723a8a9803SLisandro Dalcin 473f68a32c8SEmil Constantinescu if (bembed) { 4749566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s, &t->bembed)); 4759566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->bembed, bembed, s)); 476f68a32c8SEmil Constantinescu } 477f68a32c8SEmil Constantinescu 4789371c9d4SSatish Balay if (!binterp) { 4799371c9d4SSatish Balay p = 1; 4809371c9d4SSatish Balay binterp = t->b; 4819371c9d4SSatish Balay } 4823a8a9803SLisandro Dalcin t->p = p; 4839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s * p, &t->binterp)); 4849566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->binterp, binterp, s * p)); 4853a8a9803SLisandro Dalcin 486f68a32c8SEmil Constantinescu link->next = RKTableauList; 487f68a32c8SEmil Constantinescu RKTableauList = link; 488f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 489f68a32c8SEmil Constantinescu } 490f68a32c8SEmil Constantinescu 4919371c9d4SSatish Balay 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) { 4920f7a28cdSStefano Zampini TS_RK *rk = (TS_RK *)ts->data; 4930f7a28cdSStefano Zampini RKTableau tab = rk->tableau; 4940f7a28cdSStefano Zampini 4950f7a28cdSStefano Zampini PetscFunctionBegin; 4960f7a28cdSStefano Zampini if (s) *s = tab->s; 4970f7a28cdSStefano Zampini if (A) *A = tab->A; 4980f7a28cdSStefano Zampini if (b) *b = tab->b; 4990f7a28cdSStefano Zampini if (c) *c = tab->c; 5000f7a28cdSStefano Zampini if (bembed) *bembed = tab->bembed; 5010f7a28cdSStefano Zampini if (p) *p = tab->p; 5020f7a28cdSStefano Zampini if (binterp) *binterp = tab->binterp; 5030f7a28cdSStefano Zampini if (FSAL) *FSAL = tab->FSAL; 5040f7a28cdSStefano Zampini PetscFunctionReturn(0); 5050f7a28cdSStefano Zampini } 5060f7a28cdSStefano Zampini 5070f7a28cdSStefano Zampini /*@C 5080f7a28cdSStefano Zampini TSRKGetTableau - Get info on the RK tableau 5090f7a28cdSStefano Zampini 5100f7a28cdSStefano Zampini Not Collective 5110f7a28cdSStefano Zampini 512f899ff85SJose E. Roman Input Parameter: 5130f7a28cdSStefano Zampini . ts - timestepping context 5140f7a28cdSStefano Zampini 5150f7a28cdSStefano Zampini Output Parameters: 5160f7a28cdSStefano Zampini + s - number of stages, this is the dimension of the matrices below 5170f7a28cdSStefano Zampini . A - stage coefficients (dimension s*s, row-major) 5180f7a28cdSStefano Zampini . b - step completion table (dimension s) 5190f7a28cdSStefano Zampini . c - abscissa (dimension s) 5200f7a28cdSStefano Zampini . bembed - completion table for embedded method (dimension s; NULL if not available) 5210f7a28cdSStefano Zampini . p - Order of the interpolation scheme, equal to the number of columns of binterp 5220f7a28cdSStefano Zampini . binterp - Coefficients of the interpolation formula (dimension s*p) 5230f7a28cdSStefano Zampini - FSAL - wheather or not the scheme has the First Same As Last property 5240f7a28cdSStefano Zampini 5250f7a28cdSStefano Zampini Level: developer 5260f7a28cdSStefano Zampini 527db781477SPatrick Sanan .seealso: `TSRK` 5280f7a28cdSStefano Zampini @*/ 5299371c9d4SSatish Balay 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) { 5300f7a28cdSStefano Zampini PetscFunctionBegin; 5310f7a28cdSStefano Zampini PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 5329371c9d4SSatish 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)); 5330f7a28cdSStefano Zampini PetscFunctionReturn(0); 5340f7a28cdSStefano Zampini } 5350f7a28cdSStefano Zampini 536e4dd521cSBarry Smith /* 5372c9cb062Sluzhanghpp This is for single-step RK method 538f68a32c8SEmil Constantinescu The step completion formula is 539e4dd521cSBarry Smith 540f68a32c8SEmil Constantinescu x1 = x0 + h b^T YdotRHS 541f68a32c8SEmil Constantinescu 542f68a32c8SEmil Constantinescu This function can be called before or after ts->vec_sol has been updated. 543f68a32c8SEmil Constantinescu Suppose we have a completion formula (b) and an embedded formula (be) of different order. 544f68a32c8SEmil Constantinescu We can write 545f68a32c8SEmil Constantinescu 546f68a32c8SEmil Constantinescu x1e = x0 + h be^T YdotRHS 547f68a32c8SEmil Constantinescu = x1 - h b^T YdotRHS + h be^T YdotRHS 548f68a32c8SEmil Constantinescu = x1 + h (be - b)^T YdotRHS 549f68a32c8SEmil Constantinescu 550f68a32c8SEmil Constantinescu so we can evaluate the method with different order even after the step has been optimistically completed. 551e4dd521cSBarry Smith */ 5529371c9d4SSatish Balay static PetscErrorCode TSEvaluateStep_RK(TS ts, PetscInt order, Vec X, PetscBool *done) { 553f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK *)ts->data; 554f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 555f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 556f68a32c8SEmil Constantinescu PetscReal h; 557f68a32c8SEmil Constantinescu PetscInt s = tab->s, j; 558f68a32c8SEmil Constantinescu 559f68a32c8SEmil Constantinescu PetscFunctionBegin; 560f68a32c8SEmil Constantinescu switch (rk->status) { 561f68a32c8SEmil Constantinescu case TS_STEP_INCOMPLETE: 5629371c9d4SSatish Balay case TS_STEP_PENDING: h = ts->time_step; break; 5639371c9d4SSatish Balay case TS_STEP_COMPLETE: h = ts->ptime - ts->ptime_prev; break; 564f68a32c8SEmil Constantinescu default: SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus"); 565f68a32c8SEmil Constantinescu } 566f68a32c8SEmil Constantinescu if (order == tab->order) { 567f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { 5689566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X)); 56990b67129SHong Zhang for (j = 0; j < s; j++) w[j] = h * tab->b[j] / rk->dtratio; 5709566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, rk->YdotRHS)); 5719566063dSJacob Faibussowitsch } else PetscCall(VecCopy(ts->vec_sol, X)); 572f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 573f68a32c8SEmil Constantinescu } else if (order == tab->order - 1) { 574f68a32c8SEmil Constantinescu if (!tab->bembed) goto unavailable; 575f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { /*Complete with the embedded method (be)*/ 5769566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X)); 577f68a32c8SEmil Constantinescu for (j = 0; j < s; j++) w[j] = h * tab->bembed[j]; 5789566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, rk->YdotRHS)); 579f68a32c8SEmil Constantinescu } else { /*Rollback and re-complete using (be-b) */ 5809566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X)); 581f68a32c8SEmil Constantinescu for (j = 0; j < s; j++) w[j] = h * (tab->bembed[j] - tab->b[j]); 5829566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, rk->YdotRHS)); 583f68a32c8SEmil Constantinescu } 584f68a32c8SEmil Constantinescu if (done) *done = PETSC_TRUE; 585f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 586f68a32c8SEmil Constantinescu } 587f68a32c8SEmil Constantinescu unavailable: 588f68a32c8SEmil Constantinescu if (done) *done = PETSC_FALSE; 5899371c9d4SSatish Balay else 5909371c9d4SSatish 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); 591f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 592f68a32c8SEmil Constantinescu } 593f68a32c8SEmil Constantinescu 5949371c9d4SSatish Balay static PetscErrorCode TSForwardCostIntegral_RK(TS ts) { 5950f7a1166SHong Zhang TS_RK *rk = (TS_RK *)ts->data; 596cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 5970f7a1166SHong Zhang RKTableau tab = rk->tableau; 5980f7a1166SHong Zhang const PetscInt s = tab->s; 5990f7a1166SHong Zhang const PetscReal *b = tab->b, *c = tab->c; 6000f7a1166SHong Zhang Vec *Y = rk->Y; 6010f7a1166SHong Zhang PetscInt i; 6020f7a1166SHong Zhang 6030f7a1166SHong Zhang PetscFunctionBegin; 604cd4cee2dSHong Zhang /* No need to backup quadts->vec_sol since it can be reverted in TSRollBack_RK */ 6050f7a1166SHong Zhang for (i = s - 1; i >= 0; i--) { 606cd4cee2dSHong Zhang /* Evolve quadrature TS solution to compute integrals */ 6079566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(quadts, rk->ptime + rk->time_step * c[i], Y[i], ts->vec_costintegrand)); 6089566063dSJacob Faibussowitsch PetscCall(VecAXPY(quadts->vec_sol, rk->time_step * b[i], ts->vec_costintegrand)); 6090f7a1166SHong Zhang } 6100f7a1166SHong Zhang PetscFunctionReturn(0); 6110f7a1166SHong Zhang } 6120f7a1166SHong Zhang 6139371c9d4SSatish Balay static PetscErrorCode TSAdjointCostIntegral_RK(TS ts) { 6140f7a1166SHong Zhang TS_RK *rk = (TS_RK *)ts->data; 6150f7a1166SHong Zhang RKTableau tab = rk->tableau; 616cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 6170f7a1166SHong Zhang const PetscInt s = tab->s; 6180f7a1166SHong Zhang const PetscReal *b = tab->b, *c = tab->c; 6190f7a1166SHong Zhang Vec *Y = rk->Y; 6200f7a1166SHong Zhang PetscInt i; 6210f7a1166SHong Zhang 6220f7a1166SHong Zhang PetscFunctionBegin; 6230f7a1166SHong Zhang for (i = s - 1; i >= 0; i--) { 624cd4cee2dSHong Zhang /* Evolve quadrature TS solution to compute integrals */ 6259566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(quadts, ts->ptime + ts->time_step * (1.0 - c[i]), Y[i], ts->vec_costintegrand)); 6269566063dSJacob Faibussowitsch PetscCall(VecAXPY(quadts->vec_sol, -ts->time_step * b[i], ts->vec_costintegrand)); 6270f7a1166SHong Zhang } 6280f7a1166SHong Zhang PetscFunctionReturn(0); 6290f7a1166SHong Zhang } 6300f7a1166SHong Zhang 6319371c9d4SSatish Balay static PetscErrorCode TSRollBack_RK(TS ts) { 632fecfb714SLisandro Dalcin TS_RK *rk = (TS_RK *)ts->data; 633cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 634fecfb714SLisandro Dalcin RKTableau tab = rk->tableau; 635fecfb714SLisandro Dalcin const PetscInt s = tab->s; 636cd4cee2dSHong Zhang const PetscReal *b = tab->b, *c = tab->c; 637fecfb714SLisandro Dalcin PetscScalar *w = rk->work; 638cd4cee2dSHong Zhang Vec *Y = rk->Y, *YdotRHS = rk->YdotRHS; 639fecfb714SLisandro Dalcin PetscInt j; 640fecfb714SLisandro Dalcin PetscReal h; 641fecfb714SLisandro Dalcin 642fecfb714SLisandro Dalcin PetscFunctionBegin; 643fecfb714SLisandro Dalcin switch (rk->status) { 644fecfb714SLisandro Dalcin case TS_STEP_INCOMPLETE: 6459371c9d4SSatish Balay case TS_STEP_PENDING: h = ts->time_step; break; 6469371c9d4SSatish Balay case TS_STEP_COMPLETE: h = ts->ptime - ts->ptime_prev; break; 647fecfb714SLisandro Dalcin default: SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus"); 648fecfb714SLisandro Dalcin } 649fecfb714SLisandro Dalcin for (j = 0; j < s; j++) w[j] = -h * b[j]; 6509566063dSJacob Faibussowitsch PetscCall(VecMAXPY(ts->vec_sol, s, w, YdotRHS)); 651cd4cee2dSHong Zhang if (quadts && ts->costintegralfwd) { 652cd4cee2dSHong Zhang for (j = 0; j < s; j++) { 653cd4cee2dSHong Zhang /* Revert the quadrature TS solution */ 6549566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(quadts, rk->ptime + h * c[j], Y[j], ts->vec_costintegrand)); 6559566063dSJacob Faibussowitsch PetscCall(VecAXPY(quadts->vec_sol, -h * b[j], ts->vec_costintegrand)); 656cd4cee2dSHong Zhang } 657cd4cee2dSHong Zhang } 658fecfb714SLisandro Dalcin PetscFunctionReturn(0); 659fecfb714SLisandro Dalcin } 660fecfb714SLisandro Dalcin 6619371c9d4SSatish Balay static PetscErrorCode TSForwardStep_RK(TS ts) { 662922a638cSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 663922a638cSHong Zhang RKTableau tab = rk->tableau; 664922a638cSHong Zhang Mat J, *MatsFwdSensipTemp = rk->MatsFwdSensipTemp; 665922a638cSHong Zhang const PetscInt s = tab->s; 666922a638cSHong Zhang const PetscReal *A = tab->A, *c = tab->c, *b = tab->b; 667922a638cSHong Zhang Vec *Y = rk->Y; 668922a638cSHong Zhang PetscInt i, j; 669922a638cSHong Zhang PetscReal stage_time, h = ts->time_step; 670922a638cSHong Zhang PetscBool zero; 671922a638cSHong Zhang 672922a638cSHong Zhang PetscFunctionBegin; 6739566063dSJacob Faibussowitsch PetscCall(MatCopy(ts->mat_sensip, rk->MatFwdSensip0, SAME_NONZERO_PATTERN)); 6749566063dSJacob Faibussowitsch PetscCall(TSGetRHSJacobian(ts, &J, NULL, NULL, NULL)); 675922a638cSHong Zhang 676922a638cSHong Zhang for (i = 0; i < s; i++) { 677922a638cSHong Zhang stage_time = ts->ptime + h * c[i]; 678922a638cSHong Zhang zero = PETSC_FALSE; 679922a638cSHong Zhang if (b[i] == 0 && i == s - 1) zero = PETSC_TRUE; 680922a638cSHong Zhang /* TLM Stage values */ 681922a638cSHong Zhang if (!i) { 6829566063dSJacob Faibussowitsch PetscCall(MatCopy(ts->mat_sensip, rk->MatsFwdStageSensip[i], SAME_NONZERO_PATTERN)); 683922a638cSHong Zhang } else if (!zero) { 6849566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(rk->MatsFwdStageSensip[i])); 68548a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(MatAXPY(rk->MatsFwdStageSensip[i], h * A[i * s + j], MatsFwdSensipTemp[j], SAME_NONZERO_PATTERN)); 6869566063dSJacob Faibussowitsch PetscCall(MatAXPY(rk->MatsFwdStageSensip[i], 1., ts->mat_sensip, SAME_NONZERO_PATTERN)); 687922a638cSHong Zhang } else { 6889566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(rk->MatsFwdStageSensip[i])); 689922a638cSHong Zhang } 690922a638cSHong Zhang 6919566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobian(ts, stage_time, Y[i], J, J)); 6929566063dSJacob Faibussowitsch PetscCall(MatMatMult(J, rk->MatsFwdStageSensip[i], MAT_REUSE_MATRIX, PETSC_DEFAULT, &MatsFwdSensipTemp[i])); 69313af1a74SHong Zhang if (ts->Jacprhs) { 6949566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobianP(ts, stage_time, Y[i], ts->Jacprhs)); /* get f_p */ 69513af1a74SHong Zhang if (ts->vecs_sensi2p) { /* TLM used for 2nd-order adjoint */ 69613af1a74SHong Zhang PetscScalar *xarr; 6979566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(MatsFwdSensipTemp[i], 0, &xarr)); 6989566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(rk->VecDeltaFwdSensipCol, xarr)); 6999566063dSJacob Faibussowitsch PetscCall(MatMultAdd(ts->Jacprhs, ts->vec_dir, rk->VecDeltaFwdSensipCol, rk->VecDeltaFwdSensipCol)); 7009566063dSJacob Faibussowitsch PetscCall(VecResetArray(rk->VecDeltaFwdSensipCol)); 7019566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(MatsFwdSensipTemp[i], &xarr)); 70213af1a74SHong Zhang } else { 7039566063dSJacob Faibussowitsch PetscCall(MatAXPY(MatsFwdSensipTemp[i], 1., ts->Jacprhs, SUBSET_NONZERO_PATTERN)); 70413af1a74SHong Zhang } 705922a638cSHong Zhang } 706922a638cSHong Zhang } 707922a638cSHong Zhang 70848a46eb9SPierre Jolivet for (i = 0; i < s; i++) PetscCall(MatAXPY(ts->mat_sensip, h * b[i], rk->MatsFwdSensipTemp[i], SAME_NONZERO_PATTERN)); 709922a638cSHong Zhang rk->status = TS_STEP_COMPLETE; 710922a638cSHong Zhang PetscFunctionReturn(0); 711922a638cSHong Zhang } 712922a638cSHong Zhang 7139371c9d4SSatish Balay static PetscErrorCode TSForwardGetStages_RK(TS ts, PetscInt *ns, Mat **stagesensip) { 714922a638cSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 715922a638cSHong Zhang RKTableau tab = rk->tableau; 716922a638cSHong Zhang 717922a638cSHong Zhang PetscFunctionBegin; 718922a638cSHong Zhang if (ns) *ns = tab->s; 719922a638cSHong Zhang if (stagesensip) *stagesensip = rk->MatsFwdStageSensip; 720922a638cSHong Zhang PetscFunctionReturn(0); 721922a638cSHong Zhang } 722922a638cSHong Zhang 7239371c9d4SSatish Balay static PetscErrorCode TSForwardSetUp_RK(TS ts) { 724922a638cSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 725922a638cSHong Zhang RKTableau tab = rk->tableau; 726922a638cSHong Zhang PetscInt i; 727922a638cSHong Zhang 728922a638cSHong Zhang PetscFunctionBegin; 729922a638cSHong Zhang /* backup sensitivity results for roll-backs */ 7309566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &rk->MatFwdSensip0)); 731922a638cSHong Zhang 7329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &rk->MatsFwdStageSensip)); 7339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &rk->MatsFwdSensipTemp)); 734922a638cSHong Zhang for (i = 0; i < tab->s; i++) { 7359566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &rk->MatsFwdStageSensip[i])); 7369566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &rk->MatsFwdSensipTemp[i])); 737922a638cSHong Zhang } 7389566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &rk->VecDeltaFwdSensipCol)); 739922a638cSHong Zhang PetscFunctionReturn(0); 740922a638cSHong Zhang } 741922a638cSHong Zhang 7429371c9d4SSatish Balay static PetscErrorCode TSForwardReset_RK(TS ts) { 743922a638cSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 744922a638cSHong Zhang RKTableau tab = rk->tableau; 745922a638cSHong Zhang PetscInt i; 746922a638cSHong Zhang 747922a638cSHong Zhang PetscFunctionBegin; 7489566063dSJacob Faibussowitsch PetscCall(MatDestroy(&rk->MatFwdSensip0)); 749922a638cSHong Zhang if (rk->MatsFwdStageSensip) { 75048a46eb9SPierre Jolivet for (i = 0; i < tab->s; i++) PetscCall(MatDestroy(&rk->MatsFwdStageSensip[i])); 7519566063dSJacob Faibussowitsch PetscCall(PetscFree(rk->MatsFwdStageSensip)); 752922a638cSHong Zhang } 753922a638cSHong Zhang if (rk->MatsFwdSensipTemp) { 75448a46eb9SPierre Jolivet for (i = 0; i < tab->s; i++) PetscCall(MatDestroy(&rk->MatsFwdSensipTemp[i])); 7559566063dSJacob Faibussowitsch PetscCall(PetscFree(rk->MatsFwdSensipTemp)); 756922a638cSHong Zhang } 7579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rk->VecDeltaFwdSensipCol)); 758922a638cSHong Zhang PetscFunctionReturn(0); 759922a638cSHong Zhang } 760922a638cSHong Zhang 7619371c9d4SSatish Balay static PetscErrorCode TSStep_RK(TS ts) { 762f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK *)ts->data; 763f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 764f68a32c8SEmil Constantinescu const PetscInt s = tab->s; 765fecfb714SLisandro Dalcin const PetscReal *A = tab->A, *c = tab->c; 766f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 767f68a32c8SEmil Constantinescu Vec *Y = rk->Y, *YdotRHS = rk->YdotRHS; 768d1905564SLisandro Dalcin PetscBool FSAL = tab->FSAL; 769f68a32c8SEmil Constantinescu TSAdapt adapt; 770fecfb714SLisandro Dalcin PetscInt i, j; 771be5899b3SLisandro Dalcin PetscInt rejections = 0; 772be5899b3SLisandro Dalcin PetscBool stageok, accept = PETSC_TRUE; 773be5899b3SLisandro Dalcin PetscReal next_time_step = ts->time_step; 774f68a32c8SEmil Constantinescu 775f68a32c8SEmil Constantinescu PetscFunctionBegin; 776d1905564SLisandro Dalcin if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE; 7779566063dSJacob Faibussowitsch if (FSAL) PetscCall(VecCopy(YdotRHS[s - 1], YdotRHS[0])); 778d1905564SLisandro Dalcin 779f68a32c8SEmil Constantinescu rk->status = TS_STEP_INCOMPLETE; 780be5899b3SLisandro Dalcin while (!ts->reason && rk->status != TS_STEP_COMPLETE) { 781be5899b3SLisandro Dalcin PetscReal t = ts->ptime; 782f68a32c8SEmil Constantinescu PetscReal h = ts->time_step; 783f68a32c8SEmil Constantinescu for (i = 0; i < s; i++) { 7849be3e283SDebojyoti Ghosh rk->stage_time = t + h * c[i]; 7859566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, rk->stage_time)); 7869566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, Y[i])); 787f68a32c8SEmil Constantinescu for (j = 0; j < i; j++) w[j] = h * A[i * s + j]; 7889566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Y[i], i, w, YdotRHS)); 7899566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, rk->stage_time, i, Y)); 7909566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 7919566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(adapt, ts, rk->stage_time, Y[i], &stageok)); 792fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 7938f738a7cSLisandro Dalcin if (FSAL && !i) continue; 7949566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, t + h * c[i], Y[i], YdotRHS[i])); 795f68a32c8SEmil Constantinescu } 796be5899b3SLisandro Dalcin 797be5899b3SLisandro Dalcin rk->status = TS_STEP_INCOMPLETE; 7989566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep(ts, tab->order, ts->vec_sol, NULL)); 799f68a32c8SEmil Constantinescu rk->status = TS_STEP_PENDING; 8009566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 8019566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(adapt)); 8029566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidateAdd(adapt, tab->name, tab->order, 1, tab->ccfl, (PetscReal)tab->s, PETSC_TRUE)); 8039566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 804be5899b3SLisandro Dalcin rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 805be5899b3SLisandro Dalcin if (!accept) { /* Roll back the current step */ 8069566063dSJacob Faibussowitsch PetscCall(TSRollBack_RK(ts)); 807be5899b3SLisandro Dalcin ts->time_step = next_time_step; 808be5899b3SLisandro Dalcin goto reject_step; 809be5899b3SLisandro Dalcin } 810be5899b3SLisandro Dalcin 8110f7a1166SHong Zhang if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 8120f7a1166SHong Zhang rk->ptime = ts->ptime; 8130f7a1166SHong Zhang rk->time_step = ts->time_step; 814493ed95fSHong Zhang } 815be5899b3SLisandro Dalcin 816f68a32c8SEmil Constantinescu ts->ptime += ts->time_step; 817f68a32c8SEmil Constantinescu ts->time_step = next_time_step; 818f68a32c8SEmil Constantinescu break; 819be5899b3SLisandro Dalcin 820be5899b3SLisandro Dalcin reject_step: 8219371c9d4SSatish Balay ts->reject++; 8229371c9d4SSatish Balay accept = PETSC_FALSE; 823be5899b3SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 824be5899b3SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 82563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections)); 826f68a32c8SEmil Constantinescu } 827f68a32c8SEmil Constantinescu } 828f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 829e4dd521cSBarry Smith } 830e4dd521cSBarry Smith 8319371c9d4SSatish Balay static PetscErrorCode TSAdjointSetUp_RK(TS ts) { 832f6a906c0SBarry Smith TS_RK *rk = (TS_RK *)ts->data; 833be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 834be5899b3SLisandro Dalcin PetscInt s = tab->s; 835f6a906c0SBarry Smith 836f6a906c0SBarry Smith PetscFunctionBegin; 837f6a906c0SBarry Smith if (ts->adjointsetupcalled++) PetscFunctionReturn(0); 8389566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], s * ts->numcost, &rk->VecsDeltaLam)); 8399566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &rk->VecsSensiTemp)); 84048a46eb9SPierre Jolivet if (ts->vecs_sensip) PetscCall(VecDuplicate(ts->vecs_sensip[0], &rk->VecDeltaMu)); 84113af1a74SHong Zhang if (ts->vecs_sensi2) { 8429566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], s * ts->numcost, &rk->VecsDeltaLam2)); 8439566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vecs_sensi2[0], ts->numcost, &rk->VecsSensi2Temp)); 84413af1a74SHong Zhang } 84548a46eb9SPierre Jolivet if (ts->vecs_sensi2p) PetscCall(VecDuplicate(ts->vecs_sensi2p[0], &rk->VecDeltaMu2)); 846f6a906c0SBarry Smith PetscFunctionReturn(0); 847f6a906c0SBarry Smith } 848f6a906c0SBarry Smith 84935f5172aSHong Zhang /* 85035f5172aSHong Zhang Assumptions: 85135f5172aSHong Zhang - TSStep_RK() always evaluates the step with b, not bembed. 85235f5172aSHong Zhang */ 8539371c9d4SSatish Balay static PetscErrorCode TSAdjointStep_RK(TS ts) { 854c235aa8dSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 855cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 856c235aa8dSHong Zhang RKTableau tab = rk->tableau; 8573ca0882eSHong Zhang Mat J, Jpre, Jquad; 858c235aa8dSHong Zhang const PetscInt s = tab->s; 859c235aa8dSHong Zhang const PetscReal *A = tab->A, *b = tab->b, *c = tab->c; 86013af1a74SHong Zhang PetscScalar *w = rk->work, *xarr; 8612e7b7f96SHong Zhang Vec *Y = rk->Y, *VecsDeltaLam = rk->VecsDeltaLam, VecDeltaMu = rk->VecDeltaMu, *VecsSensiTemp = rk->VecsSensiTemp; 86213af1a74SHong Zhang Vec *VecsDeltaLam2 = rk->VecsDeltaLam2, VecDeltaMu2 = rk->VecDeltaMu2, *VecsSensi2Temp = rk->VecsSensi2Temp; 863cd4cee2dSHong Zhang Vec VecDRDUTransCol = ts->vec_drdu_col, VecDRDPTransCol = ts->vec_drdp_col; 864b18ea86cSHong Zhang PetscInt i, j, nadj; 8653d3eaba7SBarry Smith PetscReal t = ts->ptime; 8663ca0882eSHong Zhang PetscReal h = ts->time_step; 867c235aa8dSHong Zhang 868d2daff3dSHong Zhang PetscFunctionBegin; 869c235aa8dSHong Zhang rk->status = TS_STEP_INCOMPLETE; 8703389c7d9SStefano Zampini 8719566063dSJacob Faibussowitsch PetscCall(TSGetRHSJacobian(ts, &J, &Jpre, NULL, NULL)); 87248a46eb9SPierre Jolivet if (quadts) PetscCall(TSGetRHSJacobian(quadts, &Jquad, NULL, NULL, NULL)); 8732e7b7f96SHong Zhang for (i = s - 1; i >= 0; i--) { 8746a1e4597SHong Zhang if (tab->FSAL && i == s - 1) { 8756a1e4597SHong Zhang /* VecsDeltaLam[nadj*s+s-1] are initialized with zeros and the values never change.*/ 8766a1e4597SHong Zhang continue; 8776a1e4597SHong Zhang } 8783ca0882eSHong Zhang rk->stage_time = t + h * (1.0 - c[i]); 8799566063dSJacob Faibussowitsch PetscCall(TSComputeSNESJacobian(ts, Y[i], J, Jpre)); 8809371c9d4SSatish Balay if (quadts) { PetscCall(TSComputeRHSJacobian(quadts, rk->stage_time, Y[i], Jquad, Jquad)); /* get r_u^T */ } 8813389c7d9SStefano Zampini if (ts->vecs_sensip) { 8829566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobianP(ts, rk->stage_time, Y[i], ts->Jacprhs)); /* get f_p */ 8839371c9d4SSatish Balay if (quadts) { PetscCall(TSComputeRHSJacobianP(quadts, rk->stage_time, Y[i], quadts->Jacprhs)); /* get f_p for the quadrature */ } 88435f5172aSHong Zhang } 8853389c7d9SStefano Zampini 88635f5172aSHong Zhang if (b[i]) { 88735f5172aSHong Zhang for (j = i + 1; j < s; j++) w[j - i - 1] = A[j * s + i] / b[i]; /* coefficients for computing VecsSensiTemp */ 88835f5172aSHong Zhang } else { 88935f5172aSHong Zhang for (j = i + 1; j < s; j++) w[j - i - 1] = A[j * s + i]; /* coefficients for computing VecsSensiTemp */ 89035f5172aSHong Zhang } 8912e7b7f96SHong Zhang 8922e7b7f96SHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 8933389c7d9SStefano Zampini /* Stage values of lambda */ 89435f5172aSHong Zhang if (b[i]) { 89535f5172aSHong Zhang /* lambda_{n+1} + \sum_{j=i+1}^s a_{ji}/b[i]*lambda_{s,j} */ 8969566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vecs_sensi[nadj], VecsSensiTemp[nadj])); /* VecDeltaLam is an vec array of size s by numcost */ 8979566063dSJacob Faibussowitsch PetscCall(VecMAXPY(VecsSensiTemp[nadj], s - i - 1, w, &VecsDeltaLam[nadj * s + i + 1])); 8989566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(J, VecsSensiTemp[nadj], VecsDeltaLam[nadj * s + i])); /* VecsSensiTemp will be reused by 2nd-order adjoint */ 8999566063dSJacob Faibussowitsch PetscCall(VecScale(VecsDeltaLam[nadj * s + i], -h * b[i])); 90035f5172aSHong Zhang if (quadts) { 9019566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(Jquad, nadj, &xarr)); 9029566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(VecDRDUTransCol, xarr)); 9039566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecsDeltaLam[nadj * s + i], -h * b[i], VecDRDUTransCol)); 9049566063dSJacob Faibussowitsch PetscCall(VecResetArray(VecDRDUTransCol)); 9059566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(Jquad, &xarr)); 906cd4cee2dSHong Zhang } 9073389c7d9SStefano Zampini } else { 9086a1e4597SHong Zhang /* \sum_{j=i+1}^s a_{ji}*lambda_{s,j} */ 9099566063dSJacob Faibussowitsch PetscCall(VecSet(VecsSensiTemp[nadj], 0)); 9109566063dSJacob Faibussowitsch PetscCall(VecMAXPY(VecsSensiTemp[nadj], s - i - 1, w, &VecsDeltaLam[nadj * s + i + 1])); 9119566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(J, VecsSensiTemp[nadj], VecsDeltaLam[nadj * s + i])); 9129566063dSJacob Faibussowitsch PetscCall(VecScale(VecsDeltaLam[nadj * s + i], -h)); 9133389c7d9SStefano Zampini } 9146497ce14SHong Zhang 915ad8e2604SHong Zhang /* Stage values of mu */ 9166497ce14SHong Zhang if (ts->vecs_sensip) { 91735f5172aSHong Zhang if (b[i]) { 9189566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ts->Jacprhs, VecsSensiTemp[nadj], VecDeltaMu)); 9199566063dSJacob Faibussowitsch PetscCall(VecScale(VecDeltaMu, -h * b[i])); 92035f5172aSHong Zhang if (quadts) { 9219566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(quadts->Jacprhs, nadj, &xarr)); 9229566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(VecDRDPTransCol, xarr)); 9239566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecDeltaMu, -h * b[i], VecDRDPTransCol)); 9249566063dSJacob Faibussowitsch PetscCall(VecResetArray(VecDRDPTransCol)); 9259566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(quadts->Jacprhs, &xarr)); 926493ed95fSHong Zhang } 92735f5172aSHong Zhang } else { 9289566063dSJacob Faibussowitsch PetscCall(VecScale(VecDeltaMu, -h)); 92935f5172aSHong Zhang } 9309566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensip[nadj], 1., VecDeltaMu)); /* update sensip for each stage */ 931ad8e2604SHong Zhang } 932c235aa8dSHong Zhang } 93313af1a74SHong Zhang 93413af1a74SHong Zhang if (ts->vecs_sensi2 && ts->forward_solve) { /* 2nd-order adjoint, TLM mode has to be turned on */ 93513af1a74SHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 9369566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(rk->MatsFwdStageSensip[i], 0, &xarr)); 9379566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr)); 93813af1a74SHong Zhang /* lambda_s^T F_UU w_1 */ 9399566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionUU(ts, rk->stage_time, Y[i], VecsSensiTemp, ts->vec_sensip_col, ts->vecs_guu)); 94035f5172aSHong Zhang if (quadts) { 94135f5172aSHong Zhang /* R_UU w_1 */ 9429566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionUU(quadts, rk->stage_time, Y[i], NULL, ts->vec_sensip_col, ts->vecs_guu)); 94335f5172aSHong Zhang } 94435f5172aSHong Zhang if (ts->vecs_sensip) { 94513af1a74SHong Zhang /* lambda_s^T F_UP w_2 */ 9469566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionUP(ts, rk->stage_time, Y[i], VecsSensiTemp, ts->vec_dir, ts->vecs_gup)); 94735f5172aSHong Zhang if (quadts) { 94835f5172aSHong Zhang /* R_UP w_2 */ 9499566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionUP(quadts, rk->stage_time, Y[i], NULL, ts->vec_sensip_col, ts->vecs_gup)); 95035f5172aSHong Zhang } 95135f5172aSHong Zhang } 95213af1a74SHong Zhang if (ts->vecs_sensi2p) { 95313af1a74SHong Zhang /* lambda_s^T F_PU w_1 */ 9549566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionPU(ts, rk->stage_time, Y[i], VecsSensiTemp, ts->vec_sensip_col, ts->vecs_gpu)); 95535f5172aSHong Zhang /* lambda_s^T F_PP w_2 */ 9569566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionPP(ts, rk->stage_time, Y[i], VecsSensiTemp, ts->vec_dir, ts->vecs_gpp)); 95735f5172aSHong Zhang if (b[i] && quadts) { 95835f5172aSHong Zhang /* R_PU w_1 */ 9599566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionPU(quadts, rk->stage_time, Y[i], NULL, ts->vec_sensip_col, ts->vecs_gpu)); 96035f5172aSHong Zhang /* R_PP w_2 */ 9619566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionPP(quadts, rk->stage_time, Y[i], NULL, ts->vec_dir, ts->vecs_gpp)); 96235f5172aSHong Zhang } 96313af1a74SHong Zhang } 9649566063dSJacob Faibussowitsch PetscCall(VecResetArray(ts->vec_sensip_col)); 9659566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(rk->MatsFwdStageSensip[i], &xarr)); 96613af1a74SHong Zhang 96713af1a74SHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 96813af1a74SHong Zhang /* Stage values of lambda */ 96935f5172aSHong Zhang if (b[i]) { 97035f5172aSHong Zhang /* J_i^T*(Lambda_{n+1}+\sum_{j=i+1}^s a_{ji}/b_i*Lambda_{s,j} */ 9719566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vecs_sensi2[nadj], VecsSensi2Temp[nadj])); 9729566063dSJacob Faibussowitsch PetscCall(VecMAXPY(VecsSensi2Temp[nadj], s - i - 1, w, &VecsDeltaLam2[nadj * s + i + 1])); 9739566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(J, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj * s + i])); 9749566063dSJacob Faibussowitsch PetscCall(VecScale(VecsDeltaLam2[nadj * s + i], -h * b[i])); 9759566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecsDeltaLam2[nadj * s + i], -h * b[i], ts->vecs_guu[nadj])); 97648a46eb9SPierre Jolivet if (ts->vecs_sensip) PetscCall(VecAXPY(VecsDeltaLam2[nadj * s + i], -h * b[i], ts->vecs_gup[nadj])); 97713af1a74SHong Zhang } else { 97835f5172aSHong Zhang /* \sum_{j=i+1}^s a_{ji}*Lambda_{s,j} */ 9799566063dSJacob Faibussowitsch PetscCall(VecSet(VecsDeltaLam2[nadj * s + i], 0)); 9809566063dSJacob Faibussowitsch PetscCall(VecMAXPY(VecsSensi2Temp[nadj], s - i - 1, w, &VecsDeltaLam2[nadj * s + i + 1])); 9819566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(J, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj * s + i])); 9829566063dSJacob Faibussowitsch PetscCall(VecScale(VecsDeltaLam2[nadj * s + i], -h)); 9839566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecsDeltaLam2[nadj * s + i], -h, ts->vecs_guu[nadj])); 98448a46eb9SPierre Jolivet if (ts->vecs_sensip) PetscCall(VecAXPY(VecsDeltaLam2[nadj * s + i], -h, ts->vecs_gup[nadj])); 98535f5172aSHong Zhang } 98635f5172aSHong Zhang if (ts->vecs_sensi2p) { /* 2nd-order adjoint for parameters */ 9879566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ts->Jacprhs, VecsSensi2Temp[nadj], VecDeltaMu2)); 98835f5172aSHong Zhang if (b[i]) { 9899566063dSJacob Faibussowitsch PetscCall(VecScale(VecDeltaMu2, -h * b[i])); 9909566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecDeltaMu2, -h * b[i], ts->vecs_gpu[nadj])); 9919566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecDeltaMu2, -h * b[i], ts->vecs_gpp[nadj])); 99213af1a74SHong Zhang } else { 9939566063dSJacob Faibussowitsch PetscCall(VecScale(VecDeltaMu2, -h)); 9949566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecDeltaMu2, -h, ts->vecs_gpu[nadj])); 9959566063dSJacob Faibussowitsch PetscCall(VecAXPY(VecDeltaMu2, -h, ts->vecs_gpp[nadj])); 99613af1a74SHong Zhang } 9979566063dSJacob Faibussowitsch PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], 1, VecDeltaMu2)); /* update sensi2p for each stage */ 99813af1a74SHong Zhang } 99913af1a74SHong Zhang } 100013af1a74SHong Zhang } 10016497ce14SHong Zhang } 1002c235aa8dSHong Zhang 1003c235aa8dSHong Zhang for (j = 0; j < s; j++) w[j] = 1.0; 10042e7b7f96SHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { /* no need to do this for mu's */ 10059566063dSJacob Faibussowitsch PetscCall(VecMAXPY(ts->vecs_sensi[nadj], s, w, &VecsDeltaLam[nadj * s])); 100648a46eb9SPierre Jolivet if (ts->vecs_sensi2) PetscCall(VecMAXPY(ts->vecs_sensi2[nadj], s, w, &VecsDeltaLam2[nadj * s])); 10076497ce14SHong Zhang } 1008c235aa8dSHong Zhang rk->status = TS_STEP_COMPLETE; 1009d2daff3dSHong Zhang PetscFunctionReturn(0); 1010d2daff3dSHong Zhang } 1011d2daff3dSHong Zhang 10129371c9d4SSatish Balay static PetscErrorCode TSAdjointReset_RK(TS ts) { 101313af1a74SHong Zhang TS_RK *rk = (TS_RK *)ts->data; 101413af1a74SHong Zhang RKTableau tab = rk->tableau; 101513af1a74SHong Zhang 101613af1a74SHong Zhang PetscFunctionBegin; 10179566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s * ts->numcost, &rk->VecsDeltaLam)); 10189566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(ts->numcost, &rk->VecsSensiTemp)); 10199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rk->VecDeltaMu)); 10209566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s * ts->numcost, &rk->VecsDeltaLam2)); 10219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rk->VecDeltaMu2)); 10229566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(ts->numcost, &rk->VecsSensi2Temp)); 102313af1a74SHong Zhang PetscFunctionReturn(0); 102413af1a74SHong Zhang } 102513af1a74SHong Zhang 10269371c9d4SSatish Balay static PetscErrorCode TSInterpolate_RK(TS ts, PetscReal itime, Vec X) { 102755de54fcSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 102855de54fcSHong Zhang PetscInt s = rk->tableau->s, p = rk->tableau->p, i, j; 102955de54fcSHong Zhang PetscReal h; 103055de54fcSHong Zhang PetscReal tt, t; 103155de54fcSHong Zhang PetscScalar *b; 103255de54fcSHong Zhang const PetscReal *B = rk->tableau->binterp; 103355de54fcSHong Zhang 103455de54fcSHong Zhang PetscFunctionBegin; 10353c633725SBarry Smith PetscCheck(B, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSRK %s does not have an interpolation formula", rk->tableau->name); 103655de54fcSHong Zhang 103755de54fcSHong Zhang switch (rk->status) { 103855de54fcSHong Zhang case TS_STEP_INCOMPLETE: 103955de54fcSHong Zhang case TS_STEP_PENDING: 104055de54fcSHong Zhang h = ts->time_step; 104155de54fcSHong Zhang t = (itime - ts->ptime) / h; 104255de54fcSHong Zhang break; 104355de54fcSHong Zhang case TS_STEP_COMPLETE: 104455de54fcSHong Zhang h = ts->ptime - ts->ptime_prev; 104555de54fcSHong Zhang t = (itime - ts->ptime) / h + 1; /* In the interval [0,1] */ 104655de54fcSHong Zhang break; 104755de54fcSHong Zhang default: SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus"); 104855de54fcSHong Zhang } 10499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s, &b)); 105055de54fcSHong Zhang for (i = 0; i < s; i++) b[i] = 0; 105155de54fcSHong Zhang for (j = 0, tt = t; j < p; j++, tt *= t) { 1052ad540459SPierre Jolivet for (i = 0; i < s; i++) b[i] += h * B[i * p + j] * tt; 105355de54fcSHong Zhang } 10549566063dSJacob Faibussowitsch PetscCall(VecCopy(rk->Y[0], X)); 10559566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, b, rk->YdotRHS)); 10569566063dSJacob Faibussowitsch PetscCall(PetscFree(b)); 105755de54fcSHong Zhang PetscFunctionReturn(0); 105855de54fcSHong Zhang } 105955de54fcSHong Zhang 106055de54fcSHong Zhang /*------------------------------------------------------------*/ 106155de54fcSHong Zhang 10629371c9d4SSatish Balay static PetscErrorCode TSRKTableauReset(TS ts) { 1063be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK *)ts->data; 1064be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 1065be5899b3SLisandro Dalcin 1066be5899b3SLisandro Dalcin PetscFunctionBegin; 1067be5899b3SLisandro Dalcin if (!tab) PetscFunctionReturn(0); 10689566063dSJacob Faibussowitsch PetscCall(PetscFree(rk->work)); 10699566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &rk->Y)); 10709566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &rk->YdotRHS)); 1071be5899b3SLisandro Dalcin PetscFunctionReturn(0); 1072be5899b3SLisandro Dalcin } 1073be5899b3SLisandro Dalcin 10749371c9d4SSatish Balay static PetscErrorCode TSReset_RK(TS ts) { 1075e4dd521cSBarry Smith PetscFunctionBegin; 10769566063dSJacob Faibussowitsch PetscCall(TSRKTableauReset(ts)); 10770fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 1078cac4c232SBarry Smith PetscTryMethod(ts, "TSReset_RK_MultirateSplit_C", (TS), (ts)); 10790fe4d17eSHong Zhang } else { 1080cac4c232SBarry Smith PetscTryMethod(ts, "TSReset_RK_MultirateNonsplit_C", (TS), (ts)); 10810fe4d17eSHong Zhang } 1082277b19d0SLisandro Dalcin PetscFunctionReturn(0); 1083e4dd521cSBarry Smith } 1084277b19d0SLisandro Dalcin 10859371c9d4SSatish Balay static PetscErrorCode DMCoarsenHook_TSRK(DM fine, DM coarse, void *ctx) { 1086f68a32c8SEmil Constantinescu PetscFunctionBegin; 1087f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1088f68a32c8SEmil Constantinescu } 1089f68a32c8SEmil Constantinescu 10909371c9d4SSatish Balay static PetscErrorCode DMRestrictHook_TSRK(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) { 1091f68a32c8SEmil Constantinescu PetscFunctionBegin; 1092f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1093f68a32c8SEmil Constantinescu } 1094f68a32c8SEmil Constantinescu 10959371c9d4SSatish Balay static PetscErrorCode DMSubDomainHook_TSRK(DM dm, DM subdm, void *ctx) { 1096f68a32c8SEmil Constantinescu PetscFunctionBegin; 1097f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1098f68a32c8SEmil Constantinescu } 1099f68a32c8SEmil Constantinescu 11009371c9d4SSatish Balay static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) { 1101f68a32c8SEmil Constantinescu PetscFunctionBegin; 1102f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1103f68a32c8SEmil Constantinescu } 1104be5899b3SLisandro Dalcin 11059371c9d4SSatish Balay static PetscErrorCode TSRKTableauSetUp(TS ts) { 1106be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK *)ts->data; 1107be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 1108be5899b3SLisandro Dalcin 1109be5899b3SLisandro Dalcin PetscFunctionBegin; 11109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &rk->work)); 11119566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &rk->Y)); 11129566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &rk->YdotRHS)); 1113be5899b3SLisandro Dalcin PetscFunctionReturn(0); 1114be5899b3SLisandro Dalcin } 1115be5899b3SLisandro Dalcin 11169371c9d4SSatish Balay static PetscErrorCode TSSetUp_RK(TS ts) { 1117cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 1118f68a32c8SEmil Constantinescu DM dm; 1119f68a32c8SEmil Constantinescu 1120f68a32c8SEmil Constantinescu PetscFunctionBegin; 11219566063dSJacob Faibussowitsch PetscCall(TSCheckImplicitTerm(ts)); 11229566063dSJacob Faibussowitsch PetscCall(TSRKTableauSetUp(ts)); 1123cd4cee2dSHong Zhang if (quadts && ts->costintegralfwd) { 1124cd4cee2dSHong Zhang Mat Jquad; 11259566063dSJacob Faibussowitsch PetscCall(TSGetRHSJacobian(quadts, &Jquad, NULL, NULL, NULL)); 11260f7a1166SHong Zhang } 11279566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 11289566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSRK, DMRestrictHook_TSRK, ts)); 11299566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSRK, DMSubDomainRestrictHook_TSRK, ts)); 11300fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 1131cac4c232SBarry Smith PetscTryMethod(ts, "TSSetUp_RK_MultirateSplit_C", (TS), (ts)); 11320fe4d17eSHong Zhang } else { 1133cac4c232SBarry Smith PetscTryMethod(ts, "TSSetUp_RK_MultirateNonsplit_C", (TS), (ts)); 11340fe4d17eSHong Zhang } 1135e4dd521cSBarry Smith PetscFunctionReturn(0); 1136e4dd521cSBarry Smith } 1137d2daff3dSHong Zhang 11389371c9d4SSatish Balay static PetscErrorCode TSSetFromOptions_RK(TS ts, PetscOptionItems *PetscOptionsObject) { 1139be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK *)ts->data; 1140262ff7bbSBarry Smith 1141e4dd521cSBarry Smith PetscFunctionBegin; 1142d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "RK ODE solver options"); 1143f68a32c8SEmil Constantinescu { 1144f68a32c8SEmil Constantinescu RKTableauLink link; 1145f68a32c8SEmil Constantinescu PetscInt count, choice; 11460fe4d17eSHong Zhang PetscBool flg, use_multirate = PETSC_FALSE; 1147f68a32c8SEmil Constantinescu const char **namelist; 11482c9cb062Sluzhanghpp 11499371c9d4SSatish Balay for (link = RKTableauList, count = 0; link; link = link->next, count++) 11509371c9d4SSatish Balay ; 11519566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count, (char ***)&namelist)); 1152f68a32c8SEmil Constantinescu for (link = RKTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name; 11539566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_rk_multirate", "Use interpolation-based multirate RK method", "TSRKSetMultirate", rk->use_multirate, &use_multirate, &flg)); 11541baa6e33SBarry Smith if (flg) PetscCall(TSRKSetMultirate(ts, use_multirate)); 11559566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_rk_type", "Family of RK method", "TSRKSetType", (const char *const *)namelist, count, rk->tableau->name, &choice, &flg)); 11569566063dSJacob Faibussowitsch if (flg) PetscCall(TSRKSetType(ts, namelist[choice])); 11579566063dSJacob Faibussowitsch PetscCall(PetscFree(namelist)); 1158f68a32c8SEmil Constantinescu } 1159d0609cedSBarry Smith PetscOptionsHeadEnd(); 1160d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)ts), NULL, "Multirate methods options", ""); 11619566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_rk_dtratio", "time step ratio between slow and fast", "", rk->dtratio, &rk->dtratio, NULL)); 1162d0609cedSBarry Smith PetscOptionsEnd(); 1163e4dd521cSBarry Smith PetscFunctionReturn(0); 1164e4dd521cSBarry Smith } 1165e4dd521cSBarry Smith 11669371c9d4SSatish Balay static PetscErrorCode TSView_RK(TS ts, PetscViewer viewer) { 11675f70b478SJed Brown TS_RK *rk = (TS_RK *)ts->data; 11688370ee3bSLisandro Dalcin PetscBool iascii; 11692cdf8120Sasbjorn 1170e4dd521cSBarry Smith PetscFunctionBegin; 11719566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 11728370ee3bSLisandro Dalcin if (iascii) { 11739c334d8fSLisandro Dalcin RKTableau tab = rk->tableau; 1174f68a32c8SEmil Constantinescu TSRKType rktype; 11750f7a28cdSStefano Zampini const PetscReal *c; 11760f7a28cdSStefano Zampini PetscInt s; 1177f68a32c8SEmil Constantinescu char buf[512]; 11780f7a28cdSStefano Zampini PetscBool FSAL; 11790f7a28cdSStefano Zampini 11809566063dSJacob Faibussowitsch PetscCall(TSRKGetType(ts, &rktype)); 11819566063dSJacob Faibussowitsch PetscCall(TSRKGetTableau(ts, &s, NULL, NULL, &c, NULL, NULL, NULL, &FSAL)); 11829566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " RK type %s\n", rktype)); 118363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Order: %" PetscInt_FMT "\n", tab->order)); 11849566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " FSAL property: %s\n", FSAL ? "yes" : "no")); 11859566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", s, c)); 11869566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Abscissa c = %s\n", buf)); 11878370ee3bSLisandro Dalcin } 1188f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1189f68a32c8SEmil Constantinescu } 1190f68a32c8SEmil Constantinescu 11919371c9d4SSatish Balay static PetscErrorCode TSLoad_RK(TS ts, PetscViewer viewer) { 11929c334d8fSLisandro Dalcin TSAdapt adapt; 1193f68a32c8SEmil Constantinescu 1194f68a32c8SEmil Constantinescu PetscFunctionBegin; 11959566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 11969566063dSJacob Faibussowitsch PetscCall(TSAdaptLoad(adapt, viewer)); 1197f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1198f68a32c8SEmil Constantinescu } 1199f68a32c8SEmil Constantinescu 12002ea87ba9SLisandro Dalcin /*@ 12012ea87ba9SLisandro Dalcin TSRKGetOrder - Get the order of RK scheme 12022ea87ba9SLisandro Dalcin 12032ea87ba9SLisandro Dalcin Not collective 12042ea87ba9SLisandro Dalcin 12052ea87ba9SLisandro Dalcin Input Parameter: 12062ea87ba9SLisandro Dalcin . ts - timestepping context 12072ea87ba9SLisandro Dalcin 12082ea87ba9SLisandro Dalcin Output Parameter: 12092ea87ba9SLisandro Dalcin . order - order of RK-scheme 12102ea87ba9SLisandro Dalcin 12112ea87ba9SLisandro Dalcin Level: intermediate 12122ea87ba9SLisandro Dalcin 1213db781477SPatrick Sanan .seealso: `TSRKGetType()` 12142ea87ba9SLisandro Dalcin @*/ 12159371c9d4SSatish Balay PetscErrorCode TSRKGetOrder(TS ts, PetscInt *order) { 12162ea87ba9SLisandro Dalcin PetscFunctionBegin; 12172ea87ba9SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 12182ea87ba9SLisandro Dalcin PetscValidIntPointer(order, 2); 1219cac4c232SBarry Smith PetscUseMethod(ts, "TSRKGetOrder_C", (TS, PetscInt *), (ts, order)); 12202ea87ba9SLisandro Dalcin PetscFunctionReturn(0); 12212ea87ba9SLisandro Dalcin } 12222ea87ba9SLisandro Dalcin 1223f68a32c8SEmil Constantinescu /*@C 1224f68a32c8SEmil Constantinescu TSRKSetType - Set the type of RK scheme 1225f68a32c8SEmil Constantinescu 1226f68a32c8SEmil Constantinescu Logically collective 1227f68a32c8SEmil Constantinescu 1228d8d19677SJose E. Roman Input Parameters: 1229f68a32c8SEmil Constantinescu + ts - timestepping context 1230f68a32c8SEmil Constantinescu - rktype - type of RK-scheme 1231f68a32c8SEmil Constantinescu 1232287c2655SBarry Smith Options Database: 12339bd3e852SBarry Smith . -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs> 1234287c2655SBarry Smith 1235f68a32c8SEmil Constantinescu Level: intermediate 1236f68a32c8SEmil Constantinescu 1237db781477SPatrick Sanan .seealso: `TSRKGetType()`, `TSRK`, `TSRKType`, `TSRK1FE`, `TSRK2A`, `TSRK2B`, `TSRK3`, `TSRK3BS`, `TSRK4`, `TSRK5F`, `TSRK5DP`, `TSRK5BS`, `TSRK6VR`, `TSRK7VR`, `TSRK8VR` 1238f68a32c8SEmil Constantinescu @*/ 12399371c9d4SSatish Balay PetscErrorCode TSRKSetType(TS ts, TSRKType rktype) { 1240f68a32c8SEmil Constantinescu PetscFunctionBegin; 1241f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1242b92453a8SLisandro Dalcin PetscValidCharPointer(rktype, 2); 1243cac4c232SBarry Smith PetscTryMethod(ts, "TSRKSetType_C", (TS, TSRKType), (ts, rktype)); 1244f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1245f68a32c8SEmil Constantinescu } 1246f68a32c8SEmil Constantinescu 1247f68a32c8SEmil Constantinescu /*@C 1248f68a32c8SEmil Constantinescu TSRKGetType - Get the type of RK scheme 1249f68a32c8SEmil Constantinescu 12502ea87ba9SLisandro Dalcin Not collective 1251f68a32c8SEmil Constantinescu 1252f68a32c8SEmil Constantinescu Input Parameter: 1253f68a32c8SEmil Constantinescu . ts - timestepping context 1254f68a32c8SEmil Constantinescu 1255f68a32c8SEmil Constantinescu Output Parameter: 1256f68a32c8SEmil Constantinescu . rktype - type of RK-scheme 1257f68a32c8SEmil Constantinescu 1258f68a32c8SEmil Constantinescu Level: intermediate 1259f68a32c8SEmil Constantinescu 1260db781477SPatrick Sanan .seealso: `TSRKSetType()` 1261f68a32c8SEmil Constantinescu @*/ 12629371c9d4SSatish Balay PetscErrorCode TSRKGetType(TS ts, TSRKType *rktype) { 1263f68a32c8SEmil Constantinescu PetscFunctionBegin; 1264f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1265cac4c232SBarry Smith PetscUseMethod(ts, "TSRKGetType_C", (TS, TSRKType *), (ts, rktype)); 1266f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1267f68a32c8SEmil Constantinescu } 1268f68a32c8SEmil Constantinescu 12699371c9d4SSatish Balay static PetscErrorCode TSRKGetOrder_RK(TS ts, PetscInt *order) { 12702ea87ba9SLisandro Dalcin TS_RK *rk = (TS_RK *)ts->data; 12712ea87ba9SLisandro Dalcin 12722ea87ba9SLisandro Dalcin PetscFunctionBegin; 12732ea87ba9SLisandro Dalcin *order = rk->tableau->order; 12742ea87ba9SLisandro Dalcin PetscFunctionReturn(0); 12752ea87ba9SLisandro Dalcin } 12762ea87ba9SLisandro Dalcin 12779371c9d4SSatish Balay static PetscErrorCode TSRKGetType_RK(TS ts, TSRKType *rktype) { 1278f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK *)ts->data; 1279f68a32c8SEmil Constantinescu 1280f68a32c8SEmil Constantinescu PetscFunctionBegin; 1281f68a32c8SEmil Constantinescu *rktype = rk->tableau->name; 1282f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1283f68a32c8SEmil Constantinescu } 12842c9cb062Sluzhanghpp 12859371c9d4SSatish Balay static PetscErrorCode TSRKSetType_RK(TS ts, TSRKType rktype) { 1286f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK *)ts->data; 1287f68a32c8SEmil Constantinescu PetscBool match; 1288f68a32c8SEmil Constantinescu RKTableauLink link; 1289f68a32c8SEmil Constantinescu 1290f68a32c8SEmil Constantinescu PetscFunctionBegin; 1291f68a32c8SEmil Constantinescu if (rk->tableau) { 12929566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(rk->tableau->name, rktype, &match)); 1293f68a32c8SEmil Constantinescu if (match) PetscFunctionReturn(0); 1294f68a32c8SEmil Constantinescu } 1295f68a32c8SEmil Constantinescu for (link = RKTableauList; link; link = link->next) { 12969566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(link->tab.name, rktype, &match)); 1297f68a32c8SEmil Constantinescu if (match) { 12989566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSRKTableauReset(ts)); 1299f68a32c8SEmil Constantinescu rk->tableau = &link->tab; 13009566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSRKTableauSetUp(ts)); 13012ffb9264SLisandro Dalcin ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 1302f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1303f68a32c8SEmil Constantinescu } 1304f68a32c8SEmil Constantinescu } 130598921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", rktype); 1306e4dd521cSBarry Smith } 1307e4dd521cSBarry Smith 13089371c9d4SSatish Balay static PetscErrorCode TSGetStages_RK(TS ts, PetscInt *ns, Vec **Y) { 1309ff22ae23SHong Zhang TS_RK *rk = (TS_RK *)ts->data; 1310ff22ae23SHong Zhang 1311ff22ae23SHong Zhang PetscFunctionBegin; 13120429704eSStefano Zampini if (ns) *ns = rk->tableau->s; 1313d2daff3dSHong Zhang if (Y) *Y = rk->Y; 1314ff22ae23SHong Zhang PetscFunctionReturn(0); 1315ff22ae23SHong Zhang } 1316ff22ae23SHong Zhang 13179371c9d4SSatish Balay static PetscErrorCode TSDestroy_RK(TS ts) { 1318b3a6b972SJed Brown PetscFunctionBegin; 13199566063dSJacob Faibussowitsch PetscCall(TSReset_RK(ts)); 1320b3a6b972SJed Brown if (ts->dm) { 13219566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSRK, DMRestrictHook_TSRK, ts)); 13229566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSRK, DMSubDomainRestrictHook_TSRK, ts)); 1323b3a6b972SJed Brown } 13249566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 13259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetOrder_C", NULL)); 13269566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetType_C", NULL)); 13279566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKSetType_C", NULL)); 13289566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetTableau_C", NULL)); 13299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKSetMultirate_C", NULL)); 13309566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetMultirate_C", NULL)); 13312e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSetUp_RK_MultirateSplit_C", NULL)); 13322e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSReset_RK_MultirateSplit_C", NULL)); 13332e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSSetUp_RK_MultirateNonsplit_C", NULL)); 13342e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSReset_RK_MultirateNonsplit_C", NULL)); 1335b3a6b972SJed Brown PetscFunctionReturn(0); 1336b3a6b972SJed Brown } 1337ff22ae23SHong Zhang 13383ca0882eSHong Zhang /* 13393ca0882eSHong Zhang This defines the nonlinear equation that is to be solved with SNES 13403ca0882eSHong Zhang We do not need to solve the equation; we just use SNES to approximate the Jacobian 13413ca0882eSHong Zhang */ 13429371c9d4SSatish Balay static PetscErrorCode SNESTSFormFunction_RK(SNES snes, Vec x, Vec y, TS ts) { 13433ca0882eSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 13443ca0882eSHong Zhang DM dm, dmsave; 13453ca0882eSHong Zhang 13463ca0882eSHong Zhang PetscFunctionBegin; 13479566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 13483ca0882eSHong Zhang /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 13493ca0882eSHong Zhang dmsave = ts->dm; 13503ca0882eSHong Zhang ts->dm = dm; 13519566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, rk->stage_time, x, y)); 13523ca0882eSHong Zhang ts->dm = dmsave; 13533ca0882eSHong Zhang PetscFunctionReturn(0); 13543ca0882eSHong Zhang } 13553ca0882eSHong Zhang 13569371c9d4SSatish Balay static PetscErrorCode SNESTSFormJacobian_RK(SNES snes, Vec x, Mat A, Mat B, TS ts) { 13573ca0882eSHong Zhang TS_RK *rk = (TS_RK *)ts->data; 13583ca0882eSHong Zhang DM dm, dmsave; 13593ca0882eSHong Zhang 13603ca0882eSHong Zhang PetscFunctionBegin; 13619566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 13623ca0882eSHong Zhang dmsave = ts->dm; 13633ca0882eSHong Zhang ts->dm = dm; 13649566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobian(ts, rk->stage_time, x, A, B)); 13653ca0882eSHong Zhang ts->dm = dmsave; 13663ca0882eSHong Zhang PetscFunctionReturn(0); 13673ca0882eSHong Zhang } 13683ca0882eSHong Zhang 136921052c0cSHong Zhang /*@C 137021052c0cSHong Zhang TSRKSetMultirate - Use the interpolation-based multirate RK method 137121052c0cSHong Zhang 137221052c0cSHong Zhang Logically collective 137321052c0cSHong Zhang 1374d8d19677SJose E. Roman Input Parameters: 137521052c0cSHong Zhang + ts - timestepping context 137621052c0cSHong Zhang - use_multirate - PETSC_TRUE enables the multirate RK method, sets the basic method to be RK2A and sets the ratio between slow stepsize and fast stepsize to be 2 137721052c0cSHong Zhang 137821052c0cSHong Zhang Options Database: 137921052c0cSHong Zhang . -ts_rk_multirate - <true,false> 138021052c0cSHong Zhang 138121052c0cSHong Zhang Notes: 138221052c0cSHong Zhang The multirate method requires interpolation. The default interpolation works for 1st- and 2nd- order RK, but not for high-order RKs except TSRK5DP which comes with the interpolation coeffcients (binterp). 138321052c0cSHong Zhang 138421052c0cSHong Zhang Level: intermediate 138521052c0cSHong Zhang 1386db781477SPatrick Sanan .seealso: `TSRKGetMultirate()` 138721052c0cSHong Zhang @*/ 13889371c9d4SSatish Balay PetscErrorCode TSRKSetMultirate(TS ts, PetscBool use_multirate) { 13898a4be4dbSHong Zhang PetscFunctionBegin; 1390cac4c232SBarry Smith PetscTryMethod(ts, "TSRKSetMultirate_C", (TS, PetscBool), (ts, use_multirate)); 139121052c0cSHong Zhang PetscFunctionReturn(0); 139221052c0cSHong Zhang } 139321052c0cSHong Zhang 139421052c0cSHong Zhang /*@C 139521052c0cSHong Zhang TSRKGetMultirate - Gets whether to Use the interpolation-based multirate RK method 139621052c0cSHong Zhang 139721052c0cSHong Zhang Not collective 139821052c0cSHong Zhang 139921052c0cSHong Zhang Input Parameter: 140021052c0cSHong Zhang . ts - timestepping context 140121052c0cSHong Zhang 140221052c0cSHong Zhang Output Parameter: 140321052c0cSHong Zhang . use_multirate - PETSC_TRUE if the multirate RK method is enabled, PETSC_FALSE otherwise 140421052c0cSHong Zhang 140521052c0cSHong Zhang Level: intermediate 140621052c0cSHong Zhang 1407db781477SPatrick Sanan .seealso: `TSRKSetMultirate()` 140821052c0cSHong Zhang @*/ 14099371c9d4SSatish Balay PetscErrorCode TSRKGetMultirate(TS ts, PetscBool *use_multirate) { 14107dbacdf2SHong Zhang PetscFunctionBegin; 1411cac4c232SBarry Smith PetscUseMethod(ts, "TSRKGetMultirate_C", (TS, PetscBool *), (ts, use_multirate)); 141221052c0cSHong Zhang PetscFunctionReturn(0); 141321052c0cSHong Zhang } 141421052c0cSHong Zhang 1415ebe3b25bSBarry Smith /*MC 1416f68a32c8SEmil Constantinescu TSRK - ODE and DAE solver using Runge-Kutta schemes 1417ebe3b25bSBarry Smith 1418f68a32c8SEmil Constantinescu The user should provide the right hand side of the equation 1419f68a32c8SEmil Constantinescu using TSSetRHSFunction(). 1420ebe3b25bSBarry Smith 1421f68a32c8SEmil Constantinescu Notes: 142298164e13SLisandro Dalcin The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type 1423ebe3b25bSBarry Smith 1424d41222bbSBarry Smith Level: beginner 1425d41222bbSBarry Smith 1426db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSRKSetType()`, `TSRKGetType()`, `TSRK2D`, `TTSRK2E`, `TSRK3`, 1427db781477SPatrick Sanan `TSRK4`, `TSRK5`, `TSRKPRSSP2`, `TSRKBPR3`, `TSRKType`, `TSRKRegister()`, `TSRKSetMultirate()`, `TSRKGetMultirate()` 1428ebe3b25bSBarry Smith 1429ebe3b25bSBarry Smith M*/ 14309371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts) { 14311566a47fSLisandro Dalcin TS_RK *rk; 1432e4dd521cSBarry Smith 1433e4dd521cSBarry Smith PetscFunctionBegin; 14349566063dSJacob Faibussowitsch PetscCall(TSRKInitializePackage()); 1435f68a32c8SEmil Constantinescu 1436f68a32c8SEmil Constantinescu ts->ops->reset = TSReset_RK; 14375f70b478SJed Brown ts->ops->destroy = TSDestroy_RK; 14385f70b478SJed Brown ts->ops->view = TSView_RK; 1439f68a32c8SEmil Constantinescu ts->ops->load = TSLoad_RK; 1440f68a32c8SEmil Constantinescu ts->ops->setup = TSSetUp_RK; 1441f68a32c8SEmil Constantinescu ts->ops->interpolate = TSInterpolate_RK; 14422c9cb062Sluzhanghpp ts->ops->step = TSStep_RK; 1443f68a32c8SEmil Constantinescu ts->ops->evaluatestep = TSEvaluateStep_RK; 1444fecfb714SLisandro Dalcin ts->ops->rollback = TSRollBack_RK; 1445f68a32c8SEmil Constantinescu ts->ops->setfromoptions = TSSetFromOptions_RK; 1446ff22ae23SHong Zhang ts->ops->getstages = TSGetStages_RK; 1447e4dd521cSBarry Smith 14483ca0882eSHong Zhang ts->ops->snesfunction = SNESTSFormFunction_RK; 14493ca0882eSHong Zhang ts->ops->snesjacobian = SNESTSFormJacobian_RK; 14500f7a1166SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_RK; 145113af1a74SHong Zhang ts->ops->adjointsetup = TSAdjointSetUp_RK; 145213af1a74SHong Zhang ts->ops->adjointstep = TSAdjointStep_RK; 145313af1a74SHong Zhang ts->ops->adjointreset = TSAdjointReset_RK; 14540f7a1166SHong Zhang 145513af1a74SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_RK; 1456922a638cSHong Zhang ts->ops->forwardsetup = TSForwardSetUp_RK; 1457922a638cSHong Zhang ts->ops->forwardreset = TSForwardReset_RK; 1458922a638cSHong Zhang ts->ops->forwardstep = TSForwardStep_RK; 1459922a638cSHong Zhang ts->ops->forwardgetstages = TSForwardGetStages_RK; 1460922a638cSHong Zhang 1461*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&rk)); 14621566a47fSLisandro Dalcin ts->data = (void *)rk; 1463e4dd521cSBarry Smith 14649566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetOrder_C", TSRKGetOrder_RK)); 14659566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetType_C", TSRKGetType_RK)); 14669566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKSetType_C", TSRKSetType_RK)); 14679566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetTableau_C", TSRKGetTableau_RK)); 14689566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKSetMultirate_C", TSRKSetMultirate_RK)); 14699566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRKGetMultirate_C", TSRKGetMultirate_RK)); 1470be5899b3SLisandro Dalcin 14719566063dSJacob Faibussowitsch PetscCall(TSRKSetType(ts, TSRKDefault)); 147290b67129SHong Zhang rk->dtratio = 1; 1473e4dd521cSBarry Smith PetscFunctionReturn(0); 1474e4dd521cSBarry Smith } 1475