MFEM  v4.1.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 {
90 
91 protected:
92  /// The mesh that FE space lives on (not owned).
94 
95  /// Associated FE collection (not owned).
97 
98  /// %Vector dimension (number of unknowns per degree of freedom).
99  int vdim;
100 
101  /** Type of ordering of the vector dofs when #vdim > 1.
102  - Ordering::byNODES - first nodes, then vector dimension,
103  - Ordering::byVDIM - first vector dimension, then nodes */
105 
106  /// Number of degrees of freedom. Number of unknowns is #ndofs * #vdim.
107  int ndofs;
108 
110  int *fdofs, *bdofs;
111 
112  mutable Table *elem_dof; // if NURBS FE space, not owned; otherwise, owned.
113  Table *bdrElem_dof; // used only with NURBS FE spaces; not owned.
114 
116 
118  int own_ext;
119 
120  /** Matrix representing the prolongation from the global conforming dofs to
121  a set of intermediate partially conforming dofs, e.g. the dofs associated
122  with a "cut" space on a non-conforming mesh. */
123  mutable SparseMatrix *cP; // owned
124  /// Conforming restriction matrix such that cR.cP=I.
125  mutable SparseMatrix *cR; // owned
126  mutable bool cP_is_set;
127 
128  /// Transformation to apply to GridFunctions after space Update().
130 
131  /// The element restriction operators, see GetElementRestriction().
133  /// The face restriction operators, see GetFaceRestriction().
134  using key_face = std::tuple<bool, ElementDofOrdering, FaceType, L2FaceValues>;
135  struct key_hash
136  {
137  std::size_t operator()(const key_face& k) const
138  {
139  return std::get<0>(k)
140  + 2 * (int)std::get<1>(k)
141  + 4 * (int)std::get<2>(k)
142  + 8 * (int)std::get<3>(k);
143  }
144  };
145  using map_L2F = std::unordered_map<const key_face,Operator*,key_hash>;
146  mutable map_L2F L2F;
147 
151 
152  long sequence; // should match Mesh::GetSequence
153 
154  void UpdateNURBS();
155 
156  void Construct();
157  void Destroy();
158 
159  void BuildElementToDofTable() const;
160 
161  /// Helper to remove encoded sign from a DOF
162  static inline int DecodeDof(int dof, double& sign)
163  { return (dof >= 0) ? (sign = 1, dof) : (sign = -1, (-1 - dof)); }
164 
165  /// Helper to get vertex, edge or face DOFs (entity=0,1,2 resp.).
166  void GetEntityDofs(int entity, int index, Array<int> &dofs,
167  Geometry::Type master_geom = Geometry::INVALID) const;
168  // Get degenerate face DOFs: see explanation in method implementation.
169  void GetDegenerateFaceDofs(int index, Array<int> &dofs,
170  Geometry::Type master_geom) const;
171 
172  /// Calculate the cP and cR matrices for a nonconforming mesh.
173  void BuildConformingInterpolation() const;
174 
175  static void AddDependencies(SparseMatrix& deps, Array<int>& master_dofs,
176  Array<int>& slave_dofs, DenseMatrix& I);
177 
178  static bool DofFinalizable(int dof, const Array<bool>& finalized,
179  const SparseMatrix& deps);
180 
181  /// Replicate 'mat' in the vector dimension, according to vdim ordering mode.
182  void MakeVDimMatrix(SparseMatrix &mat) const;
183 
184  /// GridFunction interpolation operator applicable after mesh refinement.
186  {
187  const FiniteElementSpace* fespace;
189  Table* old_elem_dof; // Owned.
190 
191  public:
192  /** Construct the operator based on the elem_dof table of the original
193  (coarse) space. The class takes ownership of the table. */
194  RefinementOperator(const FiniteElementSpace* fespace,
195  Table *old_elem_dof/*takes ownership*/, int old_ndofs);
196  RefinementOperator(const FiniteElementSpace *fespace,
197  const FiniteElementSpace *coarse_fes);
198  virtual void Mult(const Vector &x, Vector &y) const;
199  virtual ~RefinementOperator();
200  };
201 
202  // Derefinement operator, used by the friend class InterpolationGridTransfer.
204  {
205  const FiniteElementSpace *fine_fes; // Not owned.
207  Table *coarse_elem_dof; // Owned.
208  Table coarse_to_fine;
209  Array<int> coarse_to_ref_type;
210  Array<Geometry::Type> ref_type_to_geom;
211  Array<int> ref_type_to_fine_elem_offset;
212 
213  public:
215  const FiniteElementSpace *c_fes,
216  BilinearFormIntegrator *mass_integ);
217  virtual void Mult(const Vector &x, Vector &y) const;
218  virtual ~DerefinementOperator();
219  };
220 
221  // This method makes the same assumptions as the method:
222  // void GetLocalRefinementMatrices(
223  // const FiniteElementSpace &coarse_fes, Geometry::Type geom,
224  // DenseTensor &localP) const
225  // which is defined below. It also assumes that the coarse fes and this have
226  // the same vector dimension, vdim.
227  SparseMatrix *RefinementMatrix_main(const int coarse_ndofs,
228  const Table &coarse_elem_dof,
229  const DenseTensor localP[]) const;
230 
232  DenseTensor &localP) const;
234  DenseTensor &localR) const;
235 
236  /** Calculate explicit GridFunction interpolation matrix (after mesh
237  refinement). NOTE: consider using the RefinementOperator class instead
238  of the fully assembled matrix, which can take a lot of memory. */
239  SparseMatrix* RefinementMatrix(int old_ndofs, const Table* old_elem_dof);
240 
241  /// Calculate GridFunction restriction matrix after mesh derefinement.
242  SparseMatrix* DerefinementMatrix(int old_ndofs, const Table* old_elem_dof);
243 
244  // This method assumes that this->mesh is a refinement of coarse_fes->mesh
245  // and that the CoarseFineTransformations of this->mesh are set accordingly.
246  // Another assumption is that the FEs of this use the same MapType as the FEs
247  // of coarse_fes. Finally, it assumes that the spaces this and coarse_fes are
248  // NOT variable-order spaces.
249  void GetLocalRefinementMatrices(const FiniteElementSpace &coarse_fes,
251  DenseTensor &localP) const;
252 
253  /// Help function for constructors + Load().
254  void Constructor(Mesh *mesh, NURBSExtension *ext,
256  int vdim = 1, int ordering = Ordering::byNODES);
257 
258 public:
259  /** @brief Default constructor: the object is invalid until initialized using
260  the method Load(). */
262 
263  /** @brief Copy constructor: deep copy all data from @a orig except the Mesh,
264  the FiniteElementCollection, ans some derived data. */
265  /** If the @a mesh or @a fec pointers are NULL (default), then the new
266  FiniteElementSpace will reuse the respective pointers from @a orig. If
267  any of these pointers is not NULL, the given pointer will be used instead
268  of the one used by @a orig.
269 
270  @note The objects pointed to by the @a mesh and @a fec parameters must be
271  either the same objects as the ones used by @a orig, or copies of them.
272  Otherwise, the behavior is undefined.
273 
274  @note Derived data objects, such as the conforming prolongation and
275  restriction matrices, and the update operator, will not be copied, even
276  if they are created in the @a orig object. */
277  FiniteElementSpace(const FiniteElementSpace &orig, Mesh *mesh = NULL,
278  const FiniteElementCollection *fec = NULL);
279 
282  int vdim = 1, int ordering = Ordering::byNODES)
283  { Constructor(mesh, NULL, fec, vdim, ordering); }
284 
285  /// Construct a NURBS FE space based on the given NURBSExtension, @a ext.
286  /** @note If the pointer @a ext is NULL, this constructor is equivalent to
287  the standard constructor with the same arguments minus the
288  NURBSExtension, @a ext. */
291  int vdim = 1, int ordering = Ordering::byNODES)
292  { Constructor(mesh, ext, fec, vdim, ordering); }
293 
294  /// Returns the mesh
295  inline Mesh *GetMesh() const { return mesh; }
296 
297  const NURBSExtension *GetNURBSext() const { return NURBSext; }
300 
301  bool Conforming() const { return mesh->Conforming(); }
302  bool Nonconforming() const { return mesh->Nonconforming(); }
303 
304  /// The returned SparseMatrix is owned by the FiniteElementSpace.
306 
307  /// The returned SparseMatrix is owned by the FiniteElementSpace.
308  const SparseMatrix *GetConformingRestriction() const;
309 
310  /// The returned Operator is owned by the FiniteElementSpace.
311  virtual const Operator *GetProlongationMatrix() const
312  { return GetConformingProlongation(); }
313 
314  /// The returned SparseMatrix is owned by the FiniteElementSpace.
315  virtual const SparseMatrix *GetRestrictionMatrix() const
316  { return GetConformingRestriction(); }
317 
318  /// Return an Operator that converts L-vectors to E-vectors.
319  /** An L-vector is a vector of size GetVSize() which is the same size as a
320  GridFunction. An E-vector represents the element-wise discontinuous
321  version of the FE space.
322 
323  The layout of the E-vector is: ND x VDIM x NE, where ND is the number of
324  degrees of freedom, VDIM is the vector dimension of the FE space, and NE
325  is the number of the mesh elements.
326 
327  The parameter @a e_ordering describes how the local DOFs in each element
328  should be ordered, see ElementDofOrdering.
329 
330  For discontinuous spaces, the element restriction corresponds to a
331  permutation of the degrees of freedom, implemented by the
332  L2ElementRestriction class.
333 
334  The returned Operator is owned by the FiniteElementSpace. */
335  const Operator *GetElementRestriction(ElementDofOrdering e_ordering) const;
336 
337  /// Return an Operator that converts L-vectors to E-vectors on each face.
338  virtual const Operator *GetFaceRestriction(
339  ElementDofOrdering e_ordering, FaceType,
341 
342  /** @brief Return a QuadratureInterpolator that interpolates E-vectors to
343  quadrature point values and/or derivatives (Q-vectors). */
344  /** An E-vector represents the element-wise discontinuous version of the FE
345  space and can be obtained, for example, from a GridFunction using the
346  Operator returned by GetElementRestriction().
347 
348  All elements will use the same IntegrationRule, @a ir as the target
349  quadrature points. */
351  const IntegrationRule &ir) const;
352 
353  /** @brief Return a QuadratureInterpolator that interpolates E-vectors to
354  quadrature point values and/or derivatives (Q-vectors). */
355  /** An E-vector represents the element-wise discontinuous version of the FE
356  space and can be obtained, for example, from a GridFunction using the
357  Operator returned by GetElementRestriction().
358 
359  The target quadrature points in the elements are described by the given
360  QuadratureSpace, @a qs. */
362  const QuadratureSpace &qs) const;
363 
364  /** @brief Return a FaceQuadratureInterpolator that interpolates E-vectors to
365  quadrature point values and/or derivatives (Q-vectors). */
367  const IntegrationRule &ir, FaceType type) const;
368 
369  /// Returns vector dimension.
370  inline int GetVDim() const { return vdim; }
371 
372  /// Returns the order of the i'th finite element
373  int GetOrder(int i) const;
374  /// Returns the order of the i'th face finite element
375  int GetFaceOrder(int i) const;
376 
377  /// Returns number of degrees of freedom.
378  inline int GetNDofs() const { return ndofs; }
379 
380  /// Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
381  inline int GetVSize() const { return vdim * ndofs; }
382 
383  /// Return the number of vector true (conforming) dofs.
384  virtual int GetTrueVSize() const { return GetConformingVSize(); }
385 
386  /// Returns the number of conforming ("true") degrees of freedom
387  /// (if the space is on a nonconforming mesh with hanging nodes).
388  int GetNConformingDofs() const;
389 
390  int GetConformingVSize() const { return vdim * GetNConformingDofs(); }
391 
392  /// Return the ordering method.
393  inline Ordering::Type GetOrdering() const { return ordering; }
394 
395  const FiniteElementCollection *FEColl() const { return fec; }
396 
397  /// Number of all scalar vertex dofs
398  int GetNVDofs() const { return nvdofs; }
399  /// Number of all scalar edge-interior dofs
400  int GetNEDofs() const { return nedofs; }
401  /// Number of all scalar face-interior dofs
402  int GetNFDofs() const { return nfdofs; }
403 
404  /// Returns number of vertices in the mesh.
405  inline int GetNV() const { return mesh->GetNV(); }
406 
407  /// Returns number of elements in the mesh.
408  inline int GetNE() const { return mesh->GetNE(); }
409 
410  /// Returns number of faces (i.e. co-dimension 1 entities) in the mesh.
411  /** The co-dimension 1 entities are those that have dimension 1 less than the
412  mesh dimension, e.g. for a 2D mesh, the faces are the 1D entities, i.e.
413  the edges. */
414  inline int GetNF() const { return mesh->GetNumFaces(); }
415 
416  /// Returns number of boundary elements in the mesh.
417  inline int GetNBE() const { return mesh->GetNBE(); }
418 
419  /// Returns the number of faces according to the requested type.
420  /** If type==Boundary returns only the "true" number of boundary faces
421  contrary to GetNBE() that returns "fake" boundary faces associated to
422  visualization for GLVis.
423  Similarly, if type==Interior, the "fake" boundary faces associated to
424  visualization are counted as interior faces. */
425  inline int GetNFbyType(FaceType type) const
426  { return mesh->GetNFbyType(type); }
427 
428  /// Returns the type of element i.
429  inline int GetElementType(int i) const
430  { return mesh->GetElementType(i); }
431 
432  /// Returns the vertices of element i.
433  inline void GetElementVertices(int i, Array<int> &vertices) const
434  { mesh->GetElementVertices(i, vertices); }
435 
436  /// Returns the type of boundary element i.
437  inline int GetBdrElementType(int i) const
438  { return mesh->GetBdrElementType(i); }
439 
440  /// Returns ElementTransformation for the @a i-th element.
442  { return mesh->GetElementTransformation(i); }
443 
444  /** @brief Returns the transformation defining the @a i-th element in the
445  user-defined variable @a ElTr. */
447  { mesh->GetElementTransformation(i, ElTr); }
448 
449  /// Returns ElementTransformation for the @a i-th boundary element.
451  { return mesh->GetBdrElementTransformation(i); }
452 
453  int GetAttribute(int i) const { return mesh->GetAttribute(i); }
454 
455  int GetBdrAttribute(int i) const { return mesh->GetBdrAttribute(i); }
456 
457  /// Returns indexes of degrees of freedom in array dofs for i'th element.
458  virtual void GetElementDofs(int i, Array<int> &dofs) const;
459 
460  /// Returns indexes of degrees of freedom for i'th boundary element.
461  virtual void GetBdrElementDofs(int i, Array<int> &dofs) const;
462 
463  /** Returns the indexes of the degrees of freedom for i'th face
464  including the dofs for the edges and the vertices of the face. */
465  virtual void GetFaceDofs(int i, Array<int> &dofs) const;
466 
467  /** Returns the indexes of the degrees of freedom for i'th edge
468  including the dofs for the vertices of the edge. */
469  void GetEdgeDofs(int i, Array<int> &dofs) const;
470 
471  void GetVertexDofs(int i, Array<int> &dofs) const;
472 
473  void GetElementInteriorDofs(int i, Array<int> &dofs) const;
474 
475  void GetFaceInteriorDofs(int i, Array<int> &dofs) const;
476 
477  int GetNumElementInteriorDofs(int i) const
479 
480  void GetEdgeInteriorDofs(int i, Array<int> &dofs) const;
481 
482  void DofsToVDofs(Array<int> &dofs, int ndofs = -1) const;
483 
484  void DofsToVDofs(int vd, Array<int> &dofs, int ndofs = -1) const;
485 
486  int DofToVDof(int dof, int vd, int ndofs = -1) const;
487 
488  int VDofToDof(int vdof) const
489  { return (ordering == Ordering::byNODES) ? (vdof%ndofs) : (vdof/vdim); }
490 
491  static void AdjustVDofs(Array<int> &vdofs);
492 
493  /// Returns indexes of degrees of freedom in array dofs for i'th element.
494  void GetElementVDofs(int i, Array<int> &vdofs) const;
495 
496  /// Returns indexes of degrees of freedom for i'th boundary element.
497  void GetBdrElementVDofs(int i, Array<int> &vdofs) const;
498 
499  /// Returns indexes of degrees of freedom for i'th face element (2D and 3D).
500  void GetFaceVDofs(int i, Array<int> &vdofs) const;
501 
502  /// Returns indexes of degrees of freedom for i'th edge.
503  void GetEdgeVDofs(int i, Array<int> &vdofs) const;
504 
505  void GetVertexVDofs(int i, Array<int> &vdofs) const;
506 
507  void GetElementInteriorVDofs(int i, Array<int> &vdofs) const;
508 
509  void GetEdgeInteriorVDofs(int i, Array<int> &vdofs) const;
510 
512 
513  /** @brief Reorder the scalar DOFs based on the element ordering.
514 
515  The new ordering is constructed as follows: 1) loop over all elements as
516  ordered in the Mesh; 2) for each element, assign new indices to all of
517  its current DOFs that are still unassigned; the new indices we assign are
518  simply the sequence `0,1,2,...`; if there are any signed DOFs their sign
519  is preserved. */
521 
522  void BuildDofToArrays();
523 
524  const Table &GetElementToDofTable() const { return *elem_dof; }
525  const Table &GetBdrElementToDofTable() const { return *bdrElem_dof; }
526 
527  int GetElementForDof(int i) const { return dof_elem_array[i]; }
528  int GetLocalDofForDof(int i) const { return dof_ldof_array[i]; }
529 
530  /// Returns pointer to the FiniteElement associated with i'th element.
531  const FiniteElement *GetFE(int i) const;
532 
533  /// Returns pointer to the FiniteElement for the i'th boundary element.
534  const FiniteElement *GetBE(int i) const;
535 
536  const FiniteElement *GetFaceElement(int i) const;
537 
538  const FiniteElement *GetEdgeElement(int i) const;
539 
540  /// Return the trace element from element 'i' to the given 'geom_type'
541  const FiniteElement *GetTraceElement(int i, Geometry::Type geom_type) const;
542 
543  /** Mark degrees of freedom associated with boundary elements with
544  the specified boundary attributes (marked in 'bdr_attr_is_ess').
545  For spaces with 'vdim' > 1, the 'component' parameter can be used
546  to restricts the marked vDOFs to the specified component. */
547  virtual void GetEssentialVDofs(const Array<int> &bdr_attr_is_ess,
548  Array<int> &ess_vdofs,
549  int component = -1) const;
550 
551  /** Get a list of essential true dofs, ess_tdof_list, corresponding to the
552  boundary attributes marked in the array bdr_attr_is_ess.
553  For spaces with 'vdim' > 1, the 'component' parameter can be used
554  to restricts the marked tDOFs to the specified component. */
555  virtual void GetEssentialTrueDofs(const Array<int> &bdr_attr_is_ess,
556  Array<int> &ess_tdof_list,
557  int component = -1);
558 
559  /// Convert a Boolean marker array to a list containing all marked indices.
560  static void MarkerToList(const Array<int> &marker, Array<int> &list);
561 
562  /** Convert an array of indices (list) to a Boolean marker array where all
563  indices in the list are marked with the given value and the rest are set
564  to zero. */
565  static void ListToMarker(const Array<int> &list, int marker_size,
566  Array<int> &marker, int mark_val = -1);
567 
568  /** For a partially conforming FE space, convert a marker array (nonzero
569  entries are true) on the partially conforming dofs to a marker array on
570  the conforming dofs. A conforming dofs is marked iff at least one of its
571  dependent dofs is marked. */
572  void ConvertToConformingVDofs(const Array<int> &dofs, Array<int> &cdofs);
573 
574  /** For a partially conforming FE space, convert a marker array (nonzero
575  entries are true) on the conforming dofs to a marker array on the
576  (partially conforming) dofs. A dof is marked iff it depends on a marked
577  conforming dofs, where dependency is defined by the ConformingRestriction
578  matrix; in other words, a dof is marked iff it corresponds to a marked
579  conforming dof. */
580  void ConvertFromConformingVDofs(const Array<int> &cdofs, Array<int> &dofs);
581 
582  /** Generate the global restriction matrix from a discontinuous
583  FE space to the continuous FE space of the same polynomial degree. */
585 
586  /** Generate the global restriction matrix from a discontinuous
587  FE space to the piecewise constant FE space. */
589 
590  /** Construct the restriction matrix from the FE space given by
591  (*this) to the lower degree FE space given by (*lfes) which
592  is defined on the same mesh. */
594 
595  /** @brief Construct and return an Operator that can be used to transfer
596  GridFunction data from @a coarse_fes, defined on a coarse mesh, to @a
597  this FE space, defined on a refined mesh. */
598  /** It is assumed that the mesh of this FE space is a refinement of the mesh
599  of @a coarse_fes and the CoarseFineTransformations returned by the method
600  Mesh::GetRefinementTransforms() of the refined mesh are set accordingly.
601  The Operator::Type of @a T can be set to request an Operator of the set
602  type. Currently, only Operator::MFEM_SPARSEMAT and Operator::ANY_TYPE
603  (matrix-free) are supported. When Operator::ANY_TYPE is requested, the
604  choice of the particular Operator sub-class is left to the method. This
605  method also works in parallel because the transfer operator is local to
606  the MPI task when the input is a synchronized ParGridFunction. */
607  void GetTransferOperator(const FiniteElementSpace &coarse_fes,
608  OperatorHandle &T) const;
609 
610  /** @brief Construct and return an Operator that can be used to transfer
611  true-dof data from @a coarse_fes, defined on a coarse mesh, to @a this FE
612  space, defined on a refined mesh.
613 
614  This method calls GetTransferOperator() and multiplies the result by the
615  prolongation operator of @a coarse_fes on the right, and by the
616  restriction operator of this FE space on the left.
617 
618  The Operator::Type of @a T can be set to request an Operator of the set
619  type. In serial, the supported types are: Operator::MFEM_SPARSEMAT and
620  Operator::ANY_TYPE (matrix-free). In parallel, the supported types are:
621  Operator::Hypre_ParCSR and Operator::ANY_TYPE. Any other type is treated
622  as Operator::ANY_TYPE: the operator representation choice is made by this
623  method. */
624  virtual void GetTrueTransferOperator(const FiniteElementSpace &coarse_fes,
625  OperatorHandle &T) const;
626 
627  /** Reflect changes in the mesh: update number of DOFs, etc. Also, calculate
628  GridFunction transformation operator (unless want_transform is false).
629  Safe to call multiple times, does nothing if space already up to date. */
630  virtual void Update(bool want_transform = true);
631 
632  /// Get the GridFunction update operator.
633  const Operator* GetUpdateOperator() { Update(); return Th.Ptr(); }
634 
635  /// Return the update operator in the given OperatorHandle, @a T.
637 
638  /** @brief Set the ownership of the update operator: if set to false, the
639  Operator returned by GetUpdateOperator() must be deleted outside the
640  FiniteElementSpace. */
641  /** The update operator ownership is automatically reset to true when a new
642  update operator is created by the Update() method. */
643  void SetUpdateOperatorOwner(bool own) { Th.SetOperatorOwner(own); }
644 
645  /// Specify the Operator::Type to be used by the update operators.
646  /** The default type is Operator::ANY_TYPE which leaves the choice to this
647  class. The other currently supported option is Operator::MFEM_SPARSEMAT
648  which is only guaranteed to be honored for a refinement update operator.
649  Any other type will be treated as Operator::ANY_TYPE.
650  @note This operation destroys the current update operator (if owned). */
652 
653  /// Free the GridFunction update operator (if any), to save memory.
654  virtual void UpdatesFinished() { Th.Clear(); }
655 
656  /// Return update counter (see Mesh::sequence)
657  long GetSequence() const { return sequence; }
658 
659  void Save(std::ostream &out) const;
660 
661  /** @brief Read a FiniteElementSpace from a stream. The returned
662  FiniteElementCollection is owned by the caller. */
663  FiniteElementCollection *Load(Mesh *m, std::istream &input);
664 
665  virtual ~FiniteElementSpace();
666 };
667 
668 
669 /// Class representing the storage layout of a QuadratureFunction.
670 /** Multiple QuadratureFunction%s can share the same QuadratureSpace. */
672 {
673 protected:
674  friend class QuadratureFunction; // Uses the element_offsets.
675 
677  int order;
678  int size;
679 
681  int *element_offsets; // scalar offsets; size = number of elements + 1
682 
683  // protected functions
684 
685  // Assuming mesh and order are set, construct the members: int_rule,
686  // element_offsets, and size.
687  void Construct();
688 
689 public:
690  /// Create a QuadratureSpace based on the global rules from #IntRules.
691  QuadratureSpace(Mesh *mesh_, int order_)
692  : mesh(mesh_), order(order_) { Construct(); }
693 
694  /// Read a QuadratureSpace from the stream @a in.
695  QuadratureSpace(Mesh *mesh_, std::istream &in);
696 
697  virtual ~QuadratureSpace() { delete [] element_offsets; }
698 
699  /// Return the total number of quadrature points.
700  int GetSize() const { return size; }
701 
702  /// Get the IntegrationRule associated with mesh element @a idx.
703  const IntegrationRule &GetElementIntRule(int idx) const
704  { return *int_rule[mesh->GetElementBaseGeometry(idx)]; }
705 
706  /// Write the QuadratureSpace to the stream @a out.
707  void Save(std::ostream &out) const;
708 };
709 
710 
711 /** @brief Base class for transfer algorithms that construct transfer Operator%s
712  between two finite element (FE) spaces. */
713 /** Generally, the two FE spaces (domain and range) can be defined on different
714  meshes. */
716 {
717 protected:
718  FiniteElementSpace &dom_fes; ///< Domain FE space
719  FiniteElementSpace &ran_fes; ///< Range FE space
720 
721  /** @brief Desired Operator::Type for the construction of all operators
722  defined by the underlying transfer algorithm. It can be ignored by
723  derived classes. */
725 
726  OperatorHandle fw_t_oper; ///< Forward true-dof operator
727  OperatorHandle bw_t_oper; ///< Backward true-dof operator
728 
729 #ifdef MFEM_USE_MPI
730  bool parallel;
731 #endif
732  bool Parallel() const
733  {
734 #ifndef MFEM_USE_MPI
735  return false;
736 #else
737  return parallel;
738 #endif
739  }
740 
742  FiniteElementSpace &fes_out,
743  const Operator &oper,
744  OperatorHandle &t_oper);
745 
746 public:
747  /** Construct a transfer algorithm between the domain, @a dom_fes_, and
748  range, @a ran_fes_, FE spaces. */
749  GridTransfer(FiniteElementSpace &dom_fes_, FiniteElementSpace &ran_fes_);
750 
751  /// Virtual destructor
752  virtual ~GridTransfer() { }
753 
754  /** @brief Set the desired Operator::Type for the construction of all
755  operators defined by the underlying transfer algorithm. */
756  /** The default value is Operator::ANY_TYPE which typically corresponds to
757  a matrix-free operator representation. Note that derived classes are not
758  required to support this setting and can ignore it. */
759  void SetOperatorType(Operator::Type type) { oper_type = type; }
760 
761  /** @brief Return an Operator that transfers GridFunction%s from the domain
762  FE space to GridFunction%s in the range FE space. */
763  virtual const Operator &ForwardOperator() = 0;
764 
765  /** @brief Return an Operator that transfers GridFunction%s from the range
766  FE space back to GridFunction%s in the domain FE space. */
767  virtual const Operator &BackwardOperator() = 0;
768 
769  /** @brief Return an Operator that transfers true-dof Vector%s from the
770  domain FE space to true-dof Vector%s in the range FE space. */
771  /** This method is implemented in the base class, based on ForwardOperator(),
772  however, derived classes can overload the construction, if necessary. */
773  virtual const Operator &TrueForwardOperator()
774  {
776  }
777 
778  /** @brief Return an Operator that transfers true-dof Vector%s from the
779  range FE space back to true-dof Vector%s in the domain FE space. */
780  /** This method is implemented in the base class, based on
781  BackwardOperator(), however, derived classes can overload the
782  construction, if necessary. */
784  {
786  }
787 };
788 
789 
790 /** @brief Transfer data between a coarse mesh and an embedded refined mesh
791  using interpolation. */
792 /** The forward, coarse-to-fine, transfer uses nodal interpolation. The
793  backward, fine-to-coarse, transfer is defined locally (on a coarse element)
794  as B = (F^t M_f F)^{-1} F^t M_f, where F is the forward transfer matrix, and
795  M_f is a mass matrix on the union of all fine elements comprising the coarse
796  element. Note that the backward transfer operator, B, is a left inverse of
797  the forward transfer operator, F, i.e. B F = I. Both F and B are defined in
798  reference space and do not depend on the actual physical shape of the mesh
799  elements.
800 
801  It is assumed that both the coarse and the fine FiniteElementSpace%s use
802  compatible types of elements, e.g. finite elements with the same map-type
803  (VALUE, INTEGRAL, H_DIV, H_CURL - see class FiniteElement). Generally, the
804  FE spaces can have different orders, however, in order for the backward
805  operator to be well-defined, the (local) number of the fine dofs should not
806  be smaller than the number of coarse dofs. */
808 {
809 protected:
810  BilinearFormIntegrator *mass_integ; ///< Ownership depends on #own_mass_integ
811  bool own_mass_integ; ///< Ownership flag for #mass_integ
812 
813  OperatorHandle F; ///< Forward, coarse-to-fine, operator
814  OperatorHandle B; ///< Backward, fine-to-coarse, operator
815 
816 public:
818  FiniteElementSpace &fine_fes)
819  : GridTransfer(coarse_fes, fine_fes),
820  mass_integ(NULL), own_mass_integ(false)
821  { }
822 
823  virtual ~InterpolationGridTransfer();
824 
825  /** @brief Assign a mass integrator to be used in the construction of the
826  backward, fine-to-coarse, transfer operator. */
827  void SetMassIntegrator(BilinearFormIntegrator *mass_integ_,
828  bool own_mass_integ_ = true);
829 
830  virtual const Operator &ForwardOperator();
831 
832  virtual const Operator &BackwardOperator();
833 };
834 
835 
836 /** @brief Transfer data between a coarse mesh and an embedded refined mesh
837  using L2 projection. */
838 /** The forward, coarse-to-fine, transfer uses L2 projection. The backward,
839  fine-to-coarse, transfer is defined locally (on a coarse element) as
840  B = (F^t M_f F)^{-1} F^t M_f, where F is the forward transfer matrix, and
841  M_f is the mass matrix on the union of all fine elements comprising the
842  coarse element. Note that the backward transfer operator, B, is a left
843  inverse of the forward transfer operator, F, i.e. B F = I. Both F and B are
844  defined in physical space and, generally, vary between different mesh
845  elements.
846 
847  This class currently only fully supports L2 finite element spaces and fine
848  meshes that are a uniform refinement of the coarse mesh. Generally, the
849  coarse and fine FE spaces can have different orders, however, in order for
850  the backward operator to be well-defined, the number of the fine dofs (in a
851  coarse element) should not be smaller than the number of coarse dofs.
852 
853  If used on H1 finite element spaces, the transfer will be performed locally,
854  and the value of shared (interface) degrees of freedom will be determined by
855  the value of the last transfer to be performed (according to the element
856  numbering in the finite element space). As a consequence, the mass
857  conservation properties for this operator from the L2 case do not carry over
858  to H1 spaces. */
860 {
861 protected:
862  /** Class representing projection operator between a high-order L2 finite
863  element space on a coarse mesh, and a low-order L2 finite element space
864  on a refined mesh (LOR). We assume that the low-order space, fes_lor,
865  lives on a mesh obtained by refining the mesh of the high-order space,
866  fes_ho. */
867  class L2Projection : public Operator
868  {
869  const FiniteElementSpace &fes_ho;
870  const FiniteElementSpace &fes_lor;
871 
872  int ndof_lor, ndof_ho, nref;
873 
874  Table ho2lor;
875 
876  DenseTensor R, P;
877 
878  public:
879  L2Projection(const FiniteElementSpace &fes_ho_,
880  const FiniteElementSpace &fes_lor_);
881  /// Perform the L2 projection onto the LOR space
882  virtual void Mult(const Vector &x, Vector &y) const;
883  /// Perform the mass conservative left-inverse prolongation operation.
884  /// This functionality is also provided as an Operator by L2Prolongation.
885  void Prolongate(const Vector &x, Vector &y) const;
886  virtual ~L2Projection() { }
887  };
888 
889  /** Mass-conservative prolongation operator going in the opposite direction
890  as L2Projection. This operator is a left inverse to the L2Projection. */
891  class L2Prolongation : public Operator
892  {
893  const L2Projection &l2proj;
894 
895  public:
896  L2Prolongation(const L2Projection &l2proj_) : l2proj(l2proj_) { }
897  void Mult(const Vector &x, Vector &y) const
898  {
899  l2proj.Prolongate(x, y);
900  }
901  virtual ~L2Prolongation() { }
902  };
903 
904  L2Projection *F; ///< Forward, coarse-to-fine, operator
905  L2Prolongation *B; ///< Backward, fine-to-coarse, operator
906 
907 public:
909  FiniteElementSpace &fine_fes)
910  : GridTransfer(coarse_fes, fine_fes),
911  F(NULL), B(NULL)
912  { }
913 
914  virtual const Operator &ForwardOperator();
915 
916  virtual const Operator &BackwardOperator();
917 };
918 
919 inline bool UsesTensorBasis(const FiniteElementSpace& fes)
920 {
921  return dynamic_cast<const mfem::TensorBasisElement *>(fes.GetFE(0))!=nullptr;
922 }
923 
924 }
925 
926 #endif
Abstract class for Finite Elements.
Definition: fe.hpp:232
int GetNFbyType(FaceType type) const
Returns the number of faces according to the requested type.
Definition: fespace.hpp:425
const Operator & MakeTrueOperator(FiniteElementSpace &fes_in, FiniteElementSpace &fes_out, const Operator &oper, OperatorHandle &t_oper)
Definition: fespace.cpp:2374
Ordering::Type GetOrdering() const
Return the ordering method.
Definition: fespace.hpp:393
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:381
void GetVertexVDofs(int i, Array< int > &vdofs) const
Definition: fespace.cpp:196
int GetAttribute(int i) const
Definition: fespace.hpp:453
int GetNDofs() const
Returns number of degrees of freedom.
Definition: fespace.hpp:378
int ndofs
Number of degrees of freedom. Number of unknowns is ndofs * vdim.
Definition: fespace.hpp:107
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:1007
const Table & GetBdrElementToDofTable() const
Definition: fespace.hpp:525
void GetElementInteriorVDofs(int i, Array< int > &vdofs) const
Definition: fespace.cpp:202
int DofToVDof(int dof, int vd, int ndofs=-1) const
Definition: fespace.cpp:143
const SparseMatrix * GetConformingProlongation() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:813
virtual void Update(bool want_transform=true)
Definition: fespace.cpp:2057
std::unordered_map< const key_face, Operator *, key_hash > map_L2F
Definition: fespace.hpp:145
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:433
BilinearFormIntegrator * mass_integ
Ownership depends on own_mass_integ.
Definition: fespace.hpp:810
void BuildElementToDofTable() const
Definition: fespace.cpp:214
const IntegrationRule & GetElementIntRule(int idx) const
Get the IntegrationRule associated with mesh element idx.
Definition: fespace.hpp:703
DerefinementOperator(const FiniteElementSpace *f_fes, const FiniteElementSpace *c_fes, BilinearFormIntegrator *mass_integ)
Definition: fespace.cpp:1144
void Prolongate(const Vector &x, Vector &y) const
Definition: fespace.cpp:2711
L2Prolongation(const L2Projection &l2proj_)
Definition: fespace.hpp:896
static const int NumGeom
Definition: geom.hpp:41
FiniteElementSpace & dom_fes
Domain FE space.
Definition: fespace.hpp:718
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:1936
bool Conforming() const
Definition: mesh.hpp:1167
SparseMatrix * RefinementMatrix(int old_ndofs, const Table *old_elem_dof)
Definition: fespace.cpp:1036
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:696
SparseMatrix * DerefinementMatrix(int old_ndofs, const Table *old_elem_dof)
Calculate GridFunction restriction matrix after mesh derefinement.
Definition: fespace.cpp:1303
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:172
Ordering::Type ordering
Definition: fespace.hpp:104
const Geometry::Type geom
Definition: ex1.cpp:40
static int DecodeDof(int dof, double &sign)
Helper to remove encoded sign from a DOF.
Definition: fespace.hpp:162
SparseMatrix * RefinementMatrix_main(const int coarse_ndofs, const Table &coarse_elem_dof, const DenseTensor localP[]) const
Definition: fespace.cpp:950
int GetSize() const
Return the total number of quadrature points.
Definition: fespace.hpp:700
Element::Type GetElementType(int i) const
Returns the type of element i.
Definition: mesh.cpp:4656
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
int VDofToDof(int vdof) const
Definition: fespace.hpp:488
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: fespace.hpp:897
Element::Type GetBdrElementType(int i) const
Returns the type of boundary element i.
Definition: mesh.cpp:4661
int GetFaceOrder(int i) const
Returns the order of the i&#39;th face finite element.
Definition: fespace.cpp:101
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
Definition: fespace.cpp:297
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:99
OperatorHandle L2E_lex
Definition: fespace.hpp:132
Operator::Type oper_type
Desired Operator::Type for the construction of all operators defined by the underlying transfer algor...
Definition: fespace.hpp:724
int GetNumElementInteriorDofs(int i) const
Definition: fespace.hpp:477
int GetBdrAttribute(int i) const
Definition: fespace.hpp:455
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: fespace.cpp:366
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:693
void GetDegenerateFaceDofs(int index, Array< int > &dofs, Geometry::Type master_geom) const
Definition: fespace.cpp:570
int GetBdrElementType(int i) const
Returns the type of boundary element i.
Definition: fespace.hpp:437
int GetNV() const
Returns number of vertices in the mesh.
Definition: fespace.hpp:405
const FiniteElement * GetFaceElement(int i) const
Definition: fespace.cpp:1908
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:2016
FiniteElementCollection * Load(Mesh *m, std::istream &input)
Read a FiniteElementSpace from a stream. The returned FiniteElementCollection is owned by the caller...
Definition: fespace.cpp:2203
virtual const Operator & ForwardOperator()
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
Definition: fespace.cpp:2739
bool UsesTensorBasis(const FiniteElementSpace &fes)
Definition: fespace.hpp:919
const IntegrationRule * int_rule[Geometry::NumGeom]
Definition: fespace.hpp:680
L2ProjectionGridTransfer(FiniteElementSpace &coarse_fes, FiniteElementSpace &fine_fes)
Definition: fespace.hpp:908
const FiniteElementCollection * fec
Associated FE collection (not owned).
Definition: fespace.hpp:96
OperatorHandle Th
Transformation to apply to GridFunctions after space Update().
Definition: fespace.hpp:129
NURBSExtension * GetNURBSext()
Definition: fespace.hpp:298
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:790
void GetLocalDerefinementMatrices(Geometry::Type geom, DenseTensor &localR) const
Definition: fespace.cpp:1277
int GetConformingVSize() const
Definition: fespace.hpp:390
Array< int > dof_elem_array
Definition: fespace.hpp:115
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:385
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:783
void GetFaceInteriorDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1867
virtual int DofForGeometry(Geometry::Type GeomType) const =0
virtual void GetFaceDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1742
virtual ~GridTransfer()
Virtual destructor.
Definition: fespace.hpp:752
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:719
void ConvertFromConformingVDofs(const Array< int > &cdofs, Array< int > &dofs)
Definition: fespace.cpp:422
const Operator * GetElementRestriction(ElementDofOrdering e_ordering) const
Return an Operator that converts L-vectors to E-vectors.
Definition: fespace.cpp:833
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition: mesh.cpp:3943
const SparseMatrix * GetConformingRestriction() const
The returned SparseMatrix is owned by the FiniteElementSpace.
Definition: fespace.cpp:820
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: fespace.cpp:1246
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:408
void GetVertexDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1830
virtual const Operator & BackwardOperator()
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
Definition: fespace.cpp:2527
bool Nonconforming() const
Definition: mesh.hpp:1168
FaceType
Definition: mesh.hpp:42
const FiniteElement * GetEdgeElement(int i) const
Definition: fespace.cpp:1931
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: fespace.cpp:1095
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition: fespace.hpp:417
Data type sparse matrix.
Definition: sparsemat.hpp:40
Native ordering as defined by the FiniteElement.
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:295
OperatorHandle F
Forward, coarse-to-fine, operator.
Definition: fespace.hpp:813
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:2359
Type
Ordering methods:
Definition: fespace.hpp:32
SparseMatrix * cR
Conforming restriction matrix such that cR.cP=I.
Definition: fespace.hpp:125
SparseMatrix * D2Const_GlobalRestrictionMatrix(FiniteElementSpace *cfes)
Definition: fespace.cpp:463
SparseMatrix * D2C_GlobalRestrictionMatrix(FiniteElementSpace *cfes)
Definition: fespace.cpp:431
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: fespace.hpp:311
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:1540
void BuildConformingInterpolation() const
Calculate the cP and cR matrices for a nonconforming mesh.
Definition: fespace.cpp:621
NURBSExtension * StealNURBSext()
Definition: fespace.cpp:1448
int GetElementForDof(int i) const
Definition: fespace.hpp:527
ElementTransformation * GetBdrElementTransformation(int i) const
Returns ElementTransformation for the i-th boundary element.
Definition: fespace.hpp:450
L2Projection(const FiniteElementSpace &fes_ho_, const FiniteElementSpace &fes_lor_)
Definition: fespace.cpp:2569
static void AddDependencies(SparseMatrix &deps, Array< int > &master_dofs, Array< int > &slave_dofs, DenseMatrix &I)
Definition: fespace.cpp:532
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Returns the transformation defining the i-th element in the user-defined variable ElTr...
Definition: fespace.hpp:446
Array< int > dof_ldof_array
Definition: fespace.hpp:115
InterpolationGridTransfer(FiniteElementSpace &coarse_fes, FiniteElementSpace &fine_fes)
Definition: fespace.hpp:817
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:384
bool Conforming() const
Definition: fespace.hpp:301
OperatorHandle bw_t_oper
Backward true-dof operator.
Definition: fespace.hpp:727
SparseMatrix * cP
Definition: fespace.hpp:123
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:370
int GetNF() const
Returns number of faces (i.e. co-dimension 1 entities) in the mesh.
Definition: fespace.hpp:414
void SetOperatorType(Operator::Type type)
Set the desired Operator::Type for the construction of all operators defined by the underlying transf...
Definition: fespace.hpp:759
Array< QuadratureInterpolator * > E2Q_array
Definition: fespace.hpp:148
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
Definition: fespace.cpp:1801
Operator * Ptr() const
Access the underlying Operator pointer.
Definition: handle.hpp:82
OperatorHandle L2E_nat
The element restriction operators, see GetElementRestriction().
Definition: fespace.hpp:132
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:1991
OperatorHandle B
Backward, fine-to-coarse, operator.
Definition: fespace.hpp:814
void GetElementInteriorDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1842
void GetLocalRefinementMatrices(Geometry::Type geom, DenseTensor &localP) const
Definition: fespace.cpp:1012
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:2745
Base class for transfer algorithms that construct transfer Operators between two finite element (FE) ...
Definition: fespace.hpp:715
OperatorHandle fw_t_oper
Forward true-dof operator.
Definition: fespace.hpp:726
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:862
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:2485
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition: fespace.hpp:441
Type
Enumeration defining IDs for some classes derived from Operator.
Definition: operator.hpp:229
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:24
RefinementOperator(const FiniteElementSpace *fespace, Table *old_elem_dof, int old_ndofs)
Definition: fespace.cpp:1054
const NURBSExtension * GetNURBSext() const
Definition: fespace.hpp:297
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:807
int GetLocalDofForDof(int i) const
Definition: fespace.hpp:528
virtual ~FiniteElementSpace()
Definition: fespace.cpp:1942
SparseMatrix * H2L_GlobalRestrictionMatrix(FiniteElementSpace *lfes)
Definition: fespace.cpp:494
virtual ~QuadratureSpace()
Definition: fespace.hpp:697
void GetUpdateOperator(OperatorHandle &T)
Return the update operator in the given OperatorHandle, T.
Definition: fespace.hpp:636
L2Prolongation * B
Backward, fine-to-coarse, operator.
Definition: fespace.hpp:905
int GetOrder(int i) const
Returns the order of the i&#39;th finite element.
Definition: fespace.cpp:95
void Save(std::ostream &out) const
Write the QuadratureSpace to the stream out.
Definition: fespace.cpp:2351
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:785
bool own_mass_integ
Ownership flag for mass_integ.
Definition: fespace.hpp:811
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:208
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:608
L2Projection * F
Forward, coarse-to-fine, operator.
Definition: fespace.hpp:904
int GetElementType(int i) const
Returns the type of element i.
Definition: fespace.hpp:429
FiniteElementSpace(Mesh *mesh, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES)
Definition: fespace.hpp:280
bool Nonconforming() const
Definition: fespace.hpp:302
virtual void GetBdrElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:1659
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:690
int GetNVDofs() const
Number of all scalar vertex dofs.
Definition: fespace.hpp:398
virtual const Operator & ForwardOperator()
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
Definition: fespace.cpp:2494
void ConvertToConformingVDofs(const Array< int > &dofs, Array< int > &cdofs)
Definition: fespace.cpp:414
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
Definition: handle.hpp:112
bool Parallel() const
Definition: fespace.hpp:732
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:773
Mesh * mesh
The mesh that FE space lives on (not owned).
Definition: fespace.hpp:93
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:338
void ReorderElementToDofTable()
Reorder the scalar DOFs based on the element ordering.
Definition: fespace.cpp:243
GridFunction interpolation operator applicable after mesh refinement.
Definition: fespace.hpp:185
Array< FaceQuadratureInterpolator * > E2BFQ_array
Definition: fespace.hpp:150
int GetNConformingDofs() const
Definition: fespace.cpp:827
void GetEdgeInteriorDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1855
const QuadratureInterpolator * GetQuadratureInterpolator(const IntegrationRule &ir) const
Return a QuadratureInterpolator that interpolates E-vectors to quadrature point values and/or derivat...
Definition: fespace.cpp:890
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition: fespace.hpp:65
Array< FaceQuadratureInterpolator * > E2IFQ_array
Definition: fespace.hpp:149
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement associated with i&#39;th element.
Definition: fespace.cpp:1642
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:633
Transfer data between a coarse mesh and an embedded refined mesh using L2 projection.
Definition: fespace.hpp:859
Lexicographic ordering for tensor-product FiniteElements.
FiniteElementSpace()
Default constructor: the object is invalid until initialized using the method Load().
Definition: fespace.cpp:59
QuadratureSpace(Mesh *mesh_, int order_)
Create a QuadratureSpace based on the global rules from IntRules.
Definition: fespace.hpp:691
void Constructor(Mesh *mesh, NURBSExtension *ext, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES)
Help function for constructors + Load().
Definition: fespace.cpp:1402
void SetUpdateOperatorOwner(bool own)
Set the ownership of the update operator: if set to false, the Operator returned by GetUpdateOperator...
Definition: fespace.hpp:643
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:184
void GetEdgeVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th edge.
Definition: fespace.cpp:190
static void ListToMarker(const Array< int > &list, int marker_size, Array< int > &marker, int mark_val=-1)
Definition: fespace.cpp:403
virtual void Mult(const Vector &x, Vector &y) const
Perform the L2 projection onto the LOR space.
Definition: fespace.cpp:2685
std::size_t operator()(const key_face &k) const
Definition: fespace.hpp:137
Vector data type.
Definition: vector.hpp:48
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:395
std::tuple< bool, ElementDofOrdering, FaceType, L2FaceValues > key_face
The face restriction operators, see GetFaceRestriction().
Definition: fespace.hpp:134
int GetNFbyType(FaceType type) const
Returns the number of faces according to the requested type.
Definition: mesh.cpp:3969
virtual void UpdatesFinished()
Free the GridFunction update operator (if any), to save memory.
Definition: fespace.hpp:654
const Table & GetElementToDofTable() const
Definition: fespace.hpp:524
void SetUpdateOperatorType(Operator::Type tid)
Specify the Operator::Type to be used by the update operators.
Definition: fespace.hpp:651
static bool DofFinalizable(int dof, const Array< bool > &finalized, const SparseMatrix &deps)
Definition: fespace.cpp:556
Class representing the storage layout of a QuadratureFunction.
Definition: fespace.hpp:671
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement for the i&#39;th boundary element.
Definition: fespace.cpp:1882
NURBSExtension * NURBSext
Definition: fespace.hpp:117
void GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:178
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:315
long GetSequence() const
Return update counter (see Mesh::sequence)
Definition: fespace.hpp:657
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:725
void Save(std::ostream &out) const
Definition: fespace.cpp:2134
int GetAttribute(int i) const
Return the attribute of element i.
Definition: mesh.hpp:1001
Class representing a function through its values (scalar or vector) at quadrature points...
Definition: gridfunc.hpp:508
int GetNFDofs() const
Number of all scalar face-interior dofs.
Definition: fespace.hpp:402
int GetNEDofs() const
Number of all scalar edge-interior dofs.
Definition: fespace.hpp:400
void GetElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of element i.
Definition: mesh.hpp:834
static void AdjustVDofs(Array< int > &vdofs)
Definition: fespace.cpp:159
void SetType(Operator::Type tid)
Invoke Clear() and set a new type id.
Definition: handle.hpp:124
void DofsToVDofs(Array< int > &dofs, int ndofs=-1) const
Definition: fespace.cpp:107
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:919
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:289