OTB  10.0.0
Orfeo Toolbox
otbSarParametricMapFunction.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 otbSarParametricMapFunction_h
22 #define otbSarParametricMapFunction_h
23 
24 #include "itkImageFunction.h"
25 #include "itkPointSet.h"
26 #include "itkVariableSizeMatrix.h"
27 
28 namespace otb
29 {
30 
42 template <class TInputImage, class TCoordRep = float>
43 class ITK_EXPORT SarParametricMapFunction
44  : public itk::ImageFunction<TInputImage, typename itk::NumericTraits<typename TInputImage::PixelType>::ScalarRealType, TCoordRep>
45 {
46 public:
49  typedef itk::ImageFunction<TInputImage, typename itk::NumericTraits<typename TInputImage::PixelType>::ScalarRealType, TCoordRep> Superclass;
50  typedef itk::SmartPointer<Self> Pointer;
51  typedef itk::SmartPointer<const Self> ConstPointer;
52 
54  itkTypeMacro(SarParametricMapFunction, itk::ImageFunction);
55 
57  itkNewMacro(Self);
58 
60  typedef TInputImage InputImageType;
61  typedef typename InputImageType::PixelType InputPixelType;
62  typedef typename Superclass::OutputType OutputType;
63  typedef typename Superclass::IndexType IndexType;
64  typedef typename std::array<int, 2> ArrayIndexType;
65  typedef typename Superclass::ContinuousIndexType ContinuousIndexType;
66 
67  itkStaticConstMacro(ImageDimension, unsigned int, InputImageType::ImageDimension);
68 
69  typedef itk::PointSet<OutputType, ImageDimension> PointSetType;
70  typedef typename PointSetType::Pointer PointSetPointer;
71  typedef typename PointSetType::ConstPointer PointSetConstPointer;
72  typedef typename PointSetType::PointType PointType;
73  typedef typename PointSetType::PixelType PixelType;
74 
75  typedef itk::VariableSizeMatrix<double> MatrixType;
76 
78  typedef typename itk::NumericTraits<InputPixelType>::ScalarRealType RealType;
79 
81  RealType Evaluate(const PointType& point) const override;
82 
84  RealType EvaluateAtIndex(const IndexType& index) const override
85  {
86  PointType point;
87  point[0] = static_cast<typename PointType::ValueType>(index[0]);
88  point[1] = static_cast<typename PointType::ValueType>(index[1]);
89  return this->Evaluate(point);
90  }
92 
94  {
95  IndexType index;
96  this->ConvertContinuousIndexToNearestIndex(cindex, index);
97  return this->EvaluateAtIndex(index);
98  }
99 
102  itkGetConstObjectMacro(PointSet, PointSetType);
104  {
105  m_IsInitialize = false;
106  m_PointSet = val;
107  this->Modified();
108  }
110 
112  itkSetMacro(Coeff, MatrixType);
113  itkGetMacro(Coeff, MatrixType);
114  itkGetConstMacro(Coeff, MatrixType);
116 
118  itkGetConstReferenceMacro(IsInitialize, bool);
119 
121  void SetPolynomalSize(const IndexType PolynomalSize);
122  void SetPolynomalSize(const ArrayIndexType PolynomalSize);
124 
126  void EvaluateParametricCoefficient();
127 
129  void SetConstantValue(const RealType& value);
130 
131 protected:
134  {
135  }
136  void PrintSelf(std::ostream& os, itk::Indent indent) const override;
137 
138 private:
139  SarParametricMapFunction(const Self&) = delete;
140  void operator=(const Self&) = delete;
141 
142  double Horner(PointType point) const;
143 
147  double m_ProductWidth;
149 };
150 
151 } // end namespace otb
152 
153 #ifndef OTB_MANUAL_INSTANTIATION
155 #endif
156 
157 #endif
Evaluates a parametric bidimensionnal polynomial model from a PointSet.
RealType EvaluateAtIndex(const IndexType &index) const override
itkGetObjectMacro(PointSet, PointSetType)
itk::SmartPointer< const Self > ConstPointer
itk::PointSet< OutputType, ImageDimension > PointSetType
double m_ProductHeight
the width of the complete product (read from metadata)
bool m_IsInitialize
the width of the complete product (read from metadata)
PointSetType::ConstPointer PointSetConstPointer
RealType EvaluateAtContinuousIndex(const ContinuousIndexType &cindex) const override
the width of the complete product (read from metadata)
void SetPointSet(PointSetPointer val)
the width of the complete product (read from metadata)
Superclass::ContinuousIndexType ContinuousIndexType
itk::ImageFunction< TInputImage, typename itk::NumericTraits< typename TInputImage::PixelType >::ScalarRealType, TCoordRep > Superclass
itk::NumericTraits< InputPixelType >::ScalarRealType RealType
double m_ProductWidth
the width of the complete product (read from metadata)
itk::VariableSizeMatrix< double > MatrixType
InputImageType::PixelType InputPixelType
void operator=(const Self &)=delete
the width of the complete product (read from metadata)
MatrixType m_Coeff
the width of the complete product (read from metadata)
~SarParametricMapFunction() override
the width of the complete product (read from metadata)
PointSetPointer m_PointSet
the width of the complete product (read from metadata)
SarParametricMapFunction(const Self &)=delete
the width of the complete product (read from metadata)
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.