OTB  9.0.0
Orfeo Toolbox
otbStereorectificationDisplacementFieldSource.hxx
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2005-2022 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_hxx
22 #define otbStereoSensorModelToElevationMapFilter_hxx
23 
25 #include "itkProgressReporter.h"
26 #include "otbMath.h"
27 
28 // For partial specialization
29 #include "otbVectorImage.h"
30 
31 #include "otbDEMHandler.h"
32 
33 namespace otb
34 {
35 
36 template <class TInputImage, class TOutputImage>
38  : m_ElevationOffset(50),
39  m_Scale(1),
40  m_GridStep(1),
41  m_LeftImage(),
42  m_RightImage(),
43  m_LeftToRightTransform(),
44  m_RightToLeftTransform(),
45  m_OutputOriginInLeftImage(),
46  m_MeanBaselineRatio(0),
47  m_UseDEM(false)
48 {
49  // Set the number of outputs to 2 (one deformation field for each
50  // image)
51  this->SetNumberOfRequiredOutputs(2);
52 
53  // Create the 2nd input (not created by default)
54  this->SetNthOutput(0, OutputImageType::New());
55  this->SetNthOutput(1, OutputImageType::New());
56 
57  // Build the RS Transforms
60 }
61 
62 template <class TInputImage, class TOutputImage>
65 {
66  if (this->GetNumberOfOutputs() < 1)
67  {
68  return 0;
69  }
70  return static_cast<const OutputImageType*>(this->itk::ProcessObject::GetOutput(0));
71 }
72 
73 template <class TInputImage, class TOutputImage>
76 {
77  if (this->GetNumberOfOutputs() < 1)
78  {
79  return nullptr;
80  }
81  return static_cast<OutputImageType*>(this->itk::ProcessObject::GetOutput(0));
82 }
83 
84 template <class TInputImage, class TOutputImage>
87 {
88  if (this->GetNumberOfOutputs() < 2)
89  {
90  return 0;
91  }
92  return static_cast<const OutputImageType*>(this->itk::ProcessObject::GetOutput(1));
93 }
94 
95 template <class TInputImage, class TOutputImage>
98 {
99  if (this->GetNumberOfOutputs() < 2)
100  {
101  return nullptr;
102  }
103  return static_cast<OutputImageType*>(this->itk::ProcessObject::GetOutput(1));
104 }
105 
106 template <class TInputImage, class TOutputImage>
108 {
109  // Superclass::GenerateOutputInformation();
110 
111  // Ensure that both left image and right image are available
112  if (m_LeftImage.IsNull() || m_RightImage.IsNull())
113  {
114  itkExceptionMacro(<< "Either left image or right image pointer is null, can not perform stereo-rectification.");
115  }
116 
117  // Ensure input images have up-to-date information
118  m_LeftImage->UpdateOutputInformation();
119  m_RightImage->UpdateOutputInformation();
120 
121  // Setup the DEM handler if needed
122  auto & demHandler = DEMHandler::GetInstance();
123 
124  // Set-up a transform to use the DEMHandler
125  typedef otb::GenericRSTransform<> RSTransform2DType;
126  RSTransform2DType::Pointer leftToGroundTransform = RSTransform2DType::New();
127  leftToGroundTransform->SetInputImageMetadata(&(m_LeftImage->GetImageMetadata()));
128 
129  leftToGroundTransform->InstantiateTransform();
130 
131  // Retrieve the deformation field pointers
132  OutputImageType* leftDFPtr = this->GetLeftDisplacementFieldOutput();
133  OutputImageType* rightDFPtr = this->GetRightDisplacementFieldOutput();
134 
135  // Set up the RS transforms
136  m_LeftToRightTransform->SetInputImageMetadata(&(m_LeftImage->GetImageMetadata()));
137  m_LeftToRightTransform->SetOutputImageMetadata(&(m_RightImage->GetImageMetadata()));
138  m_LeftToRightTransform->InstantiateTransform();
139 
140  m_RightToLeftTransform->SetInputImageMetadata(&(m_RightImage->GetImageMetadata()));
141  m_RightToLeftTransform->SetOutputImageMetadata(&(m_LeftImage->GetImageMetadata()));
142  m_RightToLeftTransform->InstantiateTransform();
143 
144  // Now, we must determine the optimized size, spacing and origin of the
145  // stereo-rectified images, as well as the position of the origin in
146  // the left image
147 
148  // First, spacing : choose a square spacing,
149  SpacingType outputSpacing;
150  outputSpacing.Fill(m_Scale * m_GridStep);
151  double mean_spacing = 0.5 * (std::abs(m_LeftImage->GetSignedSpacing()[0]) + std::abs(m_LeftImage->GetSignedSpacing()[1]));
152  // double ratio_x = mean_spacing / std::abs(m_LeftImage->GetSignedSpacing()[0]);
153  // double ratio_y = mean_spacing / std::abs(m_LeftImage->GetSignedSpacing()[1]);
154 
155  outputSpacing[0] *= mean_spacing;
156  outputSpacing[1] *= mean_spacing;
157 
158  // Then, we retrieve the origin of the left input image
159  double localElevation = demHandler.GetDefaultHeightAboveEllipsoid();
160 
161  if (m_UseDEM)
162  {
163  RSTransform2DType::InputPointType tmpPoint;
164  localElevation = demHandler.GetHeightAboveEllipsoid(leftToGroundTransform->TransformPoint(m_LeftImage->GetOrigin()));
165  }
166 
167  TDPointType leftInputOrigin;
168  leftInputOrigin[0] = m_LeftImage->GetOrigin()[0];
169  leftInputOrigin[1] = m_LeftImage->GetOrigin()[1];
170  leftInputOrigin[2] = localElevation;
171 
172  // Next, we will compute the parameters of the local epipolar line
173  // at the left image origin
174  TDPointType rightEpiPoint, leftEpiLineStart, leftEpiLineEnd;
175 
176  // This point is the image of the left input image origin at the
177  // average elevation
178  rightEpiPoint = m_LeftToRightTransform->TransformPoint(leftInputOrigin);
179 
180  // The beginning of the epipolar line in the left image is the image
181  // of rightEpiPoint at a lower elevation (using the offset)
182  rightEpiPoint[2] = localElevation - m_ElevationOffset;
183  leftEpiLineStart = m_RightToLeftTransform->TransformPoint(rightEpiPoint);
184 
185  // The ending of the epipolar line in the left image is the image
186  // of rightEpiPoint at a higher elevation (using the offset)
187  rightEpiPoint[2] = localElevation + m_ElevationOffset;
188  leftEpiLineEnd = m_RightToLeftTransform->TransformPoint(rightEpiPoint);
189 
190  // Now, we can compute the equation of the epipolar line y = a*x+b
191  // epipolar angle is computed in left image physical space
192  double alpha = 0;
193  if (leftEpiLineEnd[0] == leftEpiLineStart[0])
194  {
195  if (leftEpiLineEnd[1] > leftEpiLineStart[1])
196  {
197  alpha = 0.5 * otb::CONST_PI;
198  }
199  else
200  {
201  alpha = -0.5 * otb::CONST_PI;
202  }
203  }
204  else
205  {
206  double a = (leftEpiLineEnd[1] - leftEpiLineStart[1]) / (leftEpiLineEnd[0] - leftEpiLineStart[0]);
207  if (leftEpiLineEnd[0] > leftEpiLineStart[0])
208  {
209  alpha = std::atan(a);
210  }
211  else
212  {
213  alpha = otb::CONST_PI + std::atan(a);
214  }
215  }
216 
217  // And compute the unitary vectors of the new axis (equivalent to
218  // the column of the rotation matrix)
219  // TODO: Check if we need to use the input image spacing here
220  double ux = std::cos(alpha);
221  double uy = std::sin(alpha);
222  double vx = -std::sin(alpha);
223  double vy = std::cos(alpha);
224 
225  // Now, we will compute the bounding box of the left input image in
226  // this rotated geometry
227 
228  // First we compute coordinates of the 4 corners (we omit ulx which
229  // coordinates are {0,0})
230  double urx = ux * m_LeftImage->GetLargestPossibleRegion().GetSize()[0] * m_LeftImage->GetSignedSpacing()[0];
231  double ury = vx * m_LeftImage->GetLargestPossibleRegion().GetSize()[0] * m_LeftImage->GetSignedSpacing()[0];
232  double llx = uy * m_LeftImage->GetLargestPossibleRegion().GetSize()[1] * m_LeftImage->GetSignedSpacing()[1];
233  double lly = vy * m_LeftImage->GetLargestPossibleRegion().GetSize()[1] * m_LeftImage->GetSignedSpacing()[1];
234  double lrx = ux * m_LeftImage->GetLargestPossibleRegion().GetSize()[0] * m_LeftImage->GetSignedSpacing()[0] +
235  uy * m_LeftImage->GetLargestPossibleRegion().GetSize()[1] * m_LeftImage->GetSignedSpacing()[1];
236  double lry = vx * m_LeftImage->GetLargestPossibleRegion().GetSize()[0] * m_LeftImage->GetSignedSpacing()[0] +
237  vy * m_LeftImage->GetLargestPossibleRegion().GetSize()[1] * m_LeftImage->GetSignedSpacing()[1];
238 
239  // Bounding box (this time we do not omit ulx)
240  double minx = std::min(std::min(std::min(urx, llx), lrx), 0.);
241  double miny = std::min(std::min(std::min(ury, lly), lry), 0.);
242  double maxx = std::max(std::max(std::max(urx, llx), lrx), 0.);
243  double maxy = std::max(std::max(std::max(ury, lly), lry), 0.);
244 
245  // We can now estimate the output image size, taking into account
246  // the scale parameter
247  m_RectifiedImageSize[0] = static_cast<unsigned int>((maxx - minx) / (mean_spacing * m_Scale));
248  m_RectifiedImageSize[1] = static_cast<unsigned int>((maxy - miny) / (mean_spacing * m_Scale));
249 
250  // Now, we can compute the origin of the epipolar images in the left
251  // input image geometry (we rotate back)
252  m_OutputOriginInLeftImage[0] = leftInputOrigin[0] + (ux * minx + vx * miny);
253  m_OutputOriginInLeftImage[1] = leftInputOrigin[1] + (uy * minx + vy * miny);
254  m_OutputOriginInLeftImage[2] = localElevation;
255 
256  // And also the size of the deformation field
257  SizeType outputSize;
258  outputSize[0] = (m_RectifiedImageSize[0] / m_GridStep + 2);
259  outputSize[1] = (m_RectifiedImageSize[1] / m_GridStep + 2);
260 
261  // Build the output largest region
262  RegionType outputLargestRegion;
263  outputLargestRegion.SetSize(outputSize);
264 
265  // Update the information
266  leftDFPtr->SetLargestPossibleRegion(outputLargestRegion);
267  rightDFPtr->SetLargestPossibleRegion(outputLargestRegion);
268 
269  leftDFPtr->SetSignedSpacing(outputSpacing);
270  rightDFPtr->SetSignedSpacing(outputSpacing);
271 
272  leftDFPtr->SetNumberOfComponentsPerPixel(2);
273  rightDFPtr->SetNumberOfComponentsPerPixel(2);
274 }
275 
276 template <class TInputImage, class TOutputImage>
278 {
279  // Call superclass
280  Superclass::EnlargeOutputRequestedRegion(output);
281 
282  // Retrieve the deformation field pointers
283  OutputImageType* leftDFPtr = this->GetLeftDisplacementFieldOutput();
284  OutputImageType* rightDFPtr = this->GetRightDisplacementFieldOutput();
285 
286  // Prevent from streaming
287  leftDFPtr->SetRequestedRegionToLargestPossibleRegion();
288  rightDFPtr->SetRequestedRegionToLargestPossibleRegion();
289 }
290 
291 template <class TInputImage, class TOutputImage>
293 {
294  // Allocate the output
295  this->AllocateOutputs();
296 
297  // Setup the DEM handler if needed
298  auto & demHandler = DEMHandler::GetInstance();
299  auto const& threadDemHandler = demHandler.GetHandlerForCurrentThread();
300 
301  // Set-up a transform to use the DEMHandler
302  typedef otb::GenericRSTransform<> RSTransform2DType;
303  RSTransform2DType::Pointer leftToGroundTransform = RSTransform2DType::New();
304 
305  leftToGroundTransform->SetInputImageMetadata(&(m_LeftImage->GetImageMetadata()));
306 
307  leftToGroundTransform->InstantiateTransform();
308 
309  // Retrieve the output pointers
310  OutputImageType* leftDFPtr = this->GetLeftDisplacementFieldOutput();
311  OutputImageType* rightDFPtr = this->GetRightDisplacementFieldOutput();
312 
313  // Declare all the TDPoint variables we will need
314  TDPointType currentPoint1, currentPoint2, nextLineStart1, nextLineStart2, startLine1, endLine1, startLine2, endLine2, epiPoint1, epiPoint2;
315 
316  // Then, we retrieve the origin of the left input image
317  double localElevation = demHandler.GetDefaultHeightAboveEllipsoid();
318 
319  // Use the mean spacing as before
320  double mean_spacing = 0.5 * (std::abs(m_LeftImage->GetSignedSpacing()[0]) + std::abs(m_LeftImage->GetSignedSpacing()[1]));
321 
322  // Initialize
323  currentPoint1 = m_OutputOriginInLeftImage;
324  if (m_UseDEM)
325  {
326  RSTransform2DType::InputPointType tmpPoint;
327  tmpPoint[0] = currentPoint1[0];
328  tmpPoint[1] = currentPoint1[1];
329  localElevation = demHandler.GetHeightAboveEllipsoid(leftToGroundTransform->TransformPoint(tmpPoint));
330  }
331  currentPoint1[2] = localElevation;
332  currentPoint2 = m_LeftToRightTransform->TransformPoint(currentPoint1);
333  currentPoint2[2] = currentPoint1[2];
334 
335  // These are the points were the next stereo-rectified image line starts
336  nextLineStart1 = currentPoint1;
337  nextLineStart2 = currentPoint2;
338 
339  // We define the iterators we will use
340  typedef itk::ImageRegionIteratorWithIndex<OutputImageType> IteratorType;
341 
342  IteratorType it1(leftDFPtr, leftDFPtr->GetLargestPossibleRegion());
343  IteratorType it2(rightDFPtr, rightDFPtr->GetLargestPossibleRegion());
344 
345  it1.GoToBegin();
346  it2.GoToBegin();
347 
348  // Reset the mean baseline ratio
349  m_MeanBaselineRatio = 0;
350 
351  // Set-up progress reporting
352  itk::ProgressReporter progress(this, 0, leftDFPtr->GetLargestPossibleRegion().GetNumberOfPixels());
353 
354 
355  // We loop on the deformation fields
356  while (!it1.IsAtEnd() && !it2.IsAtEnd())
357  {
358  // 1 - We start by handling the special case where we are at beginning
359  // of a new line
360  if (it1.GetIndex()[0] == 0 || it2.GetIndex()[0] == 0)
361  {
362  // In which case we reset the current points
363  currentPoint1 = nextLineStart1;
364  currentPoint2 = nextLineStart2;
365  }
366 
367  // 2 - Next, we will fill the deformation fields
368  typename OutputImageType::PixelType dFValue1 = it1.Get();
369  typename OutputImageType::PixelType dFValue2 = it2.Get();
370 
371  // We must cast iterators position to physical space
372  PointType currentDFPoint1, currentDFPoint2;
373  leftDFPtr->TransformIndexToPhysicalPoint(it1.GetIndex(), currentDFPoint1);
374  rightDFPtr->TransformIndexToPhysicalPoint(it2.GetIndex(), currentDFPoint2);
375 
376  // Now we can compute the shifts
377  dFValue1[0] = currentPoint1[0] - currentDFPoint1[0];
378  dFValue1[1] = currentPoint1[1] - currentDFPoint1[1];
379  dFValue2[0] = currentPoint2[0] - currentDFPoint2[0];
380  dFValue2[1] = currentPoint2[1] - currentDFPoint2[1];
381 
382  // And set the values
383  it1.Set(dFValue1);
384  it2.Set(dFValue2);
385 
386  // 3 - Next, we will compute epipolar lines direction in both
387  // images
388  double a1;
389 
390  // First, for image 1
391 
392  // This point is the image of the left input image origin at the
393  // average elevation
394  epiPoint2 = m_LeftToRightTransform->TransformPoint(currentPoint1);
395 
396  // The beginning of the epipolar line in the left image is the image
397  // of epiPoint2 at a lower elevation (using the offset)
398  epiPoint2[2] = localElevation - m_ElevationOffset;
399  startLine1 = m_RightToLeftTransform->TransformPoint(epiPoint2);
400 
401  // The endning of the epipolar line in the left image is the image
402  // of epiPoint2 at a higher elevation (using the offset)
403  epiPoint2[2] = localElevation + m_ElevationOffset;
404  endLine1 = m_RightToLeftTransform->TransformPoint(epiPoint2);
405 
406  // Estimate the local baseline ratio
407  double localBaselineRatio =
408  std::sqrt((endLine1[0] - startLine1[0]) * (endLine1[0] - startLine1[0]) + (endLine1[1] - startLine1[1]) * (endLine1[1] - startLine1[1])) /
409  (2 * m_ElevationOffset);
410 
411  m_MeanBaselineRatio += localBaselineRatio;
412 
413  // Now, we can compute the equation of the epipolar line y = a*x+b
414  // (compute angle in physical space)
415  double alpha1 = 0;
416  if (endLine1[0] == startLine1[0])
417  {
418  if (endLine1[1] > startLine1[1])
419  {
420  alpha1 = 0.5 * otb::CONST_PI;
421  }
422  else
423  {
424  alpha1 = -0.5 * otb::CONST_PI;
425  }
426  }
427  else
428  {
429  a1 = (endLine1[1] - startLine1[1]) / (endLine1[0] - startLine1[0]);
430  if (endLine1[0] > startLine1[0])
431  {
432  alpha1 = std::atan(a1);
433  }
434  else
435  {
436  alpha1 = otb::CONST_PI + std::atan(a1);
437  }
438  }
439 
440  // We do the same for image 2
441  currentPoint2[2] = localElevation;
442  epiPoint1 = m_RightToLeftTransform->TransformPoint(currentPoint2);
443 
444  epiPoint1[2] = localElevation - m_ElevationOffset;
445  startLine2 = m_LeftToRightTransform->TransformPoint(epiPoint1);
446 
447  epiPoint1[2] = localElevation + m_ElevationOffset;
448  endLine2 = m_LeftToRightTransform->TransformPoint(epiPoint1);
449 
450  // 4 - Determine position of next points
451 
452  // We want to move m_Scale pixels away in the epipolar line of the
453  // first image
454  // Take into account height direction
455  // double alpha1 = otb::CONST_PI - std::atan(a1);
456  double deltax1 = m_Scale * m_GridStep * mean_spacing * std::cos(alpha1);
457  double deltay1 = m_Scale * m_GridStep * mean_spacing * std::sin(alpha1);
458 
459  // Now we move currentPoint1
460  currentPoint1[0] += deltax1;
461  currentPoint1[1] += deltay1;
462  if (m_UseDEM)
463  {
464  RSTransform2DType::InputPointType tmpPoint;
465  tmpPoint[0] = currentPoint1[0];
466  tmpPoint[1] = currentPoint1[1];
467  localElevation = GetHeightAboveEllipsoid(threadDemHandler, leftToGroundTransform->TransformPoint(tmpPoint));
468  }
469  currentPoint1[2] = localElevation;
470 
471  // And we compute the equivalent displacement in right image
472  currentPoint2 = m_LeftToRightTransform->TransformPoint(currentPoint1);
473 
474  // 5 - Finally, we have to handle a special case for beginning of
475  // line, since at this position we are able to compute the
476  // position of the beginning of next line
477  if (it1.GetIndex()[0] == 0 || it2.GetIndex()[0] == 0)
478  {
479  // We want to move 1 pixel away in the direction orthogonal to
480  // epipolar line
481  double nextdeltax1 = -m_Scale * mean_spacing * m_GridStep * std::sin(alpha1);
482  double nextdeltay1 = m_Scale * mean_spacing * m_GridStep * std::cos(alpha1);
483 
484  // We can then update nextLineStart1
485  nextLineStart1[0] = currentPoint1[0] - deltax1 + nextdeltax1;
486  nextLineStart1[1] = currentPoint1[1] - deltay1 + nextdeltay1;
487  nextLineStart1[2] = localElevation;
488  if (m_UseDEM)
489  {
490  RSTransform2DType::InputPointType tmpPoint;
491  tmpPoint[0] = nextLineStart1[0];
492  tmpPoint[1] = nextLineStart1[1];
493  nextLineStart1[2] = GetHeightAboveEllipsoid(threadDemHandler, leftToGroundTransform->TransformPoint(tmpPoint));
494  }
495 
496 
497  // By construction, nextLineStart2 is always the image of
498  // nextLineStart1 by the left to right transform at the m_AverageElevation
499  nextLineStart2 = m_LeftToRightTransform->TransformPoint(nextLineStart1);
500  }
501 
502  // Last, we move forward
503  ++it1;
504  ++it2;
505 
506  // Update progress
507  progress.CompletedPixel();
508  }
509 
510  // Compute the mean baseline ratio
511  m_MeanBaselineRatio /= leftDFPtr->GetBufferedRegion().GetNumberOfPixels();
512 }
513 
514 template <class TInputImage, class TOutputImage>
516 {
517  // Call superclass implementation
518  Superclass::PrintSelf(os, indent);
519 }
520 
521 } // end namespace otb
522 
523 #endif
otb::DEMHandler::GetInstance
static DEMHandler & GetInstance()
otb::CONST_PI
constexpr double CONST_PI
Definition: otbMath.h:49
otb::StereorectificationDisplacementFieldSource::TDPointType
RSTransformType::InputPointType TDPointType
Definition: otbStereorectificationDisplacementFieldSource.h:117
otbVectorImage.h
otb::StereorectificationDisplacementFieldSource::GenerateData
void GenerateData() override
Definition: otbStereorectificationDisplacementFieldSource.hxx:292
otb::StereorectificationDisplacementFieldSource::EnlargeOutputRequestedRegion
void EnlargeOutputRequestedRegion(itk::DataObject *) override
Definition: otbStereorectificationDisplacementFieldSource.hxx:277
otbDEMHandler.h
otbMath.h
otb
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Definition: otbJoinContainer.h:32
otb::StereorectificationDisplacementFieldSource::m_RightToLeftTransform
RSTransformPointerType m_RightToLeftTransform
Definition: otbStereorectificationDisplacementFieldSource.h:214
otb::StereorectificationDisplacementFieldSource::SizeType
OutputImageType::SizeType SizeType
Definition: otbStereorectificationDisplacementFieldSource.h:106
otb::StereorectificationDisplacementFieldSource::StereorectificationDisplacementFieldSource
StereorectificationDisplacementFieldSource(void)
Definition: otbStereorectificationDisplacementFieldSource.hxx:37
otb::StereorectificationDisplacementFieldSource::RegionType
OutputImageType::RegionType RegionType
Definition: otbStereorectificationDisplacementFieldSource.h:109
otb::StereorectificationDisplacementFieldSource::PrintSelf
void PrintSelf(std::ostream &os, itk::Indent indent) const override
Definition: otbStereorectificationDisplacementFieldSource.hxx:515
otbStereorectificationDisplacementFieldSource.h
otb::StereorectificationDisplacementFieldSource::GetLeftDisplacementFieldOutput
const OutputImageType * GetLeftDisplacementFieldOutput() const
Definition: otbStereorectificationDisplacementFieldSource.hxx:64
otb::StereorectificationDisplacementFieldSource::GenerateOutputInformation
void GenerateOutputInformation() override
Definition: otbStereorectificationDisplacementFieldSource.hxx:107
otb::GenericRSTransform
This is the class to handle generic remote sensing transform.
Definition: otbGenericRSTransform.h:57
otb::StereorectificationDisplacementFieldSource::SpacingType
OutputImageType::SpacingType SpacingType
Definition: otbStereorectificationDisplacementFieldSource.h:108
otb::StereorectificationDisplacementFieldSource::m_LeftToRightTransform
RSTransformPointerType m_LeftToRightTransform
Definition: otbStereorectificationDisplacementFieldSource.h:211
otb::GenericRSTransform::New
static Pointer New()
otb::StereorectificationDisplacementFieldSource::GetRightDisplacementFieldOutput
const OutputImageType * GetRightDisplacementFieldOutput() const
Definition: otbStereorectificationDisplacementFieldSource.hxx:86
otb::GetHeightAboveEllipsoid
double GetHeightAboveEllipsoid(DEMHandlerTLS const &, double lon, double lat)
otb::StereorectificationDisplacementFieldSource::PointType
OutputImageType::PointType PointType
Definition: otbStereorectificationDisplacementFieldSource.h:107
otb::StereorectificationDisplacementFieldSource::OutputImageType
TOutputImage OutputImageType
Definition: otbStereorectificationDisplacementFieldSource.h:102