version 3.8
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// 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_IO_VTK_VTKREADER_HH
13#define DUMUX_IO_VTK_VTKREADER_HH
14
15#include <iostream>
16#include <iterator>
17#include <algorithm>
18#include <memory>
19#include <type_traits>
20#include <unordered_map>
21#include <utility>
22
23#include <dune/common/parallel/mpihelper.hh>
24#include <dune/common/exceptions.hh>
25#include <dune/grid/common/capabilities.hh>
26#include <dune/grid/io/file/vtk/common.hh>
27#include <dune/grid/common/gridfactory.hh>
28
29#include <dumux/io/container.hh>
30#include <dumux/io/xml/tinyxml2.h>
31
32namespace Dumux {
33
39{
40public:
44 enum class DataType { cellData, pointData };
45
47 using Data = std::unordered_map<std::string, std::vector<double>>;
48
52 explicit VTKReader(const std::string& fileName)
53 {
54 using namespace tinyxml2;
55 fileName_ = Dune::MPIHelper::getCommunication().size() > 1 ?
56 getProcessFileName_(fileName) : fileName;
57
58 const auto eResult = doc_.LoadFile(fileName_.c_str());
59 if (eResult != tinyxml2::XML_SUCCESS)
60 DUNE_THROW(Dune::IOError, "Couldn't open XML file " << fileName_ << ".");
61
62 const XMLElement* pieceNode = getPieceNode_();
63 if (pieceNode == nullptr)
64 DUNE_THROW(Dune::IOError, "Couldn't get 'Piece' node in " << fileName_ << ".");
65 }
66
72 bool hasData(const std::string& name, const DataType& type) const
73 {
74 using namespace tinyxml2;
75
76 const XMLElement* pieceNode = getPieceNode_();
77 const XMLElement* dataNode = getDataNode_(pieceNode, type);
78 if (dataNode == nullptr)
79 return false;
80
81 const XMLElement* dataArray = findDataArray_(dataNode, name);
82 if (dataArray == nullptr)
83 return false;
84
85 return true;
86 }
87
94 template<class Container>
95 Container readData(const std::string& name, const DataType& type) const
96 {
97 using namespace tinyxml2;
98
99 const XMLElement* pieceNode = getPieceNode_();
100 const XMLElement* dataNode = getDataNode_(pieceNode, type);
101 if (dataNode == nullptr)
102 DUNE_THROW(Dune::IOError, "Couldn't get 'PointData' or 'CellData' node in " << fileName_ << ".");
103
104 const XMLElement* dataArray = findDataArray_(dataNode, name);
105 if (dataArray == nullptr)
106 DUNE_THROW(Dune::IOError, "Couldn't find the data array " << name << ".");
107
108 return parseDataArray_<Container>(dataArray);
109 }
110
115 template<class Grid>
116 std::unique_ptr<Grid> readGrid(bool verbose = false) const
117 {
118 static_assert(!Dune::Capabilities::isCartesian<Grid>::v, "Grid reader only supports unstructured grid implementations");
119
120 if (verbose) std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << fileName_ << "." << std::endl;
121
122 // make a grid factory
123 Dune::GridFactory<Grid> factory;
124
125 readGrid_(factory, verbose);
126
127 return std::unique_ptr<Grid>(factory.createGrid());
128 }
129
136 template<class Grid>
137 std::unique_ptr<Grid> readGrid(Dune::GridFactory<Grid>& factory, bool verbose = false) const
138 {
139 static_assert(!Dune::Capabilities::isCartesian<Grid>::v, "Grid reader only supports unstructured grid implementations");
140
141 if (verbose) std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << fileName_ << "." << std::endl;
142
143 readGrid_(factory, verbose);
144
145 return std::unique_ptr<Grid>(factory.createGrid());
146 }
147
156 template<class Grid>
157 std::unique_ptr<Grid> readGrid(Dune::GridFactory<Grid>& factory, Data& cellData, Data& pointData, bool verbose = false) const
158 {
159 static_assert(!Dune::Capabilities::isCartesian<Grid>::v, "Grid reader only supports unstructured grid implementations");
160
161 if (verbose) std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << fileName_ << "." << std::endl;
162
163 readGrid_(factory, verbose);
164 readGridData_(cellData, pointData, verbose);
165
166 return std::unique_ptr<Grid>(factory.createGrid());
167 }
168
169private:
173 std::string getProcessFileName_(const std::string& pvtkFileName)
174 {
175 using namespace tinyxml2;
176
177 XMLDocument pDoc;
178 const auto eResult = pDoc.LoadFile(pvtkFileName.c_str());
179 if (eResult != XML_SUCCESS)
180 DUNE_THROW(Dune::IOError, "Couldn't open XML file " << pvtkFileName << ".");
181
182 // get the first piece node
183 const XMLElement* pieceNode = getPieceNode_(pDoc, pvtkFileName);
184 const auto myrank = Dune::MPIHelper::getCommunication().rank();
185 for (int rank = 0; rank < myrank; ++rank)
186 {
187 pieceNode = pieceNode->NextSiblingElement("Piece");
188 if (pieceNode == nullptr)
189 DUNE_THROW(Dune::IOError, "Couldn't find 'Piece' node for rank "
190 << rank << " in " << pvtkFileName << ".");
191 }
192
193 const char *vtkFileName = pieceNode->Attribute("Source");
194 if (vtkFileName == nullptr)
195 DUNE_THROW(Dune::IOError, "Couldn't get 'Source' attribute of 'Piece' node no. " << myrank << " in " << pvtkFileName);
196
197 return vtkFileName;
198 }
199
205 template<class Grid>
206 void readGrid_(Dune::GridFactory<Grid>& factory, bool verbose = false) const
207 {
208 using namespace tinyxml2;
209
210 const XMLElement* pieceNode = getPieceNode_();
211 const XMLElement* pointsNode = pieceNode->FirstChildElement("Points")->FirstChildElement("DataArray");
212 if (pointsNode == nullptr)
213 DUNE_THROW(Dune::IOError, "Couldn't get data array of points in " << fileName_ << ".");
214
215 using Point3D = Dune::FieldVector<double, 3>;
216 std::vector<Point3D> points3D;
217 std::stringstream dataStream(pointsNode->GetText());
218 std::istream_iterator<Point3D> it(dataStream);
219 std::copy(it, std::istream_iterator<Point3D>(), std::back_inserter(points3D));
220
221 // adapt point dimensions if grid dimension is smaller than 3
222 auto points = adaptPointDimension_<Grid::dimensionworld>(std::move(points3D));
223
224 if (verbose) std::cout << "Found " << points.size() << " vertices." << std::endl;
225
226 // insert vertices to the grid factory
227 for (auto&& point : points)
228 factory.insertVertex(std::move(point));
229
230 const XMLElement* cellsNode = pieceNode->FirstChildElement("Cells");
231 const XMLElement* linesNode = pieceNode->FirstChildElement("Lines");
232 if (cellsNode)
233 {
234 const XMLElement* connectivityNode = findDataArray_(cellsNode, "connectivity");
235 const XMLElement* offsetsNode = findDataArray_(cellsNode, "offsets");
236 const XMLElement* typesNode = findDataArray_(cellsNode, "types");
237
238 const auto connectivity = parseDataArray_<std::vector<unsigned int>>(connectivityNode);
239 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
240 const auto types = parseDataArray_<std::vector<unsigned int>>(typesNode);
241
242 if (verbose) std::cout << "Found " << offsets.size() << " element." << std::endl;
243
244 unsigned int lastOffset = 0;
245 for (unsigned int i = 0; i < offsets.size(); ++i)
246 {
247 const auto geomType = vtkToDuneGeomType_(types[i]);
248 unsigned int offset = offsets[i];
249 std::vector<unsigned int> corners; corners.resize(offset-lastOffset);
250 for (unsigned int j = 0; j < offset-lastOffset; ++j)
251 corners[Dune::VTK::renumber(geomType, j)] = connectivity[lastOffset+j];
252 factory.insertElement(geomType, std::move(corners));
253 lastOffset = offset;
254 }
255 }
256 // for poly data
257 else if (linesNode)
258 {
259 // sanity check
260 if (Grid::dimension != 1)
261 DUNE_THROW(Dune::IOError, "Grid expects dimension " << Grid::dimension
262 << " but " << fileName_ << " contains a 1D grid.");
263
264 const XMLElement* connectivityNode = findDataArray_(linesNode, "connectivity");
265 const XMLElement* offsetsNode = findDataArray_(linesNode, "offsets");
266
267 const auto connectivity = parseDataArray_<std::vector<unsigned int>>(connectivityNode);
268 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
269
270 if (verbose) std::cout << "Found " << offsets.size() << " polylines." << std::endl;
271
272 unsigned int lastOffset = 0;
273 for (unsigned int i = 0; i < offsets.size(); ++i)
274 {
275 // a polyline can have many points in the VTK format
276 // split the line in segments with two points
277 unsigned int offset = offsets[i];
278 for (unsigned int j = 0; j < offset-lastOffset-1; ++j)
279 factory.insertElement(Dune::GeometryTypes::line,
280 std::vector<unsigned int>({ connectivity[lastOffset+j], connectivity[lastOffset+j+1] }));
281 lastOffset = offset;
282 }
283 }
284 else
285 DUNE_THROW(Dune::IOError, "No Cells or Lines element found in " << fileName_);
286 }
287
294 void readGridData_(Data& cellData, Data& pointData, bool verbose = false) const
295 {
296 using namespace tinyxml2;
297
298 const XMLElement* pieceNode = getPieceNode_();
299 const XMLElement* cellDataNode = getDataNode_(pieceNode, DataType::cellData);
300 if (cellDataNode != nullptr)
301 {
302 const XMLElement* cellsNode = pieceNode->FirstChildElement("Cells");
303 const XMLElement* linesNode = pieceNode->FirstChildElement("Lines");
304
305 if (cellsNode)
306 {
307 const XMLElement* dataArray = cellDataNode->FirstChildElement("DataArray");
308 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
309 {
310 const char *attributeText = dataArray->Attribute("Name");
311
312 if (attributeText == nullptr)
313 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a cell data array.");
314
315 cellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
316 }
317 }
318 // for poly data
319 else if (linesNode)
320 {
321 // first parse all the cell data (each cell in this sense can be a polyline)
322 const XMLElement* dataArray = cellDataNode->FirstChildElement("DataArray");
323 if (dataArray)
324 {
325 Data polyLineCellData;
326 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
327 {
328 const char *attributeText = dataArray->Attribute("Name");
329
330 if (attributeText == nullptr)
331 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a cell data array.");
332
333 polyLineCellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
334 }
335
336 // a polyline can have many points in the VTK format
337 // we split the line in segments with two points
338 // so we also need to duplicate the cell data to fit the increased line number
339 const XMLElement* offsetsNode = findDataArray_(linesNode, "offsets");
340 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
341
342 if (offsets.size() != polyLineCellData.begin()->second.size())
343 DUNE_THROW(Dune::IOError, "Expected the same number of cell data entries (is "
344 << polyLineCellData.begin()->second.size()
345 << ") as polylines (" << offsets.size() << ")!");
346
347 // count the number of Dune cells to be able to resize the data vectors
348 unsigned int lastOffset = 0;
349 std::size_t numCells = 0;
350 for (unsigned int i = 0; i < offsets.size(); ++i)
351 {
352 unsigned int offset = offsets[i];
353 for (unsigned int j = 0; j < offset-lastOffset-1; ++j)
354 ++numCells;
355 lastOffset = offset;
356 }
357
358 // create the data arrays
359 for (const auto& dArray : polyLineCellData)
360 {
361 cellData[dArray.first] = std::vector<double>(numCells);
362 auto& cd = cellData[dArray.first];
363 const auto& pd = dArray.second;
364
365 lastOffset = 0;
366 std::size_t cellIdx = 0;
367 for (unsigned int i = 0; i < offsets.size(); ++i)
368 {
369 unsigned int offset = offsets[i];
370 for (unsigned int j = 0; j < offset-lastOffset-1; ++j)
371 cd[cellIdx++] = pd[i];
372 lastOffset = offset;
373 }
374 }
375 }
376 }
377 else
378 DUNE_THROW(Dune::IOError, "No Cells or Lines element found in " << fileName_);
379 }
380
381 const XMLElement* pointDataNode = getDataNode_(pieceNode, DataType::pointData);
382 if (pointDataNode != nullptr)
383 {
384 const XMLElement* dataArray = pointDataNode->FirstChildElement("DataArray");
385 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
386 {
387 const char *attributeText = dataArray->Attribute("Name");
388
389 if (attributeText == nullptr)
390 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a point data array.");
391
392 pointData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
393 }
394 }
395 }
396
401 const tinyxml2::XMLElement* getPieceNode_() const
402 { return getPieceNode_(doc_, fileName_); }
403
410 const tinyxml2::XMLElement* getPieceNode_(const tinyxml2::XMLDocument& doc, const std::string& fileName) const
411 {
412 using namespace tinyxml2;
413
414 const XMLElement* pieceNode = doc.FirstChildElement("VTKFile");
415 if (pieceNode == nullptr)
416 DUNE_THROW(Dune::IOError, "Couldn't get 'VTKFile' node in " << fileName << ".");
417
418 pieceNode = pieceNode->FirstChildElement("UnstructuredGrid");
419 if (pieceNode == nullptr)
420 pieceNode = doc.FirstChildElement("VTKFile")->FirstChildElement("PolyData");
421 if (pieceNode == nullptr)
422 pieceNode = doc.FirstChildElement("VTKFile")->FirstChildElement("PUnstructuredGrid");
423 if (pieceNode == nullptr)
424 pieceNode = doc.FirstChildElement("VTKFile")->FirstChildElement("PPolyData");
425 if (pieceNode == nullptr)
426 DUNE_THROW(Dune::IOError, "Couldn't get 'UnstructuredGrid', 'PUnstructuredGrid', 'PolyData', or 'PPolyData' node in " << fileName << ".");
427
428 return pieceNode->FirstChildElement("Piece");
429 }
430
437 const tinyxml2::XMLElement* getDataNode_(const tinyxml2::XMLElement* pieceNode, const DataType& type) const
438 {
439 using namespace tinyxml2;
440
441 const XMLElement* dataNode = nullptr;
442 if (type == DataType::pointData)
443 dataNode = pieceNode->FirstChildElement("PointData");
444 else if (type == DataType::cellData)
445 dataNode = pieceNode->FirstChildElement("CellData");
446 else
447 DUNE_THROW(Dune::IOError, "Only cell and point data are supported.");
448
449 return dataNode;
450 }
451
458 const tinyxml2::XMLElement* findDataArray_(const tinyxml2::XMLElement* dataNode, const std::string& name) const
459 {
460 using namespace tinyxml2;
461
462 // loop over XML node siblings to find the correct data array
463 const XMLElement* dataArray = dataNode->FirstChildElement("DataArray");
464 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
465 {
466 const char *attributeText = dataArray->Attribute("Name");
467
468 if (attributeText == nullptr)
469 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a data array.");
470
471 if (attributeText == name)
472 break;
473 }
474
475 return dataArray;
476 }
477
483 template<class Container>
484 Container parseDataArray_(const tinyxml2::XMLElement* dataArray) const
485 {
486 std::stringstream dataStream(dataArray->GetText());
487 return readStreamToContainer<Container>(dataStream);
488 }
489
494 Dune::GeometryType vtkToDuneGeomType_(unsigned int vtkCellType) const
495 {
496 switch (vtkCellType)
497 {
498 case Dune::VTK::GeometryType::vertex: return Dune::GeometryTypes::vertex;
500 case Dune::VTK::GeometryType::triangle: return Dune::GeometryTypes::triangle;
501 case Dune::VTK::GeometryType::quadrilateral: return Dune::GeometryTypes::quadrilateral;
502 case Dune::VTK::GeometryType::hexahedron: return Dune::GeometryTypes::hexahedron;
503 case Dune::VTK::GeometryType::prism: return Dune::GeometryTypes::prism;
504 case Dune::VTK::GeometryType::pyramid: return Dune::GeometryTypes::pyramid;
505 default: DUNE_THROW(Dune::NotImplemented, "VTK cell type " << vtkCellType);
506 }
507 }
508
509 template<int dim, std::enable_if_t<(dim < 3), int> = 0>
510 std::vector<Dune::FieldVector<double, dim>>
511 adaptPointDimension_(std::vector<Dune::FieldVector<double, 3>>&& points3D) const
512 {
513 std::vector<Dune::FieldVector<double, dim>> points(points3D.size());
514 for (std::size_t i = 0; i < points.size(); ++i)
515 for (int j = 0; j < dim; ++j)
516 points[i][j] = points3D[i][j];
517 return points;
518 }
519
520 template<int dim, std::enable_if_t<(dim == 3), int> = 0>
521 std::vector<Dune::FieldVector<double, dim>>
522 adaptPointDimension_(std::vector<Dune::FieldVector<double, 3>>&& points3D) const
523 { return std::move(points3D); }
524
525 std::string fileName_;
526 tinyxml2::XMLDocument doc_;
527};
528
529} // end namespace Dumux
530
531#endif
A vtk file reader using tinyxml2 as xml backend.
Definition: vtkreader.hh:39
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:157
bool hasData(const std::string &name, const DataType &type) const
Reviews data from the vtk file to check if there is a data array with a specified name.
Definition: vtkreader.hh:72
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:95
VTKReader(const std::string &fileName)
The constructor creates a tinyxml2::XMLDocument from file.
Definition: vtkreader.hh:52
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:137
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:116
DataType
The data array types.
Definition: vtkreader.hh:44
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:47
Free functions to write and read a sequence container to and from a file.
constexpr Line line
Definition: couplingmanager1d3d_line.hh:31
Definition: adapt.hh:17