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  Copyright (C) 2011-2013 by Andreas Lauser
5 
6  This file is part of the Open Porous Media project (OPM).
7 
8  OPM is free software: you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation, either version 2 of the License, or
11  (at your option) any later version.
12 
13  OPM is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with OPM. If not, see <http://www.gnu.org/licenses/>.
20 */
26 #ifndef EWOMS_FV_BASE_BOUNDARY_CONTEXT_HH
27 #define EWOMS_FV_BASE_BOUNDARY_CONTEXT_HH
28 
29 #include "fvbaseproperties.hh"
30 
31 #include <dune/common/fvector.hh>
32 
33 namespace Ewoms {
34 
40 template<class TypeTag>
42 {
43  typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
44  typedef typename GET_PROP_TYPE(TypeTag, Problem) Problem;
45  typedef typename GET_PROP_TYPE(TypeTag, Model) Model;
46  typedef typename GET_PROP_TYPE(TypeTag, Stencil) Stencil;
47  typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
48  typedef typename GET_PROP_TYPE(TypeTag, IntensiveQuantities) IntensiveQuantities;
49  typedef typename GET_PROP_TYPE(TypeTag, ExtensiveQuantities) ExtensiveQuantities;
50  typedef typename GET_PROP_TYPE(TypeTag, GradientCalculator) GradientCalculator;
51 
52  typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
53  typedef typename GridView::template Codim<0>::Entity Element;
54  typedef typename GridView::IntersectionIterator IntersectionIterator;
55  typedef typename GridView::Intersection Intersection;
56 
57  enum { dim = GridView::dimension };
58  enum { dimWorld = GridView::dimensionworld };
59 
60  typedef typename GridView::ctype CoordScalar;
61  typedef Dune::FieldVector<CoordScalar, dim> GlobalPosition;
62  typedef Dune::FieldVector<Scalar, dimWorld> Vector;
63 
64 public:
68  explicit FvBaseBoundaryContext(const ElementContext &elemCtx)
69  : elemCtx_(elemCtx)
70  , intersectionIt_(gridView().ibegin(element()))
71  { }
72 
76  const Problem &problem() const
77  { return elemCtx_.problem(); }
78 
82  const Model &model() const
83  { return elemCtx_.model(); }
84 
88  const GridView &gridView() const
89  { return elemCtx_.gridView(); }
90 
94  const Element &element() const
95  { return elemCtx_.element(); }
96 
100  const ElementContext &elementContext() const
101  { return elemCtx_; }
102 
106  const GradientCalculator &gradientCalculator() const
107  { return elemCtx_.gradientCalculator(); }
108 
112  int numDof(int timeIdx) const
113  { return elemCtx_.numDof(timeIdx); }
114 
118  int numPrimaryDof(int timeIdx) const
119  { return elemCtx_.numPrimaryDof(timeIdx); }
120 
124  int numInteriorFaces(int timeIdx) const
125  { return elemCtx_.numInteriorFaces(timeIdx); }
126 
130  int numBoundaryFaces(int timeIdx) const
131  { return elemCtx_.stencil(timeIdx).numBoundaryFaces(); }
132 
136  const Stencil &stencil(int timeIdx) const
137  { return elemCtx_.stencil(timeIdx); }
138 
145  Vector normal(int boundaryFaceIdx, int timeIdx) const
146  {
147  auto tmp = stencil(timeIdx).boundaryFace[boundaryFaceIdx].normal;
148  tmp /= tmp.two_norm();
149  return tmp;
150  }
151 
155  Scalar boundarySegmentArea(int boundaryFaceIdx, int timeIdx) const
156  { return elemCtx_.stencil(timeIdx).boundaryFace(boundaryFaceIdx).area(); }
157 
164  const GlobalPosition &pos(int boundaryFaceIdx, int timeIdx) const
165  { return stencil(timeIdx).boundaryFace(boundaryFaceIdx).integrationPos(); }
166 
173  const GlobalPosition &cvCenter(int boundaryFaceIdx, int timeIdx) const
174  {
175  int scvIdx = stencil(timeIdx).boundaryFace(boundaryFaceIdx).interiorIndex();
176  return stencil(timeIdx).subControlVolume(scvIdx).globalPos();
177  }
178 
186  short interiorScvIndex(int boundaryFaceIdx, int timeIdx) const
187  {
188  return stencil(timeIdx).boundaryFace(boundaryFaceIdx).interiorIndex();
189  }
190 
198  int globalSpaceIndex(int boundaryFaceIdx, int timeIdx) const
199  { return elemCtx_.globalSpaceIndex(interiorScvIndex(boundaryFaceIdx, timeIdx), timeIdx); }
200 
208  const IntensiveQuantities &intensiveQuantities(int boundaryFaceIdx, int timeIdx) const
209  {
210  short interiorScvIdx = this->interiorScvIndex(boundaryFaceIdx, timeIdx);
211  return elemCtx_.intensiveQuantities(interiorScvIdx, timeIdx);
212  }
213 
220  const ExtensiveQuantities &extensiveQuantities(int boundaryFaceIdx, int timeIdx) const
221  { return elemCtx_.boundaryExtensiveQuantities(boundaryFaceIdx, timeIdx); }
222 
232  const Intersection &intersection(int boundaryFaceIdx) const
233  { return *intersectionIt_; }
234 
243  IntersectionIterator &intersectionIt()
244  { return intersectionIt_; }
245 
246 protected:
247  const ElementContext &elemCtx_;
248  IntersectionIterator intersectionIt_;
249 };
250 
251 } // namespace Ewoms
252 
253 #endif
int numInteriorFaces(int timeIdx) const
Definition: fvbaseboundarycontext.hh:124
const ExtensiveQuantities & extensiveQuantities(int boundaryFaceIdx, int timeIdx) const
Return the extensive quantities for a given boundary face.
Definition: fvbaseboundarycontext.hh:220
const Intersection & intersection(int boundaryFaceIdx) const
Return the intersection for the neumann segment.
Definition: fvbaseboundarycontext.hh:232
int numPrimaryDof(int timeIdx) const
Definition: fvbaseboundarycontext.hh:118
const Model & model() const
Definition: fvbaseboundarycontext.hh:82
const ElementContext & elementContext() const
Returns a reference to the element context object.
Definition: fvbaseboundarycontext.hh:100
Declare the properties used by the infrastructure code of the finite volume discretizations.
int numBoundaryFaces(int timeIdx) const
Return the number of boundary segments of the current element.
Definition: fvbaseboundarycontext.hh:130
IntersectionIterator & intersectionIt()
Return the intersection for the neumann segment.
Definition: fvbaseboundarycontext.hh:243
Scalar boundarySegmentArea(int boundaryFaceIdx, int timeIdx) const
Returns the area [m^2] of a given boudary segment.
Definition: fvbaseboundarycontext.hh:155
const ElementContext & elemCtx_
Definition: fvbaseboundarycontext.hh:247
Vector normal(int boundaryFaceIdx, int timeIdx) const
Returns the outer unit normal of the boundary segment.
Definition: fvbaseboundarycontext.hh:145
Definition: baseauxiliarymodule.hh:35
FvBaseBoundaryContext(const ElementContext &elemCtx)
The constructor.
Definition: fvbaseboundarycontext.hh:68
const Element & element() const
Definition: fvbaseboundarycontext.hh:94
const GlobalPosition & pos(int boundaryFaceIdx, int timeIdx) const
Return the position of a local entity in global coordinates.
Definition: fvbaseboundarycontext.hh:164
const GlobalPosition & cvCenter(int boundaryFaceIdx, int timeIdx) const
Return the position of a control volume's center in global coordinates.
Definition: fvbaseboundarycontext.hh:173
int numDof(int timeIdx) const
Definition: fvbaseboundarycontext.hh:112
int globalSpaceIndex(int boundaryFaceIdx, int timeIdx) const
Return the global space index of the sub-control volume at the interior of a boundary segment...
Definition: fvbaseboundarycontext.hh:198
Represents all quantities which available on boundary segments.
Definition: fvbaseboundarycontext.hh:41
const GradientCalculator & gradientCalculator() const
Returns a reference to the current gradient calculator.
Definition: fvbaseboundarycontext.hh:106
IntersectionIterator intersectionIt_
Definition: fvbaseboundarycontext.hh:248
const GridView & gridView() const
Definition: fvbaseboundarycontext.hh:88
const IntensiveQuantities & intensiveQuantities(int boundaryFaceIdx, int timeIdx) const
Return the intensive quantities for the finite volume in the interiour of a boundary segment...
Definition: fvbaseboundarycontext.hh:208
const Problem & problem() const
Definition: fvbaseboundarycontext.hh:76
short interiorScvIndex(int boundaryFaceIdx, int timeIdx) const
Return the local sub-control volume index of the interior of a boundary segment.
Definition: fvbaseboundarycontext.hh:186
const Stencil & stencil(int timeIdx) const
Definition: fvbaseboundarycontext.hh:136