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 <dune/common/fvector.hh>
32
34
35#include <cstddef>
36
37namespace Opm {
38
44template<class TypeTag>
46{
55
57 using Element = typename GridView::template Codim<0>::Entity;
58 using IntersectionIterator = typename GridView::IntersectionIterator;
59 using Intersection = typename GridView::Intersection;
60
61 enum { dim = GridView::dimension };
62 enum { dimWorld = GridView::dimensionworld };
63
64 using CoordScalar = typename GridView::ctype;
65 using GlobalPosition = Dune::FieldVector<CoordScalar, dimWorld>;
66 using Vector = Dune::FieldVector<Scalar, dimWorld>;
67
68public:
72 explicit FvBaseBoundaryContext(const ElementContext& elemCtx)
73 : elemCtx_(elemCtx)
74 , intersectionIt_(gridView().ibegin(element()))
75 {}
76
77 void increment()
78 {
79 const auto& iend = gridView().iend(element());
80
81 if (intersectionIt_ == iend) {
82 return;
83 }
84
86 // iterate to the next boundary intersection
87 while (intersectionIt_ != iend && intersectionIt_->neighbor()) {
89 }
90 }
91
95 const Problem& problem() const
96 { return elemCtx_.problem(); }
97
101 const Model& model() const
102 { return elemCtx_.model(); }
103
107 const GridView& gridView() const
108 { return elemCtx_.gridView(); }
109
113 const Element& element() const
114 { return elemCtx_.element(); }
115
119 const ElementContext& elementContext() const
120 { return elemCtx_; }
121
125 const GradientCalculator& gradientCalculator() const
126 { return elemCtx_.gradientCalculator(); }
127
131 std::size_t numDof(unsigned timeIdx) const
132 { return elemCtx_.numDof(timeIdx); }
133
137 std::size_t numPrimaryDof(unsigned timeIdx) const
138 { return elemCtx_.numPrimaryDof(timeIdx); }
139
143 std::size_t numInteriorFaces(unsigned timeIdx) const
144 { return elemCtx_.numInteriorFaces(timeIdx); }
145
149 std::size_t numBoundaryFaces(unsigned timeIdx) const
150 { return elemCtx_.stencil(timeIdx).numBoundaryFaces(); }
151
155 const Stencil& stencil(unsigned timeIdx) const
156 { return elemCtx_.stencil(timeIdx); }
157
164 Vector normal(unsigned boundaryFaceIdx, unsigned timeIdx) const
165 {
166 auto tmp = stencil(timeIdx).boundaryFace[boundaryFaceIdx].normal;
167 tmp /= tmp.two_norm();
168 return tmp;
169 }
170
174 Scalar boundarySegmentArea(unsigned boundaryFaceIdx, unsigned timeIdx) const
175 { return elemCtx_.stencil(timeIdx).boundaryFace(boundaryFaceIdx).area(); }
176
183 const GlobalPosition& pos(unsigned boundaryFaceIdx, unsigned timeIdx) const
184 { return stencil(timeIdx).boundaryFace(boundaryFaceIdx).integrationPos(); }
185
192 const GlobalPosition& cvCenter(unsigned boundaryFaceIdx, unsigned timeIdx) const
193 {
194 const unsigned scvIdx = stencil(timeIdx).boundaryFace(boundaryFaceIdx).interiorIndex();
195 return stencil(timeIdx).subControlVolume(scvIdx).globalPos();
196 }
197
202 unsigned focusDofIndex() const
203 { return elemCtx_.focusDofIndex(); }
204
212 unsigned interiorScvIndex(unsigned boundaryFaceIdx, unsigned timeIdx) const
213 { return stencil(timeIdx).boundaryFace(boundaryFaceIdx).interiorIndex(); }
214
222 unsigned globalSpaceIndex(unsigned boundaryFaceIdx, unsigned timeIdx) const
223 { return elemCtx_.globalSpaceIndex(interiorScvIndex(boundaryFaceIdx, timeIdx), timeIdx); }
224
232 const IntensiveQuantities& intensiveQuantities(unsigned boundaryFaceIdx, unsigned timeIdx) const
233 {
234 const unsigned interiorScvIdx = this->interiorScvIndex(boundaryFaceIdx, timeIdx);
235 return elemCtx_.intensiveQuantities(interiorScvIdx, timeIdx);
236 }
237
244 const ExtensiveQuantities& extensiveQuantities(unsigned boundaryFaceIdx, unsigned timeIdx) const
245 { return elemCtx_.boundaryExtensiveQuantities(boundaryFaceIdx, timeIdx); }
246
256 const Intersection intersection(unsigned) const
257 { return *intersectionIt_; }
258
267 IntersectionIterator& intersectionIt()
268 { return intersectionIt_; }
269
270protected:
271 const ElementContext& elemCtx_;
272 IntersectionIterator intersectionIt_;
273};
274
275} // namespace Opm
276
277#endif
Represents all quantities which available on boundary segments.
Definition: fvbaseboundarycontext.hh:46
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
IntersectionIterator intersectionIt_
Definition: fvbaseboundarycontext.hh:272
std::size_t numPrimaryDof(unsigned timeIdx) const
Definition: fvbaseboundarycontext.hh:137
const ElementContext & elemCtx_
Definition: fvbaseboundarycontext.hh:271
void increment()
Definition: fvbaseboundarycontext.hh:77
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.
Definition: blackoilboundaryratevector.hh:39
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