Skip to content
Snippets Groups Projects
Commit 7ace7d31 authored by Julien Michel's avatar Julien Michel
Browse files

Un nouveau filtre pour l'OTB, qui permet de générer des quicklook en streaming...

Un nouveau filtre pour l'OTB, qui permet de générer des quicklook en streaming de manière très rapide.
parent 1346cba2
Branches
Tags
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 _otbStreamingShrinkImageFilter_h
#define _otbStreamingShrinkImageFilter_h
#include "itkImageToImageFilter.h"
#include "itkDataObject.h"
namespace otb
{
/** \class StreamingShrinkImageFilter
* \brief This class performs a streaming isotropic shrinking operation without smoothing.
*
* It is intended to be used where a fast quicklook generation is needed for huge images
* (for instance for visualization applications).
*
* It computes the size of the output according to the size of the input image, a read only
* the strip of the input image needed to build a line of the output image. In this strip,
* the pixel are directly selected and passed to the output image.
*
* For exemple, with a 6000X6000 image and a 10 shrinkFactor, it will read 600 lines of 5990 pixels
* instead of the whole image.
*/
template <class TImage>
class ITK_EXPORT StreamingShrinkImageFilter
: public itk::ImageToImageFilter<TImage,TImage>
{
public:
/** Standard typedefs */
typedef StreamingShrinkImageFilter Self;
typedef itk::ImageToImageFilter<TImage,TImage> Superclass;
typedef itk::SmartPointer<Self> Pointer;
typedef itk::SmartPointer<const Self> ConstPointer;
/** Creation through object factory macro */
itkNewMacro(Self);
/** Type macro */
itkTypeMacro(StreamingShrinkImageFilter,ImageToImageFilter);
/** Template parameter typedefs */
typedef TImage ImageType;
typedef typename ImageType::Pointer ImagePointerType;
/** Shrink factor accessor */
itkSetMacro(ShrinkFactor,unsigned int);
itkGetMacro(ShrinkFactor, unsigned int);
/**
* StreamingShrinkImageFilter produces an ouptut whose size is different from its input.
* As such, it must override the GenerateOutputInformation method in order to compute
* the output size from the input size.
*/
void GenerateOutputInformation(void);
/** Main computation method */
virtual void UpdateOutputData(itk::DataObject * output);
protected:
/** Constructor */
StreamingShrinkImageFilter();
/** Destructor */
~StreamingShrinkImageFilter();
/** PrintSelf method */
void PrintSelf(std::ostream& os, itk::Indent indent) const;
private:
/** The shrink factor */
unsigned int m_ShrinkFactor;
};
} // End namespace otb
#ifndef OTB_MANUAL_INSTANTIATION
#include "otbStreamingShrinkImageFilter.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.
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 _otbStreamingShrinkImageFilter_txx
#define _otbStreamingShrinkImageFilter_txx
#include "otbStreamingShrinkImageFilter.h"
#include "itkImageRegionIteratorWithIndex.h"
#include "otbMacro.h"
namespace otb
{
/** Constructor */
template <class TImage>
StreamingShrinkImageFilter<TImage>
::StreamingShrinkImageFilter()
{
// Default shrink factor
m_ShrinkFactor=10;
}
/** Destructor */
template <class TImage>
StreamingShrinkImageFilter<TImage>
::~StreamingShrinkImageFilter()
{}
/**
* StreamingShrinkImageFilter produces an ouptut whose size is different from its input.
* As such, it must override the GenerateOutputInformation method in order to compute
* the output size from the input size.
*/
template <class TImage>
void
StreamingShrinkImageFilter<TImage>
::GenerateOutputInformation(void)
{
// call the superclass' implementation of this method
Superclass::GenerateOutputInformation();
// get pointers to the input and output
typename Superclass::InputImageConstPointer inputPtr = this->GetInput();
typename Superclass::OutputImagePointer outputPtr = this->GetOutput();
// we need to compute the output spacing, the output image size, and the
// output image start index
const typename ImageType::SpacingType&
inputSpacing = inputPtr->GetSpacing();
const typename ImageType::SizeType& inputSize
= inputPtr->GetLargestPossibleRegion().GetSize();
const typename ImageType::IndexType& inputStartIndex
= inputPtr->GetLargestPossibleRegion().GetIndex();
otbMsgDebugMacro(<<"Input idnex "<<inputStartIndex);
otbMsgDebugMacro(<<"Input size: "<<inputSize);
typename ImageType::SpacingType outputSpacing;
typename ImageType::SizeType outputSize;
typename ImageType::IndexType outputStartIndex;
for (unsigned int i = 0; i < ImageType::ImageDimension; i++)
{
outputSpacing[i] = inputSpacing[i] * static_cast<double>( m_ShrinkFactor);
outputSize[i] = inputSize[i]/m_ShrinkFactor;
outputStartIndex[i] = inputStartIndex[i];
}
outputPtr->SetSpacing( outputSpacing );
typename ImageType::RegionType outputLargestPossibleRegion;
outputLargestPossibleRegion.SetSize( outputSize );
outputLargestPossibleRegion.SetIndex( outputStartIndex );
otbMsgDebugMacro(<<"Output largest possible region: "<<outputLargestPossibleRegion);
outputPtr->SetLargestPossibleRegion( outputLargestPossibleRegion );
}
/** Main computation method */
template <class TImage>
void
StreamingShrinkImageFilter<TImage>
::UpdateOutputData(itk::DataObject *itkNotUsed(output))
{
/**
* prevent chasing our tail
*/
if (this->m_Updating)
{
return;
}
/**
* Prepare all the outputs. This may deallocate previous bulk data.
*/
this->PrepareOutputs();
/**
* Make sure we have the necessary inputs
*/
unsigned int ninputs = this->GetNumberOfValidRequiredInputs();
if (ninputs < 1)
{
itkExceptionMacro(<< "At least 1 input is required but only " << ninputs << " are specified.");
return;
}
this->SetAbortGenerateData(0);
this->SetProgress(0.0);
this->m_Updating = true;
/**
* Tell all Observers that the filter is starting
*/
this->InvokeEvent(itk::StartEvent());
/**
* Allocate the output buffer.
*/
ImagePointerType outputPtr = this->GetOutput(0);
typename ImageType::RegionType outputRegion = outputPtr->GetRequestedRegion();
outputPtr->SetBufferedRegion( outputRegion );
outputPtr->Allocate();
/**
* Grab the input
*/
ImagePointerType inputPtr = const_cast<ImageType * >( this->GetInput(0) );
typedef itk::ImageRegionIteratorWithIndex<ImageType> IteratorType;
IteratorType it(outputPtr,outputRegion);
it.GoToBegin();
typename ImageType::SizeType size = outputRegion.GetSize();
for(unsigned int i=0;i<size[1]&&!it.IsAtEnd();++i)
{
typename ImageType::IndexType readIndex;
readIndex[0] = 0;
readIndex[1] = i*m_ShrinkFactor;
typename ImageType::SizeType readSize;
readSize[0]=size[0]*m_ShrinkFactor;
readSize[1]=1;
typename ImageType::RegionType readRegion;
readRegion.SetSize(readSize);
readRegion.SetIndex(readIndex);
inputPtr->SetRequestedRegion(readRegion);
inputPtr->PropagateRequestedRegion();
inputPtr->UpdateOutputData();
IteratorType readIt(inputPtr,readRegion);
unsigned int count=0;
for(readIt.GoToBegin();!readIt.IsAtEnd()&&!it.IsAtEnd();++readIt,++count)
{
if(count%m_ShrinkFactor==0)
{
it.Set(readIt.Get());
++it;
}
}
this->UpdateProgress(static_cast<float>(i)/ static_cast<float>(size[1]));
}
/**
* If we ended due to aborting, push the progress up to 1.0 (since
* it probably didn't end there)
*/
if ( !this->GetAbortGenerateData() )
{
this->UpdateProgress(1.0);
}
// Notify end event observers
this->InvokeEvent(itk::EndEvent() );
/**
* Now we have to mark the data as up to data.
*/
if (this->GetOutput(0))
{
this->GetOutput(0)->DataHasBeenGenerated();
}
/**
* Release any inputs if marked for release
*/
this->ReleaseInputs();
// Mark that we are no longer updating the data in this filter
this->m_Updating = false;
}
/** PrintSelf method */
template <class TImage>
void
StreamingShrinkImageFilter<TImage>
::PrintSelf(std::ostream& os, itk::Indent indent) const
{
Superclass::PrintSelf(os,indent);
os << indent << "Shrink facotr: " << m_ShrinkFactor
<< std::endl;
}
} // End namespace otb
#endif
......@@ -140,6 +140,31 @@ ADD_TEST(bfTvInverseLogPolarTransformResample ${BASICFILTERS_TESTS}
${TEMP}/bfInverseLogPolarTransformResampleOutput.hdr
)
# ------- otb::StreamingShrinkImageFilter ----------------------------
ADD_TEST(bfTuStreamingShrinkImageFilterNew ${BASICFILTERS_TESTS}
otbStreamingShrinkImageFilterNew)
ADD_TEST(bfTvStreamingShrinkImageFilterQBPAN ${BASICFILTERS_TESTS}
--compare-image ${TOL}
${BASELINE}/bfTvStreamingShrinkImageFilterQBPANOutput.hdr
${TEMP}/bfTvStreamingShrinkImageFilterQBPANOutput.hdr
otbStreamingShrinkImageFilter
${IMAGEDATA}/TOULOUSE/QuickBird/000000128955_01_P001_PAN/02APR01105228-P1BS-000000128955_01_P001.TIF
${TEMP}/bfTvStreamingShrinkImageFilterQBPANOutput.hdr
50
)
ADD_TEST(bfTvStreamingShrinkImageFilterQBMUL ${BASICFILTERS_TESTS}
--compare-image ${TOL}
${BASELINE}/bfTvStreamingShrinkImageFilterQBMULOutput.hdr
${TEMP}/bfTvStreamingShrinkImageFilterQBMULOutput.hdr
otbStreamingShrinkImageFilter
${IMAGEDATA}/TOULOUSE/QuickBird/000000128955_01_P001_MUL/02APR01105228-M1BS-000000128955_01_P001.TIF
${TEMP}/bfTvStreamingShrinkImageFilterQBMULOutput.hdr
20
)
# A enrichir
SET(BasicFilters_SRCS
......@@ -160,6 +185,8 @@ otbLogPolarTransformResample.cxx
otbInverseLogPolarTransformNew.cxx
otbInverseLogPolarTransform.cxx
otbInverseLogPolarTransformResample.cxx
otbStreamingShrinkImageFilterNew.cxx
otbStreamingShrinkImageFilter.cxx
)
......
......@@ -44,4 +44,6 @@ REGISTER_TEST(otbLogPolarTransformResample);
REGISTER_TEST(otbInverseLogPolarTransformNew);
REGISTER_TEST(otbInverseLogPolarTransform);
REGISTER_TEST(otbInverseLogPolarTransformResample);
REGISTER_TEST(otbStreamingShrinkImageFilterNew);
REGISTER_TEST(otbStreamingShrinkImageFilter);
}
/*=========================================================================
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 "otbImageFileReader.h"
#include "otbImageFileWriter.h"
#include "otbVectorImage.h"
#include "otbStreamingShrinkImageFilter.h"
int otbStreamingShrinkImageFilter( int argc, char * argv[] )
{
try
{
char * inputFilename = argv[1];
char * outputFilename = argv[2];
unsigned int shrinkFactor = atoi(argv[3]);
const unsigned int Dimension = 2;
typedef unsigned int PixelType;
typedef otb::VectorImage<PixelType,Dimension> ImageType;
typedef otb::ImageFileReader<ImageType> ReaderType;
typedef otb::ImageFileWriter<ImageType> WriterType;
typedef otb::StreamingShrinkImageFilter<ImageType> ShrinkType;
ReaderType::Pointer reader = ReaderType::New();
ShrinkType::Pointer shrink = ShrinkType::New();
WriterType::Pointer writer = WriterType::New();
reader->SetFileName(inputFilename);
writer->SetFileName(outputFilename);
shrink->SetShrinkFactor(shrinkFactor);
shrink->SetInput(reader->GetOutput());
writer->SetInput(shrink->GetOutput());
writer->Update();
}
catch( itk::ExceptionObject & err )
{
std::cout << "Exception itk::ExceptionObject levee !" << std::endl;
std::cout << err << std::endl;
return EXIT_FAILURE;
}
catch( ... )
{
std::cout << "Exception levee inconnue !" << 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 "otbVectorImage.h"
#include "otbStreamingShrinkImageFilter.h"
int otbStreamingShrinkImageFilterNew( int argc, char * argv[] )
{
try
{
const unsigned int Dimension = 2;
typedef unsigned char PixelType;
typedef otb::VectorImage<PixelType,Dimension> ImageType;
typedef otb::StreamingShrinkImageFilter<ImageType> ShrinkType;
ShrinkType::Pointer shrink = ShrinkType::New();
}
catch( itk::ExceptionObject & err )
{
std::cout << "Exception itk::ExceptionObject levee !" << std::endl;
std::cout << err << std::endl;
return EXIT_FAILURE;
}
catch( ... )
{
std::cout << "Exception levee inconnue !" << std::endl;
return EXIT_FAILURE;
}
return EXIT_SUCCESS;
}
......@@ -49,14 +49,6 @@ ADD_TEST(viTvFixedSizFullImageWidget ${VISU_TESTS}
${INPUTDATA}/poupees.png
)
ADD_TEST(viTvFullLargeImageWidget ${VISU_TESTS}
otbFullLargeImageWidget
//${IMAGEDATA}/TOULOUSE/QuickBird/000000128955_01_P001_MUL/02APR01105228-M1BS-000000128955_01_P001.TIF
${IMAGEDATA}/TOULOUSE/QuickBird/000000128955_01_P001_PAN/02APR01105228-P1BS-000000128955_01_P001.TIF
)
# ------- otb::FullResolutionImageWidget -----------------------------------
ADD_TEST(viTuFullResolutionImageWidgetNew ${VISU_TESTS}
......@@ -87,7 +79,6 @@ otbImageViewerNew.cxx
otbImageWidgetBaseNew.cxx
otbFixedSizeFullImageWidgetNew.cxx
otbFixedSizeFullImageWidget.cxx
otbFullLargeImageWidget.cxx
otbFullResolutionImageWidgetNew.cxx
otbFullResolutionImageWidget.cxx
otbZoomableImageWidgetNew.cxx
......
......@@ -31,7 +31,6 @@ REGISTER_TEST(otbImageViewerNew);
REGISTER_TEST(otbImageWidgetBaseNew);
REGISTER_TEST(otbFixedSizeFullImageWidgetNew);
REGISTER_TEST(otbFixedSizeFullImageWidget);
REGISTER_TEST(otbFullLargeImageWidget);
REGISTER_TEST(otbFullResolutionImageWidgetNew);
REGISTER_TEST(otbFullResolutionImageWidget);
REGISTER_TEST(otbZoomableImageWidgetNew);
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment