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);
50 const Eigen::Vector3d &pointOnPlane,
51 const Eigen::Vector3d &planeNormal,
52 const Eigen::Vector3d &firstPointSegment,
53 const Eigen::Vector3d &secondPointSegment,
54 Eigen::Vector3d &intersectionPoint);
67template <
typename VECTORA_T,
typename VECTORB_T,
typename VECTORC_T>
71 const VECTORC_T &toCheck);
84template <
typename Derived>
86 const Eigen::MatrixBase<Derived> &a,
87 const Eigen::MatrixBase<Derived> &b,
88 const Eigen::MatrixBase<Derived> &c);
101template <
typename Derived>
103 const Eigen::MatrixBase<Derived> &a,
104 const Eigen::MatrixBase<Derived> &b,
105 const Eigen::MatrixBase<Derived> &c,
106 const Eigen::MatrixBase<Derived> &d);
115 const Eigen::VectorXd &a,
116 const Eigen::VectorXd &b,
117 const Eigen::VectorXd &c);
121 const Eigen::Vector3d &a,
122 const Eigen::Vector3d &b,
123 const Eigen::Vector3d &c,
124 const Eigen::Vector3d &d);
128 const Eigen::Vector3d &
vector,
129 const int indexDimensionToRemove);
139template <
class Derived>
141 const Eigen::MatrixBase<Derived> &sidelengths,
142 const Eigen::MatrixBase<Derived> ¢er,
143 const Eigen::MatrixBase<Derived> &testPoint);
147template <
typename VECTORA_T,
typename VECTORB_T,
typename VECTORC_T>
151 const VECTORC_T &toCheck)
166template <
typename Derived>
168 const Eigen::MatrixBase<Derived> &a,
169 const Eigen::MatrixBase<Derived> &b,
170 const Eigen::MatrixBase<Derived> &c)
174 double triangleOutline = (b - a).norm() + (c - b).norm() + (a - c).norm();
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)
194template <
class Derived>
196 const Eigen::MatrixBase<Derived> &sidelengths,
197 const Eigen::MatrixBase<Derived> ¢er,
198 const Eigen::MatrixBase<Derived> &testPoint)
200 int dim = sidelengths.size();
201 Eigen::VectorXd toCenter = testPoint - center;
202 toCenter = toCenter.cwiseAbs();
205 bool touching =
false;
206 for (
int i = 0; i < dim; i++) {
207 diff = 0.5 * sidelengths(i) - toCenter(i);
#define PRECICE_ASSERT(...)
Provides computational geometry operations.
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.
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.
int containedInHyperrectangle(const Eigen::MatrixBase< Derived > &sidelengths, const Eigen::MatrixBase< Derived > ¢er, const Eigen::MatrixBase< Derived > &testPoint)
Tests, if a vertex is contained in a hyperrectangle.
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.
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.
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.
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)
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 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)
std::array< int, 4 > vertexOrder