OTB  10.0.0
Orfeo Toolbox
otbFineRegistrationImageFilter.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 otbFineRegistrationImageFilter_h
22 #define otbFineRegistrationImageFilter_h
23 
24 #include "itkImageToImageFilter.h"
25 #include "itkInterpolateImageFunction.h"
26 #include "itkContinuousIndex.h"
27 
28 #include "itkTranslationTransform.h"
29 #include "itkImageToImageMetric.h"
30 
31 namespace otb
32 {
33 
73 template <class TInputImage, class T0utputCorrelation, class TOutputDisplacementField>
74 class ITK_EXPORT FineRegistrationImageFilter : public itk::ImageToImageFilter<TInputImage, T0utputCorrelation>
75 {
76 public:
79  typedef itk::ImageToImageFilter<TInputImage, T0utputCorrelation> Superclass;
80  typedef itk::SmartPointer<Self> Pointer;
81  typedef itk::SmartPointer<const Self> ConstPointer;
82 
84  itkNewMacro(Self);
85 
87  itkTypeMacro(FineRegistrationImageFilter, ImageToImageFilter);
88 
90  typedef typename T0utputCorrelation::RegionType OutputImageRegionType;
91  typedef typename TOutputDisplacementField::PixelType DisplacementValueType;
92  typedef typename TInputImage::Pointer InputImagePointerType;
93  typedef typename TInputImage::RegionType InputImageRegionType;
94  typedef typename TInputImage::SizeType SizeType;
95  typedef typename TInputImage::IndexType IndexType;
96  typedef typename TInputImage::SpacingType SpacingType;
97  typedef typename TInputImage::PointType PointType;
98  typedef typename TInputImage::OffsetType OffsetType;
99  typedef itk::InterpolateImageFunction<TInputImage, double> InterpolatorType;
100  typedef typename InterpolatorType::Pointer InterpolatorPointerType;
101  typedef itk::ContinuousIndex<double, 2> ContinuousIndexType;
102  typedef itk::ImageToImageMetric<TInputImage, TInputImage> MetricType;
103  typedef typename MetricType::Pointer MetricPointerType;
104  typedef itk::TranslationTransform<double, 2> TranslationType;
105  typedef typename TranslationType::Pointer TranslationPointerType;
106  typedef typename itk::Transform<double, 2, 2> TransformType;
107  typedef typename TransformType::Pointer TransformPointerType;
108 
110  itkSetObjectMacro(Metric, MetricType);
113 
115  itkSetObjectMacro(Interpolator, InterpolatorType);
118 
120  void SetFixedInput(const TInputImage* image);
121 
123  void SetMovingInput(const TInputImage* image);
124 
126  const TInputImage* GetFixedInput();
127  const TInputImage* GetMovingInput();
129 
131  TOutputDisplacementField* GetOutputDisplacementField();
132 
134  itkSetMacro(Radius, SizeType);
135  itkGetMacro(Radius, SizeType);
137 
139  itkSetMacro(SearchRadius, SizeType);
140  itkGetMacro(SearchRadius, SizeType);
142 
144  itkSetMacro(ConvergenceAccuracy, double);
145  itkGetMacro(ConvergenceAccuracy, double);
147 
149  itkSetMacro(SubPixelAccuracy, double);
150  itkGetMacro(SubPixelAccuracy, double);
152 
154  itkSetMacro(MaxIter, int);
155  itkGetMacro(MaxIter, int);
157 
159  itkSetMacro(Minimize, bool);
160  itkBooleanMacro(Minimize);
162 
164  itkSetMacro(UseSpacing, bool);
165  itkBooleanMacro(UseSpacing);
167 
169  itkSetMacro(InitialOffset, SpacingType);
170  itkGetConstReferenceMacro(InitialOffset, SpacingType);
172 
174  itkSetMacro(GridStep, OffsetType);
175  itkGetConstReferenceMacro(GridStep, OffsetType);
177 
179  void SetRadius(unsigned int radius)
180  {
181  m_Radius.Fill(radius);
182  }
183 
185  void SetSearchRadius(unsigned int radius)
186  {
187  m_SearchRadius.Fill(radius);
188  }
189 
191  void SetGridStep(unsigned int step)
192  {
193  m_GridStep.Fill(step);
194  }
195 
197  itkSetObjectMacro(Transform, TransformType);
198  itkGetConstObjectMacro(Transform, TransformType);
200 
201 protected:
204 
207 
209  void GenerateData() override;
210 
212  void GenerateInputRequestedRegion(void) override;
213 
215  void GenerateOutputInformation(void) override;
216 
217 private:
219  void operator=(const Self&) = delete;
220 
221  inline double callMetric(double val1, double val2, double& oldRes, bool& flag);
222  inline void updateOptParams(double potBestVal, double parx, double pary, // inputs
223  double& bestVal, typename TranslationType::ParametersType& optParams); // outputs
224  inline void updatePoints(double& gn, double& in1, double& in2, double& in3, // inputs
225  double& out1, double& out2, double& out3, double& out4); // outputs
226  inline void updateMinimize(double& a, double& b);
227 
230 
233 
236 
239 
243 
246 
249 
252 
255 
258 
261 
264 };
265 
266 } // end namespace otb
267 
268 #ifndef OTB_MANUAL_INSTANTIATION
270 #endif
271 
272 #endif
Computes a displacement field between two images using a given metric.
itk::Transform< double, 2, 2 > TransformType
itk::ImageToImageMetric< TInputImage, TInputImage > MetricType
TOutputDisplacementField::PixelType DisplacementValueType
itk::TranslationTransform< double, 2 > TranslationType
itk::InterpolateImageFunction< TInputImage, double > InterpolatorType
T0utputCorrelation::RegionType OutputImageRegionType
void operator=(const Self &)=delete
FineRegistrationImageFilter(const Self &)=delete
itkGetObjectMacro(Metric, MetricType)
itkGetObjectMacro(Interpolator, InterpolatorType)
itk::ImageToImageFilter< TInputImage, T0utputCorrelation > Superclass
itk::SmartPointer< const Self > ConstPointer
itk::ContinuousIndex< double, 2 > ContinuousIndexType
Class to overload method passed to virtual pure in ITK V4.
Definition: otbTransform.h:42
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.