1e4dd521cSBarry Smith /* 2474dd773SHong Zhang Code for time stepping with the Runge-Kutta method 3f68a32c8SEmil Constantinescu 4f68a32c8SEmil Constantinescu Notes: 5474dd773SHong Zhang The general system is written as 6474dd773SHong Zhang 7474dd773SHong Zhang Udot = F(t,U) 8474dd773SHong Zhang 9e4dd521cSBarry Smith */ 102c9cb062Sluzhanghpp 11af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 12f68a32c8SEmil Constantinescu #include <petscdm.h> 13474dd773SHong Zhang #include <../src/ts/impls/explicit/rk/rk.h> 1421052c0cSHong Zhang #include <../src/ts/impls/explicit/rk/mrk.h> 15f68a32c8SEmil Constantinescu 16484bcdc7SDebojyoti Ghosh static TSRKType TSRKDefault = TSRK3BS; 17f68a32c8SEmil Constantinescu static PetscBool TSRKRegisterAllCalled; 18f68a32c8SEmil Constantinescu static PetscBool TSRKPackageInitialized; 19f68a32c8SEmil Constantinescu 20ab8985f5SHong Zhang static RKTableauLink RKTableauList; 21ab8985f5SHong Zhang 22f68a32c8SEmil Constantinescu /*MC 23287c2655SBarry Smith TSRK1FE - First order forward Euler scheme. 24262ff7bbSBarry Smith 25f68a32c8SEmil Constantinescu This method has one stage. 26f68a32c8SEmil Constantinescu 27287c2655SBarry Smith Options database: 289bd3e852SBarry Smith . -ts_rk_type 1fe 29287c2655SBarry Smith 30f68a32c8SEmil Constantinescu Level: advanced 31f68a32c8SEmil Constantinescu 32287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 33f68a32c8SEmil Constantinescu M*/ 34f68a32c8SEmil Constantinescu /*MC 35e7685601SHong Zhang TSRK2A - Second order RK scheme (Heun's method). 36f68a32c8SEmil Constantinescu 37f68a32c8SEmil Constantinescu This method has two stages. 38f68a32c8SEmil Constantinescu 39287c2655SBarry Smith Options database: 409bd3e852SBarry Smith . -ts_rk_type 2a 41287c2655SBarry Smith 42f68a32c8SEmil Constantinescu Level: advanced 43f68a32c8SEmil Constantinescu 44287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 45f68a32c8SEmil Constantinescu M*/ 46f68a32c8SEmil Constantinescu /*MC 47e7685601SHong Zhang TSRK2B - Second order RK scheme (the midpoint method). 48e7685601SHong Zhang 49e7685601SHong Zhang This method has two stages. 50e7685601SHong Zhang 51e7685601SHong Zhang Options database: 52e7685601SHong Zhang . -ts_rk_type 2b 53e7685601SHong Zhang 54e7685601SHong Zhang Level: advanced 55e7685601SHong Zhang 56e7685601SHong Zhang .seealso: TSRK, TSRKType, TSRKSetType() 57e7685601SHong Zhang M*/ 58e7685601SHong Zhang /*MC 59f68a32c8SEmil Constantinescu TSRK3 - Third order RK scheme. 60f68a32c8SEmil Constantinescu 61f68a32c8SEmil Constantinescu This method has three stages. 62f68a32c8SEmil Constantinescu 63287c2655SBarry Smith Options database: 649bd3e852SBarry Smith . -ts_rk_type 3 65287c2655SBarry Smith 66f68a32c8SEmil Constantinescu Level: advanced 67f68a32c8SEmil Constantinescu 68287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 69f68a32c8SEmil Constantinescu M*/ 70f68a32c8SEmil Constantinescu /*MC 712109b73fSDebojyoti Ghosh TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method. 722109b73fSDebojyoti Ghosh 73d1905564SLisandro Dalcin This method has four stages with the First Same As Last (FSAL) property. 742109b73fSDebojyoti Ghosh 75287c2655SBarry Smith Options database: 769bd3e852SBarry Smith . -ts_rk_type 3bs 77287c2655SBarry Smith 782109b73fSDebojyoti Ghosh Level: advanced 792109b73fSDebojyoti Ghosh 80*606c0280SSatish Balay References: 81*606c0280SSatish Balay . * - https://doi.org/10.1016/0893-9659(89)90079-7 82d1905564SLisandro Dalcin 83287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 842109b73fSDebojyoti Ghosh M*/ 852109b73fSDebojyoti Ghosh /*MC 86f68a32c8SEmil Constantinescu TSRK4 - Fourth order RK scheme. 87f68a32c8SEmil Constantinescu 882e698f8bSDebojyoti Ghosh This is the classical Runge-Kutta method with four stages. 89f68a32c8SEmil Constantinescu 90287c2655SBarry Smith Options database: 919bd3e852SBarry Smith . -ts_rk_type 4 92287c2655SBarry Smith 93f68a32c8SEmil Constantinescu Level: advanced 94f68a32c8SEmil Constantinescu 95287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 96f68a32c8SEmil Constantinescu M*/ 97f68a32c8SEmil Constantinescu /*MC 982e698f8bSDebojyoti Ghosh TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method. 99f68a32c8SEmil Constantinescu 100f68a32c8SEmil Constantinescu This method has six stages. 101f68a32c8SEmil Constantinescu 102287c2655SBarry Smith Options database: 1039bd3e852SBarry Smith . -ts_rk_type 5f 104287c2655SBarry Smith 105f68a32c8SEmil Constantinescu Level: advanced 106f68a32c8SEmil Constantinescu 107287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 108f68a32c8SEmil Constantinescu M*/ 1092109b73fSDebojyoti Ghosh /*MC 1102e698f8bSDebojyoti Ghosh TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method. 1112109b73fSDebojyoti Ghosh 112d1905564SLisandro Dalcin This method has seven stages with the First Same As Last (FSAL) property. 1132109b73fSDebojyoti Ghosh 114287c2655SBarry Smith Options database: 1159bd3e852SBarry Smith . -ts_rk_type 5dp 116287c2655SBarry Smith 1172109b73fSDebojyoti Ghosh Level: advanced 1182109b73fSDebojyoti Ghosh 119*606c0280SSatish Balay References: 120*606c0280SSatish Balay . * - https://doi.org/10.1016/0771-050X(80)90013-3 121d1905564SLisandro Dalcin 122287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 1232109b73fSDebojyoti Ghosh M*/ 12405e23783SLisandro Dalcin /*MC 12505e23783SLisandro Dalcin TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method. 12605e23783SLisandro Dalcin 12705e23783SLisandro Dalcin This method has eight stages with the First Same As Last (FSAL) property. 12805e23783SLisandro Dalcin 129287c2655SBarry Smith Options database: 1309bd3e852SBarry Smith . -ts_rk_type 5bs 131287c2655SBarry Smith 13205e23783SLisandro Dalcin Level: advanced 13305e23783SLisandro Dalcin 134*606c0280SSatish Balay References: 135*606c0280SSatish Balay . * - https://doi.org/10.1016/0898-1221(96)00141-1 13605e23783SLisandro Dalcin 137287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 13805e23783SLisandro Dalcin M*/ 13937acaa02SHendrik Ranocha /*MC 14037acaa02SHendrik Ranocha TSRK6VR - Sixth order robust Verner RK scheme with fifth order embedded method. 14137acaa02SHendrik Ranocha 14237acaa02SHendrik Ranocha This method has nine stages with the First Same As Last (FSAL) property. 14337acaa02SHendrik Ranocha 14437acaa02SHendrik Ranocha Options database: 14537acaa02SHendrik Ranocha . -ts_rk_type 6vr 14637acaa02SHendrik Ranocha 14737acaa02SHendrik Ranocha Level: advanced 14837acaa02SHendrik Ranocha 149*606c0280SSatish Balay References: 150*606c0280SSatish Balay . * - http://people.math.sfu.ca/~jverner/RKV65.IIIXb.Robust.00010102836.081204.CoeffsOnlyRAT 15137acaa02SHendrik Ranocha 15237acaa02SHendrik Ranocha .seealso: TSRK, TSRKType, TSRKSetType() 15337acaa02SHendrik Ranocha M*/ 15437acaa02SHendrik Ranocha /*MC 15537acaa02SHendrik Ranocha TSRK7VR - Seventh order robust Verner RK scheme with sixth order embedded method. 15637acaa02SHendrik Ranocha 1578f3d7ee2SCarsten Uphoff This method has ten stages. 15837acaa02SHendrik Ranocha 15937acaa02SHendrik Ranocha Options database: 16037acaa02SHendrik Ranocha . -ts_rk_type 7vr 16137acaa02SHendrik Ranocha 16237acaa02SHendrik Ranocha Level: advanced 16337acaa02SHendrik Ranocha 164*606c0280SSatish Balay References: 165*606c0280SSatish Balay . * - http://people.math.sfu.ca/~jverner/RKV76.IIa.Robust.000027015646.081206.CoeffsOnlyRAT 16637acaa02SHendrik Ranocha 16737acaa02SHendrik Ranocha .seealso: TSRK, TSRKType, TSRKSetType() 16837acaa02SHendrik Ranocha M*/ 16937acaa02SHendrik Ranocha /*MC 17037acaa02SHendrik Ranocha TSRK8VR - Eigth order robust Verner RK scheme with seventh order embedded method. 17137acaa02SHendrik Ranocha 1728f3d7ee2SCarsten Uphoff This method has thirteen stages. 17337acaa02SHendrik Ranocha 17437acaa02SHendrik Ranocha Options database: 17537acaa02SHendrik Ranocha . -ts_rk_type 8vr 17637acaa02SHendrik Ranocha 17737acaa02SHendrik Ranocha Level: advanced 17837acaa02SHendrik Ranocha 179*606c0280SSatish Balay References: 180*606c0280SSatish Balay . * - http://people.math.sfu.ca/~jverner/RKV87.IIa.Robust.00000754677.081208.CoeffsOnlyRATandFLOAT 18137acaa02SHendrik Ranocha 18237acaa02SHendrik Ranocha .seealso: TSRK, TSRKType, TSRKSetType() 18337acaa02SHendrik Ranocha M*/ 184262ff7bbSBarry Smith 185f68a32c8SEmil Constantinescu /*@C 186f68a32c8SEmil Constantinescu TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK 187262ff7bbSBarry Smith 188f68a32c8SEmil Constantinescu Not Collective, but should be called by all processes which will need the schemes to be registered 189262ff7bbSBarry Smith 190f68a32c8SEmil Constantinescu Level: advanced 191262ff7bbSBarry Smith 192f68a32c8SEmil Constantinescu .seealso: TSRKRegisterDestroy() 193262ff7bbSBarry Smith @*/ 194f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void) 195262ff7bbSBarry Smith { 1964ac538c5SBarry Smith PetscErrorCode ierr; 197262ff7bbSBarry Smith 198262ff7bbSBarry Smith PetscFunctionBegin; 199f68a32c8SEmil Constantinescu if (TSRKRegisterAllCalled) PetscFunctionReturn(0); 200f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_TRUE; 201e4dd521cSBarry Smith 2024e82626cSLisandro Dalcin #define RC PetscRealConstant 203e4dd521cSBarry Smith { 204f68a32c8SEmil Constantinescu const PetscReal 2054e82626cSLisandro Dalcin A[1][1] = {{0}}, 2064e82626cSLisandro Dalcin b[1] = {RC(1.0)}; 2073a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 208e4dd521cSBarry Smith } 209db046809SBarry Smith { 210f68a32c8SEmil Constantinescu const PetscReal 2114e82626cSLisandro Dalcin A[2][2] = {{0,0}, 2124e82626cSLisandro Dalcin {RC(1.0),0}}, 2134e82626cSLisandro Dalcin b[2] = {RC(0.5),RC(0.5)}, 2144e82626cSLisandro Dalcin bembed[2] = {RC(1.0),0}; 2153a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 216db046809SBarry Smith } 217f68a32c8SEmil Constantinescu { 218f68a32c8SEmil Constantinescu const PetscReal 219e7685601SHong Zhang A[2][2] = {{0,0}, 220e7685601SHong Zhang {RC(0.5),0}}, 221e7685601SHong Zhang b[2] = {0,RC(1.0)}, 222e7685601SHong Zhang ierr = TSRKRegister(TSRK2B,2,2,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 223e7685601SHong Zhang } 224e7685601SHong Zhang { 225e7685601SHong Zhang const PetscReal 22617f6384fSBarry Smith A[3][3] = {{0,0,0}, 2274e82626cSLisandro Dalcin {RC(2.0)/RC(3.0),0,0}, 2284e82626cSLisandro Dalcin {RC(-1.0)/RC(3.0),RC(1.0),0}}, 2294e82626cSLisandro Dalcin b[3] = {RC(0.25),RC(0.5),RC(0.25)}; 2303a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 231fdefd5e5SDebojyoti Ghosh } 232fdefd5e5SDebojyoti Ghosh { 233fdefd5e5SDebojyoti Ghosh const PetscReal 23417f6384fSBarry Smith A[4][4] = {{0,0,0,0}, 2354e82626cSLisandro Dalcin {RC(1.0)/RC(2.0),0,0,0}, 2364e82626cSLisandro Dalcin {0,RC(3.0)/RC(4.0),0,0}, 2374e82626cSLisandro Dalcin {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}}, 2384e82626cSLisandro Dalcin b[4] = {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}, 2394e82626cSLisandro 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)}; 2403a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 241db046809SBarry Smith } 242f68a32c8SEmil Constantinescu { 243f68a32c8SEmil Constantinescu const PetscReal 244f68a32c8SEmil Constantinescu A[4][4] = {{0,0,0,0}, 2454e82626cSLisandro Dalcin {RC(0.5),0,0,0}, 2464e82626cSLisandro Dalcin {0,RC(0.5),0,0}, 2474e82626cSLisandro Dalcin {0,0,RC(1.0),0}}, 2484e82626cSLisandro 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)}; 2493a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 250db046809SBarry Smith } 251f68a32c8SEmil Constantinescu { 252f68a32c8SEmil Constantinescu const PetscReal 25317f6384fSBarry Smith A[6][6] = {{0,0,0,0,0,0}, 2544e82626cSLisandro Dalcin {RC(0.25),0,0,0,0,0}, 2554e82626cSLisandro Dalcin {RC(3.0)/RC(32.0),RC(9.0)/RC(32.0),0,0,0,0}, 2564e82626cSLisandro Dalcin {RC(1932.0)/RC(2197.0),RC(-7200.0)/RC(2197.0),RC(7296.0)/RC(2197.0),0,0,0}, 2574e82626cSLisandro Dalcin {RC(439.0)/RC(216.0),RC(-8.0),RC(3680.0)/RC(513.0),RC(-845.0)/RC(4104.0),0,0}, 2584e82626cSLisandro 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}}, 2594e82626cSLisandro 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)}, 2604e82626cSLisandro 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}; 2613a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 262fdefd5e5SDebojyoti Ghosh } 263fdefd5e5SDebojyoti Ghosh { 264fdefd5e5SDebojyoti Ghosh const PetscReal 26517f6384fSBarry Smith A[7][7] = {{0,0,0,0,0,0,0}, 2664e82626cSLisandro Dalcin {RC(1.0)/RC(5.0),0,0,0,0,0,0}, 2674e82626cSLisandro Dalcin {RC(3.0)/RC(40.0),RC(9.0)/RC(40.0),0,0,0,0,0}, 2684e82626cSLisandro Dalcin {RC(44.0)/RC(45.0),RC(-56.0)/RC(15.0),RC(32.0)/RC(9.0),0,0,0,0}, 2694e82626cSLisandro 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}, 2704e82626cSLisandro 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}, 2714e82626cSLisandro 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}}, 2724e82626cSLisandro 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}, 273a685a763Sluzhanghpp 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)}, 274a685a763Sluzhanghpp 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)}, 275a685a763Sluzhanghpp {0,0,0,0,0}, 276a685a763Sluzhanghpp {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)}, 277a685a763Sluzhanghpp {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)}, 278a685a763Sluzhanghpp {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)}, 279a685a763Sluzhanghpp {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)}, 280a685a763Sluzhanghpp {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)}}; 281a685a763Sluzhanghpp ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,5,binterp[0]);CHKERRQ(ierr); 282f68a32c8SEmil Constantinescu } 28305e23783SLisandro Dalcin { 28405e23783SLisandro Dalcin const PetscReal 28517f6384fSBarry Smith A[8][8] = {{0,0,0,0,0,0,0,0}, 2864e82626cSLisandro Dalcin {RC(1.0)/RC(6.0),0,0,0,0,0,0,0}, 2874e82626cSLisandro Dalcin {RC(2.0)/RC(27.0),RC(4.0)/RC(27.0),0,0,0,0,0,0}, 2884e82626cSLisandro Dalcin {RC(183.0)/RC(1372.0),RC(-162.0)/RC(343.0),RC(1053.0)/RC(1372.0),0,0,0,0,0}, 2894e82626cSLisandro 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}, 2904e82626cSLisandro 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}, 2914e82626cSLisandro 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}, 2924e82626cSLisandro 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}}, 2934e82626cSLisandro 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}, 2944e82626cSLisandro 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)}; 29505e23783SLisandro Dalcin ierr = TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 29605e23783SLisandro Dalcin } 29737acaa02SHendrik Ranocha { 29837acaa02SHendrik Ranocha const PetscReal 29937acaa02SHendrik Ranocha A[9][9] = {{0,0,0,0,0,0,0,0,0}, 30063fe90e8SHendrik Ranocha {RC(1.8000000000000000000000000000000000000000e-01),0,0,0,0,0,0,0,0}, 30163fe90e8SHendrik Ranocha {RC(8.9506172839506172839506172839506172839506e-02),RC(7.7160493827160493827160493827160493827160e-02),0,0,0,0,0,0,0}, 30263fe90e8SHendrik Ranocha {RC(6.2500000000000000000000000000000000000000e-02),0,RC(1.8750000000000000000000000000000000000000e-01),0,0,0,0,0,0}, 30363fe90e8SHendrik Ranocha {RC(3.1651600000000000000000000000000000000000e-01),0,RC(-1.0449480000000000000000000000000000000000e+00),RC(1.2584320000000000000000000000000000000000e+00),0,0,0,0,0}, 30463fe90e8SHendrik Ranocha {RC(2.7232612736485626257225065566674305502508e-01),0,RC(-8.2513360323886639676113360323886639676113e-01),RC(1.0480917678812415654520917678812415654521e+00),RC(1.0471570799276856873679117969088177628396e-01),0,0,0,0}, 30563fe90e8SHendrik Ranocha {RC(-1.6699418599716514314329607278961797333198e-01),0,RC(6.3170850202429149797570850202429149797571e-01),RC(1.7461044552773876082146758838488161796432e-01),RC(-1.0665356459086066122525194734018680677781e+00),RC(1.2272108843537414965986394557823129251701e+00),0,0,0}, 30663fe90e8SHendrik Ranocha {RC(3.6423751686909581646423751686909581646424e-01),0,RC(-2.0404858299595141700404858299595141700405e-01),RC(-3.4883737816068643136312309244640071707741e-01),RC(3.2619323032856867443333608747142581729048e+00),RC(-2.7551020408163265306122448979591836734694e+00),RC(6.8181818181818181818181818181818181818182e-01),0,0}, 30763fe90e8SHendrik Ranocha {RC(7.6388888888888888888888888888888888888889e-02),0,0,RC(3.6940836940836940836940836940836940836941e-01),0,RC(2.4801587301587301587301587301587301587302e-01),RC(2.3674242424242424242424242424242424242424e-01),RC(6.9444444444444444444444444444444444444444e-02),0}}, 30863fe90e8SHendrik Ranocha b[9] = {RC(7.6388888888888888888888888888888888888889e-02),0,0,RC(3.6940836940836940836940836940836940836941e-01),0,RC(2.4801587301587301587301587301587301587302e-01),RC(2.3674242424242424242424242424242424242424e-01),RC(6.9444444444444444444444444444444444444444e-02),0}, 30963fe90e8SHendrik Ranocha bembed[9] = {RC(5.8700209643605870020964360587002096436059e-02),0,0,RC(4.8072562358276643990929705215419501133787e-01),RC(-8.5341242076919085578832094861228313083563e-01),RC(1.2046485260770975056689342403628117913832e+00),0,RC(-5.9242373072160306202859394348756050883710e-02),RC(1.6858043453788134639198468985703028256220e-01)}; 31037acaa02SHendrik Ranocha ierr = TSRKRegister(TSRK6VR,6,9,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 31137acaa02SHendrik Ranocha } 31237acaa02SHendrik Ranocha { 31337acaa02SHendrik Ranocha const PetscReal 31437acaa02SHendrik Ranocha A[10][10] = {{0,0,0,0,0,0,0,0,0,0}, 31563fe90e8SHendrik Ranocha {RC(5.0000000000000000000000000000000000000000e-03),0,0,0,0,0,0,0,0,0}, 31663fe90e8SHendrik Ranocha {RC(-1.0767901234567901234567901234567901234568e+00),RC(1.1856790123456790123456790123456790123457e+00),0,0,0,0,0,0,0,0}, 31763fe90e8SHendrik Ranocha {RC(4.0833333333333333333333333333333333333333e-02),0,RC(1.2250000000000000000000000000000000000000e-01),0,0,0,0,0,0,0}, 31863fe90e8SHendrik Ranocha {RC(6.3607142857142857142857142857142857142857e-01),0,RC(-2.4444642857142857142857142857142857142857e+00),RC(2.2633928571428571428571428571428571428571e+00),0,0,0,0,0,0}, 31963fe90e8SHendrik Ranocha {RC(-2.5351211079349245229256383554660215487207e+00),0,RC(1.0299374654449267920438514460756024913612e+01),RC(-7.9513032885990579949493217458266876536482e+00),RC(7.9301148923100592201226014271115261823800e-01),0,0,0,0,0}, 32063fe90e8SHendrik Ranocha {RC(1.0018765812524632961969196583094999808207e+00),0,RC(-4.1665712824423798331313938005470971453189e+00),RC(3.8343432929128642412552665218251378665197e+00),RC(-5.0233333560710847547464330228611765612403e-01),RC(6.6768474388416077115385092269857695410259e-01),0,0,0,0}, 32163fe90e8SHendrik Ranocha {RC(2.7255018354630767130333963819175005717348e+01),0,RC(-4.2004617278410638355318645443909295369611e+01),RC(-1.0535713126619489917921081600546526103722e+01),RC(8.0495536711411937147983652158926826634202e+01),RC(-6.7343882271790513468549075963212975640927e+01),RC(1.3048657610777937463471187029566964762710e+01),0,0,0}, 32263fe90e8SHendrik Ranocha {RC(-3.0397378057114965146943658658755763226883e+00),0,RC(1.0138161410329801111857946190709700150441e+01),RC(-6.4293056748647215721462825629555298064437e+00),RC(-1.5864371483408276587115312853798610579467e+00),RC(1.8921781841968424410864308909131353365021e+00),RC(1.9699335407608869061292360163336442838006e-02),RC(5.4416989827933235465102724247952572977903e-03),0,0}, 32363fe90e8SHendrik Ranocha {RC(-1.4449518916777735137351003179355712360517e+00),0,RC(8.0318913859955919224117033223019560435041e+00),RC(-7.5831741663401346820798883023671588604984e+00),RC(3.5816169353190074211247685442452878696855e+00),RC(-2.4369722632199529411183809065693752383733e+00),RC(8.5158999992326179339689766032486142173390e-01),0,0,0}}, 32463fe90e8SHendrik Ranocha b[10] = {RC(4.7425837833706756083569172717574534698932e-02),0,0,RC(2.5622361659370562659961727458274623448160e-01),RC(2.6951376833074206619473817258075952886764e-01),RC(1.2686622409092782845989138364739173247882e-01),RC(2.4887225942060071622046449427647492767292e-01),RC(3.0744837408200631335304388479099184768645e-03),RC(4.8023809989496943308189063347143123323209e-02),0}, 32563fe90e8SHendrik Ranocha bembed[10] = {RC(4.7485247699299631037531273805727961552268e-02),0,0,RC(2.5599412588690633297154918245905393870497e-01),RC(2.7058478081067688722530891099268135732387e-01),RC(1.2505618684425992913638822323746917920448e-01),RC(2.5204468723743860507184043820197442562182e-01),0,0,RC(4.8834971521418614557381971303093137592592e-02)}; 32637acaa02SHendrik Ranocha ierr = TSRKRegister(TSRK7VR,7,10,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 32737acaa02SHendrik Ranocha } 32837acaa02SHendrik Ranocha { 32937acaa02SHendrik Ranocha const PetscReal 33037acaa02SHendrik Ranocha A[13][13] = {{0,0,0,0,0,0,0,0,0,0,0,0,0}, 33163fe90e8SHendrik Ranocha {RC(2.5000000000000000000000000000000000000000e-01),0,0,0,0,0,0,0,0,0,0,0,0}, 33263fe90e8SHendrik Ranocha {RC(8.7400846504915232052686327594877411977046e-02),RC(2.5487604938654321753087950620345685135815e-02),0,0,0,0,0,0,0,0,0,0,0}, 33363fe90e8SHendrik Ranocha {RC(4.2333169291338582677165354330708661417323e-02),0,RC(1.2699950787401574803149606299212598425197e-01),0,0,0,0,0,0,0,0,0,0}, 33463fe90e8SHendrik Ranocha {RC(4.2609505888742261494881445237572274090942e-01),0,RC(-1.5987952846591523265427733230657181117089e+00),RC(1.5967002257717297115939588706899953707994e+00),0,0,0,0,0,0,0,0,0}, 33563fe90e8SHendrik Ranocha {RC(5.0719337296713929515090618138513639239329e-02),0,0,RC(2.5433377264600407582754714408877778031369e-01),RC(2.0394689005728199465736223777270858044698e-01),0,0,0,0,0,0,0,0}, 33663fe90e8SHendrik Ranocha {RC(-2.9000374717523110970388379285425896124091e-01),0,0,RC(1.3441873910260789889438681109414337003184e+00),RC(-2.8647779433614427309611103827036562829470e+00),RC(2.6775942995105948517211260646164815438695e+00),0,0,0,0,0,0,0}, 33763fe90e8SHendrik Ranocha {RC(9.8535011337993546469740402980727014284756e-02),0,0,0,RC(2.2192680630751384842024036498197387903583e-01),RC(-1.8140622911806994312690338288073952457474e-01),RC(1.0944411472562548236922614918038631254153e-02),0,0,0,0,0,0}, 33863fe90e8SHendrik Ranocha {RC(3.8711052545731144679444618165166373405645e-01),0,0,RC(-1.4424454974855277571256745553077927767173e+00),RC(2.9053981890699509317691346449233848441744e+00),RC(-1.8537710696301059290843332675811978025183e+00),RC(1.4003648098728154269497325109771241479223e-01),RC(5.7273940811495816575746774624447706488753e-01),0,0,0,0,0}, 33963fe90e8SHendrik Ranocha {RC(-1.6124403444439308100630016197913480595436e-01),0,0,RC(-1.7339602957358984083578404473962567894901e-01),RC(-1.3012892814065147406016812745172492529744e+00),RC(1.1379503751738617308558792131431003472124e+00),RC(-3.1747649663966880106923521138043024698980e-02),RC(9.3351293824933666439811064486056884856590e-01),RC(-8.3786318334733852703300855629616433201504e-02),0,0,0,0}, 34063fe90e8SHendrik Ranocha {RC(-1.9199444881589533281510804651483576073142e-02),0,0,RC(2.7330857265264284907942326254016124275617e-01),RC(-6.7534973206944372919691611210942380856240e-01),RC(3.4151849813846016071738489974728382711981e-01),RC(-6.7950064803375772478920516198524629391910e-02),RC(9.6591752247623878884265586491216376509746e-02),RC(1.3253082511182101180721038466545389951226e-01),RC(3.6854959360386113446906329951531666812946e-01),0,0,0}, 34163fe90e8SHendrik Ranocha {RC(6.0918774036452898676888412111588817784584e-01),0,0,RC(-2.2725690858980016768999800931413088399719e+00),RC(4.7578983426940290068155255881914785497547e+00),RC(-5.5161067066927584824294689667844248244842e+00),RC(2.9005963696801192709095818565946174378180e-01),RC(5.6914239633590368229109858454801849145630e-01),RC(7.9267957603321670271339916205893327579951e-01),RC(1.5473720453288822894126190771849898232047e-01),RC(1.6149708956621816247083215106334544434974e+00),0,0}, 34263fe90e8SHendrik Ranocha {RC(8.8735762208534719663211694051981022704884e-01),0,0,RC(-2.9754597821085367558513632804709301581977e+00),RC(5.6007170094881630597990392548350098923829e+00),RC(-5.9156074505366744680014930189941657351840e+00),RC(2.2029689156134927016879142540807638331238e-01),RC(1.0155097824462216666143271340902996997549e-01),RC(1.1514345647386055909780397752125850553556e+00),RC(1.9297101665271239396134361900805843653065e+00),0,0,0}}, 34363fe90e8SHendrik Ranocha b[13] = {RC(4.4729564666695714203015840429049382466467e-02),0,0,0,0,RC(1.5691033527708199813368698010726645409175e-01),RC(1.8460973408151637740702451873526277892035e-01),RC(2.2516380602086991042479419400350721970920e-01),RC(1.4794615651970234687005179885449141753736e-01),RC(7.6055542444955825269798361910336491012732e-02),RC(1.2277290235018619610824346315921437388535e-01),RC(4.1811958638991631583384842800871882376786e-02),0}, 34463fe90e8SHendrik Ranocha bembed[13] = {RC(4.5847111400495925878664730122010282095875e-02),0,0,0,0,RC(2.6231891404152387437443356584845803392392e-01),RC(1.9169372337852611904485738635688429008025e-01),RC(2.1709172327902618330978407422906448568196e-01),RC(1.2738189624833706796803169450656737867900e-01),RC(1.1510530385365326258240515750043192148894e-01),0,0,RC(4.0561327798437566841823391436583608050053e-02)}; 34537acaa02SHendrik Ranocha ierr = TSRKRegister(TSRK8VR,8,13,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 34637acaa02SHendrik Ranocha } 3474e82626cSLisandro Dalcin #undef RC 348db046809SBarry Smith PetscFunctionReturn(0); 349db046809SBarry Smith } 350db046809SBarry Smith 351f68a32c8SEmil Constantinescu /*@C 352f68a32c8SEmil Constantinescu TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister(). 353f68a32c8SEmil Constantinescu 354f68a32c8SEmil Constantinescu Not Collective 355f68a32c8SEmil Constantinescu 356f68a32c8SEmil Constantinescu Level: advanced 357f68a32c8SEmil Constantinescu 358f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll() 359f68a32c8SEmil Constantinescu @*/ 360f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void) 361e4dd521cSBarry Smith { 362f68a32c8SEmil Constantinescu PetscErrorCode ierr; 363f68a32c8SEmil Constantinescu RKTableauLink link; 364f68a32c8SEmil Constantinescu 365f68a32c8SEmil Constantinescu PetscFunctionBegin; 366f68a32c8SEmil Constantinescu while ((link = RKTableauList)) { 367f68a32c8SEmil Constantinescu RKTableau t = &link->tab; 368f68a32c8SEmil Constantinescu RKTableauList = link->next; 369f68a32c8SEmil Constantinescu ierr = PetscFree3(t->A,t->b,t->c);CHKERRQ(ierr); 370f68a32c8SEmil Constantinescu ierr = PetscFree(t->bembed);CHKERRQ(ierr); 371f68a32c8SEmil Constantinescu ierr = PetscFree(t->binterp);CHKERRQ(ierr); 372f68a32c8SEmil Constantinescu ierr = PetscFree(t->name);CHKERRQ(ierr); 373f68a32c8SEmil Constantinescu ierr = PetscFree(link);CHKERRQ(ierr); 374f68a32c8SEmil Constantinescu } 375f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_FALSE; 376f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 377f68a32c8SEmil Constantinescu } 378f68a32c8SEmil Constantinescu 379f68a32c8SEmil Constantinescu /*@C 380f68a32c8SEmil Constantinescu TSRKInitializePackage - This function initializes everything in the TSRK package. It is called 3818a690491SBarry Smith from TSInitializePackage(). 382f68a32c8SEmil Constantinescu 383f68a32c8SEmil Constantinescu Level: developer 384f68a32c8SEmil Constantinescu 385f68a32c8SEmil Constantinescu .seealso: PetscInitialize() 386f68a32c8SEmil Constantinescu @*/ 387f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void) 388f68a32c8SEmil Constantinescu { 389186e87acSLisandro Dalcin PetscErrorCode ierr; 390e4dd521cSBarry Smith 391e4dd521cSBarry Smith PetscFunctionBegin; 392f68a32c8SEmil Constantinescu if (TSRKPackageInitialized) PetscFunctionReturn(0); 393f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_TRUE; 394f68a32c8SEmil Constantinescu ierr = TSRKRegisterAll();CHKERRQ(ierr); 395f68a32c8SEmil Constantinescu ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr); 396f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 397f68a32c8SEmil Constantinescu } 398186e87acSLisandro Dalcin 399f68a32c8SEmil Constantinescu /*@C 400f68a32c8SEmil Constantinescu TSRKFinalizePackage - This function destroys everything in the TSRK package. It is 401f68a32c8SEmil Constantinescu called from PetscFinalize(). 402186e87acSLisandro Dalcin 403f68a32c8SEmil Constantinescu Level: developer 404f68a32c8SEmil Constantinescu 405f68a32c8SEmil Constantinescu .seealso: PetscFinalize() 406f68a32c8SEmil Constantinescu @*/ 407f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void) 408f68a32c8SEmil Constantinescu { 409f68a32c8SEmil Constantinescu PetscErrorCode ierr; 410f68a32c8SEmil Constantinescu 411f68a32c8SEmil Constantinescu PetscFunctionBegin; 412f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_FALSE; 413f68a32c8SEmil Constantinescu ierr = TSRKRegisterDestroy();CHKERRQ(ierr); 414f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 415f68a32c8SEmil Constantinescu } 416f68a32c8SEmil Constantinescu 417f68a32c8SEmil Constantinescu /*@C 418f68a32c8SEmil Constantinescu TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 419f68a32c8SEmil Constantinescu 420f68a32c8SEmil Constantinescu Not Collective, but the same schemes should be registered on all processes on which they will be used 421f68a32c8SEmil Constantinescu 422f68a32c8SEmil Constantinescu Input Parameters: 423f68a32c8SEmil Constantinescu + name - identifier for method 424f68a32c8SEmil Constantinescu . order - approximation order of method 425f68a32c8SEmil Constantinescu . s - number of stages, this is the dimension of the matrices below 426f68a32c8SEmil Constantinescu . A - stage coefficients (dimension s*s, row-major) 427f68a32c8SEmil Constantinescu . b - step completion table (dimension s; NULL to use last row of A) 428f68a32c8SEmil Constantinescu . c - abscissa (dimension s; NULL to use row sums of A) 429f68a32c8SEmil Constantinescu . bembed - completion table for embedded method (dimension s; NULL if not available) 4303a8a9803SLisandro Dalcin . p - Order of the interpolation scheme, equal to the number of columns of binterp 4313a8a9803SLisandro Dalcin - binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1) 432f68a32c8SEmil Constantinescu 433f68a32c8SEmil Constantinescu Notes: 434f68a32c8SEmil Constantinescu Several RK methods are provided, this function is only needed to create new methods. 435f68a32c8SEmil Constantinescu 436f68a32c8SEmil Constantinescu Level: advanced 437f68a32c8SEmil Constantinescu 438f68a32c8SEmil Constantinescu .seealso: TSRK 439f68a32c8SEmil Constantinescu @*/ 440f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s, 441f68a32c8SEmil Constantinescu const PetscReal A[],const PetscReal b[],const PetscReal c[], 4423a8a9803SLisandro Dalcin const PetscReal bembed[],PetscInt p,const PetscReal binterp[]) 443f68a32c8SEmil Constantinescu { 444f68a32c8SEmil Constantinescu PetscErrorCode ierr; 445f68a32c8SEmil Constantinescu RKTableauLink link; 446f68a32c8SEmil Constantinescu RKTableau t; 447f68a32c8SEmil Constantinescu PetscInt i,j; 448f68a32c8SEmil Constantinescu 449f68a32c8SEmil Constantinescu PetscFunctionBegin; 4503a8a9803SLisandro Dalcin PetscValidCharPointer(name,1); 4513a8a9803SLisandro Dalcin PetscValidRealPointer(A,4); 4523a8a9803SLisandro Dalcin if (b) PetscValidRealPointer(b,5); 4533a8a9803SLisandro Dalcin if (c) PetscValidRealPointer(c,6); 4543a8a9803SLisandro Dalcin if (bembed) PetscValidRealPointer(bembed,7); 4553a8a9803SLisandro Dalcin if (binterp || p > 1) PetscValidRealPointer(binterp,9); 4563a8a9803SLisandro Dalcin 4571d36bdfdSBarry Smith ierr = TSRKInitializePackage();CHKERRQ(ierr); 45895dccacaSBarry Smith ierr = PetscNew(&link);CHKERRQ(ierr); 459f68a32c8SEmil Constantinescu t = &link->tab; 4603a8a9803SLisandro Dalcin 461f68a32c8SEmil Constantinescu ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 462f68a32c8SEmil Constantinescu t->order = order; 463f68a32c8SEmil Constantinescu t->s = s; 464dcca6d9dSJed Brown ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr); 465580bdb30SBarry Smith ierr = PetscArraycpy(t->A,A,s*s);CHKERRQ(ierr); 466580bdb30SBarry Smith if (b) { ierr = PetscArraycpy(t->b,b,s);CHKERRQ(ierr); } 467f68a32c8SEmil Constantinescu else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i]; 468580bdb30SBarry Smith if (c) { ierr = PetscArraycpy(t->c,c,s);CHKERRQ(ierr); } 469f68a32c8SEmil 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]; 470d760c35bSDebojyoti Ghosh t->FSAL = PETSC_TRUE; 471d760c35bSDebojyoti Ghosh for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE; 4723a8a9803SLisandro Dalcin 473f68a32c8SEmil Constantinescu if (bembed) { 474785e854fSJed Brown ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr); 475580bdb30SBarry Smith ierr = PetscArraycpy(t->bembed,bembed,s);CHKERRQ(ierr); 476f68a32c8SEmil Constantinescu } 477f68a32c8SEmil Constantinescu 4783a8a9803SLisandro Dalcin if (!binterp) { p = 1; binterp = t->b; } 4793a8a9803SLisandro Dalcin t->p = p; 4803a8a9803SLisandro Dalcin ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr); 481580bdb30SBarry Smith ierr = PetscArraycpy(t->binterp,binterp,s*p);CHKERRQ(ierr); 4823a8a9803SLisandro Dalcin 483f68a32c8SEmil Constantinescu link->next = RKTableauList; 484f68a32c8SEmil Constantinescu RKTableauList = link; 485f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 486f68a32c8SEmil Constantinescu } 487f68a32c8SEmil Constantinescu 4880f7a28cdSStefano Zampini PetscErrorCode TSRKGetTableau_RK(TS ts, PetscInt *s, const PetscReal **A, const PetscReal **b, const PetscReal **c, const PetscReal **bembed, 4890f7a28cdSStefano Zampini PetscInt *p, const PetscReal **binterp, PetscBool *FSAL) 4900f7a28cdSStefano Zampini { 4910f7a28cdSStefano Zampini TS_RK *rk = (TS_RK*)ts->data; 4920f7a28cdSStefano Zampini RKTableau tab = rk->tableau; 4930f7a28cdSStefano Zampini 4940f7a28cdSStefano Zampini PetscFunctionBegin; 4950f7a28cdSStefano Zampini if (s) *s = tab->s; 4960f7a28cdSStefano Zampini if (A) *A = tab->A; 4970f7a28cdSStefano Zampini if (b) *b = tab->b; 4980f7a28cdSStefano Zampini if (c) *c = tab->c; 4990f7a28cdSStefano Zampini if (bembed) *bembed = tab->bembed; 5000f7a28cdSStefano Zampini if (p) *p = tab->p; 5010f7a28cdSStefano Zampini if (binterp) *binterp = tab->binterp; 5020f7a28cdSStefano Zampini if (FSAL) *FSAL = tab->FSAL; 5030f7a28cdSStefano Zampini PetscFunctionReturn(0); 5040f7a28cdSStefano Zampini } 5050f7a28cdSStefano Zampini 5060f7a28cdSStefano Zampini /*@C 5070f7a28cdSStefano Zampini TSRKGetTableau - Get info on the RK tableau 5080f7a28cdSStefano Zampini 5090f7a28cdSStefano Zampini Not Collective 5100f7a28cdSStefano Zampini 511f899ff85SJose E. Roman Input Parameter: 5120f7a28cdSStefano Zampini . ts - timestepping context 5130f7a28cdSStefano Zampini 5140f7a28cdSStefano Zampini Output Parameters: 5150f7a28cdSStefano Zampini + s - number of stages, this is the dimension of the matrices below 5160f7a28cdSStefano Zampini . A - stage coefficients (dimension s*s, row-major) 5170f7a28cdSStefano Zampini . b - step completion table (dimension s) 5180f7a28cdSStefano Zampini . c - abscissa (dimension s) 5190f7a28cdSStefano Zampini . bembed - completion table for embedded method (dimension s; NULL if not available) 5200f7a28cdSStefano Zampini . p - Order of the interpolation scheme, equal to the number of columns of binterp 5210f7a28cdSStefano Zampini . binterp - Coefficients of the interpolation formula (dimension s*p) 5220f7a28cdSStefano Zampini - FSAL - wheather or not the scheme has the First Same As Last property 5230f7a28cdSStefano Zampini 5240f7a28cdSStefano Zampini Level: developer 5250f7a28cdSStefano Zampini 5260f7a28cdSStefano Zampini .seealso: TSRK 5270f7a28cdSStefano Zampini @*/ 5280f7a28cdSStefano Zampini PetscErrorCode TSRKGetTableau(TS ts, PetscInt *s, const PetscReal **A, const PetscReal **b, const PetscReal **c, const PetscReal **bembed, 5290f7a28cdSStefano Zampini PetscInt *p, const PetscReal **binterp, PetscBool *FSAL) 5300f7a28cdSStefano Zampini { 5310f7a28cdSStefano Zampini PetscErrorCode ierr; 5320f7a28cdSStefano Zampini 5330f7a28cdSStefano Zampini PetscFunctionBegin; 5340f7a28cdSStefano Zampini PetscValidHeaderSpecific(ts,TS_CLASSID,1); 5350f7a28cdSStefano Zampini ierr = PetscUseMethod(ts,"TSRKGetTableau_C",(TS,PetscInt*,const PetscReal**,const PetscReal**,const PetscReal**,const PetscReal**, 5360f7a28cdSStefano Zampini PetscInt*,const PetscReal**,PetscBool*),(ts,s,A,b,c,bembed,p,binterp,FSAL));CHKERRQ(ierr); 5370f7a28cdSStefano Zampini PetscFunctionReturn(0); 5380f7a28cdSStefano Zampini } 5390f7a28cdSStefano Zampini 540e4dd521cSBarry Smith /* 5412c9cb062Sluzhanghpp This is for single-step RK method 542f68a32c8SEmil Constantinescu The step completion formula is 543e4dd521cSBarry Smith 544f68a32c8SEmil Constantinescu x1 = x0 + h b^T YdotRHS 545f68a32c8SEmil Constantinescu 546f68a32c8SEmil Constantinescu This function can be called before or after ts->vec_sol has been updated. 547f68a32c8SEmil Constantinescu Suppose we have a completion formula (b) and an embedded formula (be) of different order. 548f68a32c8SEmil Constantinescu We can write 549f68a32c8SEmil Constantinescu 550f68a32c8SEmil Constantinescu x1e = x0 + h be^T YdotRHS 551f68a32c8SEmil Constantinescu = x1 - h b^T YdotRHS + h be^T YdotRHS 552f68a32c8SEmil Constantinescu = x1 + h (be - b)^T YdotRHS 553f68a32c8SEmil Constantinescu 554f68a32c8SEmil Constantinescu so we can evaluate the method with different order even after the step has been optimistically completed. 555e4dd521cSBarry Smith */ 556f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done) 557f68a32c8SEmil Constantinescu { 558f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 559f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 560f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 561f68a32c8SEmil Constantinescu PetscReal h; 562f68a32c8SEmil Constantinescu PetscInt s = tab->s,j; 563f68a32c8SEmil Constantinescu PetscErrorCode ierr; 564f68a32c8SEmil Constantinescu 565f68a32c8SEmil Constantinescu PetscFunctionBegin; 566f68a32c8SEmil Constantinescu switch (rk->status) { 567f68a32c8SEmil Constantinescu case TS_STEP_INCOMPLETE: 568f68a32c8SEmil Constantinescu case TS_STEP_PENDING: 569f68a32c8SEmil Constantinescu h = ts->time_step; break; 570f68a32c8SEmil Constantinescu case TS_STEP_COMPLETE: 571be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 572f68a32c8SEmil Constantinescu default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 573f68a32c8SEmil Constantinescu } 574f68a32c8SEmil Constantinescu if (order == tab->order) { 575f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { 576f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 57790b67129SHong Zhang for (j=0; j<s; j++) w[j] = h*tab->b[j]/rk->dtratio; 578f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 579f68a32c8SEmil Constantinescu } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 580f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 581f68a32c8SEmil Constantinescu } else if (order == tab->order-1) { 582f68a32c8SEmil Constantinescu if (!tab->bembed) goto unavailable; 583f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { /*Complete with the embedded method (be)*/ 584f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 585f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 586f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 587f68a32c8SEmil Constantinescu } else { /*Rollback and re-complete using (be-b) */ 588f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 589f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 590f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 591f68a32c8SEmil Constantinescu } 592f68a32c8SEmil Constantinescu if (done) *done = PETSC_TRUE; 593f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 594f68a32c8SEmil Constantinescu } 595f68a32c8SEmil Constantinescu unavailable: 596f68a32c8SEmil Constantinescu if (done) *done = PETSC_FALSE; 59798921bdaSJacob Faibussowitsch else SETERRQ(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); 598f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 599f68a32c8SEmil Constantinescu } 600f68a32c8SEmil Constantinescu 6010f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts) 6020f7a1166SHong Zhang { 6030f7a1166SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 604cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 6050f7a1166SHong Zhang RKTableau tab = rk->tableau; 6060f7a1166SHong Zhang const PetscInt s = tab->s; 6070f7a1166SHong Zhang const PetscReal *b = tab->b,*c = tab->c; 6080f7a1166SHong Zhang Vec *Y = rk->Y; 6090f7a1166SHong Zhang PetscInt i; 6100f7a1166SHong Zhang PetscErrorCode ierr; 6110f7a1166SHong Zhang 6120f7a1166SHong Zhang PetscFunctionBegin; 613cd4cee2dSHong Zhang /* No need to backup quadts->vec_sol since it can be reverted in TSRollBack_RK */ 6140f7a1166SHong Zhang for (i=s-1; i>=0; i--) { 615cd4cee2dSHong Zhang /* Evolve quadrature TS solution to compute integrals */ 616cd4cee2dSHong Zhang ierr = TSComputeRHSFunction(quadts,rk->ptime+rk->time_step*c[i],Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 617cd4cee2dSHong Zhang ierr = VecAXPY(quadts->vec_sol,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 6180f7a1166SHong Zhang } 6190f7a1166SHong Zhang PetscFunctionReturn(0); 6200f7a1166SHong Zhang } 6210f7a1166SHong Zhang 6220f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts) 6230f7a1166SHong Zhang { 6240f7a1166SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 6250f7a1166SHong Zhang RKTableau tab = rk->tableau; 626cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 6270f7a1166SHong Zhang const PetscInt s = tab->s; 6280f7a1166SHong Zhang const PetscReal *b = tab->b,*c = tab->c; 6290f7a1166SHong Zhang Vec *Y = rk->Y; 6300f7a1166SHong Zhang PetscInt i; 6310f7a1166SHong Zhang PetscErrorCode ierr; 6320f7a1166SHong Zhang 6330f7a1166SHong Zhang PetscFunctionBegin; 6340f7a1166SHong Zhang for (i=s-1; i>=0; i--) { 635cd4cee2dSHong Zhang /* Evolve quadrature TS solution to compute integrals */ 636cd4cee2dSHong Zhang ierr = TSComputeRHSFunction(quadts,ts->ptime+ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 637cd4cee2dSHong Zhang ierr = VecAXPY(quadts->vec_sol,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 6380f7a1166SHong Zhang } 6390f7a1166SHong Zhang PetscFunctionReturn(0); 6400f7a1166SHong Zhang } 6410f7a1166SHong Zhang 642fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts) 643fecfb714SLisandro Dalcin { 644fecfb714SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 645cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 646fecfb714SLisandro Dalcin RKTableau tab = rk->tableau; 647fecfb714SLisandro Dalcin const PetscInt s = tab->s; 648cd4cee2dSHong Zhang const PetscReal *b = tab->b,*c = tab->c; 649fecfb714SLisandro Dalcin PetscScalar *w = rk->work; 650cd4cee2dSHong Zhang Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 651fecfb714SLisandro Dalcin PetscInt j; 652fecfb714SLisandro Dalcin PetscReal h; 653fecfb714SLisandro Dalcin PetscErrorCode ierr; 654fecfb714SLisandro Dalcin 655fecfb714SLisandro Dalcin PetscFunctionBegin; 656fecfb714SLisandro Dalcin switch (rk->status) { 657fecfb714SLisandro Dalcin case TS_STEP_INCOMPLETE: 658fecfb714SLisandro Dalcin case TS_STEP_PENDING: 659fecfb714SLisandro Dalcin h = ts->time_step; break; 660fecfb714SLisandro Dalcin case TS_STEP_COMPLETE: 661fecfb714SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 662fecfb714SLisandro Dalcin default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 663fecfb714SLisandro Dalcin } 664fecfb714SLisandro Dalcin for (j=0; j<s; j++) w[j] = -h*b[j]; 665fecfb714SLisandro Dalcin ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr); 666cd4cee2dSHong Zhang if (quadts && ts->costintegralfwd) { 667cd4cee2dSHong Zhang for (j=0; j<s; j++) { 668cd4cee2dSHong Zhang /* Revert the quadrature TS solution */ 669cd4cee2dSHong Zhang ierr = TSComputeRHSFunction(quadts,rk->ptime+h*c[j],Y[j],ts->vec_costintegrand);CHKERRQ(ierr); 670cd4cee2dSHong Zhang ierr = VecAXPY(quadts->vec_sol,-h*b[j],ts->vec_costintegrand);CHKERRQ(ierr); 671cd4cee2dSHong Zhang } 672cd4cee2dSHong Zhang } 673fecfb714SLisandro Dalcin PetscFunctionReturn(0); 674fecfb714SLisandro Dalcin } 675fecfb714SLisandro Dalcin 676922a638cSHong Zhang static PetscErrorCode TSForwardStep_RK(TS ts) 677922a638cSHong Zhang { 678922a638cSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 679922a638cSHong Zhang RKTableau tab = rk->tableau; 680922a638cSHong Zhang Mat J,*MatsFwdSensipTemp = rk->MatsFwdSensipTemp; 681922a638cSHong Zhang const PetscInt s = tab->s; 682922a638cSHong Zhang const PetscReal *A = tab->A,*c = tab->c,*b = tab->b; 683922a638cSHong Zhang Vec *Y = rk->Y; 684922a638cSHong Zhang PetscInt i,j; 685922a638cSHong Zhang PetscReal stage_time,h = ts->time_step; 686922a638cSHong Zhang PetscBool zero; 687922a638cSHong Zhang PetscErrorCode ierr; 688922a638cSHong Zhang 689922a638cSHong Zhang PetscFunctionBegin; 690922a638cSHong Zhang ierr = MatCopy(ts->mat_sensip,rk->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 691922a638cSHong Zhang ierr = TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr); 692922a638cSHong Zhang 693922a638cSHong Zhang for (i=0; i<s; i++) { 694922a638cSHong Zhang stage_time = ts->ptime + h*c[i]; 695922a638cSHong Zhang zero = PETSC_FALSE; 696922a638cSHong Zhang if (b[i] == 0 && i == s-1) zero = PETSC_TRUE; 697922a638cSHong Zhang /* TLM Stage values */ 698922a638cSHong Zhang if (!i) { 699922a638cSHong Zhang ierr = MatCopy(ts->mat_sensip,rk->MatsFwdStageSensip[i],SAME_NONZERO_PATTERN);CHKERRQ(ierr); 700922a638cSHong Zhang } else if (!zero) { 701922a638cSHong Zhang ierr = MatZeroEntries(rk->MatsFwdStageSensip[i]);CHKERRQ(ierr); 702922a638cSHong Zhang for (j=0; j<i; j++) { 703922a638cSHong Zhang ierr = MatAXPY(rk->MatsFwdStageSensip[i],h*A[i*s+j],MatsFwdSensipTemp[j],SAME_NONZERO_PATTERN);CHKERRQ(ierr); 704922a638cSHong Zhang } 705922a638cSHong Zhang ierr = MatAXPY(rk->MatsFwdStageSensip[i],1.,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 706922a638cSHong Zhang } else { 707922a638cSHong Zhang ierr = MatZeroEntries(rk->MatsFwdStageSensip[i]);CHKERRQ(ierr); 708922a638cSHong Zhang } 709922a638cSHong Zhang 710922a638cSHong Zhang ierr = TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);CHKERRQ(ierr); 711922a638cSHong Zhang ierr = MatMatMult(J,rk->MatsFwdStageSensip[i],MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatsFwdSensipTemp[i]);CHKERRQ(ierr); 71213af1a74SHong Zhang if (ts->Jacprhs) { 71313af1a74SHong Zhang ierr = TSComputeRHSJacobianP(ts,stage_time,Y[i],ts->Jacprhs);CHKERRQ(ierr); /* get f_p */ 71413af1a74SHong Zhang if (ts->vecs_sensi2p) { /* TLM used for 2nd-order adjoint */ 71513af1a74SHong Zhang PetscScalar *xarr; 71613af1a74SHong Zhang ierr = MatDenseGetColumn(MatsFwdSensipTemp[i],0,&xarr);CHKERRQ(ierr); 71713af1a74SHong Zhang ierr = VecPlaceArray(rk->VecDeltaFwdSensipCol,xarr);CHKERRQ(ierr); 71813af1a74SHong Zhang ierr = MatMultAdd(ts->Jacprhs,ts->vec_dir,rk->VecDeltaFwdSensipCol,rk->VecDeltaFwdSensipCol);CHKERRQ(ierr); 719c3b366b1Sprj- ierr = VecResetArray(rk->VecDeltaFwdSensipCol);CHKERRQ(ierr); 72013af1a74SHong Zhang ierr = MatDenseRestoreColumn(MatsFwdSensipTemp[i],&xarr);CHKERRQ(ierr); 72113af1a74SHong Zhang } else { 72213af1a74SHong Zhang ierr = MatAXPY(MatsFwdSensipTemp[i],1.,ts->Jacprhs,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 72313af1a74SHong Zhang } 724922a638cSHong Zhang } 725922a638cSHong Zhang } 726922a638cSHong Zhang 727922a638cSHong Zhang for (i=0; i<s; i++) { 728922a638cSHong Zhang ierr = MatAXPY(ts->mat_sensip,h*b[i],rk->MatsFwdSensipTemp[i],SAME_NONZERO_PATTERN);CHKERRQ(ierr); 729922a638cSHong Zhang } 730922a638cSHong Zhang rk->status = TS_STEP_COMPLETE; 731922a638cSHong Zhang PetscFunctionReturn(0); 732922a638cSHong Zhang } 733922a638cSHong Zhang 734922a638cSHong Zhang static PetscErrorCode TSForwardGetStages_RK(TS ts,PetscInt *ns,Mat **stagesensip) 735922a638cSHong Zhang { 736922a638cSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 737922a638cSHong Zhang RKTableau tab = rk->tableau; 738922a638cSHong Zhang 739922a638cSHong Zhang PetscFunctionBegin; 740922a638cSHong Zhang if (ns) *ns = tab->s; 741922a638cSHong Zhang if (stagesensip) *stagesensip = rk->MatsFwdStageSensip; 742922a638cSHong Zhang PetscFunctionReturn(0); 743922a638cSHong Zhang } 744922a638cSHong Zhang 745922a638cSHong Zhang static PetscErrorCode TSForwardSetUp_RK(TS ts) 746922a638cSHong Zhang { 747922a638cSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 748922a638cSHong Zhang RKTableau tab = rk->tableau; 749922a638cSHong Zhang PetscInt i; 750922a638cSHong Zhang PetscErrorCode ierr; 751922a638cSHong Zhang 752922a638cSHong Zhang PetscFunctionBegin; 753922a638cSHong Zhang /* backup sensitivity results for roll-backs */ 754922a638cSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatFwdSensip0);CHKERRQ(ierr); 755922a638cSHong Zhang 756922a638cSHong Zhang ierr = PetscMalloc1(tab->s,&rk->MatsFwdStageSensip);CHKERRQ(ierr); 757922a638cSHong Zhang ierr = PetscMalloc1(tab->s,&rk->MatsFwdSensipTemp);CHKERRQ(ierr); 758922a638cSHong Zhang for (i=0; i<tab->s; i++) { 759922a638cSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatsFwdStageSensip[i]);CHKERRQ(ierr); 760922a638cSHong Zhang ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatsFwdSensipTemp[i]);CHKERRQ(ierr); 761922a638cSHong Zhang } 762922a638cSHong Zhang ierr = VecDuplicate(ts->vec_sol,&rk->VecDeltaFwdSensipCol);CHKERRQ(ierr); 763922a638cSHong Zhang PetscFunctionReturn(0); 764922a638cSHong Zhang } 765922a638cSHong Zhang 766922a638cSHong Zhang static PetscErrorCode TSForwardReset_RK(TS ts) 767922a638cSHong Zhang { 768922a638cSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 769922a638cSHong Zhang RKTableau tab = rk->tableau; 770922a638cSHong Zhang PetscInt i; 771922a638cSHong Zhang PetscErrorCode ierr; 772922a638cSHong Zhang 773922a638cSHong Zhang PetscFunctionBegin; 774922a638cSHong Zhang ierr = MatDestroy(&rk->MatFwdSensip0);CHKERRQ(ierr); 775922a638cSHong Zhang if (rk->MatsFwdStageSensip) { 776922a638cSHong Zhang for (i=0; i<tab->s; i++) { 777922a638cSHong Zhang ierr = MatDestroy(&rk->MatsFwdStageSensip[i]);CHKERRQ(ierr); 778922a638cSHong Zhang } 779922a638cSHong Zhang ierr = PetscFree(rk->MatsFwdStageSensip);CHKERRQ(ierr); 780922a638cSHong Zhang } 781922a638cSHong Zhang if (rk->MatsFwdSensipTemp) { 782922a638cSHong Zhang for (i=0; i<tab->s; i++) { 783922a638cSHong Zhang ierr = MatDestroy(&rk->MatsFwdSensipTemp[i]);CHKERRQ(ierr); 784922a638cSHong Zhang } 785922a638cSHong Zhang ierr = PetscFree(rk->MatsFwdSensipTemp);CHKERRQ(ierr); 786922a638cSHong Zhang } 787922a638cSHong Zhang ierr = VecDestroy(&rk->VecDeltaFwdSensipCol);CHKERRQ(ierr); 788922a638cSHong Zhang PetscFunctionReturn(0); 789922a638cSHong Zhang } 790922a638cSHong Zhang 791f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts) 792f68a32c8SEmil Constantinescu { 793f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 794f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 795f68a32c8SEmil Constantinescu const PetscInt s = tab->s; 796fecfb714SLisandro Dalcin const PetscReal *A = tab->A,*c = tab->c; 797f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 798f68a32c8SEmil Constantinescu Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 799d1905564SLisandro Dalcin PetscBool FSAL = tab->FSAL; 800f68a32c8SEmil Constantinescu TSAdapt adapt; 801fecfb714SLisandro Dalcin PetscInt i,j; 802be5899b3SLisandro Dalcin PetscInt rejections = 0; 803be5899b3SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 804be5899b3SLisandro Dalcin PetscReal next_time_step = ts->time_step; 805f68a32c8SEmil Constantinescu PetscErrorCode ierr; 806f68a32c8SEmil Constantinescu 807f68a32c8SEmil Constantinescu PetscFunctionBegin; 808d1905564SLisandro Dalcin if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE; 809d1905564SLisandro Dalcin if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); } 810d1905564SLisandro Dalcin 811f68a32c8SEmil Constantinescu rk->status = TS_STEP_INCOMPLETE; 812be5899b3SLisandro Dalcin while (!ts->reason && rk->status != TS_STEP_COMPLETE) { 813be5899b3SLisandro Dalcin PetscReal t = ts->ptime; 814f68a32c8SEmil Constantinescu PetscReal h = ts->time_step; 815f68a32c8SEmil Constantinescu for (i=0; i<s; i++) { 8169be3e283SDebojyoti Ghosh rk->stage_time = t + h*c[i]; 8179be3e283SDebojyoti Ghosh ierr = TSPreStage(ts,rk->stage_time);CHKERRQ(ierr); 818f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 819f68a32c8SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 820f68a32c8SEmil Constantinescu ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 8219be3e283SDebojyoti Ghosh ierr = TSPostStage(ts,rk->stage_time,i,Y);CHKERRQ(ierr); 822f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 823be5899b3SLisandro Dalcin ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr); 824fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 8258f738a7cSLisandro Dalcin if (FSAL && !i) continue; 826f68a32c8SEmil Constantinescu ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 827f68a32c8SEmil Constantinescu } 828be5899b3SLisandro Dalcin 829be5899b3SLisandro Dalcin rk->status = TS_STEP_INCOMPLETE; 830f68a32c8SEmil Constantinescu ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 831f68a32c8SEmil Constantinescu rk->status = TS_STEP_PENDING; 832f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 833f68a32c8SEmil Constantinescu ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 8341917a363SLisandro Dalcin ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr); 835fecfb714SLisandro Dalcin ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 836be5899b3SLisandro Dalcin rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 837be5899b3SLisandro Dalcin if (!accept) { /* Roll back the current step */ 838fecfb714SLisandro Dalcin ierr = TSRollBack_RK(ts);CHKERRQ(ierr); 839be5899b3SLisandro Dalcin ts->time_step = next_time_step; 840be5899b3SLisandro Dalcin goto reject_step; 841be5899b3SLisandro Dalcin } 842be5899b3SLisandro Dalcin 8430f7a1166SHong Zhang if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 8440f7a1166SHong Zhang rk->ptime = ts->ptime; 8450f7a1166SHong Zhang rk->time_step = ts->time_step; 846493ed95fSHong Zhang } 847be5899b3SLisandro Dalcin 848f68a32c8SEmil Constantinescu ts->ptime += ts->time_step; 849f68a32c8SEmil Constantinescu ts->time_step = next_time_step; 850f68a32c8SEmil Constantinescu break; 851be5899b3SLisandro Dalcin 852be5899b3SLisandro Dalcin reject_step: 853fecfb714SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 854be5899b3SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 855be5899b3SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 8567d3de750SJacob Faibussowitsch ierr = PetscInfo(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 857f68a32c8SEmil Constantinescu } 858f68a32c8SEmil Constantinescu } 859f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 860e4dd521cSBarry Smith } 861e4dd521cSBarry Smith 862f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts) 863f6a906c0SBarry Smith { 864f6a906c0SBarry Smith TS_RK *rk = (TS_RK*)ts->data; 865be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 866be5899b3SLisandro Dalcin PetscInt s = tab->s; 867f6a906c0SBarry Smith PetscErrorCode ierr; 868f6a906c0SBarry Smith 869f6a906c0SBarry Smith PetscFunctionBegin; 870f6a906c0SBarry Smith if (ts->adjointsetupcalled++) PetscFunctionReturn(0); 8712e7b7f96SHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecsDeltaLam);CHKERRQ(ierr); 8722e7b7f96SHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&rk->VecsSensiTemp);CHKERRQ(ierr); 873f6a906c0SBarry Smith if (ts->vecs_sensip) { 8742e7b7f96SHong Zhang ierr = VecDuplicate(ts->vecs_sensip[0],&rk->VecDeltaMu);CHKERRQ(ierr); 875f6a906c0SBarry Smith } 87613af1a74SHong Zhang if (ts->vecs_sensi2) { 87713af1a74SHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecsDeltaLam2);CHKERRQ(ierr); 87813af1a74SHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&rk->VecsSensi2Temp);CHKERRQ(ierr); 87913af1a74SHong Zhang } 88013af1a74SHong Zhang if (ts->vecs_sensi2p) { 88113af1a74SHong Zhang ierr = VecDuplicate(ts->vecs_sensi2p[0],&rk->VecDeltaMu2);CHKERRQ(ierr); 88213af1a74SHong Zhang } 883f6a906c0SBarry Smith PetscFunctionReturn(0); 884f6a906c0SBarry Smith } 885f6a906c0SBarry Smith 88635f5172aSHong Zhang /* 88735f5172aSHong Zhang Assumptions: 88835f5172aSHong Zhang - TSStep_RK() always evaluates the step with b, not bembed. 88935f5172aSHong Zhang */ 89042f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts) 891d2daff3dSHong Zhang { 892c235aa8dSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 893cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 894c235aa8dSHong Zhang RKTableau tab = rk->tableau; 8953ca0882eSHong Zhang Mat J,Jpre,Jquad; 896c235aa8dSHong Zhang const PetscInt s = tab->s; 897c235aa8dSHong Zhang const PetscReal *A = tab->A,*b = tab->b,*c = tab->c; 89813af1a74SHong Zhang PetscScalar *w = rk->work,*xarr; 8992e7b7f96SHong Zhang Vec *Y = rk->Y,*VecsDeltaLam = rk->VecsDeltaLam,VecDeltaMu = rk->VecDeltaMu,*VecsSensiTemp = rk->VecsSensiTemp; 90013af1a74SHong Zhang Vec *VecsDeltaLam2 = rk->VecsDeltaLam2,VecDeltaMu2 = rk->VecDeltaMu2,*VecsSensi2Temp = rk->VecsSensi2Temp; 901cd4cee2dSHong Zhang Vec VecDRDUTransCol = ts->vec_drdu_col,VecDRDPTransCol = ts->vec_drdp_col; 902b18ea86cSHong Zhang PetscInt i,j,nadj; 9033d3eaba7SBarry Smith PetscReal t = ts->ptime; 9043ca0882eSHong Zhang PetscReal h = ts->time_step; 905d2daff3dSHong Zhang PetscErrorCode ierr; 906c235aa8dSHong Zhang 907d2daff3dSHong Zhang PetscFunctionBegin; 908c235aa8dSHong Zhang rk->status = TS_STEP_INCOMPLETE; 9093389c7d9SStefano Zampini 9103ca0882eSHong Zhang ierr = TSGetRHSJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr); 91135f5172aSHong Zhang if (quadts) { 912cd4cee2dSHong Zhang ierr = TSGetRHSJacobian(quadts,&Jquad,NULL,NULL,NULL);CHKERRQ(ierr); 91335f5172aSHong Zhang } 9142e7b7f96SHong Zhang for (i=s-1; i>=0; i--) { 9156a1e4597SHong Zhang if (tab->FSAL && i == s-1) { 9166a1e4597SHong Zhang /* VecsDeltaLam[nadj*s+s-1] are initialized with zeros and the values never change.*/ 9176a1e4597SHong Zhang continue; 9186a1e4597SHong Zhang } 9193ca0882eSHong Zhang rk->stage_time = t + h*(1.0-c[i]); 920fd850964SHong Zhang ierr = TSComputeSNESJacobian(ts,Y[i],J,Jpre);CHKERRQ(ierr); 92135f5172aSHong Zhang if (quadts) { 9223ca0882eSHong Zhang ierr = TSComputeRHSJacobian(quadts,rk->stage_time,Y[i],Jquad,Jquad);CHKERRQ(ierr); /* get r_u^T */ 92335f5172aSHong Zhang } 9243389c7d9SStefano Zampini if (ts->vecs_sensip) { 9253ca0882eSHong Zhang ierr = TSComputeRHSJacobianP(ts,rk->stage_time,Y[i],ts->Jacprhs);CHKERRQ(ierr); /* get f_p */ 92635f5172aSHong Zhang if (quadts) { 9273ca0882eSHong Zhang ierr = TSComputeRHSJacobianP(quadts,rk->stage_time,Y[i],quadts->Jacprhs);CHKERRQ(ierr); /* get f_p for the quadrature */ 9283389c7d9SStefano Zampini } 92935f5172aSHong Zhang } 9303389c7d9SStefano Zampini 93135f5172aSHong Zhang if (b[i]) { 93235f5172aSHong Zhang for (j=i+1; j<s; j++) w[j-i-1] = A[j*s+i]/b[i]; /* coefficients for computing VecsSensiTemp */ 93335f5172aSHong Zhang } else { 93435f5172aSHong Zhang for (j=i+1; j<s; j++) w[j-i-1] = A[j*s+i]; /* coefficients for computing VecsSensiTemp */ 93535f5172aSHong Zhang } 9362e7b7f96SHong Zhang 9372e7b7f96SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 9383389c7d9SStefano Zampini /* Stage values of lambda */ 93935f5172aSHong Zhang if (b[i]) { 94035f5172aSHong Zhang /* lambda_{n+1} + \sum_{j=i+1}^s a_{ji}/b[i]*lambda_{s,j} */ 9412e7b7f96SHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); /* VecDeltaLam is an vec array of size s by numcost */ 9422e7b7f96SHong Zhang ierr = VecMAXPY(VecsSensiTemp[nadj],s-i-1,w,&VecsDeltaLam[nadj*s+i+1]);CHKERRQ(ierr); 94335f5172aSHong Zhang ierr = MatMultTranspose(J,VecsSensiTemp[nadj],VecsDeltaLam[nadj*s+i]);CHKERRQ(ierr); /* VecsSensiTemp will be reused by 2nd-order adjoint */ 94435f5172aSHong Zhang ierr = VecScale(VecsDeltaLam[nadj*s+i],-h*b[i]);CHKERRQ(ierr); 94535f5172aSHong Zhang if (quadts) { 946cd4cee2dSHong Zhang ierr = MatDenseGetColumn(Jquad,nadj,&xarr);CHKERRQ(ierr); 947cd4cee2dSHong Zhang ierr = VecPlaceArray(VecDRDUTransCol,xarr);CHKERRQ(ierr); 948cd4cee2dSHong Zhang ierr = VecAXPY(VecsDeltaLam[nadj*s+i],-h*b[i],VecDRDUTransCol);CHKERRQ(ierr); 949cd4cee2dSHong Zhang ierr = VecResetArray(VecDRDUTransCol);CHKERRQ(ierr); 950cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(Jquad,&xarr);CHKERRQ(ierr); 951cd4cee2dSHong Zhang } 9523389c7d9SStefano Zampini } else { 9536a1e4597SHong Zhang /* \sum_{j=i+1}^s a_{ji}*lambda_{s,j} */ 9546a1e4597SHong Zhang ierr = VecSet(VecsSensiTemp[nadj],0);CHKERRQ(ierr); 9556a1e4597SHong Zhang ierr = VecMAXPY(VecsSensiTemp[nadj],s-i-1,w,&VecsDeltaLam[nadj*s+i+1]);CHKERRQ(ierr); 9566a1e4597SHong Zhang ierr = MatMultTranspose(J,VecsSensiTemp[nadj],VecsDeltaLam[nadj*s+i]);CHKERRQ(ierr); 9576a1e4597SHong Zhang ierr = VecScale(VecsDeltaLam[nadj*s+i],-h);CHKERRQ(ierr); 9583389c7d9SStefano Zampini } 9596497ce14SHong Zhang 960ad8e2604SHong Zhang /* Stage values of mu */ 9616497ce14SHong Zhang if (ts->vecs_sensip) { 96235f5172aSHong Zhang if (b[i]) { 963d9227288SHong Zhang ierr = MatMultTranspose(ts->Jacprhs,VecsSensiTemp[nadj],VecDeltaMu);CHKERRQ(ierr); 96435f5172aSHong Zhang ierr = VecScale(VecDeltaMu,-h*b[i]);CHKERRQ(ierr); 96535f5172aSHong Zhang if (quadts) { 966cd4cee2dSHong Zhang ierr = MatDenseGetColumn(quadts->Jacprhs,nadj,&xarr);CHKERRQ(ierr); 967cd4cee2dSHong Zhang ierr = VecPlaceArray(VecDRDPTransCol,xarr);CHKERRQ(ierr); 968cd4cee2dSHong Zhang ierr = VecAXPY(VecDeltaMu,-h*b[i],VecDRDPTransCol);CHKERRQ(ierr); 969cd4cee2dSHong Zhang ierr = VecResetArray(VecDRDPTransCol);CHKERRQ(ierr); 970cd4cee2dSHong Zhang ierr = MatDenseRestoreColumn(quadts->Jacprhs,&xarr);CHKERRQ(ierr); 971493ed95fSHong Zhang } 97235f5172aSHong Zhang } else { 97335f5172aSHong Zhang ierr = VecScale(VecDeltaMu,-h);CHKERRQ(ierr); 97435f5172aSHong Zhang } 9752e7b7f96SHong Zhang ierr = VecAXPY(ts->vecs_sensip[nadj],1.,VecDeltaMu);CHKERRQ(ierr); /* update sensip for each stage */ 976ad8e2604SHong Zhang } 977c235aa8dSHong Zhang } 97813af1a74SHong Zhang 97913af1a74SHong Zhang if (ts->vecs_sensi2 && ts->forward_solve) { /* 2nd-order adjoint, TLM mode has to be turned on */ 98013af1a74SHong Zhang /* Get w1 at t_{n+1} from TLM matrix */ 98113af1a74SHong Zhang ierr = MatDenseGetColumn(rk->MatsFwdStageSensip[i],0,&xarr);CHKERRQ(ierr); 98213af1a74SHong Zhang ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 98313af1a74SHong Zhang /* lambda_s^T F_UU w_1 */ 9843ca0882eSHong Zhang ierr = TSComputeRHSHessianProductFunctionUU(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_sensip_col,ts->vecs_guu);CHKERRQ(ierr); 98535f5172aSHong Zhang if (quadts) { 98635f5172aSHong Zhang /* R_UU w_1 */ 9873ca0882eSHong Zhang ierr = TSComputeRHSHessianProductFunctionUU(quadts,rk->stage_time,Y[i],NULL,ts->vec_sensip_col,ts->vecs_guu);CHKERRQ(ierr); 98835f5172aSHong Zhang } 98935f5172aSHong Zhang if (ts->vecs_sensip) { 99013af1a74SHong Zhang /* lambda_s^T F_UP w_2 */ 9913ca0882eSHong Zhang ierr = TSComputeRHSHessianProductFunctionUP(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_dir,ts->vecs_gup);CHKERRQ(ierr); 99235f5172aSHong Zhang if (quadts) { 99335f5172aSHong Zhang /* R_UP w_2 */ 9943ca0882eSHong Zhang ierr = TSComputeRHSHessianProductFunctionUP(quadts,rk->stage_time,Y[i],NULL,ts->vec_sensip_col,ts->vecs_gup);CHKERRQ(ierr); 99535f5172aSHong Zhang } 99635f5172aSHong Zhang } 99713af1a74SHong Zhang if (ts->vecs_sensi2p) { 99813af1a74SHong Zhang /* lambda_s^T F_PU w_1 */ 9993ca0882eSHong Zhang ierr = TSComputeRHSHessianProductFunctionPU(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_sensip_col,ts->vecs_gpu);CHKERRQ(ierr); 100035f5172aSHong Zhang /* lambda_s^T F_PP w_2 */ 10013ca0882eSHong Zhang ierr = TSComputeRHSHessianProductFunctionPP(ts,rk->stage_time,Y[i],VecsSensiTemp,ts->vec_dir,ts->vecs_gpp);CHKERRQ(ierr); 100235f5172aSHong Zhang if (b[i] && quadts) { 100335f5172aSHong Zhang /* R_PU w_1 */ 10043ca0882eSHong Zhang ierr = TSComputeRHSHessianProductFunctionPU(quadts,rk->stage_time,Y[i],NULL,ts->vec_sensip_col,ts->vecs_gpu);CHKERRQ(ierr); 100535f5172aSHong Zhang /* R_PP w_2 */ 10063ca0882eSHong Zhang ierr = TSComputeRHSHessianProductFunctionPP(quadts,rk->stage_time,Y[i],NULL,ts->vec_dir,ts->vecs_gpp);CHKERRQ(ierr); 100735f5172aSHong Zhang } 100813af1a74SHong Zhang } 100913af1a74SHong Zhang ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 101013af1a74SHong Zhang ierr = MatDenseRestoreColumn(rk->MatsFwdStageSensip[i],&xarr);CHKERRQ(ierr); 101113af1a74SHong Zhang 101213af1a74SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 101313af1a74SHong Zhang /* Stage values of lambda */ 101435f5172aSHong Zhang if (b[i]) { 101535f5172aSHong Zhang /* J_i^T*(Lambda_{n+1}+\sum_{j=i+1}^s a_{ji}/b_i*Lambda_{s,j} */ 101613af1a74SHong Zhang ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 101713af1a74SHong Zhang ierr = VecMAXPY(VecsSensi2Temp[nadj],s-i-1,w,&VecsDeltaLam2[nadj*s+i+1]);CHKERRQ(ierr); 101813af1a74SHong Zhang ierr = MatMultTranspose(J,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj*s+i]);CHKERRQ(ierr); 101935f5172aSHong Zhang ierr = VecScale(VecsDeltaLam2[nadj*s+i],-h*b[i]);CHKERRQ(ierr); 102035f5172aSHong Zhang ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],-h*b[i],ts->vecs_guu[nadj]);CHKERRQ(ierr); 102135f5172aSHong Zhang if (ts->vecs_sensip) { 102235f5172aSHong Zhang ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],-h*b[i],ts->vecs_gup[nadj]);CHKERRQ(ierr); 102313af1a74SHong Zhang } 102413af1a74SHong Zhang } else { 102535f5172aSHong Zhang /* \sum_{j=i+1}^s a_{ji}*Lambda_{s,j} */ 102613af1a74SHong Zhang ierr = VecSet(VecsDeltaLam2[nadj*s+i],0);CHKERRQ(ierr); 102735f5172aSHong Zhang ierr = VecMAXPY(VecsSensi2Temp[nadj],s-i-1,w,&VecsDeltaLam2[nadj*s+i+1]);CHKERRQ(ierr); 102835f5172aSHong Zhang ierr = MatMultTranspose(J,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj*s+i]);CHKERRQ(ierr); 102935f5172aSHong Zhang ierr = VecScale(VecsDeltaLam2[nadj*s+i],-h);CHKERRQ(ierr); 103035f5172aSHong Zhang ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],-h,ts->vecs_guu[nadj]);CHKERRQ(ierr); 103135f5172aSHong Zhang if (ts->vecs_sensip) { 103235f5172aSHong Zhang ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],-h,ts->vecs_gup[nadj]);CHKERRQ(ierr); 103313af1a74SHong Zhang } 103435f5172aSHong Zhang } 103535f5172aSHong Zhang if (ts->vecs_sensi2p) { /* 2nd-order adjoint for parameters */ 103613af1a74SHong Zhang ierr = MatMultTranspose(ts->Jacprhs,VecsSensi2Temp[nadj],VecDeltaMu2);CHKERRQ(ierr); 103735f5172aSHong Zhang if (b[i]) { 103835f5172aSHong Zhang ierr = VecScale(VecDeltaMu2,-h*b[i]);CHKERRQ(ierr); 103935f5172aSHong Zhang ierr = VecAXPY(VecDeltaMu2,-h*b[i],ts->vecs_gpu[nadj]);CHKERRQ(ierr); 104035f5172aSHong Zhang ierr = VecAXPY(VecDeltaMu2,-h*b[i],ts->vecs_gpp[nadj]);CHKERRQ(ierr); 104113af1a74SHong Zhang } else { 104235f5172aSHong Zhang ierr = VecScale(VecDeltaMu2,-h);CHKERRQ(ierr); 104335f5172aSHong Zhang ierr = VecAXPY(VecDeltaMu2,-h,ts->vecs_gpu[nadj]);CHKERRQ(ierr); 104435f5172aSHong Zhang ierr = VecAXPY(VecDeltaMu2,-h,ts->vecs_gpp[nadj]);CHKERRQ(ierr); 104513af1a74SHong Zhang } 104613af1a74SHong Zhang ierr = VecAXPY(ts->vecs_sensi2p[nadj],1,VecDeltaMu2);CHKERRQ(ierr); /* update sensi2p for each stage */ 104713af1a74SHong Zhang } 104813af1a74SHong Zhang } 104913af1a74SHong Zhang } 10506497ce14SHong Zhang } 1051c235aa8dSHong Zhang 1052c235aa8dSHong Zhang for (j=0; j<s; j++) w[j] = 1.0; 10532e7b7f96SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { /* no need to do this for mu's */ 10542e7b7f96SHong Zhang ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecsDeltaLam[nadj*s]);CHKERRQ(ierr); 105513af1a74SHong Zhang if (ts->vecs_sensi2) { 105613af1a74SHong Zhang ierr = VecMAXPY(ts->vecs_sensi2[nadj],s,w,&VecsDeltaLam2[nadj*s]);CHKERRQ(ierr); 105713af1a74SHong Zhang } 10586497ce14SHong Zhang } 1059c235aa8dSHong Zhang rk->status = TS_STEP_COMPLETE; 1060d2daff3dSHong Zhang PetscFunctionReturn(0); 1061d2daff3dSHong Zhang } 1062d2daff3dSHong Zhang 106313af1a74SHong Zhang static PetscErrorCode TSAdjointReset_RK(TS ts) 106413af1a74SHong Zhang { 106513af1a74SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 106613af1a74SHong Zhang RKTableau tab = rk->tableau; 106713af1a74SHong Zhang PetscErrorCode ierr; 106813af1a74SHong Zhang 106913af1a74SHong Zhang PetscFunctionBegin; 107013af1a74SHong Zhang ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam);CHKERRQ(ierr); 107113af1a74SHong Zhang ierr = VecDestroyVecs(ts->numcost,&rk->VecsSensiTemp);CHKERRQ(ierr); 107213af1a74SHong Zhang ierr = VecDestroy(&rk->VecDeltaMu);CHKERRQ(ierr); 107313af1a74SHong Zhang ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam2);CHKERRQ(ierr); 107413af1a74SHong Zhang ierr = VecDestroy(&rk->VecDeltaMu2);CHKERRQ(ierr); 107513af1a74SHong Zhang ierr = VecDestroyVecs(ts->numcost,&rk->VecsSensi2Temp);CHKERRQ(ierr); 107613af1a74SHong Zhang PetscFunctionReturn(0); 107713af1a74SHong Zhang } 107813af1a74SHong Zhang 107955de54fcSHong Zhang static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X) 108055de54fcSHong Zhang { 108155de54fcSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 108255de54fcSHong Zhang PetscInt s = rk->tableau->s,p = rk->tableau->p,i,j; 108355de54fcSHong Zhang PetscReal h; 108455de54fcSHong Zhang PetscReal tt,t; 108555de54fcSHong Zhang PetscScalar *b; 108655de54fcSHong Zhang const PetscReal *B = rk->tableau->binterp; 108755de54fcSHong Zhang PetscErrorCode ierr; 108855de54fcSHong Zhang 108955de54fcSHong Zhang PetscFunctionBegin; 10902c71b3e2SJacob Faibussowitsch PetscCheckFalse(!B,PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 109155de54fcSHong Zhang 109255de54fcSHong Zhang switch (rk->status) { 109355de54fcSHong Zhang case TS_STEP_INCOMPLETE: 109455de54fcSHong Zhang case TS_STEP_PENDING: 109555de54fcSHong Zhang h = ts->time_step; 109655de54fcSHong Zhang t = (itime - ts->ptime)/h; 109755de54fcSHong Zhang break; 109855de54fcSHong Zhang case TS_STEP_COMPLETE: 109955de54fcSHong Zhang h = ts->ptime - ts->ptime_prev; 110055de54fcSHong Zhang t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 110155de54fcSHong Zhang break; 110255de54fcSHong Zhang default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 110355de54fcSHong Zhang } 110455de54fcSHong Zhang ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 110555de54fcSHong Zhang for (i=0; i<s; i++) b[i] = 0; 110655de54fcSHong Zhang for (j=0,tt=t; j<p; j++,tt*=t) { 110755de54fcSHong Zhang for (i=0; i<s; i++) { 110855de54fcSHong Zhang b[i] += h * B[i*p+j] * tt; 110955de54fcSHong Zhang } 111055de54fcSHong Zhang } 111155de54fcSHong Zhang ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr); 111255de54fcSHong Zhang ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr); 111355de54fcSHong Zhang ierr = PetscFree(b);CHKERRQ(ierr); 111455de54fcSHong Zhang PetscFunctionReturn(0); 111555de54fcSHong Zhang } 111655de54fcSHong Zhang 111755de54fcSHong Zhang /*------------------------------------------------------------*/ 111855de54fcSHong Zhang 1119be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts) 1120be5899b3SLisandro Dalcin { 1121be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 1122be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 1123be5899b3SLisandro Dalcin PetscErrorCode ierr; 1124be5899b3SLisandro Dalcin 1125be5899b3SLisandro Dalcin PetscFunctionBegin; 1126be5899b3SLisandro Dalcin if (!tab) PetscFunctionReturn(0); 1127be5899b3SLisandro Dalcin ierr = PetscFree(rk->work);CHKERRQ(ierr); 1128be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr); 1129be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr); 1130be5899b3SLisandro Dalcin PetscFunctionReturn(0); 1131be5899b3SLisandro Dalcin } 1132be5899b3SLisandro Dalcin 1133277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts) 1134e4dd521cSBarry Smith { 11356849ba73SBarry Smith PetscErrorCode ierr; 1136e4dd521cSBarry Smith 1137e4dd521cSBarry Smith PetscFunctionBegin; 1138be5899b3SLisandro Dalcin ierr = TSRKTableauReset(ts);CHKERRQ(ierr); 11390fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 114002555c94SHong Zhang ierr = PetscTryMethod(ts,"TSReset_RK_MultirateSplit_C",(TS),(ts));CHKERRQ(ierr); 11410fe4d17eSHong Zhang } else { 114202555c94SHong Zhang ierr = PetscTryMethod(ts,"TSReset_RK_MultirateNonsplit_C",(TS),(ts));CHKERRQ(ierr); 11430fe4d17eSHong Zhang } 1144277b19d0SLisandro Dalcin PetscFunctionReturn(0); 1145e4dd521cSBarry Smith } 1146277b19d0SLisandro Dalcin 1147f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx) 1148f68a32c8SEmil Constantinescu { 1149f68a32c8SEmil Constantinescu PetscFunctionBegin; 1150f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1151f68a32c8SEmil Constantinescu } 1152f68a32c8SEmil Constantinescu 1153f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 1154f68a32c8SEmil Constantinescu { 1155f68a32c8SEmil Constantinescu PetscFunctionBegin; 1156f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1157f68a32c8SEmil Constantinescu } 1158f68a32c8SEmil Constantinescu 1159f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx) 1160f68a32c8SEmil Constantinescu { 1161f68a32c8SEmil Constantinescu PetscFunctionBegin; 1162f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1163f68a32c8SEmil Constantinescu } 1164f68a32c8SEmil Constantinescu 1165f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 1166f68a32c8SEmil Constantinescu { 1167f68a32c8SEmil Constantinescu 1168f68a32c8SEmil Constantinescu PetscFunctionBegin; 1169f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1170f68a32c8SEmil Constantinescu } 1171be5899b3SLisandro Dalcin 1172be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts) 1173be5899b3SLisandro Dalcin { 1174be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 1175be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 1176be5899b3SLisandro Dalcin PetscErrorCode ierr; 1177be5899b3SLisandro Dalcin 1178be5899b3SLisandro Dalcin PetscFunctionBegin; 1179be5899b3SLisandro Dalcin ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr); 1180be5899b3SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr); 1181be5899b3SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr); 1182be5899b3SLisandro Dalcin PetscFunctionReturn(0); 1183be5899b3SLisandro Dalcin } 1184be5899b3SLisandro Dalcin 1185f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts) 1186f68a32c8SEmil Constantinescu { 1187cd4cee2dSHong Zhang TS quadts = ts->quadraturets; 1188f68a32c8SEmil Constantinescu PetscErrorCode ierr; 1189f68a32c8SEmil Constantinescu DM dm; 1190f68a32c8SEmil Constantinescu 1191f68a32c8SEmil Constantinescu PetscFunctionBegin; 11923dd83b38SBarry Smith ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr); 1193be5899b3SLisandro Dalcin ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr); 1194cd4cee2dSHong Zhang if (quadts && ts->costintegralfwd) { 1195cd4cee2dSHong Zhang Mat Jquad; 1196cd4cee2dSHong Zhang ierr = TSGetRHSJacobian(quadts,&Jquad,NULL,NULL,NULL);CHKERRQ(ierr); 11970f7a1166SHong Zhang } 1198f68a32c8SEmil Constantinescu ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 1199f68a32c8SEmil Constantinescu ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 1200f68a32c8SEmil Constantinescu ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 12010fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 120202555c94SHong Zhang ierr = PetscTryMethod(ts,"TSSetUp_RK_MultirateSplit_C",(TS),(ts));CHKERRQ(ierr); 12030fe4d17eSHong Zhang } else { 120402555c94SHong Zhang ierr = PetscTryMethod(ts,"TSSetUp_RK_MultirateNonsplit_C",(TS),(ts));CHKERRQ(ierr); 12050fe4d17eSHong Zhang } 1206e4dd521cSBarry Smith PetscFunctionReturn(0); 1207e4dd521cSBarry Smith } 1208d2daff3dSHong Zhang 12094416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts) 1210e4dd521cSBarry Smith { 1211be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 1212dfbe8321SBarry Smith PetscErrorCode ierr; 1213262ff7bbSBarry Smith 1214e4dd521cSBarry Smith PetscFunctionBegin; 1215e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr); 1216f68a32c8SEmil Constantinescu { 1217f68a32c8SEmil Constantinescu RKTableauLink link; 1218f68a32c8SEmil Constantinescu PetscInt count,choice; 12190fe4d17eSHong Zhang PetscBool flg,use_multirate = PETSC_FALSE; 1220f68a32c8SEmil Constantinescu const char **namelist; 12212c9cb062Sluzhanghpp 1222f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) ; 1223f489ac74SBarry Smith ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr); 1224f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 12250fe4d17eSHong Zhang ierr = PetscOptionsBool("-ts_rk_multirate","Use interpolation-based multirate RK method","TSRKSetMultirate",rk->use_multirate,&use_multirate,&flg);CHKERRQ(ierr); 12260fe4d17eSHong Zhang if (flg) { 12270fe4d17eSHong Zhang ierr = TSRKSetMultirate(ts,use_multirate);CHKERRQ(ierr); 12280fe4d17eSHong Zhang } 1229be5899b3SLisandro Dalcin ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr); 1230be5899b3SLisandro Dalcin if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);} 1231f68a32c8SEmil Constantinescu ierr = PetscFree(namelist);CHKERRQ(ierr); 1232f68a32c8SEmil Constantinescu } 1233262ff7bbSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 1234ea13f565SStefano Zampini ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)ts),NULL,"Multirate methods options","");CHKERRQ(ierr); 12352c9cb062Sluzhanghpp ierr = PetscOptionsInt("-ts_rk_dtratio","time step ratio between slow and fast","",rk->dtratio,&rk->dtratio,NULL);CHKERRQ(ierr); 12362c9cb062Sluzhanghpp ierr = PetscOptionsEnd();CHKERRQ(ierr); 1237e4dd521cSBarry Smith PetscFunctionReturn(0); 1238e4dd521cSBarry Smith } 1239e4dd521cSBarry Smith 12405f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer) 1241e4dd521cSBarry Smith { 12425f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 12438370ee3bSLisandro Dalcin PetscBool iascii; 1244dfbe8321SBarry Smith PetscErrorCode ierr; 12452cdf8120Sasbjorn 1246e4dd521cSBarry Smith PetscFunctionBegin; 1247251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 12488370ee3bSLisandro Dalcin if (iascii) { 12499c334d8fSLisandro Dalcin RKTableau tab = rk->tableau; 1250f68a32c8SEmil Constantinescu TSRKType rktype; 12510f7a28cdSStefano Zampini const PetscReal *c; 12520f7a28cdSStefano Zampini PetscInt s; 1253f68a32c8SEmil Constantinescu char buf[512]; 12540f7a28cdSStefano Zampini PetscBool FSAL; 12550f7a28cdSStefano Zampini 1256f68a32c8SEmil Constantinescu ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr); 12570f7a28cdSStefano Zampini ierr = TSRKGetTableau(ts,&s,NULL,NULL,&c,NULL,NULL,NULL,&FSAL);CHKERRQ(ierr); 1258efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RK type %s\n",rktype);CHKERRQ(ierr); 12590c37eb8eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," Order: %D\n",tab->order);CHKERRQ(ierr); 12600f7a28cdSStefano Zampini ierr = PetscViewerASCIIPrintf(viewer," FSAL property: %s\n",FSAL ? "yes" : "no");CHKERRQ(ierr); 12610f7a28cdSStefano Zampini ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",s,c);CHKERRQ(ierr); 1262f68a32c8SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," Abscissa c = %s\n",buf);CHKERRQ(ierr); 12638370ee3bSLisandro Dalcin } 1264f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1265f68a32c8SEmil Constantinescu } 1266f68a32c8SEmil Constantinescu 1267f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer) 1268f68a32c8SEmil Constantinescu { 1269f68a32c8SEmil Constantinescu PetscErrorCode ierr; 12709c334d8fSLisandro Dalcin TSAdapt adapt; 1271f68a32c8SEmil Constantinescu 1272f68a32c8SEmil Constantinescu PetscFunctionBegin; 12739c334d8fSLisandro Dalcin ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 12749c334d8fSLisandro Dalcin ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr); 1275f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1276f68a32c8SEmil Constantinescu } 1277f68a32c8SEmil Constantinescu 12782ea87ba9SLisandro Dalcin /*@ 12792ea87ba9SLisandro Dalcin TSRKGetOrder - Get the order of RK scheme 12802ea87ba9SLisandro Dalcin 12812ea87ba9SLisandro Dalcin Not collective 12822ea87ba9SLisandro Dalcin 12832ea87ba9SLisandro Dalcin Input Parameter: 12842ea87ba9SLisandro Dalcin . ts - timestepping context 12852ea87ba9SLisandro Dalcin 12862ea87ba9SLisandro Dalcin Output Parameter: 12872ea87ba9SLisandro Dalcin . order - order of RK-scheme 12882ea87ba9SLisandro Dalcin 12892ea87ba9SLisandro Dalcin Level: intermediate 12902ea87ba9SLisandro Dalcin 12912ea87ba9SLisandro Dalcin .seealso: TSRKGetType() 12922ea87ba9SLisandro Dalcin @*/ 12932ea87ba9SLisandro Dalcin PetscErrorCode TSRKGetOrder(TS ts,PetscInt *order) 12942ea87ba9SLisandro Dalcin { 12952ea87ba9SLisandro Dalcin PetscErrorCode ierr; 12962ea87ba9SLisandro Dalcin 12972ea87ba9SLisandro Dalcin PetscFunctionBegin; 12982ea87ba9SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 12992ea87ba9SLisandro Dalcin PetscValidIntPointer(order,2); 13002ea87ba9SLisandro Dalcin ierr = PetscUseMethod(ts,"TSRKGetOrder_C",(TS,PetscInt*),(ts,order));CHKERRQ(ierr); 13012ea87ba9SLisandro Dalcin PetscFunctionReturn(0); 13022ea87ba9SLisandro Dalcin } 13032ea87ba9SLisandro Dalcin 1304f68a32c8SEmil Constantinescu /*@C 1305f68a32c8SEmil Constantinescu TSRKSetType - Set the type of RK scheme 1306f68a32c8SEmil Constantinescu 1307f68a32c8SEmil Constantinescu Logically collective 1308f68a32c8SEmil Constantinescu 1309d8d19677SJose E. Roman Input Parameters: 1310f68a32c8SEmil Constantinescu + ts - timestepping context 1311f68a32c8SEmil Constantinescu - rktype - type of RK-scheme 1312f68a32c8SEmil Constantinescu 1313287c2655SBarry Smith Options Database: 13149bd3e852SBarry Smith . -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs> 1315287c2655SBarry Smith 1316f68a32c8SEmil Constantinescu Level: intermediate 1317f68a32c8SEmil Constantinescu 1318e7685601SHong Zhang .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK2B, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS, TSRK6VR, TSRK7VR, TSRK8VR 1319f68a32c8SEmil Constantinescu @*/ 1320f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype) 1321f68a32c8SEmil Constantinescu { 1322f68a32c8SEmil Constantinescu PetscErrorCode ierr; 1323f68a32c8SEmil Constantinescu 1324f68a32c8SEmil Constantinescu PetscFunctionBegin; 1325f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1326b92453a8SLisandro Dalcin PetscValidCharPointer(rktype,2); 1327f68a32c8SEmil Constantinescu ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr); 1328f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1329f68a32c8SEmil Constantinescu } 1330f68a32c8SEmil Constantinescu 1331f68a32c8SEmil Constantinescu /*@C 1332f68a32c8SEmil Constantinescu TSRKGetType - Get the type of RK scheme 1333f68a32c8SEmil Constantinescu 13342ea87ba9SLisandro Dalcin Not collective 1335f68a32c8SEmil Constantinescu 1336f68a32c8SEmil Constantinescu Input Parameter: 1337f68a32c8SEmil Constantinescu . ts - timestepping context 1338f68a32c8SEmil Constantinescu 1339f68a32c8SEmil Constantinescu Output Parameter: 1340f68a32c8SEmil Constantinescu . rktype - type of RK-scheme 1341f68a32c8SEmil Constantinescu 1342f68a32c8SEmil Constantinescu Level: intermediate 1343f68a32c8SEmil Constantinescu 13442ea87ba9SLisandro Dalcin .seealso: TSRKSetType() 1345f68a32c8SEmil Constantinescu @*/ 1346f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype) 1347f68a32c8SEmil Constantinescu { 1348f68a32c8SEmil Constantinescu PetscErrorCode ierr; 1349f68a32c8SEmil Constantinescu 1350f68a32c8SEmil Constantinescu PetscFunctionBegin; 1351f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1352f68a32c8SEmil Constantinescu ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr); 1353f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1354f68a32c8SEmil Constantinescu } 1355f68a32c8SEmil Constantinescu 13562ea87ba9SLisandro Dalcin static PetscErrorCode TSRKGetOrder_RK(TS ts,PetscInt *order) 13572ea87ba9SLisandro Dalcin { 13582ea87ba9SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 13592ea87ba9SLisandro Dalcin 13602ea87ba9SLisandro Dalcin PetscFunctionBegin; 13612ea87ba9SLisandro Dalcin *order = rk->tableau->order; 13622ea87ba9SLisandro Dalcin PetscFunctionReturn(0); 13632ea87ba9SLisandro Dalcin } 13642ea87ba9SLisandro Dalcin 1365560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype) 1366f68a32c8SEmil Constantinescu { 1367f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 1368f68a32c8SEmil Constantinescu 1369f68a32c8SEmil Constantinescu PetscFunctionBegin; 1370f68a32c8SEmil Constantinescu *rktype = rk->tableau->name; 1371f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1372f68a32c8SEmil Constantinescu } 13732c9cb062Sluzhanghpp 1374560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype) 1375f68a32c8SEmil Constantinescu { 1376f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 1377f68a32c8SEmil Constantinescu PetscErrorCode ierr; 1378f68a32c8SEmil Constantinescu PetscBool match; 1379f68a32c8SEmil Constantinescu RKTableauLink link; 1380f68a32c8SEmil Constantinescu 1381f68a32c8SEmil Constantinescu PetscFunctionBegin; 1382f68a32c8SEmil Constantinescu if (rk->tableau) { 1383f68a32c8SEmil Constantinescu ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr); 1384f68a32c8SEmil Constantinescu if (match) PetscFunctionReturn(0); 1385f68a32c8SEmil Constantinescu } 1386f68a32c8SEmil Constantinescu for (link = RKTableauList; link; link=link->next) { 1387f68a32c8SEmil Constantinescu ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr); 1388f68a32c8SEmil Constantinescu if (match) { 1389be5899b3SLisandro Dalcin if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);} 1390f68a32c8SEmil Constantinescu rk->tableau = &link->tab; 1391be5899b3SLisandro Dalcin if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);} 13922ffb9264SLisandro Dalcin ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 1393f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1394f68a32c8SEmil Constantinescu } 1395f68a32c8SEmil Constantinescu } 139698921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype); 1397e4dd521cSBarry Smith } 1398e4dd521cSBarry Smith 1399ff22ae23SHong Zhang static PetscErrorCode TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y) 1400ff22ae23SHong Zhang { 1401ff22ae23SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 1402ff22ae23SHong Zhang 1403ff22ae23SHong Zhang PetscFunctionBegin; 14040429704eSStefano Zampini if (ns) *ns = rk->tableau->s; 1405d2daff3dSHong Zhang if (Y) *Y = rk->Y; 1406ff22ae23SHong Zhang PetscFunctionReturn(0); 1407ff22ae23SHong Zhang } 1408ff22ae23SHong Zhang 1409b3a6b972SJed Brown static PetscErrorCode TSDestroy_RK(TS ts) 1410b3a6b972SJed Brown { 1411b3a6b972SJed Brown PetscErrorCode ierr; 1412b3a6b972SJed Brown 1413b3a6b972SJed Brown PetscFunctionBegin; 1414b3a6b972SJed Brown ierr = TSReset_RK(ts);CHKERRQ(ierr); 1415b3a6b972SJed Brown if (ts->dm) { 1416b3a6b972SJed Brown ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 1417b3a6b972SJed Brown ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 1418b3a6b972SJed Brown } 1419b3a6b972SJed Brown ierr = PetscFree(ts->data);CHKERRQ(ierr); 14202ea87ba9SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetOrder_C",NULL);CHKERRQ(ierr); 1421b3a6b972SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr); 1422b3a6b972SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr); 14230f7a28cdSStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetTableau_C",NULL);CHKERRQ(ierr); 14240fe4d17eSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",NULL);CHKERRQ(ierr); 14250fe4d17eSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",NULL);CHKERRQ(ierr); 1426b3a6b972SJed Brown PetscFunctionReturn(0); 1427b3a6b972SJed Brown } 1428ff22ae23SHong Zhang 14293ca0882eSHong Zhang /* 14303ca0882eSHong Zhang This defines the nonlinear equation that is to be solved with SNES 14313ca0882eSHong Zhang We do not need to solve the equation; we just use SNES to approximate the Jacobian 14323ca0882eSHong Zhang */ 14333ca0882eSHong Zhang static PetscErrorCode SNESTSFormFunction_RK(SNES snes,Vec x,Vec y,TS ts) 14343ca0882eSHong Zhang { 14353ca0882eSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 14363ca0882eSHong Zhang DM dm,dmsave; 14373ca0882eSHong Zhang PetscErrorCode ierr; 14383ca0882eSHong Zhang 14393ca0882eSHong Zhang PetscFunctionBegin; 14403ca0882eSHong Zhang ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 14413ca0882eSHong Zhang /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 14423ca0882eSHong Zhang dmsave = ts->dm; 14433ca0882eSHong Zhang ts->dm = dm; 14443ca0882eSHong Zhang ierr = TSComputeRHSFunction(ts,rk->stage_time,x,y);CHKERRQ(ierr); 14453ca0882eSHong Zhang ts->dm = dmsave; 14463ca0882eSHong Zhang PetscFunctionReturn(0); 14473ca0882eSHong Zhang } 14483ca0882eSHong Zhang 14493ca0882eSHong Zhang static PetscErrorCode SNESTSFormJacobian_RK(SNES snes,Vec x,Mat A,Mat B,TS ts) 14503ca0882eSHong Zhang { 14513ca0882eSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 14523ca0882eSHong Zhang DM dm,dmsave; 14533ca0882eSHong Zhang PetscErrorCode ierr; 14543ca0882eSHong Zhang 14553ca0882eSHong Zhang PetscFunctionBegin; 14563ca0882eSHong Zhang ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 14573ca0882eSHong Zhang dmsave = ts->dm; 14583ca0882eSHong Zhang ts->dm = dm; 14593ca0882eSHong Zhang ierr = TSComputeRHSJacobian(ts,rk->stage_time,x,A,B);CHKERRQ(ierr); 14603ca0882eSHong Zhang ts->dm = dmsave; 14613ca0882eSHong Zhang PetscFunctionReturn(0); 14623ca0882eSHong Zhang } 14633ca0882eSHong Zhang 146421052c0cSHong Zhang /*@C 146521052c0cSHong Zhang TSRKSetMultirate - Use the interpolation-based multirate RK method 146621052c0cSHong Zhang 146721052c0cSHong Zhang Logically collective 146821052c0cSHong Zhang 1469d8d19677SJose E. Roman Input Parameters: 147021052c0cSHong Zhang + ts - timestepping context 147121052c0cSHong Zhang - use_multirate - PETSC_TRUE enables the multirate RK method, sets the basic method to be RK2A and sets the ratio between slow stepsize and fast stepsize to be 2 147221052c0cSHong Zhang 147321052c0cSHong Zhang Options Database: 147421052c0cSHong Zhang . -ts_rk_multirate - <true,false> 147521052c0cSHong Zhang 147621052c0cSHong Zhang Notes: 147721052c0cSHong Zhang The multirate method requires interpolation. The default interpolation works for 1st- and 2nd- order RK, but not for high-order RKs except TSRK5DP which comes with the interpolation coeffcients (binterp). 147821052c0cSHong Zhang 147921052c0cSHong Zhang Level: intermediate 148021052c0cSHong Zhang 148121052c0cSHong Zhang .seealso: TSRKGetMultirate() 148221052c0cSHong Zhang @*/ 148321052c0cSHong Zhang PetscErrorCode TSRKSetMultirate(TS ts,PetscBool use_multirate) 148421052c0cSHong Zhang { 1485f06fb28fSHong Zhang PetscErrorCode ierr; 14867dbacdf2SHong Zhang 14878a4be4dbSHong Zhang PetscFunctionBegin; 1488f06fb28fSHong Zhang ierr = PetscTryMethod(ts,"TSRKSetMultirate_C",(TS,PetscBool),(ts,use_multirate));CHKERRQ(ierr); 148921052c0cSHong Zhang PetscFunctionReturn(0); 149021052c0cSHong Zhang } 149121052c0cSHong Zhang 149221052c0cSHong Zhang /*@C 149321052c0cSHong Zhang TSRKGetMultirate - Gets whether to Use the interpolation-based multirate RK method 149421052c0cSHong Zhang 149521052c0cSHong Zhang Not collective 149621052c0cSHong Zhang 149721052c0cSHong Zhang Input Parameter: 149821052c0cSHong Zhang . ts - timestepping context 149921052c0cSHong Zhang 150021052c0cSHong Zhang Output Parameter: 150121052c0cSHong Zhang . use_multirate - PETSC_TRUE if the multirate RK method is enabled, PETSC_FALSE otherwise 150221052c0cSHong Zhang 150321052c0cSHong Zhang Level: intermediate 150421052c0cSHong Zhang 150521052c0cSHong Zhang .seealso: TSRKSetMultirate() 150621052c0cSHong Zhang @*/ 150721052c0cSHong Zhang PetscErrorCode TSRKGetMultirate(TS ts,PetscBool *use_multirate) 150821052c0cSHong Zhang { 1509f06fb28fSHong Zhang PetscErrorCode ierr; 15107dbacdf2SHong Zhang 15117dbacdf2SHong Zhang PetscFunctionBegin; 1512f06fb28fSHong Zhang ierr = PetscUseMethod(ts,"TSRKGetMultirate_C",(TS,PetscBool*),(ts,use_multirate));CHKERRQ(ierr); 151321052c0cSHong Zhang PetscFunctionReturn(0); 151421052c0cSHong Zhang } 151521052c0cSHong Zhang 1516ebe3b25bSBarry Smith /*MC 1517f68a32c8SEmil Constantinescu TSRK - ODE and DAE solver using Runge-Kutta schemes 1518ebe3b25bSBarry Smith 1519f68a32c8SEmil Constantinescu The user should provide the right hand side of the equation 1520f68a32c8SEmil Constantinescu using TSSetRHSFunction(). 1521ebe3b25bSBarry Smith 1522f68a32c8SEmil Constantinescu Notes: 152398164e13SLisandro Dalcin The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type 1524ebe3b25bSBarry Smith 1525d41222bbSBarry Smith Level: beginner 1526d41222bbSBarry Smith 15275d1808f1SStefano Zampini .seealso: TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRK2D, TTSRK2E, TSRK3, 15280fe4d17eSHong Zhang TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister(), TSRKSetMultirate(), TSRKGetMultirate() 1529ebe3b25bSBarry Smith 1530ebe3b25bSBarry Smith M*/ 15318cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts) 1532e4dd521cSBarry Smith { 15331566a47fSLisandro Dalcin TS_RK *rk; 1534dfbe8321SBarry Smith PetscErrorCode ierr; 1535e4dd521cSBarry Smith 1536e4dd521cSBarry Smith PetscFunctionBegin; 1537f68a32c8SEmil Constantinescu ierr = TSRKInitializePackage();CHKERRQ(ierr); 1538f68a32c8SEmil Constantinescu 1539f68a32c8SEmil Constantinescu ts->ops->reset = TSReset_RK; 15405f70b478SJed Brown ts->ops->destroy = TSDestroy_RK; 15415f70b478SJed Brown ts->ops->view = TSView_RK; 1542f68a32c8SEmil Constantinescu ts->ops->load = TSLoad_RK; 1543f68a32c8SEmil Constantinescu ts->ops->setup = TSSetUp_RK; 1544f68a32c8SEmil Constantinescu ts->ops->interpolate = TSInterpolate_RK; 15452c9cb062Sluzhanghpp ts->ops->step = TSStep_RK; 1546f68a32c8SEmil Constantinescu ts->ops->evaluatestep = TSEvaluateStep_RK; 1547fecfb714SLisandro Dalcin ts->ops->rollback = TSRollBack_RK; 1548f68a32c8SEmil Constantinescu ts->ops->setfromoptions = TSSetFromOptions_RK; 1549ff22ae23SHong Zhang ts->ops->getstages = TSGetStages_RK; 1550e4dd521cSBarry Smith 15513ca0882eSHong Zhang ts->ops->snesfunction = SNESTSFormFunction_RK; 15523ca0882eSHong Zhang ts->ops->snesjacobian = SNESTSFormJacobian_RK; 15530f7a1166SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_RK; 155413af1a74SHong Zhang ts->ops->adjointsetup = TSAdjointSetUp_RK; 155513af1a74SHong Zhang ts->ops->adjointstep = TSAdjointStep_RK; 155613af1a74SHong Zhang ts->ops->adjointreset = TSAdjointReset_RK; 15570f7a1166SHong Zhang 155813af1a74SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_RK; 1559922a638cSHong Zhang ts->ops->forwardsetup = TSForwardSetUp_RK; 1560922a638cSHong Zhang ts->ops->forwardreset = TSForwardReset_RK; 1561922a638cSHong Zhang ts->ops->forwardstep = TSForwardStep_RK; 1562922a638cSHong Zhang ts->ops->forwardgetstages= TSForwardGetStages_RK; 1563922a638cSHong Zhang 15641566a47fSLisandro Dalcin ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr); 15651566a47fSLisandro Dalcin ts->data = (void*)rk; 1566e4dd521cSBarry Smith 15672ea87ba9SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetOrder_C",TSRKGetOrder_RK);CHKERRQ(ierr); 1568f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr); 1569f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr); 15700f7a28cdSStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetTableau_C",TSRKGetTableau_RK);CHKERRQ(ierr); 157121052c0cSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",TSRKSetMultirate_RK);CHKERRQ(ierr); 157221052c0cSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",TSRKGetMultirate_RK);CHKERRQ(ierr); 1573be5899b3SLisandro Dalcin 1574be5899b3SLisandro Dalcin ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr); 157590b67129SHong Zhang rk->dtratio = 1; 1576e4dd521cSBarry Smith PetscFunctionReturn(0); 1577e4dd521cSBarry Smith } 1578