![]() |
ACloudViewer
3.9.4
A Modern Library for 3D Data Processing
|
The octree structure used throughout the library. More...
#include <DgmOctree.h>


Classes | |
| struct | BoxNeighbourhood |
| Input/output parameters structure for getPointsInBoxNeighbourhood. More... | |
| struct | CellDescriptor |
| Structure used during nearest neighbour search. More... | |
| struct | CylindricalNeighbourhood |
| 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 | NearestNeighboursSphericalSearchStruct |
| 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 |
Public Types | |
| enum | RayCastProcess { RC_NEAREST_POINT , RC_CLOSE_POINTS } |
| Ray casting processes. More... | |
| using | CellCode = unsigned |
| Type of the code of an octree cell. More... | |
| using | cellCodesContainer = std::vector< CellCode > |
| Octree cell codes container. More... | |
| using | cellIndexesContainer = std::vector< unsigned int > |
| Octree cell indexes container. More... | |
| using | NeighboursSet = std::vector< PointDescriptor > |
| A set of neighbours. More... | |
| using | NeighbourCellsSet = std::vector< CellDescriptor > |
| A set of neighbour cells descriptors. More... | |
| using | cellsContainer = std::vector< IndexAndCode > |
| Container of 'IndexAndCode' structures. More... | |
| using | octreeCellFunc = bool(*)(const octreeCell &, void **, NormalizedProgress *) |
Public Member Functions | |
| DgmOctree (GenericIndexedCloudPersist *cloud) | |
| DgmOctree constructor. More... | |
| ~DgmOctree () override=default | |
| DgmOctree destructor. More... | |
| virtual void | clear () |
| Clears the octree. More... | |
| 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 CCVector3 & | getOctreeMins () const |
| Returns the lower boundaries of the octree. More... | |
| const CCVector3 & | getOctreeMaxs () 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 |
| const int * | getMaxFillIndexes (unsigned char level) const |
| 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 |
| void | getCellDistanceFromBorders (const Tuple3i &cellPos, unsigned char level, int neighbourhoodLength, int *cellDists) const |
| 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... | |
| ReferenceCloud * | getPointsInCellsWithSortedCellCodes (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 |
| unsigned | findNearestNeighborsStartingFromCell (NearestNeighboursSearchStruct &nNSS, bool getOnlyPointsWithValidScalar=false) const |
| 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 ¶ms) const |
| Returns the points falling inside a cylinder. More... | |
| std::size_t | getPointsInCylindricalNeighbourhoodProgressive (ProgressiveCylindricalNeighbourhood ¶ms) const |
| Same as getPointsInCylindricalNeighbourhood with progressive approach. More... | |
| std::size_t | getPointsInBoxNeighbourhood (BoxNeighbourhood ¶ms) const |
| Returns the points falling inside a box. More... | |
| void | getTheCellPosWhichIncludesThePoint (const CCVector3 *thePoint, Tuple3i &cellPos) const |
| void | getTheCellPosWhichIncludesThePoint (const CCVector3 *thePoint, Tuple3i &cellPos, unsigned char level) const |
| void | getTheCellPosWhichIncludesThePoint (const CCVector3 *thePoint, Tuple3i &cellPos, unsigned char level, bool &inBounds) const |
| void | getCellPos (CellCode code, unsigned char level, Tuple3i &cellPos, bool isCodeTruncated) const |
| void | computeCellCenter (CellCode code, unsigned char level, CCVector3 ¢er, bool isCodeTruncated=false) const |
| void | computeCellCenter (const Tuple3i &cellPos, unsigned char level, CCVector3 ¢er) const |
| void | computeCellCenter (const Tuple3s &cellPos, unsigned char level, CCVector3 ¢er) const |
| Short version of computeCellCenter. More... | |
| 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 char | findBestLevelForAGivenNeighbourhoodSizeExtraction (PointCoordinateType radius) const |
| unsigned char | findBestLevelForComparisonWithOctree (const DgmOctree *theOtherOctree) const |
| unsigned char | findBestLevelForAGivenPopulationPerCell (unsigned indicativeNumberOfPointsPerCell) const |
| unsigned char | findBestLevelForAGivenCellNumber (unsigned indicativeNumberOfCells) const |
| const CellCode & | getCellCode (unsigned index) const |
| Returns the ith cell code. More... | |
| bool | getCellCodes (unsigned char level, cellCodesContainer &vec, bool truncatedCodes=false) const |
| bool | getCellIndexes (unsigned char level, cellIndexesContainer &vec) const |
| bool | getCellCodesAndIndexes (unsigned char level, cellsContainer &vec, bool truncatedCodes=false) const |
| void | diff (const cellCodesContainer &codesA, const cellCodesContainer &codesB, cellCodesContainer &diffA, cellCodesContainer &diffB) const |
| bool | diff (unsigned char octreeLevel, const cellsContainer &codesA, const cellsContainer &codesB, int &diffA, int &diffB, int &cellsA, int &cellsB) const |
| const unsigned & | getCellNumber (unsigned char level) const |
| Returns the number of cells for a given level of subdivision. More... | |
| double | computeMeanOctreeDensity (unsigned char level) const |
| int | extractCCs (const cellCodesContainer &cellCodes, unsigned char level, bool sixConnexity, GenericProgressCallback *progressCb=nullptr) const |
| int | extractCCs (unsigned char level, bool sixConnexity, GenericProgressCallback *progressCb=nullptr) const |
| 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) |
| unsigned | executeFunctionForAllCellsAtLevel (unsigned char level, octreeCellFunc func, void **additionalParameters, bool multiThread=false, GenericProgressCallback *progressCb=nullptr, const char *functionTitle=nullptr, int maxThreadCount=0) |
| bool | rayCast (const CCVector3 &rayAxis, const CCVector3 &rayOrigin, double maxRadiusOrFov, bool isFOV, RayCastProcess process, std::vector< PointDescriptor > &output) const |
| Ray casting algorithm. More... | |
| GenericIndexedCloudPersist * | associatedCloud () const |
| Returns the associated cloud. More... | |
| const cellsContainer & | pointsAndTheirCellCodes () const |
| Returns the octree 'structure'. More... | |
Public Member Functions inherited from cloudViewer::GenericOctree | |
| virtual | ~GenericOctree ()=default |
| Default destructor. More... | |
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) |
| static CellCode | GenerateTruncatedCellCode (const Tuple3i &cellPos, unsigned char level) |
| static CellCode | GenerateTruncatedCellCode (const Tuple3s &pos, unsigned char level) |
| Short version of GenerateTruncatedCellCode. More... | |
| static PointCoordinateType | ComputeMinDistanceToCellBorder (const CCVector3 &queryPoint, PointCoordinateType cs, const CCVector3 &cellCenter) |
| static bool | MultiThreadSupport () |
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 | updateMinAndMaxTables () |
| Updates the tables containing octree limits and boundaries. More... | |
| void | updateCellSizeTable () |
| void | updateCellCountTable () |
| 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 |
| 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 bitDec) const |
| Returns the index of a given cell represented by its code. More... | |
| unsigned | getCellIndex (CellCode truncatedCellCode, unsigned char bitDec, 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... | |
| GenericIndexedCloudPersist * | m_theAssociatedCloud |
| Associated cloud. More... | |
| unsigned | m_numberOfProjectedPoints |
| Number of points projected in the octree. More... | |
| unsigned | m_nearestPow2 |
| CCVector3 | m_dimMin |
| Min coordinates of the octree bounding-box. More... | |
| CCVector3 | m_dimMax |
| Max coordinates of the octree bounding-box. More... | |
| CCVector3 | m_pointsMin |
| CCVector3 | m_pointsMax |
| PointCoordinateType | m_cellSize [MAX_OCTREE_LEVEL+2] |
| Cell dimensions for all subdivision levels. More... | |
| int | m_fillIndexes [(MAX_OCTREE_LEVEL+1) *6] |
| unsigned | m_cellCount [MAX_OCTREE_LEVEL+1] |
| Number of cells per level of subdivision. More... | |
| unsigned | m_maxCellPopulation [MAX_OCTREE_LEVEL+1] |
| Max cell population per level of subdivision. More... | |
| double | m_averageCellPopulation [MAX_OCTREE_LEVEL+1] |
| Average cell population per level of subdivision. More... | |
| double | m_stdDevCellPopulation [MAX_OCTREE_LEVEL+1] |
| Std. dev. of cell population per level of subdivision. More... | |
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).
Definition at line 39 of file DgmOctree.h.
| using cloudViewer::DgmOctree::CellCode = unsigned |
Type of the code of an octree cell.
Definition at line 78 of file DgmOctree.h.
| using cloudViewer::DgmOctree::cellCodesContainer = std::vector<CellCode> |
Octree cell codes container.
Definition at line 92 of file DgmOctree.h.
| using cloudViewer::DgmOctree::cellIndexesContainer = std::vector<unsigned int> |
Octree cell indexes container.
Definition at line 95 of file DgmOctree.h.
| using cloudViewer::DgmOctree::cellsContainer = std::vector<IndexAndCode> |
Container of 'IndexAndCode' structures.
Definition at line 351 of file DgmOctree.h.
| using cloudViewer::DgmOctree::NeighbourCellsSet = std::vector<CellDescriptor> |
A set of neighbour cells descriptors.
Definition at line 150 of file DgmOctree.h.
| using cloudViewer::DgmOctree::NeighboursSet = std::vector<PointDescriptor> |
A set of neighbours.
Definition at line 133 of file DgmOctree.h.
| using cloudViewer::DgmOctree::octreeCellFunc = bool (*)(const octreeCell&, void**, NormalizedProgress*) |
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:
Definition at line 393 of file DgmOctree.h.
Ray casting processes.
| Enumerator | |
|---|---|
| RC_NEAREST_POINT | |
| RC_CLOSE_POINTS | |
Definition at line 1179 of file DgmOctree.h.
|
explicit |
DgmOctree constructor.
METHODS
| cloud | the cloud to construct the octree on |
Definition at line 174 of file DgmOctree.cpp.
References clear(), and m_theAssociatedCloud.
|
overridedefault |
DgmOctree destructor.
|
inline |
Returns the associated cloud.
Definition at line 1191 of file DgmOctree.h.
Referenced by FastMarchingForFacetExtraction::addCellToCurrentFacet(), cloudViewer::DistanceComputationTools::computeApproxCloud2CloudDistance(), cloudViewer::DistanceComputationTools::computeCloud2MeshDistancesWithOctree(), cloudViewer::DistanceComputationTools::computeCloud2MeshDistanceWithOctree(), FastMarchingForFacetExtraction::computeTCoefApprox(), FastMarchingForFacetExtraction::init(), cloudViewer::FastMarchingForPropagation::init(), cloudViewer::CloudSamplingTools::resampleCloudSpatially(), cloudViewer::FastMarchingForPropagation::setPropagationTimingsAsDistances(), FastMarchingForFacetExtraction::setSeedCell(), and FastMarchingForFacetExtraction::updateFlagsTable().
| 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.
| octreeMin | the lower limits for the octree cells along X, Y and Z |
| octreeMax | the upper limits for the octree cells along X, Y and Z |
| pointsMinFilter | the lower limits for the projected points along X, Y and Z (is specified) |
| pointsMaxFilter | the upper limits for the projected points along X, Y and Z (is specified) |
| progressCb | the client application can get some notification of the process progress through this callback mechanism (see GenericProgressCallback) |
Definition at line 204 of file DgmOctree.cpp.
References clear(), genericBuild(), m_dimMax, m_dimMin, m_pointsMax, m_pointsMin, and m_thePointsAndTheirCellCodes.
| int DgmOctree::build | ( | GenericProgressCallback * | progressCb = nullptr | ) |
Builds the structure.
Octree 3D limits are determined automatically.
| progressCb | the client application can get some notification of the process progress through this callback mechanism (see GenericProgressCallback) |
Definition at line 196 of file DgmOctree.cpp.
References clear(), genericBuild(), m_thePointsAndTheirCellCodes, and updateMinAndMaxTables().
Referenced by cloudViewer::ScalarFieldTools::applyScalarFieldGaussianFilter(), qCanupoProcess::Classify(), cloudViewer::GeometricalAnalysisTools::ComputeCharactersitic(), cloudViewer::DistanceComputationTools::computeCloud2MeshDistances(), qCanupoTools::ComputeCorePointsDescriptors(), qM3C2Normals::ComputeCorePointsNormals(), cloudViewer::DistanceComputationTools::computeGeodesicDistances(), cloudViewer::GeometricalAnalysisTools::ComputeLocalDensityApprox(), ccEntityAction::computeOctree(), cloudViewer::ScalarFieldTools::computeScalarFieldGradient(), FastMarchingForFacetExtraction::ExtractPlanarFacets(), cloudViewer::GeometricalAnalysisTools::FlagDuplicatePoints(), cloudViewer::AutoSegmentationTools::frontPropagationBasedSegmentation(), cloudViewer::AutoSegmentationTools::labelConnectedComponents(), cloudViewer::CloudSamplingTools::noiseFilter(), cloudViewer::ICPRegistrationTools::Register(), cloudViewer::CloudSamplingTools::resampleCloudSpatially(), cloudViewer::CloudSamplingTools::resampleCloudWithOctree(), cloudViewer::CloudSamplingTools::resampleCloudWithOctreeAtLevel(), cloudViewer::CloudSamplingTools::sorFilter(), cloudViewer::CloudSamplingTools::subsampleCloudWithOctree(), cloudViewer::CloudSamplingTools::subsampleCloudWithOctreeAtLevel(), cloudViewer::DistanceComputationTools::synchronizeOctrees(), and cloudViewer::StatisticalTestingTools::testCloudWithStatisticalModel().
|
virtual |
Clears the octree.
Reimplemented in ccOctree.
Definition at line 183 of file DgmOctree.cpp.
References m_cellSize, m_dimMax, m_dimMin, m_fillIndexes, m_nearestPow2, m_numberOfProjectedPoints, m_pointsMax, m_pointsMin, m_thePointsAndTheirCellCodes, MAX_OCTREE_LEVEL, and updateCellCountTable().
Referenced by build(), ccEntityAction::computeOctree(), DgmOctree(), and cloudViewer::DistanceComputationTools::synchronizeOctrees().
|
inline |
Returns the cell center for a given level of subdivision of a cell designated by its code
| code | the cell code |
| level | the level of subdivision |
| center | the computed center |
| isCodeTruncated | indicates if the given code is truncated or not |
Definition at line 862 of file DgmOctree.h.
Referenced by cloudViewer::CloudSamplingTools::applyNoiseFilterAtLevel(), cloudViewer::CloudSamplingTools::applySORFilterAtLevel(), cloudViewer::GeometricalAnalysisTools::ComputeApproxPointsDensityInACellAtLevel(), cloudViewer::ScalarFieldTools::computeCellGaussianFilter(), cloudViewer::DistanceComputationTools::computeCellHausdorffDistance(), cloudViewer::DistanceComputationTools::computeCellHausdorffDistanceWithLocalModel(), cloudViewer::DistanceComputationTools::computeCloud2MeshDistanceWithOctree(), cloudViewer::GeometricalAnalysisTools::ComputeGeomCharacteristicAtLevel(), cloudViewer::StatisticalTestingTools::computeLocalChi2DistAtLevel(), cloudViewer::ScalarFieldTools::computeMeanGradientOnPatch(), ComputeNeighborhood2MeshDistancesWithOctree(), findPointNeighbourhood(), cloudViewer::GeometricalAnalysisTools::FlagDuplicatePointsInACellAtLevel(), cloudViewer::DistanceComputationTools::intersectMeshWithOctree(), masc::Tools::PrepareFeatures(), cloudViewer::CloudSamplingTools::resampleCellAtLevel(), and cloudViewer::CloudSamplingTools::subsampleCellAtLevel().
|
inline |
Returns the cell center for a given level of subdivision of a cell designated by its position
| cellPos | the cell position |
| level | the level of subdivision |
| center | the computed center |
Definition at line 878 of file DgmOctree.h.
References Tuple3Tpl< Type >::x, Tuple3Tpl< Type >::y, and Tuple3Tpl< Type >::z.
|
inline |
Short version of computeCellCenter.
Definition at line 895 of file DgmOctree.h.
References Tuple3Tpl< Type >::x, Tuple3Tpl< Type >::y, and Tuple3Tpl< Type >::z.
| void DgmOctree::computeCellLimits | ( | CellCode | code, |
| unsigned char | level, | ||
| CCVector3 & | cellMin, | ||
| CCVector3 & | cellMax, | ||
| bool | isCodeTruncated = false |
||
| ) | const |
Returns the spatial limits of a given cell.
| code | the cell code |
| level | the level of subdivision |
| cellMin | the minimum coordinates along each dimension |
| cellMax | the maximum coordinates along each dimension |
| isCodeTruncated | indicates if the given code is truncated or not |
Definition at line 520 of file DgmOctree.cpp.
References getCellPos(), getCellSize(), m_dimMin, Tuple3Tpl< Type >::x, Tuple3Tpl< Type >::y, and Tuple3Tpl< Type >::z.
|
protected |
Computes statistics about cells for a given level of subdivision.
This method requires some computation, therefore it shouldn't be called too often.
| level | the level of subdivision |
Definition at line 423 of file DgmOctree.cpp.
References GET_BIT_SHIFT(), m_averageCellPopulation, m_cellCount, m_maxCellPopulation, m_stdDevCellPopulation, m_thePointsAndTheirCellCodes, and MAX_OCTREE_LEVEL.
Referenced by updateCellCountTable().
| double DgmOctree::computeMeanOctreeDensity | ( | unsigned char | level | ) | const |
Computes mean octree density (point/cell) at a given level of subdivision
| level | the level of subdivision |
Definition at line 2789 of file DgmOctree.cpp.
References getCellNumber(), and m_numberOfProjectedPoints.
|
inlinestatic |
Computes the minimal distance between a point and the borders (faces) of the cell (cube) in which it is included
| queryPoint | the point |
| cs | the cell size (as cells are cubical, it's the same along every dimension) |
| cellCenter | the cell center |
Definition at line 1057 of file DgmOctree.h.
References abs(), Tuple3Tpl< Type >::x, Tuple3Tpl< Type >::y, and Tuple3Tpl< Type >::z.
Referenced by cloudViewer::DistanceComputationTools::computeCloud2MeshDistanceWithOctree(), ComputeNeighborhood2MeshDistancesWithOctree(), findNearestNeighborsStartingFromCell(), findNeighborsInASphereStartingFromCell(), and findTheNearestNeighborStartingFromCell().
| 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.
| codesA | the cell codes of the first octree for the implicit level |
| codesB | the cell codes of the second octree for the implicit level |
| diffA | the cells of the first octree that are not in the second octree |
| diffB | the cells of the second octree that are not in the first octree |
Definition at line 2955 of file DgmOctree.cpp.
Referenced by findBestLevelForComparisonWithOctree().
| 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.
| octreeLevel | the octree level |
| codesA | the cell codes (and point index) of the first octree |
| codesB | the cell codes (and point index) of the second octree |
| diffA | the number of cells of the first octree that are not in the second octree |
| diffB | the number of cells of the second octree that are not in the first octree |
| cellsA | the number of cells of the first octree for the given number of subdivision |
| cellsB | the number of cells of the second octree for the given number of subdivision |
Definition at line 2980 of file DgmOctree.cpp.
References GET_BIT_SHIFT(), and octreeLevel.
| 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 QtConcurrent::map system.
\param level the level of subdivision
\param func the function to apply
\param additionalParameters the function parameters
\param multiThread whether to use parallel processing or not
\param progressCb the client application can get some notification
of the process progress through this callback mechanism (see GenericProgressCallback)
| functionTitle | function title |
| maxThreadCount | the maximum number of threads to use (0 = all). Ignored if 'multiThread' is false. |
Definition at line 3573 of file DgmOctree.cpp.
References cloudViewer::ReferenceCloud::addPointIndex(), GET_BIT_SHIFT(), getCellNumber(), cloudViewer::DgmOctree::octreeCell::index, cloudViewer::DgmOctree::octreeCell::level, m_averageCellPopulation, m_maxCellPopulation, m_stdDevCellPopulation, m_theAssociatedCloud, m_thePointsAndTheirCellCodes, max(), cloudViewer::DgmOctree::octreeCell::points, cloudViewer::ReferenceCloud::reserve(), result, cloudViewer::GenericProgressCallback::setInfo(), cloudViewer::GenericProgressCallback::setMethodTitle(), cloudViewer::GenericCloud::size(), cloudViewer::GenericProgressCallback::start(), cloudViewer::GenericProgressCallback::stop(), cloudViewer::GenericProgressCallback::textCanBeEdited(), cloudViewer::DgmOctree::octreeCell::truncatedCode, and cloudViewer::GenericProgressCallback::update().
Referenced by cloudViewer::ScalarFieldTools::applyScalarFieldGaussianFilter(), cloudViewer::GeometricalAnalysisTools::ComputeCharactersitic(), cloudViewer::DistanceComputationTools::computeCloud2CloudDistances(), cloudViewer::GeometricalAnalysisTools::ComputeLocalDensityApprox(), cloudViewer::ScalarFieldTools::computeScalarFieldGradient(), cloudViewer::GeometricalAnalysisTools::FlagDuplicatePoints(), cloudViewer::CloudSamplingTools::noiseFilter(), cloudViewer::CloudSamplingTools::resampleCloudWithOctreeAtLevel(), cloudViewer::CloudSamplingTools::sorFilter(), and cloudViewer::CloudSamplingTools::subsampleCloudWithOctreeAtLevel().
| 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 QtConcurrent::map system.
\param startingLevel the initial level of subdivision
\param func the function to apply
\param additionalParameters the function parameters
\param minNumberOfPointsPerCell minimal number of points inside
a cell (indicative)
| maxNumberOfPointsPerCell | maximum number of points inside a cell (indicative) |
| multiThread | whether to use parallel processing or not |
| progressCb | the client application can get some notification of the process progress through this callback mechanism (see GenericProgressCallback) |
| functionTitle | function title |
| maxThreadCount | the maximum number of threads to use (0 = all). Ignored if 'multiThread' is false. |
Definition at line 3820 of file DgmOctree.cpp.
References cloudViewer::ReferenceCloud::addPointIndex(), cloudViewer::ReferenceCloud::clear(), GET_BIT_SHIFT(), getCellNumber(), cloudViewer::DgmOctree::octreeCell::index, cloudViewer::GenericProgressCallback::isCancelRequested(), cloudViewer::DgmOctree::octreeCell::level, m_averageCellPopulation, m_maxCellPopulation, m_numberOfProjectedPoints, m_stdDevCellPopulation, m_theAssociatedCloud, m_thePointsAndTheirCellCodes, max(), MAX_OCTREE_LEVEL, nProgress, cloudViewer::DgmOctree::octreeCell::points, cloudViewer::ReferenceCloud::reserve(), result, cloudViewer::GenericProgressCallback::setInfo(), cloudViewer::GenericProgressCallback::setMethodTitle(), cloudViewer::GenericCloud::size(), cloudViewer::GenericProgressCallback::start(), stddev(), cloudViewer::GenericProgressCallback::stop(), cloudViewer::GenericProgressCallback::textCanBeEdited(), cloudViewer::DgmOctree::octreeCell::truncatedCode, and cloudViewer::GenericProgressCallback::update().
Referenced by cloudViewer::StatisticalTestingTools::testCloudWithStatisticalModel().
| 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).
| cellCodes | the cell codes to consider for the CC computation |
| level | the level of subdivision at which to perform the algorithm |
| sixConnexity | indicates if the CC's 3D connexity should be 6 (26 otherwise) |
| progressCb | the client application can get some notification of the process progress through this callback mechanism (see GenericProgressCallback) |
Definition at line 3093 of file DgmOctree.cpp.
References cloudViewer::ReferenceCloud::forwardIterator(), GET_BIT_SHIFT(), getCellPos(), getPointsInCell(), IndexAndCodeExt::indexComp(), m_theAssociatedCloud, cloudViewer::NormalizedProgress::oneStep(), ParallelSort, cloudViewer::ReferenceCloud::placeIteratorAtBeginning(), cloudViewer::ReferenceCloud::setCurrentPointScalarValue(), cloudViewer::GenericProgressCallback::setInfo(), cloudViewer::GenericProgressCallback::setMethodTitle(), cloudViewer::ReferenceCloud::size(), cloudViewer::GenericProgressCallback::start(), cloudViewer::GenericProgressCallback::stop(), std::swap(), cloudViewer::GenericProgressCallback::textCanBeEdited(), Tuple3Tpl< Type >::u, cloudViewer::GenericProgressCallback::update(), Tuple3Tpl< Type >::x, Tuple3Tpl< Type >::y, and Tuple3Tpl< Type >::z.
Referenced by extractCCs(), and cloudViewer::AutoSegmentationTools::labelConnectedComponents().
| 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.
| level | the level of subdivision at which to perform the algorithm |
| sixConnexity | indicates if the CC's 3D connexity should be 6 (26 otherwise) |
| progressCb | the client application can get some notification of the process progress through this callback mechanism (see GenericProgressCallback) |
Definition at line 3058 of file DgmOctree.cpp.
References extractCCs(), and getCellCodes().
| unsigned char DgmOctree::findBestLevelForAGivenCellNumber | ( | unsigned | indicativeNumberOfCells | ) | const |
Determines the best subdivision level of the octree to match a given number of cells
| indicativeNumberOfCells | 'desired' number of cells |
Definition at line 2766 of file DgmOctree.cpp.
References abs(), getCellNumber(), and MAX_OCTREE_LEVEL.
Referenced by cloudViewer::CloudSamplingTools::resampleCloudWithOctree(), and cloudViewer::CloudSamplingTools::subsampleCloudWithOctree().
| 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
| radius | the sphere radius |
Definition at line 2664 of file DgmOctree.cpp.
References getCellSize(), m_averageCellPopulation, and MAX_OCTREE_LEVEL.
Referenced by cloudViewer::ScalarFieldTools::applyScalarFieldGaussianFilter(), qCanupoProcess::Classify(), cloudViewer::GeometricalAnalysisTools::ComputeCharactersitic(), qCanupoTools::ComputeCorePointsDescriptors(), qM3C2Normals::ComputeCorePointsNormals(), cloudViewer::ScalarFieldTools::computeScalarFieldGradient(), cloudViewer::GeometricalAnalysisTools::FlagDuplicatePoints(), cloudViewer::CloudSamplingTools::noiseFilter(), masc::Tools::PrepareFeatures(), cloudViewer::CloudSamplingTools::resampleCloudSpatially(), and qCanupoTools::TrainClassifier().
| 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
| indicativeNumberOfPointsPerCell | 'desired' average number of points per cell |
Definition at line 2737 of file DgmOctree.cpp.
References m_averageCellPopulation, and MAX_OCTREE_LEVEL.
Referenced by cloudViewer::GeometricalAnalysisTools::ComputeLocalDensityApprox(), ComputeMathOpWithNearestNeighbor(), cloudViewer::ScalarFieldTools::computeScalarFieldGradient(), cloudViewer::CloudSamplingTools::noiseFilter(), rayCast(), cloudViewer::CloudSamplingTools::sorFilter(), and cloudViewer::StatisticalTestingTools::testCloudWithStatisticalModel().
| 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.
| theOtherOctree | the octree of the other cloud |
Definition at line 2691 of file DgmOctree.cpp.
References diff(), getNumberOfProjectedPoints(), m_thePointsAndTheirCellCodes, max(), MAX_OCTREE_LEVEL, and min().
Referenced by cloudViewer::DistanceComputationTools::computeCloud2CloudDistances(), and cloudViewer::ICPRegistrationTools::Register().
| 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.
| nNSS | NN search parameters |
| getOnlyPointsWithValidScalar | whether to ignore points having an invalid associated scalar value |
Definition at line 1655 of file DgmOctree.cpp.
References cloudViewer::DgmOctree::NearestNeighboursSearchStruct::alreadyVisitedNeighbourhoodSize, cloudViewer::utility::ceil(), cloudViewer::DgmOctree::NearestNeighboursSearchStruct::cellCenter, cloudViewer::DgmOctree::NearestNeighboursSearchStruct::cellPos, ComputeMinDistanceToCellBorder(), cloudViewer::DgmOctree::PointDescriptor::distComp(), GenerateTruncatedCellCode(), GET_BIT_SHIFT(), getCellIndex(), getCellSize(), cloudViewer::GenericIndexedCloudPersist::getPointPersistentPtr(), cloudViewer::GenericCloud::getPointScalarValue(), getPointsInNeighbourCellsAround(), INVALID_CELL_CODE, cloudViewer::DgmOctree::NearestNeighboursSearchStruct::level, m_fillIndexes, m_numberOfProjectedPoints, m_theAssociatedCloud, m_thePointsAndTheirCellCodes, max(), cloudViewer::DgmOctree::NearestNeighboursSearchStruct::maxSearchSquareDistd, cloudViewer::DgmOctree::NearestNeighboursSearchStruct::minNumberOfNeighbors, cloudViewer::DgmOctree::NearestNeighboursSearchStruct::pointsInNeighbourhood, cloudViewer::DgmOctree::NearestNeighboursSearchStruct::queryPoint, std::swap(), Tuple3Tpl< Type >::u, and cloudViewer::ScalarField::ValidValue().
Referenced by cloudViewer::CloudSamplingTools::applyNoiseFilterAtLevel(), cloudViewer::CloudSamplingTools::applySORFilterAtLevel(), cloudViewer::GeometricalAnalysisTools::ComputeApproxPointsDensityInACellAtLevel(), cloudViewer::DistanceComputationTools::computeCellHausdorffDistanceWithLocalModel(), cloudViewer::StatisticalTestingTools::computeLocalChi2DistAtLevel(), and findPointNeighbourhood().
| 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.
| nNSS | a pack of parameters |
| radius | the sphere radius |
| sortValues | specifies if the neighbours needs to be sorted by their distance to the query point or not |
TEST_CELLS_FOR_SPHERICAL_NN
Definition at line 2479 of file DgmOctree.cpp.
References cloudViewer::DgmOctree::NearestNeighboursSearchStruct::alreadyVisitedNeighbourhoodSize, cloudViewer::utility::ceil(), cloudViewer::DgmOctree::NearestNeighboursSearchStruct::cellCenter, ComputeMinDistanceToCellBorder(), copy, count, cloudViewer::DgmOctree::PointDescriptor::distComp(), getCellSize(), getPointsInNeighbourCellsAround(), cloudViewer::DgmOctree::NearestNeighboursSearchStruct::level, cloudViewer::DgmOctree::PointDescriptor::point, cloudViewer::DgmOctree::NearestNeighboursSearchStruct::pointsInNeighbourhood, cloudViewer::DgmOctree::NearestNeighboursSearchStruct::queryPoint, SQRT_3, cloudViewer::DgmOctree::PointDescriptor::squareDistd, and std::swap().
Referenced by cloudViewer::CloudSamplingTools::applyNoiseFilterAtLevel(), cloudViewer::ScalarFieldTools::computeCellGaussianFilter(), cloudViewer::DistanceComputationTools::computeCellHausdorffDistanceWithLocalModel(), cloudViewer::GeometricalAnalysisTools::ComputeGeomCharacteristicAtLevel(), cloudViewer::ScalarFieldTools::computeMeanGradientOnPatch(), cloudViewer::GeometricalAnalysisTools::FlagDuplicatePointsInACellAtLevel(), and masc::Tools::PrepareFeatures().
| 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.
| _queryPoint | the query point | |
| Yk | the nearest neighbours | |
| maxNumberOfNeighbors | the maximal number of points to find | |
| level | the subdivision level of the octree at which to perform the search | |
| maxSquareDist | the square distance between the farthest "nearest neighbour" and the query point | |
| maxSearchDist | the maximum search distance (ignored if <= 0) | |
| [out] | the | final neighborhood (half)size (optional) |
Definition at line 721 of file DgmOctree.cpp.
References cloudViewer::ReferenceCloud::addPointIndex(), cloudViewer::DgmOctree::NearestNeighboursSearchStruct::alreadyVisitedNeighbourhoodSize, cloudViewer::DgmOctree::NearestNeighboursSearchStruct::cellCenter, cloudViewer::DgmOctree::NearestNeighboursSearchStruct::cellPos, computeCellCenter(), findNearestNeighborsStartingFromCell(), findTheNearestNeighborStartingFromCell(), getTheCellPosWhichIncludesThePoint(), cloudViewer::DgmOctree::NearestNeighboursSearchStruct::level, cloudViewer::DgmOctree::NearestNeighboursSearchStruct::maxSearchSquareDistd, min(), cloudViewer::DgmOctree::NearestNeighboursSearchStruct::minNumberOfNeighbors, cloudViewer::DgmOctree::NearestNeighboursSearchStruct::pointsInNeighbourhood, cloudViewer::DgmOctree::NearestNeighboursSearchStruct::queryPoint, and cloudViewer::DgmOctree::NearestNeighboursSearchStruct::theNearestPointIndex.
Referenced by qCanupoProcess::Classify(), ComputeMathOpWithNearestNeighbor(), and RefinePointClassif().
| 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.
| nNSS | NN search parameters |
Definition at line 1468 of file DgmOctree.cpp.
References cloudViewer::DgmOctree::NearestNeighboursSearchStruct::alreadyVisitedNeighbourhoodSize, cloudViewer::utility::ceil(), cloudViewer::DgmOctree::NearestNeighboursSearchStruct::cellCenter, cloudViewer::DgmOctree::NearestNeighboursSearchStruct::cellPos, ComputeMinDistanceToCellBorder(), GenerateTruncatedCellCode(), GET_BIT_SHIFT(), getCellIndex(), getCellSize(), getNeighborCellsAround(), cloudViewer::GenericIndexedCloudPersist::getPointPersistentPtr(), INVALID_CELL_CODE, cloudViewer::DgmOctree::NearestNeighboursSearchStruct::level, m_fillIndexes, m_numberOfProjectedPoints, m_theAssociatedCloud, m_thePointsAndTheirCellCodes, max(), cloudViewer::DgmOctree::NearestNeighboursSearchStruct::maxSearchSquareDistd, cloudViewer::DgmOctree::NearestNeighboursSearchStruct::minimalCellsSetToVisit, cloudViewer::DgmOctree::NearestNeighboursSearchStruct::queryPoint, cloudViewer::DgmOctree::NearestNeighboursSearchStruct::theNearestPointIndex, and Tuple3Tpl< Type >::u.
Referenced by cloudViewer::DistanceComputationTools::computeCellHausdorffDistance(), cloudViewer::DistanceComputationTools::computeCellHausdorffDistanceWithLocalModel(), and findPointNeighbourhood().
|
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).
| cellPos | the cell position |
| level | the level of subdivision |
Definition at line 120 of file DgmOctree.cpp.
References GET_BIT_SHIFT(), MAX_OCTREE_LEVEL, PRE_COMPUTED_POS_CODES, MonoDimensionalCellCodes::VALUE_COUNT, MonoDimensionalCellCodes::values, Tuple3Tpl< Type >::x, Tuple3Tpl< Type >::y, and Tuple3Tpl< Type >::z.
Referenced by findNearestNeighborsStartingFromCell(), findTheNearestNeighborStartingFromCell(), genericBuild(), getPointsInBoxNeighbourhood(), getPointsInCylindricalNeighbourhood(), getPointsInCylindricalNeighbourhoodProgressive(), and getPointsInSphericalNeighbourhood().
|
static |
Short version of GenerateTruncatedCellCode.
Definition at line 139 of file DgmOctree.cpp.
References GET_BIT_SHIFT(), MAX_OCTREE_LEVEL, PRE_COMPUTED_POS_CODES, MonoDimensionalCellCodes::VALUE_COUNT, MonoDimensionalCellCodes::values, Tuple3Tpl< Type >::x, Tuple3Tpl< Type >::y, and Tuple3Tpl< Type >::z.
|
protected |
Generic method to build the octree structure.
METHODS
| progressCb | the client application can get some notification of the process progress through this callback mechanism (see GenericProgressCallback) |
Definition at line 221 of file DgmOctree.cpp.
References cloudViewer::DgmOctree::IndexAndCode::codeComp(), GenerateTruncatedCellCode(), cloudViewer::GenericIndexedCloud::getPoint(), getTheCellPosWhichIncludesThePoint(), LOG_NAT_2, m_fillIndexes, m_nearestPow2, m_numberOfProjectedPoints, m_pointsMax, m_pointsMin, m_theAssociatedCloud, m_thePointsAndTheirCellCodes, MAX_OCTREE_LENGTH, MAX_OCTREE_LEVEL, cloudViewer::NormalizedProgress::oneStep(), ParallelSort, cloudViewer::GenericProgressCallback::setInfo(), cloudViewer::GenericProgressCallback::setMethodTitle(), cloudViewer::GenericCloud::size(), cloudViewer::GenericProgressCallback::start(), cloudViewer::GenericProgressCallback::stop(), cloudViewer::GenericProgressCallback::textCanBeEdited(), cloudViewer::GenericProgressCallback::update(), updateCellCountTable(), updateCellSizeTable(), Tuple3Tpl< Type >::x, Tuple3Tpl< Type >::y, and Tuple3Tpl< Type >::z.
Referenced by build().
|
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.
| level | the level of subdivision |
Definition at line 113 of file DgmOctree.cpp.
References PRE_COMPUTED_BIT_SHIFT_VALUES, and BitShiftValues::values.
Referenced by computeCellsStatistics(), ccComparisonDlg::determineBestOctreeLevel(), diff(), executeFunctionForAllCellsAtLevel(), executeFunctionForAllCellsStartingAtLevel(), extractCCs(), findNearestNeighborsStartingFromCell(), findTheNearestNeighborStartingFromCell(), GenerateTruncatedCellCode(), getCellCodes(), getCellCodesAndIndexes(), getCellIndexes(), getCellPos(), getNeighborCellsAround(), getPointsInBoxNeighbourhood(), getPointsInCell(), getPointsInCellByCellIndex(), getPointsInCellsWithSortedCellCodes(), getPointsInCylindricalNeighbourhood(), getPointsInCylindricalNeighbourhoodProgressive(), getPointsInNeighbourCellsAround(), getPointsInSphericalNeighbourhood(), and rayCast().
Returns the octree bounding box.
Method to request the octree bounding box limits
| bbMin | lower bounding-box limits (Xmin,Ymin,Zmin) |
| bbMax | higher bounding-box limits (Xmax,Ymax,Zmax) |
Definition at line 493 of file DgmOctree.cpp.
|
inline |
Returns the ith cell code.
Definition at line 963 of file DgmOctree.h.
Referenced by cloudViewer::DistanceComputationTools::computeApproxCloud2CloudDistance().
| 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.
| level | the level of subdivision |
| vec | the list of codes |
| truncatedCodes | indicates if the resulting codes should be truncated or not |
Definition at line 2822 of file DgmOctree.cpp.
References GET_BIT_SHIFT(), m_numberOfProjectedPoints, and m_thePointsAndTheirCellCodes.
Referenced by cloudViewer::DistanceComputationTools::computeApproxCloud2CloudDistance(), extractCCs(), FastMarchingForFacetExtraction::init(), and cloudViewer::FastMarchingForPropagation::init().
| 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.
| level | the level of subdivision |
| vec | the list of codes & indexes |
| truncatedCodes | indicates if the resulting codes should be truncated or not |
Definition at line 2794 of file DgmOctree.cpp.
References GET_BIT_SHIFT(), m_numberOfProjectedPoints, and m_thePointsAndTheirCellCodes.
Referenced by cloudViewer::DistanceComputationTools::computeCloud2MeshDistancesWithOctree(), and cloudViewer::DistanceComputationTools::computeCloud2MeshDistanceWithOctree().
| 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
| cellPos | cell position |
| level | level at which octree grid is considered |
| cellDists | output |
Definition at line 781 of file DgmOctree.cpp.
References m_fillIndexes, Tuple3Tpl< Type >::x, Tuple3Tpl< Type >::y, and Tuple3Tpl< Type >::z.
Referenced by getNeighborCellsAround(), and getPointsInNeighbourCellsAround().
| 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.
| cellPos | center cell position |
| level | level at which octree grid is considered |
| neighbourhoodLength | cell neighbourhood "radius" |
| cellDists | output |
Definition at line 795 of file DgmOctree.cpp.
References m_fillIndexes, and Tuple3Tpl< Type >::u.
|
protected |
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).
| truncatedCellCode | truncated cell code (i.e. original cell code shifted of 'bitDec' bits) |
| bitDec | binary shift corresponding to the level of subdivision (see GET_BIT_SHIFT) |
Definition at line 559 of file DgmOctree.cpp.
References m_nearestPow2, m_numberOfProjectedPoints, and m_thePointsAndTheirCellCodes.
Referenced by findNearestNeighborsStartingFromCell(), findTheNearestNeighborStartingFromCell(), getNeighborCellsAround(), getPointsInBoxNeighbourhood(), getPointsInCell(), getPointsInCylindricalNeighbourhood(), getPointsInCylindricalNeighbourhoodProgressive(), getPointsInNeighbourCellsAround(), and getPointsInSphericalNeighbourhood().
|
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.
| truncatedCellCode | truncated cell code (i.e. original cell code shifted of 'bitDec' bits) |
| bitDec | binary shift corresponding to the level of subdivision (see GET_BIT_SHIFT) |
| begin | first index of the sub-list in which to perform the binary search |
| end | last index of the sub-list in which to perform the binary search |
Definition at line 669 of file DgmOctree.cpp.
References count, INVALID_CELL_CODE, LOG_NAT_2, m_numberOfProjectedPoints, and m_thePointsAndTheirCellCodes.
| 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.
| level | the level of subdivision |
| vec | the list of indexes |
Definition at line 2850 of file DgmOctree.cpp.
References GET_BIT_SHIFT(), m_cellCount, m_numberOfProjectedPoints, and m_thePointsAndTheirCellCodes.
Referenced by cloudViewer::DistanceComputationTools::computeApproxCloud2CloudDistance().
|
inline |
Returns the number of cells for a given level of subdivision.
Definition at line 1038 of file DgmOctree.h.
Referenced by computeMeanOctreeDensity(), executeFunctionForAllCellsAtLevel(), executeFunctionForAllCellsStartingAtLevel(), ccPropertiesTreeDelegate::fillWithPointOctree(), findBestLevelForAGivenCellNumber(), cloudViewer::CloudSamplingTools::resampleCloudWithOctreeAtLevel(), and cloudViewer::CloudSamplingTools::subsampleCloudWithOctreeAtLevel().
| 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
| code | the cell code |
| level | the level of subdivision |
| cellPos | the computed position |
| isCodeTruncated | indicates if the given code is truncated or not |
Definition at line 498 of file DgmOctree.cpp.
References GET_BIT_SHIFT(), Tuple3Tpl< Type >::x, Tuple3Tpl< Type >::y, and Tuple3Tpl< Type >::z.
Referenced by cloudViewer::CloudSamplingTools::applyNoiseFilterAtLevel(), cloudViewer::CloudSamplingTools::applySORFilterAtLevel(), cloudViewer::DistanceComputationTools::computeApproxCloud2CloudDistance(), cloudViewer::GeometricalAnalysisTools::ComputeApproxPointsDensityInACellAtLevel(), cloudViewer::ScalarFieldTools::computeCellGaussianFilter(), cloudViewer::DistanceComputationTools::computeCellHausdorffDistance(), cloudViewer::DistanceComputationTools::computeCellHausdorffDistanceWithLocalModel(), computeCellLimits(), cloudViewer::DistanceComputationTools::computeCloud2MeshDistancesWithOctree(), cloudViewer::DistanceComputationTools::computeCloud2MeshDistanceWithOctree(), cloudViewer::GeometricalAnalysisTools::ComputeGeomCharacteristicAtLevel(), cloudViewer::StatisticalTestingTools::computeLocalChi2DistAtLevel(), cloudViewer::ScalarFieldTools::computeMeanGradientOnPatch(), extractCCs(), cloudViewer::GeometricalAnalysisTools::FlagDuplicatePointsInACellAtLevel(), FastMarchingForFacetExtraction::init(), cloudViewer::FastMarchingForPropagation::init(), and rayCast().
|
inline |
Returns the octree cells length for a given level of subdivision.
As the octree is cubical, cells are cubical.
| level | the level of subdivision (up to MAX_OCTREE_LEVEL+1 for convenience) |
Definition at line 494 of file DgmOctree.h.
Referenced by cloudViewer::CloudSamplingTools::applyNoiseFilterAtLevel(), cloudViewer::DistanceComputationTools::computeApproxCloud2CloudDistance(), cloudViewer::DistanceComputationTools::computeCellHausdorffDistanceWithLocalModel(), computeCellLimits(), cloudViewer::DistanceComputationTools::computeCloud2MeshDistancesWithOctree(), cloudViewer::DistanceComputationTools::computeCloud2MeshDistanceWithOctree(), cloudViewer::GeometricalAnalysisTools::ComputeGeomCharacteristicAtLevel(), cloudViewer::ScalarFieldTools::computeMeanGradientOnPatch(), ComputeNeighborhood2MeshDistancesWithOctree(), cloudViewer::DistanceComputationTools::computePoint2MeshDistancesWithOctree(), cloudViewer::ScalarFieldTools::computeScalarFieldGradient(), ccPropertiesTreeDelegate::fillWithPointOctree(), findBestLevelForAGivenNeighbourhoodSizeExtraction(), findNearestNeighborsStartingFromCell(), findNeighborsInASphereStartingFromCell(), findTheNearestNeighborStartingFromCell(), cloudViewer::GeometricalAnalysisTools::FlagDuplicatePointsInACellAtLevel(), getPointsInBoxNeighbourhood(), getPointsInCylindricalNeighbourhood(), getPointsInCylindricalNeighbourhoodProgressive(), getPointsInSphericalNeighbourhood(), cloudViewer::FastMarching::initGridWithOctree(), cloudViewer::DistanceComputationTools::intersectMeshWithOctree(), and rayCast().
|
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.
| level | the level of subdivision |
Definition at line 485 of file DgmOctree.h.
Referenced by cloudViewer::DistanceComputationTools::computeApproxCloud2CloudDistance(), getPointsInBoxNeighbourhood(), getPointsInCylindricalNeighbourhood(), getPointsInCylindricalNeighbourhoodProgressive(), and cloudViewer::FastMarching::initGridWithOctree().
|
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.
| level | the level of subdivision |
Definition at line 474 of file DgmOctree.h.
Referenced by cloudViewer::DistanceComputationTools::computeApproxCloud2CloudDistance(), getPointsInBoxNeighbourhood(), getPointsInCylindricalNeighbourhood(), getPointsInCylindricalNeighbourhoodProgressive(), and cloudViewer::FastMarching::initGridWithOctree().
|
protected |
Returns the indexes of the neighbourhing (existing) cells of a given cell This function is used by the nearest neighbours search algorithms.
| cellPos | the query cell |
| neighborCellsIndexes | the found neighbourhing cells |
| neighbourhoodLength | the distance (in terms of cells) at which to look for neighbour cells |
| level | the level of subdivision |
Definition at line 825 of file DgmOctree.cpp.
References abs(), GenerateCellCodeForDim(), GET_BIT_SHIFT(), getCellDistanceFromBorders(), getCellIndex(), m_numberOfProjectedPoints, Tuple3Tpl< Type >::x, Tuple3Tpl< Type >::y, and Tuple3Tpl< Type >::z.
Referenced by findTheNearestNeighborStartingFromCell().
|
inline |
Returns the number of points projected into the octree.
Definition at line 446 of file DgmOctree.h.
Referenced by findBestLevelForComparisonWithOctree(), and cloudViewer::DistanceComputationTools::synchronizeOctrees().
|
inline |
Returns the higher boundaries of the octree.
Definition at line 458 of file DgmOctree.h.
Referenced by cloudViewer::DistanceComputationTools::computeCloud2MeshDistances(), and cloudViewer::DistanceComputationTools::synchronizeOctrees().
|
inline |
Returns the lower boundaries of the octree.
Definition at line 453 of file DgmOctree.h.
Referenced by cloudViewer::DistanceComputationTools::computeCloud2MeshDistances(), cloudViewer::DistanceComputationTools::intersectMeshWithOctree(), and cloudViewer::DistanceComputationTools::synchronizeOctrees().
| std::size_t DgmOctree::getPointsInBoxNeighbourhood | ( | BoxNeighbourhood & | params | ) | const |
Returns the points falling inside a box.
Definition at line 1951 of file DgmOctree.cpp.
References abs(), GenerateTruncatedCellCode(), GET_BIT_SHIFT(), getCellIndex(), getCellSize(), getMaxFillIndexes(), getMinFillIndexes(), cloudViewer::GenericIndexedCloud::getPoint(), getTheCellPosWhichIncludesThePoint(), m_dimMin, m_numberOfProjectedPoints, m_theAssociatedCloud, m_thePointsAndTheirCellCodes, params, Tuple3Tpl< Type >::x, Tuple3Tpl< Type >::y, and Tuple3Tpl< Type >::z.
| 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.
| cellCode | the unique cell code | |
| level | the level of subdivision | |
| [out] | subset | set of points lying in the cell (references, no duplication) |
| isCodeTruncated | specifies if the code is given in a truncated form or not | |
| clearOutputCloud | whether to clear or not the output cloud (subest) if no points lie in the specified cell |
Definition at line 537 of file DgmOctree.cpp.
References cloudViewer::ReferenceCloud::clear(), GET_BIT_SHIFT(), getCellIndex(), getPointsInCellByCellIndex(), and m_numberOfProjectedPoints.
Referenced by FastMarchingForFacetExtraction::addCellToCurrentFacet(), FastMarchingForFacetExtraction::computeTCoefApprox(), extractCCs(), cloudViewer::FastMarchingForPropagation::extractPropagatedPoints(), FastMarchingForFacetExtraction::init(), cloudViewer::FastMarchingForPropagation::init(), cloudViewer::FastMarchingForPropagation::setPropagationTimingsAsDistances(), and FastMarchingForFacetExtraction::updateFlagsTable().
| 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 ;)
| cloud | ReferenceCloud to store the points lying inside the cell |
| cellIndex | the cell index |
| level | the level of subdivision |
| clearOutputCloud | whether to clear the input cloud prior to inserting the points or not |
Definition at line 2879 of file DgmOctree.cpp.
References cloudViewer::ReferenceCloud::addPointIndex(), cloudViewer::ReferenceCloud::clear(), GET_BIT_SHIFT(), cloudViewer::ReferenceCloud::getAssociatedCloud(), m_theAssociatedCloud, and m_thePointsAndTheirCellCodes.
Referenced by cloudViewer::DistanceComputationTools::computeApproxCloud2CloudDistance(), cloudViewer::DistanceComputationTools::computeCloud2MeshDistancesWithOctree(), cloudViewer::DistanceComputationTools::computeCloud2MeshDistanceWithOctree(), and getPointsInCell().
| 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.
| cellCodes | the cells codes | |
| level | the level of subdivision | |
| [out] | subset | set of points lying in the cell (references, no duplication) |
| areCodesTruncated | specifies if the codes are given in a truncated form or not |
Definition at line 2909 of file DgmOctree.cpp.
References cloudViewer::ReferenceCloud::addPointIndex(), cloudViewer::ReferenceCloud::clear(), GET_BIT_SHIFT(), m_numberOfProjectedPoints, and m_thePointsAndTheirCellCodes.
| 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 ;).
| params | input/output parameters structure |
Definition at line 2111 of file DgmOctree.cpp.
References Vector3Tpl< Type >::dot(), dot(), GenerateTruncatedCellCode(), GET_BIT_SHIFT(), getCellIndex(), getCellSize(), getMaxFillIndexes(), getMinFillIndexes(), cloudViewer::GenericIndexedCloud::getPoint(), getTheCellPosWhichIncludesThePoint(), m_dimMin, m_numberOfProjectedPoints, m_theAssociatedCloud, m_thePointsAndTheirCellCodes, max(), min(), params, SQRT_3, Tuple3Tpl< Type >::x, Tuple3Tpl< Type >::y, and Tuple3Tpl< Type >::z.
| 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).
Definition at line 2264 of file DgmOctree.cpp.
References Vector3Tpl< Type >::dot(), dot(), GenerateTruncatedCellCode(), GET_BIT_SHIFT(), getCellIndex(), getCellSize(), getMaxFillIndexes(), getMinFillIndexes(), cloudViewer::GenericIndexedCloud::getPoint(), getTheCellPosWhichIncludesThePoint(), m_dimMin, m_numberOfProjectedPoints, m_theAssociatedCloud, m_thePointsAndTheirCellCodes, max(), min(), params, SQRT_3, std::swap(), Tuple3Tpl< Type >::x, Tuple3Tpl< Type >::y, and Tuple3Tpl< Type >::z.
|
protected |
Gets point in the neighbourhing cells of a specific cell.
| nNSS | NN search parameters (from which are used: cellPos, pointsInNeighbourCells and level) |
| neighbourhoodLength | the new distance (in terms of cells) at which to look for neighbour cells |
| getOnlyPointsWithValidScalar | whether to ignore points having an invalid associated scalar value |
Definition at line 906 of file DgmOctree.cpp.
References abs(), cloudViewer::DgmOctree::NearestNeighboursSearchStruct::alreadyVisitedNeighbourhoodSize, cloudViewer::utility::ceil(), cloudViewer::DgmOctree::NearestNeighboursSearchStruct::cellPos, GenerateCellCodeForDim(), GET_BIT_SHIFT(), getCellDistanceFromBorders(), getCellIndex(), cloudViewer::DgmOctree::NearestNeighboursSearchStruct::level, m_averageCellPopulation, m_numberOfProjectedPoints, m_theAssociatedCloud, m_thePointsAndTheirCellCodes, cloudViewer::DgmOctree::NearestNeighboursSearchStruct::pointsInNeighbourhood, cloudViewer::ScalarField::ValidValue(), Tuple3Tpl< Type >::x, Tuple3Tpl< Type >::y, and Tuple3Tpl< Type >::z.
Referenced by findNearestNeighborsStartingFromCell(), and findNeighborsInASphereStartingFromCell().
| 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 ;).
| sphereCenter | center | |
| radius | radius | |
| [out] | neighbours | points falling inside the sphere |
| level | subdivision level at which to apply the extraction process |
Definition at line 1846 of file DgmOctree.cpp.
References GenerateTruncatedCellCode(), GET_BIT_SHIFT(), getCellIndex(), getCellSize(), cloudViewer::GenericIndexedCloud::getPoint(), getTheCellPosWhichIncludesThePoint(), m_dimMin, m_numberOfProjectedPoints, m_theAssociatedCloud, m_thePointsAndTheirCellCodes, OCTREE_LENGTH(), SQRT_3, Tuple3Tpl< Type >::x, Tuple3Tpl< Type >::y, and Tuple3Tpl< Type >::z.
Referenced by RefinePointClassif(), cloudViewer::CloudSamplingTools::resampleCloudSpatially(), and qCanupoTools::TrainClassifier().
|
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.
| thePoint | the query point |
| cellPos | the computed position |
Definition at line 779 of file DgmOctree.h.
References Tuple3Tpl< Type >::x, Tuple3Tpl< Type >::y, and Tuple3Tpl< Type >::z.
Referenced by cloudViewer::DistanceComputationTools::computeCellHausdorffDistanceWithLocalModel(), cloudViewer::DistanceComputationTools::computeGeodesicDistances(), cloudViewer::DistanceComputationTools::computePoint2MeshDistancesWithOctree(), FastMarchingForFacetExtraction::ExtractPlanarFacets(), findPointNeighbourhood(), cloudViewer::AutoSegmentationTools::frontPropagationBasedSegmentation(), genericBuild(), getPointsInBoxNeighbourhood(), getPointsInCylindricalNeighbourhood(), getPointsInCylindricalNeighbourhoodProgressive(), getPointsInSphericalNeighbourhood(), cloudViewer::DistanceComputationTools::intersectMeshWithOctree(), and masc::Tools::PrepareFeatures().
|
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.
| thePoint | the query point |
| cellPos | the computed position |
| level | the level of subdivision |
Definition at line 798 of file DgmOctree.h.
References Tuple3Tpl< Type >::x, Tuple3Tpl< Type >::y, and Tuple3Tpl< Type >::z.
|
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.
| thePoint | the query point |
| cellPos | the computed position |
| level | the level of subdivision |
| inBounds | indicates if the query point is inside or outside the octree bounding-box |
Definition at line 823 of file DgmOctree.h.
References Tuple3Tpl< Type >::x, Tuple3Tpl< Type >::y, and Tuple3Tpl< Type >::z.
|
static |
Returns whether multi-threading (parallel) computation is supported or not
Definition at line 162 of file DgmOctree.cpp.
|
static |
Returns the octree length (in term of cells) for a given level of subdivision
| level | the level of subdivision |
Definition at line 118 of file DgmOctree.cpp.
Referenced by cloudViewer::AutoSegmentationTools::frontPropagationBasedSegmentation(), and getPointsInSphericalNeighbourhood().
|
inline |
Returns the octree 'structure'.
Definition at line 1196 of file DgmOctree.h.
| bool DgmOctree::rayCast | ( | const CCVector3 & | rayAxis, |
| const CCVector3 & | rayOrigin, | ||
| double | maxRadiusOrFov, | ||
| bool | isFOV, | ||
| RayCastProcess | process, | ||
| std::vector< PointDescriptor > & | output | ||
| ) | const |
Ray casting algorithm.
Definition at line 4353 of file DgmOctree.cpp.
References cloudViewer::GenericCloud::enableScalarField(), findBestLevelForAGivenPopulationPerCell(), GET_BIT_SHIFT(), getCellPos(), getCellSize(), cloudViewer::GenericIndexedCloud::getPoint(), AABB< T >::intersects(), INVALID_CELL_CODE, m_dimMax, m_dimMin, m_theAssociatedCloud, m_thePointsAndTheirCellCodes, Ray< T >::radialSquareDistance(), RC_CLOSE_POINTS, RC_NEAREST_POINT, cloudViewer::GenericCloud::setPointScalarValue(), SQRT_3, Ray< T >::squareDistances(), Ray< T >::squareDistanceToOrigin(), Tuple3Tpl< Type >::x, Tuple3Tpl< Type >::y, and Tuple3Tpl< Type >::z.
|
protected |
Updates the tables containing the number of octree cells for each level of subdivision
Definition at line 416 of file DgmOctree.cpp.
References computeCellsStatistics(), and MAX_OCTREE_LEVEL.
Referenced by clear(), and genericBuild().
|
protected |
Updates the tables containing the octree cells length for each level of subdivision
Definition at line 405 of file DgmOctree.cpp.
References m_cellSize, m_dimMax, m_dimMin, MAX_OCTREE_LEVEL, and Tuple3Tpl< Type >::x.
Referenced by genericBuild().
|
protected |
Updates the tables containing octree limits and boundaries.
Definition at line 395 of file DgmOctree.cpp.
References cloudViewer::GenericCloud::getBoundingBox(), m_dimMax, m_dimMin, m_pointsMax, m_pointsMin, m_theAssociatedCloud, and cloudViewer::CCMiscTools::MakeMinAndMaxCubical().
Referenced by build().
Invalid cell code.
Definition at line 89 of file DgmOctree.h.
Referenced by findNearestNeighborsStartingFromCell(), findTheNearestNeighborStartingFromCell(), getCellIndex(), and rayCast().
|
protected |
Average cell population per level of subdivision.
Definition at line 1264 of file DgmOctree.h.
Referenced by computeCellsStatistics(), executeFunctionForAllCellsAtLevel(), executeFunctionForAllCellsStartingAtLevel(), findBestLevelForAGivenNeighbourhoodSizeExtraction(), findBestLevelForAGivenPopulationPerCell(), and getPointsInNeighbourCellsAround().
|
protected |
Number of cells per level of subdivision.
Definition at line 1260 of file DgmOctree.h.
Referenced by computeCellsStatistics(), and getCellIndexes().
|
protected |
Cell dimensions for all subdivision levels.
Definition at line 1255 of file DgmOctree.h.
Referenced by clear(), and updateCellSizeTable().
|
protected |
Max coordinates of the octree bounding-box.
Definition at line 1245 of file DgmOctree.h.
Referenced by build(), clear(), getBoundingBox(), rayCast(), updateCellSizeTable(), and updateMinAndMaxTables().
|
protected |
Min coordinates of the octree bounding-box.
Definition at line 1243 of file DgmOctree.h.
Referenced by build(), clear(), computeCellLimits(), getBoundingBox(), getPointsInBoxNeighbourhood(), getPointsInCylindricalNeighbourhood(), getPointsInCylindricalNeighbourhoodProgressive(), getPointsInSphericalNeighbourhood(), rayCast(), updateCellSizeTable(), and updateMinAndMaxTables().
|
protected |
Min and max occupied cells indexes, for all dimensions and every subdivision level
Definition at line 1258 of file DgmOctree.h.
Referenced by clear(), findNearestNeighborsStartingFromCell(), findTheNearestNeighborStartingFromCell(), genericBuild(), and getCellDistanceFromBorders().
|
protected |
Max cell population per level of subdivision.
Definition at line 1262 of file DgmOctree.h.
Referenced by computeCellsStatistics(), executeFunctionForAllCellsAtLevel(), and executeFunctionForAllCellsStartingAtLevel().
|
protected |
Nearest power of 2 less than the number of points (used for binary search)
Definition at line 1240 of file DgmOctree.h.
Referenced by clear(), genericBuild(), and getCellIndex().
|
protected |
Number of points projected in the octree.
Definition at line 1236 of file DgmOctree.h.
Referenced by clear(), computeMeanOctreeDensity(), executeFunctionForAllCellsStartingAtLevel(), findNearestNeighborsStartingFromCell(), findTheNearestNeighborStartingFromCell(), genericBuild(), getCellCodes(), getCellCodesAndIndexes(), getCellIndex(), getCellIndexes(), getNeighborCellsAround(), getPointsInBoxNeighbourhood(), getPointsInCell(), getPointsInCellsWithSortedCellCodes(), getPointsInCylindricalNeighbourhood(), getPointsInCylindricalNeighbourhoodProgressive(), getPointsInNeighbourCellsAround(), and getPointsInSphericalNeighbourhood().
|
protected |
Max coordinates of the bounding-box of the set of points projected in the octree
Definition at line 1252 of file DgmOctree.h.
Referenced by build(), clear(), genericBuild(), and updateMinAndMaxTables().
|
protected |
Min coordinates of the bounding-box of the set of points projected in the octree
Definition at line 1249 of file DgmOctree.h.
Referenced by build(), clear(), genericBuild(), and updateMinAndMaxTables().
|
protected |
Std. dev. of cell population per level of subdivision.
Definition at line 1266 of file DgmOctree.h.
Referenced by computeCellsStatistics(), executeFunctionForAllCellsAtLevel(), and executeFunctionForAllCellsStartingAtLevel().
|
protected |
Associated cloud.
Definition at line 1233 of file DgmOctree.h.
Referenced by DgmOctree(), executeFunctionForAllCellsAtLevel(), executeFunctionForAllCellsStartingAtLevel(), extractCCs(), findNearestNeighborsStartingFromCell(), findTheNearestNeighborStartingFromCell(), genericBuild(), getPointsInBoxNeighbourhood(), getPointsInCellByCellIndex(), getPointsInCylindricalNeighbourhood(), getPointsInCylindricalNeighbourhoodProgressive(), getPointsInNeighbourCellsAround(), getPointsInSphericalNeighbourhood(), cloudViewer::DgmOctree::octreeCell::octreeCell(), rayCast(), and updateMinAndMaxTables().
|
protected |
The coded octree structure.
ATTRIBUTES
Definition at line 1230 of file DgmOctree.h.
Referenced by build(), clear(), computeCellsStatistics(), executeFunctionForAllCellsAtLevel(), executeFunctionForAllCellsStartingAtLevel(), findBestLevelForComparisonWithOctree(), findNearestNeighborsStartingFromCell(), findTheNearestNeighborStartingFromCell(), genericBuild(), getCellCodes(), getCellCodesAndIndexes(), getCellIndex(), getCellIndexes(), getPointsInBoxNeighbourhood(), getPointsInCellByCellIndex(), getPointsInCellsWithSortedCellCodes(), getPointsInCylindricalNeighbourhood(), getPointsInCylindricalNeighbourhoodProgressive(), getPointsInNeighbourCellsAround(), getPointsInSphericalNeighbourhood(), and rayCast().
|
static |
Max octree length at last level of subdivision (number of cells)
Definition at line 84 of file DgmOctree.h.
Referenced by genericBuild().
|
static |
Max octree subdivision level.
STRUCTURES
Number of bits used to code the cells position: 3*MAX_OCTREE_LEVEL
Definition at line 67 of file DgmOctree.h.
Referenced by FastMarchingForFacetExtraction::addCellToCurrentFacet(), BitShiftValues::BitShiftValues(), ccComparisonDlg::ccComparisonDlg(), ccComputeOctreeDlg::ccComputeOctreeDlg(), ccLabelingDlg::ccLabelingDlg(), ccAlignDlg::changeSamplingMethod(), ccSubsamplingDlg::changeSamplingMethod(), clear(), cloudViewer::DistanceComputationTools::computeApproxCloud2CloudDistance(), computeCellsStatistics(), ccComparisonDlg::computeDistances(), ComputeMathOpWithNearestNeighbor(), ccEntityAction::computeOctree(), ccPropertiesTreeDelegate::createEditor(), ccComparisonDlg::determineBestOctreeLevel(), cloudViewer::DistanceComputationTools::diff(), executeFunctionForAllCellsStartingAtLevel(), cloudViewer::FastMarchingForPropagation::extractPropagatedPoints(), ccPropertiesTreeDelegate::fillWithPointOctree(), findBestLevelForAGivenCellNumber(), findBestLevelForAGivenNeighbourhoodSizeExtraction(), findBestLevelForAGivenPopulationPerCell(), findBestLevelForComparisonWithOctree(), GenerateTruncatedCellCode(), genericBuild(), qM3C2Tools::GuessBestParams(), cloudViewer::FastMarching::initGridWithOctree(), MonoDimensionalCellCodes::MonoDimensionalCellCodes(), ccEntityAction::orientNormalsFM(), masc::ContextBasedFeature::prepare(), CommandSubsample::process(), CommandExtractCCs::process(), cloudViewer::FastMarchingForPropagation::setPropagationTimingsAsDistances(), updateCellCountTable(), and updateCellSizeTable().