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 */ 10b45d2f2cSJed Brown #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 */ 41c235aa8dSHong Zhang Vec *VecDeltaLam; /* Increament of the adjoint sensitivity variable at stage*/ 42*b18ea86cSHong Zhang Vec *VecSensiTemp; /* Vector to be timed with Jacobian transpose*/ 43f68a32c8SEmil Constantinescu PetscScalar *work; /* Scalar work */ 44f68a32c8SEmil Constantinescu PetscReal stage_time; 45f68a32c8SEmil Constantinescu TSStepStatus status; 465f70b478SJed Brown } TS_RK; 47e4dd521cSBarry Smith 48f68a32c8SEmil Constantinescu /*MC 49f68a32c8SEmil Constantinescu TSRK1 - First order forward Euler scheme. 50262ff7bbSBarry Smith 51f68a32c8SEmil Constantinescu This method has one stage. 52f68a32c8SEmil Constantinescu 53f68a32c8SEmil Constantinescu Level: advanced 54f68a32c8SEmil Constantinescu 55f68a32c8SEmil Constantinescu .seealso: TSRK 56f68a32c8SEmil Constantinescu M*/ 57f68a32c8SEmil Constantinescu /*MC 582109b73fSDebojyoti Ghosh TSRK2A - Second order RK scheme. 59f68a32c8SEmil Constantinescu 60f68a32c8SEmil Constantinescu This method has two stages. 61f68a32c8SEmil Constantinescu 62f68a32c8SEmil Constantinescu Level: advanced 63f68a32c8SEmil Constantinescu 64f68a32c8SEmil Constantinescu .seealso: TSRK 65f68a32c8SEmil Constantinescu M*/ 66f68a32c8SEmil Constantinescu /*MC 67f68a32c8SEmil Constantinescu TSRK3 - Third order RK scheme. 68f68a32c8SEmil Constantinescu 69f68a32c8SEmil Constantinescu This method has three stages. 70f68a32c8SEmil Constantinescu 71f68a32c8SEmil Constantinescu Level: advanced 72f68a32c8SEmil Constantinescu 73f68a32c8SEmil Constantinescu .seealso: TSRK 74f68a32c8SEmil Constantinescu M*/ 75f68a32c8SEmil Constantinescu /*MC 762109b73fSDebojyoti Ghosh TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method. 772109b73fSDebojyoti Ghosh 782109b73fSDebojyoti Ghosh This method has four stages. 792109b73fSDebojyoti Ghosh 802109b73fSDebojyoti Ghosh Level: advanced 812109b73fSDebojyoti Ghosh 822109b73fSDebojyoti Ghosh .seealso: TSRK 832109b73fSDebojyoti Ghosh M*/ 842109b73fSDebojyoti Ghosh /*MC 85f68a32c8SEmil Constantinescu TSRK4 - Fourth order RK scheme. 86f68a32c8SEmil Constantinescu 872e698f8bSDebojyoti Ghosh This is the classical Runge-Kutta method with four stages. 88f68a32c8SEmil Constantinescu 89f68a32c8SEmil Constantinescu Level: advanced 90f68a32c8SEmil Constantinescu 91f68a32c8SEmil Constantinescu .seealso: TSRK 92f68a32c8SEmil Constantinescu M*/ 93f68a32c8SEmil Constantinescu /*MC 942e698f8bSDebojyoti Ghosh TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method. 95f68a32c8SEmil Constantinescu 96f68a32c8SEmil Constantinescu This method has six stages. 97f68a32c8SEmil Constantinescu 98f68a32c8SEmil Constantinescu Level: advanced 99f68a32c8SEmil Constantinescu 100f68a32c8SEmil Constantinescu .seealso: TSRK 101f68a32c8SEmil Constantinescu M*/ 1022109b73fSDebojyoti Ghosh /*MC 1032e698f8bSDebojyoti Ghosh TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method. 1042109b73fSDebojyoti Ghosh 1052109b73fSDebojyoti Ghosh This method has seven stages. 1062109b73fSDebojyoti Ghosh 1072109b73fSDebojyoti Ghosh Level: advanced 1082109b73fSDebojyoti Ghosh 1092109b73fSDebojyoti Ghosh .seealso: TSRK 1102109b73fSDebojyoti Ghosh M*/ 111262ff7bbSBarry Smith 112262ff7bbSBarry Smith #undef __FUNCT__ 113f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKRegisterAll" 114f68a32c8SEmil Constantinescu /*@C 115f68a32c8SEmil Constantinescu TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK 116262ff7bbSBarry Smith 117f68a32c8SEmil Constantinescu Not Collective, but should be called by all processes which will need the schemes to be registered 118262ff7bbSBarry Smith 119f68a32c8SEmil Constantinescu Level: advanced 120262ff7bbSBarry Smith 121f68a32c8SEmil Constantinescu .keywords: TS, TSRK, register, all 122262ff7bbSBarry Smith 123f68a32c8SEmil Constantinescu .seealso: TSRKRegisterDestroy() 124262ff7bbSBarry Smith @*/ 125f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void) 126262ff7bbSBarry Smith { 1274ac538c5SBarry Smith PetscErrorCode ierr; 128262ff7bbSBarry Smith 129262ff7bbSBarry Smith PetscFunctionBegin; 130f68a32c8SEmil Constantinescu if (TSRKRegisterAllCalled) PetscFunctionReturn(0); 131f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_TRUE; 132e4dd521cSBarry Smith 133e4dd521cSBarry Smith { 134f68a32c8SEmil Constantinescu const PetscReal 135f68a32c8SEmil Constantinescu A[1][1] = {{0.0}}, 136f68a32c8SEmil Constantinescu b[1] = {1.0}; 137f68a32c8SEmil Constantinescu ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,1,b);CHKERRQ(ierr); 138e4dd521cSBarry Smith } 139db046809SBarry Smith { 140f68a32c8SEmil Constantinescu const PetscReal 141f68a32c8SEmil Constantinescu A[2][2] = {{0.0,0.0}, 142f68a32c8SEmil Constantinescu {1.0,0.0}}, 143f68a32c8SEmil Constantinescu b[2] = {0.5,0.5}, 144f68a32c8SEmil Constantinescu bembed[2] = {1.0,0}; 145fdefd5e5SDebojyoti Ghosh ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,2,b);CHKERRQ(ierr); 146db046809SBarry Smith } 147f68a32c8SEmil Constantinescu { 148f68a32c8SEmil Constantinescu const PetscReal 149f68a32c8SEmil Constantinescu A[3][3] = {{0,0,0}, 150f68a32c8SEmil Constantinescu {2.0/3.0,0,0}, 151f68a32c8SEmil Constantinescu {-1.0/3.0,1.0,0}}, 152f68a32c8SEmil Constantinescu b[3] = {0.25,0.5,0.25}; 153fdefd5e5SDebojyoti Ghosh ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,3,b);CHKERRQ(ierr); 154fdefd5e5SDebojyoti Ghosh } 155fdefd5e5SDebojyoti Ghosh { 156fdefd5e5SDebojyoti Ghosh const PetscReal 157fdefd5e5SDebojyoti Ghosh A[4][4] = {{0,0,0,0}, 158fdefd5e5SDebojyoti Ghosh {1.0/2.0,0,0,0}, 159fdefd5e5SDebojyoti Ghosh {0,3.0/4.0,0,0}, 160fdefd5e5SDebojyoti Ghosh {2.0/9.0,1.0/3.0,4.0/9.0,0}}, 161fdefd5e5SDebojyoti Ghosh b[4] = {2.0/9.0,1.0/3.0,4.0/9.0,0}, 162fdefd5e5SDebojyoti Ghosh bembed[4] = {7.0/24.0,1.0/4.0,1.0/3.0,1.0/8.0}; 163fdefd5e5SDebojyoti Ghosh ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,3,b);CHKERRQ(ierr); 164db046809SBarry Smith } 165f68a32c8SEmil Constantinescu { 166f68a32c8SEmil Constantinescu const PetscReal 167f68a32c8SEmil Constantinescu A[4][4] = {{0,0,0,0}, 168f68a32c8SEmil Constantinescu {0.5,0,0,0}, 169f68a32c8SEmil Constantinescu {0,0.5,0,0}, 170f68a32c8SEmil Constantinescu {0,0,1.0,0}}, 171f68a32c8SEmil Constantinescu b[4] = {1.0/6.0,1.0/3.0,1.0/3.0,1.0/6.0}; 172fdefd5e5SDebojyoti Ghosh ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,4,b);CHKERRQ(ierr); 173db046809SBarry Smith } 174f68a32c8SEmil Constantinescu { 175f68a32c8SEmil Constantinescu const PetscReal 176f68a32c8SEmil Constantinescu A[6][6] = {{0,0,0,0,0,0}, 177f68a32c8SEmil Constantinescu {0.25,0,0,0,0,0}, 178f68a32c8SEmil Constantinescu {3.0/32.0,9.0/32.0,0,0,0,0}, 179f68a32c8SEmil Constantinescu {1932.0/2197.0,-7200.0/2197.0,7296.0/2197.0,0,0,0}, 180f68a32c8SEmil Constantinescu {439.0/216.0,-8.0,3680.0/513.0,-845.0/4104.0,0,0}, 181f68a32c8SEmil Constantinescu {-8.0/27.0,2.0,-3544.0/2565.0,1859.0/4104.0,-11.0/40.0,0}}, 182f68a32c8SEmil Constantinescu b[6] = {16.0/135.0,0,6656.0/12825.0,28561.0/56430.0,-9.0/50.0,2.0/55.0}, 183f68a32c8SEmil Constantinescu bembed[6] = {25.0/216.0,0,1408.0/2565.0,2197.0/4104.0,-1.0/5.0,0}; 184fdefd5e5SDebojyoti Ghosh ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,5,b);CHKERRQ(ierr); 185fdefd5e5SDebojyoti Ghosh } 186fdefd5e5SDebojyoti Ghosh { 187fdefd5e5SDebojyoti Ghosh const PetscReal 188fdefd5e5SDebojyoti Ghosh A[7][7] = {{0,0,0,0,0,0,0}, 189fdefd5e5SDebojyoti Ghosh {1.0/5.0,0,0,0,0,0,0}, 190fdefd5e5SDebojyoti Ghosh {3.0/40.0,9.0/40.0,0,0,0,0,0}, 191fdefd5e5SDebojyoti Ghosh {44.0/45.0,-56.0/15.0,32.0/9.0,0,0,0,0}, 192fdefd5e5SDebojyoti Ghosh {19372.0/6561.0,-25360.0/2187.0,64448.0/6561.0,-212.0/729.0,0,0,0}, 193fdefd5e5SDebojyoti Ghosh {9017.0/3168.0,-355.0/33.0,46732.0/5247.0,49.0/176.0,-5103.0/18656.0,0,0}, 194fdefd5e5SDebojyoti Ghosh {35.0/384.0,0,500.0/1113.0,125.0/192.0,-2187.0/6784.0,11.0/84.0,0}}, 195fdefd5e5SDebojyoti 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}, 196fdefd5e5SDebojyoti 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}; 197fdefd5e5SDebojyoti Ghosh ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,5,b);CHKERRQ(ierr); 198f68a32c8SEmil Constantinescu } 199db046809SBarry Smith PetscFunctionReturn(0); 200db046809SBarry Smith } 201db046809SBarry Smith 202e4dd521cSBarry Smith #undef __FUNCT__ 203f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKRegisterDestroy" 204f68a32c8SEmil Constantinescu /*@C 205f68a32c8SEmil Constantinescu TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister(). 206f68a32c8SEmil Constantinescu 207f68a32c8SEmil Constantinescu Not Collective 208f68a32c8SEmil Constantinescu 209f68a32c8SEmil Constantinescu Level: advanced 210f68a32c8SEmil Constantinescu 211f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy 212f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll() 213f68a32c8SEmil Constantinescu @*/ 214f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void) 215e4dd521cSBarry Smith { 216f68a32c8SEmil Constantinescu PetscErrorCode ierr; 217f68a32c8SEmil Constantinescu RKTableauLink link; 218f68a32c8SEmil Constantinescu 219f68a32c8SEmil Constantinescu PetscFunctionBegin; 220f68a32c8SEmil Constantinescu while ((link = RKTableauList)) { 221f68a32c8SEmil Constantinescu RKTableau t = &link->tab; 222f68a32c8SEmil Constantinescu RKTableauList = link->next; 223f68a32c8SEmil Constantinescu ierr = PetscFree3(t->A,t->b,t->c); CHKERRQ(ierr); 224f68a32c8SEmil Constantinescu ierr = PetscFree (t->bembed); CHKERRQ(ierr); 225f68a32c8SEmil Constantinescu ierr = PetscFree (t->binterp); CHKERRQ(ierr); 226f68a32c8SEmil Constantinescu ierr = PetscFree (t->name); CHKERRQ(ierr); 227f68a32c8SEmil Constantinescu ierr = PetscFree (link); CHKERRQ(ierr); 228f68a32c8SEmil Constantinescu } 229f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_FALSE; 230f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 231f68a32c8SEmil Constantinescu } 232f68a32c8SEmil Constantinescu 233f68a32c8SEmil Constantinescu #undef __FUNCT__ 234f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKInitializePackage" 235f68a32c8SEmil Constantinescu /*@C 236f68a32c8SEmil Constantinescu TSRKInitializePackage - This function initializes everything in the TSRK package. It is called 237f68a32c8SEmil Constantinescu from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RK() 238f68a32c8SEmil Constantinescu when using static libraries. 239f68a32c8SEmil Constantinescu 240f68a32c8SEmil Constantinescu Level: developer 241f68a32c8SEmil Constantinescu 242f68a32c8SEmil Constantinescu .keywords: TS, TSRK, initialize, package 243f68a32c8SEmil Constantinescu .seealso: PetscInitialize() 244f68a32c8SEmil Constantinescu @*/ 245f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void) 246f68a32c8SEmil Constantinescu { 247186e87acSLisandro Dalcin PetscErrorCode ierr; 248e4dd521cSBarry Smith 249e4dd521cSBarry Smith PetscFunctionBegin; 250f68a32c8SEmil Constantinescu if (TSRKPackageInitialized) PetscFunctionReturn(0); 251f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_TRUE; 252f68a32c8SEmil Constantinescu ierr = TSRKRegisterAll();CHKERRQ(ierr); 253f68a32c8SEmil Constantinescu ierr = PetscObjectComposedDataRegister(&explicit_stage_time_id);CHKERRQ(ierr); 254f68a32c8SEmil Constantinescu ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr); 255f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 256f68a32c8SEmil Constantinescu } 257186e87acSLisandro Dalcin 258f68a32c8SEmil Constantinescu #undef __FUNCT__ 259f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKFinalizePackage" 260f68a32c8SEmil Constantinescu /*@C 261f68a32c8SEmil Constantinescu TSRKFinalizePackage - This function destroys everything in the TSRK package. It is 262f68a32c8SEmil Constantinescu called from PetscFinalize(). 263186e87acSLisandro Dalcin 264f68a32c8SEmil Constantinescu Level: developer 265f68a32c8SEmil Constantinescu 266f68a32c8SEmil Constantinescu .keywords: Petsc, destroy, package 267f68a32c8SEmil Constantinescu .seealso: PetscFinalize() 268f68a32c8SEmil Constantinescu @*/ 269f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void) 270f68a32c8SEmil Constantinescu { 271f68a32c8SEmil Constantinescu PetscErrorCode ierr; 272f68a32c8SEmil Constantinescu 273f68a32c8SEmil Constantinescu PetscFunctionBegin; 274f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_FALSE; 275f68a32c8SEmil Constantinescu ierr = TSRKRegisterDestroy();CHKERRQ(ierr); 276f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 277f68a32c8SEmil Constantinescu } 278f68a32c8SEmil Constantinescu 279f68a32c8SEmil Constantinescu #undef __FUNCT__ 280f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKRegister" 281f68a32c8SEmil Constantinescu /*@C 282f68a32c8SEmil Constantinescu TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 283f68a32c8SEmil Constantinescu 284f68a32c8SEmil Constantinescu Not Collective, but the same schemes should be registered on all processes on which they will be used 285f68a32c8SEmil Constantinescu 286f68a32c8SEmil Constantinescu Input Parameters: 287f68a32c8SEmil Constantinescu + name - identifier for method 288f68a32c8SEmil Constantinescu . order - approximation order of method 289f68a32c8SEmil Constantinescu . s - number of stages, this is the dimension of the matrices below 290f68a32c8SEmil Constantinescu . A - stage coefficients (dimension s*s, row-major) 291f68a32c8SEmil Constantinescu . b - step completion table (dimension s; NULL to use last row of A) 292f68a32c8SEmil Constantinescu . c - abscissa (dimension s; NULL to use row sums of A) 293f68a32c8SEmil Constantinescu . bembed - completion table for embedded method (dimension s; NULL if not available) 294f68a32c8SEmil Constantinescu . pinterp - Order of the interpolation scheme, equal to the number of columns of binterp 295f68a32c8SEmil Constantinescu - binterp - Coefficients of the interpolation formula (dimension s*pinterp; NULL to reuse binterpt) 296f68a32c8SEmil Constantinescu 297f68a32c8SEmil Constantinescu Notes: 298f68a32c8SEmil Constantinescu Several RK methods are provided, this function is only needed to create new methods. 299f68a32c8SEmil Constantinescu 300f68a32c8SEmil Constantinescu Level: advanced 301f68a32c8SEmil Constantinescu 302f68a32c8SEmil Constantinescu .keywords: TS, register 303f68a32c8SEmil Constantinescu 304f68a32c8SEmil Constantinescu .seealso: TSRK 305f68a32c8SEmil Constantinescu @*/ 306f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s, 307f68a32c8SEmil Constantinescu const PetscReal A[],const PetscReal b[],const PetscReal c[], 308f68a32c8SEmil Constantinescu const PetscReal bembed[], 309f68a32c8SEmil Constantinescu PetscInt pinterp,const PetscReal binterp[]) 310f68a32c8SEmil Constantinescu { 311f68a32c8SEmil Constantinescu PetscErrorCode ierr; 312f68a32c8SEmil Constantinescu RKTableauLink link; 313f68a32c8SEmil Constantinescu RKTableau t; 314f68a32c8SEmil Constantinescu PetscInt i,j; 315f68a32c8SEmil Constantinescu 316f68a32c8SEmil Constantinescu PetscFunctionBegin; 317f68a32c8SEmil Constantinescu ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr); 318f68a32c8SEmil Constantinescu ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr); 319f68a32c8SEmil Constantinescu t = &link->tab; 320f68a32c8SEmil Constantinescu ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 321f68a32c8SEmil Constantinescu t->order = order; 322f68a32c8SEmil Constantinescu t->s = s; 323dcca6d9dSJed Brown ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr); 324f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 325f68a32c8SEmil Constantinescu if (b) { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); } 326f68a32c8SEmil Constantinescu else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i]; 327f68a32c8SEmil Constantinescu if (c) { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); } 328f68a32c8SEmil 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]; 329d760c35bSDebojyoti Ghosh t->FSAL = PETSC_TRUE; 330d760c35bSDebojyoti Ghosh for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE; 331f68a32c8SEmil Constantinescu if (bembed) { 332785e854fSJed Brown ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr); 333f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr); 334f68a32c8SEmil Constantinescu } 335f68a32c8SEmil Constantinescu 336f68a32c8SEmil Constantinescu t->pinterp = pinterp; 337785e854fSJed Brown ierr = PetscMalloc1(s*pinterp,&t->binterp);CHKERRQ(ierr); 338f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->binterp,binterp,s*pinterp*sizeof(binterp[0]));CHKERRQ(ierr); 339f68a32c8SEmil Constantinescu link->next = RKTableauList; 340f68a32c8SEmil Constantinescu RKTableauList = link; 341f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 342f68a32c8SEmil Constantinescu } 343f68a32c8SEmil Constantinescu 344f68a32c8SEmil Constantinescu #undef __FUNCT__ 345f68a32c8SEmil Constantinescu #define __FUNCT__ "TSEvaluateStep_RK" 346e4dd521cSBarry Smith /* 347f68a32c8SEmil Constantinescu The step completion formula is 348e4dd521cSBarry Smith 349f68a32c8SEmil Constantinescu x1 = x0 + h b^T YdotRHS 350f68a32c8SEmil Constantinescu 351f68a32c8SEmil Constantinescu This function can be called before or after ts->vec_sol has been updated. 352f68a32c8SEmil Constantinescu Suppose we have a completion formula (b) and an embedded formula (be) of different order. 353f68a32c8SEmil Constantinescu We can write 354f68a32c8SEmil Constantinescu 355f68a32c8SEmil Constantinescu x1e = x0 + h be^T YdotRHS 356f68a32c8SEmil Constantinescu = x1 - h b^T YdotRHS + h be^T YdotRHS 357f68a32c8SEmil Constantinescu = x1 + h (be - b)^T YdotRHS 358f68a32c8SEmil Constantinescu 359f68a32c8SEmil Constantinescu so we can evaluate the method with different order even after the step has been optimistically completed. 360e4dd521cSBarry Smith */ 361f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done) 362f68a32c8SEmil Constantinescu { 363f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 364f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 365f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 366f68a32c8SEmil Constantinescu PetscReal h; 367f68a32c8SEmil Constantinescu PetscInt s = tab->s,j; 368f68a32c8SEmil Constantinescu PetscErrorCode ierr; 369f68a32c8SEmil Constantinescu 370f68a32c8SEmil Constantinescu PetscFunctionBegin; 371f68a32c8SEmil Constantinescu switch (rk->status) { 372f68a32c8SEmil Constantinescu case TS_STEP_INCOMPLETE: 373f68a32c8SEmil Constantinescu case TS_STEP_PENDING: 374f68a32c8SEmil Constantinescu h = ts->time_step; break; 375f68a32c8SEmil Constantinescu case TS_STEP_COMPLETE: 376f68a32c8SEmil Constantinescu h = ts->time_step_prev; break; 377f68a32c8SEmil Constantinescu default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 378f68a32c8SEmil Constantinescu } 379f68a32c8SEmil Constantinescu if (order == tab->order) { 380f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { 381f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 382f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->b[j]; 383f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 384f68a32c8SEmil Constantinescu } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 385f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 386f68a32c8SEmil Constantinescu } else if (order == tab->order-1) { 387f68a32c8SEmil Constantinescu if (!tab->bembed) goto unavailable; 388f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (be) */ 389f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 390f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 391f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 392f68a32c8SEmil Constantinescu } else { /* Rollback and re-complete using (be-b) */ 393f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 394f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 395f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 396f68a32c8SEmil Constantinescu } 397f68a32c8SEmil Constantinescu if (done) *done = PETSC_TRUE; 398f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 399f68a32c8SEmil Constantinescu } 400f68a32c8SEmil Constantinescu unavailable: 401f68a32c8SEmil Constantinescu if (done) *done = PETSC_FALSE; 402f68a32c8SEmil Constantinescu else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"RK '%s' of order %D cannot evaluate step at order %D",tab->name,tab->order,order); 403f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 404f68a32c8SEmil Constantinescu } 405f68a32c8SEmil Constantinescu 406f68a32c8SEmil Constantinescu #undef __FUNCT__ 407f68a32c8SEmil Constantinescu #define __FUNCT__ "TSStep_RK" 408f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts) 409f68a32c8SEmil Constantinescu { 410f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 411f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 412f68a32c8SEmil Constantinescu const PetscInt s = tab->s; 413f68a32c8SEmil Constantinescu const PetscReal *A = tab->A,*b = tab->b,*c = tab->c; 414f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 415f68a32c8SEmil Constantinescu Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 416f68a32c8SEmil Constantinescu TSAdapt adapt; 417f68a32c8SEmil Constantinescu PetscInt i,j,reject,next_scheme; 418f68a32c8SEmil Constantinescu PetscReal next_time_step; 419f68a32c8SEmil Constantinescu PetscReal t; 420f68a32c8SEmil Constantinescu PetscBool accept; 421f68a32c8SEmil Constantinescu PetscErrorCode ierr; 422f68a32c8SEmil Constantinescu 423f68a32c8SEmil Constantinescu PetscFunctionBegin; 424f68a32c8SEmil Constantinescu 425f68a32c8SEmil Constantinescu next_time_step = ts->time_step; 426f68a32c8SEmil Constantinescu t = ts->ptime; 427f68a32c8SEmil Constantinescu accept = PETSC_TRUE; 428f68a32c8SEmil Constantinescu rk->status = TS_STEP_INCOMPLETE; 429f68a32c8SEmil Constantinescu 430f68a32c8SEmil Constantinescu 431f68a32c8SEmil Constantinescu for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) { 432f68a32c8SEmil Constantinescu PetscReal h = ts->time_step; 433f68a32c8SEmil Constantinescu ierr = TSPreStep(ts);CHKERRQ(ierr); 434f68a32c8SEmil Constantinescu for (i=0; i<s; i++) { 4359be3e283SDebojyoti Ghosh rk->stage_time = t + h*c[i]; 4369be3e283SDebojyoti Ghosh ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr); 437f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 438f68a32c8SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 439f68a32c8SEmil Constantinescu ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 4409be3e283SDebojyoti Ghosh ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 441f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 442f68a32c8SEmil Constantinescu ierr = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr); 443f68a32c8SEmil Constantinescu if (!accept) goto reject_step; 444f68a32c8SEmil Constantinescu ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 445f68a32c8SEmil Constantinescu } 446f68a32c8SEmil Constantinescu ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 447f68a32c8SEmil Constantinescu rk->status = TS_STEP_PENDING; 448f68a32c8SEmil Constantinescu 449f68a32c8SEmil Constantinescu /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */ 450f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 451f68a32c8SEmil Constantinescu ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 452f68a32c8SEmil Constantinescu ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr); 453f68a32c8SEmil Constantinescu ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr); 454f68a32c8SEmil Constantinescu if (accept) { 455f68a32c8SEmil Constantinescu /* ignore next_scheme for now */ 456f68a32c8SEmil Constantinescu ts->ptime += ts->time_step; 457f68a32c8SEmil Constantinescu ts->time_step = next_time_step; 458e3caeda6SBarry Smith ts->steps++; 459f68a32c8SEmil Constantinescu rk->status = TS_STEP_COMPLETE; 460f68a32c8SEmil Constantinescu ierr = PetscObjectComposedDataSetReal((PetscObject)ts->vec_sol,explicit_stage_time_id,ts->ptime);CHKERRQ(ierr); 461f68a32c8SEmil Constantinescu break; 462f68a32c8SEmil Constantinescu } else { /* Roll back the current step */ 463f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = -h*b[j]; 464f68a32c8SEmil Constantinescu ierr = VecMAXPY(ts->vec_sol,s,w,rk->YdotRHS);CHKERRQ(ierr); 465f68a32c8SEmil Constantinescu ts->time_step = next_time_step; 466f68a32c8SEmil Constantinescu rk->status = TS_STEP_INCOMPLETE; 467f68a32c8SEmil Constantinescu } 468f68a32c8SEmil Constantinescu reject_step: continue; 469f68a32c8SEmil Constantinescu } 470f68a32c8SEmil Constantinescu if (rk->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED; 471f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 472e4dd521cSBarry Smith } 473e4dd521cSBarry Smith 474f68a32c8SEmil Constantinescu #undef __FUNCT__ 475d2daff3dSHong Zhang #define __FUNCT__ "TSStepAdj_RK" 476c235aa8dSHong Zhang static PetscErrorCode TSStepAdj_RK(TS ts) 477d2daff3dSHong Zhang { 478c235aa8dSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 479c235aa8dSHong Zhang RKTableau tab = rk->tableau; 480c235aa8dSHong Zhang const PetscInt s = tab->s; 481c235aa8dSHong Zhang const PetscReal *A = tab->A,*b = tab->b,*c = tab->c; 482c235aa8dSHong Zhang PetscScalar *w = rk->work; 483*b18ea86cSHong Zhang Vec *Y = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecSensiTemp = rk->VecSensiTemp; 484*b18ea86cSHong Zhang PetscInt i,j,nadj; 485c235aa8dSHong Zhang PetscReal t; 486d2daff3dSHong Zhang PetscErrorCode ierr; 487c235aa8dSHong Zhang 488d2daff3dSHong Zhang PetscFunctionBegin; 489d2daff3dSHong Zhang 490*b18ea86cSHong Zhang /*ierr = TSStep_RK(ts);CHKERRQ(ierr); // reuse TSStep */ 491d2daff3dSHong Zhang 492c235aa8dSHong Zhang t = ts->ptime; 493c235aa8dSHong Zhang rk->status = TS_STEP_INCOMPLETE; 494*b18ea86cSHong Zhang PetscReal h = ts->time_step; 495c235aa8dSHong Zhang ierr = TSPreStep(ts);CHKERRQ(ierr); 496c235aa8dSHong Zhang for (i=s-1; i>=0; i--) { 497c235aa8dSHong Zhang Mat J,Jp; 498c235aa8dSHong Zhang rk->stage_time = t + h*(1.0-c[i]); 499*b18ea86cSHong Zhang for (nadj=0; nadj<ts->numberadjs; nadj++) { 500*b18ea86cSHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr); 501*b18ea86cSHong Zhang ierr = VecScale(VecSensiTemp[nadj],-h*b[i]); 502c235aa8dSHong Zhang for (j=i+1; j<s; j++) { 503*b18ea86cSHong Zhang ierr = VecAXPY(VecSensiTemp[nadj],-h*A[j*s+i],VecDeltaLam[nadj*s+j]); 504*b18ea86cSHong Zhang } 505c235aa8dSHong Zhang } 506c235aa8dSHong Zhang ierr = TSGetRHSJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 507c235aa8dSHong Zhang ierr = TSComputeRHSJacobian(ts,rk->stage_time,Y[i],J,Jp);CHKERRQ(ierr); 508*b18ea86cSHong Zhang for (nadj=0; nadj<ts->numberadjs; nadj++) { 509*b18ea86cSHong Zhang ierr = MatMultTranspose(J,VecSensiTemp[nadj],VecDeltaLam[nadj*s+i]);CHKERRQ(ierr); 510*b18ea86cSHong Zhang } 511c235aa8dSHong Zhang } 512c235aa8dSHong Zhang 513c235aa8dSHong Zhang for (j=0; j<s; j++) w[j] = 1.0; 514*b18ea86cSHong Zhang for (nadj=0; nadj<ts->numberadjs; nadj++) { 515*b18ea86cSHong Zhang ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr); 516*b18ea86cSHong Zhang } 517c235aa8dSHong Zhang ts->ptime += ts->time_step; 518c235aa8dSHong Zhang ts->steps++; 519c235aa8dSHong Zhang rk->status = TS_STEP_COMPLETE; 520*b18ea86cSHong Zhang ierr = PetscObjectComposedDataSetReal((PetscObject)ts->vecs_sensi,explicit_stage_time_id,ts->ptime);CHKERRQ(ierr); 521d2daff3dSHong Zhang PetscFunctionReturn(0); 522d2daff3dSHong Zhang } 523d2daff3dSHong Zhang 524d2daff3dSHong Zhang #undef __FUNCT__ 525f68a32c8SEmil Constantinescu #define __FUNCT__ "TSInterpolate_RK" 526f68a32c8SEmil Constantinescu static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X) 527f68a32c8SEmil Constantinescu { 528f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 529f68a32c8SEmil Constantinescu PetscInt s = rk->tableau->s,pinterp = rk->tableau->pinterp,i,j; 530f68a32c8SEmil Constantinescu PetscReal h; 531f68a32c8SEmil Constantinescu PetscReal tt,t; 532f68a32c8SEmil Constantinescu PetscScalar *b; 533f68a32c8SEmil Constantinescu const PetscReal *B = rk->tableau->binterp; 534f68a32c8SEmil Constantinescu PetscErrorCode ierr; 535e4dd521cSBarry Smith 536f68a32c8SEmil Constantinescu PetscFunctionBegin; 537f68a32c8SEmil Constantinescu if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 538f68a32c8SEmil Constantinescu switch (rk->status) { 539f68a32c8SEmil Constantinescu case TS_STEP_INCOMPLETE: 540f68a32c8SEmil Constantinescu case TS_STEP_PENDING: 541f68a32c8SEmil Constantinescu h = ts->time_step; 542f68a32c8SEmil Constantinescu t = (itime - ts->ptime)/h; 543f68a32c8SEmil Constantinescu break; 544f68a32c8SEmil Constantinescu case TS_STEP_COMPLETE: 545f68a32c8SEmil Constantinescu h = ts->time_step_prev; 546f68a32c8SEmil Constantinescu t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 547f68a32c8SEmil Constantinescu break; 548f68a32c8SEmil Constantinescu default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 549e4dd521cSBarry Smith } 550785e854fSJed Brown ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 551f68a32c8SEmil Constantinescu for (i=0; i<s; i++) b[i] = 0; 552f68a32c8SEmil Constantinescu for (j=0,tt=t; j<pinterp; j++,tt*=t) { 553f68a32c8SEmil Constantinescu for (i=0; i<s; i++) { 554f68a32c8SEmil Constantinescu b[i] += h * B[i*pinterp+j] * tt; 555e4dd521cSBarry Smith } 556f68a32c8SEmil Constantinescu } 557f68a32c8SEmil Constantinescu ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr); 558f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr); 559f68a32c8SEmil Constantinescu ierr = PetscFree(b);CHKERRQ(ierr); 560e4dd521cSBarry Smith PetscFunctionReturn(0); 561e4dd521cSBarry Smith } 562e4dd521cSBarry Smith 563e4dd521cSBarry Smith /*------------------------------------------------------------*/ 564e4dd521cSBarry Smith #undef __FUNCT__ 565277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_RK" 566277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts) 567e4dd521cSBarry Smith { 5685f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 569f68a32c8SEmil Constantinescu PetscInt s; 5706849ba73SBarry Smith PetscErrorCode ierr; 571e4dd521cSBarry Smith 572e4dd521cSBarry Smith PetscFunctionBegin; 573f68a32c8SEmil Constantinescu if (!rk->tableau) PetscFunctionReturn(0); 574f68a32c8SEmil Constantinescu s = rk->tableau->s; 575f68a32c8SEmil Constantinescu ierr = VecDestroyVecs(s,&rk->Y);CHKERRQ(ierr); 576f68a32c8SEmil Constantinescu ierr = VecDestroyVecs(s,&rk->YdotRHS);CHKERRQ(ierr); 577c235aa8dSHong Zhang if(ts->reverse_mode) { 578*b18ea86cSHong Zhang ierr = VecDestroyVecs(s*ts->numberadjs,&rk->VecDeltaLam);CHKERRQ(ierr); 579*b18ea86cSHong Zhang ierr = VecDestroyVecs(ts->numberadjs,&rk->VecSensiTemp);CHKERRQ(ierr); 580c235aa8dSHong Zhang } 581f68a32c8SEmil Constantinescu ierr = PetscFree(rk->work);CHKERRQ(ierr); 582277b19d0SLisandro Dalcin PetscFunctionReturn(0); 583e4dd521cSBarry Smith } 584277b19d0SLisandro Dalcin 585277b19d0SLisandro Dalcin #undef __FUNCT__ 586277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_RK" 587277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_RK(TS ts) 588277b19d0SLisandro Dalcin { 589277b19d0SLisandro Dalcin PetscErrorCode ierr; 590277b19d0SLisandro Dalcin 591277b19d0SLisandro Dalcin PetscFunctionBegin; 592277b19d0SLisandro Dalcin ierr = TSReset_RK(ts);CHKERRQ(ierr); 593277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 594f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr); 595f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr); 596f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 597f68a32c8SEmil Constantinescu } 598f68a32c8SEmil Constantinescu 599f68a32c8SEmil Constantinescu 600f68a32c8SEmil Constantinescu #undef __FUNCT__ 601f68a32c8SEmil Constantinescu #define __FUNCT__ "DMCoarsenHook_TSRK" 602f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx) 603f68a32c8SEmil Constantinescu { 604f68a32c8SEmil Constantinescu PetscFunctionBegin; 605f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 606f68a32c8SEmil Constantinescu } 607f68a32c8SEmil Constantinescu 608f68a32c8SEmil Constantinescu #undef __FUNCT__ 609f68a32c8SEmil Constantinescu #define __FUNCT__ "DMRestrictHook_TSRK" 610f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 611f68a32c8SEmil Constantinescu { 612f68a32c8SEmil Constantinescu PetscFunctionBegin; 613f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 614f68a32c8SEmil Constantinescu } 615f68a32c8SEmil Constantinescu 616f68a32c8SEmil Constantinescu 617f68a32c8SEmil Constantinescu #undef __FUNCT__ 618f68a32c8SEmil Constantinescu #define __FUNCT__ "DMSubDomainHook_TSRK" 619f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx) 620f68a32c8SEmil Constantinescu { 621f68a32c8SEmil Constantinescu PetscFunctionBegin; 622f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 623f68a32c8SEmil Constantinescu } 624f68a32c8SEmil Constantinescu 625f68a32c8SEmil Constantinescu #undef __FUNCT__ 626f68a32c8SEmil Constantinescu #define __FUNCT__ "DMSubDomainRestrictHook_TSRK" 627f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 628f68a32c8SEmil Constantinescu { 629f68a32c8SEmil Constantinescu 630f68a32c8SEmil Constantinescu PetscFunctionBegin; 631f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 632f68a32c8SEmil Constantinescu } 633c235aa8dSHong Zhang /* 634f68a32c8SEmil Constantinescu #undef __FUNCT__ 635d2daff3dSHong Zhang #define __FUNCT__ "RKSetAdjCoe" 636d2daff3dSHong Zhang static PetscErrorCode RKSetAdjCoe(RKTableau tab) 637d2daff3dSHong Zhang { 638d2daff3dSHong Zhang PetscReal *A,*b; 639d2daff3dSHong Zhang PetscInt s,i,j; 640d2daff3dSHong Zhang PetscErrorCode ierr; 641d2daff3dSHong Zhang 642d2daff3dSHong Zhang PetscFunctionBegin; 643d2daff3dSHong Zhang s = tab->s; 644d2daff3dSHong Zhang ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr); 645d2daff3dSHong Zhang 646d2daff3dSHong Zhang for (i=0; i<s; i++) 647d2daff3dSHong Zhang for (j=0; j<s; j++) { 648d2daff3dSHong 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]; 649d2daff3dSHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr); 650d2daff3dSHong Zhang } 651d2daff3dSHong Zhang 652d2daff3dSHong Zhang for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i]; 653d2daff3dSHong Zhang 654d2daff3dSHong Zhang ierr = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 655d2daff3dSHong Zhang ierr = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr); 656d2daff3dSHong Zhang ierr = PetscFree2(A,b);CHKERRQ(ierr); 657d2daff3dSHong Zhang PetscFunctionReturn(0); 658d2daff3dSHong Zhang } 659c235aa8dSHong Zhang */ 660d2daff3dSHong Zhang #undef __FUNCT__ 661f68a32c8SEmil Constantinescu #define __FUNCT__ "TSSetUp_RK" 662f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts) 663f68a32c8SEmil Constantinescu { 664f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 665f68a32c8SEmil Constantinescu RKTableau tab; 666f68a32c8SEmil Constantinescu PetscInt s; 667f68a32c8SEmil Constantinescu PetscErrorCode ierr; 668f68a32c8SEmil Constantinescu DM dm; 669f68a32c8SEmil Constantinescu 670f68a32c8SEmil Constantinescu PetscFunctionBegin; 671f68a32c8SEmil Constantinescu if (!rk->tableau) { 672f68a32c8SEmil Constantinescu ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr); 673f68a32c8SEmil Constantinescu } 674f68a32c8SEmil Constantinescu tab = rk->tableau; 675f68a32c8SEmil Constantinescu s = tab->s; 676c235aa8dSHong Zhang // old 677c235aa8dSHong Zhang //if (ts->reverse_mode) { 678c235aa8dSHong Zhang // ierr = RKSetAdjCoe(tab);CHKERRQ(ierr); 679c235aa8dSHong Zhang //} 680f68a32c8SEmil Constantinescu ierr = VecDuplicateVecs(ts->vec_sol,s,&rk->Y);CHKERRQ(ierr); 681f68a32c8SEmil Constantinescu ierr = VecDuplicateVecs(ts->vec_sol,s,&rk->YdotRHS);CHKERRQ(ierr); 682c235aa8dSHong Zhang if (ts->reverse_mode) { 683*b18ea86cSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numberadjs,&rk->VecDeltaLam);CHKERRQ(ierr); 684*b18ea86cSHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numberadjs,&rk->VecSensiTemp);CHKERRQ(ierr); 685c235aa8dSHong Zhang } 686785e854fSJed Brown ierr = PetscMalloc1(s,&rk->work);CHKERRQ(ierr); 687f68a32c8SEmil Constantinescu ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 688f68a32c8SEmil Constantinescu if (dm) { 689f68a32c8SEmil Constantinescu ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 690f68a32c8SEmil Constantinescu ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 691f68a32c8SEmil Constantinescu } 692e4dd521cSBarry Smith PetscFunctionReturn(0); 693e4dd521cSBarry Smith } 694d2daff3dSHong Zhang 695e4dd521cSBarry Smith /*------------------------------------------------------------*/ 696e4dd521cSBarry Smith 697e4dd521cSBarry Smith #undef __FUNCT__ 6985f70b478SJed Brown #define __FUNCT__ "TSSetFromOptions_RK" 6998c34d3f5SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptions *PetscOptionsObject,TS ts) 700e4dd521cSBarry Smith { 701dfbe8321SBarry Smith PetscErrorCode ierr; 702f68a32c8SEmil Constantinescu char rktype[256]; 703262ff7bbSBarry Smith 704e4dd521cSBarry Smith PetscFunctionBegin; 705e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr); 706f68a32c8SEmil Constantinescu { 707f68a32c8SEmil Constantinescu RKTableauLink link; 708f68a32c8SEmil Constantinescu PetscInt count,choice; 709f68a32c8SEmil Constantinescu PetscBool flg; 710f68a32c8SEmil Constantinescu const char **namelist; 711f68a32c8SEmil Constantinescu ierr = PetscStrncpy(rktype,TSRKDefault,sizeof(rktype));CHKERRQ(ierr); 712f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) ; 713785e854fSJed Brown ierr = PetscMalloc1(count,&namelist);CHKERRQ(ierr); 714f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 715f68a32c8SEmil Constantinescu ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rktype,&choice,&flg);CHKERRQ(ierr); 716f68a32c8SEmil Constantinescu ierr = TSRKSetType(ts,flg ? namelist[choice] : rktype);CHKERRQ(ierr); 717f68a32c8SEmil Constantinescu ierr = PetscFree(namelist);CHKERRQ(ierr); 718f68a32c8SEmil Constantinescu } 719262ff7bbSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 720e4dd521cSBarry Smith PetscFunctionReturn(0); 721e4dd521cSBarry Smith } 722e4dd521cSBarry Smith 723e4dd521cSBarry Smith #undef __FUNCT__ 724f68a32c8SEmil Constantinescu #define __FUNCT__ "PetscFormatRealArray" 725f68a32c8SEmil Constantinescu static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[]) 726f68a32c8SEmil Constantinescu { 727f68a32c8SEmil Constantinescu PetscErrorCode ierr; 728f68a32c8SEmil Constantinescu PetscInt i; 729f68a32c8SEmil Constantinescu size_t left,count; 730f68a32c8SEmil Constantinescu char *p; 731f68a32c8SEmil Constantinescu 732f68a32c8SEmil Constantinescu PetscFunctionBegin; 733f68a32c8SEmil Constantinescu for (i=0,p=buf,left=len; i<n; i++) { 734f68a32c8SEmil Constantinescu ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr); 735f68a32c8SEmil Constantinescu if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer"); 736f68a32c8SEmil Constantinescu left -= count; 737f68a32c8SEmil Constantinescu p += count; 738f68a32c8SEmil Constantinescu *p++ = ' '; 739f68a32c8SEmil Constantinescu } 740f68a32c8SEmil Constantinescu p[i ? 0 : -1] = 0; 741f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 742f68a32c8SEmil Constantinescu } 743f68a32c8SEmil Constantinescu 744f68a32c8SEmil Constantinescu #undef __FUNCT__ 7455f70b478SJed Brown #define __FUNCT__ "TSView_RK" 7465f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer) 747e4dd521cSBarry Smith { 7485f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 749f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 7508370ee3bSLisandro Dalcin PetscBool iascii; 751dfbe8321SBarry Smith PetscErrorCode ierr; 752f68a32c8SEmil Constantinescu TSAdapt adapt; 7532cdf8120Sasbjorn 754e4dd521cSBarry Smith PetscFunctionBegin; 755251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 7568370ee3bSLisandro Dalcin if (iascii) { 757f68a32c8SEmil Constantinescu TSRKType rktype; 758f68a32c8SEmil Constantinescu char buf[512]; 759f68a32c8SEmil Constantinescu ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr); 760f68a32c8SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," RK %s\n",rktype);CHKERRQ(ierr); 761f68a32c8SEmil Constantinescu ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 762f68a32c8SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," Abscissa c = %s\n",buf);CHKERRQ(ierr); 763d760c35bSDebojyoti Ghosh ierr = PetscViewerASCIIPrintf(viewer,"FSAL: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr); 7648370ee3bSLisandro Dalcin } 765f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 766f68a32c8SEmil Constantinescu ierr = TSAdaptView(adapt,viewer);CHKERRQ(ierr); 767f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 768f68a32c8SEmil Constantinescu } 769f68a32c8SEmil Constantinescu 770f68a32c8SEmil Constantinescu #undef __FUNCT__ 771f68a32c8SEmil Constantinescu #define __FUNCT__ "TSLoad_RK" 772f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer) 773f68a32c8SEmil Constantinescu { 774f68a32c8SEmil Constantinescu PetscErrorCode ierr; 775f68a32c8SEmil Constantinescu TSAdapt tsadapt; 776f68a32c8SEmil Constantinescu 777f68a32c8SEmil Constantinescu PetscFunctionBegin; 778f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&tsadapt);CHKERRQ(ierr); 779f68a32c8SEmil Constantinescu ierr = TSAdaptLoad(tsadapt,viewer);CHKERRQ(ierr); 780f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 781f68a32c8SEmil Constantinescu } 782f68a32c8SEmil Constantinescu 783f68a32c8SEmil Constantinescu #undef __FUNCT__ 784f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKSetType" 785f68a32c8SEmil Constantinescu /*@C 786f68a32c8SEmil Constantinescu TSRKSetType - Set the type of RK scheme 787f68a32c8SEmil Constantinescu 788f68a32c8SEmil Constantinescu Logically collective 789f68a32c8SEmil Constantinescu 790f68a32c8SEmil Constantinescu Input Parameter: 791f68a32c8SEmil Constantinescu + ts - timestepping context 792f68a32c8SEmil Constantinescu - rktype - type of RK-scheme 793f68a32c8SEmil Constantinescu 794f68a32c8SEmil Constantinescu Level: intermediate 795f68a32c8SEmil Constantinescu 796f68a32c8SEmil Constantinescu .seealso: TSRKGetType(), TSRK, TSRK2, TSRK3, TSRKPRSSP2, TSRK4, TSRK5 797f68a32c8SEmil Constantinescu @*/ 798f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype) 799f68a32c8SEmil Constantinescu { 800f68a32c8SEmil Constantinescu PetscErrorCode ierr; 801f68a32c8SEmil Constantinescu 802f68a32c8SEmil Constantinescu PetscFunctionBegin; 803f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 804f68a32c8SEmil Constantinescu ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr); 805f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 806f68a32c8SEmil Constantinescu } 807f68a32c8SEmil Constantinescu 808f68a32c8SEmil Constantinescu #undef __FUNCT__ 809f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKGetType" 810f68a32c8SEmil Constantinescu /*@C 811f68a32c8SEmil Constantinescu TSRKGetType - Get the type of RK scheme 812f68a32c8SEmil Constantinescu 813f68a32c8SEmil Constantinescu Logically collective 814f68a32c8SEmil Constantinescu 815f68a32c8SEmil Constantinescu Input Parameter: 816f68a32c8SEmil Constantinescu . ts - timestepping context 817f68a32c8SEmil Constantinescu 818f68a32c8SEmil Constantinescu Output Parameter: 819f68a32c8SEmil Constantinescu . rktype - type of RK-scheme 820f68a32c8SEmil Constantinescu 821f68a32c8SEmil Constantinescu Level: intermediate 822f68a32c8SEmil Constantinescu 823f68a32c8SEmil Constantinescu .seealso: TSRKGetType() 824f68a32c8SEmil Constantinescu @*/ 825f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype) 826f68a32c8SEmil Constantinescu { 827f68a32c8SEmil Constantinescu PetscErrorCode ierr; 828f68a32c8SEmil Constantinescu 829f68a32c8SEmil Constantinescu PetscFunctionBegin; 830f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 831f68a32c8SEmil Constantinescu ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr); 832f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 833f68a32c8SEmil Constantinescu } 834f68a32c8SEmil Constantinescu 835f68a32c8SEmil Constantinescu #undef __FUNCT__ 836f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKGetType_RK" 837f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype) 838f68a32c8SEmil Constantinescu { 839f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 840f68a32c8SEmil Constantinescu PetscErrorCode ierr; 841f68a32c8SEmil Constantinescu 842f68a32c8SEmil Constantinescu PetscFunctionBegin; 843f68a32c8SEmil Constantinescu if (!rk->tableau) { 844f68a32c8SEmil Constantinescu ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr); 845f68a32c8SEmil Constantinescu } 846f68a32c8SEmil Constantinescu *rktype = rk->tableau->name; 847f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 848f68a32c8SEmil Constantinescu } 849f68a32c8SEmil Constantinescu #undef __FUNCT__ 850f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKSetType_RK" 851f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype) 852f68a32c8SEmil Constantinescu { 853f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 854f68a32c8SEmil Constantinescu PetscErrorCode ierr; 855f68a32c8SEmil Constantinescu PetscBool match; 856f68a32c8SEmil Constantinescu RKTableauLink link; 857f68a32c8SEmil Constantinescu 858f68a32c8SEmil Constantinescu PetscFunctionBegin; 859f68a32c8SEmil Constantinescu if (rk->tableau) { 860f68a32c8SEmil Constantinescu ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr); 861f68a32c8SEmil Constantinescu if (match) PetscFunctionReturn(0); 862f68a32c8SEmil Constantinescu } 863f68a32c8SEmil Constantinescu for (link = RKTableauList; link; link=link->next) { 864f68a32c8SEmil Constantinescu ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr); 865f68a32c8SEmil Constantinescu if (match) { 866f68a32c8SEmil Constantinescu ierr = TSReset_RK(ts);CHKERRQ(ierr); 867f68a32c8SEmil Constantinescu rk->tableau = &link->tab; 868f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 869f68a32c8SEmil Constantinescu } 870f68a32c8SEmil Constantinescu } 871f68a32c8SEmil Constantinescu SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype); 872e4dd521cSBarry Smith PetscFunctionReturn(0); 873e4dd521cSBarry Smith } 874e4dd521cSBarry Smith 875ff22ae23SHong Zhang #undef __FUNCT__ 876ff22ae23SHong Zhang #define __FUNCT__ "TSGetStages_RK" 877ff22ae23SHong Zhang static PetscErrorCode TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y) 878ff22ae23SHong Zhang { 879ff22ae23SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 880ff22ae23SHong Zhang 881ff22ae23SHong Zhang PetscFunctionBegin; 882ff22ae23SHong Zhang *ns = rk->tableau->s; 883d2daff3dSHong Zhang if(Y) *Y = rk->Y; 884ff22ae23SHong Zhang PetscFunctionReturn(0); 885ff22ae23SHong Zhang } 886ff22ae23SHong Zhang 887ff22ae23SHong Zhang 888e4dd521cSBarry Smith /* ------------------------------------------------------------ */ 889ebe3b25bSBarry Smith /*MC 890f68a32c8SEmil Constantinescu TSRK - ODE and DAE solver using Runge-Kutta schemes 891ebe3b25bSBarry Smith 892f68a32c8SEmil Constantinescu The user should provide the right hand side of the equation 893f68a32c8SEmil Constantinescu using TSSetRHSFunction(). 894ebe3b25bSBarry Smith 895f68a32c8SEmil Constantinescu Notes: 896f68a32c8SEmil Constantinescu The default is TSRK3, it can be changed with TSRKSetType() or -ts_rk_type 897ebe3b25bSBarry Smith 898d41222bbSBarry Smith Level: beginner 899d41222bbSBarry Smith 900f68a32c8SEmil Constantinescu .seealso: TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3, 901f68a32c8SEmil Constantinescu TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister() 902ebe3b25bSBarry Smith 903ebe3b25bSBarry Smith M*/ 904e4dd521cSBarry Smith #undef __FUNCT__ 9055f70b478SJed Brown #define __FUNCT__ "TSCreate_RK" 9068cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts) 907e4dd521cSBarry Smith { 908f68a32c8SEmil Constantinescu TS_RK *th; 909dfbe8321SBarry Smith PetscErrorCode ierr; 910e4dd521cSBarry Smith 911e4dd521cSBarry Smith PetscFunctionBegin; 912f68a32c8SEmil Constantinescu #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 913f68a32c8SEmil Constantinescu ierr = TSRKInitializePackage();CHKERRQ(ierr); 914f68a32c8SEmil Constantinescu #endif 915f68a32c8SEmil Constantinescu 916f68a32c8SEmil Constantinescu ts->ops->reset = TSReset_RK; 9175f70b478SJed Brown ts->ops->destroy = TSDestroy_RK; 9185f70b478SJed Brown ts->ops->view = TSView_RK; 919f68a32c8SEmil Constantinescu ts->ops->load = TSLoad_RK; 920f68a32c8SEmil Constantinescu ts->ops->setup = TSSetUp_RK; 921f68a32c8SEmil Constantinescu ts->ops->step = TSStep_RK; 922f68a32c8SEmil Constantinescu ts->ops->interpolate = TSInterpolate_RK; 923f68a32c8SEmil Constantinescu ts->ops->evaluatestep = TSEvaluateStep_RK; 924f68a32c8SEmil Constantinescu ts->ops->setfromoptions = TSSetFromOptions_RK; 925ff22ae23SHong Zhang ts->ops->getstages = TSGetStages_RK; 926d2daff3dSHong Zhang ts->ops->stepadj = TSStepAdj_RK; 927e4dd521cSBarry Smith 928b00a9115SJed Brown ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 929f68a32c8SEmil Constantinescu ts->data = (void*)th; 930e4dd521cSBarry Smith 931f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr); 932f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr); 933e4dd521cSBarry Smith PetscFunctionReturn(0); 934e4dd521cSBarry Smith } 935