ITK  4.6.0
Insight Segmentation and Registration Toolkit
itkPointsLocator.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * http://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 #ifndef __itkPointsLocator_h
19 #define __itkPointsLocator_h
20 
21 #include "itkObject.h"
22 
23 #include "itkPoint.h"
24 #include "itkIntTypes.h"
25 #include "itkKdTree.h"
26 #include "itkKdTreeGenerator.h"
27 #include "itkVectorContainer.h"
29 
30 namespace itk
31 {
32 
41 template<
42  typename TPointsContainer = VectorContainer< IdentifierType, Point< float, 3 > >
43  >
44 class PointsLocator : public Object
45 {
46 public:
49  typedef Object Superclass;
52 
54  itkNewMacro( Self );
55 
57  itkTypeMacro( PointsLocator, Object );
58 
60  typedef TPointsContainer PointsContainer;
61  typedef typename PointsContainer::Pointer PointsContainerPointer;
62  typedef typename PointsContainer::ConstPointer PointsContainerConstPointer;
63  typedef typename PointsContainer::ElementIdentifier PointIdentifier;
64  typedef typename PointsContainer::Element PointType;
65 
67  itkStaticConstMacro( PointDimension, unsigned int, PointType::PointDimension );
68 
70  typedef typename PointsContainer::ConstIterator PointsContainerConstIterator;
71  typedef typename PointsContainer::Iterator PointsContainerIterator;
72 
77 
84 
86  itkSetObjectMacro( Points, PointsContainer );
87 
89  itkGetModifiableObjectMacro(Points, PointsContainer );
90 
92  void Initialize();
93 
95  PointIdentifier FindClosestPoint( const PointType &query ) const;
96 
98  void Search( const PointType &, unsigned int, NeighborsIdentifierType & )
99  const;
100 
102  void FindClosestNPoints( const PointType &, unsigned int,
103  NeighborsIdentifierType & ) const;
104 
106  void Search( const PointType &, double, NeighborsIdentifierType & ) const;
107 
109  void FindPointsWithinRadius( const PointType &, double,
110  NeighborsIdentifierType & ) const;
111 
112 protected:
113  PointsLocator();
114  ~PointsLocator();
115  virtual void PrintSelf(std::ostream& os, Indent indent) const ITK_OVERRIDE;
116 
117 private:
118  PointsLocator( const Self& ); //purposely not implemented
119  void operator=( const Self& ); //purposely not implemented
120 
125 };
126 
127 } // end namespace itk
128 
129 #ifndef ITK_MANUAL_INSTANTIATION
130 #include "itkPointsLocator.hxx"
131 #endif
132 
133 #endif
SmartPointer< Self > Pointer
SmartPointer< const Self > ConstPointer
PointIdentifier FindClosestPoint(const PointType &query) const
Light weight base class for most itk classes.
PointsContainer::Iterator PointsContainerIterator
This class provides ListSample interface to ITK VectorContainer.
PointsContainer::Pointer PointsContainerPointer
SampleAdaptorPointer m_SampleAdaptor
PointsContainer::ConstIterator PointsContainerConstIterator
TreeConstPointer m_Tree
TreeGeneratorType::KdTreeType TreeType
void operator=(const Self &)
PointsContainer::ConstPointer PointsContainerConstPointer
PointsLocator Self
Statistics::KdTreeGenerator< SampleAdaptorType > TreeGeneratorType
Accelerate geometric searches for points.
void FindPointsWithinRadius(const PointType &, double, NeighborsIdentifierType &) const
TreeGeneratorPointer m_KdTreeGenerator
void Search(const PointType &, unsigned int, NeighborsIdentifierType &) const
PointsContainer::Element PointType
PointsContainer::ElementIdentifier PointIdentifier
This class generates a KdTree object without centroid information.
TreeGeneratorType::Pointer TreeGeneratorPointer
TreeType::ConstPointer TreeConstPointer
std::vector< InstanceIdentifier > InstanceIdentifierVectorType
Definition: itkKdTree.h:522
TreeType::InstanceIdentifierVectorType NeighborsIdentifierType
TPointsContainer PointsContainer
void FindClosestNPoints(const PointType &, unsigned int, NeighborsIdentifierType &) const
static const unsigned int PointDimension
Control indentation during Print() invocation.
Definition: itkIndent.h:49
PointsContainerPointer m_Points
Base class for most ITK classes.
Definition: itkObject.h:57
SampleAdaptorType::Pointer SampleAdaptorPointer
Statistics::VectorContainerToListSampleAdaptor< PointsContainer > SampleAdaptorType
This class provides methods for k-nearest neighbor search and related data structures for a k-d tree...
Definition: itkKdTree.h:483
virtual void PrintSelf(std::ostream &os, Indent indent) const ITK_OVERRIDE