OTB  10.0.0
Orfeo Toolbox
otbImageToSURFKeyPointSetFilter.h
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2005-2024 Centre National d'Etudes Spatiales (CNES)
3  *
4  * This file is part of Orfeo Toolbox
5  *
6  * https://www.orfeo-toolbox.org/
7  *
8  * Licensed under the Apache License, Version 2.0 (the "License");
9  * you may not use this file except in compliance with the License.
10  * You may obtain a copy of the License at
11  *
12  * http://www.apache.org/licenses/LICENSE-2.0
13  *
14  * Unless required by applicable law or agreed to in writing, software
15  * distributed under the License is distributed on an "AS IS" BASIS,
16  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
17  * See the License for the specific language governing permissions and
18  * limitations under the License.
19  */
20 
21 #ifndef otbImageToSURFKeyPointSetFilter_h
22 #define otbImageToSURFKeyPointSetFilter_h
23 
24 #include "itkConstNeighborhoodIterator.h"
26 #include "itkPointSet.h"
27 #include <itkRescaleIntensityImageFilter.h>
28 #include <otbImageList.h>
29 #include "itkResampleImageFilter.h"
31 #include "itkVector.h"
32 
33 namespace otb
34 {
35 
61 template <class TInputImage, class TOutputPointSet>
62 class ITK_EXPORT ImageToSURFKeyPointSetFilter : public ImageToPointSetFilter<TInputImage, TOutputPointSet>
63 {
64 
65 public:
69  typedef itk::SmartPointer<Self> Pointer;
70  typedef itk::SmartPointer<const Self> ConstPointer;
71 
73  itkNewMacro(Self);
74 
77 
79  typedef TInputImage InputImageType;
80  typedef typename InputImageType::Pointer InputImagePointerType;
81  typedef typename InputImageType::IndexType PixelIndex;
82  typedef typename InputImageType::IndexType IndexType;
83  typedef typename InputImageType::PixelType PixelValue;
84  typedef typename InputImageType::SpacingType SpacingType;
85  typedef typename InputImageType::SizeType SizeType;
86  typedef typename InputImageType::PointType PointImageType;
87 
88  typedef TOutputPointSet OutputPointSetType;
89  typedef typename TOutputPointSet::Pointer OutputPointSetPointerType;
90  typedef typename TOutputPointSet::PixelType OutputPixelType;
91  typedef typename TOutputPointSet::PointType OutputPointType;
92  typedef typename TOutputPointSet::PointIdentifier OutputPointIdentifierType;
93 
95  itkSetMacro(OctavesNumber, int);
96  itkGetMacro(OctavesNumber, int);
98 
100  itkSetMacro(ScalesNumber, int);
101  itkGetMacro(ScalesNumber, int);
103 
105  itkGetMacro(NumberOfPoints, int);
106 
108  typedef itk::ConstNeighborhoodIterator<InputImageType> NeighborhoodIteratorType;
109  typedef typename NeighborhoodIteratorType::NeighborhoodType NeighborhoodType;
110  typedef typename NeighborhoodType::OffsetType OffsetType;
111 
114 
116  typedef itk::ResampleImageFilter<InputImageType, InputImageType> ResampleFilterType;
117  typedef typename ResampleFilterType::Pointer ResampleFilterPointerType;
118 
122 
124  typedef std::vector<double> VectorType;
125  typedef itk::Vector<PixelValue, 3> VectorPointType;
126 
127 protected:
132 
137 
141  void PrintSelf(std::ostream& os, itk::Indent indent) const override;
142 
146  void GenerateData() override;
147 
154  virtual bool IsLocalExtremum(const NeighborhoodType& neigh);
155 
163  virtual bool IsLocalExtremumAround(const NeighborhoodType& neigh, double CenterValue);
164 
177  bool RefineLocationKeyPoint(const NeighborhoodIteratorType& currentScale, const NeighborhoodIteratorType& previousScale,
178  const NeighborhoodIteratorType& nextScale, VectorPointType& solution);
179 
187  virtual double AssignOrientation(const NeighborhoodType& neigh, double S);
188 
197  virtual VectorType ComputeDescriptor(const NeighborhoodType& neigh, double O, double S);
198 
202  virtual int GetMin(int a, int b, int c);
203 
204 private:
206  void operator=(const Self&) = delete;
207 
210 
213 
216 
219 
225 
228 
229  /*Resample Filter*/
231 
234 
237 
239  OffsetType m_Offsets[8];
240 };
241 }
242 #ifndef OTB_MANUAL_INSTANTIATION
244 #endif
245 
246 #endif
This class represent a list of images.
Definition: otbImageList.h:40
itk::SmartPointer< Self > Pointer
Definition: otbImageList.h:45
This class compute the Hessian determinant of each pixel of an input image.
Base class to output PointSet data with image data as input.
This class extracts key points from an image through a pyramidal gaussian based decomposition.
TOutputPointSet::PointIdentifier OutputPointIdentifierType
void PrintSelf(std::ostream &os, itk::Indent indent) const override
ImageToPointSetFilter< TInputImage, TOutputPointSet > Superclass
ResampleFilterType::Pointer ResampleFilterPointerType
virtual VectorType ComputeDescriptor(const NeighborhoodType &neigh, double O, double S)
bool RefineLocationKeyPoint(const NeighborhoodIteratorType &currentScale, const NeighborhoodIteratorType &previousScale, const NeighborhoodIteratorType &nextScale, VectorPointType &solution)
virtual bool IsLocalExtremumAround(const NeighborhoodType &neigh, double CenterValue)
itk::ConstNeighborhoodIterator< InputImageType > NeighborhoodIteratorType
virtual bool IsLocalExtremum(const NeighborhoodType &neigh)
otb::ImageList< InputImageType > ImageListType
otb::ImageToHessianDeterminantImageFilter< InputImageType, InputImageType > ImageToDetHessianImageType
ImageToDetHessianImageType::Pointer DetHessianPointerFilter
void operator=(const Self &)=delete
itk::ResampleImageFilter< InputImageType, InputImageType > ResampleFilterType
ImageToSURFKeyPointSetFilter(const Self &)=delete
virtual int GetMin(int a, int b, int c)
NeighborhoodIteratorType::NeighborhoodType NeighborhoodType
virtual double AssignOrientation(const NeighborhoodType &neigh, double S)
itk::SmartPointer< const Self > ConstPointer
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.