3.3.0
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/reservedvector.hh>
33
37
38namespace Dumux {
39
50template<class GridView>
52{
53 using Grid = typename GridView::Grid;
54 static constexpr int dim = Grid::dimension;
55 static constexpr int dimWorld = Grid::dimensionworld;
56
57 // we use geometry traits that use static corner vectors to and a fixed geometry type
58 template <class ct>
59 struct ScvfMLGTraits : public Dune::MultiLinearGeometryTraits<ct>
60 {
61 // we use static vectors to store the corners as we know
62 // the number of corners in advance (2^(dim-1) corners (1<<(dim-1))
63 // However, on fracture scvs the number might be smaller (use ReservedVector)
64 template< int mydim, int cdim >
66 {
67 using Type = Dune::ReservedVector< Dune::FieldVector< ct, cdim >, (1<<(dim-1)) >;
68 };
69
70 // we know all scvfs will have the same geometry type
71 template< int mydim >
73 {
74 static const bool v = true;
75 static const unsigned int topologyId = Dune::Impl::CubeTopology< mydim >::type::id;
76 };
77 };
78
79 using GridIndexType = typename Grid::LeafGridView::IndexSet::IndexType;
80 using LocalIndexType = unsigned int;
81 using Scalar = typename Grid::ctype;
82 using Geometry = Dune::MultiLinearGeometry<Scalar, dim-1, dimWorld, ScvfMLGTraits<Scalar>>;
84 using GlobalPosition = typename CornerStorage::value_type;
86};
87
95template<class GV,
98: public SubControlVolumeFaceBase<BoxDfmSubControlVolumeFace<GV, T>, T>
99{
102 using GridIndexType = typename T::GridIndexType;
103 using LocalIndexType = typename T::LocalIndexType;
104 using Scalar = typename T::Scalar;
105 using GlobalPosition = typename T::GlobalPosition;
106 using CornerStorage = typename T::CornerStorage;
107 using Geometry = typename T::Geometry;
108 using BoundaryFlag = typename T::BoundaryFlag;
109
110 static_assert(T::dim == 2 || T::dim == 3, "Box-Dfm sub-control volume face only implemented in 2d or 3d");
111
112public:
114 using Traits = T;
115
118
120 template<class GeometryHelper, class Element>
121 BoxDfmSubControlVolumeFace(const GeometryHelper& geometryHelper,
122 const Element& element,
123 const typename Element::Geometry& elemGeometry,
124 GridIndexType scvfIndex,
125 std::vector<LocalIndexType>&& scvIndices)
126 : corners_(geometryHelper.getScvfCorners(scvfIndex))
127 , center_(0.0)
128 , unitOuterNormal_(geometryHelper.normal(corners_, scvIndices))
129 , area_(geometryHelper.scvfArea(corners_))
130 , scvfIndex_(scvfIndex)
131 , scvIndices_(std::move(scvIndices))
132 , boundary_(false)
133 , isFractureScvf_(false)
134 , boundaryFlag_{}
135 , facetIdx_(0)
136 {
137 for (const auto& corner : corners_)
138 center_ += corner;
139 center_ /= corners_.size();
140 }
141
143 template<class GeometryHelper, class Intersection>
144 BoxDfmSubControlVolumeFace(const GeometryHelper& geometryHelper,
145 const Intersection& intersection,
146 const typename Intersection::Geometry& isGeometry,
147 LocalIndexType indexInIntersection,
148 GridIndexType scvfIndex,
149 std::vector<LocalIndexType>&& scvIndices)
150 : corners_(geometryHelper.getBoundaryScvfCorners(intersection, isGeometry, indexInIntersection))
151 , center_(0.0)
152 , unitOuterNormal_(intersection.centerUnitOuterNormal())
153 , area_(geometryHelper.scvfArea(corners_))
154 , scvfIndex_(scvfIndex)
155 , scvIndices_(std::move(scvIndices))
156 , boundary_(true)
157 , isFractureScvf_(false)
158 , boundaryFlag_{intersection}
159 , facetIdx_(0)
160 {
161 for (const auto& corner : corners_)
162 center_ += corner;
163 center_ /= corners_.size();
164 }
165
167 template<class GeometryHelper, class Intersection>
168 BoxDfmSubControlVolumeFace(const GeometryHelper& geometryHelper,
169 const Intersection& intersection,
170 const typename Intersection::Geometry& isGeometry,
171 LocalIndexType indexInIntersection,
172 GridIndexType scvfIndex,
173 std::vector<LocalIndexType>&& scvIndices,
174 bool boundary)
175 : corners_(geometryHelper.getFractureScvfCorners(intersection, isGeometry, indexInIntersection))
176 , center_(0.0)
177 , scvfIndex_(scvfIndex)
178 , scvIndices_(std::move(scvIndices))
179 , boundary_(boundary)
180 , isFractureScvf_(true)
181 , boundaryFlag_{intersection}
182 , facetIdx_(intersection.indexInInside())
183 {
184 // The area here is given in meters. In order to
185 // get the right dimensions, the user has to provide
186 // the appropriate aperture in the problem (via an extrusion factor)
187 if (T::dim == 3)
188 area_ = (corners_[1]-corners_[0]).two_norm();
189 else if (T::dim == 2)
190 area_ = 1.0;
191
192 // obtain the unit normal vector
193 unitOuterNormal_ = geometryHelper.fractureNormal(corners_, intersection, indexInIntersection);
194
195 // compute the scvf center
196 for (const auto& corner : corners_)
197 center_ += corner;
198 center_ /= corners_.size();
199 }
200
202 const GlobalPosition& center() const
203 { return center_; }
204
206 const GlobalPosition& ipGlobal() const
207 { return center_; }
208
210 Scalar area() const
211 { return area_; }
212
214 bool boundary() const
215 { return boundary_; }
216
218 const GlobalPosition& unitOuterNormal() const
219 { return unitOuterNormal_; }
220
222 GridIndexType index() const
223 { return scvfIndex_; }
224
226 bool isOnFracture() const
227 { return isFractureScvf_; }
228
230 LocalIndexType facetIndexInElement() const
231 { assert(isFractureScvf_); return facetIdx_; }
232
235 { return boundaryFlag_.get(); }
236
238 LocalIndexType insideScvIdx() const
239 { return scvIndices_[0]; }
240
242 // This results in undefined behaviour if boundary is true
243 LocalIndexType outsideScvIdx() const
244 {
245 assert(!boundary());
246 return scvIndices_[1];
247 }
248
249 const GlobalPosition& corner(unsigned int localIdx) const
250 {
251 assert(localIdx < corners_.size() && "provided index exceeds the number of corners");
252 return corners_[localIdx];
253 }
254
256 Geometry geometry() const
257 {
258 if (isFractureScvf_)
259 DUNE_THROW(Dune::InvalidStateException, "The geometry object cannot be defined for fracture scvs "
260 "because the number of known corners is insufficient. "
261 "You can do this manually by extract the corners from this scv "
262 "and extrude them by the corresponding aperture. ");
263
264 return Geometry(Dune::GeometryTypes::cube(Geometry::mydimension), corners_);
265 }
266
267private:
268 CornerStorage corners_;
269 GlobalPosition center_;
270 GlobalPosition unitOuterNormal_;
271 Scalar area_;
272 GridIndexType scvfIndex_;
273 std::vector<LocalIndexType> scvIndices_;
274 bool boundary_;
275 bool isFractureScvf_;
276 BoundaryFlag boundaryFlag_;
277 LocalIndexType facetIdx_;
278};
279
280} // end namespace Dumux
281
282#endif
Boundary flag to store e.g. in sub control volume faces.
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
Base class for a sub control volume face.
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: geometry/normal.hh:36
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:52
static constexpr int dim
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:54
Dune::MultiLinearGeometry< Scalar, dim-1, dimWorld, ScvfMLGTraits< Scalar > > Geometry
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:82
typename Grid::LeafGridView::IndexSet::IndexType GridIndexType
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:79
static constexpr int dimWorld
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:55
typename CornerStorage::value_type GlobalPosition
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:84
typename Grid::ctype Scalar
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:81
typename GridView::Grid Grid
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:53
unsigned int LocalIndexType
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:80
typename ScvfMLGTraits< Scalar >::template CornerStorage< dim-1, dimWorld >::Type CornerStorage
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:83
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:60
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:66
Dune::ReservedVector< Dune::FieldVector< ct, cdim >,(1<<(dim-1)) > Type
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:67
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:73
static const unsigned int topologyId
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:75
static const bool v
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:74
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:99
const GlobalPosition & corner(unsigned int localIdx) const
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:249
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:168
const GlobalPosition & unitOuterNormal() const
returns the unit normal vector pointing outwards
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:218
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:144
BoundaryFlag::value_type boundaryFlag() const
Returns the boundary flag.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:234
LocalIndexType insideScvIdx() const
Index of the inside sub control volume for spatial param evaluation.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:238
const GlobalPosition & ipGlobal() const
The integration point for flux evaluations in global coordinates.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:206
const GlobalPosition & center() const
The center of the sub control volume face.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:202
LocalIndexType outsideScvIdx() const
Index of the outside sub control volume for spatial param evaluation.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:243
LocalIndexType facetIndexInElement() const
The element-local facet index for which a fracture scv was created.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:230
BoxDfmSubControlVolumeFace()=default
The default constructor.
bool isOnFracture() const
Return if this is a fracture scvf.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:226
bool boundary() const
returns bolean if the sub control volume face is on the boundary
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:214
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:121
T Traits
State the traits public and thus export all types.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:114
Geometry geometry() const
The geometry of the sub control volume face.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:256
GridIndexType index() const
The global index of this sub control volume face.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:222
Scalar area() const
The area of the sub control volume face.
Definition: porousmediumflow/boxdfm/subcontrolvolumeface.hh:210