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 165f68a32c8SEmil Constantinescu A[3][3] = {{0,0,0}, 166f68a32c8SEmil Constantinescu {2.0/3.0,0,0}, 167f68a32c8SEmil Constantinescu {-1.0/3.0,1.0,0}}, 168f68a32c8SEmil Constantinescu b[3] = {0.25,0.5,0.25}; 1693a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 170fdefd5e5SDebojyoti Ghosh } 171fdefd5e5SDebojyoti Ghosh { 172fdefd5e5SDebojyoti Ghosh const PetscReal 173fdefd5e5SDebojyoti Ghosh A[4][4] = {{0,0,0,0}, 174fdefd5e5SDebojyoti Ghosh {1.0/2.0,0,0,0}, 175fdefd5e5SDebojyoti Ghosh {0,3.0/4.0,0,0}, 176fdefd5e5SDebojyoti Ghosh {2.0/9.0,1.0/3.0,4.0/9.0,0}}, 177fdefd5e5SDebojyoti Ghosh b[4] = {2.0/9.0,1.0/3.0,4.0/9.0,0}, 178fdefd5e5SDebojyoti Ghosh bembed[4] = {7.0/24.0,1.0/4.0,1.0/3.0,1.0/8.0}; 1793a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 180db046809SBarry Smith } 181f68a32c8SEmil Constantinescu { 182f68a32c8SEmil Constantinescu const PetscReal 183f68a32c8SEmil Constantinescu A[4][4] = {{0,0,0,0}, 184f68a32c8SEmil Constantinescu {0.5,0,0,0}, 185f68a32c8SEmil Constantinescu {0,0.5,0,0}, 186f68a32c8SEmil Constantinescu {0,0,1.0,0}}, 187f68a32c8SEmil Constantinescu b[4] = {1.0/6.0,1.0/3.0,1.0/3.0,1.0/6.0}; 1883a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 189db046809SBarry Smith } 190f68a32c8SEmil Constantinescu { 191f68a32c8SEmil Constantinescu const PetscReal 192f68a32c8SEmil Constantinescu A[6][6] = {{0,0,0,0,0,0}, 193f68a32c8SEmil Constantinescu {0.25,0,0,0,0,0}, 194f68a32c8SEmil Constantinescu {3.0/32.0,9.0/32.0,0,0,0,0}, 195f68a32c8SEmil Constantinescu {1932.0/2197.0,-7200.0/2197.0,7296.0/2197.0,0,0,0}, 196f68a32c8SEmil Constantinescu {439.0/216.0,-8.0,3680.0/513.0,-845.0/4104.0,0,0}, 197f68a32c8SEmil Constantinescu {-8.0/27.0,2.0,-3544.0/2565.0,1859.0/4104.0,-11.0/40.0,0}}, 198f68a32c8SEmil Constantinescu b[6] = {16.0/135.0,0,6656.0/12825.0,28561.0/56430.0,-9.0/50.0,2.0/55.0}, 199f68a32c8SEmil Constantinescu bembed[6] = {25.0/216.0,0,1408.0/2565.0,2197.0/4104.0,-1.0/5.0,0}; 2003a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 201fdefd5e5SDebojyoti Ghosh } 202fdefd5e5SDebojyoti Ghosh { 203fdefd5e5SDebojyoti Ghosh const PetscReal 204fdefd5e5SDebojyoti Ghosh A[7][7] = {{0,0,0,0,0,0,0}, 205fdefd5e5SDebojyoti Ghosh {1.0/5.0,0,0,0,0,0,0}, 206fdefd5e5SDebojyoti Ghosh {3.0/40.0,9.0/40.0,0,0,0,0,0}, 207fdefd5e5SDebojyoti Ghosh {44.0/45.0,-56.0/15.0,32.0/9.0,0,0,0,0}, 208fdefd5e5SDebojyoti Ghosh {19372.0/6561.0,-25360.0/2187.0,64448.0/6561.0,-212.0/729.0,0,0,0}, 209fdefd5e5SDebojyoti Ghosh {9017.0/3168.0,-355.0/33.0,46732.0/5247.0,49.0/176.0,-5103.0/18656.0,0,0}, 210fdefd5e5SDebojyoti Ghosh {35.0/384.0,0,500.0/1113.0,125.0/192.0,-2187.0/6784.0,11.0/84.0,0}}, 211fdefd5e5SDebojyoti 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}, 212fdefd5e5SDebojyoti 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}; 2133a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 214f68a32c8SEmil Constantinescu } 21505e23783SLisandro Dalcin { 21605e23783SLisandro Dalcin const PetscReal 21705e23783SLisandro Dalcin A[8][8] = {{0,0,0,0,0,0,0,0}, 21805e23783SLisandro Dalcin {1.0/6.0,0,0,0,0,0,0,0}, 21905e23783SLisandro Dalcin {2.0/27.0,4.0/27.0,0,0,0,0,0,0}, 22005e23783SLisandro Dalcin {183.0/1372.0,-162.0/343.0,1053.0/1372.0,0,0,0,0,0}, 22105e23783SLisandro Dalcin {68.0/297.0,-4.0/11.0,42.0/143.0,1960.0/3861.0,0,0,0,0}, 22205e23783SLisandro Dalcin {597.0/22528.0,81.0/352.0,63099.0/585728.0,58653.0/366080.0,4617.0/20480.0,0,0,0}, 22305e23783SLisandro 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}, 22405e23783SLisandro 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}}, 22505e23783SLisandro 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}, 22605e23783SLisandro 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}; 22705e23783SLisandro Dalcin ierr = TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 22805e23783SLisandro Dalcin } 229db046809SBarry Smith PetscFunctionReturn(0); 230db046809SBarry Smith } 231db046809SBarry Smith 232f68a32c8SEmil Constantinescu /*@C 233f68a32c8SEmil Constantinescu TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister(). 234f68a32c8SEmil Constantinescu 235f68a32c8SEmil Constantinescu Not Collective 236f68a32c8SEmil Constantinescu 237f68a32c8SEmil Constantinescu Level: advanced 238f68a32c8SEmil Constantinescu 239f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy 240f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll() 241f68a32c8SEmil Constantinescu @*/ 242f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void) 243e4dd521cSBarry Smith { 244f68a32c8SEmil Constantinescu PetscErrorCode ierr; 245f68a32c8SEmil Constantinescu RKTableauLink link; 246f68a32c8SEmil Constantinescu 247f68a32c8SEmil Constantinescu PetscFunctionBegin; 248f68a32c8SEmil Constantinescu while ((link = RKTableauList)) { 249f68a32c8SEmil Constantinescu RKTableau t = &link->tab; 250f68a32c8SEmil Constantinescu RKTableauList = link->next; 251f68a32c8SEmil Constantinescu ierr = PetscFree3(t->A,t->b,t->c); CHKERRQ(ierr); 252f68a32c8SEmil Constantinescu ierr = PetscFree (t->bembed); CHKERRQ(ierr); 253f68a32c8SEmil Constantinescu ierr = PetscFree (t->binterp); CHKERRQ(ierr); 254f68a32c8SEmil Constantinescu ierr = PetscFree (t->name); CHKERRQ(ierr); 255f68a32c8SEmil Constantinescu ierr = PetscFree (link); CHKERRQ(ierr); 256f68a32c8SEmil Constantinescu } 257f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_FALSE; 258f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 259f68a32c8SEmil Constantinescu } 260f68a32c8SEmil Constantinescu 261f68a32c8SEmil Constantinescu /*@C 262f68a32c8SEmil Constantinescu TSRKInitializePackage - This function initializes everything in the TSRK package. It is called 263f68a32c8SEmil Constantinescu from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RK() 264f68a32c8SEmil Constantinescu when using static libraries. 265f68a32c8SEmil Constantinescu 266f68a32c8SEmil Constantinescu Level: developer 267f68a32c8SEmil Constantinescu 268f68a32c8SEmil Constantinescu .keywords: TS, TSRK, initialize, package 269f68a32c8SEmil Constantinescu .seealso: PetscInitialize() 270f68a32c8SEmil Constantinescu @*/ 271f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void) 272f68a32c8SEmil Constantinescu { 273186e87acSLisandro Dalcin PetscErrorCode ierr; 274e4dd521cSBarry Smith 275e4dd521cSBarry Smith PetscFunctionBegin; 276f68a32c8SEmil Constantinescu if (TSRKPackageInitialized) PetscFunctionReturn(0); 277f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_TRUE; 278f68a32c8SEmil Constantinescu ierr = TSRKRegisterAll();CHKERRQ(ierr); 279f68a32c8SEmil Constantinescu ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr); 280f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 281f68a32c8SEmil Constantinescu } 282186e87acSLisandro Dalcin 283f68a32c8SEmil Constantinescu /*@C 284f68a32c8SEmil Constantinescu TSRKFinalizePackage - This function destroys everything in the TSRK package. It is 285f68a32c8SEmil Constantinescu called from PetscFinalize(). 286186e87acSLisandro Dalcin 287f68a32c8SEmil Constantinescu Level: developer 288f68a32c8SEmil Constantinescu 289f68a32c8SEmil Constantinescu .keywords: Petsc, destroy, package 290f68a32c8SEmil Constantinescu .seealso: PetscFinalize() 291f68a32c8SEmil Constantinescu @*/ 292f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void) 293f68a32c8SEmil Constantinescu { 294f68a32c8SEmil Constantinescu PetscErrorCode ierr; 295f68a32c8SEmil Constantinescu 296f68a32c8SEmil Constantinescu PetscFunctionBegin; 297f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_FALSE; 298f68a32c8SEmil Constantinescu ierr = TSRKRegisterDestroy();CHKERRQ(ierr); 299f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 300f68a32c8SEmil Constantinescu } 301f68a32c8SEmil Constantinescu 302f68a32c8SEmil Constantinescu /*@C 303f68a32c8SEmil Constantinescu TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 304f68a32c8SEmil Constantinescu 305f68a32c8SEmil Constantinescu Not Collective, but the same schemes should be registered on all processes on which they will be used 306f68a32c8SEmil Constantinescu 307f68a32c8SEmil Constantinescu Input Parameters: 308f68a32c8SEmil Constantinescu + name - identifier for method 309f68a32c8SEmil Constantinescu . order - approximation order of method 310f68a32c8SEmil Constantinescu . s - number of stages, this is the dimension of the matrices below 311f68a32c8SEmil Constantinescu . A - stage coefficients (dimension s*s, row-major) 312f68a32c8SEmil Constantinescu . b - step completion table (dimension s; NULL to use last row of A) 313f68a32c8SEmil Constantinescu . c - abscissa (dimension s; NULL to use row sums of A) 314f68a32c8SEmil Constantinescu . bembed - completion table for embedded method (dimension s; NULL if not available) 3153a8a9803SLisandro Dalcin . p - Order of the interpolation scheme, equal to the number of columns of binterp 3163a8a9803SLisandro Dalcin - binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1) 317f68a32c8SEmil Constantinescu 318f68a32c8SEmil Constantinescu Notes: 319f68a32c8SEmil Constantinescu Several RK methods are provided, this function is only needed to create new methods. 320f68a32c8SEmil Constantinescu 321f68a32c8SEmil Constantinescu Level: advanced 322f68a32c8SEmil Constantinescu 323f68a32c8SEmil Constantinescu .keywords: TS, register 324f68a32c8SEmil Constantinescu 325f68a32c8SEmil Constantinescu .seealso: TSRK 326f68a32c8SEmil Constantinescu @*/ 327f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s, 328f68a32c8SEmil Constantinescu const PetscReal A[],const PetscReal b[],const PetscReal c[], 3293a8a9803SLisandro Dalcin const PetscReal bembed[],PetscInt p,const PetscReal binterp[]) 330f68a32c8SEmil Constantinescu { 331f68a32c8SEmil Constantinescu PetscErrorCode ierr; 332f68a32c8SEmil Constantinescu RKTableauLink link; 333f68a32c8SEmil Constantinescu RKTableau t; 334f68a32c8SEmil Constantinescu PetscInt i,j; 335f68a32c8SEmil Constantinescu 336f68a32c8SEmil Constantinescu PetscFunctionBegin; 3373a8a9803SLisandro Dalcin PetscValidCharPointer(name,1); 3383a8a9803SLisandro Dalcin PetscValidRealPointer(A,4); 3393a8a9803SLisandro Dalcin if (b) PetscValidRealPointer(b,5); 3403a8a9803SLisandro Dalcin if (c) PetscValidRealPointer(c,6); 3413a8a9803SLisandro Dalcin if (bembed) PetscValidRealPointer(bembed,7); 3423a8a9803SLisandro Dalcin if (binterp || p > 1) PetscValidRealPointer(binterp,9); 3433a8a9803SLisandro Dalcin 34495dccacaSBarry Smith ierr = PetscNew(&link);CHKERRQ(ierr); 345f68a32c8SEmil Constantinescu t = &link->tab; 3463a8a9803SLisandro Dalcin 347f68a32c8SEmil Constantinescu ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 348f68a32c8SEmil Constantinescu t->order = order; 349f68a32c8SEmil Constantinescu t->s = s; 350dcca6d9dSJed Brown ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr); 351f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 352f68a32c8SEmil Constantinescu if (b) { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); } 353f68a32c8SEmil Constantinescu else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i]; 354f68a32c8SEmil Constantinescu if (c) { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); } 355f68a32c8SEmil 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]; 356d760c35bSDebojyoti Ghosh t->FSAL = PETSC_TRUE; 357d760c35bSDebojyoti Ghosh for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE; 3583a8a9803SLisandro Dalcin 359f68a32c8SEmil Constantinescu if (bembed) { 360785e854fSJed Brown ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr); 361f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr); 362f68a32c8SEmil Constantinescu } 363f68a32c8SEmil Constantinescu 3643a8a9803SLisandro Dalcin if (!binterp) { p = 1; binterp = t->b; } 3653a8a9803SLisandro Dalcin t->p = p; 3663a8a9803SLisandro Dalcin ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr); 3673a8a9803SLisandro Dalcin ierr = PetscMemcpy(t->binterp,binterp,s*p*sizeof(binterp[0]));CHKERRQ(ierr); 3683a8a9803SLisandro Dalcin 369f68a32c8SEmil Constantinescu link->next = RKTableauList; 370f68a32c8SEmil Constantinescu RKTableauList = link; 371f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 372f68a32c8SEmil Constantinescu } 373f68a32c8SEmil Constantinescu 374e4dd521cSBarry Smith /* 375f68a32c8SEmil Constantinescu The step completion formula is 376e4dd521cSBarry Smith 377f68a32c8SEmil Constantinescu x1 = x0 + h b^T YdotRHS 378f68a32c8SEmil Constantinescu 379f68a32c8SEmil Constantinescu This function can be called before or after ts->vec_sol has been updated. 380f68a32c8SEmil Constantinescu Suppose we have a completion formula (b) and an embedded formula (be) of different order. 381f68a32c8SEmil Constantinescu We can write 382f68a32c8SEmil Constantinescu 383f68a32c8SEmil Constantinescu x1e = x0 + h be^T YdotRHS 384f68a32c8SEmil Constantinescu = x1 - h b^T YdotRHS + h be^T YdotRHS 385f68a32c8SEmil Constantinescu = x1 + h (be - b)^T YdotRHS 386f68a32c8SEmil Constantinescu 387f68a32c8SEmil Constantinescu so we can evaluate the method with different order even after the step has been optimistically completed. 388e4dd521cSBarry Smith */ 389f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done) 390f68a32c8SEmil Constantinescu { 391f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 392f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 393f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 394f68a32c8SEmil Constantinescu PetscReal h; 395f68a32c8SEmil Constantinescu PetscInt s = tab->s,j; 396f68a32c8SEmil Constantinescu PetscErrorCode ierr; 397f68a32c8SEmil Constantinescu 398f68a32c8SEmil Constantinescu PetscFunctionBegin; 399f68a32c8SEmil Constantinescu switch (rk->status) { 400f68a32c8SEmil Constantinescu case TS_STEP_INCOMPLETE: 401f68a32c8SEmil Constantinescu case TS_STEP_PENDING: 402f68a32c8SEmil Constantinescu h = ts->time_step; break; 403f68a32c8SEmil Constantinescu case TS_STEP_COMPLETE: 404be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 405f68a32c8SEmil Constantinescu default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 406f68a32c8SEmil Constantinescu } 407f68a32c8SEmil Constantinescu if (order == tab->order) { 408f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { 409f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 410f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->b[j]; 411f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 412f68a32c8SEmil Constantinescu } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 413f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 414f68a32c8SEmil Constantinescu } else if (order == tab->order-1) { 415f68a32c8SEmil Constantinescu if (!tab->bembed) goto unavailable; 416f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (be) */ 417f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 418f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 419f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 420f68a32c8SEmil Constantinescu } else { /* Rollback and re-complete using (be-b) */ 421f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 422f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 423f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 4240f7a1166SHong Zhang if (ts->vec_costintegral && ts->costintegralfwd) { 4250f7a1166SHong Zhang ierr = VecCopy(rk->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr); 4260f7a1166SHong Zhang } 427f68a32c8SEmil Constantinescu } 428f68a32c8SEmil Constantinescu if (done) *done = PETSC_TRUE; 429f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 430f68a32c8SEmil Constantinescu } 431f68a32c8SEmil Constantinescu unavailable: 432f68a32c8SEmil Constantinescu if (done) *done = PETSC_FALSE; 433a7fac7c2SEmil Constantinescu else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"RK '%s' of order %D cannot evaluate step at order %D. Consider using -ts_adapt_type none or a different method that has an embedded estimate.",tab->name,tab->order,order); 434f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 435f68a32c8SEmil Constantinescu } 436f68a32c8SEmil Constantinescu 4370f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts) 4380f7a1166SHong Zhang { 4390f7a1166SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 4400f7a1166SHong Zhang RKTableau tab = rk->tableau; 4410f7a1166SHong Zhang const PetscInt s = tab->s; 4420f7a1166SHong Zhang const PetscReal *b = tab->b,*c = tab->c; 4430f7a1166SHong Zhang Vec *Y = rk->Y; 4440f7a1166SHong Zhang PetscInt i; 4450f7a1166SHong Zhang PetscErrorCode ierr; 4460f7a1166SHong Zhang 4470f7a1166SHong Zhang PetscFunctionBegin; 4480f7a1166SHong Zhang /* backup cost integral */ 4490f7a1166SHong Zhang ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr); 4500f7a1166SHong Zhang for (i=s-1; i>=0; i--) { 4510f7a1166SHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 4520f7a1166SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,rk->ptime+rk->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 4530f7a1166SHong Zhang ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 4540f7a1166SHong Zhang } 4550f7a1166SHong Zhang PetscFunctionReturn(0); 4560f7a1166SHong Zhang } 4570f7a1166SHong Zhang 4580f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts) 4590f7a1166SHong Zhang { 4600f7a1166SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 4610f7a1166SHong Zhang RKTableau tab = rk->tableau; 4620f7a1166SHong Zhang const PetscInt s = tab->s; 4630f7a1166SHong Zhang const PetscReal *b = tab->b,*c = tab->c; 4640f7a1166SHong Zhang Vec *Y = rk->Y; 4650f7a1166SHong Zhang PetscInt i; 4660f7a1166SHong Zhang PetscErrorCode ierr; 4670f7a1166SHong Zhang 4680f7a1166SHong Zhang PetscFunctionBegin; 4690f7a1166SHong Zhang for (i=s-1; i>=0; i--) { 4700f7a1166SHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 4710f7a1166SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime-ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 4720f7a1166SHong Zhang ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 4730f7a1166SHong Zhang } 4740f7a1166SHong Zhang PetscFunctionReturn(0); 4750f7a1166SHong Zhang } 4760f7a1166SHong Zhang 477fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts) 478fecfb714SLisandro Dalcin { 479fecfb714SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 480fecfb714SLisandro Dalcin RKTableau tab = rk->tableau; 481fecfb714SLisandro Dalcin const PetscInt s = tab->s; 482fecfb714SLisandro Dalcin const PetscReal *b = tab->b; 483fecfb714SLisandro Dalcin PetscScalar *w = rk->work; 484fecfb714SLisandro Dalcin Vec *YdotRHS = rk->YdotRHS; 485fecfb714SLisandro Dalcin PetscInt j; 486fecfb714SLisandro Dalcin PetscReal h; 487fecfb714SLisandro Dalcin PetscErrorCode ierr; 488fecfb714SLisandro Dalcin 489fecfb714SLisandro Dalcin PetscFunctionBegin; 490fecfb714SLisandro Dalcin switch (rk->status) { 491fecfb714SLisandro Dalcin case TS_STEP_INCOMPLETE: 492fecfb714SLisandro Dalcin case TS_STEP_PENDING: 493fecfb714SLisandro Dalcin h = ts->time_step; break; 494fecfb714SLisandro Dalcin case TS_STEP_COMPLETE: 495fecfb714SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 496fecfb714SLisandro Dalcin default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 497fecfb714SLisandro Dalcin } 498fecfb714SLisandro Dalcin for (j=0; j<s; j++) w[j] = -h*b[j]; 499fecfb714SLisandro Dalcin ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr); 500fecfb714SLisandro Dalcin PetscFunctionReturn(0); 501fecfb714SLisandro Dalcin } 502fecfb714SLisandro Dalcin 503f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts) 504f68a32c8SEmil Constantinescu { 505f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 506f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 507f68a32c8SEmil Constantinescu const PetscInt s = tab->s; 508fecfb714SLisandro Dalcin const PetscReal *A = tab->A,*c = tab->c; 509f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 510f68a32c8SEmil Constantinescu Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 511d1905564SLisandro Dalcin PetscBool FSAL = tab->FSAL; 512f68a32c8SEmil Constantinescu TSAdapt adapt; 513fecfb714SLisandro Dalcin PetscInt i,j; 514be5899b3SLisandro Dalcin PetscInt rejections = 0; 515be5899b3SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 516be5899b3SLisandro Dalcin PetscReal next_time_step = ts->time_step; 517f68a32c8SEmil Constantinescu PetscErrorCode ierr; 518f68a32c8SEmil Constantinescu 519f68a32c8SEmil Constantinescu PetscFunctionBegin; 520d1905564SLisandro Dalcin if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE; 521d1905564SLisandro Dalcin if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); } 522d1905564SLisandro Dalcin 523f68a32c8SEmil Constantinescu rk->status = TS_STEP_INCOMPLETE; 524be5899b3SLisandro Dalcin while (!ts->reason && rk->status != TS_STEP_COMPLETE) { 525be5899b3SLisandro Dalcin PetscReal t = ts->ptime; 526f68a32c8SEmil Constantinescu PetscReal h = ts->time_step; 527f68a32c8SEmil Constantinescu for (i=0; i<s; i++) { 5289be3e283SDebojyoti Ghosh rk->stage_time = t + h*c[i]; 5299be3e283SDebojyoti Ghosh ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr); 530f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 531f68a32c8SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 532f68a32c8SEmil Constantinescu ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 5339be3e283SDebojyoti Ghosh ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 534f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 535be5899b3SLisandro Dalcin ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr); 536fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 5378f738a7cSLisandro Dalcin if (FSAL && !i) continue; 538f68a32c8SEmil Constantinescu ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 539f68a32c8SEmil Constantinescu } 540be5899b3SLisandro Dalcin 541be5899b3SLisandro Dalcin rk->status = TS_STEP_INCOMPLETE; 542f68a32c8SEmil Constantinescu ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 543f68a32c8SEmil Constantinescu rk->status = TS_STEP_PENDING; 544f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 545f68a32c8SEmil Constantinescu ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 546*1917a363SLisandro Dalcin ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr); 547fecfb714SLisandro Dalcin ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 548be5899b3SLisandro Dalcin rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 549be5899b3SLisandro Dalcin if (!accept) { /* Roll back the current step */ 550fecfb714SLisandro Dalcin ierr = TSRollBack_RK(ts);CHKERRQ(ierr); 551be5899b3SLisandro Dalcin ts->time_step = next_time_step; 552be5899b3SLisandro Dalcin goto reject_step; 553be5899b3SLisandro Dalcin } 554be5899b3SLisandro Dalcin 5550f7a1166SHong Zhang if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation*/ 5560f7a1166SHong Zhang rk->ptime = ts->ptime; 5570f7a1166SHong Zhang rk->time_step = ts->time_step; 558493ed95fSHong Zhang } 559be5899b3SLisandro Dalcin 560f68a32c8SEmil Constantinescu ts->ptime += ts->time_step; 561f68a32c8SEmil Constantinescu ts->time_step = next_time_step; 562f68a32c8SEmil Constantinescu break; 563be5899b3SLisandro Dalcin 564be5899b3SLisandro Dalcin reject_step: 565fecfb714SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 566be5899b3SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 567be5899b3SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 568be5899b3SLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 569f68a32c8SEmil Constantinescu } 570f68a32c8SEmil Constantinescu } 571f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 572e4dd521cSBarry Smith } 573e4dd521cSBarry Smith 574f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts) 575f6a906c0SBarry Smith { 576f6a906c0SBarry Smith TS_RK *rk = (TS_RK*)ts->data; 577be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 578be5899b3SLisandro Dalcin PetscInt s = tab->s; 579f6a906c0SBarry Smith PetscErrorCode ierr; 580f6a906c0SBarry Smith 581f6a906c0SBarry Smith PetscFunctionBegin; 582f6a906c0SBarry Smith if (ts->adjointsetupcalled++) PetscFunctionReturn(0); 583abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr); 584f6a906c0SBarry Smith if(ts->vecs_sensip) { 585abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr); 586f6a906c0SBarry Smith } 587abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr); 588f6a906c0SBarry Smith PetscFunctionReturn(0); 589f6a906c0SBarry Smith } 590f6a906c0SBarry Smith 59142f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts) 592d2daff3dSHong Zhang { 593c235aa8dSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 594c235aa8dSHong Zhang RKTableau tab = rk->tableau; 595c235aa8dSHong Zhang const PetscInt s = tab->s; 596c235aa8dSHong Zhang const PetscReal *A = tab->A,*b = tab->b,*c = tab->c; 597c235aa8dSHong Zhang PetscScalar *w = rk->work; 598ad8e2604SHong Zhang Vec *Y = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu,*VecSensiTemp = rk->VecSensiTemp; 599b18ea86cSHong Zhang PetscInt i,j,nadj; 6003d3eaba7SBarry Smith PetscReal t = ts->ptime; 601d2daff3dSHong Zhang PetscErrorCode ierr; 602cef76868SBarry Smith PetscReal h = ts->time_step; 603cef76868SBarry Smith Mat J,Jp; 604c235aa8dSHong Zhang 605d2daff3dSHong Zhang PetscFunctionBegin; 606c235aa8dSHong Zhang rk->status = TS_STEP_INCOMPLETE; 607c235aa8dSHong Zhang for (i=s-1; i>=0; i--) { 608c235aa8dSHong Zhang rk->stage_time = t + h*(1.0-c[i]); 609abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 610b18ea86cSHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr); 611302440fdSBarry Smith ierr = VecScale(VecSensiTemp[nadj],-h*b[i]);CHKERRQ(ierr); 612c235aa8dSHong Zhang for (j=i+1; j<s; j++) { 613302440fdSBarry Smith ierr = VecAXPY(VecSensiTemp[nadj],-h*A[j*s+i],VecDeltaLam[nadj*s+j]);CHKERRQ(ierr); 614b18ea86cSHong Zhang } 615c235aa8dSHong Zhang } 616ad8e2604SHong Zhang /* Stage values of lambda */ 617c235aa8dSHong Zhang ierr = TSGetRHSJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 618c235aa8dSHong Zhang ierr = TSComputeRHSJacobian(ts,rk->stage_time,Y[i],J,Jp);CHKERRQ(ierr); 619493ed95fSHong Zhang if (ts->vec_costintegral) { 620493ed95fSHong Zhang ierr = TSAdjointComputeDRDYFunction(ts,rk->stage_time,Y[i],ts->vecs_drdy);CHKERRQ(ierr); 621493ed95fSHong Zhang } 622abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 623b18ea86cSHong Zhang ierr = MatMultTranspose(J,VecSensiTemp[nadj],VecDeltaLam[nadj*s+i]);CHKERRQ(ierr); 624493ed95fSHong Zhang if (ts->vec_costintegral) { 625493ed95fSHong Zhang ierr = VecAXPY(VecDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdy[nadj]);CHKERRQ(ierr); 626493ed95fSHong Zhang } 627b18ea86cSHong Zhang } 6286497ce14SHong Zhang 629ad8e2604SHong Zhang /* Stage values of mu */ 6306497ce14SHong Zhang if(ts->vecs_sensip) { 6315bf8c567SBarry Smith ierr = TSAdjointComputeRHSJacobian(ts,rk->stage_time,Y[i],ts->Jacp);CHKERRQ(ierr); 632493ed95fSHong Zhang if (ts->vec_costintegral) { 633493ed95fSHong Zhang ierr = TSAdjointComputeDRDPFunction(ts,rk->stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr); 634493ed95fSHong Zhang } 635493ed95fSHong Zhang 636abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 637ad8e2604SHong Zhang ierr = MatMultTranspose(ts->Jacp,VecSensiTemp[nadj],VecDeltaMu[nadj*s+i]);CHKERRQ(ierr); 638493ed95fSHong Zhang if (ts->vec_costintegral) { 639493ed95fSHong Zhang ierr = VecAXPY(VecDeltaMu[nadj*s+i],-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr); 640493ed95fSHong Zhang } 641ad8e2604SHong Zhang } 642c235aa8dSHong Zhang } 6436497ce14SHong Zhang } 644c235aa8dSHong Zhang 645c235aa8dSHong Zhang for (j=0; j<s; j++) w[j] = 1.0; 646abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 647b18ea86cSHong Zhang ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr); 6486497ce14SHong Zhang if(ts->vecs_sensip) { 649ad8e2604SHong Zhang ierr = VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);CHKERRQ(ierr); 650b18ea86cSHong Zhang } 6516497ce14SHong Zhang } 652c235aa8dSHong Zhang rk->status = TS_STEP_COMPLETE; 653d2daff3dSHong Zhang PetscFunctionReturn(0); 654d2daff3dSHong Zhang } 655d2daff3dSHong Zhang 656f68a32c8SEmil Constantinescu static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X) 657f68a32c8SEmil Constantinescu { 658f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 6593a8a9803SLisandro Dalcin PetscInt s = rk->tableau->s,p = rk->tableau->p,i,j; 660f68a32c8SEmil Constantinescu PetscReal h; 661f68a32c8SEmil Constantinescu PetscReal tt,t; 662f68a32c8SEmil Constantinescu PetscScalar *b; 663f68a32c8SEmil Constantinescu const PetscReal *B = rk->tableau->binterp; 664f68a32c8SEmil Constantinescu PetscErrorCode ierr; 665e4dd521cSBarry Smith 666f68a32c8SEmil Constantinescu PetscFunctionBegin; 667f68a32c8SEmil Constantinescu if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 668be5899b3SLisandro Dalcin 669f68a32c8SEmil Constantinescu switch (rk->status) { 670f68a32c8SEmil Constantinescu case TS_STEP_INCOMPLETE: 671f68a32c8SEmil Constantinescu case TS_STEP_PENDING: 672f68a32c8SEmil Constantinescu h = ts->time_step; 673f68a32c8SEmil Constantinescu t = (itime - ts->ptime)/h; 674f68a32c8SEmil Constantinescu break; 675f68a32c8SEmil Constantinescu case TS_STEP_COMPLETE: 676be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; 677f68a32c8SEmil Constantinescu t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 678f68a32c8SEmil Constantinescu break; 679f68a32c8SEmil Constantinescu default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 680e4dd521cSBarry Smith } 681785e854fSJed Brown ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 682f68a32c8SEmil Constantinescu for (i=0; i<s; i++) b[i] = 0; 6833a8a9803SLisandro Dalcin for (j=0,tt=t; j<p; j++,tt*=t) { 684f68a32c8SEmil Constantinescu for (i=0; i<s; i++) { 6853a8a9803SLisandro Dalcin b[i] += h * B[i*p+j] * tt; 686e4dd521cSBarry Smith } 687f68a32c8SEmil Constantinescu } 688be5899b3SLisandro Dalcin 689f68a32c8SEmil Constantinescu ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr); 690f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr); 691be5899b3SLisandro Dalcin 692f68a32c8SEmil Constantinescu ierr = PetscFree(b);CHKERRQ(ierr); 693e4dd521cSBarry Smith PetscFunctionReturn(0); 694e4dd521cSBarry Smith } 695e4dd521cSBarry Smith 696e4dd521cSBarry Smith /*------------------------------------------------------------*/ 697be5899b3SLisandro Dalcin 698be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts) 699be5899b3SLisandro Dalcin { 700be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 701be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 702be5899b3SLisandro Dalcin PetscErrorCode ierr; 703be5899b3SLisandro Dalcin 704be5899b3SLisandro Dalcin PetscFunctionBegin; 705be5899b3SLisandro Dalcin if (!tab) PetscFunctionReturn(0); 706be5899b3SLisandro Dalcin ierr = PetscFree(rk->work);CHKERRQ(ierr); 707be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr); 708be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr); 709be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr); 710be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr); 711be5899b3SLisandro Dalcin PetscFunctionReturn(0); 712be5899b3SLisandro Dalcin } 713be5899b3SLisandro Dalcin 714277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts) 715e4dd521cSBarry Smith { 7165f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 7176849ba73SBarry Smith PetscErrorCode ierr; 718e4dd521cSBarry Smith 719e4dd521cSBarry Smith PetscFunctionBegin; 720be5899b3SLisandro Dalcin ierr = TSRKTableauReset(ts);CHKERRQ(ierr); 7210f7a1166SHong Zhang ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr); 722abc2977eSBarry Smith ierr = VecDestroyVecs(ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr); 723277b19d0SLisandro Dalcin PetscFunctionReturn(0); 724e4dd521cSBarry Smith } 725277b19d0SLisandro Dalcin 726277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_RK(TS ts) 727277b19d0SLisandro Dalcin { 728277b19d0SLisandro Dalcin PetscErrorCode ierr; 729277b19d0SLisandro Dalcin 730277b19d0SLisandro Dalcin PetscFunctionBegin; 731277b19d0SLisandro Dalcin ierr = TSReset_RK(ts);CHKERRQ(ierr); 732277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 733f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr); 734f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr); 735f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 736f68a32c8SEmil Constantinescu } 737f68a32c8SEmil Constantinescu 738f68a32c8SEmil Constantinescu 739f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx) 740f68a32c8SEmil Constantinescu { 741f68a32c8SEmil Constantinescu PetscFunctionBegin; 742f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 743f68a32c8SEmil Constantinescu } 744f68a32c8SEmil Constantinescu 745f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 746f68a32c8SEmil Constantinescu { 747f68a32c8SEmil Constantinescu PetscFunctionBegin; 748f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 749f68a32c8SEmil Constantinescu } 750f68a32c8SEmil Constantinescu 751f68a32c8SEmil Constantinescu 752f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx) 753f68a32c8SEmil Constantinescu { 754f68a32c8SEmil Constantinescu PetscFunctionBegin; 755f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 756f68a32c8SEmil Constantinescu } 757f68a32c8SEmil Constantinescu 758f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 759f68a32c8SEmil Constantinescu { 760f68a32c8SEmil Constantinescu 761f68a32c8SEmil Constantinescu PetscFunctionBegin; 762f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 763f68a32c8SEmil Constantinescu } 764c235aa8dSHong Zhang /* 765d2daff3dSHong Zhang static PetscErrorCode RKSetAdjCoe(RKTableau tab) 766d2daff3dSHong Zhang { 767d2daff3dSHong Zhang PetscReal *A,*b; 768d2daff3dSHong Zhang PetscInt s,i,j; 769d2daff3dSHong Zhang PetscErrorCode ierr; 770d2daff3dSHong Zhang 771d2daff3dSHong Zhang PetscFunctionBegin; 772d2daff3dSHong Zhang s = tab->s; 773d2daff3dSHong Zhang ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr); 774d2daff3dSHong Zhang 775d2daff3dSHong Zhang for (i=0; i<s; i++) 776d2daff3dSHong Zhang for (j=0; j<s; j++) { 777d2daff3dSHong 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]; 778d2daff3dSHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr); 779d2daff3dSHong Zhang } 780d2daff3dSHong Zhang 781d2daff3dSHong Zhang for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i]; 782d2daff3dSHong Zhang 783d2daff3dSHong Zhang ierr = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 784d2daff3dSHong Zhang ierr = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr); 785d2daff3dSHong Zhang ierr = PetscFree2(A,b);CHKERRQ(ierr); 786d2daff3dSHong Zhang PetscFunctionReturn(0); 787d2daff3dSHong Zhang } 788c235aa8dSHong Zhang */ 789be5899b3SLisandro Dalcin 790be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts) 791be5899b3SLisandro Dalcin { 792be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 793be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 794be5899b3SLisandro Dalcin PetscErrorCode ierr; 795be5899b3SLisandro Dalcin 796be5899b3SLisandro Dalcin PetscFunctionBegin; 797be5899b3SLisandro Dalcin ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr); 798be5899b3SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr); 799be5899b3SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr); 800be5899b3SLisandro Dalcin PetscFunctionReturn(0); 801be5899b3SLisandro Dalcin } 802be5899b3SLisandro Dalcin 803be5899b3SLisandro Dalcin 804f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts) 805f68a32c8SEmil Constantinescu { 806f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 807f68a32c8SEmil Constantinescu PetscErrorCode ierr; 808f68a32c8SEmil Constantinescu DM dm; 809f68a32c8SEmil Constantinescu 810f68a32c8SEmil Constantinescu PetscFunctionBegin; 8113dd83b38SBarry Smith ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr); 812be5899b3SLisandro Dalcin ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr); 8130f7a1166SHong Zhang if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 8140f7a1166SHong Zhang ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr); 8150f7a1166SHong Zhang } 816f68a32c8SEmil Constantinescu ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 817f68a32c8SEmil Constantinescu if (dm) { 818f68a32c8SEmil Constantinescu ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 819f68a32c8SEmil Constantinescu ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 820f68a32c8SEmil Constantinescu } 821e4dd521cSBarry Smith PetscFunctionReturn(0); 822e4dd521cSBarry Smith } 823d2daff3dSHong Zhang 824f6a906c0SBarry Smith 825e4dd521cSBarry Smith /*------------------------------------------------------------*/ 826e4dd521cSBarry Smith 8274416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts) 828e4dd521cSBarry Smith { 829be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 830dfbe8321SBarry Smith PetscErrorCode ierr; 831262ff7bbSBarry Smith 832e4dd521cSBarry Smith PetscFunctionBegin; 833e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr); 834f68a32c8SEmil Constantinescu { 835f68a32c8SEmil Constantinescu RKTableauLink link; 836f68a32c8SEmil Constantinescu PetscInt count,choice; 837f68a32c8SEmil Constantinescu PetscBool flg; 838f68a32c8SEmil Constantinescu const char **namelist; 839f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) ; 840785e854fSJed Brown ierr = PetscMalloc1(count,&namelist);CHKERRQ(ierr); 841f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 842be5899b3SLisandro Dalcin ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr); 843be5899b3SLisandro Dalcin if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);} 844f68a32c8SEmil Constantinescu ierr = PetscFree(namelist);CHKERRQ(ierr); 845f68a32c8SEmil Constantinescu } 846262ff7bbSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 847e4dd521cSBarry Smith PetscFunctionReturn(0); 848e4dd521cSBarry Smith } 849e4dd521cSBarry Smith 850f68a32c8SEmil Constantinescu static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[]) 851f68a32c8SEmil Constantinescu { 852f68a32c8SEmil Constantinescu PetscErrorCode ierr; 853f68a32c8SEmil Constantinescu PetscInt i; 854f68a32c8SEmil Constantinescu size_t left,count; 855f68a32c8SEmil Constantinescu char *p; 856f68a32c8SEmil Constantinescu 857f68a32c8SEmil Constantinescu PetscFunctionBegin; 858f68a32c8SEmil Constantinescu for (i=0,p=buf,left=len; i<n; i++) { 859f68a32c8SEmil Constantinescu ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr); 860f68a32c8SEmil Constantinescu if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer"); 861f68a32c8SEmil Constantinescu left -= count; 862f68a32c8SEmil Constantinescu p += count; 863f68a32c8SEmil Constantinescu *p++ = ' '; 864f68a32c8SEmil Constantinescu } 865f68a32c8SEmil Constantinescu p[i ? 0 : -1] = 0; 866f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 867f68a32c8SEmil Constantinescu } 868f68a32c8SEmil Constantinescu 8695f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer) 870e4dd521cSBarry Smith { 8715f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 8728370ee3bSLisandro Dalcin PetscBool iascii; 873dfbe8321SBarry Smith PetscErrorCode ierr; 8742cdf8120Sasbjorn 875e4dd521cSBarry Smith PetscFunctionBegin; 876251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 8778370ee3bSLisandro Dalcin if (iascii) { 8789c334d8fSLisandro Dalcin RKTableau tab = rk->tableau; 879f68a32c8SEmil Constantinescu TSRKType rktype; 880f68a32c8SEmil Constantinescu char buf[512]; 881f68a32c8SEmil Constantinescu ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr); 8820c37eb8eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," RK: %s\n",rktype);CHKERRQ(ierr); 8830c37eb8eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," Order: %D\n",tab->order);CHKERRQ(ierr); 8840c37eb8eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr); 885f68a32c8SEmil Constantinescu ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 886f68a32c8SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," Abscissa c = %s\n",buf);CHKERRQ(ierr); 8878370ee3bSLisandro Dalcin } 8889c334d8fSLisandro Dalcin if (ts->adapt) {ierr = TSAdaptView(ts->adapt,viewer);CHKERRQ(ierr);} 889f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 890f68a32c8SEmil Constantinescu } 891f68a32c8SEmil Constantinescu 892f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer) 893f68a32c8SEmil Constantinescu { 894f68a32c8SEmil Constantinescu PetscErrorCode ierr; 8959c334d8fSLisandro Dalcin TSAdapt adapt; 896f68a32c8SEmil Constantinescu 897f68a32c8SEmil Constantinescu PetscFunctionBegin; 8989c334d8fSLisandro Dalcin ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 8999c334d8fSLisandro Dalcin ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr); 900f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 901f68a32c8SEmil Constantinescu } 902f68a32c8SEmil Constantinescu 903f68a32c8SEmil Constantinescu /*@C 904f68a32c8SEmil Constantinescu TSRKSetType - Set the type of RK scheme 905f68a32c8SEmil Constantinescu 906f68a32c8SEmil Constantinescu Logically collective 907f68a32c8SEmil Constantinescu 908f68a32c8SEmil Constantinescu Input Parameter: 909f68a32c8SEmil Constantinescu + ts - timestepping context 910f68a32c8SEmil Constantinescu - rktype - type of RK-scheme 911f68a32c8SEmil Constantinescu 912f68a32c8SEmil Constantinescu Level: intermediate 913f68a32c8SEmil Constantinescu 914f68a32c8SEmil Constantinescu .seealso: TSRKGetType(), TSRK, TSRK2, TSRK3, TSRKPRSSP2, TSRK4, TSRK5 915f68a32c8SEmil Constantinescu @*/ 916f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype) 917f68a32c8SEmil Constantinescu { 918f68a32c8SEmil Constantinescu PetscErrorCode ierr; 919f68a32c8SEmil Constantinescu 920f68a32c8SEmil Constantinescu PetscFunctionBegin; 921f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 922f68a32c8SEmil Constantinescu ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr); 923f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 924f68a32c8SEmil Constantinescu } 925f68a32c8SEmil Constantinescu 926f68a32c8SEmil Constantinescu /*@C 927f68a32c8SEmil Constantinescu TSRKGetType - Get the type of RK scheme 928f68a32c8SEmil Constantinescu 929f68a32c8SEmil Constantinescu Logically collective 930f68a32c8SEmil Constantinescu 931f68a32c8SEmil Constantinescu Input Parameter: 932f68a32c8SEmil Constantinescu . ts - timestepping context 933f68a32c8SEmil Constantinescu 934f68a32c8SEmil Constantinescu Output Parameter: 935f68a32c8SEmil Constantinescu . rktype - type of RK-scheme 936f68a32c8SEmil Constantinescu 937f68a32c8SEmil Constantinescu Level: intermediate 938f68a32c8SEmil Constantinescu 939f68a32c8SEmil Constantinescu .seealso: TSRKGetType() 940f68a32c8SEmil Constantinescu @*/ 941f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype) 942f68a32c8SEmil Constantinescu { 943f68a32c8SEmil Constantinescu PetscErrorCode ierr; 944f68a32c8SEmil Constantinescu 945f68a32c8SEmil Constantinescu PetscFunctionBegin; 946f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 947f68a32c8SEmil Constantinescu ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr); 948f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 949f68a32c8SEmil Constantinescu } 950f68a32c8SEmil Constantinescu 951560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype) 952f68a32c8SEmil Constantinescu { 953f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 954f68a32c8SEmil Constantinescu 955f68a32c8SEmil Constantinescu PetscFunctionBegin; 956f68a32c8SEmil Constantinescu *rktype = rk->tableau->name; 957f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 958f68a32c8SEmil Constantinescu } 959560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype) 960f68a32c8SEmil Constantinescu { 961f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 962f68a32c8SEmil Constantinescu PetscErrorCode ierr; 963f68a32c8SEmil Constantinescu PetscBool match; 964f68a32c8SEmil Constantinescu RKTableauLink link; 965f68a32c8SEmil Constantinescu 966f68a32c8SEmil Constantinescu PetscFunctionBegin; 967f68a32c8SEmil Constantinescu if (rk->tableau) { 968f68a32c8SEmil Constantinescu ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr); 969f68a32c8SEmil Constantinescu if (match) PetscFunctionReturn(0); 970f68a32c8SEmil Constantinescu } 971f68a32c8SEmil Constantinescu for (link = RKTableauList; link; link=link->next) { 972f68a32c8SEmil Constantinescu ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr); 973f68a32c8SEmil Constantinescu if (match) { 974be5899b3SLisandro Dalcin if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);} 975f68a32c8SEmil Constantinescu rk->tableau = &link->tab; 976be5899b3SLisandro Dalcin if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);} 977f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 978f68a32c8SEmil Constantinescu } 979f68a32c8SEmil Constantinescu } 980f68a32c8SEmil Constantinescu SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype); 981e4dd521cSBarry Smith PetscFunctionReturn(0); 982e4dd521cSBarry Smith } 983e4dd521cSBarry Smith 984ff22ae23SHong Zhang static PetscErrorCode TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y) 985ff22ae23SHong Zhang { 986ff22ae23SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 987ff22ae23SHong Zhang 988ff22ae23SHong Zhang PetscFunctionBegin; 989ff22ae23SHong Zhang *ns = rk->tableau->s; 990d2daff3dSHong Zhang if(Y) *Y = rk->Y; 991ff22ae23SHong Zhang PetscFunctionReturn(0); 992ff22ae23SHong Zhang } 993ff22ae23SHong Zhang 994ff22ae23SHong Zhang 995e4dd521cSBarry Smith /* ------------------------------------------------------------ */ 996ebe3b25bSBarry Smith /*MC 997f68a32c8SEmil Constantinescu TSRK - ODE and DAE solver using Runge-Kutta schemes 998ebe3b25bSBarry Smith 999f68a32c8SEmil Constantinescu The user should provide the right hand side of the equation 1000f68a32c8SEmil Constantinescu using TSSetRHSFunction(). 1001ebe3b25bSBarry Smith 1002f68a32c8SEmil Constantinescu Notes: 100398164e13SLisandro Dalcin The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type 1004ebe3b25bSBarry Smith 1005d41222bbSBarry Smith Level: beginner 1006d41222bbSBarry Smith 1007f68a32c8SEmil Constantinescu .seealso: TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3, 1008f68a32c8SEmil Constantinescu TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister() 1009ebe3b25bSBarry Smith 1010ebe3b25bSBarry Smith M*/ 10118cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts) 1012e4dd521cSBarry Smith { 10131566a47fSLisandro Dalcin TS_RK *rk; 1014dfbe8321SBarry Smith PetscErrorCode ierr; 1015e4dd521cSBarry Smith 1016e4dd521cSBarry Smith PetscFunctionBegin; 1017f68a32c8SEmil Constantinescu ierr = TSRKInitializePackage();CHKERRQ(ierr); 1018f68a32c8SEmil Constantinescu 1019f68a32c8SEmil Constantinescu ts->ops->reset = TSReset_RK; 10205f70b478SJed Brown ts->ops->destroy = TSDestroy_RK; 10215f70b478SJed Brown ts->ops->view = TSView_RK; 1022f68a32c8SEmil Constantinescu ts->ops->load = TSLoad_RK; 1023f68a32c8SEmil Constantinescu ts->ops->setup = TSSetUp_RK; 102442f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_RK; 1025f68a32c8SEmil Constantinescu ts->ops->step = TSStep_RK; 1026f68a32c8SEmil Constantinescu ts->ops->interpolate = TSInterpolate_RK; 1027f68a32c8SEmil Constantinescu ts->ops->evaluatestep = TSEvaluateStep_RK; 1028fecfb714SLisandro Dalcin ts->ops->rollback = TSRollBack_RK; 1029f68a32c8SEmil Constantinescu ts->ops->setfromoptions = TSSetFromOptions_RK; 1030ff22ae23SHong Zhang ts->ops->getstages = TSGetStages_RK; 103142f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_RK; 1032e4dd521cSBarry Smith 10330f7a1166SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_RK; 10340f7a1166SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_RK; 10350f7a1166SHong Zhang 10361566a47fSLisandro Dalcin ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr); 10371566a47fSLisandro Dalcin ts->data = (void*)rk; 1038e4dd521cSBarry Smith 1039f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr); 1040f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr); 1041be5899b3SLisandro Dalcin 1042be5899b3SLisandro Dalcin ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr); 1043e4dd521cSBarry Smith PetscFunctionReturn(0); 1044e4dd521cSBarry Smith } 1045