MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
fespace.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2021, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_FESPACE
13 #define MFEM_FESPACE
14 
15 #include "../config/config.hpp"
16 #include "../linalg/sparsemat.hpp"
17 #include "../mesh/mesh.hpp"
18 #include "fe_coll.hpp"
19 #include "restriction.hpp"
20 #include <iostream>
21 #include <unordered_map>
22 
23 namespace mfem
24 {
25 
26 /** @brief The ordering method used when the number of unknowns per mesh node
27  (vector dimension) is bigger than 1. */
28 class Ordering
29 {
30 public:
31  /// %Ordering methods:
32  enum Type
33  {
34  byNODES, /**< loop first over the nodes (inner loop) then over the vector
35  dimension (outer loop); symbolically it can be represented
36  as: XXX...,YYY...,ZZZ... */
37  byVDIM /**< loop first over the vector dimension (inner loop) then over
38  the nodes (outer loop); symbolically it can be represented
39  as: XYZ,XYZ,XYZ,... */
40  };
41 
42  template <Type Ord>
43  static inline int Map(int ndofs, int vdim, int dof, int vd);
44 
45  template <Type Ord>
46  static void DofsToVDofs(int ndofs, int vdim, Array<int> &dofs);
47 };
48 
49 template <> inline int
50 Ordering::Map<Ordering::byNODES>(int ndofs, int vdim, int dof, int vd)
51 {
52  MFEM_ASSERT(dof < ndofs && -1-dof < ndofs && 0 <= vd && vd < vdim, "");
53  return (dof >= 0) ? dof+ndofs*vd : dof-ndofs*vd;
54 }
55 
56 template <> inline int
57 Ordering::Map<Ordering::byVDIM>(int ndofs, int vdim, int dof, int vd)
58 {
59  MFEM_ASSERT(dof < ndofs && -1-dof < ndofs && 0 <= vd && vd < vdim, "");
60  return (dof >= 0) ? vd+vdim*dof : -1-(vd+vdim*(-1-dof));
61 }
62 
63 
64 /// Constants describing the possible orderings of the DOFs in one element.
66 {
67  /// Native ordering as defined by the FiniteElement.
68  /** This ordering can be used by tensor-product elements when the
69  interpolation from the DOFs to quadrature points does not use the
70  tensor-product structure. */
71  NATIVE,
72  /// Lexicographic ordering for tensor-product FiniteElements.
73  /** This ordering can be used only with tensor-product elements. */
75 };
76 
77 // Forward declarations
78 class NURBSExtension;
79 class BilinearFormIntegrator;
80 class QuadratureSpace;
81 class QuadratureInterpolator;
82 class FaceQuadratureInterpolator;
83 
84 
85 /** @brief Class FiniteElementSpace - responsible for providing FEM view of the
86  mesh, mainly managing the set of degrees of freedom. */
88 {
91  friend void Mesh::Swap(Mesh &, bool);
92  friend class LORBase;
93 
94 protected:
95  /// The mesh that FE space lives on (not owned).
97 
98  /// Associated FE collection (not owned).
100 
101  /// %Vector dimension (number of unknowns per degree of freedom).
102  int vdim;
103 
104  /** Type of ordering of the vector dofs when #vdim > 1.
105  - Ordering::byNODES - first nodes, then vector dimension,
106  - Ordering::byVDIM - first vector dimension, then nodes */
108 
109  /// Number of degrees of freedom. Number of unknowns is #ndofs * #vdim.
110  int ndofs;
111 
112  /** Polynomial order for each element. If empty, all elements are assumed
113  to be of the default order (fec->GetOrder()). */
115 
117  int uni_fdof; ///< # of single face DOFs if all faces uniform; -1 otherwise
118  int *bdofs; ///< internal DOFs of elements if mixed/var-order; NULL otherwise
119 
120  /** Variable order spaces only: DOF assignments for edges and faces, see
121  docs in MakeDofTable. For constant order spaces the tables are empty. */
123  Table var_face_dofs; ///< NOTE: also used for spaces with mixed faces
124 
125  /** Additional data for the var_*_dofs tables: individual variant orders
126  (these are basically alternate J arrays for var_edge/face_dofs). */
128 
129  // precalculated DOFs for each element, boundary element, and face
130  mutable Table *elem_dof; // owned (except in NURBS FE space)
131  mutable Table *bdr_elem_dof; // owned (except in NURBS FE space)
132  mutable Table *face_dof; // owned; in var-order space contains variant 0 DOFs
133 
135 
137  int own_ext;
138  mutable Array<int> face_to_be; // NURBS FE space only
139 
140  /** Matrix representing the prolongation from the global conforming dofs to
141  a set of intermediate partially conforming dofs, e.g. the dofs associated
142  with a "cut" space on a non-conforming mesh. */
143  mutable SparseMatrix *cP; // owned
144  /// Conforming restriction matrix such that cR.cP=I.
145  mutable SparseMatrix *cR; // owned
146  /// A version of the conforming restriction matrix for variable-order spaces.
147  mutable SparseMatrix *cR_hp; // owned
148  mutable bool cP_is_set;
149 
150  /// Transformation to apply to GridFunctions after space Update().
152 
153  /// The element restriction operators, see GetElementRestriction().
155  /// The face restriction operators, see GetFaceRestriction().
156  using key_face = std::tuple<bool, ElementDofOrdering, FaceType, L2FaceValues>;
157  struct key_hash
158  {
159  std::size_t operator()(const key_face& k) const
160  {
161  return std::get<0>(k)
162  + 2 * (int)std::get<1>(k)
163  + 4 * (int)std::get<2>(k)
164  + 8 * (int)std::get<3>(k);
165  }
166  };
167  using map_L2F = std::unordered_map<const key_face,FaceRestriction*,key_hash>;
168  mutable map_L2F L2F;
169 
173 
174  /** Update counter, incremented every time the space is constructed/updated.
175  Used by GridFunctions to check if they are up to date with the space. */
176  long sequence;
177 
178  /** Mesh sequence number last seen when constructing the space. The space
179  needs updating if Mesh::GetSequence() is larger than this. */
181 
182  /// True if at least one element order changed (variable-order space only).
184 
185  bool relaxed_hp; // see SetRelaxedHpConformity()
186 
187  void UpdateNURBS();
188 
189  void Construct();
190  void Destroy();
191 
192  void BuildElementToDofTable() const;
193  void BuildBdrElementToDofTable() const;
194  void BuildFaceToDofTable() const;
195 
196  /** @brief Generates partial face_dof table for a NURBS space.
197 
198  The table is only defined for exterior faces that coincide with a
199  boundary. */
200  void BuildNURBSFaceToDofTable() const;
201 
202  /// Bit-mask representing a set of orders needed by an edge/face.
203  typedef std::uint64_t VarOrderBits;
204  static constexpr int MaxVarOrder = 8*sizeof(VarOrderBits) - 1;
205 
206  /// Return the minimum order (least significant bit set) in the bit mask.
207  static int MinOrder(VarOrderBits bits);
208 
209  /// Return element order: internal version of GetElementOrder without checks.
210  int GetElementOrderImpl(int i) const;
211 
212  /** In a variable order space, calculate a bitmask of polynomial orders that
213  need to be represented on each edge and face. */
214  void CalcEdgeFaceVarOrders(Array<VarOrderBits> &edge_orders,
215  Array<VarOrderBits> &face_orders) const;
216 
217  /** Build the table var_edge_dofs (or var_face_dofs) in a variable order
218  space; return total edge/face DOFs. */
219  int MakeDofTable(int ent_dim, const Array<int> &entity_orders,
220  Table &entity_dofs, Array<char> *var_ent_order);
221 
222  /// Search row of a DOF table for a DOF set of size 'ndof', return first DOF.
223  int FindDofs(const Table &var_dof_table, int row, int ndof) const;
224 
225  /** In a variable order space, return edge DOFs associated with a polynomial
226  order that has 'ndof' degrees of freedom. */
227  int FindEdgeDof(int edge, int ndof) const
228  { return FindDofs(var_edge_dofs, edge, ndof); }
229 
230  /// Similar to FindEdgeDof, but used for mixed meshes too.
231  int FindFaceDof(int face, int ndof) const
232  { return FindDofs(var_face_dofs, face, ndof); }
233 
234  int FirstFaceDof(int face, int variant = 0) const
235  { return uni_fdof >= 0 ? face*uni_fdof : var_face_dofs.GetRow(face)[variant];}
236 
237  /// Return number of possible DOF variants for edge/face (var. order spaces).
238  int GetNVariants(int entity, int index) const;
239 
240  /// Helper to encode a sign flip into a DOF index (for Hcurl/Hdiv shapes).
241  static inline int EncodeDof(int entity_base, int idx)
242  { return (idx >= 0) ? (entity_base + idx) : (-1-(entity_base + (-1-idx))); }
243 
244  /// Helpers to remove encoded sign from a DOF
245  static inline int DecodeDof(int dof)
246  { return (dof >= 0) ? dof : (-1 - dof); }
247 
248  static inline int DecodeDof(int dof, double& sign)
249  { return (dof >= 0) ? (sign = 1, dof) : (sign = -1, (-1 - dof)); }
250 
251  /// Helper to get vertex, edge or face DOFs (entity=0,1,2 resp.).
252  int GetEntityDofs(int entity, int index, Array<int> &dofs,
253  Geometry::Type master_geom = Geometry::INVALID,
254  int variant = 0) const;
255 
256  // Get degenerate face DOFs: see explanation in method implementation.
257  int GetDegenerateFaceDofs(int index, Array<int> &dofs,
258  Geometry::Type master_geom, int variant) const;
259 
260  int GetNumBorderDofs(Geometry::Type geom, int order) const;
261 
262  /// Calculate the cP and cR matrices for a nonconforming mesh.
263  void BuildConformingInterpolation() const;
264 
265  static void AddDependencies(SparseMatrix& deps, Array<int>& master_dofs,
266  Array<int>& slave_dofs, DenseMatrix& I,
267  int skipfirst = 0);
268 
269  static bool DofFinalizable(int dof, const Array<bool>& finalized,
270  const SparseMatrix& deps);
271 
272  void AddEdgeFaceDependencies(SparseMatrix &deps, Array<int>& master_dofs,
273  const FiniteElement *master_fe,
274  Array<int> &slave_dofs, int slave_face,
275  const DenseMatrix *pm) const;
276 
277  /// Replicate 'mat' in the vector dimension, according to vdim ordering mode.
278  void MakeVDimMatrix(SparseMatrix &mat) const;
279 
280  /// GridFunction interpolation operator applicable after mesh refinement.
282  {
283  const FiniteElementSpace* fespace;
285  Table* old_elem_dof; // Owned.
286 
287  public:
288  /** Construct the operator based on the elem_dof table of the original
289  (coarse) space. The class takes ownership of the table. */
290  RefinementOperator(const FiniteElementSpace* fespace,
291  Table *old_elem_dof/*takes ownership*/, int old_ndofs);
292  RefinementOperator(const FiniteElementSpace *fespace,
293  const FiniteElementSpace *coarse_fes);
294  virtual void Mult(const Vector &x, Vector &y) const;
295  virtual void MultTranspose(const Vector &x, Vector &y) const;
296  virtual ~RefinementOperator();
297  };
298 
299  /// Derefinement operator, used by the friend class InterpolationGridTransfer.
301  {
302  const FiniteElementSpace *fine_fes; // Not owned.
304  Table *coarse_elem_dof; // Owned.
305  Table coarse_to_fine;
306  Array<int> coarse_to_ref_type;
307  Array<Geometry::Type> ref_type_to_geom;
308  Array<int> ref_type_to_fine_elem_offset;
309 
310  public:
312  const FiniteElementSpace *c_fes,
313  BilinearFormIntegrator *mass_integ);
314  virtual void Mult(const Vector &x, Vector &y) const;
315  virtual ~DerefinementOperator();
316  };
317 
318  /** This method makes the same assumptions as the method:
319  void GetLocalRefinementMatrices(
320  const FiniteElementSpace &coarse_fes, Geometry::Type geom,
321  DenseTensor &localP) const
322  which is defined below. It also assumes that the coarse fes and this have
323  the same vector dimension, vdim. */
324  SparseMatrix *RefinementMatrix_main(const int coarse_ndofs,
325  const Table &coarse_elem_dof,
326  const DenseTensor localP[]) const;
327 
329  DenseTensor &localP) const;
331  DenseTensor &localR) const;
332 
333  /** Calculate explicit GridFunction interpolation matrix (after mesh
334  refinement). NOTE: consider using the RefinementOperator class instead
335  of the fully assembled matrix, which can take a lot of memory. */
336  SparseMatrix* RefinementMatrix(int old_ndofs, const Table* old_elem_dof);
337 
338  /// Calculate GridFunction restriction matrix after mesh derefinement.
339  SparseMatrix* DerefinementMatrix(int old_ndofs, const Table* old_elem_dof);
340 
341  /** @brief Return in @a localP the local refinement matrices that map
342  between fespaces after mesh refinement. */
343  /** This method assumes that this->mesh is a refinement of coarse_fes->mesh
344  and that the CoarseFineTransformations of this->mesh are set accordingly.
345  Another assumption is that the FEs of this use the same MapType as the FEs
346  of coarse_fes. Finally, it assumes that the spaces this and coarse_fes are
347  NOT variable-order spaces. */
348  void GetLocalRefinementMatrices(const FiniteElementSpace &coarse_fes,
350  DenseTensor &localP) const;
351 
352  /// Help function for constructors + Load().
353  void Constructor(Mesh *mesh, NURBSExtension *ext,
355  int vdim = 1, int ordering = Ordering::byNODES);
356 
357  /// Updates the internal mesh pointer. @warning @a new_mesh must be
358  /// <b>topologically identical</b> to the existing mesh. Used if the address
359  /// of the Mesh object has changed, e.g. in @a Mesh::Swap.
360  virtual void UpdateMeshPointer(Mesh *new_mesh);
361 
362  /// Resize the elem_order array on mesh change.
363  void UpdateElementOrders();
364 
365  /// @brief Copies the prolongation and restriction matrices from @a fes.
366  ///
367  /// Used for low order preconditioning on non-conforming meshes. If the DOFs
368  /// require a permutation, it will be supplied by non-NULL @a perm. NULL @a
369  /// perm indicates that no permutation is required.
370  virtual void CopyProlongationAndRestriction(const FiniteElementSpace &fes,
371  const Array<int> *perm);
372 
373 public:
374  /** @brief Default constructor: the object is invalid until initialized using
375  the method Load(). */
377 
378  /** @brief Copy constructor: deep copy all data from @a orig except the Mesh,
379  the FiniteElementCollection, ans some derived data. */
380  /** If the @a mesh or @a fec pointers are NULL (default), then the new
381  FiniteElementSpace will reuse the respective pointers from @a orig. If
382  any of these pointers is not NULL, the given pointer will be used instead
383  of the one used by @a orig.
384 
385  @note The objects pointed to by the @a mesh and @a fec parameters must be
386  either the same objects as the ones used by @a orig, or copies of them.
387  Otherwise, the behavior is undefined.
388 
389  @note Derived data objects, such as the conforming prolongation and
390  restriction matrices, and the update operator, will not be copied, even
391  if they are created in the @a orig object. */
392  FiniteElementSpace(const FiniteElementSpace &orig, Mesh *mesh = NULL,
393  const FiniteElementCollection *fec = NULL);
394 
397  int vdim = 1, int ordering = Ordering::byNODES)
398  { Constructor(mesh, NULL, fec, vdim, ordering); }
399 
400  /// Construct a NURBS FE space based on the given NURBSExtension, @a ext.
401  /** @note If the pointer @a ext is NULL, this constructor is equivalent to
402  the standard constructor with the same arguments minus the
403  NURBSExtension, @a ext. */
406  int vdim = 1, int ordering = Ordering::byNODES)
407  { Constructor(mesh, ext, fec, vdim, ordering); }
408 
409  /// Returns the mesh
410  inline Mesh *GetMesh() const { return mesh; }
411 
412  const NURBSExtension *GetNURBSext() const { return NURBSext; }
415 
416  bool Conforming() const { return mesh->Conforming() && cP == NULL; }
417  bool Nonconforming() const { return mesh->Nonconforming() || cP != NULL; }
418 
419  /// Sets the order of the i'th finite element.
420  /** By default, all elements are assumed to be of fec->GetOrder(). Once
421  SetElementOrder is called, the space becomes a variable order space. */
422  void SetElementOrder(int i, int p);
423 
424  /// Returns the order of the i'th finite element.
425  int GetElementOrder(int i) const;
426 
427  /// Return the maximum polynomial order.
428  int GetMaxElementOrder() const
429  { return IsVariableOrder() ? elem_order.Max() : fec->GetOrder(); }
430 
431  /// Returns true if the space contains elements of varying polynomial orders.
432  bool IsVariableOrder() const { return elem_order.Size(); }
433 
434  /// The returned SparseMatrix is owned by the FiniteElementSpace.
436 
437  /// The returned SparseMatrix is owned by the FiniteElementSpace.
438  const SparseMatrix *GetConformingRestriction() const;
439 
440  /** Return a version of the conforming restriction matrix for variable-order
441  spaces with complex hp interfaces, where some true DOFs are not owned by
442  any elements and need to be interpolated from higher order edge/face
443  variants (see also @a SetRelaxedHpConformity()). */
444  /// The returned SparseMatrix is owned by the FiniteElementSpace.
446 
447  /// The returned Operator is owned by the FiniteElementSpace.
448  virtual const Operator *GetProlongationMatrix() const
449  { return GetConformingProlongation(); }
450 
451  /// Return an operator that performs the transpose of GetRestrictionOperator
452  /** The returned operator is owned by the FiniteElementSpace. In serial this
453  is the same as GetProlongationMatrix() */
455  { return GetConformingProlongation(); }
456 
457  /// An abstract operator that performs the same action as GetRestrictionMatrix
458  /** In some cases this is an optimized matrix-free implementation. The
459  returned operator is owned by the FiniteElementSpace. */
460  virtual const Operator *GetRestrictionOperator() const
461  { return GetConformingRestriction(); }
462 
463  /// The returned SparseMatrix is owned by the FiniteElementSpace.
464  virtual const SparseMatrix *GetRestrictionMatrix() const
465  { return GetConformingRestriction(); }
466 
467  /// The returned SparseMatrix is owned by the FiniteElementSpace.
468  virtual const SparseMatrix *GetHpRestrictionMatrix() const
469  { return GetHpConformingRestriction(); }
470 
471  /// Return an Operator that converts L-vectors to E-vectors.
472  /** An L-vector is a vector of size GetVSize() which is the same size as a
473  GridFunction. An E-vector represents the element-wise discontinuous
474  version of the FE space.
475 
476  The layout of the E-vector is: ND x VDIM x NE, where ND is the number of
477  degrees of freedom, VDIM is the vector dimension of the FE space, and NE
478  is the number of the mesh elements.
479 
480  The parameter @a e_ordering describes how the local DOFs in each element
481  should be ordered, see ElementDofOrdering.
482 
483  For discontinuous spaces, the element restriction corresponds to a
484  permutation of the degrees of freedom, implemented by the
485  L2ElementRestriction class.
486 
487  The returned Operator is owned by the FiniteElementSpace. */
488  const Operator *GetElementRestriction(ElementDofOrdering e_ordering) const;
489 
490  /// Return an Operator that converts L-vectors to E-vectors on each face.
491  virtual const FaceRestriction *GetFaceRestriction(
492  ElementDofOrdering e_ordering, FaceType,
494 
495  /** @brief Return a QuadratureInterpolator that interpolates E-vectors to
496  quadrature point values and/or derivatives (Q-vectors). */
497  /** An E-vector represents the element-wise discontinuous version of the FE
498  space and can be obtained, for example, from a GridFunction using the
499  Operator returned by GetElementRestriction().
500 
501  All elements will use the same IntegrationRule, @a ir as the target
502  quadrature points. */
504  const IntegrationRule &ir) const;
505 
506  /** @brief Return a QuadratureInterpolator that interpolates E-vectors to
507  quadrature point values and/or derivatives (Q-vectors). */
508  /** An E-vector represents the element-wise discontinuous version of the FE
509  space and can be obtained, for example, from a GridFunction using the
510  Operator returned by GetElementRestriction().
511 
512  The target quadrature points in the elements are described by the given
513  QuadratureSpace, @a qs. */
515  const QuadratureSpace &qs) const;
516 
517  /** @brief Return a FaceQuadratureInterpolator that interpolates E-vectors to
518  quadrature point values and/or derivatives (Q-vectors). */
520  const IntegrationRule &ir, FaceType type) const;
521 
522  /// Returns the polynomial degree of the i'th finite element.
523  /** NOTE: it is recommended to use GetElementOrder in new code. */
524  int GetOrder(int i) const { return GetElementOrder(i); }
525 
526  /** Return the order of an edge. In a variable order space, return the order
527  of a specific variant, or -1 if there are no more variants. */
528  int GetEdgeOrder(int edge, int variant = 0) const;
529 
530  /// Returns the polynomial degree of the i'th face finite element
531  int GetFaceOrder(int face, int variant = 0) const;
532 
533  /// Returns vector dimension.
534  inline int GetVDim() const { return vdim; }
535 
536  /// Returns number of degrees of freedom.
537  inline int GetNDofs() const { return ndofs; }
538 
539  /// Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
540  inline int GetVSize() const { return vdim * ndofs; }
541 
542  /// Return the number of vector true (conforming) dofs.
543  virtual int GetTrueVSize() const { return GetConformingVSize(); }
544 
545  /// Returns the number of conforming ("true") degrees of freedom
546  /// (if the space is on a nonconforming mesh with hanging nodes).
547  int GetNConformingDofs() const;
548 
549  int GetConformingVSize() const { return vdim * GetNConformingDofs(); }
550 
551  /// Return the ordering method.
552  inline Ordering::Type GetOrdering() const { return ordering; }
553 
554  const FiniteElementCollection *FEColl() const { return fec; }
555 
556  /// Number of all scalar vertex dofs
557  int GetNVDofs() const { return nvdofs; }
558  /// Number of all scalar edge-interior dofs
559  int GetNEDofs() const { return nedofs; }
560  /// Number of all scalar face-interior dofs
561  int GetNFDofs() const { return nfdofs; }
562 
563  /// Returns number of vertices in the mesh.
564  inline int GetNV() const { return mesh->GetNV(); }
565 
566  /// Returns number of elements in the mesh.
567  inline int GetNE() const { return mesh->GetNE(); }
568 
569  /// Returns number of faces (i.e. co-dimension 1 entities) in the mesh.
570  /** The co-dimension 1 entities are those that have dimension 1 less than the
571  mesh dimension, e.g. for a 2D mesh, the faces are the 1D entities, i.e.
572  the edges. */
573  inline int GetNF() const { return mesh->GetNumFaces(); }
574 
575  /// Returns number of boundary elements in the mesh.
576  inline int GetNBE() const { return mesh->GetNBE(); }
577 
578  /// Returns the number of faces according to the requested type.
579  /** If type==Boundary returns only the "true" number of boundary faces
580  contrary to GetNBE() that returns "fake" boundary faces associated to
581  visualization for GLVis.
582  Similarly, if type==Interior, the "fake" boundary faces associated to
583  visualization are counted as interior faces. */
584  inline int GetNFbyType(FaceType type) const
585  { return mesh->GetNFbyType(type); }
586 
587  /// Returns the type of element i.
588  inline int GetElementType(int i) const
589  { return mesh->GetElementType(i); }
590 
591  /// Returns the vertices of element i.
592  inline void GetElementVertices(int i, Array<int> &vertices) const
593  { mesh->GetElementVertices(i, vertices); }
594 
595  /// Returns the type of boundary element i.
596  inline int GetBdrElementType(int i) const
597  { return mesh->GetBdrElementType(i); }
598 
599  /// Returns ElementTransformation for the @a i-th element.
601  { return mesh->GetElementTransformation(i); }
602 
603  /** @brief Returns the transformation defining the @a i-th element in the
604  user-defined variable @a ElTr. */
606  { mesh->GetElementTransformation(i, ElTr); }
607 
608  /// Returns ElementTransformation for the @a i-th boundary element.
610  { return mesh->GetBdrElementTransformation(i); }
611 
612  int GetAttribute(int i) const { return mesh->GetAttribute(i); }
613 
614  int GetBdrAttribute(int i) const { return mesh->GetBdrAttribute(i); }
615 
616  /// Returns indices of degrees of freedom of element 'elem'.
617  virtual void GetElementDofs(int elem, Array<int> &dofs) const;
618 
619  /// Returns indices of degrees of freedom for boundary element 'bel'.
620  virtual void GetBdrElementDofs(int bel, Array<int> &dofs) const;
621 
622  /** @brief Returns the indices of the degrees of freedom for the specified
623  face, including the DOFs for the edges and the vertices of the face. */
624  /** In variable order spaces, multiple variants of DOFs can be returned.
625  See @a GetEdgeDofs for more details.
626  @return Order of the selected variant, or -1 if there are no more
627  variants.*/
628  virtual int GetFaceDofs(int face, Array<int> &dofs, int variant = 0) const;
629 
630  /** @brief Returns the indices of the degrees of freedom for the specified
631  edge, including the DOFs for the vertices of the edge. */
632  /** In variable order spaces, multiple sets of DOFs may exist on an edge,
633  corresponding to the different polynomial orders of incident elements.
634  The 'variant' parameter is the zero-based index of the desired DOF set.
635  The variants are ordered from lowest polynomial degree to the highest.
636  @return Order of the selected variant, or -1 if there are no more
637  variants. */
638  int GetEdgeDofs(int edge, Array<int> &dofs, int variant = 0) const;
639 
640  void GetVertexDofs(int i, Array<int> &dofs) const;
641 
642  void GetElementInteriorDofs(int i, Array<int> &dofs) const;
643 
644  void GetFaceInteriorDofs(int i, Array<int> &dofs) const;
645 
646  int GetNumElementInteriorDofs(int i) const;
647 
648  void GetEdgeInteriorDofs(int i, Array<int> &dofs) const;
649 
650  /** @brief Returns the indices of all of the VDofs for the specified
651  dimension 'vd'. */
652  /** The 'ndofs' parameter defines the number of Dofs in the
653  FiniteElementSpace. If 'ndofs' is -1 (the default value), then the
654  number of Dofs is determined by the FiniteElementSpace. */
655  void GetVDofs(int vd, Array<int> &dofs, int ndofs = -1) const;
656 
657  void DofsToVDofs(Array<int> &dofs, int ndofs = -1) const;
658 
659  void DofsToVDofs(int vd, Array<int> &dofs, int ndofs = -1) const;
660 
661  int DofToVDof(int dof, int vd, int ndofs = -1) const;
662 
663  int VDofToDof(int vdof) const
664  { return (ordering == Ordering::byNODES) ? (vdof%ndofs) : (vdof/vdim); }
665 
666  static void AdjustVDofs(Array<int> &vdofs);
667 
668  /// Returns indexes of degrees of freedom in array dofs for i'th element.
669  void GetElementVDofs(int i, Array<int> &vdofs) const;
670 
671  /// Returns indexes of degrees of freedom for i'th boundary element.
672  void GetBdrElementVDofs(int i, Array<int> &vdofs) const;
673 
674  /// Returns indexes of degrees of freedom for i'th face element (2D and 3D).
675  void GetFaceVDofs(int i, Array<int> &vdofs) const;
676 
677  /// Returns indexes of degrees of freedom for i'th edge.
678  void GetEdgeVDofs(int i, Array<int> &vdofs) const;
679 
680  void GetVertexVDofs(int i, Array<int> &vdofs) const;
681 
682  void GetElementInteriorVDofs(int i, Array<int> &vdofs) const;
683 
684  void GetEdgeInteriorVDofs(int i, Array<int> &vdofs) const;
685 
686  /// (@deprecated) Use the Update() method if the space or mesh changed.
687  MFEM_DEPRECATED void RebuildElementToDofTable();
688 
689  /** @brief Reorder the scalar DOFs based on the element ordering.
690 
691  The new ordering is constructed as follows: 1) loop over all elements as
692  ordered in the Mesh; 2) for each element, assign new indices to all of
693  its current DOFs that are still unassigned; the new indices we assign are
694  simply the sequence `0,1,2,...`; if there are any signed DOFs their sign
695  is preserved. */
697 
698  /** @brief Return a reference to the internal Table that stores the lists of
699  scalar dofs, for each mesh element, as returned by GetElementDofs(). */
700  const Table &GetElementToDofTable() const { return *elem_dof; }
701 
702  /** @brief Return a reference to the internal Table that stores the lists of
703  scalar dofs, for each boundary mesh element, as returned by
704  GetBdrElementDofs(). */
706  { if (!bdr_elem_dof) { BuildBdrElementToDofTable(); } return *bdr_elem_dof; }
707 
708  /** @brief Return a reference to the internal Table that stores the lists of
709  scalar dofs, for each face in the mesh, as returned by GetFaceDofs(). In
710  this context, "face" refers to a (dim-1)-dimensional mesh entity. */
711  /** @note In the case of a NURBS space, the rows corresponding to interior
712  faces will be empty. */
713  const Table &GetFaceToDofTable() const
714  { if (!face_dof) { BuildFaceToDofTable(); } return *face_dof; }
715 
716  /** @brief Initialize internal data that enables the use of the methods
717  GetElementForDof() and GetLocalDofForDof(). */
718  void BuildDofToArrays();
719 
720  /// Return the index of the first element that contains dof @a i.
721  /** This method can be called only after setup is performed using the method
722  BuildDofToArrays(). */
723  int GetElementForDof(int i) const { return dof_elem_array[i]; }
724  /// Return the local dof index in the first element that contains dof @a i.
725  /** This method can be called only after setup is performed using the method
726  BuildDofToArrays(). */
727  int GetLocalDofForDof(int i) const { return dof_ldof_array[i]; }
728 
729  /** @brief Returns pointer to the FiniteElement in the FiniteElementCollection
730  associated with i'th element in the mesh object. */
731  virtual const FiniteElement *GetFE(int i) const;
732 
733  /** @brief Returns pointer to the FiniteElement in the FiniteElementCollection
734  associated with i'th boundary face in the mesh object. */
735  const FiniteElement *GetBE(int i) const;
736 
737  /** @brief Returns pointer to the FiniteElement in the FiniteElementCollection
738  associated with i'th face in the mesh object. Faces in this case refer
739  to the MESHDIM-1 primitive so in 2D they are segments and in 1D they are
740  points.*/
741  const FiniteElement *GetFaceElement(int i) const;
742 
743  /** @brief Returns pointer to the FiniteElement in the FiniteElementCollection
744  associated with i'th edge in the mesh object. */
745  const FiniteElement *GetEdgeElement(int i, int variant = 0) const;
746 
747  /// Return the trace element from element 'i' to the given 'geom_type'
748  const FiniteElement *GetTraceElement(int i, Geometry::Type geom_type) const;
749 
750  /** @brief Mark degrees of freedom associated with boundary elements with
751  the specified boundary attributes (marked in 'bdr_attr_is_ess').
752  For spaces with 'vdim' > 1, the 'component' parameter can be used
753  to restricts the marked vDOFs to the specified component. */
754  virtual void GetEssentialVDofs(const Array<int> &bdr_attr_is_ess,
755  Array<int> &ess_vdofs,
756  int component = -1) const;
757 
758  /** @brief Get a list of essential true dofs, ess_tdof_list, corresponding to the
759  boundary attributes marked in the array bdr_attr_is_ess.
760  For spaces with 'vdim' > 1, the 'component' parameter can be used
761  to restricts the marked tDOFs to the specified component. */
762  virtual void GetEssentialTrueDofs(const Array<int> &bdr_attr_is_ess,
763  Array<int> &ess_tdof_list,
764  int component = -1);
765 
766  /** @brief Get a list of all boundary true dofs, @a boundary_dofs. For spaces
767  with 'vdim' > 1, the 'component' parameter can be used to restricts the
768  marked tDOFs to the specified component. Equivalent to
769  FiniteElementSpace::GetEssentialTrueDofs with all boundary attributes
770  marked as essential. */
771  void GetBoundaryTrueDofs(Array<int> &boundary_dofs, int component = -1);
772 
773  /// Convert a Boolean marker array to a list containing all marked indices.
774  static void MarkerToList(const Array<int> &marker, Array<int> &list);
775 
776  /** @brief Convert an array of indices (list) to a Boolean marker array where all
777  indices in the list are marked with the given value and the rest are set
778  to zero. */
779  static void ListToMarker(const Array<int> &list, int marker_size,
780  Array<int> &marker, int mark_val = -1);
781 
782  /** @brief For a partially conforming FE space, convert a marker array (nonzero
783  entries are true) on the partially conforming dofs to a marker array on
784  the conforming dofs. A conforming dofs is marked iff at least one of its
785  dependent dofs is marked. */
786  void ConvertToConformingVDofs(const Array<int> &dofs, Array<int> &cdofs);
787 
788  /** @brief For a partially conforming FE space, convert a marker array (nonzero
789  entries are true) on the conforming dofs to a marker array on the
790  (partially conforming) dofs. A dof is marked iff it depends on a marked
791  conforming dofs, where dependency is defined by the ConformingRestriction
792  matrix; in other words, a dof is marked iff it corresponds to a marked
793  conforming dof. */
794  void ConvertFromConformingVDofs(const Array<int> &cdofs, Array<int> &dofs);
795 
796  /** @brief Generate the global restriction matrix from a discontinuous
797  FE space to the continuous FE space of the same polynomial degree. */
799 
800  /** @brief Generate the global restriction matrix from a discontinuous
801  FE space to the piecewise constant FE space. */
803 
804  /** @brief Construct the restriction matrix from the FE space given by
805  (*this) to the lower degree FE space given by (*lfes) which
806  is defined on the same mesh. */
808 
809  /** @brief Construct and return an Operator that can be used to transfer
810  GridFunction data from @a coarse_fes, defined on a coarse mesh, to @a
811  this FE space, defined on a refined mesh. */
812  /** It is assumed that the mesh of this FE space is a refinement of the mesh
813  of @a coarse_fes and the CoarseFineTransformations returned by the method
814  Mesh::GetRefinementTransforms() of the refined mesh are set accordingly.
815  The Operator::Type of @a T can be set to request an Operator of the set
816  type. Currently, only Operator::MFEM_SPARSEMAT and Operator::ANY_TYPE
817  (matrix-free) are supported. When Operator::ANY_TYPE is requested, the
818  choice of the particular Operator sub-class is left to the method. This
819  method also works in parallel because the transfer operator is local to
820  the MPI task when the input is a synchronized ParGridFunction. */
821  void GetTransferOperator(const FiniteElementSpace &coarse_fes,
822  OperatorHandle &T) const;
823 
824  /** @brief Construct and return an Operator that can be used to transfer
825  true-dof data from @a coarse_fes, defined on a coarse mesh, to @a this FE
826  space, defined on a refined mesh.
827 
828  This method calls GetTransferOperator() and multiplies the result by the
829  prolongation operator of @a coarse_fes on the right, and by the
830  restriction operator of this FE space on the left.
831 
832  The Operator::Type of @a T can be set to request an Operator of the set
833  type. In serial, the supported types are: Operator::MFEM_SPARSEMAT and
834  Operator::ANY_TYPE (matrix-free). In parallel, the supported types are:
835  Operator::Hypre_ParCSR and Operator::ANY_TYPE. Any other type is treated
836  as Operator::ANY_TYPE: the operator representation choice is made by this
837  method. */
838  virtual void GetTrueTransferOperator(const FiniteElementSpace &coarse_fes,
839  OperatorHandle &T) const;
840 
841  /** @brief Reflect changes in the mesh: update number of DOFs, etc. Also, calculate
842  GridFunction transformation operator (unless want_transform is false).
843  Safe to call multiple times, does nothing if space already up to date. */
844  virtual void Update(bool want_transform = true);
845 
846  /// Get the GridFunction update operator.
847  const Operator* GetUpdateOperator() { Update(); return Th.Ptr(); }
848 
849  /// Return the update operator in the given OperatorHandle, @a T.
851 
852  /** @brief Set the ownership of the update operator: if set to false, the
853  Operator returned by GetUpdateOperator() must be deleted outside the
854  FiniteElementSpace. */
855  /** The update operator ownership is automatically reset to true when a new
856  update operator is created by the Update() method. */
857  void SetUpdateOperatorOwner(bool own) { Th.SetOperatorOwner(own); }
858 
859  /// Specify the Operator::Type to be used by the update operators.
860  /** The default type is Operator::ANY_TYPE which leaves the choice to this
861  class. The other currently supported option is Operator::MFEM_SPARSEMAT
862  which is only guaranteed to be honored for a refinement update operator.
863  Any other type will be treated as Operator::ANY_TYPE.
864  @note This operation destroys the current update operator (if owned). */
866 
867  /// Free the GridFunction update operator (if any), to save memory.
868  virtual void UpdatesFinished() { Th.Clear(); }
869 
870  /** Return update counter, similar to Mesh::GetSequence(). Used by
871  GridFunction to check if it is up to date with the space. */
872  long GetSequence() const { return sequence; }
873 
874  /// Return whether or not the space is discontinuous (L2)
875  bool IsDGSpace() const
876  {
877  return dynamic_cast<const L2_FECollection*>(fec) != NULL;
878  }
879 
880  /** In variable order spaces on nonconforming (NC) meshes, this function
881  controls whether strict conformity is enforced in cases where coarse
882  edges/faces have higher polynomial order than their fine NC neighbors.
883  In the default (strict) case, the coarse side polynomial order is
884  reduced to that of the lowest order fine edge/face, so all fine
885  neighbors can interpolate the coarse side exactly. If relaxed == true,
886  some discontinuities in the solution in such cases are allowed and the
887  coarse side is not restricted. For an example, see
888  https://github.com/mfem/mfem/pull/1423#issuecomment-621340392 */
889  void SetRelaxedHpConformity(bool relaxed = true)
890  {
891  relaxed_hp = relaxed;
892  orders_changed = true; // force update
893  Update(false);
894  }
895 
896  /// Save finite element space to output stream @a out.
897  void Save(std::ostream &out) const;
898 
899  /** @brief Read a FiniteElementSpace from a stream. The returned
900  FiniteElementCollection is owned by the caller. */
901  FiniteElementCollection *Load(Mesh *m, std::istream &input);
902 
903  virtual ~FiniteElementSpace();
904 };
905 
906 
907 /// Class representing the storage layout of a QuadratureFunction.
908 /** Multiple QuadratureFunction%s can share the same QuadratureSpace. */
910 {
911 protected:
912  friend class QuadratureFunction; // Uses the element_offsets.
913 
915  int order;
916  int size;
917 
919  int *element_offsets; // scalar offsets; size = number of elements + 1
920 
921  // protected functions
922 
923  // Assuming mesh and order are set, construct the members: int_rule,
924  // element_offsets, and size.
925  void Construct();
926 
927 public:
928  /// Create a QuadratureSpace based on the global rules from #IntRules.
929  QuadratureSpace(Mesh *mesh_, int order_)
930  : mesh(mesh_), order(order_) { Construct(); }
931 
932  /// Read a QuadratureSpace from the stream @a in.
933  QuadratureSpace(Mesh *mesh_, std::istream &in);
934 
935  virtual ~QuadratureSpace() { delete [] element_offsets; }
936 
937  /// Return the total number of quadrature points.
938  int GetSize() const { return size; }
939 
940  /// Return the order of the quadrature rule(s) used by all elements.
941  int GetOrder() const { return order; }
942 
943  /// Returns the mesh
944  inline Mesh *GetMesh() const { return mesh; }
945 
946  /// Returns number of elements in the mesh.
947  inline int GetNE() const { return mesh->GetNE(); }
948 
949  /// Get the IntegrationRule associated with mesh element @a idx.
950  const IntegrationRule &GetElementIntRule(int idx) const
951  { return *int_rule[mesh->GetElementBaseGeometry(idx)]; }
952 
953  /// Write the QuadratureSpace to the stream @a out.
954  void Save(std::ostream &out) const;
955 };
956 
957 inline bool UsesTensorBasis(const FiniteElementSpace& fes)
958 {
959  // TODO: mixed meshes: return true if there is at least one tensor-product
960  // Geometry in the global mesh and the FE collection returns a
961  // TensorBasisElement for that Geometry?
962 
963  // Potential issue: empty local mesh --> no element 0.
964  return dynamic_cast<const mfem::TensorBasisElement *>(fes.GetFE(0))!=nullptr;
965 }
966 
967 }
968 
969 #endif
Abstract class for all finite elements.
Definition: fe.hpp:243
int GetNFbyType(FaceType type) const
Returns the number of faces according to the requested type.
Definition: fespace.hpp:584
Ordering::Type GetOrdering() const
Return the ordering method.
Definition: fespace.hpp:552
int Size() const
Return the logical size of the array.
Definition: array.hpp:134
Abstract base class for LORDiscretization and ParLORDiscretization classes, which construct low-order...
Definition: lor.hpp:22
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:540
void GetVertexVDofs(int i, Array< int > &vdofs) const
Definition: fespace.cpp:286
int GetAttribute(int i) const
Definition: fespace.hpp:612
int GetNDofs() const
Returns number of degrees of freedom.
Definition: fespace.hpp:537
int ndofs
Number of degrees of freedom. Number of unknowns is ndofs * vdim.
Definition: fespace.hpp:110
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
Definition: mesh.hpp:1203
int * bdofs
internal DOFs of elements if mixed/var-order; NULL otherwise
Definition: fespace.hpp:118
const Table & GetBdrElementToDofTable() const
Return a reference to the internal Table that stores the lists of scalar dofs, for each boundary mesh...
Definition: fespace.hpp:705
void GetElementInteriorVDofs(int i, Array< int > &vdofs) const
Definition: fespace.cpp:292
int DofToVDof(int dof, int vd, int ndofs=-1) const
Definition: fespace.cpp:233
const SparseMatrix * GetConformingProlongation() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:1168
static void AddDependencies(SparseMatrix &deps, Array< int > &master_dofs, Array< int > &slave_dofs, DenseMatrix &I, int skipfirst=0)
Definition: fespace.cpp:700
virtual const Operator * GetRestrictionTransposeOperator() const
Return an operator that performs the transpose of GetRestrictionOperator.
Definition: fespace.hpp:454
virtual void Update(bool want_transform=true)
Reflect changes in the mesh: update number of DOFs, etc. Also, calculate GridFunction transformation ...
Definition: fespace.cpp:2905
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
Definition: fespace.hpp:432
static int Map(int ndofs, int vdim, int dof, int vd)
void GetElementVertices(int i, Array< int > &vertices) const
Returns the vertices of element i.
Definition: fespace.hpp:592
virtual const SparseMatrix * GetHpRestrictionMatrix() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.hpp:468
void BuildElementToDofTable() const
Definition: fespace.cpp:304
const IntegrationRule & GetElementIntRule(int idx) const
Get the IntegrationRule associated with mesh element idx.
Definition: fespace.hpp:950
int MakeDofTable(int ent_dim, const Array< int > &entity_orders, Table &entity_dofs, Array< char > *var_ent_order)
Definition: fespace.cpp:2187
void BuildNURBSFaceToDofTable() const
Generates partial face_dof table for a NURBS space.
Definition: fespace.cpp:1889
const FiniteElement * GetEdgeElement(int i, int variant=0) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th edge in the ...
Definition: fespace.cpp:2749
DerefinementOperator(const FiniteElementSpace *f_fes, const FiniteElementSpace *c_fes, BilinearFormIntegrator *mass_integ)
Definition: fespace.cpp:1553
static const int NumGeom
Definition: geom.hpp:41
const FiniteElement * GetTraceElement(int i, Geometry::Type geom_type) const
Return the trace element from element &#39;i&#39; to the given &#39;geom_type&#39;.
Definition: fespace.cpp:2759
Matrix-free transfer operator between finite element spaces on the same mesh.
Definition: transfer.hpp:385
bool Conforming() const
Definition: mesh.hpp:1366
SparseMatrix * RefinementMatrix(int old_ndofs, const Table *old_elem_dof)
Definition: fespace.cpp:1401
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:849
SparseMatrix * DerefinementMatrix(int old_ndofs, const Table *old_elem_dof)
Calculate GridFunction restriction matrix after mesh derefinement.
Definition: fespace.cpp:1709
void GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:262
Ordering::Type ordering
Definition: fespace.hpp:107
const Geometry::Type geom
Definition: ex1.cpp:40
static int DecodeDof(int dof, double &sign)
Definition: fespace.hpp:248
void SetElementOrder(int i, int p)
Sets the order of the i&#39;th finite element.
Definition: fespace.cpp:133
static int DecodeDof(int dof)
Helpers to remove encoded sign from a DOF.
Definition: fespace.hpp:245
SparseMatrix * RefinementMatrix_main(const int coarse_ndofs, const Table &coarse_elem_dof, const DenseTensor localP[]) const
Definition: fespace.cpp:1316
int GetSize() const
Return the total number of quadrature points.
Definition: fespace.hpp:938
static int MinOrder(VarOrderBits bits)
Return the minimum order (least significant bit set) in the bit mask.
Definition: fespace.cpp:2060
Element::Type GetElementType(int i) const
Returns the type of element i.
Definition: mesh.cpp:5748
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
int VDofToDof(int vdof) const
Definition: fespace.hpp:663
static constexpr int MaxVarOrder
Definition: fespace.hpp:204
Element::Type GetBdrElementType(int i) const
Returns the type of boundary element i.
Definition: mesh.cpp:5753
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...
Definition: fespace.cpp:437
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int vdim
Vector dimension (number of unknowns per degree of freedom).
Definition: fespace.hpp:102
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
Definition: table.cpp:191
OperatorHandle L2E_lex
Definition: fespace.hpp:154
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:944
int GetNumElementInteriorDofs(int i) const
Definition: fespace.cpp:2644
int GetElementOrder(int i) const
Returns the order of the i&#39;th finite element.
Definition: fespace.cpp:160
int GetBdrAttribute(int i) const
Definition: fespace.hpp:614
const SparseMatrix * GetHpConformingRestriction() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:1182
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
Definition: fespace.cpp:506
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:846
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...
Definition: fespace.cpp:2571
int GetBdrElementType(int i) const
Returns the type of boundary element i.
Definition: fespace.hpp:596
MFEM_DEPRECATED void RebuildElementToDofTable()
(
Definition: fespace.cpp:375
int GetNV() const
Returns number of vertices in the mesh.
Definition: fespace.hpp:564
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 ...
Definition: fespace.cpp:1497
const FiniteElement * GetFaceElement(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th face in the ...
Definition: fespace.cpp:2719
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, defined on a coarse mesh, to this FE space, defined on a refined mesh.
Definition: fespace.cpp:2842
void AddEdgeFaceDependencies(SparseMatrix &deps, Array< int > &master_dofs, const FiniteElement *master_fe, Array< int > &slave_dofs, int slave_face, const DenseMatrix *pm) const
Definition: fespace.cpp:725
FiniteElementCollection * Load(Mesh *m, std::istream &input)
Read a FiniteElementSpace from a stream. The returned FiniteElementCollection is owned by the caller...
Definition: fespace.cpp:3077
virtual void GetElementDofs(int elem, Array< int > &dofs) const
Returns indices of degrees of freedom of element &#39;elem&#39;.
Definition: fespace.cpp:2298
bool UsesTensorBasis(const FiniteElementSpace &fes)
Definition: fespace.hpp:957
const IntegrationRule * int_rule[Geometry::NumGeom]
Definition: fespace.hpp:918
virtual void CopyProlongationAndRestriction(const FiniteElementSpace &fes, const Array< int > *perm)
Copies the prolongation and restriction matrices from fes.
Definition: fespace.cpp:96
const FiniteElementCollection * fec
Associated FE collection (not owned).
Definition: fespace.hpp:99
OperatorHandle Th
Transformation to apply to GridFunctions after space Update().
Definition: fespace.hpp:151
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...
Definition: fespace.cpp:2490
NURBSExtension * GetNURBSext()
Definition: fespace.hpp:413
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:971
void GetLocalDerefinementMatrices(Geometry::Type geom, DenseTensor &localR) const
Definition: fespace.cpp:1685
int GetConformingVSize() const
Definition: fespace.hpp:549
std::unordered_map< const key_face, FaceRestriction *, key_hash > map_L2F
Definition: fespace.hpp:167
int GetOrder() const
Return the order (polynomial degree) of the FE collection, corresponding to the order/degree returned...
Definition: fe_coll.hpp:195
std::uint64_t VarOrderBits
Bit-mask representing a set of orders needed by an edge/face.
Definition: fespace.hpp:203
const Table & GetFaceToDofTable() const
Return a reference to the internal Table that stores the lists of scalar dofs, for each face in the m...
Definition: fespace.hpp:713
Array< int > dof_elem_array
Definition: fespace.hpp:134
static void MarkerToList(const Array< int > &marker, Array< int > &list)
Convert a Boolean marker array to a list containing all marked indices.
Definition: fespace.cpp:540
void GetFaceInteriorDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:2662
Derefinement operator, used by the friend class InterpolationGridTransfer.
Definition: fespace.hpp:300
ElementTransformation * GetBdrElementTransformation(int i)
Returns the transformation defining the i-th boundary element.
Definition: mesh.cpp:428
virtual void UpdateMeshPointer(Mesh *new_mesh)
Definition: fespace.cpp:3003
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...
Definition: fespace.cpp:580
const Operator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
Definition: fespace.cpp:1195
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition: mesh.cpp:4915
int uni_fdof
of single face DOFs if all faces uniform; -1 otherwise
Definition: fespace.hpp:117
const SparseMatrix * GetConformingRestriction() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:1175
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: fespace.cpp:1654
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:567
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:947
void GetVertexDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:2619
bool Nonconforming() const
Definition: mesh.hpp:1367
FaceType
Definition: mesh.hpp:45
void BuildBdrElementToDofTable() const
Definition: fespace.cpp:327
int FindFaceDof(int face, int ndof) const
Similar to FindEdgeDof, but used for mixed meshes too.
Definition: fespace.hpp:231
int GetFaceOrder(int face, int variant=0) const
Returns the polynomial degree of the i&#39;th face finite element.
Definition: fespace.cpp:2270
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: fespace.cpp:1460
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition: fespace.hpp:576
bool orders_changed
True if at least one element order changed (variable-order space only).
Definition: fespace.hpp:183
Data type sparse matrix.
Definition: sparsemat.hpp:41
Native ordering as defined by the FiniteElement.
int FindEdgeDof(int edge, int ndof) const
Definition: fespace.hpp:227
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:410
Array< char > elem_order
Definition: fespace.hpp:114
void GetVDofs(int vd, Array< int > &dofs, int ndofs=-1) const
Returns the indices of all of the VDofs for the specified dimension &#39;vd&#39;.
Definition: fespace.cpp:177
A class that performs interpolation from an E-vector to quadrature point values and/or derivatives (Q...
Type
Ordering methods:
Definition: fespace.hpp:32
SparseMatrix * cR
Conforming restriction matrix such that cR.cP=I.
Definition: fespace.hpp:145
SparseMatrix * D2Const_GlobalRestrictionMatrix(FiniteElementSpace *cfes)
Generate the global restriction matrix from a discontinuous FE space to the piecewise constant FE spa...
Definition: fespace.cpp:621
SparseMatrix * D2C_GlobalRestrictionMatrix(FiniteElementSpace *cfes)
Generate the global restriction matrix from a discontinuous FE space to the continuous FE space of th...
Definition: fespace.cpp:589
void GetBoundaryTrueDofs(Array< int > &boundary_dofs, int component=-1)
Get a list of all boundary true dofs, boundary_dofs. For spaces with &#39;vdim&#39; &gt; 1, the &#39;component&#39; param...
Definition: fespace.cpp:524
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: fespace.hpp:448
void BuildConformingInterpolation() const
Calculate the cP and cR matrices for a nonconforming mesh.
Definition: fespace.cpp:869
NURBSExtension * StealNURBSext()
Definition: fespace.cpp:1854
int GetElementForDof(int i) const
Return the index of the first element that contains dof i.
Definition: fespace.hpp:723
int GetMaxElementOrder() const
Return the maximum polynomial order.
Definition: fespace.hpp:428
ElementTransformation * GetBdrElementTransformation(int i) const
Returns ElementTransformation for the i-th boundary element.
Definition: fespace.hpp:609
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Returns the transformation defining the i-th element in the user-defined variable ElTr...
Definition: fespace.hpp:605
bool IsDGSpace() const
Return whether or not the space is discontinuous (L2)
Definition: fespace.hpp:875
Array< int > dof_ldof_array
Definition: fespace.hpp:134
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:543
bool Conforming() const
Definition: fespace.hpp:416
SparseMatrix * cP
Definition: fespace.hpp:143
virtual void GetBdrElementDofs(int bel, Array< int > &dofs) const
Returns indices of degrees of freedom for boundary element &#39;bel&#39;.
Definition: fespace.cpp:2418
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:534
int GetNF() const
Returns number of faces (i.e. co-dimension 1 entities) in the mesh.
Definition: fespace.hpp:573
static int EncodeDof(int entity_base, int idx)
Helper to encode a sign flip into a DOF index (for Hcurl/Hdiv shapes).
Definition: fespace.hpp:241
virtual const FaceRestriction * GetFaceRestriction(ElementDofOrdering e_ordering, FaceType, L2FaceValues mul=L2FaceValues::DoubleValued) const
Return an Operator that converts L-vectors to E-vectors on each face.
Definition: fespace.cpp:1228
Array< QuadratureInterpolator * > E2Q_array
Definition: fespace.hpp:170
A class that performs interpolation from a face E-vector to quadrature point values and/or derivative...
int GetNumBorderDofs(Geometry::Type geom, int order) const
Definition: fespace.cpp:836
Operator * Ptr() const
Access the underlying Operator pointer.
Definition: handle.hpp:82
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.).
Definition: fespace.cpp:844
int FirstFaceDof(int face, int variant=0) const
Definition: fespace.hpp:234
OperatorHandle L2E_nat
The element restriction operators, see GetElementRestriction().
Definition: fespace.hpp:154
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...
Definition: fespace.cpp:2817
void GetElementInteriorDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:2629
void GetLocalRefinementMatrices(Geometry::Type geom, DenseTensor &localP) const
Definition: fespace.cpp:1378
static void DofsToVDofs(int ndofs, int vdim, Array< int > &dofs)
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition: fespace.hpp:600
Type
Enumeration defining IDs for some classes derived from Operator.
Definition: operator.hpp:252
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:34
RefinementOperator(const FiniteElementSpace *fespace, Table *old_elem_dof, int old_ndofs)
Definition: fespace.cpp:1419
const NURBSExtension * GetNURBSext() const
Definition: fespace.hpp:412
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
Transfer data between a coarse mesh and an embedded refined mesh using interpolation.
Definition: transfer.hpp:121
int GetLocalDofForDof(int i) const
Return the local dof index in the first element that contains dof i.
Definition: fespace.hpp:727
virtual ~FiniteElementSpace()
Definition: fespace.cpp:2764
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
SparseMatrix * H2L_GlobalRestrictionMatrix(FiniteElementSpace *lfes)
Construct the restriction matrix from the FE space given by (*this) to the lower degree FE space give...
Definition: fespace.cpp:652
virtual ~QuadratureSpace()
Definition: fespace.hpp:935
void GetUpdateOperator(OperatorHandle &T)
Return the update operator in the given OperatorHandle, T.
Definition: fespace.hpp:850
int GetOrder(int i) const
Returns the polynomial degree of the i&#39;th finite element.
Definition: fespace.hpp:524
void Save(std::ostream &out) const
Write the QuadratureSpace to the stream out.
Definition: fespace.cpp:3225
Array< char > var_face_orders
Definition: fespace.hpp:127
void Clear()
Clear the OperatorHandle, deleting the held Operator (if owned), while leaving the type id unchanged...
Definition: handle.hpp:119
void MakeVDimMatrix(SparseMatrix &mat) const
Replicate &#39;mat&#39; in the vector dimension, according to vdim ordering mode.
Definition: fespace.cpp:1140
The ordering method used when the number of unknowns per mesh node (vector dimension) is bigger than ...
Definition: fespace.hpp:28
void GetEdgeInteriorVDofs(int i, Array< int > &vdofs) const
Definition: fespace.cpp:298
void SetRelaxedHpConformity(bool relaxed=true)
Definition: fespace.hpp:889
int GetNVariants(int entity, int index) const
Return number of possible DOF variants for edge/face (var. order spaces).
Definition: fespace.cpp:2286
int GetOrder() const
Return the order of the quadrature rule(s) used by all elements.
Definition: fespace.hpp:941
int GetElementType(int i) const
Returns the type of element i.
Definition: fespace.hpp:588
FiniteElementSpace(Mesh *mesh, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES)
Definition: fespace.hpp:395
bool Nonconforming() const
Definition: fespace.hpp:417
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:843
int GetNVDofs() const
Number of all scalar vertex dofs.
Definition: fespace.hpp:557
int GetDegenerateFaceDofs(int index, Array< int > &dofs, Geometry::Type master_geom, int variant) const
Definition: fespace.cpp:794
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...
Definition: fespace.cpp:572
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
Definition: handle.hpp:115
void Swap(Mesh &other, bool non_geometry)
Definition: mesh.cpp:8612
A standard isoparametric element transformation.
Definition: eltrans.hpp:348
Mesh * mesh
The mesh that FE space lives on (not owned).
Definition: fespace.hpp:96
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:347
void ReorderElementToDofTable()
Reorder the scalar DOFs based on the element ordering.
Definition: fespace.cpp:382
GridFunction interpolation operator applicable after mesh refinement.
Definition: fespace.hpp:281
void CalcEdgeFaceVarOrders(Array< VarOrderBits > &edge_orders, Array< VarOrderBits > &face_orders) const
Definition: fespace.cpp:2071
Array< FaceQuadratureInterpolator * > E2BFQ_array
Definition: fespace.hpp:172
int GetNConformingDofs() const
Definition: fespace.cpp:1189
Table var_face_dofs
NOTE: also used for spaces with mixed faces.
Definition: fespace.hpp:123
void GetEdgeInteriorDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:2650
const QuadratureInterpolator * GetQuadratureInterpolator(const IntegrationRule &ir) const
Return a QuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivat...
Definition: fespace.cpp:1256
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition: fespace.hpp:65
Array< FaceQuadratureInterpolator * > E2IFQ_array
Definition: fespace.hpp:171
void BuildFaceToDofTable() const
Definition: fespace.cpp:349
void UpdateElementOrders()
Resize the elem_order array on mesh change.
Definition: fespace.cpp:2883
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th element in t...
Definition: fespace.cpp:2388
int index(int i, int j, int nx, int ny)
Definition: life.cpp:237
const Operator * GetUpdateOperator()
Get the GridFunction update operator.
Definition: fespace.hpp:847
Array< char > var_edge_orders
Definition: fespace.hpp:127
int GetEdgeOrder(int edge, int variant=0) const
Definition: fespace.cpp:2259
Lexicographic ordering for tensor-product FiniteElements.
FiniteElementSpace()
Default constructor: the object is invalid until initialized using the method Load().
Definition: fespace.cpp:59
QuadratureSpace(Mesh *mesh_, int order_)
Create a QuadratureSpace based on the global rules from IntRules.
Definition: fespace.hpp:929
void Constructor(Mesh *mesh, NURBSExtension *ext, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES)
Help function for constructors + Load().
Definition: fespace.cpp:1807
void SetUpdateOperatorOwner(bool own)
Set the ownership of the update operator: if set to false, the Operator returned by GetUpdateOperator...
Definition: fespace.hpp:857
void GetFaceVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th face element (2D and 3D).
Definition: fespace.cpp:274
void GetEdgeVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th edge.
Definition: fespace.cpp:280
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...
Definition: fespace.cpp:559
std::size_t operator()(const key_face &k) const
Definition: fespace.hpp:159
Vector data type.
Definition: vector.hpp:60
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:554
std::tuple< bool, ElementDofOrdering, FaceType, L2FaceValues > key_face
The face restriction operators, see GetFaceRestriction().
Definition: fespace.hpp:156
int GetNFbyType(FaceType type) const
Returns the number of faces according to the requested type.
Definition: mesh.cpp:4941
Array< int > face_to_be
Definition: fespace.hpp:138
virtual void UpdatesFinished()
Free the GridFunction update operator (if any), to save memory.
Definition: fespace.hpp:868
const Table & GetElementToDofTable() const
Return a reference to the internal Table that stores the lists of scalar dofs, for each mesh element...
Definition: fespace.hpp:700
SparseMatrix * cR_hp
A version of the conforming restriction matrix for variable-order spaces.
Definition: fespace.hpp:147
void SetUpdateOperatorType(Operator::Type tid)
Specify the Operator::Type to be used by the update operators.
Definition: fespace.hpp:865
static bool DofFinalizable(int dof, const Array< bool > &finalized, const SparseMatrix &deps)
Definition: fespace.cpp:780
Class representing the storage layout of a QuadratureFunction.
Definition: fespace.hpp:909
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th boundary fac...
Definition: fespace.cpp:2686
NURBSExtension * NURBSext
Definition: fespace.hpp:136
void GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:268
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:66
Base class for operators that extracts Face degrees of freedom.
Abstract operator.
Definition: operator.hpp:24
virtual const SparseMatrix * GetRestrictionMatrix() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.hpp:464
long GetSequence() const
Definition: fespace.hpp:872
int GetElementOrderImpl(int i) const
Return element order: internal version of GetElementOrder without checks.
Definition: fespace.cpp:171
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:745
virtual const Operator * GetRestrictionOperator() const
An abstract operator that performs the same action as GetRestrictionMatrix.
Definition: fespace.hpp:460
void Save(std::ostream &out) const
Save finite element space to output stream out.
Definition: fespace.cpp:3008
int GetAttribute(int i) const
Return the attribute of element i.
Definition: mesh.hpp:1197
void BuildDofToArrays()
Initialize internal data that enables the use of the methods GetElementForDof() and GetLocalDofForDof...
Definition: fespace.cpp:402
Class representing a function through its values (scalar or vector) at quadrature points...
Definition: gridfunc.hpp:734
int GetNFDofs() const
Number of all scalar face-interior dofs.
Definition: fespace.hpp:561
int GetNEDofs() const
Number of all scalar edge-interior dofs.
Definition: fespace.hpp:559
void GetElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of element i.
Definition: mesh.hpp:1011
int FindDofs(const Table &var_dof_table, int row, int ndof) const
Search row of a DOF table for a DOF set of size &#39;ndof&#39;, return first DOF.
Definition: fespace.cpp:2242
static void AdjustVDofs(Array< int > &vdofs)
Definition: fespace.cpp:249
void SetType(Operator::Type tid)
Invoke Clear() and set a new type id.
Definition: handle.hpp:127
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:285
void DofsToVDofs(Array< int > &dofs, int ndofs=-1) const
Definition: fespace.cpp:197
const FaceQuadratureInterpolator * GetFaceQuadratureInterpolator(const IntegrationRule &ir, FaceType type) const
Return a FaceQuadratureInterpolator that interpolates E-vectors to quadrature point values and/or der...
Definition: fespace.cpp:1285
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.
Definition: fespace.hpp:404