|
preCICE v3.3.0
|
Provides computational geometry operations. More...
Classes | |
| struct | ConvexityResult |
Enumerations | |
| enum | ResultConstants { NO_INTERSECTION , INTERSECTION , TOUCHING , NOT_CONTAINED , CONTAINED } |
Functions | |
| 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. | |
| 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, const int indexDimensionToRemove) |
| Projects a 3D vector to a 2D one by removing one dimension. | |
| ConvexityResult | isConvexQuad (std::array< Eigen::VectorXd, 4 > coords) |
| template<typename VECTORA_T, typename VECTORB_T, typename VECTORC_T> | |
| 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. | |
| template<typename Derived> | |
| 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) | |
| template<typename Derived> | |
| 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. | |
| template<class Derived> | |
| 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. | |
| template<class Derived> | |
| 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. | |
Variables | |
| logging::Logger | _log ("math::geometry") |
Provides computational geometry operations.
| Enumerator | |
|---|---|
| NO_INTERSECTION | |
| INTERSECTION | |
| TOUCHING | |
| NOT_CONTAINED | |
| CONTAINED | |
Definition at line 11 of file geometry.hpp.
| bool precice::math::geometry::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.
Works for 2D and 3D.
| [in] | a | First point of line segment |
| [in] | b | Second point of line segment |
| [in] | toCheck | Point to be checked for "betweenness" |
Definition at line 148 of file geometry.hpp.
| bool precice::math::geometry::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)
Works for Dim2 and Dim3.
| [in] | a | First point to check |
| [in] | b | Second point to check |
| [in] | c | Third point to check |
Definition at line 167 of file geometry.hpp.
| int precice::math::geometry::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 at line 195 of file geometry.hpp.
| ConvexityResult precice::math::geometry::isConvexQuad | ( | std::array< Eigen::VectorXd, 4 > | coords | ) |
| bool precice::math::geometry::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.
Works only for Dim2.
| [in] | a | First point on first line. |
| [in] | b | Second point on first line. |
| [in] | c | First point on second line. |
| [in] | d | Second point on second line. |
| [out] | intersectionPoint | Point of intersection. |
Definition at line 16 of file geometry.cpp.
|
static |
Determines, if two lines are parallel to each other.
Works for Dim2 and Dim3.
| [in] | a | First point on first line. |
| [in] | b | Second point on first line. |
| [in] | c | First point on second line. |
| [in] | d | Second point on second line. |
Definition at line 182 of file geometry.hpp.
| bool precice::math::geometry::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.
Works for Dim2 and Dim3.
| [in] | a | First point on first line. |
| [in] | b | Second point on first line. |
| [in] | c | First point on second line. |
| [in] | d | Second point on second line. |
Definition at line 182 of file geometry.hpp.
| Eigen::Vector2d precice::math::geometry::projectVector | ( | const Eigen::Vector3d & | vector, |
| int | indexDimensionToRemove ) |
Projects a 3D vector to a 2D one by removing one dimension.
Definition at line 125 of file geometry.cpp.
| ResultConstants precice::math::geometry::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.
| [out] | intersectionPoint | Contains coordinates of intersection point, if return value equals TOUCHING or INTERSECTION. |
Definition at line 42 of file geometry.cpp.
| double precice::math::geometry::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 at line 116 of file geometry.cpp.
| double precice::math::geometry::triangleArea | ( | const Eigen::VectorXd & | a, |
| const Eigen::VectorXd & | b, | ||
| const Eigen::VectorXd & | c ) |
Computes the signed area of a triangle in 2D.
The area is negative, if the points a, b, c are given in clockwise ordering, otherwise positive.
Definition at line 93 of file geometry.cpp.
| logging::Logger precice::math::geometry::_log("math::geometry") | ( | "math::geometry" | ) |