12#ifndef MFEM_PRESTRICTION
13#define MFEM_PRESTRICTION
24class ParFiniteElementSpace;
71 const real_t a = 1.0)
const override;
196 const bool keep_nbr_block =
false)
const override;
229 const bool keep_nbr_block =
false)
const override;
235 void ComputeScatterIndicesAndOffsets();
241 void ComputeGatherIndices();
304 const real_t a = 1.0)
const override;
333 const bool keep_nbr_block =
false)
const override;
372 const bool keep_nbr_block =
false)
const override;
379 void ComputeScatterIndicesAndOffsets();
385 void ComputeGatherIndices();
This class manages the storage and computation of the interpolations from master (coarse) face to sla...
Operator that extracts Face degrees of freedom for L2 spaces.
const FiniteElementSpace & fes
Operator that extracts face degrees of freedom for L2 nonconforming spaces.
Abstract parallel finite element space.
Operator that extracts Face degrees of freedom in parallel.
const ParFiniteElementSpace & pfes
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 Mult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector.
void FillJAndData(const Vector &ea_data, SparseMatrix &mat, SparseMatrix &face_mat) const
void FillI(SparseMatrix &mat, const bool keep_nbr_block=false) const override
ParL2FaceRestriction(const ParFiniteElementSpace &pfes_, ElementDofOrdering f_ordering, FaceType type, L2FaceValues m, bool build)
Constructs an ParL2FaceRestriction.
Operator that extracts Face degrees of freedom for NCMesh in parallel.
void NonconformingTransposeInterpolationInPlace(Vector &x) const
Apply a change of basis from fine element basis to coarse element basis for the coarse face dofs.
void AddMultTransposeInPlace(Vector &x, Vector &y) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
ParNCH1FaceRestriction(const ParFiniteElementSpace &fes, ElementDofOrdering f_ordering, FaceType type)
Constructs an ParNCH1FaceRestriction.
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 Mult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector.
void NonconformingInterpolation(Vector &x) const
Apply a change of basis from coarse element basis to fine element basis for the coarse face dofs.
InterpolationManager interpolations
void NonconformingTransposeInterpolation(const Vector &x) const
Apply a change of basis from fine element basis to coarse 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 real_t a=1.0) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
void Mult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector.
ParNCL2FaceRestriction(const ParFiniteElementSpace &fes, ElementDofOrdering f_ordering, FaceType type, L2FaceValues m=L2FaceValues::DoubleValued)
Constructs an ParNCL2FaceRestriction.
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 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 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...
void AddMultTransposeInPlace(Vector &x, Vector &y) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
void FillJAndData(const Vector &fea_data, SparseMatrix &mat, SparseMatrix &face_mat) const
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.