version 3.8
pointcloudvtkwriter.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_POINTCLOUD_VTK_WRITER_HH
13#define DUMUX_POINTCLOUD_VTK_WRITER_HH
14
15#include <string>
16#include <vector>
17#include <list>
18#include <fstream>
19#include <iomanip>
20
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>
25
26namespace Dumux {
27
37template<class Scalar, class GlobalPosition>
39{
40 // GlobalPosition is used for the point coordinates, DimWorldVector for the actual data.
41 // GlobalPosition's ctype does not necessarily equal Scalar.
42 using DimWorldVector = Dune::FieldVector<Scalar, GlobalPosition::size()>;
43
44 static constexpr unsigned int precision = 6;
45 static constexpr unsigned int numBeforeLineBreak = 15;
46
50 template<class ContainerType>
51 class VTKFunction
52 {
53 public:
61 VTKFunction(const ContainerType& data, const std::string& name, const int numComponents) : data_(data), name_(name), numComponents_(numComponents)
62 {}
63
67 const std::string& name() const
68 {
69 return name_;
70 }
71
75 int numComponents() const
76 {
77 return numComponents_;
78 }
79
85 auto& operator() (const int idx) const { return data_[idx]; }
86
87 decltype(auto) begin() const
88 {
89 return data_.begin();
90 }
91
92 decltype(auto) end() const
93 {
94 return data_.end();
95 }
96
100 int size() const
101 {
102 return data_.size();
103 }
104
105 private:
106 const ContainerType& data_;
107 const std::string name_;
108 const int numComponents_;
109 };
110
111
112public:
113 using ScalarFunction = VTKFunction<std::vector<Scalar>>;
114 using VectorFunction = VTKFunction<std::vector<DimWorldVector>>;
115
116
117 PointCloudVtkWriter(const std::vector<GlobalPosition>& coordinates) : coordinates_(coordinates)
118 {}
119
126 void write(const std::string& name, Dune::VTK::OutputType type = Dune::VTK::ascii)
127 {
128 auto filename = getSerialPieceName(name, "");
129
130 std::ofstream file;
131 file.open(filename);
132 writeHeader_(file);
133 writeCoordinates_(file, coordinates_);
134 writeDataInfo_(file);
135
136 for (auto&& data : scalarPointData_)
137 writeData_(file, data);
138
139 for (auto&& data :vectorPointData_)
140 writeData_(file, data);
141
142 if (!scalarPointData_.empty() || !vectorPointData_.empty())
143 file << "</PointData>\n";
144
145 file << "</Piece>\n";
146 file << "</PolyData>\n";
147 file << "</VTKFile>";
148
149 clear();
150 file.close();
151 }
152
168 void pwrite(const std::string & name, const std::string & path, const std::string & extendpath,
169 Dune::VTK::OutputType type = Dune::VTK::ascii)
170 {
171 DUNE_THROW(Dune::NotImplemented, "parallel point cloud vtk output not supported yet");
172 }
173
180 void addPointData(const std::vector<Scalar>& v, const std::string &name)
181 {
182 assert(v.size() == coordinates_.size());
183 scalarPointData_.push_back(ScalarFunction(v, name, 1));
184 }
185
192 void addPointData(const std::vector<DimWorldVector>& v, const std::string &name)
193 {
194 assert(v.size() == coordinates_.size());
195 vectorPointData_.push_back(VectorFunction(v, name, 3));
196 }
197
201 void clear()
202 {
203 scalarPointData_.clear();
204 vectorPointData_.clear();
205 }
206
218 std::string getSerialPieceName(const std::string& name,
219 const std::string& path) const
220 {
221 static const std::string extension = ".vtp";
222
223 return Dune::concatPaths(path, name+extension);
224 }
225
238 std::string getParallelHeaderName(const std::string& name,
239 const std::string& path,
240 int commSize) const
241 {
242 std::ostringstream s;
243 if(path.size() > 0) {
244 s << path;
245 if(path[path.size()-1] != '/')
246 s << '/';
247 }
248 s << 's' << std::setw(4) << std::setfill('0') << commSize << '-';
249 s << name;
250 s << ".pvtp";
251 return s.str();
252 }
253
254
255private:
259 void writeHeader_(std::ostream& file)
260 {
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";
265 file << header;
266 }
267
271 void writeDataInfo_(std::ostream& file)
272 {
273 std::string scalarName;
274 std::string vectorName;
275 bool foundScalar = false;
276 bool foundVector = false;
277
278 for(auto&& data : scalarPointData_)
279 {
280 if(data.numComponents() == 1 && !foundScalar)
281 {
282 scalarName = data.name();
283 foundScalar = true;
284 continue;
285 }
286
287 if(data.numComponents() > 1 && !foundVector)
288 {
289 vectorName = data.name();
290 foundVector = true;
291 }
292 }
293
294 for(auto&& data : vectorPointData_)
295 {
296 if(data.numComponents() > 1 && !foundVector)
297 {
298 vectorName = data.name();
299 foundVector = true;
300 }
301 }
302
303 if(foundScalar)
304 if(foundVector)
305 file << "<PointData Scalars=\"" << scalarName << "\" Vectors=\"" << vectorName <<"\">\n";
306 else
307 file << "<PointData Scalars=\"" << scalarName << "\">\n";
308 else if(foundVector)
309 file << "<PointData Vectors=\"" << vectorName << "\">\n";
310 else
311 return;
312 }
313
320 void writeCoordinates_(std::ostream& file, const std::vector<GlobalPosition>& positions)
321 {
322 // write the positions to the file
323 file << "<Points>\n";
324 file << "<DataArray type=\"Float32\" Name=\"Coordinates\" NumberOfComponents=\"3\" format=\"ascii\">\n";
325 int counter = 0;
326 for(auto&& x : positions)
327 {
328 file << x ;
329
330 if(x.size() == 1)
331 file << " 0 0 ";
332 if(x.size() == 2)
333 file << " 0 ";
334 if(x.size() == 3)
335 file << " ";
336
337 // introduce a line break after a certain time
338 if((++counter) > numBeforeLineBreak)
339 {
340 file << std::endl;
341 counter = 0;
342 }
343 }
344 file << "\n</DataArray>\n";
345 file << "</Points>\n";
346 }
347
354 template<class T>
355 void writeData_(std::ostream& file, const T& data)
356 {
357 file << "<DataArray type=\"Float32\" Name=\"" << data.name() << "\" NumberOfComponents=\"" << data.numComponents() << "\" format=\"ascii\">\n";
358 int counter = 0;
359 for(auto&& value : data)
360 {
361 // forward to specialized function
362 writeToFile_(file, value);
363
364 // introduce a line break after a certain time
365 if((++counter) > numBeforeLineBreak)
366 {
367 file << std::endl;
368 counter = 0;
369 }
370 }
371 file << "\n</DataArray>\n";
372 }
373
380 void writeToFile_(std::ostream& file, const Scalar& s)
381 {
382 file << s << " ";
383 }
384
391 void writeToFile_(std::ostream& file, const DimWorldVector& g)
392 {
393 assert(g.size() > 1 && g.size() < 4);
394 if(g.size() < 3)
395 file << g << " 0 ";
396 else
397 file << g << " ";
398 }
399
400 const std::vector<GlobalPosition>& coordinates_;
401 std::list<ScalarFunction> scalarPointData_;
402 std::list<VectorFunction> vectorPointData_;
403};
404} // end namespace Dumux
405
406#endif
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
Definition: adapt.hh:17