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

Ajout d'une nouvelle class de base : ImageListToImageFilter.

parent 6cd21bd0
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 _otbImageListToImageFilter_h
#define _otbImageListToImageFilter_h
#include "itkImageSource.h"
#include "otbImageList.h"
namespace otb
{
/** \class ImageListToImageFilter
*
* \brief Base class for all the filters taking an images list as input to
* produce an image.
*
* \ingroup Images
* \ingroup Lists
*/
template <class TInputImage, class TOutputImage>
class ITK_EXPORT ImageListToImageFilter
: public itk::ImageSource<TOutputImage>
{
public:
/** Standard typedefs */
typedef ImageListToImageFilter Self;
typedef itk::ImageSource<TOutputImage> Superclass;
typedef itk::SmartPointer<Self> Pointer;
typedef itk::SmartPointer<const Self> ConstPointer;
/** Type macro */
itkNewMacro(Self);
/** Creation through object factory macro */
itkTypeMacro(ImageListToImageFilter, ImageSource);
/** Template parameters typedefs */
typedef TInputImage InputImageType;
typedef typename InputImageType::ConstPointer InputImagePointer;
typedef typename InputImageType::RegionType InputImageRegionType;
typedef typename InputImageType::PixelType InputImagePixelType;
typedef typename InputImageType::SizeType SizeType;
typedef typename InputImageType::ValueType ValueType;
typedef ImageList<InputImageType> InputImageListType;
typedef typename InputImageListType::ConstPointer InputImageListConstPointerType;
/** Derived typedefs */
typedef typename Superclass::OutputImageType OutputImageType;
typedef typename OutputImageType::Pointer OutputImagePointerType;
/** InputImage dimension constant */
itkStaticConstMacro(InputImageDimension, unsigned int,TInputImage::ImageDimension);
/** Overiding the SetInput() and GetInput() methods */
virtual void SetInput( const InputImageListType * image);
InputImageListType * GetInput(void);
protected:
/** Constructor */
ImageListToImageFilter();
/** Destructor */
virtual ~ImageListToImageFilter() {};
/**PrintSelf method */
virtual void PrintSelf(std::ostream& os, itk::Indent indent) const;
private:
ImageListToImageFilter(const Self&); //purposely not implemented
void operator=(const Self&); //purposely not implemented
};
}// End namespace otb
#ifndef OTB_MANUAL_INSTANTIATION
#include "otbImageListToImageFilter.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 _otbImageListToImageFilter_txx
#define _otbImageListToImageFilter_txx
#include "otbImageListToImageFilter.h"
namespace otb
{
/**
* Constructor
*/
template <class TInputImage, class TOutputImage>
ImageListToImageFilter<TInputImage,TOutputImage>
::ImageListToImageFilter()
{
this->SetNumberOfRequiredInputs(1);
}
/**
* Input Connection
* \param image The input image.
*/
template <class TInputImage, class TOutputImage>
void
ImageListToImageFilter<TInputImage,TOutputImage>
::SetInput(const InputImageListType *image)
{
// A single input image
this->itk::ProcessObject::SetNthInput(0,const_cast<InputImageListType*>(image));
}
/**
* Input image retrieval
* \return The input image.
*/
template <class TInputImage, class TOutputImage>
typename ImageListToImageFilter<TInputImage,TOutputImage>::InputImageListType *
ImageListToImageFilter<TInputImage,TOutputImage>
::GetInput(void)
{
// If there is no input
if (this->GetNumberOfInputs()<1)
{
// exit
return 0;
}
// else return the first input
return static_cast<InputImageListType * >
(this->itk::ProcessObject::GetInput(0) );
}
/**
* PrintSelf Method
*/
template <class TInputImage, class TOutputImage>
void
ImageListToImageFilter<TInputImage,TOutputImage>
::PrintSelf(std::ostream& os, itk::Indent indent) const
{
Superclass::PrintSelf(os, indent);
}
} // End namespace otb
#endif
......@@ -270,6 +270,11 @@ ADD_TEST(coTuImageListSourceNew ${COMMON_TESTS}
ADD_TEST(coTuImageToImageListFilterNew ${COMMON_TESTS}
otbImageToImageListFilterNew)
# ------- otb::ImageListToImageFilter -------------------------------------------
ADD_TEST(coTuImageListToImageFilterNew ${COMMON_TESTS}
otbImageListToImageFilterNew)
# ------- otb::ImageListToImageListFilter -------------------------------------------
ADD_TEST(coTvImageListToImageListFilterNew ${COMMON_TESTS}
......@@ -410,6 +415,7 @@ otbImageListNew.cxx
otbImageList.cxx
otbImageListSourceNew.cxx
otbImageToImageListFilterNew.cxx
otbImageListToImageFilterNew.cxx
otbImageListToImageListFilterNew.cxx
otbListNew.cxx
otbList.cxx
......
......@@ -56,6 +56,7 @@ REGISTER_TEST(otbImageListNew);
REGISTER_TEST(otbImageList);
REGISTER_TEST(otbImageListSourceNew);
REGISTER_TEST(otbImageToImageListFilterNew);
REGISTER_TEST(otbImageListToImageFilterNew);
REGISTER_TEST(otbImageListToImageListFilterNew);
REGISTER_TEST(otbConcatenateVectorImageFilterNew);
REGISTER_TEST(otbConcatenateVectorImageFilter);
......
/*=========================================================================
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 "otbImageListToImageFilter.h"
#include "otbImage.h"
int otbImageListToImageFilterNew(int argc, char * argv[])
{
try
{
const unsigned int Dimension = 2;
typedef unsigned char InputPixelType;
typedef unsigned char OutputPixelType;
typedef otb::Image<InputPixelType,Dimension> InputImageType;
typedef otb::Image<OutputPixelType,Dimension> OutputImageType;
typedef otb::ImageListToImageFilter<InputImageType,OutputImageType> ImageListToImageFilterType;
// Instantiating ImageListSource object
ImageListToImageFilterType::Pointer imageList = ImageListToImageFilterType::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