24#ifndef DUMUX_IO_VTK_VTKREADER_HH
25#define DUMUX_IO_VTK_VTKREADER_HH
32#include <unordered_map>
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>
42#include <dumux/io/xml/tinyxml2.h>
59 using Data = std::unordered_map<std::string, std::vector<double>>;
66 using namespace tinyxml2;
67 fileName_ = Dune::MPIHelper::getCollectiveCommunication().size() > 1 ?
68 getProcessFileName_(fileName) : fileName;
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_ <<
".");
74 const XMLElement* pieceNode = getPieceNode_();
75 if (pieceNode ==
nullptr)
76 DUNE_THROW(Dune::IOError,
"Couldn't get 'Piece' node in " << fileName_ <<
".");
86 using namespace tinyxml2;
88 const XMLElement* pieceNode = getPieceNode_();
89 const XMLElement* dataNode = getDataNode_(pieceNode, type);
90 if (dataNode ==
nullptr)
93 const XMLElement* dataArray = findDataArray_(dataNode, name);
94 if (dataArray ==
nullptr)
106 template<
class Container>
109 using namespace tinyxml2;
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_ <<
".");
116 const XMLElement* dataArray = findDataArray_(dataNode, name);
117 if (dataArray ==
nullptr)
118 DUNE_THROW(Dune::IOError,
"Couldn't find the data array " << name <<
".");
120 return parseDataArray_<Container>(dataArray);
128 std::unique_ptr<Grid>
readGrid(
bool verbose =
false)
const
130 static_assert(!Dune::Capabilities::isCartesian<Grid>::v,
"Grid reader only supports unstructured grid implementations");
132 if (verbose) std::cout <<
"Reading " << Grid::dimension <<
"d grid from vtk file " << fileName_ <<
"." << std::endl;
135 Dune::GridFactory<Grid> factory;
137 readGrid_(factory, verbose);
139 return std::unique_ptr<Grid>(factory.createGrid());
149 std::unique_ptr<Grid>
readGrid(Dune::GridFactory<Grid>& factory,
bool verbose =
false)
const
151 static_assert(!Dune::Capabilities::isCartesian<Grid>::v,
"Grid reader only supports unstructured grid implementations");
153 if (verbose) std::cout <<
"Reading " << Grid::dimension <<
"d grid from vtk file " << fileName_ <<
"." << std::endl;
155 readGrid_(factory, verbose);
157 return std::unique_ptr<Grid>(factory.createGrid());
169 std::unique_ptr<Grid>
readGrid(Dune::GridFactory<Grid>& factory,
Data& cellData,
Data& pointData,
bool verbose =
false)
const
171 static_assert(!Dune::Capabilities::isCartesian<Grid>::v,
"Grid reader only supports unstructured grid implementations");
173 if (verbose) std::cout <<
"Reading " << Grid::dimension <<
"d grid from vtk file " << fileName_ <<
"." << std::endl;
175 readGrid_(factory, verbose);
176 readGridData_(cellData, pointData, verbose);
178 return std::unique_ptr<Grid>(factory.createGrid());
185 std::string getProcessFileName_(
const std::string& pvtkFileName)
187 using namespace tinyxml2;
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 <<
".");
195 const XMLElement* pieceNode = getPieceNode_(pDoc, pvtkFileName);
196 const auto myrank = Dune::MPIHelper::getCollectiveCommunication().rank();
197 for (
int rank = 0; rank < myrank; ++rank)
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 <<
".");
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);
218 void readGrid_(Dune::GridFactory<Grid>& factory,
bool verbose =
false)
const
220 using namespace tinyxml2;
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_ <<
".");
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));
234 auto points = adaptPointDimension_<Grid::dimensionworld>(std::move(points3D));
236 if (verbose) std::cout <<
"Found " << points.size() <<
" vertices." << std::endl;
239 for (
auto&& point : points)
240 factory.insertVertex(std::move(point));
242 const XMLElement* cellsNode = pieceNode->FirstChildElement(
"Cells");
243 const XMLElement* linesNode = pieceNode->FirstChildElement(
"Lines");
246 const XMLElement* connectivityNode = findDataArray_(cellsNode,
"connectivity");
247 const XMLElement* offsetsNode = findDataArray_(cellsNode,
"offsets");
248 const XMLElement* typesNode = findDataArray_(cellsNode,
"types");
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);
254 if (verbose) std::cout <<
"Found " << offsets.size() <<
" element." << std::endl;
256 unsigned int lastOffset = 0;
257 for (
unsigned int i = 0; i < offsets.size(); ++i)
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));
272 if (Grid::dimension != 1)
273 DUNE_THROW(Dune::IOError,
"Grid expects dimension " << Grid::dimension
274 <<
" but " << fileName_ <<
" contains a 1D grid.");
276 const XMLElement* connectivityNode = findDataArray_(linesNode,
"connectivity");
277 const XMLElement* offsetsNode = findDataArray_(linesNode,
"offsets");
279 const auto connectivity = parseDataArray_<std::vector<unsigned int>>(connectivityNode);
280 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
282 if (verbose) std::cout <<
"Found " << offsets.size() <<
" polylines." << std::endl;
284 unsigned int lastOffset = 0;
285 for (
unsigned int i = 0; i < offsets.size(); ++i)
289 unsigned int offset = offsets[i];
290 for (
unsigned int j = 0; j < offset-lastOffset-1; ++j)
292 std::vector<unsigned int>({ connectivity[lastOffset+j], connectivity[lastOffset+j+1] }));
297 DUNE_THROW(Dune::IOError,
"No Cells or Lines element found in " << fileName_);
306 void readGridData_(Data& cellData, Data& pointData,
bool verbose =
false)
const
308 using namespace tinyxml2;
310 const XMLElement* pieceNode = getPieceNode_();
311 const XMLElement* cellDataNode = getDataNode_(pieceNode, DataType::cellData);
312 if (cellDataNode !=
nullptr)
314 const XMLElement* cellsNode = pieceNode->FirstChildElement(
"Cells");
315 const XMLElement* linesNode = pieceNode->FirstChildElement(
"Lines");
319 const XMLElement* dataArray = cellDataNode->FirstChildElement(
"DataArray");
320 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
322 const char *attributeText = dataArray->Attribute(
"Name");
324 if (attributeText ==
nullptr)
325 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a cell data array.");
327 cellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
334 const XMLElement* dataArray = cellDataNode->FirstChildElement(
"DataArray");
337 Data polyLineCellData;
338 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
340 const char *attributeText = dataArray->Attribute(
"Name");
342 if (attributeText ==
nullptr)
343 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a cell data array.");
345 polyLineCellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
351 const XMLElement* offsetsNode = findDataArray_(linesNode,
"offsets");
352 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
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() <<
")!");
360 unsigned int lastOffset = 0;
361 std::size_t numCells = 0;
362 for (
unsigned int i = 0; i < offsets.size(); ++i)
364 unsigned int offset = offsets[i];
365 for (
unsigned int j = 0; j < offset-lastOffset-1; ++j)
371 for (
const auto& dArray : polyLineCellData)
373 cellData[dArray.first] = std::vector<double>(numCells);
374 auto& cd = cellData[dArray.first];
375 const auto& pd = dArray.second;
378 std::size_t cellIdx = 0;
379 for (
unsigned int i = 0; i < offsets.size(); ++i)
381 unsigned int offset = offsets[i];
382 for (
unsigned int j = 0; j < offset-lastOffset-1; ++j)
383 cd[cellIdx++] = pd[i];
390 DUNE_THROW(Dune::IOError,
"No Cells or Lines element found in " << fileName_);
393 const XMLElement* pointDataNode = getDataNode_(pieceNode, DataType::pointData);
394 if (pointDataNode !=
nullptr)
396 const XMLElement* dataArray = pointDataNode->FirstChildElement(
"DataArray");
397 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
399 const char *attributeText = dataArray->Attribute(
"Name");
401 if (attributeText ==
nullptr)
402 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a point data array.");
404 pointData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
413 const tinyxml2::XMLElement* getPieceNode_()
const
414 {
return getPieceNode_(doc_, fileName_); }
422 const tinyxml2::XMLElement* getPieceNode_(
const tinyxml2::XMLDocument& doc,
const std::string& fileName)
const
424 using namespace tinyxml2;
426 const XMLElement* pieceNode = doc.FirstChildElement(
"VTKFile");
427 if (pieceNode ==
nullptr)
428 DUNE_THROW(Dune::IOError,
"Couldn't get 'VTKFile' node in " << fileName <<
".");
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 <<
".");
440 return pieceNode->FirstChildElement(
"Piece");
449 const tinyxml2::XMLElement* getDataNode_(
const tinyxml2::XMLElement* pieceNode,
const DataType& type)
const
451 using namespace tinyxml2;
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");
459 DUNE_THROW(Dune::IOError,
"Only cell and point data are supported.");
470 const tinyxml2::XMLElement* findDataArray_(
const tinyxml2::XMLElement* dataNode,
const std::string& name)
const
472 using namespace tinyxml2;
475 const XMLElement* dataArray = dataNode->FirstChildElement(
"DataArray");
476 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
478 const char *attributeText = dataArray->Attribute(
"Name");
480 if (attributeText ==
nullptr)
481 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a data array.");
483 if (attributeText == name)
495 template<
class Container>
496 Container parseDataArray_(
const tinyxml2::XMLElement* dataArray)
const
498 std::stringstream dataStream(dataArray->GetText());
499 return readStreamToContainer<Container>(dataStream);
506 Dune::GeometryType vtkToDuneGeomType_(
unsigned int vtkCellType)
const
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);
521 template<
int dim, std::enable_if_t<(dim < 3),
int> = 0>
522 std::vector<Dune::FieldVector<
double, dim>>
523 adaptPo
intDimension_(std::vector<Dune::FieldVector<
double, 3>>&& po
ints3D) const
525 std::vector<Dune::FieldVector<
double, dim>> po
ints(po
ints3D.size());
526 for (std::
size_t i = 0; i < po
ints.size(); ++i)
527 for (
int j = 0; j < dim; ++j)
528 po
ints[i][j] = po
ints3D[i][j];
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); }
537 std::string fileName_;
538 tinyxml2::XMLDocument doc_;
Free functions to write and read a sequence container to and from a file.
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