version 3.8
porousmediumflow/boxdfm/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// 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_GEOMETRY_HELPER_HH
14#define DUMUX_POROUSMEDIUMFLOW_BOXDFM_GEOMETRY_HELPER_HH
15
16#include <dune/common/fvector.hh>
17#include <dune/common/reservedvector.hh>
18#include <dune/geometry/multilineargeometry.hh>
19
21
22namespace Dumux {
23
24template <class ct>
25struct BoxDfmMLGeometryTraits : public Dune::MultiLinearGeometryTraits<ct>
26{
27 // we use static vectors to store the corners as we know
28 // the number of corners in advance (2^(mydim) corners (1<<(mydim))
29 // However, on fracture scvs the number might be smaller (use ReservedVector)
30 template< int mydim, int cdim >
32 {
33 using Type = Dune::ReservedVector< Dune::FieldVector< ct, cdim >, (1<<(mydim)) >;
34 };
35
36 // we know all scvfs will have the same geometry type
37 template< int mydim >
39 {
40 static const bool v = true;
41 static const unsigned int topologyId = Dune::GeometryTypes::cube(mydim).id();
42 };
43};
44
46template<class GridView, int dim, class ScvType, class ScvfType>
48
50template <class GridView, class ScvType, class ScvfType>
51class BoxDfmGeometryHelper<GridView, 2, ScvType, ScvfType> : public BoxGeometryHelper<GridView, 2, ScvType, ScvfType>
52{
54
55 using Intersection = typename GridView::Intersection;
56 using ScvfCornerStorage = typename ScvfType::Traits::CornerStorage;
57
58 static constexpr auto dim = GridView::dimension;
59 using Scalar = typename GridView::ctype;
60
61public:
62
64 using ParentType::ParentType;
65
67 ScvfCornerStorage getFractureScvfCorners(unsigned int localFacetIndex,
68 unsigned int) const
69 {
70 const auto& geo = this->elementGeometry();
71 const auto ref = referenceElement(geo);
72 return ScvfCornerStorage({ geo.global(ref.position(localFacetIndex, 1)) });
73 }
74
77 typename ScvfType::Traits::GlobalPosition
78 fractureNormal(const ScvfCornerStorage& p,
79 const Intersection& is,
80 unsigned int edgeIndexInIntersection) const
81 {
82 const auto& geo = this->elementGeometry();
83 const auto ref = referenceElement(geo);
84 const auto v0 = ref.subEntity(is.indexInInside(), 1, 0, dim);
85 const auto v1 = ref.subEntity(is.indexInInside(), 1, 1, dim);
86 auto normal = geo.corner(v1) - geo.corner(v0);
87 normal /= normal.two_norm();
88 return normal;
89 }
90};
91
93template <class GridView, class ScvType, class ScvfType>
94class BoxDfmGeometryHelper<GridView, 3, ScvType, ScvfType> : public BoxGeometryHelper<GridView, 3, ScvType, ScvfType>
95{
97
98 using Intersection = typename GridView::Intersection;
99 using ScvfCornerStorage = typename ScvfType::Traits::CornerStorage;
100
101 static constexpr auto dim = GridView::dimension;
102 static constexpr auto dimWorld = GridView::dimensionworld;
103 using Scalar = typename GridView::ctype;
104
105public:
106
108 using ParentType::ParentType;
109
111 ScvfCornerStorage getFractureScvfCorners(unsigned int localFacetIndex,
112 unsigned int indexInFacet) const
113 {
114 constexpr int facetCodim = 1;
115
116 // we have to use the corresponding facet geometry as the intersection geometry
117 // might be rotated or flipped. This makes sure that the corners (dof location)
118 // and corresponding scvfs are sorted in the same way
119 using Dune::referenceElement;
120 const auto& geo = this->elementGeometry();
121 const auto type = referenceElement(geo).type(localFacetIndex, facetCodim);
122 if (type == Dune::GeometryTypes::triangle)
123 {
125 return Detail::Box::subEntityKeyToCornerStorage<ScvfCornerStorage>(geo, localFacetIndex, facetCodim, Corners::keys[indexInFacet]);
126 }
127 else if (type == Dune::GeometryTypes::quadrilateral)
128 {
130 return Detail::Box::subEntityKeyToCornerStorage<ScvfCornerStorage>(geo, localFacetIndex, facetCodim, Corners::keys[indexInFacet]);
131 }
132 else
133 DUNE_THROW(Dune::NotImplemented, "Box fracture scvf geometries for dim=" << dim
134 << " dimWorld=" << dimWorld
135 << " type=" << type);
136 }
137
139 typename ScvfType::Traits::GlobalPosition
140 fractureNormal(const ScvfCornerStorage& scvfCorners,
141 const Intersection& is,
142 unsigned int edgeIndexInIntersection) const
143 {
144 const auto& geo = this->elementGeometry();
145 const auto refElement = referenceElement(geo);
146
147 // first get the intersection corners (maximum "4" is for quadrilateral face)
148 typename ScvfType::Traits::GlobalPosition c[4];
149
150 const auto corners = refElement.size(is.indexInInside(), 1, dim);
151 for (int i = 0; i < corners; ++i)
152 {
153 const auto vIdxLocal = refElement.subEntity(is.indexInInside(), 1, i, dim);
154 c[i] = geo.corner(vIdxLocal);
155 }
156
157 // compute edge vector depending on number of corners
158 const auto gridEdge = [&] ()
159 {
160 // triangles
161 if (corners == 3)
162 {
163 if (edgeIndexInIntersection == 0) return c[1]-c[0];
164 else if (edgeIndexInIntersection == 1) return c[2]-c[0];
165 else if (edgeIndexInIntersection == 2) return c[2]-c[1];
166 else DUNE_THROW(Dune::InvalidStateException, "Invalid edge index");
167 }
168 else if (corners == 4)
169 {
170 if (edgeIndexInIntersection == 0) return c[2]-c[0];
171 else if (edgeIndexInIntersection == 1) return c[3]-c[1];
172 else if (edgeIndexInIntersection == 2) return c[1]-c[0];
173 else if (edgeIndexInIntersection == 3) return c[3]-c[2];
174 else DUNE_THROW(Dune::InvalidStateException, "Invalid edge index");
175 }
176 else
177 DUNE_THROW(Dune::InvalidStateException, "Invalid face geometry");
178 } ();
179
180 // compute lower edge of the scvf
181 assert(scvfCorners.size() == 2);
182 const auto scvfEdge = scvfCorners[1]-scvfCorners[0];
183
184 // compute scvf normal via 2 cross products
185 const auto faceN = crossProduct(gridEdge, scvfEdge);
186 auto n = crossProduct(scvfEdge, faceN);
187 n /= n.two_norm();
188 return n;
189 }
190};
191
192} // end namespace Dumux
193
194#endif
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
ScvfCornerStorage getFractureScvfCorners(unsigned int localFacetIndex, unsigned int) const
Get the corners of the (d-1)-dimensional fracture scvf.
Definition: porousmediumflow/boxdfm/geometryhelper.hh:67
ScvfType::Traits::GlobalPosition fractureNormal(const ScvfCornerStorage &p, const Intersection &is, unsigned int edgeIndexInIntersection) const
Definition: porousmediumflow/boxdfm/geometryhelper.hh:78
ScvfCornerStorage getFractureScvfCorners(unsigned int localFacetIndex, unsigned int indexInFacet) const
Create the sub control volume face geometries on an intersection marked as fracture.
Definition: porousmediumflow/boxdfm/geometryhelper.hh:111
ScvfType::Traits::GlobalPosition fractureNormal(const ScvfCornerStorage &scvfCorners, const Intersection &is, unsigned int edgeIndexInIntersection) const
get fracture scvf normal vector
Definition: porousmediumflow/boxdfm/geometryhelper.hh:140
Create sub control volumes and sub control volume face geometries.
Definition: porousmediumflow/boxdfm/geometryhelper.hh:47
A class to create sub control volume and sub control volume face geometries per element.
Definition: boxgeometryhelper.hh:339
A class to create sub control volume and sub control volume face geometries per element.
Definition: boxgeometryhelper.hh:475
Create sub control volumes and sub control volume face geometries.
Definition: boxgeometryhelper.hh:257
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:642
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.hh:26
Definition: adapt.hh:17
Definition: porousmediumflow/boxdfm/geometryhelper.hh:32
Dune::ReservedVector< Dune::FieldVector< ct, cdim >,(1<<(mydim)) > Type
Definition: porousmediumflow/boxdfm/geometryhelper.hh:33
Definition: porousmediumflow/boxdfm/geometryhelper.hh:39
static const bool v
Definition: porousmediumflow/boxdfm/geometryhelper.hh:40
static const unsigned int topologyId
Definition: porousmediumflow/boxdfm/geometryhelper.hh:41
Definition: porousmediumflow/boxdfm/geometryhelper.hh:26