277 using Data = std::unordered_map<std::string, std::vector<double>>;
284 using namespace tinyxml2;
285 const auto ext = std::filesystem::path(fileName).extension().string();
289 fileName_ = Dune::MPIHelper::instance().size() > 1 && ext[1] ==
'p' ?
292 const auto eResult = doc_.LoadFile(fileName_.c_str());
293 if (eResult != tinyxml2::XML_SUCCESS)
294 DUNE_THROW(Dune::IOError,
"Couldn't open XML file " << fileName_ <<
".");
296 const XMLElement* pieceNode = getPieceNode_();
297 if (pieceNode ==
nullptr)
298 DUNE_THROW(Dune::IOError,
"Couldn't get 'Piece' node in " << fileName_ <<
".");
308 using namespace tinyxml2;
310 const XMLElement* pieceNode = getPieceNode_();
311 const XMLElement* dataNode = getDataNode_(pieceNode, type);
312 if (dataNode ==
nullptr)
315 const XMLElement* dataArray = findDataArray_(dataNode, name);
316 if (dataArray ==
nullptr)
328 template<
class Container>
331 using namespace tinyxml2;
333 const XMLElement* pieceNode = getPieceNode_();
334 const XMLElement* dataNode = getDataNode_(pieceNode, type);
335 if (dataNode ==
nullptr)
336 DUNE_THROW(Dune::IOError,
"Couldn't get 'PointData' or 'CellData' node in " << fileName_ <<
".");
338 const XMLElement* dataArray = findDataArray_(dataNode, name);
339 if (dataArray ==
nullptr)
340 DUNE_THROW(Dune::IOError,
"Couldn't find the data array " << name <<
".");
342 return parseDataArray_<Container>(dataArray);
350 std::unique_ptr<Grid>
readGrid(
bool verbose =
false)
const
352 static_assert(!Dune::Capabilities::isCartesian<Grid>::v,
"Grid reader only supports unstructured grid implementations");
355 Dune::GridFactory<Grid> factory;
358 if (Dune::MPIHelper::instance().rank() == 0)
361 std::cout <<
"Reading " << Grid::dimension <<
"d grid from vtk file " << fileName_ <<
"." << std::endl;
362 readGrid_(factory, verbose);
365 return std::unique_ptr<Grid>(factory.createGrid());
375 std::unique_ptr<Grid>
readGrid(Dune::GridFactory<Grid>& factory,
bool verbose =
false)
const
377 static_assert(!Dune::Capabilities::isCartesian<Grid>::v,
"Grid reader only supports unstructured grid implementations");
379 if (Dune::MPIHelper::instance().rank() == 0)
382 std::cout <<
"Reading " << Grid::dimension <<
"d grid from vtk file " << fileName_ <<
"." << std::endl;
383 readGrid_(factory, verbose);
386 return std::unique_ptr<Grid>(factory.createGrid());
400 static_assert(!Dune::Capabilities::isCartesian<Grid>::v,
"Grid reader only supports unstructured grid implementations");
402 if (Dune::MPIHelper::instance().rank() == 0)
405 std::cout <<
"Reading " << Grid::dimension <<
"d grid from vtk file " << fileName_ <<
"." << std::endl;
406 readGrid_(factory, verbose);
410 return std::unique_ptr<Grid>(factory.createGrid());
420 void readGrid_(Dune::GridFactory<Grid>& factory,
bool verbose =
false)
const
422 using namespace tinyxml2;
424 const XMLElement* pieceNode = getPieceNode_();
425 const XMLElement* pointsNode = pieceNode->FirstChildElement(
"Points")->FirstChildElement(
"DataArray");
426 if (pointsNode ==
nullptr)
427 DUNE_THROW(Dune::IOError,
"Couldn't get data array of points in " << fileName_ <<
".");
429 using Point3D = Dune::FieldVector<double, 3>;
430 std::vector<Point3D> points3D;
431 std::stringstream dataStream(pointsNode->GetText());
432 std::istream_iterator<Point3D> it(dataStream);
433 std::copy(it, std::istream_iterator<Point3D>(), std::back_inserter(points3D));
436 auto points = adaptPointDimension_<Grid::dimensionworld>(std::move(points3D));
438 if (verbose) std::cout <<
"Found " << points.size() <<
" vertices." << std::endl;
441 for (
auto&& point : points)
442 factory.insertVertex(std::move(point));
444 const XMLElement* cellsNode = pieceNode->FirstChildElement(
"Cells");
445 const XMLElement* linesNode = pieceNode->FirstChildElement(
"Lines");
448 const XMLElement* connectivityNode = findDataArray_(cellsNode,
"connectivity");
449 const XMLElement* offsetsNode = findDataArray_(cellsNode,
"offsets");
450 const XMLElement* typesNode = findDataArray_(cellsNode,
"types");
452 const auto connectivity = parseDataArray_<std::vector<unsigned int>>(connectivityNode);
453 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
454 const auto types = parseDataArray_<std::vector<unsigned int>>(typesNode);
456 if (verbose) std::cout <<
"Found " << offsets.size() <<
" elements." << std::endl;
458 unsigned int lastOffset = 0;
459 for (
unsigned int i = 0; i < offsets.size(); ++i)
461 const auto geomType = vtkToDuneGeomType_(types[i]);
462 unsigned int offset = offsets[i];
463 std::vector<unsigned int> corners; corners.resize(offset-lastOffset);
464 for (
unsigned int j = 0; j < offset-lastOffset; ++j)
465 corners[Dune::VTK::renumber(geomType, j)] = connectivity[lastOffset+j];
466 factory.insertElement(geomType, std::move(corners));
474 if (Grid::dimension != 1)
475 DUNE_THROW(Dune::IOError,
"Grid expects dimension " << Grid::dimension
476 <<
" but " << fileName_ <<
" contains a 1D grid.");
478 const XMLElement* connectivityNode = findDataArray_(linesNode,
"connectivity");
479 const XMLElement* offsetsNode = findDataArray_(linesNode,
"offsets");
481 const auto connectivity = parseDataArray_<std::vector<unsigned int>>(connectivityNode);
482 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
484 if (verbose) std::cout <<
"Found " << offsets.size() <<
" polylines." << std::endl;
486 unsigned int lastOffset = 0;
487 for (
unsigned int i = 0; i < offsets.size(); ++i)
491 unsigned int offset = offsets[i];
492 for (
unsigned int j = 0; j < offset-lastOffset-1; ++j)
493 factory.insertElement(Dune::GeometryTypes::line,
494 std::vector<unsigned int>({ connectivity[lastOffset+j], connectivity[lastOffset+j+1] }));
499 DUNE_THROW(Dune::IOError,
"No Cells or Lines element found in " << fileName_);
510 using namespace tinyxml2;
512 const XMLElement* pieceNode = getPieceNode_();
514 if (cellDataNode !=
nullptr)
516 const XMLElement* cellsNode = pieceNode->FirstChildElement(
"Cells");
517 const XMLElement* linesNode = pieceNode->FirstChildElement(
"Lines");
521 const XMLElement* dataArray = cellDataNode->FirstChildElement(
"DataArray");
522 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
524 const char *attributeText = dataArray->Attribute(
"Name");
526 if (attributeText ==
nullptr)
527 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a cell data array.");
529 cellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
532 std::cout <<
"Read cell data field " << attributeText << std::endl;
539 const XMLElement* dataArray = cellDataNode->FirstChildElement(
"DataArray");
542 Data polyLineCellData;
543 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
545 const char *attributeText = dataArray->Attribute(
"Name");
547 if (attributeText ==
nullptr)
548 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a cell data array.");
550 polyLineCellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
553 std::cout <<
"Read cell data field " << attributeText << std::endl;
559 const XMLElement* offsetsNode = findDataArray_(linesNode,
"offsets");
560 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
562 if (offsets.size() != polyLineCellData.begin()->second.size())
563 DUNE_THROW(Dune::IOError,
"Expected the same number of cell data entries (is "
564 << polyLineCellData.begin()->second.size()
565 <<
") as polylines (" << offsets.size() <<
")!");
568 unsigned int lastOffset = 0;
569 std::size_t numCells = 0;
570 for (
unsigned int i = 0; i < offsets.size(); ++i)
572 unsigned int offset = offsets[i];
573 for (
unsigned int j = 0; j < offset-lastOffset-1; ++j)
579 for (
const auto& dArray : polyLineCellData)
581 cellData[dArray.first] = std::vector<double>(numCells);
583 const auto& pd = dArray.second;
586 std::size_t cellIdx = 0;
587 for (
unsigned int i = 0; i < offsets.size(); ++i)
589 unsigned int offset = offsets[i];
590 for (
unsigned int j = 0; j < offset-lastOffset-1; ++j)
591 cd[cellIdx++] = pd[i];
598 DUNE_THROW(Dune::IOError,
"No Cells or Lines element found in " << fileName_);
602 if (pointDataNode !=
nullptr)
604 const XMLElement* dataArray = pointDataNode->FirstChildElement(
"DataArray");
605 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
607 const char *attributeText = dataArray->Attribute(
"Name");
609 if (attributeText ==
nullptr)
610 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a point data array.");
612 pointData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
615 std::cout <<
"Read point data field " << attributeText << std::endl;
624 const tinyxml2::XMLElement* getPieceNode_()
const
633 const tinyxml2::XMLElement* getDataNode_(
const tinyxml2::XMLElement* pieceNode,
const DataType& type)
const
635 using namespace tinyxml2;
637 const XMLElement* dataNode =
nullptr;
639 dataNode = pieceNode->FirstChildElement(
"PointData");
641 dataNode = pieceNode->FirstChildElement(
"CellData");
643 DUNE_THROW(Dune::IOError,
"Only cell and point data are supported.");
654 const tinyxml2::XMLElement* findDataArray_(
const tinyxml2::XMLElement* dataNode,
const std::string& name)
const
656 using namespace tinyxml2;
659 const XMLElement* dataArray = dataNode->FirstChildElement(
"DataArray");
660 for (; dataArray !=
nullptr; dataArray = dataArray->NextSiblingElement(
"DataArray"))
662 const char *attributeText = dataArray->Attribute(
"Name");
664 if (attributeText ==
nullptr)
665 DUNE_THROW(Dune::IOError,
"Couldn't get Name attribute of a data array.");
667 if (attributeText == name)
679 template<
class Container>
680 Container parseDataArray_(
const tinyxml2::XMLElement* dataArray)
const
682 std::stringstream dataStream(dataArray->GetText());
690 Dune::GeometryType vtkToDuneGeomType_(
unsigned int vtkCellType)
const
694 case Dune::VTK::GeometryType::vertex:
return Dune::GeometryTypes::vertex;
695 case Dune::VTK::GeometryType::line:
return Dune::GeometryTypes::line;
696 case Dune::VTK::GeometryType::triangle:
return Dune::GeometryTypes::triangle;
697 case Dune::VTK::GeometryType::quadrilateral:
return Dune::GeometryTypes::quadrilateral;
698 case Dune::VTK::GeometryType::tetrahedron:
return Dune::GeometryTypes::tetrahedron;
699 case Dune::VTK::GeometryType::hexahedron:
return Dune::GeometryTypes::hexahedron;
700 case Dune::VTK::GeometryType::prism:
return Dune::GeometryTypes::prism;
701 case Dune::VTK::GeometryType::pyramid:
return Dune::GeometryTypes::pyramid;
702 default: DUNE_THROW(Dune::NotImplemented,
"VTK cell type " << vtkCellType);
706 template<
int dim, std::enable_if_t<(dim < 3),
int> = 0>
707 std::vector<Dune::FieldVector<
double, dim>>
708 adaptPo
intDimension_(std::vector<Dune::FieldVector<
double, 3>>&& po
ints3D) const
710 std::vector<Dune::FieldVector<
double, dim>> po
ints(po
ints3D.size());
711 for (std::
size_t i = 0; i < po
ints.size(); ++i)
712 for (
int j = 0; j < dim; ++j)
713 po
ints[i][j] = po
ints3D[i][j];
717 template<
int dim, std::enable_if_t<(dim == 3),
int> = 0>
718 std::vector<Dune::FieldVector<double, dim>>
719 adaptPointDimension_(std::vector<Dune::FieldVector<double, 3>>&& points3D)
const
720 {
return std::move(points3D); }
722 std::string fileName_;
723 tinyxml2::XMLDocument doc_;