184ff8c8bSHong Zhang /*
284ff8c8bSHong Zhang Note:
36aad120cSJose E. Roman -hratio is the ratio between mesh size of coarse grids and fine grids
484ff8c8bSHong Zhang */
584ff8c8bSHong Zhang
684ff8c8bSHong Zhang static const char help[] = "1D periodic Finite Volume solver in slope-limiter form with semidiscrete time stepping.\n"
784ff8c8bSHong Zhang " advect - Constant coefficient scalar advection\n"
884ff8c8bSHong Zhang " u_t + (a*u)_x = 0\n"
984ff8c8bSHong Zhang " shallow - 1D Shallow water equations (Saint Venant System)\n"
1084ff8c8bSHong Zhang " h_t + (q)_x = 0 \n"
1184ff8c8bSHong Zhang " q_t + (\frac{q^2}{h} + g/2*h^2)_x = 0 \n"
1284ff8c8bSHong Zhang " where, h(x,t) denotes the height of the water and q(x,t) the momentum.\n"
1384ff8c8bSHong Zhang " for this toy problem, we choose different meshsizes for different sub-domains (slow-fast-slow), say\n"
1484ff8c8bSHong Zhang " hxs = hratio*hxf \n"
1584ff8c8bSHong Zhang " where hxs and hxf are the grid spacings for coarse and fine grids respectively.\n"
1684ff8c8bSHong Zhang " exact - Exact Riemann solver which usually needs to perform a Newton iteration to connect\n"
1784ff8c8bSHong Zhang " the states across shocks and rarefactions\n"
1884ff8c8bSHong Zhang " simulation - use reference solution which is generated by smaller time step size to be true solution,\n"
1984ff8c8bSHong Zhang " also the reference solution should be generated by user and stored in a binary file.\n"
2084ff8c8bSHong Zhang " characteristic - Limit the characteristic variables, this is usually preferred (default)\n"
2184ff8c8bSHong Zhang " bc_type - Boundary condition for the problem, options are: periodic, outflow, inflow "
2284ff8c8bSHong Zhang "Several problem descriptions (initial data, physics specific features, boundary data) can be chosen with -initial N\n\n"
2384ff8c8bSHong Zhang "The problem size should be set with -da_grid_x M\n\n";
2484ff8c8bSHong Zhang
2584ff8c8bSHong Zhang /*
2684ff8c8bSHong Zhang Example:
27*188af4bfSBarry Smith ./ex4 -da_grid_x 40 -initial 1 -hratio 1 -limit mc -ts_time_step 0.01 -ts_max_time 7.0 -ts_type mprk -ts_mprk_type 2a22 -ts_monitor_draw_solution -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 0
28*188af4bfSBarry Smith ./ex4 -da_grid_x 40 -initial 2 -hratio 1 -limit mc -ts_time_step 0.01 -ts_max_time 2.5 -ts_type mprk -ts_mprk_type 2a22 -ts_monitor_draw_solution -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 0
29*188af4bfSBarry Smith ./ex4 -da_grid_x 40 -initial 3 -hratio 1 -limit mc -ts_time_step 0.01 -ts_max_time 4.0 -ts_type mprk -ts_mprk_type 2a22 -ts_monitor_draw_solution -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 0
30*188af4bfSBarry Smith ./ex4 -da_grid_x 40 -initial 4 -hratio 1 -limit koren3 -ts_time_step 0.01 -ts_max_time 4.0 -ts_type mprk -ts_mprk_type 2a22 -ts_monitor_draw_solution -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 0
31*188af4bfSBarry Smith ./ex4 -da_grid_x 40 -initial 5 -hratio 1 -limit mc -ts_time_step 0.01 -ts_max_time 5.0 -ts_type mprk -ts_mprk_type 2a22 -ts_monitor_draw_solution -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 0
3284ff8c8bSHong Zhang
3384ff8c8bSHong Zhang Contributed by: Aidan Hamilton <aidan@udel.edu>
3484ff8c8bSHong Zhang */
3584ff8c8bSHong Zhang
3684ff8c8bSHong Zhang #include <petscts.h>
3784ff8c8bSHong Zhang #include <petscdm.h>
3884ff8c8bSHong Zhang #include <petscdmda.h>
3984ff8c8bSHong Zhang #include <petscdraw.h>
4084ff8c8bSHong Zhang #include "finitevolume1d.h"
4184ff8c8bSHong Zhang #include <petsc/private/kernels/blockinvert.h>
4284ff8c8bSHong Zhang
RangeMod(PetscReal a,PetscReal xmin,PetscReal xmax)43d71ae5a4SJacob Faibussowitsch static inline PetscReal RangeMod(PetscReal a, PetscReal xmin, PetscReal xmax)
44d71ae5a4SJacob Faibussowitsch {
459371c9d4SSatish Balay PetscReal range = xmax - xmin;
469371c9d4SSatish Balay return xmin + PetscFmodReal(range + PetscFmodReal(a, range), range);
479371c9d4SSatish Balay }
MaxAbs(PetscReal a,PetscReal b)48d71ae5a4SJacob Faibussowitsch static inline PetscReal MaxAbs(PetscReal a, PetscReal b)
49d71ae5a4SJacob Faibussowitsch {
509371c9d4SSatish Balay return (PetscAbs(a) > PetscAbs(b)) ? a : b;
519371c9d4SSatish Balay }
5284ff8c8bSHong Zhang /* --------------------------------- Advection ----------------------------------- */
5384ff8c8bSHong Zhang typedef struct {
5484ff8c8bSHong Zhang PetscReal a; /* advective velocity */
5584ff8c8bSHong Zhang } AdvectCtx;
5684ff8c8bSHong Zhang
PhysicsRiemann_Advect(void * vctx,PetscInt m,const PetscScalar * uL,const PetscScalar * uR,PetscScalar * flux,PetscReal * maxspeed)57d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsRiemann_Advect(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed)
58d71ae5a4SJacob Faibussowitsch {
5984ff8c8bSHong Zhang AdvectCtx *ctx = (AdvectCtx *)vctx;
6084ff8c8bSHong Zhang PetscReal speed;
6184ff8c8bSHong Zhang
6284ff8c8bSHong Zhang PetscFunctionBeginUser;
6384ff8c8bSHong Zhang speed = ctx->a;
6484ff8c8bSHong Zhang flux[0] = PetscMax(0, speed) * uL[0] + PetscMin(0, speed) * uR[0];
6584ff8c8bSHong Zhang *maxspeed = speed;
663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6784ff8c8bSHong Zhang }
6884ff8c8bSHong Zhang
PhysicsCharacteristic_Advect(void * vctx,PetscInt m,const PetscScalar * u,PetscScalar * X,PetscScalar * Xi,PetscReal * speeds)69d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx, PetscInt m, const PetscScalar *u, PetscScalar *X, PetscScalar *Xi, PetscReal *speeds)
70d71ae5a4SJacob Faibussowitsch {
7184ff8c8bSHong Zhang AdvectCtx *ctx = (AdvectCtx *)vctx;
7284ff8c8bSHong Zhang
7384ff8c8bSHong Zhang PetscFunctionBeginUser;
7484ff8c8bSHong Zhang X[0] = 1.;
7584ff8c8bSHong Zhang Xi[0] = 1.;
7684ff8c8bSHong Zhang speeds[0] = ctx->a;
773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7884ff8c8bSHong Zhang }
7984ff8c8bSHong Zhang
PhysicsSample_Advect(void * vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal * u)80d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsSample_Advect(void *vctx, PetscInt initial, FVBCType bctype, PetscReal xmin, PetscReal xmax, PetscReal t, PetscReal x, PetscReal *u)
81d71ae5a4SJacob Faibussowitsch {
8284ff8c8bSHong Zhang AdvectCtx *ctx = (AdvectCtx *)vctx;
8384ff8c8bSHong Zhang PetscReal a = ctx->a, x0;
8484ff8c8bSHong Zhang
8584ff8c8bSHong Zhang PetscFunctionBeginUser;
8684ff8c8bSHong Zhang switch (bctype) {
87d71ae5a4SJacob Faibussowitsch case FVBC_OUTFLOW:
88d71ae5a4SJacob Faibussowitsch x0 = x - a * t;
89d71ae5a4SJacob Faibussowitsch break;
90d71ae5a4SJacob Faibussowitsch case FVBC_PERIODIC:
91d71ae5a4SJacob Faibussowitsch x0 = RangeMod(x - a * t, xmin, xmax);
92d71ae5a4SJacob Faibussowitsch break;
93d71ae5a4SJacob Faibussowitsch default:
94d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown BCType");
9584ff8c8bSHong Zhang }
9684ff8c8bSHong Zhang switch (initial) {
97d71ae5a4SJacob Faibussowitsch case 0:
98d71ae5a4SJacob Faibussowitsch u[0] = (x0 < 0) ? 1 : -1;
99d71ae5a4SJacob Faibussowitsch break;
100d71ae5a4SJacob Faibussowitsch case 1:
101d71ae5a4SJacob Faibussowitsch u[0] = (x0 < 0) ? -1 : 1;
102d71ae5a4SJacob Faibussowitsch break;
103d71ae5a4SJacob Faibussowitsch case 2:
104d71ae5a4SJacob Faibussowitsch u[0] = (0 < x0 && x0 < 1) ? 1 : 0;
105d71ae5a4SJacob Faibussowitsch break;
106d71ae5a4SJacob Faibussowitsch case 3:
107d71ae5a4SJacob Faibussowitsch u[0] = PetscSinReal(2 * PETSC_PI * x0);
108d71ae5a4SJacob Faibussowitsch break;
109d71ae5a4SJacob Faibussowitsch case 4:
110d71ae5a4SJacob Faibussowitsch u[0] = PetscAbs(x0);
111d71ae5a4SJacob Faibussowitsch break;
112d71ae5a4SJacob Faibussowitsch case 5:
113d71ae5a4SJacob Faibussowitsch u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2 * PETSC_PI * x0));
114d71ae5a4SJacob Faibussowitsch break;
115d71ae5a4SJacob Faibussowitsch case 6:
116d71ae5a4SJacob Faibussowitsch u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2 - x0 : 0));
117d71ae5a4SJacob Faibussowitsch break;
118d71ae5a4SJacob Faibussowitsch case 7:
119d71ae5a4SJacob Faibussowitsch u[0] = PetscPowReal(PetscSinReal(PETSC_PI * x0), 10.0);
120d71ae5a4SJacob Faibussowitsch break;
121d71ae5a4SJacob Faibussowitsch default:
122d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition");
12384ff8c8bSHong Zhang }
1243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
12584ff8c8bSHong Zhang }
12684ff8c8bSHong Zhang
PhysicsCreate_Advect(FVCtx * ctx)127d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx)
128d71ae5a4SJacob Faibussowitsch {
12984ff8c8bSHong Zhang AdvectCtx *user;
13084ff8c8bSHong Zhang
13184ff8c8bSHong Zhang PetscFunctionBeginUser;
1329566063dSJacob Faibussowitsch PetscCall(PetscNew(&user));
13384ff8c8bSHong Zhang ctx->physics2.sample2 = PhysicsSample_Advect;
13484ff8c8bSHong Zhang ctx->physics2.riemann2 = PhysicsRiemann_Advect;
13584ff8c8bSHong Zhang ctx->physics2.characteristic2 = PhysicsCharacteristic_Advect;
13684ff8c8bSHong Zhang ctx->physics2.destroy = PhysicsDestroy_SimpleFree;
13784ff8c8bSHong Zhang ctx->physics2.user = user;
13884ff8c8bSHong Zhang ctx->physics2.dof = 1;
13984ff8c8bSHong Zhang
1409566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy("u", &ctx->physics2.fieldname[0]));
14184ff8c8bSHong Zhang user->a = 1;
142d0609cedSBarry Smith PetscOptionsBegin(ctx->comm, ctx->prefix, "Options for advection", "");
1439566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-physics_advect_a", "Speed", "", user->a, &user->a, NULL));
144d0609cedSBarry Smith PetscOptionsEnd();
1453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
14684ff8c8bSHong Zhang }
14784ff8c8bSHong Zhang /* --------------------------------- Shallow Water ----------------------------------- */
14884ff8c8bSHong Zhang
14984ff8c8bSHong Zhang typedef struct {
15084ff8c8bSHong Zhang PetscReal gravity;
15184ff8c8bSHong Zhang } ShallowCtx;
15284ff8c8bSHong Zhang
ShallowFlux(ShallowCtx * phys,const PetscScalar * u,PetscScalar * f)153d71ae5a4SJacob Faibussowitsch static inline void ShallowFlux(ShallowCtx *phys, const PetscScalar *u, PetscScalar *f)
154d71ae5a4SJacob Faibussowitsch {
15584ff8c8bSHong Zhang f[0] = u[1];
15684ff8c8bSHong Zhang f[1] = PetscSqr(u[1]) / u[0] + 0.5 * phys->gravity * PetscSqr(u[0]);
15784ff8c8bSHong Zhang }
15884ff8c8bSHong Zhang
PhysicsRiemann_Shallow_Exact(void * vctx,PetscInt m,const PetscScalar * uL,const PetscScalar * uR,PetscScalar * flux,PetscReal * maxspeed)159d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsRiemann_Shallow_Exact(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed)
160d71ae5a4SJacob Faibussowitsch {
16184ff8c8bSHong Zhang ShallowCtx *phys = (ShallowCtx *)vctx;
16284ff8c8bSHong Zhang PetscScalar g = phys->gravity, ustar[2], cL, cR, c, cstar;
1639371c9d4SSatish Balay struct {
1649371c9d4SSatish Balay PetscScalar h, u;
1659371c9d4SSatish Balay } L = {uL[0], uL[1] / uL[0]}, R = {uR[0], uR[1] / uR[0]}, star;
16684ff8c8bSHong Zhang PetscInt i;
16784ff8c8bSHong Zhang
16884ff8c8bSHong Zhang PetscFunctionBeginUser;
16957508eceSPierre Jolivet PetscCheck(L.h > 0 && R.h > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Reconstructed thickness is negative");
17084ff8c8bSHong Zhang cL = PetscSqrtScalar(g * L.h);
17184ff8c8bSHong Zhang cR = PetscSqrtScalar(g * R.h);
17284ff8c8bSHong Zhang c = PetscMax(cL, cR);
17384ff8c8bSHong Zhang {
17484ff8c8bSHong Zhang /* Solve for star state */
17584ff8c8bSHong Zhang const PetscInt maxits = 50;
17684ff8c8bSHong Zhang PetscScalar tmp, res, res0 = 0, h0, h = 0.5 * (L.h + R.h); /* initial guess */
17784ff8c8bSHong Zhang h0 = h;
17884ff8c8bSHong Zhang for (i = 0; i < maxits; i++) {
17984ff8c8bSHong Zhang PetscScalar fr, fl, dfr, dfl;
1809371c9d4SSatish Balay fl = (L.h < h) ? PetscSqrtScalar(0.5 * g * (h * h - L.h * L.h) * (1 / L.h - 1 / h)) /* shock */
18184ff8c8bSHong Zhang : 2 * PetscSqrtScalar(g * h) - 2 * PetscSqrtScalar(g * L.h); /* rarefaction */
1829371c9d4SSatish Balay fr = (R.h < h) ? PetscSqrtScalar(0.5 * g * (h * h - R.h * R.h) * (1 / R.h - 1 / h)) /* shock */
18384ff8c8bSHong Zhang : 2 * PetscSqrtScalar(g * h) - 2 * PetscSqrtScalar(g * R.h); /* rarefaction */
18484ff8c8bSHong Zhang res = R.u - L.u + fr + fl;
185cad9d221SBarry Smith PetscCheck(!PetscIsInfOrNanScalar(res), PETSC_COMM_SELF, PETSC_ERR_FP, "Infinity or Not-a-Number generated in computation");
18684ff8c8bSHong Zhang if (PetscAbsScalar(res) < PETSC_SQRT_MACHINE_EPSILON || (i > 0 && PetscAbsScalar(h - h0) < PETSC_SQRT_MACHINE_EPSILON)) {
18784ff8c8bSHong Zhang star.h = h;
18884ff8c8bSHong Zhang star.u = L.u - fl;
18984ff8c8bSHong Zhang goto converged;
19084ff8c8bSHong Zhang } else if (i > 0 && PetscAbsScalar(res) >= PetscAbsScalar(res0)) { /* Line search */
19184ff8c8bSHong Zhang h = 0.8 * h0 + 0.2 * h;
19284ff8c8bSHong Zhang continue;
19384ff8c8bSHong Zhang }
19484ff8c8bSHong Zhang /* Accept the last step and take another */
19584ff8c8bSHong Zhang res0 = res;
19684ff8c8bSHong Zhang h0 = h;
19784ff8c8bSHong Zhang dfl = (L.h < h) ? 0.5 / fl * 0.5 * g * (-L.h * L.h / (h * h) - 1 + 2 * h / L.h) : PetscSqrtScalar(g / h);
19884ff8c8bSHong Zhang dfr = (R.h < h) ? 0.5 / fr * 0.5 * g * (-R.h * R.h / (h * h) - 1 + 2 * h / R.h) : PetscSqrtScalar(g / h);
19984ff8c8bSHong Zhang tmp = h - res / (dfr + dfl);
20084ff8c8bSHong Zhang if (tmp <= 0) h /= 2; /* Guard against Newton shooting off to a negative thickness */
20184ff8c8bSHong Zhang else h = tmp;
20257508eceSPierre Jolivet PetscCheck((h > 0) && PetscIsNormalScalar(h), PETSC_COMM_SELF, PETSC_ERR_FP, "non-normal iterate h=%g", (double)h);
20384ff8c8bSHong Zhang }
2043ba16761SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_NOT_CONVERGED, "Newton iteration for star.h diverged after %" PetscInt_FMT " iterations", i);
20584ff8c8bSHong Zhang }
20684ff8c8bSHong Zhang converged:
20784ff8c8bSHong Zhang cstar = PetscSqrtScalar(g * star.h);
20884ff8c8bSHong Zhang if (L.u - cL < 0 && 0 < star.u - cstar) { /* 1-wave is sonic rarefaction */
20984ff8c8bSHong Zhang PetscScalar ufan[2];
21084ff8c8bSHong Zhang ufan[0] = 1 / g * PetscSqr(L.u / 3 + 2. / 3 * cL);
21184ff8c8bSHong Zhang ufan[1] = PetscSqrtScalar(g * ufan[0]) * ufan[0];
21284ff8c8bSHong Zhang ShallowFlux(phys, ufan, flux);
21384ff8c8bSHong Zhang } else if (star.u + cstar < 0 && 0 < R.u + cR) { /* 2-wave is sonic rarefaction */
21484ff8c8bSHong Zhang PetscScalar ufan[2];
21584ff8c8bSHong Zhang ufan[0] = 1 / g * PetscSqr(R.u / 3 - 2. / 3 * cR);
21684ff8c8bSHong Zhang ufan[1] = -PetscSqrtScalar(g * ufan[0]) * ufan[0];
21784ff8c8bSHong Zhang ShallowFlux(phys, ufan, flux);
21884ff8c8bSHong Zhang } else if ((L.h >= star.h && L.u - c >= 0) || (L.h < star.h && (star.h * star.u - L.h * L.u) / (star.h - L.h) > 0)) {
21984ff8c8bSHong Zhang /* 1-wave is right-travelling shock (supersonic) */
22084ff8c8bSHong Zhang ShallowFlux(phys, uL, flux);
22184ff8c8bSHong Zhang } else if ((star.h <= R.h && R.u + c <= 0) || (star.h > R.h && (R.h * R.u - star.h * star.h) / (R.h - star.h) < 0)) {
22284ff8c8bSHong Zhang /* 2-wave is left-travelling shock (supersonic) */
22384ff8c8bSHong Zhang ShallowFlux(phys, uR, flux);
22484ff8c8bSHong Zhang } else {
22584ff8c8bSHong Zhang ustar[0] = star.h;
22684ff8c8bSHong Zhang ustar[1] = star.h * star.u;
22784ff8c8bSHong Zhang ShallowFlux(phys, ustar, flux);
22884ff8c8bSHong Zhang }
22984ff8c8bSHong Zhang *maxspeed = MaxAbs(MaxAbs(star.u - cstar, star.u + cstar), MaxAbs(L.u - cL, R.u + cR));
2303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
23184ff8c8bSHong Zhang }
23284ff8c8bSHong Zhang
PhysicsRiemann_Shallow_Rusanov(void * vctx,PetscInt m,const PetscScalar * uL,const PetscScalar * uR,PetscScalar * flux,PetscReal * maxspeed)233d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsRiemann_Shallow_Rusanov(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed)
234d71ae5a4SJacob Faibussowitsch {
23584ff8c8bSHong Zhang ShallowCtx *phys = (ShallowCtx *)vctx;
23684ff8c8bSHong Zhang PetscScalar g = phys->gravity, fL[2], fR[2], s;
2379371c9d4SSatish Balay struct {
2389371c9d4SSatish Balay PetscScalar h, u;
2399371c9d4SSatish Balay } L = {uL[0], uL[1] / uL[0]}, R = {uR[0], uR[1] / uR[0]};
24084ff8c8bSHong Zhang PetscReal tol = 1e-6;
24184ff8c8bSHong Zhang
24284ff8c8bSHong Zhang PetscFunctionBeginUser;
24384ff8c8bSHong Zhang /* Positivity preserving modification*/
24484ff8c8bSHong Zhang if (L.h < tol) L.u = 0.0;
24584ff8c8bSHong Zhang if (R.h < tol) R.u = 0.0;
24684ff8c8bSHong Zhang
24784ff8c8bSHong Zhang /*simple pos preserve limiter*/
24884ff8c8bSHong Zhang if (L.h < 0) L.h = 0;
24984ff8c8bSHong Zhang if (R.h < 0) R.h = 0;
25084ff8c8bSHong Zhang
25184ff8c8bSHong Zhang ShallowFlux(phys, uL, fL);
25284ff8c8bSHong Zhang ShallowFlux(phys, uR, fR);
25384ff8c8bSHong Zhang
25484ff8c8bSHong Zhang s = PetscMax(PetscAbs(L.u) + PetscSqrtScalar(g * L.h), PetscAbs(R.u) + PetscSqrtScalar(g * R.h));
25584ff8c8bSHong Zhang flux[0] = 0.5 * (fL[0] + fR[0]) + 0.5 * s * (L.h - R.h);
25684ff8c8bSHong Zhang flux[1] = 0.5 * (fL[1] + fR[1]) + 0.5 * s * (uL[1] - uR[1]);
25784ff8c8bSHong Zhang *maxspeed = s;
2583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
25984ff8c8bSHong Zhang }
26084ff8c8bSHong Zhang
PhysicsCharacteristic_Conservative(void * vctx,PetscInt m,const PetscScalar * u,PetscScalar * X,PetscScalar * Xi,PetscReal * speeds)261d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsCharacteristic_Conservative(void *vctx, PetscInt m, const PetscScalar *u, PetscScalar *X, PetscScalar *Xi, PetscReal *speeds)
262d71ae5a4SJacob Faibussowitsch {
26384ff8c8bSHong Zhang PetscInt i, j;
26484ff8c8bSHong Zhang
26584ff8c8bSHong Zhang PetscFunctionBeginUser;
26684ff8c8bSHong Zhang for (i = 0; i < m; i++) {
26784ff8c8bSHong Zhang for (j = 0; j < m; j++) Xi[i * m + j] = X[i * m + j] = (PetscScalar)(i == j);
26884ff8c8bSHong Zhang speeds[i] = PETSC_MAX_REAL; /* Indicates invalid */
26984ff8c8bSHong Zhang }
2703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
27184ff8c8bSHong Zhang }
27284ff8c8bSHong Zhang
PhysicsCharacteristic_Shallow(void * vctx,PetscInt m,const PetscScalar * u,PetscScalar * X,PetscScalar * Xi,PetscReal * speeds)273d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsCharacteristic_Shallow(void *vctx, PetscInt m, const PetscScalar *u, PetscScalar *X, PetscScalar *Xi, PetscReal *speeds)
274d71ae5a4SJacob Faibussowitsch {
27584ff8c8bSHong Zhang ShallowCtx *phys = (ShallowCtx *)vctx;
27684ff8c8bSHong Zhang PetscReal c;
27784ff8c8bSHong Zhang PetscReal tol = 1e-6;
27884ff8c8bSHong Zhang
27984ff8c8bSHong Zhang PetscFunctionBeginUser;
28084ff8c8bSHong Zhang c = PetscSqrtScalar(u[0] * phys->gravity);
28184ff8c8bSHong Zhang
28284ff8c8bSHong Zhang if (u[0] < tol) { /*Use conservative variables*/
28384ff8c8bSHong Zhang X[0 * 2 + 0] = 1;
28484ff8c8bSHong Zhang X[0 * 2 + 1] = 0;
28584ff8c8bSHong Zhang X[1 * 2 + 0] = 0;
28684ff8c8bSHong Zhang X[1 * 2 + 1] = 1;
28784ff8c8bSHong Zhang } else {
28884ff8c8bSHong Zhang speeds[0] = u[1] / u[0] - c;
28984ff8c8bSHong Zhang speeds[1] = u[1] / u[0] + c;
29084ff8c8bSHong Zhang X[0 * 2 + 0] = 1;
29184ff8c8bSHong Zhang X[0 * 2 + 1] = speeds[0];
29284ff8c8bSHong Zhang X[1 * 2 + 0] = 1;
29384ff8c8bSHong Zhang X[1 * 2 + 1] = speeds[1];
29484ff8c8bSHong Zhang }
29584ff8c8bSHong Zhang
2969566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(Xi, X, 4));
2979566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_2(Xi, 0, PETSC_FALSE, NULL));
2983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
29984ff8c8bSHong Zhang }
30084ff8c8bSHong Zhang
PhysicsSample_Shallow(void * vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal * u)301d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsSample_Shallow(void *vctx, PetscInt initial, FVBCType bctype, PetscReal xmin, PetscReal xmax, PetscReal t, PetscReal x, PetscReal *u)
302d71ae5a4SJacob Faibussowitsch {
30384ff8c8bSHong Zhang PetscFunctionBeginUser;
30408401ef6SPierre Jolivet PetscCheck(t <= 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "Exact solutions not implemented for t > 0");
30584ff8c8bSHong Zhang switch (initial) {
30684ff8c8bSHong Zhang case 0:
30784ff8c8bSHong Zhang u[0] = (x < 0) ? 2 : 1; /* Standard Dam Break Problem */
30884ff8c8bSHong Zhang u[1] = (x < 0) ? 0 : 0;
30984ff8c8bSHong Zhang break;
31084ff8c8bSHong Zhang case 1:
31184ff8c8bSHong Zhang u[0] = (x < 10) ? 1 : 0.1; /*The Next 5 problems are standard Riemann problem tests */
31284ff8c8bSHong Zhang u[1] = (x < 10) ? 2.5 : 0;
31384ff8c8bSHong Zhang break;
31484ff8c8bSHong Zhang case 2:
31584ff8c8bSHong Zhang u[0] = (x < 25) ? 1 : 1;
31684ff8c8bSHong Zhang u[1] = (x < 25) ? -5 : 5;
31784ff8c8bSHong Zhang break;
31884ff8c8bSHong Zhang case 3:
31984ff8c8bSHong Zhang u[0] = (x < 20) ? 1 : 1e-12;
32084ff8c8bSHong Zhang u[1] = (x < 20) ? 0 : 0;
32184ff8c8bSHong Zhang break;
32284ff8c8bSHong Zhang case 4:
32384ff8c8bSHong Zhang u[0] = (x < 30) ? 1e-12 : 1;
32484ff8c8bSHong Zhang u[1] = (x < 30) ? 0 : 0;
32584ff8c8bSHong Zhang break;
32684ff8c8bSHong Zhang case 5:
32784ff8c8bSHong Zhang u[0] = (x < 25) ? 0.1 : 0.1;
32884ff8c8bSHong Zhang u[1] = (x < 25) ? -0.3 : 0.3;
32984ff8c8bSHong Zhang break;
33084ff8c8bSHong Zhang case 6:
33184ff8c8bSHong Zhang u[0] = 1 + 0.5 * PetscSinReal(2 * PETSC_PI * x);
33284ff8c8bSHong Zhang u[1] = 1 * u[0];
33384ff8c8bSHong Zhang break;
33484ff8c8bSHong Zhang case 7:
33584ff8c8bSHong Zhang if (x < -0.1) {
33684ff8c8bSHong Zhang u[0] = 1e-9;
33784ff8c8bSHong Zhang u[1] = 0.0;
33884ff8c8bSHong Zhang } else if (x < 0.1) {
33984ff8c8bSHong Zhang u[0] = 1.0;
34084ff8c8bSHong Zhang u[1] = 0.0;
34184ff8c8bSHong Zhang } else {
34284ff8c8bSHong Zhang u[0] = 1e-9;
34384ff8c8bSHong Zhang u[1] = 0.0;
34484ff8c8bSHong Zhang }
34584ff8c8bSHong Zhang break;
34684ff8c8bSHong Zhang case 8:
34784ff8c8bSHong Zhang if (x < -0.1) {
34884ff8c8bSHong Zhang u[0] = 2;
34984ff8c8bSHong Zhang u[1] = 0.0;
35084ff8c8bSHong Zhang } else if (x < 0.1) {
35184ff8c8bSHong Zhang u[0] = 3.0;
35284ff8c8bSHong Zhang u[1] = 0.0;
35384ff8c8bSHong Zhang } else {
35484ff8c8bSHong Zhang u[0] = 2;
35584ff8c8bSHong Zhang u[1] = 0.0;
35684ff8c8bSHong Zhang }
35784ff8c8bSHong Zhang break;
358d71ae5a4SJacob Faibussowitsch default:
359d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition");
36084ff8c8bSHong Zhang }
3613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
36284ff8c8bSHong Zhang }
36384ff8c8bSHong Zhang
36484ff8c8bSHong Zhang /* Implements inflow conditions for the given initial conditions. Which conditions are actually enforced depends on
36584ff8c8bSHong Zhang on the results of PhysicsSetInflowType_Shallow. */
PhysicsInflow_Shallow(void * vctx,PetscReal t,PetscReal x,PetscReal * u)366d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsInflow_Shallow(void *vctx, PetscReal t, PetscReal x, PetscReal *u)
367d71ae5a4SJacob Faibussowitsch {
36884ff8c8bSHong Zhang FVCtx *ctx = (FVCtx *)vctx;
36984ff8c8bSHong Zhang
37084ff8c8bSHong Zhang PetscFunctionBeginUser;
37184ff8c8bSHong Zhang if (ctx->bctype == FVBC_INFLOW) {
37284ff8c8bSHong Zhang switch (ctx->initial) {
37384ff8c8bSHong Zhang case 0:
37484ff8c8bSHong Zhang case 1:
37584ff8c8bSHong Zhang case 2:
37684ff8c8bSHong Zhang case 3:
37784ff8c8bSHong Zhang case 4:
37884ff8c8bSHong Zhang case 5:
3799371c9d4SSatish Balay u[0] = 0;
3809371c9d4SSatish Balay u[1] = 0.0; /* Left boundary conditions */
3819371c9d4SSatish Balay u[2] = 0;
3829371c9d4SSatish Balay u[3] = 0.0; /* Right boundary conditions */
38384ff8c8bSHong Zhang break;
38484ff8c8bSHong Zhang case 6:
3859371c9d4SSatish Balay u[0] = 0;
3869371c9d4SSatish Balay u[1] = 0.0; /* Left boundary conditions */
3879371c9d4SSatish Balay u[2] = 0;
3889371c9d4SSatish Balay u[3] = 0.0; /* Right boundary conditions */
38984ff8c8bSHong Zhang break;
39084ff8c8bSHong Zhang case 7:
3919371c9d4SSatish Balay u[0] = 0;
3929371c9d4SSatish Balay u[1] = 0.0; /* Left boundary conditions */
3939371c9d4SSatish Balay u[2] = 0;
3949371c9d4SSatish Balay u[3] = 0.0; /* Right boundary conditions */
39584ff8c8bSHong Zhang break;
39684ff8c8bSHong Zhang case 8:
3979371c9d4SSatish Balay u[0] = 0;
3989371c9d4SSatish Balay u[1] = 1.0; /* Left boundary conditions */
3999371c9d4SSatish Balay u[2] = 0;
4009371c9d4SSatish Balay u[3] = -1.0; /* Right boundary conditions */
40184ff8c8bSHong Zhang break;
402d71ae5a4SJacob Faibussowitsch default:
403d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition");
40484ff8c8bSHong Zhang }
40584ff8c8bSHong Zhang }
4063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
40784ff8c8bSHong Zhang }
40884ff8c8bSHong Zhang
40984ff8c8bSHong Zhang /* Selects which boundary conditions are marked as inflow and which as outflow when FVBC_INFLOW is selected. */
PhysicsSetInflowType_Shallow(FVCtx * ctx)410d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsSetInflowType_Shallow(FVCtx *ctx)
411d71ae5a4SJacob Faibussowitsch {
41284ff8c8bSHong Zhang PetscFunctionBeginUser;
41384ff8c8bSHong Zhang switch (ctx->initial) {
41484ff8c8bSHong Zhang case 0:
41584ff8c8bSHong Zhang case 1:
41684ff8c8bSHong Zhang case 2:
41784ff8c8bSHong Zhang case 3:
41884ff8c8bSHong Zhang case 4:
41984ff8c8bSHong Zhang case 5:
42084ff8c8bSHong Zhang case 6:
42184ff8c8bSHong Zhang case 7:
42284ff8c8bSHong Zhang case 8: /* Fix left and right momentum, height is outflow */
42384ff8c8bSHong Zhang ctx->physics2.bcinflowindex[0] = PETSC_FALSE;
42484ff8c8bSHong Zhang ctx->physics2.bcinflowindex[1] = PETSC_TRUE;
42584ff8c8bSHong Zhang ctx->physics2.bcinflowindex[2] = PETSC_FALSE;
42684ff8c8bSHong Zhang ctx->physics2.bcinflowindex[3] = PETSC_TRUE;
42784ff8c8bSHong Zhang break;
428d71ae5a4SJacob Faibussowitsch default:
429d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition");
43084ff8c8bSHong Zhang }
4313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
43284ff8c8bSHong Zhang }
43384ff8c8bSHong Zhang
PhysicsCreate_Shallow(FVCtx * ctx)434d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsCreate_Shallow(FVCtx *ctx)
435d71ae5a4SJacob Faibussowitsch {
43684ff8c8bSHong Zhang ShallowCtx *user;
43784ff8c8bSHong Zhang PetscFunctionList rlist = 0, rclist = 0;
43884ff8c8bSHong Zhang char rname[256] = "rusanov", rcname[256] = "characteristic";
43984ff8c8bSHong Zhang
44084ff8c8bSHong Zhang PetscFunctionBeginUser;
4419566063dSJacob Faibussowitsch PetscCall(PetscNew(&user));
44284ff8c8bSHong Zhang ctx->physics2.sample2 = PhysicsSample_Shallow;
44384ff8c8bSHong Zhang ctx->physics2.inflow = PhysicsInflow_Shallow;
44484ff8c8bSHong Zhang ctx->physics2.destroy = PhysicsDestroy_SimpleFree;
44584ff8c8bSHong Zhang ctx->physics2.riemann2 = PhysicsRiemann_Shallow_Rusanov;
44684ff8c8bSHong Zhang ctx->physics2.characteristic2 = PhysicsCharacteristic_Shallow;
44784ff8c8bSHong Zhang ctx->physics2.user = user;
44884ff8c8bSHong Zhang ctx->physics2.dof = 2;
44984ff8c8bSHong Zhang
4503ba16761SJacob Faibussowitsch PetscCall(PetscMalloc1(2 * (ctx->physics2.dof), &ctx->physics2.bcinflowindex));
4513ba16761SJacob Faibussowitsch PetscCall(PhysicsSetInflowType_Shallow(ctx));
45284ff8c8bSHong Zhang
4539566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy("height", &ctx->physics2.fieldname[0]));
4549566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy("momentum", &ctx->physics2.fieldname[1]));
45584ff8c8bSHong Zhang
45684ff8c8bSHong Zhang user->gravity = 9.81;
45784ff8c8bSHong Zhang
4589566063dSJacob Faibussowitsch PetscCall(RiemannListAdd_2WaySplit(&rlist, "exact", PhysicsRiemann_Shallow_Exact));
4599566063dSJacob Faibussowitsch PetscCall(RiemannListAdd_2WaySplit(&rlist, "rusanov", PhysicsRiemann_Shallow_Rusanov));
4609566063dSJacob Faibussowitsch PetscCall(ReconstructListAdd_2WaySplit(&rclist, "characteristic", PhysicsCharacteristic_Shallow));
4619566063dSJacob Faibussowitsch PetscCall(ReconstructListAdd_2WaySplit(&rclist, "conservative", PhysicsCharacteristic_Conservative));
462d0609cedSBarry Smith PetscOptionsBegin(ctx->comm, ctx->prefix, "Options for Shallow", "");
4639566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-physics_shallow_gravity", "Gravity", "", user->gravity, &user->gravity, NULL));
4649566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-physics_shallow_riemann", "Riemann solver", "", rlist, rname, rname, sizeof(rname), NULL));
4659566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-physics_shallow_reconstruct", "Reconstruction", "", rclist, rcname, rcname, sizeof(rcname), NULL));
466d0609cedSBarry Smith PetscOptionsEnd();
4679566063dSJacob Faibussowitsch PetscCall(RiemannListFind_2WaySplit(rlist, rname, &ctx->physics2.riemann2));
4689566063dSJacob Faibussowitsch PetscCall(ReconstructListFind_2WaySplit(rclist, rcname, &ctx->physics2.characteristic2));
4699566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&rlist));
4709566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&rclist));
4713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
47284ff8c8bSHong Zhang }
47384ff8c8bSHong Zhang
FVSample_2WaySplit(FVCtx * ctx,DM da,PetscReal time,Vec U)474d71ae5a4SJacob Faibussowitsch PetscErrorCode FVSample_2WaySplit(FVCtx *ctx, DM da, PetscReal time, Vec U)
475d71ae5a4SJacob Faibussowitsch {
47684ff8c8bSHong Zhang PetscScalar *u, *uj, xj, xi;
47784ff8c8bSHong Zhang PetscInt i, j, k, dof, xs, xm, Mx;
47884ff8c8bSHong Zhang const PetscInt N = 200;
47984ff8c8bSHong Zhang PetscReal hs, hf;
48084ff8c8bSHong Zhang
48184ff8c8bSHong Zhang PetscFunctionBeginUser;
4823c633725SBarry Smith PetscCheck(ctx->physics2.sample2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Physics has not provided a sampling function");
4839566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
4849566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
4859566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u));
4869566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof, &uj));
48784ff8c8bSHong Zhang hs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
48884ff8c8bSHong Zhang hf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
48984ff8c8bSHong Zhang for (i = xs; i < xs + xm; i++) {
49084ff8c8bSHong Zhang if (i < ctx->sf) {
49184ff8c8bSHong Zhang xi = ctx->xmin + 0.5 * hs + i * hs;
49284ff8c8bSHong Zhang /* Integrate over cell i using trapezoid rule with N points. */
49384ff8c8bSHong Zhang for (k = 0; k < dof; k++) u[i * dof + k] = 0;
49484ff8c8bSHong Zhang for (j = 0; j < N + 1; j++) {
49584ff8c8bSHong Zhang xj = xi + hs * (j - N / 2) / (PetscReal)N;
4969566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
49784ff8c8bSHong Zhang for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
49884ff8c8bSHong Zhang }
49984ff8c8bSHong Zhang } else if (i < ctx->fs) {
50084ff8c8bSHong Zhang xi = ctx->xmin + ctx->sf * hs + 0.5 * hf + (i - ctx->sf) * hf;
50184ff8c8bSHong Zhang /* Integrate over cell i using trapezoid rule with N points. */
50284ff8c8bSHong Zhang for (k = 0; k < dof; k++) u[i * dof + k] = 0;
50384ff8c8bSHong Zhang for (j = 0; j < N + 1; j++) {
50484ff8c8bSHong Zhang xj = xi + hf * (j - N / 2) / (PetscReal)N;
5059566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
50684ff8c8bSHong Zhang for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
50784ff8c8bSHong Zhang }
50884ff8c8bSHong Zhang } else {
50984ff8c8bSHong Zhang xi = ctx->xmin + ctx->sf * hs + (ctx->fs - ctx->sf) * hf + 0.5 * hs + (i - ctx->fs) * hs;
51084ff8c8bSHong Zhang /* Integrate over cell i using trapezoid rule with N points. */
51184ff8c8bSHong Zhang for (k = 0; k < dof; k++) u[i * dof + k] = 0;
51284ff8c8bSHong Zhang for (j = 0; j < N + 1; j++) {
51384ff8c8bSHong Zhang xj = xi + hs * (j - N / 2) / (PetscReal)N;
5149566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
51584ff8c8bSHong Zhang for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
51684ff8c8bSHong Zhang }
51784ff8c8bSHong Zhang }
51884ff8c8bSHong Zhang }
5199566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u));
5209566063dSJacob Faibussowitsch PetscCall(PetscFree(uj));
5213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
52284ff8c8bSHong Zhang }
52384ff8c8bSHong Zhang
SolutionErrorNorms_2WaySplit(FVCtx * ctx,DM da,PetscReal t,Vec X,PetscReal * nrm1)524d71ae5a4SJacob Faibussowitsch static PetscErrorCode SolutionErrorNorms_2WaySplit(FVCtx *ctx, DM da, PetscReal t, Vec X, PetscReal *nrm1)
525d71ae5a4SJacob Faibussowitsch {
52684ff8c8bSHong Zhang Vec Y;
52784ff8c8bSHong Zhang PetscInt i, Mx;
52884ff8c8bSHong Zhang const PetscScalar *ptr_X, *ptr_Y;
52984ff8c8bSHong Zhang PetscReal hs, hf;
53084ff8c8bSHong Zhang
53184ff8c8bSHong Zhang PetscFunctionBeginUser;
5329566063dSJacob Faibussowitsch PetscCall(VecGetSize(X, &Mx));
5339566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &Y));
5349566063dSJacob Faibussowitsch PetscCall(FVSample_2WaySplit(ctx, da, t, Y));
53584ff8c8bSHong Zhang hs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
53684ff8c8bSHong Zhang hf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
5379566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &ptr_X));
5389566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &ptr_Y));
53984ff8c8bSHong Zhang for (i = 0; i < Mx; i++) {
54084ff8c8bSHong Zhang if (i < ctx->sf || i > ctx->fs - 1) *nrm1 += hs * PetscAbs(ptr_X[i] - ptr_Y[i]);
54184ff8c8bSHong Zhang else *nrm1 += hf * PetscAbs(ptr_X[i] - ptr_Y[i]);
54284ff8c8bSHong Zhang }
5439566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &ptr_X));
5449566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &ptr_Y));
5459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y));
5463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
54784ff8c8bSHong Zhang }
54884ff8c8bSHong Zhang
FVRHSFunction_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void * vctx)549d71ae5a4SJacob Faibussowitsch PetscErrorCode FVRHSFunction_2WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx)
550d71ae5a4SJacob Faibussowitsch {
55184ff8c8bSHong Zhang FVCtx *ctx = (FVCtx *)vctx;
55284ff8c8bSHong Zhang PetscInt i, j, k, Mx, dof, xs, xm, sf = ctx->sf, fs = ctx->fs;
55384ff8c8bSHong Zhang PetscReal hxf, hxs, cfl_idt = 0;
55484ff8c8bSHong Zhang PetscScalar *x, *f, *slope;
55584ff8c8bSHong Zhang Vec Xloc;
55684ff8c8bSHong Zhang DM da;
55784ff8c8bSHong Zhang
55884ff8c8bSHong Zhang PetscFunctionBeginUser;
5599566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
5609566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc)); /* Xloc contains ghost points */
5619566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); /* Mx is the number of center points */
56284ff8c8bSHong Zhang hxs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
56384ff8c8bSHong Zhang hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
5649566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); /* X is solution vector which does not contain ghost points */
5659566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
56684ff8c8bSHong Zhang
567dd8e379bSPierre Jolivet PetscCall(VecZeroEntries(F)); /* F is the right-hand side function corresponds to center points */
56884ff8c8bSHong Zhang
5699566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x));
5709566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f));
5719566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); /* contains ghost points */
5729566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
57384ff8c8bSHong Zhang
57484ff8c8bSHong Zhang if (ctx->bctype == FVBC_OUTFLOW) {
57584ff8c8bSHong Zhang for (i = xs - 2; i < 0; i++) {
57684ff8c8bSHong Zhang for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
57784ff8c8bSHong Zhang }
57884ff8c8bSHong Zhang for (i = Mx; i < xs + xm + 2; i++) {
57984ff8c8bSHong Zhang for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
58084ff8c8bSHong Zhang }
58184ff8c8bSHong Zhang }
58284ff8c8bSHong Zhang
58384ff8c8bSHong Zhang if (ctx->bctype == FVBC_INFLOW) {
58484ff8c8bSHong Zhang /* See LeVeque, R. (2002). Finite Volume Methods for Hyperbolic Problems. doi:10.1017/CBO9780511791253
58584ff8c8bSHong Zhang pages 137-138 for the scheme. */
58684ff8c8bSHong Zhang if (xs == 0) { /* Left Boundary */
5873ba16761SJacob Faibussowitsch PetscCall(ctx->physics2.inflow(ctx, time, ctx->xmin, ctx->ub));
58884ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
58984ff8c8bSHong Zhang if (ctx->physics2.bcinflowindex[j]) {
59084ff8c8bSHong Zhang for (i = -2; i < 0; i++) x[i * dof + j] = 2.0 * ctx->ub[j] - x[-(i + 1) * dof + j];
59184ff8c8bSHong Zhang } else {
59284ff8c8bSHong Zhang for (i = -2; i < 0; i++) x[i * dof + j] = x[j]; /* Outflow */
59384ff8c8bSHong Zhang }
59484ff8c8bSHong Zhang }
59584ff8c8bSHong Zhang }
59684ff8c8bSHong Zhang if (xs + xm == Mx) { /* Right Boundary */
5973ba16761SJacob Faibussowitsch PetscCall(ctx->physics2.inflow(ctx, time, ctx->xmax, ctx->ub));
59884ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
59984ff8c8bSHong Zhang if (ctx->physics2.bcinflowindex[dof + j]) {
60084ff8c8bSHong Zhang for (i = Mx; i < Mx + 2; i++) x[i * dof + j] = 2.0 * ctx->ub[dof + j] - x[(2 * Mx - (i + 1)) * dof + j];
60184ff8c8bSHong Zhang } else {
60284ff8c8bSHong Zhang for (i = Mx; i < Mx + 2; i++) x[i * dof + j] = x[(Mx - 1) * dof + j]; /* Outflow */
60384ff8c8bSHong Zhang }
60484ff8c8bSHong Zhang }
60584ff8c8bSHong Zhang }
60684ff8c8bSHong Zhang }
60784ff8c8bSHong Zhang
60884ff8c8bSHong Zhang for (i = xs - 1; i < xs + xm + 1; i++) {
60984ff8c8bSHong Zhang struct _LimitInfo info;
61084ff8c8bSHong Zhang PetscScalar *cjmpL, *cjmpR;
61184ff8c8bSHong Zhang /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
6129566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
61384ff8c8bSHong Zhang /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
6149566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
61584ff8c8bSHong Zhang cjmpL = &ctx->cjmpLR[0];
61684ff8c8bSHong Zhang cjmpR = &ctx->cjmpLR[dof];
61784ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
61884ff8c8bSHong Zhang PetscScalar jmpL, jmpR;
61984ff8c8bSHong Zhang jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
62084ff8c8bSHong Zhang jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
62184ff8c8bSHong Zhang for (k = 0; k < dof; k++) {
62284ff8c8bSHong Zhang cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
62384ff8c8bSHong Zhang cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
62484ff8c8bSHong Zhang }
62584ff8c8bSHong Zhang }
62684ff8c8bSHong Zhang /* Apply limiter to the left and right characteristic jumps */
62784ff8c8bSHong Zhang info.m = dof;
62884ff8c8bSHong Zhang info.hxs = hxs;
62984ff8c8bSHong Zhang info.hxf = hxf;
63084ff8c8bSHong Zhang (*ctx->limit2)(&info, cjmpL, cjmpR, ctx->sf, ctx->fs, i, ctx->cslope);
63184ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
63284ff8c8bSHong Zhang PetscScalar tmp = 0;
63384ff8c8bSHong Zhang for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
63484ff8c8bSHong Zhang slope[i * dof + j] = tmp;
63584ff8c8bSHong Zhang }
63684ff8c8bSHong Zhang }
63784ff8c8bSHong Zhang
63884ff8c8bSHong Zhang for (i = xs; i < xs + xm + 1; i++) {
63984ff8c8bSHong Zhang PetscReal maxspeed;
64084ff8c8bSHong Zhang PetscScalar *uL, *uR;
64184ff8c8bSHong Zhang uL = &ctx->uLR[0];
64284ff8c8bSHong Zhang uR = &ctx->uLR[dof];
64384ff8c8bSHong Zhang if (i < sf) { /* slow region */
64484ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
64584ff8c8bSHong Zhang uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
64684ff8c8bSHong Zhang uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
64784ff8c8bSHong Zhang }
6489566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
64984ff8c8bSHong Zhang if (i > xs) {
65084ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxs;
65184ff8c8bSHong Zhang }
65284ff8c8bSHong Zhang if (i < xs + xm) {
65384ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxs;
65484ff8c8bSHong Zhang }
65584ff8c8bSHong Zhang } else if (i == sf) { /* interface between the slow region and the fast region */
65684ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
65784ff8c8bSHong Zhang uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
65884ff8c8bSHong Zhang uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
65984ff8c8bSHong Zhang }
6609566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
66184ff8c8bSHong Zhang if (i > xs) {
66284ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxs;
66384ff8c8bSHong Zhang }
66484ff8c8bSHong Zhang if (i < xs + xm) {
66584ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxf;
66684ff8c8bSHong Zhang }
66784ff8c8bSHong Zhang } else if (i > sf && i < fs) { /* fast region */
66884ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
66984ff8c8bSHong Zhang uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
67084ff8c8bSHong Zhang uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
67184ff8c8bSHong Zhang }
6729566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
67384ff8c8bSHong Zhang if (i > xs) {
67484ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxf;
67584ff8c8bSHong Zhang }
67684ff8c8bSHong Zhang if (i < xs + xm) {
67784ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxf;
67884ff8c8bSHong Zhang }
67984ff8c8bSHong Zhang } else if (i == fs) { /* interface between the fast region and the slow region */
68084ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
68184ff8c8bSHong Zhang uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
68284ff8c8bSHong Zhang uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
68384ff8c8bSHong Zhang }
6849566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
68584ff8c8bSHong Zhang if (i > xs) {
68684ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxf;
68784ff8c8bSHong Zhang }
68884ff8c8bSHong Zhang if (i < xs + xm) {
68984ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxs;
69084ff8c8bSHong Zhang }
69184ff8c8bSHong Zhang } else { /* slow region */
69284ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
69384ff8c8bSHong Zhang uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
69484ff8c8bSHong Zhang uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
69584ff8c8bSHong Zhang }
6969566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
69784ff8c8bSHong Zhang cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hxs)); /* Max allowable value of 1/Delta t */
69884ff8c8bSHong Zhang if (i > xs) {
69984ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxs;
70084ff8c8bSHong Zhang }
70184ff8c8bSHong Zhang if (i < xs + xm) {
70284ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxs;
70384ff8c8bSHong Zhang }
70484ff8c8bSHong Zhang }
70584ff8c8bSHong Zhang }
7069566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
7079566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f));
7089566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
7099566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc));
710462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&cfl_idt, &ctx->cfl_idt, 1, MPIU_SCALAR, MPIU_MAX, PetscObjectComm((PetscObject)da)));
71184ff8c8bSHong Zhang if (0) {
71284ff8c8bSHong Zhang /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */
71384ff8c8bSHong Zhang PetscReal dt, tnow;
7149566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt));
7159566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &tnow));
71684ff8c8bSHong Zhang if (dt > 0.5 / ctx->cfl_idt) {
717ac530a7eSPierre Jolivet if (1) PetscCall(PetscPrintf(ctx->comm, "Stability constraint exceeded at t=%g, dt %g > %g\n", (double)tnow, (double)dt, (double)(0.5 / ctx->cfl_idt)));
718ac530a7eSPierre Jolivet else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Stability constraint exceeded, %g > %g", (double)dt, (double)(ctx->cfl / ctx->cfl_idt));
71984ff8c8bSHong Zhang }
72084ff8c8bSHong Zhang }
7213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
72284ff8c8bSHong Zhang }
72384ff8c8bSHong Zhang
72484ff8c8bSHong Zhang /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */
FVRHSFunctionslow_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void * vctx)725d71ae5a4SJacob Faibussowitsch PetscErrorCode FVRHSFunctionslow_2WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx)
726d71ae5a4SJacob Faibussowitsch {
72784ff8c8bSHong Zhang FVCtx *ctx = (FVCtx *)vctx;
72884ff8c8bSHong Zhang PetscInt i, j, k, Mx, dof, xs, xm, islow = 0, sf = ctx->sf, fs = ctx->fs, lsbwidth = ctx->lsbwidth, rsbwidth = ctx->rsbwidth;
72984ff8c8bSHong Zhang PetscReal hxs, hxf, cfl_idt = 0;
73084ff8c8bSHong Zhang PetscScalar *x, *f, *slope;
73184ff8c8bSHong Zhang Vec Xloc;
73284ff8c8bSHong Zhang DM da;
73384ff8c8bSHong Zhang
73484ff8c8bSHong Zhang PetscFunctionBeginUser;
7359566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
7369566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc));
7379566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
73884ff8c8bSHong Zhang hxs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
73984ff8c8bSHong Zhang hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
7409566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
7419566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
7429566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F));
7439566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x));
7449566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
7459566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
7469566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
74784ff8c8bSHong Zhang
74884ff8c8bSHong Zhang if (ctx->bctype == FVBC_OUTFLOW) {
74984ff8c8bSHong Zhang for (i = xs - 2; i < 0; i++) {
75084ff8c8bSHong Zhang for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
75184ff8c8bSHong Zhang }
75284ff8c8bSHong Zhang for (i = Mx; i < xs + xm + 2; i++) {
75384ff8c8bSHong Zhang for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
75484ff8c8bSHong Zhang }
75584ff8c8bSHong Zhang }
75684ff8c8bSHong Zhang if (ctx->bctype == FVBC_INFLOW) {
75784ff8c8bSHong Zhang /* See LeVeque, R. (2002). Finite Volume Methods for Hyperbolic Problems. doi:10.1017/CBO9780511791253
75884ff8c8bSHong Zhang pages 137-138 for the scheme. */
75984ff8c8bSHong Zhang if (xs == 0) { /* Left Boundary */
7603ba16761SJacob Faibussowitsch PetscCall(ctx->physics2.inflow(ctx, time, ctx->xmin, ctx->ub));
76184ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
76284ff8c8bSHong Zhang if (ctx->physics2.bcinflowindex[j] == PETSC_TRUE) {
76384ff8c8bSHong Zhang for (i = -2; i < 0; i++) x[i * dof + j] = 2.0 * ctx->ub[j] - x[-(i + 1) * dof + j];
76484ff8c8bSHong Zhang } else {
76584ff8c8bSHong Zhang for (i = -2; i < 0; i++) x[i * dof + j] = x[j]; /* Outflow */
76684ff8c8bSHong Zhang }
76784ff8c8bSHong Zhang }
76884ff8c8bSHong Zhang }
76984ff8c8bSHong Zhang if (xs + xm == Mx) { /* Right Boundary */
7703ba16761SJacob Faibussowitsch PetscCall(ctx->physics2.inflow(ctx, time, ctx->xmax, ctx->ub));
77184ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
77284ff8c8bSHong Zhang if (ctx->physics2.bcinflowindex[dof + j] == PETSC_TRUE) {
77384ff8c8bSHong Zhang for (i = Mx; i < Mx + 2; i++) x[i * dof + j] = 2.0 * ctx->ub[dof + j] - x[(2 * Mx - (i + 1)) * dof + j];
77484ff8c8bSHong Zhang } else {
77584ff8c8bSHong Zhang for (i = Mx; i < Mx + 2; i++) x[i * dof + j] = x[(Mx - 1) * dof + j]; /* Outflow */
77684ff8c8bSHong Zhang }
77784ff8c8bSHong Zhang }
77884ff8c8bSHong Zhang }
77984ff8c8bSHong Zhang }
78084ff8c8bSHong Zhang for (i = xs - 1; i < xs + xm + 1; i++) {
78184ff8c8bSHong Zhang struct _LimitInfo info;
78284ff8c8bSHong Zhang PetscScalar *cjmpL, *cjmpR;
78384ff8c8bSHong Zhang if (i < sf - lsbwidth + 1 || i > fs + rsbwidth - 2) { /* slow components and the first and last fast components */
78484ff8c8bSHong Zhang /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
7859566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
78684ff8c8bSHong Zhang /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
7879566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
78884ff8c8bSHong Zhang cjmpL = &ctx->cjmpLR[0];
78984ff8c8bSHong Zhang cjmpR = &ctx->cjmpLR[dof];
79084ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
79184ff8c8bSHong Zhang PetscScalar jmpL, jmpR;
79284ff8c8bSHong Zhang jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
79384ff8c8bSHong Zhang jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
79484ff8c8bSHong Zhang for (k = 0; k < dof; k++) {
79584ff8c8bSHong Zhang cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
79684ff8c8bSHong Zhang cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
79784ff8c8bSHong Zhang }
79884ff8c8bSHong Zhang }
79984ff8c8bSHong Zhang /* Apply limiter to the left and right characteristic jumps */
80084ff8c8bSHong Zhang info.m = dof;
80184ff8c8bSHong Zhang info.hxs = hxs;
80284ff8c8bSHong Zhang info.hxf = hxf;
80384ff8c8bSHong Zhang (*ctx->limit2)(&info, cjmpL, cjmpR, ctx->sf, ctx->fs, i, ctx->cslope);
80484ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
80584ff8c8bSHong Zhang PetscScalar tmp = 0;
80684ff8c8bSHong Zhang for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
80784ff8c8bSHong Zhang slope[i * dof + j] = tmp;
80884ff8c8bSHong Zhang }
80984ff8c8bSHong Zhang }
81084ff8c8bSHong Zhang }
81184ff8c8bSHong Zhang
81284ff8c8bSHong Zhang for (i = xs; i < xs + xm + 1; i++) {
81384ff8c8bSHong Zhang PetscReal maxspeed;
81484ff8c8bSHong Zhang PetscScalar *uL, *uR;
81584ff8c8bSHong Zhang uL = &ctx->uLR[0];
81684ff8c8bSHong Zhang uR = &ctx->uLR[dof];
81784ff8c8bSHong Zhang if (i < sf - lsbwidth) { /* slow region */
81884ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
81984ff8c8bSHong Zhang uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
82084ff8c8bSHong Zhang uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
82184ff8c8bSHong Zhang }
8229566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
82384ff8c8bSHong Zhang cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hxs)); /* Max allowable value of 1/Delta t */
82484ff8c8bSHong Zhang if (i > xs) {
82584ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
82684ff8c8bSHong Zhang }
82784ff8c8bSHong Zhang if (i < xs + xm) {
82884ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
82984ff8c8bSHong Zhang islow++;
83084ff8c8bSHong Zhang }
83184ff8c8bSHong Zhang }
83284ff8c8bSHong Zhang if (i == sf - lsbwidth) { /* interface between the slow region and the fast region */
83384ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
83484ff8c8bSHong Zhang uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
83584ff8c8bSHong Zhang uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
83684ff8c8bSHong Zhang }
8379566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
83884ff8c8bSHong Zhang if (i > xs) {
83984ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
84084ff8c8bSHong Zhang }
84184ff8c8bSHong Zhang }
84284ff8c8bSHong Zhang if (i == fs + rsbwidth) { /* slow region */
84384ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
84484ff8c8bSHong Zhang uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
84584ff8c8bSHong Zhang uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
84684ff8c8bSHong Zhang }
8479566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
84884ff8c8bSHong Zhang if (i < xs + xm) {
84984ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
85084ff8c8bSHong Zhang islow++;
85184ff8c8bSHong Zhang }
85284ff8c8bSHong Zhang }
85384ff8c8bSHong Zhang if (i > fs + rsbwidth) { /* slow region */
85484ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
85584ff8c8bSHong Zhang uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
85684ff8c8bSHong Zhang uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
85784ff8c8bSHong Zhang }
8589566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
85984ff8c8bSHong Zhang if (i > xs) {
86084ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
86184ff8c8bSHong Zhang }
86284ff8c8bSHong Zhang if (i < xs + xm) {
86384ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
86484ff8c8bSHong Zhang islow++;
86584ff8c8bSHong Zhang }
86684ff8c8bSHong Zhang }
86784ff8c8bSHong Zhang }
8689566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
8699566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
8709566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
8719566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc));
872462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&cfl_idt, &ctx->cfl_idt, 1, MPIU_SCALAR, MPIU_MAX, PetscObjectComm((PetscObject)da)));
8733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
87484ff8c8bSHong Zhang }
87584ff8c8bSHong Zhang
FVRHSFunctionslowbuffer_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void * vctx)876d71ae5a4SJacob Faibussowitsch PetscErrorCode FVRHSFunctionslowbuffer_2WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx)
877d71ae5a4SJacob Faibussowitsch {
87884ff8c8bSHong Zhang FVCtx *ctx = (FVCtx *)vctx;
87984ff8c8bSHong Zhang PetscInt i, j, k, Mx, dof, xs, xm, islow = 0, sf = ctx->sf, fs = ctx->fs, lsbwidth = ctx->lsbwidth, rsbwidth = ctx->rsbwidth;
88084ff8c8bSHong Zhang PetscReal hxs, hxf;
88184ff8c8bSHong Zhang PetscScalar *x, *f, *slope;
88284ff8c8bSHong Zhang Vec Xloc;
88384ff8c8bSHong Zhang DM da;
88484ff8c8bSHong Zhang
88584ff8c8bSHong Zhang PetscFunctionBeginUser;
8869566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
8879566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc));
8889566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
88984ff8c8bSHong Zhang hxs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
89084ff8c8bSHong Zhang hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
8919566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
8929566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
8939566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F));
8949566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x));
8959566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
8969566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
8979566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
89884ff8c8bSHong Zhang
89984ff8c8bSHong Zhang if (ctx->bctype == FVBC_OUTFLOW) {
90084ff8c8bSHong Zhang for (i = xs - 2; i < 0; i++) {
90184ff8c8bSHong Zhang for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
90284ff8c8bSHong Zhang }
90384ff8c8bSHong Zhang for (i = Mx; i < xs + xm + 2; i++) {
90484ff8c8bSHong Zhang for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
90584ff8c8bSHong Zhang }
90684ff8c8bSHong Zhang }
90784ff8c8bSHong Zhang if (ctx->bctype == FVBC_INFLOW) {
90884ff8c8bSHong Zhang /* See LeVeque, R. (2002). Finite Volume Methods for Hyperbolic Problems. doi:10.1017/CBO9780511791253
90984ff8c8bSHong Zhang pages 137-138 for the scheme. */
91084ff8c8bSHong Zhang if (xs == 0) { /* Left Boundary */
9113ba16761SJacob Faibussowitsch PetscCall(ctx->physics2.inflow(ctx, time, ctx->xmin, ctx->ub));
91284ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
91384ff8c8bSHong Zhang if (ctx->physics2.bcinflowindex[j] == PETSC_TRUE) {
91484ff8c8bSHong Zhang for (i = -2; i < 0; i++) x[i * dof + j] = 2.0 * ctx->ub[j] - x[-(i + 1) * dof + j];
91584ff8c8bSHong Zhang } else {
91684ff8c8bSHong Zhang for (i = -2; i < 0; i++) x[i * dof + j] = x[j]; /* Outflow */
91784ff8c8bSHong Zhang }
91884ff8c8bSHong Zhang }
91984ff8c8bSHong Zhang }
92084ff8c8bSHong Zhang if (xs + xm == Mx) { /* Right Boundary */
9213ba16761SJacob Faibussowitsch PetscCall(ctx->physics2.inflow(ctx, time, ctx->xmax, ctx->ub));
92284ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
92384ff8c8bSHong Zhang if (ctx->physics2.bcinflowindex[dof + j] == PETSC_TRUE) {
92484ff8c8bSHong Zhang for (i = Mx; i < Mx + 2; i++) x[i * dof + j] = 2.0 * ctx->ub[dof + j] - x[(2 * Mx - (i + 1)) * dof + j];
92584ff8c8bSHong Zhang } else {
92684ff8c8bSHong Zhang for (i = Mx; i < Mx + 2; i++) x[i * dof + j] = x[(Mx - 1) * dof + j]; /* Outflow */
92784ff8c8bSHong Zhang }
92884ff8c8bSHong Zhang }
92984ff8c8bSHong Zhang }
93084ff8c8bSHong Zhang }
93184ff8c8bSHong Zhang for (i = xs - 1; i < xs + xm + 1; i++) {
93284ff8c8bSHong Zhang struct _LimitInfo info;
93384ff8c8bSHong Zhang PetscScalar *cjmpL, *cjmpR;
93484ff8c8bSHong Zhang if ((i > sf - lsbwidth - 2 && i < sf + 1) || (i > fs - 2 && i < fs + rsbwidth + 1)) {
93584ff8c8bSHong Zhang /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
9369566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
93784ff8c8bSHong Zhang /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
9389566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
93984ff8c8bSHong Zhang cjmpL = &ctx->cjmpLR[0];
94084ff8c8bSHong Zhang cjmpR = &ctx->cjmpLR[dof];
94184ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
94284ff8c8bSHong Zhang PetscScalar jmpL, jmpR;
94384ff8c8bSHong Zhang jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
94484ff8c8bSHong Zhang jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
94584ff8c8bSHong Zhang for (k = 0; k < dof; k++) {
94684ff8c8bSHong Zhang cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
94784ff8c8bSHong Zhang cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
94884ff8c8bSHong Zhang }
94984ff8c8bSHong Zhang }
95084ff8c8bSHong Zhang /* Apply limiter to the left and right characteristic jumps */
95184ff8c8bSHong Zhang info.m = dof;
95284ff8c8bSHong Zhang info.hxs = hxs;
95384ff8c8bSHong Zhang info.hxf = hxf;
95484ff8c8bSHong Zhang (*ctx->limit2)(&info, cjmpL, cjmpR, ctx->sf, ctx->fs, i, ctx->cslope);
95584ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
95684ff8c8bSHong Zhang PetscScalar tmp = 0;
95784ff8c8bSHong Zhang for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
95884ff8c8bSHong Zhang slope[i * dof + j] = tmp;
95984ff8c8bSHong Zhang }
96084ff8c8bSHong Zhang }
96184ff8c8bSHong Zhang }
96284ff8c8bSHong Zhang
96384ff8c8bSHong Zhang for (i = xs; i < xs + xm + 1; i++) {
96484ff8c8bSHong Zhang PetscReal maxspeed;
96584ff8c8bSHong Zhang PetscScalar *uL, *uR;
96684ff8c8bSHong Zhang uL = &ctx->uLR[0];
96784ff8c8bSHong Zhang uR = &ctx->uLR[dof];
96884ff8c8bSHong Zhang if (i == sf - lsbwidth) {
96984ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
97084ff8c8bSHong Zhang uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
97184ff8c8bSHong Zhang uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
97284ff8c8bSHong Zhang }
9739566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
97484ff8c8bSHong Zhang if (i < xs + xm) {
97584ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
97684ff8c8bSHong Zhang islow++;
97784ff8c8bSHong Zhang }
97884ff8c8bSHong Zhang }
97984ff8c8bSHong Zhang if (i > sf - lsbwidth && i < sf) {
98084ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
98184ff8c8bSHong Zhang uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
98284ff8c8bSHong Zhang uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
98384ff8c8bSHong Zhang }
9849566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
98584ff8c8bSHong Zhang if (i > xs) {
98684ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
98784ff8c8bSHong Zhang }
98884ff8c8bSHong Zhang if (i < xs + xm) {
98984ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
99084ff8c8bSHong Zhang islow++;
99184ff8c8bSHong Zhang }
99284ff8c8bSHong Zhang }
99384ff8c8bSHong Zhang if (i == sf) { /* interface between the slow region and the fast region */
99484ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
99584ff8c8bSHong Zhang uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
99684ff8c8bSHong Zhang uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
99784ff8c8bSHong Zhang }
9989566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
99984ff8c8bSHong Zhang if (i > xs) {
100084ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
100184ff8c8bSHong Zhang }
100284ff8c8bSHong Zhang }
100384ff8c8bSHong Zhang if (i == fs) { /* interface between the fast region and the slow region */
100484ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
100584ff8c8bSHong Zhang uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
100684ff8c8bSHong Zhang uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
100784ff8c8bSHong Zhang }
10089566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
100984ff8c8bSHong Zhang if (i < xs + xm) {
101084ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
101184ff8c8bSHong Zhang islow++;
101284ff8c8bSHong Zhang }
101384ff8c8bSHong Zhang }
101484ff8c8bSHong Zhang if (i > fs && i < fs + rsbwidth) {
101584ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
101684ff8c8bSHong Zhang uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
101784ff8c8bSHong Zhang uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
101884ff8c8bSHong Zhang }
10199566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
102084ff8c8bSHong Zhang if (i > xs) {
102184ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
102284ff8c8bSHong Zhang }
102384ff8c8bSHong Zhang if (i < xs + xm) {
102484ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
102584ff8c8bSHong Zhang islow++;
102684ff8c8bSHong Zhang }
102784ff8c8bSHong Zhang }
102884ff8c8bSHong Zhang if (i == fs + rsbwidth) {
102984ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
103084ff8c8bSHong Zhang uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
103184ff8c8bSHong Zhang uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
103284ff8c8bSHong Zhang }
10339566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
103484ff8c8bSHong Zhang if (i > xs) {
103584ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
103684ff8c8bSHong Zhang }
103784ff8c8bSHong Zhang }
103884ff8c8bSHong Zhang }
10399566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
10409566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
10419566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
10429566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc));
10433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
104484ff8c8bSHong Zhang }
104584ff8c8bSHong Zhang
104684ff8c8bSHong Zhang /* --------------------------------- Finite Volume Solver for fast parts ----------------------------------- */
FVRHSFunctionfast_2WaySplit(TS ts,PetscReal time,Vec X,Vec F,void * vctx)1047d71ae5a4SJacob Faibussowitsch PetscErrorCode FVRHSFunctionfast_2WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx)
1048d71ae5a4SJacob Faibussowitsch {
104984ff8c8bSHong Zhang FVCtx *ctx = (FVCtx *)vctx;
105084ff8c8bSHong Zhang PetscInt i, j, k, Mx, dof, xs, xm, ifast = 0, sf = ctx->sf, fs = ctx->fs;
105184ff8c8bSHong Zhang PetscReal hxs, hxf;
105284ff8c8bSHong Zhang PetscScalar *x, *f, *slope;
105384ff8c8bSHong Zhang Vec Xloc;
105484ff8c8bSHong Zhang DM da;
105584ff8c8bSHong Zhang
105684ff8c8bSHong Zhang PetscFunctionBeginUser;
10579566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
10589566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc));
10599566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
106084ff8c8bSHong Zhang hxs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
106184ff8c8bSHong Zhang hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
10629566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
10639566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
10649566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F));
10659566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x));
10669566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
10679566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
10689566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
106984ff8c8bSHong Zhang
107084ff8c8bSHong Zhang if (ctx->bctype == FVBC_OUTFLOW) {
107184ff8c8bSHong Zhang for (i = xs - 2; i < 0; i++) {
107284ff8c8bSHong Zhang for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
107384ff8c8bSHong Zhang }
107484ff8c8bSHong Zhang for (i = Mx; i < xs + xm + 2; i++) {
107584ff8c8bSHong Zhang for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
107684ff8c8bSHong Zhang }
107784ff8c8bSHong Zhang }
107884ff8c8bSHong Zhang if (ctx->bctype == FVBC_INFLOW) {
107984ff8c8bSHong Zhang /* See LeVeque, R. (2002). Finite Volume Methods for Hyperbolic Problems. doi:10.1017/CBO9780511791253
108084ff8c8bSHong Zhang pages 137-138 for the scheme. */
108184ff8c8bSHong Zhang if (xs == 0) { /* Left Boundary */
10823ba16761SJacob Faibussowitsch PetscCall(ctx->physics2.inflow(ctx, time, ctx->xmin, ctx->ub));
108384ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
108484ff8c8bSHong Zhang if (ctx->physics2.bcinflowindex[j] == PETSC_TRUE) {
108584ff8c8bSHong Zhang for (i = -2; i < 0; i++) x[i * dof + j] = 2.0 * ctx->ub[j] - x[-(i + 1) * dof + j];
108684ff8c8bSHong Zhang } else {
108784ff8c8bSHong Zhang for (i = -2; i < 0; i++) x[i * dof + j] = x[j]; /* Outflow */
108884ff8c8bSHong Zhang }
108984ff8c8bSHong Zhang }
109084ff8c8bSHong Zhang }
109184ff8c8bSHong Zhang if (xs + xm == Mx) { /* Right Boundary */
10923ba16761SJacob Faibussowitsch PetscCall(ctx->physics2.inflow(ctx, time, ctx->xmax, ctx->ub));
109384ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
109484ff8c8bSHong Zhang if (ctx->physics2.bcinflowindex[dof + j] == PETSC_TRUE) {
109584ff8c8bSHong Zhang for (i = Mx; i < Mx + 2; i++) x[i * dof + j] = 2.0 * ctx->ub[dof + j] - x[(2 * Mx - (i + 1)) * dof + j];
109684ff8c8bSHong Zhang } else {
109784ff8c8bSHong Zhang for (i = Mx; i < Mx + 2; i++) x[i * dof + j] = x[(Mx - 1) * dof + j]; /* Outflow */
109884ff8c8bSHong Zhang }
109984ff8c8bSHong Zhang }
110084ff8c8bSHong Zhang }
110184ff8c8bSHong Zhang }
110284ff8c8bSHong Zhang for (i = xs - 1; i < xs + xm + 1; i++) {
110384ff8c8bSHong Zhang struct _LimitInfo info;
110484ff8c8bSHong Zhang PetscScalar *cjmpL, *cjmpR;
110584ff8c8bSHong Zhang if (i > sf - 2 && i < fs + 1) {
11069566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
11079566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
110884ff8c8bSHong Zhang cjmpL = &ctx->cjmpLR[0];
110984ff8c8bSHong Zhang cjmpR = &ctx->cjmpLR[dof];
111084ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
111184ff8c8bSHong Zhang PetscScalar jmpL, jmpR;
111284ff8c8bSHong Zhang jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
111384ff8c8bSHong Zhang jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
111484ff8c8bSHong Zhang for (k = 0; k < dof; k++) {
111584ff8c8bSHong Zhang cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
111684ff8c8bSHong Zhang cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
111784ff8c8bSHong Zhang }
111884ff8c8bSHong Zhang }
111984ff8c8bSHong Zhang /* Apply limiter to the left and right characteristic jumps */
112084ff8c8bSHong Zhang info.m = dof;
112184ff8c8bSHong Zhang info.hxs = hxs;
112284ff8c8bSHong Zhang info.hxf = hxf;
112384ff8c8bSHong Zhang (*ctx->limit2)(&info, cjmpL, cjmpR, ctx->sf, ctx->fs, i, ctx->cslope);
112484ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
112584ff8c8bSHong Zhang PetscScalar tmp = 0;
112684ff8c8bSHong Zhang for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
112784ff8c8bSHong Zhang slope[i * dof + j] = tmp;
112884ff8c8bSHong Zhang }
112984ff8c8bSHong Zhang }
113084ff8c8bSHong Zhang }
113184ff8c8bSHong Zhang
113284ff8c8bSHong Zhang for (i = xs; i < xs + xm + 1; i++) {
113384ff8c8bSHong Zhang PetscReal maxspeed;
113484ff8c8bSHong Zhang PetscScalar *uL, *uR;
113584ff8c8bSHong Zhang uL = &ctx->uLR[0];
113684ff8c8bSHong Zhang uR = &ctx->uLR[dof];
113784ff8c8bSHong Zhang if (i == sf) { /* interface between the slow region and the fast region */
113884ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
113984ff8c8bSHong Zhang uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
114084ff8c8bSHong Zhang uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
114184ff8c8bSHong Zhang }
11429566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
114384ff8c8bSHong Zhang if (i < xs + xm) {
114484ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[ifast * dof + j] += ctx->flux[j] / hxf;
114584ff8c8bSHong Zhang ifast++;
114684ff8c8bSHong Zhang }
114784ff8c8bSHong Zhang }
114884ff8c8bSHong Zhang if (i > sf && i < fs) { /* fast region */
114984ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
115084ff8c8bSHong Zhang uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
115184ff8c8bSHong Zhang uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
115284ff8c8bSHong Zhang }
11539566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
115484ff8c8bSHong Zhang if (i > xs) {
115584ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[(ifast - 1) * dof + j] -= ctx->flux[j] / hxf;
115684ff8c8bSHong Zhang }
115784ff8c8bSHong Zhang if (i < xs + xm) {
115884ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[ifast * dof + j] += ctx->flux[j] / hxf;
115984ff8c8bSHong Zhang ifast++;
116084ff8c8bSHong Zhang }
116184ff8c8bSHong Zhang }
116284ff8c8bSHong Zhang if (i == fs) { /* interface between the fast region and the slow region */
116384ff8c8bSHong Zhang for (j = 0; j < dof; j++) {
116484ff8c8bSHong Zhang uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
116584ff8c8bSHong Zhang uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
116684ff8c8bSHong Zhang }
11679566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
116884ff8c8bSHong Zhang if (i > xs) {
116984ff8c8bSHong Zhang for (j = 0; j < dof; j++) f[(ifast - 1) * dof + j] -= ctx->flux[j] / hxf;
117084ff8c8bSHong Zhang }
117184ff8c8bSHong Zhang }
117284ff8c8bSHong Zhang }
11739566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
11749566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
11759566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
11769566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc));
11773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
117884ff8c8bSHong Zhang }
117984ff8c8bSHong Zhang
main(int argc,char * argv[])1180d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[])
1181d71ae5a4SJacob Faibussowitsch {
118284ff8c8bSHong Zhang char lname[256] = "mc", physname[256] = "advect", final_fname[256] = "solution.m";
118384ff8c8bSHong Zhang PetscFunctionList limiters = 0, physics = 0;
118484ff8c8bSHong Zhang MPI_Comm comm;
118584ff8c8bSHong Zhang TS ts;
118684ff8c8bSHong Zhang DM da;
118784ff8c8bSHong Zhang Vec X, X0, R;
118884ff8c8bSHong Zhang FVCtx ctx;
118984ff8c8bSHong Zhang PetscInt i, k, dof, xs, xm, Mx, draw = 0, count_slow, count_fast, islow = 0, ifast = 0, islowbuffer = 0, *index_slow, *index_fast, *index_slowbuffer;
119084ff8c8bSHong Zhang PetscBool view_final = PETSC_FALSE;
119184ff8c8bSHong Zhang PetscReal ptime, maxtime;
119284ff8c8bSHong Zhang
1193327415f7SBarry Smith PetscFunctionBeginUser;
11949566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, 0, help));
119584ff8c8bSHong Zhang comm = PETSC_COMM_WORLD;
11969566063dSJacob Faibussowitsch PetscCall(PetscMemzero(&ctx, sizeof(ctx)));
119784ff8c8bSHong Zhang
119884ff8c8bSHong Zhang /* Register limiters to be available on the command line */
11999566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "upwind", Limit2_Upwind));
12009566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "lax-wendroff", Limit2_LaxWendroff));
12019566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "beam-warming", Limit2_BeamWarming));
12029566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "fromm", Limit2_Fromm));
12039566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "minmod", Limit2_Minmod));
12049566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "superbee", Limit2_Superbee));
12059566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "mc", Limit2_MC));
12069566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "koren3", Limit2_Koren3));
120784ff8c8bSHong Zhang
120884ff8c8bSHong Zhang /* Register physical models to be available on the command line */
12099566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&physics, "shallow", PhysicsCreate_Shallow));
12109566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&physics, "advect", PhysicsCreate_Advect));
121184ff8c8bSHong Zhang
121284ff8c8bSHong Zhang ctx.comm = comm;
121384ff8c8bSHong Zhang ctx.cfl = 0.9;
121484ff8c8bSHong Zhang ctx.bctype = FVBC_PERIODIC;
121584ff8c8bSHong Zhang ctx.xmin = -1.0;
121684ff8c8bSHong Zhang ctx.xmax = 1.0;
121784ff8c8bSHong Zhang ctx.initial = 1;
121884ff8c8bSHong Zhang ctx.hratio = 2;
121984ff8c8bSHong Zhang maxtime = 10.0;
122084ff8c8bSHong Zhang ctx.simulation = PETSC_FALSE;
1221d0609cedSBarry Smith PetscOptionsBegin(comm, NULL, "Finite Volume solver options", "");
12229566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-xmin", "X min", "", ctx.xmin, &ctx.xmin, NULL));
12239566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-xmax", "X max", "", ctx.xmax, &ctx.xmax, NULL));
12249566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-limit", "Name of flux imiter to use", "", limiters, lname, lname, sizeof(lname), NULL));
12259566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-physics", "Name of physics model to use", "", physics, physname, physname, sizeof(physname), NULL));
12269566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-draw", "Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)", "", draw, &draw, NULL));
12279566063dSJacob 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));
12289566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-initial", "Initial condition (depends on the physics)", "", ctx.initial, &ctx.initial, NULL));
12299566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-exact", "Compare errors with exact solution", "", ctx.exact, &ctx.exact, NULL));
12309566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-simulation", "Compare errors with reference solution", "", ctx.simulation, &ctx.simulation, NULL));
12319566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-cfl", "CFL number to time step at", "", ctx.cfl, &ctx.cfl, NULL));
12329566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-bc_type", "Boundary condition", "", FVBCTypes, (PetscEnum)ctx.bctype, (PetscEnum *)&ctx.bctype, NULL));
12339566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-hratio", "Spacing ratio", "", ctx.hratio, &ctx.hratio, NULL));
1234d0609cedSBarry Smith PetscOptionsEnd();
123584ff8c8bSHong Zhang
123684ff8c8bSHong Zhang /* Choose the limiter from the list of registered limiters */
12379566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(limiters, lname, &ctx.limit2));
12383c633725SBarry Smith PetscCheck(ctx.limit2, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Limiter '%s' not found", lname);
123984ff8c8bSHong Zhang
124084ff8c8bSHong Zhang /* Choose the physics from the list of registered models */
124184ff8c8bSHong Zhang {
124284ff8c8bSHong Zhang PetscErrorCode (*r)(FVCtx *);
12439566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(physics, physname, &r));
12443c633725SBarry Smith PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Physics '%s' not found", physname);
124584ff8c8bSHong Zhang /* Create the physics, will set the number of fields and their names */
12469566063dSJacob Faibussowitsch PetscCall((*r)(&ctx));
124784ff8c8bSHong Zhang }
124884ff8c8bSHong Zhang
124984ff8c8bSHong Zhang /* Create a DMDA to manage the parallel grid */
12509566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(comm, DM_BOUNDARY_PERIODIC, 50, ctx.physics2.dof, 2, NULL, &da));
12519566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da));
12529566063dSJacob Faibussowitsch PetscCall(DMSetUp(da));
125384ff8c8bSHong Zhang /* Inform the DMDA of the field names provided by the physics. */
125484ff8c8bSHong Zhang /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */
125548a46eb9SPierre Jolivet for (i = 0; i < ctx.physics2.dof; i++) PetscCall(DMDASetFieldName(da, i, ctx.physics2.fieldname[i]));
12569566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
12579566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
125884ff8c8bSHong Zhang
125984ff8c8bSHong Zhang /* Set coordinates of cell centers */
12609566063dSJacob 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));
126184ff8c8bSHong Zhang
126284ff8c8bSHong Zhang /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */
12639566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(dof * dof, &ctx.R, dof * dof, &ctx.Rinv, 2 * dof, &ctx.cjmpLR, 1 * dof, &ctx.cslope));
12649566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(2 * dof, &ctx.uLR, dof, &ctx.flux, dof, &ctx.speeds));
12659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * dof, &ctx.ub));
126684ff8c8bSHong Zhang
126784ff8c8bSHong Zhang /* Create a vector to store the solution and to save the initial state */
12689566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &X));
12699566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &X0));
12709566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &R));
127184ff8c8bSHong Zhang
127284ff8c8bSHong Zhang /* create index for slow parts and fast parts,
127384ff8c8bSHong Zhang count_slow + count_fast = Mx, counts_slow*hs = 0.5, counts_fast*hf = 0.5 */
12749fd2a872SSatish Balay count_slow = Mx * 3 / (3 + ctx.hratio); // compute Mx / (1.0 + ctx.hratio / 3.0);
127508401ef6SPierre Jolivet PetscCheck(count_slow % 2 == 0, PETSC_COMM_WORLD, PETSC_ERR_USER, "Please adjust grid size Mx (-da_grid_x) and hratio (-hratio) so that Mx/(1+hratio/3) is even");
127684ff8c8bSHong Zhang count_fast = Mx - count_slow;
127784ff8c8bSHong Zhang ctx.sf = count_slow / 2;
127884ff8c8bSHong Zhang ctx.fs = ctx.sf + count_fast;
12799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(xm * dof, &index_slow));
12809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(xm * dof, &index_fast));
12819566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(8 * dof, &index_slowbuffer));
128284ff8c8bSHong Zhang ctx.lsbwidth = 4;
128384ff8c8bSHong Zhang ctx.rsbwidth = 4;
128484ff8c8bSHong Zhang
128584ff8c8bSHong Zhang for (i = xs; i < xs + xm; i++) {
128684ff8c8bSHong Zhang if (i < ctx.sf - ctx.lsbwidth || i > ctx.fs + ctx.rsbwidth - 1)
128784ff8c8bSHong Zhang for (k = 0; k < dof; k++) index_slow[islow++] = i * dof + k;
128884ff8c8bSHong Zhang else if ((i >= ctx.sf - ctx.lsbwidth && i < ctx.sf) || (i > ctx.fs - 1 && i <= ctx.fs + ctx.rsbwidth - 1))
128984ff8c8bSHong Zhang for (k = 0; k < dof; k++) index_slowbuffer[islowbuffer++] = i * dof + k;
129084ff8c8bSHong Zhang else
129184ff8c8bSHong Zhang for (k = 0; k < dof; k++) index_fast[ifast++] = i * dof + k;
129284ff8c8bSHong Zhang }
12939566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islow, index_slow, PETSC_COPY_VALUES, &ctx.iss));
12949566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, ifast, index_fast, PETSC_COPY_VALUES, &ctx.isf));
12959566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islowbuffer, index_slowbuffer, PETSC_COPY_VALUES, &ctx.issb));
129684ff8c8bSHong Zhang
129784ff8c8bSHong Zhang /* Create a time-stepping object */
12989566063dSJacob Faibussowitsch PetscCall(TSCreate(comm, &ts));
12999566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da));
13009566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, R, FVRHSFunction_2WaySplit, &ctx));
13019566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts, "slow", ctx.iss));
13029566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts, "slowbuffer", ctx.issb));
13039566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts, "fast", ctx.isf));
13049566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts, "slow", NULL, FVRHSFunctionslow_2WaySplit, &ctx));
13059566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts, "fast", NULL, FVRHSFunctionfast_2WaySplit, &ctx));
13069566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts, "slowbuffer", NULL, FVRHSFunctionslowbuffer_2WaySplit, &ctx));
130784ff8c8bSHong Zhang
13089566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSMPRK));
13099566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, maxtime));
13109566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
131184ff8c8bSHong Zhang
131284ff8c8bSHong Zhang /* Compute initial conditions and starting time step */
13139566063dSJacob Faibussowitsch PetscCall(FVSample_2WaySplit(&ctx, da, 0, X0));
13149566063dSJacob Faibussowitsch PetscCall(FVRHSFunction_2WaySplit(ts, 0, X0, X, (void *)&ctx)); /* Initial function evaluation, only used to determine max speed */
13159566063dSJacob Faibussowitsch PetscCall(VecCopy(X0, X)); /* The function value was not used so we set X=X0 again */
13169566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, ctx.cfl / ctx.cfl_idt));
13179566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); /* Take runtime options */
13189566063dSJacob Faibussowitsch PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD));
131984ff8c8bSHong Zhang {
132084ff8c8bSHong Zhang PetscInt steps;
132184ff8c8bSHong Zhang PetscScalar mass_initial, mass_final, mass_difference, mass_differenceg;
132284ff8c8bSHong Zhang const PetscScalar *ptr_X, *ptr_X0;
132384ff8c8bSHong Zhang const PetscReal hs = (ctx.xmax - ctx.xmin) * 3.0 / 4.0 / count_slow;
132484ff8c8bSHong Zhang const PetscReal hf = (ctx.xmax - ctx.xmin) / 4.0 / count_fast;
132584ff8c8bSHong Zhang
13269566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, X));
13279566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ptime));
13289566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps));
132984ff8c8bSHong Zhang /* calculate the total mass at initial time and final time */
133084ff8c8bSHong Zhang mass_initial = 0.0;
133184ff8c8bSHong Zhang mass_final = 0.0;
13329566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, X0, (void *)&ptr_X0));
13339566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, X, (void *)&ptr_X));
133484ff8c8bSHong Zhang for (i = xs; i < xs + xm; i++) {
133584ff8c8bSHong Zhang if (i < ctx.sf || i > ctx.fs - 1) {
133684ff8c8bSHong Zhang for (k = 0; k < dof; k++) {
133784ff8c8bSHong Zhang mass_initial = mass_initial + hs * ptr_X0[i * dof + k];
133884ff8c8bSHong Zhang mass_final = mass_final + hs * ptr_X[i * dof + k];
133984ff8c8bSHong Zhang }
134084ff8c8bSHong Zhang } else {
134184ff8c8bSHong Zhang for (k = 0; k < dof; k++) {
134284ff8c8bSHong Zhang mass_initial = mass_initial + hf * ptr_X0[i * dof + k];
134384ff8c8bSHong Zhang mass_final = mass_final + hf * ptr_X[i * dof + k];
134484ff8c8bSHong Zhang }
134584ff8c8bSHong Zhang }
134684ff8c8bSHong Zhang }
13479566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, X0, (void *)&ptr_X0));
13489566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, X, (void *)&ptr_X));
134984ff8c8bSHong Zhang mass_difference = mass_final - mass_initial;
1350462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&mass_difference, &mass_differenceg, 1, MPIU_SCALAR, MPIU_SUM, comm));
13519566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Mass difference %g\n", (double)mass_differenceg));
135263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Final time %g, steps %" PetscInt_FMT "\n", (double)ptime, steps));
13539566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Maximum allowable stepsize according to CFL %g\n", (double)(1.0 / ctx.cfl_idt)));
135484ff8c8bSHong Zhang if (ctx.exact) {
135584ff8c8bSHong Zhang PetscReal nrm1 = 0;
13569566063dSJacob Faibussowitsch PetscCall(SolutionErrorNorms_2WaySplit(&ctx, da, ptime, X, &nrm1));
13579566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %g\n", (double)nrm1));
135884ff8c8bSHong Zhang }
135984ff8c8bSHong Zhang if (ctx.simulation) {
136084ff8c8bSHong Zhang PetscReal nrm1 = 0;
136184ff8c8bSHong Zhang PetscViewer fd;
136284ff8c8bSHong Zhang char filename[PETSC_MAX_PATH_LEN] = "binaryoutput";
136384ff8c8bSHong Zhang Vec XR;
136484ff8c8bSHong Zhang PetscBool flg;
136584ff8c8bSHong Zhang const PetscScalar *ptr_XR;
13669566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f", filename, sizeof(filename), &flg));
13673c633725SBarry Smith PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate binary file with the -f option");
13689566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, filename, FILE_MODE_READ, &fd));
13699566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X0, &XR));
13709566063dSJacob Faibussowitsch PetscCall(VecLoad(XR, fd));
13719566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd));
13729566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &ptr_X));
13739566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XR, &ptr_XR));
137484ff8c8bSHong Zhang for (i = xs; i < xs + xm; i++) {
137584ff8c8bSHong Zhang if (i < ctx.sf || i > ctx.fs - 1)
137684ff8c8bSHong Zhang for (k = 0; k < dof; k++) nrm1 = nrm1 + hs * PetscAbs(ptr_X[i * dof + k] - ptr_XR[i * dof + k]);
137784ff8c8bSHong Zhang else
137884ff8c8bSHong Zhang for (k = 0; k < dof; k++) nrm1 = nrm1 + hf * PetscAbs(ptr_X[i * dof + k] - ptr_XR[i * dof + k]);
137984ff8c8bSHong Zhang }
13809566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &ptr_X));
13819566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XR, &ptr_XR));
13829566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %g\n", (double)nrm1));
13839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XR));
138484ff8c8bSHong Zhang }
138584ff8c8bSHong Zhang }
138684ff8c8bSHong Zhang
13879566063dSJacob Faibussowitsch PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD));
13889566063dSJacob Faibussowitsch if (draw & 0x1) PetscCall(VecView(X0, PETSC_VIEWER_DRAW_WORLD));
13899566063dSJacob Faibussowitsch if (draw & 0x2) PetscCall(VecView(X, PETSC_VIEWER_DRAW_WORLD));
139084ff8c8bSHong Zhang if (draw & 0x4) {
139184ff8c8bSHong Zhang Vec Y;
13929566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &Y));
13939566063dSJacob Faibussowitsch PetscCall(FVSample_2WaySplit(&ctx, da, ptime, Y));
13949566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y, -1, X));
13959566063dSJacob Faibussowitsch PetscCall(VecView(Y, PETSC_VIEWER_DRAW_WORLD));
13969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y));
139784ff8c8bSHong Zhang }
139884ff8c8bSHong Zhang
139984ff8c8bSHong Zhang if (view_final) {
140084ff8c8bSHong Zhang PetscViewer viewer;
14019566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, final_fname, &viewer));
14029566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB));
14039566063dSJacob Faibussowitsch PetscCall(VecView(X, viewer));
14049566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer));
14059566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer));
140684ff8c8bSHong Zhang }
140784ff8c8bSHong Zhang
140884ff8c8bSHong Zhang /* Clean up */
14099566063dSJacob Faibussowitsch PetscCall((*ctx.physics2.destroy)(ctx.physics2.user));
14109566063dSJacob Faibussowitsch for (i = 0; i < ctx.physics2.dof; i++) PetscCall(PetscFree(ctx.physics2.fieldname[i]));
14119566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx.physics2.bcinflowindex));
1412d0609cedSBarry Smith PetscCall(PetscFree(ctx.ub));
14139566063dSJacob Faibussowitsch PetscCall(PetscFree4(ctx.R, ctx.Rinv, ctx.cjmpLR, ctx.cslope));
14149566063dSJacob Faibussowitsch PetscCall(PetscFree3(ctx.uLR, ctx.flux, ctx.speeds));
14159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X));
14169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X0));
14179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&R));
14189566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da));
14199566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
14209566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.iss));
14219566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.isf));
14229566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.issb));
14239566063dSJacob Faibussowitsch PetscCall(PetscFree(index_slow));
14249566063dSJacob Faibussowitsch PetscCall(PetscFree(index_fast));
14259566063dSJacob Faibussowitsch PetscCall(PetscFree(index_slowbuffer));
14269566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&limiters));
14279566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&physics));
14289566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
1429b122ec5aSJacob Faibussowitsch return 0;
143084ff8c8bSHong Zhang }
143184ff8c8bSHong Zhang
143284ff8c8bSHong Zhang /*TEST
143384ff8c8bSHong Zhang
143484ff8c8bSHong Zhang build:
143584ff8c8bSHong Zhang requires: !complex !single
143684ff8c8bSHong Zhang depends: finitevolume1d.c
143784ff8c8bSHong Zhang
143884ff8c8bSHong Zhang test:
143984ff8c8bSHong Zhang suffix: 1
1440*188af4bfSBarry Smith args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -limit mc -ts_time_step 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22
144184ff8c8bSHong Zhang output_file: output/ex4_1.out
144284ff8c8bSHong Zhang
144384ff8c8bSHong Zhang test:
144484ff8c8bSHong Zhang suffix: 2
1445*188af4bfSBarry Smith args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -limit mc -ts_time_step 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 1
144684ff8c8bSHong Zhang output_file: output/ex4_1.out
144784ff8c8bSHong Zhang
144884ff8c8bSHong Zhang test:
144984ff8c8bSHong Zhang suffix: 3
1450*188af4bfSBarry Smith args: -da_grid_x 40 -initial 1 -hratio 2 -limit mc -ts_time_step 0.1 -ts_max_steps 24 -ts_max_time 7.0 -ts_type mprk -ts_mprk_type 2a22 -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 0
145184ff8c8bSHong Zhang output_file: output/ex4_3.out
145284ff8c8bSHong Zhang
145384ff8c8bSHong Zhang test:
145484ff8c8bSHong Zhang suffix: 4
145584ff8c8bSHong Zhang nsize: 2
1456*188af4bfSBarry Smith args: -da_grid_x 40 -initial 1 -hratio 2 -limit mc -ts_time_step 0.1 -ts_max_steps 24 -ts_max_time 7.0 -ts_type mprk -ts_mprk_type 2a22 -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 1
145784ff8c8bSHong Zhang output_file: output/ex4_3.out
145884ff8c8bSHong Zhang
145984ff8c8bSHong Zhang test:
146084ff8c8bSHong Zhang suffix: 5
146184ff8c8bSHong Zhang nsize: 4
1462*188af4bfSBarry Smith args: -da_grid_x 40 -initial 1 -hratio 2 -limit mc -ts_time_step 0.1 -ts_max_steps 24 -ts_max_time 7.0 -ts_type mprk -ts_mprk_type 2a22 -physics shallow -bc_type outflow -xmin 0 -xmax 50 -ts_use_splitrhsfunction 1
146384ff8c8bSHong Zhang output_file: output/ex4_3.out
146484ff8c8bSHong Zhang TEST*/
1465