Actual source code: ex9.c

  1: static const char help[] = "1D periodic Finite Volume solver in slope-limiter form with semidiscrete time stepping.\n"
  2:                            "Solves scalar and vector problems, choose the physical model with -physics\n"
  3:                            "  advection   - Constant coefficient scalar advection\n"
  4:                            "                u_t       + (a*u)_x               = 0\n"
  5:                            "  burgers     - Burgers equation\n"
  6:                            "                u_t       + (u^2/2)_x             = 0\n"
  7:                            "  traffic     - Traffic equation\n"
  8:                            "                u_t       + (u*(1-u))_x           = 0\n"
  9:                            "  acoustics   - Acoustic wave propagation\n"
 10:                            "                u_t       + (c*z*v)_x             = 0\n"
 11:                            "                v_t       + (c/z*u)_x             = 0\n"
 12:                            "  isogas      - Isothermal gas dynamics\n"
 13:                            "                rho_t     + (rho*u)_x             = 0\n"
 14:                            "                (rho*u)_t + (rho*u^2 + c^2*rho)_x = 0\n"
 15:                            "  shallow     - Shallow water equations\n"
 16:                            "                h_t       + (h*u)_x               = 0\n"
 17:                            "                (h*u)_t   + (h*u^2 + g*h^2/2)_x   = 0\n"
 18:                            "Some of these physical models have multiple Riemann solvers, select these with -physics_xxx_riemann\n"
 19:                            "  exact       - Exact Riemann solver which usually needs to perform a Newton iteration to connect\n"
 20:                            "                the states across shocks and rarefactions\n"
 21:                            "  roe         - Linearized scheme, usually with an entropy fix inside sonic rarefactions\n"
 22:                            "The systems provide a choice of reconstructions with -physics_xxx_reconstruct\n"
 23:                            "  characteristic - Limit the characteristic variables, this is usually preferred (default)\n"
 24:                            "  conservative   - Limit the conservative variables directly, can cause undesired interaction of waves\n\n"
 25:                            "A variety of limiters for high-resolution TVD limiters are available with -limit\n"
 26:                            "  upwind,minmod,superbee,mc,vanleer,vanalbada,koren,cada-torillhon (last two are nominally third order)\n"
 27:                            "  and non-TVD schemes lax-wendroff,beam-warming,fromm\n\n"
 28:                            "To preserve the TVD property, one should time step with a strong stability preserving method.\n"
 29:                            "The optimal high order explicit Runge-Kutta methods in TSSSP are recommended for non-stiff problems.\n\n"
 30:                            "Several initial conditions can be chosen with -initial N\n\n"
 31:                            "The problem size should be set with -da_grid_x M\n\n";

 33: #include <petscts.h>
 34: #include <petscdm.h>
 35: #include <petscdmda.h>
 36: #include <petscdraw.h>

 38: #include <petsc/private/kernels/blockinvert.h>

 40: static inline PetscReal Sgn(PetscReal a)
 41: {
 42:   return (a < 0) ? -1 : 1;
 43: }
 44: static inline PetscReal Abs(PetscReal a)
 45: {
 46:   return (a < 0) ? 0 : a;
 47: }
 48: static inline PetscReal Sqr(PetscReal a)
 49: {
 50:   return a * a;
 51: }
 52: static inline PetscReal MaxAbs(PetscReal a, PetscReal b)
 53: {
 54:   return (PetscAbs(a) > PetscAbs(b)) ? a : b;
 55: }
 56: PETSC_UNUSED static inline PetscReal MinAbs(PetscReal a, PetscReal b)
 57: {
 58:   return (PetscAbs(a) < PetscAbs(b)) ? a : b;
 59: }
 60: static inline PetscReal MinMod2(PetscReal a, PetscReal b)
 61: {
 62:   return (a * b < 0) ? 0 : Sgn(a) * PetscMin(PetscAbs(a), PetscAbs(b));
 63: }
 64: static inline PetscReal MaxMod2(PetscReal a, PetscReal b)
 65: {
 66:   return (a * b < 0) ? 0 : Sgn(a) * PetscMax(PetscAbs(a), PetscAbs(b));
 67: }
 68: static inline PetscReal MinMod3(PetscReal a, PetscReal b, PetscReal c)
 69: {
 70:   return (a * b < 0 || a * c < 0) ? 0 : Sgn(a) * PetscMin(PetscAbs(a), PetscMin(PetscAbs(b), PetscAbs(c)));
 71: }

 73: static inline PetscReal RangeMod(PetscReal a, PetscReal xmin, PetscReal xmax)
 74: {
 75:   PetscReal range = xmax - xmin;
 76:   return xmin + PetscFmodReal(range + PetscFmodReal(a, range), range);
 77: }

 79: /* ----------------------- Lots of limiters, these could go in a separate library ------------------------- */
 80: typedef struct _LimitInfo {
 81:   PetscReal hx;
 82:   PetscInt  m;
 83: }          *LimitInfo;
 84: static void Limit_Upwind(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt)
 85: {
 86:   PetscInt i;
 87:   for (i = 0; i < info->m; i++) lmt[i] = 0;
 88: }
 89: static void Limit_LaxWendroff(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt)
 90: {
 91:   PetscInt i;
 92:   for (i = 0; i < info->m; i++) lmt[i] = jR[i];
 93: }
 94: static void Limit_BeamWarming(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt)
 95: {
 96:   PetscInt i;
 97:   for (i = 0; i < info->m; i++) lmt[i] = jL[i];
 98: }
 99: static void Limit_Fromm(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt)
100: {
101:   PetscInt i;
102:   for (i = 0; i < info->m; i++) lmt[i] = 0.5 * (jL[i] + jR[i]);
103: }
104: static void Limit_Minmod(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt)
105: {
106:   PetscInt i;
107:   for (i = 0; i < info->m; i++) lmt[i] = MinMod2(jL[i], jR[i]);
108: }
109: static void Limit_Superbee(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt)
110: {
111:   PetscInt i;
112:   for (i = 0; i < info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i], 2 * jR[i]), MinMod2(2 * jL[i], jR[i]));
113: }
114: static void Limit_MC(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt)
115: {
116:   PetscInt i;
117:   for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i], 0.5 * (jL[i] + jR[i]), 2 * jR[i]);
118: }
119: static void Limit_VanLeer(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt)
120: { /* phi = (t + abs(t)) / (1 + abs(t)) */
121:   PetscInt i;
122:   for (i = 0; i < info->m; i++) lmt[i] = (jL[i] * Abs(jR[i]) + Abs(jL[i]) * jR[i]) / (Abs(jL[i]) + Abs(jR[i]) + 1e-15);
123: }
124: static void Limit_VanAlbada(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt) /* differentiable */
125: {                                                                                                           /* phi = (t + t^2) / (1 + t^2) */
126:   PetscInt i;
127:   for (i = 0; i < info->m; i++) lmt[i] = (jL[i] * Sqr(jR[i]) + Sqr(jL[i]) * jR[i]) / (Sqr(jL[i]) + Sqr(jR[i]) + 1e-15);
128: }
129: static void Limit_VanAlbadaTVD(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt)
130: { /* phi = (t + t^2) / (1 + t^2) */
131:   PetscInt i;
132:   for (i = 0; i < info->m; i++) lmt[i] = (jL[i] * jR[i] < 0) ? 0 : (jL[i] * Sqr(jR[i]) + Sqr(jL[i]) * jR[i]) / (Sqr(jL[i]) + Sqr(jR[i]) + 1e-15);
133: }
134: static void Limit_Koren(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt) /* differentiable */
135: {                                                                                                       /* phi = (t + 2*t^2) / (2 - t + 2*t^2) */
136:   PetscInt i;
137:   for (i = 0; i < info->m; i++) lmt[i] = ((jL[i] * Sqr(jR[i]) + 2 * Sqr(jL[i]) * jR[i]) / (2 * Sqr(jL[i]) - jL[i] * jR[i] + 2 * Sqr(jR[i]) + 1e-15));
138: }
139: static void Limit_KorenSym(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt) /* differentiable */
140: {                                                                                                          /* Symmetric version of above */
141:   PetscInt i;
142:   for (i = 0; i < info->m; i++) lmt[i] = (1.5 * (jL[i] * Sqr(jR[i]) + Sqr(jL[i]) * jR[i]) / (2 * Sqr(jL[i]) - jL[i] * jR[i] + 2 * Sqr(jR[i]) + 1e-15));
143: }
144: static void Limit_Koren3(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt)
145: { /* Eq 11 of Cada-Torrilhon 2009 */
146:   PetscInt i;
147:   for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i], (jL[i] + 2 * jR[i]) / 3, 2 * jR[i]);
148: }
149: static PetscReal CadaTorrilhonPhiHatR_Eq13(PetscReal L, PetscReal R)
150: {
151:   return PetscMax(0, PetscMin((L + 2 * R) / 3, PetscMax(-0.5 * L, PetscMin(2 * L, PetscMin((L + 2 * R) / 3, 1.6 * R)))));
152: }
153: static void Limit_CadaTorrilhon2(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt)
154: { /* Cada-Torrilhon 2009, Eq 13 */
155:   PetscInt i;
156:   for (i = 0; i < info->m; i++) lmt[i] = CadaTorrilhonPhiHatR_Eq13(jL[i], jR[i]);
157: }
158: static void Limit_CadaTorrilhon3R(PetscReal r, LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt)
159: { /* Cada-Torrilhon 2009, Eq 22 */
160:   /* They recommend 0.001 < r < 1, but larger values are more accurate in smooth regions */
161:   const PetscReal eps = 1e-7, hx = info->hx;
162:   for (PetscInt i = 0; i < info->m; i++) {
163:     const PetscReal eta = (Sqr(jL[i]) + Sqr(jR[i])) / Sqr(r * hx);
164:     lmt[i] = ((eta < 1 - eps) ? (jL[i] + 2 * jR[i]) / 3 : ((eta > 1 + eps) ? CadaTorrilhonPhiHatR_Eq13(jL[i], jR[i]) : 0.5 * ((1 - (eta - 1) / eps) * (jL[i] + 2 * jR[i]) / 3 + (1 + (eta + 1) / eps) * CadaTorrilhonPhiHatR_Eq13(jL[i], jR[i]))));
165:   }
166: }
167: static void Limit_CadaTorrilhon3R0p1(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt)
168: {
169:   Limit_CadaTorrilhon3R(0.1, info, jL, jR, lmt);
170: }
171: static void Limit_CadaTorrilhon3R1(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt)
172: {
173:   Limit_CadaTorrilhon3R(1, info, jL, jR, lmt);
174: }
175: static void Limit_CadaTorrilhon3R10(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt)
176: {
177:   Limit_CadaTorrilhon3R(10, info, jL, jR, lmt);
178: }
179: static void Limit_CadaTorrilhon3R100(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt)
180: {
181:   Limit_CadaTorrilhon3R(100, info, jL, jR, lmt);
182: }

184: /* --------------------------------- Finite Volume data structures ----------------------------------- */

186: typedef enum {
187:   FVBC_PERIODIC,
188:   FVBC_OUTFLOW
189: } FVBCType;
190: static const char *FVBCTypes[] = {"PERIODIC", "OUTFLOW", "FVBCType", "FVBC_", 0};
191: typedef PetscErrorCode (*RiemannFunction)(void *, PetscInt, const PetscScalar *, const PetscScalar *, PetscScalar *, PetscReal *);
192: typedef PetscErrorCode (*ReconstructFunction)(void *, PetscInt, const PetscScalar *, PetscScalar *, PetscScalar *, PetscReal *);

194: typedef struct {
195:   PetscErrorCode (*sample)(void *, PetscInt, FVBCType, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal *);
196:   RiemannFunction     riemann;
197:   ReconstructFunction characteristic;
198:   PetscErrorCode (*destroy)(void *);
199:   void    *user;
200:   PetscInt dof;
201:   char    *fieldname[16];
202: } PhysicsCtx;

204: typedef struct {
205:   void (*limit)(LimitInfo, const PetscScalar *, const PetscScalar *, PetscScalar *);
206:   PhysicsCtx physics;
207:   MPI_Comm   comm;
208:   char       prefix[256];

210:   /* Local work arrays */
211:   PetscScalar *R, *Rinv; /* Characteristic basis, and it's inverse.  COLUMN-MAJOR */
212:   PetscScalar *cjmpLR;   /* Jumps at left and right edge of cell, in characteristic basis, len=2*dof */
213:   PetscScalar *cslope;   /* Limited slope, written in characteristic basis */
214:   PetscScalar *uLR;      /* Solution at left and right of interface, conservative variables, len=2*dof */
215:   PetscScalar *flux;     /* Flux across interface */
216:   PetscReal   *speeds;   /* Speeds of each wave */

218:   PetscReal cfl_idt; /* Max allowable value of 1/Delta t */
219:   PetscReal cfl;
220:   PetscReal xmin, xmax;
221:   PetscInt  initial;
222:   PetscBool exact;
223:   FVBCType  bctype;
224: } FVCtx;

226: PetscErrorCode RiemannListAdd(PetscFunctionList *flist, const char *name, RiemannFunction rsolve)
227: {
228:   PetscFunctionBeginUser;
229:   PetscCall(PetscFunctionListAdd(flist, name, rsolve));
230:   PetscFunctionReturn(PETSC_SUCCESS);
231: }

233: PetscErrorCode RiemannListFind(PetscFunctionList flist, const char *name, RiemannFunction *rsolve)
234: {
235:   PetscFunctionBeginUser;
236:   PetscCall(PetscFunctionListFind(flist, name, rsolve));
237:   PetscCheck(*rsolve, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Riemann solver \"%s\" could not be found", name);
238:   PetscFunctionReturn(PETSC_SUCCESS);
239: }

241: PetscErrorCode ReconstructListAdd(PetscFunctionList *flist, const char *name, ReconstructFunction r)
242: {
243:   PetscFunctionBeginUser;
244:   PetscCall(PetscFunctionListAdd(flist, name, r));
245:   PetscFunctionReturn(PETSC_SUCCESS);
246: }

248: PetscErrorCode ReconstructListFind(PetscFunctionList flist, const char *name, ReconstructFunction *r)
249: {
250:   PetscFunctionBeginUser;
251:   PetscCall(PetscFunctionListFind(flist, name, r));
252:   PetscCheck(*r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Reconstruction \"%s\" could not be found", name);
253:   PetscFunctionReturn(PETSC_SUCCESS);
254: }

256: /* --------------------------------- Physics ----------------------------------- */
257: /*
258:   Each physical model consists of Riemann solver and a function to determine the basis to use for reconstruction.  These
259:   are set with the PhysicsCreate_XXX function which allocates private storage and sets these methods as well as the
260:   number of fields and their names, and a function to deallocate private storage.
261: */

263: /* First a few functions useful to several different physics */
264: static PetscErrorCode PhysicsCharacteristic_Conservative(void *vctx, PetscInt m, const PetscScalar *u, PetscScalar *X, PetscScalar *Xi, PetscReal *speeds)
265: {
266:   PetscFunctionBeginUser;
267:   for (PetscInt i = 0; i < m; i++) {
268:     for (PetscInt j = 0; j < m; j++) Xi[i * m + j] = X[i * m + j] = (PetscScalar)(i == j);
269:     speeds[i] = PETSC_MAX_REAL; /* Indicates invalid */
270:   }
271:   PetscFunctionReturn(PETSC_SUCCESS);
272: }

274: static PetscErrorCode PhysicsDestroy_SimpleFree(void *vctx)
275: {
276:   PetscFunctionBeginUser;
277:   PetscCall(PetscFree(vctx));
278:   PetscFunctionReturn(PETSC_SUCCESS);
279: }

281: /* --------------------------------- Advection ----------------------------------- */

283: typedef struct {
284:   PetscReal a; /* advective velocity */
285: } AdvectCtx;

287: static PetscErrorCode PhysicsRiemann_Advect(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed)
288: {
289:   AdvectCtx *ctx = (AdvectCtx *)vctx;
290:   PetscReal  speed;

292:   PetscFunctionBeginUser;
293:   speed     = ctx->a;
294:   flux[0]   = PetscMax(0, speed) * uL[0] + PetscMin(0, speed) * uR[0];
295:   *maxspeed = speed;
296:   PetscFunctionReturn(PETSC_SUCCESS);
297: }

299: static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx, PetscInt m, const PetscScalar *u, PetscScalar *X, PetscScalar *Xi, PetscReal *speeds)
300: {
301:   AdvectCtx *ctx = (AdvectCtx *)vctx;

303:   PetscFunctionBeginUser;
304:   X[0]      = 1.;
305:   Xi[0]     = 1.;
306:   speeds[0] = ctx->a;
307:   PetscFunctionReturn(PETSC_SUCCESS);
308: }

310: static PetscErrorCode PhysicsSample_Advect(void *vctx, PetscInt initial, FVBCType bctype, PetscReal xmin, PetscReal xmax, PetscReal t, PetscReal x, PetscReal *u)
311: {
312:   AdvectCtx *ctx = (AdvectCtx *)vctx;
313:   PetscReal  a   = ctx->a, x0;

315:   PetscFunctionBeginUser;
316:   switch (bctype) {
317:   case FVBC_OUTFLOW:
318:     x0 = x - a * t;
319:     break;
320:   case FVBC_PERIODIC:
321:     x0 = RangeMod(x - a * t, xmin, xmax);
322:     break;
323:   default:
324:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown BCType");
325:   }
326:   switch (initial) {
327:   case 0:
328:     u[0] = (x0 < 0) ? 1 : -1;
329:     break;
330:   case 1:
331:     u[0] = (x0 < 0) ? -1 : 1;
332:     break;
333:   case 2:
334:     u[0] = (0 < x0 && x0 < 1) ? 1 : 0;
335:     break;
336:   case 3:
337:     u[0] = PetscSinReal(2 * PETSC_PI * x0);
338:     break;
339:   case 4:
340:     u[0] = PetscAbs(x0);
341:     break;
342:   case 5:
343:     u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2 * PETSC_PI * x0));
344:     break;
345:   case 6:
346:     u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2 - x0 : 0));
347:     break;
348:   case 7:
349:     u[0] = PetscPowReal(PetscSinReal(PETSC_PI * x0), 10.0);
350:     break;
351:   default:
352:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition");
353:   }
354:   PetscFunctionReturn(PETSC_SUCCESS);
355: }

357: static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx)
358: {
359:   AdvectCtx *user;

361:   PetscFunctionBeginUser;
362:   PetscCall(PetscNew(&user));
363:   ctx->physics.sample         = PhysicsSample_Advect;
364:   ctx->physics.riemann        = PhysicsRiemann_Advect;
365:   ctx->physics.characteristic = PhysicsCharacteristic_Advect;
366:   ctx->physics.destroy        = PhysicsDestroy_SimpleFree;
367:   ctx->physics.user           = user;
368:   ctx->physics.dof            = 1;
369:   PetscCall(PetscStrallocpy("u", &ctx->physics.fieldname[0]));
370:   user->a = 1;
371:   PetscOptionsBegin(ctx->comm, ctx->prefix, "Options for advection", "");
372:   {
373:     PetscCall(PetscOptionsReal("-physics_advect_a", "Speed", "", user->a, &user->a, NULL));
374:   }
375:   PetscOptionsEnd();
376:   PetscFunctionReturn(PETSC_SUCCESS);
377: }

379: /* --------------------------------- Burgers ----------------------------------- */

381: typedef struct {
382:   PetscReal lxf_speed;
383: } BurgersCtx;

385: static PetscErrorCode PhysicsSample_Burgers(void *vctx, PetscInt initial, FVBCType bctype, PetscReal xmin, PetscReal xmax, PetscReal t, PetscReal x, PetscReal *u)
386: {
387:   PetscFunctionBeginUser;
388:   PetscCheck(bctype != FVBC_PERIODIC || t <= 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "Exact solution not implemented for periodic");
389:   switch (initial) {
390:   case 0:
391:     u[0] = (x < 0) ? 1 : -1;
392:     break;
393:   case 1:
394:     if (x < -t) u[0] = -1;
395:     else if (x < t) u[0] = x / t;
396:     else u[0] = 1;
397:     break;
398:   case 2:
399:     if (x <= 0) u[0] = 0;
400:     else if (x < t) u[0] = x / t;
401:     else if (x < 1 + 0.5 * t) u[0] = 1;
402:     else u[0] = 0;
403:     break;
404:   case 3:
405:     if (x < 0.2 * t) u[0] = 0.2;
406:     else if (x < t) u[0] = x / t;
407:     else u[0] = 1;
408:     break;
409:   case 4:
410:     PetscCheck(t <= 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only initial condition available");
411:     u[0] = 0.7 + 0.3 * PetscSinReal(2 * PETSC_PI * ((x - xmin) / (xmax - xmin)));
412:     break;
413:   case 5: /* Pure shock solution */
414:     if (x < 0.5 * t) u[0] = 1;
415:     else u[0] = 0;
416:     break;
417:   default:
418:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition");
419:   }
420:   PetscFunctionReturn(PETSC_SUCCESS);
421: }

423: static PetscErrorCode PhysicsRiemann_Burgers_Exact(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed)
424: {
425:   PetscFunctionBeginUser;
426:   if (uL[0] < uR[0]) {                /* rarefaction */
427:     flux[0] = (uL[0] * uR[0] < 0) ? 0 /* sonic rarefaction */
428:                                   : 0.5 * PetscMin(PetscSqr(uL[0]), PetscSqr(uR[0]));
429:   } else { /* shock */
430:     flux[0] = 0.5 * PetscMax(PetscSqr(uL[0]), PetscSqr(uR[0]));
431:   }
432:   *maxspeed = (PetscAbs(uL[0]) > PetscAbs(uR[0])) ? uL[0] : uR[0];
433:   PetscFunctionReturn(PETSC_SUCCESS);
434: }

436: static PetscErrorCode PhysicsRiemann_Burgers_Roe(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed)
437: {
438:   PetscReal speed;

440:   PetscFunctionBeginUser;
441:   speed   = 0.5 * (uL[0] + uR[0]);
442:   flux[0] = 0.25 * (PetscSqr(uL[0]) + PetscSqr(uR[0])) - 0.5 * PetscAbs(speed) * (uR[0] - uL[0]);
443:   if (uL[0] <= 0 && 0 <= uR[0]) flux[0] = 0; /* Entropy fix for sonic rarefaction */
444:   *maxspeed = speed;
445:   PetscFunctionReturn(PETSC_SUCCESS);
446: }

448: static PetscErrorCode PhysicsRiemann_Burgers_LxF(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed)
449: {
450:   PetscReal   c;
451:   PetscScalar fL, fR;

453:   PetscFunctionBeginUser;
454:   c         = ((BurgersCtx *)vctx)->lxf_speed;
455:   fL        = 0.5 * PetscSqr(uL[0]);
456:   fR        = 0.5 * PetscSqr(uR[0]);
457:   flux[0]   = 0.5 * (fL + fR) - 0.5 * c * (uR[0] - uL[0]);
458:   *maxspeed = c;
459:   PetscFunctionReturn(PETSC_SUCCESS);
460: }

462: static PetscErrorCode PhysicsRiemann_Burgers_Rusanov(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed)
463: {
464:   PetscReal   c;
465:   PetscScalar fL, fR;

467:   PetscFunctionBeginUser;
468:   c         = PetscMax(PetscAbs(uL[0]), PetscAbs(uR[0]));
469:   fL        = 0.5 * PetscSqr(uL[0]);
470:   fR        = 0.5 * PetscSqr(uR[0]);
471:   flux[0]   = 0.5 * (fL + fR) - 0.5 * c * (uR[0] - uL[0]);
472:   *maxspeed = c;
473:   PetscFunctionReturn(PETSC_SUCCESS);
474: }

476: static PetscErrorCode PhysicsCreate_Burgers(FVCtx *ctx)
477: {
478:   BurgersCtx       *user;
479:   RiemannFunction   r;
480:   PetscFunctionList rlist      = 0;
481:   char              rname[256] = "exact";

483:   PetscFunctionBeginUser;
484:   PetscCall(PetscNew(&user));

486:   ctx->physics.sample         = PhysicsSample_Burgers;
487:   ctx->physics.characteristic = PhysicsCharacteristic_Conservative;
488:   ctx->physics.destroy        = PhysicsDestroy_SimpleFree;
489:   ctx->physics.user           = user;
490:   ctx->physics.dof            = 1;

492:   PetscCall(PetscStrallocpy("u", &ctx->physics.fieldname[0]));
493:   PetscCall(RiemannListAdd(&rlist, "exact", PhysicsRiemann_Burgers_Exact));
494:   PetscCall(RiemannListAdd(&rlist, "roe", PhysicsRiemann_Burgers_Roe));
495:   PetscCall(RiemannListAdd(&rlist, "lxf", PhysicsRiemann_Burgers_LxF));
496:   PetscCall(RiemannListAdd(&rlist, "rusanov", PhysicsRiemann_Burgers_Rusanov));
497:   PetscOptionsBegin(ctx->comm, ctx->prefix, "Options for advection", "");
498:   {
499:     PetscCall(PetscOptionsFList("-physics_burgers_riemann", "Riemann solver", "", rlist, rname, rname, sizeof(rname), NULL));
500:   }
501:   PetscOptionsEnd();
502:   PetscCall(RiemannListFind(rlist, rname, &r));
503:   PetscCall(PetscFunctionListDestroy(&rlist));
504:   ctx->physics.riemann = r;

506:   /* *
507:   * Hack to deal with LxF in semi-discrete form
508:   * max speed is 1 for the basic initial conditions (where |u| <= 1)
509:   * */
510:   if (r == PhysicsRiemann_Burgers_LxF) user->lxf_speed = 1;
511:   PetscFunctionReturn(PETSC_SUCCESS);
512: }

514: /* --------------------------------- Traffic ----------------------------------- */

516: typedef struct {
517:   PetscReal lxf_speed;
518:   PetscReal a;
519: } TrafficCtx;

521: static inline PetscScalar TrafficFlux(PetscScalar a, PetscScalar u)
522: {
523:   return a * u * (1 - u);
524: }

526: static PetscErrorCode PhysicsSample_Traffic(void *vctx, PetscInt initial, FVBCType bctype, PetscReal xmin, PetscReal xmax, PetscReal t, PetscReal x, PetscReal *u)
527: {
528:   PetscReal a = ((TrafficCtx *)vctx)->a;

530:   PetscFunctionBeginUser;
531:   PetscCheck(bctype != FVBC_PERIODIC || t <= 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "Exact solution not implemented for periodic");
532:   switch (initial) {
533:   case 0:
534:     u[0] = (-a * t < x) ? 2 : 0;
535:     break;
536:   case 1:
537:     if (x < PetscMin(2 * a * t, 0.5 + a * t)) u[0] = -1;
538:     else if (x < 1) u[0] = 0;
539:     else u[0] = 1;
540:     break;
541:   case 2:
542:     PetscCheck(t <= 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only initial condition available");
543:     u[0] = 0.7 + 0.3 * PetscSinReal(2 * PETSC_PI * ((x - xmin) / (xmax - xmin)));
544:     break;
545:   default:
546:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition");
547:   }
548:   PetscFunctionReturn(PETSC_SUCCESS);
549: }

551: static PetscErrorCode PhysicsRiemann_Traffic_Exact(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed)
552: {
553:   PetscReal a = ((TrafficCtx *)vctx)->a;

555:   PetscFunctionBeginUser;
556:   if (uL[0] < uR[0]) {
557:     flux[0] = PetscMin(TrafficFlux(a, uL[0]), TrafficFlux(a, uR[0]));
558:   } else {
559:     flux[0] = (uR[0] < 0.5 && 0.5 < uL[0]) ? TrafficFlux(a, 0.5) : PetscMax(TrafficFlux(a, uL[0]), TrafficFlux(a, uR[0]));
560:   }
561:   *maxspeed = a * MaxAbs(1 - 2 * uL[0], 1 - 2 * uR[0]);
562:   PetscFunctionReturn(PETSC_SUCCESS);
563: }

565: static PetscErrorCode PhysicsRiemann_Traffic_Roe(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed)
566: {
567:   PetscReal a = ((TrafficCtx *)vctx)->a;
568:   PetscReal speed;

570:   PetscFunctionBeginUser;
571:   speed     = a * (1 - (uL[0] + uR[0]));
572:   flux[0]   = 0.5 * (TrafficFlux(a, uL[0]) + TrafficFlux(a, uR[0])) - 0.5 * PetscAbs(speed) * (uR[0] - uL[0]);
573:   *maxspeed = speed;
574:   PetscFunctionReturn(PETSC_SUCCESS);
575: }

577: static PetscErrorCode PhysicsRiemann_Traffic_LxF(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed)
578: {
579:   TrafficCtx *phys = (TrafficCtx *)vctx;
580:   PetscReal   a    = phys->a;
581:   PetscReal   speed;

583:   PetscFunctionBeginUser;
584:   speed     = a * (1 - (uL[0] + uR[0]));
585:   flux[0]   = 0.5 * (TrafficFlux(a, uL[0]) + TrafficFlux(a, uR[0])) - 0.5 * phys->lxf_speed * (uR[0] - uL[0]);
586:   *maxspeed = speed;
587:   PetscFunctionReturn(PETSC_SUCCESS);
588: }

590: static PetscErrorCode PhysicsRiemann_Traffic_Rusanov(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed)
591: {
592:   PetscReal a = ((TrafficCtx *)vctx)->a;
593:   PetscReal speed;

595:   PetscFunctionBeginUser;
596:   speed     = a * PetscMax(PetscAbs(1 - 2 * uL[0]), PetscAbs(1 - 2 * uR[0]));
597:   flux[0]   = 0.5 * (TrafficFlux(a, uL[0]) + TrafficFlux(a, uR[0])) - 0.5 * speed * (uR[0] - uL[0]);
598:   *maxspeed = speed;
599:   PetscFunctionReturn(PETSC_SUCCESS);
600: }

602: static PetscErrorCode PhysicsCreate_Traffic(FVCtx *ctx)
603: {
604:   TrafficCtx       *user;
605:   RiemannFunction   r;
606:   PetscFunctionList rlist      = 0;
607:   char              rname[256] = "exact";

609:   PetscFunctionBeginUser;
610:   PetscCall(PetscNew(&user));
611:   ctx->physics.sample         = PhysicsSample_Traffic;
612:   ctx->physics.characteristic = PhysicsCharacteristic_Conservative;
613:   ctx->physics.destroy        = PhysicsDestroy_SimpleFree;
614:   ctx->physics.user           = user;
615:   ctx->physics.dof            = 1;

617:   PetscCall(PetscStrallocpy("density", &ctx->physics.fieldname[0]));
618:   user->a = 0.5;
619:   PetscCall(RiemannListAdd(&rlist, "exact", PhysicsRiemann_Traffic_Exact));
620:   PetscCall(RiemannListAdd(&rlist, "roe", PhysicsRiemann_Traffic_Roe));
621:   PetscCall(RiemannListAdd(&rlist, "lxf", PhysicsRiemann_Traffic_LxF));
622:   PetscCall(RiemannListAdd(&rlist, "rusanov", PhysicsRiemann_Traffic_Rusanov));
623:   PetscOptionsBegin(ctx->comm, ctx->prefix, "Options for Traffic", "");
624:   PetscCall(PetscOptionsReal("-physics_traffic_a", "Flux = a*u*(1-u)", "", user->a, &user->a, NULL));
625:   PetscCall(PetscOptionsFList("-physics_traffic_riemann", "Riemann solver", "", rlist, rname, rname, sizeof(rname), NULL));
626:   PetscOptionsEnd();

628:   PetscCall(RiemannListFind(rlist, rname, &r));
629:   PetscCall(PetscFunctionListDestroy(&rlist));

631:   ctx->physics.riemann = r;

633:   /* *
634:   * Hack to deal with LxF in semi-discrete form
635:   * max speed is 3*a for the basic initial conditions (-1 <= u <= 2)
636:   * */
637:   if (r == PhysicsRiemann_Traffic_LxF) user->lxf_speed = 3 * user->a;
638:   PetscFunctionReturn(PETSC_SUCCESS);
639: }

641: /* --------------------------------- Linear Acoustics ----------------------------------- */

643: /* Flux: u_t + (A u)_x
644:  * z = sqrt(rho*bulk), c = sqrt(rho/bulk)
645:  * Spectral decomposition: A = R * D * Rinv
646:  * [    cz] = [-z   z] [-c    ] [-1/2z  1/2]
647:  * [c/z   ] = [ 1   1] [     c] [ 1/2z  1/2]
648:  *
649:  * We decompose this into the left-traveling waves Al = R * D^- Rinv
650:  * and the right-traveling waves Ar = R * D^+ * Rinv
651:  * Multiplying out these expressions produces the following two matrices
652:  */

654: typedef struct {
655:   PetscReal c; /* speed of sound: c = sqrt(bulk/rho) */
656:   PetscReal z; /* impedance: z = sqrt(rho*bulk) */
657: } AcousticsCtx;

659: PETSC_UNUSED static inline void AcousticsFlux(AcousticsCtx *ctx, const PetscScalar *u, PetscScalar *f)
660: {
661:   f[0] = ctx->c * ctx->z * u[1];
662:   f[1] = ctx->c / ctx->z * u[0];
663: }

665: static PetscErrorCode PhysicsCharacteristic_Acoustics(void *vctx, PetscInt m, const PetscScalar *u, PetscScalar *X, PetscScalar *Xi, PetscReal *speeds)
666: {
667:   AcousticsCtx *phys = (AcousticsCtx *)vctx;
668:   PetscReal     z = phys->z, c = phys->c;

670:   PetscFunctionBeginUser;
671:   X[0 * 2 + 0]  = -z;
672:   X[0 * 2 + 1]  = z;
673:   X[1 * 2 + 0]  = 1;
674:   X[1 * 2 + 1]  = 1;
675:   Xi[0 * 2 + 0] = -1. / (2 * z);
676:   Xi[0 * 2 + 1] = 1. / 2;
677:   Xi[1 * 2 + 0] = 1. / (2 * z);
678:   Xi[1 * 2 + 1] = 1. / 2;
679:   speeds[0]     = -c;
680:   speeds[1]     = c;
681:   PetscFunctionReturn(PETSC_SUCCESS);
682: }

684: static PetscErrorCode PhysicsSample_Acoustics_Initial(AcousticsCtx *phys, PetscInt initial, PetscReal xmin, PetscReal xmax, PetscReal x, PetscReal *u)
685: {
686:   PetscFunctionBeginUser;
687:   switch (initial) {
688:   case 0:
689:     u[0] = (PetscAbs((x - xmin) / (xmax - xmin) - 0.2) < 0.1) ? 1 : 0.5;
690:     u[1] = (PetscAbs((x - xmin) / (xmax - xmin) - 0.7) < 0.1) ? 1 : -0.5;
691:     break;
692:   case 1:
693:     u[0] = PetscCosReal(3 * 2 * PETSC_PI * x / (xmax - xmin));
694:     u[1] = PetscExpReal(-PetscSqr(x - (xmax + xmin) / 2) / (2 * PetscSqr(0.2 * (xmax - xmin)))) - 0.5;
695:     break;
696:   default:
697:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition");
698:   }
699:   PetscFunctionReturn(PETSC_SUCCESS);
700: }

702: static PetscErrorCode PhysicsSample_Acoustics(void *vctx, PetscInt initial, FVBCType bctype, PetscReal xmin, PetscReal xmax, PetscReal t, PetscReal x, PetscReal *u)
703: {
704:   AcousticsCtx *phys = (AcousticsCtx *)vctx;
705:   PetscReal     c    = phys->c;
706:   PetscReal     x0a, x0b, u0a[2], u0b[2], tmp[2];
707:   PetscReal     X[2][2], Xi[2][2], dummy[2];

709:   PetscFunctionBeginUser;
710:   switch (bctype) {
711:   case FVBC_OUTFLOW:
712:     x0a = x + c * t;
713:     x0b = x - c * t;
714:     break;
715:   case FVBC_PERIODIC:
716:     x0a = RangeMod(x + c * t, xmin, xmax);
717:     x0b = RangeMod(x - c * t, xmin, xmax);
718:     break;
719:   default:
720:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown BCType");
721:   }
722:   PetscCall(PhysicsSample_Acoustics_Initial(phys, initial, xmin, xmax, x0a, u0a));
723:   PetscCall(PhysicsSample_Acoustics_Initial(phys, initial, xmin, xmax, x0b, u0b));
724:   PetscCall(PhysicsCharacteristic_Acoustics(vctx, 2, u, &X[0][0], &Xi[0][0], dummy));
725:   tmp[0] = Xi[0][0] * u0a[0] + Xi[0][1] * u0a[1];
726:   tmp[1] = Xi[1][0] * u0b[0] + Xi[1][1] * u0b[1];
727:   u[0]   = X[0][0] * tmp[0] + X[0][1] * tmp[1];
728:   u[1]   = X[1][0] * tmp[0] + X[1][1] * tmp[1];
729:   PetscFunctionReturn(PETSC_SUCCESS);
730: }

732: static PetscErrorCode PhysicsRiemann_Acoustics_Exact(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed)
733: {
734:   AcousticsCtx *phys = (AcousticsCtx *)vctx;
735:   PetscReal     c = phys->c, z = phys->z;
736:   PetscReal     Al[2][2] =
737:     {
738:       {-c / 2,      c * z / 2},
739:       {c / (2 * z), -c / 2   }
740:   },         /* Left traveling waves */
741:     Ar[2][2] = {{c / 2, c * z / 2}, {c / (2 * z), c / 2}}; /* Right traveling waves */

743:   PetscFunctionBeginUser;
744:   flux[0]   = Al[0][0] * uR[0] + Al[0][1] * uR[1] + Ar[0][0] * uL[0] + Ar[0][1] * uL[1];
745:   flux[1]   = Al[1][0] * uR[0] + Al[1][1] * uR[1] + Ar[1][0] * uL[0] + Ar[1][1] * uL[1];
746:   *maxspeed = c;
747:   PetscFunctionReturn(PETSC_SUCCESS);
748: }

750: static PetscErrorCode PhysicsCreate_Acoustics(FVCtx *ctx)
751: {
752:   AcousticsCtx     *user;
753:   PetscFunctionList rlist = 0, rclist = 0;
754:   char              rname[256] = "exact", rcname[256] = "characteristic";

756:   PetscFunctionBeginUser;
757:   PetscCall(PetscNew(&user));
758:   ctx->physics.sample  = PhysicsSample_Acoustics;
759:   ctx->physics.destroy = PhysicsDestroy_SimpleFree;
760:   ctx->physics.user    = user;
761:   ctx->physics.dof     = 2;

763:   PetscCall(PetscStrallocpy("u", &ctx->physics.fieldname[0]));
764:   PetscCall(PetscStrallocpy("v", &ctx->physics.fieldname[1]));

766:   user->c = 1;
767:   user->z = 1;

769:   PetscCall(RiemannListAdd(&rlist, "exact", PhysicsRiemann_Acoustics_Exact));
770:   PetscCall(ReconstructListAdd(&rclist, "characteristic", PhysicsCharacteristic_Acoustics));
771:   PetscCall(ReconstructListAdd(&rclist, "conservative", PhysicsCharacteristic_Conservative));
772:   PetscOptionsBegin(ctx->comm, ctx->prefix, "Options for linear Acoustics", "");
773:   {
774:     PetscCall(PetscOptionsReal("-physics_acoustics_c", "c = sqrt(bulk/rho)", "", user->c, &user->c, NULL));
775:     PetscCall(PetscOptionsReal("-physics_acoustics_z", "z = sqrt(bulk*rho)", "", user->z, &user->z, NULL));
776:     PetscCall(PetscOptionsFList("-physics_acoustics_riemann", "Riemann solver", "", rlist, rname, rname, sizeof(rname), NULL));
777:     PetscCall(PetscOptionsFList("-physics_acoustics_reconstruct", "Reconstruction", "", rclist, rcname, rcname, sizeof(rcname), NULL));
778:   }
779:   PetscOptionsEnd();
780:   PetscCall(RiemannListFind(rlist, rname, &ctx->physics.riemann));
781:   PetscCall(ReconstructListFind(rclist, rcname, &ctx->physics.characteristic));
782:   PetscCall(PetscFunctionListDestroy(&rlist));
783:   PetscCall(PetscFunctionListDestroy(&rclist));
784:   PetscFunctionReturn(PETSC_SUCCESS);
785: }

787: /* --------------------------------- Isothermal Gas Dynamics ----------------------------------- */

789: typedef struct {
790:   PetscReal acoustic_speed;
791: } IsoGasCtx;

793: static inline void IsoGasFlux(PetscReal c, const PetscScalar *u, PetscScalar *f)
794: {
795:   f[0] = u[1];
796:   f[1] = PetscSqr(u[1]) / u[0] + c * c * u[0];
797: }

799: static PetscErrorCode PhysicsSample_IsoGas(void *vctx, PetscInt initial, FVBCType bctype, PetscReal xmin, PetscReal xmax, PetscReal t, PetscReal x, PetscReal *u)
800: {
801:   PetscFunctionBeginUser;
802:   PetscCheck(t <= 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "Exact solutions not implemented for t > 0");
803:   switch (initial) {
804:   case 0:
805:     u[0] = (x < 0) ? 1 : 0.5;
806:     u[1] = (x < 0) ? 1 : 0.7;
807:     break;
808:   case 1:
809:     u[0] = 1 + 0.5 * PetscSinReal(2 * PETSC_PI * x);
810:     u[1] = 1 * u[0];
811:     break;
812:   default:
813:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "unknown initial condition");
814:   }
815:   PetscFunctionReturn(PETSC_SUCCESS);
816: }

818: static PetscErrorCode PhysicsRiemann_IsoGas_Roe(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed)
819: {
820:   IsoGasCtx  *phys = (IsoGasCtx *)vctx;
821:   PetscReal   c    = phys->acoustic_speed;
822:   PetscScalar ubar, du[2], a[2], fL[2], fR[2], lam[2], ustar[2], R[2][2];
823:   PetscInt    i;

825:   PetscFunctionBeginUser;
826:   ubar = (uL[1] / PetscSqrtScalar(uL[0]) + uR[1] / PetscSqrtScalar(uR[0])) / (PetscSqrtScalar(uL[0]) + PetscSqrtScalar(uR[0]));
827:   /* write fluxuations in characteristic basis */
828:   du[0] = uR[0] - uL[0];
829:   du[1] = uR[1] - uL[1];
830:   a[0]  = (1 / (2 * c)) * ((ubar + c) * du[0] - du[1]);
831:   a[1]  = (1 / (2 * c)) * ((-ubar + c) * du[0] + du[1]);
832:   /* wave speeds */
833:   lam[0] = ubar - c;
834:   lam[1] = ubar + c;
835:   /* Right eigenvectors */
836:   R[0][0] = 1;
837:   R[0][1] = ubar - c;
838:   R[1][0] = 1;
839:   R[1][1] = ubar + c;
840:   /* Compute state in star region (between the 1-wave and 2-wave) */
841:   for (i = 0; i < 2; i++) ustar[i] = uL[i] + a[0] * R[0][i];
842:   if (uL[1] / uL[0] < c && c < ustar[1] / ustar[0]) { /* 1-wave is sonic rarefaction */
843:     PetscScalar ufan[2];
844:     ufan[0] = uL[0] * PetscExpScalar(uL[1] / (uL[0] * c) - 1);
845:     ufan[1] = c * ufan[0];
846:     IsoGasFlux(c, ufan, flux);
847:   } else if (ustar[1] / ustar[0] < -c && -c < uR[1] / uR[0]) { /* 2-wave is sonic rarefaction */
848:     PetscScalar ufan[2];
849:     ufan[0] = uR[0] * PetscExpScalar(-uR[1] / (uR[0] * c) - 1);
850:     ufan[1] = -c * ufan[0];
851:     IsoGasFlux(c, ufan, flux);
852:   } else { /* Centered form */
853:     IsoGasFlux(c, uL, fL);
854:     IsoGasFlux(c, uR, fR);
855:     for (i = 0; i < 2; i++) {
856:       PetscScalar absdu = PetscAbsScalar(lam[0]) * a[0] * R[0][i] + PetscAbsScalar(lam[1]) * a[1] * R[1][i];
857:       flux[i]           = 0.5 * (fL[i] + fR[i]) - 0.5 * absdu;
858:     }
859:   }
860:   *maxspeed = MaxAbs(lam[0], lam[1]);
861:   PetscFunctionReturn(PETSC_SUCCESS);
862: }

864: static PetscErrorCode PhysicsRiemann_IsoGas_Exact(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed)
865: {
866:   IsoGasCtx  *phys = (IsoGasCtx *)vctx;
867:   PetscReal   c    = phys->acoustic_speed;
868:   PetscScalar ustar[2];
869:   struct {
870:     PetscScalar rho, u;
871:   } L = {uL[0], uL[1] / uL[0]}, R = {uR[0], uR[1] / uR[0]}, star;
872:   PetscInt i;

874:   PetscFunctionBeginUser;
875:   PetscCheck(L.rho > 0 && R.rho > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Reconstructed density is negative");
876:   {
877:     /* Solve for star state */
878:     PetscScalar res, tmp, rho = 0.5 * (L.rho + R.rho); /* initial guess */
879:     for (i = 0; i < 20; i++) {
880:       PetscScalar fr, fl, dfr, dfl;
881:       fl  = (L.rho < rho) ? (rho - L.rho) / PetscSqrtScalar(L.rho * rho) /* shock */
882:                           : PetscLogScalar(rho) - PetscLogScalar(L.rho); /* rarefaction */
883:       fr  = (R.rho < rho) ? (rho - R.rho) / PetscSqrtScalar(R.rho * rho) /* shock */
884:                           : PetscLogScalar(rho) - PetscLogScalar(R.rho); /* rarefaction */
885:       res = R.u - L.u + c * (fr + fl);
886:       PetscCheck(!PetscIsInfOrNanScalar(res), PETSC_COMM_SELF, PETSC_ERR_FP, "Infinity or Not-a-Number generated in computation");
887:       if (PetscAbsScalar(res) < 1e-10) {
888:         star.rho = rho;
889:         star.u   = L.u - c * fl;
890:         goto converged;
891:       }
892:       dfl = (L.rho < rho) ? 1 / PetscSqrtScalar(L.rho * rho) * (1 - 0.5 * (rho - L.rho) / rho) : 1 / rho;
893:       dfr = (R.rho < rho) ? 1 / PetscSqrtScalar(R.rho * rho) * (1 - 0.5 * (rho - R.rho) / rho) : 1 / rho;
894:       tmp = rho - res / (c * (dfr + dfl));
895:       if (tmp <= 0) rho /= 2; /* Guard against Newton shooting off to a negative density */
896:       else rho = tmp;
897:       PetscCheck((rho > 0) && PetscIsNormalScalar(rho), PETSC_COMM_SELF, PETSC_ERR_FP, "non-normal iterate rho=%g", (double)PetscRealPart(rho));
898:     }
899:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_CONV_FAILED, "Newton iteration for star.rho diverged after %" PetscInt_FMT " iterations", i);
900:   }
901: converged:
902:   if (L.u - c < 0 && 0 < star.u - c) { /* 1-wave is sonic rarefaction */
903:     PetscScalar ufan[2];
904:     ufan[0] = L.rho * PetscExpScalar(L.u / c - 1);
905:     ufan[1] = c * ufan[0];
906:     IsoGasFlux(c, ufan, flux);
907:   } else if (star.u + c < 0 && 0 < R.u + c) { /* 2-wave is sonic rarefaction */
908:     PetscScalar ufan[2];
909:     ufan[0] = R.rho * PetscExpScalar(-R.u / c - 1);
910:     ufan[1] = -c * ufan[0];
911:     IsoGasFlux(c, ufan, flux);
912:   } else if ((L.rho >= star.rho && L.u - c >= 0) || (L.rho < star.rho && (star.rho * star.u - L.rho * L.u) / (star.rho - L.rho) > 0)) {
913:     /* 1-wave is supersonic rarefaction, or supersonic shock */
914:     IsoGasFlux(c, uL, flux);
915:   } else if ((star.rho <= R.rho && R.u + c <= 0) || (star.rho > R.rho && (R.rho * R.u - star.rho * star.u) / (R.rho - star.rho) < 0)) {
916:     /* 2-wave is supersonic rarefaction or supersonic shock */
917:     IsoGasFlux(c, uR, flux);
918:   } else {
919:     ustar[0] = star.rho;
920:     ustar[1] = star.rho * star.u;
921:     IsoGasFlux(c, ustar, flux);
922:   }
923:   *maxspeed = MaxAbs(MaxAbs(star.u - c, star.u + c), MaxAbs(L.u - c, R.u + c));
924:   PetscFunctionReturn(PETSC_SUCCESS);
925: }

927: static PetscErrorCode PhysicsRiemann_IsoGas_Rusanov(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed)
928: {
929:   IsoGasCtx  *phys = (IsoGasCtx *)vctx;
930:   PetscScalar c    = phys->acoustic_speed, fL[2], fR[2], s;
931:   struct {
932:     PetscScalar rho, u;
933:   } L = {uL[0], uL[1] / uL[0]}, R = {uR[0], uR[1] / uR[0]};

935:   PetscFunctionBeginUser;
936:   PetscCheck(L.rho > 0 && R.rho > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Reconstructed density is negative");
937:   IsoGasFlux(c, uL, fL);
938:   IsoGasFlux(c, uR, fR);
939:   s         = PetscMax(PetscAbs(L.u), PetscAbs(R.u)) + c;
940:   flux[0]   = 0.5 * (fL[0] + fR[0]) + 0.5 * s * (uL[0] - uR[0]);
941:   flux[1]   = 0.5 * (fL[1] + fR[1]) + 0.5 * s * (uL[1] - uR[1]);
942:   *maxspeed = s;
943:   PetscFunctionReturn(PETSC_SUCCESS);
944: }

946: static PetscErrorCode PhysicsCharacteristic_IsoGas(void *vctx, PetscInt m, const PetscScalar *u, PetscScalar *X, PetscScalar *Xi, PetscReal *speeds)
947: {
948:   IsoGasCtx *phys = (IsoGasCtx *)vctx;
949:   PetscReal  c    = phys->acoustic_speed;

951:   PetscFunctionBeginUser;
952:   speeds[0]    = u[1] / u[0] - c;
953:   speeds[1]    = u[1] / u[0] + c;
954:   X[0 * 2 + 0] = 1;
955:   X[0 * 2 + 1] = speeds[0];
956:   X[1 * 2 + 0] = 1;
957:   X[1 * 2 + 1] = speeds[1];
958:   PetscCall(PetscArraycpy(Xi, X, 4));
959:   PetscCall(PetscKernel_A_gets_inverse_A_2(Xi, 0, PETSC_FALSE, NULL));
960:   PetscFunctionReturn(PETSC_SUCCESS);
961: }

963: static PetscErrorCode PhysicsCreate_IsoGas(FVCtx *ctx)
964: {
965:   IsoGasCtx        *user;
966:   PetscFunctionList rlist = 0, rclist = 0;
967:   char              rname[256] = "exact", rcname[256] = "characteristic";

969:   PetscFunctionBeginUser;
970:   PetscCall(PetscNew(&user));
971:   ctx->physics.sample  = PhysicsSample_IsoGas;
972:   ctx->physics.destroy = PhysicsDestroy_SimpleFree;
973:   ctx->physics.user    = user;
974:   ctx->physics.dof     = 2;

976:   PetscCall(PetscStrallocpy("density", &ctx->physics.fieldname[0]));
977:   PetscCall(PetscStrallocpy("momentum", &ctx->physics.fieldname[1]));

979:   user->acoustic_speed = 1;

981:   PetscCall(RiemannListAdd(&rlist, "exact", PhysicsRiemann_IsoGas_Exact));
982:   PetscCall(RiemannListAdd(&rlist, "roe", PhysicsRiemann_IsoGas_Roe));
983:   PetscCall(RiemannListAdd(&rlist, "rusanov", PhysicsRiemann_IsoGas_Rusanov));
984:   PetscCall(ReconstructListAdd(&rclist, "characteristic", PhysicsCharacteristic_IsoGas));
985:   PetscCall(ReconstructListAdd(&rclist, "conservative", PhysicsCharacteristic_Conservative));
986:   PetscOptionsBegin(ctx->comm, ctx->prefix, "Options for IsoGas", "");
987:   PetscCall(PetscOptionsReal("-physics_isogas_acoustic_speed", "Acoustic speed", "", user->acoustic_speed, &user->acoustic_speed, NULL));
988:   PetscCall(PetscOptionsFList("-physics_isogas_riemann", "Riemann solver", "", rlist, rname, rname, sizeof(rname), NULL));
989:   PetscCall(PetscOptionsFList("-physics_isogas_reconstruct", "Reconstruction", "", rclist, rcname, rcname, sizeof(rcname), NULL));
990:   PetscOptionsEnd();
991:   PetscCall(RiemannListFind(rlist, rname, &ctx->physics.riemann));
992:   PetscCall(ReconstructListFind(rclist, rcname, &ctx->physics.characteristic));
993:   PetscCall(PetscFunctionListDestroy(&rlist));
994:   PetscCall(PetscFunctionListDestroy(&rclist));
995:   PetscFunctionReturn(PETSC_SUCCESS);
996: }

998: /* --------------------------------- Shallow Water ----------------------------------- */
999: typedef struct {
1000:   PetscReal gravity;
1001: } ShallowCtx;

1003: static inline void ShallowFlux(ShallowCtx *phys, const PetscScalar *u, PetscScalar *f)
1004: {
1005:   f[0] = u[1];
1006:   f[1] = PetscSqr(u[1]) / u[0] + 0.5 * phys->gravity * PetscSqr(u[0]);
1007: }

1009: static PetscErrorCode PhysicsRiemann_Shallow_Exact(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed)
1010: {
1011:   ShallowCtx *phys = (ShallowCtx *)vctx;
1012:   PetscScalar g    = phys->gravity, ustar[2], cL, cR, c, cstar;
1013:   struct {
1014:     PetscScalar h, u;
1015:   } L = {uL[0], uL[1] / uL[0]}, R = {uR[0], uR[1] / uR[0]}, star;
1016:   PetscInt i;

1018:   PetscFunctionBeginUser;
1019:   PetscCheck(L.h > 0 && R.h > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Reconstructed thickness is negative");
1020:   cL = PetscSqrtScalar(g * L.h);
1021:   cR = PetscSqrtScalar(g * R.h);
1022:   c  = PetscMax(cL, cR);
1023:   {
1024:     /* Solve for star state */
1025:     const PetscInt maxits = 50;
1026:     PetscScalar    tmp, res, res0 = 0, h0, h = 0.5 * (L.h + R.h); /* initial guess */
1027:     h0 = h;
1028:     for (i = 0; i < maxits; i++) {
1029:       PetscScalar fr, fl, dfr, dfl;
1030:       fl  = (L.h < h) ? PetscSqrtScalar(0.5 * g * (h * h - L.h * L.h) * (1 / L.h - 1 / h)) /* shock */
1031:                       : 2 * PetscSqrtScalar(g * h) - 2 * PetscSqrtScalar(g * L.h);         /* rarefaction */
1032:       fr  = (R.h < h) ? PetscSqrtScalar(0.5 * g * (h * h - R.h * R.h) * (1 / R.h - 1 / h)) /* shock */
1033:                       : 2 * PetscSqrtScalar(g * h) - 2 * PetscSqrtScalar(g * R.h);         /* rarefaction */
1034:       res = R.u - L.u + fr + fl;
1035:       PetscCheck(!PetscIsInfOrNanScalar(res), PETSC_COMM_SELF, PETSC_ERR_FP, "Infinity or Not-a-Number generated in computation");
1036:       if (PetscAbsScalar(res) < 1e-8 || (i > 0 && PetscAbsScalar(h - h0) < 1e-8)) {
1037:         star.h = h;
1038:         star.u = L.u - fl;
1039:         goto converged;
1040:       } else if (i > 0 && PetscAbsScalar(res) >= PetscAbsScalar(res0)) { /* Line search */
1041:         h = 0.8 * h0 + 0.2 * h;
1042:         continue;
1043:       }
1044:       /* Accept the last step and take another */
1045:       res0 = res;
1046:       h0   = h;
1047:       dfl  = (L.h < h) ? 0.5 / fl * 0.5 * g * (-L.h * L.h / (h * h) - 1 + 2 * h / L.h) : PetscSqrtScalar(g / h);
1048:       dfr  = (R.h < h) ? 0.5 / fr * 0.5 * g * (-R.h * R.h / (h * h) - 1 + 2 * h / R.h) : PetscSqrtScalar(g / h);
1049:       tmp  = h - res / (dfr + dfl);
1050:       if (tmp <= 0) h /= 2; /* Guard against Newton shooting off to a negative thickness */
1051:       else h = tmp;
1052:       PetscCheck((h > 0) && PetscIsNormalScalar(h), PETSC_COMM_SELF, PETSC_ERR_FP, "non-normal iterate h=%g", (double)h);
1053:     }
1054:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_CONV_FAILED, "Newton iteration for star.h diverged after %" PetscInt_FMT " iterations", i);
1055:   }
1056: converged:
1057:   cstar = PetscSqrtScalar(g * star.h);
1058:   if (L.u - cL < 0 && 0 < star.u - cstar) { /* 1-wave is sonic rarefaction */
1059:     PetscScalar ufan[2];
1060:     ufan[0] = 1 / g * PetscSqr(L.u / 3 + 2. / 3 * cL);
1061:     ufan[1] = PetscSqrtScalar(g * ufan[0]) * ufan[0];
1062:     ShallowFlux(phys, ufan, flux);
1063:   } else if (star.u + cstar < 0 && 0 < R.u + cR) { /* 2-wave is sonic rarefaction */
1064:     PetscScalar ufan[2];
1065:     ufan[0] = 1 / g * PetscSqr(R.u / 3 - 2. / 3 * cR);
1066:     ufan[1] = -PetscSqrtScalar(g * ufan[0]) * ufan[0];
1067:     ShallowFlux(phys, ufan, flux);
1068:   } else if ((L.h >= star.h && L.u - c >= 0) || (L.h < star.h && (star.h * star.u - L.h * L.u) / (star.h - L.h) > 0)) {
1069:     /* 1-wave is right-travelling shock (supersonic) */
1070:     ShallowFlux(phys, uL, flux);
1071:   } else if ((star.h <= R.h && R.u + c <= 0) || (star.h > R.h && (R.h * R.u - star.h * star.h) / (R.h - star.h) < 0)) {
1072:     /* 2-wave is left-travelling shock (supersonic) */
1073:     ShallowFlux(phys, uR, flux);
1074:   } else {
1075:     ustar[0] = star.h;
1076:     ustar[1] = star.h * star.u;
1077:     ShallowFlux(phys, ustar, flux);
1078:   }
1079:   *maxspeed = MaxAbs(MaxAbs(star.u - cstar, star.u + cstar), MaxAbs(L.u - cL, R.u + cR));
1080:   PetscFunctionReturn(PETSC_SUCCESS);
1081: }

1083: static PetscErrorCode PhysicsRiemann_Shallow_Rusanov(void *vctx, PetscInt m, const PetscScalar *uL, const PetscScalar *uR, PetscScalar *flux, PetscReal *maxspeed)
1084: {
1085:   ShallowCtx *phys = (ShallowCtx *)vctx;
1086:   PetscScalar g    = phys->gravity, fL[2], fR[2], s;
1087:   struct {
1088:     PetscScalar h, u;
1089:   } L = {uL[0], uL[1] / uL[0]}, R = {uR[0], uR[1] / uR[0]};

1091:   PetscFunctionBeginUser;
1092:   PetscCheck(L.h > 0 && R.h > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Reconstructed thickness is negative");
1093:   ShallowFlux(phys, uL, fL);
1094:   ShallowFlux(phys, uR, fR);
1095:   s         = PetscMax(PetscAbs(L.u) + PetscSqrtScalar(g * L.h), PetscAbs(R.u) + PetscSqrtScalar(g * R.h));
1096:   flux[0]   = 0.5 * (fL[0] + fR[0]) + 0.5 * s * (uL[0] - uR[0]);
1097:   flux[1]   = 0.5 * (fL[1] + fR[1]) + 0.5 * s * (uL[1] - uR[1]);
1098:   *maxspeed = s;
1099:   PetscFunctionReturn(PETSC_SUCCESS);
1100: }

1102: static PetscErrorCode PhysicsCharacteristic_Shallow(void *vctx, PetscInt m, const PetscScalar *u, PetscScalar *X, PetscScalar *Xi, PetscReal *speeds)
1103: {
1104:   ShallowCtx *phys = (ShallowCtx *)vctx;
1105:   PetscReal   c;

1107:   PetscFunctionBeginUser;
1108:   c            = PetscSqrtScalar(u[0] * phys->gravity);
1109:   speeds[0]    = u[1] / u[0] - c;
1110:   speeds[1]    = u[1] / u[0] + c;
1111:   X[0 * 2 + 0] = 1;
1112:   X[0 * 2 + 1] = speeds[0];
1113:   X[1 * 2 + 0] = 1;
1114:   X[1 * 2 + 1] = speeds[1];
1115:   PetscCall(PetscArraycpy(Xi, X, 4));
1116:   PetscCall(PetscKernel_A_gets_inverse_A_2(Xi, 0, PETSC_FALSE, NULL));
1117:   PetscFunctionReturn(PETSC_SUCCESS);
1118: }

1120: static PetscErrorCode PhysicsCreate_Shallow(FVCtx *ctx)
1121: {
1122:   ShallowCtx       *user;
1123:   PetscFunctionList rlist = 0, rclist = 0;
1124:   char              rname[256] = "exact", rcname[256] = "characteristic";

1126:   PetscFunctionBeginUser;
1127:   PetscCall(PetscNew(&user));
1128:   /* Shallow water and Isothermal Gas dynamics are similar so we reuse initial conditions for now */
1129:   ctx->physics.sample  = PhysicsSample_IsoGas;
1130:   ctx->physics.destroy = PhysicsDestroy_SimpleFree;
1131:   ctx->physics.user    = user;
1132:   ctx->physics.dof     = 2;

1134:   PetscCall(PetscStrallocpy("density", &ctx->physics.fieldname[0]));
1135:   PetscCall(PetscStrallocpy("momentum", &ctx->physics.fieldname[1]));

1137:   user->gravity = 1;

1139:   PetscCall(RiemannListAdd(&rlist, "exact", PhysicsRiemann_Shallow_Exact));
1140:   PetscCall(RiemannListAdd(&rlist, "rusanov", PhysicsRiemann_Shallow_Rusanov));
1141:   PetscCall(ReconstructListAdd(&rclist, "characteristic", PhysicsCharacteristic_Shallow));
1142:   PetscCall(ReconstructListAdd(&rclist, "conservative", PhysicsCharacteristic_Conservative));
1143:   PetscOptionsBegin(ctx->comm, ctx->prefix, "Options for Shallow", "");
1144:   PetscCall(PetscOptionsReal("-physics_shallow_gravity", "Gravity", "", user->gravity, &user->gravity, NULL));
1145:   PetscCall(PetscOptionsFList("-physics_shallow_riemann", "Riemann solver", "", rlist, rname, rname, sizeof(rname), NULL));
1146:   PetscCall(PetscOptionsFList("-physics_shallow_reconstruct", "Reconstruction", "", rclist, rcname, rcname, sizeof(rcname), NULL));
1147:   PetscOptionsEnd();
1148:   PetscCall(RiemannListFind(rlist, rname, &ctx->physics.riemann));
1149:   PetscCall(ReconstructListFind(rclist, rcname, &ctx->physics.characteristic));
1150:   PetscCall(PetscFunctionListDestroy(&rlist));
1151:   PetscCall(PetscFunctionListDestroy(&rclist));
1152:   PetscFunctionReturn(PETSC_SUCCESS);
1153: }

1155: /* --------------------------------- Finite Volume Solver ----------------------------------- */

1157: static PetscErrorCode FVRHSFunction(TS ts, PetscReal time, Vec X, Vec F, void *vctx)
1158: {
1159:   FVCtx       *ctx = (FVCtx *)vctx;
1160:   PetscInt     i, j, k, Mx, dof, xs, xm;
1161:   PetscReal    hx, cfl_idt = 0;
1162:   PetscScalar *x, *f, *slope;
1163:   Vec          Xloc;
1164:   DM           da;

1166:   PetscFunctionBeginUser;
1167:   PetscCall(TSGetDM(ts, &da));
1168:   PetscCall(DMGetLocalVector(da, &Xloc));
1169:   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
1170:   hx = (ctx->xmax - ctx->xmin) / Mx;
1171:   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
1172:   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));

1174:   PetscCall(VecZeroEntries(F));

1176:   PetscCall(DMDAVecGetArray(da, Xloc, &x));
1177:   PetscCall(DMDAVecGetArray(da, F, &f));
1178:   PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope));

1180:   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));

1182:   if (ctx->bctype == FVBC_OUTFLOW) {
1183:     for (i = xs - 2; i < 0; i++) {
1184:       for (j = 0; j < dof; j++) x[i * dof + j] = x[j];
1185:     }
1186:     for (i = Mx; i < xs + xm + 2; i++) {
1187:       for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j];
1188:     }
1189:   }
1190:   for (i = xs - 1; i < xs + xm + 1; i++) {
1191:     struct _LimitInfo info;
1192:     PetscScalar      *cjmpL, *cjmpR;
1193:     /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
1194:     PetscCall((*ctx->physics.characteristic)(ctx->physics.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
1195:     /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
1196:     PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof));
1197:     cjmpL = &ctx->cjmpLR[0];
1198:     cjmpR = &ctx->cjmpLR[dof];
1199:     for (j = 0; j < dof; j++) {
1200:       PetscScalar jmpL, jmpR;
1201:       jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j];
1202:       jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j];
1203:       for (k = 0; k < dof; k++) {
1204:         cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL;
1205:         cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR;
1206:       }
1207:     }
1208:     /* Apply limiter to the left and right characteristic jumps */
1209:     info.m  = dof;
1210:     info.hx = hx;
1211:     (*ctx->limit)(&info, cjmpL, cjmpR, ctx->cslope);
1212:     for (j = 0; j < dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */
1213:     for (j = 0; j < dof; j++) {
1214:       PetscScalar tmp = 0;
1215:       for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k];
1216:       slope[i * dof + j] = tmp;
1217:     }
1218:   }

1220:   for (i = xs; i < xs + xm + 1; i++) {
1221:     PetscReal    maxspeed;
1222:     PetscScalar *uL, *uR;
1223:     uL = &ctx->uLR[0];
1224:     uR = &ctx->uLR[dof];
1225:     for (j = 0; j < dof; j++) {
1226:       uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2;
1227:       uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2;
1228:     }
1229:     PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed));
1230:     cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */

1232:     if (i > xs) {
1233:       for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hx;
1234:     }
1235:     if (i < xs + xm) {
1236:       for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hx;
1237:     }
1238:   }

1240:   PetscCall(DMDAVecRestoreArray(da, Xloc, &x));
1241:   PetscCall(DMDAVecRestoreArray(da, F, &f));
1242:   PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope));
1243:   PetscCall(DMRestoreLocalVector(da, &Xloc));

1245:   PetscCallMPI(MPIU_Allreduce(&cfl_idt, &ctx->cfl_idt, 1, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)da)));
1246:   if (0) {
1247:     /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */
1248:     PetscReal dt, tnow;
1249:     PetscCall(TSGetTimeStep(ts, &dt));
1250:     PetscCall(TSGetTime(ts, &tnow));
1251:     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)));
1252:   }
1253:   PetscFunctionReturn(PETSC_SUCCESS);
1254: }

1256: static PetscErrorCode SmallMatMultADB(PetscScalar *C, PetscInt bs, const PetscScalar *A, const PetscReal *D, const PetscScalar *B)
1257: {
1258:   PetscInt i, j, k;

1260:   PetscFunctionBeginUser;
1261:   for (i = 0; i < bs; i++) {
1262:     for (j = 0; j < bs; j++) {
1263:       PetscScalar tmp = 0;
1264:       for (k = 0; k < bs; k++) tmp += A[i * bs + k] * D[k] * B[k * bs + j];
1265:       C[i * bs + j] = tmp;
1266:     }
1267:   }
1268:   PetscFunctionReturn(PETSC_SUCCESS);
1269: }

1271: static PetscErrorCode FVIJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal shift, Mat A, Mat B, void *vctx)
1272: {
1273:   FVCtx             *ctx = (FVCtx *)vctx;
1274:   PetscInt           i, j, dof = ctx->physics.dof;
1275:   PetscScalar       *J;
1276:   const PetscScalar *x;
1277:   PetscReal          hx;
1278:   DM                 da;
1279:   DMDALocalInfo      dainfo;

1281:   PetscFunctionBeginUser;
1282:   PetscCall(TSGetDM(ts, &da));
1283:   PetscCall(DMDAVecGetArrayRead(da, X, (void *)&x));
1284:   PetscCall(DMDAGetLocalInfo(da, &dainfo));
1285:   hx = (ctx->xmax - ctx->xmin) / dainfo.mx;
1286:   PetscCall(PetscMalloc1(dof * dof, &J));
1287:   for (i = dainfo.xs; i < dainfo.xs + dainfo.xm; i++) {
1288:     PetscCall((*ctx->physics.characteristic)(ctx->physics.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds));
1289:     for (j = 0; j < dof; j++) ctx->speeds[j] = PetscAbs(ctx->speeds[j]);
1290:     PetscCall(SmallMatMultADB(J, dof, ctx->R, ctx->speeds, ctx->Rinv));
1291:     for (j = 0; j < dof * dof; j++) J[j] = J[j] / hx + shift * (j / dof == j % dof);
1292:     PetscCall(MatSetValuesBlocked(B, 1, &i, 1, &i, J, INSERT_VALUES));
1293:   }
1294:   PetscCall(PetscFree(J));
1295:   PetscCall(DMDAVecRestoreArrayRead(da, X, (void *)&x));

1297:   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1298:   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
1299:   if (A != B) {
1300:     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1301:     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1302:   }
1303:   PetscFunctionReturn(PETSC_SUCCESS);
1304: }

1306: static PetscErrorCode FVSample(FVCtx *ctx, DM da, PetscReal time, Vec U)
1307: {
1308:   PetscScalar *u, *uj;
1309:   PetscInt     i, j, k, dof, xs, xm, Mx;

1311:   PetscFunctionBeginUser;
1312:   PetscCheck(ctx->physics.sample, PETSC_COMM_SELF, PETSC_ERR_SUP, "Physics has not provided a sampling function");
1313:   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
1314:   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
1315:   PetscCall(DMDAVecGetArray(da, U, &u));
1316:   PetscCall(PetscMalloc1(dof, &uj));
1317:   for (i = xs; i < xs + xm; i++) {
1318:     const PetscReal h = (ctx->xmax - ctx->xmin) / Mx, xi = ctx->xmin + h / 2 + i * h;
1319:     const PetscInt  N = 200;
1320:     /* Integrate over cell i using trapezoid rule with N points. */
1321:     for (k = 0; k < dof; k++) u[i * dof + k] = 0;
1322:     for (j = 0; j < N + 1; j++) {
1323:       PetscScalar xj = xi + h * (j - N / 2) / (PetscReal)N;
1324:       PetscCall((*ctx->physics.sample)(ctx->physics.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj));
1325:       for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N;
1326:     }
1327:   }
1328:   PetscCall(DMDAVecRestoreArray(da, U, &u));
1329:   PetscCall(PetscFree(uj));
1330:   PetscFunctionReturn(PETSC_SUCCESS);
1331: }

1333: static PetscErrorCode SolutionStatsView(DM da, Vec X, PetscViewer viewer)
1334: {
1335:   PetscReal          xmin, xmax;
1336:   PetscScalar        sum, tvsum, tvgsum;
1337:   const PetscScalar *x;
1338:   PetscInt           imin, imax, Mx, i, j, xs, xm, dof;
1339:   Vec                Xloc;
1340:   PetscBool          isascii;

1342:   PetscFunctionBeginUser;
1343:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1344:   PetscCheck(isascii, PETSC_COMM_SELF, PETSC_ERR_SUP, "Viewer type not supported");
1345:   /* PETSc lacks a function to compute total variation norm (difficult in multiple dimensions), we do it here */
1346:   PetscCall(DMGetLocalVector(da, &Xloc));
1347:   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
1348:   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
1349:   PetscCall(DMDAVecGetArrayRead(da, Xloc, (void *)&x));
1350:   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));
1351:   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
1352:   tvsum = 0;
1353:   for (i = xs; i < xs + xm; i++) {
1354:     for (j = 0; j < dof; j++) tvsum += PetscAbsScalar(x[i * dof + j] - x[(i - 1) * dof + j]);
1355:   }
1356:   PetscCallMPI(MPIU_Allreduce(&tvsum, &tvgsum, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)da)));
1357:   PetscCall(DMDAVecRestoreArrayRead(da, Xloc, (void *)&x));
1358:   PetscCall(DMRestoreLocalVector(da, &Xloc));

1360:   PetscCall(VecMin(X, &imin, &xmin));
1361:   PetscCall(VecMax(X, &imax, &xmax));
1362:   PetscCall(VecSum(X, &sum));
1363:   PetscCall(PetscViewerASCIIPrintf(viewer, "Solution range [%8.5f,%8.5f] with extrema at %" PetscInt_FMT " and %" PetscInt_FMT ", mean %8.5f, ||x||_TV %8.5f\n", (double)xmin, (double)xmax, imin, imax, (double)(sum / Mx), (double)(tvgsum / Mx)));
1364:   PetscFunctionReturn(PETSC_SUCCESS);
1365: }

1367: static PetscErrorCode SolutionErrorNorms(FVCtx *ctx, DM da, PetscReal t, Vec X, PetscReal *nrm1, PetscReal *nrmsup)
1368: {
1369:   Vec      Y;
1370:   PetscInt Mx;

1372:   PetscFunctionBeginUser;
1373:   PetscCall(VecGetSize(X, &Mx));
1374:   PetscCall(VecDuplicate(X, &Y));
1375:   PetscCall(FVSample(ctx, da, t, Y));
1376:   PetscCall(VecAYPX(Y, -1, X));
1377:   PetscCall(VecNorm(Y, NORM_1, nrm1));
1378:   PetscCall(VecNorm(Y, NORM_INFINITY, nrmsup));
1379:   *nrm1 /= Mx;
1380:   PetscCall(VecDestroy(&Y));
1381:   PetscFunctionReturn(PETSC_SUCCESS);
1382: }

1384: int main(int argc, char *argv[])
1385: {
1386:   char              lname[256] = "mc", physname[256] = "advect", final_fname[256] = "solution.m";
1387:   PetscFunctionList limiters = 0, physics = 0;
1388:   MPI_Comm          comm;
1389:   TS                ts;
1390:   DM                da;
1391:   Vec               X, X0, R;
1392:   Mat               B;
1393:   FVCtx             ctx;
1394:   PetscInt          i, dof, xs, xm, Mx, draw = 0;
1395:   PetscBool         view_final = PETSC_FALSE;
1396:   PetscReal         ptime;

1398:   PetscFunctionBeginUser;
1399:   PetscCall(PetscInitialize(&argc, &argv, 0, help));
1400:   comm = PETSC_COMM_WORLD;
1401:   PetscCall(PetscMemzero(&ctx, sizeof(ctx)));

1403:   /* Register limiters to be available on the command line */
1404:   PetscCall(PetscFunctionListAdd(&limiters, "upwind", Limit_Upwind));
1405:   PetscCall(PetscFunctionListAdd(&limiters, "lax-wendroff", Limit_LaxWendroff));
1406:   PetscCall(PetscFunctionListAdd(&limiters, "beam-warming", Limit_BeamWarming));
1407:   PetscCall(PetscFunctionListAdd(&limiters, "fromm", Limit_Fromm));
1408:   PetscCall(PetscFunctionListAdd(&limiters, "minmod", Limit_Minmod));
1409:   PetscCall(PetscFunctionListAdd(&limiters, "superbee", Limit_Superbee));
1410:   PetscCall(PetscFunctionListAdd(&limiters, "mc", Limit_MC));
1411:   PetscCall(PetscFunctionListAdd(&limiters, "vanleer", Limit_VanLeer));
1412:   PetscCall(PetscFunctionListAdd(&limiters, "vanalbada", Limit_VanAlbada));
1413:   PetscCall(PetscFunctionListAdd(&limiters, "vanalbadatvd", Limit_VanAlbadaTVD));
1414:   PetscCall(PetscFunctionListAdd(&limiters, "koren", Limit_Koren));
1415:   PetscCall(PetscFunctionListAdd(&limiters, "korensym", Limit_KorenSym));
1416:   PetscCall(PetscFunctionListAdd(&limiters, "koren3", Limit_Koren3));
1417:   PetscCall(PetscFunctionListAdd(&limiters, "cada-torrilhon2", Limit_CadaTorrilhon2));
1418:   PetscCall(PetscFunctionListAdd(&limiters, "cada-torrilhon3-r0p1", Limit_CadaTorrilhon3R0p1));
1419:   PetscCall(PetscFunctionListAdd(&limiters, "cada-torrilhon3-r1", Limit_CadaTorrilhon3R1));
1420:   PetscCall(PetscFunctionListAdd(&limiters, "cada-torrilhon3-r10", Limit_CadaTorrilhon3R10));
1421:   PetscCall(PetscFunctionListAdd(&limiters, "cada-torrilhon3-r100", Limit_CadaTorrilhon3R100));

1423:   /* Register physical models to be available on the command line */
1424:   PetscCall(PetscFunctionListAdd(&physics, "advect", PhysicsCreate_Advect));
1425:   PetscCall(PetscFunctionListAdd(&physics, "burgers", PhysicsCreate_Burgers));
1426:   PetscCall(PetscFunctionListAdd(&physics, "traffic", PhysicsCreate_Traffic));
1427:   PetscCall(PetscFunctionListAdd(&physics, "acoustics", PhysicsCreate_Acoustics));
1428:   PetscCall(PetscFunctionListAdd(&physics, "isogas", PhysicsCreate_IsoGas));
1429:   PetscCall(PetscFunctionListAdd(&physics, "shallow", PhysicsCreate_Shallow));

1431:   ctx.comm   = comm;
1432:   ctx.cfl    = 0.9;
1433:   ctx.bctype = FVBC_PERIODIC;
1434:   ctx.xmin   = -1;
1435:   ctx.xmax   = 1;
1436:   PetscOptionsBegin(comm, NULL, "Finite Volume solver options", "");
1437:   PetscCall(PetscOptionsReal("-xmin", "X min", "", ctx.xmin, &ctx.xmin, NULL));
1438:   PetscCall(PetscOptionsReal("-xmax", "X max", "", ctx.xmax, &ctx.xmax, NULL));
1439:   PetscCall(PetscOptionsFList("-limit", "Name of flux limiter to use", "", limiters, lname, lname, sizeof(lname), NULL));
1440:   PetscCall(PetscOptionsFList("-physics", "Name of physics (Riemann solver and characteristics) to use", "", physics, physname, physname, sizeof(physname), NULL));
1441:   PetscCall(PetscOptionsInt("-draw", "Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)", "", draw, &draw, NULL));
1442:   PetscCall(PetscOptionsString("-view_final", "Write final solution in ASCII MATLAB format to given file name", "", final_fname, final_fname, sizeof(final_fname), &view_final));
1443:   PetscCall(PetscOptionsInt("-initial", "Initial condition (depends on the physics)", "", ctx.initial, &ctx.initial, NULL));
1444:   PetscCall(PetscOptionsBool("-exact", "Compare errors with exact solution", "", ctx.exact, &ctx.exact, NULL));
1445:   PetscCall(PetscOptionsReal("-cfl", "CFL number to time step at", "", ctx.cfl, &ctx.cfl, NULL));
1446:   PetscCall(PetscOptionsEnum("-bc_type", "Boundary condition", "", FVBCTypes, (PetscEnum)ctx.bctype, (PetscEnum *)&ctx.bctype, NULL));
1447:   PetscOptionsEnd();

1449:   /* Choose the limiter from the list of registered limiters */
1450:   PetscCall(PetscFunctionListFind(limiters, lname, &ctx.limit));
1451:   PetscCheck(ctx.limit, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Limiter '%s' not found", lname);

1453:   /* Choose the physics from the list of registered models */
1454:   {
1455:     PetscErrorCode (*r)(FVCtx *);
1456:     PetscCall(PetscFunctionListFind(physics, physname, &r));
1457:     PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Physics '%s' not found", physname);
1458:     /* Create the physics, will set the number of fields and their names */
1459:     PetscCall((*r)(&ctx));
1460:   }

1462:   /* Create a DMDA to manage the parallel grid */
1463:   PetscCall(DMDACreate1d(comm, DM_BOUNDARY_PERIODIC, 50, ctx.physics.dof, 2, NULL, &da));
1464:   PetscCall(DMSetFromOptions(da));
1465:   PetscCall(DMSetUp(da));
1466:   /* Inform the DMDA of the field names provided by the physics. */
1467:   /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */
1468:   for (i = 0; i < ctx.physics.dof; i++) PetscCall(DMDASetFieldName(da, i, ctx.physics.fieldname[i]));
1469:   PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
1470:   PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0));

1472:   /* Set coordinates of cell centers */
1473:   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));

1475:   /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */
1476:   PetscCall(PetscMalloc4(dof * dof, &ctx.R, dof * dof, &ctx.Rinv, 2 * dof, &ctx.cjmpLR, 1 * dof, &ctx.cslope));
1477:   PetscCall(PetscMalloc3(2 * dof, &ctx.uLR, dof, &ctx.flux, dof, &ctx.speeds));

1479:   /* Create a vector to store the solution and to save the initial state */
1480:   PetscCall(DMCreateGlobalVector(da, &X));
1481:   PetscCall(VecDuplicate(X, &X0));
1482:   PetscCall(VecDuplicate(X, &R));

1484:   PetscCall(DMCreateMatrix(da, &B));

1486:   /* Create a time-stepping object */
1487:   PetscCall(TSCreate(comm, &ts));
1488:   PetscCall(TSSetDM(ts, da));
1489:   PetscCall(TSSetRHSFunction(ts, R, FVRHSFunction, &ctx));
1490:   PetscCall(TSSetIJacobian(ts, B, B, FVIJacobian, &ctx));
1491:   PetscCall(TSSetType(ts, TSSSP));
1492:   PetscCall(TSSetMaxTime(ts, 10));
1493:   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));

1495:   /* Compute initial conditions and starting time step */
1496:   PetscCall(FVSample(&ctx, da, 0, X0));
1497:   PetscCall(FVRHSFunction(ts, 0, X0, X, (void *)&ctx)); /* Initial function evaluation, only used to determine max speed */
1498:   PetscCall(VecCopy(X0, X));                            /* The function value was not used so we set X=X0 again */
1499:   PetscCall(TSSetTimeStep(ts, ctx.cfl / ctx.cfl_idt));
1500:   PetscCall(TSSetFromOptions(ts)); /* Take runtime options */
1501:   PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD));
1502:   {
1503:     PetscReal nrm1, nrmsup;
1504:     PetscInt  steps;

1506:     PetscCall(TSSolve(ts, X));
1507:     PetscCall(TSGetSolveTime(ts, &ptime));
1508:     PetscCall(TSGetStepNumber(ts, &steps));

1510:     PetscCall(PetscPrintf(comm, "Final time %8.5f, steps %" PetscInt_FMT "\n", (double)ptime, steps));
1511:     if (ctx.exact) {
1512:       PetscCall(SolutionErrorNorms(&ctx, da, ptime, X, &nrm1, &nrmsup));
1513:       PetscCall(PetscPrintf(comm, "Error ||x-x_e||_1 %8.4e  ||x-x_e||_sup %8.4e\n", (double)nrm1, (double)nrmsup));
1514:     }
1515:   }

1517:   PetscCall(SolutionStatsView(da, X, PETSC_VIEWER_STDOUT_WORLD));
1518:   if (draw & 0x1) PetscCall(VecView(X0, PETSC_VIEWER_DRAW_WORLD));
1519:   if (draw & 0x2) PetscCall(VecView(X, PETSC_VIEWER_DRAW_WORLD));
1520:   if (draw & 0x4) {
1521:     Vec Y;
1522:     PetscCall(VecDuplicate(X, &Y));
1523:     PetscCall(FVSample(&ctx, da, ptime, Y));
1524:     PetscCall(VecAYPX(Y, -1, X));
1525:     PetscCall(VecView(Y, PETSC_VIEWER_DRAW_WORLD));
1526:     PetscCall(VecDestroy(&Y));
1527:   }

1529:   if (view_final) {
1530:     PetscViewer viewer;
1531:     PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, final_fname, &viewer));
1532:     PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB));
1533:     PetscCall(VecView(X, viewer));
1534:     PetscCall(PetscViewerPopFormat(viewer));
1535:     PetscCall(PetscViewerDestroy(&viewer));
1536:   }

1538:   /* Clean up */
1539:   PetscCall((*ctx.physics.destroy)(ctx.physics.user));
1540:   for (i = 0; i < ctx.physics.dof; i++) PetscCall(PetscFree(ctx.physics.fieldname[i]));
1541:   PetscCall(PetscFree4(ctx.R, ctx.Rinv, ctx.cjmpLR, ctx.cslope));
1542:   PetscCall(PetscFree3(ctx.uLR, ctx.flux, ctx.speeds));
1543:   PetscCall(VecDestroy(&X));
1544:   PetscCall(VecDestroy(&X0));
1545:   PetscCall(VecDestroy(&R));
1546:   PetscCall(MatDestroy(&B));
1547:   PetscCall(DMDestroy(&da));
1548:   PetscCall(TSDestroy(&ts));
1549:   PetscCall(PetscFunctionListDestroy(&limiters));
1550:   PetscCall(PetscFunctionListDestroy(&physics));
1551:   PetscCall(PetscFinalize());
1552:   return 0;
1553: }

1555: /*TEST

1557:     build:
1558:       requires: !complex

1560:     test:
1561:       args: -da_grid_x 100 -initial 1 -xmin -2 -xmax 5 -exact -limit mc
1562:       requires: !complex !single

1564:     test:
1565:       suffix: 2
1566:       args: -da_grid_x 100 -initial 2 -xmin -2 -xmax 2 -exact -limit mc -physics burgers -bc_type outflow -ts_max_time 1
1567:       filter:  sed "s/at 48/at 0/g"
1568:       requires: !complex !single

1570:     test:
1571:       suffix: 3
1572:       args: -da_grid_x 100 -initial 2 -xmin -2 -xmax 2 -exact -limit mc -physics burgers -bc_type outflow -ts_max_time 1
1573:       nsize: 3
1574:       filter:  sed "s/at 48/at 0/g"
1575:       requires: !complex !single

1577: TEST*/