3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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 *****************************************************************************/
24#ifndef DUMUX_GEOMETRY_INTERSECTS_POINT_GEOMETRY_HH
25#define DUMUX_GEOMETRY_INTERSECTS_POINT_GEOMETRY_HH
26
27#include <cmath>
28#include <dune/common/exceptions.hh>
29#include <dune/common/fvector.hh>
30
32
33namespace Dumux {
34
39template <class ctype, int dimworld, class Geometry, typename std::enable_if_t<(Geometry::mydimension == 3), int> = 0>
40bool intersectsPointGeometry(const Dune::FieldVector<ctype, dimworld>& point, const Geometry& g)
41{
42 // get the g type
43 const auto type = g.type();
44
45 // if it's a tetrahedron we can check directly
46 if (type.isTetrahedron())
47 return intersectsPointSimplex(point, g.corner(0), g.corner(1), g.corner(2), g.corner(3));
48
49 // split hexahedrons into five tetrahedrons
50 else if (type.isHexahedron())
51 {
52 if (intersectsPointSimplex(point, g.corner(0), g.corner(1), g.corner(3), g.corner(5))) return true;
53 if (intersectsPointSimplex(point, g.corner(0), g.corner(5), g.corner(6), g.corner(4))) return true;
54 if (intersectsPointSimplex(point, g.corner(5), g.corner(3), g.corner(6), g.corner(7))) return true;
55 if (intersectsPointSimplex(point, g.corner(0), g.corner(3), g.corner(2), g.corner(6))) return true;
56 if (intersectsPointSimplex(point, g.corner(5), g.corner(3), g.corner(0), g.corner(6))) return true;
57 return false;
58 }
59
60 // split pyramids into two tetrahedrons
61 else if (type.isPyramid())
62 {
63 if (intersectsPointSimplex(point, g.corner(0), g.corner(1), g.corner(2), g.corner(4))) return true;
64 if (intersectsPointSimplex(point, g.corner(1), g.corner(3), g.corner(2), g.corner(4))) return true;
65 return false;
66 }
67
68 // split prisms into three tetrahedrons
69 else if (type.isPrism())
70 {
71 if (intersectsPointSimplex(point, g.corner(0), g.corner(1), g.corner(2), g.corner(4))) return true;
72 if (intersectsPointSimplex(point, g.corner(3), g.corner(0), g.corner(2), g.corner(4))) return true;
73 if (intersectsPointSimplex(point, g.corner(2), g.corner(5), g.corner(3), g.corner(4))) return true;
74 return false;
75 }
76
77 else
78 DUNE_THROW(Dune::NotImplemented,
79 "Intersection for point and geometry type "
80 << type << " in " << dimworld << "-dimensional world.");
81}
82
87template <class ctype, int dimworld, class Geometry, typename std::enable_if_t<(Geometry::mydimension == 2), int> = 0>
88bool intersectsPointGeometry(const Dune::FieldVector<ctype, dimworld>& point, const Geometry& g)
89{
90 // get the g type
91 const auto type = g.type();
92
93 // if it's a triangle we can check directly
94 if (type.isTriangle())
95 return intersectsPointSimplex(point, g.corner(0), g.corner(1), g.corner(2));
96
97 // split quadrilaterals into two triangles
98 else if (type.isQuadrilateral())
99 {
100 if (intersectsPointSimplex(point, g.corner(0), g.corner(1), g.corner(3))) return true;
101 if (intersectsPointSimplex(point, g.corner(0), g.corner(3), g.corner(2))) return true;
102 return false;
103 }
104
105 else
106 DUNE_THROW(Dune::NotImplemented,
107 "Intersection for point and geometry type "
108 << type << " in " << dimworld << "-dimensional world.");
109}
110
115template <class ctype, int dimworld, class Geometry, typename std::enable_if_t<(Geometry::mydimension == 1), int> = 0>
116bool intersectsPointGeometry(const Dune::FieldVector<ctype, dimworld>& point, const Geometry& g)
117{
118 return intersectsPointSimplex(point, g.corner(0), g.corner(1));
119}
120
121} // end namespace Dumux
122
123#endif
Detect if a point intersects a simplex (including boundary)
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:40
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:38
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29