]> Creatis software - cpPlugins.git/blob - lib/cpPlugins/Interface/Image.hxx
...
[cpPlugins.git] / lib / cpPlugins / Interface / Image.hxx
1 #ifndef __CPPLUGINS__INTERFACE__IMAGE__HXX__
2 #define __CPPLUGINS__INTERFACE__IMAGE__HXX__
3
4 #include <cpPlugins/Interface/Macros.h>
5
6 #include <itkImage.h>
7 #include <itkImageToVTKImageFilter.h>
8
9 #include <itkCovariantVector.h>
10 #include <itkDiffusionTensor3D.h>
11 #include <itkPoint.h>
12 #include <itkRGBPixel.h>
13 #include <itkRGBAPixel.h>
14 #include <itkSymmetricSecondRankTensor.h>
15 #include <itkVector.h>
16
17 // -------------------------------------------------------------------------
18 #define cpPlugins_Image_Demangle( TI, T, D, o )                 \
19   if( typeid( typename TI::PixelType ) == typeid( T ) )         \
20     this->_ITK_2_VTK< T, D >( o )
21
22 // -------------------------------------------------------------------------
23 #define cpPlugins_ImageArray_Demangle( TI, T, P, DP, DI, o )           \
24   if( typeid( typename TI::PixelType ) == typeid( T< P, DP > ) )       \
25     this->_ITK_2_VTK< T< P, DP >, DI >( o )
26
27 // -------------------------------------------------------------------------
28 #define cpPlugins_Image_Demangle_Dimension( TI, D, o )                  \
29   cpPlugins_Image_Demangle( TI, char, D, o );                           \
30   else cpPlugins_Image_Demangle( TI, short, D, o );                     \
31   else cpPlugins_Image_Demangle( TI, int, D, o );                       \
32   else cpPlugins_Image_Demangle( TI, long, D, o );                      \
33   else cpPlugins_Image_Demangle( TI, float, D, o );                     \
34   else cpPlugins_Image_Demangle( TI, double, D, o );                    \
35   else cpPlugins_Image_Demangle( TI, unsigned char, D, o );             \
36   else cpPlugins_Image_Demangle( TI, unsigned short, D, o );            \
37   else cpPlugins_Image_Demangle( TI, unsigned int, D, o );              \
38   else cpPlugins_Image_Demangle( TI, unsigned long, D, o );             \
39   else cpPlugins_Image_Demangle( TI, itk::RGBPixel< char >, D, o );     \
40   else cpPlugins_Image_Demangle( TI, itk::RGBPixel< short >, D, o );    \
41   else cpPlugins_Image_Demangle(                                        \
42     TI, itk::RGBPixel< unsigned char >, D, o                            \
43     );                                                                  \
44   else cpPlugins_Image_Demangle(                                        \
45     TI, itk::RGBPixel< unsigned short >, D, o                           \
46     );                                                                  \
47   else cpPlugins_Image_Demangle( TI, itk::RGBAPixel< char >, D, o );    \
48   else cpPlugins_Image_Demangle( TI, itk::RGBAPixel< short >, D, o );   \
49   else cpPlugins_Image_Demangle(                                        \
50     TI, itk::RGBAPixel< unsigned char >, D, o                           \
51     );                                                                  \
52   else cpPlugins_Image_Demangle(                                        \
53     TI, itk::RGBAPixel< unsigned short >, D, o                          \
54     );                                                                  \
55   else cpPlugins_ImageArray_Demangle(                                   \
56     TI, itk::Vector, float, D, D, o                                     \
57     );                                                                  \
58   else cpPlugins_ImageArray_Demangle(                                   \
59     TI, itk::Vector, double, D, D, o                                    \
60     );                                                                  \
61   else cpPlugins_ImageArray_Demangle(                                   \
62     TI, itk::CovariantVector, float, D, D, o                            \
63     );                                                                  \
64   else cpPlugins_ImageArray_Demangle(                                   \
65     TI, itk::CovariantVector, double, D, D, o                           \
66     );                                                                  \
67   else cpPlugins_ImageArray_Demangle( TI, itk::Point, float, D, D, o ); \
68   else cpPlugins_ImageArray_Demangle(                                   \
69     TI, itk::Point, double, D, D, o                                     \
70     );                                                                  \
71   else cpPlugins_ImageArray_Demangle(                                   \
72     TI, itk::SymmetricSecondRankTensor, float, D, D, o                  \
73     );                                                                  \
74   else cpPlugins_ImageArray_Demangle(                                   \
75     TI, itk::SymmetricSecondRankTensor, double, D, D, o                 \
76     )
77
78 // -------------------------------------------------------------------------
79 #define cpPlugins_Image_Import( T, D )                                  \
80   cpPlugins_TEMPLATE_IMPORT(                                            \
81     2(class cpPlugins_Interface_EXPORT itk::Image< T, D >)              \
82     )
83
84 // -------------------------------------------------------------------------
85 #define cpPlugins_VTKImage_Import( T, D )                               \
86   cpPlugins_TEMPLATE_IMPORT(                                            \
87     2(                                                                  \
88       class cpPlugins_Interface_EXPORT                                  \
89       itk::ImageToVTKImageFilter< itk::Image< T, D > >                  \
90       )                                                                 \
91     )
92
93 // -------------------------------------------------------------------------
94 #define cpPlugins_ImageArray_Import( A, T, DA, DI )                     \
95   cpPlugins_TEMPLATE_IMPORT(                                            \
96     3(                                                                  \
97       class cpPlugins_Interface_EXPORT                                  \
98       itk::Image< itk::A< T, DA >, DI >                                 \
99       )                                                                 \
100     )
101
102 // -------------------------------------------------------------------------
103 #define cpPlugins_VTKImageArray_Import( A, T, DA, DI )                  \
104   cpPlugins_TEMPLATE_IMPORT(                                            \
105     3(                                                                  \
106       class cpPlugins_Interface_EXPORT                                  \
107       itk::ImageToVTKImageFilter< itk::Image< itk::A< T, DA >, DI > >   \
108       )                                                                 \
109     )
110
111 // -------------------------------------------------------------------------
112 #define cpPlugins_Image_Import_AllDimensions( T )       \
113   cpPlugins_Image_Import( T, 2 );                       \
114   cpPlugins_Image_Import( T, 3 );                       \
115   cpPlugins_Image_Import( T, 4 )
116
117 // -------------------------------------------------------------------------
118 #define cpPlugins_VTKImage_Import_AllDimensions( T )       \
119   cpPlugins_VTKImage_Import( T, 2 );                       \
120   cpPlugins_VTKImage_Import( T, 3 )
121
122 // -------------------------------------------------------------------------
123 #define cpPlugins_ImageArray_Import_AllDimensions( A, T )             \
124   cpPlugins_ImageArray_Import( A, T, 2, 2 );                          \
125   cpPlugins_ImageArray_Import( A, T, 3, 3 );                          \
126   cpPlugins_ImageArray_Import( A, T, 4, 4 )
127
128 // -------------------------------------------------------------------------
129 #define cpPlugins_VTKImageArray_Import_AllDimensions( A, T )            \
130   cpPlugins_VTKImageArray_Import( A, T, 2, 2 );                         \
131   cpPlugins_VTKImageArray_Import( A, T, 3, 3 )
132
133 // -------------------------------------------------------------------------
134 // ITK base clases
135
136 #ifndef cpPlugins_Interface_EXPORTS
137
138 cpPlugins_Image_Import_AllDimensions( char );
139 cpPlugins_Image_Import_AllDimensions( short );
140 cpPlugins_Image_Import_AllDimensions( int );
141 cpPlugins_Image_Import_AllDimensions( long );
142 cpPlugins_Image_Import_AllDimensions( unsigned char );
143 cpPlugins_Image_Import_AllDimensions( unsigned short );
144 cpPlugins_Image_Import_AllDimensions( unsigned int );
145 cpPlugins_Image_Import_AllDimensions( unsigned long );
146 cpPlugins_Image_Import_AllDimensions( float );
147 cpPlugins_Image_Import_AllDimensions( double );
148
149 cpPlugins_Image_Import_AllDimensions( std::complex< float > );
150 cpPlugins_Image_Import_AllDimensions( std::complex< double > );
151
152 cpPlugins_Image_Import_AllDimensions( itk::RGBPixel< char > );
153 cpPlugins_Image_Import_AllDimensions( itk::RGBPixel< short > );
154 cpPlugins_Image_Import_AllDimensions( itk::RGBPixel< int > );
155 cpPlugins_Image_Import_AllDimensions( itk::RGBPixel< long > );
156 cpPlugins_Image_Import_AllDimensions( itk::RGBPixel< unsigned char > );
157 cpPlugins_Image_Import_AllDimensions( itk::RGBPixel< unsigned short > );
158 cpPlugins_Image_Import_AllDimensions( itk::RGBPixel< unsigned int > );
159 cpPlugins_Image_Import_AllDimensions( itk::RGBPixel< unsigned long > );
160
161 cpPlugins_Image_Import_AllDimensions( itk::RGBAPixel< char > );
162 cpPlugins_Image_Import_AllDimensions( itk::RGBAPixel< short > );
163 cpPlugins_Image_Import_AllDimensions( itk::RGBAPixel< int > );
164 cpPlugins_Image_Import_AllDimensions( itk::RGBAPixel< long > );
165 cpPlugins_Image_Import_AllDimensions( itk::RGBAPixel< unsigned char > );
166 cpPlugins_Image_Import_AllDimensions( itk::RGBAPixel< unsigned short > );
167 cpPlugins_Image_Import_AllDimensions( itk::RGBAPixel< unsigned int > );
168 cpPlugins_Image_Import_AllDimensions( itk::RGBAPixel< unsigned long > );
169
170 cpPlugins_ImageArray_Import_AllDimensions( Vector, float );
171 cpPlugins_ImageArray_Import_AllDimensions( Vector, double );
172
173 cpPlugins_ImageArray_Import_AllDimensions( CovariantVector, float );
174 cpPlugins_ImageArray_Import_AllDimensions( CovariantVector, double );
175
176 cpPlugins_ImageArray_Import_AllDimensions( Point, float );
177 cpPlugins_ImageArray_Import_AllDimensions( Point, double );
178
179 cpPlugins_ImageArray_Import_AllDimensions(
180   SymmetricSecondRankTensor, float
181   );
182 cpPlugins_ImageArray_Import_AllDimensions(
183   SymmetricSecondRankTensor, double
184   );
185
186 cpPlugins_Image_Import( itk::DiffusionTensor3D< float >, 3 );
187 cpPlugins_Image_Import( itk::DiffusionTensor3D< double >, 3 );
188 cpPlugins_Image_Import( itk::Offset< 2 >, 2 );
189 cpPlugins_Image_Import( itk::Offset< 3 >, 3 );
190 cpPlugins_Image_Import( itk::Offset< 4 >, 4 );
191
192 // -------------------------------------------------------------------------
193 // ITK<->VTK base clases
194 cpPlugins_VTKImage_Import_AllDimensions( char );
195 cpPlugins_VTKImage_Import_AllDimensions( short );
196 cpPlugins_VTKImage_Import_AllDimensions( int );
197 cpPlugins_VTKImage_Import_AllDimensions( long );
198 cpPlugins_VTKImage_Import_AllDimensions( unsigned char );
199 cpPlugins_VTKImage_Import_AllDimensions( unsigned short );
200 cpPlugins_VTKImage_Import_AllDimensions( unsigned int );
201 cpPlugins_VTKImage_Import_AllDimensions( unsigned long );
202 cpPlugins_VTKImage_Import_AllDimensions( float );
203 cpPlugins_VTKImage_Import_AllDimensions( double );
204
205 cpPlugins_VTKImage_Import_AllDimensions( itk::RGBPixel< char > );
206 cpPlugins_VTKImage_Import_AllDimensions( itk::RGBPixel< short > );
207 cpPlugins_VTKImage_Import_AllDimensions( itk::RGBPixel< int > );
208 cpPlugins_VTKImage_Import_AllDimensions( itk::RGBPixel< long > );
209 cpPlugins_VTKImage_Import_AllDimensions( itk::RGBPixel< unsigned char > );
210 cpPlugins_VTKImage_Import_AllDimensions( itk::RGBPixel< unsigned short > );
211 cpPlugins_VTKImage_Import_AllDimensions( itk::RGBPixel< unsigned int > );
212 cpPlugins_VTKImage_Import_AllDimensions( itk::RGBPixel< unsigned long > );
213 cpPlugins_VTKImage_Import_AllDimensions( itk::RGBPixel< float > );
214 cpPlugins_VTKImage_Import_AllDimensions( itk::RGBPixel< double > );
215
216 cpPlugins_VTKImage_Import_AllDimensions( itk::RGBAPixel< char > );
217 cpPlugins_VTKImage_Import_AllDimensions( itk::RGBAPixel< short > );
218 cpPlugins_VTKImage_Import_AllDimensions( itk::RGBAPixel< int > );
219 cpPlugins_VTKImage_Import_AllDimensions( itk::RGBAPixel< long > );
220 cpPlugins_VTKImage_Import_AllDimensions( itk::RGBAPixel< unsigned char > );
221 cpPlugins_VTKImage_Import_AllDimensions( itk::RGBAPixel< unsigned short > );
222 cpPlugins_VTKImage_Import_AllDimensions( itk::RGBAPixel< unsigned int > );
223 cpPlugins_VTKImage_Import_AllDimensions( itk::RGBAPixel< unsigned long > );
224 cpPlugins_VTKImage_Import_AllDimensions( itk::RGBAPixel< float > );
225 cpPlugins_VTKImage_Import_AllDimensions( itk::RGBAPixel< double > );
226
227 cpPlugins_VTKImageArray_Import_AllDimensions( Vector, float );
228 cpPlugins_VTKImageArray_Import_AllDimensions( Vector, double );
229
230 cpPlugins_VTKImageArray_Import_AllDimensions( CovariantVector, float );
231 cpPlugins_VTKImageArray_Import_AllDimensions( CovariantVector, double );
232
233 cpPlugins_VTKImageArray_Import_AllDimensions( Point, float );
234 cpPlugins_VTKImageArray_Import_AllDimensions( Point, double );
235
236 cpPlugins_VTKImageArray_Import_AllDimensions(
237   SymmetricSecondRankTensor, float
238   );
239 cpPlugins_VTKImageArray_Import_AllDimensions(
240   SymmetricSecondRankTensor, double
241   );
242
243 cpPlugins_VTKImage_Import( itk::DiffusionTensor3D< float >, 3 );
244 cpPlugins_VTKImage_Import( itk::DiffusionTensor3D< double >, 3 );
245
246 #endif // cpPlugins_Interface_EXPORTS
247
248 // -------------------------------------------------------------------------
249 template< class I >
250 void cpPlugins::Interface::Image::
251 SetITKImage( itk::DataObject* object )
252 {
253   // Check if input object has the desired type
254   I* image = dynamic_cast< I* >( object );
255   if( image == NULL )
256   {
257     this->m_ITKObject = NULL;
258     this->m_VTKObject = NULL;
259     this->m_ITKvVTKConnection = NULL;
260     this->Modified( );
261
262   } // fi
263
264   // Connect it to VTK
265   if( I::ImageDimension == 2 )
266   {
267     cpPlugins_Image_Demangle_Dimension( I, 2, object );
268     else
269     {
270       this->m_VTKObject = NULL;
271       this->m_ITKvVTKConnection = NULL;
272
273     } // fi
274   }
275   else if( I::ImageDimension == 3 )
276   {
277     cpPlugins_Image_Demangle_Dimension( I, 3, object );
278     else cpPlugins_Image_Demangle(
279       I, itk::DiffusionTensor3D< float >, 3, object
280       );
281     else cpPlugins_Image_Demangle(
282       I, itk::DiffusionTensor3D< double >, 3, object
283       );
284     else
285     {
286       this->m_VTKObject = NULL;
287       this->m_ITKvVTKConnection = NULL;
288
289     } // fi
290   }
291   else
292   {
293     this->m_VTKObject = NULL;
294     this->m_ITKvVTKConnection = NULL;
295
296   } // fi
297
298   // Keep objects
299   this->m_ITKObject = object;
300   this->Modified( );
301 }
302
303 // -------------------------------------------------------------------------
304 template< class I >
305 I* cpPlugins::Interface::Image::
306 GetITKImage( )
307 {
308   return( dynamic_cast< I* >( this->m_ITKObject.GetPointer( ) ) );
309 }
310
311 // -------------------------------------------------------------------------
312 template< class I >
313 const I* cpPlugins::Interface::Image::
314 GetITKImage( ) const
315 {
316   return( dynamic_cast< const I* >( this->m_ITKObject.GetPointer( ) ) );
317 }
318
319 // -------------------------------------------------------------------------
320 template< class P, unsigned int D >
321 void cpPlugins::Interface::Image::
322 _ITK_2_VTK( itk::DataObject* object )
323 {
324   typedef itk::Image< P, D > _I;
325
326   // Check if input object has the desired type
327   _I* image = dynamic_cast< _I* >( object );
328   if( image == NULL )
329     return;
330
331   // Connect it to VTK
332   typename itk::ImageToVTKImageFilter< _I >::Pointer f =
333     itk::ImageToVTKImageFilter< _I >::New( );
334   f->SetInput( image );
335   f->Update( );
336
337   // Keep objects
338   this->m_VTKObject = f->GetOutput( );
339   this->m_ITKvVTKConnection = f;
340 }
341
342 #endif // __CPPLUGINS__INTERFACE__IMAGE__HXX__
343
344 // eof - $RCSfile$