1c4762a1bSJed Brown static const char help[] = "1D periodic Finite Volume solver in slope-limiter form 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 (slow-medium-fast-medium-slow), \n" 5c4762a1bSJed Brown " the meshsize ratio between two adjacient sub-domains is controlled with -hratio,\n" 6c4762a1bSJed Brown " exact - Exact Riemann solver which usually needs to perform a Newton iteration to connect\n" 7c4762a1bSJed Brown " the states across shocks and rarefactions\n" 8c4762a1bSJed Brown " simulation - use reference solution which is generated by smaller time step size to be true solution,\n" 9c4762a1bSJed Brown " also the reference solution should be generated by user and stored in a binary file.\n" 10c4762a1bSJed Brown " characteristic - Limit the characteristic variables, this is usually preferred (default)\n" 11c4762a1bSJed Brown "Several initial conditions can be chosen with -initial N\n\n" 12c4762a1bSJed Brown "The problem size should be set with -da_grid_x M\n\n"; 13c4762a1bSJed Brown 14c4762a1bSJed Brown #include <petscts.h> 15c4762a1bSJed Brown #include <petscdm.h> 16c4762a1bSJed Brown #include <petscdmda.h> 17c4762a1bSJed Brown #include <petscdraw.h> 18c4762a1bSJed Brown #include "finitevolume1d.h" 19c4762a1bSJed Brown 209371c9d4SSatish Balay static inline PetscReal RangeMod(PetscReal a, PetscReal xmin, PetscReal xmax) { 219371c9d4SSatish Balay PetscReal range = xmax - xmin; 229371c9d4SSatish Balay return xmin + PetscFmodReal(range + PetscFmodReal(a, range), range); 239371c9d4SSatish Balay } 24c4762a1bSJed Brown 25c4762a1bSJed Brown /* --------------------------------- Advection ----------------------------------- */ 26c4762a1bSJed Brown typedef struct { 27c4762a1bSJed Brown PetscReal a; /* advective velocity */ 28c4762a1bSJed Brown } AdvectCtx; 29c4762a1bSJed Brown 309371c9d4SSatish Balay static PetscErrorCode PhysicsRiemann_Advect(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed) { 31c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx *)vctx; 32c4762a1bSJed Brown PetscReal speed; 33c4762a1bSJed Brown 34c4762a1bSJed Brown PetscFunctionBeginUser; 35c4762a1bSJed Brown speed = ctx->a; 36c4762a1bSJed Brown flux[0] = PetscMax(0, speed) * uL[0] + PetscMin(0, speed) * uR[0]; 37c4762a1bSJed Brown *maxspeed = speed; 38c4762a1bSJed Brown PetscFunctionReturn(0); 39c4762a1bSJed Brown } 40c4762a1bSJed Brown 419371c9d4SSatish Balay static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx, PetscInt m, const PetscScalar *u, PetscScalar *X, PetscScalar *Xi, PetscReal *speeds) { 42c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx *)vctx; 43c4762a1bSJed Brown 44c4762a1bSJed Brown PetscFunctionBeginUser; 45c4762a1bSJed Brown X[0] = 1.; 46c4762a1bSJed Brown Xi[0] = 1.; 47c4762a1bSJed Brown speeds[0] = ctx->a; 48c4762a1bSJed Brown PetscFunctionReturn(0); 49c4762a1bSJed Brown } 50c4762a1bSJed Brown 519371c9d4SSatish Balay static PetscErrorCode PhysicsSample_Advect(void *vctx, PetscInt initial, FVBCType bctype, PetscReal xmin, PetscReal xmax, PetscReal t, PetscReal x, PetscReal *u) { 52c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx *)vctx; 53c4762a1bSJed Brown PetscReal a = ctx->a, x0; 54c4762a1bSJed Brown 55c4762a1bSJed Brown PetscFunctionBeginUser; 56c4762a1bSJed Brown switch (bctype) { 57c4762a1bSJed Brown case FVBC_OUTFLOW: x0 = x - a * t; break; 58c4762a1bSJed Brown case FVBC_PERIODIC: x0 = RangeMod(x - a * t, xmin, xmax); break; 59c4762a1bSJed Brown default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown BCType"); 60c4762a1bSJed Brown } 61c4762a1bSJed Brown switch (initial) { 62c4762a1bSJed Brown case 0: u[0] = (x0 < 0) ? 1 : -1; break; 63c4762a1bSJed Brown case 1: u[0] = (x0 < 0) ? -1 : 1; break; 64c4762a1bSJed Brown case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break; 65c4762a1bSJed Brown case 3: u[0] = PetscSinReal(2 * PETSC_PI * x0); break; 66c4762a1bSJed Brown case 4: u[0] = PetscAbs(x0); break; 67c4762a1bSJed Brown case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2 * PETSC_PI * x0)); break; 68c4762a1bSJed Brown case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2 - x0 : 0)); break; 69c4762a1bSJed Brown case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI * x0), 10.0); break; 70c4762a1bSJed Brown default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition"); 71c4762a1bSJed Brown } 72c4762a1bSJed Brown PetscFunctionReturn(0); 73c4762a1bSJed Brown } 74c4762a1bSJed Brown 759371c9d4SSatish Balay static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx) { 76c4762a1bSJed Brown AdvectCtx *user; 77c4762a1bSJed Brown 78c4762a1bSJed Brown PetscFunctionBeginUser; 799566063dSJacob Faibussowitsch PetscCall(PetscNew(&user)); 80c4762a1bSJed Brown ctx->physics2.sample2 = PhysicsSample_Advect; 81c4762a1bSJed Brown ctx->physics2.riemann2 = PhysicsRiemann_Advect; 82c4762a1bSJed Brown ctx->physics2.characteristic2 = PhysicsCharacteristic_Advect; 83c4762a1bSJed Brown ctx->physics2.destroy = PhysicsDestroy_SimpleFree; 84c4762a1bSJed Brown ctx->physics2.user = user; 85c4762a1bSJed Brown ctx->physics2.dof = 1; 869566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy("u", &ctx->physics2.fieldname[0])); 87c4762a1bSJed Brown user->a = 1; 88d0609cedSBarry Smith PetscOptionsBegin(ctx->comm, ctx->prefix, "Options for advection", ""); 899371c9d4SSatish Balay { PetscCall(PetscOptionsReal("-physics_advect_a", "Speed", "", user->a, &user->a, NULL)); } 90d0609cedSBarry Smith PetscOptionsEnd(); 91c4762a1bSJed Brown PetscFunctionReturn(0); 92c4762a1bSJed Brown } 93c4762a1bSJed Brown 949371c9d4SSatish Balay PetscErrorCode FVSample_3WaySplit(FVCtx *ctx, DM da, PetscReal time, Vec U) { 95c4762a1bSJed Brown PetscScalar *u, *uj, xj, xi; 96c4762a1bSJed Brown PetscInt i, j, k, dof, xs, xm, Mx; 97c4762a1bSJed Brown const PetscInt N = 200; 98c4762a1bSJed Brown PetscReal hs, hm, hf; 99c4762a1bSJed Brown 100c4762a1bSJed Brown PetscFunctionBeginUser; 1013c633725SBarry Smith PetscCheck(ctx->physics2.sample2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Physics has not provided a sampling function"); 1029566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); 1039566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 1049566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u)); 1059566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof, &uj)); 106c4762a1bSJed Brown 107c4762a1bSJed Brown hs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm; 108c4762a1bSJed Brown hm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm); 109c4762a1bSJed Brown hf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf); 110c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 111c4762a1bSJed Brown if (i < ctx->sm) { 112c4762a1bSJed Brown xi = ctx->xmin + 0.5 * hs + i * hs; 113c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 114c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] = 0; 115c4762a1bSJed Brown for (j = 0; j < N + 1; j++) { 116c4762a1bSJed Brown xj = xi + hs * (j - N / 2) / (PetscReal)N; 1179566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj)); 118c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N; 119c4762a1bSJed Brown } 120c4762a1bSJed Brown } else if (i < ctx->mf) { 121c4762a1bSJed Brown xi = ctx->xmin + ctx->sm * hs + 0.5 * hm + (i - ctx->sm) * hm; 122c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 123c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] = 0; 124c4762a1bSJed Brown for (j = 0; j < N + 1; j++) { 125c4762a1bSJed Brown xj = xi + hm * (j - N / 2) / (PetscReal)N; 1269566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj)); 127c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N; 128c4762a1bSJed Brown } 129c4762a1bSJed Brown } else if (i < ctx->fm) { 130c4762a1bSJed Brown xi = ctx->xmin + ctx->sm * hs + (ctx->mf - ctx->sm) * hm + 0.5 * hf + (i - ctx->mf) * hf; 131c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 132c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] = 0; 133c4762a1bSJed Brown for (j = 0; j < N + 1; j++) { 134c4762a1bSJed Brown xj = xi + hf * (j - N / 2) / (PetscReal)N; 1359566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj)); 136c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N; 137c4762a1bSJed Brown } 138c4762a1bSJed Brown } else if (i < ctx->ms) { 139c4762a1bSJed Brown xi = ctx->xmin + ctx->sm * hs + (ctx->mf - ctx->sm) * hm + (ctx->fm - ctx->mf) * hf + 0.5 * hm + (i - ctx->fm) * hm; 140c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 141c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] = 0; 142c4762a1bSJed Brown for (j = 0; j < N + 1; j++) { 143c4762a1bSJed Brown xj = xi + hm * (j - N / 2) / (PetscReal)N; 1449566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj)); 145c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N; 146c4762a1bSJed Brown } 147c4762a1bSJed Brown } else { 148c4762a1bSJed Brown 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; 149c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 150c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] = 0; 151c4762a1bSJed Brown for (j = 0; j < N + 1; j++) { 152c4762a1bSJed Brown xj = xi + hs * (j - N / 2) / (PetscReal)N; 1539566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj)); 154c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N; 155c4762a1bSJed Brown } 156c4762a1bSJed Brown } 157c4762a1bSJed Brown } 1589566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u)); 1599566063dSJacob Faibussowitsch PetscCall(PetscFree(uj)); 160c4762a1bSJed Brown PetscFunctionReturn(0); 161c4762a1bSJed Brown } 162c4762a1bSJed Brown 1639371c9d4SSatish Balay static PetscErrorCode SolutionErrorNorms_3WaySplit(FVCtx *ctx, DM da, PetscReal t, Vec X, PetscReal *nrm1) { 164c4762a1bSJed Brown Vec Y; 165c4762a1bSJed Brown PetscInt i, Mx; 166c4762a1bSJed Brown const PetscScalar *ptr_X, *ptr_Y; 167c4762a1bSJed Brown PetscReal hs, hm, hf; 168c4762a1bSJed Brown 169c4762a1bSJed Brown PetscFunctionBeginUser; 1709566063dSJacob Faibussowitsch PetscCall(VecGetSize(X, &Mx)); 171c4762a1bSJed Brown hs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm; 172c4762a1bSJed Brown hm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm); 173c4762a1bSJed Brown hf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf); 1749566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &Y)); 1759566063dSJacob Faibussowitsch PetscCall(FVSample_3WaySplit(ctx, da, t, Y)); 1769566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &ptr_X)); 1779566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &ptr_Y)); 178c4762a1bSJed Brown for (i = 0; i < Mx; i++) { 179c4762a1bSJed Brown if (i < ctx->sm || i > ctx->ms - 1) *nrm1 += hs * PetscAbs(ptr_X[i] - ptr_Y[i]); 180c4762a1bSJed Brown else if (i < ctx->mf || i > ctx->fm - 1) *nrm1 += hm * PetscAbs(ptr_X[i] - ptr_Y[i]); 181c4762a1bSJed Brown else *nrm1 += hf * PetscAbs(ptr_X[i] - ptr_Y[i]); 182c4762a1bSJed Brown } 1839566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &ptr_X)); 1849566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &ptr_Y)); 1859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 186c4762a1bSJed Brown PetscFunctionReturn(0); 187c4762a1bSJed Brown } 188c4762a1bSJed Brown 1899371c9d4SSatish Balay PetscErrorCode FVRHSFunction_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) { 190c4762a1bSJed Brown FVCtx *ctx = (FVCtx *)vctx; 191c4762a1bSJed Brown PetscInt i, j, k, Mx, dof, xs, xm, sm = ctx->sm, mf = ctx->mf, fm = ctx->fm, ms = ctx->ms; 192c4762a1bSJed Brown PetscReal hxf, hxm, hxs, cfl_idt = 0; 193c4762a1bSJed Brown PetscScalar *x, *f, *slope; 194c4762a1bSJed Brown Vec Xloc; 195c4762a1bSJed Brown DM da; 196c4762a1bSJed Brown 197c4762a1bSJed Brown PetscFunctionBeginUser; 1989566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 1999566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc)); /* Xloc contains ghost points */ 2009566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); /* Mx is the number of center points */ 201c4762a1bSJed Brown hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm; 202c4762a1bSJed Brown hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm); 203c4762a1bSJed Brown hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf); 2049566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); /* X is solution vector which does not contain ghost points */ 2059566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); 206c4762a1bSJed Brown 2079566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); /* F is the right hand side function corresponds to center points */ 208c4762a1bSJed Brown 2099566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x)); 2109566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f)); 2119566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); /* contains ghost points */ 212c4762a1bSJed Brown 2139566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 214c4762a1bSJed Brown 215c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 216c4762a1bSJed Brown for (i = xs - 2; i < 0; i++) { 217c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; 218c4762a1bSJed Brown } 219c4762a1bSJed Brown for (i = Mx; i < xs + xm + 2; i++) { 220c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; 221c4762a1bSJed Brown } 222c4762a1bSJed Brown } 223c4762a1bSJed Brown for (i = xs - 1; i < xs + xm + 1; i++) { 224c4762a1bSJed Brown struct _LimitInfo info; 225c4762a1bSJed Brown PetscScalar *cjmpL, *cjmpR; 226c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 2279566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds)); 228c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 2299566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof)); 230c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 231c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 232c4762a1bSJed Brown for (j = 0; j < dof; j++) { 233c4762a1bSJed Brown PetscScalar jmpL, jmpR; 234c4762a1bSJed Brown jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j]; 235c4762a1bSJed Brown jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j]; 236c4762a1bSJed Brown for (k = 0; k < dof; k++) { 237c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL; 238c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR; 239c4762a1bSJed Brown } 240c4762a1bSJed Brown } 241c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 242c4762a1bSJed Brown info.m = dof; 243c4762a1bSJed Brown info.hxs = hxs; 244c4762a1bSJed Brown info.hxm = hxm; 245c4762a1bSJed Brown info.hxf = hxf; 246c4762a1bSJed Brown (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope); 247c4762a1bSJed Brown for (j = 0; j < dof; j++) { 248c4762a1bSJed Brown PetscScalar tmp = 0; 249c4762a1bSJed Brown for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k]; 250c4762a1bSJed Brown slope[i * dof + j] = tmp; 251c4762a1bSJed Brown } 252c4762a1bSJed Brown } 253c4762a1bSJed Brown 254c4762a1bSJed Brown for (i = xs; i < xs + xm + 1; i++) { 255c4762a1bSJed Brown PetscReal maxspeed; 256c4762a1bSJed Brown PetscScalar *uL, *uR; 257c4762a1bSJed Brown uL = &ctx->uLR[0]; 258c4762a1bSJed Brown uR = &ctx->uLR[dof]; 259c4762a1bSJed Brown if (i < sm || i > ms) { /* slow region */ 260c4762a1bSJed Brown for (j = 0; j < dof; j++) { 261c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; 262c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; 263c4762a1bSJed Brown } 2649566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 265c4762a1bSJed Brown if (i > xs) { 266c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxs; 267c4762a1bSJed Brown } 268c4762a1bSJed Brown if (i < xs + xm) { 269c4762a1bSJed Brown for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxs; 270c4762a1bSJed Brown } 271c4762a1bSJed Brown } else if (i == sm) { /* interface between slow and medium component */ 272c4762a1bSJed Brown for (j = 0; j < dof; j++) { 273c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; 274c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; 275c4762a1bSJed Brown } 2769566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 277c4762a1bSJed Brown if (i > xs) { 278c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxs; 279c4762a1bSJed Brown } 280c4762a1bSJed Brown if (i < xs + xm) { 281c4762a1bSJed Brown for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxm; 282c4762a1bSJed Brown } 283c4762a1bSJed Brown } else if (i == ms) { /* interface between medium and slow regions */ 284c4762a1bSJed Brown for (j = 0; j < dof; j++) { 285c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; 286c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; 287c4762a1bSJed Brown } 2889566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 289c4762a1bSJed Brown if (i > xs) { 290c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxm; 291c4762a1bSJed Brown } 292c4762a1bSJed Brown if (i < xs + xm) { 293c4762a1bSJed Brown for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxs; 294c4762a1bSJed Brown } 295c4762a1bSJed Brown } else if (i < mf || i > fm) { /* medium region */ 296c4762a1bSJed Brown for (j = 0; j < dof; j++) { 297c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; 298c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; 299c4762a1bSJed Brown } 3009566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 301c4762a1bSJed Brown if (i > xs) { 302c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxm; 303c4762a1bSJed Brown } 304c4762a1bSJed Brown if (i < xs + xm) { 305c4762a1bSJed Brown for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxm; 306c4762a1bSJed Brown } 307c4762a1bSJed Brown } else if (i == mf) { /* interface between medium and fast regions */ 308c4762a1bSJed Brown for (j = 0; j < dof; j++) { 309c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; 310c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2; 311c4762a1bSJed Brown } 3129566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 313c4762a1bSJed Brown if (i > xs) { 314c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxm; 315c4762a1bSJed Brown } 316c4762a1bSJed Brown if (i < xs + xm) { 317c4762a1bSJed Brown for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxf; 318c4762a1bSJed Brown } 319c4762a1bSJed Brown } else if (i == fm) { /* interface between fast and medium regions */ 320c4762a1bSJed Brown for (j = 0; j < dof; j++) { 321c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2; 322c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; 323c4762a1bSJed Brown } 3249566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 325c4762a1bSJed Brown if (i > xs) { 326c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxf; 327c4762a1bSJed Brown } 328c4762a1bSJed Brown if (i < xs + xm) { 329c4762a1bSJed Brown for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxm; 330c4762a1bSJed Brown } 331c4762a1bSJed Brown } else { /* fast region */ 332c4762a1bSJed Brown for (j = 0; j < dof; j++) { 333c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2; 334c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2; 335c4762a1bSJed Brown } 3369566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 337c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hxs)); /* Max allowable value of 1/Delta t */ 338c4762a1bSJed Brown if (i > xs) { 339c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxf; 340c4762a1bSJed Brown } 341c4762a1bSJed Brown if (i < xs + xm) { 342c4762a1bSJed Brown for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxf; 343c4762a1bSJed Brown } 344c4762a1bSJed Brown } 345c4762a1bSJed Brown } 3469566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); 3479566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f)); 3489566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope)); 3499566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc)); 3509566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&cfl_idt, &ctx->cfl_idt, 1, MPIU_SCALAR, MPIU_MAX, PetscObjectComm((PetscObject)da))); 351c4762a1bSJed Brown if (0) { 352c4762a1bSJed Brown /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */ 353c4762a1bSJed Brown PetscReal dt, tnow; 3549566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt)); 3559566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &tnow)); 356*48a46eb9SPierre 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))); 357c4762a1bSJed Brown } 358c4762a1bSJed Brown PetscFunctionReturn(0); 359c4762a1bSJed Brown } 360c4762a1bSJed Brown 361c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */ 3629371c9d4SSatish Balay PetscErrorCode FVRHSFunctionslow_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) { 363c4762a1bSJed Brown FVCtx *ctx = (FVCtx *)vctx; 364c4762a1bSJed Brown PetscInt i, j, k, Mx, dof, xs, xm, islow = 0, sm = ctx->sm, ms = ctx->ms, lsbwidth = ctx->lsbwidth, rsbwidth = ctx->rsbwidth; 365c4762a1bSJed Brown PetscReal hxs, hxm, hxf, cfl_idt = 0; 366c4762a1bSJed Brown PetscScalar *x, *f, *slope; 367c4762a1bSJed Brown Vec Xloc; 368c4762a1bSJed Brown DM da; 369c4762a1bSJed Brown 370c4762a1bSJed Brown PetscFunctionBeginUser; 3719566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 3729566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc)); 3739566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); 374c4762a1bSJed Brown hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm; 375c4762a1bSJed Brown hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm); 376c4762a1bSJed Brown hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf); 3779566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); 3789566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); 3799566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 3809566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x)); 3819566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 3829566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); 3839566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 384c4762a1bSJed Brown 385c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 386c4762a1bSJed Brown for (i = xs - 2; i < 0; i++) { 387c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; 388c4762a1bSJed Brown } 389c4762a1bSJed Brown for (i = Mx; i < xs + xm + 2; i++) { 390c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; 391c4762a1bSJed Brown } 392c4762a1bSJed Brown } 393c4762a1bSJed Brown for (i = xs - 1; i < xs + xm + 1; i++) { 394c4762a1bSJed Brown struct _LimitInfo info; 395c4762a1bSJed Brown PetscScalar *cjmpL, *cjmpR; 396c4762a1bSJed Brown if (i < sm - lsbwidth + 1 || i > ms + rsbwidth - 2) { /* slow components and the first and last fast components */ 397c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 3989566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds)); 399c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 4009566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof)); 401c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 402c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 403c4762a1bSJed Brown for (j = 0; j < dof; j++) { 404c4762a1bSJed Brown PetscScalar jmpL, jmpR; 405c4762a1bSJed Brown jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j]; 406c4762a1bSJed Brown jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j]; 407c4762a1bSJed Brown for (k = 0; k < dof; k++) { 408c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL; 409c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR; 410c4762a1bSJed Brown } 411c4762a1bSJed Brown } 412c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 413c4762a1bSJed Brown info.m = dof; 414c4762a1bSJed Brown info.hxs = hxs; 415c4762a1bSJed Brown info.hxm = hxm; 416c4762a1bSJed Brown info.hxf = hxf; 417c4762a1bSJed Brown (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope); 418c4762a1bSJed Brown for (j = 0; j < dof; j++) { 419c4762a1bSJed Brown PetscScalar tmp = 0; 420c4762a1bSJed Brown for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k]; 421c4762a1bSJed Brown slope[i * dof + j] = tmp; 422c4762a1bSJed Brown } 423c4762a1bSJed Brown } 424c4762a1bSJed Brown } 425c4762a1bSJed Brown 426c4762a1bSJed Brown for (i = xs; i < xs + xm + 1; i++) { 427c4762a1bSJed Brown PetscReal maxspeed; 428c4762a1bSJed Brown PetscScalar *uL, *uR; 429c4762a1bSJed Brown uL = &ctx->uLR[0]; 430c4762a1bSJed Brown uR = &ctx->uLR[dof]; 431c4762a1bSJed Brown if (i < sm - lsbwidth) { /* slow region */ 432c4762a1bSJed Brown for (j = 0; j < dof; j++) { 433c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; 434c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; 435c4762a1bSJed Brown } 4369566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 437c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hxs)); /* Max allowable value of 1/Delta t */ 438c4762a1bSJed Brown if (i > xs) { 439c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs; 440c4762a1bSJed Brown } 441c4762a1bSJed Brown if (i < xs + xm) { 442c4762a1bSJed Brown for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs; 443c4762a1bSJed Brown islow++; 444c4762a1bSJed Brown } 445c4762a1bSJed Brown } 446c4762a1bSJed Brown if (i == sm - lsbwidth) { /* interface between slow and medium regions */ 447c4762a1bSJed Brown for (j = 0; j < dof; j++) { 448c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; 449c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; 450c4762a1bSJed Brown } 4519566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 452c4762a1bSJed Brown if (i > xs) { 453c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs; 454c4762a1bSJed Brown } 455c4762a1bSJed Brown } 456c4762a1bSJed Brown if (i == ms + rsbwidth) { /* interface between medium and slow regions */ 457c4762a1bSJed Brown for (j = 0; j < dof; j++) { 458c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; 459c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; 460c4762a1bSJed Brown } 4619566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 462c4762a1bSJed Brown if (i < xs + xm) { 463c4762a1bSJed Brown for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs; 464c4762a1bSJed Brown islow++; 465c4762a1bSJed Brown } 466c4762a1bSJed Brown } 467c4762a1bSJed Brown if (i > ms + rsbwidth) { /* slow region */ 468c4762a1bSJed Brown for (j = 0; j < dof; j++) { 469c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; 470c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; 471c4762a1bSJed Brown } 4729566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 473c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hxs)); /* Max allowable value of 1/Delta t */ 474c4762a1bSJed Brown if (i > xs) { 475c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs; 476c4762a1bSJed Brown } 477c4762a1bSJed Brown if (i < xs + xm) { 478c4762a1bSJed Brown for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs; 479c4762a1bSJed Brown islow++; 480c4762a1bSJed Brown } 481c4762a1bSJed Brown } 482c4762a1bSJed Brown } 4839566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); 4849566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 4859566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope)); 4869566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc)); 4879566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&cfl_idt, &ctx->cfl_idt, 1, MPIU_SCALAR, MPIU_MAX, PetscObjectComm((PetscObject)da))); 488c4762a1bSJed Brown PetscFunctionReturn(0); 489c4762a1bSJed Brown } 490c4762a1bSJed Brown 4919371c9d4SSatish Balay PetscErrorCode FVRHSFunctionslowbuffer_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) { 492c4762a1bSJed Brown FVCtx *ctx = (FVCtx *)vctx; 493c4762a1bSJed Brown PetscInt i, j, k, Mx, dof, xs, xm, islowbuffer = 0, sm = ctx->sm, ms = ctx->ms, lsbwidth = ctx->lsbwidth, rsbwidth = ctx->rsbwidth; 494c4762a1bSJed Brown PetscReal hxs, hxm, hxf; 495c4762a1bSJed Brown PetscScalar *x, *f, *slope; 496c4762a1bSJed Brown Vec Xloc; 497c4762a1bSJed Brown DM da; 498c4762a1bSJed Brown 499c4762a1bSJed Brown PetscFunctionBeginUser; 5009566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 5019566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc)); 5029566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); 503c4762a1bSJed Brown hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm; 504c4762a1bSJed Brown hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm); 505c4762a1bSJed Brown hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf); 5069566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); 5079566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); 5089566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 5099566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x)); 5109566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 5119566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); 5129566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 513c4762a1bSJed Brown 514c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 515c4762a1bSJed Brown for (i = xs - 2; i < 0; i++) { 516c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; 517c4762a1bSJed Brown } 518c4762a1bSJed Brown for (i = Mx; i < xs + xm + 2; i++) { 519c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; 520c4762a1bSJed Brown } 521c4762a1bSJed Brown } 522c4762a1bSJed Brown for (i = xs - 1; i < xs + xm + 1; i++) { 523c4762a1bSJed Brown struct _LimitInfo info; 524c4762a1bSJed Brown PetscScalar *cjmpL, *cjmpR; 525c4762a1bSJed Brown if ((i > sm - lsbwidth - 2 && i < sm + 1) || (i > ms - 2 && i < ms + rsbwidth + 1)) { 526c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 5279566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds)); 528c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 5299566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof)); 530c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 531c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 532c4762a1bSJed Brown for (j = 0; j < dof; j++) { 533c4762a1bSJed Brown PetscScalar jmpL, jmpR; 534c4762a1bSJed Brown jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j]; 535c4762a1bSJed Brown jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j]; 536c4762a1bSJed Brown for (k = 0; k < dof; k++) { 537c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL; 538c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR; 539c4762a1bSJed Brown } 540c4762a1bSJed Brown } 541c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 542c4762a1bSJed Brown info.m = dof; 543c4762a1bSJed Brown info.hxs = hxs; 544c4762a1bSJed Brown info.hxm = hxm; 545c4762a1bSJed Brown info.hxf = hxf; 546c4762a1bSJed Brown (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope); 547c4762a1bSJed Brown for (j = 0; j < dof; j++) { 548c4762a1bSJed Brown PetscScalar tmp = 0; 549c4762a1bSJed Brown for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k]; 550c4762a1bSJed Brown slope[i * dof + j] = tmp; 551c4762a1bSJed Brown } 552c4762a1bSJed Brown } 553c4762a1bSJed Brown } 554c4762a1bSJed Brown 555c4762a1bSJed Brown for (i = xs; i < xs + xm + 1; i++) { 556c4762a1bSJed Brown PetscReal maxspeed; 557c4762a1bSJed Brown PetscScalar *uL, *uR; 558c4762a1bSJed Brown uL = &ctx->uLR[0]; 559c4762a1bSJed Brown uR = &ctx->uLR[dof]; 560c4762a1bSJed Brown if (i == sm - lsbwidth) { 561c4762a1bSJed Brown for (j = 0; j < dof; j++) { 562c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; 563c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; 564c4762a1bSJed Brown } 5659566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 566c4762a1bSJed Brown if (i < xs + xm) { 567c4762a1bSJed Brown for (j = 0; j < dof; j++) f[islowbuffer * dof + j] += ctx->flux[j] / hxs; 568c4762a1bSJed Brown islowbuffer++; 569c4762a1bSJed Brown } 570c4762a1bSJed Brown } 571c4762a1bSJed Brown if (i > sm - lsbwidth && i < sm) { 572c4762a1bSJed Brown for (j = 0; j < dof; j++) { 573c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; 574c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; 575c4762a1bSJed Brown } 5769566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 577c4762a1bSJed Brown if (i > xs) { 578c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(islowbuffer - 1) * dof + j] -= ctx->flux[j] / hxs; 579c4762a1bSJed Brown } 580c4762a1bSJed Brown if (i < xs + xm) { 581c4762a1bSJed Brown for (j = 0; j < dof; j++) f[islowbuffer * dof + j] += ctx->flux[j] / hxs; 582c4762a1bSJed Brown islowbuffer++; 583c4762a1bSJed Brown } 584c4762a1bSJed Brown } 585c4762a1bSJed Brown if (i == sm) { /* interface between the slow region and the medium region */ 586c4762a1bSJed Brown for (j = 0; j < dof; j++) { 587c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; 588c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; 589c4762a1bSJed Brown } 5909566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 591c4762a1bSJed Brown if (i > xs) { 592c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(islowbuffer - 1) * dof + j] -= ctx->flux[j] / hxs; 593c4762a1bSJed Brown } 594c4762a1bSJed Brown } 595c4762a1bSJed Brown if (i == ms) { /* interface between the medium region and the slow region */ 596c4762a1bSJed Brown for (j = 0; j < dof; j++) { 597c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; 598c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; 599c4762a1bSJed Brown } 6009566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 601c4762a1bSJed Brown if (i < xs + xm) { 602c4762a1bSJed Brown for (j = 0; j < dof; j++) f[islowbuffer * dof + j] += ctx->flux[j] / hxs; 603c4762a1bSJed Brown islowbuffer++; 604c4762a1bSJed Brown } 605c4762a1bSJed Brown } 606c4762a1bSJed Brown if (i > ms && i < ms + rsbwidth) { 607c4762a1bSJed Brown for (j = 0; j < dof; j++) { 608c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; 609c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; 610c4762a1bSJed Brown } 6119566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 612c4762a1bSJed Brown if (i > xs) { 613c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(islowbuffer - 1) * dof + j] -= ctx->flux[j] / hxs; 614c4762a1bSJed Brown } 615c4762a1bSJed Brown if (i < xs + xm) { 616c4762a1bSJed Brown for (j = 0; j < dof; j++) f[islowbuffer * dof + j] += ctx->flux[j] / hxs; 617c4762a1bSJed Brown islowbuffer++; 618c4762a1bSJed Brown } 619c4762a1bSJed Brown } 620c4762a1bSJed Brown if (i == ms + rsbwidth) { 621c4762a1bSJed Brown for (j = 0; j < dof; j++) { 622c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; 623c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; 624c4762a1bSJed Brown } 6259566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 626c4762a1bSJed Brown if (i > xs) { 627c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(islowbuffer - 1) * dof + j] -= ctx->flux[j] / hxs; 628c4762a1bSJed Brown } 629c4762a1bSJed Brown } 630c4762a1bSJed Brown } 6319566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); 6329566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 6339566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope)); 6349566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc)); 635c4762a1bSJed Brown PetscFunctionReturn(0); 636c4762a1bSJed Brown } 637c4762a1bSJed Brown 638c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for medium components ----------------------------------- */ 6399371c9d4SSatish Balay PetscErrorCode FVRHSFunctionmedium_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) { 640c4762a1bSJed Brown FVCtx *ctx = (FVCtx *)vctx; 641c4762a1bSJed Brown 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; 642c4762a1bSJed Brown PetscReal hxs, hxm, hxf; 643c4762a1bSJed Brown PetscScalar *x, *f, *slope; 644c4762a1bSJed Brown Vec Xloc; 645c4762a1bSJed Brown DM da; 646c4762a1bSJed Brown 647c4762a1bSJed Brown PetscFunctionBeginUser; 6489566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 6499566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc)); 6509566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); 651c4762a1bSJed Brown hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm; 652c4762a1bSJed Brown hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm); 653c4762a1bSJed Brown hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf); 6549566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); 6559566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); 6569566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 6579566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x)); 6589566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 6599566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); 6609566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 661c4762a1bSJed Brown 662c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 663c4762a1bSJed Brown for (i = xs - 2; i < 0; i++) { 664c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; 665c4762a1bSJed Brown } 666c4762a1bSJed Brown for (i = Mx; i < xs + xm + 2; i++) { 667c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; 668c4762a1bSJed Brown } 669c4762a1bSJed Brown } 670c4762a1bSJed Brown for (i = xs - 1; i < xs + xm + 1; i++) { 671c4762a1bSJed Brown struct _LimitInfo info; 672c4762a1bSJed Brown PetscScalar *cjmpL, *cjmpR; 673c4762a1bSJed Brown if ((i > sm - 2 && i < mf - lmbwidth + 1) || (i > fm + rmbwidth - 2 && i < ms + 1)) { /* slow components and the first and last fast components */ 674c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 6759566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds)); 676c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 6779566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof)); 678c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 679c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 680c4762a1bSJed Brown for (j = 0; j < dof; j++) { 681c4762a1bSJed Brown PetscScalar jmpL, jmpR; 682c4762a1bSJed Brown jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j]; 683c4762a1bSJed Brown jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j]; 684c4762a1bSJed Brown for (k = 0; k < dof; k++) { 685c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL; 686c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR; 687c4762a1bSJed Brown } 688c4762a1bSJed Brown } 689c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 690c4762a1bSJed Brown info.m = dof; 691c4762a1bSJed Brown info.hxs = hxs; 692c4762a1bSJed Brown info.hxm = hxm; 693c4762a1bSJed Brown info.hxf = hxf; 694c4762a1bSJed Brown (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope); 695c4762a1bSJed Brown for (j = 0; j < dof; j++) { 696c4762a1bSJed Brown PetscScalar tmp = 0; 697c4762a1bSJed Brown for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k]; 698c4762a1bSJed Brown slope[i * dof + j] = tmp; 699c4762a1bSJed Brown } 700c4762a1bSJed Brown } 701c4762a1bSJed Brown } 702c4762a1bSJed Brown 703c4762a1bSJed Brown for (i = xs; i < xs + xm + 1; i++) { 704c4762a1bSJed Brown PetscReal maxspeed; 705c4762a1bSJed Brown PetscScalar *uL, *uR; 706c4762a1bSJed Brown uL = &ctx->uLR[0]; 707c4762a1bSJed Brown uR = &ctx->uLR[dof]; 708c4762a1bSJed Brown if (i == sm) { /* interface between slow and medium regions */ 709c4762a1bSJed Brown for (j = 0; j < dof; j++) { 710c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2; 711c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; 712c4762a1bSJed Brown } 7139566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 714c4762a1bSJed Brown if (i < xs + xm) { 715c4762a1bSJed Brown for (j = 0; j < dof; j++) f[imedium * dof + j] += ctx->flux[j] / hxm; 716c4762a1bSJed Brown imedium++; 717c4762a1bSJed Brown } 718c4762a1bSJed Brown } 719c4762a1bSJed Brown if (i > sm && i < mf - lmbwidth) { /* medium region */ 720c4762a1bSJed Brown for (j = 0; j < dof; j++) { 721c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; 722c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; 723c4762a1bSJed Brown } 7249566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 725c4762a1bSJed Brown if (i > xs) { 726c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(imedium - 1) * dof + j] -= ctx->flux[j] / hxm; 727c4762a1bSJed Brown } 728c4762a1bSJed Brown if (i < xs + xm) { 729c4762a1bSJed Brown for (j = 0; j < dof; j++) f[imedium * dof + j] += ctx->flux[j] / hxm; 730c4762a1bSJed Brown imedium++; 731c4762a1bSJed Brown } 732c4762a1bSJed Brown } 733c4762a1bSJed Brown if (i == mf - lmbwidth) { /* interface between medium and fast regions */ 734c4762a1bSJed Brown for (j = 0; j < dof; j++) { 735c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; 736c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; 737c4762a1bSJed Brown } 7389566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 739c4762a1bSJed Brown if (i > xs) { 740c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(imedium - 1) * dof + j] -= ctx->flux[j] / hxm; 741c4762a1bSJed Brown } 742c4762a1bSJed Brown } 743c4762a1bSJed Brown if (i == fm + rmbwidth) { /* interface between fast and medium regions */ 744c4762a1bSJed Brown for (j = 0; j < dof; j++) { 745c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; 746c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; 747c4762a1bSJed Brown } 7489566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 749c4762a1bSJed Brown if (i < xs + xm) { 750c4762a1bSJed Brown for (j = 0; j < dof; j++) f[imedium * dof + j] += ctx->flux[j] / hxm; 751c4762a1bSJed Brown imedium++; 752c4762a1bSJed Brown } 753c4762a1bSJed Brown } 754c4762a1bSJed Brown if (i > fm + rmbwidth && i < ms) { /* medium region */ 755c4762a1bSJed Brown for (j = 0; j < dof; j++) { 756c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; 757c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; 758c4762a1bSJed Brown } 7599566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 760c4762a1bSJed Brown if (i > xs) { 761c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(imedium - 1) * dof + j] -= ctx->flux[j] / hxm; 762c4762a1bSJed Brown } 763c4762a1bSJed Brown if (i < xs + xm) { 764c4762a1bSJed Brown for (j = 0; j < dof; j++) f[imedium * dof + j] += ctx->flux[j] / hxm; 765c4762a1bSJed Brown imedium++; 766c4762a1bSJed Brown } 767c4762a1bSJed Brown } 768c4762a1bSJed Brown if (i == ms) { /* interface between medium and slow regions */ 769c4762a1bSJed Brown for (j = 0; j < dof; j++) { 770c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; 771c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2; 772c4762a1bSJed Brown } 7739566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 774c4762a1bSJed Brown if (i > xs) { 775c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(imedium - 1) * dof + j] -= ctx->flux[j] / hxm; 776c4762a1bSJed Brown } 777c4762a1bSJed Brown } 778c4762a1bSJed Brown } 7799566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); 7809566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 7819566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope)); 7829566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc)); 783c4762a1bSJed Brown PetscFunctionReturn(0); 784c4762a1bSJed Brown } 785c4762a1bSJed Brown 7869371c9d4SSatish Balay PetscErrorCode FVRHSFunctionmediumbuffer_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) { 787c4762a1bSJed Brown FVCtx *ctx = (FVCtx *)vctx; 788c4762a1bSJed Brown PetscInt i, j, k, Mx, dof, xs, xm, imediumbuffer = 0, mf = ctx->mf, fm = ctx->fm, lmbwidth = ctx->lmbwidth, rmbwidth = ctx->rmbwidth; 789c4762a1bSJed Brown PetscReal hxs, hxm, hxf; 790c4762a1bSJed Brown PetscScalar *x, *f, *slope; 791c4762a1bSJed Brown Vec Xloc; 792c4762a1bSJed Brown DM da; 793c4762a1bSJed Brown 794c4762a1bSJed Brown PetscFunctionBeginUser; 7959566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 7969566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc)); 7979566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); 798c4762a1bSJed Brown hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm; 799c4762a1bSJed Brown hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm); 800c4762a1bSJed Brown hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf); 8019566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); 8029566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); 8039566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 8049566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x)); 8059566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 8069566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); 8079566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 808c4762a1bSJed Brown 809c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 810c4762a1bSJed Brown for (i = xs - 2; i < 0; i++) { 811c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; 812c4762a1bSJed Brown } 813c4762a1bSJed Brown for (i = Mx; i < xs + xm + 2; i++) { 814c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; 815c4762a1bSJed Brown } 816c4762a1bSJed Brown } 817c4762a1bSJed Brown for (i = xs - 1; i < xs + xm + 1; i++) { 818c4762a1bSJed Brown struct _LimitInfo info; 819c4762a1bSJed Brown PetscScalar *cjmpL, *cjmpR; 820c4762a1bSJed Brown if ((i > mf - lmbwidth - 2 && i < mf + 1) || (i > fm - 2 && i < fm + rmbwidth + 1)) { 821c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 8229566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds)); 823c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 8249566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof)); 825c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 826c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 827c4762a1bSJed Brown for (j = 0; j < dof; j++) { 828c4762a1bSJed Brown PetscScalar jmpL, jmpR; 829c4762a1bSJed Brown jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j]; 830c4762a1bSJed Brown jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j]; 831c4762a1bSJed Brown for (k = 0; k < dof; k++) { 832c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL; 833c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR; 834c4762a1bSJed Brown } 835c4762a1bSJed Brown } 836c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 837c4762a1bSJed Brown info.m = dof; 838c4762a1bSJed Brown info.hxs = hxs; 839c4762a1bSJed Brown info.hxm = hxm; 840c4762a1bSJed Brown info.hxf = hxf; 841c4762a1bSJed Brown (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope); 842c4762a1bSJed Brown for (j = 0; j < dof; j++) { 843c4762a1bSJed Brown PetscScalar tmp = 0; 844c4762a1bSJed Brown for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k]; 845c4762a1bSJed Brown slope[i * dof + j] = tmp; 846c4762a1bSJed Brown } 847c4762a1bSJed Brown } 848c4762a1bSJed Brown } 849c4762a1bSJed Brown 850c4762a1bSJed Brown for (i = xs; i < xs + xm + 1; i++) { 851c4762a1bSJed Brown PetscReal maxspeed; 852c4762a1bSJed Brown PetscScalar *uL, *uR; 853c4762a1bSJed Brown uL = &ctx->uLR[0]; 854c4762a1bSJed Brown uR = &ctx->uLR[dof]; 855c4762a1bSJed Brown if (i == mf - lmbwidth) { 856c4762a1bSJed Brown for (j = 0; j < dof; j++) { 857c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; 858c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; 859c4762a1bSJed Brown } 8609566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 861c4762a1bSJed Brown if (i < xs + xm) { 862c4762a1bSJed Brown for (j = 0; j < dof; j++) f[imediumbuffer * dof + j] += ctx->flux[j] / hxm; 863c4762a1bSJed Brown imediumbuffer++; 864c4762a1bSJed Brown } 865c4762a1bSJed Brown } 866c4762a1bSJed Brown if (i > mf - lmbwidth && i < mf) { 867c4762a1bSJed Brown for (j = 0; j < dof; j++) { 868c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; 869c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; 870c4762a1bSJed Brown } 8719566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 872c4762a1bSJed Brown if (i > xs) { 873c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(imediumbuffer - 1) * dof + j] -= ctx->flux[j] / hxm; 874c4762a1bSJed Brown } 875c4762a1bSJed Brown if (i < xs + xm) { 876c4762a1bSJed Brown for (j = 0; j < dof; j++) f[imediumbuffer * dof + j] += ctx->flux[j] / hxm; 877c4762a1bSJed Brown imediumbuffer++; 878c4762a1bSJed Brown } 879c4762a1bSJed Brown } 880c4762a1bSJed Brown if (i == mf) { /* interface between the medium region and the fast region */ 881c4762a1bSJed Brown for (j = 0; j < dof; j++) { 882c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; 883c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2; 884c4762a1bSJed Brown } 8859566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 886c4762a1bSJed Brown if (i > xs) { 887c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(imediumbuffer - 1) * dof + j] -= ctx->flux[j] / hxm; 888c4762a1bSJed Brown } 889c4762a1bSJed Brown } 890c4762a1bSJed Brown if (i == fm) { /* interface between the fast region and the medium region */ 891c4762a1bSJed Brown for (j = 0; j < dof; j++) { 892c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2; 893c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; 894c4762a1bSJed Brown } 8959566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 896c4762a1bSJed Brown if (i < xs + xm) { 897c4762a1bSJed Brown for (j = 0; j < dof; j++) f[imediumbuffer * dof + j] += ctx->flux[j] / hxm; 898c4762a1bSJed Brown imediumbuffer++; 899c4762a1bSJed Brown } 900c4762a1bSJed Brown } 901c4762a1bSJed Brown if (i > fm && i < fm + rmbwidth) { 902c4762a1bSJed Brown for (j = 0; j < dof; j++) { 903c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; 904c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; 905c4762a1bSJed Brown } 9069566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 907c4762a1bSJed Brown if (i > xs) { 908c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(imediumbuffer - 1) * dof + j] -= ctx->flux[j] / hxm; 909c4762a1bSJed Brown } 910c4762a1bSJed Brown if (i < xs + xm) { 911c4762a1bSJed Brown for (j = 0; j < dof; j++) f[imediumbuffer * dof + j] += ctx->flux[j] / hxm; 912c4762a1bSJed Brown imediumbuffer++; 913c4762a1bSJed Brown } 914c4762a1bSJed Brown } 915c4762a1bSJed Brown if (i == fm + rmbwidth) { 916c4762a1bSJed Brown for (j = 0; j < dof; j++) { 917c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; 918c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; 919c4762a1bSJed Brown } 9209566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 921c4762a1bSJed Brown if (i > xs) { 922c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(imediumbuffer - 1) * dof + j] -= ctx->flux[j] / hxm; 923c4762a1bSJed Brown } 924c4762a1bSJed Brown } 925c4762a1bSJed Brown } 9269566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); 9279566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 9289566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope)); 9299566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc)); 930c4762a1bSJed Brown PetscFunctionReturn(0); 931c4762a1bSJed Brown } 932c4762a1bSJed Brown 933c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for fast parts ----------------------------------- */ 9349371c9d4SSatish Balay PetscErrorCode FVRHSFunctionfast_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) { 935c4762a1bSJed Brown FVCtx *ctx = (FVCtx *)vctx; 936c4762a1bSJed Brown PetscInt i, j, k, Mx, dof, xs, xm, ifast = 0, mf = ctx->mf, fm = ctx->fm; 937c4762a1bSJed Brown PetscReal hxs, hxm, hxf; 938c4762a1bSJed Brown PetscScalar *x, *f, *slope; 939c4762a1bSJed Brown Vec Xloc; 940c4762a1bSJed Brown DM da; 941c4762a1bSJed Brown 942c4762a1bSJed Brown PetscFunctionBeginUser; 9439566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 9449566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc)); 9459566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); 946c4762a1bSJed Brown hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm; 947c4762a1bSJed Brown hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm); 948c4762a1bSJed Brown hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf); 9499566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); 9509566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); 9519566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 9529566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x)); 9539566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 9549566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); 9559566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 956c4762a1bSJed Brown 957c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 958c4762a1bSJed Brown for (i = xs - 2; i < 0; i++) { 959c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; 960c4762a1bSJed Brown } 961c4762a1bSJed Brown for (i = Mx; i < xs + xm + 2; i++) { 962c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; 963c4762a1bSJed Brown } 964c4762a1bSJed Brown } 9656aad120cSJose E. Roman 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 */ 966c4762a1bSJed Brown struct _LimitInfo info; 967c4762a1bSJed Brown PetscScalar *cjmpL, *cjmpR; 968c4762a1bSJed Brown if (i > mf - 2 && i < fm + 1) { 9699566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds)); 9709566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof)); 971c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 972c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 973c4762a1bSJed Brown for (j = 0; j < dof; j++) { 974c4762a1bSJed Brown PetscScalar jmpL, jmpR; 975c4762a1bSJed Brown jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j]; 976c4762a1bSJed Brown jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j]; 977c4762a1bSJed Brown for (k = 0; k < dof; k++) { 978c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL; 979c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR; 980c4762a1bSJed Brown } 981c4762a1bSJed Brown } 982c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 983c4762a1bSJed Brown info.m = dof; 984c4762a1bSJed Brown info.hxs = hxs; 985c4762a1bSJed Brown info.hxm = hxm; 986c4762a1bSJed Brown info.hxf = hxf; 987c4762a1bSJed Brown (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope); 988c4762a1bSJed Brown for (j = 0; j < dof; j++) { 989c4762a1bSJed Brown PetscScalar tmp = 0; 990c4762a1bSJed Brown for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k]; 991c4762a1bSJed Brown slope[i * dof + j] = tmp; 992c4762a1bSJed Brown } 993c4762a1bSJed Brown } 994c4762a1bSJed Brown } 995c4762a1bSJed Brown 996c4762a1bSJed Brown for (i = xs; i < xs + xm + 1; i++) { 997c4762a1bSJed Brown PetscReal maxspeed; 998c4762a1bSJed Brown PetscScalar *uL, *uR; 999c4762a1bSJed Brown uL = &ctx->uLR[0]; 1000c4762a1bSJed Brown uR = &ctx->uLR[dof]; 1001c4762a1bSJed Brown if (i == mf) { /* interface between medium and fast regions */ 1002c4762a1bSJed Brown for (j = 0; j < dof; j++) { 1003c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2; 1004c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2; 1005c4762a1bSJed Brown } 10069566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 1007c4762a1bSJed Brown if (i < xs + xm) { 1008c4762a1bSJed Brown for (j = 0; j < dof; j++) f[ifast * dof + j] += ctx->flux[j] / hxf; 1009c4762a1bSJed Brown ifast++; 1010c4762a1bSJed Brown } 1011c4762a1bSJed Brown } 1012c4762a1bSJed Brown if (i > mf && i < fm) { /* fast region */ 1013c4762a1bSJed Brown for (j = 0; j < dof; j++) { 1014c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2; 1015c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2; 1016c4762a1bSJed Brown } 10179566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 1018c4762a1bSJed Brown if (i > xs) { 1019c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(ifast - 1) * dof + j] -= ctx->flux[j] / hxf; 1020c4762a1bSJed Brown } 1021c4762a1bSJed Brown if (i < xs + xm) { 1022c4762a1bSJed Brown for (j = 0; j < dof; j++) f[ifast * dof + j] += ctx->flux[j] / hxf; 1023c4762a1bSJed Brown ifast++; 1024c4762a1bSJed Brown } 1025c4762a1bSJed Brown } 1026c4762a1bSJed Brown if (i == fm) { /* interface between fast and medium regions */ 1027c4762a1bSJed Brown for (j = 0; j < dof; j++) { 1028c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2; 1029c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2; 1030c4762a1bSJed Brown } 10319566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed)); 1032c4762a1bSJed Brown if (i > xs) { 1033c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(ifast - 1) * dof + j] -= ctx->flux[j] / hxf; 1034c4762a1bSJed Brown } 1035c4762a1bSJed Brown } 1036c4762a1bSJed Brown } 10379566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); 10389566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 10399566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope)); 10409566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc)); 1041c4762a1bSJed Brown PetscFunctionReturn(0); 1042c4762a1bSJed Brown } 1043c4762a1bSJed Brown 10449371c9d4SSatish Balay int main(int argc, char *argv[]) { 1045c4762a1bSJed Brown char lname[256] = "mc", physname[256] = "advect", final_fname[256] = "solution.m"; 1046c4762a1bSJed Brown PetscFunctionList limiters = 0, physics = 0; 1047c4762a1bSJed Brown MPI_Comm comm; 1048c4762a1bSJed Brown TS ts; 1049c4762a1bSJed Brown DM da; 1050c4762a1bSJed Brown Vec X, X0, R; 1051c4762a1bSJed Brown FVCtx ctx; 1052c4762a1bSJed Brown 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; 1053c4762a1bSJed Brown PetscBool view_final = PETSC_FALSE; 1054c4762a1bSJed Brown PetscReal ptime; 1055c4762a1bSJed Brown 1056327415f7SBarry Smith PetscFunctionBeginUser; 10579566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, 0, help)); 1058c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 10599566063dSJacob Faibussowitsch PetscCall(PetscMemzero(&ctx, sizeof(ctx))); 1060c4762a1bSJed Brown 1061c4762a1bSJed Brown /* Register limiters to be available on the command line */ 10629566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "upwind", Limit3_Upwind)); 10639566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "lax-wendroff", Limit3_LaxWendroff)); 10649566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "beam-warming", Limit3_BeamWarming)); 10659566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "fromm", Limit3_Fromm)); 10669566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "minmod", Limit3_Minmod)); 10679566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "superbee", Limit3_Superbee)); 10689566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "mc", Limit3_MC)); 10699566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "koren3", Limit3_Koren3)); 1070c4762a1bSJed Brown 1071c4762a1bSJed Brown /* Register physical models to be available on the command line */ 10729566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&physics, "advect", PhysicsCreate_Advect)); 1073c4762a1bSJed Brown 1074c4762a1bSJed Brown ctx.comm = comm; 1075c4762a1bSJed Brown ctx.cfl = 0.9; 1076c4762a1bSJed Brown ctx.bctype = FVBC_PERIODIC; 1077c4762a1bSJed Brown ctx.xmin = -1.0; 1078c4762a1bSJed Brown ctx.xmax = 1.0; 1079d0609cedSBarry Smith PetscOptionsBegin(comm, NULL, "Finite Volume solver options", ""); 10809566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-xmin", "X min", "", ctx.xmin, &ctx.xmin, NULL)); 10819566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-xmax", "X max", "", ctx.xmax, &ctx.xmax, NULL)); 10829566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-limit", "Name of flux imiter to use", "", limiters, lname, lname, sizeof(lname), NULL)); 10839566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-draw", "Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)", "", draw, &draw, NULL)); 10849566063dSJacob 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)); 10859566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-initial", "Initial condition (depends on the physics)", "", ctx.initial, &ctx.initial, NULL)); 10869566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-exact", "Compare errors with exact solution", "", ctx.exact, &ctx.exact, NULL)); 10879566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-simulation", "Compare errors with reference solution", "", ctx.simulation, &ctx.simulation, NULL)); 10889566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-cfl", "CFL number to time step at", "", ctx.cfl, &ctx.cfl, NULL)); 10899566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-bc_type", "Boundary condition", "", FVBCTypes, (PetscEnum)ctx.bctype, (PetscEnum *)&ctx.bctype, NULL)); 10909566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-hratio", "Spacing ratio", "", ctx.hratio, &ctx.hratio, NULL)); 1091d0609cedSBarry Smith PetscOptionsEnd(); 1092c4762a1bSJed Brown 1093c4762a1bSJed Brown /* Choose the limiter from the list of registered limiters */ 10949566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(limiters, lname, &ctx.limit3)); 10953c633725SBarry Smith PetscCheck(ctx.limit3, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Limiter '%s' not found", lname); 1096c4762a1bSJed Brown 1097c4762a1bSJed Brown /* Choose the physics from the list of registered models */ 1098c4762a1bSJed Brown { 1099c4762a1bSJed Brown PetscErrorCode (*r)(FVCtx *); 11009566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(physics, physname, &r)); 11013c633725SBarry Smith PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Physics '%s' not found", physname); 1102c4762a1bSJed Brown /* Create the physics, will set the number of fields and their names */ 11039566063dSJacob Faibussowitsch PetscCall((*r)(&ctx)); 1104c4762a1bSJed Brown } 1105c4762a1bSJed Brown 1106c4762a1bSJed Brown /* Create a DMDA to manage the parallel grid */ 11079566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(comm, DM_BOUNDARY_PERIODIC, 50, ctx.physics2.dof, 2, NULL, &da)); 11089566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 11099566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 1110c4762a1bSJed Brown /* Inform the DMDA of the field names provided by the physics. */ 1111c4762a1bSJed Brown /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */ 1112*48a46eb9SPierre Jolivet for (i = 0; i < ctx.physics2.dof; i++) PetscCall(DMDASetFieldName(da, i, ctx.physics2.fieldname[i])); 11139566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); 11149566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 1115c4762a1bSJed Brown 1116c4762a1bSJed Brown /* Set coordinates of cell centers */ 11179566063dSJacob 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)); 1118c4762a1bSJed Brown 1119c4762a1bSJed Brown /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */ 11209566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(dof * dof, &ctx.R, dof * dof, &ctx.Rinv, 2 * dof, &ctx.cjmpLR, 1 * dof, &ctx.cslope)); 11219566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(2 * dof, &ctx.uLR, dof, &ctx.flux, dof, &ctx.speeds)); 1122c4762a1bSJed Brown 1123c4762a1bSJed Brown /* Create a vector to store the solution and to save the initial state */ 11249566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &X)); 11259566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &X0)); 11269566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &R)); 1127c4762a1bSJed Brown 1128c4762a1bSJed Brown /* create index for slow parts and fast parts, 1129c4762a1bSJed Brown count_slow + count_fast = Mx, counts_slow*hs = 0.5, counts_fast*hf = 0.5 */ 1130c4762a1bSJed Brown count_slow = Mx / (1 + ctx.hratio) / (1 + ctx.hratio); 1131c4762a1bSJed Brown count_medium = 2 * ctx.hratio * count_slow; 1132cad9d221SBarry Smith 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"); 1133c4762a1bSJed Brown count_fast = ctx.hratio * ctx.hratio * count_slow; 1134c4762a1bSJed Brown ctx.sm = count_slow / 2; 1135c4762a1bSJed Brown ctx.mf = ctx.sm + count_medium / 2; 1136c4762a1bSJed Brown ctx.fm = ctx.mf + count_fast; 1137c4762a1bSJed Brown ctx.ms = ctx.fm + count_medium / 2; 11389566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(xm * dof, &index_slow)); 11399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(xm * dof, &index_medium)); 11409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(xm * dof, &index_fast)); 11419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(6 * dof, &index_slowbuffer)); 11429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(6 * dof, &index_mediumbuffer)); 1143c4762a1bSJed Brown if (((AdvectCtx *)ctx.physics2.user)->a > 0) { 1144c4762a1bSJed Brown ctx.lsbwidth = 2; 1145c4762a1bSJed Brown ctx.rsbwidth = 4; 1146c4762a1bSJed Brown ctx.lmbwidth = 2; 1147c4762a1bSJed Brown ctx.rmbwidth = 4; 1148c4762a1bSJed Brown } else { 1149c4762a1bSJed Brown ctx.lsbwidth = 4; 1150c4762a1bSJed Brown ctx.rsbwidth = 2; 1151c4762a1bSJed Brown ctx.lmbwidth = 4; 1152c4762a1bSJed Brown ctx.rmbwidth = 2; 1153c4762a1bSJed Brown } 1154c4762a1bSJed Brown 1155c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 1156c4762a1bSJed Brown if (i < ctx.sm - ctx.lsbwidth || i > ctx.ms + ctx.rsbwidth - 1) 1157c4762a1bSJed Brown for (k = 0; k < dof; k++) index_slow[islow++] = i * dof + k; 1158c4762a1bSJed Brown else if ((i >= ctx.sm - ctx.lsbwidth && i < ctx.sm) || (i > ctx.ms - 1 && i <= ctx.ms + ctx.rsbwidth - 1)) 1159c4762a1bSJed Brown for (k = 0; k < dof; k++) index_slowbuffer[islowbuffer++] = i * dof + k; 1160c4762a1bSJed Brown else if (i < ctx.mf - ctx.lmbwidth || i > ctx.fm + ctx.rmbwidth - 1) 1161c4762a1bSJed Brown for (k = 0; k < dof; k++) index_medium[imedium++] = i * dof + k; 1162c4762a1bSJed Brown else if ((i >= ctx.mf - ctx.lmbwidth && i < ctx.mf) || (i > ctx.fm - 1 && i <= ctx.fm + ctx.rmbwidth - 1)) 1163c4762a1bSJed Brown for (k = 0; k < dof; k++) index_mediumbuffer[imediumbuffer++] = i * dof + k; 1164c4762a1bSJed Brown else 1165c4762a1bSJed Brown for (k = 0; k < dof; k++) index_fast[ifast++] = i * dof + k; 1166c4762a1bSJed Brown } 11679566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islow, index_slow, PETSC_COPY_VALUES, &ctx.iss)); 11689566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, imedium, index_medium, PETSC_COPY_VALUES, &ctx.ism)); 11699566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, ifast, index_fast, PETSC_COPY_VALUES, &ctx.isf)); 11709566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islowbuffer, index_slowbuffer, PETSC_COPY_VALUES, &ctx.issb)); 11719566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, imediumbuffer, index_mediumbuffer, PETSC_COPY_VALUES, &ctx.ismb)); 1172c4762a1bSJed Brown 1173c4762a1bSJed Brown /* Create a time-stepping object */ 11749566063dSJacob Faibussowitsch PetscCall(TSCreate(comm, &ts)); 11759566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da)); 11769566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, R, FVRHSFunction_3WaySplit, &ctx)); 11779566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts, "slow", ctx.iss)); 11789566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts, "medium", ctx.ism)); 11799566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts, "fast", ctx.isf)); 11809566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts, "slowbuffer", ctx.issb)); 11819566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts, "mediumbuffer", ctx.ismb)); 11829566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts, "slow", NULL, FVRHSFunctionslow_3WaySplit, &ctx)); 11839566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts, "medium", NULL, FVRHSFunctionmedium_3WaySplit, &ctx)); 11849566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts, "fast", NULL, FVRHSFunctionfast_3WaySplit, &ctx)); 11859566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts, "slowbuffer", NULL, FVRHSFunctionslowbuffer_3WaySplit, &ctx)); 11869566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts, "mediumbuffer", NULL, FVRHSFunctionmediumbuffer_3WaySplit, &ctx)); 1187c4762a1bSJed Brown 11889566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSSSP)); 11899566063dSJacob Faibussowitsch /*PetscCall(TSSetType(ts,TSMPRK));*/ 11909566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 10)); 11919566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 1192c4762a1bSJed Brown 1193c4762a1bSJed Brown /* Compute initial conditions and starting time step */ 11949566063dSJacob Faibussowitsch PetscCall(FVSample_3WaySplit(&ctx, da, 0, X0)); 11959566063dSJacob Faibussowitsch PetscCall(FVRHSFunction_3WaySplit(ts, 0, X0, X, (void *)&ctx)); /* Initial function evaluation, only used to determine max speed */ 11969566063dSJacob Faibussowitsch PetscCall(VecCopy(X0, X)); /* The function value was not used so we set X=X0 again */ 11979566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, ctx.cfl / ctx.cfl_idt)); 11989566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); /* Take runtime options */ 11999566063dSJacob Faibussowitsch PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD)); 1200c4762a1bSJed Brown { 1201c4762a1bSJed Brown PetscInt steps; 1202c4762a1bSJed Brown PetscScalar mass_initial, mass_final, mass_difference, mass_differenceg; 1203c4762a1bSJed Brown const PetscScalar *ptr_X, *ptr_X0; 1204c4762a1bSJed Brown const PetscReal hs = (ctx.xmax - ctx.xmin) / 4.0 / count_slow; 1205c4762a1bSJed Brown const PetscReal hm = (ctx.xmax - ctx.xmin) / 2.0 / count_medium; 1206c4762a1bSJed Brown const PetscReal hf = (ctx.xmax - ctx.xmin) / 4.0 / count_fast; 1207c4762a1bSJed Brown 12089566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, X)); 12099566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ptime)); 12109566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 1211c4762a1bSJed Brown /* calculate the total mass at initial time and final time */ 1212c4762a1bSJed Brown mass_initial = 0.0; 1213c4762a1bSJed Brown mass_final = 0.0; 12149566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, X0, (void *)&ptr_X0)); 12159566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, X, (void *)&ptr_X)); 1216c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 1217c4762a1bSJed Brown if (i < ctx.sm || i > ctx.ms - 1) 1218c4762a1bSJed Brown for (k = 0; k < dof; k++) { 1219c4762a1bSJed Brown mass_initial = mass_initial + hs * ptr_X0[i * dof + k]; 1220c4762a1bSJed Brown mass_final = mass_final + hs * ptr_X[i * dof + k]; 1221c4762a1bSJed Brown } 1222c4762a1bSJed Brown else if (i < ctx.mf || i > ctx.fm - 1) 1223c4762a1bSJed Brown for (k = 0; k < dof; k++) { 1224c4762a1bSJed Brown mass_initial = mass_initial + hm * ptr_X0[i * dof + k]; 1225c4762a1bSJed Brown mass_final = mass_final + hm * ptr_X[i * dof + k]; 1226c4762a1bSJed Brown } 1227c4762a1bSJed Brown else { 1228c4762a1bSJed Brown for (k = 0; k < dof; k++) { 1229c4762a1bSJed Brown mass_initial = mass_initial + hf * ptr_X0[i * dof + k]; 1230c4762a1bSJed Brown mass_final = mass_final + hf * ptr_X[i * dof + k]; 1231c4762a1bSJed Brown } 1232c4762a1bSJed Brown } 1233c4762a1bSJed Brown } 12349566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, X0, (void *)&ptr_X0)); 12359566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, X, (void *)&ptr_X)); 1236c4762a1bSJed Brown mass_difference = mass_final - mass_initial; 12379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&mass_difference, &mass_differenceg, 1, MPIU_SCALAR, MPIU_SUM, comm)); 12389566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Mass difference %g\n", (double)mass_differenceg)); 123963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Final time %g, steps %" PetscInt_FMT "\n", (double)ptime, steps)); 12409566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Maximum allowable stepsize according to CFL %g\n", (double)(1.0 / ctx.cfl_idt))); 1241c4762a1bSJed Brown if (ctx.exact) { 1242c4762a1bSJed Brown PetscReal nrm1 = 0; 12439566063dSJacob Faibussowitsch PetscCall(SolutionErrorNorms_3WaySplit(&ctx, da, ptime, X, &nrm1)); 12449566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %g\n", (double)nrm1)); 1245c4762a1bSJed Brown } 1246c4762a1bSJed Brown if (ctx.simulation) { 1247c4762a1bSJed Brown PetscReal nrm1 = 0; 1248c4762a1bSJed Brown PetscViewer fd; 1249c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN] = "binaryoutput"; 1250c4762a1bSJed Brown Vec XR; 1251c4762a1bSJed Brown PetscBool flg; 1252c4762a1bSJed Brown const PetscScalar *ptr_XR; 12539566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f", filename, sizeof(filename), &flg)); 12543c633725SBarry Smith PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate binary file with the -f option"); 12559566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, filename, FILE_MODE_READ, &fd)); 12569566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X0, &XR)); 12579566063dSJacob Faibussowitsch PetscCall(VecLoad(XR, fd)); 12589566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 12599566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &ptr_X)); 12609566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XR, &ptr_XR)); 1261c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 1262c4762a1bSJed Brown if (i < ctx.sm || i > ctx.ms - 1) 1263c4762a1bSJed Brown for (k = 0; k < dof; k++) nrm1 = nrm1 + hs * PetscAbs(ptr_X[i * dof + k] - ptr_XR[i * dof + k]); 1264c4762a1bSJed Brown else if (i < ctx.mf || i < ctx.fm - 1) 1265c4762a1bSJed Brown for (k = 0; k < dof; k++) nrm1 = nrm1 + hm * PetscAbs(ptr_X[i * dof + k] - ptr_XR[i * dof + k]); 1266c4762a1bSJed Brown else 1267c4762a1bSJed Brown for (k = 0; k < dof; k++) nrm1 = nrm1 + hf * PetscAbs(ptr_X[i * dof + k] - ptr_XR[i * dof + k]); 1268c4762a1bSJed Brown } 12699566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &ptr_X)); 12709566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XR, &ptr_XR)); 12719566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %g\n", (double)nrm1)); 12729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XR)); 1273c4762a1bSJed Brown } 1274c4762a1bSJed Brown } 1275c4762a1bSJed Brown 12769566063dSJacob Faibussowitsch PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD)); 12779566063dSJacob Faibussowitsch if (draw & 0x1) PetscCall(VecView(X0, PETSC_VIEWER_DRAW_WORLD)); 12789566063dSJacob Faibussowitsch if (draw & 0x2) PetscCall(VecView(X, PETSC_VIEWER_DRAW_WORLD)); 1279c4762a1bSJed Brown if (draw & 0x4) { 1280c4762a1bSJed Brown Vec Y; 12819566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &Y)); 12829566063dSJacob Faibussowitsch PetscCall(FVSample_3WaySplit(&ctx, da, ptime, Y)); 12839566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y, -1, X)); 12849566063dSJacob Faibussowitsch PetscCall(VecView(Y, PETSC_VIEWER_DRAW_WORLD)); 12859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 1286c4762a1bSJed Brown } 1287c4762a1bSJed Brown 1288c4762a1bSJed Brown if (view_final) { 1289c4762a1bSJed Brown PetscViewer viewer; 12909566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, final_fname, &viewer)); 12919566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB)); 12929566063dSJacob Faibussowitsch PetscCall(VecView(X, viewer)); 12939566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 12949566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 1295c4762a1bSJed Brown } 1296c4762a1bSJed Brown 1297c4762a1bSJed Brown /* Clean up */ 12989566063dSJacob Faibussowitsch PetscCall((*ctx.physics2.destroy)(ctx.physics2.user)); 12999566063dSJacob Faibussowitsch for (i = 0; i < ctx.physics2.dof; i++) PetscCall(PetscFree(ctx.physics2.fieldname[i])); 13009566063dSJacob Faibussowitsch PetscCall(PetscFree4(ctx.R, ctx.Rinv, ctx.cjmpLR, ctx.cslope)); 13019566063dSJacob Faibussowitsch PetscCall(PetscFree3(ctx.uLR, ctx.flux, ctx.speeds)); 13029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 13039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X0)); 13049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&R)); 13059566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 13069566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 13079566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.iss)); 13089566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.ism)); 13099566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.isf)); 13109566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.issb)); 13119566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.ismb)); 13129566063dSJacob Faibussowitsch PetscCall(PetscFree(index_slow)); 13139566063dSJacob Faibussowitsch PetscCall(PetscFree(index_medium)); 13149566063dSJacob Faibussowitsch PetscCall(PetscFree(index_fast)); 13159566063dSJacob Faibussowitsch PetscCall(PetscFree(index_slowbuffer)); 13169566063dSJacob Faibussowitsch PetscCall(PetscFree(index_mediumbuffer)); 13179566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&limiters)); 13189566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&physics)); 13199566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 1320b122ec5aSJacob Faibussowitsch return 0; 1321c4762a1bSJed Brown } 1322c4762a1bSJed Brown 1323c4762a1bSJed Brown /*TEST 1324c4762a1bSJed Brown 1325c4762a1bSJed Brown build: 1326f56ea12dSJed Brown requires: !complex 1327c4762a1bSJed Brown depends: finitevolume1d.c 1328c4762a1bSJed Brown 1329c4762a1bSJed Brown test: 1330c4762a1bSJed Brown suffix: 1 1331c4762a1bSJed Brown 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 1332c4762a1bSJed Brown 1333c4762a1bSJed Brown test: 1334c4762a1bSJed Brown suffix: 2 1335c4762a1bSJed Brown 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 1336c4762a1bSJed Brown 1337c4762a1bSJed Brown TEST*/ 1338