OTB  10.0.0
Orfeo Toolbox
otbStereoSensorModelToElevationMapFilter.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 otbStereoSensorModelToElevationMapFilter_h
22 #define otbStereoSensorModelToElevationMapFilter_h
23 
24 #include "itkImageToImageFilter.h"
25 #include "itkInterpolateImageFunction.h"
26 #include "otbGenericRSTransform.h"
27 #include "itkNormalizedCorrelationImageToImageMetric.h"
28 
29 namespace otb
30 {
82 template <class TInputImage, class TOutputHeight>
83 class ITK_EXPORT StereoSensorModelToElevationFilter : public itk::ImageToImageFilter<TInputImage, TOutputHeight>
84 {
85 public:
86 
89  typedef itk::ImageToImageFilter<TInputImage, TOutputHeight> Superclass;
90  typedef itk::SmartPointer<Self> Pointer;
91  typedef itk::SmartPointer<const Self> ConstPointer;
92 
94  itkNewMacro(Self);
95 
97  itkTypeMacro(StereoSensorModelToElevationFilter, ImageToImageFilter);
98 
100  typedef TInputImage InputImageType;
101  typedef typename InputImageType::PixelType InputPixelType;
102  typedef TOutputHeight OutputImageType;
103  typedef typename OutputImageType::RegionType OutputRegionType;
104  typedef typename OutputImageType::PixelType HeightType;
105  typedef typename OutputImageType::PointType OutputPointType;
106 
107  typedef itk::InterpolateImageFunction<TInputImage, double> InterpolatorType;
108  typedef typename InterpolatorType::Pointer InterpolatorPointerType;
111  typedef itk::NormalizedCorrelationImageToImageMetric<TInputImage, TInputImage> CorrelationMetricType;
112  typedef typename InputImageType::SizeType RadiusType;
113 
115  void SetMasterInput(const TInputImage* image);
116 
118  void SetSlaveInput(const TInputImage* image);
119 
121  const TInputImage* GetMasterInput() const;
122 
124  const TInputImage* GetSlaveInput() const;
125 
127  TOutputHeight* GetCorrelationOutput();
128 
132  itkSetMacro(LowerElevation, double);
133 
137  itkSetMacro(HigherElevation, double);
138 
140  itkSetMacro(ElevationStep, double);
141 
145  itkSetMacro(CorrelationThreshold, double);
146 
150  itkGetMacro(CorrelationThreshold, double);
151 
155  itkSetMacro(VarianceThreshold, double);
156 
160  itkGetMacro(VarianceThreshold, double);
161 
166  itkSetMacro(SubtractInitialElevation, bool);
167 
172  itkGetMacro(SubtractInitialElevation, bool);
173 
178  itkBooleanMacro(SubtractInitialElevation);
179 
181  void SetRadius(unsigned int radius)
182  {
183  m_Radius.Fill(radius);
184  this->Modified();
185  }
187 
188 protected:
191 
194 
196  void DynamicThreadedGenerateData(const OutputRegionType& outputRegionForThread) override;
197 
199  void GenerateInputRequestedRegion(void) override;
200 
202  void BeforeThreadedGenerateData() override;
203 
209  void VerifyInputInformation() const override
210  {
211  }
212 
213 
214 private:
216  void operator=(const Self&) = delete;
217 
218  inline double Correlation(const std::vector<double>& master, const std::vector<double>& slave) const;
219 
222 
225 
228 
231 
234 
237 
241 
246 
249 };
250 } // End namespace otb
251 
252 #ifndef OTB_MANUAL_INSTANTIATION
254 #endif
255 
256 #endif
This is the class to handle generic remote sensing transform.
itk::SmartPointer< Self > Pointer
Using sensor models to perform correlation along epipolars.
void operator=(const Self &)=delete
itk::InterpolateImageFunction< TInputImage, double > InterpolatorType
itk::NormalizedCorrelationImageToImageMetric< TInputImage, TInputImage > CorrelationMetricType
itk::ImageToImageFilter< TInputImage, TOutputHeight > Superclass
StereoSensorModelToElevationFilter(const Self &)=delete
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.