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::ProjectGridFunctionFallback(
const GridFunction &gf)
61 if (gf.VectorDim() == 1)
63 GridFunctionCoefficient coeff(&gf);
64 coeff.Coefficient::Project(*
this);
68 VectorGridFunctionCoefficient coeff(&gf);
69 coeff.VectorCoefficient::Project(*
this);
73void QuadratureFunction::ProjectGridFunction(
const GridFunction &gf)
75 SetVDim(gf.VectorDim());
77 if (
auto *qs_elem =
dynamic_cast<QuadratureSpace*
>(qspace))
79 const FiniteElementSpace &gf_fes = *gf.FESpace();
82 ElementDofOrdering::LEXICOGRAPHIC :
83 ElementDofOrdering::NATIVE;
86 const QuadratureInterpolator *qi =
87 gf_fes.GetQuadratureInterpolator(*qs_elem);
93 ProjectGridFunctionFallback(gf);
98 const Operator *R = gf_fes.GetElementRestriction(ordering);
99 Vector e_vec(R->Height());
102 qi->SetOutputLayout(QVectorLayout::byVDIM);
103 qi->DisableTensorProducts(!use_tensor_products);
104 qi->PhysValues(e_vec, *
this);
106 else if (
auto *qs_face =
dynamic_cast<FaceQuadratureSpace*
>(qspace))
108 const FiniteElementSpace &gf_fes = *gf.FESpace();
109 const FaceType face_type = qs_face->GetFaceType();
112 ElementDofOrdering::LEXICOGRAPHIC :
113 ElementDofOrdering::NATIVE;
116 const FaceQuadratureInterpolator *qi =
117 gf_fes.GetFaceQuadratureInterpolator(qspace->GetIntRule(0), face_type);
123 if (qi ==
nullptr || ordering == ElementDofOrdering::NATIVE)
125 ProjectGridFunctionFallback(gf);
130 const Operator *R = gf_fes.GetFaceRestriction(
131 ordering, face_type, L2FaceValues::SingleValued);
132 Vector e_vec(R->Height());
135 qi->SetOutputLayout(QVectorLayout::byVDIM);
136 qi->DisableTensorProducts(!use_tensor_products);
137 qi->Values(e_vec, *
this);
142 MFEM_ABORT(
"Unsupported case.");
146std::ostream &
operator<<(std::ostream &os,
const QuadratureFunction &qf)
152void QuadratureFunction::SaveVTU(std::ostream &os, VTKFormat format,
153 int compression_level,
154 const std::string &field_name)
const
156 os << R
"(<VTKFile type="UnstructuredGrid" version="2.2")";
157 if (compression_level != 0)
159 os << R
"( compressor="vtkZLibDataCompressor")";
162 os <<
"<UnstructuredGrid>\n";
164 const char *fmt_str = (format == VTKFormat::ASCII) ?
"ascii" :
"binary";
165 const char *type_str = (format != VTKFormat::BINARY32) ?
"Float64" :
"Float32";
166 std::vector<char> buf;
168 const int np = qspace->GetSize();
169 const int ne = qspace->GetNE();
170 const int sdim = qspace->GetMesh()->SpaceDimension();
174 os <<
"<Piece NumberOfPoints=\"" << np
175 <<
"\" NumberOfCells=\"" << np <<
"\">\n";
179 os <<
"<DataArray type=\"" << type_str
180 <<
"\" NumberOfComponents=\"3\" format=\"" << fmt_str <<
"\">\n";
183 for (
int i = 0; i < ne; i++)
185 ElementTransformation &T = *qspace->GetTransformation(i);
186 const IntegrationRule &ir = GetIntRule(i);
187 for (
int j = 0; j < ir.Size(); j++)
189 T.Transform(ir[j], pt);
195 if (format == VTKFormat::ASCII) { os <<
'\n'; }
198 if (format != VTKFormat::ASCII)
202 os <<
"</DataArray>\n";
208 os << R
"(<DataArray type="Int32" Name="connectivity" format=")"
209 << fmt_str << "\">\n";
212 if (format != VTKFormat::ASCII)
216 os <<
"</DataArray>\n";
218 os << R
"(<DataArray type="Int32" Name="offsets" format=")"
219 << fmt_str << "\">\n";
221 if (format != VTKFormat::ASCII)
225 os <<
"</DataArray>\n";
227 os << R
"(<DataArray type="UInt8" Name="types" format=")"
228 << fmt_str << "\">\n";
229 for (
int i = 0; i < np; i++)
231 uint8_t vtk_cell_type = VTKGeometry::POINT;
234 if (format != VTKFormat::ASCII)
238 os <<
"</DataArray>\n";
241 os <<
"<PointData>\n";
242 os <<
"<DataArray type=\"" << type_str <<
"\" Name=\"" << field_name
243 <<
"\" format=\"" << fmt_str <<
"\" NumberOfComponents=\"" << vdim <<
"\" "
246 for (
int i = 0; i < ne; i++)
250 for (
int j = 0; j < vals.Size(); ++j)
252 for (
int vd = 0; vd < vdim; ++vd)
256 if (format == VTKFormat::ASCII) { os <<
'\n'; }
259 if (format != VTKFormat::ASCII)
263 os <<
"</DataArray>\n";
264 os <<
"</PointData>\n";
267 os <<
"</UnstructuredGrid>\n";
268 os <<
"</VTKFile>" << std::endl;
271void QuadratureFunction::SaveVTU(
const std::string &filename, VTKFormat format,
272 int compression_level,
273 const std::string &field_name)
const
275 std::ofstream
f(filename +
".vtu");
276 SaveVTU(
f, format, compression_level, field_name);
282 if (
auto *pmesh =
dynamic_cast<const ParMesh*
>(mesh))
284 MPI_Comm comm = pmesh->GetComm();
285 MPI_Allreduce(MPI_IN_PLACE, &value, 1, MPITypeMap<real_t>::mpi_type,
292real_t QuadratureFunction::Integrate()
const
294 MFEM_VERIFY(vdim == 1,
"Only scalar functions are supported.")
295 const
real_t local_integral = InnerProduct(*this, qspace->GetWeights());
296 return ReduceReal(qspace->
GetMesh(), local_integral);
299void QuadratureFunction::Integrate(Vector &integrals)
const
301 integrals.SetSize(vdim);
303 const Vector &weights = qspace->GetWeights();
304 QuadratureFunction component(qspace);
305 const int N = component.Size();
306 const int VDIM = vdim;
309 for (
int vd = 0; vd < vdim; ++vd)
312 real_t *d_c = component.Write();
315 d_c[i] = d_v[vd + i*VDIM];
317 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)
std::string VTKComponentLabels(int vdim)
Returns a string defining the component labels for vector-valued data arrays for use in XML VTU files...