version 3.11-dev
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// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_POROUSMEDIUMFLOW_BOXDFM_SUBCONTROLVOLUME_HH
14#define DUMUX_POROUSMEDIUMFLOW_BOXDFM_SUBCONTROLVOLUME_HH
15
16#include <dune/common/reservedvector.hh>
17#include <dune/geometry/type.hh>
18#include <dune/geometry/multilineargeometry.hh>
19
21#include <dumux/common/math.hh>
23
24namespace Dumux {
25
36template<class GridView>
38{
39 using Grid = typename GridView::Grid;
40
41 static const int dim = Grid::dimension;
42 static const int dimWorld = Grid::dimensionworld;
43 using GridIndexType = typename Grid::LeafGridView::IndexSet::IndexType;
44 using LocalIndexType = unsigned int;
45 using Scalar = typename Grid::ctype;
47 using Geometry = Dune::MultiLinearGeometry<Scalar, dim, dimWorld, GeometryTraits>;
48 using CornerStorage = typename GeometryTraits::template CornerStorage<dim, dimWorld>::Type;
49 using GlobalPosition = typename CornerStorage::value_type;
50};
51
59template<class GV,
62{
64 using Geometry = typename T::Geometry;
65 using GridIndexType = typename T::GridIndexType;
66 using LocalIndexType = typename T::LocalIndexType;
67 using Scalar = typename T::Scalar;
68 static constexpr int dim = Geometry::mydimension;
69
70 static_assert(dim == 2 || dim == 3, "Box-Dfm sub-control volume only implemented in 2d or 3d");
71
72public:
74 using GlobalPosition = typename T::GlobalPosition;
76 using Traits = T;
77
80
81 // the constructor for standard scvs
82 template<class GeometryHelper>
83 BoxDfmSubControlVolume(const GeometryHelper& geometryHelper,
84 LocalIndexType scvIdx,
85 GridIndexType elementIndex,
86 GridIndexType dofIndex)
87 : isFractureScv_(false)
88 , center_(0.0)
89 , elementIndex_(elementIndex)
90 , vIdxLocal_(scvIdx)
91 , elemLocalScvIdx_(scvIdx)
92 , dofIndex_(dofIndex)
93 , facetIdx_(0)
94 , indexInIntersection_(0)
95 {
96 const auto corners = geometryHelper.getScvCorners(scvIdx);
97 dofPosition_ = corners[0];
98 volume_ = Dumux::convexPolytopeVolume<T::dim>(
99 Dune::GeometryTypes::cube(T::dim),
100 [&](unsigned int i){ return corners[i]; });
101 // compute center point
102 for (const auto& corner : corners)
103 center_ += corner;
104 center_ /= corners.size();
105 }
106
116 template<class GeometryHelper, class Intersection>
117 BoxDfmSubControlVolume(const GeometryHelper& geometryHelper,
118 const Intersection& intersection,
119 const typename Intersection::Geometry& isGeometry,
120 LocalIndexType indexInIntersection,
121 LocalIndexType vIdxLocal,
122 LocalIndexType elemLocalScvIdx,
123 LocalIndexType elemLocalFacetIdx,
124 GridIndexType elementIndex,
125 GridIndexType dofIndex)
126 : isFractureScv_(true)
127 , center_(0.0)
128 , volume_(0.0)
129 , elementIndex_(elementIndex)
130 , vIdxLocal_(vIdxLocal)
131 , elemLocalScvIdx_(elemLocalScvIdx)
132 , dofIndex_(dofIndex)
133 , facetIdx_(elemLocalFacetIdx)
134 , indexInIntersection_(indexInIntersection)
135 {
136 const auto corners = geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), indexInIntersection);
137 dofPosition_ = corners[0];
138
139 // compute volume and scv center
140 volume_ = Dumux::convexPolytopeVolume<T::dim-1>(
141 Dune::GeometryTypes::cube(T::dim-1),
142 [&](unsigned int i){ return corners[i]; }
143 );
144 for (const auto& corner : corners)
145 center_ += corner;
146 center_ /= corners.size();
147 }
148
150 const GlobalPosition& center() const
151 { return center_; }
152
154 Scalar volume() const
155 { return volume_; }
156
158 LocalIndexType localDofIndex() const
159 { return vIdxLocal_; }
160
162 LocalIndexType indexInElement() const
163 { return elemLocalScvIdx_; }
164
166 LocalIndexType facetIndexInElement() const
167 { assert(isFractureScv_); return facetIdx_; }
168
170 LocalIndexType indexInsideIntersection() const
171 { assert(isFractureScv_); return indexInIntersection_; }
172
174 GridIndexType dofIndex() const
175 { return dofIndex_; }
176
177 // The position of the dof this scv is embedded in (list is defined such that first entry is vertex itself)
179 { return dofPosition_; }
180
182 GridIndexType elementIndex() const
183 { return elementIndex_; }
184
186 bool isOnFracture() const
187 { return isFractureScv_; }
188
189private:
190 bool isFractureScv_;
191 GlobalPosition dofPosition_;
192 GlobalPosition center_;
193 Scalar volume_;
194 GridIndexType elementIndex_;
195 LocalIndexType vIdxLocal_;
196 LocalIndexType elemLocalScvIdx_;
197 GridIndexType dofIndex_;
198
199 // for fracture scvs only!
200 LocalIndexType facetIdx_;
201 LocalIndexType indexInIntersection_;
202};
203
204} // end namespace
205
206#endif
the sub control volume for the box discrete fracture scheme
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:62
BoxDfmSubControlVolume()=default
The default constructor.
const GlobalPosition & center() const
The center of the sub control volume.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:150
LocalIndexType facetIndexInElement() const
The element-local facet index for which a fracture scv was created.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:166
LocalIndexType indexInElement() const
The element-local index of this scv.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:162
bool isOnFracture() const
Return true if this scv is part of the fracture domain.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:186
typename T::GlobalPosition GlobalPosition
export the type used for global coordinates
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:74
Scalar volume() const
The volume of the sub control volume.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:154
const GlobalPosition & dofPosition() const
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:178
LocalIndexType indexInsideIntersection() const
The local vertex index in the intersection.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:170
T Traits
State the traits public and thus export all types.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:76
GridIndexType elementIndex() const
The global index of the element this scv is embedded in.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:182
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:117
LocalIndexType localDofIndex() const
The element-local vertex index this scv is connected to.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:158
BoxDfmSubControlVolume(const GeometryHelper &geometryHelper, LocalIndexType scvIdx, GridIndexType elementIndex, GridIndexType dofIndex)
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:83
GridIndexType dofIndex() const
The index of the dof this scv is embedded in.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:174
auto convexPolytopeVolume(Dune::GeometryType type, const CornerF &c)
Compute the volume of several common geometry types.
Definition: volume.hh:41
Define some often used mathematical functions.
Definition: adapt.hh:17
Helper class constructing the dual grid finite volume geometries for the box discrete fracture model.
Default traits class to be used for the sub-control volumes for the box discrete fracture scheme.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:38
static const int dimWorld
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:42
static const int dim
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:41
typename GridView::Grid Grid
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:39
typename Grid::LeafGridView::IndexSet::IndexType GridIndexType
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:43
typename CornerStorage::value_type GlobalPosition
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:49
Dune::MultiLinearGeometry< Scalar, dim, dimWorld, GeometryTraits > Geometry
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:47
typename GeometryTraits::template CornerStorage< dim, dimWorld >::Type CornerStorage
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:48
unsigned int LocalIndexType
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:44
typename Grid::ctype Scalar
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:45
Definition: porousmediumflow/boxdfm/geometryhelper.hh:26
Compute the volume of several common geometry types.