]> Creatis software - FrontAlgorithms.git/commitdiff
...
authorLeonardo Flórez-Valencia <florez-l@javeriana.edu.co>
Tue, 7 Feb 2017 13:57:23 +0000 (08:57 -0500)
committerLeonardo Flórez-Valencia <florez-l@javeriana.edu.co>
Tue, 7 Feb 2017 13:57:23 +0000 (08:57 -0500)
examples/MoriRegionGrow_00.cxx
lib/fpa/Image/MoriRegionGrow.h
lib/fpa/Image/MoriRegionGrow.hxx

index d0dabff9eca404d19327bce07a350a6007b71153..1d7049c36745e9b66de3b0396e6001e89ced8f03 100644 (file)
@@ -1,4 +1,3 @@
-#include <itkCommand.h>
 #include <itkImage.h>
 #include <itkImageFileReader.h>
 #include <itkImageFileWriter.h>
@@ -10,71 +9,20 @@ typedef itk::ImageFileReader< TImage >               TReader;
 typedef itk::ImageFileWriter< TImage >               TWriter;
 typedef fpa::Image::MoriRegionGrow< TImage, TImage > TFilter;
 
-// -------------------------------------------------------------------------
-/**
- */
-class MyObserver
-  : public itk::Command
-{
-public:
-  typedef TFilter::TStartEvent     TStartEvent;
-  typedef TFilter::TEndEvent       TEndEvent;
-  typedef TFilter::TStartLoopEvent TStartLoopEvent;
-  typedef TFilter::TEndLoopEvent   TEndLoopEvent;
-  typedef TFilter::TPushEvent      TPushEvent;
-  typedef TFilter::TPopEvent       TPopEvent;
-  typedef TFilter::TMarkEvent      TMarkEvent;
-
-public:
-  itkNewMacro( MyObserver );
-
-public:
-  virtual void Execute(
-    itk::Object* caller, const itk::EventObject& event
-    ) override
-    {
-      this->Execute( const_cast< const itk::Object* >( caller ), event );
-    }
-  virtual void Execute(
-    const itk::Object* object, const itk::EventObject& event
-    ) override
-    {
-      /* TODO
-         if( TStartEvent( ).CheckEvent( &event ) )
-         std::cout << "Start" << std::endl;
-         else if( TEndEvent( ).CheckEvent( &event ) )
-         std::cout << "End" << std::endl;
-         else if( TStartLoopEvent( ).CheckEvent( &event ) )
-         std::cout << "StartLoop" << std::endl;
-         else if( TEndLoopEvent( ).CheckEvent( &event ) )
-         std::cout << "EndLoop" << std::endl;
-         else if( TMarkEvent( ).CheckEvent( &event ) )
-         {
-         const TMarkEvent* mark = dynamic_cast< const TMarkEvent* >( &event );
-         std::cout << "Mark: " << mark->Vertex << std::endl;
-
-         } // fi
-      */
-      /* TODO
-         TPushEvent;
-         TPopEvent;
-      */
-    }
-};
-
 // -------------------------------------------------------------------------
 int main( int argc, char* argv[] )
 {
-  if( argc < 3 )
+  if( argc < 4 )
   {
     std::cerr
       << "Usage: " << argv[ 0 ]
-      << " input_filename output_filename" << std::endl;
+      << " input_filename output_filename sensitivity" << std::endl;
     return( 1 );
 
   } // fi
   std::string in_fname = argv[ 1 ];
   std::string out_fname = argv[ 2 ];
+  double sensitivity = std::atof( argv[ 3 ] );
   int seed_x = 111;
   int seed_y = 91;
 
@@ -92,13 +40,9 @@ int main( int argc, char* argv[] )
   filter->SetStep( 1 );
   filter->SetInsideValue( 255 );
   filter->SetOutsideValue( 0 );
+  filter->SetSensitivity( sensitivity );
   filter->AddSeed( seed, filter->GetInsideValue( ) );
 
-  /* TODO
-     MyObserver::Pointer obs = MyObserver::New( );
-     filter->AddObserver( itk::AnyEvent( ), obs );
-  */
-
   TWriter::Pointer writer = TWriter::New( );
   writer->SetInput( filter->GetOutput( ) );
   writer->SetFileName( out_fname );
index 50e8bb1e3cc1769b260bfbfa4fbe4a24c00ac5bc..7a509ec9ac9339fb27f1f44d52aef2cc247a4b52 100644 (file)
@@ -39,10 +39,12 @@ namespace fpa
       itkGetConstMacro( Lower, TPixel );
       itkGetConstMacro( Upper, TPixel );
       itkGetConstMacro( Step, TPixel );
+      itkGetConstMacro( Sensitivity, double );
 
       itkSetMacro( Lower, TPixel );
       itkSetMacro( Upper, TPixel );
       itkSetMacro( Step, TPixel );
+      itkSetMacro( Sensitivity, double );
 
     protected:
       MoriRegionGrow( );
@@ -65,9 +67,17 @@ namespace fpa
       TPixel m_Lower;
       TPixel m_Upper;
       TPixel m_Step;
+      double m_Sensitivity;
 
       _TQueue m_NextQueue;
       unsigned long m_ActualCount;
+      unsigned long m_PrevCount;
+
+      // Standard deviation
+      double m_N;
+      double m_S1;
+      double m_S2;
+      double m_STD;
     };
 
   } // ecapseman
index 3b06e41cdd463dce08115f3e6cdf7eb8011442a8..c414f25d92d8580d1d0008fb298cacdc4ed242f6 100644 (file)
@@ -6,7 +6,8 @@ template< class _TInputImage, class _TOutputImage >
 fpa::Image::MoriRegionGrow< _TInputImage, _TOutputImage >::
 MoriRegionGrow( )
   : Superclass( ),
-    m_Step( TPixel( 1 ) )
+    m_Step( TPixel( 1 ) ),
+    m_Sensitivity( double( 1 ) )
 {
   this->m_Upper = std::numeric_limits< TPixel >::max( );
   if( std::numeric_limits< TPixel >::is_integer )
@@ -34,7 +35,32 @@ _ContinueGenerateData( )
     dynamic_cast< TBinThresholdFunction* >( this->GetGrowFunction( ) );
   TPixel u = functor->GetUpper( );
 
-  std::cout << long( u ) << " " << this->m_ActualCount << std::endl;
+  // Analyze flooding
+  bool stop = false;
+  if( this->m_N > double( 1 ) )
+  {
+    double Q = double( this->m_PrevCount );
+    Q += ( this->m_STD * this->m_Sensitivity );
+    if( double( this->m_ActualCount ) > Q )
+      stop = true;
+
+  } // fi
+  double v;
+  if( this->m_PrevCount > 0 )
+    v = double( this->m_ActualCount ) - double( this->m_PrevCount );
+  else
+    v = double( 0 );
+  this->m_PrevCount = this->m_ActualCount;
+  this->m_N += double( 1 );
+  this->m_S1 += v;
+  this->m_S2 += v * v;
+  if( this->m_N > double( 1 ) )
+    this->m_STD =
+      ( this->m_S2 - ( ( this->m_S1 * this->m_S1 ) / this->m_N ) ) /
+      ( this->m_N - double( 1 ) );
+  else
+    this->m_STD = double( 0 );
+  this->m_STD = std::sqrt( this->m_STD );
 
   if( u < this->m_Upper )
   {
@@ -45,7 +71,7 @@ _ContinueGenerateData( )
     this->m_Queue = this->m_NextQueue;
     while( this->m_NextQueue.size( ) > 0 )
       this->m_NextQueue.pop( );
-    return( true );
+    return( !stop );
   }
   else
     return( false );
@@ -61,6 +87,11 @@ _BeforeGenerateData( )
   while( this->m_NextQueue.size( ) > 0 )
     this->m_NextQueue.pop( );
   this->m_ActualCount = 0;
+  this->m_PrevCount = 0;
+  this->m_N   = double( 0 );
+  this->m_S1  = double( 0 );
+  this->m_S2  = double( 0 );
+  this->m_STD = double( 0 );
   TBinThresholdFunction* functor =
     dynamic_cast< TBinThresholdFunction* >( this->GetGrowFunction( ) );
   functor->SetLower( this->m_Lower );