3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
geometryhelper.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_GEOMETRY_HELPER_HH
26#define DUMUX_POROUSMEDIUMFLOW_BOXDFM_GEOMETRY_HELPER_HH
27
29
30namespace Dumux
31{
32
34template<class GridView, int dim, class ScvType, class ScvfType>
36
38template <class GridView, class ScvType, class ScvfType>
39class BoxDfmGeometryHelper<GridView, 2, ScvType, ScvfType> : public BoxGeometryHelper<GridView, 2, ScvType, ScvfType>
40{
42
43 using Intersection = typename GridView::Intersection;
44 using ScvfCornerStorage = typename ScvfType::Traits::CornerStorage;
45
46 static constexpr auto dim = GridView::dimension;
47 using Scalar = typename GridView::ctype;
48 using ReferenceElements = typename Dune::ReferenceElements<Scalar, dim>;
49
50public:
51
53 using ParentType::ParentType;
54
57 ScvfCornerStorage getFractureScvfCorners(const Intersection& is,
58 const typename Intersection::Geometry& isGeom,
59 unsigned int idxOnIntersection = 0) const
60 { return ScvfCornerStorage({isGeom.center()}); }
61
62
65 typename ScvfType::Traits::GlobalPosition
66 fractureNormal(const ScvfCornerStorage& p,
67 const Intersection& is,
68 unsigned int edgeIndexInIntersection = 0) const
69 {
70 const auto referenceElement = ReferenceElements::general(this->elementGeometry_.type());
71 const auto vIdxLocal0 = referenceElement.subEntity(is.indexInInside(), 1, 0, dim);
72 const auto vIdxLocal1 = referenceElement.subEntity(is.indexInInside(), 1, 1, dim);
73 auto n = this->elementGeometry_.corner(vIdxLocal1) - this->elementGeometry_.corner(vIdxLocal0);
74 n /= n.two_norm();
75 return n;
76 }
77};
78
80template <class GridView, class ScvType, class ScvfType>
81class BoxDfmGeometryHelper<GridView, 3, ScvType, ScvfType> : public BoxGeometryHelper<GridView, 3, ScvType, ScvfType>
82{
84
85 using Intersection = typename GridView::Intersection;
86 using ScvfCornerStorage = typename ScvfType::Traits::CornerStorage;
87
88 static constexpr auto dim = GridView::dimension;
89 static constexpr auto dimWorld = GridView::dimensionworld;
90 using Scalar = typename GridView::ctype;
91 using ReferenceElements = typename Dune::ReferenceElements<Scalar, dim>;
92 using FaceReferenceElements = typename Dune::ReferenceElements<Scalar, dim-1>;
93
94public:
95
97 using ParentType::ParentType;
98
100 ScvfCornerStorage getFractureScvfCorners(const Intersection& is,
101 const typename Intersection::Geometry& isGeom,
102 unsigned int edgeIndexInIntersection) const
103 {
104 const auto referenceElement = ReferenceElements::general(this->elementGeometry_.type());
105 const auto faceRefElem = FaceReferenceElements::general(isGeom.type());
106
107 // create point vector for this geometry
108 typename ScvfType::Traits::GlobalPosition pi[9];
109
110 // the facet center
111 pi[0] = isGeom.center();
112
113 // facet edge midpoints
114 const auto idxInInside = is.indexInInside();
115 for (int i = 0; i < faceRefElem.size(1); ++i)
116 {
117 const auto edgeIdxLocal = referenceElement.subEntity(idxInInside, 1, i, dim-1);
118 pi[i+1] = this->p_[edgeIdxLocal+this->corners_+1];
119 }
120
121 // proceed according to number of corners
122 const auto corners = isGeom.corners();
123 switch (corners)
124 {
125 case 3: // triangle
126 {
128 static const std::uint8_t fo = 1;
129 static const std::uint8_t map[3][2] =
130 {
131 {fo+0, 0},
132 {0, fo+1},
133 {fo+2, 0}
134 };
135
136 return ScvfCornerStorage{ {pi[map[edgeIndexInIntersection][0]],
137 pi[map[edgeIndexInIntersection][1]]} };
138 }
139 case 4: // quadrilateral
140 {
142 static const std::uint8_t fo = 1;
143 static const std::uint8_t map[4][2] =
144 {
145 {0, fo+0},
146 {fo+1, 0},
147 {fo+2, 0},
148 {0, fo+3}
149 };
150
151 return ScvfCornerStorage{ {pi[map[edgeIndexInIntersection][0]],
152 pi[map[edgeIndexInIntersection][1]]} };
153 }
154 default:
155 DUNE_THROW(Dune::NotImplemented, "Box fracture scvf geometries for dim=" << dim
156 << " dimWorld=" << dimWorld
157 << " corners=" << corners);
158 }
159 }
160
162 typename ScvfType::Traits::GlobalPosition
163 fractureNormal(const ScvfCornerStorage& p,
164 const Intersection& is,
165 unsigned int edgeIndexInIntersection) const
166 {
167 const auto referenceElement = ReferenceElements::general(this->elementGeometry_.type());
168
169 // first get the intersection corners (maximum "4" is for quadrilateral face)
170 typename ScvfType::Traits::GlobalPosition c[4];
171
172 const auto corners = referenceElement.size(is.indexInInside(), 1, dim);
173 for (int i = 0; i < corners; ++i)
174 {
175 const auto vIdxLocal = referenceElement.subEntity(is.indexInInside(), 1, i, dim);
176 c[i] = this->elementGeometry_.corner(vIdxLocal);
177 }
178
179 // compute edge vector depending on number of corners
180 const auto gridEdge = [&] ()
181 {
182 // triangles
183 if (corners == 3)
184 {
185 if (edgeIndexInIntersection == 0) return c[1]-c[0];
186 else if (edgeIndexInIntersection == 1) return c[2]-c[0];
187 else if (edgeIndexInIntersection == 2) return c[2]-c[1];
188 else DUNE_THROW(Dune::InvalidStateException, "Invalid edge index");
189 }
190 else if (corners == 4)
191 {
192 if (edgeIndexInIntersection == 0) return c[2]-c[0];
193 else if (edgeIndexInIntersection == 1) return c[3]-c[1];
194 else if (edgeIndexInIntersection == 2) return c[1]-c[0];
195 else if (edgeIndexInIntersection == 3) return c[3]-c[2];
196 else DUNE_THROW(Dune::InvalidStateException, "Invalid edge index");
197 }
198 else
199 DUNE_THROW(Dune::InvalidStateException, "Invalid face geometry");
200 } ();
201
202 // compute lower edge of the scvf
203 assert(p.size() == 2);
204 const auto scvfEdge = p[1]-p[0];
205
206 // compute scvf normal via 2 cross products
207 const auto faceN = crossProduct(gridEdge, scvfEdge);
208 auto n = crossProduct(scvfEdge, faceN);
209 n /= n.two_norm();
210 return n;
211 }
212};
213
214} // end namespace Dumux
215
216#endif
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
Dune::FieldVector< Scalar, 3 > crossProduct(const Dune::FieldVector< Scalar, 3 > &vec1, const Dune::FieldVector< Scalar, 3 > &vec2)
Cross product of two vectors in three-dimensional Euclidean space.
Definition: math.hh:631
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Create sub control volumes and sub control volume face geometries.
Definition: boxgeometryhelper.hh:37
A class to create sub control volume and sub control volume face geometries per element.
Definition: boxgeometryhelper.hh:120
A class to create sub control volume and sub control volume face geometries per element.
Definition: boxgeometryhelper.hh:340
Create sub control volumes and sub control volume face geometries.
Definition: geometryhelper.hh:35
ScvfType::Traits::GlobalPosition fractureNormal(const ScvfCornerStorage &p, const Intersection &is, unsigned int edgeIndexInIntersection=0) const
Definition: geometryhelper.hh:66
ScvfCornerStorage getFractureScvfCorners(const Intersection &is, const typename Intersection::Geometry &isGeom, unsigned int idxOnIntersection=0) const
Definition: geometryhelper.hh:57
ScvfCornerStorage getFractureScvfCorners(const Intersection &is, const typename Intersection::Geometry &isGeom, unsigned int edgeIndexInIntersection) const
Create the sub control volume face geometries on an intersection marked as fracture.
Definition: geometryhelper.hh:100
ScvfType::Traits::GlobalPosition fractureNormal(const ScvfCornerStorage &p, const Intersection &is, unsigned int edgeIndexInIntersection) const
get fracture scvf normal vector
Definition: geometryhelper.hh:163