]> Creatis software - FrontAlgorithms.git/blob - lib/fpa/Image/SkeletonToPolyDataFilter.hxx
c9d1136e99f5d3aa8ad78321d2d5cb4a038d0ecf
[FrontAlgorithms.git] / lib / fpa / Image / SkeletonToPolyDataFilter.hxx
1 // =========================================================================
2 // @author Leonardo Florez Valencia
3 // @email florez-l@javeriana.edu.co
4 // =========================================================================
5
6 #ifndef __fpa__Image__SkeletonToPolyDataFilter__hxx__
7 #define __fpa__Image__SkeletonToPolyDataFilter__hxx__
8
9 #ifdef USE_VTK
10 #  include <vtkCellArray.h>
11 #  include <vtkInformation.h>
12 #  include <vtkInformationVector.h>
13 #  include <vtkPointData.h>
14 #  include <vtkUnsignedIntArray.h>
15 #  include <vtkSmartPointer.h>
16 #endif // USE_VTK
17
18 // -------------------------------------------------------------------------
19 template< class _TSkeleton >
20 typename fpa::Image::SkeletonToPolyDataFilter< _TSkeleton >::
21 Self* fpa::Image::SkeletonToPolyDataFilter< _TSkeleton >::
22 New( )
23 {
24   return( new Self( ) );
25 }
26
27 // -------------------------------------------------------------------------
28 template< class _TSkeleton >
29 const typename
30 fpa::Image::SkeletonToPolyDataFilter< _TSkeleton >::
31 TSkeleton* fpa::Image::SkeletonToPolyDataFilter< _TSkeleton >::
32 GetInput( ) const
33 {
34   return( this->m_Skeleton );
35 }
36
37 // -------------------------------------------------------------------------
38 template< class _TSkeleton >
39 void fpa::Image::SkeletonToPolyDataFilter< _TSkeleton >::
40 SetInput( const TSkeleton* sk )
41 {
42   if( this->m_Skeleton != sk )
43   {
44     this->m_Skeleton = sk;
45 #ifdef USE_VTK
46     this->Modified( );
47 #endif // USE_VTK
48
49   } // fi
50 }
51
52 // -------------------------------------------------------------------------
53 template< class _TSkeleton >
54 fpa::Image::SkeletonToPolyDataFilter< _TSkeleton >::
55 SkeletonToPolyDataFilter( )
56 #ifdef USE_VTK
57   : vtkPolyDataAlgorithm( ),
58     m_Skeleton( NULL )
59 #endif // USE_VTK
60 {
61 #ifdef USE_VTK
62   this->SetNumberOfInputPorts( 0 );
63 #endif // USE_VTK
64 }
65
66 // -------------------------------------------------------------------------
67 template< class _TSkeleton >
68 fpa::Image::SkeletonToPolyDataFilter< _TSkeleton >::
69 ~SkeletonToPolyDataFilter( )
70 {
71 }
72
73 #ifdef USE_VTK
74 // -------------------------------------------------------------------------
75 template< class _TSkeleton >
76 int fpa::Image::SkeletonToPolyDataFilter< _TSkeleton >::
77 RequestData(
78   vtkInformation* information,
79   vtkInformationVector** input,
80   vtkInformationVector* output
81   )
82 {
83   typedef typename _TSkeleton::TPath _TPath;
84   static const unsigned int dim = _TPath::PathDimension;
85
86   if( this->m_Skeleton == NULL )
87     return( 0 );
88
89   // Get output
90   vtkInformation* info = output->GetInformationObject( 0 );
91   vtkPolyData* out = vtkPolyData::SafeDownCast(
92     info->Get( vtkDataObject::DATA_OBJECT( ) )
93     );
94
95   // Prepare data
96   out->SetPoints( vtkSmartPointer< vtkPoints >::New( ) );
97   out->SetVerts( vtkSmartPointer< vtkCellArray >::New( ) );
98   out->SetLines( vtkSmartPointer< vtkCellArray >::New( ) );
99   out->SetPolys( vtkSmartPointer< vtkCellArray >::New( ) );
100   out->SetStrips( vtkSmartPointer< vtkCellArray >::New( ) );
101   vtkSmartPointer< vtkUnsignedIntArray > darray =
102     vtkSmartPointer< vtkUnsignedIntArray >::New( );
103   darray->SetNumberOfComponents( 1 );
104   out->GetPointData( )->SetScalars( darray );
105   vtkPoints* points = out->GetPoints( );
106   vtkCellArray* lines = out->GetLines( );
107
108   // Assign all data
109   unsigned int dcount = 0;
110   typename TSkeleton::TMatrix::const_iterator  mIt = this->m_Skeleton->BeginEdgesRows( );
111   for( ; mIt != this->m_Skeleton->EndEdgesRows( ); ++mIt )
112   {
113     // TODO: mIt->first; --> this is the row index. <--
114     typename TSkeleton::TMatrixRow::const_iterator rIt = mIt->second.begin( );
115     for( ; rIt != mIt->second.end( ); ++rIt )
116     {
117       // TODO: rIt->first;  --> this is the column index.
118       typename TSkeleton::TEdges::const_iterator eIt = rIt->second.begin( );
119       for( ; eIt != rIt->second.end( ); ++eIt )
120       {
121         _TPath* path = *eIt;
122         for( unsigned long i = 0; i < path->GetSize( ); ++i )
123         {
124           auto pnt = path->GetPoint( i );
125           if( dim == 1 )
126             points->InsertNextPoint( pnt[ 0 ], 0, 0 );
127           else if( dim == 2 )
128             points->InsertNextPoint( pnt[ 0 ], pnt[ 1 ], 0 );
129           else
130             points->InsertNextPoint( pnt[ 0 ], pnt[ 1 ], pnt[ 2 ] );
131           darray->InsertNextTuple1( double( dcount ) );
132           if( i > 0 )
133           {
134             lines->InsertNextCell( 2 );
135             lines->InsertCellPoint( points->GetNumberOfPoints( ) - 2 );
136             lines->InsertCellPoint( points->GetNumberOfPoints( ) - 1 );
137
138           } // fi
139
140         } // rof
141         dcount++;
142
143       } // rof
144
145     } // rof
146
147   } // rof
148   return( 1 );
149 }
150
151 // -------------------------------------------------------------------------
152 template< class _TSkeleton >
153 int fpa::Image::SkeletonToPolyDataFilter< _TSkeleton >::
154 RequestInformation(
155   vtkInformation* information,
156   vtkInformationVector** input,
157   vtkInformationVector* output
158   )
159 {
160   return( 1 );
161 }
162 #endif // USE_VTK
163
164 #endif // __fpa__Image__SkeletonToPolyDataFilterFilter__hxx__
165
166 // eof - $RCSfile$