OptoInspect3D Inline - alg3Dlib  3.4.0
Library for measured 3D data processing
Data Structures | Enumerations | Functions
Best-fit approximations

This module contains best-fit algorithms of point clouds with geometric primitives. More...

Data Structures

struct  ALG3D_Plane
 Data type for a plane. More...
 
struct  ALG3D_Line
 Data type for a 3D-line. More...
 
struct  ALG3D_Circle
 Data type for a 3D-circle. More...
 
struct  ALG3D_Ellipse
 Data type for a 3D-Ellipse. More...
 
struct  ALG3D_Superellipse
 Data type for a 3D-superellipse (Lamé-curve). More...
 
struct  ALG3D_Sphere
 Data type for a sphere. More...
 
struct  ALG3D_Cylinder
 Data type for a cylinder. More...
 
struct  ALG3D_Cone
 Data type for a cone. More...
 
struct  ALG3D_Torus
 Data type for a torus. More...
 

Enumerations

enum  ALG3D_ApproximateFlags {
  enFitNone =1L<< 0, enFixCenter =1L<< 1, enFixDirection =1L<< 2, enFixRadius =1L<< 3,
  enFixRadius2 =1L<< 4, enFixRotAngle =1L<< 5, enFixA =1L<< 6, enFixB =1L<< 7,
  enFixConeAngle =1L<< 8, enFixEpsilon =1L<< 9, enVarCenter =1L<<10, enVarDirection =1L<<11,
  enVarRadius =1L<<12, enVarRadius2 = 1L<<13, enVarRotAngle =1L<<14, enVarA =1L<<15,
  enVarB =1L<<16, enVarConeAngle =1L<<17, enVarEpsilon =1L<<18, enUseWeights =1L<<31
}
 

Functions

int ALG3D_fitLine (const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Line *line, int iterations, int flags)
 Computes best-fit line in 3D with least-squares. More...
 
int ALG3D_fitPlane (const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Plane *plane, int iterations, int flags)
 Computes best-fit plane with least-squares. More...
 
int ALG3D_fitCircle (const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Circle *circle, int iterations, int flags)
 Computes best-fit circle in 3D with least-squares. More...
 
int ALG3D_fitEllipse (const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Ellipse *ellipse, int iterations, int flags)
 Computes best-fit ellipse in 3D with least-squares. More...
 
int ALG3D_fitSuperellipse (const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Superellipse *superellipse, int iterations, int flags)
 Computes best-fit superellipse with least-squares. More...
 
int ALG3D_fitSphere (const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Sphere *sphere, int iterations, int flags)
 Computes best-fit sphere with least-squares. More...
 
int ALG3D_fitCylinder (const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Cylinder *cylinder, int iterations, int flags)
 Computes best-fit cylinder with least-squares. More...
 
int ALG3D_fitCone (const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Cone *cone, int iterations, int flags)
 Computes best-fit cone with least-squares. Note that the cone radius is measured at the height of the cone center. The cone radius is an output parameter only and is ignored when given as input parameter. More...
 
int ALG3D_fitTorus (const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Torus *torus, int iterations, int flags)
 Computes best-fit torus with least-squares. More...
 
int ALG3D_fitLineMZ (const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Line *line, double *absMax, int iterations)
 Computes best-fit 3D line using the Chebyshev norm. More...
 
int ALG3D_fitPlaneMZ (const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Plane *plane, double *absMax, int iterations)
 Computes best-fit plane using the Chebyshev norm. More...
 
int ALG3D_fitSphereMZ (const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Sphere *sphere, double *absMax, int iterations, int flags)
 Computes best-fit sphere using the Chebyshev norm. More...
 
int ALG3D_fitCylinderMZ (const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Cylinder *cylinder, double *absMax, int iterations, int flags)
 Computes best-fit cylinder using the Chebyshev norm. More...
 
int ALG3D_fitConeMZ (const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Cone *cone, double *absMax, int iterations, int flags)
 Computes best-fit cone using the Chebyshev norm. More...
 
int ALG3D_fitTorusMZ (const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Torus *torus, double *absMax, int iterations, int flags)
 Computes best-fit torus using the Chebyshev norm. More...
 
int ALG3D_fitSphereMC (const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Sphere *sphere, int iterations)
 Computes the minimum enclosing sphere. More...
 
int ALG3D_fitCylinderMC (const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Cylinder *cylinder, int iterations)
 Computes the minimum enclosing cylinder. More...
 
int ALG3D_fitPlanesMC (const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Plane *plane1, ALG3D_Plane *plane2, double *thickness, int iterations)
 Computes the minimum enclosing planes. More...
 
int ALG3D_fitSphereMI (const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Sphere *sphere, int iterations)
 Computes the maximum inscribed sphere (Pferch sphere). More...
 
int ALG3D_fitCylinderMI (const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Cylinder *cylinder, int iterations)
 Computes the maximum inscribed cylinder (Pferch cylinder). More...
 
int ALG3D_distancesPoints2Line (const ALG3D_Point3d *ptsIn, size_t numIn, const ALG3D_Line *line, double *distances)
 Computes the distances of points to a 3D line. More...
 
int ALG3D_distancesPoints2Plane (const ALG3D_Point3d *ptsIn, size_t numIn, const ALG3D_Plane *plane, double *distances)
 Computes the distances of points to a plane. More...
 
int ALG3D_distancesPoints2Circle (const ALG3D_Point3d *ptsIn, size_t numIn, const ALG3D_Circle *circle, double *distances)
 Computes the distances of points to a 3D circle. More...
 
int ALG3D_distancesPoints2Sphere (const ALG3D_Point3d *ptsIn, size_t numIn, const ALG3D_Sphere *sphere, double *distances)
 Computes the distances of points to a sphere. More...
 
int ALG3D_distancesPoints2Cylinder (const ALG3D_Point3d *ptsIn, size_t numIn, const ALG3D_Cylinder *cylinder, double *distances)
 Computes the distances of points to a cylinder. More...
 
int ALG3D_distancesPoints2Cone (const ALG3D_Point3d *ptsIn, size_t numIn, const ALG3D_Cone *cone, double *distances)
 Computes the distances of points to a cone. More...
 
int ALG3D_distancesPoints2Torus (const ALG3D_Point3d *ptsIn, size_t numIn, const ALG3D_Torus *torus, double *distances)
 Computes the distances of points to a torus. More...
 
int ALG3D_distancesPoints2Ellipse (const ALG3D_Point3d *ptsIn, size_t numIn, const ALG3D_Ellipse *ellipse, double *distances)
 Computes the distances of points to a 3D ellipse. More...
 
int ALG3D_distancesPoints2Superellipse (const ALG3D_Point3d *ptsIn, size_t numIn, const ALG3D_Superellipse *superellipse, double *distances)
 Computes the distances of points to a 3D superellipse. More...
 
int ALG3D_projectPoints2Line (const ALG3D_Point3d *ptsIn, size_t numIn, const ALG3D_Line *line, ALG3D_Point3d *ptsOut)
 Computes the projections of points to a 3D line. More...
 
int ALG3D_projectPoints2Plane (const ALG3D_Point3d *ptsIn, size_t numIn, const ALG3D_Plane *plane, ALG3D_Point3d *ptsOut)
 Computes the projections of points to a plane. More...
 
int ALG3D_projectPoints2Sphere (const ALG3D_Point3d *ptsIn, size_t numIn, const ALG3D_Sphere *sphere, ALG3D_Point3d *ptsOut)
 Computes the projections of points to a sphere. More...
 
int ALG3D_projectPoints2Cylinder (const ALG3D_Point3d *ptsIn, size_t numIn, const ALG3D_Cylinder *cylinder, ALG3D_Point3d *ptsOut)
 Computes the projections of points to a cylinder. More...
 
int ALG3D_projectPoints2Cone (const ALG3D_Point3d *ptsIn, size_t numIn, const ALG3D_Cone *cone, ALG3D_Point3d *ptsOut)
 Computes the projections of points to a cone. More...
 
int ALG3D_projectPoints2Torus (const ALG3D_Point3d *ptsIn, size_t numIn, const ALG3D_Torus *torus, ALG3D_Point3d *ptsOut)
 Computes the projections of points to a torus. More...
 
int ALG3D_intersectLineLine (const ALG3D_Line *line1, const ALG3D_Line *line2, ALG3D_Point3d *point)
 Computes the intersection between two 3D lines. If the two lines are not parallel an intersection point will be computed. This point is computed from the position where the lines are closest. The middle point between the closest points on the lines is returned as the intersection point. More...
 
int ALG3D_intersectLinePlane (const ALG3D_Line *line, const ALG3D_Plane *plane, ALG3D_Point3d *point)
 Computes the intersection between 3D line and a plane. More...
 
int ALG3D_intersectLineSphere (const ALG3D_Line *line, const ALG3D_Sphere *sphere, ALG3D_Point3d *ptsOut, size_t *numOut)
 Computes the intersection between 3D line and a sphere. More...
 
int ALG3D_intersectLineCylinder (const ALG3D_Line *line, const ALG3D_Cylinder *cylinder, ALG3D_Point3d *ptsOut, size_t *numOut)
 Computes the intersection between a 3D line and a cylinder. More...
 
int ALG3D_intersectLineCone (const ALG3D_Line *line, const ALG3D_Cone *cone, ALG3D_Point3d *ptsOut, size_t *numOut)
 Computes the intersection between a 3D line and a cone. More...
 
int ALG3D_intersectPlanePlane (const ALG3D_Plane *plane1, const ALG3D_Plane *plane2, ALG3D_Line *line)
 Computes the intersection between two planes. More...
 
int ALG3D_intersectPlaneSphere (const ALG3D_Plane *plane, const ALG3D_Sphere *sphere, ALG3D_Circle *circle)
 Computes the intersection between a plane and a sphere. More...
 
int ALG3D_computePrincipleComponents (const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Point3d components[3])
 Computes the principal components of a 3D point set. The principal components are the directions of the main axes of a point sets. The three directions (normalized vectors) are perpendicular to each other and thus form a orthogonal system. Principal components are useful to analyze an unknown point set, e.g. in order to get to largest and smallest extents. The functions returns these components as three normalized vectors. More...
 
int ALG3D_computeMinimumBoundingCircle (const ALG3D_Point3d *v1, const ALG3D_Point3d *v2, const ALG3D_Point3d *v3, ALG3D_Point3d *center, double *radius)
 Computes the minimum bounding circle for a triangle. Computes the smallest circle that contains the 3 given vertices, but the circle must not pass through all of the vertices. The minimum bounding circle differs from the circumcircle, which requires to pass through all vertices and is equal or larger to the minimum bounding circle. More...
 

Functions for Triangles

int __cdecl ALG3D_computeConvexHull (const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Point3d **vertices, size_t *numVertices, size_t **indices, size_t *numIndices)
 Computes the convex hull of a 3D point set. More...
 

Detailed Description

This module contains best-fit algorithms of point clouds with geometric primitives.

The most important optimization strategies are supported: classical Least-Squares (Gauss) and Minimum-Zone (Chebyshev). The least-squares algorithms utilize a combination of a singular-value-decomposition and a Levenberg-Marquardt strategy to solve the equation system with high numerical robustness. The optimization criterion is the minimization of the squared-distance perpendicular to the primitive's surface. The Minimum Zone algorithms solve a much more difficult problem that is the numerical minimization of the absolute maximum distance measured perpendicular to the primitiv's surface. Numerical methods from linear and quadratic programming are implemented in a highly efficient manner that are numerically robust even for large points clouds.

For each of the best-fit function constrains may be applied, eg. in order to keep the radius of a cylinder fixed and just solve for the remaining parameters. Therefore, the ALG3D_ApproximateFlags can be used. Most of the functions require to specify a maximum number of iterations for solving the optimization problem. Since this is a maximum number one may define for example 100. Typically, less than 10 iterations are required, only in ill-conditioned cases more iterations may be required.


Data Structure Documentation

◆ ALG3D_Plane

struct ALG3D_Plane

Data type for a plane.

Data Fields
ALG3D_Point3d center Point on plane.
ALG3D_Point3d direction Direction vector of the plane.

◆ ALG3D_Line

struct ALG3D_Line

Data type for a 3D-line.

Data Fields
ALG3D_Point3d center Point on the line.
ALG3D_Point3d direction Direction vector of the line.

◆ ALG3D_Circle

struct ALG3D_Circle

Data type for a 3D-circle.

Data Fields
ALG3D_Point3d center Center of the circle.
ALG3D_Point3d direction Direction vector of the plane the circle lies in.
double radius Radius of the circle.

◆ ALG3D_Ellipse

struct ALG3D_Ellipse

Data type for a 3D-Ellipse.

Data Fields
double a
double angle Angular position within the plane with respect to the longer axis a.
double b Radii of the ellipse, with a >= b.
ALG3D_Point3d center Center of the ellipse.
ALG3D_Point3d direction Direction vector of the plane the ellipse lies in.

◆ ALG3D_Superellipse

struct ALG3D_Superellipse

Data type for a 3D-superellipse (Lamé-curve).

Data Fields
double a
double angle Angular position within the plane with respect to the longer axis a.
double b Radii of the superellipse, with a >= b.
ALG3D_Point3d center Center of the superellipse.
ALG3D_Point3d direction Direction vector of the plane the superellipse lies in.
double epsilon Shape factor of the superellipse: rectangle (e<<1), ellipse/circle (e==1), diamond (e==2), star (e>2).

◆ ALG3D_Sphere

struct ALG3D_Sphere

Data type for a sphere.

Data Fields
ALG3D_Point3d center Center of the sphere.
double radius Radius of the sphere.

◆ ALG3D_Cylinder

struct ALG3D_Cylinder

Data type for a cylinder.

Data Fields
ALG3D_Point3d center Center of the cylinder (projection of the mass center onto the cylinder axis)
ALG3D_Point3d direction Direction vector of the cylinder axis.
double radius Radius of the cylinder.

◆ ALG3D_Cone

struct ALG3D_Cone

Data type for a cone.

Data Fields
double angle Angle of the cone.
ALG3D_Point3d center Center of the cone (projection of the mass center onto the cone axis)
ALG3D_Point3d direction Direction vector of the cone axis.
double radius Radius of the cone (at the height of the center point)

◆ ALG3D_Torus

struct ALG3D_Torus

Data type for a torus.

Data Fields
ALG3D_Point3d center Center of the torus.
ALG3D_Point3d direction Direction vector of the torus plane.
double radius1 Radius of the ring.
double radius2 Radius of the tube.

Enumeration Type Documentation

◆ ALG3D_ApproximateFlags

#include <approximate.h>

Enum that provides flags to control the behaviour of the best-fit algoithms. Parameters of a geometric primitive may be initialized in order to support the optimization functions. Parameters may also be set fixed, which forces the algorithm to keep this value unchanged during optimization.

Enumerator
enFitNone 

set nothing, reset

enFixCenter 

center position

enFixDirection 

direction vector of the axis

enFixRadius 

radius (used by circle, sphere, cylinder, torus)

enFixRadius2 

second radius (used by torus)

enFixRotAngle 

rotation angle of 2D-shapes

enFixA 

length of the great ellipse/superellipse half-axis

enFixB 

length of the small ellipse/superellipse half-axis

enFixConeAngle 

cone angle

enFixEpsilon 

shape parameter of a superellipse

enVarCenter 

center position

enVarDirection 

direction vector of the axis

enVarRadius 

radius (used by circle, sphere, cylinder, torus)

enVarRadius2 

second radius (used by torus)

enVarRotAngle 

rotation angle of 2D-shapes

enVarA 

length of the great ellipse/superellipse half-axis

enVarB 

length of the small ellipse/superellipse half-axis

enVarConeAngle 

cone angle

enVarEpsilon 

shape parameter of a superellipse

enUseWeights 

enable automated large outlier weighting (use with care)

Function Documentation

◆ ALG3D_computeConvexHull()

int __cdecl ALG3D_computeConvexHull ( const ALG3D_Point3d ptsIn,
size_t  numIn,
ALG3D_Point3d **  vertices,
size_t *  numVertices,
size_t **  indices,
size_t *  numIndices 
)

#include <geometry.h>

Computes the convex hull of a 3D point set.

Computes the convex hull of a 3D point set as a triangle mesh. The function returns a list of vertices and indices that build the triangle mesh For example: the first triangle contains 3 vertices: v1=vertices[indices[0]], v2=vertices[indices[1]], v3=vertices[indices[2]]

code example
const char* fname = "./sample_data/line.xyz";
ALG3D_Point3d *pts = 0; //3d point data
size_t numPts = 0; //number of 3d points
ALG3D_Point3d *vertices = 0; //resulting hull mesh vertices
size_t numVertices = 0; //number of resulting hull mesh vertices
size_t* indices = 0; //resulting triangle indices of the hull mesh
size_t numIndices = 0; //number of resulting triangle indices
/***** read points from file *****/
if (!ALG3D_checkResult(ALG3D_readXYZ(fname, &pts, &numPts))) return 0; //read file
/***** compute convex hull *****/
if (!ALG3D_checkResult(ALG3D_computeConvexHull(pts, numPts, &vertices, &numVertices, &indices, &numIndices))) return 0; //hull
const size_t numTriangles = numIndices / 3; //number of resulting triangles
printf("triangles: %I64d\n", numTriangles);
/*printf("[convex hull mesh]\ntriangles: %I64d, vertices: %I64d\n", numTriangles, numVertices);
for (size_t i = 0; i < numTriangles; i+=3)
{
ALG3D_Point3d p1 = vertices[indices[i + 0]];
ALG3D_Point3d p2 = vertices[indices[i + 1]];
ALG3D_Point3d p3 = vertices[indices[i + 2]];
printf("[%3I64d] P1: %lf,%lf,%lf\n",i, p1.vec[0], p1.vec[1], p1.vec[2]);
printf("[%3I64d] P2: %lf,%lf,%lf\n",i, p2.vec[0], p2.vec[1], p2.vec[2]);
printf("[%3I64d] P3: %lf,%lf,%lf\n",i, p3.vec[0], p3.vec[1], p3.vec[2]);
}*/
/********************************/
Parameters
[in]ptsInpointer to begin of point data
[in]numInnumber of points
[out]verticespointer to begin of resulting point data for the vertices of the convex hull triangles (a sub of the input points)
[out]numVerticesresulting number of vertices
[out]indicespointer to begin of resulting indices for the triangles pointing to vertices (in counterclockwise order).
[out]numIndicesnumber of points
Returns
ERR_NONE on success. ERR_NUM_POINTS if less than 3 points are provided.

◆ ALG3D_computeMinimumBoundingCircle()

int ALG3D_computeMinimumBoundingCircle ( const ALG3D_Point3d v1,
const ALG3D_Point3d v2,
const ALG3D_Point3d v3,
ALG3D_Point3d center,
double *  radius 
)

#include <geometry.cpp>

Computes the minimum bounding circle for a triangle. Computes the smallest circle that contains the 3 given vertices, but the circle must not pass through all of the vertices. The minimum bounding circle differs from the circumcircle, which requires to pass through all vertices and is equal or larger to the minimum bounding circle.

code example
// ALG3D_Point3d pt; //3d point
// ALG3D_Point3d v1, v2, v3; //vertices of a triangle
// ALG3D_Point3d nn, center; //points
// double dist, radius;
//
// //some test data
// v1.vec[0] = 1.33398; v1.vec[1] = -0.467975; v1.vec[2] = -0.0388361;
// v2.vec[0] = 0.67865; v2.vec[1] = 0.842713; v2.vec[2] = 0.281547;
// v3.vec[0] = -1.11553; v3.vec[1] = 0.0310788; v3.vec[2] = -0.0679581;
// pt.vec[0] = -0.0174734; pt.vec[1] = 0.134954; pt.vec[2] = 0.348543;
//
// //printf("\n*** triangle functions ***\n");
//
// //compute closest point to triangle
// if (!ALG3D_checkResult(ALG3D_computeClosestPoint2Triangle(&pt, &v1, &v2, &v3, &nn))) return 0;
//
// dist = distance(nn, pt);
// //printf("closest point: %lf,%lf,%lf distance: %lf\n", nn.vec[0], nn.vec[1], nn.vec[2], dist);
//
// //minimum bounding circle of a triangle
// if (!ALG3D_checkResult(ALG3D_computeMinimumBoundingCircle(&v1, &v2, &v3, &center, &radius))) return 0;
// //printf("minimum bounding circle: pt: %lf,%lf,%lf radius: %lf\n", center.vec[0], center.vec[1], center.vec[2], radius);
//
Parameters
v1,v2,v3points/vertices of the triangle
centercenter of the minimum bounding circle
radiusradius of the minimum bounding circle
Returns
ERR_NONE on success

◆ ALG3D_computePrincipleComponents()

int ALG3D_computePrincipleComponents ( const ALG3D_Point3d ptsIn,
size_t  numIn,
ALG3D_Point3d  components[3] 
)

#include <geometry.cpp>

Computes the principal components of a 3D point set. The principal components are the directions of the main axes of a point sets. The three directions (normalized vectors) are perpendicular to each other and thus form a orthogonal system. Principal components are useful to analyze an unknown point set, e.g. in order to get to largest and smallest extents. The functions returns these components as three normalized vectors.

The components are sorted in decreasing order, that means, the first vector contains the direction of the largest extent of the point set.

code example
const char* fname ="./sample_data/cone.xyz";
ALG3D_Point3d *pts=0; //3d point data
size_t numPts=0; //number of 3d points
ALG3D_Point3d components[3];//3 vectors for the principal components
int i;
if(!ALG3D_checkResult(ALG3D_readXYZ(fname, &pts, &numPts))) return 0; //read file
/***** principal components *****/
if(!ALG3D_checkResult(ALG3D_computePrincipleComponents(pts, numPts, components) )) return 0; //compute the principle components
/********************************/
printf("[principal components]\n");
for(i=0; i<3; ++i){ //print out the principle components
printf("[%d]: %lf,%lf,%lf\n", i, components[i].vec[0], components[i].vec[1], components[i].vec[2]);
}
Parameters
ptsInpointer to the begin of the point data
numInnumber of points
componentsthree normalized vectors containing the principal components
Note
This function always uses all available threads.
Returns
ERR_NONE on success

◆ ALG3D_distancesPoints2Circle()

int ALG3D_distancesPoints2Circle ( const ALG3D_Point3d ptsIn,
size_t  numIn,
const ALG3D_Circle circle,
double *  distances 
)

#include <geometry.cpp>

Computes the distances of points to a 3D circle.

Parameters
ptsInpointer to begin of point data
numInnumber of points
circlecircle object
distancespointer to existing storage for the resulting distance values (must have the size of numIn). A negative distance indicates that the point lies inside the geometry.
Note
This function always uses all available threads.
Returns
ERR_NONE on success

◆ ALG3D_distancesPoints2Cone()

int ALG3D_distancesPoints2Cone ( const ALG3D_Point3d ptsIn,
size_t  numIn,
const ALG3D_Cone cone,
double *  distances 
)

#include <geometry.cpp>

Computes the distances of points to a cone.

Parameters
ptsInpointer to begin of point data
numInnumber of points
conecone object
distancespointer to existing storage for the resulting distance values (must have the size of numIn). A negative distance indicates that the point lies inside the geometry.
Note
This function always uses all available threads.
Returns
ERR_NONE on success

◆ ALG3D_distancesPoints2Cylinder()

int ALG3D_distancesPoints2Cylinder ( const ALG3D_Point3d ptsIn,
size_t  numIn,
const ALG3D_Cylinder cylinder,
double *  distances 
)

#include <geometry.cpp>

Computes the distances of points to a cylinder.

Parameters
ptsInpointer to begin of point data
numInnumber of points
cylindercylinder object
distancespointer to existing storage for the resulting distance values (must have the size of numIn). A negative distance indicates that the point lies inside the geometry.
Note
This function always uses all available threads.
Returns
ERR_NONE on success

◆ ALG3D_distancesPoints2Ellipse()

int ALG3D_distancesPoints2Ellipse ( const ALG3D_Point3d ptsIn,
size_t  numIn,
const ALG3D_Ellipse ellipse,
double *  distances 
)

#include <geometry.cpp>

Computes the distances of points to a 3D ellipse.

Parameters
ptsInpointer to begin of point data
numInnumber of points
ellipseellipse object
distancespointer to existing storage for the resulting distance values (must have the size of numIn). A negative distance indicates that the point lies inside the geometry.
Note
This function always uses all available threads.
Returns
ERR_NONE on success

◆ ALG3D_distancesPoints2Line()

int ALG3D_distancesPoints2Line ( const ALG3D_Point3d ptsIn,
size_t  numIn,
const ALG3D_Line line,
double *  distances 
)

#include <geometry.cpp>

Computes the distances of points to a 3D line.

Parameters
ptsInpointer to begin of point data
numInnumber of points
lineline object
distancespointer to existing storage for the resulting distance values (must have the size of numIn)
Note
This function always uses all available threads.
Returns
ERR_NONE on success

◆ ALG3D_distancesPoints2Plane()

int ALG3D_distancesPoints2Plane ( const ALG3D_Point3d ptsIn,
size_t  numIn,
const ALG3D_Plane plane,
double *  distances 
)

#include <geometry.cpp>

Computes the distances of points to a plane.

Parameters
ptsInpointer to begin of point data
numInnumber of points
planeplane object
distancespointer to existing storage for the resulting distance values (must have the size of numIn). A negative distance indicates that the point lies on the opposite side than the planes normal directs.
Note
This function always uses all available threads.
Returns
ERR_NONE on success

◆ ALG3D_distancesPoints2Sphere()

int ALG3D_distancesPoints2Sphere ( const ALG3D_Point3d ptsIn,
size_t  numIn,
const ALG3D_Sphere sphere,
double *  distances 
)

#include <geometry.cpp>

Computes the distances of points to a sphere.

Parameters
ptsInpointer to begin of point data
numInnumber of points
spheresphere object
distancespointer to existing storage for the resulting distance values (must have the size of numIn). A negative distance indicates that the point lies inside the geometry.
Note
This function always uses all available threads.
Returns
ERR_NONE on success

◆ ALG3D_distancesPoints2Superellipse()

int ALG3D_distancesPoints2Superellipse ( const ALG3D_Point3d ptsIn,
size_t  numIn,
const ALG3D_Superellipse superellipse,
double *  distances 
)

#include <geometry.cpp>

Computes the distances of points to a 3D superellipse.

Parameters
ptsInpointer to begin of point data
numInnumber of points
superellipsesuperellipse object
distancespointer to existing storage for the resulting distance values (must have the size of numIn). A negative distance indicates that the point lies inside the geometry.
Note
This function always uses all available threads.
Returns
ERR_NONE on success

◆ ALG3D_distancesPoints2Torus()

int ALG3D_distancesPoints2Torus ( const ALG3D_Point3d ptsIn,
size_t  numIn,
const ALG3D_Torus torus,
double *  distances 
)

#include <geometry.cpp>

Computes the distances of points to a torus.

Parameters
ptsInpointer to begin of point data
numInnumber of points
torustorus object
distancespointer to existing storage for the resulting distance values (must have the size of numIn). A negative distance indicates that the point lies inside the tube of the torus.
Note
This function always uses all available threads.
Returns
ERR_NONE on success

◆ ALG3D_fitCircle()

int ALG3D_fitCircle ( const ALG3D_Point3d ptsIn,
size_t  numIn,
ALG3D_Circle circle,
int  iterations,
int  flags 
)

#include <approximate.cpp>

Computes best-fit circle in 3D with least-squares.

code example
const char* fname = "./sample_data/circle.xyz";
ALG3D_Point3d *vec = 0; //3d point data
size_t ptsVec = 0; //number of 3d points
ALG3D_Circle circle; //circle object
int flags = 0; //parameter flags
//read file
if (!ALG3D_checkResult(ALG3D_readXYZ(fname, &vec, &ptsVec))) return 0;
//compute least-squares 3d circle
if (!ALG3D_checkResult(ALG3D_fitCircle(vec, ptsVec, &circle, 100, flags))) return 0;
printf("pt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr: %lf\n", circle.center.vec[0], circle.center.vec[1], circle.center.vec[2],
circle.direction.vec[0], circle.direction.vec[1], circle.direction.vec[2], circle.radius);
Parameters
ptsInpointer to begin of point data
numInnumber of points
circlecircle object
iterationsmaximum number of iterations (e.g. 100 is sufficient in most cases)
flagsinitialization flags
Returns
ERR_NONE on success. ERR_NUM_POINTS if less than 4 points are provided.

◆ ALG3D_fitCone()

int ALG3D_fitCone ( const ALG3D_Point3d ptsIn,
size_t  numIn,
ALG3D_Cone cone,
int  iterations,
int  flags 
)

#include <approximate.cpp>

Computes best-fit cone with least-squares. Note that the cone radius is measured at the height of the cone center. The cone radius is an output parameter only and is ignored when given as input parameter.

code example
const double PI = 3.1415926535;
const char* fname = "./sample_data/cone.xyz";
ALG3D_Point3d *vec = 0; //3d point data
size_t ptsVec = 0; //number of 3d points
double absMax = 0; //absolute maximum deviaton (for minimum zone-approximation)
ALG3D_Cone cone; //cone object
int flags = 0; //parameter flags
//read file
if (!ALG3D_checkResult(ALG3D_readXYZ(fname, &vec, &ptsVec))) return 0;
//compute least-squares cone
if (!ALG3D_checkResult(ALG3D_fitCone(vec, ptsVec, &cone, 100, flags))) return 0;
printf("[*Least-squares cone*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr:%lf\tangle: %lf deg\n", cone.center.vec[0], cone.center.vec[1], cone.center.vec[2],
cone.direction.vec[0], cone.direction.vec[1], cone.direction.vec[2], cone.radius, cone.angle * 180 / PI);
//compute minimum-zone cone
if (!ALG3D_checkResult(ALG3D_fitConeMZ(vec, ptsVec, &cone, &absMax, 100, flags))) return 0;
printf("[*Minimum zone cone*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr:%lf\tangle: %lf deg\nabsMax: %lf\n", cone.center.vec[0], cone.center.vec[1], cone.center.vec[2],
cone.direction.vec[0], cone.direction.vec[1], cone.direction.vec[2], cone.radius, cone.angle * 180 / PI, absMax);
Note
This function returns a best fit result, that is not yet compliant with the PTB. The qualitity of the fit is the same, but the PTB has a different definition of the radius and the center of a cone. To transform the result of this function into a PTB compliant result, you must apply the following steps:
ALG3D_Line xline;
xline.center = cone.center;
xline.direction = cone.direction;
// compute center of points from which you approximated the cylinder ...
ALG3D_projectPoints2Line(&ctr, 1, &xline, &ctr);
const double pdx = ALG3D_dotProduct(&cone.center, &cone.direction)
- ALG3D_dotProduct(&ctr, &cone.direction);
const double pdiff = pdx / tan(0.5 * (M_PI - cone.angle));
const double pradius = cone.radius + pdiff;
// ctr and pradius are now the PTB compliant versions of the center and
// radius of the cone
Parameters
ptsInpointer to begin of point data
numInnumber of points
conecone object
iterationsmaximum number of iterations (e.g. 100 is sufficient in most cases)
flagsinitialization flags. Note, that the radius cannot be initialized, because it is a return value
Returns
ERR_NONE on success. ERR_NUM_POINTS if less than 5 points are provided.

◆ ALG3D_fitConeMZ()

int ALG3D_fitConeMZ ( const ALG3D_Point3d ptsIn,
size_t  numIn,
ALG3D_Cone cone,
double *  absMax,
int  iterations,
int  flags 
)

#include <approximate.cpp>

Computes best-fit cone using the Chebyshev norm.

This function computes the best-fit cone by minimizing the Chebyshev norm (absolute maximum deviation). The resulting cone can be used to analyzed the conicity of the data points. The maximum absolute deviation is written to the variable absMax.

code example
const double PI = 3.1415926535;
const char* fname = "./sample_data/cone.xyz";
ALG3D_Point3d *vec = 0; //3d point data
size_t ptsVec = 0; //number of 3d points
double absMax = 0; //absolute maximum deviaton (for minimum zone-approximation)
ALG3D_Cone cone; //cone object
int flags = 0; //parameter flags
//read file
if (!ALG3D_checkResult(ALG3D_readXYZ(fname, &vec, &ptsVec))) return 0;
//compute least-squares cone
if (!ALG3D_checkResult(ALG3D_fitCone(vec, ptsVec, &cone, 100, flags))) return 0;
printf("[*Least-squares cone*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr:%lf\tangle: %lf deg\n", cone.center.vec[0], cone.center.vec[1], cone.center.vec[2],
cone.direction.vec[0], cone.direction.vec[1], cone.direction.vec[2], cone.radius, cone.angle * 180 / PI);
//compute minimum-zone cone
if (!ALG3D_checkResult(ALG3D_fitConeMZ(vec, ptsVec, &cone, &absMax, 100, flags))) return 0;
printf("[*Minimum zone cone*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr:%lf\tangle: %lf deg\nabsMax: %lf\n", cone.center.vec[0], cone.center.vec[1], cone.center.vec[2],
cone.direction.vec[0], cone.direction.vec[1], cone.direction.vec[2], cone.radius, cone.angle * 180 / PI, absMax);
Parameters
[in]ptsInpointer to begin of point data
[in]numInnumber of points
[in,out]conecone object
[out]absMaxcontains the absolute maximum deviation afterwards
[in]iterationsmaximum number of iterations (e.g. 100 is sufficient in most cases)
[in]flagsinitialization flags
Returns
ERR_NONE on success. ERR_NUM_POINTS if less than 6 points are provided.

◆ ALG3D_fitCylinder()

int ALG3D_fitCylinder ( const ALG3D_Point3d ptsIn,
size_t  numIn,
ALG3D_Cylinder cylinder,
int  iterations,
int  flags 
)

#include <approximate.cpp>

Computes best-fit cylinder with least-squares.

code example
const char* fname = "./sample_data/cylinder.xyz";
const char* fname2 = "./sample_data/pipe.xyz";
ALG3D_Point3d *vec = 0, *vec2 = 0; //3d point data
size_t ptsVec = 0, ptsVec2 = 0; //number of 3d points
ALG3D_Cylinder cylinder; //cylinder object
double absMax = 0; //absolute maximum deviaton (for minimum zone-approximation)
int flags = 0; //parameter flags
//read files
if (!ALG3D_checkResult(ALG3D_readXYZ(fname, &vec, &ptsVec))) return 0;
if (!ALG3D_checkResult(ALG3D_readXYZ(fname2, &vec2, &ptsVec2))) return 0;
//initialize cylinder with roughly known values
cylinder.radius = 17.0; //expected value or guess
flags |= enVarRadius; //set expected value to be changeable (if the expected value should not be changed, use 'enFixRadius' instead)
//compute least-squares cylinder
if (!ALG3D_checkResult(ALG3D_fitCylinder(vec, ptsVec, &cylinder, 100, flags))) return 0;
printf("[*Least-squares cylinder*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr: %lf\n", cylinder.center.vec[0], cylinder.center.vec[1], cylinder.center.vec[2],
cylinder.direction.vec[0], cylinder.direction.vec[1], cylinder.direction.vec[2], cylinder.radius);
//compute minimum-zone cylinder
if (!ALG3D_checkResult(ALG3D_fitCylinderMZ(vec, ptsVec, &cylinder, &absMax, 100, flags))) return 0;
printf("[*Minimum zone cylinder*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr: %lf\nabsMax: %lf\n", cylinder.center.vec[0], cylinder.center.vec[1], cylinder.center.vec[2],
cylinder.direction.vec[0], cylinder.direction.vec[1], cylinder.direction.vec[2], cylinder.radius, absMax);
//compute minimum enclosing cylinder
if (!ALG3D_checkResult(ALG3D_fitCylinderMC(vec, ptsVec, &cylinder, 100))) return 0;
printf("[*Minimum enclosing cylinder*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr: %lf\n", cylinder.center.vec[0], cylinder.center.vec[1], cylinder.center.vec[2],
cylinder.direction.vec[0], cylinder.direction.vec[1], cylinder.direction.vec[2], cylinder.radius);
//compute maximum inscribed cylinder
if (!ALG3D_checkResult(ALG3D_fitCylinderMI(vec2, ptsVec2, &cylinder, 100))) return 0;
printf("[*Maximum inscribed cylinder*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr: %lf\n", cylinder.center.vec[0], cylinder.center.vec[1], cylinder.center.vec[2],
cylinder.direction.vec[0], cylinder.direction.vec[1], cylinder.direction.vec[2], cylinder.radius);
Note
This function returns a best fit result, that is not yet compliant with the PTB. The qualitity of the fit is the same, but the PTB has a different definition of the center of a cylinder. To transform the result of this function into a PTB compliant result, you must apply the following steps:
ALG3D_Line xline;
xline.center = cylinder.center;
xline.direction = cylinder.direction;
// compute center of points from which you approximated the cylinder ...
ALG3D_projectPoints2Line(&ctr, 1, &xline, &ctr);
// ctr is now the PTB compliant center of the cylinder
Parameters
ptsInpointer to begin of point data
numInnumber of points
cylindercylinder object
iterationsmaximum number of iterations (e.g. 100 is sufficient in most cases)
flagsinitialization flags
Returns
ERR_NONE on success. ERR_NUM_POINTS if less than 4 points are provided.

◆ ALG3D_fitCylinderMC()

int ALG3D_fitCylinderMC ( const ALG3D_Point3d ptsIn,
size_t  numIn,
ALG3D_Cylinder cylinder,
int  iterations 
)

#include <approximate.cpp>

Computes the minimum enclosing cylinder.

This function computes the minimum enclosing cylinder, i.e. the smallest cylinder that contains all points. It is equivalent to MZ line.

code example
const char* fname = "./sample_data/cylinder.xyz";
const char* fname2 = "./sample_data/pipe.xyz";
ALG3D_Point3d *vec = 0, *vec2 = 0; //3d point data
size_t ptsVec = 0, ptsVec2 = 0; //number of 3d points
ALG3D_Cylinder cylinder; //cylinder object
double absMax = 0; //absolute maximum deviaton (for minimum zone-approximation)
int flags = 0; //parameter flags
//read files
if (!ALG3D_checkResult(ALG3D_readXYZ(fname, &vec, &ptsVec))) return 0;
if (!ALG3D_checkResult(ALG3D_readXYZ(fname2, &vec2, &ptsVec2))) return 0;
//initialize cylinder with roughly known values
cylinder.radius = 17.0; //expected value or guess
flags |= enVarRadius; //set expected value to be changeable (if the expected value should not be changed, use 'enFixRadius' instead)
//compute least-squares cylinder
if (!ALG3D_checkResult(ALG3D_fitCylinder(vec, ptsVec, &cylinder, 100, flags))) return 0;
printf("[*Least-squares cylinder*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr: %lf\n", cylinder.center.vec[0], cylinder.center.vec[1], cylinder.center.vec[2],
cylinder.direction.vec[0], cylinder.direction.vec[1], cylinder.direction.vec[2], cylinder.radius);
//compute minimum-zone cylinder
if (!ALG3D_checkResult(ALG3D_fitCylinderMZ(vec, ptsVec, &cylinder, &absMax, 100, flags))) return 0;
printf("[*Minimum zone cylinder*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr: %lf\nabsMax: %lf\n", cylinder.center.vec[0], cylinder.center.vec[1], cylinder.center.vec[2],
cylinder.direction.vec[0], cylinder.direction.vec[1], cylinder.direction.vec[2], cylinder.radius, absMax);
//compute minimum enclosing cylinder
if (!ALG3D_checkResult(ALG3D_fitCylinderMC(vec, ptsVec, &cylinder, 100))) return 0;
printf("[*Minimum enclosing cylinder*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr: %lf\n", cylinder.center.vec[0], cylinder.center.vec[1], cylinder.center.vec[2],
cylinder.direction.vec[0], cylinder.direction.vec[1], cylinder.direction.vec[2], cylinder.radius);
//compute maximum inscribed cylinder
if (!ALG3D_checkResult(ALG3D_fitCylinderMI(vec2, ptsVec2, &cylinder, 100))) return 0;
printf("[*Maximum inscribed cylinder*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr: %lf\n", cylinder.center.vec[0], cylinder.center.vec[1], cylinder.center.vec[2],
cylinder.direction.vec[0], cylinder.direction.vec[1], cylinder.direction.vec[2], cylinder.radius);
Parameters
[in]ptsInpointer to begin of point data
[in]numInnumber of points
[out]cylindercylinder object
[in]iterationsmaximum number of iterations (e.g. 100 is sufficient in most cases)
Returns
ERR_NONE on success.

◆ ALG3D_fitCylinderMI()

int ALG3D_fitCylinderMI ( const ALG3D_Point3d ptsIn,
size_t  numIn,
ALG3D_Cylinder cylinder,
int  iterations 
)

#include <approximate.cpp>

Computes the maximum inscribed cylinder (Pferch cylinder).

This function computes the maximum inscribed cylinder, i.e. the largest cylinder that fits inside the given points.

code example
const char* fname = "./sample_data/cylinder.xyz";
const char* fname2 = "./sample_data/pipe.xyz";
ALG3D_Point3d *vec = 0, *vec2 = 0; //3d point data
size_t ptsVec = 0, ptsVec2 = 0; //number of 3d points
ALG3D_Cylinder cylinder; //cylinder object
double absMax = 0; //absolute maximum deviaton (for minimum zone-approximation)
int flags = 0; //parameter flags
//read files
if (!ALG3D_checkResult(ALG3D_readXYZ(fname, &vec, &ptsVec))) return 0;
if (!ALG3D_checkResult(ALG3D_readXYZ(fname2, &vec2, &ptsVec2))) return 0;
//initialize cylinder with roughly known values
cylinder.radius = 17.0; //expected value or guess
flags |= enVarRadius; //set expected value to be changeable (if the expected value should not be changed, use 'enFixRadius' instead)
//compute least-squares cylinder
if (!ALG3D_checkResult(ALG3D_fitCylinder(vec, ptsVec, &cylinder, 100, flags))) return 0;
printf("[*Least-squares cylinder*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr: %lf\n", cylinder.center.vec[0], cylinder.center.vec[1], cylinder.center.vec[2],
cylinder.direction.vec[0], cylinder.direction.vec[1], cylinder.direction.vec[2], cylinder.radius);
//compute minimum-zone cylinder
if (!ALG3D_checkResult(ALG3D_fitCylinderMZ(vec, ptsVec, &cylinder, &absMax, 100, flags))) return 0;
printf("[*Minimum zone cylinder*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr: %lf\nabsMax: %lf\n", cylinder.center.vec[0], cylinder.center.vec[1], cylinder.center.vec[2],
cylinder.direction.vec[0], cylinder.direction.vec[1], cylinder.direction.vec[2], cylinder.radius, absMax);
//compute minimum enclosing cylinder
if (!ALG3D_checkResult(ALG3D_fitCylinderMC(vec, ptsVec, &cylinder, 100))) return 0;
printf("[*Minimum enclosing cylinder*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr: %lf\n", cylinder.center.vec[0], cylinder.center.vec[1], cylinder.center.vec[2],
cylinder.direction.vec[0], cylinder.direction.vec[1], cylinder.direction.vec[2], cylinder.radius);
//compute maximum inscribed cylinder
if (!ALG3D_checkResult(ALG3D_fitCylinderMI(vec2, ptsVec2, &cylinder, 100))) return 0;
printf("[*Maximum inscribed cylinder*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr: %lf\n", cylinder.center.vec[0], cylinder.center.vec[1], cylinder.center.vec[2],
cylinder.direction.vec[0], cylinder.direction.vec[1], cylinder.direction.vec[2], cylinder.radius);
Parameters
[in]ptsInpointer to begin of point data
[in]numInnumber of points
[out]cylindercylinder object
[in]iterationsmaximum number of iterations (e.g. 100 is sufficient in most cases)
Returns
ERR_NONE on success.
Note
The point set should represent much more than half of the cylinder's surface, otherwise there may be infinite solutions that are maximal and inside.

◆ ALG3D_fitCylinderMZ()

int ALG3D_fitCylinderMZ ( const ALG3D_Point3d ptsIn,
size_t  numIn,
ALG3D_Cylinder cylinder,
double *  absMax,
int  iterations,
int  flags 
)

#include <approximate.cpp>

Computes best-fit cylinder using the Chebyshev norm.

This function computes the best-fit cylinder by minimizing the Chebyshev norm (absolute maximum deviation). The resulting cylinder can be used to analyzed the cylindricity of the data points, which is measured as the minimal width between two concentric cylinders that enclose the data points. The maximum absolute deviation which is half this width is written to the variable absMax.

code example
const char* fname = "./sample_data/cylinder.xyz";
const char* fname2 = "./sample_data/pipe.xyz";
ALG3D_Point3d *vec = 0, *vec2 = 0; //3d point data
size_t ptsVec = 0, ptsVec2 = 0; //number of 3d points
ALG3D_Cylinder cylinder; //cylinder object
double absMax = 0; //absolute maximum deviaton (for minimum zone-approximation)
int flags = 0; //parameter flags
//read files
if (!ALG3D_checkResult(ALG3D_readXYZ(fname, &vec, &ptsVec))) return 0;
if (!ALG3D_checkResult(ALG3D_readXYZ(fname2, &vec2, &ptsVec2))) return 0;
//initialize cylinder with roughly known values
cylinder.radius = 17.0; //expected value or guess
flags |= enVarRadius; //set expected value to be changeable (if the expected value should not be changed, use 'enFixRadius' instead)
//compute least-squares cylinder
if (!ALG3D_checkResult(ALG3D_fitCylinder(vec, ptsVec, &cylinder, 100, flags))) return 0;
printf("[*Least-squares cylinder*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr: %lf\n", cylinder.center.vec[0], cylinder.center.vec[1], cylinder.center.vec[2],
cylinder.direction.vec[0], cylinder.direction.vec[1], cylinder.direction.vec[2], cylinder.radius);
//compute minimum-zone cylinder
if (!ALG3D_checkResult(ALG3D_fitCylinderMZ(vec, ptsVec, &cylinder, &absMax, 100, flags))) return 0;
printf("[*Minimum zone cylinder*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr: %lf\nabsMax: %lf\n", cylinder.center.vec[0], cylinder.center.vec[1], cylinder.center.vec[2],
cylinder.direction.vec[0], cylinder.direction.vec[1], cylinder.direction.vec[2], cylinder.radius, absMax);
//compute minimum enclosing cylinder
if (!ALG3D_checkResult(ALG3D_fitCylinderMC(vec, ptsVec, &cylinder, 100))) return 0;
printf("[*Minimum enclosing cylinder*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr: %lf\n", cylinder.center.vec[0], cylinder.center.vec[1], cylinder.center.vec[2],
cylinder.direction.vec[0], cylinder.direction.vec[1], cylinder.direction.vec[2], cylinder.radius);
//compute maximum inscribed cylinder
if (!ALG3D_checkResult(ALG3D_fitCylinderMI(vec2, ptsVec2, &cylinder, 100))) return 0;
printf("[*Maximum inscribed cylinder*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr: %lf\n", cylinder.center.vec[0], cylinder.center.vec[1], cylinder.center.vec[2],
cylinder.direction.vec[0], cylinder.direction.vec[1], cylinder.direction.vec[2], cylinder.radius);
Parameters
[in]ptsInpointer to begin of point data
[in]numInnumber of points
[in,out]cylindercylinder object
[out]absMaxcontains the absolute maximum deviation afterwards
[in]iterationsmaximum number of iterations (e.g. 100 is sufficient in most cases)
[in]flagsinitialization flags
Returns
ERR_NONE on success. ERR_NUM_POINTS if less than 6 points are provided.

◆ ALG3D_fitEllipse()

int ALG3D_fitEllipse ( const ALG3D_Point3d ptsIn,
size_t  numIn,
ALG3D_Ellipse ellipse,
int  iterations,
int  flags 
)

#include <approximate.cpp>

Computes best-fit ellipse in 3D with least-squares.

code example
const double PI = 3.1415926535;
const char* fname = "./sample_data/ellipse3d.xyz";
ALG3D_Point3d *vec = 0; //3d point data
size_t ptsVec = 0; //number of 3d points
ALG3D_Ellipse ellipse; //ellipse object
int flags = 0; //parameter flags
//read file
if (!ALG3D_checkResult(ALG3D_readXYZ(fname, &vec, &ptsVec))) return 0;
//compute least-squares 3d ellipse
if (!ALG3D_checkResult(ALG3D_fitEllipse(vec, ptsVec, &ellipse, 100, flags))) return 0;
printf("pt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\ta: %lf\tb: %lf\tphi: %lf deg\n", ellipse.center.vec[0], ellipse.center.vec[1], ellipse.center.vec[2],
ellipse.direction.vec[0], ellipse.direction.vec[1], ellipse.direction.vec[2],
ellipse.a, ellipse.b, ellipse.angle * 180 / PI);
Parameters
ptsInpointer to begin of point data
numInnumber of points
ellipseellipse object
iterationsmaximum number of iterations (e.g. 100 is sufficient in most cases)
flagsinitialization flags
Returns
ERR_NONE on success. ERR_NUM_POINTS if less than 6 points are provided.

◆ ALG3D_fitLine()

int ALG3D_fitLine ( const ALG3D_Point3d ptsIn,
size_t  numIn,
ALG3D_Line line,
int  iterations,
int  flags 
)

#include <approximate.cpp>

Computes best-fit line in 3D with least-squares.

code example
const char* fname = "./sample_data/line.xyz";
ALG3D_Point3d *vec = 0; //3d point data
size_t ptsVec = 0; //number of 3d points
ALG3D_Line line; //line object
double absMax = 0; //absolute maximum deviaton (for minimum zone-approximation)
int flags = 0; //parameter flags
//read file
if (!ALG3D_checkResult(ALG3D_readXYZ(fname, &vec, &ptsVec))) return 0;
//compute least-squares line
if (!ALG3D_checkResult(ALG3D_fitLine(vec, ptsVec, &line, 100, flags))) return 0;
printf("[*Least-squares line*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\n", line.center.vec[0], line.center.vec[1], line.center.vec[2],
line.direction.vec[0], line.direction.vec[1], line.direction.vec[2]);
//compute minimum-zone line
if (!ALG3D_checkResult(ALG3D_fitLineMZ(vec, ptsVec, &line, &absMax, 100))) return 0;
printf("[*Minimum zone line*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\nabsMax: %lf\n", line.center.vec[0], line.center.vec[1], line.center.vec[2],
line.direction.vec[0], line.direction.vec[1], line.direction.vec[2], absMax);
Parameters
ptsInpointer to begin of point data
numInnumber of points
lineLine object
iterationsmaximum number of iterations (e.g. 100 is sufficient in most cases)
flagsflags. This parameter only exists for interface compatibility and is ignored in the case of a line, because no parameter can be fixed.
Returns
ERR_NONE on success. ERR_NUM_POINTS if less than 2 points are provided.

◆ ALG3D_fitLineMZ()

int ALG3D_fitLineMZ ( const ALG3D_Point3d ptsIn,
size_t  numIn,
ALG3D_Line line,
double *  absMax,
int  iterations 
)

#include <approximate.cpp>

Computes best-fit 3D line using the Chebyshev norm.

This function computes the best-fit 3D line by minimizing the Chebyshev norm (absolute maximum deviation). The resulting line can be used to analyzed the straightness / linearity of the data points, which is measured as the radius of minimal cylinder that encloses the data points. The maximum absolute deviation is written to the variable absMax.

code example
const char* fname = "./sample_data/line.xyz";
ALG3D_Point3d *vec = 0; //3d point data
size_t ptsVec = 0; //number of 3d points
ALG3D_Line line; //line object
double absMax = 0; //absolute maximum deviaton (for minimum zone-approximation)
int flags = 0; //parameter flags
//read file
if (!ALG3D_checkResult(ALG3D_readXYZ(fname, &vec, &ptsVec))) return 0;
//compute least-squares line
if (!ALG3D_checkResult(ALG3D_fitLine(vec, ptsVec, &line, 100, flags))) return 0;
printf("[*Least-squares line*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\n", line.center.vec[0], line.center.vec[1], line.center.vec[2],
line.direction.vec[0], line.direction.vec[1], line.direction.vec[2]);
//compute minimum-zone line
if (!ALG3D_checkResult(ALG3D_fitLineMZ(vec, ptsVec, &line, &absMax, 100))) return 0;
printf("[*Minimum zone line*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\nabsMax: %lf\n", line.center.vec[0], line.center.vec[1], line.center.vec[2],
line.direction.vec[0], line.direction.vec[1], line.direction.vec[2], absMax);
Parameters
[in]ptsInpointer to begin of point data
[in]numInnumber of points
[out]lineline object
[out]absMaxcontains the absolute maximum deviation afterwards
[in]iterationsmaximum number of iterations (e.g. 100 is sufficient in most cases)
Returns
ERR_NONE on success. ERR_NUM_POINTS if less than 2 points are provided.

◆ ALG3D_fitPlane()

int ALG3D_fitPlane ( const ALG3D_Point3d ptsIn,
size_t  numIn,
ALG3D_Plane plane,
int  iterations,
int  flags 
)

#include <approximate.cpp>

Computes best-fit plane with least-squares.

code example
const char* fname = "./sample_data/plane.xyz";
ALG3D_Point3d *vec = 0; //pointer to array of 3d points
size_t ptsVec = 0; //number of 3d points in array
ALG3D_Plane plane, plane2; //plane objects
double absMax = 0; //absolute maximum deviaton (for minimum zone-approximation)
double thickness = 0; //thickness (distance between MC planes)
int flags = 0; //parameter flags
//read file
if (!ALG3D_checkResult(ALG3D_readXYZ(fname, &vec, &ptsVec))) return 0;
//compute least-squares plane
if (!ALG3D_checkResult(ALG3D_fitPlane(vec, ptsVec, &plane, 100, flags))) return 0;
printf("[*Least-squares plane*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\n", plane.center.vec[0], plane.center.vec[1], plane.center.vec[2],
plane.direction.vec[0], plane.direction.vec[1], plane.direction.vec[2]);
//compute minimum-zone plane
if (!ALG3D_checkResult(ALG3D_fitPlaneMZ(vec, ptsVec, &plane, &absMax, 100))) return 0;
printf("[*Minimum zone plane*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\nabsMax: %lf", plane.center.vec[0], plane.center.vec[1], plane.center.vec[2],
plane.direction.vec[0], plane.direction.vec[1], plane.direction.vec[2], absMax);
//compute minimum enclosing planes
if (!ALG3D_checkResult(ALG3D_fitPlanesMC(vec, ptsVec, &plane, &plane2, &thickness, 100))) return 0;
printf("[*Minimum enclosing planes*]\n");
printf(" [MC plane 1]pt: %lf, %lf, %lf\tdir : %lf, %lf, %lf\n", plane.center.vec[0], plane.center.vec[1], plane.center.vec[2],
plane.direction.vec[0], plane.direction.vec[1], plane.direction.vec[2]);
printf(" [MC plane 2]pt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\n", plane2.center.vec[0], plane2.center.vec[1], plane2.center.vec[2],
plane2.direction.vec[0], plane2.direction.vec[1], plane2.direction.vec[2]);
printf(" thickness: %lf\n", thickness);
Parameters
ptsInpointer to begin of point data
numInnumber of points
planeplane object
iterationsmaximum number of iterations (e.g. 100 is sufficient in most cases)
flagsinitialization flags
Returns
ERR_NONE on success. ERR_NUM_POINTS if less than 3 points are provided.

◆ ALG3D_fitPlaneMZ()

int ALG3D_fitPlaneMZ ( const ALG3D_Point3d ptsIn,
size_t  numIn,
ALG3D_Plane plane,
double *  absMax,
int  iterations 
)

#include <approximate.cpp>

Computes best-fit plane using the Chebyshev norm.

This function computes the best-fit plane by minimizing the Chebyshev norm (absolute maximum deviation). The resulting plane can be used to analyzed the planarity of the data points, which is measured as the minimal width between two parallel planes that enclose the data points. The maximum absolute deviation which is half this width is written to the variable absMax.

code example
const char* fname = "./sample_data/plane.xyz";
ALG3D_Point3d *vec = 0; //pointer to array of 3d points
size_t ptsVec = 0; //number of 3d points in array
ALG3D_Plane plane, plane2; //plane objects
double absMax = 0; //absolute maximum deviaton (for minimum zone-approximation)
double thickness = 0; //thickness (distance between MC planes)
int flags = 0; //parameter flags
//read file
if (!ALG3D_checkResult(ALG3D_readXYZ(fname, &vec, &ptsVec))) return 0;
//compute least-squares plane
if (!ALG3D_checkResult(ALG3D_fitPlane(vec, ptsVec, &plane, 100, flags))) return 0;
printf("[*Least-squares plane*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\n", plane.center.vec[0], plane.center.vec[1], plane.center.vec[2],
plane.direction.vec[0], plane.direction.vec[1], plane.direction.vec[2]);
//compute minimum-zone plane
if (!ALG3D_checkResult(ALG3D_fitPlaneMZ(vec, ptsVec, &plane, &absMax, 100))) return 0;
printf("[*Minimum zone plane*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\nabsMax: %lf", plane.center.vec[0], plane.center.vec[1], plane.center.vec[2],
plane.direction.vec[0], plane.direction.vec[1], plane.direction.vec[2], absMax);
//compute minimum enclosing planes
if (!ALG3D_checkResult(ALG3D_fitPlanesMC(vec, ptsVec, &plane, &plane2, &thickness, 100))) return 0;
printf("[*Minimum enclosing planes*]\n");
printf(" [MC plane 1]pt: %lf, %lf, %lf\tdir : %lf, %lf, %lf\n", plane.center.vec[0], plane.center.vec[1], plane.center.vec[2],
plane.direction.vec[0], plane.direction.vec[1], plane.direction.vec[2]);
printf(" [MC plane 2]pt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\n", plane2.center.vec[0], plane2.center.vec[1], plane2.center.vec[2],
plane2.direction.vec[0], plane2.direction.vec[1], plane2.direction.vec[2]);
printf(" thickness: %lf\n", thickness);
Parameters
[in]ptsInpointer to begin of point data
[in]numInnumber of points
[out]planeplane object
[out]absMaxcontains the absolute maximum deviation afterwards
[in]iterationsmaximum number of iterations (e.g. 100 is sufficient in most cases)
Returns
ERR_NONE on success. ERR_NUM_POINTS if less than 3 points are provided.

◆ ALG3D_fitPlanesMC()

int ALG3D_fitPlanesMC ( const ALG3D_Point3d ptsIn,
size_t  numIn,
ALG3D_Plane plane1,
ALG3D_Plane plane2,
double *  thickness,
int  iterations 
)

#include <approximate.cpp>

Computes the minimum enclosing planes.

This function computes the minimum enclosing planes, i.e. two minimum distance parallel planes that contains all points. This function can be used to compute the thickness of a point set. It is equivalent to MZ plane.

code example
const char* fname = "./sample_data/plane.xyz";
ALG3D_Point3d *vec = 0; //pointer to array of 3d points
size_t ptsVec = 0; //number of 3d points in array
ALG3D_Plane plane, plane2; //plane objects
double absMax = 0; //absolute maximum deviaton (for minimum zone-approximation)
double thickness = 0; //thickness (distance between MC planes)
int flags = 0; //parameter flags
//read file
if (!ALG3D_checkResult(ALG3D_readXYZ(fname, &vec, &ptsVec))) return 0;
//compute least-squares plane
if (!ALG3D_checkResult(ALG3D_fitPlane(vec, ptsVec, &plane, 100, flags))) return 0;
printf("[*Least-squares plane*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\n", plane.center.vec[0], plane.center.vec[1], plane.center.vec[2],
plane.direction.vec[0], plane.direction.vec[1], plane.direction.vec[2]);
//compute minimum-zone plane
if (!ALG3D_checkResult(ALG3D_fitPlaneMZ(vec, ptsVec, &plane, &absMax, 100))) return 0;
printf("[*Minimum zone plane*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\nabsMax: %lf", plane.center.vec[0], plane.center.vec[1], plane.center.vec[2],
plane.direction.vec[0], plane.direction.vec[1], plane.direction.vec[2], absMax);
//compute minimum enclosing planes
if (!ALG3D_checkResult(ALG3D_fitPlanesMC(vec, ptsVec, &plane, &plane2, &thickness, 100))) return 0;
printf("[*Minimum enclosing planes*]\n");
printf(" [MC plane 1]pt: %lf, %lf, %lf\tdir : %lf, %lf, %lf\n", plane.center.vec[0], plane.center.vec[1], plane.center.vec[2],
plane.direction.vec[0], plane.direction.vec[1], plane.direction.vec[2]);
printf(" [MC plane 2]pt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\n", plane2.center.vec[0], plane2.center.vec[1], plane2.center.vec[2],
plane2.direction.vec[0], plane2.direction.vec[1], plane2.direction.vec[2]);
printf(" thickness: %lf\n", thickness);
Parameters
[in]ptsInpointer to begin of point data
[in]numInnumber of points
[out]plane1,plane2two plane objects, their normal vectors will point towards each other
[out]thicknessthickness as the distance between the two planes
[in]iterationsmaximum number of iterations (e.g. 100 is sufficient in most cases)
Returns
ERR_NONE on success.

◆ ALG3D_fitSphere()

int ALG3D_fitSphere ( const ALG3D_Point3d ptsIn,
size_t  numIn,
ALG3D_Sphere sphere,
int  iterations,
int  flags 
)

#include <approximate.cpp>

Computes best-fit sphere with least-squares.

code example
const char* fname = "./sample_data/sphere.xyz";
const char* fname2 = "./sample_data/sphere_full.xyz";
ALG3D_Point3d *vec = 0, *vec2 = 0; //3d point data
size_t ptsVec = 0, ptsVec2 = 0; //number of 3d points
ALG3D_Sphere sphere; //sphere object
double absMax = 0; //absolute maximum deviaton (for minimum zone-approximation)
int flags = 0; //parameter flags
//read files
if (!ALG3D_checkResult(ALG3D_readXYZ(fname, &vec, &ptsVec))) return 0;
if (!ALG3D_checkResult(ALG3D_readXYZ(fname2, &vec2, &ptsVec2))) return 0;
//compute least-squares sphere
if (!ALG3D_checkResult(ALG3D_fitSphere(vec, ptsVec, &sphere, 100, flags))) return 0;
printf("[*Least-squares sphere*]\npt: %lf,%lf,%lf\tr: %lf\n", sphere.center.vec[0], sphere.center.vec[1], sphere.center.vec[2], sphere.radius);
//compute minimum-zone sphere
if (!ALG3D_checkResult(ALG3D_fitSphereMZ(vec, ptsVec, &sphere, &absMax, 100, flags))) return 0;
printf("[*Minimum zone sphere*]\npt: %lf,%lf,%lf\tr: %lf\nabsMax: %lf\n", sphere.center.vec[0], sphere.center.vec[1], sphere.center.vec[2], sphere.radius, absMax);
//compute minimum enclosing sphere
if (!ALG3D_checkResult(ALG3D_fitSphereMC(vec, ptsVec, &sphere, 100))) return 0;
printf("[*Minimum enclosing sphere*]\npt: %lf,%lf,%lf\tr: %lf\n", sphere.center.vec[0], sphere.center.vec[1], sphere.center.vec[2], sphere.radius);
//compute maximum inscribed sphere
if (!ALG3D_checkResult(ALG3D_fitSphereMI(vec2, ptsVec2, &sphere, 100))) return 0;
printf("[*Maximum inscribed sphere*]\npt: %lf,%lf,%lf\tr: %lf\n", sphere.center.vec[0], sphere.center.vec[1], sphere.center.vec[2], sphere.radius);
Parameters
ptsInpointer to begin of point data
numInnumber of points
spheresphere object
iterationsmaximum number of iterations (e.g. 100 is sufficient in most cases)
flagsinitialization flags
Returns
ERR_NONE on success. ERR_NUM_POINTS if less than 4 points are provided.

◆ ALG3D_fitSphereMC()

int ALG3D_fitSphereMC ( const ALG3D_Point3d ptsIn,
size_t  numIn,
ALG3D_Sphere sphere,
int  iterations 
)

#include <approximate.cpp>

Computes the minimum enclosing sphere.

This function computes the minimum enclosing sphere, i.e. the smallest sphere that contains all points.

code example
const char* fname = "./sample_data/sphere.xyz";
const char* fname2 = "./sample_data/sphere_full.xyz";
ALG3D_Point3d *vec = 0, *vec2 = 0; //3d point data
size_t ptsVec = 0, ptsVec2 = 0; //number of 3d points
ALG3D_Sphere sphere; //sphere object
double absMax = 0; //absolute maximum deviaton (for minimum zone-approximation)
int flags = 0; //parameter flags
//read files
if (!ALG3D_checkResult(ALG3D_readXYZ(fname, &vec, &ptsVec))) return 0;
if (!ALG3D_checkResult(ALG3D_readXYZ(fname2, &vec2, &ptsVec2))) return 0;
//compute least-squares sphere
if (!ALG3D_checkResult(ALG3D_fitSphere(vec, ptsVec, &sphere, 100, flags))) return 0;
printf("[*Least-squares sphere*]\npt: %lf,%lf,%lf\tr: %lf\n", sphere.center.vec[0], sphere.center.vec[1], sphere.center.vec[2], sphere.radius);
//compute minimum-zone sphere
if (!ALG3D_checkResult(ALG3D_fitSphereMZ(vec, ptsVec, &sphere, &absMax, 100, flags))) return 0;
printf("[*Minimum zone sphere*]\npt: %lf,%lf,%lf\tr: %lf\nabsMax: %lf\n", sphere.center.vec[0], sphere.center.vec[1], sphere.center.vec[2], sphere.radius, absMax);
//compute minimum enclosing sphere
if (!ALG3D_checkResult(ALG3D_fitSphereMC(vec, ptsVec, &sphere, 100))) return 0;
printf("[*Minimum enclosing sphere*]\npt: %lf,%lf,%lf\tr: %lf\n", sphere.center.vec[0], sphere.center.vec[1], sphere.center.vec[2], sphere.radius);
//compute maximum inscribed sphere
if (!ALG3D_checkResult(ALG3D_fitSphereMI(vec2, ptsVec2, &sphere, 100))) return 0;
printf("[*Maximum inscribed sphere*]\npt: %lf,%lf,%lf\tr: %lf\n", sphere.center.vec[0], sphere.center.vec[1], sphere.center.vec[2], sphere.radius);
Parameters
[in]ptsInpointer to begin of point data
[in]numInnumber of points
[out]spheresphere object
[in]iterationsmaximum number of iterations (e.g. 100 is sufficient in most cases)
Returns
ERR_NONE on success.

◆ ALG3D_fitSphereMI()

int ALG3D_fitSphereMI ( const ALG3D_Point3d ptsIn,
size_t  numIn,
ALG3D_Sphere sphere,
int  iterations 
)

#include <approximate.cpp>

Computes the maximum inscribed sphere (Pferch sphere).

This function computes the maximum inscribed sphere, i.e. the largest sphere that fits inside the given points.

code example
const char* fname = "./sample_data/sphere.xyz";
const char* fname2 = "./sample_data/sphere_full.xyz";
ALG3D_Point3d *vec = 0, *vec2 = 0; //3d point data
size_t ptsVec = 0, ptsVec2 = 0; //number of 3d points
ALG3D_Sphere sphere; //sphere object
double absMax = 0; //absolute maximum deviaton (for minimum zone-approximation)
int flags = 0; //parameter flags
//read files
if (!ALG3D_checkResult(ALG3D_readXYZ(fname, &vec, &ptsVec))) return 0;
if (!ALG3D_checkResult(ALG3D_readXYZ(fname2, &vec2, &ptsVec2))) return 0;
//compute least-squares sphere
if (!ALG3D_checkResult(ALG3D_fitSphere(vec, ptsVec, &sphere, 100, flags))) return 0;
printf("[*Least-squares sphere*]\npt: %lf,%lf,%lf\tr: %lf\n", sphere.center.vec[0], sphere.center.vec[1], sphere.center.vec[2], sphere.radius);
//compute minimum-zone sphere
if (!ALG3D_checkResult(ALG3D_fitSphereMZ(vec, ptsVec, &sphere, &absMax, 100, flags))) return 0;
printf("[*Minimum zone sphere*]\npt: %lf,%lf,%lf\tr: %lf\nabsMax: %lf\n", sphere.center.vec[0], sphere.center.vec[1], sphere.center.vec[2], sphere.radius, absMax);
//compute minimum enclosing sphere
if (!ALG3D_checkResult(ALG3D_fitSphereMC(vec, ptsVec, &sphere, 100))) return 0;
printf("[*Minimum enclosing sphere*]\npt: %lf,%lf,%lf\tr: %lf\n", sphere.center.vec[0], sphere.center.vec[1], sphere.center.vec[2], sphere.radius);
//compute maximum inscribed sphere
if (!ALG3D_checkResult(ALG3D_fitSphereMI(vec2, ptsVec2, &sphere, 100))) return 0;
printf("[*Maximum inscribed sphere*]\npt: %lf,%lf,%lf\tr: %lf\n", sphere.center.vec[0], sphere.center.vec[1], sphere.center.vec[2], sphere.radius);
Parameters
[in]ptsInpointer to begin of point data
[in]numInnumber of points
[out]spheresphere object
[in]iterationsmaximum number of iterations (e.g. 100 is sufficient in most cases)
Returns
ERR_NONE on success.
Note
The point set should represent much more than half of the sphere's surface, otherwise there may be infinite solutions that are maximal and inside.

◆ ALG3D_fitSphereMZ()

int ALG3D_fitSphereMZ ( const ALG3D_Point3d ptsIn,
size_t  numIn,
ALG3D_Sphere sphere,
double *  absMax,
int  iterations,
int  flags 
)

#include <approximate.cpp>

Computes best-fit sphere using the Chebyshev norm.

This function computes the best-fit sphere by minimizing the Chebyshev norm (absolute maximum deviation). The maximum absolute deviation is written to the variable absMax.

code example
const char* fname = "./sample_data/sphere.xyz";
const char* fname2 = "./sample_data/sphere_full.xyz";
ALG3D_Point3d *vec = 0, *vec2 = 0; //3d point data
size_t ptsVec = 0, ptsVec2 = 0; //number of 3d points
ALG3D_Sphere sphere; //sphere object
double absMax = 0; //absolute maximum deviaton (for minimum zone-approximation)
int flags = 0; //parameter flags
//read files
if (!ALG3D_checkResult(ALG3D_readXYZ(fname, &vec, &ptsVec))) return 0;
if (!ALG3D_checkResult(ALG3D_readXYZ(fname2, &vec2, &ptsVec2))) return 0;
//compute least-squares sphere
if (!ALG3D_checkResult(ALG3D_fitSphere(vec, ptsVec, &sphere, 100, flags))) return 0;
printf("[*Least-squares sphere*]\npt: %lf,%lf,%lf\tr: %lf\n", sphere.center.vec[0], sphere.center.vec[1], sphere.center.vec[2], sphere.radius);
//compute minimum-zone sphere
if (!ALG3D_checkResult(ALG3D_fitSphereMZ(vec, ptsVec, &sphere, &absMax, 100, flags))) return 0;
printf("[*Minimum zone sphere*]\npt: %lf,%lf,%lf\tr: %lf\nabsMax: %lf\n", sphere.center.vec[0], sphere.center.vec[1], sphere.center.vec[2], sphere.radius, absMax);
//compute minimum enclosing sphere
if (!ALG3D_checkResult(ALG3D_fitSphereMC(vec, ptsVec, &sphere, 100))) return 0;
printf("[*Minimum enclosing sphere*]\npt: %lf,%lf,%lf\tr: %lf\n", sphere.center.vec[0], sphere.center.vec[1], sphere.center.vec[2], sphere.radius);
//compute maximum inscribed sphere
if (!ALG3D_checkResult(ALG3D_fitSphereMI(vec2, ptsVec2, &sphere, 100))) return 0;
printf("[*Maximum inscribed sphere*]\npt: %lf,%lf,%lf\tr: %lf\n", sphere.center.vec[0], sphere.center.vec[1], sphere.center.vec[2], sphere.radius);
Parameters
[in]ptsInpointer to begin of point data
[in]numInnumber of points
[in,out]spheresphere object
[out]absMaxcontains the absolute maximum deviation afterwards
[in]iterationsmaximum number of iterations (e.g. 100 is sufficient in most cases)
[in]flagsinitialization flags
Returns
ERR_NONE on success. ERR_NUM_POINTS if less than 4 points are provided.

◆ ALG3D_fitSuperellipse()

int ALG3D_fitSuperellipse ( const ALG3D_Point3d ptsIn,
size_t  numIn,
ALG3D_Superellipse superellipse,
int  iterations,
int  flags 
)

#include <approximate.cpp>

Computes best-fit superellipse with least-squares.

code example
const double PI = 3.1415926535;
const char* fname = "./sample_data/superellipse.xyz";
ALG3D_Point3d *vec = 0; //3d point data
size_t ptsVec = 0; //number of 3d points
ALG3D_Superellipse ellipse; //superellipse object
int flags = 0; //parameter flags
//read file
if (!ALG3D_checkResult(ALG3D_readXYZ(fname, &vec, &ptsVec))) return 0;
//compute least-squares superellipse
if (!ALG3D_checkResult(ALG3D_fitSuperellipse(vec, ptsVec, &ellipse, 100, flags))) return 0;
printf("pt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\na: %lf\tb: %lf\tepsilon: %lf\tangle: %lf deg\n", ellipse.center.vec[0], ellipse.center.vec[1], ellipse.center.vec[2],
ellipse.direction.vec[0], ellipse.direction.vec[1], ellipse.direction.vec[2],
ellipse.a, ellipse.b, ellipse.epsilon, ellipse.angle * 180 / PI);
Parameters
ptsInpointer to begin of point data
numInnumber of points
superellipsesuperellipse object
iterationsmaximum number of iterations (e.g. 100 is sufficient in most cases)
flagsinitialization flags
Returns
ERR_NONE on success. ERR_NUM_POINTS if less than 7 points are provided.

◆ ALG3D_fitTorus()

int ALG3D_fitTorus ( const ALG3D_Point3d ptsIn,
size_t  numIn,
ALG3D_Torus torus,
int  iterations,
int  flags 
)

#include <approximate.cpp>

Computes best-fit torus with least-squares.

code example
const char* fname = "./sample_data/torus.xyz";
ALG3D_Point3d *vec = 0; //3d point data
size_t ptsVec = 0; //number of 3d points
double absMax = 0; //absolute maximum deviaton (for minimum zone-approximation)
ALG3D_Torus torus; //cone object
int flags = 0; //parameter flags
//read file
if (!ALG3D_checkResult(ALG3D_readXYZ(fname, &vec, &ptsVec))) return 0;
//compute least-squares torus
if (!ALG3D_checkResult(ALG3D_fitTorus(vec, ptsVec, &torus, 100, flags))) return 0;
printf("[*Least-squares torus*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr1: %lf\tr2: %lf\n", torus.center.vec[0], torus.center.vec[1], torus.center.vec[2],
torus.direction.vec[0], torus.direction.vec[1], torus.direction.vec[2], torus.radius1, torus.radius2);
//compute minimum zone torus
if (!ALG3D_checkResult(ALG3D_fitTorusMZ(vec, ptsVec, &torus, &absMax, 100, flags))) return 0;
printf("[*Minimum zone torus*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr1: %lf\tr2: %lf\nabsMax: %lf\n", torus.center.vec[0], torus.center.vec[1], torus.center.vec[2],
torus.direction.vec[0], torus.direction.vec[1], torus.direction.vec[2], torus.radius1, torus.radius2, absMax);
Parameters
ptsInpointer to begin of point data
numInnumber of points
torustorus object
iterationsmaximum number of iterations (e.g. 100 is sufficient in most cases)
flagsinitialization flags.
Returns
ERR_NONE on success. ERR_NUM_POINTS if less than 7 points are provided.

◆ ALG3D_fitTorusMZ()

int ALG3D_fitTorusMZ ( const ALG3D_Point3d ptsIn,
size_t  numIn,
ALG3D_Torus torus,
double *  absMax,
int  iterations,
int  flags 
)

#include <approximate.cpp>

Computes best-fit torus using the Chebyshev norm.

This function computes the best-fit torus by minimizing the Chebyshev norm (absolute maximum deviation). The maximum absolute deviation is written to the variable absMax.

code example
const char* fname = "./sample_data/torus.xyz";
ALG3D_Point3d *vec = 0; //3d point data
size_t ptsVec = 0; //number of 3d points
double absMax = 0; //absolute maximum deviaton (for minimum zone-approximation)
ALG3D_Torus torus; //cone object
int flags = 0; //parameter flags
//read file
if (!ALG3D_checkResult(ALG3D_readXYZ(fname, &vec, &ptsVec))) return 0;
//compute least-squares torus
if (!ALG3D_checkResult(ALG3D_fitTorus(vec, ptsVec, &torus, 100, flags))) return 0;
printf("[*Least-squares torus*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr1: %lf\tr2: %lf\n", torus.center.vec[0], torus.center.vec[1], torus.center.vec[2],
torus.direction.vec[0], torus.direction.vec[1], torus.direction.vec[2], torus.radius1, torus.radius2);
//compute minimum zone torus
if (!ALG3D_checkResult(ALG3D_fitTorusMZ(vec, ptsVec, &torus, &absMax, 100, flags))) return 0;
printf("[*Minimum zone torus*]\npt: %lf,%lf,%lf\tdir: %lf,%lf,%lf\tr1: %lf\tr2: %lf\nabsMax: %lf\n", torus.center.vec[0], torus.center.vec[1], torus.center.vec[2],
torus.direction.vec[0], torus.direction.vec[1], torus.direction.vec[2], torus.radius1, torus.radius2, absMax);
Parameters
[in]ptsInpointer to begin of point data
[in]numInnumber of points
[in,out]torustorus object
[out]absMaxcontains the absolute maximum deviation afterwards
[in]iterationsmaximum number of iterations (e.g. 100 is sufficient in most cases)
[in]flagsinitialization flags
Returns
ERR_NONE on success. ERR_NUM_POINTS if less than 7 points are provided.

◆ ALG3D_intersectLineCone()

int ALG3D_intersectLineCone ( const ALG3D_Line line,
const ALG3D_Cone cone,
ALG3D_Point3d ptsOut,
size_t *  numOut 
)

#include <geometry.cpp>

Computes the intersection between a 3D line and a cone.

Parameters
lineline
conecone
ptsOutexisting storage for max. 2 intersection points
numOutnumber of intersection points (0, 1 or 2)
Returns
ERR_NONE on success. If there are no or infinite intersections ERR_CONDITION is returned.

◆ ALG3D_intersectLineCylinder()

int ALG3D_intersectLineCylinder ( const ALG3D_Line line,
const ALG3D_Cylinder cylinder,
ALG3D_Point3d ptsOut,
size_t *  numOut 
)

#include <geometry.cpp>

Computes the intersection between a 3D line and a cylinder.

Parameters
lineline
cylindercylinder
ptsOutexisting storage for max. 2 intersection points
numOutnumber of intersection points (0, 1 or 2)
Returns
ERR_NONE on success. If there are no or infinite intersections ERR_CONDITION is returned.

◆ ALG3D_intersectLineLine()

int ALG3D_intersectLineLine ( const ALG3D_Line line1,
const ALG3D_Line line2,
ALG3D_Point3d point 
)

#include <geometry.cpp>

Computes the intersection between two 3D lines. If the two lines are not parallel an intersection point will be computed. This point is computed from the position where the lines are closest. The middle point between the closest points on the lines is returned as the intersection point.

Parameters
line1,line2two lines
pointintersection point (if there is no intersection, e.g. parallel lines, ERR_CONDITION is returned)
Returns
ERR_NONE on success

◆ ALG3D_intersectLinePlane()

int ALG3D_intersectLinePlane ( const ALG3D_Line line,
const ALG3D_Plane plane,
ALG3D_Point3d point 
)

#include <geometry.cpp>

Computes the intersection between 3D line and a plane.

Parameters
lineline
planeplane
pointintersection point (if there is no or infinite intersections, e.g. line is perpendicular to plane normal, ERR_CONDITION is returned)
Returns
ERR_NONE on success

◆ ALG3D_intersectLineSphere()

int ALG3D_intersectLineSphere ( const ALG3D_Line line,
const ALG3D_Sphere sphere,
ALG3D_Point3d ptsOut,
size_t *  numOut 
)

#include <geometry.cpp>

Computes the intersection between 3D line and a sphere.

Parameters
lineline
spheresphere
ptsOutexisting storage for max. 2 intersection points
numOutnumber of intersection points (0, 1 or 2)
Returns
ERR_NONE on success. If there are no or infinite intersections ERR_CONDITION is returned.

◆ ALG3D_intersectPlanePlane()

int ALG3D_intersectPlanePlane ( const ALG3D_Plane plane1,
const ALG3D_Plane plane2,
ALG3D_Line line 
)

#include <geometry.cpp>

Computes the intersection between two planes.

Parameters
plane1,plane2two planes
lineintersection line. If there are no or infinite intersections, e.g. parallel planes, ERR_CONDITION is returned.
Returns
ERR_NONE on success

◆ ALG3D_intersectPlaneSphere()

int ALG3D_intersectPlaneSphere ( const ALG3D_Plane plane,
const ALG3D_Sphere sphere,
ALG3D_Circle circle 
)

#include <geometry.cpp>

Computes the intersection between a plane and a sphere.

Parameters
planeplane
spheresphere
circleresulting intersection circle, if one exists. If the plane is tangent to sphere the resulting circle radius is 0.
Returns
ERR_NONE on success. If there is no intersection, ERR_CONDITION is returned.

◆ ALG3D_projectPoints2Cone()

int ALG3D_projectPoints2Cone ( const ALG3D_Point3d ptsIn,
size_t  numIn,
const ALG3D_Cone cone,
ALG3D_Point3d ptsOut 
)

#include <geometry.cpp>

Computes the projections of points to a cone.

Parameters
ptsInpointer to begin of point data
numInnumber of points
conecone object
ptsOutpointer to existing storage for the resulting projected points (must have the size of numIn).
Note
This function always uses all available threads.
Returns
ERR_NONE on success

◆ ALG3D_projectPoints2Cylinder()

int ALG3D_projectPoints2Cylinder ( const ALG3D_Point3d ptsIn,
size_t  numIn,
const ALG3D_Cylinder cylinder,
ALG3D_Point3d ptsOut 
)

#include <geometry.cpp>

Computes the projections of points to a cylinder.

Parameters
ptsInpointer to begin of point data
numInnumber of points
cylindercylinder object
ptsOutpointer to existing storage for the resulting projected points (must have the size of numIn).
Note
This function always uses all available threads.
Returns
ERR_NONE on success

◆ ALG3D_projectPoints2Line()

int ALG3D_projectPoints2Line ( const ALG3D_Point3d ptsIn,
size_t  numIn,
const ALG3D_Line line,
ALG3D_Point3d ptsOut 
)

#include <geometry.cpp>

Computes the projections of points to a 3D line.

Parameters
ptsInpointer to begin of point data
numInnumber of points
lineline object
ptsOutpointer to existing storage for the resulting projected points (must have the size of numIn).
Note
This function always uses all available threads.
Returns
ERR_NONE on success

◆ ALG3D_projectPoints2Plane()

int ALG3D_projectPoints2Plane ( const ALG3D_Point3d ptsIn,
size_t  numIn,
const ALG3D_Plane plane,
ALG3D_Point3d ptsOut 
)

#include <geometry.cpp>

Computes the projections of points to a plane.

Parameters
ptsInpointer to begin of point data
numInnumber of points
planeplane object
ptsOutpointer to existing storage for the resulting projected points (must have the size of numIn).
Note
This function always uses all available threads.
Returns
ERR_NONE on success

◆ ALG3D_projectPoints2Sphere()

int ALG3D_projectPoints2Sphere ( const ALG3D_Point3d ptsIn,
size_t  numIn,
const ALG3D_Sphere sphere,
ALG3D_Point3d ptsOut 
)

#include <geometry.cpp>

Computes the projections of points to a sphere.

Parameters
ptsInpointer to begin of point data
numInnumber of points
spheresphere object
ptsOutpointer to existing storage for the resulting projected points (must have the size of numIn).
Note
This function always uses all available threads.
Returns
ERR_NONE on success

◆ ALG3D_projectPoints2Torus()

int ALG3D_projectPoints2Torus ( const ALG3D_Point3d ptsIn,
size_t  numIn,
const ALG3D_Torus torus,
ALG3D_Point3d ptsOut 
)

#include <geometry.cpp>

Computes the projections of points to a torus.

Parameters
ptsInpointer to begin of point data
numInnumber of points
torustorus object
ptsOutpointer to existing storage for the resulting projected points (must have the size of numIn).
Note
This function always uses all available threads.
Returns
ERR_NONE on success
ALG3D_fitCone
int __cdecl ALG3D_fitCone(const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Cone *cone, int iterations, int flags)
Computes best-fit cone with least-squares. Note that the cone radius is measured at the height of the...
Definition: approximate.cpp:842
ALG3D_fitSuperellipse
int __cdecl ALG3D_fitSuperellipse(const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Superellipse *superellipse, int iterations, int flags)
Computes best-fit superellipse with least-squares.
Definition: approximate.cpp:640
ALG3D_fitCylinderMI
int __cdecl ALG3D_fitCylinderMI(const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Cylinder *cylinder, int iterations)
Computes the maximum inscribed cylinder (Pferch cylinder).
Definition: approximate.cpp:1347
ALG3D_Plane
Data type for a plane.
Definition: container.h:22
ALG3D_freeMemory
int __cdecl ALG3D_freeMemory(void *data)
Releases the memory of an array that was allocated internally. It is up to the application to assure ...
Definition: interface.cpp:127
ALG3D_Cone::angle
double angle
Angle of the cone.
Definition: container.h:103
ALG3D_Circle::center
ALG3D_Point3d center
Center of the circle.
Definition: container.h:43
ALG3D_Torus::radius1
double radius1
Radius of the ring.
Definition: container.h:115
ALG3D_dotProduct
double __cdecl ALG3D_dotProduct(const ALG3D_Point3d *vector1, const ALG3D_Point3d *vector2)
Returns the scalar/dot product of two vectors.
Definition: geometry.cpp:69
ALG3D_Ellipse::b
double b
Radii of the ellipse, with a >= b.
Definition: container.h:56
ALG3D_Point3d::vec
double vec[3]
point data, e.g. 0=x,1=y,2=z
Definition: interface.h:41
ALG3D_fitLine
int __cdecl ALG3D_fitLine(const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Line *line, int iterations, int flags)
Computes best-fit line in 3D with least-squares.
Definition: approximate.cpp:422
ALG3D_Circle::direction
ALG3D_Point3d direction
Direction vector of the plane the circle lies in.
Definition: container.h:44
ALG3D_checkResult
int __cdecl ALG3D_checkResult(int code)
Function that checks for error codes and prints out the result. To activate console print out use ALG...
Definition: interface.cpp:38
ALG3D_fitCylinderMZ
int __cdecl ALG3D_fitCylinderMZ(const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Cylinder *cylinder, double *absMax, int iterations, int flags)
Computes best-fit cylinder using the Chebyshev norm.
Definition: approximate.cpp:1034
ALG3D_Point3d
Data type for 3d points/vectors with a user-defined pointer.
Definition: interface.h:32
ALG3D_Ellipse
Data type for a 3D-Ellipse.
Definition: container.h:53
ALG3D_Ellipse::a
double a
Definition: container.h:56
ALG3D_Torus::radius2
double radius2
Radius of the tube.
Definition: container.h:116
ALG3D_fitPlaneMZ
int __cdecl ALG3D_fitPlaneMZ(const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Plane *plane, double *absMax, int iterations)
Computes best-fit plane using the Chebyshev norm.
Definition: approximate.cpp:945
ALG3D_Line
Data type for a 3D-line.
Definition: container.h:32
ALG3D_computePrincipleComponents
int __cdecl ALG3D_computePrincipleComponents(const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Point3d components[3])
Computes the principal components of a 3D point set. The principal components are the directions of t...
Definition: geometry.cpp:992
ALG3D_projectPoints2Line
int __cdecl ALG3D_projectPoints2Line(const ALG3D_Point3d *ptsIn, size_t numIn, const ALG3D_Line *line, ALG3D_Point3d *ptsOut)
Computes the projections of points to a 3D line.
Definition: geometry.cpp:490
enVarRadius
@ enVarRadius
radius (used by circle, sphere, cylinder, torus)
Definition: approximate.h:53
ALG3D_fitSphereMZ
int __cdecl ALG3D_fitSphereMZ(const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Sphere *sphere, double *absMax, int iterations, int flags)
Computes best-fit sphere using the Chebyshev norm.
Definition: approximate.cpp:989
ALG3D_readXYZ
int __cdecl ALG3D_readXYZ(const char *fname, ALG3D_Point3d **data, size_t *numPts)
Reads XYZ-data from a file.
Definition: io.cpp:67
ALG3D_computeConvexHull
int __cdecl ALG3D_computeConvexHull(const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Point3d **vertices, size_t *numVertices, size_t **indices, size_t *numIndices)
Computes the convex hull of a 3D point set.
Definition: geometry.cpp:1070
ALG3D_Torus
Data type for a torus.
Definition: container.h:111
ALG3D_fitLineMZ
int __cdecl ALG3D_fitLineMZ(const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Line *line, double *absMax, int iterations)
Computes best-fit 3D line using the Chebyshev norm.
Definition: approximate.cpp:900
ALG3D_fitTorusMZ
int __cdecl ALG3D_fitTorusMZ(const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Torus *torus, double *absMax, int iterations, int flags)
Computes best-fit torus using the Chebyshev norm.
Definition: approximate.cpp:1125
ALG3D_Cone
Data type for a cone.
Definition: container.h:99
ALG3D_fitSphere
int __cdecl ALG3D_fitSphere(const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Sphere *sphere, int iterations, int flags)
Computes best-fit sphere with least-squares.
Definition: approximate.cpp:748
ALG3D_fitSphereMI
int __cdecl ALG3D_fitSphereMI(const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Sphere *sphere, int iterations)
Computes the maximum inscribed sphere (Pferch sphere).
Definition: approximate.cpp:1306
ALG3D_fitConeMZ
int __cdecl ALG3D_fitConeMZ(const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Cone *cone, double *absMax, int iterations, int flags)
Computes best-fit cone using the Chebyshev norm.
Definition: approximate.cpp:1079
ALG3D_Superellipse
Data type for a 3D-superellipse (Lamé-curve).
Definition: container.h:65
ALG3D_fitSphereMC
int __cdecl ALG3D_fitSphereMC(const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Sphere *sphere, int iterations)
Computes the minimum enclosing sphere.
Definition: approximate.cpp:1172
ALG3D_Cylinder
Data type for a cylinder.
Definition: container.h:88
ALG3D_Circle::radius
double radius
Radius of the circle.
Definition: container.h:45
ALG3D_fitTorus
int __cdecl ALG3D_fitTorus(const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Torus *torus, int iterations, int flags)
Computes best-fit torus with least-squares.
Definition: approximate.cpp:866
ALG3D_fitPlane
int __cdecl ALG3D_fitPlane(const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Plane *plane, int iterations, int flags)
Computes best-fit plane with least-squares.
Definition: approximate.cpp:446
ALG3D_Circle
Data type for a 3D-circle.
Definition: container.h:42
ALG3D_fitPlanesMC
int __cdecl ALG3D_fitPlanesMC(const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Plane *plane1, ALG3D_Plane *plane2, double *thickness, int iterations)
Computes the minimum enclosing planes.
Definition: approximate.cpp:1256
ALG3D_fitCylinderMC
int __cdecl ALG3D_fitCylinderMC(const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Cylinder *cylinder, int iterations)
Computes the minimum enclosing cylinder.
Definition: approximate.cpp:1213
ALG3D_Superellipse::epsilon
double epsilon
Shape factor of the superellipse: rectangle (e<<1), ellipse/circle (e==1), diamond (e==2),...
Definition: container.h:69
ALG3D_Sphere
Data type for a sphere.
Definition: container.h:78
ALG3D_fitCircle
int __cdecl ALG3D_fitCircle(const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Circle *circle, int iterations, int flags)
Computes best-fit circle in 3D with least-squares.
Definition: approximate.cpp:470
ALG3D_fitEllipse
int __cdecl ALG3D_fitEllipse(const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Ellipse *ellipse, int iterations, int flags)
Computes best-fit ellipse in 3D with least-squares.
Definition: approximate.cpp:536
ALG3D_fitCylinder
int __cdecl ALG3D_fitCylinder(const ALG3D_Point3d *ptsIn, size_t numIn, ALG3D_Cylinder *cylinder, int iterations, int flags)
Computes best-fit cylinder with least-squares.
Definition: approximate.cpp:791