preCICE v3.2.0
Loading...
Searching...
No Matches
geometry.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <Eigen/Core>
4#include <array>
6#include "utils/assertion.hpp"
7
10
18
31 const Eigen::Ref<const Eigen::Vector2d> &a,
32 const Eigen::Ref<const Eigen::Vector2d> &b,
33 const Eigen::Ref<const Eigen::Vector2d> &c,
34 const Eigen::Ref<const Eigen::Vector2d> &d,
35 Eigen::Ref<Eigen::Vector2d> &intersectionPoint);
36
50 const Eigen::Vector3d &pointOnPlane,
51 const Eigen::Vector3d &planeNormal,
52 const Eigen::Vector3d &firstPointSegment,
53 const Eigen::Vector3d &secondPointSegment,
54 Eigen::Vector3d &intersectionPoint);
55
67template <typename VECTORA_T, typename VECTORB_T, typename VECTORC_T>
68bool between(
69 const VECTORA_T &a,
70 const VECTORB_T &b,
71 const VECTORC_T &toCheck);
72
84template <typename Derived>
85bool collinear(
86 const Eigen::MatrixBase<Derived> &a,
87 const Eigen::MatrixBase<Derived> &b,
88 const Eigen::MatrixBase<Derived> &c);
89
101template <typename Derived>
102static bool parallel(
103 const Eigen::MatrixBase<Derived> &a,
104 const Eigen::MatrixBase<Derived> &b,
105 const Eigen::MatrixBase<Derived> &c,
106 const Eigen::MatrixBase<Derived> &d);
107
114double triangleArea(
115 const Eigen::VectorXd &a,
116 const Eigen::VectorXd &b,
117 const Eigen::VectorXd &c);
118
120double tetraVolume(
121 const Eigen::Vector3d &a,
122 const Eigen::Vector3d &b,
123 const Eigen::Vector3d &c,
124 const Eigen::Vector3d &d);
125
127Eigen::Vector2d projectVector(
128 const Eigen::Vector3d &vector,
129 const int indexDimensionToRemove);
130
139template <class Derived>
141 const Eigen::MatrixBase<Derived> &sidelengths,
142 const Eigen::MatrixBase<Derived> &center,
143 const Eigen::MatrixBase<Derived> &testPoint);
144
145// --------------------------------------------------------- HEADER DEFINITIONS
146
147template <typename VECTORA_T, typename VECTORB_T, typename VECTORC_T>
149 const VECTORA_T &a,
150 const VECTORB_T &b,
151 const VECTORC_T &toCheck)
152{
153 if (not collinear(a, b, toCheck)) {
154 return false;
155 }
156
157 if (a(0) != b(0)) {
158 return (math::greaterEquals(toCheck(0), a(0)) && math::greaterEquals(b(0), toCheck(0))) ||
159 (math::greaterEquals(a(0), toCheck(0)) && math::greaterEquals(toCheck(0), b(0)));
160 } else {
161 return (math::greaterEquals(toCheck(1), a(1)) && math::greaterEquals(b(1), toCheck(1))) ||
162 (math::greaterEquals(a(1), toCheck(1)) && math::greaterEquals(toCheck(1), b(1)));
163 }
164}
165
166template <typename Derived>
168 const Eigen::MatrixBase<Derived> &a,
169 const Eigen::MatrixBase<Derived> &b,
170 const Eigen::MatrixBase<Derived> &c)
171{
172 PRECICE_ASSERT(a.size() == b.size(), a.size(), b.size());
173 PRECICE_ASSERT(a.size() == c.size(), a.size(), c.size());
174 double triangleOutline = (b - a).norm() + (c - b).norm() + (a - c).norm();
175 if (math::equals(triangleArea(a, b, c) / triangleOutline, 0.0)) {
176 return true;
177 }
178 return false;
179}
180
181template <class Derived>
183 const Eigen::MatrixBase<Derived> &a,
184 const Eigen::MatrixBase<Derived> &b,
185 const Eigen::MatrixBase<Derived> &c,
186 const Eigen::MatrixBase<Derived> &d)
187{
188 if (math::equals(triangleArea(a, b, c), 0.0) and math::equals(triangleArea(a, b, d), 0.0))
189 return true;
190
191 return false;
192}
193
194template <class Derived>
196 const Eigen::MatrixBase<Derived> &sidelengths,
197 const Eigen::MatrixBase<Derived> &center,
198 const Eigen::MatrixBase<Derived> &testPoint)
199{
200 int dim = sidelengths.size();
201 Eigen::VectorXd toCenter = testPoint - center;
202 toCenter = toCenter.cwiseAbs();
203
204 double diff = 0.0;
205 bool touching = false;
206 for (int i = 0; i < dim; i++) {
207 diff = 0.5 * sidelengths(i) - toCenter(i);
208 if (math::greater(0.0, diff)) {
209 return NOT_CONTAINED;
210 }
211 if (math::equals(diff, 0.0)) {
212 touching = true;
213 }
214 }
215 if (touching) {
216 return TOUCHING;
217 }
218 return CONTAINED;
219}
220
225
227
228} // namespace precice::math::geometry
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
STL class.
Provides computational geometry operations.
Definition geometry.cpp:12
bool between(const VECTORA_T &a, const VECTORB_T &b, const VECTORC_T &toCheck)
Determines, if a point lies on the line segment defined by a, b.
Definition geometry.hpp:148
double triangleArea(const Eigen::VectorXd &a, const Eigen::VectorXd &b, const Eigen::VectorXd &c)
Computes the signed area of a triangle in 2D.
Definition geometry.cpp:93
double tetraVolume(const Eigen::Vector3d &a, const Eigen::Vector3d &b, const Eigen::Vector3d &c, const Eigen::Vector3d &d)
Computes the (unsigned) area of a triangle in 3D.
Definition geometry.cpp:116
Eigen::Vector2d projectVector(const Eigen::Vector3d &vector, int indexDimensionToRemove)
Projects a 3D vector to a 2D one by removing one dimension.
Definition geometry.cpp:125
int containedInHyperrectangle(const Eigen::MatrixBase< Derived > &sidelengths, const Eigen::MatrixBase< Derived > &center, const Eigen::MatrixBase< Derived > &testPoint)
Tests, if a vertex is contained in a hyperrectangle.
Definition geometry.hpp:195
bool lineIntersection(const Eigen::Ref< const Eigen::Vector2d > &a, const Eigen::Ref< const Eigen::Vector2d > &b, const Eigen::Ref< const Eigen::Vector2d > &c, const Eigen::Ref< const Eigen::Vector2d > &d, Eigen::Ref< Eigen::Vector2d > &intersectionPoint)
Determines the intersection point of two lines, if one exists.
Definition geometry.cpp:16
static bool parallel(const Eigen::MatrixBase< Derived > &a, const Eigen::MatrixBase< Derived > &b, const Eigen::MatrixBase< Derived > &c, const Eigen::MatrixBase< Derived > &d)
Determines, if two lines are parallel to each other.
Definition geometry.hpp:182
ResultConstants segmentPlaneIntersection(const Eigen::Vector3d &pointOnPlane, const Eigen::Vector3d &planeNormal, const Eigen::Vector3d &firstPointSegment, const Eigen::Vector3d &secondPointSegment, Eigen::Vector3d &intersectionPoint)
Determines the intersection point of a segment with a plane in 3D.
Definition geometry.cpp:42
bool collinear(const Eigen::MatrixBase< Derived > &a, const Eigen::MatrixBase< Derived > &b, const Eigen::MatrixBase< Derived > &c)
Determines, if three points are collinear (on one line)
Definition geometry.hpp:167
ConvexityResult isConvexQuad(std::array< Eigen::VectorXd, 4 > coords)
Definition geometry.cpp:143
constexpr bool equals(const Eigen::MatrixBase< DerivedA > &A, const Eigen::MatrixBase< DerivedB > &B, double tolerance=NUMERICAL_ZERO_DIFFERENCE)
Compares two Eigen::MatrixBase for equality up to tolerance.
std::enable_if< std::is_arithmetic< Scalar >::value, bool >::type greaterEquals(Scalar A, Scalar B, Scalar tolerance=NUMERICAL_ZERO_DIFFERENCE)
std::enable_if< std::is_arithmetic< Scalar >::value, bool >::type greater(Scalar A, Scalar B, Scalar tolerance=NUMERICAL_ZERO_DIFFERENCE)