29 MFEM_ASSERT(
qspace && v.Size() == this->Size(),
"");
36 : QuadratureFunction()
38 const char *msg =
"invalid input stream";
41 qspace =
new QuadratureSpace(mesh, in);
44 in >> ident; MFEM_VERIFY(ident ==
"VDim:", msg);
47 Load(in, vdim*qspace->GetSize());
50void QuadratureFunction::Save(std::ostream &os)
const
53 os <<
"VDim: " << vdim <<
'\n'
55 Vector::Print(os, vdim);
59void QuadratureFunction::ProjectGridFunction(
const GridFunction &gf)
61 SetVDim(gf.VectorDim());
63 if (
auto *qs_elem =
dynamic_cast<QuadratureSpace*
>(qspace))
65 const FiniteElementSpace &gf_fes = *gf.FESpace();
68 ElementDofOrdering::LEXICOGRAPHIC :
69 ElementDofOrdering::NATIVE;
72 const Operator *R = gf_fes.GetElementRestriction(ordering);
73 Vector e_vec(R->Height());
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);
83 else if (
auto *qs_face =
dynamic_cast<FaceQuadratureSpace*
>(qspace))
85 const FiniteElementSpace &gf_fes = *gf.FESpace();
88 ElementDofOrdering::LEXICOGRAPHIC :
89 ElementDofOrdering::NATIVE;
91 const FaceType face_type = qs_face->GetFaceType();
94 const Operator *R = gf_fes.GetFaceRestriction(
95 ordering, face_type, L2FaceValues::SingleValued);
96 Vector e_vec(R->Height());
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);
109 MFEM_ABORT(
"Unsupported case.");
113std::ostream &
operator<<(std::ostream &os,
const QuadratureFunction &qf)
119void QuadratureFunction::SaveVTU(std::ostream &os, VTKFormat format,
120 int compression_level,
121 const std::string &field_name)
const
123 os << R
"(<VTKFile type="UnstructuredGrid" version="0.1")";
124 if (compression_level != 0)
126 os << R
"( compressor="vtkZLibDataCompressor")";
129 os <<
"<UnstructuredGrid>\n";
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;
135 const int np = qspace->GetSize();
136 const int ne = qspace->GetNE();
137 const int sdim = qspace->GetMesh()->SpaceDimension();
141 os <<
"<Piece NumberOfPoints=\"" << np
142 <<
"\" NumberOfCells=\"" << np <<
"\">\n";
146 os <<
"<DataArray type=\"" << type_str
147 <<
"\" NumberOfComponents=\"3\" format=\"" << fmt_str <<
"\">\n";
150 for (
int i = 0; i < ne; i++)
152 ElementTransformation &T = *qspace->GetTransformation(i);
153 const IntegrationRule &ir = GetIntRule(i);
154 for (
int j = 0; j < ir.Size(); j++)
156 T.Transform(ir[j], pt);
162 if (format == VTKFormat::ASCII) { os <<
'\n'; }
165 if (format != VTKFormat::ASCII)
169 os <<
"</DataArray>\n";
175 os << R
"(<DataArray type="Int32" Name="connectivity" format=")"
176 << fmt_str << "\">\n";
179 if (format != VTKFormat::ASCII)
183 os <<
"</DataArray>\n";
185 os << R
"(<DataArray type="Int32" Name="offsets" format=")"
186 << fmt_str << "\">\n";
188 if (format != VTKFormat::ASCII)
192 os <<
"</DataArray>\n";
194 os << R
"(<DataArray type="UInt8" Name="types" format=")"
195 << fmt_str << "\">\n";
196 for (
int i = 0; i < np; i++)
198 uint8_t vtk_cell_type = VTKGeometry::POINT;
201 if (format != VTKFormat::ASCII)
205 os <<
"</DataArray>\n";
208 os <<
"<PointData>\n";
209 os <<
"<DataArray type=\"" << type_str <<
"\" Name=\"" << field_name
210 <<
"\" format=\"" << fmt_str <<
"\" NumberOfComponents=\"" << vdim
212 for (
int i = 0; i < ne; i++)
216 for (
int j = 0; j < vals.Size(); ++j)
218 for (
int vd = 0; vd < vdim; ++vd)
222 if (format == VTKFormat::ASCII) { os <<
'\n'; }
225 if (format != VTKFormat::ASCII)
229 os <<
"</DataArray>\n";
230 os <<
"</PointData>\n";
233 os <<
"</UnstructuredGrid>\n";
234 os <<
"</VTKFile>" << std::endl;
237void QuadratureFunction::SaveVTU(
const std::string &filename, VTKFormat format,
238 int compression_level,
239 const std::string &field_name)
const
241 std::ofstream
f(filename +
".vtu");
242 SaveVTU(
f, format, compression_level, field_name);
245static real_t ReduceReal(
const Mesh *mesh, real_t value)
248 if (
auto *pmesh =
dynamic_cast<const ParMesh*
>(mesh))
250 MPI_Comm comm = pmesh->GetComm();
251 MPI_Allreduce(MPI_IN_PLACE, &value, 1, MPITypeMap<real_t>::mpi_type,
258real_t QuadratureFunction::Integrate()
const
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);
265void QuadratureFunction::Integrate(Vector &integrals)
const
267 integrals.SetSize(vdim);
269 const Vector &weights = qspace->GetWeights();
270 QuadratureFunction component(qspace);
271 const int N = component.Size();
272 const int VDIM = vdim;
275 for (
int vd = 0; vd < vdim; ++vd)
278 real_t *d_c = component.Write();
281 d_c[i] = d_v[vd + i*VDIM];
283 integrals[vd] = ReduceReal(qspace->GetMesh(),
QuadratureFunction()
Default constructor, results in an empty vector.
QuadratureFunction & operator=(real_t value)
Set this equal to a constant value.
QuadratureSpaceBase * qspace
Associated QuadratureSpaceBase object.
Vector & operator=(const real_t *v)
Copy Size() entries from v.
ostream & operator<<(ostream &v, void(*f)(VisMan &))
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,...
real_t InnerProduct(HypreParVector *x, HypreParVector *y)
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...
bool UsesTensorBasis(const FiniteElementSpace &fes)
Return true if the mesh contains only one topology and the elements are tensor elements.
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.
const char * VTKByteOrder()
Determine the byte order and return either "BigEndian" or "LittleEndian".
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
void forall(int N, lambda &&body)