3#include <Eigen/Geometry>
17 const Eigen::Ref<const Eigen::Vector2d> &a,
18 const Eigen::Ref<const Eigen::Vector2d> &b,
19 const Eigen::Ref<const Eigen::Vector2d> &c,
20 const Eigen::Ref<const Eigen::Vector2d> &d,
21 Eigen::Ref<Eigen::Vector2d> & intersectionPoint)
24 double D = a(0) * (d(1) - c(1)) +
25 b(0) * (c(1) - d(1)) +
26 d(0) * (b(1) - a(1)) -
36 double s = (a(0) * (d(1) - c(1)) + c(0) * (a(1) - d(1)) + d(0) * (c(1) - a(1))) / D;
38 intersectionPoint = a + s * (b - a);
43 const Eigen::Vector3d &pointOnPlane,
44 const Eigen::Vector3d &planeNormal,
45 const Eigen::Vector3d &firstPointSegment,
46 const Eigen::Vector3d &secondPointSegment,
47 Eigen::Vector3d & intersectionPoint)
53 double d = pointOnPlane.dot(planeNormal);
64 Eigen::Vector3d segmentVec(secondPointSegment - firstPointSegment);
65 double nominator = d - firstPointSegment.dot(planeNormal);
66 double denominator = segmentVec.dot(planeNormal);
74 double t = nominator / denominator;
82 intersectionPoint = firstPointSegment + segmentVec * t;
94 const Eigen::VectorXd &a,
95 const Eigen::VectorXd &b,
96 const Eigen::VectorXd &c)
101 Eigen::Vector2d A = b;
103 Eigen::Vector2d B = c;
105 return 0.5 *
std::fabs(A(0) * B(1) - A(1) * B(0));
108 Eigen::Vector3d A = b;
110 Eigen::Vector3d B = c;
112 return 0.5 * A.cross(B).norm();
117 const Eigen::Vector3d &a,
118 const Eigen::Vector3d &b,
119 const Eigen::Vector3d &c,
120 const Eigen::Vector3d &d)
122 return std::abs((a - d).dot((b - d).cross(c - d))) / 6.0;
126 const Eigen::Vector3d &vector,
127 int indexDimensionToRemove)
131 Eigen::Vector2d projectedVector;
133 for (
int dim = 0; dim < 3; dim++) {
134 if (indexDimensionToRemove == dim) {
137 projectedVector(reducedDim) = vector(dim);
140 return projectedVector;
150 [](
const auto &v) { return v.size() == 3; }),
151 "This only works in 3D.");
153 Eigen::Vector3d coordOrigin;
154 coordOrigin = coords[0];
157 Eigen::Vector3d e_1 = coords[1] - coordOrigin;
158 Eigen::Vector3d e_2 = coords[2] - coordOrigin;
159 Eigen::Vector3d normalVector = e_1.cross(e_2);
162 "Non-planar quads are not supported. The vertex coordinates are: {}.", coords);
165 for (
int i = 0; i < 4; i++) {
166 Eigen::Vector3d coordinateDifference = coords[i] - coordOrigin;
167 coords[i][0] = e_1.dot(coordinateDifference);
168 coords[i][1] = e_2.dot(coordinateDifference);
169 coords[i][2] = normalVector.dot(coordinateDifference);
179 int idLowestPoint = 0;
180 for (
int i = 1; i < 4; i++) {
181 if (coords[i][0] < coords[idLowestPoint][0]) {
188 int validVertexIDCounter = 0;
189 int currentVertex = idLowestPoint;
193 result.
vertexOrder[validVertexIDCounter] = currentVertex;
196 int nextVertex = (currentVertex + 1) % 4;
197 for (
int i = 0; i < 4; i++) {
198 double y1 = coords[currentVertex][1] - coords[nextVertex][1];
199 double y2 = coords[currentVertex][1] - coords[i][1];
200 double x1 = coords[currentVertex][0] - coords[nextVertex][0];
201 double x2 = coords[currentVertex][0] - coords[i][0];
202 double val = y2 * x1 - y1 * x2;
211 currentVertex = nextVertex;
212 validVertexIDCounter++;
213 }
while (currentVertex != idLowestPoint);
216 result.convex = (validVertexIDCounter == 4);
#define PRECICE_CHECK(check,...)
#define PRECICE_ASSERT(...)
This class provides a lightweight logger.
Provides computational geometry operations.
double triangleArea(const Eigen::VectorXd &a, const Eigen::VectorXd &b, const Eigen::VectorXd &c)
Computes the signed area of a triangle in 2D.
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.
Eigen::Vector2d projectVector(const Eigen::Vector3d &vector, int indexDimensionToRemove)
Projects a 3D vector to a 2D one by removing one dimension.
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.
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.
logging::Logger _log("math::geometry")
ConvexityResult isConvexQuad(std::array< Eigen::VectorXd, 4 > coords)
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 greater(Scalar A, Scalar B, Scalar tolerance=NUMERICAL_ZERO_DIFFERENCE)
std::array< int, 4 > vertexOrder