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());
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);
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)
584 factory.insertElement(Dune::GeometryTypes::line,
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_);
601 using namespace tinyxml2;
603 const XMLElement* pieceNode = getPieceNode_();
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);
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_);
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;
730 dataNode = pieceNode->FirstChildElement(
"PointData");
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());
781 Dune::GeometryType vtkToDuneGeomType_(
unsigned int vtkCellType)
const
785 case Dune::VTK::GeometryType::vertex:
return Dune::GeometryTypes::vertex;
786 case Dune::VTK::GeometryType::line:
return Dune::GeometryTypes::line;
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_;