]> Creatis software - FrontAlgorithms.git/blob - lib/fpa/VTK/Image3DObserver.hxx
Almost there...
[FrontAlgorithms.git] / lib / fpa / VTK / Image3DObserver.hxx
1 #ifndef __FPA__VTK__IMAGE3DOBSERVER__HXX__
2 #define __FPA__VTK__IMAGE3DOBSERVER__HXX__
3
4 #include <itkImageBase.h>
5
6 #include <vtkCellArray.h>
7 #include <vtkFloatArray.h>
8 #include <vtkPointData.h>
9 #include <vtkPoints.h>
10 #include <vtkRenderer.h>
11 #include <vtkRendererCollection.h>
12
13 /*
14   #include <vtkCellArray.h>
15   #include <vtkProperty.h>
16 */
17
18 // -------------------------------------------------------------------------
19 template< class F, class R >
20 void fpa::VTK::Image3DObserver< F, R >::
21 SetRenderWindow( R* rw )
22 {
23   this->m_RenderWindow = rw;
24 }
25
26 // -------------------------------------------------------------------------
27 template< class F, class R >
28 void fpa::VTK::Image3DObserver< F, R >::
29 Execute( const itk::Object* c, const itk::EventObject& e )
30 {
31   typedef itk::ImageBase< 3 >         _TImage;
32   typedef typename F::TStartEvent     _TStartEvent;
33   typedef typename F::TStartLoopEvent _TStartLoopEvent;
34   typedef typename F::TEndEvent       _TEndEvent;
35   typedef typename F::TEndLoopEvent   _TEndLoopEvent;
36   typedef typename F::TAliveEvent     _TAliveEvent;
37   typedef typename F::TFrontEvent     _TFrontEvent;
38   typedef typename F::TFreezeEvent    _TFreezeEvent;
39
40   typedef typename F::TStartBacktrackingEvent _TStartBacktrackingEvent;
41   typedef typename F::TEndBacktrackingEvent _TEndBacktrackingEvent;
42   typedef typename F::TBacktrackingEvent _TBacktrackingEvent;
43
44   // Check inputs
45   if( this->m_RenderWindow == NULL )
46     return;
47   vtkRenderer* ren =
48     this->m_RenderWindow->GetRenderers( )->GetFirstRenderer( );
49   if( ren == NULL )
50     return;
51   const F* filter = dynamic_cast< const F* >( c );
52   if( filter == NULL )
53     return;
54   const _TImage* image = filter->GetInput( );
55   if( image == NULL )
56     return;
57
58   const _TStartEvent* startEvt = dynamic_cast< const _TStartEvent* >( &e );
59   const _TStartBacktrackingEvent* startBackEvt =
60     dynamic_cast< const _TStartBacktrackingEvent* >( &e );
61   if( startEvt != NULL || startBackEvt != NULL )
62   {
63     // Create actor
64     _TImage::RegionType reg = image->GetLargestPossibleRegion( );
65     _TImage::SizeType siz = reg.GetSize( );
66     if( this->m_Data.GetPointer( ) == NULL )
67     {
68       this->m_Data   = vtkSmartPointer< vtkPolyData >::New( );
69       this->m_Mapper = vtkSmartPointer< vtkPolyDataMapper >::New( );
70       this->m_Actor  = vtkSmartPointer< vtkActor >::New( );
71
72       vtkSmartPointer< vtkPoints > points =
73         vtkSmartPointer< vtkPoints >::New( );
74       vtkSmartPointer< vtkCellArray > vertices =
75         vtkSmartPointer< vtkCellArray >::New( );
76       vtkSmartPointer< vtkFloatArray > scalars =
77         vtkSmartPointer< vtkFloatArray >::New( );
78       this->m_Data->SetPoints( points );
79       this->m_Data->SetVerts( vertices );
80       this->m_Data->GetPointData( )->SetScalars( scalars );
81
82       this->m_Mapper->SetInputData( this->m_Data );
83       this->m_Actor->SetMapper( this->m_Mapper );
84       ren->AddActor( this->m_Actor );
85
86       this->m_Marks = TMarks::New( );
87       this->m_Marks->SetLargestPossibleRegion(
88         image->GetLargestPossibleRegion( )
89         );
90       this->m_Marks->SetRequestedRegion( image->GetRequestedRegion( ) );
91       this->m_Marks->SetBufferedRegion( image->GetBufferedRegion( ) );
92       this->m_Marks->SetOrigin( image->GetOrigin( ) );
93       this->m_Marks->SetSpacing( image->GetSpacing( ) );
94       this->m_Marks->SetDirection( image->GetDirection( ) );
95       this->m_Marks->Allocate( );
96       this->m_Marks->FillBuffer( -1 );
97       this->m_Count = 0;
98       this->m_RenderCount = reg.GetNumberOfPixels( );
99
100     } // fi
101     return;
102
103   } // fi
104
105   const _TAliveEvent* aliveEvt = dynamic_cast< const _TAliveEvent* >( &e );
106   const _TFrontEvent* frontEvt = dynamic_cast< const _TFrontEvent* >( &e );
107   if( aliveEvt != NULL || frontEvt != NULL )
108   {
109     _TImage::IndexType vertex;
110     if( aliveEvt != NULL )
111       vertex = aliveEvt->Vertex;
112     else if( frontEvt != NULL )
113       vertex = frontEvt->Vertex;
114
115     if( this->m_Marks->GetPixel( vertex ) == -1 )
116     {
117       typename _TImage::PointType pnt;
118       image->TransformIndexToPhysicalPoint( vertex, pnt );
119
120       long pId = this->m_Data->GetNumberOfPoints( );
121       this->m_Data->GetVerts( )->InsertNextCell( 1 );
122       this->m_Data->GetVerts( )->InsertCellPoint( pId );
123       this->m_Data->GetPoints( )->
124         InsertNextPoint( pnt[ 0 ], pnt[ 1 ], pnt[ 2 ] );
125       this->m_Data->GetPointData( )->
126         GetScalars( )->InsertNextTuple1( 0.5 );
127       this->m_Data->Modified( );
128
129       this->m_Marks->SetPixel( vertex, pId );
130       this->m_Count++;
131
132       // Render visual debug
133       double per = double( this->m_RenderCount ) * this->m_RenderPercentage;
134       if( double( this->m_Count ) >= per )
135         this->m_RenderWindow->Render( );
136       if( double( this->m_Count ) >= per )
137         this->m_Count = 0;
138
139     } // fi
140     return;
141
142   } // fi
143
144   const _TEndEvent* endEvt = dynamic_cast< const _TEndEvent* >( &e );
145   if( endEvt != NULL )
146   {
147     this->m_RenderWindow->Render( );
148     ren->RemoveActor( this->m_Actor );
149     this->m_RenderWindow->Render( );
150     this->m_Marks = NULL;
151     this->m_Data = NULL;
152     this->m_Mapper = NULL;
153     this->m_Actor = NULL;
154     return;
155
156   } // fi
157
158
159   const _TBacktrackingEvent* backEvt =
160     dynamic_cast< const _TBacktrackingEvent* >( &e );
161   const _TEndBacktrackingEvent* endBackEvt =
162     dynamic_cast< const _TEndBacktrackingEvent* >( &e );
163   if( backEvt != NULL )
164   {
165       static const unsigned long nColors = 10;
166       double back_id =
167         double( backEvt->FrontId % nColors ) / double( nColors );
168       typename _TImage::PointType pnt;
169       image->TransformIndexToPhysicalPoint( backEvt->Vertex, pnt );
170
171       long pId = this->m_Data->GetNumberOfPoints( );
172       this->m_Data->GetVerts( )->InsertNextCell( 1 );
173       this->m_Data->GetVerts( )->InsertCellPoint( pId );
174       this->m_Data->GetPoints( )->
175         InsertNextPoint( pnt[ 0 ], pnt[ 1 ], pnt[ 2 ] );
176       this->m_Data->GetPointData( )->
177         GetScalars( )->InsertNextTuple1( back_id );
178       this->m_Data->Modified( );
179       return;
180
181     } // fi
182
183     if( endBackEvt != NULL )
184     {
185       this->m_RenderWindow->Render( );
186       /* TODO: DEBUG
187          std::cout << "Press enter: " << std::ends;
188          int aux;
189          std::cin >> aux;
190       */
191       return;
192
193     } // fi
194 }
195
196 // -------------------------------------------------------------------------
197 template< class F, class R >
198 fpa::VTK::Image3DObserver< F, R >::
199 Image3DObserver( )
200   : Superclass( ),
201     m_RenderWindow( NULL ),
202     m_RenderPercentage( double( 0.0001 ) )
203 {
204 }
205
206 #endif // __FPA__VTK__IMAGE3DOBSERVER__HXX__
207
208 // eof - $RCSfile$