1c4762a1bSJed Brown /* 2c4762a1bSJed Brown Note: 3c4762a1bSJed Brown To performance a convergence study, a reference solution can be generated with a small stepsize and stored into a binary file, e.g. -ts_type ssp -ts_dt 1e-5 -ts_max_steps 30000 -ts_view_solution binary:reference.bin 4c4762a1bSJed Brown Errors can be computed in the following runs with -simulation -f reference.bin 5c4762a1bSJed Brown 6c4762a1bSJed Brown Multirate RK options: 7c4762a1bSJed Brown -ts_rk_dtratio is the ratio between larger time step size and small time step size 8c4762a1bSJed Brown -ts_rk_multirate_type has three choices: none (for single step RK) 9c4762a1bSJed Brown combined (for multirate method and user just need to provide the same RHS with the single step RK) 10c4762a1bSJed Brown partitioned (for multiraet method and user need to provide two RHS, one is for fast components and the orther is for slow components 11c4762a1bSJed Brown */ 12c4762a1bSJed Brown 13c4762a1bSJed Brown static const char help[] = "1D periodic Finite Volume solver in slope-limiter form with semidiscrete time stepping.\n" 14c4762a1bSJed Brown " advection - Variable coefficient scalar advection\n" 15c4762a1bSJed Brown " u_t + (a*u)_x = 0\n" 16c4762a1bSJed Brown " a is a piecewise function with two values say a0 and a1 (both a0 and a1 are constant).\n" 17c4762a1bSJed Brown " in this toy problem we assume a=a0 when x<0 and a=a1 when x>0 with a0<a1 which has the same discontinuous point with initial condition.\n" 18c4762a1bSJed Brown " we don't provide an exact solution, so you need to generate reference solution to do convergnece staudy,\n" 19c4762a1bSJed Brown " more precisely, you need firstly generate a reference to a binary file say file.bin, then on commend line,\n" 20c4762a1bSJed Brown " you should type -simulation -f file.bin.\n" 21c4762a1bSJed Brown " you can choose the number of grids by -da_grid_x.\n" 22c4762a1bSJed Brown " you can choose the value of a by -physics_advect_a1 and -physics_advect_a2.\n"; 23c4762a1bSJed Brown 24c4762a1bSJed Brown #include <petscts.h> 25c4762a1bSJed Brown #include <petscdm.h> 26c4762a1bSJed Brown #include <petscdmda.h> 27c4762a1bSJed Brown #include <petscdraw.h> 28c4762a1bSJed Brown #include <petsc/private/tsimpl.h> 29c4762a1bSJed Brown 30c4762a1bSJed Brown #include <petsc/private/kernels/blockinvert.h> /* For the Kernel_*_gets_* stuff for BAIJ */ 31c4762a1bSJed Brown 32c4762a1bSJed Brown #include "finitevolume1d.h" 33c4762a1bSJed Brown 349371c9d4SSatish Balay static inline PetscReal RangeMod(PetscReal a, PetscReal xmin, PetscReal xmax) { 359371c9d4SSatish Balay PetscReal range = xmax - xmin; 369371c9d4SSatish Balay return xmin + PetscFmodReal(range + PetscFmodReal(a, range), range); 379371c9d4SSatish Balay } 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* --------------------------------- Advection ----------------------------------- */ 40c4762a1bSJed Brown 41c4762a1bSJed Brown typedef struct { 42c4762a1bSJed Brown PetscReal a[2]; /* advective velocity */ 43c4762a1bSJed Brown } AdvectCtx; 44c4762a1bSJed Brown 459371c9d4SSatish Balay static PetscErrorCode PhysicsRiemann_Advect(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed, PetscReal x, PetscReal xmin, PetscReal xmax) { 46c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx *)vctx; 47c4762a1bSJed Brown PetscReal *speed; 48c4762a1bSJed Brown 49c4762a1bSJed Brown PetscFunctionBeginUser; 50c4762a1bSJed Brown speed = ctx->a; 519371c9d4SSatish Balay if (x == 0 || x == xmin || x == xmax) 529371c9d4SSatish Balay flux[0] = PetscMax(0, (speed[0] + speed[1]) / 2.0) * uL[0] + PetscMin(0, (speed[0] + speed[1]) / 2.0) * uR[0]; /* if condition need to be changed base on different problem, '0' is the discontinuous point of a */ 53c4762a1bSJed Brown else if (x < 0) flux[0] = PetscMax(0, speed[0]) * uL[0] + PetscMin(0, speed[0]) * uR[0]; /* else if condition need to be changed based on diferent problem, 'x = 0' is discontinuous point of a */ 54c4762a1bSJed Brown else flux[0] = PetscMax(0, speed[1]) * uL[0] + PetscMin(0, speed[1]) * uR[0]; 55c4762a1bSJed Brown *maxspeed = *speed; 56c4762a1bSJed Brown PetscFunctionReturn(0); 57c4762a1bSJed Brown } 58c4762a1bSJed Brown 599371c9d4SSatish Balay static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx, PetscInt m, const PetscScalar *u, PetscScalar *X, PetscScalar *Xi, PetscReal *speeds, PetscReal x) { 60c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx *)vctx; 61c4762a1bSJed Brown 62c4762a1bSJed Brown PetscFunctionBeginUser; 63c4762a1bSJed Brown X[0] = 1.; 64c4762a1bSJed Brown Xi[0] = 1.; 65c4762a1bSJed Brown if (x < 0) speeds[0] = ctx->a[0]; /* x is discontinuous point of a */ 66c4762a1bSJed Brown else speeds[0] = ctx->a[1]; 67c4762a1bSJed Brown PetscFunctionReturn(0); 68c4762a1bSJed Brown } 69c4762a1bSJed Brown 709371c9d4SSatish Balay static PetscErrorCode PhysicsSample_Advect(void *vctx, PetscInt initial, FVBCType bctype, PetscReal xmin, PetscReal xmax, PetscReal t, PetscReal x, PetscReal *u) { 71c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx *)vctx; 72c4762a1bSJed Brown PetscReal *a = ctx->a, x0; 73c4762a1bSJed Brown 74c4762a1bSJed Brown PetscFunctionBeginUser; 75c4762a1bSJed Brown if (x < 0) { /* x is cell center */ 76c4762a1bSJed Brown switch (bctype) { 77c4762a1bSJed Brown case FVBC_OUTFLOW: x0 = x - a[0] * t; break; 78c4762a1bSJed Brown case FVBC_PERIODIC: x0 = RangeMod(x - a[0] * t, xmin, xmax); break; 79c4762a1bSJed Brown default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown BCType"); 80c4762a1bSJed Brown } 81c4762a1bSJed Brown switch (initial) { 82c4762a1bSJed Brown case 0: u[0] = (x0 < 0) ? 1 : -1; break; 83c4762a1bSJed Brown case 1: u[0] = (x0 < 0) ? -1 : 1; break; 84c4762a1bSJed Brown case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break; 85c4762a1bSJed Brown case 3: u[0] = PetscSinReal(2 * PETSC_PI * x0); break; 86c4762a1bSJed Brown case 4: u[0] = PetscAbs(x0); break; 87c4762a1bSJed Brown case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2 * PETSC_PI * x0)); break; 88c4762a1bSJed Brown case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2 - x0 : 0)); break; 89c4762a1bSJed Brown case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI * x0), 10.0); break; 90c4762a1bSJed Brown default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition"); 91c4762a1bSJed Brown } 929371c9d4SSatish Balay } else { 93c4762a1bSJed Brown switch (bctype) { 94c4762a1bSJed Brown case FVBC_OUTFLOW: x0 = x - a[1] * t; break; 95c4762a1bSJed Brown case FVBC_PERIODIC: x0 = RangeMod(x - a[1] * t, xmin, xmax); break; 96c4762a1bSJed Brown default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown BCType"); 97c4762a1bSJed Brown } 98c4762a1bSJed Brown switch (initial) { 99c4762a1bSJed Brown case 0: u[0] = (x0 < 0) ? 1 : -1; break; 100c4762a1bSJed Brown case 1: u[0] = (x0 < 0) ? -1 : 1; break; 101c4762a1bSJed Brown case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break; 102c4762a1bSJed Brown case 3: u[0] = PetscSinReal(2 * PETSC_PI * x0); break; 103c4762a1bSJed Brown case 4: u[0] = PetscAbs(x0); break; 104c4762a1bSJed Brown case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2 * PETSC_PI * x0)); break; 105c4762a1bSJed Brown case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2 - x0 : 0)); break; 106c4762a1bSJed Brown case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI * x0), 10.0); break; 107c4762a1bSJed Brown default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition"); 108c4762a1bSJed Brown } 109c4762a1bSJed Brown } 110c4762a1bSJed Brown PetscFunctionReturn(0); 111c4762a1bSJed Brown } 112c4762a1bSJed Brown 1139371c9d4SSatish Balay static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx) { 114c4762a1bSJed Brown AdvectCtx *user; 115c4762a1bSJed Brown 116c4762a1bSJed Brown PetscFunctionBeginUser; 1179566063dSJacob Faibussowitsch PetscCall(PetscNew(&user)); 118c4762a1bSJed Brown ctx->physics.sample = PhysicsSample_Advect; 119c4762a1bSJed Brown ctx->physics.riemann = PhysicsRiemann_Advect; 120c4762a1bSJed Brown ctx->physics.characteristic = PhysicsCharacteristic_Advect; 121c4762a1bSJed Brown ctx->physics.destroy = PhysicsDestroy_SimpleFree; 122c4762a1bSJed Brown ctx->physics.user = user; 123c4762a1bSJed Brown ctx->physics.dof = 1; 1249566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy("u", &ctx->physics.fieldname[0])); 125c4762a1bSJed Brown user->a[0] = 1; 126c4762a1bSJed Brown user->a[1] = 1; 127d0609cedSBarry Smith PetscOptionsBegin(ctx->comm, ctx->prefix, "Options for advection", ""); 128c4762a1bSJed Brown { 1299566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-physics_advect_a1", "Speed1", "", user->a[0], &user->a[0], NULL)); 1309566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-physics_advect_a2", "Speed2", "", user->a[1], &user->a[1], NULL)); 131c4762a1bSJed Brown } 132d0609cedSBarry Smith PetscOptionsEnd(); 133c4762a1bSJed Brown PetscFunctionReturn(0); 134c4762a1bSJed Brown } 135c4762a1bSJed Brown 136c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for slow parts ----------------------------------- */ 137c4762a1bSJed Brown 1389371c9d4SSatish Balay PetscErrorCode FVRHSFunctionslow(TS ts, PetscReal time, Vec X, Vec F, void *vctx) { 139c4762a1bSJed Brown FVCtx *ctx = (FVCtx *)vctx; 140c4762a1bSJed Brown PetscInt i, j, k, Mx, dof, xs, xm, len_slow; 141c4762a1bSJed Brown PetscReal hx, cfl_idt = 0; 142c4762a1bSJed Brown PetscScalar *x, *f, *slope; 143c4762a1bSJed Brown Vec Xloc; 144c4762a1bSJed Brown DM da; 145c4762a1bSJed Brown 146c4762a1bSJed Brown PetscFunctionBeginUser; 1479566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 1489566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc)); 1499566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); 150c4762a1bSJed Brown hx = (ctx->xmax - ctx->xmin) / Mx; 1519566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); 1529566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); 1539566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 1549566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x)); 1559566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 1569566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); 1579566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 1589566063dSJacob Faibussowitsch PetscCall(ISGetSize(ctx->iss, &len_slow)); 159c4762a1bSJed Brown 160c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 161c4762a1bSJed Brown for (i = xs - 2; i < 0; i++) { 162c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; 163c4762a1bSJed Brown } 164c4762a1bSJed Brown for (i = Mx; i < xs + xm + 2; i++) { 165c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; 166c4762a1bSJed Brown } 167c4762a1bSJed Brown } 168c4762a1bSJed Brown for (i = xs - 1; i < xs + xm + 1; i++) { 169c4762a1bSJed Brown struct _LimitInfo info; 170c4762a1bSJed Brown PetscScalar *cjmpL, *cjmpR; 171c4762a1bSJed Brown if (i < len_slow + 1) { 172c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 1739566063dSJacob Faibussowitsch PetscCall((*ctx->physics.characteristic)(ctx->physics.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds, ctx->xmin + hx * i)); 174c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 1759566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof)); 176c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 177c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 178c4762a1bSJed Brown for (j = 0; j < dof; j++) { 179c4762a1bSJed Brown PetscScalar jmpL, jmpR; 180c4762a1bSJed Brown jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j]; 181c4762a1bSJed Brown jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j]; 182c4762a1bSJed Brown for (k = 0; k < dof; k++) { 183c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL; 184c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR; 185c4762a1bSJed Brown } 186c4762a1bSJed Brown } 187c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 188c4762a1bSJed Brown info.m = dof; 189c4762a1bSJed Brown info.hx = hx; 190c4762a1bSJed Brown (*ctx->limit)(&info, cjmpL, cjmpR, ctx->cslope); 191c4762a1bSJed Brown for (j = 0; j < dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */ 192c4762a1bSJed Brown for (j = 0; j < dof; j++) { 193c4762a1bSJed Brown PetscScalar tmp = 0; 194c4762a1bSJed Brown for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k]; 195c4762a1bSJed Brown slope[i * dof + j] = tmp; 196c4762a1bSJed Brown } 197c4762a1bSJed Brown } 198c4762a1bSJed Brown } 199c4762a1bSJed Brown 200c4762a1bSJed Brown for (i = xs; i < xs + xm + 1; i++) { 201c4762a1bSJed Brown PetscReal maxspeed; 202c4762a1bSJed Brown PetscScalar *uL, *uR; 203c4762a1bSJed Brown if (i < len_slow) { /* slow parts can be changed based on a */ 204c4762a1bSJed Brown uL = &ctx->uLR[0]; 205c4762a1bSJed Brown uR = &ctx->uLR[dof]; 206c4762a1bSJed Brown for (j = 0; j < dof; j++) { 207c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2; 208c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2; 209c4762a1bSJed Brown } 2109566063dSJacob Faibussowitsch PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax)); 211c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */ 212c4762a1bSJed Brown if (i > xs) { 213c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hx; 214c4762a1bSJed Brown } 215c4762a1bSJed Brown if (i < xs + xm) { 216c4762a1bSJed Brown for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hx; 217c4762a1bSJed Brown } 218c4762a1bSJed Brown } 219c4762a1bSJed Brown if (i == len_slow) { /* slow parts can be changed based on a */ 220c4762a1bSJed Brown uL = &ctx->uLR[0]; 221c4762a1bSJed Brown uR = &ctx->uLR[dof]; 222c4762a1bSJed Brown for (j = 0; j < dof; j++) { 223c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2; 224c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2; 225c4762a1bSJed Brown } 2269566063dSJacob Faibussowitsch PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax)); 227c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */ 228c4762a1bSJed Brown if (i > xs) { 229c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hx; 230c4762a1bSJed Brown } 231c4762a1bSJed Brown } 232c4762a1bSJed Brown } 2339566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); 2349566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 2359566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope)); 2369566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc)); 237c4762a1bSJed Brown PetscFunctionReturn(0); 238c4762a1bSJed Brown } 239c4762a1bSJed Brown 240c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for fast parts ----------------------------------- */ 241c4762a1bSJed Brown 2429371c9d4SSatish Balay PetscErrorCode FVRHSFunctionfast(TS ts, PetscReal time, Vec X, Vec F, void *vctx) { 243c4762a1bSJed Brown FVCtx *ctx = (FVCtx *)vctx; 244c4762a1bSJed Brown PetscInt i, j, k, Mx, dof, xs, xm, len_slow; 245c4762a1bSJed Brown PetscReal hx, cfl_idt = 0; 246c4762a1bSJed Brown PetscScalar *x, *f, *slope; 247c4762a1bSJed Brown Vec Xloc; 248c4762a1bSJed Brown DM da; 249c4762a1bSJed Brown 250c4762a1bSJed Brown PetscFunctionBeginUser; 2519566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 2529566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc)); 2539566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); 254c4762a1bSJed Brown hx = (ctx->xmax - ctx->xmin) / Mx; 2559566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); 2569566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); 2579566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 2589566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x)); 2599566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 2609566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); 2619566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 2629566063dSJacob Faibussowitsch PetscCall(ISGetSize(ctx->iss, &len_slow)); 263c4762a1bSJed Brown 264c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 265c4762a1bSJed Brown for (i = xs - 2; i < 0; i++) { 266c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; 267c4762a1bSJed Brown } 268c4762a1bSJed Brown for (i = Mx; i < xs + xm + 2; i++) { 269c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; 270c4762a1bSJed Brown } 271c4762a1bSJed Brown } 272c4762a1bSJed Brown for (i = xs - 1; i < xs + xm + 1; i++) { 273c4762a1bSJed Brown struct _LimitInfo info; 274c4762a1bSJed Brown PetscScalar *cjmpL, *cjmpR; 275c4762a1bSJed Brown if (i > len_slow - 2) { 2769566063dSJacob Faibussowitsch PetscCall((*ctx->physics.characteristic)(ctx->physics.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds, ctx->xmin + hx * i)); 2779566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof)); 278c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 279c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 280c4762a1bSJed Brown for (j = 0; j < dof; j++) { 281c4762a1bSJed Brown PetscScalar jmpL, jmpR; 282c4762a1bSJed Brown jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j]; 283c4762a1bSJed Brown jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j]; 284c4762a1bSJed Brown for (k = 0; k < dof; k++) { 285c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL; 286c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR; 287c4762a1bSJed Brown } 288c4762a1bSJed Brown } 289c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 290c4762a1bSJed Brown info.m = dof; 291c4762a1bSJed Brown info.hx = hx; 292c4762a1bSJed Brown (*ctx->limit)(&info, cjmpL, cjmpR, ctx->cslope); 293c4762a1bSJed Brown for (j = 0; j < dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */ 294c4762a1bSJed Brown for (j = 0; j < dof; j++) { 295c4762a1bSJed Brown PetscScalar tmp = 0; 296c4762a1bSJed Brown for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k]; 297c4762a1bSJed Brown slope[i * dof + j] = tmp; 298c4762a1bSJed Brown } 299c4762a1bSJed Brown } 300c4762a1bSJed Brown } 301c4762a1bSJed Brown 302c4762a1bSJed Brown for (i = xs; i < xs + xm + 1; i++) { 303c4762a1bSJed Brown PetscReal maxspeed; 304c4762a1bSJed Brown PetscScalar *uL, *uR; 305c4762a1bSJed Brown if (i > len_slow) { 306c4762a1bSJed Brown uL = &ctx->uLR[0]; 307c4762a1bSJed Brown uR = &ctx->uLR[dof]; 308c4762a1bSJed Brown for (j = 0; j < dof; j++) { 309c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2; 310c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2; 311c4762a1bSJed Brown } 3129566063dSJacob Faibussowitsch PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax)); 313c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */ 314c4762a1bSJed Brown if (i > xs) { 315c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - len_slow - 1) * dof + j] -= ctx->flux[j] / hx; 316c4762a1bSJed Brown } 317c4762a1bSJed Brown if (i < xs + xm) { 318c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - len_slow) * dof + j] += ctx->flux[j] / hx; 319c4762a1bSJed Brown } 320c4762a1bSJed Brown } 321c4762a1bSJed Brown if (i == len_slow) { 322c4762a1bSJed Brown uL = &ctx->uLR[0]; 323c4762a1bSJed Brown uR = &ctx->uLR[dof]; 324c4762a1bSJed Brown for (j = 0; j < dof; j++) { 325c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2; 326c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2; 327c4762a1bSJed Brown } 3289566063dSJacob Faibussowitsch PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax)); 329c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */ 330c4762a1bSJed Brown if (i < xs + xm) { 331c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - len_slow) * dof + j] += ctx->flux[j] / hx; 332c4762a1bSJed Brown } 333c4762a1bSJed Brown } 334c4762a1bSJed Brown } 3359566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); 3369566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 3379566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope)); 3389566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc)); 339c4762a1bSJed Brown PetscFunctionReturn(0); 340c4762a1bSJed Brown } 341c4762a1bSJed Brown 3429371c9d4SSatish Balay PetscErrorCode FVRHSFunctionslow2(TS ts, PetscReal time, Vec X, Vec F, void *vctx) { 343c4762a1bSJed Brown FVCtx *ctx = (FVCtx *)vctx; 344c4762a1bSJed Brown PetscInt i, j, k, Mx, dof, xs, xm, len_slow1, len_slow2; 345c4762a1bSJed Brown PetscReal hx, cfl_idt = 0; 346c4762a1bSJed Brown PetscScalar *x, *f, *slope; 347c4762a1bSJed Brown Vec Xloc; 348c4762a1bSJed Brown DM da; 349c4762a1bSJed Brown 350c4762a1bSJed Brown PetscFunctionBeginUser; 3519566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 3529566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc)); 3539566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); 354c4762a1bSJed Brown hx = (ctx->xmax - ctx->xmin) / Mx; 3559566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); 3569566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); 3579566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 3589566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x)); 3599566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 3609566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); 3619566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 3629566063dSJacob Faibussowitsch PetscCall(ISGetSize(ctx->iss, &len_slow1)); 3639566063dSJacob Faibussowitsch PetscCall(ISGetSize(ctx->iss2, &len_slow2)); 364c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 365c4762a1bSJed Brown for (i = xs - 2; i < 0; i++) { 366c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; 367c4762a1bSJed Brown } 368c4762a1bSJed Brown for (i = Mx; i < xs + xm + 2; i++) { 369c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; 370c4762a1bSJed Brown } 371c4762a1bSJed Brown } 372c4762a1bSJed Brown for (i = xs - 1; i < xs + xm + 1; i++) { 373c4762a1bSJed Brown struct _LimitInfo info; 374c4762a1bSJed Brown PetscScalar *cjmpL, *cjmpR; 375c4762a1bSJed Brown if (i < len_slow1 + len_slow2 + 1 && i > len_slow1 - 2) { 3769566063dSJacob Faibussowitsch PetscCall((*ctx->physics.characteristic)(ctx->physics.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds, ctx->xmin + hx * i)); 3779566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof)); 378c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 379c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 380c4762a1bSJed Brown for (j = 0; j < dof; j++) { 381c4762a1bSJed Brown PetscScalar jmpL, jmpR; 382c4762a1bSJed Brown jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j]; 383c4762a1bSJed Brown jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j]; 384c4762a1bSJed Brown for (k = 0; k < dof; k++) { 385c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL; 386c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR; 387c4762a1bSJed Brown } 388c4762a1bSJed Brown } 389c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 390c4762a1bSJed Brown info.m = dof; 391c4762a1bSJed Brown info.hx = hx; 392c4762a1bSJed Brown (*ctx->limit)(&info, cjmpL, cjmpR, ctx->cslope); 393c4762a1bSJed Brown for (j = 0; j < dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */ 394c4762a1bSJed Brown for (j = 0; j < dof; j++) { 395c4762a1bSJed Brown PetscScalar tmp = 0; 396c4762a1bSJed Brown for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k]; 397c4762a1bSJed Brown slope[i * dof + j] = tmp; 398c4762a1bSJed Brown } 399c4762a1bSJed Brown } 400c4762a1bSJed Brown } 401c4762a1bSJed Brown 402c4762a1bSJed Brown for (i = xs; i < xs + xm + 1; i++) { 403c4762a1bSJed Brown PetscReal maxspeed; 404c4762a1bSJed Brown PetscScalar *uL, *uR; 405c4762a1bSJed Brown if (i < len_slow1 + len_slow2 && i > len_slow1) { /* slow parts can be changed based on a */ 406c4762a1bSJed Brown uL = &ctx->uLR[0]; 407c4762a1bSJed Brown uR = &ctx->uLR[dof]; 408c4762a1bSJed Brown for (j = 0; j < dof; j++) { 409c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2; 410c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2; 411c4762a1bSJed Brown } 4129566063dSJacob Faibussowitsch PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax)); 413c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */ 414c4762a1bSJed Brown if (i > xs) { 415c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - len_slow1 - 1) * dof + j] -= ctx->flux[j] / hx; 416c4762a1bSJed Brown } 417c4762a1bSJed Brown if (i < xs + xm) { 418c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - len_slow1) * dof + j] += ctx->flux[j] / hx; 419c4762a1bSJed Brown } 420c4762a1bSJed Brown } 421c4762a1bSJed Brown if (i == len_slow1 + len_slow2) { /* slow parts can be changed based on a */ 422c4762a1bSJed Brown uL = &ctx->uLR[0]; 423c4762a1bSJed Brown uR = &ctx->uLR[dof]; 424c4762a1bSJed Brown for (j = 0; j < dof; j++) { 425c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2; 426c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2; 427c4762a1bSJed Brown } 4289566063dSJacob Faibussowitsch PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax)); 429c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */ 430c4762a1bSJed Brown if (i > xs) { 431c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - len_slow1 - 1) * dof + j] -= ctx->flux[j] / hx; 432c4762a1bSJed Brown } 433c4762a1bSJed Brown } 434c4762a1bSJed Brown if (i == len_slow1) { /* slow parts can be changed based on a */ 435c4762a1bSJed Brown uL = &ctx->uLR[0]; 436c4762a1bSJed Brown uR = &ctx->uLR[dof]; 437c4762a1bSJed Brown for (j = 0; j < dof; j++) { 438c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2; 439c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2; 440c4762a1bSJed Brown } 4419566063dSJacob Faibussowitsch PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax)); 442c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */ 443c4762a1bSJed Brown if (i < xs + xm) { 444c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - len_slow1) * dof + j] += ctx->flux[j] / hx; 445c4762a1bSJed Brown } 446c4762a1bSJed Brown } 447c4762a1bSJed Brown } 448c4762a1bSJed Brown 4499566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); 4509566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 4519566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope)); 4529566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc)); 453c4762a1bSJed Brown PetscFunctionReturn(0); 454c4762a1bSJed Brown } 455c4762a1bSJed Brown 4569371c9d4SSatish Balay PetscErrorCode FVRHSFunctionfast2(TS ts, PetscReal time, Vec X, Vec F, void *vctx) { 457c4762a1bSJed Brown FVCtx *ctx = (FVCtx *)vctx; 458c4762a1bSJed Brown PetscInt i, j, k, Mx, dof, xs, xm, len_slow1, len_slow2; 459c4762a1bSJed Brown PetscReal hx, cfl_idt = 0; 460c4762a1bSJed Brown PetscScalar *x, *f, *slope; 461c4762a1bSJed Brown Vec Xloc; 462c4762a1bSJed Brown DM da; 463c4762a1bSJed Brown 464c4762a1bSJed Brown PetscFunctionBeginUser; 4659566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 4669566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc)); 4679566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); 468c4762a1bSJed Brown hx = (ctx->xmax - ctx->xmin) / Mx; 4699566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); 4709566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); 4719566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 4729566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x)); 4739566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 4749566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); 4759566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 4769566063dSJacob Faibussowitsch PetscCall(ISGetSize(ctx->iss, &len_slow1)); 4779566063dSJacob Faibussowitsch PetscCall(ISGetSize(ctx->iss2, &len_slow2)); 478c4762a1bSJed Brown 479c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 480c4762a1bSJed Brown for (i = xs - 2; i < 0; i++) { 481c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; 482c4762a1bSJed Brown } 483c4762a1bSJed Brown for (i = Mx; i < xs + xm + 2; i++) { 484c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; 485c4762a1bSJed Brown } 486c4762a1bSJed Brown } 487c4762a1bSJed Brown for (i = xs - 1; i < xs + xm + 1; i++) { 488c4762a1bSJed Brown struct _LimitInfo info; 489c4762a1bSJed Brown PetscScalar *cjmpL, *cjmpR; 490c4762a1bSJed Brown if (i > len_slow1 + len_slow2 - 2) { 4919566063dSJacob Faibussowitsch PetscCall((*ctx->physics.characteristic)(ctx->physics.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds, ctx->xmin + hx * i)); 4929566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof)); 493c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 494c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 495c4762a1bSJed Brown for (j = 0; j < dof; j++) { 496c4762a1bSJed Brown PetscScalar jmpL, jmpR; 497c4762a1bSJed Brown jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j]; 498c4762a1bSJed Brown jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j]; 499c4762a1bSJed Brown for (k = 0; k < dof; k++) { 500c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL; 501c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR; 502c4762a1bSJed Brown } 503c4762a1bSJed Brown } 504c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 505c4762a1bSJed Brown info.m = dof; 506c4762a1bSJed Brown info.hx = hx; 507c4762a1bSJed Brown (*ctx->limit)(&info, cjmpL, cjmpR, ctx->cslope); 508c4762a1bSJed Brown for (j = 0; j < dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */ 509c4762a1bSJed Brown for (j = 0; j < dof; j++) { 510c4762a1bSJed Brown PetscScalar tmp = 0; 511c4762a1bSJed Brown for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k]; 512c4762a1bSJed Brown slope[i * dof + j] = tmp; 513c4762a1bSJed Brown } 514c4762a1bSJed Brown } 515c4762a1bSJed Brown } 516c4762a1bSJed Brown 517c4762a1bSJed Brown for (i = xs; i < xs + xm + 1; i++) { 518c4762a1bSJed Brown PetscReal maxspeed; 519c4762a1bSJed Brown PetscScalar *uL, *uR; 520c4762a1bSJed Brown if (i > len_slow1 + len_slow2) { 521c4762a1bSJed Brown uL = &ctx->uLR[0]; 522c4762a1bSJed Brown uR = &ctx->uLR[dof]; 523c4762a1bSJed Brown for (j = 0; j < dof; j++) { 524c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2; 525c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2; 526c4762a1bSJed Brown } 5279566063dSJacob Faibussowitsch PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax)); 528c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */ 529c4762a1bSJed Brown if (i > xs) { 530c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - len_slow1 - len_slow2 - 1) * dof + j] -= ctx->flux[j] / hx; 531c4762a1bSJed Brown } 532c4762a1bSJed Brown if (i < xs + xm) { 533c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - len_slow1 - len_slow2) * dof + j] += ctx->flux[j] / hx; 534c4762a1bSJed Brown } 535c4762a1bSJed Brown } 536c4762a1bSJed Brown if (i == len_slow1 + len_slow2) { 537c4762a1bSJed Brown uL = &ctx->uLR[0]; 538c4762a1bSJed Brown uR = &ctx->uLR[dof]; 539c4762a1bSJed Brown for (j = 0; j < dof; j++) { 540c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2; 541c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2; 542c4762a1bSJed Brown } 5439566063dSJacob Faibussowitsch PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax)); 544c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */ 545c4762a1bSJed Brown if (i < xs + xm) { 546c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - len_slow1 - len_slow2) * dof + j] += ctx->flux[j] / hx; 547c4762a1bSJed Brown } 548c4762a1bSJed Brown } 549c4762a1bSJed Brown } 5509566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); 5519566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 5529566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope)); 5539566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc)); 554c4762a1bSJed Brown PetscFunctionReturn(0); 555c4762a1bSJed Brown } 556c4762a1bSJed Brown 5579371c9d4SSatish Balay int main(int argc, char *argv[]) { 558c4762a1bSJed Brown char lname[256] = "mc", physname[256] = "advect", final_fname[256] = "solution.m"; 559c4762a1bSJed Brown PetscFunctionList limiters = 0, physics = 0; 560c4762a1bSJed Brown MPI_Comm comm; 561c4762a1bSJed Brown TS ts; 562c4762a1bSJed Brown DM da; 563c4762a1bSJed Brown Vec X, X0, R; 564c4762a1bSJed Brown FVCtx ctx; 565c4762a1bSJed Brown PetscInt i, k, dof, xs, xm, Mx, draw = 0, *index_slow, *index_fast, islow = 0, ifast = 0; 566c4762a1bSJed Brown PetscBool view_final = PETSC_FALSE; 567c4762a1bSJed Brown PetscReal ptime; 568c4762a1bSJed Brown 569327415f7SBarry Smith PetscFunctionBeginUser; 5709566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, 0, help)); 571c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 5729566063dSJacob Faibussowitsch PetscCall(PetscMemzero(&ctx, sizeof(ctx))); 573c4762a1bSJed Brown 574c4762a1bSJed Brown /* Register limiters to be available on the command line */ 5759566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "upwind", Limit_Upwind)); 5769566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "lax-wendroff", Limit_LaxWendroff)); 5779566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "beam-warming", Limit_BeamWarming)); 5789566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "fromm", Limit_Fromm)); 5799566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "minmod", Limit_Minmod)); 5809566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "superbee", Limit_Superbee)); 5819566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "mc", Limit_MC)); 5829566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "vanleer", Limit_VanLeer)); 5839566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "vanalbada", Limit_VanAlbada)); 5849566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "vanalbadatvd", Limit_VanAlbadaTVD)); 5859566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "koren", Limit_Koren)); 5869566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "korensym", Limit_KorenSym)); 5879566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "koren3", Limit_Koren3)); 5889566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "cada-torrilhon2", Limit_CadaTorrilhon2)); 5899566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "cada-torrilhon3-r0p1", Limit_CadaTorrilhon3R0p1)); 5909566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "cada-torrilhon3-r1", Limit_CadaTorrilhon3R1)); 5919566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "cada-torrilhon3-r10", Limit_CadaTorrilhon3R10)); 5929566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "cada-torrilhon3-r100", Limit_CadaTorrilhon3R100)); 593c4762a1bSJed Brown 594c4762a1bSJed Brown /* Register physical models to be available on the command line */ 5959566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&physics, "advect", PhysicsCreate_Advect)); 596c4762a1bSJed Brown 597c4762a1bSJed Brown ctx.comm = comm; 598c4762a1bSJed Brown ctx.cfl = 0.9; 599c4762a1bSJed Brown ctx.bctype = FVBC_PERIODIC; 600c4762a1bSJed Brown ctx.xmin = -1.0; 601c4762a1bSJed Brown ctx.xmax = 1.0; 602d0609cedSBarry Smith PetscOptionsBegin(comm, NULL, "Finite Volume solver options", ""); 6039566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-xmin", "X min", "", ctx.xmin, &ctx.xmin, NULL)); 6049566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-xmax", "X max", "", ctx.xmax, &ctx.xmax, NULL)); 6059566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-limit", "Name of flux limiter to use", "", limiters, lname, lname, sizeof(lname), NULL)); 6069566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-draw", "Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)", "", draw, &draw, NULL)); 6079566063dSJacob 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)); 6089566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-initial", "Initial condition (depends on the physics)", "", ctx.initial, &ctx.initial, NULL)); 6099566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-simulation", "Compare errors with reference solution", "", ctx.simulation, &ctx.simulation, NULL)); 6109566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-cfl", "CFL number to time step at", "", ctx.cfl, &ctx.cfl, NULL)); 6119566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-bc_type", "Boundary condition", "", FVBCTypes, (PetscEnum)ctx.bctype, (PetscEnum *)&ctx.bctype, NULL)); 6129566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-recursive_split", "Split the domain recursively", "", ctx.recursive, &ctx.recursive, NULL)); 613d0609cedSBarry Smith PetscOptionsEnd(); 614c4762a1bSJed Brown 615c4762a1bSJed Brown /* Choose the limiter from the list of registered limiters */ 6169566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(limiters, lname, &ctx.limit)); 6173c633725SBarry Smith PetscCheck(ctx.limit, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Limiter '%s' not found", lname); 618c4762a1bSJed Brown 619c4762a1bSJed Brown /* Choose the physics from the list of registered models */ 620c4762a1bSJed Brown { 621c4762a1bSJed Brown PetscErrorCode (*r)(FVCtx *); 6229566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(physics, physname, &r)); 6233c633725SBarry Smith PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Physics '%s' not found", physname); 624c4762a1bSJed Brown /* Create the physics, will set the number of fields and their names */ 6259566063dSJacob Faibussowitsch PetscCall((*r)(&ctx)); 626c4762a1bSJed Brown } 627c4762a1bSJed Brown 628c4762a1bSJed Brown /* Create a DMDA to manage the parallel grid */ 6299566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(comm, DM_BOUNDARY_PERIODIC, 50, ctx.physics.dof, 2, NULL, &da)); 6309566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 6319566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 632c4762a1bSJed Brown /* Inform the DMDA of the field names provided by the physics. */ 633c4762a1bSJed Brown /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */ 634*48a46eb9SPierre Jolivet for (i = 0; i < ctx.physics.dof; i++) PetscCall(DMDASetFieldName(da, i, ctx.physics.fieldname[i])); 6359566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); 6369566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 637c4762a1bSJed Brown 638c4762a1bSJed Brown /* Set coordinates of cell centers */ 6399566063dSJacob 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)); 640c4762a1bSJed Brown 641c4762a1bSJed Brown /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */ 6429566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(dof * dof, &ctx.R, dof * dof, &ctx.Rinv, 2 * dof, &ctx.cjmpLR, 1 * dof, &ctx.cslope)); 6439566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(2 * dof, &ctx.uLR, dof, &ctx.flux, dof, &ctx.speeds)); 644c4762a1bSJed Brown 645c4762a1bSJed Brown /* Create a vector to store the solution and to save the initial state */ 6469566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &X)); 6479566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &X0)); 6489566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &R)); 649c4762a1bSJed Brown 650c4762a1bSJed Brown /*-------------------------------- create index for slow parts and fast parts ----------------------------------------*/ 6519566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(xm * dof, &index_slow)); 6529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(xm * dof, &index_fast)); 653c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 654c4762a1bSJed Brown if (ctx.xmin + i * (ctx.xmax - ctx.xmin) / (PetscReal)Mx + 0.5 * (ctx.xmax - ctx.xmin) / (PetscReal)Mx < 0) 655c4762a1bSJed Brown for (k = 0; k < dof; k++) index_slow[islow++] = i * dof + k; 656c4762a1bSJed Brown else 657c4762a1bSJed Brown for (k = 0; k < dof; k++) index_fast[ifast++] = i * dof + k; 658c4762a1bSJed Brown } /* this step need to be changed based on discontinuous point of a */ 6599566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islow, index_slow, PETSC_COPY_VALUES, &ctx.iss)); 6609566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, ifast, index_fast, PETSC_COPY_VALUES, &ctx.isf)); 661c4762a1bSJed Brown 662c4762a1bSJed Brown /* Create a time-stepping object */ 6639566063dSJacob Faibussowitsch PetscCall(TSCreate(comm, &ts)); 6649566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da)); 6659566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, R, FVRHSFunction, &ctx)); 6669566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts, "slow", ctx.iss)); 6679566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts, "fast", ctx.isf)); 6689566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts, "slow", NULL, FVRHSFunctionslow, &ctx)); 6699566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts, "fast", NULL, FVRHSFunctionfast, &ctx)); 670c4762a1bSJed Brown 671c4762a1bSJed Brown if (ctx.recursive) { 672c4762a1bSJed Brown TS subts; 673c4762a1bSJed Brown islow = 0; 674c4762a1bSJed Brown ifast = 0; 675c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 676c4762a1bSJed Brown PetscReal coord = ctx.xmin + i * (ctx.xmax - ctx.xmin) / (PetscReal)Mx + 0.5 * (ctx.xmax - ctx.xmin) / (PetscReal)Mx; 677c4762a1bSJed Brown if (coord >= 0 && coord < ctx.xmin + (ctx.xmax - ctx.xmin) * 3 / 4.) 678c4762a1bSJed Brown for (k = 0; k < dof; k++) index_slow[islow++] = i * dof + k; 679c4762a1bSJed Brown if (coord >= ctx.xmin + (ctx.xmax - ctx.xmin) * 3 / 4.) 680c4762a1bSJed Brown for (k = 0; k < dof; k++) index_fast[ifast++] = i * dof + k; 681c4762a1bSJed Brown } /* this step need to be changed based on discontinuous point of a */ 6829566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islow, index_slow, PETSC_COPY_VALUES, &ctx.iss2)); 6839566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, ifast, index_fast, PETSC_COPY_VALUES, &ctx.isf2)); 684c4762a1bSJed Brown 6859566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts, "fast", &subts)); 6869566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(subts, "slow", ctx.iss2)); 6879566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(subts, "fast", ctx.isf2)); 6889566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(subts, "slow", NULL, FVRHSFunctionslow2, &ctx)); 6899566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(subts, "fast", NULL, FVRHSFunctionfast2, &ctx)); 690c4762a1bSJed Brown } 691c4762a1bSJed Brown 6929566063dSJacob Faibussowitsch /*PetscCall(TSSetType(ts,TSSSP));*/ 6939566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSMPRK)); 6949566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 10)); 6959566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 696c4762a1bSJed Brown 697c4762a1bSJed Brown /* Compute initial conditions and starting time step */ 6989566063dSJacob Faibussowitsch PetscCall(FVSample(&ctx, da, 0, X0)); 6999566063dSJacob Faibussowitsch PetscCall(FVRHSFunction(ts, 0, X0, X, (void *)&ctx)); /* Initial function evaluation, only used to determine max speed */ 7009566063dSJacob Faibussowitsch PetscCall(VecCopy(X0, X)); /* The function value was not used so we set X=X0 again */ 7019566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, ctx.cfl / ctx.cfl_idt)); 7029566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); /* Take runtime options */ 7039566063dSJacob Faibussowitsch PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD)); 704c4762a1bSJed Brown { 705c4762a1bSJed Brown PetscInt steps; 706c4762a1bSJed Brown PetscScalar mass_initial, mass_final, mass_difference; 707c4762a1bSJed Brown 7089566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, X)); 7099566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ptime)); 7109566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 71163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Final time %g, steps %" PetscInt_FMT "\n", (double)ptime, steps)); 712c4762a1bSJed Brown /* calculate the total mass at initial time and final time */ 713c4762a1bSJed Brown mass_initial = 0.0; 714c4762a1bSJed Brown mass_final = 0.0; 7159566063dSJacob Faibussowitsch PetscCall(VecSum(X0, &mass_initial)); 7169566063dSJacob Faibussowitsch PetscCall(VecSum(X, &mass_final)); 717c4762a1bSJed Brown mass_difference = (ctx.xmax - ctx.xmin) / (PetscScalar)Mx * (mass_final - mass_initial); 7189566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Mass difference %g\n", (double)mass_difference)); 719c4762a1bSJed Brown if (ctx.simulation) { 720c4762a1bSJed Brown PetscViewer fd; 721c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN] = "binaryoutput"; 722c4762a1bSJed Brown Vec XR; 723c4762a1bSJed Brown PetscReal nrm1, nrmsup; 724c4762a1bSJed Brown PetscBool flg; 725d8185827SBarry Smith 7269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f", filename, sizeof(filename), &flg)); 7273c633725SBarry Smith PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate binary file with the -f option"); 7289566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, filename, FILE_MODE_READ, &fd)); 7299566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X0, &XR)); 7309566063dSJacob Faibussowitsch PetscCall(VecLoad(XR, fd)); 7319566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 7329566063dSJacob Faibussowitsch PetscCall(VecAYPX(XR, -1, X)); 7339566063dSJacob Faibussowitsch PetscCall(VecNorm(XR, NORM_1, &nrm1)); 7349566063dSJacob Faibussowitsch PetscCall(VecNorm(XR, NORM_INFINITY, &nrmsup)); 735c4762a1bSJed Brown nrm1 /= Mx; 7369566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %g ||x-x_e||_sup %g\n", (double)nrm1, (double)nrmsup)); 7379566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XR)); 738c4762a1bSJed Brown } 739c4762a1bSJed Brown } 740c4762a1bSJed Brown 7419566063dSJacob Faibussowitsch PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD)); 7429566063dSJacob Faibussowitsch if (draw & 0x1) PetscCall(VecView(X0, PETSC_VIEWER_DRAW_WORLD)); 7439566063dSJacob Faibussowitsch if (draw & 0x2) PetscCall(VecView(X, PETSC_VIEWER_DRAW_WORLD)); 744c4762a1bSJed Brown if (draw & 0x4) { 745c4762a1bSJed Brown Vec Y; 7469566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &Y)); 7479566063dSJacob Faibussowitsch PetscCall(FVSample(&ctx, da, ptime, Y)); 7489566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y, -1, X)); 7499566063dSJacob Faibussowitsch PetscCall(VecView(Y, PETSC_VIEWER_DRAW_WORLD)); 7509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 751c4762a1bSJed Brown } 752c4762a1bSJed Brown 753c4762a1bSJed Brown if (view_final) { 754c4762a1bSJed Brown PetscViewer viewer; 7559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, final_fname, &viewer)); 7569566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB)); 7579566063dSJacob Faibussowitsch PetscCall(VecView(X, viewer)); 7589566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 7599566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 760c4762a1bSJed Brown } 761c4762a1bSJed Brown 762c4762a1bSJed Brown /* Clean up */ 7639566063dSJacob Faibussowitsch PetscCall((*ctx.physics.destroy)(ctx.physics.user)); 7649566063dSJacob Faibussowitsch for (i = 0; i < ctx.physics.dof; i++) PetscCall(PetscFree(ctx.physics.fieldname[i])); 7659566063dSJacob Faibussowitsch PetscCall(PetscFree4(ctx.R, ctx.Rinv, ctx.cjmpLR, ctx.cslope)); 7669566063dSJacob Faibussowitsch PetscCall(PetscFree3(ctx.uLR, ctx.flux, ctx.speeds)); 7679566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.iss)); 7689566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.isf)); 7699566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.iss2)); 7709566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.isf2)); 7719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 7729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X0)); 7739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&R)); 7749566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 7759566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 7769566063dSJacob Faibussowitsch PetscCall(PetscFree(index_slow)); 7779566063dSJacob Faibussowitsch PetscCall(PetscFree(index_fast)); 7789566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&limiters)); 7799566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&physics)); 7809566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 781b122ec5aSJacob Faibussowitsch return 0; 782c4762a1bSJed Brown } 783c4762a1bSJed Brown 784c4762a1bSJed Brown /*TEST 785c4762a1bSJed Brown build: 786f56ea12dSJed Brown requires: !complex 787c4762a1bSJed Brown depends: finitevolume1d.c 788c4762a1bSJed Brown 789c4762a1bSJed Brown test: 790c4762a1bSJed Brown args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 0 791c4762a1bSJed Brown 792c4762a1bSJed Brown test: 793c4762a1bSJed Brown suffix: 2 794c4762a1bSJed Brown args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 1 795c4762a1bSJed Brown output_file: output/ex5_1.out 796c4762a1bSJed Brown 797c4762a1bSJed Brown test: 798c4762a1bSJed Brown suffix: 3 799c4762a1bSJed Brown args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 0 800c4762a1bSJed Brown 801c4762a1bSJed Brown test: 802c4762a1bSJed Brown suffix: 4 803c4762a1bSJed Brown args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 1 804c4762a1bSJed Brown output_file: output/ex5_3.out 805c4762a1bSJed Brown 806c4762a1bSJed Brown test: 807c4762a1bSJed Brown suffix: 5 808c4762a1bSJed Brown args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 0 -recursive_split 809c4762a1bSJed Brown 810c4762a1bSJed Brown test: 811c4762a1bSJed Brown suffix: 6 812c4762a1bSJed Brown args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 1 -recursive_split 813c4762a1bSJed Brown TEST*/ 814