version 3.8
discretization/staggered/subcontrolvolumeface.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_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
23
24#include <typeinfo>
25
26namespace Dumux {
27
32template<class GridView>
34{
35 using Element = typename GridView::template Codim<0>::Entity;
36 using Intersection = typename GridView::Intersection;
37 static constexpr int codimIntersection = 1;
38
39public:
40
41 BaseStaggeredGeometryHelper(const Element& element, const GridView& gridView)
42 : element_(element)
43 , gridView_(gridView)
44 { }
45
49 template<class IntersectionMapper>
50 void updateLocalFace(const IntersectionMapper& intersectionMapper, const Intersection& intersection)
51 {
52 intersection_ = intersection;
53 }
54
58 int dofIndex() const
59 {
60 //TODO: use proper intersection mapper!
61 const auto inIdx = intersection_.indexInInside();
62 return gridView_.indexSet().subIndex(intersection_.inside(), inIdx, codimIntersection);
63 }
64
68 int localFaceIndex() const
69 {
70 return intersection_.indexInInside();
71 }
72
73private:
74 Intersection intersection_;
75 const Element element_;
76 const GridView gridView_;
77};
78
79
86template<class GridView>
88{
91 using Scalar = typename GridView::ctype;
92 using Element = typename GridView::template Codim<0>::Entity;
93 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
94
95 static constexpr int dim = GridView::Grid::dimension;
96 static constexpr int dimWorld = GridView::Grid::dimensionworld;
97 using Geometry = Dune::AxisAlignedCubeGeometry<Scalar, dim-1, dimWorld>;
98};
99
105template<class GV,
108: public SubControlVolumeFaceBase<StaggeredSubControlVolumeFace<GV, T>, T>
109{
112 using Geometry = typename T::Geometry;
113 using GridIndexType = typename T::GridIndexType;
114 using LocalIndexType = typename T::LocalIndexType;
115
116 using Scalar = typename T::Scalar;
117 static const int dim = Geometry::mydimension;
118 static const int dimworld = Geometry::coorddimension;
119
120public:
121 using GlobalPosition = typename T::GlobalPosition;
122
124 using Traits = T;
125
126 // the default constructor
128
130 template <class Intersection, class GeometryHelper>
131 StaggeredSubControlVolumeFace(const Intersection& is,
132 const typename Intersection::Geometry& isGeometry,
133 GridIndexType scvfIndex,
134 const std::vector<GridIndexType>& scvIndices,
135 const GeometryHelper& geometryHelper)
136 : ParentType()
137 , area_(isGeometry.volume())
138 , center_(isGeometry.center())
139 , unitOuterNormal_(is.centerUnitOuterNormal())
140 , scvfIndex_(scvfIndex)
141 , scvIndices_(scvIndices)
142 , boundary_(is.boundary())
143 {
144 dofIdx_ = geometryHelper.dofIndex();
145 localFaceIdx_ = geometryHelper.localFaceIndex();
146 }
147
149 const GlobalPosition& center() const
150 {
151 return center_;
152 }
153
156 {
157 return center_;
158 }
159
162 {
163 // Return center for now
164 return center_;
165 }
166
168 Scalar area() const
169 {
170 return area_;
171 }
172
174 bool boundary() const
175 {
176 return boundary_;
177 }
178
181 {
182 return unitOuterNormal_;
183 }
184
186 GridIndexType insideScvIdx() const
187 {
188 return scvIndices_[0];
189 }
190
192 GridIndexType outsideScvIdx() const
193 {
194 return scvIndices_[1];
195 }
196
198 GridIndexType index() const
199 {
200 return scvfIndex_;
201 }
202
204 GridIndexType dofIndex() const
205 {
206 return dofIdx_;
207 }
208
210 LocalIndexType localFaceIdx() const
211 {
212 return localFaceIdx_;
213 }
214
215private:
216 Scalar area_;
217 GlobalPosition center_;
218 GlobalPosition unitOuterNormal_;
219 GridIndexType scvfIndex_;
220 std::vector<GridIndexType> scvIndices_;
221 bool boundary_;
222
223 GridIndexType dofIdx_;
224 LocalIndexType localFaceIdx_;
225};
226
227} // end namespace Dumux
228
229#endif
Base class for a staggered grid geometry helper.
Definition: discretization/staggered/subcontrolvolumeface.hh:34
void updateLocalFace(const IntersectionMapper &intersectionMapper, const Intersection &intersection)
Updates the current face, i.e. sets the correct intersection.
Definition: discretization/staggered/subcontrolvolumeface.hh:50
BaseStaggeredGeometryHelper(const Element &element, const GridView &gridView)
Definition: discretization/staggered/subcontrolvolumeface.hh:41
int dofIndex() const
Returns the global dofIdx of the intersection itself.
Definition: discretization/staggered/subcontrolvolumeface.hh:58
int localFaceIndex() const
Returns the local index of the face (i.e. the intersection)
Definition: discretization/staggered/subcontrolvolumeface.hh:68
defines an intersection mapper for mapping of global DOFs assigned to faces which also works for adap...
Definition: intersectionmapper.hh:200
Class for a sub control volume face in the staggered method, i.e a part of the boundary of a sub cont...
Definition: discretization/staggered/subcontrolvolumeface.hh:109
const GlobalPosition & ipGlobal() const
The integration point for flux evaluations in global coordinates.
Definition: discretization/staggered/subcontrolvolumeface.hh:161
GridIndexType index() const
The global index of this sub control volume face.
Definition: discretization/staggered/subcontrolvolumeface.hh:198
LocalIndexType localFaceIdx() const
The local index of this sub control volume face.
Definition: discretization/staggered/subcontrolvolumeface.hh:210
Scalar area() const
The area of the sub control volume face.
Definition: discretization/staggered/subcontrolvolumeface.hh:168
T Traits
state the traits public and thus export all types
Definition: discretization/staggered/subcontrolvolumeface.hh:124
StaggeredSubControlVolumeFace(const Intersection &is, const typename Intersection::Geometry &isGeometry, GridIndexType scvfIndex, const std::vector< GridIndexType > &scvIndices, const GeometryHelper &geometryHelper)
Constructor with intersection.
Definition: discretization/staggered/subcontrolvolumeface.hh:131
GridIndexType insideScvIdx() const
Index of the inside sub control volume.
Definition: discretization/staggered/subcontrolvolumeface.hh:186
const GlobalPosition & dofPosition() const
The position of the dof living on the face.
Definition: discretization/staggered/subcontrolvolumeface.hh:155
const GlobalPosition & unitOuterNormal() const
The unit outer normal vector.
Definition: discretization/staggered/subcontrolvolumeface.hh:180
const GlobalPosition & center() const
The center of the sub control volume face.
Definition: discretization/staggered/subcontrolvolumeface.hh:149
bool boundary() const
Returns boolean if the sub control volume face is on the boundary.
Definition: discretization/staggered/subcontrolvolumeface.hh:174
typename T::GlobalPosition GlobalPosition
Definition: discretization/staggered/subcontrolvolumeface.hh:121
GridIndexType dofIndex() const
The global index of the dof living on this face.
Definition: discretization/staggered/subcontrolvolumeface.hh:204
GridIndexType outsideScvIdx() const
Index of the outside sub control volume.
Definition: discretization/staggered/subcontrolvolumeface.hh:192
Base class for a sub control volume face, i.e a part of the boundary of a sub control volume we compu...
Definition: subcontrolvolumefacebase.hh:29
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:159
Defines the index types used for grid and local indices.
Definition: adapt.hh:17
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:27
unsigned int LocalIndex
Definition: indextraits.hh:28
Default traits class to be used for the sub-control volume faces for the staggered finite volume sche...
Definition: discretization/staggered/subcontrolvolumeface.hh:88
typename GridView::template Codim< 0 >::Entity Element
Definition: discretization/staggered/subcontrolvolumeface.hh:92
static constexpr int dim
Definition: discretization/staggered/subcontrolvolumeface.hh:95
typename IndexTraits< GridView >::GridIndex GridIndexType
Definition: discretization/staggered/subcontrolvolumeface.hh:89
typename GridView::ctype Scalar
Definition: discretization/staggered/subcontrolvolumeface.hh:91
typename IndexTraits< GridView >::LocalIndex LocalIndexType
Definition: discretization/staggered/subcontrolvolumeface.hh:90
Dune::AxisAlignedCubeGeometry< Scalar, dim-1, dimWorld > Geometry
Definition: discretization/staggered/subcontrolvolumeface.hh:97
static constexpr int dimWorld
Definition: discretization/staggered/subcontrolvolumeface.hh:96
typename Element::Geometry::GlobalCoordinate GlobalPosition
Definition: discretization/staggered/subcontrolvolumeface.hh:93
Base class for a sub control volume face.