1e4dd521cSBarry Smith /* 2f68a32c8SEmil Constantinescu Code for timestepping with Runge-Kutta method 3f68a32c8SEmil Constantinescu 4f68a32c8SEmil Constantinescu Notes: 5f68a32c8SEmil Constantinescu The general system is written as 6f68a32c8SEmil Constantinescu 7f68a32c8SEmil Constantinescu F(t,U,Udot) = G(t,U) 8f68a32c8SEmil Constantinescu 9f68a32c8SEmil Constantinescu where F represents the stiff part of the physics and G represents the non-stiff part. 10f68a32c8SEmil Constantinescu 11e4dd521cSBarry Smith */ 12b45d2f2cSJed Brown #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 13f68a32c8SEmil Constantinescu #include <petscdm.h> 14f68a32c8SEmil Constantinescu 15f68a32c8SEmil Constantinescu static TSRKType TSRKDefault = TSRK3; 16f68a32c8SEmil Constantinescu static PetscBool TSRKRegisterAllCalled; 17f68a32c8SEmil Constantinescu static PetscBool TSRKPackageInitialized; 18f68a32c8SEmil Constantinescu static PetscInt explicit_stage_time_id; 19f68a32c8SEmil Constantinescu 20f68a32c8SEmil Constantinescu typedef struct _RKTableau *RKTableau; 21f68a32c8SEmil Constantinescu struct _RKTableau { 22f68a32c8SEmil Constantinescu char *name; 23f68a32c8SEmil Constantinescu PetscInt order; /* Classical approximation order of the method */ 24f68a32c8SEmil Constantinescu PetscInt s; /* Number of stages */ 25f68a32c8SEmil Constantinescu PetscInt pinterp; /* Interpolation order */ 26f68a32c8SEmil Constantinescu PetscReal *A,*b,*c; /* Tableau */ 27f68a32c8SEmil Constantinescu PetscReal *bembed; /* Embedded formula of order one less (order-1) */ 28f68a32c8SEmil Constantinescu PetscReal *binterp; /* Dense output formula */ 29f68a32c8SEmil Constantinescu PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */ 30f68a32c8SEmil Constantinescu }; 31f68a32c8SEmil Constantinescu typedef struct _RKTableauLink *RKTableauLink; 32f68a32c8SEmil Constantinescu struct _RKTableauLink { 33f68a32c8SEmil Constantinescu struct _RKTableau tab; 34f68a32c8SEmil Constantinescu RKTableauLink next; 35f68a32c8SEmil Constantinescu }; 36f68a32c8SEmil Constantinescu static RKTableauLink RKTableauList; 37e4dd521cSBarry Smith 38e4dd521cSBarry Smith typedef struct { 39f68a32c8SEmil Constantinescu RKTableau tableau; 40f68a32c8SEmil Constantinescu Vec *Y; /* States computed during the step */ 41f68a32c8SEmil Constantinescu Vec *YdotRHS; /* Function evaluations for the non-stiff part */ 42f68a32c8SEmil Constantinescu PetscScalar *work; /* Scalar work */ 43f68a32c8SEmil Constantinescu PetscReal stage_time; 44f68a32c8SEmil Constantinescu TSStepStatus status; 455f70b478SJed Brown } TS_RK; 46e4dd521cSBarry Smith 47f68a32c8SEmil Constantinescu /*MC 48f68a32c8SEmil Constantinescu TSRK1 - First order forward Euler scheme. 49262ff7bbSBarry Smith 50f68a32c8SEmil Constantinescu This method has one stage. 51f68a32c8SEmil Constantinescu 52f68a32c8SEmil Constantinescu Level: advanced 53f68a32c8SEmil Constantinescu 54f68a32c8SEmil Constantinescu .seealso: TSRK 55f68a32c8SEmil Constantinescu M*/ 56f68a32c8SEmil Constantinescu /*MC 57f68a32c8SEmil Constantinescu TSRK2 - Second order RK scheme. 58f68a32c8SEmil Constantinescu 59f68a32c8SEmil Constantinescu This method has two stages. 60f68a32c8SEmil Constantinescu 61f68a32c8SEmil Constantinescu Level: advanced 62f68a32c8SEmil Constantinescu 63f68a32c8SEmil Constantinescu .seealso: TSRK 64f68a32c8SEmil Constantinescu M*/ 65f68a32c8SEmil Constantinescu /*MC 66f68a32c8SEmil Constantinescu TSRK3 - Third order RK scheme. 67f68a32c8SEmil Constantinescu 68f68a32c8SEmil Constantinescu This method has three stages. 69f68a32c8SEmil Constantinescu 70f68a32c8SEmil Constantinescu Level: advanced 71f68a32c8SEmil Constantinescu 72f68a32c8SEmil Constantinescu .seealso: TSRK 73f68a32c8SEmil Constantinescu M*/ 74f68a32c8SEmil Constantinescu /*MC 75f68a32c8SEmil Constantinescu TSRK4 - Fourth order RK scheme. 76f68a32c8SEmil Constantinescu 77f68a32c8SEmil Constantinescu This method has four stages. 78f68a32c8SEmil Constantinescu 79f68a32c8SEmil Constantinescu Level: advanced 80f68a32c8SEmil Constantinescu 81f68a32c8SEmil Constantinescu .seealso: TSRK 82f68a32c8SEmil Constantinescu M*/ 83f68a32c8SEmil Constantinescu /*MC 84f68a32c8SEmil Constantinescu TSRK5F - Fifth order Fehlberg RK scheme with 4th order embedded method. 85f68a32c8SEmil Constantinescu 86f68a32c8SEmil Constantinescu This method has six stages. 87f68a32c8SEmil Constantinescu 88f68a32c8SEmil Constantinescu Level: advanced 89f68a32c8SEmil Constantinescu 90f68a32c8SEmil Constantinescu .seealso: TSRK 91f68a32c8SEmil Constantinescu M*/ 92262ff7bbSBarry Smith 93262ff7bbSBarry Smith #undef __FUNCT__ 94f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKRegisterAll" 95f68a32c8SEmil Constantinescu /*@C 96f68a32c8SEmil Constantinescu TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK 97262ff7bbSBarry Smith 98f68a32c8SEmil Constantinescu Not Collective, but should be called by all processes which will need the schemes to be registered 99262ff7bbSBarry Smith 100f68a32c8SEmil Constantinescu Level: advanced 101262ff7bbSBarry Smith 102f68a32c8SEmil Constantinescu .keywords: TS, TSRK, register, all 103262ff7bbSBarry Smith 104f68a32c8SEmil Constantinescu .seealso: TSRKRegisterDestroy() 105262ff7bbSBarry Smith @*/ 106f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void) 107262ff7bbSBarry Smith { 1084ac538c5SBarry Smith PetscErrorCode ierr; 109262ff7bbSBarry Smith 110262ff7bbSBarry Smith PetscFunctionBegin; 111f68a32c8SEmil Constantinescu if (TSRKRegisterAllCalled) PetscFunctionReturn(0); 112f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_TRUE; 113e4dd521cSBarry Smith 114e4dd521cSBarry Smith { 115f68a32c8SEmil Constantinescu const PetscReal 116f68a32c8SEmil Constantinescu A[1][1] = {{0.0}}, 117f68a32c8SEmil Constantinescu b[1] = {1.0}; 118f68a32c8SEmil Constantinescu ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,1,b);CHKERRQ(ierr); 119e4dd521cSBarry Smith } 120db046809SBarry Smith { 121f68a32c8SEmil Constantinescu const PetscReal 122f68a32c8SEmil Constantinescu A[2][2] = {{0.0,0.0}, 123f68a32c8SEmil Constantinescu {1.0,0.0}}, 124f68a32c8SEmil Constantinescu b[2] = {0.5,0.5}, 125f68a32c8SEmil Constantinescu bembed[2] = {1.0,0}; 126*fdefd5e5SDebojyoti Ghosh ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,2,b);CHKERRQ(ierr); 127db046809SBarry Smith } 128f68a32c8SEmil Constantinescu { 129f68a32c8SEmil Constantinescu const PetscReal 130f68a32c8SEmil Constantinescu A[3][3] = {{0,0,0}, 131f68a32c8SEmil Constantinescu {2.0/3.0,0,0}, 132f68a32c8SEmil Constantinescu {-1.0/3.0,1.0,0}}, 133f68a32c8SEmil Constantinescu b[3] = {0.25,0.5,0.25}; 134*fdefd5e5SDebojyoti Ghosh ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,3,b);CHKERRQ(ierr); 135*fdefd5e5SDebojyoti Ghosh } 136*fdefd5e5SDebojyoti Ghosh { 137*fdefd5e5SDebojyoti Ghosh const PetscReal 138*fdefd5e5SDebojyoti Ghosh A[4][4] = {{0,0,0,0}, 139*fdefd5e5SDebojyoti Ghosh {1.0/2.0,0,0,0}, 140*fdefd5e5SDebojyoti Ghosh {0,3.0/4.0,0,0}, 141*fdefd5e5SDebojyoti Ghosh {2.0/9.0,1.0/3.0,4.0/9.0,0}}, 142*fdefd5e5SDebojyoti Ghosh b[4] = {2.0/9.0,1.0/3.0,4.0/9.0,0}, 143*fdefd5e5SDebojyoti Ghosh bembed[4] = {7.0/24.0,1.0/4.0,1.0/3.0,1.0/8.0}; 144*fdefd5e5SDebojyoti Ghosh ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,3,b);CHKERRQ(ierr); 145db046809SBarry Smith } 146f68a32c8SEmil Constantinescu { 147f68a32c8SEmil Constantinescu const PetscReal 148f68a32c8SEmil Constantinescu A[4][4] = {{0,0,0,0}, 149f68a32c8SEmil Constantinescu {0.5,0,0,0}, 150f68a32c8SEmil Constantinescu {0,0.5,0,0}, 151f68a32c8SEmil Constantinescu {0,0,1.0,0}}, 152f68a32c8SEmil Constantinescu b[4] = {1.0/6.0,1.0/3.0,1.0/3.0,1.0/6.0}; 153*fdefd5e5SDebojyoti Ghosh ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,4,b);CHKERRQ(ierr); 154db046809SBarry Smith } 155f68a32c8SEmil Constantinescu { 156f68a32c8SEmil Constantinescu const PetscReal 157f68a32c8SEmil Constantinescu A[6][6] = {{0,0,0,0,0,0}, 158f68a32c8SEmil Constantinescu {0.25,0,0,0,0,0}, 159f68a32c8SEmil Constantinescu {3.0/32.0,9.0/32.0,0,0,0,0}, 160f68a32c8SEmil Constantinescu {1932.0/2197.0,-7200.0/2197.0,7296.0/2197.0,0,0,0}, 161f68a32c8SEmil Constantinescu {439.0/216.0,-8.0,3680.0/513.0,-845.0/4104.0,0,0}, 162f68a32c8SEmil Constantinescu {-8.0/27.0,2.0,-3544.0/2565.0,1859.0/4104.0,-11.0/40.0,0}}, 163f68a32c8SEmil Constantinescu b[6] = {16.0/135.0,0,6656.0/12825.0,28561.0/56430.0,-9.0/50.0,2.0/55.0}, 164f68a32c8SEmil Constantinescu bembed[6] = {25.0/216.0,0,1408.0/2565.0,2197.0/4104.0,-1.0/5.0,0}; 165*fdefd5e5SDebojyoti Ghosh ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,5,b);CHKERRQ(ierr); 166*fdefd5e5SDebojyoti Ghosh } 167*fdefd5e5SDebojyoti Ghosh { 168*fdefd5e5SDebojyoti Ghosh const PetscReal 169*fdefd5e5SDebojyoti Ghosh A[7][7] = {{0,0,0,0,0,0,0}, 170*fdefd5e5SDebojyoti Ghosh {1.0/5.0,0,0,0,0,0,0}, 171*fdefd5e5SDebojyoti Ghosh {3.0/40.0,9.0/40.0,0,0,0,0,0}, 172*fdefd5e5SDebojyoti Ghosh {44.0/45.0,-56.0/15.0,32.0/9.0,0,0,0,0}, 173*fdefd5e5SDebojyoti Ghosh {19372.0/6561.0,-25360.0/2187.0,64448.0/6561.0,-212.0/729.0,0,0,0}, 174*fdefd5e5SDebojyoti Ghosh {9017.0/3168.0,-355.0/33.0,46732.0/5247.0,49.0/176.0,-5103.0/18656.0,0,0}, 175*fdefd5e5SDebojyoti Ghosh {35.0/384.0,0,500.0/1113.0,125.0/192.0,-2187.0/6784.0,11.0/84.0,0}}, 176*fdefd5e5SDebojyoti 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}, 177*fdefd5e5SDebojyoti 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}; 178*fdefd5e5SDebojyoti Ghosh ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,5,b);CHKERRQ(ierr); 179f68a32c8SEmil Constantinescu } 180db046809SBarry Smith PetscFunctionReturn(0); 181db046809SBarry Smith } 182db046809SBarry Smith 183e4dd521cSBarry Smith #undef __FUNCT__ 184f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKRegisterDestroy" 185f68a32c8SEmil Constantinescu /*@C 186f68a32c8SEmil Constantinescu TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister(). 187f68a32c8SEmil Constantinescu 188f68a32c8SEmil Constantinescu Not Collective 189f68a32c8SEmil Constantinescu 190f68a32c8SEmil Constantinescu Level: advanced 191f68a32c8SEmil Constantinescu 192f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy 193f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll() 194f68a32c8SEmil Constantinescu @*/ 195f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void) 196e4dd521cSBarry Smith { 197f68a32c8SEmil Constantinescu PetscErrorCode ierr; 198f68a32c8SEmil Constantinescu RKTableauLink link; 199f68a32c8SEmil Constantinescu 200f68a32c8SEmil Constantinescu PetscFunctionBegin; 201f68a32c8SEmil Constantinescu while ((link = RKTableauList)) { 202f68a32c8SEmil Constantinescu RKTableau t = &link->tab; 203f68a32c8SEmil Constantinescu RKTableauList = link->next; 204f68a32c8SEmil Constantinescu ierr = PetscFree3(t->A,t->b,t->c); CHKERRQ(ierr); 205f68a32c8SEmil Constantinescu ierr = PetscFree (t->bembed); CHKERRQ(ierr); 206f68a32c8SEmil Constantinescu ierr = PetscFree (t->binterp); CHKERRQ(ierr); 207f68a32c8SEmil Constantinescu ierr = PetscFree (t->name); CHKERRQ(ierr); 208f68a32c8SEmil Constantinescu ierr = PetscFree (link); CHKERRQ(ierr); 209f68a32c8SEmil Constantinescu } 210f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_FALSE; 211f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 212f68a32c8SEmil Constantinescu } 213f68a32c8SEmil Constantinescu 214f68a32c8SEmil Constantinescu #undef __FUNCT__ 215f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKInitializePackage" 216f68a32c8SEmil Constantinescu /*@C 217f68a32c8SEmil Constantinescu TSRKInitializePackage - This function initializes everything in the TSRK package. It is called 218f68a32c8SEmil Constantinescu from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RK() 219f68a32c8SEmil Constantinescu when using static libraries. 220f68a32c8SEmil Constantinescu 221f68a32c8SEmil Constantinescu Level: developer 222f68a32c8SEmil Constantinescu 223f68a32c8SEmil Constantinescu .keywords: TS, TSRK, initialize, package 224f68a32c8SEmil Constantinescu .seealso: PetscInitialize() 225f68a32c8SEmil Constantinescu @*/ 226f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void) 227f68a32c8SEmil Constantinescu { 228186e87acSLisandro Dalcin PetscErrorCode ierr; 229e4dd521cSBarry Smith 230e4dd521cSBarry Smith PetscFunctionBegin; 231f68a32c8SEmil Constantinescu if (TSRKPackageInitialized) PetscFunctionReturn(0); 232f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_TRUE; 233f68a32c8SEmil Constantinescu ierr = TSRKRegisterAll();CHKERRQ(ierr); 234f68a32c8SEmil Constantinescu ierr = PetscObjectComposedDataRegister(&explicit_stage_time_id);CHKERRQ(ierr); 235f68a32c8SEmil Constantinescu ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr); 236f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 237f68a32c8SEmil Constantinescu } 238186e87acSLisandro Dalcin 239f68a32c8SEmil Constantinescu #undef __FUNCT__ 240f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKFinalizePackage" 241f68a32c8SEmil Constantinescu /*@C 242f68a32c8SEmil Constantinescu TSRKFinalizePackage - This function destroys everything in the TSRK package. It is 243f68a32c8SEmil Constantinescu called from PetscFinalize(). 244186e87acSLisandro Dalcin 245f68a32c8SEmil Constantinescu Level: developer 246f68a32c8SEmil Constantinescu 247f68a32c8SEmil Constantinescu .keywords: Petsc, destroy, package 248f68a32c8SEmil Constantinescu .seealso: PetscFinalize() 249f68a32c8SEmil Constantinescu @*/ 250f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void) 251f68a32c8SEmil Constantinescu { 252f68a32c8SEmil Constantinescu PetscErrorCode ierr; 253f68a32c8SEmil Constantinescu 254f68a32c8SEmil Constantinescu PetscFunctionBegin; 255f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_FALSE; 256f68a32c8SEmil Constantinescu ierr = TSRKRegisterDestroy();CHKERRQ(ierr); 257f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 258f68a32c8SEmil Constantinescu } 259f68a32c8SEmil Constantinescu 260f68a32c8SEmil Constantinescu #undef __FUNCT__ 261f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKRegister" 262f68a32c8SEmil Constantinescu /*@C 263f68a32c8SEmil Constantinescu TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 264f68a32c8SEmil Constantinescu 265f68a32c8SEmil Constantinescu Not Collective, but the same schemes should be registered on all processes on which they will be used 266f68a32c8SEmil Constantinescu 267f68a32c8SEmil Constantinescu Input Parameters: 268f68a32c8SEmil Constantinescu + name - identifier for method 269f68a32c8SEmil Constantinescu . order - approximation order of method 270f68a32c8SEmil Constantinescu . s - number of stages, this is the dimension of the matrices below 271f68a32c8SEmil Constantinescu . A - stage coefficients (dimension s*s, row-major) 272f68a32c8SEmil Constantinescu . b - step completion table (dimension s; NULL to use last row of A) 273f68a32c8SEmil Constantinescu . c - abscissa (dimension s; NULL to use row sums of A) 274f68a32c8SEmil Constantinescu . bembed - completion table for embedded method (dimension s; NULL if not available) 275f68a32c8SEmil Constantinescu . pinterp - Order of the interpolation scheme, equal to the number of columns of binterp 276f68a32c8SEmil Constantinescu - binterp - Coefficients of the interpolation formula (dimension s*pinterp; NULL to reuse binterpt) 277f68a32c8SEmil Constantinescu 278f68a32c8SEmil Constantinescu Notes: 279f68a32c8SEmil Constantinescu Several RK methods are provided, this function is only needed to create new methods. 280f68a32c8SEmil Constantinescu 281f68a32c8SEmil Constantinescu Level: advanced 282f68a32c8SEmil Constantinescu 283f68a32c8SEmil Constantinescu .keywords: TS, register 284f68a32c8SEmil Constantinescu 285f68a32c8SEmil Constantinescu .seealso: TSRK 286f68a32c8SEmil Constantinescu @*/ 287f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s, 288f68a32c8SEmil Constantinescu const PetscReal A[],const PetscReal b[],const PetscReal c[], 289f68a32c8SEmil Constantinescu const PetscReal bembed[], 290f68a32c8SEmil Constantinescu PetscInt pinterp,const PetscReal binterp[]) 291f68a32c8SEmil Constantinescu { 292f68a32c8SEmil Constantinescu PetscErrorCode ierr; 293f68a32c8SEmil Constantinescu RKTableauLink link; 294f68a32c8SEmil Constantinescu RKTableau t; 295f68a32c8SEmil Constantinescu PetscInt i,j; 296f68a32c8SEmil Constantinescu 297f68a32c8SEmil Constantinescu PetscFunctionBegin; 298f68a32c8SEmil Constantinescu ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr); 299f68a32c8SEmil Constantinescu ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr); 300f68a32c8SEmil Constantinescu t = &link->tab; 301f68a32c8SEmil Constantinescu ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 302f68a32c8SEmil Constantinescu t->order = order; 303f68a32c8SEmil Constantinescu t->s = s; 304f68a32c8SEmil Constantinescu ierr = PetscMalloc3(s*s,PetscReal,&t->A,s,PetscReal,&t->b,s,PetscReal,&t->c);CHKERRQ(ierr); 305f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 306f68a32c8SEmil Constantinescu if (b) { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); } 307f68a32c8SEmil Constantinescu else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i]; 308f68a32c8SEmil Constantinescu if (c) { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); } 309f68a32c8SEmil 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]; 310f68a32c8SEmil Constantinescu 311f68a32c8SEmil Constantinescu if (bembed) { 312f68a32c8SEmil Constantinescu ierr = PetscMalloc(s*sizeof(PetscReal),&t->bembed);CHKERRQ(ierr); 313f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr); 314f68a32c8SEmil Constantinescu } 315f68a32c8SEmil Constantinescu 316f68a32c8SEmil Constantinescu t->pinterp = pinterp; 317f68a32c8SEmil Constantinescu ierr = PetscMalloc(s*pinterp*sizeof(PetscReal),&t->binterp);CHKERRQ(ierr); 318f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->binterp,binterp,s*pinterp*sizeof(binterp[0]));CHKERRQ(ierr); 319f68a32c8SEmil Constantinescu link->next = RKTableauList; 320f68a32c8SEmil Constantinescu RKTableauList = link; 321f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 322f68a32c8SEmil Constantinescu } 323f68a32c8SEmil Constantinescu 324f68a32c8SEmil Constantinescu #undef __FUNCT__ 325f68a32c8SEmil Constantinescu #define __FUNCT__ "TSEvaluateStep_RK" 326e4dd521cSBarry Smith /* 327f68a32c8SEmil Constantinescu The step completion formula is 328e4dd521cSBarry Smith 329f68a32c8SEmil Constantinescu x1 = x0 + h b^T YdotRHS 330f68a32c8SEmil Constantinescu 331f68a32c8SEmil Constantinescu This function can be called before or after ts->vec_sol has been updated. 332f68a32c8SEmil Constantinescu Suppose we have a completion formula (b) and an embedded formula (be) of different order. 333f68a32c8SEmil Constantinescu We can write 334f68a32c8SEmil Constantinescu 335f68a32c8SEmil Constantinescu x1e = x0 + h be^T YdotRHS 336f68a32c8SEmil Constantinescu = x1 - h b^T YdotRHS + h be^T YdotRHS 337f68a32c8SEmil Constantinescu = x1 + h (be - b)^T YdotRHS 338f68a32c8SEmil Constantinescu 339f68a32c8SEmil Constantinescu so we can evaluate the method with different order even after the step has been optimistically completed. 340e4dd521cSBarry Smith */ 341f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done) 342f68a32c8SEmil Constantinescu { 343f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 344f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 345f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 346f68a32c8SEmil Constantinescu PetscReal h; 347f68a32c8SEmil Constantinescu PetscInt s = tab->s,j; 348f68a32c8SEmil Constantinescu PetscErrorCode ierr; 349f68a32c8SEmil Constantinescu 350f68a32c8SEmil Constantinescu PetscFunctionBegin; 351f68a32c8SEmil Constantinescu switch (rk->status) { 352f68a32c8SEmil Constantinescu case TS_STEP_INCOMPLETE: 353f68a32c8SEmil Constantinescu case TS_STEP_PENDING: 354f68a32c8SEmil Constantinescu h = ts->time_step; break; 355f68a32c8SEmil Constantinescu case TS_STEP_COMPLETE: 356f68a32c8SEmil Constantinescu h = ts->time_step_prev; break; 357f68a32c8SEmil Constantinescu default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 358f68a32c8SEmil Constantinescu } 359f68a32c8SEmil Constantinescu if (order == tab->order) { 360f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { 361f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 362f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->b[j]; 363f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 364f68a32c8SEmil Constantinescu } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 365f68a32c8SEmil Constantinescu if (done) *done = PETSC_TRUE; 366f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 367f68a32c8SEmil Constantinescu } else if (order == tab->order-1) { 368f68a32c8SEmil Constantinescu if (!tab->bembed) goto unavailable; 369f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (be) */ 370f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 371f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 372f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 373f68a32c8SEmil Constantinescu } else { /* Rollback and re-complete using (be-b) */ 374f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 375f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 376f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 377f68a32c8SEmil Constantinescu } 378f68a32c8SEmil Constantinescu if (done) *done = PETSC_TRUE; 379f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 380f68a32c8SEmil Constantinescu } 381f68a32c8SEmil Constantinescu unavailable: 382f68a32c8SEmil Constantinescu if (done) *done = PETSC_FALSE; 383f68a32c8SEmil 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); 384f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 385f68a32c8SEmil Constantinescu } 386f68a32c8SEmil Constantinescu 387f68a32c8SEmil Constantinescu #undef __FUNCT__ 388f68a32c8SEmil Constantinescu #define __FUNCT__ "TSStep_RK" 389f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts) 390f68a32c8SEmil Constantinescu { 391f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 392f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 393f68a32c8SEmil Constantinescu const PetscInt s = tab->s; 394f68a32c8SEmil Constantinescu const PetscReal *A = tab->A,*b = tab->b,*c = tab->c; 395f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 396f68a32c8SEmil Constantinescu Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 397f68a32c8SEmil Constantinescu TSAdapt adapt; 398f68a32c8SEmil Constantinescu PetscInt i,j,reject,next_scheme; 399f68a32c8SEmil Constantinescu PetscReal next_time_step; 400f68a32c8SEmil Constantinescu PetscReal t; 401f68a32c8SEmil Constantinescu PetscBool accept; 402f68a32c8SEmil Constantinescu PetscErrorCode ierr; 403f68a32c8SEmil Constantinescu 404f68a32c8SEmil Constantinescu PetscFunctionBegin; 405f68a32c8SEmil Constantinescu 406f68a32c8SEmil Constantinescu next_time_step = ts->time_step; 407f68a32c8SEmil Constantinescu t = ts->ptime; 408f68a32c8SEmil Constantinescu accept = PETSC_TRUE; 409f68a32c8SEmil Constantinescu rk->status = TS_STEP_INCOMPLETE; 410f68a32c8SEmil Constantinescu 411f68a32c8SEmil Constantinescu 412f68a32c8SEmil Constantinescu for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) { 413f68a32c8SEmil Constantinescu PetscReal h = ts->time_step; 414f68a32c8SEmil Constantinescu ierr = TSPreStep(ts);CHKERRQ(ierr); 415f68a32c8SEmil Constantinescu for (i=0; i<s; i++) { 416f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 417f68a32c8SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 418f68a32c8SEmil Constantinescu ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 419f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 420f68a32c8SEmil Constantinescu ierr = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr); 421f68a32c8SEmil Constantinescu if (!accept) goto reject_step; 422f68a32c8SEmil Constantinescu ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 423f68a32c8SEmil Constantinescu } 424f68a32c8SEmil Constantinescu ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 425f68a32c8SEmil Constantinescu rk->status = TS_STEP_PENDING; 426f68a32c8SEmil Constantinescu 427f68a32c8SEmil Constantinescu /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */ 428f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 429f68a32c8SEmil Constantinescu ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 430f68a32c8SEmil Constantinescu ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr); 431f68a32c8SEmil Constantinescu ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr); 432f68a32c8SEmil Constantinescu if (accept) { 433f68a32c8SEmil Constantinescu /* ignore next_scheme for now */ 434f68a32c8SEmil Constantinescu ts->ptime += ts->time_step; 435f68a32c8SEmil Constantinescu ts->time_step = next_time_step; 436e3caeda6SBarry Smith ts->steps++; 437f68a32c8SEmil Constantinescu rk->status = TS_STEP_COMPLETE; 438f68a32c8SEmil Constantinescu ierr = PetscObjectComposedDataSetReal((PetscObject)ts->vec_sol,explicit_stage_time_id,ts->ptime);CHKERRQ(ierr); 439e4dd521cSBarry Smith 440f68a32c8SEmil Constantinescu break; 441f68a32c8SEmil Constantinescu } else { /* Roll back the current step */ 442f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = -h*b[j]; 443f68a32c8SEmil Constantinescu ierr = VecMAXPY(ts->vec_sol,s,w,rk->YdotRHS);CHKERRQ(ierr); 444f68a32c8SEmil Constantinescu ts->time_step = next_time_step; 445f68a32c8SEmil Constantinescu rk->status = TS_STEP_INCOMPLETE; 446f68a32c8SEmil Constantinescu } 447f68a32c8SEmil Constantinescu reject_step: continue; 448f68a32c8SEmil Constantinescu } 449f68a32c8SEmil Constantinescu if (rk->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED; 450f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 451e4dd521cSBarry Smith } 452e4dd521cSBarry Smith 453f68a32c8SEmil Constantinescu #undef __FUNCT__ 454f68a32c8SEmil Constantinescu #define __FUNCT__ "TSInterpolate_RK" 455f68a32c8SEmil Constantinescu static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X) 456f68a32c8SEmil Constantinescu { 457f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 458f68a32c8SEmil Constantinescu PetscInt s = rk->tableau->s,pinterp = rk->tableau->pinterp,i,j; 459f68a32c8SEmil Constantinescu PetscReal h; 460f68a32c8SEmil Constantinescu PetscReal tt,t; 461f68a32c8SEmil Constantinescu PetscScalar *b; 462f68a32c8SEmil Constantinescu const PetscReal *B = rk->tableau->binterp; 463f68a32c8SEmil Constantinescu PetscErrorCode ierr; 464e4dd521cSBarry Smith 465f68a32c8SEmil Constantinescu PetscFunctionBegin; 466f68a32c8SEmil Constantinescu if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 467f68a32c8SEmil Constantinescu switch (rk->status) { 468f68a32c8SEmil Constantinescu case TS_STEP_INCOMPLETE: 469f68a32c8SEmil Constantinescu case TS_STEP_PENDING: 470f68a32c8SEmil Constantinescu h = ts->time_step; 471f68a32c8SEmil Constantinescu t = (itime - ts->ptime)/h; 472f68a32c8SEmil Constantinescu break; 473f68a32c8SEmil Constantinescu case TS_STEP_COMPLETE: 474f68a32c8SEmil Constantinescu h = ts->time_step_prev; 475f68a32c8SEmil Constantinescu t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 476f68a32c8SEmil Constantinescu break; 477f68a32c8SEmil Constantinescu default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 478e4dd521cSBarry Smith } 479f68a32c8SEmil Constantinescu ierr = PetscMalloc(s*sizeof(PetscScalar),&b);CHKERRQ(ierr); 480f68a32c8SEmil Constantinescu for (i=0; i<s; i++) b[i] = 0; 481f68a32c8SEmil Constantinescu for (j=0,tt=t; j<pinterp; j++,tt*=t) { 482f68a32c8SEmil Constantinescu for (i=0; i<s; i++) { 483f68a32c8SEmil Constantinescu b[i] += h * B[i*pinterp+j] * tt; 484e4dd521cSBarry Smith } 485f68a32c8SEmil Constantinescu } 486f68a32c8SEmil Constantinescu ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr); 487f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr); 488f68a32c8SEmil Constantinescu ierr = PetscFree(b);CHKERRQ(ierr); 489e4dd521cSBarry Smith PetscFunctionReturn(0); 490e4dd521cSBarry Smith } 491e4dd521cSBarry Smith 492e4dd521cSBarry Smith /*------------------------------------------------------------*/ 493e4dd521cSBarry Smith #undef __FUNCT__ 494277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_RK" 495277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts) 496e4dd521cSBarry Smith { 4975f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 498f68a32c8SEmil Constantinescu PetscInt s; 4996849ba73SBarry Smith PetscErrorCode ierr; 500e4dd521cSBarry Smith 501e4dd521cSBarry Smith PetscFunctionBegin; 502f68a32c8SEmil Constantinescu if (!rk->tableau) PetscFunctionReturn(0); 503f68a32c8SEmil Constantinescu s = rk->tableau->s; 504f68a32c8SEmil Constantinescu ierr = VecDestroyVecs(s,&rk->Y);CHKERRQ(ierr); 505f68a32c8SEmil Constantinescu ierr = VecDestroyVecs(s,&rk->YdotRHS);CHKERRQ(ierr); 506f68a32c8SEmil Constantinescu ierr = PetscFree(rk->work);CHKERRQ(ierr); 507277b19d0SLisandro Dalcin PetscFunctionReturn(0); 508e4dd521cSBarry Smith } 509277b19d0SLisandro Dalcin 510277b19d0SLisandro Dalcin #undef __FUNCT__ 511277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_RK" 512277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_RK(TS ts) 513277b19d0SLisandro Dalcin { 514277b19d0SLisandro Dalcin PetscErrorCode ierr; 515277b19d0SLisandro Dalcin 516277b19d0SLisandro Dalcin PetscFunctionBegin; 517277b19d0SLisandro Dalcin ierr = TSReset_RK(ts);CHKERRQ(ierr); 518277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 519f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr); 520f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr); 521f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 522f68a32c8SEmil Constantinescu } 523f68a32c8SEmil Constantinescu 524f68a32c8SEmil Constantinescu 525f68a32c8SEmil Constantinescu #undef __FUNCT__ 526f68a32c8SEmil Constantinescu #define __FUNCT__ "DMCoarsenHook_TSRK" 527f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx) 528f68a32c8SEmil Constantinescu { 529f68a32c8SEmil Constantinescu PetscFunctionBegin; 530f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 531f68a32c8SEmil Constantinescu } 532f68a32c8SEmil Constantinescu 533f68a32c8SEmil Constantinescu #undef __FUNCT__ 534f68a32c8SEmil Constantinescu #define __FUNCT__ "DMRestrictHook_TSRK" 535f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 536f68a32c8SEmil Constantinescu { 537f68a32c8SEmil Constantinescu PetscFunctionBegin; 538f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 539f68a32c8SEmil Constantinescu } 540f68a32c8SEmil Constantinescu 541f68a32c8SEmil Constantinescu 542f68a32c8SEmil Constantinescu #undef __FUNCT__ 543f68a32c8SEmil Constantinescu #define __FUNCT__ "DMSubDomainHook_TSRK" 544f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx) 545f68a32c8SEmil Constantinescu { 546f68a32c8SEmil Constantinescu PetscFunctionBegin; 547f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 548f68a32c8SEmil Constantinescu } 549f68a32c8SEmil Constantinescu 550f68a32c8SEmil Constantinescu #undef __FUNCT__ 551f68a32c8SEmil Constantinescu #define __FUNCT__ "DMSubDomainRestrictHook_TSRK" 552f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 553f68a32c8SEmil Constantinescu { 554f68a32c8SEmil Constantinescu 555f68a32c8SEmil Constantinescu PetscFunctionBegin; 556f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 557f68a32c8SEmil Constantinescu } 558f68a32c8SEmil Constantinescu 559f68a32c8SEmil Constantinescu #undef __FUNCT__ 560f68a32c8SEmil Constantinescu #define __FUNCT__ "TSSetUp_RK" 561f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts) 562f68a32c8SEmil Constantinescu { 563f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 564f68a32c8SEmil Constantinescu RKTableau tab; 565f68a32c8SEmil Constantinescu PetscInt s; 566f68a32c8SEmil Constantinescu PetscErrorCode ierr; 567f68a32c8SEmil Constantinescu DM dm; 568f68a32c8SEmil Constantinescu 569f68a32c8SEmil Constantinescu PetscFunctionBegin; 570f68a32c8SEmil Constantinescu if (!rk->tableau) { 571f68a32c8SEmil Constantinescu ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr); 572f68a32c8SEmil Constantinescu } 573f68a32c8SEmil Constantinescu tab = rk->tableau; 574f68a32c8SEmil Constantinescu s = tab->s; 575f68a32c8SEmil Constantinescu ierr = VecDuplicateVecs(ts->vec_sol,s,&rk->Y);CHKERRQ(ierr); 576f68a32c8SEmil Constantinescu ierr = VecDuplicateVecs(ts->vec_sol,s,&rk->YdotRHS);CHKERRQ(ierr); 577f68a32c8SEmil Constantinescu ierr = PetscMalloc(s*sizeof(rk->work[0]),&rk->work);CHKERRQ(ierr); 578f68a32c8SEmil Constantinescu ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 579f68a32c8SEmil Constantinescu if (dm) { 580f68a32c8SEmil Constantinescu ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 581f68a32c8SEmil Constantinescu ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 582f68a32c8SEmil Constantinescu } 583e4dd521cSBarry Smith PetscFunctionReturn(0); 584e4dd521cSBarry Smith } 585e4dd521cSBarry Smith /*------------------------------------------------------------*/ 586e4dd521cSBarry Smith 587e4dd521cSBarry Smith #undef __FUNCT__ 5885f70b478SJed Brown #define __FUNCT__ "TSSetFromOptions_RK" 5895f70b478SJed Brown static PetscErrorCode TSSetFromOptions_RK(TS ts) 590e4dd521cSBarry Smith { 591dfbe8321SBarry Smith PetscErrorCode ierr; 592f68a32c8SEmil Constantinescu char rktype[256]; 593262ff7bbSBarry Smith 594e4dd521cSBarry Smith PetscFunctionBegin; 595262ff7bbSBarry Smith ierr = PetscOptionsHead("RK ODE solver options");CHKERRQ(ierr); 596f68a32c8SEmil Constantinescu { 597f68a32c8SEmil Constantinescu RKTableauLink link; 598f68a32c8SEmil Constantinescu PetscInt count,choice; 599f68a32c8SEmil Constantinescu PetscBool flg; 600f68a32c8SEmil Constantinescu const char **namelist; 601f68a32c8SEmil Constantinescu ierr = PetscStrncpy(rktype,TSRKDefault,sizeof(rktype));CHKERRQ(ierr); 602f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) ; 603f68a32c8SEmil Constantinescu ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr); 604f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 605f68a32c8SEmil Constantinescu ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rktype,&choice,&flg);CHKERRQ(ierr); 606f68a32c8SEmil Constantinescu ierr = TSRKSetType(ts,flg ? namelist[choice] : rktype);CHKERRQ(ierr); 607f68a32c8SEmil Constantinescu ierr = PetscFree(namelist);CHKERRQ(ierr); 608f68a32c8SEmil Constantinescu } 609262ff7bbSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 610e4dd521cSBarry Smith PetscFunctionReturn(0); 611e4dd521cSBarry Smith } 612e4dd521cSBarry Smith 613e4dd521cSBarry Smith #undef __FUNCT__ 614f68a32c8SEmil Constantinescu #define __FUNCT__ "PetscFormatRealArray" 615f68a32c8SEmil Constantinescu static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[]) 616f68a32c8SEmil Constantinescu { 617f68a32c8SEmil Constantinescu PetscErrorCode ierr; 618f68a32c8SEmil Constantinescu PetscInt i; 619f68a32c8SEmil Constantinescu size_t left,count; 620f68a32c8SEmil Constantinescu char *p; 621f68a32c8SEmil Constantinescu 622f68a32c8SEmil Constantinescu PetscFunctionBegin; 623f68a32c8SEmil Constantinescu for (i=0,p=buf,left=len; i<n; i++) { 624f68a32c8SEmil Constantinescu ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr); 625f68a32c8SEmil Constantinescu if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer"); 626f68a32c8SEmil Constantinescu left -= count; 627f68a32c8SEmil Constantinescu p += count; 628f68a32c8SEmil Constantinescu *p++ = ' '; 629f68a32c8SEmil Constantinescu } 630f68a32c8SEmil Constantinescu p[i ? 0 : -1] = 0; 631f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 632f68a32c8SEmil Constantinescu } 633f68a32c8SEmil Constantinescu 634f68a32c8SEmil Constantinescu #undef __FUNCT__ 6355f70b478SJed Brown #define __FUNCT__ "TSView_RK" 6365f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer) 637e4dd521cSBarry Smith { 6385f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 639f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 6408370ee3bSLisandro Dalcin PetscBool iascii; 641dfbe8321SBarry Smith PetscErrorCode ierr; 642f68a32c8SEmil Constantinescu TSAdapt adapt; 6432cdf8120Sasbjorn 644e4dd521cSBarry Smith PetscFunctionBegin; 645251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 6468370ee3bSLisandro Dalcin if (iascii) { 647f68a32c8SEmil Constantinescu TSRKType rktype; 648f68a32c8SEmil Constantinescu char buf[512]; 649f68a32c8SEmil Constantinescu ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr); 650f68a32c8SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," RK %s\n",rktype);CHKERRQ(ierr); 651f68a32c8SEmil Constantinescu ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 652f68a32c8SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," Abscissa c = %s\n",buf);CHKERRQ(ierr); 6538370ee3bSLisandro Dalcin } 654f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 655f68a32c8SEmil Constantinescu ierr = TSAdaptView(adapt,viewer);CHKERRQ(ierr); 656f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 657f68a32c8SEmil Constantinescu } 658f68a32c8SEmil Constantinescu 659f68a32c8SEmil Constantinescu #undef __FUNCT__ 660f68a32c8SEmil Constantinescu #define __FUNCT__ "TSLoad_RK" 661f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer) 662f68a32c8SEmil Constantinescu { 663f68a32c8SEmil Constantinescu PetscErrorCode ierr; 664f68a32c8SEmil Constantinescu TSAdapt tsadapt; 665f68a32c8SEmil Constantinescu 666f68a32c8SEmil Constantinescu PetscFunctionBegin; 667f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&tsadapt);CHKERRQ(ierr); 668f68a32c8SEmil Constantinescu ierr = TSAdaptLoad(tsadapt,viewer);CHKERRQ(ierr); 669f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 670f68a32c8SEmil Constantinescu } 671f68a32c8SEmil Constantinescu 672f68a32c8SEmil Constantinescu #undef __FUNCT__ 673f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKSetType" 674f68a32c8SEmil Constantinescu /*@C 675f68a32c8SEmil Constantinescu TSRKSetType - Set the type of RK scheme 676f68a32c8SEmil Constantinescu 677f68a32c8SEmil Constantinescu Logically collective 678f68a32c8SEmil Constantinescu 679f68a32c8SEmil Constantinescu Input Parameter: 680f68a32c8SEmil Constantinescu + ts - timestepping context 681f68a32c8SEmil Constantinescu - rktype - type of RK-scheme 682f68a32c8SEmil Constantinescu 683f68a32c8SEmil Constantinescu Level: intermediate 684f68a32c8SEmil Constantinescu 685f68a32c8SEmil Constantinescu .seealso: TSRKGetType(), TSRK, TSRK2, TSRK3, TSRKPRSSP2, TSRK4, TSRK5 686f68a32c8SEmil Constantinescu @*/ 687f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype) 688f68a32c8SEmil Constantinescu { 689f68a32c8SEmil Constantinescu PetscErrorCode ierr; 690f68a32c8SEmil Constantinescu 691f68a32c8SEmil Constantinescu PetscFunctionBegin; 692f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 693f68a32c8SEmil Constantinescu ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr); 694f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 695f68a32c8SEmil Constantinescu } 696f68a32c8SEmil Constantinescu 697f68a32c8SEmil Constantinescu #undef __FUNCT__ 698f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKGetType" 699f68a32c8SEmil Constantinescu /*@C 700f68a32c8SEmil Constantinescu TSRKGetType - Get the type of RK scheme 701f68a32c8SEmil Constantinescu 702f68a32c8SEmil Constantinescu Logically collective 703f68a32c8SEmil Constantinescu 704f68a32c8SEmil Constantinescu Input Parameter: 705f68a32c8SEmil Constantinescu . ts - timestepping context 706f68a32c8SEmil Constantinescu 707f68a32c8SEmil Constantinescu Output Parameter: 708f68a32c8SEmil Constantinescu . rktype - type of RK-scheme 709f68a32c8SEmil Constantinescu 710f68a32c8SEmil Constantinescu Level: intermediate 711f68a32c8SEmil Constantinescu 712f68a32c8SEmil Constantinescu .seealso: TSRKGetType() 713f68a32c8SEmil Constantinescu @*/ 714f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype) 715f68a32c8SEmil Constantinescu { 716f68a32c8SEmil Constantinescu PetscErrorCode ierr; 717f68a32c8SEmil Constantinescu 718f68a32c8SEmil Constantinescu PetscFunctionBegin; 719f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 720f68a32c8SEmil Constantinescu ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr); 721f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 722f68a32c8SEmil Constantinescu } 723f68a32c8SEmil Constantinescu 724f68a32c8SEmil Constantinescu #undef __FUNCT__ 725f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKGetType_RK" 726f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype) 727f68a32c8SEmil Constantinescu { 728f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 729f68a32c8SEmil Constantinescu PetscErrorCode ierr; 730f68a32c8SEmil Constantinescu 731f68a32c8SEmil Constantinescu PetscFunctionBegin; 732f68a32c8SEmil Constantinescu if (!rk->tableau) { 733f68a32c8SEmil Constantinescu ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr); 734f68a32c8SEmil Constantinescu } 735f68a32c8SEmil Constantinescu *rktype = rk->tableau->name; 736f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 737f68a32c8SEmil Constantinescu } 738f68a32c8SEmil Constantinescu #undef __FUNCT__ 739f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKSetType_RK" 740f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype) 741f68a32c8SEmil Constantinescu { 742f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 743f68a32c8SEmil Constantinescu PetscErrorCode ierr; 744f68a32c8SEmil Constantinescu PetscBool match; 745f68a32c8SEmil Constantinescu RKTableauLink link; 746f68a32c8SEmil Constantinescu 747f68a32c8SEmil Constantinescu PetscFunctionBegin; 748f68a32c8SEmil Constantinescu if (rk->tableau) { 749f68a32c8SEmil Constantinescu ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr); 750f68a32c8SEmil Constantinescu if (match) PetscFunctionReturn(0); 751f68a32c8SEmil Constantinescu } 752f68a32c8SEmil Constantinescu for (link = RKTableauList; link; link=link->next) { 753f68a32c8SEmil Constantinescu ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr); 754f68a32c8SEmil Constantinescu if (match) { 755f68a32c8SEmil Constantinescu ierr = TSReset_RK(ts);CHKERRQ(ierr); 756f68a32c8SEmil Constantinescu rk->tableau = &link->tab; 757f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 758f68a32c8SEmil Constantinescu } 759f68a32c8SEmil Constantinescu } 760f68a32c8SEmil Constantinescu SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype); 761e4dd521cSBarry Smith PetscFunctionReturn(0); 762e4dd521cSBarry Smith } 763e4dd521cSBarry Smith 764e4dd521cSBarry Smith /* ------------------------------------------------------------ */ 765ebe3b25bSBarry Smith /*MC 766f68a32c8SEmil Constantinescu TSRK - ODE and DAE solver using Runge-Kutta schemes 767ebe3b25bSBarry Smith 768f68a32c8SEmil Constantinescu The user should provide the right hand side of the equation 769f68a32c8SEmil Constantinescu using TSSetRHSFunction(). 770ebe3b25bSBarry Smith 771f68a32c8SEmil Constantinescu Notes: 772f68a32c8SEmil Constantinescu The default is TSRK3, it can be changed with TSRKSetType() or -ts_rk_type 773ebe3b25bSBarry Smith 774d41222bbSBarry Smith Level: beginner 775d41222bbSBarry Smith 776f68a32c8SEmil Constantinescu .seealso: TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3, 777f68a32c8SEmil Constantinescu TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister() 778ebe3b25bSBarry Smith 779ebe3b25bSBarry Smith M*/ 780e4dd521cSBarry Smith #undef __FUNCT__ 7815f70b478SJed Brown #define __FUNCT__ "TSCreate_RK" 7828cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts) 783e4dd521cSBarry Smith { 784f68a32c8SEmil Constantinescu TS_RK *th; 785dfbe8321SBarry Smith PetscErrorCode ierr; 786e4dd521cSBarry Smith 787e4dd521cSBarry Smith PetscFunctionBegin; 788f68a32c8SEmil Constantinescu #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 789f68a32c8SEmil Constantinescu ierr = TSRKInitializePackage();CHKERRQ(ierr); 790f68a32c8SEmil Constantinescu #endif 791f68a32c8SEmil Constantinescu 792f68a32c8SEmil Constantinescu ts->ops->reset = TSReset_RK; 7935f70b478SJed Brown ts->ops->destroy = TSDestroy_RK; 7945f70b478SJed Brown ts->ops->view = TSView_RK; 795f68a32c8SEmil Constantinescu ts->ops->load = TSLoad_RK; 796f68a32c8SEmil Constantinescu ts->ops->setup = TSSetUp_RK; 797f68a32c8SEmil Constantinescu ts->ops->step = TSStep_RK; 798f68a32c8SEmil Constantinescu ts->ops->interpolate = TSInterpolate_RK; 799f68a32c8SEmil Constantinescu ts->ops->evaluatestep = TSEvaluateStep_RK; 800f68a32c8SEmil Constantinescu ts->ops->setfromoptions = TSSetFromOptions_RK; 801e4dd521cSBarry Smith 802f68a32c8SEmil Constantinescu ierr = PetscNewLog(ts,TS_RK,&th);CHKERRQ(ierr); 803f68a32c8SEmil Constantinescu ts->data = (void*)th; 804e4dd521cSBarry Smith 805f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr); 806f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr); 807e4dd521cSBarry Smith PetscFunctionReturn(0); 808e4dd521cSBarry Smith } 809