Line | Branch | Exec | Source |
---|---|---|---|
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 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup BoxDFMModel | ||
10 | * \brief The sub control volume face class for the box discrete fracture model. | ||
11 | */ | ||
12 | |||
13 | #ifndef DUMUX_POROUSMEDIUMFLOW_BOXDFM_SUBCONTROLVOLUMEFACE_HH | ||
14 | #define DUMUX_POROUSMEDIUMFLOW_BOXDFM_SUBCONTROLVOLUMEFACE_HH | ||
15 | |||
16 | #include <utility> | ||
17 | |||
18 | #include <dune/geometry/type.hh> | ||
19 | #include <dune/geometry/multilineargeometry.hh> | ||
20 | |||
21 | #include <dumux/common/boundaryflag.hh> | ||
22 | #include <dumux/discretization/subcontrolvolumefacebase.hh> | ||
23 | #include <dumux/porousmediumflow/boxdfm/geometryhelper.hh> | ||
24 | #include <dumux/geometry/volume.hh> | ||
25 | |||
26 | namespace Dumux { | ||
27 | |||
28 | /*! | ||
29 | * \ingroup BoxDFMModel | ||
30 | * \brief Default traits class to be used for the sub-control volume faces | ||
31 | * for the box discrete fracture scheme | ||
32 | * | ||
33 | * \tparam GV the type of the grid view | ||
34 | * | ||
35 | * \note We define new traits for the box-dfm sub-control volume face | ||
36 | * as we use a different type of container for storing the scvf corners! | ||
37 | */ | ||
38 | template<class GridView> | ||
39 | struct BoxDfmDefaultScvfGeometryTraits | ||
40 | { | ||
41 | using Grid = typename GridView::Grid; | ||
42 | static constexpr int dim = Grid::dimension; | ||
43 | static constexpr int dimWorld = Grid::dimensionworld; | ||
44 | using GridIndexType = typename Grid::LeafGridView::IndexSet::IndexType; | ||
45 | using LocalIndexType = unsigned int; | ||
46 | using Scalar = typename Grid::ctype; | ||
47 | using GeometryTraits = BoxDfmMLGeometryTraits<Scalar>; | ||
48 | using Geometry = Dune::MultiLinearGeometry<Scalar, dim-1, dimWorld, GeometryTraits>; | ||
49 | using CornerStorage = typename GeometryTraits::template CornerStorage<dim-1, dimWorld>::Type; | ||
50 | using GlobalPosition = typename CornerStorage::value_type; | ||
51 | using BoundaryFlag = Dumux::BoundaryFlag<Grid>; | ||
52 | }; | ||
53 | |||
54 | /*! | ||
55 | * \ingroup BoxDFMModel | ||
56 | * \brief Class for a sub control volume face in the box discrete fracture method, i.e a | ||
57 | * part of the boundary of a sub control volume we compute fluxes on. | ||
58 | * \tparam GV the type of the grid view | ||
59 | * \tparam T the scvf geometry traits | ||
60 | */ | ||
61 | template<class GV, | ||
62 | class T = BoxDfmDefaultScvfGeometryTraits<GV> > | ||
63 | 20578442 | class BoxDfmSubControlVolumeFace | |
64 | : public SubControlVolumeFaceBase<BoxDfmSubControlVolumeFace<GV, T>, T> | ||
65 | { | ||
66 | using ThisType = BoxDfmSubControlVolumeFace<GV, T>; | ||
67 | using ParentType = SubControlVolumeFaceBase<ThisType, T>; | ||
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 Geometry = typename T::Geometry; | ||
73 | using BoundaryFlag = typename T::BoundaryFlag; | ||
74 | |||
75 | static_assert(T::dim == 2 || T::dim == 3, "Box-Dfm sub-control volume face only implemented in 2d or 3d"); | ||
76 | |||
77 | public: | ||
78 | //! State the traits public and thus export all types | ||
79 | using Traits = T; | ||
80 | |||
81 | //! The default constructor | ||
82 | 8710636 | BoxDfmSubControlVolumeFace() = default; | |
83 | |||
84 | //! Constructor for inner scvfs | ||
85 | template<class GeometryHelper, class Element> | ||
86 | 8791080 | BoxDfmSubControlVolumeFace(const GeometryHelper& geometryHelper, | |
87 | const Element& element, | ||
88 | const typename Element::Geometry& elemGeometry, | ||
89 | GridIndexType scvfIndex, | ||
90 | std::vector<LocalIndexType>&& scvIndices) | ||
91 | 8791080 | : center_(0.0) | |
92 | 8791080 | , scvfIndex_(scvfIndex) | |
93 |
1/2✓ Branch 1 taken 8791080 times.
✗ Branch 2 not taken.
|
8791080 | , scvIndices_(std::move(scvIndices)) |
94 | 8791080 | , boundary_(false) | |
95 | 8791080 | , isFractureScvf_(false) | |
96 |
1/2✓ Branch 1 taken 8791080 times.
✗ Branch 2 not taken.
|
8791080 | , boundaryFlag_{} |
97 | 8791080 | , facetIdx_(0) | |
98 |
1/2✓ Branch 1 taken 8791080 times.
✗ Branch 2 not taken.
|
8791080 | , indexInIntersection_(0) |
99 | { | ||
100 |
1/2✓ Branch 1 taken 8791080 times.
✗ Branch 2 not taken.
|
8791080 | const auto corners = geometryHelper.getScvfCorners(scvfIndex); |
101 |
1/2✓ Branch 1 taken 4410000 times.
✗ Branch 2 not taken.
|
8791080 | unitOuterNormal_ = geometryHelper.normal(corners, scvIndices_); |
102 | 17582160 | area_ = Dumux::convexPolytopeVolume<T::dim-1>( | |
103 | Dune::GeometryTypes::cube(T::dim-1), | ||
104 | 20347488 | [&](unsigned int i){ return corners[i]; }); | |
105 | |||
106 |
2/2✓ Branch 0 taken 23112816 times.
✓ Branch 1 taken 8791080 times.
|
31903896 | for (const auto& corner : corners) |
107 |
2/2✓ Branch 0 taken 57286944 times.
✓ Branch 1 taken 23112816 times.
|
80399760 | center_ += corner; |
108 | 8791080 | center_ /= corners.size(); | |
109 | 8791080 | } | |
110 | |||
111 | //! Constructor for boundary scvfs | ||
112 | template<class GeometryHelper, class Intersection> | ||
113 | 705888 | BoxDfmSubControlVolumeFace(const GeometryHelper& geometryHelper, | |
114 | const Intersection& intersection, | ||
115 | const typename Intersection::Geometry& isGeometry, | ||
116 | LocalIndexType indexInIntersection, | ||
117 | GridIndexType scvfIndex, | ||
118 | std::vector<LocalIndexType>&& scvIndices) | ||
119 | 1411776 | : center_(0.0) | |
120 | 705888 | , unitOuterNormal_(intersection.centerUnitOuterNormal()) | |
121 | 705888 | , scvfIndex_(scvfIndex) | |
122 |
1/2✓ Branch 1 taken 705888 times.
✗ Branch 2 not taken.
|
705888 | , scvIndices_(std::move(scvIndices)) |
123 | 705888 | , boundary_(true) | |
124 | 705888 | , isFractureScvf_(false) | |
125 |
1/2✓ Branch 1 taken 705888 times.
✗ Branch 2 not taken.
|
705888 | , boundaryFlag_{intersection} |
126 | 705888 | , facetIdx_(0) | |
127 | 705888 | , indexInIntersection_(0) | |
128 | { | ||
129 |
2/4✓ Branch 1 taken 705888 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 352604 times.
✗ Branch 5 not taken.
|
705888 | const auto corners = geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), indexInIntersection); |
130 | 1411776 | area_ = Dumux::convexPolytopeVolume<T::dim-1>( | |
131 | Dune::GeometryTypes::cube(T::dim-1), | ||
132 | 1926936 | [&](unsigned int i){ return corners[i]; }); | |
133 |
2/2✓ Branch 0 taken 2442096 times.
✓ Branch 1 taken 705888 times.
|
3147984 | for (const auto& corner : corners) |
134 |
2/2✓ Branch 0 taken 6944832 times.
✓ Branch 1 taken 2442096 times.
|
9386928 | center_ += corner; |
135 | 705888 | center_ /= corners.size(); | |
136 | 705888 | } | |
137 | |||
138 | //! Constructor for inner fracture scvfs | ||
139 | template<class GeometryHelper, class Intersection> | ||
140 | 336696 | BoxDfmSubControlVolumeFace(const GeometryHelper& geometryHelper, | |
141 | const Intersection& intersection, | ||
142 | const typename Intersection::Geometry& isGeometry, | ||
143 | LocalIndexType indexInIntersection, | ||
144 | GridIndexType scvfIndex, | ||
145 | std::vector<LocalIndexType>&& scvIndices, | ||
146 | bool boundary) | ||
147 | 336696 | : center_(0.0) | |
148 | 336696 | , scvfIndex_(scvfIndex) | |
149 |
1/2✓ Branch 1 taken 336696 times.
✗ Branch 2 not taken.
|
336696 | , scvIndices_(std::move(scvIndices)) |
150 | 336696 | , boundary_(boundary) | |
151 | 336696 | , isFractureScvf_(true) | |
152 |
2/4✓ Branch 1 taken 336696 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 336696 times.
✗ Branch 5 not taken.
|
336696 | , boundaryFlag_{intersection} |
153 | 336696 | , facetIdx_(intersection.indexInInside()) | |
154 |
2/4✓ Branch 1 taken 336696 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 167848 times.
✗ Branch 5 not taken.
|
504544 | , indexInIntersection_(indexInIntersection) |
155 | { | ||
156 |
1/2✓ Branch 1 taken 336696 times.
✗ Branch 2 not taken.
|
336696 | const auto corners = geometryHelper.getFractureScvfCorners(intersection.indexInInside(), indexInIntersection); |
157 | // The area here is given in meters. In order to | ||
158 | // get the right dimensions, the user has to provide | ||
159 | // the appropriate aperture in the problem (via an extrusion factor) | ||
160 | if (T::dim == 3) | ||
161 | 279840 | area_ = (corners[1]-corners[0]).two_norm(); | |
162 | else if (T::dim == 2) | ||
163 | 196776 | area_ = 1.0; | |
164 | |||
165 | // obtain the unit normal vector | ||
166 |
1/2✓ Branch 1 taken 336696 times.
✗ Branch 2 not taken.
|
336696 | unitOuterNormal_ = geometryHelper.fractureNormal(corners, intersection, indexInIntersection); |
167 | |||
168 | // compute the scvf center | ||
169 |
2/2✓ Branch 0 taken 476616 times.
✓ Branch 1 taken 336696 times.
|
813312 | for (const auto& corner : corners) |
170 |
2/2✓ Branch 0 taken 1233072 times.
✓ Branch 1 taken 476616 times.
|
1709688 | center_ += corner; |
171 | 336696 | center_ /= corners.size(); | |
172 | 336696 | } | |
173 | |||
174 | //! The center of the sub control volume face | ||
175 | const GlobalPosition& center() const | ||
176 | ✗ | { return center_; } | |
177 | |||
178 | //! The integration point for flux evaluations in global coordinates | ||
179 | 8874360 | const GlobalPosition& ipGlobal() const | |
180 |
1/2✓ Branch 1 taken 4452480 times.
✗ Branch 2 not taken.
|
10587410 | { return center_; } |
181 | |||
182 | //! The area of the sub control volume face | ||
183 | 74290634 | Scalar area() const | |
184 | 74290634 | { return area_; } | |
185 | |||
186 | //! returns true if the sub control volume face is on the boundary | ||
187 | 223668992 | bool boundary() const | |
188 |
4/4✓ Branch 0 taken 36288792 times.
✓ Branch 1 taken 2968120 times.
✓ Branch 2 taken 36288792 times.
✓ Branch 3 taken 2968120 times.
|
78513824 | { return boundary_; } |
189 | |||
190 | //! returns the unit normal vector pointing outwards | ||
191 | const GlobalPosition& unitOuterNormal() const | ||
192 |
1/2✓ Branch 0 taken 72577584 times.
✗ Branch 1 not taken.
|
145155168 | { return unitOuterNormal_; } |
193 | |||
194 | //! The global index of this sub control volume face | ||
195 | 81451944 | GridIndexType index() const | |
196 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 72577584 times.
|
81451944 | { return scvfIndex_; } |
197 | |||
198 | //! Return if this is a fracture scvf | ||
199 | 17748720 | bool isOnFracture() const | |
200 |
4/4✓ Branch 0 taken 8571720 times.
✓ Branch 1 taken 302640 times.
✓ Branch 2 taken 8571720 times.
✓ Branch 3 taken 302640 times.
|
17748720 | { return isFractureScvf_; } |
201 | |||
202 | //! The element-local facet index for which a fracture scv was created | ||
203 | 605280 | LocalIndexType facetIndexInElement() const | |
204 | 605280 | { assert(isFractureScvf_); return facetIdx_; } | |
205 | |||
206 | //! The local edge index inside the intersection | ||
207 | LocalIndexType indexInIntersection() const | ||
208 | { assert(isFractureScvf_); return indexInIntersection_; } | ||
209 | |||
210 | //! Returns the boundary flag | ||
211 | typename BoundaryFlag::value_type boundaryFlag() const | ||
212 | { return boundaryFlag_.get(); } | ||
213 | |||
214 | //! index of the inside sub control volume | ||
215 | 187380200 | LocalIndexType insideScvIdx() const | |
216 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 72577584 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 72577584 times.
|
187380200 | { return scvIndices_[0]; } |
217 | |||
218 | //! Index of the i-th outside sub control volume or boundary scv index. | ||
219 | // Results in undefined behaviour if i >= numOutsideScvs() | ||
220 | 181443960 | LocalIndexType outsideScvIdx(int i = 0) const | |
221 | { | ||
222 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 72577584 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 72577584 times.
|
145155168 | assert(!boundary()); |
223 |
3/4✓ Branch 0 taken 34864405 times.
✓ Branch 1 taken 37713179 times.
✓ Branch 2 taken 72577584 times.
✗ Branch 3 not taken.
|
181443960 | return scvIndices_[1]; |
224 | } | ||
225 | |||
226 | //! The number of scvs on the outside of this face | ||
227 | std::size_t numOutsideScvs() const | ||
228 | { | ||
229 | return static_cast<std::size_t>(!boundary()); | ||
230 | } | ||
231 | |||
232 | private: | ||
233 | GlobalPosition center_; | ||
234 | GlobalPosition unitOuterNormal_; | ||
235 | Scalar area_; | ||
236 | GridIndexType scvfIndex_; | ||
237 | std::vector<LocalIndexType> scvIndices_; | ||
238 | bool boundary_; | ||
239 | bool isFractureScvf_; | ||
240 | BoundaryFlag boundaryFlag_; | ||
241 | LocalIndexType facetIdx_; | ||
242 | LocalIndexType indexInIntersection_; | ||
243 | }; | ||
244 | |||
245 | } // end namespace Dumux | ||
246 | |||
247 | #endif | ||
248 |