3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
discretization/box/subcontrolvolume.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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program 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 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program 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 this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_DISCRETIZATION_BOX_SUBCONTROLVOLUME_HH
25#define DUMUX_DISCRETIZATION_BOX_SUBCONTROLVOLUME_HH
26
27#include <dune/geometry/type.hh>
28#include <dune/geometry/multilineargeometry.hh>
29
30#include <dumux/common/math.hh>
34
35namespace Dumux {
36
43template<class GridView>
45{
46 using Grid = typename GridView::Grid;
47
48 static const int dim = Grid::dimension;
49 static const int dimWorld = Grid::dimensionworld;
50
51 template <class ct>
52 struct ScvMLGTraits : public Dune::MultiLinearGeometryTraits<ct>
53 {
54 // we use static vectors to store the corners as we know
55 // the number of corners in advance (2^(dim) corners (1<<(dim))
56 template< int mydim, int cdim >
58 {
59 using Type = std::array< Dune::FieldVector< ct, cdim >, (1<<(dim)) >;
60 };
61
62 // we know all scvfs will have the same geometry type
63 template< int mydim >
65 {
66 static const bool v = true;
67 static const unsigned int topologyId = Dune::GeometryTypes::cube(mydim).id();
68 };
69 };
70
73 using Scalar = typename Grid::ctype;
74 using Geometry = Dune::MultiLinearGeometry<Scalar, dim, dimWorld, ScvMLGTraits<Scalar>>;
76 using GlobalPosition = typename CornerStorage::value_type;
77};
78
85template<class GV,
88: public SubControlVolumeBase<BoxSubControlVolume<GV, T>, T>
89{
92 using Geometry = typename T::Geometry;
93 using GridIndexType = typename T::GridIndexType;
94 using LocalIndexType = typename T::LocalIndexType;
95 using Scalar = typename T::Scalar;
96 using CornerStorage = typename T::CornerStorage;
97 enum { dim = Geometry::mydimension };
98
99public:
101 using GlobalPosition = typename T::GlobalPosition;
103 using Traits = T;
104
107
108 // the contructor in the box case
109 template<class GeometryHelper>
110 BoxSubControlVolume(const GeometryHelper& geometryHelper,
111 LocalIndexType scvIdx,
112 GridIndexType elementIndex,
113 GridIndexType dofIndex)
114 : corners_(geometryHelper.getScvCorners(scvIdx)),
115 center_(0.0),
116 volume_(geometryHelper.scvVolume(corners_)),
117 elementIndex_(elementIndex),
118 localDofIdx_(scvIdx),
119 dofIndex_(dofIndex)
120 {
121 // compute center point
122 for (const auto& corner : corners_)
123 center_ += corner;
124 center_ /= corners_.size();
125 }
126
128 const GlobalPosition& center() const
129 {
130 return center_;
131 }
132
134 Scalar volume() const
135 {
136 return volume_;
137 }
138
140 // e.g. for integration
141 Geometry geometry() const
142 {
143 return Geometry(Dune::GeometryTypes::cube(dim), corners_);
144 }
145
147 LocalIndexType localDofIndex() const
148 {
149 return localDofIdx_;
150 }
151
154 LocalIndexType indexInElement() const
155 {
156 return localDofIdx_;
157 }
158
160 GridIndexType dofIndex() const
161 {
162 return dofIndex_;
163 }
164
165 // The position of the dof this scv is embedded in
167 {
168 // The corner list is defined such that the first entry is the vertex itself
169 return corners_[0];
170 }
171
173 GridIndexType elementIndex() const
174 {
175 return elementIndex_;
176 }
177
179 const GlobalPosition& corner(LocalIndexType localIdx) const
180 {
181 assert(localIdx < corners_.size() && "provided index exceeds the number of corners");
182 return corners_[localIdx];
183 }
184
185private:
186 CornerStorage corners_;
187 GlobalPosition center_;
188 Scalar volume_;
189 GridIndexType elementIndex_;
190 LocalIndexType localDofIdx_;
191 GridIndexType dofIndex_;
192};
193
194} // end namespace Dumux
195
196#endif
Defines the index types used for grid and local indices.
Define some often used mathematical functions.
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
Base class for a sub control volume.
Definition: adapt.hh:29
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:39
unsigned int LocalIndex
Definition: indextraits.hh:40
Default traits class to be used for the sub-control volumes for the box scheme.
Definition: discretization/box/subcontrolvolume.hh:45
typename GridView::Grid Grid
Definition: discretization/box/subcontrolvolume.hh:46
Dune::MultiLinearGeometry< Scalar, dim, dimWorld, ScvMLGTraits< Scalar > > Geometry
Definition: discretization/box/subcontrolvolume.hh:74
typename IndexTraits< GridView >::GridIndex GridIndexType
Definition: discretization/box/subcontrolvolume.hh:71
typename ScvMLGTraits< Scalar >::template CornerStorage< dim, dimWorld >::Type CornerStorage
Definition: discretization/box/subcontrolvolume.hh:75
typename IndexTraits< GridView >::LocalIndex LocalIndexType
Definition: discretization/box/subcontrolvolume.hh:72
static const int dim
Definition: discretization/box/subcontrolvolume.hh:48
static const int dimWorld
Definition: discretization/box/subcontrolvolume.hh:49
typename Grid::ctype Scalar
Definition: discretization/box/subcontrolvolume.hh:73
typename CornerStorage::value_type GlobalPosition
Definition: discretization/box/subcontrolvolume.hh:76
Definition: discretization/box/subcontrolvolume.hh:53
Definition: discretization/box/subcontrolvolume.hh:58
std::array< Dune::FieldVector< ct, cdim >,(1<<(dim)) > Type
Definition: discretization/box/subcontrolvolume.hh:59
Definition: discretization/box/subcontrolvolume.hh:65
static const bool v
Definition: discretization/box/subcontrolvolume.hh:66
static const unsigned int topologyId
Definition: discretization/box/subcontrolvolume.hh:67
the sub control volume for the box scheme
Definition: discretization/box/subcontrolvolume.hh:89
Geometry geometry() const
The geometry of the sub control volume.
Definition: discretization/box/subcontrolvolume.hh:141
const GlobalPosition & corner(LocalIndexType localIdx) const
Return the corner for the given local index.
Definition: discretization/box/subcontrolvolume.hh:179
T Traits
state the traits public and thus export all types
Definition: discretization/box/subcontrolvolume.hh:103
LocalIndexType localDofIndex() const
The element-local index of the dof this scv is embedded in.
Definition: discretization/box/subcontrolvolume.hh:147
const GlobalPosition & center() const
The center of the sub control volume.
Definition: discretization/box/subcontrolvolume.hh:128
BoxSubControlVolume()=default
The default constructor.
typename T::GlobalPosition GlobalPosition
export the type used for global coordinates
Definition: discretization/box/subcontrolvolume.hh:101
LocalIndexType indexInElement() const
Definition: discretization/box/subcontrolvolume.hh:154
GridIndexType elementIndex() const
The global index of the element this scv is embedded in.
Definition: discretization/box/subcontrolvolume.hh:173
GridIndexType dofIndex() const
The index of the dof this scv is embedded in.
Definition: discretization/box/subcontrolvolume.hh:160
Scalar volume() const
The volume of the sub control volume.
Definition: discretization/box/subcontrolvolume.hh:134
const GlobalPosition & dofPosition() const
Definition: discretization/box/subcontrolvolume.hh:166
BoxSubControlVolume(const GeometryHelper &geometryHelper, LocalIndexType scvIdx, GridIndexType elementIndex, GridIndexType dofIndex)
Definition: discretization/box/subcontrolvolume.hh:110
Base class for a sub control volume, i.e a part of the control volume we are making the balance for....
Definition: subcontrolvolumebase.hh:38