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-FileCopyrightInfo: 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 Geometry | ||
10 | * \brief Detect if a point intersects a geometry | ||
11 | */ | ||
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 | |||
19 | #include <dumux/geometry/intersectspointsimplex.hh> | ||
20 | |||
21 | namespace Dumux { | ||
22 | |||
23 | /*! | ||
24 | * \ingroup Geometry | ||
25 | * \brief Find out whether a point is inside a three-dimensional geometry | ||
26 | */ | ||
27 | template <class ctype, int dimworld, class Geometry, typename std::enable_if_t<(Geometry::mydimension == 3), int> = 0> | ||
28 | 546200 | bool intersectsPointGeometry(const Dune::FieldVector<ctype, dimworld>& point, const Geometry& g) | |
29 | { | ||
30 | // get the g type | ||
31 | 546200 | const auto type = g.type(); | |
32 | |||
33 | // if it's a tetrahedron we can check directly | ||
34 |
2/2✓ Branch 0 taken 24751 times.
✓ Branch 1 taken 5100 times.
|
546200 | if (type.isTetrahedron()) |
35 | 24878 | return intersectsPointSimplex(point, g.corner(0), g.corner(1), g.corner(2), g.corner(3)); | |
36 | |||
37 | // split hexahedrons into five tetrahedrons | ||
38 |
2/2✓ Branch 0 taken 27131 times.
✓ Branch 1 taken 2000 times.
|
521322 | else if (type.isHexahedron()) |
39 | { | ||
40 |
4/4✓ Branch 1 taken 267278 times.
✓ Branch 2 taken 39170 times.
✓ Branch 5 taken 29734 times.
✓ Branch 6 taken 7001 times.
|
2377518 | if (intersectsPointSimplex(point, g.corner(0), g.corner(1), g.corner(3), g.corner(5))) return true; |
41 |
4/4✓ Branch 1 taken 228304 times.
✓ Branch 2 taken 38974 times.
✓ Branch 5 taken 24086 times.
✓ Branch 6 taken 5648 times.
|
2064069 | if (intersectsPointSimplex(point, g.corner(0), g.corner(5), g.corner(6), g.corner(4))) return true; |
42 |
4/4✓ Branch 1 taken 188816 times.
✓ Branch 2 taken 39488 times.
✓ Branch 5 taken 19147 times.
✓ Branch 6 taken 4939 times.
|
1773265 | if (intersectsPointSimplex(point, g.corner(5), g.corner(3), g.corner(6), g.corner(7))) return true; |
43 |
4/4✓ Branch 1 taken 150384 times.
✓ Branch 2 taken 38432 times.
✓ Branch 5 taken 14207 times.
✓ Branch 6 taken 4940 times.
|
1490216 | if (intersectsPointSimplex(point, g.corner(0), g.corner(3), g.corner(2), g.corner(6))) return true; |
44 |
4/4✓ Branch 1 taken 90967 times.
✓ Branch 2 taken 59417 times.
✓ Branch 5 taken 6899 times.
✓ Branch 6 taken 7308 times.
|
1221535 | if (intersectsPointSimplex(point, g.corner(5), g.corner(3), g.corner(0), g.corner(6))) return true; |
45 | 180462 | return false; | |
46 | } | ||
47 | |||
48 | // split pyramids into two tetrahedrons | ||
49 |
2/2✓ Branch 0 taken 1052 times.
✓ Branch 1 taken 960 times.
|
2024 | else if (type.isPyramid()) |
50 | { | ||
51 |
2/2✓ Branch 5 taken 400 times.
✓ Branch 6 taken 640 times.
|
1040 | if (intersectsPointSimplex(point, g.corner(0), g.corner(1), g.corner(2), g.corner(4))) return true; |
52 |
2/2✓ Branch 5 taken 240 times.
✓ Branch 6 taken 160 times.
|
400 | if (intersectsPointSimplex(point, g.corner(1), g.corner(3), g.corner(2), g.corner(4))) return true; |
53 | 240 | return false; | |
54 | } | ||
55 | |||
56 | // split prisms into three tetrahedrons | ||
57 |
2/2✓ Branch 0 taken 12 times.
✓ Branch 1 taken 960 times.
|
984 | else if (type.isPrism()) |
58 | { | ||
59 |
2/2✓ Branch 5 taken 486 times.
✓ Branch 6 taken 486 times.
|
984 | if (intersectsPointSimplex(point, g.corner(0), g.corner(1), g.corner(2), g.corner(4))) return true; |
60 |
2/2✓ Branch 5 taken 320 times.
✓ Branch 6 taken 166 times.
|
492 | if (intersectsPointSimplex(point, g.corner(3), g.corner(0), g.corner(2), g.corner(4))) return true; |
61 |
2/2✓ Branch 5 taken 160 times.
✓ Branch 6 taken 160 times.
|
320 | if (intersectsPointSimplex(point, g.corner(2), g.corner(5), g.corner(3), g.corner(4))) return true; |
62 | 160 | 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 | |||
71 | /*! | ||
72 | * \ingroup Geometry | ||
73 | * \brief Find out whether a point is inside a two-dimensional geometry | ||
74 | */ | ||
75 | template <class ctype, int dimworld, class Geometry, typename std::enable_if_t<(Geometry::mydimension == 2), int> = 0> | ||
76 | 1952430 | bool intersectsPointGeometry(const Dune::FieldVector<ctype, dimworld>& point, const Geometry& g) | |
77 | { | ||
78 | // get the g type | ||
79 | 1952430 | const auto type = g.type(); | |
80 | |||
81 | // if it's a triangle we can check directly | ||
82 |
2/2✓ Branch 0 taken 286225 times.
✓ Branch 1 taken 10164 times.
|
1952430 | if (type.isTriangle()) |
83 | 1293800 | return intersectsPointSimplex(point, g.corner(0), g.corner(1), g.corner(2)); | |
84 | |||
85 | // split quadrilaterals into two triangles | ||
86 |
2/2✓ Branch 0 taken 4 times.
✓ Branch 1 taken 10164 times.
|
1604941 | else if (type.isQuadrilateral()) |
87 | { | ||
88 |
6/6✓ Branch 0 taken 1148919 times.
✓ Branch 1 taken 436879 times.
✓ Branch 2 taken 1324 times.
✓ Branch 3 taken 77 times.
✓ Branch 4 taken 7162 times.
✓ Branch 5 taken 3036 times.
|
6383743 | if (intersectsPointSimplex(point, g.corner(0), g.corner(1), g.corner(3))) return true; |
89 |
6/6✓ Branch 0 taken 918704 times.
✓ Branch 1 taken 249241 times.
✓ Branch 2 taken 2136 times.
✓ Branch 3 taken 20 times.
✓ Branch 4 taken 5085 times.
✓ Branch 5 taken 2052 times.
|
4705101 | if (intersectsPointSimplex(point, g.corner(0), g.corner(3), g.corner(2))) return true; |
90 | 946262 | 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 | |||
99 | /*! | ||
100 | * \ingroup Geometry | ||
101 | * \brief Find out whether a point is inside a one-dimensional geometry | ||
102 | */ | ||
103 | template <class ctype, int dimworld, class Geometry, typename std::enable_if_t<(Geometry::mydimension == 1), int> = 0> | ||
104 | 261800 | bool intersectsPointGeometry(const Dune::FieldVector<ctype, dimworld>& point, const Geometry& g) | |
105 | { | ||
106 |
3/12✓ Branch 0 taken 20 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 20 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 20 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
525896 | return intersectsPointSimplex(point, g.corner(0), g.corner(1)); |
107 | } | ||
108 | |||
109 | } // end namespace Dumux | ||
110 | |||
111 | #endif | ||
112 |