3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_IO_VTK_VTKREADER_HH
25#define DUMUX_IO_VTK_VTKREADER_HH
26
27#include <iostream>
28#include <iterator>
29#include <algorithm>
30#include <memory>
31#include <type_traits>
32#include <unordered_map>
33
34#include <dune/common/parallel/mpihelper.hh>
35#include <dune/common/exceptions.hh>
36#include <dune/grid/common/capabilities.hh>
37#include <dune/grid/io/file/vtk/common.hh>
38#include <dumux/io/xml/tinyxml2.h>
39#include <dune/grid/common/gridfactory.hh>
40
41namespace Dumux {
42
48{
49public:
53 enum class DataType { cellData, pointData };
54
56 using Data = std::unordered_map<std::string, std::vector<double>>;
57
61 explicit VTKReader(const std::string& fileName)
62 {
63 fileName_ = Dune::MPIHelper::getCollectiveCommunication().size() > 1 ?
64 getProcessFileName_(fileName) : fileName;
65
66 const auto eResult = doc_.LoadFile(fileName_.c_str());
67 if (eResult != tinyxml2::XML_SUCCESS)
68 DUNE_THROW(Dune::IOError, "Couldn't open XML file " << fileName_ << ".");
69 }
70
77 template<class Container>
78 Container readData(const std::string& name, const DataType& type) const
79 {
80 using namespace tinyxml2;
81
82 const XMLElement* pieceNode = getPieceNode_();
83 if (pieceNode == nullptr)
84 DUNE_THROW(Dune::IOError, "Couldn't get 'Piece' node in " << fileName_ << ".");
85
86 const XMLElement* dataNode = getDataNode_(pieceNode, type);
87 if (dataNode == nullptr)
88 DUNE_THROW(Dune::IOError, "Couldn't get 'PointData' or 'CellData' node in " << fileName_ << ".");
89
90 const XMLElement* dataArray = findDataArray_(dataNode, name);
91 if (dataArray == nullptr)
92 DUNE_THROW(Dune::IOError, "Couldn't find the data array " << name << ".");
93
94 return parseDataArray_<Container>(dataArray);
95 }
96
101 template<class Grid>
102 std::unique_ptr<Grid> readGrid(bool verbose = false) const
103 {
104 static_assert(!Dune::Capabilities::isCartesian<Grid>::v, "Grid reader only supports unstructured grid implementations");
105
106 if (verbose) std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << fileName_ << "." << std::endl;
107
108 // make a grid factory
109 Dune::GridFactory<Grid> factory;
110
111 readGrid_(factory, verbose);
112
113 return std::unique_ptr<Grid>(factory.createGrid());
114 }
115
122 template<class Grid>
123 std::unique_ptr<Grid> readGrid(Dune::GridFactory<Grid>& factory, bool verbose = false) const
124 {
125 static_assert(!Dune::Capabilities::isCartesian<Grid>::v, "Grid reader only supports unstructured grid implementations");
126
127 if (verbose) std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << fileName_ << "." << std::endl;
128
129 readGrid_(factory, verbose);
130
131 return std::unique_ptr<Grid>(factory.createGrid());
132 }
133
142 template<class Grid>
143 std::unique_ptr<Grid> readGrid(Dune::GridFactory<Grid>& factory, Data& cellData, Data& pointData, bool verbose = false) const
144 {
145 static_assert(!Dune::Capabilities::isCartesian<Grid>::v, "Grid reader only supports unstructured grid implementations");
146
147 if (verbose) std::cout << "Reading " << Grid::dimension << "d grid from vtk file " << fileName_ << "." << std::endl;
148
149 readGrid_(factory, verbose);
150 readGridData_(cellData, pointData, verbose);
151
152 return std::unique_ptr<Grid>(factory.createGrid());
153 }
154
155private:
159 std::string getProcessFileName_(const std::string& pvtkFileName)
160 {
161 using namespace tinyxml2;
162
163 XMLDocument pDoc;
164 const auto eResult = pDoc.LoadFile(pvtkFileName.c_str());
165 if (eResult != XML_SUCCESS)
166 DUNE_THROW(Dune::IOError, "Couldn't open XML file " << pvtkFileName << ".");
167
168 // get the first piece node
169 const XMLElement* pieceNode = getPieceNode_(pDoc, pvtkFileName);
170 if (pieceNode == nullptr)
171 DUNE_THROW(Dune::IOError, "Couldn't get 'Piece' node in " << pvtkFileName << ".");
172
173 const auto myrank = Dune::MPIHelper::getCollectiveCommunication().rank();
174 for (int rank = 0; rank < myrank; ++rank)
175 {
176 pieceNode = pieceNode->NextSiblingElement("Piece");
177 if (pieceNode == nullptr)
178 DUNE_THROW(Dune::IOError, "Couldn't find 'Piece' node for rank "
179 << rank << " in " << pvtkFileName << ".");
180 }
181
182 const char *vtkFileName = pieceNode->Attribute("Source");
183 if (vtkFileName == nullptr)
184 DUNE_THROW(Dune::IOError, "Couldn't get 'Source' attribute of 'Piece' node no. " << myrank << " in " << pvtkFileName);
185
186 return vtkFileName;
187 }
188
194 template<class Grid>
195 void readGrid_(Dune::GridFactory<Grid>& factory, bool verbose = false) const
196 {
197 using namespace tinyxml2;
198
199 const XMLElement* pieceNode = getPieceNode_();
200 if (pieceNode == nullptr)
201 DUNE_THROW(Dune::IOError, "Couldn't get 'Piece' node in " << fileName_ << ".");
202
203 const XMLElement* pointsNode = pieceNode->FirstChildElement("Points")->FirstChildElement("DataArray");
204 if (pointsNode == nullptr)
205 DUNE_THROW(Dune::IOError, "Couldn't get data array of points in " << fileName_ << ".");
206
207 using Point3D = Dune::FieldVector<double, 3>;
208 std::vector<Point3D> points3D;
209 std::stringstream dataStream(pointsNode->GetText());
210 std::istream_iterator<Point3D> it(dataStream);
211 std::copy(it, std::istream_iterator<Point3D>(), std::back_inserter(points3D));
212
213 // adapt point dimensions if grid dimension is smaller than 3
214 auto points = adaptPointDimension_<Grid::dimensionworld>(std::move(points3D));
215
216 if (verbose) std::cout << "Found " << points.size() << " vertices." << std::endl;
217
218 // insert vertices to the grid factory
219 for (auto&& point : points)
220 factory.insertVertex(std::move(point));
221
222 const XMLElement* cellsNode = pieceNode->FirstChildElement("Cells");
223 const XMLElement* linesNode = pieceNode->FirstChildElement("Lines");
224 if (cellsNode)
225 {
226 const XMLElement* connectivityNode = findDataArray_(cellsNode, "connectivity");
227 const XMLElement* offsetsNode = findDataArray_(cellsNode, "offsets");
228 const XMLElement* typesNode = findDataArray_(cellsNode, "types");
229
230 const auto connectivity = parseDataArray_<std::vector<unsigned int>>(connectivityNode);
231 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
232 const auto types = parseDataArray_<std::vector<unsigned int>>(typesNode);
233
234 if (verbose) std::cout << "Found " << offsets.size() << " element." << std::endl;
235
236 unsigned int lastOffset = 0;
237 for (unsigned int i = 0; i < offsets.size(); ++i)
238 {
239 const auto geomType = vtkToDuneGeomType_(types[i]);
240 unsigned int offset = offsets[i];
241 std::vector<unsigned int> corners; corners.resize(offset-lastOffset);
242 for (unsigned int j = 0; j < offset-lastOffset; ++j)
243 corners[Dune::VTK::renumber(geomType, j)] = connectivity[lastOffset+j];
244 factory.insertElement(geomType, std::move(corners));
245 lastOffset = offset;
246 }
247 }
248 // for poly data
249 else if (linesNode)
250 {
251 // sanity check
252 if (Grid::dimension != 1)
253 DUNE_THROW(Dune::IOError, "Grid expects dimension " << Grid::dimension
254 << " but " << fileName_ << " contains a 1D grid.");
255
256 const XMLElement* connectivityNode = findDataArray_(linesNode, "connectivity");
257 const XMLElement* offsetsNode = findDataArray_(linesNode, "offsets");
258
259 const auto connectivity = parseDataArray_<std::vector<unsigned int>>(connectivityNode);
260 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
261
262 if (verbose) std::cout << "Found " << offsets.size() << " polylines." << std::endl;
263
264 unsigned int lastOffset = 0;
265 for (unsigned int i = 0; i < offsets.size(); ++i)
266 {
267 // a polyline can have many points in the VTK format
268 // split the line in segments with two points
269 unsigned int offset = offsets[i];
270 for (unsigned int j = 0; j < offset-lastOffset-1; ++j)
271 factory.insertElement(Dune::GeometryTypes::line,
272 std::vector<unsigned int>({ connectivity[lastOffset+j], connectivity[lastOffset+j+1] }));
273 lastOffset = offset;
274 }
275 }
276 else
277 DUNE_THROW(Dune::IOError, "No Cells or Lines element found in " << fileName_);
278 }
279
286 void readGridData_(Data& cellData, Data& pointData, bool verbose = false) const
287 {
288 using namespace tinyxml2;
289
290 const XMLElement* pieceNode = getPieceNode_();
291 if (pieceNode == nullptr)
292 DUNE_THROW(Dune::IOError, "Couldn't get 'Piece' node in " << fileName_ << ".");
293
294 const XMLElement* cellDataNode = getDataNode_(pieceNode, DataType::cellData);
295 if (cellDataNode != nullptr)
296 {
297 const XMLElement* cellsNode = pieceNode->FirstChildElement("Cells");
298 const XMLElement* linesNode = pieceNode->FirstChildElement("Lines");
299
300 if (cellsNode)
301 {
302 const XMLElement* dataArray = cellDataNode->FirstChildElement("DataArray");
303 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
304 {
305 const char *attributeText = dataArray->Attribute("Name");
306
307 if (attributeText == nullptr)
308 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a cell data array.");
309
310 cellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
311 }
312 }
313 // for poly data
314 else if (linesNode)
315 {
316 // first parse all the cell data (each cell in this sense can be a polyline)
317 const XMLElement* dataArray = cellDataNode->FirstChildElement("DataArray");
318 if (dataArray)
319 {
320 Data polyLineCellData;
321 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
322 {
323 const char *attributeText = dataArray->Attribute("Name");
324
325 if (attributeText == nullptr)
326 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a cell data array.");
327
328 polyLineCellData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
329 }
330
331 // a polyline can have many points in the VTK format
332 // we split the line in segments with two points
333 // so we also need to duplicate the cell data to fit the increased line number
334 const XMLElement* offsetsNode = findDataArray_(linesNode, "offsets");
335 const auto offsets = parseDataArray_<std::vector<unsigned int>>(offsetsNode);
336
337 if (offsets.size() != polyLineCellData.begin()->second.size())
338 DUNE_THROW(Dune::IOError, "Expected the same number of cell data entries (is "
339 << polyLineCellData.begin()->second.size()
340 << ") as polylines (" << offsets.size() << ")!");
341
342 // count the number of Dune cells to be able to resize the data vectors
343 unsigned int lastOffset = 0;
344 std::size_t numCells = 0;
345 for (unsigned int i = 0; i < offsets.size(); ++i)
346 {
347 unsigned int offset = offsets[i];
348 for (unsigned int j = 0; j < offset-lastOffset-1; ++j)
349 ++numCells;
350 lastOffset = offset;
351 }
352
353 // create the data arrays
354 for (const auto& dArray : polyLineCellData)
355 {
356 cellData[dArray.first] = std::vector<double>(numCells);
357 auto& cd = cellData[dArray.first];
358 const auto& pd = dArray.second;
359
360 lastOffset = 0;
361 std::size_t cellIdx = 0;
362 for (unsigned int i = 0; i < offsets.size(); ++i)
363 {
364 unsigned int offset = offsets[i];
365 for (unsigned int j = 0; j < offset-lastOffset-1; ++j)
366 cd[cellIdx++] = pd[i];
367 lastOffset = offset;
368 }
369 }
370 }
371 }
372 else
373 DUNE_THROW(Dune::IOError, "No Cells or Lines element found in " << fileName_);
374 }
375
376 const XMLElement* pointDataNode = getDataNode_(pieceNode, DataType::pointData);
377 if (pointDataNode != nullptr)
378 {
379 const XMLElement* dataArray = pointDataNode->FirstChildElement("DataArray");
380 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
381 {
382 const char *attributeText = dataArray->Attribute("Name");
383
384 if (attributeText == nullptr)
385 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a point data array.");
386
387 pointData[std::string(attributeText)] = parseDataArray_<std::vector<double>>(dataArray);
388 }
389 }
390 }
391
396 const tinyxml2::XMLElement* getPieceNode_() const
397 { return getPieceNode_(doc_, fileName_); }
398
405 const tinyxml2::XMLElement* getPieceNode_(const tinyxml2::XMLDocument& doc, const std::string& fileName) const
406 {
407 using namespace tinyxml2;
408
409 const XMLElement* pieceNode = doc.FirstChildElement("VTKFile");
410 if (pieceNode == nullptr)
411 DUNE_THROW(Dune::IOError, "Couldn't get 'VTKFile' node in " << fileName << ".");
412
413 pieceNode = pieceNode->FirstChildElement("UnstructuredGrid");
414 if (pieceNode == nullptr)
415 pieceNode = doc.FirstChildElement("VTKFile")->FirstChildElement("PolyData");
416 if (pieceNode == nullptr)
417 pieceNode = doc.FirstChildElement("VTKFile")->FirstChildElement("PUnstructuredGrid");
418 if (pieceNode == nullptr)
419 pieceNode = doc.FirstChildElement("VTKFile")->FirstChildElement("PPolyData");
420 if (pieceNode == nullptr)
421 DUNE_THROW(Dune::IOError, "Couldn't get 'UnstructuredGrid', 'PUnstructuredGrid', 'PolyData', or 'PPolyData' node in " << fileName << ".");
422
423 return pieceNode->FirstChildElement("Piece");
424 }
425
432 const tinyxml2::XMLElement* getDataNode_(const tinyxml2::XMLElement* pieceNode, const DataType& type) const
433 {
434 using namespace tinyxml2;
435
436 const XMLElement* dataNode = nullptr;
437 if (type == DataType::pointData)
438 dataNode = pieceNode->FirstChildElement("PointData");
439 else if (type == DataType::cellData)
440 dataNode = pieceNode->FirstChildElement("CellData");
441 else
442 DUNE_THROW(Dune::IOError, "Only cell and point data are supported.");
443
444 return dataNode;
445 }
446
453 const tinyxml2::XMLElement* findDataArray_(const tinyxml2::XMLElement* dataNode, const std::string& name) const
454 {
455 using namespace tinyxml2;
456
457 // loop over XML node siblings to find the correct data array
458 const XMLElement* dataArray = dataNode->FirstChildElement("DataArray");
459 for (; dataArray != nullptr; dataArray = dataArray->NextSiblingElement("DataArray"))
460 {
461 const char *attributeText = dataArray->Attribute("Name");
462
463 if (attributeText == nullptr)
464 DUNE_THROW(Dune::IOError, "Couldn't get Name attribute of a data array.");
465
466 if (attributeText == name)
467 break;
468 }
469
470 return dataArray;
471 }
472
478 template<class Container>
479 Container parseDataArray_(const tinyxml2::XMLElement* dataArray) const
480 {
481 Container data;
482 std::stringstream dataStream(dataArray->GetText());
483 std::istream_iterator<typename Container::value_type> it(dataStream);
484 std::copy(it, std::istream_iterator<typename Container::value_type>(), std::back_inserter(data));
485 return data;
486 }
487
492 Dune::GeometryType vtkToDuneGeomType_(unsigned int vtkCellType) const
493 {
494 switch (vtkCellType)
495 {
496 case Dune::VTK::GeometryType::vertex: return Dune::GeometryTypes::vertex;
497 case Dune::VTK::GeometryType::line: return Dune::GeometryTypes::line;
498 case Dune::VTK::GeometryType::triangle: return Dune::GeometryTypes::triangle;
499 case Dune::VTK::GeometryType::quadrilateral: return Dune::GeometryTypes::quadrilateral;
500 case Dune::VTK::GeometryType::hexahedron: return Dune::GeometryTypes::hexahedron;
501 case Dune::VTK::GeometryType::prism: return Dune::GeometryTypes::prism;
502 case Dune::VTK::GeometryType::pyramid: return Dune::GeometryTypes::pyramid;
503 default: DUNE_THROW(Dune::NotImplemented, "VTK cell type " << vtkCellType);
504 }
505 }
506
507 template<int dim, std::enable_if_t<(dim < 3), int> = 0>
508 std::vector<Dune::FieldVector<double, dim>>
509 adaptPointDimension_(std::vector<Dune::FieldVector<double, 3>>&& points3D) const
510 {
511 std::vector<Dune::FieldVector<double, dim>> points(points3D.size());
512 for (std::size_t i = 0; i < points.size(); ++i)
513 for (int j = 0; j < dim; ++j)
514 points[i][j] = points3D[i][j];
515 return points;
516 }
517
518 template<int dim, std::enable_if_t<(dim == 3), int> = 0>
519 std::vector<Dune::FieldVector<double, dim>>
520 adaptPointDimension_(std::vector<Dune::FieldVector<double, 3>>&& points3D) const
521 { return points3D; }
522
523 std::string fileName_;
524 tinyxml2::XMLDocument doc_;
525};
526
527} // end namespace Dumux
528
529#endif
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
A vtk file reader using tinyxml2 as xml backend.
Definition: vtkreader.hh:48
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:143
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:78
VTKReader(const std::string &fileName)
The contructor creates a tinyxml2::XMLDocument from file.
Definition: vtkreader.hh:61
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:123
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:102
DataType
The data array types.
Definition: vtkreader.hh:53
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:56