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