24#ifndef DUMUX_DISCRETIZATION_STAGGERED_FV_ELEMENT_GEOMETRY_HH
25#define DUMUX_DISCRETIZATION_STAGGERED_FV_ELEMENT_GEOMETRY_HH
41template<
class GG,
bool enableGr
idGeometryCache>
56 using GridView =
typename GG::GridView;
63 using ParentType::ParentType;
67 template<
class CellCenterOrFaceFVGr
idGeometry>
69 :
ParentType(gridGeometry.actualGridGeometry()) {}
72 using ParentType::scvf;
75 return this->gridGeometry().scvf(eIdx, localScvfIdx);
91 using GridView =
typename GG::GridView;
94 using Element =
typename GridView::template Codim<0>::Entity;
105 template<
class CellCenterOrFaceFVGr
idGeometry>
107 : gridGeometryPtr_(&gridGeometry.actualGridGeometry()) {}
111 : gridGeometryPtr_(&gridGeometry) {}
116 return scvf(this->gridGeometry().localToGlobalScvfIndex(eIdx, localScvfIdx));
123 if (scvIdx == scvIndices_[0])
126 return neighborScvs_[findLocalIndex_(scvIdx, neighborScvIndices_)];
133 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
134 if (it != scvfIndices_.end())
135 return scvfs_[std::distance(scvfIndices_.begin(), it)];
137 return neighborScvfs_[findLocalIndex_(scvfIdx, neighborScvfIndices_)];
145 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
148 using IteratorType =
typename std::array<SubControlVolume, 1>::const_iterator;
149 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
157 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
160 using IteratorType =
typename std::vector<SubControlVolumeFace>::const_iterator;
161 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
166 {
return scvs_.size(); }
170 {
return scvfs_.size(); }
174 void bind(
const Element& element)
176 bindElement(element);
178 neighborScvs_.reserve(element.subEntities(1));
179 neighborScvfIndices_.reserve(element.subEntities(1));
180 neighborScvfs_.reserve(element.subEntities(1));
182 std::vector<GridIndexType> handledNeighbors;
183 handledNeighbors.reserve(element.subEntities(1));
184 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
186 if (intersection.neighbor())
188 const auto outside = intersection.outside();
189 const auto outsideIdx = gridGeometry().elementMapper().index(outside);
192 if ( std::find(handledNeighbors.begin(), handledNeighbors.end(), outsideIdx) == handledNeighbors.end() )
194 makeNeighborGeometries_(outside, outsideIdx);
195 handledNeighbors.push_back(outsideIdx);
205 elementPtr_ = &element;
206 scvfs_.reserve(element.subEntities(1));
207 scvfIndices_.reserve(element.subEntities(1));
208 makeElementGeometries_(element);
213 {
return *gridGeometryPtr_; }
217 {
return hasBoundaryScvf_; }
222 void makeElementGeometries_(
const Element& element)
224 const auto eIdx = gridGeometry().elementMapper().index(element);
225 scvs_[0] = SubControlVolume(element.geometry(), eIdx);
226 scvIndices_[0] = eIdx;
228 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
229 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
231 typename GridGeometry::GeometryHelper geometryHelper(element, gridGeometry().gridView());
234 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
236 const auto& scvfNeighborVolVarIndex = neighborVolVarIndices[scvfCounter];
238 if (intersection.neighbor() || intersection.boundary())
240 geometryHelper.updateLocalFace(gridGeometry().intersectionMapper(), intersection);
241 std::vector<GridIndexType> scvIndices{eIdx, scvfNeighborVolVarIndex};
242 scvfs_.emplace_back(intersection,
243 intersection.geometry(),
244 scvFaceIndices[scvfCounter],
247 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
250 if (intersection.boundary())
251 hasBoundaryScvf_ =
true;
257 void makeNeighborGeometries_(
const Element& element,
const GridIndexType eIdx)
262 neighborScvs_.emplace_back(element.geometry(), eIdx);
263 neighborScvIndices_.push_back(eIdx);
265 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
266 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
268 typename GridGeometry::GeometryHelper geometryHelper(element, gridGeometry().gridView());
271 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
273 const auto& scvfNeighborVolVarIndex = neighborVolVarIndices[scvfCounter];
274 geometryHelper.updateLocalFace(gridGeometry().intersectionMapper(), intersection);
276 if (intersection.neighbor())
279 if (intersection.outside() == *elementPtr_)
281 std::vector<GridIndexType> scvIndices{eIdx, scvfNeighborVolVarIndex};
282 neighborScvfs_.emplace_back(intersection,
283 intersection.geometry(),
284 scvFaceIndices[scvfCounter],
288 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
292 else if (intersection.boundary())
297 const LocalIndexType findLocalIndex_(
const GridIndexType idx,
298 const std::vector<GridIndexType>& indices)
const
300 auto it = std::find(indices.begin(), indices.end(), idx);
301 assert(it != indices.end() &&
"Could not find the scv/scvf! Make sure to properly bind this class!");
302 return std::distance(indices.begin(), it);
308 scvfIndices_.clear();
311 neighborScvIndices_.clear();
312 neighborScvfIndices_.clear();
313 neighborScvs_.clear();
314 neighborScvfs_.clear();
316 hasBoundaryScvf_ =
false;
319 const Element* elementPtr_;
320 const GridGeometry* gridGeometryPtr_;
323 std::array<GridIndexType, 1> scvIndices_;
324 std::array<SubControlVolume, 1> scvs_;
326 std::vector<GridIndexType> scvfIndices_;
327 std::vector<SubControlVolumeFace> scvfs_;
329 std::vector<GridIndexType> neighborScvIndices_;
330 std::vector<SubControlVolume> neighborScvs_;
332 std::vector<GridIndexType> neighborScvfIndices_;
333 std::vector<SubControlVolumeFace> neighborScvfs_;
335 bool hasBoundaryScvf_ =
false;
Defines the index types used for grid and local indices.
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:39
unsigned int LocalIndex
Definition: indextraits.hh:40
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models This builds up th...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:50
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:60
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:70
Stencil-local finite volume geometry (scvs and scvfs) for staggered models This builds up the sub con...
Definition: discretization/staggered/fvelementgeometry.hh:42
const SubControlVolumeFace & scvf(GridIndexType eIdx, LocalIndexType localScvfIdx) const
Definition: discretization/staggered/fvelementgeometry.hh:73
StaggeredFVElementGeometry(const CellCenterOrFaceFVGridGeometry &gridGeometry)
Definition: discretization/staggered/fvelementgeometry.hh:68
Base class for the finite volume geometry vector for staggered models This locally builds up the sub ...
Definition: discretization/staggered/fvelementgeometry.hh:89
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/staggered/fvelementgeometry.hh:216
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/staggered/fvelementgeometry.hh:158
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/staggered/fvelementgeometry.hh:97
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/staggered/fvelementgeometry.hh:121
StaggeredFVElementGeometry(const CellCenterOrFaceFVGridGeometry &gridGeometry)
Definition: discretization/staggered/fvelementgeometry.hh:106
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/staggered/fvelementgeometry.hh:165
const SubControlVolumeFace & scvf(GridIndexType eIdx, LocalIndexType localScvfIdx) const
Get a sub control volume face with an element index and a local scvf index.
Definition: discretization/staggered/fvelementgeometry.hh:114
const GridGeometry & gridGeometry() const
The grid finite volume geometry we are a restriction of.
Definition: discretization/staggered/fvelementgeometry.hh:212
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/staggered/fvelementgeometry.hh:169
StaggeredFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/staggered/fvelementgeometry.hh:110
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/staggered/fvelementgeometry.hh:99
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/staggered/fvelementgeometry.hh:146
void bindElement(const Element &element)
Binding of an element preparing the geometries only inside the element.
Definition: discretization/staggered/fvelementgeometry.hh:202
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/staggered/fvelementgeometry.hh:131
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/staggered/fvelementgeometry.hh:101
void bind(const Element &element)
Definition: discretization/staggered/fvelementgeometry.hh:174
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models This builds up th...