]> Creatis software - FrontAlgorithms.git/blob - appli/examples/example_ImageAlgorithmDijkstra_03.cxx
Skeletonization debugged
[FrontAlgorithms.git] / appli / examples / example_ImageAlgorithmDijkstra_03.cxx
1 #include <cmath>
2 #include <ctime>
3 #include <deque>
4 #include <iostream>
5 #include <limits>
6 #include <map>
7 #include <string>
8
9 #include <itkConstNeighborhoodIterator.h>
10 #include <itkNeighborhoodIterator.h>
11 #include <itkDanielssonDistanceMapImageFilter.h>
12 #include <itkImage.h>
13 #include <itkImageFileReader.h>
14 #include <itkImageFileWriter.h>
15 #include <itkImageToVTKImageFilter.h>
16
17 #include <vtkPoints.h>
18 #include <vtkCellArray.h>
19 #include <vtkFloatArray.h>
20 #include <vtkPolyData.h>
21 #include <vtkSmartPointer.h>
22
23 #include <fpa/Image/DijkstraWithSphereBacktracking.h>
24 #include <fpa/VTK/ImageMPR.h>
25 #include <fpa/VTK/Image3DObserver.h>
26
27 // -------------------------------------------------------------------------
28 const unsigned int Dim = 3;
29 typedef double                               TPixel;
30 typedef double                               TScalar;
31 typedef itk::Image< TPixel, Dim >            TImage;
32 typedef itk::ImageToVTKImageFilter< TImage > TVTKImage;
33
34 typedef itk::ImageFileReader< TImage >          TImageReader;
35 typedef itk::ImageFileWriter< TImage >          TImageWriter;
36 typedef fpa::Image::DijkstraWithSphereBacktracking< TImage, TScalar > TDijkstra;
37
38 typedef fpa::VTK::ImageMPR TMPR;
39 typedef fpa::VTK::Image3DObserver< TDijkstra, vtkRenderWindow > TDijkstraObs;
40
41 // -------------------------------------------------------------------------
42 int main( int argc, char* argv[] )
43 {
44   if( argc < 6 )
45   {
46     std::cerr
47       << "Usage: " << argv[ 0 ]
48       << " input_image x y z output_image"
49       << " visual_debug"
50       << std::endl;
51     return( 1 );
52
53   } // fi
54   std::string input_image_fn = argv[ 1 ];
55   TImage::IndexType seed_idx;
56   seed_idx[ 0 ] = std::atoi( argv[ 2 ] );
57   seed_idx[ 1 ] = std::atoi( argv[ 3 ] );
58   seed_idx[ 2 ] = std::atoi( argv[ 4 ] );
59   std::string output_image_fn = argv[ 5 ];
60   bool visual_debug = false;
61   if( argc > 6 )
62     visual_debug = ( std::atoi( argv[ 6 ] ) == 1 );
63
64   // Read image
65   TImageReader::Pointer input_image_reader = TImageReader::New( );
66   input_image_reader->SetFileName( input_image_fn );
67   try
68   {
69     input_image_reader->Update( );
70   }
71   catch( itk::ExceptionObject& err )
72   {
73     std::cerr << "Error caught: " << err << std::endl;
74     return( 1 );
75
76   } // yrt
77   TImage::ConstPointer input_image = input_image_reader->GetOutput( );
78
79   // Show image
80   TVTKImage::Pointer vtk_image = TVTKImage::New( );
81   vtk_image->SetInput( input_image );
82   vtk_image->Update( );
83
84   TMPR view;
85   view.SetBackground( 0.3, 0.2, 0.8 );
86   view.SetSize( 800, 800 );
87   view.SetImage( vtk_image->GetOutput( ) );
88
89   // Allow some interaction
90   view.Start( );
91
92   // Extract paths
93   TDijkstra::Pointer paths = TDijkstra::New( );
94   paths->AddSeed( seed_idx, TScalar( 0 ) );
95   paths->SetInput( input_image );
96   paths->SetNeighborhoodOrder( 1 );
97
98   if( visual_debug )
99   {
100     // Configure observer
101     TDijkstraObs::Pointer obs = TDijkstraObs::New( );
102     obs->SetRenderWindow( view.GetWindow( ) );
103     paths->AddObserver( itk::AnyEvent( ), obs );
104     paths->ThrowEventsOn( );
105   }
106   else
107     paths->ThrowEventsOff( );
108   std::clock_t start = std::clock( );
109   paths->Update( );
110   std::clock_t end = std::clock( );
111   double seconds = double( end - start ) / double( CLOCKS_PER_SEC );
112   std::cout << "Paths extraction time = " << seconds << std::endl;
113
114   // Create polydata
115   vtkSmartPointer< vtkPoints > points =
116     vtkSmartPointer< vtkPoints >::New( );
117   vtkSmartPointer< vtkCellArray > cells =
118     vtkSmartPointer< vtkCellArray >::New( );
119   vtkSmartPointer< vtkFloatArray > scalars =
120     vtkSmartPointer< vtkFloatArray >::New( );
121
122   const TDijkstra::TVertices& endpoints = paths->GetEndPoints( );
123   const TDijkstra::TTree& tree = paths->GetFinalTree( );
124   TDijkstra::TVertices::const_iterator epIt = endpoints.begin( );
125   for( unsigned int epId = 0; epIt != endpoints.end( ); ++epIt, ++epId )
126   {
127     double pd = double( epId ) / double( endpoints.size( ) - 1 );
128
129     TDijkstra::TVertex idx = *epIt;
130     do
131     {
132       TImage::PointType pnt;
133       input_image->TransformIndexToPhysicalPoint( idx, pnt );
134
135       points->InsertNextPoint( pnt[ 0 ], pnt[ 1 ], pnt[ 2 ] );
136       scalars->InsertNextTuple1( pd );
137       if( idx != *epIt )
138       {
139         cells->InsertNextCell( 2 );
140         cells->InsertCellPoint( points->GetNumberOfPoints( ) - 2 );
141         cells->InsertCellPoint( points->GetNumberOfPoints( ) - 1 );
142
143       } // fi
144       idx = tree.find( idx )->second;
145
146     } while( idx != tree.find( idx )->second );
147
148   } // rof
149
150   vtkSmartPointer< vtkPolyData > vtk_tree =
151     vtkSmartPointer< vtkPolyData >::New( );
152   vtk_tree->SetPoints( points );
153   vtk_tree->SetLines( cells );
154   vtk_tree->GetPointData( )->SetScalars( scalars );
155
156   view.AddPolyData( vtk_tree );
157   view.Render( );
158   view.Start( );
159
160   /* TODO
161      TImageWriter::Pointer dijkstra_writer =
162      TImageWriter::New( );
163      dijkstra_writer->SetInput( paths->GetOutput( ) );
164      dijkstra_writer->SetFileName( "dijkstra.mhd" );
165      dijkstra_writer->Update( );
166   */
167
168   return( 0 );
169 }
170
171 // eof - $RCSfile$