21 #ifndef otbStereoSensorModelToElevationMapFilter_hxx
22 #define otbStereoSensorModelToElevationMapFilter_hxx
26 #include "itkImageRegionIteratorWithIndex.h"
27 #include "itkProgressReporter.h"
28 #include "itkLinearInterpolateImageFunction.h"
30 #include "itkConstNeighborhoodIterator.h"
35 template <
class TInputImage,
class TOutputHeight>
39 this->SetNumberOfRequiredInputs(2);
40 this->SetNumberOfRequiredOutputs(2);
41 this->SetNthOutput(1, OutputImageType::New());
44 m_Interpolator = itk::LinearInterpolateImageFunction<TInputImage>::New();
50 m_LowerElevation = -20;
51 m_HigherElevation = 20;
53 m_CorrelationThreshold = 0.7;
54 m_VarianceThreshold = 4;
55 m_SubtractInitialElevation =
false;
58 template <
class TInputImage,
class TOutputHeight>
63 template <
class TInputImage,
class TOutputHeight>
66 this->SetNthInput(0,
const_cast<TInputImage*
>(image));
69 template <
class TInputImage,
class TOutputHeight>
72 this->SetNthInput(1,
const_cast<TInputImage*
>(image));
75 template <
class TInputImage,
class TOutputHeight>
78 if (this->GetNumberOfInputs() < 1)
82 return static_cast<const TInputImage*
>(this->itk::ProcessObject::GetInput(0));
85 template <
class TInputImage,
class TOutputHeight>
88 if (this->GetNumberOfInputs() < 2)
92 return static_cast<const TInputImage*
>(this->itk::ProcessObject::GetInput(1));
95 template <
class TInputImage,
class TOutputHeight>
98 if (this->GetNumberOfOutputs() < 2)
102 return static_cast<TOutputHeight*
>(this->itk::ProcessObject::GetOutput(1));
105 template <
class TInputImage,
class TOutputHeight>
109 Superclass::GenerateInputRequestedRegion();
116 if (!masterPtr || !slavePtr || !outputPtr)
123 typename InputImageType::RegionType masterRequestedRegion, slaveRequestedRegion;
124 masterRequestedRegion = outputPtr->GetRequestedRegion();
127 masterRequestedRegion.PadByRadius(m_Radius);
130 typename InputImageType::IndexType mul, mur, mll, mlr;
131 mul = masterRequestedRegion.GetIndex();
132 mur = masterRequestedRegion.GetIndex();
133 mur[0] += masterRequestedRegion.GetSize()[0] - 1;
134 mll = masterRequestedRegion.GetIndex();
135 mll[1] += masterRequestedRegion.GetSize()[1] - 1;
136 mlr = masterRequestedRegion.GetIndex();
137 mlr[0] += masterRequestedRegion.GetSize()[0] - 1;
138 mlr[1] += masterRequestedRegion.GetSize()[1] - 1;
141 typename InputImageType::PointType mpul, mpur, mpll, mplr;
142 masterPtr->TransformIndexToPhysicalPoint(mul, mpul);
143 masterPtr->TransformIndexToPhysicalPoint(mur, mpur);
144 masterPtr->TransformIndexToPhysicalPoint(mll, mpll);
145 masterPtr->TransformIndexToPhysicalPoint(mlr, mplr);
149 transform->SetInputImageMetadata(&(masterPtr->GetImageMetadata()));
150 transform->SetOutputImageMetadata(&(slavePtr->GetImageMetadata()));
152 transform->InstantiateTransform();
158 typename InputImageType::PointType sMinPul, sMinPur, sMinPll, sMinPlr, sMaxPul, sMaxPur, sMaxPll, sMaxPlr;
161 params[0] = m_LowerElevation;
162 transform->SetParameters(params);
164 sMinPul = transform->TransformPoint(mpul);
165 sMinPur = transform->TransformPoint(mpur);
166 sMinPll = transform->TransformPoint(mpll);
167 sMinPlr = transform->TransformPoint(mplr);
170 params[0] = m_HigherElevation;
171 transform->SetParameters(params);
173 sMaxPul = transform->TransformPoint(mpul);
174 sMaxPur = transform->TransformPoint(mpur);
175 sMaxPll = transform->TransformPoint(mpll);
176 sMaxPlr = transform->TransformPoint(mplr);
179 typename InputImageType::IndexType sMinul, sMinur, sMinll, sMinlr, sMaxul, sMaxur, sMaxll, sMaxlr;
181 slavePtr->TransformPhysicalPointToIndex(sMinPul, sMinul);
182 slavePtr->TransformPhysicalPointToIndex(sMaxPul, sMaxul);
183 slavePtr->TransformPhysicalPointToIndex(sMinPur, sMinur);
184 slavePtr->TransformPhysicalPointToIndex(sMaxPur, sMaxur);
185 slavePtr->TransformPhysicalPointToIndex(sMinPll, sMinll);
186 slavePtr->TransformPhysicalPointToIndex(sMaxPll, sMaxll);
187 slavePtr->TransformPhysicalPointToIndex(sMinPlr, sMinlr);
188 slavePtr->TransformPhysicalPointToIndex(sMaxPlr, sMaxlr);
191 typename InputImageType::IndexType sul, slr;
193 sul[0] = std::min(std::min(std::min(sMinul[0], sMaxul[0]), std::min(sMinur[0], sMaxur[0])),
194 std::min(std::min(sMinll[0], sMaxll[0]), std::min(sMinlr[0], sMaxlr[0])));
195 sul[1] = std::min(std::min(std::min(sMinul[1], sMaxul[1]), std::min(sMinur[1], sMaxur[1])),
196 std::min(std::min(sMinll[1], sMaxll[1]), std::min(sMinlr[1], sMaxlr[1])));
197 slr[0] = std::max(std::max(std::max(sMinul[0], sMaxul[0]), std::max(sMinur[0], sMaxur[0])),
198 std::max(std::max(sMinll[0], sMaxll[0]), std::max(sMinlr[0], sMaxlr[0])));
199 slr[1] = std::max(std::max(std::max(sMinul[1], sMaxul[1]), std::max(sMinur[1], sMaxur[1])),
200 std::max(std::max(sMinll[1], sMaxll[1]), std::max(sMinlr[1], sMaxlr[1])));
203 slaveRequestedRegion.SetIndex(sul);
204 typename InputImageType::SizeType ssize;
205 ssize[0] = slr[0] - sul[0] + 1;
206 ssize[1] = slr[1] - sul[1] + 1;
207 slaveRequestedRegion.SetSize(ssize);
210 if (masterRequestedRegion.Crop(masterPtr->GetLargestPossibleRegion()))
212 masterPtr->SetRequestedRegion(masterRequestedRegion);
219 masterPtr->SetRequestedRegion(masterRequestedRegion);
222 itk::InvalidRequestedRegionError e(__FILE__, __LINE__);
223 std::ostringstream msg;
224 msg << this->GetNameOfClass() <<
"::GenerateInputRequestedRegion()";
225 e.SetLocation(msg.str());
226 e.SetDescription(
"Requested region is (at least partially) outside the largest possible region of image 1.");
227 e.SetDataObject(masterPtr);
232 if (slaveRequestedRegion.Crop(slavePtr->GetLargestPossibleRegion()))
234 slavePtr->SetRequestedRegion(slaveRequestedRegion);
241 typename InputImageType::SizeType slaveSize;
243 slaveRequestedRegion.SetSize(slaveSize);
244 typename InputImageType::IndexType slaveIndex;
246 slaveRequestedRegion.SetIndex(slaveIndex);
249 slavePtr->SetRequestedRegion(slaveRequestedRegion);
254 template <
class TInputImage,
class TOutputHeight>
258 m_Interpolator->SetInputImage(this->GetSlaveInput());
259 this->GetCorrelationOutput()->FillBuffer(0.);
269 rsTransform->SetInputImageMetadata(&(outputPtr->GetImageMetadata()));
270 rsTransform->InstantiateTransform();
273 itk::ImageRegionIteratorWithIndex<OutputImageType> outputIt(outputPtr, outputPtr->GetBufferedRegion());
275 for (outputIt.GoToBegin(); !outputIt.IsAtEnd(); ++outputIt)
279 outputPtr->TransformIndexToPhysicalPoint(outputIt.GetIndex(), outputPoint);
282 geoPoint = rsTransform->TransformPoint(outputPoint);
290 m_MasterToSlave = GenericRSTransform3DType::New();
291 m_MasterToSlave->SetInputImageMetadata(&(masterPtr->GetImageMetadata()));
292 m_MasterToSlave->SetOutputImageMetadata(&(slavePtr->GetImageMetadata()));
293 m_MasterToSlave->InstantiateTransform();
296 template <
class TInputImage,
class TOutputHeight>
298 itk::ThreadIdType threadId)
356 itk::ProgressReporter progress(
this, threadId, outputRegionForThread.GetNumberOfPixels());
359 itk::ConstNeighborhoodIterator<InputImageType> inputIt(m_Radius, masterPtr, outputRegionForThread);
360 itk::ImageRegionIterator<OutputImageType> outputIt(outputPtr, outputRegionForThread);
361 itk::ImageRegionIteratorWithIndex<OutputImageType> correlIt(correlPtr, outputRegionForThread);
364 outputIt.GoToBegin();
365 correlIt.GoToBegin();
370 std::vector<double> master;
371 master.reserve(inputIt.Size());
372 master = std::vector<double>(inputIt.Size(), 0);
375 std::vector<double> slave;
376 slave.reserve(inputIt.Size());
377 slave = std::vector<double>(inputIt.Size(), 0);
381 while (!outputIt.IsAtEnd() && !inputIt.IsAtEnd())
384 typename InputImageType::PointType inPoint, outPoint, currentPoint, optimalPoint;
386 typename InputImageType::IndexType index;
389 double initHeight = outputIt.Get();
390 double optimalHeight = initHeight;
391 double optimalCorrelation = 0;
394 if (initHeight != -32768)
397 double masterSum = 0;
398 double masterVariance = 0;
401 for (
unsigned int i = 0; i < inputIt.Size(); ++i)
404 double value = inputIt.GetPixel(i);
412 masterSum /= inputIt.Size();
415 for (
unsigned int i = 0; i < inputIt.Size(); ++i)
417 masterVariance += (master[i] - masterSum) * (master[i] - masterSum);
419 masterVariance /= (inputIt.Size() - 1);
422 if (masterVariance > m_VarianceThreshold)
425 for (
double height = initHeight + m_LowerElevation; height < initHeight + m_HigherElevation; height += m_ElevationStep)
429 for (
unsigned int i = 0; i < inputIt.Size(); ++i)
432 index = inputIt.GetIndex(i);
435 masterPtr->TransformIndexToPhysicalPoint(index, inPoint);
436 in3DPoint[0] = inPoint[0];
437 in3DPoint[1] = inPoint[1];
438 in3DPoint[2] = height;
441 out3DPoint = m_MasterToSlave->TransformPoint(in3DPoint);
442 outPoint[0] = out3DPoint[0];
443 outPoint[1] = out3DPoint[1];
446 if (m_Interpolator->IsInsideBuffer(outPoint))
449 slave[i] = m_Interpolator->Evaluate(outPoint);
460 double correlationValue = this->Correlation(master, slave);
463 if (correlationValue > optimalCorrelation)
466 optimalCorrelation = correlationValue;
467 optimalHeight = height;
475 double finalOffset = 0.;
478 if (m_SubtractInitialElevation)
480 finalOffset = initHeight;
483 if (optimalCorrelation > m_CorrelationThreshold)
485 outputIt.Set(optimalHeight - finalOffset);
486 correlIt.Set(optimalCorrelation);
490 outputIt.Set(initHeight - finalOffset);
495 progress.CompletedPixel();
504 template <
class TInputImage,
class TOutputHeight>
507 double meanSlave = 0;
508 double meanMaster = 0;
509 double correlationValue = 0;
512 for (
unsigned int i = 0; i < master.size(); ++i)
514 meanSlave += slave[i];
515 meanMaster += master[i];
517 meanSlave /= slave.size();
518 meanMaster /= master.size();
520 double crossProd = 0.;
521 double squareSumMaster = 0.;
522 double squareSumSlave = 0.;
526 for (
unsigned int i = 0; i < master.size(); ++i)
528 crossProd += (slave[i] - meanSlave) * (master[i] - meanMaster);
529 squareSumSlave += (slave[i] - meanSlave) * (slave[i] - meanSlave);
530 squareSumMaster += (master[i] - meanMaster) * (master[i] - meanMaster);
533 correlationValue = std::abs(crossProd / (std::sqrt(squareSumSlave) * std::sqrt(squareSumMaster)));
535 return correlationValue;