Skip to content

Namespace SFCGAL::algorithm#

Namespace List > algorithm

Main algorithm namespace.

Classes#

Type Name
class BoundaryVisitor
class Buffer3D
Computes a 3D buffer around a Point orLineString .
struct Circle
class ConsistentOrientationBuilder
struct LoopDetector
struct NoValidityCheck
struct Plane3DInexactUnsafe
struct Sphere
class SurfaceGraph
class WeightedCentroid
struct found_an_intersection
struct intersects_cb <Dim>

Public Types#

Type Name
enum PartitionAlgorithm

Public Functions#

Type Name
auto alphaShapes (const Geometry & g, double alpha=1, bool allow_holes=false)
auto alphaWrapping3D (const Geometry & geom, size_t relativeAlpha, size_t relativeOffset)
auto approximate (const Offset_polygon_2 & polygon, const int & n=0)
approximate an Offset_polygon_2 (filter null segments)
auto approximate (const Offset_polygon_with_holes_2 & polygon, const int & n=0)
approximate an Offset
auto approximateMedialAxis (const Geometry & geom)
build an approximate medial axis for a Polygon __
auto area (const Geometry & g, NoValidityCheck)
Compute the 2D area for a Geometry.
auto area (const Geometry & g)
Compute the 2D area for a Geometry.
auto area (const Triangle & g)
auto area (const Polygon & g)
auto area (const GeometryCollection & g)
auto area (const TriangulatedSurface & g)
auto area (const PolyhedralSurface & g)
auto area3D (const Geometry & g, NoValidityCheck)
auto area3D (const Geometry & g)
auto area3D (const Polygon & g)
auto area3D (const Triangle & g)
auto area3D (const GeometryCollection & g)
auto area3D (const PolyhedralSurface & g)
auto area3D (const TriangulatedSurface & g)
auto boundingCircle (const Geometry & geom)
auto boundingSphere (const Geometry & geom)
auto centroid (const Geometry & g)
Returns the 2D centroid for a Geometry.
auto centroid3D (const Geometry & g)
Returns the 3D centroid for a Geometry.
auto circleToPolygon (const Kernel::Circle_2 & circle)
helper to create a polygon from a circle
auto collect (const Geometry & ga, const Geometry & gb)
std::unique_ptr< Geometry > collect (GeometryIterator begin, GeometryIterator end)
auto collectionExtractPolygons (std::unique_ptr< Geometry > coll)
auto collectionHomogenize (std::unique_ptr< Geometry > coll)
auto collectionToMulti (std::unique_ptr< Geometry > coll)
auto convexHull (const Geometry & g)
auto convexHull3D (const Geometry & g)
auto covers (const Geometry & ga, const Geometry & gb)
bool covers (const detail::GeometrySet< Dim > & a, const detail::GeometrySet< Dim > & b)
bool covers (const detail::PrimitiveHandle< Dim > & a, const detail::PrimitiveHandle< Dim > & b)
auto covers3D (const Geometry & ga, const Geometry & gb)
auto difference (const Geometry & ga, const Geometry & gb, NoValidityCheck)
auto difference (const Geometry & ga, const Geometry & gb)
auto difference3D (const Geometry & ga, const Geometry & gb, NoValidityCheck)
auto difference3D (const Geometry & ga, const Geometry & gb)
auto distance (const Geometry & gA, const Geometry & gB, NoValidityCheck)
auto distance (const Geometry & gA, const Geometry & gB)
auto distance3D (const Geometry & gA, const Geometry & gB, NoValidityCheck)
auto distance3D (const Geometry & gA, const Geometry & gB)
auto distanceGeometryCollectionToGeometry (const Geometry & gA, const Geometry & gB)
auto distanceGeometryCollectionToGeometry3D (const Geometry & gA, const Geometry & gB)
auto distanceLineStringGeometry (const LineString & gA, const Geometry & gB)
auto distanceLineStringGeometry3D (const LineString & gA, const Geometry & gB)
auto distanceLineStringLineString (const LineString & gA, const LineString & gB)
auto distanceLineStringLineString3D (const LineString & gA, const LineString & gB)
auto distanceLineStringPolygon (const LineString & gA, const Polygon & gB)
auto distanceLineStringPolygon3D (const LineString & gA, const Polygon & gB)
auto distanceLineStringPolyhedralSurface3D (const LineString & lineA, const PolyhedralSurface & polySurfaceB)
auto distanceLineStringSolid3D (const LineString & gA, const Solid & gB)
auto distanceLineStringTriangle (const LineString & gA, const Triangle & gB)
auto distanceLineStringTriangle3D (const LineString & gA, const Triangle & gB)
auto distanceLineStringTriangulatedSurface3D (const LineString & lineA, const TriangulatedSurface & triangulatedSurfaceB)
auto distancePointGeometry (const Point & gA, const Geometry & gB)
auto distancePointGeometry3D (const Point & gA, const Geometry & gB)
auto distancePointLineString (const Point & gA, const LineString & gB)
auto distancePointLineString3D (const Point & gA, const LineString & gB)
auto distancePointPoint (const Point & gA, const Point & gB)
auto distancePointPoint3D (const Point & gA, const Point & gB)
auto distancePointPolygon (const Point & gA, const Polygon & gB)
auto distancePointPolygon3D (const Point & gA, const Polygon & gB)
auto distancePointPolyhedralSurface3D (const Point & pointA, const PolyhedralSurface & polySurfaceB)
auto distancePointSegment (const Point & p, const Point & a, const Point & b)
auto distancePointSegment3D (const Point & p, const Point & a, const Point & b)
auto distancePointSolid3D (const Point & gA, const Solid & gB)
auto distancePointTriangle (const Point & gA, const Triangle & gB)
auto distancePointTriangle3D (const Point & gA, const Triangle & gB)
auto distancePointTriangle3D (const Point & p_, const Point & a_, const Point & b_, const Point & c_)
auto distancePointTriangulatedSurface3D (const Point & pointA, const TriangulatedSurface & triangulatedSurfaceB)
auto distancePolygonGeometry (const Polygon & gA, const Geometry & gB)
auto distancePolygonGeometry3D (const Polygon & gA, const Geometry & gB)
auto distancePolygonPolygon (const Polygon & gA, const Polygon & gB)
auto distancePolygonTriangle (const Polygon & gA, const Triangle & gB)
auto distancePolyhedralSurfaceGeometry3D (const PolyhedralSurface & polySurfaceA, const Geometry & geomB)
auto distanceSegmentSegment (const Point & a, const Point & b, const Point & c, const Point & d)
auto distanceSegmentSegment3D (const Point & a, const Point & b, const Point & c, const Point & d)
auto distanceSegmentTriangle3D (const Point & sA_, const Point & sB_, const Point & tA_, const Point & tB_, const Point & tC_)
auto distanceSolidGeometry3D (const Solid & gA, const Geometry & gB)
auto distanceSolidSolid3D (const Solid & gA, const Solid & gB)
auto distanceTriangleGeometry (const Triangle & gA, const Geometry & gB)
auto distanceTriangleGeometry3D (const Triangle & gA, const Geometry & gB)
auto distanceTrianglePolyhedralSurface3D (const Triangle & triangleA, const PolyhedralSurface & polySurfaceB)
auto distanceTriangleSolid3D (const Triangle & gA, const Solid & gB)
auto distanceTriangleTriangle3D (const Triangle & gA, const Triangle & gB)
auto distanceTriangleTriangulatedSurface3D (const Triangle & triangleA, const TriangulatedSurface & triangulatedSurfaceB)
auto distanceTriangulatedSurfaceGeometry3D (const TriangulatedSurface & triangulatedSurfaceA, const Geometry & geomB)
auto extrude (const Geometry & g, const Kernel::Vector_3 & v)
Returns a Geometry equal to the specified Geometry, extruded by the specified displacement vector.
auto extrude (const Geometry & g, const Kernel::FT & dx, const Kernel::FT & dy, const Kernel::FT & dz, NoValidityCheck)
auto extrude (const Geometry & g, const Kernel::FT & dx, const Kernel::FT & dy, const Kernel::FT & dz)
Returns a Geometry equal to the specified Geometry, extruded by the specified displacement.
SFCGAL_API auto extrude (const Geometry & g, const double & dx, const double & dy, const double & dz)
Returns a Geometry equal to the specified Geometry, extruded by the specified displacement.
SFCGAL_API auto extrude (const Polygon & g, const double & height)
SFCGAL_API std::unique_ptr< Geometry > extrude (const Geometry & g, const Kernel::FT & dx, const Kernel::FT & dy, const Kernel::FT & dz, NoValidityCheck & nvc)
Returns a Geometry equal to the specified Geometry, extruded by the specified displacement.
auto extrudeStraightSkeleton (const Polygon & geom, double height)
auto extrudeStraightSkeleton (const Geometry & geom, double height)
auto extrudeStraightSkeleton (const Geometry & geom, double building_height, double roof_height)
SFCGAL_API auto extrudedStraightSkeleton (const Polygon & geom, double height)
build a 3D straight skeleton extruded for a Polygon __
void force2D (Geometry & g)
force a geometry to be 2D (project on O,x,y)
void force3D (Geometry & g, const Kernel::FT & defaultZ=0)
force a 2D geometry to be 3D (replace undefined Z by defaultZ, existing Z values remains unchanged)
void forceMeasured (Geometry & geometry, const double & defaultM=0)
force a 2D or M geometry to be M (replace undefined M by defaultM, existing M values remains unchanged)
auto hasConsistentOrientation3D (const TriangulatedSurface & g)
auto hasConsistentOrientation3D (const PolyhedralSurface & g)
bool hasPlane3D (const Polygon & polygon, CGAL::Point_3< Kernel > & a, CGAL::Point_3< Kernel > & b, CGAL::Point_3< Kernel > & c)
Test if a 3D plane can be extracted from a Polygon .
bool hasPlane3D (const Polygon & polygon)
auto intersection (const CGAL::Triangle_3< Kernel > & a, const CGAL::Triangle_3< Kernel > & b)
void intersection (const GeometrySet< Dim > & a, const GeometrySet< Dim > & b, GeometrySet< Dim > & output)
auto intersection (const Geometry & ga, const Geometry & gb, NoValidityCheck)
auto intersection (const Geometry & ga, const Geometry & gb)
auto intersection3D (const Geometry & ga, const Geometry & gb, NoValidityCheck)
auto intersection3D (const Geometry & ga, const Geometry & gb)
auto intersects (const detail::PrimitiveHandle< Dim > & a, const detail::PrimitiveHandle< Dim > & b)
auto intersects (const detail::GeometrySet< Dim > & a, const detail::GeometrySet< Dim > & b)
auto intersects (const Geometry & ga, const Geometry & gb)
auto intersects (const Geometry & ga, const Geometry & gb, NoValidityCheck)
bool intersects (const detail::GeometrySet< Dim > & a, const detail::GeometrySet< Dim > & b)
bool intersects (const detail::PrimitiveHandle< Dim > & a, const detail::PrimitiveHandle< Dim > & b)
auto intersects3D (const Geometry & ga, const Geometry & gb)
auto intersects3D (const Geometry & ga, const Geometry & gb, NoValidityCheck)
auto isClosed (const SurfaceGraph & graph)
auto isConnected (const SurfaceGraph & graph)
auto isCounterClockWiseOriented (const LineString & ls)
auto isCounterClockWiseOriented (const Triangle & tri)
auto isCounterClockWiseOriented (const Polygon & poly)
bool isPlane3D (const Geometry & geom, const double & toleranceAbs)
auto isSimple (const Geometry & g, const double & toleranceAbs=1e-9)
auto isValid (const Point & p)
auto isValid (const LineString & l, const double & toleranceAbs)
auto isValid (const Polygon & p, const double & toleranceAbs)
auto isValid (const Triangle & t, const double & toleranceAbs)
auto isValid (const MultiLineString & ml, const double & toleranceAbs)
auto isValid (const MultiPolygon & mp, const double & toleranceAbs)
auto isValid (const GeometryCollection & gc, const double & toleranceAbs)
auto isValid (const TriangulatedSurface & tin, const SurfaceGraph & graph, const double & toleranceAbs)
auto isValid (const TriangulatedSurface & tin, const double & toleranceAbs)
auto isValid (const PolyhedralSurface & s, const SurfaceGraph & graph, const double & toleranceAbs)
auto isValid (const PolyhedralSurface & s, const double & toleranceAbs)
auto isValid (const Solid & solid, const double & toleranceAbs)
auto isValid (const MultiSolid & ms, const double & toleranceAbs)
auto isValid (const Geometry & g, const double & toleranceAbs)
Check validity of a geometry.
auto length (const LineString & g)
Compute the 2D length for a LineString .
auto length (const GeometryCollection & g)
Compute the 2D length for a GeometryCollection .
auto length (const Geometry & g)
Compute the 2D length for a Geometry (0 for incompatible types)
auto length3D (const LineString & g)
Compute the 3D length for a LineString .
auto length3D (const GeometryCollection & g)
Compute the 3D length for a GeometryCollection .
auto length3D (const Geometry & g)
Compute the 2D length for a geometry.
auto lineSubstring (const LineString & ls, double start, double end)
Retrieve a substring of a specified LineString , between the specified fractional distances from the start of the specifiedLineString .
void makeConsistentOrientation3D (TriangulatedSurface & g)
void makeValidOrientation (CGAL::Polygon_2< Kernel > & polygon)
void makeValidOrientation (CGAL::Polygon_with_holes_2< Kernel > & polygon)
void makeValidOrientation (Polygon & polygon)
auto minkowskiSum (const Geometry & gA, const Polygon & gB, NoValidityCheck)
2D minkowski sum (p+q)
auto minkowskiSum (const Geometry & gA, const Polygon & gB)
2D minkowski sum (p+q)
auto minkowskiSum3D (const Geometry & gA, const Geometry & gB, NoValidityCheck)
3D Minkowski sum (p+q)
auto minkowskiSum3D (const Geometry & gA, const Geometry & gB)
3D Minkowski sum (p+q)
CGAL::Vector_3< Kernel > normal3D (const CGAL::Point_3< Kernel > & a, const CGAL::Point_3< Kernel > & b, const CGAL::Point_3< Kernel > & c)
CGAL::Vector_3< Kernel > normal3D (const LineString & ls, bool exact=true)
CGAL::Vector_3< Kernel > normal3D (const Polygon & polygon, bool exact=true)
void offset (const Geometry & g, const double & radius, Offset_polygon_set_2 & polygonSet)
dispatch a geometry
void offset (const Point & g, const double & radius, Offset_polygon_set_2 & polygonSet)
offset for a Point __
void offset (const LineString & g, const double & radius, Offset_polygon_set_2 & polygonSet)
offset for a LineString __
void offset (const Polygon & g, const double & radius, Offset_polygon_set_2 & polygonSet)
offset for a Polygon __
auto offset (const Geometry & g, const double & r, NoValidityCheck)
[experimental]compute polygon offset
auto offset (const Geometry & g, const double & r)
[experimental]compute polygon offset
void offsetCollection (const Geometry & g, const double & radius, Offset_polygon_set_2 & polygonSet)
offset for MultiPoint ,MultiLineString ,MultiPolygon ,TriangulatedSurface ,PolyhedralSurface __
auto optimal_alpha_shapes (const Geometry & g, bool allow_holes=false, size_t nb_components=1)
auto partition_2 (const Geometry & g, PartitionAlgorithm alg, NoValidityCheck)
auto partition_2 (const Geometry & g, PartitionAlgorithm alg=y_monotone)
void plane3D (const Polygon & polygon, CGAL::Point_3< Kernel > & a, CGAL::Point_3< Kernel > & b, CGAL::Point_3< Kernel > & c)
CGAL::Plane_3< Kernel > plane3D (const Polygon & polygon)
CGAL::Plane_3< Kernel > plane3D (const Polygon & polygon, const Plane3DInexactUnsafe &)
CGAL::Plane_3< Kernel > plane3D (const Polygon & polygon, bool exact)
auto polygonSetToMultiPolygon (const Offset_polygon_set_2 & polygonSet, const int & n)
convert Offset_polygon_set_2 to MultiPolygon __
void propagateValidityFlag (Geometry & g, bool valid)
void rotate (Geometry & g, const Kernel::FT & angle)
Rotate a geometry in 2D around the origin (0,0)
void rotate (Geometry & g, const Kernel::FT & angle, const Point & origin)
Rotate a geometry in 2D around a specified point.
void rotate (Geometry & g, const Kernel::FT & angle, const Kernel::Vector_3 & axis, const Point & origin=Point(0, 0, 0))
Rotate a geometry in 3D around a specified axis and origin.
void rotateX (Geometry & g, const Kernel::FT & angle)
Rotate a geometry around the X axis.
void rotateY (Geometry & g, const Kernel::FT & angle)
Rotate a geometry around the Y axis.
void rotateZ (Geometry & g, const Kernel::FT & angle)
Rotate a geometry around the Z axis.
void scale (Geometry & g, double s)
void scale (Geometry & g, double sx, double sy, double sz=0.0)
void scale (Geometry & g, double sx, double sy, double sz, double cx, double cy, double cz)
auto selfIntersects (const LineString & l)
auto selfIntersects (const PolyhedralSurface & s, const SurfaceGraph & g)
auto selfIntersects (const TriangulatedSurface & s, const SurfaceGraph & g)
auto selfIntersects3D (const LineString & l)
auto selfIntersects3D (const PolyhedralSurface & s, const SurfaceGraph & g)
auto selfIntersects3D (const TriangulatedSurface & s, const SurfaceGraph & g)
auto signedArea (const Triangle & g)
Compute the 2D signed area for a Triangle .
auto signedArea (const LineString & g)
Compute the 2D signed area for a closed LineString .
auto simplify (const Geometry & geometry, double threshold, bool preserveTopology)
Main entry point for geometry simplification.
auto squaredDistancePointTriangle3D (const Point_3 & p, const Triangle_3 & abc)
auto squaredDistanceSegmentTriangle3D (const Segment_3 & sAB, const Triangle_3 & tABC)
auto squaredDistanceTriangleTriangle3D (const Triangle_3 & triangleA, const Triangle_3 & triangleB)
auto straightSkeleton (const Geometry & geom, bool autoOrientation, NoValidityCheck, bool innerOnly=false, bool outputDistanceInM=false, const double & toleranceAbs=EPSILON)
build a 2D straight skeleton for a Polygon __
auto straightSkeleton (const Geometry & geom, bool autoOrientation=true, bool innerOnly=false, bool outputDistanceInM=false, const double & toleranceAbs=EPSILON)
build a 2D straight skeleton for a Polygon __
auto straightSkeleton (const Polygon & geom, bool autoOrientation=true, bool innerOnly=false, bool outputDistanceInM=false, const double & toleranceAbs=EPSILON)
build a 2D straight skeleton for a Polygon __
auto straightSkeleton (const MultiPolygon & geom, bool autoOrientation=true, bool innerOnly=false, bool outputDistanceInM=false, const double & toleranceAbs=EPSILON)
build a 2D straight skeleton for a Polygon __
auto straightSkeletonPartition (const Geometry & geom, bool autoOrientation=true)
Build a 2D straight skeleton partition for a Geometry.
auto straightSkeletonPartition (const MultiPolygon & geom, bool autoOrientation=true)
Build a 2D straight skeleton partition for a MultiPolygon .
auto straightSkeletonPartition (const Polygon & geom, bool autoOrientation=true)
Build a 2D straight skeleton partition for a Polygon .
auto tesselate (const Geometry & g, NoValidityCheck)
auto tesselate (const Geometry & g)
void translate (Geometry & g, const Kernel::Vector_3 & v)
translate a geometry from a given vector
void translate (Geometry & g, const Kernel::Vector_2 & v)
translate a geometry from a given vector
void translate (Geometry & g, const Kernel::FT & dx, const Kernel::FT & dy, const Kernel::FT & dz)
translate a geometry from a given vector
void translate (Geometry & g, const double & dx, const double & dy, const double & dz)
translate a geometry from double-coordinates
auto union3D (const Geometry & ga, const Geometry & gb, NoValidityCheck)
auto union3D (const Geometry & ga, const Geometry & gb)
void union_ (const detail::GeometrySet< Dim > & a, const detail::GeometrySet< Dim > & b, detail::GeometrySet< Dim > &)
void union_ (const detail::PrimitiveHandle< Dim > & a, const detail::PrimitiveHandle< Dim > & b, detail::GeometrySet< Dim > &)
auto visibility (const Geometry & polygon, const Geometry & point)
build the visibility polygon of a Point inside aPolygon __
auto visibility (const Geometry & polygon, const Geometry & point, NoValidityCheck)
build the visibility polygon of a Point inside aPolygon __
auto visibility (const Geometry & polygon, const Geometry & pointA, const Geometry & pointB)
build the visibility polygon of the segment [pointA ; pointB] on a Polygon __
auto visibility (const Geometry & polygon, const Geometry & pointA, const Geometry & pointB, NoValidityCheck)
build the visibility polygon of a Point inside aPolygon __
auto volume (const Solid & g, NoValidityCheck)
auto volume (const Geometry & g)
auto weightedCentroid (const Geometry & g, bool enable3DComputation=false)
auto weightedCentroid (const Triangle & g, bool enable3DComputation=false)
auto weightedCentroid (const Point & a, const Point & b, const Point & c, bool enable3DComputation=false)
auto weightedCentroid (const LineString & g, bool enable3DComputation=false)
auto weightedCentroid (const Polygon & g, bool enable3DComputation=false)
auto weightedCentroid (const GeometryCollection & g, bool enable3DComputation=false)
auto weightedCentroid (const TriangulatedSurface & g, bool enable3DComputation=false)
auto weightedCentroid (const PolyhedralSurface & g, bool enable3DComputation=false)
auto weightedCentroid (const Solid & g, bool enable3DComputation=false)

Public Types Documentation#

enum PartitionAlgorithm#

enum algorithm::PartitionAlgorithm {
    y_monotone,
    approx_convex,
    greene_approx_convex,
    optimal_convex
};

Partition algorithm available

Since:

1.4.2


Public Functions Documentation#

function alphaShapes#

auto algorithm::alphaShapes (
    const Geometry & g,
    double alpha=1,
    bool allow_holes=false
) 

Compute the 2D alpha shapes for a geometry https://doc.cgal.org/latest/Alpha_shapes_2/index.html#Chapter_2D_Alpha_Shapes

Since:

1.4.1


function alphaWrapping3D#

auto algorithm::alphaWrapping3D (
    const Geometry & geom,
    size_t relativeAlpha,
    size_t relativeOffset
) 

end of private section

Computes the 3D alpha wrapping of a geometry https://doc.cgal.org/latest/Alpha_wrap_3/index.html

Since:

2.1

Parameters:

  • geom input geometry
  • relativeAlpha This parameter is used to determine which features will appear in the output A small relativeAlpha will produce an output less complex but less faithful to the input
  • relativeOffset This parameter controls the tightness of the result A large relativeOffset parameter will tend to better preserve sharp features as projection If this parameter is equal to 0, it is computed from the alpha parameter

Returns:

A PolyhedralSurface representing the 3D alpha wrapping of the geometry


function approximate#

approximate an Offset_polygon_2 (filter null segments)

auto algorithm::approximate (
    const Offset_polygon_2 & polygon,
    const int & n=0
) 

function approximate#

approximate an Offset

auto algorithm::approximate (
    const Offset_polygon_with_holes_2 & polygon,
    const int & n=0
) 

function approximateMedialAxis#

build an approximate medial axis for a Polygon __

auto algorithm::approximateMedialAxis (
    const Geometry & geom
) 

Parameters:

  • geom input geometry

Precondition:

geom is a valid geometry

Exception:


function area#

Compute the 2D area for a Geometry.

auto algorithm::area (
    const Geometry & g,
    NoValidityCheck
) 

Warning:

Z component is ignored, there is no 2D projection for 3D geometries

Precondition:

g is a valid geometry

Warning:

No actual validity check is done


function area#

Compute the 2D area for a Geometry.

auto algorithm::area (
    const Geometry & g
) 

Warning:

Z component is ignored, there is no 2D projection for 3D geometries

Precondition:

g is a valid geometry


function area#

auto algorithm::area (
    const Triangle & g
) 

Returns Compute the 2D area for a Triangle


function area#

auto algorithm::area (
    const Polygon & g
) 

Returns Compute the 2D area for a Polygon


function area#

auto algorithm::area (
    const GeometryCollection & g
) 

Returns the 2D area for a GeometryCollection


function area#

auto algorithm::area (
    const TriangulatedSurface & g
) 

Returns the 2D area for a TriangulatedSurface


function area#

auto algorithm::area (
    const PolyhedralSurface & g
) 

Returns the 2D area for a TriangulatedSurface


function area3D#

auto algorithm::area3D (
    const Geometry & g,
    NoValidityCheck
) 

Returns 3D area for a Geometry

Warning:

Solid area is set to 0 (might be defined as the area of the surface)

Precondition:

g is a valid geometry

Warning:

No actual validity check is done


function area3D#

auto algorithm::area3D (
    const Geometry & g
) 

Returns 3D area for a Geometry

Warning:

Solid area is set to 0 (might be defined as the area of the surface)

Precondition:

g is a valid geometry


function area3D#

auto algorithm::area3D (
    const Polygon & g
) 

Returns 3D area for a Polygon


function area3D#

auto algorithm::area3D (
    const Triangle & g
) 

Returns the 3D area for a Triangle


function area3D#

auto algorithm::area3D (
    const GeometryCollection & g
) 

Returns the 3D area for a MultiPolygon


function area3D#

auto algorithm::area3D (
    const PolyhedralSurface & g
) 

Returns the 3D area for a PolyhedralSurface


function area3D#

auto algorithm::area3D (
    const TriangulatedSurface & g
) 

Returns the 3D area for a TriangulatedSurface


function boundingCircle#

auto algorithm::boundingCircle (
    const Geometry & geom
) 

function boundingSphere#

auto algorithm::boundingSphere (
    const Geometry & geom
) 

function centroid#

Returns the 2D centroid for a Geometry.

auto algorithm::centroid (
    const Geometry & g
) 

The result is the weighted centroid of a geometry. The implementation follows PostGIS one (https://postgis.net/docs/ST_Centroid.html). The weigth is computed in the XY space.

Warning:

Z component is ignored, geometries must be valid when projected in the XY plane. Vertical geometries will generate an error.

Precondition:

g is a valid geometry in 2D


function centroid3D#

Returns the 3D centroid for a Geometry.

auto algorithm::centroid3D (
    const Geometry & g
) 

The result is the weighted centroid of a geometry. The implementation follows PostGIS one (https://postgis.net/docs/ST_Centroid.html). The weigth is computed in the 3D space.

Precondition:

g is a valid geometry


function circleToPolygon#

helper to create a polygon from a circle

auto algorithm::circleToPolygon (
    const Kernel::Circle_2 & circle
) 

function collect#

auto algorithm::collect (
    const Geometry & ga,
    const Geometry & gb
) 

Returns an aggregate of ga and gb


function collect#

template<typename GeometryIterator>
std::unique_ptr< Geometry > algorithm::collect (
    GeometryIterator begin,
    GeometryIterator end
) 

Returns an aggregate of a list of geometries


function collectionExtractPolygons#

auto algorithm::collectionExtractPolygons (
    std::unique_ptr< Geometry > coll
) 

Given a geometry collection returns a MultiPolygon from triangles, polygons, polyhedral and polygons

Warning:

Ownership is taken from the parameter


function collectionHomogenize#

auto algorithm::collectionHomogenize (
    std::unique_ptr< Geometry > coll
) 

Given a geometry collection, returns the "simplest" representation of the contents. Singletons will be returned as singletons. Collections that are homogeneous will be returned as the appropriate multi-type.

Warning:

Ownership is taken from the parameter


function collectionToMulti#

auto algorithm::collectionToMulti (
    std::unique_ptr< Geometry > coll
) 

Given a geometry collection of triangles, TINs and polygons returns a MultiPolygon

Warning:

Ownership is taken from the parameter


function convexHull#

auto algorithm::convexHull (
    const Geometry & g
) 

Compute the 2D convex hull for a geometry


function convexHull3D#

auto algorithm::convexHull3D (
    const Geometry & g
) 

Compute the 3D convex hull for a geometry

Todo

improve to handle collinear points and coplanar points


function covers#

auto algorithm::covers (
    const Geometry & ga,
    const Geometry & gb
) 

Cover test on 2D geometries. Checks if gA covers gB. Force projection to z=0 if needed


function covers#

template<int Dim>
bool algorithm::covers (
    const detail::GeometrySet < Dim > & a,
    const detail::GeometrySet < Dim > & b
) 

function covers#

template<int Dim>
bool algorithm::covers (
    const detail::PrimitiveHandle < Dim > & a,
    const detail::PrimitiveHandle < Dim > & b
) 

function covers3D#

auto algorithm::covers3D (
    const Geometry & ga,
    const Geometry & gb
) 

Cover test on 3D geometries. Checks if gA covers gB. Assume z = 0 if needed


function difference#

auto algorithm::difference (
    const Geometry & ga,
    const Geometry & gb,
    NoValidityCheck
) 

Diffrence on 2D geometries. No validity check variant

Precondition:

ga and gb are valid geometries

Warning:

No actual validity check is done.


function difference#

auto algorithm::difference (
    const Geometry & ga,
    const Geometry & gb
) 

Difference on 2D geometries.

Precondition:

ga and gb are valid geometries


function difference3D#

auto algorithm::difference3D (
    const Geometry & ga,
    const Geometry & gb,
    NoValidityCheck
) 

Difference on 3D geometries. Assume z = 0 if needed

Precondition:

ga and gb are valid geometries

Warning:

@ No actual validity check is done


function difference3D#

auto algorithm::difference3D (
    const Geometry & ga,
    const Geometry & gb
) 

Difference on 3D geometries. Assume z = 0 if needed

Precondition:

ga and gb are valid geometries


function distance#

auto algorithm::distance (
    const Geometry & gA,
    const Geometry & gB,
    NoValidityCheck
) 

Compute the distance between two Geometries

Precondition:

gA is a valid geometry

Precondition:

gB is a valid geometry

Warning:

No actual validity check is done


function distance#

auto algorithm::distance (
    const Geometry & gA,
    const Geometry & gB
) 

Compute the distance between two Geometries

Precondition:

gA is a valid geometry

Precondition:

gB is a valid geometry


function distance3D#

auto algorithm::distance3D (
    const Geometry & gA,
    const Geometry & gB,
    NoValidityCheck
) 

Compute distance between two 3D Geometries

Precondition:

gA is a valid geometry

Precondition:

gB is a valid geometry

Warning:

No actual validity check is done


function distance3D#

auto algorithm::distance3D (
    const Geometry & gA,
    const Geometry & gB
) 

dispatch distance between two Geometries

Todo

complete with solid Compute distance between two 3D Geometries

Precondition:

gA is a valid geometry

Precondition:

gB is a valid geometry


function distanceGeometryCollectionToGeometry#

auto algorithm::distanceGeometryCollectionToGeometry (
    const Geometry & gA,
    const Geometry & gB
) 

dispatch distance from a collection of geometry (gA) to a Geometry (gB)


function distanceGeometryCollectionToGeometry3D#

auto algorithm::distanceGeometryCollectionToGeometry3D (
    const Geometry & gA,
    const Geometry & gB
) 

dispatch distance from a collection of geometry (gA) to a Geometry (gB)


function distanceLineStringGeometry#

auto algorithm::distanceLineStringGeometry (
    const LineString & gA,
    const Geometry & gB
) 

dispatch distance from LineString to Geometry


function distanceLineStringGeometry3D#

auto algorithm::distanceLineStringGeometry3D (
    const LineString & gA,
    const Geometry & gB
) 

dispatch distance between a LineString and a Geometry


function distanceLineStringLineString#

auto algorithm::distanceLineStringLineString (
    const LineString & gA,
    const LineString & gB
) 

distance between two LineStrings


function distanceLineStringLineString3D#

auto algorithm::distanceLineStringLineString3D (
    const LineString & gA,
    const LineString & gB
) 

distance between two LineStrings


function distanceLineStringPolygon#

auto algorithm::distanceLineStringPolygon (
    const LineString & gA,
    const Polygon & gB
) 

distance between a LineString and a Polygon


function distanceLineStringPolygon3D#

auto algorithm::distanceLineStringPolygon3D (
    const LineString & gA,
    const Polygon & gB
) 

distance between a LineString and a Polygon

Todo

same method than distancePointPolygon3D (unify if triangulate is available)


function distanceLineStringPolyhedralSurface3D#

auto algorithm::distanceLineStringPolyhedralSurface3D (
    const LineString & lineA,
    const PolyhedralSurface & polySurfaceB
) 

distance between a LineString and a PolyhedralSurface


function distanceLineStringSolid3D#

auto algorithm::distanceLineStringSolid3D (
    const LineString & gA,
    const Solid & gB
) 

distance between a LineString and a Solid


function distanceLineStringTriangle#

auto algorithm::distanceLineStringTriangle (
    const LineString & gA,
    const Triangle & gB
) 

distance between a LineString and a Triangle


function distanceLineStringTriangle3D#

auto algorithm::distanceLineStringTriangle3D (
    const LineString & gA,
    const Triangle & gB
) 

distance between a LineString and a Triangle


function distanceLineStringTriangulatedSurface3D#

auto algorithm::distanceLineStringTriangulatedSurface3D (
    const LineString & lineA,
    const TriangulatedSurface & triangulatedSurfaceB
) 

distance between a LineString and a TriangulatedSurface


function distancePointGeometry#

auto algorithm::distancePointGeometry (
    const Point & gA,
    const Geometry & gB
) 

dispatch distance from Point to Geometry


function distancePointGeometry3D#

auto algorithm::distancePointGeometry3D (
    const Point & gA,
    const Geometry & gB
) 

dispatch distance from Point to Geometry


function distancePointLineString#

auto algorithm::distancePointLineString (
    const Point & gA,
    const LineString & gB
) 

distance between a Point and a LineString


function distancePointLineString3D#

auto algorithm::distancePointLineString3D (
    const Point & gA,
    const LineString & gB
) 

distance between a Point and a LineString


function distancePointPoint#

auto algorithm::distancePointPoint (
    const Point & gA,
    const Point & gB
) 

distance between two Points


function distancePointPoint3D#

auto algorithm::distancePointPoint3D (
    const Point & gA,
    const Point & gB
) 

distance between two Points


function distancePointPolygon#

auto algorithm::distancePointPolygon (
    const Point & gA,
    const Polygon & gB
) 

distance between a Point and a Polygon


function distancePointPolygon3D#

auto algorithm::distancePointPolygon3D (
    const Point & gA,
    const Polygon & gB
) 

distance between a Point and a Triangle


function distancePointPolyhedralSurface3D#

auto algorithm::distancePointPolyhedralSurface3D (
    const Point & pointA,
    const PolyhedralSurface & polySurfaceB
) 

distance between a Point and a PolyhedralSurface


function distancePointSegment#

auto algorithm::distancePointSegment (
    const Point & p,
    const Point & a,
    const Point & b
) 

function distancePointSegment3D#

auto algorithm::distancePointSegment3D (
    const Point & p,
    const Point & a,
    const Point & b
) 

function distancePointSolid3D#

auto algorithm::distancePointSolid3D (
    const Point & gA,
    const Solid & gB
) 

distance between a Point and a Solid


function distancePointTriangle#

auto algorithm::distancePointTriangle (
    const Point & gA,
    const Triangle & gB
) 

distance between a Point and a Triangle


function distancePointTriangle3D#

auto algorithm::distancePointTriangle3D (
    const Point & gA,
    const Triangle & gB
) 

distance between a Point and a Triangle


function distancePointTriangle3D#

auto algorithm::distancePointTriangle3D (
    const Point & p_,
    const Point & a_,
    const Point & b_,
    const Point & c_
) 

function distancePointTriangulatedSurface3D#

auto algorithm::distancePointTriangulatedSurface3D (
    const Point & pointA,
    const TriangulatedSurface & triangulatedSurfaceB
) 

distance between a Point and a TriangulatedSurface


function distancePolygonGeometry#

auto algorithm::distancePolygonGeometry (
    const Polygon & gA,
    const Geometry & gB
) 

dispatch distance from Polygon to Geometry


function distancePolygonGeometry3D#

auto algorithm::distancePolygonGeometry3D (
    const Polygon & gA,
    const Geometry & gB
) 

dispatch distance between a Polygon and a Geometry


function distancePolygonPolygon#

auto algorithm::distancePolygonPolygon (
    const Polygon & gA,
    const Polygon & gB
) 

distance between two Polygons


function distancePolygonTriangle#

auto algorithm::distancePolygonTriangle (
    const Polygon & gA,
    const Triangle & gB
) 

distance between a Polygon and a Triangle


function distancePolyhedralSurfaceGeometry3D#

auto algorithm::distancePolyhedralSurfaceGeometry3D (
    const PolyhedralSurface & polySurfaceA,
    const Geometry & geomB
) 

dispatch distance between a PolyhedralSurface and a Geometry


function distanceSegmentSegment#

auto algorithm::distanceSegmentSegment (
    const Point & a,
    const Point & b,
    const Point & c,
    const Point & d
) 

function distanceSegmentSegment3D#

auto algorithm::distanceSegmentSegment3D (
    const Point & a,
    const Point & b,
    const Point & c,
    const Point & d
) 

function distanceSegmentTriangle3D#

auto algorithm::distanceSegmentTriangle3D (
    const Point & sA_,
    const Point & sB_,
    const Point & tA_,
    const Point & tB_,
    const Point & tC_
) 

function distanceSolidGeometry3D#

auto algorithm::distanceSolidGeometry3D (
    const Solid & gA,
    const Geometry & gB
) 

dispatch distance between a Solid and a Geometry


function distanceSolidSolid3D#

auto algorithm::distanceSolidSolid3D (
    const Solid & gA,
    const Solid & gB
) 

distance between two Solids


function distanceTriangleGeometry#

auto algorithm::distanceTriangleGeometry (
    const Triangle & gA,
    const Geometry & gB
) 

dispatch distance from a Triangle to a Geometry


function distanceTriangleGeometry3D#

auto algorithm::distanceTriangleGeometry3D (
    const Triangle & gA,
    const Geometry & gB
) 

dispatch distance between a Triangle and a Geometry


function distanceTrianglePolyhedralSurface3D#

auto algorithm::distanceTrianglePolyhedralSurface3D (
    const Triangle & triangleA,
    const PolyhedralSurface & polySurfaceB
) 

distance between a Triangle and a PolyhedralSurface


function distanceTriangleSolid3D#

auto algorithm::distanceTriangleSolid3D (
    const Triangle & gA,
    const Solid & gB
) 

distance between a Triangle and a Solid


function distanceTriangleTriangle3D#

auto algorithm::distanceTriangleTriangle3D (
    const Triangle & gA,
    const Triangle & gB
) 

distance between two Triangles


function distanceTriangleTriangulatedSurface3D#

auto algorithm::distanceTriangleTriangulatedSurface3D (
    const Triangle & triangleA,
    const TriangulatedSurface & triangulatedSurfaceB
) 

function distanceTriangulatedSurfaceGeometry3D#

auto algorithm::distanceTriangulatedSurfaceGeometry3D (
    const TriangulatedSurface & triangulatedSurfaceA,
    const Geometry & geomB
) 

dispatch distance between a TriangulatedSurface and a Geometry


function extrude#

Returns a Geometry equal to the specified Geometry, extruded by the specified displacement vector.

auto algorithm::extrude (
    const Geometry & g,
    const Kernel::Vector_3 & v
) 

Parameters:

  • g The specified Geometry.
  • v The specified displacement vector.

Returns:

A Geometry equal to g extruded by the displacement vector v.

Precondition:

g must be a valid geometry.

Note:

If g is such that g.isMeasured() is true, then, since there is no common expectation of the values of the measures on the returned Geometry, all measures from the result are removed.

Todo

Improve extrude for 3D surfaces - Extrude only faces whose scalar_product(v,normal) > 0 and use Polyhedron union to get output geometries with a clean topology.


function extrude#

auto algorithm::extrude (
    const Geometry & g,
    const Kernel::FT & dx,
    const Kernel::FT & dy,
    const Kernel::FT & dz,
    NoValidityCheck
) 

function extrude#

Returns a Geometry equal to the specified Geometry, extruded by the specified displacement.

auto algorithm::extrude (
    const Geometry & g,
    const Kernel::FT & dx,
    const Kernel::FT & dy,
    const Kernel::FT & dz
) 

Parameters:

  • g The specified Geometry.
  • dx The component of the specified displacement in the x-direction.
  • dy The component of the specified displacement in the y-direction.
  • dz The component of the specified displacement in the z-direction.

Returns:

A Geometry equal to g extruded by the displacement vector {dx, dy, dz}.

Precondition:

g must be a valid geometry.

Precondition:

dx, dy and dz must all be finite.

Note:

If g is such that g.isMeasured() is true, then, since there is no common expectation of the values of the measures on the returned Geometry, all measures from the result are removed.


function extrude#

Returns a Geometry equal to the specified Geometry, extruded by the specified displacement.

SFCGAL_API auto algorithm::extrude (
    const Geometry & g,
    const double & dx,
    const double & dy,
    const double & dz
) 

Parameters:

  • g The specified Geometry.
  • dx The component of the specified displacement in the x-direction.
  • dy The component of the specified displacement in the y-direction.
  • dz The component of the specified displacement in the z-direction.

Returns:

A Geometry equal to g extruded by the displacement vector {dx, dy, dz}.

Precondition:

g must be a valid geometry.

Precondition:

dx, dy and dz must all be finite.

Note:

If g is such that g.isMeasured() is true, then, since there is no common expectation of the values of the measures on the returned Geometry, all measures from the result are removed.

Warning:

No actual validity check is conducted.


function extrude#

SFCGAL_API auto algorithm::extrude (
    const Polygon & g,
    const double & height
) 

function extrude#

Returns a Geometry equal to the specified Geometry, extruded by the specified displacement.

SFCGAL_API std::unique_ptr< Geometry > algorithm::extrude (
    const Geometry & g,
    const Kernel::FT & dx,
    const Kernel::FT & dy,
    const Kernel::FT & dz,
    NoValidityCheck & nvc
) 

Parameters:

  • g The specified Geometry.
  • dx The component of the specified displacement in the x-direction.
  • dy The component of the specified displacement in the y-direction.
  • dz The component of the specified displacement in the z-direction.
  • nvc A NoValidityCheck object.

Returns:

A Geometry equal to g extruded by the displacement vector {dx, dy, dz}.

Precondition:

g must be a valid geometry.

Precondition:

dx, dy and dz must all be finite.

Note:

If g is such that g.isMeasured() is true, then, since there is no common expectation of the values of the measures on the returned Geometry, all measures from the result are removed.

Warning:

No actual validity check is conducted.


function extrudeStraightSkeleton#

auto algorithm::extrudeStraightSkeleton (
    const Polygon & geom,
    double height
) 

function extrudeStraightSkeleton#

auto algorithm::extrudeStraightSkeleton (
    const Geometry & geom,
    double height
) 

function extrudeStraightSkeleton#

auto algorithm::extrudeStraightSkeleton (
    const Geometry & geom,
    double building_height,
    double roof_height
) 

function extrudedStraightSkeleton#

build a 3D straight skeleton extruded for a Polygon __

SFCGAL_API auto algorithm::extrudedStraightSkeleton (
    const Polygon & geom,
    double height
) 

Exception:


function force2D#

force a geometry to be 2D (project on O,x,y)

void algorithm::force2D (
    Geometry & g
) 

Warning:

ignore empty geometries


function force3D#

force a 2D geometry to be 3D (replace undefined Z by defaultZ, existing Z values remains unchanged)

void algorithm::force3D (
    Geometry & g,
    const Kernel::FT & defaultZ=0
) 

Warning:

ignore empty geometries


function forceMeasured#

force a 2D or M geometry to be M (replace undefined M by defaultM, existing M values remains unchanged)

void algorithm::forceMeasured (
    Geometry & geometry,
    const double & defaultM=0
) 

Warning:

ignore empty geometries


function hasConsistentOrientation3D#

auto algorithm::hasConsistentOrientation3D (
    const TriangulatedSurface & g
) 

Test if a Geometry has a consistent orientation


function hasConsistentOrientation3D#

auto algorithm::hasConsistentOrientation3D (
    const PolyhedralSurface & g
) 

Test if a PolyhedralSurface has a consistent orientation


function hasPlane3D#

Test if a 3D plane can be extracted from a Polygon .

template<typename Kernel>
bool algorithm::hasPlane3D (
    const Polygon & polygon,
    CGAL::Point_3< Kernel > & a,
    CGAL::Point_3< Kernel > & b,
    CGAL::Point_3< Kernel > & c
) 

function hasPlane3D#

template<typename Kernel>
bool algorithm::hasPlane3D (
    const Polygon & polygon
) 

Test if a 3D plane can be extracted from a Polygon


function intersection#

auto algorithm::intersection (
    const CGAL::Triangle_3< Kernel > & a,
    const CGAL::Triangle_3< Kernel > & b
) 

function intersection#

template<int Dim>
void algorithm::intersection (
    const GeometrySet < Dim > & a,
    const GeometrySet < Dim > & b,
    GeometrySet < Dim > & output
) 

function intersection#

auto algorithm::intersection (
    const Geometry & ga,
    const Geometry & gb,
    NoValidityCheck
) 

Intersection on 2D geometries. No validity check variant

Precondition:

ga and gb are valid geometries

Warning:

No actual validity check is done.


function intersection#

auto algorithm::intersection (
    const Geometry & ga,
    const Geometry & gb
) 

Intersection on 2D geometries.

Precondition:

ga and gb are valid geometries


function intersection3D#

auto algorithm::intersection3D (
    const Geometry & ga,
    const Geometry & gb,
    NoValidityCheck
) 

Intersection on 3D geometries. Assume z = 0 if needed

Precondition:

ga and gb are valid geometries

Warning:

@ No actual validity check is done


function intersection3D#

auto algorithm::intersection3D (
    const Geometry & ga,
    const Geometry & gb
) 

Intersection on 3D geometries. Assume z = 0 if needed

Precondition:

ga and gb are valid geometries


function intersects#

template<int Dim>
auto algorithm::intersects (
    const detail::PrimitiveHandle < Dim > & a,
    const detail::PrimitiveHandle < Dim > & b
) 

Intersection test on a PrimitiveHandle


function intersects#

template<int Dim>
auto algorithm::intersects (
    const detail::GeometrySet < Dim > & a,
    const detail::GeometrySet < Dim > & b
) 

Intersection test on GeometrySet


function intersects#

auto algorithm::intersects (
    const Geometry & ga,
    const Geometry & gb
) 

Robust intersection test on 2D geometries. Force projection to z=0 if needed

Precondition:

ga and gb are valid geometries


function intersects#

auto algorithm::intersects (
    const Geometry & ga,
    const Geometry & gb,
    NoValidityCheck
) 

Intersection test on 2D geometries. Force projection to z=0 if needed

Precondition:

ga and gb are valid geometries

Warning:

the validity is assumed, no actual check is done


function intersects#

template<int Dim>
bool algorithm::intersects (
    const detail::GeometrySet < Dim > & a,
    const detail::GeometrySet < Dim > & b
) 

Intersection test on GeometrySet


function intersects#

template<int Dim>
bool algorithm::intersects (
    const detail::PrimitiveHandle < Dim > & a,
    const detail::PrimitiveHandle < Dim > & b
) 

Intersection test on a PrimitiveHandle


function intersects3D#

auto algorithm::intersects3D (
    const Geometry & ga,
    const Geometry & gb
) 

Robust intersection test on 3D geometries. Assume z = 0 if needed

Precondition:

ga and gb are valid geometries


function intersects3D#

auto algorithm::intersects3D (
    const Geometry & ga,
    const Geometry & gb,
    NoValidityCheck
) 

Intersection test on 3D geometries. Assume z = 0 if needed

Precondition:

ga and gb are valid geometries

Warning:

the validity is assumed, no actual check is done


function isClosed#

auto algorithm::isClosed (
    const SurfaceGraph & graph
) 

test if a surface is closed, the graph should be build beforehand

Note:

the surface may not be connected, eg. two spheres will yield a true result


function isConnected#

auto algorithm::isConnected (
    const SurfaceGraph & graph
) 

test if a surface is connected, the graph should be build beforehand


function isCounterClockWiseOriented#

auto algorithm::isCounterClockWiseOriented (
    const LineString & ls
) 

Test if a 2D surface is oriented counter clockwise


function isCounterClockWiseOriented#

auto algorithm::isCounterClockWiseOriented (
    const Triangle & tri
) 

Test if a 2D surface is oriented counter clockwise


function isCounterClockWiseOriented#

auto algorithm::isCounterClockWiseOriented (
    const Polygon & poly
) 

Test if a 2D surface is oriented counter clockwise


function isPlane3D#

template<typename Kernel>
bool algorithm::isPlane3D (
    const Geometry & geom,
    const double & toleranceAbs
) 

Test if all points of a geometry lie in the same plane


function isSimple#

auto algorithm::isSimple (
    const Geometry & g,
    const double & toleranceAbs=1e-9
) 

Check simplicity of a geometry


function isValid#

auto algorithm::isValid (
    const Point & p
) 

Note:

empty geometries are valid, but the test is only performed in the interface function in individual functions for implementation, an assertion !empty is present for this reason


function isValid#

auto algorithm::isValid (
    const LineString & l,
    const double & toleranceAbs
) 

function isValid#

auto algorithm::isValid (
    const Polygon & p,
    const double & toleranceAbs
) 

function isValid#

auto algorithm::isValid (
    const Triangle & t,
    const double & toleranceAbs
) 

function isValid#

auto algorithm::isValid (
    const MultiLineString & ml,
    const double & toleranceAbs
) 

function isValid#

auto algorithm::isValid (
    const MultiPolygon & mp,
    const double & toleranceAbs
) 

function isValid#

auto algorithm::isValid (
    const GeometryCollection & gc,
    const double & toleranceAbs
) 

function isValid#

auto algorithm::isValid (
    const TriangulatedSurface & tin,
    const SurfaceGraph & graph,
    const double & toleranceAbs
) 

function isValid#

auto algorithm::isValid (
    const TriangulatedSurface & tin,
    const double & toleranceAbs
) 

function isValid#

auto algorithm::isValid (
    const PolyhedralSurface & s,
    const SurfaceGraph & graph,
    const double & toleranceAbs
) 

function isValid#

auto algorithm::isValid (
    const PolyhedralSurface & s,
    const double & toleranceAbs
) 

function isValid#

auto algorithm::isValid (
    const Solid & solid,
    const double & toleranceAbs
) 

function isValid#

auto algorithm::isValid (
    const MultiSolid & ms,
    const double & toleranceAbs
) 

function isValid#

Check validity of a geometry.

auto algorithm::isValid (
    const Geometry & g,
    const double & toleranceAbs
) 

function length#

Compute the 2D length for a LineString .

auto algorithm::length (
    const LineString & g
) 

function length#

Compute the 2D length for a GeometryCollection .

auto algorithm::length (
    const GeometryCollection & g
) 

function length#

Compute the 2D length for a Geometry (0 for incompatible types)

auto algorithm::length (
    const Geometry & g
) 

function length3D#

Compute the 3D length for a LineString .

auto algorithm::length3D (
    const LineString & g
) 

function length3D#

Compute the 3D length for a GeometryCollection .

auto algorithm::length3D (
    const GeometryCollection & g
) 

function length3D#

Compute the 2D length for a geometry.

auto algorithm::length3D (
    const Geometry & g
) 

Returns:

the length of the Geometry, 0 for incompatible types


function lineSubstring#

Retrieve a substring of a specified LineString , between the specified fractional distances from the start of the specifiedLineString .

auto algorithm::lineSubstring (
    const LineString & ls,
    double start,
    double end
) 

Parameters:

  • ls The specified LineString.
  • start The fraction along the specified LineString defining the start of the desired substring.
  • end The fraction along the specified LineString defining the end of the desired substring.

Note:

Negative values of and/or will be interpreted as a fractional distance taken from the end of the specified LineString. +/-0 will always be interpreted as the start of .

Note:

For open lines, a negative length range will result in a line substring terminating at the specified points, but with an orientation reversed relative to . For closed lines the a negative range corresponds to the complentary section of with an orientation equal to that of it.

Returns:

The specified line substring.

Exception:

  • If either or have an absolute value greater than 1.

function makeConsistentOrientation3D#

void algorithm::makeConsistentOrientation3D (
    TriangulatedSurface & g
) 

Try to make consistent orientation in a TriangulatedSurface


function makeValidOrientation#

void algorithm::makeValidOrientation (
    CGAL::Polygon_2< Kernel > & polygon
) 

Make valid 2D orientation


function makeValidOrientation#

void algorithm::makeValidOrientation (
    CGAL::Polygon_with_holes_2< Kernel > & polygon
) 

Make valid 2D orientation


function makeValidOrientation#

void algorithm::makeValidOrientation (
    Polygon & polygon
) 

Make valid 2D orientation


function minkowskiSum#

2D minkowski sum (p+q)

auto algorithm::minkowskiSum (
    const Geometry & gA,
    const Polygon & gB,
    NoValidityCheck
) 

Warning:

If gA is a polygon, its orientation is taken into account. A "reversed" polygon (with a clockwise-oriented exterior ring) will involve a minkowski difference rather than a sum.

Todo

missing cases (union)

Precondition:

gA and gB are valid geometries

Warning:

@ No actual validity check is done.


function minkowskiSum#

2D minkowski sum (p+q)

auto algorithm::minkowskiSum (
    const Geometry & gA,
    const Polygon & gB
) 

Warning:

If gA is a polygon, its orientation is taken into account. A "reversed" polygon (with a clockwise-oriented exterior ring) will involve a minkowski difference rather than a sum.

Todo

missing cases (union)

Precondition:

gA and gB are valid geometries


function minkowskiSum3D#

3D Minkowski sum (p+q)

auto algorithm::minkowskiSum3D (
    const Geometry & gA,
    const Geometry & gB,
    NoValidityCheck
) 

Precondition:

gA and gB are valid 3D geometries

Warning:

No actual validity check is done.


function minkowskiSum3D#

3D Minkowski sum (p+q)

auto algorithm::minkowskiSum3D (
    const Geometry & gA,
    const Geometry & gB
) 

Precondition:

gA and gB are valid 3D geometries


function normal3D#

template<typename Kernel>
CGAL::Vector_3< Kernel > algorithm::normal3D (
    const CGAL::Point_3< Kernel > & a,
    const CGAL::Point_3< Kernel > & b,
    const CGAL::Point_3< Kernel > & c
) 

Returns the 3D normal to 3 consecutive points.


function normal3D#

template<typename Kernel>
CGAL::Vector_3< Kernel > algorithm::normal3D (
    const LineString & ls,
    bool exact=true
) 

Returns the 3D normal to a ring (supposed to be planar and closed).

Warning:

exact allows to avoid double rounding at the end of the computation


function normal3D#

template<typename Kernel>
CGAL::Vector_3< Kernel > algorithm::normal3D (
    const Polygon & polygon,
    bool exact=true
) 

Returns the 3D normal to a polygon (supposed to be planar).

Warning:

exact allows to avoid double rounding at the end of the computation


function offset#

dispatch a geometry

void algorithm::offset (
    const Geometry & g,
    const double & radius,
    Offset_polygon_set_2 & polygonSet
) 

function offset#

offset for a Point __

void algorithm::offset (
    const Point & g,
    const double & radius,
    Offset_polygon_set_2 & polygonSet
) 

build Point offset


function offset#

offset for a LineString __

void algorithm::offset (
    const LineString & g,
    const double & radius,
    Offset_polygon_set_2 & polygonSet
) 

build LineString offset


function offset#

offset for a Polygon __

void algorithm::offset (
    const Polygon & g,
    const double & radius,
    Offset_polygon_set_2 & polygonSet
) 

function offset#

[experimental]compute polygon offset

auto algorithm::offset (
    const Geometry & g,
    const double & r,
    NoValidityCheck
) 

Warning:

test in order to compare with minkowski sum

Precondition:

g is a valid Geometry

Warning:

No actual validity check is done.


function offset#

[experimental]compute polygon offset

auto algorithm::offset (
    const Geometry & g,
    const double & r
) 

Warning:

test in order to compare with minkowski sum

Precondition:

g is a valid Geometry


function offsetCollection#

offset for MultiPoint ,MultiLineString ,MultiPolygon ,TriangulatedSurface ,PolyhedralSurface __

void algorithm::offsetCollection (
    const Geometry & g,
    const double & radius,
    Offset_polygon_set_2 & polygonSet
) 

function optimal_alpha_shapes#

auto algorithm::optimal_alpha_shapes (
    const Geometry & g,
    bool allow_holes=false,
    size_t nb_components=1
) 

Compute the optimal 2D alpha shapes for a geometry https://doc.cgal.org/latest/Alpha_shapes_2/index.html#Chapter_2D_Alpha_Shapes

Since:

1.4.1


function partition_2#

auto algorithm::partition_2 (
    const Geometry & g,
    PartitionAlgorithm alg,
    NoValidityCheck
) 

Compute the partition of a 2D polygon https://doc.cgal.org/latest/Partition_2/index.html#Chapter_2D_Polygon_Partitioning

Precondition:

g is a valid geometry

Warning:

No actual validity check is done


function partition_2#

auto algorithm::partition_2 (
    const Geometry & g,
    PartitionAlgorithm alg=y_monotone
) 

Compute the partition of a 2D polygon https://doc.cgal.org/latest/Partition_2/index.html#Chapter_2D_Polygon_Partitioning

Since:

1.4.2


function plane3D#

template<typename Kernel>
void algorithm::plane3D (
    const Polygon & polygon,
    CGAL::Point_3< Kernel > & a,
    CGAL::Point_3< Kernel > & b,
    CGAL::Point_3< Kernel > & c
) 

Get 3 non collinear points from a Polygon


function plane3D#

template<typename Kernel>
CGAL::Plane_3< Kernel > algorithm::plane3D (
    const Polygon & polygon
) 

Returns the oriented 3D plane of a polygon (supposed to be planar). May return degenerate plane.


function plane3D#

template<typename Kernel>
CGAL::Plane_3< Kernel > algorithm::plane3D (
    const Polygon & polygon,
    const Plane3DInexactUnsafe &
) 

Returns the oriented 3D plane of a polygon (supposed to be planar) - inexact version.

Warning:

Will divide by zero if polygon is degenerate.

Warning:

result is rounded to double (avoid huge expression tree).


function plane3D#

template<typename Kernel>
CGAL::Plane_3< Kernel > algorithm::plane3D (
    const Polygon & polygon,
    bool exact
) 

Returns the oriented 3D plane of a polygon (supposed to be planar). This is legacy code for SFCGAL users and should be deprecated.

Warning:

result is rounded to double if exact is false (avoid huge expression tree).

Warning:

Will divide by zero if polygon is degenerate. This maintains the previous behaviour.


function polygonSetToMultiPolygon#

convert Offset_polygon_set_2 to MultiPolygon __

auto algorithm::polygonSetToMultiPolygon (
    const Offset_polygon_set_2 & polygonSet,
    const int & n
) 

function propagateValidityFlag#

void algorithm::propagateValidityFlag (
    Geometry & g,
    bool valid
) 

Sets the geometry flag on a geometry and propagate to every internal geometries


function rotate#

Rotate a geometry in 2D around the origin (0,0)

void algorithm::rotate (
    Geometry & g,
    const Kernel::FT & angle
) 

Parameters:

  • g The geometry to rotate
  • angle Rotation angle in radians

function rotate#

Rotate a geometry in 2D around a specified point.

void algorithm::rotate (
    Geometry & g,
    const Kernel::FT & angle,
    const Point & origin
) 

Parameters:

  • g The geometry to rotate
  • angle Rotation angle in radians
  • origin Point of origin for the rotation

function rotate#

Rotate a geometry in 3D around a specified axis and origin.

void algorithm::rotate (
    Geometry & g,
    const Kernel::FT & angle,
    const Kernel::Vector_3 & axis,
    const Point & origin=Point (0, 0, 0)
) 

Parameters:

  • g The geometry to rotate
  • angle Rotation angle in radians
  • axis The axis of rotation
  • origin Point of origin for the rotation

function rotateX#

Rotate a geometry around the X axis.

void algorithm::rotateX (
    Geometry & g,
    const Kernel::FT & angle
) 

Parameters:

  • g The geometry to rotate
  • angle Rotation angle in radians

function rotateY#

Rotate a geometry around the Y axis.

void algorithm::rotateY (
    Geometry & g,
    const Kernel::FT & angle
) 

Parameters:

  • g The geometry to rotate
  • angle Rotation angle in radians

function rotateZ#

Rotate a geometry around the Z axis.

void algorithm::rotateZ (
    Geometry & g,
    const Kernel::FT & angle
) 

Parameters:

  • g The geometry to rotate
  • angle Rotation angle in radians

function scale#

void algorithm::scale (
    Geometry & g,
    double s
) 

Scale a geometry by a given factor

Parameters:

  • g input geometry
  • s scale factor

function scale#

void algorithm::scale (
    Geometry & g,
    double sx,
    double sy,
    double sz=0.0
) 

Scale a geometry by different factors for each dimension

Parameters:

  • g input geometry
  • sx scale factor for x dimension
  • sy scale factor for y dimension
  • sz scale factor for z dimension

function scale#

void algorithm::scale (
    Geometry & g,
    double sx,
    double sy,
    double sz,
    double cx,
    double cy,
    double cz
) 

Scale a geometry by different factors for each dimension around a center point

Parameters:

  • g input geometry
  • sx scale factor for x dimension
  • sy scale factor for y dimension
  • sz scale factor for z dimension
  • cx x-coordinate of the center point
  • cy y-coordinate of the center point
  • cz z-coordinate of the center point

function selfIntersects#

auto algorithm::selfIntersects (
    const LineString & l
) 

Self intersection test for 2D LineString (false if only endpoint touch)


function selfIntersects#

auto algorithm::selfIntersects (
    const PolyhedralSurface & s,
    const SurfaceGraph & g
) 

Self intersection test for 2D PolyhedralSurface (false if only point touch)


function selfIntersects#

auto algorithm::selfIntersects (
    const TriangulatedSurface & s,
    const SurfaceGraph & g
) 

Self intersection test for 2D TriangulatedSurface (false if only point touch)


function selfIntersects3D#

auto algorithm::selfIntersects3D (
    const LineString & l
) 

Self intersection test for 3D LineString (false if only endpoints touch)


function selfIntersects3D#

auto algorithm::selfIntersects3D (
    const PolyhedralSurface & s,
    const SurfaceGraph & g
) 

Self intersection test for 3D PolyhedralSurface (false if only point touch)


function selfIntersects3D#

auto algorithm::selfIntersects3D (
    const TriangulatedSurface & s,
    const SurfaceGraph & g
) 

Self intersection test for 3D TriangulatedSurface (false if only point touch)


function signedArea#

Compute the 2D signed area for a Triangle .

auto algorithm::signedArea (
    const Triangle & g
) 

function signedArea#

Compute the 2D signed area for a closed LineString .

auto algorithm::signedArea (
    const LineString & g
) 

function simplify#

Main entry point for geometry simplification.

auto algorithm::simplify (
    const Geometry & geometry,
    double threshold,
    bool preserveTopology
) 

Simplifies a geometry using the CGAL algorithm https://doc.cgal.org/latest/Polyline_simplification_2/index.html.

Parameters:

  • geometry The geometry to simplify
  • threshold The distance (in geometry unit) threshold for simplification
  • preserveTopology Whether to preserve topology during simplification

Returns:

A simplified copy of the input geometry


function squaredDistancePointTriangle3D#

auto algorithm::squaredDistancePointTriangle3D (
    const Point_3 & p,
    const Triangle_3 & abc
) 

function squaredDistanceSegmentTriangle3D#

auto algorithm::squaredDistanceSegmentTriangle3D (
    const Segment_3 & sAB,
    const Triangle_3 & tABC
) 

function squaredDistanceTriangleTriangle3D#

auto algorithm::squaredDistanceTriangleTriangle3D (
    const Triangle_3 & triangleA,
    const Triangle_3 & triangleB
) 

function straightSkeleton#

build a 2D straight skeleton for a Polygon __

auto algorithm::straightSkeleton (
    const Geometry & geom,
    bool autoOrientation,
    NoValidityCheck,
    bool innerOnly=false,
    bool outputDistanceInM=false,
    const double & toleranceAbs=EPSILON
) 

Parameters:

  • geom input geometry
  • autoOrientation check and fix polygon orientation
  • innerOnly Skip non-inner edges if requested
  • outputDistanceInM whether to output the distance to border as M
  • toleranceAbs Distance tolerance between returned points. A line must have a maximum distance of toleranceAbs.

Precondition:

geom is a valid geometry

Warning:

No actual validity check is done

Exception:


function straightSkeleton#

build a 2D straight skeleton for a Polygon __

auto algorithm::straightSkeleton (
    const Geometry & geom,
    bool autoOrientation=true,
    bool innerOnly=false,
    bool outputDistanceInM=false,
    const double & toleranceAbs=EPSILON
) 

Todo

add supports for TriangulatedSurface and PolyhedralSurface

Parameters:

  • geom input geometry
  • autoOrientation check and fix polygon orientation
  • innerOnly Skip non-inner edges if requested
  • outputDistanceInM whether to output the distance to border as M
  • toleranceAbs Distance tolerance between returned points. A line must have a maximum distance of toleranceAbs.

Precondition:

geom is a valid geometry

Exception:


function straightSkeleton#

build a 2D straight skeleton for a Polygon __

auto algorithm::straightSkeleton (
    const Polygon & geom,
    bool autoOrientation=true,
    bool innerOnly=false,
    bool outputDistanceInM=false,
    const double & toleranceAbs=EPSILON
) 

Exception:


function straightSkeleton#

build a 2D straight skeleton for a Polygon __

auto algorithm::straightSkeleton (
    const MultiPolygon & geom,
    bool autoOrientation=true,
    bool innerOnly=false,
    bool outputDistanceInM=false,
    const double & toleranceAbs=EPSILON
) 

Exception:


function straightSkeletonPartition#

Build a 2D straight skeleton partition for a Geometry.

auto algorithm::straightSkeletonPartition (
    const Geometry & geom,
    bool autoOrientation=true
) 

Parameters:

  • geom The input geometry
  • autoOrientation Check and fix polygon orientation

Returns:

A unique pointer to a MultiPolygon representing the partitioned geometry

Exception:

  • Exception If CGAL fails to create the straight skeleton

Note:

Only Triangle, Polygon, and MultiPolygon geometries are supported

This function creates a partition of the input geometry based on its straight skeleton. For unsupported geometry types, an empty MultiPolygon is returned.


function straightSkeletonPartition#

Build a 2D straight skeleton partition for a MultiPolygon .

auto algorithm::straightSkeletonPartition (
    const MultiPolygon & geom,
    bool autoOrientation=true
) 

Parameters:

  • geom The input multi-polygon
  • autoOrientation Check and fix polygon orientation

Returns:

A unique pointer to a MultiPolygon representing the partitioned multi-polygon

This function applies the straight skeleton partition to each polygon in the input multi-polygon and combines the results into a single MultiPolygon.


function straightSkeletonPartition#

Build a 2D straight skeleton partition for a Polygon .

auto algorithm::straightSkeletonPartition (
    const Polygon & geom,
    bool autoOrientation=true
) 

Parameters:

  • geom The input polygon
  • autoOrientation Check and fix polygon orientation (not used in this implementation)

Returns:

A unique pointer to a MultiPolygon representing the partitioned polygon

Exception:

  • Exception If CGAL fails to create the straight skeleton

This function creates a partition of the input polygon based on its straight skeleton. It uses CGAL's Arrangement_2 to handle the intersection of skeleton and polygon edges.


function tesselate#

auto algorithm::tesselate (
    const Geometry & g,
    NoValidityCheck
) 

Tesselate a geometry: this will triangulate surfaces (including polyhedral and solid's surfaces) and keep untouched points, lines, etc.

Precondition:

g is a valid geometry

Warning:

No actual validity check is done.


function tesselate#

auto algorithm::tesselate (
    const Geometry & g
) 

Tesselate a geometry: this will triangulate surfaces (including polyhedral and solid's surfaces) and keep untouched points, lines, etc.

Precondition:

g is a valid geometry


function translate#

translate a geometry from a given vector

void algorithm::translate (
    Geometry & g,
    const Kernel::Vector_3 & v
) 

function translate#

translate a geometry from a given vector

void algorithm::translate (
    Geometry & g,
    const Kernel::Vector_2 & v
) 

function translate#

translate a geometry from a given vector

void algorithm::translate (
    Geometry & g,
    const Kernel::FT & dx,
    const Kernel::FT & dy,
    const Kernel::FT & dz
) 

function translate#

translate a geometry from double-coordinates

void algorithm::translate (
    Geometry & g,
    const double & dx,
    const double & dy,
    const double & dz
) 

function union3D#

auto algorithm::union3D (
    const Geometry & ga,
    const Geometry & gb,
    NoValidityCheck
) 

Union on 3D geometries. Assume z = 0 if needed

Precondition:

ga and gb are valid geometries

Warning:

@ No actual validity check is done


function union3D#

auto algorithm::union3D (
    const Geometry & ga,
    const Geometry & gb
) 

Union on 3D geometries. Assume z = 0 if needed

Precondition:

ga and gb are valid geometries


function union_#

template<int Dim>
void algorithm::union_ (
    const detail::GeometrySet < Dim > & a,
    const detail::GeometrySet < Dim > & b,
    detail::GeometrySet < Dim > &
) 

function union_#

template<int Dim>
void algorithm::union_ (
    const detail::PrimitiveHandle < Dim > & a,
    const detail::PrimitiveHandle < Dim > & b,
    detail::GeometrySet < Dim > &
) 

function visibility#

build the visibility polygon of a Point inside aPolygon __

auto algorithm::visibility (
    const Geometry & polygon,
    const Geometry & point
) 

Parameters:

  • polygon input geometry
  • point input geometry

Precondition:

polygon is a valid geometry

Precondition:

point must be inside polygon or on the boundary


function visibility#

build the visibility polygon of a Point inside aPolygon __

auto algorithm::visibility (
    const Geometry & polygon,
    const Geometry & point,
    NoValidityCheck
) 

Parameters:

  • polygon input geometry
  • point input geometry

Precondition:

polygon is a valid geometry

Precondition:

point must be inside polygon or on the boundary

Warning:

No actual validity check is done


function visibility#

build the visibility polygon of the segment [pointA ; pointB] on a Polygon __

auto algorithm::visibility (
    const Geometry & polygon,
    const Geometry & pointA,
    const Geometry & pointB
) 

Parameters:

  • polygon input geometry
  • pointA input geometry
  • pointB input geometry

Precondition:

polygon is a valid geometry

Precondition:

pointA and pointB must be vertices of poly, adjacents and respect the direction


function visibility#

build the visibility polygon of a Point inside aPolygon __

auto algorithm::visibility (
    const Geometry & polygon,
    const Geometry & pointA,
    const Geometry & pointB,
    NoValidityCheck
) 

Parameters:

  • polygon input geometry
  • pointA input geometry
  • pointB input geometry

Precondition:

polygon is a valid geometry

Warning:

No actual validity check is done

Precondition:

pointA and pointB must be vertices of poly, adjacents and respect the direction


function volume#

auto algorithm::volume (
    const Solid & g,
    NoValidityCheck
) 

Computes the volume of a Solid

Precondition:

(not checked) volume is closed and consistently oriented


function volume#

auto algorithm::volume (
    const Geometry & g
) 

Computes the volume of a geometry

Precondition:

g is a valid Geometry


function weightedCentroid#

auto algorithm::weightedCentroid (
    const Geometry & g,
    bool enable3DComputation=false
) 

Returns the weighted centroid for a Geometry


function weightedCentroid#

auto algorithm::weightedCentroid (
    const Triangle & g,
    bool enable3DComputation=false
) 

Returns the weighted centroid for a Triangle


function weightedCentroid#

auto algorithm::weightedCentroid (
    const Point & a,
    const Point & b,
    const Point & c,
    bool enable3DComputation=false
) 

Returns the weighted centroid for a Triangle


function weightedCentroid#

auto algorithm::weightedCentroid (
    const LineString & g,
    bool enable3DComputation=false
) 

Returns the weighted centroid for a LineString


function weightedCentroid#

auto algorithm::weightedCentroid (
    const Polygon & g,
    bool enable3DComputation=false
) 

Returns the weighted centroid for a Polygon


function weightedCentroid#

auto algorithm::weightedCentroid (
    const GeometryCollection & g,
    bool enable3DComputation=false
) 

Returns the weighted centroid for a GeometryCollection


function weightedCentroid#

auto algorithm::weightedCentroid (
    const TriangulatedSurface & g,
    bool enable3DComputation=false
) 

Returns the weighted centroid for a TriangulatedSurface


function weightedCentroid#

auto algorithm::weightedCentroid (
    const PolyhedralSurface & g,
    bool enable3DComputation=false
) 

Returns the weighted centroid for a PolyhedralSurface


function weightedCentroid#

auto algorithm::weightedCentroid (
    const Solid & g,
    bool enable3DComputation=false
) 

Returns the weighted centroid for a Solid



The documentation for this class was generated from the following file /builds/sfcgal/SFCGAL/src/algorithm/alphaShapes.cpp