]> Creatis software - FrontAlgorithms.git/commitdiff
...
authorLeonardo Flórez-Valencia <florez-l@javeriana.edu.co>
Fri, 26 May 2017 19:25:10 +0000 (14:25 -0500)
committerLeonardo Flórez-Valencia <florez-l@javeriana.edu.co>
Fri, 26 May 2017 19:25:10 +0000 (14:25 -0500)
examples/CMakeLists.txt
examples/Dijkstra_Gaussian.cxx [new file with mode: 0644]
lib/fpa/Base/Dijkstra.hxx
lib/fpa/Base/MarksInterface.hxx
lib/fpa/Image/Functors/GaussianWeight.h

index d36f2550a5e420bead78ef2e2148c237fcf78565..30954b7e5b11e283afd3170c70012b9230cf294a 100644 (file)
@@ -5,6 +5,7 @@ if(BUILD_EXAMPLES)
     RegionGrow_Tautology
     RegionGrow_BinaryThreshold
     RegionGrow_Mori
+    Dijkstra_Gaussian
     Dijkstra_Maurer
     SkeletonFilter
     #CreateMoriInputImage
diff --git a/examples/Dijkstra_Gaussian.cxx b/examples/Dijkstra_Gaussian.cxx
new file mode 100644 (file)
index 0000000..7befb2b
--- /dev/null
@@ -0,0 +1,79 @@
+#include <itkImage.h>
+#include <itkImageFileReader.h>
+#include <itkImageFileWriter.h>
+
+#include <fpa/Image/Dijkstra.h>
+#include <fpa/Image/Functors/GaussianWeight.h>
+
+// -------------------------------------------------------------------------
+static const unsigned int VDim = 2;
+typedef double                                                 TPixel;
+typedef itk::Image< TPixel, VDim >                             TImage;
+typedef itk::ImageFileReader< TImage >                         TReader;
+typedef itk::ImageFileWriter< TImage >                         TWriter;
+typedef fpa::Image::Dijkstra< TImage, TImage >                 TFilter;
+typedef fpa::Image::Functors::GaussianWeight< TImage, TPixel > TVertexFunc;
+
+// -------------------------------------------------------------------------
+int main( int argc, char* argv[] )
+{
+  // Get arguments
+  if( argc < 5 + VDim )
+  {
+    std::cerr
+      << "Usage: " << argv[ 0 ]
+      << " input_image output_image stop_at_one_front beta";
+    for( unsigned int i = 0; i < VDim; ++i )
+      std::cerr << " s_" << i;
+    std::cerr << " ..." << std::endl;
+    return( 1 );
+
+  } // fi
+  std::string input_image_filename = argv[ 1 ];
+  std::string output_image_filename = argv[ 2 ];
+  bool stop_at_one_front = ( std::atoi( argv[ 3 ] ) == 1 );
+  double beta = std::atof( argv[ 4 ] );
+
+  TReader::Pointer reader = TReader::New( );
+  reader->SetFileName( input_image_filename );
+
+  TFilter::Pointer filter = TFilter::New( );
+  filter->SetInput( reader->GetOutput( ) );
+  filter->SetStopAtOneFront( stop_at_one_front );
+
+  for( int i = 5; i < argc; i += VDim )
+  {
+    if( i + VDim <= argc )
+    {
+      TImage::IndexType seed;
+      for( int j = 0; j < VDim; ++j )
+        seed[ j ] = std::atoi( argv[ i + j ] );
+      filter->AddSeed( seed );
+
+    } // fi
+
+  } // rof
+
+  TVertexFunc::Pointer vertex_func = TVertexFunc::New( );
+  vertex_func->SetBeta( beta );
+  vertex_func->InvertOn( );
+  filter->SetFunctor( vertex_func );
+
+  TWriter::Pointer writer = TWriter::New( );
+  writer->SetInput( filter->GetOutput( ) );
+  writer->SetFileName( output_image_filename );
+
+  try
+  {
+    writer->Update( );
+  }
+  catch( std::exception& err )
+  {
+    std::cerr << "ERROR: " << err.what( ) << std::endl;
+    return( 1 );
+
+  } // yrt
+  return( 0 );
+}
+
+// eof - $RCSfile$
index 765b9a8beaf2c984511f17fc7b28812e3bb28278..3584b1ab406a57a282c2b491345840208101522c 100644 (file)
@@ -108,7 +108,7 @@ void fpa::Base::Dijkstra< _TFilter, _TMarksInterface, _TSeedsInterface, _TMST >:
 GenerateData( )
 {
   // Init objects
-  this->_ConfigureOutputs( TOutputValue( 0 ) );
+  this->_ConfigureOutputs( std::numeric_limits< TOutputValue >::max( ) );
   this->_InitMarks( this->GetNumberOfSeeds( ) );
   TMST* mst = this->GetMinimumSpanningTree( );
 
@@ -136,8 +136,10 @@ GenerateData( )
 
     // Add neighborhood
     TVertices neighbors = this->_GetNeighbors( node.Vertex );
-    for( TVertex neigh: neighbors )
+    typename TVertices::const_iterator neighIt = neighbors.begin( );
+    for( ; neighIt != neighbors.end( ); ++neighIt )
     {
+      TVertex neigh = *neighIt;
       if( this->_IsMarked( neigh ) )
       {
         // Invoke stop at collisions
@@ -178,8 +180,8 @@ GenerateData( )
   // Complete data into minimum spanning tree
   mst->ClearSeeds( );
   mst->SetCollisions( this->m_Collisions );
-  for( TVertex seed: this->GetSeeds( ) )
-    mst->AddSeed( seed );
+  for( sIt = this->BeginSeeds( ); sIt != this->EndSeeds( ); ++sIt )
+    mst->AddSeed( *sIt );
 }
 
 #endif // __fpa__Base__Dijkstra__hxx__
index 5e7a3bc56bfc9c2065b6c7774f0c3772b366de20..8f00eb35a37d95a9cc5581dcdd780f6ef23ea919 100644 (file)
@@ -118,7 +118,7 @@ _Collisions( const TVertex& a, const TVertex& b )
           q.push( n );
 
     } // elihw
-    this->m_NumberOfFronts = this->m_NumberOfSeeds - count;
+    this->m_NumberOfFronts = this->m_NumberOfSeeds - count + 1;
 
   } // fi
   return( this->m_NumberOfFronts );
index 3b2212a75ce15ab2005a04f50956ae3c44ba2b67..aaf6209c4fc534fcb95d3aa262d21d19a09530bd 100644 (file)
@@ -38,8 +38,11 @@ namespace fpa
           fpa::Image::Functors::VertexParentBase
           );
 
+        itkBooleanMacro( Invert );
         itkGetConstMacro( Beta, double );
+        itkGetConstMacro( Invert, bool );
         itkSetMacro( Beta, double );
+        itkSetMacro( Invert, bool );
 
       public:
         virtual TOutputValue Evaluate(
@@ -49,12 +52,17 @@ namespace fpa
             double va = double( this->m_Image->GetPixel( a ) );
             double vp = double( this->m_Image->GetPixel( p ) );
             double d = va - vp;
-            return( TOutputValue( std::exp( -d * d * this->m_Beta ) ) );
+            d = ( d * d ) / this->m_Beta;
+            if( this->m_Invert )
+              return( TOutputValue( double( 1 ) - std::exp( -d ) ) );
+            else
+              return( TOutputValue( std::exp( -d ) ) );
           }
 
       protected:
         GaussianWeight( )
           : Superclass( ),
+            m_Invert( false ),
             m_Beta( double( 1 ) )
           { }
         virtual ~GaussianWeight( ) { }
@@ -64,6 +72,7 @@ namespace fpa
         Self& operator=( const Self& other );
 
       protected:
+        bool   m_Invert;
         double m_Beta;
       };