12#ifndef DUMUX_IO_VTK_VTKREADER_HH
13#define DUMUX_IO_VTK_VTKREADER_HH
20#include <unordered_map>
24#include <dune/common/parallel/mpihelper.hh>
25#include <dune/common/exceptions.hh>
26#include <dune/grid/common/capabilities.hh>
27#include <dune/grid/io/file/vtk/common.hh>
28#include <dune/grid/common/gridfactory.hh>
31#include <dumux/io/xml/tinyxml2.h>
48 using Data = std::unordered_map<std::string, std::vector<double>>;
55 using namespace tinyxml2;
56 const auto ext = std::filesystem::path(fileName).extension().string();
60 fileName_ = Dune::MPIHelper::instance().size() > 1 && ext[1] ==
'p' ?
61 getProcessFileName_(fileName) : fileName;
63 const auto eResult = doc_.LoadFile(fileName_.c_str());
64 if (eResult != tinyxml2::XML_SUCCESS)
65 DUNE_THROW(Dune::IOError,
"Couldn't open XML file " << fileName_ <<
".");
67 const XMLElement* pieceNode = getPieceNode_();
68 if (pieceNode ==
nullptr)
69 DUNE_THROW(Dune::IOError,
"Couldn't get 'Piece' node in " << fileName_ <<
".");
79 using namespace tinyxml2;
81 const XMLElement* pieceNode = getPieceNode_();
82 const XMLElement* dataNode = getDataNode_(pieceNode, type);
83 if (dataNode ==
nullptr)
86 const XMLElement* dataArray = findDataArray_(dataNode, name);
87 if (dataArray ==
nullptr)
99 template<
class Container>
102 using namespace tinyxml2;
104 const XMLElement* pieceNode = getPieceNode_();
105 const XMLElement* dataNode = getDataNode_(pieceNode, type);
106 if (dataNode ==
nullptr)
107 DUNE_THROW(Dune::IOError,
"Couldn't get 'PointData' or 'CellData' node in " << fileName_ <<
".");
109 const XMLElement* dataArray = findDataArray_(dataNode, name);
110 if (dataArray ==
nullptr)
111 DUNE_THROW(Dune::IOError,
"Couldn't find the data array " << name <<
".");
113 return parseDataArray_<Container>(dataArray);
121 std::unique_ptr<Grid>
readGrid(
bool verbose =
false)
const
123 static_assert(!Dune::Capabilities::isCartesian<Grid>::v,
"Grid reader only supports unstructured grid implementations");
126 Dune::GridFactory<Grid> factory;
129 if (Dune::MPIHelper::instance().rank() == 0)
132 std::cout <<
"Reading " << Grid::dimension <<
"d grid from vtk file " << fileName_ <<
"." << std::endl;
133 readGrid_(factory, verbose);
136 return std::unique_ptr<Grid>(factory.createGrid());
146 std::unique_ptr<Grid>
readGrid(Dune::GridFactory<Grid>& factory,
bool verbose =
false)
const
148 static_assert(!Dune::Capabilities::isCartesian<Grid>::v,
"Grid reader only supports unstructured grid implementations");
150 if (Dune::MPIHelper::instance().rank() == 0)
153 std::cout <<
"Reading " << Grid::dimension <<
"d grid from vtk file " << fileName_ <<
"." << std::endl;
154 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 (Dune::MPIHelper::instance().rank() == 0)
176 std::cout <<
"Reading " << Grid::dimension <<
"d grid from vtk file " << fileName_ <<
"." << std::endl;
177 readGrid_(factory, verbose);
178 readGridData_(cellData, pointData, verbose);
181 return std::unique_ptr<Grid>(factory.createGrid());
188 std::string getProcessFileName_(
const std::string& pvtkFileName)
190 using namespace tinyxml2;
193 const auto eResult = pDoc.LoadFile(pvtkFileName.c_str());
194 if (eResult != XML_SUCCESS)
195 DUNE_THROW(Dune::IOError,
"Couldn't open XML file " << pvtkFileName <<
".");
198 const XMLElement* pieceNode = getPieceNode_(pDoc, pvtkFileName);
199 const auto myrank = Dune::MPIHelper::instance().rank();
200 for (
int rank = 0; rank < myrank; ++rank)
202 pieceNode = pieceNode->NextSiblingElement(
"Piece");
203 if (pieceNode ==
nullptr)
204 DUNE_THROW(Dune::IOError,
"Couldn't find 'Piece' node for rank "
205 << rank <<
" in " << pvtkFileName <<
".");
208 const char *vtkFileName = pieceNode->Attribute(
"Source");
209 if (vtkFileName ==
nullptr)
210 DUNE_THROW(Dune::IOError,
"Couldn't get 'Source' attribute of 'Piece' node no. " << myrank <<
" in " << pvtkFileName);
221 void readGrid_(Dune::GridFactory<Grid>& factory,
bool verbose =
false)
const
223 using namespace tinyxml2;
225 const XMLElement* pieceNode = getPieceNode_();
226 const XMLElement* pointsNode = pieceNode->FirstChildElement(
"Points")->FirstChildElement(
"DataArray");
227 if (pointsNode ==
nullptr)
228 DUNE_THROW(Dune::IOError,
"Couldn't get data array of points in " << fileName_ <<
".");
230 using Point3D = Dune::FieldVector<double, 3>;
231 std::vector<Point3D> points3D;
232 std::stringstream dataStream(pointsNode->GetText());
233 std::istream_iterator<Point3D> it(dataStream);
234 std::copy(it, std::istream_iterator<Point3D>(), std::back_inserter(points3D));
237 auto points = adaptPointDimension_<Grid::dimensionworld>(std::move(points3D));
239 if (verbose) std::cout <<
"Found " << points.size() <<
" vertices." << std::endl;
242 for (
auto&& point : points)
243 factory.insertVertex(std::move(point));
245 const XMLElement* cellsNode = pieceNode->FirstChildElement(
"Cells");
246 const XMLElement* linesNode = pieceNode->FirstChildElement(
"Lines");
249 const XMLElement* connectivityNode = findDataArray_(cellsNode,
"connectivity");
250 const XMLElement* offsetsNode = findDataArray_(cellsNode,
"offsets");
251 const XMLElement* typesNode = findDataArray_(cellsNode,
"types");
253 const auto connectivity = parseDataArray_<std::vector<unsigned int>>(connectivityNode);
254 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
255 const auto types = parseDataArray_<std::vector<unsigned int>>(typesNode);
257 if (verbose) std::cout <<
"Found " << offsets.size() <<
" elements." << std::endl;
259 unsigned int lastOffset = 0;
260 for (
unsigned int i = 0; i < offsets.size(); ++i)
262 const auto geomType = vtkToDuneGeomType_(types[i]);
263 unsigned int offset = offsets[i];
264 std::vector<unsigned int> corners; corners.resize(offset-lastOffset);
265 for (
unsigned int j = 0; j < offset-lastOffset; ++j)
266 corners[Dune::VTK::renumber(geomType, j)] = connectivity[lastOffset+j];
267 factory.insertElement(geomType, std::move(corners));
275 if (Grid::dimension != 1)
276 DUNE_THROW(Dune::IOError,
"Grid expects dimension " << Grid::dimension
277 <<
" but " << fileName_ <<
" contains a 1D grid.");
279 const XMLElement* connectivityNode = findDataArray_(linesNode,
"connectivity");
280 const XMLElement* offsetsNode = findDataArray_(linesNode,
"offsets");
282 const auto connectivity = parseDataArray_<std::vector<unsigned int>>(connectivityNode);
283 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
285 if (verbose) std::cout <<
"Found " << offsets.size() <<
" polylines." << std::endl;
287 unsigned int lastOffset = 0;
288 for (
unsigned int i = 0; i < offsets.size(); ++i)
292 unsigned int offset = offsets[i];
293 for (
unsigned int j = 0; j < offset-lastOffset-1; ++j)
295 std::vector<unsigned int>({ connectivity[lastOffset+j], connectivity[lastOffset+j+1] }));
300 DUNE_THROW(Dune::IOError,
"No Cells or Lines element found in " << fileName_);
309 void readGridData_(Data& cellData, Data& pointData,
bool verbose =
false)
const
311 using namespace tinyxml2;
313 const XMLElement* pieceNode = getPieceNode_();
314 const XMLElement* cellDataNode = getDataNode_(pieceNode, DataType::cellData);
315 if (cellDataNode !=
nullptr)
317 const XMLElement* cellsNode = pieceNode->FirstChildElement(
"Cells");
318 const XMLElement* linesNode = pieceNode->FirstChildElement(
"Lines");
322 const XMLElement* dataArray = cellDataNode->FirstChildElement(
"DataArray");
323 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
325 const char *attributeText = dataArray->Attribute(
"Name");
327 if (attributeText ==
nullptr)
328 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a cell data array.");
330 cellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
333 std::cout <<
"Read cell data field " << attributeText << std::endl;
340 const XMLElement* dataArray = cellDataNode->FirstChildElement(
"DataArray");
343 Data polyLineCellData;
344 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
346 const char *attributeText = dataArray->Attribute(
"Name");
348 if (attributeText ==
nullptr)
349 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a cell data array.");
351 polyLineCellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
354 std::cout <<
"Read cell data field " << attributeText << std::endl;
360 const XMLElement* offsetsNode = findDataArray_(linesNode,
"offsets");
361 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
363 if (offsets.size() != polyLineCellData.begin()->second.size())
364 DUNE_THROW(Dune::IOError,
"Expected the same number of cell data entries (is "
365 << polyLineCellData.begin()->second.size()
366 <<
") as polylines (" << offsets.size() <<
")!");
369 unsigned int lastOffset = 0;
370 std::size_t numCells = 0;
371 for (
unsigned int i = 0; i < offsets.size(); ++i)
373 unsigned int offset = offsets[i];
374 for (
unsigned int j = 0; j < offset-lastOffset-1; ++j)
380 for (
const auto& dArray : polyLineCellData)
382 cellData[dArray.first] = std::vector<double>(numCells);
383 auto& cd = cellData[dArray.first];
384 const auto& pd = dArray.second;
387 std::size_t cellIdx = 0;
388 for (
unsigned int i = 0; i < offsets.size(); ++i)
390 unsigned int offset = offsets[i];
391 for (
unsigned int j = 0; j < offset-lastOffset-1; ++j)
392 cd[cellIdx++] = pd[i];
399 DUNE_THROW(Dune::IOError,
"No Cells or Lines element found in " << fileName_);
402 const XMLElement* pointDataNode = getDataNode_(pieceNode, DataType::pointData);
403 if (pointDataNode !=
nullptr)
405 const XMLElement* dataArray = pointDataNode->FirstChildElement(
"DataArray");
406 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
408 const char *attributeText = dataArray->Attribute(
"Name");
410 if (attributeText ==
nullptr)
411 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a point data array.");
413 pointData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
416 std::cout <<
"Read point data field " << attributeText << std::endl;
425 const tinyxml2::XMLElement* getPieceNode_()
const
426 {
return getPieceNode_(doc_, fileName_); }
434 const tinyxml2::XMLElement* getPieceNode_(
const tinyxml2::XMLDocument& doc,
const std::string& fileName)
const
436 using namespace tinyxml2;
438 const XMLElement* pieceNode = doc.FirstChildElement(
"VTKFile");
439 if (pieceNode ==
nullptr)
440 DUNE_THROW(Dune::IOError,
"Couldn't get 'VTKFile' node in " << fileName <<
".");
442 pieceNode = pieceNode->FirstChildElement(
"UnstructuredGrid");
443 if (pieceNode ==
nullptr)
444 pieceNode = doc.FirstChildElement(
"VTKFile")->FirstChildElement(
"PolyData");
445 if (pieceNode ==
nullptr)
446 pieceNode = doc.FirstChildElement(
"VTKFile")->FirstChildElement(
"PUnstructuredGrid");
447 if (pieceNode ==
nullptr)
448 pieceNode = doc.FirstChildElement(
"VTKFile")->FirstChildElement(
"PPolyData");
449 if (pieceNode ==
nullptr)
450 DUNE_THROW(Dune::IOError,
"Couldn't get 'UnstructuredGrid', 'PUnstructuredGrid', 'PolyData', or 'PPolyData' node in " << fileName <<
".");
452 return pieceNode->FirstChildElement(
"Piece");
461 const tinyxml2::XMLElement* getDataNode_(
const tinyxml2::XMLElement* pieceNode,
const DataType& type)
const
463 using namespace tinyxml2;
465 const XMLElement* dataNode =
nullptr;
466 if (type == DataType::pointData)
467 dataNode = pieceNode->FirstChildElement(
"PointData");
468 else if (type == DataType::cellData)
469 dataNode = pieceNode->FirstChildElement(
"CellData");
471 DUNE_THROW(Dune::IOError,
"Only cell and point data are supported.");
482 const tinyxml2::XMLElement* findDataArray_(
const tinyxml2::XMLElement* dataNode,
const std::string& name)
const
484 using namespace tinyxml2;
487 const XMLElement* dataArray = dataNode->FirstChildElement(
"DataArray");
488 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
490 const char *attributeText = dataArray->Attribute(
"Name");
492 if (attributeText ==
nullptr)
493 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a data array.");
495 if (attributeText == name)
507 template<
class Container>
508 Container parseDataArray_(
const tinyxml2::XMLElement* dataArray)
const
510 std::stringstream dataStream(dataArray->GetText());
511 return readStreamToContainer<Container>(dataStream);
518 Dune::GeometryType vtkToDuneGeomType_(
unsigned int vtkCellType)
const
522 case Dune::VTK::GeometryType::vertex:
return Dune::GeometryTypes::vertex;
524 case Dune::VTK::GeometryType::triangle:
return Dune::GeometryTypes::triangle;
525 case Dune::VTK::GeometryType::quadrilateral:
return Dune::GeometryTypes::quadrilateral;
526 case Dune::VTK::GeometryType::tetrahedron:
return Dune::GeometryTypes::tetrahedron;
527 case Dune::VTK::GeometryType::hexahedron:
return Dune::GeometryTypes::hexahedron;
528 case Dune::VTK::GeometryType::prism:
return Dune::GeometryTypes::prism;
529 case Dune::VTK::GeometryType::pyramid:
return Dune::GeometryTypes::pyramid;
530 default: DUNE_THROW(Dune::NotImplemented,
"VTK cell type " << vtkCellType);
534 template<
int dim, std::enable_if_t<(dim < 3),
int> = 0>
535 std::vector<Dune::FieldVector<
double, dim>>
536 adaptPo
intDimension_(std::vector<Dune::FieldVector<
double, 3>>&& po
ints3D) const
538 std::vector<Dune::FieldVector<
double, dim>> po
ints(po
ints3D.size());
539 for (std::
size_t i = 0; i < po
ints.size(); ++i)
540 for (
int j = 0; j < dim; ++j)
541 po
ints[i][j] = po
ints3D[i][j];
545 template<
int dim, std::enable_if_t<(dim == 3),
int> = 0>
546 std::vector<Dune::FieldVector<double, dim>>
547 adaptPointDimension_(std::vector<Dune::FieldVector<double, 3>>&& points3D)
const
548 {
return std::move(points3D); }
550 std::string fileName_;
551 tinyxml2::XMLDocument doc_;
A vtk file reader using tinyxml2 as xml backend.
Definition: vtkreader.hh:40
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:77
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:100
VTKReader(const std::string &fileName)
The constructor creates a tinyxml2::XMLDocument from file.
Definition: vtkreader.hh:53
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:146
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:121
DataType
The data array types.
Definition: vtkreader.hh:45
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:48
Free functions to write and read a sequence container to and from a file.
constexpr Line line
Definition: couplingmanager1d3d_line.hh:31