]> Creatis software - FrontAlgorithms.git/blob - lib/fpa/VTK/Image3DObserver.hxx
First commit
[FrontAlgorithms.git] / lib / fpa / VTK / Image3DObserver.hxx
1 #ifndef __FPA__VTK__IMAGE3DOBSERVER__HXX__
2 #define __FPA__VTK__IMAGE3DOBSERVER__HXX__
3
4 #include <vtkCellArray.h>
5 #include <vtkPoints.h>
6 #include <vtkProperty.h>
7 #include <vtkRenderer.h>
8 #include <vtkRendererCollection.h>
9 /*
10   #include <vtkCellArray.h>
11   #include <vtkPolyDataMapper.h>
12   #include <vtkPoints.h>
13   #include <vtkPointData.h>
14 */
15
16 // -------------------------------------------------------------------------
17 template< class F, class R >
18 void fpa::VTK::Image3DObserver< F, R >::
19 SetImage( const TImage* img, R* rw )
20 {
21   this->m_Image = img;
22   this->m_RenderWindow = rw;
23   this->m_Vertices.clear( );
24
25   /*
26     unsigned int minD = TImage::ImageDimension;
27     minD = ( minD < 3 )? minD: 3;
28
29     int e[ 6 ] = { 0 };
30     typename TImage::RegionType reg = this->m_Image->GetRequestedRegion( );
31     for( unsigned int i = 0; i < minD; i++ )
32     {
33     e[ ( i << 1 ) + 0 ] = reg.GetIndex( )[ i ];
34     e[ ( i << 1 ) + 1 ] = reg.GetIndex( )[ i ] + reg.GetSize( )[ i ] - 1;
35
36     } // rof
37
38     typename TImage::SpacingType spac = this->m_Image->GetSpacing( );
39     double s[ 3 ] = { 1, 1, 1 };
40     for( unsigned int i = 0; i < minD; i++ )
41     s[ i ] = double( spac[ i ] );
42
43     typename TImage::PointType orig = this->m_Image->GetOrigin( );
44     double o[ 3 ] = { 0 };
45     for( unsigned int i = 0; i < minD; i++ )
46     o[ i ] = double( orig[ i ] );
47
48     this->m_Stencil->SetExtent( e );
49     this->m_Stencil->SetSpacing( s );
50     this->m_Stencil->SetOrigin( o );
51     this->m_Stencil->AllocateScalars( VTK_UNSIGNED_CHAR, 4 );
52     for( unsigned int i = 0; i < 3; i++ )
53     this->m_Stencil->GetPointData( )->
54     GetScalars( )->FillComponent( i, 255 );
55     this->m_Stencil->GetPointData( )->GetScalars( )->FillComponent( 3, 0 );
56
57     this->m_StencilActor->SetInputData( this->m_Stencil );
58     this->m_StencilActor->InterpolateOff( );
59   */
60
61   double nPix =
62     double( this->m_Image->GetRequestedRegion( ).GetNumberOfPixels( ) );
63   this->m_Percentage = ( unsigned long )( nPix * 0.001 );
64   this->m_Number = 0;
65
66   if( this->m_RenderWindow != NULL )
67   {
68     vtkRenderer* ren =
69       this->m_RenderWindow->GetRenderers( )->GetFirstRenderer( );
70     ren->AddActor( this->m_PolyDataActor );
71     this->Render( );
72
73   } // fi
74 }
75
76 // -------------------------------------------------------------------------
77 template< class F, class R >
78 void fpa::VTK::Image3DObserver< F, R >::
79 Render( )
80 {
81   if( this->m_RenderWindow != NULL )
82     this->m_RenderWindow->Render( );
83 }
84
85 // -------------------------------------------------------------------------
86 template< class F, class R >
87 void fpa::VTK::Image3DObserver< F, R >::
88 Execute( const itk::Object* c, const itk::EventObject& e )
89 {
90   typedef typename TImage::IndexType        TIndex;
91   typedef typename TImage::PointType        TPoint;
92   typedef typename TFilter::TEvent          TEvent;
93   typedef typename TFilter::TFrontEvent     TFrontEvent;
94   typedef typename TFilter::TMarkEvent      TMarkEvent;
95   typedef typename TFilter::TCollisionEvent TCollisionEvent;
96   typedef typename TFilter::TEndEvent       TEndEvent;
97
98   if( this->m_RenderWindow == NULL )
99     return;
100
101   const TEvent* evt = dynamic_cast< const TEvent* >( &e );
102   TFrontEvent fevt;
103   TMarkEvent mevt;
104   TCollisionEvent cevt;
105   TEndEvent eevt;
106   if( evt != NULL )
107   {
108     typename TImage::PointType pnt;
109     this->m_Image->TransformIndexToPhysicalPoint( evt->Node.Vertex, pnt );
110
111     if( fevt.CheckEvent( evt ) )
112     {
113       if(
114         this->m_Vertices.find( evt->Node.Vertex ) == this->m_Vertices.end( )
115         )
116       {
117         unsigned long pId = this->m_PolyData->GetPoints( )->InsertNextPoint(
118           double( pnt[ 0 ] ), double( pnt[ 1 ] ), double( pnt[ 2 ] )
119           );
120         unsigned long cId = this->m_PolyData->GetVerts( )->InsertNextCell( 1 );
121         this->m_PolyData->GetVerts( )->InsertCellPoint( pId );
122         this->m_PolyData->Modified( );
123         this->m_Vertices[ evt->Node.Vertex ] = TVertexIds( pId, cId );
124
125       } // rof
126     }
127     else if( mevt.CheckEvent( evt ) )
128     {
129       /*
130         typename TVertices::iterator vIt =
131         this->m_Vertices.find( evt->Node.Vertex );
132         if( vIt != this->m_Vertices.end( ) )
133         {
134         std::cout << "Erase " << evt->Node.Vertex << " " << vIt->second.second << std::endl;
135         } // fi
136       */
137     } // fi
138
139     if( cevt.CheckEvent( evt ) )
140     {
141       // std::cout << "Collision";
142     } // fi
143     /*
144       this->SetPixel(
145       evt->Node.Vertex,
146       Colors[ evt->Node.FrontId ][ 0 ],
147       Colors[ evt->Node.FrontId ][ 1 ],
148       Colors[ evt->Node.FrontId ][ 2 ],
149       Colors[ evt->Node.FrontId ][ 3 ]
150       );
151       else if( mevt.CheckEvent( evt ) )
152       this->SetPixel( evt->Node.Vertex, 255, 0, 0, 255 );
153
154       if( cevt.CheckEvent( evt ) )
155       this->SetPixel( evt->Node.Vertex, 100, 50, 25, 255 );
156     */
157     
158     // Update visualization
159     if( this->m_Number == 0 || this->m_Percentage < 0 )
160       this->Render( );
161     this->m_Number++;
162     this->m_Number %= this->m_Percentage;
163
164     if( eevt.CheckEvent( evt ) )
165     {
166       if( this->m_RenderWindow != NULL )
167       {
168         vtkRenderer* ren =
169           this->m_RenderWindow->GetRenderers( )->GetFirstRenderer( );
170         ren->RemoveActor( this->m_PolyDataActor );
171         this->Render( );
172
173       } // fi
174
175     } // fi
176
177   } // fi
178 }
179
180 // -------------------------------------------------------------------------
181 template< class F, class R >
182 fpa::VTK::Image3DObserver< F, R >::
183 Image3DObserver( )
184   : Superclass( ),
185     m_RenderWindow( NULL ),
186     m_Number( 0 ),
187     m_Percentage( 0 )
188 {
189   this->m_PolyData = vtkSmartPointer< vtkPolyData >::New( );
190   this->m_PolyDataMapper = vtkSmartPointer< vtkPolyDataMapper >::New( );
191   this->m_PolyDataActor = vtkSmartPointer< vtkActor >::New( );
192
193   this->m_PolyData->SetPoints( vtkPoints::New( ) );
194   this->m_PolyData->SetVerts( vtkCellArray::New( ) );
195
196   this->m_PolyDataMapper->SetInputData( this->m_PolyData );
197   this->m_PolyDataActor->SetMapper( this->m_PolyDataMapper );
198   this->m_PolyDataActor->GetProperty( )->SetColor( 0, 0, 1 );
199 }
200
201 #endif // __FPA__VTK__IMAGE3DOBSERVER__HXX__
202
203 // eof - $RCSfile$