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/parallel/traits.hpp>
36#include <gridformat/traits/dune.hpp>
37#include <gridformat/reader.hpp>
38#include <gridformat/decorators/reader_polylines_subdivider.hpp>
40#include "imagegrid_.hh"
45auto makeGridformatReaderFactory(
const C& c)
49 using Comm = std::remove_cvref_t<C>;
50 using GridFormatComm = std::conditional_t<
51 std::is_same_v<Comm, Dune::No_Comm>,
52 GridFormat::NullCommunicator,
56 const auto comm = [&]() -> GridFormatComm
58 if constexpr (std::is_same_v<GridFormatComm, GridFormat::NullCommunicator>)
64 return [fac = GridFormat::AnyReaderFactory<GridFormatComm>{comm}] (
const std::string& f)
65 -> std::unique_ptr<GridFormat::GridReader> {
67 if (f.ends_with(
"vtp"))
68 return std::make_unique<GridFormat::ReaderDecorators::PolylinesSubdivider>(fac.make_for(f));
69 return fac.make_for(f);
84 enum class DataType { cellData, pointData, fieldData };
87 using Data = std::unordered_map<std::string, std::vector<double>>;
89 explicit VTKReader(
const std::string& fileName)
90 : reader_(std::make_shared<GridFormat::Reader>(GridFormat::Reader::from(fileName,
91 Detail::
VTKReader::makeGridformatReaderFactory(
92 Dune::MPIHelper::instance().getCommunicator()
105 return std::ranges::any_of(
106 cell_field_names(*reader_),
107 [&] (
const auto& n) {
return name == n; }
110 return std::ranges::any_of(
111 point_field_names(*reader_),
112 [&] (
const auto& n) {
return name == n; }
114 else if (type == DataType::fieldData)
115 return std::ranges::any_of(
116 meta_data_field_names(*reader_),
117 [&] (
const auto& n) {
return name == n; }
120 DUNE_THROW(Dune::IOError,
"Unknown data type for field '" << name <<
"'");
129 template<
class Container>
134 Container values(reader_->number_of_cells());
135 reader_->cell_field(name)->export_to(values);
140 Container values(reader_->number_of_points());
141 reader_->point_field(name)->export_to(values);
144 else if (type == DataType::fieldData)
146 DUNE_THROW(Dune::NotImplemented,
"Reading meta data not yet implemented");
149 DUNE_THROW(Dune::IOError,
"Unknown data type for field '" << name <<
"'");
157 std::unique_ptr<Grid>
readGrid(
bool verbose =
false)
const
159 Dune::GridFactory<Grid> factory;
170 std::unique_ptr<Grid>
readGrid(Dune::GridFactory<Grid>& factory,
bool verbose =
false)
const
172 if (Dune::MPIHelper::instance().rank() == 0)
175 std::cout <<
"Reading " << Grid::dimension <<
"d grid from vtk file " << reader_->filename() <<
"." << std::endl;
178 GridFormat::Dune::GridFactoryAdapter<Grid> adapter{ factory };
179 reader_->export_grid(adapter);
183 return std::unique_ptr<Grid>(factory.createGrid());
195 std::unique_ptr<Grid>
readGrid(Dune::GridFactory<Grid>& factory,
Data& cellData,
Data& pointData,
bool verbose =
false)
const
197 auto grid =
readGrid(factory, verbose);
200 for (
const auto& name : cell_field_names(*reader_))
201 cellData[name] = Data::mapped_type{};
202 for (
const auto& name : point_field_names(*reader_))
203 pointData[name] = Data::mapped_type{};
206 if (Dune::MPIHelper::instance().rank() == 0)
208 for (
const auto& [name, field_ptr] : cell_fields(*reader_))
210 if (verbose) std::cout <<
"-- reading cell data '" << name <<
"'." << std::endl;
211 field_ptr->export_to(cellData[name]);
214 for (
const auto& [name, field_ptr] : point_fields(*reader_))
216 if (verbose) std::cout <<
"-- reading point data '" << name <<
"'." << std::endl;
217 field_ptr->export_to(pointData[name]);
221 return std::unique_ptr<Grid>(std::move(grid));
234 template<
class ct,
int dim>
235 std::unique_ptr<Detail::VTKReader::ImageGrid<ct, dim>> readStructuredGrid(
bool verbose =
false)
const
238 std::cout <<
"Reading " << dim <<
"d grid from vtk file " << reader_->filename() <<
"." << std::endl;
240 return std::make_unique<Detail::VTKReader::ImageGrid<ct, dim>>(reader_);
251 template<
class ct,
int dim>
252 std::unique_ptr<Detail::VTKReader::ImageGrid<ct, dim>> readStructuredGrid(
Data& cellData,
Data& pointData,
bool verbose =
false)
const
254 auto grid = readStructuredGrid<ct, dim>(verbose);
257 for (
const auto& name : cell_field_names(*reader_))
258 cellData[name] = Data::mapped_type{};
259 for (
const auto& name : point_field_names(*reader_))
260 pointData[name] = Data::mapped_type{};
263 if (Dune::MPIHelper::instance().rank() == 0)
265 for (
const auto& [name, field_ptr] : cell_fields(*reader_))
267 if (verbose) std::cout <<
"-- reading cell data '" << name <<
"'." << std::endl;
268 field_ptr->export_to(cellData[name]);
271 for (
const auto& [name, field_ptr] : point_fields(*reader_))
273 if (verbose) std::cout <<
"-- reading point data '" << name <<
"'." << std::endl;
274 field_ptr->export_to(pointData[name]);
277 return { std::move(grid) };
281 std::shared_ptr<GridFormat::Reader> reader_;
288#include <dumux/io/xml/tinyxml2.h>
300const tinyxml2::XMLElement*
getPieceNode(
const tinyxml2::XMLDocument& doc,
const std::string& fileName)
302 using namespace tinyxml2;
304 const XMLElement* pieceNode = doc.FirstChildElement(
"VTKFile");
305 if (pieceNode ==
nullptr)
306 DUNE_THROW(Dune::IOError,
"Couldn't get 'VTKFile' node in " << fileName <<
".");
308 pieceNode = pieceNode->FirstChildElement(
"UnstructuredGrid");
309 if (pieceNode ==
nullptr)
310 pieceNode = doc.FirstChildElement(
"VTKFile")->FirstChildElement(
"PolyData");
311 if (pieceNode ==
nullptr)
312 pieceNode = doc.FirstChildElement(
"VTKFile")->FirstChildElement(
"PUnstructuredGrid");
313 if (pieceNode ==
nullptr)
314 pieceNode = doc.FirstChildElement(
"VTKFile")->FirstChildElement(
"PPolyData");
315 if (pieceNode ==
nullptr)
316 DUNE_THROW(Dune::IOError,
"Couldn't get 'UnstructuredGrid', 'PUnstructuredGrid', 'PolyData', or 'PPolyData' node in " << fileName <<
".");
318 return pieceNode->FirstChildElement(
"Piece");
326 using namespace tinyxml2;
329 const auto eResult = pDoc.LoadFile(pvtkFileName.c_str());
330 if (eResult != XML_SUCCESS)
331 DUNE_THROW(Dune::IOError,
"Couldn't open XML file " << pvtkFileName <<
".");
334 const XMLElement* pieceNode =
getPieceNode(pDoc, pvtkFileName);
335 const auto myrank = Dune::MPIHelper::instance().rank();
336 for (
int rank = 0; rank < myrank; ++rank)
338 pieceNode = pieceNode->NextSiblingElement(
"Piece");
339 if (pieceNode ==
nullptr)
340 DUNE_THROW(Dune::IOError,
"Couldn't find 'Piece' node for rank "
341 << rank <<
" in " << pvtkFileName <<
".");
344 const char *vtkFileName = pieceNode->Attribute(
"Source");
345 if (vtkFileName ==
nullptr)
346 DUNE_THROW(Dune::IOError,
"Couldn't get 'Source' attribute of 'Piece' node no. " << myrank <<
" in " << pvtkFileName);
368 using Data = std::unordered_map<std::string, std::vector<double>>;
375 using namespace tinyxml2;
376 const auto ext = std::filesystem::path(fileName).extension().string();
380 fileName_ = Dune::MPIHelper::instance().size() > 1 && ext[1] ==
'p' ?
383 const auto eResult = doc_.LoadFile(fileName_.c_str());
384 if (eResult != tinyxml2::XML_SUCCESS)
385 DUNE_THROW(Dune::IOError,
"Couldn't open XML file " << fileName_ <<
".");
387 const XMLElement* pieceNode = getPieceNode_();
388 if (pieceNode ==
nullptr)
389 DUNE_THROW(Dune::IOError,
"Couldn't get 'Piece' node in " << fileName_ <<
".");
399 using namespace tinyxml2;
401 const XMLElement* pieceNode = getPieceNode_();
402 const XMLElement* dataNode = getDataNode_(pieceNode, type);
403 if (dataNode ==
nullptr)
406 const XMLElement* dataArray = findDataArray_(dataNode, name);
407 if (dataArray ==
nullptr)
419 template<
class Container>
422 using namespace tinyxml2;
424 const XMLElement* pieceNode = getPieceNode_();
425 const XMLElement* dataNode = getDataNode_(pieceNode, type);
426 if (dataNode ==
nullptr)
427 DUNE_THROW(Dune::IOError,
"Couldn't get 'PointData' or 'CellData' node in " << fileName_ <<
".");
429 const XMLElement* dataArray = findDataArray_(dataNode, name);
430 if (dataArray ==
nullptr)
431 DUNE_THROW(Dune::IOError,
"Couldn't find the data array " << name <<
".");
433 return parseDataArray_<Container>(dataArray);
441 std::unique_ptr<Grid>
readGrid(
bool verbose =
false)
const
443 static_assert(!Dune::Capabilities::isCartesian<Grid>::v,
"Grid reader only supports unstructured grid implementations");
446 Dune::GridFactory<Grid> factory;
449 if (Dune::MPIHelper::instance().rank() == 0)
452 std::cout <<
"Reading " << Grid::dimension <<
"d grid from vtk file " << fileName_ <<
"." << std::endl;
453 readGrid_(factory, verbose);
456 return std::unique_ptr<Grid>(factory.createGrid());
466 std::unique_ptr<Grid>
readGrid(Dune::GridFactory<Grid>& factory,
bool verbose =
false)
const
468 static_assert(!Dune::Capabilities::isCartesian<Grid>::v,
"Grid reader only supports unstructured grid implementations");
470 if (Dune::MPIHelper::instance().rank() == 0)
473 std::cout <<
"Reading " << Grid::dimension <<
"d grid from vtk file " << fileName_ <<
"." << std::endl;
474 readGrid_(factory, verbose);
477 return std::unique_ptr<Grid>(factory.createGrid());
489 std::unique_ptr<Grid>
readGrid(Dune::GridFactory<Grid>& factory,
Data& cellData,
Data& pointData,
bool verbose =
false)
const
491 static_assert(!Dune::Capabilities::isCartesian<Grid>::v,
"Grid reader only supports unstructured grid implementations");
493 if (Dune::MPIHelper::instance().rank() == 0)
496 std::cout <<
"Reading " << Grid::dimension <<
"d grid from vtk file " << fileName_ <<
"." << std::endl;
497 readGrid_(factory, verbose);
498 readGridData_(cellData, pointData, verbose);
501 return std::unique_ptr<Grid>(factory.createGrid());
511 void readGrid_(Dune::GridFactory<Grid>& factory,
bool verbose =
false)
const
513 using namespace tinyxml2;
515 const XMLElement* pieceNode = getPieceNode_();
516 const XMLElement* pointsNode = pieceNode->FirstChildElement(
"Points")->FirstChildElement(
"DataArray");
517 if (pointsNode ==
nullptr)
518 DUNE_THROW(Dune::IOError,
"Couldn't get data array of points in " << fileName_ <<
".");
520 using Point3D = Dune::FieldVector<double, 3>;
521 std::vector<Point3D> points3D;
522 std::stringstream dataStream(pointsNode->GetText());
523 std::istream_iterator<Point3D> it(dataStream);
524 std::copy(it, std::istream_iterator<Point3D>(), std::back_inserter(points3D));
527 auto points = adaptPointDimension_<Grid::dimensionworld>(std::move(points3D));
529 if (verbose) std::cout <<
"Found " << points.size() <<
" vertices." << std::endl;
532 for (
auto&& point : points)
533 factory.insertVertex(std::move(point));
535 const XMLElement* cellsNode = pieceNode->FirstChildElement(
"Cells");
536 const XMLElement* linesNode = pieceNode->FirstChildElement(
"Lines");
539 const XMLElement* connectivityNode = findDataArray_(cellsNode,
"connectivity");
540 const XMLElement* offsetsNode = findDataArray_(cellsNode,
"offsets");
541 const XMLElement* typesNode = findDataArray_(cellsNode,
"types");
543 const auto connectivity = parseDataArray_<std::vector<unsigned int>>(connectivityNode);
544 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
545 const auto types = parseDataArray_<std::vector<unsigned int>>(typesNode);
547 if (verbose) std::cout <<
"Found " << offsets.size() <<
" elements." << std::endl;
549 unsigned int lastOffset = 0;
550 for (
unsigned int i = 0; i < offsets.size(); ++i)
552 const auto geomType = vtkToDuneGeomType_(types[i]);
553 unsigned int offset = offsets[i];
554 std::vector<unsigned int> corners; corners.resize(offset-lastOffset);
555 for (
unsigned int j = 0; j < offset-lastOffset; ++j)
556 corners[Dune::VTK::renumber(geomType, j)] = connectivity[lastOffset+j];
557 factory.insertElement(geomType, std::move(corners));
565 if (Grid::dimension != 1)
566 DUNE_THROW(Dune::IOError,
"Grid expects dimension " << Grid::dimension
567 <<
" but " << fileName_ <<
" contains a 1D grid.");
569 const XMLElement* connectivityNode = findDataArray_(linesNode,
"connectivity");
570 const XMLElement* offsetsNode = findDataArray_(linesNode,
"offsets");
572 const auto connectivity = parseDataArray_<std::vector<unsigned int>>(connectivityNode);
573 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
575 if (verbose) std::cout <<
"Found " << offsets.size() <<
" polylines." << std::endl;
577 unsigned int lastOffset = 0;
578 for (
unsigned int i = 0; i < offsets.size(); ++i)
582 unsigned int offset = offsets[i];
583 for (
unsigned int j = 0; j < offset-lastOffset-1; ++j)
585 std::vector<unsigned int>({ connectivity[lastOffset+j], connectivity[lastOffset+j+1] }));
590 DUNE_THROW(Dune::IOError,
"No Cells or Lines element found in " << fileName_);
599 void readGridData_(Data& cellData, Data& pointData,
bool verbose =
false)
const
601 using namespace tinyxml2;
603 const XMLElement* pieceNode = getPieceNode_();
604 const XMLElement* cellDataNode = getDataNode_(pieceNode, DataType::cellData);
605 if (cellDataNode !=
nullptr)
607 const XMLElement* cellsNode = pieceNode->FirstChildElement(
"Cells");
608 const XMLElement* linesNode = pieceNode->FirstChildElement(
"Lines");
612 const XMLElement* dataArray = cellDataNode->FirstChildElement(
"DataArray");
613 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
615 const char *attributeText = dataArray->Attribute(
"Name");
617 if (attributeText ==
nullptr)
618 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a cell data array.");
620 cellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
623 std::cout <<
"Read cell data field " << attributeText << std::endl;
630 const XMLElement* dataArray = cellDataNode->FirstChildElement(
"DataArray");
633 Data polyLineCellData;
634 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
636 const char *attributeText = dataArray->Attribute(
"Name");
638 if (attributeText ==
nullptr)
639 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a cell data array.");
641 polyLineCellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
644 std::cout <<
"Read cell data field " << attributeText << std::endl;
650 const XMLElement* offsetsNode = findDataArray_(linesNode,
"offsets");
651 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
653 if (offsets.size() != polyLineCellData.begin()->second.size())
654 DUNE_THROW(Dune::IOError,
"Expected the same number of cell data entries (is "
655 << polyLineCellData.begin()->second.size()
656 <<
") as polylines (" << offsets.size() <<
")!");
659 unsigned int lastOffset = 0;
660 std::size_t numCells = 0;
661 for (
unsigned int i = 0; i < offsets.size(); ++i)
663 unsigned int offset = offsets[i];
664 for (
unsigned int j = 0; j < offset-lastOffset-1; ++j)
670 for (
const auto& dArray : polyLineCellData)
672 cellData[dArray.first] = std::vector<double>(numCells);
673 auto& cd = cellData[dArray.first];
674 const auto& pd = dArray.second;
677 std::size_t cellIdx = 0;
678 for (
unsigned int i = 0; i < offsets.size(); ++i)
680 unsigned int offset = offsets[i];
681 for (
unsigned int j = 0; j < offset-lastOffset-1; ++j)
682 cd[cellIdx++] = pd[i];
689 DUNE_THROW(Dune::IOError,
"No Cells or Lines element found in " << fileName_);
692 const XMLElement* pointDataNode = getDataNode_(pieceNode, DataType::pointData);
693 if (pointDataNode !=
nullptr)
695 const XMLElement* dataArray = pointDataNode->FirstChildElement(
"DataArray");
696 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
698 const char *attributeText = dataArray->Attribute(
"Name");
700 if (attributeText ==
nullptr)
701 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a point data array.");
703 pointData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
706 std::cout <<
"Read point data field " << attributeText << std::endl;
715 const tinyxml2::XMLElement* getPieceNode_()
const
724 const tinyxml2::XMLElement* getDataNode_(
const tinyxml2::XMLElement* pieceNode,
const DataType& type)
const
726 using namespace tinyxml2;
728 const XMLElement* dataNode =
nullptr;
729 if (type == DataType::pointData)
730 dataNode = pieceNode->FirstChildElement(
"PointData");
731 else if (type == DataType::cellData)
732 dataNode = pieceNode->FirstChildElement(
"CellData");
734 DUNE_THROW(Dune::IOError,
"Only cell and point data are supported.");
745 const tinyxml2::XMLElement* findDataArray_(
const tinyxml2::XMLElement* dataNode,
const std::string& name)
const
747 using namespace tinyxml2;
750 const XMLElement* dataArray = dataNode->FirstChildElement(
"DataArray");
751 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
753 const char *attributeText = dataArray->Attribute(
"Name");
755 if (attributeText ==
nullptr)
756 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a data array.");
758 if (attributeText == name)
770 template<
class Container>
771 Container parseDataArray_(
const tinyxml2::XMLElement* dataArray)
const
773 std::stringstream dataStream(dataArray->GetText());
774 return readStreamToContainer<Container>(dataStream);
781 Dune::GeometryType vtkToDuneGeomType_(
unsigned int vtkCellType)
const
785 case Dune::VTK::GeometryType::vertex:
return Dune::GeometryTypes::vertex;
787 case Dune::VTK::GeometryType::triangle:
return Dune::GeometryTypes::triangle;
788 case Dune::VTK::GeometryType::quadrilateral:
return Dune::GeometryTypes::quadrilateral;
789 case Dune::VTK::GeometryType::tetrahedron:
return Dune::GeometryTypes::tetrahedron;
790 case Dune::VTK::GeometryType::hexahedron:
return Dune::GeometryTypes::hexahedron;
791 case Dune::VTK::GeometryType::prism:
return Dune::GeometryTypes::prism;
792 case Dune::VTK::GeometryType::pyramid:
return Dune::GeometryTypes::pyramid;
793 default: DUNE_THROW(Dune::NotImplemented,
"VTK cell type " << vtkCellType);
797 template<
int dim, std::enable_if_t<(dim < 3),
int> = 0>
798 std::vector<Dune::FieldVector<
double, dim>>
799 adaptPo
intDimension_(std::vector<Dune::FieldVector<
double, 3>>&& po
ints3D) const
801 std::vector<Dune::FieldVector<
double, dim>> po
ints(po
ints3D.size());
802 for (std::
size_t i = 0; i < po
ints.size(); ++i)
803 for (
int j = 0; j < dim; ++j)
804 po
ints[i][j] = po
ints3D[i][j];
808 template<
int dim, std::enable_if_t<(dim == 3),
int> = 0>
809 std::vector<Dune::FieldVector<double, dim>>
810 adaptPointDimension_(std::vector<Dune::FieldVector<double, 3>>&& points3D)
const
811 {
return std::move(points3D); }
813 std::string fileName_;
814 tinyxml2::XMLDocument doc_;
A vtk file reader using tinyxml2 as xml backend.
Definition: vtkreader.hh:360
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:489
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:397
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:420
VTKReader(const std::string &fileName)
The constructor creates a tinyxml2::XMLDocument from file.
Definition: vtkreader.hh:373
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:466
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:441
DataType
The data array types.
Definition: vtkreader.hh:365
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:368
Free functions to write and read a sequence container to and from a file.
Definition: vtkreader.hh:292
const tinyxml2::XMLElement * getPieceNode(const tinyxml2::XMLDocument &doc, const std::string &fileName)
Get the piece node an xml document.
Definition: vtkreader.hh:300
std::string getProcessPieceFileName(const std::string &pvtkFileName)
get the vtk filename for the current processor
Definition: vtkreader.hh:324
constexpr Line line
Definition: couplingmanager1d3d_line.hh:31
Definition: common/pdesolver.hh:24