12#ifndef DUMUX_IO_GRID_DATA_HH
13#define DUMUX_IO_GRID_DATA_HH
19#include <dune/common/exceptions.hh>
20#include <dune/common/parallel/communication.hh>
21#include <dune/common/parallel/mpihelper.hh>
22#include <dune/grid/common/gridfactory.hh>
23#include <dune/grid/io/file/dgfparser/parser.hh>
24#include <dune/grid/io/file/dgfparser/gridptr.hh>
29#include <dune/grid/uggrid.hh>
39struct isUG :
public std::false_type {};
43struct isUG<
Dune::UGGrid<dim>> :
public std::true_type {};
55 using Intersection =
typename Grid::LeafIntersection;
56 using Element =
typename Grid::template Codim<0>::Entity;
57 using Vertex =
typename Grid::template Codim<Grid::dimension>::Entity;
60 enum DataSourceType { dgf, gmsh, vtk };
64 GridData(std::shared_ptr<Grid> grid, std::shared_ptr<Dune::GridFactory<Grid>> factory,
65 std::vector<int>&& elementMarkers, std::vector<int>&& boundaryMarkers, std::vector<int>&& faceMarkers = std::vector<int>{})
67 , gridFactory_(factory)
68 , elementMarkers_(elementMarkers)
69 , boundaryMarkers_(boundaryMarkers)
70 , faceMarkers_(faceMarkers)
71 , dataSourceType_(DataSourceType::gmsh)
77 , dataSourceType_(DataSourceType::dgf)
81 GridData(std::shared_ptr<Grid> grid, std::shared_ptr<Dune::GridFactory<Grid>> factory,
84 , gridFactory_(factory)
86 , pointData_(pointData)
87 , dataSourceType_(DataSourceType::vtk)
100 const std::vector<double>&
parameters(
const Vertex& vertex)
const
102 if (dataSourceType_ == DataSourceType::dgf)
104 if (vertex.level() != 0)
105 DUNE_THROW(Dune::IOError,
"You can only obtain parameters for level 0 vertices!");
107 return dgfGrid_.parameters(vertex);
110 DUNE_THROW(Dune::InvalidStateException,
"The parameters method is only available if the grid was constructed with a DGF file.");
116 const std::vector<double>&
parameters(
const Element& element)
const
118 if (dataSourceType_ == DataSourceType::dgf)
120 if (element.hasFather())
122 auto level0Element = element;
123 while(level0Element.hasFather())
124 level0Element = level0Element.father();
126 return dgfGrid_.parameters(level0Element);
130 return dgfGrid_.parameters(element);
134 DUNE_THROW(Dune::InvalidStateException,
"The parameters method is only available if the grid was constructed with a DGF file.");
140 template <
class Gr
idImp,
class IntersectionImp>
141 const Dune::DGFBoundaryParameter::type&
parameters(
const Dune::Intersection<GridImp, IntersectionImp>& intersection)
const
143 if (dataSourceType_ == DataSourceType::dgf)
144 return dgfGrid_.parameters(intersection);
146 DUNE_THROW(Dune::InvalidStateException,
"The parameters method is only available if the grid was constructed with a DGF file.");
163 if (dataSourceType_ != DataSourceType::gmsh)
164 DUNE_THROW(Dune::InvalidStateException,
"Domain markers are only available for gmsh grids.");
165 if (boundarySegmentIndex >= boundaryMarkers_.size())
166 DUNE_THROW(Dune::RangeError,
"Boundary segment index "<< boundarySegmentIndex <<
" bigger than number of boundary segments in grid.\n"
167 "Make sure to call this function only for boundaries that were defined as physical entities in gmsh.");
168 return boundaryMarkers_[boundarySegmentIndex];
183 {
return gridFactory_->wasInserted(intersection); }
192 if (dataSourceType_ != DataSourceType::gmsh)
193 DUNE_THROW(Dune::InvalidStateException,
"Domain markers are only available for gmsh grids.");
196 auto level0element = element;
197 while (level0element.hasFather())
198 level0element = level0element.father();
203 return elementMarkers_[gridPtr_->levelGridView(0).indexSet().index(level0element)];
205 return elementMarkers_[gridFactory_->insertionIndex(level0element)];
212 template<bool ug = Detail::isUG<Grid>::value,
typename std::enable_if_t<!ug, int> = 0>
215 return DataHandle(*gridPtr_, *gridFactory_, elementMarkers_, boundaryMarkers_, faceMarkers_);
222 template<bool ug = Detail::isUG<Grid>::value,
typename std::enable_if_t<ug, int> = 0>
225 return DataHandle(*gridPtr_, *gridFactory_, elementMarkers_, boundaryMarkers_);
240 double getParameter(
const Element& element,
const std::string& fieldName)
const
242 if (dataSourceType_ != DataSourceType::vtk)
243 DUNE_THROW(Dune::InvalidStateException,
"This access function is only available for data from VTK files.");
245 if (cellData_.count(fieldName) == 0)
246 DUNE_THROW(Dune::IOError,
"No field with the name " << fieldName <<
" found in cell data");
249 auto level0element = element;
250 while (level0element.hasFather())
251 level0element = level0element.father();
253 return cellData_.at(fieldName)[gridFactory_->insertionIndex(level0element)];
262 double getParameter(
const Vertex& vertex,
const std::string& fieldName)
const
264 if (dataSourceType_ != DataSourceType::vtk)
265 DUNE_THROW(Dune::InvalidStateException,
"This access function is only available for data from VTK files.");
267 if (vertex.level() != 0)
268 DUNE_THROW(Dune::IOError,
"You can only obtain parameters for level 0 vertices!");
270 if (pointData_.count(fieldName) == 0)
271 DUNE_THROW(Dune::IOError,
"No field with the name " << fieldName <<
" found in point data");
273 return pointData_.at(fieldName)[gridFactory_->insertionIndex(vertex)];
280 std::shared_ptr<Grid> gridPtr_;
281 std::shared_ptr<Dune::GridFactory<Grid>> gridFactory_;
287 std::vector<int> elementMarkers_;
288 std::vector<int> boundaryMarkers_;
289 std::vector<int> faceMarkers_;
297 Dune::GridPtr<Grid> dgfGrid_;
301 DataSourceType dataSourceType_;
Class for grid data attached to dgf or gmsh grid files.
Definition: griddata.hh:54
GridData(std::shared_ptr< Grid > grid, std::shared_ptr< Dune::GridFactory< Grid > > factory, std::vector< int > &&elementMarkers, std::vector< int > &&boundaryMarkers, std::vector< int > &&faceMarkers=std::vector< int >{})
constructor for gmsh grid data
Definition: griddata.hh:64
double getParameter(const Vertex &vertex, const std::string &fieldName) const
Call the parameters function of the DGF grid pointer if available for vertex data.
Definition: griddata.hh:262
int getBoundaryDomainMarker(const Intersection &intersection) const
Return the boundary domain marker (Gmsh physical entity number) of an intersection Only available whe...
Definition: griddata.hh:176
double getParameter(const Element &element, const std::string &fieldName) const
Get a element parameter.
Definition: griddata.hh:240
GridData(std::shared_ptr< Grid > grid, std::shared_ptr< Dune::GridFactory< Grid > > factory, VTKReader::Data &&cellData, VTKReader::Data &&pointData)
constructor for gmsh grid data
Definition: griddata.hh:81
const Dune::DGFBoundaryParameter::type & parameters(const Dune::Intersection< GridImp, IntersectionImp > &intersection) const
Call the parameters function of the DGF grid pointer if available.
Definition: griddata.hh:141
GridData(Dune::GridPtr< Grid > grid)
constructor for dgf grid data
Definition: griddata.hh:75
DataHandle createGmshDataHandle()
Create a data handle for communication of the data in parallel simulations.
Definition: griddata.hh:213
const std::vector< double > & parameters(const Vertex &vertex) const
Call the parameters function of the DGF grid pointer if available for vertex data.
Definition: griddata.hh:100
int getBoundaryDomainMarker(int boundarySegmentIndex) const
Return the boundary domain marker (Gmsh physical entity number) of an intersection Only available whe...
Definition: griddata.hh:161
int getElementDomainMarker(const Element &element) const
Return the element domain marker (Gmsh physical entity number) of an element. Only available when usi...
Definition: griddata.hh:190
const std::vector< double > & parameters(const Element &element) const
Call the parameters function of the DGF grid pointer if available for element data.
Definition: griddata.hh:116
bool wasInserted(const Intersection &intersection) const
Returns true if an intersection was inserted during grid creation.
Definition: griddata.hh:182
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:47
A data handle for commucating grid data for gmsh grids.
Definition: common/pdesolver.hh:24
Definition: griddata.hh:39
A data handle for commucating grid data for gmsh grids.
Definition: gmshgriddatahandle.hh:37
A vtk file reader using tinyxml2 as xml backend.