OTB  10.0.0
Orfeo Toolbox
otbSFSTexturesImageFilter.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 otbSFSTexturesImageFilter_h
22 #define otbSFSTexturesImageFilter_h
23 
24 #include "otbSFSTexturesFunctor.h"
25 #include "itkImageToImageFilter.h"
26 #include "itkImageRegionIteratorWithIndex.h"
27 #include "itkConstNeighborhoodIterator.h"
28 
29 namespace otb
30 {
31 
57 template <class TInputImage, class TOutputImage>
58 class ITK_EXPORT SFSTexturesImageFilter : public itk::ImageToImageFilter<TInputImage, TOutputImage>
59 {
60 public:
63  typedef TInputImage InputImageType;
64  typedef TOutputImage OutputImageType;
65  typedef itk::ImageToImageFilter<TInputImage, TOutputImage> Superclass;
66  typedef itk::SmartPointer<Self> Pointer;
67  typedef itk::SmartPointer<const Self> ConstPointer;
68 
70  itkNewMacro(Self);
71 
73  itkTypeMacro(UnaryFunctorNeighborhoodImageFilter, ImageToImageFilter);
74 
76  typedef typename InputImageType::ConstPointer InputImagePointerType;
77  typedef typename InputImageType::RegionType InputImageRegionType;
78  typedef typename InputImageType::PixelType InputImagePixelType;
79  typedef typename InputImageType::SizeType InputImageSizeType;
80  typedef typename OutputImageType::Pointer OutputImagePointerType;
81  typedef typename OutputImageType::RegionType OutputImageRegionType;
82  typedef typename OutputImageType::PixelType OutputImagePixelType;
83  typedef itk::ConstNeighborhoodIterator<TInputImage> NeighborhoodIteratorType;
84  typedef typename NeighborhoodIteratorType::RadiusType RadiusType;
87  typedef itk::ProcessObject ProcessObjectType;
88 
90  itkGetMacro(Radius, unsigned int);
91 
94  {
95  return m_Functor;
96  }
97 
98  const FunctorType& GetFunctor() const
99  {
100  return m_Functor;
101  }
102  void SetFunctor(const FunctorType& functor)
103  {
104  m_Functor = functor;
105  this->Modified();
106  }
107 
109  void SetSpatialThreshold(unsigned int thresh)
110  {
111  this->GetFunctor().SetSpatialThreshold(thresh);
112  m_Radius = thresh;
113  this->Modified();
114  }
115  unsigned int GetSpatialThreshold()
116  {
117  return this->GetFunctor().GetSpatialThreshold();
118  }
119 
122  {
123  this->GetFunctor().SetSpectralThreshold(thresh);
124  }
126  {
127  return this->GetFunctor().GetSpectralThreshold();
128  }
129 
131  void SetRatioMaxConsiderationNumber(unsigned int value)
132  {
133  this->GetFunctor().SetRatioMaxConsiderationNumber(value);
134  }
136  {
137  return this->GetFunctor().GetRatioMaxConsiderationNumber();
138  }
140 
142  void SetAlpha(double alpha)
143  {
144  this->GetFunctor().SetAlpha(alpha);
145  }
146  double GetAlpha()
147  {
148  return this->GetFunctor().GetAlpha();
149  }
151 
153  void SetNumberOfDirections(unsigned int D)
154  {
155  this->GetFunctor().SetNumberOfDirections(D);
156  double step = CONST_PI / static_cast<double>(D);
157  this->GetFunctor().SetDirectionStep(step);
158  }
159  unsigned int GetNumberOfDirections()
160  {
161  return this->GetFunctor().GetNumberOfDirections();
162  }
164 
174  typedef enum { LENGTH = 1, WIDTH, PSI, WMEAN, RATIO, SD } FeatureType;
175 
176  void SetFeatureStatus(FeatureType id, bool isSelected)
177  {
178  if (static_cast<unsigned int>(id) > this->GetTexturesStatus().size() || id == 0)
179  {
180  itkExceptionMacro(<< "Invalid texture index " << id << ", must be in [1;" << this->GetTexturesStatus().size() << "]");
181  }
182  else
183  {
184  this->GetFunctor().SetTextureStatus(id - 1, isSelected);
185  }
186  }
187 
188  std::vector<bool> GetTexturesStatus()
189  {
190  return this->GetFunctor().GetTexturesStatus();
191  }
192 
193  void InitFeatureStatus(bool status);
194 
196  const OutputImageType* GetLengthOutput() const;
197  OutputImageType* GetLengthOutput();
199 
201  const OutputImageType* GetWidthOutput() const;
202  OutputImageType* GetWidthOutput();
204 
206  const OutputImageType* GetPSIOutput() const;
207  OutputImageType* GetPSIOutput();
209 
211  const OutputImageType* GetWMeanOutput() const;
212  OutputImageType* GetWMeanOutput();
214 
216  const OutputImageType* GetRatioOutput() const;
217  OutputImageType* GetRatioOutput();
219 
221  const OutputImageType* GetSDOutput() const;
222  OutputImageType* GetSDOutput();
224 
225  void GenerateOutputInformation() override;
226  std::vector<FunctorType> m_FunctorList;
227 
228 protected:
231  {
232  }
233  void PrintSelf(std::ostream& os, itk::Indent indent) const override;
234 
235  void BeforeThreadedGenerateData() override;
236  void ThreadedGenerateData(const OutputImageRegionType& outputRegionForThread, itk::ThreadIdType threadId) override;
238  void GenerateInputRequestedRegion(void) override;
239 
240 private:
241  SFSTexturesImageFilter(const Self&) = delete;
242  void operator=(const Self&) = delete;
243 
244  unsigned int m_Radius;
246 };
247 
248 } // end namespace otb
249 
250 #ifndef OTB_MANUAL_INSTANTIATION
252 #endif
253 
254 #endif
This functor computes the texture describes in the following publication It is based on line directio...
InputImageType::RegionType InputImageRegionType
itk::SmartPointer< Self > Pointer
InputImageType::PixelType InputImagePixelType
void operator=(const Self &)=delete
itk::ImageToImageFilter< TInputImage, TOutputImage > Superclass
FunctorType::OutputType FunctorOutputType
void SetFeatureStatus(FeatureType id, bool isSelected)
Functor::SFSTexturesFunctor< NeighborhoodIteratorType, OutputImagePixelType > FunctorType
InputImageType::SizeType InputImageSizeType
NeighborhoodIteratorType::RadiusType RadiusType
void SetFunctor(const FunctorType &functor)
void SetSpatialThreshold(unsigned int thresh)
InputImagePixelType GetSpectralThreshold()
OutputImageType::RegionType OutputImageRegionType
const FunctorType & GetFunctor() const
itk::SmartPointer< const Self > ConstPointer
SFSTexturesImageFilter(const Self &)=delete
std::vector< FunctorType > m_FunctorList
itk::ConstNeighborhoodIterator< TInputImage > NeighborhoodIteratorType
InputImageType::ConstPointer InputImagePointerType
void SetSpectralThreshold(InputImagePixelType thresh)
OutputImageType::Pointer OutputImagePointerType
void SetNumberOfDirections(unsigned int D)
OutputImageType::PixelType OutputImagePixelType
void SetRatioMaxConsiderationNumber(unsigned int value)
Implements neighborhood-wise generic operation on image.
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
constexpr double CONST_PI
Definition: otbMath.h:49