12#ifndef DUMUX_IO_VTK_VTKREADER_HH
13#define DUMUX_IO_VTK_VTKREADER_HH
20#include <unordered_map>
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>
30#include <dumux/io/xml/tinyxml2.h>
47 using Data = std::unordered_map<std::string, std::vector<double>>;
54 using namespace tinyxml2;
55 fileName_ = Dune::MPIHelper::getCommunication().size() > 1 ?
56 getProcessFileName_(fileName) : fileName;
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_ <<
".");
62 const XMLElement* pieceNode = getPieceNode_();
63 if (pieceNode ==
nullptr)
64 DUNE_THROW(Dune::IOError,
"Couldn't get 'Piece' node in " << fileName_ <<
".");
74 using namespace tinyxml2;
76 const XMLElement* pieceNode = getPieceNode_();
77 const XMLElement* dataNode = getDataNode_(pieceNode, type);
78 if (dataNode ==
nullptr)
81 const XMLElement* dataArray = findDataArray_(dataNode, name);
82 if (dataArray ==
nullptr)
94 template<
class Container>
97 using namespace tinyxml2;
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_ <<
".");
104 const XMLElement* dataArray = findDataArray_(dataNode, name);
105 if (dataArray ==
nullptr)
106 DUNE_THROW(Dune::IOError,
"Couldn't find the data array " << name <<
".");
108 return parseDataArray_<Container>(dataArray);
116 std::unique_ptr<Grid>
readGrid(
bool verbose =
false)
const
118 static_assert(!Dune::Capabilities::isCartesian<Grid>::v,
"Grid reader only supports unstructured grid implementations");
120 if (verbose) std::cout <<
"Reading " << Grid::dimension <<
"d grid from vtk file " << fileName_ <<
"." << std::endl;
123 Dune::GridFactory<Grid> factory;
125 readGrid_(factory, verbose);
127 return std::unique_ptr<Grid>(factory.createGrid());
137 std::unique_ptr<Grid>
readGrid(Dune::GridFactory<Grid>& factory,
bool verbose =
false)
const
139 static_assert(!Dune::Capabilities::isCartesian<Grid>::v,
"Grid reader only supports unstructured grid implementations");
141 if (verbose) std::cout <<
"Reading " << Grid::dimension <<
"d grid from vtk file " << fileName_ <<
"." << std::endl;
143 readGrid_(factory, verbose);
145 return std::unique_ptr<Grid>(factory.createGrid());
157 std::unique_ptr<Grid>
readGrid(Dune::GridFactory<Grid>& factory,
Data& cellData,
Data& pointData,
bool verbose =
false)
const
159 static_assert(!Dune::Capabilities::isCartesian<Grid>::v,
"Grid reader only supports unstructured grid implementations");
161 if (verbose) std::cout <<
"Reading " << Grid::dimension <<
"d grid from vtk file " << fileName_ <<
"." << std::endl;
163 readGrid_(factory, verbose);
164 readGridData_(cellData, pointData, verbose);
166 return std::unique_ptr<Grid>(factory.createGrid());
173 std::string getProcessFileName_(
const std::string& pvtkFileName)
175 using namespace tinyxml2;
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 <<
".");
183 const XMLElement* pieceNode = getPieceNode_(pDoc, pvtkFileName);
184 const auto myrank = Dune::MPIHelper::getCommunication().rank();
185 for (
int rank = 0; rank < myrank; ++rank)
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 <<
".");
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);
206 void readGrid_(Dune::GridFactory<Grid>& factory,
bool verbose =
false)
const
208 using namespace tinyxml2;
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_ <<
".");
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));
222 auto points = adaptPointDimension_<Grid::dimensionworld>(std::move(points3D));
224 if (verbose) std::cout <<
"Found " << points.size() <<
" vertices." << std::endl;
227 for (
auto&& point : points)
228 factory.insertVertex(std::move(point));
230 const XMLElement* cellsNode = pieceNode->FirstChildElement(
"Cells");
231 const XMLElement* linesNode = pieceNode->FirstChildElement(
"Lines");
234 const XMLElement* connectivityNode = findDataArray_(cellsNode,
"connectivity");
235 const XMLElement* offsetsNode = findDataArray_(cellsNode,
"offsets");
236 const XMLElement* typesNode = findDataArray_(cellsNode,
"types");
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);
242 if (verbose) std::cout <<
"Found " << offsets.size() <<
" element." << std::endl;
244 unsigned int lastOffset = 0;
245 for (
unsigned int i = 0; i < offsets.size(); ++i)
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));
260 if (Grid::dimension != 1)
261 DUNE_THROW(Dune::IOError,
"Grid expects dimension " << Grid::dimension
262 <<
" but " << fileName_ <<
" contains a 1D grid.");
264 const XMLElement* connectivityNode = findDataArray_(linesNode,
"connectivity");
265 const XMLElement* offsetsNode = findDataArray_(linesNode,
"offsets");
267 const auto connectivity = parseDataArray_<std::vector<unsigned int>>(connectivityNode);
268 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
270 if (verbose) std::cout <<
"Found " << offsets.size() <<
" polylines." << std::endl;
272 unsigned int lastOffset = 0;
273 for (
unsigned int i = 0; i < offsets.size(); ++i)
277 unsigned int offset = offsets[i];
278 for (
unsigned int j = 0; j < offset-lastOffset-1; ++j)
280 std::vector<unsigned int>({ connectivity[lastOffset+j], connectivity[lastOffset+j+1] }));
285 DUNE_THROW(Dune::IOError,
"No Cells or Lines element found in " << fileName_);
294 void readGridData_(Data& cellData, Data& pointData,
bool verbose =
false)
const
296 using namespace tinyxml2;
298 const XMLElement* pieceNode = getPieceNode_();
299 const XMLElement* cellDataNode = getDataNode_(pieceNode, DataType::cellData);
300 if (cellDataNode !=
nullptr)
302 const XMLElement* cellsNode = pieceNode->FirstChildElement(
"Cells");
303 const XMLElement* linesNode = pieceNode->FirstChildElement(
"Lines");
307 const XMLElement* dataArray = cellDataNode->FirstChildElement(
"DataArray");
308 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
310 const char *attributeText = dataArray->Attribute(
"Name");
312 if (attributeText ==
nullptr)
313 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a cell data array.");
315 cellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
322 const XMLElement* dataArray = cellDataNode->FirstChildElement(
"DataArray");
325 Data polyLineCellData;
326 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
328 const char *attributeText = dataArray->Attribute(
"Name");
330 if (attributeText ==
nullptr)
331 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a cell data array.");
333 polyLineCellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
339 const XMLElement* offsetsNode = findDataArray_(linesNode,
"offsets");
340 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
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() <<
")!");
348 unsigned int lastOffset = 0;
349 std::size_t numCells = 0;
350 for (
unsigned int i = 0; i < offsets.size(); ++i)
352 unsigned int offset = offsets[i];
353 for (
unsigned int j = 0; j < offset-lastOffset-1; ++j)
359 for (
const auto& dArray : polyLineCellData)
361 cellData[dArray.first] = std::vector<double>(numCells);
362 auto& cd = cellData[dArray.first];
363 const auto& pd = dArray.second;
366 std::size_t cellIdx = 0;
367 for (
unsigned int i = 0; i < offsets.size(); ++i)
369 unsigned int offset = offsets[i];
370 for (
unsigned int j = 0; j < offset-lastOffset-1; ++j)
371 cd[cellIdx++] = pd[i];
378 DUNE_THROW(Dune::IOError,
"No Cells or Lines element found in " << fileName_);
381 const XMLElement* pointDataNode = getDataNode_(pieceNode, DataType::pointData);
382 if (pointDataNode !=
nullptr)
384 const XMLElement* dataArray = pointDataNode->FirstChildElement(
"DataArray");
385 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
387 const char *attributeText = dataArray->Attribute(
"Name");
389 if (attributeText ==
nullptr)
390 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a point data array.");
392 pointData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
401 const tinyxml2::XMLElement* getPieceNode_()
const
402 {
return getPieceNode_(doc_, fileName_); }
410 const tinyxml2::XMLElement* getPieceNode_(
const tinyxml2::XMLDocument& doc,
const std::string& fileName)
const
412 using namespace tinyxml2;
414 const XMLElement* pieceNode = doc.FirstChildElement(
"VTKFile");
415 if (pieceNode ==
nullptr)
416 DUNE_THROW(Dune::IOError,
"Couldn't get 'VTKFile' node in " << fileName <<
".");
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 <<
".");
428 return pieceNode->FirstChildElement(
"Piece");
437 const tinyxml2::XMLElement* getDataNode_(
const tinyxml2::XMLElement* pieceNode,
const DataType& type)
const
439 using namespace tinyxml2;
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");
447 DUNE_THROW(Dune::IOError,
"Only cell and point data are supported.");
458 const tinyxml2::XMLElement* findDataArray_(
const tinyxml2::XMLElement* dataNode,
const std::string& name)
const
460 using namespace tinyxml2;
463 const XMLElement* dataArray = dataNode->FirstChildElement(
"DataArray");
464 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
466 const char *attributeText = dataArray->Attribute(
"Name");
468 if (attributeText ==
nullptr)
469 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a data array.");
471 if (attributeText == name)
483 template<
class Container>
484 Container parseDataArray_(
const tinyxml2::XMLElement* dataArray)
const
486 std::stringstream dataStream(dataArray->GetText());
487 return readStreamToContainer<Container>(dataStream);
494 Dune::GeometryType vtkToDuneGeomType_(
unsigned int vtkCellType)
const
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);
509 template<
int dim, std::enable_if_t<(dim < 3),
int> = 0>
510 std::vector<Dune::FieldVector<
double, dim>>
511 adaptPo
intDimension_(std::vector<Dune::FieldVector<
double, 3>>&& po
ints3D) const
513 std::vector<Dune::FieldVector<
double, dim>> po
ints(po
ints3D.size());
514 for (std::
size_t i = 0; i < po
ints.size(); ++i)
515 for (
int j = 0; j < dim; ++j)
516 po
ints[i][j] = po
ints3D[i][j];
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); }
525 std::string fileName_;
526 tinyxml2::XMLDocument doc_;
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