version 3.8
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 using CornerStorage = typename T::CornerStorage;
73 enum { dim = Geometry::mydimension };
74
75 static_assert(dim == 2 || dim == 3, "Box-Dfm sub-control volume only implemented in 2d or 3d");
76
77public:
79 using Traits = T;
80
83
84 // the constructor for standard scvs
85 template<class GeometryHelper>
86 BoxDfmSubControlVolume(const GeometryHelper& geometryHelper,
87 LocalIndexType scvIdx,
88 GridIndexType elementIndex,
89 GridIndexType dofIndex)
90 : isFractureScv_(false)
91 , corners_(geometryHelper.getScvCorners(scvIdx))
92 , center_(0.0)
93 , volume_(Dumux::convexPolytopeVolume<T::dim>(
94 Dune::GeometryTypes::cube(T::dim),
95 [&](unsigned int i){ return corners_[i]; })
96 )
97 , elementIndex_(elementIndex)
98 , vIdxLocal_(scvIdx)
99 , elemLocalScvIdx_(scvIdx)
100 , dofIndex_(dofIndex)
101 , facetIdx_(0)
102 {
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 , corners_()
130 , center_(0.0)
131 , volume_(0.0)
132 , elementIndex_(elementIndex)
133 , vIdxLocal_(vIdxLocal)
134 , elemLocalScvIdx_(elemLocalScvIdx)
135 , dofIndex_(dofIndex)
136 , facetIdx_(elemLocalFacetIdx)
137 {
138 // copy corners
139 auto corners = geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), indexInIntersection);
140 const auto numCorners = corners.size();
141 corners_.resize(numCorners);
142 for (unsigned int i = 0; i < numCorners; ++i)
143 corners_[i] = corners[i];
144
145 // compute volume and scv center
146 volume_ = Dumux::convexPolytopeVolume<T::dim-1>(
147 Dune::GeometryTypes::cube(T::dim-1),
148 [&](unsigned int i){ return corners_[i]; }
149 );
150 for (const auto& corner : corners_)
151 center_ += corner;
152 center_ /= corners_.size();
153 }
154
156 const GlobalPosition& center() const
157 { return center_; }
158
160 Scalar volume() const
161 { return volume_; }
162
164 LocalIndexType localDofIndex() const
165 { return vIdxLocal_; }
166
168 LocalIndexType indexInElement() const
169 { return elemLocalScvIdx_; }
170
172 LocalIndexType facetIndexInElement() const
173 { assert(isFractureScv_); return facetIdx_; }
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 corners_[0]; }
182
184 GridIndexType elementIndex() const
185 { return elementIndex_; }
186
188 bool isOnFracture() const
189 { return isFractureScv_; }
190
192 // e.g. for integration
193 [[deprecated("Will be removed after 3.7. Use fvGeometry.geometry(scv).")]]
194 Geometry geometry() const
195 {
196 if (isFractureScv_)
197 DUNE_THROW(Dune::InvalidStateException, "The geometry object cannot be defined for fracture scvs "
198 "because the number of known corners is insufficient. "
199 "You can do this manually by extract the corners from this scv "
200 "and extrude them by the corresponding aperture. ");
201
202 return Geometry(Dune::GeometryTypes::cube(dim), corners_);
203 }
204
205private:
206 bool isFractureScv_;
207 CornerStorage corners_;
208 GlobalPosition center_;
209 Scalar volume_;
210 GridIndexType elementIndex_;
211 LocalIndexType vIdxLocal_;
212 LocalIndexType elemLocalScvIdx_;
213 GridIndexType dofIndex_;
214
215 // for fracture scvs only!
216 LocalIndexType facetIdx_;
217};
218
219} // end namespace
220
221#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:156
LocalIndexType facetIndexInElement() const
The element-local facet index for which a fracture scv was created.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:172
LocalIndexType indexInElement() const
The element-local index of this scv.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:168
Geometry geometry() const
The geometry of the sub control volume.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:194
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:160
const GlobalPosition & dofPosition() const
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:180
T Traits
State the traits public and thus export all types.
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:79
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:164
BoxDfmSubControlVolume(const GeometryHelper &geometryHelper, LocalIndexType scvIdx, GridIndexType elementIndex, GridIndexType dofIndex)
Definition: porousmediumflow/boxdfm/subcontrolvolume.hh:86
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.
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:220
Definition: adapt.hh:17
Definition: common/pdesolver.hh:24
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.