preCICE v3.1.1
Loading...
Searching...
No Matches
Classes | Enumerations | Functions | Variables
precice::math::geometry Namespace Reference

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 > &center, 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")
 

Detailed Description

Provides computational geometry operations.

Enumeration Type Documentation

◆ ResultConstants

Enumerator
NO_INTERSECTION 
INTERSECTION 
TOUCHING 
NOT_CONTAINED 
CONTAINED 

Definition at line 14 of file geometry.hpp.

Function Documentation

◆ between()

template<typename VECTORA_T , typename VECTORB_T , typename VECTORC_T >
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.

Parameters
[in]aFirst point of line segment
[in]bSecond point of line segment
[in]toCheckPoint to be checked for "betweenness"
Returns
True, if toCheck lies between a and b. False, otherwise

Definition at line 151 of file geometry.hpp.

Here is the call graph for this function:

◆ collinear()

template<typename Derived >
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.

Parameters
[in]aFirst point to check
[in]bSecond point to check
[in]cThird point to check
Returns
True, if collinear. False, otherwise

Definition at line 170 of file geometry.hpp.

Here is the call graph for this function:

◆ containedInHyperrectangle()

template<class Derived >
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.

Returns
One of:
  • CONTAINED
  • NOT_CONTAINED
  • TOUCHING

Definition at line 198 of file geometry.hpp.

Here is the call graph for this function:

◆ isConvexQuad()

ConvexityResult precice::math::geometry::isConvexQuad ( std::array< Eigen::VectorXd, 4 > coords)

Definition at line 143 of file geometry.cpp.

Here is the call graph for this function:

◆ lineIntersection()

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.

Parameters
[in]aFirst point on first line.
[in]bSecond point on first line.
[in]cFirst point on second line.
[in]dSecond point on second line.
[out]intersectionPointPoint of intersection.
Returns
True, if an intersection point exists and could be determined.

Definition at line 16 of file geometry.cpp.

Here is the call graph for this function:

◆ parallel() [1/2]

template<typename Derived >
static 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 )
static

Determines, if two lines are parallel to each other.

Works for Dim2 and Dim3.

Parameters
[in]aFirst point on first line.
[in]bSecond point on first line.
[in]cFirst point on second line.
[in]dSecond point on second line.
Returns
True, if the two lines are parallel to each other.

Definition at line 185 of file geometry.hpp.

Here is the call graph for this function:

◆ parallel() [2/2]

template<class Derived >
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.

Parameters
[in]aFirst point on first line.
[in]bSecond point on first line.
[in]cFirst point on second line.
[in]dSecond point on second line.
Returns
True, if the two lines are parallel to each other.

Definition at line 185 of file geometry.hpp.

Here is the call graph for this function:

◆ projectVector()

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.

◆ segmentPlaneIntersection()

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.

Parameters
[out]intersectionPointContains coordinates of intersection point, if return value equals TOUCHING or INTERSECTION.
Returns
Returns type of intersection. One of
  • NO_INTERSECTION
  • INTERSECTION
  • TOUCHING
  • CONTAINED

Definition at line 42 of file geometry.cpp.

Here is the call graph for this function:

◆ tetraVolume()

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.

◆ triangleArea()

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.

Here is the call graph for this function:

Variable Documentation

◆ _log

logging::Logger precice::math::geometry::_log("math::geometry") ( "math::geometry" )