MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
qfunction.hpp
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#ifndef MFEM_QFUNCTION
13#define MFEM_QFUNCTION
14
15#include "../config/config.hpp"
16#include "qspace.hpp"
17#include "gridfunc.hpp"
18
19namespace mfem
20{
21
22/// Represents values or vectors of values at quadrature points on a mesh.
24{
25protected:
26 QuadratureSpaceBase *qspace; ///< Associated QuadratureSpaceBase object.
27 bool own_qspace; ///< Does this own the associated QuadratureSpaceBase?
28 int vdim; ///< Vector dimension.
29
31
32public:
33 /// Default constructor, results in an empty vector.
34 QuadratureFunction() : qspace(nullptr), own_qspace(false), vdim(0)
35 { UseDevice(true); }
36
37 /// Create a QuadratureFunction based on the given QuadratureSpaceBase.
38 /** The QuadratureFunction does not assume ownership of the
39 QuadratureSpaceBase.
40 @note The Vector data is not initialized. */
41 QuadratureFunction(QuadratureSpaceBase &qspace_, int vdim_ = 1)
42 : Vector(vdim_*qspace_.GetSize()),
43 qspace(&qspace_), own_qspace(false), vdim(vdim_)
44 { UseDevice(true); }
45
46 /// Same as above but specify the device memory type
48 : Vector(vdim_*qspace_.GetSize(), mt),
49 qspace(&qspace_), own_qspace(false), vdim(vdim_)
50 { UseDevice(true); }
51
52 /// Create a QuadratureFunction based on the given QuadratureSpaceBase.
53 /** The QuadratureFunction does not assume ownership of the
54 QuadratureSpaceBase.
55 @warning @a qspace_ may not be NULL. */
56 QuadratureFunction(QuadratureSpaceBase *qspace_, int vdim_ = 1)
57 : QuadratureFunction(*qspace_, vdim_) { }
58
59 /// Same as above but specify the device memory type
61 : QuadratureFunction(*qspace_, mt, vdim_) { }
62
63 /** @brief Create a QuadratureFunction based on the given QuadratureSpaceBase,
64 using the external (host) data, @a qf_data. */
65 /** The QuadratureFunction does not assume ownership of the
66 QuadratureSpaceBase or the external data.
67 @warning @a qspace_ may not be NULL.
68 @note @a qf_data must be a valid **host** pointer (see the constructor
69 Vector::Vector(double *, int)). */
70 QuadratureFunction(QuadratureSpaceBase *qspace_, real_t *qf_data, int vdim_ = 1)
71 : Vector(qf_data, vdim_*qspace_->GetSize()),
72 qspace(qspace_), own_qspace(false), vdim(vdim_) { UseDevice(true); }
73
74 /** @brief Copy constructor. The QuadratureSpace ownership flag, #own_qspace,
75 in the new object is set to false. */
77 : QuadratureFunction(*orig.qspace, orig.vdim)
78 {
80 }
81
82 /// Read a QuadratureFunction from the stream @a in.
83 /** The QuadratureFunction assumes ownership of the read QuadratureSpace. */
84 QuadratureFunction(Mesh *mesh, std::istream &in);
85
86 /// Get the vector dimension.
87 int GetVDim() const { return vdim; }
88
89 /// Set the vector dimension, updating the size by calling Vector::SetSize().
90 void SetVDim(int vdim_)
91 { vdim = vdim_; SetSize(vdim*qspace->GetSize()); }
92
93 /// Get the associated QuadratureSpaceBase object.
95
96 /// Get the associated QuadratureSpaceBase object (const version).
97 const QuadratureSpaceBase *GetSpace() const { return qspace; }
98
99 /// Change the QuadratureSpaceBase and optionally the vector dimension.
100 /** If the new QuadratureSpaceBase is different from the current one, the
101 QuadratureFunction will not assume ownership of the new space; otherwise,
102 the ownership flag remains the same.
103
104 If the new vector dimension @a vdim_ < 0, the vector dimension remains
105 the same.
106
107 The data size is updated by calling Vector::SetSize(). */
108 inline void SetSpace(QuadratureSpaceBase *qspace_, int vdim_ = -1);
109
110 /** @brief Change the QuadratureSpaceBase, the data array, and optionally the
111 vector dimension. */
112 /** If the new QuadratureSpaceBase is different from the current one, the
113 QuadratureFunction will not assume ownership of the new space; otherwise,
114 the ownership flag remains the same.
115
116 If the new vector dimension @a vdim_ < 0, the vector dimension remains
117 the same.
118
119 The data array is replaced by calling Vector::NewDataAndSize(). */
120 inline void SetSpace(QuadratureSpaceBase *qspace_, real_t *qf_data,
121 int vdim_ = -1);
122
123 /// Get the QuadratureSpaceBase ownership flag.
124 bool OwnsSpace() { return own_qspace; }
125
126 /// Set the QuadratureSpaceBase ownership flag.
127 void SetOwnsSpace(bool own) { own_qspace = own; }
128
129 /// Set this equal to a constant value.
131
132 /// Copy the data from @a v.
133 /** The size of @a v must be equal to the size of the associated
134 QuadratureSpaceBase #qspace times the QuadratureFunction vector
135 dimension i.e. QuadratureFunction::Size(). */
137
138 /// Evaluate a grid function at each quadrature point.
140
141 /// Return all values associated with mesh element @a idx in a Vector.
142 /** The result is stored in the Vector @a values as a reference to the
143 global values.
144
145 Inside the Vector @a values, the index `i+vdim*j` corresponds to the
146 `i`-th vector component at the `j`-th quadrature point.
147 */
148 inline void GetValues(int idx, Vector &values);
149
150 /// Return all values associated with mesh element @a idx in a Vector.
151 /** The result is stored in the Vector @a values as a copy of the
152 global values.
153
154 Inside the Vector @a values, the index `i+vdim*j` corresponds to the
155 `i`-th vector component at the `j`-th quadrature point.
156 */
157 inline void GetValues(int idx, Vector &values) const;
158
159 /// Return the quadrature function values at an integration point.
160 /** The result is stored in the Vector @a values as a reference to the
161 global values. */
162 inline void GetValues(int idx, const int ip_num, Vector &values);
163
164 /// Return the quadrature function values at an integration point.
165 /** The result is stored in the Vector @a values as a copy to the
166 global values. */
167 inline void GetValues(int idx, const int ip_num, Vector &values) const;
168
169 /// Return all values associated with mesh element @a idx in a DenseMatrix.
170 /** The result is stored in the DenseMatrix @a values as a reference to the
171 global values.
172
173 Inside the DenseMatrix @a values, the `(i,j)` entry corresponds to the
174 `i`-th vector component at the `j`-th quadrature point.
175 */
176 inline void GetValues(int idx, DenseMatrix &values);
177
178 /// Return all values associated with mesh element @a idx in a const DenseMatrix.
179 /** The result is stored in the DenseMatrix @a values as a copy of the
180 global values.
181
182 Inside the DenseMatrix @a values, the `(i,j)` entry corresponds to the
183 `i`-th vector component at the `j`-th quadrature point.
184 */
185 inline void GetValues(int idx, DenseMatrix &values) const;
186
187 /// Get the IntegrationRule associated with entity (element or face) @a idx.
188 const IntegrationRule &GetIntRule(int idx) const
189 { return GetSpace()->GetIntRule(idx); }
190
191 /// Write the QuadratureFunction to the stream @a out.
192 void Save(std::ostream &out) const;
193
194 /// @brief Write the QuadratureFunction to @a out in VTU (ParaView) format.
195 ///
196 /// The data will be uncompressed if @a compression_level is zero, or if the
197 /// format is VTKFormat::ASCII. Otherwise, zlib compression will be used for
198 /// binary data.
199 void SaveVTU(std::ostream &out, VTKFormat format=VTKFormat::ASCII,
200 int compression_level=0, const std::string &field_name="u") const;
201
202 /// @brief Save the QuadratureFunction to a VTU (ParaView) file.
203 ///
204 /// The extension ".vtu" will be appended to @a filename.
205 /// @sa SaveVTU(std::ostream &out, VTKFormat format=VTKFormat::ASCII,
206 /// int compression_level=0)
207 void SaveVTU(const std::string &filename, VTKFormat format=VTKFormat::ASCII,
208 int compression_level=0, const std::string &field_name="u") const;
209
210
211 /// Return the integral of the quadrature function (vdim = 1 only).
213
214 /// @brief Integrate the (potentially vector-valued) quadrature function,
215 /// storing the results in @a integrals (length @a vdim).
216 void Integrate(Vector &integrals) const;
217
219 {
220 if (own_qspace) { delete qspace; }
221 }
222};
223
224// Inline methods
225
227 int idx, Vector &values)
228{
229 const int s_offset = qspace->Offset(idx);
230 const int sl_size = qspace->Offset(idx + 1) - s_offset;
231 values.MakeRef(*this, vdim*s_offset, vdim*sl_size);
232}
233
235 int idx, Vector &values) const
236{
237 const int s_offset = qspace->Offset(idx);
238 const int sl_size = qspace->Offset(idx + 1) - s_offset;
239 values.SetSize(vdim*sl_size);
240 values.HostWrite();
241 const real_t *q = HostRead() + vdim*s_offset;
242 for (int i = 0; i<values.Size(); i++)
243 {
244 values(i) = *(q++);
245 }
246}
247
249 int idx, const int ip_num, Vector &values)
250{
251 const int s_offset = qspace->Offset(idx) * vdim + ip_num * vdim;
252 values.MakeRef(*this, s_offset, vdim);
253}
254
256 int idx, const int ip_num, Vector &values) const
257{
258 const int s_offset = qspace->Offset(idx) * vdim + ip_num * vdim;
259 values.SetSize(vdim);
260 values.HostWrite();
261 const real_t *q = HostRead() + s_offset;
262 for (int i = 0; i < values.Size(); i++)
263 {
264 values(i) = *(q++);
265 }
266}
267
269 int idx, DenseMatrix &values)
270{
271 const int s_offset = qspace->Offset(idx);
272 const int sl_size = qspace->Offset(idx + 1) - s_offset;
273 // Make the values matrix memory an alias of the quadrature function memory
274 Memory<real_t> &values_mem = values.GetMemory();
275 values_mem.Delete();
276 values_mem.MakeAlias(GetMemory(), vdim*s_offset, vdim*sl_size);
277 values.SetSize(vdim, sl_size);
278}
279
281 int idx, DenseMatrix &values) const
282{
283 const int s_offset = qspace->Offset(idx);
284 const int sl_size = qspace->Offset(idx + 1) - s_offset;
285 values.SetSize(vdim, sl_size);
286 values.HostWrite();
287 const real_t *q = HostRead() + vdim*s_offset;
288 for (int j = 0; j<sl_size; j++)
289 {
290 for (int i = 0; i<vdim; i++)
291 {
292 values(i,j) = *(q++);
293 }
294 }
295}
296
297
299 int vdim_)
300{
301 if (qspace_ != qspace)
302 {
303 if (own_qspace) { delete qspace; }
304 qspace = qspace_;
305 own_qspace = false;
306 }
307 vdim = (vdim_ < 0) ? vdim : vdim_;
309}
310
312 QuadratureSpaceBase *qspace_, real_t *qf_data, int vdim_)
313{
314 if (qspace_ != qspace)
315 {
316 if (own_qspace) { delete qspace; }
317 qspace = qspace_;
318 own_qspace = false;
319 }
320 vdim = (vdim_ < 0) ? vdim : vdim_;
321 NewDataAndSize(qf_data, vdim*qspace->GetSize());
322}
323
324
325} // namespace mfem
326
327#endif
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
real_t * HostWrite()
Shortcut for mfem::Write(GetMemory(), TotalSize(), false).
Definition densemat.hpp:495
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition densemat.hpp:116
Memory< real_t > & GetMemory()
Definition densemat.hpp:128
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:32
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
Class used by MFEM to store pointers to host and/or device memory.
void MakeAlias(const Memory &base, int offset, int size)
Create a memory object that points inside the memory object base.
void Delete()
Delete the owned pointers and reset the Memory object.
Mesh data type.
Definition mesh.hpp:65
Represents values or vectors of values at quadrature points on a mesh.
Definition qfunction.hpp:24
void SetVDim(int vdim_)
Set the vector dimension, updating the size by calling Vector::SetSize().
Definition qfunction.hpp:90
QuadratureSpaceBase * GetSpace()
Get the associated QuadratureSpaceBase object.
Definition qfunction.hpp:94
void Integrate(Vector &integrals) const
Integrate the (potentially vector-valued) quadrature function, storing the results in integrals (leng...
QuadratureFunction(const QuadratureFunction &orig)
Copy constructor. The QuadratureSpace ownership flag, own_qspace, in the new object is set to false.
Definition qfunction.hpp:76
real_t Integrate() const
Return the integral of the quadrature function (vdim = 1 only).
void ProjectGridFunction(const GridFunction &gf)
Evaluate a grid function at each quadrature point.
QuadratureFunction()
Default constructor, results in an empty vector.
Definition qfunction.hpp:34
QuadratureFunction(QuadratureSpaceBase *qspace_, int vdim_=1)
Create a QuadratureFunction based on the given QuadratureSpaceBase.
Definition qfunction.hpp:56
void ProjectGridFunctionFallback(const GridFunction &gf)
void GetValues(int idx, Vector &values)
Return all values associated with mesh element idx in a Vector.
QuadratureFunction & operator=(real_t value)
Set this equal to a constant value.
void SaveVTU(const std::string &filename, VTKFormat format=VTKFormat::ASCII, int compression_level=0, const std::string &field_name="u") const
Save the QuadratureFunction to a VTU (ParaView) file.
QuadratureFunction(QuadratureSpaceBase *qspace_, MemoryType mt, int vdim_=1)
Same as above but specify the device memory type.
Definition qfunction.hpp:60
QuadratureFunction(QuadratureSpaceBase *qspace_, real_t *qf_data, int vdim_=1)
Create a QuadratureFunction based on the given QuadratureSpaceBase, using the external (host) data,...
Definition qfunction.hpp:70
const QuadratureSpaceBase * GetSpace() const
Get the associated QuadratureSpaceBase object (const version).
Definition qfunction.hpp:97
QuadratureFunction & operator=(const Vector &v)
Copy the data from v.
void SetOwnsSpace(bool own)
Set the QuadratureSpaceBase ownership flag.
void SaveVTU(std::ostream &out, VTKFormat format=VTKFormat::ASCII, int compression_level=0, const std::string &field_name="u") const
Write the QuadratureFunction to out in VTU (ParaView) format.
int GetVDim() const
Get the vector dimension.
Definition qfunction.hpp:87
QuadratureFunction(Mesh *mesh, std::istream &in)
Read a QuadratureFunction from the stream in.
QuadratureFunction(QuadratureSpaceBase &qspace_, int vdim_=1)
Create a QuadratureFunction based on the given QuadratureSpaceBase.
Definition qfunction.hpp:41
int vdim
Vector dimension.
Definition qfunction.hpp:28
bool OwnsSpace()
Get the QuadratureSpaceBase ownership flag.
QuadratureSpaceBase * qspace
Associated QuadratureSpaceBase object.
Definition qfunction.hpp:26
void SetSpace(QuadratureSpaceBase *qspace_, int vdim_=-1)
Change the QuadratureSpaceBase and optionally the vector dimension.
void Save(std::ostream &out) const
Write the QuadratureFunction to the stream out.
const IntegrationRule & GetIntRule(int idx) const
Get the IntegrationRule associated with entity (element or face) idx.
bool own_qspace
Does this own the associated QuadratureSpaceBase?
Definition qfunction.hpp:27
QuadratureFunction(QuadratureSpaceBase &qspace_, MemoryType mt, int vdim_=1)
Same as above but specify the device memory type.
Definition qfunction.hpp:47
Abstract base class for QuadratureSpace and FaceQuadratureSpace.
Definition qspace.hpp:32
int GetSize() const
Return the total number of quadrature points.
Definition qspace.hpp:107
const IntegrationRule & GetIntRule(int idx) const
Return the IntegrationRule associated with entity idx.
Definition qspace.hpp:125
int Offset(int idx) const
Gets the offset for a given entity idx.
Definition qspace.hpp:86
Vector data type.
Definition vector.hpp:82
virtual const real_t * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
Definition vector.hpp:524
Memory< real_t > & GetMemory()
Return a reference to the Memory object used by the Vector.
Definition vector.hpp:265
virtual bool UseDevice() const
Return the device flag of the Memory object used by the Vector.
Definition vector.hpp:148
int Size() const
Returns the size of the vector.
Definition vector.hpp:234
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:584
virtual real_t * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
Definition vector.hpp:532
void NewDataAndSize(real_t *d, int s)
Set the Vector data and size, deleting the old data, if owned.
Definition vector.hpp:197
Vector & operator=(const real_t *v)
Copy Size() entries from v.
Definition vector.cpp:197
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
Definition vector.hpp:660
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
VTKFormat
Data array format for VTK and VTU files.
Definition vtk.hpp:100
@ ASCII
Data arrays will be written in ASCII format.
float real_t
Definition config.hpp:46
MemoryType
Memory types supported by MFEM.