3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
porousmediumflow/boxdfm/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 *****************************************************************************/
25#ifndef DUMUX_POROUSMEDIUMFLOW_BOXDFM_SUBCONTROLVOLUME_HH
26#define DUMUX_POROUSMEDIUMFLOW_BOXDFM_SUBCONTROLVOLUME_HH
27
28#include <dune/common/reservedvector.hh>
29#include <dune/geometry/multilineargeometry.hh>
30
33#include <dumux/common/math.hh>
34
35namespace Dumux {
36
47template<class GridView>
49{
50 using Grid = typename GridView::Grid;
51
52 static const int dim = Grid::dimension;
53 static const int dimWorld = Grid::dimensionworld;
54
55 template <class ct>
56 struct ScvMLGTraits : public Dune::MultiLinearGeometryTraits<ct>
57 {
58 // we use static vectors to store the corners as we know
59 // the number of corners in advance (2^(dim) corners (1<<(dim))
60 // However, on fracture scvs the number might be smaller (use ReservedVector)
61 template< int mydim, int cdim >
63 {
64 using Type = Dune::ReservedVector< Dune::FieldVector< ct, cdim >, (1<<(dim)) >;
65 };
66
67 // we know all scvfs will have the same geometry type
68 template< int mydim >
70 {
71 static const bool v = true;
72 static const unsigned int topologyId = Dune::Impl::CubeTopology< mydim >::type::id;
73 };
74 };
75
76 using GridIndexType = typename Grid::LeafGridView::IndexSet::IndexType;
77 using LocalIndexType = unsigned int;
78 using Scalar = typename Grid::ctype;
79 using Geometry = Dune::MultiLinearGeometry<Scalar, dim, dimWorld, ScvMLGTraits<Scalar>>;
81 using GlobalPosition = typename CornerStorage::value_type;
82};
83
91template<class GV,
94: public SubControlVolumeBase<BoxDfmSubControlVolume<GV, T>, T>
95{
98 using Geometry = typename T::Geometry;
99 using GridIndexType = typename T::GridIndexType;
100 using LocalIndexType = typename T::LocalIndexType;
101 using Scalar = typename T::Scalar;
102 using GlobalPosition = typename T::GlobalPosition;
103 using CornerStorage = typename T::CornerStorage;
104 enum { dim = Geometry::mydimension };
105
106 static_assert(dim == 2 || dim == 3, "Box-Dfm sub-control volume only implemented in 2d or 3d");
107
108public:
110 using Traits = T;
111
114
115 // the constructor for standard scvs
116 template<class GeometryHelper>
117 BoxDfmSubControlVolume(const GeometryHelper& geometryHelper,
118 LocalIndexType scvIdx,
119 GridIndexType elementIndex,
120 GridIndexType dofIndex)
121 : isFractureScv_(false)
122 , corners_(geometryHelper.getScvCorners(scvIdx))
123 , center_(0.0)
124 , volume_(geometryHelper.scvVolume(corners_))
125 , elementIndex_(elementIndex)
126 , vIdxLocal_(scvIdx)
127 , elemLocalScvIdx_(scvIdx)
128 , dofIndex_(dofIndex)
129 , facetIdx_(0)
130 {
131 // compute center point
132 for (const auto& corner : corners_)
133 center_ += corner;
134 center_ /= corners_.size();
135 }
136
146 template<class GeometryHelper, class Intersection>
147 BoxDfmSubControlVolume(const GeometryHelper& geometryHelper,
148 const Intersection& intersection,
149 const typename Intersection::Geometry& isGeometry,
150 LocalIndexType indexInIntersection,
151 LocalIndexType vIdxLocal,
152 LocalIndexType elemLocalScvIdx,
153 LocalIndexType elemLocalFacetIdx,
154 GridIndexType elementIndex,
155 GridIndexType dofIndex)
156 : isFractureScv_(true)
157 , corners_()
158 , center_(0.0)
159 , volume_(0.0)
160 , elementIndex_(elementIndex)
161 , vIdxLocal_(vIdxLocal)
162 , elemLocalScvIdx_(elemLocalScvIdx)
163 , dofIndex_(dofIndex)
164 , facetIdx_(elemLocalFacetIdx)
165 {
166 // copy corners
167 auto corners = geometryHelper.getBoundaryScvfCorners(intersection, isGeometry, indexInIntersection);
168 const auto numCorners = corners.size();
169 corners_.resize(numCorners);
170 for (unsigned int i = 0; i < numCorners; ++i)
171 corners_[i] = corners[i];
172
173 // compute volume and scv center
174 volume_ = geometryHelper.scvfArea(corners);
175 for (const auto& corner : corners_)
176 center_ += corner;
177 center_ /= corners_.size();
178 }
179
181 const GlobalPosition& center() const
182 { return center_; }
183
185 Scalar volume() const
186 { return volume_; }
187
189 LocalIndexType localDofIndex() const
190 { return vIdxLocal_; }
191
193 LocalIndexType indexInElement() const
194 { return elemLocalScvIdx_; }
195
197 LocalIndexType facetIndexInElement() const
198 { assert(isFractureScv_); return facetIdx_; }
199
201 GridIndexType dofIndex() const
202 { return dofIndex_; }
203
204 // The position of the dof this scv is embedded in (list is defined such that first entry is vertex itself)
205 const GlobalPosition& dofPosition() const
206 { return corners_[0]; }
207
209 GridIndexType elementIndex() const
210 { return elementIndex_; }
211
213 bool isOnFracture() const
214 { return isFractureScv_; }
215
217 // e.g. for integration
218 Geometry geometry() const
219 {
220 if (isFractureScv_)
221 DUNE_THROW(Dune::InvalidStateException, "The geometry object cannot be defined for fracture scvs "
222 "because the number of known corners is insufficient. "
223 "You can do this manually by extract the corners from this scv "
224 "and extrude them by the corresponding aperture. ");
225
226 return Geometry(Dune::GeometryTypes::cube(dim), corners_);
227 }
228
230 const GlobalPosition& corner(LocalIndexType localIdx) const
231 {
232 assert(localIdx < corners_.size() && "provided index exceeds the number of corners");
233 return corners_[localIdx];
234 }
235
236private:
237 bool isFractureScv_;
238 CornerStorage corners_;
239 GlobalPosition center_;
240 Scalar volume_;
241 GridIndexType elementIndex_;
242 LocalIndexType vIdxLocal_;
243 LocalIndexType elemLocalScvIdx_;
244 GridIndexType dofIndex_;
245
246 // for fracture scvs only!
247 LocalIndexType facetIdx_;
248};
249
250} // end namespace
251
252#endif
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
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
Default traits class to be used for the sub-control volumes for the box discrete fracture scheme.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:49
static const int dimWorld
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:53
static const int dim
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:52
typename GridView::Grid Grid
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:50
Dune::MultiLinearGeometry< Scalar, dim, dimWorld, ScvMLGTraits< Scalar > > Geometry
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:79
typename Grid::LeafGridView::IndexSet::IndexType GridIndexType
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:76
typename ScvMLGTraits< Scalar >::template CornerStorage< dim, dimWorld >::Type CornerStorage
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:80
typename CornerStorage::value_type GlobalPosition
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:81
unsigned int LocalIndexType
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:77
typename Grid::ctype Scalar
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:78
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:57
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:63
Dune::ReservedVector< Dune::FieldVector< ct, cdim >,(1<<(dim)) > Type
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:64
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:70
static const unsigned int topologyId
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:72
static const bool v
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:71
the sub control volume for the box discrete fracture scheme
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:95
BoxDfmSubControlVolume()=default
The default constructor.
const GlobalPosition & center() const
The center of the sub control volume.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:181
LocalIndexType facetIndexInElement() const
The element-local facet index for which a fracture scv was created.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:197
LocalIndexType indexInElement() const
The element-local index of this scv.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:193
Geometry geometry() const
The geometry of the sub control volume.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:218
bool isOnFracture() const
Return true if this scv is part of the fracture domain.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:213
const GlobalPosition & corner(LocalIndexType localIdx) const
Return the corner for the given local index.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:230
Scalar volume() const
The volume of the sub control volume.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:185
const GlobalPosition & dofPosition() const
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:205
T Traits
State the traits public and thus export all types.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:110
GridIndexType elementIndex() const
The global index of the element this scv is embedded in.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:209
BoxDfmSubControlVolume(const GeometryHelper &geometryHelper, const Intersection &intersection, const typename Intersection::Geometry &isGeometry, LocalIndexType indexInIntersection, LocalIndexType vIdxLocal, LocalIndexType elemLocalScvIdx, LocalIndexType elemLocalFacetIdx, GridIndexType elementIndex, GridIndexType dofIndex)
Constructor for fracture scvs.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:147
LocalIndexType localDofIndex() const
The element-local vertex index this scv is connected to.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:189
BoxDfmSubControlVolume(const GeometryHelper &geometryHelper, LocalIndexType scvIdx, GridIndexType elementIndex, GridIndexType dofIndex)
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:117
GridIndexType dofIndex() const
The index of the dof this scv is embedded in.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:201