MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
qfunction.cpp
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#include "qfunction.hpp"
13#include "quadinterpolator.hpp"
15#include "../general/forall.hpp"
16#include "../mesh/pmesh.hpp"
17
18namespace mfem
19{
20
21QuadratureFunction &QuadratureFunction::operator=(real_t value)
22{
23 Vector::operator=(value);
24 return *this;
25}
26
27QuadratureFunction &QuadratureFunction::operator=(const Vector &v)
28{
29 MFEM_ASSERT(qspace && v.Size() == this->Size(), "");
31 return *this;
32}
33
34
35QuadratureFunction::QuadratureFunction(Mesh *mesh, std::istream &in)
36 : QuadratureFunction()
37{
38 const char *msg = "invalid input stream";
39 std::string ident;
40
41 qspace = new QuadratureSpace(mesh, in);
42 own_qspace = true;
43
44 in >> ident; MFEM_VERIFY(ident == "VDim:", msg);
45 in >> vdim;
46
47 Load(in, vdim*qspace->GetSize());
48}
49
50void QuadratureFunction::Save(std::ostream &os) const
51{
52 GetSpace()->Save(os);
53 os << "VDim: " << vdim << '\n'
54 << '\n';
55 Vector::Print(os, vdim);
56 os.flush();
57}
58
59void QuadratureFunction::ProjectGridFunction(const GridFunction &gf)
60{
61 SetVDim(gf.VectorDim());
62
63 if (auto *qs_elem = dynamic_cast<QuadratureSpace*>(qspace))
64 {
65 const FiniteElementSpace &gf_fes = *gf.FESpace();
66 const bool use_tensor_products = UsesTensorBasis(gf_fes);
67 const ElementDofOrdering ordering = use_tensor_products ?
68 ElementDofOrdering::LEXICOGRAPHIC :
69 ElementDofOrdering::NATIVE;
70
71 // Use element restriction to go from L-vector to E-vector
72 const Operator *R = gf_fes.GetElementRestriction(ordering);
73 Vector e_vec(R->Height());
74 R->Mult(gf, e_vec);
75
76 // Use quadrature interpolator to go from E-vector to Q-vector
77 const QuadratureInterpolator *qi =
78 gf_fes.GetQuadratureInterpolator(*qs_elem);
79 qi->SetOutputLayout(QVectorLayout::byVDIM);
80 qi->DisableTensorProducts(!use_tensor_products);
81 qi->Values(e_vec, *this);
82 }
83 else if (auto *qs_face = dynamic_cast<FaceQuadratureSpace*>(qspace))
84 {
85 const FiniteElementSpace &gf_fes = *gf.FESpace();
86 const bool use_tensor_products = UsesTensorBasis(gf_fes);
87 const ElementDofOrdering ordering = use_tensor_products ?
88 ElementDofOrdering::LEXICOGRAPHIC :
89 ElementDofOrdering::NATIVE;
90
91 const FaceType face_type = qs_face->GetFaceType();
92
93 // Use element restriction to go from L-vector to E-vector
94 const Operator *R = gf_fes.GetFaceRestriction(
95 ordering, face_type, L2FaceValues::SingleValued);
96 Vector e_vec(R->Height());
97 R->Mult(gf, e_vec);
98
99 // Use quadrature interpolator to go from E-vector to Q-vector
100 const FaceQuadratureInterpolator *qi =
101 gf_fes.GetFaceQuadratureInterpolator(qspace->GetIntRule(0), face_type);
102 qi->SetOutputLayout(QVectorLayout::byVDIM);
103 qi->DisableTensorProducts(!use_tensor_products);
104 qi->Values(e_vec, *this);
105 }
106 else
107 {
108 // This branch should be unreachable
109 MFEM_ABORT("Unsupported case.");
110 }
111}
112
113std::ostream &operator<<(std::ostream &os, const QuadratureFunction &qf)
114{
115 qf.Save(os);
116 return os;
117}
118
119void QuadratureFunction::SaveVTU(std::ostream &os, VTKFormat format,
120 int compression_level,
121 const std::string &field_name) const
122{
123 os << R"(<VTKFile type="UnstructuredGrid" version="0.1")";
124 if (compression_level != 0)
125 {
126 os << R"( compressor="vtkZLibDataCompressor")";
127 }
128 os << " byte_order=\"" << VTKByteOrder() << "\">\n";
129 os << "<UnstructuredGrid>\n";
130
131 const char *fmt_str = (format == VTKFormat::ASCII) ? "ascii" : "binary";
132 const char *type_str = (format != VTKFormat::BINARY32) ? "Float64" : "Float32";
133 std::vector<char> buf;
134
135 const int np = qspace->GetSize();
136 const int ne = qspace->GetNE();
137 const int sdim = qspace->GetMesh()->SpaceDimension();
138
139 // For quadrature functions, each point is a vertex cell, so number of cells
140 // is equal to number of points
141 os << "<Piece NumberOfPoints=\"" << np
142 << "\" NumberOfCells=\"" << np << "\">\n";
143
144 // print out the points
145 os << "<Points>\n";
146 os << "<DataArray type=\"" << type_str
147 << "\" NumberOfComponents=\"3\" format=\"" << fmt_str << "\">\n";
148
149 Vector pt(sdim);
150 for (int i = 0; i < ne; i++)
151 {
152 ElementTransformation &T = *qspace->GetTransformation(i);
153 const IntegrationRule &ir = GetIntRule(i);
154 for (int j = 0; j < ir.Size(); j++)
155 {
156 T.Transform(ir[j], pt);
157 WriteBinaryOrASCII(os, buf, pt[0], " ", format);
158 if (sdim > 1) { WriteBinaryOrASCII(os, buf, pt[1], " ", format); }
159 else { WriteBinaryOrASCII(os, buf, 0.0, " ", format); }
160 if (sdim > 2) { WriteBinaryOrASCII(os, buf, pt[2], "", format); }
161 else { WriteBinaryOrASCII(os, buf, 0.0, "", format); }
162 if (format == VTKFormat::ASCII) { os << '\n'; }
163 }
164 }
165 if (format != VTKFormat::ASCII)
166 {
167 WriteBase64WithSizeAndClear(os, buf, compression_level);
168 }
169 os << "</DataArray>\n";
170 os << "</Points>\n";
171
172 // Write cells (each cell is just a vertex)
173 os << "<Cells>\n";
174 // Connectivity
175 os << R"(<DataArray type="Int32" Name="connectivity" format=")"
176 << fmt_str << "\">\n";
177
178 for (int i=0; i<np; ++i) { WriteBinaryOrASCII(os, buf, i, "\n", format); }
179 if (format != VTKFormat::ASCII)
180 {
181 WriteBase64WithSizeAndClear(os, buf, compression_level);
182 }
183 os << "</DataArray>\n";
184 // Offsets
185 os << R"(<DataArray type="Int32" Name="offsets" format=")"
186 << fmt_str << "\">\n";
187 for (int i=0; i<np; ++i) { WriteBinaryOrASCII(os, buf, i, "\n", format); }
188 if (format != VTKFormat::ASCII)
189 {
190 WriteBase64WithSizeAndClear(os, buf, compression_level);
191 }
192 os << "</DataArray>\n";
193 // Types
194 os << R"(<DataArray type="UInt8" Name="types" format=")"
195 << fmt_str << "\">\n";
196 for (int i = 0; i < np; i++)
197 {
198 uint8_t vtk_cell_type = VTKGeometry::POINT;
199 WriteBinaryOrASCII(os, buf, vtk_cell_type, "\n", format);
200 }
201 if (format != VTKFormat::ASCII)
202 {
203 WriteBase64WithSizeAndClear(os, buf, compression_level);
204 }
205 os << "</DataArray>\n";
206 os << "</Cells>\n";
207
208 os << "<PointData>\n";
209 os << "<DataArray type=\"" << type_str << "\" Name=\"" << field_name
210 << "\" format=\"" << fmt_str << "\" NumberOfComponents=\"" << vdim
211 << "\">\n";
212 for (int i = 0; i < ne; i++)
213 {
214 DenseMatrix vals;
215 GetValues(i, vals);
216 for (int j = 0; j < vals.Size(); ++j)
217 {
218 for (int vd = 0; vd < vdim; ++vd)
219 {
220 WriteBinaryOrASCII(os, buf, vals(vd, j), " ", format);
221 }
222 if (format == VTKFormat::ASCII) { os << '\n'; }
223 }
224 }
225 if (format != VTKFormat::ASCII)
226 {
227 WriteBase64WithSizeAndClear(os, buf, compression_level);
228 }
229 os << "</DataArray>\n";
230 os << "</PointData>\n";
231
232 os << "</Piece>\n";
233 os << "</UnstructuredGrid>\n";
234 os << "</VTKFile>" << std::endl;
235}
236
237void QuadratureFunction::SaveVTU(const std::string &filename, VTKFormat format,
238 int compression_level,
239 const std::string &field_name) const
240{
241 std::ofstream f(filename + ".vtu");
242 SaveVTU(f, format, compression_level, field_name);
243}
244
245static real_t ReduceReal(const Mesh *mesh, real_t value)
246{
247#ifdef MFEM_USE_MPI
248 if (auto *pmesh = dynamic_cast<const ParMesh*>(mesh))
249 {
250 MPI_Comm comm = pmesh->GetComm();
251 MPI_Allreduce(MPI_IN_PLACE, &value, 1, MPITypeMap<real_t>::mpi_type,
252 MPI_SUM, comm);
253 }
254#endif
255 return value;
256}
257
258real_t QuadratureFunction::Integrate() const
259{
260 MFEM_VERIFY(vdim == 1, "Only scalar functions are supported.")
261 const real_t local_integral = InnerProduct(*this, qspace->GetWeights());
262 return ReduceReal(qspace->GetMesh(), local_integral);
263}
264
265void QuadratureFunction::Integrate(Vector &integrals) const
266{
267 integrals.SetSize(vdim);
268
269 const Vector &weights = qspace->GetWeights();
270 QuadratureFunction component(qspace);
271 const int N = component.Size();
272 const int VDIM = vdim; // avoid capturing 'this' in lambda body
273 const real_t *d_v = Read();
274
275 for (int vd = 0; vd < vdim; ++vd)
276 {
277 // Extract the component 'vd' into component.
278 real_t *d_c = component.Write();
279 mfem::forall(N, [=] MFEM_HOST_DEVICE (int i)
280 {
281 d_c[i] = d_v[vd + i*VDIM];
282 });
283 integrals[vd] = ReduceReal(qspace->GetMesh(),
284 InnerProduct(component, weights));
285 }
286}
287
288}
QuadratureFunction()
Default constructor, results in an empty vector.
Definition qfunction.hpp:32
QuadratureFunction & operator=(real_t value)
Set this equal to a constant value.
QuadratureSpaceBase * qspace
Associated QuadratureSpaceBase object.
Definition qfunction.hpp:26
Vector & operator=(const real_t *v)
Copy Size() entries from v.
Definition vector.cpp:118
ostream & operator<<(ostream &v, void(*f)(VisMan &))
Definition ex17.cpp:406
Mesh * GetMesh(int type)
Definition ex29.cpp:218
real_t f(const Vector &p)
const T * Read(const Memory< T > &mem, int size, bool on_dev=true)
Get a pointer for read access to mem with the mfem::Device's DeviceMemoryClass, if on_dev = true,...
Definition device.hpp:320
real_t InnerProduct(HypreParVector *x, HypreParVector *y)
Definition hypre.cpp:439
void WriteBase64WithSizeAndClear(std::ostream &os, std::vector< char > &buf, int compression_level)
Encode in base 64 (and potentially compress) the given data, write it to the output stream (with a he...
Definition vtk.cpp:654
bool UsesTensorBasis(const FiniteElementSpace &fes)
Return true if the mesh contains only one topology and the elements are tensor elements.
Definition fespace.hpp:1345
void WriteBinaryOrASCII(std::ostream &os, std::vector< char > &buf, const T &val, const char *suffix, VTKFormat format)
Write either ASCII data to the stream or binary data to the buffer depending on the given format.
Definition vtk.hpp:147
float real_t
Definition config.hpp:43
const char * VTKByteOrder()
Determine the byte order and return either "BigEndian" or "LittleEndian".
Definition vtk.cpp:602
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition fespace.hpp:75
void forall(int N, lambda &&body)
Definition forall.hpp:754
FaceType
Definition mesh.hpp:47