25#ifndef DUMUX_POROUSMEDIUMFLOW_BOXDFM_GEOMETRY_HELPER_HH
26#define DUMUX_POROUSMEDIUMFLOW_BOXDFM_GEOMETRY_HELPER_HH
34template<
class Gr
idView,
int dim,
class ScvType,
class ScvfType>
38template <
class Gr
idView,
class ScvType,
class ScvfType>
43 using Intersection =
typename GridView::Intersection;
44 using ScvfCornerStorage =
typename ScvfType::Traits::CornerStorage;
46 static constexpr auto dim = GridView::dimension;
47 using Scalar =
typename GridView::ctype;
48 using ReferenceElements =
typename Dune::ReferenceElements<Scalar, dim>;
53 using ParentType::ParentType;
58 const typename Intersection::Geometry& isGeom,
59 unsigned int idxOnIntersection = 0)
const
60 {
return ScvfCornerStorage({isGeom.center()}); }
65 typename ScvfType::Traits::GlobalPosition
67 const Intersection& is,
68 unsigned int edgeIndexInIntersection = 0)
const
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);
80template <
class Gr
idView,
class ScvType,
class ScvfType>
85 using Intersection =
typename GridView::Intersection;
86 using ScvfCornerStorage =
typename ScvfType::Traits::CornerStorage;
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>;
97 using ParentType::ParentType;
101 const typename Intersection::Geometry& isGeom,
102 unsigned int edgeIndexInIntersection)
const
104 const auto referenceElement = ReferenceElements::general(this->elementGeometry_.type());
105 const auto faceRefElem = FaceReferenceElements::general(isGeom.type());
108 typename ScvfType::Traits::GlobalPosition pi[9];
111 pi[0] = isGeom.center();
114 const auto idxInInside = is.indexInInside();
115 for (
int i = 0; i < faceRefElem.size(1); ++i)
117 const auto edgeIdxLocal = referenceElement.subEntity(idxInInside, 1, i, dim-1);
118 pi[i+1] = this->p_[edgeIdxLocal+this->corners_+1];
122 const auto corners = isGeom.corners();
128 static const std::uint8_t fo = 1;
129 static const std::uint8_t map[3][2] =
136 return ScvfCornerStorage{ {pi[map[edgeIndexInIntersection][0]],
137 pi[map[edgeIndexInIntersection][1]]} };
142 static const std::uint8_t fo = 1;
143 static const std::uint8_t map[4][2] =
151 return ScvfCornerStorage{ {pi[map[edgeIndexInIntersection][0]],
152 pi[map[edgeIndexInIntersection][1]]} };
155 DUNE_THROW(Dune::NotImplemented,
"Box fracture scvf geometries for dim=" << dim
156 <<
" dimWorld=" << dimWorld
157 <<
" corners=" << corners);
162 typename ScvfType::Traits::GlobalPosition
164 const Intersection& is,
165 unsigned int edgeIndexInIntersection)
const
167 const auto referenceElement = ReferenceElements::general(this->elementGeometry_.type());
170 typename ScvfType::Traits::GlobalPosition c[4];
172 const auto corners = referenceElement.size(is.indexInInside(), 1, dim);
173 for (
int i = 0; i < corners; ++i)
175 const auto vIdxLocal = referenceElement.subEntity(is.indexInInside(), 1, i, dim);
176 c[i] = this->elementGeometry_.corner(vIdxLocal);
180 const auto gridEdge = [&] ()
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");
190 else if (corners == 4)
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");
199 DUNE_THROW(Dune::InvalidStateException,
"Invalid face geometry");
203 assert(p.size() == 2);
204 const auto scvfEdge = p[1]-p[0];
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