subcontrolvolume.hh 4.28 KB
Newer Older
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
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
 *   See the file COPYING for full copying permissions.                      *
 *                                                                           *
 *   This program 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.                                     *
 *                                                                           *
 *   This program 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 this program.  If not, see <http://www.gnu.org/licenses/>.   *
 *****************************************************************************/
/*!
 * \file
 * \brief Base class for a sub control volume
 */
#ifndef DUMUX_DISCRETIZATION_BOX_SUBCONTROLVOLUME_HH
#define DUMUX_DISCRETIZATION_BOX_SUBCONTROLVOLUME_HH

#include <dune/common/fvector.hh>
#include <dumux/discretization/subcontrolvolumebase.hh>
Timo Koch's avatar
Timo Koch committed
28
#include <dumux/discretization/box/boxgeometryhelper.hh>
29
#include <dumux/common/math.hh>
30
31
32
33

namespace Dumux
{
template<class G, typename I>
34
class BoxSubControlVolume : public SubControlVolumeBase<BoxSubControlVolume<G, I>, G, I>
35
{
36
    using ParentType = SubControlVolumeBase<BoxSubControlVolume<G, I>, G, I>;
37
38
39
40
    using Geometry = G;
    using IndexType = I;

    using Scalar = typename Geometry::ctype;
41
    enum { dim = Geometry::mydimension };
42
43
44
    enum { dimworld = Geometry::coorddimension };
    using GlobalPosition = Dune::FieldVector<Scalar, dimworld>;

Timo Koch's avatar
Timo Koch committed
45

46
public:
Timo Koch's avatar
Timo Koch committed
47
48
49
    //! The default constructor
    BoxSubControlVolume() = default;

50
    // the contructor in the box case
Timo Koch's avatar
Timo Koch committed
51
52
    template<class GeometryHelper>
    BoxSubControlVolume(const GeometryHelper& geometryHelper,
53
54
55
                        IndexType scvIdx,
                        IndexType elementIndex,
                        IndexType dofIndex)
Timo Koch's avatar
Timo Koch committed
56
    : corners_(geometryHelper.getScvCorners(scvIdx)),
57
      center_(0.0),
Timo Koch's avatar
Timo Koch committed
58
59
60
61
62
      volume_(geometryHelper.scvVolume(corners_)),
      elementIndex_(elementIndex),
      scvIdx_(scvIdx),
      dofIndex_(dofIndex)
    {
63
64
65
66
        // compute center point
        for (const auto& corner : corners_)
            center_ += corner;
        center_ /= corners_.size();
Timo Koch's avatar
Timo Koch committed
67
    }
68
69
70
71
72
73

    //! The center of the sub control volume
    GlobalPosition center() const
    {
        return center_;
    }
74

75
76
77
78
79
80
81
82
83
84
85
86
87
88
    //! The volume of the sub control volume
    Scalar volume() const
    {
        return volume_;
    }

    //! The geometry of the sub control volume
    // e.g. for integration
    Geometry geometry() const
    {
        return Geometry(Dune::GeometryType(Dune::GeometryType::cube, dim), corners_);
    }

    //! The global index of this scv
89
    IndexType indexInElement() const
90
91
92
93
    {
        return scvIdx_;
    }

94
    //! The index of the dof this scv is embedded in
95
96
97
98
99
    IndexType dofIndex() const
    {
        return dofIndex_;
    }

100
    // The position of the dof this scv is embedded in
101
102
    GlobalPosition dofPosition() const
    {
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
        // The corner list is defined such that the first entry is the vertex itself
        return corner(0);
    }

    //! The global index of the element this scv is embedded in
    IndexType elementIndex() const
    {
        return elementIndex_;
    }

    //! Return the corner for the given local index
    GlobalPosition corner(unsigned int localIdx) const
    {
        assert(localIdx < corners_.size() && "provided index exceeds the number of corners");
        return corners_[localIdx];
118
119
120
    }

private:
121
122
123
124
    std::vector<GlobalPosition> corners_;
    GlobalPosition center_;
    Scalar volume_;
    IndexType elementIndex_;
125
126
127
128
129
130
131
    IndexType scvIdx_;
    IndexType dofIndex_;
};

} // end namespace

#endif