3.1-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 <dune/common/fvector.hh>
31#include <dune/common/exceptions.hh>
32#include <dune/common/path.hh>
33#include <dune/grid/io/file/vtk/common.hh>
34
35namespace Dumux {
36
46template<class Scalar, class GlobalPosition>
48{
49 // GlobalPosition is used for the point coordinates, DimWorldVector for the actual data.
50 // GlobalPosition's ctype does not necessarily equal Scalar.
51 using DimWorldVector = Dune::FieldVector<Scalar, GlobalPosition::size()>;
52
53 static constexpr unsigned int precision = 6;
54 static constexpr unsigned int numBeforeLineBreak = 15;
55
59 template<class ContainerType>
60 class VTKFunction
61 {
62 public:
70 VTKFunction(const ContainerType& data, const std::string& name, const int numComponents) : data_(data), name_(name), numComponents_(numComponents)
71 {}
72
76 const std::string& name() const
77 {
78 return name_;
79 }
80
84 int numComponents() const
85 {
86 return numComponents_;
87 }
88
94 auto& operator() (const int idx) const { return data_[idx]; }
95
96 decltype(auto) begin() const
97 {
98 return data_.begin();
99 }
100
101 decltype(auto) end() const
102 {
103 return data_.end();
104 }
105
109 int size() const
110 {
111 return data_.size();
112 }
113
114 private:
115 const ContainerType& data_;
116 const std::string name_;
117 const int numComponents_;
118 };
119
120
121public:
122 using ScalarFunction = VTKFunction<std::vector<Scalar>>;
123 using VectorFunction = VTKFunction<std::vector<DimWorldVector>>;
124
125
126 PointCloudVtkWriter(const std::vector<GlobalPosition>& coordinates) : coordinates_(coordinates)
127 {}
128
135 void write(const std::string& name, Dune::VTK::OutputType type = Dune::VTK::ascii)
136 {
137 auto filename = getSerialPieceName(name, "");
138
139 std::ofstream file;
140 file.open(filename);
141 writeHeader_(file);
142 writeCoordinates_(file, coordinates_);
143 writeDataInfo_(file);
144
145 for (auto&& data : scalarPointData_)
146 writeData_(file, data);
147
148 for (auto&& data :vectorPointData_)
149 writeData_(file, data);
150
151 if (!scalarPointData_.empty() || !vectorPointData_.empty())
152 file << "</PointData>\n";
153
154 file << "</Piece>\n";
155 file << "</PolyData>\n";
156 file << "</VTKFile>";
157
158 clear();
159 file.close();
160 }
161
177 void pwrite(const std::string & name, const std::string & path, const std::string & extendpath,
178 Dune::VTK::OutputType type = Dune::VTK::ascii)
179 {
180 DUNE_THROW(Dune::NotImplemented, "parallel point cloud vtk output not supported yet");
181 }
182
189 void addPointData(const std::vector<Scalar>& v, const std::string &name)
190 {
191 assert(v.size() == coordinates_.size());
192 scalarPointData_.push_back(ScalarFunction(v, name, 1));
193 }
194
201 void addPointData(const std::vector<DimWorldVector>& v, const std::string &name)
202 {
203 assert(v.size() == coordinates_.size());
204 vectorPointData_.push_back(VectorFunction(v, name, 3));
205 }
206
210 void clear()
211 {
212 scalarPointData_.clear();
213 vectorPointData_.clear();
214 }
215
227 std::string getSerialPieceName(const std::string& name,
228 const std::string& path) const
229 {
230 static const std::string extension = ".vtp";
231
232 return Dune::concatPaths(path, name+extension);
233 }
234
247 std::string getParallelHeaderName(const std::string& name,
248 const std::string& path,
249 int commSize) const
250 {
251 std::ostringstream s;
252 if(path.size() > 0) {
253 s << path;
254 if(path[path.size()-1] != '/')
255 s << '/';
256 }
257 s << 's' << std::setw(4) << std::setfill('0') << commSize << '-';
258 s << name;
259 s << ".pvtp";
260 return s.str();
261 }
262
263
264private:
268 void writeHeader_(std::ostream& file)
269 {
270 std::string header = "<?xml version=\"1.0\"?>\n";
271 header += "<VTKFile type=\"PolyData\" version=\"0.1\" byte_order=\"LittleEndian\">\n";
272 header += "<PolyData>\n";
273 header += "<Piece NumberOfLines=\"0\" NumberOfPoints=\"" + std::to_string(coordinates_.size()) + "\">\n";
274 file << header;
275 }
276
280 void writeDataInfo_(std::ostream& file)
281 {
282 std::string scalarName;
283 std::string vectorName;
284 bool foundScalar = false;
285 bool foundVector = false;
286
287 for(auto&& data : scalarPointData_)
288 {
289 if(data.numComponents() == 1 && !foundScalar)
290 {
291 scalarName = data.name();
292 foundScalar = true;
293 continue;
294 }
295
296 if(data.numComponents() > 1 && !foundVector)
297 {
298 vectorName = data.name();
299 foundVector = true;
300 }
301 }
302
303 for(auto&& data : vectorPointData_)
304 {
305 if(data.numComponents() > 1 && !foundVector)
306 {
307 vectorName = data.name();
308 foundVector = true;
309 }
310 }
311
312 if(foundScalar)
313 if(foundVector)
314 file << "<PointData Scalars=\"" << scalarName << "\" Vectors=\"" << vectorName <<"\">\n";
315 else
316 file << "<PointData Scalars=\"" << scalarName << "\">\n";
317 else if(foundVector)
318 file << "<PointData Vectors=\"" << vectorName << "\">\n";
319 else
320 return;
321 }
322
329 void writeCoordinates_(std::ostream& file, const std::vector<GlobalPosition>& positions)
330 {
331 // write the positions to the file
332 file << "<Points>\n";
333 file << "<DataArray type=\"Float32\" Name=\"Coordinates\" NumberOfComponents=\"3\" format=\"ascii\">\n";
334 int counter = 0;
335 for(auto&& x : positions)
336 {
337 file << x ;
338
339 if(x.size() == 1)
340 file << " 0 0 ";
341 if(x.size() == 2)
342 file << " 0 ";
343 if(x.size() == 3)
344 file << " ";
345
346 // introduce a line break after a certain time
347 if((++counter) > numBeforeLineBreak)
348 {
349 file << std::endl;
350 counter = 0;
351 }
352 }
353 file << "\n</DataArray>\n";
354 file << "</Points>\n";
355 }
356
363 template<class T>
364 void writeData_(std::ostream& file, const T& data)
365 {
366 file << "<DataArray type=\"Float32\" Name=\"" << data.name() << "\" NumberOfComponents=\"" << data.numComponents() << "\" format=\"ascii\">\n";
367 int counter = 0;
368 for(auto&& value : data)
369 {
370 // forward to specialized function
371 writeToFile_(file, value);
372
373 // introduce a line break after a certain time
374 if((++counter) > numBeforeLineBreak)
375 {
376 file << std::endl;
377 counter = 0;
378 }
379 }
380 file << "\n</DataArray>\n";
381 }
382
389 void writeToFile_(std::ostream& file, const Scalar& s)
390 {
391 file << s << " ";
392 }
393
400 void writeToFile_(std::ostream& file, const DimWorldVector& g)
401 {
402 assert(g.size() > 1 && g.size() < 4);
403 if(g.size() < 3)
404 file << g << " 0 ";
405 else
406 file << g << " ";
407 }
408
409 const std::vector<GlobalPosition>& coordinates_;
410 std::list<ScalarFunction> scalarPointData_;
411 std::list<VectorFunction> vectorPointData_;
412};
413} // end namespace Dumux
414
415#endif
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: pointcloudvtkwriter.hh:48
std::string getSerialPieceName(const std::string &name, const std::string &path) const
Return name of a serial header file.
Definition: pointcloudvtkwriter.hh:227
VTKFunction< std::vector< Scalar > > ScalarFunction
Definition: pointcloudvtkwriter.hh:122
std::string getParallelHeaderName(const std::string &name, const std::string &path, int commSize) const
Return name of a parallel header file.
Definition: pointcloudvtkwriter.hh:247
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:201
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:189
VTKFunction< std::vector< DimWorldVector > > VectorFunction
Definition: pointcloudvtkwriter.hh:123
PointCloudVtkWriter(const std::vector< GlobalPosition > &coordinates)
Definition: pointcloudvtkwriter.hh:126
void write(const std::string &name, Dune::VTK::OutputType type=Dune::VTK::ascii)
Create an output file.
Definition: pointcloudvtkwriter.hh:135
void clear()
Clears the data lists.
Definition: pointcloudvtkwriter.hh:210
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:177