xref: /petsc/src/ts/tutorials/multirate/ex8.c (revision a69119a591a03a9d906b29c0a4e9802e4d7c9795)
1 static const char help[] = "1D periodic Finite Volume solver in slope-limiter form with semidiscrete time stepping.\n"
2                            "  advection   - Constant coefficient scalar advection\n"
3                            "                u_t       + (a*u)_x               = 0\n"
4                            "  for this toy problem, we choose different meshsizes for different sub-domains (slow-medium-fast-medium-slow), \n"
5                            "  the meshsize ratio between two adjacient sub-domains is controlled with -hratio,\n"
6                            "  exact       - Exact Riemann solver which usually needs to perform a Newton iteration to connect\n"
7                            "                the states across shocks and rarefactions\n"
8                            "  simulation  - use reference solution which is generated by smaller time step size to be true solution,\n"
9                            "                also the reference solution should be generated by user and stored in a binary file.\n"
10                            "  characteristic - Limit the characteristic variables, this is usually preferred (default)\n"
11                            "Several initial conditions can be chosen with -initial N\n\n"
12                            "The problem size should be set with -da_grid_x M\n\n";
13 
14 #include <petscts.h>
15 #include <petscdm.h>
16 #include <petscdmda.h>
17 #include <petscdraw.h>
18 #include "finitevolume1d.h"
19 
20 static inline PetscReal RangeMod(PetscReal a, PetscReal xmin, PetscReal xmax) {
21   PetscReal range = xmax - xmin;
22   return xmin + PetscFmodReal(range + PetscFmodReal(a, range), range);
23 }
24 
25 /* --------------------------------- Advection ----------------------------------- */
26 typedef struct {
27   PetscReal a; /* advective velocity */
28 } AdvectCtx;
29 
30 static PetscErrorCode PhysicsRiemann_Advect(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed) {
31   AdvectCtx *ctx = (AdvectCtx *)vctx;
32   PetscReal  speed;
33 
34   PetscFunctionBeginUser;
35   speed     = ctx->a;
36   flux[0]   = PetscMax(0, speed) * uL[0] + PetscMin(0, speed) * uR[0];
37   *maxspeed = speed;
38   PetscFunctionReturn(0);
39 }
40 
41 static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx, PetscInt m, const PetscScalar *u, PetscScalar *X, PetscScalar *Xi, PetscReal *speeds) {
42   AdvectCtx *ctx = (AdvectCtx *)vctx;
43 
44   PetscFunctionBeginUser;
45   X[0]      = 1.;
46   Xi[0]     = 1.;
47   speeds[0] = ctx->a;
48   PetscFunctionReturn(0);
49 }
50 
51 static PetscErrorCode PhysicsSample_Advect(void *vctx, PetscInt initial, FVBCType bctype, PetscReal xmin, PetscReal xmax, PetscReal t, PetscReal x, PetscReal *u) {
52   AdvectCtx *ctx = (AdvectCtx *)vctx;
53   PetscReal  a   = ctx->a, x0;
54 
55   PetscFunctionBeginUser;
56   switch (bctype) {
57   case FVBC_OUTFLOW: x0 = x - a * t; break;
58   case FVBC_PERIODIC: x0 = RangeMod(x - a * t, xmin, xmax); break;
59   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown BCType");
60   }
61   switch (initial) {
62   case 0: u[0] = (x0 < 0) ? 1 : -1; break;
63   case 1: u[0] = (x0 < 0) ? -1 : 1; break;
64   case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break;
65   case 3: u[0] = PetscSinReal(2 * PETSC_PI * x0); break;
66   case 4: u[0] = PetscAbs(x0); break;
67   case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2 * PETSC_PI * x0)); break;
68   case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2 - x0 : 0)); break;
69   case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI * x0), 10.0); break;
70   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition");
71   }
72   PetscFunctionReturn(0);
73 }
74 
75 static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx) {
76   AdvectCtx *user;
77 
78   PetscFunctionBeginUser;
79   PetscCall(PetscNew(&user));
80   ctx->physics2.sample2         = PhysicsSample_Advect;
81   ctx->physics2.riemann2        = PhysicsRiemann_Advect;
82   ctx->physics2.characteristic2 = PhysicsCharacteristic_Advect;
83   ctx->physics2.destroy         = PhysicsDestroy_SimpleFree;
84   ctx->physics2.user            = user;
85   ctx->physics2.dof             = 1;
86   PetscCall(PetscStrallocpy("u", &ctx->physics2.fieldname[0]));
87   user->a = 1;
88   PetscOptionsBegin(ctx->comm, ctx->prefix, "Options for advection", "");
89   { PetscCall(PetscOptionsReal("-physics_advect_a", "Speed", "", user->a, &user->a, NULL)); }
90   PetscOptionsEnd();
91   PetscFunctionReturn(0);
92 }
93 
94 PetscErrorCode FVSample_3WaySplit(FVCtx *ctx, DM da, PetscReal time, Vec U) {
95   PetscScalar   *u, *uj, xj, xi;
96   PetscInt       i, j, k, dof, xs, xm, Mx;
97   const PetscInt N = 200;
98   PetscReal      hs, hm, hf;
99 
100   PetscFunctionBeginUser;
101   PetscCheck(ctx->physics2.sample2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Physics has not provided a sampling function");
102   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
103   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
104   PetscCall(DMDAVecGetArray(da, U, &u));
105   PetscCall(PetscMalloc1(dof, &uj));
106 
107   hs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm;
108   hm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm);
109   hf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf);
110   for (i = xs; i < xs + xm; i++) {
111     if (i < ctx->sm) {
112       xi = ctx->xmin + 0.5 * hs + i * hs;
113       /* Integrate over cell i using trapezoid rule with N points. */
114       for (k = 0; k < dof; k++) u[i * dof + k] = 0;
115       for (j = 0; j < N + 1; j++) {
116         xj = xi + hs * (j - N / 2) / (PetscReal)N;
117         PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
118         for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
119       }
120     } else if (i < ctx->mf) {
121       xi = ctx->xmin + ctx->sm * hs + 0.5 * hm + (i - ctx->sm) * hm;
122       /* Integrate over cell i using trapezoid rule with N points. */
123       for (k = 0; k < dof; k++) u[i * dof + k] = 0;
124       for (j = 0; j < N + 1; j++) {
125         xj = xi + hm * (j - N / 2) / (PetscReal)N;
126         PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
127         for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
128       }
129     } else if (i < ctx->fm) {
130       xi = ctx->xmin + ctx->sm * hs + (ctx->mf - ctx->sm) * hm + 0.5 * hf + (i - ctx->mf) * hf;
131       /* Integrate over cell i using trapezoid rule with N points. */
132       for (k = 0; k < dof; k++) u[i * dof + k] = 0;
133       for (j = 0; j < N + 1; j++) {
134         xj = xi + hf * (j - N / 2) / (PetscReal)N;
135         PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
136         for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
137       }
138     } else if (i < ctx->ms) {
139       xi = ctx->xmin + ctx->sm * hs + (ctx->mf - ctx->sm) * hm + (ctx->fm - ctx->mf) * hf + 0.5 * hm + (i - ctx->fm) * hm;
140       /* Integrate over cell i using trapezoid rule with N points. */
141       for (k = 0; k < dof; k++) u[i * dof + k] = 0;
142       for (j = 0; j < N + 1; j++) {
143         xj = xi + hm * (j - N / 2) / (PetscReal)N;
144         PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
145         for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
146       }
147     } else {
148       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;
149       /* Integrate over cell i using trapezoid rule with N points. */
150       for (k = 0; k < dof; k++) u[i * dof + k] = 0;
151       for (j = 0; j < N + 1; j++) {
152         xj = xi + hs * (j - N / 2) / (PetscReal)N;
153         PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
154         for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
155       }
156     }
157   }
158   PetscCall(DMDAVecRestoreArray(da, U, &u));
159   PetscCall(PetscFree(uj));
160   PetscFunctionReturn(0);
161 }
162 
163 static PetscErrorCode SolutionErrorNorms_3WaySplit(FVCtx *ctx, DM da, PetscReal t, Vec X, PetscReal *nrm1) {
164   Vec                Y;
165   PetscInt           i, Mx;
166   const PetscScalar *ptr_X, *ptr_Y;
167   PetscReal          hs, hm, hf;
168 
169   PetscFunctionBeginUser;
170   PetscCall(VecGetSize(X, &Mx));
171   hs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm;
172   hm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm);
173   hf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf);
174   PetscCall(VecDuplicate(X, &Y));
175   PetscCall(FVSample_3WaySplit(ctx, da, t, Y));
176   PetscCall(VecGetArrayRead(X, &ptr_X));
177   PetscCall(VecGetArrayRead(Y, &ptr_Y));
178   for (i = 0; i < Mx; i++) {
179     if (i < ctx->sm || i > ctx->ms - 1) *nrm1 += hs * PetscAbs(ptr_X[i] - ptr_Y[i]);
180     else if (i < ctx->mf || i > ctx->fm - 1) *nrm1 += hm * PetscAbs(ptr_X[i] - ptr_Y[i]);
181     else *nrm1 += hf * PetscAbs(ptr_X[i] - ptr_Y[i]);
182   }
183   PetscCall(VecRestoreArrayRead(X, &ptr_X));
184   PetscCall(VecRestoreArrayRead(Y, &ptr_Y));
185   PetscCall(VecDestroy(&Y));
186   PetscFunctionReturn(0);
187 }
188 
189 PetscErrorCode FVRHSFunction_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) {
190   FVCtx       *ctx = (FVCtx *)vctx;
191   PetscInt     i, j, k, Mx, dof, xs, xm, sm = ctx->sm, mf = ctx->mf, fm = ctx->fm, ms = ctx->ms;
192   PetscReal    hxf, hxm, hxs, cfl_idt = 0;
193   PetscScalar *x, *f, *slope;
194   Vec          Xloc;
195   DM           da;
196 
197   PetscFunctionBeginUser;
198   PetscCall(TSGetDM(ts, &da));
199   PetscCall(DMGetLocalVector(da, &Xloc));                                 /* Xloc contains ghost points                                     */
200   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); /* Mx is the number of center points                              */
201   hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm;
202   hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm);
203   hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf);
204   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); /* X is solution vector which does not contain ghost points       */
205   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
206 
207   PetscCall(VecZeroEntries(F)); /* F is the right hand side function corresponds to center points */
208 
209   PetscCall(DMDAVecGetArray(da, Xloc, &x));
210   PetscCall(DMDAVecGetArray(da, F, &f));
211   PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); /* contains ghost points                                           */
212 
213   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
214 
215   if (ctx->bctype == FVBC_OUTFLOW) {
216     for (i = xs - 2; i < 0; i++) {
217       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
218     }
219     for (i = Mx; i < xs + xm + 2; i++) {
220       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
221     }
222   }
223   for (i = xs - 1; i < xs + xm + 1; i++) {
224     struct _LimitInfo info;
225     PetscScalar      *cjmpL, *cjmpR;
226     /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
227     PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
228     /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
229     PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
230     cjmpL = &ctx->cjmpLR[0];
231     cjmpR = &ctx->cjmpLR[dof];
232     for (j = 0; j < dof; j++) {
233       PetscScalar jmpL, jmpR;
234       jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
235       jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
236       for (k = 0; k < dof; k++) {
237         cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
238         cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
239       }
240     }
241     /* Apply limiter to the left and right characteristic jumps */
242     info.m   = dof;
243     info.hxs = hxs;
244     info.hxm = hxm;
245     info.hxf = hxf;
246     (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope);
247     for (j = 0; j < dof; j++) {
248       PetscScalar tmp = 0;
249       for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
250       slope[i * dof + j] = tmp;
251     }
252   }
253 
254   for (i = xs; i < xs + xm + 1; i++) {
255     PetscReal    maxspeed;
256     PetscScalar *uL, *uR;
257     uL = &ctx->uLR[0];
258     uR = &ctx->uLR[dof];
259     if (i < sm || i > ms) { /* slow region */
260       for (j = 0; j < dof; j++) {
261         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
262         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
263       }
264       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
265       if (i > xs) {
266         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxs;
267       }
268       if (i < xs + xm) {
269         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxs;
270       }
271     } else if (i == sm) { /* interface between slow and medium component */
272       for (j = 0; j < dof; j++) {
273         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
274         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
275       }
276       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
277       if (i > xs) {
278         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxs;
279       }
280       if (i < xs + xm) {
281         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxm;
282       }
283     } else if (i == ms) { /* interface between medium and slow regions */
284       for (j = 0; j < dof; j++) {
285         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
286         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
287       }
288       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
289       if (i > xs) {
290         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxm;
291       }
292       if (i < xs + xm) {
293         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxs;
294       }
295     } else if (i < mf || i > fm) { /* medium region */
296       for (j = 0; j < dof; j++) {
297         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
298         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
299       }
300       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
301       if (i > xs) {
302         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxm;
303       }
304       if (i < xs + xm) {
305         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxm;
306       }
307     } else if (i == mf) { /* interface between medium and fast regions */
308       for (j = 0; j < dof; j++) {
309         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
310         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
311       }
312       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
313       if (i > xs) {
314         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxm;
315       }
316       if (i < xs + xm) {
317         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxf;
318       }
319     } else if (i == fm) { /* interface between fast and medium regions */
320       for (j = 0; j < dof; j++) {
321         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
322         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
323       }
324       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
325       if (i > xs) {
326         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxf;
327       }
328       if (i < xs + xm) {
329         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxm;
330       }
331     } else { /* fast region */
332       for (j = 0; j < dof; j++) {
333         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
334         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
335       }
336       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
337       cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hxs)); /* Max allowable value of 1/Delta t */
338       if (i > xs) {
339         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxf;
340       }
341       if (i < xs + xm) {
342         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxf;
343       }
344     }
345   }
346   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
347   PetscCall(DMDAVecRestoreArray(da, F, &f));
348   PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
349   PetscCall(DMRestoreLocalVector(da, &Xloc));
350   PetscCallMPI(MPI_Allreduce(&cfl_idt, &ctx->cfl_idt, 1, MPIU_SCALAR, MPIU_MAX, PetscObjectComm((PetscObject)da)));
351   if (0) {
352     /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */
353     PetscReal dt, tnow;
354     PetscCall(TSGetTimeStep(ts, &dt));
355     PetscCall(TSGetTime(ts, &tnow));
356     if (dt > 0.5 / ctx->cfl_idt) PetscCall(PetscPrintf(ctx->comm, "Stability constraint exceeded at t=%g, dt %g > %g\n", (double)tnow, (double)dt, (double)(0.5 / ctx->cfl_idt)));
357   }
358   PetscFunctionReturn(0);
359 }
360 
361 /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */
362 PetscErrorCode FVRHSFunctionslow_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) {
363   FVCtx       *ctx = (FVCtx *)vctx;
364   PetscInt     i, j, k, Mx, dof, xs, xm, islow = 0, sm = ctx->sm, ms = ctx->ms, lsbwidth = ctx->lsbwidth, rsbwidth = ctx->rsbwidth;
365   PetscReal    hxs, hxm, hxf, cfl_idt = 0;
366   PetscScalar *x, *f, *slope;
367   Vec          Xloc;
368   DM           da;
369 
370   PetscFunctionBeginUser;
371   PetscCall(TSGetDM(ts, &da));
372   PetscCall(DMGetLocalVector(da, &Xloc));
373   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
374   hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm;
375   hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm);
376   hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf);
377   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
378   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
379   PetscCall(VecZeroEntries(F));
380   PetscCall(DMDAVecGetArray(da, Xloc, &x));
381   PetscCall(VecGetArray(F, &f));
382   PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
383   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
384 
385   if (ctx->bctype == FVBC_OUTFLOW) {
386     for (i = xs - 2; i < 0; i++) {
387       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
388     }
389     for (i = Mx; i < xs + xm + 2; i++) {
390       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
391     }
392   }
393   for (i = xs - 1; i < xs + xm + 1; i++) {
394     struct _LimitInfo info;
395     PetscScalar      *cjmpL, *cjmpR;
396     if (i < sm - lsbwidth + 1 || i > ms + rsbwidth - 2) { /* slow components and the first and last fast components */
397       /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
398       PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
399       /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
400       PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
401       cjmpL = &ctx->cjmpLR[0];
402       cjmpR = &ctx->cjmpLR[dof];
403       for (j = 0; j < dof; j++) {
404         PetscScalar jmpL, jmpR;
405         jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
406         jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
407         for (k = 0; k < dof; k++) {
408           cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
409           cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
410         }
411       }
412       /* Apply limiter to the left and right characteristic jumps */
413       info.m   = dof;
414       info.hxs = hxs;
415       info.hxm = hxm;
416       info.hxf = hxf;
417       (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope);
418       for (j = 0; j < dof; j++) {
419         PetscScalar tmp = 0;
420         for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
421         slope[i * dof + j] = tmp;
422       }
423     }
424   }
425 
426   for (i = xs; i < xs + xm + 1; i++) {
427     PetscReal    maxspeed;
428     PetscScalar *uL, *uR;
429     uL = &ctx->uLR[0];
430     uR = &ctx->uLR[dof];
431     if (i < sm - lsbwidth) { /* slow region */
432       for (j = 0; j < dof; j++) {
433         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
434         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
435       }
436       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
437       cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hxs)); /* Max allowable value of 1/Delta t */
438       if (i > xs) {
439         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
440       }
441       if (i < xs + xm) {
442         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
443         islow++;
444       }
445     }
446     if (i == sm - lsbwidth) { /* interface between slow and medium regions */
447       for (j = 0; j < dof; j++) {
448         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
449         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
450       }
451       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
452       if (i > xs) {
453         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
454       }
455     }
456     if (i == ms + rsbwidth) { /* interface between medium and slow regions */
457       for (j = 0; j < dof; j++) {
458         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
459         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
460       }
461       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
462       if (i < xs + xm) {
463         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
464         islow++;
465       }
466     }
467     if (i > ms + rsbwidth) { /* slow region */
468       for (j = 0; j < dof; j++) {
469         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
470         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
471       }
472       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
473       cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hxs)); /* Max allowable value of 1/Delta t */
474       if (i > xs) {
475         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
476       }
477       if (i < xs + xm) {
478         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
479         islow++;
480       }
481     }
482   }
483   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
484   PetscCall(VecRestoreArray(F, &f));
485   PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
486   PetscCall(DMRestoreLocalVector(da, &Xloc));
487   PetscCallMPI(MPI_Allreduce(&cfl_idt, &ctx->cfl_idt, 1, MPIU_SCALAR, MPIU_MAX, PetscObjectComm((PetscObject)da)));
488   PetscFunctionReturn(0);
489 }
490 
491 PetscErrorCode FVRHSFunctionslowbuffer_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) {
492   FVCtx       *ctx = (FVCtx *)vctx;
493   PetscInt     i, j, k, Mx, dof, xs, xm, islowbuffer = 0, sm = ctx->sm, ms = ctx->ms, lsbwidth = ctx->lsbwidth, rsbwidth = ctx->rsbwidth;
494   PetscReal    hxs, hxm, hxf;
495   PetscScalar *x, *f, *slope;
496   Vec          Xloc;
497   DM           da;
498 
499   PetscFunctionBeginUser;
500   PetscCall(TSGetDM(ts, &da));
501   PetscCall(DMGetLocalVector(da, &Xloc));
502   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
503   hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm;
504   hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm);
505   hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf);
506   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
507   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
508   PetscCall(VecZeroEntries(F));
509   PetscCall(DMDAVecGetArray(da, Xloc, &x));
510   PetscCall(VecGetArray(F, &f));
511   PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
512   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
513 
514   if (ctx->bctype == FVBC_OUTFLOW) {
515     for (i = xs - 2; i < 0; i++) {
516       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
517     }
518     for (i = Mx; i < xs + xm + 2; i++) {
519       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
520     }
521   }
522   for (i = xs - 1; i < xs + xm + 1; i++) {
523     struct _LimitInfo info;
524     PetscScalar      *cjmpL, *cjmpR;
525     if ((i > sm - lsbwidth - 2 && i < sm + 1) || (i > ms - 2 && i < ms + rsbwidth + 1)) {
526       /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
527       PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
528       /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
529       PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
530       cjmpL = &ctx->cjmpLR[0];
531       cjmpR = &ctx->cjmpLR[dof];
532       for (j = 0; j < dof; j++) {
533         PetscScalar jmpL, jmpR;
534         jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
535         jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
536         for (k = 0; k < dof; k++) {
537           cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
538           cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
539         }
540       }
541       /* Apply limiter to the left and right characteristic jumps */
542       info.m   = dof;
543       info.hxs = hxs;
544       info.hxm = hxm;
545       info.hxf = hxf;
546       (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope);
547       for (j = 0; j < dof; j++) {
548         PetscScalar tmp = 0;
549         for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
550         slope[i * dof + j] = tmp;
551       }
552     }
553   }
554 
555   for (i = xs; i < xs + xm + 1; i++) {
556     PetscReal    maxspeed;
557     PetscScalar *uL, *uR;
558     uL = &ctx->uLR[0];
559     uR = &ctx->uLR[dof];
560     if (i == sm - lsbwidth) {
561       for (j = 0; j < dof; j++) {
562         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
563         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
564       }
565       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
566       if (i < xs + xm) {
567         for (j = 0; j < dof; j++) f[islowbuffer * dof + j] += ctx->flux[j] / hxs;
568         islowbuffer++;
569       }
570     }
571     if (i > sm - lsbwidth && i < sm) {
572       for (j = 0; j < dof; j++) {
573         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
574         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
575       }
576       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
577       if (i > xs) {
578         for (j = 0; j < dof; j++) f[(islowbuffer - 1) * dof + j] -= ctx->flux[j] / hxs;
579       }
580       if (i < xs + xm) {
581         for (j = 0; j < dof; j++) f[islowbuffer * dof + j] += ctx->flux[j] / hxs;
582         islowbuffer++;
583       }
584     }
585     if (i == sm) { /* interface between the slow region and the medium region */
586       for (j = 0; j < dof; j++) {
587         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
588         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
589       }
590       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
591       if (i > xs) {
592         for (j = 0; j < dof; j++) f[(islowbuffer - 1) * dof + j] -= ctx->flux[j] / hxs;
593       }
594     }
595     if (i == ms) { /* interface between the medium region and the slow region */
596       for (j = 0; j < dof; j++) {
597         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
598         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
599       }
600       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
601       if (i < xs + xm) {
602         for (j = 0; j < dof; j++) f[islowbuffer * dof + j] += ctx->flux[j] / hxs;
603         islowbuffer++;
604       }
605     }
606     if (i > ms && i < ms + rsbwidth) {
607       for (j = 0; j < dof; j++) {
608         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
609         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
610       }
611       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
612       if (i > xs) {
613         for (j = 0; j < dof; j++) f[(islowbuffer - 1) * dof + j] -= ctx->flux[j] / hxs;
614       }
615       if (i < xs + xm) {
616         for (j = 0; j < dof; j++) f[islowbuffer * dof + j] += ctx->flux[j] / hxs;
617         islowbuffer++;
618       }
619     }
620     if (i == ms + rsbwidth) {
621       for (j = 0; j < dof; j++) {
622         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
623         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
624       }
625       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
626       if (i > xs) {
627         for (j = 0; j < dof; j++) f[(islowbuffer - 1) * dof + j] -= ctx->flux[j] / hxs;
628       }
629     }
630   }
631   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
632   PetscCall(VecRestoreArray(F, &f));
633   PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
634   PetscCall(DMRestoreLocalVector(da, &Xloc));
635   PetscFunctionReturn(0);
636 }
637 
638 /* --------------------------------- Finite Volume Solver for medium components ----------------------------------- */
639 PetscErrorCode FVRHSFunctionmedium_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) {
640   FVCtx       *ctx = (FVCtx *)vctx;
641   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;
642   PetscReal    hxs, hxm, hxf;
643   PetscScalar *x, *f, *slope;
644   Vec          Xloc;
645   DM           da;
646 
647   PetscFunctionBeginUser;
648   PetscCall(TSGetDM(ts, &da));
649   PetscCall(DMGetLocalVector(da, &Xloc));
650   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
651   hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm;
652   hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm);
653   hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf);
654   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
655   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
656   PetscCall(VecZeroEntries(F));
657   PetscCall(DMDAVecGetArray(da, Xloc, &x));
658   PetscCall(VecGetArray(F, &f));
659   PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
660   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
661 
662   if (ctx->bctype == FVBC_OUTFLOW) {
663     for (i = xs - 2; i < 0; i++) {
664       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
665     }
666     for (i = Mx; i < xs + xm + 2; i++) {
667       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
668     }
669   }
670   for (i = xs - 1; i < xs + xm + 1; i++) {
671     struct _LimitInfo info;
672     PetscScalar      *cjmpL, *cjmpR;
673     if ((i > sm - 2 && i < mf - lmbwidth + 1) || (i > fm + rmbwidth - 2 && i < ms + 1)) { /* slow components and the first and last fast components */
674       /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
675       PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
676       /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
677       PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
678       cjmpL = &ctx->cjmpLR[0];
679       cjmpR = &ctx->cjmpLR[dof];
680       for (j = 0; j < dof; j++) {
681         PetscScalar jmpL, jmpR;
682         jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
683         jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
684         for (k = 0; k < dof; k++) {
685           cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
686           cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
687         }
688       }
689       /* Apply limiter to the left and right characteristic jumps */
690       info.m   = dof;
691       info.hxs = hxs;
692       info.hxm = hxm;
693       info.hxf = hxf;
694       (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope);
695       for (j = 0; j < dof; j++) {
696         PetscScalar tmp = 0;
697         for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
698         slope[i * dof + j] = tmp;
699       }
700     }
701   }
702 
703   for (i = xs; i < xs + xm + 1; i++) {
704     PetscReal    maxspeed;
705     PetscScalar *uL, *uR;
706     uL = &ctx->uLR[0];
707     uR = &ctx->uLR[dof];
708     if (i == sm) { /* interface between slow and medium regions */
709       for (j = 0; j < dof; j++) {
710         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
711         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
712       }
713       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
714       if (i < xs + xm) {
715         for (j = 0; j < dof; j++) f[imedium * dof + j] += ctx->flux[j] / hxm;
716         imedium++;
717       }
718     }
719     if (i > sm && i < mf - lmbwidth) { /* medium region */
720       for (j = 0; j < dof; j++) {
721         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
722         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
723       }
724       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
725       if (i > xs) {
726         for (j = 0; j < dof; j++) f[(imedium - 1) * dof + j] -= ctx->flux[j] / hxm;
727       }
728       if (i < xs + xm) {
729         for (j = 0; j < dof; j++) f[imedium * dof + j] += ctx->flux[j] / hxm;
730         imedium++;
731       }
732     }
733     if (i == mf - lmbwidth) { /* interface between medium and fast regions */
734       for (j = 0; j < dof; j++) {
735         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
736         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
737       }
738       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
739       if (i > xs) {
740         for (j = 0; j < dof; j++) f[(imedium - 1) * dof + j] -= ctx->flux[j] / hxm;
741       }
742     }
743     if (i == fm + rmbwidth) { /* interface between fast and medium regions */
744       for (j = 0; j < dof; j++) {
745         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
746         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
747       }
748       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
749       if (i < xs + xm) {
750         for (j = 0; j < dof; j++) f[imedium * dof + j] += ctx->flux[j] / hxm;
751         imedium++;
752       }
753     }
754     if (i > fm + rmbwidth && i < ms) { /* medium region */
755       for (j = 0; j < dof; j++) {
756         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
757         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
758       }
759       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
760       if (i > xs) {
761         for (j = 0; j < dof; j++) f[(imedium - 1) * dof + j] -= ctx->flux[j] / hxm;
762       }
763       if (i < xs + xm) {
764         for (j = 0; j < dof; j++) f[imedium * dof + j] += ctx->flux[j] / hxm;
765         imedium++;
766       }
767     }
768     if (i == ms) { /* interface between medium and slow regions */
769       for (j = 0; j < dof; j++) {
770         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
771         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
772       }
773       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
774       if (i > xs) {
775         for (j = 0; j < dof; j++) f[(imedium - 1) * dof + j] -= ctx->flux[j] / hxm;
776       }
777     }
778   }
779   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
780   PetscCall(VecRestoreArray(F, &f));
781   PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
782   PetscCall(DMRestoreLocalVector(da, &Xloc));
783   PetscFunctionReturn(0);
784 }
785 
786 PetscErrorCode FVRHSFunctionmediumbuffer_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) {
787   FVCtx       *ctx = (FVCtx *)vctx;
788   PetscInt     i, j, k, Mx, dof, xs, xm, imediumbuffer = 0, mf = ctx->mf, fm = ctx->fm, lmbwidth = ctx->lmbwidth, rmbwidth = ctx->rmbwidth;
789   PetscReal    hxs, hxm, hxf;
790   PetscScalar *x, *f, *slope;
791   Vec          Xloc;
792   DM           da;
793 
794   PetscFunctionBeginUser;
795   PetscCall(TSGetDM(ts, &da));
796   PetscCall(DMGetLocalVector(da, &Xloc));
797   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
798   hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm;
799   hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm);
800   hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf);
801   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
802   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
803   PetscCall(VecZeroEntries(F));
804   PetscCall(DMDAVecGetArray(da, Xloc, &x));
805   PetscCall(VecGetArray(F, &f));
806   PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
807   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
808 
809   if (ctx->bctype == FVBC_OUTFLOW) {
810     for (i = xs - 2; i < 0; i++) {
811       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
812     }
813     for (i = Mx; i < xs + xm + 2; i++) {
814       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
815     }
816   }
817   for (i = xs - 1; i < xs + xm + 1; i++) {
818     struct _LimitInfo info;
819     PetscScalar      *cjmpL, *cjmpR;
820     if ((i > mf - lmbwidth - 2 && i < mf + 1) || (i > fm - 2 && i < fm + rmbwidth + 1)) {
821       /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
822       PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
823       /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
824       PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
825       cjmpL = &ctx->cjmpLR[0];
826       cjmpR = &ctx->cjmpLR[dof];
827       for (j = 0; j < dof; j++) {
828         PetscScalar jmpL, jmpR;
829         jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
830         jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
831         for (k = 0; k < dof; k++) {
832           cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
833           cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
834         }
835       }
836       /* Apply limiter to the left and right characteristic jumps */
837       info.m   = dof;
838       info.hxs = hxs;
839       info.hxm = hxm;
840       info.hxf = hxf;
841       (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope);
842       for (j = 0; j < dof; j++) {
843         PetscScalar tmp = 0;
844         for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
845         slope[i * dof + j] = tmp;
846       }
847     }
848   }
849 
850   for (i = xs; i < xs + xm + 1; i++) {
851     PetscReal    maxspeed;
852     PetscScalar *uL, *uR;
853     uL = &ctx->uLR[0];
854     uR = &ctx->uLR[dof];
855     if (i == mf - lmbwidth) {
856       for (j = 0; j < dof; j++) {
857         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
858         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
859       }
860       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
861       if (i < xs + xm) {
862         for (j = 0; j < dof; j++) f[imediumbuffer * dof + j] += ctx->flux[j] / hxm;
863         imediumbuffer++;
864       }
865     }
866     if (i > mf - lmbwidth && i < mf) {
867       for (j = 0; j < dof; j++) {
868         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
869         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
870       }
871       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
872       if (i > xs) {
873         for (j = 0; j < dof; j++) f[(imediumbuffer - 1) * dof + j] -= ctx->flux[j] / hxm;
874       }
875       if (i < xs + xm) {
876         for (j = 0; j < dof; j++) f[imediumbuffer * dof + j] += ctx->flux[j] / hxm;
877         imediumbuffer++;
878       }
879     }
880     if (i == mf) { /* interface between the medium region and the fast region */
881       for (j = 0; j < dof; j++) {
882         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
883         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
884       }
885       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
886       if (i > xs) {
887         for (j = 0; j < dof; j++) f[(imediumbuffer - 1) * dof + j] -= ctx->flux[j] / hxm;
888       }
889     }
890     if (i == fm) { /* interface between the fast region and the medium region */
891       for (j = 0; j < dof; j++) {
892         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
893         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
894       }
895       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
896       if (i < xs + xm) {
897         for (j = 0; j < dof; j++) f[imediumbuffer * dof + j] += ctx->flux[j] / hxm;
898         imediumbuffer++;
899       }
900     }
901     if (i > fm && i < fm + rmbwidth) {
902       for (j = 0; j < dof; j++) {
903         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
904         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
905       }
906       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
907       if (i > xs) {
908         for (j = 0; j < dof; j++) f[(imediumbuffer - 1) * dof + j] -= ctx->flux[j] / hxm;
909       }
910       if (i < xs + xm) {
911         for (j = 0; j < dof; j++) f[imediumbuffer * dof + j] += ctx->flux[j] / hxm;
912         imediumbuffer++;
913       }
914     }
915     if (i == fm + rmbwidth) {
916       for (j = 0; j < dof; j++) {
917         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
918         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
919       }
920       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
921       if (i > xs) {
922         for (j = 0; j < dof; j++) f[(imediumbuffer - 1) * dof + j] -= ctx->flux[j] / hxm;
923       }
924     }
925   }
926   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
927   PetscCall(VecRestoreArray(F, &f));
928   PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
929   PetscCall(DMRestoreLocalVector(da, &Xloc));
930   PetscFunctionReturn(0);
931 }
932 
933 /* --------------------------------- Finite Volume Solver for fast  parts ----------------------------------- */
934 PetscErrorCode FVRHSFunctionfast_3WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) {
935   FVCtx       *ctx = (FVCtx *)vctx;
936   PetscInt     i, j, k, Mx, dof, xs, xm, ifast = 0, mf = ctx->mf, fm = ctx->fm;
937   PetscReal    hxs, hxm, hxf;
938   PetscScalar *x, *f, *slope;
939   Vec          Xloc;
940   DM           da;
941 
942   PetscFunctionBeginUser;
943   PetscCall(TSGetDM(ts, &da));
944   PetscCall(DMGetLocalVector(da, &Xloc));
945   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
946   hxs = (ctx->xmax - ctx->xmin) / 8.0 / ctx->sm;
947   hxm = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->mf - ctx->sm);
948   hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fm - ctx->mf);
949   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
950   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
951   PetscCall(VecZeroEntries(F));
952   PetscCall(DMDAVecGetArray(da, Xloc, &x));
953   PetscCall(VecGetArray(F, &f));
954   PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
955   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
956 
957   if (ctx->bctype == FVBC_OUTFLOW) {
958     for (i = xs - 2; i < 0; i++) {
959       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
960     }
961     for (i = Mx; i < xs + xm + 2; i++) {
962       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
963     }
964   }
965   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 */
966     struct _LimitInfo info;
967     PetscScalar      *cjmpL, *cjmpR;
968     if (i > mf - 2 && i < fm + 1) {
969       PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
970       PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
971       cjmpL = &ctx->cjmpLR[0];
972       cjmpR = &ctx->cjmpLR[dof];
973       for (j = 0; j < dof; j++) {
974         PetscScalar jmpL, jmpR;
975         jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
976         jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
977         for (k = 0; k < dof; k++) {
978           cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
979           cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
980         }
981       }
982       /* Apply limiter to the left and right characteristic jumps */
983       info.m   = dof;
984       info.hxs = hxs;
985       info.hxm = hxm;
986       info.hxf = hxf;
987       (*ctx->limit3)(&info, cjmpL, cjmpR, ctx->sm, ctx->mf, ctx->fm, ctx->ms, i, ctx->cslope);
988       for (j = 0; j < dof; j++) {
989         PetscScalar tmp = 0;
990         for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
991         slope[i * dof + j] = tmp;
992       }
993     }
994   }
995 
996   for (i = xs; i < xs + xm + 1; i++) {
997     PetscReal    maxspeed;
998     PetscScalar *uL, *uR;
999     uL = &ctx->uLR[0];
1000     uR = &ctx->uLR[dof];
1001     if (i == mf) { /* interface between medium and fast regions */
1002       for (j = 0; j < dof; j++) {
1003         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxm / 2;
1004         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
1005       }
1006       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
1007       if (i < xs + xm) {
1008         for (j = 0; j < dof; j++) f[ifast * dof + j] += ctx->flux[j] / hxf;
1009         ifast++;
1010       }
1011     }
1012     if (i > mf && i < fm) { /* fast region */
1013       for (j = 0; j < dof; j++) {
1014         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
1015         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
1016       }
1017       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
1018       if (i > xs) {
1019         for (j = 0; j < dof; j++) f[(ifast - 1) * dof + j] -= ctx->flux[j] / hxf;
1020       }
1021       if (i < xs + xm) {
1022         for (j = 0; j < dof; j++) f[ifast * dof + j] += ctx->flux[j] / hxf;
1023         ifast++;
1024       }
1025     }
1026     if (i == fm) { /* interface between fast and medium regions */
1027       for (j = 0; j < dof; j++) {
1028         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
1029         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxm / 2;
1030       }
1031       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
1032       if (i > xs) {
1033         for (j = 0; j < dof; j++) f[(ifast - 1) * dof + j] -= ctx->flux[j] / hxf;
1034       }
1035     }
1036   }
1037   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
1038   PetscCall(VecRestoreArray(F, &f));
1039   PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
1040   PetscCall(DMRestoreLocalVector(da, &Xloc));
1041   PetscFunctionReturn(0);
1042 }
1043 
1044 int main(int argc, char *argv[]) {
1045   char              lname[256] = "mc", physname[256] = "advect", final_fname[256] = "solution.m";
1046   PetscFunctionList limiters = 0, physics = 0;
1047   MPI_Comm          comm;
1048   TS                ts;
1049   DM                da;
1050   Vec               X, X0, R;
1051   FVCtx             ctx;
1052   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;
1053   PetscBool         view_final = PETSC_FALSE;
1054   PetscReal         ptime;
1055 
1056   PetscFunctionBeginUser;
1057   PetscCall(PetscInitialize(&argc, &argv, 0, help));
1058   comm = PETSC_COMM_WORLD;
1059   PetscCall(PetscMemzero(&ctx, sizeof(ctx)));
1060 
1061   /* Register limiters to be available on the command line */
1062   PetscCall(PetscFunctionListAdd(&limiters, "upwind", Limit3_Upwind));
1063   PetscCall(PetscFunctionListAdd(&limiters, "lax-wendroff", Limit3_LaxWendroff));
1064   PetscCall(PetscFunctionListAdd(&limiters, "beam-warming", Limit3_BeamWarming));
1065   PetscCall(PetscFunctionListAdd(&limiters, "fromm", Limit3_Fromm));
1066   PetscCall(PetscFunctionListAdd(&limiters, "minmod", Limit3_Minmod));
1067   PetscCall(PetscFunctionListAdd(&limiters, "superbee", Limit3_Superbee));
1068   PetscCall(PetscFunctionListAdd(&limiters, "mc", Limit3_MC));
1069   PetscCall(PetscFunctionListAdd(&limiters, "koren3", Limit3_Koren3));
1070 
1071   /* Register physical models to be available on the command line */
1072   PetscCall(PetscFunctionListAdd(&physics, "advect", PhysicsCreate_Advect));
1073 
1074   ctx.comm   = comm;
1075   ctx.cfl    = 0.9;
1076   ctx.bctype = FVBC_PERIODIC;
1077   ctx.xmin   = -1.0;
1078   ctx.xmax   = 1.0;
1079   PetscOptionsBegin(comm, NULL, "Finite Volume solver options", "");
1080   PetscCall(PetscOptionsReal("-xmin", "X min", "", ctx.xmin, &ctx.xmin, NULL));
1081   PetscCall(PetscOptionsReal("-xmax", "X max", "", ctx.xmax, &ctx.xmax, NULL));
1082   PetscCall(PetscOptionsFList("-limit", "Name of flux imiter to use", "", limiters, lname, lname, sizeof(lname), NULL));
1083   PetscCall(PetscOptionsInt("-draw", "Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)", "", draw, &draw, NULL));
1084   PetscCall(PetscOptionsString("-view_final", "Write final solution in ASCII MATLAB format to given file name", "", final_fname, final_fname, sizeof(final_fname), &view_final));
1085   PetscCall(PetscOptionsInt("-initial", "Initial condition (depends on the physics)", "", ctx.initial, &ctx.initial, NULL));
1086   PetscCall(PetscOptionsBool("-exact", "Compare errors with exact solution", "", ctx.exact, &ctx.exact, NULL));
1087   PetscCall(PetscOptionsBool("-simulation", "Compare errors with reference solution", "", ctx.simulation, &ctx.simulation, NULL));
1088   PetscCall(PetscOptionsReal("-cfl", "CFL number to time step at", "", ctx.cfl, &ctx.cfl, NULL));
1089   PetscCall(PetscOptionsEnum("-bc_type", "Boundary condition", "", FVBCTypes, (PetscEnum)ctx.bctype, (PetscEnum *)&ctx.bctype, NULL));
1090   PetscCall(PetscOptionsInt("-hratio", "Spacing ratio", "", ctx.hratio, &ctx.hratio, NULL));
1091   PetscOptionsEnd();
1092 
1093   /* Choose the limiter from the list of registered limiters */
1094   PetscCall(PetscFunctionListFind(limiters, lname, &ctx.limit3));
1095   PetscCheck(ctx.limit3, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Limiter '%s' not found", lname);
1096 
1097   /* Choose the physics from the list of registered models */
1098   {
1099     PetscErrorCode (*r)(FVCtx *);
1100     PetscCall(PetscFunctionListFind(physics, physname, &r));
1101     PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Physics '%s' not found", physname);
1102     /* Create the physics, will set the number of fields and their names */
1103     PetscCall((*r)(&ctx));
1104   }
1105 
1106   /* Create a DMDA to manage the parallel grid */
1107   PetscCall(DMDACreate1d(comm, DM_BOUNDARY_PERIODIC, 50, ctx.physics2.dof, 2, NULL, &da));
1108   PetscCall(DMSetFromOptions(da));
1109   PetscCall(DMSetUp(da));
1110   /* Inform the DMDA of the field names provided by the physics. */
1111   /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */
1112   for (i = 0; i < ctx.physics2.dof; i++) PetscCall(DMDASetFieldName(da, i, ctx.physics2.fieldname[i]));
1113   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
1114   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
1115 
1116   /* Set coordinates of cell centers */
1117   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));
1118 
1119   /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */
1120   PetscCall(PetscMalloc4(dof * dof, &ctx.R, dof * dof, &ctx.Rinv, 2 * dof, &ctx.cjmpLR, 1 * dof, &ctx.cslope));
1121   PetscCall(PetscMalloc3(2 * dof, &ctx.uLR, dof, &ctx.flux, dof, &ctx.speeds));
1122 
1123   /* Create a vector to store the solution and to save the initial state */
1124   PetscCall(DMCreateGlobalVector(da, &X));
1125   PetscCall(VecDuplicate(X, &X0));
1126   PetscCall(VecDuplicate(X, &R));
1127 
1128   /* create index for slow parts and fast parts,
1129      count_slow + count_fast = Mx, counts_slow*hs = 0.5, counts_fast*hf = 0.5 */
1130   count_slow   = Mx / (1 + ctx.hratio) / (1 + ctx.hratio);
1131   count_medium = 2 * ctx.hratio * count_slow;
1132   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");
1133   count_fast = ctx.hratio * ctx.hratio * count_slow;
1134   ctx.sm     = count_slow / 2;
1135   ctx.mf     = ctx.sm + count_medium / 2;
1136   ctx.fm     = ctx.mf + count_fast;
1137   ctx.ms     = ctx.fm + count_medium / 2;
1138   PetscCall(PetscMalloc1(xm * dof, &index_slow));
1139   PetscCall(PetscMalloc1(xm * dof, &index_medium));
1140   PetscCall(PetscMalloc1(xm * dof, &index_fast));
1141   PetscCall(PetscMalloc1(6 * dof, &index_slowbuffer));
1142   PetscCall(PetscMalloc1(6 * dof, &index_mediumbuffer));
1143   if (((AdvectCtx *)ctx.physics2.user)->a > 0) {
1144     ctx.lsbwidth = 2;
1145     ctx.rsbwidth = 4;
1146     ctx.lmbwidth = 2;
1147     ctx.rmbwidth = 4;
1148   } else {
1149     ctx.lsbwidth = 4;
1150     ctx.rsbwidth = 2;
1151     ctx.lmbwidth = 4;
1152     ctx.rmbwidth = 2;
1153   }
1154 
1155   for (i = xs; i < xs + xm; i++) {
1156     if (i < ctx.sm - ctx.lsbwidth || i > ctx.ms + ctx.rsbwidth - 1)
1157       for (k = 0; k < dof; k++) index_slow[islow++] = i * dof + k;
1158     else if ((i >= ctx.sm - ctx.lsbwidth && i < ctx.sm) || (i > ctx.ms - 1 && i <= ctx.ms + ctx.rsbwidth - 1))
1159       for (k = 0; k < dof; k++) index_slowbuffer[islowbuffer++] = i * dof + k;
1160     else if (i < ctx.mf - ctx.lmbwidth || i > ctx.fm + ctx.rmbwidth - 1)
1161       for (k = 0; k < dof; k++) index_medium[imedium++] = i * dof + k;
1162     else if ((i >= ctx.mf - ctx.lmbwidth && i < ctx.mf) || (i > ctx.fm - 1 && i <= ctx.fm + ctx.rmbwidth - 1))
1163       for (k = 0; k < dof; k++) index_mediumbuffer[imediumbuffer++] = i * dof + k;
1164     else
1165       for (k = 0; k < dof; k++) index_fast[ifast++] = i * dof + k;
1166   }
1167   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islow, index_slow, PETSC_COPY_VALUES, &ctx.iss));
1168   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, imedium, index_medium, PETSC_COPY_VALUES, &ctx.ism));
1169   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, ifast, index_fast, PETSC_COPY_VALUES, &ctx.isf));
1170   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islowbuffer, index_slowbuffer, PETSC_COPY_VALUES, &ctx.issb));
1171   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, imediumbuffer, index_mediumbuffer, PETSC_COPY_VALUES, &ctx.ismb));
1172 
1173   /* Create a time-stepping object */
1174   PetscCall(TSCreate(comm, &ts));
1175   PetscCall(TSSetDM(ts, da));
1176   PetscCall(TSSetRHSFunction(ts, R, FVRHSFunction_3WaySplit, &ctx));
1177   PetscCall(TSRHSSplitSetIS(ts, "slow", ctx.iss));
1178   PetscCall(TSRHSSplitSetIS(ts, "medium", ctx.ism));
1179   PetscCall(TSRHSSplitSetIS(ts, "fast", ctx.isf));
1180   PetscCall(TSRHSSplitSetIS(ts, "slowbuffer", ctx.issb));
1181   PetscCall(TSRHSSplitSetIS(ts, "mediumbuffer", ctx.ismb));
1182   PetscCall(TSRHSSplitSetRHSFunction(ts, "slow", NULL, FVRHSFunctionslow_3WaySplit, &ctx));
1183   PetscCall(TSRHSSplitSetRHSFunction(ts, "medium", NULL, FVRHSFunctionmedium_3WaySplit, &ctx));
1184   PetscCall(TSRHSSplitSetRHSFunction(ts, "fast", NULL, FVRHSFunctionfast_3WaySplit, &ctx));
1185   PetscCall(TSRHSSplitSetRHSFunction(ts, "slowbuffer", NULL, FVRHSFunctionslowbuffer_3WaySplit, &ctx));
1186   PetscCall(TSRHSSplitSetRHSFunction(ts, "mediumbuffer", NULL, FVRHSFunctionmediumbuffer_3WaySplit, &ctx));
1187 
1188   PetscCall(TSSetType(ts, TSSSP));
1189   /*PetscCall(TSSetType(ts,TSMPRK));*/
1190   PetscCall(TSSetMaxTime(ts, 10));
1191   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1192 
1193   /* Compute initial conditions and starting time step */
1194   PetscCall(FVSample_3WaySplit(&ctx, da, 0, X0));
1195   PetscCall(FVRHSFunction_3WaySplit(ts, 0, X0, X, (void *)&ctx)); /* Initial function evaluation, only used to determine max speed */
1196   PetscCall(VecCopy(X0, X));                                      /* The function value was not used so we set X=X0 again */
1197   PetscCall(TSSetTimeStep(ts, ctx.cfl / ctx.cfl_idt));
1198   PetscCall(TSSetFromOptions(ts)); /* Take runtime options */
1199   PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD));
1200   {
1201     PetscInt           steps;
1202     PetscScalar        mass_initial, mass_final, mass_difference, mass_differenceg;
1203     const PetscScalar *ptr_X, *ptr_X0;
1204     const PetscReal    hs = (ctx.xmax - ctx.xmin) / 4.0 / count_slow;
1205     const PetscReal    hm = (ctx.xmax - ctx.xmin) / 2.0 / count_medium;
1206     const PetscReal    hf = (ctx.xmax - ctx.xmin) / 4.0 / count_fast;
1207 
1208     PetscCall(TSSolve(ts, X));
1209     PetscCall(TSGetSolveTime(ts, &ptime));
1210     PetscCall(TSGetStepNumber(ts, &steps));
1211     /* calculate the total mass at initial time and final time */
1212     mass_initial = 0.0;
1213     mass_final   = 0.0;
1214     PetscCall(DMDAVecGetArrayRead(da, X0, (void *)&ptr_X0));
1215     PetscCall(DMDAVecGetArrayRead(da, X, (void *)&ptr_X));
1216     for (i = xs; i < xs + xm; i++) {
1217       if (i < ctx.sm || i > ctx.ms - 1)
1218         for (k = 0; k < dof; k++) {
1219           mass_initial = mass_initial + hs * ptr_X0[i * dof + k];
1220           mass_final   = mass_final + hs * ptr_X[i * dof + k];
1221         }
1222       else if (i < ctx.mf || i > ctx.fm - 1)
1223         for (k = 0; k < dof; k++) {
1224           mass_initial = mass_initial + hm * ptr_X0[i * dof + k];
1225           mass_final   = mass_final + hm * ptr_X[i * dof + k];
1226         }
1227       else {
1228         for (k = 0; k < dof; k++) {
1229           mass_initial = mass_initial + hf * ptr_X0[i * dof + k];
1230           mass_final   = mass_final + hf * ptr_X[i * dof + k];
1231         }
1232       }
1233     }
1234     PetscCall(DMDAVecRestoreArrayRead(da, X0, (void *)&ptr_X0));
1235     PetscCall(DMDAVecRestoreArrayRead(da, X, (void *)&ptr_X));
1236     mass_difference = mass_final - mass_initial;
1237     PetscCallMPI(MPI_Allreduce(&mass_difference, &mass_differenceg, 1, MPIU_SCALAR, MPIU_SUM, comm));
1238     PetscCall(PetscPrintf(comm, "Mass difference %g\n", (double)mass_differenceg));
1239     PetscCall(PetscPrintf(comm, "Final time %g, steps %" PetscInt_FMT "\n", (double)ptime, steps));
1240     PetscCall(PetscPrintf(comm, "Maximum allowable stepsize according to CFL %g\n", (double)(1.0 / ctx.cfl_idt)));
1241     if (ctx.exact) {
1242       PetscReal nrm1 = 0;
1243       PetscCall(SolutionErrorNorms_3WaySplit(&ctx, da, ptime, X, &nrm1));
1244       PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %g\n", (double)nrm1));
1245     }
1246     if (ctx.simulation) {
1247       PetscReal          nrm1 = 0;
1248       PetscViewer        fd;
1249       char               filename[PETSC_MAX_PATH_LEN] = "binaryoutput";
1250       Vec                XR;
1251       PetscBool          flg;
1252       const PetscScalar *ptr_XR;
1253       PetscCall(PetscOptionsGetString(NULL, NULL, "-f", filename, sizeof(filename), &flg));
1254       PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate binary file with the -f option");
1255       PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, filename, FILE_MODE_READ, &fd));
1256       PetscCall(VecDuplicate(X0, &XR));
1257       PetscCall(VecLoad(XR, fd));
1258       PetscCall(PetscViewerDestroy(&fd));
1259       PetscCall(VecGetArrayRead(X, &ptr_X));
1260       PetscCall(VecGetArrayRead(XR, &ptr_XR));
1261       for (i = xs; i < xs + xm; i++) {
1262         if (i < ctx.sm || i > ctx.ms - 1)
1263           for (k = 0; k < dof; k++) nrm1 = nrm1 + hs * PetscAbs(ptr_X[i * dof + k] - ptr_XR[i * dof + k]);
1264         else if (i < ctx.mf || i < ctx.fm - 1)
1265           for (k = 0; k < dof; k++) nrm1 = nrm1 + hm * PetscAbs(ptr_X[i * dof + k] - ptr_XR[i * dof + k]);
1266         else
1267           for (k = 0; k < dof; k++) nrm1 = nrm1 + hf * PetscAbs(ptr_X[i * dof + k] - ptr_XR[i * dof + k]);
1268       }
1269       PetscCall(VecRestoreArrayRead(X, &ptr_X));
1270       PetscCall(VecRestoreArrayRead(XR, &ptr_XR));
1271       PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %g\n", (double)nrm1));
1272       PetscCall(VecDestroy(&XR));
1273     }
1274   }
1275 
1276   PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD));
1277   if (draw & 0x1) PetscCall(VecView(X0, PETSC_VIEWER_DRAW_WORLD));
1278   if (draw & 0x2) PetscCall(VecView(X, PETSC_VIEWER_DRAW_WORLD));
1279   if (draw & 0x4) {
1280     Vec Y;
1281     PetscCall(VecDuplicate(X, &Y));
1282     PetscCall(FVSample_3WaySplit(&ctx, da, ptime, Y));
1283     PetscCall(VecAYPX(Y, -1, X));
1284     PetscCall(VecView(Y, PETSC_VIEWER_DRAW_WORLD));
1285     PetscCall(VecDestroy(&Y));
1286   }
1287 
1288   if (view_final) {
1289     PetscViewer viewer;
1290     PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, final_fname, &viewer));
1291     PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB));
1292     PetscCall(VecView(X, viewer));
1293     PetscCall(PetscViewerPopFormat(viewer));
1294     PetscCall(PetscViewerDestroy(&viewer));
1295   }
1296 
1297   /* Clean up */
1298   PetscCall((*ctx.physics2.destroy)(ctx.physics2.user));
1299   for (i = 0; i < ctx.physics2.dof; i++) PetscCall(PetscFree(ctx.physics2.fieldname[i]));
1300   PetscCall(PetscFree4(ctx.R, ctx.Rinv, ctx.cjmpLR, ctx.cslope));
1301   PetscCall(PetscFree3(ctx.uLR, ctx.flux, ctx.speeds));
1302   PetscCall(VecDestroy(&X));
1303   PetscCall(VecDestroy(&X0));
1304   PetscCall(VecDestroy(&R));
1305   PetscCall(DMDestroy(&da));
1306   PetscCall(TSDestroy(&ts));
1307   PetscCall(ISDestroy(&ctx.iss));
1308   PetscCall(ISDestroy(&ctx.ism));
1309   PetscCall(ISDestroy(&ctx.isf));
1310   PetscCall(ISDestroy(&ctx.issb));
1311   PetscCall(ISDestroy(&ctx.ismb));
1312   PetscCall(PetscFree(index_slow));
1313   PetscCall(PetscFree(index_medium));
1314   PetscCall(PetscFree(index_fast));
1315   PetscCall(PetscFree(index_slowbuffer));
1316   PetscCall(PetscFree(index_mediumbuffer));
1317   PetscCall(PetscFunctionListDestroy(&limiters));
1318   PetscCall(PetscFunctionListDestroy(&physics));
1319   PetscCall(PetscFinalize());
1320   return 0;
1321 }
1322 
1323 /*TEST
1324 
1325     build:
1326       requires: !complex
1327       depends: finitevolume1d.c
1328 
1329     test:
1330       suffix: 1
1331       args: -da_grid_x 90 -initial 7 -xmin -1 -xmax 1 -hratio 2 -limit mc -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a23 -ts_use_splitrhsfunction 0
1332 
1333     test:
1334       suffix: 2
1335       args: -da_grid_x 90 -initial 7 -xmin -1 -xmax 1 -hratio 2 -limit mc -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a23 -ts_use_splitrhsfunction 1
1336 
1337 TEST*/
1338