MFEM  v3.3.2
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 /** The ordering method used when the number of unknowns per mesh
25  node (vector dimension) is bigger than 1. */
26 class Ordering
27 {
28 public:
29  /** Ordering methods:
30  byNODES - loop first over the nodes (inner loop) then over the vector dimension (outer loop),
31  byVDIM - loop first over the vector dimension (inner loop) then over the nodes (outer loop) */
32  enum Type { byNODES, byVDIM };
33 
34  template <Type Ord>
35  static inline int Map(int ndofs, int vdim, int dof, int vd);
36 
37  template <Type Ord>
38  static void DofsToVDofs(int ndofs, int vdim, Array<int> &dofs);
39 };
40 
41 template <> inline int
42 Ordering::Map<Ordering::byNODES>(int ndofs, int vdim, int dof, int vd)
43 {
44  MFEM_ASSERT(dof < ndofs && -1-dof < ndofs && 0 <= vd && vd < vdim, "");
45  return (dof >= 0) ? dof+ndofs*vd : dof-ndofs*vd;
46 }
47 
48 template <> inline int
49 Ordering::Map<Ordering::byVDIM>(int ndofs, int vdim, int dof, int vd)
50 {
51  MFEM_ASSERT(dof < ndofs && -1-dof < ndofs && 0 <= vd && vd < vdim, "");
52  return (dof >= 0) ? vd+vdim*dof : -1-(vd+vdim*(-1-dof));
53 }
54 
55 
56 class NURBSExtension;
57 
58 /** Class FiniteElementSpace - responsible for providing FEM view of the mesh
59  (mainly managing the set of degrees of freedom). */
61 {
62 protected:
63  /// The mesh that FE space lives on.
65 
67 
68  /// Vector dimension (number of unknowns per degree of freedom).
69  int vdim;
70 
71  /** Type of ordering of dofs.
72  Ordering::byNODES - first nodes, then vector dimension,
73  Ordering::byVDIM - first vector dimension, then nodes */
75 
76  /// Number of degrees of freedom. Number of unknowns are ndofs*vdim.
77  int ndofs;
78 
80  int *fdofs, *bdofs;
81 
82  mutable Table *elem_dof;
84 
86 
88  int own_ext;
89 
90  /** Matrix representing the prolongation from the global conforming dofs to
91  a set of intermediate partially conforming dofs, e.g. the dofs associated
92  with a "cut" space on a non-conforming mesh. */
93  mutable SparseMatrix *cP;
94  /// Conforming restriction matrix such that cR.cP=I.
95  mutable SparseMatrix *cR;
96  mutable bool cP_is_set;
97 
98  /// Transformation to apply to GridFunctions after space Update().
100  bool own_T;
101 
102  long sequence; // should match Mesh::GetSequence
103 
104  void UpdateNURBS();
105 
106  void Construct();
107  void Destroy();
108 
109  void BuildElementToDofTable() const;
110 
111  /** This is a helper function to get edge (type == 0) or face (type == 1)
112  DOFs. The function is aware of ghost edges/faces in parallel, for which
113  an empty DOF list is returned. */
114  void GetEdgeFaceDofs(int type, int index, Array<int> &dofs) const;
115 
116  /// Calculate the cP and cR matrices for a nonconforming mesh.
117  void GetConformingInterpolation() const;
118 
119  void MakeVDimMatrix(SparseMatrix &mat) const;
120 
121  /// Calculate GridFunction interpolation matrix after mesh refinement.
122  SparseMatrix* RefinementMatrix(int old_ndofs, const Table* old_elem_dof);
123 
125  DenseTensor &localR);
126 
127  /// Calculate GridFunction restriction matrix after mesh derefinement.
128  SparseMatrix* DerefinementMatrix(int old_ndofs, const Table* old_elem_dof);
129 
130 
131 public:
133  int vdim = 1, int ordering = Ordering::byNODES);
134 
135  /// Returns the mesh
136  inline Mesh *GetMesh() const { return mesh; }
137 
140 
141  bool Conforming() const { return mesh->Conforming(); }
142  bool Nonconforming() const { return mesh->Nonconforming(); }
143 
145  const SparseMatrix *GetConformingRestriction() const;
146 
148  { return GetConformingProlongation(); }
150  { return GetConformingRestriction(); }
151 
152  /// Returns vector dimension.
153  inline int GetVDim() const { return vdim; }
154 
155  /// Returns the order of the i'th finite element
156  int GetOrder(int i) const;
157  /// Returns the order of the i'th face finite element
158  int GetFaceOrder(int i) const;
159 
160  /// Returns number of degrees of freedom.
161  inline int GetNDofs() const { return ndofs; }
162 
163  inline int GetVSize() const { return vdim * ndofs; }
164 
165  /// Return the number of vector true (conforming) dofs.
166  virtual int GetTrueVSize() const { return GetConformingVSize(); }
167 
168  /// Returns the number of conforming ("true") degrees of freedom
169  /// (if the space is on a nonconforming mesh with hanging nodes).
170  int GetNConformingDofs() const;
171 
172  int GetConformingVSize() const { return vdim * GetNConformingDofs(); }
173 
174  /// Return the ordering method.
175  inline Ordering::Type GetOrdering() const { return ordering; }
176 
177  const FiniteElementCollection *FEColl() const { return fec; }
178 
179  int GetNVDofs() const { return nvdofs; }
180  int GetNEDofs() const { return nedofs; }
181  int GetNFDofs() const { return nfdofs; }
182 
183  /// Returns number of elements in the mesh.
184  inline int GetNE() const { return mesh->GetNE(); }
185 
186  /// Returns number of nodes in the mesh.
187  inline int GetNV() const { return mesh->GetNV(); }
188 
189  /// Returns number of boundary elements in the mesh.
190  inline int GetNBE() const { return mesh->GetNBE(); }
191 
192  /// Returns the type of element i.
193  inline int GetElementType(int i) const
194  { return mesh->GetElementType(i); }
195 
196  /// Returns the vertices of element i.
197  inline void GetElementVertices(int i, Array<int> &vertices) const
198  { mesh->GetElementVertices(i, vertices); }
199 
200  /// Returns the type of boundary element i.
201  inline int GetBdrElementType(int i) const
202  { return mesh->GetBdrElementType(i); }
203 
204  /// Returns ElementTransformation for the i'th element.
206  { return mesh->GetElementTransformation(i); }
207 
208  /** Returns the transformation defining the i-th element in the user-defined
209  variable. */
211  { mesh->GetElementTransformation(i, ElTr); }
212 
213  /// Returns ElementTransformation for the i'th boundary element.
215  { return mesh->GetBdrElementTransformation(i); }
216 
217  int GetAttribute(int i) const { return mesh->GetAttribute(i); }
218 
219  int GetBdrAttribute(int i) const { return mesh->GetBdrAttribute(i); }
220 
221  /// Returns indexes of degrees of freedom in array dofs for i'th element.
222  virtual void GetElementDofs(int i, Array<int> &dofs) const;
223 
224  /// Returns indexes of degrees of freedom for i'th boundary element.
225  virtual void GetBdrElementDofs(int i, Array<int> &dofs) const;
226 
227  /** Returns the indexes of the degrees of freedom for i'th face
228  including the dofs for the edges and the vertices of the face. */
229  virtual void GetFaceDofs(int i, Array<int> &dofs) const;
230 
231  /** Returns the indexes of the degrees of freedom for i'th edge
232  including the dofs for the vertices of the edge. */
233  void GetEdgeDofs(int i, Array<int> &dofs) const;
234 
235  void GetVertexDofs(int i, Array<int> &dofs) const;
236 
237  void GetElementInteriorDofs(int i, Array<int> &dofs) const;
238 
239  void GetFaceInteriorDofs(int i, Array<int> &dofs) const;
240 
241  int GetNumElementInteriorDofs(int i) const
243 
244  void GetEdgeInteriorDofs(int i, Array<int> &dofs) const;
245 
246  void DofsToVDofs(Array<int> &dofs, int ndofs = -1) const;
247 
248  void DofsToVDofs(int vd, Array<int> &dofs, int ndofs = -1) const;
249 
250  int DofToVDof(int dof, int vd, int ndofs = -1) const;
251 
252  int VDofToDof(int vdof) const
253  { return (ordering == Ordering::byNODES) ? (vdof%ndofs) : (vdof/vdim); }
254 
255  static void AdjustVDofs(Array<int> &vdofs);
256 
257  /// Returns indexes of degrees of freedom in array dofs for i'th element.
258  void GetElementVDofs(int i, Array<int> &vdofs) const;
259 
260  /// Returns indexes of degrees of freedom for i'th boundary element.
261  void GetBdrElementVDofs(int i, Array<int> &vdofs) const;
262 
263  /// Returns indexes of degrees of freedom for i'th face element (2D and 3D).
264  void GetFaceVDofs(int i, Array<int> &vdofs) const;
265 
266  /// Returns indexes of degrees of freedom for i'th edge.
267  void GetEdgeVDofs(int i, Array<int> &vdofs) const;
268 
269  void GetVertexVDofs(int i, Array<int> &vdofs) const;
270 
271  void GetElementInteriorVDofs(int i, Array<int> &vdofs) const;
272 
273  void GetEdgeInteriorVDofs(int i, Array<int> &vdofs) const;
274 
276 
277  /** @brief Reorder the scalar DOFs based on the element ordering.
278 
279  The new ordering is constructed as follows: 1) loop over all elements as
280  ordered in the Mesh; 2) for each element, assign new indices to all of
281  its current DOFs that are still unassigned; the new indices we assign are
282  simply the sequence `0,1,2,...`; if there are any signed DOFs their sign
283  is preserved. */
285 
286  void BuildDofToArrays();
287 
288  const Table &GetElementToDofTable() const { return *elem_dof; }
289  const Table &GetBdrElementToDofTable() const { return *bdrElem_dof; }
290 
291  int GetElementForDof(int i) const { return dof_elem_array[i]; }
292  int GetLocalDofForDof(int i) const { return dof_ldof_array[i]; }
293 
294  /// Returns pointer to the FiniteElement associated with i'th element.
295  const FiniteElement *GetFE(int i) const;
296 
297  /// Returns pointer to the FiniteElement for the i'th boundary element.
298  const FiniteElement *GetBE(int i) const;
299 
300  const FiniteElement *GetFaceElement(int i) const;
301 
302  const FiniteElement *GetEdgeElement(int i) const;
303 
304  /// Return the trace element from element 'i' to the given 'geom_type'
305  const FiniteElement *GetTraceElement(int i, int geom_type) const;
306 
307  /** Mark degrees of freedom associated with boundary elements with
308  the specified boundary attributes (marked in 'bdr_attr_is_ess').
309  For spaces with 'vdim' > 1, the 'component' parameter can be used
310  to restricts the marked vDOFs to the specified component. */
311  virtual void GetEssentialVDofs(const Array<int> &bdr_attr_is_ess,
312  Array<int> &ess_vdofs,
313  int component = -1) const;
314 
315  /** Get a list of essential true dofs, ess_tdof_list, corresponding to the
316  boundary attributes marked in the array bdr_attr_is_ess.
317  For spaces with 'vdim' > 1, the 'component' parameter can be used
318  to restricts the marked tDOFs to the specified component. */
319  virtual void GetEssentialTrueDofs(const Array<int> &bdr_attr_is_ess,
320  Array<int> &ess_tdof_list,
321  int component = -1);
322 
323  /// Convert a Boolean marker array to a list containing all marked indices.
324  static void MarkerToList(const Array<int> &marker, Array<int> &list);
325 
326  /** Convert an array of indices (list) to a Boolean marker array where all
327  indices in the list are marked with the given value and the rest are set
328  to zero. */
329  static void ListToMarker(const Array<int> &list, int marker_size,
330  Array<int> &marker, int mark_val = -1);
331 
332  /** For a partially conforming FE space, convert a marker array (nonzero
333  entries are true) on the partially conforming dofs to a marker array on
334  the conforming dofs. A conforming dofs is marked iff at least one of its
335  dependent dofs is marked. */
336  void ConvertToConformingVDofs(const Array<int> &dofs, Array<int> &cdofs);
337 
338  /** For a partially conforming FE space, convert a marker array (nonzero
339  entries are true) on the conforming dofs to a marker array on the
340  (partially conforming) dofs. A dof is marked iff it depends on a marked
341  conforming dofs, where dependency is defined by the ConformingRestriction
342  matrix; in other words, a dof is marked iff it corresponds to a marked
343  conforming dof. */
344  void ConvertFromConformingVDofs(const Array<int> &cdofs, Array<int> &dofs);
345 
346  /** Generate the global restriction matrix from a discontinuous
347  FE space to the continuous FE space of the same polynomial degree. */
349 
350  /** Generate the global restriction matrix from a discontinuous
351  FE space to the piecewise constant FE space. */
353 
354  /** Construct the restriction matrix from the FE space given by
355  (*this) to the lower degree FE space given by (*lfes) which
356  is defined on the same mesh. */
358 
359  /** Reflect changes in the mesh: update number of DOFs, etc. Also, calculate
360  GridFunction transformation matrix (unless want_transform is false).
361  Safe to call multiple times, does nothing if space already up to date. */
362  virtual void Update(bool want_transform = true);
363 
364  /// Get the GridFunction update matrix.
365  const Operator* GetUpdateOperator() { Update(); return T; }
366 
367  /** @brief Set the ownership of the update operator: if set to false, the
368  Operator returned by GetUpdateOperator() must be deleted outside the
369  FiniteElementSpace. */
370  void SetUpdateOperatorOwner(bool own) { own_T = own; }
371 
372  /// Free GridFunction transformation matrix (if any), to save memory.
373  virtual void UpdatesFinished() { if (own_T) { delete T; } T = NULL; }
374 
375  /// Return update counter (see Mesh::sequence)
376  long GetSequence() const { return sequence; }
377 
378  void Save(std::ostream &out) const;
379 
380  virtual ~FiniteElementSpace();
381 };
382 
383 
384 /// Class representing the storage layout of a QuadratureFunction.
385 /** Multiple QuadratureFunction%s can share the same QuadratureSpace. */
387 {
388 protected:
389  friend class QuadratureFunction; // Uses the element_offsets.
390 
392  int order;
393  int size;
394 
396  int *element_offsets; // scalar offsets; size = number of elements + 1
397 
398  // protected functions
399 
400  // Assuming mesh and order are set, construct the members: int_rule,
401  // element_offsets, and size.
402  void Construct();
403 
404 public:
405  /// Create a QuadratureSpace based on the global rules from #IntRules.
406  QuadratureSpace(Mesh *mesh_, int order_)
407  : mesh(mesh_), order(order_) { Construct(); }
408 
409  /// Read a QuadratureSpace from the stream @a in.
410  QuadratureSpace(Mesh *mesh_, std::istream &in);
411 
412  virtual ~QuadratureSpace() { delete [] element_offsets; }
413 
414  /// Return the total number of quadrature points.
415  int GetSize() { return size; }
416 
417  /// Get the IntegrationRule associated with mesh element @a idx.
419  { return *int_rule[mesh->GetElementBaseGeometry(idx)]; }
420 
421  /// Write the QuadratureSpace to the stream @a out.
422  void Save(std::ostream &out) const;
423 };
424 
425 }
426 
427 #endif
Abstract class for Finite Elements.
Definition: fe.hpp:47
Ordering::Type GetOrdering() const
Return the ordering method.
Definition: fespace.hpp:175
int GetVSize() const
Definition: fespace.hpp:163
void GetVertexVDofs(int i, Array< int > &vdofs) const
Definition: fespace.cpp:157
int GetAttribute(int i) const
Definition: fespace.hpp:217
int GetNDofs() const
Returns number of degrees of freedom.
Definition: fespace.hpp:161
int ndofs
Number of degrees of freedom. Number of unknowns are ndofs*vdim.
Definition: fespace.hpp:77
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:855
const Table & GetBdrElementToDofTable() const
Definition: fespace.hpp:289
void GetElementInteriorVDofs(int i, Array< int > &vdofs) const
Definition: fespace.cpp:163
int DofToVDof(int dof, int vd, int ndofs=-1) const
Definition: fespace.cpp:104
const SparseMatrix * GetConformingProlongation() const
Definition: fespace.cpp:725
virtual void Update(bool want_transform=true)
Definition: fespace.cpp:1494
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:197
void BuildElementToDofTable() const
Definition: fespace.cpp:175
static const int NumGeom
Definition: geom.hpp:34
bool Conforming() const
Definition: mesh.hpp:998
SparseMatrix * RefinementMatrix(int old_ndofs, const Table *old_elem_dof)
Calculate GridFunction interpolation matrix after mesh refinement.
Definition: fespace.cpp:745
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:617
SparseMatrix * DerefinementMatrix(int old_ndofs, const Table *old_elem_dof)
Calculate GridFunction restriction matrix after mesh derefinement.
Definition: fespace.cpp:873
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:133
Ordering::Type ordering
Definition: fespace.hpp:74
int GetSize()
Return the total number of quadrature points.
Definition: fespace.hpp:415
const Geometry::Type geom
int VDofToDof(int vdof) const
Definition: fespace.hpp:252
int GetFaceOrder(int i) const
Returns the order of the i&#39;th face finite element.
Definition: fespace.cpp:62
virtual void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_vdofs, int component=-1) const
Definition: fespace.cpp:258
int vdim
Vector dimension (number of unknowns per degree of freedom).
Definition: fespace.hpp:69
int GetNumElementInteriorDofs(int i) const
Definition: fespace.hpp:241
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:1460
int GetBdrAttribute(int i) const
Definition: fespace.hpp:219
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: fespace.cpp:327
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:614
int GetBdrElementType(int i) const
Returns the type of boundary element i.
Definition: fespace.hpp:201
int GetNV() const
Returns number of nodes in the mesh.
Definition: fespace.hpp:187
const FiniteElement * GetFaceElement(int i) const
Definition: fespace.cpp:1432
const IntegrationRule * int_rule[Geometry::NumGeom]
Definition: fespace.hpp:395
const FiniteElementCollection * fec
Definition: fespace.hpp:66
NURBSExtension * GetNURBSext()
Definition: fespace.hpp:138
int GetConformingVSize() const
Definition: fespace.hpp:172
void GetEdgeFaceDofs(int type, int index, Array< int > &dofs) const
Definition: fespace.cpp:528
Array< int > dof_elem_array
Definition: fespace.hpp:85
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:346
void GetFaceInteriorDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1391
virtual void GetFaceDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1267
ElementTransformation * GetBdrElementTransformation(int i)
Returns the transformation defining the i-th boundary element.
Definition: mesh.cpp:332
void ConvertFromConformingVDofs(const Array< int > &cdofs, Array< int > &dofs)
Definition: fespace.cpp:382
const SparseMatrix * GetConformingRestriction() const
Definition: fespace.cpp:732
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:184
void GetVertexDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1355
bool Nonconforming() const
Definition: mesh.hpp:999
const FiniteElement * GetEdgeElement(int i) const
Definition: fespace.cpp:1455
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition: fespace.hpp:190
Data type sparse matrix.
Definition: sparsemat.hpp:38
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:136
Operator * T
Transformation to apply to GridFunctions after space Update().
Definition: fespace.hpp:99
SparseMatrix * cR
Conforming restriction matrix such that cR.cP=I.
Definition: fespace.hpp:95
SparseMatrix * D2Const_GlobalRestrictionMatrix(FiniteElementSpace *cfes)
Definition: fespace.cpp:423
SparseMatrix * D2C_GlobalRestrictionMatrix(FiniteElementSpace *cfes)
Definition: fespace.cpp:391
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:1072
NURBSExtension * StealNURBSext()
Definition: fespace.cpp:981
int GetElementForDof(int i) const
Definition: fespace.hpp:291
ElementTransformation * GetBdrElementTransformation(int i) const
Returns ElementTransformation for the i&#39;th boundary element.
Definition: fespace.hpp:214
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: fespace.hpp:210
Array< int > dof_ldof_array
Definition: fespace.hpp:85
int GetElementType(int i) const
Returns the type of element i.
Definition: mesh.cpp:3989
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:166
bool Conforming() const
Definition: fespace.hpp:141
SparseMatrix * cP
Definition: fespace.hpp:93
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:153
virtual const SparseMatrix * GetRestrictionMatrix()
Definition: fespace.hpp:149
const IntegrationRule & GetElementIntRule(int idx)
Get the IntegrationRule associated with mesh element idx.
Definition: fespace.hpp:418
void GetEdgeDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1326
int GetBdrElementType(int i) const
Returns the type of boundary element i.
Definition: mesh.cpp:3994
void GetElementInteriorDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1367
static void DofsToVDofs(int ndofs, int vdim, Array< int > &dofs)
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i&#39;th element.
Definition: fespace.hpp:205
int GetLocalDofForDof(int i) const
Definition: fespace.hpp:292
virtual ~FiniteElementSpace()
Definition: fespace.cpp:1466
SparseMatrix * H2L_GlobalRestrictionMatrix(FiniteElementSpace *lfes)
Definition: fespace.cpp:454
virtual ~QuadratureSpace()
Definition: fespace.hpp:412
int GetOrder(int i) const
Returns the order of the i&#39;th finite element.
Definition: fespace.cpp:56
void Save(std::ostream &out) const
Write the QuadratureSpace to the stream out.
Definition: fespace.cpp:1613
void MakeVDimMatrix(SparseMatrix &mat) const
Definition: fespace.cpp:698
void GetEdgeInteriorVDofs(int i, Array< int > &vdofs) const
Definition: fespace.cpp:169
virtual const Operator * GetProlongationMatrix()
Definition: fespace.hpp:147
int GetElementBaseGeometry(int i=0) const
Definition: mesh.hpp:684
int GetElementType(int i) const
Returns the type of element i.
Definition: fespace.hpp:193
FiniteElementSpace(Mesh *mesh, const FiniteElementCollection *fec, int vdim=1, int ordering=Ordering::byNODES)
Definition: fespace.cpp:930
bool Nonconforming() const
Definition: fespace.hpp:142
virtual void GetBdrElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:1184
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:611
int GetNVDofs() const
Definition: fespace.hpp:179
void ConvertToConformingVDofs(const Array< int > &dofs, Array< int > &cdofs)
Definition: fespace.cpp:374
Mesh * mesh
The mesh that FE space lives on.
Definition: fespace.hpp:64
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:255
void ReorderElementToDofTable()
Reorder the scalar DOFs based on the element ordering.
Definition: fespace.cpp:204
int GetNConformingDofs() const
Definition: fespace.cpp:739
void GetEdgeInteriorDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1379
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement associated with i&#39;th element.
Definition: fespace.cpp:1171
const Operator * GetUpdateOperator()
Get the GridFunction update matrix.
Definition: fespace.hpp:365
QuadratureSpace(Mesh *mesh_, int order_)
Create a QuadratureSpace based on the global rules from IntRules.
Definition: fespace.hpp:406
void SetUpdateOperatorOwner(bool own)
Set the ownership of the update operator: if set to false, the Operator returned by GetUpdateOperator...
Definition: fespace.hpp:370
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:145
void GetEdgeVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th edge.
Definition: fespace.cpp:151
static void ListToMarker(const Array< int > &list, int marker_size, Array< int > &marker, int mark_val=-1)
Definition: fespace.cpp:363
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:177
virtual void UpdatesFinished()
Free GridFunction transformation matrix (if any), to save memory.
Definition: fespace.hpp:373
const Table & GetElementToDofTable() const
Definition: fespace.hpp:288
void GetConformingInterpolation() const
Calculate the cP and cR matrices for a nonconforming mesh.
Definition: fespace.cpp:542
Defines the coarse-fine transformations of all fine elements.
Definition: ncmesh.hpp:52
Class representing the storage layout of a QuadratureFunction.
Definition: fespace.hpp:386
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement for the i&#39;th boundary element.
Definition: fespace.cpp:1406
NURBSExtension * NURBSext
Definition: fespace.hpp:87
void GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:139
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
long GetSequence() const
Return update counter (see Mesh::sequence)
Definition: fespace.hpp:376
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:634
void Save(std::ostream &out) const
Definition: fespace.cpp:1558
int GetAttribute(int i) const
Return the attribute of element i.
Definition: mesh.hpp:849
Class representing a function through its values (scalar or vector) at quadrature points...
Definition: gridfunc.hpp:354
int GetNFDofs() const
Definition: fespace.hpp:181
int GetNEDofs() const
Definition: fespace.hpp:180
void GetElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of element i.
Definition: mesh.hpp:691
static void AdjustVDofs(Array< int > &vdofs)
Definition: fespace.cpp:120
void GetLocalDerefinementMatrices(int geom, const CoarseFineTransformations &dt, DenseTensor &localR)
Definition: fespace.cpp:827
void DofsToVDofs(Array< int > &dofs, int ndofs=-1) const
Definition: fespace.cpp:68