1c4762a1bSJed Brown static const char help[] = "1D periodic Finite Volume solver in slope-limiter form with semidiscrete time stepping.\n"
2c4762a1bSJed Brown " advection - Constant coefficient scalar advection\n"
3c4762a1bSJed Brown " u_t + (a*u)_x = 0\n"
4c4762a1bSJed Brown " for this toy problem, we choose different meshsizes for different sub-domains (slow-medium-fast-medium-slow), \n"
5c4762a1bSJed Brown " the meshsize ratio between two adjacient sub-domains is controlled with -hratio,\n"
6c4762a1bSJed Brown " exact - Exact Riemann solver which usually needs to perform a Newton iteration to connect\n"
7c4762a1bSJed Brown " the states across shocks and rarefactions\n"
8c4762a1bSJed Brown " simulation - use reference solution which is generated by smaller time step size to be true solution,\n"
9c4762a1bSJed Brown " also the reference solution should be generated by user and stored in a binary file.\n"
10c4762a1bSJed Brown " characteristic - Limit the characteristic variables, this is usually preferred (default)\n"
11c4762a1bSJed Brown "Several initial conditions can be chosen with -initial N\n\n"
12c4762a1bSJed Brown "The problem size should be set with -da_grid_x M\n\n";
13c4762a1bSJed Brown
14c4762a1bSJed Brown #include <petscts.h>
15c4762a1bSJed Brown #include <petscdm.h>
16c4762a1bSJed Brown #include <petscdmda.h>
17c4762a1bSJed Brown #include <petscdraw.h>
18c4762a1bSJed Brown #include "finitevolume1d.h"
19c4762a1bSJed Brown
RangeMod(PetscReal a,PetscReal xmin,PetscReal xmax)20d71ae5a4SJacob Faibussowitsch static inline PetscReal RangeMod(PetscReal a, PetscReal xmin, PetscReal xmax)
21d71ae5a4SJacob Faibussowitsch {
229371c9d4SSatish Balay PetscReal range = xmax - xmin;
239371c9d4SSatish Balay return xmin + PetscFmodReal(range + PetscFmodReal(a, range), range);
249371c9d4SSatish Balay }
25c4762a1bSJed Brown
26c4762a1bSJed Brown /* --------------------------------- Advection ----------------------------------- */
27c4762a1bSJed Brown typedef struct {
28c4762a1bSJed Brown PetscReal a; /* advective velocity */
29c4762a1bSJed Brown } AdvectCtx;
30c4762a1bSJed Brown
PhysicsRiemann_Advect(void * vctx,PetscInt m,const PetscScalar * uL,const PetscScalar * uR,PetscScalar * flux,PetscReal * maxspeed)31d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsRiemann_Advect(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed)
32d71ae5a4SJacob Faibussowitsch {
33c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx *)vctx;
34c4762a1bSJed Brown PetscReal speed;
35c4762a1bSJed Brown
36c4762a1bSJed Brown PetscFunctionBeginUser;
37c4762a1bSJed Brown speed = ctx->a;
38c4762a1bSJed Brown flux[0] = PetscMax(0, speed) * uL[0] + PetscMin(0, speed) * uR[0];
39c4762a1bSJed Brown *maxspeed = speed;
403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
41c4762a1bSJed Brown }
42c4762a1bSJed Brown
PhysicsCharacteristic_Advect(void * vctx,PetscInt m,const PetscScalar * u,PetscScalar * X,PetscScalar * Xi,PetscReal * speeds)43d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx, PetscInt m, const PetscScalar *u, PetscScalar *X, PetscScalar *Xi, PetscReal *speeds)
44d71ae5a4SJacob Faibussowitsch {
45c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx *)vctx;
46c4762a1bSJed Brown
47c4762a1bSJed Brown PetscFunctionBeginUser;
48c4762a1bSJed Brown X[0] = 1.;
49c4762a1bSJed Brown Xi[0] = 1.;
50c4762a1bSJed Brown speeds[0] = ctx->a;
513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
52c4762a1bSJed Brown }
53c4762a1bSJed Brown
PhysicsSample_Advect(void * vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal * u)54d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsSample_Advect(void *vctx, PetscInt initial, FVBCType bctype, PetscReal xmin, PetscReal xmax, PetscReal t, PetscReal x, PetscReal *u)
55d71ae5a4SJacob Faibussowitsch {
56c4762a1bSJed Brown AdvectCtx *ctx = (AdvectCtx *)vctx;
57c4762a1bSJed Brown PetscReal a = ctx->a, x0;
58c4762a1bSJed Brown
59c4762a1bSJed Brown PetscFunctionBeginUser;
60c4762a1bSJed Brown switch (bctype) {
61d71ae5a4SJacob Faibussowitsch case FVBC_OUTFLOW:
62d71ae5a4SJacob Faibussowitsch x0 = x - a * t;
63d71ae5a4SJacob Faibussowitsch break;
64d71ae5a4SJacob Faibussowitsch case FVBC_PERIODIC:
65d71ae5a4SJacob Faibussowitsch x0 = RangeMod(x - a * t, xmin, xmax);
66d71ae5a4SJacob Faibussowitsch break;
67d71ae5a4SJacob Faibussowitsch default:
68d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown BCType");
69c4762a1bSJed Brown }
70c4762a1bSJed Brown switch (initial) {
71d71ae5a4SJacob Faibussowitsch case 0:
72d71ae5a4SJacob Faibussowitsch u[0] = (x0 < 0) ? 1 : -1;
73d71ae5a4SJacob Faibussowitsch break;
74d71ae5a4SJacob Faibussowitsch case 1:
75d71ae5a4SJacob Faibussowitsch u[0] = (x0 < 0) ? -1 : 1;
76d71ae5a4SJacob Faibussowitsch break;
77d71ae5a4SJacob Faibussowitsch case 2:
78d71ae5a4SJacob Faibussowitsch u[0] = (0 < x0 && x0 < 1) ? 1 : 0;
79d71ae5a4SJacob Faibussowitsch break;
80d71ae5a4SJacob Faibussowitsch case 3:
81d71ae5a4SJacob Faibussowitsch u[0] = PetscSinReal(2 * PETSC_PI * x0);
82d71ae5a4SJacob Faibussowitsch break;
83d71ae5a4SJacob Faibussowitsch case 4:
84d71ae5a4SJacob Faibussowitsch u[0] = PetscAbs(x0);
85d71ae5a4SJacob Faibussowitsch break;
86d71ae5a4SJacob Faibussowitsch case 5:
87d71ae5a4SJacob Faibussowitsch u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2 * PETSC_PI * x0));
88d71ae5a4SJacob Faibussowitsch break;
89d71ae5a4SJacob Faibussowitsch case 6:
90d71ae5a4SJacob Faibussowitsch u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2 - x0 : 0));
91d71ae5a4SJacob Faibussowitsch break;
92d71ae5a4SJacob Faibussowitsch case 7:
93d71ae5a4SJacob Faibussowitsch u[0] = PetscPowReal(PetscSinReal(PETSC_PI * x0), 10.0);
94d71ae5a4SJacob Faibussowitsch break;
95d71ae5a4SJacob Faibussowitsch default:
96d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition");
97c4762a1bSJed Brown }
983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
99c4762a1bSJed Brown }
100c4762a1bSJed Brown
PhysicsCreate_Advect(FVCtx * ctx)101d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx)
102d71ae5a4SJacob Faibussowitsch {
103c4762a1bSJed Brown AdvectCtx *user;
104c4762a1bSJed Brown
105c4762a1bSJed Brown PetscFunctionBeginUser;
1069566063dSJacob Faibussowitsch PetscCall(PetscNew(&user));
107c4762a1bSJed Brown ctx->physics2.sample2 = PhysicsSample_Advect;
108c4762a1bSJed Brown ctx->physics2.riemann2 = PhysicsRiemann_Advect;
109c4762a1bSJed Brown ctx->physics2.characteristic2 = PhysicsCharacteristic_Advect;
110c4762a1bSJed Brown ctx->physics2.destroy = PhysicsDestroy_SimpleFree;
111c4762a1bSJed Brown ctx->physics2.user = user;
112c4762a1bSJed Brown ctx->physics2.dof = 1;
1139566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy("u", &ctx->physics2.fieldname[0]));
114c4762a1bSJed Brown user->a = 1;
115d0609cedSBarry Smith PetscOptionsBegin(ctx->comm, ctx->prefix, "Options for advection", "");
116d71ae5a4SJacob Faibussowitsch {
117d71ae5a4SJacob Faibussowitsch PetscCall(PetscOptionsReal("-physics_advect_a", "Speed", "", user->a, &user->a, NULL));
118d71ae5a4SJacob Faibussowitsch }
119d0609cedSBarry Smith PetscOptionsEnd();
1203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
121c4762a1bSJed Brown }
122c4762a1bSJed Brown
FVSample_3WaySplit(FVCtx * ctx,DM da,PetscReal time,Vec U)123d71ae5a4SJacob Faibussowitsch PetscErrorCode FVSample_3WaySplit(FVCtx *ctx, DM da, PetscReal time, Vec U)
124d71ae5a4SJacob Faibussowitsch {
125c4762a1bSJed Brown PetscScalar *u, *uj, xj, xi;
126c4762a1bSJed Brown PetscInt i, j, k, dof, xs, xm, Mx;
127c4762a1bSJed Brown const PetscInt N = 200;
128c4762a1bSJed Brown PetscReal hs, hm, hf;
129c4762a1bSJed Brown
130c4762a1bSJed Brown PetscFunctionBeginUser;
1313c633725SBarry Smith PetscCheck(ctx->physics2.sample2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Physics has not provided a sampling function");
1329566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
1339566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
1349566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u));
1359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof, &uj));
136c4762a1bSJed Brown
137c4762a1bSJed Brown hs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm;
138c4762a1bSJed Brown hm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm);
139c4762a1bSJed Brown hf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf);
140c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
141c4762a1bSJed Brown if (i < ctx->sm) {
142c4762a1bSJed Brown xi = ctx->xmin + 0.5 * hs + i * hs;
143c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */
144c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] = 0;
145c4762a1bSJed Brown for (j = 0; j < N + 1; j++) {
146c4762a1bSJed Brown xj = xi + hs * (j - N / 2) / (PetscReal)N;
1479566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
148c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
149c4762a1bSJed Brown }
150c4762a1bSJed Brown } else if (i < ctx->mf) {
151c4762a1bSJed Brown xi = ctx->xmin + ctx->sm * hs + 0.5 * hm + (i - ctx->sm) * hm;
152c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */
153c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] = 0;
154c4762a1bSJed Brown for (j = 0; j < N + 1; j++) {
155c4762a1bSJed Brown xj = xi + hm * (j - N / 2) / (PetscReal)N;
1569566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
157c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
158c4762a1bSJed Brown }
159c4762a1bSJed Brown } else if (i < ctx->fm) {
160c4762a1bSJed Brown xi = ctx->xmin + ctx->sm * hs + (ctx->mf - ctx->sm) * hm + 0.5 * hf + (i - ctx->mf) * hf;
161c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */
162c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] = 0;
163c4762a1bSJed Brown for (j = 0; j < N + 1; j++) {
164c4762a1bSJed Brown xj = xi + hf * (j - N / 2) / (PetscReal)N;
1659566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
166c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
167c4762a1bSJed Brown }
168c4762a1bSJed Brown } else if (i < ctx->ms) {
169c4762a1bSJed Brown xi = ctx->xmin + ctx->sm * hs + (ctx->mf - ctx->sm) * hm + (ctx->fm - ctx->mf) * hf + 0.5 * hm + (i - ctx->fm) * hm;
170c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */
171c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] = 0;
172c4762a1bSJed Brown for (j = 0; j < N + 1; j++) {
173c4762a1bSJed Brown xj = xi + hm * (j - N / 2) / (PetscReal)N;
1749566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
175c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
176c4762a1bSJed Brown }
177c4762a1bSJed Brown } else {
178c4762a1bSJed Brown xi = ctx->xmin + ctx->sm * hs + (ctx->mf - ctx->sm) * hm + (ctx->fm - ctx->mf) * hf + (ctx->ms - ctx->fm) * hm + 0.5 * hs + (i - ctx->ms) * hs;
179c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */
180c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] = 0;
181c4762a1bSJed Brown for (j = 0; j < N + 1; j++) {
182c4762a1bSJed Brown xj = xi + hs * (j - N / 2) / (PetscReal)N;
1839566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
184c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
185c4762a1bSJed Brown }
186c4762a1bSJed Brown }
187c4762a1bSJed Brown }
1889566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u));
1899566063dSJacob Faibussowitsch PetscCall(PetscFree(uj));
1903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
191c4762a1bSJed Brown }
192c4762a1bSJed Brown
SolutionErrorNorms_3WaySplit(FVCtx * ctx,DM da,PetscReal t,Vec X,PetscReal * nrm1)193d71ae5a4SJacob Faibussowitsch static PetscErrorCode SolutionErrorNorms_3WaySplit(FVCtx *ctx, DM da, PetscReal t, Vec X, PetscReal *nrm1)
194d71ae5a4SJacob Faibussowitsch {
195c4762a1bSJed Brown Vec Y;
196c4762a1bSJed Brown PetscInt i, Mx;
197c4762a1bSJed Brown const PetscScalar *ptr_X, *ptr_Y;
198c4762a1bSJed Brown PetscReal hs, hm, hf;
199c4762a1bSJed Brown
200c4762a1bSJed Brown PetscFunctionBeginUser;
2019566063dSJacob Faibussowitsch PetscCall(VecGetSize(X, &Mx));
202c4762a1bSJed Brown hs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm;
203c4762a1bSJed Brown hm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm);
204c4762a1bSJed Brown hf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf);
2059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &Y));
2069566063dSJacob Faibussowitsch PetscCall(FVSample_3WaySplit(ctx, da, t, Y));
2079566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &ptr_X));
2089566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &ptr_Y));
209c4762a1bSJed Brown for (i = 0; i < Mx; i++) {
210c4762a1bSJed Brown if (i < ctx->sm || i > ctx->ms - 1) *nrm1 += hs * PetscAbs(ptr_X[i] - ptr_Y[i]);
211c4762a1bSJed Brown else if (i < ctx->mf || i > ctx->fm - 1) *nrm1 += hm * PetscAbs(ptr_X[i] - ptr_Y[i]);
212c4762a1bSJed Brown else *nrm1 += hf * PetscAbs(ptr_X[i] - ptr_Y[i]);
213c4762a1bSJed Brown }
2149566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &ptr_X));
2159566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &ptr_Y));
2169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y));
2173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
218c4762a1bSJed Brown }
219c4762a1bSJed Brown
FVRHSFunction_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void * vctx)220d71ae5a4SJacob Faibussowitsch PetscErrorCode FVRHSFunction_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx)
221d71ae5a4SJacob Faibussowitsch {
222c4762a1bSJed Brown FVCtx *ctx = (FVCtx *)vctx;
223c4762a1bSJed Brown PetscInt i, j, k, Mx, dof, xs, xm, sm = ctx->sm, mf = ctx->mf, fm = ctx->fm, ms = ctx->ms;
224c4762a1bSJed Brown PetscReal hxf, hxm, hxs, cfl_idt = 0;
225c4762a1bSJed Brown PetscScalar *x, *f, *slope;
226c4762a1bSJed Brown Vec Xloc;
227c4762a1bSJed Brown DM da;
228c4762a1bSJed Brown
229c4762a1bSJed Brown PetscFunctionBeginUser;
2309566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
2319566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc)); /* Xloc contains ghost points */
2329566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); /* Mx is the number of center points */
233c4762a1bSJed Brown hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm;
234c4762a1bSJed Brown hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm);
235c4762a1bSJed Brown hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf);
2369566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); /* X is solution vector which does not contain ghost points */
2379566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
238c4762a1bSJed Brown
239dd8e379bSPierre Jolivet PetscCall(VecZeroEntries(F)); /* F is the right-hand side function corresponds to center points */
240c4762a1bSJed Brown
2419566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x));
2429566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f));
2439566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); /* contains ghost points */
244c4762a1bSJed Brown
2459566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
246c4762a1bSJed Brown
247c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) {
248c4762a1bSJed Brown for (i = xs - 2; i < 0; i++) {
249c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
250c4762a1bSJed Brown }
251c4762a1bSJed Brown for (i = Mx; i < xs + xm + 2; i++) {
252c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
253c4762a1bSJed Brown }
254c4762a1bSJed Brown }
255c4762a1bSJed Brown for (i = xs - 1; i < xs + xm + 1; i++) {
256c4762a1bSJed Brown struct _LimitInfo info;
257c4762a1bSJed Brown PetscScalar *cjmpL, *cjmpR;
258c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
2599566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
260c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
2619566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
262c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0];
263c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof];
264c4762a1bSJed Brown for (j = 0; j < dof; j++) {
265c4762a1bSJed Brown PetscScalar jmpL, jmpR;
266c4762a1bSJed Brown jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
267c4762a1bSJed Brown jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
268c4762a1bSJed Brown for (k = 0; k < dof; k++) {
269c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
270c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
271c4762a1bSJed Brown }
272c4762a1bSJed Brown }
273c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */
274c4762a1bSJed Brown info.m = dof;
275c4762a1bSJed Brown info.hxs = hxs;
276c4762a1bSJed Brown info.hxm = hxm;
277c4762a1bSJed Brown info.hxf = hxf;
278c4762a1bSJed Brown (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope);
279c4762a1bSJed Brown for (j = 0; j < dof; j++) {
280c4762a1bSJed Brown PetscScalar tmp = 0;
281c4762a1bSJed Brown for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
282c4762a1bSJed Brown slope[i * dof + j] = tmp;
283c4762a1bSJed Brown }
284c4762a1bSJed Brown }
285c4762a1bSJed Brown
286c4762a1bSJed Brown for (i = xs; i < xs + xm + 1; i++) {
287c4762a1bSJed Brown PetscReal maxspeed;
288c4762a1bSJed Brown PetscScalar *uL, *uR;
289c4762a1bSJed Brown uL = &ctx->uLR[0];
290c4762a1bSJed Brown uR = &ctx->uLR[dof];
291c4762a1bSJed Brown if (i < sm || i > ms) { /* slow region */
292c4762a1bSJed Brown for (j = 0; j < dof; j++) {
293c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
294c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
295c4762a1bSJed Brown }
2969566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
297c4762a1bSJed Brown if (i > xs) {
298c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxs;
299c4762a1bSJed Brown }
300c4762a1bSJed Brown if (i < xs + xm) {
301c4762a1bSJed Brown for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxs;
302c4762a1bSJed Brown }
303c4762a1bSJed Brown } else if (i == sm) { /* interface between slow and medium component */
304c4762a1bSJed Brown for (j = 0; j < dof; j++) {
305c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
306c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
307c4762a1bSJed Brown }
3089566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
309c4762a1bSJed Brown if (i > xs) {
310c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxs;
311c4762a1bSJed Brown }
312c4762a1bSJed Brown if (i < xs + xm) {
313c4762a1bSJed Brown for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxm;
314c4762a1bSJed Brown }
315c4762a1bSJed Brown } else if (i == ms) { /* interface between medium and slow regions */
316c4762a1bSJed Brown for (j = 0; j < dof; j++) {
317c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
318c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
319c4762a1bSJed Brown }
3209566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
321c4762a1bSJed Brown if (i > xs) {
322c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxm;
323c4762a1bSJed Brown }
324c4762a1bSJed Brown if (i < xs + xm) {
325c4762a1bSJed Brown for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxs;
326c4762a1bSJed Brown }
327c4762a1bSJed Brown } else if (i < mf || i > fm) { /* medium region */
328c4762a1bSJed Brown for (j = 0; j < dof; j++) {
329c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
330c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
331c4762a1bSJed Brown }
3329566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
333c4762a1bSJed Brown if (i > xs) {
334c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxm;
335c4762a1bSJed Brown }
336c4762a1bSJed Brown if (i < xs + xm) {
337c4762a1bSJed Brown for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxm;
338c4762a1bSJed Brown }
339c4762a1bSJed Brown } else if (i == mf) { /* interface between medium and fast regions */
340c4762a1bSJed Brown for (j = 0; j < dof; j++) {
341c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
342c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
343c4762a1bSJed Brown }
3449566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
345c4762a1bSJed Brown if (i > xs) {
346c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxm;
347c4762a1bSJed Brown }
348c4762a1bSJed Brown if (i < xs + xm) {
349c4762a1bSJed Brown for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxf;
350c4762a1bSJed Brown }
351c4762a1bSJed Brown } else if (i == fm) { /* interface between fast and medium regions */
352c4762a1bSJed Brown for (j = 0; j < dof; j++) {
353c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
354c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
355c4762a1bSJed Brown }
3569566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
357c4762a1bSJed Brown if (i > xs) {
358c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxf;
359c4762a1bSJed Brown }
360c4762a1bSJed Brown if (i < xs + xm) {
361c4762a1bSJed Brown for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxm;
362c4762a1bSJed Brown }
363c4762a1bSJed Brown } else { /* fast region */
364c4762a1bSJed Brown for (j = 0; j < dof; j++) {
365c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
366c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
367c4762a1bSJed Brown }
3689566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
369c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hxs)); /* Max allowable value of 1/Delta t */
370c4762a1bSJed Brown if (i > xs) {
371c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxf;
372c4762a1bSJed Brown }
373c4762a1bSJed Brown if (i < xs + xm) {
374c4762a1bSJed Brown for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxf;
375c4762a1bSJed Brown }
376c4762a1bSJed Brown }
377c4762a1bSJed Brown }
3789566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
3799566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f));
3809566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
3819566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc));
382462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&cfl_idt, &ctx->cfl_idt, 1, MPIU_SCALAR, MPIU_MAX, PetscObjectComm((PetscObject)da)));
383c4762a1bSJed Brown if (0) {
384c4762a1bSJed Brown /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */
385c4762a1bSJed Brown PetscReal dt, tnow;
3869566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt));
3879566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &tnow));
388300f1712SStefano Zampini if (dt > 0.5 / ctx->cfl_idt) PetscCall(PetscPrintf(ctx->comm, "Stability constraint exceeded at t=%g, dt %g > %g\n", (double)tnow, (double)dt, (double)(1 / (2 * ctx->cfl_idt))));
389c4762a1bSJed Brown }
3903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
391c4762a1bSJed Brown }
392c4762a1bSJed Brown
393c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */
FVRHSFunctionslow_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void * vctx)394d71ae5a4SJacob Faibussowitsch PetscErrorCode FVRHSFunctionslow_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx)
395d71ae5a4SJacob Faibussowitsch {
396c4762a1bSJed Brown FVCtx *ctx = (FVCtx *)vctx;
397c4762a1bSJed Brown PetscInt i, j, k, Mx, dof, xs, xm, islow = 0, sm = ctx->sm, ms = ctx->ms, lsbwidth = ctx->lsbwidth, rsbwidth = ctx->rsbwidth;
398c4762a1bSJed Brown PetscReal hxs, hxm, hxf, cfl_idt = 0;
399c4762a1bSJed Brown PetscScalar *x, *f, *slope;
400c4762a1bSJed Brown Vec Xloc;
401c4762a1bSJed Brown DM da;
402c4762a1bSJed Brown
403c4762a1bSJed Brown PetscFunctionBeginUser;
4049566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
4059566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc));
4069566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
407c4762a1bSJed Brown hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm;
408c4762a1bSJed Brown hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm);
409c4762a1bSJed Brown hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf);
4109566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
4119566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
4129566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F));
4139566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x));
4149566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
4159566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
4169566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
417c4762a1bSJed Brown
418c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) {
419c4762a1bSJed Brown for (i = xs - 2; i < 0; i++) {
420c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
421c4762a1bSJed Brown }
422c4762a1bSJed Brown for (i = Mx; i < xs + xm + 2; i++) {
423c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
424c4762a1bSJed Brown }
425c4762a1bSJed Brown }
426c4762a1bSJed Brown for (i = xs - 1; i < xs + xm + 1; i++) {
427c4762a1bSJed Brown struct _LimitInfo info;
428c4762a1bSJed Brown PetscScalar *cjmpL, *cjmpR;
429c4762a1bSJed Brown if (i < sm - lsbwidth + 1 || i > ms + rsbwidth - 2) { /* slow components and the first and last fast components */
430c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
4319566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
432c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
4339566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
434c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0];
435c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof];
436c4762a1bSJed Brown for (j = 0; j < dof; j++) {
437c4762a1bSJed Brown PetscScalar jmpL, jmpR;
438c4762a1bSJed Brown jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
439c4762a1bSJed Brown jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
440c4762a1bSJed Brown for (k = 0; k < dof; k++) {
441c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
442c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
443c4762a1bSJed Brown }
444c4762a1bSJed Brown }
445c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */
446c4762a1bSJed Brown info.m = dof;
447c4762a1bSJed Brown info.hxs = hxs;
448c4762a1bSJed Brown info.hxm = hxm;
449c4762a1bSJed Brown info.hxf = hxf;
450c4762a1bSJed Brown (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope);
451c4762a1bSJed Brown for (j = 0; j < dof; j++) {
452c4762a1bSJed Brown PetscScalar tmp = 0;
453c4762a1bSJed Brown for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
454c4762a1bSJed Brown slope[i * dof + j] = tmp;
455c4762a1bSJed Brown }
456c4762a1bSJed Brown }
457c4762a1bSJed Brown }
458c4762a1bSJed Brown
459c4762a1bSJed Brown for (i = xs; i < xs + xm + 1; i++) {
460c4762a1bSJed Brown PetscReal maxspeed;
461c4762a1bSJed Brown PetscScalar *uL, *uR;
462c4762a1bSJed Brown uL = &ctx->uLR[0];
463c4762a1bSJed Brown uR = &ctx->uLR[dof];
464c4762a1bSJed Brown if (i < sm - lsbwidth) { /* slow region */
465c4762a1bSJed Brown for (j = 0; j < dof; j++) {
466c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
467c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
468c4762a1bSJed Brown }
4699566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
470c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hxs)); /* Max allowable value of 1/Delta t */
471c4762a1bSJed Brown if (i > xs) {
472c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
473c4762a1bSJed Brown }
474c4762a1bSJed Brown if (i < xs + xm) {
475c4762a1bSJed Brown for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
476c4762a1bSJed Brown islow++;
477c4762a1bSJed Brown }
478c4762a1bSJed Brown }
479c4762a1bSJed Brown if (i == sm - lsbwidth) { /* interface between slow and medium regions */
480c4762a1bSJed Brown for (j = 0; j < dof; j++) {
481c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
482c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
483c4762a1bSJed Brown }
4849566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
485c4762a1bSJed Brown if (i > xs) {
486c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
487c4762a1bSJed Brown }
488c4762a1bSJed Brown }
489c4762a1bSJed Brown if (i == ms + rsbwidth) { /* interface between medium and slow regions */
490c4762a1bSJed Brown for (j = 0; j < dof; j++) {
491c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
492c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
493c4762a1bSJed Brown }
4949566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
495c4762a1bSJed Brown if (i < xs + xm) {
496c4762a1bSJed Brown for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
497c4762a1bSJed Brown islow++;
498c4762a1bSJed Brown }
499c4762a1bSJed Brown }
500c4762a1bSJed Brown if (i > ms + rsbwidth) { /* slow region */
501c4762a1bSJed Brown for (j = 0; j < dof; j++) {
502c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
503c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
504c4762a1bSJed Brown }
5059566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
506c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hxs)); /* Max allowable value of 1/Delta t */
507c4762a1bSJed Brown if (i > xs) {
508c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
509c4762a1bSJed Brown }
510c4762a1bSJed Brown if (i < xs + xm) {
511c4762a1bSJed Brown for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
512c4762a1bSJed Brown islow++;
513c4762a1bSJed Brown }
514c4762a1bSJed Brown }
515c4762a1bSJed Brown }
5169566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
5179566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
5189566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
5199566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc));
520462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&cfl_idt, &ctx->cfl_idt, 1, MPIU_SCALAR, MPIU_MAX, PetscObjectComm((PetscObject)da)));
5213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
522c4762a1bSJed Brown }
523c4762a1bSJed Brown
FVRHSFunctionslowbuffer_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void * vctx)524d71ae5a4SJacob Faibussowitsch PetscErrorCode FVRHSFunctionslowbuffer_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx)
525d71ae5a4SJacob Faibussowitsch {
526c4762a1bSJed Brown FVCtx *ctx = (FVCtx *)vctx;
527c4762a1bSJed Brown PetscInt i, j, k, Mx, dof, xs, xm, islowbuffer = 0, sm = ctx->sm, ms = ctx->ms, lsbwidth = ctx->lsbwidth, rsbwidth = ctx->rsbwidth;
528c4762a1bSJed Brown PetscReal hxs, hxm, hxf;
529c4762a1bSJed Brown PetscScalar *x, *f, *slope;
530c4762a1bSJed Brown Vec Xloc;
531c4762a1bSJed Brown DM da;
532c4762a1bSJed Brown
533c4762a1bSJed Brown PetscFunctionBeginUser;
5349566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
5359566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc));
5369566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
537c4762a1bSJed Brown hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm;
538c4762a1bSJed Brown hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm);
539c4762a1bSJed Brown hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf);
5409566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
5419566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
5429566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F));
5439566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x));
5449566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
5459566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
5469566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
547c4762a1bSJed Brown
548c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) {
549c4762a1bSJed Brown for (i = xs - 2; i < 0; i++) {
550c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
551c4762a1bSJed Brown }
552c4762a1bSJed Brown for (i = Mx; i < xs + xm + 2; i++) {
553c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
554c4762a1bSJed Brown }
555c4762a1bSJed Brown }
556c4762a1bSJed Brown for (i = xs - 1; i < xs + xm + 1; i++) {
557c4762a1bSJed Brown struct _LimitInfo info;
558c4762a1bSJed Brown PetscScalar *cjmpL, *cjmpR;
559c4762a1bSJed Brown if ((i > sm - lsbwidth - 2 && i < sm + 1) || (i > ms - 2 && i < ms + rsbwidth + 1)) {
560c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
5619566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
562c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
5639566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
564c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0];
565c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof];
566c4762a1bSJed Brown for (j = 0; j < dof; j++) {
567c4762a1bSJed Brown PetscScalar jmpL, jmpR;
568c4762a1bSJed Brown jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
569c4762a1bSJed Brown jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
570c4762a1bSJed Brown for (k = 0; k < dof; k++) {
571c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
572c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
573c4762a1bSJed Brown }
574c4762a1bSJed Brown }
575c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */
576c4762a1bSJed Brown info.m = dof;
577c4762a1bSJed Brown info.hxs = hxs;
578c4762a1bSJed Brown info.hxm = hxm;
579c4762a1bSJed Brown info.hxf = hxf;
580c4762a1bSJed Brown (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope);
581c4762a1bSJed Brown for (j = 0; j < dof; j++) {
582c4762a1bSJed Brown PetscScalar tmp = 0;
583c4762a1bSJed Brown for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
584c4762a1bSJed Brown slope[i * dof + j] = tmp;
585c4762a1bSJed Brown }
586c4762a1bSJed Brown }
587c4762a1bSJed Brown }
588c4762a1bSJed Brown
589c4762a1bSJed Brown for (i = xs; i < xs + xm + 1; i++) {
590c4762a1bSJed Brown PetscReal maxspeed;
591c4762a1bSJed Brown PetscScalar *uL, *uR;
592c4762a1bSJed Brown uL = &ctx->uLR[0];
593c4762a1bSJed Brown uR = &ctx->uLR[dof];
594c4762a1bSJed Brown if (i == sm - lsbwidth) {
595c4762a1bSJed Brown for (j = 0; j < dof; j++) {
596c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
597c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
598c4762a1bSJed Brown }
5999566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
600c4762a1bSJed Brown if (i < xs + xm) {
601c4762a1bSJed Brown for (j = 0; j < dof; j++) f[islowbuffer * dof + j] += ctx->flux[j] / hxs;
602c4762a1bSJed Brown islowbuffer++;
603c4762a1bSJed Brown }
604c4762a1bSJed Brown }
605c4762a1bSJed Brown if (i > sm - lsbwidth && i < sm) {
606c4762a1bSJed Brown for (j = 0; j < dof; j++) {
607c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
608c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
609c4762a1bSJed Brown }
6109566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
611c4762a1bSJed Brown if (i > xs) {
612c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(islowbuffer - 1) * dof + j] -= ctx->flux[j] / hxs;
613c4762a1bSJed Brown }
614c4762a1bSJed Brown if (i < xs + xm) {
615c4762a1bSJed Brown for (j = 0; j < dof; j++) f[islowbuffer * dof + j] += ctx->flux[j] / hxs;
616c4762a1bSJed Brown islowbuffer++;
617c4762a1bSJed Brown }
618c4762a1bSJed Brown }
619c4762a1bSJed Brown if (i == sm) { /* interface between the slow region and the medium region */
620c4762a1bSJed Brown for (j = 0; j < dof; j++) {
621c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
622c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
623c4762a1bSJed Brown }
6249566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
625c4762a1bSJed Brown if (i > xs) {
626c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(islowbuffer - 1) * dof + j] -= ctx->flux[j] / hxs;
627c4762a1bSJed Brown }
628c4762a1bSJed Brown }
629c4762a1bSJed Brown if (i == ms) { /* interface between the medium region and the slow region */
630c4762a1bSJed Brown for (j = 0; j < dof; j++) {
631c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
632c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
633c4762a1bSJed Brown }
6349566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
635c4762a1bSJed Brown if (i < xs + xm) {
636c4762a1bSJed Brown for (j = 0; j < dof; j++) f[islowbuffer * dof + j] += ctx->flux[j] / hxs;
637c4762a1bSJed Brown islowbuffer++;
638c4762a1bSJed Brown }
639c4762a1bSJed Brown }
640c4762a1bSJed Brown if (i > ms && i < ms + rsbwidth) {
641c4762a1bSJed Brown for (j = 0; j < dof; j++) {
642c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
643c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
644c4762a1bSJed Brown }
6459566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
646c4762a1bSJed Brown if (i > xs) {
647c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(islowbuffer - 1) * dof + j] -= ctx->flux[j] / hxs;
648c4762a1bSJed Brown }
649c4762a1bSJed Brown if (i < xs + xm) {
650c4762a1bSJed Brown for (j = 0; j < dof; j++) f[islowbuffer * dof + j] += ctx->flux[j] / hxs;
651c4762a1bSJed Brown islowbuffer++;
652c4762a1bSJed Brown }
653c4762a1bSJed Brown }
654c4762a1bSJed Brown if (i == ms + rsbwidth) {
655c4762a1bSJed Brown for (j = 0; j < dof; j++) {
656c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
657c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
658c4762a1bSJed Brown }
6599566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
660c4762a1bSJed Brown if (i > xs) {
661c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(islowbuffer - 1) * dof + j] -= ctx->flux[j] / hxs;
662c4762a1bSJed Brown }
663c4762a1bSJed Brown }
664c4762a1bSJed Brown }
6659566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
6669566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
6679566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
6689566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc));
6693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
670c4762a1bSJed Brown }
671c4762a1bSJed Brown
672c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for medium components ----------------------------------- */
FVRHSFunctionmedium_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void * vctx)673d71ae5a4SJacob Faibussowitsch PetscErrorCode FVRHSFunctionmedium_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx)
674d71ae5a4SJacob Faibussowitsch {
675c4762a1bSJed Brown FVCtx *ctx = (FVCtx *)vctx;
676c4762a1bSJed Brown PetscInt i, j, k, Mx, dof, xs, xm, imedium = 0, sm = ctx->sm, mf = ctx->mf, fm = ctx->fm, ms = ctx->ms, lmbwidth = ctx->lmbwidth, rmbwidth = ctx->rmbwidth;
677c4762a1bSJed Brown PetscReal hxs, hxm, hxf;
678c4762a1bSJed Brown PetscScalar *x, *f, *slope;
679c4762a1bSJed Brown Vec Xloc;
680c4762a1bSJed Brown DM da;
681c4762a1bSJed Brown
682c4762a1bSJed Brown PetscFunctionBeginUser;
6839566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
6849566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc));
6859566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
686c4762a1bSJed Brown hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm;
687c4762a1bSJed Brown hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm);
688c4762a1bSJed Brown hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf);
6899566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
6909566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
6919566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F));
6929566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x));
6939566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
6949566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
6959566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
696c4762a1bSJed Brown
697c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) {
698c4762a1bSJed Brown for (i = xs - 2; i < 0; i++) {
699c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
700c4762a1bSJed Brown }
701c4762a1bSJed Brown for (i = Mx; i < xs + xm + 2; i++) {
702c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
703c4762a1bSJed Brown }
704c4762a1bSJed Brown }
705c4762a1bSJed Brown for (i = xs - 1; i < xs + xm + 1; i++) {
706c4762a1bSJed Brown struct _LimitInfo info;
707c4762a1bSJed Brown PetscScalar *cjmpL, *cjmpR;
708c4762a1bSJed Brown if ((i > sm - 2 && i < mf - lmbwidth + 1) || (i > fm + rmbwidth - 2 && i < ms + 1)) { /* slow components and the first and last fast components */
709c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
7109566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
711c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
7129566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
713c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0];
714c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof];
715c4762a1bSJed Brown for (j = 0; j < dof; j++) {
716c4762a1bSJed Brown PetscScalar jmpL, jmpR;
717c4762a1bSJed Brown jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
718c4762a1bSJed Brown jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
719c4762a1bSJed Brown for (k = 0; k < dof; k++) {
720c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
721c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
722c4762a1bSJed Brown }
723c4762a1bSJed Brown }
724c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */
725c4762a1bSJed Brown info.m = dof;
726c4762a1bSJed Brown info.hxs = hxs;
727c4762a1bSJed Brown info.hxm = hxm;
728c4762a1bSJed Brown info.hxf = hxf;
729c4762a1bSJed Brown (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope);
730c4762a1bSJed Brown for (j = 0; j < dof; j++) {
731c4762a1bSJed Brown PetscScalar tmp = 0;
732c4762a1bSJed Brown for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
733c4762a1bSJed Brown slope[i * dof + j] = tmp;
734c4762a1bSJed Brown }
735c4762a1bSJed Brown }
736c4762a1bSJed Brown }
737c4762a1bSJed Brown
738c4762a1bSJed Brown for (i = xs; i < xs + xm + 1; i++) {
739c4762a1bSJed Brown PetscReal maxspeed;
740c4762a1bSJed Brown PetscScalar *uL, *uR;
741c4762a1bSJed Brown uL = &ctx->uLR[0];
742c4762a1bSJed Brown uR = &ctx->uLR[dof];
743c4762a1bSJed Brown if (i == sm) { /* interface between slow and medium regions */
744c4762a1bSJed Brown for (j = 0; j < dof; j++) {
745c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
746c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
747c4762a1bSJed Brown }
7489566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
749c4762a1bSJed Brown if (i < xs + xm) {
750c4762a1bSJed Brown for (j = 0; j < dof; j++) f[imedium * dof + j] += ctx->flux[j] / hxm;
751c4762a1bSJed Brown imedium++;
752c4762a1bSJed Brown }
753c4762a1bSJed Brown }
754c4762a1bSJed Brown if (i > sm && i < mf - lmbwidth) { /* medium region */
755c4762a1bSJed Brown for (j = 0; j < dof; j++) {
756c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
757c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
758c4762a1bSJed Brown }
7599566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
760c4762a1bSJed Brown if (i > xs) {
761c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(imedium - 1) * dof + j] -= ctx->flux[j] / hxm;
762c4762a1bSJed Brown }
763c4762a1bSJed Brown if (i < xs + xm) {
764c4762a1bSJed Brown for (j = 0; j < dof; j++) f[imedium * dof + j] += ctx->flux[j] / hxm;
765c4762a1bSJed Brown imedium++;
766c4762a1bSJed Brown }
767c4762a1bSJed Brown }
768c4762a1bSJed Brown if (i == mf - lmbwidth) { /* interface between medium and fast regions */
769c4762a1bSJed Brown for (j = 0; j < dof; j++) {
770c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
771c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
772c4762a1bSJed Brown }
7739566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
774c4762a1bSJed Brown if (i > xs) {
775c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(imedium - 1) * dof + j] -= ctx->flux[j] / hxm;
776c4762a1bSJed Brown }
777c4762a1bSJed Brown }
778c4762a1bSJed Brown if (i == fm + rmbwidth) { /* interface between fast and medium regions */
779c4762a1bSJed Brown for (j = 0; j < dof; j++) {
780c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
781c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
782c4762a1bSJed Brown }
7839566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
784c4762a1bSJed Brown if (i < xs + xm) {
785c4762a1bSJed Brown for (j = 0; j < dof; j++) f[imedium * dof + j] += ctx->flux[j] / hxm;
786c4762a1bSJed Brown imedium++;
787c4762a1bSJed Brown }
788c4762a1bSJed Brown }
789c4762a1bSJed Brown if (i > fm + rmbwidth && i < ms) { /* medium region */
790c4762a1bSJed Brown for (j = 0; j < dof; j++) {
791c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
792c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
793c4762a1bSJed Brown }
7949566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
795c4762a1bSJed Brown if (i > xs) {
796c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(imedium - 1) * dof + j] -= ctx->flux[j] / hxm;
797c4762a1bSJed Brown }
798c4762a1bSJed Brown if (i < xs + xm) {
799c4762a1bSJed Brown for (j = 0; j < dof; j++) f[imedium * dof + j] += ctx->flux[j] / hxm;
800c4762a1bSJed Brown imedium++;
801c4762a1bSJed Brown }
802c4762a1bSJed Brown }
803c4762a1bSJed Brown if (i == ms) { /* interface between medium and slow regions */
804c4762a1bSJed Brown for (j = 0; j < dof; j++) {
805c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
806c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
807c4762a1bSJed Brown }
8089566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
809c4762a1bSJed Brown if (i > xs) {
810c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(imedium - 1) * dof + j] -= ctx->flux[j] / hxm;
811c4762a1bSJed Brown }
812c4762a1bSJed Brown }
813c4762a1bSJed Brown }
8149566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
8159566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
8169566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
8179566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc));
8183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
819c4762a1bSJed Brown }
820c4762a1bSJed Brown
FVRHSFunctionmediumbuffer_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void * vctx)821d71ae5a4SJacob Faibussowitsch PetscErrorCode FVRHSFunctionmediumbuffer_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx)
822d71ae5a4SJacob Faibussowitsch {
823c4762a1bSJed Brown FVCtx *ctx = (FVCtx *)vctx;
824c4762a1bSJed Brown PetscInt i, j, k, Mx, dof, xs, xm, imediumbuffer = 0, mf = ctx->mf, fm = ctx->fm, lmbwidth = ctx->lmbwidth, rmbwidth = ctx->rmbwidth;
825c4762a1bSJed Brown PetscReal hxs, hxm, hxf;
826c4762a1bSJed Brown PetscScalar *x, *f, *slope;
827c4762a1bSJed Brown Vec Xloc;
828c4762a1bSJed Brown DM da;
829c4762a1bSJed Brown
830c4762a1bSJed Brown PetscFunctionBeginUser;
8319566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
8329566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc));
8339566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
834c4762a1bSJed Brown hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm;
835c4762a1bSJed Brown hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm);
836c4762a1bSJed Brown hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf);
8379566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
8389566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
8399566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F));
8409566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x));
8419566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
8429566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
8439566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
844c4762a1bSJed Brown
845c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) {
846c4762a1bSJed Brown for (i = xs - 2; i < 0; i++) {
847c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
848c4762a1bSJed Brown }
849c4762a1bSJed Brown for (i = Mx; i < xs + xm + 2; i++) {
850c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
851c4762a1bSJed Brown }
852c4762a1bSJed Brown }
853c4762a1bSJed Brown for (i = xs - 1; i < xs + xm + 1; i++) {
854c4762a1bSJed Brown struct _LimitInfo info;
855c4762a1bSJed Brown PetscScalar *cjmpL, *cjmpR;
856c4762a1bSJed Brown if ((i > mf - lmbwidth - 2 && i < mf + 1) || (i > fm - 2 && i < fm + rmbwidth + 1)) {
857c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
8589566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
859c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
8609566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
861c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0];
862c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof];
863c4762a1bSJed Brown for (j = 0; j < dof; j++) {
864c4762a1bSJed Brown PetscScalar jmpL, jmpR;
865c4762a1bSJed Brown jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
866c4762a1bSJed Brown jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
867c4762a1bSJed Brown for (k = 0; k < dof; k++) {
868c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
869c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
870c4762a1bSJed Brown }
871c4762a1bSJed Brown }
872c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */
873c4762a1bSJed Brown info.m = dof;
874c4762a1bSJed Brown info.hxs = hxs;
875c4762a1bSJed Brown info.hxm = hxm;
876c4762a1bSJed Brown info.hxf = hxf;
877c4762a1bSJed Brown (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope);
878c4762a1bSJed Brown for (j = 0; j < dof; j++) {
879c4762a1bSJed Brown PetscScalar tmp = 0;
880c4762a1bSJed Brown for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
881c4762a1bSJed Brown slope[i * dof + j] = tmp;
882c4762a1bSJed Brown }
883c4762a1bSJed Brown }
884c4762a1bSJed Brown }
885c4762a1bSJed Brown
886c4762a1bSJed Brown for (i = xs; i < xs + xm + 1; i++) {
887c4762a1bSJed Brown PetscReal maxspeed;
888c4762a1bSJed Brown PetscScalar *uL, *uR;
889c4762a1bSJed Brown uL = &ctx->uLR[0];
890c4762a1bSJed Brown uR = &ctx->uLR[dof];
891c4762a1bSJed Brown if (i == mf - lmbwidth) {
892c4762a1bSJed Brown for (j = 0; j < dof; j++) {
893c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
894c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
895c4762a1bSJed Brown }
8969566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
897c4762a1bSJed Brown if (i < xs + xm) {
898c4762a1bSJed Brown for (j = 0; j < dof; j++) f[imediumbuffer * dof + j] += ctx->flux[j] / hxm;
899c4762a1bSJed Brown imediumbuffer++;
900c4762a1bSJed Brown }
901c4762a1bSJed Brown }
902c4762a1bSJed Brown if (i > mf - lmbwidth && i < mf) {
903c4762a1bSJed Brown for (j = 0; j < dof; j++) {
904c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
905c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
906c4762a1bSJed Brown }
9079566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
908c4762a1bSJed Brown if (i > xs) {
909c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(imediumbuffer - 1) * dof + j] -= ctx->flux[j] / hxm;
910c4762a1bSJed Brown }
911c4762a1bSJed Brown if (i < xs + xm) {
912c4762a1bSJed Brown for (j = 0; j < dof; j++) f[imediumbuffer * dof + j] += ctx->flux[j] / hxm;
913c4762a1bSJed Brown imediumbuffer++;
914c4762a1bSJed Brown }
915c4762a1bSJed Brown }
916c4762a1bSJed Brown if (i == mf) { /* interface between the medium region and the fast region */
917c4762a1bSJed Brown for (j = 0; j < dof; j++) {
918c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
919c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
920c4762a1bSJed Brown }
9219566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
922c4762a1bSJed Brown if (i > xs) {
923c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(imediumbuffer - 1) * dof + j] -= ctx->flux[j] / hxm;
924c4762a1bSJed Brown }
925c4762a1bSJed Brown }
926c4762a1bSJed Brown if (i == fm) { /* interface between the fast region and the medium region */
927c4762a1bSJed Brown for (j = 0; j < dof; j++) {
928c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
929c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
930c4762a1bSJed Brown }
9319566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
932c4762a1bSJed Brown if (i < xs + xm) {
933c4762a1bSJed Brown for (j = 0; j < dof; j++) f[imediumbuffer * dof + j] += ctx->flux[j] / hxm;
934c4762a1bSJed Brown imediumbuffer++;
935c4762a1bSJed Brown }
936c4762a1bSJed Brown }
937c4762a1bSJed Brown if (i > fm && i < fm + rmbwidth) {
938c4762a1bSJed Brown for (j = 0; j < dof; j++) {
939c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
940c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
941c4762a1bSJed Brown }
9429566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
943c4762a1bSJed Brown if (i > xs) {
944c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(imediumbuffer - 1) * dof + j] -= ctx->flux[j] / hxm;
945c4762a1bSJed Brown }
946c4762a1bSJed Brown if (i < xs + xm) {
947c4762a1bSJed Brown for (j = 0; j < dof; j++) f[imediumbuffer * dof + j] += ctx->flux[j] / hxm;
948c4762a1bSJed Brown imediumbuffer++;
949c4762a1bSJed Brown }
950c4762a1bSJed Brown }
951c4762a1bSJed Brown if (i == fm + rmbwidth) {
952c4762a1bSJed Brown for (j = 0; j < dof; j++) {
953c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
954c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
955c4762a1bSJed Brown }
9569566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
957c4762a1bSJed Brown if (i > xs) {
958c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(imediumbuffer - 1) * dof + j] -= ctx->flux[j] / hxm;
959c4762a1bSJed Brown }
960c4762a1bSJed Brown }
961c4762a1bSJed Brown }
9629566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
9639566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
9649566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
9659566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc));
9663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
967c4762a1bSJed Brown }
968c4762a1bSJed Brown
969c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for fast parts ----------------------------------- */
FVRHSFunctionfast_3WaySplit(TS ts,PetscReal time,Vec X,Vec F,void * vctx)970d71ae5a4SJacob Faibussowitsch PetscErrorCode FVRHSFunctionfast_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx)
971d71ae5a4SJacob Faibussowitsch {
972c4762a1bSJed Brown FVCtx *ctx = (FVCtx *)vctx;
973c4762a1bSJed Brown PetscInt i, j, k, Mx, dof, xs, xm, ifast = 0, mf = ctx->mf, fm = ctx->fm;
974c4762a1bSJed Brown PetscReal hxs, hxm, hxf;
975c4762a1bSJed Brown PetscScalar *x, *f, *slope;
976c4762a1bSJed Brown Vec Xloc;
977c4762a1bSJed Brown DM da;
978c4762a1bSJed Brown
979c4762a1bSJed Brown PetscFunctionBeginUser;
9809566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
9819566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc));
9829566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
983c4762a1bSJed Brown hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm;
984c4762a1bSJed Brown hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm);
985c4762a1bSJed Brown hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf);
9869566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
9879566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
9889566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F));
9899566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x));
9909566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
9919566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
9929566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
993c4762a1bSJed Brown
994c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) {
995c4762a1bSJed Brown for (i = xs - 2; i < 0; i++) {
996c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
997c4762a1bSJed Brown }
998c4762a1bSJed Brown for (i = Mx; i < xs + xm + 2; i++) {
999c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
1000c4762a1bSJed Brown }
1001c4762a1bSJed Brown }
10026aad120cSJose E. Roman for (i = xs - 1; i < xs + xm + 1; i++) { /* fast components and the last slow components before fast components and the first slow component after fast components */
1003c4762a1bSJed Brown struct _LimitInfo info;
1004c4762a1bSJed Brown PetscScalar *cjmpL, *cjmpR;
1005c4762a1bSJed Brown if (i > mf - 2 && i < fm + 1) {
10069566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
10079566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
1008c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0];
1009c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof];
1010c4762a1bSJed Brown for (j = 0; j < dof; j++) {
1011c4762a1bSJed Brown PetscScalar jmpL, jmpR;
1012c4762a1bSJed Brown jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
1013c4762a1bSJed Brown jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
1014c4762a1bSJed Brown for (k = 0; k < dof; k++) {
1015c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
1016c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
1017c4762a1bSJed Brown }
1018c4762a1bSJed Brown }
1019c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */
1020c4762a1bSJed Brown info.m = dof;
1021c4762a1bSJed Brown info.hxs = hxs;
1022c4762a1bSJed Brown info.hxm = hxm;
1023c4762a1bSJed Brown info.hxf = hxf;
1024c4762a1bSJed Brown (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope);
1025c4762a1bSJed Brown for (j = 0; j < dof; j++) {
1026c4762a1bSJed Brown PetscScalar tmp = 0;
1027c4762a1bSJed Brown for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
1028c4762a1bSJed Brown slope[i * dof + j] = tmp;
1029c4762a1bSJed Brown }
1030c4762a1bSJed Brown }
1031c4762a1bSJed Brown }
1032c4762a1bSJed Brown
1033c4762a1bSJed Brown for (i = xs; i < xs + xm + 1; i++) {
1034c4762a1bSJed Brown PetscReal maxspeed;
1035c4762a1bSJed Brown PetscScalar *uL, *uR;
1036c4762a1bSJed Brown uL = &ctx->uLR[0];
1037c4762a1bSJed Brown uR = &ctx->uLR[dof];
1038c4762a1bSJed Brown if (i == mf) { /* interface between medium and fast regions */
1039c4762a1bSJed Brown for (j = 0; j < dof; j++) {
1040c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
1041c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
1042c4762a1bSJed Brown }
10439566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
1044c4762a1bSJed Brown if (i < xs + xm) {
1045c4762a1bSJed Brown for (j = 0; j < dof; j++) f[ifast * dof + j] += ctx->flux[j] / hxf;
1046c4762a1bSJed Brown ifast++;
1047c4762a1bSJed Brown }
1048c4762a1bSJed Brown }
1049c4762a1bSJed Brown if (i > mf && i < fm) { /* fast region */
1050c4762a1bSJed Brown for (j = 0; j < dof; j++) {
1051c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
1052c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
1053c4762a1bSJed Brown }
10549566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
1055c4762a1bSJed Brown if (i > xs) {
1056c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(ifast - 1) * dof + j] -= ctx->flux[j] / hxf;
1057c4762a1bSJed Brown }
1058c4762a1bSJed Brown if (i < xs + xm) {
1059c4762a1bSJed Brown for (j = 0; j < dof; j++) f[ifast * dof + j] += ctx->flux[j] / hxf;
1060c4762a1bSJed Brown ifast++;
1061c4762a1bSJed Brown }
1062c4762a1bSJed Brown }
1063c4762a1bSJed Brown if (i == fm) { /* interface between fast and medium regions */
1064c4762a1bSJed Brown for (j = 0; j < dof; j++) {
1065c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
1066c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
1067c4762a1bSJed Brown }
10689566063dSJacob Faibussowitsch PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
1069c4762a1bSJed Brown if (i > xs) {
1070c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(ifast - 1) * dof + j] -= ctx->flux[j] / hxf;
1071c4762a1bSJed Brown }
1072c4762a1bSJed Brown }
1073c4762a1bSJed Brown }
10749566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
10759566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
10769566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
10779566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc));
10783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1079c4762a1bSJed Brown }
1080c4762a1bSJed Brown
main(int argc,char * argv[])1081d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[])
1082d71ae5a4SJacob Faibussowitsch {
1083c4762a1bSJed Brown char lname[256] = "mc", physname[256] = "advect", final_fname[256] = "solution.m";
1084c4762a1bSJed Brown PetscFunctionList limiters = 0, physics = 0;
1085c4762a1bSJed Brown MPI_Comm comm;
1086c4762a1bSJed Brown TS ts;
1087c4762a1bSJed Brown DM da;
1088c4762a1bSJed Brown Vec X, X0, R;
1089c4762a1bSJed Brown FVCtx ctx;
1090c4762a1bSJed Brown PetscInt i, k, dof, xs, xm, Mx, draw = 0, count_slow, count_medium, count_fast, islow = 0, imedium = 0, ifast = 0, *index_slow, *index_medium, *index_fast, islowbuffer = 0, *index_slowbuffer, imediumbuffer = 0, *index_mediumbuffer;
1091c4762a1bSJed Brown PetscBool view_final = PETSC_FALSE;
1092c4762a1bSJed Brown PetscReal ptime;
1093c4762a1bSJed Brown
1094327415f7SBarry Smith PetscFunctionBeginUser;
10959566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, 0, help));
1096c4762a1bSJed Brown comm = PETSC_COMM_WORLD;
10979566063dSJacob Faibussowitsch PetscCall(PetscMemzero(&ctx, sizeof(ctx)));
1098c4762a1bSJed Brown
1099c4762a1bSJed Brown /* Register limiters to be available on the command line */
11009566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "upwind", Limit3_Upwind));
11019566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "lax-wendroff", Limit3_LaxWendroff));
11029566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "beam-warming", Limit3_BeamWarming));
11039566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "fromm", Limit3_Fromm));
11049566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "minmod", Limit3_Minmod));
11059566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "superbee", Limit3_Superbee));
11069566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "mc", Limit3_MC));
11079566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&limiters, "koren3", Limit3_Koren3));
1108c4762a1bSJed Brown
1109c4762a1bSJed Brown /* Register physical models to be available on the command line */
11109566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&physics, "advect", PhysicsCreate_Advect));
1111c4762a1bSJed Brown
1112c4762a1bSJed Brown ctx.comm = comm;
1113c4762a1bSJed Brown ctx.cfl = 0.9;
1114c4762a1bSJed Brown ctx.bctype = FVBC_PERIODIC;
1115c4762a1bSJed Brown ctx.xmin = -1.0;
1116c4762a1bSJed Brown ctx.xmax = 1.0;
1117d0609cedSBarry Smith PetscOptionsBegin(comm, NULL, "Finite Volume solver options", "");
11189566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-xmin", "X min", "", ctx.xmin, &ctx.xmin, NULL));
11199566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-xmax", "X max", "", ctx.xmax, &ctx.xmax, NULL));
11209566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-limit", "Name of flux imiter to use", "", limiters, lname, lname, sizeof(lname), NULL));
11219566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-draw", "Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)", "", draw, &draw, NULL));
11229566063dSJacob 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));
11239566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-initial", "Initial condition (depends on the physics)", "", ctx.initial, &ctx.initial, NULL));
11249566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-exact", "Compare errors with exact solution", "", ctx.exact, &ctx.exact, NULL));
11259566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-simulation", "Compare errors with reference solution", "", ctx.simulation, &ctx.simulation, NULL));
11269566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-cfl", "CFL number to time step at", "", ctx.cfl, &ctx.cfl, NULL));
11279566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-bc_type", "Boundary condition", "", FVBCTypes, (PetscEnum)ctx.bctype, (PetscEnum *)&ctx.bctype, NULL));
11289566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-hratio", "Spacing ratio", "", ctx.hratio, &ctx.hratio, NULL));
1129d0609cedSBarry Smith PetscOptionsEnd();
1130c4762a1bSJed Brown
1131c4762a1bSJed Brown /* Choose the limiter from the list of registered limiters */
11329566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(limiters, lname, &ctx.limit3));
11333c633725SBarry Smith PetscCheck(ctx.limit3, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Limiter '%s' not found", lname);
1134c4762a1bSJed Brown
1135c4762a1bSJed Brown /* Choose the physics from the list of registered models */
1136c4762a1bSJed Brown {
1137c4762a1bSJed Brown PetscErrorCode (*r)(FVCtx *);
11389566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(physics, physname, &r));
11393c633725SBarry Smith PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Physics '%s' not found", physname);
1140c4762a1bSJed Brown /* Create the physics, will set the number of fields and their names */
11419566063dSJacob Faibussowitsch PetscCall((*r)(&ctx));
1142c4762a1bSJed Brown }
1143c4762a1bSJed Brown
1144c4762a1bSJed Brown /* Create a DMDA to manage the parallel grid */
11459566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(comm, DM_BOUNDARY_PERIODIC, 50, ctx.physics2.dof, 2, NULL, &da));
11469566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da));
11479566063dSJacob Faibussowitsch PetscCall(DMSetUp(da));
1148c4762a1bSJed Brown /* Inform the DMDA of the field names provided by the physics. */
1149c4762a1bSJed Brown /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */
115048a46eb9SPierre Jolivet for (i = 0; i < ctx.physics2.dof; i++) PetscCall(DMDASetFieldName(da, i, ctx.physics2.fieldname[i]));
11519566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
11529566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
1153c4762a1bSJed Brown
1154c4762a1bSJed Brown /* Set coordinates of cell centers */
11559566063dSJacob 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));
1156c4762a1bSJed Brown
1157c4762a1bSJed Brown /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */
11589566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(dof * dof, &ctx.R, dof * dof, &ctx.Rinv, 2 * dof, &ctx.cjmpLR, 1 * dof, &ctx.cslope));
11599566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(2 * dof, &ctx.uLR, dof, &ctx.flux, dof, &ctx.speeds));
1160c4762a1bSJed Brown
1161c4762a1bSJed Brown /* Create a vector to store the solution and to save the initial state */
11629566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &X));
11639566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &X0));
11649566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &R));
1165c4762a1bSJed Brown
1166c4762a1bSJed Brown /* create index for slow parts and fast parts,
1167c4762a1bSJed Brown count_slow + count_fast = Mx, counts_slow*hs = 0.5, counts_fast*hf = 0.5 */
1168c4762a1bSJed Brown count_slow = Mx / (1 + ctx.hratio) / (1 + ctx.hratio);
1169c4762a1bSJed Brown count_medium = 2 * ctx.hratio * count_slow;
1170cad9d221SBarry Smith PetscCheck((count_slow % 2) == 0 && (count_medium % 2) == 0, PETSC_COMM_WORLD, PETSC_ERR_USER, "Please adjust grid size Mx (-da_grid_x) and hratio (-hratio) so that Mx/(1+hartio)^2 and Mx*2*hratio/(1+hratio)^2 is even");
1171c4762a1bSJed Brown count_fast = ctx.hratio * ctx.hratio * count_slow;
1172c4762a1bSJed Brown ctx.sm = count_slow / 2;
1173c4762a1bSJed Brown ctx.mf = ctx.sm + count_medium / 2;
1174c4762a1bSJed Brown ctx.fm = ctx.mf + count_fast;
1175c4762a1bSJed Brown ctx.ms = ctx.fm + count_medium / 2;
11769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(xm * dof, &index_slow));
11779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(xm * dof, &index_medium));
11789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(xm * dof, &index_fast));
11799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(6 * dof, &index_slowbuffer));
11809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(6 * dof, &index_mediumbuffer));
1181c4762a1bSJed Brown if (((AdvectCtx *)ctx.physics2.user)->a > 0) {
1182c4762a1bSJed Brown ctx.lsbwidth = 2;
1183c4762a1bSJed Brown ctx.rsbwidth = 4;
1184c4762a1bSJed Brown ctx.lmbwidth = 2;
1185c4762a1bSJed Brown ctx.rmbwidth = 4;
1186c4762a1bSJed Brown } else {
1187c4762a1bSJed Brown ctx.lsbwidth = 4;
1188c4762a1bSJed Brown ctx.rsbwidth = 2;
1189c4762a1bSJed Brown ctx.lmbwidth = 4;
1190c4762a1bSJed Brown ctx.rmbwidth = 2;
1191c4762a1bSJed Brown }
1192c4762a1bSJed Brown
1193c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
1194c4762a1bSJed Brown if (i < ctx.sm - ctx.lsbwidth || i > ctx.ms + ctx.rsbwidth - 1)
1195c4762a1bSJed Brown for (k = 0; k < dof; k++) index_slow[islow++] = i * dof + k;
1196c4762a1bSJed Brown else if ((i >= ctx.sm - ctx.lsbwidth && i < ctx.sm) || (i > ctx.ms - 1 && i <= ctx.ms + ctx.rsbwidth - 1))
1197c4762a1bSJed Brown for (k = 0; k < dof; k++) index_slowbuffer[islowbuffer++] = i * dof + k;
1198c4762a1bSJed Brown else if (i < ctx.mf - ctx.lmbwidth || i > ctx.fm + ctx.rmbwidth - 1)
1199c4762a1bSJed Brown for (k = 0; k < dof; k++) index_medium[imedium++] = i * dof + k;
1200c4762a1bSJed Brown else if ((i >= ctx.mf - ctx.lmbwidth && i < ctx.mf) || (i > ctx.fm - 1 && i <= ctx.fm + ctx.rmbwidth - 1))
1201c4762a1bSJed Brown for (k = 0; k < dof; k++) index_mediumbuffer[imediumbuffer++] = i * dof + k;
1202c4762a1bSJed Brown else
1203c4762a1bSJed Brown for (k = 0; k < dof; k++) index_fast[ifast++] = i * dof + k;
1204c4762a1bSJed Brown }
12059566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islow, index_slow, PETSC_COPY_VALUES, &ctx.iss));
12069566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, imedium, index_medium, PETSC_COPY_VALUES, &ctx.ism));
12079566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, ifast, index_fast, PETSC_COPY_VALUES, &ctx.isf));
12089566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islowbuffer, index_slowbuffer, PETSC_COPY_VALUES, &ctx.issb));
12099566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, imediumbuffer, index_mediumbuffer, PETSC_COPY_VALUES, &ctx.ismb));
1210c4762a1bSJed Brown
1211c4762a1bSJed Brown /* Create a time-stepping object */
12129566063dSJacob Faibussowitsch PetscCall(TSCreate(comm, &ts));
12139566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da));
12149566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, R, FVRHSFunction_3WaySplit, &ctx));
12159566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts, "slow", ctx.iss));
12169566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts, "medium", ctx.ism));
12179566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts, "fast", ctx.isf));
12189566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts, "slowbuffer", ctx.issb));
12199566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts, "mediumbuffer", ctx.ismb));
12209566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts, "slow", NULL, FVRHSFunctionslow_3WaySplit, &ctx));
12219566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts, "medium", NULL, FVRHSFunctionmedium_3WaySplit, &ctx));
12229566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts, "fast", NULL, FVRHSFunctionfast_3WaySplit, &ctx));
12239566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts, "slowbuffer", NULL, FVRHSFunctionslowbuffer_3WaySplit, &ctx));
12249566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts, "mediumbuffer", NULL, FVRHSFunctionmediumbuffer_3WaySplit, &ctx));
1225c4762a1bSJed Brown
12269566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSSSP));
12279566063dSJacob Faibussowitsch /*PetscCall(TSSetType(ts,TSMPRK));*/
12289566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 10));
12299566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1230c4762a1bSJed Brown
1231c4762a1bSJed Brown /* Compute initial conditions and starting time step */
12329566063dSJacob Faibussowitsch PetscCall(FVSample_3WaySplit(&ctx, da, 0, X0));
12339566063dSJacob Faibussowitsch PetscCall(FVRHSFunction_3WaySplit(ts, 0, X0, X, (void *)&ctx)); /* Initial function evaluation, only used to determine max speed */
12349566063dSJacob Faibussowitsch PetscCall(VecCopy(X0, X)); /* The function value was not used so we set X=X0 again */
12359566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, ctx.cfl / ctx.cfl_idt));
12369566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); /* Take runtime options */
12379566063dSJacob Faibussowitsch PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD));
1238c4762a1bSJed Brown {
1239c4762a1bSJed Brown PetscInt steps;
1240c4762a1bSJed Brown PetscScalar mass_initial, mass_final, mass_difference, mass_differenceg;
1241c4762a1bSJed Brown const PetscScalar *ptr_X, *ptr_X0;
1242c4762a1bSJed Brown const PetscReal hs = (ctx.xmax - ctx.xmin) / 4.0 / count_slow;
1243c4762a1bSJed Brown const PetscReal hm = (ctx.xmax - ctx.xmin) / 2.0 / count_medium;
1244c4762a1bSJed Brown const PetscReal hf = (ctx.xmax - ctx.xmin) / 4.0 / count_fast;
1245c4762a1bSJed Brown
12469566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, X));
12479566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ptime));
12489566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps));
1249c4762a1bSJed Brown /* calculate the total mass at initial time and final time */
1250c4762a1bSJed Brown mass_initial = 0.0;
1251c4762a1bSJed Brown mass_final = 0.0;
12529566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, X0, (void *)&ptr_X0));
12539566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, X, (void *)&ptr_X));
1254c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
1255c4762a1bSJed Brown if (i < ctx.sm || i > ctx.ms - 1)
1256c4762a1bSJed Brown for (k = 0; k < dof; k++) {
1257c4762a1bSJed Brown mass_initial = mass_initial + hs * ptr_X0[i * dof + k];
1258c4762a1bSJed Brown mass_final = mass_final + hs * ptr_X[i * dof + k];
1259c4762a1bSJed Brown }
1260c4762a1bSJed Brown else if (i < ctx.mf || i > ctx.fm - 1)
1261c4762a1bSJed Brown for (k = 0; k < dof; k++) {
1262c4762a1bSJed Brown mass_initial = mass_initial + hm * ptr_X0[i * dof + k];
1263c4762a1bSJed Brown mass_final = mass_final + hm * ptr_X[i * dof + k];
1264c4762a1bSJed Brown }
1265c4762a1bSJed Brown else {
1266c4762a1bSJed Brown for (k = 0; k < dof; k++) {
1267c4762a1bSJed Brown mass_initial = mass_initial + hf * ptr_X0[i * dof + k];
1268c4762a1bSJed Brown mass_final = mass_final + hf * ptr_X[i * dof + k];
1269c4762a1bSJed Brown }
1270c4762a1bSJed Brown }
1271c4762a1bSJed Brown }
12729566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, X0, (void *)&ptr_X0));
12739566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, X, (void *)&ptr_X));
1274c4762a1bSJed Brown mass_difference = mass_final - mass_initial;
1275462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&mass_difference, &mass_differenceg, 1, MPIU_SCALAR, MPIU_SUM, comm));
12769566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Mass difference %g\n", (double)mass_differenceg));
127763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Final time %g, steps %" PetscInt_FMT "\n", (double)ptime, steps));
1278300f1712SStefano Zampini PetscCall(PetscPrintf(comm, "Maximum allowable stepsize according to CFL %g\n", (double)(1 / ctx.cfl_idt)));
1279c4762a1bSJed Brown if (ctx.exact) {
1280c4762a1bSJed Brown PetscReal nrm1 = 0;
12819566063dSJacob Faibussowitsch PetscCall(SolutionErrorNorms_3WaySplit(&ctx, da, ptime, X, &nrm1));
12829566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %g\n", (double)nrm1));
1283c4762a1bSJed Brown }
1284c4762a1bSJed Brown if (ctx.simulation) {
1285c4762a1bSJed Brown PetscReal nrm1 = 0;
1286c4762a1bSJed Brown PetscViewer fd;
1287c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN] = "binaryoutput";
1288c4762a1bSJed Brown Vec XR;
1289c4762a1bSJed Brown PetscBool flg;
1290c4762a1bSJed Brown const PetscScalar *ptr_XR;
12919566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f", filename, sizeof(filename), &flg));
12923c633725SBarry Smith PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate binary file with the -f option");
12939566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, filename, FILE_MODE_READ, &fd));
12949566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X0, &XR));
12959566063dSJacob Faibussowitsch PetscCall(VecLoad(XR, fd));
12969566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd));
12979566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &ptr_X));
12989566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(XR, &ptr_XR));
1299c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
1300c4762a1bSJed Brown if (i < ctx.sm || i > ctx.ms - 1)
1301c4762a1bSJed Brown for (k = 0; k < dof; k++) nrm1 = nrm1 + hs * PetscAbs(ptr_X[i * dof + k] - ptr_XR[i * dof + k]);
1302c4762a1bSJed Brown else if (i < ctx.mf || i < ctx.fm - 1)
1303c4762a1bSJed Brown for (k = 0; k < dof; k++) nrm1 = nrm1 + hm * PetscAbs(ptr_X[i * dof + k] - ptr_XR[i * dof + k]);
1304c4762a1bSJed Brown else
1305c4762a1bSJed Brown for (k = 0; k < dof; k++) nrm1 = nrm1 + hf * PetscAbs(ptr_X[i * dof + k] - ptr_XR[i * dof + k]);
1306c4762a1bSJed Brown }
13079566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &ptr_X));
13089566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(XR, &ptr_XR));
13099566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %g\n", (double)nrm1));
13109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&XR));
1311c4762a1bSJed Brown }
1312c4762a1bSJed Brown }
1313c4762a1bSJed Brown
13149566063dSJacob Faibussowitsch PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD));
13159566063dSJacob Faibussowitsch if (draw & 0x1) PetscCall(VecView(X0, PETSC_VIEWER_DRAW_WORLD));
13169566063dSJacob Faibussowitsch if (draw & 0x2) PetscCall(VecView(X, PETSC_VIEWER_DRAW_WORLD));
1317c4762a1bSJed Brown if (draw & 0x4) {
1318c4762a1bSJed Brown Vec Y;
13199566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &Y));
13209566063dSJacob Faibussowitsch PetscCall(FVSample_3WaySplit(&ctx, da, ptime, Y));
13219566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y, -1, X));
13229566063dSJacob Faibussowitsch PetscCall(VecView(Y, PETSC_VIEWER_DRAW_WORLD));
13239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y));
1324c4762a1bSJed Brown }
1325c4762a1bSJed Brown
1326c4762a1bSJed Brown if (view_final) {
1327c4762a1bSJed Brown PetscViewer viewer;
13289566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, final_fname, &viewer));
13299566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB));
13309566063dSJacob Faibussowitsch PetscCall(VecView(X, viewer));
13319566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer));
13329566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer));
1333c4762a1bSJed Brown }
1334c4762a1bSJed Brown
1335c4762a1bSJed Brown /* Clean up */
13369566063dSJacob Faibussowitsch PetscCall((*ctx.physics2.destroy)(ctx.physics2.user));
13379566063dSJacob Faibussowitsch for (i = 0; i < ctx.physics2.dof; i++) PetscCall(PetscFree(ctx.physics2.fieldname[i]));
13389566063dSJacob Faibussowitsch PetscCall(PetscFree4(ctx.R, ctx.Rinv, ctx.cjmpLR, ctx.cslope));
13399566063dSJacob Faibussowitsch PetscCall(PetscFree3(ctx.uLR, ctx.flux, ctx.speeds));
13409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X));
13419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X0));
13429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&R));
13439566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da));
13449566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
13459566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.iss));
13469566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.ism));
13479566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.isf));
13489566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.issb));
13499566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ctx.ismb));
13509566063dSJacob Faibussowitsch PetscCall(PetscFree(index_slow));
13519566063dSJacob Faibussowitsch PetscCall(PetscFree(index_medium));
13529566063dSJacob Faibussowitsch PetscCall(PetscFree(index_fast));
13539566063dSJacob Faibussowitsch PetscCall(PetscFree(index_slowbuffer));
13549566063dSJacob Faibussowitsch PetscCall(PetscFree(index_mediumbuffer));
13559566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&limiters));
13569566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&physics));
13579566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
1358b122ec5aSJacob Faibussowitsch return 0;
1359c4762a1bSJed Brown }
1360c4762a1bSJed Brown
1361c4762a1bSJed Brown /*TEST
1362c4762a1bSJed Brown
1363c4762a1bSJed Brown build:
1364f56ea12dSJed Brown requires: !complex
1365c4762a1bSJed Brown depends: finitevolume1d.c
1366c4762a1bSJed Brown
1367c4762a1bSJed Brown test:
1368c4762a1bSJed Brown suffix: 1
1369*188af4bfSBarry Smith args: -da_grid_x 90 -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 2a23 -ts_use_splitrhsfunction 0
1370c4762a1bSJed Brown
1371c4762a1bSJed Brown test:
1372c4762a1bSJed Brown suffix: 2
1373*188af4bfSBarry Smith args: -da_grid_x 90 -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 2a23 -ts_use_splitrhsfunction 1
1374c4762a1bSJed Brown
1375c4762a1bSJed Brown TEST*/
1376