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

Ajout du filtre de segmentation de la pyramide morphologique.

parent f6ec1a44
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 _otbMorphologicalPyramidSegmenter_h
#define _otbMorphologicalPyramidSegmenter_h
#include "itkImageToImageFilter.h"
namespace otb
{
namespace morphologicalPyramid
{
/**
* \class Segmenter
* \brief This class performs the segmentation of a detail image extracted from a
* morphological pyramid analysis.
*
* This class performs the segmentation of a detail image extracted from a
* morphological pyramid analysis.
*
* The Segmentation is perfomed using the ConnectedThresholdImageFilter. The seeds
* are extracted from the image using the ImageToPointSetFilter. The thresolds are set
* by using quantiles computed with the HistogramGenerator.
*
* \sa MorphologicalPyramidAnalyseFilter
* \sa MorphologicalPyramidSynthesisFilter
* \sa ResampleImageFilter
* \sa LinearInterpolateImageFunction
* \sa ScaleTransform
*/
template <class TInputImage, class TOutputImage>
class Segmenter
: public itk::ImageToImageFilter<TInputImage,TOutputImage>
{
public :
/** Standard typedefs */
typedef Segmenter Self;
typedef itk::ImageToImageFilter<TInputImage,TOutputImage> Superclass;
typedef itk::SmartPointer<Self> Pointer;
typedef itk::SmartPointer<const Self> ConstPointer;
/** Creation through object factory macro */
itkNewMacro(Self);
/** Type macro */
itkTypeMacro(Segmenter,ImageToImageFilter);
/** Template parameters typedefs */
typedef TInputImage InputImageType;
typedef typename InputImageType::PixelType InputPixelType;
typedef typename InputImageType::Pointer InputImagePointerType;
typedef TOutputImage OutputImageType;
typedef typename OutputImageType::Pointer OutputImagePointerType;
typedef typename OutputImageType::PixelType OutputPixelType;
typedef typename OutputImageType::RegionType OutputImageRegionType;
typedef typename OutputImageType::SizeType SizeType;
typedef typename OutputImageType::SpacingType SpacingType;
/** ImageDimension constants */
itkStaticConstMacro(DetailsImageDimension, unsigned int,
TInputImage::ImageDimension);
itkStaticConstMacro(OriginalImageDimension, unsigned int,
TInputImage::ImageDimension);
itkStaticConstMacro(OutputImageDimension, unsigned int,
TOutputImage::ImageDimension);
/**
* Set the details image.
* \param detailsImage The details image from the morphological pyramid
*/
void SetDetailsImage(const InputImageType * detailsImage);
/**
* Set the details image.
* \return detailsImage The input details image.
*/
InputImageType * GetDetailsImage(void);
/**
* Set the original image.
* \param originalImage The original image to segment.
*/
void SetOriginalImage(const InputImageType * originalImage);
/**
* Get the original image.
* \return originalImage The original image to segment.
*/
InputImageType * GetOriginalImage(void);
/** Min object size parameter accessor */
itkSetMacro(MinimumObjectSize,unsigned long);
itkGetMacro(MinimumObjectSize,unsigned long);
itkSetMacro(SeedsQuantile,float);
itkGetMacro(SeedsQuantile,float);
itkSetMacro(ConnectedThresholdQuantile,float);
itkGetMacro(ConnectedThresholdQuantile,float);
itkSetMacro(SegmentDarkDetailsBool,bool);
itkGetMacro(SegmentDarkDetailsBool,bool);
itkBooleanMacro(SegmentDarkDetailsBool);
itkGetMacro(NumberOfObjects,OutputPixelType);
protected:
/** Constructor */
Segmenter();
/** Destructor */
~Segmenter() {};
/** Main computation method */
void GenerateData(void);
/** PrintSelf method */
void PrintSelf(std::ostream& os, itk::Indent indent) const;
/**
* Configure the input datas.
*/
void GenerateInputRequestedRegion(void);
/**
* Configure the output data.
*/
void EnlargeOutputRequestedRegion(void);
private :
Segmenter(const Self&); // purposely not implemented
void operator=(const Self&); // purposely not implemented
/** Minimum size for the segmented object */
unsigned long m_MinimumObjectSize;
/** Quantile for seeds determination */
float m_SeedsQuantile;
/** Quantile to set the connectedThresholdFilter threshold */
float m_ConnectedThresholdQuantile;
/** Set to true if the details to segment are darker than background */
bool m_SegmentDarkDetailsBool;
/** Number of segmented objects */
OutputPixelType m_NumberOfObjects;
};
} // End namespace morphologicalPyramid
}// End namespace otb
#ifndef OTB_MANUAL_INSTANTIATION
#include "otbMorphologicalPyramidSegmenter.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 _otbMorphologicalPyramidSegmenter_txx
#define _otbMorphologicalPyramidSegmenter_txx
#include "otbMorphologicalPyramidSegmenter.h"
#include "otbImage.h"
#include "itkConnectedThresholdImageFilter.h"
#include "itkConnectedComponentImageFilter.h"
#include "itkRelabelComponentImageFilter.h"
#include "itkThresholdImageFilter.h"
#include "itkMultiplyImageFilter.h"
#include "itkCastImageFilter.h"
#include "itkInvertIntensityImageFilter.h"
#include "itkRescaleIntensityImageFilter.h"
#include "otbThresholdImageToPointSetFilter.h"
#include "itkScalarImageToHistogramGenerator.h"
namespace otb
{
namespace morphologicalPyramid
{
/**
* Constructor
*/
template <class TInputImage,class TOuputImage>
Segmenter<TInputImage, TOuputImage>
::Segmenter()
{
this->SetNumberOfRequiredInputs(2);
}
/**
* Set the details image.
* \param detailsImage The details image from the morphological pyramid
*/
template <class TInputImage,class TOuputImage>
void
Segmenter<TInputImage, TOuputImage>
::SetDetailsImage(const InputImageType * detailsImage)
{
this->SetNthInput(0,const_cast<TInputImage *>(detailsImage));
}
/**
* Set the details image.
* \return detailsImage The input details image.
*/
template <class TInputImage,class TOuputImage>
typename Segmenter<TInputImage, TOuputImage>::InputImageType *
Segmenter<TInputImage, TOuputImage>
::GetDetailsImage(void)
{
return const_cast<InputImageType *>(this->GetInput(0));
}
/**
* Set the original image.
* \param originalImage The original image to segment.
*/
template <class TInputImage,class TOuputImage>
void
Segmenter<TInputImage, TOuputImage>
::SetOriginalImage(const InputImageType * originalImage)
{
this->SetNthInput(1,const_cast<TInputImage *>(originalImage));
}
/**
* Get the original image.
* \return originalImage The original image to segment.
*/
template <class TInputImage,class TOuputImage>
typename Segmenter<TInputImage, TOuputImage>::InputImageType *
Segmenter<TInputImage, TOuputImage>
::GetOriginalImage(void)
{
return const_cast<InputImageType *>(this->GetInput(1));
}
/**
* Configure the input datas.
*/
template <class TInputImage,class TOuputImage>
void
Segmenter<TInputImage, TOuputImage>
::GenerateInputRequestedRegion()
{
// call the superclass' implementation of this method
Superclass::GenerateInputRequestedRegion();
// get pointers to the inputs
InputImagePointerType detailsPtr =
const_cast< InputImageType * >( this->GetInput(0) );
InputImagePointerType origPtr =
const_cast< InputImageType * >( this->GetInput(1) );
if ( !detailsPtr || !origPtr )
{
return;
}
// We need to
// configure the inputs such that all the data is available.
detailsPtr->SetRequestedRegion(detailsPtr->GetLargestPossibleRegion());
origPtr->SetRequestedRegion(origPtr->GetLargestPossibleRegion());
}
/**
* Configure the output data
*/
template <class TInputImage,class TOuputImage>
void
Segmenter<TInputImage, TOuputImage>
::EnlargeOutputRequestedRegion(void)
{
this->GetOutput()
->SetRequestedRegion( this->GetOutput()->GetLargestPossibleRegion() );
}
/**
* Main computation method
*/
template <class TInputImage,class TOuputImage>
void
Segmenter<TInputImage, TOuputImage>
::GenerateData()
{
// Input images pointers
InputImagePointerType details = this->GetDetailsImage();
InputImagePointerType original = this->GetOriginalImage();
// Typedefs for details image enhancement
typedef float FloatPixelType;
typedef otb::Image<FloatPixelType,InputImageType::ImageDimension> FloatImageType;
typedef itk::InvertIntensityImageFilter<InputImageType,InputImageType> InvertFilterType;
typedef itk::MultiplyImageFilter<FloatImageType,FloatImageType,FloatImageType> MultiplyFilterType;
typedef itk::CastImageFilter<InputImageType,FloatImageType> CastImageFilterType;
typedef itk::RescaleIntensityImageFilter<FloatImageType,InputImageType> RescaleFilterType;
//Typedefs for seeds extraction
typedef itk::PointSet<InputPixelType,InputImageType::ImageDimension> PointSetType;
typedef otb::ThresholdImageToPointSetFilter<InputImageType,PointSetType> PointSetFilterType;
typedef typename PointSetType::PointsContainer::Iterator PointSetIteratorType;
typedef typename PointSetType::PointType PointType;
// Typedefs for segmentation
typedef itk::ConnectedThresholdImageFilter<InputImageType,InputImageType> ConnectedFilterType;
typedef itk::ConnectedComponentImageFilter<InputImageType,OutputImageType> LabelFilterType;
typedef itk::RelabelComponentImageFilter<OutputImageType,OutputImageType> RelabelFilterType;
typedef itk::ThresholdImageFilter<OutputImageType> ThresholdFilterType;
// Typedefs for statistics computation
typedef itk::Statistics::ScalarImageToHistogramGenerator<InputImageType> HistGeneratorType;
typedef typename HistGeneratorType::HistogramType HistogramType;
/////////////////////////////////////
//// Details image enhancement //////
/////////////////////////////////////
// Filters instantiation
typename InvertFilterType::Pointer invert;
typename CastImageFilterType::Pointer cast1 = CastImageFilterType::New();
typename CastImageFilterType::Pointer cast2 = CastImageFilterType::New();
typename MultiplyFilterType::Pointer mult = MultiplyFilterType::New();
typename RescaleFilterType::Pointer rescaler = RescaleFilterType::New();
// Pipeline connection
cast1->SetInput(details);
// If we want to segment darker detail, the original image must have its itensity inverted
if(m_SegmentDarkDetailsBool)
{
invert = InvertFilterType::New();
invert->SetInput(original);
cast2->SetInput(invert->GetOutput());
}
else
{
cast2->SetInput(original);
}
mult->SetInput1(cast1->GetOutput());
mult->SetInput2(cast2->GetOutput());
rescaler->SetInput(mult->GetOutput());
rescaler->SetOutputMinimum(0);
rescaler->SetOutputMaximum(255);
rescaler->Update();
/////////////////////////////////////
//// Thresholds computation /////////
/////////////////////////////////////
// Filter instantiation
typename HistGeneratorType::Pointer histogram = HistGeneratorType::New();
// Thresholds are computed from the quantile
histogram->SetInput(rescaler->GetOutput());
histogram->SetNumberOfBins(255);
histogram->SetMarginalScale(10.0);
histogram->Compute();
InputPixelType pointSetThreshold =
static_cast<InputPixelType>(histogram->GetOutput()->Quantile(0,m_SeedsQuantile));
InputPixelType connectedThresholdValue =
static_cast<InputPixelType>(histogram->GetOutput()->Quantile(0,m_ConnectedThresholdQuantile));
/////////////////////////////////////
//// Seeds extraction ///////////////
/////////////////////////////////////
typename PointSetFilterType::Pointer pointSetFilter = PointSetFilterType::New();
pointSetFilter->SetInput(0,details);
pointSetFilter->SetThreshold(pointSetThreshold);
pointSetFilter->Update();
/////////////////////////////////////
//// Segmentation ///////////////////
/////////////////////////////////////
// Filters instantiation
typename ConnectedFilterType::Pointer connectedThreshold = ConnectedFilterType::New();
typename LabelFilterType::Pointer labeler = LabelFilterType::New();
typename RelabelFilterType::Pointer relabeler = RelabelFilterType::New();
typename ThresholdFilterType::Pointer threshold = ThresholdFilterType::New();
//Passing seeds to the connected filter
connectedThreshold = ConnectedFilterType::New();
connectedThreshold->ClearSeeds();
connectedThreshold->SetInput(rescaler->GetOutput());
PointSetIteratorType it = pointSetFilter->GetOutput()->GetPoints()->Begin();
while(it!=pointSetFilter->GetOutput()->GetPoints()->End())
{
typename OutputImageType::IndexType index;
index[0]=static_cast<long int>(it.Value()[0]*(rescaler->GetOutput()->GetSpacing()[0]));
index[1]=static_cast<long int>(it.Value()[1]*(rescaler->GetOutput()->GetSpacing()[1]));
connectedThreshold->AddSeed(index);
it++;
}
// segmentation
connectedThreshold->SetLower(connectedThresholdValue);
// labelling
labeler = LabelFilterType::New();
relabeler = RelabelFilterType::New();
labeler->SetInput(connectedThreshold->GetOutput());
relabeler->SetInput(labeler->GetOutput());
relabeler->SetMinimumObjectSize(m_MinimumObjectSize);
relabeler->Update();
// In some cases it might happen that the whole extent of the image is segmented as a single region.
// Since this is not desirable, we test this case here to avoid it.
threshold = ThresholdFilterType::New();
threshold->SetInput(relabeler->GetOutput());
OutputPixelType num = 0;
if(relabeler->GetNumberOfObjects()==1)
{
int surface = rescaler->GetOutput()->GetLargestPossibleRegion().GetSize()[0]
*rescaler->GetOutput()->GetLargestPossibleRegion().GetSize()[1];
if(relabeler->GetSizeOfObjectsInPixels()[0]==surface)
{
num = 0;
}
else
{
num=1;
}
}
else
{
num= static_cast<OutputPixelType>(relabeler->GetNumberOfObjects());
}
threshold->ThresholdOutside(0,num);
// Output connection
threshold->GraftOutput(this->GetOutput());
threshold->Update();
this->GraftOutput(threshold->GetOutput());
m_NumberOfObjects = num;
}
/**
* PrintSelf method
*/
template <class TInputImage,class TOuputImage>
void
Segmenter<TInputImage, TOuputImage>
::PrintSelf(std::ostream& os, itk::Indent indent) const
{
Superclass::PrintSelf(os,indent);
}
} // End namespace morphologicalPyramid
} // End namespace otb
#endif
...@@ -95,6 +95,39 @@ ADD_TEST(msTvMorphoPyrMRToMSConverter ${MULTISCALE_TESTS} ...@@ -95,6 +95,39 @@ ADD_TEST(msTvMorphoPyrMRToMSConverter ${MULTISCALE_TESTS}
4 4
2.0) 2.0)
# ------- otb::MorphologicalPyramidSegmenter ----------
ADD_TEST(msTuMorphoPyrSegmenterNew ${MULTISCALE_TESTS}
otbMorphologicalPyramidSegmenterNew)
ADD_TEST(msTvMorphoPyrSegmentBright ${MULTISCALE_TESTS}
--compare-image ${TOL}
${BASELINE}/msPyrSegmenter_IKO_Halles_4_2_sf.tif
${TEMP}/msPyrSegmenter_IKO_Halles_4_2_sf.tif
otbMorphologicalPyramidSegmenter
${TEMP}/msPyrMRToMS_IKO_Halles_4_2_sf_full.tif
${INPUTDATA}/ROI_IKO_PAN_LesHalles.tif
${TEMP}/msPyrSegmenter_IKO_Halles_4_2_sf.tif
0
0.9
0.9
10
)
ADD_TEST(msTvMorphoPyrSegmentDark ${MULTISCALE_TESTS}
--compare-image ${TOL}
${BASELINE}/msPyrSegmenter_IKO_Halles_4_2_if.tif
${TEMP}/msPyrSegmenter_IKO_Halles_4_2_if.tif
otbMorphologicalPyramidSegmenter
${TEMP}/msPyrMRToMS_IKO_Halles_4_2_if_full.tif
${INPUTDATA}/ROI_IKO_PAN_LesHalles.tif
${TEMP}/msPyrSegmenter_IKO_Halles_4_2_if.tif
1
0.9
0.9
10
)
# ------- Fichiers sources CXX ----------------------------------- # ------- Fichiers sources CXX -----------------------------------
SET(BasicMultiScale_SRCS SET(BasicMultiScale_SRCS
otbMorphologicalPyramidResamplerNew.cxx otbMorphologicalPyramidResamplerNew.cxx
...@@ -105,6 +138,8 @@ otbMorphologicalPyramidSynthesisFilterNew.cxx ...@@ -105,6 +138,8 @@ otbMorphologicalPyramidSynthesisFilterNew.cxx
otbMorphologicalPyramidSynthesisFilter.cxx otbMorphologicalPyramidSynthesisFilter.cxx
otbMorphologicalPyramidMRToMSConverterNew.cxx otbMorphologicalPyramidMRToMSConverterNew.cxx
otbMorphologicalPyramidMRToMSConverter.cxx otbMorphologicalPyramidMRToMSConverter.cxx
otbMorphologicalPyramidSegmenterNew.cxx
otbMorphologicalPyramidSegmenter.cxx
) )
......
/*=========================================================================
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 "otbMorphologicalPyramidSegmenter.h"
#include "otbImageFileReader.h"
#include "otbImageFileWriter.h"
#include "otbImage.h"
int otbMorphologicalPyramidSegmenter(int argc, char * argv[])
{
try
{
const char* inputFilename = argv[1];
const char* originalFilename = argv[2];
const char* outputFilename1 = argv[3];
const bool segmentDark = atoi(argv[4]);
const float seedsQuantile = atof(argv[5]);
const float segmentationQuantile = atof(argv[6]);
const unsigned int minObjectSize = atoi(argv[7]);
const unsigned int Dimension = 2;
typedef unsigned char InputPixelType;
typedef unsigned short OutputPixelType;
typedef otb::Image<InputPixelType,Dimension> InputImageType;
typedef otb::Image<OutputPixelType,Dimension> OutputImageType;
typedef otb::ImageFileReader<InputImageType> ReaderType;
typedef otb::ImageFileWriter<OutputImageType> WriterType;
typedef otb::morphologicalPyramid::Segmenter<InputImageType,OutputImageType>
SegmenterType;
// Input images reading
ReaderType::Pointer reader = ReaderType::New();
reader->SetFileName(inputFilename);
ReaderType::Pointer reader2 = ReaderType::New();
reader2->SetFileName(originalFilename);
// Instantiation
SegmenterType::Pointer segmenter = SegmenterType::New();
segmenter->SetDetailsImage(reader->GetOutput());
segmenter->SetOriginalImage(reader2->GetOutput());
segmenter->SetSegmentDarkDetailsBool(segmentDark);
segmenter->SetSeedsQuantile(seedsQuantile);
segmenter->SetConnectedThresholdQuantile(segmentationQuantile);
segmenter->SetMinimumObjectSize(minObjectSize);
// File writing
WriterType::Pointer writer = WriterType::New();
writer->SetInput(segmenter->GetOutput());
writer->SetFileName(outputFilename1);
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 "otbMorphologicalPyramidSegmenter.h"
#include "otbImage.h"
int otbMorphologicalPyramidSegmenterNew(int argc, char * argv[])
{
try
{
const unsigned int Dimension = 2;
typedef unsigned char InputPixelType;
typedef unsigned short OutputPixelType;
typedef otb::Image<InputPixelType,Dimension> InputImageType;
typedef otb::Image<OutputPixelType,Dimension> OutputImageType;
typedef otb::morphologicalPyramid::Segmenter<InputImageType,OutputImageType>
SegmenterType;
// Instantiation
SegmenterType::Pointer segmenter = SegmenterType::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;
}
...@@ -34,4 +34,6 @@ REGISTER_TEST(otbMorphologicalPyramidSynthesisFilterNew); ...@@ -34,4 +34,6 @@ REGISTER_TEST(otbMorphologicalPyramidSynthesisFilterNew);
REGISTER_TEST(otbMorphologicalPyramidSynthesisFilter); REGISTER_TEST(otbMorphologicalPyramidSynthesisFilter);
REGISTER_TEST(otbMorphologicalPyramidMRToMSConverterNew); REGISTER_TEST(otbMorphologicalPyramidMRToMSConverterNew);
REGISTER_TEST(otbMorphologicalPyramidMRToMSConverter); REGISTER_TEST(otbMorphologicalPyramidMRToMSConverter);
REGISTER_TEST(otbMorphologicalPyramidSegmenterNew);
REGISTER_TEST(otbMorphologicalPyramidSegmenter);
} }
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment