static const char help[] = "1D periodic Finite Volume solver in slope-limiter form with semidiscrete time stepping.\n" " advection - Constant coefficient scalar advection\n" " u_t + (a*u)_x = 0\n" " for this toy problem, we choose different meshsizes for different sub-domains (slow-medium-fast-medium-slow), \n" " the meshsize ratio between two adjacient sub-domains is controlled with -hratio,\n" " exact - Exact Riemann solver which usually needs to perform a Newton iteration to connect\n" " the states across shocks and rarefactions\n" " simulation - use reference solution which is generated by smaller time step size to be true solution,\n" " also the reference solution should be generated by user and stored in a binary file.\n" " characteristic - Limit the characteristic variables, this is usually preferred (default)\n" "Several initial conditions can be chosen with -initial N\n\n" "The problem size should be set with -da_grid_x M\n\n"; #include #include #include #include #include "finitevolume1d.h" static inline PetscReal RangeMod(PetscReal a, PetscReal xmin, PetscReal xmax) { PetscReal range = xmax - xmin; return xmin + PetscFmodReal(range + PetscFmodReal(a, range), range); } /* --------------------------------- Advection ----------------------------------- */ typedef struct { PetscReal a; /* advective velocity */ } AdvectCtx; static PetscErrorCode PhysicsRiemann_Advect(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed) { AdvectCtx *ctx = (AdvectCtx *)vctx; PetscReal speed; PetscFunctionBeginUser; speed = ctx->a; flux[0] = PetscMax(0, speed) * uL[0] + PetscMin(0, speed) * uR[0]; *maxspeed = speed; PetscFunctionReturn(PETSC_SUCCESS); } static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx, PetscInt m, const PetscScalar *u, PetscScalar *X, PetscScalar *Xi, PetscReal *speeds) { AdvectCtx *ctx = (AdvectCtx *)vctx; PetscFunctionBeginUser; X[0] = 1.; Xi[0] = 1.; speeds[0] = ctx->a; PetscFunctionReturn(PETSC_SUCCESS); } static PetscErrorCode PhysicsSample_Advect(void *vctx, PetscInt initial, FVBCType bctype, PetscReal xmin, PetscReal xmax, PetscReal t, PetscReal x, PetscReal *u) { AdvectCtx *ctx = (AdvectCtx *)vctx; PetscReal a = ctx->a, x0; PetscFunctionBeginUser; switch (bctype) { case FVBC_OUTFLOW: x0 = x - a * t; break; case FVBC_PERIODIC: x0 = RangeMod(x - a * t, xmin, xmax); break; default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown BCType"); } switch (initial) { case 0: u[0] = (x0 < 0) ? 1 : -1; break; case 1: u[0] = (x0 < 0) ? -1 : 1; break; case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break; case 3: u[0] = PetscSinReal(2 * PETSC_PI * x0); break; case 4: u[0] = PetscAbs(x0); break; case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2 * PETSC_PI * x0)); break; case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2 - x0 : 0)); break; case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI * x0), 10.0); break; default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition"); } PetscFunctionReturn(PETSC_SUCCESS); } static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx) { AdvectCtx *user; PetscFunctionBeginUser; PetscCall(PetscNew(&user)); ctx->physics2.sample2 = PhysicsSample_Advect; ctx->physics2.riemann2 = PhysicsRiemann_Advect; ctx->physics2.characteristic2 = PhysicsCharacteristic_Advect; ctx->physics2.destroy = PhysicsDestroy_SimpleFree; ctx->physics2.user = user; ctx->physics2.dof = 1; PetscCall(PetscStrallocpy("u", &ctx->physics2.fieldname[0])); user->a = 1; PetscOptionsBegin(ctx->comm, ctx->prefix, "Options for advection", ""); { PetscCall(PetscOptionsReal("-physics_advect_a", "Speed", "", user->a, &user->a, NULL)); } PetscOptionsEnd(); PetscFunctionReturn(PETSC_SUCCESS); } PetscErrorCode FVSample_3WaySplit(FVCtx *ctx, DM da, PetscReal time, Vec U) { PetscScalar *u, *uj, xj, xi; PetscInt i, j, k, dof, xs, xm, Mx; const PetscInt N = 200; PetscReal hs, hm, hf; PetscFunctionBeginUser; PetscCheck(ctx->physics2.sample2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Physics has not provided a sampling function"); PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); PetscCall(DMDAVecGetArray(da, U, &u)); PetscCall(PetscMalloc1(dof, &uj)); hs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm; hm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm); hf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf); for (i = xs; i < xs + xm; i++) { if (i < ctx->sm) { xi = ctx->xmin + 0.5 * hs + i * hs; /* Integrate over cell i using trapezoid rule with N points. */ for (k = 0; k < dof; k++) u[i * dof + k] = 0; for (j = 0; j < N + 1; j++) { xj = xi + hs * (j - N / 2) / (PetscReal)N; PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj)); for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N; } } else if (i < ctx->mf) { xi = ctx->xmin + ctx->sm * hs + 0.5 * hm + (i - ctx->sm) * hm; /* Integrate over cell i using trapezoid rule with N points. */ for (k = 0; k < dof; k++) u[i * dof + k] = 0; for (j = 0; j < N + 1; j++) { xj = xi + hm * (j - N / 2) / (PetscReal)N; PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj)); for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N; } } else if (i < ctx->fm) { xi = ctx->xmin + ctx->sm * hs + (ctx->mf - ctx->sm) * hm + 0.5 * hf + (i - ctx->mf) * hf; /* Integrate over cell i using trapezoid rule with N points. */ for (k = 0; k < dof; k++) u[i * dof + k] = 0; for (j = 0; j < N + 1; j++) { xj = xi + hf * (j - N / 2) / (PetscReal)N; PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj)); for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N; } } else if (i < ctx->ms) { xi = ctx->xmin + ctx->sm * hs + (ctx->mf - ctx->sm) * hm + (ctx->fm - ctx->mf) * hf + 0.5 * hm + (i - ctx->fm) * hm; /* Integrate over cell i using trapezoid rule with N points. */ for (k = 0; k < dof; k++) u[i * dof + k] = 0; for (j = 0; j < N + 1; j++) { xj = xi + hm * (j - N / 2) / (PetscReal)N; PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj)); for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N; } } else { xi = ctx->xmin + ctx->sm * hs + (ctx->mf - ctx->sm) * hm + (ctx->fm - ctx->mf) * hf + (ctx->ms - ctx->fm) * hm + 0.5 * hs + (i - ctx->ms) * hs; /* Integrate over cell i using trapezoid rule with N points. */ for (k = 0; k < dof; k++) u[i * dof + k] = 0; for (j = 0; j < N + 1; j++) { xj = xi + hs * (j - N / 2) / (PetscReal)N; PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj)); for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N; } } } PetscCall(DMDAVecRestoreArray(da, U, &u)); PetscCall(PetscFree(uj)); PetscFunctionReturn(PETSC_SUCCESS); } static PetscErrorCode SolutionErrorNorms_3WaySplit(FVCtx *ctx, DM da, PetscReal t, Vec X, PetscReal *nrm1) { Vec Y; PetscInt i, Mx; const PetscScalar *ptr_X, *ptr_Y; PetscReal hs, hm, hf; PetscFunctionBeginUser; PetscCall(VecGetSize(X, &Mx)); hs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm; hm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm); hf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf); PetscCall(VecDuplicate(X, &Y)); PetscCall(FVSample_3WaySplit(ctx, da, t, Y)); PetscCall(VecGetArrayRead(X, &ptr_X)); PetscCall(VecGetArrayRead(Y, &ptr_Y)); for (i = 0; i < Mx; i++) { if (i < ctx->sm || i > ctx->ms - 1) *nrm1 += hs * PetscAbs(ptr_X[i] - ptr_Y[i]); else if (i < ctx->mf || i > ctx->fm - 1) *nrm1 += hm * PetscAbs(ptr_X[i] - ptr_Y[i]); else *nrm1 += hf * PetscAbs(ptr_X[i] - ptr_Y[i]); } PetscCall(VecRestoreArrayRead(X, &ptr_X)); PetscCall(VecRestoreArrayRead(Y, &ptr_Y)); PetscCall(VecDestroy(&Y)); PetscFunctionReturn(PETSC_SUCCESS); } PetscErrorCode FVRHSFunction_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) { FVCtx *ctx = (FVCtx *)vctx; PetscInt i, j, k, Mx, dof, xs, xm, sm = ctx->sm, mf = ctx->mf, fm = ctx->fm, ms = ctx->ms; PetscReal hxf, hxm, hxs, cfl_idt = 0; PetscScalar *x, *f, *slope; Vec Xloc; DM da; PetscFunctionBeginUser; PetscCall(TSGetDM(ts, &da)); PetscCall(DMGetLocalVector(da, &Xloc)); /* Xloc contains ghost points */ PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); /* Mx is the number of center points */ hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm; hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm); hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf); PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); /* X is solution vector which does not contain ghost points */ PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); PetscCall(VecZeroEntries(F)); /* F is the right-hand side function corresponds to center points */ PetscCall(DMDAVecGetArray(da, Xloc, &x)); PetscCall(DMDAVecGetArray(da, F, &f)); PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); /* contains ghost points */ PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); if (ctx->bctype == FVBC_OUTFLOW) { for (i = xs - 2; i < 0; i++) { for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; } for (i = Mx; i < xs + xm + 2; i++) { for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; } } for (i = xs - 1; i < xs + xm + 1; i++) { struct _LimitInfo info; PetscScalar *cjmpL, *cjmpR; /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds)); /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof)); cjmpL = &ctx->cjmpLR[0]; cjmpR = &ctx->cjmpLR[dof]; for (j = 0; j < dof; j++) { PetscScalar jmpL, jmpR; jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j]; jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j]; for (k = 0; k < dof; k++) { cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL; cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR; } } /* Apply limiter to the left and right characteristic jumps */ info.m = dof; info.hxs = hxs; info.hxm = hxm; info.hxf = hxf; (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope); for (j = 0; j < dof; j++) { PetscScalar tmp = 0; for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k]; slope[i * dof + j] = tmp; } } for (i = xs; i < xs + xm + 1; i++) { PetscReal maxspeed; PetscScalar *uL, *uR; uL = &ctx->uLR[0]; uR = &ctx->uLR[dof]; if (i < sm || i > ms) { /* slow region */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i > xs) { for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxs; } if (i < xs + xm) { for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxs; } } else if (i == sm) { /* interface between slow and medium component */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i > xs) { for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxs; } if (i < xs + xm) { for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxm; } } else if (i == ms) { /* interface between medium and slow regions */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i > xs) { for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxm; } if (i < xs + xm) { for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxs; } } else if (i < mf || i > fm) { /* medium region */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i > xs) { for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxm; } if (i < xs + xm) { for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxm; } } else if (i == mf) { /* interface between medium and fast regions */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i > xs) { for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxm; } if (i < xs + xm) { for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxf; } } else if (i == fm) { /* interface between fast and medium regions */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i > xs) { for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxf; } if (i < xs + xm) { for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxm; } } else { /* fast region */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hxs)); /* Max allowable value of 1/Delta t */ if (i > xs) { for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxf; } if (i < xs + xm) { for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxf; } } } PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); PetscCall(DMDAVecRestoreArray(da, F, &f)); PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope)); PetscCall(DMRestoreLocalVector(da, &Xloc)); PetscCall(MPIU_Allreduce(&cfl_idt, &ctx->cfl_idt, 1, MPIU_SCALAR, MPIU_MAX, PetscObjectComm((PetscObject)da))); if (0) { /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */ PetscReal dt, tnow; PetscCall(TSGetTimeStep(ts, &dt)); PetscCall(TSGetTime(ts, &tnow)); 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))); } PetscFunctionReturn(PETSC_SUCCESS); } /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */ PetscErrorCode FVRHSFunctionslow_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) { FVCtx *ctx = (FVCtx *)vctx; PetscInt i, j, k, Mx, dof, xs, xm, islow = 0, sm = ctx->sm, ms = ctx->ms, lsbwidth = ctx->lsbwidth, rsbwidth = ctx->rsbwidth; PetscReal hxs, hxm, hxf, cfl_idt = 0; PetscScalar *x, *f, *slope; Vec Xloc; DM da; PetscFunctionBeginUser; PetscCall(TSGetDM(ts, &da)); PetscCall(DMGetLocalVector(da, &Xloc)); PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm; hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm); hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf); PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); PetscCall(VecZeroEntries(F)); PetscCall(DMDAVecGetArray(da, Xloc, &x)); PetscCall(VecGetArray(F, &f)); PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); if (ctx->bctype == FVBC_OUTFLOW) { for (i = xs - 2; i < 0; i++) { for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; } for (i = Mx; i < xs + xm + 2; i++) { for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; } } for (i = xs - 1; i < xs + xm + 1; i++) { struct _LimitInfo info; PetscScalar *cjmpL, *cjmpR; if (i < sm - lsbwidth + 1 || i > ms + rsbwidth - 2) { /* slow components and the first and last fast components */ /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds)); /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof)); cjmpL = &ctx->cjmpLR[0]; cjmpR = &ctx->cjmpLR[dof]; for (j = 0; j < dof; j++) { PetscScalar jmpL, jmpR; jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j]; jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j]; for (k = 0; k < dof; k++) { cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL; cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR; } } /* Apply limiter to the left and right characteristic jumps */ info.m = dof; info.hxs = hxs; info.hxm = hxm; info.hxf = hxf; (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope); for (j = 0; j < dof; j++) { PetscScalar tmp = 0; for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k]; slope[i * dof + j] = tmp; } } } for (i = xs; i < xs + xm + 1; i++) { PetscReal maxspeed; PetscScalar *uL, *uR; uL = &ctx->uLR[0]; uR = &ctx->uLR[dof]; if (i < sm - lsbwidth) { /* slow region */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hxs)); /* Max allowable value of 1/Delta t */ if (i > xs) { for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs; } if (i < xs + xm) { for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs; islow++; } } if (i == sm - lsbwidth) { /* interface between slow and medium regions */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i > xs) { for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs; } } if (i == ms + rsbwidth) { /* interface between medium and slow regions */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i < xs + xm) { for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs; islow++; } } if (i > ms + rsbwidth) { /* slow region */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hxs)); /* Max allowable value of 1/Delta t */ if (i > xs) { for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs; } if (i < xs + xm) { for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs; islow++; } } } PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); PetscCall(VecRestoreArray(F, &f)); PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope)); PetscCall(DMRestoreLocalVector(da, &Xloc)); PetscCall(MPIU_Allreduce(&cfl_idt, &ctx->cfl_idt, 1, MPIU_SCALAR, MPIU_MAX, PetscObjectComm((PetscObject)da))); PetscFunctionReturn(PETSC_SUCCESS); } PetscErrorCode FVRHSFunctionslowbuffer_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) { FVCtx *ctx = (FVCtx *)vctx; PetscInt i, j, k, Mx, dof, xs, xm, islowbuffer = 0, sm = ctx->sm, ms = ctx->ms, lsbwidth = ctx->lsbwidth, rsbwidth = ctx->rsbwidth; PetscReal hxs, hxm, hxf; PetscScalar *x, *f, *slope; Vec Xloc; DM da; PetscFunctionBeginUser; PetscCall(TSGetDM(ts, &da)); PetscCall(DMGetLocalVector(da, &Xloc)); PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm; hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm); hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf); PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); PetscCall(VecZeroEntries(F)); PetscCall(DMDAVecGetArray(da, Xloc, &x)); PetscCall(VecGetArray(F, &f)); PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); if (ctx->bctype == FVBC_OUTFLOW) { for (i = xs - 2; i < 0; i++) { for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; } for (i = Mx; i < xs + xm + 2; i++) { for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; } } for (i = xs - 1; i < xs + xm + 1; i++) { struct _LimitInfo info; PetscScalar *cjmpL, *cjmpR; if ((i > sm - lsbwidth - 2 && i < sm + 1) || (i > ms - 2 && i < ms + rsbwidth + 1)) { /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds)); /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof)); cjmpL = &ctx->cjmpLR[0]; cjmpR = &ctx->cjmpLR[dof]; for (j = 0; j < dof; j++) { PetscScalar jmpL, jmpR; jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j]; jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j]; for (k = 0; k < dof; k++) { cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL; cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR; } } /* Apply limiter to the left and right characteristic jumps */ info.m = dof; info.hxs = hxs; info.hxm = hxm; info.hxf = hxf; (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope); for (j = 0; j < dof; j++) { PetscScalar tmp = 0; for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k]; slope[i * dof + j] = tmp; } } } for (i = xs; i < xs + xm + 1; i++) { PetscReal maxspeed; PetscScalar *uL, *uR; uL = &ctx->uLR[0]; uR = &ctx->uLR[dof]; if (i == sm - lsbwidth) { for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i < xs + xm) { for (j = 0; j < dof; j++) f[islowbuffer * dof + j] += ctx->flux[j] / hxs; islowbuffer++; } } if (i > sm - lsbwidth && i < sm) { for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i > xs) { for (j = 0; j < dof; j++) f[(islowbuffer - 1) * dof + j] -= ctx->flux[j] / hxs; } if (i < xs + xm) { for (j = 0; j < dof; j++) f[islowbuffer * dof + j] += ctx->flux[j] / hxs; islowbuffer++; } } if (i == sm) { /* interface between the slow region and the medium region */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i > xs) { for (j = 0; j < dof; j++) f[(islowbuffer - 1) * dof + j] -= ctx->flux[j] / hxs; } } if (i == ms) { /* interface between the medium region and the slow region */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i < xs + xm) { for (j = 0; j < dof; j++) f[islowbuffer * dof + j] += ctx->flux[j] / hxs; islowbuffer++; } } if (i > ms && i < ms + rsbwidth) { for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i > xs) { for (j = 0; j < dof; j++) f[(islowbuffer - 1) * dof + j] -= ctx->flux[j] / hxs; } if (i < xs + xm) { for (j = 0; j < dof; j++) f[islowbuffer * dof + j] += ctx->flux[j] / hxs; islowbuffer++; } } if (i == ms + rsbwidth) { for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i > xs) { for (j = 0; j < dof; j++) f[(islowbuffer - 1) * dof + j] -= ctx->flux[j] / hxs; } } } PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); PetscCall(VecRestoreArray(F, &f)); PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope)); PetscCall(DMRestoreLocalVector(da, &Xloc)); PetscFunctionReturn(PETSC_SUCCESS); } /* --------------------------------- Finite Volume Solver for medium components ----------------------------------- */ PetscErrorCode FVRHSFunctionmedium_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) { FVCtx *ctx = (FVCtx *)vctx; PetscInt i, j, k, Mx, dof, xs, xm, imedium = 0, sm = ctx->sm, mf = ctx->mf, fm = ctx->fm, ms = ctx->ms, lmbwidth = ctx->lmbwidth, rmbwidth = ctx->rmbwidth; PetscReal hxs, hxm, hxf; PetscScalar *x, *f, *slope; Vec Xloc; DM da; PetscFunctionBeginUser; PetscCall(TSGetDM(ts, &da)); PetscCall(DMGetLocalVector(da, &Xloc)); PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm; hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm); hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf); PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); PetscCall(VecZeroEntries(F)); PetscCall(DMDAVecGetArray(da, Xloc, &x)); PetscCall(VecGetArray(F, &f)); PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); if (ctx->bctype == FVBC_OUTFLOW) { for (i = xs - 2; i < 0; i++) { for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; } for (i = Mx; i < xs + xm + 2; i++) { for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; } } for (i = xs - 1; i < xs + xm + 1; i++) { struct _LimitInfo info; PetscScalar *cjmpL, *cjmpR; if ((i > sm - 2 && i < mf - lmbwidth + 1) || (i > fm + rmbwidth - 2 && i < ms + 1)) { /* slow components and the first and last fast components */ /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds)); /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof)); cjmpL = &ctx->cjmpLR[0]; cjmpR = &ctx->cjmpLR[dof]; for (j = 0; j < dof; j++) { PetscScalar jmpL, jmpR; jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j]; jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j]; for (k = 0; k < dof; k++) { cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL; cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR; } } /* Apply limiter to the left and right characteristic jumps */ info.m = dof; info.hxs = hxs; info.hxm = hxm; info.hxf = hxf; (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope); for (j = 0; j < dof; j++) { PetscScalar tmp = 0; for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k]; slope[i * dof + j] = tmp; } } } for (i = xs; i < xs + xm + 1; i++) { PetscReal maxspeed; PetscScalar *uL, *uR; uL = &ctx->uLR[0]; uR = &ctx->uLR[dof]; if (i == sm) { /* interface between slow and medium regions */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i < xs + xm) { for (j = 0; j < dof; j++) f[imedium * dof + j] += ctx->flux[j] / hxm; imedium++; } } if (i > sm && i < mf - lmbwidth) { /* medium region */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i > xs) { for (j = 0; j < dof; j++) f[(imedium - 1) * dof + j] -= ctx->flux[j] / hxm; } if (i < xs + xm) { for (j = 0; j < dof; j++) f[imedium * dof + j] += ctx->flux[j] / hxm; imedium++; } } if (i == mf - lmbwidth) { /* interface between medium and fast regions */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i > xs) { for (j = 0; j < dof; j++) f[(imedium - 1) * dof + j] -= ctx->flux[j] / hxm; } } if (i == fm + rmbwidth) { /* interface between fast and medium regions */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i < xs + xm) { for (j = 0; j < dof; j++) f[imedium * dof + j] += ctx->flux[j] / hxm; imedium++; } } if (i > fm + rmbwidth && i < ms) { /* medium region */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i > xs) { for (j = 0; j < dof; j++) f[(imedium - 1) * dof + j] -= ctx->flux[j] / hxm; } if (i < xs + xm) { for (j = 0; j < dof; j++) f[imedium * dof + j] += ctx->flux[j] / hxm; imedium++; } } if (i == ms) { /* interface between medium and slow regions */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i > xs) { for (j = 0; j < dof; j++) f[(imedium - 1) * dof + j] -= ctx->flux[j] / hxm; } } } PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); PetscCall(VecRestoreArray(F, &f)); PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope)); PetscCall(DMRestoreLocalVector(da, &Xloc)); PetscFunctionReturn(PETSC_SUCCESS); } PetscErrorCode FVRHSFunctionmediumbuffer_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) { FVCtx *ctx = (FVCtx *)vctx; PetscInt i, j, k, Mx, dof, xs, xm, imediumbuffer = 0, mf = ctx->mf, fm = ctx->fm, lmbwidth = ctx->lmbwidth, rmbwidth = ctx->rmbwidth; PetscReal hxs, hxm, hxf; PetscScalar *x, *f, *slope; Vec Xloc; DM da; PetscFunctionBeginUser; PetscCall(TSGetDM(ts, &da)); PetscCall(DMGetLocalVector(da, &Xloc)); PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm; hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm); hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf); PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); PetscCall(VecZeroEntries(F)); PetscCall(DMDAVecGetArray(da, Xloc, &x)); PetscCall(VecGetArray(F, &f)); PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); if (ctx->bctype == FVBC_OUTFLOW) { for (i = xs - 2; i < 0; i++) { for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; } for (i = Mx; i < xs + xm + 2; i++) { for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; } } for (i = xs - 1; i < xs + xm + 1; i++) { struct _LimitInfo info; PetscScalar *cjmpL, *cjmpR; if ((i > mf - lmbwidth - 2 && i < mf + 1) || (i > fm - 2 && i < fm + rmbwidth + 1)) { /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds)); /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof)); cjmpL = &ctx->cjmpLR[0]; cjmpR = &ctx->cjmpLR[dof]; for (j = 0; j < dof; j++) { PetscScalar jmpL, jmpR; jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j]; jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j]; for (k = 0; k < dof; k++) { cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL; cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR; } } /* Apply limiter to the left and right characteristic jumps */ info.m = dof; info.hxs = hxs; info.hxm = hxm; info.hxf = hxf; (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope); for (j = 0; j < dof; j++) { PetscScalar tmp = 0; for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k]; slope[i * dof + j] = tmp; } } } for (i = xs; i < xs + xm + 1; i++) { PetscReal maxspeed; PetscScalar *uL, *uR; uL = &ctx->uLR[0]; uR = &ctx->uLR[dof]; if (i == mf - lmbwidth) { for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i < xs + xm) { for (j = 0; j < dof; j++) f[imediumbuffer * dof + j] += ctx->flux[j] / hxm; imediumbuffer++; } } if (i > mf - lmbwidth && i < mf) { for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i > xs) { for (j = 0; j < dof; j++) f[(imediumbuffer - 1) * dof + j] -= ctx->flux[j] / hxm; } if (i < xs + xm) { for (j = 0; j < dof; j++) f[imediumbuffer * dof + j] += ctx->flux[j] / hxm; imediumbuffer++; } } if (i == mf) { /* interface between the medium region and the fast region */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i > xs) { for (j = 0; j < dof; j++) f[(imediumbuffer - 1) * dof + j] -= ctx->flux[j] / hxm; } } if (i == fm) { /* interface between the fast region and the medium region */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i < xs + xm) { for (j = 0; j < dof; j++) f[imediumbuffer * dof + j] += ctx->flux[j] / hxm; imediumbuffer++; } } if (i > fm && i < fm + rmbwidth) { for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i > xs) { for (j = 0; j < dof; j++) f[(imediumbuffer - 1) * dof + j] -= ctx->flux[j] / hxm; } if (i < xs + xm) { for (j = 0; j < dof; j++) f[imediumbuffer * dof + j] += ctx->flux[j] / hxm; imediumbuffer++; } } if (i == fm + rmbwidth) { for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i > xs) { for (j = 0; j < dof; j++) f[(imediumbuffer - 1) * dof + j] -= ctx->flux[j] / hxm; } } } PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); PetscCall(VecRestoreArray(F, &f)); PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope)); PetscCall(DMRestoreLocalVector(da, &Xloc)); PetscFunctionReturn(PETSC_SUCCESS); } /* --------------------------------- Finite Volume Solver for fast parts ----------------------------------- */ PetscErrorCode FVRHSFunctionfast_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) { FVCtx *ctx = (FVCtx *)vctx; PetscInt i, j, k, Mx, dof, xs, xm, ifast = 0, mf = ctx->mf, fm = ctx->fm; PetscReal hxs, hxm, hxf; PetscScalar *x, *f, *slope; Vec Xloc; DM da; PetscFunctionBeginUser; PetscCall(TSGetDM(ts, &da)); PetscCall(DMGetLocalVector(da, &Xloc)); PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm; hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm); hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf); PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); PetscCall(VecZeroEntries(F)); PetscCall(DMDAVecGetArray(da, Xloc, &x)); PetscCall(VecGetArray(F, &f)); PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); if (ctx->bctype == FVBC_OUTFLOW) { for (i = xs - 2; i < 0; i++) { for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; } for (i = Mx; i < xs + xm + 2; i++) { for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; } } for (i = xs - 1; i < xs + xm + 1; i++) { /* fast components and the last slow components before fast components and the first slow component after fast components */ struct _LimitInfo info; PetscScalar *cjmpL, *cjmpR; if (i > mf - 2 && i < fm + 1) { PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds)); PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof)); cjmpL = &ctx->cjmpLR[0]; cjmpR = &ctx->cjmpLR[dof]; for (j = 0; j < dof; j++) { PetscScalar jmpL, jmpR; jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j]; jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j]; for (k = 0; k < dof; k++) { cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL; cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR; } } /* Apply limiter to the left and right characteristic jumps */ info.m = dof; info.hxs = hxs; info.hxm = hxm; info.hxf = hxf; (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope); for (j = 0; j < dof; j++) { PetscScalar tmp = 0; for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k]; slope[i * dof + j] = tmp; } } } for (i = xs; i < xs + xm + 1; i++) { PetscReal maxspeed; PetscScalar *uL, *uR; uL = &ctx->uLR[0]; uR = &ctx->uLR[dof]; if (i == mf) { /* interface between medium and fast regions */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i < xs + xm) { for (j = 0; j < dof; j++) f[ifast * dof + j] += ctx->flux[j] / hxf; ifast++; } } if (i > mf && i < fm) { /* fast region */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i > xs) { for (j = 0; j < dof; j++) f[(ifast - 1) * dof + j] -= ctx->flux[j] / hxf; } if (i < xs + xm) { for (j = 0; j < dof; j++) f[ifast * dof + j] += ctx->flux[j] / hxf; ifast++; } } if (i == fm) { /* interface between fast and medium regions */ for (j = 0; j < dof; j++) { uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2; uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; } PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); if (i > xs) { for (j = 0; j < dof; j++) f[(ifast - 1) * dof + j] -= ctx->flux[j] / hxf; } } } PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); PetscCall(VecRestoreArray(F, &f)); PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope)); PetscCall(DMRestoreLocalVector(da, &Xloc)); PetscFunctionReturn(PETSC_SUCCESS); } int main(int argc, char *argv[]) { char lname[256] = "mc", physname[256] = "advect", final_fname[256] = "solution.m"; PetscFunctionList limiters = 0, physics = 0; MPI_Comm comm; TS ts; DM da; Vec X, X0, R; FVCtx ctx; PetscInt i, k, dof, xs, xm, Mx, draw = 0, count_slow, count_medium, count_fast, islow = 0, imedium = 0, ifast = 0, *index_slow, *index_medium, *index_fast, islowbuffer = 0, *index_slowbuffer, imediumbuffer = 0, *index_mediumbuffer; PetscBool view_final = PETSC_FALSE; PetscReal ptime; PetscFunctionBeginUser; PetscCall(PetscInitialize(&argc, &argv, 0, help)); comm = PETSC_COMM_WORLD; PetscCall(PetscMemzero(&ctx, sizeof(ctx))); /* Register limiters to be available on the command line */ PetscCall(PetscFunctionListAdd(&limiters, "upwind", Limit3_Upwind)); PetscCall(PetscFunctionListAdd(&limiters, "lax-wendroff", Limit3_LaxWendroff)); PetscCall(PetscFunctionListAdd(&limiters, "beam-warming", Limit3_BeamWarming)); PetscCall(PetscFunctionListAdd(&limiters, "fromm", Limit3_Fromm)); PetscCall(PetscFunctionListAdd(&limiters, "minmod", Limit3_Minmod)); PetscCall(PetscFunctionListAdd(&limiters, "superbee", Limit3_Superbee)); PetscCall(PetscFunctionListAdd(&limiters, "mc", Limit3_MC)); PetscCall(PetscFunctionListAdd(&limiters, "koren3", Limit3_Koren3)); /* Register physical models to be available on the command line */ PetscCall(PetscFunctionListAdd(&physics, "advect", PhysicsCreate_Advect)); ctx.comm = comm; ctx.cfl = 0.9; ctx.bctype = FVBC_PERIODIC; ctx.xmin = -1.0; ctx.xmax = 1.0; PetscOptionsBegin(comm, NULL, "Finite Volume solver options", ""); PetscCall(PetscOptionsReal("-xmin", "X min", "", ctx.xmin, &ctx.xmin, NULL)); PetscCall(PetscOptionsReal("-xmax", "X max", "", ctx.xmax, &ctx.xmax, NULL)); PetscCall(PetscOptionsFList("-limit", "Name of flux imiter to use", "", limiters, lname, lname, sizeof(lname), NULL)); PetscCall(PetscOptionsInt("-draw", "Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)", "", draw, &draw, NULL)); PetscCall(PetscOptionsString("-view_final", "Write final solution in ASCII MATLAB format to given file name", "", final_fname, final_fname, sizeof(final_fname), &view_final)); PetscCall(PetscOptionsInt("-initial", "Initial condition (depends on the physics)", "", ctx.initial, &ctx.initial, NULL)); PetscCall(PetscOptionsBool("-exact", "Compare errors with exact solution", "", ctx.exact, &ctx.exact, NULL)); PetscCall(PetscOptionsBool("-simulation", "Compare errors with reference solution", "", ctx.simulation, &ctx.simulation, NULL)); PetscCall(PetscOptionsReal("-cfl", "CFL number to time step at", "", ctx.cfl, &ctx.cfl, NULL)); PetscCall(PetscOptionsEnum("-bc_type", "Boundary condition", "", FVBCTypes, (PetscEnum)ctx.bctype, (PetscEnum *)&ctx.bctype, NULL)); PetscCall(PetscOptionsInt("-hratio", "Spacing ratio", "", ctx.hratio, &ctx.hratio, NULL)); PetscOptionsEnd(); /* Choose the limiter from the list of registered limiters */ PetscCall(PetscFunctionListFind(limiters, lname, &ctx.limit3)); PetscCheck(ctx.limit3, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Limiter '%s' not found", lname); /* Choose the physics from the list of registered models */ { PetscErrorCode (*r)(FVCtx *); PetscCall(PetscFunctionListFind(physics, physname, &r)); PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Physics '%s' not found", physname); /* Create the physics, will set the number of fields and their names */ PetscCall((*r)(&ctx)); } /* Create a DMDA to manage the parallel grid */ PetscCall(DMDACreate1d(comm, DM_BOUNDARY_PERIODIC, 50, ctx.physics2.dof, 2, NULL, &da)); PetscCall(DMSetFromOptions(da)); PetscCall(DMSetUp(da)); /* Inform the DMDA of the field names provided by the physics. */ /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */ for (i = 0; i < ctx.physics2.dof; i++) PetscCall(DMDASetFieldName(da, i, ctx.physics2.fieldname[i])); PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); /* Set coordinates of cell centers */ 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)); /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */ PetscCall(PetscMalloc4(dof * dof, &ctx.R, dof * dof, &ctx.Rinv, 2 * dof, &ctx.cjmpLR, 1 * dof, &ctx.cslope)); PetscCall(PetscMalloc3(2 * dof, &ctx.uLR, dof, &ctx.flux, dof, &ctx.speeds)); /* Create a vector to store the solution and to save the initial state */ PetscCall(DMCreateGlobalVector(da, &X)); PetscCall(VecDuplicate(X, &X0)); PetscCall(VecDuplicate(X, &R)); /* create index for slow parts and fast parts, count_slow + count_fast = Mx, counts_slow*hs = 0.5, counts_fast*hf = 0.5 */ count_slow = Mx / (1 + ctx.hratio) / (1 + ctx.hratio); count_medium = 2 * ctx.hratio * count_slow; PetscCheck((count_slow % 2) == 0 && (count_medium % 2) == 0, PETSC_COMM_WORLD, PETSC_ERR_USER, "Please adjust grid size Mx (-da_grid_x) and hratio (-hratio) so that Mx/(1+hartio)^2 and Mx*2*hratio/(1+hratio)^2 is even"); count_fast = ctx.hratio * ctx.hratio * count_slow; ctx.sm = count_slow / 2; ctx.mf = ctx.sm + count_medium / 2; ctx.fm = ctx.mf + count_fast; ctx.ms = ctx.fm + count_medium / 2; PetscCall(PetscMalloc1(xm * dof, &index_slow)); PetscCall(PetscMalloc1(xm * dof, &index_medium)); PetscCall(PetscMalloc1(xm * dof, &index_fast)); PetscCall(PetscMalloc1(6 * dof, &index_slowbuffer)); PetscCall(PetscMalloc1(6 * dof, &index_mediumbuffer)); if (((AdvectCtx *)ctx.physics2.user)->a > 0) { ctx.lsbwidth = 2; ctx.rsbwidth = 4; ctx.lmbwidth = 2; ctx.rmbwidth = 4; } else { ctx.lsbwidth = 4; ctx.rsbwidth = 2; ctx.lmbwidth = 4; ctx.rmbwidth = 2; } for (i = xs; i < xs + xm; i++) { if (i < ctx.sm - ctx.lsbwidth || i > ctx.ms + ctx.rsbwidth - 1) for (k = 0; k < dof; k++) index_slow[islow++] = i * dof + k; else if ((i >= ctx.sm - ctx.lsbwidth && i < ctx.sm) || (i > ctx.ms - 1 && i <= ctx.ms + ctx.rsbwidth - 1)) for (k = 0; k < dof; k++) index_slowbuffer[islowbuffer++] = i * dof + k; else if (i < ctx.mf - ctx.lmbwidth || i > ctx.fm + ctx.rmbwidth - 1) for (k = 0; k < dof; k++) index_medium[imedium++] = i * dof + k; else if ((i >= ctx.mf - ctx.lmbwidth && i < ctx.mf) || (i > ctx.fm - 1 && i <= ctx.fm + ctx.rmbwidth - 1)) for (k = 0; k < dof; k++) index_mediumbuffer[imediumbuffer++] = i * dof + k; else for (k = 0; k < dof; k++) index_fast[ifast++] = i * dof + k; } PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islow, index_slow, PETSC_COPY_VALUES, &ctx.iss)); PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, imedium, index_medium, PETSC_COPY_VALUES, &ctx.ism)); PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, ifast, index_fast, PETSC_COPY_VALUES, &ctx.isf)); PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islowbuffer, index_slowbuffer, PETSC_COPY_VALUES, &ctx.issb)); PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, imediumbuffer, index_mediumbuffer, PETSC_COPY_VALUES, &ctx.ismb)); /* Create a time-stepping object */ PetscCall(TSCreate(comm, &ts)); PetscCall(TSSetDM(ts, da)); PetscCall(TSSetRHSFunction(ts, R, FVRHSFunction_3WaySplit, &ctx)); PetscCall(TSRHSSplitSetIS(ts, "slow", ctx.iss)); PetscCall(TSRHSSplitSetIS(ts, "medium", ctx.ism)); PetscCall(TSRHSSplitSetIS(ts, "fast", ctx.isf)); PetscCall(TSRHSSplitSetIS(ts, "slowbuffer", ctx.issb)); PetscCall(TSRHSSplitSetIS(ts, "mediumbuffer", ctx.ismb)); PetscCall(TSRHSSplitSetRHSFunction(ts, "slow", NULL, FVRHSFunctionslow_3WaySplit, &ctx)); PetscCall(TSRHSSplitSetRHSFunction(ts, "medium", NULL, FVRHSFunctionmedium_3WaySplit, &ctx)); PetscCall(TSRHSSplitSetRHSFunction(ts, "fast", NULL, FVRHSFunctionfast_3WaySplit, &ctx)); PetscCall(TSRHSSplitSetRHSFunction(ts, "slowbuffer", NULL, FVRHSFunctionslowbuffer_3WaySplit, &ctx)); PetscCall(TSRHSSplitSetRHSFunction(ts, "mediumbuffer", NULL, FVRHSFunctionmediumbuffer_3WaySplit, &ctx)); PetscCall(TSSetType(ts, TSSSP)); /*PetscCall(TSSetType(ts,TSMPRK));*/ PetscCall(TSSetMaxTime(ts, 10)); PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); /* Compute initial conditions and starting time step */ PetscCall(FVSample_3WaySplit(&ctx, da, 0, X0)); PetscCall(FVRHSFunction_3WaySplit(ts, 0, X0, X, (void *)&ctx)); /* Initial function evaluation, only used to determine max speed */ PetscCall(VecCopy(X0, X)); /* The function value was not used so we set X=X0 again */ PetscCall(TSSetTimeStep(ts, ctx.cfl / ctx.cfl_idt)); PetscCall(TSSetFromOptions(ts)); /* Take runtime options */ PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD)); { PetscInt steps; PetscScalar mass_initial, mass_final, mass_difference, mass_differenceg; const PetscScalar *ptr_X, *ptr_X0; const PetscReal hs = (ctx.xmax - ctx.xmin) / 4.0 / count_slow; const PetscReal hm = (ctx.xmax - ctx.xmin) / 2.0 / count_medium; const PetscReal hf = (ctx.xmax - ctx.xmin) / 4.0 / count_fast; PetscCall(TSSolve(ts, X)); PetscCall(TSGetSolveTime(ts, &ptime)); PetscCall(TSGetStepNumber(ts, &steps)); /* calculate the total mass at initial time and final time */ mass_initial = 0.0; mass_final = 0.0; PetscCall(DMDAVecGetArrayRead(da, X0, (void *)&ptr_X0)); PetscCall(DMDAVecGetArrayRead(da, X, (void *)&ptr_X)); for (i = xs; i < xs + xm; i++) { if (i < ctx.sm || i > ctx.ms - 1) for (k = 0; k < dof; k++) { mass_initial = mass_initial + hs * ptr_X0[i * dof + k]; mass_final = mass_final + hs * ptr_X[i * dof + k]; } else if (i < ctx.mf || i > ctx.fm - 1) for (k = 0; k < dof; k++) { mass_initial = mass_initial + hm * ptr_X0[i * dof + k]; mass_final = mass_final + hm * ptr_X[i * dof + k]; } else { for (k = 0; k < dof; k++) { mass_initial = mass_initial + hf * ptr_X0[i * dof + k]; mass_final = mass_final + hf * ptr_X[i * dof + k]; } } } PetscCall(DMDAVecRestoreArrayRead(da, X0, (void *)&ptr_X0)); PetscCall(DMDAVecRestoreArrayRead(da, X, (void *)&ptr_X)); mass_difference = mass_final - mass_initial; PetscCall(MPIU_Allreduce(&mass_difference, &mass_differenceg, 1, MPIU_SCALAR, MPIU_SUM, comm)); PetscCall(PetscPrintf(comm, "Mass difference %g\n", (double)mass_differenceg)); PetscCall(PetscPrintf(comm, "Final time %g, steps %" PetscInt_FMT "\n", (double)ptime, steps)); PetscCall(PetscPrintf(comm, "Maximum allowable stepsize according to CFL %g\n", (double)(1.0 / ctx.cfl_idt))); if (ctx.exact) { PetscReal nrm1 = 0; PetscCall(SolutionErrorNorms_3WaySplit(&ctx, da, ptime, X, &nrm1)); PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %g\n", (double)nrm1)); } if (ctx.simulation) { PetscReal nrm1 = 0; PetscViewer fd; char filename[PETSC_MAX_PATH_LEN] = "binaryoutput"; Vec XR; PetscBool flg; const PetscScalar *ptr_XR; PetscCall(PetscOptionsGetString(NULL, NULL, "-f", filename, sizeof(filename), &flg)); PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate binary file with the -f option"); PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, filename, FILE_MODE_READ, &fd)); PetscCall(VecDuplicate(X0, &XR)); PetscCall(VecLoad(XR, fd)); PetscCall(PetscViewerDestroy(&fd)); PetscCall(VecGetArrayRead(X, &ptr_X)); PetscCall(VecGetArrayRead(XR, &ptr_XR)); for (i = xs; i < xs + xm; i++) { if (i < ctx.sm || i > ctx.ms - 1) for (k = 0; k < dof; k++) nrm1 = nrm1 + hs * PetscAbs(ptr_X[i * dof + k] - ptr_XR[i * dof + k]); else if (i < ctx.mf || i < ctx.fm - 1) for (k = 0; k < dof; k++) nrm1 = nrm1 + hm * PetscAbs(ptr_X[i * dof + k] - ptr_XR[i * dof + k]); else for (k = 0; k < dof; k++) nrm1 = nrm1 + hf * PetscAbs(ptr_X[i * dof + k] - ptr_XR[i * dof + k]); } PetscCall(VecRestoreArrayRead(X, &ptr_X)); PetscCall(VecRestoreArrayRead(XR, &ptr_XR)); PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %g\n", (double)nrm1)); PetscCall(VecDestroy(&XR)); } } PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD)); if (draw & 0x1) PetscCall(VecView(X0, PETSC_VIEWER_DRAW_WORLD)); if (draw & 0x2) PetscCall(VecView(X, PETSC_VIEWER_DRAW_WORLD)); if (draw & 0x4) { Vec Y; PetscCall(VecDuplicate(X, &Y)); PetscCall(FVSample_3WaySplit(&ctx, da, ptime, Y)); PetscCall(VecAYPX(Y, -1, X)); PetscCall(VecView(Y, PETSC_VIEWER_DRAW_WORLD)); PetscCall(VecDestroy(&Y)); } if (view_final) { PetscViewer viewer; PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, final_fname, &viewer)); PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB)); PetscCall(VecView(X, viewer)); PetscCall(PetscViewerPopFormat(viewer)); PetscCall(PetscViewerDestroy(&viewer)); } /* Clean up */ PetscCall((*ctx.physics2.destroy)(ctx.physics2.user)); for (i = 0; i < ctx.physics2.dof; i++) PetscCall(PetscFree(ctx.physics2.fieldname[i])); PetscCall(PetscFree4(ctx.R, ctx.Rinv, ctx.cjmpLR, ctx.cslope)); PetscCall(PetscFree3(ctx.uLR, ctx.flux, ctx.speeds)); PetscCall(VecDestroy(&X)); PetscCall(VecDestroy(&X0)); PetscCall(VecDestroy(&R)); PetscCall(DMDestroy(&da)); PetscCall(TSDestroy(&ts)); PetscCall(ISDestroy(&ctx.iss)); PetscCall(ISDestroy(&ctx.ism)); PetscCall(ISDestroy(&ctx.isf)); PetscCall(ISDestroy(&ctx.issb)); PetscCall(ISDestroy(&ctx.ismb)); PetscCall(PetscFree(index_slow)); PetscCall(PetscFree(index_medium)); PetscCall(PetscFree(index_fast)); PetscCall(PetscFree(index_slowbuffer)); PetscCall(PetscFree(index_mediumbuffer)); PetscCall(PetscFunctionListDestroy(&limiters)); PetscCall(PetscFunctionListDestroy(&physics)); PetscCall(PetscFinalize()); return 0; } /*TEST build: requires: !complex depends: finitevolume1d.c test: suffix: 1 args: -da_grid_x 90 -initial 7 -xmin -1 -xmax 1 -hratio 2 -limit mc -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a23 -ts_use_splitrhsfunction 0 test: suffix: 2 args: -da_grid_x 90 -initial 7 -xmin -1 -xmax 1 -hratio 2 -limit mc -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a23 -ts_use_splitrhsfunction 1 TEST*/