Skip to content
Snippets Groups Projects
Commit ec2af780 authored by Romain Garrigues's avatar Romain Garrigues
Browse files

Ajout de otb::StreamingResampleImageFilter, basé sur itk::ResampleImageFilter

parent 778816eb
No related branches found
No related tags found
No related merge requests found
/*=========================================================================
Program: ORFEO Toolbox
Language: C++
Date: $Date$
Version: $Revision$
Copyright (c) Centre National d'Etudes Spatiales. All rights reserved.
See OTBCopyright.txt for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notices for more information.
=========================================================================*/
#ifndef __otbStreamingResampleImageFilter_h
#define __otbStreamingResampleImageFilter_h
#include "itkResampleImageFilter.h"
#include "otbStreamingTraits.h"
namespace otb
{
/** \class StreamingResampleImageFilter
* \brief Resample image filter
*
* This class add streaming aspect on it::ResampleImageFilter
*/
template <class TInputImage, class TOutputImage, class TInterpolatorPrecisionType>
class ITK_EXPORT StreamingResampleImageFilter : public itk::ResampleImageFilter< TInputImage, TOutputImage, TInterpolatorPrecisionType >
{
public:
/** Standard class typedefs. */
typedef StreamingResampleImageFilter Self;
typedef itk::ResampleImageFilter< TInputImage, TOutputImage, TInterpolatorPrecisionType> Superclass;
typedef itk::SmartPointer<Self> Pointer;
typedef itk::SmartPointer<const Self> ConstPointer;
/** Image related typedefs. */
typedef typename TInputImage::Pointer InputImagePointer;
typedef typename TOutputImage::Pointer OutputImagePointer;
typedef typename TInputImage::IndexType IndexType;
typedef typename TInputImage::SizeType SizeType;
/** Method for creation through the object factory. */
itkNewMacro(Self);
/** Run-time type information (and related methods). */
itkTypeMacro(StreamingResampleImageFilter, itk::ResampleImageFilter);
/** Type definitions */
typedef typename TInputImage::PixelType InputPixelType;
typedef typename TOutputImage::PixelType OutputPixelType;
typedef typename TInputImage::RegionType InputImageRegionType;
typedef typename TOutputImage::RegionType OutputImageRegionType;
typedef typename Superclass::InterpolatorType InterpolatorType;
typedef typename InterpolatorType::PointType PointType;
typedef typename TInputImage::SizeType SizeType;
/** Set size of neighborhood needed to interpolate points */
itkSetMacro(InterpolatorNeighborhoodRadius,unsigned int);
/** Get size of neighborhood needed to interpolate points */
itkGetMacro(InterpolatorNeighborhoodRadius,unsigned int);
/** ResampleImageFilter needs a different input requested region than
* the output requested region. As such, ResampleImageFilter needs
* to provide an implementation for GenerateInputRequestedRegion()
* in order to inform the pipeline execution model.
* \sa ProcessObject::GenerateInputRequestedRegion()
*/
virtual void GenerateInputRequestedRegion() ;
protected:
StreamingResampleImageFilter();
virtual ~StreamingResampleImageFilter() {};
void PrintSelf(std::ostream& os, itk::Indent indent) const;
private:
StreamingResampleImageFilter(const Self&); //purposely not implemented
void operator=(const Self&); //purposely not implemented
// Determine size of pad needed for interpolators neighborhood
unsigned int m_InterpolatorNeighborhoodRadius;
// Determine if interpolator radius is determined by class user
// bool m_RadiusIsDeterminedByUser;
};
} // end namespace otb
#ifndef OTB_MANUAL_INSTANTIATION
#include "otbStreamingResampleImageFilter.txx"
#endif
#endif
/*=========================================================================
Program: ORFEO Toolbox
Language: C++
Date: $Date$
Version: $Revision$
Copyright (c) Centre National d'Etudes Spatiales. All rights reserved.
See OTBCopyright.txt for details.
Some parts of this code are derived from ITK. See ITKCopyright.txt
for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notices for more information.
=========================================================================*/
#ifndef _otbStreamingResampleImageFilter_txx
#define _otbStreamingResampleImageFilter_txx
#include "otbStreamingResampleImageFilter.h"
#include "otbStreamingTraits.h"
#include "itkImageRegionIterator.h"
#include "itkImageRegionConstIterator.h"
#include "itkNumericTraits.h"
#include "itkProgressReporter.h"
#include "itkContinuousIndex.h"
namespace otb {
template<class TInputImage, class TOutputImage, class TInterpolatorPrecisionType>
StreamingResampleImageFilter<TInputImage,TOutputImage,TInterpolatorPrecisionType>
::StreamingResampleImageFilter()
{
// Default neighborhood interpolation radius is one pixel
m_InterpolatorNeighborhoodRadius = 1 ;
// m_RadiusIsDeterminedByUser = false;
}
template<class TInputImage, class TOutputImage, class TInterpolatorPrecisionType>
void
StreamingResampleImageFilter<TInputImage,TOutputImage,TInterpolatorPrecisionType>
::GenerateInputRequestedRegion()
{
Superclass::GenerateInputRequestedRegion();
if ( this->GetInput() )
{
otbDebugMacro(<< "-------------- GenerateInputRequestedRegion ---------------" << std::endl);
InputImagePointer inputImage = const_cast< typename Superclass::InputImageType *>( this->GetInput() );
OutputImagePointer outputImage = const_cast< typename Superclass::OutputImageType *>( this->GetOutput() );
IndexType index = outputImage->GetRequestedRegion().GetIndex();
SizeType size = outputImage->GetRequestedRegion().GetSize();
// Obtain coordinates of upperleft, upperright, lowerleft and lowerright points in the image
IndexType indexTmp ;
std::vector<IndexType> vPoints;
typename std::vector<IndexType>::iterator it;
otbDebugMacro(<< "Size : " << size[0] << " " << size[1]);
indexTmp[0]=index[0];
indexTmp[1]=index[1];
vPoints.push_back(indexTmp);
otbDebugMacro(<< "indexUL : (" << indexTmp[0] << "," << indexTmp[1] << ")");
indexTmp[0]=index[0]+size[0];
indexTmp[1]=index[1];
vPoints.push_back(indexTmp);
otbDebugMacro(<< "indexUR : (" << indexTmp[0] << "," << indexTmp[1] << ")");
indexTmp[0]=index[0]+size[0];
indexTmp[1]=index[1]+size[1];
vPoints.push_back(indexTmp);
otbDebugMacro(<< "indexLR : (" << indexTmp[0] << "," << indexTmp[1] << ")");
indexTmp[0]=index[0];
indexTmp[1]=index[1]+size[1];
vPoints.push_back(indexTmp);
otbDebugMacro(<< "indexLL : (" << indexTmp[0] << "," << indexTmp[1] << ")");
typedef itk::ContinuousIndex<TInterpolatorPrecisionType, 2> ContinuousIndexType;
typename ContinuousIndexType::ValueType minX = itk::NumericTraits<typename ContinuousIndexType::ValueType>::max();
typename ContinuousIndexType::ValueType maxX = 0;
typename ContinuousIndexType::ValueType minY = itk::NumericTraits<typename ContinuousIndexType::ValueType>::max();
typename ContinuousIndexType::ValueType maxY = 0;
// Coordinates of current output pixel
PointType outputPoint;
PointType inputPoint;
// Transform each "corner" point
for (it = vPoints.begin(); it != vPoints.end(); it++)
{
ContinuousIndexType indexTmpTr;
// Calculate transformed points needed for previous filter in the pipeline
outputImage->TransformIndexToPhysicalPoint( *it, outputPoint );
// Compute corresponding input pixel continuous index
inputPoint = this->GetTransform()->TransformPoint(outputPoint);
inputImage->TransformPhysicalPointToContinuousIndex(inputPoint, indexTmpTr);
if (indexTmpTr[0]>maxX)
maxX = indexTmpTr[0];
else if (indexTmpTr[0]<minX)
minX = indexTmpTr[0];
if (indexTmpTr[1]>maxY)
maxY = indexTmpTr[1];
else if (indexTmpTr[1]<minY)
minY = indexTmpTr[1];
otbDebugMacro(<< "indexTr : (" << indexTmpTr[0] << "," << indexTmpTr[1] << ")");
}
otbDebugMacro(<< "MinX : " << minX << " MinY : " << minY << " MaxX : " << maxX << " MaxY " << maxY);
// Create region needed in previous filter in the pipeline, which is the bounding box of previous transformed points
InputImageRegionType region;
index[0] = static_cast<long int>(minX);
index[1] = static_cast<long int>(minY);
size[0] = static_cast<long unsigned int>(maxX - minX);
size[1] = static_cast<long unsigned int>(maxY - minY);
otbDebugMacro(<< "Index : (" << index[0] << "," << index[1] << ") Size : (" << size[0] << "," << size[1] << ")");
region.SetSize(size);
region.SetIndex(index);
// Grow region to be sure that interpolator can found needed point on image borders
unsigned int neededRadius = otb::StreamingTraits<typename Superclass::InputImageType>::CalculateNeededRadiusForInterpolator(this->GetInterpolator());
if (neededRadius == 0)
{
otbDebugMacro(<< "Warning : if you haven't fixed interpolator radius, default is 1");
neededRadius = m_InterpolatorNeighborhoodRadius;
}
otbDebugMacro(<< "Interpolation needed radius : " << neededRadius);
region.PadByRadius(neededRadius);
otbDebugMacro(<< "Initial Region : Index(" << inputImage->GetLargestPossibleRegion().GetIndex()[0] << "," << inputImage->GetLargestPossibleRegion().GetIndex()[1] << ") Size(" << inputImage->GetLargestPossibleRegion().GetSize()[0] << "," << inputImage->GetLargestPossibleRegion().GetSize()[1] << ")");
// To be sure that requested region in pipeline is not largest than real input image
otbDebugMacro(<< "Final Region (Before Crop) : Index(" << region.GetIndex()[0] << "," << region.GetIndex()[1] << ") Size(" << region.GetSize()[0] << "," << region.GetSize()[1] << ")");
// If requested region is not contained in input image, then result region is null
if (!region.Crop(inputImage->GetLargestPossibleRegion()))
{
index[0]=0;
index[1]=0;
size[0]=0;
size[1]=0;
region.SetIndex(index);
region.SetSize(size);
}
inputImage->SetRequestedRegion(region);
otbDebugMacro(<< "Final Region (After Crop) : Index(" << region.GetIndex()[0] << "," << region.GetIndex()[1] << ") Size(" << region.GetSize()[0] << "," << region.GetSize()[1] << ")");
}
}
/*template<class TInputImage, class TOutputImage, class TInterpolatorPrecisionType>
void
StreamingResampleImageFilter<TInputImage,TOutputImage,TInterpolatorPrecisionType>
::SetInterpolatorNeighborhoodRadius(unsigned int radius)
{
m_InterpolatorNeighborhoodRadius = radius;
// m_RadiusIsDeterminedByUser = true;
}*/
template<class TInputImage, class TOutputImage, class TInterpolatorPrecisionType>
void
StreamingResampleImageFilter<TInputImage,TOutputImage,TInterpolatorPrecisionType>
::PrintSelf(std::ostream& os, itk::Indent indent) const
{
Superclass::PrintSelf(os,indent);
os << indent << "StreamingResampleImageFilter " << std::endl;
}
}// end namespace otb
#endif
......@@ -21,6 +21,12 @@ PURPOSE. See the above copyright notices for more information.
#include "otbMacro.h"
#include "otbConfigure.h"
#include "itkInterpolateImageFunction.h"
#include "itkBSplineInterpolateImageFunction.h"
#include "itkLinearInterpolateImageFunction.h"
#include "itkNearestNeighborInterpolateImageFunction.h"
//#include "itkWindowedSincInterpolateImageFunction.h"
namespace otb
{
......@@ -50,12 +56,21 @@ class ITK_EXPORT StreamingTraits
typedef typename ImageType::InternalPixelType PixelType;
typedef StreamingMode StreamingModeType;
typedef itk::InterpolateImageFunction<TImage,double> InterpolationType;
typedef itk::BSplineInterpolateImageFunction<TImage,double> BSplineInterpolationType;
typedef itk::LinearInterpolateImageFunction<TImage,double> LinearInterpolationType;
typedef itk::NearestNeighborInterpolateImageFunction<TImage,double> NearestNeighborInterpolationType;
//typedef typename InterpolationType::Pointer InterpolationPointerType;
/**
* This method computes the number of streaming divisions, based on
* the streaming mode and the image attributes. The last three parameters are ignored if
* the mode does not need them.
* \param image The image to stream,
* \param image The image to stream,
* \param region the region to stream,
* \param mode the streaming mode,
* \param numberOfStreamDivision the number of stream division if streaming mode is SET_NUMBER_OF_STREAM_DIVISIONS,
......@@ -69,23 +84,23 @@ class ITK_EXPORT StreamingTraits
unsigned long numberOfStreamDivision,
unsigned long bufferMemorySize,
unsigned long bufferNumberOfLinesDivisions)
{
{
unsigned long numDivisions(0);
switch(mode)
{
case SET_NUMBER_OF_STREAM_DIVISIONS :
{
numDivisions = numberOfStreamDivision;
}
break;
case SET_BUFFER_MEMORY_SIZE :
{
{
case SET_NUMBER_OF_STREAM_DIVISIONS :
{
numDivisions = numberOfStreamDivision;
}
break;
case SET_BUFFER_MEMORY_SIZE :
{
const unsigned long bufferMemorySize = bufferMemorySize/8;
unsigned long numberColumnsOfRegion = region.GetSize()[0]; // X dimension
unsigned long numberColumnsOfRegion = region.GetSize()[0]; // X dimension
const unsigned long sizeLine = numberColumnsOfRegion * \
image->GetNumberOfComponentsPerPixel() * \
sizeof(PixelType);
image->GetNumberOfComponentsPerPixel() * \
sizeof(PixelType);
unsigned long regionSize = region.GetSize()[1] * sizeLine;
otbMsgDevMacro(<<"image->GetNumberOfComponentsPerPixel() = "<<image->GetNumberOfComponentsPerPixel());
otbMsgDevMacro(<<"sizeof(PixelType) = "<<sizeof(PixelType));
......@@ -105,43 +120,43 @@ class ITK_EXPORT StreamingTraits
regionSize = sizeLine;
}
//Calculate NumberOfStreamDivisions
numDivisions = static_cast<unsigned long>(vcl_ceil(static_cast<double>(regionSize)/static_cast<double>(bufferMemorySize)));
numDivisions = static_cast<unsigned long>(vcl_ceil(static_cast<double>(regionSize)/static_cast<double>(bufferMemorySize)));
}
else
{
//Non streaming
numDivisions = 1;
}
}
break;
case SET_BUFFER_NUMBER_OF_LINES :
{
if( bufferNumberOfLinesDivisions < 1 )
{
itkGenericExceptionMacro(<<"Buffer number of lines division must be greater than 0 !");
}
/* Calculate number of split */
unsigned long numberLinesOfRegion = region.GetSize()[1]; // Y dimension
if( numberLinesOfRegion > bufferNumberOfLinesDivisions )
{
numDivisions = static_cast<unsigned long>(vcl_ceil(static_cast<double>(numberLinesOfRegion)/static_cast<double>(bufferNumberOfLinesDivisions)));
}
else
{
}
break;
case SET_BUFFER_NUMBER_OF_LINES :
{
if( bufferNumberOfLinesDivisions < 1 )
{
itkGenericExceptionMacro(<<"Buffer number of lines division must be greater than 0 !");
}
/* Calculate number of split */
unsigned long numberLinesOfRegion = region.GetSize()[1]; // Y dimension
if( numberLinesOfRegion > bufferNumberOfLinesDivisions )
{
numDivisions = static_cast<unsigned long>(vcl_ceil(static_cast<double>(numberLinesOfRegion)/static_cast<double>(bufferNumberOfLinesDivisions)));
}
else
{
//Non streaming
numDivisions = 1;
}
}
break;
case SET_AUTOMATIC_NUMBER_OF_STREAM_DIVISIONS :
{
numDivisions = 1;
}
}
break;
case SET_AUTOMATIC_NUMBER_OF_STREAM_DIVISIONS :
{
const unsigned long streamMaxSizeBufferForStreamingBytes = OTB_STREAM_MAX_SIZE_BUFFER_FOR_STREAMING;
const unsigned long streamImageSizeToActivateStreamingBytes = OTB_STREAM_IMAGE_SIZE_TO_ACTIVATE_STREAMING;
const unsigned long streamImageSizeToActivateStreamingBytes = OTB_STREAM_IMAGE_SIZE_TO_ACTIVATE_STREAMING;
//Convert in octet unit
unsigned long streamMaxSizeBufferForStreaming = streamMaxSizeBufferForStreamingBytes/8;
const unsigned long streamImageSizeToActivateStreaming = streamImageSizeToActivateStreamingBytes/8;
const unsigned long streamImageSizeToActivateStreaming = streamImageSizeToActivateStreamingBytes/8;
unsigned long numberColumnsOfRegion = region.GetSize()[0]; // X dimension
unsigned long numberColumnsOfRegion = region.GetSize()[0]; // X dimension
const unsigned long sizeLine = numberColumnsOfRegion * \
image->GetNumberOfComponentsPerPixel() * \
sizeof(PixelType);
......@@ -167,26 +182,62 @@ class ITK_EXPORT StreamingTraits
}
otbMsgDevMacro(<<"Buffer size : "<<streamMaxSizeBufferForStreaming);
//Calculate NumberOfStreamDivisions
numDivisions = static_cast<unsigned long>(vcl_ceil(static_cast<double>(regionSize)/static_cast<double>(streamMaxSizeBufferForStreaming)));
numDivisions = static_cast<unsigned long>(vcl_ceil(static_cast<double>(regionSize)/static_cast<double>(streamMaxSizeBufferForStreaming)));
}
else
{
//Non streaming
numDivisions = 1;
}
}
break;
default :
itkGenericExceptionMacro(<<"Method use to calculate number of stream divisions is not set !");
break;
}
if( numDivisions == 0) numDivisions = 1;
otbMsgDevMacro(<<" -> Resume : method : "<<mode<<"\n -> Number of divisions = "<<numDivisions);
return(numDivisions);
}
static unsigned int CalculateNeededRadiusForInterpolator(const InterpolationType* interpolator)
{
// unsigned int dimension = TImage::ImageDimension;
unsigned int neededRadius = 0;
std::string className;
className = interpolator->GetNameOfClass();
if (className == "LinearInterpolateImageFunction")
{
otbMsgDevMacro(<<"Linear Interpolator");
neededRadius = 1;
}
else if (className == "NearestNeighborInterpolateImageFunction")
{
otbMsgDevMacro(<<"Nearest Neighbor Interpolator");
neededRadius = 1;
}
else if (className == "BSplineInterpolateImageFunction")
{
otbMsgDevMacro(<<"Nearest Neighbor Interpolator");
neededRadius = 2;
}
break;
default :
itkGenericExceptionMacro(<<"Method use to calculate number of stream divisions is not set !");
break;
/*else if (className == "WindowedSincInterpolateImageFunction")
{
itkGenericExceptionMacro(<< "Windowed Sinc Interpolator not supported yet in resample");
otbMsgDevMacro(<<"Windowed Sinc Interpolator not supported yet in resample");
// dynamic_cast<typename itk::WindowedSincInterpolateImageFunction*>(interpolator);
}
else
{
itkGenericExceptionMacro(<< "Interpolator not recognized, please choose another type !");
} */
return neededRadius;
}
if( numDivisions == 0) numDivisions = 1;
otbMsgDevMacro(<<" -> Resume : method : "<<mode<<"\n -> Number of divisions = "<<numDivisions);
return(numDivisions);
}
};
}// End namespace otb
......
......@@ -345,6 +345,40 @@ ADD_TEST(bfTvVectorImageTo3DScalarImageFilter ${BASICFILTERS_TESTS}
${TEMP}/bfTvVectorImageTo3DScalarImageFilterOutput.png
)
# ------- otb::StreamingResampleImageFilter ----------------------------
ADD_TEST(bfTuStreamingResampleImageFilterNew ${BASICFILTERS_TESTS}
otbStreamingResampleImageFilterNew)
ADD_TEST(bfTvStreamingResampleImageFilter ${BASICFILTERS_TESTS}
otbStreamingResampleImageFilter
${INPUTDATA}/poupees.tif
${TEMP}/bfTvStreamingResamplePoupeesTest.tif
)
ADD_TEST(StreamingResampleImageFilterCompareWithITK ${BASICFILTERS_TESTS}
--compare-image ${TOL}
${TEMP}/bfTvStreamingResamplePoupeesTestITK.tif
${TEMP}/bfTvStreamingResamplePoupeesTestOTB.tif
otbStreamingResampleImageFilterCompareWithITK
${INPUTDATA}/poupees.tif
1100 1100
${TEMP}/bfTvStreamingResamplePoupeesTestITK.tif
${TEMP}/bfTvStreamingResamplePoupeesTestOTB.tif
)
#ADD_TEST(StreamingResampleImageFilterCompareWithITK ${BASICFILTERS_TESTS}
# --compare-image ${TOL}
# ${TEMP}/bfTvStreamingResampleSPOT5TestITK.tif
# ${TEMP}/bfTvStreamingResampleSPOT5TestOTB.tif
# otbStreamingResampleImageFilterCompareWithITK
# ${INPUTDATA}/spot5SubWithGcps.tif
# 1100 1100
# ${TEMP}/bfTvStreamingResampleSPOT5TestITK.tif
# ${TEMP}/bfTvStreamingResampleSPOT5TestOTB.tif
#)
# A enrichir
SET(BasicFilters_SRCS
otbLeeFilter.cxx
......@@ -387,6 +421,9 @@ otbStreamingStatisticsImageFilterNew.cxx
otbStreamingStatisticsImageFilter.cxx
otbVectorImageTo3DScalarImageFilterNew.cxx
otbVectorImageTo3DScalarImageFilter.cxx
otbStreamingResampleImageFilterNew.cxx
otbStreamingResampleImageFilter.cxx
otbStreamingResampleImageFilterCompareWithITK.cxx
)
......
......@@ -67,4 +67,7 @@ REGISTER_TEST(otbStreamingStatisticsImageFilterNew);
REGISTER_TEST(otbStreamingStatisticsImageFilter);
REGISTER_TEST(otbVectorImageTo3DScalarImageFilterNew);
REGISTER_TEST(otbVectorImageTo3DScalarImageFilter);
REGISTER_TEST(otbStreamingResampleImageFilterNew);
REGISTER_TEST(otbStreamingResampleImageFilter);
REGISTER_TEST(otbStreamingResampleImageFilterCompareWithITK);
}
/*=========================================================================
Program: ORFEO Toolbox
Language: C++
Date: $Date$
Version: $Revision$
Copyright (c) Centre National d'Etudes Spatiales. All rights reserved.
See OTBCopyright.txt for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notices for more information.
=========================================================================*/
#include "itkExceptionObject.h"
#include "otbStreamingResampleImageFilter.h"
#include "otbImage.h"
#include "otbImageFileReader.h"
#include "otbStreamingImageFileWriter.h"
#include "itkTranslationTransform.h"
int otbStreamingResampleImageFilter(int argc, char * argv[])
{
try
{
const char* inputFilename = argv[1];
const char* outputFilename = argv[2];
const unsigned int Dimension = 2;
typedef double InputPixelType;
typedef unsigned char OutputPixelType;
typedef double InterpolatorPrecisionType;
typedef otb::Image<InputPixelType,Dimension> InputImageType;
typedef otb::Image<OutputPixelType,Dimension> OutputImageType;
typedef otb::ImageFileReader<InputImageType> ReaderType;
typedef otb::StreamingImageFileWriter<OutputImageType> WriterType;
typedef itk::TranslationTransform<InputPixelType,Dimension> TransformType;
typedef otb::StreamingResampleImageFilter<InputImageType,OutputImageType,InterpolatorPrecisionType> StreamingResampleImageFilterType;
// Instantiating object
ReaderType::Pointer reader = ReaderType::New();
WriterType::Pointer writer = WriterType::New();
StreamingResampleImageFilterType::Pointer resampler = StreamingResampleImageFilterType::New();
TransformType::Pointer transform = TransformType::New();
// Input Image
reader->SetFileName(inputFilename);
// Resampler connected to input image
resampler->SetInput(reader->GetOutput());
// Size of output resampler result
StreamingResampleImageFilterType::SizeType size;
size[0]=600;
size[1]=600;
resampler->SetSize(size);
// Transformation creation
TransformType::OutputVectorType translation;
translation[0] = 10;
translation[1] = 20;
transform->SetOffset(translation);
// Resampler is updated with new transformation (default is identity)
resampler->SetTransform(transform);
// Result of resampler is written
writer->SetInput(resampler->GetOutput());
writer->SetNumberOfStreamDivisions(1);
writer->SetFileName(outputFilename);
writer->Update();
}
catch( itk::ExceptionObject & err )
{
std::cout << "Exception itk::ExceptionObject thrown !" << std::endl;
std::cout << err << std::endl;
return EXIT_FAILURE;
}
catch( ... )
{
std::cout << "Unknown exception thrown !" << std::endl;
return EXIT_FAILURE;
}
return EXIT_SUCCESS;
}
/*=========================================================================
Program: ORFEO Toolbox
Language: C++
Date: $Date$
Version: $Revision$
Copyright (c) Centre National d'Etudes Spatiales. All rights reserved.
See OTBCopyright.txt for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notices for more information.
=========================================================================*/
#include "itkExceptionObject.h"
#include "otbStreamingResampleImageFilter.h"
#include "otbImage.h"
#include "otbVectorImage.h"
#include "otbImageFileReader.h"
#include "otbImageFileWriter.h"
#include "otbStreamingImageFileWriter.h"
#include "itkTranslationTransform.h"
#include "itkNearestNeighborInterpolateImageFunction.h"
#include "otbStreamingTraits.h"
int otbStreamingResampleImageFilterCompareWithITK(int argc, char * argv[])
{
try
{
const char* inputFilename = argv[1];
unsigned int sizeXOutputImage = atoi(argv[2]);
unsigned int sizeYOutputImage = atoi(argv[3]);
const char* outputITKFilename = argv[4];
const char* outputOTBFilename = argv[5];
const unsigned int Dimension = 2;
typedef double InputPixelType;
typedef double OutputPixelType;
typedef double InterpolatorPrecisionType;
typedef otb::Image<InputPixelType,Dimension> InputImageType;
typedef otb::Image<OutputPixelType,Dimension> OutputImageType;
// typedef otb::VectorImage<OutputPixelType,Dimension> InputVectorImageType;
typedef otb::ImageFileReader<InputImageType> ReaderType;
typedef otb::ImageFileWriter<OutputImageType> WriterType;
typedef otb::StreamingImageFileWriter<OutputImageType> StreamingWriterType;
typedef itk::TranslationTransform<InputPixelType,Dimension> TransformType;
typedef itk::NearestNeighborInterpolateImageFunction<InputImageType, InterpolatorPrecisionType> NNInterpolatorType;
typedef itk::ResampleImageFilter<InputImageType,OutputImageType,InterpolatorPrecisionType> ITKResampleImageFilterType;
typedef otb::StreamingResampleImageFilter<InputImageType,OutputImageType,InterpolatorPrecisionType> OTBResampleImageFilterType;
// Instantiating object
ReaderType::Pointer reader = ReaderType::New();
WriterType::Pointer writerITKFilter = WriterType::New();
StreamingWriterType::Pointer writerOTBFilter = StreamingWriterType::New();
ITKResampleImageFilterType::Pointer resamplerITK = ITKResampleImageFilterType::New();
OTBResampleImageFilterType::Pointer resamplerOTB = OTBResampleImageFilterType::New();
TransformType::Pointer transform = TransformType::New();
// Transformation creation
TransformType::OutputVectorType translation;
translation[0] = 30;
translation[1] = 50;
transform->SetOffset(translation);
// Input Image
reader->SetFileName(inputFilename);
// Resamplers connected to input image
resamplerITK->SetInput(reader->GetOutput());
resamplerOTB->SetInput(reader->GetOutput());
// Size of output resamplers result
OTBResampleImageFilterType::SizeType size;
size[0]=sizeXOutputImage;
size[1]=sizeYOutputImage;
resamplerITK->SetSize(size);
resamplerOTB->SetSize(size);
// Set Interpolation
NNInterpolatorType::Pointer interpolator = NNInterpolatorType::New();
// Resamplers are updated with new interpolation (default is linear interpolation)
resamplerITK->SetInterpolator(interpolator);
resamplerOTB->SetInterpolator(interpolator);
// Resamplers are updated with new transformation (default is identity)
resamplerITK->SetTransform(transform);
resamplerOTB->SetTransform(transform);
// Result of resamplers is written
writerITKFilter->SetInput(resamplerITK->GetOutput());
writerITKFilter->SetFileName(outputITKFilename);
writerITKFilter->Update();
writerOTBFilter->SetInput(resamplerOTB->GetOutput());
writerOTBFilter->SetNumberOfStreamDivisions(10);
writerOTBFilter->SetFileName(outputOTBFilename);
writerOTBFilter->Update();
}
catch( itk::ExceptionObject & err )
{
std::cout << "Exception itk::ExceptionObject thrown !" << std::endl;
std::cout << err << std::endl;
return EXIT_FAILURE;
}
catch( ... )
{
std::cout << "Unknown exception thrown !" << std::endl;
return EXIT_FAILURE;
}
return EXIT_SUCCESS;
}
/*=========================================================================
Program: ORFEO Toolbox
Language: C++
Date: $Date$
Version: $Revision$
Copyright (c) Centre National d'Etudes Spatiales. All rights reserved.
See OTBCopyright.txt for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notices for more information.
=========================================================================*/
#include "itkExceptionObject.h"
#include "otbStreamingResampleImageFilter.h"
#include "otbImage.h"
int otbStreamingResampleImageFilterNew(int argc, char * argv[])
{
try
{
const unsigned int Dimension = 2;
typedef double InputPixelType;
typedef unsigned char OutputPixelType;
typedef otb::Image<InputPixelType,Dimension> InputImageType;
typedef otb::Image<OutputPixelType,Dimension> OutputImageType;
typedef double InterpolatorPrecisionType;
typedef otb::StreamingResampleImageFilter<InputImageType,OutputImageType,InterpolatorPrecisionType> StreamingResampleImageFilterType;
// Instantiating object
StreamingResampleImageFilterType::Pointer object = StreamingResampleImageFilterType::New();
}
catch( itk::ExceptionObject & err )
{
std::cout << "Exception itk::ExceptionObject thrown !" << std::endl;
std::cout << err << std::endl;
return EXIT_FAILURE;
}
catch( ... )
{
std::cout << "Unknown exception thrown !" << std::endl;
return EXIT_FAILURE;
}
return EXIT_SUCCESS;
}
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment