SpECTRE  v2023.09.07
dg Namespace Reference

Functionality related to discontinuous Galerkin schemes. More...

Classes

struct  SimpleBoundaryData
 Distinguishes between field data, which can be projected to a mortar, and extra data, which will not be projected. More...
 
class  SimpleMortarData
 Storage of boundary data on two sides of a mortar. More...
 

Typedefs

template<size_t VolumeDim>
using MortarId = std::pair<::Direction< VolumeDim >, ElementId< VolumeDim > >
 
template<size_t MortarDim>
using MortarSize = std::array< Spectral::MortarSize, MortarDim >
 
template<size_t VolumeDim, typename ValueType >
using MortarMap = std::unordered_map< MortarId< VolumeDim >, ValueType, boost::hash< MortarId< VolumeDim > > >
 

Enumerations

enum class  Formulation { StrongInertial , WeakInertial }
 The DG formulation to use. More...
 

Functions

std::ostreamoperator<< (std::ostream &os, Formulation t)
 
template<typename InboxTag , size_t Dim, typename TemporalIdType , typename... InboxTags>
bool has_received_from_all_mortars (const TemporalIdType &temporal_id, const Element< Dim > &element, const tuples::TaggedTuple< InboxTags... > &inboxes)
 Determines if data on all mortars has been received for the InboxTag at time temporal_id. More...
 
template<size_t Dim>
void metric_identity_det_jac_times_inv_jac (gsl::not_null< InverseJacobian< DataVector, Dim, Frame::ElementLogical, Frame::Inertial > * > det_jac_times_inverse_jacobian, const Mesh< Dim > &mesh, const tnsr::I< DataVector, Dim, Frame::Inertial > &inertial_coords, const Jacobian< DataVector, Dim, Frame::ElementLogical, Frame::Inertial > &jacobian)
 Compute the Jacobian determinant times the inverse Jacobian so that the result is divergence-free. More...
 
template<size_t Dim>
void metric_identity_jacobian_quantities (gsl::not_null< InverseJacobian< DataVector, Dim, Frame::ElementLogical, Frame::Inertial > * > det_jac_times_inverse_jacobian, gsl::not_null< InverseJacobian< DataVector, Dim, Frame::ElementLogical, Frame::Inertial > * > inverse_jacobian, gsl::not_null< Jacobian< DataVector, Dim, Frame::ElementLogical, Frame::Inertial > * > jacobian, gsl::not_null< Scalar< DataVector > * > det_jacobian, const Mesh< Dim > &mesh, const tnsr::I< DataVector, Dim, Frame::Inertial > &inertial_coords)
 Compute the Jacobian, inverse Jacobian, and determinant of the Jacobian so that they satisfy the metric identities. More...
 
template<size_t Dim>
Mesh< Dim > mortar_mesh (const Mesh< Dim > &face_mesh1, const Mesh< Dim > &face_mesh2)
 Find a mesh for a mortar capable of representing data from either of two faces. More...
 
template<size_t Dim>
MortarSize< Dim - 1 > mortar_size (const ElementId< Dim > &self, const ElementId< Dim > &neighbor, size_t dimension, const OrientationMap< Dim > &orientation)
 Determine the size of the mortar (i.e., the part of the face it covers) for communicating with a neighbor. This is the size relative to the size of self, and will not generally agree with that determined by neighbor.
 
template<size_t Dim>
void apply_mass_matrix (const gsl::not_null< DataVector * > data, const Mesh< Dim > &mesh)
 Apply the DG mass matrix to the data, in the diagonal mass-matrix approximation ("mass-lumping") More...
 
template<size_t Dim, typename TagsList >
void apply_mass_matrix (const gsl::not_null< Variables< TagsList > * > data, const Mesh< Dim > &mesh)
 Apply the DG mass matrix to the data, in the diagonal mass-matrix approximation ("mass-lumping") More...
 
template<size_t Dim>
void apply_inverse_mass_matrix (const gsl::not_null< DataVector * > data, const Mesh< Dim > &mesh)
 Apply the inverse DG mass matrix to the data. More...
 
template<size_t Dim, typename TagsList >
void apply_inverse_mass_matrix (const gsl::not_null< Variables< TagsList > * > data, const Mesh< Dim > &mesh)
 Apply the inverse DG mass matrix to the data. More...
 
template<typename... BoundaryCorrectionTags>
void lift_flux (const gsl::not_null< Variables< tmpl::list< BoundaryCorrectionTags... > > * > boundary_correction_terms, const size_t extent_perpendicular_to_boundary, const Scalar< DataVector > &magnitude_of_face_normal)
 Lifts the flux contribution from an interface to the volume. More...
 
template<typename... FluxTags>
auto lift_flux (Variables< tmpl::list< FluxTags... > > flux, const size_t extent_perpendicular_to_boundary, const Scalar< DataVector > &magnitude_of_face_normal) -> Variables< tmpl::list< db::remove_tag_prefix< FluxTags >... > >
 Lifts the flux contribution from an interface to the volume. More...
 
template<typename Tags , size_t Dim>
void project_to_mortar (const gsl::not_null< Variables< Tags > * > result, const Variables< Tags > &vars, const Mesh< Dim > &face_mesh, const Mesh< Dim > &mortar_mesh, const MortarSize< Dim > &mortar_size)
 Project variables from a face to a mortar.
 
template<typename Tags , size_t Dim>
Variables< Tags > project_to_mortar (const Variables< Tags > &vars, const Mesh< Dim > &face_mesh, const Mesh< Dim > &mortar_mesh, const MortarSize< Dim > &mortar_size)
 Project variables from a face to a mortar.
 
template<typename Tags , size_t Dim>
void project_from_mortar (const gsl::not_null< Variables< Tags > * > result, const Variables< Tags > &vars, const Mesh< Dim > &face_mesh, const Mesh< Dim > &mortar_mesh, const MortarSize< Dim > &mortar_size)
 Project variables from a mortar to a face.
 
template<typename Tags , size_t Dim>
Variables< Tags > project_from_mortar (const Variables< Tags > &vars, const Mesh< Dim > &face_mesh, const Mesh< Dim > &mortar_mesh, const MortarSize< Dim > &mortar_size)
 Project variables from a mortar to a face.
 

Detailed Description

Functionality related to discontinuous Galerkin schemes.

Function Documentation

◆ apply_inverse_mass_matrix() [1/2]

template<size_t Dim>
void dg::apply_inverse_mass_matrix ( const gsl::not_null< DataVector * >  data,
const Mesh< Dim > &  mesh 
)

Apply the inverse DG mass matrix to the data.

See also
dg::apply_mass_matrix

◆ apply_inverse_mass_matrix() [2/2]

template<size_t Dim, typename TagsList >
void dg::apply_inverse_mass_matrix ( const gsl::not_null< Variables< TagsList > * >  data,
const Mesh< Dim > &  mesh 
)

Apply the inverse DG mass matrix to the data.

See also
dg::apply_mass_matrix

◆ apply_mass_matrix() [1/2]

template<size_t Dim>
void dg::apply_mass_matrix ( const gsl::not_null< DataVector * >  data,
const Mesh< Dim > &  mesh 
)

Apply the DG mass matrix to the data, in the diagonal mass-matrix approximation ("mass-lumping")

The DG mass matrix is:

\begin{equation} M_{pq} = \int_{\Omega_k} \psi_p(\xi) \psi_q(\xi) \mathrm{d}V \end{equation}

where \(\psi_p(\xi)\) are the basis functions on the element \(\Omega_k\). In the diagonal mass-matrix approximation ("mass-lumping") we evaluate the integral directly on the collocation points, i.e. with a Gauss or Gauss-Lobatto quadrature determined by the element mesh. Then it reduces to:

\begin{equation} M_{pq} \approx \delta_{pq} \prod_{i=1}^d w_{p_i} \end{equation}

where \(d\) is the spatial dimension and \(w_{p_i}\) are the quadrature weights in dimension \(i\). To apply the mass matrix in coordinates different than logical, or to account for a curved background metric, the data can be pre-multiplied with the Jacobian determinant and/or a metric determinant.

Note
The mass-lumping is exact on Legendre-Gauss meshes, but omits a correction term on Legendre-Gauss-Lobatto meshes.

◆ apply_mass_matrix() [2/2]

template<size_t Dim, typename TagsList >
void dg::apply_mass_matrix ( const gsl::not_null< Variables< TagsList > * >  data,
const Mesh< Dim > &  mesh 
)

Apply the DG mass matrix to the data, in the diagonal mass-matrix approximation ("mass-lumping")

The DG mass matrix is:

\begin{equation} M_{pq} = \int_{\Omega_k} \psi_p(\xi) \psi_q(\xi) \mathrm{d}V \end{equation}

where \(\psi_p(\xi)\) are the basis functions on the element \(\Omega_k\). In the diagonal mass-matrix approximation ("mass-lumping") we evaluate the integral directly on the collocation points, i.e. with a Gauss or Gauss-Lobatto quadrature determined by the element mesh. Then it reduces to:

\begin{equation} M_{pq} \approx \delta_{pq} \prod_{i=1}^d w_{p_i} \end{equation}

where \(d\) is the spatial dimension and \(w_{p_i}\) are the quadrature weights in dimension \(i\). To apply the mass matrix in coordinates different than logical, or to account for a curved background metric, the data can be pre-multiplied with the Jacobian determinant and/or a metric determinant.

Note
The mass-lumping is exact on Legendre-Gauss meshes, but omits a correction term on Legendre-Gauss-Lobatto meshes.