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: 289bd3e852SBarry Smith . -ts_rk_type 1fe 29287c2655SBarry Smith 30f68a32c8SEmil Constantinescu Level: advanced 31f68a32c8SEmil Constantinescu 32287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 33f68a32c8SEmil Constantinescu M*/ 34f68a32c8SEmil Constantinescu /*MC 352109b73fSDebojyoti Ghosh TSRK2A - Second order RK scheme. 36f68a32c8SEmil Constantinescu 37f68a32c8SEmil Constantinescu This method has two stages. 38f68a32c8SEmil Constantinescu 39287c2655SBarry Smith Options database: 409bd3e852SBarry Smith . -ts_rk_type 2a 41287c2655SBarry Smith 42f68a32c8SEmil Constantinescu Level: advanced 43f68a32c8SEmil Constantinescu 44287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 45f68a32c8SEmil Constantinescu M*/ 46f68a32c8SEmil Constantinescu /*MC 47f68a32c8SEmil Constantinescu TSRK3 - Third order RK scheme. 48f68a32c8SEmil Constantinescu 49f68a32c8SEmil Constantinescu This method has three stages. 50f68a32c8SEmil Constantinescu 51287c2655SBarry Smith Options database: 529bd3e852SBarry Smith . -ts_rk_type 3 53287c2655SBarry Smith 54f68a32c8SEmil Constantinescu Level: advanced 55f68a32c8SEmil Constantinescu 56287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 57f68a32c8SEmil Constantinescu M*/ 58f68a32c8SEmil Constantinescu /*MC 592109b73fSDebojyoti Ghosh TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method. 602109b73fSDebojyoti Ghosh 61d1905564SLisandro Dalcin This method has four stages with the First Same As Last (FSAL) property. 622109b73fSDebojyoti Ghosh 63287c2655SBarry Smith Options database: 649bd3e852SBarry Smith . -ts_rk_type 3bs 65287c2655SBarry Smith 662109b73fSDebojyoti Ghosh Level: advanced 672109b73fSDebojyoti Ghosh 6898164e13SLisandro Dalcin References: https://doi.org/10.1016/0893-9659(89)90079-7 69d1905564SLisandro Dalcin 70287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 712109b73fSDebojyoti Ghosh M*/ 722109b73fSDebojyoti Ghosh /*MC 73f68a32c8SEmil Constantinescu TSRK4 - Fourth order RK scheme. 74f68a32c8SEmil Constantinescu 752e698f8bSDebojyoti Ghosh This is the classical Runge-Kutta method with four stages. 76f68a32c8SEmil Constantinescu 77287c2655SBarry Smith Options database: 789bd3e852SBarry Smith . -ts_rk_type 4 79287c2655SBarry Smith 80f68a32c8SEmil Constantinescu Level: advanced 81f68a32c8SEmil Constantinescu 82287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 83f68a32c8SEmil Constantinescu M*/ 84f68a32c8SEmil Constantinescu /*MC 852e698f8bSDebojyoti Ghosh TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method. 86f68a32c8SEmil Constantinescu 87f68a32c8SEmil Constantinescu This method has six stages. 88f68a32c8SEmil Constantinescu 89287c2655SBarry Smith Options database: 909bd3e852SBarry Smith . -ts_rk_type 5f 91287c2655SBarry Smith 92f68a32c8SEmil Constantinescu Level: advanced 93f68a32c8SEmil Constantinescu 94287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 95f68a32c8SEmil Constantinescu M*/ 962109b73fSDebojyoti Ghosh /*MC 972e698f8bSDebojyoti Ghosh TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method. 982109b73fSDebojyoti Ghosh 99d1905564SLisandro Dalcin This method has seven stages with the First Same As Last (FSAL) property. 1002109b73fSDebojyoti Ghosh 101287c2655SBarry Smith Options database: 1029bd3e852SBarry Smith . -ts_rk_type 5dp 103287c2655SBarry Smith 1042109b73fSDebojyoti Ghosh Level: advanced 1052109b73fSDebojyoti Ghosh 10698164e13SLisandro Dalcin References: https://doi.org/10.1016/0771-050X(80)90013-3 107d1905564SLisandro Dalcin 108287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 1092109b73fSDebojyoti Ghosh M*/ 11005e23783SLisandro Dalcin /*MC 11105e23783SLisandro Dalcin TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method. 11205e23783SLisandro Dalcin 11305e23783SLisandro Dalcin This method has eight stages with the First Same As Last (FSAL) property. 11405e23783SLisandro Dalcin 115287c2655SBarry Smith Options database: 1169bd3e852SBarry Smith . -ts_rk_type 5bs 117287c2655SBarry Smith 11805e23783SLisandro Dalcin Level: advanced 11905e23783SLisandro Dalcin 12098164e13SLisandro Dalcin References: https://doi.org/10.1016/0898-1221(96)00141-1 12105e23783SLisandro Dalcin 122287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 12305e23783SLisandro Dalcin M*/ 12437acaa02SHendrik Ranocha /*MC 12537acaa02SHendrik Ranocha TSRK6VR - Sixth order robust Verner RK scheme with fifth order embedded method. 12637acaa02SHendrik Ranocha 12737acaa02SHendrik Ranocha This method has nine stages with the First Same As Last (FSAL) property. 12837acaa02SHendrik Ranocha 12937acaa02SHendrik Ranocha Options database: 13037acaa02SHendrik Ranocha . -ts_rk_type 6vr 13137acaa02SHendrik Ranocha 13237acaa02SHendrik Ranocha Level: advanced 13337acaa02SHendrik Ranocha 13437acaa02SHendrik Ranocha References: http://people.math.sfu.ca/~jverner/RKV65.IIIXb.Robust.00010102836.081204.CoeffsOnlyRAT 13537acaa02SHendrik Ranocha 13637acaa02SHendrik Ranocha .seealso: TSRK, TSRKType, TSRKSetType() 13737acaa02SHendrik Ranocha M*/ 13837acaa02SHendrik Ranocha /*MC 13937acaa02SHendrik Ranocha TSRK7VR - Seventh order robust Verner RK scheme with sixth order embedded method. 14037acaa02SHendrik Ranocha 14137acaa02SHendrik Ranocha This method has ten stages with the First Same As Last (FSAL) property. 14237acaa02SHendrik Ranocha 14337acaa02SHendrik Ranocha Options database: 14437acaa02SHendrik Ranocha . -ts_rk_type 7vr 14537acaa02SHendrik Ranocha 14637acaa02SHendrik Ranocha Level: advanced 14737acaa02SHendrik Ranocha 14837acaa02SHendrik Ranocha References: http://people.math.sfu.ca/~jverner/RKV76.IIa.Robust.000027015646.081206.CoeffsOnlyRAT 14937acaa02SHendrik Ranocha 15037acaa02SHendrik Ranocha .seealso: TSRK, TSRKType, TSRKSetType() 15137acaa02SHendrik Ranocha M*/ 15237acaa02SHendrik Ranocha /*MC 15337acaa02SHendrik Ranocha TSRK8VR - Eigth order robust Verner RK scheme with seventh order embedded method. 15437acaa02SHendrik Ranocha 15537acaa02SHendrik Ranocha This method has thirteen stages with the First Same As Last (FSAL) property. 15637acaa02SHendrik Ranocha 15737acaa02SHendrik Ranocha Options database: 15837acaa02SHendrik Ranocha . -ts_rk_type 8vr 15937acaa02SHendrik Ranocha 16037acaa02SHendrik Ranocha Level: advanced 16137acaa02SHendrik Ranocha 16237acaa02SHendrik Ranocha References: http://people.math.sfu.ca/~jverner/RKV87.IIa.Robust.00000754677.081208.CoeffsOnlyRATandFLOAT 16337acaa02SHendrik Ranocha 16437acaa02SHendrik Ranocha .seealso: TSRK, TSRKType, TSRKSetType() 16537acaa02SHendrik Ranocha M*/ 166262ff7bbSBarry Smith 167f68a32c8SEmil Constantinescu /*@C 168f68a32c8SEmil Constantinescu TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK 169262ff7bbSBarry Smith 170f68a32c8SEmil Constantinescu Not Collective, but should be called by all processes which will need the schemes to be registered 171262ff7bbSBarry Smith 172f68a32c8SEmil Constantinescu Level: advanced 173262ff7bbSBarry Smith 174f68a32c8SEmil Constantinescu .keywords: TS, TSRK, register, all 175262ff7bbSBarry Smith 176f68a32c8SEmil Constantinescu .seealso: TSRKRegisterDestroy() 177262ff7bbSBarry Smith @*/ 178f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void) 179262ff7bbSBarry Smith { 1804ac538c5SBarry Smith PetscErrorCode ierr; 181262ff7bbSBarry Smith 182262ff7bbSBarry Smith PetscFunctionBegin; 183f68a32c8SEmil Constantinescu if (TSRKRegisterAllCalled) PetscFunctionReturn(0); 184f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_TRUE; 185e4dd521cSBarry Smith 1864e82626cSLisandro Dalcin #define RC PetscRealConstant 187e4dd521cSBarry Smith { 188f68a32c8SEmil Constantinescu const PetscReal 1894e82626cSLisandro Dalcin A[1][1] = {{0}}, 1904e82626cSLisandro Dalcin b[1] = {RC(1.0)}; 1913a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 192e4dd521cSBarry Smith } 193db046809SBarry Smith { 194f68a32c8SEmil Constantinescu const PetscReal 1954e82626cSLisandro Dalcin A[2][2] = {{0,0}, 1964e82626cSLisandro Dalcin {RC(1.0),0}}, 1974e82626cSLisandro Dalcin b[2] = {RC(0.5),RC(0.5)}, 1984e82626cSLisandro Dalcin bembed[2] = {RC(1.0),0}; 1993a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 200db046809SBarry Smith } 201f68a32c8SEmil Constantinescu { 202f68a32c8SEmil Constantinescu const PetscReal 20317f6384fSBarry Smith A[3][3] = {{0,0,0}, 2044e82626cSLisandro Dalcin {RC(2.0)/RC(3.0),0,0}, 2054e82626cSLisandro Dalcin {RC(-1.0)/RC(3.0),RC(1.0),0}}, 2064e82626cSLisandro Dalcin b[3] = {RC(0.25),RC(0.5),RC(0.25)}; 2073a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 208fdefd5e5SDebojyoti Ghosh } 209fdefd5e5SDebojyoti Ghosh { 210fdefd5e5SDebojyoti Ghosh const PetscReal 21117f6384fSBarry Smith A[4][4] = {{0,0,0,0}, 2124e82626cSLisandro Dalcin {RC(1.0)/RC(2.0),0,0,0}, 2134e82626cSLisandro Dalcin {0,RC(3.0)/RC(4.0),0,0}, 2144e82626cSLisandro Dalcin {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}}, 2154e82626cSLisandro Dalcin b[4] = {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}, 2164e82626cSLisandro Dalcin bembed[4] = {RC(7.0)/RC(24.0),RC(1.0)/RC(4.0),RC(1.0)/RC(3.0),RC(1.0)/RC(8.0)}; 2173a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 218db046809SBarry Smith } 219f68a32c8SEmil Constantinescu { 220f68a32c8SEmil Constantinescu const PetscReal 221f68a32c8SEmil Constantinescu A[4][4] = {{0,0,0,0}, 2224e82626cSLisandro Dalcin {RC(0.5),0,0,0}, 2234e82626cSLisandro Dalcin {0,RC(0.5),0,0}, 2244e82626cSLisandro Dalcin {0,0,RC(1.0),0}}, 2254e82626cSLisandro 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)}; 2263a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 227db046809SBarry Smith } 228f68a32c8SEmil Constantinescu { 229f68a32c8SEmil Constantinescu const PetscReal 23017f6384fSBarry Smith A[6][6] = {{0,0,0,0,0,0}, 2314e82626cSLisandro Dalcin {RC(0.25),0,0,0,0,0}, 2324e82626cSLisandro Dalcin {RC(3.0)/RC(32.0),RC(9.0)/RC(32.0),0,0,0,0}, 2334e82626cSLisandro Dalcin {RC(1932.0)/RC(2197.0),RC(-7200.0)/RC(2197.0),RC(7296.0)/RC(2197.0),0,0,0}, 2344e82626cSLisandro Dalcin {RC(439.0)/RC(216.0),RC(-8.0),RC(3680.0)/RC(513.0),RC(-845.0)/RC(4104.0),0,0}, 2354e82626cSLisandro Dalcin {RC(-8.0)/RC(27.0),RC(2.0),RC(-3544.0)/RC(2565.0),RC(1859.0)/RC(4104.0),RC(-11.0)/RC(40.0),0}}, 2364e82626cSLisandro 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)}, 2374e82626cSLisandro 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}; 2383a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 239fdefd5e5SDebojyoti Ghosh } 240fdefd5e5SDebojyoti Ghosh { 241fdefd5e5SDebojyoti Ghosh const PetscReal 24217f6384fSBarry Smith A[7][7] = {{0,0,0,0,0,0,0}, 2434e82626cSLisandro Dalcin {RC(1.0)/RC(5.0),0,0,0,0,0,0}, 2444e82626cSLisandro Dalcin {RC(3.0)/RC(40.0),RC(9.0)/RC(40.0),0,0,0,0,0}, 2454e82626cSLisandro Dalcin {RC(44.0)/RC(45.0),RC(-56.0)/RC(15.0),RC(32.0)/RC(9.0),0,0,0,0}, 2464e82626cSLisandro 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}, 2474e82626cSLisandro 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}, 2484e82626cSLisandro Dalcin {RC(35.0)/RC(384.0),0,RC(500.0)/RC(1113.0),RC(125.0)/RC(192.0),RC(-2187.0)/RC(6784.0),RC(11.0)/RC(84.0),0}}, 2494e82626cSLisandro 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}, 250a685a763Sluzhanghpp 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)}, 251a685a763Sluzhanghpp 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)}, 252a685a763Sluzhanghpp {0,0,0,0,0}, 253a685a763Sluzhanghpp {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)}, 254a685a763Sluzhanghpp {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)}, 255a685a763Sluzhanghpp {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)}, 256a685a763Sluzhanghpp {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)}, 257a685a763Sluzhanghpp {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)}}; 258a685a763Sluzhanghpp 259a685a763Sluzhanghpp ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,5,binterp[0]);CHKERRQ(ierr); 260f68a32c8SEmil Constantinescu } 26105e23783SLisandro Dalcin { 26205e23783SLisandro Dalcin const PetscReal 26317f6384fSBarry Smith A[8][8] = {{0,0,0,0,0,0,0,0}, 2644e82626cSLisandro Dalcin {RC(1.0)/RC(6.0),0,0,0,0,0,0,0}, 2654e82626cSLisandro Dalcin {RC(2.0)/RC(27.0),RC(4.0)/RC(27.0),0,0,0,0,0,0}, 2664e82626cSLisandro Dalcin {RC(183.0)/RC(1372.0),RC(-162.0)/RC(343.0),RC(1053.0)/RC(1372.0),0,0,0,0,0}, 2674e82626cSLisandro 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}, 2684e82626cSLisandro 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}, 2694e82626cSLisandro 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}, 2704e82626cSLisandro Dalcin {RC(587.0)/RC(8064.0),0,RC(4440339.0)/RC(15491840.0),RC(24353.0)/RC(124800.0),RC(387.0)/RC(44800.0),RC(2152.0)/RC(5985.0),RC(7267.0)/RC(94080.0),0}}, 2714e82626cSLisandro 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}, 2724e82626cSLisandro 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)}; 27305e23783SLisandro Dalcin ierr = TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 27405e23783SLisandro Dalcin } 27537acaa02SHendrik Ranocha { 27637acaa02SHendrik Ranocha const PetscReal 27737acaa02SHendrik Ranocha A[9][9] = {{0,0,0,0,0,0,0,0,0}, 278*63fe90e8SHendrik Ranocha {RC(1.8000000000000000000000000000000000000000e-01),0,0,0,0,0,0,0,0}, 279*63fe90e8SHendrik Ranocha {RC(8.9506172839506172839506172839506172839506e-02),RC(7.7160493827160493827160493827160493827160e-02),0,0,0,0,0,0,0}, 280*63fe90e8SHendrik Ranocha {RC(6.2500000000000000000000000000000000000000e-02),0,RC(1.8750000000000000000000000000000000000000e-01),0,0,0,0,0,0}, 281*63fe90e8SHendrik Ranocha {RC(3.1651600000000000000000000000000000000000e-01),0,RC(-1.0449480000000000000000000000000000000000e+00),RC(1.2584320000000000000000000000000000000000e+00),0,0,0,0,0}, 282*63fe90e8SHendrik Ranocha {RC(2.7232612736485626257225065566674305502508e-01),0,RC(-8.2513360323886639676113360323886639676113e-01),RC(1.0480917678812415654520917678812415654521e+00),RC(1.0471570799276856873679117969088177628396e-01),0,0,0,0}, 283*63fe90e8SHendrik 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}, 284*63fe90e8SHendrik 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}, 285*63fe90e8SHendrik Ranocha {RC(7.6388888888888888888888888888888888888889e-02),0,0,RC(3.6940836940836940836940836940836940836941e-01),0,RC(2.4801587301587301587301587301587301587302e-01),RC(2.3674242424242424242424242424242424242424e-01),RC(6.9444444444444444444444444444444444444444e-02),0}}, 286*63fe90e8SHendrik Ranocha b[9] = {RC(7.6388888888888888888888888888888888888889e-02),0,0,RC(3.6940836940836940836940836940836940836941e-01),0,RC(2.4801587301587301587301587301587301587302e-01),RC(2.3674242424242424242424242424242424242424e-01),RC(6.9444444444444444444444444444444444444444e-02),0}, 287*63fe90e8SHendrik 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)}; 28837acaa02SHendrik Ranocha ierr = TSRKRegister(TSRK6VR,6,9,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 28937acaa02SHendrik Ranocha } 29037acaa02SHendrik Ranocha { 29137acaa02SHendrik Ranocha const PetscReal 29237acaa02SHendrik Ranocha A[10][10] = {{0,0,0,0,0,0,0,0,0,0}, 293*63fe90e8SHendrik Ranocha {RC(5.0000000000000000000000000000000000000000e-03),0,0,0,0,0,0,0,0,0}, 294*63fe90e8SHendrik Ranocha {RC(-1.0767901234567901234567901234567901234568e+00),RC(1.1856790123456790123456790123456790123457e+00),0,0,0,0,0,0,0,0}, 295*63fe90e8SHendrik Ranocha {RC(4.0833333333333333333333333333333333333333e-02),0,RC(1.2250000000000000000000000000000000000000e-01),0,0,0,0,0,0,0}, 296*63fe90e8SHendrik Ranocha {RC(6.3607142857142857142857142857142857142857e-01),0,RC(-2.4444642857142857142857142857142857142857e+00),RC(2.2633928571428571428571428571428571428571e+00),0,0,0,0,0,0}, 297*63fe90e8SHendrik Ranocha {RC(-2.5351211079349245229256383554660215487207e+00),0,RC(1.0299374654449267920438514460756024913612e+01),RC(-7.9513032885990579949493217458266876536482e+00),RC(7.9301148923100592201226014271115261823800e-01),0,0,0,0,0}, 298*63fe90e8SHendrik 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}, 299*63fe90e8SHendrik 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}, 300*63fe90e8SHendrik 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}, 301*63fe90e8SHendrik Ranocha {RC(-1.4449518916777735137351003179355712360517e+00),0,RC(8.0318913859955919224117033223019560435041e+00),RC(-7.5831741663401346820798883023671588604984e+00),RC(3.5816169353190074211247685442452878696855e+00),RC(-2.4369722632199529411183809065693752383733e+00),RC(8.5158999992326179339689766032486142173390e-01),0,0,0}}, 302*63fe90e8SHendrik Ranocha b[10] = {RC(4.7425837833706756083569172717574534698932e-02),0,0,RC(2.5622361659370562659961727458274623448160e-01),RC(2.6951376833074206619473817258075952886764e-01),RC(1.2686622409092782845989138364739173247882e-01),RC(2.4887225942060071622046449427647492767292e-01),RC(3.0744837408200631335304388479099184768645e-03),RC(4.8023809989496943308189063347143123323209e-02),0}, 303*63fe90e8SHendrik Ranocha bembed[10] = {RC(4.7485247699299631037531273805727961552268e-02),0,0,RC(2.5599412588690633297154918245905393870497e-01),RC(2.7058478081067688722530891099268135732387e-01),RC(1.2505618684425992913638822323746917920448e-01),RC(2.5204468723743860507184043820197442562182e-01),0,0,RC(4.8834971521418614557381971303093137592592e-02)}; 30437acaa02SHendrik Ranocha ierr = TSRKRegister(TSRK7VR,7,10,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 30537acaa02SHendrik Ranocha } 30637acaa02SHendrik Ranocha { 30737acaa02SHendrik Ranocha const PetscReal 30837acaa02SHendrik Ranocha A[13][13] = {{0,0,0,0,0,0,0,0,0,0,0,0,0}, 309*63fe90e8SHendrik Ranocha {RC(2.5000000000000000000000000000000000000000e-01),0,0,0,0,0,0,0,0,0,0,0,0}, 310*63fe90e8SHendrik Ranocha {RC(8.7400846504915232052686327594877411977046e-02),RC(2.5487604938654321753087950620345685135815e-02),0,0,0,0,0,0,0,0,0,0,0}, 311*63fe90e8SHendrik Ranocha {RC(4.2333169291338582677165354330708661417323e-02),0,RC(1.2699950787401574803149606299212598425197e-01),0,0,0,0,0,0,0,0,0,0}, 312*63fe90e8SHendrik Ranocha {RC(4.2609505888742261494881445237572274090942e-01),0,RC(-1.5987952846591523265427733230657181117089e+00),RC(1.5967002257717297115939588706899953707994e+00),0,0,0,0,0,0,0,0,0}, 313*63fe90e8SHendrik Ranocha {RC(5.0719337296713929515090618138513639239329e-02),0,0,RC(2.5433377264600407582754714408877778031369e-01),RC(2.0394689005728199465736223777270858044698e-01),0,0,0,0,0,0,0,0}, 314*63fe90e8SHendrik 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}, 315*63fe90e8SHendrik 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}, 316*63fe90e8SHendrik 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}, 317*63fe90e8SHendrik 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}, 318*63fe90e8SHendrik 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}, 319*63fe90e8SHendrik 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}, 320*63fe90e8SHendrik Ranocha {RC(8.8735762208534719663211694051981022704884e-01),0,0,RC(-2.9754597821085367558513632804709301581977e+00),RC(5.6007170094881630597990392548350098923829e+00),RC(-5.9156074505366744680014930189941657351840e+00),RC(2.2029689156134927016879142540807638331238e-01),RC(1.0155097824462216666143271340902996997549e-01),RC(1.1514345647386055909780397752125850553556e+00),RC(1.9297101665271239396134361900805843653065e+00),0,0,0}}, 321*63fe90e8SHendrik Ranocha b[13] = {RC(4.4729564666695714203015840429049382466467e-02),0,0,0,0,RC(1.5691033527708199813368698010726645409175e-01),RC(1.8460973408151637740702451873526277892035e-01),RC(2.2516380602086991042479419400350721970920e-01),RC(1.4794615651970234687005179885449141753736e-01),RC(7.6055542444955825269798361910336491012732e-02),RC(1.2277290235018619610824346315921437388535e-01),RC(4.1811958638991631583384842800871882376786e-02),0}, 322*63fe90e8SHendrik Ranocha bembed[13] = {RC(4.5847111400495925878664730122010282095875e-02),0,0,0,0,RC(2.6231891404152387437443356584845803392392e-01),RC(1.9169372337852611904485738635688429008025e-01),RC(2.1709172327902618330978407422906448568196e-01),RC(1.2738189624833706796803169450656737867900e-01),RC(1.1510530385365326258240515750043192148894e-01),0,0,RC(4.0561327798437566841823391436583608050053e-02)}; 32337acaa02SHendrik Ranocha ierr = TSRKRegister(TSRK8VR,8,13,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 32437acaa02SHendrik Ranocha } 3254e82626cSLisandro Dalcin #undef RC 326db046809SBarry Smith PetscFunctionReturn(0); 327db046809SBarry Smith } 328db046809SBarry Smith 329f68a32c8SEmil Constantinescu /*@C 330f68a32c8SEmil Constantinescu TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister(). 331f68a32c8SEmil Constantinescu 332f68a32c8SEmil Constantinescu Not Collective 333f68a32c8SEmil Constantinescu 334f68a32c8SEmil Constantinescu Level: advanced 335f68a32c8SEmil Constantinescu 336f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy 337f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll() 338f68a32c8SEmil Constantinescu @*/ 339f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void) 340e4dd521cSBarry Smith { 341f68a32c8SEmil Constantinescu PetscErrorCode ierr; 342f68a32c8SEmil Constantinescu RKTableauLink link; 343f68a32c8SEmil Constantinescu 344f68a32c8SEmil Constantinescu PetscFunctionBegin; 345f68a32c8SEmil Constantinescu while ((link = RKTableauList)) { 346f68a32c8SEmil Constantinescu RKTableau t = &link->tab; 347f68a32c8SEmil Constantinescu RKTableauList = link->next; 348f68a32c8SEmil Constantinescu ierr = PetscFree3(t->A,t->b,t->c); CHKERRQ(ierr); 349f68a32c8SEmil Constantinescu ierr = PetscFree (t->bembed); CHKERRQ(ierr); 350f68a32c8SEmil Constantinescu ierr = PetscFree (t->binterp); CHKERRQ(ierr); 351f68a32c8SEmil Constantinescu ierr = PetscFree (t->name); CHKERRQ(ierr); 352f68a32c8SEmil Constantinescu ierr = PetscFree (link); CHKERRQ(ierr); 353f68a32c8SEmil Constantinescu } 354f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_FALSE; 355f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 356f68a32c8SEmil Constantinescu } 357f68a32c8SEmil Constantinescu 358f68a32c8SEmil Constantinescu /*@C 359f68a32c8SEmil Constantinescu TSRKInitializePackage - This function initializes everything in the TSRK package. It is called 3608a690491SBarry Smith from TSInitializePackage(). 361f68a32c8SEmil Constantinescu 362f68a32c8SEmil Constantinescu Level: developer 363f68a32c8SEmil Constantinescu 364f68a32c8SEmil Constantinescu .keywords: TS, TSRK, initialize, package 365f68a32c8SEmil Constantinescu .seealso: PetscInitialize() 366f68a32c8SEmil Constantinescu @*/ 367f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void) 368f68a32c8SEmil Constantinescu { 369186e87acSLisandro Dalcin PetscErrorCode ierr; 370e4dd521cSBarry Smith 371e4dd521cSBarry Smith PetscFunctionBegin; 372f68a32c8SEmil Constantinescu if (TSRKPackageInitialized) PetscFunctionReturn(0); 373f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_TRUE; 374f68a32c8SEmil Constantinescu ierr = TSRKRegisterAll();CHKERRQ(ierr); 375f68a32c8SEmil Constantinescu ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr); 376f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 377f68a32c8SEmil Constantinescu } 378186e87acSLisandro Dalcin 379f68a32c8SEmil Constantinescu /*@C 380f68a32c8SEmil Constantinescu TSRKFinalizePackage - This function destroys everything in the TSRK package. It is 381f68a32c8SEmil Constantinescu called from PetscFinalize(). 382186e87acSLisandro Dalcin 383f68a32c8SEmil Constantinescu Level: developer 384f68a32c8SEmil Constantinescu 385f68a32c8SEmil Constantinescu .keywords: Petsc, destroy, package 386f68a32c8SEmil Constantinescu .seealso: PetscFinalize() 387f68a32c8SEmil Constantinescu @*/ 388f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void) 389f68a32c8SEmil Constantinescu { 390f68a32c8SEmil Constantinescu PetscErrorCode ierr; 391f68a32c8SEmil Constantinescu 392f68a32c8SEmil Constantinescu PetscFunctionBegin; 393f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_FALSE; 394f68a32c8SEmil Constantinescu ierr = TSRKRegisterDestroy();CHKERRQ(ierr); 395f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 396f68a32c8SEmil Constantinescu } 397f68a32c8SEmil Constantinescu 398f68a32c8SEmil Constantinescu /*@C 399f68a32c8SEmil Constantinescu TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 400f68a32c8SEmil Constantinescu 401f68a32c8SEmil Constantinescu Not Collective, but the same schemes should be registered on all processes on which they will be used 402f68a32c8SEmil Constantinescu 403f68a32c8SEmil Constantinescu Input Parameters: 404f68a32c8SEmil Constantinescu + name - identifier for method 405f68a32c8SEmil Constantinescu . order - approximation order of method 406f68a32c8SEmil Constantinescu . s - number of stages, this is the dimension of the matrices below 407f68a32c8SEmil Constantinescu . A - stage coefficients (dimension s*s, row-major) 408f68a32c8SEmil Constantinescu . b - step completion table (dimension s; NULL to use last row of A) 409f68a32c8SEmil Constantinescu . c - abscissa (dimension s; NULL to use row sums of A) 410f68a32c8SEmil Constantinescu . bembed - completion table for embedded method (dimension s; NULL if not available) 4113a8a9803SLisandro Dalcin . p - Order of the interpolation scheme, equal to the number of columns of binterp 4123a8a9803SLisandro Dalcin - binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1) 413f68a32c8SEmil Constantinescu 414f68a32c8SEmil Constantinescu Notes: 415f68a32c8SEmil Constantinescu Several RK methods are provided, this function is only needed to create new methods. 416f68a32c8SEmil Constantinescu 417f68a32c8SEmil Constantinescu Level: advanced 418f68a32c8SEmil Constantinescu 419f68a32c8SEmil Constantinescu .keywords: TS, register 420f68a32c8SEmil Constantinescu 421f68a32c8SEmil Constantinescu .seealso: TSRK 422f68a32c8SEmil Constantinescu @*/ 423f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s, 424f68a32c8SEmil Constantinescu const PetscReal A[],const PetscReal b[],const PetscReal c[], 4253a8a9803SLisandro Dalcin const PetscReal bembed[],PetscInt p,const PetscReal binterp[]) 426f68a32c8SEmil Constantinescu { 427f68a32c8SEmil Constantinescu PetscErrorCode ierr; 428f68a32c8SEmil Constantinescu RKTableauLink link; 429f68a32c8SEmil Constantinescu RKTableau t; 430f68a32c8SEmil Constantinescu PetscInt i,j; 431f68a32c8SEmil Constantinescu 432f68a32c8SEmil Constantinescu PetscFunctionBegin; 4333a8a9803SLisandro Dalcin PetscValidCharPointer(name,1); 4343a8a9803SLisandro Dalcin PetscValidRealPointer(A,4); 4353a8a9803SLisandro Dalcin if (b) PetscValidRealPointer(b,5); 4363a8a9803SLisandro Dalcin if (c) PetscValidRealPointer(c,6); 4373a8a9803SLisandro Dalcin if (bembed) PetscValidRealPointer(bembed,7); 4383a8a9803SLisandro Dalcin if (binterp || p > 1) PetscValidRealPointer(binterp,9); 4393a8a9803SLisandro Dalcin 4401d36bdfdSBarry Smith ierr = TSRKInitializePackage();CHKERRQ(ierr); 44195dccacaSBarry Smith ierr = PetscNew(&link);CHKERRQ(ierr); 442f68a32c8SEmil Constantinescu t = &link->tab; 4433a8a9803SLisandro Dalcin 444f68a32c8SEmil Constantinescu ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 445f68a32c8SEmil Constantinescu t->order = order; 446f68a32c8SEmil Constantinescu t->s = s; 447dcca6d9dSJed Brown ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr); 448f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 449f68a32c8SEmil Constantinescu if (b) { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); } 450f68a32c8SEmil Constantinescu else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i]; 451f68a32c8SEmil Constantinescu if (c) { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); } 452f68a32c8SEmil Constantinescu else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j]; 453d760c35bSDebojyoti Ghosh t->FSAL = PETSC_TRUE; 454d760c35bSDebojyoti Ghosh for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE; 4553a8a9803SLisandro Dalcin 456f68a32c8SEmil Constantinescu if (bembed) { 457785e854fSJed Brown ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr); 458f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr); 459f68a32c8SEmil Constantinescu } 460f68a32c8SEmil Constantinescu 4613a8a9803SLisandro Dalcin if (!binterp) { p = 1; binterp = t->b; } 4623a8a9803SLisandro Dalcin t->p = p; 4633a8a9803SLisandro Dalcin ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr); 4643a8a9803SLisandro Dalcin ierr = PetscMemcpy(t->binterp,binterp,s*p*sizeof(binterp[0]));CHKERRQ(ierr); 4653a8a9803SLisandro Dalcin 466f68a32c8SEmil Constantinescu link->next = RKTableauList; 467f68a32c8SEmil Constantinescu RKTableauList = link; 468f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 469f68a32c8SEmil Constantinescu } 470f68a32c8SEmil Constantinescu 471e4dd521cSBarry Smith /* 4722c9cb062Sluzhanghpp This is for single-step RK method 473f68a32c8SEmil Constantinescu The step completion formula is 474e4dd521cSBarry Smith 475f68a32c8SEmil Constantinescu x1 = x0 + h b^T YdotRHS 476f68a32c8SEmil Constantinescu 477f68a32c8SEmil Constantinescu This function can be called before or after ts->vec_sol has been updated. 478f68a32c8SEmil Constantinescu Suppose we have a completion formula (b) and an embedded formula (be) of different order. 479f68a32c8SEmil Constantinescu We can write 480f68a32c8SEmil Constantinescu 481f68a32c8SEmil Constantinescu x1e = x0 + h be^T YdotRHS 482f68a32c8SEmil Constantinescu = x1 - h b^T YdotRHS + h be^T YdotRHS 483f68a32c8SEmil Constantinescu = x1 + h (be - b)^T YdotRHS 484f68a32c8SEmil Constantinescu 485f68a32c8SEmil Constantinescu so we can evaluate the method with different order even after the step has been optimistically completed. 486e4dd521cSBarry Smith */ 487f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done) 488f68a32c8SEmil Constantinescu { 489f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 490f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 491f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 492f68a32c8SEmil Constantinescu PetscReal h; 493f68a32c8SEmil Constantinescu PetscInt s = tab->s,j; 494f68a32c8SEmil Constantinescu PetscErrorCode ierr; 495f68a32c8SEmil Constantinescu 496f68a32c8SEmil Constantinescu PetscFunctionBegin; 497f68a32c8SEmil Constantinescu switch (rk->status) { 498f68a32c8SEmil Constantinescu case TS_STEP_INCOMPLETE: 499f68a32c8SEmil Constantinescu case TS_STEP_PENDING: 500f68a32c8SEmil Constantinescu h = ts->time_step; break; 501f68a32c8SEmil Constantinescu case TS_STEP_COMPLETE: 502be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 503f68a32c8SEmil Constantinescu default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 504f68a32c8SEmil Constantinescu } 505f68a32c8SEmil Constantinescu if (order == tab->order) { 506f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { 507f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 50890b67129SHong Zhang for (j=0; j<s; j++) w[j] = h*tab->b[j]/rk->dtratio; 509f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 510f68a32c8SEmil Constantinescu } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 511f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 512f68a32c8SEmil Constantinescu } else if (order == tab->order-1) { 513f68a32c8SEmil Constantinescu if (!tab->bembed) goto unavailable; 514f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { /*Complete with the embedded method (be)*/ 515f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 516f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 517f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 518f68a32c8SEmil Constantinescu } else { /*Rollback and re-complete using (be-b) */ 519f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 520f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 521f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 522f68a32c8SEmil Constantinescu } 523f68a32c8SEmil Constantinescu if (done) *done = PETSC_TRUE; 524f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 525f68a32c8SEmil Constantinescu } 526f68a32c8SEmil Constantinescu unavailable: 527f68a32c8SEmil Constantinescu if (done) *done = PETSC_FALSE; 528a7fac7c2SEmil Constantinescu else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"RK '%s' of order %D cannot evaluate step at order %D. Consider using -ts_adapt_type none or a different method that has an embedded estimate.",tab->name,tab->order,order); 529f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 530f68a32c8SEmil Constantinescu } 531f68a32c8SEmil Constantinescu 5320f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts) 5330f7a1166SHong Zhang { 5340f7a1166SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 535cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 5360f7a1166SHong Zhang RKTableau tab = rk->tableau; 5370f7a1166SHong Zhang const PetscInt s = tab->s; 5380f7a1166SHong Zhang const PetscReal *b = tab->b,*c = tab->c; 5390f7a1166SHong Zhang Vec *Y = rk->Y; 5400f7a1166SHong Zhang PetscInt i; 5410f7a1166SHong Zhang PetscErrorCode ierr; 5420f7a1166SHong Zhang 5430f7a1166SHong Zhang PetscFunctionBegin; 544cd4cee2dSHong Zhang /* No need to backup quadts->vec_sol since it can be reverted in TSRollBack_RK */ 5450f7a1166SHong Zhang for (i=s-1; i>=0; i--) { 546cd4cee2dSHong Zhang /* Evolve quadrature TS solution to compute integrals */ 547cd4cee2dSHong Zhang ierr = TSComputeRHSFunction(quadts,rk->ptime+rk->time_step*c[i],Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 548cd4cee2dSHong Zhang ierr = VecAXPY(quadts->vec_sol,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 5490f7a1166SHong Zhang } 5500f7a1166SHong Zhang PetscFunctionReturn(0); 5510f7a1166SHong Zhang } 5520f7a1166SHong Zhang 5530f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts) 5540f7a1166SHong Zhang { 5550f7a1166SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 5560f7a1166SHong Zhang RKTableau tab = rk->tableau; 557cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 5580f7a1166SHong Zhang const PetscInt s = tab->s; 5590f7a1166SHong Zhang const PetscReal *b = tab->b,*c = tab->c; 5600f7a1166SHong Zhang Vec *Y = rk->Y; 5610f7a1166SHong Zhang PetscInt i; 5620f7a1166SHong Zhang PetscErrorCode ierr; 5630f7a1166SHong Zhang 5640f7a1166SHong Zhang PetscFunctionBegin; 5650f7a1166SHong Zhang for (i=s-1; i>=0; i--) { 566cd4cee2dSHong Zhang /* Evolve quadrature TS solution to compute integrals */ 567cd4cee2dSHong Zhang ierr = TSComputeRHSFunction(quadts,ts->ptime+ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 568cd4cee2dSHong Zhang ierr = VecAXPY(quadts->vec_sol,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 5690f7a1166SHong Zhang } 5700f7a1166SHong Zhang PetscFunctionReturn(0); 5710f7a1166SHong Zhang } 5720f7a1166SHong Zhang 573fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts) 574fecfb714SLisandro Dalcin { 575fecfb714SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 576cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 577fecfb714SLisandro Dalcin RKTableau tab = rk->tableau; 578fecfb714SLisandro Dalcin const PetscInt s = tab->s; 579cd4cee2dSHong Zhang const PetscReal *b = tab->b,*c = tab->c; 580fecfb714SLisandro Dalcin PetscScalar *w = rk->work; 581cd4cee2dSHong Zhang Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 582fecfb714SLisandro Dalcin PetscInt j; 583fecfb714SLisandro Dalcin PetscReal h; 584fecfb714SLisandro Dalcin PetscErrorCode ierr; 585fecfb714SLisandro Dalcin 586fecfb714SLisandro Dalcin PetscFunctionBegin; 587fecfb714SLisandro Dalcin switch (rk->status) { 588fecfb714SLisandro Dalcin case TS_STEP_INCOMPLETE: 589fecfb714SLisandro Dalcin case TS_STEP_PENDING: 590fecfb714SLisandro Dalcin h = ts->time_step; break; 591fecfb714SLisandro Dalcin case TS_STEP_COMPLETE: 592fecfb714SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 593fecfb714SLisandro Dalcin default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 594fecfb714SLisandro Dalcin } 595fecfb714SLisandro Dalcin for (j=0; j<s; j++) w[j] = -h*b[j]; 596fecfb714SLisandro Dalcin ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr); 597cd4cee2dSHong Zhang if (quadts && ts->costintegralfwd) { 598cd4cee2dSHong Zhang for (j=0; j<s; j++) { 599cd4cee2dSHong Zhang /* Revert the quadrature TS solution */ 600cd4cee2dSHong Zhang ierr = TSComputeRHSFunction(quadts,rk->ptime+h*c[j],Y[j],ts->vec_costintegrand);CHKERRQ(ierr); 601cd4cee2dSHong Zhang ierr = VecAXPY(quadts->vec_sol,-h*b[j],ts->vec_costintegrand);CHKERRQ(ierr); 602cd4cee2dSHong Zhang } 603cd4cee2dSHong Zhang } 604fecfb714SLisandro Dalcin PetscFunctionReturn(0); 605fecfb714SLisandro Dalcin } 606fecfb714SLisandro Dalcin 607922a638cSHong Zhang static PetscErrorCode TSForwardStep_RK(TS ts) 608922a638cSHong Zhang { 609922a638cSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 610922a638cSHong Zhang RKTableau tab = rk->tableau; 611922a638cSHong Zhang Mat J,*MatsFwdSensipTemp = rk->MatsFwdSensipTemp; 612922a638cSHong Zhang const PetscInt s = tab->s; 613922a638cSHong Zhang const PetscReal *A = tab->A,*c = tab->c,*b = tab->b; 614922a638cSHong Zhang Vec *Y = rk->Y; 615922a638cSHong Zhang PetscInt i,j; 616922a638cSHong Zhang PetscReal stage_time,h = ts->time_step; 617922a638cSHong Zhang PetscBool zero; 618922a638cSHong Zhang PetscErrorCode ierr; 619922a638cSHong Zhang 620922a638cSHong Zhang PetscFunctionBegin; 621922a638cSHong Zhang ierr = MatCopy(ts->mat_sensip,rk->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 622922a638cSHong Zhang ierr = TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr); 623922a638cSHong Zhang 624922a638cSHong Zhang for (i=0; i<s; i++) { 625922a638cSHong Zhang stage_time = ts->ptime + h*c[i]; 626922a638cSHong Zhang zero = PETSC_FALSE; 627922a638cSHong Zhang if (b[i] == 0 && i == s-1) zero = PETSC_TRUE; 628922a638cSHong Zhang /* TLM Stage values */ 629922a638cSHong Zhang if(!i) { 630922a638cSHong Zhang ierr = MatCopy(ts->mat_sensip,rk->MatsFwdStageSensip[i],SAME_NONZERO_PATTERN);CHKERRQ(ierr); 631922a638cSHong Zhang } else if (!zero) { 632922a638cSHong Zhang ierr = MatZeroEntries(rk->MatsFwdStageSensip[i]);CHKERRQ(ierr); 633922a638cSHong Zhang for (j=0; j<i; j++) { 634922a638cSHong Zhang ierr = MatAXPY(rk->MatsFwdStageSensip[i],h*A[i*s+j],MatsFwdSensipTemp[j],SAME_NONZERO_PATTERN);CHKERRQ(ierr); 635922a638cSHong Zhang } 636922a638cSHong Zhang ierr = MatAXPY(rk->MatsFwdStageSensip[i],1.,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 637922a638cSHong Zhang } else { 638922a638cSHong Zhang ierr = MatZeroEntries(rk->MatsFwdStageSensip[i]);CHKERRQ(ierr); 639922a638cSHong Zhang } 640922a638cSHong Zhang 641922a638cSHong Zhang ierr = TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);CHKERRQ(ierr); 642922a638cSHong Zhang ierr = MatMatMult(J,rk->MatsFwdStageSensip[i],MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatsFwdSensipTemp[i]);CHKERRQ(ierr); 64313af1a74SHong Zhang if (ts->Jacprhs) { 64413af1a74SHong Zhang ierr = TSComputeRHSJacobianP(ts,stage_time,Y[i],ts->Jacprhs);CHKERRQ(ierr); /* get f_p */ 64513af1a74SHong Zhang if (ts->vecs_sensi2p) { /* TLM used for 2nd-order adjoint */ 64613af1a74SHong Zhang PetscScalar *xarr; 64713af1a74SHong Zhang ierr = MatDenseGetColumn(MatsFwdSensipTemp[i],0,&xarr);CHKERRQ(ierr); 64813af1a74SHong Zhang ierr = VecPlaceArray(rk->VecDeltaFwdSensipCol,xarr);CHKERRQ(ierr); 64913af1a74SHong Zhang ierr = MatMultAdd(ts->Jacprhs,ts->vec_dir,rk->VecDeltaFwdSensipCol,rk->VecDeltaFwdSensipCol);CHKERRQ(ierr); 65013af1a74SHong Zhang ierr = VecResetArray(rk->VecDeltaFwdSensipCol);CHKERRQ(ierr);CHKERRQ(ierr); 65113af1a74SHong Zhang ierr = MatDenseRestoreColumn(MatsFwdSensipTemp[i],&xarr);CHKERRQ(ierr); 65213af1a74SHong Zhang } else { 65313af1a74SHong Zhang ierr = MatAXPY(MatsFwdSensipTemp[i],1.,ts->Jacprhs,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 65413af1a74SHong Zhang } 655922a638cSHong Zhang } 656922a638cSHong Zhang } 657922a638cSHong Zhang 658922a638cSHong Zhang for (i=0; i<s; i++) { 659922a638cSHong Zhang ierr = MatAXPY(ts->mat_sensip,h*b[i],rk->MatsFwdSensipTemp[i],SAME_NONZERO_PATTERN);CHKERRQ(ierr); 660922a638cSHong Zhang } 661922a638cSHong Zhang rk->status = TS_STEP_COMPLETE; 662922a638cSHong Zhang PetscFunctionReturn(0); 663922a638cSHong Zhang } 664922a638cSHong Zhang 665922a638cSHong Zhang static PetscErrorCode TSForwardGetStages_RK(TS ts,PetscInt *ns,Mat **stagesensip) 666922a638cSHong Zhang { 667922a638cSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 668922a638cSHong Zhang RKTableau tab = rk->tableau; 669922a638cSHong Zhang 670922a638cSHong Zhang PetscFunctionBegin; 671922a638cSHong Zhang if (ns) *ns = tab->s; 672922a638cSHong Zhang if (stagesensip) *stagesensip = rk->MatsFwdStageSensip; 673922a638cSHong Zhang PetscFunctionReturn(0); 674922a638cSHong Zhang } 675922a638cSHong Zhang 676922a638cSHong Zhang static PetscErrorCode TSForwardSetUp_RK(TS ts) 677922a638cSHong Zhang { 678922a638cSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 679922a638cSHong Zhang RKTableau tab = rk->tableau; 680922a638cSHong Zhang PetscInt i; 681922a638cSHong Zhang PetscErrorCode ierr; 682922a638cSHong Zhang 683922a638cSHong Zhang PetscFunctionBegin; 684922a638cSHong Zhang /* backup sensitivity results for roll-backs */ 685922a638cSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatFwdSensip0);CHKERRQ(ierr); 686922a638cSHong Zhang 687922a638cSHong Zhang ierr = PetscMalloc1(tab->s,&rk->MatsFwdStageSensip);CHKERRQ(ierr); 688922a638cSHong Zhang ierr = PetscMalloc1(tab->s,&rk->MatsFwdSensipTemp);CHKERRQ(ierr); 689922a638cSHong Zhang for(i=0; i<tab->s; i++) { 690922a638cSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatsFwdStageSensip[i]);CHKERRQ(ierr); 691922a638cSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatsFwdSensipTemp[i]);CHKERRQ(ierr); 692922a638cSHong Zhang } 693922a638cSHong Zhang ierr = VecDuplicate(ts->vec_sol,&rk->VecDeltaFwdSensipCol);CHKERRQ(ierr); 694922a638cSHong Zhang PetscFunctionReturn(0); 695922a638cSHong Zhang } 696922a638cSHong Zhang 697922a638cSHong Zhang static PetscErrorCode TSForwardReset_RK(TS ts) 698922a638cSHong Zhang { 699922a638cSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 700922a638cSHong Zhang RKTableau tab = rk->tableau; 701922a638cSHong Zhang PetscInt i; 702922a638cSHong Zhang PetscErrorCode ierr; 703922a638cSHong Zhang 704922a638cSHong Zhang PetscFunctionBegin; 705922a638cSHong Zhang ierr = MatDestroy(&rk->MatFwdSensip0);CHKERRQ(ierr); 706922a638cSHong Zhang if (rk->MatsFwdStageSensip) { 707922a638cSHong Zhang for (i=0; i<tab->s; i++) { 708922a638cSHong Zhang ierr = MatDestroy(&rk->MatsFwdStageSensip[i]);CHKERRQ(ierr); 709922a638cSHong Zhang } 710922a638cSHong Zhang ierr = PetscFree(rk->MatsFwdStageSensip);CHKERRQ(ierr); 711922a638cSHong Zhang } 712922a638cSHong Zhang if (rk->MatsFwdSensipTemp) { 713922a638cSHong Zhang for (i=0; i<tab->s; i++) { 714922a638cSHong Zhang ierr = MatDestroy(&rk->MatsFwdSensipTemp[i]);CHKERRQ(ierr); 715922a638cSHong Zhang } 716922a638cSHong Zhang ierr = PetscFree(rk->MatsFwdSensipTemp);CHKERRQ(ierr); 717922a638cSHong Zhang } 718922a638cSHong Zhang ierr = VecDestroy(&rk->VecDeltaFwdSensipCol);CHKERRQ(ierr); 719922a638cSHong Zhang PetscFunctionReturn(0); 720922a638cSHong Zhang } 721922a638cSHong Zhang 722f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts) 723f68a32c8SEmil Constantinescu { 724f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 725f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 726f68a32c8SEmil Constantinescu const PetscInt s = tab->s; 727fecfb714SLisandro Dalcin const PetscReal *A = tab->A,*c = tab->c; 728f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 729f68a32c8SEmil Constantinescu Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 730d1905564SLisandro Dalcin PetscBool FSAL = tab->FSAL; 731f68a32c8SEmil Constantinescu TSAdapt adapt; 732fecfb714SLisandro Dalcin PetscInt i,j; 733be5899b3SLisandro Dalcin PetscInt rejections = 0; 734be5899b3SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 735be5899b3SLisandro Dalcin PetscReal next_time_step = ts->time_step; 736f68a32c8SEmil Constantinescu PetscErrorCode ierr; 737f68a32c8SEmil Constantinescu 738f68a32c8SEmil Constantinescu PetscFunctionBegin; 739d1905564SLisandro Dalcin if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE; 740d1905564SLisandro Dalcin if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); } 741d1905564SLisandro Dalcin 742f68a32c8SEmil Constantinescu rk->status = TS_STEP_INCOMPLETE; 743be5899b3SLisandro Dalcin while (!ts->reason && rk->status != TS_STEP_COMPLETE) { 744be5899b3SLisandro Dalcin PetscReal t = ts->ptime; 745f68a32c8SEmil Constantinescu PetscReal h = ts->time_step; 746f68a32c8SEmil Constantinescu for (i=0; i<s; i++) { 7479be3e283SDebojyoti Ghosh rk->stage_time = t + h*c[i]; 7489be3e283SDebojyoti Ghosh ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr); 749f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 750f68a32c8SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 751f68a32c8SEmil Constantinescu ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 7529be3e283SDebojyoti Ghosh ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 753f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 754be5899b3SLisandro Dalcin ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr); 755fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 7568f738a7cSLisandro Dalcin if (FSAL && !i) continue; 757f68a32c8SEmil Constantinescu ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 758f68a32c8SEmil Constantinescu } 759be5899b3SLisandro Dalcin 760be5899b3SLisandro Dalcin rk->status = TS_STEP_INCOMPLETE; 761f68a32c8SEmil Constantinescu ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 762f68a32c8SEmil Constantinescu rk->status = TS_STEP_PENDING; 763f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 764f68a32c8SEmil Constantinescu ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 7651917a363SLisandro Dalcin ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr); 766fecfb714SLisandro Dalcin ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 767be5899b3SLisandro Dalcin rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 768be5899b3SLisandro Dalcin if (!accept) { /* Roll back the current step */ 769fecfb714SLisandro Dalcin ierr = TSRollBack_RK(ts);CHKERRQ(ierr); 770be5899b3SLisandro Dalcin ts->time_step = next_time_step; 771be5899b3SLisandro Dalcin goto reject_step; 772be5899b3SLisandro Dalcin } 773be5899b3SLisandro Dalcin 7740f7a1166SHong Zhang if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 7750f7a1166SHong Zhang rk->ptime = ts->ptime; 7760f7a1166SHong Zhang rk->time_step = ts->time_step; 777493ed95fSHong Zhang } 778be5899b3SLisandro Dalcin 779f68a32c8SEmil Constantinescu ts->ptime += ts->time_step; 780f68a32c8SEmil Constantinescu ts->time_step = next_time_step; 781f68a32c8SEmil Constantinescu break; 782be5899b3SLisandro Dalcin 783be5899b3SLisandro Dalcin reject_step: 784fecfb714SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 785be5899b3SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 786be5899b3SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 787be5899b3SLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 788f68a32c8SEmil Constantinescu } 789f68a32c8SEmil Constantinescu } 790f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 791e4dd521cSBarry Smith } 792e4dd521cSBarry Smith 793f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts) 794f6a906c0SBarry Smith { 795f6a906c0SBarry Smith TS_RK *rk = (TS_RK*)ts->data; 796be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 797be5899b3SLisandro Dalcin PetscInt s = tab->s; 798f6a906c0SBarry Smith PetscErrorCode ierr; 799f6a906c0SBarry Smith 800f6a906c0SBarry Smith PetscFunctionBegin; 801f6a906c0SBarry Smith if (ts->adjointsetupcalled++) PetscFunctionReturn(0); 8022e7b7f96SHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecsDeltaLam);CHKERRQ(ierr); 8032e7b7f96SHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&rk->VecsSensiTemp);CHKERRQ(ierr); 804f6a906c0SBarry Smith if(ts->vecs_sensip) { 8052e7b7f96SHong Zhang ierr = VecDuplicate(ts->vecs_sensip[0],&rk->VecDeltaMu);CHKERRQ(ierr); 806f6a906c0SBarry Smith } 80713af1a74SHong Zhang if (ts->vecs_sensi2) { 80813af1a74SHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecsDeltaLam2);CHKERRQ(ierr); 80913af1a74SHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&rk->VecsSensi2Temp);CHKERRQ(ierr); 81013af1a74SHong Zhang } 81113af1a74SHong Zhang if (ts->vecs_sensi2p) { 81213af1a74SHong Zhang ierr = VecDuplicate(ts->vecs_sensi2p[0],&rk->VecDeltaMu2);CHKERRQ(ierr); 81313af1a74SHong Zhang } 814f6a906c0SBarry Smith PetscFunctionReturn(0); 815f6a906c0SBarry Smith } 816f6a906c0SBarry Smith 81735f5172aSHong Zhang /* 81835f5172aSHong Zhang Assumptions: 81935f5172aSHong Zhang - TSStep_RK() always evaluates the step with b, not bembed. 82035f5172aSHong Zhang */ 82142f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts) 822d2daff3dSHong Zhang { 823c235aa8dSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 824cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 825c235aa8dSHong Zhang RKTableau tab = rk->tableau; 826cd4cee2dSHong Zhang Mat J,Jquad; 827c235aa8dSHong Zhang const PetscInt s = tab->s; 828c235aa8dSHong Zhang const PetscReal *A = tab->A,*b = tab->b,*c = tab->c; 82913af1a74SHong Zhang PetscScalar *w = rk->work,*xarr; 8302e7b7f96SHong Zhang Vec *Y = rk->Y,*VecsDeltaLam = rk->VecsDeltaLam,VecDeltaMu = rk->VecDeltaMu,*VecsSensiTemp = rk->VecsSensiTemp; 83113af1a74SHong Zhang Vec *VecsDeltaLam2 = rk->VecsDeltaLam2,VecDeltaMu2 = rk->VecDeltaMu2,*VecsSensi2Temp = rk->VecsSensi2Temp; 832cd4cee2dSHong Zhang Vec VecDRDUTransCol = ts->vec_drdu_col,VecDRDPTransCol = ts->vec_drdp_col; 833b18ea86cSHong Zhang PetscInt i,j,nadj; 8343d3eaba7SBarry Smith PetscReal t = ts->ptime; 8352e7b7f96SHong Zhang PetscReal h = ts->time_step,stage_time; 836d2daff3dSHong Zhang PetscErrorCode ierr; 837c235aa8dSHong Zhang 838d2daff3dSHong Zhang PetscFunctionBegin; 839c235aa8dSHong Zhang rk->status = TS_STEP_INCOMPLETE; 8403389c7d9SStefano Zampini 841cd4cee2dSHong Zhang ierr = TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr); 84235f5172aSHong Zhang if (quadts) { 843cd4cee2dSHong Zhang ierr = TSGetRHSJacobian(quadts,&Jquad,NULL,NULL,NULL);CHKERRQ(ierr); 84435f5172aSHong Zhang } 8452e7b7f96SHong Zhang for (i=s-1; i>=0; i--) { 8466a1e4597SHong Zhang if (tab->FSAL && i == s-1) { 8476a1e4597SHong Zhang /* VecsDeltaLam[nadj*s+s-1] are initialized with zeros and the values never change.*/ 8486a1e4597SHong Zhang continue; 8496a1e4597SHong Zhang } 8502e7b7f96SHong Zhang stage_time = t + h*(1.0-c[i]); 8513389c7d9SStefano Zampini ierr = TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);CHKERRQ(ierr); 85235f5172aSHong Zhang if (quadts) { 853cd4cee2dSHong Zhang ierr = TSComputeRHSJacobian(quadts,stage_time,Y[i],Jquad,Jquad);CHKERRQ(ierr); /* get r_u^T */ 85435f5172aSHong Zhang } 8553389c7d9SStefano Zampini if (ts->vecs_sensip) { 85613af1a74SHong Zhang ierr = TSComputeRHSJacobianP(ts,stage_time,Y[i],ts->Jacprhs);CHKERRQ(ierr); /* get f_p */ 85735f5172aSHong Zhang if (quadts) { 858cd4cee2dSHong Zhang ierr = TSComputeRHSJacobianP(quadts,stage_time,Y[i],quadts->Jacprhs);CHKERRQ(ierr); /* get f_p for the quadrature */ 8593389c7d9SStefano Zampini } 86035f5172aSHong Zhang } 8613389c7d9SStefano Zampini 86235f5172aSHong Zhang if (b[i]) { 86335f5172aSHong Zhang for (j=i+1; j<s; j++) w[j-i-1] = A[j*s+i]/b[i]; /* coefficients for computing VecsSensiTemp */ 86435f5172aSHong Zhang } else { 86535f5172aSHong Zhang for (j=i+1; j<s; j++) w[j-i-1] = A[j*s+i]; /* coefficients for computing VecsSensiTemp */ 86635f5172aSHong Zhang } 8672e7b7f96SHong Zhang 8682e7b7f96SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 8693389c7d9SStefano Zampini /* Stage values of lambda */ 87035f5172aSHong Zhang if (b[i]) { 87135f5172aSHong Zhang /* lambda_{n+1} + \sum_{j=i+1}^s a_{ji}/b[i]*lambda_{s,j} */ 8722e7b7f96SHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); /* VecDeltaLam is an vec array of size s by numcost */ 8732e7b7f96SHong Zhang ierr = VecMAXPY(VecsSensiTemp[nadj],s-i-1,w,&VecsDeltaLam[nadj*s+i+1]);CHKERRQ(ierr); 87435f5172aSHong Zhang ierr = MatMultTranspose(J,VecsSensiTemp[nadj],VecsDeltaLam[nadj*s+i]);CHKERRQ(ierr); /* VecsSensiTemp will be reused by 2nd-order adjoint */ 87535f5172aSHong Zhang ierr = VecScale(VecsDeltaLam[nadj*s+i],-h*b[i]);CHKERRQ(ierr); 87635f5172aSHong Zhang if (quadts) { 877cd4cee2dSHong Zhang ierr = MatDenseGetColumn(Jquad,nadj,&xarr);CHKERRQ(ierr); 878cd4cee2dSHong Zhang ierr = VecPlaceArray(VecDRDUTransCol,xarr);CHKERRQ(ierr); 879cd4cee2dSHong Zhang ierr = VecAXPY(VecsDeltaLam[nadj*s+i],-h*b[i],VecDRDUTransCol);CHKERRQ(ierr); 880cd4cee2dSHong Zhang ierr = VecResetArray(VecDRDUTransCol);CHKERRQ(ierr); 881cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(Jquad,&xarr);CHKERRQ(ierr); 882cd4cee2dSHong Zhang } 8833389c7d9SStefano Zampini } else { 8846a1e4597SHong Zhang /* \sum_{j=i+1}^s a_{ji}*lambda_{s,j} */ 8856a1e4597SHong Zhang ierr = VecSet(VecsSensiTemp[nadj],0);CHKERRQ(ierr); 8866a1e4597SHong Zhang ierr = VecMAXPY(VecsSensiTemp[nadj],s-i-1,w,&VecsDeltaLam[nadj*s+i+1]);CHKERRQ(ierr); 8876a1e4597SHong Zhang ierr = MatMultTranspose(J,VecsSensiTemp[nadj],VecsDeltaLam[nadj*s+i]);CHKERRQ(ierr); 8886a1e4597SHong Zhang ierr = VecScale(VecsDeltaLam[nadj*s+i],-h);CHKERRQ(ierr); 8893389c7d9SStefano Zampini } 8906497ce14SHong Zhang 891ad8e2604SHong Zhang /* Stage values of mu */ 8926497ce14SHong Zhang if (ts->vecs_sensip) { 89313af1a74SHong Zhang ierr = MatMultTranspose(ts->Jacprhs,VecsSensiTemp[nadj],VecDeltaMu);CHKERRQ(ierr); 89435f5172aSHong Zhang if (b[i]) { 89535f5172aSHong Zhang ierr = VecScale(VecDeltaMu,-h*b[i]);CHKERRQ(ierr); 89635f5172aSHong Zhang if (quadts) { 897cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadts->Jacprhs,nadj,&xarr);CHKERRQ(ierr); 898cd4cee2dSHong Zhang ierr = VecPlaceArray(VecDRDPTransCol,xarr);CHKERRQ(ierr); 899cd4cee2dSHong Zhang ierr = VecAXPY(VecDeltaMu,-h*b[i],VecDRDPTransCol);CHKERRQ(ierr); 900cd4cee2dSHong Zhang ierr = VecResetArray(VecDRDPTransCol);CHKERRQ(ierr); 901cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadts->Jacprhs,&xarr);CHKERRQ(ierr); 902493ed95fSHong Zhang } 90335f5172aSHong Zhang } else { 90435f5172aSHong Zhang ierr = VecScale(VecDeltaMu,-h);CHKERRQ(ierr); 90535f5172aSHong Zhang } 9062e7b7f96SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],1.,VecDeltaMu);CHKERRQ(ierr); /* update sensip for each stage */ 907ad8e2604SHong Zhang } 908c235aa8dSHong Zhang } 90913af1a74SHong Zhang 91013af1a74SHong Zhang if (ts->vecs_sensi2 && ts->forward_solve) { /* 2nd-order adjoint, TLM mode has to be turned on */ 91113af1a74SHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 91213af1a74SHong Zhang ierr = MatDenseGetColumn(rk->MatsFwdStageSensip[i],0,&xarr);CHKERRQ(ierr); 91313af1a74SHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 91413af1a74SHong Zhang /* lambda_s^T F_UU w_1 */ 915b1e111ebSHong Zhang ierr = TSComputeRHSHessianProductFunctionUU(ts,stage_time,Y[i],VecsSensiTemp,ts->vec_sensip_col,ts->vecs_guu);CHKERRQ(ierr); 91635f5172aSHong Zhang if (quadts) { 91735f5172aSHong Zhang /* R_UU w_1 */ 918b1e111ebSHong Zhang ierr = TSComputeRHSHessianProductFunctionUU(quadts,stage_time,Y[i],NULL,ts->vec_sensip_col,ts->vecs_guu);CHKERRQ(ierr); 91935f5172aSHong Zhang } 92035f5172aSHong Zhang if (ts->vecs_sensip) { 92113af1a74SHong Zhang /* lambda_s^T F_UP w_2 */ 922b1e111ebSHong Zhang ierr = TSComputeRHSHessianProductFunctionUP(ts,stage_time,Y[i],VecsSensiTemp,ts->vec_dir,ts->vecs_gup);CHKERRQ(ierr); 92335f5172aSHong Zhang if (quadts) { 92435f5172aSHong Zhang /* R_UP w_2 */ 925b1e111ebSHong Zhang ierr = TSComputeRHSHessianProductFunctionUP(quadts,stage_time,Y[i],NULL,ts->vec_sensip_col,ts->vecs_gup);CHKERRQ(ierr); 92635f5172aSHong Zhang } 92735f5172aSHong Zhang } 92813af1a74SHong Zhang if (ts->vecs_sensi2p) { 92913af1a74SHong Zhang /* lambda_s^T F_PU w_1 */ 930b1e111ebSHong Zhang ierr = TSComputeRHSHessianProductFunctionPU(ts,stage_time,Y[i],VecsSensiTemp,ts->vec_sensip_col,ts->vecs_gpu);CHKERRQ(ierr); 93135f5172aSHong Zhang /* lambda_s^T F_PP w_2 */ 932b1e111ebSHong Zhang ierr = TSComputeRHSHessianProductFunctionPP(ts,stage_time,Y[i],VecsSensiTemp,ts->vec_dir,ts->vecs_gpp);CHKERRQ(ierr); 93335f5172aSHong Zhang if (b[i] && quadts) { 93435f5172aSHong Zhang /* R_PU w_1 */ 935b1e111ebSHong Zhang ierr = TSComputeRHSHessianProductFunctionPU(quadts,stage_time,Y[i],NULL,ts->vec_sensip_col,ts->vecs_gpu);CHKERRQ(ierr); 93635f5172aSHong Zhang /* R_PP w_2 */ 937b1e111ebSHong Zhang ierr = TSComputeRHSHessianProductFunctionPP(quadts,stage_time,Y[i],NULL,ts->vec_dir,ts->vecs_gpp);CHKERRQ(ierr); 93835f5172aSHong Zhang } 93913af1a74SHong Zhang } 94013af1a74SHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 94113af1a74SHong Zhang ierr = MatDenseRestoreColumn(rk->MatsFwdStageSensip[i],&xarr);CHKERRQ(ierr); 94213af1a74SHong Zhang 94313af1a74SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 94413af1a74SHong Zhang /* Stage values of lambda */ 94535f5172aSHong Zhang if (b[i]) { 94635f5172aSHong Zhang /* J_i^T*(Lambda_{n+1}+\sum_{j=i+1}^s a_{ji}/b_i*Lambda_{s,j} */ 94713af1a74SHong Zhang ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 94813af1a74SHong Zhang ierr = VecMAXPY(VecsSensi2Temp[nadj],s-i-1,w,&VecsDeltaLam2[nadj*s+i+1]);CHKERRQ(ierr); 94913af1a74SHong Zhang ierr = MatMultTranspose(J,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj*s+i]);CHKERRQ(ierr); 95035f5172aSHong Zhang ierr = VecScale(VecsDeltaLam2[nadj*s+i],-h*b[i]);CHKERRQ(ierr); 95135f5172aSHong Zhang ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],-h*b[i],ts->vecs_guu[nadj]);CHKERRQ(ierr); 95235f5172aSHong Zhang if (ts->vecs_sensip) { 95335f5172aSHong Zhang ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],-h*b[i],ts->vecs_gup[nadj]);CHKERRQ(ierr); 95413af1a74SHong Zhang } 95513af1a74SHong Zhang } else { 95635f5172aSHong Zhang /* \sum_{j=i+1}^s a_{ji}*Lambda_{s,j} */ 95713af1a74SHong Zhang ierr = VecSet(VecsDeltaLam2[nadj*s+i],0);CHKERRQ(ierr); 95835f5172aSHong Zhang ierr = VecMAXPY(VecsSensi2Temp[nadj],s-i-1,w,&VecsDeltaLam2[nadj*s+i+1]);CHKERRQ(ierr); 95935f5172aSHong Zhang ierr = MatMultTranspose(J,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj*s+i]);CHKERRQ(ierr); 96035f5172aSHong Zhang ierr = VecScale(VecsDeltaLam2[nadj*s+i],-h);CHKERRQ(ierr); 96135f5172aSHong Zhang ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],-h,ts->vecs_guu[nadj]);CHKERRQ(ierr); 96235f5172aSHong Zhang if (ts->vecs_sensip) { 96335f5172aSHong Zhang ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],-h,ts->vecs_gup[nadj]);CHKERRQ(ierr); 96413af1a74SHong Zhang } 96535f5172aSHong Zhang } 96635f5172aSHong Zhang if (ts->vecs_sensi2p) { /* 2nd-order adjoint for parameters */ 96713af1a74SHong Zhang ierr = MatMultTranspose(ts->Jacprhs,VecsSensi2Temp[nadj],VecDeltaMu2);CHKERRQ(ierr); 96835f5172aSHong Zhang if (b[i]) { 96935f5172aSHong Zhang ierr = VecScale(VecDeltaMu2,-h*b[i]);CHKERRQ(ierr); 97035f5172aSHong Zhang ierr = VecAXPY(VecDeltaMu2,-h*b[i],ts->vecs_gpu[nadj]);CHKERRQ(ierr); 97135f5172aSHong Zhang ierr = VecAXPY(VecDeltaMu2,-h*b[i],ts->vecs_gpp[nadj]);CHKERRQ(ierr); 97213af1a74SHong Zhang } else { 97335f5172aSHong Zhang ierr = VecScale(VecDeltaMu2,-h);CHKERRQ(ierr); 97435f5172aSHong Zhang ierr = VecAXPY(VecDeltaMu2,-h,ts->vecs_gpu[nadj]);CHKERRQ(ierr); 97535f5172aSHong Zhang ierr = VecAXPY(VecDeltaMu2,-h,ts->vecs_gpp[nadj]);CHKERRQ(ierr); 97613af1a74SHong Zhang } 97713af1a74SHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],1,VecDeltaMu2);CHKERRQ(ierr); /* update sensi2p for each stage */ 97813af1a74SHong Zhang } 97913af1a74SHong Zhang } 98013af1a74SHong Zhang } 9816497ce14SHong Zhang } 982c235aa8dSHong Zhang 983c235aa8dSHong Zhang for (j=0; j<s; j++) w[j] = 1.0; 9842e7b7f96SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { /* no need to do this for mu's */ 9852e7b7f96SHong Zhang ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecsDeltaLam[nadj*s]);CHKERRQ(ierr); 98613af1a74SHong Zhang if (ts->vecs_sensi2) { 98713af1a74SHong Zhang ierr = VecMAXPY(ts->vecs_sensi2[nadj],s,w,&VecsDeltaLam2[nadj*s]);CHKERRQ(ierr); 98813af1a74SHong Zhang } 9896497ce14SHong Zhang } 990c235aa8dSHong Zhang rk->status = TS_STEP_COMPLETE; 991d2daff3dSHong Zhang PetscFunctionReturn(0); 992d2daff3dSHong Zhang } 993d2daff3dSHong Zhang 99413af1a74SHong Zhang static PetscErrorCode TSAdjointReset_RK(TS ts) 99513af1a74SHong Zhang { 99613af1a74SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 99713af1a74SHong Zhang RKTableau tab = rk->tableau; 99813af1a74SHong Zhang PetscErrorCode ierr; 99913af1a74SHong Zhang 100013af1a74SHong Zhang PetscFunctionBegin; 100113af1a74SHong Zhang ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam);CHKERRQ(ierr); 100213af1a74SHong Zhang ierr = VecDestroyVecs(ts->numcost,&rk->VecsSensiTemp);CHKERRQ(ierr); 100313af1a74SHong Zhang ierr = VecDestroy(&rk->VecDeltaMu);CHKERRQ(ierr); 100413af1a74SHong Zhang ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam2);CHKERRQ(ierr); 100513af1a74SHong Zhang ierr = VecDestroy(&rk->VecDeltaMu2);CHKERRQ(ierr); 100613af1a74SHong Zhang ierr = VecDestroyVecs(ts->numcost,&rk->VecsSensi2Temp);CHKERRQ(ierr); 100713af1a74SHong Zhang PetscFunctionReturn(0); 100813af1a74SHong Zhang } 100913af1a74SHong Zhang 101055de54fcSHong Zhang static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X) 101155de54fcSHong Zhang { 101255de54fcSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 101355de54fcSHong Zhang PetscInt s = rk->tableau->s,p = rk->tableau->p,i,j; 101455de54fcSHong Zhang PetscReal h; 101555de54fcSHong Zhang PetscReal tt,t; 101655de54fcSHong Zhang PetscScalar *b; 101755de54fcSHong Zhang const PetscReal *B = rk->tableau->binterp; 101855de54fcSHong Zhang PetscErrorCode ierr; 101955de54fcSHong Zhang 102055de54fcSHong Zhang PetscFunctionBegin; 102155de54fcSHong Zhang if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 102255de54fcSHong Zhang 102355de54fcSHong Zhang switch (rk->status) { 102455de54fcSHong Zhang case TS_STEP_INCOMPLETE: 102555de54fcSHong Zhang case TS_STEP_PENDING: 102655de54fcSHong Zhang h = ts->time_step; 102755de54fcSHong Zhang t = (itime - ts->ptime)/h; 102855de54fcSHong Zhang break; 102955de54fcSHong Zhang case TS_STEP_COMPLETE: 103055de54fcSHong Zhang h = ts->ptime - ts->ptime_prev; 103155de54fcSHong Zhang t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 103255de54fcSHong Zhang break; 103355de54fcSHong Zhang default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 103455de54fcSHong Zhang } 103555de54fcSHong Zhang ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 103655de54fcSHong Zhang for (i=0; i<s; i++) b[i] = 0; 103755de54fcSHong Zhang for (j=0,tt=t; j<p; j++,tt*=t) { 103855de54fcSHong Zhang for (i=0; i<s; i++) { 103955de54fcSHong Zhang b[i] += h * B[i*p+j] * tt; 104055de54fcSHong Zhang } 104155de54fcSHong Zhang } 104255de54fcSHong Zhang ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr); 104355de54fcSHong Zhang ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr); 104455de54fcSHong Zhang ierr = PetscFree(b);CHKERRQ(ierr); 104555de54fcSHong Zhang PetscFunctionReturn(0); 104655de54fcSHong Zhang } 104755de54fcSHong Zhang 104855de54fcSHong Zhang /*------------------------------------------------------------*/ 104955de54fcSHong Zhang 1050be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts) 1051be5899b3SLisandro Dalcin { 1052be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 1053be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 1054be5899b3SLisandro Dalcin PetscErrorCode ierr; 1055be5899b3SLisandro Dalcin 1056be5899b3SLisandro Dalcin PetscFunctionBegin; 1057be5899b3SLisandro Dalcin if (!tab) PetscFunctionReturn(0); 1058be5899b3SLisandro Dalcin ierr = PetscFree(rk->work);CHKERRQ(ierr); 1059be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr); 1060be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr); 10612e7b7f96SHong Zhang ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam);CHKERRQ(ierr); 10622e7b7f96SHong Zhang ierr = VecDestroyVecs(ts->numcost,&rk->VecsSensiTemp);CHKERRQ(ierr); 10632e7b7f96SHong Zhang ierr = VecDestroy(&rk->VecDeltaMu);CHKERRQ(ierr); 1064be5899b3SLisandro Dalcin PetscFunctionReturn(0); 1065be5899b3SLisandro Dalcin } 1066be5899b3SLisandro Dalcin 1067277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts) 1068e4dd521cSBarry Smith { 10696849ba73SBarry Smith PetscErrorCode ierr; 1070e4dd521cSBarry Smith 1071e4dd521cSBarry Smith PetscFunctionBegin; 1072be5899b3SLisandro Dalcin ierr = TSRKTableauReset(ts);CHKERRQ(ierr); 10730fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 107402555c94SHong Zhang ierr = PetscTryMethod(ts,"TSReset_RK_MultirateSplit_C",(TS),(ts));CHKERRQ(ierr); 10750fe4d17eSHong Zhang } else { 107602555c94SHong Zhang ierr = PetscTryMethod(ts,"TSReset_RK_MultirateNonsplit_C",(TS),(ts));CHKERRQ(ierr); 10770fe4d17eSHong Zhang } 1078277b19d0SLisandro Dalcin PetscFunctionReturn(0); 1079e4dd521cSBarry Smith } 1080277b19d0SLisandro Dalcin 1081f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx) 1082f68a32c8SEmil Constantinescu { 1083f68a32c8SEmil Constantinescu PetscFunctionBegin; 1084f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1085f68a32c8SEmil Constantinescu } 1086f68a32c8SEmil Constantinescu 1087f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 1088f68a32c8SEmil Constantinescu { 1089f68a32c8SEmil Constantinescu PetscFunctionBegin; 1090f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1091f68a32c8SEmil Constantinescu } 1092f68a32c8SEmil Constantinescu 1093f68a32c8SEmil Constantinescu 1094f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx) 1095f68a32c8SEmil Constantinescu { 1096f68a32c8SEmil Constantinescu PetscFunctionBegin; 1097f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1098f68a32c8SEmil Constantinescu } 1099f68a32c8SEmil Constantinescu 1100f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 1101f68a32c8SEmil Constantinescu { 1102f68a32c8SEmil Constantinescu 1103f68a32c8SEmil Constantinescu PetscFunctionBegin; 1104f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1105f68a32c8SEmil Constantinescu } 1106c235aa8dSHong Zhang /* 1107d2daff3dSHong Zhang static PetscErrorCode RKSetAdjCoe(RKTableau tab) 1108d2daff3dSHong Zhang { 1109d2daff3dSHong Zhang PetscReal *A,*b; 1110d2daff3dSHong Zhang PetscInt s,i,j; 1111d2daff3dSHong Zhang PetscErrorCode ierr; 1112d2daff3dSHong Zhang 1113d2daff3dSHong Zhang PetscFunctionBegin; 1114d2daff3dSHong Zhang s = tab->s; 1115d2daff3dSHong Zhang ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr); 1116d2daff3dSHong Zhang 1117d2daff3dSHong Zhang for (i=0; i<s; i++) 1118d2daff3dSHong Zhang for (j=0; j<s; j++) { 1119d2daff3dSHong Zhang A[i*s+j] = (tab->b[s-1-i]==0) ? 0: -tab->A[s-1-i+(s-1-j)*s] * tab->b[s-1-j] / tab->b[s-1-i]; 1120d2daff3dSHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr); 1121d2daff3dSHong Zhang } 1122d2daff3dSHong Zhang 1123d2daff3dSHong Zhang for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i]; 1124d2daff3dSHong Zhang 1125d2daff3dSHong Zhang ierr = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 1126d2daff3dSHong Zhang ierr = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr); 1127d2daff3dSHong Zhang ierr = PetscFree2(A,b);CHKERRQ(ierr); 1128d2daff3dSHong Zhang PetscFunctionReturn(0); 1129d2daff3dSHong Zhang } 1130c235aa8dSHong Zhang */ 1131be5899b3SLisandro Dalcin 1132be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts) 1133be5899b3SLisandro Dalcin { 1134be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 1135be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 1136be5899b3SLisandro Dalcin PetscErrorCode ierr; 1137be5899b3SLisandro Dalcin 1138be5899b3SLisandro Dalcin PetscFunctionBegin; 1139be5899b3SLisandro Dalcin ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr); 1140be5899b3SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr); 1141be5899b3SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr); 1142be5899b3SLisandro Dalcin PetscFunctionReturn(0); 1143be5899b3SLisandro Dalcin } 1144be5899b3SLisandro Dalcin 1145f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts) 1146f68a32c8SEmil Constantinescu { 1147cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 1148f68a32c8SEmil Constantinescu PetscErrorCode ierr; 1149f68a32c8SEmil Constantinescu DM dm; 1150f68a32c8SEmil Constantinescu 1151f68a32c8SEmil Constantinescu PetscFunctionBegin; 11523dd83b38SBarry Smith ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr); 1153be5899b3SLisandro Dalcin ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr); 1154cd4cee2dSHong Zhang if (quadts && ts->costintegralfwd) { 1155cd4cee2dSHong Zhang Mat Jquad; 1156cd4cee2dSHong Zhang ierr = TSGetRHSJacobian(quadts,&Jquad,NULL,NULL,NULL);CHKERRQ(ierr); 11570f7a1166SHong Zhang } 1158f68a32c8SEmil Constantinescu ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 1159f68a32c8SEmil Constantinescu ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 1160f68a32c8SEmil Constantinescu ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 11610fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 116202555c94SHong Zhang ierr = PetscTryMethod(ts,"TSSetUp_RK_MultirateSplit_C",(TS),(ts));CHKERRQ(ierr); 11630fe4d17eSHong Zhang } else { 116402555c94SHong Zhang ierr = PetscTryMethod(ts,"TSSetUp_RK_MultirateNonsplit_C",(TS),(ts));CHKERRQ(ierr); 11650fe4d17eSHong Zhang } 1166e4dd521cSBarry Smith PetscFunctionReturn(0); 1167e4dd521cSBarry Smith } 1168d2daff3dSHong Zhang 11694416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts) 1170e4dd521cSBarry Smith { 1171be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 1172dfbe8321SBarry Smith PetscErrorCode ierr; 1173262ff7bbSBarry Smith 1174e4dd521cSBarry Smith PetscFunctionBegin; 1175e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr); 1176f68a32c8SEmil Constantinescu { 1177f68a32c8SEmil Constantinescu RKTableauLink link; 1178f68a32c8SEmil Constantinescu PetscInt count,choice; 11790fe4d17eSHong Zhang PetscBool flg,use_multirate = PETSC_FALSE; 1180f68a32c8SEmil Constantinescu const char **namelist; 11812c9cb062Sluzhanghpp 1182f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) ; 1183f489ac74SBarry Smith ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr); 1184f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 11850fe4d17eSHong Zhang ierr = PetscOptionsBool("-ts_rk_multirate","Use interpolation-based multirate RK method","TSRKSetMultirate",rk->use_multirate,&use_multirate,&flg);CHKERRQ(ierr); 11860fe4d17eSHong Zhang if (flg) { 11870fe4d17eSHong Zhang ierr = TSRKSetMultirate(ts,use_multirate);CHKERRQ(ierr); 11880fe4d17eSHong Zhang } 1189be5899b3SLisandro Dalcin ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr); 1190be5899b3SLisandro Dalcin if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);} 1191f68a32c8SEmil Constantinescu ierr = PetscFree(namelist);CHKERRQ(ierr); 1192f68a32c8SEmil Constantinescu } 1193262ff7bbSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 11942c9cb062Sluzhanghpp ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Multirate methods options","");CHKERRQ(ierr); 11952c9cb062Sluzhanghpp ierr = PetscOptionsInt("-ts_rk_dtratio","time step ratio between slow and fast","",rk->dtratio,&rk->dtratio,NULL);CHKERRQ(ierr); 11962c9cb062Sluzhanghpp ierr = PetscOptionsEnd();CHKERRQ(ierr); 1197e4dd521cSBarry Smith PetscFunctionReturn(0); 1198e4dd521cSBarry Smith } 1199e4dd521cSBarry Smith 12005f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer) 1201e4dd521cSBarry Smith { 12025f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 12038370ee3bSLisandro Dalcin PetscBool iascii; 1204dfbe8321SBarry Smith PetscErrorCode ierr; 12052cdf8120Sasbjorn 1206e4dd521cSBarry Smith PetscFunctionBegin; 1207251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 12088370ee3bSLisandro Dalcin if (iascii) { 12099c334d8fSLisandro Dalcin RKTableau tab = rk->tableau; 1210f68a32c8SEmil Constantinescu TSRKType rktype; 1211f68a32c8SEmil Constantinescu char buf[512]; 1212f68a32c8SEmil Constantinescu ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr); 1213efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RK type %s\n",rktype);CHKERRQ(ierr); 12140c37eb8eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," Order: %D\n",tab->order);CHKERRQ(ierr); 12150c37eb8eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr); 1216f68a32c8SEmil Constantinescu ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 1217f68a32c8SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," Abscissa c = %s\n",buf);CHKERRQ(ierr); 12188370ee3bSLisandro Dalcin } 1219f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1220f68a32c8SEmil Constantinescu } 1221f68a32c8SEmil Constantinescu 1222f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer) 1223f68a32c8SEmil Constantinescu { 1224f68a32c8SEmil Constantinescu PetscErrorCode ierr; 12259c334d8fSLisandro Dalcin TSAdapt adapt; 1226f68a32c8SEmil Constantinescu 1227f68a32c8SEmil Constantinescu PetscFunctionBegin; 12289c334d8fSLisandro Dalcin ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 12299c334d8fSLisandro Dalcin ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr); 1230f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1231f68a32c8SEmil Constantinescu } 1232f68a32c8SEmil Constantinescu 1233f68a32c8SEmil Constantinescu /*@C 1234f68a32c8SEmil Constantinescu TSRKSetType - Set the type of RK scheme 1235f68a32c8SEmil Constantinescu 1236f68a32c8SEmil Constantinescu Logically collective 1237f68a32c8SEmil Constantinescu 1238f68a32c8SEmil Constantinescu Input Parameter: 1239f68a32c8SEmil Constantinescu + ts - timestepping context 1240f68a32c8SEmil Constantinescu - rktype - type of RK-scheme 1241f68a32c8SEmil Constantinescu 1242287c2655SBarry Smith Options Database: 12439bd3e852SBarry Smith . -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs> 1244287c2655SBarry Smith 1245f68a32c8SEmil Constantinescu Level: intermediate 1246f68a32c8SEmil Constantinescu 124737acaa02SHendrik Ranocha .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS, TSRK6VR, TSRK7VR, TSRK8VR 1248f68a32c8SEmil Constantinescu @*/ 1249f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype) 1250f68a32c8SEmil Constantinescu { 1251f68a32c8SEmil Constantinescu PetscErrorCode ierr; 1252f68a32c8SEmil Constantinescu 1253f68a32c8SEmil Constantinescu PetscFunctionBegin; 1254f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1255b92453a8SLisandro Dalcin PetscValidCharPointer(rktype,2); 1256f68a32c8SEmil Constantinescu ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr); 1257f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1258f68a32c8SEmil Constantinescu } 1259f68a32c8SEmil Constantinescu 1260f68a32c8SEmil Constantinescu /*@C 1261f68a32c8SEmil Constantinescu TSRKGetType - Get the type of RK scheme 1262f68a32c8SEmil Constantinescu 1263f68a32c8SEmil Constantinescu Logically collective 1264f68a32c8SEmil Constantinescu 1265f68a32c8SEmil Constantinescu Input Parameter: 1266f68a32c8SEmil Constantinescu . ts - timestepping context 1267f68a32c8SEmil Constantinescu 1268f68a32c8SEmil Constantinescu Output Parameter: 1269f68a32c8SEmil Constantinescu . rktype - type of RK-scheme 1270f68a32c8SEmil Constantinescu 1271f68a32c8SEmil Constantinescu Level: intermediate 1272f68a32c8SEmil Constantinescu 1273f68a32c8SEmil Constantinescu .seealso: TSRKGetType() 1274f68a32c8SEmil Constantinescu @*/ 1275f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype) 1276f68a32c8SEmil Constantinescu { 1277f68a32c8SEmil Constantinescu PetscErrorCode ierr; 1278f68a32c8SEmil Constantinescu 1279f68a32c8SEmil Constantinescu PetscFunctionBegin; 1280f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1281f68a32c8SEmil Constantinescu ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr); 1282f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1283f68a32c8SEmil Constantinescu } 1284f68a32c8SEmil Constantinescu 1285560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype) 1286f68a32c8SEmil Constantinescu { 1287f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 1288f68a32c8SEmil Constantinescu 1289f68a32c8SEmil Constantinescu PetscFunctionBegin; 1290f68a32c8SEmil Constantinescu *rktype = rk->tableau->name; 1291f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1292f68a32c8SEmil Constantinescu } 12932c9cb062Sluzhanghpp 1294560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype) 1295f68a32c8SEmil Constantinescu { 1296f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 1297f68a32c8SEmil Constantinescu PetscErrorCode ierr; 1298f68a32c8SEmil Constantinescu PetscBool match; 1299f68a32c8SEmil Constantinescu RKTableauLink link; 1300f68a32c8SEmil Constantinescu 1301f68a32c8SEmil Constantinescu PetscFunctionBegin; 1302f68a32c8SEmil Constantinescu if (rk->tableau) { 1303f68a32c8SEmil Constantinescu ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr); 1304f68a32c8SEmil Constantinescu if (match) PetscFunctionReturn(0); 1305f68a32c8SEmil Constantinescu } 1306f68a32c8SEmil Constantinescu for (link = RKTableauList; link; link=link->next) { 1307f68a32c8SEmil Constantinescu ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr); 1308f68a32c8SEmil Constantinescu if (match) { 1309be5899b3SLisandro Dalcin if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);} 1310f68a32c8SEmil Constantinescu rk->tableau = &link->tab; 1311be5899b3SLisandro Dalcin if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);} 13122ffb9264SLisandro Dalcin ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 1313f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1314f68a32c8SEmil Constantinescu } 1315f68a32c8SEmil Constantinescu } 1316f68a32c8SEmil Constantinescu SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype); 1317e4dd521cSBarry Smith PetscFunctionReturn(0); 1318e4dd521cSBarry Smith } 1319e4dd521cSBarry Smith 1320ff22ae23SHong Zhang static PetscErrorCode TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y) 1321ff22ae23SHong Zhang { 1322ff22ae23SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 1323ff22ae23SHong Zhang 1324ff22ae23SHong Zhang PetscFunctionBegin; 13250429704eSStefano Zampini if (ns) *ns = rk->tableau->s; 1326d2daff3dSHong Zhang if (Y) *Y = rk->Y; 1327ff22ae23SHong Zhang PetscFunctionReturn(0); 1328ff22ae23SHong Zhang } 1329ff22ae23SHong Zhang 1330b3a6b972SJed Brown static PetscErrorCode TSDestroy_RK(TS ts) 1331b3a6b972SJed Brown { 1332b3a6b972SJed Brown PetscErrorCode ierr; 1333b3a6b972SJed Brown 1334b3a6b972SJed Brown PetscFunctionBegin; 1335b3a6b972SJed Brown ierr = TSReset_RK(ts);CHKERRQ(ierr); 1336b3a6b972SJed Brown if (ts->dm) { 1337b3a6b972SJed Brown ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 1338b3a6b972SJed Brown ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 1339b3a6b972SJed Brown } 1340b3a6b972SJed Brown ierr = PetscFree(ts->data);CHKERRQ(ierr); 1341b3a6b972SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr); 1342b3a6b972SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr); 13430fe4d17eSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",NULL);CHKERRQ(ierr); 13440fe4d17eSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",NULL);CHKERRQ(ierr); 1345b3a6b972SJed Brown PetscFunctionReturn(0); 1346b3a6b972SJed Brown } 1347ff22ae23SHong Zhang 134821052c0cSHong Zhang /*@C 134921052c0cSHong Zhang TSRKSetMultirate - Use the interpolation-based multirate RK method 135021052c0cSHong Zhang 135121052c0cSHong Zhang Logically collective 135221052c0cSHong Zhang 135321052c0cSHong Zhang Input Parameter: 135421052c0cSHong Zhang + ts - timestepping context 135521052c0cSHong 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 135621052c0cSHong Zhang 135721052c0cSHong Zhang Options Database: 135821052c0cSHong Zhang . -ts_rk_multirate - <true,false> 135921052c0cSHong Zhang 136021052c0cSHong Zhang Notes: 136121052c0cSHong 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). 136221052c0cSHong Zhang 136321052c0cSHong Zhang Level: intermediate 136421052c0cSHong Zhang 136521052c0cSHong Zhang .seealso: TSRKGetMultirate() 136621052c0cSHong Zhang @*/ 136721052c0cSHong Zhang PetscErrorCode TSRKSetMultirate(TS ts,PetscBool use_multirate) 136821052c0cSHong Zhang { 1369f06fb28fSHong Zhang PetscErrorCode ierr; 13707dbacdf2SHong Zhang 13718a4be4dbSHong Zhang PetscFunctionBegin; 1372f06fb28fSHong Zhang ierr = PetscTryMethod(ts,"TSRKSetMultirate_C",(TS,PetscBool),(ts,use_multirate));CHKERRQ(ierr); 137321052c0cSHong Zhang PetscFunctionReturn(0); 137421052c0cSHong Zhang } 137521052c0cSHong Zhang 137621052c0cSHong Zhang /*@C 137721052c0cSHong Zhang TSRKGetMultirate - Gets whether to Use the interpolation-based multirate RK method 137821052c0cSHong Zhang 137921052c0cSHong Zhang Not collective 138021052c0cSHong Zhang 138121052c0cSHong Zhang Input Parameter: 138221052c0cSHong Zhang . ts - timestepping context 138321052c0cSHong Zhang 138421052c0cSHong Zhang Output Parameter: 138521052c0cSHong Zhang . use_multirate - PETSC_TRUE if the multirate RK method is enabled, PETSC_FALSE otherwise 138621052c0cSHong Zhang 138721052c0cSHong Zhang Level: intermediate 138821052c0cSHong Zhang 138921052c0cSHong Zhang .seealso: TSRKSetMultirate() 139021052c0cSHong Zhang @*/ 139121052c0cSHong Zhang PetscErrorCode TSRKGetMultirate(TS ts,PetscBool *use_multirate) 139221052c0cSHong Zhang { 1393f06fb28fSHong Zhang PetscErrorCode ierr; 13947dbacdf2SHong Zhang 13957dbacdf2SHong Zhang PetscFunctionBegin; 1396f06fb28fSHong Zhang ierr = PetscUseMethod(ts,"TSRKGetMultirate_C",(TS,PetscBool*),(ts,use_multirate));CHKERRQ(ierr); 139721052c0cSHong Zhang PetscFunctionReturn(0); 139821052c0cSHong Zhang } 139921052c0cSHong Zhang 1400ebe3b25bSBarry Smith /*MC 1401f68a32c8SEmil Constantinescu TSRK - ODE and DAE solver using Runge-Kutta schemes 1402ebe3b25bSBarry Smith 1403f68a32c8SEmil Constantinescu The user should provide the right hand side of the equation 1404f68a32c8SEmil Constantinescu using TSSetRHSFunction(). 1405ebe3b25bSBarry Smith 1406f68a32c8SEmil Constantinescu Notes: 140798164e13SLisandro Dalcin The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type 1408ebe3b25bSBarry Smith 1409d41222bbSBarry Smith Level: beginner 1410d41222bbSBarry Smith 1411f68a32c8SEmil Constantinescu .seealso: TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3, 14120fe4d17eSHong Zhang TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister(), TSRKSetMultirate(), TSRKGetMultirate() 1413ebe3b25bSBarry Smith 1414ebe3b25bSBarry Smith M*/ 14158cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts) 1416e4dd521cSBarry Smith { 14171566a47fSLisandro Dalcin TS_RK *rk; 1418dfbe8321SBarry Smith PetscErrorCode ierr; 1419e4dd521cSBarry Smith 1420e4dd521cSBarry Smith PetscFunctionBegin; 1421f68a32c8SEmil Constantinescu ierr = TSRKInitializePackage();CHKERRQ(ierr); 1422f68a32c8SEmil Constantinescu 1423f68a32c8SEmil Constantinescu ts->ops->reset = TSReset_RK; 14245f70b478SJed Brown ts->ops->destroy = TSDestroy_RK; 14255f70b478SJed Brown ts->ops->view = TSView_RK; 1426f68a32c8SEmil Constantinescu ts->ops->load = TSLoad_RK; 1427f68a32c8SEmil Constantinescu ts->ops->setup = TSSetUp_RK; 1428f68a32c8SEmil Constantinescu ts->ops->interpolate = TSInterpolate_RK; 14292c9cb062Sluzhanghpp ts->ops->step = TSStep_RK; 1430f68a32c8SEmil Constantinescu ts->ops->evaluatestep = TSEvaluateStep_RK; 1431fecfb714SLisandro Dalcin ts->ops->rollback = TSRollBack_RK; 1432f68a32c8SEmil Constantinescu ts->ops->setfromoptions = TSSetFromOptions_RK; 1433ff22ae23SHong Zhang ts->ops->getstages = TSGetStages_RK; 1434e4dd521cSBarry Smith 14350f7a1166SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_RK; 143613af1a74SHong Zhang ts->ops->adjointsetup = TSAdjointSetUp_RK; 143713af1a74SHong Zhang ts->ops->adjointstep = TSAdjointStep_RK; 143813af1a74SHong Zhang ts->ops->adjointreset = TSAdjointReset_RK; 14390f7a1166SHong Zhang 144013af1a74SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_RK; 1441922a638cSHong Zhang ts->ops->forwardsetup = TSForwardSetUp_RK; 1442922a638cSHong Zhang ts->ops->forwardreset = TSForwardReset_RK; 1443922a638cSHong Zhang ts->ops->forwardstep = TSForwardStep_RK; 1444922a638cSHong Zhang ts->ops->forwardgetstages= TSForwardGetStages_RK; 1445922a638cSHong Zhang 14461566a47fSLisandro Dalcin ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr); 14471566a47fSLisandro Dalcin ts->data = (void*)rk; 1448e4dd521cSBarry Smith 1449f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr); 1450f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr); 145121052c0cSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",TSRKSetMultirate_RK);CHKERRQ(ierr); 145221052c0cSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",TSRKGetMultirate_RK);CHKERRQ(ierr); 1453be5899b3SLisandro Dalcin 1454be5899b3SLisandro Dalcin ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr); 145590b67129SHong Zhang rk->dtratio = 1; 1456e4dd521cSBarry Smith PetscFunctionReturn(0); 1457e4dd521cSBarry Smith } 1458