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 static PetscInt explicit_stage_time_id; 17f68a32c8SEmil Constantinescu 18f68a32c8SEmil Constantinescu typedef struct _RKTableau *RKTableau; 19f68a32c8SEmil Constantinescu struct _RKTableau { 20f68a32c8SEmil Constantinescu char *name; 21d760c35bSDebojyoti Ghosh PetscInt order; /* Classical approximation order of the method i */ 22f68a32c8SEmil Constantinescu PetscInt s; /* Number of stages */ 23d760c35bSDebojyoti Ghosh PetscBool FSAL; /* flag to indicate if tableau is FSAL */ 24f68a32c8SEmil Constantinescu PetscInt pinterp; /* Interpolation order */ 25f68a32c8SEmil Constantinescu PetscReal *A,*b,*c; /* Tableau */ 26f68a32c8SEmil Constantinescu PetscReal *bembed; /* Embedded formula of order one less (order-1) */ 27f68a32c8SEmil Constantinescu PetscReal *binterp; /* Dense output formula */ 28f68a32c8SEmil Constantinescu PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */ 29f68a32c8SEmil Constantinescu }; 30f68a32c8SEmil Constantinescu typedef struct _RKTableauLink *RKTableauLink; 31f68a32c8SEmil Constantinescu struct _RKTableauLink { 32f68a32c8SEmil Constantinescu struct _RKTableau tab; 33f68a32c8SEmil Constantinescu RKTableauLink next; 34f68a32c8SEmil Constantinescu }; 35f68a32c8SEmil Constantinescu static RKTableauLink RKTableauList; 36e4dd521cSBarry Smith 37e4dd521cSBarry Smith typedef struct { 38f68a32c8SEmil Constantinescu RKTableau tableau; 39f68a32c8SEmil Constantinescu Vec *Y; /* States computed during the step */ 40f68a32c8SEmil Constantinescu Vec *YdotRHS; /* Function evaluations for the non-stiff part */ 41ad8e2604SHong Zhang Vec *VecDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage */ 42ad8e2604SHong Zhang Vec *VecDeltaMu; /* Increment of the adjoint sensitivity w.r.t P at stage */ 43b18ea86cSHong Zhang Vec *VecSensiTemp; /* Vector to be timed with Jacobian transpose */ 440f7a1166SHong Zhang Vec VecCostIntegral0; /* backup for roll-backs due to events */ 45f68a32c8SEmil Constantinescu PetscScalar *work; /* Scalar work */ 46f68a32c8SEmil Constantinescu PetscReal stage_time; 47f68a32c8SEmil Constantinescu TSStepStatus status; 480f7a1166SHong Zhang PetscReal ptime; 490f7a1166SHong Zhang PetscReal time_step; 505f70b478SJed Brown } TS_RK; 51e4dd521cSBarry Smith 52f68a32c8SEmil Constantinescu /*MC 53f68a32c8SEmil Constantinescu TSRK1 - First order forward Euler scheme. 54262ff7bbSBarry Smith 55f68a32c8SEmil Constantinescu This method has one stage. 56f68a32c8SEmil Constantinescu 57f68a32c8SEmil Constantinescu Level: advanced 58f68a32c8SEmil Constantinescu 59f68a32c8SEmil Constantinescu .seealso: TSRK 60f68a32c8SEmil Constantinescu M*/ 61f68a32c8SEmil Constantinescu /*MC 622109b73fSDebojyoti Ghosh TSRK2A - Second order RK scheme. 63f68a32c8SEmil Constantinescu 64f68a32c8SEmil Constantinescu This method has two stages. 65f68a32c8SEmil Constantinescu 66f68a32c8SEmil Constantinescu Level: advanced 67f68a32c8SEmil Constantinescu 68f68a32c8SEmil Constantinescu .seealso: TSRK 69f68a32c8SEmil Constantinescu M*/ 70f68a32c8SEmil Constantinescu /*MC 71f68a32c8SEmil Constantinescu TSRK3 - Third order RK scheme. 72f68a32c8SEmil Constantinescu 73f68a32c8SEmil Constantinescu This method has three stages. 74f68a32c8SEmil Constantinescu 75f68a32c8SEmil Constantinescu Level: advanced 76f68a32c8SEmil Constantinescu 77f68a32c8SEmil Constantinescu .seealso: TSRK 78f68a32c8SEmil Constantinescu M*/ 79f68a32c8SEmil Constantinescu /*MC 802109b73fSDebojyoti Ghosh TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method. 812109b73fSDebojyoti Ghosh 822109b73fSDebojyoti Ghosh This method has four stages. 832109b73fSDebojyoti Ghosh 842109b73fSDebojyoti Ghosh Level: advanced 852109b73fSDebojyoti Ghosh 862109b73fSDebojyoti Ghosh .seealso: TSRK 872109b73fSDebojyoti Ghosh M*/ 882109b73fSDebojyoti Ghosh /*MC 89f68a32c8SEmil Constantinescu TSRK4 - Fourth order RK scheme. 90f68a32c8SEmil Constantinescu 912e698f8bSDebojyoti Ghosh This is the classical Runge-Kutta method with four stages. 92f68a32c8SEmil Constantinescu 93f68a32c8SEmil Constantinescu Level: advanced 94f68a32c8SEmil Constantinescu 95f68a32c8SEmil Constantinescu .seealso: TSRK 96f68a32c8SEmil Constantinescu M*/ 97f68a32c8SEmil Constantinescu /*MC 982e698f8bSDebojyoti Ghosh TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method. 99f68a32c8SEmil Constantinescu 100f68a32c8SEmil Constantinescu This method has six stages. 101f68a32c8SEmil Constantinescu 102f68a32c8SEmil Constantinescu Level: advanced 103f68a32c8SEmil Constantinescu 104f68a32c8SEmil Constantinescu .seealso: TSRK 105f68a32c8SEmil Constantinescu M*/ 1062109b73fSDebojyoti Ghosh /*MC 1072e698f8bSDebojyoti Ghosh TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method. 1082109b73fSDebojyoti Ghosh 1092109b73fSDebojyoti Ghosh This method has seven stages. 1102109b73fSDebojyoti Ghosh 1112109b73fSDebojyoti Ghosh Level: advanced 1122109b73fSDebojyoti Ghosh 1132109b73fSDebojyoti Ghosh .seealso: TSRK 1142109b73fSDebojyoti Ghosh M*/ 115262ff7bbSBarry Smith 116262ff7bbSBarry Smith #undef __FUNCT__ 117f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKRegisterAll" 118f68a32c8SEmil Constantinescu /*@C 119f68a32c8SEmil Constantinescu TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK 120262ff7bbSBarry Smith 121f68a32c8SEmil Constantinescu Not Collective, but should be called by all processes which will need the schemes to be registered 122262ff7bbSBarry Smith 123f68a32c8SEmil Constantinescu Level: advanced 124262ff7bbSBarry Smith 125f68a32c8SEmil Constantinescu .keywords: TS, TSRK, register, all 126262ff7bbSBarry Smith 127f68a32c8SEmil Constantinescu .seealso: TSRKRegisterDestroy() 128262ff7bbSBarry Smith @*/ 129f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void) 130262ff7bbSBarry Smith { 1314ac538c5SBarry Smith PetscErrorCode ierr; 132262ff7bbSBarry Smith 133262ff7bbSBarry Smith PetscFunctionBegin; 134f68a32c8SEmil Constantinescu if (TSRKRegisterAllCalled) PetscFunctionReturn(0); 135f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_TRUE; 136e4dd521cSBarry Smith 137e4dd521cSBarry Smith { 138f68a32c8SEmil Constantinescu const PetscReal 139f68a32c8SEmil Constantinescu A[1][1] = {{0.0}}, 140f68a32c8SEmil Constantinescu b[1] = {1.0}; 141f68a32c8SEmil Constantinescu ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,1,b);CHKERRQ(ierr); 142e4dd521cSBarry Smith } 143db046809SBarry Smith { 144f68a32c8SEmil Constantinescu const PetscReal 145f68a32c8SEmil Constantinescu A[2][2] = {{0.0,0.0}, 146f68a32c8SEmil Constantinescu {1.0,0.0}}, 147f68a32c8SEmil Constantinescu b[2] = {0.5,0.5}, 148f68a32c8SEmil Constantinescu bembed[2] = {1.0,0}; 149fdefd5e5SDebojyoti Ghosh ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,2,b);CHKERRQ(ierr); 150db046809SBarry Smith } 151f68a32c8SEmil Constantinescu { 152f68a32c8SEmil Constantinescu const PetscReal 153f68a32c8SEmil Constantinescu A[3][3] = {{0,0,0}, 154f68a32c8SEmil Constantinescu {2.0/3.0,0,0}, 155f68a32c8SEmil Constantinescu {-1.0/3.0,1.0,0}}, 156f68a32c8SEmil Constantinescu b[3] = {0.25,0.5,0.25}; 157fdefd5e5SDebojyoti Ghosh ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,3,b);CHKERRQ(ierr); 158fdefd5e5SDebojyoti Ghosh } 159fdefd5e5SDebojyoti Ghosh { 160fdefd5e5SDebojyoti Ghosh const PetscReal 161fdefd5e5SDebojyoti Ghosh A[4][4] = {{0,0,0,0}, 162fdefd5e5SDebojyoti Ghosh {1.0/2.0,0,0,0}, 163fdefd5e5SDebojyoti Ghosh {0,3.0/4.0,0,0}, 164fdefd5e5SDebojyoti Ghosh {2.0/9.0,1.0/3.0,4.0/9.0,0}}, 165fdefd5e5SDebojyoti Ghosh b[4] = {2.0/9.0,1.0/3.0,4.0/9.0,0}, 166fdefd5e5SDebojyoti Ghosh bembed[4] = {7.0/24.0,1.0/4.0,1.0/3.0,1.0/8.0}; 167fdefd5e5SDebojyoti Ghosh ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,3,b);CHKERRQ(ierr); 168db046809SBarry Smith } 169f68a32c8SEmil Constantinescu { 170f68a32c8SEmil Constantinescu const PetscReal 171f68a32c8SEmil Constantinescu A[4][4] = {{0,0,0,0}, 172f68a32c8SEmil Constantinescu {0.5,0,0,0}, 173f68a32c8SEmil Constantinescu {0,0.5,0,0}, 174f68a32c8SEmil Constantinescu {0,0,1.0,0}}, 175f68a32c8SEmil Constantinescu b[4] = {1.0/6.0,1.0/3.0,1.0/3.0,1.0/6.0}; 176fdefd5e5SDebojyoti Ghosh ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,4,b);CHKERRQ(ierr); 177db046809SBarry Smith } 178f68a32c8SEmil Constantinescu { 179f68a32c8SEmil Constantinescu const PetscReal 180f68a32c8SEmil Constantinescu A[6][6] = {{0,0,0,0,0,0}, 181f68a32c8SEmil Constantinescu {0.25,0,0,0,0,0}, 182f68a32c8SEmil Constantinescu {3.0/32.0,9.0/32.0,0,0,0,0}, 183f68a32c8SEmil Constantinescu {1932.0/2197.0,-7200.0/2197.0,7296.0/2197.0,0,0,0}, 184f68a32c8SEmil Constantinescu {439.0/216.0,-8.0,3680.0/513.0,-845.0/4104.0,0,0}, 185f68a32c8SEmil Constantinescu {-8.0/27.0,2.0,-3544.0/2565.0,1859.0/4104.0,-11.0/40.0,0}}, 186f68a32c8SEmil Constantinescu b[6] = {16.0/135.0,0,6656.0/12825.0,28561.0/56430.0,-9.0/50.0,2.0/55.0}, 187f68a32c8SEmil Constantinescu bembed[6] = {25.0/216.0,0,1408.0/2565.0,2197.0/4104.0,-1.0/5.0,0}; 188fdefd5e5SDebojyoti Ghosh ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,5,b);CHKERRQ(ierr); 189fdefd5e5SDebojyoti Ghosh } 190fdefd5e5SDebojyoti Ghosh { 191fdefd5e5SDebojyoti Ghosh const PetscReal 192fdefd5e5SDebojyoti Ghosh A[7][7] = {{0,0,0,0,0,0,0}, 193fdefd5e5SDebojyoti Ghosh {1.0/5.0,0,0,0,0,0,0}, 194fdefd5e5SDebojyoti Ghosh {3.0/40.0,9.0/40.0,0,0,0,0,0}, 195fdefd5e5SDebojyoti Ghosh {44.0/45.0,-56.0/15.0,32.0/9.0,0,0,0,0}, 196fdefd5e5SDebojyoti Ghosh {19372.0/6561.0,-25360.0/2187.0,64448.0/6561.0,-212.0/729.0,0,0,0}, 197fdefd5e5SDebojyoti Ghosh {9017.0/3168.0,-355.0/33.0,46732.0/5247.0,49.0/176.0,-5103.0/18656.0,0,0}, 198fdefd5e5SDebojyoti Ghosh {35.0/384.0,0,500.0/1113.0,125.0/192.0,-2187.0/6784.0,11.0/84.0,0}}, 199fdefd5e5SDebojyoti 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}, 200fdefd5e5SDebojyoti 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}; 201fdefd5e5SDebojyoti Ghosh ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,5,b);CHKERRQ(ierr); 202f68a32c8SEmil Constantinescu } 203db046809SBarry Smith PetscFunctionReturn(0); 204db046809SBarry Smith } 205db046809SBarry Smith 206e4dd521cSBarry Smith #undef __FUNCT__ 207f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKRegisterDestroy" 208f68a32c8SEmil Constantinescu /*@C 209f68a32c8SEmil Constantinescu TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister(). 210f68a32c8SEmil Constantinescu 211f68a32c8SEmil Constantinescu Not Collective 212f68a32c8SEmil Constantinescu 213f68a32c8SEmil Constantinescu Level: advanced 214f68a32c8SEmil Constantinescu 215f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy 216f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll() 217f68a32c8SEmil Constantinescu @*/ 218f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void) 219e4dd521cSBarry Smith { 220f68a32c8SEmil Constantinescu PetscErrorCode ierr; 221f68a32c8SEmil Constantinescu RKTableauLink link; 222f68a32c8SEmil Constantinescu 223f68a32c8SEmil Constantinescu PetscFunctionBegin; 224f68a32c8SEmil Constantinescu while ((link = RKTableauList)) { 225f68a32c8SEmil Constantinescu RKTableau t = &link->tab; 226f68a32c8SEmil Constantinescu RKTableauList = link->next; 227f68a32c8SEmil Constantinescu ierr = PetscFree3(t->A,t->b,t->c); CHKERRQ(ierr); 228f68a32c8SEmil Constantinescu ierr = PetscFree (t->bembed); CHKERRQ(ierr); 229f68a32c8SEmil Constantinescu ierr = PetscFree (t->binterp); CHKERRQ(ierr); 230f68a32c8SEmil Constantinescu ierr = PetscFree (t->name); CHKERRQ(ierr); 231f68a32c8SEmil Constantinescu ierr = PetscFree (link); CHKERRQ(ierr); 232f68a32c8SEmil Constantinescu } 233f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_FALSE; 234f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 235f68a32c8SEmil Constantinescu } 236f68a32c8SEmil Constantinescu 237f68a32c8SEmil Constantinescu #undef __FUNCT__ 238f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKInitializePackage" 239f68a32c8SEmil Constantinescu /*@C 240f68a32c8SEmil Constantinescu TSRKInitializePackage - This function initializes everything in the TSRK package. It is called 241f68a32c8SEmil Constantinescu from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RK() 242f68a32c8SEmil Constantinescu when using static libraries. 243f68a32c8SEmil Constantinescu 244f68a32c8SEmil Constantinescu Level: developer 245f68a32c8SEmil Constantinescu 246f68a32c8SEmil Constantinescu .keywords: TS, TSRK, initialize, package 247f68a32c8SEmil Constantinescu .seealso: PetscInitialize() 248f68a32c8SEmil Constantinescu @*/ 249f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void) 250f68a32c8SEmil Constantinescu { 251186e87acSLisandro Dalcin PetscErrorCode ierr; 252e4dd521cSBarry Smith 253e4dd521cSBarry Smith PetscFunctionBegin; 254f68a32c8SEmil Constantinescu if (TSRKPackageInitialized) PetscFunctionReturn(0); 255f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_TRUE; 256f68a32c8SEmil Constantinescu ierr = TSRKRegisterAll();CHKERRQ(ierr); 257f68a32c8SEmil Constantinescu ierr = PetscObjectComposedDataRegister(&explicit_stage_time_id);CHKERRQ(ierr); 258f68a32c8SEmil Constantinescu ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr); 259f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 260f68a32c8SEmil Constantinescu } 261186e87acSLisandro Dalcin 262f68a32c8SEmil Constantinescu #undef __FUNCT__ 263f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKFinalizePackage" 264f68a32c8SEmil Constantinescu /*@C 265f68a32c8SEmil Constantinescu TSRKFinalizePackage - This function destroys everything in the TSRK package. It is 266f68a32c8SEmil Constantinescu called from PetscFinalize(). 267186e87acSLisandro Dalcin 268f68a32c8SEmil Constantinescu Level: developer 269f68a32c8SEmil Constantinescu 270f68a32c8SEmil Constantinescu .keywords: Petsc, destroy, package 271f68a32c8SEmil Constantinescu .seealso: PetscFinalize() 272f68a32c8SEmil Constantinescu @*/ 273f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void) 274f68a32c8SEmil Constantinescu { 275f68a32c8SEmil Constantinescu PetscErrorCode ierr; 276f68a32c8SEmil Constantinescu 277f68a32c8SEmil Constantinescu PetscFunctionBegin; 278f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_FALSE; 279f68a32c8SEmil Constantinescu ierr = TSRKRegisterDestroy();CHKERRQ(ierr); 280f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 281f68a32c8SEmil Constantinescu } 282f68a32c8SEmil Constantinescu 283f68a32c8SEmil Constantinescu #undef __FUNCT__ 284f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKRegister" 285f68a32c8SEmil Constantinescu /*@C 286f68a32c8SEmil Constantinescu TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 287f68a32c8SEmil Constantinescu 288f68a32c8SEmil Constantinescu Not Collective, but the same schemes should be registered on all processes on which they will be used 289f68a32c8SEmil Constantinescu 290f68a32c8SEmil Constantinescu Input Parameters: 291f68a32c8SEmil Constantinescu + name - identifier for method 292f68a32c8SEmil Constantinescu . order - approximation order of method 293f68a32c8SEmil Constantinescu . s - number of stages, this is the dimension of the matrices below 294f68a32c8SEmil Constantinescu . A - stage coefficients (dimension s*s, row-major) 295f68a32c8SEmil Constantinescu . b - step completion table (dimension s; NULL to use last row of A) 296f68a32c8SEmil Constantinescu . c - abscissa (dimension s; NULL to use row sums of A) 297f68a32c8SEmil Constantinescu . bembed - completion table for embedded method (dimension s; NULL if not available) 298f68a32c8SEmil Constantinescu . pinterp - Order of the interpolation scheme, equal to the number of columns of binterp 299f68a32c8SEmil Constantinescu - binterp - Coefficients of the interpolation formula (dimension s*pinterp; NULL to reuse binterpt) 300f68a32c8SEmil Constantinescu 301f68a32c8SEmil Constantinescu Notes: 302f68a32c8SEmil Constantinescu Several RK methods are provided, this function is only needed to create new methods. 303f68a32c8SEmil Constantinescu 304f68a32c8SEmil Constantinescu Level: advanced 305f68a32c8SEmil Constantinescu 306f68a32c8SEmil Constantinescu .keywords: TS, register 307f68a32c8SEmil Constantinescu 308f68a32c8SEmil Constantinescu .seealso: TSRK 309f68a32c8SEmil Constantinescu @*/ 310f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s, 311f68a32c8SEmil Constantinescu const PetscReal A[],const PetscReal b[],const PetscReal c[], 312f68a32c8SEmil Constantinescu const PetscReal bembed[], 313f68a32c8SEmil Constantinescu PetscInt pinterp,const PetscReal binterp[]) 314f68a32c8SEmil Constantinescu { 315f68a32c8SEmil Constantinescu PetscErrorCode ierr; 316f68a32c8SEmil Constantinescu RKTableauLink link; 317f68a32c8SEmil Constantinescu RKTableau t; 318f68a32c8SEmil Constantinescu PetscInt i,j; 319f68a32c8SEmil Constantinescu 320f68a32c8SEmil Constantinescu PetscFunctionBegin; 321f68a32c8SEmil Constantinescu ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr); 322f68a32c8SEmil Constantinescu ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr); 323f68a32c8SEmil Constantinescu t = &link->tab; 324f68a32c8SEmil Constantinescu ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 325f68a32c8SEmil Constantinescu t->order = order; 326f68a32c8SEmil Constantinescu t->s = s; 327dcca6d9dSJed Brown ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr); 328f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 329f68a32c8SEmil Constantinescu if (b) { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); } 330f68a32c8SEmil Constantinescu else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i]; 331f68a32c8SEmil Constantinescu if (c) { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); } 332f68a32c8SEmil 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]; 333d760c35bSDebojyoti Ghosh t->FSAL = PETSC_TRUE; 334d760c35bSDebojyoti Ghosh for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE; 335f68a32c8SEmil Constantinescu if (bembed) { 336785e854fSJed Brown ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr); 337f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr); 338f68a32c8SEmil Constantinescu } 339f68a32c8SEmil Constantinescu 340f68a32c8SEmil Constantinescu t->pinterp = pinterp; 341785e854fSJed Brown ierr = PetscMalloc1(s*pinterp,&t->binterp);CHKERRQ(ierr); 342f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->binterp,binterp,s*pinterp*sizeof(binterp[0]));CHKERRQ(ierr); 343f68a32c8SEmil Constantinescu link->next = RKTableauList; 344f68a32c8SEmil Constantinescu RKTableauList = link; 345f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 346f68a32c8SEmil Constantinescu } 347f68a32c8SEmil Constantinescu 348f68a32c8SEmil Constantinescu #undef __FUNCT__ 349f68a32c8SEmil Constantinescu #define __FUNCT__ "TSEvaluateStep_RK" 350e4dd521cSBarry Smith /* 351f68a32c8SEmil Constantinescu The step completion formula is 352e4dd521cSBarry Smith 353f68a32c8SEmil Constantinescu x1 = x0 + h b^T YdotRHS 354f68a32c8SEmil Constantinescu 355f68a32c8SEmil Constantinescu This function can be called before or after ts->vec_sol has been updated. 356f68a32c8SEmil Constantinescu Suppose we have a completion formula (b) and an embedded formula (be) of different order. 357f68a32c8SEmil Constantinescu We can write 358f68a32c8SEmil Constantinescu 359f68a32c8SEmil Constantinescu x1e = x0 + h be^T YdotRHS 360f68a32c8SEmil Constantinescu = x1 - h b^T YdotRHS + h be^T YdotRHS 361f68a32c8SEmil Constantinescu = x1 + h (be - b)^T YdotRHS 362f68a32c8SEmil Constantinescu 363f68a32c8SEmil Constantinescu so we can evaluate the method with different order even after the step has been optimistically completed. 364e4dd521cSBarry Smith */ 365f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done) 366f68a32c8SEmil Constantinescu { 367f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 368f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 369f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 370f68a32c8SEmil Constantinescu PetscReal h; 371f68a32c8SEmil Constantinescu PetscInt s = tab->s,j; 372f68a32c8SEmil Constantinescu PetscErrorCode ierr; 373f68a32c8SEmil Constantinescu 374f68a32c8SEmil Constantinescu PetscFunctionBegin; 375f68a32c8SEmil Constantinescu switch (rk->status) { 376f68a32c8SEmil Constantinescu case TS_STEP_INCOMPLETE: 377f68a32c8SEmil Constantinescu case TS_STEP_PENDING: 378f68a32c8SEmil Constantinescu h = ts->time_step; break; 379f68a32c8SEmil Constantinescu case TS_STEP_COMPLETE: 380*be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 381f68a32c8SEmil Constantinescu default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 382f68a32c8SEmil Constantinescu } 383f68a32c8SEmil Constantinescu if (order == tab->order) { 384f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { 385f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 386f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->b[j]; 387f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 388f68a32c8SEmil Constantinescu } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 389f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 390f68a32c8SEmil Constantinescu } else if (order == tab->order-1) { 391f68a32c8SEmil Constantinescu if (!tab->bembed) goto unavailable; 392f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (be) */ 393f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 394f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 395f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 396f68a32c8SEmil Constantinescu } else { /* Rollback and re-complete using (be-b) */ 397f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 398f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 399f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 4000f7a1166SHong Zhang if (ts->vec_costintegral && ts->costintegralfwd) { 4010f7a1166SHong Zhang ierr = VecCopy(rk->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr); 4020f7a1166SHong Zhang } 403f68a32c8SEmil Constantinescu } 404f68a32c8SEmil Constantinescu if (done) *done = PETSC_TRUE; 405f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 406f68a32c8SEmil Constantinescu } 407f68a32c8SEmil Constantinescu unavailable: 408f68a32c8SEmil Constantinescu if (done) *done = PETSC_FALSE; 409a7fac7c2SEmil 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); 410f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 411f68a32c8SEmil Constantinescu } 412f68a32c8SEmil Constantinescu 413f68a32c8SEmil Constantinescu #undef __FUNCT__ 4140f7a1166SHong Zhang #define __FUNCT__ "TSForwardCostIntegral_RK" 4150f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts) 4160f7a1166SHong Zhang { 4170f7a1166SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 4180f7a1166SHong Zhang RKTableau tab = rk->tableau; 4190f7a1166SHong Zhang const PetscInt s = tab->s; 4200f7a1166SHong Zhang const PetscReal *b = tab->b,*c = tab->c; 4210f7a1166SHong Zhang Vec *Y = rk->Y; 4220f7a1166SHong Zhang PetscInt i; 4230f7a1166SHong Zhang PetscErrorCode ierr; 4240f7a1166SHong Zhang 4250f7a1166SHong Zhang PetscFunctionBegin; 4260f7a1166SHong Zhang /* backup cost integral */ 4270f7a1166SHong Zhang ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr); 4280f7a1166SHong Zhang for (i=s-1; i>=0; i--) { 4290f7a1166SHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 4300f7a1166SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,rk->ptime+rk->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 4310f7a1166SHong Zhang ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 4320f7a1166SHong Zhang } 4330f7a1166SHong Zhang PetscFunctionReturn(0); 4340f7a1166SHong Zhang } 4350f7a1166SHong Zhang 4360f7a1166SHong Zhang #undef __FUNCT__ 4370f7a1166SHong Zhang #define __FUNCT__ "TSAdjointCostIntegral_RK" 4380f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts) 4390f7a1166SHong Zhang { 4400f7a1166SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 4410f7a1166SHong Zhang RKTableau tab = rk->tableau; 4420f7a1166SHong Zhang const PetscInt s = tab->s; 4430f7a1166SHong Zhang const PetscReal *b = tab->b,*c = tab->c; 4440f7a1166SHong Zhang Vec *Y = rk->Y; 4450f7a1166SHong Zhang PetscInt i; 4460f7a1166SHong Zhang PetscErrorCode ierr; 4470f7a1166SHong Zhang 4480f7a1166SHong Zhang PetscFunctionBegin; 4490f7a1166SHong Zhang for (i=s-1; i>=0; i--) { 4500f7a1166SHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 4510f7a1166SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime-ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 4520f7a1166SHong Zhang ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 4530f7a1166SHong Zhang } 4540f7a1166SHong Zhang PetscFunctionReturn(0); 4550f7a1166SHong Zhang } 4560f7a1166SHong Zhang 4570f7a1166SHong Zhang #undef __FUNCT__ 458f68a32c8SEmil Constantinescu #define __FUNCT__ "TSStep_RK" 459f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts) 460f68a32c8SEmil Constantinescu { 461f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 462f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 463f68a32c8SEmil Constantinescu const PetscInt s = tab->s; 464f68a32c8SEmil Constantinescu const PetscReal *A = tab->A,*b = tab->b,*c = tab->c; 465f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 466f68a32c8SEmil Constantinescu Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 467f68a32c8SEmil Constantinescu TSAdapt adapt; 468*be5899b3SLisandro Dalcin PetscInt i,j,next_scheme; 469*be5899b3SLisandro Dalcin PetscInt rejections = 0; 470*be5899b3SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 471*be5899b3SLisandro Dalcin PetscReal next_time_step = ts->time_step; 472f68a32c8SEmil Constantinescu PetscErrorCode ierr; 473f68a32c8SEmil Constantinescu 474f68a32c8SEmil Constantinescu PetscFunctionBegin; 475f68a32c8SEmil Constantinescu 476f68a32c8SEmil Constantinescu rk->status = TS_STEP_INCOMPLETE; 477*be5899b3SLisandro Dalcin while (!ts->reason && rk->status != TS_STEP_COMPLETE) { 478*be5899b3SLisandro Dalcin PetscReal t = ts->ptime; 479f68a32c8SEmil Constantinescu PetscReal h = ts->time_step; 480f68a32c8SEmil Constantinescu for (i=0; i<s; i++) { 4819be3e283SDebojyoti Ghosh rk->stage_time = t + h*c[i]; 4829be3e283SDebojyoti Ghosh ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr); 483f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 484f68a32c8SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 485f68a32c8SEmil Constantinescu ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 4869be3e283SDebojyoti Ghosh ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 487f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 488*be5899b3SLisandro Dalcin ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr); 489*be5899b3SLisandro Dalcin if (!stageok) {accept = PETSC_FALSE; goto reject_step;} 490f68a32c8SEmil Constantinescu ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 491f68a32c8SEmil Constantinescu } 492*be5899b3SLisandro Dalcin 493*be5899b3SLisandro Dalcin rk->status = TS_STEP_INCOMPLETE; 494f68a32c8SEmil Constantinescu ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 495f68a32c8SEmil Constantinescu rk->status = TS_STEP_PENDING; 496f68a32c8SEmil Constantinescu /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */ 497f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 498f68a32c8SEmil Constantinescu ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 499f68a32c8SEmil Constantinescu ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr); 500f68a32c8SEmil Constantinescu ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr); 501*be5899b3SLisandro Dalcin rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 502*be5899b3SLisandro Dalcin if (!accept) { /* Roll back the current step */ 503*be5899b3SLisandro Dalcin for (j=0; j<s; j++) w[j] = -h*b[j]; 504*be5899b3SLisandro Dalcin ierr = VecMAXPY(ts->vec_sol,s,w,rk->YdotRHS);CHKERRQ(ierr); 505*be5899b3SLisandro Dalcin ts->time_step = next_time_step; 506*be5899b3SLisandro Dalcin goto reject_step; 507*be5899b3SLisandro Dalcin } 508*be5899b3SLisandro Dalcin 5090f7a1166SHong Zhang if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation*/ 5100f7a1166SHong Zhang rk->ptime = ts->ptime; 5110f7a1166SHong Zhang rk->time_step = ts->time_step; 512493ed95fSHong Zhang } 513*be5899b3SLisandro Dalcin 514*be5899b3SLisandro Dalcin /* Ignore next_scheme for now */ 515f68a32c8SEmil Constantinescu ts->ptime += ts->time_step; 516f68a32c8SEmil Constantinescu ts->time_step = next_time_step; 517f68a32c8SEmil Constantinescu ierr = PetscObjectComposedDataSetReal((PetscObject)ts->vec_sol,explicit_stage_time_id,ts->ptime);CHKERRQ(ierr); 518f68a32c8SEmil Constantinescu break; 519*be5899b3SLisandro Dalcin 520*be5899b3SLisandro Dalcin reject_step: 521*be5899b3SLisandro Dalcin ts->reject++; 522*be5899b3SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 523*be5899b3SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 524*be5899b3SLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 525f68a32c8SEmil Constantinescu } 526f68a32c8SEmil Constantinescu } 527f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 528e4dd521cSBarry Smith } 529e4dd521cSBarry Smith 530f68a32c8SEmil Constantinescu #undef __FUNCT__ 531f6a906c0SBarry Smith #define __FUNCT__ "TSAdjointSetUp_RK" 532f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts) 533f6a906c0SBarry Smith { 534f6a906c0SBarry Smith TS_RK *rk = (TS_RK*)ts->data; 535*be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 536*be5899b3SLisandro Dalcin PetscInt s = tab->s; 537f6a906c0SBarry Smith PetscErrorCode ierr; 538f6a906c0SBarry Smith 539f6a906c0SBarry Smith PetscFunctionBegin; 540f6a906c0SBarry Smith if (ts->adjointsetupcalled++) PetscFunctionReturn(0); 541abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr); 542f6a906c0SBarry Smith if(ts->vecs_sensip) { 543abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr); 544f6a906c0SBarry Smith } 545abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr); 546f6a906c0SBarry Smith PetscFunctionReturn(0); 547f6a906c0SBarry Smith } 548f6a906c0SBarry Smith 549f6a906c0SBarry Smith #undef __FUNCT__ 55042f2b339SBarry Smith #define __FUNCT__ "TSAdjointStep_RK" 55142f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts) 552d2daff3dSHong Zhang { 553c235aa8dSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 554c235aa8dSHong Zhang RKTableau tab = rk->tableau; 555c235aa8dSHong Zhang const PetscInt s = tab->s; 556c235aa8dSHong Zhang const PetscReal *A = tab->A,*b = tab->b,*c = tab->c; 557c235aa8dSHong Zhang PetscScalar *w = rk->work; 558ad8e2604SHong Zhang Vec *Y = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu,*VecSensiTemp = rk->VecSensiTemp; 559b18ea86cSHong Zhang PetscInt i,j,nadj; 5603d3eaba7SBarry Smith PetscReal t = ts->ptime; 561d2daff3dSHong Zhang PetscErrorCode ierr; 562cef76868SBarry Smith PetscReal h = ts->time_step; 563cef76868SBarry Smith Mat J,Jp; 564c235aa8dSHong Zhang 565d2daff3dSHong Zhang PetscFunctionBegin; 566c235aa8dSHong Zhang rk->status = TS_STEP_INCOMPLETE; 567c235aa8dSHong Zhang for (i=s-1; i>=0; i--) { 568c235aa8dSHong Zhang rk->stage_time = t + h*(1.0-c[i]); 569abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 570b18ea86cSHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr); 571302440fdSBarry Smith ierr = VecScale(VecSensiTemp[nadj],-h*b[i]);CHKERRQ(ierr); 572c235aa8dSHong Zhang for (j=i+1; j<s; j++) { 573302440fdSBarry Smith ierr = VecAXPY(VecSensiTemp[nadj],-h*A[j*s+i],VecDeltaLam[nadj*s+j]);CHKERRQ(ierr); 574b18ea86cSHong Zhang } 575c235aa8dSHong Zhang } 576ad8e2604SHong Zhang /* Stage values of lambda */ 577c235aa8dSHong Zhang ierr = TSGetRHSJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 578c235aa8dSHong Zhang ierr = TSComputeRHSJacobian(ts,rk->stage_time,Y[i],J,Jp);CHKERRQ(ierr); 579493ed95fSHong Zhang if (ts->vec_costintegral) { 580493ed95fSHong Zhang ierr = TSAdjointComputeDRDYFunction(ts,rk->stage_time,Y[i],ts->vecs_drdy);CHKERRQ(ierr); 581493ed95fSHong Zhang } 582abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 583b18ea86cSHong Zhang ierr = MatMultTranspose(J,VecSensiTemp[nadj],VecDeltaLam[nadj*s+i]);CHKERRQ(ierr); 584493ed95fSHong Zhang if (ts->vec_costintegral) { 585493ed95fSHong Zhang ierr = VecAXPY(VecDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdy[nadj]);CHKERRQ(ierr); 586493ed95fSHong Zhang } 587b18ea86cSHong Zhang } 5886497ce14SHong Zhang 589ad8e2604SHong Zhang /* Stage values of mu */ 5906497ce14SHong Zhang if(ts->vecs_sensip) { 5915bf8c567SBarry Smith ierr = TSAdjointComputeRHSJacobian(ts,rk->stage_time,Y[i],ts->Jacp);CHKERRQ(ierr); 592493ed95fSHong Zhang if (ts->vec_costintegral) { 593493ed95fSHong Zhang ierr = TSAdjointComputeDRDPFunction(ts,rk->stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr); 594493ed95fSHong Zhang } 595493ed95fSHong Zhang 596abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 597ad8e2604SHong Zhang ierr = MatMultTranspose(ts->Jacp,VecSensiTemp[nadj],VecDeltaMu[nadj*s+i]);CHKERRQ(ierr); 598493ed95fSHong Zhang if (ts->vec_costintegral) { 599493ed95fSHong Zhang ierr = VecAXPY(VecDeltaMu[nadj*s+i],-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr); 600493ed95fSHong Zhang } 601ad8e2604SHong Zhang } 602c235aa8dSHong Zhang } 6036497ce14SHong Zhang } 604c235aa8dSHong Zhang 605c235aa8dSHong Zhang for (j=0; j<s; j++) w[j] = 1.0; 606abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 607b18ea86cSHong Zhang ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr); 6086497ce14SHong Zhang if(ts->vecs_sensip) { 609ad8e2604SHong Zhang ierr = VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);CHKERRQ(ierr); 610b18ea86cSHong Zhang } 6116497ce14SHong Zhang } 612c235aa8dSHong Zhang rk->status = TS_STEP_COMPLETE; 613d2daff3dSHong Zhang PetscFunctionReturn(0); 614d2daff3dSHong Zhang } 615d2daff3dSHong Zhang 616d2daff3dSHong Zhang #undef __FUNCT__ 617f68a32c8SEmil Constantinescu #define __FUNCT__ "TSInterpolate_RK" 618f68a32c8SEmil Constantinescu static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X) 619f68a32c8SEmil Constantinescu { 620f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 621f68a32c8SEmil Constantinescu PetscInt s = rk->tableau->s,pinterp = rk->tableau->pinterp,i,j; 622f68a32c8SEmil Constantinescu PetscReal h; 623f68a32c8SEmil Constantinescu PetscReal tt,t; 624f68a32c8SEmil Constantinescu PetscScalar *b; 625f68a32c8SEmil Constantinescu const PetscReal *B = rk->tableau->binterp; 626f68a32c8SEmil Constantinescu PetscErrorCode ierr; 627e4dd521cSBarry Smith 628f68a32c8SEmil Constantinescu PetscFunctionBegin; 629f68a32c8SEmil Constantinescu if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 630*be5899b3SLisandro Dalcin 631f68a32c8SEmil Constantinescu switch (rk->status) { 632f68a32c8SEmil Constantinescu case TS_STEP_INCOMPLETE: 633f68a32c8SEmil Constantinescu case TS_STEP_PENDING: 634f68a32c8SEmil Constantinescu h = ts->time_step; 635f68a32c8SEmil Constantinescu t = (itime - ts->ptime)/h; 636f68a32c8SEmil Constantinescu break; 637f68a32c8SEmil Constantinescu case TS_STEP_COMPLETE: 638*be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; 639f68a32c8SEmil Constantinescu t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 640f68a32c8SEmil Constantinescu break; 641f68a32c8SEmil Constantinescu default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 642e4dd521cSBarry Smith } 643785e854fSJed Brown ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 644f68a32c8SEmil Constantinescu for (i=0; i<s; i++) b[i] = 0; 645f68a32c8SEmil Constantinescu for (j=0,tt=t; j<pinterp; j++,tt*=t) { 646f68a32c8SEmil Constantinescu for (i=0; i<s; i++) { 647f68a32c8SEmil Constantinescu b[i] += h * B[i*pinterp+j] * tt; 648e4dd521cSBarry Smith } 649f68a32c8SEmil Constantinescu } 650*be5899b3SLisandro Dalcin 651f68a32c8SEmil Constantinescu ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr); 652f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr); 653*be5899b3SLisandro Dalcin 654f68a32c8SEmil Constantinescu ierr = PetscFree(b);CHKERRQ(ierr); 655e4dd521cSBarry Smith PetscFunctionReturn(0); 656e4dd521cSBarry Smith } 657e4dd521cSBarry Smith 658e4dd521cSBarry Smith /*------------------------------------------------------------*/ 659*be5899b3SLisandro Dalcin 660*be5899b3SLisandro Dalcin #undef __FUNCT__ 661*be5899b3SLisandro Dalcin #define __FUNCT__ "TSRKTableauReset" 662*be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts) 663*be5899b3SLisandro Dalcin { 664*be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 665*be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 666*be5899b3SLisandro Dalcin PetscErrorCode ierr; 667*be5899b3SLisandro Dalcin 668*be5899b3SLisandro Dalcin PetscFunctionBegin; 669*be5899b3SLisandro Dalcin if (!tab) PetscFunctionReturn(0); 670*be5899b3SLisandro Dalcin ierr = PetscFree(rk->work);CHKERRQ(ierr); 671*be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr); 672*be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr); 673*be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr); 674*be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr); 675*be5899b3SLisandro Dalcin PetscFunctionReturn(0); 676*be5899b3SLisandro Dalcin } 677*be5899b3SLisandro Dalcin 678e4dd521cSBarry Smith #undef __FUNCT__ 679277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_RK" 680277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts) 681e4dd521cSBarry Smith { 6825f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 6836849ba73SBarry Smith PetscErrorCode ierr; 684e4dd521cSBarry Smith 685e4dd521cSBarry Smith PetscFunctionBegin; 686*be5899b3SLisandro Dalcin ierr = TSRKTableauReset(ts);CHKERRQ(ierr); 6870f7a1166SHong Zhang ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr); 688abc2977eSBarry Smith ierr = VecDestroyVecs(ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr); 689277b19d0SLisandro Dalcin PetscFunctionReturn(0); 690e4dd521cSBarry Smith } 691277b19d0SLisandro Dalcin 692277b19d0SLisandro Dalcin #undef __FUNCT__ 693277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_RK" 694277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_RK(TS ts) 695277b19d0SLisandro Dalcin { 696277b19d0SLisandro Dalcin PetscErrorCode ierr; 697277b19d0SLisandro Dalcin 698277b19d0SLisandro Dalcin PetscFunctionBegin; 699277b19d0SLisandro Dalcin ierr = TSReset_RK(ts);CHKERRQ(ierr); 700277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 701f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr); 702f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr); 703f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 704f68a32c8SEmil Constantinescu } 705f68a32c8SEmil Constantinescu 706f68a32c8SEmil Constantinescu 707f68a32c8SEmil Constantinescu #undef __FUNCT__ 708f68a32c8SEmil Constantinescu #define __FUNCT__ "DMCoarsenHook_TSRK" 709f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx) 710f68a32c8SEmil Constantinescu { 711f68a32c8SEmil Constantinescu PetscFunctionBegin; 712f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 713f68a32c8SEmil Constantinescu } 714f68a32c8SEmil Constantinescu 715f68a32c8SEmil Constantinescu #undef __FUNCT__ 716f68a32c8SEmil Constantinescu #define __FUNCT__ "DMRestrictHook_TSRK" 717f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 718f68a32c8SEmil Constantinescu { 719f68a32c8SEmil Constantinescu PetscFunctionBegin; 720f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 721f68a32c8SEmil Constantinescu } 722f68a32c8SEmil Constantinescu 723f68a32c8SEmil Constantinescu 724f68a32c8SEmil Constantinescu #undef __FUNCT__ 725f68a32c8SEmil Constantinescu #define __FUNCT__ "DMSubDomainHook_TSRK" 726f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx) 727f68a32c8SEmil Constantinescu { 728f68a32c8SEmil Constantinescu PetscFunctionBegin; 729f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 730f68a32c8SEmil Constantinescu } 731f68a32c8SEmil Constantinescu 732f68a32c8SEmil Constantinescu #undef __FUNCT__ 733f68a32c8SEmil Constantinescu #define __FUNCT__ "DMSubDomainRestrictHook_TSRK" 734f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 735f68a32c8SEmil Constantinescu { 736f68a32c8SEmil Constantinescu 737f68a32c8SEmil Constantinescu PetscFunctionBegin; 738f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 739f68a32c8SEmil Constantinescu } 740c235aa8dSHong Zhang /* 741f68a32c8SEmil Constantinescu #undef __FUNCT__ 742d2daff3dSHong Zhang #define __FUNCT__ "RKSetAdjCoe" 743d2daff3dSHong Zhang static PetscErrorCode RKSetAdjCoe(RKTableau tab) 744d2daff3dSHong Zhang { 745d2daff3dSHong Zhang PetscReal *A,*b; 746d2daff3dSHong Zhang PetscInt s,i,j; 747d2daff3dSHong Zhang PetscErrorCode ierr; 748d2daff3dSHong Zhang 749d2daff3dSHong Zhang PetscFunctionBegin; 750d2daff3dSHong Zhang s = tab->s; 751d2daff3dSHong Zhang ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr); 752d2daff3dSHong Zhang 753d2daff3dSHong Zhang for (i=0; i<s; i++) 754d2daff3dSHong Zhang for (j=0; j<s; j++) { 755d2daff3dSHong 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]; 756d2daff3dSHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr); 757d2daff3dSHong Zhang } 758d2daff3dSHong Zhang 759d2daff3dSHong Zhang for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i]; 760d2daff3dSHong Zhang 761d2daff3dSHong Zhang ierr = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 762d2daff3dSHong Zhang ierr = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr); 763d2daff3dSHong Zhang ierr = PetscFree2(A,b);CHKERRQ(ierr); 764d2daff3dSHong Zhang PetscFunctionReturn(0); 765d2daff3dSHong Zhang } 766c235aa8dSHong Zhang */ 767*be5899b3SLisandro Dalcin 768*be5899b3SLisandro Dalcin #undef __FUNCT__ 769*be5899b3SLisandro Dalcin #define __FUNCT__ "TSRKTableauSetUp" 770*be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts) 771*be5899b3SLisandro Dalcin { 772*be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 773*be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 774*be5899b3SLisandro Dalcin PetscErrorCode ierr; 775*be5899b3SLisandro Dalcin 776*be5899b3SLisandro Dalcin PetscFunctionBegin; 777*be5899b3SLisandro Dalcin ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr); 778*be5899b3SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr); 779*be5899b3SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr); 780*be5899b3SLisandro Dalcin PetscFunctionReturn(0); 781*be5899b3SLisandro Dalcin } 782*be5899b3SLisandro Dalcin 783*be5899b3SLisandro Dalcin 784d2daff3dSHong Zhang #undef __FUNCT__ 785f68a32c8SEmil Constantinescu #define __FUNCT__ "TSSetUp_RK" 786f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts) 787f68a32c8SEmil Constantinescu { 788f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 789f68a32c8SEmil Constantinescu PetscErrorCode ierr; 790f68a32c8SEmil Constantinescu DM dm; 791f68a32c8SEmil Constantinescu 792f68a32c8SEmil Constantinescu PetscFunctionBegin; 793*be5899b3SLisandro Dalcin ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr); 7940f7a1166SHong Zhang if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 7950f7a1166SHong Zhang ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr); 7960f7a1166SHong Zhang } 797f68a32c8SEmil Constantinescu ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 798f68a32c8SEmil Constantinescu if (dm) { 799f68a32c8SEmil Constantinescu ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 800f68a32c8SEmil Constantinescu ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 801f68a32c8SEmil Constantinescu } 802e4dd521cSBarry Smith PetscFunctionReturn(0); 803e4dd521cSBarry Smith } 804d2daff3dSHong Zhang 805f6a906c0SBarry Smith 806e4dd521cSBarry Smith /*------------------------------------------------------------*/ 807e4dd521cSBarry Smith 808e4dd521cSBarry Smith #undef __FUNCT__ 8095f70b478SJed Brown #define __FUNCT__ "TSSetFromOptions_RK" 8104416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts) 811e4dd521cSBarry Smith { 812*be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 813dfbe8321SBarry Smith PetscErrorCode ierr; 814262ff7bbSBarry Smith 815e4dd521cSBarry Smith PetscFunctionBegin; 816e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr); 817f68a32c8SEmil Constantinescu { 818f68a32c8SEmil Constantinescu RKTableauLink link; 819f68a32c8SEmil Constantinescu PetscInt count,choice; 820f68a32c8SEmil Constantinescu PetscBool flg; 821f68a32c8SEmil Constantinescu const char **namelist; 822f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) ; 823785e854fSJed Brown ierr = PetscMalloc1(count,&namelist);CHKERRQ(ierr); 824f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 825*be5899b3SLisandro Dalcin ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr); 826*be5899b3SLisandro Dalcin if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);} 827f68a32c8SEmil Constantinescu ierr = PetscFree(namelist);CHKERRQ(ierr); 828f68a32c8SEmil Constantinescu } 829262ff7bbSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 830e4dd521cSBarry Smith PetscFunctionReturn(0); 831e4dd521cSBarry Smith } 832e4dd521cSBarry Smith 833e4dd521cSBarry Smith #undef __FUNCT__ 834f68a32c8SEmil Constantinescu #define __FUNCT__ "PetscFormatRealArray" 835f68a32c8SEmil Constantinescu static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[]) 836f68a32c8SEmil Constantinescu { 837f68a32c8SEmil Constantinescu PetscErrorCode ierr; 838f68a32c8SEmil Constantinescu PetscInt i; 839f68a32c8SEmil Constantinescu size_t left,count; 840f68a32c8SEmil Constantinescu char *p; 841f68a32c8SEmil Constantinescu 842f68a32c8SEmil Constantinescu PetscFunctionBegin; 843f68a32c8SEmil Constantinescu for (i=0,p=buf,left=len; i<n; i++) { 844f68a32c8SEmil Constantinescu ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr); 845f68a32c8SEmil Constantinescu if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer"); 846f68a32c8SEmil Constantinescu left -= count; 847f68a32c8SEmil Constantinescu p += count; 848f68a32c8SEmil Constantinescu *p++ = ' '; 849f68a32c8SEmil Constantinescu } 850f68a32c8SEmil Constantinescu p[i ? 0 : -1] = 0; 851f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 852f68a32c8SEmil Constantinescu } 853f68a32c8SEmil Constantinescu 854f68a32c8SEmil Constantinescu #undef __FUNCT__ 8555f70b478SJed Brown #define __FUNCT__ "TSView_RK" 8565f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer) 857e4dd521cSBarry Smith { 8585f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 8598370ee3bSLisandro Dalcin PetscBool iascii; 860dfbe8321SBarry Smith PetscErrorCode ierr; 8612cdf8120Sasbjorn 862e4dd521cSBarry Smith PetscFunctionBegin; 863251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 8648370ee3bSLisandro Dalcin if (iascii) { 8659c334d8fSLisandro Dalcin RKTableau tab = rk->tableau; 866f68a32c8SEmil Constantinescu TSRKType rktype; 867f68a32c8SEmil Constantinescu char buf[512]; 868f68a32c8SEmil Constantinescu ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr); 869f68a32c8SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," RK %s\n",rktype);CHKERRQ(ierr); 870f68a32c8SEmil Constantinescu ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 871f68a32c8SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," Abscissa c = %s\n",buf);CHKERRQ(ierr); 872d760c35bSDebojyoti Ghosh ierr = PetscViewerASCIIPrintf(viewer,"FSAL: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr); 8738370ee3bSLisandro Dalcin } 8749c334d8fSLisandro Dalcin if (ts->adapt) {ierr = TSAdaptView(ts->adapt,viewer);CHKERRQ(ierr);} 875f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 876f68a32c8SEmil Constantinescu } 877f68a32c8SEmil Constantinescu 878f68a32c8SEmil Constantinescu #undef __FUNCT__ 879f68a32c8SEmil Constantinescu #define __FUNCT__ "TSLoad_RK" 880f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer) 881f68a32c8SEmil Constantinescu { 882f68a32c8SEmil Constantinescu PetscErrorCode ierr; 8839c334d8fSLisandro Dalcin TSAdapt adapt; 884f68a32c8SEmil Constantinescu 885f68a32c8SEmil Constantinescu PetscFunctionBegin; 8869c334d8fSLisandro Dalcin ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 8879c334d8fSLisandro Dalcin ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr); 888f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 889f68a32c8SEmil Constantinescu } 890f68a32c8SEmil Constantinescu 891f68a32c8SEmil Constantinescu #undef __FUNCT__ 892f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKSetType" 893f68a32c8SEmil Constantinescu /*@C 894f68a32c8SEmil Constantinescu TSRKSetType - Set the type of RK scheme 895f68a32c8SEmil Constantinescu 896f68a32c8SEmil Constantinescu Logically collective 897f68a32c8SEmil Constantinescu 898f68a32c8SEmil Constantinescu Input Parameter: 899f68a32c8SEmil Constantinescu + ts - timestepping context 900f68a32c8SEmil Constantinescu - rktype - type of RK-scheme 901f68a32c8SEmil Constantinescu 902f68a32c8SEmil Constantinescu Level: intermediate 903f68a32c8SEmil Constantinescu 904f68a32c8SEmil Constantinescu .seealso: TSRKGetType(), TSRK, TSRK2, TSRK3, TSRKPRSSP2, TSRK4, TSRK5 905f68a32c8SEmil Constantinescu @*/ 906f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype) 907f68a32c8SEmil Constantinescu { 908f68a32c8SEmil Constantinescu PetscErrorCode ierr; 909f68a32c8SEmil Constantinescu 910f68a32c8SEmil Constantinescu PetscFunctionBegin; 911f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 912f68a32c8SEmil Constantinescu ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr); 913f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 914f68a32c8SEmil Constantinescu } 915f68a32c8SEmil Constantinescu 916f68a32c8SEmil Constantinescu #undef __FUNCT__ 917f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKGetType" 918f68a32c8SEmil Constantinescu /*@C 919f68a32c8SEmil Constantinescu TSRKGetType - Get the type of RK scheme 920f68a32c8SEmil Constantinescu 921f68a32c8SEmil Constantinescu Logically collective 922f68a32c8SEmil Constantinescu 923f68a32c8SEmil Constantinescu Input Parameter: 924f68a32c8SEmil Constantinescu . ts - timestepping context 925f68a32c8SEmil Constantinescu 926f68a32c8SEmil Constantinescu Output Parameter: 927f68a32c8SEmil Constantinescu . rktype - type of RK-scheme 928f68a32c8SEmil Constantinescu 929f68a32c8SEmil Constantinescu Level: intermediate 930f68a32c8SEmil Constantinescu 931f68a32c8SEmil Constantinescu .seealso: TSRKGetType() 932f68a32c8SEmil Constantinescu @*/ 933f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype) 934f68a32c8SEmil Constantinescu { 935f68a32c8SEmil Constantinescu PetscErrorCode ierr; 936f68a32c8SEmil Constantinescu 937f68a32c8SEmil Constantinescu PetscFunctionBegin; 938f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 939f68a32c8SEmil Constantinescu ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr); 940f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 941f68a32c8SEmil Constantinescu } 942f68a32c8SEmil Constantinescu 943f68a32c8SEmil Constantinescu #undef __FUNCT__ 944f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKGetType_RK" 945560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype) 946f68a32c8SEmil Constantinescu { 947f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 948f68a32c8SEmil Constantinescu 949f68a32c8SEmil Constantinescu PetscFunctionBegin; 950f68a32c8SEmil Constantinescu *rktype = rk->tableau->name; 951f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 952f68a32c8SEmil Constantinescu } 953f68a32c8SEmil Constantinescu #undef __FUNCT__ 954f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKSetType_RK" 955560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype) 956f68a32c8SEmil Constantinescu { 957f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 958f68a32c8SEmil Constantinescu PetscErrorCode ierr; 959f68a32c8SEmil Constantinescu PetscBool match; 960f68a32c8SEmil Constantinescu RKTableauLink link; 961f68a32c8SEmil Constantinescu 962f68a32c8SEmil Constantinescu PetscFunctionBegin; 963f68a32c8SEmil Constantinescu if (rk->tableau) { 964f68a32c8SEmil Constantinescu ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr); 965f68a32c8SEmil Constantinescu if (match) PetscFunctionReturn(0); 966f68a32c8SEmil Constantinescu } 967f68a32c8SEmil Constantinescu for (link = RKTableauList; link; link=link->next) { 968f68a32c8SEmil Constantinescu ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr); 969f68a32c8SEmil Constantinescu if (match) { 970*be5899b3SLisandro Dalcin if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);} 971f68a32c8SEmil Constantinescu rk->tableau = &link->tab; 972*be5899b3SLisandro Dalcin if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);} 973f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 974f68a32c8SEmil Constantinescu } 975f68a32c8SEmil Constantinescu } 976f68a32c8SEmil Constantinescu SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype); 977e4dd521cSBarry Smith PetscFunctionReturn(0); 978e4dd521cSBarry Smith } 979e4dd521cSBarry Smith 980ff22ae23SHong Zhang #undef __FUNCT__ 981ff22ae23SHong Zhang #define __FUNCT__ "TSGetStages_RK" 982ff22ae23SHong Zhang static PetscErrorCode TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y) 983ff22ae23SHong Zhang { 984ff22ae23SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 985ff22ae23SHong Zhang 986ff22ae23SHong Zhang PetscFunctionBegin; 987ff22ae23SHong Zhang *ns = rk->tableau->s; 988d2daff3dSHong Zhang if(Y) *Y = rk->Y; 989ff22ae23SHong Zhang PetscFunctionReturn(0); 990ff22ae23SHong Zhang } 991ff22ae23SHong Zhang 992ff22ae23SHong Zhang 993e4dd521cSBarry Smith /* ------------------------------------------------------------ */ 994ebe3b25bSBarry Smith /*MC 995f68a32c8SEmil Constantinescu TSRK - ODE and DAE solver using Runge-Kutta schemes 996ebe3b25bSBarry Smith 997f68a32c8SEmil Constantinescu The user should provide the right hand side of the equation 998f68a32c8SEmil Constantinescu using TSSetRHSFunction(). 999ebe3b25bSBarry Smith 1000f68a32c8SEmil Constantinescu Notes: 1001f68a32c8SEmil Constantinescu The default is TSRK3, it can be changed with TSRKSetType() or -ts_rk_type 1002ebe3b25bSBarry Smith 1003d41222bbSBarry Smith Level: beginner 1004d41222bbSBarry Smith 1005f68a32c8SEmil Constantinescu .seealso: TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3, 1006f68a32c8SEmil Constantinescu TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister() 1007ebe3b25bSBarry Smith 1008ebe3b25bSBarry Smith M*/ 1009e4dd521cSBarry Smith #undef __FUNCT__ 10105f70b478SJed Brown #define __FUNCT__ "TSCreate_RK" 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; 1028f68a32c8SEmil Constantinescu ts->ops->setfromoptions = TSSetFromOptions_RK; 1029ff22ae23SHong Zhang ts->ops->getstages = TSGetStages_RK; 103042f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_RK; 1031e4dd521cSBarry Smith 10320f7a1166SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_RK; 10330f7a1166SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_RK; 10340f7a1166SHong Zhang 10351566a47fSLisandro Dalcin ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr); 10361566a47fSLisandro Dalcin ts->data = (void*)rk; 1037e4dd521cSBarry Smith 1038f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr); 1039f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr); 1040*be5899b3SLisandro Dalcin 1041*be5899b3SLisandro Dalcin ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr); 1042e4dd521cSBarry Smith PetscFunctionReturn(0); 1043e4dd521cSBarry Smith } 1044