MFEM v4.7.0
Finite element discretization library
|
Abstract parallel finite element space. More...
#include <pfespace.hpp>
Public Member Functions | |
ParFiniteElementSpace (const ParFiniteElementSpace &orig, ParMesh *pmesh=NULL, const FiniteElementCollection *fec=NULL) | |
Copy constructor: deep copy all data from orig except the ParMesh, the FiniteElementCollection, and some derived data. | |
ParFiniteElementSpace (const FiniteElementSpace &orig, ParMesh &pmesh, const FiniteElementCollection *fec=NULL) | |
Convert/copy the local (Par)FiniteElementSpace orig to ParFiniteElementSpace: deep copy all data from orig except the Mesh, the FiniteElementCollection, and some derived data. | |
ParFiniteElementSpace (ParMesh *pm, const FiniteElementSpace *global_fes, const int *partitioning, const FiniteElementCollection *f=NULL) | |
Construct the local ParFiniteElementSpace corresponding to the global FE space, global_fes. | |
ParFiniteElementSpace (ParMesh *pm, const FiniteElementCollection *f, int dim=1, int ordering=Ordering::byNODES) | |
ParFiniteElementSpace (ParMesh *pm, NURBSExtension *ext, const FiniteElementCollection *f, int dim=1, int ordering=Ordering::byNODES) | |
Construct a NURBS FE space based on the given NURBSExtension, ext. | |
MPI_Comm | GetComm () const |
int | GetNRanks () const |
int | GetMyRank () const |
ParMesh * | GetParMesh () const |
int | GetDofSign (int i) |
HYPRE_BigInt * | GetDofOffsets () const |
HYPRE_BigInt * | GetTrueDofOffsets () const |
HYPRE_BigInt | GlobalVSize () const |
HYPRE_BigInt | GlobalTrueVSize () const |
int | GetTrueVSize () const override |
Return the number of local vector true dofs. | |
void | GetElementDofs (int i, Array< int > &dofs, DofTransformation &doftrans) const override |
The same as GetElementDofs(), but with a user-allocated DofTransformation object. doftrans must be allocated in advance and will be owned by the caller. The user can use the DofTransformation::GetDofTransformation method on the returned doftrans object to detect if the DofTransformation should actually be used. | |
void | GetBdrElementDofs (int i, Array< int > &dofs, DofTransformation &doftrans) const override |
The same as GetBdrElementDofs(), but with a user-allocated DofTransformation object. doftrans must be allocated in advance and will be owned by the caller. The user can use the DofTransformation::GetDofTransformation method on the returned doftrans object to detect if the DofTransformation should actually be used. | |
int | GetFaceDofs (int i, Array< int > &dofs, int variant=0) const override |
const FiniteElement * | GetFE (int i) const override |
const FaceRestriction * | GetFaceRestriction (ElementDofOrdering f_ordering, FaceType type, L2FaceValues mul=L2FaceValues::DoubleValued) const override |
void | GetSharedEdgeDofs (int group, int ei, Array< int > &dofs) const |
void | GetSharedTriangleDofs (int group, int fi, Array< int > &dofs) const |
void | GetSharedQuadrilateralDofs (int group, int fi, Array< int > &dofs) const |
HypreParMatrix * | Dof_TrueDof_Matrix () const |
The true dof-to-dof interpolation matrix. | |
HypreParMatrix * | GetPartialConformingInterpolation () |
For a non-conforming mesh, construct and return the interpolation matrix from the partially conforming true dofs to the local dofs. | |
HypreParVector * | NewTrueDofVector () |
void | DivideByGroupSize (real_t *vec) |
Scale a vector of true dofs. | |
GroupCommunicator & | GroupComm () |
Return a reference to the internal GroupCommunicator (on VDofs) | |
const GroupCommunicator & | GroupComm () const |
Return a const reference to the internal GroupCommunicator (on VDofs) | |
GroupCommunicator * | ScalarGroupComm () |
Return a new GroupCommunicator on scalar dofs, i.e. for VDim = 1. | |
void | Synchronize (Array< int > &ldof_marker) const |
Given an integer array on the local degrees of freedom, perform a bitwise OR between the shared dofs. | |
void | GetEssentialVDofs (const Array< int > &bdr_attr_is_ess, Array< int > &ess_dofs, int component=-1) const override |
Determine the boundary degrees of freedom. | |
void | GetEssentialTrueDofs (const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const override |
int | GetLocalTDofNumber (int ldof) const |
HYPRE_BigInt | GetGlobalTDofNumber (int ldof) const |
Returns the global tdof number of the given local degree of freedom. | |
HYPRE_BigInt | GetGlobalScalarTDofNumber (int sldof) |
HYPRE_BigInt | GetMyDofOffset () const |
HYPRE_BigInt | GetMyTDofOffset () const |
const Operator * | GetProlongationMatrix () const override |
The returned Operator is owned by the FiniteElementSpace. | |
const Operator * | GetRestrictionOperator () const override |
const SparseMatrix * | GetRestrictionMatrix () const override |
Get the R matrix which restricts a local dof vector to true dof vector. | |
void | ExchangeFaceNbrData () |
int | GetFaceNbrVSize () const |
void | GetFaceNbrElementVDofs (int i, Array< int > &vdofs, DofTransformation &doftrans) const |
DofTransformation * | GetFaceNbrElementVDofs (int i, Array< int > &vdofs) const |
void | GetFaceNbrFaceVDofs (int i, Array< int > &vdofs) const |
const FiniteElement * | GetFaceNbrFE (int i) const |
const FiniteElement * | GetFaceNbrFaceFE (int i) const |
const HYPRE_BigInt * | GetFaceNbrGlobalDofMap () |
ElementTransformation * | GetFaceNbrElementTransformation (int i) const |
void | Lose_Dof_TrueDof_Matrix () |
void | LoseDofOffsets () |
void | LoseTrueDofOffsets () |
bool | Conforming () const |
bool | Nonconforming () const |
bool | SharedNDTriangleDofs () const |
void | GetTrueTransferOperator (const FiniteElementSpace &coarse_fes, OperatorHandle &T) const override |
Construct and return an Operator that can be used to transfer true-dof data from coarse_fes, defined on a coarse mesh, to this FE space, defined on a refined mesh. | |
void | Update (bool want_transform=true) override |
void | UpdatesFinished () override |
Free ParGridFunction transformation matrix (if any), to save memory. | |
virtual | ~ParFiniteElementSpace () |
void | PrintPartitionStats () |
int | TrueVSize () const |
Obsolete, kept for backward compatibility. | |
DofTransformation * | GetElementDofs (int elem, Array< int > &dofs) const |
DofTransformation * | GetBdrElementDofs (int bel, Array< int > &dofs) const |
Public Member Functions inherited from mfem::FiniteElementSpace | |
FiniteElementSpace () | |
Default constructor: the object is invalid until initialized using the method Load(). | |
FiniteElementSpace (const FiniteElementSpace &orig, Mesh *mesh=NULL, const FiniteElementCollection *fec=NULL) | |
Copy constructor: deep copy all data from orig except the Mesh, the FiniteElementCollection, and some derived data. | |
FiniteElementSpace (Mesh *mesh, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES) | |
FiniteElementSpace (Mesh *mesh, NURBSExtension *ext, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES) | |
Construct a NURBS FE space based on the given NURBSExtension, ext. | |
FiniteElementSpace & | operator= (const FiniteElementSpace &)=delete |
Copy assignment not supported. | |
Mesh * | GetMesh () const |
Returns the mesh. | |
const NURBSExtension * | GetNURBSext () const |
NURBSExtension * | GetNURBSext () |
NURBSExtension * | StealNURBSext () |
bool | Conforming () const |
bool | Nonconforming () const |
void | SetElementOrder (int i, int p) |
Sets the order of the i'th finite element. | |
int | GetElementOrder (int i) const |
Returns the order of the i'th finite element. | |
int | GetMaxElementOrder () const |
Return the maximum polynomial order. | |
bool | IsVariableOrder () const |
Returns true if the space contains elements of varying polynomial orders. | |
const SparseMatrix * | GetConformingProlongation () const |
The returned SparseMatrix is owned by the FiniteElementSpace. | |
const SparseMatrix * | GetConformingRestriction () const |
The returned SparseMatrix is owned by the FiniteElementSpace. | |
const SparseMatrix * | GetHpConformingRestriction () const |
The returned SparseMatrix is owned by the FiniteElementSpace. | |
const Operator * | GetRestrictionTransposeOperator () const |
Return an operator that performs the transpose of GetRestrictionOperator. | |
virtual const SparseMatrix * | GetHpRestrictionMatrix () const |
The returned SparseMatrix is owned by the FiniteElementSpace. | |
const ElementRestrictionOperator * | GetElementRestriction (ElementDofOrdering e_ordering) const |
Return an Operator that converts L-vectors to E-vectors. | |
const QuadratureInterpolator * | GetQuadratureInterpolator (const IntegrationRule &ir) const |
Return a QuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivatives (Q-vectors). | |
const QuadratureInterpolator * | GetQuadratureInterpolator (const QuadratureSpace &qs) const |
Return a QuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivatives (Q-vectors). | |
const FaceQuadratureInterpolator * | GetFaceQuadratureInterpolator (const IntegrationRule &ir, FaceType type) const |
Return a FaceQuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivatives (Q-vectors). | |
int | GetOrder (int i) const |
Returns the polynomial degree of the i'th finite element. | |
int | GetEdgeOrder (int edge, int variant=0) const |
int | GetFaceOrder (int face, int variant=0) const |
Returns the polynomial degree of the i'th face finite element. | |
int | GetVDim () const |
Returns vector dimension. | |
int | GetNDofs () const |
Returns number of degrees of freedom. This is the number of Local Degrees of Freedom. | |
int | GetVSize () const |
Return the number of vector dofs, i.e. GetNDofs() x GetVDim(). | |
int | GetNConformingDofs () const |
int | GetConformingVSize () const |
Ordering::Type | GetOrdering () const |
Return the ordering method. | |
const FiniteElementCollection * | FEColl () const |
int | GetNVDofs () const |
Number of all scalar vertex dofs. | |
int | GetNEDofs () const |
Number of all scalar edge-interior dofs. | |
int | GetNFDofs () const |
Number of all scalar face-interior dofs. | |
int | GetNV () const |
Returns number of vertices in the mesh. | |
int | GetNE () const |
Returns number of elements in the mesh. | |
int | GetNF () const |
Returns number of faces (i.e. co-dimension 1 entities) in the mesh. | |
int | GetNBE () const |
Returns number of boundary elements in the mesh. | |
int | GetNFbyType (FaceType type) const |
Returns the number of faces according to the requested type. | |
int | GetElementType (int i) const |
Returns the type of element i. | |
void | GetElementVertices (int i, Array< int > &vertices) const |
Returns the vertices of element i. | |
int | GetBdrElementType (int i) const |
Returns the type of boundary element i. | |
ElementTransformation * | GetElementTransformation (int i) const |
Returns ElementTransformation for the i-th element. | |
void | GetElementTransformation (int i, IsoparametricTransformation *ElTr) |
Returns the transformation defining the i-th element in the user-defined variable ElTr. | |
ElementTransformation * | GetBdrElementTransformation (int i) const |
Returns ElementTransformation for the i-th boundary element. | |
int | GetAttribute (int i) const |
int | GetBdrAttribute (int i) const |
void | GetPatchDofs (int patch, Array< int > &dofs) const |
Returns indices of degrees of freedom for NURBS patch index patch. Cartesian ordering is used, for the tensor-product degrees of freedom. | |
MFEM_DEPRECATED void | RebuildElementToDofTable () |
( | |
void | ReorderElementToDofTable () |
Reorder the scalar DOFs based on the element ordering. | |
const Table * | GetElementToFaceOrientationTable () const |
const Table & | GetElementToDofTable () const |
Return a reference to the internal Table that stores the lists of scalar dofs, for each mesh element, as returned by GetElementDofs(). | |
const Table & | GetBdrElementToDofTable () const |
Return a reference to the internal Table that stores the lists of scalar dofs, for each boundary mesh element, as returned by GetBdrElementDofs(). | |
const Table & | GetFaceToDofTable () const |
Return a reference to the internal Table that stores the lists of scalar dofs, for each face in the mesh, as returned by GetFaceDofs(). In this context, "face" refers to a (dim-1)-dimensional mesh entity. | |
void | BuildDofToArrays () |
Initialize internal data that enables the use of the methods GetElementForDof() and GetLocalDofForDof(). | |
int | GetElementForDof (int i) const |
Return the index of the first element that contains dof i. | |
int | GetLocalDofForDof (int i) const |
Return the local dof index in the first element that contains dof i. | |
const FiniteElement * | GetBE (int i) const |
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th boundary face in the mesh object. | |
const FiniteElement * | GetFaceElement (int i) const |
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th face in the mesh object. Faces in this case refer to the MESHDIM-1 primitive so in 2D they are segments and in 1D they are points. | |
const FiniteElement * | GetEdgeElement (int i, int variant=0) const |
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th edge in the mesh object. | |
const FiniteElement * | GetTraceElement (int i, Geometry::Type geom_type) const |
Return the trace element from element 'i' to the given 'geom_type'. | |
void | GetBoundaryTrueDofs (Array< int > &boundary_dofs, int component=-1) |
Get a list of all boundary true dofs, boundary_dofs. For spaces with 'vdim' > 1, the 'component' parameter can be used to restricts the marked tDOFs to the specified component. Equivalent to FiniteElementSpace::GetEssentialTrueDofs with all boundary attributes marked as essential. | |
void | ConvertToConformingVDofs (const Array< int > &dofs, Array< int > &cdofs) |
For a partially conforming FE space, convert a marker array (nonzero entries are true) on the partially conforming dofs to a marker array on the conforming dofs. A conforming dofs is marked iff at least one of its dependent dofs is marked. | |
void | ConvertFromConformingVDofs (const Array< int > &cdofs, Array< int > &dofs) |
For a partially conforming FE space, convert a marker array (nonzero entries are true) on the conforming dofs to a marker array on the (partially conforming) dofs. A dof is marked iff it depends on a marked conforming dofs, where dependency is defined by the ConformingRestriction matrix; in other words, a dof is marked iff it corresponds to a marked conforming dof. | |
SparseMatrix * | D2C_GlobalRestrictionMatrix (FiniteElementSpace *cfes) |
Generate the global restriction matrix from a discontinuous FE space to the continuous FE space of the same polynomial degree. | |
SparseMatrix * | D2Const_GlobalRestrictionMatrix (FiniteElementSpace *cfes) |
Generate the global restriction matrix from a discontinuous FE space to the piecewise constant FE space. | |
SparseMatrix * | H2L_GlobalRestrictionMatrix (FiniteElementSpace *lfes) |
Construct the restriction matrix from the FE space given by (*this) to the lower degree FE space given by (*lfes) which is defined on the same mesh. | |
void | GetTransferOperator (const FiniteElementSpace &coarse_fes, OperatorHandle &T) const |
Construct and return an Operator that can be used to transfer GridFunction data from coarse_fes, defined on a coarse mesh, to this FE space, defined on a refined mesh. | |
const Operator * | GetUpdateOperator () |
Get the GridFunction update operator. | |
void | GetUpdateOperator (OperatorHandle &T) |
Return the update operator in the given OperatorHandle, T. | |
void | SetUpdateOperatorOwner (bool own) |
Set the ownership of the update operator: if set to false, the Operator returned by GetUpdateOperator() must be deleted outside the FiniteElementSpace. | |
void | SetUpdateOperatorType (Operator::Type tid) |
Specify the Operator::Type to be used by the update operators. | |
long | GetSequence () const |
bool | IsDGSpace () const |
Return whether or not the space is discontinuous (L2) | |
void | SetRelaxedHpConformity (bool relaxed=true) |
void | Save (std::ostream &out) const |
Save finite element space to output stream out. | |
FiniteElementCollection * | Load (Mesh *m, std::istream &input) |
Read a FiniteElementSpace from a stream. The returned FiniteElementCollection is owned by the caller. | |
virtual | ~FiniteElementSpace () |
DofTransformation * | GetElementDofs (int elem, Array< int > &dofs) const |
Returns indices of degrees of freedom of element 'elem'. The returned indices are offsets into an ldof vector. See also GetElementVDofs(). | |
DofTransformation * | GetBdrElementDofs (int bel, Array< int > &dofs) const |
Returns indices of degrees of freedom for boundary element 'bel'. The returned indices are offsets into an ldof vector. See also GetBdrElementVDofs(). | |
int | GetEdgeDofs (int edge, Array< int > &dofs, int variant=0) const |
Returns the indices of the degrees of freedom for the specified edge, including the DOFs for the vertices of the edge. | |
void | GetVertexDofs (int i, Array< int > &dofs) const |
Returns the indices of the degrees of freedom for the specified vertices. | |
void | GetElementInteriorDofs (int i, Array< int > &dofs) const |
Returns the indices of the degrees of freedom for the interior of the specified element. | |
int | GetNumElementInteriorDofs (int i) const |
Returns the number of degrees of freedom associated with the interior of the specified element. | |
void | GetFaceInteriorDofs (int i, Array< int > &dofs) const |
Returns the indices of the degrees of freedom for the interior of the specified face. | |
void | GetEdgeInteriorDofs (int i, Array< int > &dofs) const |
Returns the indices of the degrees of freedom for the interior of the specified edge. | |
void | GetVDofs (int vd, Array< int > &dofs, int ndofs=-1) const |
Returns the indices of all of the VDofs for the specified dimension 'vd'. | |
void | DofsToVDofs (Array< int > &dofs, int ndofs=-1) const |
Compute the full set of vdofs corresponding to each entry in dofs. | |
void | DofsToVDofs (int vd, Array< int > &dofs, int ndofs=-1) const |
Compute the set of vdofs corresponding to each entry in dofs for the given vector index vd. | |
int | DofToVDof (int dof, int vd, int ndofs=-1) const |
Compute a single vdof corresponding to the index dof and the vector index vd. | |
int | VDofToDof (int vdof) const |
Compute the inverse of the Dof to VDof mapping for a single index vdof. | |
DofTransformation * | GetElementVDofs (int i, Array< int > &vdofs) const |
Returns indices of degrees of freedom for the i'th element. The returned indices are offsets into an ldof vector with vdim not necessarily equal to 1. The returned indices are always ordered byNODES, irrespective of whether the space is byNODES or byVDIM. See also GetElementDofs(). | |
void | GetElementVDofs (int i, Array< int > &vdofs, DofTransformation &doftrans) const |
The same as GetElementVDofs(), but with a user-allocated DofTransformation object. doftrans must be allocated in advance and will be owned by the caller. The user can use the DofTransformation::GetDofTransformation method on the returned doftrans object to detect if the DofTransformation should actually be used. | |
DofTransformation * | GetBdrElementVDofs (int i, Array< int > &vdofs) const |
Returns indices of degrees of freedom for i'th boundary element. The returned indices are offsets into an ldof vector with vdim not necessarily equal to 1. See also GetBdrElementDofs(). | |
void | GetBdrElementVDofs (int i, Array< int > &vdofs, DofTransformation &doftrans) const |
The same as GetBdrElementVDofs(), but with a user-allocated DofTransformation object. doftrans must be allocated in advance and will be owned by the caller. The user can use the DofTransformation::GetDofTransformation method on the returned doftrans object to detect if the DofTransformation should actually be used. | |
void | GetPatchVDofs (int i, Array< int > &vdofs) const |
Returns indices of degrees of freedom in vdofs for NURBS patch i. | |
void | GetFaceVDofs (int i, Array< int > &vdofs) const |
Returns the indices of the degrees of freedom for the specified face, including the DOFs for the edges and the vertices of the face. | |
void | GetEdgeVDofs (int i, Array< int > &vdofs) const |
Returns the indices of the degrees of freedom for the specified edge, including the DOFs for the vertices of the edge. | |
void | GetVertexVDofs (int i, Array< int > &vdofs) const |
Returns the indices of the degrees of freedom for the specified vertices. | |
void | GetElementInteriorVDofs (int i, Array< int > &vdofs) const |
Returns the indices of the degrees of freedom for the interior of the specified element. | |
void | GetEdgeInteriorVDofs (int i, Array< int > &vdofs) const |
Returns the indices of the degrees of freedom for the interior of the specified edge. | |
Public Attributes | |
int | num_face_nbr_dofs |
Table | face_nbr_element_dof |
Table | face_nbr_element_fos |
Table | face_nbr_ldof |
Array< HYPRE_BigInt > | face_nbr_glob_dof_map |
Table | send_face_nbr_ldof |
Additional Inherited Members | |
Static Public Member Functions inherited from mfem::FiniteElementSpace | |
static void | AdjustVDofs (Array< int > &vdofs) |
Remove the orientation information encoded into an array of dofs Some basis function types have a relative orientation associated with degrees of freedom shared between neighboring elements, see ldof for more information. An orientation mismatch is indicated in the dof indices by a negative index value. This method replaces such negative indices with the corresponding positive offsets. | |
static int | EncodeDof (int entity_base, int idx) |
Helper to encode a sign flip into a DOF index (for Hcurl/Hdiv shapes). | |
static int | DecodeDof (int dof) |
Helper to return the DOF associated with a sign encoded DOF. | |
static int | DecodeDof (int dof, real_t &sign) |
Helper to determine the DOF and sign of a sign encoded DOF. | |
static void | MarkerToList (const Array< int > &marker, Array< int > &list) |
Convert a Boolean marker array to a list containing all marked indices. | |
static void | ListToMarker (const Array< int > &list, int marker_size, Array< int > &marker, int mark_val=-1) |
Convert an array of indices (list) to a Boolean marker array where all indices in the list are marked with the given value and the rest are set to zero. | |
Protected Types inherited from mfem::FiniteElementSpace | |
using | key_face = std::tuple<bool, ElementDofOrdering, FaceType, L2FaceValues> |
The face restriction operators, see GetFaceRestriction(). | |
using | map_L2F = std::unordered_map<const key_face,FaceRestriction*,key_hash> |
typedef std::uint64_t | VarOrderBits |
Bit-mask representing a set of orders needed by an edge/face. | |
Protected Member Functions inherited from mfem::FiniteElementSpace | |
void | UpdateNURBS () |
void | Construct () |
void | Destroy () |
void | ConstructDoFTransArray () |
void | DestroyDoFTransArray () |
void | BuildElementToDofTable () const |
void | BuildBdrElementToDofTable () const |
void | BuildFaceToDofTable () const |
void | BuildNURBSFaceToDofTable () const |
Generates partial face_dof table for a NURBS space. | |
int | GetElementOrderImpl (int i) const |
Return element order: internal version of GetElementOrder without checks. | |
void | CalcEdgeFaceVarOrders (Array< VarOrderBits > &edge_orders, Array< VarOrderBits > &face_orders) const |
int | MakeDofTable (int ent_dim, const Array< int > &entity_orders, Table &entity_dofs, Array< char > *var_ent_order) |
int | FindDofs (const Table &var_dof_table, int row, int ndof) const |
Search row of a DOF table for a DOF set of size 'ndof', return first DOF. | |
int | FindEdgeDof (int edge, int ndof) const |
int | FindFaceDof (int face, int ndof) const |
Similar to FindEdgeDof, but used for mixed meshes too. | |
int | FirstFaceDof (int face, int variant=0) const |
int | GetNVariants (int entity, int index) const |
Return number of possible DOF variants for edge/face (var. order spaces). | |
int | GetEntityDofs (int entity, int index, Array< int > &dofs, Geometry::Type master_geom=Geometry::INVALID, int variant=0) const |
Helper to get vertex, edge or face DOFs (entity=0,1,2 resp.). | |
int | GetEntityVDofs (int entity, int index, Array< int > &dofs, Geometry::Type master_geom=Geometry::INVALID, int variant=0) const |
Helper to get vertex, edge or face VDOFs (entity=0,1,2 resp.). | |
int | GetDegenerateFaceDofs (int index, Array< int > &dofs, Geometry::Type master_geom, int variant) const |
int | GetNumBorderDofs (Geometry::Type geom, int order) const |
void | BuildConformingInterpolation () const |
Calculate the cP and cR matrices for a nonconforming mesh. | |
void | AddEdgeFaceDependencies (SparseMatrix &deps, Array< int > &master_dofs, const FiniteElement *master_fe, Array< int > &slave_dofs, int slave_face, const DenseMatrix *pm) const |
void | MakeVDimMatrix (SparseMatrix &mat) const |
Replicate 'mat' in the vector dimension, according to vdim ordering mode. | |
SparseMatrix * | RefinementMatrix_main (const int coarse_ndofs, const Table &coarse_elem_dof, const Table *coarse_elem_fos, const DenseTensor localP[]) const |
void | GetLocalRefinementMatrices (Geometry::Type geom, DenseTensor &localP) const |
void | GetLocalDerefinementMatrices (Geometry::Type geom, DenseTensor &localR) const |
SparseMatrix * | RefinementMatrix (int old_ndofs, const Table *old_elem_dof, const Table *old_elem_fos) |
SparseMatrix * | DerefinementMatrix (int old_ndofs, const Table *old_elem_dof, const Table *old_elem_fos) |
Calculate GridFunction restriction matrix after mesh derefinement. | |
void | GetLocalRefinementMatrices (const FiniteElementSpace &coarse_fes, Geometry::Type geom, DenseTensor &localP) const |
Return in localP the local refinement matrices that map between fespaces after mesh refinement. | |
void | Constructor (Mesh *mesh, NURBSExtension *ext, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES) |
Help function for constructors + Load(). | |
void | UpdateElementOrders () |
Resize the elem_order array on mesh change. | |
Static Protected Member Functions inherited from mfem::FiniteElementSpace | |
static int | MinOrder (VarOrderBits bits) |
Return the minimum order (least significant bit set) in the bit mask. | |
static void | AddDependencies (SparseMatrix &deps, Array< int > &master_dofs, Array< int > &slave_dofs, DenseMatrix &I, int skipfirst=0) |
static bool | DofFinalizable (int dof, const Array< bool > &finalized, const SparseMatrix &deps) |
Protected Attributes inherited from mfem::FiniteElementSpace | |
Mesh * | mesh |
The mesh that FE space lives on (not owned). | |
const FiniteElementCollection * | fec |
Associated FE collection (not owned). | |
int | vdim |
Vector dimension (number of unknowns per degree of freedom). | |
Ordering::Type | ordering |
int | ndofs |
Number of degrees of freedom. Number of unknowns is ndofs * vdim. | |
Array< char > | elem_order |
int | nvdofs |
int | nedofs |
int | nfdofs |
int | nbdofs |
int | uni_fdof |
int * | bdofs |
internal DOFs of elements if mixed/var-order; NULL otherwise | |
Table | var_edge_dofs |
Table | var_face_dofs |
NOTE: also used for spaces with mixed faces. | |
Array< char > | var_edge_orders |
Array< char > | var_face_orders |
Table * | elem_dof |
Table * | elem_fos |
Table * | bdr_elem_dof |
Table * | bdr_elem_fos |
Table * | face_dof |
Array< int > | dof_elem_array |
Array< int > | dof_ldof_array |
NURBSExtension * | NURBSext |
int | own_ext |
Array< int > | face_to_be |
Array< StatelessDofTransformation * > | DoFTransArray |
DofTransformation | DoFTrans |
std::unique_ptr< SparseMatrix > | cP |
std::unique_ptr< SparseMatrix > | cR |
Conforming restriction matrix such that cR.cP=I. | |
std::unique_ptr< SparseMatrix > | cR_hp |
A version of the conforming restriction matrix for variable-order spaces. | |
bool | cP_is_set |
std::unique_ptr< Operator > | R_transpose |
Operator computing the action of the transpose of the restriction. | |
OperatorHandle | Th |
Transformation to apply to GridFunctions after space Update(). | |
OperatorHandle | L2E_nat |
The element restriction operators, see GetElementRestriction(). | |
OperatorHandle | L2E_lex |
map_L2F | L2F |
Array< QuadratureInterpolator * > | E2Q_array |
Array< FaceQuadratureInterpolator * > | E2IFQ_array |
Array< FaceQuadratureInterpolator * > | E2BFQ_array |
long | sequence |
long | mesh_sequence |
bool | orders_changed |
True if at least one element order changed (variable-order space only). | |
bool | relaxed_hp |
Static Protected Attributes inherited from mfem::FiniteElementSpace | |
static constexpr int | MaxVarOrder = 8*sizeof(VarOrderBits) - 1 |
Abstract parallel finite element space.
Definition at line 28 of file pfespace.hpp.
mfem::ParFiniteElementSpace::ParFiniteElementSpace | ( | const ParFiniteElementSpace & | orig, |
ParMesh * | pmesh = NULL, | ||
const FiniteElementCollection * | fec = NULL ) |
Copy constructor: deep copy all data from orig except the ParMesh, the FiniteElementCollection, and some derived data.
If the pmesh or fec pointers are NULL (default), then the new ParFiniteElementSpace will reuse the respective pointers from orig. If any of these pointers is not NULL, the given pointer will be used instead of the one used by orig.
Definition at line 29 of file pfespace.cpp.
mfem::ParFiniteElementSpace::ParFiniteElementSpace | ( | const FiniteElementSpace & | orig, |
ParMesh & | pmesh, | ||
const FiniteElementCollection * | fec = NULL ) |
Convert/copy the local (Par)FiniteElementSpace orig to ParFiniteElementSpace: deep copy all data from orig except the Mesh, the FiniteElementCollection, and some derived data.
Definition at line 37 of file pfespace.cpp.
mfem::ParFiniteElementSpace::ParFiniteElementSpace | ( | ParMesh * | pm, |
const FiniteElementSpace * | global_fes, | ||
const int * | partitioning, | ||
const FiniteElementCollection * | f = NULL ) |
Construct the local ParFiniteElementSpace corresponding to the global FE space, global_fes.
The parameter pm is the local ParMesh obtained by decomposing the global Mesh used by global_fes. The array partitioning represents the parallel decomposition - it maps global element ids to MPI ranks. If the FiniteElementCollection, f, is NULL (default), the FE collection used by global_fes will be reused. If f is not NULL, it must be the same as, or a copy of, the FE collection used by global_fes.
Definition at line 45 of file pfespace.cpp.
mfem::ParFiniteElementSpace::ParFiniteElementSpace | ( | ParMesh * | pm, |
const FiniteElementCollection * | f, | ||
int | dim = 1, | ||
int | ordering = Ordering::byNODES ) |
Definition at line 61 of file pfespace.cpp.
mfem::ParFiniteElementSpace::ParFiniteElementSpace | ( | ParMesh * | pm, |
NURBSExtension * | ext, | ||
const FiniteElementCollection * | f, | ||
int | dim = 1, | ||
int | ordering = Ordering::byNODES ) |
Construct a NURBS FE space based on the given NURBSExtension, ext.
The parameter ext will be deleted by this constructor, replaced by a ParNURBSExtension owned by the ParFiniteElementSpace.
Definition at line 68 of file pfespace.cpp.
|
inlinevirtual |
Definition at line 431 of file pfespace.hpp.
|
inline |
Definition at line 409 of file pfespace.hpp.
void mfem::ParFiniteElementSpace::DivideByGroupSize | ( | real_t * | vec | ) |
Scale a vector of true dofs.
Definition at line 974 of file pfespace.cpp.
|
inline |
The true dof-to-dof interpolation matrix.
Definition at line 327 of file pfespace.hpp.
void mfem::ParFiniteElementSpace::ExchangeFaceNbrData | ( | ) |
Definition at line 1226 of file pfespace.cpp.
DofTransformation * mfem::FiniteElementSpace::GetBdrElementDofs | ( | int | bel, |
Array< int > & | dofs ) const |
Returns indexes of degrees of freedom for i'th boundary element and returns the DofTransformation data in a user-provided object.
Definition at line 833 of file fespace.cpp.
|
overridevirtual |
The same as GetBdrElementDofs(), but with a user-allocated DofTransformation object. doftrans must be allocated in advance and will be owned by the caller. The user can use the DofTransformation::GetDofTransformation method on the returned doftrans object to detect if the DofTransformation should actually be used.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 493 of file pfespace.cpp.
|
inline |
Definition at line 273 of file pfespace.hpp.
|
inline |
Definition at line 281 of file pfespace.hpp.
|
inline |
Definition at line 279 of file pfespace.hpp.
DofTransformation * mfem::FiniteElementSpace::GetElementDofs | ( | int | elem, |
Array< int > & | dofs ) const |
Returns indexes of degrees of freedom in array dofs for i'th element and returns the DofTransformation data in a user-provided object.
Definition at line 810 of file fespace.cpp.
|
overridevirtual |
The same as GetElementDofs(), but with a user-allocated DofTransformation object. doftrans must be allocated in advance and will be owned by the caller. The user can use the DofTransformation::GetDofTransformation method on the returned doftrans object to detect if the DofTransformation should actually be used.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 468 of file pfespace.cpp.
|
overridevirtual |
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in the array bdr_attr_is_ess.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 1031 of file pfespace.cpp.
|
overridevirtual |
Determine the boundary degrees of freedom.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 1020 of file pfespace.cpp.
|
overridevirtual |
Returns the indexes of the degrees of freedom for i'th face including the dofs for the edges and the vertices of the face.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 518 of file pfespace.cpp.
|
inline |
Definition at line 402 of file pfespace.hpp.
DofTransformation * mfem::ParFiniteElementSpace::GetFaceNbrElementVDofs | ( | int | i, |
Array< int > & | vdofs ) const |
Definition at line 1483 of file pfespace.cpp.
void mfem::ParFiniteElementSpace::GetFaceNbrElementVDofs | ( | int | i, |
Array< int > & | vdofs, | ||
DofTransformation & | doftrans ) const |
Definition at line 1467 of file pfespace.cpp.
const FiniteElement * mfem::ParFiniteElementSpace::GetFaceNbrFaceFE | ( | int | i | ) | const |
Definition at line 1533 of file pfespace.cpp.
void mfem::ParFiniteElementSpace::GetFaceNbrFaceVDofs | ( | int | i, |
Array< int > & | vdofs ) const |
Definition at line 1491 of file pfespace.cpp.
const FiniteElement * mfem::ParFiniteElementSpace::GetFaceNbrFE | ( | int | i | ) | const |
Definition at line 1519 of file pfespace.cpp.
|
inline |
Definition at line 401 of file pfespace.hpp.
|
inline |
Definition at line 394 of file pfespace.hpp.
|
overridevirtual |
Returns an Operator that converts L-vectors to E-vectors on each face. The parallel version is different from the serial one because of the presence of shared faces. Shared faces are treated as interior faces, the returned operator handles the communication needed to get the shared face values from other MPI ranks
Reimplemented from mfem::FiniteElementSpace.
Definition at line 541 of file pfespace.cpp.
|
overridevirtual |
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in the mesh object. If i is greater than or equal to the number of local mesh elements, i will be interpreted as a shifted index of a face neighbor element.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 534 of file pfespace.cpp.
HYPRE_BigInt mfem::ParFiniteElementSpace::GetGlobalScalarTDofNumber | ( | int | sldof | ) |
Returns the global tdof number of the given local degree of freedom in the scalar version of the current finite element space. The input should be a scalar local dof.
Definition at line 1115 of file pfespace.cpp.
HYPRE_BigInt mfem::ParFiniteElementSpace::GetGlobalTDofNumber | ( | int | ldof | ) | const |
Returns the global tdof number of the given local degree of freedom.
Definition at line 1092 of file pfespace.cpp.
int mfem::ParFiniteElementSpace::GetLocalTDofNumber | ( | int | ldof | ) | const |
If the given ldof is owned by the current processor, return its local tdof number, otherwise return -1
Definition at line 1071 of file pfespace.cpp.
HYPRE_BigInt mfem::ParFiniteElementSpace::GetMyDofOffset | ( | ) | const |
Definition at line 1152 of file pfespace.cpp.
|
inline |
Definition at line 275 of file pfespace.hpp.
HYPRE_BigInt mfem::ParFiniteElementSpace::GetMyTDofOffset | ( | ) | const |
Definition at line 1157 of file pfespace.cpp.
|
inline |
Definition at line 274 of file pfespace.hpp.
|
inline |
Definition at line 277 of file pfespace.hpp.
HypreParMatrix * mfem::ParFiniteElementSpace::GetPartialConformingInterpolation | ( | ) |
For a non-conforming mesh, construct and return the interpolation matrix from the partially conforming true dofs to the local dofs.
Definition at line 963 of file pfespace.cpp.
|
overridevirtual |
The returned Operator is owned by the FiniteElementSpace.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 1162 of file pfespace.cpp.
|
inlineoverridevirtual |
Get the R matrix which restricts a local dof vector to true dof vector.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 389 of file pfespace.hpp.
|
overridevirtual |
Get an Operator that performs the action of GetRestrictionMatrix(), but potentially with a non-assembled optimized matrix-free implementation.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 1193 of file pfespace.cpp.
void mfem::ParFiniteElementSpace::GetSharedEdgeDofs | ( | int | group, |
int | ei, | ||
Array< int > & | dofs ) const |
Definition at line 583 of file pfespace.cpp.
void mfem::ParFiniteElementSpace::GetSharedQuadrilateralDofs | ( | int | group, |
int | fi, | ||
Array< int > & | dofs ) const |
Definition at line 630 of file pfespace.cpp.
void mfem::ParFiniteElementSpace::GetSharedTriangleDofs | ( | int | group, |
int | fi, | ||
Array< int > & | dofs ) const |
Definition at line 606 of file pfespace.cpp.
|
inline |
Definition at line 282 of file pfespace.hpp.
|
overridevirtual |
Construct and return an Operator that can be used to transfer true-dof data from coarse_fes, defined on a coarse mesh, to this FE space, defined on a refined mesh.
This method calls GetTransferOperator() and multiplies the result by the prolongation operator of coarse_fes on the right, and by the restriction operator of this FE space on the left.
The Operator::Type of T can be set to request an Operator of the set type. In serial, the supported types are: Operator::MFEM_SPARSEMAT and Operator::ANY_TYPE (matrix-free). In parallel, the supported types are: Operator::Hypre_ParCSR and Operator::ANY_TYPE. Any other type is treated as Operator::ANY_TYPE: the operator representation choice is made by this method.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 3387 of file pfespace.cpp.
|
inlineoverridevirtual |
Return the number of local vector true dofs.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 289 of file pfespace.hpp.
|
inline |
Definition at line 285 of file pfespace.hpp.
|
inline |
Definition at line 283 of file pfespace.hpp.
|
inline |
Return a reference to the internal GroupCommunicator (on VDofs)
Definition at line 344 of file pfespace.hpp.
|
inline |
Return a const reference to the internal GroupCommunicator (on VDofs)
Definition at line 347 of file pfespace.hpp.
void mfem::ParFiniteElementSpace::Lose_Dof_TrueDof_Matrix | ( | ) |
Definition at line 1545 of file pfespace.cpp.
|
inline |
Definition at line 406 of file pfespace.hpp.
|
inline |
Definition at line 407 of file pfespace.hpp.
|
inline |
Create and return a new HypreParVector on the true dofs, which is owned by (i.e. it must be destroyed by) the calling function.
Definition at line 337 of file pfespace.hpp.
|
inline |
Definition at line 410 of file pfespace.hpp.
void mfem::ParFiniteElementSpace::PrintPartitionStats | ( | ) |
Definition at line 193 of file pfespace.cpp.
GroupCommunicator * mfem::ParFiniteElementSpace::ScalarGroupComm | ( | ) |
Return a new GroupCommunicator on scalar dofs, i.e. for VDim = 1.
Definition at line 994 of file pfespace.cpp.
|
inline |
Definition at line 412 of file pfespace.hpp.
void mfem::ParFiniteElementSpace::Synchronize | ( | Array< int > & | ldof_marker | ) | const |
Given an integer array on the local degrees of freedom, perform a bitwise OR between the shared dofs.
For non-conforming mesh, synchronization is performed on the cut (aka "partially conforming") space.
Definition at line 1008 of file pfespace.cpp.
|
inline |
Obsolete, kept for backward compatibility.
Definition at line 436 of file pfespace.hpp.
|
overridevirtual |
Reflect changes in the mesh. Calculate one of the refinement/derefinement /rebalance matrices, unless want_transform is false.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 3414 of file pfespace.cpp.
|
inlineoverridevirtual |
Free ParGridFunction transformation matrix (if any), to save memory.
Reimplemented from mfem::FiniteElementSpace.
Definition at line 425 of file pfespace.hpp.
Table mfem::ParFiniteElementSpace::face_nbr_element_dof |
Definition at line 208 of file pfespace.hpp.
Table mfem::ParFiniteElementSpace::face_nbr_element_fos |
Definition at line 210 of file pfespace.hpp.
Array<HYPRE_BigInt> mfem::ParFiniteElementSpace::face_nbr_glob_dof_map |
Definition at line 214 of file pfespace.hpp.
Table mfem::ParFiniteElementSpace::face_nbr_ldof |
Definition at line 212 of file pfespace.hpp.
int mfem::ParFiniteElementSpace::num_face_nbr_dofs |
Definition at line 206 of file pfespace.hpp.
Table mfem::ParFiniteElementSpace::send_face_nbr_ldof |
Definition at line 216 of file pfespace.hpp.