Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017 #ifndef __itkResampleImageFilter_h
00018 #define __itkResampleImageFilter_h
00019
00020
00021
00022
00023 #include "itkConfigure.h"
00024
00025 #ifdef ITK_USE_OPTIMIZED_REGISTRATION_METHODS
00026 #include "itkOptResampleImageFilter.h"
00027 #else
00028
00029 #include "itkFixedArray.h"
00030 #include "itkTransform.h"
00031 #include "itkImageFunction.h"
00032 #include "itkImageRegionIterator.h"
00033 #include "itkImageToImageFilter.h"
00034 #include "itkInterpolateImageFunction.h"
00035 #include "itkSize.h"
00036
00037 namespace itk
00038 {
00039
00077 template <class TInputImage, class TOutputImage,
00078 class TInterpolatorPrecisionType=double>
00079 class ITK_EXPORT ResampleImageFilter:
00080 public ImageToImageFilter<TInputImage, TOutputImage>
00081 {
00082 public:
00084 typedef ResampleImageFilter Self;
00085 typedef ImageToImageFilter<TInputImage,TOutputImage> Superclass;
00086 typedef SmartPointer<Self> Pointer;
00087 typedef SmartPointer<const Self> ConstPointer;
00088
00089 typedef TInputImage InputImageType;
00090 typedef TOutputImage OutputImageType;
00091 typedef typename InputImageType::Pointer InputImagePointer;
00092 typedef typename InputImageType::ConstPointer InputImageConstPointer;
00093 typedef typename OutputImageType::Pointer OutputImagePointer;
00094 typedef typename InputImageType::RegionType InputImageRegionType;
00095
00097 itkNewMacro(Self);
00098
00100 itkTypeMacro(ResampleImageFilter, ImageToImageFilter);
00101
00103 itkStaticConstMacro(ImageDimension, unsigned int,
00104 TOutputImage::ImageDimension);
00105 itkStaticConstMacro(InputImageDimension, unsigned int,
00106 TInputImage::ImageDimension);
00108
00109
00111 typedef Transform<TInterpolatorPrecisionType,
00112 itkGetStaticConstMacro(ImageDimension),
00113 itkGetStaticConstMacro(ImageDimension)> TransformType;
00114 typedef typename TransformType::ConstPointer TransformPointerType;
00116
00118 typedef InterpolateImageFunction<InputImageType, TInterpolatorPrecisionType> InterpolatorType;
00119 typedef typename InterpolatorType::Pointer InterpolatorPointerType;
00120
00122 typedef Size<itkGetStaticConstMacro(ImageDimension)> SizeType;
00123
00125 typedef typename TOutputImage::IndexType IndexType;
00126
00128 typedef typename InterpolatorType::PointType PointType;
00129
00130
00132 typedef typename TOutputImage::PixelType PixelType;
00133 typedef typename TInputImage::PixelType InputPixelType;
00134
00136 typedef typename TOutputImage::RegionType OutputImageRegionType;
00137
00139 typedef typename TOutputImage::SpacingType SpacingType;
00140 typedef typename TOutputImage::PointType OriginPointType;
00141 typedef typename TOutputImage::DirectionType DirectionType;
00142
00144 typedef ImageBase<itkGetStaticConstMacro(ImageDimension)> ImageBaseType;
00145
00153 itkSetConstObjectMacro( Transform, TransformType );
00154
00156 itkGetConstObjectMacro( Transform, TransformType );
00157
00164 itkSetObjectMacro( Interpolator, InterpolatorType );
00165
00167 itkGetConstObjectMacro( Interpolator, InterpolatorType );
00168
00170 itkSetMacro( Size, SizeType );
00171
00173 itkGetConstReferenceMacro( Size, SizeType );
00174
00177 itkSetMacro( DefaultPixelValue, PixelType );
00178
00180 itkGetConstReferenceMacro( DefaultPixelValue, PixelType );
00181
00183 itkSetMacro( OutputSpacing, SpacingType );
00184 virtual void SetOutputSpacing( const double* values );
00186
00188 itkGetConstReferenceMacro( OutputSpacing, SpacingType );
00189
00191 itkSetMacro( OutputOrigin, OriginPointType );
00192 virtual void SetOutputOrigin( const double* values);
00194
00196 itkGetConstReferenceMacro( OutputOrigin, OriginPointType );
00197
00199 itkSetMacro( OutputDirection, DirectionType );
00200 itkGetConstReferenceMacro( OutputDirection, DirectionType );
00202
00204 void SetOutputParametersFromImage ( const ImageBaseType * image );
00205
00208 itkLegacyMacro( void SetOutputParametersFromConstImage ( const ImageBaseType * image ) );
00209
00212 itkSetMacro( OutputStartIndex, IndexType );
00213
00215 itkGetConstReferenceMacro( OutputStartIndex, IndexType );
00216
00223 void SetReferenceImage ( const TOutputImage *image );
00224 const TOutputImage * GetReferenceImage( void ) const;
00226
00227 itkSetMacro( UseReferenceImage, bool );
00228 itkBooleanMacro( UseReferenceImage );
00229 itkGetConstMacro( UseReferenceImage, bool );
00230
00236 virtual void GenerateOutputInformation( void );
00237
00243 virtual void GenerateInputRequestedRegion( void );
00244
00247 virtual void BeforeThreadedGenerateData( void );
00248
00251 virtual void AfterThreadedGenerateData( void );
00252
00254 unsigned long GetMTime( void ) const;
00255
00256 #ifdef ITK_USE_CONCEPT_CHECKING
00257
00258 itkConceptMacro(OutputHasNumericTraitsCheck,
00259 (Concept::HasNumericTraits<PixelType>));
00260
00262 #endif
00263
00264 protected:
00265 ResampleImageFilter( void );
00266 ~ResampleImageFilter( void ) {};
00267
00268 void PrintSelf( std::ostream& os, Indent indent ) const;
00269
00278 void ThreadedGenerateData( const OutputImageRegionType& outputRegionForThread,
00279 int threadId );
00280
00283 void NonlinearThreadedGenerateData( const OutputImageRegionType& outputRegionForThread,
00284 int threadId );
00285
00289 void LinearThreadedGenerateData( const OutputImageRegionType&
00290 outputRegionForThread,
00291 int threadId );
00292
00293
00294 private:
00295 ResampleImageFilter( const Self& );
00296 void operator=( const Self& );
00297
00298 SizeType m_Size;
00299 TransformPointerType m_Transform;
00300 InterpolatorPointerType m_Interpolator;
00301
00302 PixelType m_DefaultPixelValue;
00303
00304
00305 SpacingType m_OutputSpacing;
00306 OriginPointType m_OutputOrigin;
00307 DirectionType m_OutputDirection;
00308 IndexType m_OutputStartIndex;
00309 bool m_UseReferenceImage;
00310
00311 };
00312
00313
00314 }
00315
00316 #ifndef ITK_MANUAL_INSTANTIATION
00317 #include "itkResampleImageFilter.txx"
00318 #endif
00319
00320 #endif
00321
00322 #endif
00323