File: fvbaseboundarycontext.hh

package info (click to toggle)
opm-simulators 2024.10%2Bds-6
  • links: PTS, VCS
  • area: main
  • in suites: trixie
  • size: 19,416 kB
  • sloc: cpp: 165,337; sh: 1,285; lisp: 1,108; python: 355; makefile: 24; awk: 10
file content (274 lines) | stat: -rw-r--r-- 9,427 bytes parent folder | download | duplicates (3)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*
  This file is part of the Open Porous Media project (OPM).

  OPM is free software: you can redistribute it and/or modify
  it under the terms of the GNU General Public License as published by
  the Free Software Foundation, either version 2 of the License, or
  (at your option) any later version.

  OPM is distributed in the hope that it will be useful,
  but WITHOUT ANY WARRANTY; without even the implied warranty of
  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
  GNU General Public License for more details.

  You should have received a copy of the GNU General Public License
  along with OPM.  If not, see <http://www.gnu.org/licenses/>.

  Consult the COPYING file in the top-level source directory of this
  module for the precise wording of the license and the list of
  copyright holders.
*/
/*!
 * \file
 *
 * \copydoc Opm::FvBaseBoundaryContext
 */
#ifndef EWOMS_FV_BASE_BOUNDARY_CONTEXT_HH
#define EWOMS_FV_BASE_BOUNDARY_CONTEXT_HH

#include "fvbaseproperties.hh"

#include <dune/common/fvector.hh>

namespace Opm {

/*!
 * \ingroup FiniteVolumeDiscretizations
 *
 * \brief Represents all quantities which available on boundary segments
 */
template<class TypeTag>
class FvBaseBoundaryContext
{
    using Scalar = GetPropType<TypeTag, Properties::Scalar>;
    using Problem = GetPropType<TypeTag, Properties::Problem>;
    using Model = GetPropType<TypeTag, Properties::Model>;
    using Stencil = GetPropType<TypeTag, Properties::Stencil>;
    using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
    using IntensiveQuantities = GetPropType<TypeTag, Properties::IntensiveQuantities>;
    using ExtensiveQuantities = GetPropType<TypeTag, Properties::ExtensiveQuantities>;
    using GradientCalculator = GetPropType<TypeTag, Properties::GradientCalculator>;

    using GridView = GetPropType<TypeTag, Properties::GridView>;
    using Element = typename GridView::template Codim<0>::Entity;
    using IntersectionIterator = typename GridView::IntersectionIterator;
    using Intersection = typename GridView::Intersection;

    enum { dim = GridView::dimension };
    enum { dimWorld = GridView::dimensionworld };

    using CoordScalar = typename GridView::ctype;
    using GlobalPosition = Dune::FieldVector<CoordScalar, dimWorld>;
    using Vector = Dune::FieldVector<Scalar, dimWorld>;

public:
    /*!
     * \brief The constructor.
     */
    explicit FvBaseBoundaryContext(const ElementContext& elemCtx)
        : elemCtx_(elemCtx)
        , intersectionIt_(gridView().ibegin(element()))
    { }

    void increment()
    {
        const auto& iend = gridView().iend(element());

        if(intersectionIt_ == iend)
          return;

        ++intersectionIt_;
        // iterate to the next boundary intersection
        while (intersectionIt_ != iend && intersectionIt_->neighbor()) {
            ++intersectionIt_;
        }
    }

    /*!
     * \copydoc Opm::ElementContext::problem()
     */
    const Problem& problem() const
    { return elemCtx_.problem(); }

    /*!
     * \copydoc Opm::ElementContext::model()
     */
    const Model& model() const
    { return elemCtx_.model(); }

    /*!
     * \copydoc Opm::ElementContext::gridView()
     */
    const GridView& gridView() const
    { return elemCtx_.gridView(); }

    /*!
     * \copydoc Opm::ElementContext::element()
     */
    const Element& element() const
    { return elemCtx_.element(); }

    /*!
     * \brief Returns a reference to the element context object.
     */
    const ElementContext& elementContext() const
    { return elemCtx_; }

    /*!
     * \brief Returns a reference to the current gradient calculator.
     */
    const GradientCalculator& gradientCalculator() const
    { return elemCtx_.gradientCalculator(); }

    /*!
     * \copydoc Opm::ElementContext::numDof()
     */
    size_t numDof(unsigned timeIdx) const
    { return elemCtx_.numDof(timeIdx); }

    /*!
     * \copydoc Opm::ElementContext::numPrimaryDof()
     */
    size_t numPrimaryDof(unsigned timeIdx) const
    { return elemCtx_.numPrimaryDof(timeIdx); }

    /*!
     * \copydoc Opm::ElementContext::numInteriorFaces()
     */
    size_t numInteriorFaces(unsigned timeIdx) const
    { return elemCtx_.numInteriorFaces(timeIdx); }

    /*!
     * \brief Return the number of boundary segments of the current element
     */
    size_t numBoundaryFaces(unsigned timeIdx) const
    { return elemCtx_.stencil(timeIdx).numBoundaryFaces(); }

    /*!
     * \copydoc Opm::ElementContext::stencil()
     */
    const Stencil& stencil(unsigned timeIdx) const
    { return elemCtx_.stencil(timeIdx); }

    /*!
     * \brief Returns the outer unit normal of the boundary segment
     *
     * \param boundaryFaceIdx The local index of the boundary segment
     * \param timeIdx The index of the solution used by the time discretization
     */
    Vector normal(unsigned boundaryFaceIdx, unsigned timeIdx) const
    {
        auto tmp = stencil(timeIdx).boundaryFace[boundaryFaceIdx].normal;
        tmp /= tmp.two_norm();
        return tmp;
    }

    /*!
     * \brief Returns the area [m^2] of a given boudary segment.
     */
    Scalar boundarySegmentArea(unsigned boundaryFaceIdx, unsigned timeIdx) const
    { return elemCtx_.stencil(timeIdx).boundaryFace(boundaryFaceIdx).area(); }

    /*!
     * \brief Return the position of a local entity in global coordinates.
     *
     * \param boundaryFaceIdx The local index of the boundary segment
     * \param timeIdx The index of the solution used by the time discretization
     */
    const GlobalPosition& pos(unsigned boundaryFaceIdx, unsigned timeIdx) const
    { return stencil(timeIdx).boundaryFace(boundaryFaceIdx).integrationPos(); }

    /*!
     * \brief Return the position of a control volume's center in global coordinates.
     *
     * \param boundaryFaceIdx The local index of the boundary segment
     * \param timeIdx The index of the solution used by the time discretization
     */
    const GlobalPosition& cvCenter(unsigned boundaryFaceIdx, unsigned timeIdx) const
    {
        unsigned scvIdx = stencil(timeIdx).boundaryFace(boundaryFaceIdx).interiorIndex();
        return stencil(timeIdx).subControlVolume(scvIdx).globalPos();
    }

    /*!
     * \brief Return the local sub-control volume index upon which the linearization is
     *        currently focused.
     */
    unsigned focusDofIndex() const
    { return elemCtx_.focusDofIndex(); }

    /*!
     * \brief Return the local sub-control volume index of the
     *        interior of a boundary segment
     *
     * \param boundaryFaceIdx The local index of the boundary segment
     * \param timeIdx The index of the solution used by the time discretization
     */
    unsigned interiorScvIndex(unsigned boundaryFaceIdx, unsigned timeIdx) const
    { return stencil(timeIdx).boundaryFace(boundaryFaceIdx).interiorIndex(); }

    /*!
     * \brief Return the global space index of the sub-control volume
     *        at the interior of a boundary segment
     *
     * \param boundaryFaceIdx The local index of the boundary segment
     * \param timeIdx The index of the solution used by the time discretization
     */
    unsigned globalSpaceIndex(unsigned boundaryFaceIdx, unsigned timeIdx) const
    { return elemCtx_.globalSpaceIndex(interiorScvIndex(boundaryFaceIdx, timeIdx), timeIdx); }

    /*!
     * \brief Return the intensive quantities for the finite volume in the
     *        interiour of a boundary segment
     *
     * \param boundaryFaceIdx The local index of the boundary segment
     * \param timeIdx The index of the solution used by the time discretization
     */
    const IntensiveQuantities& intensiveQuantities(unsigned boundaryFaceIdx, unsigned timeIdx) const
    {
        unsigned interiorScvIdx = this->interiorScvIndex(boundaryFaceIdx, timeIdx);
        return elemCtx_.intensiveQuantities(interiorScvIdx, timeIdx);
    }

    /*!
     * \brief Return the extensive quantities for a given boundary face.
     *
     * \param boundaryFaceIdx The local index of the boundary segment
     * \param timeIdx The index of the solution used by the time discretization
     */
    const ExtensiveQuantities& extensiveQuantities(unsigned boundaryFaceIdx, unsigned timeIdx) const
    { return elemCtx_.boundaryExtensiveQuantities(boundaryFaceIdx, timeIdx); }

    /*!
     * \brief Return the intersection for the neumann segment
     *
     * TODO/HACK: The intersection should take a local index as an
     * argument. since that's not supported efficiently by the DUNE
     * grid interface, we just ignore the index argument here!
     *
     * \param boundaryFaceIdx The local index of the boundary segment
     */
    const Intersection intersection(unsigned) const
    { return *intersectionIt_; }

    /*!
     * \brief Return the intersection for the neumann segment
     *
     * TODO/HACK: the intersection iterator can basically be
     * considered as an index which is manipulated externally, but
     * context classes should not store any indices. it is done this
     * way for performance reasons
     */
    IntersectionIterator& intersectionIt()
    { return intersectionIt_; }

protected:
    const ElementContext& elemCtx_;
    IntersectionIterator intersectionIt_;
};

} // namespace Opm

#endif