Skip to content
Snippets Groups Projects
Commit 840dbf34 authored by Thomas Feuvrier's avatar Thomas Feuvrier
Browse files

Finalisation des classes d'extraction de ROI et de canaux.

parent ab1db447
No related branches found
No related tags found
No related merge requests found
......@@ -20,7 +20,8 @@ namespace otb
{
/** \class ExtractROI
* \brief Extrait une partie d'une image mono-canal.
* \brief Extrait une partie d'une image (mono-canal ou multi-canal).
* Le Pixel peut etre de type simple ou alors un itk::RGBPixel, etc.
* Cette classe s'appuie sur la classe "otb::ExtractROIBase"
*
*/
......
......@@ -20,31 +20,8 @@ namespace otb
{
/** \class ExtractROIBase
* \brief Decrease the image size by cropping the image to the selected
* region bounds.
*
* ExtractROIBase changes the image boundary of an image by removing
* pixels outside the target region. The target region must be specified.
*
* ExtractROIBase also collapses dimensions so that the input image
* may have more dimensions than the output image (i.e. 4-D input image
* to a 3-D output image). To specify what dimensions to collapse,
* the ExtractionRegion must be specified. For any dimension dim where
* ExtractionRegion.Size[dim] = 0, that dimension is collapsed. The
* index to collapse on is specified by ExtractionRegion.Index[dim].
* For example, we have a image 4D = a 4x4x4x4 image, and we want
* to get a 3D image, 3D = a 4x4x4 image, specified as [x,y,z,2] from 4D
* (i.e. the 3rd "time" slice from 4D). The ExtractionRegion.Size =
* [4,4,4,0] and ExtractionRegion.Index = [0,0,0,2].
*
* The number of dimension in ExtractionRegion.Size and Index must =
* InputImageDimension. The number of non-zero dimensions in
* ExtractionRegion.Size must = OutputImageDimension.
*
* This filter is implemented as a multithreaded filter. It provides a
* ThreadedGenerateData() method for its implementation.
* \brief Classe de base, permettant d'extraire une partie d'une image (mono-canal ou multi-canal).
*
* \ingroup GeometricTransforms
*/
template <class TInputImage, class TOutputImage>
class ITK_EXPORT ExtractROIBase:
......
#ifndef __otbVectorToRGBImageFilter_h
#define __otbVectorToRGBImageFilter_h
#include "itkImageToImageFilter.h"
#include "itkImage.h"
#include "itkVectorImage.h"
#include "itkRGBPixel.h"
namespace otb
{
/** \class VectorToRGBImageFilter
*
* \brief Creer une image RGB a partir d'une image multi-canal de type VectorImage
* Il est possible de slectionner les trois canaux RGB de l'image source que l'on souhaite insrer dans l'image destination.
*
* \ingroup Multithreaded
*/
template <class TInputPixel, class TOutputPixel, unsigned int VImageDimension=2>
class ITK_EXPORT VectorToRGBImageFilter :
public itk::ImageToImageFilter< itk::VectorImage < TInputPixel, VImageDimension> , itk::Image<itk::RGBPixel<TOutputPixel> , VImageDimension> >
{
public:
/** Standard class typedefs. */
typedef VectorToRGBImageFilter Self;
typedef itk::ImageToImageFilter< itk::VectorImage < TInputPixel, VImageDimension> , itk::Image<itk::RGBPixel<TOutputPixel> , VImageDimension> > Superclass;
typedef itk::SmartPointer<Self> Pointer;
typedef itk::SmartPointer<const Self> ConstPointer;
/** Standard class macros */
itkNewMacro(Self);
itkTypeMacro(VectorToRGBImageFilter, itk::ImageToImageFilter) ;
typedef typename Superclass::OutputImageRegionType OutputImageRegionType ;
typedef itk::VectorImage<TInputPixel,VImageDimension> InputImageType;
typedef itk::Image<itk::RGBPixel<TOutputPixel>,VImageDimension> OutputImageType;
typedef typename OutputImageType::PixelType OutputImagePixelType ;
/** Set/Get le numero de canal traiter pour le canal rouge */
itkSetMacro(RedChannel,unsigned int);
itkGetConstMacro(RedChannel,unsigned int);
/** Set/Get le numero de canal traiter pour le canal vert */
itkSetMacro(GreenChannel,unsigned int);
itkGetConstMacro(GreenChannel,unsigned int);
/** Set/Get le numero de canal traiter pour le canal bleu */
itkSetMacro(BlueChannel,unsigned int);
itkGetConstMacro(BlueChannel,unsigned int);
protected:
VectorToRGBImageFilter() ;
virtual ~VectorToRGBImageFilter() {}
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread,
int threadId) ;
private:
VectorToRGBImageFilter(const Self&); //purposely not implemented
void operator=(const Self&); //purposely not implemented
/** Contient le numero de canal rouge a extraire */
unsigned int m_RedChannel;
/** Contient le numero de canal vert a extraire */
unsigned int m_GreenChannel;
/** Contient le numero de canal bleu a extraire */
unsigned int m_BlueChannel;
};
} // end namespace otb
#ifndef OTB_MANUAL_INSTANTIATION
#include "otbVectorToRGBImageFilter.txx"
#endif
#endif
#ifndef _otbVectorToRGBImageFilter_txx
#define _otbVectorToRGBImageFilter_txx
#include "otbVectorToRGBImageFilter.h"
#include "itkImageRegionIterator.h"
#include "itkImageRegionConstIterator.h"
#include "itkProgressReporter.h"
#include "itkPixelTraits.h"
namespace otb
{
/**
* Constructor
*/
template <class TInputPixel, class TOutputPixel, unsigned int VImageDimension>
VectorToRGBImageFilter< TInputPixel, TOutputPixel, VImageDimension >
::VectorToRGBImageFilter() : m_RedChannel(1),
m_GreenChannel(2),
m_BlueChannel(3)
{
}
/**
* ThreadedGenerateData Performs the pixel-wise addition
*/
template <class TInputPixel, class TOutputPixel, unsigned int VImageDimension>
void
VectorToRGBImageFilter< TInputPixel, TOutputPixel, VImageDimension >
::ThreadedGenerateData( const OutputImageRegionType &outputRegionForThread,
int threadId)
{
typename OutputImageType::Pointer outputPtr = this->GetOutput();
typename InputImageType::Pointer inputPtr = (InputImageType*)(this->GetInput(0));
itk::ImageRegionIterator<OutputImageType> outputIt(outputPtr, outputRegionForThread);
itk::ProgressReporter progress(this,
threadId,
outputRegionForThread.GetNumberOfPixels());
outputIt.GoToBegin();
typename OutputImageType::PixelType arrayPixel ;
itk::ImageRegionConstIterator< InputImageType > inputIt(inputPtr, outputRegionForThread);
inputIt.GoToBegin();
while( !outputIt.IsAtEnd() )
{
typename OutputImageType::PixelType pix;
pix.SetRed(inputIt.Get().GetElement(m_RedChannel-1));
pix.SetGreen(inputIt.Get().GetElement(m_GreenChannel-1));
pix.SetBlue(inputIt.Get().GetElement(m_BlueChannel-1));
outputIt.Set(pix);
++inputIt;
++outputIt;
progress.CompletedPixel();
}
}
} // end namespace otb
#endif
#ifndef __otbVectorToScalarImageFilter_h
#define __otbVectorToScalarImageFilter_h
#include "itkImageToImageFilter.h"
#include "itkImage.h"
#include "itkVectorImage.h"
namespace otb
{
/** \class VectorToScalarImageFilter
*
* \brief Creer une image mono-canal a partir d'une image multi-canal (type VectorImage)
* Il est possible de slectionner le canal de l'image source que l'on souhaite extraire dans l'image destination.
*
* \ingroup Multithreaded
*/
template <class TInputPixel, class TOutputPixel, unsigned int VImageDimension=2>
class ITK_EXPORT VectorToScalarImageFilter :
public itk::ImageToImageFilter< itk::VectorImage < TInputPixel, VImageDimension> , itk::Image<TOutputPixel , VImageDimension> >
{
public:
/** Standard class typedefs. */
typedef VectorToScalarImageFilter Self;
typedef itk::ImageToImageFilter< itk::VectorImage < TInputPixel, VImageDimension> , itk::Image<TOutputPixel , VImageDimension> > Superclass;
typedef itk::SmartPointer<Self> Pointer;
typedef itk::SmartPointer<const Self> ConstPointer;
/** Standard class macros */
itkNewMacro(Self);
itkTypeMacro(VectorToScalarImageFilter, itk::ImageToImageFilter) ;
typedef typename Superclass::OutputImageRegionType OutputImageRegionType ;
typedef itk::VectorImage<TInputPixel,VImageDimension> InputImageType;
typedef itk::Image<TOutputPixel,VImageDimension> OutputImageType;
typedef typename OutputImageType::PixelType OutputImagePixelType ;
/** Set/Get le numero de canal traiter */
itkSetMacro(Channel,unsigned int);
itkGetConstMacro(Channel,unsigned int);
protected:
VectorToScalarImageFilter() ;
virtual ~VectorToScalarImageFilter() {}
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread,
int threadId) ;
private:
VectorToScalarImageFilter(const Self&); //purposely not implemented
void operator=(const Self&); //purposely not implemented
/** Contient le numero de canal a extraire */
unsigned int m_Channel; // [1...]
};
} // end namespace otb
#ifndef OTB_MANUAL_INSTANTIATION
#include "otbVectorToScalarImageFilter.txx"
#endif
#endif
#ifndef _otbVectorToScalarImageFilter_txx
#define _otbVectorToScalarImageFilter_txx
#include "otbVectorToScalarImageFilter.h"
#include "itkImageRegionIterator.h"
#include "itkImageRegionConstIterator.h"
#include "itkProgressReporter.h"
#include "itkPixelTraits.h"
namespace otb
{
/**
* Constructor
*/
template <class TInputPixel, class TOutputPixel, unsigned int VImageDimension>
VectorToScalarImageFilter< TInputPixel, TOutputPixel, VImageDimension >
::VectorToScalarImageFilter() : m_Channel(1)
{
}
/**
* ThreadedGenerateData Performs the pixel-wise addition
*/
template <class TInputPixel, class TOutputPixel, unsigned int VImageDimension>
void
VectorToScalarImageFilter< TInputPixel, TOutputPixel, VImageDimension >
::ThreadedGenerateData( const OutputImageRegionType &outputRegionForThread,
int threadId)
{
typename OutputImageType::Pointer outputPtr = this->GetOutput();
typename InputImageType::Pointer inputPtr = (InputImageType*)(this->GetInput(0));
itk::ImageRegionIterator<OutputImageType> outputIt(outputPtr, outputRegionForThread);
itk::ProgressReporter progress(this,
threadId,
outputRegionForThread.GetNumberOfPixels());
outputIt.GoToBegin();
typename OutputImageType::PixelType arrayPixel ;
itk::ImageRegionConstIterator< InputImageType > inputIt(inputPtr, outputRegionForThread);
inputIt.GoToBegin();
while( !outputIt.IsAtEnd() )
{
// typename OutputImageType::PixelType pix;
// pix.Set();
outputIt.Set(inputIt.Get().GetElement(m_Channel-1));
++inputIt;
++outputIt;
progress.CompletedPixel();
}
}
} // end namespace otb
#endif
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