12#ifndef DUMUX_POINTCLOUD_VTK_WRITER_HH
13#define DUMUX_POINTCLOUD_VTK_WRITER_HH
21#include <dune/common/fvector.hh>
22#include <dune/common/exceptions.hh>
23#include <dune/common/path.hh>
24#include <dune/grid/io/file/vtk/common.hh>
37template<
class Scalar,
class GlobalPosition>
42 using DimWorldVector = Dune::FieldVector<Scalar, GlobalPosition::size()>;
44 static constexpr unsigned int precision = 6;
45 static constexpr unsigned int numBeforeLineBreak = 15;
50 template<
class ContainerType>
61 VTKFunction(
const ContainerType& data,
const std::string& name,
const int numComponents) : data_(data), name_(name), numComponents_(numComponents)
67 const std::string& name()
const
75 int numComponents()
const
77 return numComponents_;
85 auto& operator() (
const int idx)
const {
return data_[idx]; }
87 decltype(
auto) begin()
const
92 decltype(
auto) end()
const
106 const ContainerType& data_;
107 const std::string name_;
108 const int numComponents_;
126 void write(
const std::string& name, Dune::VTK::OutputType type = Dune::VTK::ascii)
133 writeCoordinates_(file, coordinates_);
134 writeDataInfo_(file);
136 for (
auto&& data : scalarPointData_)
137 writeData_(file, data);
139 for (
auto&& data :vectorPointData_)
140 writeData_(file, data);
142 if (!scalarPointData_.empty() || !vectorPointData_.empty())
143 file <<
"</PointData>\n";
145 file <<
"</Piece>\n";
146 file <<
"</PolyData>\n";
147 file <<
"</VTKFile>";
168 void pwrite(
const std::string & name,
const std::string & path,
const std::string & extendpath,
169 Dune::VTK::OutputType type = Dune::VTK::ascii)
171 DUNE_THROW(Dune::NotImplemented,
"parallel point cloud vtk output not supported yet");
180 void addPointData(
const std::vector<Scalar>& v,
const std::string &name)
182 assert(v.size() == coordinates_.size());
192 void addPointData(
const std::vector<DimWorldVector>& v,
const std::string &name)
194 assert(v.size() == coordinates_.size());
203 scalarPointData_.clear();
204 vectorPointData_.clear();
219 const std::string& path)
const
221 static const std::string extension =
".vtp";
223 return Dune::concatPaths(path, name+extension);
239 const std::string& path,
242 std::ostringstream s;
243 if(path.size() > 0) {
245 if(path[path.size()-1] !=
'/')
248 s <<
's' << std::setw(4) << std::setfill(
'0') << commSize <<
'-';
259 void writeHeader_(std::ostream& file)
261 std::string header =
"<?xml version=\"1.0\"?>\n";
262 header +=
"<VTKFile type=\"PolyData\" version=\"0.1\" byte_order=\"LittleEndian\">\n";
263 header +=
"<PolyData>\n";
264 header +=
"<Piece NumberOfLines=\"0\" NumberOfPoints=\"" + std::to_string(coordinates_.size()) +
"\">\n";
271 void writeDataInfo_(std::ostream& file)
273 std::string scalarName;
274 std::string vectorName;
275 bool foundScalar =
false;
276 bool foundVector =
false;
278 for(
auto&& data : scalarPointData_)
280 if(data.numComponents() == 1 && !foundScalar)
282 scalarName = data.name();
287 if(data.numComponents() > 1 && !foundVector)
289 vectorName = data.name();
294 for(
auto&& data : vectorPointData_)
296 if(data.numComponents() > 1 && !foundVector)
298 vectorName = data.name();
305 file <<
"<PointData Scalars=\"" << scalarName <<
"\" Vectors=\"" << vectorName <<
"\">\n";
307 file <<
"<PointData Scalars=\"" << scalarName <<
"\">\n";
309 file <<
"<PointData Vectors=\"" << vectorName <<
"\">\n";
320 void writeCoordinates_(std::ostream& file,
const std::vector<GlobalPosition>& positions)
323 file <<
"<Points>\n";
324 file <<
"<DataArray type=\"Float32\" Name=\"Coordinates\" NumberOfComponents=\"3\" format=\"ascii\">\n";
326 for(
auto&& x : positions)
338 if((++counter) > numBeforeLineBreak)
344 file <<
"\n</DataArray>\n";
345 file <<
"</Points>\n";
355 void writeData_(std::ostream& file,
const T& data)
357 file <<
"<DataArray type=\"Float32\" Name=\"" << data.name() <<
"\" NumberOfComponents=\"" << data.numComponents() <<
"\" format=\"ascii\">\n";
359 for(
auto&& value : data)
362 writeToFile_(file, value);
365 if((++counter) > numBeforeLineBreak)
371 file <<
"\n</DataArray>\n";
380 void writeToFile_(std::ostream& file,
const Scalar& s)
391 void writeToFile_(std::ostream& file,
const DimWorldVector& g)
393 assert(g.size() > 1 && g.size() < 4);
400 const std::vector<GlobalPosition>& coordinates_;
401 std::list<ScalarFunction> scalarPointData_;
402 std::list<VectorFunction> vectorPointData_;
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: pointcloudvtkwriter.hh:39
std::string getSerialPieceName(const std::string &name, const std::string &path) const
Return name of a serial header file.
Definition: pointcloudvtkwriter.hh:218
VTKFunction< std::vector< Scalar > > ScalarFunction
Definition: pointcloudvtkwriter.hh:113
std::string getParallelHeaderName(const std::string &name, const std::string &path, int commSize) const
Return name of a parallel header file.
Definition: pointcloudvtkwriter.hh:238
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:192
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:180
VTKFunction< std::vector< DimWorldVector > > VectorFunction
Definition: pointcloudvtkwriter.hh:114
PointCloudVtkWriter(const std::vector< GlobalPosition > &coordinates)
Definition: pointcloudvtkwriter.hh:117
void write(const std::string &name, Dune::VTK::OutputType type=Dune::VTK::ascii)
Create an output file.
Definition: pointcloudvtkwriter.hh:126
void clear()
Clears the data lists.
Definition: pointcloudvtkwriter.hh:201
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:168