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

CompFlow 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/CompFlow/Physics.h
  • Problem - problem configuration, see PDE/CompFlow/Problem.h

Constructors, destructors, conversion operators

CompFlow() 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
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>>&, tk::Fields& unk, tk::real t, const std::size_t nielem) const
void lhs(const tk::Fields& geoElem, tk::Fields& l) const
void updateInterfaceCells(tk::Fields&, std::size_t, std::vector<std::size_t>&) const
void updatePrimitives(const tk::Fields&, const tk::Fields&, const tk::Fields&, tk::Fields&, std::size_t) const
void cleanTraceMaterial(tk::real, const tk::Fields&, tk::Fields&, tk::Fields&, std::size_t) const
void reconstruct(tk::real t, const tk::Fields& geoFace, const tk::Fields& geoElem, const inciter::FaceData& fd, const std::map<std::size_t, std::vector<std::size_t>>&, const std::vector<std::size_t>& inpoel, const tk::UnsMesh::Coords& coord, tk::Fields& U, tk::Fields& P) const
void limit(] tk::real t, ] 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>>&, const std::vector<std::vector<tk::real>>& mtInv, tk::Fields& U, tk::Fields&, std::vector<std::size_t>& shockmarker) const
void CPL(const tk::Fields&, const tk::Fields&, const std::vector<std::size_t>&, const tk::UnsMesh::Coords&, tk::Fields&, std::size_t) const
auto cellAvgDeformGrad(const tk::Fields&, std::size_t) const -> std::array<std::vector<tk::real>, 9>
void rhs(tk::real t, 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>>& boxelems, 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>& coord, const std::vector<std::size_t>& inpoel, const inciter::FaceData& fd, const tk::Fields& geoFace, const tk::Fields& geoElem, const std::vector<std::size_t>& ndofel, const tk::Fields& U, const tk::Fields&, const std::size_t) const -> tk::real
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&) 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::CompFlow<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::CompFlow<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::CompFlow<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::CompFlow<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>
void inciter::dg::CompFlow<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>>&, 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
unk in/out Array of unknowns
in Physical time
nielem in Number of internal elements

Initalize the compressible flow equations, prepare for time integration

template<class Physics, class Problem>
void inciter::dg::CompFlow<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::CompFlow<Physics, Problem>::updateInterfaceCells(tk::Fields&, std::size_t, std::vector<std::size_t>&) const

Update the interface cells to first order dofs

This function resets the high-order terms in interface cells, and is currently not used in compflow.

template<class Physics, class Problem>
void inciter::dg::CompFlow<Physics, Problem>::updatePrimitives(const tk::Fields&, const tk::Fields&, const tk::Fields&, tk::Fields&, std::size_t) const

Update the primitives for this PDE system

This function computes and stores the dofs for primitive quantities, which is currently unused for compflow. But if a limiter requires primitive variables for example, this would be the place to add the computation of the primitive variables.

template<class Physics, class Problem>
void inciter::dg::CompFlow<Physics, Problem>::cleanTraceMaterial(tk::real, const tk::Fields&, tk::Fields&, tk::Fields&, std::size_t) const

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. This is unused for compflow.

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

Parameters
in Physical time
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/out Solution vector at recent time step
in/out Primitive vector at recent time step

Reconstruct second-order solution from first-order using least-squares

template<class Physics, class Problem>
void inciter::dg::CompFlow<Physics, Problem>::limit(] tk::real t, ] 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>>&, const std::vector<std::vector<tk::real>>& mtInv, tk::Fields& U, tk::Fields&, std::vector<std::size_t>& shockmarker) const

Parameters
in Physical time
geoFace in Face geometry array
geoElem in Element geometry array
fd in Face connectivity and boundary conditions object
esup in Elements surrounding points
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
mtInv in Inverse of Taylor mass matrix
in/out Solution vector at recent time step
shockmarker in/out Vector of shock-marker values

Limit second-order solution

template<class Physics, class Problem>
void inciter::dg::CompFlow<Physics, Problem>::CPL(const tk::Fields&, const tk::Fields&, const std::vector<std::size_t>&, const tk::UnsMesh::Coords&, tk::Fields&, std::size_t) const

Update the conservative variable solution for this PDE system

This function computes the updated dofs for conservative quantities based on the limited solution and is currently not used in compflow.

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

Return cell-average deformation gradient tensor (no-op for compflow)

This function is a no-op in compflow.

template<class Physics, class Problem>
void inciter::dg::CompFlow<Physics, Problem>::rhs(tk::real t, 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>>& boxelems, 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
geoFace in Face geometry array
geoElem in Element geometry array
fd in Face connectivity and boundary conditions object
inpoel in Element-node connectivity
boxelems in Mesh node ids within user-defined IC boxes
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::CompFlow<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::CompFlow<Physics, Problem>::dt(const std::array<std::vector<tk::real>, 3>& coord, const std::vector<std::size_t>& inpoel, const inciter::FaceData& fd, const tk::Fields& geoFace, const tk::Fields& geoElem, const std::vector<std::size_t>& ndofel, const tk::Fields& U, const tk::Fields&, const std::size_t) const

Parameters
coord in Mesh node coordinates
inpoel in Mesh element connectivity
fd in Face connectivity and boundary conditions object
geoFace in Face geometry array
geoElem in Element geometry array
ndofel in Vector of local number of degrees of freedom
in Solution vector at recent time step
Returns Minimum time step size

Compute the minimum time step size

template<class Physics, class Problem>
std::array<std::array<tk::real, 4>, 3> inciter::dg::CompFlow<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::CompFlow<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::CompFlow<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::CompFlow<Physics, Problem>::histNames() const

Returns Vector of strings labeling 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::CompFlow<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&) const

Parameters
in History point data
inpoel in Element-node connectivity
coord in Array of nodal coordinates
in Array of unknowns

Return time history field output evaluated at time history points

template<class Physics, class Problem>
std::vector<std::string> inciter::dg::CompFlow<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::CompFlow<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::CompFlow<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::CompFlow<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