]> Creatis software - cpPlugins.git/blob - lib/ivq/VTK/ImageToFourierSeriesFilter.cxx
...
[cpPlugins.git] / lib / ivq / VTK / ImageToFourierSeriesFilter.cxx
1 // =========================================================================
2 // @author Leonardo Florez Valencia
3 // @email florez-l@javeriana.edu.co
4 // =========================================================================
5
6 #include <ivq/VTK/ImageToFourierSeriesFilter.h>
7 #include <ivq/VTK/PolyDataToFourierSeriesFilter.h>
8 #include <vtkImageData.h>
9 #include <vtkInformation.h>
10 #include <vtkInformationVector.h>
11 #include <vtkMarchingSquares.h>
12 #include <vtkPolyDataConnectivityFilter.h>
13 #include <vtkSmartPointer.h>
14
15 // -------------------------------------------------------------------------
16 template< class _TFourierSeries >
17 typename ivq::VTK::ImageToFourierSeriesFilter< _TFourierSeries >::
18 Self* ivq::VTK::ImageToFourierSeriesFilter< _TFourierSeries >::
19 New( )
20 {
21   return( new Self( ) );
22 }
23
24 // -------------------------------------------------------------------------
25 template< class _TFourierSeries >
26 const _TFourierSeries&
27 ivq::VTK::ImageToFourierSeriesFilter< _TFourierSeries >::
28 GetOutput( ) const
29 {
30   return( this->m_FourierSeries );
31 }
32
33 // -------------------------------------------------------------------------
34 template< class _TFourierSeries >
35 void ivq::VTK::ImageToFourierSeriesFilter< _TFourierSeries >::
36 SetNumberOfHarmonics( unsigned int q )
37 {
38   if( this->m_FourierSeries.GetNumberOfHarmonics( ) != q )
39   {
40     this->m_FourierSeries.SetNumberOfHarmonics( q );
41     this->Modified( );
42
43   } // fi
44 }
45
46 // -------------------------------------------------------------------------
47 template< class _TFourierSeries >
48 unsigned int ivq::VTK::ImageToFourierSeriesFilter< _TFourierSeries >::
49 GetNumberOfHarmonics( ) const
50 {
51   return( this->m_FourierSeries.GetNumberOfHarmonics( ) );
52 }
53
54 // -------------------------------------------------------------------------
55 template< class _TFourierSeries >
56 const double& ivq::VTK::ImageToFourierSeriesFilter< _TFourierSeries >::
57 GetContourValue( ) const
58 {
59   return( this->m_ContourValue );
60 }
61
62 // -------------------------------------------------------------------------
63 template< class _TFourierSeries >
64 void ivq::VTK::ImageToFourierSeriesFilter< _TFourierSeries >::
65 SetContourValue( const double& v )
66 {
67   if( this->m_ContourValue != v )
68   {
69     this->m_ContourValue = v;
70     this->Modified( );
71
72   } // fi
73 }
74
75 // -------------------------------------------------------------------------
76 template< class _TFourierSeries >
77 ivq::VTK::ImageToFourierSeriesFilter< _TFourierSeries >::
78 ImageToFourierSeriesFilter( )
79   : vtkImageAlgorithm( ),
80     m_ContourValue( double( 0 ) )
81 {
82   this->SetNumberOfInputPorts( 1 );
83   this->SetNumberOfOutputPorts( 0 );
84   this->m_FourierSeries.SetNumberOfHarmonics( 1 );
85 }
86
87 // -------------------------------------------------------------------------
88 template< class _TFourierSeries >
89 ivq::VTK::ImageToFourierSeriesFilter< _TFourierSeries >::
90 ~ImageToFourierSeriesFilter( )
91 {
92 }
93
94 // -------------------------------------------------------------------------
95 template< class _TFourierSeries >
96 int ivq::VTK::ImageToFourierSeriesFilter< _TFourierSeries >::
97 RequestData(
98   vtkInformation* information,
99   vtkInformationVector** input,
100   vtkInformationVector* output
101   )
102 {
103   // Get input polydata
104   vtkInformation* inInfo = input[ 0 ]->GetInformationObject( 0 );
105   vtkImageData* image =
106     vtkImageData::SafeDownCast(
107       inInfo->Get( vtkDataObject::DATA_OBJECT( ) )
108       );
109
110   // Extract contour
111   vtkSmartPointer< vtkMarchingSquares > ms =
112     vtkSmartPointer< vtkMarchingSquares >::New( );
113   ms->SetInputData( image );
114   ms->SetValue( 0, this->m_ContourValue );
115
116   vtkSmartPointer< vtkPolyDataConnectivityFilter > conn =
117     vtkSmartPointer< vtkPolyDataConnectivityFilter >::New( );
118   conn->SetInputConnection( ms->GetOutputPort( ) );
119   conn->SetExtractionModeToClosestPointRegion( );
120   conn->SetClosestPoint( 0, 0, 0 );
121
122   // Compute fourier series
123   typedef ivq::VTK::PolyDataToFourierSeriesFilter< TFourierSeries > _TFilter;
124   vtkSmartPointer< _TFilter > filter = vtkSmartPointer< _TFilter >::New( );
125   filter->SetInputConnection( conn->GetOutputPort( ) );
126   filter->SetNumberOfHarmonics( this->m_FourierSeries.GetNumberOfHarmonics( ) );
127   filter->Update( );
128   this->m_FourierSeries = filter->GetOutput( );
129
130   return( 1 );
131 }
132
133 // -------------------------------------------------------------------------
134 template< class _TFourierSeries >
135 int ivq::VTK::ImageToFourierSeriesFilter< _TFourierSeries >::
136 RequestInformation(
137   vtkInformation* information,
138   vtkInformationVector** input,
139   vtkInformationVector* output
140   )
141 {
142   return( 1 );
143 }
144
145 // -------------------------------------------------------------------------
146 #include <ivq/ITK/FourierSeries.h>
147
148 template class IVQ_EXPORT ivq::VTK::ImageToFourierSeriesFilter< ivq::ITK::FourierSeries< float, 2 > >;
149 template class IVQ_EXPORT ivq::VTK::ImageToFourierSeriesFilter< ivq::ITK::FourierSeries< float, 3 > >;
150 template class IVQ_EXPORT ivq::VTK::ImageToFourierSeriesFilter< ivq::ITK::FourierSeries< double, 2 > >;
151 template class IVQ_EXPORT ivq::VTK::ImageToFourierSeriesFilter< ivq::ITK::FourierSeries< double, 3 > >;
152
153 // eof - $RCSfile$