28#ifndef EWOMS_FV_BASE_BOUNDARY_CONTEXT_HH
29#define EWOMS_FV_BASE_BOUNDARY_CONTEXT_HH
31#include <dune/common/fvector.hh>
44template<
class TypeTag>
57 using Element =
typename GridView::template Codim<0>::Entity;
58 using IntersectionIterator =
typename GridView::IntersectionIterator;
59 using Intersection =
typename GridView::Intersection;
61 enum { dim = GridView::dimension };
62 enum { dimWorld = GridView::dimensionworld };
64 using CoordScalar =
typename GridView::ctype;
65 using GlobalPosition = Dune::FieldVector<CoordScalar, dimWorld>;
66 using Vector = Dune::FieldVector<Scalar, dimWorld>;
81 if (intersectionIt_ == iend) {
87 while (intersectionIt_ != iend && intersectionIt_->neighbor()) {
96 {
return elemCtx_.problem(); }
102 {
return elemCtx_.model(); }
108 {
return elemCtx_.gridView(); }
114 {
return elemCtx_.element(); }
126 {
return elemCtx_.gradientCalculator(); }
131 std::size_t
numDof(
unsigned timeIdx)
const
132 {
return elemCtx_.numDof(timeIdx); }
138 {
return elemCtx_.numPrimaryDof(timeIdx); }
144 {
return elemCtx_.numInteriorFaces(timeIdx); }
150 {
return elemCtx_.stencil(timeIdx).numBoundaryFaces(); }
155 const Stencil&
stencil(
unsigned timeIdx)
const
156 {
return elemCtx_.stencil(timeIdx); }
164 Vector
normal(
unsigned boundaryFaceIdx,
unsigned timeIdx)
const
166 auto tmp =
stencil(timeIdx).boundaryFace[boundaryFaceIdx].normal;
167 tmp /= tmp.two_norm();
175 {
return elemCtx_.stencil(timeIdx).boundaryFace(boundaryFaceIdx).area(); }
183 const GlobalPosition&
pos(
unsigned boundaryFaceIdx,
unsigned timeIdx)
const
184 {
return stencil(timeIdx).boundaryFace(boundaryFaceIdx).integrationPos(); }
192 const GlobalPosition&
cvCenter(
unsigned boundaryFaceIdx,
unsigned timeIdx)
const
194 const unsigned scvIdx =
stencil(timeIdx).boundaryFace(boundaryFaceIdx).interiorIndex();
195 return stencil(timeIdx).subControlVolume(scvIdx).globalPos();
203 {
return elemCtx_.focusDofIndex(); }
213 {
return stencil(timeIdx).boundaryFace(boundaryFaceIdx).interiorIndex(); }
223 {
return elemCtx_.globalSpaceIndex(
interiorScvIndex(boundaryFaceIdx, timeIdx), timeIdx); }
234 const unsigned interiorScvIdx = this->
interiorScvIndex(boundaryFaceIdx, timeIdx);
235 return elemCtx_.intensiveQuantities(interiorScvIdx, timeIdx);
245 {
return elemCtx_.boundaryExtensiveQuantities(boundaryFaceIdx, timeIdx); }
257 {
return *intersectionIt_; }
268 {
return intersectionIt_; }
271 const ElementContext& elemCtx_;
272 IntersectionIterator intersectionIt_;
const GlobalPosition & pos(unsigned boundaryFaceIdx, unsigned timeIdx) const
Return the position of a local entity in global coordinates.
Definition fvbaseboundarycontext.hh:183
const GlobalPosition & cvCenter(unsigned boundaryFaceIdx, unsigned timeIdx) const
Return the position of a control volume's center in global coordinates.
Definition fvbaseboundarycontext.hh:192
Vector normal(unsigned boundaryFaceIdx, unsigned timeIdx) const
Returns the outer unit normal of the boundary segment.
Definition fvbaseboundarycontext.hh:164
unsigned globalSpaceIndex(unsigned boundaryFaceIdx, unsigned timeIdx) const
Return the global space index of the sub-control volume at the interior of a boundary segment.
Definition fvbaseboundarycontext.hh:222
const Stencil & stencil(unsigned timeIdx) const
Definition fvbaseboundarycontext.hh:155
FvBaseBoundaryContext(const ElementContext &elemCtx)
The constructor.
Definition fvbaseboundarycontext.hh:72
Scalar boundarySegmentArea(unsigned boundaryFaceIdx, unsigned timeIdx) const
Returns the area [m^2] of a given boudary segment.
Definition fvbaseboundarycontext.hh:174
unsigned focusDofIndex() const
Return the local sub-control volume index upon which the linearization is currently focused.
Definition fvbaseboundarycontext.hh:202
const Intersection intersection(unsigned) const
Return the intersection for the neumann segment.
Definition fvbaseboundarycontext.hh:256
const GridView & gridView() const
Definition fvbaseboundarycontext.hh:107
const ElementContext & elementContext() const
Returns a reference to the element context object.
Definition fvbaseboundarycontext.hh:119
const IntensiveQuantities & intensiveQuantities(unsigned boundaryFaceIdx, unsigned timeIdx) const
Return the intensive quantities for the finite volume in the interiour of a boundary segment.
Definition fvbaseboundarycontext.hh:232
const Problem & problem() const
Definition fvbaseboundarycontext.hh:95
std::size_t numDof(unsigned timeIdx) const
Definition fvbaseboundarycontext.hh:131
const Element & element() const
Definition fvbaseboundarycontext.hh:113
const GradientCalculator & gradientCalculator() const
Returns a reference to the current gradient calculator.
Definition fvbaseboundarycontext.hh:125
std::size_t numBoundaryFaces(unsigned timeIdx) const
Return the number of boundary segments of the current element.
Definition fvbaseboundarycontext.hh:149
std::size_t numInteriorFaces(unsigned timeIdx) const
Definition fvbaseboundarycontext.hh:143
unsigned interiorScvIndex(unsigned boundaryFaceIdx, unsigned timeIdx) const
Return the local sub-control volume index of the interior of a boundary segment.
Definition fvbaseboundarycontext.hh:212
std::size_t numPrimaryDof(unsigned timeIdx) const
Definition fvbaseboundarycontext.hh:137
const Model & model() const
Definition fvbaseboundarycontext.hh:101
const ExtensiveQuantities & extensiveQuantities(unsigned boundaryFaceIdx, unsigned timeIdx) const
Return the extensive quantities for a given boundary face.
Definition fvbaseboundarycontext.hh:244
IntersectionIterator & intersectionIt()
Return the intersection for the neumann segment.
Definition fvbaseboundarycontext.hh:267
Declare the properties used by the infrastructure code of the finite volume discretizations.
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition propertysystem.hh:233