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

This module contains functions to align point clouds and triangle meshes. More...

Data Structures

struct  icpParam
 ICP input parameters. More...
 
struct  icpMultiViewParam
 MultiView-ICP input parameters. More...
 
struct  icpResult
 ICP result parameters. More...
 
struct  icpMultiViewTransformations
 MultiView-ICP transformations. More...
 

Enumerations

enum  icpOptions {
  icpRXaxis = 1L << 0, icpRYaxis = 1L << 1, icpRZaxis = 1L << 2, icpTXaxis = 1L << 3,
  icpTYaxis = 1L << 4, icpTZaxis = 1L << 5, icpScaling = 1L << 6, icpStandardDoF = 63,
  icpAllDoF = 127, icpMatchCoG = 1L << 7, icpOptMiniMax = 1L << 8, icpRotCenterIsOrigin = 1L << 9
}
 ICP option flags. More...
 

Functions

int __cdecl ALG3D_computeRegistration (const ALG3D_Point3d *ptsStat, size_t nStat, const ALG3D_Point3d *ptsDyn, size_t nDyn, icpParam prm, icpResult *res, double R[9], double T[3], double S[1])
 Computes the registration/alignment of two 3d point clouds. More...
 
int __cdecl ALG3D_computeRegistrationTree (const ALG3D_KDTREE tree, const ALG3D_Point3d *ptsDyn, size_t nDyn, icpParam prm, icpResult *res, double R[9], double T[3], double S[1])
 Computes the registration/alignment of two 3d point clouds. More...
 
int __cdecl ALG3D_computeRegistrationMesh (const ALG3D_MESHTREE meshTree, ALG3D_Point3d *ptsDyn, size_t nDyn, icpParam prm, icpResult *res, double R[9], double T[3], double S[1])
 Computes the registration/alignment of a3d point clouds to a mesh (e.g. CAD model). More...
 
int __cdecl ALG3D_computeRegistrationCorrespondences (const ALG3D_Point3d *ptsStat, const ALG3D_Point3d *ptsDyn, size_t numPts, int options, double R[9], double T[3], double S[1])
 Computes the registration/alignment of two corresponding sets of 3D points. More...
 
int __cdecl ALG3D_computeRegistrationPlane (const ALG3D_Point3d *ptsStat, const ALG3D_Point3d *normalsStat, size_t nStat, const ALG3D_Point3d *ptsDyn, size_t nDyn, icpParam prm, icpResult *res, double R[9], double T[3], double S[1])
 Computes the registration/alignment of two 3d point clouds. More...
 
int __cdecl ALG3D_computeRegistrationPlaneTree (const ALG3D_KDTREE tree, const ALG3D_Point3d *normalsStat, const ALG3D_Point3d *ptsDyn, size_t nDyn, icpParam prm, icpResult *res, double R[9], double T[3], double S[1])
 Computes the registration/alignment of two 3d point clouds. More...
 
int __cdecl ALG3D_computeMultiViewRegistration (const ALG3D_Point3d **pointclouds, size_t *pointcloudsSizes, size_t numPointclouds, icpMultiViewParam prm, icpResult *res, icpMultiViewTransformations *transformations)
 Computes the registration/alignment of multiple 3d point clouds. More...
 

Detailed Description

This module contains functions to align point clouds and triangle meshes.

This module also provides the functions:


Data Structure Documentation

◆ icpParam

struct icpParam

ICP input parameters.

Examples
testKdTree.c.
Data Fields
size_t iterations maximum number of iterations
int options sets options (see icpOptions). If set to 0 icpStandardDoF is used.
double radius maximum search radius (0=auto)
size_t samples number of sample points that should be used for the alignment (0=all)
double tolerance required minimum distance (iterations stop when this value has been reached)

◆ icpMultiViewParam

struct icpMultiViewParam

MultiView-ICP input parameters.

Data Fields
int fixedIndex index of the dataset that should stay fixed (-1 = none)
size_t iterations maximum number of iterations
double radius maximum search radius (0=auto)
double thinningRadius internal thinning radius
double tolerance required minimum distance (iterations stop when this value has been reached)

◆ icpResult

struct icpResult

ICP result parameters.

Examples
testKdTree.c.
Data Fields
size_t iterations number of iterations needed
double maxDist maximal point distance
double meanDev mean distance between the data sets
double minDist minimal point distance
size_t pairs number of correlating data pairs
double stdDev standard deviation of the point distances

◆ icpMultiViewTransformations

struct icpMultiViewTransformations

MultiView-ICP transformations.

Data Fields
double R[9]
double T[3]

Enumeration Type Documentation

◆ icpOptions

enum icpOptions

#include <algorithm.h>

ICP option flags.

Enumerator
icpRXaxis 

X-Axis (Rotation)

icpRYaxis 

Y-Axis (Rotation)

icpRZaxis 

Z-Axis (Rotation)

icpTXaxis 

X-Axis (Translation)

icpTYaxis 

Y-Axis (Translation)

icpTZaxis 

Z-Axis (Translation)

icpScaling 

Scaling.

icpStandardDoF 

Enables all rotation and translation degrees of freedom (without scaling)

icpAllDoF 

Enables all degrees of freedom.

icpMatchCoG 

enable intial alignment of the center of gravities

icpOptMiniMax 

enable MiniMax/Chebyshev optimization (minimize the absolute maximum distance) instead of least-squares. This option is valid for triangle meshes only.

icpRotCenterIsOrigin 

force the rotation center to stay at the origin (0,0,0)

Function Documentation

◆ ALG3D_computeMultiViewRegistration()

int __cdecl ALG3D_computeMultiViewRegistration ( const ALG3D_Point3d **  pointclouds,
size_t *  pointcloudsSizes,
size_t  numPointclouds,
icpMultiViewParam  prm,
icpResult res,
icpMultiViewTransformations transformations 
)

#include <algorithm.h>

Computes the registration/alignment of multiple 3d point clouds.

Parameters
[in]pointcloudsArray of point cloud pointers (begin of each point cloud)
[in]pointcloudsSizesArray with the sizes of each point cloud
[in]numPointcloudsNumber of point clouds
[in]prmSet of parameters to configure the algorithm
[out]resResult statistics
[out]transformationsArray of result transformations (must be preallocated, size = numPointclouds)
Returns
ERR_NONE on success
Note
Also check the result parameters to find out detailed status information
code example
const char* filenames[] = { "./sample_data/schiller05.xyz",
"./sample_data/schiller06.xyz", "./sample_data/schiller07.xyz",
"./sample_data/schiller08.xyz", "./sample_data/schiller11.xyz",
"./sample_data/schiller12.xyz" };
size_t sizes[6];
size_t numPointClouds = 6;
ALG3D_Point3d* pointclouds[6];
for (size_t i = 0; i < numPointClouds; ++i)
{
if (!ALG3D_checkResult(ALG3D_readXYZ(filenames[i], &pointclouds[i], &sizes[i])))
return 0;
printf("Points in pointcloud %s: %I64d\n", filenames[i], sizes[i]);
}
icpMultiViewTransformations transformations[6];
icpResult res;
prm.fixedIndex = -1; // no fixed dataset
prm.iterations = 200; // max iterations
prm.tolerance = 0.001; // tolerance
prm.radius = 2.5; // max search radius
prm.thinningRadius = 0; // no additional thinning
&sizes[0], numPointClouds, prm, &res, &transformations[0])))
return 0;
printf("minDist: %lf, maxDist: %lf, meanDev: %lf, stdDev: %3.10lf, pairs: %I64d, iterations: %I64d\n",
res.minDist, res.maxDist, res.meanDev, res.stdDev, res.pairs, res.iterations);
printf("** Rs **\n");
for (size_t j = 0; j < numPointClouds; ++j)
{
printf("** %I64d **\n", j);
for (int i = 0; i < 9; ++i)
{
if (i % 3 == 0)
printf("\n");
printf("%lf\t", transformations[j].R[i]);
}
printf("\n");
}
printf("\n** Ts **\n");
for (size_t j = 0; j < numPointClouds; ++j)
{
printf("** %I64d **\n", j);
for (int i = 0; i < 3; ++i)
{
printf("%lf\t", transformations[j].T[i]);
}
printf("\n");
}
for (size_t i = 0; i < numPointClouds; ++i)
{
ALG3D_freeMemory(pointclouds[i]);
}

◆ ALG3D_computeRegistration()

int ALG3D_computeRegistration ( const ALG3D_Point3d ptsStat,
size_t  nStat,
const ALG3D_Point3d ptsDyn,
size_t  nDyn,
icpParam  prm,
icpResult res,
double  R[9],
double  T[3],
double  S[1] 
)

#include <algorithm.h>

Computes the registration/alignment of two 3d point clouds.

A given 'static' point cloud, i.e. a kdTree, is used as a reference. For the second, 'dynamic' point cloud the translation, rotation and scaling are computed iteratively . The iterations stop when a (local) distance minimum is found, i.e. until there is no change anymore or one of the termination criteria has been reached (min. tolerance, max. iterations).

code example
ALG3D_Point3d *vecStat = 0, *vecDyn = 0;
size_t ptsStat = 0, ptsDyn = 0;
double R[9] = { 0, 0, 0, 0, 0, 0, 0, 0, 0 }, T[3] = { 0, 0, 0 }, S[1] = { 0 };
icpResult res;
icpParam prm;
prm.radius = 100; //max. search radius (0=auto)
prm.iterations = 150; //max. number of iterations
prm.tolerance = 0.001;//tolerance
prm.samples = 2000; //maximum number of points used for the registraton
prm.options = icpStandardDoF; //enable translation and rotation only
printf("\n***test_Registration***\n");
//read 3d point data from file
if (!ALG3D_checkResult(ALG3D_readXYZ(fnameStat, &vecStat, &ptsStat))) return 0;
if (!ALG3D_checkResult(ALG3D_readXYZ(fnameDyn, &vecDyn, &ptsDyn))) return 0;
//compute registration
if (!ALG3D_checkResult(ALG3D_computeRegistration(vecStat, ptsStat, vecDyn, ptsDyn, prm, &res, R, T, S))) return 0;
//apply transformation
if (!ALG3D_checkResult(ALG3D_transformPoints(vecDyn, ptsDyn, R, T, S))) return 0;
//write tranformed 3d point data to file
if (!ALG3D_checkResult(ALG3D_writeXYZ(fnameOut, vecDyn, ptsDyn))) return 0;
Parameters
[in]ptsStatpointer to the begin of the static data
[in]nStatnumber of points in static data
[in]ptsDynpointer to the begin of the dynamic data
[in]nDynnumber of points in dynamic data
[in]prmset of parameters to configure the algorithm
[out]resresult statistics
[out]R3x3 rotation matrix (row-major)
[out]T3x1 translation matrix (vector)
[out]S1x1 scaling value
Returns
ERR_NONE on success (at least one termination criteria has been reached), ERR_CONDITION if the computation stopped early (if not enough correlating pairs are found)
Note
Also check the result parameters to find out detailed status information

◆ ALG3D_computeRegistrationCorrespondences()

int ALG3D_computeRegistrationCorrespondences ( const ALG3D_Point3d ptsStat,
const ALG3D_Point3d ptsDyn,
size_t  numPts,
int  options,
double  R[9],
double  T[3],
double  S[1] 
)

#include <algorithm.h>

Computes the registration/alignment of two corresponding sets of 3D points.

Computes the transformation that aligns the points from the ptsDyn array with the points from the ptsStat array. The points represent correspondences (i.e. the first points from ptsStat corresponds to the first points in ptsDyn, etc.). Use this function if you already know the correspondences (e.g. the user selected 3 corresponding points in both point clouds). If the correspondences are unknown you have to use the ALG3D_computeRegistration function.

Precondition
ptsStat and ptsDyn need to have the same size and contain at least 3 points.
Warning
Only some combinations of rotation axes in options are allowed in this function:
  • Enabling all 3 rotation axes is allowed
  • Enabling only 1 rotation axis is allowed
  • Enabling no rotation axis is allowed (translation only)
Enabling just two rotation axes is not allowed. ERR_INVALID_ARG is returned in this case.
Be aware that you can get wrong or surprising results if you have deactivated some translation axes. If the centers of the point clouds do not match the origin, a translation is usually required to achieve a correct result (since rotations are always performed around the global origin).
Examples:
  • All rotations allowed + only X and Y translation allowed ==> this can lead to wrong results
  • Z rotation allowed + only X translation allowed ==> this can lead to wrong results
  • Y rotation allowed + X and Z translation allowed ==> this is valid
code example
double R[9] = { 0, 0, 0, 0, 0, 0, 0, 0, 0 }, T[3] = { 0, 0, 0 }, S[1] = { 0 };
ALG3D_Point3d *vecStat = 0, *vecDyn = 0;
vecStat = (ALG3D_Point3d*)malloc(3 * sizeof(ALG3D_Point3d));
vecDyn = (ALG3D_Point3d*)malloc(3 * sizeof(ALG3D_Point3d));
//Fill with known correspondences
vecStat[0].vec[0] = 1; vecStat[0].vec[1] = 4; vecStat[0].vec[2] = 2;
vecStat[1].vec[0] = 2; vecStat[1].vec[1] = 1; vecStat[1].vec[2] = 2;
vecStat[2].vec[0] = 4; vecStat[2].vec[1] = 3; vecStat[2].vec[2] = 2;
vecDyn[0].vec[0] = 1; vecDyn[0].vec[1] = -1; vecDyn[0].vec[2] = 4;
vecDyn[1].vec[0] = 1; vecDyn[1].vec[1] = 0; vecDyn[1].vec[2] = 1;
vecDyn[2].vec[0] = 1; vecDyn[2].vec[1] = 2; vecDyn[2].vec[2] = 3;
if (!ALG3D_checkResult(ALG3D_computeRegistrationCorrespondences(vecStat, vecDyn, 3, icpStandardDoF, R, T, S))) return 0;
free(vecDyn);
free(vecStat);
Parameters
[in]ptsStatpointer to the begin of the static data
[in]ptsDynpointer to the begin of the dynamic data
[in]numPtsnumber of correspondences
[in]optionsdegrees of freedom flags
[out]R3x3 rotation matrix (row-major)
[out]T3x1 translation matrix (vector)
[out]S1x1 scaling value
Returns
ERR_NONE on success

◆ ALG3D_computeRegistrationMesh()

int ALG3D_computeRegistrationMesh ( const ALG3D_MESHTREE  meshTree,
ALG3D_Point3d ptsDyn,
size_t  nDyn,
icpParam  prm,
icpResult res,
double  R[9],
double  T[3],
double  S[1] 
)

#include <algorithm.h>

Computes the registration/alignment of a3d point clouds to a mesh (e.g. CAD model).

A given mesh, i.e. a CAD model, is used as a reference. For the 'dynamic' point cloud the translation, rotation and scaling are computed iteratively . The iterations stop when a (local) distance minimum is found, i.e. until there is no change anymore or one of the termination criteria has been reached (min. tolerance, max. iterations).

code example
//setting ICP paramers
prm.radius = 0; //max. search radius (0=auto)
prm.iterations= 100; //max. number of iterations
prm.tolerance = 0.001; //tolerance
prm.samples = 2000; //number of points used for registration
prm.options = icpStandardDoF | icpMatchCoG; //enable translation and rotation and matching the centers of gravity
//compute ICP registration
if (!ALG3D_checkResult(ALG3D_computeRegistrationMesh(tree, dataXYZ, ptsXYZ, prm, &res, R, T, S))) return 0;
Parameters
[in]meshTreekdTree object containing the mesh
[in]ptsDynpointer to the begin of the dynamic data
[in]nDynnumber of points in dynamic data
[in]prmset of parameters to configure the algorithm
[out]resresult statistics
[out]R3x3 rotation matrix (row-major)
[out]T3x1 translation matrix (vector)
[out]S1x1 scaling value
Returns
ERR_NONE on success (at least one termination criteria has been reached), ERR_CONDITION if the computation stopped early (if not enough correlating pairs are found)
Note
Also check the result parameters to find out detailed status information

◆ ALG3D_computeRegistrationPlane()

int ALG3D_computeRegistrationPlane ( const ALG3D_Point3d ptsStat,
const ALG3D_Point3d normalsStat,
size_t  nStat,
const ALG3D_Point3d ptsDyn,
size_t  nDyn,
icpParam  prm,
icpResult res,
double  R[9],
double  T[3],
double  S[1] 
)

#include <algorithm.h>

Computes the registration/alignment of two 3d point clouds.

A given 'static' point cloud, i.e. a kdTree, is used as a reference. For the second, 'dynamic' point cloud the translation, rotation and scaling are computed iteratively. The iterations stop when a (local) distance minimum is found, i.e. until there is no change anymore or one of the termination criteria has been reached (min. tolerance, max. iterations).

This algorithm uses the "point-to-plane" metric to compute correspondences between the static and dynamic point cloud. This metric requires the static point cloud to have normal vectors. These must be provided by the user. You can use one of the ALG3D_computeNormalVectors functions to perform this operation or use your own one.

Note
The elements of ptsStat and normalsStat must correspond to each other (e.g. normalsStat[0] = normal of ptsStat[0] etc.).
code example
ALG3D_Point3d *vecStat = 0, *vecDyn = 0;
ALG3D_Point3d *normalsStat = 0;
size_t ptsStat = 0, ptsDyn = 0;
double R[9] = { 0, 0, 0, 0, 0, 0, 0, 0, 0 }, T[3] = { 0, 0, 0 }, S[1] = { 0 };
icpResult res;
icpParam prm;
prm.radius = 100; //max. search radius (0=auto)
prm.iterations = 150; //max. number of iterations
prm.tolerance = 0.001;//tolerance
prm.samples = 2000; //maximum number of points used for the registraton
prm.options = icpStandardDoF; //enable translation and rotation only
printf("\n***test_RegistrationPlane***\n");
//read 3d point data from file
if (!ALG3D_checkResult(ALG3D_readXYZ(fnameStat, &vecStat, &ptsStat))) return 0;
if (!ALG3D_checkResult(ALG3D_readXYZ(fnameDyn, &vecDyn, &ptsDyn))) return 0;
normalsStat = (ALG3D_Point3d*)malloc(ptsStat * sizeof(ALG3D_Point3d));
//compute normal vectors
if (!ALG3D_checkResult(ALG3D_computeNormalVectorsKNN(vecStat, ptsStat, normalsStat, 30))) return 0;
//compute registration
if (!ALG3D_checkResult(ALG3D_computeRegistrationPlane(vecStat, normalsStat, ptsStat, vecDyn, ptsDyn, prm, &res, R, T, S))) return 0;
//apply transformation
if (!ALG3D_checkResult(ALG3D_transformPoints(vecDyn, ptsDyn, R, T, S))) return 0;
//write tranformed 3d point data to file
if (!ALG3D_checkResult(ALG3D_writeXYZ(fnameOut, vecDyn, ptsDyn))) return 0;
Parameters
[in]ptsStatpointer to the begin of the static data
[in]normalsStatpointer to the begin of the normal vectors of ptsStat
[in]nStatnumber of points in static data
[in]ptsDynpointer to the begin of the dynamic data
[in]nDynnumber of points in dynamic data
[in]prmset of parameters to configure the algorithm
[out]resresult statistics
[out]R3x3 rotation matrix (row-major)
[out]T3x1 translation matrix (vector)
[out]S1x1 scaling value
Returns
ERR_NONE on success (at least one termination criteria has been reached), ERR_CONDITION if the computation stopped early (if not enough correlating pairs are found)
Note
Also check the result parameters to find out detailed status information

◆ ALG3D_computeRegistrationPlaneTree()

int ALG3D_computeRegistrationPlaneTree ( const ALG3D_KDTREE  tree,
const ALG3D_Point3d normalsStat,
const ALG3D_Point3d ptsDyn,
size_t  nDyn,
icpParam  prm,
icpResult res,
double  R[9],
double  T[3],
double  S[1] 
)

#include <algorithm.h>

Computes the registration/alignment of two 3d point clouds.

A given 'static' point cloud, i.e. a kdTree, is used as a reference. For the second, 'dynamic' point cloud the translation, rotation and scaling are computed iteratively. The iterations stop when a (local) distance minimum is found, i.e. until there is no change anymore or one of the termination criteria has been reached (min. tolerance, max. iterations).

This algorithm uses the "point-to-plane" metric to compute correspondences between the static and dynamic point cloud. This metric requires the static point cloud to have normal vectors. These must be provided by the user. You can use one of the ALG3D_computeNormalVectors functions to perform this operation or use your own one.

Note
The elements of normalsStat must correspond to the points from which tree was built (e.g. normalsStat[0] = normal of first point of point cloud from which tree was built etc.).
code example
//read static point data
if(!ALG3D_checkResult(ALG3D_readXYZ(fnameStat,&dataStat,&ptsStat))) return 0;
//read dynamic point data (will be transformed)
if(!ALG3D_checkResult(ALG3D_readXYZ(fnameDyn ,&dataDyn ,&ptsDyn))) return 0;
//create kdTree from static data
if(!ALG3D_checkResult(ALG3D_createKdTree(dataStat, ptsStat, &tree) )) return 0;
normalsStat = (ALG3D_Point3d*)malloc(ptsStat * sizeof(ALG3D_Point3d));
//compute normal vectors using the KdTree
if(!ALG3D_checkResult(ALG3D_computeNormalVectorsKNNTree(tree, normalsStat, 30) )) return 0;
//compute registration, i.e. transformation matrix
if (!ALG3D_checkResult(ALG3D_computeRegistrationPlaneTree(tree, normalsStat, dataDyn, ptsDyn, prm, &res, R, T, S))) return 0;
Parameters
[in]treekdTree object containing the static point set
[in]normalsStatpointer to the begin of the normal vectors of tree
[in]ptsDynpointer to the begin of the dynamic data
[in]nDynnumber of points in dynamic data
[in]prmset of parameters to configure the algorithm
[out]resresult statistics
[out]R3x3 rotation matrix (row-major)
[out]T3x1 translation matrix (vector)
[out]S1x1 scaling value
Returns
ERR_NONE on success (at least one termination criteria has been reached), ERR_CONDITION if the computation stopped early (if not enough correlating pairs are found)
Note
Also check the result parameters to find out detailed status information
Examples
testKdTree.c.

◆ ALG3D_computeRegistrationTree()

int ALG3D_computeRegistrationTree ( const ALG3D_KDTREE  tree,
const ALG3D_Point3d ptsDyn,
size_t  nDyn,
icpParam  prm,
icpResult res,
double  R[9],
double  T[3],
double  S[1] 
)

#include <algorithm.h>

Computes the registration/alignment of two 3d point clouds.

A given 'static' point cloud, i.e. a kdTree, is used as a reference. For the second, 'dynamic' point cloud the translation, rotation and scaling are computed iteratively . The iterations stop when a (local) distance minimum is found, i.e. until there is no change anymore or one of the termination criteria has been reached (min. tolerance, max. iterations).

code example
//read static point data
if(!ALG3D_checkResult(ALG3D_readXYZ(fnameStat,&dataStat,&ptsStat))) return 0;
//read dynamic point data (will be transformed)
if(!ALG3D_checkResult(ALG3D_readXYZ(fnameDyn ,&dataDyn ,&ptsDyn))) return 0;
//create kdTree from static data
if(!ALG3D_checkResult(ALG3D_createKdTree(dataStat, ptsStat, &tree) )) return 0;
//compute registration, i.e. transformation matrix
if (!ALG3D_checkResult(ALG3D_computeRegistrationTree(tree, dataDyn, ptsDyn, prm, &res, R, T, S))) return 0;
Parameters
[in]treekdTree object containing the static point set
[in]ptsDynpointer to the begin of the dynamic data
[in]nDynnumber of points in dynamic data
[in]prmset of parameters to configure the algorithm
[out]resresult statistics
[out]R3x3 rotation matrix (row-major)
[out]T3x1 translation matrix (vector)
[out]S1x1 scaling value
Returns
ERR_NONE on success (at least one termination criteria has been reached), ERR_CONDITION if the computation stopped early (if not enough correlating pairs are found)
Note
Also check the result parameters to find out detailed status information
Examples
testKdTree.c.
ALG3D_computeNormalVectorsKNNTree
int __cdecl ALG3D_computeNormalVectorsKNNTree(const ALG3D_KDTREE tree, ALG3D_Point3d *normals, size_t maxPoints)
Computes surface normals vector for the given pointcloud (represented by its KdTree) from K-Nearest-N...
Definition: algorithm.cpp:701
icpMultiViewParam
MultiView-ICP input parameters.
Definition: algorithm.h:126
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
icpResult::meanDev
double meanDev
mean distance between the data sets
Definition: algorithm.h:141
icpParam::samples
size_t samples
number of sample points that should be used for the alignment (0=all)
Definition: algorithm.h:121
ALG3D_createKdTree
int __cdecl ALG3D_createKdTree(const ALG3D_Point3d *data, size_t numPts, ALG3D_KDTREE *tree)
Builds a kdTree from given 3D data.
Definition: algorithm.cpp:1475
ALG3D_Point3d::vec
double vec[3]
point data, e.g. 0=x,1=y,2=z
Definition: interface.h:41
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
icpResult::pairs
size_t pairs
number of correlating data pairs
Definition: algorithm.h:143
ALG3D_Point3d
Data type for 3d points/vectors with a user-defined pointer.
Definition: interface.h:32
ALG3D_computeNormalVectorsKNN
int __cdecl ALG3D_computeNormalVectorsKNN(const ALG3D_Point3d *pts, size_t numPts, ALG3D_Point3d *normals, size_t maxPoints)
Computes surface normals vector for the given pointcloud from K-Nearest-Neighbours.
Definition: algorithm.cpp:603
icpMultiViewParam::fixedIndex
int fixedIndex
index of the dataset that should stay fixed (-1 = none)
Definition: algorithm.h:133
icpResult::stdDev
double stdDev
standard deviation of the point distances
Definition: algorithm.h:142
ALG3D_computeRegistrationMesh
int __cdecl ALG3D_computeRegistrationMesh(const ALG3D_MESHTREE meshTree, ALG3D_Point3d *ptsDyn, size_t nDyn, icpParam prm, icpResult *res, double R[9], double T[3], double S[1])
Computes the registration/alignment of a3d point clouds to a mesh (e.g. CAD model).
Definition: algorithm.cpp:3044
icpStandardDoF
@ icpStandardDoF
Enables all rotation and translation degrees of freedom (without scaling)
Definition: algorithm.h:106
icpParam
ICP input parameters.
Definition: algorithm.h:116
ALG3D_computeMultiViewRegistration
int __cdecl ALG3D_computeMultiViewRegistration(const ALG3D_Point3d **pointclouds, size_t *pointcloudsSizes, size_t numPointclouds, icpMultiViewParam prm, icpResult *res, icpMultiViewTransformations *transformations)
Computes the registration/alignment of multiple 3d point clouds.
Definition: algorithm.cpp:2616
ALG3D_computeRegistrationPlaneTree
int __cdecl ALG3D_computeRegistrationPlaneTree(const ALG3D_KDTREE tree, const ALG3D_Point3d *normalsStat, const ALG3D_Point3d *ptsDyn, size_t nDyn, icpParam prm, icpResult *res, double R[9], double T[3], double S[1])
Computes the registration/alignment of two 3d point clouds.
Definition: algorithm.cpp:2491
icpResult
ICP result parameters.
Definition: algorithm.h:137
ALG3D_transformPoints
int __cdecl ALG3D_transformPoints(ALG3D_Point3d *pts, size_t numPts, const double R[9], const double T[3], const double S[1])
Transforms points by the given tranformation.
Definition: geometry.cpp:102
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
icpMultiViewParam::radius
double radius
maximum search radius (0=auto)
Definition: algorithm.h:128
ALG3D_computeRegistrationCorrespondences
int __cdecl ALG3D_computeRegistrationCorrespondences(const ALG3D_Point3d *ptsStat, const ALG3D_Point3d *ptsDyn, size_t numPts, int options, double R[9], double T[3], double S[1])
Computes the registration/alignment of two corresponding sets of 3D points.
Definition: algorithm.cpp:2311
icpResult::minDist
double minDist
minimal point distance
Definition: algorithm.h:139
icpMultiViewParam::iterations
size_t iterations
maximum number of iterations
Definition: algorithm.h:129
ALG3D_writeXYZ
int __cdecl ALG3D_writeXYZ(const char *fname, const ALG3D_Point3d *data, size_t numPts)
Writes 3D point data to a XYZ-ASCII-file.
Definition: io.cpp:123
ALG3D_computeRegistration
int __cdecl ALG3D_computeRegistration(const ALG3D_Point3d *ptsStat, size_t nStat, const ALG3D_Point3d *ptsDyn, size_t nDyn, icpParam prm, icpResult *res, double R[9], double T[3], double S[1])
Computes the registration/alignment of two 3d point clouds.
Definition: algorithm.cpp:109
ALG3D_computeRegistrationPlane
int __cdecl ALG3D_computeRegistrationPlane(const ALG3D_Point3d *ptsStat, const ALG3D_Point3d *normalsStat, size_t nStat, const ALG3D_Point3d *ptsDyn, size_t nDyn, icpParam prm, icpResult *res, double R[9], double T[3], double S[1])
Computes the registration/alignment of two 3d point clouds.
Definition: algorithm.cpp:2429
icpResult::maxDist
double maxDist
maximal point distance
Definition: algorithm.h:140
icpMultiViewTransformations
MultiView-ICP transformations.
Definition: algorithm.h:148
icpMatchCoG
@ icpMatchCoG
enable intial alignment of the center of gravities
Definition: algorithm.h:110
icpMultiViewParam::tolerance
double tolerance
required minimum distance (iterations stop when this value has been reached)
Definition: algorithm.h:130
icpMultiViewParam::thinningRadius
double thinningRadius
internal thinning radius
Definition: algorithm.h:132
ALG3D_computeRegistrationTree
int __cdecl ALG3D_computeRegistrationTree(const ALG3D_KDTREE tree, const ALG3D_Point3d *ptsDyn, size_t nDyn, icpParam prm, icpResult *res, double R[9], double T[3], double S[1])
Computes the registration/alignment of two 3d point clouds.
Definition: algorithm.cpp:2130
icpParam::options
int options
sets options (see icpOptions). If set to 0 icpStandardDoF is used.
Definition: algorithm.h:122