version 3.11-dev
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages Concepts
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/traits/dune.hpp>
36#include <gridformat/reader.hpp>
37#include <gridformat/decorators/reader_polylines_subdivider.hpp>
38
39#include "imagegrid_.hh"
40
42
43template<GridFormat::Concepts::Communicator C>
44auto makeGridformatReaderFactory(const C& c)
45{
46 return [fac = GridFormat::AnyReaderFactory<C>{c}] (const std::string& f)
47 -> std::unique_ptr<GridFormat::GridReader> {
48 // use adapter for poly data that splits polylines into segments
49 if (f.ends_with("vtp"))
50 return std::make_unique<GridFormat::ReaderDecorators::PolylinesSubdivider>(fac.make_for(f));
51 return fac.make_for(f);
52 };
53}
54
55} // end namespace Dumux::Detail::VTKReader
56
57namespace Dumux {
58
63class VTKReader
64{
65public:
66 enum class DataType { cellData, pointData, fieldData };
67
69 using Data = std::unordered_map<std::string, std::vector<double>>;
70
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()
75 )
76 )))
77 {}
78
84 bool hasData(const std::string& name, const DataType& type) const
85 {
86 if (type == DataType::cellData)
87 return std::ranges::any_of(
88 cell_field_names(*reader_),
89 [&] (const auto& n) { return name == n; }
90 );
91 else if (type == DataType::pointData)
92 return std::ranges::any_of(
93 point_field_names(*reader_),
94 [&] (const auto& n) { return name == n; }
95 );
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; }
100 );
101 else
102 DUNE_THROW(Dune::IOError, "Unknown data type for field '" << name << "'");
103 }
104
111 template<class Container>
112 Container readData(const std::string& name, const DataType& type) const
113 {
114 if (type == DataType::cellData)
115 {
116 Container values(reader_->number_of_cells());
117 reader_->cell_field(name)->export_to(values);
118 return values;
119 }
120 else if (type == DataType::pointData)
121 {
122 Container values(reader_->number_of_points());
123 reader_->point_field(name)->export_to(values);
124 return values;
125 }
126 else if (type == DataType::fieldData)
127 {
128 DUNE_THROW(Dune::NotImplemented, "Reading meta data not yet implemented");
129 }
130 else
131 DUNE_THROW(Dune::IOError, "Unknown data type for field '" << name << "'");
132 }
133
138 template<class Grid>
139 std::unique_ptr<Grid> readGrid(bool verbose = false) const
140 {
141 Dune::GridFactory<Grid> factory;
142 return readGrid(factory, verbose);
143 }
144
151 template<class Grid>
152 std::unique_ptr<Grid> readGrid(Dune::GridFactory<Grid>& factory, bool verbose = false) const
153 {
154 if (Dune::MPIHelper::instance().rank() == 0)
155 {
156 if (verbose)
157 std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << reader_->filename() << "." << std::endl;
158
159 {
160 GridFormat::Dune::GridFactoryAdapter<Grid> adapter{ factory };
161 reader_->export_grid(adapter);
162 }
163 }
164
165 return std::unique_ptr<Grid>(factory.createGrid());
166 }
167
176 template<class Grid>
177 std::unique_ptr<Grid> readGrid(Dune::GridFactory<Grid>& factory, Data& cellData, Data& pointData, bool verbose = false) const
178 {
179 auto grid = readGrid(factory, verbose);
180
181 // read field names on all processes
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{};
186
187 // read data arrays only on rank 0
188 if (Dune::MPIHelper::instance().rank() == 0)
189 {
190 for (const auto& [name, field_ptr] : cell_fields(*reader_))
191 {
192 if (verbose) std::cout << "-- reading cell data '" << name << "'." << std::endl;
193 field_ptr->export_to(cellData[name]);
194 }
195
196 for (const auto& [name, field_ptr] : point_fields(*reader_))
197 {
198 if (verbose) std::cout << "-- reading point data '" << name << "'." << std::endl;
199 field_ptr->export_to(pointData[name]);
200 }
201 }
202
203 return std::unique_ptr<Grid>(std::move(grid));
204 }
205
216 template<class ct, int dim>
217 std::unique_ptr<Detail::VTKReader::ImageGrid<ct, dim>> readStructuredGrid(bool verbose = false) const
218 {
219 if (verbose)
220 std::cout << "Reading " << dim << "d grid from vtk file " << reader_->filename() << "." << std::endl;
221
222 return std::make_unique<Detail::VTKReader::ImageGrid<ct, dim>>(reader_);
223 }
224
233 template<class ct, int dim>
234 std::unique_ptr<Detail::VTKReader::ImageGrid<ct, dim>> readStructuredGrid(Data& cellData, Data& pointData, bool verbose = false) const
235 {
236 auto grid = readStructuredGrid<ct, dim>(verbose);
237
238 // read field names on all processes
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{};
243
244 // read data arrays only on rank 0
245 if (Dune::MPIHelper::instance().rank() == 0)
246 {
247 for (const auto& [name, field_ptr] : cell_fields(*reader_))
248 {
249 if (verbose) std::cout << "-- reading cell data '" << name << "'." << std::endl;
250 field_ptr->export_to(cellData[name]);
251 }
252
253 for (const auto& [name, field_ptr] : point_fields(*reader_))
254 {
255 if (verbose) std::cout << "-- reading point data '" << name << "'." << std::endl;
256 field_ptr->export_to(pointData[name]);
257 }
258 }
259 return { std::move(grid) };
260 }
261
262private:
263 std::shared_ptr<GridFormat::Reader> reader_;
264};
265
266} // namespace Dumux
267
268#else // DUMUX_HAVE_GRIDFORMAT
269
270#include <dumux/io/xml/tinyxml2.h>
271// fallback to simple vtk reader using tinyxml2
272// this reader only support ascii files and limited VTK file formats
273
275
282const tinyxml2::XMLElement* getPieceNode(const tinyxml2::XMLDocument& doc, const std::string& fileName)
283{
284 using namespace tinyxml2;
285
286 const XMLElement* pieceNode = doc.FirstChildElement("VTKFile");
287 if (pieceNode == nullptr)
288 DUNE_THROW(Dune::IOError, "Couldn't get 'VTKFile' node in " << fileName << ".");
289
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 << ".");
299
300 return pieceNode->FirstChildElement("Piece");
301}
302
306std::string getProcessPieceFileName(const std::string& pvtkFileName)
307{
308 using namespace tinyxml2;
309
310 XMLDocument pDoc;
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 << ".");
314
315 // get the first piece node
316 const XMLElement* pieceNode = getPieceNode(pDoc, pvtkFileName);
317 const auto myrank = Dune::MPIHelper::instance().rank();
318 for (int rank = 0; rank < myrank; ++rank)
319 {
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 << ".");
324 }
325
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);
329
330 return vtkFileName;
331}
332
333} // end namespace Dumux::Detail::VTKReader
334
335namespace Dumux {
336
342{
343public:
347 enum class DataType { cellData, pointData };
348
350 using Data = std::unordered_map<std::string, std::vector<double>>;
351
355 explicit VTKReader(const std::string& fileName)
356 {
357 using namespace tinyxml2;
358 const auto ext = std::filesystem::path(fileName).extension().string();
359 // If in parallel and the file to read is a parallel piece collection (pvtu/pvtp)
360 // read only the piece belonging to the own process. For this to work, the files
361 // have to have exactly the same amount of pieces than processes.
362 fileName_ = Dune::MPIHelper::instance().size() > 1 && ext[1] == 'p' ?
364
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_ << ".");
368
369 const XMLElement* pieceNode = getPieceNode_();
370 if (pieceNode == nullptr)
371 DUNE_THROW(Dune::IOError, "Couldn't get 'Piece' node in " << fileName_ << ".");
372 }
373
379 bool hasData(const std::string& name, const DataType& type) const
380 {
381 using namespace tinyxml2;
382
383 const XMLElement* pieceNode = getPieceNode_();
384 const XMLElement* dataNode = getDataNode_(pieceNode, type);
385 if (dataNode == nullptr)
386 return false;
387
388 const XMLElement* dataArray = findDataArray_(dataNode, name);
389 if (dataArray == nullptr)
390 return false;
391
392 return true;
393 }
394
401 template<class Container>
402 Container readData(const std::string& name, const DataType& type) const
403 {
404 using namespace tinyxml2;
405
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_ << ".");
410
411 const XMLElement* dataArray = findDataArray_(dataNode, name);
412 if (dataArray == nullptr)
413 DUNE_THROW(Dune::IOError, "Couldn't find the data array " << name << ".");
414
415 return parseDataArray_<Container>(dataArray);
416 }
417
422 template<class Grid>
423 std::unique_ptr<Grid> readGrid(bool verbose = false) const
424 {
425 static_assert(!Dune::Capabilities::isCartesian<Grid>::v, "Grid reader only supports unstructured grid implementations");
426
427 // make a grid factory
428 Dune::GridFactory<Grid> factory;
429
430 // only read on rank 0
431 if (Dune::MPIHelper::instance().rank() == 0)
432 {
433 if (verbose)
434 std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << fileName_ << "." << std::endl;
435 readGrid_(factory, verbose);
436 }
437
438 return std::unique_ptr<Grid>(factory.createGrid());
439 }
440
447 template<class Grid>
448 std::unique_ptr<Grid> readGrid(Dune::GridFactory<Grid>& factory, bool verbose = false) const
449 {
450 static_assert(!Dune::Capabilities::isCartesian<Grid>::v, "Grid reader only supports unstructured grid implementations");
451
452 if (Dune::MPIHelper::instance().rank() == 0)
453 {
454 if (verbose)
455 std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << fileName_ << "." << std::endl;
456 readGrid_(factory, verbose);
457 }
458
459 return std::unique_ptr<Grid>(factory.createGrid());
460 }
461
470 template<class Grid>
471 std::unique_ptr<Grid> readGrid(Dune::GridFactory<Grid>& factory, Data& cellData, Data& pointData, bool verbose = false) const
472 {
473 static_assert(!Dune::Capabilities::isCartesian<Grid>::v, "Grid reader only supports unstructured grid implementations");
474
475 if (Dune::MPIHelper::instance().rank() == 0)
476 {
477 if (verbose)
478 std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << fileName_ << "." << std::endl;
479 readGrid_(factory, verbose);
480 readGridData_(cellData, pointData, verbose);
481 }
482
483 return std::unique_ptr<Grid>(factory.createGrid());
484 }
485
486private:
492 template<class Grid>
493 void readGrid_(Dune::GridFactory<Grid>& factory, bool verbose = false) const
494 {
495 using namespace tinyxml2;
496
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_ << ".");
501
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));
507
508 // adapt point dimensions if grid dimension is smaller than 3
509 auto points = adaptPointDimension_<Grid::dimensionworld>(std::move(points3D));
510
511 if (verbose) std::cout << "Found " << points.size() << " vertices." << std::endl;
512
513 // insert vertices to the grid factory
514 for (auto&& point : points)
515 factory.insertVertex(std::move(point));
516
517 const XMLElement* cellsNode = pieceNode->FirstChildElement("Cells");
518 const XMLElement* linesNode = pieceNode->FirstChildElement("Lines");
519 if (cellsNode)
520 {
521 const XMLElement* connectivityNode = findDataArray_(cellsNode, "connectivity");
522 const XMLElement* offsetsNode = findDataArray_(cellsNode, "offsets");
523 const XMLElement* typesNode = findDataArray_(cellsNode, "types");
524
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);
528
529 if (verbose) std::cout << "Found " << offsets.size() << " elements." << std::endl;
530
531 unsigned int lastOffset = 0;
532 for (unsigned int i = 0; i < offsets.size(); ++i)
533 {
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));
540 lastOffset = offset;
541 }
542 }
543 // for poly data
544 else if (linesNode)
545 {
546 // sanity check
547 if (Grid::dimension != 1)
548 DUNE_THROW(Dune::IOError, "Grid expects dimension " << Grid::dimension
549 << " but " << fileName_ << " contains a 1D grid.");
550
551 const XMLElement* connectivityNode = findDataArray_(linesNode, "connectivity");
552 const XMLElement* offsetsNode = findDataArray_(linesNode, "offsets");
553
554 const auto connectivity = parseDataArray_<std::vector<unsigned int>>(connectivityNode);
555 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
556
557 if (verbose) std::cout << "Found " << offsets.size() << " polylines." << std::endl;
558
559 unsigned int lastOffset = 0;
560 for (unsigned int i = 0; i < offsets.size(); ++i)
561 {
562 // a polyline can have many points in the VTK format
563 // split the line in segments with two points
564 unsigned int offset = offsets[i];
565 for (unsigned int j = 0; j < offset-lastOffset-1; ++j)
566 factory.insertElement(Dune::GeometryTypes::line,
567 std::vector<unsigned int>({ connectivity[lastOffset+j], connectivity[lastOffset+j+1] }));
568 lastOffset = offset;
569 }
570 }
571 else
572 DUNE_THROW(Dune::IOError, "No Cells or Lines element found in " << fileName_);
573 }
574
581 void readGridData_(Data& cellData, Data& pointData, bool verbose = false) const
582 {
583 using namespace tinyxml2;
584
585 const XMLElement* pieceNode = getPieceNode_();
586 const XMLElement* cellDataNode = getDataNode_(pieceNode, DataType::cellData);
587 if (cellDataNode != nullptr)
588 {
589 const XMLElement* cellsNode = pieceNode->FirstChildElement("Cells");
590 const XMLElement* linesNode = pieceNode->FirstChildElement("Lines");
591
592 if (cellsNode)
593 {
594 const XMLElement* dataArray = cellDataNode->FirstChildElement("DataArray");
595 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
596 {
597 const char *attributeText = dataArray->Attribute("Name");
598
599 if (attributeText == nullptr)
600 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a cell data array.");
601
602 cellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
603
604 if (verbose)
605 std::cout << "Read cell data field " << attributeText << std::endl;
606 }
607 }
608 // for poly data
609 else if (linesNode)
610 {
611 // first parse all the cell data (each cell in this sense can be a polyline)
612 const XMLElement* dataArray = cellDataNode->FirstChildElement("DataArray");
613 if (dataArray)
614 {
615 Data polyLineCellData;
616 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
617 {
618 const char *attributeText = dataArray->Attribute("Name");
619
620 if (attributeText == nullptr)
621 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a cell data array.");
622
623 polyLineCellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
624
625 if (verbose)
626 std::cout << "Read cell data field " << attributeText << std::endl;
627 }
628
629 // a polyline can have many points in the VTK format
630 // we split the line in segments with two points
631 // so we also need to duplicate the cell data to fit the increased line number
632 const XMLElement* offsetsNode = findDataArray_(linesNode, "offsets");
633 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
634
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() << ")!");
639
640 // count the number of Dune cells to be able to resize the data vectors
641 unsigned int lastOffset = 0;
642 std::size_t numCells = 0;
643 for (unsigned int i = 0; i < offsets.size(); ++i)
644 {
645 unsigned int offset = offsets[i];
646 for (unsigned int j = 0; j < offset-lastOffset-1; ++j)
647 ++numCells;
648 lastOffset = offset;
649 }
650
651 // create the data arrays
652 for (const auto& dArray : polyLineCellData)
653 {
654 cellData[dArray.first] = std::vector<double>(numCells);
655 auto& cd = cellData[dArray.first];
656 const auto& pd = dArray.second;
657
658 lastOffset = 0;
659 std::size_t cellIdx = 0;
660 for (unsigned int i = 0; i < offsets.size(); ++i)
661 {
662 unsigned int offset = offsets[i];
663 for (unsigned int j = 0; j < offset-lastOffset-1; ++j)
664 cd[cellIdx++] = pd[i];
665 lastOffset = offset;
666 }
667 }
668 }
669 }
670 else
671 DUNE_THROW(Dune::IOError, "No Cells or Lines element found in " << fileName_);
672 }
673
674 const XMLElement* pointDataNode = getDataNode_(pieceNode, DataType::pointData);
675 if (pointDataNode != nullptr)
676 {
677 const XMLElement* dataArray = pointDataNode->FirstChildElement("DataArray");
678 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
679 {
680 const char *attributeText = dataArray->Attribute("Name");
681
682 if (attributeText == nullptr)
683 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a point data array.");
684
685 pointData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
686
687 if (verbose)
688 std::cout << "Read point data field " << attributeText << std::endl;
689 }
690 }
691 }
692
697 const tinyxml2::XMLElement* getPieceNode_() const
698 { return Detail::VTKReader::getPieceNode(doc_, fileName_); }
699
706 const tinyxml2::XMLElement* getDataNode_(const tinyxml2::XMLElement* pieceNode, const DataType& type) const
707 {
708 using namespace tinyxml2;
709
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");
715 else
716 DUNE_THROW(Dune::IOError, "Only cell and point data are supported.");
717
718 return dataNode;
719 }
720
727 const tinyxml2::XMLElement* findDataArray_(const tinyxml2::XMLElement* dataNode, const std::string& name) const
728 {
729 using namespace tinyxml2;
730
731 // loop over XML node siblings to find the correct data array
732 const XMLElement* dataArray = dataNode->FirstChildElement("DataArray");
733 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
734 {
735 const char *attributeText = dataArray->Attribute("Name");
736
737 if (attributeText == nullptr)
738 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a data array.");
739
740 if (attributeText == name)
741 break;
742 }
743
744 return dataArray;
745 }
746
752 template<class Container>
753 Container parseDataArray_(const tinyxml2::XMLElement* dataArray) const
754 {
755 std::stringstream dataStream(dataArray->GetText());
756 return readStreamToContainer<Container>(dataStream);
757 }
758
763 Dune::GeometryType vtkToDuneGeomType_(unsigned int vtkCellType) const
764 {
765 switch (vtkCellType)
766 {
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);
776 }
777 }
778
779 template<int dim, std::enable_if_t<(dim < 3), int> = 0>
780 std::vector<Dune::FieldVector<double, dim>>
781 adaptPointDimension_(std::vector<Dune::FieldVector<double, 3>>&& points3D) const
782 {
783 std::vector<Dune::FieldVector<double, dim>> points(points3D.size());
784 for (std::size_t i = 0; i < points.size(); ++i)
785 for (int j = 0; j < dim; ++j)
786 points[i][j] = points3D[i][j];
787 return points;
788 }
789
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); }
794
795 std::string fileName_;
796 tinyxml2::XMLDocument doc_;
797};
798
799} // end namespace Dumux
800
801#endif // DUMUX_HAVE_GRIDFORMAT
802
803#endif
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: adapt.hh:17
Definition: common/pdesolver.hh:24