Line | Branch | Exec | Source |
---|---|---|---|
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-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder | ||
5 | // SPDX-License-Identifier: GPL-3.0-or-later | ||
6 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup StaggeredDiscretization | ||
10 | * \copydoc Dumux::StaggeredSubControlVolumeFace | ||
11 | */ | ||
12 | #ifndef DUMUX_DISCRETIZATION_STAGGERED_SUBCONTROLVOLUMEFACE_HH | ||
13 | #define DUMUX_DISCRETIZATION_STAGGERED_SUBCONTROLVOLUMEFACE_HH | ||
14 | |||
15 | #include <utility> | ||
16 | |||
17 | #include <dune/geometry/axisalignedcubegeometry.hh> | ||
18 | #include <dune/common/fvector.hh> | ||
19 | #include <dune/geometry/type.hh> | ||
20 | |||
21 | #include <dumux/common/indextraits.hh> | ||
22 | |||
23 | #include <typeinfo> | ||
24 | |||
25 | namespace Dumux { | ||
26 | |||
27 | /*! | ||
28 | * \ingroup StaggeredDiscretization | ||
29 | * \brief Base class for a staggered grid geometry helper | ||
30 | */ | ||
31 | template<class GridView> | ||
32 | class BaseStaggeredGeometryHelper | ||
33 | { | ||
34 | using Element = typename GridView::template Codim<0>::Entity; | ||
35 | using Intersection = typename GridView::Intersection; | ||
36 | static constexpr int codimIntersection = 1; | ||
37 | |||
38 | public: | ||
39 | |||
40 | 8 | BaseStaggeredGeometryHelper(const Element& element, const GridView& gridView) | |
41 | 8 | : element_(element) | |
42 | 16 | , gridView_(gridView) | |
43 | { } | ||
44 | |||
45 | /*! | ||
46 | * \brief Updates the current face, i.e. sets the correct intersection | ||
47 | */ | ||
48 | template<class IntersectionMapper> | ||
49 | 32 | void updateLocalFace(const IntersectionMapper& intersectionMapper, const Intersection& intersection) | |
50 | { | ||
51 |
2/2✓ Branch 0 taken 26 times.
✓ Branch 1 taken 6 times.
|
32 | intersection_ = intersection; |
52 | } | ||
53 | |||
54 | /*! | ||
55 | * \brief Returns the global dofIdx of the intersection itself | ||
56 | */ | ||
57 | 32 | int dofIndex() const | |
58 | { | ||
59 | //TODO: use proper intersection mapper! | ||
60 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 32 times.
|
32 | const auto inIdx = intersection_.indexInInside(); |
61 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 32 times.
|
32 | return gridView_.indexSet().subIndex(intersection_.inside(), inIdx, codimIntersection); |
62 | } | ||
63 | |||
64 | /*! | ||
65 | * \brief Returns the local index of the face (i.e. the intersection) | ||
66 | */ | ||
67 | 64 | int localFaceIndex() const | |
68 | { | ||
69 |
2/2✓ Branch 0 taken 26 times.
✓ Branch 1 taken 6 times.
|
32 | return intersection_.indexInInside(); |
70 | } | ||
71 | |||
72 | private: | ||
73 | Intersection intersection_; //!< The intersection of interest | ||
74 | const Element element_; //!< The respective element | ||
75 | const GridView gridView_; | ||
76 | }; | ||
77 | |||
78 | |||
79 | /*! | ||
80 | * \ingroup StaggeredDiscretization | ||
81 | * \brief Default traits class to be used for the sub-control volume faces | ||
82 | * for the staggered finite volume scheme | ||
83 | * \tparam GV the type of the grid view | ||
84 | */ | ||
85 | template<class GridView> | ||
86 | struct StaggeredDefaultScvfGeometryTraits | ||
87 | { | ||
88 | using GridIndexType = typename IndexTraits<GridView>::GridIndex; | ||
89 | using LocalIndexType = typename IndexTraits<GridView>::LocalIndex; | ||
90 | using Scalar = typename GridView::ctype; | ||
91 | using Element = typename GridView::template Codim<0>::Entity; | ||
92 | using GlobalPosition = typename Element::Geometry::GlobalCoordinate; | ||
93 | |||
94 | static constexpr int dim = GridView::Grid::dimension; | ||
95 | static constexpr int dimWorld = GridView::Grid::dimensionworld; | ||
96 | using Geometry = Dune::AxisAlignedCubeGeometry<Scalar, dim-1, dimWorld>; | ||
97 | }; | ||
98 | |||
99 | /*! | ||
100 | * \ingroup StaggeredDiscretization | ||
101 | * \brief Class for a sub control volume face in the staggered method, i.e a part of the boundary | ||
102 | * of a sub control volume we compute fluxes on. | ||
103 | */ | ||
104 | template<class GV, | ||
105 | class T = StaggeredDefaultScvfGeometryTraits<GV> > | ||
106 | 32 | class StaggeredSubControlVolumeFace | |
107 | { | ||
108 | using ThisType = StaggeredSubControlVolumeFace<GV, T>; | ||
109 | using Geometry = typename T::Geometry; | ||
110 | using GridIndexType = typename T::GridIndexType; | ||
111 | using LocalIndexType = typename T::LocalIndexType; | ||
112 | |||
113 | using Scalar = typename T::Scalar; | ||
114 | static const int dim = Geometry::mydimension; | ||
115 | static const int dimworld = Geometry::coorddimension; | ||
116 | |||
117 | public: | ||
118 | using GlobalPosition = typename T::GlobalPosition; | ||
119 | |||
120 | //! state the traits public and thus export all types | ||
121 | using Traits = T; | ||
122 | |||
123 | // the default constructor | ||
124 | StaggeredSubControlVolumeFace() = default; | ||
125 | |||
126 | //! Constructor with intersection | ||
127 | template <class Intersection, class GeometryHelper> | ||
128 | 32 | StaggeredSubControlVolumeFace(const Intersection& is, | |
129 | const typename Intersection::Geometry& isGeometry, | ||
130 | GridIndexType scvfIndex, | ||
131 | const std::vector<GridIndexType>& scvIndices, | ||
132 | const GeometryHelper& geometryHelper) | ||
133 | 32 | : area_(isGeometry.volume()) | |
134 |
2/2✓ Branch 0 taken 16 times.
✓ Branch 1 taken 16 times.
|
32 | , center_(isGeometry.center()) |
135 | 32 | , unitOuterNormal_(is.centerUnitOuterNormal()) | |
136 | 32 | , scvfIndex_(scvfIndex) | |
137 | 32 | , scvIndices_(scvIndices) | |
138 | 32 | , boundary_(is.boundary()) | |
139 | { | ||
140 | 32 | dofIdx_ = geometryHelper.dofIndex(); | |
141 | 32 | localFaceIdx_ = geometryHelper.localFaceIndex(); | |
142 | 32 | } | |
143 | |||
144 | //! The center of the sub control volume face | ||
145 | const GlobalPosition& center() const | ||
146 | { | ||
147 | return center_; | ||
148 | } | ||
149 | |||
150 | //! The position of the dof living on the face | ||
151 | const GlobalPosition& dofPosition() const | ||
152 | { | ||
153 | return center_; | ||
154 | } | ||
155 | |||
156 | //! The integration point for flux evaluations in global coordinates | ||
157 | const GlobalPosition& ipGlobal() const | ||
158 | { | ||
159 | // Return center for now | ||
160 |
1/2✓ Branch 1 taken 32 times.
✗ Branch 2 not taken.
|
32 | return center_; |
161 | } | ||
162 | |||
163 | //! The area of the sub control volume face | ||
164 | Scalar area() const | ||
165 | { | ||
166 | return area_; | ||
167 | } | ||
168 | |||
169 | //! Returns boolean if the sub control volume face is on the boundary | ||
170 | 32 | bool boundary() const | |
171 | { | ||
172 |
2/2✓ Branch 0 taken 12 times.
✓ Branch 1 taken 20 times.
|
32 | return boundary_; |
173 | } | ||
174 | |||
175 | //! The unit outer normal vector | ||
176 | 32 | const GlobalPosition& unitOuterNormal() const | |
177 | { | ||
178 |
1/2✓ Branch 1 taken 32 times.
✗ Branch 2 not taken.
|
32 | return unitOuterNormal_; |
179 | } | ||
180 | |||
181 | //! Index of the inside sub control volume | ||
182 | GridIndexType insideScvIdx() const | ||
183 | { | ||
184 | return scvIndices_[0]; | ||
185 | } | ||
186 | |||
187 | //! Index of the outside sub control volume | ||
188 | GridIndexType outsideScvIdx() const | ||
189 | { | ||
190 | return scvIndices_[1]; | ||
191 | } | ||
192 | |||
193 | //! The global index of this sub control volume face | ||
194 | 32 | GridIndexType index() const | |
195 | { | ||
196 |
1/2✓ Branch 1 taken 32 times.
✗ Branch 2 not taken.
|
32 | return scvfIndex_; |
197 | } | ||
198 | |||
199 | //! The global index of the dof living on this face | ||
200 | GridIndexType dofIndex() const | ||
201 | { | ||
202 | return dofIdx_; | ||
203 | } | ||
204 | |||
205 | //! The local index of this sub control volume face | ||
206 | LocalIndexType localFaceIdx() const | ||
207 | { | ||
208 | return localFaceIdx_; | ||
209 | } | ||
210 | |||
211 | private: | ||
212 | Scalar area_; | ||
213 | GlobalPosition center_; | ||
214 | GlobalPosition unitOuterNormal_; | ||
215 | GridIndexType scvfIndex_; | ||
216 | std::vector<GridIndexType> scvIndices_; | ||
217 | bool boundary_; | ||
218 | |||
219 | GridIndexType dofIdx_; | ||
220 | LocalIndexType localFaceIdx_; | ||
221 | }; | ||
222 | |||
223 | } // end namespace Dumux | ||
224 | |||
225 | #endif | ||
226 |