24#ifndef DUMUX_DISCRETIZATION_STAGGERED_FV_ELEMENT_GEOMETRY_HH
25#define DUMUX_DISCRETIZATION_STAGGERED_FV_ELEMENT_GEOMETRY_HH
42template<
class GG,
bool enableGr
idGeometryCache>
57 using GridView =
typename GG::GridView;
62 using Element =
typename GridView::template Codim<0>::Entity;
66 using ParentType::ParentType;
70 template<
class CellCenterOrFaceFVGr
idGeometry>
72 :
ParentType(gridGeometry.actualGridGeometry()) {}
75 using ParentType::scvf;
78 return this->gridGeometry().scvf(eIdx, localScvfIdx);
94 using GridView =
typename GG::GridView;
99 using Element =
typename GridView::template Codim<0>::Entity;
109 template<
class CellCenterOrFaceFVGr
idGeometry>
111 : gridGeometryPtr_(&gridGeometry.actualGridGeometry()) {}
115 : gridGeometryPtr_(&gridGeometry) {}
120 return scvf(this->gridGeometry().localToGlobalScvfIndex(eIdx, localScvfIdx));
127 if (scvIdx == scvIndices_[0])
130 return neighborScvs_[findLocalIndex_(scvIdx, neighborScvIndices_)];
137 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
138 if (it != scvfIndices_.end())
141 return neighborScvfs_[findLocalIndex_(scvfIdx, neighborScvfIndices_)];
149 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
152 using IteratorType =
typename std::array<SubControlVolume, 1>::const_iterator;
153 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
161 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
164 using IteratorType =
typename std::vector<SubControlVolumeFace>::const_iterator;
165 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
170 {
return scvs_.size(); }
174 {
return scvfs_.size(); }
183 this->bind_(element);
184 return std::move(*
this);
189 { this->bind_(element); }
198 this->bindElement_(element);
199 return std::move(*
this);
204 { this->bindElement_(element); }
208 {
return static_cast<bool>(element_); }
212 {
return *element_; }
216 {
return *gridGeometryPtr_; }
220 {
return hasBoundaryScvf_; }
225 return scvf.geometry();
230 void bindElement_(
const Element& element)
234 scvfs_.reserve(element.subEntities(1));
235 scvfIndices_.reserve(element.subEntities(1));
236 makeElementGeometries_();
241 void bind_(
const Element& element)
243 bindElement_(element);
245 neighborScvs_.reserve(element.subEntities(1));
246 neighborScvfIndices_.reserve(element.subEntities(1));
247 neighborScvfs_.reserve(element.subEntities(1));
249 std::vector<GridIndexType> handledNeighbors;
250 handledNeighbors.reserve(element.subEntities(1));
251 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
253 if (intersection.neighbor())
255 const auto outside = intersection.outside();
256 const auto outsideIdx = gridGeometry().elementMapper().index(outside);
259 if ( std::find(handledNeighbors.begin(), handledNeighbors.end(), outsideIdx) == handledNeighbors.end() )
261 makeNeighborGeometries_(outside, outsideIdx);
262 handledNeighbors.push_back(outsideIdx);
269 void makeElementGeometries_()
271 const auto&
element = *element_;
272 const auto eIdx = gridGeometry().elementMapper().index(element);
273 scvs_[0] = SubControlVolume(
element.geometry(), eIdx);
274 scvIndices_[0] = eIdx;
276 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
277 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
279 typename GridGeometry::GeometryHelper geometryHelper(element, gridGeometry().gridView());
282 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
284 const auto& scvfNeighborVolVarIndex = neighborVolVarIndices[scvfCounter];
286 if (intersection.neighbor() || intersection.boundary())
288 geometryHelper.updateLocalFace(gridGeometry().intersectionMapper(), intersection);
289 std::vector<GridIndexType> scvIndices{eIdx, scvfNeighborVolVarIndex};
290 scvfs_.emplace_back(intersection,
291 intersection.geometry(),
292 scvFaceIndices[scvfCounter],
295 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
298 if (intersection.boundary())
299 hasBoundaryScvf_ =
true;
305 void makeNeighborGeometries_(
const Element& element,
const GridIndexType eIdx)
310 neighborScvs_.emplace_back(
element.geometry(), eIdx);
311 neighborScvIndices_.push_back(eIdx);
313 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
314 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
316 typename GridGeometry::GeometryHelper geometryHelper(element, gridGeometry().gridView());
319 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
321 const auto& scvfNeighborVolVarIndex = neighborVolVarIndices[scvfCounter];
322 geometryHelper.updateLocalFace(gridGeometry().intersectionMapper(), intersection);
324 if (intersection.neighbor())
327 if (intersection.outside() == *element_)
329 std::vector<GridIndexType> scvIndices{eIdx, scvfNeighborVolVarIndex};
330 neighborScvfs_.emplace_back(intersection,
331 intersection.geometry(),
332 scvFaceIndices[scvfCounter],
336 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
340 else if (intersection.boundary())
345 const LocalIndexType findLocalIndex_(
const GridIndexType idx,
346 const std::vector<GridIndexType>& indices)
const
348 auto it = std::find(indices.begin(), indices.end(), idx);
349 assert(it != indices.end() &&
"Could not find the scv/scvf! Make sure to properly bind this class!");
356 scvfIndices_.clear();
359 neighborScvIndices_.clear();
360 neighborScvfIndices_.clear();
361 neighborScvs_.clear();
362 neighborScvfs_.clear();
364 hasBoundaryScvf_ =
false;
367 std::optional<Element> element_;
368 const GridGeometry* gridGeometryPtr_;
371 std::array<GridIndexType, 1> scvIndices_;
372 std::array<SubControlVolume, 1> scvs_;
374 std::vector<GridIndexType> scvfIndices_;
375 std::vector<SubControlVolumeFace> scvfs_;
377 std::vector<GridIndexType> neighborScvIndices_;
378 std::vector<SubControlVolume> neighborScvs_;
380 std::vector<GridIndexType> neighborScvfIndices_;
381 std::vector<SubControlVolumeFace> neighborScvfs_;
383 bool hasBoundaryScvf_ =
false;
Defines the index types used for grid and local indices.
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:294
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
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:52
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:62
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:69
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:73
Stencil-local finite volume geometry (scvs and scvfs) for staggered models This builds up the sub con...
Definition: discretization/staggered/fvelementgeometry.hh:43
const SubControlVolumeFace & scvf(GridIndexType eIdx, LocalIndexType localScvfIdx) const
Definition: discretization/staggered/fvelementgeometry.hh:76
StaggeredFVElementGeometry(const CellCenterOrFaceFVGridGeometry &gridGeometry)
Definition: discretization/staggered/fvelementgeometry.hh:71
Base class for the finite volume geometry vector for staggered models This locally builds up the sub ...
Definition: discretization/staggered/fvelementgeometry.hh:92
StaggeredFVElementGeometry bind(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/staggered/fvelementgeometry.hh:181
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/staggered/fvelementgeometry.hh:99
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/staggered/fvelementgeometry.hh:219
void bindElement(const Element &element) &
bind this local view to a specific element
Definition: discretization/staggered/fvelementgeometry.hh:203
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/staggered/fvelementgeometry.hh:162
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/staggered/fvelementgeometry.hh:101
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Definition: discretization/staggered/fvelementgeometry.hh:222
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/staggered/fvelementgeometry.hh:125
StaggeredFVElementGeometry(const CellCenterOrFaceFVGridGeometry &gridGeometry)
Definition: discretization/staggered/fvelementgeometry.hh:110
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/staggered/fvelementgeometry.hh:169
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:118
const GridGeometry & gridGeometry() const
The grid finite volume geometry we are a restriction of.
Definition: discretization/staggered/fvelementgeometry.hh:215
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/staggered/fvelementgeometry.hh:173
StaggeredFVElementGeometry bindElement(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/staggered/fvelementgeometry.hh:196
StaggeredFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/staggered/fvelementgeometry.hh:114
const Element & element() const
The bound element.
Definition: discretization/staggered/fvelementgeometry.hh:211
void bind(const Element &element) &
bind this local view to a specific element (full stencil)
Definition: discretization/staggered/fvelementgeometry.hh:188
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/staggered/fvelementgeometry.hh:207
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/staggered/fvelementgeometry.hh:103
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/staggered/fvelementgeometry.hh:150
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/staggered/fvelementgeometry.hh:135
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/staggered/fvelementgeometry.hh:105
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models This builds up th...