Lagrange
ExactPredicates Class Referenceabstract

Inherited by ExactPredicatesShewchuk.

Public Member Functions

virtual short collinear3D (double p1[3], double p2[3], double p3[3]) const
 Tests whether p1, p2, and p3 are collinear in 3D. More...
 
virtual short orient2D (double p1[2], double p2[2], double p3[2]) const =0
 Exact 2D orientation test. More...
 
virtual short orient3D (double p1[3], double p2[3], double p3[3], double p4[3]) const =0
 Exact 3D orientation test. More...
 
virtual short incircle (double p1[2], double p2[2], double p3[2], double p4[2]) const =0
 Exact 2D incircle test. More...
 
virtual short insphere (double p1[3], double p2[3], double p3[3], double p4[3], double p5[3]) const =0
 Exact 3D insphere test. More...
 

Static Public Member Functions

static std::unique_ptr< ExactPredicatescreate (const std::string &engine)
 Factory method to create an exact predicate engine. More...
 

Member Function Documentation

◆ create()

std::unique_ptr< ExactPredicates > create ( const std::string &  engine)
static

Factory method to create an exact predicate engine.

Currently only the only engine supported is "shewchuk", so you might as well create it directly.

Parameters
[in]engineEngine name.
Returns
The created engine.

◆ collinear3D()

short collinear3D ( double  p1[3],
double  p2[3],
double  p3[3] 
) const
virtual

Tests whether p1, p2, and p3 are collinear in 3D.

This works by calling orient2D successively on the xy, yz and zx projections of p1, p2 p3.

Parameters
p1First 3D point.
p2Second 3D point.
p3Third 3D point.
Returns
1 if the points are collinear, 0 otherwise.

◆ orient2D()

virtual short orient2D ( double  p1[2],
double  p2[2],
double  p3[2] 
) const
pure virtual

Exact 2D orientation test.

Parameters
p1First 2D point.
p2Second 2D point.
p3Third 2D point.
Returns
Return a positive value if the points p1, p2, and p3 occur in counterclockwise order; a negative value if they occur in clockwise order; and zero if they are collinear.

Implemented in ExactPredicatesShewchuk.

◆ orient3D()

virtual short orient3D ( double  p1[3],
double  p2[3],
double  p3[3],
double  p4[3] 
) const
pure virtual

Exact 3D orientation test.

Parameters
p1First 3D point.
p2Second 3D point.
p3Third 3D point.
p4Fourth 3D point.
Returns
Return a positive value if the point pd lies below the plane passing through p1, p2, and p3; "below" is defined so that p1, p2, and p3 appear in counterclockwise order when viewed from above the plane. Returns a negative value if pd lies above the plane. Returns zero if the points are coplanar.

Implemented in ExactPredicatesShewchuk.

◆ incircle()

virtual short incircle ( double  p1[2],
double  p2[2],
double  p3[2],
double  p4[2] 
) const
pure virtual

Exact 2D incircle test.

Parameters
p1First 2D point.
p2Second 2D point.
p3Third 2D point.
p4Fourth 3D point.
Returns
Return a positive value if the point pd lies inside the circle passing through p1, p2, and p3; a negative value if it lies outside; and zero if the four points are cocircular. The points p1, p2, and p3 must be in counterclockwise order, or the sign of the result will be reversed.

Implemented in ExactPredicatesShewchuk.

◆ insphere()

virtual short insphere ( double  p1[3],
double  p2[3],
double  p3[3],
double  p4[3],
double  p5[3] 
) const
pure virtual

Exact 3D insphere test.

Parameters
p1First 3D point.
p2Second 3D point.
p3Third 3D point.
p4Fourth 3D point.
p5Fifth 3D point.
Returns
Return a positive value if the point p5 lies inside the sphere passing through p1, p2, p3, and p4; a negative value if it lies outside; and zero if the five points are cospherical. The points p1, p2, p3, and p4 must be ordered so that they have a positive orientation (as defined by orient3d()), or the sign of the result will be reversed.

Implemented in ExactPredicatesShewchuk.


The documentation for this class was generated from the following files: