12#ifndef DUMUX_DISCRETIZATION_STAGGERED_SUBCONTROLVOLUMEFACE_HH 
   13#define DUMUX_DISCRETIZATION_STAGGERED_SUBCONTROLVOLUMEFACE_HH 
   17#include <dune/geometry/axisalignedcubegeometry.hh> 
   18#include <dune/common/fvector.hh> 
   19#include <dune/geometry/type.hh> 
   32template<
class Gr
idView>
 
   35    using Element = 
typename GridView::template Codim<0>::Entity;
 
   36    using Intersection = 
typename GridView::Intersection;
 
   37    static constexpr int codimIntersection =  1;
 
   49    template<
class IntersectionMapper>
 
   52        intersection_ = intersection;
 
 
   61       const auto inIdx = intersection_.indexInInside();
 
   62       return gridView_.indexSet().subIndex(intersection_.inside(), inIdx, codimIntersection);
 
 
   70       return intersection_.indexInInside();
 
 
   74   Intersection intersection_; 
 
   75   const Element element_; 
 
   76   const GridView gridView_;
 
 
   86template<
class Gr
idView>
 
   91    using Scalar = 
typename GridView::ctype;
 
   92    using Element = 
typename GridView::template Codim<0>::Entity;
 
   95    static constexpr int dim = GridView::Grid::dimension;
 
   96    static constexpr int dimWorld = GridView::Grid::dimensionworld;
 
 
  112    using Geometry = 
typename T::Geometry;
 
  113    using GridIndexType = 
typename T::GridIndexType;
 
  114    using LocalIndexType = 
typename T::LocalIndexType;
 
  116    using Scalar = 
typename T::Scalar;
 
  117    static const int dim = Geometry::mydimension;
 
  118    static const int dimworld = Geometry::coorddimension;
 
  130    template <
class Intersection, 
class GeometryHelper>
 
  132                                  const typename Intersection::Geometry& isGeometry,
 
  133                                  GridIndexType scvfIndex,
 
  134                                  const std::vector<GridIndexType>& scvIndices,
 
  135                                  const GeometryHelper& geometryHelper)
 
  137    , area_(isGeometry.
volume())
 
  138    , center_(isGeometry.
center())
 
  139    , unitOuterNormal_(is.centerUnitOuterNormal())
 
  140    , scvfIndex_(scvfIndex)
 
  141    , scvIndices_(scvIndices)
 
  144        dofIdx_ = geometryHelper.dofIndex();
 
  145        localFaceIdx_ = geometryHelper.localFaceIndex();
 
 
  182        return unitOuterNormal_;
 
 
  188        return scvIndices_[0];
 
 
  194        return scvIndices_[1];
 
 
  212        return localFaceIdx_;
 
 
  219    GridIndexType scvfIndex_;
 
  220    std::vector<GridIndexType> scvIndices_;
 
  223    GridIndexType dofIdx_;
 
  224    LocalIndexType localFaceIdx_;
 
 
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
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
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
Definition discretization/staggered/subcontrolvolumeface.hh:149
bool boundary() const
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
StaggeredSubControlVolumeFace()=default
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
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.