24#ifndef DUMUX_POINTCLOUD_VTK_WRITER_HH
25#define DUMUX_POINTCLOUD_VTK_WRITER_HH
33#include <dune/common/fvector.hh>
34#include <dune/common/exceptions.hh>
35#include <dune/common/path.hh>
36#include <dune/grid/io/file/vtk/common.hh>
49template<
class Scalar,
class GlobalPosition>
54 using DimWorldVector = Dune::FieldVector<Scalar, GlobalPosition::size()>;
56 static constexpr unsigned int precision = 6;
57 static constexpr unsigned int numBeforeLineBreak = 15;
62 template<
class ContainerType>
73 VTKFunction(
const ContainerType& data,
const std::string& name,
const int numComponents) : data_(data), name_(name), numComponents_(numComponents)
79 const std::string& name()
const
87 int numComponents()
const
89 return numComponents_;
97 auto& operator() (
const int idx)
const {
return data_[idx]; }
99 decltype(
auto) begin()
const
101 return data_.begin();
104 decltype(
auto) end()
const
118 const ContainerType& data_;
119 const std::string name_;
120 const int numComponents_;
138 void write(
const std::string& name, Dune::VTK::OutputType type = Dune::VTK::ascii)
145 writeCoordinates_(file, coordinates_);
146 writeDataInfo_(file);
148 for (
auto&& data : scalarPointData_)
149 writeData_(file, data);
151 for (
auto&& data :vectorPointData_)
152 writeData_(file, data);
154 if (!scalarPointData_.empty() || !vectorPointData_.empty())
155 file <<
"</PointData>\n";
157 file <<
"</Piece>\n";
158 file <<
"</PolyData>\n";
159 file <<
"</VTKFile>";
180 void pwrite(
const std::string & name,
const std::string & path,
const std::string & extendpath,
181 Dune::VTK::OutputType type = Dune::VTK::ascii)
183 DUNE_THROW(Dune::NotImplemented,
"parallel point cloud vtk output not supported yet");
192 void addPointData(
const std::vector<Scalar>& v,
const std::string &name)
194 assert(v.size() == coordinates_.size());
204 void addPointData(
const std::vector<DimWorldVector>& v,
const std::string &name)
206 assert(v.size() == coordinates_.size());
215 scalarPointData_.clear();
216 vectorPointData_.clear();
231 const std::string& path)
const
233 static const std::string extension =
".vtp";
235 return Dune::concatPaths(path, name+extension);
251 const std::string& path,
254 std::ostringstream s;
255 if(path.size() > 0) {
257 if(path[path.size()-1] !=
'/')
260 s <<
's' << std::setw(4) << std::setfill(
'0') << commSize <<
'-';
271 void writeHeader_(std::ostream& file)
273 std::string header =
"<?xml version=\"1.0\"?>\n";
274 header +=
"<VTKFile type=\"PolyData\" version=\"0.1\" byte_order=\"LittleEndian\">\n";
275 header +=
"<PolyData>\n";
276 header +=
"<Piece NumberOfLines=\"0\" NumberOfPoints=\"" + std::to_string(coordinates_.size()) +
"\">\n";
283 void writeDataInfo_(std::ostream& file)
285 std::string scalarName;
286 std::string vectorName;
287 bool foundScalar =
false;
288 bool foundVector =
false;
290 for(
auto&& data : scalarPointData_)
292 if(data.numComponents() == 1 && !foundScalar)
294 scalarName = data.name();
299 if(data.numComponents() > 1 && !foundVector)
301 vectorName = data.name();
306 for(
auto&& data : vectorPointData_)
308 if(data.numComponents() > 1 && !foundVector)
310 vectorName = data.name();
317 file <<
"<PointData Scalars=\"" << scalarName <<
"\" Vectors=\"" << vectorName <<
"\">\n";
319 file <<
"<PointData Scalars=\"" << scalarName <<
"\">\n";
321 file <<
"<PointData Vectors=\"" << vectorName <<
"\">\n";
332 void writeCoordinates_(std::ostream& file,
const std::vector<GlobalPosition>& positions)
335 file <<
"<Points>\n";
336 file <<
"<DataArray type=\"Float32\" Name=\"Coordinates\" NumberOfComponents=\"3\" format=\"ascii\">\n";
338 for(
auto&& x : positions)
350 if((++counter) > numBeforeLineBreak)
356 file <<
"\n</DataArray>\n";
357 file <<
"</Points>\n";
367 void writeData_(std::ostream& file,
const T& data)
369 file <<
"<DataArray type=\"Float32\" Name=\"" << data.name() <<
"\" NumberOfComponents=\"" << data.numComponents() <<
"\" format=\"ascii\">\n";
371 for(
auto&& value : data)
374 writeToFile_(file, value);
377 if((++counter) > numBeforeLineBreak)
383 file <<
"\n</DataArray>\n";
392 void writeToFile_(std::ostream& file,
const Scalar& s)
403 void writeToFile_(std::ostream& file,
const DimWorldVector& g)
405 assert(g.size() > 1 && g.size() < 4);
412 const std::vector<GlobalPosition>& coordinates_;
413 std::list<ScalarFunction> scalarPointData_;
414 std::list<VectorFunction> vectorPointData_;
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: pointcloudvtkwriter.hh:51
std::string getSerialPieceName(const std::string &name, const std::string &path) const
Return name of a serial header file.
Definition: pointcloudvtkwriter.hh:230
VTKFunction< std::vector< Scalar > > ScalarFunction
Definition: pointcloudvtkwriter.hh:125
std::string getParallelHeaderName(const std::string &name, const std::string &path, int commSize) const
Return name of a parallel header file.
Definition: pointcloudvtkwriter.hh:250
void addPointData(const std::vector< DimWorldVector > &v, const std::string &name)
Add a vector of vector data that live on arbitrary points to the visualization.
Definition: pointcloudvtkwriter.hh:204
void addPointData(const std::vector< Scalar > &v, const std::string &name)
Add a vector of scalar data that live on arbitrary points to the visualization.
Definition: pointcloudvtkwriter.hh:192
VTKFunction< std::vector< DimWorldVector > > VectorFunction
Definition: pointcloudvtkwriter.hh:126
PointCloudVtkWriter(const std::vector< GlobalPosition > &coordinates)
Definition: pointcloudvtkwriter.hh:129
void write(const std::string &name, Dune::VTK::OutputType type=Dune::VTK::ascii)
Create an output file.
Definition: pointcloudvtkwriter.hh:138
void clear()
Clears the data lists.
Definition: pointcloudvtkwriter.hh:213
void pwrite(const std::string &name, const std::string &path, const std::string &extendpath, Dune::VTK::OutputType type=Dune::VTK::ascii)
Create an output file in parallel.
Definition: pointcloudvtkwriter.hh:180