MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
fespace.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_FESPACE
13 #define MFEM_FESPACE
14 
15 #include "../config/config.hpp"
16 #include "../linalg/sparsemat.hpp"
17 #include "../mesh/mesh.hpp"
18 #include "fe_coll.hpp"
19 #include "restriction.hpp"
20 #include <iostream>
21 #include <unordered_map>
22 
23 namespace mfem
24 {
25 
26 /** @brief The ordering method used when the number of unknowns per mesh node
27  (vector dimension) is bigger than 1. */
28 class Ordering
29 {
30 public:
31  /// %Ordering methods:
32  enum Type
33  {
34  byNODES, /**< loop first over the nodes (inner loop) then over the vector
35  dimension (outer loop); symbolically it can be represented
36  as: XXX...,YYY...,ZZZ... */
37  byVDIM /**< loop first over the vector dimension (inner loop) then over
38  the nodes (outer loop); symbolically it can be represented
39  as: XYZ,XYZ,XYZ,... */
40  };
41 
42  template <Type Ord>
43  static inline int Map(int ndofs, int vdim, int dof, int vd);
44 
45  template <Type Ord>
46  static void DofsToVDofs(int ndofs, int vdim, Array<int> &dofs);
47 };
48 
49 template <> inline int
50 Ordering::Map<Ordering::byNODES>(int ndofs, int vdim, int dof, int vd)
51 {
52  MFEM_ASSERT(dof < ndofs && -1-dof < ndofs && 0 <= vd && vd < vdim, "");
53  return (dof >= 0) ? dof+ndofs*vd : dof-ndofs*vd;
54 }
55 
56 template <> inline int
57 Ordering::Map<Ordering::byVDIM>(int ndofs, int vdim, int dof, int vd)
58 {
59  MFEM_ASSERT(dof < ndofs && -1-dof < ndofs && 0 <= vd && vd < vdim, "");
60  return (dof >= 0) ? vd+vdim*dof : -1-(vd+vdim*(-1-dof));
61 }
62 
63 
64 /// Constants describing the possible orderings of the DOFs in one element.
66 {
67  /// Native ordering as defined by the FiniteElement.
68  /** This ordering can be used by tensor-product elements when the
69  interpolation from the DOFs to quadrature points does not use the
70  tensor-product structure. */
71  NATIVE,
72  /// Lexicographic ordering for tensor-product FiniteElements.
73  /** This ordering can be used only with tensor-product elements. */
75 };
76 
77 // Forward declarations
78 class NURBSExtension;
79 class BilinearFormIntegrator;
80 class QuadratureSpace;
81 class QuadratureInterpolator;
82 class FaceQuadratureInterpolator;
83 
84 
85 /** @brief Class FiniteElementSpace - responsible for providing FEM view of the
86  mesh, mainly managing the set of degrees of freedom. */
88 {
91 
92 protected:
93  /// The mesh that FE space lives on (not owned).
95 
96  /// Associated FE collection (not owned).
98 
99  /// %Vector dimension (number of unknowns per degree of freedom).
100  int vdim;
101 
102  /** Type of ordering of the vector dofs when #vdim > 1.
103  - Ordering::byNODES - first nodes, then vector dimension,
104  - Ordering::byVDIM - first vector dimension, then nodes */
106 
107  /// Number of degrees of freedom. Number of unknowns is #ndofs * #vdim.
108  int ndofs;
109 
111  int *fdofs, *bdofs;
112 
113  mutable Table *elem_dof; // if NURBS FE space, not owned; otherwise, owned.
114  mutable Table *bdrElem_dof; // not owned only if NURBS FE space.
115  mutable Table *face_dof; // owned
116  mutable Array<int> face_to_be; // used only with NURBS FE spaces; owned.
117 
119 
121  int own_ext;
122 
123  /** Matrix representing the prolongation from the global conforming dofs to
124  a set of intermediate partially conforming dofs, e.g. the dofs associated
125  with a "cut" space on a non-conforming mesh. */
126  mutable SparseMatrix *cP; // owned
127  /// Conforming restriction matrix such that cR.cP=I.
128  mutable SparseMatrix *cR; // owned
129  mutable bool cP_is_set;
130 
131  /// Transformation to apply to GridFunctions after space Update().
133 
134  /// The element restriction operators, see GetElementRestriction().
136  /// The face restriction operators, see GetFaceRestriction().
137  using key_face = std::tuple<bool, ElementDofOrdering, FaceType, L2FaceValues>;
138  struct key_hash
139  {
140  std::size_t operator()(const key_face& k) const
141  {
142  return std::get<0>(k)
143  + 2 * (int)std::get<1>(k)
144  + 4 * (int)std::get<2>(k)
145  + 8 * (int)std::get<3>(k);
146  }
147  };
148  using map_L2F = std::unordered_map<const key_face,Operator*,key_hash>;
149  mutable map_L2F L2F;
150 
154 
155  long sequence; // should match Mesh::GetSequence
156 
157  void UpdateNURBS();
158 
159  void Construct();
160  void Destroy();
161 
162  void BuildElementToDofTable() const;
163  void BuildBdrElementToDofTable() const;
164  void BuildFaceToDofTable() const;
165 
166  /** @brief Generates partial face_dof table for a NURBS space.
167 
168  The table is only defined for exterior faces that coincide with a
169  boundary. */
170  void BuildNURBSFaceToDofTable() const;
171 
172  /// Helpers to remove encoded sign from a DOF
173  static inline int DecodeDof(int dof)
174  {
175  return (dof >= 0) ? dof : (-1 - dof);
176  }
177 
178  static inline int DecodeDof(int dof, double& sign)
179  { return (dof >= 0) ? (sign = 1, dof) : (sign = -1, (-1 - dof)); }
180 
181  /// Helper to get vertex, edge or face DOFs (entity=0,1,2 resp.).
182  void GetEntityDofs(int entity, int index, Array<int> &dofs,
183  Geometry::Type master_geom = Geometry::INVALID) const;
184  // Get degenerate face DOFs: see explanation in method implementation.
185  void GetDegenerateFaceDofs(int index, Array<int> &dofs,
186  Geometry::Type master_geom) const;
187 
188  /// Calculate the cP and cR matrices for a nonconforming mesh.
189  void BuildConformingInterpolation() const;
190 
191  static void AddDependencies(SparseMatrix& deps, Array<int>& master_dofs,
192  Array<int>& slave_dofs, DenseMatrix& I);
193 
194  static bool DofFinalizable(int dof, const Array<bool>& finalized,
195  const SparseMatrix& deps);
196 
197  /// Replicate 'mat' in the vector dimension, according to vdim ordering mode.
198  void MakeVDimMatrix(SparseMatrix &mat) const;
199 
200  /// GridFunction interpolation operator applicable after mesh refinement.
202  {
203  const FiniteElementSpace* fespace;
205  Table* old_elem_dof; // Owned.
206 
207  public:
208  /** Construct the operator based on the elem_dof table of the original
209  (coarse) space. The class takes ownership of the table. */
210  RefinementOperator(const FiniteElementSpace* fespace,
211  Table *old_elem_dof/*takes ownership*/, int old_ndofs);
212  RefinementOperator(const FiniteElementSpace *fespace,
213  const FiniteElementSpace *coarse_fes);
214  virtual void Mult(const Vector &x, Vector &y) const;
215  virtual void MultTranspose(const Vector &x, Vector &y) const;
216  virtual ~RefinementOperator();
217  };
218 
219  /// Derefinement operator, used by the friend class InterpolationGridTransfer.
221  {
222  const FiniteElementSpace *fine_fes; // Not owned.
224  Table *coarse_elem_dof; // Owned.
225  Table coarse_to_fine;
226  Array<int> coarse_to_ref_type;
227  Array<Geometry::Type> ref_type_to_geom;
228  Array<int> ref_type_to_fine_elem_offset;
229 
230  public:
232  const FiniteElementSpace *c_fes,
233  BilinearFormIntegrator *mass_integ);
234  virtual void Mult(const Vector &x, Vector &y) const;
235  virtual ~DerefinementOperator();
236  };
237 
238  /** This method makes the same assumptions as the method:
239  void GetLocalRefinementMatrices(
240  const FiniteElementSpace &coarse_fes, Geometry::Type geom,
241  DenseTensor &localP) const
242  which is defined below. It also assumes that the coarse fes and this have
243  the same vector dimension, vdim. */
244  SparseMatrix *RefinementMatrix_main(const int coarse_ndofs,
245  const Table &coarse_elem_dof,
246  const DenseTensor localP[]) const;
247 
249  DenseTensor &localP) const;
251  DenseTensor &localR) const;
252 
253  /** Calculate explicit GridFunction interpolation matrix (after mesh
254  refinement). NOTE: consider using the RefinementOperator class instead
255  of the fully assembled matrix, which can take a lot of memory. */
256  SparseMatrix* RefinementMatrix(int old_ndofs, const Table* old_elem_dof);
257 
258  /// Calculate GridFunction restriction matrix after mesh derefinement.
259  SparseMatrix* DerefinementMatrix(int old_ndofs, const Table* old_elem_dof);
260 
261  /** @brief Return in @a localP the local refinement matrices that map
262  between fespaces after mesh refinement. */
263  /** This method assumes that this->mesh is a refinement of coarse_fes->mesh
264  and that the CoarseFineTransformations of this->mesh are set accordingly.
265  Another assumption is that the FEs of this use the same MapType as the FEs
266  of coarse_fes. Finally, it assumes that the spaces this and coarse_fes are
267  NOT variable-order spaces. */
268  void GetLocalRefinementMatrices(const FiniteElementSpace &coarse_fes,
270  DenseTensor &localP) const;
271 
272  /// Help function for constructors + Load().
273  void Constructor(Mesh *mesh, NURBSExtension *ext,
275  int vdim = 1, int ordering = Ordering::byNODES);
276 
277 public:
278  /** @brief Default constructor: the object is invalid until initialized using
279  the method Load(). */
281 
282  /** @brief Copy constructor: deep copy all data from @a orig except the Mesh,
283  the FiniteElementCollection, ans some derived data. */
284  /** If the @a mesh or @a fec pointers are NULL (default), then the new
285  FiniteElementSpace will reuse the respective pointers from @a orig. If
286  any of these pointers is not NULL, the given pointer will be used instead
287  of the one used by @a orig.
288 
289  @note The objects pointed to by the @a mesh and @a fec parameters must be
290  either the same objects as the ones used by @a orig, or copies of them.
291  Otherwise, the behavior is undefined.
292 
293  @note Derived data objects, such as the conforming prolongation and
294  restriction matrices, and the update operator, will not be copied, even
295  if they are created in the @a orig object. */
296  FiniteElementSpace(const FiniteElementSpace &orig, Mesh *mesh = NULL,
297  const FiniteElementCollection *fec = NULL);
298 
301  int vdim = 1, int ordering = Ordering::byNODES)
302  { Constructor(mesh, NULL, fec, vdim, ordering); }
303 
304  /// Construct a NURBS FE space based on the given NURBSExtension, @a ext.
305  /** @note If the pointer @a ext is NULL, this constructor is equivalent to
306  the standard constructor with the same arguments minus the
307  NURBSExtension, @a ext. */
310  int vdim = 1, int ordering = Ordering::byNODES)
311  { Constructor(mesh, ext, fec, vdim, ordering); }
312 
313  /// Returns the mesh
314  inline Mesh *GetMesh() const { return mesh; }
315 
316  const NURBSExtension *GetNURBSext() const { return NURBSext; }
319 
320  bool Conforming() const { return mesh->Conforming(); }
321  bool Nonconforming() const { return mesh->Nonconforming(); }
322 
323  /// The returned SparseMatrix is owned by the FiniteElementSpace.
325 
326  /// The returned SparseMatrix is owned by the FiniteElementSpace.
327  const SparseMatrix *GetConformingRestriction() const;
328 
329  /// The returned Operator is owned by the FiniteElementSpace.
330  virtual const Operator *GetProlongationMatrix() const
331  { return GetConformingProlongation(); }
332 
333  /// The returned SparseMatrix is owned by the FiniteElementSpace.
334  virtual const SparseMatrix *GetRestrictionMatrix() const
335  { return GetConformingRestriction(); }
336 
337  /// Return an Operator that converts L-vectors to E-vectors.
338  /** An L-vector is a vector of size GetVSize() which is the same size as a
339  GridFunction. An E-vector represents the element-wise discontinuous
340  version of the FE space.
341 
342  The layout of the E-vector is: ND x VDIM x NE, where ND is the number of
343  degrees of freedom, VDIM is the vector dimension of the FE space, and NE
344  is the number of the mesh elements.
345 
346  The parameter @a e_ordering describes how the local DOFs in each element
347  should be ordered, see ElementDofOrdering.
348 
349  For discontinuous spaces, the element restriction corresponds to a
350  permutation of the degrees of freedom, implemented by the
351  L2ElementRestriction class.
352 
353  The returned Operator is owned by the FiniteElementSpace. */
354  const Operator *GetElementRestriction(ElementDofOrdering e_ordering) const;
355 
356  /// Return an Operator that converts L-vectors to E-vectors on each face.
357  virtual const Operator *GetFaceRestriction(
358  ElementDofOrdering e_ordering, FaceType,
360 
361  /** @brief Return a QuadratureInterpolator that interpolates E-vectors to
362  quadrature point values and/or derivatives (Q-vectors). */
363  /** An E-vector represents the element-wise discontinuous version of the FE
364  space and can be obtained, for example, from a GridFunction using the
365  Operator returned by GetElementRestriction().
366 
367  All elements will use the same IntegrationRule, @a ir as the target
368  quadrature points. */
370  const IntegrationRule &ir) const;
371 
372  /** @brief Return a QuadratureInterpolator that interpolates E-vectors to
373  quadrature point values and/or derivatives (Q-vectors). */
374  /** An E-vector represents the element-wise discontinuous version of the FE
375  space and can be obtained, for example, from a GridFunction using the
376  Operator returned by GetElementRestriction().
377 
378  The target quadrature points in the elements are described by the given
379  QuadratureSpace, @a qs. */
381  const QuadratureSpace &qs) const;
382 
383  /** @brief Return a FaceQuadratureInterpolator that interpolates E-vectors to
384  quadrature point values and/or derivatives (Q-vectors). */
386  const IntegrationRule &ir, FaceType type) const;
387 
388  /// Returns vector dimension.
389  inline int GetVDim() const { return vdim; }
390 
391  /// Returns the order of the i'th finite element
392  int GetOrder(int i) const;
393  /// Returns the order of the i'th face finite element
394  int GetFaceOrder(int i) const;
395 
396  /// Returns number of degrees of freedom.
397  inline int GetNDofs() const { return ndofs; }
398 
399  /// Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
400  inline int GetVSize() const { return vdim * ndofs; }
401 
402  /// Return the number of vector true (conforming) dofs.
403  virtual int GetTrueVSize() const { return GetConformingVSize(); }
404 
405  /// Returns the number of conforming ("true") degrees of freedom
406  /// (if the space is on a nonconforming mesh with hanging nodes).
407  int GetNConformingDofs() const;
408 
409  int GetConformingVSize() const { return vdim * GetNConformingDofs(); }
410 
411  /// Return the ordering method.
412  inline Ordering::Type GetOrdering() const { return ordering; }
413 
414  const FiniteElementCollection *FEColl() const { return fec; }
415 
416  /// Number of all scalar vertex dofs
417  int GetNVDofs() const { return nvdofs; }
418  /// Number of all scalar edge-interior dofs
419  int GetNEDofs() const { return nedofs; }
420  /// Number of all scalar face-interior dofs
421  int GetNFDofs() const { return nfdofs; }
422 
423  /// Returns number of vertices in the mesh.
424  inline int GetNV() const { return mesh->GetNV(); }
425 
426  /// Returns number of elements in the mesh.
427  inline int GetNE() const { return mesh->GetNE(); }
428 
429  /// Returns number of faces (i.e. co-dimension 1 entities) in the mesh.
430  /** The co-dimension 1 entities are those that have dimension 1 less than the
431  mesh dimension, e.g. for a 2D mesh, the faces are the 1D entities, i.e.
432  the edges. */
433  inline int GetNF() const { return mesh->GetNumFaces(); }
434 
435  /// Returns number of boundary elements in the mesh.
436  inline int GetNBE() const { return mesh->GetNBE(); }
437 
438  /// Returns the number of faces according to the requested type.
439  /** If type==Boundary returns only the "true" number of boundary faces
440  contrary to GetNBE() that returns "fake" boundary faces associated to
441  visualization for GLVis.
442  Similarly, if type==Interior, the "fake" boundary faces associated to
443  visualization are counted as interior faces. */
444  inline int GetNFbyType(FaceType type) const
445  { return mesh->GetNFbyType(type); }
446 
447  /// Returns the type of element i.
448  inline int GetElementType(int i) const
449  { return mesh->GetElementType(i); }
450 
451  /// Returns the vertices of element i.
452  inline void GetElementVertices(int i, Array<int> &vertices) const
453  { mesh->GetElementVertices(i, vertices); }
454 
455  /// Returns the type of boundary element i.
456  inline int GetBdrElementType(int i) const
457  { return mesh->GetBdrElementType(i); }
458 
459  /// Returns ElementTransformation for the @a i-th element.
461  { return mesh->GetElementTransformation(i); }
462 
463  /** @brief Returns the transformation defining the @a i-th element in the
464  user-defined variable @a ElTr. */
466  { mesh->GetElementTransformation(i, ElTr); }
467 
468  /// Returns ElementTransformation for the @a i-th boundary element.
470  { return mesh->GetBdrElementTransformation(i); }
471 
472  int GetAttribute(int i) const { return mesh->GetAttribute(i); }
473 
474  int GetBdrAttribute(int i) const { return mesh->GetBdrAttribute(i); }
475 
476  /// Returns indexes of degrees of freedom in array dofs for i'th element.
477  virtual void GetElementDofs(int i, Array<int> &dofs) const;
478 
479  /// Returns indexes of degrees of freedom for i'th boundary element.
480  virtual void GetBdrElementDofs(int i, Array<int> &dofs) const;
481 
482  /** @brief Returns the indexes of the degrees of freedom for i'th face
483  including the dofs for the edges and the vertices of the face. */
484  virtual void GetFaceDofs(int i, Array<int> &dofs) const;
485 
486  /** @brief Returns the indexes of the degrees of freedom for i'th edge
487  including the dofs for the vertices of the edge. */
488  void GetEdgeDofs(int i, Array<int> &dofs) const;
489 
490  void GetVertexDofs(int i, Array<int> &dofs) const;
491 
492  void GetElementInteriorDofs(int i, Array<int> &dofs) const;
493 
494  void GetFaceInteriorDofs(int i, Array<int> &dofs) const;
495 
496  int GetNumElementInteriorDofs(int i) const
498 
499  void GetEdgeInteriorDofs(int i, Array<int> &dofs) const;
500 
501  void DofsToVDofs(Array<int> &dofs, int ndofs = -1) const;
502 
503  void DofsToVDofs(int vd, Array<int> &dofs, int ndofs = -1) const;
504 
505  int DofToVDof(int dof, int vd, int ndofs = -1) const;
506 
507  int VDofToDof(int vdof) const
508  { return (ordering == Ordering::byNODES) ? (vdof%ndofs) : (vdof/vdim); }
509 
510  static void AdjustVDofs(Array<int> &vdofs);
511 
512  /// Returns indexes of degrees of freedom in array dofs for i'th element.
513  void GetElementVDofs(int i, Array<int> &vdofs) const;
514 
515  /// Returns indexes of degrees of freedom for i'th boundary element.
516  void GetBdrElementVDofs(int i, Array<int> &vdofs) const;
517 
518  /// Returns indexes of degrees of freedom for i'th face element (2D and 3D).
519  void GetFaceVDofs(int i, Array<int> &vdofs) const;
520 
521  /// Returns indexes of degrees of freedom for i'th edge.
522  void GetEdgeVDofs(int i, Array<int> &vdofs) const;
523 
524  void GetVertexVDofs(int i, Array<int> &vdofs) const;
525 
526  void GetElementInteriorVDofs(int i, Array<int> &vdofs) const;
527 
528  void GetEdgeInteriorVDofs(int i, Array<int> &vdofs) const;
529 
531 
532  /** @brief Reorder the scalar DOFs based on the element ordering.
533 
534  The new ordering is constructed as follows: 1) loop over all elements as
535  ordered in the Mesh; 2) for each element, assign new indices to all of
536  its current DOFs that are still unassigned; the new indices we assign are
537  simply the sequence `0,1,2,...`; if there are any signed DOFs their sign
538  is preserved. */
540 
541  /** @brief Return a reference to the internal Table that stores the lists of
542  scalar dofs, for each mesh element, as returned by GetElementDofs(). */
543  const Table &GetElementToDofTable() const { return *elem_dof; }
544 
545  /** @brief Return a reference to the internal Table that stores the lists of
546  scalar dofs, for each boundary mesh element, as returned by
547  GetBdrElementDofs(). */
549  { if (!bdrElem_dof) { BuildBdrElementToDofTable(); } return *bdrElem_dof; }
550 
551  /** @brief Return a reference to the internal Table that stores the lists of
552  scalar dofs, for each face in the mesh, as returned by GetFaceDofs(). In
553  this context, "face" refers to a (dim-1)-dimensional mesh entity. */
554  /** @note In the case of a NURBS space, the rows corresponding to interior
555  faces will be empty. */
556  const Table &GetFaceToDofTable() const
557  { if (!face_dof) { BuildFaceToDofTable(); } return *face_dof; }
558 
559  /** @brief Initialize internal data that enables the use of the methods
560  GetElementForDof() and GetLocalDofForDof(). */
561  void BuildDofToArrays();
562 
563  /// Return the index of the first element that contains dof @a i.
564  /** This method can be called only after setup is performed using the method
565  BuildDofToArrays(). */
566  int GetElementForDof(int i) const { return dof_elem_array[i]; }
567  /// Return the local dof index in the first element that contains dof @a i.
568  /** This method can be called only after setup is performed using the method
569  BuildDofToArrays(). */
570  int GetLocalDofForDof(int i) const { return dof_ldof_array[i]; }
571 
572  /** @brief Returns pointer to the FiniteElement in the FiniteElementCollection
573  associated with i'th element in the mesh object. */
574  const FiniteElement *GetFE(int i) const;
575 
576  /** @brief Returns pointer to the FiniteElement in the FiniteElementCollection
577  associated with i'th boundary face in the mesh object. */
578  const FiniteElement *GetBE(int i) const;
579 
580  /** @brief Returns pointer to the FiniteElement in the FiniteElementCollection
581  associated with i'th face in the mesh object. Faces in this case refer
582  to the MESHDIM-1 primitive so in 2D they are segments and in 1D they are
583  points.*/
584  const FiniteElement *GetFaceElement(int i) const;
585 
586  /** @brief Returns pointer to the FiniteElement in the FiniteElementCollection
587  associated with i'th edge in the mesh object. */
588  const FiniteElement *GetEdgeElement(int i) const;
589 
590  /// Return the trace element from element 'i' to the given 'geom_type'
591  const FiniteElement *GetTraceElement(int i, Geometry::Type geom_type) const;
592 
593  /** @brief Mark degrees of freedom associated with boundary elements with
594  the specified boundary attributes (marked in 'bdr_attr_is_ess').
595  For spaces with 'vdim' > 1, the 'component' parameter can be used
596  to restricts the marked vDOFs to the specified component. */
597  virtual void GetEssentialVDofs(const Array<int> &bdr_attr_is_ess,
598  Array<int> &ess_vdofs,
599  int component = -1) const;
600 
601  /** @brief Get a list of essential true dofs, ess_tdof_list, corresponding to the
602  boundary attributes marked in the array bdr_attr_is_ess.
603  For spaces with 'vdim' > 1, the 'component' parameter can be used
604  to restricts the marked tDOFs to the specified component. */
605  virtual void GetEssentialTrueDofs(const Array<int> &bdr_attr_is_ess,
606  Array<int> &ess_tdof_list,
607  int component = -1);
608 
609  /// Convert a Boolean marker array to a list containing all marked indices.
610  static void MarkerToList(const Array<int> &marker, Array<int> &list);
611 
612  /** @brief Convert an array of indices (list) to a Boolean marker array where all
613  indices in the list are marked with the given value and the rest are set
614  to zero. */
615  static void ListToMarker(const Array<int> &list, int marker_size,
616  Array<int> &marker, int mark_val = -1);
617 
618  /** @brief For a partially conforming FE space, convert a marker array (nonzero
619  entries are true) on the partially conforming dofs to a marker array on
620  the conforming dofs. A conforming dofs is marked iff at least one of its
621  dependent dofs is marked. */
622  void ConvertToConformingVDofs(const Array<int> &dofs, Array<int> &cdofs);
623 
624  /** @brief For a partially conforming FE space, convert a marker array (nonzero
625  entries are true) on the conforming dofs to a marker array on the
626  (partially conforming) dofs. A dof is marked iff it depends on a marked
627  conforming dofs, where dependency is defined by the ConformingRestriction
628  matrix; in other words, a dof is marked iff it corresponds to a marked
629  conforming dof. */
630  void ConvertFromConformingVDofs(const Array<int> &cdofs, Array<int> &dofs);
631 
632  /** @brief Generate the global restriction matrix from a discontinuous
633  FE space to the continuous FE space of the same polynomial degree. */
635 
636  /** @brief Generate the global restriction matrix from a discontinuous
637  FE space to the piecewise constant FE space. */
639 
640  /** @brief Construct the restriction matrix from the FE space given by
641  (*this) to the lower degree FE space given by (*lfes) which
642  is defined on the same mesh. */
644 
645  /** @brief Construct and return an Operator that can be used to transfer
646  GridFunction data from @a coarse_fes, defined on a coarse mesh, to @a
647  this FE space, defined on a refined mesh. */
648  /** It is assumed that the mesh of this FE space is a refinement of the mesh
649  of @a coarse_fes and the CoarseFineTransformations returned by the method
650  Mesh::GetRefinementTransforms() of the refined mesh are set accordingly.
651  The Operator::Type of @a T can be set to request an Operator of the set
652  type. Currently, only Operator::MFEM_SPARSEMAT and Operator::ANY_TYPE
653  (matrix-free) are supported. When Operator::ANY_TYPE is requested, the
654  choice of the particular Operator sub-class is left to the method. This
655  method also works in parallel because the transfer operator is local to
656  the MPI task when the input is a synchronized ParGridFunction. */
657  void GetTransferOperator(const FiniteElementSpace &coarse_fes,
658  OperatorHandle &T) const;
659 
660  /** @brief Construct and return an Operator that can be used to transfer
661  true-dof data from @a coarse_fes, defined on a coarse mesh, to @a this FE
662  space, defined on a refined mesh.
663 
664  This method calls GetTransferOperator() and multiplies the result by the
665  prolongation operator of @a coarse_fes on the right, and by the
666  restriction operator of this FE space on the left.
667 
668  The Operator::Type of @a T can be set to request an Operator of the set
669  type. In serial, the supported types are: Operator::MFEM_SPARSEMAT and
670  Operator::ANY_TYPE (matrix-free). In parallel, the supported types are:
671  Operator::Hypre_ParCSR and Operator::ANY_TYPE. Any other type is treated
672  as Operator::ANY_TYPE: the operator representation choice is made by this
673  method. */
674  virtual void GetTrueTransferOperator(const FiniteElementSpace &coarse_fes,
675  OperatorHandle &T) const;
676 
677  /** @brief Reflect changes in the mesh: update number of DOFs, etc. Also, calculate
678  GridFunction transformation operator (unless want_transform is false).
679  Safe to call multiple times, does nothing if space already up to date. */
680  virtual void Update(bool want_transform = true);
681 
682  /// Get the GridFunction update operator.
683  const Operator* GetUpdateOperator() { Update(); return Th.Ptr(); }
684 
685  /// Return the update operator in the given OperatorHandle, @a T.
687 
688  /** @brief Set the ownership of the update operator: if set to false, the
689  Operator returned by GetUpdateOperator() must be deleted outside the
690  FiniteElementSpace. */
691  /** The update operator ownership is automatically reset to true when a new
692  update operator is created by the Update() method. */
693  void SetUpdateOperatorOwner(bool own) { Th.SetOperatorOwner(own); }
694 
695  /// Specify the Operator::Type to be used by the update operators.
696  /** The default type is Operator::ANY_TYPE which leaves the choice to this
697  class. The other currently supported option is Operator::MFEM_SPARSEMAT
698  which is only guaranteed to be honored for a refinement update operator.
699  Any other type will be treated as Operator::ANY_TYPE.
700  @note This operation destroys the current update operator (if owned). */
702 
703  /// Free the GridFunction update operator (if any), to save memory.
704  virtual void UpdatesFinished() { Th.Clear(); }
705 
706  /// Return update counter (see Mesh::sequence)
707  long GetSequence() const { return sequence; }
708 
709  /// Return whether or not the space is discontinuous (L2)
710  bool IsDGSpace() const
711  {
712  return dynamic_cast<const L2_FECollection*>(fec) != NULL;
713  }
714 
715  /// Save finite element space to output stream @a out.
716  void Save(std::ostream &out) const;
717 
718  /** @brief Read a FiniteElementSpace from a stream. The returned
719  FiniteElementCollection is owned by the caller. */
720  FiniteElementCollection *Load(Mesh *m, std::istream &input);
721 
722  virtual ~FiniteElementSpace();
723 };
724 
725 
726 /// Class representing the storage layout of a QuadratureFunction.
727 /** Multiple QuadratureFunction%s can share the same QuadratureSpace. */
729 {
730 protected:
731  friend class QuadratureFunction; // Uses the element_offsets.
732 
734  int order;
735  int size;
736 
738  int *element_offsets; // scalar offsets; size = number of elements + 1
739 
740  // protected functions
741 
742  // Assuming mesh and order are set, construct the members: int_rule,
743  // element_offsets, and size.
744  void Construct();
745 
746 public:
747  /// Create a QuadratureSpace based on the global rules from #IntRules.
748  QuadratureSpace(Mesh *mesh_, int order_)
749  : mesh(mesh_), order(order_) { Construct(); }
750 
751  /// Read a QuadratureSpace from the stream @a in.
752  QuadratureSpace(Mesh *mesh_, std::istream &in);
753 
754  virtual ~QuadratureSpace() { delete [] element_offsets; }
755 
756  /// Return the total number of quadrature points.
757  int GetSize() const { return size; }
758 
759  /// Returns the mesh
760  inline Mesh *GetMesh() const { return mesh; }
761 
762  /// Returns number of elements in the mesh.
763  inline int GetNE() const { return mesh->GetNE(); }
764 
765  /// Get the IntegrationRule associated with mesh element @a idx.
766  const IntegrationRule &GetElementIntRule(int idx) const
767  { return *int_rule[mesh->GetElementBaseGeometry(idx)]; }
768 
769  /// Write the QuadratureSpace to the stream @a out.
770  void Save(std::ostream &out) const;
771 };
772 
773 
774 /** @brief Base class for transfer algorithms that construct transfer Operator%s
775  between two finite element (FE) spaces. */
776 /** Generally, the two FE spaces (domain and range) can be defined on different
777  meshes. */
779 {
780 protected:
781  FiniteElementSpace &dom_fes; ///< Domain FE space
782  FiniteElementSpace &ran_fes; ///< Range FE space
783 
784  /** @brief Desired Operator::Type for the construction of all operators
785  defined by the underlying transfer algorithm. It can be ignored by
786  derived classes. */
788 
789  OperatorHandle fw_t_oper; ///< Forward true-dof operator
790  OperatorHandle bw_t_oper; ///< Backward true-dof operator
791 
792 #ifdef MFEM_USE_MPI
793  bool parallel;
794 #endif
795  bool Parallel() const
796  {
797 #ifndef MFEM_USE_MPI
798  return false;
799 #else
800  return parallel;
801 #endif
802  }
803 
805  FiniteElementSpace &fes_out,
806  const Operator &oper,
807  OperatorHandle &t_oper);
808 
809 public:
810  /** Construct a transfer algorithm between the domain, @a dom_fes_, and
811  range, @a ran_fes_, FE spaces. */
812  GridTransfer(FiniteElementSpace &dom_fes_, FiniteElementSpace &ran_fes_);
813 
814  /// Virtual destructor
815  virtual ~GridTransfer() { }
816 
817  /** @brief Set the desired Operator::Type for the construction of all
818  operators defined by the underlying transfer algorithm. */
819  /** The default value is Operator::ANY_TYPE which typically corresponds to
820  a matrix-free operator representation. Note that derived classes are not
821  required to support this setting and can ignore it. */
822  void SetOperatorType(Operator::Type type) { oper_type = type; }
823 
824  /** @brief Return an Operator that transfers GridFunction%s from the domain
825  FE space to GridFunction%s in the range FE space. */
826  virtual const Operator &ForwardOperator() = 0;
827 
828  /** @brief Return an Operator that transfers GridFunction%s from the range
829  FE space back to GridFunction%s in the domain FE space. */
830  virtual const Operator &BackwardOperator() = 0;
831 
832  /** @brief Return an Operator that transfers true-dof Vector%s from the
833  domain FE space to true-dof Vector%s in the range FE space. */
834  /** This method is implemented in the base class, based on ForwardOperator(),
835  however, derived classes can overload the construction, if necessary. */
836  virtual const Operator &TrueForwardOperator()
837  {
839  }
840 
841  /** @brief Return an Operator that transfers true-dof Vector%s from the
842  range FE space back to true-dof Vector%s in the domain FE space. */
843  /** This method is implemented in the base class, based on
844  BackwardOperator(), however, derived classes can overload the
845  construction, if necessary. */
847  {
849  }
850 };
851 
852 
853 /** @brief Transfer data between a coarse mesh and an embedded refined mesh
854  using interpolation. */
855 /** The forward, coarse-to-fine, transfer uses nodal interpolation. The
856  backward, fine-to-coarse, transfer is defined locally (on a coarse element)
857  as B = (F^t M_f F)^{-1} F^t M_f, where F is the forward transfer matrix, and
858  M_f is a mass matrix on the union of all fine elements comprising the coarse
859  element. Note that the backward transfer operator, B, is a left inverse of
860  the forward transfer operator, F, i.e. B F = I. Both F and B are defined in
861  reference space and do not depend on the actual physical shape of the mesh
862  elements.
863 
864  It is assumed that both the coarse and the fine FiniteElementSpace%s use
865  compatible types of elements, e.g. finite elements with the same map-type
866  (VALUE, INTEGRAL, H_DIV, H_CURL - see class FiniteElement). Generally, the
867  FE spaces can have different orders, however, in order for the backward
868  operator to be well-defined, the (local) number of the fine dofs should not
869  be smaller than the number of coarse dofs. */
871 {
872 protected:
873  BilinearFormIntegrator *mass_integ; ///< Ownership depends on #own_mass_integ
874  bool own_mass_integ; ///< Ownership flag for #mass_integ
875 
876  OperatorHandle F; ///< Forward, coarse-to-fine, operator
877  OperatorHandle B; ///< Backward, fine-to-coarse, operator
878 
879 public:
881  FiniteElementSpace &fine_fes)
882  : GridTransfer(coarse_fes, fine_fes),
883  mass_integ(NULL), own_mass_integ(false)
884  { }
885 
886  virtual ~InterpolationGridTransfer();
887 
888  /** @brief Assign a mass integrator to be used in the construction of the
889  backward, fine-to-coarse, transfer operator. */
890  void SetMassIntegrator(BilinearFormIntegrator *mass_integ_,
891  bool own_mass_integ_ = true);
892 
893  virtual const Operator &ForwardOperator();
894 
895  virtual const Operator &BackwardOperator();
896 };
897 
898 
899 /** @brief Transfer data between a coarse mesh and an embedded refined mesh
900  using L2 projection. */
901 /** The forward, coarse-to-fine, transfer uses L2 projection. The backward,
902  fine-to-coarse, transfer is defined locally (on a coarse element) as
903  B = (F^t M_f F)^{-1} F^t M_f, where F is the forward transfer matrix, and
904  M_f is the mass matrix on the union of all fine elements comprising the
905  coarse element. Note that the backward transfer operator, B, is a left
906  inverse of the forward transfer operator, F, i.e. B F = I. Both F and B are
907  defined in physical space and, generally, vary between different mesh
908  elements.
909 
910  This class currently only fully supports L2 finite element spaces and fine
911  meshes that are a uniform refinement of the coarse mesh. Generally, the
912  coarse and fine FE spaces can have different orders, however, in order for
913  the backward operator to be well-defined, the number of the fine dofs (in a
914  coarse element) should not be smaller than the number of coarse dofs.
915 
916  If used on H1 finite element spaces, the transfer will be performed locally,
917  and the value of shared (interface) degrees of freedom will be determined by
918  the value of the last transfer to be performed (according to the element
919  numbering in the finite element space). As a consequence, the mass
920  conservation properties for this operator from the L2 case do not carry over
921  to H1 spaces. */
923 {
924 protected:
925  /** Class representing projection operator between a high-order L2 finite
926  element space on a coarse mesh, and a low-order L2 finite element space
927  on a refined mesh (LOR). We assume that the low-order space, fes_lor,
928  lives on a mesh obtained by refining the mesh of the high-order space,
929  fes_ho. */
930  class L2Projection : public Operator
931  {
932  const FiniteElementSpace &fes_ho;
933  const FiniteElementSpace &fes_lor;
934 
935  int ndof_lor, ndof_ho, nref;
936 
937  Table ho2lor;
938 
939  DenseTensor R, P;
940 
941  public:
942  L2Projection(const FiniteElementSpace &fes_ho_,
943  const FiniteElementSpace &fes_lor_);
944  /// Perform the L2 projection onto the LOR space
945  virtual void Mult(const Vector &x, Vector &y) const;
946  /// Perform the mass conservative left-inverse prolongation operation.
947  /// This functionality is also provided as an Operator by L2Prolongation.
948  void Prolongate(const Vector &x, Vector &y) const;
949  virtual ~L2Projection() { }
950  };
951 
952  /** Mass-conservative prolongation operator going in the opposite direction
953  as L2Projection. This operator is a left inverse to the L2Projection. */
954  class L2Prolongation : public Operator
955  {
956  const L2Projection &l2proj;
957 
958  public:
959  L2Prolongation(const L2Projection &l2proj_)
960  : Operator(l2proj_.Width(), l2proj_.Height()), l2proj(l2proj_) { }
961  void Mult(const Vector &x, Vector &y) const
962  {
963  l2proj.Prolongate(x, y);
964  }
965  virtual ~L2Prolongation() { }
966  };
967 
968  L2Projection *F; ///< Forward, coarse-to-fine, operator
969  L2Prolongation *B; ///< Backward, fine-to-coarse, operator
970 
971 public:
973  FiniteElementSpace &fine_fes)
974  : GridTransfer(coarse_fes, fine_fes),
975  F(NULL), B(NULL)
976  { }
977 
978  virtual const Operator &ForwardOperator();
979 
980  virtual const Operator &BackwardOperator();
981 };
982 
983 inline bool UsesTensorBasis(const FiniteElementSpace& fes)
984 {
985  return dynamic_cast<const mfem::TensorBasisElement *>(fes.GetFE(0))!=nullptr;
986 }
987 
988 }
989 
990 #endif
Abstract class for all finite elements.
Definition: fe.hpp:235
int GetNFbyType(FaceType type) const
Returns the number of faces according to the requested type.
Definition: fespace.hpp:444
const Operator & MakeTrueOperator(FiniteElementSpace &fes_in, FiniteElementSpace &fes_out, const Operator &oper, OperatorHandle &t_oper)
Definition: fespace.cpp:2550
Ordering::Type GetOrdering() const
Return the ordering method.
Definition: fespace.hpp:412
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:400
void GetVertexVDofs(int i, Array< int > &vdofs) const
Definition: fespace.cpp:197
int GetAttribute(int i) const
Definition: fespace.hpp:472
int GetNDofs() const
Returns number of degrees of freedom.
Definition: fespace.hpp:397
int ndofs
Number of degrees of freedom. Number of unknowns is ndofs * vdim.
Definition: fespace.hpp:108
L2FaceValues
Definition: restriction.hpp:26
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:1053
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:548
void GetElementInteriorVDofs(int i, Array< int > &vdofs) const
Definition: fespace.cpp:203
int DofToVDof(int dof, int vd, int ndofs=-1) const
Definition: fespace.cpp:144
const SparseMatrix * GetConformingProlongation() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:874
virtual void Update(bool want_transform=true)
Reflect changes in the mesh: update number of DOFs, etc. Also, calculate GridFunction transformation ...
Definition: fespace.cpp:2233
std::unordered_map< const key_face, Operator *, key_hash > map_L2F
Definition: fespace.hpp:148
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:452
BilinearFormIntegrator * mass_integ
Ownership depends on own_mass_integ.
Definition: fespace.hpp:873
void BuildElementToDofTable() const
Definition: fespace.cpp:215
const IntegrationRule & GetElementIntRule(int idx) const
Get the IntegrationRule associated with mesh element idx.
Definition: fespace.hpp:766
void BuildNURBSFaceToDofTable() const
Generates partial face_dof table for a NURBS space.
Definition: fespace.cpp:1581
DerefinementOperator(const FiniteElementSpace *f_fes, const FiniteElementSpace *c_fes, BilinearFormIntegrator *mass_integ)
Definition: fespace.cpp:1247
void Prolongate(const Vector &x, Vector &y) const
Definition: fespace.cpp:2888
L2Prolongation(const L2Projection &l2proj_)
Definition: fespace.hpp:959
static const int NumGeom
Definition: geom.hpp:41
FiniteElementSpace & dom_fes
Domain FE space.
Definition: fespace.hpp:781
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:2108
Matrix-free transfer operator between finite element spaces on the same mesh.
Definition: transfer.hpp:57
bool Conforming() const
Definition: mesh.hpp:1213
SparseMatrix * RefinementMatrix(int old_ndofs, const Table *old_elem_dof)
Definition: fespace.cpp:1095
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:740
SparseMatrix * DerefinementMatrix(int old_ndofs, const Table *old_elem_dof)
Calculate GridFunction restriction matrix after mesh derefinement.
Definition: fespace.cpp:1403
void GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:173
Ordering::Type ordering
Definition: fespace.hpp:105
const Geometry::Type geom
Definition: ex1.cpp:40
static int DecodeDof(int dof, double &sign)
Definition: fespace.hpp:178
static int DecodeDof(int dof)
Helpers to remove encoded sign from a DOF.
Definition: fespace.hpp:173
SparseMatrix * RefinementMatrix_main(const int coarse_ndofs, const Table &coarse_elem_dof, const DenseTensor localP[]) const
Definition: fespace.cpp:1010
int GetSize() const
Return the total number of quadrature points.
Definition: fespace.hpp:757
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:71
Element::Type GetElementType(int i) const
Returns the type of element i.
Definition: mesh.cpp:5030
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
int VDofToDof(int vdof) const
Definition: fespace.hpp:507
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: fespace.hpp:961
Element::Type GetBdrElementType(int i) const
Returns the type of boundary element i.
Definition: mesh.cpp:5035
int GetFaceOrder(int i) const
Returns the order of the i&#39;th face finite element.
Definition: fespace.cpp:102
virtual const Operator & BackwardOperator()=0
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
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:346
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:100
OperatorHandle L2E_lex
Definition: fespace.hpp:135
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:760
Operator::Type oper_type
Desired Operator::Type for the construction of all operators defined by the underlying transfer algor...
Definition: fespace.hpp:787
int GetNumElementInteriorDofs(int i) const
Definition: fespace.hpp:496
int GetBdrAttribute(int i) const
Definition: fespace.hpp:474
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:415
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:737
void GetDegenerateFaceDofs(int index, Array< int > &dofs, Geometry::Type master_geom) const
Definition: fespace.cpp:632
int GetBdrElementType(int i) const
Returns the type of boundary element i.
Definition: fespace.hpp:456
int GetNV() const
Returns number of vertices in the mesh.
Definition: fespace.hpp:424
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:1191
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:2073
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:2192
FiniteElementCollection * Load(Mesh *m, std::istream &input)
Read a FiniteElementSpace from a stream. The returned FiniteElementCollection is owned by the caller...
Definition: fespace.cpp:2379
virtual const Operator & ForwardOperator()
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
Definition: fespace.cpp:2916
bool UsesTensorBasis(const FiniteElementSpace &fes)
Definition: fespace.hpp:983
const IntegrationRule * int_rule[Geometry::NumGeom]
Definition: fespace.hpp:737
L2ProjectionGridTransfer(FiniteElementSpace &coarse_fes, FiniteElementSpace &fine_fes)
Definition: fespace.hpp:972
const FiniteElementCollection * fec
Associated FE collection (not owned).
Definition: fespace.hpp:97
OperatorHandle Th
Transformation to apply to GridFunctions after space Update().
Definition: fespace.hpp:132
NURBSExtension * GetNURBSext()
Definition: fespace.hpp:317
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:834
void GetLocalDerefinementMatrices(Geometry::Type geom, DenseTensor &localR) const
Definition: fespace.cpp:1379
int GetConformingVSize() const
Definition: fespace.hpp:409
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:556
Array< int > dof_elem_array
Definition: fespace.hpp:118
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:434
virtual const Operator & TrueBackwardOperator()
Return an Operator that transfers true-dof Vectors from the range FE space back to true-dof Vectors i...
Definition: fespace.hpp:846
void GetFaceInteriorDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:2032
Derefinement operator, used by the friend class InterpolationGridTransfer.
Definition: fespace.hpp:220
virtual int DofForGeometry(Geometry::Type GeomType) const =0
virtual void GetFaceDofs(int i, Array< int > &dofs) const
Returns the indexes of the degrees of freedom for i&#39;th face including the dofs for the edges and the ...
Definition: fespace.cpp:1898
virtual ~GridTransfer()
Virtual destructor.
Definition: fespace.hpp:815
ElementTransformation * GetBdrElementTransformation(int i)
Returns the transformation defining the i-th boundary element.
Definition: mesh.cpp:417
FiniteElementSpace & ran_fes
Range FE space.
Definition: fespace.hpp:782
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:474
const Operator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
Definition: fespace.cpp:894
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition: mesh.cpp:4184
const SparseMatrix * GetConformingRestriction() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:881
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: fespace.cpp:1348
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:427
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:763
void GetVertexDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1995
virtual const Operator & BackwardOperator()
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
Definition: fespace.cpp:2703
bool Nonconforming() const
Definition: mesh.hpp:1214
FaceType
Definition: mesh.hpp:45
void BuildBdrElementToDofTable() const
Definition: fespace.cpp:237
const FiniteElement * GetEdgeElement(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th edge in the ...
Definition: fespace.cpp:2102
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: fespace.cpp:1154
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition: fespace.hpp:436
Data type sparse matrix.
Definition: sparsemat.hpp:46
Native ordering as defined by the FiniteElement.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:65
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:314
OperatorHandle F
Forward, coarse-to-fine, operator.
Definition: fespace.hpp:876
virtual const Operator & ForwardOperator()=0
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
A class that performs interpolation from an E-vector to quadrature point values and/or derivatives (Q...
GridTransfer(FiniteElementSpace &dom_fes_, FiniteElementSpace &ran_fes_)
Definition: fespace.cpp:2535
Type
Ordering methods:
Definition: fespace.hpp:32
SparseMatrix * cR
Conforming restriction matrix such that cR.cP=I.
Definition: fespace.hpp:128
SparseMatrix * D2Const_GlobalRestrictionMatrix(FiniteElementSpace *cfes)
Generate the global restriction matrix from a discontinuous FE space to the piecewise constant FE spa...
Definition: fespace.cpp:515
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:483
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: fespace.hpp:330
virtual void GetElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:1696
void BuildConformingInterpolation() const
Calculate the cP and cR matrices for a nonconforming mesh.
Definition: fespace.cpp:683
NURBSExtension * StealNURBSext()
Definition: fespace.cpp:1548
int GetElementForDof(int i) const
Return the index of the first element that contains dof i.
Definition: fespace.hpp:566
ElementTransformation * GetBdrElementTransformation(int i) const
Returns ElementTransformation for the i-th boundary element.
Definition: fespace.hpp:469
L2Projection(const FiniteElementSpace &fes_ho_, const FiniteElementSpace &fes_lor_)
Definition: fespace.cpp:2745
static void AddDependencies(SparseMatrix &deps, Array< int > &master_dofs, Array< int > &slave_dofs, DenseMatrix &I)
Definition: fespace.cpp:594
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Returns the transformation defining the i-th element in the user-defined variable ElTr...
Definition: fespace.hpp:465
bool IsDGSpace() const
Return whether or not the space is discontinuous (L2)
Definition: fespace.hpp:710
Array< int > dof_ldof_array
Definition: fespace.hpp:118
InterpolationGridTransfer(FiniteElementSpace &coarse_fes, FiniteElementSpace &fine_fes)
Definition: fespace.hpp:880
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:403
bool Conforming() const
Definition: fespace.hpp:320
OperatorHandle bw_t_oper
Backward true-dof operator.
Definition: fespace.hpp:790
SparseMatrix * cP
Definition: fespace.hpp:126
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:389
int GetNF() const
Returns number of faces (i.e. co-dimension 1 entities) in the mesh.
Definition: fespace.hpp:433
void SetOperatorType(Operator::Type type)
Set the desired Operator::Type for the construction of all operators defined by the underlying transf...
Definition: fespace.hpp:822
Array< QuadratureInterpolator * > E2Q_array
Definition: fespace.hpp:151
A class that performs interpolation from a face E-vector to quadrature point values and/or derivative...
void GetEdgeDofs(int i, Array< int > &dofs) const
Returns the indexes of the degrees of freedom for i&#39;th edge including the dofs for the vertices of th...
Definition: fespace.cpp:1966
Operator * Ptr() const
Access the underlying Operator pointer.
Definition: handle.hpp:82
OperatorHandle L2E_nat
The element restriction operators, see GetElementRestriction().
Definition: fespace.hpp:135
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:2167
OperatorHandle B
Backward, fine-to-coarse, operator.
Definition: fespace.hpp:877
void GetElementInteriorDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:2007
void GetLocalRefinementMatrices(Geometry::Type geom, DenseTensor &localP) const
Definition: fespace.cpp:1072
static void DofsToVDofs(int ndofs, int vdim, Array< int > &dofs)
virtual const Operator & BackwardOperator()
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
Definition: fespace.cpp:2922
Base class for transfer algorithms that construct transfer Operators between two finite element (FE) ...
Definition: fespace.hpp:778
OperatorHandle fw_t_oper
Forward true-dof operator.
Definition: fespace.hpp:789
virtual const Operator * 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:922
void SetMassIntegrator(BilinearFormIntegrator *mass_integ_, bool own_mass_integ_=true)
Assign a mass integrator to be used in the construction of the backward, fine-to-coarse, transfer operator.
Definition: fespace.cpp:2661
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition: fespace.hpp:460
Type
Enumeration defining IDs for some classes derived from Operator.
Definition: operator.hpp:237
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:31
RefinementOperator(const FiniteElementSpace *fespace, Table *old_elem_dof, int old_ndofs)
Definition: fespace.cpp:1113
const NURBSExtension * GetNURBSext() const
Definition: fespace.hpp:316
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
Transfer data between a coarse mesh and an embedded refined mesh using interpolation.
Definition: fespace.hpp:870
int GetLocalDofForDof(int i) const
Return the local dof index in the first element that contains dof i.
Definition: fespace.hpp:570
virtual ~FiniteElementSpace()
Definition: fespace.cpp:2114
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:546
virtual ~QuadratureSpace()
Definition: fespace.hpp:754
void GetUpdateOperator(OperatorHandle &T)
Return the update operator in the given OperatorHandle, T.
Definition: fespace.hpp:686
L2Prolongation * B
Backward, fine-to-coarse, operator.
Definition: fespace.hpp:969
int GetOrder(int i) const
Returns the order of the i&#39;th finite element.
Definition: fespace.cpp:96
void Save(std::ostream &out) const
Write the QuadratureSpace to the stream out.
Definition: fespace.cpp:2527
void Clear()
Clear the OperatorHandle, deleting the held Operator (if owned), while leaving the type id unchanged...
Definition: handle.hpp:116
void MakeVDimMatrix(SparseMatrix &mat) const
Replicate &#39;mat&#39; in the vector dimension, according to vdim ordering mode.
Definition: fespace.cpp:846
bool own_mass_integ
Ownership flag for mass_integ.
Definition: fespace.hpp:874
The ordering method used when the number of unknowns per mesh node (vector dimension) is bigger than ...
Definition: fespace.hpp:28
void GetEdgeInteriorVDofs(int i, Array< int > &vdofs) const
Definition: fespace.cpp:209
void GetEntityDofs(int entity, int index, Array< int > &dofs, Geometry::Type master_geom=Geometry::INVALID) const
Helper to get vertex, edge or face DOFs (entity=0,1,2 resp.).
Definition: fespace.cpp:670
L2Projection * F
Forward, coarse-to-fine, operator.
Definition: fespace.hpp:968
int GetElementType(int i) const
Returns the type of element i.
Definition: fespace.hpp:448
FiniteElementSpace(Mesh *mesh, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES)
Definition: fespace.hpp:299
bool Nonconforming() const
Definition: fespace.hpp:321
virtual void GetBdrElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:1815
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:734
int GetNVDofs() const
Number of all scalar vertex dofs.
Definition: fespace.hpp:417
virtual const Operator & ForwardOperator()
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
Definition: fespace.cpp:2670
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:466
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
Definition: handle.hpp:112
bool Parallel() const
Definition: fespace.hpp:795
virtual const Operator & TrueForwardOperator()
Return an Operator that transfers true-dof Vectors from the domain FE space to true-dof Vectors in th...
Definition: fespace.hpp:836
A standard isoparametric element transformation.
Definition: eltrans.hpp:348
Mesh * mesh
The mesh that FE space lives on (not owned).
Definition: fespace.hpp:94
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:336
void ReorderElementToDofTable()
Reorder the scalar DOFs based on the element ordering.
Definition: fespace.cpp:292
GridFunction interpolation operator applicable after mesh refinement.
Definition: fespace.hpp:201
Array< FaceQuadratureInterpolator * > E2BFQ_array
Definition: fespace.hpp:153
int GetNConformingDofs() const
Definition: fespace.cpp:888
void GetEdgeInteriorDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:2020
const QuadratureInterpolator * GetQuadratureInterpolator(const IntegrationRule &ir) const
Return a QuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivat...
Definition: fespace.cpp:950
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition: fespace.hpp:65
Array< FaceQuadratureInterpolator * > E2IFQ_array
Definition: fespace.hpp:152
void BuildFaceToDofTable() const
Definition: fespace.cpp:259
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:1798
int index(int i, int j, int nx, int ny)
Definition: life.cpp:241
const Operator * GetUpdateOperator()
Get the GridFunction update operator.
Definition: fespace.hpp:683
Transfer data between a coarse mesh and an embedded refined mesh using L2 projection.
Definition: fespace.hpp:922
Lexicographic ordering for tensor-product FiniteElements.
FiniteElementSpace()
Default constructor: the object is invalid until initialized using the method Load().
Definition: fespace.cpp:60
QuadratureSpace(Mesh *mesh_, int order_)
Create a QuadratureSpace based on the global rules from IntRules.
Definition: fespace.hpp:748
void Constructor(Mesh *mesh, NURBSExtension *ext, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES)
Help function for constructors + Load().
Definition: fespace.cpp:1501
void SetUpdateOperatorOwner(bool own)
Set the ownership of the update operator: if set to false, the Operator returned by GetUpdateOperator...
Definition: fespace.hpp:693
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:185
void GetEdgeVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th edge.
Definition: fespace.cpp:191
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:453
virtual void Mult(const Vector &x, Vector &y) const
Perform the L2 projection onto the LOR space.
Definition: fespace.cpp:2862
std::size_t operator()(const key_face &k) const
Definition: fespace.hpp:140
Vector data type.
Definition: vector.hpp:51
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:414
std::tuple< bool, ElementDofOrdering, FaceType, L2FaceValues > key_face
The face restriction operators, see GetFaceRestriction().
Definition: fespace.hpp:137
int GetNFbyType(FaceType type) const
Returns the number of faces according to the requested type.
Definition: mesh.cpp:4210
Array< int > face_to_be
Definition: fespace.hpp:116
virtual void UpdatesFinished()
Free the GridFunction update operator (if any), to save memory.
Definition: fespace.hpp:704
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:543
void SetUpdateOperatorType(Operator::Type tid)
Specify the Operator::Type to be used by the update operators.
Definition: fespace.hpp:701
static bool DofFinalizable(int dof, const Array< bool > &finalized, const SparseMatrix &deps)
Definition: fespace.cpp:618
Class representing the storage layout of a QuadratureFunction.
Definition: fespace.hpp:728
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th boundary fac...
Definition: fespace.cpp:2047
NURBSExtension * NURBSext
Definition: fespace.hpp:120
void GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:179
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
Abstract operator.
Definition: operator.hpp:24
virtual const SparseMatrix * GetRestrictionMatrix() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.hpp:334
long GetSequence() const
Return update counter (see Mesh::sequence)
Definition: fespace.hpp:707
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:728
void Save(std::ostream &out) const
Save finite element space to output stream out.
Definition: fespace.cpp:2310
int GetAttribute(int i) const
Return the attribute of element i.
Definition: mesh.hpp:1047
void BuildDofToArrays()
Initialize internal data that enables the use of the methods GetElementForDof() and GetLocalDofForDof...
Definition: fespace.cpp:312
Class representing a function through its values (scalar or vector) at quadrature points...
Definition: gridfunc.hpp:670
int GetNFDofs() const
Number of all scalar face-interior dofs.
Definition: fespace.hpp:421
int GetNEDofs() const
Number of all scalar edge-interior dofs.
Definition: fespace.hpp:419
void GetElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of element i.
Definition: mesh.hpp:878
static void AdjustVDofs(Array< int > &vdofs)
Definition: fespace.cpp:160
void SetType(Operator::Type tid)
Invoke Clear() and set a new type id.
Definition: handle.hpp:124
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:221
void DofsToVDofs(Array< int > &dofs, int ndofs=-1) const
Definition: fespace.cpp:108
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:979
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:308