3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_POINTCLOUD_VTK_WRITER_HH
25#define DUMUX_POINTCLOUD_VTK_WRITER_HH
26
27#include <string>
28#include <vector>
29#include <list>
30#include <fstream>
31#include <iomanip>
32
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>
37
38namespace Dumux {
39
49template<class Scalar, class GlobalPosition>
51{
52 // GlobalPosition is used for the point coordinates, DimWorldVector for the actual data.
53 // GlobalPosition's ctype does not necessarily equal Scalar.
54 using DimWorldVector = Dune::FieldVector<Scalar, GlobalPosition::size()>;
55
56 static constexpr unsigned int precision = 6;
57 static constexpr unsigned int numBeforeLineBreak = 15;
58
62 template<class ContainerType>
63 class VTKFunction
64 {
65 public:
73 VTKFunction(const ContainerType& data, const std::string& name, const int numComponents) : data_(data), name_(name), numComponents_(numComponents)
74 {}
75
79 const std::string& name() const
80 {
81 return name_;
82 }
83
87 int numComponents() const
88 {
89 return numComponents_;
90 }
91
97 auto& operator() (const int idx) const { return data_[idx]; }
98
99 decltype(auto) begin() const
100 {
101 return data_.begin();
102 }
103
104 decltype(auto) end() const
105 {
106 return data_.end();
107 }
108
112 int size() const
113 {
114 return data_.size();
115 }
116
117 private:
118 const ContainerType& data_;
119 const std::string name_;
120 const int numComponents_;
121 };
122
123
124public:
125 using ScalarFunction = VTKFunction<std::vector<Scalar>>;
126 using VectorFunction = VTKFunction<std::vector<DimWorldVector>>;
127
128
129 PointCloudVtkWriter(const std::vector<GlobalPosition>& coordinates) : coordinates_(coordinates)
130 {}
131
138 void write(const std::string& name, Dune::VTK::OutputType type = Dune::VTK::ascii)
139 {
140 auto filename = getSerialPieceName(name, "");
141
142 std::ofstream file;
143 file.open(filename);
144 writeHeader_(file);
145 writeCoordinates_(file, coordinates_);
146 writeDataInfo_(file);
147
148 for (auto&& data : scalarPointData_)
149 writeData_(file, data);
150
151 for (auto&& data :vectorPointData_)
152 writeData_(file, data);
153
154 if (!scalarPointData_.empty() || !vectorPointData_.empty())
155 file << "</PointData>\n";
156
157 file << "</Piece>\n";
158 file << "</PolyData>\n";
159 file << "</VTKFile>";
160
161 clear();
162 file.close();
163 }
164
180 void pwrite(const std::string & name, const std::string & path, const std::string & extendpath,
181 Dune::VTK::OutputType type = Dune::VTK::ascii)
182 {
183 DUNE_THROW(Dune::NotImplemented, "parallel point cloud vtk output not supported yet");
184 }
185
192 void addPointData(const std::vector<Scalar>& v, const std::string &name)
193 {
194 assert(v.size() == coordinates_.size());
195 scalarPointData_.push_back(ScalarFunction(v, name, 1));
196 }
197
204 void addPointData(const std::vector<DimWorldVector>& v, const std::string &name)
205 {
206 assert(v.size() == coordinates_.size());
207 vectorPointData_.push_back(VectorFunction(v, name, 3));
208 }
209
213 void clear()
214 {
215 scalarPointData_.clear();
216 vectorPointData_.clear();
217 }
218
230 std::string getSerialPieceName(const std::string& name,
231 const std::string& path) const
232 {
233 static const std::string extension = ".vtp";
234
235 return Dune::concatPaths(path, name+extension);
236 }
237
250 std::string getParallelHeaderName(const std::string& name,
251 const std::string& path,
252 int commSize) const
253 {
254 std::ostringstream s;
255 if(path.size() > 0) {
256 s << path;
257 if(path[path.size()-1] != '/')
258 s << '/';
259 }
260 s << 's' << std::setw(4) << std::setfill('0') << commSize << '-';
261 s << name;
262 s << ".pvtp";
263 return s.str();
264 }
265
266
267private:
271 void writeHeader_(std::ostream& file)
272 {
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";
277 file << header;
278 }
279
283 void writeDataInfo_(std::ostream& file)
284 {
285 std::string scalarName;
286 std::string vectorName;
287 bool foundScalar = false;
288 bool foundVector = false;
289
290 for(auto&& data : scalarPointData_)
291 {
292 if(data.numComponents() == 1 && !foundScalar)
293 {
294 scalarName = data.name();
295 foundScalar = true;
296 continue;
297 }
298
299 if(data.numComponents() > 1 && !foundVector)
300 {
301 vectorName = data.name();
302 foundVector = true;
303 }
304 }
305
306 for(auto&& data : vectorPointData_)
307 {
308 if(data.numComponents() > 1 && !foundVector)
309 {
310 vectorName = data.name();
311 foundVector = true;
312 }
313 }
314
315 if(foundScalar)
316 if(foundVector)
317 file << "<PointData Scalars=\"" << scalarName << "\" Vectors=\"" << vectorName <<"\">\n";
318 else
319 file << "<PointData Scalars=\"" << scalarName << "\">\n";
320 else if(foundVector)
321 file << "<PointData Vectors=\"" << vectorName << "\">\n";
322 else
323 return;
324 }
325
332 void writeCoordinates_(std::ostream& file, const std::vector<GlobalPosition>& positions)
333 {
334 // write the positions to the file
335 file << "<Points>\n";
336 file << "<DataArray type=\"Float32\" Name=\"Coordinates\" NumberOfComponents=\"3\" format=\"ascii\">\n";
337 int counter = 0;
338 for(auto&& x : positions)
339 {
340 file << x ;
341
342 if(x.size() == 1)
343 file << " 0 0 ";
344 if(x.size() == 2)
345 file << " 0 ";
346 if(x.size() == 3)
347 file << " ";
348
349 // introduce a line break after a certain time
350 if((++counter) > numBeforeLineBreak)
351 {
352 file << std::endl;
353 counter = 0;
354 }
355 }
356 file << "\n</DataArray>\n";
357 file << "</Points>\n";
358 }
359
366 template<class T>
367 void writeData_(std::ostream& file, const T& data)
368 {
369 file << "<DataArray type=\"Float32\" Name=\"" << data.name() << "\" NumberOfComponents=\"" << data.numComponents() << "\" format=\"ascii\">\n";
370 int counter = 0;
371 for(auto&& value : data)
372 {
373 // forward to specialized function
374 writeToFile_(file, value);
375
376 // introduce a line break after a certain time
377 if((++counter) > numBeforeLineBreak)
378 {
379 file << std::endl;
380 counter = 0;
381 }
382 }
383 file << "\n</DataArray>\n";
384 }
385
392 void writeToFile_(std::ostream& file, const Scalar& s)
393 {
394 file << s << " ";
395 }
396
403 void writeToFile_(std::ostream& file, const DimWorldVector& g)
404 {
405 assert(g.size() > 1 && g.size() < 4);
406 if(g.size() < 3)
407 file << g << " 0 ";
408 else
409 file << g << " ";
410 }
411
412 const std::vector<GlobalPosition>& coordinates_;
413 std::list<ScalarFunction> scalarPointData_;
414 std::list<VectorFunction> vectorPointData_;
415};
416} // end namespace Dumux
417
418#endif
Definition: adapt.hh:29
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