24#ifndef DUMUX_DISCRETIZATION_STAGGERED_SUBCONTROLVOLUMEFACE_HH
25#define DUMUX_DISCRETIZATION_STAGGERED_SUBCONTROLVOLUMEFACE_HH
29#include <dune/common/fvector.hh>
30#include <dune/geometry/type.hh>
43template<
class Gr
idView>
46 using Element =
typename GridView::template Codim<0>::Entity;
47 using Intersection =
typename GridView::Intersection;
48 static constexpr int codimIntersection = 1;
60 template<
class IntersectionMapper>
63 intersection_ = intersection;
72 const auto inIdx = intersection_.indexInInside();
73 return gridView_.indexSet().subIndex(intersection_.inside(), inIdx, codimIntersection);
81 return intersection_.indexInInside();
85 Intersection intersection_;
86 const Element element_;
87 const GridView gridView_;
97template<
class Gr
idView>
100 using Geometry =
typename GridView::template Codim<1>::Geometry;
119 using Geometry =
typename T::Geometry;
120 using GridIndexType =
typename T::GridIndexType;
121 using LocalIndexType =
typename T::LocalIndexType;
123 using Scalar =
typename T::Scalar;
124 static const int dim = Geometry::mydimension;
125 static const int dimworld = Geometry::coorddimension;
137 template <
class Intersection,
class GeometryHelper>
139 const typename Intersection::Geometry& isGeometry,
140 GridIndexType scvfIndex,
141 const std::vector<GridIndexType>& scvIndices,
142 const GeometryHelper& geometryHelper)
144 , geomType_(isGeometry.type())
145 , area_(isGeometry.volume())
146 , center_(isGeometry.
center())
147 , unitOuterNormal_(is.centerUnitOuterNormal())
148 , scvfIndex_(scvfIndex)
149 , scvIndices_(scvIndices)
152 corners_.resize(isGeometry.corners());
153 for (
int i = 0; i < isGeometry.corners(); ++i)
154 corners_[i] = isGeometry.corner(i);
156 dofIdx_ = geometryHelper.dofIndex();
157 localFaceIdx_ = geometryHelper.localFaceIndex();
194 return unitOuterNormal_;
200 return scvIndices_[0];
206 return scvIndices_[1];
218 assert(localIdx < corners_.size() &&
"provided index exceeds the number of corners");
219 return corners_[localIdx];
225 return Geometry(geomType_, corners_);
237 return localFaceIdx_;
241 Dune::GeometryType geomType_;
242 std::vector<GlobalPosition> corners_;
246 GridIndexType scvfIndex_;
247 std::vector<GridIndexType> scvIndices_;
250 GridIndexType dofIdx_;
251 LocalIndexType localFaceIdx_;
Defines the index types used for grid and local indices.
Base class for a sub control volume face.
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:39
unsigned int LocalIndex
Definition: indextraits.hh:40
defines an intersection mapper for mapping of global DOFs assigned to faces which also works for adap...
Definition: intersectionmapper.hh:185
Base class for a staggered grid geometry helper.
Definition: discretization/staggered/subcontrolvolumeface.hh:45
void updateLocalFace(const IntersectionMapper &intersectionMapper, const Intersection &intersection)
Updates the current face, i.e. sets the correct intersection.
Definition: discretization/staggered/subcontrolvolumeface.hh:61
BaseStaggeredGeometryHelper(const Element &element, const GridView &gridView)
Definition: discretization/staggered/subcontrolvolumeface.hh:52
int dofIndex() const
Returns the global dofIdx of the intersection itself.
Definition: discretization/staggered/subcontrolvolumeface.hh:69
int localFaceIndex() const
Returns the local index of the face (i.e. the intersection)
Definition: discretization/staggered/subcontrolvolumeface.hh:79
Default traits class to be used for the sub-control volume faces for the staggered finite volume sche...
Definition: discretization/staggered/subcontrolvolumeface.hh:99
Dune::FieldVector< Scalar, GridView::dimensionworld > GlobalPosition
Definition: discretization/staggered/subcontrolvolumeface.hh:104
typename GridView::template Codim< 1 >::Geometry Geometry
Definition: discretization/staggered/subcontrolvolumeface.hh:100
typename IndexTraits< GridView >::GridIndex GridIndexType
Definition: discretization/staggered/subcontrolvolumeface.hh:101
typename GridView::ctype Scalar
Definition: discretization/staggered/subcontrolvolumeface.hh:103
typename IndexTraits< GridView >::LocalIndex LocalIndexType
Definition: discretization/staggered/subcontrolvolumeface.hh:102
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:116
const GlobalPosition & ipGlobal() const
The integration point for flux evaluations in global coordinates.
Definition: discretization/staggered/subcontrolvolumeface.hh:173
GridIndexType index() const
The global index of this sub control volume face.
Definition: discretization/staggered/subcontrolvolumeface.hh:210
LocalIndexType localFaceIdx() const
The local index of this sub control volume face.
Definition: discretization/staggered/subcontrolvolumeface.hh:235
Scalar area() const
The area of the sub control volume face.
Definition: discretization/staggered/subcontrolvolumeface.hh:180
T Traits
state the traits public and thus export all types
Definition: discretization/staggered/subcontrolvolumeface.hh:131
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:138
GridIndexType insideScvIdx() const
Index of the inside sub control volume.
Definition: discretization/staggered/subcontrolvolumeface.hh:198
const GlobalPosition & dofPosition() const
The position of the dof living on the face.
Definition: discretization/staggered/subcontrolvolumeface.hh:167
const GlobalPosition & unitOuterNormal() const
The unit outer normal vector.
Definition: discretization/staggered/subcontrolvolumeface.hh:192
const GlobalPosition & center() const
The center of the sub control volume face.
Definition: discretization/staggered/subcontrolvolumeface.hh:161
bool boundary() const
Returns bolean if the sub control volume face is on the boundary.
Definition: discretization/staggered/subcontrolvolumeface.hh:186
const Geometry geometry() const
The geometry of the sub control volume face.
Definition: discretization/staggered/subcontrolvolumeface.hh:223
typename T::GlobalPosition GlobalPosition
Definition: discretization/staggered/subcontrolvolumeface.hh:128
GridIndexType dofIndex() const
The global index of the dof living on this face.
Definition: discretization/staggered/subcontrolvolumeface.hh:229
GridIndexType outsideScvIdx() const
Index of the outside sub control volume.
Definition: discretization/staggered/subcontrolvolumeface.hh:204
StaggeredSubControlVolumeFace()=default
const GlobalPosition & corner(unsigned int localIdx) const
The positions of the corners.
Definition: discretization/staggered/subcontrolvolumeface.hh:216
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:41