OTB  10.0.0
Orfeo Toolbox
otbHaralickTexturesImageFunction.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 otbHaralickTexturesImageFunction_h
22 #define otbHaralickTexturesImageFunction_h
23 
25 #include "itkImageToImageFilter.h"
26 #include "itkImageFunction.h"
27 #include "itkFixedArray.h"
28 
29 namespace otb
30 {
92 
93 
94 template <class TInputImage, class TCoordRep = double>
96  : public itk::ImageFunction<TInputImage, itk::FixedArray<typename itk::NumericTraits<typename TInputImage::PixelType>::RealType, 8>, TCoordRep>
97 {
98 public:
101  typedef itk::ImageFunction<TInputImage, itk::FixedArray<typename itk::NumericTraits<typename TInputImage::PixelType>::RealType, 8>, TCoordRep> Superclass;
102  typedef itk::SmartPointer<Self> Pointer;
103  typedef itk::SmartPointer<const Self> ConstPointer;
104 
106  itkTypeMacro(HaralickTexturesImageFunction, ImageFunction);
107 
109  itkNewMacro(Self);
110 
112  typedef TInputImage InputImageType;
113  typedef typename Superclass::IndexType IndexType;
114  typedef typename Superclass::ContinuousIndexType ContinuousIndexType;
115  typedef typename Superclass::PointType PointType;
116  typedef typename InputImageType::Pointer InputImagePointerType;
117  typedef typename InputImageType::PixelType InputPixelType;
118  typedef typename InputImageType::RegionType InputRegionType;
119  typedef typename InputImageType::OffsetType OffsetType;
120  typedef typename InputRegionType::SizeType SizeType;
121 
122  // Output typedef support
123  typedef typename Superclass::OutputType OutputType;
124  typedef typename OutputType::ValueType ScalarRealType;
125  // Textures tools
133 
134  typedef typename VectorType::iterator VectorIteratorType;
135  typedef typename VectorType::const_iterator VectorConstIteratorType;
136  // Coord rep
137  typedef TCoordRep CoordRepType;
138 
140  itkStaticConstMacro(ImageDimension, unsigned int, InputImageType::ImageDimension);
141 
143  OutputType EvaluateAtIndex(const IndexType& index) const override;
144 
146  OutputType Evaluate(const PointType& point) const override
147  {
148  IndexType index;
149  this->ConvertPointToNearestIndex(point, index);
150  return this->EvaluateAtIndex(index);
151  }
153  {
154  IndexType index;
155  this->ConvertContinuousIndexToNearestIndex(cindex, index);
156  return this->EvaluateAtIndex(index);
157  }
159 
164  itkSetMacro(NeighborhoodRadius, unsigned int);
165 
170  itkGetConstReferenceMacro(NeighborhoodRadius, unsigned int);
171 
173  itkSetMacro(Offset, OffsetType);
174 
176  itkGetMacro(Offset, OffsetType);
177 
179  itkSetMacro(NumberOfBinsPerAxis, unsigned int);
180 
182  itkGetMacro(NumberOfBinsPerAxis, unsigned int);
183 
185  itkSetMacro(InputImageMinimum, InputPixelType);
186 
188  itkGetMacro(InputImageMinimum, InputPixelType);
189 
191  itkSetMacro(InputImageMaximum, InputPixelType);
192 
194  itkGetMacro(InputImageMaximum, InputPixelType);
195 
196 protected:
199  {
200  }
201  void PrintSelf(std::ostream& os, itk::Indent indent) const override;
202 
203 private:
205  void operator=(const Self&) = delete;
206 
208  unsigned int m_NeighborhoodRadius;
209 
212 
214  unsigned int m_NumberOfBinsPerAxis;
215 
218 
221 
222  // TODO: should we use constexpr? only c++11 and problem for msvc
223  inline double GetPixelValueTolerance() const
224  {
225  return 0.0001;
226  }
227 };
228 
229 } // namespace otb
230 
231 #ifndef OTB_MANUAL_INSTANTIATION
233 #endif
234 
235 #endif
This class holds a VectorType of CooccurrencePairType with each pair is a combination of pixel index ...
itk::NumericTraits< FrequencyType >::RealType RelativeFrequencyType
itk::NumericTraits< PixelType >::RealType PixelValueType
Compute the 8 Haralick texture indices on the neighborhood of the point.
HaralickTexturesImageFunction(const Self &)=delete
CooccurrenceIndexedListType::RelativeFrequencyType RelativeFrequencyType
OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &cindex) const override
CooccurrenceIndexedListType::IndexType CooccurrenceIndexType
OutputType Evaluate(const PointType &point) const override
GreyLevelCooccurrenceIndexedList< InputPixelType > CooccurrenceIndexedListType
void operator=(const Self &)=delete
CooccurrenceIndexedListType::Pointer CooccurrenceIndexedListPointerType
CooccurrenceIndexedListType::ConstPointer CooccurrenceIndexedListConstPointerType
itk::ImageFunction< TInputImage, itk::FixedArray< typename itk::NumericTraits< typename TInputImage::PixelType >::RealType, 8 >, TCoordRep > Superclass
CooccurrenceIndexedListType::VectorType VectorType
Superclass::ContinuousIndexType ContinuousIndexType
CooccurrenceIndexedListType::PixelValueType PixelValueType
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.