xref: /petsc/src/ts/tutorials/multirate/ex5.c (revision 2a1887a77e7b2c6e00dd0ba96d1387c839460237)
1c4762a1bSJed Brown /*
2c4762a1bSJed Brown   Note:
3*188af4bfSBarry Smith   To performance a convergence study, a reference solution can be generated with a small stepsize and stored into a binary file, e.g. -ts_type ssp -ts_time_step 1e-5 -ts_max_steps 30000 -ts_view_solution binary:reference.bin
4c4762a1bSJed Brown   Errors can be computed in the following runs with -simulation -f reference.bin
5c4762a1bSJed Brown 
6c4762a1bSJed Brown   Multirate RK options:
7c4762a1bSJed Brown   -ts_rk_dtratio is the ratio between larger time step size and small time step size
8c4762a1bSJed Brown   -ts_rk_multirate_type has three choices: none (for single step RK)
9c4762a1bSJed Brown                                            combined (for multirate method and user just need to provide the same RHS with the single step RK)
10c4762a1bSJed Brown                                            partitioned (for multiraet method and user need to provide two RHS, one is for fast components and the orther is for slow components
11c4762a1bSJed Brown */
12c4762a1bSJed Brown 
13c4762a1bSJed Brown static const char help[] = "1D periodic Finite Volume solver in slope-limiter form with semidiscrete time stepping.\n"
14c4762a1bSJed Brown                            " advection   - Variable coefficient scalar advection\n"
15c4762a1bSJed Brown                            "                u_t       + (a*u)_x               = 0\n"
16c4762a1bSJed Brown                            " a is a piecewise function with two values say a0 and a1 (both a0 and a1 are constant).\n"
17c4762a1bSJed Brown                            " in this toy problem we assume a=a0 when x<0 and a=a1 when x>0 with a0<a1 which has the same discontinuous point with initial condition.\n"
18c4762a1bSJed Brown                            " we don't provide an exact solution, so you need to generate reference solution to do convergnece staudy,\n"
19aaa8cc7dSPierre Jolivet                            " more precisely, you need firstly generate a reference to a binary file say file.bin, then on command line,\n"
20c4762a1bSJed Brown                            " you should type -simulation -f file.bin.\n"
21c4762a1bSJed Brown                            " you can choose the number of grids by -da_grid_x.\n"
22c4762a1bSJed Brown                            " you can choose the value of a by -physics_advect_a1 and -physics_advect_a2.\n";
23c4762a1bSJed Brown 
24c4762a1bSJed Brown #include <petscts.h>
25c4762a1bSJed Brown #include <petscdm.h>
26c4762a1bSJed Brown #include <petscdmda.h>
27c4762a1bSJed Brown #include <petscdraw.h>
28c4762a1bSJed Brown #include <petsc/private/tsimpl.h>
29c4762a1bSJed Brown 
30c4762a1bSJed Brown #include <petsc/private/kernels/blockinvert.h> /* For the Kernel_*_gets_* stuff for BAIJ */
31c4762a1bSJed Brown 
32c4762a1bSJed Brown #include "finitevolume1d.h"
33c4762a1bSJed Brown 
RangeMod(PetscReal a,PetscReal xmin,PetscReal xmax)34d71ae5a4SJacob Faibussowitsch static inline PetscReal RangeMod(PetscReal a, PetscReal xmin, PetscReal xmax)
35d71ae5a4SJacob Faibussowitsch {
369371c9d4SSatish Balay   PetscReal range = xmax - xmin;
379371c9d4SSatish Balay   return xmin + PetscFmodReal(range + PetscFmodReal(a, range), range);
389371c9d4SSatish Balay }
39c4762a1bSJed Brown 
40c4762a1bSJed Brown /* --------------------------------- Advection ----------------------------------- */
41c4762a1bSJed Brown 
42c4762a1bSJed Brown typedef struct {
43c4762a1bSJed Brown   PetscReal a[2]; /* advective velocity */
44c4762a1bSJed Brown } AdvectCtx;
45c4762a1bSJed Brown 
PhysicsRiemann_Advect(void * vctx,PetscInt m,const PetscScalar * uL,const PetscScalar * uR,PetscScalar * flux,PetscReal * maxspeed,PetscReal x,PetscReal xmin,PetscReal xmax)46d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsRiemann_Advect(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed, PetscReal x, PetscReal xmin, PetscReal xmax)
47d71ae5a4SJacob Faibussowitsch {
48c4762a1bSJed Brown   AdvectCtx *ctx = (AdvectCtx *)vctx;
49c4762a1bSJed Brown   PetscReal *speed;
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   PetscFunctionBeginUser;
52c4762a1bSJed Brown   speed = ctx->a;
539371c9d4SSatish Balay   if (x == 0 || x == xmin || x == xmax)
549371c9d4SSatish Balay     flux[0] = PetscMax(0, (speed[0] + speed[1]) / 2.0) * uL[0] + PetscMin(0, (speed[0] + speed[1]) / 2.0) * uR[0]; /* if condition need to be changed base on different problem, '0' is the discontinuous point of a */
55da81f932SPierre Jolivet   else if (x < 0) flux[0] = PetscMax(0, speed[0]) * uL[0] + PetscMin(0, speed[0]) * uR[0];                         /* else if condition need to be changed based on different problem, 'x = 0' is discontinuous point of a */
56c4762a1bSJed Brown   else flux[0] = PetscMax(0, speed[1]) * uL[0] + PetscMin(0, speed[1]) * uR[0];
57c4762a1bSJed Brown   *maxspeed = *speed;
583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
59c4762a1bSJed Brown }
60c4762a1bSJed Brown 
PhysicsCharacteristic_Advect(void * vctx,PetscInt m,const PetscScalar * u,PetscScalar * X,PetscScalar * Xi,PetscReal * speeds,PetscReal x)61d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx, PetscInt m, const PetscScalar *u, PetscScalar *X, PetscScalar *Xi, PetscReal *speeds, PetscReal x)
62d71ae5a4SJacob Faibussowitsch {
63c4762a1bSJed Brown   AdvectCtx *ctx = (AdvectCtx *)vctx;
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   PetscFunctionBeginUser;
66c4762a1bSJed Brown   X[0]  = 1.;
67c4762a1bSJed Brown   Xi[0] = 1.;
68c4762a1bSJed Brown   if (x < 0) speeds[0] = ctx->a[0]; /* x is discontinuous point of a */
69c4762a1bSJed Brown   else speeds[0] = ctx->a[1];
703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
71c4762a1bSJed Brown }
72c4762a1bSJed Brown 
PhysicsSample_Advect(void * vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal * u)73d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsSample_Advect(void *vctx, PetscInt initial, FVBCType bctype, PetscReal xmin, PetscReal xmax, PetscReal t, PetscReal x, PetscReal *u)
74d71ae5a4SJacob Faibussowitsch {
75c4762a1bSJed Brown   AdvectCtx *ctx = (AdvectCtx *)vctx;
76c4762a1bSJed Brown   PetscReal *a   = ctx->a, x0;
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   PetscFunctionBeginUser;
79c4762a1bSJed Brown   if (x < 0) { /* x is cell center */
80c4762a1bSJed Brown     switch (bctype) {
81d71ae5a4SJacob Faibussowitsch     case FVBC_OUTFLOW:
82d71ae5a4SJacob Faibussowitsch       x0 = x - a[0] * t;
83d71ae5a4SJacob Faibussowitsch       break;
84d71ae5a4SJacob Faibussowitsch     case FVBC_PERIODIC:
85d71ae5a4SJacob Faibussowitsch       x0 = RangeMod(x - a[0] * t, xmin, xmax);
86d71ae5a4SJacob Faibussowitsch       break;
87d71ae5a4SJacob Faibussowitsch     default:
88d71ae5a4SJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown BCType");
89c4762a1bSJed Brown     }
90c4762a1bSJed Brown     switch (initial) {
91d71ae5a4SJacob Faibussowitsch     case 0:
92d71ae5a4SJacob Faibussowitsch       u[0] = (x0 < 0) ? 1 : -1;
93d71ae5a4SJacob Faibussowitsch       break;
94d71ae5a4SJacob Faibussowitsch     case 1:
95d71ae5a4SJacob Faibussowitsch       u[0] = (x0 < 0) ? -1 : 1;
96d71ae5a4SJacob Faibussowitsch       break;
97d71ae5a4SJacob Faibussowitsch     case 2:
98d71ae5a4SJacob Faibussowitsch       u[0] = (0 < x0 && x0 < 1) ? 1 : 0;
99d71ae5a4SJacob Faibussowitsch       break;
100d71ae5a4SJacob Faibussowitsch     case 3:
101d71ae5a4SJacob Faibussowitsch       u[0] = PetscSinReal(2 * PETSC_PI * x0);
102d71ae5a4SJacob Faibussowitsch       break;
103d71ae5a4SJacob Faibussowitsch     case 4:
104d71ae5a4SJacob Faibussowitsch       u[0] = PetscAbs(x0);
105d71ae5a4SJacob Faibussowitsch       break;
106d71ae5a4SJacob Faibussowitsch     case 5:
107d71ae5a4SJacob Faibussowitsch       u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2 * PETSC_PI * x0));
108d71ae5a4SJacob Faibussowitsch       break;
109d71ae5a4SJacob Faibussowitsch     case 6:
110d71ae5a4SJacob Faibussowitsch       u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2 - x0 : 0));
111d71ae5a4SJacob Faibussowitsch       break;
112d71ae5a4SJacob Faibussowitsch     case 7:
113d71ae5a4SJacob Faibussowitsch       u[0] = PetscPowReal(PetscSinReal(PETSC_PI * x0), 10.0);
114d71ae5a4SJacob Faibussowitsch       break;
115d71ae5a4SJacob Faibussowitsch     default:
116d71ae5a4SJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition");
117c4762a1bSJed Brown     }
1189371c9d4SSatish Balay   } else {
119c4762a1bSJed Brown     switch (bctype) {
120d71ae5a4SJacob Faibussowitsch     case FVBC_OUTFLOW:
121d71ae5a4SJacob Faibussowitsch       x0 = x - a[1] * t;
122d71ae5a4SJacob Faibussowitsch       break;
123d71ae5a4SJacob Faibussowitsch     case FVBC_PERIODIC:
124d71ae5a4SJacob Faibussowitsch       x0 = RangeMod(x - a[1] * t, xmin, xmax);
125d71ae5a4SJacob Faibussowitsch       break;
126d71ae5a4SJacob Faibussowitsch     default:
127d71ae5a4SJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown BCType");
128c4762a1bSJed Brown     }
129c4762a1bSJed Brown     switch (initial) {
130d71ae5a4SJacob Faibussowitsch     case 0:
131d71ae5a4SJacob Faibussowitsch       u[0] = (x0 < 0) ? 1 : -1;
132d71ae5a4SJacob Faibussowitsch       break;
133d71ae5a4SJacob Faibussowitsch     case 1:
134d71ae5a4SJacob Faibussowitsch       u[0] = (x0 < 0) ? -1 : 1;
135d71ae5a4SJacob Faibussowitsch       break;
136d71ae5a4SJacob Faibussowitsch     case 2:
137d71ae5a4SJacob Faibussowitsch       u[0] = (0 < x0 && x0 < 1) ? 1 : 0;
138d71ae5a4SJacob Faibussowitsch       break;
139d71ae5a4SJacob Faibussowitsch     case 3:
140d71ae5a4SJacob Faibussowitsch       u[0] = PetscSinReal(2 * PETSC_PI * x0);
141d71ae5a4SJacob Faibussowitsch       break;
142d71ae5a4SJacob Faibussowitsch     case 4:
143d71ae5a4SJacob Faibussowitsch       u[0] = PetscAbs(x0);
144d71ae5a4SJacob Faibussowitsch       break;
145d71ae5a4SJacob Faibussowitsch     case 5:
146d71ae5a4SJacob Faibussowitsch       u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2 * PETSC_PI * x0));
147d71ae5a4SJacob Faibussowitsch       break;
148d71ae5a4SJacob Faibussowitsch     case 6:
149d71ae5a4SJacob Faibussowitsch       u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2 - x0 : 0));
150d71ae5a4SJacob Faibussowitsch       break;
151d71ae5a4SJacob Faibussowitsch     case 7:
152d71ae5a4SJacob Faibussowitsch       u[0] = PetscPowReal(PetscSinReal(PETSC_PI * x0), 10.0);
153d71ae5a4SJacob Faibussowitsch       break;
154d71ae5a4SJacob Faibussowitsch     default:
155d71ae5a4SJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition");
156c4762a1bSJed Brown     }
157c4762a1bSJed Brown   }
1583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
159c4762a1bSJed Brown }
160c4762a1bSJed Brown 
PhysicsCreate_Advect(FVCtx * ctx)161d71ae5a4SJacob Faibussowitsch static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx)
162d71ae5a4SJacob Faibussowitsch {
163c4762a1bSJed Brown   AdvectCtx *user;
164c4762a1bSJed Brown 
165c4762a1bSJed Brown   PetscFunctionBeginUser;
1669566063dSJacob Faibussowitsch   PetscCall(PetscNew(&user));
167c4762a1bSJed Brown   ctx->physics.sample         = PhysicsSample_Advect;
168c4762a1bSJed Brown   ctx->physics.riemann        = PhysicsRiemann_Advect;
169c4762a1bSJed Brown   ctx->physics.characteristic = PhysicsCharacteristic_Advect;
170c4762a1bSJed Brown   ctx->physics.destroy        = PhysicsDestroy_SimpleFree;
171c4762a1bSJed Brown   ctx->physics.user           = user;
172c4762a1bSJed Brown   ctx->physics.dof            = 1;
1739566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("u", &ctx->physics.fieldname[0]));
174c4762a1bSJed Brown   user->a[0] = 1;
175c4762a1bSJed Brown   user->a[1] = 1;
176d0609cedSBarry Smith   PetscOptionsBegin(ctx->comm, ctx->prefix, "Options for advection", "");
177c4762a1bSJed Brown   {
1789566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-physics_advect_a1", "Speed1", "", user->a[0], &user->a[0], NULL));
1799566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-physics_advect_a2", "Speed2", "", user->a[1], &user->a[1], NULL));
180c4762a1bSJed Brown   }
181d0609cedSBarry Smith   PetscOptionsEnd();
1823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
183c4762a1bSJed Brown }
184c4762a1bSJed Brown 
185c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for slow parts ----------------------------------- */
186c4762a1bSJed Brown 
FVRHSFunctionslow(TS ts,PetscReal time,Vec X,Vec F,void * vctx)187d71ae5a4SJacob Faibussowitsch PetscErrorCode FVRHSFunctionslow(TS ts, PetscReal time, Vec X, Vec F, void *vctx)
188d71ae5a4SJacob Faibussowitsch {
189c4762a1bSJed Brown   FVCtx       *ctx = (FVCtx *)vctx;
190c4762a1bSJed Brown   PetscInt     i, j, k, Mx, dof, xs, xm, len_slow;
191c4762a1bSJed Brown   PetscReal    hx, cfl_idt = 0;
192c4762a1bSJed Brown   PetscScalar *x, *f, *slope;
193c4762a1bSJed Brown   Vec          Xloc;
194c4762a1bSJed Brown   DM           da;
195c4762a1bSJed Brown 
196c4762a1bSJed Brown   PetscFunctionBeginUser;
1979566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
1989566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &Xloc));
1999566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
200c4762a1bSJed Brown   hx = (ctx->xmax - ctx->xmin) / Mx;
2019566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
2029566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
2039566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
2049566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Xloc, &x));
2059566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
2069566063dSJacob Faibussowitsch   PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
2079566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
2089566063dSJacob Faibussowitsch   PetscCall(ISGetSize(ctx->iss, &len_slow));
209c4762a1bSJed Brown 
210c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
211c4762a1bSJed Brown     for (i = xs - 2; i < 0; i++) {
212c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
213c4762a1bSJed Brown     }
214c4762a1bSJed Brown     for (i = Mx; i < xs + xm + 2; i++) {
215c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
216c4762a1bSJed Brown     }
217c4762a1bSJed Brown   }
218c4762a1bSJed Brown   for (i = xs - 1; i < xs + xm + 1; i++) {
219c4762a1bSJed Brown     struct _LimitInfo info;
220c4762a1bSJed Brown     PetscScalar      *cjmpL, *cjmpR;
221c4762a1bSJed Brown     if (i < len_slow + 1) {
222c4762a1bSJed Brown       /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
2239566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.characteristic)(ctx->physics.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds, ctx->xmin + hx * i));
224c4762a1bSJed Brown       /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
2259566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
226c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
227c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
228c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
229c4762a1bSJed Brown         PetscScalar jmpL, jmpR;
230c4762a1bSJed Brown         jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
231c4762a1bSJed Brown         jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
232c4762a1bSJed Brown         for (k = 0; k < dof; k++) {
233c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
234c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
235c4762a1bSJed Brown         }
236c4762a1bSJed Brown       }
237c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
238c4762a1bSJed Brown       info.m  = dof;
239c4762a1bSJed Brown       info.hx = hx;
240c4762a1bSJed Brown       (*ctx->limit)(&info, cjmpL, cjmpR, ctx->cslope);
241c4762a1bSJed Brown       for (j = 0; j < dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */
242c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
243c4762a1bSJed Brown         PetscScalar tmp = 0;
244c4762a1bSJed Brown         for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
245c4762a1bSJed Brown         slope[i * dof + j] = tmp;
246c4762a1bSJed Brown       }
247c4762a1bSJed Brown     }
248c4762a1bSJed Brown   }
249c4762a1bSJed Brown 
250c4762a1bSJed Brown   for (i = xs; i < xs + xm + 1; i++) {
251c4762a1bSJed Brown     PetscReal    maxspeed;
252c4762a1bSJed Brown     PetscScalar *uL, *uR;
253c4762a1bSJed Brown     if (i < len_slow) { /* slow parts can be changed based on a */
254c4762a1bSJed Brown       uL = &ctx->uLR[0];
255c4762a1bSJed Brown       uR = &ctx->uLR[dof];
256c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
257c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2;
258c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2;
259c4762a1bSJed Brown       }
2609566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax));
261c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */
262c4762a1bSJed Brown       if (i > xs) {
263c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hx;
264c4762a1bSJed Brown       }
265c4762a1bSJed Brown       if (i < xs + xm) {
266c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hx;
267c4762a1bSJed Brown       }
268c4762a1bSJed Brown     }
269c4762a1bSJed Brown     if (i == len_slow) { /* slow parts can be changed based on a */
270c4762a1bSJed Brown       uL = &ctx->uLR[0];
271c4762a1bSJed Brown       uR = &ctx->uLR[dof];
272c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
273c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2;
274c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2;
275c4762a1bSJed Brown       }
2769566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax));
277c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */
278c4762a1bSJed Brown       if (i > xs) {
279c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hx;
280c4762a1bSJed Brown       }
281c4762a1bSJed Brown     }
282c4762a1bSJed Brown   }
2839566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
2849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
2859566063dSJacob Faibussowitsch   PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
2869566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &Xloc));
2873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
288c4762a1bSJed Brown }
289c4762a1bSJed Brown 
290c4762a1bSJed Brown /* --------------------------------- Finite Volume Solver for fast  parts ----------------------------------- */
291c4762a1bSJed Brown 
FVRHSFunctionfast(TS ts,PetscReal time,Vec X,Vec F,void * vctx)292d71ae5a4SJacob Faibussowitsch PetscErrorCode FVRHSFunctionfast(TS ts, PetscReal time, Vec X, Vec F, void *vctx)
293d71ae5a4SJacob Faibussowitsch {
294c4762a1bSJed Brown   FVCtx       *ctx = (FVCtx *)vctx;
295c4762a1bSJed Brown   PetscInt     i, j, k, Mx, dof, xs, xm, len_slow;
296c4762a1bSJed Brown   PetscReal    hx, cfl_idt = 0;
297c4762a1bSJed Brown   PetscScalar *x, *f, *slope;
298c4762a1bSJed Brown   Vec          Xloc;
299c4762a1bSJed Brown   DM           da;
300c4762a1bSJed Brown 
301c4762a1bSJed Brown   PetscFunctionBeginUser;
3029566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
3039566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &Xloc));
3049566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
305c4762a1bSJed Brown   hx = (ctx->xmax - ctx->xmin) / Mx;
3069566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
3079566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
3089566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
3099566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Xloc, &x));
3109566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
3119566063dSJacob Faibussowitsch   PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
3129566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
3139566063dSJacob Faibussowitsch   PetscCall(ISGetSize(ctx->iss, &len_slow));
314c4762a1bSJed Brown 
315c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
316c4762a1bSJed Brown     for (i = xs - 2; i < 0; i++) {
317c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
318c4762a1bSJed Brown     }
319c4762a1bSJed Brown     for (i = Mx; i < xs + xm + 2; i++) {
320c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
321c4762a1bSJed Brown     }
322c4762a1bSJed Brown   }
323c4762a1bSJed Brown   for (i = xs - 1; i < xs + xm + 1; i++) {
324c4762a1bSJed Brown     struct _LimitInfo info;
325c4762a1bSJed Brown     PetscScalar      *cjmpL, *cjmpR;
326c4762a1bSJed Brown     if (i > len_slow - 2) {
3279566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.characteristic)(ctx->physics.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds, ctx->xmin + hx * i));
3289566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
329c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
330c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
331c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
332c4762a1bSJed Brown         PetscScalar jmpL, jmpR;
333c4762a1bSJed Brown         jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
334c4762a1bSJed Brown         jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
335c4762a1bSJed Brown         for (k = 0; k < dof; k++) {
336c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
337c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
338c4762a1bSJed Brown         }
339c4762a1bSJed Brown       }
340c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
341c4762a1bSJed Brown       info.m  = dof;
342c4762a1bSJed Brown       info.hx = hx;
343c4762a1bSJed Brown       (*ctx->limit)(&info, cjmpL, cjmpR, ctx->cslope);
344c4762a1bSJed Brown       for (j = 0; j < dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */
345c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
346c4762a1bSJed Brown         PetscScalar tmp = 0;
347c4762a1bSJed Brown         for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
348c4762a1bSJed Brown         slope[i * dof + j] = tmp;
349c4762a1bSJed Brown       }
350c4762a1bSJed Brown     }
351c4762a1bSJed Brown   }
352c4762a1bSJed Brown 
353c4762a1bSJed Brown   for (i = xs; i < xs + xm + 1; i++) {
354c4762a1bSJed Brown     PetscReal    maxspeed;
355c4762a1bSJed Brown     PetscScalar *uL, *uR;
356c4762a1bSJed Brown     if (i > len_slow) {
357c4762a1bSJed Brown       uL = &ctx->uLR[0];
358c4762a1bSJed Brown       uR = &ctx->uLR[dof];
359c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
360c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2;
361c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2;
362c4762a1bSJed Brown       }
3639566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax));
364c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */
365c4762a1bSJed Brown       if (i > xs) {
366c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - len_slow - 1) * dof + j] -= ctx->flux[j] / hx;
367c4762a1bSJed Brown       }
368c4762a1bSJed Brown       if (i < xs + xm) {
369c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - len_slow) * dof + j] += ctx->flux[j] / hx;
370c4762a1bSJed Brown       }
371c4762a1bSJed Brown     }
372c4762a1bSJed Brown     if (i == len_slow) {
373c4762a1bSJed Brown       uL = &ctx->uLR[0];
374c4762a1bSJed Brown       uR = &ctx->uLR[dof];
375c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
376c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2;
377c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2;
378c4762a1bSJed Brown       }
3799566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax));
380c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */
381c4762a1bSJed Brown       if (i < xs + xm) {
382c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - len_slow) * dof + j] += ctx->flux[j] / hx;
383c4762a1bSJed Brown       }
384c4762a1bSJed Brown     }
385c4762a1bSJed Brown   }
3869566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
3879566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
3889566063dSJacob Faibussowitsch   PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
3899566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &Xloc));
3903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
391c4762a1bSJed Brown }
392c4762a1bSJed Brown 
FVRHSFunctionslow2(TS ts,PetscReal time,Vec X,Vec F,void * vctx)393d71ae5a4SJacob Faibussowitsch PetscErrorCode FVRHSFunctionslow2(TS ts, PetscReal time, Vec X, Vec F, void *vctx)
394d71ae5a4SJacob Faibussowitsch {
395c4762a1bSJed Brown   FVCtx       *ctx = (FVCtx *)vctx;
396c4762a1bSJed Brown   PetscInt     i, j, k, Mx, dof, xs, xm, len_slow1, len_slow2;
397c4762a1bSJed Brown   PetscReal    hx, cfl_idt = 0;
398c4762a1bSJed Brown   PetscScalar *x, *f, *slope;
399c4762a1bSJed Brown   Vec          Xloc;
400c4762a1bSJed Brown   DM           da;
401c4762a1bSJed Brown 
402c4762a1bSJed Brown   PetscFunctionBeginUser;
4039566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
4049566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &Xloc));
4059566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
406c4762a1bSJed Brown   hx = (ctx->xmax - ctx->xmin) / Mx;
4079566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
4089566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
4099566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
4109566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Xloc, &x));
4119566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
4129566063dSJacob Faibussowitsch   PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
4139566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
4149566063dSJacob Faibussowitsch   PetscCall(ISGetSize(ctx->iss, &len_slow1));
4159566063dSJacob Faibussowitsch   PetscCall(ISGetSize(ctx->iss2, &len_slow2));
416c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
417c4762a1bSJed Brown     for (i = xs - 2; i < 0; i++) {
418c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
419c4762a1bSJed Brown     }
420c4762a1bSJed Brown     for (i = Mx; i < xs + xm + 2; i++) {
421c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
422c4762a1bSJed Brown     }
423c4762a1bSJed Brown   }
424c4762a1bSJed Brown   for (i = xs - 1; i < xs + xm + 1; i++) {
425c4762a1bSJed Brown     struct _LimitInfo info;
426c4762a1bSJed Brown     PetscScalar      *cjmpL, *cjmpR;
427c4762a1bSJed Brown     if (i < len_slow1 + len_slow2 + 1 && i > len_slow1 - 2) {
4289566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.characteristic)(ctx->physics.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds, ctx->xmin + hx * i));
4299566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
430c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
431c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
432c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
433c4762a1bSJed Brown         PetscScalar jmpL, jmpR;
434c4762a1bSJed Brown         jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
435c4762a1bSJed Brown         jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
436c4762a1bSJed Brown         for (k = 0; k < dof; k++) {
437c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
438c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
439c4762a1bSJed Brown         }
440c4762a1bSJed Brown       }
441c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
442c4762a1bSJed Brown       info.m  = dof;
443c4762a1bSJed Brown       info.hx = hx;
444c4762a1bSJed Brown       (*ctx->limit)(&info, cjmpL, cjmpR, ctx->cslope);
445c4762a1bSJed Brown       for (j = 0; j < dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */
446c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
447c4762a1bSJed Brown         PetscScalar tmp = 0;
448c4762a1bSJed Brown         for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
449c4762a1bSJed Brown         slope[i * dof + j] = tmp;
450c4762a1bSJed Brown       }
451c4762a1bSJed Brown     }
452c4762a1bSJed Brown   }
453c4762a1bSJed Brown 
454c4762a1bSJed Brown   for (i = xs; i < xs + xm + 1; i++) {
455c4762a1bSJed Brown     PetscReal    maxspeed;
456c4762a1bSJed Brown     PetscScalar *uL, *uR;
457c4762a1bSJed Brown     if (i < len_slow1 + len_slow2 && i > len_slow1) { /* slow parts can be changed based on a */
458c4762a1bSJed Brown       uL = &ctx->uLR[0];
459c4762a1bSJed Brown       uR = &ctx->uLR[dof];
460c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
461c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2;
462c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2;
463c4762a1bSJed Brown       }
4649566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax));
465c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */
466c4762a1bSJed Brown       if (i > xs) {
467c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - len_slow1 - 1) * dof + j] -= ctx->flux[j] / hx;
468c4762a1bSJed Brown       }
469c4762a1bSJed Brown       if (i < xs + xm) {
470c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - len_slow1) * dof + j] += ctx->flux[j] / hx;
471c4762a1bSJed Brown       }
472c4762a1bSJed Brown     }
473c4762a1bSJed Brown     if (i == len_slow1 + len_slow2) { /* slow parts can be changed based on a */
474c4762a1bSJed Brown       uL = &ctx->uLR[0];
475c4762a1bSJed Brown       uR = &ctx->uLR[dof];
476c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
477c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2;
478c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2;
479c4762a1bSJed Brown       }
4809566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax));
481c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */
482c4762a1bSJed Brown       if (i > xs) {
483c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - len_slow1 - 1) * dof + j] -= ctx->flux[j] / hx;
484c4762a1bSJed Brown       }
485c4762a1bSJed Brown     }
486c4762a1bSJed Brown     if (i == len_slow1) { /* slow parts can be changed based on a */
487c4762a1bSJed Brown       uL = &ctx->uLR[0];
488c4762a1bSJed Brown       uR = &ctx->uLR[dof];
489c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
490c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2;
491c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2;
492c4762a1bSJed Brown       }
4939566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax));
494c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */
495c4762a1bSJed Brown       if (i < xs + xm) {
496c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - len_slow1) * dof + j] += ctx->flux[j] / hx;
497c4762a1bSJed Brown       }
498c4762a1bSJed Brown     }
499c4762a1bSJed Brown   }
500c4762a1bSJed Brown 
5019566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
5029566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
5039566063dSJacob Faibussowitsch   PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
5049566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &Xloc));
5053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
506c4762a1bSJed Brown }
507c4762a1bSJed Brown 
FVRHSFunctionfast2(TS ts,PetscReal time,Vec X,Vec F,void * vctx)508d71ae5a4SJacob Faibussowitsch PetscErrorCode FVRHSFunctionfast2(TS ts, PetscReal time, Vec X, Vec F, void *vctx)
509d71ae5a4SJacob Faibussowitsch {
510c4762a1bSJed Brown   FVCtx       *ctx = (FVCtx *)vctx;
511c4762a1bSJed Brown   PetscInt     i, j, k, Mx, dof, xs, xm, len_slow1, len_slow2;
512c4762a1bSJed Brown   PetscReal    hx, cfl_idt = 0;
513c4762a1bSJed Brown   PetscScalar *x, *f, *slope;
514c4762a1bSJed Brown   Vec          Xloc;
515c4762a1bSJed Brown   DM           da;
516c4762a1bSJed Brown 
517c4762a1bSJed Brown   PetscFunctionBeginUser;
5189566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
5199566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &Xloc));
5209566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
521c4762a1bSJed Brown   hx = (ctx->xmax - ctx->xmin) / Mx;
5229566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
5239566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
5249566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
5259566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Xloc, &x));
5269566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
5279566063dSJacob Faibussowitsch   PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));
5289566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
5299566063dSJacob Faibussowitsch   PetscCall(ISGetSize(ctx->iss, &len_slow1));
5309566063dSJacob Faibussowitsch   PetscCall(ISGetSize(ctx->iss2, &len_slow2));
531c4762a1bSJed Brown 
532c4762a1bSJed Brown   if (ctx->bctype == FVBC_OUTFLOW) {
533c4762a1bSJed Brown     for (i = xs - 2; i < 0; i++) {
534c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
535c4762a1bSJed Brown     }
536c4762a1bSJed Brown     for (i = Mx; i < xs + xm + 2; i++) {
537c4762a1bSJed Brown       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
538c4762a1bSJed Brown     }
539c4762a1bSJed Brown   }
540c4762a1bSJed Brown   for (i = xs - 1; i < xs + xm + 1; i++) {
541c4762a1bSJed Brown     struct _LimitInfo info;
542c4762a1bSJed Brown     PetscScalar      *cjmpL, *cjmpR;
543c4762a1bSJed Brown     if (i > len_slow1 + len_slow2 - 2) {
5449566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.characteristic)(ctx->physics.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds, ctx->xmin + hx * i));
5459566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
546c4762a1bSJed Brown       cjmpL = &ctx->cjmpLR[0];
547c4762a1bSJed Brown       cjmpR = &ctx->cjmpLR[dof];
548c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
549c4762a1bSJed Brown         PetscScalar jmpL, jmpR;
550c4762a1bSJed Brown         jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
551c4762a1bSJed Brown         jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
552c4762a1bSJed Brown         for (k = 0; k < dof; k++) {
553c4762a1bSJed Brown           cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
554c4762a1bSJed Brown           cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
555c4762a1bSJed Brown         }
556c4762a1bSJed Brown       }
557c4762a1bSJed Brown       /* Apply limiter to the left and right characteristic jumps */
558c4762a1bSJed Brown       info.m  = dof;
559c4762a1bSJed Brown       info.hx = hx;
560c4762a1bSJed Brown       (*ctx->limit)(&info, cjmpL, cjmpR, ctx->cslope);
561c4762a1bSJed Brown       for (j = 0; j < dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */
562c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
563c4762a1bSJed Brown         PetscScalar tmp = 0;
564c4762a1bSJed Brown         for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
565c4762a1bSJed Brown         slope[i * dof + j] = tmp;
566c4762a1bSJed Brown       }
567c4762a1bSJed Brown     }
568c4762a1bSJed Brown   }
569c4762a1bSJed Brown 
570c4762a1bSJed Brown   for (i = xs; i < xs + xm + 1; i++) {
571c4762a1bSJed Brown     PetscReal    maxspeed;
572c4762a1bSJed Brown     PetscScalar *uL, *uR;
573c4762a1bSJed Brown     if (i > len_slow1 + len_slow2) {
574c4762a1bSJed Brown       uL = &ctx->uLR[0];
575c4762a1bSJed Brown       uR = &ctx->uLR[dof];
576c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
577c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2;
578c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2;
579c4762a1bSJed Brown       }
5809566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax));
581c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */
582c4762a1bSJed Brown       if (i > xs) {
583c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - len_slow1 - len_slow2 - 1) * dof + j] -= ctx->flux[j] / hx;
584c4762a1bSJed Brown       }
585c4762a1bSJed Brown       if (i < xs + xm) {
586c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - len_slow1 - len_slow2) * dof + j] += ctx->flux[j] / hx;
587c4762a1bSJed Brown       }
588c4762a1bSJed Brown     }
589c4762a1bSJed Brown     if (i == len_slow1 + len_slow2) {
590c4762a1bSJed Brown       uL = &ctx->uLR[0];
591c4762a1bSJed Brown       uR = &ctx->uLR[dof];
592c4762a1bSJed Brown       for (j = 0; j < dof; j++) {
593c4762a1bSJed Brown         uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2;
594c4762a1bSJed Brown         uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2;
595c4762a1bSJed Brown       }
5969566063dSJacob Faibussowitsch       PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax));
597c4762a1bSJed Brown       cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */
598c4762a1bSJed Brown       if (i < xs + xm) {
599c4762a1bSJed Brown         for (j = 0; j < dof; j++) f[(i - len_slow1 - len_slow2) * dof + j] += ctx->flux[j] / hx;
600c4762a1bSJed Brown       }
601c4762a1bSJed Brown     }
602c4762a1bSJed Brown   }
6039566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
6049566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
6059566063dSJacob Faibussowitsch   PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
6069566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &Xloc));
6073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
608c4762a1bSJed Brown }
609c4762a1bSJed Brown 
main(int argc,char * argv[])610d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[])
611d71ae5a4SJacob Faibussowitsch {
612c4762a1bSJed Brown   char              lname[256] = "mc", physname[256] = "advect", final_fname[256] = "solution.m";
613c4762a1bSJed Brown   PetscFunctionList limiters = 0, physics = 0;
614c4762a1bSJed Brown   MPI_Comm          comm;
615c4762a1bSJed Brown   TS                ts;
616c4762a1bSJed Brown   DM                da;
617c4762a1bSJed Brown   Vec               X, X0, R;
618c4762a1bSJed Brown   FVCtx             ctx;
619c4762a1bSJed Brown   PetscInt          i, k, dof, xs, xm, Mx, draw = 0, *index_slow, *index_fast, islow = 0, ifast = 0;
620c4762a1bSJed Brown   PetscBool         view_final = PETSC_FALSE;
621c4762a1bSJed Brown   PetscReal         ptime;
622c4762a1bSJed Brown 
623327415f7SBarry Smith   PetscFunctionBeginUser;
6249566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, 0, help));
625c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
6269566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&ctx, sizeof(ctx)));
627c4762a1bSJed Brown 
628c4762a1bSJed Brown   /* Register limiters to be available on the command line */
6299566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "upwind", Limit_Upwind));
6309566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "lax-wendroff", Limit_LaxWendroff));
6319566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "beam-warming", Limit_BeamWarming));
6329566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "fromm", Limit_Fromm));
6339566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "minmod", Limit_Minmod));
6349566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "superbee", Limit_Superbee));
6359566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "mc", Limit_MC));
6369566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "vanleer", Limit_VanLeer));
6379566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "vanalbada", Limit_VanAlbada));
6389566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "vanalbadatvd", Limit_VanAlbadaTVD));
6399566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "koren", Limit_Koren));
6409566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "korensym", Limit_KorenSym));
6419566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "koren3", Limit_Koren3));
6429566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "cada-torrilhon2", Limit_CadaTorrilhon2));
6439566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "cada-torrilhon3-r0p1", Limit_CadaTorrilhon3R0p1));
6449566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "cada-torrilhon3-r1", Limit_CadaTorrilhon3R1));
6459566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "cada-torrilhon3-r10", Limit_CadaTorrilhon3R10));
6469566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&limiters, "cada-torrilhon3-r100", Limit_CadaTorrilhon3R100));
647c4762a1bSJed Brown 
648c4762a1bSJed Brown   /* Register physical models to be available on the command line */
6499566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&physics, "advect", PhysicsCreate_Advect));
650c4762a1bSJed Brown 
651c4762a1bSJed Brown   ctx.comm   = comm;
652c4762a1bSJed Brown   ctx.cfl    = 0.9;
653c4762a1bSJed Brown   ctx.bctype = FVBC_PERIODIC;
654c4762a1bSJed Brown   ctx.xmin   = -1.0;
655c4762a1bSJed Brown   ctx.xmax   = 1.0;
656d0609cedSBarry Smith   PetscOptionsBegin(comm, NULL, "Finite Volume solver options", "");
6579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-xmin", "X min", "", ctx.xmin, &ctx.xmin, NULL));
6589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-xmax", "X max", "", ctx.xmax, &ctx.xmax, NULL));
6599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-limit", "Name of flux limiter to use", "", limiters, lname, lname, sizeof(lname), NULL));
6609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-draw", "Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)", "", draw, &draw, NULL));
6619566063dSJacob 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));
6629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-initial", "Initial condition (depends on the physics)", "", ctx.initial, &ctx.initial, NULL));
6639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-simulation", "Compare errors with reference solution", "", ctx.simulation, &ctx.simulation, NULL));
6649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-cfl", "CFL number to time step at", "", ctx.cfl, &ctx.cfl, NULL));
6659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-bc_type", "Boundary condition", "", FVBCTypes, (PetscEnum)ctx.bctype, (PetscEnum *)&ctx.bctype, NULL));
6669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-recursive_split", "Split the domain recursively", "", ctx.recursive, &ctx.recursive, NULL));
667d0609cedSBarry Smith   PetscOptionsEnd();
668c4762a1bSJed Brown 
669c4762a1bSJed Brown   /* Choose the limiter from the list of registered limiters */
6709566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(limiters, lname, &ctx.limit));
6713c633725SBarry Smith   PetscCheck(ctx.limit, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Limiter '%s' not found", lname);
672c4762a1bSJed Brown 
673c4762a1bSJed Brown   /* Choose the physics from the list of registered models */
674c4762a1bSJed Brown   {
675c4762a1bSJed Brown     PetscErrorCode (*r)(FVCtx *);
6769566063dSJacob Faibussowitsch     PetscCall(PetscFunctionListFind(physics, physname, &r));
6773c633725SBarry Smith     PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Physics '%s' not found", physname);
678c4762a1bSJed Brown     /* Create the physics, will set the number of fields and their names */
6799566063dSJacob Faibussowitsch     PetscCall((*r)(&ctx));
680c4762a1bSJed Brown   }
681c4762a1bSJed Brown 
682c4762a1bSJed Brown   /* Create a DMDA to manage the parallel grid */
6839566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(comm, DM_BOUNDARY_PERIODIC, 50, ctx.physics.dof, 2, NULL, &da));
6849566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
6859566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
686c4762a1bSJed Brown   /* Inform the DMDA of the field names provided by the physics. */
687c4762a1bSJed Brown   /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */
68848a46eb9SPierre Jolivet   for (i = 0; i < ctx.physics.dof; i++) PetscCall(DMDASetFieldName(da, i, ctx.physics.fieldname[i]));
6899566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
6909566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
691c4762a1bSJed Brown 
692c4762a1bSJed Brown   /* Set coordinates of cell centers */
6939566063dSJacob 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));
694c4762a1bSJed Brown 
695c4762a1bSJed Brown   /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */
6969566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(dof * dof, &ctx.R, dof * dof, &ctx.Rinv, 2 * dof, &ctx.cjmpLR, 1 * dof, &ctx.cslope));
6979566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(2 * dof, &ctx.uLR, dof, &ctx.flux, dof, &ctx.speeds));
698c4762a1bSJed Brown 
699c4762a1bSJed Brown   /* Create a vector to store the solution and to save the initial state */
7009566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &X));
7019566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X, &X0));
7029566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X, &R));
703c4762a1bSJed Brown 
704c4762a1bSJed Brown   /*-------------------------------- create index for slow parts and fast parts ----------------------------------------*/
7059566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(xm * dof, &index_slow));
7069566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(xm * dof, &index_fast));
707c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
708c4762a1bSJed Brown     if (ctx.xmin + i * (ctx.xmax - ctx.xmin) / (PetscReal)Mx + 0.5 * (ctx.xmax - ctx.xmin) / (PetscReal)Mx < 0)
709c4762a1bSJed Brown       for (k = 0; k < dof; k++) index_slow[islow++] = i * dof + k;
710c4762a1bSJed Brown     else
711c4762a1bSJed Brown       for (k = 0; k < dof; k++) index_fast[ifast++] = i * dof + k;
712c4762a1bSJed Brown   } /* this step need to be changed based on discontinuous point of a */
7139566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islow, index_slow, PETSC_COPY_VALUES, &ctx.iss));
7149566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, ifast, index_fast, PETSC_COPY_VALUES, &ctx.isf));
715c4762a1bSJed Brown 
716c4762a1bSJed Brown   /* Create a time-stepping object */
7179566063dSJacob Faibussowitsch   PetscCall(TSCreate(comm, &ts));
7189566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, da));
7199566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, R, FVRHSFunction, &ctx));
7209566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetIS(ts, "slow", ctx.iss));
7219566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetIS(ts, "fast", ctx.isf));
7229566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetRHSFunction(ts, "slow", NULL, FVRHSFunctionslow, &ctx));
7239566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetRHSFunction(ts, "fast", NULL, FVRHSFunctionfast, &ctx));
724c4762a1bSJed Brown 
725c4762a1bSJed Brown   if (ctx.recursive) {
726c4762a1bSJed Brown     TS subts;
727c4762a1bSJed Brown     islow = 0;
728c4762a1bSJed Brown     ifast = 0;
729c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
730c4762a1bSJed Brown       PetscReal coord = ctx.xmin + i * (ctx.xmax - ctx.xmin) / (PetscReal)Mx + 0.5 * (ctx.xmax - ctx.xmin) / (PetscReal)Mx;
731c4762a1bSJed Brown       if (coord >= 0 && coord < ctx.xmin + (ctx.xmax - ctx.xmin) * 3 / 4.)
732c4762a1bSJed Brown         for (k = 0; k < dof; k++) index_slow[islow++] = i * dof + k;
733c4762a1bSJed Brown       if (coord >= ctx.xmin + (ctx.xmax - ctx.xmin) * 3 / 4.)
734c4762a1bSJed Brown         for (k = 0; k < dof; k++) index_fast[ifast++] = i * dof + k;
735c4762a1bSJed Brown     } /* this step need to be changed based on discontinuous point of a */
7369566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, islow, index_slow, PETSC_COPY_VALUES, &ctx.iss2));
7379566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, ifast, index_fast, PETSC_COPY_VALUES, &ctx.isf2));
738c4762a1bSJed Brown 
7399566063dSJacob Faibussowitsch     PetscCall(TSRHSSplitGetSubTS(ts, "fast", &subts));
7409566063dSJacob Faibussowitsch     PetscCall(TSRHSSplitSetIS(subts, "slow", ctx.iss2));
7419566063dSJacob Faibussowitsch     PetscCall(TSRHSSplitSetIS(subts, "fast", ctx.isf2));
7429566063dSJacob Faibussowitsch     PetscCall(TSRHSSplitSetRHSFunction(subts, "slow", NULL, FVRHSFunctionslow2, &ctx));
7439566063dSJacob Faibussowitsch     PetscCall(TSRHSSplitSetRHSFunction(subts, "fast", NULL, FVRHSFunctionfast2, &ctx));
744c4762a1bSJed Brown   }
745c4762a1bSJed Brown 
7469566063dSJacob Faibussowitsch   /*PetscCall(TSSetType(ts,TSSSP));*/
7479566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSMPRK));
7489566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 10));
7499566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
750c4762a1bSJed Brown 
751c4762a1bSJed Brown   /* Compute initial conditions and starting time step */
7529566063dSJacob Faibussowitsch   PetscCall(FVSample(&ctx, da, 0, X0));
7539566063dSJacob Faibussowitsch   PetscCall(FVRHSFunction(ts, 0, X0, X, (void *)&ctx)); /* Initial function evaluation, only used to determine max speed */
7549566063dSJacob Faibussowitsch   PetscCall(VecCopy(X0, X));                            /* The function value was not used so we set X=X0 again */
7559566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, ctx.cfl / ctx.cfl_idt));
7569566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts)); /* Take runtime options */
7579566063dSJacob Faibussowitsch   PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD));
758c4762a1bSJed Brown   {
759c4762a1bSJed Brown     PetscInt    steps;
760c4762a1bSJed Brown     PetscScalar mass_initial, mass_final, mass_difference;
761c4762a1bSJed Brown 
7629566063dSJacob Faibussowitsch     PetscCall(TSSolve(ts, X));
7639566063dSJacob Faibussowitsch     PetscCall(TSGetSolveTime(ts, &ptime));
7649566063dSJacob Faibussowitsch     PetscCall(TSGetStepNumber(ts, &steps));
76563a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Final time %g, steps %" PetscInt_FMT "\n", (double)ptime, steps));
766c4762a1bSJed Brown     /* calculate the total mass at initial time and final time */
767c4762a1bSJed Brown     mass_initial = 0.0;
768c4762a1bSJed Brown     mass_final   = 0.0;
7699566063dSJacob Faibussowitsch     PetscCall(VecSum(X0, &mass_initial));
7709566063dSJacob Faibussowitsch     PetscCall(VecSum(X, &mass_final));
771c4762a1bSJed Brown     mass_difference = (ctx.xmax - ctx.xmin) / (PetscScalar)Mx * (mass_final - mass_initial);
7729566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Mass difference %g\n", (double)mass_difference));
773c4762a1bSJed Brown     if (ctx.simulation) {
774c4762a1bSJed Brown       PetscViewer fd;
775c4762a1bSJed Brown       char        filename[PETSC_MAX_PATH_LEN] = "binaryoutput";
776c4762a1bSJed Brown       Vec         XR;
777c4762a1bSJed Brown       PetscReal   nrm1, nrmsup;
778c4762a1bSJed Brown       PetscBool   flg;
779d8185827SBarry Smith 
7809566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetString(NULL, NULL, "-f", filename, sizeof(filename), &flg));
7813c633725SBarry Smith       PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate binary file with the -f option");
7829566063dSJacob Faibussowitsch       PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, filename, FILE_MODE_READ, &fd));
7839566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(X0, &XR));
7849566063dSJacob Faibussowitsch       PetscCall(VecLoad(XR, fd));
7859566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&fd));
7869566063dSJacob Faibussowitsch       PetscCall(VecAYPX(XR, -1, X));
7879566063dSJacob Faibussowitsch       PetscCall(VecNorm(XR, NORM_1, &nrm1));
7889566063dSJacob Faibussowitsch       PetscCall(VecNorm(XR, NORM_INFINITY, &nrmsup));
789c4762a1bSJed Brown       nrm1 /= Mx;
7909566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %g  ||x-x_e||_sup %g\n", (double)nrm1, (double)nrmsup));
7919566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&XR));
792c4762a1bSJed Brown     }
793c4762a1bSJed Brown   }
794c4762a1bSJed Brown 
7959566063dSJacob Faibussowitsch   PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD));
7969566063dSJacob Faibussowitsch   if (draw & 0x1) PetscCall(VecView(X0, PETSC_VIEWER_DRAW_WORLD));
7979566063dSJacob Faibussowitsch   if (draw & 0x2) PetscCall(VecView(X, PETSC_VIEWER_DRAW_WORLD));
798c4762a1bSJed Brown   if (draw & 0x4) {
799c4762a1bSJed Brown     Vec Y;
8009566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(X, &Y));
8019566063dSJacob Faibussowitsch     PetscCall(FVSample(&ctx, da, ptime, Y));
8029566063dSJacob Faibussowitsch     PetscCall(VecAYPX(Y, -1, X));
8039566063dSJacob Faibussowitsch     PetscCall(VecView(Y, PETSC_VIEWER_DRAW_WORLD));
8049566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&Y));
805c4762a1bSJed Brown   }
806c4762a1bSJed Brown 
807c4762a1bSJed Brown   if (view_final) {
808c4762a1bSJed Brown     PetscViewer viewer;
8099566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, final_fname, &viewer));
8109566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB));
8119566063dSJacob Faibussowitsch     PetscCall(VecView(X, viewer));
8129566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(viewer));
8139566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
814c4762a1bSJed Brown   }
815c4762a1bSJed Brown 
816c4762a1bSJed Brown   /* Clean up */
8179566063dSJacob Faibussowitsch   PetscCall((*ctx.physics.destroy)(ctx.physics.user));
8189566063dSJacob Faibussowitsch   for (i = 0; i < ctx.physics.dof; i++) PetscCall(PetscFree(ctx.physics.fieldname[i]));
8199566063dSJacob Faibussowitsch   PetscCall(PetscFree4(ctx.R, ctx.Rinv, ctx.cjmpLR, ctx.cslope));
8209566063dSJacob Faibussowitsch   PetscCall(PetscFree3(ctx.uLR, ctx.flux, ctx.speeds));
8219566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ctx.iss));
8229566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ctx.isf));
8239566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ctx.iss2));
8249566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ctx.isf2));
8259566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
8269566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X0));
8279566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&R));
8289566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
8299566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
8309566063dSJacob Faibussowitsch   PetscCall(PetscFree(index_slow));
8319566063dSJacob Faibussowitsch   PetscCall(PetscFree(index_fast));
8329566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&limiters));
8339566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&physics));
8349566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
835b122ec5aSJacob Faibussowitsch   return 0;
836c4762a1bSJed Brown }
837c4762a1bSJed Brown 
838c4762a1bSJed Brown /*TEST
839c4762a1bSJed Brown     build:
840f56ea12dSJed Brown       requires: !complex
841c4762a1bSJed Brown       depends: finitevolume1d.c
842c4762a1bSJed Brown 
843c4762a1bSJed Brown     test:
844*188af4bfSBarry Smith       args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_time_step 0.025 -ts_max_steps 24 -ts_type rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 0
845c4762a1bSJed Brown 
846c4762a1bSJed Brown     test:
847c4762a1bSJed Brown       suffix: 2
848*188af4bfSBarry Smith       args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_time_step 0.025 -ts_max_steps 24 -ts_type rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 1
849c4762a1bSJed Brown       output_file: output/ex5_1.out
850c4762a1bSJed Brown 
851c4762a1bSJed Brown     test:
852c4762a1bSJed Brown       suffix: 3
853*188af4bfSBarry Smith       args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_time_step 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 0
854c4762a1bSJed Brown 
855c4762a1bSJed Brown     test:
856c4762a1bSJed Brown       suffix: 4
857*188af4bfSBarry Smith       args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_time_step 0.025 -ts_max_steps 24 -ts_type mprk -ts_mprk_type 2a22 -ts_use_splitrhsfunction 1
858c4762a1bSJed Brown       output_file: output/ex5_3.out
859c4762a1bSJed Brown 
860c4762a1bSJed Brown     test:
861c4762a1bSJed Brown       suffix: 5
862*188af4bfSBarry Smith       args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_time_step 0.025 -ts_max_steps 24 -ts_type rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 0 -recursive_split
863c4762a1bSJed Brown 
864c4762a1bSJed Brown     test:
865c4762a1bSJed Brown       suffix: 6
866*188af4bfSBarry Smith       args: -da_grid_x 60 -initial 1 -xmin -1 -xmax 1 -limit mc -physics_advect_a1 1 -physics_advect_a2 2 -ts_time_step 0.025 -ts_max_steps 24 -ts_type rk -ts_rk_type 2a -ts_rk_dtratio 2 -ts_rk_multirate -ts_use_splitrhsfunction 1 -recursive_split
867c4762a1bSJed Brown TEST*/
868