12#ifndef MFEM_RESTRICTION
13#define MFEM_RESTRICTION
22class FiniteElementSpace;
25class FaceQuadratureSpace;
34 const real_t a = 1.0)
const override = 0;
45 static const int MaxNbNbr = 16;
64 const real_t a = 1.0)
const override;
96 template <
bool ADD>
void TAddMultTranspose(
const Vector &x,
Vector &y)
const;
123 const real_t a = 1.0)
const override;
134 template <
bool ADD>
void TAddMultTranspose(
const Vector &x,
Vector &y)
const;
188 const real_t a = 1.0)
const override = 0;
249 MFEM_ABORT(
"Not implemented for this restriction operator.");
263 MFEM_ABORT(
"Not implemented for this restriction operator.");
337 const real_t a = 1.0)
const override;
344 const real_t a = 1.0)
const override;
384 const int face_index,
395 const int face_index,
485 const real_t a = 1.0)
const override;
515 const bool keep_nbr_block =
false)
const;
555 Vector &y)
const override;
560 void ComputeScatterIndicesAndOffsets();
566 void ComputeGatherIndices();
569 void EnsureNormalDerivativeRestriction()
const;
586 const int face_index);
596 const int face_index);
607 const int face_index);
616 const int face_index);
627 const int face_index);
639 const int face_index);
764 using Key = std::pair<const DenseMatrix*,int>;
766 using Map = std::map<Key, std::pair<int,const DenseMatrix*>>;
924 const real_t a = 1.0)
const override;
954 const bool keep_nbr_block =
false)
const override;
975 const bool keep_nbr_block =
false)
const override;
994 Vector &ea_data)
const override;
1001 void ComputeScatterIndicesAndOffsets();
1007 void ComputeGatherIndices();
1094 const int face_id2,
const int orientation,
1095 const int size1d,
const int index);
Data type dense matrix using column-major storage.
Abstract base class that defines an interface for element restrictions.
void AddMultTranspose(const Vector &x, Vector &y, const real_t a=1.0) const override=0
Add the E-vector degrees of freedom x to the L-vector degrees of freedom y.
Operator that converts FiniteElementSpace L-vectors to E-vectors.
ElementRestriction(const FiniteElementSpace &, ElementDofOrdering)
const FiniteElementSpace & fes
void Mult(const Vector &x, Vector &y) const override
Operator application: y=A(x).
const Array< int > & GatherMap() const
void FillSparseMatrix(const Vector &mat_ea, SparseMatrix &mat) const
Fill a Sparse Matrix with Element Matrices.
void AddMultTranspose(const Vector &x, Vector &y, const real_t a=1.0) const override
Add the E-vector degrees of freedom x to the L-vector degrees of freedom y.
void FillJAndData(const Vector &ea_data, SparseMatrix &mat) const
void MultTranspose(const Vector &x, Vector &y) const override
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
void MultTransposeUnsigned(const Vector &x, Vector &y) const
Compute MultTranspose without applying signs based on DOF orientations.
const Array< int > & Offsets() const
void MultLeftInverse(const Vector &x, Vector &y) const
const Array< int > & Indices() const
void BooleanMask(Vector &y) const
Fills the E-vector y with boolean values 0.0 and 1.0 such that each each entry of the L-vector is uni...
void MultUnsigned(const Vector &x, Vector &y) const
Compute Mult without applying signs based on DOF orientations.
int FillI(SparseMatrix &mat) const
Base class for operators that extracts Face degrees of freedom.
virtual void AddMultTransposeInPlace(Vector &x, Vector &y) const
Add the face degrees of freedom x to the element degrees of freedom y. Perform the same computation a...
virtual void AddMultTranspose(const Vector &x, Vector &y, const real_t a=1.0) const override=0
Add the face degrees of freedom x to the element degrees of freedom y.
virtual void NormalDerivativeAddMultTranspose(const Vector &x, Vector &y) const
Add the face reference-normal derivative degrees of freedom in x to the element degrees of freedom in...
virtual void AddMultTransposeUnsigned(const Vector &x, Vector &y, const real_t a=1.0) const
Add the face degrees of freedom x to the element degrees of freedom y ignoring the signs from DOF ori...
FaceRestriction(int h, int w)
virtual void NormalDerivativeMult(const Vector &x, Vector &y) const
For each face, sets y to the partial derivative of x with respect to the reference coordinate whose d...
void MultTranspose(const Vector &x, Vector &y) const override
Set the face degrees of freedom in the element degrees of freedom y to the values given in x.
virtual ~FaceRestriction()
void Mult(const Vector &x, Vector &y) const override=0
Extract the face degrees of freedom from x into y.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
This class manages the storage and computation of the interpolations from master (coarse) face to sla...
std::map< Key, std::pair< int, const DenseMatrix * > > Map
The temporary map used to store the different interpolators.
std::pair< const DenseMatrix *, int > Key
const Array< NCInterpConfig > & GetNCFaceInterpConfig() const
Return an array containing the interpolation configuration for each face registered with RegisterFace...
void RegisterFaceCoarseToFineInterpolation(const Mesh::FaceInformation &face, int face_index)
Register the face with face and index face_index as a conforming face for the interpolation of the de...
int GetNumInterpolators() const
Return the total number of interpolators.
Array< NCInterpConfig > nc_interp_config
InterpolationManager()=delete
const Vector & GetInterpolators() const
Return an mfem::Vector containing the interpolators in the following format: face_dofs x face_dofs x ...
void InitializeNCInterpConfig()
void LinearizeInterpolatorMapIntoVector()
Transform the interpolation matrix map into a contiguous memory structure.
const FiniteElementSpace & fes
const ElementDofOrdering ordering
Array< InterpConfig > interp_config
void RegisterFaceConformingInterpolation(const Mesh::FaceInformation &face, int face_index)
Register the face with face and index face_index as a conforming face for the interpolation of the de...
const Array< InterpConfig > & GetFaceInterpConfig() const
Return an array containing the interpolation configuration for each face registered with RegisterFace...
Operator that converts L2 FiniteElementSpace L-vectors to E-vectors.
L2ElementRestriction(const FiniteElementSpace &)
void AddMultTranspose(const Vector &x, Vector &y, const real_t a=1.0) const override
Add the E-vector degrees of freedom x to the L-vector degrees of freedom y.
void FillJAndData(const Vector &ea_data, SparseMatrix &mat) const
void MultTranspose(const Vector &x, Vector &y) const override
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
void FillI(SparseMatrix &mat) const
void Mult(const Vector &x, Vector &y) const override
Operator application: y=A(x).
Operator that extracts Face degrees of freedom for L2 spaces.
void PermuteAndSetFaceDofsGatherIndices2(const Mesh::FaceInformation &face, const int face_index)
Permute and set the gathering indices of elem2 for the interior face described by the face....
virtual void FillI(SparseMatrix &mat, const bool keep_nbr_block=false) const
Fill the I array of SparseMatrix corresponding to the sparsity pattern given by this L2FaceRestrictio...
Array< int > scatter_indices2
void SingleValuedConformingMult(const Vector &x, Vector &y) const
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector. Should only be used with co...
L2FaceRestriction(const FiniteElementSpace &fes, const ElementDofOrdering f_ordering, const FaceType type, const L2FaceValues m, bool build)
Constructs an L2FaceRestriction.
void NormalDerivativeMult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector.
void AddMultTranspose(const Vector &x, Vector &y, const real_t a=1.0) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
void NormalDerivativeAddMultTranspose(const Vector &x, Vector &y) const override
Add the face reference-normal derivative degrees of freedom in x to the element degrees of freedom in...
void Mult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector.
void SingleValuedConformingAddMultTranspose(const Vector &x, Vector &y) const
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector. Should only be used with con...
std::unique_ptr< L2NormalDerivativeFaceRestriction > normal_deriv_restr
void PermuteAndSetSharedFaceDofsScatterIndices2(const Mesh::FaceInformation &face, const int face_index)
Permute and set the scattering indices of elem2 for the shared face described by the face....
void CheckFESpace()
Verify that L2FaceRestriction is built from an L2 FESpace.
Array< int > gather_offsets
virtual void DoubleValuedConformingMult(const Vector &x, Vector &y) const
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector. Should only be used with co...
void PermuteAndSetFaceDofsScatterIndices2(const Mesh::FaceInformation &face, const int face_index)
Permute and set the scattering indices of elem2, and increment the offsets for the face described by ...
Array< int > scatter_indices1
const FiniteElementSpace & fes
void SetBoundaryDofsScatterIndices2(const Mesh::FaceInformation &face, const int face_index)
Set the scattering indices of elem2 for the boundary face described by the face.
void SetFaceDofsScatterIndices1(const Mesh::FaceInformation &face, const int face_index)
Set the scattering indices of elem1, and increment the offsets for the face described by the face....
void SetFaceDofsGatherIndices1(const Mesh::FaceInformation &face, const int face_index)
Set the gathering indices of elem1 for the interior face described by the face.
const ElementDofOrdering ordering
void DoubleValuedConformingAddMultTranspose(const Vector &x, Vector &y) const
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector. Should only be used with con...
virtual void FillJAndData(const Vector &fea_data, SparseMatrix &mat, const bool keep_nbr_block=false) const
Fill the J and Data arrays of the SparseMatrix corresponding to the sparsity pattern given by this L2...
Array< int > gather_indices
virtual void AddFaceMatricesToElementMatrices(const Vector &fea_data, Vector &ea_data) const
This methods adds the DG face matrices to the element matrices.
Operator that extracts face degrees of freedom for L2 nonconforming spaces.
void FillI(SparseMatrix &mat, const bool keep_nbr_block=false) const override
Fill the I array of SparseMatrix corresponding to the sparsity pattern given by this NCL2FaceRestrict...
void AddFaceMatricesToElementMatrices(const Vector &fea_data, Vector &ea_data) const override
This methods adds the DG face matrices to the element matrices.
void FillJAndData(const Vector &fea_data, SparseMatrix &mat, const bool keep_nbr_block=false) const override
Fill the J and Data arrays of the SparseMatrix corresponding to the sparsity pattern given by this NC...
void AddMultTransposeInPlace(Vector &x, Vector &y) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
void SingleValuedNonconformingTransposeInterpolation(const Vector &x) const
Apply a change of basis from fine element basis to coarse element basis for the coarse face dofs....
void SingleValuedNonconformingTransposeInterpolationInPlace(Vector &x) const
Apply a change of basis from fine element basis to coarse element basis for the coarse face dofs....
InterpolationManager interpolations
virtual void DoubleValuedNonconformingMult(const Vector &x, Vector &y) const
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector. Should only be used with no...
void DoubleValuedNonconformingTransposeInterpolation(const Vector &x) const
Apply a change of basis from fine element basis to coarse element basis for the coarse face dofs....
void DoubleValuedNonconformingInterpolation(Vector &x) const
Apply a change of basis from coarse element basis to fine element basis for the coarse face dofs.
void Mult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector.
void AddMultTranspose(const Vector &x, Vector &y, const real_t a=1.0) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
void DoubleValuedNonconformingTransposeInterpolationInPlace(Vector &x) const
Apply a change of basis from fine element basis to coarse element basis for the coarse face dofs....
NCL2FaceRestriction(const FiniteElementSpace &fes, const ElementDofOrdering f_ordering, const FaceType type, const L2FaceValues m, bool build)
Constructs an NCL2FaceRestriction, this is a specialization of a L2FaceRestriction for nonconforming ...
int index(int i, int j, int nx, int ny)
Vector GetLVectorFaceNbrData(const FiniteElementSpace &fes, const Vector &x, FaceType ftype)
Return the face-neighbor data given the L-vector x.
int ToLexOrdering(const int dim, const int face_id, const int size1d, const int index)
Convert a dof face index from Native ordering to lexicographic ordering for quads and hexes.
int PermuteFaceL2(const int dim, const int face_id1, const int face_id2, const int orientation, const int size1d, const int index)
Compute the dof face index of elem2 corresponding to the given dof face index.
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
InterpConfig & operator=(const InterpConfig &rhs)=default
InterpConfig(int master_side, int nc_index)
uint32_t is_non_conforming
InterpConfig(const InterpConfig &)=default
NCInterpConfig & operator=(const NCInterpConfig &rhs)=default
NCInterpConfig(const NCInterpConfig &)=default
NCInterpConfig(int face_index, InterpConfig &config)
uint32_t is_non_conforming
NCInterpConfig(int face_index, int master_side, int nc_index)