static char help[] = "Second Order TVD Finite Volume Example.\n"; /*F We use a second order TVD finite volume method to evolve a system of PDEs. Our simple upwinded residual evaluation loops over all mesh faces and uses a Riemann solver to produce the flux given the face geometry and cell values, \begin{equation} f_i = \mathrm{riemann}(\mathrm{phys}, p_\mathrm{centroid}, \hat n, x^L, x^R) \end{equation} and then update the cell values given the cell volume. \begin{eqnarray} f^L_i &-=& \frac{f_i}{vol^L} \\ f^R_i &+=& \frac{f_i}{vol^R} \end{eqnarray} As an example, we can consider the shallow water wave equation, \begin{eqnarray} h_t + \nabla\cdot \left( uh \right) &=& 0 \\ (uh)_t + \nabla\cdot \left( u\otimes uh + \frac{g h^2}{2} I \right) &=& 0 \end{eqnarray} where $h$ is wave height, $u$ is wave velocity, and $g$ is the acceleration due to gravity. A representative Riemann solver for the shallow water equations is given in the PhysicsRiemann_SW() function, \begin{eqnarray} f^{L,R}_h &=& uh^{L,R} \cdot \hat n \\ f^{L,R}_{uh} &=& \frac{f^{L,R}_h}{h^{L,R}} uh^{L,R} + g (h^{L,R})^2 \hat n \\ c^{L,R} &=& \sqrt{g h^{L,R}} \\ s &=& \max\left( \left|\frac{uh^L \cdot \hat n}{h^L}\right| + c^L, \left|\frac{uh^R \cdot \hat n}{h^R}\right| + c^R \right) \\ f_i &=& \frac{A_\mathrm{face}}{2} \left( f^L_i + f^R_i + s \left( x^L_i - x^R_i \right) \right) \end{eqnarray} where $c$ is the local gravity wave speed and $f_i$ is a Rusanov flux. The more sophisticated residual evaluation in RHSFunctionLocal_LS() uses a least-squares fit to a quadratic polynomial over a neighborhood of the given element. The mesh is read in from an ExodusII file, usually generated by Cubit. The example also shows how to handle AMR in a time-dependent TS solver. F*/ #include #include #include #include #include "ex11.h" static PetscFunctionList PhysicsList, PhysicsRiemannList_SW, PhysicsRiemannList_Euler; /* 'User' implements a discretization of a continuous model. */ typedef struct _n_User *User; typedef PetscErrorCode (*SolutionFunction)(Model, PetscReal, const PetscReal *, PetscScalar *, void *); typedef PetscErrorCode (*SetUpBCFunction)(DM, PetscDS, Physics); typedef PetscErrorCode (*FunctionalFunction)(Model, PetscReal, const PetscReal *, const PetscScalar *, PetscReal *, void *); typedef PetscErrorCode (*SetupFields)(Physics, PetscSection); static PetscErrorCode ModelSolutionSetDefault(Model, SolutionFunction, void *); static PetscErrorCode ModelFunctionalRegister(Model, const char *, PetscInt *, FunctionalFunction, void *); static PetscErrorCode OutputVTK(DM, const char *, PetscViewer *); typedef struct _n_FunctionalLink *FunctionalLink; struct _n_FunctionalLink { char *name; FunctionalFunction func; void *ctx; PetscInt offset; FunctionalLink next; }; struct _n_Model { MPI_Comm comm; /* Does not do collective communication, but some error conditions can be collective */ Physics physics; FunctionalLink functionalRegistry; PetscInt maxComputed; PetscInt numMonitored; FunctionalLink *functionalMonitored; PetscInt numCall; FunctionalLink *functionalCall; SolutionFunction solution; SetUpBCFunction setupbc; void *solutionctx; PetscReal maxspeed; /* estimate of global maximum speed (for CFL calculation) */ PetscReal bounds[2 * DIM]; PetscErrorCode (*errorIndicator)(PetscInt, PetscReal, PetscInt, const PetscScalar[], const PetscScalar[], PetscReal *, void *); void *errorCtx; PetscErrorCode (*setupCEED)(DM, Physics); }; struct _n_User { PetscInt vtkInterval; /* For monitor */ char outputBasename[PETSC_MAX_PATH_LEN]; /* Basename for output files */ PetscInt monitorStepOffset; Model model; PetscBool vtkmon; }; #ifdef PETSC_HAVE_LIBCEED // Free a plain data context that was allocated using PETSc; returning libCEED error codes static int FreeContextPetsc(void *data) { if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed"); return CEED_ERROR_SUCCESS; } #endif /******************* Advect ********************/ typedef enum { ADVECT_SOL_TILTED, ADVECT_SOL_BUMP, ADVECT_SOL_BUMP_CAVITY } AdvectSolType; static const char *const AdvectSolTypes[] = {"TILTED", "BUMP", "BUMP_CAVITY", "AdvectSolType", "ADVECT_SOL_", 0}; typedef enum { ADVECT_SOL_BUMP_CONE, ADVECT_SOL_BUMP_COS } AdvectSolBumpType; static const char *const AdvectSolBumpTypes[] = {"CONE", "COS", "AdvectSolBumpType", "ADVECT_SOL_BUMP_", 0}; typedef struct { PetscReal wind[DIM]; } Physics_Advect_Tilted; typedef struct { PetscReal center[DIM]; PetscReal radius; AdvectSolBumpType type; } Physics_Advect_Bump; typedef struct { PetscReal inflowState; AdvectSolType soltype; union { Physics_Advect_Tilted tilted; Physics_Advect_Bump bump; } sol; struct { PetscInt Solution; PetscInt Error; } functional; } Physics_Advect; static const struct FieldDescription PhysicsFields_Advect[] = { {"U", 1}, {NULL, 0} }; static PetscErrorCode PhysicsBoundary_Advect_Inflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, PetscCtx ctx) { Physics phys = (Physics)ctx; Physics_Advect *advect = (Physics_Advect *)phys->data; PetscFunctionBeginUser; xG[0] = advect->inflowState; PetscFunctionReturn(PETSC_SUCCESS); } static PetscErrorCode PhysicsBoundary_Advect_Outflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, PetscCtx ctx) { PetscFunctionBeginUser; xG[0] = xI[0]; PetscFunctionReturn(PETSC_SUCCESS); } static void PhysicsRiemann_Advect(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) { Physics_Advect *advect = (Physics_Advect *)phys->data; PetscReal wind[DIM], wn; switch (advect->soltype) { case ADVECT_SOL_TILTED: { Physics_Advect_Tilted *tilted = &advect->sol.tilted; wind[0] = tilted->wind[0]; wind[1] = tilted->wind[1]; } break; case ADVECT_SOL_BUMP: wind[0] = -qp[1]; wind[1] = qp[0]; break; case ADVECT_SOL_BUMP_CAVITY: { PetscInt i; PetscReal comp2[3] = {0., 0., 0.}, rad2; rad2 = 0.; for (i = 0; i < dim; i++) { comp2[i] = qp[i] * qp[i]; rad2 += comp2[i]; } wind[0] = -qp[1]; wind[1] = qp[0]; if (rad2 > 1.) { PetscInt maxI = 0; PetscReal maxComp2 = comp2[0]; for (i = 1; i < dim; i++) { if (comp2[i] > maxComp2) { maxI = i; maxComp2 = comp2[i]; } } wind[maxI] = 0.; } } break; default: { PetscInt i; for (i = 0; i < DIM; ++i) wind[i] = 0.0; } /* default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for solution type %s",AdvectSolBumpTypes[advect->soltype]); */ } wn = Dot2Real(wind, n); flux[0] = (wn > 0 ? xL[0] : xR[0]) * wn; } static PetscErrorCode PhysicsSolution_Advect(Model mod, PetscReal time, const PetscReal *x, PetscScalar *u, PetscCtx ctx) { Physics phys = (Physics)ctx; Physics_Advect *advect = (Physics_Advect *)phys->data; PetscFunctionBeginUser; switch (advect->soltype) { case ADVECT_SOL_TILTED: { PetscReal x0[DIM]; Physics_Advect_Tilted *tilted = &advect->sol.tilted; Waxpy2Real(-time, tilted->wind, x, x0); if (x0[1] > 0) u[0] = 1. * x[0] + 3. * x[1]; else u[0] = advect->inflowState; } break; case ADVECT_SOL_BUMP_CAVITY: case ADVECT_SOL_BUMP: { Physics_Advect_Bump *bump = &advect->sol.bump; PetscReal x0[DIM], v[DIM], r, cost, sint; cost = PetscCosReal(time); sint = PetscSinReal(time); x0[0] = cost * x[0] + sint * x[1]; x0[1] = -sint * x[0] + cost * x[1]; Waxpy2Real(-1, bump->center, x0, v); r = Norm2Real(v); switch (bump->type) { case ADVECT_SOL_BUMP_CONE: u[0] = PetscMax(1 - r / bump->radius, 0); break; case ADVECT_SOL_BUMP_COS: u[0] = 0.5 + 0.5 * PetscCosReal(PetscMin(r / bump->radius, 1) * PETSC_PI); break; } } break; default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unknown solution type"); } PetscFunctionReturn(PETSC_SUCCESS); } static PetscErrorCode PhysicsFunctional_Advect(Model mod, PetscReal time, const PetscReal *x, const PetscScalar *y, PetscReal *f, PetscCtx ctx) { Physics phys = (Physics)ctx; Physics_Advect *advect = (Physics_Advect *)phys->data; PetscScalar yexact[1] = {0.0}; PetscFunctionBeginUser; PetscCall(PhysicsSolution_Advect(mod, time, x, yexact, phys)); f[advect->functional.Solution] = PetscRealPart(y[0]); f[advect->functional.Error] = PetscAbsScalar(y[0] - yexact[0]); PetscFunctionReturn(PETSC_SUCCESS); } static PetscErrorCode SetUpBC_Advect(DM dm, PetscDS prob, Physics phys) { const PetscInt inflowids[] = {100, 200, 300}, outflowids[] = {101}; DMLabel label; PetscFunctionBeginUser; /* Register "canned" boundary conditions and defaults for where to apply. */ PetscCall(DMGetLabel(dm, "Face Sets", &label)); PetscCall(PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "inflow", label, PETSC_STATIC_ARRAY_LENGTH(inflowids), inflowids, 0, 0, NULL, (PetscVoidFn *)PhysicsBoundary_Advect_Inflow, NULL, phys, NULL)); PetscCall(PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "outflow", label, PETSC_STATIC_ARRAY_LENGTH(outflowids), outflowids, 0, 0, NULL, (PetscVoidFn *)PhysicsBoundary_Advect_Outflow, NULL, phys, NULL)); PetscFunctionReturn(PETSC_SUCCESS); } static PetscErrorCode PhysicsCreate_Advect(Model mod, Physics phys, PetscOptionItems PetscOptionsObject) { Physics_Advect *advect; PetscFunctionBeginUser; phys->field_desc = PhysicsFields_Advect; phys->riemann = (PetscRiemannFn *)(PetscVoidFn *)PhysicsRiemann_Advect; PetscCall(PetscNew(&advect)); phys->data = advect; mod->setupbc = SetUpBC_Advect; PetscOptionsHeadBegin(PetscOptionsObject, "Advect options"); { PetscInt two = 2, dof = 1; advect->soltype = ADVECT_SOL_TILTED; PetscCall(PetscOptionsEnum("-advect_sol_type", "solution type", "", AdvectSolTypes, (PetscEnum)advect->soltype, (PetscEnum *)&advect->soltype, NULL)); switch (advect->soltype) { case ADVECT_SOL_TILTED: { Physics_Advect_Tilted *tilted = &advect->sol.tilted; two = 2; tilted->wind[0] = 0.0; tilted->wind[1] = 1.0; PetscCall(PetscOptionsRealArray("-advect_tilted_wind", "background wind vx,vy", "", tilted->wind, &two, NULL)); advect->inflowState = -2.0; PetscCall(PetscOptionsRealArray("-advect_tilted_inflow", "Inflow state", "", &advect->inflowState, &dof, NULL)); phys->maxspeed = Norm2Real(tilted->wind); } break; case ADVECT_SOL_BUMP_CAVITY: case ADVECT_SOL_BUMP: { Physics_Advect_Bump *bump = &advect->sol.bump; two = 2; bump->center[0] = 2.; bump->center[1] = 0.; PetscCall(PetscOptionsRealArray("-advect_bump_center", "location of center of bump x,y", "", bump->center, &two, NULL)); bump->radius = 0.9; PetscCall(PetscOptionsReal("-advect_bump_radius", "radius of bump", "", bump->radius, &bump->radius, NULL)); bump->type = ADVECT_SOL_BUMP_CONE; PetscCall(PetscOptionsEnum("-advect_bump_type", "type of bump", "", AdvectSolBumpTypes, (PetscEnum)bump->type, (PetscEnum *)&bump->type, NULL)); phys->maxspeed = 3.; /* radius of mesh, kludge */ } break; } } PetscOptionsHeadEnd(); /* Initial/transient solution with default boundary conditions */ PetscCall(ModelSolutionSetDefault(mod, PhysicsSolution_Advect, phys)); /* Register "canned" functionals */ PetscCall(ModelFunctionalRegister(mod, "Solution", &advect->functional.Solution, PhysicsFunctional_Advect, phys)); PetscCall(ModelFunctionalRegister(mod, "Error", &advect->functional.Error, PhysicsFunctional_Advect, phys)); PetscFunctionReturn(PETSC_SUCCESS); } /******************* Shallow Water ********************/ static const struct FieldDescription PhysicsFields_SW[] = { {"Height", 1 }, {"Momentum", DIM}, {NULL, 0 } }; static PetscErrorCode PhysicsBoundary_SW_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, PetscCtx ctx) { PetscFunctionBeginUser; xG[0] = xI[0]; xG[1] = -xI[1]; xG[2] = -xI[2]; PetscFunctionReturn(PETSC_SUCCESS); } static PetscErrorCode PhysicsSolution_SW(Model mod, PetscReal time, const PetscReal *x, PetscScalar *u, PetscCtx ctx) { PetscReal dx[2], r, sigma; PetscFunctionBeginUser; PetscCheck(time == 0.0, mod->comm, PETSC_ERR_SUP, "No solution known for time %g", (double)time); dx[0] = x[0] - 1.5; dx[1] = x[1] - 1.0; r = Norm2Real(dx); sigma = 0.5; u[0] = 1 + 2 * PetscExpReal(-PetscSqr(r) / (2 * PetscSqr(sigma))); u[1] = 0.0; u[2] = 0.0; PetscFunctionReturn(PETSC_SUCCESS); } static PetscErrorCode PhysicsFunctional_SW(Model mod, PetscReal time, const PetscReal *coord, const PetscScalar *xx, PetscReal *f, PetscCtx ctx) { Physics phys = (Physics)ctx; Physics_SW *sw = (Physics_SW *)phys->data; const SWNode *x = (const SWNode *)xx; PetscReal u[2]; PetscReal h; PetscFunctionBeginUser; h = x->h; Scale2Real(1. / x->h, x->uh, u); f[sw->functional.Height] = h; f[sw->functional.Speed] = Norm2Real(u) + PetscSqrtReal(sw->gravity * h); f[sw->functional.Energy] = 0.5 * (Dot2Real(x->uh, u) + sw->gravity * PetscSqr(h)); PetscFunctionReturn(PETSC_SUCCESS); } static PetscErrorCode SetUpBC_SW(DM dm, PetscDS prob, Physics phys) { const PetscInt wallids[] = {100, 101, 200, 300}; DMLabel label; PetscFunctionBeginUser; PetscCall(DMGetLabel(dm, "Face Sets", &label)); PetscCall(PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", label, PETSC_STATIC_ARRAY_LENGTH(wallids), wallids, 0, 0, NULL, (PetscVoidFn *)PhysicsBoundary_SW_Wall, NULL, phys, NULL)); PetscFunctionReturn(PETSC_SUCCESS); } #ifdef PETSC_HAVE_LIBCEED static PetscErrorCode CreateQFunctionContext_SW(Physics phys, Ceed ceed, CeedQFunctionContext *qfCtx) { Physics_SW *in = (Physics_SW *)phys->data; Physics_SW *sw; PetscFunctionBeginUser; PetscCall(PetscCalloc1(1, &sw)); sw->gravity = in->gravity; PetscCallCEED(CeedQFunctionContextCreate(ceed, qfCtx)); PetscCallCEED(CeedQFunctionContextSetData(*qfCtx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*sw), sw)); PetscCallCEED(CeedQFunctionContextSetDataDestroy(*qfCtx, CEED_MEM_HOST, FreeContextPetsc)); PetscCallCEED(CeedQFunctionContextRegisterDouble(*qfCtx, "gravity", offsetof(Physics_SW, gravity), 1, "Acceleration due to gravity")); PetscFunctionReturn(PETSC_SUCCESS); } #endif static PetscErrorCode SetupCEED_SW(DM dm, Physics physics) { #ifdef PETSC_HAVE_LIBCEED Ceed ceed; CeedQFunctionContext qfCtx; #endif PetscFunctionBegin; #ifdef PETSC_HAVE_LIBCEED PetscCall(DMGetCeed(dm, &ceed)); PetscCall(CreateQFunctionContext_SW(physics, ceed, &qfCtx)); PetscCall(DMCeedCreateFVM(dm, PETSC_TRUE, PhysicsRiemann_SW_Rusanov_CEED, PhysicsRiemann_SW_Rusanov_CEED_loc, qfCtx)); PetscCallCEED(CeedQFunctionContextDestroy(&qfCtx)); #endif PetscFunctionReturn(PETSC_SUCCESS); } static PetscErrorCode PhysicsCreate_SW(Model mod, Physics phys, PetscOptionItems PetscOptionsObject) { Physics_SW *sw; char sw_riemann[64] = "rusanov"; PetscFunctionBeginUser; phys->field_desc = PhysicsFields_SW; PetscCall(PetscNew(&sw)); phys->data = sw; mod->setupbc = SetUpBC_SW; mod->setupCEED = SetupCEED_SW; PetscCall(PetscFunctionListAdd(&PhysicsRiemannList_SW, "rusanov", PhysicsRiemann_SW_Rusanov)); PetscCall(PetscFunctionListAdd(&PhysicsRiemannList_SW, "hll", PhysicsRiemann_SW_HLL)); #ifdef PETSC_HAVE_LIBCEED PetscCall(PetscFunctionListAdd(&PhysicsRiemannList_SW, "rusanov_ceed", PhysicsRiemann_SW_Rusanov_CEED)); #endif PetscOptionsHeadBegin(PetscOptionsObject, "SW options"); { void (*PhysicsRiemann_SW)(PetscInt, PetscInt, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, PetscInt, const PetscScalar, PetscScalar *, Physics); sw->gravity = 1.0; PetscCall(PetscOptionsReal("-sw_gravity", "Gravitational constant", "", sw->gravity, &sw->gravity, NULL)); PetscCall(PetscOptionsFList("-sw_riemann", "Riemann solver", "", PhysicsRiemannList_SW, sw_riemann, sw_riemann, sizeof sw_riemann, NULL)); PetscCall(PetscFunctionListFind(PhysicsRiemannList_SW, sw_riemann, &PhysicsRiemann_SW)); phys->riemann = (PetscRiemannFn *)(PetscVoidFn *)PhysicsRiemann_SW; } PetscOptionsHeadEnd(); phys->maxspeed = PetscSqrtReal(2.0 * sw->gravity); /* Mach 1 for depth of 2 */ PetscCall(ModelSolutionSetDefault(mod, PhysicsSolution_SW, phys)); PetscCall(ModelFunctionalRegister(mod, "Height", &sw->functional.Height, PhysicsFunctional_SW, phys)); PetscCall(ModelFunctionalRegister(mod, "Speed", &sw->functional.Speed, PhysicsFunctional_SW, phys)); PetscCall(ModelFunctionalRegister(mod, "Energy", &sw->functional.Energy, PhysicsFunctional_SW, phys)); PetscFunctionReturn(PETSC_SUCCESS); } /******************* Euler Density Shock (EULER_IV_SHOCK,EULER_SS_SHOCK) ********************/ /* An initial-value and self-similar solutions of the compressible Euler equations */ /* Ravi Samtaney and D. I. Pullin */ /* Phys. Fluids 8, 2650 (1996); http://dx.doi.org/10.1063/1.869050 */ static const struct FieldDescription PhysicsFields_Euler[] = { {"Density", 1 }, {"Momentum", DIM}, {"Energy", 1 }, {NULL, 0 } }; /* initial condition */ int initLinearWave(EulerNode *ux, const PetscReal gamma, const PetscReal coord[], const PetscReal Lx); static PetscErrorCode PhysicsSolution_Euler(Model mod, PetscReal time, const PetscReal *x, PetscScalar *u, PetscCtx ctx) { PetscInt i; Physics phys = (Physics)ctx; Physics_Euler *eu = (Physics_Euler *)phys->data; EulerNode *uu = (EulerNode *)u; PetscReal p0, gamma, c = 0.; PetscFunctionBeginUser; PetscCheck(time == 0.0, mod->comm, PETSC_ERR_SUP, "No solution known for time %g", (double)time); for (i = 0; i < DIM; i++) uu->ru[i] = 0.0; /* zero out initial velocity */ /* set E and rho */ gamma = eu->gamma; if (eu->type == EULER_IV_SHOCK || eu->type == EULER_SS_SHOCK) { /******************* Euler Density Shock ********************/ /* On initial-value and self-similar solutions of the compressible Euler equations */ /* Ravi Samtaney and D. I. Pullin */ /* Phys. Fluids 8, 2650 (1996); http://dx.doi.org/10.1063/1.869050 */ /* initial conditions 1: left of shock, 0: left of discontinuity 2: right of discontinuity, */ p0 = 1.; if (x[0] < 0.0 + x[1] * eu->itana) { if (x[0] < mod->bounds[0] * 0.5) { /* left of shock (1) */ PetscReal amach, rho, press, gas1, p1; amach = eu->amach; rho = 1.; press = p0; p1 = press * (1.0 + 2.0 * gamma / (gamma + 1.0) * (amach * amach - 1.0)); gas1 = (gamma - 1.0) / (gamma + 1.0); uu->r = rho * (p1 / press + gas1) / (gas1 * p1 / press + 1.0); uu->ru[0] = ((uu->r - rho) * PetscSqrtReal(gamma * press / rho) * amach); uu->E = p1 / (gamma - 1.0) + .5 / uu->r * uu->ru[0] * uu->ru[0]; } else { /* left of discontinuity (0) */ uu->r = 1.; /* rho = 1 */ uu->E = p0 / (gamma - 1.0); } } else { /* right of discontinuity (2) */ uu->r = eu->rhoR; uu->E = p0 / (gamma - 1.0); } } else if (eu->type == EULER_SHOCK_TUBE) { /* For (xx0) set (rho,u,p)=(1,0,1). Choose x0 to the midpoint of the domain in the x-direction. */ if (x[0] < 0.0) { uu->r = 8.; uu->E = 10. / (gamma - 1.); } else { uu->r = 1.; uu->E = 1. / (gamma - 1.); } } else if (eu->type == EULER_LINEAR_WAVE) { initLinearWave(uu, gamma, x, mod->bounds[1] - mod->bounds[0]); } else SETERRQ(mod->comm, PETSC_ERR_SUP, "Unknown type %d", eu->type); /* set phys->maxspeed: (mod->maxspeed = phys->maxspeed) in main; */ PetscCall(SpeedOfSound_PG(gamma, uu, &c)); c = (uu->ru[0] / uu->r) + c; if (c > phys->maxspeed) phys->maxspeed = c; PetscFunctionReturn(PETSC_SUCCESS); } /* PetscReal* => EulerNode* conversion */ static PetscErrorCode PhysicsBoundary_Euler_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *a_xI, PetscScalar *a_xG, PetscCtx ctx) { PetscInt i; const EulerNode *xI = (const EulerNode *)a_xI; EulerNode *xG = (EulerNode *)a_xG; Physics phys = (Physics)ctx; Physics_Euler *eu = (Physics_Euler *)phys->data; PetscFunctionBeginUser; xG->r = xI->r; /* ghost cell density - same */ xG->E = xI->E; /* ghost cell energy - same */ if (n[1] != 0.) { /* top and bottom */ xG->ru[0] = xI->ru[0]; /* copy tang to wall */ xG->ru[1] = -xI->ru[1]; /* reflect perp to t/b wall */ } else { /* sides */ for (i = 0; i < DIM; i++) xG->ru[i] = xI->ru[i]; /* copy */ } if (eu->type == EULER_LINEAR_WAVE) { /* debug */ #if 0 PetscPrintf(PETSC_COMM_WORLD,"%s coord=%g,%g\n",PETSC_FUNCTION_NAME,(double)c[0],(double)c[1]); #endif } PetscFunctionReturn(PETSC_SUCCESS); } static PetscErrorCode PhysicsFunctional_Euler(Model mod, PetscReal time, const PetscReal *coord, const PetscScalar *xx, PetscReal *f, PetscCtx ctx) { Physics phys = (Physics)ctx; Physics_Euler *eu = (Physics_Euler *)phys->data; const EulerNode *x = (const EulerNode *)xx; PetscReal p; PetscFunctionBeginUser; f[eu->monitor.Density] = x->r; f[eu->monitor.Momentum] = NormDIM(x->ru); f[eu->monitor.Energy] = x->E; f[eu->monitor.Speed] = NormDIM(x->ru) / x->r; PetscCall(Pressure_PG(eu->gamma, x, &p)); f[eu->monitor.Pressure] = p; PetscFunctionReturn(PETSC_SUCCESS); } static PetscErrorCode SetUpBC_Euler(DM dm, PetscDS prob, Physics phys) { Physics_Euler *eu = (Physics_Euler *)phys->data; DMLabel label; PetscFunctionBeginUser; PetscCall(DMGetLabel(dm, "Face Sets", &label)); if (eu->type == EULER_LINEAR_WAVE) { const PetscInt wallids[] = {100, 101}; PetscCall(PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", label, PETSC_STATIC_ARRAY_LENGTH(wallids), wallids, 0, 0, NULL, (PetscVoidFn *)PhysicsBoundary_Euler_Wall, NULL, phys, NULL)); } else { const PetscInt wallids[] = {100, 101, 200, 300}; PetscCall(PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", label, PETSC_STATIC_ARRAY_LENGTH(wallids), wallids, 0, 0, NULL, (PetscVoidFn *)PhysicsBoundary_Euler_Wall, NULL, phys, NULL)); } PetscFunctionReturn(PETSC_SUCCESS); } #ifdef PETSC_HAVE_LIBCEED static PetscErrorCode CreateQFunctionContext_Euler(Physics phys, Ceed ceed, CeedQFunctionContext *qfCtx) { Physics_Euler *in = (Physics_Euler *)phys->data; Physics_Euler *eu; PetscFunctionBeginUser; PetscCall(PetscCalloc1(1, &eu)); eu->gamma = in->gamma; PetscCallCEED(CeedQFunctionContextCreate(ceed, qfCtx)); PetscCallCEED(CeedQFunctionContextSetData(*qfCtx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*eu), eu)); PetscCallCEED(CeedQFunctionContextSetDataDestroy(*qfCtx, CEED_MEM_HOST, FreeContextPetsc)); PetscCallCEED(CeedQFunctionContextRegisterDouble(*qfCtx, "gamma", offsetof(Physics_Euler, gamma), 1, "Heat capacity ratio")); PetscFunctionReturn(PETSC_SUCCESS); } #endif static PetscErrorCode SetupCEED_Euler(DM dm, Physics physics) { #ifdef PETSC_HAVE_LIBCEED Ceed ceed; CeedQFunctionContext qfCtx; #endif PetscFunctionBegin; #ifdef PETSC_HAVE_LIBCEED PetscCall(DMGetCeed(dm, &ceed)); PetscCall(CreateQFunctionContext_Euler(physics, ceed, &qfCtx)); PetscCall(DMCeedCreateFVM(dm, PETSC_TRUE, PhysicsRiemann_Euler_Godunov_CEED, PhysicsRiemann_Euler_Godunov_CEED_loc, qfCtx)); PetscCallCEED(CeedQFunctionContextDestroy(&qfCtx)); #endif PetscFunctionReturn(PETSC_SUCCESS); } static PetscErrorCode PhysicsCreate_Euler(Model mod, Physics phys, PetscOptionItems PetscOptionsObject) { Physics_Euler *eu; PetscFunctionBeginUser; phys->field_desc = PhysicsFields_Euler; phys->riemann = (PetscRiemannFn *)(PetscVoidFn *)PhysicsRiemann_Euler_Godunov; PetscCall(PetscNew(&eu)); phys->data = eu; mod->setupbc = SetUpBC_Euler; mod->setupCEED = SetupCEED_Euler; PetscCall(PetscFunctionListAdd(&PhysicsRiemannList_Euler, "godunov", PhysicsRiemann_Euler_Godunov)); #ifdef PETSC_HAVE_LIBCEED PetscCall(PetscFunctionListAdd(&PhysicsRiemannList_Euler, "godunov_ceed", PhysicsRiemann_Euler_Godunov_CEED)); #endif PetscOptionsHeadBegin(PetscOptionsObject, "Euler options"); { void (*PhysicsRiemann_Euler)(PetscInt, PetscInt, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, PetscInt, const PetscScalar, PetscScalar *, Physics); PetscReal alpha; char type[64] = "linear_wave"; char eu_riemann[64] = "godunov"; PetscBool is; eu->gamma = 1.4; eu->amach = 2.02; eu->rhoR = 3.0; eu->itana = 0.57735026918963; /* angle of Euler self similar (SS) shock */ PetscCall(PetscOptionsFList("-eu_riemann", "Riemann solver", "", PhysicsRiemannList_Euler, eu_riemann, eu_riemann, sizeof eu_riemann, NULL)); PetscCall(PetscFunctionListFind(PhysicsRiemannList_Euler, eu_riemann, &PhysicsRiemann_Euler)); phys->riemann = (PetscRiemannFn *)(PetscVoidFn *)PhysicsRiemann_Euler; PetscCall(PetscOptionsReal("-eu_gamma", "Heat capacity ratio", "", eu->gamma, &eu->gamma, NULL)); PetscCall(PetscOptionsReal("-eu_amach", "Shock speed (Mach)", "", eu->amach, &eu->amach, NULL)); PetscCall(PetscOptionsReal("-eu_rho2", "Density right of discontinuity", "", eu->rhoR, &eu->rhoR, NULL)); alpha = 60.; PetscCall(PetscOptionsRangeReal("-eu_alpha", "Angle of discontinuity", "", alpha, &alpha, NULL, 0.0, 90.0)); eu->itana = 1. / PetscTanReal(alpha * PETSC_PI / 180.0); PetscCall(PetscOptionsString("-eu_type", "Type of Euler test", "", type, type, sizeof(type), NULL)); PetscCall(PetscStrcmp(type, "linear_wave", &is)); if (is) { /* Remember this should be periodic */ eu->type = EULER_LINEAR_WAVE; PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s set Euler type: %s\n", PETSC_FUNCTION_NAME, "linear_wave")); } else { PetscCheck(DIM == 2, PETSC_COMM_WORLD, PETSC_ERR_SUP, "DIM must be 2 unless linear wave test %s", type); PetscCall(PetscStrcmp(type, "iv_shock", &is)); if (is) { eu->type = EULER_IV_SHOCK; PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s set Euler type: %s\n", PETSC_FUNCTION_NAME, "iv_shock")); } else { PetscCall(PetscStrcmp(type, "ss_shock", &is)); if (is) { eu->type = EULER_SS_SHOCK; PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s set Euler type: %s\n", PETSC_FUNCTION_NAME, "ss_shock")); } else { PetscCall(PetscStrcmp(type, "shock_tube", &is)); PetscCheck(is, PETSC_COMM_WORLD, PETSC_ERR_SUP, "Unknown Euler type %s", type); eu->type = EULER_SHOCK_TUBE; PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s set Euler type: %s\n", PETSC_FUNCTION_NAME, "shock_tube")); } } } } PetscOptionsHeadEnd(); phys->maxspeed = 0.; /* will get set in solution */ PetscCall(ModelSolutionSetDefault(mod, PhysicsSolution_Euler, phys)); PetscCall(ModelFunctionalRegister(mod, "Speed", &eu->monitor.Speed, PhysicsFunctional_Euler, phys)); PetscCall(ModelFunctionalRegister(mod, "Energy", &eu->monitor.Energy, PhysicsFunctional_Euler, phys)); PetscCall(ModelFunctionalRegister(mod, "Density", &eu->monitor.Density, PhysicsFunctional_Euler, phys)); PetscCall(ModelFunctionalRegister(mod, "Momentum", &eu->monitor.Momentum, PhysicsFunctional_Euler, phys)); PetscCall(ModelFunctionalRegister(mod, "Pressure", &eu->monitor.Pressure, PhysicsFunctional_Euler, phys)); PetscFunctionReturn(PETSC_SUCCESS); } static PetscErrorCode ErrorIndicator_Simple(PetscInt dim, PetscReal volume, PetscInt numComps, const PetscScalar u[], const PetscScalar grad[], PetscReal *error, PetscCtx ctx) { PetscReal err = 0.; PetscInt i, j; PetscFunctionBeginUser; for (i = 0; i < numComps; i++) { for (j = 0; j < dim; j++) err += PetscSqr(PetscRealPart(grad[i * dim + j])); } *error = volume * err; PetscFunctionReturn(PETSC_SUCCESS); } PetscErrorCode CreatePartitionVec(DM dm, DM *dmCell, Vec *partition) { PetscSF sfPoint; PetscSection coordSection; Vec coordinates; PetscSection sectionCell; PetscScalar *part; PetscInt cStart, cEnd, c; PetscMPIInt rank; PetscFunctionBeginUser; PetscCall(DMGetCoordinateSection(dm, &coordSection)); PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); PetscCall(DMClone(dm, dmCell)); PetscCall(DMGetPointSF(dm, &sfPoint)); PetscCall(DMSetPointSF(*dmCell, sfPoint)); PetscCall(DMSetCoordinateSection(*dmCell, PETSC_DETERMINE, coordSection)); PetscCall(DMSetCoordinatesLocal(*dmCell, coordinates)); PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ionCell)); PetscCall(DMPlexGetHeightStratum(*dmCell, 0, &cStart, &cEnd)); PetscCall(PetscSectionSetChart(sectionCell, cStart, cEnd)); for (c = cStart; c < cEnd; ++c) PetscCall(PetscSectionSetDof(sectionCell, c, 1)); PetscCall(PetscSectionSetUp(sectionCell)); PetscCall(DMSetLocalSection(*dmCell, sectionCell)); PetscCall(PetscSectionDestroy(§ionCell)); PetscCall(DMCreateLocalVector(*dmCell, partition)); PetscCall(PetscObjectSetName((PetscObject)*partition, "partition")); PetscCall(VecGetArray(*partition, &part)); for (c = cStart; c < cEnd; ++c) { PetscScalar *p; PetscCall(DMPlexPointLocalRef(*dmCell, c, part, &p)); p[0] = rank; } PetscCall(VecRestoreArray(*partition, &part)); PetscFunctionReturn(PETSC_SUCCESS); } PetscErrorCode CreateMassMatrix(DM dm, Vec *massMatrix, User user) { DM plex, dmMass, dmFace, dmCell, dmCoord; PetscSection coordSection; Vec coordinates, facegeom, cellgeom; PetscSection sectionMass; PetscScalar *m; const PetscScalar *fgeom, *cgeom, *coords; PetscInt vStart, vEnd, v; PetscFunctionBeginUser; PetscCall(DMConvert(dm, DMPLEX, &plex)); PetscCall(DMGetCoordinateSection(dm, &coordSection)); PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); PetscCall(DMClone(dm, &dmMass)); PetscCall(DMSetCoordinateSection(dmMass, PETSC_DETERMINE, coordSection)); PetscCall(DMSetCoordinatesLocal(dmMass, coordinates)); PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ionMass)); PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); PetscCall(PetscSectionSetChart(sectionMass, vStart, vEnd)); for (v = vStart; v < vEnd; ++v) { PetscInt numFaces; PetscCall(DMPlexGetSupportSize(dmMass, v, &numFaces)); PetscCall(PetscSectionSetDof(sectionMass, v, numFaces * numFaces)); } PetscCall(PetscSectionSetUp(sectionMass)); PetscCall(DMSetLocalSection(dmMass, sectionMass)); PetscCall(PetscSectionDestroy(§ionMass)); PetscCall(DMGetLocalVector(dmMass, massMatrix)); PetscCall(VecGetArray(*massMatrix, &m)); PetscCall(DMPlexGetGeometryFVM(plex, &facegeom, &cellgeom, NULL)); PetscCall(VecGetDM(facegeom, &dmFace)); PetscCall(VecGetArrayRead(facegeom, &fgeom)); PetscCall(VecGetDM(cellgeom, &dmCell)); PetscCall(VecGetArrayRead(cellgeom, &cgeom)); PetscCall(DMGetCoordinateDM(dm, &dmCoord)); PetscCall(VecGetArrayRead(coordinates, &coords)); for (v = vStart; v < vEnd; ++v) { const PetscInt *faces; PetscFVFaceGeom *fgA, *fgB, *cg; PetscScalar *vertex; PetscInt numFaces, sides[2], f, g; PetscCall(DMPlexPointLocalRead(dmCoord, v, coords, &vertex)); PetscCall(DMPlexGetSupportSize(dmMass, v, &numFaces)); PetscCall(DMPlexGetSupport(dmMass, v, &faces)); for (f = 0; f < numFaces; ++f) { sides[0] = faces[f]; PetscCall(DMPlexPointLocalRead(dmFace, faces[f], fgeom, &fgA)); for (g = 0; g < numFaces; ++g) { const PetscInt *cells = NULL; PetscReal area = 0.0; PetscInt numCells; sides[1] = faces[g]; PetscCall(DMPlexPointLocalRead(dmFace, faces[g], fgeom, &fgB)); PetscCall(DMPlexGetJoin(dmMass, 2, sides, &numCells, &cells)); PetscCheck(numCells == 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "Invalid join for faces"); PetscCall(DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cg)); area += PetscAbsScalar((vertex[0] - cg->centroid[0]) * (fgA->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1]) * (fgA->centroid[0] - cg->centroid[0])); area += PetscAbsScalar((vertex[0] - cg->centroid[0]) * (fgB->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1]) * (fgB->centroid[0] - cg->centroid[0])); m[f * numFaces + g] = Dot2Real(fgA->normal, fgB->normal) * area * 0.5; PetscCall(DMPlexRestoreJoin(dmMass, 2, sides, &numCells, &cells)); } } } PetscCall(VecRestoreArrayRead(facegeom, &fgeom)); PetscCall(VecRestoreArrayRead(cellgeom, &cgeom)); PetscCall(VecRestoreArrayRead(coordinates, &coords)); PetscCall(VecRestoreArray(*massMatrix, &m)); PetscCall(DMDestroy(&dmMass)); PetscCall(DMDestroy(&plex)); PetscFunctionReturn(PETSC_SUCCESS); } /* Behavior will be different for multi-physics or when using non-default boundary conditions */ static PetscErrorCode ModelSolutionSetDefault(Model mod, SolutionFunction func, PetscCtx ctx) { PetscFunctionBeginUser; mod->solution = func; mod->solutionctx = ctx; PetscFunctionReturn(PETSC_SUCCESS); } static PetscErrorCode ModelFunctionalRegister(Model mod, const char *name, PetscInt *offset, FunctionalFunction func, PetscCtx ctx) { FunctionalLink link, *ptr; PetscInt lastoffset = -1; PetscFunctionBeginUser; for (ptr = &mod->functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset; PetscCall(PetscNew(&link)); PetscCall(PetscStrallocpy(name, &link->name)); link->offset = lastoffset + 1; link->func = func; link->ctx = ctx; link->next = NULL; *ptr = link; *offset = link->offset; PetscFunctionReturn(PETSC_SUCCESS); } static PetscErrorCode ModelFunctionalSetFromOptions(Model mod, PetscOptionItems PetscOptionsObject) { PetscInt i, j; FunctionalLink link; char *names[256]; PetscFunctionBeginUser; mod->numMonitored = PETSC_STATIC_ARRAY_LENGTH(names); PetscCall(PetscOptionsStringArray("-monitor", "list of functionals to monitor", "", names, &mod->numMonitored, NULL)); /* Create list of functionals that will be computed somehow */ PetscCall(PetscMalloc1(mod->numMonitored, &mod->functionalMonitored)); /* Create index of calls that we will have to make to compute these functionals (over-allocation in general). */ PetscCall(PetscMalloc1(mod->numMonitored, &mod->functionalCall)); mod->numCall = 0; for (i = 0; i < mod->numMonitored; i++) { for (link = mod->functionalRegistry; link; link = link->next) { PetscBool match; PetscCall(PetscStrcasecmp(names[i], link->name, &match)); if (match) break; } PetscCheck(link, mod->comm, PETSC_ERR_USER, "No known functional '%s'", names[i]); mod->functionalMonitored[i] = link; for (j = 0; j < i; j++) { if (mod->functionalCall[j]->func == link->func && mod->functionalCall[j]->ctx == link->ctx) goto next_name; } mod->functionalCall[mod->numCall++] = link; /* Just points to the first link using the result. There may be more results. */ next_name: PetscCall(PetscFree(names[i])); } /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */ mod->maxComputed = -1; for (link = mod->functionalRegistry; link; link = link->next) { for (i = 0; i < mod->numCall; i++) { FunctionalLink call = mod->functionalCall[i]; if (link->func == call->func && link->ctx == call->ctx) mod->maxComputed = PetscMax(mod->maxComputed, link->offset); } } PetscFunctionReturn(PETSC_SUCCESS); } static PetscErrorCode FunctionalLinkDestroy(FunctionalLink *link) { FunctionalLink l, next; PetscFunctionBeginUser; if (!link) PetscFunctionReturn(PETSC_SUCCESS); l = *link; *link = NULL; for (; l; l = next) { next = l->next; PetscCall(PetscFree(l->name)); PetscCall(PetscFree(l)); } PetscFunctionReturn(PETSC_SUCCESS); } /* put the solution callback into a functional callback */ static PetscErrorCode SolutionFunctional(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *modctx) { Model mod; PetscFunctionBeginUser; mod = (Model)modctx; PetscCall((*mod->solution)(mod, time, x, u, mod->solutionctx)); PetscFunctionReturn(PETSC_SUCCESS); } PetscErrorCode SetInitialCondition(DM dm, Vec X, User user) { PetscErrorCode (*func[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx); PetscCtx ctx[1]; Model mod = user->model; PetscFunctionBeginUser; func[0] = SolutionFunctional; ctx[0] = (void *)mod; PetscCall(DMProjectFunction(dm, 0.0, func, ctx, INSERT_ALL_VALUES, X)); PetscFunctionReturn(PETSC_SUCCESS); } static PetscErrorCode OutputVTK(DM dm, const char *filename, PetscViewer *viewer) { PetscFunctionBeginUser; PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)dm), viewer)); PetscCall(PetscViewerSetType(*viewer, PETSCVIEWERVTK)); PetscCall(PetscViewerFileSetName(*viewer, filename)); PetscFunctionReturn(PETSC_SUCCESS); } static PetscErrorCode MonitorVTK(TS ts, PetscInt stepnum, PetscReal time, Vec X, PetscCtx ctx) { User user = (User)ctx; DM dm, plex; PetscViewer viewer; char filename[PETSC_MAX_PATH_LEN], *ftable = NULL; PetscReal xnorm; PetscBool rollback; PetscFunctionBeginUser; PetscCall(TSGetStepRollBack(ts, &rollback)); if (rollback) PetscFunctionReturn(PETSC_SUCCESS); PetscCall(PetscObjectSetName((PetscObject)X, "u")); PetscCall(VecGetDM(X, &dm)); PetscCall(VecNorm(X, NORM_INFINITY, &xnorm)); if (stepnum >= 0) stepnum += user->monitorStepOffset; if (stepnum >= 0) { /* No summary for final time */ Model mod = user->model; Vec cellgeom; PetscInt c, cStart, cEnd, fcount, i; size_t ftableused, ftablealloc; const PetscScalar *cgeom, *x; DM dmCell; DMLabel vtkLabel; PetscReal *fmin, *fmax, *fintegral, *ftmp; PetscCall(DMConvert(dm, DMPLEX, &plex)); PetscCall(DMPlexGetGeometryFVM(plex, NULL, &cellgeom, NULL)); fcount = mod->maxComputed + 1; PetscCall(PetscMalloc4(fcount, &fmin, fcount, &fmax, fcount, &fintegral, fcount, &ftmp)); for (i = 0; i < fcount; i++) { fmin[i] = PETSC_MAX_REAL; fmax[i] = PETSC_MIN_REAL; fintegral[i] = 0; } PetscCall(VecGetDM(cellgeom, &dmCell)); PetscCall(DMPlexGetSimplexOrBoxCells(dmCell, 0, &cStart, &cEnd)); PetscCall(VecGetArrayRead(cellgeom, &cgeom)); PetscCall(VecGetArrayRead(X, &x)); PetscCall(DMGetLabel(dm, "vtk", &vtkLabel)); for (c = cStart; c < cEnd; ++c) { PetscFVCellGeom *cg; const PetscScalar *cx = NULL; PetscInt vtkVal = 0; /* not that these two routines as currently implemented work for any dm with a * localSection/globalSection */ PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); PetscCall(DMPlexPointGlobalRead(dm, c, x, &cx)); if (vtkLabel) PetscCall(DMLabelGetValue(vtkLabel, c, &vtkVal)); if (!vtkVal || !cx) continue; /* ghost, or not a global cell */ for (i = 0; i < mod->numCall; i++) { FunctionalLink flink = mod->functionalCall[i]; PetscCall((*flink->func)(mod, time, cg->centroid, cx, ftmp, flink->ctx)); } for (i = 0; i < fcount; i++) { fmin[i] = PetscMin(fmin[i], ftmp[i]); fmax[i] = PetscMax(fmax[i], ftmp[i]); fintegral[i] += cg->volume * ftmp[i]; } } PetscCall(VecRestoreArrayRead(cellgeom, &cgeom)); PetscCall(VecRestoreArrayRead(X, &x)); PetscCall(DMDestroy(&plex)); PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, fmin, (PetscMPIInt)fcount, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)ts))); PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, fmax, (PetscMPIInt)fcount, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts))); PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, fintegral, (PetscMPIInt)fcount, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)ts))); ftablealloc = fcount * 100; ftableused = 0; PetscCall(PetscMalloc1(ftablealloc, &ftable)); for (i = 0; i < mod->numMonitored; i++) { size_t countused; char buffer[256], *p; FunctionalLink flink = mod->functionalMonitored[i]; PetscInt id = flink->offset; if (i % 3) { PetscCall(PetscArraycpy(buffer, " ", 2)); p = buffer + 2; } else if (i) { char newline[] = "\n"; PetscCall(PetscMemcpy(buffer, newline, sizeof(newline) - 1)); p = buffer + sizeof(newline) - 1; } else { p = buffer; } PetscCall(PetscSNPrintfCount(p, sizeof buffer - (p - buffer), "%12s [%10.7g,%10.7g] int %10.7g", &countused, flink->name, (double)fmin[id], (double)fmax[id], (double)fintegral[id])); countused--; countused += p - buffer; if (countused > ftablealloc - ftableused - 1) { /* reallocate */ char *ftablenew; ftablealloc = 2 * ftablealloc + countused; PetscCall(PetscMalloc(ftablealloc, &ftablenew)); PetscCall(PetscArraycpy(ftablenew, ftable, ftableused)); PetscCall(PetscFree(ftable)); ftable = ftablenew; } PetscCall(PetscArraycpy(ftable + ftableused, buffer, countused)); ftableused += countused; ftable[ftableused] = 0; } PetscCall(PetscFree4(fmin, fmax, fintegral, ftmp)); PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "% 3" PetscInt_FMT " time %8.4g |x| %8.4g %s\n", stepnum, (double)time, (double)xnorm, ftable ? ftable : "")); PetscCall(PetscFree(ftable)); } if (user->vtkInterval < 1) PetscFunctionReturn(PETSC_SUCCESS); if ((stepnum == -1) ^ (stepnum % user->vtkInterval == 0)) { if (stepnum == -1) { /* Final time is not multiple of normal time interval, write it anyway */ PetscCall(TSGetStepNumber(ts, &stepnum)); } PetscCall(PetscSNPrintf(filename, sizeof filename, "%s-%03" PetscInt_FMT ".vtu", user->outputBasename, stepnum)); PetscCall(OutputVTK(dm, filename, &viewer)); PetscCall(VecView(X, viewer)); PetscCall(PetscViewerDestroy(&viewer)); } PetscFunctionReturn(PETSC_SUCCESS); } static PetscErrorCode initializeTS(DM dm, User user, TS *ts) { #ifdef PETSC_HAVE_LIBCEED PetscBool useCeed; #endif PetscFunctionBeginUser; PetscCall(TSCreate(PetscObjectComm((PetscObject)dm), ts)); PetscCall(TSSetType(*ts, TSSSP)); PetscCall(TSSetDM(*ts, dm)); if (user->vtkmon) PetscCall(TSMonitorSet(*ts, MonitorVTK, user, NULL)); PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, user)); #ifdef PETSC_HAVE_LIBCEED PetscCall(DMPlexGetUseCeed(dm, &useCeed)); if (useCeed) PetscCall(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVMCEED, user)); else #endif PetscCall(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, user)); PetscCall(TSSetMaxTime(*ts, 2.0)); PetscCall(TSSetExactFinalTime(*ts, TS_EXACTFINALTIME_STEPOVER)); PetscFunctionReturn(PETSC_SUCCESS); } typedef struct { PetscFV fvm; VecTagger refineTag; VecTagger coarsenTag; DM adaptedDM; User user; PetscReal cfl; PetscLimiter limiter; PetscLimiter noneLimiter; } TransferCtx; static PetscErrorCode adaptToleranceFVMSetUp(TS ts, PetscInt nstep, PetscReal time, Vec sol, PetscBool *resize, PetscCtx ctx) { TransferCtx *tctx = (TransferCtx *)ctx; PetscFV fvm = tctx->fvm; VecTagger refineTag = tctx->refineTag; VecTagger coarsenTag = tctx->coarsenTag; User user = tctx->user; DM dm, gradDM, plex, cellDM, adaptedDM = NULL; Vec cellGeom, faceGeom; PetscBool computeGradient; Vec grad, locGrad, locX, errVec; PetscInt cStart, cEnd, c, dim, nRefine, nCoarsen; PetscReal minMaxInd[2] = {PETSC_MAX_REAL, PETSC_MIN_REAL}, minMaxIndGlobal[2]; PetscScalar *errArray; const PetscScalar *pointVals; const PetscScalar *pointGrads; const PetscScalar *pointGeom; DMLabel adaptLabel = NULL; IS refineIS, coarsenIS; PetscFunctionBeginUser; *resize = PETSC_FALSE; PetscCall(VecGetDM(sol, &dm)); PetscCall(DMGetDimension(dm, &dim)); PetscCall(PetscFVSetLimiter(fvm, tctx->noneLimiter)); PetscCall(PetscFVGetComputeGradients(fvm, &computeGradient)); PetscCall(PetscFVSetComputeGradients(fvm, PETSC_TRUE)); PetscCall(DMConvert(dm, DMPLEX, &plex)); PetscCall(DMPlexGetDataFVM(plex, fvm, &cellGeom, &faceGeom, &gradDM)); PetscCall(DMCreateLocalVector(plex, &locX)); PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, faceGeom, cellGeom, NULL)); PetscCall(DMGlobalToLocalBegin(plex, sol, INSERT_VALUES, locX)); PetscCall(DMGlobalToLocalEnd(plex, sol, INSERT_VALUES, locX)); PetscCall(DMCreateGlobalVector(gradDM, &grad)); PetscCall(DMPlexReconstructGradientsFVM(plex, locX, grad)); PetscCall(DMCreateLocalVector(gradDM, &locGrad)); PetscCall(DMGlobalToLocalBegin(gradDM, grad, INSERT_VALUES, locGrad)); PetscCall(DMGlobalToLocalEnd(gradDM, grad, INSERT_VALUES, locGrad)); PetscCall(VecDestroy(&grad)); PetscCall(DMPlexGetSimplexOrBoxCells(plex, 0, &cStart, &cEnd)); PetscCall(VecGetArrayRead(locGrad, &pointGrads)); PetscCall(VecGetArrayRead(cellGeom, &pointGeom)); PetscCall(VecGetArrayRead(locX, &pointVals)); PetscCall(VecGetDM(cellGeom, &cellDM)); PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel)); PetscCall(VecCreateFromOptions(PetscObjectComm((PetscObject)plex), NULL, 1, cEnd - cStart, PETSC_DETERMINE, &errVec)); PetscCall(VecSetUp(errVec)); PetscCall(VecGetArray(errVec, &errArray)); for (c = cStart; c < cEnd; c++) { PetscReal errInd = 0.; PetscScalar *pointGrad; PetscScalar *pointVal; PetscFVCellGeom *cg; PetscCall(DMPlexPointLocalRead(gradDM, c, pointGrads, &pointGrad)); PetscCall(DMPlexPointLocalRead(cellDM, c, pointGeom, &cg)); PetscCall(DMPlexPointLocalRead(plex, c, pointVals, &pointVal)); PetscCall((*user->model->errorIndicator)(dim, cg->volume, user->model->physics->dof, pointVal, pointGrad, &errInd, user->model->errorCtx)); errArray[c - cStart] = errInd; minMaxInd[0] = PetscMin(minMaxInd[0], errInd); minMaxInd[1] = PetscMax(minMaxInd[1], errInd); } PetscCall(VecRestoreArray(errVec, &errArray)); PetscCall(VecRestoreArrayRead(locX, &pointVals)); PetscCall(VecRestoreArrayRead(cellGeom, &pointGeom)); PetscCall(VecRestoreArrayRead(locGrad, &pointGrads)); PetscCall(VecDestroy(&locGrad)); PetscCall(VecDestroy(&locX)); PetscCall(DMDestroy(&plex)); PetscCall(VecTaggerComputeIS(refineTag, errVec, &refineIS, NULL)); PetscCall(VecTaggerComputeIS(coarsenTag, errVec, &coarsenIS, NULL)); PetscCall(ISGetSize(refineIS, &nRefine)); PetscCall(ISGetSize(coarsenIS, &nCoarsen)); if (nRefine) PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_REFINE, refineIS)); if (nCoarsen) PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_COARSEN, coarsenIS)); PetscCall(ISDestroy(&coarsenIS)); PetscCall(ISDestroy(&refineIS)); PetscCall(VecDestroy(&errVec)); PetscCall(PetscFVSetComputeGradients(fvm, computeGradient)); PetscCall(PetscFVSetLimiter(fvm, tctx->limiter)); minMaxInd[1] = -minMaxInd[1]; PetscCallMPI(MPIU_Allreduce(minMaxInd, minMaxIndGlobal, 2, MPIU_REAL, MPI_MIN, PetscObjectComm((PetscObject)dm))); PetscCall(PetscInfo(ts, "error indicator range (%E, %E)\n", (double)minMaxIndGlobal[0], (double)(-minMaxIndGlobal[1]))); if (nRefine || nCoarsen) { /* at least one cell is over the refinement threshold */ PetscCall(DMAdaptLabel(dm, adaptLabel, &adaptedDM)); } PetscCall(DMLabelDestroy(&adaptLabel)); if (adaptedDM) { tctx->adaptedDM = adaptedDM; PetscCall(PetscInfo(ts, "Step %" PetscInt_FMT " adapted mesh, marking %" PetscInt_FMT " cells for refinement, and %" PetscInt_FMT " cells for coarsening\n", nstep, nRefine, nCoarsen)); *resize = PETSC_TRUE; } else { PetscCall(PetscInfo(ts, "Step %" PetscInt_FMT " no adaptation\n", nstep)); } PetscFunctionReturn(PETSC_SUCCESS); } static PetscErrorCode Transfer(TS ts, PetscInt nv, Vec vecsin[], Vec vecsout[], PetscCtx ctx) { TransferCtx *tctx = (TransferCtx *)ctx; DM dm; PetscReal time; PetscFunctionBeginUser; PetscCall(TSGetDM(ts, &dm)); PetscCall(TSGetTime(ts, &time)); PetscCheck(tctx->adaptedDM, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "Missing adaptedDM"); for (PetscInt i = 0; i < nv; i++) { PetscCall(DMCreateGlobalVector(tctx->adaptedDM, &vecsout[i])); PetscCall(DMForestTransferVec(dm, vecsin[i], tctx->adaptedDM, vecsout[i], PETSC_TRUE, time)); } PetscCall(DMForestSetAdaptivityForest(tctx->adaptedDM, NULL)); /* clear internal references to the previous dm */ Model mod = tctx->user->model; Physics phys = mod->physics; PetscReal minRadius; PetscCall(DMPlexGetGeometryFVM(tctx->adaptedDM, NULL, NULL, &minRadius)); PetscCallMPI(MPIU_Allreduce(&phys->maxspeed, &mod->maxspeed, 1, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts))); PetscCheck(mod->maxspeed > 0, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "Physics did not set maxspeed"); PetscReal dt = tctx->cfl * minRadius / mod->maxspeed; PetscCall(TSSetTimeStep(ts, dt)); PetscCall(TSSetDM(ts, tctx->adaptedDM)); PetscCall(DMDestroy(&tctx->adaptedDM)); PetscFunctionReturn(PETSC_SUCCESS); } int main(int argc, char **argv) { MPI_Comm comm; PetscDS prob; PetscFV fvm; PetscLimiter limiter = NULL, noneLimiter = NULL; User user; Model mod; Physics phys; DM dm, plex; PetscReal ftime, cfl, dt, minRadius; PetscInt dim, nsteps; TS ts; TSConvergedReason reason; Vec X; PetscViewer viewer; PetscBool vtkCellGeom, useAMR; PetscInt adaptInterval; char physname[256] = "advect"; VecTagger refineTag = NULL, coarsenTag = NULL; TransferCtx tctx; PetscFunctionBeginUser; PetscCall(PetscInitialize(&argc, &argv, NULL, help)); comm = PETSC_COMM_WORLD; PetscCall(PetscNew(&user)); PetscCall(PetscNew(&user->model)); PetscCall(PetscNew(&user->model->physics)); mod = user->model; phys = mod->physics; mod->comm = comm; useAMR = PETSC_FALSE; adaptInterval = 1; /* Register physical models to be available on the command line */ PetscCall(PetscFunctionListAdd(&PhysicsList, "advect", PhysicsCreate_Advect)); PetscCall(PetscFunctionListAdd(&PhysicsList, "sw", PhysicsCreate_SW)); PetscCall(PetscFunctionListAdd(&PhysicsList, "euler", PhysicsCreate_Euler)); PetscOptionsBegin(comm, NULL, "Unstructured Finite Volume Mesh Options", ""); { cfl = 0.9 * 4; /* default SSPRKS2 with s=5 stages is stable for CFL number s-1 */ PetscCall(PetscOptionsReal("-ufv_cfl", "CFL number per step", "", cfl, &cfl, NULL)); user->vtkInterval = 1; PetscCall(PetscOptionsInt("-ufv_vtk_interval", "VTK output interval (0 to disable)", "", user->vtkInterval, &user->vtkInterval, NULL)); user->vtkmon = PETSC_TRUE; PetscCall(PetscOptionsBool("-ufv_vtk_monitor", "Use VTKMonitor routine", "", user->vtkmon, &user->vtkmon, NULL)); vtkCellGeom = PETSC_FALSE; PetscCall(PetscStrncpy(user->outputBasename, "ex11", sizeof(user->outputBasename))); PetscCall(PetscOptionsString("-ufv_vtk_basename", "VTK output basename", "", user->outputBasename, user->outputBasename, sizeof(user->outputBasename), NULL)); PetscCall(PetscOptionsBool("-ufv_vtk_cellgeom", "Write cell geometry (for debugging)", "", vtkCellGeom, &vtkCellGeom, NULL)); PetscCall(PetscOptionsBool("-ufv_use_amr", "use local adaptive mesh refinement", "", useAMR, &useAMR, NULL)); PetscCall(PetscOptionsInt("-ufv_adapt_interval", "time steps between AMR", "", adaptInterval, &adaptInterval, NULL)); } PetscOptionsEnd(); if (useAMR) { VecTaggerBox refineBox, coarsenBox; refineBox.min = refineBox.max = PETSC_MAX_REAL; coarsenBox.min = coarsenBox.max = PETSC_MIN_REAL; PetscCall(VecTaggerCreate(comm, &refineTag)); PetscCall(PetscObjectSetOptionsPrefix((PetscObject)refineTag, "refine_")); PetscCall(VecTaggerSetType(refineTag, VECTAGGERABSOLUTE)); PetscCall(VecTaggerAbsoluteSetBox(refineTag, &refineBox)); PetscCall(VecTaggerSetFromOptions(refineTag)); PetscCall(VecTaggerSetUp(refineTag)); PetscCall(PetscObjectViewFromOptions((PetscObject)refineTag, NULL, "-tag_view")); PetscCall(VecTaggerCreate(comm, &coarsenTag)); PetscCall(PetscObjectSetOptionsPrefix((PetscObject)coarsenTag, "coarsen_")); PetscCall(VecTaggerSetType(coarsenTag, VECTAGGERABSOLUTE)); PetscCall(VecTaggerAbsoluteSetBox(coarsenTag, &coarsenBox)); PetscCall(VecTaggerSetFromOptions(coarsenTag)); PetscCall(VecTaggerSetUp(coarsenTag)); PetscCall(PetscObjectViewFromOptions((PetscObject)coarsenTag, NULL, "-tag_view")); } PetscOptionsBegin(comm, NULL, "Unstructured Finite Volume Physics Options", ""); { PetscErrorCode (*physcreate)(Model, Physics, PetscOptionItems); PetscCall(PetscOptionsFList("-physics", "Physics module to solve", "", PhysicsList, physname, physname, sizeof physname, NULL)); PetscCall(PetscFunctionListFind(PhysicsList, physname, &physcreate)); PetscCall(PetscMemzero(phys, sizeof(struct _n_Physics))); PetscCall((*physcreate)(mod, phys, PetscOptionsObject)); /* Count number of fields and dofs */ for (phys->nfields = 0, phys->dof = 0; phys->field_desc[phys->nfields].name; phys->nfields++) phys->dof += phys->field_desc[phys->nfields].dof; PetscCheck(phys->dof > 0, comm, PETSC_ERR_ARG_WRONGSTATE, "Physics '%s' did not set dof", physname); PetscCall(ModelFunctionalSetFromOptions(mod, PetscOptionsObject)); } PetscOptionsEnd(); /* Create mesh */ { PetscInt i; PetscCall(DMCreate(comm, &dm)); PetscCall(DMSetType(dm, DMPLEX)); PetscCall(DMSetFromOptions(dm)); for (i = 0; i < DIM; i++) { mod->bounds[2 * i] = 0.; mod->bounds[2 * i + 1] = 1.; } dim = DIM; { /* a null name means just do a hex box */ PetscInt cells[3] = {1, 1, 1}, n = 3; PetscBool flg2, skew = PETSC_FALSE; PetscInt nret2 = 2 * DIM; PetscOptionsBegin(comm, NULL, "Rectangular mesh options", ""); PetscCall(PetscOptionsRealArray("-grid_bounds", "bounds of the mesh in each direction (i.e., x_min,x_max,y_min,y_max", "", mod->bounds, &nret2, &flg2)); PetscCall(PetscOptionsBool("-grid_skew_60", "Skew grid for 60 degree shock mesh", "", skew, &skew, NULL)); PetscCall(PetscOptionsIntArray("-dm_plex_box_faces", "Number of faces along each dimension", "", cells, &n, NULL)); PetscOptionsEnd(); /* TODO Rewrite this with Mark, and remove grid_bounds at that time */ if (flg2) { PetscInt dimEmbed, i; PetscInt nCoords; PetscScalar *coords; Vec coordinates; PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); PetscCall(DMGetCoordinateDim(dm, &dimEmbed)); PetscCall(VecGetLocalSize(coordinates, &nCoords)); PetscCheck(!(nCoords % dimEmbed), PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Coordinate vector the wrong size"); PetscCall(VecGetArray(coordinates, &coords)); for (i = 0; i < nCoords; i += dimEmbed) { PetscInt j; PetscScalar *coord = &coords[i]; for (j = 0; j < dimEmbed; j++) { coord[j] = mod->bounds[2 * j] + coord[j] * (mod->bounds[2 * j + 1] - mod->bounds[2 * j]); if (dim == 2 && cells[1] == 1 && j == 0 && skew) { if (cells[0] == 2 && i == 8) { coord[j] = .57735026918963; /* hack to get 60 deg skewed mesh */ } else if (cells[0] == 3) { if (i == 2 || i == 10) coord[j] = mod->bounds[1] / 4.; else if (i == 4) coord[j] = mod->bounds[1] / 2.; else if (i == 12) coord[j] = 1.57735026918963 * mod->bounds[1] / 2.; } } } } PetscCall(VecRestoreArray(coordinates, &coords)); PetscCall(DMSetCoordinatesLocal(dm, coordinates)); } } } PetscCall(DMViewFromOptions(dm, NULL, "-orig_dm_view")); PetscCall(DMGetDimension(dm, &dim)); /* set up BCs, functions, tags */ PetscCall(DMCreateLabel(dm, "Face Sets")); mod->errorIndicator = ErrorIndicator_Simple; { DM gdm; PetscCall(DMPlexConstructGhostCells(dm, NULL, NULL, &gdm)); PetscCall(DMDestroy(&dm)); dm = gdm; PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); } PetscCall(PetscFVCreate(comm, &fvm)); PetscCall(PetscFVSetFromOptions(fvm)); PetscCall(PetscFVSetNumComponents(fvm, phys->dof)); PetscCall(PetscFVSetSpatialDimension(fvm, dim)); PetscCall(PetscObjectSetName((PetscObject)fvm, "")); { PetscInt f, dof; for (f = 0, dof = 0; f < phys->nfields; f++) { PetscInt newDof = phys->field_desc[f].dof; if (newDof == 1) { PetscCall(PetscFVSetComponentName(fvm, dof, phys->field_desc[f].name)); } else { PetscInt j; for (j = 0; j < newDof; j++) { char compName[256] = "Unknown"; PetscCall(PetscSNPrintf(compName, sizeof(compName), "%s_%" PetscInt_FMT, phys->field_desc[f].name, j)); PetscCall(PetscFVSetComponentName(fvm, dof + j, compName)); } } dof += newDof; } } /* FV is now structured with one field having all physics as components */ PetscCall(DMAddField(dm, NULL, (PetscObject)fvm)); PetscCall(DMCreateDS(dm)); PetscCall(DMGetDS(dm, &prob)); PetscCall(PetscDSSetRiemannSolver(prob, 0, user->model->physics->riemann)); PetscCall(PetscDSSetContext(prob, 0, user->model->physics)); PetscCall((*mod->setupbc)(dm, prob, phys)); PetscCall(PetscDSSetFromOptions(prob)); { char convType[256]; PetscBool flg; PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX"); PetscCall(PetscOptionsFList("-dm_type", "Convert DMPlex to another format", "ex12", DMList, DMPLEX, convType, 256, &flg)); PetscOptionsEnd(); if (flg) { DM dmConv; PetscCall(DMConvert(dm, convType, &dmConv)); if (dmConv) { PetscCall(DMViewFromOptions(dmConv, NULL, "-dm_conv_view")); PetscCall(DMDestroy(&dm)); dm = dmConv; PetscCall(DMSetFromOptions(dm)); } } } #ifdef PETSC_HAVE_LIBCEED { PetscBool useCeed; PetscCall(DMPlexGetUseCeed(dm, &useCeed)); if (useCeed) PetscCall((*user->model->setupCEED)(dm, user->model->physics)); } #endif PetscCall(initializeTS(dm, user, &ts)); PetscCall(DMCreateGlobalVector(dm, &X)); PetscCall(PetscObjectSetName((PetscObject)X, "solution")); PetscCall(SetInitialCondition(dm, X, user)); if (useAMR) { PetscInt adaptIter; /* use no limiting when reconstructing gradients for adaptivity */ PetscCall(PetscFVGetLimiter(fvm, &limiter)); PetscCall(PetscObjectReference((PetscObject)limiter)); PetscCall(PetscLimiterCreate(PetscObjectComm((PetscObject)fvm), &noneLimiter)); PetscCall(PetscLimiterSetType(noneLimiter, PETSCLIMITERNONE)); /* Refinement context */ tctx.fvm = fvm; tctx.refineTag = refineTag; tctx.coarsenTag = coarsenTag; tctx.adaptedDM = NULL; tctx.user = user; tctx.noneLimiter = noneLimiter; tctx.limiter = limiter; tctx.cfl = cfl; /* Do some initial refinement steps */ for (adaptIter = 0;; ++adaptIter) { PetscLogDouble bytes; PetscBool resize; PetscCall(PetscMemoryGetCurrentUsage(&bytes)); PetscCall(PetscInfo(ts, "refinement loop %" PetscInt_FMT ": memory used %g\n", adaptIter, bytes)); PetscCall(DMViewFromOptions(dm, NULL, "-initial_dm_view")); PetscCall(VecViewFromOptions(X, NULL, "-initial_vec_view")); PetscCall(adaptToleranceFVMSetUp(ts, -1, 0.0, X, &resize, &tctx)); if (!resize) break; PetscCall(DMDestroy(&dm)); PetscCall(VecDestroy(&X)); dm = tctx.adaptedDM; tctx.adaptedDM = NULL; PetscCall(TSSetDM(ts, dm)); PetscCall(DMCreateGlobalVector(dm, &X)); PetscCall(PetscObjectSetName((PetscObject)X, "solution")); PetscCall(SetInitialCondition(dm, X, user)); } } PetscCall(DMConvert(dm, DMPLEX, &plex)); if (vtkCellGeom) { DM dmCell; Vec cellgeom, partition; PetscCall(DMPlexGetGeometryFVM(plex, NULL, &cellgeom, NULL)); PetscCall(OutputVTK(dm, "ex11-cellgeom.vtk", &viewer)); PetscCall(VecView(cellgeom, viewer)); PetscCall(PetscViewerDestroy(&viewer)); PetscCall(CreatePartitionVec(dm, &dmCell, &partition)); PetscCall(OutputVTK(dmCell, "ex11-partition.vtk", &viewer)); PetscCall(VecView(partition, viewer)); PetscCall(PetscViewerDestroy(&viewer)); PetscCall(VecDestroy(&partition)); PetscCall(DMDestroy(&dmCell)); } /* collect max maxspeed from all processes -- todo */ PetscCall(DMPlexGetGeometryFVM(plex, NULL, NULL, &minRadius)); PetscCall(DMDestroy(&plex)); PetscCallMPI(MPIU_Allreduce(&phys->maxspeed, &mod->maxspeed, 1, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)ts))); PetscCheck(mod->maxspeed > 0, comm, PETSC_ERR_ARG_WRONGSTATE, "Physics '%s' did not set maxspeed", physname); dt = cfl * minRadius / mod->maxspeed; PetscCall(TSSetTimeStep(ts, dt)); PetscCall(TSSetFromOptions(ts)); /* When using adaptive mesh refinement specify callbacks to refine the solution and interpolate data from old to new mesh When mesh adaption is requested, the step will be rolled back */ if (useAMR) PetscCall(TSSetResize(ts, PETSC_TRUE, adaptToleranceFVMSetUp, Transfer, &tctx)); PetscCall(TSSetSolution(ts, X)); PetscCall(VecDestroy(&X)); PetscCall(TSSolve(ts, NULL)); PetscCall(TSGetSolveTime(ts, &ftime)); PetscCall(TSGetStepNumber(ts, &nsteps)); PetscCall(TSGetConvergedReason(ts, &reason)); PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %" PetscInt_FMT " steps\n", TSConvergedReasons[reason], (double)ftime, nsteps)); PetscCall(TSDestroy(&ts)); PetscCall(VecTaggerDestroy(&refineTag)); PetscCall(VecTaggerDestroy(&coarsenTag)); PetscCall(PetscFunctionListDestroy(&PhysicsList)); PetscCall(PetscFunctionListDestroy(&PhysicsRiemannList_SW)); PetscCall(PetscFunctionListDestroy(&PhysicsRiemannList_Euler)); PetscCall(FunctionalLinkDestroy(&user->model->functionalRegistry)); PetscCall(PetscFree(user->model->functionalMonitored)); PetscCall(PetscFree(user->model->functionalCall)); PetscCall(PetscFree(user->model->physics->data)); PetscCall(PetscFree(user->model->physics)); PetscCall(PetscFree(user->model)); PetscCall(PetscFree(user)); PetscCall(PetscLimiterDestroy(&limiter)); PetscCall(PetscLimiterDestroy(&noneLimiter)); PetscCall(PetscFVDestroy(&fvm)); PetscCall(DMDestroy(&dm)); PetscCall(PetscFinalize()); return 0; } /* Subroutine to set up the initial conditions for the */ /* Shock Interface interaction or linear wave (Ravi Samtaney,Mark Adams). */ /* ----------------------------------------------------------------------- */ int projecteqstate(PetscReal wc[], const PetscReal ueq[], PetscReal lv[][3]) { int j, k; /* Wc=matmul(lv,Ueq) 3 vars */ for (k = 0; k < 3; ++k) { wc[k] = 0.; for (j = 0; j < 3; ++j) wc[k] += lv[k][j] * ueq[j]; } return 0; } /* ----------------------------------------------------------------------- */ int projecttoprim(PetscReal v[], const PetscReal wc[], PetscReal rv[][3]) { int k, j; /* V=matmul(rv,WC) 3 vars */ for (k = 0; k < 3; ++k) { v[k] = 0.; for (j = 0; j < 3; ++j) v[k] += rv[k][j] * wc[j]; } return 0; } /* ---------------------------------------------------------------------- */ int eigenvectors(PetscReal rv[][3], PetscReal lv[][3], const PetscReal ueq[], PetscReal gamma) { int j, k; PetscReal rho, csnd, p0; /* PetscScalar u; */ for (k = 0; k < 3; ++k) for (j = 0; j < 3; ++j) { lv[k][j] = 0.; rv[k][j] = 0.; } rho = ueq[0]; /* u = ueq[1]; */ p0 = ueq[2]; csnd = PetscSqrtReal(gamma * p0 / rho); lv[0][1] = rho * .5; lv[0][2] = -.5 / csnd; lv[1][0] = csnd; lv[1][2] = -1. / csnd; lv[2][1] = rho * .5; lv[2][2] = .5 / csnd; rv[0][0] = -1. / csnd; rv[1][0] = 1. / rho; rv[2][0] = -csnd; rv[0][1] = 1. / csnd; rv[0][2] = 1. / csnd; rv[1][2] = 1. / rho; rv[2][2] = csnd; return 0; } int initLinearWave(EulerNode *ux, const PetscReal gamma, const PetscReal coord[], const PetscReal Lx) { PetscReal p0, u0, wcp[3], wc[3]; PetscReal lv[3][3]; PetscReal vp[3]; PetscReal rv[3][3]; PetscReal eps, ueq[3], rho0, twopi; /* Function Body */ twopi = 2. * PETSC_PI; eps = 1e-4; /* perturbation */ rho0 = 1e3; /* density of water */ p0 = 101325.; /* init pressure of 1 atm (?) */ u0 = 0.; ueq[0] = rho0; ueq[1] = u0; ueq[2] = p0; /* Project initial state to characteristic variables */ eigenvectors(rv, lv, ueq, gamma); projecteqstate(wc, ueq, lv); wcp[0] = wc[0]; wcp[1] = wc[1]; wcp[2] = wc[2] + eps * PetscCosReal(coord[0] * 2. * twopi / Lx); projecttoprim(vp, wcp, rv); ux->r = vp[0]; /* density */ ux->ru[0] = vp[0] * vp[1]; /* x momentum */ ux->ru[1] = 0.; /* E = rho * e + rho * v^2/2 = p/(gam-1) + rho*v^2/2 */ ux->E = vp[2] / (gamma - 1.) + 0.5 * vp[0] * vp[1] * vp[1]; return 0; } /*TEST testset: args: -dm_plex_adj_cone -dm_plex_adj_closure 0 test: suffix: adv_2d_tri_0 requires: triangle TODO: how did this ever get in main when there is no support for this args: -ufv_vtk_interval 0 -simplex -dm_refine 3 -dm_plex_faces 1,1 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3 test: suffix: adv_2d_tri_1 requires: triangle TODO: how did this ever get in main when there is no support for this args: -ufv_vtk_interval 0 -simplex -dm_refine 5 -dm_plex_faces 1,1 -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1 test: suffix: tut_1 requires: exodusii nsize: 1 args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo - test: suffix: tut_2 requires: exodusii nsize: 1 args: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw test: suffix: tut_3 requires: exodusii nsize: 4 args: -dm_distribute_overlap 1 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -monitor Error -advect_sol_type bump -petscfv_type leastsquares -petsclimiter_type sin test: suffix: tut_4 requires: exodusii !single nsize: 4 args: -dm_distribute_overlap 1 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -physics sw -monitor Height,Energy -petscfv_type leastsquares -petsclimiter_type minmod testset: args: -dm_plex_adj_cone -dm_plex_adj_closure 0 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1 # 2D Advection 0-10 test: suffix: 0 requires: exodusii args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo test: suffix: 1 requires: exodusii args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo test: suffix: 2 requires: exodusii nsize: 2 args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo test: suffix: 3 requires: exodusii nsize: 2 args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo test: suffix: 4 requires: exodusii nsize: 4 args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo -petscpartitioner_type simple test: suffix: 5 requires: exodusii args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw -ts_adapt_reject_safety 1 test: suffix: 7 requires: exodusii args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1 test: suffix: 8 requires: exodusii nsize: 2 args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1 test: suffix: 9 requires: exodusii nsize: 8 args: -dm_distribute_overlap 1 -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1 test: suffix: 10 requires: exodusii args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo # 2D Shallow water testset: args: -physics sw -ufv_vtk_interval 0 -dm_plex_adj_cone -dm_plex_adj_closure 0 test: suffix: sw_0 requires: exodusii args: -bc_wall 100,101 -ufv_cfl 5 -petscfv_type leastsquares -petsclimiter_type sin \ -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo \ -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10 \ -monitor height,energy test: suffix: sw_ceed requires: exodusii libceed args: -sw_riemann rusanov_ceed -bc_wall 100,101 -ufv_cfl 5 -petsclimiter_type sin \ -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -dm_plex_use_ceed \ -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10 \ -monitor height,energy test: suffix: sw_ceed_small requires: exodusii libceed args: -sw_riemann rusanov_ceed -bc_wall 1,3 -ufv_cfl 5 -petsclimiter_type sin -dm_plex_use_ceed \ -dm_plex_shape annulus -dm_plex_simplex 0 -dm_plex_box_lower 0,1 -dm_plex_box_upper 6.28,3 -dm_plex_box_faces 8,2 \ -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10 \ -monitor height,energy test: suffix: sw_1 nsize: 2 args: -bc_wall 1,3 -ufv_cfl 5 -petsclimiter_type sin \ -dm_plex_shape annulus -dm_plex_simplex 0 -dm_plex_box_faces 24,12 -dm_plex_box_lower 0,1 -dm_plex_box_upper 1,3 -dm_distribute_overlap 1 \ -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10 \ -monitor height,energy test: suffix: sw_hll args: -sw_riemann hll -bc_wall 1,2,3,4 -ufv_cfl 3 -petscfv_type leastsquares -petsclimiter_type sin \ -grid_bounds 0,5,0,5 -dm_plex_simplex 0 -dm_plex_box_faces 25,25 \ -ts_max_steps 5 -ts_ssp_type rks2 -ts_ssp_nstages 10 \ -monitor height,energy # 2D Euler testset: args: -physics euler -eu_type linear_wave -eu_gamma 1.4 -dm_plex_adj_cone -dm_plex_adj_closure 0 \ -ufv_vtk_interval 0 -ufv_vtk_basename ${wPETSC_DIR}/ex11 -monitor density,energy test: suffix: euler_0 requires: exodusii !complex !single args: -eu_riemann godunov -bc_wall 100,101 -ufv_cfl 5 -petsclimiter_type sin \ -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo \ -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10 test: suffix: euler_ceed requires: exodusii libceed args: -eu_riemann godunov_ceed -bc_wall 100,101 -ufv_cfl 5 -petsclimiter_type sin \ -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -dm_plex_use_ceed \ -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10 testset: args: -dm_plex_adj_cone -dm_plex_adj_closure 0 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1 # 2D Advection: p4est test: suffix: p4est_advec_2d requires: p4est args: -ufv_vtk_interval 0 -dm_type p4est -dm_forest_minimum_refinement 1 -dm_forest_initial_refinement 2 -dm_p4est_refine_pattern hash -dm_forest_maximum_refinement 5 # Advection in a box test: suffix: adv_2d_quad_0 args: -ufv_vtk_interval 0 -dm_refine 3 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3 test: suffix: adv_2d_quad_1 args: -ufv_vtk_interval 0 -dm_refine 3 -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1 timeoutfactor: 3 test: suffix: adv_2d_quad_p4est_0 requires: p4est args: -ufv_vtk_interval 0 -dm_refine 5 -dm_type p4est -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3 test: suffix: adv_2d_quad_p4est_1 requires: p4est args: -ufv_vtk_interval 0 -dm_refine 5 -dm_type p4est -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1 timeoutfactor: 3 test: # broken for quad precision suffix: adv_2d_quad_p4est_adapt_0 requires: p4est !__float128 args: -ufv_vtk_interval 0 -dm_refine 3 -dm_type p4est -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1 -ufv_use_amr -refine_vec_tagger_box 0.005,inf -coarsen_vec_tagger_box 0,1.e-5 -petscfv_type leastsquares -ts_max_time 0.01 timeoutfactor: 3 test: suffix: adv_0 requires: exodusii args: -ufv_vtk_interval 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -bc_inflow 100,101,200 -bc_outflow 201 # Run with -dm_forest_maximum_refinement 6 -ts_max_time 0.5 instead to get the full movie test: suffix: shock_0 requires: p4est !single !complex args: -dm_plex_box_faces 2,1 -grid_bounds -1,1.,0.,1 -grid_skew_60 \ -dm_type p4est -dm_forest_partition_overlap 1 -dm_forest_maximum_refinement 2 -dm_forest_minimum_refinement 2 -dm_forest_initial_refinement 2 \ -ufv_use_amr -refine_vec_tagger_box 0.5,inf -coarsen_vec_tagger_box 0,1.e-2 -refine_tag_view -coarsen_tag_view \ -bc_wall 1,2,3,4 -physics euler -eu_type iv_shock -ufv_cfl 10 -eu_alpha 60. -eu_gamma 1.4 -eu_amach 2.02 -eu_rho2 3. \ -petscfv_type leastsquares -petsclimiter_type minmod -petscfv_compute_gradients 0 \ -ts_max_steps 3 -ts_ssp_type rks2 -ts_ssp_nstages 10 \ -ufv_vtk_basename ${wPETSC_DIR}/ex11 -ufv_vtk_interval 0 -monitor density,energy # Test GLVis visualization of PetscFV fields test: suffix: glvis_adv_2d_tri args: -ufv_vtk_interval 0 -ufv_vtk_monitor 0 \ -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_gmsh_periodic 0 \ -ts_monitor_solution glvis: -ts_max_steps 0 test: suffix: glvis_adv_2d_quad args: -ufv_vtk_interval 0 -ufv_vtk_monitor 0 -bc_inflow 1,2,4 -bc_outflow 3 \ -dm_refine 5 -dm_plex_separate_marker \ -ts_monitor_solution glvis: -ts_max_steps 0 # Test CGNS file writing for PetscFV fields test: suffix: cgns_adv_2d_tri requires: cgns args: -ufv_vtk_interval 0 -ufv_vtk_monitor 0 \ -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_gmsh_periodic 0 \ -ts_monitor_solution cgns:sol.cgns -ts_max_steps 0 test: suffix: cgns_adv_2d_quad requires: cgns args: -ufv_vtk_interval 0 -ufv_vtk_monitor 0 -bc_inflow 1,2,4 -bc_outflow 3 \ -dm_refine 5 -dm_plex_separate_marker \ -ts_monitor_solution cgns:sol.cgns -ts_max_steps 0 # Test CGNS file writing, cgns_batch_size, ts_run_steps, ts_monitor_skip_initial test: suffix: cgns_adv_2d_tri_monitor requires: cgns args: -ufv_vtk_interval 0 -ufv_vtk_monitor 0 \ -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_gmsh_periodic 0 \ -ts_monitor_solution cgns:sol-%d.cgns -ts_run_steps 4 -ts_monitor_solution_interval 2 \ -viewer_cgns_batch_size 1 -ts_monitor_solution_skip_initial # Test VTK file writing for PetscFV fields with -ts_monitor_solution_vtk test: suffix: vtk_adv_2d_tri args: -ufv_vtk_interval 0 -ufv_vtk_monitor 0 \ -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_gmsh_periodic 0 \ -ts_monitor_solution_vtk 'bar-%03d.vtu' -ts_monitor_solution_vtk_interval 2 -ts_max_steps 5 TEST*/