1e4dd521cSBarry Smith /* 22c9cb062Sluzhanghpp Code for time stepping with the Runge-Kutta method(single-rate RK and multi-rate RK) 3f68a32c8SEmil Constantinescu 4f68a32c8SEmil Constantinescu Notes: 52c9cb062Sluzhanghpp 1) The general system is written as 62c9cb062Sluzhanghpp Udot = F(t,U) for single-rate RK; 72c9cb062Sluzhanghpp 2) The general system is written as 855de54fcSHong Zhang Udot = F(t,U) for the nonsplit version of multi-rate RK, 92c9cb062Sluzhanghpp user should give the indexes for both slow and fast components; 102c9cb062Sluzhanghpp 3) The general system is written as 112c9cb062Sluzhanghpp Usdot = Fs(t,Us,Uf) 1255de54fcSHong Zhang Ufdot = Ff(t,Us,Uf) for multi-rate RK with RHS splits, 132c9cb062Sluzhanghpp user should partioned RHS by themselves and also provide the indexes for both slow and fast components. 14e4dd521cSBarry Smith */ 152c9cb062Sluzhanghpp 16af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 17f68a32c8SEmil Constantinescu #include <petscdm.h> 1855de54fcSHong Zhang #include <petscdmshell.h> 19f68a32c8SEmil Constantinescu 20484bcdc7SDebojyoti Ghosh static TSRKType TSRKDefault = TSRK3BS; 212c9cb062Sluzhanghpp static TSRKType TSMRKDefault = TSRK2A; 22f68a32c8SEmil Constantinescu static PetscBool TSRKRegisterAllCalled; 23f68a32c8SEmil Constantinescu static PetscBool TSRKPackageInitialized; 24f68a32c8SEmil Constantinescu 25f68a32c8SEmil Constantinescu typedef struct _RKTableau *RKTableau; 26f68a32c8SEmil Constantinescu struct _RKTableau { 27f68a32c8SEmil Constantinescu char *name; 28d760c35bSDebojyoti Ghosh PetscInt order; /* Classical approximation order of the method i */ 29f68a32c8SEmil Constantinescu PetscInt s; /* Number of stages */ 303a8a9803SLisandro Dalcin PetscInt p; /* Interpolation order */ 31d760c35bSDebojyoti Ghosh PetscBool FSAL; /* flag to indicate if tableau is FSAL */ 32f68a32c8SEmil Constantinescu PetscReal *A,*b,*c; /* Tableau */ 33f68a32c8SEmil Constantinescu PetscReal *bembed; /* Embedded formula of order one less (order-1) */ 34f68a32c8SEmil Constantinescu PetscReal *binterp; /* Dense output formula */ 35f68a32c8SEmil Constantinescu PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */ 36f68a32c8SEmil Constantinescu }; 37f68a32c8SEmil Constantinescu typedef struct _RKTableauLink *RKTableauLink; 38f68a32c8SEmil Constantinescu struct _RKTableauLink { 39f68a32c8SEmil Constantinescu struct _RKTableau tab; 40f68a32c8SEmil Constantinescu RKTableauLink next; 41f68a32c8SEmil Constantinescu }; 42f68a32c8SEmil Constantinescu static RKTableauLink RKTableauList; 43e4dd521cSBarry Smith 44e4dd521cSBarry Smith typedef struct { 45f68a32c8SEmil Constantinescu RKTableau tableau; 4655de54fcSHong Zhang Vec X0; 47f68a32c8SEmil Constantinescu Vec *Y; /* States computed during the step */ 482c9cb062Sluzhanghpp Vec *YdotRHS; /* Function evaluations for the non-stiff part and contains all components */ 492c9cb062Sluzhanghpp Vec *YdotRHS_fast; /* Function evaluations for the non-stiff part and contains fast components */ 502c9cb062Sluzhanghpp Vec *YdotRHS_slow; /* Function evaluations for the non-stiff part and contains slow components */ 51ad8e2604SHong Zhang Vec *VecDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage */ 52ad8e2604SHong Zhang Vec *VecDeltaMu; /* Increment of the adjoint sensitivity w.r.t P at stage */ 530f7a1166SHong Zhang Vec VecCostIntegral0; /* backup for roll-backs due to events */ 54f68a32c8SEmil Constantinescu PetscScalar *work; /* Scalar work */ 552c9cb062Sluzhanghpp PetscInt slow; /* flag indicates call slow components solver (0) or fast components solver (1) */ 56f68a32c8SEmil Constantinescu PetscReal stage_time; 57f68a32c8SEmil Constantinescu TSStepStatus status; 580f7a1166SHong Zhang PetscReal ptime; 590f7a1166SHong Zhang PetscReal time_step; 602c9cb062Sluzhanghpp PetscInt dtratio; /* ratio between slow time step size and fast step size */ 6155de54fcSHong Zhang IS is_fast,is_slow; 6255de54fcSHong Zhang TS subts_fast,subts_slow; 635f70b478SJed Brown } TS_RK; 64e4dd521cSBarry Smith 652c9cb062Sluzhanghpp 66f68a32c8SEmil Constantinescu /*MC 67287c2655SBarry Smith TSRK1FE - First order forward Euler scheme. 68262ff7bbSBarry Smith 69f68a32c8SEmil Constantinescu This method has one stage. 70f68a32c8SEmil Constantinescu 71287c2655SBarry Smith Options database: 729bd3e852SBarry Smith . -ts_rk_type 1fe 73287c2655SBarry Smith 74f68a32c8SEmil Constantinescu Level: advanced 75f68a32c8SEmil Constantinescu 76287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 77f68a32c8SEmil Constantinescu M*/ 78f68a32c8SEmil Constantinescu /*MC 792109b73fSDebojyoti Ghosh TSRK2A - Second order RK scheme. 80f68a32c8SEmil Constantinescu 81f68a32c8SEmil Constantinescu This method has two stages. 82f68a32c8SEmil Constantinescu 83287c2655SBarry Smith Options database: 849bd3e852SBarry Smith . -ts_rk_type 2a 85287c2655SBarry Smith 86f68a32c8SEmil Constantinescu Level: advanced 87f68a32c8SEmil Constantinescu 88287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 89f68a32c8SEmil Constantinescu M*/ 90f68a32c8SEmil Constantinescu /*MC 91f68a32c8SEmil Constantinescu TSRK3 - Third order RK scheme. 92f68a32c8SEmil Constantinescu 93f68a32c8SEmil Constantinescu This method has three stages. 94f68a32c8SEmil Constantinescu 95287c2655SBarry Smith Options database: 969bd3e852SBarry Smith . -ts_rk_type 3 97287c2655SBarry Smith 98f68a32c8SEmil Constantinescu Level: advanced 99f68a32c8SEmil Constantinescu 100287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 101f68a32c8SEmil Constantinescu M*/ 102f68a32c8SEmil Constantinescu /*MC 1032109b73fSDebojyoti Ghosh TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method. 1042109b73fSDebojyoti Ghosh 105d1905564SLisandro Dalcin This method has four stages with the First Same As Last (FSAL) property. 1062109b73fSDebojyoti Ghosh 107287c2655SBarry Smith Options database: 1089bd3e852SBarry Smith . -ts_rk_type 3bs 109287c2655SBarry Smith 1102109b73fSDebojyoti Ghosh Level: advanced 1112109b73fSDebojyoti Ghosh 11298164e13SLisandro Dalcin References: https://doi.org/10.1016/0893-9659(89)90079-7 113d1905564SLisandro Dalcin 114287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 1152109b73fSDebojyoti Ghosh M*/ 1162109b73fSDebojyoti Ghosh /*MC 117f68a32c8SEmil Constantinescu TSRK4 - Fourth order RK scheme. 118f68a32c8SEmil Constantinescu 1192e698f8bSDebojyoti Ghosh This is the classical Runge-Kutta method with four stages. 120f68a32c8SEmil Constantinescu 121287c2655SBarry Smith Options database: 1229bd3e852SBarry Smith . -ts_rk_type 4 123287c2655SBarry Smith 124f68a32c8SEmil Constantinescu Level: advanced 125f68a32c8SEmil Constantinescu 126287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 127f68a32c8SEmil Constantinescu M*/ 128f68a32c8SEmil Constantinescu /*MC 1292e698f8bSDebojyoti Ghosh TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method. 130f68a32c8SEmil Constantinescu 131f68a32c8SEmil Constantinescu This method has six stages. 132f68a32c8SEmil Constantinescu 133287c2655SBarry Smith Options database: 1349bd3e852SBarry Smith . -ts_rk_type 5f 135287c2655SBarry Smith 136f68a32c8SEmil Constantinescu Level: advanced 137f68a32c8SEmil Constantinescu 138287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 139f68a32c8SEmil Constantinescu M*/ 1402109b73fSDebojyoti Ghosh /*MC 1412e698f8bSDebojyoti Ghosh TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method. 1422109b73fSDebojyoti Ghosh 143d1905564SLisandro Dalcin This method has seven stages with the First Same As Last (FSAL) property. 1442109b73fSDebojyoti Ghosh 145287c2655SBarry Smith Options database: 1469bd3e852SBarry Smith . -ts_rk_type 5dp 147287c2655SBarry Smith 1482109b73fSDebojyoti Ghosh Level: advanced 1492109b73fSDebojyoti Ghosh 15098164e13SLisandro Dalcin References: https://doi.org/10.1016/0771-050X(80)90013-3 151d1905564SLisandro Dalcin 152287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 1532109b73fSDebojyoti Ghosh M*/ 15405e23783SLisandro Dalcin /*MC 15505e23783SLisandro Dalcin TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method. 15605e23783SLisandro Dalcin 15705e23783SLisandro Dalcin This method has eight stages with the First Same As Last (FSAL) property. 15805e23783SLisandro Dalcin 159287c2655SBarry Smith Options database: 1609bd3e852SBarry Smith . -ts_rk_type 5bs 161287c2655SBarry Smith 16205e23783SLisandro Dalcin Level: advanced 16305e23783SLisandro Dalcin 16498164e13SLisandro Dalcin References: https://doi.org/10.1016/0898-1221(96)00141-1 16505e23783SLisandro Dalcin 166287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 16705e23783SLisandro Dalcin M*/ 168262ff7bbSBarry Smith 169f68a32c8SEmil Constantinescu /*@C 170f68a32c8SEmil Constantinescu TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK 171262ff7bbSBarry Smith 172f68a32c8SEmil Constantinescu Not Collective, but should be called by all processes which will need the schemes to be registered 173262ff7bbSBarry Smith 174f68a32c8SEmil Constantinescu Level: advanced 175262ff7bbSBarry Smith 176f68a32c8SEmil Constantinescu .keywords: TS, TSRK, register, all 177262ff7bbSBarry Smith 178f68a32c8SEmil Constantinescu .seealso: TSRKRegisterDestroy() 179262ff7bbSBarry Smith @*/ 180f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void) 181262ff7bbSBarry Smith { 1824ac538c5SBarry Smith PetscErrorCode ierr; 183262ff7bbSBarry Smith 184262ff7bbSBarry Smith PetscFunctionBegin; 185f68a32c8SEmil Constantinescu if (TSRKRegisterAllCalled) PetscFunctionReturn(0); 186f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_TRUE; 187e4dd521cSBarry Smith 1884e82626cSLisandro Dalcin #define RC PetscRealConstant 189e4dd521cSBarry Smith { 190f68a32c8SEmil Constantinescu const PetscReal 1914e82626cSLisandro Dalcin A[1][1] = {{0}}, 1924e82626cSLisandro Dalcin b[1] = {RC(1.0)}; 1933a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 194e4dd521cSBarry Smith } 195db046809SBarry Smith { 196f68a32c8SEmil Constantinescu const PetscReal 1974e82626cSLisandro Dalcin A[2][2] = {{0,0}, 1984e82626cSLisandro Dalcin {RC(1.0),0}}, 1994e82626cSLisandro Dalcin b[2] = {RC(0.5),RC(0.5)}, 2004e82626cSLisandro Dalcin bembed[2] = {RC(1.0),0}; 2013a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 202db046809SBarry Smith } 203f68a32c8SEmil Constantinescu { 204f68a32c8SEmil Constantinescu const PetscReal 20517f6384fSBarry Smith A[3][3] = {{0,0,0}, 2064e82626cSLisandro Dalcin {RC(2.0)/RC(3.0),0,0}, 2074e82626cSLisandro Dalcin {RC(-1.0)/RC(3.0),RC(1.0),0}}, 2084e82626cSLisandro Dalcin b[3] = {RC(0.25),RC(0.5),RC(0.25)}; 2093a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 210fdefd5e5SDebojyoti Ghosh } 211fdefd5e5SDebojyoti Ghosh { 212fdefd5e5SDebojyoti Ghosh const PetscReal 21317f6384fSBarry Smith A[4][4] = {{0,0,0,0}, 2144e82626cSLisandro Dalcin {RC(1.0)/RC(2.0),0,0,0}, 2154e82626cSLisandro Dalcin {0,RC(3.0)/RC(4.0),0,0}, 2164e82626cSLisandro Dalcin {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}}, 2174e82626cSLisandro Dalcin b[4] = {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}, 2184e82626cSLisandro 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)}; 2193a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 220db046809SBarry Smith } 221f68a32c8SEmil Constantinescu { 222f68a32c8SEmil Constantinescu const PetscReal 223f68a32c8SEmil Constantinescu A[4][4] = {{0,0,0,0}, 2244e82626cSLisandro Dalcin {RC(0.5),0,0,0}, 2254e82626cSLisandro Dalcin {0,RC(0.5),0,0}, 2264e82626cSLisandro Dalcin {0,0,RC(1.0),0}}, 2274e82626cSLisandro 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)}; 2283a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 229db046809SBarry Smith } 230f68a32c8SEmil Constantinescu { 231f68a32c8SEmil Constantinescu const PetscReal 23217f6384fSBarry Smith A[6][6] = {{0,0,0,0,0,0}, 2334e82626cSLisandro Dalcin {RC(0.25),0,0,0,0,0}, 2344e82626cSLisandro Dalcin {RC(3.0)/RC(32.0),RC(9.0)/RC(32.0),0,0,0,0}, 2354e82626cSLisandro Dalcin {RC(1932.0)/RC(2197.0),RC(-7200.0)/RC(2197.0),RC(7296.0)/RC(2197.0),0,0,0}, 2364e82626cSLisandro Dalcin {RC(439.0)/RC(216.0),RC(-8.0),RC(3680.0)/RC(513.0),RC(-845.0)/RC(4104.0),0,0}, 2374e82626cSLisandro 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}}, 2384e82626cSLisandro 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)}, 2394e82626cSLisandro 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}; 2403a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 241fdefd5e5SDebojyoti Ghosh } 242fdefd5e5SDebojyoti Ghosh { 243fdefd5e5SDebojyoti Ghosh const PetscReal 24417f6384fSBarry Smith A[7][7] = {{0,0,0,0,0,0,0}, 2454e82626cSLisandro Dalcin {RC(1.0)/RC(5.0),0,0,0,0,0,0}, 2464e82626cSLisandro Dalcin {RC(3.0)/RC(40.0),RC(9.0)/RC(40.0),0,0,0,0,0}, 2474e82626cSLisandro Dalcin {RC(44.0)/RC(45.0),RC(-56.0)/RC(15.0),RC(32.0)/RC(9.0),0,0,0,0}, 2484e82626cSLisandro 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}, 2494e82626cSLisandro 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}, 2504e82626cSLisandro 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}}, 2514e82626cSLisandro 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}, 252a685a763Sluzhanghpp 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)}, 253a685a763Sluzhanghpp 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)}, 254a685a763Sluzhanghpp {0,0,0,0,0}, 255a685a763Sluzhanghpp {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)}, 256a685a763Sluzhanghpp {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)}, 257a685a763Sluzhanghpp {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)}, 258a685a763Sluzhanghpp {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)}, 259a685a763Sluzhanghpp {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)}}; 260a685a763Sluzhanghpp 261a685a763Sluzhanghpp ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,5,binterp[0]);CHKERRQ(ierr); 262f68a32c8SEmil Constantinescu } 26305e23783SLisandro Dalcin { 26405e23783SLisandro Dalcin const PetscReal 26517f6384fSBarry Smith A[8][8] = {{0,0,0,0,0,0,0,0}, 2664e82626cSLisandro Dalcin {RC(1.0)/RC(6.0),0,0,0,0,0,0,0}, 2674e82626cSLisandro Dalcin {RC(2.0)/RC(27.0),RC(4.0)/RC(27.0),0,0,0,0,0,0}, 2684e82626cSLisandro Dalcin {RC(183.0)/RC(1372.0),RC(-162.0)/RC(343.0),RC(1053.0)/RC(1372.0),0,0,0,0,0}, 2694e82626cSLisandro 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}, 2704e82626cSLisandro 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}, 2714e82626cSLisandro 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}, 2724e82626cSLisandro 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}}, 2734e82626cSLisandro 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}, 2744e82626cSLisandro 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)}; 27505e23783SLisandro Dalcin ierr = TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 27605e23783SLisandro Dalcin } 2774e82626cSLisandro Dalcin #undef RC 278db046809SBarry Smith PetscFunctionReturn(0); 279db046809SBarry Smith } 280db046809SBarry Smith 281f68a32c8SEmil Constantinescu /*@C 282f68a32c8SEmil Constantinescu TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister(). 283f68a32c8SEmil Constantinescu 284f68a32c8SEmil Constantinescu Not Collective 285f68a32c8SEmil Constantinescu 286f68a32c8SEmil Constantinescu Level: advanced 287f68a32c8SEmil Constantinescu 288f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy 289f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll() 290f68a32c8SEmil Constantinescu @*/ 291f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void) 292e4dd521cSBarry Smith { 293f68a32c8SEmil Constantinescu PetscErrorCode ierr; 294f68a32c8SEmil Constantinescu RKTableauLink link; 295f68a32c8SEmil Constantinescu 296f68a32c8SEmil Constantinescu PetscFunctionBegin; 297f68a32c8SEmil Constantinescu while ((link = RKTableauList)) { 298f68a32c8SEmil Constantinescu RKTableau t = &link->tab; 299f68a32c8SEmil Constantinescu RKTableauList = link->next; 300f68a32c8SEmil Constantinescu ierr = PetscFree3(t->A,t->b,t->c); CHKERRQ(ierr); 301f68a32c8SEmil Constantinescu ierr = PetscFree (t->bembed); CHKERRQ(ierr); 302f68a32c8SEmil Constantinescu ierr = PetscFree (t->binterp); CHKERRQ(ierr); 303f68a32c8SEmil Constantinescu ierr = PetscFree (t->name); CHKERRQ(ierr); 304f68a32c8SEmil Constantinescu ierr = PetscFree (link); CHKERRQ(ierr); 305f68a32c8SEmil Constantinescu } 306f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_FALSE; 307f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 308f68a32c8SEmil Constantinescu } 309f68a32c8SEmil Constantinescu 310f68a32c8SEmil Constantinescu /*@C 311f68a32c8SEmil Constantinescu TSRKInitializePackage - This function initializes everything in the TSRK package. It is called 3128a690491SBarry Smith from TSInitializePackage(). 313f68a32c8SEmil Constantinescu 314f68a32c8SEmil Constantinescu Level: developer 315f68a32c8SEmil Constantinescu 316f68a32c8SEmil Constantinescu .keywords: TS, TSRK, initialize, package 317f68a32c8SEmil Constantinescu .seealso: PetscInitialize() 318f68a32c8SEmil Constantinescu @*/ 319f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void) 320f68a32c8SEmil Constantinescu { 321186e87acSLisandro Dalcin PetscErrorCode ierr; 322e4dd521cSBarry Smith 323e4dd521cSBarry Smith PetscFunctionBegin; 324f68a32c8SEmil Constantinescu if (TSRKPackageInitialized) PetscFunctionReturn(0); 325f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_TRUE; 326f68a32c8SEmil Constantinescu ierr = TSRKRegisterAll();CHKERRQ(ierr); 327f68a32c8SEmil Constantinescu ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr); 328f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 329f68a32c8SEmil Constantinescu } 330186e87acSLisandro Dalcin 331f68a32c8SEmil Constantinescu /*@C 332f68a32c8SEmil Constantinescu TSRKFinalizePackage - This function destroys everything in the TSRK package. It is 333f68a32c8SEmil Constantinescu called from PetscFinalize(). 334186e87acSLisandro Dalcin 335f68a32c8SEmil Constantinescu Level: developer 336f68a32c8SEmil Constantinescu 337f68a32c8SEmil Constantinescu .keywords: Petsc, destroy, package 338f68a32c8SEmil Constantinescu .seealso: PetscFinalize() 339f68a32c8SEmil Constantinescu @*/ 340f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void) 341f68a32c8SEmil Constantinescu { 342f68a32c8SEmil Constantinescu PetscErrorCode ierr; 343f68a32c8SEmil Constantinescu 344f68a32c8SEmil Constantinescu PetscFunctionBegin; 345f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_FALSE; 346f68a32c8SEmil Constantinescu ierr = TSRKRegisterDestroy();CHKERRQ(ierr); 347f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 348f68a32c8SEmil Constantinescu } 349f68a32c8SEmil Constantinescu 350f68a32c8SEmil Constantinescu /*@C 351f68a32c8SEmil Constantinescu TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 352f68a32c8SEmil Constantinescu 353f68a32c8SEmil Constantinescu Not Collective, but the same schemes should be registered on all processes on which they will be used 354f68a32c8SEmil Constantinescu 355f68a32c8SEmil Constantinescu Input Parameters: 356f68a32c8SEmil Constantinescu + name - identifier for method 357f68a32c8SEmil Constantinescu . order - approximation order of method 358f68a32c8SEmil Constantinescu . s - number of stages, this is the dimension of the matrices below 359f68a32c8SEmil Constantinescu . A - stage coefficients (dimension s*s, row-major) 360f68a32c8SEmil Constantinescu . b - step completion table (dimension s; NULL to use last row of A) 361f68a32c8SEmil Constantinescu . c - abscissa (dimension s; NULL to use row sums of A) 362f68a32c8SEmil Constantinescu . bembed - completion table for embedded method (dimension s; NULL if not available) 3633a8a9803SLisandro Dalcin . p - Order of the interpolation scheme, equal to the number of columns of binterp 3643a8a9803SLisandro Dalcin - binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1) 365f68a32c8SEmil Constantinescu 366f68a32c8SEmil Constantinescu Notes: 367f68a32c8SEmil Constantinescu Several RK methods are provided, this function is only needed to create new methods. 368f68a32c8SEmil Constantinescu 369f68a32c8SEmil Constantinescu Level: advanced 370f68a32c8SEmil Constantinescu 371f68a32c8SEmil Constantinescu .keywords: TS, register 372f68a32c8SEmil Constantinescu 373f68a32c8SEmil Constantinescu .seealso: TSRK 374f68a32c8SEmil Constantinescu @*/ 375f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s, 376f68a32c8SEmil Constantinescu const PetscReal A[],const PetscReal b[],const PetscReal c[], 3773a8a9803SLisandro Dalcin const PetscReal bembed[],PetscInt p,const PetscReal binterp[]) 378f68a32c8SEmil Constantinescu { 379f68a32c8SEmil Constantinescu PetscErrorCode ierr; 380f68a32c8SEmil Constantinescu RKTableauLink link; 381f68a32c8SEmil Constantinescu RKTableau t; 382f68a32c8SEmil Constantinescu PetscInt i,j; 383f68a32c8SEmil Constantinescu 384f68a32c8SEmil Constantinescu PetscFunctionBegin; 3853a8a9803SLisandro Dalcin PetscValidCharPointer(name,1); 3863a8a9803SLisandro Dalcin PetscValidRealPointer(A,4); 3873a8a9803SLisandro Dalcin if (b) PetscValidRealPointer(b,5); 3883a8a9803SLisandro Dalcin if (c) PetscValidRealPointer(c,6); 3893a8a9803SLisandro Dalcin if (bembed) PetscValidRealPointer(bembed,7); 3903a8a9803SLisandro Dalcin if (binterp || p > 1) PetscValidRealPointer(binterp,9); 3913a8a9803SLisandro Dalcin 3921d36bdfdSBarry Smith ierr = TSRKInitializePackage();CHKERRQ(ierr); 39395dccacaSBarry Smith ierr = PetscNew(&link);CHKERRQ(ierr); 394f68a32c8SEmil Constantinescu t = &link->tab; 3953a8a9803SLisandro Dalcin 396f68a32c8SEmil Constantinescu ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 397f68a32c8SEmil Constantinescu t->order = order; 398f68a32c8SEmil Constantinescu t->s = s; 399dcca6d9dSJed Brown ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr); 400f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 401f68a32c8SEmil Constantinescu if (b) { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); } 402f68a32c8SEmil Constantinescu else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i]; 403f68a32c8SEmil Constantinescu if (c) { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); } 404f68a32c8SEmil 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]; 405d760c35bSDebojyoti Ghosh t->FSAL = PETSC_TRUE; 406d760c35bSDebojyoti Ghosh for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE; 4073a8a9803SLisandro Dalcin 408f68a32c8SEmil Constantinescu if (bembed) { 409785e854fSJed Brown ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr); 410f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr); 411f68a32c8SEmil Constantinescu } 412f68a32c8SEmil Constantinescu 4133a8a9803SLisandro Dalcin if (!binterp) { p = 1; binterp = t->b; } 4143a8a9803SLisandro Dalcin t->p = p; 4153a8a9803SLisandro Dalcin ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr); 4163a8a9803SLisandro Dalcin ierr = PetscMemcpy(t->binterp,binterp,s*p*sizeof(binterp[0]));CHKERRQ(ierr); 4173a8a9803SLisandro Dalcin 418f68a32c8SEmil Constantinescu link->next = RKTableauList; 419f68a32c8SEmil Constantinescu RKTableauList = link; 420f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 421f68a32c8SEmil Constantinescu } 422f68a32c8SEmil Constantinescu 423e4dd521cSBarry Smith /* 4242c9cb062Sluzhanghpp This is for single-step RK method 425f68a32c8SEmil Constantinescu The step completion formula is 426e4dd521cSBarry Smith 427f68a32c8SEmil Constantinescu x1 = x0 + h b^T YdotRHS 428f68a32c8SEmil Constantinescu 429f68a32c8SEmil Constantinescu This function can be called before or after ts->vec_sol has been updated. 430f68a32c8SEmil Constantinescu Suppose we have a completion formula (b) and an embedded formula (be) of different order. 431f68a32c8SEmil Constantinescu We can write 432f68a32c8SEmil Constantinescu 433f68a32c8SEmil Constantinescu x1e = x0 + h be^T YdotRHS 434f68a32c8SEmil Constantinescu = x1 - h b^T YdotRHS + h be^T YdotRHS 435f68a32c8SEmil Constantinescu = x1 + h (be - b)^T YdotRHS 436f68a32c8SEmil Constantinescu 437f68a32c8SEmil Constantinescu so we can evaluate the method with different order even after the step has been optimistically completed. 438e4dd521cSBarry Smith */ 439f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done) 440f68a32c8SEmil Constantinescu { 441f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 442f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 443f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 444f68a32c8SEmil Constantinescu PetscReal h; 445f68a32c8SEmil Constantinescu PetscInt s = tab->s,j; 446f68a32c8SEmil Constantinescu PetscErrorCode ierr; 447f68a32c8SEmil Constantinescu 448f68a32c8SEmil Constantinescu PetscFunctionBegin; 449f68a32c8SEmil Constantinescu switch (rk->status) { 450f68a32c8SEmil Constantinescu case TS_STEP_INCOMPLETE: 451f68a32c8SEmil Constantinescu case TS_STEP_PENDING: 452f68a32c8SEmil Constantinescu h = ts->time_step; break; 453f68a32c8SEmil Constantinescu case TS_STEP_COMPLETE: 454be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 455f68a32c8SEmil Constantinescu default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 456f68a32c8SEmil Constantinescu } 457f68a32c8SEmil Constantinescu if (order == tab->order) { 458f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { 459f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 460*90b67129SHong Zhang for (j=0; j<s; j++) w[j] = h*tab->b[j]/rk->dtratio; 461f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 462f68a32c8SEmil Constantinescu } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 463f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 464f68a32c8SEmil Constantinescu } else if (order == tab->order-1) { 465f68a32c8SEmil Constantinescu if (!tab->bembed) goto unavailable; 466f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { /*Complete with the embedded method (be)*/ 467f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 468f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 469f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 470f68a32c8SEmil Constantinescu } else { /*Rollback and re-complete using (be-b) */ 471f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 472f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 473f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 474f68a32c8SEmil Constantinescu } 475f68a32c8SEmil Constantinescu if (done) *done = PETSC_TRUE; 476f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 477f68a32c8SEmil Constantinescu } 478f68a32c8SEmil Constantinescu unavailable: 479f68a32c8SEmil Constantinescu if (done) *done = PETSC_FALSE; 480a7fac7c2SEmil 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); 481f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 482f68a32c8SEmil Constantinescu } 483f68a32c8SEmil Constantinescu 4840f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts) 4850f7a1166SHong Zhang { 4860f7a1166SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 4870f7a1166SHong Zhang RKTableau tab = rk->tableau; 4880f7a1166SHong Zhang const PetscInt s = tab->s; 4890f7a1166SHong Zhang const PetscReal *b = tab->b,*c = tab->c; 4900f7a1166SHong Zhang Vec *Y = rk->Y; 4910f7a1166SHong Zhang PetscInt i; 4920f7a1166SHong Zhang PetscErrorCode ierr; 4930f7a1166SHong Zhang 4940f7a1166SHong Zhang PetscFunctionBegin; 4950f7a1166SHong Zhang /* backup cost integral */ 4960f7a1166SHong Zhang ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr); 4970f7a1166SHong Zhang for (i=s-1; i>=0; i--) { 4980f7a1166SHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 49951a7f651SHong Zhang ierr = TSComputeCostIntegrand(ts,rk->ptime+rk->time_step*c[i],Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 5000f7a1166SHong Zhang ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 5010f7a1166SHong Zhang } 5020f7a1166SHong Zhang PetscFunctionReturn(0); 5030f7a1166SHong Zhang } 5040f7a1166SHong Zhang 5050f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts) 5060f7a1166SHong Zhang { 5070f7a1166SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 5080f7a1166SHong Zhang RKTableau tab = rk->tableau; 5090f7a1166SHong Zhang const PetscInt s = tab->s; 5100f7a1166SHong Zhang const PetscReal *b = tab->b,*c = tab->c; 5110f7a1166SHong Zhang Vec *Y = rk->Y; 5120f7a1166SHong Zhang PetscInt i; 5130f7a1166SHong Zhang PetscErrorCode ierr; 5140f7a1166SHong Zhang 5150f7a1166SHong Zhang PetscFunctionBegin; 5160f7a1166SHong Zhang for (i=s-1; i>=0; i--) { 5170f7a1166SHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 51851a7f651SHong Zhang ierr = TSComputeCostIntegrand(ts,ts->ptime+ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 5190f7a1166SHong Zhang ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 5200f7a1166SHong Zhang } 5210f7a1166SHong Zhang PetscFunctionReturn(0); 5220f7a1166SHong Zhang } 5230f7a1166SHong Zhang 524fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts) 525fecfb714SLisandro Dalcin { 526fecfb714SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 527fecfb714SLisandro Dalcin RKTableau tab = rk->tableau; 528fecfb714SLisandro Dalcin const PetscInt s = tab->s; 529fecfb714SLisandro Dalcin const PetscReal *b = tab->b; 530fecfb714SLisandro Dalcin PetscScalar *w = rk->work; 531fecfb714SLisandro Dalcin Vec *YdotRHS = rk->YdotRHS; 532fecfb714SLisandro Dalcin PetscInt j; 533fecfb714SLisandro Dalcin PetscReal h; 534fecfb714SLisandro Dalcin PetscErrorCode ierr; 535fecfb714SLisandro Dalcin 536fecfb714SLisandro Dalcin PetscFunctionBegin; 537fecfb714SLisandro Dalcin switch (rk->status) { 538fecfb714SLisandro Dalcin case TS_STEP_INCOMPLETE: 539fecfb714SLisandro Dalcin case TS_STEP_PENDING: 540fecfb714SLisandro Dalcin h = ts->time_step; break; 541fecfb714SLisandro Dalcin case TS_STEP_COMPLETE: 542fecfb714SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 543fecfb714SLisandro Dalcin default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 544fecfb714SLisandro Dalcin } 545fecfb714SLisandro Dalcin for (j=0; j<s; j++) w[j] = -h*b[j]; 546fecfb714SLisandro Dalcin ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr); 547fecfb714SLisandro Dalcin PetscFunctionReturn(0); 548fecfb714SLisandro Dalcin } 549fecfb714SLisandro Dalcin 550f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts) 551f68a32c8SEmil Constantinescu { 552f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 553f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 554f68a32c8SEmil Constantinescu const PetscInt s = tab->s; 555fecfb714SLisandro Dalcin const PetscReal *A = tab->A,*c = tab->c; 556f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 557f68a32c8SEmil Constantinescu Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 558d1905564SLisandro Dalcin PetscBool FSAL = tab->FSAL; 559f68a32c8SEmil Constantinescu TSAdapt adapt; 560fecfb714SLisandro Dalcin PetscInt i,j; 561be5899b3SLisandro Dalcin PetscInt rejections = 0; 562be5899b3SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 563be5899b3SLisandro Dalcin PetscReal next_time_step = ts->time_step; 564f68a32c8SEmil Constantinescu PetscErrorCode ierr; 565f68a32c8SEmil Constantinescu 566f68a32c8SEmil Constantinescu PetscFunctionBegin; 567d1905564SLisandro Dalcin if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE; 568d1905564SLisandro Dalcin if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); } 569d1905564SLisandro Dalcin 570f68a32c8SEmil Constantinescu rk->status = TS_STEP_INCOMPLETE; 571be5899b3SLisandro Dalcin while (!ts->reason && rk->status != TS_STEP_COMPLETE) { 572be5899b3SLisandro Dalcin PetscReal t = ts->ptime; 573f68a32c8SEmil Constantinescu PetscReal h = ts->time_step; 574f68a32c8SEmil Constantinescu for (i=0; i<s; i++) { 5759be3e283SDebojyoti Ghosh rk->stage_time = t + h*c[i]; 5769be3e283SDebojyoti Ghosh ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr); 577f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 578f68a32c8SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 579f68a32c8SEmil Constantinescu ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 5809be3e283SDebojyoti Ghosh ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 581f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 582be5899b3SLisandro Dalcin ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr); 583fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 5848f738a7cSLisandro Dalcin if (FSAL && !i) continue; 585f68a32c8SEmil Constantinescu ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 586f68a32c8SEmil Constantinescu } 587be5899b3SLisandro Dalcin 588be5899b3SLisandro Dalcin rk->status = TS_STEP_INCOMPLETE; 589f68a32c8SEmil Constantinescu ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 590f68a32c8SEmil Constantinescu rk->status = TS_STEP_PENDING; 591f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 592f68a32c8SEmil Constantinescu ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 5931917a363SLisandro Dalcin ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr); 594fecfb714SLisandro Dalcin ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 595be5899b3SLisandro Dalcin rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 596be5899b3SLisandro Dalcin if (!accept) { /*Roll back the current step */ 597fecfb714SLisandro Dalcin ierr = TSRollBack_RK(ts);CHKERRQ(ierr); 598be5899b3SLisandro Dalcin ts->time_step = next_time_step; 599be5899b3SLisandro Dalcin goto reject_step; 600be5899b3SLisandro Dalcin } 601be5899b3SLisandro Dalcin 6020f7a1166SHong Zhang if (ts->costintegralfwd) { /*Save the info for the later use in cost integral evaluation*/ 6030f7a1166SHong Zhang rk->ptime = ts->ptime; 6040f7a1166SHong Zhang rk->time_step = ts->time_step; 605493ed95fSHong Zhang } 606be5899b3SLisandro Dalcin 607f68a32c8SEmil Constantinescu ts->ptime += ts->time_step; 608f68a32c8SEmil Constantinescu ts->time_step = next_time_step; 609f68a32c8SEmil Constantinescu break; 610be5899b3SLisandro Dalcin 611be5899b3SLisandro Dalcin reject_step: 612fecfb714SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 613be5899b3SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 614be5899b3SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 615be5899b3SLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 616f68a32c8SEmil Constantinescu } 617f68a32c8SEmil Constantinescu } 618f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 619e4dd521cSBarry Smith } 620e4dd521cSBarry Smith 621f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts) 622f6a906c0SBarry Smith { 623f6a906c0SBarry Smith TS_RK *rk = (TS_RK*)ts->data; 624be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 625be5899b3SLisandro Dalcin PetscInt s = tab->s; 626f6a906c0SBarry Smith PetscErrorCode ierr; 627f6a906c0SBarry Smith 628f6a906c0SBarry Smith PetscFunctionBegin; 629f6a906c0SBarry Smith if (ts->adjointsetupcalled++) PetscFunctionReturn(0); 630abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr); 631f6a906c0SBarry Smith if(ts->vecs_sensip) { 632abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr); 633f6a906c0SBarry Smith } 634f6a906c0SBarry Smith PetscFunctionReturn(0); 635f6a906c0SBarry Smith } 636f6a906c0SBarry Smith 63742f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts) 638d2daff3dSHong Zhang { 639c235aa8dSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 640c235aa8dSHong Zhang RKTableau tab = rk->tableau; 641c235aa8dSHong Zhang const PetscInt s = tab->s; 642c235aa8dSHong Zhang const PetscReal *A = tab->A,*b = tab->b,*c = tab->c; 643c235aa8dSHong Zhang PetscScalar *w = rk->work; 6443389c7d9SStefano Zampini Vec *Y = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu; 645b18ea86cSHong Zhang PetscInt i,j,nadj; 6463d3eaba7SBarry Smith PetscReal t = ts->ptime; 647d2daff3dSHong Zhang PetscErrorCode ierr; 648cef76868SBarry Smith PetscReal h = ts->time_step; 649c235aa8dSHong Zhang 650d2daff3dSHong Zhang PetscFunctionBegin; 651c235aa8dSHong Zhang rk->status = TS_STEP_INCOMPLETE; 652c235aa8dSHong Zhang for (i=s-1; i>=0; i--) { 6533389c7d9SStefano Zampini Mat J; 6543389c7d9SStefano Zampini PetscReal stage_time = t + h*(1.0-c[i]); 6553389c7d9SStefano Zampini PetscBool zero = PETSC_FALSE; 6563389c7d9SStefano Zampini 6573389c7d9SStefano Zampini ierr = TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr); 6583389c7d9SStefano Zampini ierr = TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);CHKERRQ(ierr); 659493ed95fSHong Zhang if (ts->vec_costintegral) { 6603389c7d9SStefano Zampini ierr = TSAdjointComputeDRDYFunction(ts,stage_time,Y[i],ts->vecs_drdy);CHKERRQ(ierr); 661493ed95fSHong Zhang } 6623389c7d9SStefano Zampini /* Stage values of mu */ 6633389c7d9SStefano Zampini if (ts->vecs_sensip) { 6643389c7d9SStefano Zampini ierr = TSAdjointComputeRHSJacobian(ts,stage_time,Y[i],ts->Jacp);CHKERRQ(ierr); 6653389c7d9SStefano Zampini if (ts->vec_costintegral) { 6663389c7d9SStefano Zampini ierr = TSAdjointComputeDRDPFunction(ts,stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr); 6673389c7d9SStefano Zampini } 6683389c7d9SStefano Zampini } 6693389c7d9SStefano Zampini 6703389c7d9SStefano Zampini if (b[i] == 0 && i == s-1) zero = PETSC_TRUE; 671abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 6723389c7d9SStefano Zampini DM dm; 6733389c7d9SStefano Zampini Vec VecSensiTemp; 6743389c7d9SStefano Zampini 6753389c7d9SStefano Zampini ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 6763389c7d9SStefano Zampini ierr = DMGetGlobalVector(dm,&VecSensiTemp);CHKERRQ(ierr); 6773389c7d9SStefano Zampini /* Stage values of lambda */ 6783389c7d9SStefano Zampini if (!zero) { 6793389c7d9SStefano Zampini ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp);CHKERRQ(ierr); 6803389c7d9SStefano Zampini ierr = VecScale(VecSensiTemp,-h*b[i]);CHKERRQ(ierr); 6813389c7d9SStefano Zampini for (j=i+1; j<s; j++) w[j-i-1] = -h*A[j*s+i]; 6823389c7d9SStefano Zampini ierr = VecMAXPY(VecSensiTemp,s-i-1,w,&VecDeltaLam[nadj*s+i+1]);CHKERRQ(ierr); 6833389c7d9SStefano Zampini ierr = MatMultTranspose(J,VecSensiTemp,VecDeltaLam[nadj*s+i]);CHKERRQ(ierr); 6843389c7d9SStefano Zampini } else { 6853389c7d9SStefano Zampini ierr = VecSet(VecDeltaLam[nadj*s+i],0);CHKERRQ(ierr); 6863389c7d9SStefano Zampini } 687493ed95fSHong Zhang if (ts->vec_costintegral) { 688493ed95fSHong Zhang ierr = VecAXPY(VecDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdy[nadj]);CHKERRQ(ierr); 689493ed95fSHong Zhang } 6906497ce14SHong Zhang 691ad8e2604SHong Zhang /* Stage values of mu */ 6926497ce14SHong Zhang if (ts->vecs_sensip) { 6933389c7d9SStefano Zampini if (!zero) { 6943389c7d9SStefano Zampini ierr = MatMultTranspose(ts->Jacp,VecSensiTemp,VecDeltaMu[nadj*s+i]);CHKERRQ(ierr); 6953389c7d9SStefano Zampini } else { 6963389c7d9SStefano Zampini ierr = VecSet(VecDeltaMu[nadj*s+i],0);CHKERRQ(ierr); 697493ed95fSHong Zhang } 698493ed95fSHong Zhang if (ts->vec_costintegral) { 699493ed95fSHong Zhang ierr = VecAXPY(VecDeltaMu[nadj*s+i],-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr); 700493ed95fSHong Zhang } 701ad8e2604SHong Zhang } 7023389c7d9SStefano Zampini ierr = DMRestoreGlobalVector(dm,&VecSensiTemp);CHKERRQ(ierr); 703c235aa8dSHong Zhang } 7046497ce14SHong Zhang } 705c235aa8dSHong Zhang 706c235aa8dSHong Zhang for (j=0; j<s; j++) w[j] = 1.0; 707abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 708b18ea86cSHong Zhang ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr); 7096497ce14SHong Zhang if (ts->vecs_sensip) { 710ad8e2604SHong Zhang ierr = VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);CHKERRQ(ierr); 711b18ea86cSHong Zhang } 7126497ce14SHong Zhang } 713c235aa8dSHong Zhang rk->status = TS_STEP_COMPLETE; 714d2daff3dSHong Zhang PetscFunctionReturn(0); 715d2daff3dSHong Zhang } 716d2daff3dSHong Zhang 71755de54fcSHong Zhang static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X) 71855de54fcSHong Zhang { 71955de54fcSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 72055de54fcSHong Zhang PetscInt s = rk->tableau->s,p = rk->tableau->p,i,j; 72155de54fcSHong Zhang PetscReal h; 72255de54fcSHong Zhang PetscReal tt,t; 72355de54fcSHong Zhang PetscScalar *b; 72455de54fcSHong Zhang const PetscReal *B = rk->tableau->binterp; 72555de54fcSHong Zhang PetscErrorCode ierr; 72655de54fcSHong Zhang 72755de54fcSHong Zhang PetscFunctionBegin; 72855de54fcSHong Zhang if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 72955de54fcSHong Zhang 73055de54fcSHong Zhang switch (rk->status) { 73155de54fcSHong Zhang case TS_STEP_INCOMPLETE: 73255de54fcSHong Zhang case TS_STEP_PENDING: 73355de54fcSHong Zhang h = ts->time_step; 73455de54fcSHong Zhang t = (itime - ts->ptime)/h; 73555de54fcSHong Zhang break; 73655de54fcSHong Zhang case TS_STEP_COMPLETE: 73755de54fcSHong Zhang h = ts->ptime - ts->ptime_prev; 73855de54fcSHong Zhang t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 73955de54fcSHong Zhang break; 74055de54fcSHong Zhang default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 74155de54fcSHong Zhang } 74255de54fcSHong Zhang ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 74355de54fcSHong Zhang for (i=0; i<s; i++) b[i] = 0; 74455de54fcSHong Zhang for (j=0,tt=t; j<p; j++,tt*=t) { 74555de54fcSHong Zhang for (i=0; i<s; i++) { 74655de54fcSHong Zhang b[i] += h * B[i*p+j] * tt; 74755de54fcSHong Zhang } 74855de54fcSHong Zhang } 74955de54fcSHong Zhang ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr); 75055de54fcSHong Zhang ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr); 75155de54fcSHong Zhang ierr = PetscFree(b);CHKERRQ(ierr); 75255de54fcSHong Zhang PetscFunctionReturn(0); 75355de54fcSHong Zhang } 75455de54fcSHong Zhang 75555de54fcSHong Zhang /*------------------------------------------------------------*/ 75655de54fcSHong Zhang 757be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts) 758be5899b3SLisandro Dalcin { 759be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 760be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 761be5899b3SLisandro Dalcin PetscErrorCode ierr; 762be5899b3SLisandro Dalcin 763be5899b3SLisandro Dalcin PetscFunctionBegin; 764be5899b3SLisandro Dalcin if (!tab) PetscFunctionReturn(0); 765be5899b3SLisandro Dalcin ierr = PetscFree(rk->work);CHKERRQ(ierr); 766be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr); 767be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr); 768be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr); 769be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr); 770be5899b3SLisandro Dalcin PetscFunctionReturn(0); 771be5899b3SLisandro Dalcin } 772be5899b3SLisandro Dalcin 773277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts) 774e4dd521cSBarry Smith { 7755f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 7766849ba73SBarry Smith PetscErrorCode ierr; 777e4dd521cSBarry Smith 778e4dd521cSBarry Smith PetscFunctionBegin; 779be5899b3SLisandro Dalcin ierr = TSRKTableauReset(ts);CHKERRQ(ierr); 7800f7a1166SHong Zhang ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr); 78155de54fcSHong Zhang ierr = PetscTryMethod(ts,"TSReset_MRKSPLIT_C",(TS),(ts));CHKERRQ(ierr); 78255de54fcSHong Zhang ierr = PetscTryMethod(ts,"TSReset_MRKNONSPLIT_C",(TS),(ts));CHKERRQ(ierr); 783277b19d0SLisandro Dalcin PetscFunctionReturn(0); 784e4dd521cSBarry Smith } 785277b19d0SLisandro Dalcin 786f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx) 787f68a32c8SEmil Constantinescu { 788f68a32c8SEmil Constantinescu PetscFunctionBegin; 789f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 790f68a32c8SEmil Constantinescu } 791f68a32c8SEmil Constantinescu 792f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 793f68a32c8SEmil Constantinescu { 794f68a32c8SEmil Constantinescu PetscFunctionBegin; 795f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 796f68a32c8SEmil Constantinescu } 797f68a32c8SEmil Constantinescu 798f68a32c8SEmil Constantinescu 799f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx) 800f68a32c8SEmil Constantinescu { 801f68a32c8SEmil Constantinescu PetscFunctionBegin; 802f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 803f68a32c8SEmil Constantinescu } 804f68a32c8SEmil Constantinescu 805f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 806f68a32c8SEmil Constantinescu { 807f68a32c8SEmil Constantinescu 808f68a32c8SEmil Constantinescu PetscFunctionBegin; 809f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 810f68a32c8SEmil Constantinescu } 811c235aa8dSHong Zhang /* 812d2daff3dSHong Zhang static PetscErrorCode RKSetAdjCoe(RKTableau tab) 813d2daff3dSHong Zhang { 814d2daff3dSHong Zhang PetscReal *A,*b; 815d2daff3dSHong Zhang PetscInt s,i,j; 816d2daff3dSHong Zhang PetscErrorCode ierr; 817d2daff3dSHong Zhang 818d2daff3dSHong Zhang PetscFunctionBegin; 819d2daff3dSHong Zhang s = tab->s; 820d2daff3dSHong Zhang ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr); 821d2daff3dSHong Zhang 822d2daff3dSHong Zhang for (i=0; i<s; i++) 823d2daff3dSHong Zhang for (j=0; j<s; j++) { 824d2daff3dSHong 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]; 825d2daff3dSHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr); 826d2daff3dSHong Zhang } 827d2daff3dSHong Zhang 828d2daff3dSHong Zhang for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i]; 829d2daff3dSHong Zhang 830d2daff3dSHong Zhang ierr = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 831d2daff3dSHong Zhang ierr = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr); 832d2daff3dSHong Zhang ierr = PetscFree2(A,b);CHKERRQ(ierr); 833d2daff3dSHong Zhang PetscFunctionReturn(0); 834d2daff3dSHong Zhang } 835c235aa8dSHong Zhang */ 836be5899b3SLisandro Dalcin 837be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts) 838be5899b3SLisandro Dalcin { 839be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 840be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 841be5899b3SLisandro Dalcin PetscErrorCode ierr; 842be5899b3SLisandro Dalcin 843be5899b3SLisandro Dalcin PetscFunctionBegin; 844be5899b3SLisandro Dalcin ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr); 845be5899b3SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr); 846be5899b3SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr); 847be5899b3SLisandro Dalcin PetscFunctionReturn(0); 848be5899b3SLisandro Dalcin } 849be5899b3SLisandro Dalcin 850f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts) 851f68a32c8SEmil Constantinescu { 852f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 853f68a32c8SEmil Constantinescu PetscErrorCode ierr; 854f68a32c8SEmil Constantinescu DM dm; 855f68a32c8SEmil Constantinescu 856f68a32c8SEmil Constantinescu PetscFunctionBegin; 8573dd83b38SBarry Smith ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr); 858be5899b3SLisandro Dalcin ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr); 8590f7a1166SHong Zhang if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 8600f7a1166SHong Zhang ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr); 8610f7a1166SHong Zhang } 862f68a32c8SEmil Constantinescu ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 863f68a32c8SEmil Constantinescu ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 864f68a32c8SEmil Constantinescu ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 86555de54fcSHong Zhang ierr = PetscTryMethod(ts,"TSSetUp_MRKSPLIT_C",(TS),(ts));CHKERRQ(ierr); 86655de54fcSHong Zhang ierr = PetscTryMethod(ts,"TSSetUp_MRKNONSPLIT_C",(TS),(ts));CHKERRQ(ierr); 867e4dd521cSBarry Smith PetscFunctionReturn(0); 868e4dd521cSBarry Smith } 869d2daff3dSHong Zhang 8702c9cb062Sluzhanghpp /* 87155de54fcSHong Zhang construct a database to choose single-step RK, or nonsplit version of mutirate RK or multirate RK with RHS splits 8722c9cb062Sluzhanghpp */ 87355de54fcSHong Zhang const char *const TSMRKTypes[] = {"NONE","NONSPLIT","SPLIT","TSMRKType","TSMRK",0}; 874e4dd521cSBarry Smith 8754416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts) 876e4dd521cSBarry Smith { 877be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 878dfbe8321SBarry Smith PetscErrorCode ierr; 879262ff7bbSBarry Smith 880e4dd521cSBarry Smith PetscFunctionBegin; 881e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr); 882f68a32c8SEmil Constantinescu { 883f68a32c8SEmil Constantinescu RKTableauLink link; 884f68a32c8SEmil Constantinescu PetscInt count,choice; 885f68a32c8SEmil Constantinescu PetscBool flg; 886f68a32c8SEmil Constantinescu const char **namelist; 88755de54fcSHong Zhang PetscInt mrktype = 0; 8882c9cb062Sluzhanghpp 889f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) ; 890f489ac74SBarry Smith ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr); 891f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 89255de54fcSHong Zhang ierr = PetscOptionsEList("-ts_rk_multirate_type","Use Multirate or partioned RHS Multirate or single rate RK method","TSRKSetMultirateType",TSMRKTypes,3,TSMRKTypes[0],&mrktype,&flg);CHKERRQ(ierr); 89355de54fcSHong Zhang if (flg) {ierr = TSRKSetMultirateType(ts,mrktype);CHKERRQ(ierr);} 894be5899b3SLisandro Dalcin ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr); 895be5899b3SLisandro Dalcin if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);} 896f68a32c8SEmil Constantinescu ierr = PetscFree(namelist);CHKERRQ(ierr); 897f68a32c8SEmil Constantinescu } 898262ff7bbSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 8992c9cb062Sluzhanghpp ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Multirate methods options","");CHKERRQ(ierr); 9002c9cb062Sluzhanghpp ierr = PetscOptionsInt("-ts_rk_dtratio","time step ratio between slow and fast","",rk->dtratio,&rk->dtratio,NULL);CHKERRQ(ierr); 9012c9cb062Sluzhanghpp ierr = PetscOptionsEnd();CHKERRQ(ierr); 902e4dd521cSBarry Smith PetscFunctionReturn(0); 903e4dd521cSBarry Smith } 904e4dd521cSBarry Smith 9055f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer) 906e4dd521cSBarry Smith { 9075f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 9088370ee3bSLisandro Dalcin PetscBool iascii; 909dfbe8321SBarry Smith PetscErrorCode ierr; 9102cdf8120Sasbjorn 911e4dd521cSBarry Smith PetscFunctionBegin; 912251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 9138370ee3bSLisandro Dalcin if (iascii) { 9149c334d8fSLisandro Dalcin RKTableau tab = rk->tableau; 915f68a32c8SEmil Constantinescu TSRKType rktype; 916f68a32c8SEmil Constantinescu char buf[512]; 917f68a32c8SEmil Constantinescu ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr); 918efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RK type %s\n",rktype);CHKERRQ(ierr); 9190c37eb8eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," Order: %D\n",tab->order);CHKERRQ(ierr); 9200c37eb8eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr); 921f68a32c8SEmil Constantinescu ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 922f68a32c8SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," Abscissa c = %s\n",buf);CHKERRQ(ierr); 9238370ee3bSLisandro Dalcin } 924f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 925f68a32c8SEmil Constantinescu } 926f68a32c8SEmil Constantinescu 927f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer) 928f68a32c8SEmil Constantinescu { 929f68a32c8SEmil Constantinescu PetscErrorCode ierr; 9309c334d8fSLisandro Dalcin TSAdapt adapt; 931f68a32c8SEmil Constantinescu 932f68a32c8SEmil Constantinescu PetscFunctionBegin; 9339c334d8fSLisandro Dalcin ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 9349c334d8fSLisandro Dalcin ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr); 935f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 936f68a32c8SEmil Constantinescu } 937f68a32c8SEmil Constantinescu 938f68a32c8SEmil Constantinescu /*@C 939f68a32c8SEmil Constantinescu TSRKSetType - Set the type of RK scheme 940f68a32c8SEmil Constantinescu 941f68a32c8SEmil Constantinescu Logically collective 942f68a32c8SEmil Constantinescu 943f68a32c8SEmil Constantinescu Input Parameter: 944f68a32c8SEmil Constantinescu + ts - timestepping context 945f68a32c8SEmil Constantinescu - rktype - type of RK-scheme 946f68a32c8SEmil Constantinescu 947287c2655SBarry Smith Options Database: 9489bd3e852SBarry Smith . -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs> 949287c2655SBarry Smith 950f68a32c8SEmil Constantinescu Level: intermediate 951f68a32c8SEmil Constantinescu 952287c2655SBarry Smith .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS 953f68a32c8SEmil Constantinescu @*/ 954f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype) 955f68a32c8SEmil Constantinescu { 956f68a32c8SEmil Constantinescu PetscErrorCode ierr; 957f68a32c8SEmil Constantinescu 958f68a32c8SEmil Constantinescu PetscFunctionBegin; 959f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 960b92453a8SLisandro Dalcin PetscValidCharPointer(rktype,2); 961f68a32c8SEmil Constantinescu ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr); 962f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 963f68a32c8SEmil Constantinescu } 964f68a32c8SEmil Constantinescu 965f68a32c8SEmil Constantinescu /*@C 966f68a32c8SEmil Constantinescu TSRKGetType - Get the type of RK scheme 967f68a32c8SEmil Constantinescu 968f68a32c8SEmil Constantinescu Logically collective 969f68a32c8SEmil Constantinescu 970f68a32c8SEmil Constantinescu Input Parameter: 971f68a32c8SEmil Constantinescu . ts - timestepping context 972f68a32c8SEmil Constantinescu 973f68a32c8SEmil Constantinescu Output Parameter: 974f68a32c8SEmil Constantinescu . rktype - type of RK-scheme 975f68a32c8SEmil Constantinescu 976f68a32c8SEmil Constantinescu Level: intermediate 977f68a32c8SEmil Constantinescu 978f68a32c8SEmil Constantinescu .seealso: TSRKGetType() 979f68a32c8SEmil Constantinescu @*/ 980f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype) 981f68a32c8SEmil Constantinescu { 982f68a32c8SEmil Constantinescu PetscErrorCode ierr; 983f68a32c8SEmil Constantinescu 984f68a32c8SEmil Constantinescu PetscFunctionBegin; 985f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 986f68a32c8SEmil Constantinescu ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr); 987f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 988f68a32c8SEmil Constantinescu } 989f68a32c8SEmil Constantinescu 990560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype) 991f68a32c8SEmil Constantinescu { 992f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 993f68a32c8SEmil Constantinescu 994f68a32c8SEmil Constantinescu PetscFunctionBegin; 995f68a32c8SEmil Constantinescu *rktype = rk->tableau->name; 996f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 997f68a32c8SEmil Constantinescu } 9982c9cb062Sluzhanghpp 999560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype) 1000f68a32c8SEmil Constantinescu { 1001f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 1002f68a32c8SEmil Constantinescu PetscErrorCode ierr; 1003f68a32c8SEmil Constantinescu PetscBool match; 1004f68a32c8SEmil Constantinescu RKTableauLink link; 1005f68a32c8SEmil Constantinescu 1006f68a32c8SEmil Constantinescu PetscFunctionBegin; 1007f68a32c8SEmil Constantinescu if (rk->tableau) { 1008f68a32c8SEmil Constantinescu ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr); 1009f68a32c8SEmil Constantinescu if (match) PetscFunctionReturn(0); 1010f68a32c8SEmil Constantinescu } 1011f68a32c8SEmil Constantinescu for (link = RKTableauList; link; link=link->next) { 1012f68a32c8SEmil Constantinescu ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr); 1013f68a32c8SEmil Constantinescu if (match) { 1014be5899b3SLisandro Dalcin if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);} 1015f68a32c8SEmil Constantinescu rk->tableau = &link->tab; 1016be5899b3SLisandro Dalcin if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);} 10172ffb9264SLisandro Dalcin ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 1018f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1019f68a32c8SEmil Constantinescu } 1020f68a32c8SEmil Constantinescu } 1021f68a32c8SEmil Constantinescu SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype); 1022e4dd521cSBarry Smith PetscFunctionReturn(0); 1023e4dd521cSBarry Smith } 1024e4dd521cSBarry Smith 1025ff22ae23SHong Zhang static PetscErrorCode TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y) 1026ff22ae23SHong Zhang { 1027ff22ae23SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 1028ff22ae23SHong Zhang 1029ff22ae23SHong Zhang PetscFunctionBegin; 10300429704eSStefano Zampini if (ns) *ns = rk->tableau->s; 1031d2daff3dSHong Zhang if (Y) *Y = rk->Y; 1032ff22ae23SHong Zhang PetscFunctionReturn(0); 1033ff22ae23SHong Zhang } 1034ff22ae23SHong Zhang 1035b3a6b972SJed Brown static PetscErrorCode TSDestroy_RK(TS ts) 1036b3a6b972SJed Brown { 1037b3a6b972SJed Brown PetscErrorCode ierr; 1038b3a6b972SJed Brown 1039b3a6b972SJed Brown PetscFunctionBegin; 1040b3a6b972SJed Brown ierr = TSReset_RK(ts);CHKERRQ(ierr); 1041b3a6b972SJed Brown if (ts->dm) { 1042b3a6b972SJed Brown ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 1043b3a6b972SJed Brown ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 1044b3a6b972SJed Brown } 1045b3a6b972SJed Brown ierr = PetscFree(ts->data);CHKERRQ(ierr); 1046b3a6b972SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr); 1047b3a6b972SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr); 10482c9cb062Sluzhanghpp ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirateType_C",NULL);CHKERRQ(ierr); 1049b3a6b972SJed Brown PetscFunctionReturn(0); 1050b3a6b972SJed Brown } 1051ff22ae23SHong Zhang 1052ebe3b25bSBarry Smith /*MC 1053f68a32c8SEmil Constantinescu TSRK - ODE and DAE solver using Runge-Kutta schemes 1054ebe3b25bSBarry Smith 1055f68a32c8SEmil Constantinescu The user should provide the right hand side of the equation 1056f68a32c8SEmil Constantinescu using TSSetRHSFunction(). 1057ebe3b25bSBarry Smith 1058f68a32c8SEmil Constantinescu Notes: 105998164e13SLisandro Dalcin The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type 1060ebe3b25bSBarry Smith 1061d41222bbSBarry Smith Level: beginner 1062d41222bbSBarry Smith 1063f68a32c8SEmil Constantinescu .seealso: TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3, 10642c9cb062Sluzhanghpp TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister(), TSRKSetMultirateType() 1065ebe3b25bSBarry Smith 1066ebe3b25bSBarry Smith M*/ 10678cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts) 1068e4dd521cSBarry Smith { 10691566a47fSLisandro Dalcin TS_RK *rk; 1070dfbe8321SBarry Smith PetscErrorCode ierr; 1071e4dd521cSBarry Smith 1072e4dd521cSBarry Smith PetscFunctionBegin; 1073f68a32c8SEmil Constantinescu ierr = TSRKInitializePackage();CHKERRQ(ierr); 1074f68a32c8SEmil Constantinescu 1075f68a32c8SEmil Constantinescu ts->ops->reset = TSReset_RK; 10765f70b478SJed Brown ts->ops->destroy = TSDestroy_RK; 10775f70b478SJed Brown ts->ops->view = TSView_RK; 1078f68a32c8SEmil Constantinescu ts->ops->load = TSLoad_RK; 1079f68a32c8SEmil Constantinescu ts->ops->setup = TSSetUp_RK; 108042f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_RK; 1081f68a32c8SEmil Constantinescu ts->ops->interpolate = TSInterpolate_RK; 10822c9cb062Sluzhanghpp ts->ops->step = TSStep_RK; 1083f68a32c8SEmil Constantinescu ts->ops->evaluatestep = TSEvaluateStep_RK; 1084fecfb714SLisandro Dalcin ts->ops->rollback = TSRollBack_RK; 1085f68a32c8SEmil Constantinescu ts->ops->setfromoptions = TSSetFromOptions_RK; 1086ff22ae23SHong Zhang ts->ops->getstages = TSGetStages_RK; 108742f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_RK; 1088e4dd521cSBarry Smith 10890f7a1166SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_RK; 10900f7a1166SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_RK; 10910f7a1166SHong Zhang 10921566a47fSLisandro Dalcin ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr); 10931566a47fSLisandro Dalcin ts->data = (void*)rk; 1094e4dd521cSBarry Smith 1095f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr); 1096f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr); 10972c9cb062Sluzhanghpp ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirateType_C",TSRKSetMultirateType);CHKERRQ(ierr); 1098be5899b3SLisandro Dalcin 1099be5899b3SLisandro Dalcin ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr); 1100*90b67129SHong Zhang rk->dtratio = 1; 1101e4dd521cSBarry Smith PetscFunctionReturn(0); 1102e4dd521cSBarry Smith } 110355de54fcSHong Zhang 110455de54fcSHong Zhang 110555de54fcSHong Zhang /*------------ Multirate support for partitioned systems --------------*/ 110655de54fcSHong Zhang 110755de54fcSHong Zhang static PetscErrorCode TSSetUp_MRKSPLIT(TS ts) 110855de54fcSHong Zhang { 110955de54fcSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 111055de54fcSHong Zhang RKTableau tab = rk->tableau; 111155de54fcSHong Zhang DM dm,subdm,newdm; 111255de54fcSHong Zhang PetscErrorCode ierr; 111355de54fcSHong Zhang 111455de54fcSHong Zhang PetscFunctionBegin; 111555de54fcSHong Zhang ierr = TSRHSSplitGetIS(ts,"slow",&rk->is_slow);CHKERRQ(ierr); 111655de54fcSHong Zhang ierr = TSRHSSplitGetIS(ts,"fast",&rk->is_fast);CHKERRQ(ierr); 111755de54fcSHong Zhang if (!rk->is_slow || !rk->is_fast) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set up RHSSplits with TSRHSSplitSetIS() using split names 'slow' and 'fast' respectively in order to use -ts_type bsi"); 111855de54fcSHong Zhang ierr = TSRHSSplitGetSubTS(ts,"slow",&rk->subts_slow);CHKERRQ(ierr); 111955de54fcSHong Zhang ierr = TSRHSSplitGetSubTS(ts,"fast",&rk->subts_fast);CHKERRQ(ierr); 112055de54fcSHong Zhang if (!rk->subts_slow || !rk->subts_fast) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set up the RHSFunctions for 'slow' and 'fast' components using TSRHSSplitSetRHSFunction() or calling TSSetRHSFunction() for each sub-TS"); 112155de54fcSHong Zhang 112255de54fcSHong Zhang /* Only copy */ 112355de54fcSHong Zhang ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 112455de54fcSHong Zhang ierr = DMClone(dm,&newdm);CHKERRQ(ierr); 112555de54fcSHong Zhang ierr = TSGetDM(rk->subts_fast,&subdm);CHKERRQ(ierr); 112655de54fcSHong Zhang ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr); 112755de54fcSHong Zhang ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr); 112855de54fcSHong Zhang ierr = TSSetDM(rk->subts_fast,newdm);CHKERRQ(ierr); 112955de54fcSHong Zhang ierr = DMDestroy(&newdm);CHKERRQ(ierr); 113055de54fcSHong Zhang ierr = DMClone(dm,&newdm);CHKERRQ(ierr); 113155de54fcSHong Zhang ierr = TSGetDM(rk->subts_slow,&subdm);CHKERRQ(ierr); 113255de54fcSHong Zhang ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr); 113355de54fcSHong Zhang ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr); 113455de54fcSHong Zhang ierr = TSSetDM(rk->subts_slow,newdm);CHKERRQ(ierr); 113555de54fcSHong Zhang ierr = DMDestroy(&newdm);CHKERRQ(ierr); 113655de54fcSHong Zhang 113755de54fcSHong Zhang ierr = PetscMalloc1(tab->s,&rk->YdotRHS_fast);CHKERRQ(ierr); 113855de54fcSHong Zhang ierr = PetscMalloc1(tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr); 113955de54fcSHong Zhang ierr = VecDuplicate(ts->vec_sol,&rk->X0);CHKERRQ(ierr); 114055de54fcSHong Zhang PetscFunctionReturn(0); 114155de54fcSHong Zhang } 114255de54fcSHong Zhang 114355de54fcSHong Zhang static PetscErrorCode TSReset_MRKSPLIT(TS ts) 114455de54fcSHong Zhang { 114555de54fcSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 114655de54fcSHong Zhang PetscErrorCode ierr; 114755de54fcSHong Zhang 114855de54fcSHong Zhang PetscFunctionBegin; 1149*90b67129SHong Zhang ierr = PetscFree(rk->YdotRHS_fast);CHKERRQ(ierr); 1150*90b67129SHong Zhang ierr = PetscFree(rk->YdotRHS_slow);CHKERRQ(ierr); 115155de54fcSHong Zhang ierr = VecDestroy(&rk->X0);CHKERRQ(ierr); 115255de54fcSHong Zhang PetscFunctionReturn(0); 115355de54fcSHong Zhang } 115455de54fcSHong Zhang 115555de54fcSHong Zhang static PetscErrorCode TSSetUp_MRKNONSPLIT(TS ts) 115655de54fcSHong Zhang { 115755de54fcSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 115855de54fcSHong Zhang RKTableau tab = rk->tableau; 115955de54fcSHong Zhang PetscErrorCode ierr; 116055de54fcSHong Zhang 116155de54fcSHong Zhang PetscFunctionBegin; 116255de54fcSHong Zhang ierr = TSRHSSplitGetIS(ts,"slow",&rk->is_slow);CHKERRQ(ierr); 116355de54fcSHong Zhang ierr = TSRHSSplitGetIS(ts,"fast",&rk->is_fast);CHKERRQ(ierr); 116455de54fcSHong Zhang if (!rk->is_slow || !rk->is_fast) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set up RHSSplits with TSRHSSplitSetIS() using split names 'slow' and 'fast' respectively in order to use -ts_type bsi"); 116555de54fcSHong Zhang ierr = VecDuplicate(ts->vec_sol,&rk->X0);CHKERRQ(ierr); 116655de54fcSHong Zhang ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr); 116755de54fcSHong Zhang PetscFunctionReturn(0); 116855de54fcSHong Zhang } 116955de54fcSHong Zhang 117055de54fcSHong Zhang static PetscErrorCode TSReset_MRKNONSPLIT(TS ts) 117155de54fcSHong Zhang { 117255de54fcSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 117355de54fcSHong Zhang RKTableau tab = rk->tableau; 117455de54fcSHong Zhang PetscErrorCode ierr; 117555de54fcSHong Zhang 117655de54fcSHong Zhang PetscFunctionBegin; 117755de54fcSHong Zhang ierr = VecDestroy(&rk->X0);CHKERRQ(ierr); 117855de54fcSHong Zhang ierr = VecDestroyVecs(tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr); 117955de54fcSHong Zhang PetscFunctionReturn(0); 118055de54fcSHong Zhang } 118155de54fcSHong Zhang 118255de54fcSHong Zhang static PetscErrorCode TSInterpolate_MRKNONSPLIT(TS ts,PetscReal itime,Vec X) 118355de54fcSHong Zhang { 118455de54fcSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 118555de54fcSHong Zhang PetscInt s = rk->tableau->s,p = rk->tableau->p,i,j; 118655de54fcSHong Zhang PetscReal h; 118755de54fcSHong Zhang PetscReal tt,t; 118855de54fcSHong Zhang PetscScalar *b; 118955de54fcSHong Zhang const PetscReal *B = rk->tableau->binterp; 119055de54fcSHong Zhang PetscErrorCode ierr; 119155de54fcSHong Zhang 119255de54fcSHong Zhang PetscFunctionBegin; 119355de54fcSHong Zhang if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 119455de54fcSHong Zhang 119555de54fcSHong Zhang switch (rk->status) { 119655de54fcSHong Zhang case TS_STEP_INCOMPLETE: 119755de54fcSHong Zhang case TS_STEP_PENDING: 119855de54fcSHong Zhang h = ts->time_step; 119955de54fcSHong Zhang t = (itime - ts->ptime)/h; 120055de54fcSHong Zhang break; 120155de54fcSHong Zhang case TS_STEP_COMPLETE: 120255de54fcSHong Zhang h = ts->ptime - ts->ptime_prev; 120355de54fcSHong Zhang t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 120455de54fcSHong Zhang break; 120555de54fcSHong Zhang default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 120655de54fcSHong Zhang } 120755de54fcSHong Zhang ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 120855de54fcSHong Zhang for (i=0; i<s; i++) b[i] = 0; 120955de54fcSHong Zhang for (j=0,tt=t; j<p; j++,tt*=t) { 121055de54fcSHong Zhang for (i=0; i<s; i++) { 121155de54fcSHong Zhang b[i] += h * B[i*p+j] * tt; 121255de54fcSHong Zhang } 121355de54fcSHong Zhang } 121455de54fcSHong Zhang ierr = VecCopy(rk->X0,X);CHKERRQ(ierr); 121555de54fcSHong Zhang ierr = VecMAXPY(X,s,b,rk->YdotRHS_slow);CHKERRQ(ierr); 121655de54fcSHong Zhang ierr = PetscFree(b);CHKERRQ(ierr); 121755de54fcSHong Zhang PetscFunctionReturn(0); 121855de54fcSHong Zhang } 121955de54fcSHong Zhang 122055de54fcSHong Zhang 122155de54fcSHong Zhang static PetscErrorCode TSInterpolate_MRKSPLIT(TS ts,PetscReal itime,Vec X) 122255de54fcSHong Zhang { 122355de54fcSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 122455de54fcSHong Zhang PetscInt s = rk->tableau->s,p = rk->tableau->p,i,j; 122555de54fcSHong Zhang Vec Yslow; /* vector holds the slow components in Y[0] */ 122655de54fcSHong Zhang PetscReal h; 122755de54fcSHong Zhang PetscReal tt,t; 122855de54fcSHong Zhang PetscScalar *b; 122955de54fcSHong Zhang const PetscReal *B = rk->tableau->binterp; 123055de54fcSHong Zhang PetscErrorCode ierr; 123155de54fcSHong Zhang 123255de54fcSHong Zhang PetscFunctionBegin; 123355de54fcSHong Zhang if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 123455de54fcSHong Zhang 123555de54fcSHong Zhang switch (rk->status) { 123655de54fcSHong Zhang case TS_STEP_INCOMPLETE: 123755de54fcSHong Zhang case TS_STEP_PENDING: 123855de54fcSHong Zhang h = ts->time_step; 123955de54fcSHong Zhang t = (itime - ts->ptime)/h; 124055de54fcSHong Zhang break; 124155de54fcSHong Zhang case TS_STEP_COMPLETE: 124255de54fcSHong Zhang h = ts->ptime - ts->ptime_prev; 124355de54fcSHong Zhang t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 124455de54fcSHong Zhang break; 124555de54fcSHong Zhang default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 124655de54fcSHong Zhang } 124755de54fcSHong Zhang ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 124855de54fcSHong Zhang for (i=0; i<s; i++) b[i] = 0; 124955de54fcSHong Zhang for (j=0,tt=t; j<p; j++,tt*=t) { 125055de54fcSHong Zhang for (i=0; i<s; i++) { 125155de54fcSHong Zhang b[i] += h * B[i*p+j] * tt; 125255de54fcSHong Zhang } 125355de54fcSHong Zhang } 125455de54fcSHong Zhang ierr = VecGetSubVector(rk->X0,rk->is_slow,&Yslow);CHKERRQ(ierr); 125555de54fcSHong Zhang ierr = VecCopy(Yslow,X);CHKERRQ(ierr); 125655de54fcSHong Zhang ierr = VecMAXPY(X,s,b,rk->YdotRHS_slow);CHKERRQ(ierr); 125755de54fcSHong Zhang ierr = VecRestoreSubVector(rk->X0,rk->is_slow,&Yslow);CHKERRQ(ierr); 125855de54fcSHong Zhang ierr = PetscFree(b);CHKERRQ(ierr); 125955de54fcSHong Zhang PetscFunctionReturn(0); 126055de54fcSHong Zhang } 126155de54fcSHong Zhang 126255de54fcSHong Zhang static PetscErrorCode TSStep_MRKNONSPLIT(TS ts) 126355de54fcSHong Zhang { 126455de54fcSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 126555de54fcSHong Zhang RKTableau tab = rk->tableau; 126655de54fcSHong Zhang Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS,*YdotRHS_slow = rk->YdotRHS_slow; 1267*90b67129SHong Zhang Vec stage_slow,sol_slow; /* vectors store the slow components */ 1268*90b67129SHong Zhang Vec subvec_slow; /* sub vector to store the slow components */ 126955de54fcSHong Zhang const PetscInt s = tab->s; 127055de54fcSHong Zhang const PetscReal *A = tab->A,*c = tab->c; 127155de54fcSHong Zhang PetscScalar *w = rk->work; 1272*90b67129SHong Zhang PetscInt i,j,k,dtratio = rk->dtratio; 127355de54fcSHong Zhang PetscReal next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step; 127455de54fcSHong Zhang PetscErrorCode ierr; 127555de54fcSHong Zhang 127655de54fcSHong Zhang PetscFunctionBegin; 127755de54fcSHong Zhang rk->status = TS_STEP_INCOMPLETE; 127855de54fcSHong Zhang ierr = VecDuplicate(ts->vec_sol,&stage_slow);CHKERRQ(ierr); 1279*90b67129SHong Zhang ierr = VecDuplicate(ts->vec_sol,&sol_slow);CHKERRQ(ierr); 128055de54fcSHong Zhang ierr = VecCopy(ts->vec_sol,rk->X0);CHKERRQ(ierr); 128155de54fcSHong Zhang for (i=0; i<s; i++) { 128255de54fcSHong Zhang rk->stage_time = t + h*c[i]; 128355de54fcSHong Zhang ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr); 128455de54fcSHong Zhang ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 128555de54fcSHong Zhang for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 128655de54fcSHong Zhang ierr = VecMAXPY(Y[i],i,w,YdotRHS_slow);CHKERRQ(ierr); 128755de54fcSHong Zhang ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 128855de54fcSHong Zhang /* compute the stage RHS */ 128955de54fcSHong Zhang ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr); 129055de54fcSHong Zhang } 1291*90b67129SHong Zhang /* update the slow components in the solution */ 1292*90b67129SHong Zhang rk->YdotRHS = YdotRHS_slow; 1293*90b67129SHong Zhang rk->dtratio = 1; 1294*90b67129SHong Zhang ierr = TSEvaluateStep(ts,tab->order,sol_slow,NULL);CHKERRQ(ierr); 1295*90b67129SHong Zhang rk->dtratio = dtratio; 1296*90b67129SHong Zhang rk->YdotRHS = YdotRHS; 129755de54fcSHong Zhang for (k=0; k<rk->dtratio; k++) { 129855de54fcSHong Zhang for (i=0; i<s; i++) { 129955de54fcSHong Zhang rk->stage_time = t + h/rk->dtratio*c[i]; 130055de54fcSHong Zhang ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr); 1301*90b67129SHong Zhang /* update the fast components in the stage value, the slow components will be overwritten, so it is ok to have garbage in the slow components */ 1302*90b67129SHong Zhang ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 130355de54fcSHong Zhang for (j=0; j<i; j++) w[j] = h/rk->dtratio*A[i*s+j]; 130455de54fcSHong Zhang ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 130555de54fcSHong Zhang ierr = TSInterpolate_MRKNONSPLIT(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],stage_slow);CHKERRQ(ierr); 1306*90b67129SHong Zhang /* update the slow components in the stage value */ 1307*90b67129SHong Zhang ierr = VecGetSubVector(stage_slow,rk->is_slow,&subvec_slow);CHKERRQ(ierr); 1308*90b67129SHong Zhang ierr = VecISCopy(Y[i],rk->is_slow,SCATTER_FORWARD,subvec_slow);CHKERRQ(ierr); 1309*90b67129SHong Zhang ierr = VecRestoreSubVector(stage_slow,rk->is_slow,&subvec_slow);CHKERRQ(ierr); 131055de54fcSHong Zhang ierr = TSPostStage(ts,rk->stage_time,i,Y);CHKERRQ(ierr); 131155de54fcSHong Zhang /* compute the stage RHS */ 131255de54fcSHong Zhang ierr = TSComputeRHSFunction(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 131355de54fcSHong Zhang } 1314*90b67129SHong Zhang ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 131555de54fcSHong Zhang } 1316*90b67129SHong Zhang /* update the slow components in the solution */ 1317*90b67129SHong Zhang ierr = VecGetSubVector(sol_slow,rk->is_slow,&subvec_slow);CHKERRQ(ierr); 1318*90b67129SHong Zhang ierr = VecISCopy(ts->vec_sol,rk->is_slow,SCATTER_FORWARD,subvec_slow);CHKERRQ(ierr); 1319*90b67129SHong Zhang ierr = VecRestoreSubVector(sol_slow,rk->is_slow,&subvec_slow);CHKERRQ(ierr); 1320*90b67129SHong Zhang 132155de54fcSHong Zhang ts->ptime += ts->time_step; 132255de54fcSHong Zhang ts->time_step = next_time_step; 132355de54fcSHong Zhang rk->status = TS_STEP_COMPLETE; 132455de54fcSHong Zhang /* free memory of work vectors */ 1325*90b67129SHong Zhang ierr = VecDestroy(&stage_slow);CHKERRQ(ierr); 1326*90b67129SHong Zhang ierr = VecDestroy(&sol_slow);CHKERRQ(ierr); 132755de54fcSHong Zhang PetscFunctionReturn(0); 132855de54fcSHong Zhang } 132955de54fcSHong Zhang 133055de54fcSHong Zhang /* 133155de54fcSHong Zhang This is for partitioned RHS multirate RK method 133255de54fcSHong Zhang The step completion formula is 133355de54fcSHong Zhang 133455de54fcSHong Zhang x1 = x0 + h b^T YdotRHS 133555de54fcSHong Zhang 133655de54fcSHong Zhang */ 133755de54fcSHong Zhang static PetscErrorCode TSEvaluateStep_MRKSPLIT(TS ts,PetscInt order,Vec X,PetscBool *done) 133855de54fcSHong Zhang { 133955de54fcSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 134055de54fcSHong Zhang RKTableau tab = rk->tableau; 134155de54fcSHong Zhang Vec Xslow,Xfast; /* subvectors of X which store slow components and fast components respectively */ 134255de54fcSHong Zhang PetscScalar *w = rk->work; 134355de54fcSHong Zhang PetscReal h = ts->time_step; 134455de54fcSHong Zhang PetscInt s = tab->s,j; 134555de54fcSHong Zhang PetscErrorCode ierr; 134655de54fcSHong Zhang 134755de54fcSHong Zhang PetscFunctionBegin; 134855de54fcSHong Zhang ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 134955de54fcSHong Zhang if (rk->slow) { 135055de54fcSHong Zhang for (j=0; j<s; j++) w[j] = h*tab->b[j]; 135155de54fcSHong Zhang ierr = VecGetSubVector(ts->vec_sol,rk->is_slow,&Xslow);CHKERRQ(ierr); 135255de54fcSHong Zhang ierr = VecMAXPY(Xslow,s,w,rk->YdotRHS_slow);CHKERRQ(ierr); 135355de54fcSHong Zhang ierr = VecRestoreSubVector(ts->vec_sol,rk->is_slow,&Xslow);CHKERRQ(ierr);; 135455de54fcSHong Zhang } else { 135555de54fcSHong Zhang for (j=0; j<s; j++) w[j] = h/rk->dtratio*tab->b[j]; 135655de54fcSHong Zhang ierr = VecGetSubVector(X,rk->is_fast,&Xfast);CHKERRQ(ierr); 135755de54fcSHong Zhang ierr = VecMAXPY(Xfast,s,w,rk->YdotRHS_fast);CHKERRQ(ierr); 135855de54fcSHong Zhang ierr = VecRestoreSubVector(X,rk->is_fast,&Xfast);CHKERRQ(ierr); 135955de54fcSHong Zhang } 136055de54fcSHong Zhang PetscFunctionReturn(0); 136155de54fcSHong Zhang } 136255de54fcSHong Zhang 136355de54fcSHong Zhang static PetscErrorCode TSStep_MRKSPLIT(TS ts) 136455de54fcSHong Zhang { 136555de54fcSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 136655de54fcSHong Zhang RKTableau tab = rk->tableau; 136755de54fcSHong Zhang Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 136855de54fcSHong Zhang Vec *YdotRHS_fast = rk->YdotRHS_fast,*YdotRHS_slow = rk->YdotRHS_slow; 136955de54fcSHong Zhang Vec Yslow,Yfast; /* subvectors store the stges of slow components and fast components respectively */ 137055de54fcSHong Zhang const PetscInt s = tab->s; 137155de54fcSHong Zhang const PetscReal *A = tab->A,*c = tab->c; 137255de54fcSHong Zhang PetscScalar *w = rk->work; 137355de54fcSHong Zhang PetscInt i,j,k; 137455de54fcSHong Zhang PetscReal next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step; 137555de54fcSHong Zhang PetscErrorCode ierr; 137655de54fcSHong Zhang 137755de54fcSHong Zhang PetscFunctionBegin; 137855de54fcSHong Zhang rk->status = TS_STEP_INCOMPLETE; 137955de54fcSHong Zhang for (i=0; i<s; i++) { 138055de54fcSHong Zhang ierr = VecGetSubVector(YdotRHS[i],rk->is_slow,&YdotRHS_slow[i]);CHKERRQ(ierr); 138155de54fcSHong Zhang ierr = VecGetSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr); 138255de54fcSHong Zhang } 138355de54fcSHong Zhang ierr = VecCopy(ts->vec_sol,rk->X0);CHKERRQ(ierr); 138455de54fcSHong Zhang /* propogate both slow and fast components using large time steps */ 138555de54fcSHong Zhang for (i=0; i<s; i++) { 138655de54fcSHong Zhang rk->stage_time = t + h*c[i]; 138755de54fcSHong Zhang ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr); 138855de54fcSHong Zhang ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 138955de54fcSHong Zhang /*ierr = VecCopy(ts->vec_sol,Yc[i]);CHKERRQ(ierr);*/ 139055de54fcSHong Zhang ierr = VecGetSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr); 139155de54fcSHong Zhang ierr = VecGetSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr); 139255de54fcSHong Zhang for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 139355de54fcSHong Zhang ierr = VecMAXPY(Yslow,i,w,YdotRHS_slow);CHKERRQ(ierr); 139455de54fcSHong Zhang ierr = VecMAXPY(Yfast,i,w,YdotRHS_fast);CHKERRQ(ierr); 139555de54fcSHong Zhang ierr = VecRestoreSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr); 139655de54fcSHong Zhang ierr = VecRestoreSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr); 139755de54fcSHong Zhang ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 139855de54fcSHong Zhang ierr = TSComputeRHSFunction(rk->subts_slow,t+h*c[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr); 139955de54fcSHong Zhang ierr = TSComputeRHSFunction(rk->subts_fast,t+h*c[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr); 140055de54fcSHong Zhang } 140155de54fcSHong Zhang rk->slow = PETSC_TRUE; 140255de54fcSHong Zhang ierr = TSEvaluateStep_MRKSPLIT(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 140355de54fcSHong Zhang for (k=0; k<rk->dtratio; k++) { 140455de54fcSHong Zhang /* propogate fast component using small time steps */ 140555de54fcSHong Zhang for (i=0; i<s; i++) { 140655de54fcSHong Zhang rk->stage_time = t + h/rk->dtratio*c[i]; 140755de54fcSHong Zhang ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr); 140855de54fcSHong Zhang ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 140955de54fcSHong Zhang /* stage value for fast components */ 141055de54fcSHong Zhang for (j=0; j<i; j++) w[j] = h/rk->dtratio*A[i*s+j]; 141155de54fcSHong Zhang ierr = VecGetSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr); 141255de54fcSHong Zhang ierr = VecMAXPY(Yfast,i,w,YdotRHS_fast);CHKERRQ(ierr); 141355de54fcSHong Zhang ierr = VecRestoreSubVector(Y[i],rk->is_fast,&Yfast);CHKERRQ(ierr); 141455de54fcSHong Zhang /* stage value for slow components */ 141555de54fcSHong Zhang ierr = VecGetSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr); 141655de54fcSHong Zhang ierr = TSInterpolate_MRKSPLIT(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Yslow);CHKERRQ(ierr); 141755de54fcSHong Zhang ierr = VecRestoreSubVector(Y[i],rk->is_slow,&Yslow);CHKERRQ(ierr); 141855de54fcSHong Zhang ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 141955de54fcSHong Zhang /* compute the stage RHS for fast components */ 142055de54fcSHong Zhang ierr = TSComputeRHSFunction(rk->subts_fast,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr); 142155de54fcSHong Zhang } 142255de54fcSHong Zhang /* update the value of fast components when using fast time step */ 142355de54fcSHong Zhang rk->slow = PETSC_FALSE; 142455de54fcSHong Zhang ierr = TSEvaluateStep_MRKSPLIT(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 142555de54fcSHong Zhang } 142655de54fcSHong Zhang for (i=0; i<s; i++) { 142755de54fcSHong Zhang ierr = VecRestoreSubVector(YdotRHS[i],rk->is_slow,&YdotRHS_slow[i]);CHKERRQ(ierr); 142855de54fcSHong Zhang ierr = VecRestoreSubVector(YdotRHS[i],rk->is_fast,&YdotRHS_fast[i]);CHKERRQ(ierr); 142955de54fcSHong Zhang } 143055de54fcSHong Zhang ts->ptime += ts->time_step; 143155de54fcSHong Zhang ts->time_step = next_time_step; 143255de54fcSHong Zhang rk->status = TS_STEP_COMPLETE; 143355de54fcSHong Zhang PetscFunctionReturn(0); 143455de54fcSHong Zhang } 143555de54fcSHong Zhang 143655de54fcSHong Zhang /*@C 143755de54fcSHong Zhang TSRKSetMultirateType - Set the type of RK Multirate scheme 143855de54fcSHong Zhang 143955de54fcSHong Zhang Logically collective 144055de54fcSHong Zhang 144155de54fcSHong Zhang Input Parameter: 144255de54fcSHong Zhang + ts - timestepping context 144355de54fcSHong Zhang - mrktype - type of MRK-scheme 144455de54fcSHong Zhang 144555de54fcSHong Zhang Options Database: 144655de54fcSHong Zhang . -ts_rk_multiarte_type - <none,nonsplit,split> 144755de54fcSHong Zhang 144855de54fcSHong Zhang Level: intermediate 144955de54fcSHong Zhang @*/ 145055de54fcSHong Zhang PetscErrorCode TSRKSetMultirateType(TS ts, TSMRKType mrktype) 145155de54fcSHong Zhang { 145255de54fcSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 145355de54fcSHong Zhang PetscErrorCode ierr; 145455de54fcSHong Zhang 145555de54fcSHong Zhang PetscFunctionBegin; 145655de54fcSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 145755de54fcSHong Zhang switch(mrktype){ 145855de54fcSHong Zhang case TSMRKNONE: 145955de54fcSHong Zhang break; 146055de54fcSHong Zhang case TSMRKNONSPLIT: 146155de54fcSHong Zhang ts->ops->step = TSStep_MRKNONSPLIT; 146255de54fcSHong Zhang ts->ops->interpolate = TSInterpolate_MRKNONSPLIT; 146355de54fcSHong Zhang rk->dtratio = 2; 146455de54fcSHong Zhang ierr = TSRKSetType(ts,TSMRKDefault);CHKERRQ(ierr); 146555de54fcSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_MRKNONSPLIT_C",TSSetUp_MRKNONSPLIT);CHKERRQ(ierr); 146655de54fcSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_MRKNONSPLIT_C",TSReset_MRKNONSPLIT);CHKERRQ(ierr); 146755de54fcSHong Zhang break; 146855de54fcSHong Zhang case TSMRKSPLIT: 146955de54fcSHong Zhang ts->ops->step = TSStep_MRKSPLIT; 147055de54fcSHong Zhang ts->ops->evaluatestep = TSEvaluateStep_MRKSPLIT; 147155de54fcSHong Zhang ts->ops->interpolate = TSInterpolate_MRKSPLIT; 147255de54fcSHong Zhang rk->dtratio = 2; 147355de54fcSHong Zhang ierr = TSRKSetType(ts,TSMRKDefault);CHKERRQ(ierr); 147455de54fcSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSetUp_MRKSPLIT_C",TSSetUp_MRKSPLIT);CHKERRQ(ierr); 147555de54fcSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSReset_MRKSPLIT_C",TSReset_MRKSPLIT);CHKERRQ(ierr); 147655de54fcSHong Zhang break; 147755de54fcSHong Zhang default : 147855de54fcSHong Zhang SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown type '%s'",mrktype); 147955de54fcSHong Zhang } 148055de54fcSHong Zhang PetscFunctionReturn(0); 148155de54fcSHong Zhang } 1482