OTB  9.0.0
Orfeo Toolbox
otbRealMomentsImageFunction.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 otbRealMomentsImageFunction_h
22 #define otbRealMomentsImageFunction_h
23 
24 #include "itkImageFunction.h"
25 
26 namespace otb
27 {
28 
41 template <class TInputImage, class TCoordRep = double>
42 class ITK_EXPORT RealMomentsImageFunction
43  : public itk::ImageFunction<TInputImage, std::vector<std::vector<typename itk::NumericTraits<typename TInputImage::PixelType>::RealType>>, TCoordRep>
44 {
45 public:
48  typedef itk::ImageFunction<TInputImage, std::vector<std::vector<typename itk::NumericTraits<typename TInputImage::PixelType>::RealType>>, TCoordRep>
50  typedef itk::SmartPointer<Self> Pointer;
51  typedef itk::SmartPointer<const Self> ConstPointer;
52 
54  itkTypeMacro(RealMomentsImageFunction, ImageFunction);
55 
57  itkNewMacro(Self);
58 
60  typedef TInputImage InputImageType;
61  typedef typename Superclass::IndexType IndexType;
62  typedef typename Superclass::ContinuousIndexType ContinuousIndexType;
63  typedef typename Superclass::PointType PointType;
64 
65  typedef typename Superclass::OutputType OutputType;
66  typedef float ScalarRealType;
67 
68  typedef TCoordRep CoordRepType;
69 
71  itkStaticConstMacro(ImageDimension, unsigned int, InputImageType::ImageDimension);
72 
74  OutputType EvaluateAtIndex(const IndexType& index) const override;
75 
77  OutputType Evaluate(const PointType& point) const override
78  {
79  IndexType index;
80  this->ConvertPointToNearestIndex(point, index);
81  return this->EvaluateAtIndex(index);
82  }
84  {
85  IndexType index;
86  this->ConvertContinuousIndexToNearestIndex(cindex, index);
87  return this->EvaluateAtIndex(index);
88  }
90 
94  itkSetMacro(NeighborhoodRadius, unsigned int);
95  itkGetConstReferenceMacro(NeighborhoodRadius, unsigned int);
97 
98  itkSetMacro(Pmax, unsigned int);
99  itkGetConstReferenceMacro(Pmax, unsigned int);
100  itkSetMacro(Qmax, unsigned int);
101  itkGetConstReferenceMacro(Qmax, unsigned int);
102 
103 protected:
106  {
107  }
108  void PrintSelf(std::ostream& os, itk::Indent indent) const override;
109 
110 private:
111  RealMomentsImageFunction(const Self&) = delete;
112  void operator=(const Self&) = delete;
113 
114  unsigned int m_Pmax;
115  unsigned int m_Qmax;
116  unsigned int m_NeighborhoodRadius;
117 };
118 
119 } // namespace otb
120 
121 #ifndef OTB_MANUAL_INSTANTIATION
123 #endif
124 
125 #endif
otb::RealMomentsImageFunction::m_Pmax
unsigned int m_Pmax
Definition: otbRealMomentsImageFunction.h:114
otb::RealMomentsImageFunction::OutputType
Superclass::OutputType OutputType
Definition: otbRealMomentsImageFunction.h:65
otb
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Definition: otbJoinContainer.h:32
otb::RealMomentsImageFunction::InputImageType
TInputImage InputImageType
Definition: otbRealMomentsImageFunction.h:57
otb::RealMomentsImageFunction::Superclass
itk::ImageFunction< TInputImage, std::vector< std::vector< typename itk::NumericTraits< typename TInputImage::PixelType >::RealType > >, TCoordRep > Superclass
Definition: otbRealMomentsImageFunction.h:49
otbRealMomentsImageFunction.hxx
otb::RealMomentsImageFunction::IndexType
Superclass::IndexType IndexType
Definition: otbRealMomentsImageFunction.h:61
otb::RealMomentsImageFunction::m_NeighborhoodRadius
unsigned int m_NeighborhoodRadius
Definition: otbRealMomentsImageFunction.h:116
otb::RealMomentsImageFunction::Self
RealMomentsImageFunction Self
Definition: otbRealMomentsImageFunction.h:47
otb::RealMomentsImageFunction::m_Qmax
unsigned int m_Qmax
Definition: otbRealMomentsImageFunction.h:115
otb::RealMomentsImageFunction::ContinuousIndexType
Superclass::ContinuousIndexType ContinuousIndexType
Definition: otbRealMomentsImageFunction.h:62
otb::RealMomentsImageFunction::Evaluate
OutputType Evaluate(const PointType &point) const override
Definition: otbRealMomentsImageFunction.h:77
otb::RealMomentsImageFunction::CoordRepType
TCoordRep CoordRepType
Definition: otbRealMomentsImageFunction.h:68
otb::RealMomentsImageFunction::PointType
Superclass::PointType PointType
Definition: otbRealMomentsImageFunction.h:63
otb::RealMomentsImageFunction::Pointer
itk::SmartPointer< Self > Pointer
Definition: otbRealMomentsImageFunction.h:50
otb::RealMomentsImageFunction::EvaluateAtContinuousIndex
OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &cindex) const override
Definition: otbRealMomentsImageFunction.h:83
otb::RealMomentsImageFunction
Calculate the moment values in the specified neighborhood.
Definition: otbRealMomentsImageFunction.h:42
otb::RealMomentsImageFunction::ConstPointer
itk::SmartPointer< const Self > ConstPointer
Definition: otbRealMomentsImageFunction.h:51
otb::RealMomentsImageFunction::ScalarRealType
float ScalarRealType
Definition: otbRealMomentsImageFunction.h:66
otb::RealMomentsImageFunction::~RealMomentsImageFunction
~RealMomentsImageFunction() override
Definition: otbRealMomentsImageFunction.h:105