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>
32#if DUMUX_HAVE_GRIDFORMAT
34#include <gridformat/gridformat.hpp>
35#include <gridformat/traits/dune.hpp>
36#include <gridformat/reader.hpp>
37#include <gridformat/decorators/reader_polylines_subdivider.hpp>
39#include "imagegrid_.hh"
43template<Gr
idFormat::Concepts::Communicator C>
44auto makeGridformatReaderFactory(
const C& c)
46 return [fac = GridFormat::AnyReaderFactory<C>{c}] (
const std::string& f)
47 -> std::unique_ptr<GridFormat::GridReader> {
49 if (f.ends_with(
"vtp"))
50 return std::make_unique<GridFormat::ReaderDecorators::PolylinesSubdivider>(fac.make_for(f));
51 return fac.make_for(f);
66 enum class DataType { cellData, pointData, fieldData };
69 using Data = std::unordered_map<std::string, std::vector<double>>;
71 explicit VTKReader(
const std::string& fileName)
72 : reader_(std::make_shared<GridFormat::Reader>(GridFormat::Reader::from(fileName,
73 Detail::
VTKReader::makeGridformatReaderFactory(
74 Dune::MPIHelper::instance().getCommunicator()
87 return std::ranges::any_of(
88 cell_field_names(*reader_),
89 [&] (
const auto& n) {
return name == n; }
92 return std::ranges::any_of(
93 point_field_names(*reader_),
94 [&] (
const auto& n) {
return name == n; }
96 else if (type == DataType::fieldData)
97 return std::ranges::any_of(
98 meta_data_field_names(*reader_),
99 [&] (
const auto& n) {
return name == n; }
102 DUNE_THROW(Dune::IOError,
"Unknown data type for field '" << name <<
"'");
111 template<
class Container>
116 Container values(reader_->number_of_cells());
117 reader_->cell_field(name)->export_to(values);
122 Container values(reader_->number_of_points());
123 reader_->point_field(name)->export_to(values);
126 else if (type == DataType::fieldData)
128 DUNE_THROW(Dune::NotImplemented,
"Reading meta data not yet implemented");
131 DUNE_THROW(Dune::IOError,
"Unknown data type for field '" << name <<
"'");
139 std::unique_ptr<Grid>
readGrid(
bool verbose =
false)
const
141 Dune::GridFactory<Grid> factory;
152 std::unique_ptr<Grid>
readGrid(Dune::GridFactory<Grid>& factory,
bool verbose =
false)
const
154 if (Dune::MPIHelper::instance().rank() == 0)
157 std::cout <<
"Reading " << Grid::dimension <<
"d grid from vtk file " << reader_->filename() <<
"." << std::endl;
160 GridFormat::Dune::GridFactoryAdapter<Grid> adapter{ factory };
161 reader_->export_grid(adapter);
165 return std::unique_ptr<Grid>(factory.createGrid());
177 std::unique_ptr<Grid>
readGrid(Dune::GridFactory<Grid>& factory,
Data& cellData,
Data& pointData,
bool verbose =
false)
const
179 auto grid =
readGrid(factory, verbose);
182 for (
const auto& name : cell_field_names(*reader_))
183 cellData[name] = Data::mapped_type{};
184 for (
const auto& name : point_field_names(*reader_))
185 pointData[name] = Data::mapped_type{};
188 if (Dune::MPIHelper::instance().rank() == 0)
190 for (
const auto& [name, field_ptr] : cell_fields(*reader_))
192 if (verbose) std::cout <<
"-- reading cell data '" << name <<
"'." << std::endl;
193 field_ptr->export_to(cellData[name]);
196 for (
const auto& [name, field_ptr] : point_fields(*reader_))
198 if (verbose) std::cout <<
"-- reading point data '" << name <<
"'." << std::endl;
199 field_ptr->export_to(pointData[name]);
203 return std::unique_ptr<Grid>(std::move(grid));
216 template<
class ct,
int dim>
217 std::unique_ptr<Detail::VTKReader::ImageGrid<ct, dim>> readStructuredGrid(
bool verbose =
false)
const
220 std::cout <<
"Reading " << dim <<
"d grid from vtk file " << reader_->filename() <<
"." << std::endl;
222 return std::make_unique<Detail::VTKReader::ImageGrid<ct, dim>>(reader_);
233 template<
class ct,
int dim>
234 std::unique_ptr<Detail::VTKReader::ImageGrid<ct, dim>> readStructuredGrid(
Data& cellData,
Data& pointData,
bool verbose =
false)
const
236 auto grid = readStructuredGrid<ct, dim>(verbose);
239 for (
const auto& name : cell_field_names(*reader_))
240 cellData[name] = Data::mapped_type{};
241 for (
const auto& name : point_field_names(*reader_))
242 pointData[name] = Data::mapped_type{};
245 if (Dune::MPIHelper::instance().rank() == 0)
247 for (
const auto& [name, field_ptr] : cell_fields(*reader_))
249 if (verbose) std::cout <<
"-- reading cell data '" << name <<
"'." << std::endl;
250 field_ptr->export_to(cellData[name]);
253 for (
const auto& [name, field_ptr] : point_fields(*reader_))
255 if (verbose) std::cout <<
"-- reading point data '" << name <<
"'." << std::endl;
256 field_ptr->export_to(pointData[name]);
259 return { std::move(grid) };
263 std::shared_ptr<GridFormat::Reader> reader_;
270#include <dumux/io/xml/tinyxml2.h>
282const tinyxml2::XMLElement*
getPieceNode(
const tinyxml2::XMLDocument& doc,
const std::string& fileName)
284 using namespace tinyxml2;
286 const XMLElement* pieceNode = doc.FirstChildElement(
"VTKFile");
287 if (pieceNode ==
nullptr)
288 DUNE_THROW(Dune::IOError,
"Couldn't get 'VTKFile' node in " << fileName <<
".");
290 pieceNode = pieceNode->FirstChildElement(
"UnstructuredGrid");
291 if (pieceNode ==
nullptr)
292 pieceNode = doc.FirstChildElement(
"VTKFile")->FirstChildElement(
"PolyData");
293 if (pieceNode ==
nullptr)
294 pieceNode = doc.FirstChildElement(
"VTKFile")->FirstChildElement(
"PUnstructuredGrid");
295 if (pieceNode ==
nullptr)
296 pieceNode = doc.FirstChildElement(
"VTKFile")->FirstChildElement(
"PPolyData");
297 if (pieceNode ==
nullptr)
298 DUNE_THROW(Dune::IOError,
"Couldn't get 'UnstructuredGrid', 'PUnstructuredGrid', 'PolyData', or 'PPolyData' node in " << fileName <<
".");
300 return pieceNode->FirstChildElement(
"Piece");
308 using namespace tinyxml2;
311 const auto eResult = pDoc.LoadFile(pvtkFileName.c_str());
312 if (eResult != XML_SUCCESS)
313 DUNE_THROW(Dune::IOError,
"Couldn't open XML file " << pvtkFileName <<
".");
316 const XMLElement* pieceNode =
getPieceNode(pDoc, pvtkFileName);
317 const auto myrank = Dune::MPIHelper::instance().rank();
318 for (
int rank = 0; rank < myrank; ++rank)
320 pieceNode = pieceNode->NextSiblingElement(
"Piece");
321 if (pieceNode ==
nullptr)
322 DUNE_THROW(Dune::IOError,
"Couldn't find 'Piece' node for rank "
323 << rank <<
" in " << pvtkFileName <<
".");
326 const char *vtkFileName = pieceNode->Attribute(
"Source");
327 if (vtkFileName ==
nullptr)
328 DUNE_THROW(Dune::IOError,
"Couldn't get 'Source' attribute of 'Piece' node no. " << myrank <<
" in " << pvtkFileName);
350 using Data = std::unordered_map<std::string, std::vector<double>>;
357 using namespace tinyxml2;
358 const auto ext = std::filesystem::path(fileName).extension().string();
362 fileName_ = Dune::MPIHelper::instance().size() > 1 && ext[1] ==
'p' ?
365 const auto eResult = doc_.LoadFile(fileName_.c_str());
366 if (eResult != tinyxml2::XML_SUCCESS)
367 DUNE_THROW(Dune::IOError,
"Couldn't open XML file " << fileName_ <<
".");
369 const XMLElement* pieceNode = getPieceNode_();
370 if (pieceNode ==
nullptr)
371 DUNE_THROW(Dune::IOError,
"Couldn't get 'Piece' node in " << fileName_ <<
".");
381 using namespace tinyxml2;
383 const XMLElement* pieceNode = getPieceNode_();
384 const XMLElement* dataNode = getDataNode_(pieceNode, type);
385 if (dataNode ==
nullptr)
388 const XMLElement* dataArray = findDataArray_(dataNode, name);
389 if (dataArray ==
nullptr)
401 template<
class Container>
404 using namespace tinyxml2;
406 const XMLElement* pieceNode = getPieceNode_();
407 const XMLElement* dataNode = getDataNode_(pieceNode, type);
408 if (dataNode ==
nullptr)
409 DUNE_THROW(Dune::IOError,
"Couldn't get 'PointData' or 'CellData' node in " << fileName_ <<
".");
411 const XMLElement* dataArray = findDataArray_(dataNode, name);
412 if (dataArray ==
nullptr)
413 DUNE_THROW(Dune::IOError,
"Couldn't find the data array " << name <<
".");
415 return parseDataArray_<Container>(dataArray);
423 std::unique_ptr<Grid>
readGrid(
bool verbose =
false)
const
425 static_assert(!Dune::Capabilities::isCartesian<Grid>::v,
"Grid reader only supports unstructured grid implementations");
428 Dune::GridFactory<Grid> factory;
431 if (Dune::MPIHelper::instance().rank() == 0)
434 std::cout <<
"Reading " << Grid::dimension <<
"d grid from vtk file " << fileName_ <<
"." << std::endl;
435 readGrid_(factory, verbose);
438 return std::unique_ptr<Grid>(factory.createGrid());
448 std::unique_ptr<Grid>
readGrid(Dune::GridFactory<Grid>& factory,
bool verbose =
false)
const
450 static_assert(!Dune::Capabilities::isCartesian<Grid>::v,
"Grid reader only supports unstructured grid implementations");
452 if (Dune::MPIHelper::instance().rank() == 0)
455 std::cout <<
"Reading " << Grid::dimension <<
"d grid from vtk file " << fileName_ <<
"." << std::endl;
456 readGrid_(factory, verbose);
459 return std::unique_ptr<Grid>(factory.createGrid());
471 std::unique_ptr<Grid>
readGrid(Dune::GridFactory<Grid>& factory,
Data& cellData,
Data& pointData,
bool verbose =
false)
const
473 static_assert(!Dune::Capabilities::isCartesian<Grid>::v,
"Grid reader only supports unstructured grid implementations");
475 if (Dune::MPIHelper::instance().rank() == 0)
478 std::cout <<
"Reading " << Grid::dimension <<
"d grid from vtk file " << fileName_ <<
"." << std::endl;
479 readGrid_(factory, verbose);
480 readGridData_(cellData, pointData, verbose);
483 return std::unique_ptr<Grid>(factory.createGrid());
493 void readGrid_(Dune::GridFactory<Grid>& factory,
bool verbose =
false)
const
495 using namespace tinyxml2;
497 const XMLElement* pieceNode = getPieceNode_();
498 const XMLElement* pointsNode = pieceNode->FirstChildElement(
"Points")->FirstChildElement(
"DataArray");
499 if (pointsNode ==
nullptr)
500 DUNE_THROW(Dune::IOError,
"Couldn't get data array of points in " << fileName_ <<
".");
502 using Point3D = Dune::FieldVector<double, 3>;
503 std::vector<Point3D> points3D;
504 std::stringstream dataStream(pointsNode->GetText());
505 std::istream_iterator<Point3D> it(dataStream);
506 std::copy(it, std::istream_iterator<Point3D>(), std::back_inserter(points3D));
509 auto points = adaptPointDimension_<Grid::dimensionworld>(std::move(points3D));
511 if (verbose) std::cout <<
"Found " << points.size() <<
" vertices." << std::endl;
514 for (
auto&& point : points)
515 factory.insertVertex(std::move(point));
517 const XMLElement* cellsNode = pieceNode->FirstChildElement(
"Cells");
518 const XMLElement* linesNode = pieceNode->FirstChildElement(
"Lines");
521 const XMLElement* connectivityNode = findDataArray_(cellsNode,
"connectivity");
522 const XMLElement* offsetsNode = findDataArray_(cellsNode,
"offsets");
523 const XMLElement* typesNode = findDataArray_(cellsNode,
"types");
525 const auto connectivity = parseDataArray_<std::vector<unsigned int>>(connectivityNode);
526 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
527 const auto types = parseDataArray_<std::vector<unsigned int>>(typesNode);
529 if (verbose) std::cout <<
"Found " << offsets.size() <<
" elements." << std::endl;
531 unsigned int lastOffset = 0;
532 for (
unsigned int i = 0; i < offsets.size(); ++i)
534 const auto geomType = vtkToDuneGeomType_(types[i]);
535 unsigned int offset = offsets[i];
536 std::vector<unsigned int> corners; corners.resize(offset-lastOffset);
537 for (
unsigned int j = 0; j < offset-lastOffset; ++j)
538 corners[Dune::VTK::renumber(geomType, j)] = connectivity[lastOffset+j];
539 factory.insertElement(geomType, std::move(corners));
547 if (Grid::dimension != 1)
548 DUNE_THROW(Dune::IOError,
"Grid expects dimension " << Grid::dimension
549 <<
" but " << fileName_ <<
" contains a 1D grid.");
551 const XMLElement* connectivityNode = findDataArray_(linesNode,
"connectivity");
552 const XMLElement* offsetsNode = findDataArray_(linesNode,
"offsets");
554 const auto connectivity = parseDataArray_<std::vector<unsigned int>>(connectivityNode);
555 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
557 if (verbose) std::cout <<
"Found " << offsets.size() <<
" polylines." << std::endl;
559 unsigned int lastOffset = 0;
560 for (
unsigned int i = 0; i < offsets.size(); ++i)
564 unsigned int offset = offsets[i];
565 for (
unsigned int j = 0; j < offset-lastOffset-1; ++j)
567 std::vector<unsigned int>({ connectivity[lastOffset+j], connectivity[lastOffset+j+1] }));
572 DUNE_THROW(Dune::IOError,
"No Cells or Lines element found in " << fileName_);
581 void readGridData_(Data& cellData, Data& pointData,
bool verbose =
false)
const
583 using namespace tinyxml2;
585 const XMLElement* pieceNode = getPieceNode_();
586 const XMLElement* cellDataNode = getDataNode_(pieceNode, DataType::cellData);
587 if (cellDataNode !=
nullptr)
589 const XMLElement* cellsNode = pieceNode->FirstChildElement(
"Cells");
590 const XMLElement* linesNode = pieceNode->FirstChildElement(
"Lines");
594 const XMLElement* dataArray = cellDataNode->FirstChildElement(
"DataArray");
595 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
597 const char *attributeText = dataArray->Attribute(
"Name");
599 if (attributeText ==
nullptr)
600 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a cell data array.");
602 cellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
605 std::cout <<
"Read cell data field " << attributeText << std::endl;
612 const XMLElement* dataArray = cellDataNode->FirstChildElement(
"DataArray");
615 Data polyLineCellData;
616 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
618 const char *attributeText = dataArray->Attribute(
"Name");
620 if (attributeText ==
nullptr)
621 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a cell data array.");
623 polyLineCellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
626 std::cout <<
"Read cell data field " << attributeText << std::endl;
632 const XMLElement* offsetsNode = findDataArray_(linesNode,
"offsets");
633 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
635 if (offsets.size() != polyLineCellData.begin()->second.size())
636 DUNE_THROW(Dune::IOError,
"Expected the same number of cell data entries (is "
637 << polyLineCellData.begin()->second.size()
638 <<
") as polylines (" << offsets.size() <<
")!");
641 unsigned int lastOffset = 0;
642 std::size_t numCells = 0;
643 for (
unsigned int i = 0; i < offsets.size(); ++i)
645 unsigned int offset = offsets[i];
646 for (
unsigned int j = 0; j < offset-lastOffset-1; ++j)
652 for (
const auto& dArray : polyLineCellData)
654 cellData[dArray.first] = std::vector<double>(numCells);
655 auto& cd = cellData[dArray.first];
656 const auto& pd = dArray.second;
659 std::size_t cellIdx = 0;
660 for (
unsigned int i = 0; i < offsets.size(); ++i)
662 unsigned int offset = offsets[i];
663 for (
unsigned int j = 0; j < offset-lastOffset-1; ++j)
664 cd[cellIdx++] = pd[i];
671 DUNE_THROW(Dune::IOError,
"No Cells or Lines element found in " << fileName_);
674 const XMLElement* pointDataNode = getDataNode_(pieceNode, DataType::pointData);
675 if (pointDataNode !=
nullptr)
677 const XMLElement* dataArray = pointDataNode->FirstChildElement(
"DataArray");
678 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
680 const char *attributeText = dataArray->Attribute(
"Name");
682 if (attributeText ==
nullptr)
683 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a point data array.");
685 pointData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
688 std::cout <<
"Read point data field " << attributeText << std::endl;
697 const tinyxml2::XMLElement* getPieceNode_()
const
706 const tinyxml2::XMLElement* getDataNode_(
const tinyxml2::XMLElement* pieceNode,
const DataType& type)
const
708 using namespace tinyxml2;
710 const XMLElement* dataNode =
nullptr;
711 if (type == DataType::pointData)
712 dataNode = pieceNode->FirstChildElement(
"PointData");
713 else if (type == DataType::cellData)
714 dataNode = pieceNode->FirstChildElement(
"CellData");
716 DUNE_THROW(Dune::IOError,
"Only cell and point data are supported.");
727 const tinyxml2::XMLElement* findDataArray_(
const tinyxml2::XMLElement* dataNode,
const std::string& name)
const
729 using namespace tinyxml2;
732 const XMLElement* dataArray = dataNode->FirstChildElement(
"DataArray");
733 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
735 const char *attributeText = dataArray->Attribute(
"Name");
737 if (attributeText ==
nullptr)
738 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a data array.");
740 if (attributeText == name)
752 template<
class Container>
753 Container parseDataArray_(
const tinyxml2::XMLElement* dataArray)
const
755 std::stringstream dataStream(dataArray->GetText());
756 return readStreamToContainer<Container>(dataStream);
763 Dune::GeometryType vtkToDuneGeomType_(
unsigned int vtkCellType)
const
767 case Dune::VTK::GeometryType::vertex:
return Dune::GeometryTypes::vertex;
769 case Dune::VTK::GeometryType::triangle:
return Dune::GeometryTypes::triangle;
770 case Dune::VTK::GeometryType::quadrilateral:
return Dune::GeometryTypes::quadrilateral;
771 case Dune::VTK::GeometryType::tetrahedron:
return Dune::GeometryTypes::tetrahedron;
772 case Dune::VTK::GeometryType::hexahedron:
return Dune::GeometryTypes::hexahedron;
773 case Dune::VTK::GeometryType::prism:
return Dune::GeometryTypes::prism;
774 case Dune::VTK::GeometryType::pyramid:
return Dune::GeometryTypes::pyramid;
775 default: DUNE_THROW(Dune::NotImplemented,
"VTK cell type " << vtkCellType);
779 template<
int dim, std::enable_if_t<(dim < 3),
int> = 0>
780 std::vector<Dune::FieldVector<
double, dim>>
781 adaptPo
intDimension_(std::vector<Dune::FieldVector<
double, 3>>&& po
ints3D) const
783 std::vector<Dune::FieldVector<
double, dim>> po
ints(po
ints3D.size());
784 for (std::
size_t i = 0; i < po
ints.size(); ++i)
785 for (
int j = 0; j < dim; ++j)
786 po
ints[i][j] = po
ints3D[i][j];
790 template<
int dim, std::enable_if_t<(dim == 3),
int> = 0>
791 std::vector<Dune::FieldVector<double, dim>>
792 adaptPointDimension_(std::vector<Dune::FieldVector<double, 3>>&& points3D)
const
793 {
return std::move(points3D); }
795 std::string fileName_;
796 tinyxml2::XMLDocument doc_;
A vtk file reader using tinyxml2 as xml backend.
Definition: vtkreader.hh:342
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:471
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:379
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:402
VTKReader(const std::string &fileName)
The constructor creates a tinyxml2::XMLDocument from file.
Definition: vtkreader.hh:355
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:448
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:423
DataType
The data array types.
Definition: vtkreader.hh:347
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:350
Free functions to write and read a sequence container to and from a file.
Definition: vtkreader.hh:274
const tinyxml2::XMLElement * getPieceNode(const tinyxml2::XMLDocument &doc, const std::string &fileName)
Get the piece node an xml document.
Definition: vtkreader.hh:282
std::string getProcessPieceFileName(const std::string &pvtkFileName)
get the vtk filename for the current processor
Definition: vtkreader.hh:306
constexpr Line line
Definition: couplingmanager1d3d_line.hh:31
Definition: common/pdesolver.hh:24