1c4762a1bSJed Brown static const char help[] = "1D periodic Finite Volume solver by a particular slope limiter with semidiscrete time stepping.\n" 2c4762a1bSJed Brown " advection - Constant coefficient scalar advection\n" 3c4762a1bSJed Brown " u_t + (a*u)_x = 0\n" 4c4762a1bSJed Brown " for this toy problem, we choose different meshsizes for different sub-domains, say\n" 5c4762a1bSJed Brown " hxs = (xmax - xmin)/2.0*(hratio+1.0)/Mx, \n" 6c4762a1bSJed Brown " hxf = (xmax - xmin)/2.0*(1.0+1.0/hratio)/Mx, \n" 76aad120cSJose E. Roman " with x belongs to (xmin,xmax), the number of total mesh points is Mx and the ratio between the meshsize of coarse\n\n" 8c4762a1bSJed Brown " grids and fine grids is hratio.\n" 9c4762a1bSJed Brown " exact - Exact Riemann solver which usually needs to perform a Newton iteration to connect\n" 10c4762a1bSJed Brown " the states across shocks and rarefactions\n" 11c4762a1bSJed Brown " simulation - use reference solution which is generated by smaller time step size to be true solution,\n" 12c4762a1bSJed Brown " also the reference solution should be generated by user and stored in a binary file.\n" 13c4762a1bSJed Brown " characteristic - Limit the characteristic variables, this is usually preferred (default)\n" 14c4762a1bSJed Brown "Several initial conditions can be chosen with -initial N\n\n" 15c4762a1bSJed Brown "The problem size should be set with -da_grid_x M\n\n" 16c4762a1bSJed Brown "This script choose the slope limiter by biased second-order upwind procedure which is proposed by Van Leer in 1994\n" 17c4762a1bSJed Brown " u(x_(k+1/2),t) = u(x_k,t) + phi(x_(k+1/2),t)*(u(x_k,t)-u(x_(k-1),t)) \n" 18c4762a1bSJed Brown " limiter phi(x_(k+1/2),t) = max(0,min(r(k+1/2),min(2,gamma(k+1/2)*r(k+1/2)+alpha(k+1/2)))) \n" 19c4762a1bSJed Brown " r(k+1/2) = (u(x_(k+1))-u(x_k))/(u(x_k)-u(x_(k-1))) \n" 20c4762a1bSJed Brown " alpha(k+1/2) = (h_k*h_(k+1))/(h_(k-1)+h_k)/(h_(k-1)+h_k+h_(k+1)) \n" 21c4762a1bSJed Brown " gamma(k+1/2) = h_k*(h_(k-1)+h_k)/(h_k+h_(k+1))/(h_(k-1)+h_k+h_(k+1)) \n"; 22c4762a1bSJed Brown 23c4762a1bSJed Brown #include <petscts.h> 24c4762a1bSJed Brown #include <petscdm.h> 25c4762a1bSJed Brown #include <petscdmda.h> 26c4762a1bSJed Brown #include <petscdraw.h> 27c4762a1bSJed Brown #include <petscmath.h> 28c4762a1bSJed Brown 29d71ae5a4SJacob Faibussowitsch static inline PetscReal RangeMod(PetscReal a, PetscReal xmin, PetscReal xmax) 30d71ae5a4SJacob Faibussowitsch { 319371c9d4SSatish Balay PetscReal range = xmax - xmin; 329371c9d4SSatish Balay return xmin + PetscFmodReal(range + PetscFmodReal(a, range), range); 339371c9d4SSatish Balay } 34c4762a1bSJed Brown 35c4762a1bSJed Brown /* --------------------------------- Finite Volume data structures ----------------------------------- */ 36c4762a1bSJed Brown 379371c9d4SSatish Balay typedef enum { 389371c9d4SSatish Balay FVBC_PERIODIC, 399371c9d4SSatish Balay FVBC_OUTFLOW 409371c9d4SSatish Balay } FVBCType; 41c4762a1bSJed Brown static const char *FVBCTypes[] = {"PERIODIC", "OUTFLOW", "FVBCType", "FVBC_", 0}; 42c4762a1bSJed Brown 43c4762a1bSJed Brown typedef struct { 44c4762a1bSJed Brown PetscErrorCode (*sample)(void *, PetscInt, FVBCType, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal *); 45c4762a1bSJed Brown PetscErrorCode (*flux)(void *, const PetscScalar *, PetscScalar *, PetscReal *); 46c4762a1bSJed Brown PetscErrorCode (*destroy)(void *); 47c4762a1bSJed Brown void *user; 48c4762a1bSJed Brown PetscInt dof; 49c4762a1bSJed Brown char *fieldname[16]; 50c4762a1bSJed Brown } PhysicsCtx; 51c4762a1bSJed Brown 52c4762a1bSJed Brown typedef struct { 53c4762a1bSJed Brown PhysicsCtx physics; 54c4762a1bSJed Brown MPI_Comm comm; 55c4762a1bSJed Brown char prefix[256]; 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* Local work arrays */ 58c4762a1bSJed Brown PetscScalar *flux; /* Flux across interface */ 59c4762a1bSJed Brown PetscReal *speeds; /* Speeds of each wave */ 60c4762a1bSJed Brown PetscScalar *u; /* value at face */ 61c4762a1bSJed Brown 62c4762a1bSJed Brown PetscReal cfl_idt; /* Max allowable value of 1/Delta t */ 63c4762a1bSJed Brown PetscReal cfl; 64c4762a1bSJed Brown PetscReal xmin, xmax; 65c4762a1bSJed Brown PetscInt initial; 66c4762a1bSJed Brown PetscBool exact; 67c4762a1bSJed Brown PetscBool simulation; 68c4762a1bSJed Brown FVBCType bctype; 69c4762a1bSJed Brown PetscInt hratio; /* hratio = hslow/hfast */ 70c4762a1bSJed Brown IS isf, iss; 71c4762a1bSJed Brown PetscInt sf, fs; /* slow-fast and fast-slow interfaces */ 72c4762a1bSJed Brown } FVCtx; 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* --------------------------------- Physics ----------------------------------- */ 75d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsDestroy_SimpleFree(void *vctx) 76d71ae5a4SJacob Faibussowitsch { 77c4762a1bSJed Brown PetscFunctionBeginUser; 789566063dSJacob Faibussowitsch PetscCall(PetscFree(vctx)); 793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 80c4762a1bSJed Brown } 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* --------------------------------- Advection ----------------------------------- */ 83c4762a1bSJed Brown typedef struct { 84c4762a1bSJed Brown PetscReal a; /* advective velocity */ 85c4762a1bSJed Brown } AdvectCtx; 86c4762a1bSJed Brown 87d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsFlux_Advect(void *vctx, const PetscScalar *u, PetscScalar *flux, PetscReal *maxspeed) 88d71ae5a4SJacob Faibussowitsch { 89c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx *)vctx; 90c4762a1bSJed Brown PetscReal speed; 91c4762a1bSJed Brown 92c4762a1bSJed Brown PetscFunctionBeginUser; 93c4762a1bSJed Brown speed = ctx->a; 94c4762a1bSJed Brown flux[0] = speed * u[0]; 95c4762a1bSJed Brown *maxspeed = speed; 963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 97c4762a1bSJed Brown } 98c4762a1bSJed Brown 99d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsSample_Advect(void *vctx, PetscInt initial, FVBCType bctype, PetscReal xmin, PetscReal xmax, PetscReal t, PetscReal x, PetscReal *u) 100d71ae5a4SJacob Faibussowitsch { 101c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx *)vctx; 102c4762a1bSJed Brown PetscReal a = ctx->a, x0; 103c4762a1bSJed Brown 104c4762a1bSJed Brown PetscFunctionBeginUser; 105c4762a1bSJed Brown switch (bctype) { 106d71ae5a4SJacob Faibussowitsch case FVBC_OUTFLOW: 107d71ae5a4SJacob Faibussowitsch x0 = x - a * t; 108d71ae5a4SJacob Faibussowitsch break; 109d71ae5a4SJacob Faibussowitsch case FVBC_PERIODIC: 110d71ae5a4SJacob Faibussowitsch x0 = RangeMod(x - a * t, xmin, xmax); 111d71ae5a4SJacob Faibussowitsch break; 112d71ae5a4SJacob Faibussowitsch default: 113d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown BCType"); 114c4762a1bSJed Brown } 115c4762a1bSJed Brown switch (initial) { 116d71ae5a4SJacob Faibussowitsch case 0: 117d71ae5a4SJacob Faibussowitsch u[0] = (x0 < 0) ? 1 : -1; 118d71ae5a4SJacob Faibussowitsch break; 119d71ae5a4SJacob Faibussowitsch case 1: 120d71ae5a4SJacob Faibussowitsch u[0] = (x0 < 0) ? -1 : 1; 121d71ae5a4SJacob Faibussowitsch break; 122d71ae5a4SJacob Faibussowitsch case 2: 123d71ae5a4SJacob Faibussowitsch u[0] = (0 < x0 && x0 < 1) ? 1 : 0; 124d71ae5a4SJacob Faibussowitsch break; 125d71ae5a4SJacob Faibussowitsch case 3: 126d71ae5a4SJacob Faibussowitsch u[0] = PetscSinReal(2 * PETSC_PI * x0); 127d71ae5a4SJacob Faibussowitsch break; 128d71ae5a4SJacob Faibussowitsch case 4: 129d71ae5a4SJacob Faibussowitsch u[0] = PetscAbs(x0); 130d71ae5a4SJacob Faibussowitsch break; 131d71ae5a4SJacob Faibussowitsch case 5: 132d71ae5a4SJacob Faibussowitsch u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2 * PETSC_PI * x0)); 133d71ae5a4SJacob Faibussowitsch break; 134d71ae5a4SJacob Faibussowitsch case 6: 135d71ae5a4SJacob Faibussowitsch u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2 - x0 : 0)); 136d71ae5a4SJacob Faibussowitsch break; 137d71ae5a4SJacob Faibussowitsch case 7: 138d71ae5a4SJacob Faibussowitsch u[0] = PetscPowReal(PetscSinReal(PETSC_PI * x0), 10.0); 139d71ae5a4SJacob Faibussowitsch break; 140d71ae5a4SJacob Faibussowitsch default: 141d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition"); 142c4762a1bSJed Brown } 1433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 144c4762a1bSJed Brown } 145c4762a1bSJed Brown 146d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx) 147d71ae5a4SJacob Faibussowitsch { 148c4762a1bSJed Brown AdvectCtx *user; 149c4762a1bSJed Brown 150c4762a1bSJed Brown PetscFunctionBeginUser; 1519566063dSJacob Faibussowitsch PetscCall(PetscNew(&user)); 152c4762a1bSJed Brown ctx->physics.sample = PhysicsSample_Advect; 153c4762a1bSJed Brown ctx->physics.flux = PhysicsFlux_Advect; 154c4762a1bSJed Brown ctx->physics.destroy = PhysicsDestroy_SimpleFree; 155c4762a1bSJed Brown ctx->physics.user = user; 156c4762a1bSJed Brown ctx->physics.dof = 1; 1579566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy("u", &ctx->physics.fieldname[0])); 158c4762a1bSJed Brown user->a = 1; 159d0609cedSBarry Smith PetscOptionsBegin(ctx->comm, ctx->prefix, "Options for advection", ""); 160d71ae5a4SJacob Faibussowitsch { 161d71ae5a4SJacob Faibussowitsch PetscCall(PetscOptionsReal("-physics_advect_a", "Speed", "", user->a, &user->a, NULL)); 162d71ae5a4SJacob Faibussowitsch } 163d0609cedSBarry Smith PetscOptionsEnd(); 1643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 165c4762a1bSJed Brown } 166c4762a1bSJed Brown 167c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver ----------------------------------- */ 168c4762a1bSJed Brown 169d71ae5a4SJacob Faibussowitsch static PetscErrorCode FVRHSFunction(TS ts, PetscReal time, Vec X, Vec F, void *vctx) 170d71ae5a4SJacob Faibussowitsch { 171c4762a1bSJed Brown FVCtx *ctx = (FVCtx *)vctx; 172c4762a1bSJed Brown PetscInt i, j, Mx, dof, xs, xm, sf = ctx->sf, fs = ctx->fs; 173c4762a1bSJed Brown PetscReal hf, hs, cfl_idt = 0; 174c4762a1bSJed Brown PetscScalar *x, *f, *r, *min, *alpha, *gamma; 175c4762a1bSJed Brown Vec Xloc; 176c4762a1bSJed Brown DM da; 177c4762a1bSJed Brown 178c4762a1bSJed Brown PetscFunctionBeginUser; 1799566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 1809566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc)); /* Xloc contains ghost points */ 1819566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); /* Mx is the number of center points */ 182c4762a1bSJed Brown hs = (ctx->xmax - ctx->xmin) / 2.0 * (ctx->hratio + 1.0) / Mx; 183c4762a1bSJed Brown hf = (ctx->xmax - ctx->xmin) / 2.0 * (1.0 + 1.0 / ctx->hratio) / Mx; 1849566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); /* X is solution vector which does not contain ghost points */ 1859566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); 186*dd8e379bSPierre Jolivet PetscCall(VecZeroEntries(F)); /* F is the right-hand side function corresponds to center points */ 1879566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x)); 1889566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f)); 1899566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 1909566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(dof, &r, dof, &min, dof, &alpha, dof, &gamma)); 191c4762a1bSJed Brown 192c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 193c4762a1bSJed Brown for (i = xs - 2; i < 0; i++) { 194c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; 195c4762a1bSJed Brown } 196c4762a1bSJed Brown for (i = Mx; i < xs + xm + 2; i++) { 197c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; 198c4762a1bSJed Brown } 199c4762a1bSJed Brown } 200c4762a1bSJed Brown 201c4762a1bSJed Brown for (i = xs; i < xs + xm + 1; i++) { 202c4762a1bSJed Brown PetscReal maxspeed; 203c4762a1bSJed Brown PetscScalar *u; 204c4762a1bSJed Brown if (i < sf || i > fs + 1) { 205c4762a1bSJed Brown u = &ctx->u[0]; 206c4762a1bSJed Brown alpha[0] = 1.0 / 6.0; 207c4762a1bSJed Brown gamma[0] = 1.0 / 3.0; 208c4762a1bSJed Brown for (j = 0; j < dof; j++) { 209c4762a1bSJed Brown r[j] = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 210c4762a1bSJed Brown min[j] = PetscMin(r[j], 2.0); 211c4762a1bSJed Brown u[j] = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 212c4762a1bSJed Brown } 2139566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed)); 214c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hs)); 215c4762a1bSJed Brown if (i > xs) { 216c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hs; 217c4762a1bSJed Brown } 218c4762a1bSJed Brown if (i < xs + xm) { 219c4762a1bSJed Brown for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hs; 220c4762a1bSJed Brown } 221c4762a1bSJed Brown } else if (i == sf) { 222c4762a1bSJed Brown u = &ctx->u[0]; 223c4762a1bSJed Brown alpha[0] = hs * hf / (hs + hs) / (hs + hs + hf); 224c4762a1bSJed Brown gamma[0] = hs * (hs + hs) / (hs + hf) / (hs + hs + hf); 225c4762a1bSJed Brown for (j = 0; j < dof; j++) { 226c4762a1bSJed Brown r[j] = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 227c4762a1bSJed Brown min[j] = PetscMin(r[j], 2.0); 228c4762a1bSJed Brown u[j] = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 229c4762a1bSJed Brown } 2309566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed)); 231c4762a1bSJed Brown if (i > xs) { 232c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hs; 233c4762a1bSJed Brown } 234c4762a1bSJed Brown if (i < xs + xm) { 235c4762a1bSJed Brown for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hf; 236c4762a1bSJed Brown } 237c4762a1bSJed Brown } else if (i == sf + 1) { 238c4762a1bSJed Brown u = &ctx->u[0]; 239c4762a1bSJed Brown alpha[0] = hf * hf / (hs + hf) / (hs + hf + hf); 240c4762a1bSJed Brown gamma[0] = hf * (hs + hf) / (hf + hf) / (hs + hf + hf); 241c4762a1bSJed Brown for (j = 0; j < dof; j++) { 242c4762a1bSJed Brown r[j] = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 243c4762a1bSJed Brown min[j] = PetscMin(r[j], 2.0); 244c4762a1bSJed Brown u[j] = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 245c4762a1bSJed Brown } 2469566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed)); 247c4762a1bSJed Brown if (i > xs) { 248c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hf; 249c4762a1bSJed Brown } 250c4762a1bSJed Brown if (i < xs + xm) { 251c4762a1bSJed Brown for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hf; 252c4762a1bSJed Brown } 253c4762a1bSJed Brown } else if (i > sf + 1 && i < fs) { 254c4762a1bSJed Brown u = &ctx->u[0]; 255c4762a1bSJed Brown alpha[0] = 1.0 / 6.0; 256c4762a1bSJed Brown gamma[0] = 1.0 / 3.0; 257c4762a1bSJed Brown for (j = 0; j < dof; j++) { 258c4762a1bSJed Brown r[j] = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 259c4762a1bSJed Brown min[j] = PetscMin(r[j], 2.0); 260c4762a1bSJed Brown u[j] = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 261c4762a1bSJed Brown } 2629566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed)); 263c4762a1bSJed Brown if (i > xs) { 264c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hf; 265c4762a1bSJed Brown } 266c4762a1bSJed Brown if (i < xs + xm) { 267c4762a1bSJed Brown for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hf; 268c4762a1bSJed Brown } 269c4762a1bSJed Brown } else if (i == fs) { 270c4762a1bSJed Brown u = &ctx->u[0]; 271c4762a1bSJed Brown alpha[0] = hf * hs / (hf + hf) / (hf + hf + hs); 272c4762a1bSJed Brown gamma[0] = hf * (hf + hf) / (hf + hs) / (hf + hf + hs); 273c4762a1bSJed Brown for (j = 0; j < dof; j++) { 274c4762a1bSJed Brown r[j] = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 275c4762a1bSJed Brown min[j] = PetscMin(r[j], 2.0); 276c4762a1bSJed Brown u[j] = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 277c4762a1bSJed Brown } 2789566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed)); 279c4762a1bSJed Brown if (i > xs) { 280c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hf; 281c4762a1bSJed Brown } 282c4762a1bSJed Brown if (i < xs + xm) { 283c4762a1bSJed Brown for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hs; 284c4762a1bSJed Brown } 285c4762a1bSJed Brown } else if (i == fs + 1) { 286c4762a1bSJed Brown u = &ctx->u[0]; 287c4762a1bSJed Brown alpha[0] = hs * hs / (hf + hs) / (hf + hs + hs); 288c4762a1bSJed Brown gamma[0] = hs * (hf + hs) / (hs + hs) / (hf + hs + hs); 289c4762a1bSJed Brown for (j = 0; j < dof; j++) { 290c4762a1bSJed Brown r[j] = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 291c4762a1bSJed Brown min[j] = PetscMin(r[j], 2.0); 292c4762a1bSJed Brown u[j] = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 293c4762a1bSJed Brown } 2949566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed)); 295c4762a1bSJed Brown if (i > xs) { 296c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hs; 297c4762a1bSJed Brown } 298c4762a1bSJed Brown if (i < xs + xm) { 299c4762a1bSJed Brown for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hs; 300c4762a1bSJed Brown } 301c4762a1bSJed Brown } 302c4762a1bSJed Brown } 3039566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); 3049566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f)); 3059566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc)); 306712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(&cfl_idt, &ctx->cfl_idt, 1, MPIU_SCALAR, MPIU_MAX, PetscObjectComm((PetscObject)da))); 307c4762a1bSJed Brown if (0) { 308c2a16db8SHong Zhang /* We need a way to inform the TS of a CFL constraint, this is a debugging fragment */ 309c4762a1bSJed Brown PetscReal dt, tnow; 3109566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt)); 3119566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &tnow)); 31248a46eb9SPierre Jolivet if (dt > 0.5 / ctx->cfl_idt) PetscCall(PetscPrintf(ctx->comm, "Stability constraint exceeded at t=%g, dt %g > %g\n", (double)tnow, (double)dt, (double)(0.5 / ctx->cfl_idt))); 313c4762a1bSJed Brown } 3149566063dSJacob Faibussowitsch PetscCall(PetscFree4(r, min, alpha, gamma)); 3153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 316c4762a1bSJed Brown } 317c4762a1bSJed Brown 318d71ae5a4SJacob Faibussowitsch static PetscErrorCode FVRHSFunctionslow(TS ts, PetscReal time, Vec X, Vec F, void *vctx) 319d71ae5a4SJacob Faibussowitsch { 320c4762a1bSJed Brown FVCtx *ctx = (FVCtx *)vctx; 321c4762a1bSJed Brown PetscInt i, j, Mx, dof, xs, xm, islow = 0, sf = ctx->sf, fs = ctx->fs; 322c4762a1bSJed Brown PetscReal hf, hs; 323c4762a1bSJed Brown PetscScalar *x, *f, *r, *min, *alpha, *gamma; 324c4762a1bSJed Brown Vec Xloc; 325c4762a1bSJed Brown DM da; 326c4762a1bSJed Brown 327c4762a1bSJed Brown PetscFunctionBeginUser; 3289566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 3299566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc)); /* Xloc contains ghost points */ 3309566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); /* Mx is the number of center points */ 331c4762a1bSJed Brown hs = (ctx->xmax - ctx->xmin) / 2.0 * (ctx->hratio + 1.0) / Mx; 332c4762a1bSJed Brown hf = (ctx->xmax - ctx->xmin) / 2.0 * (1.0 + 1.0 / ctx->hratio) / Mx; 3339566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); /* X is solution vector which does not contain ghost points */ 3349566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); 335*dd8e379bSPierre Jolivet PetscCall(VecZeroEntries(F)); /* F is the right-hand side function corresponds to center points */ 3369566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x)); 3379566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 3389566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 3399566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(dof, &r, dof, &min, dof, &alpha, dof, &gamma)); 340c4762a1bSJed Brown 341c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 342c4762a1bSJed Brown for (i = xs - 2; i < 0; i++) { 343c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; 344c4762a1bSJed Brown } 345c4762a1bSJed Brown for (i = Mx; i < xs + xm + 2; i++) { 346c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; 347c4762a1bSJed Brown } 348c4762a1bSJed Brown } 349c4762a1bSJed Brown 350c4762a1bSJed Brown for (i = xs; i < xs + xm + 1; i++) { 351c4762a1bSJed Brown PetscReal maxspeed; 352c4762a1bSJed Brown PetscScalar *u; 353c4762a1bSJed Brown if (i < sf) { 354c4762a1bSJed Brown u = &ctx->u[0]; 355c4762a1bSJed Brown alpha[0] = 1.0 / 6.0; 356c4762a1bSJed Brown gamma[0] = 1.0 / 3.0; 357c4762a1bSJed Brown for (j = 0; j < dof; j++) { 358c4762a1bSJed Brown r[j] = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 359c4762a1bSJed Brown min[j] = PetscMin(r[j], 2.0); 360c4762a1bSJed Brown u[j] = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 361c4762a1bSJed Brown } 3629566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed)); 363c4762a1bSJed Brown if (i > xs) { 364c2a16db8SHong Zhang for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hs; 365c4762a1bSJed Brown } 366c4762a1bSJed Brown if (i < xs + xm) { 367c2a16db8SHong Zhang for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hs; 368c4762a1bSJed Brown islow++; 369c4762a1bSJed Brown } 370c4762a1bSJed Brown } else if (i == sf) { 371c4762a1bSJed Brown u = &ctx->u[0]; 372c4762a1bSJed Brown alpha[0] = hs * hf / (hs + hs) / (hs + hs + hf); 373c4762a1bSJed Brown gamma[0] = hs * (hs + hs) / (hs + hf) / (hs + hs + hf); 374c4762a1bSJed Brown for (j = 0; j < dof; j++) { 375c4762a1bSJed Brown r[j] = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 376c4762a1bSJed Brown min[j] = PetscMin(r[j], 2.0); 377c4762a1bSJed Brown u[j] = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 378c4762a1bSJed Brown } 3799566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed)); 380c2a16db8SHong Zhang if (i > xs) { 381c2a16db8SHong Zhang for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hs; 382c4762a1bSJed Brown } 383c4762a1bSJed Brown } else if (i == fs) { 384c4762a1bSJed Brown u = &ctx->u[0]; 385c4762a1bSJed Brown alpha[0] = hf * hs / (hf + hf) / (hf + hf + hs); 386c4762a1bSJed Brown gamma[0] = hf * (hf + hf) / (hf + hs) / (hf + hf + hs); 387c4762a1bSJed Brown for (j = 0; j < dof; j++) { 388c4762a1bSJed Brown r[j] = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 389c4762a1bSJed Brown min[j] = PetscMin(r[j], 2.0); 390c4762a1bSJed Brown u[j] = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 391c4762a1bSJed Brown } 3929566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed)); 393c4762a1bSJed Brown if (i < xs + xm) { 394c2a16db8SHong Zhang for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hs; 395c4762a1bSJed Brown islow++; 396c4762a1bSJed Brown } 397c4762a1bSJed Brown } else if (i == fs + 1) { 398c4762a1bSJed Brown u = &ctx->u[0]; 399c4762a1bSJed Brown alpha[0] = hs * hs / (hf + hs) / (hf + hs + hs); 400c4762a1bSJed Brown gamma[0] = hs * (hf + hs) / (hs + hs) / (hf + hs + hs); 401c4762a1bSJed Brown for (j = 0; j < dof; j++) { 402c4762a1bSJed Brown r[j] = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 403c4762a1bSJed Brown min[j] = PetscMin(r[j], 2.0); 404c4762a1bSJed Brown u[j] = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 405c4762a1bSJed Brown } 4069566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed)); 407c4762a1bSJed Brown if (i > xs) { 408c2a16db8SHong Zhang for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hs; 409c4762a1bSJed Brown } 410c4762a1bSJed Brown if (i < xs + xm) { 411c2a16db8SHong Zhang for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hs; 412c2a16db8SHong Zhang islow++; 413c4762a1bSJed Brown } 414c4762a1bSJed Brown } else if (i > fs + 1) { 415c4762a1bSJed Brown u = &ctx->u[0]; 416c4762a1bSJed Brown alpha[0] = 1.0 / 6.0; 417c4762a1bSJed Brown gamma[0] = 1.0 / 3.0; 418c4762a1bSJed Brown for (j = 0; j < dof; j++) { 419c4762a1bSJed Brown r[j] = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 420c4762a1bSJed Brown min[j] = PetscMin(r[j], 2.0); 421c4762a1bSJed Brown u[j] = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 422c4762a1bSJed Brown } 4239566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed)); 424c4762a1bSJed Brown if (i > xs) { 425c2a16db8SHong Zhang for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hs; 426c4762a1bSJed Brown } 427c4762a1bSJed Brown if (i < xs + xm) { 428c2a16db8SHong Zhang for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hs; 429c2a16db8SHong Zhang islow++; 430c4762a1bSJed Brown } 431c4762a1bSJed Brown } 432c4762a1bSJed Brown } 4339566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); 4349566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 4359566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc)); 4369566063dSJacob Faibussowitsch PetscCall(PetscFree4(r, min, alpha, gamma)); 4373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 438c4762a1bSJed Brown } 439c4762a1bSJed Brown 440d71ae5a4SJacob Faibussowitsch static PetscErrorCode FVRHSFunctionfast(TS ts, PetscReal time, Vec X, Vec F, void *vctx) 441d71ae5a4SJacob Faibussowitsch { 442c4762a1bSJed Brown FVCtx *ctx = (FVCtx *)vctx; 443c4762a1bSJed Brown PetscInt i, j, Mx, dof, xs, xm, ifast = 0, sf = ctx->sf, fs = ctx->fs; 444c4762a1bSJed Brown PetscReal hf, hs; 445c4762a1bSJed Brown PetscScalar *x, *f, *r, *min, *alpha, *gamma; 446c4762a1bSJed Brown Vec Xloc; 447c4762a1bSJed Brown DM da; 448c4762a1bSJed Brown 449c4762a1bSJed Brown PetscFunctionBeginUser; 4509566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 4519566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc)); /* Xloc contains ghost points */ 4529566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); /* Mx is the number of center points */ 453c4762a1bSJed Brown hs = (ctx->xmax - ctx->xmin) / 2.0 * (ctx->hratio + 1.0) / Mx; 454c4762a1bSJed Brown hf = (ctx->xmax - ctx->xmin) / 2.0 * (1.0 + 1.0 / ctx->hratio) / Mx; 4559566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); /* X is solution vector which does not contain ghost points */ 4569566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); 457*dd8e379bSPierre Jolivet PetscCall(VecZeroEntries(F)); /* F is the right-hand side function corresponds to center points */ 4589566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x)); 4599566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 4609566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 4619566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(dof, &r, dof, &min, dof, &alpha, dof, &gamma)); 462c4762a1bSJed Brown 463c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 464c4762a1bSJed Brown for (i = xs - 2; i < 0; i++) { 465c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; 466c4762a1bSJed Brown } 467c4762a1bSJed Brown for (i = Mx; i < xs + xm + 2; i++) { 468c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; 469c4762a1bSJed Brown } 470c4762a1bSJed Brown } 471c4762a1bSJed Brown 472c4762a1bSJed Brown for (i = xs; i < xs + xm + 1; i++) { 473c4762a1bSJed Brown PetscReal maxspeed; 474c4762a1bSJed Brown PetscScalar *u; 475c4762a1bSJed Brown if (i == sf) { 476c4762a1bSJed Brown u = &ctx->u[0]; 477c4762a1bSJed Brown alpha[0] = hs * hf / (hs + hs) / (hs + hs + hf); 478c4762a1bSJed Brown gamma[0] = hs * (hs + hs) / (hs + hf) / (hs + hs + hf); 479c4762a1bSJed Brown for (j = 0; j < dof; j++) { 480c4762a1bSJed Brown r[j] = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 481c4762a1bSJed Brown min[j] = PetscMin(r[j], 2.0); 482c4762a1bSJed Brown u[j] = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 483c4762a1bSJed Brown } 4849566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed)); 485c4762a1bSJed Brown if (i < xs + xm) { 486c2a16db8SHong Zhang for (j = 0; j < dof; j++) f[ifast * dof + j] += ctx->flux[j] / hf; 487c4762a1bSJed Brown ifast++; 488c4762a1bSJed Brown } 489c4762a1bSJed Brown } else if (i == sf + 1) { 490c4762a1bSJed Brown u = &ctx->u[0]; 491c4762a1bSJed Brown alpha[0] = hf * hf / (hs + hf) / (hs + hf + hf); 492c4762a1bSJed Brown gamma[0] = hf * (hs + hf) / (hf + hf) / (hs + hf + hf); 493c4762a1bSJed Brown for (j = 0; j < dof; j++) { 494c4762a1bSJed Brown r[j] = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 495c4762a1bSJed Brown min[j] = PetscMin(r[j], 2.0); 496c4762a1bSJed Brown u[j] = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 497c4762a1bSJed Brown } 4989566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed)); 499c4762a1bSJed Brown if (i > xs) { 500c2a16db8SHong Zhang for (j = 0; j < dof; j++) f[(ifast - 1) * dof + j] -= ctx->flux[j] / hf; 501c4762a1bSJed Brown } 502c4762a1bSJed Brown if (i < xs + xm) { 503c2a16db8SHong Zhang for (j = 0; j < dof; j++) f[ifast * dof + j] += ctx->flux[j] / hf; 504c4762a1bSJed Brown ifast++; 505c4762a1bSJed Brown } 506c4762a1bSJed Brown } else if (i > sf + 1 && i < fs) { 507c4762a1bSJed Brown u = &ctx->u[0]; 508c4762a1bSJed Brown alpha[0] = 1.0 / 6.0; 509c4762a1bSJed Brown gamma[0] = 1.0 / 3.0; 510c4762a1bSJed Brown for (j = 0; j < dof; j++) { 511c4762a1bSJed Brown r[j] = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 512c4762a1bSJed Brown min[j] = PetscMin(r[j], 2.0); 513c4762a1bSJed Brown u[j] = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 514c4762a1bSJed Brown } 5159566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed)); 516c4762a1bSJed Brown if (i > xs) { 517c2a16db8SHong Zhang for (j = 0; j < dof; j++) f[(ifast - 1) * dof + j] -= ctx->flux[j] / hf; 518c4762a1bSJed Brown } 519c4762a1bSJed Brown if (i < xs + xm) { 520c2a16db8SHong Zhang for (j = 0; j < dof; j++) f[ifast * dof + j] += ctx->flux[j] / hf; 521c4762a1bSJed Brown ifast++; 522c4762a1bSJed Brown } 523c4762a1bSJed Brown } else if (i == fs) { 524c4762a1bSJed Brown u = &ctx->u[0]; 525c4762a1bSJed Brown alpha[0] = hf * hs / (hf + hf) / (hf + hf + hs); 526c4762a1bSJed Brown gamma[0] = hf * (hf + hf) / (hf + hs) / (hf + hf + hs); 527c4762a1bSJed Brown for (j = 0; j < dof; j++) { 528c4762a1bSJed Brown r[j] = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 529c4762a1bSJed Brown min[j] = PetscMin(r[j], 2.0); 530c4762a1bSJed Brown u[j] = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]); 531c4762a1bSJed Brown } 5329566063dSJacob Faibussowitsch PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed)); 533c4762a1bSJed Brown if (i > xs) { 534c2a16db8SHong Zhang for (j = 0; j < dof; j++) f[(ifast - 1) * dof + j] -= ctx->flux[j] / hf; 535c4762a1bSJed Brown } 536c4762a1bSJed Brown } 537c4762a1bSJed Brown } 5389566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); 5399566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 5409566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc)); 5419566063dSJacob Faibussowitsch PetscCall(PetscFree4(r, min, alpha, gamma)); 5423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 543c4762a1bSJed Brown } 544c4762a1bSJed Brown 545c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */ 546c4762a1bSJed Brown 547d71ae5a4SJacob Faibussowitsch PetscErrorCode FVSample(FVCtx *ctx, DM da, PetscReal time, Vec U) 548d71ae5a4SJacob Faibussowitsch { 549c4762a1bSJed Brown PetscScalar *u, *uj, xj, xi; 550c4762a1bSJed Brown PetscInt i, j, k, dof, xs, xm, Mx, count_slow, count_fast; 551c4762a1bSJed Brown const PetscInt N = 200; 552c4762a1bSJed Brown 553c4762a1bSJed Brown PetscFunctionBeginUser; 5543c633725SBarry Smith PetscCheck(ctx->physics.sample, PETSC_COMM_SELF, PETSC_ERR_SUP, "Physics has not provided a sampling function"); 5559566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); 5569566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 5579566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u)); 5589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof, &uj)); 559c4762a1bSJed Brown const PetscReal hs = (ctx->xmax - ctx->xmin) / 2.0 * (ctx->hratio + 1.0) / Mx; 560c4762a1bSJed Brown const PetscReal hf = (ctx->xmax - ctx->xmin) / 2.0 * (1.0 + 1.0 / ctx->hratio) / Mx; 561c4762a1bSJed Brown count_slow = Mx / (1 + ctx->hratio); 562c4762a1bSJed Brown count_fast = Mx - count_slow; 563c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 564c4762a1bSJed Brown if (i * hs + 0.5 * hs < (ctx->xmax - ctx->xmin) * 0.25) { 565c4762a1bSJed Brown xi = ctx->xmin + 0.5 * hs + i * hs; 566c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 567c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] = 0; 568c4762a1bSJed Brown for (j = 0; j < N + 1; j++) { 569c4762a1bSJed Brown xj = xi + hs * (j - N / 2) / (PetscReal)N; 5709566063dSJacob Faibussowitsch PetscCall((*ctx->physics.sample)(ctx->physics.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj)); 571c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N; 572c4762a1bSJed Brown } 573c4762a1bSJed Brown } else if ((ctx->xmax - ctx->xmin) * 0.25 + (i - count_slow / 2) * hf + 0.5 * hf < (ctx->xmax - ctx->xmin) * 0.75) { 574c4762a1bSJed Brown xi = ctx->xmin + (ctx->xmax - ctx->xmin) * 0.25 + 0.5 * hf + (i - count_slow / 2) * hf; 575c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 576c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] = 0; 577c4762a1bSJed Brown for (j = 0; j < N + 1; j++) { 578c4762a1bSJed Brown xj = xi + hf * (j - N / 2) / (PetscReal)N; 5799566063dSJacob Faibussowitsch PetscCall((*ctx->physics.sample)(ctx->physics.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj)); 580c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N; 581c4762a1bSJed Brown } 582c4762a1bSJed Brown } else { 583c4762a1bSJed Brown xi = ctx->xmin + (ctx->xmax - ctx->xmin) * 0.75 + 0.5 * hs + (i - count_slow / 2 - count_fast) * hs; 584c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 585c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] = 0; 586c4762a1bSJed Brown for (j = 0; j < N + 1; j++) { 587c4762a1bSJed Brown xj = xi + hs * (j - N / 2) / (PetscReal)N; 5889566063dSJacob Faibussowitsch PetscCall((*ctx->physics.sample)(ctx->physics.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj)); 589c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N; 590c4762a1bSJed Brown } 591c4762a1bSJed Brown } 592c4762a1bSJed Brown } 5939566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u)); 5949566063dSJacob Faibussowitsch PetscCall(PetscFree(uj)); 5953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 596c4762a1bSJed Brown } 597c4762a1bSJed Brown 598d71ae5a4SJacob Faibussowitsch static PetscErrorCode SolutionStatsView(DM da, Vec X, PetscViewer viewer) 599d71ae5a4SJacob Faibussowitsch { 600c4762a1bSJed Brown PetscReal xmin, xmax; 601c4762a1bSJed Brown PetscScalar sum, tvsum, tvgsum; 602c4762a1bSJed Brown const PetscScalar *x; 603c4762a1bSJed Brown PetscInt imin, imax, Mx, i, j, xs, xm, dof; 604c4762a1bSJed Brown Vec Xloc; 605c4762a1bSJed Brown PetscBool iascii; 606c4762a1bSJed Brown 607c4762a1bSJed Brown PetscFunctionBeginUser; 6089566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 609c4762a1bSJed Brown if (iascii) { 610c4762a1bSJed Brown /* PETSc lacks a function to compute total variation norm (difficult in multiple dimensions), we do it here */ 6119566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc)); 6129566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); 6139566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); 6149566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, Xloc, (void *)&x)); 6159566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 6169566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); 617c4762a1bSJed Brown tvsum = 0; 618c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 619c4762a1bSJed Brown for (j = 0; j < dof; j++) tvsum += PetscAbsScalar(x[i * dof + j] - x[(i - 1) * dof + j]); 620c4762a1bSJed Brown } 621712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(&tvsum, &tvgsum, 1, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)da))); 6229566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, Xloc, (void *)&x)); 6239566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc)); 624c4762a1bSJed Brown 6259566063dSJacob Faibussowitsch PetscCall(VecMin(X, &imin, &xmin)); 6269566063dSJacob Faibussowitsch PetscCall(VecMax(X, &imax, &xmax)); 6279566063dSJacob Faibussowitsch PetscCall(VecSum(X, &sum)); 62863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Solution range [%g,%g] with minimum at %" PetscInt_FMT ", mean %g, ||x||_TV %g\n", (double)xmin, (double)xmax, imin, (double)(sum / Mx), (double)(tvgsum / Mx))); 629c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Viewer type not supported"); 6303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 631c4762a1bSJed Brown } 632c4762a1bSJed Brown 633d71ae5a4SJacob Faibussowitsch static PetscErrorCode SolutionErrorNorms(FVCtx *ctx, DM da, PetscReal t, Vec X, PetscReal *nrm1) 634d71ae5a4SJacob Faibussowitsch { 635c4762a1bSJed Brown Vec Y; 636c4762a1bSJed Brown PetscInt i, Mx, count_slow = 0, count_fast = 0; 637c4762a1bSJed Brown const PetscScalar *ptr_X, *ptr_Y; 638c4762a1bSJed Brown 639c4762a1bSJed Brown PetscFunctionBeginUser; 6409566063dSJacob Faibussowitsch PetscCall(VecGetSize(X, &Mx)); 6419566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &Y)); 6429566063dSJacob Faibussowitsch PetscCall(FVSample(ctx, da, t, Y)); 643c4762a1bSJed Brown const PetscReal hs = (ctx->xmax - ctx->xmin) / 2.0 * (ctx->hratio + 1.0) / Mx; 644c4762a1bSJed Brown const PetscReal hf = (ctx->xmax - ctx->xmin) / 2.0 * (1.0 + 1.0 / ctx->hratio) / Mx; 645c4762a1bSJed Brown count_slow = (PetscReal)Mx / (1.0 + ctx->hratio); 646c4762a1bSJed Brown count_fast = Mx - count_slow; 6479566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &ptr_X)); 6489566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &ptr_Y)); 649c4762a1bSJed Brown for (i = 0; i < Mx; i++) { 650c4762a1bSJed Brown if (i < count_slow / 2 || i > count_slow / 2 + count_fast - 1) *nrm1 += hs * PetscAbs(ptr_X[i] - ptr_Y[i]); 651c4762a1bSJed Brown else *nrm1 += hf * PetscAbs(ptr_X[i] - ptr_Y[i]); 652c4762a1bSJed Brown } 6539566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &ptr_X)); 6549566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &ptr_Y)); 6559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 6563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 657c4762a1bSJed Brown } 658c4762a1bSJed Brown 659d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[]) 660d71ae5a4SJacob Faibussowitsch { 661c4762a1bSJed Brown char physname[256] = "advect", final_fname[256] = "solution.m"; 662c4762a1bSJed Brown PetscFunctionList physics = 0; 663c4762a1bSJed Brown MPI_Comm comm; 664c4762a1bSJed Brown TS ts; 665c4762a1bSJed Brown DM da; 666c4762a1bSJed Brown Vec X, X0, R; 667c4762a1bSJed Brown FVCtx ctx; 668c4762a1bSJed Brown PetscInt i, k, dof, xs, xm, Mx, draw = 0, count_slow, count_fast, islow = 0, ifast = 0, *index_slow, *index_fast; 669c4762a1bSJed Brown PetscBool view_final = PETSC_FALSE; 670c4762a1bSJed Brown PetscReal ptime; 671c4762a1bSJed Brown 672327415f7SBarry Smith PetscFunctionBeginUser; 6739566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, 0, help)); 674c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 6759566063dSJacob Faibussowitsch PetscCall(PetscMemzero(&ctx, sizeof(ctx))); 676c4762a1bSJed Brown 677c4762a1bSJed Brown /* Register physical models to be available on the command line */ 6789566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&physics, "advect", PhysicsCreate_Advect)); 679c4762a1bSJed Brown 680c4762a1bSJed Brown ctx.comm = comm; 681c4762a1bSJed Brown ctx.cfl = 0.9; 682c4762a1bSJed Brown ctx.bctype = FVBC_PERIODIC; 683c4762a1bSJed Brown ctx.xmin = -1.0; 684c4762a1bSJed Brown ctx.xmax = 1.0; 685d0609cedSBarry Smith PetscOptionsBegin(comm, NULL, "Finite Volume solver options", ""); 6869566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-xmin", "X min", "", ctx.xmin, &ctx.xmin, NULL)); 6879566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-xmax", "X max", "", ctx.xmax, &ctx.xmax, NULL)); 6889566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-draw", "Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)", "", draw, &draw, NULL)); 6899566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-view_final", "Write final solution in ASCII MATLAB format to given file name", "", final_fname, final_fname, sizeof(final_fname), &view_final)); 6909566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-initial", "Initial condition (depends on the physics)", "", ctx.initial, &ctx.initial, NULL)); 6919566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-exact", "Compare errors with exact solution", "", ctx.exact, &ctx.exact, NULL)); 6929566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-simulation", "Compare errors with reference solution", "", ctx.simulation, &ctx.simulation, NULL)); 6939566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-cfl", "CFL number to time step at", "", ctx.cfl, &ctx.cfl, NULL)); 6949566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-bc_type", "Boundary condition", "", FVBCTypes, (PetscEnum)ctx.bctype, (PetscEnum *)&ctx.bctype, NULL)); 6959566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-hratio", "Spacing ratio", "", ctx.hratio, &ctx.hratio, NULL)); 696d0609cedSBarry Smith PetscOptionsEnd(); 697c4762a1bSJed Brown 698c4762a1bSJed Brown /* Choose the physics from the list of registered models */ 699c4762a1bSJed Brown { 700c4762a1bSJed Brown PetscErrorCode (*r)(FVCtx *); 7019566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(physics, physname, &r)); 7023c633725SBarry Smith PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Physics '%s' not found", physname); 703c4762a1bSJed Brown /* Create the physics, will set the number of fields and their names */ 7049566063dSJacob Faibussowitsch PetscCall((*r)(&ctx)); 705c4762a1bSJed Brown } 706c4762a1bSJed Brown 707c4762a1bSJed Brown /* Create a DMDA to manage the parallel grid */ 7089566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(comm, DM_BOUNDARY_PERIODIC, 50, ctx.physics.dof, 2, NULL, &da)); 7099566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 7109566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 711c4762a1bSJed Brown /* Inform the DMDA of the field names provided by the physics. */ 712c4762a1bSJed Brown /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */ 71348a46eb9SPierre Jolivet for (i = 0; i < ctx.physics.dof; i++) PetscCall(DMDASetFieldName(da, i, ctx.physics.fieldname[i])); 7149566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); 7159566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 716c4762a1bSJed Brown 717c4762a1bSJed Brown /* Set coordinates of cell centers */ 7189566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, ctx.xmin + 0.5 * (ctx.xmax - ctx.xmin) / Mx, ctx.xmax + 0.5 * (ctx.xmax - ctx.xmin) / Mx, 0, 0, 0, 0)); 719c4762a1bSJed Brown 720c4762a1bSJed Brown /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */ 7219566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(dof, &ctx.u, dof, &ctx.flux, dof, &ctx.speeds)); 722c4762a1bSJed Brown 723c4762a1bSJed Brown /* Create a vector to store the solution and to save the initial state */ 7249566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &X)); 7259566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &X0)); 7269566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &R)); 727c4762a1bSJed Brown 728c4762a1bSJed Brown /* create index for slow parts and fast parts*/ 729c4762a1bSJed Brown count_slow = Mx / (1 + ctx.hratio); 73008401ef6SPierre Jolivet PetscCheck(count_slow % 2 == 0, PETSC_COMM_WORLD, PETSC_ERR_USER, "Please adjust grid size Mx (-da_grid_x) and hratio (-hratio) so that Mx/(1+hartio) is even"); 731c4762a1bSJed Brown count_fast = Mx - count_slow; 732c4762a1bSJed Brown ctx.sf = count_slow / 2; 733c4762a1bSJed Brown ctx.fs = ctx.sf + count_fast; 7349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(xm * dof, &index_slow)); 7359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(xm * dof, &index_fast)); 736c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 737c4762a1bSJed Brown if (i < count_slow / 2 || i > count_slow / 2 + count_fast - 1) 738c4762a1bSJed Brown for (k = 0; k < dof; k++) index_slow[islow++] = i * dof + k; 739c4762a1bSJed Brown else 740c4762a1bSJed Brown for (k = 0; k < dof; k++) index_fast[ifast++] = i * dof + k; 741c4762a1bSJed Brown } 7429566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islow, index_slow, PETSC_COPY_VALUES, &ctx.iss)); 7439566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, ifast, index_fast, PETSC_COPY_VALUES, &ctx.isf)); 744c4762a1bSJed Brown 745c4762a1bSJed Brown /* Create a time-stepping object */ 7469566063dSJacob Faibussowitsch PetscCall(TSCreate(comm, &ts)); 7479566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da)); 7489566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, R, FVRHSFunction, &ctx)); 7499566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts, "slow", ctx.iss)); 7509566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts, "fast", ctx.isf)); 7519566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts, "slow", NULL, FVRHSFunctionslow, &ctx)); 7529566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts, "fast", NULL, FVRHSFunctionfast, &ctx)); 753c4762a1bSJed Brown 7549566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSMPRK)); 7559566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 10)); 7569566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 757c4762a1bSJed Brown 758c4762a1bSJed Brown /* Compute initial conditions and starting time step */ 7599566063dSJacob Faibussowitsch PetscCall(FVSample(&ctx, da, 0, X0)); 7609566063dSJacob Faibussowitsch PetscCall(FVRHSFunction(ts, 0, X0, X, (void *)&ctx)); /* Initial function evaluation, only used to determine max speed */ 7619566063dSJacob Faibussowitsch PetscCall(VecCopy(X0, X)); /* The function value was not used so we set X=X0 again */ 7629566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, ctx.cfl / ctx.cfl_idt)); 7639566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); /* Take runtime options */ 7649566063dSJacob Faibussowitsch PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD)); 765c4762a1bSJed Brown { 766c4762a1bSJed Brown PetscInt steps; 767c4762a1bSJed Brown PetscScalar mass_initial, mass_final, mass_difference, mass_differenceg; 768c4762a1bSJed Brown const PetscScalar *ptr_X, *ptr_X0; 769c2a16db8SHong Zhang const PetscReal hs = (ctx.xmax - ctx.xmin) / 2.0 / count_slow; 770c2a16db8SHong Zhang const PetscReal hf = (ctx.xmax - ctx.xmin) / 2.0 / count_fast; 7719566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, X)); 7729566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ptime)); 7739566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 774c4762a1bSJed Brown /* calculate the total mass at initial time and final time */ 775c4762a1bSJed Brown mass_initial = 0.0; 776c4762a1bSJed Brown mass_final = 0.0; 7779566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, X0, (void *)&ptr_X0)); 7789566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, X, (void *)&ptr_X)); 779c2a16db8SHong Zhang for (i = xs; i < xs + xm; i++) { 780c2a16db8SHong Zhang if (i < ctx.sf || i > ctx.fs - 1) { 781c2a16db8SHong Zhang for (k = 0; k < dof; k++) { 782c2a16db8SHong Zhang mass_initial = mass_initial + hs * ptr_X0[i * dof + k]; 783c2a16db8SHong Zhang mass_final = mass_final + hs * ptr_X[i * dof + k]; 784c2a16db8SHong Zhang } 785c4762a1bSJed Brown } else { 786c2a16db8SHong Zhang for (k = 0; k < dof; k++) { 787c2a16db8SHong Zhang mass_initial = mass_initial + hf * ptr_X0[i * dof + k]; 788c2a16db8SHong Zhang mass_final = mass_final + hf * ptr_X[i * dof + k]; 789c2a16db8SHong Zhang } 790c4762a1bSJed Brown } 791c4762a1bSJed Brown } 7929566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, X0, (void *)&ptr_X0)); 7939566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, X, (void *)&ptr_X)); 794c4762a1bSJed Brown mass_difference = mass_final - mass_initial; 795712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(&mass_difference, &mass_differenceg, 1, MPIU_SCALAR, MPIU_SUM, comm)); 7969566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Mass difference %g\n", (double)mass_differenceg)); 79763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Final time %g, steps %" PetscInt_FMT "\n", (double)ptime, steps)); 798c4762a1bSJed Brown if (ctx.exact) { 799c4762a1bSJed Brown PetscReal nrm1 = 0; 8009566063dSJacob Faibussowitsch PetscCall(SolutionErrorNorms(&ctx, da, ptime, X, &nrm1)); 8019566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %g\n", (double)nrm1)); 802c4762a1bSJed Brown } 803c4762a1bSJed Brown if (ctx.simulation) { 804c4762a1bSJed Brown PetscReal nrm1 = 0; 805c4762a1bSJed Brown PetscViewer fd; 806c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN] = "binaryoutput"; 807c4762a1bSJed Brown Vec XR; 808c4762a1bSJed Brown PetscBool flg; 809c4762a1bSJed Brown const PetscScalar *ptr_XR; 8109566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f", filename, sizeof(filename), &flg)); 8113c633725SBarry Smith PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate binary file with the -f option"); 8129566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, filename, FILE_MODE_READ, &fd)); 8139566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X0, &XR)); 8149566063dSJacob Faibussowitsch PetscCall(VecLoad(XR, fd)); 8159566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 8169566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &ptr_X)); 8179566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XR, &ptr_XR)); 818c4762a1bSJed Brown for (i = 0; i < Mx; i++) { 819c4762a1bSJed Brown if (i < count_slow / 2 || i > count_slow / 2 + count_fast - 1) nrm1 = nrm1 + hs * PetscAbs(ptr_X[i] - ptr_XR[i]); 820c4762a1bSJed Brown else nrm1 = nrm1 + hf * PetscAbs(ptr_X[i] - ptr_XR[i]); 821c4762a1bSJed Brown } 8229566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &ptr_X)); 8239566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XR, &ptr_XR)); 8249566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %g\n", (double)nrm1)); 8259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XR)); 826c4762a1bSJed Brown } 827c4762a1bSJed Brown } 828c4762a1bSJed Brown 8299566063dSJacob Faibussowitsch PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD)); 8309566063dSJacob Faibussowitsch if (draw & 0x1) PetscCall(VecView(X0, PETSC_VIEWER_DRAW_WORLD)); 8319566063dSJacob Faibussowitsch if (draw & 0x2) PetscCall(VecView(X, PETSC_VIEWER_DRAW_WORLD)); 832c4762a1bSJed Brown if (draw & 0x4) { 833c4762a1bSJed Brown Vec Y; 8349566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &Y)); 8359566063dSJacob Faibussowitsch PetscCall(FVSample(&ctx, da, ptime, Y)); 8369566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y, -1, X)); 8379566063dSJacob Faibussowitsch PetscCall(VecView(Y, PETSC_VIEWER_DRAW_WORLD)); 8389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 839c4762a1bSJed Brown } 840c4762a1bSJed Brown 841c4762a1bSJed Brown if (view_final) { 842c4762a1bSJed Brown PetscViewer viewer; 8439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, final_fname, &viewer)); 8449566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB)); 8459566063dSJacob Faibussowitsch PetscCall(VecView(X, viewer)); 8469566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 8479566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 848c4762a1bSJed Brown } 849c4762a1bSJed Brown 850c4762a1bSJed Brown /* Clean up */ 8519566063dSJacob Faibussowitsch PetscCall((*ctx.physics.destroy)(ctx.physics.user)); 8529566063dSJacob Faibussowitsch for (i = 0; i < ctx.physics.dof; i++) PetscCall(PetscFree(ctx.physics.fieldname[i])); 8539566063dSJacob Faibussowitsch PetscCall(PetscFree3(ctx.u, ctx.flux, ctx.speeds)); 8549566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.iss)); 8559566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.isf)); 8569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 8579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X0)); 8589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&R)); 8599566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 8609566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 8619566063dSJacob Faibussowitsch PetscCall(PetscFree(index_slow)); 8629566063dSJacob Faibussowitsch PetscCall(PetscFree(index_fast)); 8639566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&physics)); 8649566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 865b122ec5aSJacob Faibussowitsch return 0; 866c4762a1bSJed Brown } 867c4762a1bSJed Brown 868c4762a1bSJed Brown /*TEST 869c4762a1bSJed Brown 870c4762a1bSJed Brown build: 871f56ea12dSJed Brown requires: !complex 872c4762a1bSJed Brown 873c4762a1bSJed Brown test: 874c4762a1bSJed Brown args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 0 875c4762a1bSJed Brown 876c4762a1bSJed Brown test: 877c4762a1bSJed Brown suffix: 2 878c4762a1bSJed Brown args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 1 879c4762a1bSJed Brown output_file: output/ex7_1.out 880c4762a1bSJed Brown 881c4762a1bSJed Brown test: 882c4762a1bSJed Brown suffix: 3 883c4762a1bSJed Brown args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 0 884c4762a1bSJed Brown 885c4762a1bSJed Brown test: 886c4762a1bSJed Brown suffix: 4 887c4762a1bSJed Brown args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 1 888c4762a1bSJed Brown output_file: output/ex7_3.out 889c2a16db8SHong Zhang 890c2a16db8SHong Zhang test: 891c2a16db8SHong Zhang suffix: 5 892c2a16db8SHong Zhang nsize: 2 893c2a16db8SHong Zhang args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 1 894c2a16db8SHong Zhang output_file: output/ex7_3.out 895c4762a1bSJed Brown TEST*/ 896