version 3.10-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-FileCopyrightInfo: 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
22#include <dumux/common/math.hh>
24
25namespace Dumux {
26
37template<class GridView>
39{
40 using Grid = typename GridView::Grid;
41
42 static const int dim = Grid::dimension;
43 static const int dimWorld = Grid::dimensionworld;
44 using GridIndexType = typename Grid::LeafGridView::IndexSet::IndexType;
45 using LocalIndexType = unsigned int;
46 using Scalar = typename Grid::ctype;
48 using Geometry = Dune::MultiLinearGeometry<Scalar, dim, dimWorld, GeometryTraits>;
49 using CornerStorage = typename GeometryTraits::template CornerStorage<dim, dimWorld>::Type;
50 using GlobalPosition = typename CornerStorage::value_type;
51};
52
60template<class GV,
63: public SubControlVolumeBase<BoxDfmSubControlVolume<GV, T>, T>
64{
67 using Geometry = typename T::Geometry;
68 using GridIndexType = typename T::GridIndexType;
69 using LocalIndexType = typename T::LocalIndexType;
70 using Scalar = typename T::Scalar;
71 using GlobalPosition = typename T::GlobalPosition;
72 enum { dim = Geometry::mydimension };
73
74 static_assert(dim == 2 || dim == 3, "Box-Dfm sub-control volume only implemented in 2d or 3d");
75
76public:
78 using Traits = T;
79
82
83 // the constructor for standard scvs
84 template<class GeometryHelper>
85 BoxDfmSubControlVolume(const GeometryHelper& geometryHelper,
86 LocalIndexType scvIdx,
87 GridIndexType elementIndex,
88 GridIndexType dofIndex)
89 : isFractureScv_(false)
90 , center_(0.0)
91 , elementIndex_(elementIndex)
92 , vIdxLocal_(scvIdx)
93 , elemLocalScvIdx_(scvIdx)
94 , dofIndex_(dofIndex)
95 , facetIdx_(0)
96 , indexInIntersection_(0)
97 {
98 const auto corners = geometryHelper.getScvCorners(scvIdx);
99 dofPosition_ = corners[0];
100 volume_ = Dumux::convexPolytopeVolume<T::dim>(
101 Dune::GeometryTypes::cube(T::dim),
102 [&](unsigned int i){ return corners[i]; });
103 // compute center point
104 for (const auto& corner : corners)
105 center_ += corner;
106 center_ /= corners.size();
107 }
108
118 template<class GeometryHelper, class Intersection>
119 BoxDfmSubControlVolume(const GeometryHelper& geometryHelper,
120 const Intersection& intersection,
121 const typename Intersection::Geometry& isGeometry,
122 LocalIndexType indexInIntersection,
123 LocalIndexType vIdxLocal,
124 LocalIndexType elemLocalScvIdx,
125 LocalIndexType elemLocalFacetIdx,
126 GridIndexType elementIndex,
127 GridIndexType dofIndex)
128 : isFractureScv_(true)
129 , center_(0.0)
130 , volume_(0.0)
131 , elementIndex_(elementIndex)
132 , vIdxLocal_(vIdxLocal)
133 , elemLocalScvIdx_(elemLocalScvIdx)
134 , dofIndex_(dofIndex)
135 , facetIdx_(elemLocalFacetIdx)
136 , indexInIntersection_(indexInIntersection)
137 {
138 const auto corners = geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), indexInIntersection);
139 dofPosition_ = corners[0];
140
141 // compute volume and scv center
142 volume_ = Dumux::convexPolytopeVolume<T::dim-1>(
143 Dune::GeometryTypes::cube(T::dim-1),
144 [&](unsigned int i){ return corners[i]; }
145 );
146 for (const auto& corner : corners)
147 center_ += corner;
148 center_ /= corners.size();
149 }
150
152 const GlobalPosition& center() const
153 { return center_; }
154
156 Scalar volume() const
157 { return volume_; }
158
160 LocalIndexType localDofIndex() const
161 { return vIdxLocal_; }
162
164 LocalIndexType indexInElement() const
165 { return elemLocalScvIdx_; }
166
168 LocalIndexType facetIndexInElement() const
169 { assert(isFractureScv_); return facetIdx_; }
170
172 LocalIndexType indexInsideIntersection() const
173 { assert(isFractureScv_); return indexInIntersection_; }
174
176 GridIndexType dofIndex() const
177 { return dofIndex_; }
178
179 // The position of the dof this scv is embedded in (list is defined such that first entry is vertex itself)
180 const GlobalPosition& dofPosition() const
181 { return dofPosition_; }
182
184 GridIndexType elementIndex() const
185 { return elementIndex_; }
186
188 bool isOnFracture() const
189 { return isFractureScv_; }
190
191private:
192 bool isFractureScv_;
193 GlobalPosition dofPosition_;
194 GlobalPosition center_;
195 Scalar volume_;
196 GridIndexType elementIndex_;
197 LocalIndexType vIdxLocal_;
198 LocalIndexType elemLocalScvIdx_;
199 GridIndexType dofIndex_;
200
201 // for fracture scvs only!
202 LocalIndexType facetIdx_;
203 LocalIndexType indexInIntersection_;
204};
205
206} // end namespace
207
208#endif
the sub control volume for the box discrete fracture scheme
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:64
BoxDfmSubControlVolume()=default
The default constructor.
const GlobalPosition & center() const
The center of the sub control volume.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:152
LocalIndexType facetIndexInElement() const
The element-local facet index for which a fracture scv was created.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:168
LocalIndexType indexInElement() const
The element-local index of this scv.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:164
bool isOnFracture() const
Return true if this scv is part of the fracture domain.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:188
Scalar volume() const
The volume of the sub control volume.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:156
const GlobalPosition & dofPosition() const
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:180
LocalIndexType indexInsideIntersection() const
The local vertex index in the intersection.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:172
T Traits
State the traits public and thus export all types.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:78
GridIndexType elementIndex() const
The global index of the element this scv is embedded in.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:184
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:119
LocalIndexType localDofIndex() const
The element-local vertex index this scv is connected to.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:160
BoxDfmSubControlVolume(const GeometryHelper &geometryHelper, LocalIndexType scvIdx, GridIndexType elementIndex, GridIndexType dofIndex)
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:85
GridIndexType dofIndex() const
The index of the dof this scv is embedded in.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:176
Base class for a sub control volume, i.e a part of the control volume we are making the balance for....
Definition: subcontrolvolumebase.hh:26
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:39
static const int dimWorld
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:43
static const int dim
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:42
typename GridView::Grid Grid
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:40
typename Grid::LeafGridView::IndexSet::IndexType GridIndexType
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:44
typename CornerStorage::value_type GlobalPosition
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:50
Dune::MultiLinearGeometry< Scalar, dim, dimWorld, GeometryTraits > Geometry
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:48
typename GeometryTraits::template CornerStorage< dim, dimWorld >::Type CornerStorage
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:49
unsigned int LocalIndexType
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:45
typename Grid::ctype Scalar
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:46
Definition: porousmediumflow/boxdfm/geometryhelper.hh:26
Base class for a sub control volume.
Compute the volume of several common geometry types.