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