3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
subgriddata.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_GRID_PORENETWORK_SUBGRID_DATA_HH
25#define DUMUX_IO_GRID_PORENETWORK_SUBGRID_DATA_HH
26
27#include <vector>
28#include <memory>
29#include "griddata.hh"
30
31namespace Dumux::PoreNetwork {
32
34template<class HostGrid, class SubGrid>
36{
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;
44 using SmallLocalIndex = typename IndexTraits<GridView>::SmallLocalIndex;
45
46public:
47 SubGridData(const SubGrid& subGrid,
48 std::shared_ptr<const GridData<HostGrid>> hostGridData)
49 : subGrid_(subGrid)
50 , hostGridData_(hostGridData)
51 {}
52
57 const std::vector<double>& parameters(const Vertex& vertex) const
58 { return hostGridData_->parameters(vertex.impl().hostEntity()); }
59
63 const std::vector<double>& parameters(const Element& element) const
64 { return hostGridData_->parameters(element.impl().hostEntity()); }
65
69 auto getParameter(const Element& element, const std::string& param) const
70 { return hostGridData_->getParameter(element.impl().hostEntity(), param); }
71
75 auto getParameter(const Vertex& vertex, const std::string& param) const
76 { return hostGridData_->getParameter(vertex.impl().hostEntity(), param); }
77
81 template <class GridImp, class IntersectionImp>
82 const Dune::DGFBoundaryParameter::type& parameters(const Dune::Intersection<GridImp, IntersectionImp>& intersection) const
83 { return hostGridData_->paramters(intersection); }
84
90 auto throatLabel(const Element& element) const
91 { return hostGridData_->throatLabel(element.impl().hostEntity()); }
92
98 int boundaryFaceMarkerAtPos(const GlobalPosition& pos) const
99 { return hostGridData_->boundaryFaceMarkerAtPos(pos); }
100
104 std::vector<SmallLocalIndex> getCoordinationNumbers() const
105 {
106 const auto gridView = subGrid_.leafGridView();
107
108 std::vector<SmallLocalIndex> coordinationNumbers(gridView.size(dim), 0);
109
110 for (const auto &element : elements(gridView))
111 {
112 for (SmallLocalIndex vIdxLocal = 0; vIdxLocal < 2; ++vIdxLocal)
113 {
114 const auto vIdxGlobal = gridView.indexSet().subIndex(element, vIdxLocal, dim);
115 coordinationNumbers[vIdxGlobal] += 1;
116 }
117 }
118
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");
121
122 return coordinationNumbers;
123 }
124
128 int parameterIndex(const std::string& paramName) const
129 { return hostGridData_->parameterIndex(paramName); }
130
135 const std::string& paramGroup() const
136 { return hostGridData_->paramGroup(); }
137
141 bool gridHasElementParameter(const std::string& param) const
142 { return hostGridData_->gridHasElementParameter(param); }
143
147 bool gridHasVertexParameter(const std::string& param) const
148 { return hostGridData_->gridHasVertexParameter(param); }
149
153 const std::vector<std::string>& vertexParameterNames() const
154 { return hostGridData_->vertexParameterNames() ; }
155
159 const std::vector<std::string>& elementParameterNames() const
160 { return hostGridData_->elementParameterNames() ; }
161
162private:
163 const SubGrid& subGrid_;
164 std::shared_ptr<const GridData<HostGrid>> hostGridData_;
165};
166
167} // end namespace Dumux::PoreNetwork
168
169#endif
Definition: discretization/porenetwork/fvelementgeometry.hh:34
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 &param) 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 &paramName) 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 &param) 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 &param) const
Returns the value of an vertex parameter.
Definition: subgriddata.hh:75
bool gridHasVertexParameter(const std::string &param) 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.