12 #ifndef MFEM_PRESTRICTION 13 #define MFEM_PRESTRICTION 15 #include "../config/config.hpp" 19 #include "restriction.hpp" 24 class ParFiniteElementSpace;
71 const double a = 1.0)
const override;
194 const bool keep_nbr_block =
false)
const override;
227 const bool keep_nbr_block =
false)
const override;
310 const double a = 1.0)
const override;
339 const bool keep_nbr_block =
false)
const override;
378 const bool keep_nbr_block =
false)
const override;
431 #endif // MFEM_USE_MPI 433 #endif // MFEM_PRESTRICTION void Mult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector.
Operator that extracts Face degrees of freedom in parallel.
void DoubleValuedNonconformingMult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector. Should only be used with no...
void FillJAndData(const Vector &fea_data, SparseMatrix &mat, SparseMatrix &face_mat) const
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 NonconformingTransposeInterpolation(const Vector &x) const
Apply a change of basis from fine element basis to coarse element basis for the coarse face dofs...
ParNCH1FaceRestriction(const ParFiniteElementSpace &fes, ElementDofOrdering f_ordering, FaceType type)
Constructs an ParNCH1FaceRestriction.
InterpolationManager interpolations
Operator that extracts face degrees of freedom for L2 nonconforming spaces.
Abstract parallel finite element space.
void Mult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector.
ParL2FaceRestriction(const ParFiniteElementSpace &fes, ElementDofOrdering f_ordering, FaceType type, L2FaceValues m, bool build)
Constructs an ParL2FaceRestriction.
void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
void SingleValuedNonconformingMult(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 AddMultTransposeInPlace(Vector &x, Vector &y) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
const FiniteElementSpace & fes
void AddMultTransposeInPlace(Vector &x, Vector &y) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
Operator that extracts Face degrees of freedom for L2 spaces.
Operator that extracts Face degrees of freedom for NCMesh in parallel.
This class manages the storage and computation of the interpolations from master (coarse) face to sla...
void DoubleValuedConformingMult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector. Should only be used with co...
void FillJAndData(const Vector &ea_data, SparseMatrix &mat, SparseMatrix &face_mat) const
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
void NonconformingTransposeInterpolationInPlace(Vector &x) const
Apply a change of basis from fine element basis to coarse element basis for the coarse face dofs...
void NonconformingInterpolation(Vector &x) const
Apply a change of basis from coarse element basis to fine element basis for the coarse face dofs...
Operator that extracts Face degrees of freedom for NCMesh in parallel.
void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
ParNCL2FaceRestriction(const ParFiniteElementSpace &fes, ElementDofOrdering f_ordering, FaceType type, L2FaceValues m=L2FaceValues::DoubleValued)
Constructs an ParNCL2FaceRestriction.
void FillI(SparseMatrix &mat, const bool keep_nbr_block=false) const override
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 ParNCL2FaceRestr...