3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
vtkreader.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_IO_VTK_VTKREADER_HH
25#define DUMUX_IO_VTK_VTKREADER_HH
26
27#include <iostream>
28#include <iterator>
29#include <algorithm>
30#include <memory>
31#include <type_traits>
32#include <unordered_map>
33#include <utility>
34
35#include <dune/common/parallel/mpihelper.hh>
36#include <dune/common/exceptions.hh>
37#include <dune/grid/common/capabilities.hh>
38#include <dune/grid/io/file/vtk/common.hh>
39#include <dumux/io/xml/tinyxml2.h>
40#include <dune/grid/common/gridfactory.hh>
41
42namespace Dumux {
43
49{
50public:
54 enum class DataType { cellData, pointData };
55
57 using Data = std::unordered_map<std::string, std::vector<double>>;
58
62 explicit VTKReader(const std::string& fileName)
63 {
64 fileName_ = Dune::MPIHelper::getCollectiveCommunication().size() > 1 ?
65 getProcessFileName_(fileName) : fileName;
66
67 const auto eResult = doc_.LoadFile(fileName_.c_str());
68 if (eResult != tinyxml2::XML_SUCCESS)
69 DUNE_THROW(Dune::IOError, "Couldn't open XML file " << fileName_ << ".");
70 }
71
78 template<class Container>
79 Container readData(const std::string& name, const DataType& type) const
80 {
81 using namespace tinyxml2;
82
83 const XMLElement* pieceNode = getPieceNode_();
84 if (pieceNode == nullptr)
85 DUNE_THROW(Dune::IOError, "Couldn't get 'Piece' node in " << fileName_ << ".");
86
87 const XMLElement* dataNode = getDataNode_(pieceNode, type);
88 if (dataNode == nullptr)
89 DUNE_THROW(Dune::IOError, "Couldn't get 'PointData' or 'CellData' node in " << fileName_ << ".");
90
91 const XMLElement* dataArray = findDataArray_(dataNode, name);
92 if (dataArray == nullptr)
93 DUNE_THROW(Dune::IOError, "Couldn't find the data array " << name << ".");
94
95 return parseDataArray_<Container>(dataArray);
96 }
97
102 template<class Grid>
103 std::unique_ptr<Grid> readGrid(bool verbose = false) const
104 {
105 static_assert(!Dune::Capabilities::isCartesian<Grid>::v, "Grid reader only supports unstructured grid implementations");
106
107 if (verbose) std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << fileName_ << "." << std::endl;
108
109 // make a grid factory
110 Dune::GridFactory<Grid> factory;
111
112 readGrid_(factory, verbose);
113
114 return std::unique_ptr<Grid>(factory.createGrid());
115 }
116
123 template<class Grid>
124 std::unique_ptr<Grid> readGrid(Dune::GridFactory<Grid>& factory, bool verbose = false) const
125 {
126 static_assert(!Dune::Capabilities::isCartesian<Grid>::v, "Grid reader only supports unstructured grid implementations");
127
128 if (verbose) std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << fileName_ << "." << std::endl;
129
130 readGrid_(factory, verbose);
131
132 return std::unique_ptr<Grid>(factory.createGrid());
133 }
134
143 template<class Grid>
144 std::unique_ptr<Grid> readGrid(Dune::GridFactory<Grid>& factory, Data& cellData, Data& pointData, bool verbose = false) const
145 {
146 static_assert(!Dune::Capabilities::isCartesian<Grid>::v, "Grid reader only supports unstructured grid implementations");
147
148 if (verbose) std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << fileName_ << "." << std::endl;
149
150 readGrid_(factory, verbose);
151 readGridData_(cellData, pointData, verbose);
152
153 return std::unique_ptr<Grid>(factory.createGrid());
154 }
155
156private:
160 std::string getProcessFileName_(const std::string& pvtkFileName)
161 {
162 using namespace tinyxml2;
163
164 XMLDocument pDoc;
165 const auto eResult = pDoc.LoadFile(pvtkFileName.c_str());
166 if (eResult != XML_SUCCESS)
167 DUNE_THROW(Dune::IOError, "Couldn't open XML file " << pvtkFileName << ".");
168
169 // get the first piece node
170 const XMLElement* pieceNode = getPieceNode_(pDoc, pvtkFileName);
171 if (pieceNode == nullptr)
172 DUNE_THROW(Dune::IOError, "Couldn't get 'Piece' node in " << pvtkFileName << ".");
173
174 const auto myrank = Dune::MPIHelper::getCollectiveCommunication().rank();
175 for (int rank = 0; rank < myrank; ++rank)
176 {
177 pieceNode = pieceNode->NextSiblingElement("Piece");
178 if (pieceNode == nullptr)
179 DUNE_THROW(Dune::IOError, "Couldn't find 'Piece' node for rank "
180 << rank << " in " << pvtkFileName << ".");
181 }
182
183 const char *vtkFileName = pieceNode->Attribute("Source");
184 if (vtkFileName == nullptr)
185 DUNE_THROW(Dune::IOError, "Couldn't get 'Source' attribute of 'Piece' node no. " << myrank << " in " << pvtkFileName);
186
187 return vtkFileName;
188 }
189
195 template<class Grid>
196 void readGrid_(Dune::GridFactory<Grid>& factory, bool verbose = false) const
197 {
198 using namespace tinyxml2;
199
200 const XMLElement* pieceNode = getPieceNode_();
201 if (pieceNode == nullptr)
202 DUNE_THROW(Dune::IOError, "Couldn't get 'Piece' node in " << fileName_ << ".");
203
204 const XMLElement* pointsNode = pieceNode->FirstChildElement("Points")->FirstChildElement("DataArray");
205 if (pointsNode == nullptr)
206 DUNE_THROW(Dune::IOError, "Couldn't get data array of points in " << fileName_ << ".");
207
208 using Point3D = Dune::FieldVector<double, 3>;
209 std::vector<Point3D> points3D;
210 std::stringstream dataStream(pointsNode->GetText());
211 std::istream_iterator<Point3D> it(dataStream);
212 std::copy(it, std::istream_iterator<Point3D>(), std::back_inserter(points3D));
213
214 // adapt point dimensions if grid dimension is smaller than 3
215 auto points = adaptPointDimension_<Grid::dimensionworld>(std::move(points3D));
216
217 if (verbose) std::cout << "Found " << points.size() << " vertices." << std::endl;
218
219 // insert vertices to the grid factory
220 for (auto&& point : points)
221 factory.insertVertex(std::move(point));
222
223 const XMLElement* cellsNode = pieceNode->FirstChildElement("Cells");
224 const XMLElement* linesNode = pieceNode->FirstChildElement("Lines");
225 if (cellsNode)
226 {
227 const XMLElement* connectivityNode = findDataArray_(cellsNode, "connectivity");
228 const XMLElement* offsetsNode = findDataArray_(cellsNode, "offsets");
229 const XMLElement* typesNode = findDataArray_(cellsNode, "types");
230
231 const auto connectivity = parseDataArray_<std::vector<unsigned int>>(connectivityNode);
232 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
233 const auto types = parseDataArray_<std::vector<unsigned int>>(typesNode);
234
235 if (verbose) std::cout << "Found " << offsets.size() << " element." << std::endl;
236
237 unsigned int lastOffset = 0;
238 for (unsigned int i = 0; i < offsets.size(); ++i)
239 {
240 const auto geomType = vtkToDuneGeomType_(types[i]);
241 unsigned int offset = offsets[i];
242 std::vector<unsigned int> corners; corners.resize(offset-lastOffset);
243 for (unsigned int j = 0; j < offset-lastOffset; ++j)
244 corners[Dune::VTK::renumber(geomType, j)] = connectivity[lastOffset+j];
245 factory.insertElement(geomType, std::move(corners));
246 lastOffset = offset;
247 }
248 }
249 // for poly data
250 else if (linesNode)
251 {
252 // sanity check
253 if (Grid::dimension != 1)
254 DUNE_THROW(Dune::IOError, "Grid expects dimension " << Grid::dimension
255 << " but " << fileName_ << " contains a 1D grid.");
256
257 const XMLElement* connectivityNode = findDataArray_(linesNode, "connectivity");
258 const XMLElement* offsetsNode = findDataArray_(linesNode, "offsets");
259
260 const auto connectivity = parseDataArray_<std::vector<unsigned int>>(connectivityNode);
261 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
262
263 if (verbose) std::cout << "Found " << offsets.size() << " polylines." << std::endl;
264
265 unsigned int lastOffset = 0;
266 for (unsigned int i = 0; i < offsets.size(); ++i)
267 {
268 // a polyline can have many points in the VTK format
269 // split the line in segments with two points
270 unsigned int offset = offsets[i];
271 for (unsigned int j = 0; j < offset-lastOffset-1; ++j)
272 factory.insertElement(Dune::GeometryTypes::line,
273 std::vector<unsigned int>({ connectivity[lastOffset+j], connectivity[lastOffset+j+1] }));
274 lastOffset = offset;
275 }
276 }
277 else
278 DUNE_THROW(Dune::IOError, "No Cells or Lines element found in " << fileName_);
279 }
280
287 void readGridData_(Data& cellData, Data& pointData, bool verbose = false) const
288 {
289 using namespace tinyxml2;
290
291 const XMLElement* pieceNode = getPieceNode_();
292 if (pieceNode == nullptr)
293 DUNE_THROW(Dune::IOError, "Couldn't get 'Piece' node in " << fileName_ << ".");
294
295 const XMLElement* cellDataNode = getDataNode_(pieceNode, DataType::cellData);
296 if (cellDataNode != nullptr)
297 {
298 const XMLElement* cellsNode = pieceNode->FirstChildElement("Cells");
299 const XMLElement* linesNode = pieceNode->FirstChildElement("Lines");
300
301 if (cellsNode)
302 {
303 const XMLElement* dataArray = cellDataNode->FirstChildElement("DataArray");
304 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
305 {
306 const char *attributeText = dataArray->Attribute("Name");
307
308 if (attributeText == nullptr)
309 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a cell data array.");
310
311 cellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
312 }
313 }
314 // for poly data
315 else if (linesNode)
316 {
317 // first parse all the cell data (each cell in this sense can be a polyline)
318 const XMLElement* dataArray = cellDataNode->FirstChildElement("DataArray");
319 if (dataArray)
320 {
321 Data polyLineCellData;
322 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
323 {
324 const char *attributeText = dataArray->Attribute("Name");
325
326 if (attributeText == nullptr)
327 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a cell data array.");
328
329 polyLineCellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
330 }
331
332 // a polyline can have many points in the VTK format
333 // we split the line in segments with two points
334 // so we also need to duplicate the cell data to fit the increased line number
335 const XMLElement* offsetsNode = findDataArray_(linesNode, "offsets");
336 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
337
338 if (offsets.size() != polyLineCellData.begin()->second.size())
339 DUNE_THROW(Dune::IOError, "Expected the same number of cell data entries (is "
340 << polyLineCellData.begin()->second.size()
341 << ") as polylines (" << offsets.size() << ")!");
342
343 // count the number of Dune cells to be able to resize the data vectors
344 unsigned int lastOffset = 0;
345 std::size_t numCells = 0;
346 for (unsigned int i = 0; i < offsets.size(); ++i)
347 {
348 unsigned int offset = offsets[i];
349 for (unsigned int j = 0; j < offset-lastOffset-1; ++j)
350 ++numCells;
351 lastOffset = offset;
352 }
353
354 // create the data arrays
355 for (const auto& dArray : polyLineCellData)
356 {
357 cellData[dArray.first] = std::vector<double>(numCells);
358 auto& cd = cellData[dArray.first];
359 const auto& pd = dArray.second;
360
361 lastOffset = 0;
362 std::size_t cellIdx = 0;
363 for (unsigned int i = 0; i < offsets.size(); ++i)
364 {
365 unsigned int offset = offsets[i];
366 for (unsigned int j = 0; j < offset-lastOffset-1; ++j)
367 cd[cellIdx++] = pd[i];
368 lastOffset = offset;
369 }
370 }
371 }
372 }
373 else
374 DUNE_THROW(Dune::IOError, "No Cells or Lines element found in " << fileName_);
375 }
376
377 const XMLElement* pointDataNode = getDataNode_(pieceNode, DataType::pointData);
378 if (pointDataNode != nullptr)
379 {
380 const XMLElement* dataArray = pointDataNode->FirstChildElement("DataArray");
381 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
382 {
383 const char *attributeText = dataArray->Attribute("Name");
384
385 if (attributeText == nullptr)
386 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a point data array.");
387
388 pointData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
389 }
390 }
391 }
392
397 const tinyxml2::XMLElement* getPieceNode_() const
398 { return getPieceNode_(doc_, fileName_); }
399
406 const tinyxml2::XMLElement* getPieceNode_(const tinyxml2::XMLDocument& doc, const std::string& fileName) const
407 {
408 using namespace tinyxml2;
409
410 const XMLElement* pieceNode = doc.FirstChildElement("VTKFile");
411 if (pieceNode == nullptr)
412 DUNE_THROW(Dune::IOError, "Couldn't get 'VTKFile' node in " << fileName << ".");
413
414 pieceNode = pieceNode->FirstChildElement("UnstructuredGrid");
415 if (pieceNode == nullptr)
416 pieceNode = doc.FirstChildElement("VTKFile")->FirstChildElement("PolyData");
417 if (pieceNode == nullptr)
418 pieceNode = doc.FirstChildElement("VTKFile")->FirstChildElement("PUnstructuredGrid");
419 if (pieceNode == nullptr)
420 pieceNode = doc.FirstChildElement("VTKFile")->FirstChildElement("PPolyData");
421 if (pieceNode == nullptr)
422 DUNE_THROW(Dune::IOError, "Couldn't get 'UnstructuredGrid', 'PUnstructuredGrid', 'PolyData', or 'PPolyData' node in " << fileName << ".");
423
424 return pieceNode->FirstChildElement("Piece");
425 }
426
433 const tinyxml2::XMLElement* getDataNode_(const tinyxml2::XMLElement* pieceNode, const DataType& type) const
434 {
435 using namespace tinyxml2;
436
437 const XMLElement* dataNode = nullptr;
438 if (type == DataType::pointData)
439 dataNode = pieceNode->FirstChildElement("PointData");
440 else if (type == DataType::cellData)
441 dataNode = pieceNode->FirstChildElement("CellData");
442 else
443 DUNE_THROW(Dune::IOError, "Only cell and point data are supported.");
444
445 return dataNode;
446 }
447
454 const tinyxml2::XMLElement* findDataArray_(const tinyxml2::XMLElement* dataNode, const std::string& name) const
455 {
456 using namespace tinyxml2;
457
458 // loop over XML node siblings to find the correct data array
459 const XMLElement* dataArray = dataNode->FirstChildElement("DataArray");
460 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
461 {
462 const char *attributeText = dataArray->Attribute("Name");
463
464 if (attributeText == nullptr)
465 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a data array.");
466
467 if (attributeText == name)
468 break;
469 }
470
471 return dataArray;
472 }
473
479 template<class Container>
480 Container parseDataArray_(const tinyxml2::XMLElement* dataArray) const
481 {
482 Container data;
483 std::stringstream dataStream(dataArray->GetText());
484 std::istream_iterator<typename Container::value_type> it(dataStream);
485 std::copy(it, std::istream_iterator<typename Container::value_type>(), std::back_inserter(data));
486 return data;
487 }
488
493 Dune::GeometryType vtkToDuneGeomType_(unsigned int vtkCellType) const
494 {
495 switch (vtkCellType)
496 {
497 case Dune::VTK::GeometryType::vertex: return Dune::GeometryTypes::vertex;
498 case Dune::VTK::GeometryType::line: return Dune::GeometryTypes::line;
499 case Dune::VTK::GeometryType::triangle: return Dune::GeometryTypes::triangle;
500 case Dune::VTK::GeometryType::quadrilateral: return Dune::GeometryTypes::quadrilateral;
501 case Dune::VTK::GeometryType::hexahedron: return Dune::GeometryTypes::hexahedron;
502 case Dune::VTK::GeometryType::prism: return Dune::GeometryTypes::prism;
503 case Dune::VTK::GeometryType::pyramid: return Dune::GeometryTypes::pyramid;
504 default: DUNE_THROW(Dune::NotImplemented, "VTK cell type " << vtkCellType);
505 }
506 }
507
508 template<int dim, std::enable_if_t<(dim < 3), int> = 0>
509 std::vector<Dune::FieldVector<double, dim>>
510 adaptPointDimension_(std::vector<Dune::FieldVector<double, 3>>&& points3D) const
511 {
512 std::vector<Dune::FieldVector<double, dim>> points(points3D.size());
513 for (std::size_t i = 0; i < points.size(); ++i)
514 for (int j = 0; j < dim; ++j)
515 points[i][j] = points3D[i][j];
516 return points;
517 }
518
519 template<int dim, std::enable_if_t<(dim == 3), int> = 0>
520 std::vector<Dune::FieldVector<double, dim>>
521 adaptPointDimension_(std::vector<Dune::FieldVector<double, 3>>&& points3D) const
522 { return std::move(points3D); }
523
524 std::string fileName_;
525 tinyxml2::XMLDocument doc_;
526};
527
528} // end namespace Dumux
529
530#endif
Definition: adapt.hh:29
A vtk file reader using tinyxml2 as xml backend.
Definition: vtkreader.hh:49
std::unique_ptr< Grid > readGrid(Dune::GridFactory< Grid > &factory, Data &cellData, Data &pointData, bool verbose=false) const
Read a grid from a vtk/vtu/vtp file, reading all cell and point data.
Definition: vtkreader.hh:144
Container readData(const std::string &name, const DataType &type) const
read data from the vtk file to a container, e.g. std::vector<double>
Definition: vtkreader.hh:79
VTKReader(const std::string &fileName)
The contructor creates a tinyxml2::XMLDocument from file.
Definition: vtkreader.hh:62
std::unique_ptr< Grid > readGrid(Dune::GridFactory< Grid > &factory, bool verbose=false) const
Read a grid from a vtk/vtu/vtp file, ignoring cell and point data.
Definition: vtkreader.hh:124
std::unique_ptr< Grid > readGrid(bool verbose=false) const
Read a grid from a vtk/vtu/vtp file, ignoring cell and point data.
Definition: vtkreader.hh:103
DataType
The data array types.
Definition: vtkreader.hh:54
std::unordered_map< std::string, std::vector< double > > Data
the cell / point data type for point data read from a grid file
Definition: vtkreader.hh:57