24#ifndef DUMUX_IO_GRID_PORENETWORK_SUBGRID_DATA_HH
25#define DUMUX_IO_GRID_PORENETWORK_SUBGRID_DATA_HH
34template<
class HostGr
id,
class SubGr
id>
37 static constexpr int dim = SubGrid::dimension;
38 static constexpr int dimWorld = SubGrid::dimensionworld;
39 using Intersection =
typename SubGrid::LeafIntersection;
40 using Element =
typename SubGrid::template Codim<0>::Entity;
41 using Vertex =
typename SubGrid::template Codim<dim>::Entity;
42 using GridView =
typename SubGrid::LeafGridView;
43 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
50 , hostGridData_(hostGridData)
57 const std::vector<double>&
parameters(
const Vertex& vertex)
const
58 {
return hostGridData_->parameters(vertex.impl().hostEntity()); }
63 const std::vector<double>&
parameters(
const Element& element)
const
64 {
return hostGridData_->parameters(element.impl().hostEntity()); }
69 auto getParameter(
const Element& element,
const std::string& param)
const
70 {
return hostGridData_->getParameter(element.impl().hostEntity(), param); }
75 auto getParameter(
const Vertex& vertex,
const std::string& param)
const
76 {
return hostGridData_->getParameter(vertex.impl().hostEntity(), param); }
81 template <
class Gr
idImp,
class IntersectionImp>
82 const Dune::DGFBoundaryParameter::type&
parameters(
const Dune::Intersection<GridImp, IntersectionImp>& intersection)
const
83 {
return hostGridData_->paramters(intersection); }
91 {
return hostGridData_->throatLabel(element.impl().hostEntity()); }
99 {
return hostGridData_->boundaryFaceMarkerAtPos(pos); }
106 const auto gridView = subGrid_.leafGridView();
108 std::vector<SmallLocalIndex> coordinationNumbers(gridView.size(dim), 0);
110 for (
const auto &element : elements(gridView))
112 for (SmallLocalIndex vIdxLocal = 0; vIdxLocal < 2; ++vIdxLocal)
114 const auto vIdxGlobal = gridView.indexSet().subIndex(element, vIdxLocal, dim);
115 coordinationNumbers[vIdxGlobal] += 1;
119 if (std::any_of(coordinationNumbers.begin(), coordinationNumbers.end(), [](
auto i){ return i == 0; }))
120 DUNE_THROW(Dune::InvalidStateException,
"One of the pores is not connected to another pore. SanitizeGrid will not help in this case. Check your grid file");
122 return coordinationNumbers;
129 {
return hostGridData_->parameterIndex(paramName); }
136 {
return hostGridData_->paramGroup(); }
142 {
return hostGridData_->gridHasElementParameter(param); }
148 {
return hostGridData_->gridHasVertexParameter(param); }
154 {
return hostGridData_->vertexParameterNames() ; }
160 {
return hostGridData_->elementParameterNames() ; }
163 const SubGrid& subGrid_;
164 std::shared_ptr<const GridData<HostGrid>> hostGridData_;
Definition: discretization/porenetwork/fvelementgeometry.hh:33
std::uint_least8_t SmallLocalIndex
Definition: indextraits.hh:41
wrapper for subgrid data
Definition: subgriddata.hh:36
auto getParameter(const Element &element, const std::string ¶m) const
Returns the value of an element parameter.
Definition: subgriddata.hh:69
std::vector< SmallLocalIndex > getCoordinationNumbers() const
Returns the coordination numbers for all pore bodies.
Definition: subgriddata.hh:104
int parameterIndex(const std::string ¶mName) const
Return the index for a given parameter name.
Definition: subgriddata.hh:128
SubGridData(const SubGrid &subGrid, std::shared_ptr< const GridData< HostGrid > > hostGridData)
Definition: subgriddata.hh:47
bool gridHasElementParameter(const std::string ¶m) const
Return if a given element parameter is provided by the grid.
Definition: subgriddata.hh:141
const std::string & paramGroup() const
Return the parameter group.
Definition: subgriddata.hh:135
const std::vector< double > & parameters(const Element &element) const
Call the parameters function of the DGF grid pointer if available for element data.
Definition: subgriddata.hh:63
auto throatLabel(const Element &element) const
Computes and returns the label of a given throat.
Definition: subgriddata.hh:90
const std::vector< std::string > & vertexParameterNames() const
Returns the names of the vertex parameters.
Definition: subgriddata.hh:153
const std::vector< double > & parameters(const Vertex &vertex) const
Call the parameters function of the DGF grid pointer if available for vertex data.
Definition: subgriddata.hh:57
auto getParameter(const Vertex &vertex, const std::string ¶m) const
Returns the value of an vertex parameter.
Definition: subgriddata.hh:75
bool gridHasVertexParameter(const std::string ¶m) const
Return if a given vertex parameter is provided by the grid.
Definition: subgriddata.hh:147
const Dune::DGFBoundaryParameter::type & parameters(const Dune::Intersection< GridImp, IntersectionImp > &intersection) const
Call the parameters function of the DGF grid pointer if available.
Definition: subgriddata.hh:82
const std::vector< std::string > & elementParameterNames() const
Returns the names of the element parameters.
Definition: subgriddata.hh:159
int boundaryFaceMarkerAtPos(const GlobalPosition &pos) const
Returns the boundary face marker index at given position.
Definition: subgriddata.hh:98
Class for grid data attached to dgf or gmsh grid files.