CCCoreLib 31 May 2022
CloudCompare Core algorithms
Classes | Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | List of all members
CCCoreLib::DgmOctree Class Reference

The octree structure used throughout the library. More...

#include <DgmOctree.h>

Inheritance diagram for CCCoreLib::DgmOctree:
Inheritance graph
[legend]
Collaboration diagram for CCCoreLib::DgmOctree:
Collaboration graph
[legend]

Classes

struct  BoxNeighbourhood
 Input/output parameters structure for getPointsInBoxNeighbourhood. More...
 
struct  CellDescriptor
 Structure used during nearest neighbour search. More...
 
struct  CylindricalNeighbourhood
 Input/output parameters structure for getPointsInCylindricalNeighbourhood. More...
 
struct  IndexAndCode
 Association between an index and the code of an octree cell. More...
 
struct  NearestNeighboursSearchStruct
 Container of in/out parameters for nearest neighbour(s) search. More...
 
struct  octreeCell
 Octree cell descriptor. More...
 
struct  octreeTopDownScanStruct
 Internal structure used to perform a top-down scan of the octree. More...
 
struct  PointDescriptor
 Structure used during nearest neighbour search. More...
 
struct  ProgressiveCylindricalNeighbourhood
 Input/output parameters structure for getPointsInCylindricalNeighbourhoodProgressive. More...
 

Public Types

enum  RayCastProcess { RC_NEAREST_POINT , RC_CLOSE_POINTS }
 Ray casting processes.
 
using CellCode = unsigned
 Type of the code of an octree cell. More...
 
using cellCodesContainer = std::vector< CellCode >
 Octree cell codes container.
 
using cellIndexesContainer = std::vector< unsigned int >
 Octree cell indexes container.
 
using NeighboursSet = std::vector< PointDescriptor >
 A set of neighbours.
 
using NeighbourCellsSet = std::vector< CellDescriptor >
 A set of neighbour cells descriptors.
 
using cellsContainer = std::vector< IndexAndCode >
 Container of 'IndexAndCode' structures.
 
using octreeCellFunc = bool(*)(const octreeCell &, void **, NormalizedProgress *)
 Generic form of a function that can be applied automatically to all cells of the octree. More...
 

Public Member Functions

 DgmOctree (GenericIndexedCloudPersist *cloud)
 DgmOctree constructor. More...
 
 ~DgmOctree () override=default
 DgmOctree destructor.
 
virtual void clear ()
 Clears the octree.
 
int build (GenericProgressCallback *progressCb=nullptr)
 Builds the structure. More...
 
int build (const CCVector3 &octreeMin, const CCVector3 &octreeMax, const CCVector3 *pointsMinFilter=nullptr, const CCVector3 *pointsMaxFilter=nullptr, GenericProgressCallback *progressCb=nullptr)
 Builds the structure with constraints. More...
 
unsigned getNumberOfProjectedPoints () const
 Returns the number of points projected into the octree. More...
 
const CCVector3getOctreeMins () const
 Returns the lower boundaries of the octree. More...
 
const CCVector3getOctreeMaxs () const
 Returns the higher boundaries of the octree. More...
 
void getBoundingBox (CCVector3 &bbMin, CCVector3 &bbMax) const
 Returns the octree bounding box. More...
 
const int * getMinFillIndexes (unsigned char level) const
 Returns the lowest cell positions in the octree along all dimensions and for a given level of subdivision. More...
 
const int * getMaxFillIndexes (unsigned char level) const
 Returns the highest cell positions in the octree along all dimensions and for a given level of subdivision. More...
 
const PointCoordinateType & getCellSize (unsigned char level) const
 Returns the octree cells length for a given level of subdivision. More...
 
void getCellDistanceFromBorders (const Tuple3i &cellPos, unsigned char level, int *cellDists) const
 Returns distance form a cell to the filled octree borders in all directions. More...
 
void getCellDistanceFromBorders (const Tuple3i &cellPos, unsigned char level, int neighbourhoodLength, int *cellDists) const
 Returns distance from cell center to cell neighbourhood INSIDE filled octree. More...
 
bool getPointsInCellByCellIndex (ReferenceCloud *cloud, unsigned cellIndex, unsigned char level, bool clearOutputCloud=true) const
 Returns the points lying in a specific cell. More...
 
bool getPointsInCell (CellCode cellCode, unsigned char level, ReferenceCloud *subset, bool isCodeTruncated=false, bool clearOutputCloud=true) const
 Returns the points lying in a specific cell. More...
 
ReferenceCloudgetPointsInCellsWithSortedCellCodes (cellCodesContainer &cellCodes, unsigned char level, ReferenceCloud *subset, bool areCodesTruncated=false) const
 Returns the points lying in multiple cells. More...
 
unsigned findPointNeighbourhood (const CCVector3 *_queryPoint, ReferenceCloud *Yk, unsigned maxNumberOfNeighbors, unsigned char level, double &maxSquareDist, double maxSearchDist=0, int *finalNeighbourhoodSize=nullptr) const
 Finds the nearest neighbours around a query point. More...
 
double findTheNearestNeighborStartingFromCell (NearestNeighboursSearchStruct &nNSS) const
 Advanced form of the nearest neighbour search algorithm (unique neighbour) More...
 
unsigned findNearestNeighborsStartingFromCell (NearestNeighboursSearchStruct &nNSS, bool getOnlyPointsWithValidScalar=false) const
 Advanced form of the nearest neighbours search algorithm (multiple neighbours) More...
 
int findNeighborsInASphereStartingFromCell (NearestNeighboursSearchStruct &nNSS, double radius, bool sortValues=true) const
 Advanced form of the nearest neighbours search algorithm (in a sphere) More...
 
int getPointsInSphericalNeighbourhood (const CCVector3 &sphereCenter, PointCoordinateType radius, NeighboursSet &neighbours, unsigned char level) const
 Returns the points falling inside a sphere. More...
 
std::size_t getPointsInCylindricalNeighbourhood (CylindricalNeighbourhood &params) const
 Returns the points falling inside a cylinder. More...
 
std::size_t getPointsInCylindricalNeighbourhoodProgressive (ProgressiveCylindricalNeighbourhood &params) const
 Same as getPointsInCylindricalNeighbourhood with progressive approach. More...
 
std::size_t getPointsInBoxNeighbourhood (BoxNeighbourhood &params) const
 Returns the points falling inside a box. More...
 
void getTheCellPosWhichIncludesThePoint (const CCVector3 *thePoint, Tuple3i &cellPos) const
 Returns the position FOR THE DEEPEST LEVEL OF SUBDIVISION of the cell that includes a given point. More...
 
void getTheCellPosWhichIncludesThePoint (const CCVector3 *thePoint, Tuple3i &cellPos, unsigned char level) const
 Returns the position for a given level of subdivision of the cell that includes a given point. More...
 
void getTheCellPosWhichIncludesThePoint (const CCVector3 *thePoint, Tuple3i &cellPos, unsigned char level, bool &inBounds) const
 Returns the position for a given level of subdivision of the cell that includes a given point. More...
 
void getCellPos (CellCode code, unsigned char level, Tuple3i &cellPos, bool isCodeTruncated) const
 Returns the cell position for a given level of subdivision of a cell designated by its code. More...
 
void computeCellCenter (CellCode code, unsigned char level, CCVector3 &center, bool isCodeTruncated=false) const
 Returns the cell center for a given level of subdivision of a cell designated by its code. More...
 
void computeCellCenter (const Tuple3i &cellPos, unsigned char level, CCVector3 &center) const
 Returns the cell center for a given level of subdivision of a cell designated by its position. More...
 
void computeCellCenter (const Tuple3s &cellPos, unsigned char level, CCVector3 &center) const
 Short version of computeCellCenter.
 
void computeCellLimits (CellCode code, unsigned char level, CCVector3 &cellMin, CCVector3 &cellMax, bool isCodeTruncated=false) const
 Returns the spatial limits of a given cell. More...
 
unsigned getCellIndex (CellCode truncatedCellCode, unsigned char bitShift) const
 Returns the index of a given cell represented by its code. More...
 
unsigned char findBestLevelForAGivenNeighbourhoodSizeExtraction (PointCoordinateType radius) const
 Determines the best level of subdivision of the octree at which to apply the nearest neighbours search algorithm (inside a sphere) depending on the sphere radius. More...
 
unsigned char findBestLevelForComparisonWithOctree (const DgmOctree *theOtherOctree) const
 Determines the best level of subdivision of the octree at which to apply a cloud-2-cloud distance computation algorithm. More...
 
unsigned char findBestLevelForAGivenPopulationPerCell (unsigned indicativeNumberOfPointsPerCell) const
 Determines the best subdivision level of the octree that gives the average population per cell closest to the input value. More...
 
unsigned char findBestLevelForAGivenCellNumber (unsigned indicativeNumberOfCells) const
 Determines the best subdivision level of the octree to match a given number of cells. More...
 
const CellCodegetCellCode (unsigned index) const
 Returns the ith cell code.
 
bool getCellCodes (unsigned char level, cellCodesContainer &vec, bool truncatedCodes=false) const
 Returns the list of codes corresponding to the octree cells for a given level of subdivision. More...
 
bool getCellIndexes (unsigned char level, cellIndexesContainer &vec) const
 Returns the list of indexes corresponding to the octree cells for a given level of subdivision. More...
 
bool getCellCodesAndIndexes (unsigned char level, cellsContainer &vec, bool truncatedCodes=false) const
 Returns the list of indexes and codes corresponding to the octree cells for a given level of subdivision. More...
 
void diff (const cellCodesContainer &codesA, const cellCodesContainer &codesB, cellCodesContainer &diffA, cellCodesContainer &diffB) const
 Returns the cells that differ between two octrees (for a same implicit level of subdivision) More...
 
bool diff (unsigned char octreeLevel, const cellsContainer &codesA, const cellsContainer &codesB, int &diffA, int &diffB, int &cellsA, int &cellsB) const
 Returns the differences (in terms of number of cells) between two octrees for a given level of subdivision. More...
 
const unsigned & getCellNumber (unsigned char level) const
 Returns the number of cells for a given level of subdivision.
 
double computeMeanOctreeDensity (unsigned char level) const
 Computes mean octree density (point/cell) at a given level of subdivision. More...
 
int extractCCs (const cellCodesContainer &cellCodes, unsigned char level, bool sixConnexity, GenericProgressCallback *progressCb=nullptr) const
 Computes the connected components (considering the octree cells only) for a given level of subdivision (partial) More...
 
int extractCCs (unsigned char level, bool sixConnexity, GenericProgressCallback *progressCb=nullptr) const
 Computes the connected components (considering the octree cells only) for a given level of subdivision (complete) More...
 
unsigned executeFunctionForAllCellsStartingAtLevel (unsigned char startingLevel, octreeCellFunc func, void **additionalParameters, unsigned minNumberOfPointsPerCell, unsigned maxNumberOfPointsPerCell, bool multiThread=true, GenericProgressCallback *progressCb=nullptr, const char *functionTitle=nullptr, int maxThreadCount=0)
 Method to apply automatically a specific function to each cell of the octree. More...
 
unsigned executeFunctionForAllCellsAtLevel (unsigned char level, octreeCellFunc func, void **additionalParameters, bool multiThread=false, GenericProgressCallback *progressCb=nullptr, const char *functionTitle=nullptr, int maxThreadCount=0)
 Method to apply automatically a specific function to each cell of the octree. More...
 
bool rayCast (const CCVector3 &rayAxis, const CCVector3 &rayOrigin, double maxRadiusOrFov, bool isFOV, RayCastProcess process, std::vector< PointDescriptor > &output) const
 Ray casting algorithm.
 
GenericIndexedCloudPersistassociatedCloud () const
 Returns the associated cloud.
 
const cellsContainerpointsAndTheirCellCodes () const
 Returns the octree 'structure'.
 
- Public Member Functions inherited from CCCoreLib::GenericOctree
virtual ~GenericOctree ()=default
 Default destructor.
 

Static Public Member Functions

static unsigned char GET_BIT_SHIFT (unsigned char level)
 Returns the binary shift for a given level of subdivision. More...
 
static int OCTREE_LENGTH (int level)
 Returns the octree length (in term of cells) for a given level of subdivision. More...
 
static CellCode GenerateTruncatedCellCode (const Tuple3i &cellPos, unsigned char level)
 Generates the truncated cell code of a cell given its position at a given level of subdivision. More...
 
static CellCode GenerateTruncatedCellCode (const Tuple3s &pos, unsigned char level)
 Short version of GenerateTruncatedCellCode.
 
static PointCoordinateType ComputeMinDistanceToCellBorder (const CCVector3 &queryPoint, PointCoordinateType cs, const CCVector3 &cellCenter)
 Computes the minimal distance between a point and the borders (faces) of the cell (cube) in which it is included. More...
 
static bool MultiThreadSupport ()
 Returns whether multi-threading (parallel) computation is supported or not.
 

Static Public Attributes

static const int MAX_OCTREE_LEVEL = 10
 Max octree subdivision level. More...
 
static const int MAX_OCTREE_LENGTH = (1 << MAX_OCTREE_LEVEL)
 Max octree length at last level of subdivision (number of cells) More...
 
static const CellCode INVALID_CELL_CODE = (~static_cast<CellCode>(0))
 Invalid cell code. More...
 

Protected Member Functions

int genericBuild (GenericProgressCallback *progressCb=nullptr)
 Generic method to build the octree structure. More...
 
void updateCellSizeTable ()
 Updates the tables containing the octree cells length for each level of subdivision.
 
void updateCellCountTable ()
 Updates the tables containing the number of octree cells for each level of subdivision.
 
void computeCellsStatistics (unsigned char level)
 Computes statistics about cells for a given level of subdivision. More...
 
void getNeighborCellsAround (const Tuple3i &cellPos, cellIndexesContainer &neighborCellsIndexes, int neighbourhoodLength, unsigned char level) const
 Returns the indexes of the neighbourhing (existing) cells of a given cell. More...
 
void getPointsInNeighbourCellsAround (NearestNeighboursSearchStruct &nNSS, int neighbourhoodLength, bool getOnlyPointsWithValidScalar=false) const
 Gets point in the neighbourhing cells of a specific cell. More...
 
unsigned getCellIndex (CellCode truncatedCellCode, unsigned char bitShift, unsigned begin, unsigned end) const
 Returns the index of a given cell represented by its code. More...
 

Protected Attributes

cellsContainer m_thePointsAndTheirCellCodes
 The coded octree structure. More...
 
GenericIndexedCloudPersistm_theAssociatedCloud
 Associated cloud.
 
unsigned m_numberOfProjectedPoints
 Number of points projected in the octree.
 
unsigned m_nearestPow2
 Nearest power of 2 smaller than the number of points (used for binary search)
 
CCVector3 m_dimMin
 Min coordinates of the octree bounding-box.
 
CCVector3 m_dimMax
 Max coordinates of the octree bounding-box.
 
CCVector3 m_pointsMin
 Min coordinates of the bounding-box of the set of points projected in the octree.
 
CCVector3 m_pointsMax
 Max coordinates of the bounding-box of the set of points projected in the octree.
 
PointCoordinateType m_cellSize [MAX_OCTREE_LEVEL+2]
 Cell dimensions for all subdivision levels.
 
int m_fillIndexes [(MAX_OCTREE_LEVEL+1) *6]
 Min and max occupied cells indexes, for all dimensions and every subdivision level.
 
unsigned m_cellCount [MAX_OCTREE_LEVEL+1]
 Number of cells per level of subdivision.
 
unsigned m_maxCellPopulation [MAX_OCTREE_LEVEL+1]
 Max cell population per level of subdivision.
 
double m_averageCellPopulation [MAX_OCTREE_LEVEL+1]
 Average cell population per level of subdivision.
 
double m_stdDevCellPopulation [MAX_OCTREE_LEVEL+1]
 Std. dev. of cell population per level of subdivision.
 

Detailed Description

The octree structure used throughout the library.

Implements the GenericOctree interface. Corresponds to the octree structure developed during Daniel Girardeau-Montaut's PhD (see PhD manuscript, Chapter 4).

Member Typedef Documentation

◆ CellCode

Type of the code of an octree cell.

Warning
3 bits per level are required.
Never pass a 'constant initializer' by reference

◆ octreeCellFunc

Generic form of a function that can be applied automatically to all cells of the octree.

See DgmOctree::executeFunctionForAllCellsAtLevel and DgmOctree::executeFunctionForAllCellsStartingAtLevel. The parameters of such a function are:

  • (octreeCell) cell descriptor
  • (void**) table of user parameters for the function (maybe void)
  • (NormalizedProgress*) optional (normalized) progress callback
  • return success

Constructor & Destructor Documentation

◆ DgmOctree()

DgmOctree::DgmOctree ( GenericIndexedCloudPersist cloud)
explicit

DgmOctree constructor.

METHODS

Parameters
cloudthe cloud to construct the octree on

Member Function Documentation

◆ build() [1/2]

int DgmOctree::build ( const CCVector3 octreeMin,
const CCVector3 octreeMax,
const CCVector3 pointsMinFilter = nullptr,
const CCVector3 pointsMaxFilter = nullptr,
GenericProgressCallback progressCb = nullptr 
)

Builds the structure with constraints.

Octree spatial limits must be specified. Also, if specified, points falling outside the "pointsFilter" limits won't be projected into the octree structure. Otherwise, all points will be taken into account. Octree 3D limits in space should be cubical.

Parameters
octreeMinthe lower limits for the octree cells along X, Y and Z
octreeMaxthe upper limits for the octree cells along X, Y and Z
pointsMinFilterthe lower limits for the projected points along X, Y and Z (if specified)
pointsMaxFilterthe upper limits for the projected points along X, Y and Z (if specified)
progressCbthe client application can get some notification of the process progress through this callback mechanism (see GenericProgressCallback)
Returns
the number of points projected in the octree

◆ build() [2/2]

int DgmOctree::build ( GenericProgressCallback progressCb = nullptr)

Builds the structure.

Octree 3D limits are determined automatically.

Parameters
progressCbthe client application can get some notification of the process progress through this callback mechanism (see GenericProgressCallback)
Returns
the number of points projected in the octree

◆ computeCellCenter() [1/2]

void CCCoreLib::DgmOctree::computeCellCenter ( CellCode  code,
unsigned char  level,
CCVector3 center,
bool  isCodeTruncated = false 
) const
inline

Returns the cell center for a given level of subdivision of a cell designated by its code.

Parameters
codethe cell code
levelthe level of subdivision
centerthe computed center
isCodeTruncatedindicates if the given code is truncated or not

◆ computeCellCenter() [2/2]

void CCCoreLib::DgmOctree::computeCellCenter ( const Tuple3i cellPos,
unsigned char  level,
CCVector3 center 
) const
inline

Returns the cell center for a given level of subdivision of a cell designated by its position.

Parameters
cellPosthe cell position
levelthe level of subdivision
centerthe computed center

◆ computeCellLimits()

void DgmOctree::computeCellLimits ( CellCode  code,
unsigned char  level,
CCVector3 cellMin,
CCVector3 cellMax,
bool  isCodeTruncated = false 
) const

Returns the spatial limits of a given cell.

Parameters
codethe cell code
levelthe level of subdivision
cellMinthe minimum coordinates along each dimension
cellMaxthe maximum coordinates along each dimension
isCodeTruncatedindicates if the given code is truncated or not

◆ computeCellsStatistics()

void DgmOctree::computeCellsStatistics ( unsigned char  level)
protected

Computes statistics about cells for a given level of subdivision.

This method requires some computation, therefore it shouldn't be called too often.

Parameters
levelthe level of subdivision

◆ computeMeanOctreeDensity()

double DgmOctree::computeMeanOctreeDensity ( unsigned char  level) const

Computes mean octree density (point/cell) at a given level of subdivision.

Parameters
levelthe level of subdivision
Returns
mean density (point/cell)

◆ ComputeMinDistanceToCellBorder()

static PointCoordinateType CCCoreLib::DgmOctree::ComputeMinDistanceToCellBorder ( const CCVector3 queryPoint,
PointCoordinateType  cs,
const CCVector3 cellCenter 
)
inlinestatic

Computes the minimal distance between a point and the borders (faces) of the cell (cube) in which it is included.

Parameters
queryPointthe point
csthe cell size (as cells are cubical, it's the same along every dimension)
cellCenterthe cell center
Returns
the minimal distance

◆ diff() [1/2]

void DgmOctree::diff ( const cellCodesContainer codesA,
const cellCodesContainer codesB,
cellCodesContainer diffA,
cellCodesContainer diffB 
) const

Returns the cells that differ between two octrees (for a same implicit level of subdivision)

Warning: the two octrees should have been computed with the same bounding-box.

Parameters
codesAthe cell codes of the first octree for the implicit level
codesBthe cell codes of the second octree for the implicit level
diffAthe cells of the first octree that are not in the second octree
diffBthe cells of the second octree that are not in the first octree

◆ diff() [2/2]

bool DgmOctree::diff ( unsigned char  octreeLevel,
const cellsContainer codesA,
const cellsContainer codesB,
int &  diffA,
int &  diffB,
int &  cellsA,
int &  cellsB 
) const

Returns the differences (in terms of number of cells) between two octrees for a given level of subdivision.

Warning: the two octrees should have been computed with the same bounding-box.

Parameters
octreeLevelthe octree level
codesAthe cell codes (and point index) of the first octree
codesBthe cell codes (and point index) of the second octree
diffAthe number of cells of the first octree that are not in the second octree
diffBthe number of cells of the second octree that are not in the first octree
cellsAthe number of cells of the first octree for the given number of subdivision
cellsBthe number of cells of the second octree for the given number of subdivision
Returns
false if it could not calculate the differences

◆ executeFunctionForAllCellsAtLevel()

unsigned DgmOctree::executeFunctionForAllCellsAtLevel ( unsigned char  level,
octreeCellFunc  func,
void **  additionalParameters,
bool  multiThread = false,
GenericProgressCallback progressCb = nullptr,
const char *  functionTitle = nullptr,
int  maxThreadCount = 0 
)

Method to apply automatically a specific function to each cell of the octree.

The function to apply should be of the form DgmOctree::octreeCellFunc. In this case the octree cells are scanned one by one at the same level of subdivision.

Parallel processing is based on tbb::parallel_for or QtConcurrent.

Parameters
levelthe level of subdivision
functhe function to apply
additionalParametersthe function parameters
multiThreadwhether to use parallel processing or not
progressCbthe client application can get some notification of the process progress through this callback mechanism (see GenericProgressCallback)
functionTitlefunction title
maxThreadCountthe maximum number of threads to use (0 = all). Ignored if 'multiThread' is false or with tbb.
Returns
the number of processed cells (or 0 is something went wrong)

◆ executeFunctionForAllCellsStartingAtLevel()

unsigned DgmOctree::executeFunctionForAllCellsStartingAtLevel ( unsigned char  startingLevel,
octreeCellFunc  func,
void **  additionalParameters,
unsigned  minNumberOfPointsPerCell,
unsigned  maxNumberOfPointsPerCell,
bool  multiThread = true,
GenericProgressCallback progressCb = nullptr,
const char *  functionTitle = nullptr,
int  maxThreadCount = 0 
)

Method to apply automatically a specific function to each cell of the octree.

The function to apply should be of the form DgmOctree::octreeCellFunc. In this case the octree cells are scanned one by one at the same level of subdivision, but the scan can also be sometimes done (inside the cell) at deeper levels, in order to limit the number of points in a cell. Thanks to this, the function is applied on a limited number of points, avoiding great loss of performances. The only limitation is when the level of subdivision is deepest level. In this case no more splitting is possible.

Parallel processing is based on tbb::parallel_for or QtCocncurrent.

Parameters
startingLevelthe initial level of subdivision
functhe function to apply
additionalParametersthe function parameters
minNumberOfPointsPerCellminimal number of points inside a cell (indicative)
maxNumberOfPointsPerCellmaximum number of points inside a cell (indicative)
multiThreadwhether to use parallel processing or not
progressCbthe client application can get some notification of the process progress through this callback mechanism (see GenericProgressCallback)
maxThreadCountthe maximum number of threads to use (0 = all). Ignored if 'multiThread' is false or with tbb.
functionTitlefunction title
Returns
the number of processed cells (or 0 is something went wrong)

◆ extractCCs() [1/2]

int DgmOctree::extractCCs ( const cellCodesContainer cellCodes,
unsigned char  level,
bool  sixConnexity,
GenericProgressCallback progressCb = nullptr 
) const

Computes the connected components (considering the octree cells only) for a given level of subdivision (partial)

The octree is seen as a regular 3D grid, and each cell of this grid is either set to 0 (if no points lies in it) or to 1 (if some points lie in it, e.g. if it is indeed a cell of this octree). This version of the algorithm can be applied by considering only a specified list of octree cells (ignoring the others).

Parameters
cellCodesthe cell codes to consider for the CC computation
levelthe level of subdivision at which to perform the algorithm
sixConnexityindicates if the CC's 3D connexity should be 6 (26 otherwise)
progressCbthe client application can get some notification of the process progress through this callback mechanism (see GenericProgressCallback)
Returns
error code:
  • '>= 0' = number of components
  • '-1' = no cells (input)
  • '-2' = not enough memory
  • '-3' = no CC found

◆ extractCCs() [2/2]

int DgmOctree::extractCCs ( unsigned char  level,
bool  sixConnexity,
GenericProgressCallback progressCb = nullptr 
) const

Computes the connected components (considering the octree cells only) for a given level of subdivision (complete)

The octree is seen as a regular 3D grid, and each cell of this grid is either set to 0 (if no points lies in it) or to 1 (if some points lie in it, e.g. if it is indeed a cell of this octree). This version of the algorithm is directly applied on the whole octree.

Parameters
levelthe level of subdivision at which to perform the algorithm
sixConnexityindicates if the CC's 3D connexity should be 6 (26 otherwise)
progressCbthe client application can get some notification of the process progress through this callback mechanism (see GenericProgressCallback)
Returns
error code:
  • '>= 0' = number of components
  • '-1' = no cells (input)
  • '-2' = not enough memory
  • '-3' = no CC found

◆ findBestLevelForAGivenCellNumber()

unsigned char DgmOctree::findBestLevelForAGivenCellNumber ( unsigned  indicativeNumberOfCells) const

Determines the best subdivision level of the octree to match a given number of cells.

Parameters
indicativeNumberOfCells'desired' number of cells
Returns
the 'best' level

◆ findBestLevelForAGivenNeighbourhoodSizeExtraction()

unsigned char DgmOctree::findBestLevelForAGivenNeighbourhoodSizeExtraction ( PointCoordinateType  radius) const

Determines the best level of subdivision of the octree at which to apply the nearest neighbours search algorithm (inside a sphere) depending on the sphere radius.

Parameters
radiusthe sphere radius
Returns
the 'best' level

◆ findBestLevelForAGivenPopulationPerCell()

unsigned char DgmOctree::findBestLevelForAGivenPopulationPerCell ( unsigned  indicativeNumberOfPointsPerCell) const

Determines the best subdivision level of the octree that gives the average population per cell closest to the input value.

Parameters
indicativeNumberOfPointsPerCell'desired' average number of points per cell
Returns
the 'best' level

◆ findBestLevelForComparisonWithOctree()

unsigned char DgmOctree::findBestLevelForComparisonWithOctree ( const DgmOctree theOtherOctree) const

Determines the best level of subdivision of the octree at which to apply a cloud-2-cloud distance computation algorithm.

The octree instance on which is "applied" this method should be the compared cloud's one. "theOtherOctree" should be the reference cloud's octree.

Parameters
theOtherOctreethe octree of the other cloud
Returns
the 'best' level

◆ findNearestNeighborsStartingFromCell()

unsigned DgmOctree::findNearestNeighborsStartingFromCell ( NearestNeighboursSearchStruct nNSS,
bool  getOnlyPointsWithValidScalar = false 
) const

Advanced form of the nearest neighbours search algorithm (multiple neighbours)

This version is optimized for a multiple nearest neighbours search that is applied around several query points included in the same octree cell. See DgmOctree::NearestNeighboursSearchStruct for more details.

Parameters
nNSSNN search parameters
getOnlyPointsWithValidScalarwhether to ignore points having an invalid associated scalar value
Returns
the number of neighbours found

◆ findNeighborsInASphereStartingFromCell()

int DgmOctree::findNeighborsInASphereStartingFromCell ( NearestNeighboursSearchStruct nNSS,
double  radius,
bool  sortValues = true 
) const

Advanced form of the nearest neighbours search algorithm (in a sphere)

This version is optimized for a spatially bounded search instead of a search bounded by a number of neighbours.

Warning
the number of points in the output buffer (nNSS.pointsInNeighbourhood) may be greater than the actual count of closest points inside the sphere! (which is returned by the method). Only the 'k' first points are actually inside the sphere (the others are not removed for the sake of performance).
Parameters
nNSSa pack of parameters
radiusthe sphere radius
sortValuesspecifies if the neighbours needs to be sorted by their distance to the query point or not
Returns
the number of neighbours found

◆ findPointNeighbourhood()

unsigned DgmOctree::findPointNeighbourhood ( const CCVector3 _queryPoint,
ReferenceCloud Yk,
unsigned  maxNumberOfNeighbors,
unsigned char  level,
double &  maxSquareDist,
double  maxSearchDist = 0,
int *  finalNeighbourhoodSize = nullptr 
) const

Finds the nearest neighbours around a query point.

This is the simplest form of the nearest neighbour search algorithm. It should only be used for unique/few requests as it is not optimized for repetitive search around points lying in the same octree cell (see DgmOctree::findNearestNeighborsStartingFromCell for example). Moreover, distances between each neighbour and the query aren't stored in this version of the algorithm.

Parameters
_queryPointthe query point
Ykthe nearest neighbours
maxNumberOfNeighborsthe maximal number of points to find
levelthe subdivision level of the octree at which to perform the search
maxSquareDistthe square distance between the farthest "nearest neighbour" and the query point
maxSearchDistthe maximum search distance (ignored if <= 0)
[out]finalNeighbourhoodSizethe final neighborhood (half)size (optional)
Returns
the number of neighbours found

◆ findTheNearestNeighborStartingFromCell()

double DgmOctree::findTheNearestNeighborStartingFromCell ( NearestNeighboursSearchStruct nNSS) const

Advanced form of the nearest neighbour search algorithm (unique neighbour)

This version is optimized for a unique nearest-neighbour search. See DgmOctree::NearestNeighboursSearchStruct for more details.

Parameters
nNSSNN search parameters
Returns
the square distance between the query point and its nearest neighbour (or -1 if none was found - i.e. maxSearchDist was reached)

◆ GenerateTruncatedCellCode()

DgmOctree::CellCode DgmOctree::GenerateTruncatedCellCode ( const Tuple3i cellPos,
unsigned char  level 
)
static

Generates the truncated cell code of a cell given its position at a given level of subdivision.

For a given level of subdivision (lets call it N), the cell position can be expressed as 3 integer coordinates between 0 and 2^N-1 (the number of cells along each dimension). This method computes the corresponding cell code, truncated at the level N (meaning that it is only valid for the Nth level, not for other levels).

Parameters
cellPosthe cell position
levelthe level of subdivision
Returns
the truncated cell code

◆ genericBuild()

int DgmOctree::genericBuild ( GenericProgressCallback progressCb = nullptr)
protected

Generic method to build the octree structure.

METHODS

Parameters
progressCbthe client application can get some notification of the process progress through this callback mechanism (see GenericProgressCallback)
Returns
the number of points projected in the octree

◆ GET_BIT_SHIFT()

unsigned char DgmOctree::GET_BIT_SHIFT ( unsigned char  level)
static

Returns the binary shift for a given level of subdivision.

This binary shift is used to truncate a full cell code in order to deduce the cell code for a given level of subdivision.

Parameters
levelthe level of subdivision
Returns
the binary shift

◆ getBoundingBox()

void DgmOctree::getBoundingBox ( CCVector3 bbMin,
CCVector3 bbMax 
) const

Returns the octree bounding box.

Method to request the octree bounding box limits

Parameters
bbMinlower bounding-box limits (Xmin,Ymin,Zmin)
bbMaxhigher bounding-box limits (Xmax,Ymax,Zmax)

◆ getCellCodes()

bool DgmOctree::getCellCodes ( unsigned char  level,
cellCodesContainer vec,
bool  truncatedCodes = false 
) const

Returns the list of codes corresponding to the octree cells for a given level of subdivision.

Only the non empty cells are represented in the octree structure.

Parameters
levelthe level of subdivision
vecthe list of codes
truncatedCodesindicates if the resulting codes should be truncated or not
Returns
false if an error occurred (e.g. not enough memory)

◆ getCellCodesAndIndexes()

bool DgmOctree::getCellCodesAndIndexes ( unsigned char  level,
cellsContainer vec,
bool  truncatedCodes = false 
) const

Returns the list of indexes and codes corresponding to the octree cells for a given level of subdivision.

Only the non empty cells are represented in the octree structure. Cell indexes are expressed relatively to the DgmOctree structure. They correspond to the indexes of the first points of each cell.

Parameters
levelthe level of subdivision
vecthe list of codes & indexes
truncatedCodesindicates if the resulting codes should be truncated or not
Returns
false if an error occurred (e.g. not enough memory)

◆ getCellDistanceFromBorders() [1/2]

void DgmOctree::getCellDistanceFromBorders ( const Tuple3i cellPos,
unsigned char  level,
int *  cellDists 
) const

Returns distance form a cell to the filled octree borders in all directions.

WARNING: distance values may be negative! (if cell is outside

Parameters
cellPoscell position
levellevel at which octree grid is considered
cellDistsoutput

◆ getCellDistanceFromBorders() [2/2]

void DgmOctree::getCellDistanceFromBorders ( const Tuple3i cellPos,
unsigned char  level,
int  neighbourhoodLength,
int *  cellDists 
) const

Returns distance from cell center to cell neighbourhood INSIDE filled octree.

WARNING: if cell neighbourhood is totally outside filled octree, the method returns false and cellDists is invalid.

Parameters
cellPoscenter cell position
levellevel at which octree grid is considered
neighbourhoodLengthcell neighbourhood "radius"
cellDistsoutput

◆ getCellIndex() [1/2]

unsigned DgmOctree::getCellIndex ( CellCode  truncatedCellCode,
unsigned char  bitShift 
) const

Returns the index of a given cell represented by its code.

The index is found thanks to a binary search. The index of an existing cell is between 0 and the number of points projected in the octree minus 1. If the cell code cannot be found in the octree structure, then the method returns an index equal to the number of projected points (m_numberOfProjectedPoints).

Parameters
truncatedCellCodetruncated cell code (i.e. original cell code shifted of 'bitShift' bits)
bitShiftbinary shift corresponding to the level of subdivision (see GET_BIT_SHIFT)
Returns
the index of the cell (or 'm_numberOfProjectedPoints' if none found)

◆ getCellIndex() [2/2]

unsigned DgmOctree::getCellIndex ( CellCode  truncatedCellCode,
unsigned char  bitShift,
unsigned  begin,
unsigned  end 
) const
protected

Returns the index of a given cell represented by its code.

Same algorithm as the other "getCellIndex" method, but in an optimized form. The binary search can be performed on a sub-part of the DgmOctree structure.

Parameters
truncatedCellCodetruncated cell code (i.e. original cell code shifted of 'bitShift' bits)
bitShiftbinary shift corresponding to the level of subdivision (see GET_BIT_SHIFT)
beginfirst index of the sub-list in which to perform the binary search
endlast index of the sub-list in which to perform the binary search
Returns
the index of the cell (or 'm_numberOfProjectedPoints' if none found)

◆ getCellIndexes()

bool DgmOctree::getCellIndexes ( unsigned char  level,
cellIndexesContainer vec 
) const

Returns the list of indexes corresponding to the octree cells for a given level of subdivision.

Only the non empty cells are represented in the octree structure. Cell indexes are expressed relatively to the DgmOctree structure. They correspond to the indexes of the first points of each cell.

Parameters
levelthe level of subdivision
vecthe list of indexes
Returns
false if an error occurred (e.g. not enough memory)

◆ getCellPos()

void DgmOctree::getCellPos ( CellCode  code,
unsigned char  level,
Tuple3i cellPos,
bool  isCodeTruncated 
) const

Returns the cell position for a given level of subdivision of a cell designated by its code.

Parameters
codethe cell code
levelthe level of subdivision
cellPosthe computed position
isCodeTruncatedindicates if the given code is truncated or not

◆ getCellSize()

const PointCoordinateType & CCCoreLib::DgmOctree::getCellSize ( unsigned char  level) const
inline

Returns the octree cells length for a given level of subdivision.

As the octree is cubical, cells are cubical.

Parameters
levelthe level of subdivision (up to MAX_OCTREE_LEVEL+1 for convenience)
Returns
the cell size

◆ getMaxFillIndexes()

const int * CCCoreLib::DgmOctree::getMaxFillIndexes ( unsigned char  level) const
inline

Returns the highest cell positions in the octree along all dimensions and for a given level of subdivision.

For example, at a level n, the octree length is 2^n cells along each dimension. The highest cell position along each dimension will be expressed between 0 and 2^n-1.

Parameters
levelthe level of subdivision
Returns
the highest cell position along X,Y and Z for a given level of subdivision

◆ getMinFillIndexes()

const int * CCCoreLib::DgmOctree::getMinFillIndexes ( unsigned char  level) const
inline

Returns the lowest cell positions in the octree along all dimensions and for a given level of subdivision.

For example, at a level n, the octree length is 2^n cells along each dimension. The lowest cell position along each dimension will be expressed between 0 and 2^n-1.

Parameters
levelthe level of subdivision
Returns
the lowest cell position along X,Y and Z for a given level of subdivision

◆ getNeighborCellsAround()

void DgmOctree::getNeighborCellsAround ( const Tuple3i cellPos,
cellIndexesContainer neighborCellsIndexes,
int  neighbourhoodLength,
unsigned char  level 
) const
protected

Returns the indexes of the neighbourhing (existing) cells of a given cell.

This function is used by the nearest neighbours search algorithms.

Parameters
cellPosthe query cell
neighborCellsIndexesthe found neighbourhing cells
neighbourhoodLengththe distance (in terms of cells) at which to look for neighbour cells
levelthe level of subdivision

◆ getNumberOfProjectedPoints()

unsigned CCCoreLib::DgmOctree::getNumberOfProjectedPoints ( ) const
inline

Returns the number of points projected into the octree.

Returns
the number of projected points

◆ getOctreeMaxs()

const CCVector3 & CCCoreLib::DgmOctree::getOctreeMaxs ( ) const
inline

Returns the higher boundaries of the octree.

Returns
the higher coordinates along X,Y and Z

◆ getOctreeMins()

const CCVector3 & CCCoreLib::DgmOctree::getOctreeMins ( ) const
inline

Returns the lower boundaries of the octree.

Returns
the lower coordinates along X,Y and Z

◆ getPointsInBoxNeighbourhood()

std::size_t DgmOctree::getPointsInBoxNeighbourhood ( BoxNeighbourhood params) const

Returns the points falling inside a box.

Warning
the 'squareDistd' field of each neighbour in the NeighboursSet structure is not used/set
Returns
the number of extracted points

◆ getPointsInCell()

bool DgmOctree::getPointsInCell ( CellCode  cellCode,
unsigned char  level,
ReferenceCloud subset,
bool  isCodeTruncated = false,
bool  clearOutputCloud = true 
) const

Returns the points lying in a specific cell.

Parameters
cellCodethe unique cell code
levelthe level of subdivision
[out]subsetset of points lying in the cell (references, no duplication)
isCodeTruncatedspecifies if the code is given in a truncated form or not
clearOutputCloudwhether to clear or not the output cloud (subest) if no points lie in the specified cell
Returns
success

◆ getPointsInCellByCellIndex()

bool DgmOctree::getPointsInCellByCellIndex ( ReferenceCloud cloud,
unsigned  cellIndex,
unsigned char  level,
bool  clearOutputCloud = true 
) const

Returns the points lying in a specific cell.

Each cell at a given level of subdivision can be recognized by the index in the DgmOctree structure of the first point that lies inside it. By construction, we are assured that every point lying in the same cell for a given level of subdivision are next to each others in the octree structure (which is the vector "m_thePointsAndTheirCellCodes" in practical). This is the quickest way to access the points inside a given cell (but its kind of hard to know directly what is the index of a given cell ;)

Parameters
cloudReferenceCloud to store the points lying inside the cell
cellIndexthe cell index
levelthe level of subdivision
clearOutputCloudwhether to clear the input cloud prior to inserting the points or not
Returns
success

◆ getPointsInCellsWithSortedCellCodes()

ReferenceCloud * DgmOctree::getPointsInCellsWithSortedCellCodes ( cellCodesContainer cellCodes,
unsigned char  level,
ReferenceCloud subset,
bool  areCodesTruncated = false 
) const

Returns the points lying in multiple cells.

Cells are recognized here by their unique "code". They should be sorted along by codes, with an ascendant order. See DgmOctree::getPointsInCellByCellIndex for more information.

Parameters
cellCodesthe cells codes
levelthe level of subdivision
[out]subsetset of points lying in the cell (references, no duplication)
areCodesTruncatedspecifies if the codes are given in a truncated form or not
Returns
the set of points lying in the cell (references, no duplication)

◆ getPointsInCylindricalNeighbourhood()

std::size_t DgmOctree::getPointsInCylindricalNeighbourhood ( CylindricalNeighbourhood params) const

Returns the points falling inside a cylinder.

Use findBestLevelForAGivenNeighbourhoodSizeExtraction to get the right value for 'level' (only once as it only depends on the radius value ;).

Warning
the 'squareDistd' field of each neighbour in the NeighboursSet structure is in fact the signed distance (not squared) of the point relatively to the cylinder's center and projected along its axis.
Parameters
paramsinput/output parameters structure
Returns
the number of extracted points

◆ getPointsInCylindricalNeighbourhoodProgressive()

std::size_t DgmOctree::getPointsInCylindricalNeighbourhoodProgressive ( ProgressiveCylindricalNeighbourhood params) const

Same as getPointsInCylindricalNeighbourhood with progressive approach.

Can be called multiple times (the 'currentHalfLength' parameter will increase each time until 'maxHalfLength' is reached).

◆ getPointsInNeighbourCellsAround()

void DgmOctree::getPointsInNeighbourCellsAround ( NearestNeighboursSearchStruct nNSS,
int  neighbourhoodLength,
bool  getOnlyPointsWithValidScalar = false 
) const
protected

Gets point in the neighbourhing cells of a specific cell.

Parameters
nNSSNN search parameters (from which are used: cellPos, pointsInNeighbourCells and level)
neighbourhoodLengththe new distance (in terms of cells) at which to look for neighbour cells
getOnlyPointsWithValidScalarwhether to ignore points having an invalid associated scalar value

◆ getPointsInSphericalNeighbourhood()

int DgmOctree::getPointsInSphericalNeighbourhood ( const CCVector3 sphereCenter,
PointCoordinateType  radius,
NeighboursSet neighbours,
unsigned char  level 
) const

Returns the points falling inside a sphere.

Use findBestLevelForAGivenNeighbourhoodSizeExtraction to get the right value for 'level' (only once as it only depends on the radius value ;).

Parameters
sphereCentercenter
radiusradius
[out]neighbourspoints falling inside the sphere
levelsubdivision level at which to apply the extraction process
Returns
the number of extracted points

◆ getTheCellPosWhichIncludesThePoint() [1/3]

void CCCoreLib::DgmOctree::getTheCellPosWhichIncludesThePoint ( const CCVector3 thePoint,
Tuple3i cellPos 
) const
inline

Returns the position FOR THE DEEPEST LEVEL OF SUBDIVISION of the cell that includes a given point.

The cell coordinates can be negative or greater than 2^MAX_OCTREE_LEVEL-1 as the point can lie outside the octree bounding-box.

Parameters
thePointthe query point
cellPosthe computed position

◆ getTheCellPosWhichIncludesThePoint() [2/3]

void CCCoreLib::DgmOctree::getTheCellPosWhichIncludesThePoint ( const CCVector3 thePoint,
Tuple3i cellPos,
unsigned char  level 
) const
inline

Returns the position for a given level of subdivision of the cell that includes a given point.

The cell coordinates can be negative or greater than 2^N-1 (where N is the level of subdivision) as the point can lie outside the octree bounding-box.

Parameters
thePointthe query point
cellPosthe computed position
levelthe level of subdivision

◆ getTheCellPosWhichIncludesThePoint() [3/3]

void CCCoreLib::DgmOctree::getTheCellPosWhichIncludesThePoint ( const CCVector3 thePoint,
Tuple3i cellPos,
unsigned char  level,
bool &  inBounds 
) const
inline

Returns the position for a given level of subdivision of the cell that includes a given point.

The cell coordinates can be negative or greater than 2^N-1 (where N is the level of subdivision) as the point can lie outside the octree bounding-box. In this version, method indicates if the query point is inside ("inbounds") or outside the octree bounding-box.

Parameters
thePointthe query point
cellPosthe computed position
levelthe level of subdivision
inBoundsindicates if the query point is inside or outside the octree bounding-box

◆ OCTREE_LENGTH()

int DgmOctree::OCTREE_LENGTH ( int  level)
static

Returns the octree length (in term of cells) for a given level of subdivision.

Parameters
levelthe level of subdivision
Returns
2^level

Member Data Documentation

◆ INVALID_CELL_CODE

const CellCode CCCoreLib::DgmOctree::INVALID_CELL_CODE = (~static_cast<CellCode>(0))
static

Invalid cell code.

Warning
Never pass a 'constant initializer' by reference

◆ m_thePointsAndTheirCellCodes

cellsContainer CCCoreLib::DgmOctree::m_thePointsAndTheirCellCodes
protected

The coded octree structure.

ATTRIBUTES

◆ MAX_OCTREE_LENGTH

const int CCCoreLib::DgmOctree::MAX_OCTREE_LENGTH = (1 << MAX_OCTREE_LEVEL)
static

Max octree length at last level of subdivision (number of cells)

Warning
Never pass a 'constant initializer' by reference

◆ MAX_OCTREE_LEVEL

const int CCCoreLib::DgmOctree::MAX_OCTREE_LEVEL = 10
static

Max octree subdivision level.

STRUCTURES
Number of bits used to code the cells position: 3*MAX_OCTREE_LEVEL

Warning
Never pass a 'constant initializer' by reference

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