version 3.10.0
Loading...
Searching...
No Matches
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
40
41template<GridFormat::Concepts::Communicator C>
42auto makeGridformatReaderFactory(const C& c)
43{
44 return [fac = GridFormat::AnyReaderFactory<C>{c}] (const std::string& f)
45 -> std::unique_ptr<GridFormat::GridReader> {
46 // use adapter for poly data that splits polylines into segments
47 if (f.ends_with("vtp"))
48 return std::make_unique<GridFormat::ReaderDecorators::PolylinesSubdivider>(fac.make_for(f));
49 return fac.make_for(f);
50 };
51}
52
53} // end namespace Dumux::Detail::VTKReader
54
55namespace Dumux {
56
61class VTKReader
62{
63public:
64 enum class DataType { cellData, pointData, fieldData };
65
67 using Data = std::unordered_map<std::string, std::vector<double>>;
68
69 explicit VTKReader(const std::string& fileName)
70 : reader_(GridFormat::Reader::from(fileName,
71 Detail::VTKReader::makeGridformatReaderFactory(
72 Dune::MPIHelper::instance().getCommunicator()
73 )
74 ))
75 {}
76
82 bool hasData(const std::string& name, const DataType& type) const
83 {
84 if (type == DataType::cellData)
85 return std::ranges::any_of(
86 cell_field_names(reader_),
87 [&] (const auto& n) { return name == n; }
88 );
89 else if (type == DataType::pointData)
90 return std::ranges::any_of(
91 point_field_names(reader_),
92 [&] (const auto& n) { return name == n; }
93 );
94 else if (type == DataType::fieldData)
95 return std::ranges::any_of(
96 meta_data_field_names(reader_),
97 [&] (const auto& n) { return name == n; }
98 );
99 else
100 DUNE_THROW(Dune::IOError, "Unknown data type for field '" << name << "'");
101 }
102
109 template<class Container>
110 Container readData(const std::string& name, const DataType& type) const
111 {
112 if (type == DataType::cellData)
113 {
114 Container values(reader_.number_of_cells());
115 reader_.cell_field(name)->export_to(values);
116 return values;
117 }
118 else if (type == DataType::pointData)
119 {
120 Container values(reader_.number_of_points());
121 reader_.point_field(name)->export_to(values);
122 return values;
123 }
124 else if (type == DataType::fieldData)
125 {
126 DUNE_THROW(Dune::NotImplemented, "Reading meta data not yet implemented");
127 }
128 else
129 DUNE_THROW(Dune::IOError, "Unknown data type for field '" << name << "'");
130 }
131
136 template<class Grid>
137 std::unique_ptr<Grid> readGrid(bool verbose = false) const
138 {
139 Dune::GridFactory<Grid> factory;
140 return readGrid(factory, verbose);
141 }
142
149 template<class Grid>
150 std::unique_ptr<Grid> readGrid(Dune::GridFactory<Grid>& factory, bool verbose = false) const
151 {
152 if (Dune::MPIHelper::instance().rank() == 0)
153 {
154 if (verbose)
155 std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << reader_.filename() << "." << std::endl;
156
157 {
158 GridFormat::Dune::GridFactoryAdapter<Grid> adapter{ factory };
159 reader_.export_grid(adapter);
160 }
161 }
162
163 return std::unique_ptr<Grid>(factory.createGrid());
164 }
165
174 template<class Grid>
175 std::unique_ptr<Grid> readGrid(Dune::GridFactory<Grid>& factory, Data& cellData, Data& pointData, bool verbose = false) const
176 {
177 auto grid = readGrid(factory, verbose);
178 if (Dune::MPIHelper::instance().rank() == 0)
179 {
180 for (const auto& [name, field_ptr] : cell_fields(reader_))
181 field_ptr->export_to(cellData[name]);
182
183 for (const auto& [name, field_ptr] : point_fields(reader_))
184 field_ptr->export_to(pointData[name]);
185 }
186 return std::unique_ptr<Grid>(std::move(grid));
187 }
188
189private:
190 GridFormat::Reader reader_;
191};
192
193} // namespace Dumux
194
195#else // DUMUX_HAVE_GRIDFORMAT
196
197#include <dumux/io/xml/tinyxml2.h>
198// fallback to simple vtk reader using tinyxml2
199// this reader only support ascii files and limited VTK file formats
200
202
209const tinyxml2::XMLElement* getPieceNode(const tinyxml2::XMLDocument& doc, const std::string& fileName)
210{
211 using namespace tinyxml2;
212
213 const XMLElement* pieceNode = doc.FirstChildElement("VTKFile");
214 if (pieceNode == nullptr)
215 DUNE_THROW(Dune::IOError, "Couldn't get 'VTKFile' node in " << fileName << ".");
216
217 pieceNode = pieceNode->FirstChildElement("UnstructuredGrid");
218 if (pieceNode == nullptr)
219 pieceNode = doc.FirstChildElement("VTKFile")->FirstChildElement("PolyData");
220 if (pieceNode == nullptr)
221 pieceNode = doc.FirstChildElement("VTKFile")->FirstChildElement("PUnstructuredGrid");
222 if (pieceNode == nullptr)
223 pieceNode = doc.FirstChildElement("VTKFile")->FirstChildElement("PPolyData");
224 if (pieceNode == nullptr)
225 DUNE_THROW(Dune::IOError, "Couldn't get 'UnstructuredGrid', 'PUnstructuredGrid', 'PolyData', or 'PPolyData' node in " << fileName << ".");
226
227 return pieceNode->FirstChildElement("Piece");
228}
229
233std::string getProcessPieceFileName(const std::string& pvtkFileName)
234{
235 using namespace tinyxml2;
236
237 XMLDocument pDoc;
238 const auto eResult = pDoc.LoadFile(pvtkFileName.c_str());
239 if (eResult != XML_SUCCESS)
240 DUNE_THROW(Dune::IOError, "Couldn't open XML file " << pvtkFileName << ".");
241
242 // get the first piece node
243 const XMLElement* pieceNode = getPieceNode(pDoc, pvtkFileName);
244 const auto myrank = Dune::MPIHelper::instance().rank();
245 for (int rank = 0; rank < myrank; ++rank)
246 {
247 pieceNode = pieceNode->NextSiblingElement("Piece");
248 if (pieceNode == nullptr)
249 DUNE_THROW(Dune::IOError, "Couldn't find 'Piece' node for rank "
250 << rank << " in " << pvtkFileName << ".");
251 }
252
253 const char *vtkFileName = pieceNode->Attribute("Source");
254 if (vtkFileName == nullptr)
255 DUNE_THROW(Dune::IOError, "Couldn't get 'Source' attribute of 'Piece' node no. " << myrank << " in " << pvtkFileName);
256
257 return vtkFileName;
258}
259
260} // end namespace Dumux::Detail::VTKReader
261
262namespace Dumux {
263
269{
270public:
274 enum class DataType { cellData, pointData };
275
277 using Data = std::unordered_map<std::string, std::vector<double>>;
278
282 explicit VTKReader(const std::string& fileName)
283 {
284 using namespace tinyxml2;
285 const auto ext = std::filesystem::path(fileName).extension().string();
286 // If in parallel and the file to read is a parallel piece collection (pvtu/pvtp)
287 // read only the piece belonging to the own process. For this to work, the files
288 // have to have exactly the same amount of pieces than processes.
289 fileName_ = Dune::MPIHelper::instance().size() > 1 && ext[1] == 'p' ?
291
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_ << ".");
295
296 const XMLElement* pieceNode = getPieceNode_();
297 if (pieceNode == nullptr)
298 DUNE_THROW(Dune::IOError, "Couldn't get 'Piece' node in " << fileName_ << ".");
299 }
300
306 bool hasData(const std::string& name, const DataType& type) const
307 {
308 using namespace tinyxml2;
309
310 const XMLElement* pieceNode = getPieceNode_();
311 const XMLElement* dataNode = getDataNode_(pieceNode, type);
312 if (dataNode == nullptr)
313 return false;
314
315 const XMLElement* dataArray = findDataArray_(dataNode, name);
316 if (dataArray == nullptr)
317 return false;
318
319 return true;
320 }
321
328 template<class Container>
329 Container readData(const std::string& name, const DataType& type) const
330 {
331 using namespace tinyxml2;
332
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_ << ".");
337
338 const XMLElement* dataArray = findDataArray_(dataNode, name);
339 if (dataArray == nullptr)
340 DUNE_THROW(Dune::IOError, "Couldn't find the data array " << name << ".");
341
342 return parseDataArray_<Container>(dataArray);
343 }
344
349 template<class Grid>
350 std::unique_ptr<Grid> readGrid(bool verbose = false) const
351 {
352 static_assert(!Dune::Capabilities::isCartesian<Grid>::v, "Grid reader only supports unstructured grid implementations");
353
354 // make a grid factory
355 Dune::GridFactory<Grid> factory;
356
357 // only read on rank 0
358 if (Dune::MPIHelper::instance().rank() == 0)
359 {
360 if (verbose)
361 std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << fileName_ << "." << std::endl;
362 readGrid_(factory, verbose);
363 }
364
365 return std::unique_ptr<Grid>(factory.createGrid());
366 }
367
374 template<class Grid>
375 std::unique_ptr<Grid> readGrid(Dune::GridFactory<Grid>& factory, bool verbose = false) const
376 {
377 static_assert(!Dune::Capabilities::isCartesian<Grid>::v, "Grid reader only supports unstructured grid implementations");
378
379 if (Dune::MPIHelper::instance().rank() == 0)
380 {
381 if (verbose)
382 std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << fileName_ << "." << std::endl;
383 readGrid_(factory, verbose);
384 }
385
386 return std::unique_ptr<Grid>(factory.createGrid());
387 }
388
397 template<class Grid>
398 std::unique_ptr<Grid> readGrid(Dune::GridFactory<Grid>& factory, Data& cellData, Data& pointData, bool verbose = false) const
399 {
400 static_assert(!Dune::Capabilities::isCartesian<Grid>::v, "Grid reader only supports unstructured grid implementations");
401
402 if (Dune::MPIHelper::instance().rank() == 0)
403 {
404 if (verbose)
405 std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << fileName_ << "." << std::endl;
406 readGrid_(factory, verbose);
407 readGridData_(cellData, pointData, verbose);
408 }
409
410 return std::unique_ptr<Grid>(factory.createGrid());
411 }
412
413private:
419 template<class Grid>
420 void readGrid_(Dune::GridFactory<Grid>& factory, bool verbose = false) const
421 {
422 using namespace tinyxml2;
423
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_ << ".");
428
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));
434
435 // adapt point dimensions if grid dimension is smaller than 3
436 auto points = adaptPointDimension_<Grid::dimensionworld>(std::move(points3D));
437
438 if (verbose) std::cout << "Found " << points.size() << " vertices." << std::endl;
439
440 // insert vertices to the grid factory
441 for (auto&& point : points)
442 factory.insertVertex(std::move(point));
443
444 const XMLElement* cellsNode = pieceNode->FirstChildElement("Cells");
445 const XMLElement* linesNode = pieceNode->FirstChildElement("Lines");
446 if (cellsNode)
447 {
448 const XMLElement* connectivityNode = findDataArray_(cellsNode, "connectivity");
449 const XMLElement* offsetsNode = findDataArray_(cellsNode, "offsets");
450 const XMLElement* typesNode = findDataArray_(cellsNode, "types");
451
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);
455
456 if (verbose) std::cout << "Found " << offsets.size() << " elements." << std::endl;
457
458 unsigned int lastOffset = 0;
459 for (unsigned int i = 0; i < offsets.size(); ++i)
460 {
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));
467 lastOffset = offset;
468 }
469 }
470 // for poly data
471 else if (linesNode)
472 {
473 // sanity check
474 if (Grid::dimension != 1)
475 DUNE_THROW(Dune::IOError, "Grid expects dimension " << Grid::dimension
476 << " but " << fileName_ << " contains a 1D grid.");
477
478 const XMLElement* connectivityNode = findDataArray_(linesNode, "connectivity");
479 const XMLElement* offsetsNode = findDataArray_(linesNode, "offsets");
480
481 const auto connectivity = parseDataArray_<std::vector<unsigned int>>(connectivityNode);
482 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
483
484 if (verbose) std::cout << "Found " << offsets.size() << " polylines." << std::endl;
485
486 unsigned int lastOffset = 0;
487 for (unsigned int i = 0; i < offsets.size(); ++i)
488 {
489 // a polyline can have many points in the VTK format
490 // split the line in segments with two points
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] }));
495 lastOffset = offset;
496 }
497 }
498 else
499 DUNE_THROW(Dune::IOError, "No Cells or Lines element found in " << fileName_);
500 }
501
508 void readGridData_(Data& cellData, Data& pointData, bool verbose = false) const
509 {
510 using namespace tinyxml2;
511
512 const XMLElement* pieceNode = getPieceNode_();
513 const XMLElement* cellDataNode = getDataNode_(pieceNode, DataType::cellData);
514 if (cellDataNode != nullptr)
515 {
516 const XMLElement* cellsNode = pieceNode->FirstChildElement("Cells");
517 const XMLElement* linesNode = pieceNode->FirstChildElement("Lines");
518
519 if (cellsNode)
520 {
521 const XMLElement* dataArray = cellDataNode->FirstChildElement("DataArray");
522 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
523 {
524 const char *attributeText = dataArray->Attribute("Name");
525
526 if (attributeText == nullptr)
527 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a cell data array.");
528
529 cellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
530
531 if (verbose)
532 std::cout << "Read cell data field " << attributeText << std::endl;
533 }
534 }
535 // for poly data
536 else if (linesNode)
537 {
538 // first parse all the cell data (each cell in this sense can be a polyline)
539 const XMLElement* dataArray = cellDataNode->FirstChildElement("DataArray");
540 if (dataArray)
541 {
542 Data polyLineCellData;
543 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
544 {
545 const char *attributeText = dataArray->Attribute("Name");
546
547 if (attributeText == nullptr)
548 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a cell data array.");
549
550 polyLineCellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
551
552 if (verbose)
553 std::cout << "Read cell data field " << attributeText << std::endl;
554 }
555
556 // a polyline can have many points in the VTK format
557 // we split the line in segments with two points
558 // so we also need to duplicate the cell data to fit the increased line number
559 const XMLElement* offsetsNode = findDataArray_(linesNode, "offsets");
560 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
561
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() << ")!");
566
567 // count the number of Dune cells to be able to resize the data vectors
568 unsigned int lastOffset = 0;
569 std::size_t numCells = 0;
570 for (unsigned int i = 0; i < offsets.size(); ++i)
571 {
572 unsigned int offset = offsets[i];
573 for (unsigned int j = 0; j < offset-lastOffset-1; ++j)
574 ++numCells;
575 lastOffset = offset;
576 }
577
578 // create the data arrays
579 for (const auto& dArray : polyLineCellData)
580 {
581 cellData[dArray.first] = std::vector<double>(numCells);
582 auto& cd = cellData[dArray.first];
583 const auto& pd = dArray.second;
584
585 lastOffset = 0;
586 std::size_t cellIdx = 0;
587 for (unsigned int i = 0; i < offsets.size(); ++i)
588 {
589 unsigned int offset = offsets[i];
590 for (unsigned int j = 0; j < offset-lastOffset-1; ++j)
591 cd[cellIdx++] = pd[i];
592 lastOffset = offset;
593 }
594 }
595 }
596 }
597 else
598 DUNE_THROW(Dune::IOError, "No Cells or Lines element found in " << fileName_);
599 }
600
601 const XMLElement* pointDataNode = getDataNode_(pieceNode, DataType::pointData);
602 if (pointDataNode != nullptr)
603 {
604 const XMLElement* dataArray = pointDataNode->FirstChildElement("DataArray");
605 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
606 {
607 const char *attributeText = dataArray->Attribute("Name");
608
609 if (attributeText == nullptr)
610 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a point data array.");
611
612 pointData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
613
614 if (verbose)
615 std::cout << "Read point data field " << attributeText << std::endl;
616 }
617 }
618 }
619
624 const tinyxml2::XMLElement* getPieceNode_() const
625 { return Detail::VTKReader::getPieceNode(doc_, fileName_); }
626
633 const tinyxml2::XMLElement* getDataNode_(const tinyxml2::XMLElement* pieceNode, const DataType& type) const
634 {
635 using namespace tinyxml2;
636
637 const XMLElement* dataNode = nullptr;
638 if (type == DataType::pointData)
639 dataNode = pieceNode->FirstChildElement("PointData");
640 else if (type == DataType::cellData)
641 dataNode = pieceNode->FirstChildElement("CellData");
642 else
643 DUNE_THROW(Dune::IOError, "Only cell and point data are supported.");
644
645 return dataNode;
646 }
647
654 const tinyxml2::XMLElement* findDataArray_(const tinyxml2::XMLElement* dataNode, const std::string& name) const
655 {
656 using namespace tinyxml2;
657
658 // loop over XML node siblings to find the correct data array
659 const XMLElement* dataArray = dataNode->FirstChildElement("DataArray");
660 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
661 {
662 const char *attributeText = dataArray->Attribute("Name");
663
664 if (attributeText == nullptr)
665 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a data array.");
666
667 if (attributeText == name)
668 break;
669 }
670
671 return dataArray;
672 }
673
679 template<class Container>
680 Container parseDataArray_(const tinyxml2::XMLElement* dataArray) const
681 {
682 std::stringstream dataStream(dataArray->GetText());
683 return readStreamToContainer<Container>(dataStream);
684 }
685
690 Dune::GeometryType vtkToDuneGeomType_(unsigned int vtkCellType) const
691 {
692 switch (vtkCellType)
693 {
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);
703 }
704 }
705
706 template<int dim, std::enable_if_t<(dim < 3), int> = 0>
707 std::vector<Dune::FieldVector<double, dim>>
708 adaptPointDimension_(std::vector<Dune::FieldVector<double, 3>>&& points3D) const
709 {
710 std::vector<Dune::FieldVector<double, dim>> points(points3D.size());
711 for (std::size_t i = 0; i < points.size(); ++i)
712 for (int j = 0; j < dim; ++j)
713 points[i][j] = points3D[i][j];
714 return points;
715 }
716
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); }
721
722 std::string fileName_;
723 tinyxml2::XMLDocument doc_;
724};
725
726} // end namespace Dumux
727
728#endif // DUMUX_HAVE_GRIDFORMAT
729
730#endif
A vtk file reader using tinyxml2 as xml backend.
Definition vtkreader.hh:269
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:398
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:306
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:329
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:277
VTKReader(const std::string &fileName)
The constructor creates a tinyxml2::XMLDocument from file.
Definition vtkreader.hh:282
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:375
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:350
DataType
The data array types.
Definition vtkreader.hh:274
@ pointData
Definition vtkreader.hh:274
@ cellData
Definition vtkreader.hh:274
Free functions to write and read a sequence container to and from a file.
Definition vtkreader.hh:201
const tinyxml2::XMLElement * getPieceNode(const tinyxml2::XMLDocument &doc, const std::string &fileName)
Get the piece node an xml document.
Definition vtkreader.hh:209
std::string getProcessPieceFileName(const std::string &pvtkFileName)
get the vtk filename for the current processor
Definition vtkreader.hh:233
Definition adapt.hh:17
Container readStreamToContainer(std::istream &stream)
Read an input stream into a container.
Definition container.hh:55