MFEM  v3.4
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, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
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 <iostream>
20 
21 namespace mfem
22 {
23 
24 /** @brief The ordering method used when the number of unknowns per mesh node
25  (vector dimension) is bigger than 1. */
26 class Ordering
27 {
28 public:
29  /// %Ordering methods:
30  enum Type
31  {
32  byNODES, /**< loop first over the nodes (inner loop) then over the vector
33  dimension (outer loop); symbolically it can be represented
34  as: XXX...,YYY...,ZZZ... */
35  byVDIM /**< loop first over the vector dimension (inner loop) then over
36  the nodes (outer loop); symbolically it can be represented
37  as: XYZ,XYZ,XYZ,... */
38  };
39 
40  template <Type Ord>
41  static inline int Map(int ndofs, int vdim, int dof, int vd);
42 
43  template <Type Ord>
44  static void DofsToVDofs(int ndofs, int vdim, Array<int> &dofs);
45 };
46 
47 template <> inline int
48 Ordering::Map<Ordering::byNODES>(int ndofs, int vdim, int dof, int vd)
49 {
50  MFEM_ASSERT(dof < ndofs && -1-dof < ndofs && 0 <= vd && vd < vdim, "");
51  return (dof >= 0) ? dof+ndofs*vd : dof-ndofs*vd;
52 }
53 
54 template <> inline int
55 Ordering::Map<Ordering::byVDIM>(int ndofs, int vdim, int dof, int vd)
56 {
57  MFEM_ASSERT(dof < ndofs && -1-dof < ndofs && 0 <= vd && vd < vdim, "");
58  return (dof >= 0) ? vd+vdim*dof : -1-(vd+vdim*(-1-dof));
59 }
60 
61 
62 class NURBSExtension;
63 
64 /** @brief Class FiniteElementSpace - responsible for providing FEM view of the
65  mesh, mainly managing the set of degrees of freedom. */
67 {
68 protected:
69  /// The mesh that FE space lives on (not owned).
71 
72  /// Associated FE collection (not owned).
74 
75  /// %Vector dimension (number of unknowns per degree of freedom).
76  int vdim;
77 
78  /** Type of ordering of the vector dofs when #vdim > 1.
79  - Ordering::byNODES - first nodes, then vector dimension,
80  - Ordering::byVDIM - first vector dimension, then nodes */
82 
83  /// Number of degrees of freedom. Number of unknowns is #ndofs * #vdim.
84  int ndofs;
85 
87  int *fdofs, *bdofs;
88 
89  mutable Table *elem_dof; // if NURBS FE space, not owned; otherwise, owned.
90  Table *bdrElem_dof; // used only with NURBS FE spaces; not owned.
91 
93 
95  int own_ext;
96 
97  /** Matrix representing the prolongation from the global conforming dofs to
98  a set of intermediate partially conforming dofs, e.g. the dofs associated
99  with a "cut" space on a non-conforming mesh. */
100  mutable SparseMatrix *cP; // owned
101  /// Conforming restriction matrix such that cR.cP=I.
102  mutable SparseMatrix *cR; // owned
103  mutable bool cP_is_set;
104 
105  /// Transformation to apply to GridFunctions after space Update().
107 
108  long sequence; // should match Mesh::GetSequence
109 
110  void UpdateNURBS();
111 
112  void Construct();
113  void Destroy();
114 
115  void BuildElementToDofTable() const;
116 
117  /// Helper to remove encoded sign from a DOF
118  static inline int DecodeDof(int dof, double& sign)
119  { return (dof >= 0) ? (sign = 1, dof) : (sign = -1, (-1 - dof)); }
120 
121  /// Helper to get vertex, edge or face DOFs (entity=0,1,2 resp.).
122  void GetEntityDofs(int entity, int index, Array<int> &dofs) const;
123 
124  /// Calculate the cP and cR matrices for a nonconforming mesh.
125  void BuildConformingInterpolation() const;
126 
127  static void AddDependencies(SparseMatrix& deps, Array<int>& master_dofs,
128  Array<int>& slave_dofs, DenseMatrix& I);
129 
130  static bool DofFinalizable(int dof, const Array<bool>& finalized,
131  const SparseMatrix& deps);
132 
133  void MakeVDimMatrix(SparseMatrix &mat) const;
134 
135  /// GridFunction interpolation operator applicable after mesh refinement.
137  {
138  const FiniteElementSpace* fespace;
139  DenseTensor localP;
140  Table* old_elem_dof; // Owned.
141 
142  public:
143  /** Construct the operator based on the elem_dof table of the original
144  (coarse) space. The class takes ownership of the table. */
145  RefinementOperator(const FiniteElementSpace* fespace,
146  Table *old_elem_dof/*takes ownership*/, int old_ndofs);
147  RefinementOperator(const FiniteElementSpace *fespace,
148  const FiniteElementSpace *coarse_fes);
149  virtual void Mult(const Vector &x, Vector &y) const;
150  virtual ~RefinementOperator();
151  };
152 
153  // This method makes the same assumptions as the method:
154  // void GetLocalRefinementMatrices(
155  // const FiniteElementSpace &coarse_fes, DenseTensor &localP) const
156  // which is defined below. It also assumes that the coarse fes and this have
157  // the same vector dimension, vdim.
158  SparseMatrix *RefinementMatrix_main(const int coarse_ndofs,
159  const Table &coarse_elem_dof,
160  const DenseTensor &localP) const;
161 
162  void GetLocalRefinementMatrices(DenseTensor &localP) const;
163  void GetLocalDerefinementMatrices(DenseTensor &localR) const;
164 
165  /** Calculate explicit GridFunction interpolation matrix (after mesh
166  refinement). NOTE: consider using the RefinementOperator class instead
167  of the fully assembled matrix, which can take a lot of memory. */
168  SparseMatrix* RefinementMatrix(int old_ndofs, const Table* old_elem_dof);
169 
170  /// Calculate GridFunction restriction matrix after mesh derefinement.
171  SparseMatrix* DerefinementMatrix(int old_ndofs, const Table* old_elem_dof);
172 
173  // This method assumes that this->mesh is a refinement of coarse_fes->mesh
174  // and that the CoarseFineTransformations of this->mesh are set accordingly.
175  // Another assumption is that the FEs of this use the same MapType as the FEs
176  // of coarse_fes. Finally, it assumes that this->mesh and coarse_fes->mesh
177  // are NOT mixed meshes, and that the spaces this and coarse_fes are NOT
178  // variable-order spaces.
179  void GetLocalRefinementMatrices(const FiniteElementSpace &coarse_fes,
180  DenseTensor &localP) const;
181 
182  /// Help function for constructors + Load().
183  void Constructor(Mesh *mesh, NURBSExtension *ext,
185  int vdim = 1, int ordering = Ordering::byNODES);
186 
187 public:
188  /** @brief Default constructor: the object is invalid until initialized using
189  the method Load(). */
191 
192  /** @brief Copy constructor: deep copy all data from @a orig except the Mesh,
193  the FiniteElementCollection, ans some derived data. */
194  /** If the @a mesh or @a fec pointers are NULL (default), then the new
195  FiniteElementSpace will reuse the respective pointers from @a orig. If
196  any of these pointers is not NULL, the given pointer will be used instead
197  of the one used by @a orig.
198 
199  @note The objects pointed to by the @a mesh and @a fec parameters must be
200  either the same objects as the ones used by @a orig, or copies of them.
201  Otherwise, the behavior is undefined.
202 
203  @note Derived data objects, such as the conforming prolongation and
204  restriction matrices, and the update operator, will not be copied, even
205  if they are created in the @a orig object. */
206  FiniteElementSpace(const FiniteElementSpace &orig, Mesh *mesh = NULL,
207  const FiniteElementCollection *fec = NULL);
208 
211  int vdim = 1, int ordering = Ordering::byNODES)
212  { Constructor(mesh, NULL, fec, vdim, ordering); }
213 
214  /// Construct a NURBS FE space based on the given NURBSExtension, @a ext.
215  /** @note If the pointer @a ext is NULL, this constructor is equivalent to
216  the standard constructor with the same arguments minus the
217  NURBSExtension, @a ext. */
220  int vdim = 1, int ordering = Ordering::byNODES)
221  { Constructor(mesh, ext, fec, vdim, ordering); }
222 
223  /// Returns the mesh
224  inline Mesh *GetMesh() const { return mesh; }
225 
226  const NURBSExtension *GetNURBSext() const { return NURBSext; }
229 
230  bool Conforming() const { return mesh->Conforming(); }
231  bool Nonconforming() const { return mesh->Nonconforming(); }
232 
234  const SparseMatrix *GetConformingRestriction() const;
235 
236  virtual const Operator *GetProlongationMatrix() const
237  { return GetConformingProlongation(); }
238  virtual const SparseMatrix *GetRestrictionMatrix() const
239  { return GetConformingRestriction(); }
240 
241  /// Returns vector dimension.
242  inline int GetVDim() const { return vdim; }
243 
244  /// Returns the order of the i'th finite element
245  int GetOrder(int i) const;
246  /// Returns the order of the i'th face finite element
247  int GetFaceOrder(int i) const;
248 
249  /// Returns number of degrees of freedom.
250  inline int GetNDofs() const { return ndofs; }
251 
252  /// Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
253  inline int GetVSize() const { return vdim * ndofs; }
254 
255  /// Return the number of vector true (conforming) dofs.
256  virtual int GetTrueVSize() const { return GetConformingVSize(); }
257 
258  /// Returns the number of conforming ("true") degrees of freedom
259  /// (if the space is on a nonconforming mesh with hanging nodes).
260  int GetNConformingDofs() const;
261 
262  int GetConformingVSize() const { return vdim * GetNConformingDofs(); }
263 
264  /// Return the ordering method.
265  inline Ordering::Type GetOrdering() const { return ordering; }
266 
267  const FiniteElementCollection *FEColl() const { return fec; }
268 
269  int GetNVDofs() const { return nvdofs; }
270  int GetNEDofs() const { return nedofs; }
271  int GetNFDofs() const { return nfdofs; }
272 
273  /// Returns number of vertices in the mesh.
274  inline int GetNV() const { return mesh->GetNV(); }
275 
276  /// Returns number of elements in the mesh.
277  inline int GetNE() const { return mesh->GetNE(); }
278 
279  /// Returns number of faces (i.e. co-dimension 1 entities) in the mesh.
280  /** The co-dimension 1 entities are those that have dimension 1 less than the
281  mesh dimension, e.g. for a 2D mesh, the faces are the 1D entities, i.e.
282  the edges. */
283  inline int GetNF() const { return mesh->GetNumFaces(); }
284 
285  /// Returns number of boundary elements in the mesh.
286  inline int GetNBE() const { return mesh->GetNBE(); }
287 
288  /// Returns the type of element i.
289  inline int GetElementType(int i) const
290  { return mesh->GetElementType(i); }
291 
292  /// Returns the vertices of element i.
293  inline void GetElementVertices(int i, Array<int> &vertices) const
294  { mesh->GetElementVertices(i, vertices); }
295 
296  /// Returns the type of boundary element i.
297  inline int GetBdrElementType(int i) const
298  { return mesh->GetBdrElementType(i); }
299 
300  /// Returns ElementTransformation for the @a i-th element.
302  { return mesh->GetElementTransformation(i); }
303 
304  /** @brief Returns the transformation defining the @a i-th element in the
305  user-defined variable @a ElTr. */
307  { mesh->GetElementTransformation(i, ElTr); }
308 
309  /// Returns ElementTransformation for the @a i-th boundary element.
311  { return mesh->GetBdrElementTransformation(i); }
312 
313  int GetAttribute(int i) const { return mesh->GetAttribute(i); }
314 
315  int GetBdrAttribute(int i) const { return mesh->GetBdrAttribute(i); }
316 
317  /// Returns indexes of degrees of freedom in array dofs for i'th element.
318  virtual void GetElementDofs(int i, Array<int> &dofs) const;
319 
320  /// Returns indexes of degrees of freedom for i'th boundary element.
321  virtual void GetBdrElementDofs(int i, Array<int> &dofs) const;
322 
323  /** Returns the indexes of the degrees of freedom for i'th face
324  including the dofs for the edges and the vertices of the face. */
325  virtual void GetFaceDofs(int i, Array<int> &dofs) const;
326 
327  /** Returns the indexes of the degrees of freedom for i'th edge
328  including the dofs for the vertices of the edge. */
329  void GetEdgeDofs(int i, Array<int> &dofs) const;
330 
331  void GetVertexDofs(int i, Array<int> &dofs) const;
332 
333  void GetElementInteriorDofs(int i, Array<int> &dofs) const;
334 
335  void GetFaceInteriorDofs(int i, Array<int> &dofs) const;
336 
337  int GetNumElementInteriorDofs(int i) const
339 
340  void GetEdgeInteriorDofs(int i, Array<int> &dofs) const;
341 
342  void DofsToVDofs(Array<int> &dofs, int ndofs = -1) const;
343 
344  void DofsToVDofs(int vd, Array<int> &dofs, int ndofs = -1) const;
345 
346  int DofToVDof(int dof, int vd, int ndofs = -1) const;
347 
348  int VDofToDof(int vdof) const
349  { return (ordering == Ordering::byNODES) ? (vdof%ndofs) : (vdof/vdim); }
350 
351  static void AdjustVDofs(Array<int> &vdofs);
352 
353  /// Returns indexes of degrees of freedom in array dofs for i'th element.
354  void GetElementVDofs(int i, Array<int> &vdofs) const;
355 
356  /// Returns indexes of degrees of freedom for i'th boundary element.
357  void GetBdrElementVDofs(int i, Array<int> &vdofs) const;
358 
359  /// Returns indexes of degrees of freedom for i'th face element (2D and 3D).
360  void GetFaceVDofs(int i, Array<int> &vdofs) const;
361 
362  /// Returns indexes of degrees of freedom for i'th edge.
363  void GetEdgeVDofs(int i, Array<int> &vdofs) const;
364 
365  void GetVertexVDofs(int i, Array<int> &vdofs) const;
366 
367  void GetElementInteriorVDofs(int i, Array<int> &vdofs) const;
368 
369  void GetEdgeInteriorVDofs(int i, Array<int> &vdofs) const;
370 
372 
373  /** @brief Reorder the scalar DOFs based on the element ordering.
374 
375  The new ordering is constructed as follows: 1) loop over all elements as
376  ordered in the Mesh; 2) for each element, assign new indices to all of
377  its current DOFs that are still unassigned; the new indices we assign are
378  simply the sequence `0,1,2,...`; if there are any signed DOFs their sign
379  is preserved. */
381 
382  void BuildDofToArrays();
383 
384  const Table &GetElementToDofTable() const { return *elem_dof; }
385  const Table &GetBdrElementToDofTable() const { return *bdrElem_dof; }
386 
387  int GetElementForDof(int i) const { return dof_elem_array[i]; }
388  int GetLocalDofForDof(int i) const { return dof_ldof_array[i]; }
389 
390  /// Returns pointer to the FiniteElement associated with i'th element.
391  const FiniteElement *GetFE(int i) const;
392 
393  /// Returns pointer to the FiniteElement for the i'th boundary element.
394  const FiniteElement *GetBE(int i) const;
395 
396  const FiniteElement *GetFaceElement(int i) const;
397 
398  const FiniteElement *GetEdgeElement(int i) const;
399 
400  /// Return the trace element from element 'i' to the given 'geom_type'
401  const FiniteElement *GetTraceElement(int i, int geom_type) const;
402 
403  /** Mark degrees of freedom associated with boundary elements with
404  the specified boundary attributes (marked in 'bdr_attr_is_ess').
405  For spaces with 'vdim' > 1, the 'component' parameter can be used
406  to restricts the marked vDOFs to the specified component. */
407  virtual void GetEssentialVDofs(const Array<int> &bdr_attr_is_ess,
408  Array<int> &ess_vdofs,
409  int component = -1) const;
410 
411  /** Get a list of essential true dofs, ess_tdof_list, corresponding to the
412  boundary attributes marked in the array bdr_attr_is_ess.
413  For spaces with 'vdim' > 1, the 'component' parameter can be used
414  to restricts the marked tDOFs to the specified component. */
415  virtual void GetEssentialTrueDofs(const Array<int> &bdr_attr_is_ess,
416  Array<int> &ess_tdof_list,
417  int component = -1);
418 
419  /// Convert a Boolean marker array to a list containing all marked indices.
420  static void MarkerToList(const Array<int> &marker, Array<int> &list);
421 
422  /** Convert an array of indices (list) to a Boolean marker array where all
423  indices in the list are marked with the given value and the rest are set
424  to zero. */
425  static void ListToMarker(const Array<int> &list, int marker_size,
426  Array<int> &marker, int mark_val = -1);
427 
428  /** For a partially conforming FE space, convert a marker array (nonzero
429  entries are true) on the partially conforming dofs to a marker array on
430  the conforming dofs. A conforming dofs is marked iff at least one of its
431  dependent dofs is marked. */
432  void ConvertToConformingVDofs(const Array<int> &dofs, Array<int> &cdofs);
433 
434  /** For a partially conforming FE space, convert a marker array (nonzero
435  entries are true) on the conforming dofs to a marker array on the
436  (partially conforming) dofs. A dof is marked iff it depends on a marked
437  conforming dofs, where dependency is defined by the ConformingRestriction
438  matrix; in other words, a dof is marked iff it corresponds to a marked
439  conforming dof. */
440  void ConvertFromConformingVDofs(const Array<int> &cdofs, Array<int> &dofs);
441 
442  /** Generate the global restriction matrix from a discontinuous
443  FE space to the continuous FE space of the same polynomial degree. */
445 
446  /** Generate the global restriction matrix from a discontinuous
447  FE space to the piecewise constant FE space. */
449 
450  /** Construct the restriction matrix from the FE space given by
451  (*this) to the lower degree FE space given by (*lfes) which
452  is defined on the same mesh. */
454 
455  /** @brief Construct and return an Operator that can be used to transfer
456  GridFunction data from @a coarse_fes, defined on a coarse mesh, to @a
457  this FE space, defined on a refined mesh. */
458  /** It is assumed that the mesh of this FE space is a refinement of the mesh
459  of @a coarse_fes and the CoarseFineTransformations returned by the method
460  Mesh::GetRefinementTransforms() of the refined mesh are set accordingly.
461  The Operator::Type of @a T can be set to request an Operator of the set
462  type. Currently, only Operator::MFEM_SPARSEMAT and Operator::ANY_TYPE
463  (matrix-free) are supported. When Operator::ANY_TYPE is requested, the
464  choice of the particular Operator sub-class is left to the method. This
465  method also works in parallel because the transfer operator is local to
466  the MPI task when the input is a synchronized ParGridFunction. */
467  void GetTransferOperator(const FiniteElementSpace &coarse_fes,
468  OperatorHandle &T) const;
469 
470  /** @brief Construct and return an Operator that can be used to transfer
471  true-dof data from @a coarse_fes, defined on a coarse mesh, to @a this FE
472  space, defined on a refined mesh.
473 
474  This method calls GetTransferOperator() and multiplies the result by the
475  prolongation operator of @a coarse_fes on the right, and by the
476  restriction operator of this FE space on the left.
477 
478  The Operator::Type of @a T can be set to request an Operator of the set
479  type. In serial, the supported types are: Operator::MFEM_SPARSEMAT and
480  Operator::ANY_TYPE (matrix-free). In parallel, the supported types are:
481  Operator::Hypre_ParCSR and Operator::ANY_TYPE. Any other type is treated
482  as Operator::ANY_TYPE: the operator representation choice is made by this
483  method. */
484  virtual void GetTrueTransferOperator(const FiniteElementSpace &coarse_fes,
485  OperatorHandle &T) const;
486 
487  /** Reflect changes in the mesh: update number of DOFs, etc. Also, calculate
488  GridFunction transformation operator (unless want_transform is false).
489  Safe to call multiple times, does nothing if space already up to date. */
490  virtual void Update(bool want_transform = true);
491 
492  /// Get the GridFunction update operator.
493  const Operator* GetUpdateOperator() { Update(); return Th.Ptr(); }
494 
495  /// Return the update operator in the given OperatorHandle, @a T.
497 
498  /** @brief Set the ownership of the update operator: if set to false, the
499  Operator returned by GetUpdateOperator() must be deleted outside the
500  FiniteElementSpace. */
501  /** The update operator ownership is automatically reset to true when a new
502  update operator is created by the Update() method. */
503  void SetUpdateOperatorOwner(bool own) { Th.SetOperatorOwner(own); }
504 
505  /// Specify the Operator::Type to be used by the update operators.
506  /** The default type is Operator::ANY_TYPE which leaves the choice to this
507  class. The other currently supported option is Operator::MFEM_SPARSEMAT
508  which is only guaranteed to be honored for a refinement update operator.
509  Any other type will be treated as Operator::ANY_TYPE.
510  @note This operation destroys the current update operator (if owned). */
512 
513  /// Free the GridFunction update operator (if any), to save memory.
514  virtual void UpdatesFinished() { Th.Clear(); }
515 
516  /// Return update counter (see Mesh::sequence)
517  long GetSequence() const { return sequence; }
518 
519  void Save(std::ostream &out) const;
520 
521  /** @brief Read a FiniteElementSpace from a stream. The returned
522  FiniteElementCollection is owned by the caller. */
523  FiniteElementCollection *Load(Mesh *m, std::istream &input);
524 
525  virtual ~FiniteElementSpace();
526 };
527 
528 
529 /// Class representing the storage layout of a QuadratureFunction.
530 /** Multiple QuadratureFunction%s can share the same QuadratureSpace. */
532 {
533 protected:
534  friend class QuadratureFunction; // Uses the element_offsets.
535 
537  int order;
538  int size;
539 
541  int *element_offsets; // scalar offsets; size = number of elements + 1
542 
543  // protected functions
544 
545  // Assuming mesh and order are set, construct the members: int_rule,
546  // element_offsets, and size.
547  void Construct();
548 
549 public:
550  /// Create a QuadratureSpace based on the global rules from #IntRules.
551  QuadratureSpace(Mesh *mesh_, int order_)
552  : mesh(mesh_), order(order_) { Construct(); }
553 
554  /// Read a QuadratureSpace from the stream @a in.
555  QuadratureSpace(Mesh *mesh_, std::istream &in);
556 
557  virtual ~QuadratureSpace() { delete [] element_offsets; }
558 
559  /// Return the total number of quadrature points.
560  int GetSize() const { return size; }
561 
562  /// Get the IntegrationRule associated with mesh element @a idx.
563  const IntegrationRule &GetElementIntRule(int idx) const
564  { return *int_rule[mesh->GetElementBaseGeometry(idx)]; }
565 
566  /// Write the QuadratureSpace to the stream @a out.
567  void Save(std::ostream &out) const;
568 };
569 
570 }
571 
572 #endif
Abstract class for Finite Elements.
Definition: fe.hpp:140
Ordering::Type GetOrdering() const
Return the ordering method.
Definition: fespace.hpp:265
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:253
void GetVertexVDofs(int i, Array< int > &vdofs) const
Definition: fespace.cpp:195
int GetAttribute(int i) const
Definition: fespace.hpp:313
int GetNDofs() const
Returns number of degrees of freedom.
Definition: fespace.hpp:250
int ndofs
Number of degrees of freedom. Number of unknowns is ndofs * vdim.
Definition: fespace.hpp:84
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:83
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
Definition: mesh.hpp:868
const Table & GetBdrElementToDofTable() const
Definition: fespace.hpp:385
void GetElementInteriorVDofs(int i, Array< int > &vdofs) const
Definition: fespace.cpp:201
int DofToVDof(int dof, int vd, int ndofs=-1) const
Definition: fespace.cpp:142
void GetLocalRefinementMatrices(DenseTensor &localP) const
Definition: fespace.cpp:832
const SparseMatrix * GetConformingProlongation() const
Definition: fespace.cpp:761
virtual void Update(bool want_transform=true)
Definition: fespace.cpp:1719
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:293
void BuildElementToDofTable() const
Definition: fespace.cpp:213
const IntegrationRule & GetElementIntRule(int idx) const
Get the IntegrationRule associated with mesh element idx.
Definition: fespace.hpp:563
static const int NumGeom
Definition: geom.hpp:34
bool Conforming() const
Definition: mesh.hpp:1013
SparseMatrix * RefinementMatrix(int old_ndofs, const Table *old_elem_dof)
Definition: fespace.cpp:855
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:621
SparseMatrix * DerefinementMatrix(int old_ndofs, const Table *old_elem_dof)
Calculate GridFunction restriction matrix after mesh derefinement.
Definition: fespace.cpp:1011
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:171
Ordering::Type ordering
Definition: fespace.hpp:81
static int DecodeDof(int dof, double &sign)
Helper to remove encoded sign from a DOF.
Definition: fespace.hpp:118
int GetSize() const
Return the total number of quadrature points.
Definition: fespace.hpp:560
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
int VDofToDof(int vdof) const
Definition: fespace.hpp:348
int GetFaceOrder(int i) const
Returns the order of the i&#39;th face finite element.
Definition: fespace.cpp:100
virtual void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_vdofs, int component=-1) const
Definition: fespace.cpp:296
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:76
void GetEntityDofs(int entity, int index, Array< int > &dofs) const
Helper to get vertex, edge or face DOFs (entity=0,1,2 resp.).
Definition: fespace.cpp:565
int GetNumElementInteriorDofs(int i) const
Definition: fespace.hpp:337
const FiniteElement * GetTraceElement(int i, int geom_type) const
Return the trace element from element &#39;i&#39; to the given &#39;geom_type&#39;.
Definition: fespace.cpp:1625
int GetBdrAttribute(int i) const
Definition: fespace.hpp:315
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: fespace.cpp:365
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:618
int GetBdrElementType(int i) const
Returns the type of boundary element i.
Definition: fespace.hpp:297
int GetNV() const
Returns number of vertices in the mesh.
Definition: fespace.hpp:274
const FiniteElement * GetFaceElement(int i) const
Definition: fespace.cpp:1597
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:1678
FiniteElementCollection * Load(Mesh *m, std::istream &input)
Read a FiniteElementSpace from a stream. The returned FiniteElementCollection is owned by the caller...
Definition: fespace.cpp:1857
const IntegrationRule * int_rule[Geometry::NumGeom]
Definition: fespace.hpp:540
const FiniteElementCollection * fec
Associated FE collection (not owned).
Definition: fespace.hpp:73
OperatorHandle Th
Transformation to apply to GridFunctions after space Update().
Definition: fespace.hpp:106
NURBSExtension * GetNURBSext()
Definition: fespace.hpp:227
int GetConformingVSize() const
Definition: fespace.hpp:262
Array< int > dof_elem_array
Definition: fespace.hpp:92
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:384
void GetFaceInteriorDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1556
virtual void GetFaceDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1431
ElementTransformation * GetBdrElementTransformation(int i)
Returns the transformation defining the i-th boundary element.
Definition: mesh.cpp:356
void ConvertFromConformingVDofs(const Array< int > &cdofs, Array< int > &dofs)
Definition: fespace.cpp:420
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition: mesh.cpp:3374
const SparseMatrix * GetConformingRestriction() const
Definition: fespace.cpp:768
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:277
void GetVertexDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1519
bool Nonconforming() const
Definition: mesh.hpp:1014
const FiniteElement * GetEdgeElement(int i) const
Definition: fespace.cpp:1620
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: fespace.cpp:896
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition: fespace.hpp:286
Data type sparse matrix.
Definition: sparsemat.hpp:38
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:224
Type
Ordering methods:
Definition: fespace.hpp:30
SparseMatrix * cR
Conforming restriction matrix such that cR.cP=I.
Definition: fespace.hpp:102
SparseMatrix * D2Const_GlobalRestrictionMatrix(FiniteElementSpace *cfes)
Definition: fespace.cpp:461
SparseMatrix * D2C_GlobalRestrictionMatrix(FiniteElementSpace *cfes)
Definition: fespace.cpp:429
virtual const Operator * GetProlongationMatrix() const
Definition: fespace.hpp:236
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:1233
void BuildConformingInterpolation() const
Calculate the cP and cR matrices for a nonconforming mesh.
Definition: fespace.cpp:576
NURBSExtension * StealNURBSext()
Definition: fespace.cpp:1138
int GetElementForDof(int i) const
Definition: fespace.hpp:387
ElementTransformation * GetBdrElementTransformation(int i) const
Returns ElementTransformation for the i-th boundary element.
Definition: fespace.hpp:310
static void AddDependencies(SparseMatrix &deps, Array< int > &master_dofs, Array< int > &slave_dofs, DenseMatrix &I)
Definition: fespace.cpp:526
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Returns the transformation defining the i-th element in the user-defined variable ElTr...
Definition: fespace.hpp:306
Array< int > dof_ldof_array
Definition: fespace.hpp:92
int GetElementType(int i) const
Returns the type of element i.
Definition: mesh.cpp:4029
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:256
bool Conforming() const
Definition: fespace.hpp:230
SparseMatrix * cP
Definition: fespace.hpp:100
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:242
int GetNF() const
Returns number of faces (i.e. co-dimension 1 entities) in the mesh.
Definition: fespace.hpp:283
void GetEdgeDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1490
Operator * Ptr() const
Access the underlying Operator pointer.
Definition: handle.hpp:82
int GetBdrElementType(int i) const
Returns the type of boundary element i.
Definition: mesh.cpp:4034
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:1659
void GetElementInteriorDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1531
static void DofsToVDofs(int ndofs, int vdim, Array< int > &dofs)
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition: fespace.hpp:301
Type
Enumeration defining IDs for some classes derived from Operator.
Definition: operator.hpp:124
RefinementOperator(const FiniteElementSpace *fespace, Table *old_elem_dof, int old_ndofs)
Definition: fespace.cpp:867
const NURBSExtension * GetNURBSext() const
Definition: fespace.hpp:226
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:66
int GetLocalDofForDof(int i) const
Definition: fespace.hpp:388
virtual ~FiniteElementSpace()
Definition: fespace.cpp:1631
SparseMatrix * H2L_GlobalRestrictionMatrix(FiniteElementSpace *lfes)
Definition: fespace.cpp:492
virtual ~QuadratureSpace()
Definition: fespace.hpp:557
void GetUpdateOperator(OperatorHandle &T)
Return the update operator in the given OperatorHandle, T.
Definition: fespace.hpp:496
int GetOrder(int i) const
Returns the order of the i&#39;th finite element.
Definition: fespace.cpp:94
void Save(std::ostream &out) const
Write the QuadratureSpace to the stream out.
Definition: fespace.cpp:1998
void Clear()
Clear the OperatorHandle, deleting the held Operator (if owned), while leaving the type id unchanged...
Definition: handle.hpp:110
void GetLocalDerefinementMatrices(DenseTensor &localR) const
Definition: fespace.cpp:961
void MakeVDimMatrix(SparseMatrix &mat) const
Definition: fespace.cpp:733
The ordering method used when the number of unknowns per mesh node (vector dimension) is bigger than ...
Definition: fespace.hpp:26
void GetEdgeInteriorVDofs(int i, Array< int > &vdofs) const
Definition: fespace.cpp:207
int GetElementBaseGeometry(int i=0) const
Definition: mesh.hpp:688
int GetElementType(int i) const
Returns the type of element i.
Definition: fespace.hpp:289
FiniteElementSpace(Mesh *mesh, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES)
Definition: fespace.hpp:209
bool Nonconforming() const
Definition: fespace.hpp:231
virtual void GetBdrElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:1348
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:615
int GetNVDofs() const
Definition: fespace.hpp:269
void ConvertToConformingVDofs(const Array< int > &dofs, Array< int > &cdofs)
Definition: fespace.cpp:412
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
Definition: handle.hpp:106
Mesh * mesh
The mesh that FE space lives on (not owned).
Definition: fespace.hpp:70
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:279
void ReorderElementToDofTable()
Reorder the scalar DOFs based on the element ordering.
Definition: fespace.cpp:242
GridFunction interpolation operator applicable after mesh refinement.
Definition: fespace.hpp:136
int GetNConformingDofs() const
Definition: fespace.cpp:775
void GetEdgeInteriorDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1544
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement associated with i&#39;th element.
Definition: fespace.cpp:1335
const Operator * GetUpdateOperator()
Get the GridFunction update operator.
Definition: fespace.hpp:493
FiniteElementSpace()
Default constructor: the object is invalid until initialized using the method Load().
Definition: fespace.cpp:58
QuadratureSpace(Mesh *mesh_, int order_)
Create a QuadratureSpace based on the global rules from IntRules.
Definition: fespace.hpp:551
void Constructor(Mesh *mesh, NURBSExtension *ext, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES)
Help function for constructors + Load().
Definition: fespace.cpp:1092
void SetUpdateOperatorOwner(bool own)
Set the ownership of the update operator: if set to false, the Operator returned by GetUpdateOperator...
Definition: fespace.hpp:503
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:183
void GetEdgeVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th edge.
Definition: fespace.cpp:189
static void ListToMarker(const Array< int > &list, int marker_size, Array< int > &marker, int mark_val=-1)
Definition: fespace.cpp:401
Vector data type.
Definition: vector.hpp:48
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:267
virtual void UpdatesFinished()
Free the GridFunction update operator (if any), to save memory.
Definition: fespace.hpp:514
const Table & GetElementToDofTable() const
Definition: fespace.hpp:384
void SetUpdateOperatorType(Operator::Type tid)
Specify the Operator::Type to be used by the update operators.
Definition: fespace.hpp:511
SparseMatrix * RefinementMatrix_main(const int coarse_ndofs, const Table &coarse_elem_dof, const DenseTensor &localP) const
Definition: fespace.cpp:781
static bool DofFinalizable(int dof, const Array< bool > &finalized, const SparseMatrix &deps)
Definition: fespace.cpp:550
Class representing the storage layout of a QuadratureFunction.
Definition: fespace.hpp:531
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement for the i&#39;th boundary element.
Definition: fespace.cpp:1571
NURBSExtension * NURBSext
Definition: fespace.hpp:94
void GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:177
virtual int DofForGeometry(int GeomType) const =0
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:64
Abstract operator.
Definition: operator.hpp:21
virtual const SparseMatrix * GetRestrictionMatrix() const
Definition: fespace.hpp:238
long GetSequence() const
Return update counter (see Mesh::sequence)
Definition: fespace.hpp:517
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:638
void Save(std::ostream &out) const
Definition: fespace.cpp:1796
int GetAttribute(int i) const
Return the attribute of element i.
Definition: mesh.hpp:862
Class representing a function through its values (scalar or vector) at quadrature points...
Definition: gridfunc.hpp:401
int GetNFDofs() const
Definition: fespace.hpp:271
int GetNEDofs() const
Definition: fespace.hpp:270
void GetElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of element i.
Definition: mesh.hpp:695
static void AdjustVDofs(Array< int > &vdofs)
Definition: fespace.cpp:158
void SetType(Operator::Type tid)
Invoke Clear() and set a new type id.
Definition: handle.hpp:118
void DofsToVDofs(Array< int > &dofs, int ndofs=-1) const
Definition: fespace.cpp:106
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:218