1e4dd521cSBarry Smith /* 22e698f8bSDebojyoti Ghosh Code for time stepping with the Runge-Kutta method 3f68a32c8SEmil Constantinescu 4f68a32c8SEmil Constantinescu Notes: 5f68a32c8SEmil Constantinescu The general system is written as 6f68a32c8SEmil Constantinescu 72e698f8bSDebojyoti Ghosh Udot = F(t,U) 8f68a32c8SEmil Constantinescu 9e4dd521cSBarry Smith */ 10af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 11f68a32c8SEmil Constantinescu #include <petscdm.h> 12f68a32c8SEmil Constantinescu 13484bcdc7SDebojyoti Ghosh static TSRKType TSRKDefault = TSRK3BS; 14f68a32c8SEmil Constantinescu static PetscBool TSRKRegisterAllCalled; 15f68a32c8SEmil Constantinescu static PetscBool TSRKPackageInitialized; 16f68a32c8SEmil Constantinescu 17f68a32c8SEmil Constantinescu typedef struct _RKTableau *RKTableau; 18f68a32c8SEmil Constantinescu struct _RKTableau { 19f68a32c8SEmil Constantinescu char *name; 20d760c35bSDebojyoti Ghosh PetscInt order; /* Classical approximation order of the method i */ 21f68a32c8SEmil Constantinescu PetscInt s; /* Number of stages */ 22d760c35bSDebojyoti Ghosh PetscBool FSAL; /* flag to indicate if tableau is FSAL */ 23f68a32c8SEmil Constantinescu PetscInt pinterp; /* Interpolation order */ 24f68a32c8SEmil Constantinescu PetscReal *A,*b,*c; /* Tableau */ 25f68a32c8SEmil Constantinescu PetscReal *bembed; /* Embedded formula of order one less (order-1) */ 26f68a32c8SEmil Constantinescu PetscReal *binterp; /* Dense output formula */ 27f68a32c8SEmil Constantinescu PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */ 28f68a32c8SEmil Constantinescu }; 29f68a32c8SEmil Constantinescu typedef struct _RKTableauLink *RKTableauLink; 30f68a32c8SEmil Constantinescu struct _RKTableauLink { 31f68a32c8SEmil Constantinescu struct _RKTableau tab; 32f68a32c8SEmil Constantinescu RKTableauLink next; 33f68a32c8SEmil Constantinescu }; 34f68a32c8SEmil Constantinescu static RKTableauLink RKTableauList; 35e4dd521cSBarry Smith 36e4dd521cSBarry Smith typedef struct { 37f68a32c8SEmil Constantinescu RKTableau tableau; 38f68a32c8SEmil Constantinescu Vec *Y; /* States computed during the step */ 39f68a32c8SEmil Constantinescu Vec *YdotRHS; /* Function evaluations for the non-stiff part */ 40ad8e2604SHong Zhang Vec *VecDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage */ 41ad8e2604SHong Zhang Vec *VecDeltaMu; /* Increment of the adjoint sensitivity w.r.t P at stage */ 42b18ea86cSHong Zhang Vec *VecSensiTemp; /* Vector to be timed with Jacobian transpose */ 430f7a1166SHong Zhang Vec VecCostIntegral0; /* backup for roll-backs due to events */ 44f68a32c8SEmil Constantinescu PetscScalar *work; /* Scalar work */ 45f68a32c8SEmil Constantinescu PetscReal stage_time; 46f68a32c8SEmil Constantinescu TSStepStatus status; 470f7a1166SHong Zhang PetscReal ptime; 480f7a1166SHong Zhang PetscReal time_step; 495f70b478SJed Brown } TS_RK; 50e4dd521cSBarry Smith 51f68a32c8SEmil Constantinescu /*MC 52f68a32c8SEmil Constantinescu TSRK1 - First order forward Euler scheme. 53262ff7bbSBarry Smith 54f68a32c8SEmil Constantinescu This method has one stage. 55f68a32c8SEmil Constantinescu 56f68a32c8SEmil Constantinescu Level: advanced 57f68a32c8SEmil Constantinescu 58f68a32c8SEmil Constantinescu .seealso: TSRK 59f68a32c8SEmil Constantinescu M*/ 60f68a32c8SEmil Constantinescu /*MC 612109b73fSDebojyoti Ghosh TSRK2A - Second order RK scheme. 62f68a32c8SEmil Constantinescu 63f68a32c8SEmil Constantinescu This method has two stages. 64f68a32c8SEmil Constantinescu 65f68a32c8SEmil Constantinescu Level: advanced 66f68a32c8SEmil Constantinescu 67f68a32c8SEmil Constantinescu .seealso: TSRK 68f68a32c8SEmil Constantinescu M*/ 69f68a32c8SEmil Constantinescu /*MC 70f68a32c8SEmil Constantinescu TSRK3 - Third order RK scheme. 71f68a32c8SEmil Constantinescu 72f68a32c8SEmil Constantinescu This method has three stages. 73f68a32c8SEmil Constantinescu 74f68a32c8SEmil Constantinescu Level: advanced 75f68a32c8SEmil Constantinescu 76f68a32c8SEmil Constantinescu .seealso: TSRK 77f68a32c8SEmil Constantinescu M*/ 78f68a32c8SEmil Constantinescu /*MC 792109b73fSDebojyoti Ghosh TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method. 802109b73fSDebojyoti Ghosh 812109b73fSDebojyoti Ghosh This method has four stages. 822109b73fSDebojyoti Ghosh 832109b73fSDebojyoti Ghosh Level: advanced 842109b73fSDebojyoti Ghosh 852109b73fSDebojyoti Ghosh .seealso: TSRK 862109b73fSDebojyoti Ghosh M*/ 872109b73fSDebojyoti Ghosh /*MC 88f68a32c8SEmil Constantinescu TSRK4 - Fourth order RK scheme. 89f68a32c8SEmil Constantinescu 902e698f8bSDebojyoti Ghosh This is the classical Runge-Kutta method with four stages. 91f68a32c8SEmil Constantinescu 92f68a32c8SEmil Constantinescu Level: advanced 93f68a32c8SEmil Constantinescu 94f68a32c8SEmil Constantinescu .seealso: TSRK 95f68a32c8SEmil Constantinescu M*/ 96f68a32c8SEmil Constantinescu /*MC 972e698f8bSDebojyoti Ghosh TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method. 98f68a32c8SEmil Constantinescu 99f68a32c8SEmil Constantinescu This method has six stages. 100f68a32c8SEmil Constantinescu 101f68a32c8SEmil Constantinescu Level: advanced 102f68a32c8SEmil Constantinescu 103f68a32c8SEmil Constantinescu .seealso: TSRK 104f68a32c8SEmil Constantinescu M*/ 1052109b73fSDebojyoti Ghosh /*MC 1062e698f8bSDebojyoti Ghosh TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method. 1072109b73fSDebojyoti Ghosh 1082109b73fSDebojyoti Ghosh This method has seven stages. 1092109b73fSDebojyoti Ghosh 1102109b73fSDebojyoti Ghosh Level: advanced 1112109b73fSDebojyoti Ghosh 1122109b73fSDebojyoti Ghosh .seealso: TSRK 1132109b73fSDebojyoti Ghosh M*/ 114262ff7bbSBarry Smith 115f68a32c8SEmil Constantinescu /*@C 116f68a32c8SEmil Constantinescu TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK 117262ff7bbSBarry Smith 118f68a32c8SEmil Constantinescu Not Collective, but should be called by all processes which will need the schemes to be registered 119262ff7bbSBarry Smith 120f68a32c8SEmil Constantinescu Level: advanced 121262ff7bbSBarry Smith 122f68a32c8SEmil Constantinescu .keywords: TS, TSRK, register, all 123262ff7bbSBarry Smith 124f68a32c8SEmil Constantinescu .seealso: TSRKRegisterDestroy() 125262ff7bbSBarry Smith @*/ 126f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void) 127262ff7bbSBarry Smith { 1284ac538c5SBarry Smith PetscErrorCode ierr; 129262ff7bbSBarry Smith 130262ff7bbSBarry Smith PetscFunctionBegin; 131f68a32c8SEmil Constantinescu if (TSRKRegisterAllCalled) PetscFunctionReturn(0); 132f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_TRUE; 133e4dd521cSBarry Smith 134e4dd521cSBarry Smith { 135f68a32c8SEmil Constantinescu const PetscReal 136f68a32c8SEmil Constantinescu A[1][1] = {{0.0}}, 137f68a32c8SEmil Constantinescu b[1] = {1.0}; 138f68a32c8SEmil Constantinescu ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,1,b);CHKERRQ(ierr); 139e4dd521cSBarry Smith } 140db046809SBarry Smith { 141f68a32c8SEmil Constantinescu const PetscReal 142f68a32c8SEmil Constantinescu A[2][2] = {{0.0,0.0}, 143f68a32c8SEmil Constantinescu {1.0,0.0}}, 144f68a32c8SEmil Constantinescu b[2] = {0.5,0.5}, 145f68a32c8SEmil Constantinescu bembed[2] = {1.0,0}; 1468886ffe3SEmil Constantinescu ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,1,b);CHKERRQ(ierr); 147db046809SBarry Smith } 148f68a32c8SEmil Constantinescu { 149f68a32c8SEmil Constantinescu const PetscReal 150f68a32c8SEmil Constantinescu A[3][3] = {{0,0,0}, 151f68a32c8SEmil Constantinescu {2.0/3.0,0,0}, 152f68a32c8SEmil Constantinescu {-1.0/3.0,1.0,0}}, 153f68a32c8SEmil Constantinescu b[3] = {0.25,0.5,0.25}; 1548886ffe3SEmil Constantinescu ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,1,b);CHKERRQ(ierr); 155fdefd5e5SDebojyoti Ghosh } 156fdefd5e5SDebojyoti Ghosh { 157fdefd5e5SDebojyoti Ghosh const PetscReal 158fdefd5e5SDebojyoti Ghosh A[4][4] = {{0,0,0,0}, 159fdefd5e5SDebojyoti Ghosh {1.0/2.0,0,0,0}, 160fdefd5e5SDebojyoti Ghosh {0,3.0/4.0,0,0}, 161fdefd5e5SDebojyoti Ghosh {2.0/9.0,1.0/3.0,4.0/9.0,0}}, 162fdefd5e5SDebojyoti Ghosh b[4] = {2.0/9.0,1.0/3.0,4.0/9.0,0}, 163fdefd5e5SDebojyoti Ghosh bembed[4] = {7.0/24.0,1.0/4.0,1.0/3.0,1.0/8.0}; 1648886ffe3SEmil Constantinescu ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,1,b);CHKERRQ(ierr); 165db046809SBarry Smith } 166f68a32c8SEmil Constantinescu { 167f68a32c8SEmil Constantinescu const PetscReal 168f68a32c8SEmil Constantinescu A[4][4] = {{0,0,0,0}, 169f68a32c8SEmil Constantinescu {0.5,0,0,0}, 170f68a32c8SEmil Constantinescu {0,0.5,0,0}, 171f68a32c8SEmil Constantinescu {0,0,1.0,0}}, 172f68a32c8SEmil Constantinescu b[4] = {1.0/6.0,1.0/3.0,1.0/3.0,1.0/6.0}; 1738886ffe3SEmil Constantinescu ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,1,b);CHKERRQ(ierr); 174db046809SBarry Smith } 175f68a32c8SEmil Constantinescu { 176f68a32c8SEmil Constantinescu const PetscReal 177f68a32c8SEmil Constantinescu A[6][6] = {{0,0,0,0,0,0}, 178f68a32c8SEmil Constantinescu {0.25,0,0,0,0,0}, 179f68a32c8SEmil Constantinescu {3.0/32.0,9.0/32.0,0,0,0,0}, 180f68a32c8SEmil Constantinescu {1932.0/2197.0,-7200.0/2197.0,7296.0/2197.0,0,0,0}, 181f68a32c8SEmil Constantinescu {439.0/216.0,-8.0,3680.0/513.0,-845.0/4104.0,0,0}, 182f68a32c8SEmil Constantinescu {-8.0/27.0,2.0,-3544.0/2565.0,1859.0/4104.0,-11.0/40.0,0}}, 183f68a32c8SEmil Constantinescu b[6] = {16.0/135.0,0,6656.0/12825.0,28561.0/56430.0,-9.0/50.0,2.0/55.0}, 184f68a32c8SEmil Constantinescu bembed[6] = {25.0/216.0,0,1408.0/2565.0,2197.0/4104.0,-1.0/5.0,0}; 1858886ffe3SEmil Constantinescu ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,1,b);CHKERRQ(ierr); 186fdefd5e5SDebojyoti Ghosh } 187fdefd5e5SDebojyoti Ghosh { 188fdefd5e5SDebojyoti Ghosh const PetscReal 189fdefd5e5SDebojyoti Ghosh A[7][7] = {{0,0,0,0,0,0,0}, 190fdefd5e5SDebojyoti Ghosh {1.0/5.0,0,0,0,0,0,0}, 191fdefd5e5SDebojyoti Ghosh {3.0/40.0,9.0/40.0,0,0,0,0,0}, 192fdefd5e5SDebojyoti Ghosh {44.0/45.0,-56.0/15.0,32.0/9.0,0,0,0,0}, 193fdefd5e5SDebojyoti Ghosh {19372.0/6561.0,-25360.0/2187.0,64448.0/6561.0,-212.0/729.0,0,0,0}, 194fdefd5e5SDebojyoti Ghosh {9017.0/3168.0,-355.0/33.0,46732.0/5247.0,49.0/176.0,-5103.0/18656.0,0,0}, 195fdefd5e5SDebojyoti Ghosh {35.0/384.0,0,500.0/1113.0,125.0/192.0,-2187.0/6784.0,11.0/84.0,0}}, 196fdefd5e5SDebojyoti 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}, 197fdefd5e5SDebojyoti 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}; 1988886ffe3SEmil Constantinescu ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,1,b);CHKERRQ(ierr); 199f68a32c8SEmil Constantinescu } 200db046809SBarry Smith PetscFunctionReturn(0); 201db046809SBarry Smith } 202db046809SBarry Smith 203f68a32c8SEmil Constantinescu /*@C 204f68a32c8SEmil Constantinescu TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister(). 205f68a32c8SEmil Constantinescu 206f68a32c8SEmil Constantinescu Not Collective 207f68a32c8SEmil Constantinescu 208f68a32c8SEmil Constantinescu Level: advanced 209f68a32c8SEmil Constantinescu 210f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy 211f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll() 212f68a32c8SEmil Constantinescu @*/ 213f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void) 214e4dd521cSBarry Smith { 215f68a32c8SEmil Constantinescu PetscErrorCode ierr; 216f68a32c8SEmil Constantinescu RKTableauLink link; 217f68a32c8SEmil Constantinescu 218f68a32c8SEmil Constantinescu PetscFunctionBegin; 219f68a32c8SEmil Constantinescu while ((link = RKTableauList)) { 220f68a32c8SEmil Constantinescu RKTableau t = &link->tab; 221f68a32c8SEmil Constantinescu RKTableauList = link->next; 222f68a32c8SEmil Constantinescu ierr = PetscFree3(t->A,t->b,t->c); CHKERRQ(ierr); 223f68a32c8SEmil Constantinescu ierr = PetscFree (t->bembed); CHKERRQ(ierr); 224f68a32c8SEmil Constantinescu ierr = PetscFree (t->binterp); CHKERRQ(ierr); 225f68a32c8SEmil Constantinescu ierr = PetscFree (t->name); CHKERRQ(ierr); 226f68a32c8SEmil Constantinescu ierr = PetscFree (link); CHKERRQ(ierr); 227f68a32c8SEmil Constantinescu } 228f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_FALSE; 229f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 230f68a32c8SEmil Constantinescu } 231f68a32c8SEmil Constantinescu 232f68a32c8SEmil Constantinescu /*@C 233f68a32c8SEmil Constantinescu TSRKInitializePackage - This function initializes everything in the TSRK package. It is called 234f68a32c8SEmil Constantinescu from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RK() 235f68a32c8SEmil Constantinescu when using static libraries. 236f68a32c8SEmil Constantinescu 237f68a32c8SEmil Constantinescu Level: developer 238f68a32c8SEmil Constantinescu 239f68a32c8SEmil Constantinescu .keywords: TS, TSRK, initialize, package 240f68a32c8SEmil Constantinescu .seealso: PetscInitialize() 241f68a32c8SEmil Constantinescu @*/ 242f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void) 243f68a32c8SEmil Constantinescu { 244186e87acSLisandro Dalcin PetscErrorCode ierr; 245e4dd521cSBarry Smith 246e4dd521cSBarry Smith PetscFunctionBegin; 247f68a32c8SEmil Constantinescu if (TSRKPackageInitialized) PetscFunctionReturn(0); 248f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_TRUE; 249f68a32c8SEmil Constantinescu ierr = TSRKRegisterAll();CHKERRQ(ierr); 250f68a32c8SEmil Constantinescu ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr); 251f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 252f68a32c8SEmil Constantinescu } 253186e87acSLisandro Dalcin 254f68a32c8SEmil Constantinescu /*@C 255f68a32c8SEmil Constantinescu TSRKFinalizePackage - This function destroys everything in the TSRK package. It is 256f68a32c8SEmil Constantinescu called from PetscFinalize(). 257186e87acSLisandro Dalcin 258f68a32c8SEmil Constantinescu Level: developer 259f68a32c8SEmil Constantinescu 260f68a32c8SEmil Constantinescu .keywords: Petsc, destroy, package 261f68a32c8SEmil Constantinescu .seealso: PetscFinalize() 262f68a32c8SEmil Constantinescu @*/ 263f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void) 264f68a32c8SEmil Constantinescu { 265f68a32c8SEmil Constantinescu PetscErrorCode ierr; 266f68a32c8SEmil Constantinescu 267f68a32c8SEmil Constantinescu PetscFunctionBegin; 268f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_FALSE; 269f68a32c8SEmil Constantinescu ierr = TSRKRegisterDestroy();CHKERRQ(ierr); 270f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 271f68a32c8SEmil Constantinescu } 272f68a32c8SEmil Constantinescu 273f68a32c8SEmil Constantinescu /*@C 274f68a32c8SEmil Constantinescu TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 275f68a32c8SEmil Constantinescu 276f68a32c8SEmil Constantinescu Not Collective, but the same schemes should be registered on all processes on which they will be used 277f68a32c8SEmil Constantinescu 278f68a32c8SEmil Constantinescu Input Parameters: 279f68a32c8SEmil Constantinescu + name - identifier for method 280f68a32c8SEmil Constantinescu . order - approximation order of method 281f68a32c8SEmil Constantinescu . s - number of stages, this is the dimension of the matrices below 282f68a32c8SEmil Constantinescu . A - stage coefficients (dimension s*s, row-major) 283f68a32c8SEmil Constantinescu . b - step completion table (dimension s; NULL to use last row of A) 284f68a32c8SEmil Constantinescu . c - abscissa (dimension s; NULL to use row sums of A) 285f68a32c8SEmil Constantinescu . bembed - completion table for embedded method (dimension s; NULL if not available) 286f68a32c8SEmil Constantinescu . pinterp - Order of the interpolation scheme, equal to the number of columns of binterp 287f68a32c8SEmil Constantinescu - binterp - Coefficients of the interpolation formula (dimension s*pinterp; NULL to reuse binterpt) 288f68a32c8SEmil Constantinescu 289f68a32c8SEmil Constantinescu Notes: 290f68a32c8SEmil Constantinescu Several RK methods are provided, this function is only needed to create new methods. 291f68a32c8SEmil Constantinescu 292f68a32c8SEmil Constantinescu Level: advanced 293f68a32c8SEmil Constantinescu 294f68a32c8SEmil Constantinescu .keywords: TS, register 295f68a32c8SEmil Constantinescu 296f68a32c8SEmil Constantinescu .seealso: TSRK 297f68a32c8SEmil Constantinescu @*/ 298f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s, 299f68a32c8SEmil Constantinescu const PetscReal A[],const PetscReal b[],const PetscReal c[], 300f68a32c8SEmil Constantinescu const PetscReal bembed[], 301f68a32c8SEmil Constantinescu PetscInt pinterp,const PetscReal binterp[]) 302f68a32c8SEmil Constantinescu { 303f68a32c8SEmil Constantinescu PetscErrorCode ierr; 304f68a32c8SEmil Constantinescu RKTableauLink link; 305f68a32c8SEmil Constantinescu RKTableau t; 306f68a32c8SEmil Constantinescu PetscInt i,j; 307f68a32c8SEmil Constantinescu 308f68a32c8SEmil Constantinescu PetscFunctionBegin; 30995dccacaSBarry Smith ierr = PetscNew(&link);CHKERRQ(ierr); 310f68a32c8SEmil Constantinescu t = &link->tab; 311f68a32c8SEmil Constantinescu ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 312f68a32c8SEmil Constantinescu t->order = order; 313f68a32c8SEmil Constantinescu t->s = s; 314dcca6d9dSJed Brown ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr); 315f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 316f68a32c8SEmil Constantinescu if (b) { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); } 317f68a32c8SEmil Constantinescu else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i]; 318f68a32c8SEmil Constantinescu if (c) { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); } 319f68a32c8SEmil 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]; 320d760c35bSDebojyoti Ghosh t->FSAL = PETSC_TRUE; 321d760c35bSDebojyoti Ghosh for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE; 322f68a32c8SEmil Constantinescu if (bembed) { 323785e854fSJed Brown ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr); 324f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr); 325f68a32c8SEmil Constantinescu } 326f68a32c8SEmil Constantinescu 327f68a32c8SEmil Constantinescu t->pinterp = pinterp; 328785e854fSJed Brown ierr = PetscMalloc1(s*pinterp,&t->binterp);CHKERRQ(ierr); 329f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->binterp,binterp,s*pinterp*sizeof(binterp[0]));CHKERRQ(ierr); 330f68a32c8SEmil Constantinescu link->next = RKTableauList; 331f68a32c8SEmil Constantinescu RKTableauList = link; 332f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 333f68a32c8SEmil Constantinescu } 334f68a32c8SEmil Constantinescu 335e4dd521cSBarry Smith /* 336f68a32c8SEmil Constantinescu The step completion formula is 337e4dd521cSBarry Smith 338f68a32c8SEmil Constantinescu x1 = x0 + h b^T YdotRHS 339f68a32c8SEmil Constantinescu 340f68a32c8SEmil Constantinescu This function can be called before or after ts->vec_sol has been updated. 341f68a32c8SEmil Constantinescu Suppose we have a completion formula (b) and an embedded formula (be) of different order. 342f68a32c8SEmil Constantinescu We can write 343f68a32c8SEmil Constantinescu 344f68a32c8SEmil Constantinescu x1e = x0 + h be^T YdotRHS 345f68a32c8SEmil Constantinescu = x1 - h b^T YdotRHS + h be^T YdotRHS 346f68a32c8SEmil Constantinescu = x1 + h (be - b)^T YdotRHS 347f68a32c8SEmil Constantinescu 348f68a32c8SEmil Constantinescu so we can evaluate the method with different order even after the step has been optimistically completed. 349e4dd521cSBarry Smith */ 350f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done) 351f68a32c8SEmil Constantinescu { 352f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 353f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 354f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 355f68a32c8SEmil Constantinescu PetscReal h; 356f68a32c8SEmil Constantinescu PetscInt s = tab->s,j; 357f68a32c8SEmil Constantinescu PetscErrorCode ierr; 358f68a32c8SEmil Constantinescu 359f68a32c8SEmil Constantinescu PetscFunctionBegin; 360f68a32c8SEmil Constantinescu switch (rk->status) { 361f68a32c8SEmil Constantinescu case TS_STEP_INCOMPLETE: 362f68a32c8SEmil Constantinescu case TS_STEP_PENDING: 363f68a32c8SEmil Constantinescu h = ts->time_step; break; 364f68a32c8SEmil Constantinescu case TS_STEP_COMPLETE: 365be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 366f68a32c8SEmil Constantinescu default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 367f68a32c8SEmil Constantinescu } 368f68a32c8SEmil Constantinescu if (order == tab->order) { 369f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { 370f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 371f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->b[j]; 372f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 373f68a32c8SEmil Constantinescu } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 374f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 375f68a32c8SEmil Constantinescu } else if (order == tab->order-1) { 376f68a32c8SEmil Constantinescu if (!tab->bembed) goto unavailable; 377f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (be) */ 378f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 379f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 380f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 381f68a32c8SEmil Constantinescu } else { /* Rollback and re-complete using (be-b) */ 382f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 383f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 384f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 3850f7a1166SHong Zhang if (ts->vec_costintegral && ts->costintegralfwd) { 3860f7a1166SHong Zhang ierr = VecCopy(rk->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr); 3870f7a1166SHong Zhang } 388f68a32c8SEmil Constantinescu } 389f68a32c8SEmil Constantinescu if (done) *done = PETSC_TRUE; 390f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 391f68a32c8SEmil Constantinescu } 392f68a32c8SEmil Constantinescu unavailable: 393f68a32c8SEmil Constantinescu if (done) *done = PETSC_FALSE; 394a7fac7c2SEmil 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); 395f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 396f68a32c8SEmil Constantinescu } 397f68a32c8SEmil Constantinescu 3980f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts) 3990f7a1166SHong Zhang { 4000f7a1166SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 4010f7a1166SHong Zhang RKTableau tab = rk->tableau; 4020f7a1166SHong Zhang const PetscInt s = tab->s; 4030f7a1166SHong Zhang const PetscReal *b = tab->b,*c = tab->c; 4040f7a1166SHong Zhang Vec *Y = rk->Y; 4050f7a1166SHong Zhang PetscInt i; 4060f7a1166SHong Zhang PetscErrorCode ierr; 4070f7a1166SHong Zhang 4080f7a1166SHong Zhang PetscFunctionBegin; 4090f7a1166SHong Zhang /* backup cost integral */ 4100f7a1166SHong Zhang ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr); 4110f7a1166SHong Zhang for (i=s-1; i>=0; i--) { 4120f7a1166SHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 4130f7a1166SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,rk->ptime+rk->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 4140f7a1166SHong Zhang ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 4150f7a1166SHong Zhang } 4160f7a1166SHong Zhang PetscFunctionReturn(0); 4170f7a1166SHong Zhang } 4180f7a1166SHong Zhang 4190f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts) 4200f7a1166SHong Zhang { 4210f7a1166SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 4220f7a1166SHong Zhang RKTableau tab = rk->tableau; 4230f7a1166SHong Zhang const PetscInt s = tab->s; 4240f7a1166SHong Zhang const PetscReal *b = tab->b,*c = tab->c; 4250f7a1166SHong Zhang Vec *Y = rk->Y; 4260f7a1166SHong Zhang PetscInt i; 4270f7a1166SHong Zhang PetscErrorCode ierr; 4280f7a1166SHong Zhang 4290f7a1166SHong Zhang PetscFunctionBegin; 4300f7a1166SHong Zhang for (i=s-1; i>=0; i--) { 4310f7a1166SHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 4320f7a1166SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime-ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 4330f7a1166SHong Zhang ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 4340f7a1166SHong Zhang } 4350f7a1166SHong Zhang PetscFunctionReturn(0); 4360f7a1166SHong Zhang } 4370f7a1166SHong Zhang 438fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts) 439fecfb714SLisandro Dalcin { 440fecfb714SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 441fecfb714SLisandro Dalcin RKTableau tab = rk->tableau; 442fecfb714SLisandro Dalcin const PetscInt s = tab->s; 443fecfb714SLisandro Dalcin const PetscReal *b = tab->b; 444fecfb714SLisandro Dalcin PetscScalar *w = rk->work; 445fecfb714SLisandro Dalcin Vec *YdotRHS = rk->YdotRHS; 446fecfb714SLisandro Dalcin PetscInt j; 447fecfb714SLisandro Dalcin PetscReal h; 448fecfb714SLisandro Dalcin PetscErrorCode ierr; 449fecfb714SLisandro Dalcin 450fecfb714SLisandro Dalcin PetscFunctionBegin; 451fecfb714SLisandro Dalcin switch (rk->status) { 452fecfb714SLisandro Dalcin case TS_STEP_INCOMPLETE: 453fecfb714SLisandro Dalcin case TS_STEP_PENDING: 454fecfb714SLisandro Dalcin h = ts->time_step; break; 455fecfb714SLisandro Dalcin case TS_STEP_COMPLETE: 456fecfb714SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 457fecfb714SLisandro Dalcin default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 458fecfb714SLisandro Dalcin } 459fecfb714SLisandro Dalcin for (j=0; j<s; j++) w[j] = -h*b[j]; 460fecfb714SLisandro Dalcin ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr); 461fecfb714SLisandro Dalcin PetscFunctionReturn(0); 462fecfb714SLisandro Dalcin } 463fecfb714SLisandro Dalcin 464f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts) 465f68a32c8SEmil Constantinescu { 466f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 467f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 468f68a32c8SEmil Constantinescu const PetscInt s = tab->s; 469fecfb714SLisandro Dalcin const PetscReal *A = tab->A,*c = tab->c; 470f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 471f68a32c8SEmil Constantinescu Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 472f68a32c8SEmil Constantinescu TSAdapt adapt; 473fecfb714SLisandro Dalcin PetscInt i,j; 474be5899b3SLisandro Dalcin PetscInt rejections = 0; 475be5899b3SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 476be5899b3SLisandro Dalcin PetscReal next_time_step = ts->time_step; 477f68a32c8SEmil Constantinescu PetscErrorCode ierr; 478f68a32c8SEmil Constantinescu 479f68a32c8SEmil Constantinescu PetscFunctionBegin; 480f68a32c8SEmil Constantinescu rk->status = TS_STEP_INCOMPLETE; 481be5899b3SLisandro Dalcin while (!ts->reason && rk->status != TS_STEP_COMPLETE) { 482be5899b3SLisandro Dalcin PetscReal t = ts->ptime; 483f68a32c8SEmil Constantinescu PetscReal h = ts->time_step; 484f68a32c8SEmil Constantinescu for (i=0; i<s; i++) { 4859be3e283SDebojyoti Ghosh rk->stage_time = t + h*c[i]; 4869be3e283SDebojyoti Ghosh ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr); 487f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 488f68a32c8SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 489f68a32c8SEmil Constantinescu ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 4909be3e283SDebojyoti Ghosh ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 491f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 492be5899b3SLisandro Dalcin ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr); 493fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 494f68a32c8SEmil Constantinescu ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 495f68a32c8SEmil Constantinescu } 496be5899b3SLisandro Dalcin 497be5899b3SLisandro Dalcin rk->status = TS_STEP_INCOMPLETE; 498f68a32c8SEmil Constantinescu ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 499f68a32c8SEmil Constantinescu rk->status = TS_STEP_PENDING; 500f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 501f68a32c8SEmil Constantinescu ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 502f68a32c8SEmil Constantinescu ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr); 503fecfb714SLisandro Dalcin ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 504be5899b3SLisandro Dalcin rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 505be5899b3SLisandro Dalcin if (!accept) { /* Roll back the current step */ 506fecfb714SLisandro Dalcin ierr = TSRollBack_RK(ts);CHKERRQ(ierr); 507be5899b3SLisandro Dalcin ts->time_step = next_time_step; 508be5899b3SLisandro Dalcin goto reject_step; 509be5899b3SLisandro Dalcin } 510be5899b3SLisandro Dalcin 5110f7a1166SHong Zhang if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation*/ 5120f7a1166SHong Zhang rk->ptime = ts->ptime; 5130f7a1166SHong Zhang rk->time_step = ts->time_step; 514493ed95fSHong Zhang } 515be5899b3SLisandro Dalcin 516f68a32c8SEmil Constantinescu ts->ptime += ts->time_step; 517f68a32c8SEmil Constantinescu ts->time_step = next_time_step; 518f68a32c8SEmil Constantinescu break; 519be5899b3SLisandro Dalcin 520be5899b3SLisandro Dalcin reject_step: 521fecfb714SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 522be5899b3SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 523be5899b3SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 524be5899b3SLisandro 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 530f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts) 531f6a906c0SBarry Smith { 532f6a906c0SBarry Smith TS_RK *rk = (TS_RK*)ts->data; 533be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 534be5899b3SLisandro Dalcin PetscInt s = tab->s; 535f6a906c0SBarry Smith PetscErrorCode ierr; 536f6a906c0SBarry Smith 537f6a906c0SBarry Smith PetscFunctionBegin; 538f6a906c0SBarry Smith if (ts->adjointsetupcalled++) PetscFunctionReturn(0); 539abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr); 540f6a906c0SBarry Smith if(ts->vecs_sensip) { 541abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr); 542f6a906c0SBarry Smith } 543abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr); 544f6a906c0SBarry Smith PetscFunctionReturn(0); 545f6a906c0SBarry Smith } 546f6a906c0SBarry Smith 54742f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts) 548d2daff3dSHong Zhang { 549c235aa8dSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 550c235aa8dSHong Zhang RKTableau tab = rk->tableau; 551c235aa8dSHong Zhang const PetscInt s = tab->s; 552c235aa8dSHong Zhang const PetscReal *A = tab->A,*b = tab->b,*c = tab->c; 553c235aa8dSHong Zhang PetscScalar *w = rk->work; 554ad8e2604SHong Zhang Vec *Y = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu,*VecSensiTemp = rk->VecSensiTemp; 555b18ea86cSHong Zhang PetscInt i,j,nadj; 5563d3eaba7SBarry Smith PetscReal t = ts->ptime; 557d2daff3dSHong Zhang PetscErrorCode ierr; 558cef76868SBarry Smith PetscReal h = ts->time_step; 559cef76868SBarry Smith Mat J,Jp; 560c235aa8dSHong Zhang 561d2daff3dSHong Zhang PetscFunctionBegin; 562c235aa8dSHong Zhang rk->status = TS_STEP_INCOMPLETE; 563c235aa8dSHong Zhang for (i=s-1; i>=0; i--) { 564c235aa8dSHong Zhang rk->stage_time = t + h*(1.0-c[i]); 565abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 566b18ea86cSHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr); 567302440fdSBarry Smith ierr = VecScale(VecSensiTemp[nadj],-h*b[i]);CHKERRQ(ierr); 568c235aa8dSHong Zhang for (j=i+1; j<s; j++) { 569302440fdSBarry Smith ierr = VecAXPY(VecSensiTemp[nadj],-h*A[j*s+i],VecDeltaLam[nadj*s+j]);CHKERRQ(ierr); 570b18ea86cSHong Zhang } 571c235aa8dSHong Zhang } 572ad8e2604SHong Zhang /* Stage values of lambda */ 573c235aa8dSHong Zhang ierr = TSGetRHSJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 574c235aa8dSHong Zhang ierr = TSComputeRHSJacobian(ts,rk->stage_time,Y[i],J,Jp);CHKERRQ(ierr); 575493ed95fSHong Zhang if (ts->vec_costintegral) { 576493ed95fSHong Zhang ierr = TSAdjointComputeDRDYFunction(ts,rk->stage_time,Y[i],ts->vecs_drdy);CHKERRQ(ierr); 577493ed95fSHong Zhang } 578abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 579b18ea86cSHong Zhang ierr = MatMultTranspose(J,VecSensiTemp[nadj],VecDeltaLam[nadj*s+i]);CHKERRQ(ierr); 580493ed95fSHong Zhang if (ts->vec_costintegral) { 581493ed95fSHong Zhang ierr = VecAXPY(VecDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdy[nadj]);CHKERRQ(ierr); 582493ed95fSHong Zhang } 583b18ea86cSHong Zhang } 5846497ce14SHong Zhang 585ad8e2604SHong Zhang /* Stage values of mu */ 5866497ce14SHong Zhang if(ts->vecs_sensip) { 5875bf8c567SBarry Smith ierr = TSAdjointComputeRHSJacobian(ts,rk->stage_time,Y[i],ts->Jacp);CHKERRQ(ierr); 588493ed95fSHong Zhang if (ts->vec_costintegral) { 589493ed95fSHong Zhang ierr = TSAdjointComputeDRDPFunction(ts,rk->stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr); 590493ed95fSHong Zhang } 591493ed95fSHong Zhang 592abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 593ad8e2604SHong Zhang ierr = MatMultTranspose(ts->Jacp,VecSensiTemp[nadj],VecDeltaMu[nadj*s+i]);CHKERRQ(ierr); 594493ed95fSHong Zhang if (ts->vec_costintegral) { 595493ed95fSHong Zhang ierr = VecAXPY(VecDeltaMu[nadj*s+i],-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr); 596493ed95fSHong Zhang } 597ad8e2604SHong Zhang } 598c235aa8dSHong Zhang } 5996497ce14SHong Zhang } 600c235aa8dSHong Zhang 601c235aa8dSHong Zhang for (j=0; j<s; j++) w[j] = 1.0; 602abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 603b18ea86cSHong Zhang ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr); 6046497ce14SHong Zhang if(ts->vecs_sensip) { 605ad8e2604SHong Zhang ierr = VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);CHKERRQ(ierr); 606b18ea86cSHong Zhang } 6076497ce14SHong Zhang } 608c235aa8dSHong Zhang rk->status = TS_STEP_COMPLETE; 609d2daff3dSHong Zhang PetscFunctionReturn(0); 610d2daff3dSHong Zhang } 611d2daff3dSHong Zhang 612f68a32c8SEmil Constantinescu static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X) 613f68a32c8SEmil Constantinescu { 614f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 615f68a32c8SEmil Constantinescu PetscInt s = rk->tableau->s,pinterp = rk->tableau->pinterp,i,j; 616f68a32c8SEmil Constantinescu PetscReal h; 617f68a32c8SEmil Constantinescu PetscReal tt,t; 618f68a32c8SEmil Constantinescu PetscScalar *b; 619f68a32c8SEmil Constantinescu const PetscReal *B = rk->tableau->binterp; 620f68a32c8SEmil Constantinescu PetscErrorCode ierr; 621e4dd521cSBarry Smith 622f68a32c8SEmil Constantinescu PetscFunctionBegin; 623f68a32c8SEmil Constantinescu if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 624be5899b3SLisandro Dalcin 625f68a32c8SEmil Constantinescu switch (rk->status) { 626f68a32c8SEmil Constantinescu case TS_STEP_INCOMPLETE: 627f68a32c8SEmil Constantinescu case TS_STEP_PENDING: 628f68a32c8SEmil Constantinescu h = ts->time_step; 629f68a32c8SEmil Constantinescu t = (itime - ts->ptime)/h; 630f68a32c8SEmil Constantinescu break; 631f68a32c8SEmil Constantinescu case TS_STEP_COMPLETE: 632be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; 633f68a32c8SEmil Constantinescu t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 634f68a32c8SEmil Constantinescu break; 635f68a32c8SEmil Constantinescu default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 636e4dd521cSBarry Smith } 637785e854fSJed Brown ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 638f68a32c8SEmil Constantinescu for (i=0; i<s; i++) b[i] = 0; 639f68a32c8SEmil Constantinescu for (j=0,tt=t; j<pinterp; j++,tt*=t) { 640f68a32c8SEmil Constantinescu for (i=0; i<s; i++) { 641f68a32c8SEmil Constantinescu b[i] += h * B[i*pinterp+j] * tt; 642e4dd521cSBarry Smith } 643f68a32c8SEmil Constantinescu } 644be5899b3SLisandro Dalcin 645f68a32c8SEmil Constantinescu ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr); 646f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr); 647be5899b3SLisandro Dalcin 648f68a32c8SEmil Constantinescu ierr = PetscFree(b);CHKERRQ(ierr); 649e4dd521cSBarry Smith PetscFunctionReturn(0); 650e4dd521cSBarry Smith } 651e4dd521cSBarry Smith 652e4dd521cSBarry Smith /*------------------------------------------------------------*/ 653be5899b3SLisandro Dalcin 654be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts) 655be5899b3SLisandro Dalcin { 656be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 657be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 658be5899b3SLisandro Dalcin PetscErrorCode ierr; 659be5899b3SLisandro Dalcin 660be5899b3SLisandro Dalcin PetscFunctionBegin; 661be5899b3SLisandro Dalcin if (!tab) PetscFunctionReturn(0); 662be5899b3SLisandro Dalcin ierr = PetscFree(rk->work);CHKERRQ(ierr); 663be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr); 664be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr); 665be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr); 666be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr); 667be5899b3SLisandro Dalcin PetscFunctionReturn(0); 668be5899b3SLisandro Dalcin } 669be5899b3SLisandro Dalcin 670277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts) 671e4dd521cSBarry Smith { 6725f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 6736849ba73SBarry Smith PetscErrorCode ierr; 674e4dd521cSBarry Smith 675e4dd521cSBarry Smith PetscFunctionBegin; 676be5899b3SLisandro Dalcin ierr = TSRKTableauReset(ts);CHKERRQ(ierr); 6770f7a1166SHong Zhang ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr); 678abc2977eSBarry Smith ierr = VecDestroyVecs(ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr); 679277b19d0SLisandro Dalcin PetscFunctionReturn(0); 680e4dd521cSBarry Smith } 681277b19d0SLisandro Dalcin 682277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_RK(TS ts) 683277b19d0SLisandro Dalcin { 684277b19d0SLisandro Dalcin PetscErrorCode ierr; 685277b19d0SLisandro Dalcin 686277b19d0SLisandro Dalcin PetscFunctionBegin; 687277b19d0SLisandro Dalcin ierr = TSReset_RK(ts);CHKERRQ(ierr); 688277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 689f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr); 690f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr); 691f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 692f68a32c8SEmil Constantinescu } 693f68a32c8SEmil Constantinescu 694f68a32c8SEmil Constantinescu 695f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx) 696f68a32c8SEmil Constantinescu { 697f68a32c8SEmil Constantinescu PetscFunctionBegin; 698f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 699f68a32c8SEmil Constantinescu } 700f68a32c8SEmil Constantinescu 701f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 702f68a32c8SEmil Constantinescu { 703f68a32c8SEmil Constantinescu PetscFunctionBegin; 704f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 705f68a32c8SEmil Constantinescu } 706f68a32c8SEmil Constantinescu 707f68a32c8SEmil Constantinescu 708f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx) 709f68a32c8SEmil Constantinescu { 710f68a32c8SEmil Constantinescu PetscFunctionBegin; 711f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 712f68a32c8SEmil Constantinescu } 713f68a32c8SEmil Constantinescu 714f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 715f68a32c8SEmil Constantinescu { 716f68a32c8SEmil Constantinescu 717f68a32c8SEmil Constantinescu PetscFunctionBegin; 718f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 719f68a32c8SEmil Constantinescu } 720c235aa8dSHong Zhang /* 721d2daff3dSHong Zhang static PetscErrorCode RKSetAdjCoe(RKTableau tab) 722d2daff3dSHong Zhang { 723d2daff3dSHong Zhang PetscReal *A,*b; 724d2daff3dSHong Zhang PetscInt s,i,j; 725d2daff3dSHong Zhang PetscErrorCode ierr; 726d2daff3dSHong Zhang 727d2daff3dSHong Zhang PetscFunctionBegin; 728d2daff3dSHong Zhang s = tab->s; 729d2daff3dSHong Zhang ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr); 730d2daff3dSHong Zhang 731d2daff3dSHong Zhang for (i=0; i<s; i++) 732d2daff3dSHong Zhang for (j=0; j<s; j++) { 733d2daff3dSHong 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]; 734d2daff3dSHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr); 735d2daff3dSHong Zhang } 736d2daff3dSHong Zhang 737d2daff3dSHong Zhang for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i]; 738d2daff3dSHong Zhang 739d2daff3dSHong Zhang ierr = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 740d2daff3dSHong Zhang ierr = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr); 741d2daff3dSHong Zhang ierr = PetscFree2(A,b);CHKERRQ(ierr); 742d2daff3dSHong Zhang PetscFunctionReturn(0); 743d2daff3dSHong Zhang } 744c235aa8dSHong Zhang */ 745be5899b3SLisandro Dalcin 746be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts) 747be5899b3SLisandro Dalcin { 748be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 749be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 750be5899b3SLisandro Dalcin PetscErrorCode ierr; 751be5899b3SLisandro Dalcin 752be5899b3SLisandro Dalcin PetscFunctionBegin; 753be5899b3SLisandro Dalcin ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr); 754be5899b3SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr); 755be5899b3SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr); 756be5899b3SLisandro Dalcin PetscFunctionReturn(0); 757be5899b3SLisandro Dalcin } 758be5899b3SLisandro Dalcin 759be5899b3SLisandro Dalcin 760f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts) 761f68a32c8SEmil Constantinescu { 762f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 763f68a32c8SEmil Constantinescu PetscErrorCode ierr; 764f68a32c8SEmil Constantinescu DM dm; 765f68a32c8SEmil Constantinescu 766f68a32c8SEmil Constantinescu PetscFunctionBegin; 7673dd83b38SBarry Smith ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr); 768be5899b3SLisandro Dalcin ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr); 7690f7a1166SHong Zhang if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 7700f7a1166SHong Zhang ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr); 7710f7a1166SHong Zhang } 772f68a32c8SEmil Constantinescu ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 773f68a32c8SEmil Constantinescu if (dm) { 774f68a32c8SEmil Constantinescu ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 775f68a32c8SEmil Constantinescu ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 776f68a32c8SEmil Constantinescu } 777e4dd521cSBarry Smith PetscFunctionReturn(0); 778e4dd521cSBarry Smith } 779d2daff3dSHong Zhang 780f6a906c0SBarry Smith 781e4dd521cSBarry Smith /*------------------------------------------------------------*/ 782e4dd521cSBarry Smith 7834416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts) 784e4dd521cSBarry Smith { 785be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 786dfbe8321SBarry Smith PetscErrorCode ierr; 787262ff7bbSBarry Smith 788e4dd521cSBarry Smith PetscFunctionBegin; 789e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr); 790f68a32c8SEmil Constantinescu { 791f68a32c8SEmil Constantinescu RKTableauLink link; 792f68a32c8SEmil Constantinescu PetscInt count,choice; 793f68a32c8SEmil Constantinescu PetscBool flg; 794f68a32c8SEmil Constantinescu const char **namelist; 795f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) ; 796785e854fSJed Brown ierr = PetscMalloc1(count,&namelist);CHKERRQ(ierr); 797f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 798be5899b3SLisandro Dalcin ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr); 799be5899b3SLisandro Dalcin if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);} 800f68a32c8SEmil Constantinescu ierr = PetscFree(namelist);CHKERRQ(ierr); 801f68a32c8SEmil Constantinescu } 802262ff7bbSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 803e4dd521cSBarry Smith PetscFunctionReturn(0); 804e4dd521cSBarry Smith } 805e4dd521cSBarry Smith 806f68a32c8SEmil Constantinescu static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[]) 807f68a32c8SEmil Constantinescu { 808f68a32c8SEmil Constantinescu PetscErrorCode ierr; 809f68a32c8SEmil Constantinescu PetscInt i; 810f68a32c8SEmil Constantinescu size_t left,count; 811f68a32c8SEmil Constantinescu char *p; 812f68a32c8SEmil Constantinescu 813f68a32c8SEmil Constantinescu PetscFunctionBegin; 814f68a32c8SEmil Constantinescu for (i=0,p=buf,left=len; i<n; i++) { 815f68a32c8SEmil Constantinescu ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr); 816f68a32c8SEmil Constantinescu if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer"); 817f68a32c8SEmil Constantinescu left -= count; 818f68a32c8SEmil Constantinescu p += count; 819f68a32c8SEmil Constantinescu *p++ = ' '; 820f68a32c8SEmil Constantinescu } 821f68a32c8SEmil Constantinescu p[i ? 0 : -1] = 0; 822f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 823f68a32c8SEmil Constantinescu } 824f68a32c8SEmil Constantinescu 8255f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer) 826e4dd521cSBarry Smith { 8275f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 8288370ee3bSLisandro Dalcin PetscBool iascii; 829dfbe8321SBarry Smith PetscErrorCode ierr; 8302cdf8120Sasbjorn 831e4dd521cSBarry Smith PetscFunctionBegin; 832251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 8338370ee3bSLisandro Dalcin if (iascii) { 8349c334d8fSLisandro Dalcin RKTableau tab = rk->tableau; 835f68a32c8SEmil Constantinescu TSRKType rktype; 836f68a32c8SEmil Constantinescu char buf[512]; 837f68a32c8SEmil Constantinescu ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr); 838*0c37eb8eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," RK: %s\n",rktype);CHKERRQ(ierr); 839*0c37eb8eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," Order: %D\n",tab->order);CHKERRQ(ierr); 840*0c37eb8eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr); 841f68a32c8SEmil Constantinescu ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 842f68a32c8SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," Abscissa c = %s\n",buf);CHKERRQ(ierr); 8438370ee3bSLisandro Dalcin } 8449c334d8fSLisandro Dalcin if (ts->adapt) {ierr = TSAdaptView(ts->adapt,viewer);CHKERRQ(ierr);} 845f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 846f68a32c8SEmil Constantinescu } 847f68a32c8SEmil Constantinescu 848f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer) 849f68a32c8SEmil Constantinescu { 850f68a32c8SEmil Constantinescu PetscErrorCode ierr; 8519c334d8fSLisandro Dalcin TSAdapt adapt; 852f68a32c8SEmil Constantinescu 853f68a32c8SEmil Constantinescu PetscFunctionBegin; 8549c334d8fSLisandro Dalcin ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 8559c334d8fSLisandro Dalcin ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr); 856f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 857f68a32c8SEmil Constantinescu } 858f68a32c8SEmil Constantinescu 859f68a32c8SEmil Constantinescu /*@C 860f68a32c8SEmil Constantinescu TSRKSetType - Set the type of RK scheme 861f68a32c8SEmil Constantinescu 862f68a32c8SEmil Constantinescu Logically collective 863f68a32c8SEmil Constantinescu 864f68a32c8SEmil Constantinescu Input Parameter: 865f68a32c8SEmil Constantinescu + ts - timestepping context 866f68a32c8SEmil Constantinescu - rktype - type of RK-scheme 867f68a32c8SEmil Constantinescu 868f68a32c8SEmil Constantinescu Level: intermediate 869f68a32c8SEmil Constantinescu 870f68a32c8SEmil Constantinescu .seealso: TSRKGetType(), TSRK, TSRK2, TSRK3, TSRKPRSSP2, TSRK4, TSRK5 871f68a32c8SEmil Constantinescu @*/ 872f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype) 873f68a32c8SEmil Constantinescu { 874f68a32c8SEmil Constantinescu PetscErrorCode ierr; 875f68a32c8SEmil Constantinescu 876f68a32c8SEmil Constantinescu PetscFunctionBegin; 877f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 878f68a32c8SEmil Constantinescu ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr); 879f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 880f68a32c8SEmil Constantinescu } 881f68a32c8SEmil Constantinescu 882f68a32c8SEmil Constantinescu /*@C 883f68a32c8SEmil Constantinescu TSRKGetType - Get the type of RK scheme 884f68a32c8SEmil Constantinescu 885f68a32c8SEmil Constantinescu Logically collective 886f68a32c8SEmil Constantinescu 887f68a32c8SEmil Constantinescu Input Parameter: 888f68a32c8SEmil Constantinescu . ts - timestepping context 889f68a32c8SEmil Constantinescu 890f68a32c8SEmil Constantinescu Output Parameter: 891f68a32c8SEmil Constantinescu . rktype - type of RK-scheme 892f68a32c8SEmil Constantinescu 893f68a32c8SEmil Constantinescu Level: intermediate 894f68a32c8SEmil Constantinescu 895f68a32c8SEmil Constantinescu .seealso: TSRKGetType() 896f68a32c8SEmil Constantinescu @*/ 897f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype) 898f68a32c8SEmil Constantinescu { 899f68a32c8SEmil Constantinescu PetscErrorCode ierr; 900f68a32c8SEmil Constantinescu 901f68a32c8SEmil Constantinescu PetscFunctionBegin; 902f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 903f68a32c8SEmil Constantinescu ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr); 904f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 905f68a32c8SEmil Constantinescu } 906f68a32c8SEmil Constantinescu 907560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype) 908f68a32c8SEmil Constantinescu { 909f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 910f68a32c8SEmil Constantinescu 911f68a32c8SEmil Constantinescu PetscFunctionBegin; 912f68a32c8SEmil Constantinescu *rktype = rk->tableau->name; 913f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 914f68a32c8SEmil Constantinescu } 915560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype) 916f68a32c8SEmil Constantinescu { 917f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 918f68a32c8SEmil Constantinescu PetscErrorCode ierr; 919f68a32c8SEmil Constantinescu PetscBool match; 920f68a32c8SEmil Constantinescu RKTableauLink link; 921f68a32c8SEmil Constantinescu 922f68a32c8SEmil Constantinescu PetscFunctionBegin; 923f68a32c8SEmil Constantinescu if (rk->tableau) { 924f68a32c8SEmil Constantinescu ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr); 925f68a32c8SEmil Constantinescu if (match) PetscFunctionReturn(0); 926f68a32c8SEmil Constantinescu } 927f68a32c8SEmil Constantinescu for (link = RKTableauList; link; link=link->next) { 928f68a32c8SEmil Constantinescu ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr); 929f68a32c8SEmil Constantinescu if (match) { 930be5899b3SLisandro Dalcin if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);} 931f68a32c8SEmil Constantinescu rk->tableau = &link->tab; 932be5899b3SLisandro Dalcin if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);} 933f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 934f68a32c8SEmil Constantinescu } 935f68a32c8SEmil Constantinescu } 936f68a32c8SEmil Constantinescu SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype); 937e4dd521cSBarry Smith PetscFunctionReturn(0); 938e4dd521cSBarry Smith } 939e4dd521cSBarry Smith 940ff22ae23SHong Zhang static PetscErrorCode TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y) 941ff22ae23SHong Zhang { 942ff22ae23SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 943ff22ae23SHong Zhang 944ff22ae23SHong Zhang PetscFunctionBegin; 945ff22ae23SHong Zhang *ns = rk->tableau->s; 946d2daff3dSHong Zhang if(Y) *Y = rk->Y; 947ff22ae23SHong Zhang PetscFunctionReturn(0); 948ff22ae23SHong Zhang } 949ff22ae23SHong Zhang 950ff22ae23SHong Zhang 951e4dd521cSBarry Smith /* ------------------------------------------------------------ */ 952ebe3b25bSBarry Smith /*MC 953f68a32c8SEmil Constantinescu TSRK - ODE and DAE solver using Runge-Kutta schemes 954ebe3b25bSBarry Smith 955f68a32c8SEmil Constantinescu The user should provide the right hand side of the equation 956f68a32c8SEmil Constantinescu using TSSetRHSFunction(). 957ebe3b25bSBarry Smith 958f68a32c8SEmil Constantinescu Notes: 959f68a32c8SEmil Constantinescu The default is TSRK3, it can be changed with TSRKSetType() or -ts_rk_type 960ebe3b25bSBarry Smith 961d41222bbSBarry Smith Level: beginner 962d41222bbSBarry Smith 963f68a32c8SEmil Constantinescu .seealso: TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3, 964f68a32c8SEmil Constantinescu TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister() 965ebe3b25bSBarry Smith 966ebe3b25bSBarry Smith M*/ 9678cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts) 968e4dd521cSBarry Smith { 9691566a47fSLisandro Dalcin TS_RK *rk; 970dfbe8321SBarry Smith PetscErrorCode ierr; 971e4dd521cSBarry Smith 972e4dd521cSBarry Smith PetscFunctionBegin; 973f68a32c8SEmil Constantinescu ierr = TSRKInitializePackage();CHKERRQ(ierr); 974f68a32c8SEmil Constantinescu 975f68a32c8SEmil Constantinescu ts->ops->reset = TSReset_RK; 9765f70b478SJed Brown ts->ops->destroy = TSDestroy_RK; 9775f70b478SJed Brown ts->ops->view = TSView_RK; 978f68a32c8SEmil Constantinescu ts->ops->load = TSLoad_RK; 979f68a32c8SEmil Constantinescu ts->ops->setup = TSSetUp_RK; 98042f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_RK; 981f68a32c8SEmil Constantinescu ts->ops->step = TSStep_RK; 982f68a32c8SEmil Constantinescu ts->ops->interpolate = TSInterpolate_RK; 983f68a32c8SEmil Constantinescu ts->ops->evaluatestep = TSEvaluateStep_RK; 984fecfb714SLisandro Dalcin ts->ops->rollback = TSRollBack_RK; 985f68a32c8SEmil Constantinescu ts->ops->setfromoptions = TSSetFromOptions_RK; 986ff22ae23SHong Zhang ts->ops->getstages = TSGetStages_RK; 98742f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_RK; 988e4dd521cSBarry Smith 9890f7a1166SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_RK; 9900f7a1166SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_RK; 9910f7a1166SHong Zhang 9921566a47fSLisandro Dalcin ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr); 9931566a47fSLisandro Dalcin ts->data = (void*)rk; 994e4dd521cSBarry Smith 995f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr); 996f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr); 997be5899b3SLisandro Dalcin 998be5899b3SLisandro Dalcin ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr); 999e4dd521cSBarry Smith PetscFunctionReturn(0); 1000e4dd521cSBarry Smith } 1001