xref: /petsc/src/ts/tutorials/multirate/ex7.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown static const char help[] = "1D periodic Finite Volume solver by a particular slope limiter with semidiscrete time stepping.\n"
2c4762a1bSJed Brown                            "  advection   - Constant coefficient scalar advection\n"
3c4762a1bSJed Brown                            "                u_t       + (a*u)_x               = 0\n"
4c4762a1bSJed Brown                            "  for this toy problem, we choose different meshsizes for different sub-domains, say\n"
5c4762a1bSJed Brown                            "                hxs  = (xmax - xmin)/2.0*(hratio+1.0)/Mx, \n"
6c4762a1bSJed Brown                            "                hxf  = (xmax - xmin)/2.0*(1.0+1.0/hratio)/Mx, \n"
76aad120cSJose E. Roman                            "  with x belongs to (xmin,xmax), the number of total mesh points is Mx and the ratio between the meshsize of coarse\n\n"
8c4762a1bSJed Brown                            "  grids and fine grids is hratio.\n"
9c4762a1bSJed Brown                            "  exact       - Exact Riemann solver which usually needs to perform a Newton iteration to connect\n"
10c4762a1bSJed Brown                            "                the states across shocks and rarefactions\n"
11c4762a1bSJed Brown                            "  simulation  - use reference solution which is generated by smaller time step size to be true solution,\n"
12c4762a1bSJed Brown                            "                also the reference solution should be generated by user and stored in a binary file.\n"
13c4762a1bSJed Brown                            "  characteristic - Limit the characteristic variables, this is usually preferred (default)\n"
14c4762a1bSJed Brown                            "Several initial conditions can be chosen with -initial N\n\n"
15c4762a1bSJed Brown                            "The problem size should be set with -da_grid_x M\n\n"
16c4762a1bSJed Brown                            "This script choose the slope limiter by biased second-order upwind procedure which is proposed by Van Leer in 1994\n"
17c4762a1bSJed Brown                            "                             u(x_(k+1/2),t) = u(x_k,t) + phi(x_(k+1/2),t)*(u(x_k,t)-u(x_(k-1),t))                 \n"
18c4762a1bSJed Brown                            "                     limiter phi(x_(k+1/2),t) = max(0,min(r(k+1/2),min(2,gamma(k+1/2)*r(k+1/2)+alpha(k+1/2))))    \n"
19c4762a1bSJed Brown                            "                             r(k+1/2) = (u(x_(k+1))-u(x_k))/(u(x_k)-u(x_(k-1)))                                   \n"
20c4762a1bSJed Brown                            "                             alpha(k+1/2) = (h_k*h_(k+1))/(h_(k-1)+h_k)/(h_(k-1)+h_k+h_(k+1))                     \n"
21c4762a1bSJed Brown                            "                             gamma(k+1/2) = h_k*(h_(k-1)+h_k)/(h_k+h_(k+1))/(h_(k-1)+h_k+h_(k+1))                 \n";
22c4762a1bSJed Brown 
23c4762a1bSJed Brown #include <petscts.h>
24c4762a1bSJed Brown #include <petscdm.h>
25c4762a1bSJed Brown #include <petscdmda.h>
26c4762a1bSJed Brown #include <petscdraw.h>
27c4762a1bSJed Brown #include <petscmath.h>
28c4762a1bSJed Brown 
299371c9d4SSatish Balay static inline PetscReal RangeMod(PetscReal a, PetscReal xmin, PetscReal xmax) {
309371c9d4SSatish Balay   PetscReal range = xmax - xmin;
319371c9d4SSatish Balay   return xmin + PetscFmodReal(range + PetscFmodReal(a, range), range);
329371c9d4SSatish Balay }
33c4762a1bSJed Brown 
34c4762a1bSJed Brown /* --------------------------------- Finite Volume data structures ----------------------------------- */
35c4762a1bSJed Brown 
369371c9d4SSatish Balay typedef enum {
379371c9d4SSatish Balay   FVBC_PERIODIC,
389371c9d4SSatish Balay   FVBC_OUTFLOW
399371c9d4SSatish Balay } FVBCType;
40c4762a1bSJed Brown static const char *FVBCTypes[] = {"PERIODIC", "OUTFLOW", "FVBCType", "FVBC_", 0};
41c4762a1bSJed Brown 
42c4762a1bSJed Brown typedef struct {
43c4762a1bSJed Brown   PetscErrorCode (*sample)(void *, PetscInt, FVBCType, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal *);
44c4762a1bSJed Brown   PetscErrorCode (*flux)(void *, const PetscScalar *, PetscScalar *, PetscReal *);
45c4762a1bSJed Brown   PetscErrorCode (*destroy)(void *);
46c4762a1bSJed Brown   void    *user;
47c4762a1bSJed Brown   PetscInt dof;
48c4762a1bSJed Brown   char    *fieldname[16];
49c4762a1bSJed Brown } PhysicsCtx;
50c4762a1bSJed Brown 
51c4762a1bSJed Brown typedef struct {
52c4762a1bSJed Brown   PhysicsCtx physics;
53c4762a1bSJed Brown   MPI_Comm   comm;
54c4762a1bSJed Brown   char       prefix[256];
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   /* Local work arrays */
57c4762a1bSJed Brown   PetscScalar *flux;   /* Flux across interface                                                      */
58c4762a1bSJed Brown   PetscReal   *speeds; /* Speeds of each wave                                                        */
59c4762a1bSJed Brown   PetscScalar *u;      /* value at face                                                              */
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   PetscReal cfl_idt; /* Max allowable value of 1/Delta t                                           */
62c4762a1bSJed Brown   PetscReal cfl;
63c4762a1bSJed Brown   PetscReal xmin, xmax;
64c4762a1bSJed Brown   PetscInt  initial;
65c4762a1bSJed Brown   PetscBool exact;
66c4762a1bSJed Brown   PetscBool simulation;
67c4762a1bSJed Brown   FVBCType  bctype;
68c4762a1bSJed Brown   PetscInt  hratio; /* hratio = hslow/hfast */
69c4762a1bSJed Brown   IS        isf, iss;
70c4762a1bSJed Brown   PetscInt  sf, fs; /* slow-fast and fast-slow interfaces */
71c4762a1bSJed Brown } FVCtx;
72c4762a1bSJed Brown 
73c4762a1bSJed Brown /* --------------------------------- Physics ----------------------------------- */
749371c9d4SSatish Balay static PetscErrorCode PhysicsDestroy_SimpleFree(void *vctx) {
75c4762a1bSJed Brown   PetscFunctionBeginUser;
769566063dSJacob Faibussowitsch   PetscCall(PetscFree(vctx));
77c4762a1bSJed Brown   PetscFunctionReturn(0);
78c4762a1bSJed Brown }
79c4762a1bSJed Brown 
80c4762a1bSJed Brown /* --------------------------------- Advection ----------------------------------- */
81c4762a1bSJed Brown typedef struct {
82c4762a1bSJed Brown   PetscReal a; /* advective velocity */
83c4762a1bSJed Brown } AdvectCtx;
84c4762a1bSJed Brown 
859371c9d4SSatish Balay static PetscErrorCode PhysicsFlux_Advect(void *vctx, const PetscScalar *u, PetscScalar *flux, PetscReal *maxspeed) {
86c4762a1bSJed Brown   AdvectCtx *ctx = (AdvectCtx *)vctx;
87c4762a1bSJed Brown   PetscReal  speed;
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   PetscFunctionBeginUser;
90c4762a1bSJed Brown   speed     = ctx->a;
91c4762a1bSJed Brown   flux[0]   = speed * u[0];
92c4762a1bSJed Brown   *maxspeed = speed;
93c4762a1bSJed Brown   PetscFunctionReturn(0);
94c4762a1bSJed Brown }
95c4762a1bSJed Brown 
969371c9d4SSatish Balay static PetscErrorCode PhysicsSample_Advect(void *vctx, PetscInt initial, FVBCType bctype, PetscReal xmin, PetscReal xmax, PetscReal t, PetscReal x, PetscReal *u) {
97c4762a1bSJed Brown   AdvectCtx *ctx = (AdvectCtx *)vctx;
98c4762a1bSJed Brown   PetscReal  a   = ctx->a, x0;
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   PetscFunctionBeginUser;
101c4762a1bSJed Brown   switch (bctype) {
102c4762a1bSJed Brown   case FVBC_OUTFLOW: x0 = x - a * t; break;
103c4762a1bSJed Brown   case FVBC_PERIODIC: x0 = RangeMod(x - a * t, xmin, xmax); break;
104c4762a1bSJed Brown   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown BCType");
105c4762a1bSJed Brown   }
106c4762a1bSJed Brown   switch (initial) {
107c4762a1bSJed Brown   case 0: u[0] = (x0 < 0) ? 1 : -1; break;
108c4762a1bSJed Brown   case 1: u[0] = (x0 < 0) ? -1 : 1; break;
109c4762a1bSJed Brown   case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break;
110c4762a1bSJed Brown   case 3: u[0] = PetscSinReal(2 * PETSC_PI * x0); break;
111c4762a1bSJed Brown   case 4: u[0] = PetscAbs(x0); break;
112c4762a1bSJed Brown   case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2 * PETSC_PI * x0)); break;
113c4762a1bSJed Brown   case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2 - x0 : 0)); break;
114c4762a1bSJed Brown   case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI * x0), 10.0); break;
115c4762a1bSJed Brown   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition");
116c4762a1bSJed Brown   }
117c4762a1bSJed Brown   PetscFunctionReturn(0);
118c4762a1bSJed Brown }
119c4762a1bSJed Brown 
1209371c9d4SSatish Balay static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx) {
121c4762a1bSJed Brown   AdvectCtx *user;
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   PetscFunctionBeginUser;
1249566063dSJacob Faibussowitsch   PetscCall(PetscNew(&user));
125c4762a1bSJed Brown   ctx->physics.sample  = PhysicsSample_Advect;
126c4762a1bSJed Brown   ctx->physics.flux    = PhysicsFlux_Advect;
127c4762a1bSJed Brown   ctx->physics.destroy = PhysicsDestroy_SimpleFree;
128c4762a1bSJed Brown   ctx->physics.user    = user;
129c4762a1bSJed Brown   ctx->physics.dof     = 1;
1309566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("u", &ctx->physics.fieldname[0]));
131c4762a1bSJed Brown   user->a = 1;
132d0609cedSBarry Smith   PetscOptionsBegin(ctx->comm, ctx->prefix, "Options for advection", "");
1339371c9d4SSatish Balay   { PetscCall(PetscOptionsReal("-physics_advect_a", "Speed", "", user->a, &user->a, NULL)); }
134d0609cedSBarry Smith   PetscOptionsEnd();
135c4762a1bSJed Brown   PetscFunctionReturn(0);
136c4762a1bSJed Brown }
137c4762a1bSJed Brown 
138c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver ----------------------------------- */
139c4762a1bSJed Brown 
1409371c9d4SSatish Balay static PetscErrorCode FVRHSFunction(TS ts, PetscReal time, Vec X, Vec F, void *vctx) {
141c4762a1bSJed Brown   FVCtx       *ctx = (FVCtx *)vctx;
142c4762a1bSJed Brown   PetscInt     i, j, Mx, dof, xs, xm, sf = ctx->sf, fs = ctx->fs;
143c4762a1bSJed Brown   PetscReal    hf, hs, cfl_idt = 0;
144c4762a1bSJed Brown   PetscScalar *x, *f, *r, *min, *alpha, *gamma;
145c4762a1bSJed Brown   Vec          Xloc;
146c4762a1bSJed Brown   DM           da;
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   PetscFunctionBeginUser;
1499566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
1509566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &Xloc));                                 /* Xloc contains ghost points                                     */
1519566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); /* Mx is the number of center points                              */
152c4762a1bSJed Brown   hs = (ctx->xmax - ctx->xmin) / 2.0 * (ctx->hratio + 1.0) / Mx;
153c4762a1bSJed Brown   hf = (ctx->xmax - ctx->xmin) / 2.0 * (1.0 + 1.0 / ctx->hratio) / Mx;
1549566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); /* X is solution vector which does not contain ghost points       */
1559566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
1569566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F)); /* F is the right hand side function corresponds to center points */
1579566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Xloc, &x));
1589566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, F, &f));
1599566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
1609566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(dof, &r, dof, &min, dof, &alpha, dof, &gamma));
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
163c4762a1bSJed Brown     for (i = xs - 2; i < 0; i++) {
164c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
165c4762a1bSJed Brown     }
166c4762a1bSJed Brown     for (i = Mx; i < xs + xm + 2; i++) {
167c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
168c4762a1bSJed Brown     }
169c4762a1bSJed Brown   }
170c4762a1bSJed Brown 
171c4762a1bSJed Brown   for (i = xs; i < xs + xm + 1; i++) {
172c4762a1bSJed Brown     PetscReal    maxspeed;
173c4762a1bSJed Brown     PetscScalar *u;
174c4762a1bSJed Brown     if (i < sf || i > fs + 1) {
175c4762a1bSJed Brown       u        = &ctx->u[0];
176c4762a1bSJed Brown       alpha[0] = 1.0 / 6.0;
177c4762a1bSJed Brown       gamma[0] = 1.0 / 3.0;
178c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
179c4762a1bSJed Brown         r[j]   = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
180c4762a1bSJed Brown         min[j] = PetscMin(r[j], 2.0);
181c4762a1bSJed Brown         u[j]   = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
182c4762a1bSJed Brown       }
1839566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed));
184c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hs));
185c4762a1bSJed Brown       if (i > xs) {
186c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hs;
187c4762a1bSJed Brown       }
188c4762a1bSJed Brown       if (i < xs + xm) {
189c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hs;
190c4762a1bSJed Brown       }
191c4762a1bSJed Brown     } else if (i == sf) {
192c4762a1bSJed Brown       u        = &ctx->u[0];
193c4762a1bSJed Brown       alpha[0] = hs * hf / (hs + hs) / (hs + hs + hf);
194c4762a1bSJed Brown       gamma[0] = hs * (hs + hs) / (hs + hf) / (hs + hs + hf);
195c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
196c4762a1bSJed Brown         r[j]   = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
197c4762a1bSJed Brown         min[j] = PetscMin(r[j], 2.0);
198c4762a1bSJed Brown         u[j]   = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
199c4762a1bSJed Brown       }
2009566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed));
201c4762a1bSJed Brown       if (i > xs) {
202c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hs;
203c4762a1bSJed Brown       }
204c4762a1bSJed Brown       if (i < xs + xm) {
205c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hf;
206c4762a1bSJed Brown       }
207c4762a1bSJed Brown     } else if (i == sf + 1) {
208c4762a1bSJed Brown       u        = &ctx->u[0];
209c4762a1bSJed Brown       alpha[0] = hf * hf / (hs + hf) / (hs + hf + hf);
210c4762a1bSJed Brown       gamma[0] = hf * (hs + hf) / (hf + hf) / (hs + hf + hf);
211c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
212c4762a1bSJed Brown         r[j]   = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
213c4762a1bSJed Brown         min[j] = PetscMin(r[j], 2.0);
214c4762a1bSJed Brown         u[j]   = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
215c4762a1bSJed Brown       }
2169566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed));
217c4762a1bSJed Brown       if (i > xs) {
218c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hf;
219c4762a1bSJed Brown       }
220c4762a1bSJed Brown       if (i < xs + xm) {
221c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hf;
222c4762a1bSJed Brown       }
223c4762a1bSJed Brown     } else if (i > sf + 1 && i < fs) {
224c4762a1bSJed Brown       u        = &ctx->u[0];
225c4762a1bSJed Brown       alpha[0] = 1.0 / 6.0;
226c4762a1bSJed Brown       gamma[0] = 1.0 / 3.0;
227c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
228c4762a1bSJed Brown         r[j]   = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
229c4762a1bSJed Brown         min[j] = PetscMin(r[j], 2.0);
230c4762a1bSJed Brown         u[j]   = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
231c4762a1bSJed Brown       }
2329566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed));
233c4762a1bSJed Brown       if (i > xs) {
234c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hf;
235c4762a1bSJed Brown       }
236c4762a1bSJed Brown       if (i < xs + xm) {
237c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hf;
238c4762a1bSJed Brown       }
239c4762a1bSJed Brown     } else if (i == fs) {
240c4762a1bSJed Brown       u        = &ctx->u[0];
241c4762a1bSJed Brown       alpha[0] = hf * hs / (hf + hf) / (hf + hf + hs);
242c4762a1bSJed Brown       gamma[0] = hf * (hf + hf) / (hf + hs) / (hf + hf + hs);
243c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
244c4762a1bSJed Brown         r[j]   = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
245c4762a1bSJed Brown         min[j] = PetscMin(r[j], 2.0);
246c4762a1bSJed Brown         u[j]   = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
247c4762a1bSJed Brown       }
2489566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed));
249c4762a1bSJed Brown       if (i > xs) {
250c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hf;
251c4762a1bSJed Brown       }
252c4762a1bSJed Brown       if (i < xs + xm) {
253c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hs;
254c4762a1bSJed Brown       }
255c4762a1bSJed Brown     } else if (i == fs + 1) {
256c4762a1bSJed Brown       u        = &ctx->u[0];
257c4762a1bSJed Brown       alpha[0] = hs * hs / (hf + hs) / (hf + hs + hs);
258c4762a1bSJed Brown       gamma[0] = hs * (hf + hs) / (hs + hs) / (hf + hs + hs);
259c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
260c4762a1bSJed Brown         r[j]   = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
261c4762a1bSJed Brown         min[j] = PetscMin(r[j], 2.0);
262c4762a1bSJed Brown         u[j]   = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
263c4762a1bSJed Brown       }
2649566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed));
265c4762a1bSJed Brown       if (i > xs) {
266c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hs;
267c4762a1bSJed Brown       }
268c4762a1bSJed Brown       if (i < xs + xm) {
269c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hs;
270c4762a1bSJed Brown       }
271c4762a1bSJed Brown     }
272c4762a1bSJed Brown   }
2739566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
2749566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, F, &f));
2759566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &Xloc));
2769566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&cfl_idt, &ctx->cfl_idt, 1, MPIU_SCALAR, MPIU_MAX, PetscObjectComm((PetscObject)da)));
277c4762a1bSJed Brown   if (0) {
278c2a16db8SHong Zhang     /* We need a way to inform the TS of a CFL constraint, this is a debugging fragment */
279c4762a1bSJed Brown     PetscReal dt, tnow;
2809566063dSJacob Faibussowitsch     PetscCall(TSGetTimeStep(ts, &dt));
2819566063dSJacob Faibussowitsch     PetscCall(TSGetTime(ts, &tnow));
282*48a46eb9SPierre Jolivet     if (dt > 0.5 / ctx->cfl_idt) PetscCall(PetscPrintf(ctx->comm, "Stability constraint exceeded at t=%g, dt %g > %g\n", (double)tnow, (double)dt, (double)(0.5 / ctx->cfl_idt)));
283c4762a1bSJed Brown   }
2849566063dSJacob Faibussowitsch   PetscCall(PetscFree4(r, min, alpha, gamma));
285c4762a1bSJed Brown   PetscFunctionReturn(0);
286c4762a1bSJed Brown }
287c4762a1bSJed Brown 
2889371c9d4SSatish Balay static PetscErrorCode FVRHSFunctionslow(TS ts, PetscReal time, Vec X, Vec F, void *vctx) {
289c4762a1bSJed Brown   FVCtx       *ctx = (FVCtx *)vctx;
290c4762a1bSJed Brown   PetscInt     i, j, Mx, dof, xs, xm, islow = 0, sf = ctx->sf, fs = ctx->fs;
291c4762a1bSJed Brown   PetscReal    hf, hs;
292c4762a1bSJed Brown   PetscScalar *x, *f, *r, *min, *alpha, *gamma;
293c4762a1bSJed Brown   Vec          Xloc;
294c4762a1bSJed Brown   DM           da;
295c4762a1bSJed Brown 
296c4762a1bSJed Brown   PetscFunctionBeginUser;
2979566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
2989566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &Xloc));                                 /* Xloc contains ghost points                                     */
2999566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); /* Mx is the number of center points                              */
300c4762a1bSJed Brown   hs = (ctx->xmax - ctx->xmin) / 2.0 * (ctx->hratio + 1.0) / Mx;
301c4762a1bSJed Brown   hf = (ctx->xmax - ctx->xmin) / 2.0 * (1.0 + 1.0 / ctx->hratio) / Mx;
3029566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); /* X is solution vector which does not contain ghost points       */
3039566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
3049566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F)); /* F is the right hand side function corresponds to center points */
3059566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Xloc, &x));
3069566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
3079566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
3089566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(dof, &r, dof, &min, dof, &alpha, dof, &gamma));
309c4762a1bSJed Brown 
310c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
311c4762a1bSJed Brown     for (i = xs - 2; i < 0; i++) {
312c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
313c4762a1bSJed Brown     }
314c4762a1bSJed Brown     for (i = Mx; i < xs + xm + 2; i++) {
315c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
316c4762a1bSJed Brown     }
317c4762a1bSJed Brown   }
318c4762a1bSJed Brown 
319c4762a1bSJed Brown   for (i = xs; i < xs + xm + 1; i++) {
320c4762a1bSJed Brown     PetscReal    maxspeed;
321c4762a1bSJed Brown     PetscScalar *u;
322c4762a1bSJed Brown     if (i < sf) {
323c4762a1bSJed Brown       u        = &ctx->u[0];
324c4762a1bSJed Brown       alpha[0] = 1.0 / 6.0;
325c4762a1bSJed Brown       gamma[0] = 1.0 / 3.0;
326c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
327c4762a1bSJed Brown         r[j]   = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
328c4762a1bSJed Brown         min[j] = PetscMin(r[j], 2.0);
329c4762a1bSJed Brown         u[j]   = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
330c4762a1bSJed Brown       }
3319566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed));
332c4762a1bSJed Brown       if (i > xs) {
333c2a16db8SHong Zhang         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hs;
334c4762a1bSJed Brown       }
335c4762a1bSJed Brown       if (i < xs + xm) {
336c2a16db8SHong Zhang         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hs;
337c4762a1bSJed Brown         islow++;
338c4762a1bSJed Brown       }
339c4762a1bSJed Brown     } else if (i == sf) {
340c4762a1bSJed Brown       u        = &ctx->u[0];
341c4762a1bSJed Brown       alpha[0] = hs * hf / (hs + hs) / (hs + hs + hf);
342c4762a1bSJed Brown       gamma[0] = hs * (hs + hs) / (hs + hf) / (hs + hs + hf);
343c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
344c4762a1bSJed Brown         r[j]   = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
345c4762a1bSJed Brown         min[j] = PetscMin(r[j], 2.0);
346c4762a1bSJed Brown         u[j]   = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
347c4762a1bSJed Brown       }
3489566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed));
349c2a16db8SHong Zhang       if (i > xs) {
350c2a16db8SHong Zhang         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hs;
351c4762a1bSJed Brown       }
352c4762a1bSJed Brown     } else if (i == fs) {
353c4762a1bSJed Brown       u        = &ctx->u[0];
354c4762a1bSJed Brown       alpha[0] = hf * hs / (hf + hf) / (hf + hf + hs);
355c4762a1bSJed Brown       gamma[0] = hf * (hf + hf) / (hf + hs) / (hf + hf + hs);
356c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
357c4762a1bSJed Brown         r[j]   = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
358c4762a1bSJed Brown         min[j] = PetscMin(r[j], 2.0);
359c4762a1bSJed Brown         u[j]   = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
360c4762a1bSJed Brown       }
3619566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed));
362c4762a1bSJed Brown       if (i < xs + xm) {
363c2a16db8SHong Zhang         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hs;
364c4762a1bSJed Brown         islow++;
365c4762a1bSJed Brown       }
366c4762a1bSJed Brown     } else if (i == fs + 1) {
367c4762a1bSJed Brown       u        = &ctx->u[0];
368c4762a1bSJed Brown       alpha[0] = hs * hs / (hf + hs) / (hf + hs + hs);
369c4762a1bSJed Brown       gamma[0] = hs * (hf + hs) / (hs + hs) / (hf + hs + hs);
370c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
371c4762a1bSJed Brown         r[j]   = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
372c4762a1bSJed Brown         min[j] = PetscMin(r[j], 2.0);
373c4762a1bSJed Brown         u[j]   = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
374c4762a1bSJed Brown       }
3759566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed));
376c4762a1bSJed Brown       if (i > xs) {
377c2a16db8SHong Zhang         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hs;
378c4762a1bSJed Brown       }
379c4762a1bSJed Brown       if (i < xs + xm) {
380c2a16db8SHong Zhang         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hs;
381c2a16db8SHong Zhang         islow++;
382c4762a1bSJed Brown       }
383c4762a1bSJed Brown     } else if (i > fs + 1) {
384c4762a1bSJed Brown       u        = &ctx->u[0];
385c4762a1bSJed Brown       alpha[0] = 1.0 / 6.0;
386c4762a1bSJed Brown       gamma[0] = 1.0 / 3.0;
387c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
388c4762a1bSJed Brown         r[j]   = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
389c4762a1bSJed Brown         min[j] = PetscMin(r[j], 2.0);
390c4762a1bSJed Brown         u[j]   = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
391c4762a1bSJed Brown       }
3929566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed));
393c4762a1bSJed Brown       if (i > xs) {
394c2a16db8SHong Zhang         for (j = 0; j < dof; j++) f[(islow - 1) * dof + j] -= ctx->flux[j] / hs;
395c4762a1bSJed Brown       }
396c4762a1bSJed Brown       if (i < xs + xm) {
397c2a16db8SHong Zhang         for (j = 0; j < dof; j++) f[islow * dof + j] += ctx->flux[j] / hs;
398c2a16db8SHong Zhang         islow++;
399c4762a1bSJed Brown       }
400c4762a1bSJed Brown     }
401c4762a1bSJed Brown   }
4029566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
4039566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
4049566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &Xloc));
4059566063dSJacob Faibussowitsch   PetscCall(PetscFree4(r, min, alpha, gamma));
406c4762a1bSJed Brown   PetscFunctionReturn(0);
407c4762a1bSJed Brown }
408c4762a1bSJed Brown 
4099371c9d4SSatish Balay static PetscErrorCode FVRHSFunctionfast(TS ts, PetscReal time, Vec X, Vec F, void *vctx) {
410c4762a1bSJed Brown   FVCtx       *ctx = (FVCtx *)vctx;
411c4762a1bSJed Brown   PetscInt     i, j, Mx, dof, xs, xm, ifast = 0, sf = ctx->sf, fs = ctx->fs;
412c4762a1bSJed Brown   PetscReal    hf, hs;
413c4762a1bSJed Brown   PetscScalar *x, *f, *r, *min, *alpha, *gamma;
414c4762a1bSJed Brown   Vec          Xloc;
415c4762a1bSJed Brown   DM           da;
416c4762a1bSJed Brown 
417c4762a1bSJed Brown   PetscFunctionBeginUser;
4189566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
4199566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &Xloc));                                 /* Xloc contains ghost points                                     */
4209566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); /* Mx is the number of center points                              */
421c4762a1bSJed Brown   hs = (ctx->xmax - ctx->xmin) / 2.0 * (ctx->hratio + 1.0) / Mx;
422c4762a1bSJed Brown   hf = (ctx->xmax - ctx->xmin) / 2.0 * (1.0 + 1.0 / ctx->hratio) / Mx;
4239566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); /* X is solution vector which does not contain ghost points       */
4249566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
4259566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F)); /* F is the right hand side function corresponds to center points */
4269566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Xloc, &x));
4279566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
4289566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
4299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(dof, &r, dof, &min, dof, &alpha, dof, &gamma));
430c4762a1bSJed Brown 
431c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
432c4762a1bSJed Brown     for (i = xs - 2; i < 0; i++) {
433c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
434c4762a1bSJed Brown     }
435c4762a1bSJed Brown     for (i = Mx; i < xs + xm + 2; i++) {
436c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
437c4762a1bSJed Brown     }
438c4762a1bSJed Brown   }
439c4762a1bSJed Brown 
440c4762a1bSJed Brown   for (i = xs; i < xs + xm + 1; i++) {
441c4762a1bSJed Brown     PetscReal    maxspeed;
442c4762a1bSJed Brown     PetscScalar *u;
443c4762a1bSJed Brown     if (i == sf) {
444c4762a1bSJed Brown       u        = &ctx->u[0];
445c4762a1bSJed Brown       alpha[0] = hs * hf / (hs + hs) / (hs + hs + hf);
446c4762a1bSJed Brown       gamma[0] = hs * (hs + hs) / (hs + hf) / (hs + hs + hf);
447c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
448c4762a1bSJed Brown         r[j]   = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
449c4762a1bSJed Brown         min[j] = PetscMin(r[j], 2.0);
450c4762a1bSJed Brown         u[j]   = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
451c4762a1bSJed Brown       }
4529566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed));
453c4762a1bSJed Brown       if (i < xs + xm) {
454c2a16db8SHong Zhang         for (j = 0; j < dof; j++) f[ifast * dof + j] += ctx->flux[j] / hf;
455c4762a1bSJed Brown         ifast++;
456c4762a1bSJed Brown       }
457c4762a1bSJed Brown     } else if (i == sf + 1) {
458c4762a1bSJed Brown       u        = &ctx->u[0];
459c4762a1bSJed Brown       alpha[0] = hf * hf / (hs + hf) / (hs + hf + hf);
460c4762a1bSJed Brown       gamma[0] = hf * (hs + hf) / (hf + hf) / (hs + hf + hf);
461c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
462c4762a1bSJed Brown         r[j]   = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
463c4762a1bSJed Brown         min[j] = PetscMin(r[j], 2.0);
464c4762a1bSJed Brown         u[j]   = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
465c4762a1bSJed Brown       }
4669566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed));
467c4762a1bSJed Brown       if (i > xs) {
468c2a16db8SHong Zhang         for (j = 0; j < dof; j++) f[(ifast - 1) * dof + j] -= ctx->flux[j] / hf;
469c4762a1bSJed Brown       }
470c4762a1bSJed Brown       if (i < xs + xm) {
471c2a16db8SHong Zhang         for (j = 0; j < dof; j++) f[ifast * dof + j] += ctx->flux[j] / hf;
472c4762a1bSJed Brown         ifast++;
473c4762a1bSJed Brown       }
474c4762a1bSJed Brown     } else if (i > sf + 1 && i < fs) {
475c4762a1bSJed Brown       u        = &ctx->u[0];
476c4762a1bSJed Brown       alpha[0] = 1.0 / 6.0;
477c4762a1bSJed Brown       gamma[0] = 1.0 / 3.0;
478c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
479c4762a1bSJed Brown         r[j]   = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
480c4762a1bSJed Brown         min[j] = PetscMin(r[j], 2.0);
481c4762a1bSJed Brown         u[j]   = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
482c4762a1bSJed Brown       }
4839566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed));
484c4762a1bSJed Brown       if (i > xs) {
485c2a16db8SHong Zhang         for (j = 0; j < dof; j++) f[(ifast - 1) * dof + j] -= ctx->flux[j] / hf;
486c4762a1bSJed Brown       }
487c4762a1bSJed Brown       if (i < xs + xm) {
488c2a16db8SHong Zhang         for (j = 0; j < dof; j++) f[ifast * dof + j] += ctx->flux[j] / hf;
489c4762a1bSJed Brown         ifast++;
490c4762a1bSJed Brown       }
491c4762a1bSJed Brown     } else if (i == fs) {
492c4762a1bSJed Brown       u        = &ctx->u[0];
493c4762a1bSJed Brown       alpha[0] = hf * hs / (hf + hf) / (hf + hf + hs);
494c4762a1bSJed Brown       gamma[0] = hf * (hf + hf) / (hf + hs) / (hf + hf + hs);
495c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
496c4762a1bSJed Brown         r[j]   = (x[i * dof + j] - x[(i - 1) * dof + j]) / (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
497c4762a1bSJed Brown         min[j] = PetscMin(r[j], 2.0);
498c4762a1bSJed Brown         u[j]   = x[(i - 1) * dof + j] + PetscMax(0, PetscMin(min[j], alpha[0] + gamma[0] * r[j])) * (x[(i - 1) * dof + j] - x[(i - 2) * dof + j]);
499c4762a1bSJed Brown       }
5009566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.flux)(ctx->physics.user, u, ctx->flux, &maxspeed));
501c4762a1bSJed Brown       if (i > xs) {
502c2a16db8SHong Zhang         for (j = 0; j < dof; j++) f[(ifast - 1) * dof + j] -= ctx->flux[j] / hf;
503c4762a1bSJed Brown       }
504c4762a1bSJed Brown     }
505c4762a1bSJed Brown   }
5069566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
5079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
5089566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &Xloc));
5099566063dSJacob Faibussowitsch   PetscCall(PetscFree4(r, min, alpha, gamma));
510c4762a1bSJed Brown   PetscFunctionReturn(0);
511c4762a1bSJed Brown }
512c4762a1bSJed Brown 
513c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for slow components ----------------------------------- */
514c4762a1bSJed Brown 
5159371c9d4SSatish Balay PetscErrorCode FVSample(FVCtx *ctx, DM da, PetscReal time, Vec U) {
516c4762a1bSJed Brown   PetscScalar   *u, *uj, xj, xi;
517c4762a1bSJed Brown   PetscInt       i, j, k, dof, xs, xm, Mx, count_slow, count_fast;
518c4762a1bSJed Brown   const PetscInt N = 200;
519c4762a1bSJed Brown 
520c4762a1bSJed Brown   PetscFunctionBeginUser;
5213c633725SBarry Smith   PetscCheck(ctx->physics.sample, PETSC_COMM_SELF, PETSC_ERR_SUP, "Physics has not provided a sampling function");
5229566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
5239566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
5249566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, U, &u));
5259566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dof, &uj));
526c4762a1bSJed Brown   const PetscReal hs = (ctx->xmax - ctx->xmin) / 2.0 * (ctx->hratio + 1.0) / Mx;
527c4762a1bSJed Brown   const PetscReal hf = (ctx->xmax - ctx->xmin) / 2.0 * (1.0 + 1.0 / ctx->hratio) / Mx;
528c4762a1bSJed Brown   count_slow         = Mx / (1 + ctx->hratio);
529c4762a1bSJed Brown   count_fast         = Mx - count_slow;
530c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
531c4762a1bSJed Brown     if (i * hs + 0.5 * hs < (ctx->xmax - ctx->xmin) * 0.25) {
532c4762a1bSJed Brown       xi = ctx->xmin + 0.5 * hs + i * hs;
533c4762a1bSJed Brown       /* Integrate over cell i using trapezoid rule with N points. */
534c4762a1bSJed Brown       for (k = 0; k < dof; k++) u[i * dof + k] = 0;
535c4762a1bSJed Brown       for (j = 0; j < N + 1; j++) {
536c4762a1bSJed Brown         xj = xi + hs * (j - N / 2) / (PetscReal)N;
5379566063dSJacob Faibussowitsch         PetscCall((*ctx->physics.sample)(ctx->physics.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
538c4762a1bSJed Brown         for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
539c4762a1bSJed Brown       }
540c4762a1bSJed Brown     } else if ((ctx->xmax - ctx->xmin) * 0.25 + (i - count_slow / 2) * hf + 0.5 * hf < (ctx->xmax - ctx->xmin) * 0.75) {
541c4762a1bSJed Brown       xi = ctx->xmin + (ctx->xmax - ctx->xmin) * 0.25 + 0.5 * hf + (i - count_slow / 2) * hf;
542c4762a1bSJed Brown       /* Integrate over cell i using trapezoid rule with N points. */
543c4762a1bSJed Brown       for (k = 0; k < dof; k++) u[i * dof + k] = 0;
544c4762a1bSJed Brown       for (j = 0; j < N + 1; j++) {
545c4762a1bSJed Brown         xj = xi + hf * (j - N / 2) / (PetscReal)N;
5469566063dSJacob Faibussowitsch         PetscCall((*ctx->physics.sample)(ctx->physics.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
547c4762a1bSJed Brown         for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
548c4762a1bSJed Brown       }
549c4762a1bSJed Brown     } else {
550c4762a1bSJed Brown       xi = ctx->xmin + (ctx->xmax - ctx->xmin) * 0.75 + 0.5 * hs + (i - count_slow / 2 - count_fast) * hs;
551c4762a1bSJed Brown       /* Integrate over cell i using trapezoid rule with N points. */
552c4762a1bSJed Brown       for (k = 0; k < dof; k++) u[i * dof + k] = 0;
553c4762a1bSJed Brown       for (j = 0; j < N + 1; j++) {
554c4762a1bSJed Brown         xj = xi + hs * (j - N / 2) / (PetscReal)N;
5559566063dSJacob Faibussowitsch         PetscCall((*ctx->physics.sample)(ctx->physics.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
556c4762a1bSJed Brown         for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
557c4762a1bSJed Brown       }
558c4762a1bSJed Brown     }
559c4762a1bSJed Brown   }
5609566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, U, &u));
5619566063dSJacob Faibussowitsch   PetscCall(PetscFree(uj));
562c4762a1bSJed Brown   PetscFunctionReturn(0);
563c4762a1bSJed Brown }
564c4762a1bSJed Brown 
5659371c9d4SSatish Balay static PetscErrorCode SolutionStatsView(DM da, Vec X, PetscViewer viewer) {
566c4762a1bSJed Brown   PetscReal          xmin, xmax;
567c4762a1bSJed Brown   PetscScalar        sum, tvsum, tvgsum;
568c4762a1bSJed Brown   const PetscScalar *x;
569c4762a1bSJed Brown   PetscInt           imin, imax, Mx, i, j, xs, xm, dof;
570c4762a1bSJed Brown   Vec                Xloc;
571c4762a1bSJed Brown   PetscBool          iascii;
572c4762a1bSJed Brown 
573c4762a1bSJed Brown   PetscFunctionBeginUser;
5749566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
575c4762a1bSJed Brown   if (iascii) {
576c4762a1bSJed Brown     /* PETSc lacks a function to compute total variation norm (difficult in multiple dimensions), we do it here */
5779566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(da, &Xloc));
5789566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
5799566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
5809566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArrayRead(da, Xloc, (void *)&x));
5819566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
5829566063dSJacob Faibussowitsch     PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
583c4762a1bSJed Brown     tvsum = 0;
584c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
585c4762a1bSJed Brown       for (j = 0; j < dof; j++) tvsum += PetscAbsScalar(x[i * dof + j] - x[(i - 1) * dof + j]);
586c4762a1bSJed Brown     }
5879566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allreduce(&tvsum, &tvgsum, 1, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)da)));
5889566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArrayRead(da, Xloc, (void *)&x));
5899566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(da, &Xloc));
590c4762a1bSJed Brown 
5919566063dSJacob Faibussowitsch     PetscCall(VecMin(X, &imin, &xmin));
5929566063dSJacob Faibussowitsch     PetscCall(VecMax(X, &imax, &xmax));
5939566063dSJacob Faibussowitsch     PetscCall(VecSum(X, &sum));
59463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Solution range [%g,%g] with minimum at %" PetscInt_FMT ", mean %g, ||x||_TV %g\n", (double)xmin, (double)xmax, imin, (double)(sum / Mx), (double)(tvgsum / Mx)));
595c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Viewer type not supported");
596c4762a1bSJed Brown   PetscFunctionReturn(0);
597c4762a1bSJed Brown }
598c4762a1bSJed Brown 
5999371c9d4SSatish Balay static PetscErrorCode SolutionErrorNorms(FVCtx *ctx, DM da, PetscReal t, Vec X, PetscReal *nrm1) {
600c4762a1bSJed Brown   Vec                Y;
601c4762a1bSJed Brown   PetscInt           i, Mx, count_slow = 0, count_fast = 0;
602c4762a1bSJed Brown   const PetscScalar *ptr_X, *ptr_Y;
603c4762a1bSJed Brown 
604c4762a1bSJed Brown   PetscFunctionBeginUser;
6059566063dSJacob Faibussowitsch   PetscCall(VecGetSize(X, &Mx));
6069566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X, &Y));
6079566063dSJacob Faibussowitsch   PetscCall(FVSample(ctx, da, t, Y));
608c4762a1bSJed Brown   const PetscReal hs = (ctx->xmax - ctx->xmin) / 2.0 * (ctx->hratio + 1.0) / Mx;
609c4762a1bSJed Brown   const PetscReal hf = (ctx->xmax - ctx->xmin) / 2.0 * (1.0 + 1.0 / ctx->hratio) / Mx;
610c4762a1bSJed Brown   count_slow         = (PetscReal)Mx / (1.0 + ctx->hratio);
611c4762a1bSJed Brown   count_fast         = Mx - count_slow;
6129566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &ptr_X));
6139566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &ptr_Y));
614c4762a1bSJed Brown   for (i = 0; i < Mx; i++) {
615c4762a1bSJed Brown     if (i < count_slow / 2 || i > count_slow / 2 + count_fast - 1) *nrm1 += hs * PetscAbs(ptr_X[i] - ptr_Y[i]);
616c4762a1bSJed Brown     else *nrm1 += hf * PetscAbs(ptr_X[i] - ptr_Y[i]);
617c4762a1bSJed Brown   }
6189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &ptr_X));
6199566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &ptr_Y));
6209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Y));
621c4762a1bSJed Brown   PetscFunctionReturn(0);
622c4762a1bSJed Brown }
623c4762a1bSJed Brown 
6249371c9d4SSatish Balay int main(int argc, char *argv[]) {
625c4762a1bSJed Brown   char              physname[256] = "advect", final_fname[256] = "solution.m";
626c4762a1bSJed Brown   PetscFunctionList physics = 0;
627c4762a1bSJed Brown   MPI_Comm          comm;
628c4762a1bSJed Brown   TS                ts;
629c4762a1bSJed Brown   DM                da;
630c4762a1bSJed Brown   Vec               X, X0, R;
631c4762a1bSJed Brown   FVCtx             ctx;
632c4762a1bSJed Brown   PetscInt          i, k, dof, xs, xm, Mx, draw = 0, count_slow, count_fast, islow = 0, ifast = 0, *index_slow, *index_fast;
633c4762a1bSJed Brown   PetscBool         view_final = PETSC_FALSE;
634c4762a1bSJed Brown   PetscReal         ptime;
635c4762a1bSJed Brown 
636327415f7SBarry Smith   PetscFunctionBeginUser;
6379566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, 0, help));
638c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
6399566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&ctx, sizeof(ctx)));
640c4762a1bSJed Brown 
641c4762a1bSJed Brown   /* Register physical models to be available on the command line */
6429566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&physics, "advect", PhysicsCreate_Advect));
643c4762a1bSJed Brown 
644c4762a1bSJed Brown   ctx.comm   = comm;
645c4762a1bSJed Brown   ctx.cfl    = 0.9;
646c4762a1bSJed Brown   ctx.bctype = FVBC_PERIODIC;
647c4762a1bSJed Brown   ctx.xmin   = -1.0;
648c4762a1bSJed Brown   ctx.xmax   = 1.0;
649d0609cedSBarry Smith   PetscOptionsBegin(comm, NULL, "Finite Volume solver options", "");
6509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-xmin", "X min", "", ctx.xmin, &ctx.xmin, NULL));
6519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-xmax", "X max", "", ctx.xmax, &ctx.xmax, NULL));
6529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-draw", "Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)", "", draw, &draw, NULL));
6539566063dSJacob 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));
6549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-initial", "Initial condition (depends on the physics)", "", ctx.initial, &ctx.initial, NULL));
6559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-exact", "Compare errors with exact solution", "", ctx.exact, &ctx.exact, NULL));
6569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-simulation", "Compare errors with reference solution", "", ctx.simulation, &ctx.simulation, NULL));
6579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-cfl", "CFL number to time step at", "", ctx.cfl, &ctx.cfl, NULL));
6589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-bc_type", "Boundary condition", "", FVBCTypes, (PetscEnum)ctx.bctype, (PetscEnum *)&ctx.bctype, NULL));
6599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-hratio", "Spacing ratio", "", ctx.hratio, &ctx.hratio, NULL));
660d0609cedSBarry Smith   PetscOptionsEnd();
661c4762a1bSJed Brown 
662c4762a1bSJed Brown   /* Choose the physics from the list of registered models */
663c4762a1bSJed Brown   {
664c4762a1bSJed Brown     PetscErrorCode (*r)(FVCtx *);
6659566063dSJacob Faibussowitsch     PetscCall(PetscFunctionListFind(physics, physname, &r));
6663c633725SBarry Smith     PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Physics '%s' not found", physname);
667c4762a1bSJed Brown     /* Create the physics, will set the number of fields and their names */
6689566063dSJacob Faibussowitsch     PetscCall((*r)(&ctx));
669c4762a1bSJed Brown   }
670c4762a1bSJed Brown 
671c4762a1bSJed Brown   /* Create a DMDA to manage the parallel grid */
6729566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(comm, DM_BOUNDARY_PERIODIC, 50, ctx.physics.dof, 2, NULL, &da));
6739566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
6749566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
675c4762a1bSJed Brown   /* Inform the DMDA of the field names provided by the physics. */
676c4762a1bSJed Brown   /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */
677*48a46eb9SPierre Jolivet   for (i = 0; i < ctx.physics.dof; i++) PetscCall(DMDASetFieldName(da, i, ctx.physics.fieldname[i]));
6789566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
6799566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
680c4762a1bSJed Brown 
681c4762a1bSJed Brown   /* Set coordinates of cell centers */
6829566063dSJacob 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));
683c4762a1bSJed Brown 
684c4762a1bSJed Brown   /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */
6859566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(dof, &ctx.u, dof, &ctx.flux, dof, &ctx.speeds));
686c4762a1bSJed Brown 
687c4762a1bSJed Brown   /* Create a vector to store the solution and to save the initial state */
6889566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &X));
6899566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X, &X0));
6909566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X, &R));
691c4762a1bSJed Brown 
692c4762a1bSJed Brown   /* create index for slow parts and fast parts*/
693c4762a1bSJed Brown   count_slow = Mx / (1 + ctx.hratio);
69408401ef6SPierre 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) is even");
695c4762a1bSJed Brown   count_fast = Mx - count_slow;
696c4762a1bSJed Brown   ctx.sf     = count_slow / 2;
697c4762a1bSJed Brown   ctx.fs     = ctx.sf + count_fast;
6989566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(xm * dof, &index_slow));
6999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(xm * dof, &index_fast));
700c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
701c4762a1bSJed Brown     if (i < count_slow / 2 || i > count_slow / 2 + count_fast - 1)
702c4762a1bSJed Brown       for (k = 0; k < dof; k++) index_slow[islow++] = i * dof + k;
703c4762a1bSJed Brown     else
704c4762a1bSJed Brown       for (k = 0; k < dof; k++) index_fast[ifast++] = i * dof + k;
705c4762a1bSJed Brown   }
7069566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islow, index_slow, PETSC_COPY_VALUES, &ctx.iss));
7079566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, ifast, index_fast, PETSC_COPY_VALUES, &ctx.isf));
708c4762a1bSJed Brown 
709c4762a1bSJed Brown   /* Create a time-stepping object */
7109566063dSJacob Faibussowitsch   PetscCall(TSCreate(comm, &ts));
7119566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, da));
7129566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, R, FVRHSFunction, &ctx));
7139566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetIS(ts, "slow", ctx.iss));
7149566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetIS(ts, "fast", ctx.isf));
7159566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetRHSFunction(ts, "slow", NULL, FVRHSFunctionslow, &ctx));
7169566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetRHSFunction(ts, "fast", NULL, FVRHSFunctionfast, &ctx));
717c4762a1bSJed Brown 
7189566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSMPRK));
7199566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 10));
7209566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
721c4762a1bSJed Brown 
722c4762a1bSJed Brown   /* Compute initial conditions and starting time step */
7239566063dSJacob Faibussowitsch   PetscCall(FVSample(&ctx, da, 0, X0));
7249566063dSJacob Faibussowitsch   PetscCall(FVRHSFunction(ts, 0, X0, X, (void *)&ctx)); /* Initial function evaluation, only used to determine max speed */
7259566063dSJacob Faibussowitsch   PetscCall(VecCopy(X0, X));                            /* The function value was not used so we set X=X0 again */
7269566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, ctx.cfl / ctx.cfl_idt));
7279566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts)); /* Take runtime options */
7289566063dSJacob Faibussowitsch   PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD));
729c4762a1bSJed Brown   {
730c4762a1bSJed Brown     PetscInt           steps;
731c4762a1bSJed Brown     PetscScalar        mass_initial, mass_final, mass_difference, mass_differenceg;
732c4762a1bSJed Brown     const PetscScalar *ptr_X, *ptr_X0;
733c2a16db8SHong Zhang     const PetscReal    hs = (ctx.xmax - ctx.xmin) / 2.0 / count_slow;
734c2a16db8SHong Zhang     const PetscReal    hf = (ctx.xmax - ctx.xmin) / 2.0 / count_fast;
7359566063dSJacob Faibussowitsch     PetscCall(TSSolve(ts, X));
7369566063dSJacob Faibussowitsch     PetscCall(TSGetSolveTime(ts, &ptime));
7379566063dSJacob Faibussowitsch     PetscCall(TSGetStepNumber(ts, &steps));
738c4762a1bSJed Brown     /* calculate the total mass at initial time and final time */
739c4762a1bSJed Brown     mass_initial = 0.0;
740c4762a1bSJed Brown     mass_final   = 0.0;
7419566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArrayRead(da, X0, (void *)&ptr_X0));
7429566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArrayRead(da, X, (void *)&ptr_X));
743c2a16db8SHong Zhang     for (i = xs; i < xs + xm; i++) {
744c2a16db8SHong Zhang       if (i < ctx.sf || i > ctx.fs - 1) {
745c2a16db8SHong Zhang         for (k = 0; k < dof; k++) {
746c2a16db8SHong Zhang           mass_initial = mass_initial + hs * ptr_X0[i * dof + k];
747c2a16db8SHong Zhang           mass_final   = mass_final + hs * ptr_X[i * dof + k];
748c2a16db8SHong Zhang         }
749c4762a1bSJed Brown       } else {
750c2a16db8SHong Zhang         for (k = 0; k < dof; k++) {
751c2a16db8SHong Zhang           mass_initial = mass_initial + hf * ptr_X0[i * dof + k];
752c2a16db8SHong Zhang           mass_final   = mass_final + hf * ptr_X[i * dof + k];
753c2a16db8SHong Zhang         }
754c4762a1bSJed Brown       }
755c4762a1bSJed Brown     }
7569566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArrayRead(da, X0, (void *)&ptr_X0));
7579566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArrayRead(da, X, (void *)&ptr_X));
758c4762a1bSJed Brown     mass_difference = mass_final - mass_initial;
7599566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allreduce(&mass_difference, &mass_differenceg, 1, MPIU_SCALAR, MPIU_SUM, comm));
7609566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Mass difference %g\n", (double)mass_differenceg));
76163a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Final time %g, steps %" PetscInt_FMT "\n", (double)ptime, steps));
762c4762a1bSJed Brown     if (ctx.exact) {
763c4762a1bSJed Brown       PetscReal nrm1 = 0;
7649566063dSJacob Faibussowitsch       PetscCall(SolutionErrorNorms(&ctx, da, ptime, X, &nrm1));
7659566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %g\n", (double)nrm1));
766c4762a1bSJed Brown     }
767c4762a1bSJed Brown     if (ctx.simulation) {
768c4762a1bSJed Brown       PetscReal          nrm1 = 0;
769c4762a1bSJed Brown       PetscViewer        fd;
770c4762a1bSJed Brown       char               filename[PETSC_MAX_PATH_LEN] = "binaryoutput";
771c4762a1bSJed Brown       Vec                XR;
772c4762a1bSJed Brown       PetscBool          flg;
773c4762a1bSJed Brown       const PetscScalar *ptr_XR;
7749566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetString(NULL, NULL, "-f", filename, sizeof(filename), &flg));
7753c633725SBarry Smith       PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate binary file with the -f option");
7769566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, filename, FILE_MODE_READ, &fd));
7779566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(X0, &XR));
7789566063dSJacob Faibussowitsch       PetscCall(VecLoad(XR, fd));
7799566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&fd));
7809566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(X, &ptr_X));
7819566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(XR, &ptr_XR));
782c4762a1bSJed Brown       for (i = 0; i < Mx; i++) {
783c4762a1bSJed Brown         if (i < count_slow / 2 || i > count_slow / 2 + count_fast - 1) nrm1 = nrm1 + hs * PetscAbs(ptr_X[i] - ptr_XR[i]);
784c4762a1bSJed Brown         else nrm1 = nrm1 + hf * PetscAbs(ptr_X[i] - ptr_XR[i]);
785c4762a1bSJed Brown       }
7869566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(X, &ptr_X));
7879566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(XR, &ptr_XR));
7889566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %g\n", (double)nrm1));
7899566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&XR));
790c4762a1bSJed Brown     }
791c4762a1bSJed Brown   }
792c4762a1bSJed Brown 
7939566063dSJacob Faibussowitsch   PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD));
7949566063dSJacob Faibussowitsch   if (draw & 0x1) PetscCall(VecView(X0, PETSC_VIEWER_DRAW_WORLD));
7959566063dSJacob Faibussowitsch   if (draw & 0x2) PetscCall(VecView(X, PETSC_VIEWER_DRAW_WORLD));
796c4762a1bSJed Brown   if (draw & 0x4) {
797c4762a1bSJed Brown     Vec Y;
7989566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(X, &Y));
7999566063dSJacob Faibussowitsch     PetscCall(FVSample(&ctx, da, ptime, Y));
8009566063dSJacob Faibussowitsch     PetscCall(VecAYPX(Y, -1, X));
8019566063dSJacob Faibussowitsch     PetscCall(VecView(Y, PETSC_VIEWER_DRAW_WORLD));
8029566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&Y));
803c4762a1bSJed Brown   }
804c4762a1bSJed Brown 
805c4762a1bSJed Brown   if (view_final) {
806c4762a1bSJed Brown     PetscViewer viewer;
8079566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, final_fname, &viewer));
8089566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB));
8099566063dSJacob Faibussowitsch     PetscCall(VecView(X, viewer));
8109566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(viewer));
8119566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
812c4762a1bSJed Brown   }
813c4762a1bSJed Brown 
814c4762a1bSJed Brown   /* Clean up */
8159566063dSJacob Faibussowitsch   PetscCall((*ctx.physics.destroy)(ctx.physics.user));
8169566063dSJacob Faibussowitsch   for (i = 0; i < ctx.physics.dof; i++) PetscCall(PetscFree(ctx.physics.fieldname[i]));
8179566063dSJacob Faibussowitsch   PetscCall(PetscFree3(ctx.u, ctx.flux, ctx.speeds));
8189566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ctx.iss));
8199566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ctx.isf));
8209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
8219566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X0));
8229566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&R));
8239566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
8249566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
8259566063dSJacob Faibussowitsch   PetscCall(PetscFree(index_slow));
8269566063dSJacob Faibussowitsch   PetscCall(PetscFree(index_fast));
8279566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&physics));
8289566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
829b122ec5aSJacob Faibussowitsch   return 0;
830c4762a1bSJed Brown }
831c4762a1bSJed Brown 
832c4762a1bSJed Brown /*TEST
833c4762a1bSJed Brown 
834c4762a1bSJed Brown     build:
835f56ea12dSJed Brown       requires: !complex
836c4762a1bSJed Brown 
837c4762a1bSJed Brown     test:
838c4762a1bSJed Brown       args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 0
839c4762a1bSJed Brown 
840c4762a1bSJed Brown     test:
841c4762a1bSJed Brown       suffix: 2
842c4762a1bSJed Brown       args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 1
843c4762a1bSJed Brown       output_file: output/ex7_1.out
844c4762a1bSJed Brown 
845c4762a1bSJed Brown     test:
846c4762a1bSJed Brown       suffix: 3
847c4762a1bSJed Brown       args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 0
848c4762a1bSJed Brown 
849c4762a1bSJed Brown     test:
850c4762a1bSJed Brown       suffix: 4
851c4762a1bSJed Brown       args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 1
852c4762a1bSJed Brown       output_file: output/ex7_3.out
853c2a16db8SHong Zhang 
854c2a16db8SHong Zhang     test:
855c2a16db8SHong Zhang       suffix: 5
856c2a16db8SHong Zhang       nsize: 2
857c2a16db8SHong Zhang       args: -da_grid_x 60 -initial 7 -xmin -1 -xmax 1 -hratio 2 -ts_dt 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 1
858c2a16db8SHong Zhang       output_file: output/ex7_3.out
859c4762a1bSJed Brown TEST*/
860