21 #ifndef otbMultiDisparityMapTo3DFilter_hxx
22 #define otbMultiDisparityMapTo3DFilter_hxx
25 #include "itkImageRegionConstIteratorWithIndex.h"
26 #include "itkImageRegionIterator.h"
31 template <
class TDisparityImage,
class TOutputImage,
class TMaskImage,
class TRes
idueImage>
35 this->SetNumberOfRequiredInputs(3);
36 this->SetNumberOfRequiredInputs(1);
37 this->m_MovingImageMetadatas.resize(1);
40 this->SetNumberOfRequiredOutputs(2);
41 this->SetNthOutput(0, TOutputImage::New());
42 this->SetNthOutput(1, TResidueImage::New());
45 template <
class TDisparityImage,
class TOutputImage,
class TMaskImage,
class TRes
idueImage>
50 template <
class TDisparityImage,
class TOutputImage,
class TMaskImage,
class TRes
idueImage>
55 this->SetNumberOfRequiredInputs(3 * nb);
56 this->m_MovingImageMetadatas.resize(nb);
60 template <
class TDisparityImage,
class TOutputImage,
class TMaskImage,
class TRes
idueImage>
63 return this->m_MovingImageMetadatas.size();
66 template <
class TDisparityImage,
class TOutputImage,
class TMaskImage,
class TRes
idueImage>
68 const TDisparityImage* hmap)
70 if ((3 * (index + 1)) > this->GetNumberOfRequiredInputs())
72 itkExceptionMacro(<<
"Index is greater than the number of moving images");
75 this->SetNthInput(3 * index,
const_cast<TDisparityImage*
>(hmap));
78 template <
class TDisparityImage,
class TOutputImage,
class TMaskImage,
class TRes
idueImage>
80 const TDisparityImage* vmap)
82 if ((3 * (index + 1)) > this->GetNumberOfRequiredInputs())
84 itkExceptionMacro(<<
"Index is greater than the number of moving images");
87 this->SetNthInput(3 * index + 1,
const_cast<TDisparityImage*
>(vmap));
90 template <
class TDisparityImage,
class TOutputImage,
class TMaskImage,
class TRes
idueImage>
93 if ((3 * (index + 1)) > this->GetNumberOfRequiredInputs())
95 itkExceptionMacro(<<
"Index is greater than the number of moving images");
98 this->SetNthInput(3 * index + 2,
const_cast<TMaskImage*
>(mask));
101 template <
class TDisparityImage,
class TOutputImage,
class TMaskImage,
class TRes
idueImage>
102 const TDisparityImage*
105 if ((3 * (index + 1)) > this->GetNumberOfRequiredInputs())
109 return static_cast<const TDisparityImage*
>(this->itk::ProcessObject::GetInput(3 * index));
112 template <
class TDisparityImage,
class TOutputImage,
class TMaskImage,
class TRes
idueImage>
113 const TDisparityImage*
116 if ((3 * (index + 1)) > this->GetNumberOfRequiredInputs())
120 return static_cast<const TDisparityImage*
>(this->itk::ProcessObject::GetInput(3 * index + 1));
123 template <
class TDisparityImage,
class TOutputImage,
class TMaskImage,
class TRes
idueImage>
126 if ((3 * (index + 1)) > this->GetNumberOfRequiredInputs())
130 return static_cast<const TMaskImage*
>(this->itk::ProcessObject::GetInput(3 * index + 2));
133 template <
class TDisparityImage,
class TOutputImage,
class TMaskImage,
class TRes
idueImage>
136 return static_cast<const TResidueImage*
>(this->itk::ProcessObject::GetOutput(1));
139 template <
class TDisparityImage,
class TOutputImage,
class TMaskImage,
class TRes
idueImage>
142 return static_cast<TResidueImage*
>(this->itk::ProcessObject::GetOutput(1));
145 template <
class TDisparityImage,
class TOutputImage,
class TMaskImage,
class TRes
idueImage>
149 if (this->m_MovingImageMetadatas.size() > index)
151 this->m_MovingImageMetadatas[index] = imd;
155 template <
class TDisparityImage,
class TOutputImage,
class TMaskImage,
class TRes
idueImage>
159 if (this->m_MovingImageMetadatas.size() <= index)
161 itkExceptionMacro(<<
"ImageMetadata index is outside the container");
163 return this->m_MovingImageMetadatas[index];
166 template <
class TDisparityImage,
class TOutputImage,
class TMaskImage,
class TRes
idueImage>
169 const TDisparityImage* horizDisp = this->GetHorizontalDisparityMapInput(0);
170 TOutputImage* outputPtr = this->GetOutput();
171 TResidueImage* residuePtr = this->GetResidueOutput();
175 outputPtr->SetLargestPossibleRegion(horizDisp->GetLargestPossibleRegion());
176 outputPtr->SetNumberOfComponentsPerPixel(3);
178 residuePtr->SetLargestPossibleRegion(horizDisp->GetLargestPossibleRegion());
179 residuePtr->SetNumberOfComponentsPerPixel(1);
182 outputPtr->SetOrigin(horizDisp->GetOrigin());
183 outputPtr->SetSignedSpacing(horizDisp->GetSignedSpacing());
185 residuePtr->SetOrigin(horizDisp->GetOrigin());
186 residuePtr->SetSignedSpacing(horizDisp->GetSignedSpacing());
188 if (this->m_ReferenceImageMetadata)
190 auto outputmetadata = *m_ReferenceImageMetadata;
192 outputmetadata.Bands.clear();
193 outputPtr->SetImageMetadata(outputmetadata);
194 residuePtr->SetImageMetadata(outputmetadata);
198 itkExceptionMacro(<<
"Reference ImageMetadata is missing");
203 itkExceptionMacro(<<
"First horizontal disparity map is missing");
207 template <
class TDisparityImage,
class TOutputImage,
class TMaskImage,
class TRes
idueImage>
210 RegionType requested = this->GetOutput()->GetRequestedRegion();
211 RegionType largest = this->GetHorizontalDisparityMapInput(0)->GetLargestPossibleRegion();
213 for (
unsigned int i = 0; i < this->GetNumberOfRequiredInputs(); ++i)
215 unsigned int index = i % 3;
216 unsigned int pos = i / 3;
222 TDisparityImage* horizDisp =
const_cast<TDisparityImage*
>(this->GetHorizontalDisparityMapInput(pos));
225 if (horizDisp->GetLargestPossibleRegion() != largest)
227 itkExceptionMacro(<<
"Horizontal disparity map at position " << pos <<
" has a different largest region");
229 horizDisp->SetRequestedRegion(requested);
233 itkExceptionMacro(<<
"Horizontal disparity map at position " << pos <<
" is missing");
240 TDisparityImage* vertiDisp =
const_cast<TDisparityImage*
>(this->GetVerticalDisparityMapInput(pos));
243 if (vertiDisp->GetLargestPossibleRegion() != largest)
245 itkExceptionMacro(<<
"Vertical disparity map at position " << pos <<
" has a different largest region");
247 vertiDisp->SetRequestedRegion(requested);
254 TMaskImage* maskDisp =
const_cast<TMaskImage*
>(this->GetDisparityMaskInput(pos));
257 if (maskDisp->GetLargestPossibleRegion() != largest)
259 itkExceptionMacro(<<
"Disparity mask at position " << pos <<
" has a different largest region");
261 maskDisp->SetRequestedRegion(requested);
266 itkExceptionMacro(<<
"Unexpected value for (i%3) ");
271 for (
unsigned int k = 0; k < this->m_MovingImageMetadatas.size(); ++k)
273 if (!this->m_MovingImageMetadatas[k])
275 itkExceptionMacro(<<
"ImageMetadata of moving image at position " << k <<
" is empty");
280 template <
class TDisparityImage,
class TOutputImage,
class TMaskImage,
class TRes
idueImage>
283 auto referenceToGroundTransform = RSTransformType::New();
285 referenceToGroundTransform->SetInputImageMetadata(this->m_ReferenceImageMetadata);
286 referenceToGroundTransform->InstantiateTransform();
289 std::vector<RSTransformType::Pointer> movingToGroundTransform;
291 for (
unsigned int k = 0; k < this->m_MovingImageMetadatas.size(); ++k)
294 transfo->SetInputImageMetadata(this->m_MovingImageMetadatas[k]);
295 transfo->InstantiateTransform();
296 movingToGroundTransform.push_back(transfo);
300 TOutputImage* outputPtr = this->GetOutput();
301 TResidueImage* residuePtr = this->GetResidueOutput();
303 itk::ImageRegionIteratorWithIndex<OutputImageType> outIt(outputPtr, outputRegionForThread);
304 itk::ImageRegionIterator<ResidueImageType> resIt(residuePtr, outputRegionForThread);
312 for (
unsigned int k = 0; k < this->m_MovingImageMetadatas.size(); ++k)
315 hDispIts[k] = itk::ImageRegionConstIterator<DisparityMapType>(this->GetHorizontalDisparityMapInput(k), outputRegionForThread);
316 hDispIts[k].GoToBegin();
318 if (this->GetVerticalDisparityMapInput(k))
320 vDispIts[k] = itk::ImageRegionConstIterator<DisparityMapType>(this->GetVerticalDisparityMapInput(k), outputRegionForThread);
321 vDispIts[k].GoToBegin();
324 if (this->GetDisparityMaskInput(k))
326 maskIts[k] = itk::ImageRegionConstIterator<MaskImageType>(this->GetDisparityMaskInput(k), outputRegionForThread);
327 maskIts[k].GoToBegin();
337 typename OutputImageType::PointType pointRef;
340 typename OutputImageType::PixelType outPixel(3);
343 typename PointSetType::Pointer pointSetA = PointSetType::New();
344 typename PointSetType::Pointer pointSetB = PointSetType::New();
346 while (!outIt.IsAtEnd())
348 pointSetA->Initialize();
349 pointSetB->Initialize();
354 outputPtr->TransformIndexToPhysicalPoint(outIt.GetIndex(), pointRef);
356 currentPoint[0] = pointRef[0];
357 currentPoint[1] = pointRef[1];
358 currentPoint[2] = altiMax;
360 pointA = referenceToGroundTransform->TransformPoint(currentPoint);
362 currentPoint[2] = altiMin;
363 pointB = referenceToGroundTransform->TransformPoint(currentPoint);
365 pointSetA->SetPoint(0, pointA);
366 pointSetB->SetPoint(0, pointB);
367 pointSetA->SetPointData(0, 0);
368 pointSetB->SetPointData(0, 0);
370 unsigned int nbPoints = 1;
372 for (
unsigned int k = 0; k < this->m_MovingImageMetadatas.size(); ++k)
377 if (maskIts.count(k) && !(maskIts[k].Get() > 0))
382 currentPoint[0] = pointRef[0] + hDispIts[k].Get();
383 currentPoint[1] = pointRef[1];
384 if (vDispIts.count(k))
386 currentPoint[1] += vDispIts[k].Get();
389 currentPoint[2] = altiMax;
390 pointAi = movingToGroundTransform[k]->TransformPoint(currentPoint);
391 currentPoint[2] = altiMin;
392 pointBi = movingToGroundTransform[k]->TransformPoint(currentPoint);
393 pointSetA->SetPoint(nbPoints, pointAi);
394 pointSetB->SetPoint(nbPoints, pointBi);
395 pointSetA->SetPointData(nbPoints, k + 1);
396 pointSetB->SetPointData(nbPoints, k + 1);
403 TDPointType intersection = optimizer->Compute(pointSetA, pointSetB);
404 outPixel[0] = intersection[0];
405 outPixel[1] = intersection[1];
406 outPixel[2] = intersection[2];
407 globalResidue = optimizer->GetGlobalResidue();
416 resIt.Set(globalResidue);
419 for (
typename DispMapIteratorList::iterator hIt = hDispIts.begin(); hIt != hDispIts.end(); ++hIt)
423 for (
typename DispMapIteratorList::iterator vIt = vDispIts.begin(); vIt != vDispIts.end(); ++vIt)
427 for (
typename MaskIteratorList::iterator mIt = maskIts.begin(); mIt != maskIts.end(); ++mIt)
itk::SmartPointer< Self > Pointer
void SetNumberOfMovingImages(unsigned int nb)
void SetDisparityMaskInput(unsigned int index, const TMaskImage *mask)
void GenerateInputRequestedRegion() override
std::map< unsigned int, itk::ImageRegionConstIterator< DisparityMapType > > DispMapIteratorList
~MultiDisparityMapTo3DFilter() override
const TDisparityImage * GetVerticalDisparityMapInput(unsigned int index) const
void SetMovingImageMetadata(unsigned int index, const ImageMetadata *imd)
void SetVerticalDisparityMapInput(unsigned int index, const TDisparityImage *vmap)
unsigned int GetNumberOfMovingImages()
const TMaskImage * GetDisparityMaskInput(unsigned int index) const
const TResidueImage * GetResidueOutput() const
MultiDisparityMapTo3DFilter()
OutputImageType::RegionType RegionType
void GenerateOutputInformation() override
void SetHorizontalDisparityMapInput(unsigned int index, const TDisparityImage *hmap)
const ImageMetadata * GetMovingImageMetadata(unsigned int index) const
std::map< unsigned int, itk::ImageRegionConstIterator< MaskImageType > > MaskIteratorList
void DynamicThreadedGenerateData(const RegionType &outputRegionForThread) override
const TDisparityImage * GetHorizontalDisparityMapInput(unsigned int index) const
RSTransformType::InputPointType TDPointType
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.