version 3.11-dev
vtkreader.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_IO_VTK_VTKREADER_HH
13#define DUMUX_IO_VTK_VTKREADER_HH
14
15#include <iostream>
16#include <iterator>
17#include <algorithm>
18#include <memory>
19#include <type_traits>
20#include <unordered_map>
21#include <utility>
22#include <filesystem>
23
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>
29
30#include <dumux/io/container.hh>
31
32#if DUMUX_HAVE_GRIDFORMAT
33
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>
39
40#include "imagegrid_.hh"
41
43
44template<class C>
45auto makeGridformatReaderFactory(const C& c)
46{
47 // Dune::No_Comm does not satisfy the GridFormat communicator concept
48 // we use GridFormat::NullCommunicator in this case
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,
53 Comm
54 >;
55
56 const auto comm = [&]() -> GridFormatComm
57 {
58 if constexpr (std::is_same_v<GridFormatComm, GridFormat::NullCommunicator>)
59 return {};
60 else
61 return c;
62 }();
63
64 return [fac = GridFormat::AnyReaderFactory<GridFormatComm>{comm}] (const std::string& f)
65 -> std::unique_ptr<GridFormat::GridReader> {
66 // use adapter for poly data that splits polylines into segments
67 if (f.ends_with("vtp"))
68 return std::make_unique<GridFormat::ReaderDecorators::PolylinesSubdivider>(fac.make_for(f));
69 return fac.make_for(f);
70 };
71}
72
73} // end namespace Dumux::Detail::VTKReader
74
75namespace Dumux {
76
81class VTKReader
82{
83public:
84 enum class DataType { cellData, pointData, fieldData };
85
87 using Data = std::unordered_map<std::string, std::vector<double>>;
88
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()
93 )
94 )))
95 {}
96
102 bool hasData(const std::string& name, const DataType& type) const
103 {
104 if (type == DataType::cellData)
105 return std::ranges::any_of(
106 cell_field_names(*reader_),
107 [&] (const auto& n) { return name == n; }
108 );
109 else if (type == DataType::pointData)
110 return std::ranges::any_of(
111 point_field_names(*reader_),
112 [&] (const auto& n) { return name == n; }
113 );
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; }
118 );
119 else
120 DUNE_THROW(Dune::IOError, "Unknown data type for field '" << name << "'");
121 }
122
129 template<class Container>
130 Container readData(const std::string& name, const DataType& type) const
131 {
132 if (type == DataType::cellData)
133 {
134 Container values(reader_->number_of_cells());
135 reader_->cell_field(name)->export_to(values);
136 return values;
137 }
138 else if (type == DataType::pointData)
139 {
140 Container values(reader_->number_of_points());
141 reader_->point_field(name)->export_to(values);
142 return values;
143 }
144 else if (type == DataType::fieldData)
145 {
146 DUNE_THROW(Dune::NotImplemented, "Reading meta data not yet implemented");
147 }
148 else
149 DUNE_THROW(Dune::IOError, "Unknown data type for field '" << name << "'");
150 }
151
156 template<class Grid>
157 std::unique_ptr<Grid> readGrid(bool verbose = false) const
158 {
159 Dune::GridFactory<Grid> factory;
160 return readGrid(factory, verbose);
161 }
162
169 template<class Grid>
170 std::unique_ptr<Grid> readGrid(Dune::GridFactory<Grid>& factory, bool verbose = false) const
171 {
172 if (Dune::MPIHelper::instance().rank() == 0)
173 {
174 if (verbose)
175 std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << reader_->filename() << "." << std::endl;
176
177 {
178 GridFormat::Dune::GridFactoryAdapter<Grid> adapter{ factory };
179 reader_->export_grid(adapter);
180 }
181 }
182
183 return std::unique_ptr<Grid>(factory.createGrid());
184 }
185
194 template<class Grid>
195 std::unique_ptr<Grid> readGrid(Dune::GridFactory<Grid>& factory, Data& cellData, Data& pointData, bool verbose = false) const
196 {
197 auto grid = readGrid(factory, verbose);
198
199 // read field names on all processes
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{};
204
205 // read data arrays only on rank 0
206 if (Dune::MPIHelper::instance().rank() == 0)
207 {
208 for (const auto& [name, field_ptr] : cell_fields(*reader_))
209 {
210 if (verbose) std::cout << "-- reading cell data '" << name << "'." << std::endl;
211 field_ptr->export_to(cellData[name]);
212 }
213
214 for (const auto& [name, field_ptr] : point_fields(*reader_))
215 {
216 if (verbose) std::cout << "-- reading point data '" << name << "'." << std::endl;
217 field_ptr->export_to(pointData[name]);
218 }
219 }
220
221 return std::unique_ptr<Grid>(std::move(grid));
222 }
223
234 template<class ct, int dim>
235 std::unique_ptr<Detail::VTKReader::ImageGrid<ct, dim>> readStructuredGrid(bool verbose = false) const
236 {
237 if (verbose)
238 std::cout << "Reading " << dim << "d grid from vtk file " << reader_->filename() << "." << std::endl;
239
240 return std::make_unique<Detail::VTKReader::ImageGrid<ct, dim>>(reader_);
241 }
242
251 template<class ct, int dim>
252 std::unique_ptr<Detail::VTKReader::ImageGrid<ct, dim>> readStructuredGrid(Data& cellData, Data& pointData, bool verbose = false) const
253 {
254 auto grid = readStructuredGrid<ct, dim>(verbose);
255
256 // read field names on all processes
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{};
261
262 // read data arrays only on rank 0
263 if (Dune::MPIHelper::instance().rank() == 0)
264 {
265 for (const auto& [name, field_ptr] : cell_fields(*reader_))
266 {
267 if (verbose) std::cout << "-- reading cell data '" << name << "'." << std::endl;
268 field_ptr->export_to(cellData[name]);
269 }
270
271 for (const auto& [name, field_ptr] : point_fields(*reader_))
272 {
273 if (verbose) std::cout << "-- reading point data '" << name << "'." << std::endl;
274 field_ptr->export_to(pointData[name]);
275 }
276 }
277 return { std::move(grid) };
278 }
279
280private:
281 std::shared_ptr<GridFormat::Reader> reader_;
282};
283
284} // namespace Dumux
285
286#else // DUMUX_HAVE_GRIDFORMAT
287
288#include <dumux/io/xml/tinyxml2.h>
289// fallback to simple vtk reader using tinyxml2
290// this reader only support ascii files and limited VTK file formats
291
293
300const tinyxml2::XMLElement* getPieceNode(const tinyxml2::XMLDocument& doc, const std::string& fileName)
301{
302 using namespace tinyxml2;
303
304 const XMLElement* pieceNode = doc.FirstChildElement("VTKFile");
305 if (pieceNode == nullptr)
306 DUNE_THROW(Dune::IOError, "Couldn't get 'VTKFile' node in " << fileName << ".");
307
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 << ".");
317
318 return pieceNode->FirstChildElement("Piece");
319}
320
324std::string getProcessPieceFileName(const std::string& pvtkFileName)
325{
326 using namespace tinyxml2;
327
328 XMLDocument pDoc;
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 << ".");
332
333 // get the first piece node
334 const XMLElement* pieceNode = getPieceNode(pDoc, pvtkFileName);
335 const auto myrank = Dune::MPIHelper::instance().rank();
336 for (int rank = 0; rank < myrank; ++rank)
337 {
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 << ".");
342 }
343
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);
347
348 return vtkFileName;
349}
350
351} // end namespace Dumux::Detail::VTKReader
352
353namespace Dumux {
354
360{
361public:
365 enum class DataType { cellData, pointData };
366
368 using Data = std::unordered_map<std::string, std::vector<double>>;
369
373 explicit VTKReader(const std::string& fileName)
374 {
375 using namespace tinyxml2;
376 const auto ext = std::filesystem::path(fileName).extension().string();
377 // If in parallel and the file to read is a parallel piece collection (pvtu/pvtp)
378 // read only the piece belonging to the own process. For this to work, the files
379 // have to have exactly the same amount of pieces than processes.
380 fileName_ = Dune::MPIHelper::instance().size() > 1 && ext[1] == 'p' ?
382
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_ << ".");
386
387 const XMLElement* pieceNode = getPieceNode_();
388 if (pieceNode == nullptr)
389 DUNE_THROW(Dune::IOError, "Couldn't get 'Piece' node in " << fileName_ << ".");
390 }
391
397 bool hasData(const std::string& name, const DataType& type) const
398 {
399 using namespace tinyxml2;
400
401 const XMLElement* pieceNode = getPieceNode_();
402 const XMLElement* dataNode = getDataNode_(pieceNode, type);
403 if (dataNode == nullptr)
404 return false;
405
406 const XMLElement* dataArray = findDataArray_(dataNode, name);
407 if (dataArray == nullptr)
408 return false;
409
410 return true;
411 }
412
419 template<class Container>
420 Container readData(const std::string& name, const DataType& type) const
421 {
422 using namespace tinyxml2;
423
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_ << ".");
428
429 const XMLElement* dataArray = findDataArray_(dataNode, name);
430 if (dataArray == nullptr)
431 DUNE_THROW(Dune::IOError, "Couldn't find the data array " << name << ".");
432
433 return parseDataArray_<Container>(dataArray);
434 }
435
440 template<class Grid>
441 std::unique_ptr<Grid> readGrid(bool verbose = false) const
442 {
443 static_assert(!Dune::Capabilities::isCartesian<Grid>::v, "Grid reader only supports unstructured grid implementations");
444
445 // make a grid factory
446 Dune::GridFactory<Grid> factory;
447
448 // only read on rank 0
449 if (Dune::MPIHelper::instance().rank() == 0)
450 {
451 if (verbose)
452 std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << fileName_ << "." << std::endl;
453 readGrid_(factory, verbose);
454 }
455
456 return std::unique_ptr<Grid>(factory.createGrid());
457 }
458
465 template<class Grid>
466 std::unique_ptr<Grid> readGrid(Dune::GridFactory<Grid>& factory, bool verbose = false) const
467 {
468 static_assert(!Dune::Capabilities::isCartesian<Grid>::v, "Grid reader only supports unstructured grid implementations");
469
470 if (Dune::MPIHelper::instance().rank() == 0)
471 {
472 if (verbose)
473 std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << fileName_ << "." << std::endl;
474 readGrid_(factory, verbose);
475 }
476
477 return std::unique_ptr<Grid>(factory.createGrid());
478 }
479
488 template<class Grid>
489 std::unique_ptr<Grid> readGrid(Dune::GridFactory<Grid>& factory, Data& cellData, Data& pointData, bool verbose = false) const
490 {
491 static_assert(!Dune::Capabilities::isCartesian<Grid>::v, "Grid reader only supports unstructured grid implementations");
492
493 if (Dune::MPIHelper::instance().rank() == 0)
494 {
495 if (verbose)
496 std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << fileName_ << "." << std::endl;
497 readGrid_(factory, verbose);
498 readGridData_(cellData, pointData, verbose);
499 }
500
501 return std::unique_ptr<Grid>(factory.createGrid());
502 }
503
504private:
510 template<class Grid>
511 void readGrid_(Dune::GridFactory<Grid>& factory, bool verbose = false) const
512 {
513 using namespace tinyxml2;
514
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_ << ".");
519
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));
525
526 // adapt point dimensions if grid dimension is smaller than 3
527 auto points = adaptPointDimension_<Grid::dimensionworld>(std::move(points3D));
528
529 if (verbose) std::cout << "Found " << points.size() << " vertices." << std::endl;
530
531 // insert vertices to the grid factory
532 for (auto&& point : points)
533 factory.insertVertex(std::move(point));
534
535 const XMLElement* cellsNode = pieceNode->FirstChildElement("Cells");
536 const XMLElement* linesNode = pieceNode->FirstChildElement("Lines");
537 if (cellsNode)
538 {
539 const XMLElement* connectivityNode = findDataArray_(cellsNode, "connectivity");
540 const XMLElement* offsetsNode = findDataArray_(cellsNode, "offsets");
541 const XMLElement* typesNode = findDataArray_(cellsNode, "types");
542
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);
546
547 if (verbose) std::cout << "Found " << offsets.size() << " elements." << std::endl;
548
549 unsigned int lastOffset = 0;
550 for (unsigned int i = 0; i < offsets.size(); ++i)
551 {
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));
558 lastOffset = offset;
559 }
560 }
561 // for poly data
562 else if (linesNode)
563 {
564 // sanity check
565 if (Grid::dimension != 1)
566 DUNE_THROW(Dune::IOError, "Grid expects dimension " << Grid::dimension
567 << " but " << fileName_ << " contains a 1D grid.");
568
569 const XMLElement* connectivityNode = findDataArray_(linesNode, "connectivity");
570 const XMLElement* offsetsNode = findDataArray_(linesNode, "offsets");
571
572 const auto connectivity = parseDataArray_<std::vector<unsigned int>>(connectivityNode);
573 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
574
575 if (verbose) std::cout << "Found " << offsets.size() << " polylines." << std::endl;
576
577 unsigned int lastOffset = 0;
578 for (unsigned int i = 0; i < offsets.size(); ++i)
579 {
580 // a polyline can have many points in the VTK format
581 // split the line in segments with two points
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] }));
586 lastOffset = offset;
587 }
588 }
589 else
590 DUNE_THROW(Dune::IOError, "No Cells or Lines element found in " << fileName_);
591 }
592
599 void readGridData_(Data& cellData, Data& pointData, bool verbose = false) const
600 {
601 using namespace tinyxml2;
602
603 const XMLElement* pieceNode = getPieceNode_();
604 const XMLElement* cellDataNode = getDataNode_(pieceNode, DataType::cellData);
605 if (cellDataNode != nullptr)
606 {
607 const XMLElement* cellsNode = pieceNode->FirstChildElement("Cells");
608 const XMLElement* linesNode = pieceNode->FirstChildElement("Lines");
609
610 if (cellsNode)
611 {
612 const XMLElement* dataArray = cellDataNode->FirstChildElement("DataArray");
613 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
614 {
615 const char *attributeText = dataArray->Attribute("Name");
616
617 if (attributeText == nullptr)
618 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a cell data array.");
619
620 cellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
621
622 if (verbose)
623 std::cout << "Read cell data field " << attributeText << std::endl;
624 }
625 }
626 // for poly data
627 else if (linesNode)
628 {
629 // first parse all the cell data (each cell in this sense can be a polyline)
630 const XMLElement* dataArray = cellDataNode->FirstChildElement("DataArray");
631 if (dataArray)
632 {
633 Data polyLineCellData;
634 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
635 {
636 const char *attributeText = dataArray->Attribute("Name");
637
638 if (attributeText == nullptr)
639 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a cell data array.");
640
641 polyLineCellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
642
643 if (verbose)
644 std::cout << "Read cell data field " << attributeText << std::endl;
645 }
646
647 // a polyline can have many points in the VTK format
648 // we split the line in segments with two points
649 // so we also need to duplicate the cell data to fit the increased line number
650 const XMLElement* offsetsNode = findDataArray_(linesNode, "offsets");
651 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
652
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() << ")!");
657
658 // count the number of Dune cells to be able to resize the data vectors
659 unsigned int lastOffset = 0;
660 std::size_t numCells = 0;
661 for (unsigned int i = 0; i < offsets.size(); ++i)
662 {
663 unsigned int offset = offsets[i];
664 for (unsigned int j = 0; j < offset-lastOffset-1; ++j)
665 ++numCells;
666 lastOffset = offset;
667 }
668
669 // create the data arrays
670 for (const auto& dArray : polyLineCellData)
671 {
672 cellData[dArray.first] = std::vector<double>(numCells);
673 auto& cd = cellData[dArray.first];
674 const auto& pd = dArray.second;
675
676 lastOffset = 0;
677 std::size_t cellIdx = 0;
678 for (unsigned int i = 0; i < offsets.size(); ++i)
679 {
680 unsigned int offset = offsets[i];
681 for (unsigned int j = 0; j < offset-lastOffset-1; ++j)
682 cd[cellIdx++] = pd[i];
683 lastOffset = offset;
684 }
685 }
686 }
687 }
688 else
689 DUNE_THROW(Dune::IOError, "No Cells or Lines element found in " << fileName_);
690 }
691
692 const XMLElement* pointDataNode = getDataNode_(pieceNode, DataType::pointData);
693 if (pointDataNode != nullptr)
694 {
695 const XMLElement* dataArray = pointDataNode->FirstChildElement("DataArray");
696 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
697 {
698 const char *attributeText = dataArray->Attribute("Name");
699
700 if (attributeText == nullptr)
701 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a point data array.");
702
703 pointData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
704
705 if (verbose)
706 std::cout << "Read point data field " << attributeText << std::endl;
707 }
708 }
709 }
710
715 const tinyxml2::XMLElement* getPieceNode_() const
716 { return Detail::VTKReader::getPieceNode(doc_, fileName_); }
717
724 const tinyxml2::XMLElement* getDataNode_(const tinyxml2::XMLElement* pieceNode, const DataType& type) const
725 {
726 using namespace tinyxml2;
727
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");
733 else
734 DUNE_THROW(Dune::IOError, "Only cell and point data are supported.");
735
736 return dataNode;
737 }
738
745 const tinyxml2::XMLElement* findDataArray_(const tinyxml2::XMLElement* dataNode, const std::string& name) const
746 {
747 using namespace tinyxml2;
748
749 // loop over XML node siblings to find the correct data array
750 const XMLElement* dataArray = dataNode->FirstChildElement("DataArray");
751 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
752 {
753 const char *attributeText = dataArray->Attribute("Name");
754
755 if (attributeText == nullptr)
756 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a data array.");
757
758 if (attributeText == name)
759 break;
760 }
761
762 return dataArray;
763 }
764
770 template<class Container>
771 Container parseDataArray_(const tinyxml2::XMLElement* dataArray) const
772 {
773 std::stringstream dataStream(dataArray->GetText());
774 return readStreamToContainer<Container>(dataStream);
775 }
776
781 Dune::GeometryType vtkToDuneGeomType_(unsigned int vtkCellType) const
782 {
783 switch (vtkCellType)
784 {
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);
794 }
795 }
796
797 template<int dim, std::enable_if_t<(dim < 3), int> = 0>
798 std::vector<Dune::FieldVector<double, dim>>
799 adaptPointDimension_(std::vector<Dune::FieldVector<double, 3>>&& points3D) const
800 {
801 std::vector<Dune::FieldVector<double, dim>> points(points3D.size());
802 for (std::size_t i = 0; i < points.size(); ++i)
803 for (int j = 0; j < dim; ++j)
804 points[i][j] = points3D[i][j];
805 return points;
806 }
807
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); }
812
813 std::string fileName_;
814 tinyxml2::XMLDocument doc_;
815};
816
817} // end namespace Dumux
818
819#endif // DUMUX_HAVE_GRIDFORMAT
820
821#endif
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: adapt.hh:17
Definition: common/pdesolver.hh:24