3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
porousmediumflow/boxdfm/subcontrolvolumeface.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_SUBCONTROLVOLUMEFACE_HH
26#define DUMUX_POROUSMEDIUMFLOW_BOXDFM_SUBCONTROLVOLUMEFACE_HH
27
28#include <utility>
29
30#include <dune/geometry/type.hh>
31#include <dune/geometry/multilineargeometry.hh>
32#include <dune/common/version.hh>
33#include <dune/common/reservedvector.hh>
34
38
39namespace Dumux {
40
51template<class GridView>
53{
54 using Grid = typename GridView::Grid;
55 static constexpr int dim = Grid::dimension;
56 static constexpr int dimWorld = Grid::dimensionworld;
57
58 // we use geometry traits that use static corner vectors to and a fixed geometry type
59 template <class ct>
60 struct ScvfMLGTraits : public Dune::MultiLinearGeometryTraits<ct>
61 {
62 // we use static vectors to store the corners as we know
63 // the number of corners in advance (2^(dim-1) corners (1<<(dim-1))
64 // However, on fracture scvs the number might be smaller (use ReservedVector)
65 template< int mydim, int cdim >
67 {
68 using Type = Dune::ReservedVector< Dune::FieldVector< ct, cdim >, (1<<(dim-1)) >;
69 };
70
71 // we know all scvfs will have the same geometry type
72 template< int mydim >
74 {
75 static const bool v = true;
76 static const unsigned int topologyId = Dune::Impl::CubeTopology< mydim >::type::id;
77 };
78 };
79
80 using GridIndexType = typename Grid::LeafGridView::IndexSet::IndexType;
81 using LocalIndexType = unsigned int;
82 using Scalar = typename Grid::ctype;
83 using Geometry = Dune::MultiLinearGeometry<Scalar, dim-1, dimWorld, ScvfMLGTraits<Scalar>>;
85 using GlobalPosition = typename CornerStorage::value_type;
87};
88
96template<class GV,
99: public SubControlVolumeFaceBase<BoxDfmSubControlVolumeFace<GV, T>, T>
100{
103 using GridIndexType = typename T::GridIndexType;
104 using LocalIndexType = typename T::LocalIndexType;
105 using Scalar = typename T::Scalar;
106 using GlobalPosition = typename T::GlobalPosition;
107 using CornerStorage = typename T::CornerStorage;
108 using Geometry = typename T::Geometry;
109 using BoundaryFlag = typename T::BoundaryFlag;
110
111 static_assert(T::dim == 2 || T::dim == 3, "Box-Dfm sub-control volume face only implemented in 2d or 3d");
112
113public:
115 using Traits = T;
116
119
121 template<class GeometryHelper, class Element>
122 BoxDfmSubControlVolumeFace(const GeometryHelper& geometryHelper,
123 const Element& element,
124 const typename Element::Geometry& elemGeometry,
125 GridIndexType scvfIndex,
126 std::vector<LocalIndexType>&& scvIndices)
127 : corners_(geometryHelper.getScvfCorners(scvfIndex))
128 , center_(0.0)
129 , unitOuterNormal_(geometryHelper.normal(corners_, scvIndices))
130 , area_(geometryHelper.scvfArea(corners_))
131 , scvfIndex_(scvfIndex)
132 , scvIndices_(std::move(scvIndices))
133 , boundary_(false)
134 , isFractureScvf_(false)
135 , boundaryFlag_{}
136 , facetIdx_(0)
137 {
138 for (const auto& corner : corners_)
139 center_ += corner;
140 center_ /= corners_.size();
141 }
142
144 template<class GeometryHelper, class Intersection>
145 BoxDfmSubControlVolumeFace(const GeometryHelper& geometryHelper,
146 const Intersection& intersection,
147 const typename Intersection::Geometry& isGeometry,
148 LocalIndexType indexInIntersection,
149 GridIndexType scvfIndex,
150 std::vector<LocalIndexType>&& scvIndices)
151 : corners_(geometryHelper.getBoundaryScvfCorners(intersection, isGeometry, indexInIntersection))
152 , center_(0.0)
153 , unitOuterNormal_(intersection.centerUnitOuterNormal())
154 , area_(geometryHelper.scvfArea(corners_))
155 , scvfIndex_(scvfIndex)
156 , scvIndices_(std::move(scvIndices))
157 , boundary_(true)
158 , isFractureScvf_(false)
159 , boundaryFlag_{intersection}
160 , facetIdx_(0)
161 {
162 for (const auto& corner : corners_)
163 center_ += corner;
164 center_ /= corners_.size();
165 }
166
168 template<class GeometryHelper, class Intersection>
169 BoxDfmSubControlVolumeFace(const GeometryHelper& geometryHelper,
170 const Intersection& intersection,
171 const typename Intersection::Geometry& isGeometry,
172 LocalIndexType indexInIntersection,
173 GridIndexType scvfIndex,
174 std::vector<LocalIndexType>&& scvIndices,
175 bool boundary)
176 : corners_(geometryHelper.getFractureScvfCorners(intersection, isGeometry, indexInIntersection))
177 , center_(0.0)
178 , scvfIndex_(scvfIndex)
179 , scvIndices_(std::move(scvIndices))
180 , boundary_(boundary)
181 , isFractureScvf_(true)
182 , boundaryFlag_{intersection}
183 , facetIdx_(intersection.indexInInside())
184 {
185 // The area here is given in meters. In order to
186 // get the right dimensions, the user has to provide
187 // the appropriate aperture in the problem (via an extrusion factor)
188 if (T::dim == 3)
189 area_ = (corners_[1]-corners_[0]).two_norm();
190 else if (T::dim == 2)
191 area_ = 1.0;
192
193 // obtain the unit normal vector
194 unitOuterNormal_ = geometryHelper.fractureNormal(corners_, intersection, indexInIntersection);
195
196 // compute the scvf center
197 for (const auto& corner : corners_)
198 center_ += corner;
199 center_ /= corners_.size();
200 }
201
203 const GlobalPosition& center() const
204 { return center_; }
205
207 const GlobalPosition& ipGlobal() const
208 { return center_; }
209
211 Scalar area() const
212 { return area_; }
213
215 bool boundary() const
216 { return boundary_; }
217
219 const GlobalPosition& unitOuterNormal() const
220 { return unitOuterNormal_; }
221
223 GridIndexType index() const
224 { return scvfIndex_; }
225
227 bool isOnFracture() const
228 { return isFractureScvf_; }
229
231 LocalIndexType facetIndexInElement() const
232 { assert(isFractureScvf_); return facetIdx_; }
233
236 { return boundaryFlag_.get(); }
237
239 LocalIndexType insideScvIdx() const
240 { return scvIndices_[0]; }
241
243 // This results in undefined behaviour if boundary is true
244 LocalIndexType outsideScvIdx() const
245 {
246 assert(!boundary());
247 return scvIndices_[1];
248 }
249
250 const GlobalPosition& corner(unsigned int localIdx) const
251 {
252 assert(localIdx < corners_.size() && "provided index exceeds the number of corners");
253 return corners_[localIdx];
254 }
255
257 Geometry geometry() const
258 {
259 if (isFractureScvf_)
260 DUNE_THROW(Dune::InvalidStateException, "The geometry object cannot be defined for fracture scvs "
261 "because the number of known corners is insufficient. "
262 "You can do this manually by extract the corners from this scv "
263 "and extrude them by the corresponding aperture. ");
264
265#if DUNE_VERSION_NEWER(DUNE_COMMON,2,6)
266 return Geometry(Dune::GeometryTypes::cube(Geometry::mydimension), corners_);
267#else
268 return Geometry(Dune::GeometryType(Dune::GeometryType::cube, Geometry::mydimension), corners_);
269#endif
270 }
271
272private:
273 CornerStorage corners_;
274 GlobalPosition center_;
275 GlobalPosition unitOuterNormal_;
276 Scalar area_;
277 GridIndexType scvfIndex_;
278 std::vector<LocalIndexType> scvIndices_;
279 bool boundary_;
280 bool isFractureScvf_;
281 BoundaryFlag boundaryFlag_;
282 LocalIndexType facetIdx_;
283};
284
285} // end namespace Dumux
286
287#endif
Boundary flag to store e.g. in sub control volume faces.
Base class for a sub control volume face.
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
std::size_t value_type
Definition: boundaryflag.hh:51
Definition: boundaryflag.hh:68
Base class for a sub control volume face, i.e a part of the boundary of a sub control volume we compu...
Definition: subcontrolvolumefacebase.hh:41
Default traits class to be used for the sub-control volume faces for the box discrete fracture scheme...
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:53
static constexpr int dim
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:55
Dune::MultiLinearGeometry< Scalar, dim-1, dimWorld, ScvfMLGTraits< Scalar > > Geometry
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:83
typename Grid::LeafGridView::IndexSet::IndexType GridIndexType
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:80
static constexpr int dimWorld
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:56
typename CornerStorage::value_type GlobalPosition
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:85
typename Grid::ctype Scalar
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:82
typename GridView::Grid Grid
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:54
unsigned int LocalIndexType
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:81
typename ScvfMLGTraits< Scalar >::template CornerStorage< dim-1, dimWorld >::Type CornerStorage
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:84
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:61
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:67
Dune::ReservedVector< Dune::FieldVector< ct, cdim >,(1<<(dim-1)) > Type
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:68
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:74
static const unsigned int topologyId
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:76
static const bool v
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:75
Class for a sub control volume face in the box discrete fracture method, i.e a part of the boundary o...
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:100
const GlobalPosition & corner(unsigned int localIdx) const
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:250
BoxDfmSubControlVolumeFace(const GeometryHelper &geometryHelper, const Intersection &intersection, const typename Intersection::Geometry &isGeometry, LocalIndexType indexInIntersection, GridIndexType scvfIndex, std::vector< LocalIndexType > &&scvIndices, bool boundary)
Constructor for inner fracture scvfs.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:169
const GlobalPosition & unitOuterNormal() const
returns the unit normal vector pointing outwards
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:219
BoxDfmSubControlVolumeFace(const GeometryHelper &geometryHelper, const Intersection &intersection, const typename Intersection::Geometry &isGeometry, LocalIndexType indexInIntersection, GridIndexType scvfIndex, std::vector< LocalIndexType > &&scvIndices)
Constructor for boundary scvfs.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:145
BoundaryFlag::value_type boundaryFlag() const
Returns the boundary flag.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:235
LocalIndexType insideScvIdx() const
Index of the inside sub control volume for spatial param evaluation.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:239
const GlobalPosition & ipGlobal() const
The integration point for flux evaluations in global coordinates.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:207
const GlobalPosition & center() const
The center of the sub control volume face.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:203
LocalIndexType outsideScvIdx() const
Index of the outside sub control volume for spatial param evaluation.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:244
LocalIndexType facetIndexInElement() const
The element-local facet index for which a fracture scv was created.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:231
BoxDfmSubControlVolumeFace()=default
The default constructor.
bool isOnFracture() const
Return if this is a fracture scvf.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:227
bool boundary() const
returns bolean if the sub control volume face is on the boundary
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:215
BoxDfmSubControlVolumeFace(const GeometryHelper &geometryHelper, const Element &element, const typename Element::Geometry &elemGeometry, GridIndexType scvfIndex, std::vector< LocalIndexType > &&scvIndices)
Constructor for inner scvfs.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:122
T Traits
State the traits public and thus export all types.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:115
Geometry geometry() const
The geometry of the sub control volume face.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:257
GridIndexType index() const
The global index of this sub control volume face.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:223
Scalar area() const
The area of the sub control volume face.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:211