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