OTB  6.7.0
Orfeo Toolbox
otbRealMomentsImageFunction.h
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2005-2019 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,
44  std::vector< std::vector<
45  typename itk::NumericTraits<
46  typename TInputImage::PixelType>::RealType > >,
47  TCoordRep>
48 {
49 public:
52  typedef itk::ImageFunction<TInputImage,
53  std::vector< std::vector<
54  typename itk::NumericTraits<
55  typename TInputImage::PixelType>::RealType > >,
56  TCoordRep> Superclass;
59 
61  itkTypeMacro(RealMomentsImageFunction, ImageFunction);
62 
64  itkNewMacro(Self);
65 
67  typedef TInputImage InputImageType;
68  typedef typename Superclass::IndexType IndexType;
69  typedef typename Superclass::ContinuousIndexType ContinuousIndexType;
70  typedef typename Superclass::PointType PointType;
71 
73  typedef float ScalarRealType;
74 
75  typedef TCoordRep CoordRepType;
76 
78  itkStaticConstMacro(ImageDimension, unsigned int,
79  InputImageType::ImageDimension);
80 
82  OutputType EvaluateAtIndex(const IndexType& index) const override;
83 
85  OutputType Evaluate(const PointType& point) const override
86  {
87  IndexType index;
88  this->ConvertPointToNearestIndex(point, index);
89  return this->EvaluateAtIndex(index);
90  }
92  const ContinuousIndexType& cindex) const override
93  {
94  IndexType index;
95  this->ConvertContinuousIndexToNearestIndex(cindex, index);
96  return this->EvaluateAtIndex(index);
97  }
99 
103  itkSetMacro( NeighborhoodRadius, unsigned int );
104  itkGetConstReferenceMacro( NeighborhoodRadius, unsigned int );
106 
107  itkSetMacro(Pmax, unsigned int);
108  itkGetConstReferenceMacro(Pmax, unsigned int);
109  itkSetMacro(Qmax, unsigned int);
110  itkGetConstReferenceMacro(Qmax, unsigned int);
111 
112 protected:
115  void PrintSelf(std::ostream& os, itk::Indent indent) const override;
116 
117 private:
118  RealMomentsImageFunction(const Self &) = delete;
119  void operator =(const Self&) = delete;
120 
121  unsigned int m_Pmax;
122  unsigned int m_Qmax;
123  unsigned int m_NeighborhoodRadius;
124 };
125 
126 } // namespace otb
127 
128 #ifndef OTB_MANUAL_INSTANTIATION
130 #endif
131 
132 #endif
OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &cindex) const override
itk::ImageFunction< TInputImage, std::vector< std::vector< typename itk::NumericTraits< typename TInputImage::PixelType >::RealType > >, TCoordRep > Superclass
OutputType Evaluate(const PointType &point) const override
itk::Index< Monteverdi_DIMENSION > IndexType
Definition: mvdTypes.h:133
Monteverdi_FLOATING_TYPE RealType
Definition: mvdTypes.h:84
Superclass::ContinuousIndexType ContinuousIndexType
itk::SmartPointer< const Self > ConstPointer
Calculate the moment values in the specified neighborhood.
VectorImageType::PointType PointType
Definition: mvdTypes.h:189