]> Creatis software - FrontAlgorithms.git/blob - lib/fpa/VTK/Image3DObserver.hxx
Big bug stamped on wall
[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::TEvent                _TEvent;
33   typedef typename F::TFrontEvent           _TFrontEvent;
34   typedef typename F::TMarkEvent            _TMarkEvent;
35   typedef typename F::TCollisionEvent       _TCollisionEvent;
36   typedef typename F::TEndEvent             _TEndEvent;
37   typedef typename F::TBacktrackingEvent    _TBacktrackingEvent;
38   typedef typename F::TEndBacktrackingEvent _TEndBacktrackingEvent;
39
40   // Check inputs
41   if( this->m_RenderWindow == NULL )
42     return;
43   vtkRenderer* ren =
44     this->m_RenderWindow->GetRenderers( )->GetFirstRenderer( );
45   if( ren == NULL )
46     return;
47   const F* filter = dynamic_cast< const F* >( c );
48   if( filter == NULL )
49     return;
50   const _TImage* image = filter->GetInput( );
51   if( image == NULL )
52     return;
53
54   // Create actor
55   _TImage::RegionType reg = image->GetLargestPossibleRegion( );
56   _TImage::SizeType siz = reg.GetSize( );
57   if( this->m_Data.GetPointer( ) == NULL )
58   {
59     this->m_Data   = vtkSmartPointer< vtkPolyData >::New( );
60     this->m_Mapper = vtkSmartPointer< vtkPolyDataMapper >::New( );
61     this->m_Actor  = vtkSmartPointer< vtkActor >::New( );
62
63     vtkSmartPointer< vtkPoints > points =
64       vtkSmartPointer< vtkPoints >::New( );
65     vtkSmartPointer< vtkCellArray > vertices =
66       vtkSmartPointer< vtkCellArray >::New( );
67     vtkSmartPointer< vtkFloatArray > scalars =
68       vtkSmartPointer< vtkFloatArray >::New( );
69     this->m_Data->SetPoints( points );
70     this->m_Data->SetVerts( vertices );
71     this->m_Data->GetPointData( )->SetScalars( scalars );
72
73     this->m_Mapper->SetInputData( this->m_Data );
74     this->m_Actor->SetMapper( this->m_Mapper );
75     ren->AddActor( this->m_Actor );
76
77     this->m_Marks = TMarks::New( );
78     this->m_Marks->SetLargestPossibleRegion(
79       image->GetLargestPossibleRegion( )
80       );
81     this->m_Marks->SetRequestedRegion( image->GetRequestedRegion( ) );
82     this->m_Marks->SetBufferedRegion( image->GetBufferedRegion( ) );
83     this->m_Marks->SetOrigin( image->GetOrigin( ) );
84     this->m_Marks->SetSpacing( image->GetSpacing( ) );
85     this->m_Marks->SetDirection( image->GetDirection( ) );
86     this->m_Marks->Allocate( );
87     this->m_Marks->FillBuffer( -1 );
88     this->m_Count = 0;
89     this->m_RenderCount = reg.GetNumberOfPixels( );
90
91   } // fi
92
93   // Get possible events
94   const _TEvent* evt = dynamic_cast< const _TEvent* >( &e );
95   _TFrontEvent fevt;
96   _TMarkEvent mevt;
97   _TCollisionEvent cevt;
98   _TEndEvent eevt;
99
100   if( evt != NULL )
101   {
102     if( fevt.CheckEvent( evt ) )
103     {
104       if( this->m_Marks->GetPixel( evt->Node.Vertex ) == -1 )
105       {
106         typename _TImage::PointType pnt;
107         image->TransformIndexToPhysicalPoint( evt->Node.Vertex, pnt );
108
109         long pId = this->m_Data->GetNumberOfPoints( );
110         this->m_Data->GetVerts( )->InsertNextCell( 1 );
111         this->m_Data->GetVerts( )->InsertCellPoint( pId );
112         this->m_Data->GetPoints( )->
113           InsertNextPoint( pnt[ 0 ], pnt[ 1 ], pnt[ 2 ] );
114         this->m_Data->GetPointData( )->
115           GetScalars( )->InsertNextTuple1( 0.5 );
116         this->m_Data->Modified( );
117
118         this->m_Marks->SetPixel( evt->Node.Vertex, pId );
119         this->m_Count++;
120
121         // Render visual debug
122         double per = double( this->m_RenderCount ) * this->m_RenderPercentage;
123         if( double( this->m_Count ) >= per )
124           this->m_RenderWindow->Render( );
125         if( double( this->m_Count ) >= per )
126           this->m_Count = 0;
127
128       } // fi
129       return;
130     }
131     else if( mevt.CheckEvent( evt ) )
132     {
133       // TODO: std::cout << "mark" << std::endl;
134       return;
135
136     } // fi
137
138     if( cevt.CheckEvent( evt ) )
139     {
140       // TODO: std::cout << "collision" << std::endl;
141       return;
142
143     } // fi
144
145     if( eevt.CheckEvent( evt ) )
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   else
159   {
160     const _TBacktrackingEvent* bevt =
161       dynamic_cast< const _TBacktrackingEvent* >( &e );
162     const _TEndBacktrackingEvent* ebevt =
163       dynamic_cast< const _TEndBacktrackingEvent* >( &e );
164     if( bevt != NULL )
165     {
166       static const unsigned long nColors = 10;
167       double back_id = double( bevt->BackId % nColors ) / double( nColors );
168       typename _TImage::PointType pnt;
169       image->TransformIndexToPhysicalPoint( bevt->Node, 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( ebevt != NULL )
184     {
185       this->m_RenderWindow->Render( );
186       return;
187
188     } // fi
189
190   } // fi
191 }
192
193 // -------------------------------------------------------------------------
194 template< class F, class R >
195 fpa::VTK::Image3DObserver< F, R >::
196 Image3DObserver( )
197   : Superclass( ),
198     m_RenderWindow( NULL ),
199     m_RenderPercentage( double( 0.0001 ) )
200 {
201 }
202
203 #endif // __FPA__VTK__IMAGE3DOBSERVER__HXX__
204
205 // eof - $RCSfile$