3.3.0
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 <dune/grid/common/gridfactory.hh>
40
41#include <dumux/io/container.hh>
42#include <dumux/io/xml/tinyxml2.h>
43
44namespace Dumux {
45
51{
52public:
56 enum class DataType { cellData, pointData };
57
59 using Data = std::unordered_map<std::string, std::vector<double>>;
60
64 explicit VTKReader(const std::string& fileName)
65 {
66 using namespace tinyxml2;
67 fileName_ = Dune::MPIHelper::getCollectiveCommunication().size() > 1 ?
68 getProcessFileName_(fileName) : fileName;
69
70 const auto eResult = doc_.LoadFile(fileName_.c_str());
71 if (eResult != tinyxml2::XML_SUCCESS)
72 DUNE_THROW(Dune::IOError, "Couldn't open XML file " << fileName_ << ".");
73
74 const XMLElement* pieceNode = getPieceNode_();
75 if (pieceNode == nullptr)
76 DUNE_THROW(Dune::IOError, "Couldn't get 'Piece' node in " << fileName_ << ".");
77 }
78
84 bool hasData(const std::string& name, const DataType& type) const
85 {
86 using namespace tinyxml2;
87
88 const XMLElement* pieceNode = getPieceNode_();
89 const XMLElement* dataNode = getDataNode_(pieceNode, type);
90 if (dataNode == nullptr)
91 return false;
92
93 const XMLElement* dataArray = findDataArray_(dataNode, name);
94 if (dataArray == nullptr)
95 return false;
96
97 return true;
98 }
99
106 template<class Container>
107 Container readData(const std::string& name, const DataType& type) const
108 {
109 using namespace tinyxml2;
110
111 const XMLElement* pieceNode = getPieceNode_();
112 const XMLElement* dataNode = getDataNode_(pieceNode, type);
113 if (dataNode == nullptr)
114 DUNE_THROW(Dune::IOError, "Couldn't get 'PointData' or 'CellData' node in " << fileName_ << ".");
115
116 const XMLElement* dataArray = findDataArray_(dataNode, name);
117 if (dataArray == nullptr)
118 DUNE_THROW(Dune::IOError, "Couldn't find the data array " << name << ".");
119
120 return parseDataArray_<Container>(dataArray);
121 }
122
127 template<class Grid>
128 std::unique_ptr<Grid> readGrid(bool verbose = false) const
129 {
130 static_assert(!Dune::Capabilities::isCartesian<Grid>::v, "Grid reader only supports unstructured grid implementations");
131
132 if (verbose) std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << fileName_ << "." << std::endl;
133
134 // make a grid factory
135 Dune::GridFactory<Grid> factory;
136
137 readGrid_(factory, verbose);
138
139 return std::unique_ptr<Grid>(factory.createGrid());
140 }
141
148 template<class Grid>
149 std::unique_ptr<Grid> readGrid(Dune::GridFactory<Grid>& factory, bool verbose = false) const
150 {
151 static_assert(!Dune::Capabilities::isCartesian<Grid>::v, "Grid reader only supports unstructured grid implementations");
152
153 if (verbose) std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << fileName_ << "." << std::endl;
154
155 readGrid_(factory, verbose);
156
157 return std::unique_ptr<Grid>(factory.createGrid());
158 }
159
168 template<class Grid>
169 std::unique_ptr<Grid> readGrid(Dune::GridFactory<Grid>& factory, Data& cellData, Data& pointData, bool verbose = false) const
170 {
171 static_assert(!Dune::Capabilities::isCartesian<Grid>::v, "Grid reader only supports unstructured grid implementations");
172
173 if (verbose) std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << fileName_ << "." << std::endl;
174
175 readGrid_(factory, verbose);
176 readGridData_(cellData, pointData, verbose);
177
178 return std::unique_ptr<Grid>(factory.createGrid());
179 }
180
181private:
185 std::string getProcessFileName_(const std::string& pvtkFileName)
186 {
187 using namespace tinyxml2;
188
189 XMLDocument pDoc;
190 const auto eResult = pDoc.LoadFile(pvtkFileName.c_str());
191 if (eResult != XML_SUCCESS)
192 DUNE_THROW(Dune::IOError, "Couldn't open XML file " << pvtkFileName << ".");
193
194 // get the first piece node
195 const XMLElement* pieceNode = getPieceNode_(pDoc, pvtkFileName);
196 const auto myrank = Dune::MPIHelper::getCollectiveCommunication().rank();
197 for (int rank = 0; rank < myrank; ++rank)
198 {
199 pieceNode = pieceNode->NextSiblingElement("Piece");
200 if (pieceNode == nullptr)
201 DUNE_THROW(Dune::IOError, "Couldn't find 'Piece' node for rank "
202 << rank << " in " << pvtkFileName << ".");
203 }
204
205 const char *vtkFileName = pieceNode->Attribute("Source");
206 if (vtkFileName == nullptr)
207 DUNE_THROW(Dune::IOError, "Couldn't get 'Source' attribute of 'Piece' node no. " << myrank << " in " << pvtkFileName);
208
209 return vtkFileName;
210 }
211
217 template<class Grid>
218 void readGrid_(Dune::GridFactory<Grid>& factory, bool verbose = false) const
219 {
220 using namespace tinyxml2;
221
222 const XMLElement* pieceNode = getPieceNode_();
223 const XMLElement* pointsNode = pieceNode->FirstChildElement("Points")->FirstChildElement("DataArray");
224 if (pointsNode == nullptr)
225 DUNE_THROW(Dune::IOError, "Couldn't get data array of points in " << fileName_ << ".");
226
227 using Point3D = Dune::FieldVector<double, 3>;
228 std::vector<Point3D> points3D;
229 std::stringstream dataStream(pointsNode->GetText());
230 std::istream_iterator<Point3D> it(dataStream);
231 std::copy(it, std::istream_iterator<Point3D>(), std::back_inserter(points3D));
232
233 // adapt point dimensions if grid dimension is smaller than 3
234 auto points = adaptPointDimension_<Grid::dimensionworld>(std::move(points3D));
235
236 if (verbose) std::cout << "Found " << points.size() << " vertices." << std::endl;
237
238 // insert vertices to the grid factory
239 for (auto&& point : points)
240 factory.insertVertex(std::move(point));
241
242 const XMLElement* cellsNode = pieceNode->FirstChildElement("Cells");
243 const XMLElement* linesNode = pieceNode->FirstChildElement("Lines");
244 if (cellsNode)
245 {
246 const XMLElement* connectivityNode = findDataArray_(cellsNode, "connectivity");
247 const XMLElement* offsetsNode = findDataArray_(cellsNode, "offsets");
248 const XMLElement* typesNode = findDataArray_(cellsNode, "types");
249
250 const auto connectivity = parseDataArray_<std::vector<unsigned int>>(connectivityNode);
251 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
252 const auto types = parseDataArray_<std::vector<unsigned int>>(typesNode);
253
254 if (verbose) std::cout << "Found " << offsets.size() << " element." << std::endl;
255
256 unsigned int lastOffset = 0;
257 for (unsigned int i = 0; i < offsets.size(); ++i)
258 {
259 const auto geomType = vtkToDuneGeomType_(types[i]);
260 unsigned int offset = offsets[i];
261 std::vector<unsigned int> corners; corners.resize(offset-lastOffset);
262 for (unsigned int j = 0; j < offset-lastOffset; ++j)
263 corners[Dune::VTK::renumber(geomType, j)] = connectivity[lastOffset+j];
264 factory.insertElement(geomType, std::move(corners));
265 lastOffset = offset;
266 }
267 }
268 // for poly data
269 else if (linesNode)
270 {
271 // sanity check
272 if (Grid::dimension != 1)
273 DUNE_THROW(Dune::IOError, "Grid expects dimension " << Grid::dimension
274 << " but " << fileName_ << " contains a 1D grid.");
275
276 const XMLElement* connectivityNode = findDataArray_(linesNode, "connectivity");
277 const XMLElement* offsetsNode = findDataArray_(linesNode, "offsets");
278
279 const auto connectivity = parseDataArray_<std::vector<unsigned int>>(connectivityNode);
280 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
281
282 if (verbose) std::cout << "Found " << offsets.size() << " polylines." << std::endl;
283
284 unsigned int lastOffset = 0;
285 for (unsigned int i = 0; i < offsets.size(); ++i)
286 {
287 // a polyline can have many points in the VTK format
288 // split the line in segments with two points
289 unsigned int offset = offsets[i];
290 for (unsigned int j = 0; j < offset-lastOffset-1; ++j)
291 factory.insertElement(Dune::GeometryTypes::line,
292 std::vector<unsigned int>({ connectivity[lastOffset+j], connectivity[lastOffset+j+1] }));
293 lastOffset = offset;
294 }
295 }
296 else
297 DUNE_THROW(Dune::IOError, "No Cells or Lines element found in " << fileName_);
298 }
299
306 void readGridData_(Data& cellData, Data& pointData, bool verbose = false) const
307 {
308 using namespace tinyxml2;
309
310 const XMLElement* pieceNode = getPieceNode_();
311 const XMLElement* cellDataNode = getDataNode_(pieceNode, DataType::cellData);
312 if (cellDataNode != nullptr)
313 {
314 const XMLElement* cellsNode = pieceNode->FirstChildElement("Cells");
315 const XMLElement* linesNode = pieceNode->FirstChildElement("Lines");
316
317 if (cellsNode)
318 {
319 const XMLElement* dataArray = cellDataNode->FirstChildElement("DataArray");
320 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
321 {
322 const char *attributeText = dataArray->Attribute("Name");
323
324 if (attributeText == nullptr)
325 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a cell data array.");
326
327 cellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
328 }
329 }
330 // for poly data
331 else if (linesNode)
332 {
333 // first parse all the cell data (each cell in this sense can be a polyline)
334 const XMLElement* dataArray = cellDataNode->FirstChildElement("DataArray");
335 if (dataArray)
336 {
337 Data polyLineCellData;
338 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
339 {
340 const char *attributeText = dataArray->Attribute("Name");
341
342 if (attributeText == nullptr)
343 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a cell data array.");
344
345 polyLineCellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
346 }
347
348 // a polyline can have many points in the VTK format
349 // we split the line in segments with two points
350 // so we also need to duplicate the cell data to fit the increased line number
351 const XMLElement* offsetsNode = findDataArray_(linesNode, "offsets");
352 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
353
354 if (offsets.size() != polyLineCellData.begin()->second.size())
355 DUNE_THROW(Dune::IOError, "Expected the same number of cell data entries (is "
356 << polyLineCellData.begin()->second.size()
357 << ") as polylines (" << offsets.size() << ")!");
358
359 // count the number of Dune cells to be able to resize the data vectors
360 unsigned int lastOffset = 0;
361 std::size_t numCells = 0;
362 for (unsigned int i = 0; i < offsets.size(); ++i)
363 {
364 unsigned int offset = offsets[i];
365 for (unsigned int j = 0; j < offset-lastOffset-1; ++j)
366 ++numCells;
367 lastOffset = offset;
368 }
369
370 // create the data arrays
371 for (const auto& dArray : polyLineCellData)
372 {
373 cellData[dArray.first] = std::vector<double>(numCells);
374 auto& cd = cellData[dArray.first];
375 const auto& pd = dArray.second;
376
377 lastOffset = 0;
378 std::size_t cellIdx = 0;
379 for (unsigned int i = 0; i < offsets.size(); ++i)
380 {
381 unsigned int offset = offsets[i];
382 for (unsigned int j = 0; j < offset-lastOffset-1; ++j)
383 cd[cellIdx++] = pd[i];
384 lastOffset = offset;
385 }
386 }
387 }
388 }
389 else
390 DUNE_THROW(Dune::IOError, "No Cells or Lines element found in " << fileName_);
391 }
392
393 const XMLElement* pointDataNode = getDataNode_(pieceNode, DataType::pointData);
394 if (pointDataNode != nullptr)
395 {
396 const XMLElement* dataArray = pointDataNode->FirstChildElement("DataArray");
397 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
398 {
399 const char *attributeText = dataArray->Attribute("Name");
400
401 if (attributeText == nullptr)
402 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a point data array.");
403
404 pointData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
405 }
406 }
407 }
408
413 const tinyxml2::XMLElement* getPieceNode_() const
414 { return getPieceNode_(doc_, fileName_); }
415
422 const tinyxml2::XMLElement* getPieceNode_(const tinyxml2::XMLDocument& doc, const std::string& fileName) const
423 {
424 using namespace tinyxml2;
425
426 const XMLElement* pieceNode = doc.FirstChildElement("VTKFile");
427 if (pieceNode == nullptr)
428 DUNE_THROW(Dune::IOError, "Couldn't get 'VTKFile' node in " << fileName << ".");
429
430 pieceNode = pieceNode->FirstChildElement("UnstructuredGrid");
431 if (pieceNode == nullptr)
432 pieceNode = doc.FirstChildElement("VTKFile")->FirstChildElement("PolyData");
433 if (pieceNode == nullptr)
434 pieceNode = doc.FirstChildElement("VTKFile")->FirstChildElement("PUnstructuredGrid");
435 if (pieceNode == nullptr)
436 pieceNode = doc.FirstChildElement("VTKFile")->FirstChildElement("PPolyData");
437 if (pieceNode == nullptr)
438 DUNE_THROW(Dune::IOError, "Couldn't get 'UnstructuredGrid', 'PUnstructuredGrid', 'PolyData', or 'PPolyData' node in " << fileName << ".");
439
440 return pieceNode->FirstChildElement("Piece");
441 }
442
449 const tinyxml2::XMLElement* getDataNode_(const tinyxml2::XMLElement* pieceNode, const DataType& type) const
450 {
451 using namespace tinyxml2;
452
453 const XMLElement* dataNode = nullptr;
454 if (type == DataType::pointData)
455 dataNode = pieceNode->FirstChildElement("PointData");
456 else if (type == DataType::cellData)
457 dataNode = pieceNode->FirstChildElement("CellData");
458 else
459 DUNE_THROW(Dune::IOError, "Only cell and point data are supported.");
460
461 return dataNode;
462 }
463
470 const tinyxml2::XMLElement* findDataArray_(const tinyxml2::XMLElement* dataNode, const std::string& name) const
471 {
472 using namespace tinyxml2;
473
474 // loop over XML node siblings to find the correct data array
475 const XMLElement* dataArray = dataNode->FirstChildElement("DataArray");
476 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
477 {
478 const char *attributeText = dataArray->Attribute("Name");
479
480 if (attributeText == nullptr)
481 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a data array.");
482
483 if (attributeText == name)
484 break;
485 }
486
487 return dataArray;
488 }
489
495 template<class Container>
496 Container parseDataArray_(const tinyxml2::XMLElement* dataArray) const
497 {
498 std::stringstream dataStream(dataArray->GetText());
499 return readStreamToContainer<Container>(dataStream);
500 }
501
506 Dune::GeometryType vtkToDuneGeomType_(unsigned int vtkCellType) const
507 {
508 switch (vtkCellType)
509 {
510 case Dune::VTK::GeometryType::vertex: return Dune::GeometryTypes::vertex;
512 case Dune::VTK::GeometryType::triangle: return Dune::GeometryTypes::triangle;
513 case Dune::VTK::GeometryType::quadrilateral: return Dune::GeometryTypes::quadrilateral;
514 case Dune::VTK::GeometryType::hexahedron: return Dune::GeometryTypes::hexahedron;
515 case Dune::VTK::GeometryType::prism: return Dune::GeometryTypes::prism;
516 case Dune::VTK::GeometryType::pyramid: return Dune::GeometryTypes::pyramid;
517 default: DUNE_THROW(Dune::NotImplemented, "VTK cell type " << vtkCellType);
518 }
519 }
520
521 template<int dim, std::enable_if_t<(dim < 3), int> = 0>
522 std::vector<Dune::FieldVector<double, dim>>
523 adaptPointDimension_(std::vector<Dune::FieldVector<double, 3>>&& points3D) const
524 {
525 std::vector<Dune::FieldVector<double, dim>> points(points3D.size());
526 for (std::size_t i = 0; i < points.size(); ++i)
527 for (int j = 0; j < dim; ++j)
528 points[i][j] = points3D[i][j];
529 return points;
530 }
531
532 template<int dim, std::enable_if_t<(dim == 3), int> = 0>
533 std::vector<Dune::FieldVector<double, dim>>
534 adaptPointDimension_(std::vector<Dune::FieldVector<double, 3>>&& points3D) const
535 { return std::move(points3D); }
536
537 std::string fileName_;
538 tinyxml2::XMLDocument doc_;
539};
540
541} // end namespace Dumux
542
543#endif
Free functions to write and read a sequence container to and from a file.
Definition: adapt.hh:29
constexpr Line line
Definition: couplingmanager1d3d_line.hh:52
A vtk file reader using tinyxml2 as xml backend.
Definition: vtkreader.hh:51
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:169
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:84
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:107
VTKReader(const std::string &fileName)
The contructor creates a tinyxml2::XMLDocument from file.
Definition: vtkreader.hh:64
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:149
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:128
DataType
The data array types.
Definition: vtkreader.hh:56
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:59