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 82c9cb062Sluzhanghpp Udot = F(t,U) for combined RHS 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) 122c9cb062Sluzhanghpp Ufdot = Ff(t,Us,Uf) for partitioned RHS multi-rate RK, 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> 18f68a32c8SEmil Constantinescu 19484bcdc7SDebojyoti Ghosh static TSRKType TSRKDefault = TSRK3BS; 202c9cb062Sluzhanghpp static TSRKType TSMRKDefault = TSRK2A; 21f68a32c8SEmil Constantinescu static PetscBool TSRKRegisterAllCalled; 22f68a32c8SEmil Constantinescu static PetscBool TSRKPackageInitialized; 23f68a32c8SEmil Constantinescu 24f68a32c8SEmil Constantinescu typedef struct _RKTableau *RKTableau; 25f68a32c8SEmil Constantinescu struct _RKTableau { 26f68a32c8SEmil Constantinescu char *name; 27d760c35bSDebojyoti Ghosh PetscInt order; /* Classical approximation order of the method i */ 28f68a32c8SEmil Constantinescu PetscInt s; /* Number of stages */ 293a8a9803SLisandro Dalcin PetscInt p; /* Interpolation order */ 30d760c35bSDebojyoti Ghosh PetscBool FSAL; /* flag to indicate if tableau is FSAL */ 31f68a32c8SEmil Constantinescu PetscReal *A,*b,*c; /* Tableau */ 32f68a32c8SEmil Constantinescu PetscReal *bembed; /* Embedded formula of order one less (order-1) */ 33f68a32c8SEmil Constantinescu PetscReal *binterp; /* Dense output formula */ 34f68a32c8SEmil Constantinescu PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */ 35f68a32c8SEmil Constantinescu }; 36f68a32c8SEmil Constantinescu typedef struct _RKTableauLink *RKTableauLink; 37f68a32c8SEmil Constantinescu struct _RKTableauLink { 38f68a32c8SEmil Constantinescu struct _RKTableau tab; 39f68a32c8SEmil Constantinescu RKTableauLink next; 40f68a32c8SEmil Constantinescu }; 41f68a32c8SEmil Constantinescu static RKTableauLink RKTableauList; 42e4dd521cSBarry Smith 43e4dd521cSBarry Smith typedef struct { 44f68a32c8SEmil Constantinescu RKTableau tableau; 45f68a32c8SEmil Constantinescu Vec *Y; /* States computed during the step */ 462c9cb062Sluzhanghpp Vec *YdotRHS; /* Function evaluations for the non-stiff part and contains all components */ 472c9cb062Sluzhanghpp Vec *YdotRHS_fast; /* Function evaluations for the non-stiff part and contains fast components */ 482c9cb062Sluzhanghpp Vec *YdotRHS_slow; /* Function evaluations for the non-stiff part and contains slow components */ 49ad8e2604SHong Zhang Vec *VecDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage */ 50ad8e2604SHong Zhang Vec *VecDeltaMu; /* Increment of the adjoint sensitivity w.r.t P at stage */ 510f7a1166SHong Zhang Vec VecCostIntegral0; /* backup for roll-backs due to events */ 52f68a32c8SEmil Constantinescu PetscScalar *work; /* Scalar work */ 532c9cb062Sluzhanghpp PetscInt slow; /* flag indicates call slow components solver (0) or fast components solver (1) */ 54f68a32c8SEmil Constantinescu PetscReal stage_time; 55f68a32c8SEmil Constantinescu TSStepStatus status; 560f7a1166SHong Zhang PetscReal ptime; 570f7a1166SHong Zhang PetscReal time_step; 582c9cb062Sluzhanghpp PetscInt dtratio; /* ratio between slow time step size and fast step size */ 595f70b478SJed Brown } TS_RK; 60e4dd521cSBarry Smith 612c9cb062Sluzhanghpp 62f68a32c8SEmil Constantinescu /*MC 63287c2655SBarry Smith TSRK1FE - First order forward Euler scheme. 64262ff7bbSBarry Smith 65f68a32c8SEmil Constantinescu This method has one stage. 66f68a32c8SEmil Constantinescu 67287c2655SBarry Smith Options database: 689bd3e852SBarry Smith . -ts_rk_type 1fe 69287c2655SBarry Smith 70f68a32c8SEmil Constantinescu Level: advanced 71f68a32c8SEmil Constantinescu 72287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 73f68a32c8SEmil Constantinescu M*/ 74f68a32c8SEmil Constantinescu /*MC 752109b73fSDebojyoti Ghosh TSRK2A - Second order RK scheme. 76f68a32c8SEmil Constantinescu 77f68a32c8SEmil Constantinescu This method has two stages. 78f68a32c8SEmil Constantinescu 79287c2655SBarry Smith Options database: 809bd3e852SBarry Smith . -ts_rk_type 2a 81287c2655SBarry Smith 82f68a32c8SEmil Constantinescu Level: advanced 83f68a32c8SEmil Constantinescu 84287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 85f68a32c8SEmil Constantinescu M*/ 86f68a32c8SEmil Constantinescu /*MC 87f68a32c8SEmil Constantinescu TSRK3 - Third order RK scheme. 88f68a32c8SEmil Constantinescu 89f68a32c8SEmil Constantinescu This method has three stages. 90f68a32c8SEmil Constantinescu 91287c2655SBarry Smith Options database: 929bd3e852SBarry Smith . -ts_rk_type 3 93287c2655SBarry Smith 94f68a32c8SEmil Constantinescu Level: advanced 95f68a32c8SEmil Constantinescu 96287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 97f68a32c8SEmil Constantinescu M*/ 98f68a32c8SEmil Constantinescu /*MC 992109b73fSDebojyoti Ghosh TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method. 1002109b73fSDebojyoti Ghosh 101d1905564SLisandro Dalcin This method has four stages with the First Same As Last (FSAL) property. 1022109b73fSDebojyoti Ghosh 103287c2655SBarry Smith Options database: 1049bd3e852SBarry Smith . -ts_rk_type 3bs 105287c2655SBarry Smith 1062109b73fSDebojyoti Ghosh Level: advanced 1072109b73fSDebojyoti Ghosh 10898164e13SLisandro Dalcin References: https://doi.org/10.1016/0893-9659(89)90079-7 109d1905564SLisandro Dalcin 110287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 1112109b73fSDebojyoti Ghosh M*/ 1122109b73fSDebojyoti Ghosh /*MC 113f68a32c8SEmil Constantinescu TSRK4 - Fourth order RK scheme. 114f68a32c8SEmil Constantinescu 1152e698f8bSDebojyoti Ghosh This is the classical Runge-Kutta method with four stages. 116f68a32c8SEmil Constantinescu 117287c2655SBarry Smith Options database: 1189bd3e852SBarry Smith . -ts_rk_type 4 119287c2655SBarry Smith 120f68a32c8SEmil Constantinescu Level: advanced 121f68a32c8SEmil Constantinescu 122287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 123f68a32c8SEmil Constantinescu M*/ 124f68a32c8SEmil Constantinescu /*MC 1252e698f8bSDebojyoti Ghosh TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method. 126f68a32c8SEmil Constantinescu 127f68a32c8SEmil Constantinescu This method has six stages. 128f68a32c8SEmil Constantinescu 129287c2655SBarry Smith Options database: 1309bd3e852SBarry Smith . -ts_rk_type 5f 131287c2655SBarry Smith 132f68a32c8SEmil Constantinescu Level: advanced 133f68a32c8SEmil Constantinescu 134287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 135f68a32c8SEmil Constantinescu M*/ 1362109b73fSDebojyoti Ghosh /*MC 1372e698f8bSDebojyoti Ghosh TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method. 1382109b73fSDebojyoti Ghosh 139d1905564SLisandro Dalcin This method has seven stages with the First Same As Last (FSAL) property. 1402109b73fSDebojyoti Ghosh 141287c2655SBarry Smith Options database: 1429bd3e852SBarry Smith . -ts_rk_type 5dp 143287c2655SBarry Smith 1442109b73fSDebojyoti Ghosh Level: advanced 1452109b73fSDebojyoti Ghosh 14698164e13SLisandro Dalcin References: https://doi.org/10.1016/0771-050X(80)90013-3 147d1905564SLisandro Dalcin 148287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 1492109b73fSDebojyoti Ghosh M*/ 15005e23783SLisandro Dalcin /*MC 15105e23783SLisandro Dalcin TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method. 15205e23783SLisandro Dalcin 15305e23783SLisandro Dalcin This method has eight stages with the First Same As Last (FSAL) property. 15405e23783SLisandro Dalcin 155287c2655SBarry Smith Options database: 1569bd3e852SBarry Smith . -ts_rk_type 5bs 157287c2655SBarry Smith 15805e23783SLisandro Dalcin Level: advanced 15905e23783SLisandro Dalcin 16098164e13SLisandro Dalcin References: https://doi.org/10.1016/0898-1221(96)00141-1 16105e23783SLisandro Dalcin 162287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 16305e23783SLisandro Dalcin M*/ 164262ff7bbSBarry Smith 165f68a32c8SEmil Constantinescu /*@C 166f68a32c8SEmil Constantinescu TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK 167262ff7bbSBarry Smith 168f68a32c8SEmil Constantinescu Not Collective, but should be called by all processes which will need the schemes to be registered 169262ff7bbSBarry Smith 170f68a32c8SEmil Constantinescu Level: advanced 171262ff7bbSBarry Smith 172f68a32c8SEmil Constantinescu .keywords: TS, TSRK, register, all 173262ff7bbSBarry Smith 174f68a32c8SEmil Constantinescu .seealso: TSRKRegisterDestroy() 175262ff7bbSBarry Smith @*/ 176f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void) 177262ff7bbSBarry Smith { 1784ac538c5SBarry Smith PetscErrorCode ierr; 179262ff7bbSBarry Smith 180262ff7bbSBarry Smith PetscFunctionBegin; 181f68a32c8SEmil Constantinescu if (TSRKRegisterAllCalled) PetscFunctionReturn(0); 182f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_TRUE; 183e4dd521cSBarry Smith 1844e82626cSLisandro Dalcin #define RC PetscRealConstant 185e4dd521cSBarry Smith { 186f68a32c8SEmil Constantinescu const PetscReal 1874e82626cSLisandro Dalcin A[1][1] = {{0}}, 1884e82626cSLisandro Dalcin b[1] = {RC(1.0)}; 1893a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 190e4dd521cSBarry Smith } 191db046809SBarry Smith { 192f68a32c8SEmil Constantinescu const PetscReal 1934e82626cSLisandro Dalcin A[2][2] = {{0,0}, 1944e82626cSLisandro Dalcin {RC(1.0),0}}, 1954e82626cSLisandro Dalcin b[2] = {RC(0.5),RC(0.5)}, 1964e82626cSLisandro Dalcin bembed[2] = {RC(1.0),0}; 1973a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 198db046809SBarry Smith } 199f68a32c8SEmil Constantinescu { 200f68a32c8SEmil Constantinescu const PetscReal 20117f6384fSBarry Smith A[3][3] = {{0,0,0}, 2024e82626cSLisandro Dalcin {RC(2.0)/RC(3.0),0,0}, 2034e82626cSLisandro Dalcin {RC(-1.0)/RC(3.0),RC(1.0),0}}, 2044e82626cSLisandro Dalcin b[3] = {RC(0.25),RC(0.5),RC(0.25)}; 2053a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 206fdefd5e5SDebojyoti Ghosh } 207fdefd5e5SDebojyoti Ghosh { 208fdefd5e5SDebojyoti Ghosh const PetscReal 20917f6384fSBarry Smith A[4][4] = {{0,0,0,0}, 2104e82626cSLisandro Dalcin {RC(1.0)/RC(2.0),0,0,0}, 2114e82626cSLisandro Dalcin {0,RC(3.0)/RC(4.0),0,0}, 2124e82626cSLisandro Dalcin {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}}, 2134e82626cSLisandro Dalcin b[4] = {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}, 2144e82626cSLisandro 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)}; 2153a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 216db046809SBarry Smith } 217f68a32c8SEmil Constantinescu { 218f68a32c8SEmil Constantinescu const PetscReal 219f68a32c8SEmil Constantinescu A[4][4] = {{0,0,0,0}, 2204e82626cSLisandro Dalcin {RC(0.5),0,0,0}, 2214e82626cSLisandro Dalcin {0,RC(0.5),0,0}, 2224e82626cSLisandro Dalcin {0,0,RC(1.0),0}}, 2234e82626cSLisandro 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)}; 2243a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 225db046809SBarry Smith } 226f68a32c8SEmil Constantinescu { 227f68a32c8SEmil Constantinescu const PetscReal 22817f6384fSBarry Smith A[6][6] = {{0,0,0,0,0,0}, 2294e82626cSLisandro Dalcin {RC(0.25),0,0,0,0,0}, 2304e82626cSLisandro Dalcin {RC(3.0)/RC(32.0),RC(9.0)/RC(32.0),0,0,0,0}, 2314e82626cSLisandro Dalcin {RC(1932.0)/RC(2197.0),RC(-7200.0)/RC(2197.0),RC(7296.0)/RC(2197.0),0,0,0}, 2324e82626cSLisandro Dalcin {RC(439.0)/RC(216.0),RC(-8.0),RC(3680.0)/RC(513.0),RC(-845.0)/RC(4104.0),0,0}, 2334e82626cSLisandro 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}}, 2344e82626cSLisandro 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)}, 2354e82626cSLisandro 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}; 2363a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 237fdefd5e5SDebojyoti Ghosh } 238fdefd5e5SDebojyoti Ghosh { 239fdefd5e5SDebojyoti Ghosh const PetscReal 24017f6384fSBarry Smith A[7][7] = {{0,0,0,0,0,0,0}, 2414e82626cSLisandro Dalcin {RC(1.0)/RC(5.0),0,0,0,0,0,0}, 2424e82626cSLisandro Dalcin {RC(3.0)/RC(40.0),RC(9.0)/RC(40.0),0,0,0,0,0}, 2434e82626cSLisandro Dalcin {RC(44.0)/RC(45.0),RC(-56.0)/RC(15.0),RC(32.0)/RC(9.0),0,0,0,0}, 2444e82626cSLisandro 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}, 2454e82626cSLisandro 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}, 2464e82626cSLisandro 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}}, 2474e82626cSLisandro 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}, 248*a685a763Sluzhanghpp 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)}, 249*a685a763Sluzhanghpp 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)}, 250*a685a763Sluzhanghpp {0,0,0,0,0}, 251*a685a763Sluzhanghpp {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)}, 252*a685a763Sluzhanghpp {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)}, 253*a685a763Sluzhanghpp {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)}, 254*a685a763Sluzhanghpp {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)}, 255*a685a763Sluzhanghpp {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)}}; 256*a685a763Sluzhanghpp 257*a685a763Sluzhanghpp ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,5,binterp[0]);CHKERRQ(ierr); 258f68a32c8SEmil Constantinescu } 25905e23783SLisandro Dalcin { 26005e23783SLisandro Dalcin const PetscReal 26117f6384fSBarry Smith A[8][8] = {{0,0,0,0,0,0,0,0}, 2624e82626cSLisandro Dalcin {RC(1.0)/RC(6.0),0,0,0,0,0,0,0}, 2634e82626cSLisandro Dalcin {RC(2.0)/RC(27.0),RC(4.0)/RC(27.0),0,0,0,0,0,0}, 2644e82626cSLisandro Dalcin {RC(183.0)/RC(1372.0),RC(-162.0)/RC(343.0),RC(1053.0)/RC(1372.0),0,0,0,0,0}, 2654e82626cSLisandro 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}, 2664e82626cSLisandro 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}, 2674e82626cSLisandro 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}, 2684e82626cSLisandro 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}}, 2694e82626cSLisandro 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}, 2704e82626cSLisandro 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)}; 27105e23783SLisandro Dalcin ierr = TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 27205e23783SLisandro Dalcin } 2734e82626cSLisandro Dalcin #undef RC 274db046809SBarry Smith PetscFunctionReturn(0); 275db046809SBarry Smith } 276db046809SBarry Smith 277f68a32c8SEmil Constantinescu /*@C 278f68a32c8SEmil Constantinescu TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister(). 279f68a32c8SEmil Constantinescu 280f68a32c8SEmil Constantinescu Not Collective 281f68a32c8SEmil Constantinescu 282f68a32c8SEmil Constantinescu Level: advanced 283f68a32c8SEmil Constantinescu 284f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy 285f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll() 286f68a32c8SEmil Constantinescu @*/ 287f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void) 288e4dd521cSBarry Smith { 289f68a32c8SEmil Constantinescu PetscErrorCode ierr; 290f68a32c8SEmil Constantinescu RKTableauLink link; 291f68a32c8SEmil Constantinescu 292f68a32c8SEmil Constantinescu PetscFunctionBegin; 293f68a32c8SEmil Constantinescu while ((link = RKTableauList)) { 294f68a32c8SEmil Constantinescu RKTableau t = &link->tab; 295f68a32c8SEmil Constantinescu RKTableauList = link->next; 296f68a32c8SEmil Constantinescu ierr = PetscFree3(t->A,t->b,t->c); CHKERRQ(ierr); 297f68a32c8SEmil Constantinescu ierr = PetscFree (t->bembed); CHKERRQ(ierr); 298f68a32c8SEmil Constantinescu ierr = PetscFree (t->binterp); CHKERRQ(ierr); 299f68a32c8SEmil Constantinescu ierr = PetscFree (t->name); CHKERRQ(ierr); 300f68a32c8SEmil Constantinescu ierr = PetscFree (link); CHKERRQ(ierr); 301f68a32c8SEmil Constantinescu } 302f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_FALSE; 303f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 304f68a32c8SEmil Constantinescu } 305f68a32c8SEmil Constantinescu 306f68a32c8SEmil Constantinescu /*@C 307f68a32c8SEmil Constantinescu TSRKInitializePackage - This function initializes everything in the TSRK package. It is called 3088a690491SBarry Smith from TSInitializePackage(). 309f68a32c8SEmil Constantinescu 310f68a32c8SEmil Constantinescu Level: developer 311f68a32c8SEmil Constantinescu 312f68a32c8SEmil Constantinescu .keywords: TS, TSRK, initialize, package 313f68a32c8SEmil Constantinescu .seealso: PetscInitialize() 314f68a32c8SEmil Constantinescu @*/ 315f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void) 316f68a32c8SEmil Constantinescu { 317186e87acSLisandro Dalcin PetscErrorCode ierr; 318e4dd521cSBarry Smith 319e4dd521cSBarry Smith PetscFunctionBegin; 320f68a32c8SEmil Constantinescu if (TSRKPackageInitialized) PetscFunctionReturn(0); 321f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_TRUE; 322f68a32c8SEmil Constantinescu ierr = TSRKRegisterAll();CHKERRQ(ierr); 323f68a32c8SEmil Constantinescu ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr); 324f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 325f68a32c8SEmil Constantinescu } 326186e87acSLisandro Dalcin 327f68a32c8SEmil Constantinescu /*@C 328f68a32c8SEmil Constantinescu TSRKFinalizePackage - This function destroys everything in the TSRK package. It is 329f68a32c8SEmil Constantinescu called from PetscFinalize(). 330186e87acSLisandro Dalcin 331f68a32c8SEmil Constantinescu Level: developer 332f68a32c8SEmil Constantinescu 333f68a32c8SEmil Constantinescu .keywords: Petsc, destroy, package 334f68a32c8SEmil Constantinescu .seealso: PetscFinalize() 335f68a32c8SEmil Constantinescu @*/ 336f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void) 337f68a32c8SEmil Constantinescu { 338f68a32c8SEmil Constantinescu PetscErrorCode ierr; 339f68a32c8SEmil Constantinescu 340f68a32c8SEmil Constantinescu PetscFunctionBegin; 341f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_FALSE; 342f68a32c8SEmil Constantinescu ierr = TSRKRegisterDestroy();CHKERRQ(ierr); 343f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 344f68a32c8SEmil Constantinescu } 345f68a32c8SEmil Constantinescu 346f68a32c8SEmil Constantinescu /*@C 347f68a32c8SEmil Constantinescu TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 348f68a32c8SEmil Constantinescu 349f68a32c8SEmil Constantinescu Not Collective, but the same schemes should be registered on all processes on which they will be used 350f68a32c8SEmil Constantinescu 351f68a32c8SEmil Constantinescu Input Parameters: 352f68a32c8SEmil Constantinescu + name - identifier for method 353f68a32c8SEmil Constantinescu . order - approximation order of method 354f68a32c8SEmil Constantinescu . s - number of stages, this is the dimension of the matrices below 355f68a32c8SEmil Constantinescu . A - stage coefficients (dimension s*s, row-major) 356f68a32c8SEmil Constantinescu . b - step completion table (dimension s; NULL to use last row of A) 357f68a32c8SEmil Constantinescu . c - abscissa (dimension s; NULL to use row sums of A) 358f68a32c8SEmil Constantinescu . bembed - completion table for embedded method (dimension s; NULL if not available) 3593a8a9803SLisandro Dalcin . p - Order of the interpolation scheme, equal to the number of columns of binterp 3603a8a9803SLisandro Dalcin - binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1) 361f68a32c8SEmil Constantinescu 362f68a32c8SEmil Constantinescu Notes: 363f68a32c8SEmil Constantinescu Several RK methods are provided, this function is only needed to create new methods. 364f68a32c8SEmil Constantinescu 365f68a32c8SEmil Constantinescu Level: advanced 366f68a32c8SEmil Constantinescu 367f68a32c8SEmil Constantinescu .keywords: TS, register 368f68a32c8SEmil Constantinescu 369f68a32c8SEmil Constantinescu .seealso: TSRK 370f68a32c8SEmil Constantinescu @*/ 371f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s, 372f68a32c8SEmil Constantinescu const PetscReal A[],const PetscReal b[],const PetscReal c[], 3733a8a9803SLisandro Dalcin const PetscReal bembed[],PetscInt p,const PetscReal binterp[]) 374f68a32c8SEmil Constantinescu { 375f68a32c8SEmil Constantinescu PetscErrorCode ierr; 376f68a32c8SEmil Constantinescu RKTableauLink link; 377f68a32c8SEmil Constantinescu RKTableau t; 378f68a32c8SEmil Constantinescu PetscInt i,j; 379f68a32c8SEmil Constantinescu 380f68a32c8SEmil Constantinescu PetscFunctionBegin; 3813a8a9803SLisandro Dalcin PetscValidCharPointer(name,1); 3823a8a9803SLisandro Dalcin PetscValidRealPointer(A,4); 3833a8a9803SLisandro Dalcin if (b) PetscValidRealPointer(b,5); 3843a8a9803SLisandro Dalcin if (c) PetscValidRealPointer(c,6); 3853a8a9803SLisandro Dalcin if (bembed) PetscValidRealPointer(bembed,7); 3863a8a9803SLisandro Dalcin if (binterp || p > 1) PetscValidRealPointer(binterp,9); 3873a8a9803SLisandro Dalcin 3881d36bdfdSBarry Smith ierr = TSRKInitializePackage();CHKERRQ(ierr); 38995dccacaSBarry Smith ierr = PetscNew(&link);CHKERRQ(ierr); 390f68a32c8SEmil Constantinescu t = &link->tab; 3913a8a9803SLisandro Dalcin 392f68a32c8SEmil Constantinescu ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 393f68a32c8SEmil Constantinescu t->order = order; 394f68a32c8SEmil Constantinescu t->s = s; 395dcca6d9dSJed Brown ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr); 396f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 397f68a32c8SEmil Constantinescu if (b) { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); } 398f68a32c8SEmil Constantinescu else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i]; 399f68a32c8SEmil Constantinescu if (c) { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); } 400f68a32c8SEmil 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]; 401d760c35bSDebojyoti Ghosh t->FSAL = PETSC_TRUE; 402d760c35bSDebojyoti Ghosh for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE; 4033a8a9803SLisandro Dalcin 404f68a32c8SEmil Constantinescu if (bembed) { 405785e854fSJed Brown ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr); 406f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr); 407f68a32c8SEmil Constantinescu } 408f68a32c8SEmil Constantinescu 4093a8a9803SLisandro Dalcin if (!binterp) { p = 1; binterp = t->b; } 4103a8a9803SLisandro Dalcin t->p = p; 4113a8a9803SLisandro Dalcin ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr); 4123a8a9803SLisandro Dalcin ierr = PetscMemcpy(t->binterp,binterp,s*p*sizeof(binterp[0]));CHKERRQ(ierr); 4133a8a9803SLisandro Dalcin 414f68a32c8SEmil Constantinescu link->next = RKTableauList; 415f68a32c8SEmil Constantinescu RKTableauList = link; 416f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 417f68a32c8SEmil Constantinescu } 418f68a32c8SEmil Constantinescu 419e4dd521cSBarry Smith /* 4202c9cb062Sluzhanghpp This is for single-step RK method 421f68a32c8SEmil Constantinescu The step completion formula is 422e4dd521cSBarry Smith 423f68a32c8SEmil Constantinescu x1 = x0 + h b^T YdotRHS 424f68a32c8SEmil Constantinescu 425f68a32c8SEmil Constantinescu This function can be called before or after ts->vec_sol has been updated. 426f68a32c8SEmil Constantinescu Suppose we have a completion formula (b) and an embedded formula (be) of different order. 427f68a32c8SEmil Constantinescu We can write 428f68a32c8SEmil Constantinescu 429f68a32c8SEmil Constantinescu x1e = x0 + h be^T YdotRHS 430f68a32c8SEmil Constantinescu = x1 - h b^T YdotRHS + h be^T YdotRHS 431f68a32c8SEmil Constantinescu = x1 + h (be - b)^T YdotRHS 432f68a32c8SEmil Constantinescu 433f68a32c8SEmil Constantinescu so we can evaluate the method with different order even after the step has been optimistically completed. 434e4dd521cSBarry Smith */ 435f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done) 436f68a32c8SEmil Constantinescu { 437f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 438f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 439f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 440f68a32c8SEmil Constantinescu PetscReal h; 441f68a32c8SEmil Constantinescu PetscInt s = tab->s,j; 442f68a32c8SEmil Constantinescu PetscErrorCode ierr; 443f68a32c8SEmil Constantinescu 444f68a32c8SEmil Constantinescu PetscFunctionBegin; 445f68a32c8SEmil Constantinescu switch (rk->status) { 446f68a32c8SEmil Constantinescu case TS_STEP_INCOMPLETE: 447f68a32c8SEmil Constantinescu case TS_STEP_PENDING: 448f68a32c8SEmil Constantinescu h = ts->time_step; break; 449f68a32c8SEmil Constantinescu case TS_STEP_COMPLETE: 450be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 451f68a32c8SEmil Constantinescu default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 452f68a32c8SEmil Constantinescu } 453f68a32c8SEmil Constantinescu if (order == tab->order) { 454f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { 455f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 456f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->b[j]; 457f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 458f68a32c8SEmil Constantinescu } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 459f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 460f68a32c8SEmil Constantinescu } else if (order == tab->order-1) { 461f68a32c8SEmil Constantinescu if (!tab->bembed) goto unavailable; 462f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { /*Complete with the embedded method (be)*/ 463f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 464f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 465f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 466f68a32c8SEmil Constantinescu } else { /*Rollback and re-complete using (be-b) */ 467f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 468f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 469f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 470f68a32c8SEmil Constantinescu } 471f68a32c8SEmil Constantinescu if (done) *done = PETSC_TRUE; 472f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 473f68a32c8SEmil Constantinescu } 474f68a32c8SEmil Constantinescu unavailable: 475f68a32c8SEmil Constantinescu if (done) *done = PETSC_FALSE; 476a7fac7c2SEmil 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); 477f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 478f68a32c8SEmil Constantinescu } 479f68a32c8SEmil Constantinescu 4800f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts) 4810f7a1166SHong Zhang { 4820f7a1166SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 4830f7a1166SHong Zhang RKTableau tab = rk->tableau; 4840f7a1166SHong Zhang const PetscInt s = tab->s; 4850f7a1166SHong Zhang const PetscReal *b = tab->b,*c = tab->c; 4860f7a1166SHong Zhang Vec *Y = rk->Y; 4870f7a1166SHong Zhang PetscInt i; 4880f7a1166SHong Zhang PetscErrorCode ierr; 4890f7a1166SHong Zhang 4900f7a1166SHong Zhang PetscFunctionBegin; 4910f7a1166SHong Zhang /* backup cost integral */ 4920f7a1166SHong Zhang ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr); 4930f7a1166SHong Zhang for (i=s-1; i>=0; i--) { 4940f7a1166SHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 49551a7f651SHong Zhang ierr = TSComputeCostIntegrand(ts,rk->ptime+rk->time_step*c[i],Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 4960f7a1166SHong Zhang ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 4970f7a1166SHong Zhang } 4980f7a1166SHong Zhang PetscFunctionReturn(0); 4990f7a1166SHong Zhang } 5000f7a1166SHong Zhang 5010f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts) 5020f7a1166SHong Zhang { 5030f7a1166SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 5040f7a1166SHong Zhang RKTableau tab = rk->tableau; 5050f7a1166SHong Zhang const PetscInt s = tab->s; 5060f7a1166SHong Zhang const PetscReal *b = tab->b,*c = tab->c; 5070f7a1166SHong Zhang Vec *Y = rk->Y; 5080f7a1166SHong Zhang PetscInt i; 5090f7a1166SHong Zhang PetscErrorCode ierr; 5100f7a1166SHong Zhang 5110f7a1166SHong Zhang PetscFunctionBegin; 5120f7a1166SHong Zhang for (i=s-1; i>=0; i--) { 5130f7a1166SHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 51451a7f651SHong Zhang ierr = TSComputeCostIntegrand(ts,ts->ptime+ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 5150f7a1166SHong Zhang ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 5160f7a1166SHong Zhang } 5170f7a1166SHong Zhang PetscFunctionReturn(0); 5180f7a1166SHong Zhang } 5190f7a1166SHong Zhang 520fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts) 521fecfb714SLisandro Dalcin { 522fecfb714SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 523fecfb714SLisandro Dalcin RKTableau tab = rk->tableau; 524fecfb714SLisandro Dalcin const PetscInt s = tab->s; 525fecfb714SLisandro Dalcin const PetscReal *b = tab->b; 526fecfb714SLisandro Dalcin PetscScalar *w = rk->work; 527fecfb714SLisandro Dalcin Vec *YdotRHS = rk->YdotRHS; 528fecfb714SLisandro Dalcin PetscInt j; 529fecfb714SLisandro Dalcin PetscReal h; 530fecfb714SLisandro Dalcin PetscErrorCode ierr; 531fecfb714SLisandro Dalcin 532fecfb714SLisandro Dalcin PetscFunctionBegin; 533fecfb714SLisandro Dalcin switch (rk->status) { 534fecfb714SLisandro Dalcin case TS_STEP_INCOMPLETE: 535fecfb714SLisandro Dalcin case TS_STEP_PENDING: 536fecfb714SLisandro Dalcin h = ts->time_step; break; 537fecfb714SLisandro Dalcin case TS_STEP_COMPLETE: 538fecfb714SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 539fecfb714SLisandro Dalcin default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 540fecfb714SLisandro Dalcin } 541fecfb714SLisandro Dalcin for (j=0; j<s; j++) w[j] = -h*b[j]; 542fecfb714SLisandro Dalcin ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr); 543fecfb714SLisandro Dalcin PetscFunctionReturn(0); 544fecfb714SLisandro Dalcin } 545fecfb714SLisandro Dalcin 546f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts) 547f68a32c8SEmil Constantinescu { 548f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 549f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 550f68a32c8SEmil Constantinescu const PetscInt s = tab->s; 551fecfb714SLisandro Dalcin const PetscReal *A = tab->A,*c = tab->c; 552f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 553f68a32c8SEmil Constantinescu Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 554d1905564SLisandro Dalcin PetscBool FSAL = tab->FSAL; 555f68a32c8SEmil Constantinescu TSAdapt adapt; 556fecfb714SLisandro Dalcin PetscInt i,j; 557be5899b3SLisandro Dalcin PetscInt rejections = 0; 558be5899b3SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 559be5899b3SLisandro Dalcin PetscReal next_time_step = ts->time_step; 560f68a32c8SEmil Constantinescu PetscErrorCode ierr; 561f68a32c8SEmil Constantinescu 562f68a32c8SEmil Constantinescu PetscFunctionBegin; 563d1905564SLisandro Dalcin if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE; 564d1905564SLisandro Dalcin if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); } 565d1905564SLisandro Dalcin 566f68a32c8SEmil Constantinescu rk->status = TS_STEP_INCOMPLETE; 567be5899b3SLisandro Dalcin while (!ts->reason && rk->status != TS_STEP_COMPLETE) { 568be5899b3SLisandro Dalcin PetscReal t = ts->ptime; 569f68a32c8SEmil Constantinescu PetscReal h = ts->time_step; 570f68a32c8SEmil Constantinescu for (i=0; i<s; i++) { 5719be3e283SDebojyoti Ghosh rk->stage_time = t + h*c[i]; 5729be3e283SDebojyoti Ghosh ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr); 573f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 574f68a32c8SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 575f68a32c8SEmil Constantinescu ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 5769be3e283SDebojyoti Ghosh ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 577f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 578be5899b3SLisandro Dalcin ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr); 579fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 5808f738a7cSLisandro Dalcin if (FSAL && !i) continue; 581f68a32c8SEmil Constantinescu ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 582f68a32c8SEmil Constantinescu } 583be5899b3SLisandro Dalcin 584be5899b3SLisandro Dalcin rk->status = TS_STEP_INCOMPLETE; 585f68a32c8SEmil Constantinescu ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 586f68a32c8SEmil Constantinescu rk->status = TS_STEP_PENDING; 587f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 588f68a32c8SEmil Constantinescu ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 5891917a363SLisandro Dalcin ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr); 590fecfb714SLisandro Dalcin ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 591be5899b3SLisandro Dalcin rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 592be5899b3SLisandro Dalcin if (!accept) { /*Roll back the current step */ 593fecfb714SLisandro Dalcin ierr = TSRollBack_RK(ts);CHKERRQ(ierr); 594be5899b3SLisandro Dalcin ts->time_step = next_time_step; 595be5899b3SLisandro Dalcin goto reject_step; 596be5899b3SLisandro Dalcin } 597be5899b3SLisandro Dalcin 5980f7a1166SHong Zhang if (ts->costintegralfwd) { /*Save the info for the later use in cost integral evaluation*/ 5990f7a1166SHong Zhang rk->ptime = ts->ptime; 6000f7a1166SHong Zhang rk->time_step = ts->time_step; 601493ed95fSHong Zhang } 602be5899b3SLisandro Dalcin 603f68a32c8SEmil Constantinescu ts->ptime += ts->time_step; 604f68a32c8SEmil Constantinescu ts->time_step = next_time_step; 605f68a32c8SEmil Constantinescu break; 606be5899b3SLisandro Dalcin 607be5899b3SLisandro Dalcin reject_step: 608fecfb714SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 609be5899b3SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 610be5899b3SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 611be5899b3SLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 612f68a32c8SEmil Constantinescu } 613f68a32c8SEmil Constantinescu } 614f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 615e4dd521cSBarry Smith } 616e4dd521cSBarry Smith 6172c9cb062Sluzhanghpp static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X) 6182c9cb062Sluzhanghpp { 6192c9cb062Sluzhanghpp TS_RK *rk = (TS_RK*)ts->data; 6202c9cb062Sluzhanghpp PetscInt s = rk->tableau->s,p = rk->tableau->p,i,j; 6212c9cb062Sluzhanghpp PetscReal h; 6222c9cb062Sluzhanghpp PetscReal tt,t; 6232c9cb062Sluzhanghpp PetscScalar *b; 6242c9cb062Sluzhanghpp const PetscReal *B = rk->tableau->binterp; 6252c9cb062Sluzhanghpp PetscErrorCode ierr; 6262c9cb062Sluzhanghpp 6272c9cb062Sluzhanghpp PetscFunctionBegin; 6282c9cb062Sluzhanghpp if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 6292c9cb062Sluzhanghpp 6302c9cb062Sluzhanghpp switch (rk->status) { 6312c9cb062Sluzhanghpp case TS_STEP_INCOMPLETE: 6322c9cb062Sluzhanghpp case TS_STEP_PENDING: 6332c9cb062Sluzhanghpp h = ts->time_step; 6342c9cb062Sluzhanghpp t = (itime - ts->ptime)/h; 6352c9cb062Sluzhanghpp break; 6362c9cb062Sluzhanghpp case TS_STEP_COMPLETE: 6372c9cb062Sluzhanghpp h = ts->ptime - ts->ptime_prev; 6382c9cb062Sluzhanghpp t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 6392c9cb062Sluzhanghpp break; 6402c9cb062Sluzhanghpp default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 6412c9cb062Sluzhanghpp } 6422c9cb062Sluzhanghpp ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 6432c9cb062Sluzhanghpp for (i=0; i<s; i++) b[i] = 0; 6442c9cb062Sluzhanghpp for (j=0,tt=t; j<p; j++,tt*=t) { 6452c9cb062Sluzhanghpp for (i=0; i<s; i++) { 6462c9cb062Sluzhanghpp b[i] += h * B[i*p+j] * tt; 6472c9cb062Sluzhanghpp } 6482c9cb062Sluzhanghpp } 6492c9cb062Sluzhanghpp ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr); 6502c9cb062Sluzhanghpp ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr); 6512c9cb062Sluzhanghpp ierr = PetscFree(b);CHKERRQ(ierr); 6522c9cb062Sluzhanghpp PetscFunctionReturn(0); 6532c9cb062Sluzhanghpp } 6542c9cb062Sluzhanghpp 6552c9cb062Sluzhanghpp /* 6562c9cb062Sluzhanghpp This is interpolate formula for partitioned RHS multirate RK method 6572c9cb062Sluzhanghpp */ 6582c9cb062Sluzhanghpp 6592c9cb062Sluzhanghpp static PetscErrorCode TSInterpolate_PARTITIONEDMRK(TS ts,PetscReal itime,Vec X) 6602c9cb062Sluzhanghpp { 6612c9cb062Sluzhanghpp TS_RK *rk = (TS_RK*)ts->data; 6622c9cb062Sluzhanghpp PetscInt s = rk->tableau->s,p = rk->tableau->p,i,j; 6632c9cb062Sluzhanghpp Vec Yslow; /* vector holds the slow components in Y[0] */ 6642c9cb062Sluzhanghpp PetscReal h; 6652c9cb062Sluzhanghpp PetscReal tt,t; 6662c9cb062Sluzhanghpp PetscScalar *b; 6672c9cb062Sluzhanghpp const PetscReal *B = rk->tableau->binterp; 6682c9cb062Sluzhanghpp PetscErrorCode ierr; 6692c9cb062Sluzhanghpp 6702c9cb062Sluzhanghpp PetscFunctionBegin; 6712c9cb062Sluzhanghpp if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 6722c9cb062Sluzhanghpp 6732c9cb062Sluzhanghpp switch (rk->status) { 6742c9cb062Sluzhanghpp case TS_STEP_INCOMPLETE: 6752c9cb062Sluzhanghpp case TS_STEP_PENDING: 6762c9cb062Sluzhanghpp h = ts->time_step; 6772c9cb062Sluzhanghpp t = (itime - ts->ptime)/h; 6782c9cb062Sluzhanghpp break; 6792c9cb062Sluzhanghpp case TS_STEP_COMPLETE: 6802c9cb062Sluzhanghpp h = ts->ptime - ts->ptime_prev; 6812c9cb062Sluzhanghpp t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 6822c9cb062Sluzhanghpp break; 6832c9cb062Sluzhanghpp default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 6842c9cb062Sluzhanghpp } 6852c9cb062Sluzhanghpp ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 6862c9cb062Sluzhanghpp for (i=0; i<s; i++) b[i] = 0; 6872c9cb062Sluzhanghpp for (j=0,tt=t; j<p; j++,tt*=t) { 6882c9cb062Sluzhanghpp for (i=0; i<s; i++) { 6892c9cb062Sluzhanghpp b[i] += h * B[i*p+j] * tt; 6902c9cb062Sluzhanghpp } 6912c9cb062Sluzhanghpp } 6922c9cb062Sluzhanghpp ierr = VecGetSubVector(rk->Y[0],ts->iss,&Yslow);CHKERRQ(ierr); 6932c9cb062Sluzhanghpp ierr = VecCopy(Yslow,X);CHKERRQ(ierr); 6942c9cb062Sluzhanghpp ierr = VecMAXPY(X,s,b,rk->YdotRHS_slow);CHKERRQ(ierr); 6952c9cb062Sluzhanghpp ierr = VecRestoreSubVector(rk->Y[0],ts->iss,&Yslow);CHKERRQ(ierr); 6962c9cb062Sluzhanghpp ierr = PetscFree(b);CHKERRQ(ierr); 6972c9cb062Sluzhanghpp PetscFunctionReturn(0); 6982c9cb062Sluzhanghpp } 6992c9cb062Sluzhanghpp 7002c9cb062Sluzhanghpp /* 7012c9cb062Sluzhanghpp This is for combined RHS multirate RK method 7022c9cb062Sluzhanghpp The step completion formula is 7032c9cb062Sluzhanghpp 7042c9cb062Sluzhanghpp x1 = x0 + h b^T YdotRHS 7052c9cb062Sluzhanghpp 7062c9cb062Sluzhanghpp */ 7072c9cb062Sluzhanghpp static PetscErrorCode TSEvaluateStep_RKMULTIRATE(TS ts,PetscInt order,Vec X,PetscBool *done) 7082c9cb062Sluzhanghpp { 7092c9cb062Sluzhanghpp TS_RK *rk = (TS_RK*)ts->data; 7102c9cb062Sluzhanghpp RKTableau tab = rk->tableau; 7112c9cb062Sluzhanghpp Vec Xtemp; /* temporary work vector for X */ 7122c9cb062Sluzhanghpp PetscScalar *w = rk->work; 7132c9cb062Sluzhanghpp PetscScalar *x_ptr,*xtemp_ptr; /* location to put the pointer to X and Xtemp respectively */ 7142c9cb062Sluzhanghpp PetscReal h = ts->time_step; 7152c9cb062Sluzhanghpp PetscInt s = tab->s,j; 7162c9cb062Sluzhanghpp PetscInt len_slow,len_fast; /* the number of slow components and fast components respectively */ 7172c9cb062Sluzhanghpp const PetscInt *is_slow,*is_fast; /* the index for slow components and fast components respectively */ 7182c9cb062Sluzhanghpp PetscErrorCode ierr; 7192c9cb062Sluzhanghpp 7202c9cb062Sluzhanghpp PetscFunctionBegin; 7212c9cb062Sluzhanghpp ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 7222c9cb062Sluzhanghpp ierr = VecDuplicate(ts->vec_sol,&Xtemp);CHKERRQ(ierr); 7232c9cb062Sluzhanghpp ierr = VecCopy(ts->vec_sol,Xtemp);CHKERRQ(ierr); 7242c9cb062Sluzhanghpp if (!rk->slow){ 7252c9cb062Sluzhanghpp for (j=0; j<s; j++) w[j] = h*tab->b[j]; 7262c9cb062Sluzhanghpp /* update the value of slow components, and discard the updating value of fast components */ 7272c9cb062Sluzhanghpp ierr = VecMAXPY(Xtemp,s,w,rk->YdotRHS);CHKERRQ(ierr); 7282c9cb062Sluzhanghpp ierr = VecGetArray(X,&x_ptr);CHKERRQ(ierr); 7292c9cb062Sluzhanghpp ierr = VecGetArray(Xtemp,&xtemp_ptr);CHKERRQ(ierr); 7302c9cb062Sluzhanghpp ierr = ISGetSize(ts->iss,&len_slow);CHKERRQ(ierr); 7312c9cb062Sluzhanghpp ierr = ISGetIndices(ts->iss,&is_slow);CHKERRQ(ierr); 7322c9cb062Sluzhanghpp ierr = ISRestoreIndices(ts->iss,&is_slow);CHKERRQ(ierr); 7332c9cb062Sluzhanghpp for (j=0; j<len_slow; j++) x_ptr[is_slow[j]] = xtemp_ptr[is_slow[j]]; 7342c9cb062Sluzhanghpp ierr = VecRestoreArray(X,&x_ptr);CHKERRQ(ierr); 7352c9cb062Sluzhanghpp ierr = VecRestoreArray(Xtemp,&xtemp_ptr);CHKERRQ(ierr); 7362c9cb062Sluzhanghpp } else { 7372c9cb062Sluzhanghpp /* update the value of fast components, and discard the updating value of slow components */ 7382c9cb062Sluzhanghpp for (j=0; j<s; j++) w[j] = h/rk->dtratio*tab->b[j]; 7392c9cb062Sluzhanghpp ierr = VecMAXPY(Xtemp,s,w,rk->YdotRHS);CHKERRQ(ierr); 7402c9cb062Sluzhanghpp ierr = VecGetArray(X,&x_ptr);CHKERRQ(ierr); 7412c9cb062Sluzhanghpp ierr = VecGetArray(Xtemp,&xtemp_ptr);CHKERRQ(ierr); 7422c9cb062Sluzhanghpp ierr = ISGetSize(ts->isf,&len_fast);CHKERRQ(ierr); 7432c9cb062Sluzhanghpp ierr = ISGetIndices(ts->isf,&is_fast);CHKERRQ(ierr); 7442c9cb062Sluzhanghpp ierr = ISRestoreIndices(ts->isf,&is_fast);CHKERRQ(ierr); 7452c9cb062Sluzhanghpp for (j=0; j<len_fast; j++) x_ptr[is_fast[j]] = xtemp_ptr[is_fast[j]]; 7462c9cb062Sluzhanghpp ierr = VecRestoreArray(X,&x_ptr);CHKERRQ(ierr); 7472c9cb062Sluzhanghpp ierr = VecRestoreArray(Xtemp,&xtemp_ptr);CHKERRQ(ierr); 7482c9cb062Sluzhanghpp } 7492c9cb062Sluzhanghpp /* free temporary vector space */ 7502c9cb062Sluzhanghpp ierr = VecDestroy(&Xtemp);CHKERRQ(ierr); 7512c9cb062Sluzhanghpp PetscFunctionReturn(0); 7522c9cb062Sluzhanghpp } 7532c9cb062Sluzhanghpp 7542c9cb062Sluzhanghpp static PetscErrorCode TSStep_RKMULTIRATE(TS ts) 7552c9cb062Sluzhanghpp { 7562c9cb062Sluzhanghpp TS_RK *rk = (TS_RK*)ts->data; 7572c9cb062Sluzhanghpp RKTableau tab = rk->tableau; 7582c9cb062Sluzhanghpp Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 7592c9cb062Sluzhanghpp Vec stage_fast,stage_slow; /* vectors store the stage value of fast and slow components respectively */ 7602c9cb062Sluzhanghpp Vec presol_fast,newslo_fast; /* vectors store the previous and new time solution of fast components respectively */ 7612c9cb062Sluzhanghpp Vec YdotRHS_prev,prev_sol; /* vectors store the previous YdotRHS and solution respectively */ 7622c9cb062Sluzhanghpp const PetscInt s = tab->s,*is_slow,*is_fast; /* is_slow: index of slow components; is_fast: index of fast components */ 7632c9cb062Sluzhanghpp const PetscReal *A = tab->A,*c = tab->c; 7642c9cb062Sluzhanghpp PetscScalar *w = rk->work; 7652c9cb062Sluzhanghpp PetscScalar *y_ptr,*faststage_ptr,*slowstage_ptr; /* location to put the pointer to Y, stage_fast and stage_slow respectively */ 7662c9cb062Sluzhanghpp PetscScalar *ydotrhsprev_ptr,*ydotrhs_ptr; /* location to put the pointer to YdotRHS_prev and YdotRHS respectively */ 7672c9cb062Sluzhanghpp PetscInt i,j,k,len_slow,len_fast; /* len_slow: the number of slow comonents; len_fast: the number of fast components */ 7682c9cb062Sluzhanghpp PetscReal next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step; 7692c9cb062Sluzhanghpp PetscErrorCode ierr; 7702c9cb062Sluzhanghpp 7712c9cb062Sluzhanghpp PetscFunctionBegin; 7722c9cb062Sluzhanghpp rk->status = TS_STEP_INCOMPLETE; 7732c9cb062Sluzhanghpp ierr = VecDuplicate(ts->vec_sol,&stage_fast);CHKERRQ(ierr); 7742c9cb062Sluzhanghpp ierr = VecDuplicate(ts->vec_sol,&stage_slow);CHKERRQ(ierr); 7752c9cb062Sluzhanghpp ierr = VecDuplicate(ts->vec_sol,&YdotRHS_prev);CHKERRQ(ierr); 7762c9cb062Sluzhanghpp ierr = VecDuplicate(ts->vec_sol,&prev_sol);CHKERRQ(ierr); 7772c9cb062Sluzhanghpp ierr = ISGetSize(ts->iss,&len_slow);CHKERRQ(ierr); 7782c9cb062Sluzhanghpp ierr = ISGetSize(ts->isf,&len_fast);CHKERRQ(ierr); 7792c9cb062Sluzhanghpp ierr = ISGetIndices(ts->iss,&is_slow);CHKERRQ(ierr); 7802c9cb062Sluzhanghpp ierr = ISGetIndices(ts->isf,&is_fast);CHKERRQ(ierr); 7812c9cb062Sluzhanghpp ierr = ISRestoreIndices(ts->isf,&is_fast);CHKERRQ(ierr); 7822c9cb062Sluzhanghpp ierr = ISRestoreIndices(ts->iss,&is_slow);CHKERRQ(ierr); 7832c9cb062Sluzhanghpp for (i=0; i<s; i++) { 7842c9cb062Sluzhanghpp rk->stage_time = t + h*c[i]; 7852c9cb062Sluzhanghpp ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr); 7862c9cb062Sluzhanghpp ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 7872c9cb062Sluzhanghpp for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 7882c9cb062Sluzhanghpp ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 7892c9cb062Sluzhanghpp ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 7902c9cb062Sluzhanghpp /* compute the stage RHS */ 7912c9cb062Sluzhanghpp ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 7922c9cb062Sluzhanghpp } 7932c9cb062Sluzhanghpp ierr = VecCopy(ts->vec_sol,prev_sol);CHKERRQ(ierr); 7942c9cb062Sluzhanghpp rk->slow = 0; 7952c9cb062Sluzhanghpp ierr = TSEvaluateStep_RKMULTIRATE(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 7962c9cb062Sluzhanghpp for (k=0; k<rk->dtratio; k++){ 7972c9cb062Sluzhanghpp for (i=0; i<s; i++) { 7982c9cb062Sluzhanghpp rk->stage_time = t + h/rk->dtratio*c[i]; 7992c9cb062Sluzhanghpp ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr); 8002c9cb062Sluzhanghpp ierr = VecCopy(prev_sol,Y[i]);CHKERRQ(ierr); 8012c9cb062Sluzhanghpp ierr = VecCopy(prev_sol,stage_slow);CHKERRQ(ierr); 8022c9cb062Sluzhanghpp ierr = VecCopy(prev_sol,stage_fast);CHKERRQ(ierr); 8032c9cb062Sluzhanghpp /* guarantee the stage RHS for slow components do not change while updating the stage RHS for fast components */ 8042c9cb062Sluzhanghpp ierr = VecCopy(YdotRHS[i],YdotRHS_prev);CHKERRQ(ierr); 8052c9cb062Sluzhanghpp /* update the fast components stage value */ 8062c9cb062Sluzhanghpp for (j=0; j<i; j++) w[j] = h/rk->dtratio*A[i*s+j]; 8072c9cb062Sluzhanghpp ierr = VecMAXPY(stage_fast,i,w,YdotRHS);CHKERRQ(ierr); 8082c9cb062Sluzhanghpp /* update the slow components stage value */ 8092c9cb062Sluzhanghpp ierr = VecCopy(prev_sol,Y[0]);CHKERRQ(ierr); 8102c9cb062Sluzhanghpp ierr = TSInterpolate_RK(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],stage_slow);CHKERRQ(ierr); 8112c9cb062Sluzhanghpp /* combine the update fast components stage value and slow components stage value to Y[i] */ 8122c9cb062Sluzhanghpp ierr = VecGetArray(Y[i],&y_ptr);CHKERRQ(ierr); 8132c9cb062Sluzhanghpp ierr = VecGetArray(stage_slow,&slowstage_ptr);CHKERRQ(ierr); 8142c9cb062Sluzhanghpp ierr = VecGetArray(stage_fast,&faststage_ptr);CHKERRQ(ierr); 8152c9cb062Sluzhanghpp for (j=0; j<len_slow; j++) y_ptr[is_slow[j]] = slowstage_ptr[is_slow[j]]; 8162c9cb062Sluzhanghpp for (j=0; j<len_fast; j++) y_ptr[is_fast[j]] = faststage_ptr[is_fast[j]]; 8172c9cb062Sluzhanghpp ierr = VecRestoreArray(Y[i],&y_ptr);CHKERRQ(ierr); 8182c9cb062Sluzhanghpp ierr = VecRestoreArray(stage_slow,&slowstage_ptr);CHKERRQ(ierr); 8192c9cb062Sluzhanghpp ierr = VecRestoreArray(stage_fast,&faststage_ptr);CHKERRQ(ierr); 8202c9cb062Sluzhanghpp ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 8212c9cb062Sluzhanghpp /* compute the stage RHS */ 8222c9cb062Sluzhanghpp ierr = TSComputeRHSFunction(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 8232c9cb062Sluzhanghpp /* re-assign the value of slow components in YdotRHS[i] to the calculated value before working on smaller time step-size */ 8242c9cb062Sluzhanghpp ierr = VecGetArray(YdotRHS[i],&ydotrhs_ptr);CHKERRQ(ierr); 8252c9cb062Sluzhanghpp ierr = VecGetArray(YdotRHS_prev,&ydotrhsprev_ptr);CHKERRQ(ierr); 8262c9cb062Sluzhanghpp for (j=0; j<len_slow; j++) ydotrhs_ptr[is_slow[j]] = ydotrhsprev_ptr[is_slow[j]]; 8272c9cb062Sluzhanghpp ierr = VecRestoreArray(YdotRHS[i],&ydotrhs_ptr);CHKERRQ(ierr); 8282c9cb062Sluzhanghpp ierr = VecRestoreArray(YdotRHS_prev,&ydotrhsprev_ptr);CHKERRQ(ierr); 8292c9cb062Sluzhanghpp } 8302c9cb062Sluzhanghpp rk->slow = 1; 8312c9cb062Sluzhanghpp ierr = TSEvaluateStep_RKMULTIRATE(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 8322c9cb062Sluzhanghpp /* update the intial value for fast components */ 8332c9cb062Sluzhanghpp ierr = VecGetSubVector(prev_sol,ts->isf,&presol_fast);CHKERRQ(ierr); 8342c9cb062Sluzhanghpp ierr = VecGetSubVector(ts->vec_sol,ts->isf,&newslo_fast);CHKERRQ(ierr); 8352c9cb062Sluzhanghpp ierr = VecCopy(newslo_fast,presol_fast);CHKERRQ(ierr); 8362c9cb062Sluzhanghpp ierr = VecRestoreSubVector(prev_sol,ts->isf,&presol_fast);CHKERRQ(ierr); 8372c9cb062Sluzhanghpp ierr = VecRestoreSubVector(ts->vec_sol,ts->isf,&newslo_fast);CHKERRQ(ierr); 8382c9cb062Sluzhanghpp } 8392c9cb062Sluzhanghpp ts->ptime += ts->time_step; 8402c9cb062Sluzhanghpp ts->time_step = next_time_step; 8412c9cb062Sluzhanghpp rk->slow = 0; 8422c9cb062Sluzhanghpp rk->status = TS_STEP_COMPLETE; 8432c9cb062Sluzhanghpp /* free memory of work vectors */ 8442c9cb062Sluzhanghpp ierr = VecDestroy(&stage_fast);CHKERRQ(ierr); 8452c9cb062Sluzhanghpp ierr = VecDestroy(&stage_slow);CHKERRQ(ierr); 8462c9cb062Sluzhanghpp ierr = VecDestroy(&YdotRHS_prev);CHKERRQ(ierr); 8472c9cb062Sluzhanghpp ierr = VecDestroy(&prev_sol);CHKERRQ(ierr); 8482c9cb062Sluzhanghpp PetscFunctionReturn(0); 8492c9cb062Sluzhanghpp } 8502c9cb062Sluzhanghpp 8512c9cb062Sluzhanghpp /* 8522c9cb062Sluzhanghpp This is for partitioned RHS multirate RK method 8532c9cb062Sluzhanghpp The step completion formula is 8542c9cb062Sluzhanghpp 8552c9cb062Sluzhanghpp x1 = x0 + h b^T YdotRHS 8562c9cb062Sluzhanghpp 8572c9cb062Sluzhanghpp */ 8582c9cb062Sluzhanghpp static PetscErrorCode TSEvaluateStep_RKPMULTIRATE(TS ts,PetscInt order,Vec X,PetscBool *done) 8592c9cb062Sluzhanghpp { 8602c9cb062Sluzhanghpp TS_RK *rk = (TS_RK*)ts->data; 8612c9cb062Sluzhanghpp RKTableau tab = rk->tableau; 8622c9cb062Sluzhanghpp Vec Xslow,Xfast; /* subvectors of X which store slow components and fast components respectively */ 8632c9cb062Sluzhanghpp PetscScalar *w = rk->work; 8642c9cb062Sluzhanghpp PetscReal h = ts->time_step; 8652c9cb062Sluzhanghpp PetscInt s = tab->s,j; 8662c9cb062Sluzhanghpp PetscErrorCode ierr; 8672c9cb062Sluzhanghpp 8682c9cb062Sluzhanghpp PetscFunctionBegin; 8692c9cb062Sluzhanghpp ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 8702c9cb062Sluzhanghpp if (!rk->slow){ 8712c9cb062Sluzhanghpp for (j=0; j<s; j++) w[j] = h*tab->b[j]; 8722c9cb062Sluzhanghpp ierr = VecGetSubVector(ts->vec_sol,ts->iss,&Xslow);CHKERRQ(ierr); 8732c9cb062Sluzhanghpp ierr = VecMAXPY(Xslow,s,w,rk->YdotRHS_slow);CHKERRQ(ierr); 8742c9cb062Sluzhanghpp ierr = VecRestoreSubVector(ts->vec_sol,ts->iss,&Xslow);CHKERRQ(ierr);; 8752c9cb062Sluzhanghpp } else {; 8762c9cb062Sluzhanghpp for (j=0; j<s; j++) w[j] = h/rk->dtratio*tab->b[j]; 8772c9cb062Sluzhanghpp ierr = VecGetSubVector(X,ts->isf,&Xfast);CHKERRQ(ierr); 8782c9cb062Sluzhanghpp ierr = VecMAXPY(Xfast,s,w,rk->YdotRHS_fast);CHKERRQ(ierr); 8792c9cb062Sluzhanghpp ierr = VecRestoreSubVector(X,ts->isf,&Xfast);CHKERRQ(ierr); 8802c9cb062Sluzhanghpp } 8812c9cb062Sluzhanghpp PetscFunctionReturn(0); 8822c9cb062Sluzhanghpp } 8832c9cb062Sluzhanghpp 8842c9cb062Sluzhanghpp static PetscErrorCode TSStep_RKPMULTIRATE(TS ts) 8852c9cb062Sluzhanghpp { 8862c9cb062Sluzhanghpp TS_RK *rk = (TS_RK*)ts->data; 8872c9cb062Sluzhanghpp RKTableau tab = rk->tableau; 8882c9cb062Sluzhanghpp Vec *Y = rk->Y; 8892c9cb062Sluzhanghpp Vec *YdotRHS_fast = rk->YdotRHS_fast,*YdotRHS_slow = rk->YdotRHS_slow; 8902c9cb062Sluzhanghpp Vec Yslow,Yfast; /* subvectors store the stges of slow components and fast components respectively */ 8912c9cb062Sluzhanghpp Vec Yfast_prev,Yfast_new,prev_sol; /* subvectors store the previous and new solution of fast components and vector store the previous solution */ 8922c9cb062Sluzhanghpp const PetscInt s = tab->s; 8932c9cb062Sluzhanghpp const PetscReal *A = tab->A,*c = tab->c; 8942c9cb062Sluzhanghpp PetscScalar *w = rk->work; 8952c9cb062Sluzhanghpp PetscInt i,j,k; 8962c9cb062Sluzhanghpp PetscReal next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step; 8972c9cb062Sluzhanghpp PetscErrorCode ierr; 8982c9cb062Sluzhanghpp 8992c9cb062Sluzhanghpp PetscFunctionBegin; 9002c9cb062Sluzhanghpp rk->status = TS_STEP_INCOMPLETE; 9012c9cb062Sluzhanghpp for (i=0; i<s; i++) { 9022c9cb062Sluzhanghpp rk->stage_time = t + h*c[i]; 9032c9cb062Sluzhanghpp ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr); 9042c9cb062Sluzhanghpp ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 9052c9cb062Sluzhanghpp /*ierr = VecCopy(ts->vec_sol,Yc[i]);CHKERRQ(ierr);*/ 9062c9cb062Sluzhanghpp ierr = VecGetSubVector(Y[i],ts->iss,&Yslow);CHKERRQ(ierr); 9072c9cb062Sluzhanghpp ierr = VecGetSubVector(Y[i],ts->isf,&Yfast);CHKERRQ(ierr); 9082c9cb062Sluzhanghpp for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 9092c9cb062Sluzhanghpp ierr = VecMAXPY(Yslow,i,w,YdotRHS_slow);CHKERRQ(ierr); 9102c9cb062Sluzhanghpp ierr = VecMAXPY(Yfast,i,w,YdotRHS_fast);CHKERRQ(ierr); 9112c9cb062Sluzhanghpp ierr = VecRestoreSubVector(Y[i],ts->iss,&Yslow);CHKERRQ(ierr); 9122c9cb062Sluzhanghpp ierr = VecRestoreSubVector(Y[i],ts->isf,&Yfast);CHKERRQ(ierr); 9132c9cb062Sluzhanghpp ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 9142c9cb062Sluzhanghpp /* compute the stage RHS for both slow and fast components */ 9152c9cb062Sluzhanghpp ierr = TSComputeRHSFunctionslow(ts,t+h*c[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr); 9162c9cb062Sluzhanghpp ierr = TSComputeRHSFunctionfast(ts,t+h*c[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr); 9172c9cb062Sluzhanghpp } 9182c9cb062Sluzhanghpp /* store the value of slow components at previous time */ 9192c9cb062Sluzhanghpp ierr = VecDuplicate(ts->vec_sol,&prev_sol);CHKERRQ(ierr); 9202c9cb062Sluzhanghpp ierr = VecCopy(ts->vec_sol,prev_sol);CHKERRQ(ierr); 9212c9cb062Sluzhanghpp rk->slow = 0; 9222c9cb062Sluzhanghpp ierr = TSEvaluateStep_RKPMULTIRATE(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 9232c9cb062Sluzhanghpp for (k=0; k<rk->dtratio; k++){ 9242c9cb062Sluzhanghpp for (i=0; i<s; i++) { 9252c9cb062Sluzhanghpp rk->stage_time = t + h/rk->dtratio*c[i]; 9262c9cb062Sluzhanghpp ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr); 9272c9cb062Sluzhanghpp ierr = VecCopy(prev_sol,Y[i]);CHKERRQ(ierr); 9282c9cb062Sluzhanghpp /* stage value for fast components */ 9292c9cb062Sluzhanghpp for (j=0; j<i; j++) w[j] = h/rk->dtratio*A[i*s+j]; 9302c9cb062Sluzhanghpp ierr = VecGetSubVector(Y[i],ts->isf,&Yfast);CHKERRQ(ierr); 9312c9cb062Sluzhanghpp ierr = VecMAXPY(Yfast,i,w,YdotRHS_fast);CHKERRQ(ierr); 9322c9cb062Sluzhanghpp ierr = VecRestoreSubVector(Y[i],ts->isf,&Yfast);CHKERRQ(ierr); 9332c9cb062Sluzhanghpp /* stage value for slow components */ 9342c9cb062Sluzhanghpp ierr = VecCopy(prev_sol,Y[0]);CHKERRQ(ierr); 9352c9cb062Sluzhanghpp ierr = VecGetSubVector(Y[i],ts->iss,&Yslow);CHKERRQ(ierr); 9362c9cb062Sluzhanghpp ierr = TSInterpolate_PARTITIONEDMRK(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Yslow);CHKERRQ(ierr); 9372c9cb062Sluzhanghpp ierr = VecRestoreSubVector(Y[i],ts->iss,&Yslow);CHKERRQ(ierr); 9382c9cb062Sluzhanghpp ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 9392c9cb062Sluzhanghpp /* compute the stage RHS for fast components */ 9402c9cb062Sluzhanghpp ierr = TSComputeRHSFunctionfast(ts,t+k*h/rk->dtratio+h/rk->dtratio*c[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr); 9412c9cb062Sluzhanghpp } 9422c9cb062Sluzhanghpp /* update the value of fast components whenusing fast time step */ 9432c9cb062Sluzhanghpp rk->slow = 1; 9442c9cb062Sluzhanghpp ierr = TSEvaluateStep_RKPMULTIRATE(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 9452c9cb062Sluzhanghpp ierr = VecGetSubVector(prev_sol,ts->isf,&Yfast_prev);CHKERRQ(ierr); 9462c9cb062Sluzhanghpp ierr = VecGetSubVector(ts->vec_sol,ts->isf,&Yfast_new);CHKERRQ(ierr); 9472c9cb062Sluzhanghpp ierr = VecCopy(Yfast_new,Yfast_prev);CHKERRQ(ierr); 9482c9cb062Sluzhanghpp ierr = VecRestoreSubVector(prev_sol,ts->isf,&Yfast_prev);CHKERRQ(ierr); 9492c9cb062Sluzhanghpp ierr = VecRestoreSubVector(ts->vec_sol,ts->isf,&Yfast_new);CHKERRQ(ierr); 9502c9cb062Sluzhanghpp } 9512c9cb062Sluzhanghpp ts->ptime += ts->time_step; 9522c9cb062Sluzhanghpp ts->time_step = next_time_step; 9532c9cb062Sluzhanghpp rk->slow = 0; 9542c9cb062Sluzhanghpp rk->status = TS_STEP_COMPLETE; 9552c9cb062Sluzhanghpp ierr = VecDestroy(&prev_sol);CHKERRQ(ierr); 9562c9cb062Sluzhanghpp PetscFunctionReturn(0); 9572c9cb062Sluzhanghpp } 9582c9cb062Sluzhanghpp 959f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts) 960f6a906c0SBarry Smith { 961f6a906c0SBarry Smith TS_RK *rk = (TS_RK*)ts->data; 962be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 963be5899b3SLisandro Dalcin PetscInt s = tab->s; 964f6a906c0SBarry Smith PetscErrorCode ierr; 965f6a906c0SBarry Smith 966f6a906c0SBarry Smith PetscFunctionBegin; 967f6a906c0SBarry Smith if (ts->adjointsetupcalled++) PetscFunctionReturn(0); 968abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr); 969f6a906c0SBarry Smith if(ts->vecs_sensip) { 970abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr); 971f6a906c0SBarry Smith } 972f6a906c0SBarry Smith PetscFunctionReturn(0); 973f6a906c0SBarry Smith } 974f6a906c0SBarry Smith 97542f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts) 976d2daff3dSHong Zhang { 977c235aa8dSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 978c235aa8dSHong Zhang RKTableau tab = rk->tableau; 979c235aa8dSHong Zhang const PetscInt s = tab->s; 980c235aa8dSHong Zhang const PetscReal *A = tab->A,*b = tab->b,*c = tab->c; 981c235aa8dSHong Zhang PetscScalar *w = rk->work; 9823389c7d9SStefano Zampini Vec *Y = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu; 983b18ea86cSHong Zhang PetscInt i,j,nadj; 9843d3eaba7SBarry Smith PetscReal t = ts->ptime; 985d2daff3dSHong Zhang PetscErrorCode ierr; 986cef76868SBarry Smith PetscReal h = ts->time_step; 987c235aa8dSHong Zhang 988d2daff3dSHong Zhang PetscFunctionBegin; 989c235aa8dSHong Zhang rk->status = TS_STEP_INCOMPLETE; 990c235aa8dSHong Zhang for (i=s-1; i>=0; i--) { 9913389c7d9SStefano Zampini Mat J; 9923389c7d9SStefano Zampini PetscReal stage_time = t + h*(1.0-c[i]); 9933389c7d9SStefano Zampini PetscBool zero = PETSC_FALSE; 9943389c7d9SStefano Zampini 9953389c7d9SStefano Zampini ierr = TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr); 9963389c7d9SStefano Zampini ierr = TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);CHKERRQ(ierr); 997493ed95fSHong Zhang if (ts->vec_costintegral) { 9983389c7d9SStefano Zampini ierr = TSAdjointComputeDRDYFunction(ts,stage_time,Y[i],ts->vecs_drdy);CHKERRQ(ierr); 999493ed95fSHong Zhang } 10003389c7d9SStefano Zampini /* Stage values of mu */ 10013389c7d9SStefano Zampini if (ts->vecs_sensip) { 10023389c7d9SStefano Zampini ierr = TSAdjointComputeRHSJacobian(ts,stage_time,Y[i],ts->Jacp);CHKERRQ(ierr); 10033389c7d9SStefano Zampini if (ts->vec_costintegral) { 10043389c7d9SStefano Zampini ierr = TSAdjointComputeDRDPFunction(ts,stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr); 10053389c7d9SStefano Zampini } 10063389c7d9SStefano Zampini } 10073389c7d9SStefano Zampini 10083389c7d9SStefano Zampini if (b[i] == 0 && i == s-1) zero = PETSC_TRUE; 1009abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 10103389c7d9SStefano Zampini DM dm; 10113389c7d9SStefano Zampini Vec VecSensiTemp; 10123389c7d9SStefano Zampini 10133389c7d9SStefano Zampini ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 10143389c7d9SStefano Zampini ierr = DMGetGlobalVector(dm,&VecSensiTemp);CHKERRQ(ierr); 10153389c7d9SStefano Zampini /* Stage values of lambda */ 10163389c7d9SStefano Zampini if (!zero) { 10173389c7d9SStefano Zampini ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp);CHKERRQ(ierr); 10183389c7d9SStefano Zampini ierr = VecScale(VecSensiTemp,-h*b[i]);CHKERRQ(ierr); 10193389c7d9SStefano Zampini for (j=i+1; j<s; j++) w[j-i-1] = -h*A[j*s+i]; 10203389c7d9SStefano Zampini ierr = VecMAXPY(VecSensiTemp,s-i-1,w,&VecDeltaLam[nadj*s+i+1]);CHKERRQ(ierr); 10213389c7d9SStefano Zampini ierr = MatMultTranspose(J,VecSensiTemp,VecDeltaLam[nadj*s+i]);CHKERRQ(ierr); 10223389c7d9SStefano Zampini } else { 10233389c7d9SStefano Zampini ierr = VecSet(VecDeltaLam[nadj*s+i],0);CHKERRQ(ierr); 10243389c7d9SStefano Zampini } 1025493ed95fSHong Zhang if (ts->vec_costintegral) { 1026493ed95fSHong Zhang ierr = VecAXPY(VecDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdy[nadj]);CHKERRQ(ierr); 1027493ed95fSHong Zhang } 10286497ce14SHong Zhang 1029ad8e2604SHong Zhang /* Stage values of mu */ 10306497ce14SHong Zhang if (ts->vecs_sensip) { 10313389c7d9SStefano Zampini if (!zero) { 10323389c7d9SStefano Zampini ierr = MatMultTranspose(ts->Jacp,VecSensiTemp,VecDeltaMu[nadj*s+i]);CHKERRQ(ierr); 10333389c7d9SStefano Zampini } else { 10343389c7d9SStefano Zampini ierr = VecSet(VecDeltaMu[nadj*s+i],0);CHKERRQ(ierr); 1035493ed95fSHong Zhang } 1036493ed95fSHong Zhang if (ts->vec_costintegral) { 1037493ed95fSHong Zhang ierr = VecAXPY(VecDeltaMu[nadj*s+i],-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr); 1038493ed95fSHong Zhang } 1039ad8e2604SHong Zhang } 10403389c7d9SStefano Zampini ierr = DMRestoreGlobalVector(dm,&VecSensiTemp);CHKERRQ(ierr); 1041c235aa8dSHong Zhang } 10426497ce14SHong Zhang } 1043c235aa8dSHong Zhang 1044c235aa8dSHong Zhang for (j=0; j<s; j++) w[j] = 1.0; 1045abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 1046b18ea86cSHong Zhang ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr); 10476497ce14SHong Zhang if (ts->vecs_sensip) { 1048ad8e2604SHong Zhang ierr = VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);CHKERRQ(ierr); 1049b18ea86cSHong Zhang } 10506497ce14SHong Zhang } 1051c235aa8dSHong Zhang rk->status = TS_STEP_COMPLETE; 1052d2daff3dSHong Zhang PetscFunctionReturn(0); 1053d2daff3dSHong Zhang } 1054d2daff3dSHong Zhang 1055be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts) 1056be5899b3SLisandro Dalcin { 1057be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 1058be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 1059be5899b3SLisandro Dalcin PetscErrorCode ierr; 1060be5899b3SLisandro Dalcin 1061be5899b3SLisandro Dalcin PetscFunctionBegin; 1062be5899b3SLisandro Dalcin if (!tab) PetscFunctionReturn(0); 1063be5899b3SLisandro Dalcin ierr = PetscFree(rk->work);CHKERRQ(ierr); 1064be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr); 1065be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr); 10662c9cb062Sluzhanghpp ierr = VecDestroyVecs(tab->s,&rk->YdotRHS_fast);CHKERRQ(ierr); 10672c9cb062Sluzhanghpp ierr = VecDestroyVecs(tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr); 1068be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr); 1069be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr); 1070be5899b3SLisandro Dalcin PetscFunctionReturn(0); 1071be5899b3SLisandro Dalcin } 1072be5899b3SLisandro Dalcin 1073277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts) 1074e4dd521cSBarry Smith { 10755f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 10766849ba73SBarry Smith PetscErrorCode ierr; 1077e4dd521cSBarry Smith 1078e4dd521cSBarry Smith PetscFunctionBegin; 1079be5899b3SLisandro Dalcin ierr = TSRKTableauReset(ts);CHKERRQ(ierr); 10800f7a1166SHong Zhang ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr); 1081277b19d0SLisandro Dalcin PetscFunctionReturn(0); 1082e4dd521cSBarry Smith } 1083277b19d0SLisandro Dalcin 1084f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx) 1085f68a32c8SEmil Constantinescu { 1086f68a32c8SEmil Constantinescu PetscFunctionBegin; 1087f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1088f68a32c8SEmil Constantinescu } 1089f68a32c8SEmil Constantinescu 1090f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 1091f68a32c8SEmil Constantinescu { 1092f68a32c8SEmil Constantinescu PetscFunctionBegin; 1093f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1094f68a32c8SEmil Constantinescu } 1095f68a32c8SEmil Constantinescu 1096f68a32c8SEmil Constantinescu 1097f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx) 1098f68a32c8SEmil Constantinescu { 1099f68a32c8SEmil Constantinescu PetscFunctionBegin; 1100f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1101f68a32c8SEmil Constantinescu } 1102f68a32c8SEmil Constantinescu 1103f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 1104f68a32c8SEmil Constantinescu { 1105f68a32c8SEmil Constantinescu 1106f68a32c8SEmil Constantinescu PetscFunctionBegin; 1107f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1108f68a32c8SEmil Constantinescu } 1109c235aa8dSHong Zhang /* 1110d2daff3dSHong Zhang static PetscErrorCode RKSetAdjCoe(RKTableau tab) 1111d2daff3dSHong Zhang { 1112d2daff3dSHong Zhang PetscReal *A,*b; 1113d2daff3dSHong Zhang PetscInt s,i,j; 1114d2daff3dSHong Zhang PetscErrorCode ierr; 1115d2daff3dSHong Zhang 1116d2daff3dSHong Zhang PetscFunctionBegin; 1117d2daff3dSHong Zhang s = tab->s; 1118d2daff3dSHong Zhang ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr); 1119d2daff3dSHong Zhang 1120d2daff3dSHong Zhang for (i=0; i<s; i++) 1121d2daff3dSHong Zhang for (j=0; j<s; j++) { 1122d2daff3dSHong 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]; 1123d2daff3dSHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr); 1124d2daff3dSHong Zhang } 1125d2daff3dSHong Zhang 1126d2daff3dSHong Zhang for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i]; 1127d2daff3dSHong Zhang 1128d2daff3dSHong Zhang ierr = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 1129d2daff3dSHong Zhang ierr = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr); 1130d2daff3dSHong Zhang ierr = PetscFree2(A,b);CHKERRQ(ierr); 1131d2daff3dSHong Zhang PetscFunctionReturn(0); 1132d2daff3dSHong Zhang } 1133c235aa8dSHong Zhang */ 1134be5899b3SLisandro Dalcin 1135be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts) 1136be5899b3SLisandro Dalcin { 1137be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 1138be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 11392c9cb062Sluzhanghpp Vec YdotRHS_fast,YdotRHS_slow; 1140be5899b3SLisandro Dalcin PetscErrorCode ierr; 1141be5899b3SLisandro Dalcin 1142be5899b3SLisandro Dalcin PetscFunctionBegin; 1143be5899b3SLisandro Dalcin ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr); 1144be5899b3SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr); 1145be5899b3SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr); 11462c9cb062Sluzhanghpp ierr = VecGetSubVector(ts->vec_sol,ts->isf,&YdotRHS_fast);CHKERRQ(ierr); 11472c9cb062Sluzhanghpp ierr = VecDuplicateVecs(YdotRHS_fast,tab->s,&rk->YdotRHS_fast);CHKERRQ(ierr); 11482c9cb062Sluzhanghpp ierr = VecRestoreSubVector(ts->vec_sol,ts->isf,&YdotRHS_fast);CHKERRQ(ierr); 11492c9cb062Sluzhanghpp ierr = VecGetSubVector(ts->vec_sol,ts->iss,&YdotRHS_slow);CHKERRQ(ierr); 11502c9cb062Sluzhanghpp ierr = VecDuplicateVecs(YdotRHS_slow,tab->s,&rk->YdotRHS_slow);CHKERRQ(ierr); 11512c9cb062Sluzhanghpp ierr = VecRestoreSubVector(ts->vec_sol,ts->iss,&YdotRHS_slow);CHKERRQ(ierr); 11522c9cb062Sluzhanghpp 1153be5899b3SLisandro Dalcin PetscFunctionReturn(0); 1154be5899b3SLisandro Dalcin } 1155be5899b3SLisandro Dalcin 1156be5899b3SLisandro Dalcin 1157f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts) 1158f68a32c8SEmil Constantinescu { 1159f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 1160f68a32c8SEmil Constantinescu PetscErrorCode ierr; 1161f68a32c8SEmil Constantinescu DM dm; 1162f68a32c8SEmil Constantinescu 1163f68a32c8SEmil Constantinescu PetscFunctionBegin; 11643dd83b38SBarry Smith ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr); 1165be5899b3SLisandro Dalcin ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr); 11660f7a1166SHong Zhang if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 11670f7a1166SHong Zhang ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr); 11680f7a1166SHong Zhang } 1169f68a32c8SEmil Constantinescu ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 1170f68a32c8SEmil Constantinescu ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 1171f68a32c8SEmil Constantinescu ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 1172e4dd521cSBarry Smith PetscFunctionReturn(0); 1173e4dd521cSBarry Smith } 1174d2daff3dSHong Zhang 11752c9cb062Sluzhanghpp /* 11762c9cb062Sluzhanghpp construct a database to choose single-step rk method, or combined rhs mutirate rk method or partitioned rhs rk method 11772c9cb062Sluzhanghpp */ 11782c9cb062Sluzhanghpp const char *const TSRKMultirateTypes[] = {"NONE","COMBINED","PARTITIONED","TSRKMultirateType","RKM_",0}; 1179e4dd521cSBarry Smith 11804416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts) 1181e4dd521cSBarry Smith { 1182be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 1183dfbe8321SBarry Smith PetscErrorCode ierr; 1184262ff7bbSBarry Smith 1185e4dd521cSBarry Smith PetscFunctionBegin; 1186e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr); 1187f68a32c8SEmil Constantinescu { 1188f68a32c8SEmil Constantinescu RKTableauLink link; 1189f68a32c8SEmil Constantinescu PetscInt count,choice; 1190f68a32c8SEmil Constantinescu PetscBool flg; 1191f68a32c8SEmil Constantinescu const char **namelist; 11922c9cb062Sluzhanghpp PetscInt rkmtype = 0; 11932c9cb062Sluzhanghpp 1194f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) ; 1195f489ac74SBarry Smith ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr); 1196f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 11972c9cb062Sluzhanghpp ierr = PetscOptionsEList("-ts_rk_multirate_type","Use Multirate or partioned RHS Multirate or single rate RK method","TSRKSetMultirateType",TSRKMultirateTypes,3,TSRKMultirateTypes[0],&rkmtype,&flg);CHKERRQ(ierr); 11982c9cb062Sluzhanghpp if (flg) {ierr = TSRKSetMultirateType(ts,rkmtype);CHKERRQ(ierr);} 1199be5899b3SLisandro Dalcin ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr); 1200be5899b3SLisandro Dalcin if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);} 1201f68a32c8SEmil Constantinescu ierr = PetscFree(namelist);CHKERRQ(ierr); 1202f68a32c8SEmil Constantinescu } 1203262ff7bbSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 12042c9cb062Sluzhanghpp ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Multirate methods options","");CHKERRQ(ierr); 12052c9cb062Sluzhanghpp ierr = PetscOptionsInt("-ts_rk_dtratio","time step ratio between slow and fast","",rk->dtratio,&rk->dtratio,NULL);CHKERRQ(ierr); 12062c9cb062Sluzhanghpp ierr = PetscOptionsEnd();CHKERRQ(ierr); 1207e4dd521cSBarry Smith PetscFunctionReturn(0); 1208e4dd521cSBarry Smith } 1209e4dd521cSBarry Smith 12105f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer) 1211e4dd521cSBarry Smith { 12125f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 12138370ee3bSLisandro Dalcin PetscBool iascii; 1214dfbe8321SBarry Smith PetscErrorCode ierr; 12152cdf8120Sasbjorn 1216e4dd521cSBarry Smith PetscFunctionBegin; 1217251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 12188370ee3bSLisandro Dalcin if (iascii) { 12199c334d8fSLisandro Dalcin RKTableau tab = rk->tableau; 1220f68a32c8SEmil Constantinescu TSRKType rktype; 1221f68a32c8SEmil Constantinescu char buf[512]; 1222f68a32c8SEmil Constantinescu ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr); 1223efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RK type %s\n",rktype);CHKERRQ(ierr); 12240c37eb8eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," Order: %D\n",tab->order);CHKERRQ(ierr); 12250c37eb8eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr); 1226f68a32c8SEmil Constantinescu ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 1227f68a32c8SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," Abscissa c = %s\n",buf);CHKERRQ(ierr); 12288370ee3bSLisandro Dalcin } 1229f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1230f68a32c8SEmil Constantinescu } 1231f68a32c8SEmil Constantinescu 1232f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer) 1233f68a32c8SEmil Constantinescu { 1234f68a32c8SEmil Constantinescu PetscErrorCode ierr; 12359c334d8fSLisandro Dalcin TSAdapt adapt; 1236f68a32c8SEmil Constantinescu 1237f68a32c8SEmil Constantinescu PetscFunctionBegin; 12389c334d8fSLisandro Dalcin ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 12399c334d8fSLisandro Dalcin ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr); 1240f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1241f68a32c8SEmil Constantinescu } 1242f68a32c8SEmil Constantinescu 1243f68a32c8SEmil Constantinescu /*@C 1244f68a32c8SEmil Constantinescu TSRKSetType - Set the type of RK scheme 1245f68a32c8SEmil Constantinescu 1246f68a32c8SEmil Constantinescu Logically collective 1247f68a32c8SEmil Constantinescu 1248f68a32c8SEmil Constantinescu Input Parameter: 1249f68a32c8SEmil Constantinescu + ts - timestepping context 1250f68a32c8SEmil Constantinescu - rktype - type of RK-scheme 1251f68a32c8SEmil Constantinescu 1252287c2655SBarry Smith Options Database: 12539bd3e852SBarry Smith . -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs> 1254287c2655SBarry Smith 1255f68a32c8SEmil Constantinescu Level: intermediate 1256f68a32c8SEmil Constantinescu 1257287c2655SBarry Smith .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS 1258f68a32c8SEmil Constantinescu @*/ 1259f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype) 1260f68a32c8SEmil Constantinescu { 1261f68a32c8SEmil Constantinescu PetscErrorCode ierr; 1262f68a32c8SEmil Constantinescu 1263f68a32c8SEmil Constantinescu PetscFunctionBegin; 1264f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1265b92453a8SLisandro Dalcin PetscValidCharPointer(rktype,2); 1266f68a32c8SEmil Constantinescu ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr); 1267f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1268f68a32c8SEmil Constantinescu } 1269f68a32c8SEmil Constantinescu 1270f68a32c8SEmil Constantinescu /*@C 12712c9cb062Sluzhanghpp TSRKSetMultirateType - Set the type of RK Multirate scheme 12722c9cb062Sluzhanghpp 12732c9cb062Sluzhanghpp Logically collective 12742c9cb062Sluzhanghpp 12752c9cb062Sluzhanghpp Input Parameter: 12762c9cb062Sluzhanghpp + ts - timestepping context 12772c9cb062Sluzhanghpp - rkmtype - type of RKM-scheme 12782c9cb062Sluzhanghpp 12792c9cb062Sluzhanghpp Options Database: 12802c9cb062Sluzhanghpp . -ts_rk_multiarte_type - <none,combined,partitioned> 12812c9cb062Sluzhanghpp 12822c9cb062Sluzhanghpp Level: intermediate 12832c9cb062Sluzhanghpp @*/ 12842c9cb062Sluzhanghpp PetscErrorCode TSRKSetMultirateType(TS ts, TSRKMultirateType rkmtype) 12852c9cb062Sluzhanghpp { 12862c9cb062Sluzhanghpp TS_RK *rk = (TS_RK*)ts->data; 12872c9cb062Sluzhanghpp PetscErrorCode ierr; 12882c9cb062Sluzhanghpp 12892c9cb062Sluzhanghpp PetscFunctionBegin; 12902c9cb062Sluzhanghpp PetscValidHeaderSpecific(ts,TS_CLASSID,1); 12912c9cb062Sluzhanghpp switch(rkmtype){ 12922c9cb062Sluzhanghpp case RKM_NONE: 12932c9cb062Sluzhanghpp break; 12942c9cb062Sluzhanghpp case RKM_COMBINED: 12952c9cb062Sluzhanghpp ts->ops->step = TSStep_RKMULTIRATE; 12962c9cb062Sluzhanghpp ts->ops->evaluatestep = TSEvaluateStep_RKMULTIRATE; 12972c9cb062Sluzhanghpp rk->dtratio = 2; 12982c9cb062Sluzhanghpp ierr = TSRKSetType(ts,TSMRKDefault);CHKERRQ(ierr); 12992c9cb062Sluzhanghpp break; 13002c9cb062Sluzhanghpp case RKM_PARTITIONED: 13012c9cb062Sluzhanghpp ts->ops->step = TSStep_RKPMULTIRATE; 13022c9cb062Sluzhanghpp ts->ops->evaluatestep = TSEvaluateStep_RKPMULTIRATE; 13032c9cb062Sluzhanghpp ts->ops->interpolate = TSInterpolate_PARTITIONEDMRK; 13042c9cb062Sluzhanghpp rk->dtratio = 2; 13052c9cb062Sluzhanghpp ierr = TSRKSetType(ts,TSMRKDefault);CHKERRQ(ierr); 13062c9cb062Sluzhanghpp break; 13072c9cb062Sluzhanghpp default : 13082c9cb062Sluzhanghpp SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown type '%s'",rkmtype); 13092c9cb062Sluzhanghpp } 13102c9cb062Sluzhanghpp PetscFunctionReturn(0); 13112c9cb062Sluzhanghpp } 13122c9cb062Sluzhanghpp 13132c9cb062Sluzhanghpp /*@C 1314f68a32c8SEmil Constantinescu TSRKGetType - Get the type of RK scheme 1315f68a32c8SEmil Constantinescu 1316f68a32c8SEmil Constantinescu Logically collective 1317f68a32c8SEmil Constantinescu 1318f68a32c8SEmil Constantinescu Input Parameter: 1319f68a32c8SEmil Constantinescu . ts - timestepping context 1320f68a32c8SEmil Constantinescu 1321f68a32c8SEmil Constantinescu Output Parameter: 1322f68a32c8SEmil Constantinescu . rktype - type of RK-scheme 1323f68a32c8SEmil Constantinescu 1324f68a32c8SEmil Constantinescu Level: intermediate 1325f68a32c8SEmil Constantinescu 1326f68a32c8SEmil Constantinescu .seealso: TSRKGetType() 1327f68a32c8SEmil Constantinescu @*/ 1328f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype) 1329f68a32c8SEmil Constantinescu { 1330f68a32c8SEmil Constantinescu PetscErrorCode ierr; 1331f68a32c8SEmil Constantinescu 1332f68a32c8SEmil Constantinescu PetscFunctionBegin; 1333f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1334f68a32c8SEmil Constantinescu ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr); 1335f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1336f68a32c8SEmil Constantinescu } 1337f68a32c8SEmil Constantinescu 1338560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype) 1339f68a32c8SEmil Constantinescu { 1340f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 1341f68a32c8SEmil Constantinescu 1342f68a32c8SEmil Constantinescu PetscFunctionBegin; 1343f68a32c8SEmil Constantinescu *rktype = rk->tableau->name; 1344f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1345f68a32c8SEmil Constantinescu } 13462c9cb062Sluzhanghpp 1347560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype) 1348f68a32c8SEmil Constantinescu { 1349f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 1350f68a32c8SEmil Constantinescu PetscErrorCode ierr; 1351f68a32c8SEmil Constantinescu PetscBool match; 1352f68a32c8SEmil Constantinescu RKTableauLink link; 1353f68a32c8SEmil Constantinescu 1354f68a32c8SEmil Constantinescu PetscFunctionBegin; 1355f68a32c8SEmil Constantinescu if (rk->tableau) { 1356f68a32c8SEmil Constantinescu ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr); 1357f68a32c8SEmil Constantinescu if (match) PetscFunctionReturn(0); 1358f68a32c8SEmil Constantinescu } 1359f68a32c8SEmil Constantinescu for (link = RKTableauList; link; link=link->next) { 1360f68a32c8SEmil Constantinescu ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr); 1361f68a32c8SEmil Constantinescu if (match) { 1362be5899b3SLisandro Dalcin if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);} 1363f68a32c8SEmil Constantinescu rk->tableau = &link->tab; 1364be5899b3SLisandro Dalcin if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);} 13652ffb9264SLisandro Dalcin ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 1366f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1367f68a32c8SEmil Constantinescu } 1368f68a32c8SEmil Constantinescu } 1369f68a32c8SEmil Constantinescu SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype); 1370e4dd521cSBarry Smith PetscFunctionReturn(0); 1371e4dd521cSBarry Smith } 1372e4dd521cSBarry Smith 1373ff22ae23SHong Zhang static PetscErrorCode TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y) 1374ff22ae23SHong Zhang { 1375ff22ae23SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 1376ff22ae23SHong Zhang 1377ff22ae23SHong Zhang PetscFunctionBegin; 13780429704eSStefano Zampini if (ns) *ns = rk->tableau->s; 1379d2daff3dSHong Zhang if (Y) *Y = rk->Y; 1380ff22ae23SHong Zhang PetscFunctionReturn(0); 1381ff22ae23SHong Zhang } 1382ff22ae23SHong Zhang 1383b3a6b972SJed Brown static PetscErrorCode TSDestroy_RK(TS ts) 1384b3a6b972SJed Brown { 1385b3a6b972SJed Brown PetscErrorCode ierr; 1386b3a6b972SJed Brown 1387b3a6b972SJed Brown PetscFunctionBegin; 1388b3a6b972SJed Brown ierr = TSReset_RK(ts);CHKERRQ(ierr); 1389b3a6b972SJed Brown if (ts->dm) { 1390b3a6b972SJed Brown ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 1391b3a6b972SJed Brown ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 1392b3a6b972SJed Brown } 1393b3a6b972SJed Brown ierr = PetscFree(ts->data);CHKERRQ(ierr); 1394b3a6b972SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr); 1395b3a6b972SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr); 13962c9cb062Sluzhanghpp ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirateType_C",NULL);CHKERRQ(ierr); 1397b3a6b972SJed Brown PetscFunctionReturn(0); 1398b3a6b972SJed Brown } 1399ff22ae23SHong Zhang 1400ebe3b25bSBarry Smith /*MC 1401f68a32c8SEmil Constantinescu TSRK - ODE and DAE solver using Runge-Kutta schemes 1402ebe3b25bSBarry Smith 1403f68a32c8SEmil Constantinescu The user should provide the right hand side of the equation 1404f68a32c8SEmil Constantinescu using TSSetRHSFunction(). 1405ebe3b25bSBarry Smith 1406f68a32c8SEmil Constantinescu Notes: 140798164e13SLisandro Dalcin The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type 1408ebe3b25bSBarry Smith 1409d41222bbSBarry Smith Level: beginner 1410d41222bbSBarry Smith 1411f68a32c8SEmil Constantinescu .seealso: TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3, 14122c9cb062Sluzhanghpp TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister(), TSRKSetMultirateType() 1413ebe3b25bSBarry Smith 1414ebe3b25bSBarry Smith M*/ 14158cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts) 1416e4dd521cSBarry Smith { 14171566a47fSLisandro Dalcin TS_RK *rk; 1418dfbe8321SBarry Smith PetscErrorCode ierr; 1419e4dd521cSBarry Smith 1420e4dd521cSBarry Smith PetscFunctionBegin; 1421f68a32c8SEmil Constantinescu ierr = TSRKInitializePackage();CHKERRQ(ierr); 1422f68a32c8SEmil Constantinescu 1423f68a32c8SEmil Constantinescu ts->ops->reset = TSReset_RK; 14245f70b478SJed Brown ts->ops->destroy = TSDestroy_RK; 14255f70b478SJed Brown ts->ops->view = TSView_RK; 1426f68a32c8SEmil Constantinescu ts->ops->load = TSLoad_RK; 1427f68a32c8SEmil Constantinescu ts->ops->setup = TSSetUp_RK; 142842f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_RK; 1429f68a32c8SEmil Constantinescu ts->ops->interpolate = TSInterpolate_RK; 14302c9cb062Sluzhanghpp ts->ops->step = TSStep_RK; 1431f68a32c8SEmil Constantinescu ts->ops->evaluatestep = TSEvaluateStep_RK; 1432fecfb714SLisandro Dalcin ts->ops->rollback = TSRollBack_RK; 1433f68a32c8SEmil Constantinescu ts->ops->setfromoptions = TSSetFromOptions_RK; 1434ff22ae23SHong Zhang ts->ops->getstages = TSGetStages_RK; 143542f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_RK; 1436e4dd521cSBarry Smith 14370f7a1166SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_RK; 14380f7a1166SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_RK; 14390f7a1166SHong Zhang 14401566a47fSLisandro Dalcin ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr); 14411566a47fSLisandro Dalcin ts->data = (void*)rk; 1442e4dd521cSBarry Smith 1443f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr); 1444f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr); 14452c9cb062Sluzhanghpp ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirateType_C",TSRKSetMultirateType);CHKERRQ(ierr); 1446be5899b3SLisandro Dalcin 1447be5899b3SLisandro Dalcin ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr); 1448e4dd521cSBarry Smith PetscFunctionReturn(0); 1449e4dd521cSBarry Smith } 1450