OTB  10.0.0
Orfeo Toolbox
otbBSplineInterpolateImageFunction.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 otbBSplineInterpolateImageFunction_h
22 #define otbBSplineInterpolateImageFunction_h
23 
24 #include <vector>
25 
26 #include "itkInterpolateImageFunction.h"
27 #include "vnl/vnl_matrix.h"
28 
30 #include "itkConceptChecking.h"
31 #include "itkCovariantVector.h"
32 
33 namespace otb
34 {
47 template <class TImageType, class TCoordRep = double, class TCoefficientType = double>
48 class ITK_EXPORT BSplineInterpolateImageFunction : public itk::InterpolateImageFunction<TImageType, TCoordRep>
49 {
50 public:
53  typedef itk::InterpolateImageFunction<TImageType, TCoordRep> Superclass;
54  typedef itk::SmartPointer<Self> Pointer;
55  typedef itk::SmartPointer<const Self> ConstPointer;
56 
58  itkTypeMacro(BSplineInterpolateImageFunction, InterpolateImageFunction);
59 
61  itkNewMacro(Self);
62 
64  typedef typename Superclass::OutputType OutputType;
65 
67  typedef typename Superclass::InputImageType InputImageType;
68 
70  itkStaticConstMacro(ImageDimension, unsigned int, Superclass::ImageDimension);
71 
73  typedef typename Superclass::IndexType IndexType;
74 
76  typedef typename InputImageType::RegionType RegionType;
77 
79  typedef typename Superclass::ContinuousIndexType ContinuousIndexType;
80 
82  typedef typename Superclass::PointType PointType;
83 
84  typedef typename Superclass::SizeType SizeType;
85 
87  typedef itk::ImageLinearIteratorWithIndex<TImageType> Iterator;
88 
90  typedef TCoefficientType CoefficientDataType;
91  typedef itk::Image<CoefficientDataType, itkGetStaticConstMacro(ImageDimension)> CoefficientImageType;
92 
96 
105  OutputType EvaluateAtContinuousIndex(const ContinuousIndexType& index) const override;
106 
108  typedef itk::CovariantVector<OutputType, itkGetStaticConstMacro(ImageDimension)> CovariantVectorType;
109 
111  {
112  ContinuousIndexType index;
113  this->GetInputImage()->TransformPhysicalPointToContinuousIndex(point, index);
114  return (this->EvaluateDerivativeAtContinuousIndex(index));
115  }
116 
117  CovariantVectorType EvaluateDerivativeAtContinuousIndex(const ContinuousIndexType& x) const;
118 
121  void SetSplineOrder(unsigned int SplineOrder);
122  itkGetMacro(SplineOrder, int);
124 
126  void SetInputImage(const TImageType* inputData) override;
127 
130  virtual void UpdateCoefficientsFilter(void);
131 
132  virtual SizeType GetRadius() const
133  {
134  typename itk::InterpolateImageFunction<TImageType, TCoordRep>::SizeType radius({2,2});
135  return radius;
136  }
137 
138 protected:
140  ~BSplineInterpolateImageFunction() override = default;
141 
142  void operator=(const Self&) = delete;
143  void PrintSelf(std::ostream& os, itk::Indent indent) const override;
144 
145  // These are needed by the smoothing spline routine.
146  std::vector<CoefficientDataType> m_Scratch; // temp storage for processing of Coefficients
147  typename TImageType::SizeType m_DataLength; // Image size
148  unsigned int m_SplineOrder; // User specified spline order (3rd or cubic is the default)
149 
150  typename CoefficientImageType::ConstPointer m_Coefficients; // Spline coefficients
151 
152 private:
155  void SetInterpolationWeights(const ContinuousIndexType& x, const vnl_matrix<long>& EvaluateIndex, vnl_matrix<double>& weights,
156  unsigned int splineOrder) const;
157 
159  void SetDerivativeWeights(const ContinuousIndexType& x, const vnl_matrix<long>& EvaluateIndex, vnl_matrix<double>& weights, unsigned int splineOrder) const;
160 
163  void GeneratePointsToIndex();
164 
166  void DetermineRegionOfSupport(vnl_matrix<long>& evaluateIndex, const ContinuousIndexType& x, unsigned int splineOrder) const;
167 
170  void ApplyMirrorBoundaryConditions(vnl_matrix<long>& evaluateIndex, unsigned int splineOrder) const;
171 
172  Iterator m_CIterator; // Iterator for traversing spline coefficients.
173  unsigned long m_MaxNumberInterpolationPoints; // number of neighborhood points used for interpolation
174  std::vector<IndexType> m_PointsToIndex; // Preallocation of interpolation neighborhood indices
175 
177 
179 };
180 
181 } // namespace otb
182 
183 #ifndef OTB_MANUAL_INSTANTIATION
185 #endif
186 
187 #endif
This class is an evolution of the itk::BSplineDecompositionImageFilter to handle huge images with thi...
This class is an evolution of the itk::BSplineInterpolateImageFunction to handle huge images with thi...
itk::ImageLinearIteratorWithIndex< TImageType > Iterator
BSplineInterpolateImageFunction(const Self &)=delete
~BSplineInterpolateImageFunction() override=default
itk::CovariantVector< OutputType, itkGetStaticConstMacro(ImageDimension)> CovariantVectorType
otb::BSplineDecompositionImageFilter< TImageType, CoefficientImageType > CoefficientFilter
void operator=(const Self &)=delete
CovariantVectorType EvaluateDerivative(const PointType &point) const
itk::Image< CoefficientDataType, itkGetStaticConstMacro(ImageDimension)> CoefficientImageType
itk::InterpolateImageFunction< TImageType, TCoordRep > Superclass
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.