version 3.10-dev
discretization/facecentered/diamond/fvelementgeometry.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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_DISCRETIZATION_FACECENTERED_DIAMOND_FV_ELEMENT_GEOMETRY_HH
13#define DUMUX_DISCRETIZATION_FACECENTERED_DIAMOND_FV_ELEMENT_GEOMETRY_HH
14
15#include <type_traits>
16#include <optional>
17
18#include <dune/common/reservedvector.hh>
19#include <dune/common/iteratorrange.hh>
20
24
25namespace Dumux {
26
31template<class GG, bool cachingEnabled>
33
38template<class GG>
39class FaceCenteredDiamondFVElementGeometry<GG, /*cachingEnabled*/true>
40{
41 using ThisType = FaceCenteredDiamondFVElementGeometry<GG, /*cachingEnabled*/true>;
42 using GridView = typename GG::GridView;
43 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
44 using LocalIndexType = typename IndexTraits<GridView>::SmallLocalIndex;
45 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
46 using GGCache = typename GG::Cache;
47 using GeometryHelper = typename GGCache::GeometryHelper;
48
49public:
51 using SubControlVolume = typename GG::SubControlVolume;
52 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
53 using Element = typename GridView::template Codim<0>::Entity;
54 using GridGeometry = GG;
55
57 static constexpr std::size_t maxNumElementScvs = 2*GridView::dimension;
58
60 : ggCache_(&ggCache)
61 {}
62
64 const SubControlVolume& scv(LocalIndexType scvIdx) const
65 { return ggCache_->scvs(eIdx_)[scvIdx]; }
66
68 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
69 { return ggCache_->scvf(eIdx_)[scvfIdx]; }
70
76 friend inline auto
78 {
79 using Iter = typename std::vector<SubControlVolume>::const_iterator;
80 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
81 return Dune::IteratorRange<Iter>(s.begin(), s.end());
82 }
83
89 friend inline auto
91 {
92 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
93 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
94 return Dune::IteratorRange<Iter>(s.begin(), s.end());
95 }
96
98 const FeLocalBasis& feLocalBasis() const
99 {
100 return gridGeometry().feCache().get(element().type()).localBasis();
101 }
102
104 std::size_t numScv() const
105 {
106 return ggCache_->scvs(eIdx_).size();
107 }
108
110 std::size_t numScvf() const
111 {
112 return ggCache_->scvfs(eIdx_).size();
113 }
114
116 bool hasBoundaryScvf() const
117 { return ggCache_->hasBoundaryScvf(eIdx_); }
118
125 {
126 this->bindElement(element);
127 return std::move(*this);
128 }
129
130 void bind(const Element& element) &
131 {
132 this->bindElement(element);
133 }
134
141 {
142 this->bindElement(element);
143 return std::move(*this);
144 }
145
147 void bindElement(const Element& element) &
148 {
149 element_ = element;
150 eIdx_ = gridGeometry().elementMapper().index(element);
151 }
152
154 bool isBound() const
155 { return static_cast<bool>(element_); }
156
158 const Element& element() const
159 { return *element_; }
160
163 { return ggCache_->gridGeometry(); }
164
166 std::size_t elementIndex() const
167 { return eIdx_; }
168
170 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
171 {
172 assert(isBound());
173 const auto geo = element().geometry();
174 return {
175 SubControlVolume::Traits::geometryType(geo.type()),
176 GeometryHelper(geo).getScvCorners(scv.indexInElement())
177 };
178 }
179
181 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
182 {
183 assert(isBound());
184 const auto geo = element().geometry();
185 if (scvf.boundary())
186 {
187 // use the information that each boundary scvf corresponds to one scv constructed around the same facet
188 const auto localFacetIndex = scvf.insideScvIdx();
189 return {
190 referenceElement(geo).type(localFacetIndex, 1),
191 GeometryHelper(geo).getBoundaryScvfCorners(localFacetIndex)
192 };
193 }
194 else
195 {
196 return {
197 SubControlVolumeFace::Traits::interiorGeometryType(geo.type()),
198 GeometryHelper(geo).getScvfCorners(scvf.index())
199 };
200 }
201 }
202
203private:
204 std::optional<Element> element_;
205 GridIndexType eIdx_;
206 const GGCache* ggCache_;
207};
208
209} // end namespace Dumux
210
211#endif
Element-wise grid geometry (local view)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:40
friend auto scvs(const FaceCenteredDiamondFVElementGeometry &fvGeometry)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:77
void bindElement(const Element &element) &
Bind only element-local.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:147
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:162
const Element & element() const
The bound element.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:158
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume face
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:51
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:116
typename GG::SubControlVolumeFace SubControlVolumeFace
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:52
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:154
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:110
friend auto scvfs(const FaceCenteredDiamondFVElementGeometry &fvGeometry)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:90
std::size_t elementIndex() const
The bound element index.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:166
GG GridGeometry
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:54
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:68
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:64
typename GridView::template Codim< 0 >::Entity Element
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:53
FaceCenteredDiamondFVElementGeometry 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/facecentered/diamond/fvelementgeometry.hh:124
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:104
FaceCenteredDiamondFVElementGeometry(const GGCache &ggCache)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:59
void bind(const Element &element) &
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:130
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:98
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:181
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:170
FaceCenteredDiamondFVElementGeometry 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/facecentered/diamond/fvelementgeometry.hh:140
Element-wise grid geometry (local view)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:32
Defines the index types used for grid and local indices.
Definition: adapt.hh:17
Class providing iterators over sub control volumes and sub control volume faces of an element.
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:27
std::uint_least8_t SmallLocalIndex
Definition: indextraits.hh:29