xref: /petsc/src/ts/tutorials/multirate/ex6.c (revision a69119a591a03a9d906b29c0a4e9802e4d7c9795)
1 /*
2   Note:
3     -hratio is the ratio between mesh size of coarse grids and fine grids
4     -ts_rk_dtratio is the ratio between the large stepsize and the small stepsize
5 */
6 
7 static const char help[] = "1D periodic Finite Volume solver in slope-limiter form with semidiscrete time stepping.\n"
8                            "  advection   - Constant coefficient scalar advection\n"
9                            "                u_t       + (a*u)_x               = 0\n"
10                            "  for this toy problem, we choose different meshsizes for different sub-domains (slow-fast-slow), say\n"
11                            "                hxs  = hratio*hxf \n"
12                            "  where hxs and hxf are the grid spacings for coarse and fine grids respectively.\n"
13                            "  exact       - Exact Riemann solver which usually needs to perform a Newton iteration to connect\n"
14                            "                the states across shocks and rarefactions\n"
15                            "  simulation  - use reference solution which is generated by smaller time step size to be true solution,\n"
16                            "                also the reference solution should be generated by user and stored in a binary file.\n"
17                            "  characteristic - Limit the characteristic variables, this is usually preferred (default)\n"
18                            "Several initial conditions can be chosen with -initial N\n\n"
19                            "The problem size should be set with -da_grid_x M\n\n";
20 
21 #include <petscts.h>
22 #include <petscdm.h>
23 #include <petscdmda.h>
24 #include <petscdraw.h>
25 #include "finitevolume1d.h"
26 
27 static inline PetscReal RangeMod(PetscReal a, PetscReal xmin, PetscReal xmax) {
28   PetscReal range = xmax - xmin;
29   return xmin + PetscFmodReal(range + PetscFmodReal(a, range), range);
30 }
31 
32 /* --------------------------------- Advection ----------------------------------- */
33 typedef struct {
34   PetscReal a; /* advective velocity */
35 } AdvectCtx;
36 
37 static PetscErrorCode PhysicsRiemann_Advect(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed) {
38   AdvectCtx *ctx = (AdvectCtx *)vctx;
39   PetscReal  speed;
40 
41   PetscFunctionBeginUser;
42   speed     = ctx->a;
43   flux[0]   = PetscMax(0, speed) * uL[0] + PetscMin(0, speed) * uR[0];
44   *maxspeed = speed;
45   PetscFunctionReturn(0);
46 }
47 
48 static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx, PetscInt m, const PetscScalar *u, PetscScalar *X, PetscScalar *Xi, PetscReal *speeds) {
49   AdvectCtx *ctx = (AdvectCtx *)vctx;
50 
51   PetscFunctionBeginUser;
52   X[0]      = 1.;
53   Xi[0]     = 1.;
54   speeds[0] = ctx->a;
55   PetscFunctionReturn(0);
56 }
57 
58 static PetscErrorCode PhysicsSample_Advect(void *vctx, PetscInt initial, FVBCType bctype, PetscReal xmin, PetscReal xmax, PetscReal t, PetscReal x, PetscReal *u) {
59   AdvectCtx *ctx = (AdvectCtx *)vctx;
60   PetscReal  a   = ctx->a, x0;
61 
62   PetscFunctionBeginUser;
63   switch (bctype) {
64   case FVBC_OUTFLOW: x0 = x - a * t; break;
65   case FVBC_PERIODIC: x0 = RangeMod(x - a * t, xmin, xmax); break;
66   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown BCType");
67   }
68   switch (initial) {
69   case 0: u[0] = (x0 < 0) ? 1 : -1; break;
70   case 1: u[0] = (x0 < 0) ? -1 : 1; break;
71   case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break;
72   case 3: u[0] = PetscSinReal(2 * PETSC_PI * x0); break;
73   case 4: u[0] = PetscAbs(x0); break;
74   case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2 * PETSC_PI * x0)); break;
75   case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2 - x0 : 0)); break;
76   case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI * x0), 10.0); break;
77   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition");
78   }
79   PetscFunctionReturn(0);
80 }
81 
82 static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx) {
83   AdvectCtx *user;
84 
85   PetscFunctionBeginUser;
86   PetscCall(PetscNew(&user));
87   ctx->physics2.sample2         = PhysicsSample_Advect;
88   ctx->physics2.riemann2        = PhysicsRiemann_Advect;
89   ctx->physics2.characteristic2 = PhysicsCharacteristic_Advect;
90   ctx->physics2.destroy         = PhysicsDestroy_SimpleFree;
91   ctx->physics2.user            = user;
92   ctx->physics2.dof             = 1;
93   PetscCall(PetscStrallocpy("u", &ctx->physics2.fieldname[0]));
94   user->a = 1;
95   PetscOptionsBegin(ctx->comm, ctx->prefix, "Options for advection", "");
96   { PetscCall(PetscOptionsReal("-physics_advect_a", "Speed", "", user->a, &user->a, NULL)); }
97   PetscOptionsEnd();
98   PetscFunctionReturn(0);
99 }
100 
101 PetscErrorCode FVSample_2WaySplit(FVCtx *ctx, DM da, PetscReal time, Vec U) {
102   PetscScalar   *u, *uj, xj, xi;
103   PetscInt       i, j, k, dof, xs, xm, Mx;
104   const PetscInt N = 200;
105   PetscReal      hs, hf;
106 
107   PetscFunctionBeginUser;
108   PetscCheck(ctx->physics2.sample2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Physics has not provided a sampling function");
109   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
110   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
111   PetscCall(DMDAVecGetArray(da, U, &u));
112   PetscCall(PetscMalloc1(dof, &uj));
113   hs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
114   hf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
115   for (i = xs; i < xs + xm; i++) {
116     if (i < ctx->sf) {
117       xi = ctx->xmin + 0.5 * hs + i * hs;
118       /* Integrate over cell i using trapezoid rule with N points. */
119       for (k = 0; k < dof; k++) u[i * dof + k] = 0;
120       for (j = 0; j < N + 1; j++) {
121         xj = xi + hs * (j - N / 2) / (PetscReal)N;
122         PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
123         for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
124       }
125     } else if (i < ctx->fs) {
126       xi = ctx->xmin + ctx->sf * hs + 0.5 * hf + (i - ctx->sf) * hf;
127       /* Integrate over cell i using trapezoid rule with N points. */
128       for (k = 0; k < dof; k++) u[i * dof + k] = 0;
129       for (j = 0; j < N + 1; j++) {
130         xj = xi + hf * (j - N / 2) / (PetscReal)N;
131         PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
132         for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
133       }
134     } else {
135       xi = ctx->xmin + ctx->sf * hs + (ctx->fs - ctx->sf) * hf + 0.5 * hs + (i - ctx->fs) * hs;
136       /* Integrate over cell i using trapezoid rule with N points. */
137       for (k = 0; k < dof; k++) u[i * dof + k] = 0;
138       for (j = 0; j < N + 1; j++) {
139         xj = xi + hs * (j - N / 2) / (PetscReal)N;
140         PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
141         for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
142       }
143     }
144   }
145   PetscCall(DMDAVecRestoreArray(da, U, &u));
146   PetscCall(PetscFree(uj));
147   PetscFunctionReturn(0);
148 }
149 
150 static PetscErrorCode SolutionErrorNorms_2WaySplit(FVCtx *ctx, DM da, PetscReal t, Vec X, PetscReal *nrm1) {
151   Vec                Y;
152   PetscInt           i, Mx;
153   const PetscScalar *ptr_X, *ptr_Y;
154   PetscReal          hs, hf;
155 
156   PetscFunctionBeginUser;
157   PetscCall(VecGetSize(X, &Mx));
158   PetscCall(VecDuplicate(X, &Y));
159   PetscCall(FVSample_2WaySplit(ctx, da, t, Y));
160   hs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
161   hf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
162   PetscCall(VecGetArrayRead(X, &ptr_X));
163   PetscCall(VecGetArrayRead(Y, &ptr_Y));
164   for (i = 0; i < Mx; i++) {
165     if (i < ctx->sf || i > ctx->fs - 1) *nrm1 += hs * PetscAbs(ptr_X[i] - ptr_Y[i]);
166     else *nrm1 += hf * PetscAbs(ptr_X[i] - ptr_Y[i]);
167   }
168   PetscCall(VecRestoreArrayRead(X, &ptr_X));
169   PetscCall(VecRestoreArrayRead(Y, &ptr_Y));
170   PetscCall(VecDestroy(&Y));
171   PetscFunctionReturn(0);
172 }
173 
174 PetscErrorCode FVRHSFunction_2WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) {
175   FVCtx       *ctx = (FVCtx *)vctx;
176   PetscInt     i, j, k, Mx, dof, xs, xm, sf = ctx->sf, fs = ctx->fs;
177   PetscReal    hxf, hxs, cfl_idt = 0;
178   PetscScalar *x, *f, *slope;
179   Vec          Xloc;
180   DM           da;
181 
182   PetscFunctionBeginUser;
183   PetscCall(TSGetDM(ts, &da));
184   PetscCall(DMGetLocalVector(da, &Xloc));                                 /* Xloc contains ghost points                                     */
185   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); /* Mx is the number of center points                              */
186   hxs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
187   hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
188   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); /* X is solution vector which does not contain ghost points       */
189   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
190 
191   PetscCall(VecZeroEntries(F)); /* F is the right hand side function corresponds to center points */
192 
193   PetscCall(DMDAVecGetArray(da, Xloc, &x));
194   PetscCall(DMDAVecGetArray(da, F, &f));
195   PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); /* contains ghost points                                           */
196 
197   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
198 
199   if (ctx->bctype == FVBC_OUTFLOW) {
200     for (i = xs - 2; i < 0; i++) {
201       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
202     }
203     for (i = Mx; i < xs + xm + 2; i++) {
204       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
205     }
206   }
207   for (i = xs - 1; i < xs + xm + 1; i++) {
208     struct _LimitInfo info;
209     PetscScalar      *cjmpL, *cjmpR;
210     /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
211     PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
212     /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
213     PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
214     cjmpL = &ctx->cjmpLR[0];
215     cjmpR = &ctx->cjmpLR[dof];
216     for (j = 0; j < dof; j++) {
217       PetscScalar jmpL, jmpR;
218       jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
219       jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
220       for (k = 0; k < dof; k++) {
221         cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
222         cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
223       }
224     }
225     /* Apply limiter to the left and right characteristic jumps */
226     info.m   = dof;
227     info.hxs = hxs;
228     info.hxf = hxf;
229     (*ctx->limit2)(&info, cjmpL, cjmpR, ctx->sf, ctx->fs, i, ctx->cslope);
230     for (j = 0; j < dof; j++) {
231       PetscScalar tmp = 0;
232       for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
233       slope[i * dof + j] = tmp;
234     }
235   }
236 
237   for (i = xs; i < xs + xm + 1; i++) {
238     PetscReal    maxspeed;
239     PetscScalar *uL, *uR;
240     uL = &ctx->uLR[0];
241     uR = &ctx->uLR[dof];
242     if (i < sf) { /* slow region */
243       for (j = 0; j < dof; j++) {
244         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
245         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
246       }
247       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
248       if (i > xs) {
249         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxs;
250       }
251       if (i < xs + xm) {
252         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxs;
253       }
254     } else if (i == sf) { /* interface between the slow region and the fast region */
255       for (j = 0; j < dof; j++) {
256         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
257         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
258       }
259       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
260       if (i > xs) {
261         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxs;
262       }
263       if (i < xs + xm) {
264         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxf;
265       }
266     } else if (i > sf && i < fs) { /* fast region */
267       for (j = 0; j < dof; j++) {
268         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
269         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
270       }
271       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
272       if (i > xs) {
273         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxf;
274       }
275       if (i < xs + xm) {
276         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxf;
277       }
278     } else if (i == fs) { /* interface between the fast region and the slow region */
279       for (j = 0; j < dof; j++) {
280         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
281         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
282       }
283       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
284       if (i > xs) {
285         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxf;
286       }
287       if (i < xs + xm) {
288         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxs;
289       }
290     } else { /* slow region */
291       for (j = 0; j < dof; j++) {
292         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
293         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
294       }
295       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
296       cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hxs)); /* Max allowable value of 1/Delta t */
297       if (i > xs) {
298         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxs;
299       }
300       if (i < xs + xm) {
301         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxs;
302       }
303     }
304   }
305   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
306   PetscCall(DMDAVecRestoreArray(da, F, &f));
307   PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
308   PetscCall(DMRestoreLocalVector(da, &Xloc));
309   PetscCallMPI(MPI_Allreduce(&cfl_idt, &ctx->cfl_idt, 1, MPIU_SCALAR, MPIU_MAX, PetscObjectComm((PetscObject)da)));
310   if (0) {
311     /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */
312     PetscReal dt, tnow;
313     PetscCall(TSGetTimeStep(ts, &dt));
314     PetscCall(TSGetTime(ts, &tnow));
315     if (dt > 0.5 / ctx->cfl_idt) {
316       if (1) {
317         PetscCall(PetscPrintf(ctx->comm, "Stability constraint exceeded at t=%g, dt %g > %g\n", (double)tnow, (double)dt, (double)(0.5 / ctx->cfl_idt)));
318       } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Stability constraint exceeded, %g > %g", (double)dt, (double)(ctx->cfl / ctx->cfl_idt));
319     }
320   }
321   PetscFunctionReturn(0);
322 }
323 
324 /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */
325 PetscErrorCode FVRHSFunctionslow_2WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) {
326   FVCtx       *ctx = (FVCtx *)vctx;
327   PetscInt     i, j, k, Mx, dof, xs, xm, islow = 0, sf = ctx->sf, fs = ctx->fs, lsbwidth = ctx->lsbwidth, rsbwidth = ctx->rsbwidth;
328   PetscReal    hxs, hxf, cfl_idt = 0;
329   PetscScalar *x, *f, *slope;
330   Vec          Xloc;
331   DM           da;
332 
333   PetscFunctionBeginUser;
334   PetscCall(TSGetDM(ts, &da));
335   PetscCall(DMGetLocalVector(da, &Xloc));
336   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
337   hxs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
338   hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
339   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
340   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
341   PetscCall(VecZeroEntries(F));
342   PetscCall(DMDAVecGetArray(da, Xloc, &x));
343   PetscCall(VecGetArray(F, &f));
344   PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
345   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
346 
347   if (ctx->bctype == FVBC_OUTFLOW) {
348     for (i = xs - 2; i < 0; i++) {
349       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
350     }
351     for (i = Mx; i < xs + xm + 2; i++) {
352       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
353     }
354   }
355   for (i = xs - 1; i < xs + xm + 1; i++) {
356     struct _LimitInfo info;
357     PetscScalar      *cjmpL, *cjmpR;
358     if (i < sf - lsbwidth + 1 || i > fs + rsbwidth - 2) { /* slow components and the first and last fast components */
359       /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
360       PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
361       /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
362       PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
363       cjmpL = &ctx->cjmpLR[0];
364       cjmpR = &ctx->cjmpLR[dof];
365       for (j = 0; j < dof; j++) {
366         PetscScalar jmpL, jmpR;
367         jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
368         jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
369         for (k = 0; k < dof; k++) {
370           cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
371           cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
372         }
373       }
374       /* Apply limiter to the left and right characteristic jumps */
375       info.m   = dof;
376       info.hxs = hxs;
377       info.hxf = hxf;
378       (*ctx->limit2)(&info, cjmpL, cjmpR, ctx->sf, ctx->fs, i, ctx->cslope);
379       for (j = 0; j < dof; j++) {
380         PetscScalar tmp = 0;
381         for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
382         slope[i * dof + j] = tmp;
383       }
384     }
385   }
386 
387   for (i = xs; i < xs + xm + 1; i++) {
388     PetscReal    maxspeed;
389     PetscScalar *uL, *uR;
390     uL = &ctx->uLR[0];
391     uR = &ctx->uLR[dof];
392     if (i < sf - lsbwidth) { /* slow region */
393       for (j = 0; j < dof; j++) {
394         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
395         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
396       }
397       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
398       cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hxs)); /* Max allowable value of 1/Delta t */
399       if (i > xs) {
400         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
401       }
402       if (i < xs + xm) {
403         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
404         islow++;
405       }
406     }
407     if (i == sf - lsbwidth) { /* interface between the slow region and the fast region */
408       for (j = 0; j < dof; j++) {
409         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
410         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
411       }
412       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
413       if (i > xs) {
414         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
415       }
416     }
417     if (i == fs + rsbwidth) { /* slow region */
418       for (j = 0; j < dof; j++) {
419         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
420         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
421       }
422       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
423       if (i < xs + xm) {
424         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
425         islow++;
426       }
427     }
428     if (i > fs + rsbwidth) { /* slow region */
429       for (j = 0; j < dof; j++) {
430         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
431         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
432       }
433       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
434       if (i > xs) {
435         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
436       }
437       if (i < xs + xm) {
438         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
439         islow++;
440       }
441     }
442   }
443   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
444   PetscCall(VecRestoreArray(F, &f));
445   PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
446   PetscCall(DMRestoreLocalVector(da, &Xloc));
447   PetscCallMPI(MPI_Allreduce(&cfl_idt, &ctx->cfl_idt, 1, MPIU_SCALAR, MPIU_MAX, PetscObjectComm((PetscObject)da)));
448   PetscFunctionReturn(0);
449 }
450 
451 PetscErrorCode FVRHSFunctionslowbuffer_2WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) {
452   FVCtx       *ctx = (FVCtx *)vctx;
453   PetscInt     i, j, k, Mx, dof, xs, xm, islow = 0, sf = ctx->sf, fs = ctx->fs, lsbwidth = ctx->lsbwidth, rsbwidth = ctx->rsbwidth;
454   PetscReal    hxs, hxf;
455   PetscScalar *x, *f, *slope;
456   Vec          Xloc;
457   DM           da;
458 
459   PetscFunctionBeginUser;
460   PetscCall(TSGetDM(ts, &da));
461   PetscCall(DMGetLocalVector(da, &Xloc));
462   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
463   hxs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
464   hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
465   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
466   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
467   PetscCall(VecZeroEntries(F));
468   PetscCall(DMDAVecGetArray(da, Xloc, &x));
469   PetscCall(VecGetArray(F, &f));
470   PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
471   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
472 
473   if (ctx->bctype == FVBC_OUTFLOW) {
474     for (i = xs - 2; i < 0; i++) {
475       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
476     }
477     for (i = Mx; i < xs + xm + 2; i++) {
478       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
479     }
480   }
481   for (i = xs - 1; i < xs + xm + 1; i++) {
482     struct _LimitInfo info;
483     PetscScalar      *cjmpL, *cjmpR;
484     if ((i > sf - lsbwidth - 2 && i < sf + 1) || (i > fs - 2 && i < fs + rsbwidth + 1)) {
485       /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
486       PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
487       /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
488       PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
489       cjmpL = &ctx->cjmpLR[0];
490       cjmpR = &ctx->cjmpLR[dof];
491       for (j = 0; j < dof; j++) {
492         PetscScalar jmpL, jmpR;
493         jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
494         jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
495         for (k = 0; k < dof; k++) {
496           cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
497           cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
498         }
499       }
500       /* Apply limiter to the left and right characteristic jumps */
501       info.m   = dof;
502       info.hxs = hxs;
503       info.hxf = hxf;
504       (*ctx->limit2)(&info, cjmpL, cjmpR, ctx->sf, ctx->fs, i, ctx->cslope);
505       for (j = 0; j < dof; j++) {
506         PetscScalar tmp = 0;
507         for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
508         slope[i * dof + j] = tmp;
509       }
510     }
511   }
512 
513   for (i = xs; i < xs + xm + 1; i++) {
514     PetscReal    maxspeed;
515     PetscScalar *uL, *uR;
516     uL = &ctx->uLR[0];
517     uR = &ctx->uLR[dof];
518     if (i == sf - lsbwidth) {
519       for (j = 0; j < dof; j++) {
520         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
521         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
522       }
523       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
524       if (i < xs + xm) {
525         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
526         islow++;
527       }
528     }
529     if (i > sf - lsbwidth && i < sf) {
530       for (j = 0; j < dof; j++) {
531         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
532         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
533       }
534       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
535       if (i > xs) {
536         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
537       }
538       if (i < xs + xm) {
539         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
540         islow++;
541       }
542     }
543     if (i == sf) { /* interface between the slow region and the fast region */
544       for (j = 0; j < dof; j++) {
545         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
546         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
547       }
548       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
549       if (i > xs) {
550         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
551       }
552     }
553     if (i == fs) { /* interface between the fast region and the slow region */
554       for (j = 0; j < dof; j++) {
555         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
556         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
557       }
558       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
559       if (i < xs + xm) {
560         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
561         islow++;
562       }
563     }
564     if (i > fs && i < fs + rsbwidth) {
565       for (j = 0; j < dof; j++) {
566         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
567         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
568       }
569       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
570       if (i > xs) {
571         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
572       }
573       if (i < xs + xm) {
574         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
575         islow++;
576       }
577     }
578     if (i == fs + rsbwidth) {
579       for (j = 0; j < dof; j++) {
580         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
581         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
582       }
583       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
584       if (i > xs) {
585         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
586       }
587     }
588   }
589   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
590   PetscCall(VecRestoreArray(F, &f));
591   PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
592   PetscCall(DMRestoreLocalVector(da, &Xloc));
593   PetscFunctionReturn(0);
594 }
595 
596 /* --------------------------------- Finite Volume Solver for fast  parts ----------------------------------- */
597 PetscErrorCode FVRHSFunctionfast_2WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) {
598   FVCtx       *ctx = (FVCtx *)vctx;
599   PetscInt     i, j, k, Mx, dof, xs, xm, ifast = 0, sf = ctx->sf, fs = ctx->fs;
600   PetscReal    hxs, hxf;
601   PetscScalar *x, *f, *slope;
602   Vec          Xloc;
603   DM           da;
604 
605   PetscFunctionBeginUser;
606   PetscCall(TSGetDM(ts, &da));
607   PetscCall(DMGetLocalVector(da, &Xloc));
608   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
609   hxs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
610   hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
611   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
612   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
613   PetscCall(VecZeroEntries(F));
614   PetscCall(DMDAVecGetArray(da, Xloc, &x));
615   PetscCall(VecGetArray(F, &f));
616   PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
617   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
618 
619   if (ctx->bctype == FVBC_OUTFLOW) {
620     for (i = xs - 2; i < 0; i++) {
621       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
622     }
623     for (i = Mx; i < xs + xm + 2; i++) {
624       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
625     }
626   }
627   for (i = xs - 1; i < xs + xm + 1; i++) {
628     struct _LimitInfo info;
629     PetscScalar      *cjmpL, *cjmpR;
630     if (i > sf - 2 && i < fs + 1) {
631       PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
632       PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
633       cjmpL = &ctx->cjmpLR[0];
634       cjmpR = &ctx->cjmpLR[dof];
635       for (j = 0; j < dof; j++) {
636         PetscScalar jmpL, jmpR;
637         jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
638         jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
639         for (k = 0; k < dof; k++) {
640           cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
641           cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
642         }
643       }
644       /* Apply limiter to the left and right characteristic jumps */
645       info.m   = dof;
646       info.hxs = hxs;
647       info.hxf = hxf;
648       (*ctx->limit2)(&info, cjmpL, cjmpR, ctx->sf, ctx->fs, i, ctx->cslope);
649       for (j = 0; j < dof; j++) {
650         PetscScalar tmp = 0;
651         for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
652         slope[i * dof + j] = tmp;
653       }
654     }
655   }
656 
657   for (i = xs; i < xs + xm + 1; i++) {
658     PetscReal    maxspeed;
659     PetscScalar *uL, *uR;
660     uL = &ctx->uLR[0];
661     uR = &ctx->uLR[dof];
662     if (i == sf) { /* interface between the slow region and the fast region */
663       for (j = 0; j < dof; j++) {
664         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
665         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
666       }
667       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
668       if (i < xs + xm) {
669         for (j = 0; j < dof; j++) f[ifast * dof + j] += ctx->flux[j] / hxf;
670         ifast++;
671       }
672     }
673     if (i > sf && i < fs) { /* fast region */
674       for (j = 0; j < dof; j++) {
675         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
676         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
677       }
678       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
679       if (i > xs) {
680         for (j = 0; j < dof; j++) f[(ifast - 1) * dof + j] -= ctx->flux[j] / hxf;
681       }
682       if (i < xs + xm) {
683         for (j = 0; j < dof; j++) f[ifast * dof + j] += ctx->flux[j] / hxf;
684         ifast++;
685       }
686     }
687     if (i == fs) { /* interface between the fast region and the slow region */
688       for (j = 0; j < dof; j++) {
689         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
690         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
691       }
692       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
693       if (i > xs) {
694         for (j = 0; j < dof; j++) f[(ifast - 1) * dof + j] -= ctx->flux[j] / hxf;
695       }
696     }
697   }
698   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
699   PetscCall(VecRestoreArray(F, &f));
700   PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
701   PetscCall(DMRestoreLocalVector(da, &Xloc));
702   PetscFunctionReturn(0);
703 }
704 
705 int main(int argc, char *argv[]) {
706   char              lname[256] = "mc", physname[256] = "advect", final_fname[256] = "solution.m";
707   PetscFunctionList limiters = 0, physics = 0;
708   MPI_Comm          comm;
709   TS                ts;
710   DM                da;
711   Vec               X, X0, R;
712   FVCtx             ctx;
713   PetscInt          i, k, dof, xs, xm, Mx, draw = 0, count_slow, count_fast, islow = 0, ifast = 0, islowbuffer = 0, *index_slow, *index_fast, *index_slowbuffer;
714   PetscBool         view_final = PETSC_FALSE;
715   PetscReal         ptime;
716 
717   PetscFunctionBeginUser;
718   PetscCall(PetscInitialize(&argc, &argv, 0, help));
719   comm = PETSC_COMM_WORLD;
720   PetscCall(PetscMemzero(&ctx, sizeof(ctx)));
721 
722   /* Register limiters to be available on the command line */
723   PetscCall(PetscFunctionListAdd(&limiters, "upwind", Limit2_Upwind));
724   PetscCall(PetscFunctionListAdd(&limiters, "lax-wendroff", Limit2_LaxWendroff));
725   PetscCall(PetscFunctionListAdd(&limiters, "beam-warming", Limit2_BeamWarming));
726   PetscCall(PetscFunctionListAdd(&limiters, "fromm", Limit2_Fromm));
727   PetscCall(PetscFunctionListAdd(&limiters, "minmod", Limit2_Minmod));
728   PetscCall(PetscFunctionListAdd(&limiters, "superbee", Limit2_Superbee));
729   PetscCall(PetscFunctionListAdd(&limiters, "mc", Limit2_MC));
730   PetscCall(PetscFunctionListAdd(&limiters, "koren3", Limit2_Koren3));
731 
732   /* Register physical models to be available on the command line */
733   PetscCall(PetscFunctionListAdd(&physics, "advect", PhysicsCreate_Advect));
734 
735   ctx.comm   = comm;
736   ctx.cfl    = 0.9;
737   ctx.bctype = FVBC_PERIODIC;
738   ctx.xmin   = -1.0;
739   ctx.xmax   = 1.0;
740   PetscOptionsBegin(comm, NULL, "Finite Volume solver options", "");
741   PetscCall(PetscOptionsReal("-xmin", "X min", "", ctx.xmin, &ctx.xmin, NULL));
742   PetscCall(PetscOptionsReal("-xmax", "X max", "", ctx.xmax, &ctx.xmax, NULL));
743   PetscCall(PetscOptionsFList("-limit", "Name of flux imiter to use", "", limiters, lname, lname, sizeof(lname), NULL));
744   PetscCall(PetscOptionsInt("-draw", "Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)", "", draw, &draw, NULL));
745   PetscCall(PetscOptionsString("-view_final", "Write final solution in ASCII MATLAB format to given file name", "", final_fname, final_fname, sizeof(final_fname), &view_final));
746   PetscCall(PetscOptionsInt("-initial", "Initial condition (depends on the physics)", "", ctx.initial, &ctx.initial, NULL));
747   PetscCall(PetscOptionsBool("-exact", "Compare errors with exact solution", "", ctx.exact, &ctx.exact, NULL));
748   PetscCall(PetscOptionsBool("-simulation", "Compare errors with reference solution", "", ctx.simulation, &ctx.simulation, NULL));
749   PetscCall(PetscOptionsReal("-cfl", "CFL number to time step at", "", ctx.cfl, &ctx.cfl, NULL));
750   PetscCall(PetscOptionsEnum("-bc_type", "Boundary condition", "", FVBCTypes, (PetscEnum)ctx.bctype, (PetscEnum *)&ctx.bctype, NULL));
751   PetscCall(PetscOptionsInt("-hratio", "Spacing ratio", "", ctx.hratio, &ctx.hratio, NULL));
752   PetscOptionsEnd();
753 
754   /* Choose the limiter from the list of registered limiters */
755   PetscCall(PetscFunctionListFind(limiters, lname, &ctx.limit2));
756   PetscCheck(ctx.limit2, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Limiter '%s' not found", lname);
757 
758   /* Choose the physics from the list of registered models */
759   {
760     PetscErrorCode (*r)(FVCtx *);
761     PetscCall(PetscFunctionListFind(physics, physname, &r));
762     PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Physics '%s' not found", physname);
763     /* Create the physics, will set the number of fields and their names */
764     PetscCall((*r)(&ctx));
765   }
766 
767   /* Create a DMDA to manage the parallel grid */
768   PetscCall(DMDACreate1d(comm, DM_BOUNDARY_PERIODIC, 50, ctx.physics2.dof, 2, NULL, &da));
769   PetscCall(DMSetFromOptions(da));
770   PetscCall(DMSetUp(da));
771   /* Inform the DMDA of the field names provided by the physics. */
772   /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */
773   for (i = 0; i < ctx.physics2.dof; i++) PetscCall(DMDASetFieldName(da, i, ctx.physics2.fieldname[i]));
774   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
775   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
776 
777   /* Set coordinates of cell centers */
778   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));
779 
780   /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */
781   PetscCall(PetscMalloc4(dof * dof, &ctx.R, dof * dof, &ctx.Rinv, 2 * dof, &ctx.cjmpLR, 1 * dof, &ctx.cslope));
782   PetscCall(PetscMalloc3(2 * dof, &ctx.uLR, dof, &ctx.flux, dof, &ctx.speeds));
783 
784   /* Create a vector to store the solution and to save the initial state */
785   PetscCall(DMCreateGlobalVector(da, &X));
786   PetscCall(VecDuplicate(X, &X0));
787   PetscCall(VecDuplicate(X, &R));
788 
789   /* create index for slow parts and fast parts,
790      count_slow + count_fast = Mx, counts_slow*hs = 0.5, counts_fast*hf = 0.5 */
791   count_slow = Mx / (1.0 + ctx.hratio / 3.0);
792   PetscCheck(count_slow % 2 == 0, PETSC_COMM_WORLD, PETSC_ERR_USER, "Please adjust grid size Mx (-da_grid_x) and hratio (-hratio) so that Mx/(1+hartio/3) is even");
793   count_fast = Mx - count_slow;
794   ctx.sf     = count_slow / 2;
795   ctx.fs     = ctx.sf + count_fast;
796   PetscCall(PetscMalloc1(xm * dof, &index_slow));
797   PetscCall(PetscMalloc1(xm * dof, &index_fast));
798   PetscCall(PetscMalloc1(6 * dof, &index_slowbuffer));
799   if (((AdvectCtx *)ctx.physics2.user)->a > 0) {
800     ctx.lsbwidth = 2;
801     ctx.rsbwidth = 4;
802   } else {
803     ctx.lsbwidth = 4;
804     ctx.rsbwidth = 2;
805   }
806   for (i = xs; i < xs + xm; i++) {
807     if (i < ctx.sf - ctx.lsbwidth || i > ctx.fs + ctx.rsbwidth - 1)
808       for (k = 0; k < dof; k++) index_slow[islow++] = i * dof + k;
809     else if ((i >= ctx.sf - ctx.lsbwidth && i < ctx.sf) || (i > ctx.fs - 1 && i <= ctx.fs + ctx.rsbwidth - 1))
810       for (k = 0; k < dof; k++) index_slowbuffer[islowbuffer++] = i * dof + k;
811     else
812       for (k = 0; k < dof; k++) index_fast[ifast++] = i * dof + k;
813   }
814   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islow, index_slow, PETSC_COPY_VALUES, &ctx.iss));
815   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, ifast, index_fast, PETSC_COPY_VALUES, &ctx.isf));
816   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islowbuffer, index_slowbuffer, PETSC_COPY_VALUES, &ctx.issb));
817 
818   /* Create a time-stepping object */
819   PetscCall(TSCreate(comm, &ts));
820   PetscCall(TSSetDM(ts, da));
821   PetscCall(TSSetRHSFunction(ts, R, FVRHSFunction_2WaySplit, &ctx));
822   PetscCall(TSRHSSplitSetIS(ts, "slow", ctx.iss));
823   PetscCall(TSRHSSplitSetIS(ts, "slowbuffer", ctx.issb));
824   PetscCall(TSRHSSplitSetIS(ts, "fast", ctx.isf));
825   PetscCall(TSRHSSplitSetRHSFunction(ts, "slow", NULL, FVRHSFunctionslow_2WaySplit, &ctx));
826   PetscCall(TSRHSSplitSetRHSFunction(ts, "fast", NULL, FVRHSFunctionfast_2WaySplit, &ctx));
827   PetscCall(TSRHSSplitSetRHSFunction(ts, "slowbuffer", NULL, FVRHSFunctionslowbuffer_2WaySplit, &ctx));
828 
829   PetscCall(TSSetType(ts, TSSSP));
830   /*PetscCall(TSSetType(ts,TSMPRK));*/
831   PetscCall(TSSetMaxTime(ts, 10));
832   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
833 
834   /* Compute initial conditions and starting time step */
835   PetscCall(FVSample_2WaySplit(&ctx, da, 0, X0));
836   PetscCall(FVRHSFunction_2WaySplit(ts, 0, X0, X, (void *)&ctx)); /* Initial function evaluation, only used to determine max speed */
837   PetscCall(VecCopy(X0, X));                                      /* The function value was not used so we set X=X0 again */
838   PetscCall(TSSetTimeStep(ts, ctx.cfl / ctx.cfl_idt));
839   PetscCall(TSSetFromOptions(ts)); /* Take runtime options */
840   PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD));
841   {
842     PetscInt           steps;
843     PetscScalar        mass_initial, mass_final, mass_difference, mass_differenceg;
844     const PetscScalar *ptr_X, *ptr_X0;
845     const PetscReal    hs = (ctx.xmax - ctx.xmin) * 3.0 / 4.0 / count_slow;
846     const PetscReal    hf = (ctx.xmax - ctx.xmin) / 4.0 / count_fast;
847 
848     PetscCall(TSSolve(ts, X));
849     PetscCall(TSGetSolveTime(ts, &ptime));
850     PetscCall(TSGetStepNumber(ts, &steps));
851     /* calculate the total mass at initial time and final time */
852     mass_initial = 0.0;
853     mass_final   = 0.0;
854     PetscCall(DMDAVecGetArrayRead(da, X0, (void *)&ptr_X0));
855     PetscCall(DMDAVecGetArrayRead(da, X, (void *)&ptr_X));
856     for (i = xs; i < xs + xm; i++) {
857       if (i < ctx.sf || i > ctx.fs - 1) {
858         for (k = 0; k < dof; k++) {
859           mass_initial = mass_initial + hs * ptr_X0[i * dof + k];
860           mass_final   = mass_final + hs * ptr_X[i * dof + k];
861         }
862       } else {
863         for (k = 0; k < dof; k++) {
864           mass_initial = mass_initial + hf * ptr_X0[i * dof + k];
865           mass_final   = mass_final + hf * ptr_X[i * dof + k];
866         }
867       }
868     }
869     PetscCall(DMDAVecRestoreArrayRead(da, X0, (void *)&ptr_X0));
870     PetscCall(DMDAVecRestoreArrayRead(da, X, (void *)&ptr_X));
871     mass_difference = mass_final - mass_initial;
872     PetscCallMPI(MPI_Allreduce(&mass_difference, &mass_differenceg, 1, MPIU_SCALAR, MPIU_SUM, comm));
873     PetscCall(PetscPrintf(comm, "Mass difference %g\n", (double)mass_differenceg));
874     PetscCall(PetscPrintf(comm, "Final time %g, steps %" PetscInt_FMT "\n", (double)ptime, steps));
875     PetscCall(PetscPrintf(comm, "Maximum allowable stepsize according to CFL %g\n", (double)(1.0 / ctx.cfl_idt)));
876     if (ctx.exact) {
877       PetscReal nrm1 = 0;
878       PetscCall(SolutionErrorNorms_2WaySplit(&ctx, da, ptime, X, &nrm1));
879       PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %g\n", (double)nrm1));
880     }
881     if (ctx.simulation) {
882       PetscReal          nrm1 = 0;
883       PetscViewer        fd;
884       char               filename[PETSC_MAX_PATH_LEN] = "binaryoutput";
885       Vec                XR;
886       PetscBool          flg;
887       const PetscScalar *ptr_XR;
888       PetscCall(PetscOptionsGetString(NULL, NULL, "-f", filename, sizeof(filename), &flg));
889       PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate binary file with the -f option");
890       PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, filename, FILE_MODE_READ, &fd));
891       PetscCall(VecDuplicate(X0, &XR));
892       PetscCall(VecLoad(XR, fd));
893       PetscCall(PetscViewerDestroy(&fd));
894       PetscCall(VecGetArrayRead(X, &ptr_X));
895       PetscCall(VecGetArrayRead(XR, &ptr_XR));
896       for (i = xs; i < xs + xm; i++) {
897         if (i < ctx.sf || i > ctx.fs - 1)
898           for (k = 0; k < dof; k++) nrm1 = nrm1 + hs * PetscAbs(ptr_X[i * dof + k] - ptr_XR[i * dof + k]);
899         else
900           for (k = 0; k < dof; k++) nrm1 = nrm1 + hf * PetscAbs(ptr_X[i * dof + k] - ptr_XR[i * dof + k]);
901       }
902       PetscCall(VecRestoreArrayRead(X, &ptr_X));
903       PetscCall(VecRestoreArrayRead(XR, &ptr_XR));
904       PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %g\n", (double)nrm1));
905       PetscCall(VecDestroy(&XR));
906     }
907   }
908 
909   PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD));
910   if (draw & 0x1) PetscCall(VecView(X0, PETSC_VIEWER_DRAW_WORLD));
911   if (draw & 0x2) PetscCall(VecView(X, PETSC_VIEWER_DRAW_WORLD));
912   if (draw & 0x4) {
913     Vec Y;
914     PetscCall(VecDuplicate(X, &Y));
915     PetscCall(FVSample_2WaySplit(&ctx, da, ptime, Y));
916     PetscCall(VecAYPX(Y, -1, X));
917     PetscCall(VecView(Y, PETSC_VIEWER_DRAW_WORLD));
918     PetscCall(VecDestroy(&Y));
919   }
920 
921   if (view_final) {
922     PetscViewer viewer;
923     PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, final_fname, &viewer));
924     PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB));
925     PetscCall(VecView(X, viewer));
926     PetscCall(PetscViewerPopFormat(viewer));
927     PetscCall(PetscViewerDestroy(&viewer));
928   }
929 
930   /* Clean up */
931   PetscCall((*ctx.physics2.destroy)(ctx.physics2.user));
932   for (i = 0; i < ctx.physics2.dof; i++) PetscCall(PetscFree(ctx.physics2.fieldname[i]));
933   PetscCall(PetscFree4(ctx.R, ctx.Rinv, ctx.cjmpLR, ctx.cslope));
934   PetscCall(PetscFree3(ctx.uLR, ctx.flux, ctx.speeds));
935   PetscCall(VecDestroy(&X));
936   PetscCall(VecDestroy(&X0));
937   PetscCall(VecDestroy(&R));
938   PetscCall(DMDestroy(&da));
939   PetscCall(TSDestroy(&ts));
940   PetscCall(ISDestroy(&ctx.iss));
941   PetscCall(ISDestroy(&ctx.isf));
942   PetscCall(ISDestroy(&ctx.issb));
943   PetscCall(PetscFree(index_slow));
944   PetscCall(PetscFree(index_fast));
945   PetscCall(PetscFree(index_slowbuffer));
946   PetscCall(PetscFunctionListDestroy(&limiters));
947   PetscCall(PetscFunctionListDestroy(&physics));
948   PetscCall(PetscFinalize());
949   return 0;
950 }
951 
952 /*TEST
953 
954     build:
955       requires: !complex
956       depends: finitevolume1d.c
957 
958     test:
959       suffix: 1
960       args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -limit mc -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22
961 
962     test:
963       suffix: 2
964       args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -limit mc -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 0
965       output_file: output/ex6_1.out
966 
967 TEST*/
968