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> 14f68a32c8SEmil Constantinescu 15484bcdc7SDebojyoti Ghosh static TSRKType TSRKDefault = TSRK3BS; 16f68a32c8SEmil Constantinescu static PetscBool TSRKRegisterAllCalled; 17f68a32c8SEmil Constantinescu static PetscBool TSRKPackageInitialized; 18f68a32c8SEmil Constantinescu 19f68a32c8SEmil Constantinescu /*MC 20287c2655SBarry Smith TSRK1FE - First order forward Euler scheme. 21262ff7bbSBarry Smith 22f68a32c8SEmil Constantinescu This method has one stage. 23f68a32c8SEmil Constantinescu 24287c2655SBarry Smith Options database: 259bd3e852SBarry Smith . -ts_rk_type 1fe 26287c2655SBarry Smith 27f68a32c8SEmil Constantinescu Level: advanced 28f68a32c8SEmil Constantinescu 29287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 30f68a32c8SEmil Constantinescu M*/ 31f68a32c8SEmil Constantinescu /*MC 322109b73fSDebojyoti Ghosh TSRK2A - Second order RK scheme. 33f68a32c8SEmil Constantinescu 34f68a32c8SEmil Constantinescu This method has two stages. 35f68a32c8SEmil Constantinescu 36287c2655SBarry Smith Options database: 379bd3e852SBarry Smith . -ts_rk_type 2a 38287c2655SBarry Smith 39f68a32c8SEmil Constantinescu Level: advanced 40f68a32c8SEmil Constantinescu 41287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 42f68a32c8SEmil Constantinescu M*/ 43f68a32c8SEmil Constantinescu /*MC 44f68a32c8SEmil Constantinescu TSRK3 - Third order RK scheme. 45f68a32c8SEmil Constantinescu 46f68a32c8SEmil Constantinescu This method has three stages. 47f68a32c8SEmil Constantinescu 48287c2655SBarry Smith Options database: 499bd3e852SBarry Smith . -ts_rk_type 3 50287c2655SBarry Smith 51f68a32c8SEmil Constantinescu Level: advanced 52f68a32c8SEmil Constantinescu 53287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 54f68a32c8SEmil Constantinescu M*/ 55f68a32c8SEmil Constantinescu /*MC 562109b73fSDebojyoti Ghosh TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method. 572109b73fSDebojyoti Ghosh 58d1905564SLisandro Dalcin This method has four stages with the First Same As Last (FSAL) property. 592109b73fSDebojyoti Ghosh 60287c2655SBarry Smith Options database: 619bd3e852SBarry Smith . -ts_rk_type 3bs 62287c2655SBarry Smith 632109b73fSDebojyoti Ghosh Level: advanced 642109b73fSDebojyoti Ghosh 6598164e13SLisandro Dalcin References: https://doi.org/10.1016/0893-9659(89)90079-7 66d1905564SLisandro Dalcin 67287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 682109b73fSDebojyoti Ghosh M*/ 692109b73fSDebojyoti Ghosh /*MC 70f68a32c8SEmil Constantinescu TSRK4 - Fourth order RK scheme. 71f68a32c8SEmil Constantinescu 722e698f8bSDebojyoti Ghosh This is the classical Runge-Kutta method with four stages. 73f68a32c8SEmil Constantinescu 74287c2655SBarry Smith Options database: 759bd3e852SBarry Smith . -ts_rk_type 4 76287c2655SBarry Smith 77f68a32c8SEmil Constantinescu Level: advanced 78f68a32c8SEmil Constantinescu 79287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 80f68a32c8SEmil Constantinescu M*/ 81f68a32c8SEmil Constantinescu /*MC 822e698f8bSDebojyoti Ghosh TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method. 83f68a32c8SEmil Constantinescu 84f68a32c8SEmil Constantinescu This method has six stages. 85f68a32c8SEmil Constantinescu 86287c2655SBarry Smith Options database: 879bd3e852SBarry Smith . -ts_rk_type 5f 88287c2655SBarry Smith 89f68a32c8SEmil Constantinescu Level: advanced 90f68a32c8SEmil Constantinescu 91287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 92f68a32c8SEmil Constantinescu M*/ 932109b73fSDebojyoti Ghosh /*MC 942e698f8bSDebojyoti Ghosh TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method. 952109b73fSDebojyoti Ghosh 96d1905564SLisandro Dalcin This method has seven stages with the First Same As Last (FSAL) property. 972109b73fSDebojyoti Ghosh 98287c2655SBarry Smith Options database: 999bd3e852SBarry Smith . -ts_rk_type 5dp 100287c2655SBarry Smith 1012109b73fSDebojyoti Ghosh Level: advanced 1022109b73fSDebojyoti Ghosh 10398164e13SLisandro Dalcin References: https://doi.org/10.1016/0771-050X(80)90013-3 104d1905564SLisandro Dalcin 105287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 1062109b73fSDebojyoti Ghosh M*/ 10705e23783SLisandro Dalcin /*MC 10805e23783SLisandro Dalcin TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method. 10905e23783SLisandro Dalcin 11005e23783SLisandro Dalcin This method has eight stages with the First Same As Last (FSAL) property. 11105e23783SLisandro Dalcin 112287c2655SBarry Smith Options database: 1139bd3e852SBarry Smith . -ts_rk_type 5bs 114287c2655SBarry Smith 11505e23783SLisandro Dalcin Level: advanced 11605e23783SLisandro Dalcin 11798164e13SLisandro Dalcin References: https://doi.org/10.1016/0898-1221(96)00141-1 11805e23783SLisandro Dalcin 119287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 12005e23783SLisandro Dalcin M*/ 121262ff7bbSBarry Smith 122f68a32c8SEmil Constantinescu /*@C 123f68a32c8SEmil Constantinescu TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK 124262ff7bbSBarry Smith 125f68a32c8SEmil Constantinescu Not Collective, but should be called by all processes which will need the schemes to be registered 126262ff7bbSBarry Smith 127f68a32c8SEmil Constantinescu Level: advanced 128262ff7bbSBarry Smith 129f68a32c8SEmil Constantinescu .keywords: TS, TSRK, register, all 130262ff7bbSBarry Smith 131f68a32c8SEmil Constantinescu .seealso: TSRKRegisterDestroy() 132262ff7bbSBarry Smith @*/ 133f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void) 134262ff7bbSBarry Smith { 1354ac538c5SBarry Smith PetscErrorCode ierr; 136262ff7bbSBarry Smith 137262ff7bbSBarry Smith PetscFunctionBegin; 138f68a32c8SEmil Constantinescu if (TSRKRegisterAllCalled) PetscFunctionReturn(0); 139f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_TRUE; 140e4dd521cSBarry Smith 1414e82626cSLisandro Dalcin #define RC PetscRealConstant 142e4dd521cSBarry Smith { 143f68a32c8SEmil Constantinescu const PetscReal 1444e82626cSLisandro Dalcin A[1][1] = {{0}}, 1454e82626cSLisandro Dalcin b[1] = {RC(1.0)}; 1463a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 147e4dd521cSBarry Smith } 148db046809SBarry Smith { 149f68a32c8SEmil Constantinescu const PetscReal 1504e82626cSLisandro Dalcin A[2][2] = {{0,0}, 1514e82626cSLisandro Dalcin {RC(1.0),0}}, 1524e82626cSLisandro Dalcin b[2] = {RC(0.5),RC(0.5)}, 1534e82626cSLisandro Dalcin bembed[2] = {RC(1.0),0}; 1543a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 155db046809SBarry Smith } 156f68a32c8SEmil Constantinescu { 157f68a32c8SEmil Constantinescu const PetscReal 15817f6384fSBarry Smith A[3][3] = {{0,0,0}, 1594e82626cSLisandro Dalcin {RC(2.0)/RC(3.0),0,0}, 1604e82626cSLisandro Dalcin {RC(-1.0)/RC(3.0),RC(1.0),0}}, 1614e82626cSLisandro Dalcin b[3] = {RC(0.25),RC(0.5),RC(0.25)}; 1623a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 163fdefd5e5SDebojyoti Ghosh } 164fdefd5e5SDebojyoti Ghosh { 165fdefd5e5SDebojyoti Ghosh const PetscReal 16617f6384fSBarry Smith A[4][4] = {{0,0,0,0}, 1674e82626cSLisandro Dalcin {RC(1.0)/RC(2.0),0,0,0}, 1684e82626cSLisandro Dalcin {0,RC(3.0)/RC(4.0),0,0}, 1694e82626cSLisandro Dalcin {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}}, 1704e82626cSLisandro Dalcin b[4] = {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}, 1714e82626cSLisandro 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)}; 1723a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 173db046809SBarry Smith } 174f68a32c8SEmil Constantinescu { 175f68a32c8SEmil Constantinescu const PetscReal 176f68a32c8SEmil Constantinescu A[4][4] = {{0,0,0,0}, 1774e82626cSLisandro Dalcin {RC(0.5),0,0,0}, 1784e82626cSLisandro Dalcin {0,RC(0.5),0,0}, 1794e82626cSLisandro Dalcin {0,0,RC(1.0),0}}, 1804e82626cSLisandro 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)}; 1813a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 182db046809SBarry Smith } 183f68a32c8SEmil Constantinescu { 184f68a32c8SEmil Constantinescu const PetscReal 18517f6384fSBarry Smith A[6][6] = {{0,0,0,0,0,0}, 1864e82626cSLisandro Dalcin {RC(0.25),0,0,0,0,0}, 1874e82626cSLisandro Dalcin {RC(3.0)/RC(32.0),RC(9.0)/RC(32.0),0,0,0,0}, 1884e82626cSLisandro Dalcin {RC(1932.0)/RC(2197.0),RC(-7200.0)/RC(2197.0),RC(7296.0)/RC(2197.0),0,0,0}, 1894e82626cSLisandro Dalcin {RC(439.0)/RC(216.0),RC(-8.0),RC(3680.0)/RC(513.0),RC(-845.0)/RC(4104.0),0,0}, 1904e82626cSLisandro 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}}, 1914e82626cSLisandro 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)}, 1924e82626cSLisandro 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}; 1933a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 194fdefd5e5SDebojyoti Ghosh } 195fdefd5e5SDebojyoti Ghosh { 196fdefd5e5SDebojyoti Ghosh const PetscReal 19717f6384fSBarry Smith A[7][7] = {{0,0,0,0,0,0,0}, 1984e82626cSLisandro Dalcin {RC(1.0)/RC(5.0),0,0,0,0,0,0}, 1994e82626cSLisandro Dalcin {RC(3.0)/RC(40.0),RC(9.0)/RC(40.0),0,0,0,0,0}, 2004e82626cSLisandro Dalcin {RC(44.0)/RC(45.0),RC(-56.0)/RC(15.0),RC(32.0)/RC(9.0),0,0,0,0}, 2014e82626cSLisandro 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}, 2024e82626cSLisandro 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}, 2034e82626cSLisandro 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}}, 2044e82626cSLisandro 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}, 205a685a763Sluzhanghpp 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)}, 206a685a763Sluzhanghpp 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)}, 207a685a763Sluzhanghpp {0,0,0,0,0}, 208a685a763Sluzhanghpp {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)}, 209a685a763Sluzhanghpp {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)}, 210a685a763Sluzhanghpp {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)}, 211a685a763Sluzhanghpp {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)}, 212a685a763Sluzhanghpp {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)}}; 213a685a763Sluzhanghpp 214a685a763Sluzhanghpp ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,5,binterp[0]);CHKERRQ(ierr); 215f68a32c8SEmil Constantinescu } 21605e23783SLisandro Dalcin { 21705e23783SLisandro Dalcin const PetscReal 21817f6384fSBarry Smith A[8][8] = {{0,0,0,0,0,0,0,0}, 2194e82626cSLisandro Dalcin {RC(1.0)/RC(6.0),0,0,0,0,0,0,0}, 2204e82626cSLisandro Dalcin {RC(2.0)/RC(27.0),RC(4.0)/RC(27.0),0,0,0,0,0,0}, 2214e82626cSLisandro Dalcin {RC(183.0)/RC(1372.0),RC(-162.0)/RC(343.0),RC(1053.0)/RC(1372.0),0,0,0,0,0}, 2224e82626cSLisandro 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}, 2234e82626cSLisandro 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}, 2244e82626cSLisandro 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}, 2254e82626cSLisandro 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}}, 2264e82626cSLisandro 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}, 2274e82626cSLisandro 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)}; 22805e23783SLisandro Dalcin ierr = TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 22905e23783SLisandro Dalcin } 2304e82626cSLisandro Dalcin #undef RC 231db046809SBarry Smith PetscFunctionReturn(0); 232db046809SBarry Smith } 233db046809SBarry Smith 234f68a32c8SEmil Constantinescu /*@C 235f68a32c8SEmil Constantinescu TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister(). 236f68a32c8SEmil Constantinescu 237f68a32c8SEmil Constantinescu Not Collective 238f68a32c8SEmil Constantinescu 239f68a32c8SEmil Constantinescu Level: advanced 240f68a32c8SEmil Constantinescu 241f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy 242f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll() 243f68a32c8SEmil Constantinescu @*/ 244f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void) 245e4dd521cSBarry Smith { 246f68a32c8SEmil Constantinescu PetscErrorCode ierr; 247f68a32c8SEmil Constantinescu RKTableauLink link; 248f68a32c8SEmil Constantinescu 249f68a32c8SEmil Constantinescu PetscFunctionBegin; 250f68a32c8SEmil Constantinescu while ((link = RKTableauList)) { 251f68a32c8SEmil Constantinescu RKTableau t = &link->tab; 252f68a32c8SEmil Constantinescu RKTableauList = link->next; 253f68a32c8SEmil Constantinescu ierr = PetscFree3(t->A,t->b,t->c); CHKERRQ(ierr); 254f68a32c8SEmil Constantinescu ierr = PetscFree (t->bembed); CHKERRQ(ierr); 255f68a32c8SEmil Constantinescu ierr = PetscFree (t->binterp); CHKERRQ(ierr); 256f68a32c8SEmil Constantinescu ierr = PetscFree (t->name); CHKERRQ(ierr); 257f68a32c8SEmil Constantinescu ierr = PetscFree (link); CHKERRQ(ierr); 258f68a32c8SEmil Constantinescu } 259f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_FALSE; 260f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 261f68a32c8SEmil Constantinescu } 262f68a32c8SEmil Constantinescu 263f68a32c8SEmil Constantinescu /*@C 264f68a32c8SEmil Constantinescu TSRKInitializePackage - This function initializes everything in the TSRK package. It is called 2658a690491SBarry Smith from TSInitializePackage(). 266f68a32c8SEmil Constantinescu 267f68a32c8SEmil Constantinescu Level: developer 268f68a32c8SEmil Constantinescu 269f68a32c8SEmil Constantinescu .keywords: TS, TSRK, initialize, package 270f68a32c8SEmil Constantinescu .seealso: PetscInitialize() 271f68a32c8SEmil Constantinescu @*/ 272f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void) 273f68a32c8SEmil Constantinescu { 274186e87acSLisandro Dalcin PetscErrorCode ierr; 275e4dd521cSBarry Smith 276e4dd521cSBarry Smith PetscFunctionBegin; 277f68a32c8SEmil Constantinescu if (TSRKPackageInitialized) PetscFunctionReturn(0); 278f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_TRUE; 279f68a32c8SEmil Constantinescu ierr = TSRKRegisterAll();CHKERRQ(ierr); 280f68a32c8SEmil Constantinescu ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr); 281f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 282f68a32c8SEmil Constantinescu } 283186e87acSLisandro Dalcin 284f68a32c8SEmil Constantinescu /*@C 285f68a32c8SEmil Constantinescu TSRKFinalizePackage - This function destroys everything in the TSRK package. It is 286f68a32c8SEmil Constantinescu called from PetscFinalize(). 287186e87acSLisandro Dalcin 288f68a32c8SEmil Constantinescu Level: developer 289f68a32c8SEmil Constantinescu 290f68a32c8SEmil Constantinescu .keywords: Petsc, destroy, package 291f68a32c8SEmil Constantinescu .seealso: PetscFinalize() 292f68a32c8SEmil Constantinescu @*/ 293f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void) 294f68a32c8SEmil Constantinescu { 295f68a32c8SEmil Constantinescu PetscErrorCode ierr; 296f68a32c8SEmil Constantinescu 297f68a32c8SEmil Constantinescu PetscFunctionBegin; 298f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_FALSE; 299f68a32c8SEmil Constantinescu ierr = TSRKRegisterDestroy();CHKERRQ(ierr); 300f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 301f68a32c8SEmil Constantinescu } 302f68a32c8SEmil Constantinescu 303f68a32c8SEmil Constantinescu /*@C 304f68a32c8SEmil Constantinescu TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 305f68a32c8SEmil Constantinescu 306f68a32c8SEmil Constantinescu Not Collective, but the same schemes should be registered on all processes on which they will be used 307f68a32c8SEmil Constantinescu 308f68a32c8SEmil Constantinescu Input Parameters: 309f68a32c8SEmil Constantinescu + name - identifier for method 310f68a32c8SEmil Constantinescu . order - approximation order of method 311f68a32c8SEmil Constantinescu . s - number of stages, this is the dimension of the matrices below 312f68a32c8SEmil Constantinescu . A - stage coefficients (dimension s*s, row-major) 313f68a32c8SEmil Constantinescu . b - step completion table (dimension s; NULL to use last row of A) 314f68a32c8SEmil Constantinescu . c - abscissa (dimension s; NULL to use row sums of A) 315f68a32c8SEmil Constantinescu . bembed - completion table for embedded method (dimension s; NULL if not available) 3163a8a9803SLisandro Dalcin . p - Order of the interpolation scheme, equal to the number of columns of binterp 3173a8a9803SLisandro Dalcin - binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1) 318f68a32c8SEmil Constantinescu 319f68a32c8SEmil Constantinescu Notes: 320f68a32c8SEmil Constantinescu Several RK methods are provided, this function is only needed to create new methods. 321f68a32c8SEmil Constantinescu 322f68a32c8SEmil Constantinescu Level: advanced 323f68a32c8SEmil Constantinescu 324f68a32c8SEmil Constantinescu .keywords: TS, register 325f68a32c8SEmil Constantinescu 326f68a32c8SEmil Constantinescu .seealso: TSRK 327f68a32c8SEmil Constantinescu @*/ 328f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s, 329f68a32c8SEmil Constantinescu const PetscReal A[],const PetscReal b[],const PetscReal c[], 3303a8a9803SLisandro Dalcin const PetscReal bembed[],PetscInt p,const PetscReal binterp[]) 331f68a32c8SEmil Constantinescu { 332f68a32c8SEmil Constantinescu PetscErrorCode ierr; 333f68a32c8SEmil Constantinescu RKTableauLink link; 334f68a32c8SEmil Constantinescu RKTableau t; 335f68a32c8SEmil Constantinescu PetscInt i,j; 336f68a32c8SEmil Constantinescu 337f68a32c8SEmil Constantinescu PetscFunctionBegin; 3383a8a9803SLisandro Dalcin PetscValidCharPointer(name,1); 3393a8a9803SLisandro Dalcin PetscValidRealPointer(A,4); 3403a8a9803SLisandro Dalcin if (b) PetscValidRealPointer(b,5); 3413a8a9803SLisandro Dalcin if (c) PetscValidRealPointer(c,6); 3423a8a9803SLisandro Dalcin if (bembed) PetscValidRealPointer(bembed,7); 3433a8a9803SLisandro Dalcin if (binterp || p > 1) PetscValidRealPointer(binterp,9); 3443a8a9803SLisandro Dalcin 3451d36bdfdSBarry Smith ierr = TSRKInitializePackage();CHKERRQ(ierr); 34695dccacaSBarry Smith ierr = PetscNew(&link);CHKERRQ(ierr); 347f68a32c8SEmil Constantinescu t = &link->tab; 3483a8a9803SLisandro Dalcin 349f68a32c8SEmil Constantinescu ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 350f68a32c8SEmil Constantinescu t->order = order; 351f68a32c8SEmil Constantinescu t->s = s; 352dcca6d9dSJed Brown ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr); 353f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 354f68a32c8SEmil Constantinescu if (b) { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); } 355f68a32c8SEmil Constantinescu else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i]; 356f68a32c8SEmil Constantinescu if (c) { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); } 357f68a32c8SEmil 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]; 358d760c35bSDebojyoti Ghosh t->FSAL = PETSC_TRUE; 359d760c35bSDebojyoti Ghosh for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE; 3603a8a9803SLisandro Dalcin 361f68a32c8SEmil Constantinescu if (bembed) { 362785e854fSJed Brown ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr); 363f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr); 364f68a32c8SEmil Constantinescu } 365f68a32c8SEmil Constantinescu 3663a8a9803SLisandro Dalcin if (!binterp) { p = 1; binterp = t->b; } 3673a8a9803SLisandro Dalcin t->p = p; 3683a8a9803SLisandro Dalcin ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr); 3693a8a9803SLisandro Dalcin ierr = PetscMemcpy(t->binterp,binterp,s*p*sizeof(binterp[0]));CHKERRQ(ierr); 3703a8a9803SLisandro Dalcin 371f68a32c8SEmil Constantinescu link->next = RKTableauList; 372f68a32c8SEmil Constantinescu RKTableauList = link; 373f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 374f68a32c8SEmil Constantinescu } 375f68a32c8SEmil Constantinescu 376e4dd521cSBarry Smith /* 3772c9cb062Sluzhanghpp This is for single-step RK method 378f68a32c8SEmil Constantinescu The step completion formula is 379e4dd521cSBarry Smith 380f68a32c8SEmil Constantinescu x1 = x0 + h b^T YdotRHS 381f68a32c8SEmil Constantinescu 382f68a32c8SEmil Constantinescu This function can be called before or after ts->vec_sol has been updated. 383f68a32c8SEmil Constantinescu Suppose we have a completion formula (b) and an embedded formula (be) of different order. 384f68a32c8SEmil Constantinescu We can write 385f68a32c8SEmil Constantinescu 386f68a32c8SEmil Constantinescu x1e = x0 + h be^T YdotRHS 387f68a32c8SEmil Constantinescu = x1 - h b^T YdotRHS + h be^T YdotRHS 388f68a32c8SEmil Constantinescu = x1 + h (be - b)^T YdotRHS 389f68a32c8SEmil Constantinescu 390f68a32c8SEmil Constantinescu so we can evaluate the method with different order even after the step has been optimistically completed. 391e4dd521cSBarry Smith */ 392f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done) 393f68a32c8SEmil Constantinescu { 394f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 395f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 396f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 397f68a32c8SEmil Constantinescu PetscReal h; 398f68a32c8SEmil Constantinescu PetscInt s = tab->s,j; 399f68a32c8SEmil Constantinescu PetscErrorCode ierr; 400f68a32c8SEmil Constantinescu 401f68a32c8SEmil Constantinescu PetscFunctionBegin; 402f68a32c8SEmil Constantinescu switch (rk->status) { 403f68a32c8SEmil Constantinescu case TS_STEP_INCOMPLETE: 404f68a32c8SEmil Constantinescu case TS_STEP_PENDING: 405f68a32c8SEmil Constantinescu h = ts->time_step; break; 406f68a32c8SEmil Constantinescu case TS_STEP_COMPLETE: 407be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 408f68a32c8SEmil Constantinescu default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 409f68a32c8SEmil Constantinescu } 410f68a32c8SEmil Constantinescu if (order == tab->order) { 411f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { 412f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 41390b67129SHong Zhang for (j=0; j<s; j++) w[j] = h*tab->b[j]/rk->dtratio; 414f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 415f68a32c8SEmil Constantinescu } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 416f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 417f68a32c8SEmil Constantinescu } else if (order == tab->order-1) { 418f68a32c8SEmil Constantinescu if (!tab->bembed) goto unavailable; 419f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { /*Complete with the embedded method (be)*/ 420f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 421f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 422f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 423f68a32c8SEmil Constantinescu } else { /*Rollback and re-complete using (be-b) */ 424f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 425f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 426f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 427f68a32c8SEmil Constantinescu } 428f68a32c8SEmil Constantinescu if (done) *done = PETSC_TRUE; 429f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 430f68a32c8SEmil Constantinescu } 431f68a32c8SEmil Constantinescu unavailable: 432f68a32c8SEmil Constantinescu if (done) *done = PETSC_FALSE; 433a7fac7c2SEmil 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); 434f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 435f68a32c8SEmil Constantinescu } 436f68a32c8SEmil Constantinescu 4370f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts) 4380f7a1166SHong Zhang { 4390f7a1166SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 4400f7a1166SHong Zhang RKTableau tab = rk->tableau; 4410f7a1166SHong Zhang const PetscInt s = tab->s; 4420f7a1166SHong Zhang const PetscReal *b = tab->b,*c = tab->c; 4430f7a1166SHong Zhang Vec *Y = rk->Y; 4440f7a1166SHong Zhang PetscInt i; 4450f7a1166SHong Zhang PetscErrorCode ierr; 4460f7a1166SHong Zhang 4470f7a1166SHong Zhang PetscFunctionBegin; 4480f7a1166SHong Zhang /* backup cost integral */ 4490f7a1166SHong Zhang ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr); 4500f7a1166SHong Zhang for (i=s-1; i>=0; i--) { 4510f7a1166SHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 45251a7f651SHong Zhang ierr = TSComputeCostIntegrand(ts,rk->ptime+rk->time_step*c[i],Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 4530f7a1166SHong Zhang ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 4540f7a1166SHong Zhang } 4550f7a1166SHong Zhang PetscFunctionReturn(0); 4560f7a1166SHong Zhang } 4570f7a1166SHong Zhang 4580f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts) 4590f7a1166SHong Zhang { 4600f7a1166SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 4610f7a1166SHong Zhang RKTableau tab = rk->tableau; 4620f7a1166SHong Zhang const PetscInt s = tab->s; 4630f7a1166SHong Zhang const PetscReal *b = tab->b,*c = tab->c; 4640f7a1166SHong Zhang Vec *Y = rk->Y; 4650f7a1166SHong Zhang PetscInt i; 4660f7a1166SHong Zhang PetscErrorCode ierr; 4670f7a1166SHong Zhang 4680f7a1166SHong Zhang PetscFunctionBegin; 4690f7a1166SHong Zhang for (i=s-1; i>=0; i--) { 4700f7a1166SHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 47151a7f651SHong Zhang ierr = TSComputeCostIntegrand(ts,ts->ptime+ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 4720f7a1166SHong Zhang ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 4730f7a1166SHong Zhang } 4740f7a1166SHong Zhang PetscFunctionReturn(0); 4750f7a1166SHong Zhang } 4760f7a1166SHong Zhang 477fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts) 478fecfb714SLisandro Dalcin { 479fecfb714SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 480fecfb714SLisandro Dalcin RKTableau tab = rk->tableau; 481fecfb714SLisandro Dalcin const PetscInt s = tab->s; 482fecfb714SLisandro Dalcin const PetscReal *b = tab->b; 483fecfb714SLisandro Dalcin PetscScalar *w = rk->work; 484fecfb714SLisandro Dalcin Vec *YdotRHS = rk->YdotRHS; 485fecfb714SLisandro Dalcin PetscInt j; 486fecfb714SLisandro Dalcin PetscReal h; 487fecfb714SLisandro Dalcin PetscErrorCode ierr; 488fecfb714SLisandro Dalcin 489fecfb714SLisandro Dalcin PetscFunctionBegin; 490fecfb714SLisandro Dalcin switch (rk->status) { 491fecfb714SLisandro Dalcin case TS_STEP_INCOMPLETE: 492fecfb714SLisandro Dalcin case TS_STEP_PENDING: 493fecfb714SLisandro Dalcin h = ts->time_step; break; 494fecfb714SLisandro Dalcin case TS_STEP_COMPLETE: 495fecfb714SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 496fecfb714SLisandro Dalcin default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 497fecfb714SLisandro Dalcin } 498fecfb714SLisandro Dalcin for (j=0; j<s; j++) w[j] = -h*b[j]; 499fecfb714SLisandro Dalcin ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr); 500fecfb714SLisandro Dalcin PetscFunctionReturn(0); 501fecfb714SLisandro Dalcin } 502fecfb714SLisandro Dalcin 503f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts) 504f68a32c8SEmil Constantinescu { 505f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 506f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 507f68a32c8SEmil Constantinescu const PetscInt s = tab->s; 508fecfb714SLisandro Dalcin const PetscReal *A = tab->A,*c = tab->c; 509f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 510f68a32c8SEmil Constantinescu Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 511d1905564SLisandro Dalcin PetscBool FSAL = tab->FSAL; 512f68a32c8SEmil Constantinescu TSAdapt adapt; 513fecfb714SLisandro Dalcin PetscInt i,j; 514be5899b3SLisandro Dalcin PetscInt rejections = 0; 515be5899b3SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 516be5899b3SLisandro Dalcin PetscReal next_time_step = ts->time_step; 517f68a32c8SEmil Constantinescu PetscErrorCode ierr; 518f68a32c8SEmil Constantinescu 519f68a32c8SEmil Constantinescu PetscFunctionBegin; 520d1905564SLisandro Dalcin if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE; 521d1905564SLisandro Dalcin if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); } 522d1905564SLisandro Dalcin 523f68a32c8SEmil Constantinescu rk->status = TS_STEP_INCOMPLETE; 524be5899b3SLisandro Dalcin while (!ts->reason && rk->status != TS_STEP_COMPLETE) { 525be5899b3SLisandro Dalcin PetscReal t = ts->ptime; 526f68a32c8SEmil Constantinescu PetscReal h = ts->time_step; 527f68a32c8SEmil Constantinescu for (i=0; i<s; i++) { 5289be3e283SDebojyoti Ghosh rk->stage_time = t + h*c[i]; 5299be3e283SDebojyoti Ghosh ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr); 530f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 531f68a32c8SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 532f68a32c8SEmil Constantinescu ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 5339be3e283SDebojyoti Ghosh ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 534f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 535be5899b3SLisandro Dalcin ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr); 536fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 5378f738a7cSLisandro Dalcin if (FSAL && !i) continue; 538f68a32c8SEmil Constantinescu ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 539f68a32c8SEmil Constantinescu } 540be5899b3SLisandro Dalcin 541be5899b3SLisandro Dalcin rk->status = TS_STEP_INCOMPLETE; 542f68a32c8SEmil Constantinescu ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 543f68a32c8SEmil Constantinescu rk->status = TS_STEP_PENDING; 544f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 545f68a32c8SEmil Constantinescu ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 5461917a363SLisandro Dalcin ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr); 547fecfb714SLisandro Dalcin ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 548be5899b3SLisandro Dalcin rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 549be5899b3SLisandro Dalcin if (!accept) { /* Roll back the current step */ 550fecfb714SLisandro Dalcin ierr = TSRollBack_RK(ts);CHKERRQ(ierr); 551be5899b3SLisandro Dalcin ts->time_step = next_time_step; 552be5899b3SLisandro Dalcin goto reject_step; 553be5899b3SLisandro Dalcin } 554be5899b3SLisandro Dalcin 5550f7a1166SHong Zhang if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 5560f7a1166SHong Zhang rk->ptime = ts->ptime; 5570f7a1166SHong Zhang rk->time_step = ts->time_step; 558493ed95fSHong Zhang } 559be5899b3SLisandro Dalcin 560f68a32c8SEmil Constantinescu ts->ptime += ts->time_step; 561f68a32c8SEmil Constantinescu ts->time_step = next_time_step; 562f68a32c8SEmil Constantinescu break; 563be5899b3SLisandro Dalcin 564be5899b3SLisandro Dalcin reject_step: 565fecfb714SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 566be5899b3SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 567be5899b3SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 568be5899b3SLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 569f68a32c8SEmil Constantinescu } 570f68a32c8SEmil Constantinescu } 571f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 572e4dd521cSBarry Smith } 573e4dd521cSBarry Smith 574f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts) 575f6a906c0SBarry Smith { 576f6a906c0SBarry Smith TS_RK *rk = (TS_RK*)ts->data; 577be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 578be5899b3SLisandro Dalcin PetscInt s = tab->s; 579f6a906c0SBarry Smith PetscErrorCode ierr; 580f6a906c0SBarry Smith 581f6a906c0SBarry Smith PetscFunctionBegin; 582f6a906c0SBarry Smith if (ts->adjointsetupcalled++) PetscFunctionReturn(0); 583abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr); 584f6a906c0SBarry Smith if(ts->vecs_sensip) { 585abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr); 586f6a906c0SBarry Smith } 587f6a906c0SBarry Smith PetscFunctionReturn(0); 588f6a906c0SBarry Smith } 589f6a906c0SBarry Smith 59042f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts) 591d2daff3dSHong Zhang { 592c235aa8dSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 593c235aa8dSHong Zhang RKTableau tab = rk->tableau; 594c235aa8dSHong Zhang const PetscInt s = tab->s; 595c235aa8dSHong Zhang const PetscReal *A = tab->A,*b = tab->b,*c = tab->c; 596c235aa8dSHong Zhang PetscScalar *w = rk->work; 5973389c7d9SStefano Zampini Vec *Y = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu; 598b18ea86cSHong Zhang PetscInt i,j,nadj; 5993d3eaba7SBarry Smith PetscReal t = ts->ptime; 600d2daff3dSHong Zhang PetscErrorCode ierr; 601cef76868SBarry Smith PetscReal h = ts->time_step; 602c235aa8dSHong Zhang 603d2daff3dSHong Zhang PetscFunctionBegin; 604c235aa8dSHong Zhang rk->status = TS_STEP_INCOMPLETE; 605c235aa8dSHong Zhang for (i=s-1; i>=0; i--) { 6063389c7d9SStefano Zampini Mat J; 6073389c7d9SStefano Zampini PetscReal stage_time = t + h*(1.0-c[i]); 6083389c7d9SStefano Zampini PetscBool zero = PETSC_FALSE; 6093389c7d9SStefano Zampini 6103389c7d9SStefano Zampini ierr = TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr); 6113389c7d9SStefano Zampini ierr = TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);CHKERRQ(ierr); 612493ed95fSHong Zhang if (ts->vec_costintegral) { 6133389c7d9SStefano Zampini ierr = TSAdjointComputeDRDYFunction(ts,stage_time,Y[i],ts->vecs_drdy);CHKERRQ(ierr); 614493ed95fSHong Zhang } 6153389c7d9SStefano Zampini /* Stage values of mu */ 6163389c7d9SStefano Zampini if (ts->vecs_sensip) { 6173389c7d9SStefano Zampini ierr = TSAdjointComputeRHSJacobian(ts,stage_time,Y[i],ts->Jacp);CHKERRQ(ierr); 6183389c7d9SStefano Zampini if (ts->vec_costintegral) { 6193389c7d9SStefano Zampini ierr = TSAdjointComputeDRDPFunction(ts,stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr); 6203389c7d9SStefano Zampini } 6213389c7d9SStefano Zampini } 6223389c7d9SStefano Zampini 6233389c7d9SStefano Zampini if (b[i] == 0 && i == s-1) zero = PETSC_TRUE; 624abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 6253389c7d9SStefano Zampini DM dm; 6263389c7d9SStefano Zampini Vec VecSensiTemp; 6273389c7d9SStefano Zampini 6283389c7d9SStefano Zampini ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 6293389c7d9SStefano Zampini ierr = DMGetGlobalVector(dm,&VecSensiTemp);CHKERRQ(ierr); 6303389c7d9SStefano Zampini /* Stage values of lambda */ 6313389c7d9SStefano Zampini if (!zero) { 6323389c7d9SStefano Zampini ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp);CHKERRQ(ierr); 6333389c7d9SStefano Zampini ierr = VecScale(VecSensiTemp,-h*b[i]);CHKERRQ(ierr); 6343389c7d9SStefano Zampini for (j=i+1; j<s; j++) w[j-i-1] = -h*A[j*s+i]; 6353389c7d9SStefano Zampini ierr = VecMAXPY(VecSensiTemp,s-i-1,w,&VecDeltaLam[nadj*s+i+1]);CHKERRQ(ierr); 6363389c7d9SStefano Zampini ierr = MatMultTranspose(J,VecSensiTemp,VecDeltaLam[nadj*s+i]);CHKERRQ(ierr); 6373389c7d9SStefano Zampini } else { 6383389c7d9SStefano Zampini ierr = VecSet(VecDeltaLam[nadj*s+i],0);CHKERRQ(ierr); 6393389c7d9SStefano Zampini } 640493ed95fSHong Zhang if (ts->vec_costintegral) { 641493ed95fSHong Zhang ierr = VecAXPY(VecDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdy[nadj]);CHKERRQ(ierr); 642493ed95fSHong Zhang } 6436497ce14SHong Zhang 644ad8e2604SHong Zhang /* Stage values of mu */ 6456497ce14SHong Zhang if (ts->vecs_sensip) { 6463389c7d9SStefano Zampini if (!zero) { 6473389c7d9SStefano Zampini ierr = MatMultTranspose(ts->Jacp,VecSensiTemp,VecDeltaMu[nadj*s+i]);CHKERRQ(ierr); 6483389c7d9SStefano Zampini } else { 6493389c7d9SStefano Zampini ierr = VecSet(VecDeltaMu[nadj*s+i],0);CHKERRQ(ierr); 650493ed95fSHong Zhang } 651493ed95fSHong Zhang if (ts->vec_costintegral) { 652493ed95fSHong Zhang ierr = VecAXPY(VecDeltaMu[nadj*s+i],-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr); 653493ed95fSHong Zhang } 654ad8e2604SHong Zhang } 6553389c7d9SStefano Zampini ierr = DMRestoreGlobalVector(dm,&VecSensiTemp);CHKERRQ(ierr); 656c235aa8dSHong Zhang } 6576497ce14SHong Zhang } 658c235aa8dSHong Zhang 659c235aa8dSHong Zhang for (j=0; j<s; j++) w[j] = 1.0; 660abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 661b18ea86cSHong Zhang ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr); 6626497ce14SHong Zhang if (ts->vecs_sensip) { 663ad8e2604SHong Zhang ierr = VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);CHKERRQ(ierr); 664b18ea86cSHong Zhang } 6656497ce14SHong Zhang } 666c235aa8dSHong Zhang rk->status = TS_STEP_COMPLETE; 667d2daff3dSHong Zhang PetscFunctionReturn(0); 668d2daff3dSHong Zhang } 669d2daff3dSHong Zhang 67055de54fcSHong Zhang static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X) 67155de54fcSHong Zhang { 67255de54fcSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 67355de54fcSHong Zhang PetscInt s = rk->tableau->s,p = rk->tableau->p,i,j; 67455de54fcSHong Zhang PetscReal h; 67555de54fcSHong Zhang PetscReal tt,t; 67655de54fcSHong Zhang PetscScalar *b; 67755de54fcSHong Zhang const PetscReal *B = rk->tableau->binterp; 67855de54fcSHong Zhang PetscErrorCode ierr; 67955de54fcSHong Zhang 68055de54fcSHong Zhang PetscFunctionBegin; 68155de54fcSHong Zhang if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 68255de54fcSHong Zhang 68355de54fcSHong Zhang switch (rk->status) { 68455de54fcSHong Zhang case TS_STEP_INCOMPLETE: 68555de54fcSHong Zhang case TS_STEP_PENDING: 68655de54fcSHong Zhang h = ts->time_step; 68755de54fcSHong Zhang t = (itime - ts->ptime)/h; 68855de54fcSHong Zhang break; 68955de54fcSHong Zhang case TS_STEP_COMPLETE: 69055de54fcSHong Zhang h = ts->ptime - ts->ptime_prev; 69155de54fcSHong Zhang t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 69255de54fcSHong Zhang break; 69355de54fcSHong Zhang default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 69455de54fcSHong Zhang } 69555de54fcSHong Zhang ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 69655de54fcSHong Zhang for (i=0; i<s; i++) b[i] = 0; 69755de54fcSHong Zhang for (j=0,tt=t; j<p; j++,tt*=t) { 69855de54fcSHong Zhang for (i=0; i<s; i++) { 69955de54fcSHong Zhang b[i] += h * B[i*p+j] * tt; 70055de54fcSHong Zhang } 70155de54fcSHong Zhang } 70255de54fcSHong Zhang ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr); 70355de54fcSHong Zhang ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr); 70455de54fcSHong Zhang ierr = PetscFree(b);CHKERRQ(ierr); 70555de54fcSHong Zhang PetscFunctionReturn(0); 70655de54fcSHong Zhang } 70755de54fcSHong Zhang 70855de54fcSHong Zhang /*------------------------------------------------------------*/ 70955de54fcSHong Zhang 710be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts) 711be5899b3SLisandro Dalcin { 712be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 713be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 714be5899b3SLisandro Dalcin PetscErrorCode ierr; 715be5899b3SLisandro Dalcin 716be5899b3SLisandro Dalcin PetscFunctionBegin; 717be5899b3SLisandro Dalcin if (!tab) PetscFunctionReturn(0); 718be5899b3SLisandro Dalcin ierr = PetscFree(rk->work);CHKERRQ(ierr); 719be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr); 720be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr); 721be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr); 722be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr); 723be5899b3SLisandro Dalcin PetscFunctionReturn(0); 724be5899b3SLisandro Dalcin } 725be5899b3SLisandro Dalcin 726277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts) 727e4dd521cSBarry Smith { 7285f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 7296849ba73SBarry Smith PetscErrorCode ierr; 730e4dd521cSBarry Smith 731e4dd521cSBarry Smith PetscFunctionBegin; 732be5899b3SLisandro Dalcin ierr = TSRKTableauReset(ts);CHKERRQ(ierr); 7330f7a1166SHong Zhang ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr); 734*0fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 73555de54fcSHong Zhang ierr = PetscTryMethod(ts,"TSReset_MRKSPLIT_C",(TS),(ts));CHKERRQ(ierr); 736*0fe4d17eSHong Zhang } else { 73755de54fcSHong Zhang ierr = PetscTryMethod(ts,"TSReset_MRKNONSPLIT_C",(TS),(ts));CHKERRQ(ierr); 738*0fe4d17eSHong Zhang } 739277b19d0SLisandro Dalcin PetscFunctionReturn(0); 740e4dd521cSBarry Smith } 741277b19d0SLisandro Dalcin 742f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx) 743f68a32c8SEmil Constantinescu { 744f68a32c8SEmil Constantinescu PetscFunctionBegin; 745f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 746f68a32c8SEmil Constantinescu } 747f68a32c8SEmil Constantinescu 748f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 749f68a32c8SEmil Constantinescu { 750f68a32c8SEmil Constantinescu PetscFunctionBegin; 751f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 752f68a32c8SEmil Constantinescu } 753f68a32c8SEmil Constantinescu 754f68a32c8SEmil Constantinescu 755f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx) 756f68a32c8SEmil Constantinescu { 757f68a32c8SEmil Constantinescu PetscFunctionBegin; 758f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 759f68a32c8SEmil Constantinescu } 760f68a32c8SEmil Constantinescu 761f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 762f68a32c8SEmil Constantinescu { 763f68a32c8SEmil Constantinescu 764f68a32c8SEmil Constantinescu PetscFunctionBegin; 765f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 766f68a32c8SEmil Constantinescu } 767c235aa8dSHong Zhang /* 768d2daff3dSHong Zhang static PetscErrorCode RKSetAdjCoe(RKTableau tab) 769d2daff3dSHong Zhang { 770d2daff3dSHong Zhang PetscReal *A,*b; 771d2daff3dSHong Zhang PetscInt s,i,j; 772d2daff3dSHong Zhang PetscErrorCode ierr; 773d2daff3dSHong Zhang 774d2daff3dSHong Zhang PetscFunctionBegin; 775d2daff3dSHong Zhang s = tab->s; 776d2daff3dSHong Zhang ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr); 777d2daff3dSHong Zhang 778d2daff3dSHong Zhang for (i=0; i<s; i++) 779d2daff3dSHong Zhang for (j=0; j<s; j++) { 780d2daff3dSHong 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]; 781d2daff3dSHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr); 782d2daff3dSHong Zhang } 783d2daff3dSHong Zhang 784d2daff3dSHong Zhang for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i]; 785d2daff3dSHong Zhang 786d2daff3dSHong Zhang ierr = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 787d2daff3dSHong Zhang ierr = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr); 788d2daff3dSHong Zhang ierr = PetscFree2(A,b);CHKERRQ(ierr); 789d2daff3dSHong Zhang PetscFunctionReturn(0); 790d2daff3dSHong Zhang } 791c235aa8dSHong Zhang */ 792be5899b3SLisandro Dalcin 793be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts) 794be5899b3SLisandro Dalcin { 795be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 796be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 797be5899b3SLisandro Dalcin PetscErrorCode ierr; 798be5899b3SLisandro Dalcin 799be5899b3SLisandro Dalcin PetscFunctionBegin; 800be5899b3SLisandro Dalcin ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr); 801be5899b3SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr); 802be5899b3SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr); 803be5899b3SLisandro Dalcin PetscFunctionReturn(0); 804be5899b3SLisandro Dalcin } 805be5899b3SLisandro Dalcin 806f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts) 807f68a32c8SEmil Constantinescu { 808f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 809f68a32c8SEmil Constantinescu PetscErrorCode ierr; 810f68a32c8SEmil Constantinescu DM dm; 811f68a32c8SEmil Constantinescu 812f68a32c8SEmil Constantinescu PetscFunctionBegin; 8133dd83b38SBarry Smith ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr); 814be5899b3SLisandro Dalcin ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr); 8150f7a1166SHong Zhang if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 8160f7a1166SHong Zhang ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr); 8170f7a1166SHong Zhang } 818f68a32c8SEmil Constantinescu ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 819f68a32c8SEmil Constantinescu ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 820f68a32c8SEmil Constantinescu ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 821*0fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 82255de54fcSHong Zhang ierr = PetscTryMethod(ts,"TSSetUp_MRKSPLIT_C",(TS),(ts));CHKERRQ(ierr); 823*0fe4d17eSHong Zhang } else { 82455de54fcSHong Zhang ierr = PetscTryMethod(ts,"TSSetUp_MRKNONSPLIT_C",(TS),(ts));CHKERRQ(ierr); 825*0fe4d17eSHong Zhang } 826e4dd521cSBarry Smith PetscFunctionReturn(0); 827e4dd521cSBarry Smith } 828d2daff3dSHong Zhang 8294416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts) 830e4dd521cSBarry Smith { 831be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 832dfbe8321SBarry Smith PetscErrorCode ierr; 833262ff7bbSBarry Smith 834e4dd521cSBarry Smith PetscFunctionBegin; 835e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr); 836f68a32c8SEmil Constantinescu { 837f68a32c8SEmil Constantinescu RKTableauLink link; 838f68a32c8SEmil Constantinescu PetscInt count,choice; 839*0fe4d17eSHong Zhang PetscBool flg,use_multirate = PETSC_FALSE; 840f68a32c8SEmil Constantinescu const char **namelist; 8412c9cb062Sluzhanghpp 842f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) ; 843f489ac74SBarry Smith ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr); 844f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 845*0fe4d17eSHong Zhang ierr = PetscOptionsBool("-ts_rk_multirate","Use interpolation-based multirate RK method","TSRKSetMultirate",rk->use_multirate,&use_multirate,&flg);CHKERRQ(ierr); 846*0fe4d17eSHong Zhang if (flg) { 847*0fe4d17eSHong Zhang ierr = TSRKSetMultirate(ts,use_multirate);CHKERRQ(ierr); 848*0fe4d17eSHong Zhang } 849be5899b3SLisandro Dalcin ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr); 850be5899b3SLisandro Dalcin if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);} 851f68a32c8SEmil Constantinescu ierr = PetscFree(namelist);CHKERRQ(ierr); 852f68a32c8SEmil Constantinescu } 853262ff7bbSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 8542c9cb062Sluzhanghpp ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Multirate methods options","");CHKERRQ(ierr); 8552c9cb062Sluzhanghpp ierr = PetscOptionsInt("-ts_rk_dtratio","time step ratio between slow and fast","",rk->dtratio,&rk->dtratio,NULL);CHKERRQ(ierr); 8562c9cb062Sluzhanghpp ierr = PetscOptionsEnd();CHKERRQ(ierr); 857e4dd521cSBarry Smith PetscFunctionReturn(0); 858e4dd521cSBarry Smith } 859e4dd521cSBarry Smith 8605f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer) 861e4dd521cSBarry Smith { 8625f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 8638370ee3bSLisandro Dalcin PetscBool iascii; 864dfbe8321SBarry Smith PetscErrorCode ierr; 8652cdf8120Sasbjorn 866e4dd521cSBarry Smith PetscFunctionBegin; 867251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 8688370ee3bSLisandro Dalcin if (iascii) { 8699c334d8fSLisandro Dalcin RKTableau tab = rk->tableau; 870f68a32c8SEmil Constantinescu TSRKType rktype; 871f68a32c8SEmil Constantinescu char buf[512]; 872f68a32c8SEmil Constantinescu ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr); 873efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RK type %s\n",rktype);CHKERRQ(ierr); 8740c37eb8eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," Order: %D\n",tab->order);CHKERRQ(ierr); 8750c37eb8eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr); 876f68a32c8SEmil Constantinescu ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 877f68a32c8SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," Abscissa c = %s\n",buf);CHKERRQ(ierr); 8788370ee3bSLisandro Dalcin } 879f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 880f68a32c8SEmil Constantinescu } 881f68a32c8SEmil Constantinescu 882f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer) 883f68a32c8SEmil Constantinescu { 884f68a32c8SEmil Constantinescu PetscErrorCode ierr; 8859c334d8fSLisandro Dalcin TSAdapt adapt; 886f68a32c8SEmil Constantinescu 887f68a32c8SEmil Constantinescu PetscFunctionBegin; 8889c334d8fSLisandro Dalcin ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 8899c334d8fSLisandro Dalcin ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr); 890f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 891f68a32c8SEmil Constantinescu } 892f68a32c8SEmil Constantinescu 893f68a32c8SEmil Constantinescu /*@C 894f68a32c8SEmil Constantinescu TSRKSetType - Set the type of RK scheme 895f68a32c8SEmil Constantinescu 896f68a32c8SEmil Constantinescu Logically collective 897f68a32c8SEmil Constantinescu 898f68a32c8SEmil Constantinescu Input Parameter: 899f68a32c8SEmil Constantinescu + ts - timestepping context 900f68a32c8SEmil Constantinescu - rktype - type of RK-scheme 901f68a32c8SEmil Constantinescu 902287c2655SBarry Smith Options Database: 9039bd3e852SBarry Smith . -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs> 904287c2655SBarry Smith 905f68a32c8SEmil Constantinescu Level: intermediate 906f68a32c8SEmil Constantinescu 907287c2655SBarry Smith .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS 908f68a32c8SEmil Constantinescu @*/ 909f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype) 910f68a32c8SEmil Constantinescu { 911f68a32c8SEmil Constantinescu PetscErrorCode ierr; 912f68a32c8SEmil Constantinescu 913f68a32c8SEmil Constantinescu PetscFunctionBegin; 914f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 915b92453a8SLisandro Dalcin PetscValidCharPointer(rktype,2); 916f68a32c8SEmil Constantinescu ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr); 917f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 918f68a32c8SEmil Constantinescu } 919f68a32c8SEmil Constantinescu 920f68a32c8SEmil Constantinescu /*@C 921f68a32c8SEmil Constantinescu TSRKGetType - Get the type of RK scheme 922f68a32c8SEmil Constantinescu 923f68a32c8SEmil Constantinescu Logically collective 924f68a32c8SEmil Constantinescu 925f68a32c8SEmil Constantinescu Input Parameter: 926f68a32c8SEmil Constantinescu . ts - timestepping context 927f68a32c8SEmil Constantinescu 928f68a32c8SEmil Constantinescu Output Parameter: 929f68a32c8SEmil Constantinescu . rktype - type of RK-scheme 930f68a32c8SEmil Constantinescu 931f68a32c8SEmil Constantinescu Level: intermediate 932f68a32c8SEmil Constantinescu 933f68a32c8SEmil Constantinescu .seealso: TSRKGetType() 934f68a32c8SEmil Constantinescu @*/ 935f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype) 936f68a32c8SEmil Constantinescu { 937f68a32c8SEmil Constantinescu PetscErrorCode ierr; 938f68a32c8SEmil Constantinescu 939f68a32c8SEmil Constantinescu PetscFunctionBegin; 940f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 941f68a32c8SEmil Constantinescu ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr); 942f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 943f68a32c8SEmil Constantinescu } 944f68a32c8SEmil Constantinescu 945560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype) 946f68a32c8SEmil Constantinescu { 947f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 948f68a32c8SEmil Constantinescu 949f68a32c8SEmil Constantinescu PetscFunctionBegin; 950f68a32c8SEmil Constantinescu *rktype = rk->tableau->name; 951f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 952f68a32c8SEmil Constantinescu } 9532c9cb062Sluzhanghpp 954560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype) 955f68a32c8SEmil Constantinescu { 956f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 957f68a32c8SEmil Constantinescu PetscErrorCode ierr; 958f68a32c8SEmil Constantinescu PetscBool match; 959f68a32c8SEmil Constantinescu RKTableauLink link; 960f68a32c8SEmil Constantinescu 961f68a32c8SEmil Constantinescu PetscFunctionBegin; 962f68a32c8SEmil Constantinescu if (rk->tableau) { 963f68a32c8SEmil Constantinescu ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr); 964f68a32c8SEmil Constantinescu if (match) PetscFunctionReturn(0); 965f68a32c8SEmil Constantinescu } 966f68a32c8SEmil Constantinescu for (link = RKTableauList; link; link=link->next) { 967f68a32c8SEmil Constantinescu ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr); 968f68a32c8SEmil Constantinescu if (match) { 969be5899b3SLisandro Dalcin if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);} 970f68a32c8SEmil Constantinescu rk->tableau = &link->tab; 971be5899b3SLisandro Dalcin if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);} 9722ffb9264SLisandro Dalcin ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 973f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 974f68a32c8SEmil Constantinescu } 975f68a32c8SEmil Constantinescu } 976f68a32c8SEmil Constantinescu SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype); 977e4dd521cSBarry Smith PetscFunctionReturn(0); 978e4dd521cSBarry Smith } 979e4dd521cSBarry Smith 980ff22ae23SHong Zhang static PetscErrorCode TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y) 981ff22ae23SHong Zhang { 982ff22ae23SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 983ff22ae23SHong Zhang 984ff22ae23SHong Zhang PetscFunctionBegin; 9850429704eSStefano Zampini if (ns) *ns = rk->tableau->s; 986d2daff3dSHong Zhang if (Y) *Y = rk->Y; 987ff22ae23SHong Zhang PetscFunctionReturn(0); 988ff22ae23SHong Zhang } 989ff22ae23SHong Zhang 990b3a6b972SJed Brown static PetscErrorCode TSDestroy_RK(TS ts) 991b3a6b972SJed Brown { 992b3a6b972SJed Brown PetscErrorCode ierr; 993b3a6b972SJed Brown 994b3a6b972SJed Brown PetscFunctionBegin; 995b3a6b972SJed Brown ierr = TSReset_RK(ts);CHKERRQ(ierr); 996b3a6b972SJed Brown if (ts->dm) { 997b3a6b972SJed Brown ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 998b3a6b972SJed Brown ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 999b3a6b972SJed Brown } 1000b3a6b972SJed Brown ierr = PetscFree(ts->data);CHKERRQ(ierr); 1001b3a6b972SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr); 1002b3a6b972SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr); 1003*0fe4d17eSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",NULL);CHKERRQ(ierr); 1004*0fe4d17eSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",NULL);CHKERRQ(ierr); 1005b3a6b972SJed Brown PetscFunctionReturn(0); 1006b3a6b972SJed Brown } 1007ff22ae23SHong Zhang 1008ebe3b25bSBarry Smith /*MC 1009f68a32c8SEmil Constantinescu TSRK - ODE and DAE solver using Runge-Kutta schemes 1010ebe3b25bSBarry Smith 1011f68a32c8SEmil Constantinescu The user should provide the right hand side of the equation 1012f68a32c8SEmil Constantinescu using TSSetRHSFunction(). 1013ebe3b25bSBarry Smith 1014f68a32c8SEmil Constantinescu Notes: 101598164e13SLisandro Dalcin The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type 1016ebe3b25bSBarry Smith 1017d41222bbSBarry Smith Level: beginner 1018d41222bbSBarry Smith 1019f68a32c8SEmil Constantinescu .seealso: TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3, 1020*0fe4d17eSHong Zhang TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister(), TSRKSetMultirate(), TSRKGetMultirate() 1021ebe3b25bSBarry Smith 1022ebe3b25bSBarry Smith M*/ 10238cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts) 1024e4dd521cSBarry Smith { 10251566a47fSLisandro Dalcin TS_RK *rk; 1026dfbe8321SBarry Smith PetscErrorCode ierr; 1027e4dd521cSBarry Smith 1028e4dd521cSBarry Smith PetscFunctionBegin; 1029f68a32c8SEmil Constantinescu ierr = TSRKInitializePackage();CHKERRQ(ierr); 1030f68a32c8SEmil Constantinescu 1031f68a32c8SEmil Constantinescu ts->ops->reset = TSReset_RK; 10325f70b478SJed Brown ts->ops->destroy = TSDestroy_RK; 10335f70b478SJed Brown ts->ops->view = TSView_RK; 1034f68a32c8SEmil Constantinescu ts->ops->load = TSLoad_RK; 1035f68a32c8SEmil Constantinescu ts->ops->setup = TSSetUp_RK; 103642f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_RK; 1037f68a32c8SEmil Constantinescu ts->ops->interpolate = TSInterpolate_RK; 10382c9cb062Sluzhanghpp ts->ops->step = TSStep_RK; 1039f68a32c8SEmil Constantinescu ts->ops->evaluatestep = TSEvaluateStep_RK; 1040fecfb714SLisandro Dalcin ts->ops->rollback = TSRollBack_RK; 1041f68a32c8SEmil Constantinescu ts->ops->setfromoptions = TSSetFromOptions_RK; 1042ff22ae23SHong Zhang ts->ops->getstages = TSGetStages_RK; 104342f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_RK; 1044e4dd521cSBarry Smith 10450f7a1166SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_RK; 10460f7a1166SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_RK; 10470f7a1166SHong Zhang 10481566a47fSLisandro Dalcin ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr); 10491566a47fSLisandro Dalcin ts->data = (void*)rk; 1050e4dd521cSBarry Smith 1051f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr); 1052f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr); 1053*0fe4d17eSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",TSRKSetMultirate);CHKERRQ(ierr); 1054*0fe4d17eSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",TSRKGetMultirate);CHKERRQ(ierr); 1055be5899b3SLisandro Dalcin 1056be5899b3SLisandro Dalcin ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr); 105790b67129SHong Zhang rk->dtratio = 1; 1058e4dd521cSBarry Smith PetscFunctionReturn(0); 1059e4dd521cSBarry Smith } 1060