OTB  10.0.0
Orfeo Toolbox
otbRadiometricMomentsImageFunction.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 otbRadiometricMomentsImageFunction_h
22 #define otbRadiometricMomentsImageFunction_h
23 
24 #include "itkImageFunction.h"
26 #include "itkConstNeighborhoodIterator.h"
27 #include "itkFixedArray.h"
28 
29 namespace otb
30 {
31 
47 template <class TInputImage, class TCoordRep = double>
49  : public itk::ImageFunction<TInputImage, itk::FixedArray<typename itk::NumericTraits<typename TInputImage::PixelType>::RealType, 4>, TCoordRep>
50 {
51 public:
54  typedef itk::ImageFunction<TInputImage, itk::FixedArray<typename itk::NumericTraits<typename TInputImage::PixelType>::RealType, 4>, TCoordRep> Superclass;
55  typedef itk::SmartPointer<Self> Pointer;
56  typedef itk::SmartPointer<const Self> ConstPointer;
57 
59  itkTypeMacro(RadiometricMomentsImageFunction, ImageFunction);
60 
62  itkNewMacro(Self);
63 
65  typedef TInputImage InputImageType;
66  typedef typename Superclass::IndexType IndexType;
67  typedef typename Superclass::ContinuousIndexType ContinuousIndexType;
68  typedef typename Superclass::PointType PointType;
69 
70  typedef typename Superclass::OutputType OutputType;
71  typedef typename OutputType::ValueType ScalarRealType;
72 
73  typedef TCoordRep CoordRepType;
74 
76 
78  itkStaticConstMacro(ImageDimension, unsigned int, InputImageType::ImageDimension);
79 
81  OutputType EvaluateAtIndex(const IndexType& index) const override;
82 
84  OutputType Evaluate(const PointType& point) const override
85  {
86  IndexType index;
87  this->ConvertPointToNearestIndex(point, index);
88  return this->EvaluateAtIndex(index);
89  }
91  {
92  IndexType index;
93  this->ConvertContinuousIndexToNearestIndex(cindex, index);
94  return this->EvaluateAtIndex(index);
95  }
97 
101  itkSetMacro(NeighborhoodRadius, unsigned int);
102  itkGetConstReferenceMacro(NeighborhoodRadius, unsigned int);
104 
105 protected:
108  {
109  }
110  void PrintSelf(std::ostream& os, itk::Indent indent) const override;
111 
112 private:
114  void operator=(const Self&) = delete;
115 
116  unsigned int m_NeighborhoodRadius;
118 };
119 
120 } // namespace otb
121 
122 #ifndef OTB_MANUAL_INSTANTIATION
124 #endif
125 
126 #endif
OutputType Evaluate(const PointType &point) const override
Functor::RadiometricMomentsFunctor< itk::ConstNeighborhoodIterator< InputImageType >, ScalarRealType > FunctorType
itk::ImageFunction< TInputImage, itk::FixedArray< typename itk::NumericTraits< typename TInputImage::PixelType >::RealType, 4 >, TCoordRep > Superclass
OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &cindex) const override
RadiometricMomentsImageFunction(const Self &)=delete
void operator=(const Self &)=delete
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.