]> Creatis software - FrontAlgorithms.git/blob - tests/image/Dijkstra/Gaussian.cxx
...
[FrontAlgorithms.git] / tests / image / Dijkstra / Gaussian.cxx
1 // =========================================================================
2 // @author Leonardo Florez Valencia
3 // @email florez-l@javeriana.edu.co
4 // =========================================================================
5
6 #include <itkImage.h>
7 #include <itkImageFileReader.h>
8 #include <itkImageFileWriter.h>
9 #include <fpa/Filters/Image/Dijkstra.h>
10 #include <fpa/Functors/Dijkstra/Image/Gaussian.h>
11
12 // -------------------------------------------------------------------------
13 const unsigned int Dim = 2;
14 typedef unsigned char TInputPixel;
15 typedef float         TOutputPixel;
16 typedef itk::Image< TInputPixel, Dim >  TInputImage;
17 typedef itk::Image< TOutputPixel, Dim > TOutputImage;
18
19 // -------------------------------------------------------------------------
20 int main( int argc, char* argv[] )
21 {
22   // Get arguments
23   if( argc < 7 )
24   {
25     std::cerr
26       << "Usage: " << argv[ 0 ]
27       << " input_image output_image output_marks output_mst"
28       << " beta epsilon [seeds]"
29       << std::endl;
30     return( 1 );
31
32   } // fi
33   std::string input_image_filename = argv[ 1 ];
34   std::string output_image_filename = argv[ 2 ];
35   std::string output_marks_filename = argv[ 3 ];
36   std::string output_mst_filename = argv[ 4 ];
37   double beta = std::atof( argv[ 5 ] );
38   double eps = std::atof( argv[ 6 ] );
39
40   // Read image
41   typedef itk::ImageFileReader< TInputImage > TInputImageReader;
42   TInputImageReader::Pointer input_image_reader = TInputImageReader::New( );
43   input_image_reader->SetFileName( input_image_filename );
44
45   // Prepare weight functor
46   typedef fpa::Functors::Dijkstra::Image::Gaussian< TInputImage, TOutputPixel > TWeight;
47   TWeight::Pointer weight = TWeight::New( );
48   weight->SetBeta( beta );
49   weight->SetEpsilon( eps );
50
51   // Prepare filter
52   typedef fpa::Filters::Image::Dijkstra< TInputImage, TOutputImage > TFilter;
53   TFilter::Pointer filter = TFilter::New( );
54   filter->SetInput( input_image_reader->GetOutput( ) );
55   filter->SetWeightFunction( weight );
56
57   // Get all seeds
58   for( int i = 7; i < argc; i += Dim )
59   {
60     TInputImage::IndexType seed;
61     for( int j = 0; j < Dim; ++j )
62       if( i + j < argc )
63         seed[ j ] = std::atoi( argv[ i + j ] );
64     filter->AddSeed( seed );
65
66   } // rof
67
68   // Execute filter
69   filter->Update( );
70
71   // Save results
72   typedef itk::ImageFileWriter< TFilter::TOutputImage > TOutputWriter;
73   TOutputWriter::Pointer output_writer = TOutputWriter::New( );
74   output_writer->SetInput( filter->GetOutput( ) );
75   output_writer->SetFileName( output_image_filename );
76
77   typedef itk::ImageFileWriter< TFilter::TMarksImage > TMarksWriter;
78   TMarksWriter::Pointer marks_writer = TMarksWriter::New( );
79   marks_writer->SetInput( filter->GetMarks( ) );
80   marks_writer->SetFileName( output_marks_filename );
81
82   typedef itk::ImageFileWriter< TFilter::TMST > TMSTWriter;
83   TMSTWriter::Pointer mst_writer = TMSTWriter::New( );
84   mst_writer->SetInput( filter->GetMinimumSpanningTree( ) );
85   mst_writer->SetFileName( output_mst_filename );
86
87   try
88   {
89     output_writer->Update( );
90     marks_writer->Update( );
91     mst_writer->Update( );
92   }
93   catch( std::exception& err )
94   {
95     std::cerr << "Error caught: " << err.what( ) << std::endl;
96     return( 1 );
97
98   } // yrt
99
100   return( 0 );
101 }
102
103 // eof - $RCSfile$