MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
qspace.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, 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#include "qspace.hpp"
13#include "qfunction.hpp"
14#include "../general/forall.hpp"
15
16namespace mfem
17{
18
20 const IntegrationRule &ir)
21 : mesh(mesh_), order(ir.GetOrder())
22{
23 for (int g = 0; g < Geometry::NumGeom; g++)
24 {
25 int_rule[g] = nullptr;
26 }
27 int_rule[geom] = &ir;
28}
29
31{
33 mesh.GetGeometries(dim, geoms);
34 for (Geometry::Type geom : geoms)
35 {
36 int_rule[geom] = &IntRules.Get(geom, order);
37 }
38}
39
41 QSpaceOffsetStorage storage) const
42{
43 if (storage == QSpaceOffsetStorage::COMPRESSED || offsets.Size() > 1)
44 {
45 return offsets;
46 }
47 else
48 {
49 if (full_offset_cache.Size() == 0)
50 {
51 const int nq = size / ne;
53 int *d_full_offset_cache = full_offset_cache.Write();
54 mfem::forall(ne + 1, [=] MFEM_HOST_DEVICE (int e)
55 {
56 d_full_offset_cache[e] = nq * e;
57 });
58 }
59 return full_offset_cache;
60 }
61}
62
63namespace
64{
65
66void ScaleByQuadratureWeights(Vector &weights, const IntegrationRule &ir)
67{
68 const int N = weights.Size();
69 const int n = ir.Size();
70 real_t *d_weights = weights.ReadWrite();
71 const real_t *d_w = ir.GetWeights().Read();
72
73 mfem::forall(N, [=] MFEM_HOST_DEVICE (int i)
74 {
75 d_weights[i] *= d_w[i%n];
76 });
77}
78
79} // anonymous namespace
80
82{
83 // First get the Jacobian determinants (without the quadrature weight
84 // contributions). We also store the pointer to the Vector object, so that
85 // we know when the cached weights are invalidated.
88
89 // Then scale by the quadrature weights.
90 const IntegrationRule &ir = GetIntRule(0);
91 ScaleByQuadratureWeights(weights, ir);
92}
93
95{
96 if (GetNE() == 0) { return weights; }
98 {
100 }
101 return weights;
102}
103
105{
106 QuadratureFunction qf(const_cast<QuadratureSpaceBase*>(this));
107 coeff.Project(qf);
108 return qf.Integrate();
109}
110
112 Vector &integrals) const
113{
114 const int vdim = coeff.GetVDim();
115 QuadratureFunction qf(const_cast<QuadratureSpaceBase*>(this), vdim);
116 coeff.Project(qf);
117 qf.Integrate(integrals);
118}
119
121{
122 const int num_elem = mesh.GetNE();
123 ne = num_elem;
124
126 {
129 offsets.SetSize(1);
131 offsets[0] = int_rule[geoms[0]]->GetNPoints();
132 size = num_elem * offsets[0];
133 }
134 else
135 {
136 offsets.SetSize(num_elem + 1);
137 int offset = 0;
138 for (int i = 0; i < num_elem; i++)
139 {
140 offsets[i] = offset;
142 MFEM_ASSERT(int_rule[geom] != nullptr, "Missing integration rule.");
143 offset += int_rule[geom]->GetNPoints();
144 }
145 offsets[num_elem] = offset;
146 size = offsets.Last();
147 }
148}
149
155
156QuadratureSpace::QuadratureSpace(Mesh *mesh_, std::istream &in)
157 : QuadratureSpaceBase(*mesh_)
158{
159 const char *msg = "invalid input stream";
160 std::string ident;
161
162 in >> ident; MFEM_VERIFY(ident == "QuadratureSpace", msg);
163 in >> ident; MFEM_VERIFY(ident == "Type:", msg);
164 in >> ident;
165 if (ident == "default_quadrature")
166 {
167 in >> ident; MFEM_VERIFY(ident == "Order:", msg);
168 in >> order;
169 }
170 else
171 {
172 MFEM_ABORT("unknown QuadratureSpace type: " << ident);
173 return;
174 }
175
176 Construct();
177}
178
180 : QuadratureSpaceBase(mesh_, mesh_.GetTypicalElementGeometry(), ir)
181{
182 MFEM_VERIFY(mesh.GetNumGeometries(mesh.Dimension()) <= 1,
183 "Constructor not valid for mixed meshes");
185}
186
187void QuadratureSpace::Save(std::ostream &os) const
188{
189 os << "QuadratureSpace\n"
190 << "Type: default_quadrature\n"
191 << "Order: " << order << '\n';
192}
193
195{
197 // TODO: assumes only one integration rule. This should be fixed once
198 // Mesh::GetGeometricFactors acceps a QuadratureSpace instead of
199 // IntegrationRule.
200 const IntegrationRule &ir = GetIntRule(0);
201 auto *geom = mesh.GetGeometricFactors(ir, flags);
202 return geom->detJ;
203}
204
206 FaceType face_type_)
207 : QuadratureSpaceBase(mesh_, order_), face_type(face_type_),
208 face_indices(mesh.GetFaceIndices(face_type_)),
209 face_indices_inv(mesh.GetInvFaceIndices(face_type_))
210{
211 Construct();
212}
213
215 FaceType face_type_)
216 : QuadratureSpaceBase(mesh_, mesh_.GetTypicalFaceGeometry(), ir),
217 face_type(face_type_),
218 face_indices(mesh.GetFaceIndices(face_type_)),
219 face_indices_inv(mesh.GetInvFaceIndices(face_type_))
220{
221 MFEM_VERIFY(mesh.GetNumGeometries(mesh.Dimension() - 1) <= 1,
222 "Constructor not valid for mixed meshes");
223 ConstructOffsets();
224}
225
226void FaceQuadratureSpace::ConstructOffsets()
227{
228 ne = face_indices.Size();
229
230 if (mesh.GetNumGeometries(mesh.Dimension() - 1) == 1)
231 {
233 mesh.GetGeometries(mesh.Dimension() - 1, geoms);
234 offsets.SetSize(1);
236 offsets[0] = int_rule[geoms[0]]->GetNPoints();
237 size = ne * offsets[0];
238 }
239 else
240 {
241 offsets.SetSize(face_indices.Size() + 1);
242 int offset = 0;
243 for (int i = 0; i < mesh.GetNFbyType(face_type); ++i)
244 {
245 offsets[i] = offset;
246 Geometry::Type geom = mesh.GetFaceGeometry(face_indices[i]);
247 MFEM_ASSERT(int_rule[geom] != nullptr, "Missing integration rule");
248 offset += int_rule[geom]->GetNPoints();
249 }
250 offsets[face_indices.Size()] = size = offset;
251 }
252}
253
254void FaceQuadratureSpace::Construct()
255{
257 ConstructOffsets();
258}
259
261{
262 const int f_idx = face_indices[idx];
264 {
265 const int dim = mesh.Dimension();
266 const IntegrationRule &ir = GetIntRule(idx);
267 const int q1d = (int)floor(pow(ir.GetNPoints(), 1.0/(dim-1)) + 0.5);
269 return ToLexOrdering(dim, face.element[0].local_face_id, q1d, iq);
270 }
271 else
272 {
273 return iq;
274 }
275}
276
278{
279 ElementTransformation *T = mesh.GetFaceTransformation(face_indices[idx]);
280 if (face_type == FaceType::Boundary)
281 {
283 }
284 return T;
285}
286
288{
289 auto get_face_index = [this](const int idx)
290 {
291 const auto it = face_indices_inv.find(idx);
292 if (it == face_indices_inv.end()) { return -1; }
293 else { return it->second; }
294 };
295
296 switch (T.ElementType)
297 {
299 return get_face_index(T.ElementNo);
302 return get_face_index(mesh.GetBdrElementFaceIndex(T.ElementNo));
303 default:
304 MFEM_ABORT("Invalid element type.");
305 return -1;
306 }
307}
308
309void FaceQuadratureSpace::Save(std::ostream &os) const
310{
311 os << "FaceQuadratureSpace\n"
312 << "Type: default_quadrature\n"
313 << "Order: " << order << '\n';
314}
315
316const Vector &FaceQuadratureSpace::GetGeometricFactorWeights() const
317{
319 // TODO: assumes only one integration rule. This should be fixed once
320 // Mesh::GetFaceGeometricFactors acceps a QuadratureSpace instead of
321 // IntegrationRule.
322 const IntegrationRule &ir = GetIntRule(0);
323 auto *geom = mesh.GetFaceGeometricFactors(ir, flags, face_type);
324 return geom->detJ;
325}
326
327} // namespace mfem
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:840
int Size() const
Return the logical size of the array.
Definition array.hpp:166
T * Write(bool on_dev=true)
Shortcut for mfem::Write(a.GetMemory(), a.Size(), on_dev).
Definition array.hpp:389
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition array.hpp:381
T * HostWrite()
Shortcut for mfem::Write(a.GetMemory(), a.Size(), false).
Definition array.hpp:393
T & Last()
Return the last element in the array.
Definition array.hpp:945
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
virtual void Project(QuadratureFunction &qf)
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
Vector detJ
Determinants of the Jacobians at all quadrature points.
Definition mesh.hpp:3169
int GetEntityIndex(const ElementTransformation &T) const override
Returns the index associated with the face described by T.
Definition qspace.cpp:287
Geometry::Type GetGeometry(int idx) const override
Returns the geometry type of face idx.
Definition qspace.hpp:246
void Save(std::ostream &out) const override
Write the FaceQuadratureSpace to the stream out.
Definition qspace.cpp:309
int GetPermutedIndex(int idx, int iq) const override
Returns the permuted index of the iq quadrature point in entity idx.
Definition qspace.cpp:260
FaceQuadratureSpace(Mesh &mesh_, int order_, FaceType face_type_)
Create a FaceQuadratureSpace based on the global rules from IntRules.
Definition qspace.cpp:205
ElementTransformation * GetTransformation(int idx) override
Returns the face transformation of face idx.
Definition qspace.cpp:277
Vector detJ
Determinants of the Jacobians at all quadrature points.
Definition mesh.hpp:3121
static const int NumGeom
Definition geom.hpp:46
static bool IsTensorProduct(Type geom)
Definition geom.hpp:112
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition intrules.hpp:256
const Array< real_t > & GetWeights() const
Return the quadrature weights in a contiguous array.
Definition intrules.cpp:86
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Mesh data type.
Definition mesh.hpp:65
void GetGeometries(int dim, Array< Geometry::Type > &el_geoms) const
Return all element geometries of the given dimension present in the mesh.
Definition mesh.cpp:7569
Geometry::Type GetFaceGeometry(int i) const
Return the Geometry::Type associated with face i.
Definition mesh.cpp:1576
int GetBdrElementFaceIndex(int be_idx) const
Return the local face (codimension-1) index for the given boundary element index.
Definition mesh.hpp:1689
virtual int GetNFbyType(FaceType type) const
Returns the number of faces according to the requested type, does not count master nonconforming face...
Definition mesh.cpp:6862
ElementTransformation * GetFaceTransformation(int FaceNo)
Returns a pointer to the transformation defining the given face element.
Definition mesh.cpp:610
const FaceGeometricFactors * GetFaceGeometricFactors(const IntegrationRule &ir, const int flags, FaceType type, MemoryType d_mt=MemoryType::DEFAULT)
Return the mesh geometric factors for the faces corresponding to the given integration rule.
Definition mesh.cpp:903
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1377
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1306
FaceInformation GetFaceInformation(int f) const
Definition mesh.cpp:1293
const GeometricFactors * GetGeometricFactors(const IntegrationRule &ir, const int flags, MemoryType d_mt=MemoryType::DEFAULT)
Return the mesh geometric factors corresponding to the given integration rule.
Definition mesh.cpp:883
const Array< int > & GetBdrFaceAttributes() const
Returns the attributes for all boundary elements in this mesh.
Definition mesh.cpp:925
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
Definition mesh.cpp:7558
Geometry::Type GetElementBaseGeometry(int i) const
Definition mesh.hpp:1556
long GetNodesSequence() const
Return the nodes update counter.
Definition mesh.hpp:2526
Represents values or vectors of values at quadrature points on a mesh.
Definition qfunction.hpp:24
real_t Integrate() const
Return the integral of the quadrature function (vdim = 1 only).
Abstract base class for QuadratureSpace and FaceQuadratureSpace.
Definition qspace.hpp:32
Vector weights
Integration weights.
Definition qspace.hpp:40
const Vector & GetWeights() const
Return the integration weights (including geometric factors).
Definition qspace.cpp:94
void ConstructWeights() const
Compute the integration weights.
Definition qspace.cpp:81
int GetNE() const
Return the number of entities.
Definition qspace.hpp:113
void ConstructIntRules(int dim)
Fill the int_rule array for each geometry type using order.
Definition qspace.cpp:30
const IntegrationRule * int_rule[Geometry::NumGeom]
The quadrature rules used for each geometry type.
Definition qspace.hpp:62
QuadratureSpaceBase(Mesh &mesh_, int order_=0)
Protected constructor. Used by derived classes.
Definition qspace.hpp:65
long nodes_sequence
Nodes counter for cache invalidation.
Definition qspace.hpp:41
int order
The order of integration rule.
Definition qspace.hpp:37
const IntegrationRule & GetIntRule(int idx) const
Return the IntegrationRule associated with entity idx.
Definition qspace.hpp:125
virtual const Vector & GetGeometricFactorWeights() const =0
Compute the det(J) (volume or faces, depending on the type).
int ne
Number of entities.
Definition qspace.hpp:39
const Array< int > & Offsets(QSpaceOffsetStorage storage) const
Entity quadrature point offset array.
Definition qspace.cpp:40
Array< int > full_offset_cache
Cached version of the "full" offsets, returned by Offsets() when QSpaceOffsetStorage::FULL is provide...
Definition qspace.hpp:59
Mesh & mesh
The underlying mesh.
Definition qspace.hpp:36
Array< int > offsets
Entity quadrature point offset array.
Definition qspace.hpp:52
int size
Total number of quadrature points.
Definition qspace.hpp:38
real_t Integrate(Coefficient &coeff) const
Return the integral of the scalar Coefficient coeff.
Definition qspace.cpp:104
const Vector & GetGeometricFactorWeights() const override
Compute the det(J) (volume or faces, depending on the type).
Definition qspace.cpp:194
void Save(std::ostream &out) const override
Write the QuadratureSpace to the stream out.
Definition qspace.cpp:187
QuadratureSpace(Mesh *mesh_, int order_)
Create a QuadratureSpace based on the global rules from IntRules.
Definition qspace.hpp:171
Base class for vector Coefficients that optionally depend on time and space.
int GetVDim()
Returns dimension of the vector.
virtual void Project(QuadratureFunction &qf)
Fill the QuadratureFunction qf by evaluating the coefficient at the quadrature points.
Vector data type.
Definition vector.hpp:82
virtual real_t * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:536
int Size() const
Returns the size of the vector.
Definition vector.hpp:234
int dim
Definition ex24.cpp:53
int ToLexOrdering(const int dim, const int face_id, const int size1d, const int index)
Convert a dof face index from Native ordering to lexicographic ordering for quads and hexes.
QSpaceOffsetStorage
Definition qspace.hpp:23
float real_t
Definition config.hpp:46
void forall(int N, lambda &&body)
Definition forall.hpp:839
FaceType
Definition mesh.hpp:49
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition intrules.hpp:492
This structure is used as a human readable output format that deciphers the information contained in ...
Definition mesh.hpp:2081
struct mfem::Mesh::FaceInformation::@15 element[2]