]> Creatis software - FrontAlgorithms.git/blob - lib/fpa/VTK/Image/PolyLineParametricPathToPolyDataFilter.hxx
ad995aab2d08190f73b0ce6fccb4fb565f03fe03
[FrontAlgorithms.git] / lib / fpa / VTK / Image / PolyLineParametricPathToPolyDataFilter.hxx
1 // =========================================================================
2 // @author Leonardo Florez Valencia
3 // @email florez-l@javeriana.edu.co
4 // =========================================================================
5 #ifndef __fpa__VTK__Image__PolyLineParametricPathToPolyDataFilter__hxx__
6 #define __fpa__VTK__Image__PolyLineParametricPathToPolyDataFilter__hxx__
7
8 #include <vtkCellArray.h>
9 #include <vtkInformation.h>
10 #include <vtkInformationVector.h>
11 #include <vtkPointData.h>
12 #include <vtkUnsignedIntArray.h>
13 #include <vtkSmartPointer.h>
14
15 // -------------------------------------------------------------------------
16 template< class _TPath >
17 typename fpa::VTK::Image::PolyLineParametricPathToPolyDataFilter< _TPath >::
18 Self* fpa::VTK::Image::PolyLineParametricPathToPolyDataFilter< _TPath >::
19 New( )
20 {
21   return( new Self( ) );
22 }
23
24 // -------------------------------------------------------------------------
25 template< class _TPath >
26 const typename
27 fpa::VTK::Image::PolyLineParametricPathToPolyDataFilter< _TPath >::
28 TPath* fpa::VTK::Image::PolyLineParametricPathToPolyDataFilter< _TPath >::
29 GetInput( ) const
30 {
31   return( this->m_PolyLineParametricPath );
32 }
33
34 // -------------------------------------------------------------------------
35 template< class _TPath >
36 void fpa::VTK::Image::PolyLineParametricPathToPolyDataFilter< _TPath >::
37 SetInput( const TPath* path )
38 {
39   if( this->m_PolyLineParametricPath != path )
40   {
41     this->m_PolyLineParametricPath = path;
42     this->Modified( );
43
44   } // fi
45 }
46
47 // -------------------------------------------------------------------------
48 template< class _TPath >
49 fpa::VTK::Image::PolyLineParametricPathToPolyDataFilter< _TPath >::
50 PolyLineParametricPathToPolyDataFilter( )
51   : vtkPolyDataAlgorithm( ),
52     m_PolyLineParametricPath( NULL )
53 {
54   this->SetNumberOfInputPorts( 0 );
55 }
56
57 // -------------------------------------------------------------------------
58 template< class _TPath >
59 fpa::VTK::Image::PolyLineParametricPathToPolyDataFilter< _TPath >::
60 ~PolyLineParametricPathToPolyDataFilter( )
61 {
62 }
63
64 // -------------------------------------------------------------------------
65 template< class _TPath >
66 int fpa::VTK::Image::PolyLineParametricPathToPolyDataFilter< _TPath >::
67 RequestData(
68   vtkInformation* information,
69   vtkInformationVector** input,
70   vtkInformationVector* output
71   )
72 {
73   static const unsigned int dim = TPath::PathDimension;
74   if( this->m_PolyLineParametricPath == NULL )
75     return( 0 );
76
77   // Get output
78   vtkInformation* info = output->GetInformationObject( 0 );
79   vtkPolyData* out = vtkPolyData::SafeDownCast(
80     info->Get( vtkDataObject::DATA_OBJECT( ) )
81     );
82
83   // Prepare data
84   out->SetPoints( vtkSmartPointer< vtkPoints >::New( ) );
85   out->SetVerts( vtkSmartPointer< vtkCellArray >::New( ) );
86   out->SetLines( vtkSmartPointer< vtkCellArray >::New( ) );
87   out->SetPolys( vtkSmartPointer< vtkCellArray >::New( ) );
88   out->SetStrips( vtkSmartPointer< vtkCellArray >::New( ) );
89   vtkSmartPointer< vtkUnsignedIntArray > darray =
90     vtkSmartPointer< vtkUnsignedIntArray >::New( );
91   darray->SetNumberOfComponents( 1 );
92   out->GetPointData( )->SetScalars( darray );
93   vtkPoints* points = out->GetPoints( );
94   vtkCellArray* lines = out->GetLines( );
95
96   // Assign all data
97   const TPath* path = this->GetInput( );
98   for( unsigned long i = 0; i < path->GetSize( ); ++i )
99   {
100     auto pnt = path->GetPoint( i );
101     if( dim == 1 )
102       points->InsertNextPoint( pnt[ 0 ], 0, 0 );
103     else if( dim == 2 )
104       points->InsertNextPoint( pnt[ 0 ], pnt[ 1 ], 0 );
105     else
106       points->InsertNextPoint( pnt[ 0 ], pnt[ 1 ], pnt[ 2 ] );
107     darray->InsertNextTuple1( double( i ) );
108     if( i > 0 )
109     {
110       lines->InsertNextCell( 2 );
111       lines->InsertCellPoint( points->GetNumberOfPoints( ) - 2 );
112       lines->InsertCellPoint( points->GetNumberOfPoints( ) - 1 );
113
114     } // fi
115
116   } // rof
117   return( 1 );
118 }
119
120 // -------------------------------------------------------------------------
121 template< class _TPath >
122 int fpa::VTK::Image::PolyLineParametricPathToPolyDataFilter< _TPath >::
123 RequestInformation(
124   vtkInformation* information,
125   vtkInformationVector** input,
126   vtkInformationVector* output
127   )
128 {
129   return( 1 );
130 }
131
132 #endif // __fpa__VTK__Image__PolyLineParametricPathToPolyDataFilterFilter__hxx__
133 // eof - $RCSfile$