]> Creatis software - FrontAlgorithms.git/blob - examples/RegionGrow_Mori.cxx
...
[FrontAlgorithms.git] / examples / RegionGrow_Mori.cxx
1 #include <chrono>
2 #include <iomanip>
3 #include <itkCommand.h>
4 #include <itkImage.h>
5 #include <itkImageFileReader.h>
6 #include <itkImageFileWriter.h>
7 #include <itkBinaryThresholdImageFilter.h>
8
9 #include <fpa/Image/MoriRegionGrow.h>
10
11 // -------------------------------------------------------------------------
12 static const unsigned int VDim = 3;
13 typedef short                                             TPixel;
14 typedef itk::Image< TPixel, VDim >                        TImage;
15 typedef itk::ImageFileReader< TImage >                    TReader;
16 typedef itk::ImageFileWriter< TImage >                    TWriter;
17 typedef fpa::Image::MoriRegionGrow< TImage, TImage >      TFilter;
18 typedef itk::BinaryThresholdImageFilter< TImage, TImage > TThreshold;
19
20 // -------------------------------------------------------------------------
21 class ShowProgressObject
22 {
23 public:
24   ShowProgressObject( itk::ProcessObject* o )
25     {
26       this->m_Process = o;
27     }
28   void ShowProgress( )
29     {
30       std::cerr
31         << "\rProgress " << std::fixed << std::setprecision( 2 )
32         << ( this->m_Process->GetProgress( ) * 100 )
33         << " %" << std::flush;
34     }
35   itk::ProcessObject::Pointer m_Process;
36 };
37
38 // -------------------------------------------------------------------------
39 int main( int argc, char* argv[] )
40 {
41   // Get arguments
42   if( argc < 6 + VDim )
43   {
44     std::cerr
45       << "Usage: " << argv[ 0 ]
46       << " input_image output_image lower upper delta";
47     for( unsigned int i = 0; i < VDim; ++i )
48       std::cerr << " s_" << i;
49     std::cerr << std::endl;
50     return( 1 );
51
52   } // fi
53   std::string input_image_filename = argv[ 1 ];
54   std::string output_image_filename = argv[ 2 ];
55   TPixel lower = std::atof( argv[ 3 ] );
56   TPixel upper = std::atof( argv[ 4 ] );
57   TPixel delta = std::atof( argv[ 5 ] );
58
59   TReader::Pointer reader = TReader::New( );
60   reader->SetFileName( input_image_filename );
61   try
62   {
63     reader->Update( );
64   }
65   catch( std::exception& err )
66   {
67     std::cerr << "ERROR: " << err.what( ) << std::endl;
68     return( 1 );
69
70   } // yrt
71
72   TFilter::Pointer filter = TFilter::New( );
73   filter->SetInput( reader->GetOutput( ) );
74   filter->SetThresholdRange( lower, upper, delta );
75   TImage::PointType pnt;
76   for( int i = 0; i < VDim; ++i )
77     pnt[ i ] = std::atof( argv[ i + 6 ] );
78
79   TImage::IndexType seed;
80   if( !( reader->GetOutput( )->TransformPhysicalPointToIndex( pnt, seed ) ) )
81   {
82     std::cerr << "ERROR: seed outside image." << std::endl;
83     return( 1 );
84
85   } // fi
86   filter->AddSeed( seed );
87
88   // to test ProgressReporter
89   /* TODO
90      ShowProgressObject progressWatch( filter );
91      typedef itk::SimpleMemberCommand< ShowProgressObject > CommandType;
92      CommandType::Pointer command = CommandType::New();
93      command->SetCallbackFunction( &progressWatch,
94      &ShowProgressObject::ShowProgress );
95      filter->AddObserver( itk::ProgressEvent( ), command );
96   */
97   std::chrono::time_point< std::chrono::system_clock > start, end;
98   start = std::chrono::system_clock::now( );
99   filter->Update( );
100   end = std::chrono::system_clock::now( );
101   std::chrono::duration< double > elapsed_seconds = end - start;
102
103   TThreshold::Pointer threshold = TThreshold::New( );
104   threshold->SetInput( filter->GetOutput( ) );
105   threshold->SetInsideValue( 255 );
106   threshold->SetOutsideValue( 0 );
107   threshold->SetLowerThreshold( 0 );
108   threshold->SetUpperThreshold( filter->GetOptimumThreshold( ) );
109   
110   TWriter::Pointer writer = TWriter::New( );
111   writer->SetInput( threshold->GetOutput( ) );
112   writer->SetFileName( output_image_filename );
113   try
114   {
115     writer->Update( );
116   }
117   catch( std::exception& err )
118   {
119     std::cerr << "ERROR: " << err.what( ) << std::endl;
120     return( 1 );
121
122   } // yrt
123
124   // Show data
125   TFilter::TCurve curve = filter->GetCurve( );
126   TFilter::TCurve::const_iterator data = curve.begin( );
127   for( ; data != curve.end( ); ++data )
128   {
129     std::cout << data->XValue << " " << data->YValue << " " << data->Diff1 << std::endl;
130   }
131   std::cout
132     << std::endl
133     << "# Opt: "
134     << curve[ filter->GetOptimumThreshold( ) ].XValue
135     << "("
136     << filter->GetOptimumThreshold( )
137      << ")"
138     << std::endl;
139   std::cout << "Time: " << elapsed_seconds.count( ) << "s" << std::endl;
140   return( 0 );
141 }
142
143 // eof - $RCSfile$