24#ifndef DUMUX_IO_GRID_DATA_HH
25#define DUMUX_IO_GRID_DATA_HH
31#include <dune/common/exceptions.hh>
33#include <dune/common/version.hh>
34#if DUNE_VERSION_LT(DUNE_COMMON,2,7)
35#include <dune/common/parallel/collectivecommunication.hh>
37#include <dune/common/parallel/communication.hh>
40#include <dune/common/parallel/mpihelper.hh>
41#include <dune/grid/common/gridfactory.hh>
42#include <dune/grid/io/file/dgfparser/parser.hh>
43#include <dune/grid/io/file/dgfparser/gridptr.hh>
48#include <dune/grid/uggrid.hh>
58struct isUG :
public std::false_type {};
62struct isUG<
Dune::UGGrid<dim>> :
public std::true_type {};
74 using Intersection =
typename Grid::LeafIntersection;
75 using Element =
typename Grid::template Codim<0>::Entity;
76 using Vertex =
typename Grid::template Codim<Grid::dimension>::Entity;
79 enum DataSourceType { dgf, gmsh, vtk };
83 GridData(std::shared_ptr<Grid> grid, std::shared_ptr<Dune::GridFactory<Grid>> factory,
84 std::vector<int>&& elementMarkers, std::vector<int>&& boundaryMarkers, std::vector<int>&& faceMarkers = std::vector<int>{})
86 , gridFactory_(factory)
87 , elementMarkers_(elementMarkers)
88 , boundaryMarkers_(boundaryMarkers)
89 , faceMarkers_(faceMarkers)
90 , dataSourceType_(DataSourceType::gmsh)
96 , dataSourceType_(DataSourceType::dgf)
100 GridData(std::shared_ptr<Grid> grid, std::shared_ptr<Dune::GridFactory<Grid>> factory,
103 , gridFactory_(factory)
104 , cellData_(cellData)
105 , pointData_(pointData)
106 , dataSourceType_(DataSourceType::vtk)
119 const std::vector<double>&
parameters(
const Vertex& vertex)
const
121 if (dataSourceType_ == DataSourceType::dgf)
123 if (vertex.level() != 0)
124 DUNE_THROW(Dune::IOError,
"You can only obtain parameters for level 0 vertices!");
126 return dgfGrid_.parameters(vertex);
129 DUNE_THROW(Dune::InvalidStateException,
"The parameters method is only available if the grid was constructed with a DGF file.");
135 const std::vector<double>&
parameters(
const Element& element)
const
137 if (dataSourceType_ == DataSourceType::dgf)
139 if (element.hasFather())
141 auto level0Element = element;
142 while(level0Element.hasFather())
143 level0Element = level0Element.father();
145 return dgfGrid_.parameters(level0Element);
149 return dgfGrid_.parameters(element);
153 DUNE_THROW(Dune::InvalidStateException,
"The parameters method is only available if the grid was constructed with a DGF file.");
159 template <
class Gr
idImp,
class IntersectionImp>
160 const Dune::DGFBoundaryParameter::type&
parameters(
const Dune::Intersection<GridImp, IntersectionImp>& intersection)
const
162 if (dataSourceType_ == DataSourceType::dgf)
163 return dgfGrid_.parameters(intersection);
165 DUNE_THROW(Dune::InvalidStateException,
"The parameters method is only available if the grid was constructed with a DGF file.");
182 if (dataSourceType_ != DataSourceType::gmsh)
183 DUNE_THROW(Dune::InvalidStateException,
"Domain markers are only available for gmsh grids.");
184 if (boundarySegmentIndex >= boundaryMarkers_.size())
185 DUNE_THROW(Dune::RangeError,
"Boundary segment index "<< boundarySegmentIndex <<
" bigger than number of boundary segments in grid.\n"
186 "Make sure to call this function only for boundaries that were defined as physical entities in gmsh.");
187 return boundaryMarkers_[boundarySegmentIndex];
202 {
return gridFactory_->wasInserted(intersection); }
211 if (dataSourceType_ != DataSourceType::gmsh)
212 DUNE_THROW(Dune::InvalidStateException,
"Domain markers are only available for gmsh grids.");
215 auto level0element = element;
216 while (level0element.hasFather())
217 level0element = level0element.father();
222 return elementMarkers_[gridPtr_->levelGridView(0).indexSet().index(level0element)];
224 return elementMarkers_[gridFactory_->insertionIndex(level0element)];
231 template<bool ug = Detail::isUG<Grid>::value,
typename std::enable_if_t<!ug, int> = 0>
234 return DataHandle(*gridPtr_, *gridFactory_, elementMarkers_, boundaryMarkers_, faceMarkers_);
241 template<bool ug = Detail::isUG<Grid>::value,
typename std::enable_if_t<ug, int> = 0>
244 return DataHandle(*gridPtr_, *gridFactory_, elementMarkers_, boundaryMarkers_);
259 double getParameter(
const Element& element,
const std::string& fieldName)
const
261 if (dataSourceType_ != DataSourceType::vtk)
262 DUNE_THROW(Dune::InvalidStateException,
"This access function is only available for data from VTK files.");
264 if (cellData_.count(fieldName) == 0)
265 DUNE_THROW(Dune::IOError,
"No field with the name " << fieldName <<
" found in cell data");
268 auto level0element = element;
269 while (level0element.hasFather())
270 level0element = level0element.father();
272 return cellData_.at(fieldName)[gridFactory_->insertionIndex(level0element)];
281 double getParameter(
const Vertex& vertex,
const std::string& fieldName)
const
283 if (dataSourceType_ != DataSourceType::vtk)
284 DUNE_THROW(Dune::InvalidStateException,
"This access function is only available for data from VTK files.");
286 if (vertex.level() != 0)
287 DUNE_THROW(Dune::IOError,
"You can only obtain parameters for level 0 vertices!");
289 if (pointData_.count(fieldName) == 0)
290 DUNE_THROW(Dune::IOError,
"No field with the name " << fieldName <<
" found in point data");
292 return pointData_.at(fieldName)[gridFactory_->insertionIndex(vertex)];
299 std::shared_ptr<Grid> gridPtr_;
300 std::shared_ptr<Dune::GridFactory<Grid>> gridFactory_;
306 std::vector<int> elementMarkers_;
307 std::vector<int> boundaryMarkers_;
308 std::vector<int> faceMarkers_;
316 Dune::GridPtr<Grid> dgfGrid_;
320 DataSourceType dataSourceType_;
A data handle for commucating grid data for gmsh grids.
A vtk file reader using tinyxml2 as xml backend.
Definition: common/pdesolver.hh:35
A data handle for commucating grid data for gmsh grids.
Definition: gmshgriddatahandle.hh:55
Definition: griddata.hh:58
Class for grid data attached to dgf or gmsh grid files.
Definition: griddata.hh:73
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:83
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:281
int getBoundaryDomainMarker(const Intersection &intersection) const
Return the boundary domain marker (Gmsh physical entity number) of an intersection Only available whe...
Definition: griddata.hh:195
double getParameter(const Element &element, const std::string &fieldName) const
Get a element parameter.
Definition: griddata.hh:259
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:100
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:160
GridData(Dune::GridPtr< Grid > grid)
constructor for dgf grid data
Definition: griddata.hh:94
DataHandle createGmshDataHandle()
Create a data handle for communication of the data in parallel simulations.
Definition: griddata.hh:232
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:119
int getBoundaryDomainMarker(int boundarySegmentIndex) const
Return the boundary domain marker (Gmsh physical entity number) of an intersection Only available whe...
Definition: griddata.hh:180
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:209
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:135
bool wasInserted(const Intersection &intersection) const
Returns true if an intersection was inserted during grid creation.
Definition: griddata.hh:201
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:57