OTB  9.0.0
Orfeo Toolbox
otbFlusserMomentsImageFunction.h
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2005-2022 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 otbFlusserMomentsImageFunction_h
22 #define otbFlusserMomentsImageFunction_h
23 
24 #include "itkImageFunction.h"
25 #include "itkFixedArray.h"
26 
27 namespace otb
28 {
29 
69 
70 template <class TInputImage, class TCoordRep = double>
72  : public itk::ImageFunction<TInputImage, itk::FixedArray<typename itk::NumericTraits<typename TInputImage::PixelType>::RealType, 11>, TCoordRep>
73 {
74 public:
77  typedef itk::ImageFunction<TInputImage, itk::FixedArray<typename itk::NumericTraits<typename TInputImage::PixelType>::RealType, 11>, TCoordRep> Superclass;
78  typedef itk::SmartPointer<Self> Pointer;
79  typedef itk::SmartPointer<const Self> ConstPointer;
80 
82  itkTypeMacro(FlusserMomentsImageFunction, ImageFunction);
83 
85  itkNewMacro(Self);
86 
88  typedef TInputImage InputImageType;
89  typedef typename Superclass::IndexType IndexType;
90  typedef typename Superclass::ContinuousIndexType ContinuousIndexType;
91  typedef typename Superclass::PointType PointType;
92 
93  typedef typename Superclass::OutputType OutputType;
94  typedef typename OutputType::ValueType ScalarRealType;
95 
96  typedef TCoordRep CoordRepType;
97 
99  itkStaticConstMacro(ImageDimension, unsigned int, InputImageType::ImageDimension);
100 
102  OutputType EvaluateAtIndex(const IndexType& index) const override;
103 
105  OutputType Evaluate(const PointType& point) const override
106  {
107  IndexType index;
108  this->ConvertPointToNearestIndex(point, index);
109  return this->EvaluateAtIndex(index);
110  }
112  {
113  IndexType index;
114  this->ConvertContinuousIndexToNearestIndex(cindex, index);
115  return this->EvaluateAtIndex(index);
116  }
118 
122  itkSetMacro(NeighborhoodRadius, unsigned int);
123  itkGetConstReferenceMacro(NeighborhoodRadius, unsigned int);
125 
126 protected:
129  {
130  }
131  void PrintSelf(std::ostream& os, itk::Indent indent) const override;
132 
133 private:
134  FlusserMomentsImageFunction(const Self&) = delete;
135  void operator=(const Self&) = delete;
136 
137  unsigned int m_NeighborhoodRadius;
138 };
139 
140 } // namespace otb
141 
142 #ifndef OTB_MANUAL_INSTANTIATION
144 #endif
145 
146 #endif
otb::FlusserMomentsImageFunction::ContinuousIndexType
Superclass::ContinuousIndexType ContinuousIndexType
Definition: otbFlusserMomentsImageFunction.h:90
otb::FlusserMomentsImageFunction::Superclass
itk::ImageFunction< TInputImage, itk::FixedArray< typename itk::NumericTraits< typename TInputImage::PixelType >::RealType, 11 >, TCoordRep > Superclass
Definition: otbFlusserMomentsImageFunction.h:77
otb::FlusserMomentsImageFunction::EvaluateAtContinuousIndex
OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &cindex) const override
Definition: otbFlusserMomentsImageFunction.h:111
otb::FlusserMomentsImageFunction::OutputType
Superclass::OutputType OutputType
Definition: otbFlusserMomentsImageFunction.h:93
otb
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Definition: otbJoinContainer.h:32
otb::FlusserMomentsImageFunction::m_NeighborhoodRadius
unsigned int m_NeighborhoodRadius
Definition: otbFlusserMomentsImageFunction.h:137
otb::FlusserMomentsImageFunction::Self
FlusserMomentsImageFunction Self
Definition: otbFlusserMomentsImageFunction.h:76
otb::FlusserMomentsImageFunction
Calculate the Flusser's invariant parameters.
Definition: otbFlusserMomentsImageFunction.h:71
otb::FlusserMomentsImageFunction::IndexType
Superclass::IndexType IndexType
Definition: otbFlusserMomentsImageFunction.h:89
otb::FlusserMomentsImageFunction::PointType
Superclass::PointType PointType
Definition: otbFlusserMomentsImageFunction.h:91
otb::FlusserMomentsImageFunction::~FlusserMomentsImageFunction
~FlusserMomentsImageFunction() override
Definition: otbFlusserMomentsImageFunction.h:128
otb::FlusserMomentsImageFunction::ScalarRealType
OutputType::ValueType ScalarRealType
Definition: otbFlusserMomentsImageFunction.h:94
otb::FlusserMomentsImageFunction::CoordRepType
TCoordRep CoordRepType
Definition: otbFlusserMomentsImageFunction.h:96
otb::FlusserMomentsImageFunction::InputImageType
TInputImage InputImageType
Definition: otbFlusserMomentsImageFunction.h:85
otbFlusserMomentsImageFunction.hxx
otb::FlusserMomentsImageFunction::ConstPointer
itk::SmartPointer< const Self > ConstPointer
Definition: otbFlusserMomentsImageFunction.h:79
otb::FlusserMomentsImageFunction::Evaluate
OutputType Evaluate(const PointType &point) const override
Definition: otbFlusserMomentsImageFunction.h:105
otb::FlusserMomentsImageFunction::Pointer
itk::SmartPointer< Self > Pointer
Definition: otbFlusserMomentsImageFunction.h:78