fvbaseboundarycontext.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 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef EWOMS_FV_BASE_BOUNDARY_CONTEXT_HH
29#define EWOMS_FV_BASE_BOUNDARY_CONTEXT_HH
30
31#include "fvbaseproperties.hh"
32
33#include <dune/common/fvector.hh>
34
35namespace Opm {
36
42template<class TypeTag>
44{
53
55 using Element = typename GridView::template Codim<0>::Entity;
56 using IntersectionIterator = typename GridView::IntersectionIterator;
57 using Intersection = typename GridView::Intersection;
58
59 enum { dim = GridView::dimension };
60 enum { dimWorld = GridView::dimensionworld };
61
62 using CoordScalar = typename GridView::ctype;
63 using GlobalPosition = Dune::FieldVector<CoordScalar, dimWorld>;
64 using Vector = Dune::FieldVector<Scalar, dimWorld>;
65
66public:
70 explicit FvBaseBoundaryContext(const ElementContext& elemCtx)
71 : elemCtx_(elemCtx)
72 , intersectionIt_(gridView().ibegin(element()))
73 { }
74
75 void increment()
76 {
77 const auto& iend = gridView().iend(element());
78
79 if(intersectionIt_ == iend)
80 return;
81
83 // iterate to the next boundary intersection
84 while (intersectionIt_ != iend && intersectionIt_->neighbor()) {
86 }
87 }
88
92 const Problem& problem() const
93 { return elemCtx_.problem(); }
94
98 const Model& model() const
99 { return elemCtx_.model(); }
100
104 const GridView& gridView() const
105 { return elemCtx_.gridView(); }
106
110 const Element& element() const
111 { return elemCtx_.element(); }
112
116 const ElementContext& elementContext() const
117 { return elemCtx_; }
118
122 const GradientCalculator& gradientCalculator() const
123 { return elemCtx_.gradientCalculator(); }
124
128 size_t numDof(unsigned timeIdx) const
129 { return elemCtx_.numDof(timeIdx); }
130
134 size_t numPrimaryDof(unsigned timeIdx) const
135 { return elemCtx_.numPrimaryDof(timeIdx); }
136
140 size_t numInteriorFaces(unsigned timeIdx) const
141 { return elemCtx_.numInteriorFaces(timeIdx); }
142
146 size_t numBoundaryFaces(unsigned timeIdx) const
147 { return elemCtx_.stencil(timeIdx).numBoundaryFaces(); }
148
152 const Stencil& stencil(unsigned timeIdx) const
153 { return elemCtx_.stencil(timeIdx); }
154
161 Vector normal(unsigned boundaryFaceIdx, unsigned timeIdx) const
162 {
163 auto tmp = stencil(timeIdx).boundaryFace[boundaryFaceIdx].normal;
164 tmp /= tmp.two_norm();
165 return tmp;
166 }
167
171 Scalar boundarySegmentArea(unsigned boundaryFaceIdx, unsigned timeIdx) const
172 { return elemCtx_.stencil(timeIdx).boundaryFace(boundaryFaceIdx).area(); }
173
180 const GlobalPosition& pos(unsigned boundaryFaceIdx, unsigned timeIdx) const
181 { return stencil(timeIdx).boundaryFace(boundaryFaceIdx).integrationPos(); }
182
189 const GlobalPosition& cvCenter(unsigned boundaryFaceIdx, unsigned timeIdx) const
190 {
191 unsigned scvIdx = stencil(timeIdx).boundaryFace(boundaryFaceIdx).interiorIndex();
192 return stencil(timeIdx).subControlVolume(scvIdx).globalPos();
193 }
194
199 unsigned focusDofIndex() const
200 { return elemCtx_.focusDofIndex(); }
201
209 unsigned interiorScvIndex(unsigned boundaryFaceIdx, unsigned timeIdx) const
210 { return stencil(timeIdx).boundaryFace(boundaryFaceIdx).interiorIndex(); }
211
219 unsigned globalSpaceIndex(unsigned boundaryFaceIdx, unsigned timeIdx) const
220 { return elemCtx_.globalSpaceIndex(interiorScvIndex(boundaryFaceIdx, timeIdx), timeIdx); }
221
229 const IntensiveQuantities& intensiveQuantities(unsigned boundaryFaceIdx, unsigned timeIdx) const
230 {
231 unsigned interiorScvIdx = this->interiorScvIndex(boundaryFaceIdx, timeIdx);
232 return elemCtx_.intensiveQuantities(interiorScvIdx, timeIdx);
233 }
234
241 const ExtensiveQuantities& extensiveQuantities(unsigned boundaryFaceIdx, unsigned timeIdx) const
242 { return elemCtx_.boundaryExtensiveQuantities(boundaryFaceIdx, timeIdx); }
243
253 const Intersection intersection(unsigned) const
254 { return *intersectionIt_; }
255
264 IntersectionIterator& intersectionIt()
265 { return intersectionIt_; }
266
267protected:
268 const ElementContext& elemCtx_;
269 IntersectionIterator intersectionIt_;
270};
271
272} // namespace Opm
273
274#endif
Represents all quantities which available on boundary segments.
Definition: fvbaseboundarycontext.hh:44
const GlobalPosition & pos(unsigned boundaryFaceIdx, unsigned timeIdx) const
Return the position of a local entity in global coordinates.
Definition: fvbaseboundarycontext.hh:180
const GlobalPosition & cvCenter(unsigned boundaryFaceIdx, unsigned timeIdx) const
Return the position of a control volume's center in global coordinates.
Definition: fvbaseboundarycontext.hh:189
Vector normal(unsigned boundaryFaceIdx, unsigned timeIdx) const
Returns the outer unit normal of the boundary segment.
Definition: fvbaseboundarycontext.hh:161
size_t numInteriorFaces(unsigned timeIdx) const
Definition: fvbaseboundarycontext.hh:140
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:219
const Stencil & stencil(unsigned timeIdx) const
Definition: fvbaseboundarycontext.hh:152
FvBaseBoundaryContext(const ElementContext &elemCtx)
The constructor.
Definition: fvbaseboundarycontext.hh:70
Scalar boundarySegmentArea(unsigned boundaryFaceIdx, unsigned timeIdx) const
Returns the area [m^2] of a given boudary segment.
Definition: fvbaseboundarycontext.hh:171
unsigned focusDofIndex() const
Return the local sub-control volume index upon which the linearization is currently focused.
Definition: fvbaseboundarycontext.hh:199
const Intersection intersection(unsigned) const
Return the intersection for the neumann segment.
Definition: fvbaseboundarycontext.hh:253
const GridView & gridView() const
Definition: fvbaseboundarycontext.hh:104
size_t numBoundaryFaces(unsigned timeIdx) const
Return the number of boundary segments of the current element.
Definition: fvbaseboundarycontext.hh:146
const ElementContext & elementContext() const
Returns a reference to the element context object.
Definition: fvbaseboundarycontext.hh:116
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:229
const Problem & problem() const
Definition: fvbaseboundarycontext.hh:92
const Element & element() const
Definition: fvbaseboundarycontext.hh:110
const GradientCalculator & gradientCalculator() const
Returns a reference to the current gradient calculator.
Definition: fvbaseboundarycontext.hh:122
unsigned interiorScvIndex(unsigned boundaryFaceIdx, unsigned timeIdx) const
Return the local sub-control volume index of the interior of a boundary segment.
Definition: fvbaseboundarycontext.hh:209
size_t numDof(unsigned timeIdx) const
Definition: fvbaseboundarycontext.hh:128
IntersectionIterator intersectionIt_
Definition: fvbaseboundarycontext.hh:269
size_t numPrimaryDof(unsigned timeIdx) const
Definition: fvbaseboundarycontext.hh:134
const ElementContext & elemCtx_
Definition: fvbaseboundarycontext.hh:268
void increment()
Definition: fvbaseboundarycontext.hh:75
const Model & model() const
Definition: fvbaseboundarycontext.hh:98
const ExtensiveQuantities & extensiveQuantities(unsigned boundaryFaceIdx, unsigned timeIdx) const
Return the extensive quantities for a given boundary face.
Definition: fvbaseboundarycontext.hh:241
IntersectionIterator & intersectionIt()
Return the intersection for the neumann segment.
Definition: fvbaseboundarycontext.hh:264
Declare the properties used by the infrastructure code of the finite volume discretizations.
Definition: blackoilboundaryratevector.hh:37
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:242