22#include <unordered_map>
44 static inline int Map(
int ndofs,
int vdim,
int dof,
int vd);
59Ordering::Map<Ordering::byNODES>(
int ndofs,
int vdim,
int dof,
int vd)
61 MFEM_ASSERT(dof < ndofs && -1-dof < ndofs && 0 <= vd && vd < vdim,
"");
62 return (dof >= 0) ? dof+ndofs*vd : dof-ndofs*vd;
66Ordering::Map<Ordering::byVDIM>(
int ndofs,
int vdim,
int dof,
int vd)
68 MFEM_ASSERT(dof < ndofs && -1-dof < ndofs && 0 <= vd && vd < vdim,
"");
69 return (dof >= 0) ? vd+vdim*dof : -1-(vd+vdim*(-1-dof));
88class BilinearFormIntegrator;
90class QuadratureInterpolator;
91class FaceQuadratureInterpolator;
280 mutable std::unique_ptr<SparseMatrix>
cP;
282 mutable std::unique_ptr<SparseMatrix>
cR;
284 mutable std::unique_ptr<SparseMatrix>
cR_hp;
295 using key_face = std::tuple<bool, ElementDofOrdering, FaceType, L2FaceValues>;
300 return std::get<0>(k)
301 + 2 * (int)std::get<1>(k)
302 + 4 * (int)std::get<2>(k)
303 + 8 * (int)std::get<3>(k);
306 using map_L2F = std::unordered_map<const key_face,FaceRestriction*,key_hash>;
365 int FindDofs(
const Table &var_dof_table,
int row,
int ndof)
const;
385 int variant = 0)
const;
389 int variant = 0)
const;
426 void ConstructDoFTransArray();
433 Table *old_elem_fos,
int old_ndofs);
446 Table *coarse_elem_dof;
448 Table coarse_to_fine;
468 const Table &coarse_elem_dof,
469 const Table *coarse_elem_fos,
481 const Table* old_elem_fos);
485 const Table* old_elem_fos);
1014 {
return (idx >= 0) ? (entity_base + idx) : (-1-(entity_base + (-1-idx))); }
1018 {
return (dof >= 0) ? dof : (-1 - dof); }
1022 {
return (dof >= 0) ? (sign = 1, dof) : (sign = -1, (-1 - dof)); }
1194 int component = -1)
const;
1202 int component = -1)
const;
1282 virtual void Update(
bool want_transform =
true);
1335 void Save(std::ostream &
out)
const;
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
int Size() const
Return the logical size of the array.
Data type dense matrix using column-major storage.
Rank 3 tensor (array of matrices)
Abstract base class that defines an interface for element restrictions.
A class that performs interpolation from a face E-vector to quadrature point values and/or derivative...
Base class for operators that extracts Face degrees of freedom.
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
int GetOrder() const
Return the order (polynomial degree) of the FE collection, corresponding to the order/degree returned...
Derefinement operator, used by the friend class InterpolationGridTransfer.
virtual ~DerefinementOperator()
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
DerefinementOperator(const FiniteElementSpace *f_fes, const FiniteElementSpace *c_fes, BilinearFormIntegrator *mass_integ)
TODO: Implement DofTransformation support.
GridFunction interpolation operator applicable after mesh refinement.
virtual ~RefinementOperator()
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
RefinementOperator(const FiniteElementSpace *fespace, Table *old_elem_dof, Table *old_elem_fos, int old_ndofs)
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
void Save(std::ostream &out) const
Save finite element space to output stream out.
void DestroyDoFTransArray()
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.).
void GetVDofs(int vd, Array< int > &dofs, int ndofs=-1) const
Returns the indices of all of the VDofs for the specified dimension 'vd'.
DofTransformation DoFTrans
int GetNVDofs() const
Number of all scalar vertex dofs.
static int EncodeDof(int entity_base, int idx)
Helper to encode a sign flip into a DOF index (for Hcurl/Hdiv shapes).
const SparseMatrix * GetConformingRestriction() const
The returned SparseMatrix is owned by the FiniteElementSpace.
void ReorderElementToDofTable()
Reorder the scalar DOFs based on the element ordering.
Array< char > var_face_orders
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
void BuildNURBSFaceToDofTable() const
Generates partial face_dof table for a NURBS space.
static void AddDependencies(SparseMatrix &deps, Array< int > &master_dofs, Array< int > &slave_dofs, DenseMatrix &I, int skipfirst=0)
const Table & GetElementToDofTable() const
Return a reference to the internal Table that stores the lists of scalar dofs, for each mesh element,...
Array< int > dof_ldof_array
int GetElementType(int i) const
Returns the type of element i.
Array< StatelessDofTransformation * > DoFTransArray
void DofsToVDofs(Array< int > &dofs, int ndofs=-1) const
Compute the full set of vdofs corresponding to each entry in dofs.
void GetEdgeInteriorDofs(int i, Array< int > &dofs) const
Returns the indices of the degrees of freedom for the interior of the specified edge.
Array< FaceQuadratureInterpolator * > E2BFQ_array
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th boundary fac...
Array< FaceQuadratureInterpolator * > E2IFQ_array
NURBSExtension * GetNURBSext()
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
DofTransformation * GetElementDofs(int elem, Array< int > &dofs) const
Returns indices of degrees of freedom of element 'elem'. The returned indices are offsets into an ldo...
int GetNumElementInteriorDofs(int i) const
Returns the number of degrees of freedom associated with the interior of the specified element.
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 GetEdgeVDofs(int i, Array< int > &vdofs) const
Returns the indices of the degrees of freedom for the specified edge, including the DOFs for the vert...
NURBSExtension * NURBSext
void GetUpdateOperator(OperatorHandle &T)
Return the update operator in the given OperatorHandle, T.
virtual const SparseMatrix * GetRestrictionMatrix() const
The returned SparseMatrix is owned by the FiniteElementSpace.
virtual int GetFaceDofs(int face, Array< int > &dofs, int variant=0) const
Returns the indices of the degrees of freedom for the specified face, including the DOFs for the edge...
virtual void GetTrueTransferOperator(const FiniteElementSpace &coarse_fes, OperatorHandle &T) const
Construct and return an Operator that can be used to transfer true-dof data from coarse_fes,...
static void AdjustVDofs(Array< int > &vdofs)
Remove the orientation information encoded into an array of dofs Some basis function types have a rel...
const Table & GetFaceToDofTable() const
Return a reference to the internal Table that stores the lists of scalar dofs, for each face in the m...
int GetEdgeOrder(int edge, int variant=0) const
bool Nonconforming() const
void GetVertexVDofs(int i, Array< int > &vdofs) const
Returns the indices of the degrees of freedom for the specified vertices.
void GetVertexDofs(int i, Array< int > &dofs) const
Returns the indices of the degrees of freedom for the specified vertices.
void BuildFaceToDofTable() const
int GetFaceOrder(int face, int variant=0) const
Returns the polynomial degree of the i'th face finite element.
int GetNEDofs() const
Number of all scalar edge-interior dofs.
int GetNumBorderDofs(Geometry::Type geom, int order) const
FiniteElementSpace()
Default constructor: the object is invalid until initialized using the method Load().
Array< int > dof_elem_array
int FindFaceDof(int face, int ndof) const
Similar to FindEdgeDof, but used for mixed meshes too.
static int MinOrder(VarOrderBits bits)
Return the minimum order (least significant bit set) in the bit mask.
int GetAttribute(int i) const
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...
MFEM_DEPRECATED void RebuildElementToDofTable()
(
int GetNDofs() const
Returns number of degrees of freedom. This is the number of Local Degrees of Freedom.
bool orders_changed
True if at least one element order changed (variable-order space only).
SparseMatrix * DerefinementMatrix(int old_ndofs, const Table *old_elem_dof, const Table *old_elem_fos)
Calculate GridFunction restriction matrix after mesh derefinement.
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Returns the transformation defining the i-th element in the user-defined variable ElTr.
void GetLocalRefinementMatrices(Geometry::Type geom, DenseTensor &localP) const
void GetEdgeInteriorVDofs(int i, Array< int > &vdofs) const
Returns the indices of the degrees of freedom for the interior of the specified edge.
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,...
SparseMatrix * D2C_GlobalRestrictionMatrix(FiniteElementSpace *cfes)
Generate the global restriction matrix from a discontinuous FE space to the continuous FE space of th...
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
int GetNBE() const
Returns number of boundary elements in the mesh.
virtual void UpdateMeshPointer(Mesh *new_mesh)
const NURBSExtension * GetNURBSext() const
int GetDegenerateFaceDofs(int index, Array< int > &dofs, Geometry::Type master_geom, int variant) const
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
const QuadratureInterpolator * GetQuadratureInterpolator(const IntegrationRule &ir) const
Return a QuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivat...
int GetBdrAttribute(int i) const
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.).
virtual const Operator * GetRestrictionOperator() const
An abstract operator that performs the same action as GetRestrictionMatrix.
int GetLocalDofForDof(int i) const
Return the local dof index in the first element that contains dof i.
static constexpr int MaxVarOrder
int GetElementForDof(int i) const
Return the index of the first element that contains dof i.
virtual void CopyProlongationAndRestriction(const FiniteElementSpace &fes, const Array< int > *perm)
Copies the prolongation and restriction matrices from fes.
Array< char > var_edge_orders
std::unique_ptr< Operator > R_transpose
Operator computing the action of the transpose of the restriction.
void UpdateElementOrders()
Resize the elem_order array on mesh change.
void MakeVDimMatrix(SparseMatrix &mat) const
Replicate 'mat' in the vector dimension, according to vdim ordering mode.
int GetNF() const
Returns number of faces (i.e. co-dimension 1 entities) in the mesh.
SparseMatrix * H2L_GlobalRestrictionMatrix(FiniteElementSpace *lfes)
Construct the restriction matrix from the FE space given by (*this) to the lower degree FE space give...
SparseMatrix * D2Const_GlobalRestrictionMatrix(FiniteElementSpace *cfes)
Generate the global restriction matrix from a discontinuous FE space to the piecewise constant FE spa...
const FiniteElementCollection * fec
Associated FE collection (not owned).
FiniteElementCollection * Load(Mesh *m, std::istream &input)
Read a FiniteElementSpace from a stream. The returned FiniteElementCollection is owned by the caller.
int VDofToDof(int vdof) const
Compute the inverse of the Dof to VDof mapping for a single index vdof.
void BuildBdrElementToDofTable() const
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 ...
Array< QuadratureInterpolator * > E2Q_array
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
Ordering::Type GetOrdering() const
Return the ordering method.
NURBSExtension * StealNURBSext()
ElementTransformation * GetBdrElementTransformation(int i) const
Returns ElementTransformation for the i-th boundary element.
void BuildDofToArrays()
Initialize internal data that enables the use of the methods GetElementForDof() and GetLocalDofForDof...
void GetPatchDofs(int patch, Array< int > &dofs) const
Returns indices of degrees of freedom for NURBS patch index patch. Cartesian ordering is used,...
FiniteElementSpace & operator=(const FiniteElementSpace &)=delete
Copy assignment not supported.
const Operator * GetRestrictionTransposeOperator() const
Return an operator that performs the transpose of GetRestrictionOperator.
int FirstFaceDof(int face, int variant=0) const
static bool DofFinalizable(int dof, const Array< bool > &finalized, const SparseMatrix &deps)
void GetLocalDerefinementMatrices(Geometry::Type geom, DenseTensor &localR) const
void BuildConformingInterpolation() const
Calculate the cP and cR matrices for a nonconforming mesh.
int GetNE() const
Returns number of elements in the mesh.
int vdim
Vector dimension (number of unknowns per degree of freedom).
Table var_face_dofs
NOTE: also used for spaces with mixed faces.
int GetOrder(int i) const
Returns the polynomial degree of the i'th finite element.
OperatorHandle L2E_nat
The element restriction operators, see GetElementRestriction().
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.
int GetConformingVSize() const
virtual void UpdatesFinished()
Free the GridFunction update operator (if any), to save memory.
int * bdofs
internal DOFs of elements if mixed/var-order; NULL otherwise
const ElementRestrictionOperator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
const FiniteElement * GetEdgeElement(int i, int variant=0) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th edge in the ...
void SetUpdateOperatorOwner(bool own)
Set the ownership of the update operator: if set to false, the Operator returned by GetUpdateOperator...
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 vert...
std::tuple< bool, ElementDofOrdering, FaceType, L2FaceValues > key_face
The face restriction operators, see GetFaceRestriction().
std::unique_ptr< SparseMatrix > cR
Conforming restriction matrix such that cR.cP=I.
void ConstructDoFTransArray()
void GetElementInteriorVDofs(int i, Array< int > &vdofs) const
Returns the indices of the degrees of freedom for the interior of the specified element.
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 edge...
virtual void Update(bool want_transform=true)
Reflect changes in the mesh: update number of DOFs, etc. Also, calculate GridFunction transformation ...
OperatorHandle Th
Transformation to apply to GridFunctions after space Update().
void CalcEdgeFaceVarOrders(Array< VarOrderBits > &edge_orders, Array< VarOrderBits > &face_orders) const
SparseMatrix * RefinementMatrix(int old_ndofs, const Table *old_elem_dof, const Table *old_elem_fos)
std::uint64_t VarOrderBits
Bit-mask representing a set of orders needed by an edge/face.
void GetElementVertices(int i, Array< int > &vertices) const
Returns the vertices of element i.
SparseMatrix * RefinementMatrix_main(const int coarse_ndofs, const Table &coarse_elem_dof, const Table *coarse_elem_fos, const DenseTensor localP[]) const
int ndofs
Number of degrees of freedom. Number of unknowns is ndofs * vdim.
std::unique_ptr< SparseMatrix > cP
void SetUpdateOperatorType(Operator::Type tid)
Specify the Operator::Type to be used by the update operators.
void AddEdgeFaceDependencies(SparseMatrix &deps, Array< int > &master_dofs, const FiniteElement *master_fe, Array< int > &slave_dofs, int slave_face, const DenseMatrix *pm) const
const FiniteElement * GetTraceElement(int i, Geometry::Type geom_type) const
Return the trace element from element 'i' to the given 'geom_type'.
void SetRelaxedHpConformity(bool relaxed=true)
FiniteElementSpace(Mesh *mesh, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES)
int MakeDofTable(int ent_dim, const Array< int > &entity_orders, Table &entity_dofs, Array< char > *var_ent_order)
static void MarkerToList(const Array< int > &marker, Array< int > &list)
Convert a Boolean marker array to a list containing all marked indices.
int GetNFDofs() const
Number of all scalar face-interior dofs.
int GetElementOrder(int i) const
Returns the order of the i'th finite element.
Mesh * mesh
The mesh that FE space lives on (not owned).
const FiniteElement * GetFaceElement(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th face in the ...
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 partial...
void SetElementOrder(int i, int p)
Sets the order of the i'th finite element.
const SparseMatrix * GetHpConformingRestriction() const
The returned SparseMatrix is owned by the FiniteElementSpace.
int GetNFbyType(FaceType type) const
Returns the number of faces according to the requested type.
const SparseMatrix * GetConformingProlongation() const
The returned SparseMatrix is owned by the FiniteElementSpace.
const FiniteElementCollection * FEColl() const
Mesh * GetMesh() const
Returns the mesh.
int GetNV() const
Returns number of vertices in the mesh.
int GetElementOrderImpl(int i) const
Return element order: internal version of GetElementOrder without checks.
int GetNVariants(int entity, int index) const
Return number of possible DOF variants for edge/face (var. order spaces).
virtual ~FiniteElementSpace()
int GetMaxElementOrder() const
Return the maximum polynomial order.
const Table & GetBdrElementToDofTable() const
Return a reference to the internal Table that stores the lists of scalar dofs, for each boundary mesh...
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
const Table * GetElementToFaceOrientationTable() const
DofTransformation * GetBdrElementDofs(int bel, Array< int > &dofs) const
Returns indices of degrees of freedom for boundary element 'bel'. The returned indices are offsets in...
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' para...
int GetVDim() const
Returns vector dimension.
const FaceQuadratureInterpolator * GetFaceQuadratureInterpolator(const IntegrationRule &ir, FaceType type) const
Return a FaceQuadratureInterpolator that interpolates E-vectors to quadrature point values and/or der...
bool IsDGSpace() const
Return whether or not the space is discontinuous (L2)
void GetFaceInteriorDofs(int i, Array< int > &dofs) const
Returns the indices of the degrees of freedom for the interior of the specified face.
virtual const SparseMatrix * GetHpRestrictionMatrix() const
The returned SparseMatrix is owned by the FiniteElementSpace.
int GetNConformingDofs() const
std::unordered_map< const key_face, FaceRestriction *, key_hash > map_L2F
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.
void GetElementInteriorDofs(int i, Array< int > &dofs) const
Returns the indices of the degrees of freedom for the interior of the specified element.
void Constructor(Mesh *mesh, NURBSExtension *ext, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES)
Help function for constructors + Load().
DofTransformation * GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for i'th boundary element. The returned indices are offsets int...
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 conform...
virtual void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_vdofs, int component=-1) const
Mark degrees of freedom associated with boundary elements with the specified boundary attributes (mar...
int DofToVDof(int dof, int vd, int ndofs=-1) const
Compute a single vdof corresponding to the index dof and the vector index vd.
virtual const FaceRestriction * GetFaceRestriction(ElementDofOrdering f_ordering, FaceType, L2FaceValues mul=L2FaceValues::DoubleValued) const
Return an Operator that converts L-vectors to E-vectors on each face.
const Operator * GetUpdateOperator()
Get the GridFunction update operator.
int FindEdgeDof(int edge, int ndof) const
void GetPatchVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom in vdofs for NURBS patch i.
void BuildElementToDofTable() const
std::unique_ptr< SparseMatrix > cR_hp
A version of the conforming restriction matrix for variable-order spaces.
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.
Abstract class for all finite elements.
Class for an integration rule - an Array of IntegrationPoint.
Transfer data between a coarse mesh and an embedded refined mesh using interpolation.
Arbitrary order "L2-conforming" discontinuous finite elements.
Abstract base class for LORDiscretization and ParLORDiscretization classes, which construct low-order...
Element::Type GetElementType(int i) const
Returns the type of element i.
Element::Type GetBdrElementType(int i) const
Returns the type of boundary element i.
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
int GetAttribute(int i) const
Return the attribute of element i.
void GetElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of element i.
bool Nonconforming() const
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
virtual int GetNFbyType(FaceType type) const
Returns the number of faces according to the requested type, does not count master nonconforming face...
int GetNE() const
Returns number of elements.
int Dimension() const
Dimension of the reference space used within the elements.
void GetElementTransformation(int i, IsoparametricTransformation *ElTr) const
Builds the transformation defining the i-th element in ElTr. ElTr must be allocated in advance and wi...
ElementTransformation * GetBdrElementTransformation(int i)
Returns a pointer to the transformation defining the i-th boundary element.
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
void Swap(Mesh &other, bool non_geometry)
int GetNBE() const
Returns number of boundary elements.
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
Pointer to an Operator of a specified type.
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
void SetType(Operator::Type tid)
Invoke Clear() and set a new type id.
Operator * Ptr() const
Access the underlying Operator pointer.
void Clear()
Clear the OperatorHandle, deleting the held Operator (if owned), while leaving the type id unchanged.
Type
Enumeration defining IDs for some classes derived from Operator.
The ordering method used when the number of unknowns per mesh node (vector dimension) is bigger than ...
static int Map(int ndofs, int vdim, int dof, int vd)
static void DofsToVDofs(int ndofs, int vdim, Array< int > &dofs)
Matrix-free transfer operator between finite element spaces on the same mesh.
A class that performs interpolation from an E-vector to quadrature point values and/or derivatives (Q...
Class representing the storage layout of a QuadratureFunction.
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
int index(int i, int j, int nx, int ny)
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
bool UsesTensorBasis(const FiniteElementSpace &fes)
Return true if the mesh contains only one topology and the elements are tensor elements.
QVectorLayout
Type describing possible layouts for Q-vectors.
@ byNODES
NQPT x VDIM x NE (values) / NQPT x VDIM x DIM x NE (grads)
@ byVDIM
VDIM x NQPT x NE (values) / VDIM x DIM x NQPT x NE (grads)
ElementDofOrdering GetEVectorOrdering(const FiniteElementSpace &fes)
Return LEXICOGRAPHIC if mesh contains only one topology and the elements are tensor elements,...
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
@ LEXICOGRAPHIC
Lexicographic ordering for tensor-product FiniteElements.
@ NATIVE
Native ordering as defined by the FiniteElement.
real_t p(const Vector &x, real_t t)
std::size_t operator()(const key_face &k) const