template<class Physics, class Problem>
inciter::dg::MultiMat class

MultiMat used polymorphically with tk::DGPDE.

The template arguments specify policies and are used to configure the behavior of the class. The policies are:

  • Physics - physics configuration, see PDE/MultiMat/Physics.h
  • Problem - problem configuration, see PDE/MultiMat/Problem.h

Constructors, destructors, conversion operators

MultiMat() explicit
Constructor.

Public functions

auto nprim() const -> std::size_t
auto nmat() const -> std::size_t
void numEquationDofs(std::vector<std::size_t>& numEqDof) const
void IcBoxElems(const tk::Fields& geoElem, std::size_t nielem, std::vector<std::unordered_set<std::size_t>>& inbox) const
auto nstiffeq() const -> std::size_t
auto nnonstiffeq() const -> std::size_t
void setStiffEqIdx(std::vector<std::size_t>& stiffEqIdx) const
void setNonStiffEqIdx(std::vector<std::size_t>& nonStiffEqIdx) const
void initialize(const tk::Fields& L, const std::vector<std::size_t>& inpoel, const tk::UnsMesh::Coords& coord, const std::vector<std::unordered_set<std::size_t>>& inbox, const std::unordered_map<std::size_t, std::set<std::size_t>>& elemblkid, tk::Fields& unk, tk::real t, const std::size_t nielem) const
void computeDensityConstr(std::size_t nelem, tk::Fields& unk, std::vector<tk::real>& densityConstr) const
void lhs(const tk::Fields& geoElem, tk::Fields& l) const
void updateInterfaceCells(tk::Fields& unk, std::size_t nielem, std::vector<std::size_t>& ndofel, std::vector<std::size_t>& interface) const
void updatePrimitives(const tk::Fields& unk, const tk::Fields& L, const tk::Fields& geoElem, tk::Fields& prim, std::size_t nielem, const std::vector<std::size_t>& ndofel) const
void cleanTraceMaterial(tk::real t, const tk::Fields& geoElem, tk::Fields& unk, tk::Fields& prim, std::size_t nielem) const
void reconstruct(tk::real, const tk::Fields&, const tk::Fields& geoElem, const inciter::FaceData& fd, const std::map<std::size_t, std::vector<std::size_t>>& esup, const std::vector<std::size_t>& inpoel, const tk::UnsMesh::Coords& coord, tk::Fields& U, tk::Fields& P, const bool pref, const std::vector<std::size_t>& ndofel) const
void limit(] tk::real t, const bool pref, const tk::Fields& geoFace, const tk::Fields& geoElem, const inciter::FaceData& fd, const std::map<std::size_t, std::vector<std::size_t>>& esup, const std::vector<std::size_t>& inpoel, const tk::UnsMesh::Coords& coord, const std::vector<std::size_t>& ndofel, const std::vector<std::size_t>& gid, const std::unordered_map<std::size_t, std::size_t>& bid, const std::vector<std::vector<tk::real>>& uNodalExtrm, const std::vector<std::vector<tk::real>>& pNodalExtrm, const std::vector<std::vector<tk::real>>& mtInv, tk::Fields& U, tk::Fields& P, std::vector<std::size_t>& shockmarker) const
void CPL(const tk::Fields& prim, const tk::Fields& geoElem, const std::vector<std::size_t>& inpoel, const tk::UnsMesh::Coords& coord, tk::Fields& unk, std::size_t nielem) const
auto cellAvgDeformGrad(const tk::Fields& unk, std::size_t nielem) const -> std::array<std::vector<tk::real>, 9>
void resetAdapSol(const inciter::FaceData& fd, tk::Fields& unk, tk::Fields& prim, const std::vector<std::size_t>& ndofel) const
void rhs(tk::real t, const bool pref, const tk::Fields& geoFace, const tk::Fields& geoElem, const inciter::FaceData& fd, const std::vector<std::size_t>& inpoel, const std::vector<std::unordered_set<std::size_t>>&, const tk::UnsMesh::Coords& coord, const tk::Fields& U, const tk::Fields& P, const std::vector<std::size_t>& ndofel, const tk::real dt, tk::Fields& R) const
void eval_ndof(std::size_t nunk, ] const tk::UnsMesh::Coords& coord, ] const std::vector<std::size_t>& inpoel, const inciter::FaceData& fd, const tk::Fields& unk, ] const tk::Fields& prim, inciter::ctr::PrefIndicatorType indicator, std::size_t ndof, std::size_t ndofmax, tk::real tolref, std::vector<std::size_t>& ndofel) const
auto dt(const std::array<std::vector<tk::real>, 3>&, const std::vector<std::size_t>&, const inciter::FaceData& fd, const tk::Fields& geoFace, const tk::Fields& geoElem, const std::vector<std::size_t>&, const tk::Fields& U, const tk::Fields& P, const std::size_t nielem) const -> tk::real
void stiff_rhs(std::size_t e, const tk::Fields& geoElem, const std::vector<std::size_t>& inpoel, const tk::UnsMesh::Coords& coord, const tk::Fields& U, const tk::Fields& P, const std::vector<std::size_t>& ndofel, tk::Fields& R) const
auto velocity(const tk::Fields& U, const std::array<std::vector<tk::real>, 3>&, const std::array<std::size_t, 4>& N) const -> std::array<std::array<tk::real, 4>, 3>
auto OutVarFn() const -> std::map<std::string, tk::GetVarFn>
auto analyticFieldNames() const -> std::vector<std::string>
auto histNames() const -> std::vector<std::string>
auto surfOutput(const std::map<int, std::vector<std::size_t>>&, tk::Fields&) const -> std::vector<std::vector<tk::real>>
Return surface field output going to file.
auto histOutput(const std::vector<HistData>& h, const std::vector<std::size_t>& inpoel, const tk::UnsMesh::Coords& coord, const tk::Fields& U, const tk::Fields& P) const -> std::vector<std::vector<tk::real>>
auto names() const -> std::vector<std::string>
auto analyticSolution(tk::real xi, tk::real yi, tk::real zi, tk::real t) const -> std::vector<tk::real>
auto solution(tk::real xi, tk::real yi, tk::real zi, tk::real t) const -> std::vector<tk::real>
auto sp_totalenergy(std::size_t e, const tk::Fields& unk) const -> tk::real

Function documentation

template<class Physics, class Problem>
std::size_t inciter::dg::MultiMat<Physics, Problem>::nprim() const

Returns The number of primitive quantities required to be stored for this PDE system

Find the number of primitive quantities required for this PDE system

template<class Physics, class Problem>
std::size_t inciter::dg::MultiMat<Physics, Problem>::nmat() const

Returns The number of materials set up for this PDE system

Find the number of materials set up for this PDE system

template<class Physics, class Problem>
void inciter::dg::MultiMat<Physics, Problem>::numEquationDofs(std::vector<std::size_t>& numEqDof) const

Parameters
numEqDof in/out Array storing number of Dofs for each PDE equation

Assign number of DOFs per equation in the PDE system

template<class Physics, class Problem>
void inciter::dg::MultiMat<Physics, Problem>::IcBoxElems(const tk::Fields& geoElem, std::size_t nielem, std::vector<std::unordered_set<std::size_t>>& inbox) const

Parameters
geoElem in Element geometry array
nielem in Number of internal elements
inbox in/out List of nodes at which box user ICs are set for each IC box

Determine elements that lie inside the user-defined IC box

template<class Physics, class Problem>
std::size_t inciter::dg::MultiMat<Physics, Problem>::nstiffeq() const

Returns number of stiff equations

Find how many 'stiff equations', which are the inverse deformation equations because of plasticity

template<class Physics, class Problem>
std::size_t inciter::dg::MultiMat<Physics, Problem>::nnonstiffeq() const

Returns number of stiff equations

Find how many 'non-stiff equations', which are the inverse deformation equations because of plasticity

template<class Physics, class Problem>
void inciter::dg::MultiMat<Physics, Problem>::setStiffEqIdx(std::vector<std::size_t>& stiffEqIdx) const

Parameters
stiffEqIdx out list with pointers to stiff equations

Locate the stiff equations.

template<class Physics, class Problem>
void inciter::dg::MultiMat<Physics, Problem>::setNonStiffEqIdx(std::vector<std::size_t>& nonStiffEqIdx) const

Parameters
nonStiffEqIdx out list with pointers to nonstiff equations

Locate the nonstiff equations.

template<class Physics, class Problem>
void inciter::dg::MultiMat<Physics, Problem>::initialize(const tk::Fields& L, const std::vector<std::size_t>& inpoel, const tk::UnsMesh::Coords& coord, const std::vector<std::unordered_set<std::size_t>>& inbox, const std::unordered_map<std::size_t, std::set<std::size_t>>& elemblkid, tk::Fields& unk, tk::real t, const std::size_t nielem) const

Parameters
in Block diagonal mass matrix
inpoel in Element-node connectivity
coord in Array of nodal coordinates
inbox in List of elements at which box user ICs are set for each IC box
elemblkid in Element ids associated with mesh block ids where user ICs are set
unk in/out Array of unknowns
in Physical time
nielem in Number of internal elements

Initialize the compressible flow equations, prepare for time integration

template<class Physics, class Problem>
void inciter::dg::MultiMat<Physics, Problem>::computeDensityConstr(std::size_t nelem, tk::Fields& unk, std::vector<tk::real>& densityConstr) const

Parameters
nelem in Number of elements
unk in Array of unknowns
densityConstr out Density Constraint: rho/(rho0*det(g))

Compute density constraint for a given material

template<class Physics, class Problem>
void inciter::dg::MultiMat<Physics, Problem>::lhs(const tk::Fields& geoElem, tk::Fields& l) const

Parameters
geoElem in Element geometry array
in/out Block diagonal mass matrix

Compute the left hand side block-diagonal mass matrix

template<class Physics, class Problem>
void inciter::dg::MultiMat<Physics, Problem>::updateInterfaceCells(tk::Fields& unk, std::size_t nielem, std::vector<std::size_t>& ndofel, std::vector<std::size_t>& interface) const

Parameters
unk in Array of unknowns
nielem in Number of internal elements
ndofel in/out Array of dofs
interface in/out Vector of interface marker

Update the interface cells to first order dofs This function resets the high-order terms in interface cells.

template<class Physics, class Problem>
void inciter::dg::MultiMat<Physics, Problem>::updatePrimitives(const tk::Fields& unk, const tk::Fields& L, const tk::Fields& geoElem, tk::Fields& prim, std::size_t nielem, const std::vector<std::size_t>& ndofel) const

Parameters
unk in Array of unknowns
in The left hand side block-diagonal mass matrix
geoElem in Element geometry array
prim in/out Array of primitives
nielem in Number of internal elements
ndofel in Array of dofs

Update the primitives for this PDE system This function computes and stores the dofs for primitive quantities, which are required for obtaining reconstructed states used in the Riemann solver. See /PDE/Riemann/AUSM.hpp, where the normal velocity for advection is calculated from independently reconstructed velocities.

template<class Physics, class Problem>
void inciter::dg::MultiMat<Physics, Problem>::cleanTraceMaterial(tk::real t, const tk::Fields& geoElem, tk::Fields& unk, tk::Fields& prim, std::size_t nielem) const

Parameters
in Physical time
geoElem in Element geometry array
unk in/out Array of unknowns
prim in/out Array of primitives
nielem in Number of internal elements

Clean up the state of trace materials for this PDE system This function cleans up the state of materials present in trace quantities in each cell. Specifically, the state of materials with very low volume-fractions in a cell is replaced by the state of the material which is present in the largest quantity in that cell. This becomes necessary when shocks pass through cells which contain a very small amount of material. The state of that tiny material might become unphysical and cause solution to diverge; thus requiring such a "reset".

template<class Physics, class Problem>
void inciter::dg::MultiMat<Physics, Problem>::reconstruct(tk::real, const tk::Fields&, const tk::Fields& geoElem, const inciter::FaceData& fd, const std::map<std::size_t, std::vector<std::size_t>>& esup, const std::vector<std::size_t>& inpoel, const tk::UnsMesh::Coords& coord, tk::Fields& U, tk::Fields& P, const bool pref, const std::vector<std::size_t>& ndofel) const

Parameters
geoElem in Element geometry array
fd in Face connectivity and boundary conditions object
esup in Elements-surrounding-nodes connectivity
inpoel in Element-node connectivity
coord in Array of nodal coordinates
in/out Solution vector at recent time step
in/out Vector of primitives at recent time step
pref in Indicator for p-adaptive algorithm
ndofel in Vector of local number of degrees of freedome

Reconstruct second-order solution from first-order

template<class Physics, class Problem>
void inciter::dg::MultiMat<Physics, Problem>::limit(] tk::real t, const bool pref, const tk::Fields& geoFace, const tk::Fields& geoElem, const inciter::FaceData& fd, const std::map<std::size_t, std::vector<std::size_t>>& esup, const std::vector<std::size_t>& inpoel, const tk::UnsMesh::Coords& coord, const std::vector<std::size_t>& ndofel, const std::vector<std::size_t>& gid, const std::unordered_map<std::size_t, std::size_t>& bid, const std::vector<std::vector<tk::real>>& uNodalExtrm, const std::vector<std::vector<tk::real>>& pNodalExtrm, const std::vector<std::vector<tk::real>>& mtInv, tk::Fields& U, tk::Fields& P, std::vector<std::size_t>& shockmarker) const

Parameters
in Physical time
pref in Indicator for p-adaptive algorithm
geoFace in Face geometry array
geoElem in Element geometry array
fd in Face connectivity and boundary conditions object
esup in Elements-surrounding-nodes connectivity
inpoel in Element-node connectivity
coord in Array of nodal coordinates
ndofel in Vector of local number of degrees of freedome
gid in Local->global node id map
bid in Local chare-boundary node ids (value) associated to global node ids (key)
uNodalExtrm in Chare-boundary nodal extrema for conservative variables
pNodalExtrm in Chare-boundary nodal extrema for primitive variables
mtInv in Inverse of Taylor mass matrix
in/out Solution vector at recent time step
in/out Vector of primitives at recent time step
shockmarker in/out Vector of shock-marker values

Limit second-order solution, and primitive quantities separately

template<class Physics, class Problem>
void inciter::dg::MultiMat<Physics, Problem>::CPL(const tk::Fields& prim, const tk::Fields& geoElem, const std::vector<std::size_t>& inpoel, const tk::UnsMesh::Coords& coord, tk::Fields& unk, std::size_t nielem) const

Parameters
prim in Array of primitive variables
geoElem in Element geometry array
inpoel in Element-node connectivity
coord in Array of nodal coordinates
unk in/out Array of conservative variables
nielem in Number of internal elements

Apply CPL to the conservative variable solution for this PDE system This function applies CPL to obtain consistent dofs for conservative quantities based on the limited primitive quantities. See Pandare et al. (2023). On the Design of Stable, Consistent, and Conservative High-Order Methods for Multi-Material Hydrodynamics. J Comp Phys, 112313.

template<class Physics, class Problem>
std::array<std::vector<tk::real>, 9> inciter::dg::MultiMat<Physics, Problem>::cellAvgDeformGrad(const tk::Fields& unk, std::size_t nielem) const

Parameters
unk in Solution vector at recent time step
nielem in Number of internal elements

Return cell-average deformation gradient tensor This function returns the bulk cell-average inverse deformation gradient tensor

template<class Physics, class Problem>
void inciter::dg::MultiMat<Physics, Problem>::resetAdapSol(const inciter::FaceData& fd, tk::Fields& unk, tk::Fields& prim, const std::vector<std::size_t>& ndofel) const

Parameters
fd in Face connectivity and boundary conditions object
unk in/out Solution vector at recent time step
prim in/out Primitive vector at recent time step
ndofel in Vector of local number of degrees of freedome

Reset the high order solution for p-adaptive scheme This function reset the high order coefficient for p-adaptive solution polynomials. Unlike compflow class, the high order of fv solution will not be reset since p0p1 is the base scheme for multi-material p-adaptive DG method.

template<class Physics, class Problem>
void inciter::dg::MultiMat<Physics, Problem>::rhs(tk::real t, const bool pref, const tk::Fields& geoFace, const tk::Fields& geoElem, const inciter::FaceData& fd, const std::vector<std::size_t>& inpoel, const std::vector<std::unordered_set<std::size_t>>&, const tk::UnsMesh::Coords& coord, const tk::Fields& U, const tk::Fields& P, const std::vector<std::size_t>& ndofel, const tk::real dt, tk::Fields& R) const

Parameters
in Physical time
pref in Indicator for p-adaptive algorithm
geoFace in Face geometry array
geoElem in Element geometry array
fd in Face connectivity and boundary conditions object
inpoel in Element-node connectivity
coord in Array of nodal coordinates
in Solution vector at recent time step
in Primitive vector at recent time step
ndofel in Vector of local number of degrees of freedom
dt in Delta time
in/out Right-hand side vector computed

Compute right hand side

template<class Physics, class Problem>
void inciter::dg::MultiMat<Physics, Problem>::eval_ndof(std::size_t nunk, ] const tk::UnsMesh::Coords& coord, ] const std::vector<std::size_t>& inpoel, const inciter::FaceData& fd, const tk::Fields& unk, ] const tk::Fields& prim, inciter::ctr::PrefIndicatorType indicator, std::size_t ndof, std::size_t ndofmax, tk::real tolref, std::vector<std::size_t>& ndofel) const

Parameters
nunk in Number of unknowns
coord in Array of nodal coordinates
inpoel in Element-node connectivity
fd in Face connectivity and boundary conditions object
unk in Array of unknowns
prim in Array of primitive quantities
indicator in p-refinement indicator type
ndof in Number of degrees of freedom in the solution
ndofmax in Max number of degrees of freedom for p-refinement
tolref in Tolerance for p-refinement
ndofel in/out Vector of local number of degrees of freedome

Evaluate the adaptive indicator and mark the ndof for each element

template<class Physics, class Problem>
tk::real inciter::dg::MultiMat<Physics, Problem>::dt(const std::array<std::vector<tk::real>, 3>&, const std::vector<std::size_t>&, const inciter::FaceData& fd, const tk::Fields& geoFace, const tk::Fields& geoElem, const std::vector<std::size_t>&, const tk::Fields& U, const tk::Fields& P, const std::size_t nielem) const

Parameters
fd in Face connectivity and boundary conditions object
geoFace in Face geometry array
geoElem in Element geometry array
in Solution vector at recent time step
in Vector of primitive quantities at recent time step
nielem in Number of internal elements
Returns Minimum time step size

Compute the minimum time step size The allowable dt is calculated by looking at the maximum wave-speed in elements surrounding each face, times the area of that face. Once the maximum of this quantity over the mesh is determined, the volume of each cell is divided by this quantity. A minimum of this ratio is found over the entire mesh, which gives the allowable dt.

template<class Physics, class Problem>
void inciter::dg::MultiMat<Physics, Problem>::stiff_rhs(std::size_t e, const tk::Fields& geoElem, const std::vector<std::size_t>& inpoel, const tk::UnsMesh::Coords& coord, const tk::Fields& U, const tk::Fields& P, const std::vector<std::size_t>& ndofel, tk::Fields& R) const

Parameters
in Element number
geoElem in Element geometry array
inpoel in Element-node connectivity
coord in Array of nodal coordinates
in Solution vector at recent time step
in Primitive vector at recent time step
ndofel in Vector of local number of degrees of freedom
in/out Right-hand side vector computed

Compute stiff terms for a single element

template<class Physics, class Problem>
std::array<std::array<tk::real, 4>, 3> inciter::dg::MultiMat<Physics, Problem>::velocity(const tk::Fields& U, const std::array<std::vector<tk::real>, 3>&, const std::array<std::size_t, 4>& N) const

Parameters
in Solution vector at recent time step
in Element node indices
Returns Array of the four values of the velocity field

Extract the velocity field at cell nodes. Currently unused.

template<class Physics, class Problem>
std::map<std::string, tk::GetVarFn> inciter::dg::MultiMat<Physics, Problem>::OutVarFn() const

Returns Map that associates user-specified strings to functions that compute relevant quantities to be output to file

Return a map that associates user-specified strings to functions

template<class Physics, class Problem>
std::vector<std::string> inciter::dg::MultiMat<Physics, Problem>::analyticFieldNames() const

Returns Vector of strings labelling analytic fields output in file

Return analytic field names to be output to file

template<class Physics, class Problem>
std::vector<std::string> inciter::dg::MultiMat<Physics, Problem>::histNames() const

Returns Vector of strings labelling time history fields output in file

Return time history field names to be output to file

template<class Physics, class Problem>
std::vector<std::vector<tk::real>> inciter::dg::MultiMat<Physics, Problem>::histOutput(const std::vector<HistData>& h, const std::vector<std::size_t>& inpoel, const tk::UnsMesh::Coords& coord, const tk::Fields& U, const tk::Fields& P) const

Parameters
in History point data
inpoel in Element-node connectivity
coord in Array of nodal coordinates
in Array of unknowns
in Array of primitive quantities
Returns Vector of time history output of bulk flow quantities (density, velocity, total energy, and pressure) evaluated at time history points

Return time history field output evaluated at time history points

template<class Physics, class Problem>
std::vector<std::string> inciter::dg::MultiMat<Physics, Problem>::names() const

Returns Vector of strings labelling integral variables output

Return names of integral variables to be output to diagnostics file

template<class Physics, class Problem>
std::vector<tk::real> inciter::dg::MultiMat<Physics, Problem>::analyticSolution(tk::real xi, tk::real yi, tk::real zi, tk::real t) const

Parameters
xi in X-coordinate at which to evaluate the analytic solution
yi in Y-coordinate at which to evaluate the analytic solution
zi in Z-coordinate at which to evaluate the analytic solution
in Physical time at which to evaluate the analytic solution
Returns Vector of analytic solution at given location and time

Return analytic solution (if defined by Problem) at xi, yi, zi, t

template<class Physics, class Problem>
std::vector<tk::real> inciter::dg::MultiMat<Physics, Problem>::solution(tk::real xi, tk::real yi, tk::real zi, tk::real t) const

Parameters
xi in X-coordinate at which to evaluate the analytic solution
yi in Y-coordinate at which to evaluate the analytic solution
zi in Z-coordinate at which to evaluate the analytic solution
in Physical time at which to evaluate the analytic solution
Returns Vector of analytic solution at given location and time

Return analytic solution for conserved variables

template<class Physics, class Problem>
tk::real inciter::dg::MultiMat<Physics, Problem>::sp_totalenergy(std::size_t e, const tk::Fields& unk) const

Parameters
in Element id for which total energy is required
unk in Vector of conserved quantities
Returns Cell-averaged specific total energy for given element

Return cell-averaged specific total energy for an element