24#ifndef DUMUX_IO_GRID_DATA_HH
25#define DUMUX_IO_GRID_DATA_HH
31#include <dune/common/exceptions.hh>
32#include <dune/common/parallel/communication.hh>
33#include <dune/common/parallel/mpihelper.hh>
34#include <dune/grid/common/gridfactory.hh>
35#include <dune/grid/io/file/dgfparser/parser.hh>
36#include <dune/grid/io/file/dgfparser/gridptr.hh>
41#include <dune/grid/uggrid.hh>
51struct isUG :
public std::false_type {};
55struct isUG<
Dune::UGGrid<dim>> :
public std::true_type {};
67 using Intersection =
typename Grid::LeafIntersection;
68 using Element =
typename Grid::template Codim<0>::Entity;
69 using Vertex =
typename Grid::template Codim<Grid::dimension>::Entity;
72 enum DataSourceType { dgf, gmsh, vtk };
76 GridData(std::shared_ptr<Grid> grid, std::shared_ptr<Dune::GridFactory<Grid>> factory,
77 std::vector<int>&& elementMarkers, std::vector<int>&& boundaryMarkers, std::vector<int>&& faceMarkers = std::vector<int>{})
79 , gridFactory_(factory)
80 , elementMarkers_(elementMarkers)
81 , boundaryMarkers_(boundaryMarkers)
82 , faceMarkers_(faceMarkers)
83 , dataSourceType_(DataSourceType::gmsh)
89 , dataSourceType_(DataSourceType::dgf)
93 GridData(std::shared_ptr<Grid> grid, std::shared_ptr<Dune::GridFactory<Grid>> factory,
96 , gridFactory_(factory)
98 , pointData_(pointData)
99 , dataSourceType_(DataSourceType::vtk)
112 const std::vector<double>&
parameters(
const Vertex& vertex)
const
114 if (dataSourceType_ == DataSourceType::dgf)
116 if (vertex.level() != 0)
117 DUNE_THROW(Dune::IOError,
"You can only obtain parameters for level 0 vertices!");
119 return dgfGrid_.parameters(vertex);
122 DUNE_THROW(Dune::InvalidStateException,
"The parameters method is only available if the grid was constructed with a DGF file.");
128 const std::vector<double>&
parameters(
const Element& element)
const
130 if (dataSourceType_ == DataSourceType::dgf)
132 if (element.hasFather())
134 auto level0Element = element;
135 while(level0Element.hasFather())
136 level0Element = level0Element.father();
138 return dgfGrid_.parameters(level0Element);
142 return dgfGrid_.parameters(element);
146 DUNE_THROW(Dune::InvalidStateException,
"The parameters method is only available if the grid was constructed with a DGF file.");
152 template <
class Gr
idImp,
class IntersectionImp>
153 const Dune::DGFBoundaryParameter::type&
parameters(
const Dune::Intersection<GridImp, IntersectionImp>& intersection)
const
155 if (dataSourceType_ == DataSourceType::dgf)
156 return dgfGrid_.parameters(intersection);
158 DUNE_THROW(Dune::InvalidStateException,
"The parameters method is only available if the grid was constructed with a DGF file.");
175 if (dataSourceType_ != DataSourceType::gmsh)
176 DUNE_THROW(Dune::InvalidStateException,
"Domain markers are only available for gmsh grids.");
177 if (boundarySegmentIndex >= boundaryMarkers_.size())
178 DUNE_THROW(Dune::RangeError,
"Boundary segment index "<< boundarySegmentIndex <<
" bigger than number of boundary segments in grid.\n"
179 "Make sure to call this function only for boundaries that were defined as physical entities in gmsh.");
180 return boundaryMarkers_[boundarySegmentIndex];
195 {
return gridFactory_->wasInserted(intersection); }
204 if (dataSourceType_ != DataSourceType::gmsh)
205 DUNE_THROW(Dune::InvalidStateException,
"Domain markers are only available for gmsh grids.");
208 auto level0element = element;
209 while (level0element.hasFather())
210 level0element = level0element.father();
215 return elementMarkers_[gridPtr_->levelGridView(0).indexSet().index(level0element)];
217 return elementMarkers_[gridFactory_->insertionIndex(level0element)];
224 template<bool ug = Detail::isUG<Grid>::value,
typename std::enable_if_t<!ug, int> = 0>
227 return DataHandle(*gridPtr_, *gridFactory_, elementMarkers_, boundaryMarkers_, faceMarkers_);
234 template<bool ug = Detail::isUG<Grid>::value,
typename std::enable_if_t<ug, int> = 0>
237 return DataHandle(*gridPtr_, *gridFactory_, elementMarkers_, boundaryMarkers_);
252 double getParameter(
const Element& element,
const std::string& fieldName)
const
254 if (dataSourceType_ != DataSourceType::vtk)
255 DUNE_THROW(Dune::InvalidStateException,
"This access function is only available for data from VTK files.");
257 if (cellData_.count(fieldName) == 0)
258 DUNE_THROW(Dune::IOError,
"No field with the name " << fieldName <<
" found in cell data");
261 auto level0element = element;
262 while (level0element.hasFather())
263 level0element = level0element.father();
265 return cellData_.at(fieldName)[gridFactory_->insertionIndex(level0element)];
274 double getParameter(
const Vertex& vertex,
const std::string& fieldName)
const
276 if (dataSourceType_ != DataSourceType::vtk)
277 DUNE_THROW(Dune::InvalidStateException,
"This access function is only available for data from VTK files.");
279 if (vertex.level() != 0)
280 DUNE_THROW(Dune::IOError,
"You can only obtain parameters for level 0 vertices!");
282 if (pointData_.count(fieldName) == 0)
283 DUNE_THROW(Dune::IOError,
"No field with the name " << fieldName <<
" found in point data");
285 return pointData_.at(fieldName)[gridFactory_->insertionIndex(vertex)];
292 std::shared_ptr<Grid> gridPtr_;
293 std::shared_ptr<Dune::GridFactory<Grid>> gridFactory_;
299 std::vector<int> elementMarkers_;
300 std::vector<int> boundaryMarkers_;
301 std::vector<int> faceMarkers_;
309 Dune::GridPtr<Grid> dgfGrid_;
313 DataSourceType dataSourceType_;
A data handle for commucating grid data for gmsh grids.
A vtk file reader using tinyxml2 as xml backend.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
Definition: deprecated.hh:149
A data handle for commucating grid data for gmsh grids.
Definition: gmshgriddatahandle.hh:49
Definition: griddata.hh:51
Class for grid data attached to dgf or gmsh grid files.
Definition: griddata.hh:66
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:76
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:274
int getBoundaryDomainMarker(const Intersection &intersection) const
Return the boundary domain marker (Gmsh physical entity number) of an intersection Only available whe...
Definition: griddata.hh:188
double getParameter(const Element &element, const std::string &fieldName) const
Get a element parameter.
Definition: griddata.hh:252
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:93
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:153
GridData(Dune::GridPtr< Grid > grid)
constructor for dgf grid data
Definition: griddata.hh:87
DataHandle createGmshDataHandle()
Create a data handle for communication of the data in parallel simulations.
Definition: griddata.hh:225
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:112
int getBoundaryDomainMarker(int boundarySegmentIndex) const
Return the boundary domain marker (Gmsh physical entity number) of an intersection Only available whe...
Definition: griddata.hh:173
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:202
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:128
bool wasInserted(const Intersection &intersection) const
Returns true if an intersection was inserted during grid creation.
Definition: griddata.hh:194
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:59