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*/ 124262ff7bbSBarry Smith 125f68a32c8SEmil Constantinescu /*@C 126f68a32c8SEmil Constantinescu TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK 127262ff7bbSBarry Smith 128f68a32c8SEmil Constantinescu Not Collective, but should be called by all processes which will need the schemes to be registered 129262ff7bbSBarry Smith 130f68a32c8SEmil Constantinescu Level: advanced 131262ff7bbSBarry Smith 132f68a32c8SEmil Constantinescu .keywords: TS, TSRK, register, all 133262ff7bbSBarry Smith 134f68a32c8SEmil Constantinescu .seealso: TSRKRegisterDestroy() 135262ff7bbSBarry Smith @*/ 136f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void) 137262ff7bbSBarry Smith { 1384ac538c5SBarry Smith PetscErrorCode ierr; 139262ff7bbSBarry Smith 140262ff7bbSBarry Smith PetscFunctionBegin; 141f68a32c8SEmil Constantinescu if (TSRKRegisterAllCalled) PetscFunctionReturn(0); 142f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_TRUE; 143e4dd521cSBarry Smith 1444e82626cSLisandro Dalcin #define RC PetscRealConstant 145e4dd521cSBarry Smith { 146f68a32c8SEmil Constantinescu const PetscReal 1474e82626cSLisandro Dalcin A[1][1] = {{0}}, 1484e82626cSLisandro Dalcin b[1] = {RC(1.0)}; 1493a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 150e4dd521cSBarry Smith } 151db046809SBarry Smith { 152f68a32c8SEmil Constantinescu const PetscReal 1534e82626cSLisandro Dalcin A[2][2] = {{0,0}, 1544e82626cSLisandro Dalcin {RC(1.0),0}}, 1554e82626cSLisandro Dalcin b[2] = {RC(0.5),RC(0.5)}, 1564e82626cSLisandro Dalcin bembed[2] = {RC(1.0),0}; 1573a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 158db046809SBarry Smith } 159f68a32c8SEmil Constantinescu { 160f68a32c8SEmil Constantinescu const PetscReal 16117f6384fSBarry Smith A[3][3] = {{0,0,0}, 1624e82626cSLisandro Dalcin {RC(2.0)/RC(3.0),0,0}, 1634e82626cSLisandro Dalcin {RC(-1.0)/RC(3.0),RC(1.0),0}}, 1644e82626cSLisandro Dalcin b[3] = {RC(0.25),RC(0.5),RC(0.25)}; 1653a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 166fdefd5e5SDebojyoti Ghosh } 167fdefd5e5SDebojyoti Ghosh { 168fdefd5e5SDebojyoti Ghosh const PetscReal 16917f6384fSBarry Smith A[4][4] = {{0,0,0,0}, 1704e82626cSLisandro Dalcin {RC(1.0)/RC(2.0),0,0,0}, 1714e82626cSLisandro Dalcin {0,RC(3.0)/RC(4.0),0,0}, 1724e82626cSLisandro Dalcin {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}}, 1734e82626cSLisandro Dalcin b[4] = {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}, 1744e82626cSLisandro 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)}; 1753a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 176db046809SBarry Smith } 177f68a32c8SEmil Constantinescu { 178f68a32c8SEmil Constantinescu const PetscReal 179f68a32c8SEmil Constantinescu A[4][4] = {{0,0,0,0}, 1804e82626cSLisandro Dalcin {RC(0.5),0,0,0}, 1814e82626cSLisandro Dalcin {0,RC(0.5),0,0}, 1824e82626cSLisandro Dalcin {0,0,RC(1.0),0}}, 1834e82626cSLisandro 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)}; 1843a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 185db046809SBarry Smith } 186f68a32c8SEmil Constantinescu { 187f68a32c8SEmil Constantinescu const PetscReal 18817f6384fSBarry Smith A[6][6] = {{0,0,0,0,0,0}, 1894e82626cSLisandro Dalcin {RC(0.25),0,0,0,0,0}, 1904e82626cSLisandro Dalcin {RC(3.0)/RC(32.0),RC(9.0)/RC(32.0),0,0,0,0}, 1914e82626cSLisandro Dalcin {RC(1932.0)/RC(2197.0),RC(-7200.0)/RC(2197.0),RC(7296.0)/RC(2197.0),0,0,0}, 1924e82626cSLisandro Dalcin {RC(439.0)/RC(216.0),RC(-8.0),RC(3680.0)/RC(513.0),RC(-845.0)/RC(4104.0),0,0}, 1934e82626cSLisandro 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}}, 1944e82626cSLisandro 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)}, 1954e82626cSLisandro 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}; 1963a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 197fdefd5e5SDebojyoti Ghosh } 198fdefd5e5SDebojyoti Ghosh { 199fdefd5e5SDebojyoti Ghosh const PetscReal 20017f6384fSBarry Smith A[7][7] = {{0,0,0,0,0,0,0}, 2014e82626cSLisandro Dalcin {RC(1.0)/RC(5.0),0,0,0,0,0,0}, 2024e82626cSLisandro Dalcin {RC(3.0)/RC(40.0),RC(9.0)/RC(40.0),0,0,0,0,0}, 2034e82626cSLisandro Dalcin {RC(44.0)/RC(45.0),RC(-56.0)/RC(15.0),RC(32.0)/RC(9.0),0,0,0,0}, 2044e82626cSLisandro 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}, 2054e82626cSLisandro 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}, 2064e82626cSLisandro 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}}, 2074e82626cSLisandro 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}, 208a685a763Sluzhanghpp 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)}, 209a685a763Sluzhanghpp 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)}, 210a685a763Sluzhanghpp {0,0,0,0,0}, 211a685a763Sluzhanghpp {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)}, 212a685a763Sluzhanghpp {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)}, 213a685a763Sluzhanghpp {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)}, 214a685a763Sluzhanghpp {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)}, 215a685a763Sluzhanghpp {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)}}; 216a685a763Sluzhanghpp 217a685a763Sluzhanghpp ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,5,binterp[0]);CHKERRQ(ierr); 218f68a32c8SEmil Constantinescu } 21905e23783SLisandro Dalcin { 22005e23783SLisandro Dalcin const PetscReal 22117f6384fSBarry Smith A[8][8] = {{0,0,0,0,0,0,0,0}, 2224e82626cSLisandro Dalcin {RC(1.0)/RC(6.0),0,0,0,0,0,0,0}, 2234e82626cSLisandro Dalcin {RC(2.0)/RC(27.0),RC(4.0)/RC(27.0),0,0,0,0,0,0}, 2244e82626cSLisandro Dalcin {RC(183.0)/RC(1372.0),RC(-162.0)/RC(343.0),RC(1053.0)/RC(1372.0),0,0,0,0,0}, 2254e82626cSLisandro 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}, 2264e82626cSLisandro 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}, 2274e82626cSLisandro 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}, 2284e82626cSLisandro 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}}, 2294e82626cSLisandro 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}, 2304e82626cSLisandro 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)}; 23105e23783SLisandro Dalcin ierr = TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 23205e23783SLisandro Dalcin } 2334e82626cSLisandro Dalcin #undef RC 234db046809SBarry Smith PetscFunctionReturn(0); 235db046809SBarry Smith } 236db046809SBarry Smith 237f68a32c8SEmil Constantinescu /*@C 238f68a32c8SEmil Constantinescu TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister(). 239f68a32c8SEmil Constantinescu 240f68a32c8SEmil Constantinescu Not Collective 241f68a32c8SEmil Constantinescu 242f68a32c8SEmil Constantinescu Level: advanced 243f68a32c8SEmil Constantinescu 244f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy 245f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll() 246f68a32c8SEmil Constantinescu @*/ 247f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void) 248e4dd521cSBarry Smith { 249f68a32c8SEmil Constantinescu PetscErrorCode ierr; 250f68a32c8SEmil Constantinescu RKTableauLink link; 251f68a32c8SEmil Constantinescu 252f68a32c8SEmil Constantinescu PetscFunctionBegin; 253f68a32c8SEmil Constantinescu while ((link = RKTableauList)) { 254f68a32c8SEmil Constantinescu RKTableau t = &link->tab; 255f68a32c8SEmil Constantinescu RKTableauList = link->next; 256f68a32c8SEmil Constantinescu ierr = PetscFree3(t->A,t->b,t->c); CHKERRQ(ierr); 257f68a32c8SEmil Constantinescu ierr = PetscFree (t->bembed); CHKERRQ(ierr); 258f68a32c8SEmil Constantinescu ierr = PetscFree (t->binterp); CHKERRQ(ierr); 259f68a32c8SEmil Constantinescu ierr = PetscFree (t->name); CHKERRQ(ierr); 260f68a32c8SEmil Constantinescu ierr = PetscFree (link); CHKERRQ(ierr); 261f68a32c8SEmil Constantinescu } 262f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_FALSE; 263f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 264f68a32c8SEmil Constantinescu } 265f68a32c8SEmil Constantinescu 266f68a32c8SEmil Constantinescu /*@C 267f68a32c8SEmil Constantinescu TSRKInitializePackage - This function initializes everything in the TSRK package. It is called 2688a690491SBarry Smith from TSInitializePackage(). 269f68a32c8SEmil Constantinescu 270f68a32c8SEmil Constantinescu Level: developer 271f68a32c8SEmil Constantinescu 272f68a32c8SEmil Constantinescu .keywords: TS, TSRK, initialize, package 273f68a32c8SEmil Constantinescu .seealso: PetscInitialize() 274f68a32c8SEmil Constantinescu @*/ 275f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void) 276f68a32c8SEmil Constantinescu { 277186e87acSLisandro Dalcin PetscErrorCode ierr; 278e4dd521cSBarry Smith 279e4dd521cSBarry Smith PetscFunctionBegin; 280f68a32c8SEmil Constantinescu if (TSRKPackageInitialized) PetscFunctionReturn(0); 281f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_TRUE; 282f68a32c8SEmil Constantinescu ierr = TSRKRegisterAll();CHKERRQ(ierr); 283f68a32c8SEmil Constantinescu ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr); 284f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 285f68a32c8SEmil Constantinescu } 286186e87acSLisandro Dalcin 287f68a32c8SEmil Constantinescu /*@C 288f68a32c8SEmil Constantinescu TSRKFinalizePackage - This function destroys everything in the TSRK package. It is 289f68a32c8SEmil Constantinescu called from PetscFinalize(). 290186e87acSLisandro Dalcin 291f68a32c8SEmil Constantinescu Level: developer 292f68a32c8SEmil Constantinescu 293f68a32c8SEmil Constantinescu .keywords: Petsc, destroy, package 294f68a32c8SEmil Constantinescu .seealso: PetscFinalize() 295f68a32c8SEmil Constantinescu @*/ 296f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void) 297f68a32c8SEmil Constantinescu { 298f68a32c8SEmil Constantinescu PetscErrorCode ierr; 299f68a32c8SEmil Constantinescu 300f68a32c8SEmil Constantinescu PetscFunctionBegin; 301f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_FALSE; 302f68a32c8SEmil Constantinescu ierr = TSRKRegisterDestroy();CHKERRQ(ierr); 303f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 304f68a32c8SEmil Constantinescu } 305f68a32c8SEmil Constantinescu 306f68a32c8SEmil Constantinescu /*@C 307f68a32c8SEmil Constantinescu TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 308f68a32c8SEmil Constantinescu 309f68a32c8SEmil Constantinescu Not Collective, but the same schemes should be registered on all processes on which they will be used 310f68a32c8SEmil Constantinescu 311f68a32c8SEmil Constantinescu Input Parameters: 312f68a32c8SEmil Constantinescu + name - identifier for method 313f68a32c8SEmil Constantinescu . order - approximation order of method 314f68a32c8SEmil Constantinescu . s - number of stages, this is the dimension of the matrices below 315f68a32c8SEmil Constantinescu . A - stage coefficients (dimension s*s, row-major) 316f68a32c8SEmil Constantinescu . b - step completion table (dimension s; NULL to use last row of A) 317f68a32c8SEmil Constantinescu . c - abscissa (dimension s; NULL to use row sums of A) 318f68a32c8SEmil Constantinescu . bembed - completion table for embedded method (dimension s; NULL if not available) 3193a8a9803SLisandro Dalcin . p - Order of the interpolation scheme, equal to the number of columns of binterp 3203a8a9803SLisandro Dalcin - binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1) 321f68a32c8SEmil Constantinescu 322f68a32c8SEmil Constantinescu Notes: 323f68a32c8SEmil Constantinescu Several RK methods are provided, this function is only needed to create new methods. 324f68a32c8SEmil Constantinescu 325f68a32c8SEmil Constantinescu Level: advanced 326f68a32c8SEmil Constantinescu 327f68a32c8SEmil Constantinescu .keywords: TS, register 328f68a32c8SEmil Constantinescu 329f68a32c8SEmil Constantinescu .seealso: TSRK 330f68a32c8SEmil Constantinescu @*/ 331f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s, 332f68a32c8SEmil Constantinescu const PetscReal A[],const PetscReal b[],const PetscReal c[], 3333a8a9803SLisandro Dalcin const PetscReal bembed[],PetscInt p,const PetscReal binterp[]) 334f68a32c8SEmil Constantinescu { 335f68a32c8SEmil Constantinescu PetscErrorCode ierr; 336f68a32c8SEmil Constantinescu RKTableauLink link; 337f68a32c8SEmil Constantinescu RKTableau t; 338f68a32c8SEmil Constantinescu PetscInt i,j; 339f68a32c8SEmil Constantinescu 340f68a32c8SEmil Constantinescu PetscFunctionBegin; 3413a8a9803SLisandro Dalcin PetscValidCharPointer(name,1); 3423a8a9803SLisandro Dalcin PetscValidRealPointer(A,4); 3433a8a9803SLisandro Dalcin if (b) PetscValidRealPointer(b,5); 3443a8a9803SLisandro Dalcin if (c) PetscValidRealPointer(c,6); 3453a8a9803SLisandro Dalcin if (bembed) PetscValidRealPointer(bembed,7); 3463a8a9803SLisandro Dalcin if (binterp || p > 1) PetscValidRealPointer(binterp,9); 3473a8a9803SLisandro Dalcin 3481d36bdfdSBarry Smith ierr = TSRKInitializePackage();CHKERRQ(ierr); 34995dccacaSBarry Smith ierr = PetscNew(&link);CHKERRQ(ierr); 350f68a32c8SEmil Constantinescu t = &link->tab; 3513a8a9803SLisandro Dalcin 352f68a32c8SEmil Constantinescu ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 353f68a32c8SEmil Constantinescu t->order = order; 354f68a32c8SEmil Constantinescu t->s = s; 355dcca6d9dSJed Brown ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr); 356f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 357f68a32c8SEmil Constantinescu if (b) { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); } 358f68a32c8SEmil Constantinescu else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i]; 359f68a32c8SEmil Constantinescu if (c) { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); } 360f68a32c8SEmil 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]; 361d760c35bSDebojyoti Ghosh t->FSAL = PETSC_TRUE; 362d760c35bSDebojyoti Ghosh for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE; 3633a8a9803SLisandro Dalcin 364f68a32c8SEmil Constantinescu if (bembed) { 365785e854fSJed Brown ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr); 366f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr); 367f68a32c8SEmil Constantinescu } 368f68a32c8SEmil Constantinescu 3693a8a9803SLisandro Dalcin if (!binterp) { p = 1; binterp = t->b; } 3703a8a9803SLisandro Dalcin t->p = p; 3713a8a9803SLisandro Dalcin ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr); 3723a8a9803SLisandro Dalcin ierr = PetscMemcpy(t->binterp,binterp,s*p*sizeof(binterp[0]));CHKERRQ(ierr); 3733a8a9803SLisandro Dalcin 374f68a32c8SEmil Constantinescu link->next = RKTableauList; 375f68a32c8SEmil Constantinescu RKTableauList = link; 376f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 377f68a32c8SEmil Constantinescu } 378f68a32c8SEmil Constantinescu 379e4dd521cSBarry Smith /* 3802c9cb062Sluzhanghpp This is for single-step RK method 381f68a32c8SEmil Constantinescu The step completion formula is 382e4dd521cSBarry Smith 383f68a32c8SEmil Constantinescu x1 = x0 + h b^T YdotRHS 384f68a32c8SEmil Constantinescu 385f68a32c8SEmil Constantinescu This function can be called before or after ts->vec_sol has been updated. 386f68a32c8SEmil Constantinescu Suppose we have a completion formula (b) and an embedded formula (be) of different order. 387f68a32c8SEmil Constantinescu We can write 388f68a32c8SEmil Constantinescu 389f68a32c8SEmil Constantinescu x1e = x0 + h be^T YdotRHS 390f68a32c8SEmil Constantinescu = x1 - h b^T YdotRHS + h be^T YdotRHS 391f68a32c8SEmil Constantinescu = x1 + h (be - b)^T YdotRHS 392f68a32c8SEmil Constantinescu 393f68a32c8SEmil Constantinescu so we can evaluate the method with different order even after the step has been optimistically completed. 394e4dd521cSBarry Smith */ 395f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done) 396f68a32c8SEmil Constantinescu { 397f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 398f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 399f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 400f68a32c8SEmil Constantinescu PetscReal h; 401f68a32c8SEmil Constantinescu PetscInt s = tab->s,j; 402f68a32c8SEmil Constantinescu PetscErrorCode ierr; 403f68a32c8SEmil Constantinescu 404f68a32c8SEmil Constantinescu PetscFunctionBegin; 405f68a32c8SEmil Constantinescu switch (rk->status) { 406f68a32c8SEmil Constantinescu case TS_STEP_INCOMPLETE: 407f68a32c8SEmil Constantinescu case TS_STEP_PENDING: 408f68a32c8SEmil Constantinescu h = ts->time_step; break; 409f68a32c8SEmil Constantinescu case TS_STEP_COMPLETE: 410be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 411f68a32c8SEmil Constantinescu default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 412f68a32c8SEmil Constantinescu } 413f68a32c8SEmil Constantinescu if (order == tab->order) { 414f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { 415f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 41690b67129SHong Zhang for (j=0; j<s; j++) w[j] = h*tab->b[j]/rk->dtratio; 417f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 418f68a32c8SEmil Constantinescu } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 419f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 420f68a32c8SEmil Constantinescu } else if (order == tab->order-1) { 421f68a32c8SEmil Constantinescu if (!tab->bembed) goto unavailable; 422f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { /*Complete with the embedded method (be)*/ 423f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 424f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 425f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 426f68a32c8SEmil Constantinescu } else { /*Rollback and re-complete using (be-b) */ 427f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 428f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 429f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 430f68a32c8SEmil Constantinescu } 431f68a32c8SEmil Constantinescu if (done) *done = PETSC_TRUE; 432f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 433f68a32c8SEmil Constantinescu } 434f68a32c8SEmil Constantinescu unavailable: 435f68a32c8SEmil Constantinescu if (done) *done = PETSC_FALSE; 436a7fac7c2SEmil 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); 437f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 438f68a32c8SEmil Constantinescu } 439f68a32c8SEmil Constantinescu 4400f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts) 4410f7a1166SHong Zhang { 4420f7a1166SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 4430f7a1166SHong Zhang RKTableau tab = rk->tableau; 4440f7a1166SHong Zhang const PetscInt s = tab->s; 4450f7a1166SHong Zhang const PetscReal *b = tab->b,*c = tab->c; 4460f7a1166SHong Zhang Vec *Y = rk->Y; 4470f7a1166SHong Zhang PetscInt i; 4480f7a1166SHong Zhang PetscErrorCode ierr; 4490f7a1166SHong Zhang 4500f7a1166SHong Zhang PetscFunctionBegin; 4510f7a1166SHong Zhang /* backup cost integral */ 4520f7a1166SHong Zhang ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr); 4530f7a1166SHong Zhang for (i=s-1; i>=0; i--) { 4540f7a1166SHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 45551a7f651SHong Zhang ierr = TSComputeCostIntegrand(ts,rk->ptime+rk->time_step*c[i],Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 4560f7a1166SHong Zhang ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 4570f7a1166SHong Zhang } 4580f7a1166SHong Zhang PetscFunctionReturn(0); 4590f7a1166SHong Zhang } 4600f7a1166SHong Zhang 4610f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts) 4620f7a1166SHong Zhang { 4630f7a1166SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 4640f7a1166SHong Zhang RKTableau tab = rk->tableau; 4650f7a1166SHong Zhang const PetscInt s = tab->s; 4660f7a1166SHong Zhang const PetscReal *b = tab->b,*c = tab->c; 4670f7a1166SHong Zhang Vec *Y = rk->Y; 4680f7a1166SHong Zhang PetscInt i; 4690f7a1166SHong Zhang PetscErrorCode ierr; 4700f7a1166SHong Zhang 4710f7a1166SHong Zhang PetscFunctionBegin; 4720f7a1166SHong Zhang for (i=s-1; i>=0; i--) { 4730f7a1166SHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 47451a7f651SHong Zhang ierr = TSComputeCostIntegrand(ts,ts->ptime+ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 4750f7a1166SHong Zhang ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 4760f7a1166SHong Zhang } 4770f7a1166SHong Zhang PetscFunctionReturn(0); 4780f7a1166SHong Zhang } 4790f7a1166SHong Zhang 480fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts) 481fecfb714SLisandro Dalcin { 482fecfb714SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 483fecfb714SLisandro Dalcin RKTableau tab = rk->tableau; 484fecfb714SLisandro Dalcin const PetscInt s = tab->s; 485fecfb714SLisandro Dalcin const PetscReal *b = tab->b; 486fecfb714SLisandro Dalcin PetscScalar *w = rk->work; 487fecfb714SLisandro Dalcin Vec *YdotRHS = rk->YdotRHS; 488fecfb714SLisandro Dalcin PetscInt j; 489fecfb714SLisandro Dalcin PetscReal h; 490fecfb714SLisandro Dalcin PetscErrorCode ierr; 491fecfb714SLisandro Dalcin 492fecfb714SLisandro Dalcin PetscFunctionBegin; 493fecfb714SLisandro Dalcin switch (rk->status) { 494fecfb714SLisandro Dalcin case TS_STEP_INCOMPLETE: 495fecfb714SLisandro Dalcin case TS_STEP_PENDING: 496fecfb714SLisandro Dalcin h = ts->time_step; break; 497fecfb714SLisandro Dalcin case TS_STEP_COMPLETE: 498fecfb714SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 499fecfb714SLisandro Dalcin default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 500fecfb714SLisandro Dalcin } 501fecfb714SLisandro Dalcin for (j=0; j<s; j++) w[j] = -h*b[j]; 502fecfb714SLisandro Dalcin ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr); 503fecfb714SLisandro Dalcin PetscFunctionReturn(0); 504fecfb714SLisandro Dalcin } 505fecfb714SLisandro Dalcin 506f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts) 507f68a32c8SEmil Constantinescu { 508f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 509f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 510f68a32c8SEmil Constantinescu const PetscInt s = tab->s; 511fecfb714SLisandro Dalcin const PetscReal *A = tab->A,*c = tab->c; 512f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 513f68a32c8SEmil Constantinescu Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 514d1905564SLisandro Dalcin PetscBool FSAL = tab->FSAL; 515f68a32c8SEmil Constantinescu TSAdapt adapt; 516fecfb714SLisandro Dalcin PetscInt i,j; 517be5899b3SLisandro Dalcin PetscInt rejections = 0; 518be5899b3SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 519be5899b3SLisandro Dalcin PetscReal next_time_step = ts->time_step; 520f68a32c8SEmil Constantinescu PetscErrorCode ierr; 521f68a32c8SEmil Constantinescu 522f68a32c8SEmil Constantinescu PetscFunctionBegin; 523d1905564SLisandro Dalcin if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE; 524d1905564SLisandro Dalcin if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); } 525d1905564SLisandro Dalcin 526f68a32c8SEmil Constantinescu rk->status = TS_STEP_INCOMPLETE; 527be5899b3SLisandro Dalcin while (!ts->reason && rk->status != TS_STEP_COMPLETE) { 528be5899b3SLisandro Dalcin PetscReal t = ts->ptime; 529f68a32c8SEmil Constantinescu PetscReal h = ts->time_step; 530f68a32c8SEmil Constantinescu for (i=0; i<s; i++) { 5319be3e283SDebojyoti Ghosh rk->stage_time = t + h*c[i]; 5329be3e283SDebojyoti Ghosh ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr); 533f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 534f68a32c8SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 535f68a32c8SEmil Constantinescu ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 5369be3e283SDebojyoti Ghosh ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 537f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 538be5899b3SLisandro Dalcin ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr); 539fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 5408f738a7cSLisandro Dalcin if (FSAL && !i) continue; 541f68a32c8SEmil Constantinescu ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 542f68a32c8SEmil Constantinescu } 543be5899b3SLisandro Dalcin 544be5899b3SLisandro Dalcin rk->status = TS_STEP_INCOMPLETE; 545f68a32c8SEmil Constantinescu ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 546f68a32c8SEmil Constantinescu rk->status = TS_STEP_PENDING; 547f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 548f68a32c8SEmil Constantinescu ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 5491917a363SLisandro Dalcin ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr); 550fecfb714SLisandro Dalcin ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 551be5899b3SLisandro Dalcin rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 552be5899b3SLisandro Dalcin if (!accept) { /* Roll back the current step */ 553fecfb714SLisandro Dalcin ierr = TSRollBack_RK(ts);CHKERRQ(ierr); 554be5899b3SLisandro Dalcin ts->time_step = next_time_step; 555be5899b3SLisandro Dalcin goto reject_step; 556be5899b3SLisandro Dalcin } 557be5899b3SLisandro Dalcin 5580f7a1166SHong Zhang if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 5590f7a1166SHong Zhang rk->ptime = ts->ptime; 5600f7a1166SHong Zhang rk->time_step = ts->time_step; 561493ed95fSHong Zhang } 562be5899b3SLisandro Dalcin 563f68a32c8SEmil Constantinescu ts->ptime += ts->time_step; 564f68a32c8SEmil Constantinescu ts->time_step = next_time_step; 565f68a32c8SEmil Constantinescu break; 566be5899b3SLisandro Dalcin 567be5899b3SLisandro Dalcin reject_step: 568fecfb714SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 569be5899b3SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 570be5899b3SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 571be5899b3SLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 572f68a32c8SEmil Constantinescu } 573f68a32c8SEmil Constantinescu } 574f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 575e4dd521cSBarry Smith } 576e4dd521cSBarry Smith 577f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts) 578f6a906c0SBarry Smith { 579f6a906c0SBarry Smith TS_RK *rk = (TS_RK*)ts->data; 580be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 581be5899b3SLisandro Dalcin PetscInt s = tab->s; 582f6a906c0SBarry Smith PetscErrorCode ierr; 583f6a906c0SBarry Smith 584f6a906c0SBarry Smith PetscFunctionBegin; 585f6a906c0SBarry Smith if (ts->adjointsetupcalled++) PetscFunctionReturn(0); 5862e7b7f96SHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecsDeltaLam);CHKERRQ(ierr); 5872e7b7f96SHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&rk->VecsSensiTemp);CHKERRQ(ierr); 588f6a906c0SBarry Smith if(ts->vecs_sensip) { 5892e7b7f96SHong Zhang ierr = VecDuplicate(ts->vecs_sensip[0],&rk->VecDeltaMu);CHKERRQ(ierr); 590f6a906c0SBarry Smith } 591f6a906c0SBarry Smith PetscFunctionReturn(0); 592f6a906c0SBarry Smith } 593f6a906c0SBarry Smith 59442f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts) 595d2daff3dSHong Zhang { 596c235aa8dSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 597c235aa8dSHong Zhang RKTableau tab = rk->tableau; 5982e7b7f96SHong Zhang Mat J; 599c235aa8dSHong Zhang const PetscInt s = tab->s; 600c235aa8dSHong Zhang const PetscReal *A = tab->A,*b = tab->b,*c = tab->c; 601c235aa8dSHong Zhang PetscScalar *w = rk->work; 6022e7b7f96SHong Zhang Vec *Y = rk->Y,*VecsDeltaLam = rk->VecsDeltaLam,VecDeltaMu = rk->VecDeltaMu,*VecsSensiTemp = rk->VecsSensiTemp; 603b18ea86cSHong Zhang PetscInt i,j,nadj; 6043d3eaba7SBarry Smith PetscReal t = ts->ptime; 6052e7b7f96SHong Zhang PetscReal h = ts->time_step,stage_time; 6062e7b7f96SHong Zhang PetscBool zero; 607d2daff3dSHong Zhang PetscErrorCode ierr; 608c235aa8dSHong Zhang 609d2daff3dSHong Zhang PetscFunctionBegin; 610c235aa8dSHong Zhang rk->status = TS_STEP_INCOMPLETE; 6113389c7d9SStefano Zampini 6122e7b7f96SHong Zhang for (i=s-1; i>=0; i--) { 613*6a1e4597SHong Zhang if (tab->FSAL && i == s-1) { 614*6a1e4597SHong Zhang /* VecsDeltaLam[nadj*s+s-1] are initialized with zeros and the values never change.*/ 615*6a1e4597SHong Zhang continue; 616*6a1e4597SHong Zhang } 6172e7b7f96SHong Zhang stage_time = t + h*(1.0-c[i]); 6182e7b7f96SHong Zhang zero = PETSC_FALSE; 6193389c7d9SStefano Zampini ierr = TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr); 6203389c7d9SStefano Zampini ierr = TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);CHKERRQ(ierr); 621493ed95fSHong Zhang if (ts->vec_costintegral) { 622c9ad9de0SHong Zhang ierr = TSComputeDRDUFunction(ts,stage_time,Y[i],ts->vecs_drdu);CHKERRQ(ierr); 623493ed95fSHong Zhang } 6243389c7d9SStefano Zampini if (ts->vecs_sensip) { 6252e7b7f96SHong Zhang ierr = TSComputeRHSJacobianP(ts,stage_time,Y[i],ts->Jacp);CHKERRQ(ierr); /* get f_p */ 6263389c7d9SStefano Zampini if (ts->vec_costintegral) { 627c9ad9de0SHong Zhang ierr = TSComputeDRDPFunction(ts,stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr); 6283389c7d9SStefano Zampini } 6293389c7d9SStefano Zampini } 6303389c7d9SStefano Zampini 6313389c7d9SStefano Zampini if (b[i] == 0 && i == s-1) zero = PETSC_TRUE; 6323389c7d9SStefano Zampini 6332e7b7f96SHong Zhang for (j=i+1; j<s; j++) w[j-i-1] = -h*A[j*s+i]; /* coefficients for computing VecsSensiTemp */ 6342e7b7f96SHong Zhang 6352e7b7f96SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 6363389c7d9SStefano Zampini /* Stage values of lambda */ 6373389c7d9SStefano Zampini if (!zero) { 6382e7b7f96SHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); /* VecDeltaLam is an vec array of size s by numcost */ 6392e7b7f96SHong Zhang ierr = VecScale(VecsSensiTemp[nadj],-h*b[i]);CHKERRQ(ierr); 6402e7b7f96SHong Zhang ierr = VecMAXPY(VecsSensiTemp[nadj],s-i-1,w,&VecsDeltaLam[nadj*s+i+1]);CHKERRQ(ierr); 6412e7b7f96SHong Zhang ierr = MatMultTranspose(J,VecsSensiTemp[nadj],VecsDeltaLam[nadj*s+i]);CHKERRQ(ierr); 6423389c7d9SStefano Zampini } else { 643*6a1e4597SHong Zhang /* \sum_{j=i+1}^s a_{ji}*lambda_{s,j} */ 644*6a1e4597SHong Zhang ierr = VecSet(VecsSensiTemp[nadj],0);CHKERRQ(ierr); 645*6a1e4597SHong Zhang ierr = VecMAXPY(VecsSensiTemp[nadj],s-i-1,w,&VecsDeltaLam[nadj*s+i+1]);CHKERRQ(ierr); 646*6a1e4597SHong Zhang ierr = MatMultTranspose(J,VecsSensiTemp[nadj],VecsDeltaLam[nadj*s+i]);CHKERRQ(ierr); 647*6a1e4597SHong Zhang ierr = VecScale(VecsDeltaLam[nadj*s+i],-h);CHKERRQ(ierr); 6483389c7d9SStefano Zampini } 649493ed95fSHong Zhang if (ts->vec_costintegral) { 6502e7b7f96SHong Zhang ierr = VecAXPY(VecsDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdu[nadj]);CHKERRQ(ierr); 651493ed95fSHong Zhang } 6526497ce14SHong Zhang 653ad8e2604SHong Zhang /* Stage values of mu */ 6546497ce14SHong Zhang if (ts->vecs_sensip) { 6553389c7d9SStefano Zampini if (!zero) { 6562e7b7f96SHong Zhang ierr = MatMultTranspose(ts->Jacp,VecsSensiTemp[nadj],VecDeltaMu);CHKERRQ(ierr); 6573389c7d9SStefano Zampini } else { 6582e7b7f96SHong Zhang ierr = VecSet(VecDeltaMu,0);CHKERRQ(ierr); 659493ed95fSHong Zhang } 660493ed95fSHong Zhang if (ts->vec_costintegral) { 6612e7b7f96SHong Zhang ierr = VecAXPY(VecDeltaMu,-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr); 662493ed95fSHong Zhang } 6632e7b7f96SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],1.,VecDeltaMu);CHKERRQ(ierr); /* update sensip for each stage */ 664ad8e2604SHong Zhang } 665c235aa8dSHong Zhang } 6666497ce14SHong Zhang } 667c235aa8dSHong Zhang 668c235aa8dSHong Zhang for (j=0; j<s; j++) w[j] = 1.0; 6692e7b7f96SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { /* no need to do this for mu's */ 6702e7b7f96SHong Zhang ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecsDeltaLam[nadj*s]);CHKERRQ(ierr); 6716497ce14SHong Zhang } 672c235aa8dSHong Zhang rk->status = TS_STEP_COMPLETE; 673d2daff3dSHong Zhang PetscFunctionReturn(0); 674d2daff3dSHong Zhang } 675d2daff3dSHong Zhang 67655de54fcSHong Zhang static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X) 67755de54fcSHong Zhang { 67855de54fcSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 67955de54fcSHong Zhang PetscInt s = rk->tableau->s,p = rk->tableau->p,i,j; 68055de54fcSHong Zhang PetscReal h; 68155de54fcSHong Zhang PetscReal tt,t; 68255de54fcSHong Zhang PetscScalar *b; 68355de54fcSHong Zhang const PetscReal *B = rk->tableau->binterp; 68455de54fcSHong Zhang PetscErrorCode ierr; 68555de54fcSHong Zhang 68655de54fcSHong Zhang PetscFunctionBegin; 68755de54fcSHong Zhang if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 68855de54fcSHong Zhang 68955de54fcSHong Zhang switch (rk->status) { 69055de54fcSHong Zhang case TS_STEP_INCOMPLETE: 69155de54fcSHong Zhang case TS_STEP_PENDING: 69255de54fcSHong Zhang h = ts->time_step; 69355de54fcSHong Zhang t = (itime - ts->ptime)/h; 69455de54fcSHong Zhang break; 69555de54fcSHong Zhang case TS_STEP_COMPLETE: 69655de54fcSHong Zhang h = ts->ptime - ts->ptime_prev; 69755de54fcSHong Zhang t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 69855de54fcSHong Zhang break; 69955de54fcSHong Zhang default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 70055de54fcSHong Zhang } 70155de54fcSHong Zhang ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 70255de54fcSHong Zhang for (i=0; i<s; i++) b[i] = 0; 70355de54fcSHong Zhang for (j=0,tt=t; j<p; j++,tt*=t) { 70455de54fcSHong Zhang for (i=0; i<s; i++) { 70555de54fcSHong Zhang b[i] += h * B[i*p+j] * tt; 70655de54fcSHong Zhang } 70755de54fcSHong Zhang } 70855de54fcSHong Zhang ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr); 70955de54fcSHong Zhang ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr); 71055de54fcSHong Zhang ierr = PetscFree(b);CHKERRQ(ierr); 71155de54fcSHong Zhang PetscFunctionReturn(0); 71255de54fcSHong Zhang } 71355de54fcSHong Zhang 71455de54fcSHong Zhang /*------------------------------------------------------------*/ 71555de54fcSHong Zhang 716be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts) 717be5899b3SLisandro Dalcin { 718be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 719be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 720be5899b3SLisandro Dalcin PetscErrorCode ierr; 721be5899b3SLisandro Dalcin 722be5899b3SLisandro Dalcin PetscFunctionBegin; 723be5899b3SLisandro Dalcin if (!tab) PetscFunctionReturn(0); 724be5899b3SLisandro Dalcin ierr = PetscFree(rk->work);CHKERRQ(ierr); 725be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr); 726be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr); 7272e7b7f96SHong Zhang ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam);CHKERRQ(ierr); 7282e7b7f96SHong Zhang ierr = VecDestroyVecs(ts->numcost,&rk->VecsSensiTemp);CHKERRQ(ierr); 7292e7b7f96SHong Zhang ierr = VecDestroy(&rk->VecDeltaMu);CHKERRQ(ierr); 730be5899b3SLisandro Dalcin PetscFunctionReturn(0); 731be5899b3SLisandro Dalcin } 732be5899b3SLisandro Dalcin 733277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts) 734e4dd521cSBarry Smith { 7355f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 7366849ba73SBarry Smith PetscErrorCode ierr; 737e4dd521cSBarry Smith 738e4dd521cSBarry Smith PetscFunctionBegin; 739be5899b3SLisandro Dalcin ierr = TSRKTableauReset(ts);CHKERRQ(ierr); 7400f7a1166SHong Zhang ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr); 7410fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 74202555c94SHong Zhang ierr = PetscTryMethod(ts,"TSReset_RK_MultirateSplit_C",(TS),(ts));CHKERRQ(ierr); 7430fe4d17eSHong Zhang } else { 74402555c94SHong Zhang ierr = PetscTryMethod(ts,"TSReset_RK_MultirateNonsplit_C",(TS),(ts));CHKERRQ(ierr); 7450fe4d17eSHong Zhang } 746277b19d0SLisandro Dalcin PetscFunctionReturn(0); 747e4dd521cSBarry Smith } 748277b19d0SLisandro Dalcin 749f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx) 750f68a32c8SEmil Constantinescu { 751f68a32c8SEmil Constantinescu PetscFunctionBegin; 752f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 753f68a32c8SEmil Constantinescu } 754f68a32c8SEmil Constantinescu 755f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 756f68a32c8SEmil Constantinescu { 757f68a32c8SEmil Constantinescu PetscFunctionBegin; 758f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 759f68a32c8SEmil Constantinescu } 760f68a32c8SEmil Constantinescu 761f68a32c8SEmil Constantinescu 762f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx) 763f68a32c8SEmil Constantinescu { 764f68a32c8SEmil Constantinescu PetscFunctionBegin; 765f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 766f68a32c8SEmil Constantinescu } 767f68a32c8SEmil Constantinescu 768f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 769f68a32c8SEmil Constantinescu { 770f68a32c8SEmil Constantinescu 771f68a32c8SEmil Constantinescu PetscFunctionBegin; 772f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 773f68a32c8SEmil Constantinescu } 774c235aa8dSHong Zhang /* 775d2daff3dSHong Zhang static PetscErrorCode RKSetAdjCoe(RKTableau tab) 776d2daff3dSHong Zhang { 777d2daff3dSHong Zhang PetscReal *A,*b; 778d2daff3dSHong Zhang PetscInt s,i,j; 779d2daff3dSHong Zhang PetscErrorCode ierr; 780d2daff3dSHong Zhang 781d2daff3dSHong Zhang PetscFunctionBegin; 782d2daff3dSHong Zhang s = tab->s; 783d2daff3dSHong Zhang ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr); 784d2daff3dSHong Zhang 785d2daff3dSHong Zhang for (i=0; i<s; i++) 786d2daff3dSHong Zhang for (j=0; j<s; j++) { 787d2daff3dSHong 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]; 788d2daff3dSHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr); 789d2daff3dSHong Zhang } 790d2daff3dSHong Zhang 791d2daff3dSHong Zhang for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i]; 792d2daff3dSHong Zhang 793d2daff3dSHong Zhang ierr = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 794d2daff3dSHong Zhang ierr = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr); 795d2daff3dSHong Zhang ierr = PetscFree2(A,b);CHKERRQ(ierr); 796d2daff3dSHong Zhang PetscFunctionReturn(0); 797d2daff3dSHong Zhang } 798c235aa8dSHong Zhang */ 799be5899b3SLisandro Dalcin 800be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts) 801be5899b3SLisandro Dalcin { 802be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 803be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 804be5899b3SLisandro Dalcin PetscErrorCode ierr; 805be5899b3SLisandro Dalcin 806be5899b3SLisandro Dalcin PetscFunctionBegin; 807be5899b3SLisandro Dalcin ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr); 808be5899b3SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr); 809be5899b3SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr); 810be5899b3SLisandro Dalcin PetscFunctionReturn(0); 811be5899b3SLisandro Dalcin } 812be5899b3SLisandro Dalcin 813f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts) 814f68a32c8SEmil Constantinescu { 815f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 816f68a32c8SEmil Constantinescu PetscErrorCode ierr; 817f68a32c8SEmil Constantinescu DM dm; 818f68a32c8SEmil Constantinescu 819f68a32c8SEmil Constantinescu PetscFunctionBegin; 8203dd83b38SBarry Smith ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr); 821be5899b3SLisandro Dalcin ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr); 8220f7a1166SHong Zhang if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 8230f7a1166SHong Zhang ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr); 8240f7a1166SHong Zhang } 825f68a32c8SEmil Constantinescu ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 826f68a32c8SEmil Constantinescu ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 827f68a32c8SEmil Constantinescu ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 8280fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 82902555c94SHong Zhang ierr = PetscTryMethod(ts,"TSSetUp_RK_MultirateSplit_C",(TS),(ts));CHKERRQ(ierr); 8300fe4d17eSHong Zhang } else { 83102555c94SHong Zhang ierr = PetscTryMethod(ts,"TSSetUp_RK_MultirateNonsplit_C",(TS),(ts));CHKERRQ(ierr); 8320fe4d17eSHong Zhang } 833e4dd521cSBarry Smith PetscFunctionReturn(0); 834e4dd521cSBarry Smith } 835d2daff3dSHong Zhang 8364416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts) 837e4dd521cSBarry Smith { 838be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 839dfbe8321SBarry Smith PetscErrorCode ierr; 840262ff7bbSBarry Smith 841e4dd521cSBarry Smith PetscFunctionBegin; 842e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr); 843f68a32c8SEmil Constantinescu { 844f68a32c8SEmil Constantinescu RKTableauLink link; 845f68a32c8SEmil Constantinescu PetscInt count,choice; 8460fe4d17eSHong Zhang PetscBool flg,use_multirate = PETSC_FALSE; 847f68a32c8SEmil Constantinescu const char **namelist; 8482c9cb062Sluzhanghpp 849f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) ; 850f489ac74SBarry Smith ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr); 851f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 8520fe4d17eSHong Zhang ierr = PetscOptionsBool("-ts_rk_multirate","Use interpolation-based multirate RK method","TSRKSetMultirate",rk->use_multirate,&use_multirate,&flg);CHKERRQ(ierr); 8530fe4d17eSHong Zhang if (flg) { 8540fe4d17eSHong Zhang ierr = TSRKSetMultirate(ts,use_multirate);CHKERRQ(ierr); 8550fe4d17eSHong Zhang } 856be5899b3SLisandro Dalcin ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr); 857be5899b3SLisandro Dalcin if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);} 858f68a32c8SEmil Constantinescu ierr = PetscFree(namelist);CHKERRQ(ierr); 859f68a32c8SEmil Constantinescu } 860262ff7bbSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 8612c9cb062Sluzhanghpp ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Multirate methods options","");CHKERRQ(ierr); 8622c9cb062Sluzhanghpp ierr = PetscOptionsInt("-ts_rk_dtratio","time step ratio between slow and fast","",rk->dtratio,&rk->dtratio,NULL);CHKERRQ(ierr); 8632c9cb062Sluzhanghpp ierr = PetscOptionsEnd();CHKERRQ(ierr); 864e4dd521cSBarry Smith PetscFunctionReturn(0); 865e4dd521cSBarry Smith } 866e4dd521cSBarry Smith 8675f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer) 868e4dd521cSBarry Smith { 8695f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 8708370ee3bSLisandro Dalcin PetscBool iascii; 871dfbe8321SBarry Smith PetscErrorCode ierr; 8722cdf8120Sasbjorn 873e4dd521cSBarry Smith PetscFunctionBegin; 874251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 8758370ee3bSLisandro Dalcin if (iascii) { 8769c334d8fSLisandro Dalcin RKTableau tab = rk->tableau; 877f68a32c8SEmil Constantinescu TSRKType rktype; 878f68a32c8SEmil Constantinescu char buf[512]; 879f68a32c8SEmil Constantinescu ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr); 880efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RK type %s\n",rktype);CHKERRQ(ierr); 8810c37eb8eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," Order: %D\n",tab->order);CHKERRQ(ierr); 8820c37eb8eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr); 883f68a32c8SEmil Constantinescu ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 884f68a32c8SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," Abscissa c = %s\n",buf);CHKERRQ(ierr); 8858370ee3bSLisandro Dalcin } 886f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 887f68a32c8SEmil Constantinescu } 888f68a32c8SEmil Constantinescu 889f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer) 890f68a32c8SEmil Constantinescu { 891f68a32c8SEmil Constantinescu PetscErrorCode ierr; 8929c334d8fSLisandro Dalcin TSAdapt adapt; 893f68a32c8SEmil Constantinescu 894f68a32c8SEmil Constantinescu PetscFunctionBegin; 8959c334d8fSLisandro Dalcin ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 8969c334d8fSLisandro Dalcin ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr); 897f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 898f68a32c8SEmil Constantinescu } 899f68a32c8SEmil Constantinescu 900f68a32c8SEmil Constantinescu /*@C 901f68a32c8SEmil Constantinescu TSRKSetType - Set the type of RK scheme 902f68a32c8SEmil Constantinescu 903f68a32c8SEmil Constantinescu Logically collective 904f68a32c8SEmil Constantinescu 905f68a32c8SEmil Constantinescu Input Parameter: 906f68a32c8SEmil Constantinescu + ts - timestepping context 907f68a32c8SEmil Constantinescu - rktype - type of RK-scheme 908f68a32c8SEmil Constantinescu 909287c2655SBarry Smith Options Database: 9109bd3e852SBarry Smith . -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs> 911287c2655SBarry Smith 912f68a32c8SEmil Constantinescu Level: intermediate 913f68a32c8SEmil Constantinescu 914287c2655SBarry Smith .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS 915f68a32c8SEmil Constantinescu @*/ 916f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype) 917f68a32c8SEmil Constantinescu { 918f68a32c8SEmil Constantinescu PetscErrorCode ierr; 919f68a32c8SEmil Constantinescu 920f68a32c8SEmil Constantinescu PetscFunctionBegin; 921f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 922b92453a8SLisandro Dalcin PetscValidCharPointer(rktype,2); 923f68a32c8SEmil Constantinescu ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr); 924f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 925f68a32c8SEmil Constantinescu } 926f68a32c8SEmil Constantinescu 927f68a32c8SEmil Constantinescu /*@C 928f68a32c8SEmil Constantinescu TSRKGetType - Get the type of RK scheme 929f68a32c8SEmil Constantinescu 930f68a32c8SEmil Constantinescu Logically collective 931f68a32c8SEmil Constantinescu 932f68a32c8SEmil Constantinescu Input Parameter: 933f68a32c8SEmil Constantinescu . ts - timestepping context 934f68a32c8SEmil Constantinescu 935f68a32c8SEmil Constantinescu Output Parameter: 936f68a32c8SEmil Constantinescu . rktype - type of RK-scheme 937f68a32c8SEmil Constantinescu 938f68a32c8SEmil Constantinescu Level: intermediate 939f68a32c8SEmil Constantinescu 940f68a32c8SEmil Constantinescu .seealso: TSRKGetType() 941f68a32c8SEmil Constantinescu @*/ 942f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype) 943f68a32c8SEmil Constantinescu { 944f68a32c8SEmil Constantinescu PetscErrorCode ierr; 945f68a32c8SEmil Constantinescu 946f68a32c8SEmil Constantinescu PetscFunctionBegin; 947f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 948f68a32c8SEmil Constantinescu ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr); 949f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 950f68a32c8SEmil Constantinescu } 951f68a32c8SEmil Constantinescu 952560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype) 953f68a32c8SEmil Constantinescu { 954f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 955f68a32c8SEmil Constantinescu 956f68a32c8SEmil Constantinescu PetscFunctionBegin; 957f68a32c8SEmil Constantinescu *rktype = rk->tableau->name; 958f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 959f68a32c8SEmil Constantinescu } 9602c9cb062Sluzhanghpp 961560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype) 962f68a32c8SEmil Constantinescu { 963f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 964f68a32c8SEmil Constantinescu PetscErrorCode ierr; 965f68a32c8SEmil Constantinescu PetscBool match; 966f68a32c8SEmil Constantinescu RKTableauLink link; 967f68a32c8SEmil Constantinescu 968f68a32c8SEmil Constantinescu PetscFunctionBegin; 969f68a32c8SEmil Constantinescu if (rk->tableau) { 970f68a32c8SEmil Constantinescu ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr); 971f68a32c8SEmil Constantinescu if (match) PetscFunctionReturn(0); 972f68a32c8SEmil Constantinescu } 973f68a32c8SEmil Constantinescu for (link = RKTableauList; link; link=link->next) { 974f68a32c8SEmil Constantinescu ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr); 975f68a32c8SEmil Constantinescu if (match) { 976be5899b3SLisandro Dalcin if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);} 977f68a32c8SEmil Constantinescu rk->tableau = &link->tab; 978be5899b3SLisandro Dalcin if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);} 9792ffb9264SLisandro Dalcin ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 980f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 981f68a32c8SEmil Constantinescu } 982f68a32c8SEmil Constantinescu } 983f68a32c8SEmil Constantinescu SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype); 984e4dd521cSBarry Smith PetscFunctionReturn(0); 985e4dd521cSBarry Smith } 986e4dd521cSBarry Smith 987ff22ae23SHong Zhang static PetscErrorCode TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y) 988ff22ae23SHong Zhang { 989ff22ae23SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 990ff22ae23SHong Zhang 991ff22ae23SHong Zhang PetscFunctionBegin; 9920429704eSStefano Zampini if (ns) *ns = rk->tableau->s; 993d2daff3dSHong Zhang if (Y) *Y = rk->Y; 994ff22ae23SHong Zhang PetscFunctionReturn(0); 995ff22ae23SHong Zhang } 996ff22ae23SHong Zhang 997b3a6b972SJed Brown static PetscErrorCode TSDestroy_RK(TS ts) 998b3a6b972SJed Brown { 999b3a6b972SJed Brown PetscErrorCode ierr; 1000b3a6b972SJed Brown 1001b3a6b972SJed Brown PetscFunctionBegin; 1002b3a6b972SJed Brown ierr = TSReset_RK(ts);CHKERRQ(ierr); 1003b3a6b972SJed Brown if (ts->dm) { 1004b3a6b972SJed Brown ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 1005b3a6b972SJed Brown ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 1006b3a6b972SJed Brown } 1007b3a6b972SJed Brown ierr = PetscFree(ts->data);CHKERRQ(ierr); 1008b3a6b972SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr); 1009b3a6b972SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr); 10100fe4d17eSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",NULL);CHKERRQ(ierr); 10110fe4d17eSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",NULL);CHKERRQ(ierr); 1012b3a6b972SJed Brown PetscFunctionReturn(0); 1013b3a6b972SJed Brown } 1014ff22ae23SHong Zhang 101521052c0cSHong Zhang /*@C 101621052c0cSHong Zhang TSRKSetMultirate - Use the interpolation-based multirate RK method 101721052c0cSHong Zhang 101821052c0cSHong Zhang Logically collective 101921052c0cSHong Zhang 102021052c0cSHong Zhang Input Parameter: 102121052c0cSHong Zhang + ts - timestepping context 102221052c0cSHong 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 102321052c0cSHong Zhang 102421052c0cSHong Zhang Options Database: 102521052c0cSHong Zhang . -ts_rk_multirate - <true,false> 102621052c0cSHong Zhang 102721052c0cSHong Zhang Notes: 102821052c0cSHong 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). 102921052c0cSHong Zhang 103021052c0cSHong Zhang Level: intermediate 103121052c0cSHong Zhang 103221052c0cSHong Zhang .seealso: TSRKGetMultirate() 103321052c0cSHong Zhang @*/ 103421052c0cSHong Zhang PetscErrorCode TSRKSetMultirate(TS ts,PetscBool use_multirate) 103521052c0cSHong Zhang { 1036f06fb28fSHong Zhang PetscErrorCode ierr; 10377dbacdf2SHong Zhang 10388a4be4dbSHong Zhang PetscFunctionBegin; 1039f06fb28fSHong Zhang ierr = PetscTryMethod(ts,"TSRKSetMultirate_C",(TS,PetscBool),(ts,use_multirate));CHKERRQ(ierr); 104021052c0cSHong Zhang PetscFunctionReturn(0); 104121052c0cSHong Zhang } 104221052c0cSHong Zhang 104321052c0cSHong Zhang /*@C 104421052c0cSHong Zhang TSRKGetMultirate - Gets whether to Use the interpolation-based multirate RK method 104521052c0cSHong Zhang 104621052c0cSHong Zhang Not collective 104721052c0cSHong Zhang 104821052c0cSHong Zhang Input Parameter: 104921052c0cSHong Zhang . ts - timestepping context 105021052c0cSHong Zhang 105121052c0cSHong Zhang Output Parameter: 105221052c0cSHong Zhang . use_multirate - PETSC_TRUE if the multirate RK method is enabled, PETSC_FALSE otherwise 105321052c0cSHong Zhang 105421052c0cSHong Zhang Level: intermediate 105521052c0cSHong Zhang 105621052c0cSHong Zhang .seealso: TSRKSetMultirate() 105721052c0cSHong Zhang @*/ 105821052c0cSHong Zhang PetscErrorCode TSRKGetMultirate(TS ts,PetscBool *use_multirate) 105921052c0cSHong Zhang { 1060f06fb28fSHong Zhang PetscErrorCode ierr; 10617dbacdf2SHong Zhang 10627dbacdf2SHong Zhang PetscFunctionBegin; 1063f06fb28fSHong Zhang ierr = PetscUseMethod(ts,"TSRKGetMultirate_C",(TS,PetscBool*),(ts,use_multirate));CHKERRQ(ierr); 106421052c0cSHong Zhang PetscFunctionReturn(0); 106521052c0cSHong Zhang } 106621052c0cSHong Zhang 1067ebe3b25bSBarry Smith /*MC 1068f68a32c8SEmil Constantinescu TSRK - ODE and DAE solver using Runge-Kutta schemes 1069ebe3b25bSBarry Smith 1070f68a32c8SEmil Constantinescu The user should provide the right hand side of the equation 1071f68a32c8SEmil Constantinescu using TSSetRHSFunction(). 1072ebe3b25bSBarry Smith 1073f68a32c8SEmil Constantinescu Notes: 107498164e13SLisandro Dalcin The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type 1075ebe3b25bSBarry Smith 1076d41222bbSBarry Smith Level: beginner 1077d41222bbSBarry Smith 1078f68a32c8SEmil Constantinescu .seealso: TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3, 10790fe4d17eSHong Zhang TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister(), TSRKSetMultirate(), TSRKGetMultirate() 1080ebe3b25bSBarry Smith 1081ebe3b25bSBarry Smith M*/ 10828cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts) 1083e4dd521cSBarry Smith { 10841566a47fSLisandro Dalcin TS_RK *rk; 1085dfbe8321SBarry Smith PetscErrorCode ierr; 1086e4dd521cSBarry Smith 1087e4dd521cSBarry Smith PetscFunctionBegin; 1088f68a32c8SEmil Constantinescu ierr = TSRKInitializePackage();CHKERRQ(ierr); 1089f68a32c8SEmil Constantinescu 1090f68a32c8SEmil Constantinescu ts->ops->reset = TSReset_RK; 10915f70b478SJed Brown ts->ops->destroy = TSDestroy_RK; 10925f70b478SJed Brown ts->ops->view = TSView_RK; 1093f68a32c8SEmil Constantinescu ts->ops->load = TSLoad_RK; 1094f68a32c8SEmil Constantinescu ts->ops->setup = TSSetUp_RK; 109542f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_RK; 1096f68a32c8SEmil Constantinescu ts->ops->interpolate = TSInterpolate_RK; 10972c9cb062Sluzhanghpp ts->ops->step = TSStep_RK; 1098f68a32c8SEmil Constantinescu ts->ops->evaluatestep = TSEvaluateStep_RK; 1099fecfb714SLisandro Dalcin ts->ops->rollback = TSRollBack_RK; 1100f68a32c8SEmil Constantinescu ts->ops->setfromoptions = TSSetFromOptions_RK; 1101ff22ae23SHong Zhang ts->ops->getstages = TSGetStages_RK; 110242f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_RK; 1103e4dd521cSBarry Smith 11040f7a1166SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_RK; 11050f7a1166SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_RK; 11060f7a1166SHong Zhang 11071566a47fSLisandro Dalcin ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr); 11081566a47fSLisandro Dalcin ts->data = (void*)rk; 1109e4dd521cSBarry Smith 1110f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr); 1111f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr); 111221052c0cSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",TSRKSetMultirate_RK);CHKERRQ(ierr); 111321052c0cSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",TSRKGetMultirate_RK);CHKERRQ(ierr); 1114be5899b3SLisandro Dalcin 1115be5899b3SLisandro Dalcin ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr); 111690b67129SHong Zhang rk->dtratio = 1; 1117e4dd521cSBarry Smith PetscFunctionReturn(0); 1118e4dd521cSBarry Smith } 1119