MFEM  v4.6.0
Finite element discretization library
fespace.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, 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.
96 
97  @details The term "degree of freedom", or "dof" for short, can mean
98  different things in different contexts. In MFEM we use "dof" to refer to
99  four closely related types of data; @ref edof "edofs", @ref ldof "ldofs",
100  @ref tdof "tdofs", and @ref vdof "vdofs".
101 
102  @anchor edof @par Element DoF:
103  %Element dofs, sometimes referred to as @b edofs, are the expansion
104  coefficients used to build the linear combination of basis functions which
105  approximate a field within one element of the computational mesh. The
106  arrangement of the element dofs is determined by the basis function and
107  element types.
108  @par
109  %Element dofs are usually accessed one element at a time but they can be
110  concatenated together into a global vector when minimizing access time is
111  crucial. The global number of element dofs is not directly available from
112  the FiniteElementSpace. It can be determined by repeatedly calling
113  FiniteElementSpace::GetElementDofs and summing the lengths of the resulting
114  @a dofs arrays.
115 
116  @anchor ldof @par Local DoF:
117  Most basis function types share many of their element dofs with neighboring
118  elements. Consequently, the global @ref edof "edof" vector suggested above
119  would contain many redundant entries. One of the primary roles of the
120  FiniteElementSpace is to collapse out these redundancies and
121  define a unique ordering of the remaining degrees of freedom. The
122  collapsed set of dofs are called @b "local dofs" or @b ldofs in
123  the MFEM parlance.
124  @par
125  The term @b local in this context refers to the local rank in a parallel
126  processing environment. MFEM can, of course, be used in sequential
127  computing environments but it is designed with parallel processing in mind
128  and this terminology reflects that design focus.
129  @par
130  When running in parallel the set of local dofs contains all of the degrees
131  of freedom associated with locally owned elements. When running in serial
132  all elements are locally owned so all element dofs are represented in the
133  set of local dofs.
134  @par
135  There are two important caveats regarding local dofs. First, some basis
136  function types, Nedelec and Raviart-Thomas are the prime examples, have an
137  orientation associated with each basis function. The relative orientations
138  of such basis functions in neighboring elements can lead to shared degrees
139  of freedom with opposite signs from the point of view of these neighboring
140  elements. MFEM typically chooses the orientation of the first such shared
141  degree of freedom that it encounters as the default orientation for the
142  corresponding local dof. When this local dof is referenced by a neighboring
143  element which happens to require the opposite orientation the local dof
144  index will be returned (by calls to functions such as
145  FiniteElementSpace::GetElementDofs) as a negative integer. In such cases
146  the actual offset into the vector of local dofs is @b -index-1 and the
147  value expected by this element should have the opposite sign to the value
148  stored in the local dof vector.
149  @par
150  The second important caveat only pertains to high order Nedelec basis
151  functions when shared triangular faces are present in the mesh. In this
152  very particular case the relative orientation of the face with respect to
153  its two neighboring elements can lead to different definitions of the
154  degrees of freedom associated with the interior of the face which cannot
155  be handled by simply flipping the signs of the corresponding values. The
156  DofTransformation class is designed to manage the necessary @b edof to
157  @b ldof transformations in this case. In the majority of cases the
158  DofTransformation is unnecessary and a NULL pointer will be returned in
159  place of a pointer to this object. See DofTransformation for more
160  information.
161 
162  @anchor tdof @par True DoF:
163  As the name suggests "true dofs" or @b tdofs form the minimal set of data
164  values needed (along with mesh and basis function definitions) to uniquely
165  define a finite element discretization of a field. The number of true dofs
166  determines the size of the linear systems which typically need to be solved
167  in FEM simulations.
168  @par
169  Often the true dofs and the local dofs are identical, however, there are
170  important cases where they differ significantly. The first such case is
171  related to non-conforming meshes. On non-conforming meshes it is common
172  for degrees of freedom associated with "hanging" nodes, edges, or faces to
173  be constrained by degrees of freedom associated with another mesh entity.
174  In such cases the "hanging" degrees of freedom should not be considered
175  "true" degrees of freedom since their values cannot be independently
176  assigned. For this reason the FiniteElementSpace must process these
177  constraints and define a reduced set of "true" degrees of freedom which are
178  distinct from the local degrees of freedom.
179  @par
180  The second important distinction arises in parallel processing. When
181  distributing a linear system in parallel each degree of freedom must be
182  assigned to a particular processor, its owner. From the finite element
183  point of view it is convenient to distribute a computational mesh and
184  define an owning processor for each element. Since degrees of freedom may
185  be shared between neighboring elements they may also be shared between
186  neighboring processors. Another role of the FiniteElementSpace is to
187  identify the ownership of degrees of freedom which must be shared between
188  processors. Therefore the set of "true" degrees of freedom must also remove
189  redundant degrees of freedom which are owned by other processors.
190  @par
191  To summarize the set of true degrees of freedom are those degrees of
192  freedom needed to solve a linear system representing the partial
193  differential equation being modeled. True dofs differ from "local" dofs by
194  eliminating redundancies across processor boundaries and applying
195  the constraints needed to properly define fields on non-conforming meshes.
196 
197  @anchor vdof @par Vector DoF:
198  %Vector dofs or @b vdofs are related to fields which are constructed using
199  multiple copies of the same set of basis functions. A typical example would
200  be the use of three instances of the scalar H1 basis functions to
201  approximate the x, y, and z components of a displacement vector field in
202  three dimensional space as often seen in elasticity simulations.
203  @par
204  %Vector dofs do not represent a specific index space the way the three
205  previous types of dofs do. Rather they are related to modifications of
206  these other index spaces to accomodate multiple copies of the underlying
207  function spaces.
208  @par
209  When using @b vdofs, i.e. when @b vdim != 1, the FiniteElementSpace only
210  manages a single set of degrees of freedom and then uses simple rules to
211  determine the appropriate offsets into the full index spaces. Two ordering
212  rules are supported; @b byNODES and @b byVDIM. See Ordering::Type for
213  details.
214  @par
215  Clearly the notion of a @b vdof is relevant in each of the three contexts
216  mentioned above so extra care must be taken whenever @b vdim != 1 to ensure
217  that the @b edof, @b ldof, or @b tdof is being interpreted correctly.
218  */
220 {
223  friend void Mesh::Swap(Mesh &, bool);
224  friend class LORBase;
225 
226 protected:
227  /// The mesh that FE space lives on (not owned).
229 
230  /// Associated FE collection (not owned).
232 
233  /// %Vector dimension (number of unknowns per degree of freedom).
234  int vdim;
235 
236  /** Type of ordering of the vector dofs when #vdim > 1.
237  - Ordering::byNODES - first nodes, then vector dimension,
238  - Ordering::byVDIM - first vector dimension, then nodes */
240 
241  /// Number of degrees of freedom. Number of unknowns is #ndofs * #vdim.
242  int ndofs;
243 
244  /** Polynomial order for each element. If empty, all elements are assumed
245  to be of the default order (fec->GetOrder()). */
247 
249  int uni_fdof; ///< # of single face DOFs if all faces uniform; -1 otherwise
250  int *bdofs; ///< internal DOFs of elements if mixed/var-order; NULL otherwise
251 
252  /** Variable order spaces only: DOF assignments for edges and faces, see
253  docs in MakeDofTable. For constant order spaces the tables are empty. */
255  Table var_face_dofs; ///< NOTE: also used for spaces with mixed faces
256 
257  /** Additional data for the var_*_dofs tables: individual variant orders
258  (these are basically alternate J arrays for var_edge/face_dofs). */
260 
261  // precalculated DOFs for each element, boundary element, and face
262  mutable Table *elem_dof; // owned (except in NURBS FE space)
263  mutable Table *elem_fos; // face orientations by element index
264  mutable Table *bdr_elem_dof; // owned (except in NURBS FE space)
265  mutable Table *bdr_elem_fos; // bdr face orientations by bdr element index
266  mutable Table *face_dof; // owned; in var-order space contains variant 0 DOFs
267 
269 
271  int own_ext;
272  mutable Array<int> face_to_be; // NURBS FE space only
273 
276 
277  /** Matrix representing the prolongation from the global conforming dofs to
278  a set of intermediate partially conforming dofs, e.g. the dofs associated
279  with a "cut" space on a non-conforming mesh. */
280  mutable std::unique_ptr<SparseMatrix> cP;
281  /// Conforming restriction matrix such that cR.cP=I.
282  mutable std::unique_ptr<SparseMatrix> cR;
283  /// A version of the conforming restriction matrix for variable-order spaces.
284  mutable std::unique_ptr<SparseMatrix> cR_hp;
285  mutable bool cP_is_set;
286  /// Operator computing the action of the transpose of the restriction.
287  mutable std::unique_ptr<Operator> R_transpose;
288 
289  /// Transformation to apply to GridFunctions after space Update().
291 
292  /// The element restriction operators, see GetElementRestriction().
294  /// The face restriction operators, see GetFaceRestriction().
295  using key_face = std::tuple<bool, ElementDofOrdering, FaceType, L2FaceValues>;
296  struct key_hash
297  {
298  std::size_t operator()(const key_face& k) const
299  {
300  return std::get<0>(k)
301  + 2 * (int)std::get<1>(k)
302  + 4 * (int)std::get<2>(k)
303  + 8 * (int)std::get<3>(k);
304  }
305  };
306  using map_L2F = std::unordered_map<const key_face,FaceRestriction*,key_hash>;
307  mutable map_L2F L2F;
308 
312 
313  /** Update counter, incremented every time the space is constructed/updated.
314  Used by GridFunctions to check if they are up to date with the space. */
315  long sequence;
316 
317  /** Mesh sequence number last seen when constructing the space. The space
318  needs updating if Mesh::GetSequence() is larger than this. */
320 
321  /// True if at least one element order changed (variable-order space only).
323 
324  bool relaxed_hp; // see SetRelaxedHpConformity()
325 
326  void UpdateNURBS();
327 
328  void Construct();
329  void Destroy();
330 
331  void ConstructDoFTrans();
332  void DestroyDoFTrans();
333 
334  void BuildElementToDofTable() const;
335  void BuildBdrElementToDofTable() const;
336  void BuildFaceToDofTable() const;
337 
338  /** @brief Generates partial face_dof table for a NURBS space.
339 
340  The table is only defined for exterior faces that coincide with a
341  boundary. */
342  void BuildNURBSFaceToDofTable() const;
343 
344  /// Bit-mask representing a set of orders needed by an edge/face.
345  typedef std::uint64_t VarOrderBits;
346  static constexpr int MaxVarOrder = 8*sizeof(VarOrderBits) - 1;
347 
348  /// Return the minimum order (least significant bit set) in the bit mask.
349  static int MinOrder(VarOrderBits bits);
350 
351  /// Return element order: internal version of GetElementOrder without checks.
352  int GetElementOrderImpl(int i) const;
353 
354  /** In a variable order space, calculate a bitmask of polynomial orders that
355  need to be represented on each edge and face. */
356  void CalcEdgeFaceVarOrders(Array<VarOrderBits> &edge_orders,
357  Array<VarOrderBits> &face_orders) const;
358 
359  /** Build the table var_edge_dofs (or var_face_dofs) in a variable order
360  space; return total edge/face DOFs. */
361  int MakeDofTable(int ent_dim, const Array<int> &entity_orders,
362  Table &entity_dofs, Array<char> *var_ent_order);
363 
364  /// Search row of a DOF table for a DOF set of size 'ndof', return first DOF.
365  int FindDofs(const Table &var_dof_table, int row, int ndof) const;
366 
367  /** In a variable order space, return edge DOFs associated with a polynomial
368  order that has 'ndof' degrees of freedom. */
369  int FindEdgeDof(int edge, int ndof) const
370  { return FindDofs(var_edge_dofs, edge, ndof); }
371 
372  /// Similar to FindEdgeDof, but used for mixed meshes too.
373  int FindFaceDof(int face, int ndof) const
374  { return FindDofs(var_face_dofs, face, ndof); }
375 
376  int FirstFaceDof(int face, int variant = 0) const
377  { return uni_fdof >= 0 ? face*uni_fdof : var_face_dofs.GetRow(face)[variant];}
378 
379  /// Return number of possible DOF variants for edge/face (var. order spaces).
380  int GetNVariants(int entity, int index) const;
381 
382  /// Helper to get vertex, edge or face DOFs (entity=0,1,2 resp.).
383  int GetEntityDofs(int entity, int index, Array<int> &dofs,
384  Geometry::Type master_geom = Geometry::INVALID,
385  int variant = 0) const;
386 
387  // Get degenerate face DOFs: see explanation in method implementation.
388  int GetDegenerateFaceDofs(int index, Array<int> &dofs,
389  Geometry::Type master_geom, int variant) const;
390 
391  int GetNumBorderDofs(Geometry::Type geom, int order) const;
392 
393  /// Calculate the cP and cR matrices for a nonconforming mesh.
394  void BuildConformingInterpolation() const;
395 
396  static void AddDependencies(SparseMatrix& deps, Array<int>& master_dofs,
397  Array<int>& slave_dofs, DenseMatrix& I,
398  int skipfirst = 0);
399 
400  static bool DofFinalizable(int dof, const Array<bool>& finalized,
401  const SparseMatrix& deps);
402 
403  void AddEdgeFaceDependencies(SparseMatrix &deps, Array<int>& master_dofs,
404  const FiniteElement *master_fe,
405  Array<int> &slave_dofs, int slave_face,
406  const DenseMatrix *pm) const;
407 
408  /// Replicate 'mat' in the vector dimension, according to vdim ordering mode.
409  void MakeVDimMatrix(SparseMatrix &mat) const;
410 
411  /// GridFunction interpolation operator applicable after mesh refinement.
413  {
414  const FiniteElementSpace* fespace;
416  Table* old_elem_dof; // Owned.
417  Table* old_elem_fos; // Owned.
418 
419  Array<DofTransformation*> old_DoFTrans;
420  mutable VDofTransformation old_VDoFTrans;
421 
422  void ConstructDoFTrans();
423 
424  public:
425  /** Construct the operator based on the elem_dof table of the original
426  (coarse) space. The class takes ownership of the table. */
427  RefinementOperator(const FiniteElementSpace* fespace,
428  Table *old_elem_dof/*takes ownership*/,
429  Table *old_elem_fos/*takes ownership*/, int old_ndofs);
430  RefinementOperator(const FiniteElementSpace *fespace,
431  const FiniteElementSpace *coarse_fes);
432  virtual void Mult(const Vector &x, Vector &y) const;
433  virtual void MultTranspose(const Vector &x, Vector &y) const;
434  virtual ~RefinementOperator();
435  };
436 
437  /// Derefinement operator, used by the friend class InterpolationGridTransfer.
439  {
440  const FiniteElementSpace *fine_fes; // Not owned.
442  Table *coarse_elem_dof; // Owned.
443  // Table *coarse_elem_fos; // Owned.
444  Table coarse_to_fine;
445  Array<int> coarse_to_ref_type;
446  Array<Geometry::Type> ref_type_to_geom;
447  Array<int> ref_type_to_fine_elem_offset;
448 
449  public:
451  const FiniteElementSpace *c_fes,
452  BilinearFormIntegrator *mass_integ);
453  virtual void Mult(const Vector &x, Vector &y) const;
454  virtual ~DerefinementOperator();
455  };
456 
457  /** This method makes the same assumptions as the method:
458  void GetLocalRefinementMatrices(
459  const FiniteElementSpace &coarse_fes, Geometry::Type geom,
460  DenseTensor &localP) const
461  which is defined below. It also assumes that the coarse fes and this have
462  the same vector dimension, vdim. */
463  SparseMatrix *RefinementMatrix_main(const int coarse_ndofs,
464  const Table &coarse_elem_dof,
465  const Table *coarse_elem_fos,
466  const DenseTensor localP[]) const;
467 
469  DenseTensor &localP) const;
471  DenseTensor &localR) const;
472 
473  /** Calculate explicit GridFunction interpolation matrix (after mesh
474  refinement). NOTE: consider using the RefinementOperator class instead
475  of the fully assembled matrix, which can take a lot of memory. */
476  SparseMatrix* RefinementMatrix(int old_ndofs, const Table* old_elem_dof,
477  const Table* old_elem_fos);
478 
479  /// Calculate GridFunction restriction matrix after mesh derefinement.
480  SparseMatrix* DerefinementMatrix(int old_ndofs, const Table* old_elem_dof,
481  const Table* old_elem_fos);
482 
483  /** @brief Return in @a localP the local refinement matrices that map
484  between fespaces after mesh refinement. */
485  /** This method assumes that this->mesh is a refinement of coarse_fes->mesh
486  and that the CoarseFineTransformations of this->mesh are set accordingly.
487  Another assumption is that the FEs of this use the same MapType as the FEs
488  of coarse_fes. Finally, it assumes that the spaces this and coarse_fes are
489  NOT variable-order spaces. */
490  void GetLocalRefinementMatrices(const FiniteElementSpace &coarse_fes,
491  Geometry::Type geom,
492  DenseTensor &localP) const;
493 
494  /// Help function for constructors + Load().
495  void Constructor(Mesh *mesh, NURBSExtension *ext,
497  int vdim = 1, int ordering = Ordering::byNODES);
498 
499  /// Updates the internal mesh pointer. @warning @a new_mesh must be
500  /// <b>topologically identical</b> to the existing mesh. Used if the address
501  /// of the Mesh object has changed, e.g. in @a Mesh::Swap.
502  virtual void UpdateMeshPointer(Mesh *new_mesh);
503 
504  /// Resize the elem_order array on mesh change.
505  void UpdateElementOrders();
506 
507  /// @brief Copies the prolongation and restriction matrices from @a fes.
508  ///
509  /// Used for low order preconditioning on non-conforming meshes. If the DOFs
510  /// require a permutation, it will be supplied by non-NULL @a perm. NULL @a
511  /// perm indicates that no permutation is required.
512  virtual void CopyProlongationAndRestriction(const FiniteElementSpace &fes,
513  const Array<int> *perm);
514 
515 public:
516  /** @brief Default constructor: the object is invalid until initialized using
517  the method Load(). */
519 
520  /** @brief Copy constructor: deep copy all data from @a orig except the Mesh,
521  the FiniteElementCollection, and some derived data. */
522  /** If the @a mesh or @a fec pointers are NULL (default), then the new
523  FiniteElementSpace will reuse the respective pointers from @a orig. If
524  any of these pointers is not NULL, the given pointer will be used instead
525  of the one used by @a orig.
526 
527  @note The objects pointed to by the @a mesh and @a fec parameters must be
528  either the same objects as the ones used by @a orig, or copies of them.
529  Otherwise, the behavior is undefined.
530 
531  @note Derived data objects, such as the conforming prolongation and
532  restriction matrices, and the update operator, will not be copied, even
533  if they are created in the @a orig object. */
534  FiniteElementSpace(const FiniteElementSpace &orig, Mesh *mesh = NULL,
535  const FiniteElementCollection *fec = NULL);
536 
539  int vdim = 1, int ordering = Ordering::byNODES)
540  { Constructor(mesh, NULL, fec, vdim, ordering); }
541 
542  /// Construct a NURBS FE space based on the given NURBSExtension, @a ext.
543  /** @note If the pointer @a ext is NULL, this constructor is equivalent to
544  the standard constructor with the same arguments minus the
545  NURBSExtension, @a ext. */
548  int vdim = 1, int ordering = Ordering::byNODES)
549  { Constructor(mesh, ext, fec, vdim, ordering); }
550 
551  /// Copy assignment not supported
553 
554  /// Returns the mesh
555  inline Mesh *GetMesh() const { return mesh; }
556 
557  const NURBSExtension *GetNURBSext() const { return NURBSext; }
560 
561  bool Conforming() const { return mesh->Conforming() && cP == NULL; }
562  bool Nonconforming() const { return mesh->Nonconforming() || cP != NULL; }
563 
564  /// Sets the order of the i'th finite element.
565  /** By default, all elements are assumed to be of fec->GetOrder(). Once
566  SetElementOrder is called, the space becomes a variable order space. */
567  void SetElementOrder(int i, int p);
568 
569  /// Returns the order of the i'th finite element.
570  int GetElementOrder(int i) const;
571 
572  /// Return the maximum polynomial order.
573  int GetMaxElementOrder() const
574  { return IsVariableOrder() ? elem_order.Max() : fec->GetOrder(); }
575 
576  /// Returns true if the space contains elements of varying polynomial orders.
577  bool IsVariableOrder() const { return elem_order.Size(); }
578 
579  /// The returned SparseMatrix is owned by the FiniteElementSpace.
581 
582  /// The returned SparseMatrix is owned by the FiniteElementSpace.
583  const SparseMatrix *GetConformingRestriction() const;
584 
585  /** Return a version of the conforming restriction matrix for variable-order
586  spaces with complex hp interfaces, where some true DOFs are not owned by
587  any elements and need to be interpolated from higher order edge/face
588  variants (see also @a SetRelaxedHpConformity()). */
589  /// The returned SparseMatrix is owned by the FiniteElementSpace.
591 
592  /// The returned Operator is owned by the FiniteElementSpace.
593  virtual const Operator *GetProlongationMatrix() const
594  { return GetConformingProlongation(); }
595 
596  /// Return an operator that performs the transpose of GetRestrictionOperator
597  /** The returned operator is owned by the FiniteElementSpace.
598 
599  For a serial conforming space, this returns NULL, indicating the identity
600  operator.
601 
602  For a parallel conforming space, this will return a matrix-free
603  (Device)ConformingProlongationOperator.
604 
605  For a non-conforming mesh this will return a TransposeOperator wrapping
606  the restriction matrix. */
608 
609  /// An abstract operator that performs the same action as GetRestrictionMatrix
610  /** In some cases this is an optimized matrix-free implementation. The
611  returned operator is owned by the FiniteElementSpace. */
612  virtual const Operator *GetRestrictionOperator() const
613  { return GetConformingRestriction(); }
614 
615  /// The returned SparseMatrix is owned by the FiniteElementSpace.
616  virtual const SparseMatrix *GetRestrictionMatrix() const
617  { return GetConformingRestriction(); }
618 
619  /// The returned SparseMatrix is owned by the FiniteElementSpace.
620  virtual const SparseMatrix *GetHpRestrictionMatrix() const
621  { return GetHpConformingRestriction(); }
622 
623  /// Return an Operator that converts L-vectors to E-vectors.
624  /** An L-vector is a vector of size GetVSize() which is the same size as a
625  GridFunction. An E-vector represents the element-wise discontinuous
626  version of the FE space.
627 
628  The layout of the E-vector is: ND x VDIM x NE, where ND is the number of
629  degrees of freedom, VDIM is the vector dimension of the FE space, and NE
630  is the number of the mesh elements.
631 
632  The parameter @a e_ordering describes how the local DOFs in each element
633  should be ordered in the E-vector, see ElementDofOrdering.
634 
635  For discontinuous spaces, the element restriction corresponds to a
636  permutation of the degrees of freedom, implemented by the
637  L2ElementRestriction class.
638 
639  The returned Operator is owned by the FiniteElementSpace. */
641  ElementDofOrdering e_ordering) const;
642 
643  /// Return an Operator that converts L-vectors to E-vectors on each face.
644  virtual const FaceRestriction *GetFaceRestriction(
645  ElementDofOrdering f_ordering, FaceType,
647 
648  /** @brief Return a QuadratureInterpolator that interpolates E-vectors to
649  quadrature point values and/or derivatives (Q-vectors). */
650  /** An E-vector represents the element-wise discontinuous version of the FE
651  space and can be obtained, for example, from a GridFunction using the
652  Operator returned by GetElementRestriction().
653 
654  All elements will use the same IntegrationRule, @a ir as the target
655  quadrature points.
656 
657  @note The returned pointer is shared. A good practice, before using it,
658  is to set all its properties to their expected values, as other parts of
659  the code may also change them. That is, it's good to call
660  SetOutputLayout() and DisableTensorProducts() before interpolating. */
662  const IntegrationRule &ir) const;
663 
664  /** @brief Return a QuadratureInterpolator that interpolates E-vectors to
665  quadrature point values and/or derivatives (Q-vectors). */
666  /** An E-vector represents the element-wise discontinuous version of the FE
667  space and can be obtained, for example, from a GridFunction using the
668  Operator returned by GetElementRestriction().
669 
670  The target quadrature points in the elements are described by the given
671  QuadratureSpace, @a qs.
672 
673  @note The returned pointer is shared. A good practice, before using it,
674  is to set all its properties to their expected values, as other parts of
675  the code may also change them. That is, it's good to call
676  SetOutputLayout() and DisableTensorProducts() before interpolating. */
678  const QuadratureSpace &qs) const;
679 
680  /** @brief Return a FaceQuadratureInterpolator that interpolates E-vectors to
681  quadrature point values and/or derivatives (Q-vectors).
682 
683  @note The returned pointer is shared. A good practice, before using it,
684  is to set all its properties to their expected values, as other parts of
685  the code may also change them. That is, it's good to call
686  SetOutputLayout() and DisableTensorProducts() before interpolating. */
688  const IntegrationRule &ir, FaceType type) const;
689 
690  /// Returns the polynomial degree of the i'th finite element.
691  /** NOTE: it is recommended to use GetElementOrder in new code. */
692  int GetOrder(int i) const { return GetElementOrder(i); }
693 
694  /** Return the order of an edge. In a variable order space, return the order
695  of a specific variant, or -1 if there are no more variants. */
696  int GetEdgeOrder(int edge, int variant = 0) const;
697 
698  /// Returns the polynomial degree of the i'th face finite element
699  int GetFaceOrder(int face, int variant = 0) const;
700 
701  /// Returns vector dimension.
702  inline int GetVDim() const { return vdim; }
703 
704  /// @brief Returns number of degrees of freedom.
705  /// This is the number of @ref ldof "Local Degrees of Freedom"
706  inline int GetNDofs() const { return ndofs; }
707 
708  /// @brief Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
709  inline int GetVSize() const { return vdim * ndofs; }
710 
711  /// @brief Return the number of vector true (conforming) dofs.
712  virtual int GetTrueVSize() const { return GetConformingVSize(); }
713 
714  /// Returns the number of conforming ("true") degrees of freedom
715  /// (if the space is on a nonconforming mesh with hanging nodes).
716  int GetNConformingDofs() const;
717 
718  int GetConformingVSize() const { return vdim * GetNConformingDofs(); }
719 
720  /// Return the ordering method.
721  inline Ordering::Type GetOrdering() const { return ordering; }
722 
723  const FiniteElementCollection *FEColl() const { return fec; }
724 
725  /// Number of all scalar vertex dofs
726  int GetNVDofs() const { return nvdofs; }
727  /// Number of all scalar edge-interior dofs
728  int GetNEDofs() const { return nedofs; }
729  /// Number of all scalar face-interior dofs
730  int GetNFDofs() const { return nfdofs; }
731 
732  /// Returns number of vertices in the mesh.
733  inline int GetNV() const { return mesh->GetNV(); }
734 
735  /// Returns number of elements in the mesh.
736  inline int GetNE() const { return mesh->GetNE(); }
737 
738  /// Returns number of faces (i.e. co-dimension 1 entities) in the mesh.
739  /** The co-dimension 1 entities are those that have dimension 1 less than the
740  mesh dimension, e.g. for a 2D mesh, the faces are the 1D entities, i.e.
741  the edges. */
742  inline int GetNF() const { return mesh->GetNumFaces(); }
743 
744  /// Returns number of boundary elements in the mesh.
745  inline int GetNBE() const { return mesh->GetNBE(); }
746 
747  /// Returns the number of faces according to the requested type.
748  /** If type==Boundary returns only the "true" number of boundary faces
749  contrary to GetNBE() that returns "fake" boundary faces associated to
750  visualization for GLVis.
751  Similarly, if type==Interior, the "fake" boundary faces associated to
752  visualization are counted as interior faces. */
753  inline int GetNFbyType(FaceType type) const
754  { return mesh->GetNFbyType(type); }
755 
756  /// Returns the type of element i.
757  inline int GetElementType(int i) const
758  { return mesh->GetElementType(i); }
759 
760  /// Returns the vertices of element i.
761  inline void GetElementVertices(int i, Array<int> &vertices) const
762  { mesh->GetElementVertices(i, vertices); }
763 
764  /// Returns the type of boundary element i.
765  inline int GetBdrElementType(int i) const
766  { return mesh->GetBdrElementType(i); }
767 
768  /// Returns ElementTransformation for the @a i-th element.
770  { return mesh->GetElementTransformation(i); }
771 
772  /** @brief Returns the transformation defining the @a i-th element in the
773  user-defined variable @a ElTr. */
775  { mesh->GetElementTransformation(i, ElTr); }
776 
777  /// Returns ElementTransformation for the @a i-th boundary element.
779  { return mesh->GetBdrElementTransformation(i); }
780 
781  int GetAttribute(int i) const { return mesh->GetAttribute(i); }
782 
783  int GetBdrAttribute(int i) const { return mesh->GetBdrAttribute(i); }
784 
785  /// @anchor getdof @name Local DoF Access Members
786  /// These member functions produce arrays of local degree of freedom
787  /// indices, see @ref ldof. If @b vdim == 1 these indices can be used to
788  /// access entries in GridFunction, LinearForm, and BilinearForm objects.
789  /// If @b vdim != 1 the corresponding @ref getvdof "Get*VDofs" methods
790  /// should be used instead or one of the @ref dof2vdof "DofToVDof" methods
791  /// could be used to produce the appropriate offsets from these local dofs.
792  ///@{
793 
794  /// @brief Returns indices of degrees of freedom of element 'elem'.
795  /// The returned indices are offsets into an @ref ldof vector. See also
796  /// GetElementVDofs().
797  ///
798  /// @note In many cases the returned DofTransformation object will be NULL.
799  /// In other cases see the documentation of the DofTransformation class for
800  /// guidance on its role in performing @ref edof to @ref ldof transformations
801  /// on local vectors and matrices. At present the DofTransformation is only
802  /// needed for Nedelec basis functions of order 2 and above on 3D elements
803  /// with triangular faces.
804  ///
805  /// @note The returned object should NOT be deleted by the caller.
806  virtual DofTransformation *GetElementDofs(int elem, Array<int> &dofs) const;
807 
808  /// @brief Returns indices of degrees of freedom for boundary element 'bel'.
809  /// The returned indices are offsets into an @ref ldof vector. See also
810  /// GetBdrElementVDofs().
811  ///
812  /// @note In many cases the returned DofTransformation object will be NULL.
813  /// In other cases see the documentation of the DofTransformation class for
814  /// guidance on its role in performing @ref edof to @ref ldof transformations
815  /// on local vectors and matrices. At present the DofTransformation is only
816  /// needed for Nedelec basis functions of order 2 and above on 3D elements
817  /// with triangular faces.
818  ///
819  /// @note The returned object should NOT be deleted by the caller.
820  virtual DofTransformation *GetBdrElementDofs(int bel,
821  Array<int> &dofs) const;
822 
823  /** @brief Returns indices of degrees of freedom for NURBS patch index
824  @a patch. Cartesian ordering is used, for the tensor-product degrees of
825  freedom. */
826  void GetPatchDofs(int patch, Array<int> &dofs) const;
827 
828  /// @brief Returns the indices of the degrees of freedom for the specified
829  /// face, including the DOFs for the edges and the vertices of the face.
830  ///
831  /// In variable order spaces, multiple variants of DOFs can be returned.
832  /// See GetEdgeDofs() for more details.
833  /// @return Order of the selected variant, or -1 if there are no more
834  /// variants.
835  ///
836  /// The returned indices are offsets into an @ref ldof vector. See also
837  /// GetFaceVDofs().
838  virtual int GetFaceDofs(int face, Array<int> &dofs, int variant = 0) const;
839 
840  /// @brief Returns the indices of the degrees of freedom for the specified
841  /// edge, including the DOFs for the vertices of the edge.
842  ///
843  /// In variable order spaces, multiple sets of DOFs may exist on an edge,
844  /// corresponding to the different polynomial orders of incident elements.
845  /// The 'variant' parameter is the zero-based index of the desired DOF set.
846  /// The variants are ordered from lowest polynomial degree to the highest.
847  /// @return Order of the selected variant, or -1 if there are no more
848  /// variants.
849  ///
850  /// The returned indices are offsets into an @ref ldof vector. See also
851  /// GetEdgeVDofs().
852  int GetEdgeDofs(int edge, Array<int> &dofs, int variant = 0) const;
853 
854  /// @brief Returns the indices of the degrees of freedom for the specified
855  /// vertices.
856  ///
857  /// The returned indices are offsets into an @ref ldof vector. See also
858  /// GetVertexVDofs().
859  void GetVertexDofs(int i, Array<int> &dofs) const;
860 
861  /// @brief Returns the indices of the degrees of freedom for the interior
862  /// of the specified element.
863  ///
864  /// Specifically this refers to degrees of freedom which are not associated
865  /// with the vertices, edges, or faces of the mesh. This method may be
866  /// useful in conjunction with schemes which process shared and non-shared
867  /// degrees of freedom differently such as static condensation.
868  ///
869  /// The returned indices are offsets into an @ref ldof vector. See also
870  /// GetElementInteriorVDofs().
871  void GetElementInteriorDofs(int i, Array<int> &dofs) const;
872 
873  /// @brief Returns the indices of the degrees of freedom for the interior
874  /// of the specified face.
875  ///
876  /// Specifically this refers to degrees of freedom which are not associated
877  /// with the vertices, edges, or cell interiors of the mesh. This method may
878  /// be useful in conjunction with schemes which process shared and non-shared
879  /// degrees of freedom differently such as static condensation.
880  ///
881  /// The returned indices are offsets into an @ref ldof vector. See also
882  /// GetFaceInteriorVDofs().
883  void GetFaceInteriorDofs(int i, Array<int> &dofs) const;
884 
885  /// @brief Returns the number of degrees of freedom associated with the
886  /// interior of the specified element.
887  ///
888  /// See GetElementInteriorDofs() for more information or to obtain the
889  /// relevant indices.
890  int GetNumElementInteriorDofs(int i) const;
891 
892  /// @brief Returns the indices of the degrees of freedom for the interior
893  /// of the specified edge.
894  ///
895  /// The returned indices are offsets into an @ref ldof vector. See also
896  /// GetEdgeInteriorVDofs().
897  void GetEdgeInteriorDofs(int i, Array<int> &dofs) const;
898  ///@}
899 
900  /// @anchor dof2vdof @name DoF To VDoF Conversion methods
901  /// These methods convert between local dof and local vector dof using the
902  /// appropriate relationship based on the Ordering::Type defined in this
903  /// FiniteElementSpace object.
904  ///
905  /// These methods assume the index set has a range [0, GetNDofs()) which
906  /// will be mapped to the range [0, GetVSize()). This assumption can be
907  /// changed in the forward mappings by passing a value for @a ndofs which
908  /// differs from that returned by GetNDofs().
909  ///
910  /// @note These methods, with the exception of VDofToDof(), are designed to
911  /// produce the correctly encoded values when dof entries are negative,
912  /// see @ref ldof for more on negative dof indices.
913  ///
914  /// @warning When MFEM_DEBUG is enabled at build time the forward mappings
915  /// will verify that each @a dof lies in the proper range. If MFEM_DEBUG is
916  /// disabled no range checking is performed.
917  ///@{
918 
919  /// @brief Returns the indices of all of the VDofs for the specified
920  /// dimension 'vd'.
921  ///
922  /// The @a ndofs parameter can be used to indicate the total number of Dofs
923  /// associated with each component of @b vdim. If @a ndofs is -1 (the
924  /// default value), then the number of Dofs is determined by the
925  /// FiniteElementSpace::GetNDofs().
926  ///
927  /// @note This method does not resize the @a dofs array. It takes the range
928  /// of dofs [0, dofs.Size()) and converts these to @ref vdof "vdofs" and
929  /// stores the results in the @a dofs array.
930  void GetVDofs(int vd, Array<int> &dofs, int ndofs = -1) const;
931 
932  /// @brief Compute the full set of @ref vdof "vdofs" corresponding to each
933  /// entry in @a dofs.
934  ///
935  /// @details Produces a set of @ref vdof "vdofs" of
936  /// length @b vdim * dofs.Size() corresponding to the entries contained in
937  /// the @a dofs array.
938  ///
939  /// The @a ndofs parameter can be used to indicate the total number of Dofs
940  /// associated with each component of @b vdim. If @a ndofs is -1 (the
941  /// default value), then the number of Dofs is <determined by the
942  /// FiniteElementSpace::GetNDofs().
943  ///
944  /// @note The @a dofs array is overwritten and resized to accomodate the
945  /// new values.
946  void DofsToVDofs(Array<int> &dofs, int ndofs = -1) const;
947 
948  /// @brief Compute the set of @ref vdof "vdofs" corresponding to each entry
949  /// in @a dofs for the given vector index @a vd.
950  ///
951  /// The @a ndofs parameter can be used to indicate the total number of Dofs
952  /// associated with each component of @b vdim. If @a ndofs is -1 (the
953  /// default value), then the number of Dofs is <determined by the
954  /// FiniteElementSpace::GetNDofs().
955  ///
956  /// @note The @a dofs array is overwritten with the new values but its size
957  /// will not be altered.
958  void DofsToVDofs(int vd, Array<int> &dofs, int ndofs = -1) const;
959 
960  /// @brief Compute a single @ref vdof corresponding to the index @a dof and
961  /// the vector index @a vd.
962  ///
963  /// The @a ndofs parameter can be used to indicate the total number of Dofs
964  /// associated with each component of @b vdim. If @a ndofs is -1 (the
965  /// default value), then the number of Dofs is <determined by the
966  /// FiniteElementSpace::GetNDofs().
967  int DofToVDof(int dof, int vd, int ndofs = -1) const;
968 
969  /// @brief Compute the inverse of the Dof to VDof mapping for a single
970  /// index @a vdof.
971  ///
972  /// @warning This method is only intended for use with positive indices.
973  /// Passing a negative value for @a vdof will produce an invalid result.
974  int VDofToDof(int vdof) const
975  { return (ordering == Ordering::byNODES) ? (vdof%ndofs) : (vdof/vdim); }
976 
977  ///@}
978 
979  /// @brief Remove the orientation information encoded into an array of dofs
980  /// Some basis function types have a relative orientation associated with
981  /// degrees of freedom shared between neighboring elements, see @ref ldof
982  /// for more information. An orientation mismatch is indicated in the dof
983  /// indices by a negative index value. This method replaces such negative
984  /// indices with the corresponding positive offsets.
985  ///
986  /// @note The name of this method reflects the fact that it is most often
987  /// applied to sets of @ref vdof "Vector Dofs" but it would work equally
988  /// well on sets of @ref ldof "Local Dofs".
989  static void AdjustVDofs(Array<int> &vdofs);
990 
991  /// Helper to encode a sign flip into a DOF index (for Hcurl/Hdiv shapes).
992  static inline int EncodeDof(int entity_base, int idx)
993  { return (idx >= 0) ? (entity_base + idx) : (-1-(entity_base + (-1-idx))); }
994 
995  /// Helper to return the DOF associated with a sign encoded DOF
996  static inline int DecodeDof(int dof)
997  { return (dof >= 0) ? dof : (-1 - dof); }
998 
999  /// Helper to determine the DOF and sign of a sign encoded DOF
1000  static inline int DecodeDof(int dof, double& sign)
1001  { return (dof >= 0) ? (sign = 1, dof) : (sign = -1, (-1 - dof)); }
1002 
1003  /// @anchor getvdof @name Local Vector DoF Access Members
1004  /// These member functions produce arrays of local vector degree of freedom
1005  /// indices, see @ref ldof and @ref vdof. These indices can be used to
1006  /// access entries in GridFunction, LinearForm, and BilinearForm objects
1007  /// regardless of the value of @b vdim.
1008  /// @{
1009 
1010  /// @brief Returns indices of degrees of freedom for the @a i'th element.
1011  /// The returned indices are offsets into an @ref ldof vector with @b vdim
1012  /// not necessarily equal to 1. The returned indices are always ordered
1013  /// byNODES, irrespective of whether the space is byNODES or byVDIM.
1014  /// See also GetElementDofs().
1015  ///
1016  /// @note In many cases the returned DofTransformation object will be NULL.
1017  /// In other cases see the documentation of the DofTransformation class for
1018  /// guidance on its role in performing @ref edof to @ref ldof transformations
1019  /// on local vectors and matrices. At present the DofTransformation is only
1020  /// needed for Nedelec basis functions of order 2 and above on 3D elements
1021  /// with triangular faces.
1022  ///
1023  /// @note The returned object should NOT be deleted by the caller.
1024  DofTransformation *GetElementVDofs(int i, Array<int> &vdofs) const;
1025 
1026  /// @brief Returns indices of degrees of freedom for @a i'th boundary
1027  /// element.
1028  /// The returned indices are offsets into an @ref ldof vector with @b vdim
1029  /// not necessarily equal to 1. See also GetBdrElementDofs().
1030  ///
1031  /// @note In many cases the returned DofTransformation object will be NULL.
1032  /// In other cases see the documentation of the DofTransformation class for
1033  /// guidance on its role in performing @ref edof to @ref ldof transformations
1034  /// on local vectors and matrices. At present the DofTransformation is only
1035  /// needed for Nedelec basis functions of order 2 and above on 3D elements
1036  /// with triangular faces.
1037  ///
1038  /// @note The returned object should NOT be deleted by the caller.
1039  DofTransformation *GetBdrElementVDofs(int i, Array<int> &vdofs) const;
1040 
1041  /// Returns indices of degrees of freedom in @a vdofs for NURBS patch @a i.
1042  void GetPatchVDofs(int i, Array<int> &vdofs) const;
1043 
1044  /// @brief Returns the indices of the degrees of freedom for the specified
1045  /// face, including the DOFs for the edges and the vertices of the face.
1046  ///
1047  /// The returned indices are offsets into an @ref ldof vector with @b vdim
1048  /// not necessarily equal to 1. See GetFaceDofs() for more information.
1049  void GetFaceVDofs(int i, Array<int> &vdofs) const;
1050 
1051  /// @brief Returns the indices of the degrees of freedom for the specified
1052  /// edge, including the DOFs for the vertices of the edge.
1053  ///
1054  /// The returned indices are offsets into an @ref ldof vector with @b vdim
1055  /// not necessarily equal to 1. See GetEdgeDofs() for more information.
1056  void GetEdgeVDofs(int i, Array<int> &vdofs) const;
1057 
1058  /// @brief Returns the indices of the degrees of freedom for the specified
1059  /// vertices.
1060  ///
1061  /// The returned indices are offsets into an @ref ldof vector with @b vdim
1062  /// not necessarily equal to 1. See also GetVertexDofs().
1063  void GetVertexVDofs(int i, Array<int> &vdofs) const;
1064 
1065  /// @brief Returns the indices of the degrees of freedom for the interior
1066  /// of the specified element.
1067  ///
1068  /// The returned indices are offsets into an @ref ldof vector with @b vdim
1069  /// not necessarily equal to 1. See GetElementInteriorDofs() for more
1070  /// information.
1071  void GetElementInteriorVDofs(int i, Array<int> &vdofs) const;
1072 
1073  /// @brief Returns the indices of the degrees of freedom for the interior
1074  /// of the specified edge.
1075  ///
1076  /// The returned indices are offsets into an @ref ldof vector with @b vdim
1077  /// not necessarily equal to 1. See also GetEdgeInteriorDofs().
1078  void GetEdgeInteriorVDofs(int i, Array<int> &vdofs) const;
1079  /// @}
1080 
1081  /// (@deprecated) Use the Update() method if the space or mesh changed.
1082  MFEM_DEPRECATED void RebuildElementToDofTable();
1083 
1084  /** @brief Reorder the scalar DOFs based on the element ordering.
1085 
1086  The new ordering is constructed as follows: 1) loop over all elements as
1087  ordered in the Mesh; 2) for each element, assign new indices to all of
1088  its current DOFs that are still unassigned; the new indices we assign are
1089  simply the sequence `0,1,2,...`; if there are any signed DOFs their sign
1090  is preserved. */
1091  void ReorderElementToDofTable();
1092 
1094 
1095  /** @brief Return a reference to the internal Table that stores the lists of
1096  scalar dofs, for each mesh element, as returned by GetElementDofs(). */
1097  const Table &GetElementToDofTable() const { return *elem_dof; }
1098 
1099  /** @brief Return a reference to the internal Table that stores the lists of
1100  scalar dofs, for each boundary mesh element, as returned by
1101  GetBdrElementDofs(). */
1103  { if (!bdr_elem_dof) { BuildBdrElementToDofTable(); } return *bdr_elem_dof; }
1104 
1105  /** @brief Return a reference to the internal Table that stores the lists of
1106  scalar dofs, for each face in the mesh, as returned by GetFaceDofs(). In
1107  this context, "face" refers to a (dim-1)-dimensional mesh entity. */
1108  /** @note In the case of a NURBS space, the rows corresponding to interior
1109  faces will be empty. */
1110  const Table &GetFaceToDofTable() const
1111  { if (!face_dof) { BuildFaceToDofTable(); } return *face_dof; }
1112 
1113  /** @brief Initialize internal data that enables the use of the methods
1114  GetElementForDof() and GetLocalDofForDof(). */
1115  void BuildDofToArrays();
1116 
1117  /// Return the index of the first element that contains dof @a i.
1118  /** This method can be called only after setup is performed using the method
1119  BuildDofToArrays(). */
1120  int GetElementForDof(int i) const { return dof_elem_array[i]; }
1121  /// Return the local dof index in the first element that contains dof @a i.
1122  /** This method can be called only after setup is performed using the method
1123  BuildDofToArrays(). */
1124  int GetLocalDofForDof(int i) const { return dof_ldof_array[i]; }
1125 
1126  /** @brief Returns pointer to the FiniteElement in the FiniteElementCollection
1127  associated with i'th element in the mesh object.
1128  Note: The method has been updated to abort instead of returning NULL for
1129  an empty partition. */
1130  virtual const FiniteElement *GetFE(int i) const;
1131 
1132  /** @brief Returns pointer to the FiniteElement in the FiniteElementCollection
1133  associated with i'th boundary face in the mesh object. */
1134  const FiniteElement *GetBE(int i) const;
1135 
1136  /** @brief Returns pointer to the FiniteElement in the FiniteElementCollection
1137  associated with i'th face in the mesh object. Faces in this case refer
1138  to the MESHDIM-1 primitive so in 2D they are segments and in 1D they are
1139  points.*/
1140  const FiniteElement *GetFaceElement(int i) const;
1141 
1142  /** @brief Returns pointer to the FiniteElement in the FiniteElementCollection
1143  associated with i'th edge in the mesh object. */
1144  const FiniteElement *GetEdgeElement(int i, int variant = 0) const;
1145 
1146  /// Return the trace element from element 'i' to the given 'geom_type'
1147  const FiniteElement *GetTraceElement(int i, Geometry::Type geom_type) const;
1148 
1149  /** @brief Mark degrees of freedom associated with boundary elements with
1150  the specified boundary attributes (marked in 'bdr_attr_is_ess').
1151  For spaces with 'vdim' > 1, the 'component' parameter can be used
1152  to restricts the marked vDOFs to the specified component. */
1153  virtual void GetEssentialVDofs(const Array<int> &bdr_attr_is_ess,
1154  Array<int> &ess_vdofs,
1155  int component = -1) const;
1156 
1157  /** @brief Get a list of essential true dofs, ess_tdof_list, corresponding to the
1158  boundary attributes marked in the array bdr_attr_is_ess.
1159  For spaces with 'vdim' > 1, the 'component' parameter can be used
1160  to restricts the marked tDOFs to the specified component. */
1161  virtual void GetEssentialTrueDofs(const Array<int> &bdr_attr_is_ess,
1162  Array<int> &ess_tdof_list,
1163  int component = -1);
1164 
1165  /** @brief Get a list of all boundary true dofs, @a boundary_dofs. For spaces
1166  with 'vdim' > 1, the 'component' parameter can be used to restricts the
1167  marked tDOFs to the specified component. Equivalent to
1168  FiniteElementSpace::GetEssentialTrueDofs with all boundary attributes
1169  marked as essential. */
1170  void GetBoundaryTrueDofs(Array<int> &boundary_dofs, int component = -1);
1171 
1172  /// Convert a Boolean marker array to a list containing all marked indices.
1173  static void MarkerToList(const Array<int> &marker, Array<int> &list);
1174 
1175  /** @brief Convert an array of indices (list) to a Boolean marker array where all
1176  indices in the list are marked with the given value and the rest are set
1177  to zero. */
1178  static void ListToMarker(const Array<int> &list, int marker_size,
1179  Array<int> &marker, int mark_val = -1);
1180 
1181  /** @brief For a partially conforming FE space, convert a marker array (nonzero
1182  entries are true) on the partially conforming dofs to a marker array on
1183  the conforming dofs. A conforming dofs is marked iff at least one of its
1184  dependent dofs is marked. */
1185  void ConvertToConformingVDofs(const Array<int> &dofs, Array<int> &cdofs);
1186 
1187  /** @brief For a partially conforming FE space, convert a marker array (nonzero
1188  entries are true) on the conforming dofs to a marker array on the
1189  (partially conforming) dofs. A dof is marked iff it depends on a marked
1190  conforming dofs, where dependency is defined by the ConformingRestriction
1191  matrix; in other words, a dof is marked iff it corresponds to a marked
1192  conforming dof. */
1193  void ConvertFromConformingVDofs(const Array<int> &cdofs, Array<int> &dofs);
1194 
1195  /** @brief Generate the global restriction matrix from a discontinuous
1196  FE space to the continuous FE space of the same polynomial degree. */
1198 
1199  /** @brief Generate the global restriction matrix from a discontinuous
1200  FE space to the piecewise constant FE space. */
1202 
1203  /** @brief Construct the restriction matrix from the FE space given by
1204  (*this) to the lower degree FE space given by (*lfes) which
1205  is defined on the same mesh. */
1207 
1208  /** @brief Construct and return an Operator that can be used to transfer
1209  GridFunction data from @a coarse_fes, defined on a coarse mesh, to @a
1210  this FE space, defined on a refined mesh. */
1211  /** It is assumed that the mesh of this FE space is a refinement of the mesh
1212  of @a coarse_fes and the CoarseFineTransformations returned by the method
1213  Mesh::GetRefinementTransforms() of the refined mesh are set accordingly.
1214  The Operator::Type of @a T can be set to request an Operator of the set
1215  type. Currently, only Operator::MFEM_SPARSEMAT and Operator::ANY_TYPE
1216  (matrix-free) are supported. When Operator::ANY_TYPE is requested, the
1217  choice of the particular Operator sub-class is left to the method. This
1218  method also works in parallel because the transfer operator is local to
1219  the MPI task when the input is a synchronized ParGridFunction. */
1220  void GetTransferOperator(const FiniteElementSpace &coarse_fes,
1221  OperatorHandle &T) const;
1222 
1223  /** @brief Construct and return an Operator that can be used to transfer
1224  true-dof data from @a coarse_fes, defined on a coarse mesh, to @a this FE
1225  space, defined on a refined mesh.
1226 
1227  This method calls GetTransferOperator() and multiplies the result by the
1228  prolongation operator of @a coarse_fes on the right, and by the
1229  restriction operator of this FE space on the left.
1230 
1231  The Operator::Type of @a T can be set to request an Operator of the set
1232  type. In serial, the supported types are: Operator::MFEM_SPARSEMAT and
1233  Operator::ANY_TYPE (matrix-free). In parallel, the supported types are:
1234  Operator::Hypre_ParCSR and Operator::ANY_TYPE. Any other type is treated
1235  as Operator::ANY_TYPE: the operator representation choice is made by this
1236  method. */
1237  virtual void GetTrueTransferOperator(const FiniteElementSpace &coarse_fes,
1238  OperatorHandle &T) const;
1239 
1240  /** @brief Reflect changes in the mesh: update number of DOFs, etc. Also, calculate
1241  GridFunction transformation operator (unless want_transform is false).
1242  Safe to call multiple times, does nothing if space already up to date. */
1243  virtual void Update(bool want_transform = true);
1244 
1245  /// Get the GridFunction update operator.
1246  const Operator* GetUpdateOperator() { Update(); return Th.Ptr(); }
1247 
1248  /// Return the update operator in the given OperatorHandle, @a T.
1250 
1251  /** @brief Set the ownership of the update operator: if set to false, the
1252  Operator returned by GetUpdateOperator() must be deleted outside the
1253  FiniteElementSpace. */
1254  /** The update operator ownership is automatically reset to true when a new
1255  update operator is created by the Update() method. */
1256  void SetUpdateOperatorOwner(bool own) { Th.SetOperatorOwner(own); }
1257 
1258  /// Specify the Operator::Type to be used by the update operators.
1259  /** The default type is Operator::ANY_TYPE which leaves the choice to this
1260  class. The other currently supported option is Operator::MFEM_SPARSEMAT
1261  which is only guaranteed to be honored for a refinement update operator.
1262  Any other type will be treated as Operator::ANY_TYPE.
1263  @note This operation destroys the current update operator (if owned). */
1265 
1266  /// Free the GridFunction update operator (if any), to save memory.
1267  virtual void UpdatesFinished() { Th.Clear(); }
1268 
1269  /** Return update counter, similar to Mesh::GetSequence(). Used by
1270  GridFunction to check if it is up to date with the space. */
1271  long GetSequence() const { return sequence; }
1272 
1273  /// Return whether or not the space is discontinuous (L2)
1274  bool IsDGSpace() const
1275  {
1276  return dynamic_cast<const L2_FECollection*>(fec) != NULL;
1277  }
1278 
1279  /** In variable order spaces on nonconforming (NC) meshes, this function
1280  controls whether strict conformity is enforced in cases where coarse
1281  edges/faces have higher polynomial order than their fine NC neighbors.
1282  In the default (strict) case, the coarse side polynomial order is
1283  reduced to that of the lowest order fine edge/face, so all fine
1284  neighbors can interpolate the coarse side exactly. If relaxed == true,
1285  some discontinuities in the solution in such cases are allowed and the
1286  coarse side is not restricted. For an example, see
1287  https://github.com/mfem/mfem/pull/1423#issuecomment-621340392 */
1288  void SetRelaxedHpConformity(bool relaxed = true)
1289  {
1290  relaxed_hp = relaxed;
1291  orders_changed = true; // force update
1292  Update(false);
1293  }
1294 
1295  /// Save finite element space to output stream @a out.
1296  void Save(std::ostream &out) const;
1297 
1298  /** @brief Read a FiniteElementSpace from a stream. The returned
1299  FiniteElementCollection is owned by the caller. */
1300  FiniteElementCollection *Load(Mesh *m, std::istream &input);
1301 
1302  virtual ~FiniteElementSpace();
1303 };
1304 
1305 /// @brief Return true if the mesh contains only one topology and the elements are tensor elements.
1306 inline bool UsesTensorBasis(const FiniteElementSpace& fes)
1307 {
1308  Mesh & mesh = *fes.GetMesh();
1309  const bool mixed = mesh.GetNumGeometries(mesh.Dimension()) > 1;
1310  // Potential issue: empty local mesh --> no element 0.
1311  return !mixed &&
1312  dynamic_cast<const mfem::TensorBasisElement *>(fes.GetFE(0))!=nullptr;
1313 }
1314 
1315 }
1316 
1317 #endif
SparseMatrix * RefinementMatrix(int old_ndofs, const Table *old_elem_dof, const Table *old_elem_fos)
Definition: fespace.cpp:1517
Abstract class for all finite elements.
Definition: fe_base.hpp:233
VDofTransformation VDoFTrans
Definition: fespace.hpp:275
int GetBdrElementType(int i) const
Returns the type of boundary element i.
Definition: fespace.hpp:765
void GetEdgeInteriorDofs(int i, Array< int > &dofs) const
Returns the indices of the degrees of freedom for the interior of the specified edge.
Definition: fespace.cpp:3130
Abstract base class for LORDiscretization and ParLORDiscretization classes, which construct low-order...
Definition: lor.hpp:22
int ndofs
Number of degrees of freedom. Number of unknowns is ndofs * vdim.
Definition: fespace.hpp:242
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:96
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
Definition: fespace.hpp:577
int * bdofs
internal DOFs of elements if mixed/var-order; NULL otherwise
Definition: fespace.hpp:250
static void AddDependencies(SparseMatrix &deps, Array< int > &master_dofs, Array< int > &slave_dofs, DenseMatrix &I, int skipfirst=0)
Definition: fespace.cpp:781
int GetConformingVSize() const
Definition: fespace.hpp:718
virtual void Update(bool want_transform=true)
Reflect changes in the mesh: update number of DOFs, etc. Also, calculate GridFunction transformation ...
Definition: fespace.cpp:3402
static int Map(int ndofs, int vdim, int dof, int vd)
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:2085
int MakeDofTable(int ent_dim, const Array< int > &entity_orders, Table &entity_dofs, Array< char > *var_ent_order)
Definition: fespace.cpp:2619
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:1102
std::unique_ptr< SparseMatrix > cP
Definition: fespace.hpp:280
DerefinementOperator(const FiniteElementSpace *f_fes, const FiniteElementSpace *c_fes, BilinearFormIntegrator *mass_integ)
TODO: Implement DofTransformation support.
Definition: fespace.cpp:1926
static const int NumGeom
Definition: geom.hpp:42
bool Nonconforming() const
Definition: fespace.hpp:562
int Dimension() const
Dimension of the reference space used within the elements.
Definition: mesh.hpp:1020
Matrix-free transfer operator between finite element spaces on the same mesh.
Definition: transfer.hpp:428
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition: mesh.cpp:5668
Ordering::Type ordering
Definition: fespace.hpp:239
static int DecodeDof(int dof, double &sign)
Helper to determine the DOF and sign of a sign encoded DOF.
Definition: fespace.hpp:1000
int GetElementType(int i) const
Returns the type of element i.
Definition: fespace.hpp:757
void GetEdgeInteriorVDofs(int i, Array< int > &vdofs) const
Returns the indices of the degrees of freedom for the interior of the specified edge.
Definition: fespace.cpp:342
void GetVertexVDofs(int i, Array< int > &vdofs) const
Returns the indices of the degrees of freedom for the specified vertices.
Definition: fespace.cpp:330
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)
Helper to return the DOF associated with a sign encoded DOF.
Definition: fespace.hpp:996
const SparseMatrix * GetConformingRestriction() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:1275
void GetFaceVDofs(int i, Array< int > &vdofs) const
Returns the indices of the degrees of freedom for the specified face, including the DOFs for the edge...
Definition: fespace.cpp:318
Element::Type GetBdrElementType(int i) const
Returns the type of boundary element i.
Definition: mesh.cpp:6649
static int MinOrder(VarOrderBits bits)
Return the minimum order (least significant bit set) in the bit mask.
Definition: fespace.cpp:2492
int GetOrder() const
Return the order (polynomial degree) of the FE collection, corresponding to the order/degree returned...
Definition: fe_coll.hpp:225
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
virtual DofTransformation * GetBdrElementDofs(int bel, Array< int > &dofs) const
Returns indices of degrees of freedom for boundary element &#39;bel&#39;. The returned indices are offsets in...
Definition: fespace.cpp:2878
static constexpr int MaxVarOrder
Definition: fespace.hpp:346
std::unique_ptr< SparseMatrix > cR_hp
A version of the conforming restriction matrix for variable-order spaces.
Definition: fespace.hpp:284
virtual const FaceRestriction * GetFaceRestriction(ElementDofOrdering f_ordering, FaceType, L2FaceValues mul=L2FaceValues::DoubleValued) const
Return an Operator that converts L-vectors to E-vectors on each face.
Definition: fespace.cpp:1335
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:806
bool Nonconforming() const
Definition: mesh.hpp:1969
virtual const Operator * GetRestrictionOperator() const
An abstract operator that performs the same action as GetRestrictionMatrix.
Definition: fespace.hpp:612
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition: array.cpp:68
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
const NURBSExtension * GetNURBSext() const
Definition: fespace.hpp:557
int vdim
Vector dimension (number of unknowns per degree of freedom).
Definition: fespace.hpp:234
int GetDegenerateFaceDofs(int index, Array< int > &dofs, Geometry::Type master_geom, int variant) const
Definition: fespace.cpp:875
OperatorHandle L2E_lex
Definition: fespace.hpp:293
int GetNDofs() const
Returns number of degrees of freedom. This is the number of Local Degrees of Freedom.
Definition: fespace.hpp:706
void BuildBdrElementToDofTable() const
Definition: fespace.cpp:389
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:587
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:3239
const Operator * GetRestrictionTransposeOperator() const
Return an operator that performs the transpose of GetRestrictionOperator.
Definition: fespace.cpp:1290
MFEM_DEPRECATED void RebuildElementToDofTable()
(
Definition: fespace.cpp:454
void CalcEdgeFaceVarOrders(Array< VarOrderBits > &edge_orders, Array< VarOrderBits > &face_orders) const
Definition: fespace.cpp:2503
void BuildConformingInterpolation() const
Calculate the cP and cR matrices for a nonconforming mesh.
Definition: fespace.cpp:951
FiniteElementCollection * Load(Mesh *m, std::istream &input)
Read a FiniteElementSpace from a stream. The returned FiniteElementCollection is owned by the caller...
Definition: fespace.cpp:3581
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
bool UsesTensorBasis(const FiniteElementSpace &fes)
Return true if the mesh contains only one topology and the elements are tensor elements.
Definition: fespace.hpp:1306
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: fespace.cpp:1633
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:231
OperatorHandle Th
Transformation to apply to GridFunctions after space Update().
Definition: fespace.hpp:290
void GetVertexDofs(int i, Array< int > &dofs) const
Returns the indices of the degrees of freedom for the specified vertices.
Definition: fespace.cpp:3099
NURBSExtension * GetNURBSext()
Definition: fespace.hpp:558
const SparseMatrix * GetHpConformingRestriction() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:1283
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:1089
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:2841
std::unordered_map< const key_face, FaceRestriction *, key_hash > map_L2F
Definition: fespace.hpp:306
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
Definition: mesh.cpp:6274
std::uint64_t VarOrderBits
Bit-mask representing a set of orders needed by an edge/face.
Definition: fespace.hpp:345
int FindFaceDof(int face, int ndof) const
Similar to FindEdgeDof, but used for mixed meshes too.
Definition: fespace.hpp:373
Array< int > dof_elem_array
Definition: fespace.hpp:268
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:621
Derefinement operator, used by the friend class InterpolationGridTransfer.
Definition: fespace.hpp:438
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:3199
const ElementRestrictionOperator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
Definition: fespace.cpp:1302
int GetNFbyType(FaceType type) const
Returns the number of faces according to the requested type.
Definition: fespace.hpp:753
ElementTransformation * GetBdrElementTransformation(int i)
Definition: mesh.cpp:438
virtual void UpdateMeshPointer(Mesh *new_mesh)
Definition: fespace.cpp:3507
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:661
void GetElementInteriorVDofs(int i, Array< int > &vdofs) const
Returns the indices of the degrees of freedom for the interior of the specified element.
Definition: fespace.cpp:336
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:1083
void GetPatchVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom in vdofs for NURBS patch i.
Definition: fespace.cpp:312
int uni_fdof
of single face DOFs if all faces uniform; -1 otherwise
Definition: fespace.hpp:249
int GetAttribute(int i) const
Return the attribute of element i.
Definition: mesh.hpp:1190
void BuildFaceToDofTable() const
Definition: fespace.cpp:428
const Table * GetElementToFaceOrientationTable() const
Definition: fespace.hpp:1093
int GetNEDofs() const
Number of all scalar edge-interior dofs.
Definition: fespace.hpp:728
FaceType
Definition: mesh.hpp:45
void BuildElementToDofTable() const
Definition: fespace.cpp:348
Array< DofTransformation * > DoFTrans
Definition: fespace.hpp:274
int GetMaxElementOrder() const
Return the maximum polynomial order.
Definition: fespace.hpp:573
bool orders_changed
True if at least one element order changed (variable-order space only).
Definition: fespace.hpp:322
Data type sparse matrix.
Definition: sparsemat.hpp:50
Native ordering as defined by the FiniteElement.
std::unique_ptr< Operator > R_transpose
Operator computing the action of the transpose of the restriction.
Definition: fespace.hpp:287
Array< char > elem_order
Definition: fespace.hpp:246
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:723
ElementTransformation * GetBdrElementTransformation(int i) const
Returns ElementTransformation for the i-th boundary element.
Definition: fespace.hpp:778
A class that performs interpolation from an E-vector to quadrature point values and/or derivatives (Q...
virtual const SparseMatrix * GetRestrictionMatrix() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.hpp:616
Type
Ordering methods:
Definition: fespace.hpp:33
SparseMatrix * D2Const_GlobalRestrictionMatrix(FiniteElementSpace *cfes)
Generate the global restriction matrix from a discontinuous FE space to the piecewise constant FE spa...
Definition: fespace.cpp:702
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:670
bool Conforming() const
Definition: mesh.hpp:1968
int GetNumBorderDofs(Geometry::Type geom, int order) const
Definition: fespace.cpp:917
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; > 1, the &#39;component&#39; para...
Definition: fespace.cpp:605
NURBSExtension * StealNURBSext()
Definition: fespace.cpp:2281
long GetSequence() const
Definition: fespace.hpp:1271
void GetElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of element i.
Definition: mesh.hpp:1293
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:1399
const SparseMatrix * GetConformingProlongation() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:1268
FiniteElementSpace & operator=(const FiniteElementSpace &)=delete
Copy assignment not supported.
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Returns the transformation defining the i-th element in the user-defined variable ElTr...
Definition: fespace.hpp:774
Array< int > dof_ldof_array
Definition: fespace.hpp:268
int FindEdgeDof(int edge, int ndof) const
Definition: fespace.hpp:369
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:992
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
Definition: table.cpp:187
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:736
void GetFaceInteriorDofs(int i, Array< int > &dofs) const
Returns the indices of the degrees of freedom for the interior of the specified face.
Definition: fespace.cpp:3142
Array< QuadratureInterpolator * > E2Q_array
Definition: fespace.hpp:309
A class that performs interpolation from a face E-vector to quadrature point values and/or derivative...
Abstract base class that defines an interface for element restrictions.
Definition: restriction.hpp:25
OperatorHandle L2E_nat
The element restriction operators, see GetElementRestriction().
Definition: fespace.hpp:293
static void DofsToVDofs(int ndofs, int vdim, Array< int > &dofs)
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:702
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for the i&#39;th element. The returned indices are offsets into an ...
Definition: fespace.cpp:281
int GetNF() const
Returns number of faces (i.e. co-dimension 1 entities) in the mesh.
Definition: fespace.hpp:742
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
Definition: mesh.hpp:1196
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:712
Type
Enumeration defining IDs for some classes derived from Operator.
Definition: operator.hpp:283
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:555
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:24
std::size_t operator()(const key_face &k) const
Definition: fespace.hpp:298
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:219
SparseMatrix * RefinementMatrix_main(const int coarse_ndofs, const Table &coarse_elem_dof, const Table *coarse_elem_fos, const DenseTensor localP[]) const
Definition: fespace.cpp:1430
Transfer data between a coarse mesh and an embedded refined mesh using interpolation.
Definition: transfer.hpp:123
virtual ~FiniteElementSpace()
Definition: fespace.cpp:3244
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:733
int GetNumElementInteriorDofs(int i) const
Returns the number of degrees of freedom associated with the interior of the specified element...
Definition: fespace.cpp:3124
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
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:518
void GetPatchDofs(int patch, Array< int > &dofs) const
Returns indices of degrees of freedom for NURBS patch index patch. Cartesian ordering is used...
Definition: fespace.cpp:2834
int GetOrder(int i) const
Returns the polynomial degree of the i&#39;th finite element.
Definition: fespace.hpp:692
int GetElementForDof(int i) const
Return the index of the first element that contains dof i.
Definition: fespace.hpp:1120
int GetEdgeOrder(int edge, int variant=0) const
Definition: fespace.cpp:2691
void GetUpdateOperator(OperatorHandle &T)
Return the update operator in the given OperatorHandle, T.
Definition: fespace.hpp:1249
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:3051
RefinementOperator(const FiniteElementSpace *fespace, Table *old_elem_dof, Table *old_elem_fos, int old_ndofs)
Definition: fespace.cpp:1537
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:2674
Array< char > var_face_orders
Definition: fespace.hpp:259
int GetFaceOrder(int face, int variant=0) const
Returns the polynomial degree of the i&#39;th face finite element.
Definition: fespace.cpp:2702
void Clear()
Clear the OperatorHandle, deleting the held Operator (if owned), while leaving the type id unchanged...
Definition: handle.hpp:124
void GetElementVertices(int i, Array< int > &vertices) const
Returns the vertices of element i.
Definition: fespace.hpp:761
virtual DofTransformation * GetElementDofs(int elem, Array< int > &dofs) const
Returns indices of degrees of freedom of element &#39;elem&#39;. The returned indices are offsets into an ldo...
Definition: fespace.cpp:2731
int GetNVDofs() const
Number of all scalar vertex dofs.
Definition: fespace.hpp:726
The ordering method used when the number of unknowns per mesh node (vector dimension) is bigger than ...
Definition: fespace.hpp:29
void SetRelaxedHpConformity(bool relaxed=true)
Definition: fespace.hpp:1288
int GetAttribute(int i) const
Definition: fespace.hpp:781
std::unique_ptr< SparseMatrix > cR
Conforming restriction matrix such that cR.cP=I.
Definition: fespace.hpp:282
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:3339
FiniteElementSpace(Mesh *mesh, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES)
Definition: fespace.hpp:537
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:1086
QVectorLayout
Type describing possible layouts for Q-vectors.
Definition: fespace.hpp:52
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: fespace.hpp:593
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:653
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
Definition: handle.hpp:120
Element::Type GetElementType(int i) const
Returns the type of element i.
Definition: mesh.cpp:6644
const QuadratureInterpolator * GetQuadratureInterpolator(const IntegrationRule &ir) const
Return a QuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivat...
Definition: fespace.cpp:1370
void Swap(Mesh &other, bool non_geometry)
Definition: mesh.cpp:9699
int GetElementOrder(int i) const
Returns the order of the i&#39;th finite element.
Definition: fespace.cpp:178
A standard isoparametric element transformation.
Definition: eltrans.hpp:361
Mesh * mesh
The mesh that FE space lives on (not owned).
Definition: fespace.hpp:228
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:355
void ReorderElementToDofTable()
Reorder the scalar DOFs based on the element ordering.
Definition: fespace.cpp:463
GridFunction interpolation operator applicable after mesh refinement.
Definition: fespace.hpp:412
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:5684
Array< FaceQuadratureInterpolator * > E2BFQ_array
Definition: fespace.hpp:311
int GetBdrAttribute(int i) const
Definition: fespace.hpp:783
int GetElementOrderImpl(int i) const
Return element order: internal version of GetElementOrder without checks.
Definition: fespace.cpp:189
Ordering::Type GetOrdering() const
Return the ordering method.
Definition: fespace.hpp:721
DofTransformation * GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for i&#39;th boundary element. The returned indices are offsets int...
Definition: fespace.cpp:297
Table var_face_dofs
NOTE: also used for spaces with mixed faces.
Definition: fespace.hpp:255
void Save(std::ostream &out) const
Save finite element space to output stream out.
Definition: fespace.cpp:3512
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition: fespace.hpp:74
Array< FaceQuadratureInterpolator * > E2IFQ_array
Definition: fespace.hpp:310
void UpdateElementOrders()
Resize the elem_order array on mesh change.
Definition: fespace.cpp:3380
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:926
void GetLocalDerefinementMatrices(Geometry::Type geom, DenseTensor &localR) const
Definition: fespace.cpp:2061
void GetEdgeVDofs(int i, Array< int > &vdofs) const
Returns the indices of the degrees of freedom for the specified edge, including the DOFs for the vert...
Definition: fespace.cpp:324
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition: fespace.hpp:745
void MakeVDimMatrix(SparseMatrix &mat) const
Replicate &#39;mat&#39; in the vector dimension, according to vdim ordering mode.
Definition: fespace.cpp:1240
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:1246
int FirstFaceDof(int face, int variant=0) const
Definition: fespace.hpp:376
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:709
Array< char > var_edge_orders
Definition: fespace.hpp:259
Lexicographic ordering for tensor-product FiniteElements.
FiniteElementSpace()
Default constructor: the object is invalid until initialized using the method Load().
Definition: fespace.cpp:59
void Constructor(Mesh *mesh, NURBSExtension *ext, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES)
Help function for constructors + Load().
Definition: fespace.cpp:2188
Operator * Ptr() const
Access the underlying Operator pointer.
Definition: handle.hpp:87
void SetUpdateOperatorOwner(bool own)
Set the ownership of the update operator: if set to false, the Operator returned by GetUpdateOperator...
Definition: fespace.hpp:1256
int Size() const
Return the logical size of the array.
Definition: array.hpp:141
NQPT x VDIM x NE (values) / NQPT x VDIM x DIM x NE (grads)
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:640
void BuildNURBSFaceToDofTable() const
Generates partial face_dof table for a NURBS space.
Definition: fespace.cpp:2316
int DofToVDof(int dof, int vd, int ndofs=-1) const
Compute a single vdof corresponding to the index dof and the vector index vd.
Definition: fespace.cpp:251
int GetNV() const
Returns number of vertices in the mesh.
Definition: fespace.hpp:733
Vector data type.
Definition: vector.hpp:58
std::tuple< bool, ElementDofOrdering, FaceType, L2FaceValues > key_face
The face restriction operators, see GetFaceRestriction().
Definition: fespace.hpp:295
bool IsDGSpace() const
Return whether or not the space is discontinuous (L2)
Definition: fespace.hpp:1274
Array< int > face_to_be
Definition: fespace.hpp:272
virtual void UpdatesFinished()
Free the GridFunction update operator (if any), to save memory.
Definition: fespace.hpp:1267
void GetElementInteriorDofs(int i, Array< int > &dofs) const
Returns the indices of the degrees of freedom for the interior of the specified element.
Definition: fespace.cpp:3109
void SetUpdateOperatorType(Operator::Type tid)
Specify the Operator::Type to be used by the update operators.
Definition: fespace.hpp:1264
static bool DofFinalizable(int dof, const Array< bool > &finalized, const SparseMatrix &deps)
Definition: fespace.cpp:861
Class representing the storage layout of a QuadratureFunction.
Definition: qspace.hpp:101
NURBSExtension * NURBSext
Definition: fespace.hpp:270
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:3229
int GetNVariants(int entity, int index) const
Return number of possible DOF variants for edge/face (var. order spaces).
Definition: fespace.cpp:2718
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: fespace.cpp:2028
Base class for operators that extracts Face degrees of freedom.
Abstract operator.
Definition: operator.hpp:24
int GetNConformingDofs() const
Definition: fespace.cpp:1296
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:2970
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:1096
void BuildDofToArrays()
Initialize internal data that enables the use of the methods GetElementForDof() and GetLocalDofForDof...
Definition: fespace.cpp:483
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:1097
virtual const SparseMatrix * GetHpRestrictionMatrix() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.hpp:620
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th boundary fac...
Definition: fespace.cpp:3166
int GetNFDofs() const
Number of all scalar face-interior dofs.
Definition: fespace.hpp:730
void GetLocalRefinementMatrices(Geometry::Type geom, DenseTensor &localP) const
Definition: fespace.cpp:1494
int GetLocalDofForDof(int i) const
Return the local dof index in the first element that contains dof i.
Definition: fespace.hpp:1124
bool Conforming() const
Definition: fespace.hpp:561
void DofsToVDofs(Array< int > &dofs, int ndofs=-1) const
Compute the full set of vdofs corresponding to each entry in dofs.
Definition: fespace.cpp:215
static void AdjustVDofs(Array< int > &vdofs)
Remove the orientation information encoded into an array of dofs Some basis function types have a rel...
Definition: fespace.cpp:267
void SetType(Operator::Type tid)
Invoke Clear() and set a new type id.
Definition: handle.hpp:132
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:1110
int VDofToDof(int vdof) const
Compute the inverse of the Dof to VDof mapping for a single index vdof.
Definition: fespace.hpp:974
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:3312
VDIM x NQPT x NE (values) / VDIM x DIM x NQPT x NE (grads)
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:1706
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition: fe_coll.hpp:327
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition: fespace.hpp:769
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:546