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