MFEM  v4.6.0
Finite element discretization library
qfunction.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_QFUNCTION
13 #define MFEM_QFUNCTION
14 
15 #include "../config/config.hpp"
16 #include "qspace.hpp"
17 #include "gridfunc.hpp"
18 
19 namespace mfem
20 {
21 
22 /// Represents values or vectors of values at quadrature points on a mesh.
23 class QuadratureFunction : public Vector
24 {
25 protected:
26  QuadratureSpaceBase *qspace; ///< Associated QuadratureSpaceBase object.
27  bool own_qspace; ///< Does this own the associated QuadratureSpaceBase?
28  int vdim; ///< Vector dimension.
29 
30 public:
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_, double *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  {
67  Vector::operator=(orig);
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_, double *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.
118  QuadratureFunction &operator=(double 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.
127  void ProjectGridFunction(const GridFunction &gf);
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 
199  {
200  if (own_qspace) { delete qspace; }
201  }
202 };
203 
204 // Inline methods
205 
207  int idx, Vector &values)
208 {
209  const int s_offset = qspace->offsets[idx];
210  const int sl_size = qspace->offsets[idx+1] - s_offset;
211  values.MakeRef(*this, vdim*s_offset, vdim*sl_size);
212 }
213 
215  int idx, Vector &values) const
216 {
217  const int s_offset = qspace->offsets[idx];
218  const int sl_size = qspace->offsets[idx+1] - s_offset;
219  values.SetSize(vdim*sl_size);
220  values.HostWrite();
221  const double *q = HostRead() + vdim*s_offset;
222  for (int i = 0; i<values.Size(); i++)
223  {
224  values(i) = *(q++);
225  }
226 }
227 
229  int idx, const int ip_num, Vector &values)
230 {
231  const int s_offset = qspace->offsets[idx] * vdim + ip_num * vdim;
232  values.MakeRef(*this, s_offset, vdim);
233 }
234 
236  int idx, const int ip_num, Vector &values) const
237 {
238  const int s_offset = qspace->offsets[idx] * vdim + ip_num * vdim;
239  values.SetSize(vdim);
240  values.HostWrite();
241  const double *q = HostRead() + s_offset;
242  for (int i = 0; i < values.Size(); i++)
243  {
244  values(i) = *(q++);
245  }
246 }
247 
249  int idx, DenseMatrix &values)
250 {
251  const int s_offset = qspace->offsets[idx];
252  const int sl_size = qspace->offsets[idx+1] - s_offset;
253  // Make the values matrix memory an alias of the quadrature function memory
254  Memory<double> &values_mem = values.GetMemory();
255  values_mem.Delete();
256  values_mem.MakeAlias(GetMemory(), vdim*s_offset, vdim*sl_size);
257  values.SetSize(vdim, sl_size);
258 }
259 
261  int idx, DenseMatrix &values) const
262 {
263  const int s_offset = qspace->offsets[idx];
264  const int sl_size = qspace->offsets[idx+1] - s_offset;
265  values.SetSize(vdim, sl_size);
266  values.HostWrite();
267  const double *q = HostRead() + vdim*s_offset;
268  for (int j = 0; j<sl_size; j++)
269  {
270  for (int i = 0; i<vdim; i++)
271  {
272  values(i,j) = *(q++);
273  }
274  }
275 }
276 
277 
279  int vdim_)
280 {
281  if (qspace_ != qspace)
282  {
283  if (own_qspace) { delete qspace; }
284  qspace = qspace_;
285  own_qspace = false;
286  }
287  vdim = (vdim_ < 0) ? vdim : vdim_;
289 }
290 
292  QuadratureSpaceBase *qspace_, double *qf_data, int vdim_)
293 {
294  if (qspace_ != qspace)
295  {
296  if (own_qspace) { delete qspace; }
297  qspace = qspace_;
298  own_qspace = false;
299  }
300  vdim = (vdim_ < 0) ? vdim : vdim_;
301  NewDataAndSize(qf_data, vdim*qspace->GetSize());
302 }
303 
304 
305 } // namespace mfem
306 
307 #endif
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:96
void NewDataAndSize(double *d, int s)
Set the Vector data and size, deleting the old data, if owned.
Definition: vector.hpp:160
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
void ProjectGridFunction(const GridFunction &gf)
Evaluate a grid function at each quadrature point.
Definition: qfunction.cpp:57
double * HostWrite()
Shortcut for mfem::Write(GetMemory(), TotalSize(), false).
Definition: densemat.hpp:478
Memory< double > & GetMemory()
Definition: densemat.hpp:117
QuadratureFunction(QuadratureSpaceBase *qspace_, double *qf_data, int vdim_=1)
Create a QuadratureFunction based on the given QuadratureSpaceBase, using the external (host) data...
Definition: qfunction.hpp:58
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:517
void Delete()
Delete the owned pointers and reset the Memory object.
QuadratureSpaceBase * qspace
Associated QuadratureSpaceBase object.
Definition: qfunction.hpp:26
virtual const double * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:457
bool own_qspace
Does this own the associated QuadratureSpaceBase?
Definition: qfunction.hpp:27
int Size() const
Returns the size of the vector.
Definition: vector.hpp:197
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
virtual double * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:465
int vdim
Vector dimension.
Definition: qfunction.hpp:28
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.
Definition: qfunction.cpp:117
Data arrays will be written in ASCII format.
bool OwnsSpace()
Get the QuadratureSpaceBase ownership flag.
Definition: qfunction.hpp:112
Memory< double > & GetMemory()
Return a reference to the Memory object used by the Vector.
Definition: vector.hpp:228
Vector & operator=(const double *v)
Copy Size() entries from v.
Definition: vector.cpp:119
void SetOwnsSpace(bool own)
Set the QuadratureSpaceBase ownership flag.
Definition: qfunction.hpp:115
QuadratureFunction(QuadratureSpaceBase &qspace_, int vdim_=1)
Create a QuadratureFunction based on the given QuadratureSpaceBase.
Definition: qfunction.hpp:39
int GetVDim() const
Get the vector dimension.
Definition: qfunction.hpp:75
const IntegrationRule & GetIntRule(int idx) const
Get the IntegrationRule associated with entity (element or face) idx.
Definition: qfunction.hpp:176
QuadratureFunction()
Default constructor, results in an empty vector.
Definition: qfunction.hpp:32
VTKFormat
Data array format for VTK and VTU files.
Definition: vtk.hpp:98
void SetVDim(int vdim_)
Set the vector dimension, updating the size by calling Vector::SetSize().
Definition: qfunction.hpp:78
QuadratureFunction & operator=(double value)
Set this equal to a constant value.
Definition: qfunction.cpp:19
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
void Save(std::ostream &out) const
Write the QuadratureFunction to the stream out.
Definition: qfunction.cpp:48
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
void GetValues(int idx, Vector &values)
Return all values associated with mesh element idx in a Vector.
Definition: qfunction.hpp:206
void MakeAlias(const Memory &base, int offset, int size)
Create a memory object that points inside the memory object base.
QuadratureSpaceBase * GetSpace()
Get the associated QuadratureSpaceBase object.
Definition: qfunction.hpp:82
void SetSpace(QuadratureSpaceBase *qspace_, int vdim_=-1)
Change the QuadratureSpaceBase and optionally the vector dimension.
Definition: qfunction.hpp:278
QuadratureFunction(const QuadratureFunction &orig)
Copy constructor. The QuadratureSpace ownership flag, own_qspace, in the new object is set to false...
Definition: qfunction.hpp:64
virtual bool UseDevice() const
Return the device flag of the Memory object used by the Vector.
Definition: vector.hpp:118
Vector data type.
Definition: vector.hpp:58
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
Definition: vector.hpp:581
QuadratureFunction(QuadratureSpaceBase *qspace_, int vdim_=1)
Create a QuadratureFunction based on the given QuadratureSpaceBase.
Definition: qfunction.hpp:48
int GetSize() const
Return the total number of quadrature points.
Definition: qspace.hpp:55
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:105
const QuadratureSpaceBase * GetSpace() const
Get the associated QuadratureSpaceBase object (const version).
Definition: qfunction.hpp:85
Represents values or vectors of values at quadrature points on a mesh.
Definition: qfunction.hpp:23