xref: /petsc/src/ts/tutorials/multirate/ex6.c (revision dd8e379b23d2ef935d8131fb74f7eb73fef09263)
1c4762a1bSJed Brown /*
2c4762a1bSJed Brown   Note:
36aad120cSJose E. Roman     -hratio is the ratio between mesh size of coarse grids and fine grids
4c4762a1bSJed Brown     -ts_rk_dtratio is the ratio between the large stepsize and the small stepsize
5c4762a1bSJed Brown */
6c4762a1bSJed Brown 
7c4762a1bSJed Brown static const char help[] = "1D periodic Finite Volume solver in slope-limiter form with semidiscrete time stepping.\n"
8c4762a1bSJed Brown                            "  advection   - Constant coefficient scalar advection\n"
9c4762a1bSJed Brown                            "                u_t       + (a*u)_x               = 0\n"
10c4762a1bSJed Brown                            "  for this toy problem, we choose different meshsizes for different sub-domains (slow-fast-slow), say\n"
11c4762a1bSJed Brown                            "                hxs  = hratio*hxf \n"
12c4762a1bSJed Brown                            "  where hxs and hxf are the grid spacings for coarse and fine grids respectively.\n"
13c4762a1bSJed Brown                            "  exact       - Exact Riemann solver which usually needs to perform a Newton iteration to connect\n"
14c4762a1bSJed Brown                            "                the states across shocks and rarefactions\n"
15c4762a1bSJed Brown                            "  simulation  - use reference solution which is generated by smaller time step size to be true solution,\n"
16c4762a1bSJed Brown                            "                also the reference solution should be generated by user and stored in a binary file.\n"
17c4762a1bSJed Brown                            "  characteristic - Limit the characteristic variables, this is usually preferred (default)\n"
18c4762a1bSJed Brown                            "Several initial conditions can be chosen with -initial N\n\n"
19c4762a1bSJed Brown                            "The problem size should be set with -da_grid_x M\n\n";
20c4762a1bSJed Brown 
21c4762a1bSJed Brown #include <petscts.h>
22c4762a1bSJed Brown #include <petscdm.h>
23c4762a1bSJed Brown #include <petscdmda.h>
24c4762a1bSJed Brown #include <petscdraw.h>
25c4762a1bSJed Brown #include "finitevolume1d.h"
26c4762a1bSJed Brown 
27d71ae5a4SJacob Faibussowitsch static inline PetscReal RangeMod(PetscReal a, PetscReal xmin, PetscReal xmax)
28d71ae5a4SJacob Faibussowitsch {
299371c9d4SSatish Balay   PetscReal range = xmax - xmin;
309371c9d4SSatish Balay   return xmin + PetscFmodReal(range + PetscFmodReal(a, range), range);
319371c9d4SSatish Balay }
32c4762a1bSJed Brown 
33c4762a1bSJed Brown /* --------------------------------- Advection ----------------------------------- */
34c4762a1bSJed Brown typedef struct {
35c4762a1bSJed Brown   PetscReal a; /* advective velocity */
36c4762a1bSJed Brown } AdvectCtx;
37c4762a1bSJed Brown 
38d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsRiemann_Advect(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed)
39d71ae5a4SJacob Faibussowitsch {
40c4762a1bSJed Brown   AdvectCtx *ctx = (AdvectCtx *)vctx;
41c4762a1bSJed Brown   PetscReal  speed;
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   PetscFunctionBeginUser;
44c4762a1bSJed Brown   speed     = ctx->a;
45c4762a1bSJed Brown   flux[0]   = PetscMax(0, speed) * uL[0] + PetscMin(0, speed) * uR[0];
46c4762a1bSJed Brown   *maxspeed = speed;
473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
48c4762a1bSJed Brown }
49c4762a1bSJed Brown 
50d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx, PetscInt m, const PetscScalar *u, PetscScalar *X, PetscScalar *Xi, PetscReal *speeds)
51d71ae5a4SJacob Faibussowitsch {
52c4762a1bSJed Brown   AdvectCtx *ctx = (AdvectCtx *)vctx;
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   PetscFunctionBeginUser;
55c4762a1bSJed Brown   X[0]      = 1.;
56c4762a1bSJed Brown   Xi[0]     = 1.;
57c4762a1bSJed Brown   speeds[0] = ctx->a;
583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
59c4762a1bSJed Brown }
60c4762a1bSJed Brown 
61d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsSample_Advect(void *vctx, PetscInt initial, FVBCType bctype, PetscReal xmin, PetscReal xmax, PetscReal t, PetscReal x, PetscReal *u)
62d71ae5a4SJacob Faibussowitsch {
63c4762a1bSJed Brown   AdvectCtx *ctx = (AdvectCtx *)vctx;
64c4762a1bSJed Brown   PetscReal  a   = ctx->a, x0;
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   PetscFunctionBeginUser;
67c4762a1bSJed Brown   switch (bctype) {
68d71ae5a4SJacob Faibussowitsch   case FVBC_OUTFLOW:
69d71ae5a4SJacob Faibussowitsch     x0 = x - a * t;
70d71ae5a4SJacob Faibussowitsch     break;
71d71ae5a4SJacob Faibussowitsch   case FVBC_PERIODIC:
72d71ae5a4SJacob Faibussowitsch     x0 = RangeMod(x - a * t, xmin, xmax);
73d71ae5a4SJacob Faibussowitsch     break;
74d71ae5a4SJacob Faibussowitsch   default:
75d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown BCType");
76c4762a1bSJed Brown   }
77c4762a1bSJed Brown   switch (initial) {
78d71ae5a4SJacob Faibussowitsch   case 0:
79d71ae5a4SJacob Faibussowitsch     u[0] = (x0 < 0) ? 1 : -1;
80d71ae5a4SJacob Faibussowitsch     break;
81d71ae5a4SJacob Faibussowitsch   case 1:
82d71ae5a4SJacob Faibussowitsch     u[0] = (x0 < 0) ? -1 : 1;
83d71ae5a4SJacob Faibussowitsch     break;
84d71ae5a4SJacob Faibussowitsch   case 2:
85d71ae5a4SJacob Faibussowitsch     u[0] = (0 < x0 && x0 < 1) ? 1 : 0;
86d71ae5a4SJacob Faibussowitsch     break;
87d71ae5a4SJacob Faibussowitsch   case 3:
88d71ae5a4SJacob Faibussowitsch     u[0] = PetscSinReal(2 * PETSC_PI * x0);
89d71ae5a4SJacob Faibussowitsch     break;
90d71ae5a4SJacob Faibussowitsch   case 4:
91d71ae5a4SJacob Faibussowitsch     u[0] = PetscAbs(x0);
92d71ae5a4SJacob Faibussowitsch     break;
93d71ae5a4SJacob Faibussowitsch   case 5:
94d71ae5a4SJacob Faibussowitsch     u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2 * PETSC_PI * x0));
95d71ae5a4SJacob Faibussowitsch     break;
96d71ae5a4SJacob Faibussowitsch   case 6:
97d71ae5a4SJacob Faibussowitsch     u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2 - x0 : 0));
98d71ae5a4SJacob Faibussowitsch     break;
99d71ae5a4SJacob Faibussowitsch   case 7:
100d71ae5a4SJacob Faibussowitsch     u[0] = PetscPowReal(PetscSinReal(PETSC_PI * x0), 10.0);
101d71ae5a4SJacob Faibussowitsch     break;
102d71ae5a4SJacob Faibussowitsch   default:
103d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition");
104c4762a1bSJed Brown   }
1053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
106c4762a1bSJed Brown }
107c4762a1bSJed Brown 
108d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx)
109d71ae5a4SJacob Faibussowitsch {
110c4762a1bSJed Brown   AdvectCtx *user;
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   PetscFunctionBeginUser;
1139566063dSJacob Faibussowitsch   PetscCall(PetscNew(&user));
114c4762a1bSJed Brown   ctx->physics2.sample2         = PhysicsSample_Advect;
115c4762a1bSJed Brown   ctx->physics2.riemann2        = PhysicsRiemann_Advect;
116c4762a1bSJed Brown   ctx->physics2.characteristic2 = PhysicsCharacteristic_Advect;
117c4762a1bSJed Brown   ctx->physics2.destroy         = PhysicsDestroy_SimpleFree;
118c4762a1bSJed Brown   ctx->physics2.user            = user;
119c4762a1bSJed Brown   ctx->physics2.dof             = 1;
1209566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("u", &ctx->physics2.fieldname[0]));
121c4762a1bSJed Brown   user->a = 1;
122d0609cedSBarry Smith   PetscOptionsBegin(ctx->comm, ctx->prefix, "Options for advection", "");
123d71ae5a4SJacob Faibussowitsch   {
124d71ae5a4SJacob Faibussowitsch     PetscCall(PetscOptionsReal("-physics_advect_a", "Speed", "", user->a, &user->a, NULL));
125d71ae5a4SJacob Faibussowitsch   }
126d0609cedSBarry Smith   PetscOptionsEnd();
1273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
128c4762a1bSJed Brown }
129c4762a1bSJed Brown 
130d71ae5a4SJacob Faibussowitsch PetscErrorCode FVSample_2WaySplit(FVCtx *ctx, DM da, PetscReal time, Vec U)
131d71ae5a4SJacob Faibussowitsch {
132c4762a1bSJed Brown   PetscScalar   *u, *uj, xj, xi;
133c4762a1bSJed Brown   PetscInt       i, j, k, dof, xs, xm, Mx;
134c4762a1bSJed Brown   const PetscInt N = 200;
135c4762a1bSJed Brown   PetscReal      hs, hf;
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   PetscFunctionBeginUser;
1383c633725SBarry Smith   PetscCheck(ctx->physics2.sample2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Physics has not provided a sampling function");
1399566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
1409566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
1419566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, U, &u));
1429566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dof, &uj));
143c4762a1bSJed Brown   hs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
144c4762a1bSJed Brown   hf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
145c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
146c4762a1bSJed Brown     if (i < ctx->sf) {
147c4762a1bSJed Brown       xi = ctx->xmin + 0.5 * hs + i * hs;
148c4762a1bSJed Brown       /* Integrate over cell i using trapezoid rule with N points. */
149c4762a1bSJed Brown       for (k = 0; k < dof; k++) u[i * dof + k] = 0;
150c4762a1bSJed Brown       for (j = 0; j < N + 1; j++) {
151c4762a1bSJed Brown         xj = xi + hs * (j - N / 2) / (PetscReal)N;
1529566063dSJacob Faibussowitsch         PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
153c4762a1bSJed Brown         for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
154c4762a1bSJed Brown       }
155c4762a1bSJed Brown     } else if (i < ctx->fs) {
156c4762a1bSJed Brown       xi = ctx->xmin + ctx->sf * hs + 0.5 * hf + (i - ctx->sf) * hf;
157c4762a1bSJed Brown       /* Integrate over cell i using trapezoid rule with N points. */
158c4762a1bSJed Brown       for (k = 0; k < dof; k++) u[i * dof + k] = 0;
159c4762a1bSJed Brown       for (j = 0; j < N + 1; j++) {
160c4762a1bSJed Brown         xj = xi + hf * (j - N / 2) / (PetscReal)N;
1619566063dSJacob Faibussowitsch         PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
162c4762a1bSJed Brown         for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
163c4762a1bSJed Brown       }
164c4762a1bSJed Brown     } else {
165c4762a1bSJed Brown       xi = ctx->xmin + ctx->sf * hs + (ctx->fs - ctx->sf) * hf + 0.5 * hs + (i - ctx->fs) * hs;
166c4762a1bSJed Brown       /* Integrate over cell i using trapezoid rule with N points. */
167c4762a1bSJed Brown       for (k = 0; k < dof; k++) u[i * dof + k] = 0;
168c4762a1bSJed Brown       for (j = 0; j < N + 1; j++) {
169c4762a1bSJed Brown         xj = xi + hs * (j - N / 2) / (PetscReal)N;
1709566063dSJacob Faibussowitsch         PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
171c4762a1bSJed Brown         for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
172c4762a1bSJed Brown       }
173c4762a1bSJed Brown     }
174c4762a1bSJed Brown   }
1759566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, U, &u));
1769566063dSJacob Faibussowitsch   PetscCall(PetscFree(uj));
1773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
178c4762a1bSJed Brown }
179c4762a1bSJed Brown 
180d71ae5a4SJacob Faibussowitsch static PetscErrorCode SolutionErrorNorms_2WaySplit(FVCtx *ctx, DM da, PetscReal t, Vec X, PetscReal *nrm1)
181d71ae5a4SJacob Faibussowitsch {
182c4762a1bSJed Brown   Vec                Y;
183c4762a1bSJed Brown   PetscInt           i, Mx;
184c4762a1bSJed Brown   const PetscScalar *ptr_X, *ptr_Y;
185c4762a1bSJed Brown   PetscReal          hs, hf;
186c4762a1bSJed Brown 
187c4762a1bSJed Brown   PetscFunctionBeginUser;
1889566063dSJacob Faibussowitsch   PetscCall(VecGetSize(X, &Mx));
1899566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X, &Y));
1909566063dSJacob Faibussowitsch   PetscCall(FVSample_2WaySplit(ctx, da, t, Y));
191c4762a1bSJed Brown   hs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
192c4762a1bSJed Brown   hf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
1939566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &ptr_X));
1949566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &ptr_Y));
195c4762a1bSJed Brown   for (i = 0; i < Mx; i++) {
196c4762a1bSJed Brown     if (i < ctx->sf || i > ctx->fs - 1) *nrm1 += hs * PetscAbs(ptr_X[i] - ptr_Y[i]);
197c4762a1bSJed Brown     else *nrm1 += hf * PetscAbs(ptr_X[i] - ptr_Y[i]);
198c4762a1bSJed Brown   }
1999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &ptr_X));
2009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &ptr_Y));
2019566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Y));
2023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
203c4762a1bSJed Brown }
204c4762a1bSJed Brown 
205d71ae5a4SJacob Faibussowitsch PetscErrorCode FVRHSFunction_2WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx)
206d71ae5a4SJacob Faibussowitsch {
207c4762a1bSJed Brown   FVCtx       *ctx = (FVCtx *)vctx;
208c4762a1bSJed Brown   PetscInt     i, j, k, Mx, dof, xs, xm, sf = ctx->sf, fs = ctx->fs;
209c4762a1bSJed Brown   PetscReal    hxf, hxs, cfl_idt = 0;
210c4762a1bSJed Brown   PetscScalar *x, *f, *slope;
211c4762a1bSJed Brown   Vec          Xloc;
212c4762a1bSJed Brown   DM           da;
213c4762a1bSJed Brown 
214c4762a1bSJed Brown   PetscFunctionBeginUser;
2159566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
2169566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &Xloc));                                 /* Xloc contains ghost points                                     */
2179566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); /* Mx is the number of center points                              */
218c4762a1bSJed Brown   hxs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
219c4762a1bSJed Brown   hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
2209566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); /* X is solution vector which does not contain ghost points       */
2219566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
222c4762a1bSJed Brown 
223*dd8e379bSPierre Jolivet   PetscCall(VecZeroEntries(F)); /* F is the right-hand side function corresponds to center points */
224c4762a1bSJed Brown 
2259566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Xloc, &x));
2269566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, F, &f));
2279566063dSJacob Faibussowitsch   PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); /* contains ghost points                                           */
228c4762a1bSJed Brown 
2299566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
230c4762a1bSJed Brown 
231c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
232c4762a1bSJed Brown     for (i = xs - 2; i < 0; i++) {
233c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
234c4762a1bSJed Brown     }
235c4762a1bSJed Brown     for (i = Mx; i < xs + xm + 2; i++) {
236c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
237c4762a1bSJed Brown     }
238c4762a1bSJed Brown   }
239c4762a1bSJed Brown   for (i = xs - 1; i < xs + xm + 1; i++) {
240c4762a1bSJed Brown     struct _LimitInfo info;
241c4762a1bSJed Brown     PetscScalar      *cjmpL, *cjmpR;
242c4762a1bSJed Brown     /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
2439566063dSJacob Faibussowitsch     PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
244c4762a1bSJed Brown     /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
2459566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
246c4762a1bSJed Brown     cjmpL = &ctx->cjmpLR[0];
247c4762a1bSJed Brown     cjmpR = &ctx->cjmpLR[dof];
248c4762a1bSJed Brown     for (j = 0; j < dof; j++) {
249c4762a1bSJed Brown       PetscScalar jmpL, jmpR;
250c4762a1bSJed Brown       jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
251c4762a1bSJed Brown       jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
252c4762a1bSJed Brown       for (k = 0; k < dof; k++) {
253c4762a1bSJed Brown         cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
254c4762a1bSJed Brown         cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
255c4762a1bSJed Brown       }
256c4762a1bSJed Brown     }
257c4762a1bSJed Brown     /* Apply limiter to the left and right characteristic jumps */
258c4762a1bSJed Brown     info.m   = dof;
259c4762a1bSJed Brown     info.hxs = hxs;
260c4762a1bSJed Brown     info.hxf = hxf;
261c4762a1bSJed Brown     (*ctx->limit2)(&info, cjmpL, cjmpR, ctx->sf, ctx->fs, i, ctx->cslope);
262c4762a1bSJed Brown     for (j = 0; j < dof; j++) {
263c4762a1bSJed Brown       PetscScalar tmp = 0;
264c4762a1bSJed Brown       for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
265c4762a1bSJed Brown       slope[i * dof + j] = tmp;
266c4762a1bSJed Brown     }
267c4762a1bSJed Brown   }
268c4762a1bSJed Brown 
269c4762a1bSJed Brown   for (i = xs; i < xs + xm + 1; i++) {
270c4762a1bSJed Brown     PetscReal    maxspeed;
271c4762a1bSJed Brown     PetscScalar *uL, *uR;
272c4762a1bSJed Brown     uL = &ctx->uLR[0];
273c4762a1bSJed Brown     uR = &ctx->uLR[dof];
274c4762a1bSJed Brown     if (i < sf) { /* slow region */
275c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
276c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
277c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
278c4762a1bSJed Brown       }
2799566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
280c4762a1bSJed Brown       if (i > xs) {
281c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxs;
282c4762a1bSJed Brown       }
283c4762a1bSJed Brown       if (i < xs + xm) {
284c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxs;
285c4762a1bSJed Brown       }
286c4762a1bSJed Brown     } else if (i == sf) { /* interface between the slow region and the fast region */
287c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
288c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
289c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
290c4762a1bSJed Brown       }
2919566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
292c4762a1bSJed Brown       if (i > xs) {
293c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxs;
294c4762a1bSJed Brown       }
295c4762a1bSJed Brown       if (i < xs + xm) {
296c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxf;
297c4762a1bSJed Brown       }
298c4762a1bSJed Brown     } else if (i > sf && i < fs) { /* fast region */
299c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
300c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
301c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
302c4762a1bSJed Brown       }
3039566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
304c4762a1bSJed Brown       if (i > xs) {
305c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxf;
306c4762a1bSJed Brown       }
307c4762a1bSJed Brown       if (i < xs + xm) {
308c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxf;
309c4762a1bSJed Brown       }
310c4762a1bSJed Brown     } else if (i == fs) { /* interface between the fast region and the slow region */
311c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
312c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
313c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
314c4762a1bSJed Brown       }
3159566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
316c4762a1bSJed Brown       if (i > xs) {
317c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxf;
318c4762a1bSJed Brown       }
319c4762a1bSJed Brown       if (i < xs + xm) {
320c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxs;
321c4762a1bSJed Brown       }
322c4762a1bSJed Brown     } else { /* slow region */
323c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
324c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
325c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
326c4762a1bSJed Brown       }
3279566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
328c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hxs)); /* Max allowable value of 1/Delta t */
329c4762a1bSJed Brown       if (i > xs) {
330c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxs;
331c4762a1bSJed Brown       }
332c4762a1bSJed Brown       if (i < xs + xm) {
333c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxs;
334c4762a1bSJed Brown       }
335c4762a1bSJed Brown     }
336c4762a1bSJed Brown   }
3379566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
3389566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, F, &f));
3399566063dSJacob Faibussowitsch   PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
3409566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &Xloc));
341712fec58SPierre Jolivet   PetscCall(MPIU_Allreduce(&cfl_idt, &ctx->cfl_idt, 1, MPIU_SCALAR, MPIU_MAX, PetscObjectComm((PetscObject)da)));
342c4762a1bSJed Brown   if (0) {
343c4762a1bSJed Brown     /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */
344c4762a1bSJed Brown     PetscReal dt, tnow;
3459566063dSJacob Faibussowitsch     PetscCall(TSGetTimeStep(ts, &dt));
3469566063dSJacob Faibussowitsch     PetscCall(TSGetTime(ts, &tnow));
347c4762a1bSJed Brown     if (dt > 0.5 / ctx->cfl_idt) {
348c4762a1bSJed Brown       if (1) {
3499566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(ctx->comm, "Stability constraint exceeded at t=%g, dt %g > %g\n", (double)tnow, (double)dt, (double)(0.5 / ctx->cfl_idt)));
35098921bdaSJacob Faibussowitsch       } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Stability constraint exceeded, %g > %g", (double)dt, (double)(ctx->cfl / ctx->cfl_idt));
351c4762a1bSJed Brown     }
352c4762a1bSJed Brown   }
3533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
354c4762a1bSJed Brown }
355c4762a1bSJed Brown 
356c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */
357d71ae5a4SJacob Faibussowitsch PetscErrorCode FVRHSFunctionslow_2WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx)
358d71ae5a4SJacob Faibussowitsch {
359c4762a1bSJed Brown   FVCtx       *ctx = (FVCtx *)vctx;
360c4762a1bSJed Brown   PetscInt     i, j, k, Mx, dof, xs, xm, islow = 0, sf = ctx->sf, fs = ctx->fs, lsbwidth = ctx->lsbwidth, rsbwidth = ctx->rsbwidth;
361c4762a1bSJed Brown   PetscReal    hxs, hxf, cfl_idt = 0;
362c4762a1bSJed Brown   PetscScalar *x, *f, *slope;
363c4762a1bSJed Brown   Vec          Xloc;
364c4762a1bSJed Brown   DM           da;
365c4762a1bSJed Brown 
366c4762a1bSJed Brown   PetscFunctionBeginUser;
3679566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
3689566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &Xloc));
3699566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
370c4762a1bSJed Brown   hxs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
371c4762a1bSJed Brown   hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
3729566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
3739566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
3749566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
3759566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Xloc, &x));
3769566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
3779566063dSJacob Faibussowitsch   PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
3789566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
379c4762a1bSJed Brown 
380c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
381c4762a1bSJed Brown     for (i = xs - 2; i < 0; i++) {
382c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
383c4762a1bSJed Brown     }
384c4762a1bSJed Brown     for (i = Mx; i < xs + xm + 2; i++) {
385c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
386c4762a1bSJed Brown     }
387c4762a1bSJed Brown   }
388c4762a1bSJed Brown   for (i = xs - 1; i < xs + xm + 1; i++) {
389c4762a1bSJed Brown     struct _LimitInfo info;
390c4762a1bSJed Brown     PetscScalar      *cjmpL, *cjmpR;
391c4762a1bSJed Brown     if (i < sf - lsbwidth + 1 || i > fs + rsbwidth - 2) { /* slow components and the first and last fast components */
392c4762a1bSJed Brown       /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
3939566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
394c4762a1bSJed Brown       /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
3959566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
396c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
397c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
398c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
399c4762a1bSJed Brown         PetscScalar jmpL, jmpR;
400c4762a1bSJed Brown         jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
401c4762a1bSJed Brown         jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
402c4762a1bSJed Brown         for (k = 0; k < dof; k++) {
403c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
404c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
405c4762a1bSJed Brown         }
406c4762a1bSJed Brown       }
407c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
408c4762a1bSJed Brown       info.m   = dof;
409c4762a1bSJed Brown       info.hxs = hxs;
410c4762a1bSJed Brown       info.hxf = hxf;
411c4762a1bSJed Brown       (*ctx->limit2)(&info, cjmpL, cjmpR, ctx->sf, ctx->fs, i, ctx->cslope);
412c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
413c4762a1bSJed Brown         PetscScalar tmp = 0;
414c4762a1bSJed Brown         for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
415c4762a1bSJed Brown         slope[i * dof + j] = tmp;
416c4762a1bSJed Brown       }
417c4762a1bSJed Brown     }
418c4762a1bSJed Brown   }
419c4762a1bSJed Brown 
420c4762a1bSJed Brown   for (i = xs; i < xs + xm + 1; i++) {
421c4762a1bSJed Brown     PetscReal    maxspeed;
422c4762a1bSJed Brown     PetscScalar *uL, *uR;
423c4762a1bSJed Brown     uL = &ctx->uLR[0];
424c4762a1bSJed Brown     uR = &ctx->uLR[dof];
425c4762a1bSJed Brown     if (i < sf - lsbwidth) { /* slow region */
426c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
427c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
428c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
429c4762a1bSJed Brown       }
4309566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
431c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hxs)); /* Max allowable value of 1/Delta t */
432c4762a1bSJed Brown       if (i > xs) {
433c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
434c4762a1bSJed Brown       }
435c4762a1bSJed Brown       if (i < xs + xm) {
436c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
437c4762a1bSJed Brown         islow++;
438c4762a1bSJed Brown       }
439c4762a1bSJed Brown     }
440c4762a1bSJed Brown     if (i == sf - lsbwidth) { /* interface between the slow region and the fast region */
441c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
442c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
443c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
444c4762a1bSJed Brown       }
4459566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
446c4762a1bSJed Brown       if (i > xs) {
447c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
448c4762a1bSJed Brown       }
449c4762a1bSJed Brown     }
450c4762a1bSJed Brown     if (i == fs + rsbwidth) { /* slow region */
451c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
452c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
453c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
454c4762a1bSJed Brown       }
4559566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
456c4762a1bSJed Brown       if (i < xs + xm) {
457c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
458c4762a1bSJed Brown         islow++;
459c4762a1bSJed Brown       }
460c4762a1bSJed Brown     }
461c4762a1bSJed Brown     if (i > fs + rsbwidth) { /* slow region */
462c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
463c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
464c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
465c4762a1bSJed Brown       }
4669566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
467c4762a1bSJed Brown       if (i > xs) {
468c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
469c4762a1bSJed Brown       }
470c4762a1bSJed Brown       if (i < xs + xm) {
471c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
472c4762a1bSJed Brown         islow++;
473c4762a1bSJed Brown       }
474c4762a1bSJed Brown     }
475c4762a1bSJed Brown   }
4769566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
4779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
4789566063dSJacob Faibussowitsch   PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
4799566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &Xloc));
480712fec58SPierre Jolivet   PetscCall(MPIU_Allreduce(&cfl_idt, &ctx->cfl_idt, 1, MPIU_SCALAR, MPIU_MAX, PetscObjectComm((PetscObject)da)));
4813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
482c4762a1bSJed Brown }
483c4762a1bSJed Brown 
484d71ae5a4SJacob Faibussowitsch PetscErrorCode FVRHSFunctionslowbuffer_2WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx)
485d71ae5a4SJacob Faibussowitsch {
486c4762a1bSJed Brown   FVCtx       *ctx = (FVCtx *)vctx;
487c4762a1bSJed Brown   PetscInt     i, j, k, Mx, dof, xs, xm, islow = 0, sf = ctx->sf, fs = ctx->fs, lsbwidth = ctx->lsbwidth, rsbwidth = ctx->rsbwidth;
488c4762a1bSJed Brown   PetscReal    hxs, hxf;
489c4762a1bSJed Brown   PetscScalar *x, *f, *slope;
490c4762a1bSJed Brown   Vec          Xloc;
491c4762a1bSJed Brown   DM           da;
492c4762a1bSJed Brown 
493c4762a1bSJed Brown   PetscFunctionBeginUser;
4949566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
4959566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &Xloc));
4969566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
497c4762a1bSJed Brown   hxs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
498c4762a1bSJed Brown   hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
4999566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
5009566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
5019566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
5029566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Xloc, &x));
5039566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
5049566063dSJacob Faibussowitsch   PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
5059566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
506c4762a1bSJed Brown 
507c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
508c4762a1bSJed Brown     for (i = xs - 2; i < 0; i++) {
509c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
510c4762a1bSJed Brown     }
511c4762a1bSJed Brown     for (i = Mx; i < xs + xm + 2; i++) {
512c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
513c4762a1bSJed Brown     }
514c4762a1bSJed Brown   }
515c4762a1bSJed Brown   for (i = xs - 1; i < xs + xm + 1; i++) {
516c4762a1bSJed Brown     struct _LimitInfo info;
517c4762a1bSJed Brown     PetscScalar      *cjmpL, *cjmpR;
518c4762a1bSJed Brown     if ((i > sf - lsbwidth - 2 && i < sf + 1) || (i > fs - 2 && i < fs + rsbwidth + 1)) {
519c4762a1bSJed Brown       /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
5209566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
521c4762a1bSJed Brown       /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
5229566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
523c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
524c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
525c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
526c4762a1bSJed Brown         PetscScalar jmpL, jmpR;
527c4762a1bSJed Brown         jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
528c4762a1bSJed Brown         jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
529c4762a1bSJed Brown         for (k = 0; k < dof; k++) {
530c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
531c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
532c4762a1bSJed Brown         }
533c4762a1bSJed Brown       }
534c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
535c4762a1bSJed Brown       info.m   = dof;
536c4762a1bSJed Brown       info.hxs = hxs;
537c4762a1bSJed Brown       info.hxf = hxf;
538c4762a1bSJed Brown       (*ctx->limit2)(&info, cjmpL, cjmpR, ctx->sf, ctx->fs, i, ctx->cslope);
539c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
540c4762a1bSJed Brown         PetscScalar tmp = 0;
541c4762a1bSJed Brown         for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
542c4762a1bSJed Brown         slope[i * dof + j] = tmp;
543c4762a1bSJed Brown       }
544c4762a1bSJed Brown     }
545c4762a1bSJed Brown   }
546c4762a1bSJed Brown 
547c4762a1bSJed Brown   for (i = xs; i < xs + xm + 1; i++) {
548c4762a1bSJed Brown     PetscReal    maxspeed;
549c4762a1bSJed Brown     PetscScalar *uL, *uR;
550c4762a1bSJed Brown     uL = &ctx->uLR[0];
551c4762a1bSJed Brown     uR = &ctx->uLR[dof];
552c4762a1bSJed Brown     if (i == sf - lsbwidth) {
553c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
554c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
555c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
556c4762a1bSJed Brown       }
5579566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
558c4762a1bSJed Brown       if (i < xs + xm) {
559c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
560c4762a1bSJed Brown         islow++;
561c4762a1bSJed Brown       }
562c4762a1bSJed Brown     }
563c4762a1bSJed Brown     if (i > sf - lsbwidth && i < sf) {
564c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
565c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
566c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
567c4762a1bSJed Brown       }
5689566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
569c4762a1bSJed Brown       if (i > xs) {
570c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
571c4762a1bSJed Brown       }
572c4762a1bSJed Brown       if (i < xs + xm) {
573c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
574c4762a1bSJed Brown         islow++;
575c4762a1bSJed Brown       }
576c4762a1bSJed Brown     }
577c4762a1bSJed Brown     if (i == sf) { /* interface between the slow region and the fast region */
578c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
579c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
580c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
581c4762a1bSJed Brown       }
5829566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
583c4762a1bSJed Brown       if (i > xs) {
584c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
585c4762a1bSJed Brown       }
586c4762a1bSJed Brown     }
587c4762a1bSJed Brown     if (i == fs) { /* interface between the fast region and the slow region */
588c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
589c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
590c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
591c4762a1bSJed Brown       }
5929566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
593c4762a1bSJed Brown       if (i < xs + xm) {
594c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
595c4762a1bSJed Brown         islow++;
596c4762a1bSJed Brown       }
597c4762a1bSJed Brown     }
598c4762a1bSJed Brown     if (i > fs && i < fs + rsbwidth) {
599c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
600c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
601c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
602c4762a1bSJed Brown       }
6039566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
604c4762a1bSJed Brown       if (i > xs) {
605c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
606c4762a1bSJed Brown       }
607c4762a1bSJed Brown       if (i < xs + xm) {
608c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
609c4762a1bSJed Brown         islow++;
610c4762a1bSJed Brown       }
611c4762a1bSJed Brown     }
612c4762a1bSJed Brown     if (i == fs + rsbwidth) {
613c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
614c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
615c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
616c4762a1bSJed Brown       }
6179566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
618c4762a1bSJed Brown       if (i > xs) {
619c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
620c4762a1bSJed Brown       }
621c4762a1bSJed Brown     }
622c4762a1bSJed Brown   }
6239566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
6249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
6259566063dSJacob Faibussowitsch   PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
6269566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &Xloc));
6273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
628c4762a1bSJed Brown }
629c4762a1bSJed Brown 
630c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for fast  parts ----------------------------------- */
631d71ae5a4SJacob Faibussowitsch PetscErrorCode FVRHSFunctionfast_2WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx)
632d71ae5a4SJacob Faibussowitsch {
633c4762a1bSJed Brown   FVCtx       *ctx = (FVCtx *)vctx;
634c4762a1bSJed Brown   PetscInt     i, j, k, Mx, dof, xs, xm, ifast = 0, sf = ctx->sf, fs = ctx->fs;
635c4762a1bSJed Brown   PetscReal    hxs, hxf;
636c4762a1bSJed Brown   PetscScalar *x, *f, *slope;
637c4762a1bSJed Brown   Vec          Xloc;
638c4762a1bSJed Brown   DM           da;
639c4762a1bSJed Brown 
640c4762a1bSJed Brown   PetscFunctionBeginUser;
6419566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
6429566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &Xloc));
6439566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
644c4762a1bSJed Brown   hxs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
645c4762a1bSJed Brown   hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
6469566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
6479566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
6489566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
6499566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Xloc, &x));
6509566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
6519566063dSJacob Faibussowitsch   PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
6529566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
653c4762a1bSJed Brown 
654c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
655c4762a1bSJed Brown     for (i = xs - 2; i < 0; i++) {
656c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
657c4762a1bSJed Brown     }
658c4762a1bSJed Brown     for (i = Mx; i < xs + xm + 2; i++) {
659c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
660c4762a1bSJed Brown     }
661c4762a1bSJed Brown   }
662c4762a1bSJed Brown   for (i = xs - 1; i < xs + xm + 1; i++) {
663c4762a1bSJed Brown     struct _LimitInfo info;
664c4762a1bSJed Brown     PetscScalar      *cjmpL, *cjmpR;
665c4762a1bSJed Brown     if (i > sf - 2 && i < fs + 1) {
6669566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
6679566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
668c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
669c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
670c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
671c4762a1bSJed Brown         PetscScalar jmpL, jmpR;
672c4762a1bSJed Brown         jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
673c4762a1bSJed Brown         jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
674c4762a1bSJed Brown         for (k = 0; k < dof; k++) {
675c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
676c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
677c4762a1bSJed Brown         }
678c4762a1bSJed Brown       }
679c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
680c4762a1bSJed Brown       info.m   = dof;
681c4762a1bSJed Brown       info.hxs = hxs;
682c4762a1bSJed Brown       info.hxf = hxf;
683c4762a1bSJed Brown       (*ctx->limit2)(&info, cjmpL, cjmpR, ctx->sf, ctx->fs, i, ctx->cslope);
684c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
685c4762a1bSJed Brown         PetscScalar tmp = 0;
686c4762a1bSJed Brown         for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
687c4762a1bSJed Brown         slope[i * dof + j] = tmp;
688c4762a1bSJed Brown       }
689c4762a1bSJed Brown     }
690c4762a1bSJed Brown   }
691c4762a1bSJed Brown 
692c4762a1bSJed Brown   for (i = xs; i < xs + xm + 1; i++) {
693c4762a1bSJed Brown     PetscReal    maxspeed;
694c4762a1bSJed Brown     PetscScalar *uL, *uR;
695c4762a1bSJed Brown     uL = &ctx->uLR[0];
696c4762a1bSJed Brown     uR = &ctx->uLR[dof];
697c4762a1bSJed Brown     if (i == sf) { /* interface between the slow region and the fast region */
698c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
699c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
700c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
701c4762a1bSJed Brown       }
7029566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
703c4762a1bSJed Brown       if (i < xs + xm) {
704c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[ifast * dof + j] += ctx->flux[j] / hxf;
705c4762a1bSJed Brown         ifast++;
706c4762a1bSJed Brown       }
707c4762a1bSJed Brown     }
708c4762a1bSJed Brown     if (i > sf && i < fs) { /* fast region */
709c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
710c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
711c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
712c4762a1bSJed Brown       }
7139566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
714c4762a1bSJed Brown       if (i > xs) {
715c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(ifast - 1) * dof + j] -= ctx->flux[j] / hxf;
716c4762a1bSJed Brown       }
717c4762a1bSJed Brown       if (i < xs + xm) {
718c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[ifast * dof + j] += ctx->flux[j] / hxf;
719c4762a1bSJed Brown         ifast++;
720c4762a1bSJed Brown       }
721c4762a1bSJed Brown     }
722c4762a1bSJed Brown     if (i == fs) { /* interface between the fast region and the slow region */
723c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
724c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
725c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
726c4762a1bSJed Brown       }
7279566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
728c4762a1bSJed Brown       if (i > xs) {
729c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(ifast - 1) * dof + j] -= ctx->flux[j] / hxf;
730c4762a1bSJed Brown       }
731c4762a1bSJed Brown     }
732c4762a1bSJed Brown   }
7339566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
7349566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
7359566063dSJacob Faibussowitsch   PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
7369566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &Xloc));
7373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
738c4762a1bSJed Brown }
739c4762a1bSJed Brown 
740d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[])
741d71ae5a4SJacob Faibussowitsch {
742c4762a1bSJed Brown   char              lname[256] = "mc", physname[256] = "advect", final_fname[256] = "solution.m";
743c4762a1bSJed Brown   PetscFunctionList limiters = 0, physics = 0;
744c4762a1bSJed Brown   MPI_Comm          comm;
745c4762a1bSJed Brown   TS                ts;
746c4762a1bSJed Brown   DM                da;
747c4762a1bSJed Brown   Vec               X, X0, R;
748c4762a1bSJed Brown   FVCtx             ctx;
749c4762a1bSJed Brown   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;
750c4762a1bSJed Brown   PetscBool         view_final = PETSC_FALSE;
751c4762a1bSJed Brown   PetscReal         ptime;
752c4762a1bSJed Brown 
753327415f7SBarry Smith   PetscFunctionBeginUser;
7549566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, 0, help));
755c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
7569566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&ctx, sizeof(ctx)));
757c4762a1bSJed Brown 
758c4762a1bSJed Brown   /* Register limiters to be available on the command line */
7599566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "upwind", Limit2_Upwind));
7609566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "lax-wendroff", Limit2_LaxWendroff));
7619566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "beam-warming", Limit2_BeamWarming));
7629566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "fromm", Limit2_Fromm));
7639566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "minmod", Limit2_Minmod));
7649566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "superbee", Limit2_Superbee));
7659566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "mc", Limit2_MC));
7669566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "koren3", Limit2_Koren3));
767c4762a1bSJed Brown 
768c4762a1bSJed Brown   /* Register physical models to be available on the command line */
7699566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&physics, "advect", PhysicsCreate_Advect));
770c4762a1bSJed Brown 
771c4762a1bSJed Brown   ctx.comm   = comm;
772c4762a1bSJed Brown   ctx.cfl    = 0.9;
773c4762a1bSJed Brown   ctx.bctype = FVBC_PERIODIC;
774c4762a1bSJed Brown   ctx.xmin   = -1.0;
775c4762a1bSJed Brown   ctx.xmax   = 1.0;
776d0609cedSBarry Smith   PetscOptionsBegin(comm, NULL, "Finite Volume solver options", "");
7779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-xmin", "X min", "", ctx.xmin, &ctx.xmin, NULL));
7789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-xmax", "X max", "", ctx.xmax, &ctx.xmax, NULL));
7799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-limit", "Name of flux imiter to use", "", limiters, lname, lname, sizeof(lname), NULL));
7809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-draw", "Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)", "", draw, &draw, NULL));
7819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-view_final", "Write final solution in ASCII MATLAB format to given file name", "", final_fname, final_fname, sizeof(final_fname), &view_final));
7829566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-initial", "Initial condition (depends on the physics)", "", ctx.initial, &ctx.initial, NULL));
7839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-exact", "Compare errors with exact solution", "", ctx.exact, &ctx.exact, NULL));
7849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-simulation", "Compare errors with reference solution", "", ctx.simulation, &ctx.simulation, NULL));
7859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-cfl", "CFL number to time step at", "", ctx.cfl, &ctx.cfl, NULL));
7869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-bc_type", "Boundary condition", "", FVBCTypes, (PetscEnum)ctx.bctype, (PetscEnum *)&ctx.bctype, NULL));
7879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-hratio", "Spacing ratio", "", ctx.hratio, &ctx.hratio, NULL));
788d0609cedSBarry Smith   PetscOptionsEnd();
789c4762a1bSJed Brown 
790c4762a1bSJed Brown   /* Choose the limiter from the list of registered limiters */
7919566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(limiters, lname, &ctx.limit2));
7923c633725SBarry Smith   PetscCheck(ctx.limit2, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Limiter '%s' not found", lname);
793c4762a1bSJed Brown 
794c4762a1bSJed Brown   /* Choose the physics from the list of registered models */
795c4762a1bSJed Brown   {
796c4762a1bSJed Brown     PetscErrorCode (*r)(FVCtx *);
7979566063dSJacob Faibussowitsch     PetscCall(PetscFunctionListFind(physics, physname, &r));
7983c633725SBarry Smith     PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Physics '%s' not found", physname);
799c4762a1bSJed Brown     /* Create the physics, will set the number of fields and their names */
8009566063dSJacob Faibussowitsch     PetscCall((*r)(&ctx));
801c4762a1bSJed Brown   }
802c4762a1bSJed Brown 
803c4762a1bSJed Brown   /* Create a DMDA to manage the parallel grid */
8049566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(comm, DM_BOUNDARY_PERIODIC, 50, ctx.physics2.dof, 2, NULL, &da));
8059566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
8069566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
807c4762a1bSJed Brown   /* Inform the DMDA of the field names provided by the physics. */
808c4762a1bSJed Brown   /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */
80948a46eb9SPierre Jolivet   for (i = 0; i < ctx.physics2.dof; i++) PetscCall(DMDASetFieldName(da, i, ctx.physics2.fieldname[i]));
8109566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
8119566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
812c4762a1bSJed Brown 
813c4762a1bSJed Brown   /* Set coordinates of cell centers */
8149566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da, ctx.xmin + 0.5 * (ctx.xmax - ctx.xmin) / Mx, ctx.xmax + 0.5 * (ctx.xmax - ctx.xmin) / Mx, 0, 0, 0, 0));
815c4762a1bSJed Brown 
816c4762a1bSJed Brown   /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */
8179566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(dof * dof, &ctx.R, dof * dof, &ctx.Rinv, 2 * dof, &ctx.cjmpLR, 1 * dof, &ctx.cslope));
8189566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(2 * dof, &ctx.uLR, dof, &ctx.flux, dof, &ctx.speeds));
819c4762a1bSJed Brown 
820c4762a1bSJed Brown   /* Create a vector to store the solution and to save the initial state */
8219566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &X));
8229566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X, &X0));
8239566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X, &R));
824c4762a1bSJed Brown 
825c4762a1bSJed Brown   /* create index for slow parts and fast parts,
826c4762a1bSJed Brown      count_slow + count_fast = Mx, counts_slow*hs = 0.5, counts_fast*hf = 0.5 */
827c4762a1bSJed Brown   count_slow = Mx / (1.0 + ctx.hratio / 3.0);
82808401ef6SPierre Jolivet   PetscCheck(count_slow % 2 == 0, PETSC_COMM_WORLD, PETSC_ERR_USER, "Please adjust grid size Mx (-da_grid_x) and hratio (-hratio) so that Mx/(1+hartio/3) is even");
829c4762a1bSJed Brown   count_fast = Mx - count_slow;
830c4762a1bSJed Brown   ctx.sf     = count_slow / 2;
831c4762a1bSJed Brown   ctx.fs     = ctx.sf + count_fast;
8329566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(xm * dof, &index_slow));
8339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(xm * dof, &index_fast));
8349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(6 * dof, &index_slowbuffer));
835c4762a1bSJed Brown   if (((AdvectCtx *)ctx.physics2.user)->a > 0) {
836c4762a1bSJed Brown     ctx.lsbwidth = 2;
837c4762a1bSJed Brown     ctx.rsbwidth = 4;
838c4762a1bSJed Brown   } else {
839c4762a1bSJed Brown     ctx.lsbwidth = 4;
840c4762a1bSJed Brown     ctx.rsbwidth = 2;
841c4762a1bSJed Brown   }
842c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
843c4762a1bSJed Brown     if (i < ctx.sf - ctx.lsbwidth || i > ctx.fs + ctx.rsbwidth - 1)
844c4762a1bSJed Brown       for (k = 0; k < dof; k++) index_slow[islow++] = i * dof + k;
845c4762a1bSJed Brown     else if ((i >= ctx.sf - ctx.lsbwidth && i < ctx.sf) || (i > ctx.fs - 1 && i <= ctx.fs + ctx.rsbwidth - 1))
846c4762a1bSJed Brown       for (k = 0; k < dof; k++) index_slowbuffer[islowbuffer++] = i * dof + k;
847c4762a1bSJed Brown     else
848c4762a1bSJed Brown       for (k = 0; k < dof; k++) index_fast[ifast++] = i * dof + k;
849c4762a1bSJed Brown   }
8509566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islow, index_slow, PETSC_COPY_VALUES, &ctx.iss));
8519566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, ifast, index_fast, PETSC_COPY_VALUES, &ctx.isf));
8529566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islowbuffer, index_slowbuffer, PETSC_COPY_VALUES, &ctx.issb));
853c4762a1bSJed Brown 
854c4762a1bSJed Brown   /* Create a time-stepping object */
8559566063dSJacob Faibussowitsch   PetscCall(TSCreate(comm, &ts));
8569566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, da));
8579566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, R, FVRHSFunction_2WaySplit, &ctx));
8589566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetIS(ts, "slow", ctx.iss));
8599566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetIS(ts, "slowbuffer", ctx.issb));
8609566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetIS(ts, "fast", ctx.isf));
8619566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetRHSFunction(ts, "slow", NULL, FVRHSFunctionslow_2WaySplit, &ctx));
8629566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetRHSFunction(ts, "fast", NULL, FVRHSFunctionfast_2WaySplit, &ctx));
8639566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetRHSFunction(ts, "slowbuffer", NULL, FVRHSFunctionslowbuffer_2WaySplit, &ctx));
864c4762a1bSJed Brown 
8659566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSSSP));
8669566063dSJacob Faibussowitsch   /*PetscCall(TSSetType(ts,TSMPRK));*/
8679566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 10));
8689566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
869c4762a1bSJed Brown 
870c4762a1bSJed Brown   /* Compute initial conditions and starting time step */
8719566063dSJacob Faibussowitsch   PetscCall(FVSample_2WaySplit(&ctx, da, 0, X0));
8729566063dSJacob Faibussowitsch   PetscCall(FVRHSFunction_2WaySplit(ts, 0, X0, X, (void *)&ctx)); /* Initial function evaluation, only used to determine max speed */
8739566063dSJacob Faibussowitsch   PetscCall(VecCopy(X0, X));                                      /* The function value was not used so we set X=X0 again */
8749566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, ctx.cfl / ctx.cfl_idt));
8759566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts)); /* Take runtime options */
8769566063dSJacob Faibussowitsch   PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD));
877c4762a1bSJed Brown   {
878c4762a1bSJed Brown     PetscInt           steps;
879c4762a1bSJed Brown     PetscScalar        mass_initial, mass_final, mass_difference, mass_differenceg;
880c4762a1bSJed Brown     const PetscScalar *ptr_X, *ptr_X0;
881c4762a1bSJed Brown     const PetscReal    hs = (ctx.xmax - ctx.xmin) * 3.0 / 4.0 / count_slow;
882c4762a1bSJed Brown     const PetscReal    hf = (ctx.xmax - ctx.xmin) / 4.0 / count_fast;
883c4762a1bSJed Brown 
8849566063dSJacob Faibussowitsch     PetscCall(TSSolve(ts, X));
8859566063dSJacob Faibussowitsch     PetscCall(TSGetSolveTime(ts, &ptime));
8869566063dSJacob Faibussowitsch     PetscCall(TSGetStepNumber(ts, &steps));
887c4762a1bSJed Brown     /* calculate the total mass at initial time and final time */
888c4762a1bSJed Brown     mass_initial = 0.0;
889c4762a1bSJed Brown     mass_final   = 0.0;
8909566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArrayRead(da, X0, (void *)&ptr_X0));
8919566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArrayRead(da, X, (void *)&ptr_X));
892c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
893c4762a1bSJed Brown       if (i < ctx.sf || i > ctx.fs - 1) {
894c4762a1bSJed Brown         for (k = 0; k < dof; k++) {
895c4762a1bSJed Brown           mass_initial = mass_initial + hs * ptr_X0[i * dof + k];
896c4762a1bSJed Brown           mass_final   = mass_final + hs * ptr_X[i * dof + k];
897c4762a1bSJed Brown         }
898c4762a1bSJed Brown       } else {
899c4762a1bSJed Brown         for (k = 0; k < dof; k++) {
900c4762a1bSJed Brown           mass_initial = mass_initial + hf * ptr_X0[i * dof + k];
901c4762a1bSJed Brown           mass_final   = mass_final + hf * ptr_X[i * dof + k];
902c4762a1bSJed Brown         }
903c4762a1bSJed Brown       }
904c4762a1bSJed Brown     }
9059566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArrayRead(da, X0, (void *)&ptr_X0));
9069566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArrayRead(da, X, (void *)&ptr_X));
907c4762a1bSJed Brown     mass_difference = mass_final - mass_initial;
908712fec58SPierre Jolivet     PetscCall(MPIU_Allreduce(&mass_difference, &mass_differenceg, 1, MPIU_SCALAR, MPIU_SUM, comm));
9099566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Mass difference %g\n", (double)mass_differenceg));
91063a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Final time %g, steps %" PetscInt_FMT "\n", (double)ptime, steps));
9119566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Maximum allowable stepsize according to CFL %g\n", (double)(1.0 / ctx.cfl_idt)));
912c4762a1bSJed Brown     if (ctx.exact) {
913c4762a1bSJed Brown       PetscReal nrm1 = 0;
9149566063dSJacob Faibussowitsch       PetscCall(SolutionErrorNorms_2WaySplit(&ctx, da, ptime, X, &nrm1));
9159566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %g\n", (double)nrm1));
916c4762a1bSJed Brown     }
917c4762a1bSJed Brown     if (ctx.simulation) {
918c4762a1bSJed Brown       PetscReal          nrm1 = 0;
919c4762a1bSJed Brown       PetscViewer        fd;
920c4762a1bSJed Brown       char               filename[PETSC_MAX_PATH_LEN] = "binaryoutput";
921c4762a1bSJed Brown       Vec                XR;
922c4762a1bSJed Brown       PetscBool          flg;
923c4762a1bSJed Brown       const PetscScalar *ptr_XR;
9249566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetString(NULL, NULL, "-f", filename, sizeof(filename), &flg));
9253c633725SBarry Smith       PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate binary file with the -f option");
9269566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, filename, FILE_MODE_READ, &fd));
9279566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(X0, &XR));
9289566063dSJacob Faibussowitsch       PetscCall(VecLoad(XR, fd));
9299566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&fd));
9309566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(X, &ptr_X));
9319566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(XR, &ptr_XR));
932c4762a1bSJed Brown       for (i = xs; i < xs + xm; i++) {
933c4762a1bSJed Brown         if (i < ctx.sf || i > ctx.fs - 1)
934c4762a1bSJed Brown           for (k = 0; k < dof; k++) nrm1 = nrm1 + hs * PetscAbs(ptr_X[i * dof + k] - ptr_XR[i * dof + k]);
935c4762a1bSJed Brown         else
936c4762a1bSJed Brown           for (k = 0; k < dof; k++) nrm1 = nrm1 + hf * PetscAbs(ptr_X[i * dof + k] - ptr_XR[i * dof + k]);
937c4762a1bSJed Brown       }
9389566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(X, &ptr_X));
9399566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(XR, &ptr_XR));
9409566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %g\n", (double)nrm1));
9419566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&XR));
942c4762a1bSJed Brown     }
943c4762a1bSJed Brown   }
944c4762a1bSJed Brown 
9459566063dSJacob Faibussowitsch   PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD));
9469566063dSJacob Faibussowitsch   if (draw & 0x1) PetscCall(VecView(X0, PETSC_VIEWER_DRAW_WORLD));
9479566063dSJacob Faibussowitsch   if (draw & 0x2) PetscCall(VecView(X, PETSC_VIEWER_DRAW_WORLD));
948c4762a1bSJed Brown   if (draw & 0x4) {
949c4762a1bSJed Brown     Vec Y;
9509566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(X, &Y));
9519566063dSJacob Faibussowitsch     PetscCall(FVSample_2WaySplit(&ctx, da, ptime, Y));
9529566063dSJacob Faibussowitsch     PetscCall(VecAYPX(Y, -1, X));
9539566063dSJacob Faibussowitsch     PetscCall(VecView(Y, PETSC_VIEWER_DRAW_WORLD));
9549566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&Y));
955c4762a1bSJed Brown   }
956c4762a1bSJed Brown 
957c4762a1bSJed Brown   if (view_final) {
958c4762a1bSJed Brown     PetscViewer viewer;
9599566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, final_fname, &viewer));
9609566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB));
9619566063dSJacob Faibussowitsch     PetscCall(VecView(X, viewer));
9629566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(viewer));
9639566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
964c4762a1bSJed Brown   }
965c4762a1bSJed Brown 
966c4762a1bSJed Brown   /* Clean up */
9679566063dSJacob Faibussowitsch   PetscCall((*ctx.physics2.destroy)(ctx.physics2.user));
9689566063dSJacob Faibussowitsch   for (i = 0; i < ctx.physics2.dof; i++) PetscCall(PetscFree(ctx.physics2.fieldname[i]));
9699566063dSJacob Faibussowitsch   PetscCall(PetscFree4(ctx.R, ctx.Rinv, ctx.cjmpLR, ctx.cslope));
9709566063dSJacob Faibussowitsch   PetscCall(PetscFree3(ctx.uLR, ctx.flux, ctx.speeds));
9719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
9729566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X0));
9739566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&R));
9749566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
9759566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
9769566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ctx.iss));
9779566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ctx.isf));
9789566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ctx.issb));
9799566063dSJacob Faibussowitsch   PetscCall(PetscFree(index_slow));
9809566063dSJacob Faibussowitsch   PetscCall(PetscFree(index_fast));
9819566063dSJacob Faibussowitsch   PetscCall(PetscFree(index_slowbuffer));
9829566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&limiters));
9839566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&physics));
9849566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
985b122ec5aSJacob Faibussowitsch   return 0;
986c4762a1bSJed Brown }
987c4762a1bSJed Brown 
988c4762a1bSJed Brown /*TEST
989c4762a1bSJed Brown 
990c4762a1bSJed Brown     build:
991f56ea12dSJed Brown       requires: !complex
992c4762a1bSJed Brown       depends: finitevolume1d.c
993c4762a1bSJed Brown 
994c4762a1bSJed Brown     test:
995c4762a1bSJed Brown       suffix: 1
996c4762a1bSJed Brown       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
997c4762a1bSJed Brown 
998c4762a1bSJed Brown     test:
999c4762a1bSJed Brown       suffix: 2
1000c4762a1bSJed Brown       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
1001c4762a1bSJed Brown       output_file: output/ex6_1.out
1002c4762a1bSJed Brown 
1003c4762a1bSJed Brown TEST*/
1004