]> Creatis software - FrontAlgorithms.git/blob - appli/examples/example_ImageAlgorithm_Skeletonization.cxx
Second order neighbors debugged in 3D
[FrontAlgorithms.git] / appli / examples / example_ImageAlgorithm_Skeletonization.cxx
1 #include <ctime>
2 #include <iostream>
3 #include <string>
4
5 #include <itkConstNeighborhoodIterator.h>
6 #include <itkDanielssonDistanceMapImageFilter.h>
7 #include <itkImage.h>
8 #include <itkImageFileReader.h>
9 #include <itkImageFileWriter.h>
10 #include <itkImageToVTKImageFilter.h>
11
12 #include <vtkImageMarchingCubes.h>
13
14 #include <fpa/Image/RegionGrowWithMultipleThresholds.h>
15 #include <fpa/Image/DijkstraWithSphereBacktracking.h>
16 #include <fpa/VTK/ImageMPR.h>
17 #include <fpa/VTK/Image3DObserver.h>
18
19 // -------------------------------------------------------------------------
20 const unsigned int Dim = 3;
21 typedef short                                TPixel;
22 typedef float                                TScalar;
23 typedef itk::Image< TPixel, Dim >            TImage;
24 typedef itk::Image< TScalar, Dim >           TDistanceMap;
25 typedef itk::ImageToVTKImageFilter< TImage > TVTKImage;
26 typedef itk::ImageToVTKImageFilter< TDistanceMap > TVTKDistanceMap;
27
28 typedef itk::ImageFileReader< TImage >                       TImageReader;
29 typedef itk::ImageFileWriter< TImage >                       TImageWriter;
30 typedef itk::ImageFileWriter< TDistanceMap >           TDistanceMapWriter;
31 typedef fpa::Image::RegionGrowWithMultipleThresholds< TImage > TSegmentor;
32 typedef
33 itk::DanielssonDistanceMapImageFilter< TImage, TDistanceMap >
34 TDistanceFilter;
35 typedef
36 fpa::Image::DijkstraWithSphereBacktracking< TDistanceMap, TScalar >
37 TDijkstra;
38
39 typedef fpa::VTK::ImageMPR TMPR;
40 typedef fpa::VTK::Image3DObserver< TSegmentor, vtkRenderWindow > TSegmentorObs;
41 typedef fpa::VTK::Image3DObserver< TDijkstra, vtkRenderWindow > TDijkstraObs;
42
43 // -------------------------------------------------------------------------
44 int main( int argc, char* argv[] )
45 {
46   if( argc < 8 )
47   {
48     std::cerr
49       << "Usage: " << argv[ 0 ]
50       << " input_image thr_0 thr_1 step"
51       << " output_segmentation output_distancemap output_dijkstra"
52       << " visual_debug"
53       << std::endl;
54     return( 1 );
55
56   } // fi
57   std::string input_image_fn = argv[ 1 ];
58   TPixel thr_0 = TPixel( std::atof( argv[ 2 ] ) );
59   TPixel thr_1 = TPixel( std::atof( argv[ 3 ] ) );
60   unsigned int step = std::atoi( argv[ 4 ] );
61   std::string output_segmentation_fn = argv[ 5 ];
62   std::string output_distancemap_fn = argv[ 6 ];
63   std::string output_dijkstra_fn = argv[ 7 ];
64   bool visual_debug = false;
65   if( argc > 8 )
66     visual_debug = ( std::atoi( argv[ 8 ] ) == 1 );
67
68   // Read image
69   TImageReader::Pointer input_image_reader = TImageReader::New( );
70   input_image_reader->SetFileName( input_image_fn );
71   try
72   {
73     input_image_reader->Update( );
74   }
75   catch( itk::ExceptionObject& err )
76   {
77     std::cerr << "Error caught: " << err << std::endl;
78     return( 1 );
79
80   } // yrt
81   TImage::ConstPointer input_image = input_image_reader->GetOutput( );
82
83   // Show image
84   TVTKImage::Pointer vtk_image = TVTKImage::New( );
85   vtk_image->SetInput( input_image );
86   vtk_image->Update( );
87
88   TMPR view;
89   view.SetBackground( 0.3, 0.2, 0.8 );
90   view.SetSize( 800, 800 );
91   view.SetImage( vtk_image->GetOutput( ) );
92
93   // Wait for a seed to be given
94   while( view.GetNumberOfSeeds( ) == 0 )
95     view.Start( );
96
97   // Get seed from user
98   double seed[ 3 ];
99   view.GetSeed( 0, seed );
100   TImage::PointType seed_pnt;
101   seed_pnt[ 0 ] = seed[ 0 ];
102   seed_pnt[ 1 ] = seed[ 1 ];
103   seed_pnt[ 2 ] = seed[ 2 ];
104   TImage::IndexType seed_idx;
105   input_image->TransformPhysicalPointToIndex( seed_pnt, seed_idx );
106
107   // Segment input image
108   TSegmentor::Pointer segmentor = TSegmentor::New( );
109   segmentor->AddThresholds( thr_0, thr_1, step );
110   segmentor->AddSeed( seed_idx, 0 );
111   segmentor->SetInput( input_image );
112   segmentor->SetNeighborhoodOrder( 1 );
113   segmentor->SetDifferenceThreshold( double( 3 ) );
114   segmentor->SetInsideValue( TPixel( 0 ) );  // WARNING: This is to optimize
115   segmentor->SetOutsideValue( TPixel( 1 ) ); // distance map calculation
116   if( visual_debug )
117   {
118     // Configure observer
119     TSegmentorObs::Pointer obs = TSegmentorObs::New( );
120     obs->SetRenderWindow( view.GetWindow( ) );
121     segmentor->AddObserver( itk::AnyEvent( ), obs );
122     segmentor->ThrowEventsOn( );
123   }
124   else
125     segmentor->ThrowEventsOff( );
126   std::clock_t start = std::clock( );
127   segmentor->Update( );
128   std::clock_t end = std::clock( );
129   double seconds = double( end - start ) / double( CLOCKS_PER_SEC );
130   std::cout << "Segmentation time = " << seconds << std::endl;
131
132   // Extract distance map
133   TDistanceFilter::Pointer distanceMap = TDistanceFilter::New( );
134   distanceMap->SetInput( segmentor->GetOutput( ) );
135   distanceMap->InputIsBinaryOn( );
136   distanceMap->SquaredDistanceOn( );
137   distanceMap->UseImageSpacingOn( );
138   start = std::clock( );
139   distanceMap->Update( );
140   end = std::clock( );
141   seconds = double( end - start ) / double( CLOCKS_PER_SEC );
142   std::cout << "Distance map time = " << seconds << std::endl;
143
144   TVTKDistanceMap::Pointer vtk_distanceMap = TVTKDistanceMap::New( );
145   vtk_distanceMap->SetInput( distanceMap->GetOutput( ) );
146   vtk_distanceMap->Update( );
147
148   vtkSmartPointer< vtkImageMarchingCubes > mc =
149     vtkSmartPointer< vtkImageMarchingCubes >::New( );
150   mc->SetInputData( vtk_distanceMap->GetOutput( ) );
151   mc->SetValue( 0, 1e-2 );
152   mc->Update( );
153   view.AddPolyData( mc->GetOutput( ), 1, 1, 1, 0.4 );
154   view.Render( );
155
156   // Extract paths
157   TDijkstra::Pointer paths = TDijkstra::New( );
158   paths->AddSeed( segmentor->GetSeed( 0 ), TScalar( 0 ) );
159   paths->SetInput( distanceMap->GetOutput( ) );
160   paths->SetNeighborhoodOrder( 2 );
161
162   if( visual_debug )
163   {
164     // Configure observer
165     TDijkstraObs::Pointer obs = TDijkstraObs::New( );
166     obs->SetRenderWindow( view.GetWindow( ) );
167     paths->AddObserver( itk::AnyEvent( ), obs );
168     paths->ThrowEventsOn( );
169   }
170   else
171     paths->ThrowEventsOff( );
172   start = std::clock( );
173   paths->Update( );
174   end = std::clock( );
175   seconds = double( end - start ) / double( CLOCKS_PER_SEC );
176   std::cout << "Paths extraction time = " << seconds << std::endl;
177
178   std::cout << "Final seed: " << paths->GetSeed( 0 ) << std::endl;
179
180   // Create polydata
181   vtkSmartPointer< vtkPoints > points =
182     vtkSmartPointer< vtkPoints >::New( );
183   vtkSmartPointer< vtkCellArray > cells =
184     vtkSmartPointer< vtkCellArray >::New( );
185   vtkSmartPointer< vtkFloatArray > scalars =
186     vtkSmartPointer< vtkFloatArray >::New( );
187
188   const TDijkstra::TVertices& endpoints = paths->GetEndPoints( );
189   const TDijkstra::TTree& tree = paths->GetFinalTree( );
190   TDijkstra::TVertices::const_iterator epIt = endpoints.begin( );
191   for( unsigned int epId = 0; epIt != endpoints.end( ); ++epIt, ++epId )
192   {
193     double pd = double( epId ) / double( endpoints.size( ) - 1 );
194
195     TDijkstra::TVertex idx = *epIt;
196     do
197     {
198       TImage::PointType pnt;
199       input_image->TransformIndexToPhysicalPoint( idx, pnt );
200
201       points->InsertNextPoint( pnt[ 0 ], pnt[ 1 ], pnt[ 2 ] );
202       scalars->InsertNextTuple1( pd );
203       if( idx != *epIt )
204       {
205         cells->InsertNextCell( 2 );
206         cells->InsertCellPoint( points->GetNumberOfPoints( ) - 2 );
207         cells->InsertCellPoint( points->GetNumberOfPoints( ) - 1 );
208
209       } // fi
210       idx = tree.find( idx )->second;
211
212     } while( idx != tree.find( idx )->second );
213
214   } // rof
215
216   vtkSmartPointer< vtkPolyData > vtk_tree =
217     vtkSmartPointer< vtkPolyData >::New( );
218   vtk_tree->SetPoints( points );
219   vtk_tree->SetLines( cells );
220   vtk_tree->GetPointData( )->SetScalars( scalars );
221
222   view.AddPolyData( vtk_tree );
223   view.Render( );
224   view.Start( );
225
226   itk::ImageFileWriter< TImage >::Pointer segmentation_writer =
227     itk::ImageFileWriter< TImage >::New( );
228   segmentation_writer->SetInput( segmentor->GetOutput( ) );
229   segmentation_writer->SetFileName( output_segmentation_fn );
230
231   itk::ImageFileWriter< TDistanceMap >::Pointer dmap_writer =
232     itk::ImageFileWriter< TDistanceMap >::New( );
233   dmap_writer->SetInput( distanceMap->GetOutput( ) );
234   dmap_writer->SetFileName( output_distancemap_fn );
235
236   itk::ImageFileWriter< TDistanceMap >::Pointer dijk_writer =
237     itk::ImageFileWriter< TDistanceMap >::New( );
238   dijk_writer->SetInput( paths->GetOutput( ) );
239   dijk_writer->SetFileName( output_dijkstra_fn );
240
241   try
242   {
243     segmentation_writer->Update( );
244     dmap_writer->Update( );
245     dijk_writer->Update( );
246   }
247   catch( itk::ExceptionObject& err )
248   {
249     std::cerr << "Error caught: " << err << std::endl;
250     return( -1 );
251
252   } // yrt
253
254   return( 0 );
255 }
256
257 // eof - $RCSfile$