xref: /petsc/src/ts/tutorials/multirate/ex6.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
279371c9d4SSatish Balay static inline PetscReal RangeMod(PetscReal a, PetscReal xmin, PetscReal xmax) {
289371c9d4SSatish Balay   PetscReal range = xmax - xmin;
299371c9d4SSatish Balay   return xmin + PetscFmodReal(range + PetscFmodReal(a, range), range);
309371c9d4SSatish Balay }
31c4762a1bSJed Brown 
32c4762a1bSJed Brown /* --------------------------------- Advection ----------------------------------- */
33c4762a1bSJed Brown typedef struct {
34c4762a1bSJed Brown   PetscReal a; /* advective velocity */
35c4762a1bSJed Brown } AdvectCtx;
36c4762a1bSJed Brown 
379371c9d4SSatish Balay static PetscErrorCode PhysicsRiemann_Advect(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed) {
38c4762a1bSJed Brown   AdvectCtx *ctx = (AdvectCtx *)vctx;
39c4762a1bSJed Brown   PetscReal  speed;
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   PetscFunctionBeginUser;
42c4762a1bSJed Brown   speed     = ctx->a;
43c4762a1bSJed Brown   flux[0]   = PetscMax(0, speed) * uL[0] + PetscMin(0, speed) * uR[0];
44c4762a1bSJed Brown   *maxspeed = speed;
45c4762a1bSJed Brown   PetscFunctionReturn(0);
46c4762a1bSJed Brown }
47c4762a1bSJed Brown 
489371c9d4SSatish Balay static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx, PetscInt m, const PetscScalar *u, PetscScalar *X, PetscScalar *Xi, PetscReal *speeds) {
49c4762a1bSJed Brown   AdvectCtx *ctx = (AdvectCtx *)vctx;
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   PetscFunctionBeginUser;
52c4762a1bSJed Brown   X[0]      = 1.;
53c4762a1bSJed Brown   Xi[0]     = 1.;
54c4762a1bSJed Brown   speeds[0] = ctx->a;
55c4762a1bSJed Brown   PetscFunctionReturn(0);
56c4762a1bSJed Brown }
57c4762a1bSJed Brown 
589371c9d4SSatish Balay static PetscErrorCode PhysicsSample_Advect(void *vctx, PetscInt initial, FVBCType bctype, PetscReal xmin, PetscReal xmax, PetscReal t, PetscReal x, PetscReal *u) {
59c4762a1bSJed Brown   AdvectCtx *ctx = (AdvectCtx *)vctx;
60c4762a1bSJed Brown   PetscReal  a   = ctx->a, x0;
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   PetscFunctionBeginUser;
63c4762a1bSJed Brown   switch (bctype) {
64c4762a1bSJed Brown   case FVBC_OUTFLOW: x0 = x - a * t; break;
65c4762a1bSJed Brown   case FVBC_PERIODIC: x0 = RangeMod(x - a * t, xmin, xmax); break;
66c4762a1bSJed Brown   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown BCType");
67c4762a1bSJed Brown   }
68c4762a1bSJed Brown   switch (initial) {
69c4762a1bSJed Brown   case 0: u[0] = (x0 < 0) ? 1 : -1; break;
70c4762a1bSJed Brown   case 1: u[0] = (x0 < 0) ? -1 : 1; break;
71c4762a1bSJed Brown   case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break;
72c4762a1bSJed Brown   case 3: u[0] = PetscSinReal(2 * PETSC_PI * x0); break;
73c4762a1bSJed Brown   case 4: u[0] = PetscAbs(x0); break;
74c4762a1bSJed Brown   case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2 * PETSC_PI * x0)); break;
75c4762a1bSJed Brown   case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2 - x0 : 0)); break;
76c4762a1bSJed Brown   case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI * x0), 10.0); break;
77c4762a1bSJed Brown   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition");
78c4762a1bSJed Brown   }
79c4762a1bSJed Brown   PetscFunctionReturn(0);
80c4762a1bSJed Brown }
81c4762a1bSJed Brown 
829371c9d4SSatish Balay static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx) {
83c4762a1bSJed Brown   AdvectCtx *user;
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   PetscFunctionBeginUser;
869566063dSJacob Faibussowitsch   PetscCall(PetscNew(&user));
87c4762a1bSJed Brown   ctx->physics2.sample2         = PhysicsSample_Advect;
88c4762a1bSJed Brown   ctx->physics2.riemann2        = PhysicsRiemann_Advect;
89c4762a1bSJed Brown   ctx->physics2.characteristic2 = PhysicsCharacteristic_Advect;
90c4762a1bSJed Brown   ctx->physics2.destroy         = PhysicsDestroy_SimpleFree;
91c4762a1bSJed Brown   ctx->physics2.user            = user;
92c4762a1bSJed Brown   ctx->physics2.dof             = 1;
939566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("u", &ctx->physics2.fieldname[0]));
94c4762a1bSJed Brown   user->a = 1;
95d0609cedSBarry Smith   PetscOptionsBegin(ctx->comm, ctx->prefix, "Options for advection", "");
969371c9d4SSatish Balay   { PetscCall(PetscOptionsReal("-physics_advect_a", "Speed", "", user->a, &user->a, NULL)); }
97d0609cedSBarry Smith   PetscOptionsEnd();
98c4762a1bSJed Brown   PetscFunctionReturn(0);
99c4762a1bSJed Brown }
100c4762a1bSJed Brown 
1019371c9d4SSatish Balay PetscErrorCode FVSample_2WaySplit(FVCtx *ctx, DM da, PetscReal time, Vec U) {
102c4762a1bSJed Brown   PetscScalar   *u, *uj, xj, xi;
103c4762a1bSJed Brown   PetscInt       i, j, k, dof, xs, xm, Mx;
104c4762a1bSJed Brown   const PetscInt N = 200;
105c4762a1bSJed Brown   PetscReal      hs, hf;
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   PetscFunctionBeginUser;
1083c633725SBarry Smith   PetscCheck(ctx->physics2.sample2, PETSC_COMM_SELF, PETSC_ERR_SUP, "Physics has not provided a sampling function");
1099566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
1109566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
1119566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, U, &u));
1129566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dof, &uj));
113c4762a1bSJed Brown   hs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
114c4762a1bSJed Brown   hf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
115c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
116c4762a1bSJed Brown     if (i < ctx->sf) {
117c4762a1bSJed Brown       xi = ctx->xmin + 0.5 * hs + i * hs;
118c4762a1bSJed Brown       /* Integrate over cell i using trapezoid rule with N points. */
119c4762a1bSJed Brown       for (k = 0; k < dof; k++) u[i * dof + k] = 0;
120c4762a1bSJed Brown       for (j = 0; j < N + 1; j++) {
121c4762a1bSJed Brown         xj = xi + hs * (j - N / 2) / (PetscReal)N;
1229566063dSJacob Faibussowitsch         PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
123c4762a1bSJed Brown         for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
124c4762a1bSJed Brown       }
125c4762a1bSJed Brown     } else if (i < ctx->fs) {
126c4762a1bSJed Brown       xi = ctx->xmin + ctx->sf * hs + 0.5 * hf + (i - ctx->sf) * hf;
127c4762a1bSJed Brown       /* Integrate over cell i using trapezoid rule with N points. */
128c4762a1bSJed Brown       for (k = 0; k < dof; k++) u[i * dof + k] = 0;
129c4762a1bSJed Brown       for (j = 0; j < N + 1; j++) {
130c4762a1bSJed Brown         xj = xi + hf * (j - N / 2) / (PetscReal)N;
1319566063dSJacob Faibussowitsch         PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
132c4762a1bSJed Brown         for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
133c4762a1bSJed Brown       }
134c4762a1bSJed Brown     } else {
135c4762a1bSJed Brown       xi = ctx->xmin + ctx->sf * hs + (ctx->fs - ctx->sf) * hf + 0.5 * hs + (i - ctx->fs) * hs;
136c4762a1bSJed Brown       /* Integrate over cell i using trapezoid rule with N points. */
137c4762a1bSJed Brown       for (k = 0; k < dof; k++) u[i * dof + k] = 0;
138c4762a1bSJed Brown       for (j = 0; j < N + 1; j++) {
139c4762a1bSJed Brown         xj = xi + hs * (j - N / 2) / (PetscReal)N;
1409566063dSJacob Faibussowitsch         PetscCall((*ctx->physics2.sample2)(ctx->physics2.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
141c4762a1bSJed Brown         for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
142c4762a1bSJed Brown       }
143c4762a1bSJed Brown     }
144c4762a1bSJed Brown   }
1459566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, U, &u));
1469566063dSJacob Faibussowitsch   PetscCall(PetscFree(uj));
147c4762a1bSJed Brown   PetscFunctionReturn(0);
148c4762a1bSJed Brown }
149c4762a1bSJed Brown 
1509371c9d4SSatish Balay static PetscErrorCode SolutionErrorNorms_2WaySplit(FVCtx *ctx, DM da, PetscReal t, Vec X, PetscReal *nrm1) {
151c4762a1bSJed Brown   Vec                Y;
152c4762a1bSJed Brown   PetscInt           i, Mx;
153c4762a1bSJed Brown   const PetscScalar *ptr_X, *ptr_Y;
154c4762a1bSJed Brown   PetscReal          hs, hf;
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   PetscFunctionBeginUser;
1579566063dSJacob Faibussowitsch   PetscCall(VecGetSize(X, &Mx));
1589566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X, &Y));
1599566063dSJacob Faibussowitsch   PetscCall(FVSample_2WaySplit(ctx, da, t, Y));
160c4762a1bSJed Brown   hs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
161c4762a1bSJed Brown   hf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
1629566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &ptr_X));
1639566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &ptr_Y));
164c4762a1bSJed Brown   for (i = 0; i < Mx; i++) {
165c4762a1bSJed Brown     if (i < ctx->sf || i > ctx->fs - 1) *nrm1 += hs * PetscAbs(ptr_X[i] - ptr_Y[i]);
166c4762a1bSJed Brown     else *nrm1 += hf * PetscAbs(ptr_X[i] - ptr_Y[i]);
167c4762a1bSJed Brown   }
1689566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &ptr_X));
1699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &ptr_Y));
1709566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Y));
171c4762a1bSJed Brown   PetscFunctionReturn(0);
172c4762a1bSJed Brown }
173c4762a1bSJed Brown 
1749371c9d4SSatish Balay PetscErrorCode FVRHSFunction_2WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) {
175c4762a1bSJed Brown   FVCtx       *ctx = (FVCtx *)vctx;
176c4762a1bSJed Brown   PetscInt     i, j, k, Mx, dof, xs, xm, sf = ctx->sf, fs = ctx->fs;
177c4762a1bSJed Brown   PetscReal    hxf, hxs, cfl_idt = 0;
178c4762a1bSJed Brown   PetscScalar *x, *f, *slope;
179c4762a1bSJed Brown   Vec          Xloc;
180c4762a1bSJed Brown   DM           da;
181c4762a1bSJed Brown 
182c4762a1bSJed Brown   PetscFunctionBeginUser;
1839566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
1849566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &Xloc));                                 /* Xloc contains ghost points                                     */
1859566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); /* Mx is the number of center points                              */
186c4762a1bSJed Brown   hxs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
187c4762a1bSJed Brown   hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
1889566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); /* X is solution vector which does not contain ghost points       */
1899566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
190c4762a1bSJed Brown 
1919566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F)); /* F is the right hand side function corresponds to center points */
192c4762a1bSJed Brown 
1939566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Xloc, &x));
1949566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, F, &f));
1959566063dSJacob Faibussowitsch   PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); /* contains ghost points                                           */
196c4762a1bSJed Brown 
1979566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
198c4762a1bSJed Brown 
199c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
200c4762a1bSJed Brown     for (i = xs - 2; i < 0; i++) {
201c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
202c4762a1bSJed Brown     }
203c4762a1bSJed Brown     for (i = Mx; i < xs + xm + 2; i++) {
204c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
205c4762a1bSJed Brown     }
206c4762a1bSJed Brown   }
207c4762a1bSJed Brown   for (i = xs - 1; i < xs + xm + 1; i++) {
208c4762a1bSJed Brown     struct _LimitInfo info;
209c4762a1bSJed Brown     PetscScalar      *cjmpL, *cjmpR;
210c4762a1bSJed Brown     /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
2119566063dSJacob Faibussowitsch     PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
212c4762a1bSJed Brown     /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
2139566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
214c4762a1bSJed Brown     cjmpL = &ctx->cjmpLR[0];
215c4762a1bSJed Brown     cjmpR = &ctx->cjmpLR[dof];
216c4762a1bSJed Brown     for (j = 0; j < dof; j++) {
217c4762a1bSJed Brown       PetscScalar jmpL, jmpR;
218c4762a1bSJed Brown       jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
219c4762a1bSJed Brown       jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
220c4762a1bSJed Brown       for (k = 0; k < dof; k++) {
221c4762a1bSJed Brown         cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
222c4762a1bSJed Brown         cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
223c4762a1bSJed Brown       }
224c4762a1bSJed Brown     }
225c4762a1bSJed Brown     /* Apply limiter to the left and right characteristic jumps */
226c4762a1bSJed Brown     info.m   = dof;
227c4762a1bSJed Brown     info.hxs = hxs;
228c4762a1bSJed Brown     info.hxf = hxf;
229c4762a1bSJed Brown     (*ctx->limit2)(&info, cjmpL, cjmpR, ctx->sf, ctx->fs, i, ctx->cslope);
230c4762a1bSJed Brown     for (j = 0; j < dof; j++) {
231c4762a1bSJed Brown       PetscScalar tmp = 0;
232c4762a1bSJed Brown       for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
233c4762a1bSJed Brown       slope[i * dof + j] = tmp;
234c4762a1bSJed Brown     }
235c4762a1bSJed Brown   }
236c4762a1bSJed Brown 
237c4762a1bSJed Brown   for (i = xs; i < xs + xm + 1; i++) {
238c4762a1bSJed Brown     PetscReal    maxspeed;
239c4762a1bSJed Brown     PetscScalar *uL, *uR;
240c4762a1bSJed Brown     uL = &ctx->uLR[0];
241c4762a1bSJed Brown     uR = &ctx->uLR[dof];
242c4762a1bSJed Brown     if (i < sf) { /* slow region */
243c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
244c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
245c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
246c4762a1bSJed Brown       }
2479566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
248c4762a1bSJed Brown       if (i > xs) {
249c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxs;
250c4762a1bSJed Brown       }
251c4762a1bSJed Brown       if (i < xs + xm) {
252c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxs;
253c4762a1bSJed Brown       }
254c4762a1bSJed Brown     } else if (i == sf) { /* interface between the slow region and the fast region */
255c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
256c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
257c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
258c4762a1bSJed Brown       }
2599566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
260c4762a1bSJed Brown       if (i > xs) {
261c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxs;
262c4762a1bSJed Brown       }
263c4762a1bSJed Brown       if (i < xs + xm) {
264c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxf;
265c4762a1bSJed Brown       }
266c4762a1bSJed Brown     } else if (i > sf && i < fs) { /* fast region */
267c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
268c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
269c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
270c4762a1bSJed Brown       }
2719566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
272c4762a1bSJed Brown       if (i > xs) {
273c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxf;
274c4762a1bSJed Brown       }
275c4762a1bSJed Brown       if (i < xs + xm) {
276c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxf;
277c4762a1bSJed Brown       }
278c4762a1bSJed Brown     } else if (i == fs) { /* interface between the fast region and the slow region */
279c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
280c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
281c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
282c4762a1bSJed Brown       }
2839566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
284c4762a1bSJed Brown       if (i > xs) {
285c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxf;
286c4762a1bSJed Brown       }
287c4762a1bSJed Brown       if (i < xs + xm) {
288c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxs;
289c4762a1bSJed Brown       }
290c4762a1bSJed Brown     } else { /* slow region */
291c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
292c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
293c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
294c4762a1bSJed Brown       }
2959566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
296c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hxs)); /* Max allowable value of 1/Delta t */
297c4762a1bSJed Brown       if (i > xs) {
298c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hxs;
299c4762a1bSJed Brown       }
300c4762a1bSJed Brown       if (i < xs + xm) {
301c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hxs;
302c4762a1bSJed Brown       }
303c4762a1bSJed Brown     }
304c4762a1bSJed Brown   }
3059566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
3069566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, F, &f));
3079566063dSJacob Faibussowitsch   PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
3089566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &Xloc));
3099566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&cfl_idt, &ctx->cfl_idt, 1, MPIU_SCALAR, MPIU_MAX, PetscObjectComm((PetscObject)da)));
310c4762a1bSJed Brown   if (0) {
311c4762a1bSJed Brown     /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */
312c4762a1bSJed Brown     PetscReal dt, tnow;
3139566063dSJacob Faibussowitsch     PetscCall(TSGetTimeStep(ts, &dt));
3149566063dSJacob Faibussowitsch     PetscCall(TSGetTime(ts, &tnow));
315c4762a1bSJed Brown     if (dt > 0.5 / ctx->cfl_idt) {
316c4762a1bSJed Brown       if (1) {
3179566063dSJacob 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)));
31898921bdaSJacob Faibussowitsch       } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Stability constraint exceeded, %g > %g", (double)dt, (double)(ctx->cfl / ctx->cfl_idt));
319c4762a1bSJed Brown     }
320c4762a1bSJed Brown   }
321c4762a1bSJed Brown   PetscFunctionReturn(0);
322c4762a1bSJed Brown }
323c4762a1bSJed Brown 
324c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */
3259371c9d4SSatish Balay PetscErrorCode FVRHSFunctionslow_2WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) {
326c4762a1bSJed Brown   FVCtx       *ctx = (FVCtx *)vctx;
327c4762a1bSJed Brown   PetscInt     i, j, k, Mx, dof, xs, xm, islow = 0, sf = ctx->sf, fs = ctx->fs, lsbwidth = ctx->lsbwidth, rsbwidth = ctx->rsbwidth;
328c4762a1bSJed Brown   PetscReal    hxs, hxf, cfl_idt = 0;
329c4762a1bSJed Brown   PetscScalar *x, *f, *slope;
330c4762a1bSJed Brown   Vec          Xloc;
331c4762a1bSJed Brown   DM           da;
332c4762a1bSJed Brown 
333c4762a1bSJed Brown   PetscFunctionBeginUser;
3349566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
3359566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &Xloc));
3369566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
337c4762a1bSJed Brown   hxs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
338c4762a1bSJed Brown   hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
3399566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
3409566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
3419566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
3429566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Xloc, &x));
3439566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
3449566063dSJacob Faibussowitsch   PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
3459566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
346c4762a1bSJed Brown 
347c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
348c4762a1bSJed Brown     for (i = xs - 2; i < 0; i++) {
349c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
350c4762a1bSJed Brown     }
351c4762a1bSJed Brown     for (i = Mx; i < xs + xm + 2; i++) {
352c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
353c4762a1bSJed Brown     }
354c4762a1bSJed Brown   }
355c4762a1bSJed Brown   for (i = xs - 1; i < xs + xm + 1; i++) {
356c4762a1bSJed Brown     struct _LimitInfo info;
357c4762a1bSJed Brown     PetscScalar      *cjmpL, *cjmpR;
358c4762a1bSJed Brown     if (i < sf - lsbwidth + 1 || i > fs + rsbwidth - 2) { /* slow components and the first and last fast components */
359c4762a1bSJed Brown       /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
3609566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
361c4762a1bSJed Brown       /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
3629566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
363c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
364c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
365c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
366c4762a1bSJed Brown         PetscScalar jmpL, jmpR;
367c4762a1bSJed Brown         jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
368c4762a1bSJed Brown         jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
369c4762a1bSJed Brown         for (k = 0; k < dof; k++) {
370c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
371c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
372c4762a1bSJed Brown         }
373c4762a1bSJed Brown       }
374c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
375c4762a1bSJed Brown       info.m   = dof;
376c4762a1bSJed Brown       info.hxs = hxs;
377c4762a1bSJed Brown       info.hxf = hxf;
378c4762a1bSJed Brown       (*ctx->limit2)(&info, cjmpL, cjmpR, ctx->sf, ctx->fs, i, ctx->cslope);
379c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
380c4762a1bSJed Brown         PetscScalar tmp = 0;
381c4762a1bSJed Brown         for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
382c4762a1bSJed Brown         slope[i * dof + j] = tmp;
383c4762a1bSJed Brown       }
384c4762a1bSJed Brown     }
385c4762a1bSJed Brown   }
386c4762a1bSJed Brown 
387c4762a1bSJed Brown   for (i = xs; i < xs + xm + 1; i++) {
388c4762a1bSJed Brown     PetscReal    maxspeed;
389c4762a1bSJed Brown     PetscScalar *uL, *uR;
390c4762a1bSJed Brown     uL = &ctx->uLR[0];
391c4762a1bSJed Brown     uR = &ctx->uLR[dof];
392c4762a1bSJed Brown     if (i < sf - lsbwidth) { /* slow region */
393c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
394c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
395c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
396c4762a1bSJed Brown       }
3979566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
398c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hxs)); /* Max allowable value of 1/Delta t */
399c4762a1bSJed Brown       if (i > xs) {
400c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
401c4762a1bSJed Brown       }
402c4762a1bSJed Brown       if (i < xs + xm) {
403c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
404c4762a1bSJed Brown         islow++;
405c4762a1bSJed Brown       }
406c4762a1bSJed Brown     }
407c4762a1bSJed Brown     if (i == sf - lsbwidth) { /* interface between the slow region and the fast region */
408c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
409c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
410c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
411c4762a1bSJed Brown       }
4129566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
413c4762a1bSJed Brown       if (i > xs) {
414c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
415c4762a1bSJed Brown       }
416c4762a1bSJed Brown     }
417c4762a1bSJed Brown     if (i == fs + rsbwidth) { /* slow region */
418c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
419c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
420c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
421c4762a1bSJed Brown       }
4229566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
423c4762a1bSJed Brown       if (i < xs + xm) {
424c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
425c4762a1bSJed Brown         islow++;
426c4762a1bSJed Brown       }
427c4762a1bSJed Brown     }
428c4762a1bSJed Brown     if (i > fs + rsbwidth) { /* slow region */
429c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
430c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
431c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
432c4762a1bSJed Brown       }
4339566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
434c4762a1bSJed Brown       if (i > xs) {
435c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
436c4762a1bSJed Brown       }
437c4762a1bSJed Brown       if (i < xs + xm) {
438c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
439c4762a1bSJed Brown         islow++;
440c4762a1bSJed Brown       }
441c4762a1bSJed Brown     }
442c4762a1bSJed Brown   }
4439566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
4449566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
4459566063dSJacob Faibussowitsch   PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
4469566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &Xloc));
4479566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&cfl_idt, &ctx->cfl_idt, 1, MPIU_SCALAR, MPIU_MAX, PetscObjectComm((PetscObject)da)));
448c4762a1bSJed Brown   PetscFunctionReturn(0);
449c4762a1bSJed Brown }
450c4762a1bSJed Brown 
4519371c9d4SSatish Balay PetscErrorCode FVRHSFunctionslowbuffer_2WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) {
452c4762a1bSJed Brown   FVCtx       *ctx = (FVCtx *)vctx;
453c4762a1bSJed Brown   PetscInt     i, j, k, Mx, dof, xs, xm, islow = 0, sf = ctx->sf, fs = ctx->fs, lsbwidth = ctx->lsbwidth, rsbwidth = ctx->rsbwidth;
454c4762a1bSJed Brown   PetscReal    hxs, hxf;
455c4762a1bSJed Brown   PetscScalar *x, *f, *slope;
456c4762a1bSJed Brown   Vec          Xloc;
457c4762a1bSJed Brown   DM           da;
458c4762a1bSJed Brown 
459c4762a1bSJed Brown   PetscFunctionBeginUser;
4609566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
4619566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &Xloc));
4629566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
463c4762a1bSJed Brown   hxs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
464c4762a1bSJed Brown   hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
4659566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
4669566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
4679566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
4689566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Xloc, &x));
4699566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
4709566063dSJacob Faibussowitsch   PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
4719566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
472c4762a1bSJed Brown 
473c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
474c4762a1bSJed Brown     for (i = xs - 2; i < 0; i++) {
475c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
476c4762a1bSJed Brown     }
477c4762a1bSJed Brown     for (i = Mx; i < xs + xm + 2; i++) {
478c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
479c4762a1bSJed Brown     }
480c4762a1bSJed Brown   }
481c4762a1bSJed Brown   for (i = xs - 1; i < xs + xm + 1; i++) {
482c4762a1bSJed Brown     struct _LimitInfo info;
483c4762a1bSJed Brown     PetscScalar      *cjmpL, *cjmpR;
484c4762a1bSJed Brown     if ((i > sf - lsbwidth - 2 && i < sf + 1) || (i > fs - 2 && i < fs + rsbwidth + 1)) {
485c4762a1bSJed Brown       /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
4869566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
487c4762a1bSJed Brown       /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
4889566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
489c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
490c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
491c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
492c4762a1bSJed Brown         PetscScalar jmpL, jmpR;
493c4762a1bSJed Brown         jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
494c4762a1bSJed Brown         jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
495c4762a1bSJed Brown         for (k = 0; k < dof; k++) {
496c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
497c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
498c4762a1bSJed Brown         }
499c4762a1bSJed Brown       }
500c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
501c4762a1bSJed Brown       info.m   = dof;
502c4762a1bSJed Brown       info.hxs = hxs;
503c4762a1bSJed Brown       info.hxf = hxf;
504c4762a1bSJed Brown       (*ctx->limit2)(&info, cjmpL, cjmpR, ctx->sf, ctx->fs, i, ctx->cslope);
505c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
506c4762a1bSJed Brown         PetscScalar tmp = 0;
507c4762a1bSJed Brown         for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
508c4762a1bSJed Brown         slope[i * dof + j] = tmp;
509c4762a1bSJed Brown       }
510c4762a1bSJed Brown     }
511c4762a1bSJed Brown   }
512c4762a1bSJed Brown 
513c4762a1bSJed Brown   for (i = xs; i < xs + xm + 1; i++) {
514c4762a1bSJed Brown     PetscReal    maxspeed;
515c4762a1bSJed Brown     PetscScalar *uL, *uR;
516c4762a1bSJed Brown     uL = &ctx->uLR[0];
517c4762a1bSJed Brown     uR = &ctx->uLR[dof];
518c4762a1bSJed Brown     if (i == sf - lsbwidth) {
519c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
520c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
521c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
522c4762a1bSJed Brown       }
5239566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
524c4762a1bSJed Brown       if (i < xs + xm) {
525c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
526c4762a1bSJed Brown         islow++;
527c4762a1bSJed Brown       }
528c4762a1bSJed Brown     }
529c4762a1bSJed Brown     if (i > sf - lsbwidth && i < sf) {
530c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
531c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
532c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
533c4762a1bSJed Brown       }
5349566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
535c4762a1bSJed Brown       if (i > xs) {
536c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
537c4762a1bSJed Brown       }
538c4762a1bSJed Brown       if (i < xs + xm) {
539c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
540c4762a1bSJed Brown         islow++;
541c4762a1bSJed Brown       }
542c4762a1bSJed Brown     }
543c4762a1bSJed Brown     if (i == sf) { /* interface between the slow region and the fast region */
544c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
545c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
546c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
547c4762a1bSJed Brown       }
5489566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
549c4762a1bSJed Brown       if (i > xs) {
550c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
551c4762a1bSJed Brown       }
552c4762a1bSJed Brown     }
553c4762a1bSJed Brown     if (i == fs) { /* interface between the fast region and the slow region */
554c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
555c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
556c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
557c4762a1bSJed Brown       }
5589566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
559c4762a1bSJed Brown       if (i < xs + xm) {
560c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
561c4762a1bSJed Brown         islow++;
562c4762a1bSJed Brown       }
563c4762a1bSJed Brown     }
564c4762a1bSJed Brown     if (i > fs && i < fs + rsbwidth) {
565c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
566c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
567c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
568c4762a1bSJed Brown       }
5699566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
570c4762a1bSJed Brown       if (i > xs) {
571c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
572c4762a1bSJed Brown       }
573c4762a1bSJed Brown       if (i < xs + xm) {
574c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hxs;
575c4762a1bSJed Brown         islow++;
576c4762a1bSJed Brown       }
577c4762a1bSJed Brown     }
578c4762a1bSJed Brown     if (i == fs + rsbwidth) {
579c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
580c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
581c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
582c4762a1bSJed Brown       }
5839566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
584c4762a1bSJed Brown       if (i > xs) {
585c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hxs;
586c4762a1bSJed Brown       }
587c4762a1bSJed Brown     }
588c4762a1bSJed Brown   }
5899566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
5909566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
5919566063dSJacob Faibussowitsch   PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
5929566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &Xloc));
593c4762a1bSJed Brown   PetscFunctionReturn(0);
594c4762a1bSJed Brown }
595c4762a1bSJed Brown 
596c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for fast  parts ----------------------------------- */
5979371c9d4SSatish Balay PetscErrorCode FVRHSFunctionfast_2WaySplit(TS ts, PetscReal time, Vec X, Vec F, void *vctx) {
598c4762a1bSJed Brown   FVCtx       *ctx = (FVCtx *)vctx;
599c4762a1bSJed Brown   PetscInt     i, j, k, Mx, dof, xs, xm, ifast = 0, sf = ctx->sf, fs = ctx->fs;
600c4762a1bSJed Brown   PetscReal    hxs, hxf;
601c4762a1bSJed Brown   PetscScalar *x, *f, *slope;
602c4762a1bSJed Brown   Vec          Xloc;
603c4762a1bSJed Brown   DM           da;
604c4762a1bSJed Brown 
605c4762a1bSJed Brown   PetscFunctionBeginUser;
6069566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
6079566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &Xloc));
6089566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
609c4762a1bSJed Brown   hxs = (ctx->xmax - ctx->xmin) * 3.0 / 8.0 / ctx->sf;
610c4762a1bSJed Brown   hxf = (ctx->xmax - ctx->xmin) / 4.0 / (ctx->fs - ctx->sf);
6119566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
6129566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
6139566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
6149566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Xloc, &x));
6159566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
6169566063dSJacob Faibussowitsch   PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
6179566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
618c4762a1bSJed Brown 
619c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
620c4762a1bSJed Brown     for (i = xs - 2; i < 0; i++) {
621c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
622c4762a1bSJed Brown     }
623c4762a1bSJed Brown     for (i = Mx; i < xs + xm + 2; i++) {
624c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
625c4762a1bSJed Brown     }
626c4762a1bSJed Brown   }
627c4762a1bSJed Brown   for (i = xs - 1; i < xs + xm + 1; i++) {
628c4762a1bSJed Brown     struct _LimitInfo info;
629c4762a1bSJed Brown     PetscScalar      *cjmpL, *cjmpR;
630c4762a1bSJed Brown     if (i > sf - 2 && i < fs + 1) {
6319566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.characteristic2)(ctx->physics2.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
6329566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
633c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
634c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
635c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
636c4762a1bSJed Brown         PetscScalar jmpL, jmpR;
637c4762a1bSJed Brown         jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
638c4762a1bSJed Brown         jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
639c4762a1bSJed Brown         for (k = 0; k < dof; k++) {
640c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
641c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
642c4762a1bSJed Brown         }
643c4762a1bSJed Brown       }
644c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
645c4762a1bSJed Brown       info.m   = dof;
646c4762a1bSJed Brown       info.hxs = hxs;
647c4762a1bSJed Brown       info.hxf = hxf;
648c4762a1bSJed Brown       (*ctx->limit2)(&info, cjmpL, cjmpR, ctx->sf, ctx->fs, i, ctx->cslope);
649c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
650c4762a1bSJed Brown         PetscScalar tmp = 0;
651c4762a1bSJed Brown         for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
652c4762a1bSJed Brown         slope[i * dof + j] = tmp;
653c4762a1bSJed Brown       }
654c4762a1bSJed Brown     }
655c4762a1bSJed Brown   }
656c4762a1bSJed Brown 
657c4762a1bSJed Brown   for (i = xs; i < xs + xm + 1; i++) {
658c4762a1bSJed Brown     PetscReal    maxspeed;
659c4762a1bSJed Brown     PetscScalar *uL, *uR;
660c4762a1bSJed Brown     uL = &ctx->uLR[0];
661c4762a1bSJed Brown     uR = &ctx->uLR[dof];
662c4762a1bSJed Brown     if (i == sf) { /* interface between the slow region and the fast region */
663c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
664c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxs / 2;
665c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
666c4762a1bSJed Brown       }
6679566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
668c4762a1bSJed Brown       if (i < xs + xm) {
669c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[ifast * dof + j] += ctx->flux[j] / hxf;
670c4762a1bSJed Brown         ifast++;
671c4762a1bSJed Brown       }
672c4762a1bSJed Brown     }
673c4762a1bSJed Brown     if (i > sf && i < fs) { /* fast region */
674c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
675c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
676c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxf / 2;
677c4762a1bSJed Brown       }
6789566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
679c4762a1bSJed Brown       if (i > xs) {
680c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(ifast - 1) * dof + j] -= ctx->flux[j] / hxf;
681c4762a1bSJed Brown       }
682c4762a1bSJed Brown       if (i < xs + xm) {
683c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[ifast * dof + j] += ctx->flux[j] / hxf;
684c4762a1bSJed Brown         ifast++;
685c4762a1bSJed Brown       }
686c4762a1bSJed Brown     }
687c4762a1bSJed Brown     if (i == fs) { /* interface between the fast region and the slow region */
688c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
689c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hxf / 2;
690c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hxs / 2;
691c4762a1bSJed Brown       }
6929566063dSJacob Faibussowitsch       PetscCall((*ctx->physics2.riemann2)(ctx->physics2.user, dof, uL, uR, ctx->flux, &maxspeed));
693c4762a1bSJed Brown       if (i > xs) {
694c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(ifast - 1) * dof + j] -= ctx->flux[j] / hxf;
695c4762a1bSJed Brown       }
696c4762a1bSJed Brown     }
697c4762a1bSJed Brown   }
6989566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
6999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
7009566063dSJacob Faibussowitsch   PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
7019566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &Xloc));
702c4762a1bSJed Brown   PetscFunctionReturn(0);
703c4762a1bSJed Brown }
704c4762a1bSJed Brown 
7059371c9d4SSatish Balay int main(int argc, char *argv[]) {
706c4762a1bSJed Brown   char              lname[256] = "mc", physname[256] = "advect", final_fname[256] = "solution.m";
707c4762a1bSJed Brown   PetscFunctionList limiters = 0, physics = 0;
708c4762a1bSJed Brown   MPI_Comm          comm;
709c4762a1bSJed Brown   TS                ts;
710c4762a1bSJed Brown   DM                da;
711c4762a1bSJed Brown   Vec               X, X0, R;
712c4762a1bSJed Brown   FVCtx             ctx;
713c4762a1bSJed 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;
714c4762a1bSJed Brown   PetscBool         view_final = PETSC_FALSE;
715c4762a1bSJed Brown   PetscReal         ptime;
716c4762a1bSJed Brown 
717327415f7SBarry Smith   PetscFunctionBeginUser;
7189566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, 0, help));
719c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
7209566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&ctx, sizeof(ctx)));
721c4762a1bSJed Brown 
722c4762a1bSJed Brown   /* Register limiters to be available on the command line */
7239566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "upwind", Limit2_Upwind));
7249566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "lax-wendroff", Limit2_LaxWendroff));
7259566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "beam-warming", Limit2_BeamWarming));
7269566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "fromm", Limit2_Fromm));
7279566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "minmod", Limit2_Minmod));
7289566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "superbee", Limit2_Superbee));
7299566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "mc", Limit2_MC));
7309566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "koren3", Limit2_Koren3));
731c4762a1bSJed Brown 
732c4762a1bSJed Brown   /* Register physical models to be available on the command line */
7339566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&physics, "advect", PhysicsCreate_Advect));
734c4762a1bSJed Brown 
735c4762a1bSJed Brown   ctx.comm   = comm;
736c4762a1bSJed Brown   ctx.cfl    = 0.9;
737c4762a1bSJed Brown   ctx.bctype = FVBC_PERIODIC;
738c4762a1bSJed Brown   ctx.xmin   = -1.0;
739c4762a1bSJed Brown   ctx.xmax   = 1.0;
740d0609cedSBarry Smith   PetscOptionsBegin(comm, NULL, "Finite Volume solver options", "");
7419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-xmin", "X min", "", ctx.xmin, &ctx.xmin, NULL));
7429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-xmax", "X max", "", ctx.xmax, &ctx.xmax, NULL));
7439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-limit", "Name of flux imiter to use", "", limiters, lname, lname, sizeof(lname), NULL));
7449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-draw", "Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)", "", draw, &draw, NULL));
7459566063dSJacob 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));
7469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-initial", "Initial condition (depends on the physics)", "", ctx.initial, &ctx.initial, NULL));
7479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-exact", "Compare errors with exact solution", "", ctx.exact, &ctx.exact, NULL));
7489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-simulation", "Compare errors with reference solution", "", ctx.simulation, &ctx.simulation, NULL));
7499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-cfl", "CFL number to time step at", "", ctx.cfl, &ctx.cfl, NULL));
7509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-bc_type", "Boundary condition", "", FVBCTypes, (PetscEnum)ctx.bctype, (PetscEnum *)&ctx.bctype, NULL));
7519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-hratio", "Spacing ratio", "", ctx.hratio, &ctx.hratio, NULL));
752d0609cedSBarry Smith   PetscOptionsEnd();
753c4762a1bSJed Brown 
754c4762a1bSJed Brown   /* Choose the limiter from the list of registered limiters */
7559566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(limiters, lname, &ctx.limit2));
7563c633725SBarry Smith   PetscCheck(ctx.limit2, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Limiter '%s' not found", lname);
757c4762a1bSJed Brown 
758c4762a1bSJed Brown   /* Choose the physics from the list of registered models */
759c4762a1bSJed Brown   {
760c4762a1bSJed Brown     PetscErrorCode (*r)(FVCtx *);
7619566063dSJacob Faibussowitsch     PetscCall(PetscFunctionListFind(physics, physname, &r));
7623c633725SBarry Smith     PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Physics '%s' not found", physname);
763c4762a1bSJed Brown     /* Create the physics, will set the number of fields and their names */
7649566063dSJacob Faibussowitsch     PetscCall((*r)(&ctx));
765c4762a1bSJed Brown   }
766c4762a1bSJed Brown 
767c4762a1bSJed Brown   /* Create a DMDA to manage the parallel grid */
7689566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(comm, DM_BOUNDARY_PERIODIC, 50, ctx.physics2.dof, 2, NULL, &da));
7699566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
7709566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
771c4762a1bSJed Brown   /* Inform the DMDA of the field names provided by the physics. */
772c4762a1bSJed Brown   /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */
773*48a46eb9SPierre Jolivet   for (i = 0; i < ctx.physics2.dof; i++) PetscCall(DMDASetFieldName(da, i, ctx.physics2.fieldname[i]));
7749566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
7759566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
776c4762a1bSJed Brown 
777c4762a1bSJed Brown   /* Set coordinates of cell centers */
7789566063dSJacob 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));
779c4762a1bSJed Brown 
780c4762a1bSJed Brown   /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */
7819566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(dof * dof, &ctx.R, dof * dof, &ctx.Rinv, 2 * dof, &ctx.cjmpLR, 1 * dof, &ctx.cslope));
7829566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(2 * dof, &ctx.uLR, dof, &ctx.flux, dof, &ctx.speeds));
783c4762a1bSJed Brown 
784c4762a1bSJed Brown   /* Create a vector to store the solution and to save the initial state */
7859566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &X));
7869566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X, &X0));
7879566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X, &R));
788c4762a1bSJed Brown 
789c4762a1bSJed Brown   /* create index for slow parts and fast parts,
790c4762a1bSJed Brown      count_slow + count_fast = Mx, counts_slow*hs = 0.5, counts_fast*hf = 0.5 */
791c4762a1bSJed Brown   count_slow = Mx / (1.0 + ctx.hratio / 3.0);
79208401ef6SPierre 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");
793c4762a1bSJed Brown   count_fast = Mx - count_slow;
794c4762a1bSJed Brown   ctx.sf     = count_slow / 2;
795c4762a1bSJed Brown   ctx.fs     = ctx.sf + count_fast;
7969566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(xm * dof, &index_slow));
7979566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(xm * dof, &index_fast));
7989566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(6 * dof, &index_slowbuffer));
799c4762a1bSJed Brown   if (((AdvectCtx *)ctx.physics2.user)->a > 0) {
800c4762a1bSJed Brown     ctx.lsbwidth = 2;
801c4762a1bSJed Brown     ctx.rsbwidth = 4;
802c4762a1bSJed Brown   } else {
803c4762a1bSJed Brown     ctx.lsbwidth = 4;
804c4762a1bSJed Brown     ctx.rsbwidth = 2;
805c4762a1bSJed Brown   }
806c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
807c4762a1bSJed Brown     if (i < ctx.sf - ctx.lsbwidth || i > ctx.fs + ctx.rsbwidth - 1)
808c4762a1bSJed Brown       for (k = 0; k < dof; k++) index_slow[islow++] = i * dof + k;
809c4762a1bSJed Brown     else if ((i >= ctx.sf - ctx.lsbwidth && i < ctx.sf) || (i > ctx.fs - 1 && i <= ctx.fs + ctx.rsbwidth - 1))
810c4762a1bSJed Brown       for (k = 0; k < dof; k++) index_slowbuffer[islowbuffer++] = i * dof + k;
811c4762a1bSJed Brown     else
812c4762a1bSJed Brown       for (k = 0; k < dof; k++) index_fast[ifast++] = i * dof + k;
813c4762a1bSJed Brown   }
8149566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islow, index_slow, PETSC_COPY_VALUES, &ctx.iss));
8159566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, ifast, index_fast, PETSC_COPY_VALUES, &ctx.isf));
8169566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islowbuffer, index_slowbuffer, PETSC_COPY_VALUES, &ctx.issb));
817c4762a1bSJed Brown 
818c4762a1bSJed Brown   /* Create a time-stepping object */
8199566063dSJacob Faibussowitsch   PetscCall(TSCreate(comm, &ts));
8209566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, da));
8219566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, R, FVRHSFunction_2WaySplit, &ctx));
8229566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetIS(ts, "slow", ctx.iss));
8239566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetIS(ts, "slowbuffer", ctx.issb));
8249566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetIS(ts, "fast", ctx.isf));
8259566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetRHSFunction(ts, "slow", NULL, FVRHSFunctionslow_2WaySplit, &ctx));
8269566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetRHSFunction(ts, "fast", NULL, FVRHSFunctionfast_2WaySplit, &ctx));
8279566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetRHSFunction(ts, "slowbuffer", NULL, FVRHSFunctionslowbuffer_2WaySplit, &ctx));
828c4762a1bSJed Brown 
8299566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSSSP));
8309566063dSJacob Faibussowitsch   /*PetscCall(TSSetType(ts,TSMPRK));*/
8319566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 10));
8329566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
833c4762a1bSJed Brown 
834c4762a1bSJed Brown   /* Compute initial conditions and starting time step */
8359566063dSJacob Faibussowitsch   PetscCall(FVSample_2WaySplit(&ctx, da, 0, X0));
8369566063dSJacob Faibussowitsch   PetscCall(FVRHSFunction_2WaySplit(ts, 0, X0, X, (void *)&ctx)); /* Initial function evaluation, only used to determine max speed */
8379566063dSJacob Faibussowitsch   PetscCall(VecCopy(X0, X));                                      /* The function value was not used so we set X=X0 again */
8389566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, ctx.cfl / ctx.cfl_idt));
8399566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts)); /* Take runtime options */
8409566063dSJacob Faibussowitsch   PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD));
841c4762a1bSJed Brown   {
842c4762a1bSJed Brown     PetscInt           steps;
843c4762a1bSJed Brown     PetscScalar        mass_initial, mass_final, mass_difference, mass_differenceg;
844c4762a1bSJed Brown     const PetscScalar *ptr_X, *ptr_X0;
845c4762a1bSJed Brown     const PetscReal    hs = (ctx.xmax - ctx.xmin) * 3.0 / 4.0 / count_slow;
846c4762a1bSJed Brown     const PetscReal    hf = (ctx.xmax - ctx.xmin) / 4.0 / count_fast;
847c4762a1bSJed Brown 
8489566063dSJacob Faibussowitsch     PetscCall(TSSolve(ts, X));
8499566063dSJacob Faibussowitsch     PetscCall(TSGetSolveTime(ts, &ptime));
8509566063dSJacob Faibussowitsch     PetscCall(TSGetStepNumber(ts, &steps));
851c4762a1bSJed Brown     /* calculate the total mass at initial time and final time */
852c4762a1bSJed Brown     mass_initial = 0.0;
853c4762a1bSJed Brown     mass_final   = 0.0;
8549566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArrayRead(da, X0, (void *)&ptr_X0));
8559566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArrayRead(da, X, (void *)&ptr_X));
856c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
857c4762a1bSJed Brown       if (i < ctx.sf || i > ctx.fs - 1) {
858c4762a1bSJed Brown         for (k = 0; k < dof; k++) {
859c4762a1bSJed Brown           mass_initial = mass_initial + hs * ptr_X0[i * dof + k];
860c4762a1bSJed Brown           mass_final   = mass_final + hs * ptr_X[i * dof + k];
861c4762a1bSJed Brown         }
862c4762a1bSJed Brown       } else {
863c4762a1bSJed Brown         for (k = 0; k < dof; k++) {
864c4762a1bSJed Brown           mass_initial = mass_initial + hf * ptr_X0[i * dof + k];
865c4762a1bSJed Brown           mass_final   = mass_final + hf * ptr_X[i * dof + k];
866c4762a1bSJed Brown         }
867c4762a1bSJed Brown       }
868c4762a1bSJed Brown     }
8699566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArrayRead(da, X0, (void *)&ptr_X0));
8709566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArrayRead(da, X, (void *)&ptr_X));
871c4762a1bSJed Brown     mass_difference = mass_final - mass_initial;
8729566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allreduce(&mass_difference, &mass_differenceg, 1, MPIU_SCALAR, MPIU_SUM, comm));
8739566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Mass difference %g\n", (double)mass_differenceg));
87463a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Final time %g, steps %" PetscInt_FMT "\n", (double)ptime, steps));
8759566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Maximum allowable stepsize according to CFL %g\n", (double)(1.0 / ctx.cfl_idt)));
876c4762a1bSJed Brown     if (ctx.exact) {
877c4762a1bSJed Brown       PetscReal nrm1 = 0;
8789566063dSJacob Faibussowitsch       PetscCall(SolutionErrorNorms_2WaySplit(&ctx, da, ptime, X, &nrm1));
8799566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %g\n", (double)nrm1));
880c4762a1bSJed Brown     }
881c4762a1bSJed Brown     if (ctx.simulation) {
882c4762a1bSJed Brown       PetscReal          nrm1 = 0;
883c4762a1bSJed Brown       PetscViewer        fd;
884c4762a1bSJed Brown       char               filename[PETSC_MAX_PATH_LEN] = "binaryoutput";
885c4762a1bSJed Brown       Vec                XR;
886c4762a1bSJed Brown       PetscBool          flg;
887c4762a1bSJed Brown       const PetscScalar *ptr_XR;
8889566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetString(NULL, NULL, "-f", filename, sizeof(filename), &flg));
8893c633725SBarry Smith       PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate binary file with the -f option");
8909566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, filename, FILE_MODE_READ, &fd));
8919566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(X0, &XR));
8929566063dSJacob Faibussowitsch       PetscCall(VecLoad(XR, fd));
8939566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&fd));
8949566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(X, &ptr_X));
8959566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(XR, &ptr_XR));
896c4762a1bSJed Brown       for (i = xs; i < xs + xm; i++) {
897c4762a1bSJed Brown         if (i < ctx.sf || i > ctx.fs - 1)
898c4762a1bSJed Brown           for (k = 0; k < dof; k++) nrm1 = nrm1 + hs * PetscAbs(ptr_X[i * dof + k] - ptr_XR[i * dof + k]);
899c4762a1bSJed Brown         else
900c4762a1bSJed Brown           for (k = 0; k < dof; k++) nrm1 = nrm1 + hf * PetscAbs(ptr_X[i * dof + k] - ptr_XR[i * dof + k]);
901c4762a1bSJed Brown       }
9029566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(X, &ptr_X));
9039566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(XR, &ptr_XR));
9049566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %g\n", (double)nrm1));
9059566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&XR));
906c4762a1bSJed Brown     }
907c4762a1bSJed Brown   }
908c4762a1bSJed Brown 
9099566063dSJacob Faibussowitsch   PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD));
9109566063dSJacob Faibussowitsch   if (draw & 0x1) PetscCall(VecView(X0, PETSC_VIEWER_DRAW_WORLD));
9119566063dSJacob Faibussowitsch   if (draw & 0x2) PetscCall(VecView(X, PETSC_VIEWER_DRAW_WORLD));
912c4762a1bSJed Brown   if (draw & 0x4) {
913c4762a1bSJed Brown     Vec Y;
9149566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(X, &Y));
9159566063dSJacob Faibussowitsch     PetscCall(FVSample_2WaySplit(&ctx, da, ptime, Y));
9169566063dSJacob Faibussowitsch     PetscCall(VecAYPX(Y, -1, X));
9179566063dSJacob Faibussowitsch     PetscCall(VecView(Y, PETSC_VIEWER_DRAW_WORLD));
9189566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&Y));
919c4762a1bSJed Brown   }
920c4762a1bSJed Brown 
921c4762a1bSJed Brown   if (view_final) {
922c4762a1bSJed Brown     PetscViewer viewer;
9239566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, final_fname, &viewer));
9249566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB));
9259566063dSJacob Faibussowitsch     PetscCall(VecView(X, viewer));
9269566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(viewer));
9279566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
928c4762a1bSJed Brown   }
929c4762a1bSJed Brown 
930c4762a1bSJed Brown   /* Clean up */
9319566063dSJacob Faibussowitsch   PetscCall((*ctx.physics2.destroy)(ctx.physics2.user));
9329566063dSJacob Faibussowitsch   for (i = 0; i < ctx.physics2.dof; i++) PetscCall(PetscFree(ctx.physics2.fieldname[i]));
9339566063dSJacob Faibussowitsch   PetscCall(PetscFree4(ctx.R, ctx.Rinv, ctx.cjmpLR, ctx.cslope));
9349566063dSJacob Faibussowitsch   PetscCall(PetscFree3(ctx.uLR, ctx.flux, ctx.speeds));
9359566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
9369566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X0));
9379566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&R));
9389566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
9399566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
9409566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ctx.iss));
9419566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ctx.isf));
9429566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ctx.issb));
9439566063dSJacob Faibussowitsch   PetscCall(PetscFree(index_slow));
9449566063dSJacob Faibussowitsch   PetscCall(PetscFree(index_fast));
9459566063dSJacob Faibussowitsch   PetscCall(PetscFree(index_slowbuffer));
9469566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&limiters));
9479566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&physics));
9489566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
949b122ec5aSJacob Faibussowitsch   return 0;
950c4762a1bSJed Brown }
951c4762a1bSJed Brown 
952c4762a1bSJed Brown /*TEST
953c4762a1bSJed Brown 
954c4762a1bSJed Brown     build:
955f56ea12dSJed Brown       requires: !complex
956c4762a1bSJed Brown       depends: finitevolume1d.c
957c4762a1bSJed Brown 
958c4762a1bSJed Brown     test:
959c4762a1bSJed Brown       suffix: 1
960c4762a1bSJed 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
961c4762a1bSJed Brown 
962c4762a1bSJed Brown     test:
963c4762a1bSJed Brown       suffix: 2
964c4762a1bSJed 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
965c4762a1bSJed Brown       output_file: output/ex6_1.out
966c4762a1bSJed Brown 
967c4762a1bSJed Brown TEST*/
968