]> Creatis software - cpPlugins.git/blob - lib/ivq/ITK/IsoImageSlicer.h
...
[cpPlugins.git] / lib / ivq / ITK / IsoImageSlicer.h
1 // =======================================================================
2 // @author: Leonardo Florez-Valencia
3 // @email: florez-l@javeriana.edu.co
4 // =======================================================================
5 #ifndef __ivq__ITK__IsoImageSlicer__h__
6 #define __ivq__ITK__IsoImageSlicer__h__
7
8 #include <itkAffineTransform.h>
9 #include <itkExtractImageFilter.h>
10 #include <itkImage.h>
11 #include <itkImageToImageFilter.h>
12 #include <itkInterpolateImageFunction.h>
13 #include <itkResampleImageFilter.h>
14 #include <itkVectorResampleImageFilter.h>
15 #include <itkVectorInterpolateImageFunction.h>
16
17 namespace ivq
18 {
19   namespace ITK
20   {
21     /**
22      */
23     template< class _TSlicer, class _TInterpolateFunction >
24     class BaseImageSlicer
25       : public itk::ImageToImageFilter< typename _TSlicer::InputImageType, itk::Image< typename _TSlicer::InputImageType::PixelType, _TSlicer::ImageDimension - 1 > >
26     {
27     public:
28       // Basic types
29       typedef BaseImageSlicer                                Self;
30       typedef _TSlicer                                    TSlicer;
31       typedef _TInterpolateFunction          TInterpolateFunction;
32       typedef typename TSlicer::InputImageType             TImage;
33       typedef typename TInterpolateFunction::CoordRepType TScalar;
34       typedef typename TImage::PixelType                   TPixel;
35       enum
36       {
37         Dim      = TImage::ImageDimension,
38         SliceDim = TImage::ImageDimension - 1
39       };
40       typedef itk::Image< TPixel, Self::SliceDim > TSliceImage;
41
42       // itk types
43       typedef itk::ImageToImageFilter< TImage, TSliceImage > Superclass;
44       typedef itk::SmartPointer< Self >                      Pointer;
45       typedef itk::SmartPointer< const Self >                ConstPointer;
46
47       // Internal filters
48       typedef itk::ExtractImageFilter< TImage, TSliceImage > TCollapsor;
49
50       // Various types
51       typedef typename TImage::IndexType   TIndex;
52       typedef typename TImage::RegionType  TRegion;
53       typedef typename TImage::SizeType    TSize;
54       typedef typename TImage::SpacingType TSpacing;
55       typedef typename TSpacing::ValueType TSpacingValue;
56
57       typedef itk::AffineTransform< TScalar, Self::Dim > TTransform;
58       typedef typename TTransform::MatrixType            TMatrix;
59       typedef typename TTransform::OffsetType            TVector;
60
61     public:
62       itkNewMacro( Self );
63       itkTypeMacro( BaseImageSlicer, itkImageToImageFilter );
64
65       itkBooleanMacro( SizeFromMaximum );
66       itkBooleanMacro( SizeFromMinimum );
67       itkBooleanMacro( SpacingFromMaximum );
68       itkBooleanMacro( SpacingFromMinimum );
69
70       itkGetConstObjectMacro( Transform, TTransform );
71       itkGetConstMacro( DefaultValue, TPixel );
72       itkGetConstMacro( Size, TVector );
73       itkGetConstMacro( SizeFromMaximum, bool );
74       itkGetConstMacro( SizeFromMinimum, bool );
75       itkGetConstMacro( Spacing, TSpacingValue );
76       itkGetConstMacro( SpacingFromMaximum, bool );
77       itkGetConstMacro( SpacingFromMinimum, bool );
78
79       itkSetObjectMacro( Transform, TTransform );
80       itkSetMacro( Size, TVector );
81       itkSetMacro( DefaultValue, TPixel );
82       itkSetMacro( SizeFromMaximum, bool );
83       itkSetMacro( SizeFromMinimum, bool );
84       itkSetMacro( Spacing, TSpacingValue );
85       itkSetMacro( SpacingFromMaximum, bool );
86       itkSetMacro( SpacingFromMinimum, bool );
87
88     public:
89       virtual unsigned long GetMTime( ) const override;
90
91       const TInterpolateFunction* GetInterpolator( ) const;
92       const TMatrix& GetRotation( ) const;
93       const TVector& GetTranslation( ) const;
94
95       void SetInterpolator( TInterpolateFunction* f );
96
97       template< class _TOtherMatrix >
98       void SetRotation( const _TOtherMatrix& r );
99
100       template< class _TOtherVector >
101       void SetTranslation( const _TOtherVector& t );
102
103       void SetSize( TScalar s );
104
105     protected:
106       BaseImageSlicer( );
107       virtual ~BaseImageSlicer( );
108
109       virtual void GenerateOutputInformation( ) override; // TODO { }
110       virtual void GenerateInputRequestedRegion( ) override;
111       virtual void GenerateData( ) override;
112
113     private:
114       // Purposely not implemented
115       BaseImageSlicer( const Self& );
116       void operator=( const Self& );
117
118     protected:
119       typename TSlicer::Pointer    m_Slicer;
120       typename TCollapsor::Pointer m_Collapsor;
121       typename TTransform::Pointer m_Transform;
122
123       TPixel m_DefaultValue;
124
125       TVector m_Size;
126       bool    m_SizeFromMaximum;
127       bool    m_SizeFromMinimum;
128
129       TSpacingValue m_Spacing;
130       bool          m_SpacingFromMaximum;
131       bool          m_SpacingFromMinimum;
132     };
133
134   } // ecapseman
135
136 } // ecapseman
137
138 // -------------------------------------------------------------------------
139 #define ivq_ITK_DefineIsoImageSlicer( name, R, F )                      \
140   template< class _TImage, class _TScalar = double >                    \
141   class name                                                            \
142     : public ivq::ITK::BaseImageSlicer< R< _TImage, _TImage, _TScalar >, F< _TImage, _TScalar > > \
143   {                                                                     \
144   public:                                                               \
145     typedef ivq::ITK::BaseImageSlicer< R< _TImage, _TImage, _TScalar >, F< _TImage, _TScalar > > Superclass; \
146     typedef name                                       Self;            \
147     typedef itk::SmartPointer< Self >                  Pointer;         \
148     typedef itk::SmartPointer< const Self >            ConstPointer;    \
149   public:                                                               \
150     itkNewMacro( Self );                                                \
151     itkTypeMacro( ivq::ITK::name, ivq::ITK::BaseSlicer );               \
152   protected:                                                            \
153     name( ) : Superclass( ) { }                                         \
154     virtual ~name( ) { }                                                \
155   private:                                                              \
156     name( const Self& );                                                \
157     void operator=( const Self& );                                      \
158   };
159
160 namespace ivq
161 {
162   namespace ITK
163   {
164     ivq_ITK_DefineIsoImageSlicer(
165       IsoImageSlicer,
166       itk::ResampleImageFilter,
167       itk::InterpolateImageFunction
168       );
169     ivq_ITK_DefineIsoImageSlicer(
170       VectorIsoImageSlicer,
171       itk::VectorResampleImageFilter,
172       itk::VectorInterpolateImageFunction
173       );
174
175   } // ecapseman
176
177 } // ecapseman
178
179 // -------------------------------------------------------------------------
180 template< class _TSlicer, class _TInterpolateFunction >
181 template< class _TOtherMatrix >
182 void ivq::ITK::BaseImageSlicer< _TSlicer, _TInterpolateFunction >::
183 SetRotation( const _TOtherMatrix& r )
184 {
185   TMatrix rotation;
186   for( unsigned int i = 0; i < Self::Dim; ++i )
187     for( unsigned int j = 0; j < Self::Dim; ++j )
188       rotation[ i ][ j ] = r[ i ][ j ];
189   this->m_Transform->SetMatrix( rotation );
190   this->Modified( );
191 }
192
193 // -------------------------------------------------------------------------
194 template< class _TSlicer, class _TInterpolateFunction >
195 template< class _TOtherVector >
196 void ivq::ITK::BaseImageSlicer< _TSlicer, _TInterpolateFunction >::
197 SetTranslation( const _TOtherVector& t )
198 {
199   TVector off;
200   for( unsigned int i = 0; i < Self::Dim; ++i )
201     off[ i ] = t[ i ];
202   this->m_Transform->SetOffset( off );
203   this->Modified( );
204 }
205
206 #ifndef ITK_MANUAL_INSTANTIATION
207 #  include <ivq/ITK/IsoImageSlicer.hxx>
208 #endif // ITK_MANUAL_INSTANTIATION
209 #endif // __ivq__ITK__IsoImageSlicer__h__
210
211 // eof - $RCSfile$