1e4dd521cSBarry Smith /* 22e698f8bSDebojyoti Ghosh Code for time stepping with the Runge-Kutta method 3f68a32c8SEmil Constantinescu 4f68a32c8SEmil Constantinescu Notes: 5f68a32c8SEmil Constantinescu The general system is written as 6f68a32c8SEmil Constantinescu 72e698f8bSDebojyoti Ghosh Udot = F(t,U) 8f68a32c8SEmil Constantinescu 9e4dd521cSBarry Smith */ 10af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 11f68a32c8SEmil Constantinescu #include <petscdm.h> 12f68a32c8SEmil Constantinescu 13484bcdc7SDebojyoti Ghosh static TSRKType TSRKDefault = TSRK3BS; 14f68a32c8SEmil Constantinescu static PetscBool TSRKRegisterAllCalled; 15f68a32c8SEmil Constantinescu static PetscBool TSRKPackageInitialized; 16f68a32c8SEmil Constantinescu 17f68a32c8SEmil Constantinescu typedef struct _RKTableau *RKTableau; 18f68a32c8SEmil Constantinescu struct _RKTableau { 19f68a32c8SEmil Constantinescu char *name; 20d760c35bSDebojyoti Ghosh PetscInt order; /* Classical approximation order of the method i */ 21f68a32c8SEmil Constantinescu PetscInt s; /* Number of stages */ 223a8a9803SLisandro Dalcin PetscInt p; /* Interpolation order */ 23d760c35bSDebojyoti Ghosh PetscBool FSAL; /* flag to indicate if tableau is FSAL */ 24f68a32c8SEmil Constantinescu PetscReal *A,*b,*c; /* Tableau */ 25f68a32c8SEmil Constantinescu PetscReal *bembed; /* Embedded formula of order one less (order-1) */ 26f68a32c8SEmil Constantinescu PetscReal *binterp; /* Dense output formula */ 27f68a32c8SEmil Constantinescu PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */ 28f68a32c8SEmil Constantinescu }; 29f68a32c8SEmil Constantinescu typedef struct _RKTableauLink *RKTableauLink; 30f68a32c8SEmil Constantinescu struct _RKTableauLink { 31f68a32c8SEmil Constantinescu struct _RKTableau tab; 32f68a32c8SEmil Constantinescu RKTableauLink next; 33f68a32c8SEmil Constantinescu }; 34f68a32c8SEmil Constantinescu static RKTableauLink RKTableauList; 35e4dd521cSBarry Smith 36e4dd521cSBarry Smith typedef struct { 37f68a32c8SEmil Constantinescu RKTableau tableau; 38f68a32c8SEmil Constantinescu Vec *Y; /* States computed during the step */ 39f68a32c8SEmil Constantinescu Vec *YdotRHS; /* Function evaluations for the non-stiff part */ 40ad8e2604SHong Zhang Vec *VecDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage */ 41ad8e2604SHong Zhang Vec *VecDeltaMu; /* Increment of the adjoint sensitivity w.r.t P at stage */ 42b18ea86cSHong Zhang Vec *VecSensiTemp; /* Vector to be timed with Jacobian transpose */ 430f7a1166SHong Zhang Vec VecCostIntegral0; /* backup for roll-backs due to events */ 44f68a32c8SEmil Constantinescu PetscScalar *work; /* Scalar work */ 45f68a32c8SEmil Constantinescu PetscReal stage_time; 46f68a32c8SEmil Constantinescu TSStepStatus status; 470f7a1166SHong Zhang PetscReal ptime; 480f7a1166SHong Zhang PetscReal time_step; 495f70b478SJed Brown } TS_RK; 50e4dd521cSBarry Smith 51f68a32c8SEmil Constantinescu /*MC 52f68a32c8SEmil Constantinescu TSRK1 - First order forward Euler scheme. 53262ff7bbSBarry Smith 54f68a32c8SEmil Constantinescu This method has one stage. 55f68a32c8SEmil Constantinescu 56f68a32c8SEmil Constantinescu Level: advanced 57f68a32c8SEmil Constantinescu 58f68a32c8SEmil Constantinescu .seealso: TSRK 59f68a32c8SEmil Constantinescu M*/ 60f68a32c8SEmil Constantinescu /*MC 612109b73fSDebojyoti Ghosh TSRK2A - Second order RK scheme. 62f68a32c8SEmil Constantinescu 63f68a32c8SEmil Constantinescu This method has two stages. 64f68a32c8SEmil Constantinescu 65f68a32c8SEmil Constantinescu Level: advanced 66f68a32c8SEmil Constantinescu 67f68a32c8SEmil Constantinescu .seealso: TSRK 68f68a32c8SEmil Constantinescu M*/ 69f68a32c8SEmil Constantinescu /*MC 70f68a32c8SEmil Constantinescu TSRK3 - Third order RK scheme. 71f68a32c8SEmil Constantinescu 72f68a32c8SEmil Constantinescu This method has three stages. 73f68a32c8SEmil Constantinescu 74f68a32c8SEmil Constantinescu Level: advanced 75f68a32c8SEmil Constantinescu 76f68a32c8SEmil Constantinescu .seealso: TSRK 77f68a32c8SEmil Constantinescu M*/ 78f68a32c8SEmil Constantinescu /*MC 792109b73fSDebojyoti Ghosh TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method. 802109b73fSDebojyoti Ghosh 81d1905564SLisandro Dalcin This method has four stages with the First Same As Last (FSAL) property. 822109b73fSDebojyoti Ghosh 832109b73fSDebojyoti Ghosh Level: advanced 842109b73fSDebojyoti Ghosh 8598164e13SLisandro Dalcin References: https://doi.org/10.1016/0893-9659(89)90079-7 86d1905564SLisandro Dalcin 872109b73fSDebojyoti Ghosh .seealso: TSRK 882109b73fSDebojyoti Ghosh M*/ 892109b73fSDebojyoti Ghosh /*MC 90f68a32c8SEmil Constantinescu TSRK4 - Fourth order RK scheme. 91f68a32c8SEmil Constantinescu 922e698f8bSDebojyoti Ghosh This is the classical Runge-Kutta method with four stages. 93f68a32c8SEmil Constantinescu 94f68a32c8SEmil Constantinescu Level: advanced 95f68a32c8SEmil Constantinescu 96f68a32c8SEmil Constantinescu .seealso: TSRK 97f68a32c8SEmil Constantinescu M*/ 98f68a32c8SEmil Constantinescu /*MC 992e698f8bSDebojyoti Ghosh TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method. 100f68a32c8SEmil Constantinescu 101f68a32c8SEmil Constantinescu This method has six stages. 102f68a32c8SEmil Constantinescu 103f68a32c8SEmil Constantinescu Level: advanced 104f68a32c8SEmil Constantinescu 105f68a32c8SEmil Constantinescu .seealso: TSRK 106f68a32c8SEmil Constantinescu M*/ 1072109b73fSDebojyoti Ghosh /*MC 1082e698f8bSDebojyoti Ghosh TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method. 1092109b73fSDebojyoti Ghosh 110d1905564SLisandro Dalcin This method has seven stages with the First Same As Last (FSAL) property. 1112109b73fSDebojyoti Ghosh 1122109b73fSDebojyoti Ghosh Level: advanced 1132109b73fSDebojyoti Ghosh 11498164e13SLisandro Dalcin References: https://doi.org/10.1016/0771-050X(80)90013-3 115d1905564SLisandro Dalcin 1162109b73fSDebojyoti Ghosh .seealso: TSRK 1172109b73fSDebojyoti Ghosh M*/ 11805e23783SLisandro Dalcin /*MC 11905e23783SLisandro Dalcin TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method. 12005e23783SLisandro Dalcin 12105e23783SLisandro Dalcin This method has eight stages with the First Same As Last (FSAL) property. 12205e23783SLisandro Dalcin 12305e23783SLisandro Dalcin Level: advanced 12405e23783SLisandro Dalcin 12598164e13SLisandro Dalcin References: https://doi.org/10.1016/0898-1221(96)00141-1 12605e23783SLisandro Dalcin 12705e23783SLisandro Dalcin .seealso: TSRK 12805e23783SLisandro Dalcin M*/ 129262ff7bbSBarry Smith 130f68a32c8SEmil Constantinescu /*@C 131f68a32c8SEmil Constantinescu TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK 132262ff7bbSBarry Smith 133f68a32c8SEmil Constantinescu Not Collective, but should be called by all processes which will need the schemes to be registered 134262ff7bbSBarry Smith 135f68a32c8SEmil Constantinescu Level: advanced 136262ff7bbSBarry Smith 137f68a32c8SEmil Constantinescu .keywords: TS, TSRK, register, all 138262ff7bbSBarry Smith 139f68a32c8SEmil Constantinescu .seealso: TSRKRegisterDestroy() 140262ff7bbSBarry Smith @*/ 141f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void) 142262ff7bbSBarry Smith { 1434ac538c5SBarry Smith PetscErrorCode ierr; 144262ff7bbSBarry Smith 145262ff7bbSBarry Smith PetscFunctionBegin; 146f68a32c8SEmil Constantinescu if (TSRKRegisterAllCalled) PetscFunctionReturn(0); 147f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_TRUE; 148e4dd521cSBarry Smith 149e4dd521cSBarry Smith { 150f68a32c8SEmil Constantinescu const PetscReal 151f68a32c8SEmil Constantinescu A[1][1] = {{0.0}}, 152f68a32c8SEmil Constantinescu b[1] = {1.0}; 1533a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 154e4dd521cSBarry Smith } 155db046809SBarry Smith { 156f68a32c8SEmil Constantinescu const PetscReal 157f68a32c8SEmil Constantinescu A[2][2] = {{0.0,0.0}, 158f68a32c8SEmil Constantinescu {1.0,0.0}}, 159f68a32c8SEmil Constantinescu b[2] = {0.5,0.5}, 160f68a32c8SEmil Constantinescu bembed[2] = {1.0,0}; 1613a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 162db046809SBarry Smith } 163f68a32c8SEmil Constantinescu { 164f68a32c8SEmil Constantinescu const PetscReal 165*17f6384fSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) 166*17f6384fSBarry Smith A[3][3] = {{0,0,0}, 167*17f6384fSBarry Smith {2.0q/3.0q,0,0}, 168*17f6384fSBarry Smith {-1.0q/3.0q,1.0,0}}, 169*17f6384fSBarry Smith #else 170f68a32c8SEmil Constantinescu A[3][3] = {{0,0,0}, 171f68a32c8SEmil Constantinescu {2.0/3.0,0,0}, 172f68a32c8SEmil Constantinescu {-1.0/3.0,1.0,0}}, 173*17f6384fSBarry Smith #endif 174f68a32c8SEmil Constantinescu b[3] = {0.25,0.5,0.25}; 1753a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 176fdefd5e5SDebojyoti Ghosh } 177fdefd5e5SDebojyoti Ghosh { 178fdefd5e5SDebojyoti Ghosh const PetscReal 179*17f6384fSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) 180*17f6384fSBarry Smith A[4][4] = {{0,0,0,0}, 181*17f6384fSBarry Smith {1.0q/2.0q,0,0,0}, 182*17f6384fSBarry Smith {0,3.0q/4.0q,0,0}, 183*17f6384fSBarry Smith {2.0q/9.0q,1.0q/3.0q,4.0q/9.0q,0}}, 184*17f6384fSBarry Smith b[4] = {2.0q/9.0q,1.0q/3.0q,4.0q/9.0q,0}, 185*17f6384fSBarry Smith bembed[4] = {7.0q/24.0q,1.0q/4.0q,1.0q/3.0q,1.0q/8.0q}; 186*17f6384fSBarry Smith #else 187fdefd5e5SDebojyoti Ghosh A[4][4] = {{0,0,0,0}, 188fdefd5e5SDebojyoti Ghosh {1.0/2.0,0,0,0}, 189fdefd5e5SDebojyoti Ghosh {0,3.0/4.0,0,0}, 190fdefd5e5SDebojyoti Ghosh {2.0/9.0,1.0/3.0,4.0/9.0,0}}, 191fdefd5e5SDebojyoti Ghosh b[4] = {2.0/9.0,1.0/3.0,4.0/9.0,0}, 192fdefd5e5SDebojyoti Ghosh bembed[4] = {7.0/24.0,1.0/4.0,1.0/3.0,1.0/8.0}; 193*17f6384fSBarry Smith #endif 1943a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 195db046809SBarry Smith } 196f68a32c8SEmil Constantinescu { 197f68a32c8SEmil Constantinescu const PetscReal 198f68a32c8SEmil Constantinescu A[4][4] = {{0,0,0,0}, 199f68a32c8SEmil Constantinescu {0.5,0,0,0}, 200f68a32c8SEmil Constantinescu {0,0.5,0,0}, 201f68a32c8SEmil Constantinescu {0,0,1.0,0}}, 202*17f6384fSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) 203*17f6384fSBarry Smith b[4] = {1.0q/6.0q,1.0q/3.0q,1.0q/3.0q,1.0q/6.0q}; 204*17f6384fSBarry Smith #else 205f68a32c8SEmil Constantinescu b[4] = {1.0/6.0,1.0/3.0,1.0/3.0,1.0/6.0}; 206*17f6384fSBarry Smith #endif 2073a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 208db046809SBarry Smith } 209f68a32c8SEmil Constantinescu { 210f68a32c8SEmil Constantinescu const PetscReal 211*17f6384fSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) 212*17f6384fSBarry Smith A[6][6] = {{0,0,0,0,0,0}, 213*17f6384fSBarry Smith {0.25q,0,0,0,0,0}, 214*17f6384fSBarry Smith {3.0q/32.0q,9.0q/32.0q,0,0,0,0}, 215*17f6384fSBarry Smith {1932.0q/2197.0q,-7200.0q/2197.0q,7296.0q/2197.0q,0,0,0}, 216*17f6384fSBarry Smith {439.0q/216.0q,-8.0q,3680.0q/513.0q,-845.0q/4104.0q,0,0}, 217*17f6384fSBarry Smith {-8.0q/27.0q,2.0q,-3544.0q/2565.0q,1859.0q/4104.0q,-11.0q/40.0q,0}}, 218*17f6384fSBarry Smith b[6] = {16.0q/135.0q,0,6656.0q/12825.0q,28561.0q/56430.0q,-9.0q/50.0q,2.0q/55.0q}, 219*17f6384fSBarry Smith bembed[6] = {25.0q/216.0q,0,1408.0q/2565.0q,2197.0q/4104.0q,-1.0q/5.0q,0}; 220*17f6384fSBarry Smith #else 221f68a32c8SEmil Constantinescu A[6][6] = {{0,0,0,0,0,0}, 222f68a32c8SEmil Constantinescu {0.25,0,0,0,0,0}, 223f68a32c8SEmil Constantinescu {3.0/32.0,9.0/32.0,0,0,0,0}, 224f68a32c8SEmil Constantinescu {1932.0/2197.0,-7200.0/2197.0,7296.0/2197.0,0,0,0}, 225f68a32c8SEmil Constantinescu {439.0/216.0,-8.0,3680.0/513.0,-845.0/4104.0,0,0}, 226f68a32c8SEmil Constantinescu {-8.0/27.0,2.0,-3544.0/2565.0,1859.0/4104.0,-11.0/40.0,0}}, 227f68a32c8SEmil Constantinescu b[6] = {16.0/135.0,0,6656.0/12825.0,28561.0/56430.0,-9.0/50.0,2.0/55.0}, 228f68a32c8SEmil Constantinescu bembed[6] = {25.0/216.0,0,1408.0/2565.0,2197.0/4104.0,-1.0/5.0,0}; 229*17f6384fSBarry Smith #endif 2303a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 231fdefd5e5SDebojyoti Ghosh } 232fdefd5e5SDebojyoti Ghosh { 233fdefd5e5SDebojyoti Ghosh const PetscReal 234*17f6384fSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) 235*17f6384fSBarry Smith A[7][7] = {{0,0,0,0,0,0,0}, 236*17f6384fSBarry Smith {1.0q/5.0q,0,0,0,0,0,0}, 237*17f6384fSBarry Smith {3.0q/40.0q,9.0q/40.0q,0,0,0,0,0}, 238*17f6384fSBarry Smith {44.0q/45.0q,-56.0q/15.0q,32.0q/9.0q,0,0,0,0}, 239*17f6384fSBarry Smith {19372.0q/6561.0q,-25360.0q/2187.0q,64448.0q/6561.0q,-212.0q/729.0q,0,0,0}, 240*17f6384fSBarry Smith {9017.0q/3168.0q,-355.0q/33.0q,46732.0q/5247.0q,49.0q/176.0q,-5103.0q/18656.0q,0,0}, 241*17f6384fSBarry Smith {35.0q/384.0q,0,500.0q/1113.0q,125.0q/192.0q,-2187.0q/6784.0q,11.0q/84.0q,0}}, 242*17f6384fSBarry Smith b[7] = {35.0q/384.0q,0,500.0q/1113.0q,125.0q/192.0q,-2187.0q/6784.0q,11.0q/84.0q,0}, 243*17f6384fSBarry Smith bembed[7] = {5179.0q/57600.0q,0,7571.0q/16695.0q,393.0q/640.0q,-92097.0q/339200.0q,187.0q/2100.0q,1.0q/40.0q}; 244*17f6384fSBarry Smith #else 245fdefd5e5SDebojyoti Ghosh A[7][7] = {{0,0,0,0,0,0,0}, 246fdefd5e5SDebojyoti Ghosh {1.0/5.0,0,0,0,0,0,0}, 247fdefd5e5SDebojyoti Ghosh {3.0/40.0,9.0/40.0,0,0,0,0,0}, 248fdefd5e5SDebojyoti Ghosh {44.0/45.0,-56.0/15.0,32.0/9.0,0,0,0,0}, 249fdefd5e5SDebojyoti Ghosh {19372.0/6561.0,-25360.0/2187.0,64448.0/6561.0,-212.0/729.0,0,0,0}, 250fdefd5e5SDebojyoti Ghosh {9017.0/3168.0,-355.0/33.0,46732.0/5247.0,49.0/176.0,-5103.0/18656.0,0,0}, 251fdefd5e5SDebojyoti Ghosh {35.0/384.0,0,500.0/1113.0,125.0/192.0,-2187.0/6784.0,11.0/84.0,0}}, 252fdefd5e5SDebojyoti Ghosh b[7] = {35.0/384.0,0,500.0/1113.0,125.0/192.0,-2187.0/6784.0,11.0/84.0,0}, 253fdefd5e5SDebojyoti Ghosh bembed[7] = {5179.0/57600.0,0,7571.0/16695.0,393.0/640.0,-92097.0/339200.0,187.0/2100.0,1.0/40.0}; 254*17f6384fSBarry Smith #endif 2553a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 256f68a32c8SEmil Constantinescu } 25705e23783SLisandro Dalcin { 25805e23783SLisandro Dalcin const PetscReal 259*17f6384fSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) 260*17f6384fSBarry Smith A[8][8] = {{0,0,0,0,0,0,0,0}, 261*17f6384fSBarry Smith {1.0q/6.0q,0,0,0,0,0,0,0}, 262*17f6384fSBarry Smith {2.0q/27.0q,4.0q/27.0q,0,0,0,0,0,0}, 263*17f6384fSBarry Smith {183.0q/1372.0q,-162.0q/343.0q,1053.0q/1372.0q,0,0,0,0,0}, 264*17f6384fSBarry Smith {68.0q/297.0q,-4.0q/11.0q,42.0q/143.0q,1960.0q/3861.0q,0,0,0,0}, 265*17f6384fSBarry Smith {597.0q/22528.0q,81.0q/352.0q,63099.0q/585728.0q,58653.0q/366080.0q,4617.0q/20480.0q,0,0,0}, 266*17f6384fSBarry Smith {174197.0q/959244.0q,-30942.0q/79937.0q,8152137.0q/19744439.0q,666106.0q/1039181.0q,-29421.0q/29068.0q,482048.0q/414219.0q,0,0}, 267*17f6384fSBarry Smith {587.0q/8064.0q,0,4440339.0q/15491840.0q,24353.0q/124800.0q,387.0q/44800.0q,2152.0q/5985.0q,7267.0q/94080.0q,0}}, 268*17f6384fSBarry Smith b[8] = {587.0q/8064.0q,0,4440339.0q/15491840.0q,24353.0q/124800.0q,387.0q/44800.0q,2152.0q/5985.0q,7267.0q/94080.0q,0}, 269*17f6384fSBarry Smith bembed[8] = {2479.0q/34992.0q,0,123.0q/416.0q,612941.0q/3411720.0q,43.0q/1440.0q,2272.0q/6561.0q,79937.0q/1113912.0q,3293.0q/556956.0q}; 270*17f6384fSBarry Smith #else 27105e23783SLisandro Dalcin A[8][8] = {{0,0,0,0,0,0,0,0}, 27205e23783SLisandro Dalcin {1.0/6.0,0,0,0,0,0,0,0}, 27305e23783SLisandro Dalcin {2.0/27.0,4.0/27.0,0,0,0,0,0,0}, 27405e23783SLisandro Dalcin {183.0/1372.0,-162.0/343.0,1053.0/1372.0,0,0,0,0,0}, 27505e23783SLisandro Dalcin {68.0/297.0,-4.0/11.0,42.0/143.0,1960.0/3861.0,0,0,0,0}, 27605e23783SLisandro Dalcin {597.0/22528.0,81.0/352.0,63099.0/585728.0,58653.0/366080.0,4617.0/20480.0,0,0,0}, 27705e23783SLisandro Dalcin {174197.0/959244.0,-30942.0/79937.0,8152137.0/19744439.0,666106.0/1039181.0,-29421.0/29068.0,482048.0/414219.0,0,0}, 27805e23783SLisandro Dalcin {587.0/8064.0,0,4440339.0/15491840.0,24353.0/124800.0,387.0/44800.0,2152.0/5985.0,7267.0/94080.0,0}}, 27905e23783SLisandro Dalcin b[8] = {587.0/8064.0,0,4440339.0/15491840.0,24353.0/124800.0,387.0/44800.0,2152.0/5985.0,7267.0/94080.0,0}, 28005e23783SLisandro Dalcin bembed[8] = {2479.0/34992.0,0,123.0/416.0,612941.0/3411720.0,43.0/1440.0,2272.0/6561.0,79937.0/1113912.0,3293.0/556956.0}; 281*17f6384fSBarry Smith #endif 28205e23783SLisandro Dalcin ierr = TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 28305e23783SLisandro Dalcin } 284db046809SBarry Smith PetscFunctionReturn(0); 285db046809SBarry Smith } 286db046809SBarry Smith 287f68a32c8SEmil Constantinescu /*@C 288f68a32c8SEmil Constantinescu TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister(). 289f68a32c8SEmil Constantinescu 290f68a32c8SEmil Constantinescu Not Collective 291f68a32c8SEmil Constantinescu 292f68a32c8SEmil Constantinescu Level: advanced 293f68a32c8SEmil Constantinescu 294f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy 295f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll() 296f68a32c8SEmil Constantinescu @*/ 297f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void) 298e4dd521cSBarry Smith { 299f68a32c8SEmil Constantinescu PetscErrorCode ierr; 300f68a32c8SEmil Constantinescu RKTableauLink link; 301f68a32c8SEmil Constantinescu 302f68a32c8SEmil Constantinescu PetscFunctionBegin; 303f68a32c8SEmil Constantinescu while ((link = RKTableauList)) { 304f68a32c8SEmil Constantinescu RKTableau t = &link->tab; 305f68a32c8SEmil Constantinescu RKTableauList = link->next; 306f68a32c8SEmil Constantinescu ierr = PetscFree3(t->A,t->b,t->c); CHKERRQ(ierr); 307f68a32c8SEmil Constantinescu ierr = PetscFree (t->bembed); CHKERRQ(ierr); 308f68a32c8SEmil Constantinescu ierr = PetscFree (t->binterp); CHKERRQ(ierr); 309f68a32c8SEmil Constantinescu ierr = PetscFree (t->name); CHKERRQ(ierr); 310f68a32c8SEmil Constantinescu ierr = PetscFree (link); CHKERRQ(ierr); 311f68a32c8SEmil Constantinescu } 312f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_FALSE; 313f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 314f68a32c8SEmil Constantinescu } 315f68a32c8SEmil Constantinescu 316f68a32c8SEmil Constantinescu /*@C 317f68a32c8SEmil Constantinescu TSRKInitializePackage - This function initializes everything in the TSRK package. It is called 318f68a32c8SEmil Constantinescu from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RK() 319f68a32c8SEmil Constantinescu when using static libraries. 320f68a32c8SEmil Constantinescu 321f68a32c8SEmil Constantinescu Level: developer 322f68a32c8SEmil Constantinescu 323f68a32c8SEmil Constantinescu .keywords: TS, TSRK, initialize, package 324f68a32c8SEmil Constantinescu .seealso: PetscInitialize() 325f68a32c8SEmil Constantinescu @*/ 326f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void) 327f68a32c8SEmil Constantinescu { 328186e87acSLisandro Dalcin PetscErrorCode ierr; 329e4dd521cSBarry Smith 330e4dd521cSBarry Smith PetscFunctionBegin; 331f68a32c8SEmil Constantinescu if (TSRKPackageInitialized) PetscFunctionReturn(0); 332f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_TRUE; 333f68a32c8SEmil Constantinescu ierr = TSRKRegisterAll();CHKERRQ(ierr); 334f68a32c8SEmil Constantinescu ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr); 335f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 336f68a32c8SEmil Constantinescu } 337186e87acSLisandro Dalcin 338f68a32c8SEmil Constantinescu /*@C 339f68a32c8SEmil Constantinescu TSRKFinalizePackage - This function destroys everything in the TSRK package. It is 340f68a32c8SEmil Constantinescu called from PetscFinalize(). 341186e87acSLisandro Dalcin 342f68a32c8SEmil Constantinescu Level: developer 343f68a32c8SEmil Constantinescu 344f68a32c8SEmil Constantinescu .keywords: Petsc, destroy, package 345f68a32c8SEmil Constantinescu .seealso: PetscFinalize() 346f68a32c8SEmil Constantinescu @*/ 347f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void) 348f68a32c8SEmil Constantinescu { 349f68a32c8SEmil Constantinescu PetscErrorCode ierr; 350f68a32c8SEmil Constantinescu 351f68a32c8SEmil Constantinescu PetscFunctionBegin; 352f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_FALSE; 353f68a32c8SEmil Constantinescu ierr = TSRKRegisterDestroy();CHKERRQ(ierr); 354f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 355f68a32c8SEmil Constantinescu } 356f68a32c8SEmil Constantinescu 357f68a32c8SEmil Constantinescu /*@C 358f68a32c8SEmil Constantinescu TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 359f68a32c8SEmil Constantinescu 360f68a32c8SEmil Constantinescu Not Collective, but the same schemes should be registered on all processes on which they will be used 361f68a32c8SEmil Constantinescu 362f68a32c8SEmil Constantinescu Input Parameters: 363f68a32c8SEmil Constantinescu + name - identifier for method 364f68a32c8SEmil Constantinescu . order - approximation order of method 365f68a32c8SEmil Constantinescu . s - number of stages, this is the dimension of the matrices below 366f68a32c8SEmil Constantinescu . A - stage coefficients (dimension s*s, row-major) 367f68a32c8SEmil Constantinescu . b - step completion table (dimension s; NULL to use last row of A) 368f68a32c8SEmil Constantinescu . c - abscissa (dimension s; NULL to use row sums of A) 369f68a32c8SEmil Constantinescu . bembed - completion table for embedded method (dimension s; NULL if not available) 3703a8a9803SLisandro Dalcin . p - Order of the interpolation scheme, equal to the number of columns of binterp 3713a8a9803SLisandro Dalcin - binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1) 372f68a32c8SEmil Constantinescu 373f68a32c8SEmil Constantinescu Notes: 374f68a32c8SEmil Constantinescu Several RK methods are provided, this function is only needed to create new methods. 375f68a32c8SEmil Constantinescu 376f68a32c8SEmil Constantinescu Level: advanced 377f68a32c8SEmil Constantinescu 378f68a32c8SEmil Constantinescu .keywords: TS, register 379f68a32c8SEmil Constantinescu 380f68a32c8SEmil Constantinescu .seealso: TSRK 381f68a32c8SEmil Constantinescu @*/ 382f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s, 383f68a32c8SEmil Constantinescu const PetscReal A[],const PetscReal b[],const PetscReal c[], 3843a8a9803SLisandro Dalcin const PetscReal bembed[],PetscInt p,const PetscReal binterp[]) 385f68a32c8SEmil Constantinescu { 386f68a32c8SEmil Constantinescu PetscErrorCode ierr; 387f68a32c8SEmil Constantinescu RKTableauLink link; 388f68a32c8SEmil Constantinescu RKTableau t; 389f68a32c8SEmil Constantinescu PetscInt i,j; 390f68a32c8SEmil Constantinescu 391f68a32c8SEmil Constantinescu PetscFunctionBegin; 3923a8a9803SLisandro Dalcin PetscValidCharPointer(name,1); 3933a8a9803SLisandro Dalcin PetscValidRealPointer(A,4); 3943a8a9803SLisandro Dalcin if (b) PetscValidRealPointer(b,5); 3953a8a9803SLisandro Dalcin if (c) PetscValidRealPointer(c,6); 3963a8a9803SLisandro Dalcin if (bembed) PetscValidRealPointer(bembed,7); 3973a8a9803SLisandro Dalcin if (binterp || p > 1) PetscValidRealPointer(binterp,9); 3983a8a9803SLisandro Dalcin 39995dccacaSBarry Smith ierr = PetscNew(&link);CHKERRQ(ierr); 400f68a32c8SEmil Constantinescu t = &link->tab; 4013a8a9803SLisandro Dalcin 402f68a32c8SEmil Constantinescu ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 403f68a32c8SEmil Constantinescu t->order = order; 404f68a32c8SEmil Constantinescu t->s = s; 405dcca6d9dSJed Brown ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr); 406f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 407f68a32c8SEmil Constantinescu if (b) { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); } 408f68a32c8SEmil Constantinescu else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i]; 409f68a32c8SEmil Constantinescu if (c) { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); } 410f68a32c8SEmil 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]; 411d760c35bSDebojyoti Ghosh t->FSAL = PETSC_TRUE; 412d760c35bSDebojyoti Ghosh for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE; 4133a8a9803SLisandro Dalcin 414f68a32c8SEmil Constantinescu if (bembed) { 415785e854fSJed Brown ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr); 416f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr); 417f68a32c8SEmil Constantinescu } 418f68a32c8SEmil Constantinescu 4193a8a9803SLisandro Dalcin if (!binterp) { p = 1; binterp = t->b; } 4203a8a9803SLisandro Dalcin t->p = p; 4213a8a9803SLisandro Dalcin ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr); 4223a8a9803SLisandro Dalcin ierr = PetscMemcpy(t->binterp,binterp,s*p*sizeof(binterp[0]));CHKERRQ(ierr); 4233a8a9803SLisandro Dalcin 424f68a32c8SEmil Constantinescu link->next = RKTableauList; 425f68a32c8SEmil Constantinescu RKTableauList = link; 426f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 427f68a32c8SEmil Constantinescu } 428f68a32c8SEmil Constantinescu 429e4dd521cSBarry Smith /* 430f68a32c8SEmil Constantinescu The step completion formula is 431e4dd521cSBarry Smith 432f68a32c8SEmil Constantinescu x1 = x0 + h b^T YdotRHS 433f68a32c8SEmil Constantinescu 434f68a32c8SEmil Constantinescu This function can be called before or after ts->vec_sol has been updated. 435f68a32c8SEmil Constantinescu Suppose we have a completion formula (b) and an embedded formula (be) of different order. 436f68a32c8SEmil Constantinescu We can write 437f68a32c8SEmil Constantinescu 438f68a32c8SEmil Constantinescu x1e = x0 + h be^T YdotRHS 439f68a32c8SEmil Constantinescu = x1 - h b^T YdotRHS + h be^T YdotRHS 440f68a32c8SEmil Constantinescu = x1 + h (be - b)^T YdotRHS 441f68a32c8SEmil Constantinescu 442f68a32c8SEmil Constantinescu so we can evaluate the method with different order even after the step has been optimistically completed. 443e4dd521cSBarry Smith */ 444f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done) 445f68a32c8SEmil Constantinescu { 446f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 447f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 448f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 449f68a32c8SEmil Constantinescu PetscReal h; 450f68a32c8SEmil Constantinescu PetscInt s = tab->s,j; 451f68a32c8SEmil Constantinescu PetscErrorCode ierr; 452f68a32c8SEmil Constantinescu 453f68a32c8SEmil Constantinescu PetscFunctionBegin; 454f68a32c8SEmil Constantinescu switch (rk->status) { 455f68a32c8SEmil Constantinescu case TS_STEP_INCOMPLETE: 456f68a32c8SEmil Constantinescu case TS_STEP_PENDING: 457f68a32c8SEmil Constantinescu h = ts->time_step; break; 458f68a32c8SEmil Constantinescu case TS_STEP_COMPLETE: 459be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 460f68a32c8SEmil Constantinescu default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 461f68a32c8SEmil Constantinescu } 462f68a32c8SEmil Constantinescu if (order == tab->order) { 463f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { 464f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 465f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->b[j]; 466f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 467f68a32c8SEmil Constantinescu } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 468f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 469f68a32c8SEmil Constantinescu } else if (order == tab->order-1) { 470f68a32c8SEmil Constantinescu if (!tab->bembed) goto unavailable; 471f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (be) */ 472f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 473f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 474f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 475f68a32c8SEmil Constantinescu } else { /* Rollback and re-complete using (be-b) */ 476f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 477f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 478f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 4790f7a1166SHong Zhang if (ts->vec_costintegral && ts->costintegralfwd) { 4800f7a1166SHong Zhang ierr = VecCopy(rk->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr); 4810f7a1166SHong Zhang } 482f68a32c8SEmil Constantinescu } 483f68a32c8SEmil Constantinescu if (done) *done = PETSC_TRUE; 484f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 485f68a32c8SEmil Constantinescu } 486f68a32c8SEmil Constantinescu unavailable: 487f68a32c8SEmil Constantinescu if (done) *done = PETSC_FALSE; 488a7fac7c2SEmil Constantinescu else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"RK '%s' of order %D cannot evaluate step at order %D. Consider using -ts_adapt_type none or a different method that has an embedded estimate.",tab->name,tab->order,order); 489f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 490f68a32c8SEmil Constantinescu } 491f68a32c8SEmil Constantinescu 4920f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts) 4930f7a1166SHong Zhang { 4940f7a1166SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 4950f7a1166SHong Zhang RKTableau tab = rk->tableau; 4960f7a1166SHong Zhang const PetscInt s = tab->s; 4970f7a1166SHong Zhang const PetscReal *b = tab->b,*c = tab->c; 4980f7a1166SHong Zhang Vec *Y = rk->Y; 4990f7a1166SHong Zhang PetscInt i; 5000f7a1166SHong Zhang PetscErrorCode ierr; 5010f7a1166SHong Zhang 5020f7a1166SHong Zhang PetscFunctionBegin; 5030f7a1166SHong Zhang /* backup cost integral */ 5040f7a1166SHong Zhang ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr); 5050f7a1166SHong Zhang for (i=s-1; i>=0; i--) { 5060f7a1166SHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 5070f7a1166SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,rk->ptime+rk->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 5080f7a1166SHong Zhang ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 5090f7a1166SHong Zhang } 5100f7a1166SHong Zhang PetscFunctionReturn(0); 5110f7a1166SHong Zhang } 5120f7a1166SHong Zhang 5130f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts) 5140f7a1166SHong Zhang { 5150f7a1166SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 5160f7a1166SHong Zhang RKTableau tab = rk->tableau; 5170f7a1166SHong Zhang const PetscInt s = tab->s; 5180f7a1166SHong Zhang const PetscReal *b = tab->b,*c = tab->c; 5190f7a1166SHong Zhang Vec *Y = rk->Y; 5200f7a1166SHong Zhang PetscInt i; 5210f7a1166SHong Zhang PetscErrorCode ierr; 5220f7a1166SHong Zhang 5230f7a1166SHong Zhang PetscFunctionBegin; 5240f7a1166SHong Zhang for (i=s-1; i>=0; i--) { 5250f7a1166SHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 5260f7a1166SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime-ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 5270f7a1166SHong Zhang ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 5280f7a1166SHong Zhang } 5290f7a1166SHong Zhang PetscFunctionReturn(0); 5300f7a1166SHong Zhang } 5310f7a1166SHong Zhang 532fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts) 533fecfb714SLisandro Dalcin { 534fecfb714SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 535fecfb714SLisandro Dalcin RKTableau tab = rk->tableau; 536fecfb714SLisandro Dalcin const PetscInt s = tab->s; 537fecfb714SLisandro Dalcin const PetscReal *b = tab->b; 538fecfb714SLisandro Dalcin PetscScalar *w = rk->work; 539fecfb714SLisandro Dalcin Vec *YdotRHS = rk->YdotRHS; 540fecfb714SLisandro Dalcin PetscInt j; 541fecfb714SLisandro Dalcin PetscReal h; 542fecfb714SLisandro Dalcin PetscErrorCode ierr; 543fecfb714SLisandro Dalcin 544fecfb714SLisandro Dalcin PetscFunctionBegin; 545fecfb714SLisandro Dalcin switch (rk->status) { 546fecfb714SLisandro Dalcin case TS_STEP_INCOMPLETE: 547fecfb714SLisandro Dalcin case TS_STEP_PENDING: 548fecfb714SLisandro Dalcin h = ts->time_step; break; 549fecfb714SLisandro Dalcin case TS_STEP_COMPLETE: 550fecfb714SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 551fecfb714SLisandro Dalcin default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 552fecfb714SLisandro Dalcin } 553fecfb714SLisandro Dalcin for (j=0; j<s; j++) w[j] = -h*b[j]; 554fecfb714SLisandro Dalcin ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr); 555fecfb714SLisandro Dalcin PetscFunctionReturn(0); 556fecfb714SLisandro Dalcin } 557fecfb714SLisandro Dalcin 558f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts) 559f68a32c8SEmil Constantinescu { 560f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 561f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 562f68a32c8SEmil Constantinescu const PetscInt s = tab->s; 563fecfb714SLisandro Dalcin const PetscReal *A = tab->A,*c = tab->c; 564f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 565f68a32c8SEmil Constantinescu Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 566d1905564SLisandro Dalcin PetscBool FSAL = tab->FSAL; 567f68a32c8SEmil Constantinescu TSAdapt adapt; 568fecfb714SLisandro Dalcin PetscInt i,j; 569be5899b3SLisandro Dalcin PetscInt rejections = 0; 570be5899b3SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 571be5899b3SLisandro Dalcin PetscReal next_time_step = ts->time_step; 572f68a32c8SEmil Constantinescu PetscErrorCode ierr; 573f68a32c8SEmil Constantinescu 574f68a32c8SEmil Constantinescu PetscFunctionBegin; 575d1905564SLisandro Dalcin if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE; 576d1905564SLisandro Dalcin if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); } 577d1905564SLisandro Dalcin 578f68a32c8SEmil Constantinescu rk->status = TS_STEP_INCOMPLETE; 579be5899b3SLisandro Dalcin while (!ts->reason && rk->status != TS_STEP_COMPLETE) { 580be5899b3SLisandro Dalcin PetscReal t = ts->ptime; 581f68a32c8SEmil Constantinescu PetscReal h = ts->time_step; 582f68a32c8SEmil Constantinescu for (i=0; i<s; i++) { 5839be3e283SDebojyoti Ghosh rk->stage_time = t + h*c[i]; 5849be3e283SDebojyoti Ghosh ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr); 585f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 586f68a32c8SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 587f68a32c8SEmil Constantinescu ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 5889be3e283SDebojyoti Ghosh ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 589f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 590be5899b3SLisandro Dalcin ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr); 591fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 5928f738a7cSLisandro Dalcin if (FSAL && !i) continue; 593f68a32c8SEmil Constantinescu ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 594f68a32c8SEmil Constantinescu } 595be5899b3SLisandro Dalcin 596be5899b3SLisandro Dalcin rk->status = TS_STEP_INCOMPLETE; 597f68a32c8SEmil Constantinescu ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 598f68a32c8SEmil Constantinescu rk->status = TS_STEP_PENDING; 599f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 600f68a32c8SEmil Constantinescu ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 6011917a363SLisandro Dalcin ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr); 602fecfb714SLisandro Dalcin ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 603be5899b3SLisandro Dalcin rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 604be5899b3SLisandro Dalcin if (!accept) { /* Roll back the current step */ 605fecfb714SLisandro Dalcin ierr = TSRollBack_RK(ts);CHKERRQ(ierr); 606be5899b3SLisandro Dalcin ts->time_step = next_time_step; 607be5899b3SLisandro Dalcin goto reject_step; 608be5899b3SLisandro Dalcin } 609be5899b3SLisandro Dalcin 6100f7a1166SHong Zhang if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation*/ 6110f7a1166SHong Zhang rk->ptime = ts->ptime; 6120f7a1166SHong Zhang rk->time_step = ts->time_step; 613493ed95fSHong Zhang } 614be5899b3SLisandro Dalcin 615f68a32c8SEmil Constantinescu ts->ptime += ts->time_step; 616f68a32c8SEmil Constantinescu ts->time_step = next_time_step; 617f68a32c8SEmil Constantinescu break; 618be5899b3SLisandro Dalcin 619be5899b3SLisandro Dalcin reject_step: 620fecfb714SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 621be5899b3SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 622be5899b3SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 623be5899b3SLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 624f68a32c8SEmil Constantinescu } 625f68a32c8SEmil Constantinescu } 626f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 627e4dd521cSBarry Smith } 628e4dd521cSBarry Smith 629f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts) 630f6a906c0SBarry Smith { 631f6a906c0SBarry Smith TS_RK *rk = (TS_RK*)ts->data; 632be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 633be5899b3SLisandro Dalcin PetscInt s = tab->s; 634f6a906c0SBarry Smith PetscErrorCode ierr; 635f6a906c0SBarry Smith 636f6a906c0SBarry Smith PetscFunctionBegin; 637f6a906c0SBarry Smith if (ts->adjointsetupcalled++) PetscFunctionReturn(0); 638abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr); 639f6a906c0SBarry Smith if(ts->vecs_sensip) { 640abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr); 641f6a906c0SBarry Smith } 642abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr); 643f6a906c0SBarry Smith PetscFunctionReturn(0); 644f6a906c0SBarry Smith } 645f6a906c0SBarry Smith 64642f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts) 647d2daff3dSHong Zhang { 648c235aa8dSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 649c235aa8dSHong Zhang RKTableau tab = rk->tableau; 650c235aa8dSHong Zhang const PetscInt s = tab->s; 651c235aa8dSHong Zhang const PetscReal *A = tab->A,*b = tab->b,*c = tab->c; 652c235aa8dSHong Zhang PetscScalar *w = rk->work; 653ad8e2604SHong Zhang Vec *Y = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu,*VecSensiTemp = rk->VecSensiTemp; 654b18ea86cSHong Zhang PetscInt i,j,nadj; 6553d3eaba7SBarry Smith PetscReal t = ts->ptime; 656d2daff3dSHong Zhang PetscErrorCode ierr; 657cef76868SBarry Smith PetscReal h = ts->time_step; 658cef76868SBarry Smith Mat J,Jp; 659c235aa8dSHong Zhang 660d2daff3dSHong Zhang PetscFunctionBegin; 661c235aa8dSHong Zhang rk->status = TS_STEP_INCOMPLETE; 662c235aa8dSHong Zhang for (i=s-1; i>=0; i--) { 663c235aa8dSHong Zhang rk->stage_time = t + h*(1.0-c[i]); 664abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 665b18ea86cSHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr); 666302440fdSBarry Smith ierr = VecScale(VecSensiTemp[nadj],-h*b[i]);CHKERRQ(ierr); 667c235aa8dSHong Zhang for (j=i+1; j<s; j++) { 668302440fdSBarry Smith ierr = VecAXPY(VecSensiTemp[nadj],-h*A[j*s+i],VecDeltaLam[nadj*s+j]);CHKERRQ(ierr); 669b18ea86cSHong Zhang } 670c235aa8dSHong Zhang } 671ad8e2604SHong Zhang /* Stage values of lambda */ 672c235aa8dSHong Zhang ierr = TSGetRHSJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 673c235aa8dSHong Zhang ierr = TSComputeRHSJacobian(ts,rk->stage_time,Y[i],J,Jp);CHKERRQ(ierr); 674493ed95fSHong Zhang if (ts->vec_costintegral) { 675493ed95fSHong Zhang ierr = TSAdjointComputeDRDYFunction(ts,rk->stage_time,Y[i],ts->vecs_drdy);CHKERRQ(ierr); 676493ed95fSHong Zhang } 677abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 678b18ea86cSHong Zhang ierr = MatMultTranspose(J,VecSensiTemp[nadj],VecDeltaLam[nadj*s+i]);CHKERRQ(ierr); 679493ed95fSHong Zhang if (ts->vec_costintegral) { 680493ed95fSHong Zhang ierr = VecAXPY(VecDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdy[nadj]);CHKERRQ(ierr); 681493ed95fSHong Zhang } 682b18ea86cSHong Zhang } 6836497ce14SHong Zhang 684ad8e2604SHong Zhang /* Stage values of mu */ 6856497ce14SHong Zhang if(ts->vecs_sensip) { 6865bf8c567SBarry Smith ierr = TSAdjointComputeRHSJacobian(ts,rk->stage_time,Y[i],ts->Jacp);CHKERRQ(ierr); 687493ed95fSHong Zhang if (ts->vec_costintegral) { 688493ed95fSHong Zhang ierr = TSAdjointComputeDRDPFunction(ts,rk->stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr); 689493ed95fSHong Zhang } 690493ed95fSHong Zhang 691abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 692ad8e2604SHong Zhang ierr = MatMultTranspose(ts->Jacp,VecSensiTemp[nadj],VecDeltaMu[nadj*s+i]);CHKERRQ(ierr); 693493ed95fSHong Zhang if (ts->vec_costintegral) { 694493ed95fSHong Zhang ierr = VecAXPY(VecDeltaMu[nadj*s+i],-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr); 695493ed95fSHong Zhang } 696ad8e2604SHong Zhang } 697c235aa8dSHong Zhang } 6986497ce14SHong Zhang } 699c235aa8dSHong Zhang 700c235aa8dSHong Zhang for (j=0; j<s; j++) w[j] = 1.0; 701abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 702b18ea86cSHong Zhang ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr); 7036497ce14SHong Zhang if(ts->vecs_sensip) { 704ad8e2604SHong Zhang ierr = VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);CHKERRQ(ierr); 705b18ea86cSHong Zhang } 7066497ce14SHong Zhang } 707c235aa8dSHong Zhang rk->status = TS_STEP_COMPLETE; 708d2daff3dSHong Zhang PetscFunctionReturn(0); 709d2daff3dSHong Zhang } 710d2daff3dSHong Zhang 711f68a32c8SEmil Constantinescu static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X) 712f68a32c8SEmil Constantinescu { 713f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 7143a8a9803SLisandro Dalcin PetscInt s = rk->tableau->s,p = rk->tableau->p,i,j; 715f68a32c8SEmil Constantinescu PetscReal h; 716f68a32c8SEmil Constantinescu PetscReal tt,t; 717f68a32c8SEmil Constantinescu PetscScalar *b; 718f68a32c8SEmil Constantinescu const PetscReal *B = rk->tableau->binterp; 719f68a32c8SEmil Constantinescu PetscErrorCode ierr; 720e4dd521cSBarry Smith 721f68a32c8SEmil Constantinescu PetscFunctionBegin; 722f68a32c8SEmil Constantinescu if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 723be5899b3SLisandro Dalcin 724f68a32c8SEmil Constantinescu switch (rk->status) { 725f68a32c8SEmil Constantinescu case TS_STEP_INCOMPLETE: 726f68a32c8SEmil Constantinescu case TS_STEP_PENDING: 727f68a32c8SEmil Constantinescu h = ts->time_step; 728f68a32c8SEmil Constantinescu t = (itime - ts->ptime)/h; 729f68a32c8SEmil Constantinescu break; 730f68a32c8SEmil Constantinescu case TS_STEP_COMPLETE: 731be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; 732f68a32c8SEmil Constantinescu t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 733f68a32c8SEmil Constantinescu break; 734f68a32c8SEmil Constantinescu default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 735e4dd521cSBarry Smith } 736785e854fSJed Brown ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 737f68a32c8SEmil Constantinescu for (i=0; i<s; i++) b[i] = 0; 7383a8a9803SLisandro Dalcin for (j=0,tt=t; j<p; j++,tt*=t) { 739f68a32c8SEmil Constantinescu for (i=0; i<s; i++) { 7403a8a9803SLisandro Dalcin b[i] += h * B[i*p+j] * tt; 741e4dd521cSBarry Smith } 742f68a32c8SEmil Constantinescu } 743be5899b3SLisandro Dalcin 744f68a32c8SEmil Constantinescu ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr); 745f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr); 746be5899b3SLisandro Dalcin 747f68a32c8SEmil Constantinescu ierr = PetscFree(b);CHKERRQ(ierr); 748e4dd521cSBarry Smith PetscFunctionReturn(0); 749e4dd521cSBarry Smith } 750e4dd521cSBarry Smith 751e4dd521cSBarry Smith /*------------------------------------------------------------*/ 752be5899b3SLisandro Dalcin 753be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts) 754be5899b3SLisandro Dalcin { 755be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 756be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 757be5899b3SLisandro Dalcin PetscErrorCode ierr; 758be5899b3SLisandro Dalcin 759be5899b3SLisandro Dalcin PetscFunctionBegin; 760be5899b3SLisandro Dalcin if (!tab) PetscFunctionReturn(0); 761be5899b3SLisandro Dalcin ierr = PetscFree(rk->work);CHKERRQ(ierr); 762be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr); 763be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr); 764be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr); 765be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr); 766be5899b3SLisandro Dalcin PetscFunctionReturn(0); 767be5899b3SLisandro Dalcin } 768be5899b3SLisandro Dalcin 769277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts) 770e4dd521cSBarry Smith { 7715f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 7726849ba73SBarry Smith PetscErrorCode ierr; 773e4dd521cSBarry Smith 774e4dd521cSBarry Smith PetscFunctionBegin; 775be5899b3SLisandro Dalcin ierr = TSRKTableauReset(ts);CHKERRQ(ierr); 7760f7a1166SHong Zhang ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr); 777abc2977eSBarry Smith ierr = VecDestroyVecs(ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr); 778277b19d0SLisandro Dalcin PetscFunctionReturn(0); 779e4dd521cSBarry Smith } 780277b19d0SLisandro Dalcin 781277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_RK(TS ts) 782277b19d0SLisandro Dalcin { 783277b19d0SLisandro Dalcin PetscErrorCode ierr; 784277b19d0SLisandro Dalcin 785277b19d0SLisandro Dalcin PetscFunctionBegin; 786277b19d0SLisandro Dalcin ierr = TSReset_RK(ts);CHKERRQ(ierr); 787277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 788f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr); 789f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr); 790f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 791f68a32c8SEmil Constantinescu } 792f68a32c8SEmil Constantinescu 793f68a32c8SEmil Constantinescu 794f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx) 795f68a32c8SEmil Constantinescu { 796f68a32c8SEmil Constantinescu PetscFunctionBegin; 797f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 798f68a32c8SEmil Constantinescu } 799f68a32c8SEmil Constantinescu 800f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 801f68a32c8SEmil Constantinescu { 802f68a32c8SEmil Constantinescu PetscFunctionBegin; 803f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 804f68a32c8SEmil Constantinescu } 805f68a32c8SEmil Constantinescu 806f68a32c8SEmil Constantinescu 807f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx) 808f68a32c8SEmil Constantinescu { 809f68a32c8SEmil Constantinescu PetscFunctionBegin; 810f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 811f68a32c8SEmil Constantinescu } 812f68a32c8SEmil Constantinescu 813f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 814f68a32c8SEmil Constantinescu { 815f68a32c8SEmil Constantinescu 816f68a32c8SEmil Constantinescu PetscFunctionBegin; 817f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 818f68a32c8SEmil Constantinescu } 819c235aa8dSHong Zhang /* 820d2daff3dSHong Zhang static PetscErrorCode RKSetAdjCoe(RKTableau tab) 821d2daff3dSHong Zhang { 822d2daff3dSHong Zhang PetscReal *A,*b; 823d2daff3dSHong Zhang PetscInt s,i,j; 824d2daff3dSHong Zhang PetscErrorCode ierr; 825d2daff3dSHong Zhang 826d2daff3dSHong Zhang PetscFunctionBegin; 827d2daff3dSHong Zhang s = tab->s; 828d2daff3dSHong Zhang ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr); 829d2daff3dSHong Zhang 830d2daff3dSHong Zhang for (i=0; i<s; i++) 831d2daff3dSHong Zhang for (j=0; j<s; j++) { 832d2daff3dSHong Zhang A[i*s+j] = (tab->b[s-1-i]==0) ? 0: -tab->A[s-1-i+(s-1-j)*s] * tab->b[s-1-j] / tab->b[s-1-i]; 833d2daff3dSHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr); 834d2daff3dSHong Zhang } 835d2daff3dSHong Zhang 836d2daff3dSHong Zhang for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i]; 837d2daff3dSHong Zhang 838d2daff3dSHong Zhang ierr = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 839d2daff3dSHong Zhang ierr = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr); 840d2daff3dSHong Zhang ierr = PetscFree2(A,b);CHKERRQ(ierr); 841d2daff3dSHong Zhang PetscFunctionReturn(0); 842d2daff3dSHong Zhang } 843c235aa8dSHong Zhang */ 844be5899b3SLisandro Dalcin 845be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts) 846be5899b3SLisandro Dalcin { 847be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 848be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 849be5899b3SLisandro Dalcin PetscErrorCode ierr; 850be5899b3SLisandro Dalcin 851be5899b3SLisandro Dalcin PetscFunctionBegin; 852be5899b3SLisandro Dalcin ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr); 853be5899b3SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr); 854be5899b3SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr); 855be5899b3SLisandro Dalcin PetscFunctionReturn(0); 856be5899b3SLisandro Dalcin } 857be5899b3SLisandro Dalcin 858be5899b3SLisandro Dalcin 859f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts) 860f68a32c8SEmil Constantinescu { 861f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 862f68a32c8SEmil Constantinescu PetscErrorCode ierr; 863f68a32c8SEmil Constantinescu DM dm; 864f68a32c8SEmil Constantinescu 865f68a32c8SEmil Constantinescu PetscFunctionBegin; 8663dd83b38SBarry Smith ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr); 867be5899b3SLisandro Dalcin ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr); 8680f7a1166SHong Zhang if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 8690f7a1166SHong Zhang ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr); 8700f7a1166SHong Zhang } 871f68a32c8SEmil Constantinescu ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 872f68a32c8SEmil Constantinescu if (dm) { 873f68a32c8SEmil Constantinescu ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 874f68a32c8SEmil Constantinescu ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 875f68a32c8SEmil Constantinescu } 876e4dd521cSBarry Smith PetscFunctionReturn(0); 877e4dd521cSBarry Smith } 878d2daff3dSHong Zhang 879f6a906c0SBarry Smith 880e4dd521cSBarry Smith /*------------------------------------------------------------*/ 881e4dd521cSBarry Smith 8824416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts) 883e4dd521cSBarry Smith { 884be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 885dfbe8321SBarry Smith PetscErrorCode ierr; 886262ff7bbSBarry Smith 887e4dd521cSBarry Smith PetscFunctionBegin; 888e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr); 889f68a32c8SEmil Constantinescu { 890f68a32c8SEmil Constantinescu RKTableauLink link; 891f68a32c8SEmil Constantinescu PetscInt count,choice; 892f68a32c8SEmil Constantinescu PetscBool flg; 893f68a32c8SEmil Constantinescu const char **namelist; 894f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) ; 895785e854fSJed Brown ierr = PetscMalloc1(count,&namelist);CHKERRQ(ierr); 896f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 897be5899b3SLisandro Dalcin ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr); 898be5899b3SLisandro Dalcin if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);} 899f68a32c8SEmil Constantinescu ierr = PetscFree(namelist);CHKERRQ(ierr); 900f68a32c8SEmil Constantinescu } 901262ff7bbSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 902e4dd521cSBarry Smith PetscFunctionReturn(0); 903e4dd521cSBarry Smith } 904e4dd521cSBarry Smith 905f68a32c8SEmil Constantinescu static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[]) 906f68a32c8SEmil Constantinescu { 907f68a32c8SEmil Constantinescu PetscErrorCode ierr; 908f68a32c8SEmil Constantinescu PetscInt i; 909f68a32c8SEmil Constantinescu size_t left,count; 910f68a32c8SEmil Constantinescu char *p; 911f68a32c8SEmil Constantinescu 912f68a32c8SEmil Constantinescu PetscFunctionBegin; 913f68a32c8SEmil Constantinescu for (i=0,p=buf,left=len; i<n; i++) { 914f68a32c8SEmil Constantinescu ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr); 915f68a32c8SEmil Constantinescu if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer"); 916f68a32c8SEmil Constantinescu left -= count; 917f68a32c8SEmil Constantinescu p += count; 918f68a32c8SEmil Constantinescu *p++ = ' '; 919f68a32c8SEmil Constantinescu } 920f68a32c8SEmil Constantinescu p[i ? 0 : -1] = 0; 921f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 922f68a32c8SEmil Constantinescu } 923f68a32c8SEmil Constantinescu 9245f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer) 925e4dd521cSBarry Smith { 9265f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 9278370ee3bSLisandro Dalcin PetscBool iascii; 928dfbe8321SBarry Smith PetscErrorCode ierr; 9292cdf8120Sasbjorn 930e4dd521cSBarry Smith PetscFunctionBegin; 931251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 9328370ee3bSLisandro Dalcin if (iascii) { 9339c334d8fSLisandro Dalcin RKTableau tab = rk->tableau; 934f68a32c8SEmil Constantinescu TSRKType rktype; 935f68a32c8SEmil Constantinescu char buf[512]; 936f68a32c8SEmil Constantinescu ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr); 9370c37eb8eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," RK: %s\n",rktype);CHKERRQ(ierr); 9380c37eb8eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," Order: %D\n",tab->order);CHKERRQ(ierr); 9390c37eb8eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr); 940f68a32c8SEmil Constantinescu ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 941f68a32c8SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," Abscissa c = %s\n",buf);CHKERRQ(ierr); 9428370ee3bSLisandro Dalcin } 9439c334d8fSLisandro Dalcin if (ts->adapt) {ierr = TSAdaptView(ts->adapt,viewer);CHKERRQ(ierr);} 944f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 945f68a32c8SEmil Constantinescu } 946f68a32c8SEmil Constantinescu 947f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer) 948f68a32c8SEmil Constantinescu { 949f68a32c8SEmil Constantinescu PetscErrorCode ierr; 9509c334d8fSLisandro Dalcin TSAdapt adapt; 951f68a32c8SEmil Constantinescu 952f68a32c8SEmil Constantinescu PetscFunctionBegin; 9539c334d8fSLisandro Dalcin ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 9549c334d8fSLisandro Dalcin ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr); 955f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 956f68a32c8SEmil Constantinescu } 957f68a32c8SEmil Constantinescu 958f68a32c8SEmil Constantinescu /*@C 959f68a32c8SEmil Constantinescu TSRKSetType - Set the type of RK scheme 960f68a32c8SEmil Constantinescu 961f68a32c8SEmil Constantinescu Logically collective 962f68a32c8SEmil Constantinescu 963f68a32c8SEmil Constantinescu Input Parameter: 964f68a32c8SEmil Constantinescu + ts - timestepping context 965f68a32c8SEmil Constantinescu - rktype - type of RK-scheme 966f68a32c8SEmil Constantinescu 967f68a32c8SEmil Constantinescu Level: intermediate 968f68a32c8SEmil Constantinescu 969f68a32c8SEmil Constantinescu .seealso: TSRKGetType(), TSRK, TSRK2, TSRK3, TSRKPRSSP2, TSRK4, TSRK5 970f68a32c8SEmil Constantinescu @*/ 971f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype) 972f68a32c8SEmil Constantinescu { 973f68a32c8SEmil Constantinescu PetscErrorCode ierr; 974f68a32c8SEmil Constantinescu 975f68a32c8SEmil Constantinescu PetscFunctionBegin; 976f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 977b92453a8SLisandro Dalcin PetscValidCharPointer(rktype,2); 978f68a32c8SEmil Constantinescu ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr); 979f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 980f68a32c8SEmil Constantinescu } 981f68a32c8SEmil Constantinescu 982f68a32c8SEmil Constantinescu /*@C 983f68a32c8SEmil Constantinescu TSRKGetType - Get the type of RK scheme 984f68a32c8SEmil Constantinescu 985f68a32c8SEmil Constantinescu Logically collective 986f68a32c8SEmil Constantinescu 987f68a32c8SEmil Constantinescu Input Parameter: 988f68a32c8SEmil Constantinescu . ts - timestepping context 989f68a32c8SEmil Constantinescu 990f68a32c8SEmil Constantinescu Output Parameter: 991f68a32c8SEmil Constantinescu . rktype - type of RK-scheme 992f68a32c8SEmil Constantinescu 993f68a32c8SEmil Constantinescu Level: intermediate 994f68a32c8SEmil Constantinescu 995f68a32c8SEmil Constantinescu .seealso: TSRKGetType() 996f68a32c8SEmil Constantinescu @*/ 997f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype) 998f68a32c8SEmil Constantinescu { 999f68a32c8SEmil Constantinescu PetscErrorCode ierr; 1000f68a32c8SEmil Constantinescu 1001f68a32c8SEmil Constantinescu PetscFunctionBegin; 1002f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1003f68a32c8SEmil Constantinescu ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr); 1004f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1005f68a32c8SEmil Constantinescu } 1006f68a32c8SEmil Constantinescu 1007560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype) 1008f68a32c8SEmil Constantinescu { 1009f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 1010f68a32c8SEmil Constantinescu 1011f68a32c8SEmil Constantinescu PetscFunctionBegin; 1012f68a32c8SEmil Constantinescu *rktype = rk->tableau->name; 1013f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1014f68a32c8SEmil Constantinescu } 1015560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype) 1016f68a32c8SEmil Constantinescu { 1017f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 1018f68a32c8SEmil Constantinescu PetscErrorCode ierr; 1019f68a32c8SEmil Constantinescu PetscBool match; 1020f68a32c8SEmil Constantinescu RKTableauLink link; 1021f68a32c8SEmil Constantinescu 1022f68a32c8SEmil Constantinescu PetscFunctionBegin; 1023f68a32c8SEmil Constantinescu if (rk->tableau) { 1024f68a32c8SEmil Constantinescu ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr); 1025f68a32c8SEmil Constantinescu if (match) PetscFunctionReturn(0); 1026f68a32c8SEmil Constantinescu } 1027f68a32c8SEmil Constantinescu for (link = RKTableauList; link; link=link->next) { 1028f68a32c8SEmil Constantinescu ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr); 1029f68a32c8SEmil Constantinescu if (match) { 1030be5899b3SLisandro Dalcin if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);} 1031f68a32c8SEmil Constantinescu rk->tableau = &link->tab; 1032be5899b3SLisandro Dalcin if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);} 10332ffb9264SLisandro Dalcin ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 1034f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1035f68a32c8SEmil Constantinescu } 1036f68a32c8SEmil Constantinescu } 1037f68a32c8SEmil Constantinescu SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype); 1038e4dd521cSBarry Smith PetscFunctionReturn(0); 1039e4dd521cSBarry Smith } 1040e4dd521cSBarry Smith 1041ff22ae23SHong Zhang static PetscErrorCode TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y) 1042ff22ae23SHong Zhang { 1043ff22ae23SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 1044ff22ae23SHong Zhang 1045ff22ae23SHong Zhang PetscFunctionBegin; 1046ff22ae23SHong Zhang *ns = rk->tableau->s; 1047d2daff3dSHong Zhang if(Y) *Y = rk->Y; 1048ff22ae23SHong Zhang PetscFunctionReturn(0); 1049ff22ae23SHong Zhang } 1050ff22ae23SHong Zhang 1051ff22ae23SHong Zhang 1052e4dd521cSBarry Smith /* ------------------------------------------------------------ */ 1053ebe3b25bSBarry Smith /*MC 1054f68a32c8SEmil Constantinescu TSRK - ODE and DAE solver using Runge-Kutta schemes 1055ebe3b25bSBarry Smith 1056f68a32c8SEmil Constantinescu The user should provide the right hand side of the equation 1057f68a32c8SEmil Constantinescu using TSSetRHSFunction(). 1058ebe3b25bSBarry Smith 1059f68a32c8SEmil Constantinescu Notes: 106098164e13SLisandro Dalcin The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type 1061ebe3b25bSBarry Smith 1062d41222bbSBarry Smith Level: beginner 1063d41222bbSBarry Smith 1064f68a32c8SEmil Constantinescu .seealso: TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3, 1065f68a32c8SEmil Constantinescu TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister() 1066ebe3b25bSBarry Smith 1067ebe3b25bSBarry Smith M*/ 10688cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts) 1069e4dd521cSBarry Smith { 10701566a47fSLisandro Dalcin TS_RK *rk; 1071dfbe8321SBarry Smith PetscErrorCode ierr; 1072e4dd521cSBarry Smith 1073e4dd521cSBarry Smith PetscFunctionBegin; 1074f68a32c8SEmil Constantinescu ierr = TSRKInitializePackage();CHKERRQ(ierr); 1075f68a32c8SEmil Constantinescu 1076f68a32c8SEmil Constantinescu ts->ops->reset = TSReset_RK; 10775f70b478SJed Brown ts->ops->destroy = TSDestroy_RK; 10785f70b478SJed Brown ts->ops->view = TSView_RK; 1079f68a32c8SEmil Constantinescu ts->ops->load = TSLoad_RK; 1080f68a32c8SEmil Constantinescu ts->ops->setup = TSSetUp_RK; 108142f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_RK; 1082f68a32c8SEmil Constantinescu ts->ops->step = TSStep_RK; 1083f68a32c8SEmil Constantinescu ts->ops->interpolate = TSInterpolate_RK; 1084f68a32c8SEmil Constantinescu ts->ops->evaluatestep = TSEvaluateStep_RK; 1085fecfb714SLisandro Dalcin ts->ops->rollback = TSRollBack_RK; 1086f68a32c8SEmil Constantinescu ts->ops->setfromoptions = TSSetFromOptions_RK; 1087ff22ae23SHong Zhang ts->ops->getstages = TSGetStages_RK; 108842f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_RK; 1089e4dd521cSBarry Smith 10900f7a1166SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_RK; 10910f7a1166SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_RK; 10920f7a1166SHong Zhang 10931566a47fSLisandro Dalcin ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr); 10941566a47fSLisandro Dalcin ts->data = (void*)rk; 1095e4dd521cSBarry Smith 1096f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr); 1097f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr); 1098be5899b3SLisandro Dalcin 1099be5899b3SLisandro Dalcin ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr); 1100e4dd521cSBarry Smith PetscFunctionReturn(0); 1101e4dd521cSBarry Smith } 1102