1de8de29fSMatthew G. Knepley #include <petscdm.h> 2de8de29fSMatthew G. Knepley #include <petscdmceed.h> 3de8de29fSMatthew G. Knepley 4a7d931c6SMatthew G. Knepley #ifdef __CUDACC_RTC__ 5a7d931c6SMatthew G. Knepley #define PETSC_HAVE_LIBCEED 6f0b74427SPierre Jolivet // Define PETSc types to be equal to Ceed types 7a7d931c6SMatthew G. Knepley typedef CeedInt PetscInt; 8a7d931c6SMatthew G. Knepley typedef CeedScalar PetscReal; 9a7d931c6SMatthew G. Knepley typedef CeedScalar PetscScalar; 10a7d931c6SMatthew G. Knepley typedef CeedInt PetscErrorCode; 11a7d931c6SMatthew G. Knepley // Define things we are missing from PETSc headers 12a7d931c6SMatthew G. Knepley #undef PETSC_SUCCESS 13a7d931c6SMatthew G. Knepley #define PETSC_SUCCESS 0 14a7d931c6SMatthew G. Knepley #define PETSC_COMM_SELF MPI_COMM_SELF 15a7d931c6SMatthew G. Knepley #undef PetscFunctionBeginUser 16a7d931c6SMatthew G. Knepley #define PetscFunctionBeginUser 17a7d931c6SMatthew G. Knepley #undef PetscFunctionReturn 18a7d931c6SMatthew G. Knepley #define PetscFunctionReturn(x) return x 19a7d931c6SMatthew G. Knepley #undef PetscCall 20a7d931c6SMatthew G. Knepley #define PetscCall(a) a 21a7d931c6SMatthew G. Knepley #define PetscFunctionReturnVoid() return 22a7d931c6SMatthew G. Knepley // Math definitions 23a7d931c6SMatthew G. Knepley #undef PetscSqrtReal 24a7d931c6SMatthew G. Knepley #define PetscSqrtReal(x) sqrt(x) 25a7d931c6SMatthew G. Knepley #undef PetscSqrtScalar 26a7d931c6SMatthew G. Knepley #define PetscSqrtScalar(x) sqrt(x) 27a7d931c6SMatthew G. Knepley #undef PetscSqr 28a7d931c6SMatthew G. Knepley #define PetscSqr(x) (x * x) 29a7d931c6SMatthew G. Knepley #define PetscSqrReal(x) (x * x) 30a7d931c6SMatthew G. Knepley #define PetscAbsReal(x) abs(x) 31a7d931c6SMatthew G. Knepley #define PetscAbsScalar(x) abs(x) 32a7d931c6SMatthew G. Knepley #define PetscMax(x, y) x > y ? x : y 33a7d931c6SMatthew G. Knepley #define PetscMin(x, y) x < y ? x : y 34a7d931c6SMatthew G. Knepley #define PetscRealPart(a) a 35a7d931c6SMatthew G. Knepley #define PetscPowScalar(a, b) pow(a, b) 36a7d931c6SMatthew G. Knepley #endif 37a7d931c6SMatthew G. Knepley 38de8de29fSMatthew G. Knepley #define DIM 2 /* Geometric dimension */ 39de8de29fSMatthew G. Knepley 40de8de29fSMatthew G. Knepley /* Represents continuum physical equations. */ 41de8de29fSMatthew G. Knepley typedef struct _n_Physics *Physics; 42de8de29fSMatthew G. Knepley 43de8de29fSMatthew G. Knepley /* Physical model includes boundary conditions, initial conditions, and functionals of interest. It is 44de8de29fSMatthew G. Knepley * discretization-independent, but its members depend on the scenario being solved. */ 45de8de29fSMatthew G. Knepley typedef struct _n_Model *Model; 46de8de29fSMatthew G. Knepley 47de8de29fSMatthew G. Knepley struct FieldDescription { 48de8de29fSMatthew G. Knepley const char *name; 49de8de29fSMatthew G. Knepley PetscInt dof; 50de8de29fSMatthew G. Knepley }; 51de8de29fSMatthew G. Knepley 52de8de29fSMatthew G. Knepley struct _n_Physics { 53de8de29fSMatthew G. Knepley void (*riemann)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *); 54de8de29fSMatthew G. Knepley PetscInt dof; /* number of degrees of freedom per cell */ 55de8de29fSMatthew G. Knepley PetscReal maxspeed; /* kludge to pick initial time step, need to add monitoring and step control */ 56de8de29fSMatthew G. Knepley void *data; 57de8de29fSMatthew G. Knepley PetscInt nfields; 58de8de29fSMatthew G. Knepley const struct FieldDescription *field_desc; 59de8de29fSMatthew G. Knepley }; 60de8de29fSMatthew G. Knepley 61de8de29fSMatthew G. Knepley typedef struct { 62de8de29fSMatthew G. Knepley PetscReal gravity; 63de8de29fSMatthew G. Knepley struct { 64de8de29fSMatthew G. Knepley PetscInt Height; 65de8de29fSMatthew G. Knepley PetscInt Speed; 66de8de29fSMatthew G. Knepley PetscInt Energy; 67de8de29fSMatthew G. Knepley } functional; 68de8de29fSMatthew G. Knepley } Physics_SW; 69de8de29fSMatthew G. Knepley 70de8de29fSMatthew G. Knepley typedef struct { 71de8de29fSMatthew G. Knepley PetscReal h; 72de8de29fSMatthew G. Knepley PetscReal uh[DIM]; 73de8de29fSMatthew G. Knepley } SWNode; 74de8de29fSMatthew G. Knepley typedef union 75de8de29fSMatthew G. Knepley { 76de8de29fSMatthew G. Knepley SWNode swnode; 77de8de29fSMatthew G. Knepley PetscReal vals[DIM + 1]; 78de8de29fSMatthew G. Knepley } SWNodeUnion; 79de8de29fSMatthew G. Knepley 80de8de29fSMatthew G. Knepley typedef enum { 81de8de29fSMatthew G. Knepley EULER_IV_SHOCK, 82de8de29fSMatthew G. Knepley EULER_SS_SHOCK, 83de8de29fSMatthew G. Knepley EULER_SHOCK_TUBE, 84de8de29fSMatthew G. Knepley EULER_LINEAR_WAVE 85de8de29fSMatthew G. Knepley } EulerType; 86de8de29fSMatthew G. Knepley 87de8de29fSMatthew G. Knepley typedef struct { 88de8de29fSMatthew G. Knepley PetscReal gamma; 89de8de29fSMatthew G. Knepley PetscReal rhoR; 90de8de29fSMatthew G. Knepley PetscReal amach; 91de8de29fSMatthew G. Knepley PetscReal itana; 92de8de29fSMatthew G. Knepley EulerType type; 93de8de29fSMatthew G. Knepley struct { 94de8de29fSMatthew G. Knepley PetscInt Density; 95de8de29fSMatthew G. Knepley PetscInt Momentum; 96de8de29fSMatthew G. Knepley PetscInt Energy; 97de8de29fSMatthew G. Knepley PetscInt Pressure; 98de8de29fSMatthew G. Knepley PetscInt Speed; 99de8de29fSMatthew G. Knepley } monitor; 100de8de29fSMatthew G. Knepley } Physics_Euler; 101de8de29fSMatthew G. Knepley 102de8de29fSMatthew G. Knepley typedef struct { 103de8de29fSMatthew G. Knepley PetscReal r; 104de8de29fSMatthew G. Knepley PetscReal ru[DIM]; 105de8de29fSMatthew G. Knepley PetscReal E; 106de8de29fSMatthew G. Knepley } EulerNode; 107de8de29fSMatthew G. Knepley typedef union 108de8de29fSMatthew G. Knepley { 109de8de29fSMatthew G. Knepley EulerNode eulernode; 110de8de29fSMatthew G. Knepley PetscReal vals[DIM + 2]; 111de8de29fSMatthew G. Knepley } EulerNodeUnion; 112de8de29fSMatthew G. Knepley 113de8de29fSMatthew G. Knepley static inline PetscReal Dot2Real(const PetscReal *x, const PetscReal *y) 114de8de29fSMatthew G. Knepley { 115de8de29fSMatthew G. Knepley return x[0] * y[0] + x[1] * y[1]; 116de8de29fSMatthew G. Knepley } 117de8de29fSMatthew G. Knepley static inline PetscReal Norm2Real(const PetscReal *x) 118de8de29fSMatthew G. Knepley { 119de8de29fSMatthew G. Knepley return PetscSqrtReal(PetscAbsReal(Dot2Real(x, x))); 120de8de29fSMatthew G. Knepley } 121de8de29fSMatthew G. Knepley static inline void Normalize2Real(PetscReal *x) 122de8de29fSMatthew G. Knepley { 123de8de29fSMatthew G. Knepley PetscReal a = 1. / Norm2Real(x); 124de8de29fSMatthew G. Knepley x[0] *= a; 125de8de29fSMatthew G. Knepley x[1] *= a; 126de8de29fSMatthew G. Knepley } 127de8de29fSMatthew G. Knepley static inline void Scale2Real(PetscReal a, const PetscReal *x, PetscReal *y) 128de8de29fSMatthew G. Knepley { 129de8de29fSMatthew G. Knepley y[0] = a * x[0]; 130de8de29fSMatthew G. Knepley y[1] = a * x[1]; 131de8de29fSMatthew G. Knepley } 132de8de29fSMatthew G. Knepley 133de8de29fSMatthew G. Knepley static inline PetscReal DotDIMReal(const PetscReal *x, const PetscReal *y) 134de8de29fSMatthew G. Knepley { 135de8de29fSMatthew G. Knepley PetscInt i; 136de8de29fSMatthew G. Knepley PetscReal prod = 0.0; 137de8de29fSMatthew G. Knepley 138de8de29fSMatthew G. Knepley for (i = 0; i < DIM; i++) prod += x[i] * y[i]; 139de8de29fSMatthew G. Knepley return prod; 140de8de29fSMatthew G. Knepley } 141de8de29fSMatthew G. Knepley static inline PetscReal NormDIM(const PetscReal *x) 142de8de29fSMatthew G. Knepley { 143de8de29fSMatthew G. Knepley return PetscSqrtReal(PetscAbsReal(DotDIMReal(x, x))); 144de8de29fSMatthew G. Knepley } 145de8de29fSMatthew G. Knepley static inline void Waxpy2Real(PetscReal a, const PetscReal *x, const PetscReal *y, PetscReal *w) 146de8de29fSMatthew G. Knepley { 147de8de29fSMatthew G. Knepley w[0] = a * x[0] + y[0]; 148de8de29fSMatthew G. Knepley w[1] = a * x[1] + y[1]; 149de8de29fSMatthew G. Knepley } 150de8de29fSMatthew G. Knepley 151de8de29fSMatthew G. Knepley /* 152de8de29fSMatthew G. Knepley * h_t + div(uh) = 0 153de8de29fSMatthew G. Knepley * (uh)_t + div (u\otimes uh + g h^2 / 2 I) = 0 154de8de29fSMatthew G. Knepley * 155de8de29fSMatthew G. Knepley * */ 156de8de29fSMatthew G. Knepley static PetscErrorCode SWFlux(Physics phys, const PetscReal *n, const SWNode *x, SWNode *f) 157de8de29fSMatthew G. Knepley { 158de8de29fSMatthew G. Knepley Physics_SW *sw = (Physics_SW *)phys->data; 159de8de29fSMatthew G. Knepley PetscReal uhn, u[DIM]; 160de8de29fSMatthew G. Knepley PetscInt i; 161de8de29fSMatthew G. Knepley 162de8de29fSMatthew G. Knepley PetscFunctionBeginUser; 163de8de29fSMatthew G. Knepley Scale2Real(1. / x->h, x->uh, u); 164de8de29fSMatthew G. Knepley uhn = x->uh[0] * n[0] + x->uh[1] * n[1]; 165de8de29fSMatthew G. Knepley f->h = uhn; 166de8de29fSMatthew G. Knepley for (i = 0; i < DIM; i++) f->uh[i] = u[i] * uhn + sw->gravity * PetscSqr(x->h) * n[i]; 167de8de29fSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 168de8de29fSMatthew G. Knepley } 169de8de29fSMatthew G. Knepley 170de8de29fSMatthew G. Knepley static void PhysicsRiemann_SW_Rusanov(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, Physics phys) 171de8de29fSMatthew G. Knepley { 172de8de29fSMatthew G. Knepley Physics_SW *sw = (Physics_SW *)phys->data; 173de8de29fSMatthew G. Knepley PetscReal cL, cR, speed; 174de8de29fSMatthew G. Knepley PetscReal nn[DIM]; 175de8de29fSMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX) 176de8de29fSMatthew G. Knepley const SWNode *uL = (const SWNode *)xL, *uR = (const SWNode *)xR; 177de8de29fSMatthew G. Knepley #else 178de8de29fSMatthew G. Knepley SWNodeUnion uLreal, uRreal; 179de8de29fSMatthew G. Knepley const SWNode *uL = &uLreal.swnode; 180de8de29fSMatthew G. Knepley const SWNode *uR = &uRreal.swnode; 181de8de29fSMatthew G. Knepley #endif 182de8de29fSMatthew G. Knepley SWNodeUnion fL, fR; 183de8de29fSMatthew G. Knepley PetscInt i; 184de8de29fSMatthew G. Knepley PetscReal zero = 0.; 185a7d931c6SMatthew G. Knepley PetscErrorCode ierr; 186de8de29fSMatthew G. Knepley 187de8de29fSMatthew G. Knepley #if defined(PETSC_USE_COMPLEX) 188de8de29fSMatthew G. Knepley uLreal.swnode.h = 0; 189de8de29fSMatthew G. Knepley uRreal.swnode.h = 0; 190de8de29fSMatthew G. Knepley for (i = 0; i < 1 + dim; i++) uLreal.vals[i] = PetscRealPart(xL[i]); 191de8de29fSMatthew G. Knepley for (i = 0; i < 1 + dim; i++) uRreal.vals[i] = PetscRealPart(xR[i]); 192de8de29fSMatthew G. Knepley #endif 193a7d931c6SMatthew G. Knepley 194de8de29fSMatthew G. Knepley if (uL->h < 0 || uR->h < 0) { 195a7d931c6SMatthew G. Knepley // reconstructed thickness is negative 196a7d931c6SMatthew G. Knepley PetscCallVoid(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 197a7d931c6SMatthew G. Knepley for (i = 0; i < 1 + dim; ++i) flux[i] = zero / zero; 198a7d931c6SMatthew G. Knepley PetscCallVoid(PetscFPTrapPop()); 199de8de29fSMatthew G. Knepley return; 200a7d931c6SMatthew G. Knepley } 201a7d931c6SMatthew G. Knepley 202de8de29fSMatthew G. Knepley nn[0] = n[0]; 203de8de29fSMatthew G. Knepley nn[1] = n[1]; 204de8de29fSMatthew G. Knepley Normalize2Real(nn); 205a7d931c6SMatthew G. Knepley ierr = SWFlux(phys, nn, uL, &fL.swnode); 206a7d931c6SMatthew G. Knepley if (ierr) { 207a7d931c6SMatthew G. Knepley PetscCallVoid(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 208a7d931c6SMatthew G. Knepley for (i = 0; i < 1 + dim; ++i) fL.vals[i] = zero / zero; 209a7d931c6SMatthew G. Knepley PetscCallVoid(PetscFPTrapPop()); 210a7d931c6SMatthew G. Knepley } 211a7d931c6SMatthew G. Knepley ierr = SWFlux(phys, nn, uR, &fR.swnode); 212a7d931c6SMatthew G. Knepley if (ierr) { 213a7d931c6SMatthew G. Knepley PetscCallVoid(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 214a7d931c6SMatthew G. Knepley for (i = 0; i < 1 + dim; ++i) fR.vals[i] = zero / zero; 215a7d931c6SMatthew G. Knepley PetscCallVoid(PetscFPTrapPop()); 216a7d931c6SMatthew G. Knepley } 217de8de29fSMatthew G. Knepley cL = PetscSqrtReal(sw->gravity * uL->h); 218de8de29fSMatthew G. Knepley cR = PetscSqrtReal(sw->gravity * uR->h); /* gravity wave speed */ 219de8de29fSMatthew G. Knepley speed = PetscMax(PetscAbsReal(Dot2Real(uL->uh, nn) / uL->h) + cL, PetscAbsReal(Dot2Real(uR->uh, nn) / uR->h) + cR); 220de8de29fSMatthew G. Knepley for (i = 0; i < 1 + dim; i++) flux[i] = (0.5 * (fL.vals[i] + fR.vals[i]) + 0.5 * speed * (xL[i] - xR[i])) * Norm2Real(n); 221de8de29fSMatthew G. Knepley #if 0 222de8de29fSMatthew G. Knepley PetscPrintf(PETSC_COMM_SELF, "Rusanov Flux (%g)\n", sw->gravity); 2233a7d0413SPierre Jolivet for (PetscInt j = 0; j < 3; ++j) PetscPrintf(PETSC_COMM_SELF, " | %g |\n", flux[j]); 224de8de29fSMatthew G. Knepley #endif 225de8de29fSMatthew G. Knepley } 226de8de29fSMatthew G. Knepley 227de8de29fSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED 228*2a8381b2SBarry Smith CEED_QFUNCTION(PhysicsRiemann_SW_Rusanov_CEED)(PetscCtx ctx, CeedInt Q, const CeedScalar *const in[], CeedScalar *const out[]) 229de8de29fSMatthew G. Knepley { 230de8de29fSMatthew G. Knepley const CeedScalar *xL = in[0], *xR = in[1], *geom = in[2]; 231de8de29fSMatthew G. Knepley CeedScalar *cL = out[0], *cR = out[1]; 232de8de29fSMatthew G. Knepley const Physics_SW *sw = (Physics_SW *)ctx; 233de8de29fSMatthew G. Knepley struct _n_Physics phys; 234a7d931c6SMatthew G. Knepley #if 0 235a7d931c6SMatthew G. Knepley const CeedScalar *info = in[3]; 236a7d931c6SMatthew G. Knepley #endif 237de8de29fSMatthew G. Knepley 238de8de29fSMatthew G. Knepley phys.data = (void *)sw; 239de8de29fSMatthew G. Knepley CeedPragmaSIMD for (CeedInt i = 0; i < Q; ++i) 240de8de29fSMatthew G. Knepley { 241de8de29fSMatthew G. Knepley const CeedScalar qL[3] = {xL[i + Q * 0], xL[i + Q * 1], xL[i + Q * 2]}; 242de8de29fSMatthew G. Knepley const CeedScalar qR[3] = {xR[i + Q * 0], xR[i + Q * 1], xR[i + Q * 2]}; 243de8de29fSMatthew G. Knepley const CeedScalar n[2] = {geom[i + Q * 0], geom[i + Q * 1]}; 244de8de29fSMatthew G. Knepley CeedScalar flux[3]; 245de8de29fSMatthew G. Knepley 246de8de29fSMatthew G. Knepley #if 0 247a7d931c6SMatthew G. Knepley PetscPrintf(PETSC_COMM_SELF, "Face %d Normal\n", (int)info[i + Q * 0]); 2483a7d0413SPierre Jolivet for (CeedInt j = 0; j < DIM; ++j) PetscPrintf(PETSC_COMM_SELF, " | %g |\n", n[j]); 249a7d931c6SMatthew G. Knepley PetscPrintf(PETSC_COMM_SELF, "Cell %d Element Residual: left state\n", (int)info[i + Q * 1]); 2503a7d0413SPierre Jolivet for (CeedInt j = 0; j < DIM + 1; ++j) PetscPrintf(PETSC_COMM_SELF, " | %g |\n", qL[j]); 251a7d931c6SMatthew G. Knepley PetscPrintf(PETSC_COMM_SELF, "Cell %d Element Residual: right state\n", (int)info[i + Q * 2]); 2523a7d0413SPierre Jolivet for (CeedInt j = 0; j < DIM + 1; ++j) PetscPrintf(PETSC_COMM_SELF, " | %g |\n", qR[j]); 253de8de29fSMatthew G. Knepley #endif 254de8de29fSMatthew G. Knepley PhysicsRiemann_SW_Rusanov(DIM, DIM + 1, NULL, n, qL, qR, 0, NULL, flux, &phys); 255de8de29fSMatthew G. Knepley for (CeedInt j = 0; j < 3; ++j) { 256de8de29fSMatthew G. Knepley cL[i + Q * j] = -flux[j] / geom[i + Q * 2]; 257de8de29fSMatthew G. Knepley cR[i + Q * j] = flux[j] / geom[i + Q * 3]; 258de8de29fSMatthew G. Knepley } 259de8de29fSMatthew G. Knepley #if 0 260a7d931c6SMatthew G. Knepley PetscPrintf(PETSC_COMM_SELF, "Cell %d Element Residual: left flux\n", (int)info[i + Q * 1]); 2613a7d0413SPierre Jolivet for (CeedInt j = 0; j < DIM + 1; ++j) PetscPrintf(PETSC_COMM_SELF, " | %g | (%g)\n", cL[i + Q * j], geom[i + Q * 2]); 262a7d931c6SMatthew G. Knepley PetscPrintf(PETSC_COMM_SELF, "Cell %d Element Residual: right flux\n", (int)info[i + Q * 2]); 2633a7d0413SPierre Jolivet for (CeedInt j = 0; j < DIM + 1; ++j) PetscPrintf(PETSC_COMM_SELF, " | %g | (%g)\n", cR[i + Q * j], geom[i + Q * 3]); 264de8de29fSMatthew G. Knepley #endif 265de8de29fSMatthew G. Knepley } 266de8de29fSMatthew G. Knepley return CEED_ERROR_SUCCESS; 267de8de29fSMatthew G. Knepley } 268de8de29fSMatthew G. Knepley #endif 269de8de29fSMatthew G. Knepley 270de8de29fSMatthew G. Knepley static void PhysicsRiemann_SW_HLL(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, Physics phys) 271de8de29fSMatthew G. Knepley { 272de8de29fSMatthew G. Knepley Physics_SW *sw = (Physics_SW *)phys->data; 273de8de29fSMatthew G. Knepley PetscReal aL, aR; 274de8de29fSMatthew G. Knepley PetscReal nn[DIM]; 275de8de29fSMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX) 276de8de29fSMatthew G. Knepley const SWNode *uL = (const SWNode *)xL, *uR = (const SWNode *)xR; 277de8de29fSMatthew G. Knepley #else 278de8de29fSMatthew G. Knepley SWNodeUnion uLreal, uRreal; 279de8de29fSMatthew G. Knepley const SWNode *uL = &uLreal.swnode; 280de8de29fSMatthew G. Knepley const SWNode *uR = &uRreal.swnode; 281de8de29fSMatthew G. Knepley #endif 282de8de29fSMatthew G. Knepley SWNodeUnion fL, fR; 283de8de29fSMatthew G. Knepley PetscInt i; 284de8de29fSMatthew G. Knepley PetscReal zero = 0.; 285a7d931c6SMatthew G. Knepley PetscErrorCode ierr; 286de8de29fSMatthew G. Knepley 287de8de29fSMatthew G. Knepley #if defined(PETSC_USE_COMPLEX) 288de8de29fSMatthew G. Knepley uLreal.swnode.h = 0; 289de8de29fSMatthew G. Knepley uRreal.swnode.h = 0; 290de8de29fSMatthew G. Knepley for (i = 0; i < 1 + dim; i++) uLreal.vals[i] = PetscRealPart(xL[i]); 291de8de29fSMatthew G. Knepley for (i = 0; i < 1 + dim; i++) uRreal.vals[i] = PetscRealPart(xR[i]); 292de8de29fSMatthew G. Knepley #endif 293de8de29fSMatthew G. Knepley if (uL->h <= 0 || uR->h <= 0) { 294a7d931c6SMatthew G. Knepley PetscCallVoid(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 295de8de29fSMatthew G. Knepley for (i = 0; i < 1 + dim; i++) flux[i] = zero; 296a7d931c6SMatthew G. Knepley PetscCallVoid(PetscFPTrapPop()); 297de8de29fSMatthew G. Knepley return; 298de8de29fSMatthew G. Knepley } /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed thickness is negative"); */ 299de8de29fSMatthew G. Knepley nn[0] = n[0]; 300de8de29fSMatthew G. Knepley nn[1] = n[1]; 301de8de29fSMatthew G. Knepley Normalize2Real(nn); 302a7d931c6SMatthew G. Knepley ierr = SWFlux(phys, nn, uL, &fL.swnode); 303a7d931c6SMatthew G. Knepley if (ierr) { 304a7d931c6SMatthew G. Knepley PetscCallVoid(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 305a7d931c6SMatthew G. Knepley for (i = 0; i < 1 + dim; ++i) fL.vals[i] = zero / zero; 306a7d931c6SMatthew G. Knepley PetscCallVoid(PetscFPTrapPop()); 307a7d931c6SMatthew G. Knepley } 308a7d931c6SMatthew G. Knepley ierr = SWFlux(phys, nn, uR, &fR.swnode); 309a7d931c6SMatthew G. Knepley if (ierr) { 310a7d931c6SMatthew G. Knepley PetscCallVoid(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 311a7d931c6SMatthew G. Knepley for (i = 0; i < 1 + dim; ++i) fR.vals[i] = zero / zero; 312a7d931c6SMatthew G. Knepley PetscCallVoid(PetscFPTrapPop()); 313a7d931c6SMatthew G. Knepley } 314de8de29fSMatthew G. Knepley /* gravity wave speed */ 315de8de29fSMatthew G. Knepley aL = PetscSqrtReal(sw->gravity * uL->h); 316de8de29fSMatthew G. Knepley aR = PetscSqrtReal(sw->gravity * uR->h); 317de8de29fSMatthew G. Knepley // Defining u_tilda and v_tilda as u and v 318de8de29fSMatthew G. Knepley PetscReal u_L, u_R; 319de8de29fSMatthew G. Knepley u_L = Dot2Real(uL->uh, nn) / uL->h; 320de8de29fSMatthew G. Knepley u_R = Dot2Real(uR->uh, nn) / uR->h; 321de8de29fSMatthew G. Knepley PetscReal sL, sR; 322de8de29fSMatthew G. Knepley sL = PetscMin(u_L - aL, u_R - aR); 323de8de29fSMatthew G. Knepley sR = PetscMax(u_L + aL, u_R + aR); 324de8de29fSMatthew G. Knepley if (sL > zero) { 325de8de29fSMatthew G. Knepley for (i = 0; i < dim + 1; i++) flux[i] = fL.vals[i] * Norm2Real(n); 326de8de29fSMatthew G. Knepley } else if (sR < zero) { 327de8de29fSMatthew G. Knepley for (i = 0; i < dim + 1; i++) flux[i] = fR.vals[i] * Norm2Real(n); 328de8de29fSMatthew G. Knepley } else { 329de8de29fSMatthew G. Knepley for (i = 0; i < dim + 1; i++) flux[i] = ((sR * fL.vals[i] - sL * fR.vals[i] + sR * sL * (xR[i] - xL[i])) / (sR - sL)) * Norm2Real(n); 330de8de29fSMatthew G. Knepley } 331de8de29fSMatthew G. Knepley } 332de8de29fSMatthew G. Knepley 333de8de29fSMatthew G. Knepley static PetscErrorCode Pressure_PG(const PetscReal gamma, const EulerNode *x, PetscReal *p) 334de8de29fSMatthew G. Knepley { 335de8de29fSMatthew G. Knepley PetscReal ru2; 336de8de29fSMatthew G. Knepley 337de8de29fSMatthew G. Knepley PetscFunctionBeginUser; 338de8de29fSMatthew G. Knepley ru2 = DotDIMReal(x->ru, x->ru); 339de8de29fSMatthew G. Knepley (*p) = (x->E - 0.5 * ru2 / x->r) * (gamma - 1.0); /* (E - rho V^2/2)(gamma-1) = e rho (gamma-1) */ 340de8de29fSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 341de8de29fSMatthew G. Knepley } 342de8de29fSMatthew G. Knepley 343de8de29fSMatthew G. Knepley static PetscErrorCode SpeedOfSound_PG(const PetscReal gamma, const EulerNode *x, PetscReal *c) 344de8de29fSMatthew G. Knepley { 345de8de29fSMatthew G. Knepley PetscReal p; 346de8de29fSMatthew G. Knepley 347de8de29fSMatthew G. Knepley PetscFunctionBeginUser; 348de8de29fSMatthew G. Knepley PetscCall(Pressure_PG(gamma, x, &p)); 349de8de29fSMatthew G. Knepley /* gamma = heat capacity ratio */ 350de8de29fSMatthew G. Knepley (*c) = PetscSqrtReal(gamma * p / x->r); 351de8de29fSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 352de8de29fSMatthew G. Knepley } 353de8de29fSMatthew G. Knepley 354de8de29fSMatthew G. Knepley /* 355de8de29fSMatthew G. Knepley * x = (rho,rho*(u_1),...,rho*e)^T 356de8de29fSMatthew G. Knepley * x_t+div(f_1(x))+...+div(f_DIM(x)) = 0 357de8de29fSMatthew G. Knepley * 358de8de29fSMatthew G. Knepley * f_i(x) = u_i*x+(0,0,...,p,...,p*u_i)^T 359de8de29fSMatthew G. Knepley * 360de8de29fSMatthew G. Knepley */ 361de8de29fSMatthew G. Knepley static PetscErrorCode EulerFlux(Physics phys, const PetscReal *n, const EulerNode *x, EulerNode *f) 362de8de29fSMatthew G. Knepley { 363de8de29fSMatthew G. Knepley Physics_Euler *eu = (Physics_Euler *)phys->data; 364de8de29fSMatthew G. Knepley PetscReal nu, p; 365de8de29fSMatthew G. Knepley PetscInt i; 366de8de29fSMatthew G. Knepley 367de8de29fSMatthew G. Knepley PetscFunctionBeginUser; 368de8de29fSMatthew G. Knepley PetscCall(Pressure_PG(eu->gamma, x, &p)); 369de8de29fSMatthew G. Knepley nu = DotDIMReal(x->ru, n); 370de8de29fSMatthew G. Knepley f->r = nu; /* A rho u */ 371de8de29fSMatthew G. Knepley nu /= x->r; /* A u */ 372de8de29fSMatthew G. Knepley for (i = 0; i < DIM; i++) f->ru[i] = nu * x->ru[i] + n[i] * p; /* r u^2 + p */ 373de8de29fSMatthew G. Knepley f->E = nu * (x->E + p); /* u(e+p) */ 374de8de29fSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 375de8de29fSMatthew G. Knepley } 376de8de29fSMatthew G. Knepley 377de8de29fSMatthew G. Knepley /* Godunov fluxs */ 378de8de29fSMatthew G. Knepley static PetscScalar cvmgp_(PetscScalar *a, PetscScalar *b, PetscScalar *test) 379de8de29fSMatthew G. Knepley { 380de8de29fSMatthew G. Knepley /* System generated locals */ 381de8de29fSMatthew G. Knepley PetscScalar ret_val; 382de8de29fSMatthew G. Knepley 383de8de29fSMatthew G. Knepley if (PetscRealPart(*test) > 0.) goto L10; 384de8de29fSMatthew G. Knepley ret_val = *b; 385de8de29fSMatthew G. Knepley return ret_val; 386de8de29fSMatthew G. Knepley L10: 387de8de29fSMatthew G. Knepley ret_val = *a; 388de8de29fSMatthew G. Knepley return ret_val; 389de8de29fSMatthew G. Knepley } /* cvmgp_ */ 390de8de29fSMatthew G. Knepley 391de8de29fSMatthew G. Knepley static PetscScalar cvmgm_(PetscScalar *a, PetscScalar *b, PetscScalar *test) 392de8de29fSMatthew G. Knepley { 393de8de29fSMatthew G. Knepley /* System generated locals */ 394de8de29fSMatthew G. Knepley PetscScalar ret_val; 395de8de29fSMatthew G. Knepley 396de8de29fSMatthew G. Knepley if (PetscRealPart(*test) < 0.) goto L10; 397de8de29fSMatthew G. Knepley ret_val = *b; 398de8de29fSMatthew G. Knepley return ret_val; 399de8de29fSMatthew G. Knepley L10: 400de8de29fSMatthew G. Knepley ret_val = *a; 401de8de29fSMatthew G. Knepley return ret_val; 402de8de29fSMatthew G. Knepley } /* cvmgm_ */ 403de8de29fSMatthew G. Knepley 404de8de29fSMatthew G. Knepley static int riem1mdt(PetscScalar *gaml, PetscScalar *gamr, PetscScalar *rl, PetscScalar *pl, PetscScalar *uxl, PetscScalar *rr, PetscScalar *pr, PetscScalar *uxr, PetscScalar *rstarl, PetscScalar *rstarr, PetscScalar *pstar, PetscScalar *ustar) 405de8de29fSMatthew G. Knepley { 406de8de29fSMatthew G. Knepley /* Initialized data */ 407de8de29fSMatthew G. Knepley 408de8de29fSMatthew G. Knepley static PetscScalar smallp = 1e-8; 409de8de29fSMatthew G. Knepley 410de8de29fSMatthew G. Knepley /* System generated locals */ 411de8de29fSMatthew G. Knepley int i__1; 412de8de29fSMatthew G. Knepley PetscScalar d__1, d__2; 413de8de29fSMatthew G. Knepley 414de8de29fSMatthew G. Knepley /* Local variables */ 415de8de29fSMatthew G. Knepley static int i0; 416de8de29fSMatthew G. Knepley static PetscScalar cl, cr, wl, zl, wr, zr, pst, durl, skpr1, skpr2; 417de8de29fSMatthew G. Knepley static int iwave; 418de8de29fSMatthew G. Knepley static PetscScalar gascl4, gascr4, cstarl, dpstar, cstarr; 419de8de29fSMatthew G. Knepley /* static PetscScalar csqrl, csqrr, gascl1, gascl2, gascl3, gascr1, gascr2, gascr3; */ 420de8de29fSMatthew G. Knepley static int iterno; 421de8de29fSMatthew G. Knepley static PetscScalar ustarl, ustarr, rarepr1, rarepr2; 422de8de29fSMatthew G. Knepley 423de8de29fSMatthew G. Knepley /* gascl1 = *gaml - 1.; */ 424de8de29fSMatthew G. Knepley /* gascl2 = (*gaml + 1.) * .5; */ 425de8de29fSMatthew G. Knepley /* gascl3 = gascl2 / *gaml; */ 426de8de29fSMatthew G. Knepley gascl4 = 1. / (*gaml - 1.); 427de8de29fSMatthew G. Knepley 428de8de29fSMatthew G. Knepley /* gascr1 = *gamr - 1.; */ 429de8de29fSMatthew G. Knepley /* gascr2 = (*gamr + 1.) * .5; */ 430de8de29fSMatthew G. Knepley /* gascr3 = gascr2 / *gamr; */ 431de8de29fSMatthew G. Knepley gascr4 = 1. / (*gamr - 1.); 432de8de29fSMatthew G. Knepley iterno = 10; 433de8de29fSMatthew G. Knepley /* find pstar: */ 434de8de29fSMatthew G. Knepley cl = PetscSqrtScalar(*gaml * *pl / *rl); 435de8de29fSMatthew G. Knepley cr = PetscSqrtScalar(*gamr * *pr / *rr); 436de8de29fSMatthew G. Knepley wl = *rl * cl; 437de8de29fSMatthew G. Knepley wr = *rr * cr; 438de8de29fSMatthew G. Knepley /* csqrl = wl * wl; */ 439de8de29fSMatthew G. Knepley /* csqrr = wr * wr; */ 440de8de29fSMatthew G. Knepley *pstar = (wl * *pr + wr * *pl) / (wl + wr); 441de8de29fSMatthew G. Knepley *pstar = PetscMax(PetscRealPart(*pstar), PetscRealPart(smallp)); 442de8de29fSMatthew G. Knepley pst = *pl / *pr; 443de8de29fSMatthew G. Knepley skpr1 = cr * (pst - 1.) * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst))); 444de8de29fSMatthew G. Knepley d__1 = (*gamr - 1.) / (*gamr * 2.); 445de8de29fSMatthew G. Knepley rarepr2 = gascr4 * 2. * cr * (1. - PetscPowScalar(pst, d__1)); 446de8de29fSMatthew G. Knepley pst = *pr / *pl; 447de8de29fSMatthew G. Knepley skpr2 = cl * (pst - 1.) * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst))); 448de8de29fSMatthew G. Knepley d__1 = (*gaml - 1.) / (*gaml * 2.); 449de8de29fSMatthew G. Knepley rarepr1 = gascl4 * 2. * cl * (1. - PetscPowScalar(pst, d__1)); 450de8de29fSMatthew G. Knepley durl = *uxr - *uxl; 451de8de29fSMatthew G. Knepley if (PetscRealPart(*pr) < PetscRealPart(*pl)) { 452de8de29fSMatthew G. Knepley if (PetscRealPart(durl) >= PetscRealPart(rarepr1)) { 453de8de29fSMatthew G. Knepley iwave = 100; 454de8de29fSMatthew G. Knepley } else if (PetscRealPart(durl) <= PetscRealPart(-skpr1)) { 455de8de29fSMatthew G. Knepley iwave = 300; 456de8de29fSMatthew G. Knepley } else { 457de8de29fSMatthew G. Knepley iwave = 400; 458de8de29fSMatthew G. Knepley } 459de8de29fSMatthew G. Knepley } else { 460de8de29fSMatthew G. Knepley if (PetscRealPart(durl) >= PetscRealPart(rarepr2)) { 461de8de29fSMatthew G. Knepley iwave = 100; 462de8de29fSMatthew G. Knepley } else if (PetscRealPart(durl) <= PetscRealPart(-skpr2)) { 463de8de29fSMatthew G. Knepley iwave = 300; 464de8de29fSMatthew G. Knepley } else { 465de8de29fSMatthew G. Knepley iwave = 200; 466de8de29fSMatthew G. Knepley } 467de8de29fSMatthew G. Knepley } 468de8de29fSMatthew G. Knepley if (iwave == 100) { 469de8de29fSMatthew G. Knepley /* 1-wave: rarefaction wave, 3-wave: rarefaction wave */ 470de8de29fSMatthew G. Knepley /* case (100) */ 471de8de29fSMatthew G. Knepley i__1 = iterno; 472de8de29fSMatthew G. Knepley for (i0 = 1; i0 <= i__1; ++i0) { 473de8de29fSMatthew G. Knepley d__1 = *pstar / *pl; 474de8de29fSMatthew G. Knepley d__2 = 1. / *gaml; 475de8de29fSMatthew G. Knepley *rstarl = *rl * PetscPowScalar(d__1, d__2); 476de8de29fSMatthew G. Knepley cstarl = PetscSqrtScalar(*gaml * *pstar / *rstarl); 477de8de29fSMatthew G. Knepley ustarl = *uxl - gascl4 * 2. * (cstarl - cl); 478de8de29fSMatthew G. Knepley zl = *rstarl * cstarl; 479de8de29fSMatthew G. Knepley d__1 = *pstar / *pr; 480de8de29fSMatthew G. Knepley d__2 = 1. / *gamr; 481de8de29fSMatthew G. Knepley *rstarr = *rr * PetscPowScalar(d__1, d__2); 482de8de29fSMatthew G. Knepley cstarr = PetscSqrtScalar(*gamr * *pstar / *rstarr); 483de8de29fSMatthew G. Knepley ustarr = *uxr + gascr4 * 2. * (cstarr - cr); 484de8de29fSMatthew G. Knepley zr = *rstarr * cstarr; 485de8de29fSMatthew G. Knepley dpstar = zl * zr * (ustarr - ustarl) / (zl + zr); 486de8de29fSMatthew G. Knepley *pstar -= dpstar; 487de8de29fSMatthew G. Knepley *pstar = PetscMax(PetscRealPart(*pstar), PetscRealPart(smallp)); 488de8de29fSMatthew G. Knepley if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) { 489de8de29fSMatthew G. Knepley #if 0 490de8de29fSMatthew G. Knepley break; 491de8de29fSMatthew G. Knepley #endif 492de8de29fSMatthew G. Knepley } 493de8de29fSMatthew G. Knepley } 494de8de29fSMatthew G. Knepley /* 1-wave: shock wave, 3-wave: rarefaction wave */ 495de8de29fSMatthew G. Knepley } else if (iwave == 200) { 496de8de29fSMatthew G. Knepley /* case (200) */ 497de8de29fSMatthew G. Knepley i__1 = iterno; 498de8de29fSMatthew G. Knepley for (i0 = 1; i0 <= i__1; ++i0) { 499de8de29fSMatthew G. Knepley pst = *pstar / *pl; 500de8de29fSMatthew G. Knepley ustarl = *uxl - (pst - 1.) * cl * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst))); 501de8de29fSMatthew G. Knepley zl = *pl / cl * PetscSqrtScalar(*gaml * 2. * (*gaml - 1. + (*gaml + 1.) * pst)) * (*gaml - 1. + (*gaml + 1.) * pst) / (*gaml * 3. - 1. + (*gaml + 1.) * pst); 502de8de29fSMatthew G. Knepley d__1 = *pstar / *pr; 503de8de29fSMatthew G. Knepley d__2 = 1. / *gamr; 504de8de29fSMatthew G. Knepley *rstarr = *rr * PetscPowScalar(d__1, d__2); 505de8de29fSMatthew G. Knepley cstarr = PetscSqrtScalar(*gamr * *pstar / *rstarr); 506de8de29fSMatthew G. Knepley zr = *rstarr * cstarr; 507de8de29fSMatthew G. Knepley ustarr = *uxr + gascr4 * 2. * (cstarr - cr); 508de8de29fSMatthew G. Knepley dpstar = zl * zr * (ustarr - ustarl) / (zl + zr); 509de8de29fSMatthew G. Knepley *pstar -= dpstar; 510de8de29fSMatthew G. Knepley *pstar = PetscMax(PetscRealPart(*pstar), PetscRealPart(smallp)); 511de8de29fSMatthew G. Knepley if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) { 512de8de29fSMatthew G. Knepley #if 0 513de8de29fSMatthew G. Knepley break; 514de8de29fSMatthew G. Knepley #endif 515de8de29fSMatthew G. Knepley } 516de8de29fSMatthew G. Knepley } 517de8de29fSMatthew G. Knepley /* 1-wave: shock wave, 3-wave: shock */ 518de8de29fSMatthew G. Knepley } else if (iwave == 300) { 519de8de29fSMatthew G. Knepley /* case (300) */ 520de8de29fSMatthew G. Knepley i__1 = iterno; 521de8de29fSMatthew G. Knepley for (i0 = 1; i0 <= i__1; ++i0) { 522de8de29fSMatthew G. Knepley pst = *pstar / *pl; 523de8de29fSMatthew G. Knepley ustarl = *uxl - (pst - 1.) * cl * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst))); 524de8de29fSMatthew G. Knepley zl = *pl / cl * PetscSqrtScalar(*gaml * 2. * (*gaml - 1. + (*gaml + 1.) * pst)) * (*gaml - 1. + (*gaml + 1.) * pst) / (*gaml * 3. - 1. + (*gaml + 1.) * pst); 525de8de29fSMatthew G. Knepley pst = *pstar / *pr; 526de8de29fSMatthew G. Knepley ustarr = *uxr + (pst - 1.) * cr * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst))); 527de8de29fSMatthew G. Knepley zr = *pr / cr * PetscSqrtScalar(*gamr * 2. * (*gamr - 1. + (*gamr + 1.) * pst)) * (*gamr - 1. + (*gamr + 1.) * pst) / (*gamr * 3. - 1. + (*gamr + 1.) * pst); 528de8de29fSMatthew G. Knepley dpstar = zl * zr * (ustarr - ustarl) / (zl + zr); 529de8de29fSMatthew G. Knepley *pstar -= dpstar; 530de8de29fSMatthew G. Knepley *pstar = PetscMax(PetscRealPart(*pstar), PetscRealPart(smallp)); 531de8de29fSMatthew G. Knepley if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) { 532de8de29fSMatthew G. Knepley #if 0 533de8de29fSMatthew G. Knepley break; 534de8de29fSMatthew G. Knepley #endif 535de8de29fSMatthew G. Knepley } 536de8de29fSMatthew G. Knepley } 537de8de29fSMatthew G. Knepley /* 1-wave: rarefaction wave, 3-wave: shock */ 538de8de29fSMatthew G. Knepley } else if (iwave == 400) { 539de8de29fSMatthew G. Knepley /* case (400) */ 540de8de29fSMatthew G. Knepley i__1 = iterno; 541de8de29fSMatthew G. Knepley for (i0 = 1; i0 <= i__1; ++i0) { 542de8de29fSMatthew G. Knepley d__1 = *pstar / *pl; 543de8de29fSMatthew G. Knepley d__2 = 1. / *gaml; 544de8de29fSMatthew G. Knepley *rstarl = *rl * PetscPowScalar(d__1, d__2); 545de8de29fSMatthew G. Knepley cstarl = PetscSqrtScalar(*gaml * *pstar / *rstarl); 546de8de29fSMatthew G. Knepley ustarl = *uxl - gascl4 * 2. * (cstarl - cl); 547de8de29fSMatthew G. Knepley zl = *rstarl * cstarl; 548de8de29fSMatthew G. Knepley pst = *pstar / *pr; 549de8de29fSMatthew G. Knepley ustarr = *uxr + (pst - 1.) * cr * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst))); 550de8de29fSMatthew G. Knepley zr = *pr / cr * PetscSqrtScalar(*gamr * 2. * (*gamr - 1. + (*gamr + 1.) * pst)) * (*gamr - 1. + (*gamr + 1.) * pst) / (*gamr * 3. - 1. + (*gamr + 1.) * pst); 551de8de29fSMatthew G. Knepley dpstar = zl * zr * (ustarr - ustarl) / (zl + zr); 552de8de29fSMatthew G. Knepley *pstar -= dpstar; 553de8de29fSMatthew G. Knepley *pstar = PetscMax(PetscRealPart(*pstar), PetscRealPart(smallp)); 554de8de29fSMatthew G. Knepley if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) { 555de8de29fSMatthew G. Knepley #if 0 556de8de29fSMatthew G. Knepley break; 557de8de29fSMatthew G. Knepley #endif 558de8de29fSMatthew G. Knepley } 559de8de29fSMatthew G. Knepley } 560de8de29fSMatthew G. Knepley } 561de8de29fSMatthew G. Knepley 562de8de29fSMatthew G. Knepley *ustar = (zl * ustarr + zr * ustarl) / (zl + zr); 563de8de29fSMatthew G. Knepley if (PetscRealPart(*pstar) > PetscRealPart(*pl)) { 564de8de29fSMatthew G. Knepley pst = *pstar / *pl; 565de8de29fSMatthew G. Knepley *rstarl = ((*gaml + 1.) * pst + *gaml - 1.) / ((*gaml - 1.) * pst + *gaml + 1.) * *rl; 566de8de29fSMatthew G. Knepley } 567de8de29fSMatthew G. Knepley if (PetscRealPart(*pstar) > PetscRealPart(*pr)) { 568de8de29fSMatthew G. Knepley pst = *pstar / *pr; 569de8de29fSMatthew G. Knepley *rstarr = ((*gamr + 1.) * pst + *gamr - 1.) / ((*gamr - 1.) * pst + *gamr + 1.) * *rr; 570de8de29fSMatthew G. Knepley } 571de8de29fSMatthew G. Knepley return iwave; 572de8de29fSMatthew G. Knepley } 573de8de29fSMatthew G. Knepley 574de8de29fSMatthew G. Knepley static PetscScalar sign(PetscScalar x) 575de8de29fSMatthew G. Knepley { 576de8de29fSMatthew G. Knepley if (PetscRealPart(x) > 0) return 1.0; 577de8de29fSMatthew G. Knepley if (PetscRealPart(x) < 0) return -1.0; 578de8de29fSMatthew G. Knepley return 0.0; 579de8de29fSMatthew G. Knepley } 580de8de29fSMatthew G. Knepley /* Riemann Solver */ 581de8de29fSMatthew G. Knepley /* -------------------------------------------------------------------- */ 582de8de29fSMatthew G. Knepley static int riemannsolver(PetscScalar *xcen, PetscScalar *xp, PetscScalar *dtt, PetscScalar *rl, PetscScalar *uxl, PetscScalar *pl, PetscScalar *utl, PetscScalar *ubl, PetscScalar *gaml, PetscScalar *rho1l, PetscScalar *rr, PetscScalar *uxr, PetscScalar *pr, PetscScalar *utr, PetscScalar *ubr, PetscScalar *gamr, PetscScalar *rho1r, PetscScalar *rx, PetscScalar *uxm, PetscScalar *px, PetscScalar *utx, PetscScalar *ubx, PetscScalar *gam, PetscScalar *rho1) 583de8de29fSMatthew G. Knepley { 584de8de29fSMatthew G. Knepley /* System generated locals */ 585de8de29fSMatthew G. Knepley PetscScalar d__1, d__2; 586de8de29fSMatthew G. Knepley 587de8de29fSMatthew G. Knepley /* Local variables */ 588de8de29fSMatthew G. Knepley static PetscScalar s, c0, p0, r0, u0, w0, x0, x2, ri, cx, sgn0, wsp0, gasc1, gasc2, gasc3, gasc4; 589de8de29fSMatthew G. Knepley static PetscScalar cstar, pstar, rstar, ustar, xstar, wspst, ushock, streng, rstarl, rstarr, rstars; 590de8de29fSMatthew G. Knepley int iwave; 591de8de29fSMatthew G. Knepley 592de8de29fSMatthew G. Knepley if (*rl == *rr && *pr == *pl && *uxl == *uxr && *gaml == *gamr) { 593de8de29fSMatthew G. Knepley *rx = *rl; 594de8de29fSMatthew G. Knepley *px = *pl; 595de8de29fSMatthew G. Knepley *uxm = *uxl; 596de8de29fSMatthew G. Knepley *gam = *gaml; 597de8de29fSMatthew G. Knepley x2 = *xcen + *uxm * *dtt; 598de8de29fSMatthew G. Knepley 599de8de29fSMatthew G. Knepley if (PetscRealPart(*xp) >= PetscRealPart(x2)) { 600de8de29fSMatthew G. Knepley *utx = *utr; 601de8de29fSMatthew G. Knepley *ubx = *ubr; 602de8de29fSMatthew G. Knepley *rho1 = *rho1r; 603de8de29fSMatthew G. Knepley } else { 604de8de29fSMatthew G. Knepley *utx = *utl; 605de8de29fSMatthew G. Knepley *ubx = *ubl; 606de8de29fSMatthew G. Knepley *rho1 = *rho1l; 607de8de29fSMatthew G. Knepley } 608de8de29fSMatthew G. Knepley return 0; 609de8de29fSMatthew G. Knepley } 610de8de29fSMatthew G. Knepley iwave = riem1mdt(gaml, gamr, rl, pl, uxl, rr, pr, uxr, &rstarl, &rstarr, &pstar, &ustar); 611de8de29fSMatthew G. Knepley 612de8de29fSMatthew G. Knepley x2 = *xcen + ustar * *dtt; 613de8de29fSMatthew G. Knepley d__1 = *xp - x2; 614de8de29fSMatthew G. Knepley sgn0 = sign(d__1); 615de8de29fSMatthew G. Knepley /* x is in 3-wave if sgn0 = 1 */ 616de8de29fSMatthew G. Knepley /* x is in 1-wave if sgn0 = -1 */ 617de8de29fSMatthew G. Knepley r0 = cvmgm_(rl, rr, &sgn0); 618de8de29fSMatthew G. Knepley p0 = cvmgm_(pl, pr, &sgn0); 619de8de29fSMatthew G. Knepley u0 = cvmgm_(uxl, uxr, &sgn0); 620de8de29fSMatthew G. Knepley *gam = cvmgm_(gaml, gamr, &sgn0); 621de8de29fSMatthew G. Knepley gasc1 = *gam - 1.; 622de8de29fSMatthew G. Knepley gasc2 = (*gam + 1.) * .5; 623de8de29fSMatthew G. Knepley gasc3 = gasc2 / *gam; 624de8de29fSMatthew G. Knepley gasc4 = 1. / (*gam - 1.); 625de8de29fSMatthew G. Knepley c0 = PetscSqrtScalar(*gam * p0 / r0); 626de8de29fSMatthew G. Knepley streng = pstar - p0; 627de8de29fSMatthew G. Knepley w0 = *gam * r0 * p0 * (gasc3 * streng / p0 + 1.); 628de8de29fSMatthew G. Knepley rstars = r0 / (1. - r0 * streng / w0); 629de8de29fSMatthew G. Knepley d__1 = p0 / pstar; 630de8de29fSMatthew G. Knepley d__2 = -1. / *gam; 631de8de29fSMatthew G. Knepley rstarr = r0 * PetscPowScalar(d__1, d__2); 632de8de29fSMatthew G. Knepley rstar = cvmgm_(&rstarr, &rstars, &streng); 633de8de29fSMatthew G. Knepley w0 = PetscSqrtScalar(w0); 634de8de29fSMatthew G. Knepley cstar = PetscSqrtScalar(*gam * pstar / rstar); 635de8de29fSMatthew G. Knepley wsp0 = u0 + sgn0 * c0; 636de8de29fSMatthew G. Knepley wspst = ustar + sgn0 * cstar; 637de8de29fSMatthew G. Knepley ushock = ustar + sgn0 * w0 / rstar; 638de8de29fSMatthew G. Knepley wspst = cvmgp_(&ushock, &wspst, &streng); 639de8de29fSMatthew G. Knepley wsp0 = cvmgp_(&ushock, &wsp0, &streng); 640de8de29fSMatthew G. Knepley x0 = *xcen + wsp0 * *dtt; 641de8de29fSMatthew G. Knepley xstar = *xcen + wspst * *dtt; 642de8de29fSMatthew G. Knepley /* using gas formula to evaluate rarefaction wave */ 643de8de29fSMatthew G. Knepley /* ri : reiman invariant */ 644de8de29fSMatthew G. Knepley ri = u0 - sgn0 * 2. * gasc4 * c0; 645de8de29fSMatthew G. Knepley cx = sgn0 * .5 * gasc1 / gasc2 * ((*xp - *xcen) / *dtt - ri); 646de8de29fSMatthew G. Knepley *uxm = ri + sgn0 * 2. * gasc4 * cx; 647de8de29fSMatthew G. Knepley s = p0 / PetscPowScalar(r0, *gam); 648de8de29fSMatthew G. Knepley d__1 = cx * cx / (*gam * s); 649de8de29fSMatthew G. Knepley *rx = PetscPowScalar(d__1, gasc4); 650de8de29fSMatthew G. Knepley *px = cx * cx * *rx / *gam; 651de8de29fSMatthew G. Knepley d__1 = sgn0 * (x0 - *xp); 652de8de29fSMatthew G. Knepley *rx = cvmgp_(rx, &r0, &d__1); 653de8de29fSMatthew G. Knepley d__1 = sgn0 * (x0 - *xp); 654de8de29fSMatthew G. Knepley *px = cvmgp_(px, &p0, &d__1); 655de8de29fSMatthew G. Knepley d__1 = sgn0 * (x0 - *xp); 656de8de29fSMatthew G. Knepley *uxm = cvmgp_(uxm, &u0, &d__1); 657de8de29fSMatthew G. Knepley d__1 = sgn0 * (xstar - *xp); 658de8de29fSMatthew G. Knepley *rx = cvmgm_(rx, &rstar, &d__1); 659de8de29fSMatthew G. Knepley d__1 = sgn0 * (xstar - *xp); 660de8de29fSMatthew G. Knepley *px = cvmgm_(px, &pstar, &d__1); 661de8de29fSMatthew G. Knepley d__1 = sgn0 * (xstar - *xp); 662de8de29fSMatthew G. Knepley *uxm = cvmgm_(uxm, &ustar, &d__1); 663de8de29fSMatthew G. Knepley if (PetscRealPart(*xp) >= PetscRealPart(x2)) { 664de8de29fSMatthew G. Knepley *utx = *utr; 665de8de29fSMatthew G. Knepley *ubx = *ubr; 666de8de29fSMatthew G. Knepley *rho1 = *rho1r; 667de8de29fSMatthew G. Knepley } else { 668de8de29fSMatthew G. Knepley *utx = *utl; 669de8de29fSMatthew G. Knepley *ubx = *ubl; 670de8de29fSMatthew G. Knepley *rho1 = *rho1l; 671de8de29fSMatthew G. Knepley } 672de8de29fSMatthew G. Knepley return iwave; 673de8de29fSMatthew G. Knepley } 674de8de29fSMatthew G. Knepley 675de8de29fSMatthew G. Knepley static int godunovflux(const PetscScalar *ul, const PetscScalar *ur, PetscScalar *flux, const PetscReal *nn, int ndim, PetscReal gamma) 676de8de29fSMatthew G. Knepley { 677de8de29fSMatthew G. Knepley /* System generated locals */ 678de8de29fSMatthew G. Knepley int i__1, iwave; 679de8de29fSMatthew G. Knepley PetscScalar d__1, d__2, d__3; 680de8de29fSMatthew G. Knepley 681de8de29fSMatthew G. Knepley /* Local variables */ 682de8de29fSMatthew G. Knepley static int k; 683de8de29fSMatthew G. Knepley static PetscScalar bn[3], fn, ft, tg[3], pl, rl, pm, pr, rr, xp, ubl, ubm, ubr, dtt, unm, tmp, utl, utm, uxl, utr, uxr, gaml, gamm, gamr, xcen, rhom, rho1l, rho1m, rho1r; 684de8de29fSMatthew G. Knepley 685de8de29fSMatthew G. Knepley /* Function Body */ 686de8de29fSMatthew G. Knepley xcen = 0.; 687de8de29fSMatthew G. Knepley xp = 0.; 688de8de29fSMatthew G. Knepley i__1 = ndim; 689de8de29fSMatthew G. Knepley for (k = 1; k <= i__1; ++k) { 690de8de29fSMatthew G. Knepley tg[k - 1] = 0.; 691de8de29fSMatthew G. Knepley bn[k - 1] = 0.; 692de8de29fSMatthew G. Knepley } 693de8de29fSMatthew G. Knepley dtt = 1.; 694de8de29fSMatthew G. Knepley if (ndim == 3) { 695de8de29fSMatthew G. Knepley if (nn[0] == 0. && nn[1] == 0.) { 696de8de29fSMatthew G. Knepley tg[0] = 1.; 697de8de29fSMatthew G. Knepley } else { 698de8de29fSMatthew G. Knepley tg[0] = -nn[1]; 699de8de29fSMatthew G. Knepley tg[1] = nn[0]; 700de8de29fSMatthew G. Knepley } 701de8de29fSMatthew G. Knepley /* tmp=dsqrt(tg(1)**2+tg(2)**2) */ 702de8de29fSMatthew G. Knepley /* tg=tg/tmp */ 703de8de29fSMatthew G. Knepley bn[0] = -nn[2] * tg[1]; 704de8de29fSMatthew G. Knepley bn[1] = nn[2] * tg[0]; 705de8de29fSMatthew G. Knepley bn[2] = nn[0] * tg[1] - nn[1] * tg[0]; 706de8de29fSMatthew G. Knepley /* Computing 2nd power */ 707de8de29fSMatthew G. Knepley d__1 = bn[0]; 708de8de29fSMatthew G. Knepley /* Computing 2nd power */ 709de8de29fSMatthew G. Knepley d__2 = bn[1]; 710de8de29fSMatthew G. Knepley /* Computing 2nd power */ 711de8de29fSMatthew G. Knepley d__3 = bn[2]; 712de8de29fSMatthew G. Knepley tmp = PetscSqrtScalar(d__1 * d__1 + d__2 * d__2 + d__3 * d__3); 713de8de29fSMatthew G. Knepley i__1 = ndim; 714de8de29fSMatthew G. Knepley for (k = 1; k <= i__1; ++k) bn[k - 1] /= tmp; 715de8de29fSMatthew G. Knepley } else if (ndim == 2) { 716de8de29fSMatthew G. Knepley tg[0] = -nn[1]; 717de8de29fSMatthew G. Knepley tg[1] = nn[0]; 718de8de29fSMatthew G. Knepley /* tmp=dsqrt(tg(1)**2+tg(2)**2) */ 719de8de29fSMatthew G. Knepley /* tg=tg/tmp */ 720de8de29fSMatthew G. Knepley bn[0] = 0.; 721de8de29fSMatthew G. Knepley bn[1] = 0.; 722de8de29fSMatthew G. Knepley bn[2] = 1.; 723de8de29fSMatthew G. Knepley } 724de8de29fSMatthew G. Knepley rl = ul[0]; 725de8de29fSMatthew G. Knepley rr = ur[0]; 726de8de29fSMatthew G. Knepley uxl = 0.; 727de8de29fSMatthew G. Knepley uxr = 0.; 728de8de29fSMatthew G. Knepley utl = 0.; 729de8de29fSMatthew G. Knepley utr = 0.; 730de8de29fSMatthew G. Knepley ubl = 0.; 731de8de29fSMatthew G. Knepley ubr = 0.; 732de8de29fSMatthew G. Knepley i__1 = ndim; 733de8de29fSMatthew G. Knepley for (k = 1; k <= i__1; ++k) { 734de8de29fSMatthew G. Knepley uxl += ul[k] * nn[k - 1]; 735de8de29fSMatthew G. Knepley uxr += ur[k] * nn[k - 1]; 736de8de29fSMatthew G. Knepley utl += ul[k] * tg[k - 1]; 737de8de29fSMatthew G. Knepley utr += ur[k] * tg[k - 1]; 738de8de29fSMatthew G. Knepley ubl += ul[k] * bn[k - 1]; 739de8de29fSMatthew G. Knepley ubr += ur[k] * bn[k - 1]; 740de8de29fSMatthew G. Knepley } 741de8de29fSMatthew G. Knepley uxl /= rl; 742de8de29fSMatthew G. Knepley uxr /= rr; 743de8de29fSMatthew G. Knepley utl /= rl; 744de8de29fSMatthew G. Knepley utr /= rr; 745de8de29fSMatthew G. Knepley ubl /= rl; 746de8de29fSMatthew G. Knepley ubr /= rr; 747de8de29fSMatthew G. Knepley 748de8de29fSMatthew G. Knepley gaml = gamma; 749de8de29fSMatthew G. Knepley gamr = gamma; 750de8de29fSMatthew G. Knepley /* Computing 2nd power */ 751de8de29fSMatthew G. Knepley d__1 = uxl; 752de8de29fSMatthew G. Knepley /* Computing 2nd power */ 753de8de29fSMatthew G. Knepley d__2 = utl; 754de8de29fSMatthew G. Knepley /* Computing 2nd power */ 755de8de29fSMatthew G. Knepley d__3 = ubl; 756de8de29fSMatthew G. Knepley pl = (gamma - 1.) * (ul[ndim + 1] - rl * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3)); 757de8de29fSMatthew G. Knepley /* Computing 2nd power */ 758de8de29fSMatthew G. Knepley d__1 = uxr; 759de8de29fSMatthew G. Knepley /* Computing 2nd power */ 760de8de29fSMatthew G. Knepley d__2 = utr; 761de8de29fSMatthew G. Knepley /* Computing 2nd power */ 762de8de29fSMatthew G. Knepley d__3 = ubr; 763de8de29fSMatthew G. Knepley pr = (gamma - 1.) * (ur[ndim + 1] - rr * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3)); 764de8de29fSMatthew G. Knepley rho1l = rl; 765de8de29fSMatthew G. Knepley rho1r = rr; 766de8de29fSMatthew G. Knepley 767de8de29fSMatthew G. Knepley iwave = riemannsolver(&xcen, &xp, &dtt, &rl, &uxl, &pl, &utl, &ubl, &gaml, &rho1l, &rr, &uxr, &pr, &utr, &ubr, &gamr, &rho1r, &rhom, &unm, &pm, &utm, &ubm, &gamm, &rho1m); 768de8de29fSMatthew G. Knepley 769de8de29fSMatthew G. Knepley flux[0] = rhom * unm; 770de8de29fSMatthew G. Knepley fn = rhom * unm * unm + pm; 771de8de29fSMatthew G. Knepley ft = rhom * unm * utm; 772de8de29fSMatthew G. Knepley /* flux(2)=fn*nn(1)+ft*nn(2) */ 773de8de29fSMatthew G. Knepley /* flux(3)=fn*tg(1)+ft*tg(2) */ 774de8de29fSMatthew G. Knepley flux[1] = fn * nn[0] + ft * tg[0]; 775de8de29fSMatthew G. Knepley flux[2] = fn * nn[1] + ft * tg[1]; 776de8de29fSMatthew G. Knepley /* flux(2)=rhom*unm*(unm)+pm */ 777de8de29fSMatthew G. Knepley /* flux(3)=rhom*(unm)*utm */ 778de8de29fSMatthew G. Knepley if (ndim == 3) flux[3] = rhom * unm * ubm; 779de8de29fSMatthew G. Knepley flux[ndim + 1] = (rhom * .5 * (unm * unm + utm * utm + ubm * ubm) + gamm / (gamm - 1.) * pm) * unm; 780de8de29fSMatthew G. Knepley return iwave; 781de8de29fSMatthew G. Knepley } /* godunovflux_ */ 782de8de29fSMatthew G. Knepley 783de8de29fSMatthew G. Knepley /* PetscReal* => EulerNode* conversion */ 784de8de29fSMatthew G. Knepley static void PhysicsRiemann_Euler_Godunov(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, Physics phys) 785de8de29fSMatthew G. Knepley { 786de8de29fSMatthew G. Knepley Physics_Euler *eu = (Physics_Euler *)phys->data; 787de8de29fSMatthew G. Knepley const PetscReal gamma = eu->gamma; 788a7d931c6SMatthew G. Knepley PetscReal zero = 0.; 789de8de29fSMatthew G. Knepley PetscReal cL, cR, speed, velL, velR, nn[DIM], s2; 790de8de29fSMatthew G. Knepley PetscInt i; 791de8de29fSMatthew G. Knepley PetscErrorCode ierr; 792de8de29fSMatthew G. Knepley 793de8de29fSMatthew G. Knepley PetscFunctionBeginUser; 794de8de29fSMatthew G. Knepley for (i = 0, s2 = 0.; i < DIM; i++) { 795de8de29fSMatthew G. Knepley nn[i] = n[i]; 796de8de29fSMatthew G. Knepley s2 += nn[i] * nn[i]; 797de8de29fSMatthew G. Knepley } 798de8de29fSMatthew G. Knepley s2 = PetscSqrtReal(s2); /* |n|_2 = sum(n^2)^1/2 */ 799de8de29fSMatthew G. Knepley for (i = 0.; i < DIM; i++) nn[i] /= s2; 800de8de29fSMatthew G. Knepley if (0) { /* Rusanov */ 801de8de29fSMatthew G. Knepley const EulerNode *uL = (const EulerNode *)xL, *uR = (const EulerNode *)xR; 802de8de29fSMatthew G. Knepley EulerNodeUnion fL, fR; 803a7d931c6SMatthew G. Knepley ierr = EulerFlux(phys, nn, uL, &fL.eulernode); 804a7d931c6SMatthew G. Knepley if (ierr) { 805a7d931c6SMatthew G. Knepley PetscCallVoid(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 806a7d931c6SMatthew G. Knepley for (i = 0; i < 2 + dim; i++) fL.vals[i] = zero / zero; 807a7d931c6SMatthew G. Knepley PetscCallVoid(PetscFPTrapPop()); 808a7d931c6SMatthew G. Knepley } 809a7d931c6SMatthew G. Knepley ierr = EulerFlux(phys, nn, uR, &fR.eulernode); 810a7d931c6SMatthew G. Knepley if (ierr) { 811a7d931c6SMatthew G. Knepley PetscCallVoid(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 812a7d931c6SMatthew G. Knepley for (i = 0; i < 2 + dim; i++) fR.vals[i] = zero / zero; 813a7d931c6SMatthew G. Knepley PetscCallVoid(PetscFPTrapPop()); 814a7d931c6SMatthew G. Knepley } 815de8de29fSMatthew G. Knepley ierr = SpeedOfSound_PG(gamma, uL, &cL); 816a7d931c6SMatthew G. Knepley if (ierr) { 817a7d931c6SMatthew G. Knepley PetscCallVoid(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 818a7d931c6SMatthew G. Knepley cL = zero / zero; 819a7d931c6SMatthew G. Knepley PetscCallVoid(PetscFPTrapPop()); 820a7d931c6SMatthew G. Knepley } 821de8de29fSMatthew G. Knepley ierr = SpeedOfSound_PG(gamma, uR, &cR); 822a7d931c6SMatthew G. Knepley if (ierr) { 823a7d931c6SMatthew G. Knepley PetscCallVoid(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 824a7d931c6SMatthew G. Knepley cR = zero / zero; 825a7d931c6SMatthew G. Knepley PetscCallVoid(PetscFPTrapPop()); 826a7d931c6SMatthew G. Knepley } 827de8de29fSMatthew G. Knepley velL = DotDIMReal(uL->ru, nn) / uL->r; 828de8de29fSMatthew G. Knepley velR = DotDIMReal(uR->ru, nn) / uR->r; 829de8de29fSMatthew G. Knepley speed = PetscMax(velR + cR, velL + cL); 830de8de29fSMatthew G. Knepley for (i = 0; i < 2 + dim; i++) flux[i] = 0.5 * ((fL.vals[i] + fR.vals[i]) + speed * (xL[i] - xR[i])) * s2; 831de8de29fSMatthew G. Knepley } else { 832de8de29fSMatthew G. Knepley /* int iwave = */ 833de8de29fSMatthew G. Knepley godunovflux(xL, xR, flux, nn, DIM, gamma); 834de8de29fSMatthew G. Knepley for (i = 0; i < 2 + dim; i++) flux[i] *= s2; 835de8de29fSMatthew G. Knepley } 836de8de29fSMatthew G. Knepley PetscFunctionReturnVoid(); 837de8de29fSMatthew G. Knepley } 838de8de29fSMatthew G. Knepley 839de8de29fSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED 840*2a8381b2SBarry Smith CEED_QFUNCTION(PhysicsRiemann_Euler_Godunov_CEED)(PetscCtx ctx, CeedInt Q, const CeedScalar *const in[], CeedScalar *const out[]) 841de8de29fSMatthew G. Knepley { 842de8de29fSMatthew G. Knepley const CeedScalar *xL = in[0], *xR = in[1], *geom = in[2]; 843de8de29fSMatthew G. Knepley CeedScalar *cL = out[0], *cR = out[1]; 844de8de29fSMatthew G. Knepley const Physics_Euler *eu = (Physics_Euler *)ctx; 845de8de29fSMatthew G. Knepley struct _n_Physics phys; 846de8de29fSMatthew G. Knepley 847de8de29fSMatthew G. Knepley phys.data = (void *)eu; 848de8de29fSMatthew G. Knepley CeedPragmaSIMD for (CeedInt i = 0; i < Q; ++i) 849de8de29fSMatthew G. Knepley { 850de8de29fSMatthew G. Knepley const CeedScalar qL[DIM + 2] = {xL[i + Q * 0], xL[i + Q * 1], xL[i + Q * 2], xL[i + Q * 3]}; 851de8de29fSMatthew G. Knepley const CeedScalar qR[DIM + 2] = {xR[i + Q * 0], xR[i + Q * 1], xR[i + Q * 2], xR[i + Q * 3]}; 852de8de29fSMatthew G. Knepley const CeedScalar n[DIM] = {geom[i + Q * 0], geom[i + Q * 1]}; 853de8de29fSMatthew G. Knepley CeedScalar flux[DIM + 2]; 854de8de29fSMatthew G. Knepley 855de8de29fSMatthew G. Knepley #if 0 856de8de29fSMatthew G. Knepley PetscPrintf(PETSC_COMM_SELF, "Cell %d Normal\n", 0); 8573a7d0413SPierre Jolivet for (CeedInt j = 0; j < DIM; ++j) PetscPrintf(PETSC_COMM_SELF, " | %g |\n", n[j]); 858de8de29fSMatthew G. Knepley PetscPrintf(PETSC_COMM_SELF, "Cell %d Element Residual: left state\n", 0); 8593a7d0413SPierre Jolivet for (CeedInt j = 0; j < DIM + 2; ++j) PetscPrintf(PETSC_COMM_SELF, " | %g |\n", qL[j]); 860de8de29fSMatthew G. Knepley PetscPrintf(PETSC_COMM_SELF, "Cell %d Element Residual: right state\n", 0); 8613a7d0413SPierre Jolivet for (CeedInt j = 0; j < DIM + 2; ++j) PetscPrintf(PETSC_COMM_SELF, " | %g |\n", qR[j]); 862de8de29fSMatthew G. Knepley #endif 863de8de29fSMatthew G. Knepley PhysicsRiemann_Euler_Godunov(DIM, DIM + 2, NULL, n, qL, qR, 0, NULL, flux, &phys); 864de8de29fSMatthew G. Knepley for (CeedInt j = 0; j < DIM + 2; ++j) { 865de8de29fSMatthew G. Knepley cL[i + Q * j] = -flux[j] / geom[i + Q * 2]; 866de8de29fSMatthew G. Knepley cR[i + Q * j] = flux[j] / geom[i + Q * 3]; 867de8de29fSMatthew G. Knepley } 868de8de29fSMatthew G. Knepley #if 0 869de8de29fSMatthew G. Knepley PetscPrintf(PETSC_COMM_SELF, "Cell %d Element Residual: left flux\n", 0); 8703a7d0413SPierre Jolivet for (CeedInt j = 0; j < DIM + 2; ++j) PetscPrintf(PETSC_COMM_SELF, " | %g | (%g)\n", cL[i + Q * j], geom[i + Q * 2]); 871de8de29fSMatthew G. Knepley PetscPrintf(PETSC_COMM_SELF, "Cell %d Element Residual: right flux\n", 0); 8723a7d0413SPierre Jolivet for (CeedInt j = 0; j < DIM + 2; ++j) PetscPrintf(PETSC_COMM_SELF, " | %g | (%g)\n", cR[i + Q * j], geom[i + Q * 3]); 873de8de29fSMatthew G. Knepley #endif 874de8de29fSMatthew G. Knepley } 875de8de29fSMatthew G. Knepley return CEED_ERROR_SUCCESS; 876de8de29fSMatthew G. Knepley } 877de8de29fSMatthew G. Knepley #endif 878