MFEM  v4.6.0
Finite element discretization library
qspace.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_QSPACE
13 #define MFEM_QSPACE
14 
15 #include "../config/config.hpp"
16 #include "fespace.hpp"
17 #include <unordered_map>
18 
19 namespace mfem
20 {
21 
22 /// Abstract base class for QuadratureSpace and FaceQuadratureSpace.
23 /** This class represents the storage layout for QuadratureFunction%s, that may
24  be defined either on mesh elements or mesh faces. */
26 {
27 protected:
28  friend class QuadratureFunction; // Uses the offsets.
29 
30  Mesh &mesh; ///< The underlying mesh.
31  int order; ///< The order of integration rule.
32  int size; ///< Total number of quadrature points.
33 
34  /// @brief Entity quadrature point offset array, of size num_entities + 1.
35  ///
36  /// The quadrature point values for entity i are stored in the indices between
37  /// offsets[i] and offsets[i+1].
39  /// The quadrature rules used for each geometry type.
41 
42  /// Protected constructor. Used by derived classes.
43  QuadratureSpaceBase(Mesh &mesh_, int order_ = 0)
44  : mesh(mesh_), order(order_) { }
45 
46  /// Protected constructor. Used by derived classes.
48  const IntegrationRule &ir);
49 
50  /// Fill the @ref int_rule array for each geometry type using @ref order.
51  void ConstructIntRules(int dim);
52 
53 public:
54  /// Return the total number of quadrature points.
55  int GetSize() const { return size; }
56 
57  /// Return the order of the quadrature rule(s) used by all elements.
58  int GetOrder() const { return order; }
59 
60  /// Return the number of entities.
61  int GetNE() const { return offsets.Size() - 1; }
62 
63  /// Returns the mesh.
64  inline Mesh *GetMesh() const { return &mesh; }
65 
66  /// Get the (element or face) transformation of entity @a idx.
67  virtual ElementTransformation *GetTransformation(int idx) = 0;
68 
69  /// Return the geometry type of entity (element or face) @a idx.
70  virtual Geometry::Type GetGeometry(int idx) const = 0;
71 
72  /// Return the IntegrationRule associated with entity @a idx.
73  const IntegrationRule &GetIntRule(int idx) const
74  { return *int_rule[GetGeometry(idx)]; }
75 
76  /// @brief Returns the permuted index of the @a iq quadrature point in entity
77  /// @a idx.
78  ///
79  /// For tensor-product faces, returns the lexicographic index of the
80  /// quadrature point, oriented relative to "element 1". For QuadratureSpace%s
81  /// defined on elements (not faces), the permutation is trivial, and this
82  /// returns @a iq.
83  virtual int GetPermutedIndex(int idx, int iq) const = 0;
84 
85  /// @brief Returns the index in the quadrature space of the entity associated
86  /// with the transformation @a T.
87  ///
88  /// For a QuadratureSpace defined on elements, this just returns the element
89  /// index. For FaceQuadratureSpace, the returned index depends on the chosen
90  /// FaceType.
91  virtual int GetEntityIndex(const ElementTransformation &T) const = 0;
92 
93  /// Write the QuadratureSpace to the stream @a out.
94  virtual void Save(std::ostream &out) const = 0;
95 
96  virtual ~QuadratureSpaceBase() { }
97 };
98 
99 /// Class representing the storage layout of a QuadratureFunction.
100 /** Multiple QuadratureFunction%s can share the same QuadratureSpace. */
102 {
103 protected:
104  void ConstructOffsets();
105  void Construct();
106 public:
107  /// Create a QuadratureSpace based on the global rules from #IntRules.
108  QuadratureSpace(Mesh *mesh_, int order_)
109  : QuadratureSpaceBase(*mesh_, order_) { Construct(); }
110 
111  /// @brief Create a QuadratureSpace with an IntegrationRule, valid only when
112  /// the mesh has one element type.
113  QuadratureSpace(Mesh &mesh_, const IntegrationRule &ir);
114 
115  /// Read a QuadratureSpace from the stream @a in.
116  QuadratureSpace(Mesh *mesh_, std::istream &in);
117 
118  /// Returns number of elements in the mesh.
119  inline int GetNE() const { return mesh.GetNE(); }
120 
121  /// Returns the element transformation of element @a idx.
123  { return mesh.GetElementTransformation(idx); }
124 
125  /// Returns the geometry type of element @a idx.
126  Geometry::Type GetGeometry(int idx) const override
127  { return mesh.GetElementGeometry(idx); }
128 
129  /// Get the IntegrationRule associated with mesh element @a idx.
130  const IntegrationRule &GetElementIntRule(int idx) const
131  { return *int_rule[mesh.GetElementBaseGeometry(idx)]; }
132 
133  /// @brief Returns the permuted index of the @a iq quadrature point in entity
134  /// @a idx.
135  ///
136  /// The member function QuadratureSpace::GetPermutedIndex always returns @a
137  /// iq, the permutation is only nontrivial for FaceQuadratureSpace.
138  int GetPermutedIndex(int idx, int iq) const override { return iq; }
139 
140  /// Returns the element index of @a T.
141  int GetEntityIndex(const ElementTransformation &T) const override { return T.ElementNo; }
142 
143  /// Write the QuadratureSpace to the stream @a out.
144  void Save(std::ostream &out) const override;
145 };
146 
147 /// Class representing the storage layout of a FaceQuadratureFunction.
148 /** FaceQuadratureSpace is defined on either the interior or boundary faces
149  of a mesh, depending on the provided FaceType. */
151 {
152  FaceType face_type; ///< Is the space defined on interior or boundary faces?
153  const int num_faces; ///< Number of faces.
154 
155  /// Map from boundary or interior face indices to mesh face indices.
156  Array<int> face_indices;
157 
158  /// Inverse of the map @a face_indices.
159  std::unordered_map<int,int> face_indices_inv;
160 
161  void ConstructOffsets();
162  void Construct();
163 
164 public:
165  /// Create a FaceQuadratureSpace based on the global rules from #IntRules.
166  FaceQuadratureSpace(Mesh &mesh_, int order_, FaceType face_type_);
167 
168  /// @brief Create a FaceQuadratureSpace with an IntegrationRule, valid only
169  /// when the mesh has one type of face geometry.
170  FaceQuadratureSpace(Mesh &mesh_, const IntegrationRule &ir,
171  FaceType face_type_);
172 
173  /// Returns number of faces in the mesh.
174  inline int GetNumFaces() const { return num_faces; }
175 
176  /// Returns the face type (boundary or interior).
177  FaceType GetFaceType() const { return face_type; }
178 
179  /// Returns the face transformation of face @a idx.
181  { return mesh.GetFaceTransformation(face_indices[idx]); }
182 
183  /// Returns the geometry type of face @a idx.
184  Geometry::Type GetGeometry(int idx) const override
185  { return mesh.GetFaceGeometry(face_indices[idx]); }
186 
187  /// Get the IntegrationRule associated with mesh element @a idx.
188  const IntegrationRule &GetFaceIntRule(int idx) const
189  { return *int_rule[GetGeometry(idx)]; }
190 
191  /// @brief Returns the permuted index of the @a iq quadrature point in entity
192  /// @a idx.
193  ///
194  /// For tensor-product faces, returns the lexicographic index of the
195  /// quadrature point, oriented relative to "element 1".
196  int GetPermutedIndex(int idx, int iq) const override;
197 
198  /// @brief Returns the index associated with the face described by @a T.
199  ///
200  /// The index may differ from the mesh face or boundary element index
201  /// depending on the FaceType used to construct the FaceQuadratureSpace.
202  int GetEntityIndex(const ElementTransformation &T) const override;
203 
204  /// Write the FaceQuadratureSpace to the stream @a out.
205  void Save(std::ostream &out) const override;
206 };
207 
208 }
209 
210 #endif
ElementTransformation * GetTransformation(int idx) override
Returns the element transformation of element idx.
Definition: qspace.hpp:122
Geometry::Type GetGeometry(int idx) const override
Returns the geometry type of face idx.
Definition: qspace.hpp:184
int GetNE() const
Return the number of entities.
Definition: qspace.hpp:61
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:96
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:1242
static const int NumGeom
Definition: geom.hpp:42
int GetNE() const
Returns number of elements in the mesh.
Definition: qspace.hpp:119
virtual int GetPermutedIndex(int idx, int iq) const =0
Returns the permuted index of the iq quadrature point in entity idx.
Geometry::Type GetFaceGeometry(int i) const
Return the Geometry::Type associated with face i.
Definition: mesh.cpp:1421
const IntegrationRule & GetElementIntRule(int idx) const
Get the IntegrationRule associated with mesh element idx.
Definition: qspace.hpp:130
QuadratureSpaceBase(Mesh &mesh_, int order_=0)
Protected constructor. Used by derived classes.
Definition: qspace.hpp:43
int size
Total number of quadrature points.
Definition: qspace.hpp:32
int GetPermutedIndex(int idx, int iq) const override
Returns the permuted index of the iq quadrature point in entity idx.
Definition: qspace.cpp:148
virtual ElementTransformation * GetTransformation(int idx)=0
Get the (element or face) transformation of entity idx.
FaceType
Definition: mesh.hpp:45
void GetFaceTransformation(int i, IsoparametricTransformation *FTr)
Definition: mesh.cpp:504
virtual Geometry::Type GetGeometry(int idx) const =0
Return the geometry type of entity (element or face) idx.
Geometry::Type GetGeometry(int idx) const override
Returns the geometry type of element idx.
Definition: qspace.hpp:126
virtual ~QuadratureSpaceBase()
Definition: qspace.hpp:96
int GetNumFaces() const
Returns number of faces in the mesh.
Definition: qspace.hpp:174
Mesh * GetMesh() const
Returns the mesh.
Definition: qspace.hpp:64
void Save(std::ostream &out) const override
Write the QuadratureSpace to the stream out.
Definition: qspace.cpp:90
virtual void Save(std::ostream &out) const =0
Write the QuadratureSpace to the stream out.
Array< int > offsets
Entity quadrature point offset array, of size num_entities + 1.
Definition: qspace.hpp:38
Abstract base class for QuadratureSpace and FaceQuadratureSpace.
Definition: qspace.hpp:25
Geometry::Type GetElementGeometry(int i) const
Definition: mesh.hpp:1228
FaceQuadratureSpace(Mesh &mesh_, int order_, FaceType face_type_)
Create a FaceQuadratureSpace based on the global rules from IntRules.
Definition: qspace.cpp:97
const IntegrationRule & GetFaceIntRule(int idx) const
Get the IntegrationRule associated with mesh element idx.
Definition: qspace.hpp:188
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
const IntegrationRule & GetIntRule(int idx) const
Return the IntegrationRule associated with entity idx.
Definition: qspace.hpp:73
Class representing the storage layout of a FaceQuadratureFunction.
Definition: qspace.hpp:150
ElementTransformation * GetTransformation(int idx) override
Returns the face transformation of face idx.
Definition: qspace.hpp:180
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:1086
int GetEntityIndex(const ElementTransformation &T) const override
Returns the index associated with the face described by T.
Definition: qspace.cpp:165
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:355
Mesh & mesh
The underlying mesh.
Definition: qspace.hpp:30
const IntegrationRule * int_rule[Geometry::NumGeom]
The quadrature rules used for each geometry type.
Definition: qspace.hpp:40
int GetPermutedIndex(int idx, int iq) const override
Returns the permuted index of the iq quadrature point in entity idx.
Definition: qspace.hpp:138
int dim
Definition: ex24.cpp:53
void Save(std::ostream &out) const override
Write the FaceQuadratureSpace to the stream out.
Definition: qspace.cpp:180
QuadratureSpace(Mesh *mesh_, int order_)
Create a QuadratureSpace based on the global rules from IntRules.
Definition: qspace.hpp:108
int Size() const
Return the logical size of the array.
Definition: array.hpp:141
void ConstructOffsets()
Definition: qspace.cpp:38
FaceType GetFaceType() const
Returns the face type (boundary or interior).
Definition: qspace.hpp:177
void ConstructIntRules(int dim)
Fill the int_rule array for each geometry type using order.
Definition: qspace.cpp:28
virtual int GetEntityIndex(const ElementTransformation &T) const =0
Returns the index in the quadrature space of the entity associated with the transformation T...
Class representing the storage layout of a QuadratureFunction.
Definition: qspace.hpp:101
int GetSize() const
Return the total number of quadrature points.
Definition: qspace.hpp:55
int GetEntityIndex(const ElementTransformation &T) const override
Returns the element index of T.
Definition: qspace.hpp:141
Represents values or vectors of values at quadrature points on a mesh.
Definition: qfunction.hpp:23
int order
The order of integration rule.
Definition: qspace.hpp:31
int GetOrder() const
Return the order of the quadrature rule(s) used by all elements.
Definition: qspace.hpp:58