MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
qspace.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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
19namespace 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{
27protected:
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 mutable Vector weights; ///< Integration weights.
34 mutable long nodes_sequence = 0; ///< Nodes counter for cache invalidation.
35
36 /// @brief Entity quadrature point offset array, of size num_entities + 1.
37 ///
38 /// The quadrature point values for entity i are stored in the indices between
39 /// offsets[i] and offsets[i+1].
41 /// The quadrature rules used for each geometry type.
43
44 /// Protected constructor. Used by derived classes.
45 QuadratureSpaceBase(Mesh &mesh_, int order_ = 0)
46 : mesh(mesh_), order(order_) { }
47
48 /// Protected constructor. Used by derived classes.
50 const IntegrationRule &ir);
51
52 /// Fill the @ref int_rule array for each geometry type using @ref order.
53 void ConstructIntRules(int dim);
54
55 /// Compute the det(J) (volume or faces, depending on the type).
56 virtual const Vector &GetGeometricFactorWeights() const = 0;
57
58 /// Compute the integration weights.
59 void ConstructWeights() const;
60
61public:
62 /// Return the total number of quadrature points.
63 int GetSize() const { return size; }
64
65 /// Return the order of the quadrature rule(s) used by all elements.
66 int GetOrder() const { return order; }
67
68 /// Return the number of entities.
69 int GetNE() const { return offsets.Size() - 1; }
70
71 /// Returns the mesh.
72 inline Mesh *GetMesh() const { return &mesh; }
73
74 /// Get the (element or face) transformation of entity @a idx.
76
77 /// Return the geometry type of entity (element or face) @a idx.
78 virtual Geometry::Type GetGeometry(int idx) const = 0;
79
80 /// Return the IntegrationRule associated with entity @a idx.
81 const IntegrationRule &GetIntRule(int idx) const
82 { return *int_rule[GetGeometry(idx)]; }
83
84 /// @brief Returns the permuted index of the @a iq quadrature point in entity
85 /// @a idx.
86 ///
87 /// For tensor-product faces, returns the lexicographic index of the
88 /// quadrature point, oriented relative to "element 1". For QuadratureSpace%s
89 /// defined on elements (not faces), the permutation is trivial, and this
90 /// returns @a iq.
91 virtual int GetPermutedIndex(int idx, int iq) const = 0;
92
93 /// @brief Returns the index in the quadrature space of the entity associated
94 /// with the transformation @a T.
95 ///
96 /// For a QuadratureSpace defined on elements, this just returns the element
97 /// index. For FaceQuadratureSpace, the returned index depends on the chosen
98 /// FaceType. If the entity is not found (for example, if @a T represents an
99 /// interior face, and the space has FaceType::Boundary) then -1 is returned.
100 virtual int GetEntityIndex(const ElementTransformation &T) const = 0;
101
102 /// Write the QuadratureSpace to the stream @a out.
103 virtual void Save(std::ostream &out) const = 0;
104
105 /// Return the integration weights (including geometric factors).
106 const Vector &GetWeights() const;
107
108 /// Return the integral of the scalar Coefficient @a coeff.
109 real_t Integrate(Coefficient &coeff) const;
110
111 /// Return the integral of the VectorCoefficient @a coeff in @a integrals.
112 void Integrate(VectorCoefficient &coeff, Vector &integrals) const;
113
115};
116
117/// Class representing the storage layout of a QuadratureFunction.
118/** Multiple QuadratureFunction%s can share the same QuadratureSpace. */
120{
121protected:
122 const Vector &GetGeometricFactorWeights() const override;
123 void ConstructOffsets();
124 void Construct();
125public:
126 /// Create a QuadratureSpace based on the global rules from #IntRules.
127 QuadratureSpace(Mesh *mesh_, int order_)
128 : QuadratureSpaceBase(*mesh_, order_) { Construct(); }
129
130 /// @brief Create a QuadratureSpace with an IntegrationRule, valid only when
131 /// the mesh has one element type.
132 QuadratureSpace(Mesh &mesh_, const IntegrationRule &ir);
133
134 /// Read a QuadratureSpace from the stream @a in.
135 QuadratureSpace(Mesh *mesh_, std::istream &in);
136
137 /// Returns number of elements in the mesh.
138 inline int GetNE() const { return mesh.GetNE(); }
139
140 /// Returns the element transformation of element @a idx.
143
144 /// Returns the geometry type of element @a idx.
145 Geometry::Type GetGeometry(int idx) const override
146 { return mesh.GetElementGeometry(idx); }
147
148 /// Get the IntegrationRule associated with mesh element @a idx.
149 const IntegrationRule &GetElementIntRule(int idx) const
150 { return *int_rule[mesh.GetElementBaseGeometry(idx)]; }
151
152 /// @brief Returns the permuted index of the @a iq quadrature point in entity
153 /// @a idx.
154 ///
155 /// The member function QuadratureSpace::GetPermutedIndex always returns @a
156 /// iq, the permutation is only nontrivial for FaceQuadratureSpace.
157 int GetPermutedIndex(int idx, int iq) const override { return iq; }
158
159 /// Returns the element index of @a T.
160 int GetEntityIndex(const ElementTransformation &T) const override { return T.ElementNo; }
161
162 /// Write the QuadratureSpace to the stream @a out.
163 void Save(std::ostream &out) const override;
164};
165
166/// Class representing the storage layout of a FaceQuadratureFunction.
167/** FaceQuadratureSpace is defined on either the interior or boundary faces
168 of a mesh, depending on the provided FaceType. */
170{
171 FaceType face_type; ///< Is the space defined on interior or boundary faces?
172 const int num_faces; ///< Number of faces.
173
174 /// Map from boundary or interior face indices to mesh face indices.
175 Array<int> face_indices;
176
177 /// Inverse of the map @a face_indices.
178 std::unordered_map<int,int> face_indices_inv;
179
180 const Vector &GetGeometricFactorWeights() const override;
181 void ConstructOffsets();
182 void Construct();
183
184public:
185 /// Create a FaceQuadratureSpace based on the global rules from #IntRules.
186 FaceQuadratureSpace(Mesh &mesh_, int order_, FaceType face_type_);
187
188 /// @brief Create a FaceQuadratureSpace with an IntegrationRule, valid only
189 /// when the mesh has one type of face geometry.
190 FaceQuadratureSpace(Mesh &mesh_, const IntegrationRule &ir,
191 FaceType face_type_);
192
193 /// Returns number of faces in the mesh.
194 inline int GetNumFaces() const { return num_faces; }
195
196 /// Returns the face type (boundary or interior).
197 FaceType GetFaceType() const { return face_type; }
198
199 /// Returns the face transformation of face @a idx.
201 { return mesh.GetFaceTransformation(face_indices[idx]); }
202
203 /// Returns the geometry type of face @a idx.
204 Geometry::Type GetGeometry(int idx) const override
205 { return mesh.GetFaceGeometry(face_indices[idx]); }
206
207 /// Get the IntegrationRule associated with mesh element @a idx.
208 const IntegrationRule &GetFaceIntRule(int idx) const
209 { return *int_rule[GetGeometry(idx)]; }
210
211 /// @brief Returns the permuted index of the @a iq quadrature point in entity
212 /// @a idx.
213 ///
214 /// For tensor-product faces, returns the lexicographic index of the
215 /// quadrature point, oriented relative to "element 1".
216 int GetPermutedIndex(int idx, int iq) const override;
217
218 /// @brief Get the face index (in the standard Mesh numbering) associated
219 /// with face @a idx in the FaceQuadratureSpace.
220 int GetMeshFaceIndex(int idx) const { return face_indices[idx]; }
221
222 /// @brief Returns the index associated with the face described by @a T.
223 ///
224 /// The index may differ from the mesh face or boundary element index
225 /// depending on the FaceType used to construct the FaceQuadratureSpace.
226 int GetEntityIndex(const ElementTransformation &T) const override;
227
228 /// Write the FaceQuadratureSpace to the stream @a out.
229 void Save(std::ostream &out) const override;
230};
231
232}
233
234#endif
int Size() const
Return the logical size of the array.
Definition array.hpp:144
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Class representing the storage layout of a FaceQuadratureFunction.
Definition qspace.hpp:170
FaceType GetFaceType() const
Returns the face type (boundary or interior).
Definition qspace.hpp:197
int GetMeshFaceIndex(int idx) const
Get the face index (in the standard Mesh numbering) associated with face idx in the FaceQuadratureSpa...
Definition qspace.hpp:220
int GetEntityIndex(const ElementTransformation &T) const override
Returns the index associated with the face described by T.
Definition qspace.cpp:235
Geometry::Type GetGeometry(int idx) const override
Returns the geometry type of face idx.
Definition qspace.hpp:204
void Save(std::ostream &out) const override
Write the FaceQuadratureSpace to the stream out.
Definition qspace.cpp:257
int GetPermutedIndex(int idx, int iq) const override
Returns the permuted index of the iq quadrature point in entity idx.
Definition qspace.cpp:218
const IntegrationRule & GetFaceIntRule(int idx) const
Get the IntegrationRule associated with mesh element idx.
Definition qspace.hpp:208
FaceQuadratureSpace(Mesh &mesh_, int order_, FaceType face_type_)
Create a FaceQuadratureSpace based on the global rules from IntRules.
Definition qspace.cpp:167
int GetNumFaces() const
Returns number of faces in the mesh.
Definition qspace.hpp:194
ElementTransformation * GetTransformation(int idx) override
Returns the face transformation of face idx.
Definition qspace.hpp:200
static const int NumGeom
Definition geom.hpp:42
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
Mesh data type.
Definition mesh.hpp:56
Geometry::Type GetFaceGeometry(int i) const
Return the Geometry::Type associated with face i.
Definition mesh.cpp:1452
Geometry::Type GetElementGeometry(int i) const
Definition mesh.hpp:1371
ElementTransformation * GetFaceTransformation(int FaceNo)
Returns a pointer to the transformation defining the given face element.
Definition mesh.cpp:583
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1226
void GetElementTransformation(int i, IsoparametricTransformation *ElTr) const
Builds the transformation defining the i-th element in ElTr. ElTr must be allocated in advance and wi...
Definition mesh.cpp:357
Geometry::Type GetElementBaseGeometry(int i) const
Definition mesh.hpp:1385
Represents values or vectors of values at quadrature points on a mesh.
Definition qfunction.hpp:24
Abstract base class for QuadratureSpace and FaceQuadratureSpace.
Definition qspace.hpp:26
Vector weights
Integration weights.
Definition qspace.hpp:33
virtual Geometry::Type GetGeometry(int idx) const =0
Return the geometry type of entity (element or face) idx.
const Vector & GetWeights() const
Return the integration weights (including geometric factors).
Definition qspace.cpp:71
void ConstructWeights() const
Compute the integration weights.
Definition qspace.cpp:58
int GetNE() const
Return the number of entities.
Definition qspace.hpp:69
virtual void Save(std::ostream &out) const =0
Write the QuadratureSpace to the stream out.
void ConstructIntRules(int dim)
Fill the int_rule array for each geometry type using order.
Definition qspace.cpp:30
virtual ElementTransformation * GetTransformation(int idx)=0
Get the (element or face) transformation of entity idx.
int GetSize() const
Return the total number of quadrature points.
Definition qspace.hpp:63
const IntegrationRule * int_rule[Geometry::NumGeom]
The quadrature rules used for each geometry type.
Definition qspace.hpp:42
QuadratureSpaceBase(Mesh &mesh_, int order_=0)
Protected constructor. Used by derived classes.
Definition qspace.hpp:45
long nodes_sequence
Nodes counter for cache invalidation.
Definition qspace.hpp:34
int order
The order of integration rule.
Definition qspace.hpp:31
const IntegrationRule & GetIntRule(int idx) const
Return the IntegrationRule associated with entity idx.
Definition qspace.hpp:81
virtual const Vector & GetGeometricFactorWeights() const =0
Compute the det(J) (volume or faces, depending on the type).
virtual int GetEntityIndex(const ElementTransformation &T) const =0
Returns the index in the quadrature space of the entity associated with the transformation T.
virtual int GetPermutedIndex(int idx, int iq) const =0
Returns the permuted index of the iq quadrature point in entity idx.
int GetOrder() const
Return the order of the quadrature rule(s) used by all elements.
Definition qspace.hpp:66
Mesh & mesh
The underlying mesh.
Definition qspace.hpp:30
Array< int > offsets
Entity quadrature point offset array, of size num_entities + 1.
Definition qspace.hpp:40
int size
Total number of quadrature points.
Definition qspace.hpp:32
virtual ~QuadratureSpaceBase()
Definition qspace.hpp:114
Mesh * GetMesh() const
Returns the mesh.
Definition qspace.hpp:72
real_t Integrate(Coefficient &coeff) const
Return the integral of the scalar Coefficient coeff.
Definition qspace.cpp:81
Class representing the storage layout of a QuadratureFunction.
Definition qspace.hpp:120
int GetEntityIndex(const ElementTransformation &T) const override
Returns the element index of T.
Definition qspace.hpp:160
int GetPermutedIndex(int idx, int iq) const override
Returns the permuted index of the iq quadrature point in entity idx.
Definition qspace.hpp:157
const Vector & GetGeometricFactorWeights() const override
Compute the det(J) (volume or faces, depending on the type).
Definition qspace.cpp:156
int GetNE() const
Returns number of elements in the mesh.
Definition qspace.hpp:138
const IntegrationRule & GetElementIntRule(int idx) const
Get the IntegrationRule associated with mesh element idx.
Definition qspace.hpp:149
ElementTransformation * GetTransformation(int idx) override
Returns the element transformation of element idx.
Definition qspace.hpp:141
void Save(std::ostream &out) const override
Write the QuadratureSpace to the stream out.
Definition qspace.cpp:149
QuadratureSpace(Mesh *mesh_, int order_)
Create a QuadratureSpace based on the global rules from IntRules.
Definition qspace.hpp:127
Geometry::Type GetGeometry(int idx) const override
Returns the geometry type of element idx.
Definition qspace.hpp:145
Base class for vector Coefficients that optionally depend on time and space.
Vector data type.
Definition vector.hpp:80
int dim
Definition ex24.cpp:53
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
float real_t
Definition config.hpp:43
FaceType
Definition mesh.hpp:47