version 3.8
intersectspointgeometry.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//
12#ifndef DUMUX_GEOMETRY_INTERSECTS_POINT_GEOMETRY_HH
13#define DUMUX_GEOMETRY_INTERSECTS_POINT_GEOMETRY_HH
14
15#include <cmath>
16#include <dune/common/exceptions.hh>
17#include <dune/common/fvector.hh>
18
20
21namespace Dumux {
22
27template <class ctype, int dimworld, class Geometry, typename std::enable_if_t<(Geometry::mydimension == 3), int> = 0>
28bool intersectsPointGeometry(const Dune::FieldVector<ctype, dimworld>& point, const Geometry& g)
29{
30 // get the g type
31 const auto type = g.type();
32
33 // if it's a tetrahedron we can check directly
34 if (type.isTetrahedron())
35 return intersectsPointSimplex(point, g.corner(0), g.corner(1), g.corner(2), g.corner(3));
36
37 // split hexahedrons into five tetrahedrons
38 else if (type.isHexahedron())
39 {
40 if (intersectsPointSimplex(point, g.corner(0), g.corner(1), g.corner(3), g.corner(5))) return true;
41 if (intersectsPointSimplex(point, g.corner(0), g.corner(5), g.corner(6), g.corner(4))) return true;
42 if (intersectsPointSimplex(point, g.corner(5), g.corner(3), g.corner(6), g.corner(7))) return true;
43 if (intersectsPointSimplex(point, g.corner(0), g.corner(3), g.corner(2), g.corner(6))) return true;
44 if (intersectsPointSimplex(point, g.corner(5), g.corner(3), g.corner(0), g.corner(6))) return true;
45 return false;
46 }
47
48 // split pyramids into two tetrahedrons
49 else if (type.isPyramid())
50 {
51 if (intersectsPointSimplex(point, g.corner(0), g.corner(1), g.corner(2), g.corner(4))) return true;
52 if (intersectsPointSimplex(point, g.corner(1), g.corner(3), g.corner(2), g.corner(4))) return true;
53 return false;
54 }
55
56 // split prisms into three tetrahedrons
57 else if (type.isPrism())
58 {
59 if (intersectsPointSimplex(point, g.corner(0), g.corner(1), g.corner(2), g.corner(4))) return true;
60 if (intersectsPointSimplex(point, g.corner(3), g.corner(0), g.corner(2), g.corner(4))) return true;
61 if (intersectsPointSimplex(point, g.corner(2), g.corner(5), g.corner(3), g.corner(4))) return true;
62 return false;
63 }
64
65 else
66 DUNE_THROW(Dune::NotImplemented,
67 "Intersection for point and geometry type "
68 << type << " in " << dimworld << "-dimensional world.");
69}
70
75template <class ctype, int dimworld, class Geometry, typename std::enable_if_t<(Geometry::mydimension == 2), int> = 0>
76bool intersectsPointGeometry(const Dune::FieldVector<ctype, dimworld>& point, const Geometry& g)
77{
78 // get the g type
79 const auto type = g.type();
80
81 // if it's a triangle we can check directly
82 if (type.isTriangle())
83 return intersectsPointSimplex(point, g.corner(0), g.corner(1), g.corner(2));
84
85 // split quadrilaterals into two triangles
86 else if (type.isQuadrilateral())
87 {
88 if (intersectsPointSimplex(point, g.corner(0), g.corner(1), g.corner(3))) return true;
89 if (intersectsPointSimplex(point, g.corner(0), g.corner(3), g.corner(2))) return true;
90 return false;
91 }
92
93 else
94 DUNE_THROW(Dune::NotImplemented,
95 "Intersection for point and geometry type "
96 << type << " in " << dimworld << "-dimensional world.");
97}
98
103template <class ctype, int dimworld, class Geometry, typename std::enable_if_t<(Geometry::mydimension == 1), int> = 0>
104bool intersectsPointGeometry(const Dune::FieldVector<ctype, dimworld>& point, const Geometry& g)
105{
106 return intersectsPointSimplex(point, g.corner(0), g.corner(1));
107}
108
109} // end namespace Dumux
110
111#endif
bool intersectsPointGeometry(const Dune::FieldVector< ctype, dimworld > &point, const Geometry &g)
Find out whether a point is inside a three-dimensional geometry.
Definition: intersectspointgeometry.hh:28
bool intersectsPointSimplex(const Dune::FieldVector< ctype, dimworld > &point, const Dune::FieldVector< ctype, dimworld > &p0, const Dune::FieldVector< ctype, dimworld > &p1, const Dune::FieldVector< ctype, dimworld > &p2, const Dune::FieldVector< ctype, dimworld > &p3)
Find out whether a point is inside the tetrahedron (p0, p1, p2, p3) (dimworld is 3)
Definition: intersectspointsimplex.hh:26
Detect if a point intersects a simplex (including boundary)
Definition: adapt.hh:17