3.6-git
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/type.hh>
30#include <dune/geometry/multilineargeometry.hh>
31
34#include <dumux/common/math.hh>
36
37namespace Dumux {
38
49template<class GridView>
51{
52 using Grid = typename GridView::Grid;
53
54 static const int dim = Grid::dimension;
55 static const int dimWorld = Grid::dimensionworld;
56 using GridIndexType = typename Grid::LeafGridView::IndexSet::IndexType;
57 using LocalIndexType = unsigned int;
58 using Scalar = typename Grid::ctype;
60 using Geometry = Dune::MultiLinearGeometry<Scalar, dim, dimWorld, GeometryTraits>;
61 using CornerStorage = typename GeometryTraits::template CornerStorage<dim, dimWorld>::Type;
62 using GlobalPosition = typename CornerStorage::value_type;
63};
64
72template<class GV,
75: public SubControlVolumeBase<BoxDfmSubControlVolume<GV, T>, T>
76{
79 using Geometry = typename T::Geometry;
80 using GridIndexType = typename T::GridIndexType;
81 using LocalIndexType = typename T::LocalIndexType;
82 using Scalar = typename T::Scalar;
83 using GlobalPosition = typename T::GlobalPosition;
84 using CornerStorage = typename T::CornerStorage;
85 enum { dim = Geometry::mydimension };
86
87 static_assert(dim == 2 || dim == 3, "Box-Dfm sub-control volume only implemented in 2d or 3d");
88
89public:
91 using Traits = T;
92
95
96 // the constructor for standard scvs
97 template<class GeometryHelper>
98 BoxDfmSubControlVolume(const GeometryHelper& geometryHelper,
99 LocalIndexType scvIdx,
100 GridIndexType elementIndex,
101 GridIndexType dofIndex)
102 : isFractureScv_(false)
103 , corners_(geometryHelper.getScvCorners(scvIdx))
104 , center_(0.0)
105 , volume_(Dumux::convexPolytopeVolume<T::dim>(
106 Dune::GeometryTypes::cube(T::dim),
107 [&](unsigned int i){ return corners_[i]; })
108 )
109 , elementIndex_(elementIndex)
110 , vIdxLocal_(scvIdx)
111 , elemLocalScvIdx_(scvIdx)
112 , dofIndex_(dofIndex)
113 , facetIdx_(0)
114 {
115 // compute center point
116 for (const auto& corner : corners_)
117 center_ += corner;
118 center_ /= corners_.size();
119 }
120
130 template<class GeometryHelper, class Intersection>
131 BoxDfmSubControlVolume(const GeometryHelper& geometryHelper,
132 const Intersection& intersection,
133 const typename Intersection::Geometry& isGeometry,
134 LocalIndexType indexInIntersection,
135 LocalIndexType vIdxLocal,
136 LocalIndexType elemLocalScvIdx,
137 LocalIndexType elemLocalFacetIdx,
138 GridIndexType elementIndex,
139 GridIndexType dofIndex)
140 : isFractureScv_(true)
141 , corners_()
142 , center_(0.0)
143 , volume_(0.0)
144 , elementIndex_(elementIndex)
145 , vIdxLocal_(vIdxLocal)
146 , elemLocalScvIdx_(elemLocalScvIdx)
147 , dofIndex_(dofIndex)
148 , facetIdx_(elemLocalFacetIdx)
149 {
150 // copy corners
151 auto corners = geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), indexInIntersection);
152 const auto numCorners = corners.size();
153 corners_.resize(numCorners);
154 for (unsigned int i = 0; i < numCorners; ++i)
155 corners_[i] = corners[i];
156
157 // compute volume and scv center
158 volume_ = Dumux::convexPolytopeVolume<T::dim-1>(
159 Dune::GeometryTypes::cube(T::dim-1),
160 [&](unsigned int i){ return corners_[i]; }
161 );
162 for (const auto& corner : corners_)
163 center_ += corner;
164 center_ /= corners_.size();
165 }
166
168 const GlobalPosition& center() const
169 { return center_; }
170
172 Scalar volume() const
173 { return volume_; }
174
176 LocalIndexType localDofIndex() const
177 { return vIdxLocal_; }
178
180 LocalIndexType indexInElement() const
181 { return elemLocalScvIdx_; }
182
184 LocalIndexType facetIndexInElement() const
185 { assert(isFractureScv_); return facetIdx_; }
186
188 GridIndexType dofIndex() const
189 { return dofIndex_; }
190
191 // The position of the dof this scv is embedded in (list is defined such that first entry is vertex itself)
192 const GlobalPosition& dofPosition() const
193 { return corners_[0]; }
194
196 GridIndexType elementIndex() const
197 { return elementIndex_; }
198
200 bool isOnFracture() const
201 { return isFractureScv_; }
202
204 // e.g. for integration
205 Geometry geometry() const
206 {
207 if (isFractureScv_)
208 DUNE_THROW(Dune::InvalidStateException, "The geometry object cannot be defined for fracture scvs "
209 "because the number of known corners is insufficient. "
210 "You can do this manually by extract the corners from this scv "
211 "and extrude them by the corresponding aperture. ");
212
213 return Geometry(Dune::GeometryTypes::cube(dim), corners_);
214 }
215
217 const GlobalPosition& corner(LocalIndexType localIdx) const
218 {
219 assert(localIdx < corners_.size() && "provided index exceeds the number of corners");
220 return corners_[localIdx];
221 }
222
223private:
224 bool isFractureScv_;
225 CornerStorage corners_;
226 GlobalPosition center_;
227 Scalar volume_;
228 GridIndexType elementIndex_;
229 LocalIndexType vIdxLocal_;
230 LocalIndexType elemLocalScvIdx_;
231 GridIndexType dofIndex_;
232
233 // for fracture scvs only!
234 LocalIndexType facetIdx_;
235};
236
237} // end namespace
238
239#endif
Define some often used mathematical functions.
Base class for a sub control volume.
Compute the volume of several common geometry types.
auto convexPolytopeVolume(Dune::GeometryType type, const CornerF &c)
Compute the volume of several common geometry types.
Definition: volume.hh:53
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
Definition: deprecated.hh:149
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:232
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
Definition: porousmediumflow/boxdfm/geometryhelper.hh:38
Default traits class to be used for the sub-control volumes for the box discrete fracture scheme.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:51
static const int dimWorld
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:55
static const int dim
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:54
typename GridView::Grid Grid
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:52
typename Grid::LeafGridView::IndexSet::IndexType GridIndexType
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:56
typename CornerStorage::value_type GlobalPosition
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:62
Dune::MultiLinearGeometry< Scalar, dim, dimWorld, GeometryTraits > Geometry
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:60
typename GeometryTraits::template CornerStorage< dim, dimWorld >::Type CornerStorage
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:61
unsigned int LocalIndexType
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:57
typename Grid::ctype Scalar
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:58
the sub control volume for the box discrete fracture scheme
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:76
BoxDfmSubControlVolume()=default
The default constructor.
const GlobalPosition & center() const
The center of the sub control volume.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:168
LocalIndexType facetIndexInElement() const
The element-local facet index for which a fracture scv was created.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:184
LocalIndexType indexInElement() const
The element-local index of this scv.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:180
Geometry geometry() const
The geometry of the sub control volume.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:205
bool isOnFracture() const
Return true if this scv is part of the fracture domain.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:200
const GlobalPosition & corner(LocalIndexType localIdx) const
Return the corner for the given local index.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:217
Scalar volume() const
The volume of the sub control volume.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:172
const GlobalPosition & dofPosition() const
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:192
T Traits
State the traits public and thus export all types.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:91
GridIndexType elementIndex() const
The global index of the element this scv is embedded in.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:196
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:131
LocalIndexType localDofIndex() const
The element-local vertex index this scv is connected to.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:176
BoxDfmSubControlVolume(const GeometryHelper &geometryHelper, LocalIndexType scvIdx, GridIndexType elementIndex, GridIndexType dofIndex)
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:98
GridIndexType dofIndex() const
The index of the dof this scv is embedded in.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:188
Helper class constructing the dual grid finite volume geometries for the box discrete fracture model.