OTB  10.0.0
Orfeo Toolbox
otbComplexMomentsImageFunction.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 otbComplexMomentsImageFunction_h
22 #define otbComplexMomentsImageFunction_h
23 
24 #include "itkImageFunction.h"
25 
26 #include <complex>
27 
28 namespace otb
29 {
30 
55 template <class TInputImage, class TCoordRep = double>
56 class ITK_EXPORT ComplexMomentsImageFunction : public itk::ImageFunction<TInputImage, std::vector<std::vector<std::complex<double>>>, TCoordRep>
57 {
58 public:
61  typedef itk::ImageFunction<TInputImage, std::vector<std::vector<std::complex<double>>>, TCoordRep> Superclass;
62  typedef itk::SmartPointer<Self> Pointer;
63  typedef itk::SmartPointer<const Self> ConstPointer;
64 
66  itkTypeMacro(ComplexMomentsImageFunction, ImageFunction);
67 
69  itkNewMacro(Self);
70 
72  typedef TInputImage InputImageType;
73  typedef typename Superclass::IndexType IndexType;
74  typedef typename Superclass::ContinuousIndexType ContinuousIndexType;
75  typedef typename Superclass::PointType PointType;
76  typedef typename Superclass::OutputType OutputType;
77 
78  typedef double ScalarRealType;
79  typedef typename std::complex<ScalarRealType> ScalarComplexType;
80 
81  typedef TCoordRep CoordRepType;
82 
84  itkStaticConstMacro(ImageDimension, unsigned int, InputImageType::ImageDimension);
85 
87  OutputType EvaluateAtIndex(const IndexType& index) const override;
88 
90  OutputType Evaluate(const PointType& point) const override
91  {
92  IndexType index;
93  this->ConvertPointToNearestIndex(point, index);
94  return this->EvaluateAtIndex(index);
95  }
97  {
98  IndexType index;
99  this->ConvertContinuousIndexToNearestIndex(cindex, index);
100  return this->EvaluateAtIndex(index);
101  }
103 
107  itkSetMacro(NeighborhoodRadius, unsigned int);
108  itkGetConstReferenceMacro(NeighborhoodRadius, unsigned int);
110 
111  itkSetMacro(Pmax, unsigned int);
112  itkGetConstReferenceMacro(Pmax, unsigned int);
113  itkSetMacro(Qmax, unsigned int);
114  itkGetConstReferenceMacro(Qmax, unsigned int);
115 
116 protected:
119  {
120  }
121  void PrintSelf(std::ostream& os, itk::Indent indent) const override;
122 
123 private:
125  void operator=(const Self&) = delete;
126 
127  unsigned int m_Pmax;
128  unsigned int m_Qmax;
129  unsigned int m_NeighborhoodRadius;
130 };
131 
132 } // namespace otb
133 
134 #ifndef OTB_MANUAL_INSTANTIATION
136 #endif
137 
138 #endif
Calculate the complex moment values in the specified neighborhood.
std::complex< ScalarRealType > ScalarComplexType
ComplexMomentsImageFunction(const Self &)=delete
itk::ImageFunction< TInputImage, std::vector< std::vector< std::complex< double > > >, TCoordRep > Superclass
OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &cindex) const override
OutputType Evaluate(const PointType &point) const override
itk::SmartPointer< const Self > ConstPointer
Superclass::ContinuousIndexType ContinuousIndexType
void operator=(const Self &)=delete
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.