]> Creatis software - cpPlugins.git/commitdiff
Discontinuity map (Wan Yu) filter added.
authorLeonardo Flórez-Valencia <florez-l@javeriana.edu.co>
Tue, 7 Mar 2017 16:43:41 +0000 (11:43 -0500)
committerLeonardo Flórez-Valencia <florez-l@javeriana.edu.co>
Tue, 7 Mar 2017 16:43:41 +0000 (11:43 -0500)
lib/cpExtensions/Algorithms/DiscontinuityMapImageFilter.h [new file with mode: 0644]
lib/cpExtensions/Algorithms/DiscontinuityMapImageFilter.hxx [new file with mode: 0644]
lib/cpInstances/Images/ITKBoxImageFilters.i [new file with mode: 0644]
lib/cpPlugins/Pipeline/Functor.h
plugins/cpExtensions/DiscontinuityMapImageFilter.cxx [new file with mode: 0644]
plugins/cpExtensions/DiscontinuityMapImageFilter.h [new file with mode: 0644]
plugins/cpExtensions/ImageFunctorFilter.cxx
plugins/cpExtensions/cpExtensions.i

diff --git a/lib/cpExtensions/Algorithms/DiscontinuityMapImageFilter.h b/lib/cpExtensions/Algorithms/DiscontinuityMapImageFilter.h
new file mode 100644 (file)
index 0000000..43f44c3
--- /dev/null
@@ -0,0 +1,66 @@
+// -------------------------------------------------------------------------
+// @author Leonardo Florez-Valencia (florez-l@javeriana.edu.co)
+// -------------------------------------------------------------------------
+
+#ifndef __cpExtensions__Algorithms__DiscontinuityMapImageFilter__h__
+#define __cpExtensions__Algorithms__DiscontinuityMapImageFilter__h__
+
+#include <itkImage.h>
+#include <itkBoxImageFilter.h>
+
+namespace cpExtensions
+{
+  namespace Algorithms
+  {
+    /**
+     */
+    template< class _TInputImage >
+    class DiscontinuityMapImageFilter
+      : public itk::BoxImageFilter< _TInputImage, itk::Image< double, _TInputImage::ImageDimension > >
+    {
+    public:
+      typedef _TInputImage                                        TInputImage;
+      typedef itk::Image< double, _TInputImage::ImageDimension > TOutputImage;
+      typedef DiscontinuityMapImageFilter                      Self;
+      typedef itk::BoxImageFilter< TInputImage, TOutputImage > Superclass;
+      typedef itk::SmartPointer< Self >                        Pointer;
+      typedef itk::SmartPointer< const Self >                  ConstPointer;
+
+      typedef typename TInputImage::PixelType  TInputPixel;
+      typedef typename TInputImage::RegionType TRegion;
+      typedef typename TOutputImage::PixelType TOutputPixel;
+
+    public:
+      itkNewMacro( Self );
+      itkTypeMacro( DiscontinuityMapImageFilter, itkBoxImageFilter );
+
+      itkBooleanMacro( UseSquareRoot );
+      itkGetConstMacro( UseSquareRoot, bool );
+      itkSetMacro( UseSquareRoot, bool );
+
+    protected:
+      DiscontinuityMapImageFilter( );
+      virtual ~DiscontinuityMapImageFilter( );
+
+      virtual void ThreadedGenerateData( const TRegion& region, itk::ThreadIdType threadId ) override;
+
+    private:
+      // Purposely not implemented.
+      DiscontinuityMapImageFilter( const Self& );
+      void operator=( const Self& );
+
+    protected:
+      bool m_UseSquareRoot;
+    };
+
+  } // ecapseman
+
+} // ecapseman
+
+#ifndef ITK_MANUAL_INSTANTIATION
+#  include <cpExtensions/Algorithms/DiscontinuityMapImageFilter.hxx>
+#endif // ITK_MANUAL_INSTANTIATION
+
+#endif // __cpExtensions__Algorithms__DiscontinuityMapImageFilter__h__
+
+// eof - $RCSfile$
diff --git a/lib/cpExtensions/Algorithms/DiscontinuityMapImageFilter.hxx b/lib/cpExtensions/Algorithms/DiscontinuityMapImageFilter.hxx
new file mode 100644 (file)
index 0000000..877edbc
--- /dev/null
@@ -0,0 +1,92 @@
+// -------------------------------------------------------------------------
+// @author Leonardo Florez-Valencia (florez-l@javeriana.edu.co)
+// -------------------------------------------------------------------------
+
+#ifndef __cpExtensions__Algorithms__DiscontinuityMapImageFilter__hxx__
+#define __cpExtensions__Algorithms__DiscontinuityMapImageFilter__hxx__
+
+#include <cmath>
+#include <itkConstNeighborhoodIterator.h>
+#include <itkImageRegionIterator.h>
+#include <itkNeighborhoodAlgorithm.h>
+#include <itkZeroFluxNeumannBoundaryCondition.h>
+
+// -------------------------------------------------------------------------
+template< class _TInputImage >
+cpExtensions::Algorithms::DiscontinuityMapImageFilter< _TInputImage >::
+DiscontinuityMapImageFilter( )
+  : Superclass( ),
+    m_UseSquareRoot( false )
+{
+}
+
+// -------------------------------------------------------------------------
+template< class _TInputImage >
+cpExtensions::Algorithms::DiscontinuityMapImageFilter< _TInputImage >::
+~DiscontinuityMapImageFilter( )
+{
+}
+
+// -------------------------------------------------------------------------
+template< class _TInputImage >
+void cpExtensions::Algorithms::DiscontinuityMapImageFilter< _TInputImage >::
+ThreadedGenerateData( const TRegion& region, itk::ThreadIdType threadId )
+{
+  unsigned int i;
+
+  itk::ZeroFluxNeumannBoundaryCondition< TInputImage > nbc;
+  itk::ConstNeighborhoodIterator< TInputImage > bit;
+  itk::ImageRegionIterator< TOutputImage > it;
+
+  // Allocate output
+  typename TOutputImage::Pointer output = this->GetOutput( );
+  typename TInputImage::ConstPointer input = this->GetInput( );
+
+  // Find the data-set boundary "faces"
+  typename itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator< TInputImage >::FaceListType faceList;
+  itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator< TInputImage > bC;
+  faceList = bC( input, region, this->GetRadius( ) );
+
+  typename itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator< TInputImage >::FaceListType::iterator fit;
+
+  TOutputPixel sum1, sum2;
+  TOutputPixel neighborhoodSize;
+
+  // Process each of the boundary faces.  These are N-d regions which border
+  // the edge of the buffer.
+  for( fit = faceList.begin( ); fit != faceList.end( ); ++fit )
+  {
+    bit = itk::ConstNeighborhoodIterator< TInputImage >( this->GetRadius( ), input, *fit );
+    neighborhoodSize = bit.Size( );
+    it = itk::ImageRegionIterator< TOutputImage >( output, *fit );
+    bit.OverrideBoundaryCondition( &nbc );
+    bit.GoToBegin( );
+    while( !bit.IsAtEnd( ) )
+    {
+      sum1 = itk::NumericTraits< TOutputPixel >::ZeroValue( );
+      sum2 = itk::NumericTraits< TOutputPixel >::ZeroValue( );
+      for ( i = 0; i < neighborhoodSize; ++i )
+      {
+        TOutputPixel v = TOutputPixel( bit.GetPixel( i ) );
+        sum1 += v;
+        sum2 += v * v;
+      }
+      TOutputPixel mean = ( sum1 / TOutputPixel( neighborhoodSize ) );
+      TOutputPixel var =
+        ( sum2 - ( ( sum1 * sum1 ) / TOutputPixel( neighborhoodSize ) ) ) /
+        ( TOutputPixel( neighborhoodSize ) - TOutputPixel( 1 ) );
+      if( this->m_UseSquareRoot )
+        it.Set( std::sqrt( double( var ) ) );
+      else
+        it.Set( var );
+      ++bit;
+      ++it;
+
+    } // elihw
+
+  } // rof
+}
+
+#endif // __cpExtensions__Algorithms__DiscontinuityMapImageFilter__hxx__
+
+// eof - $RCSfile$
diff --git a/lib/cpInstances/Images/ITKBoxImageFilters.i b/lib/cpInstances/Images/ITKBoxImageFilters.i
new file mode 100644 (file)
index 0000000..69f40fd
--- /dev/null
@@ -0,0 +1,10 @@
+header #define ITK_MANUAL_INSTANTIATION
+
+define i_scalars=#scalar_types#
+define o_scalars=#scalar_types#
+
+tinclude itkBoxImageFilter:h|hxx
+
+instances itk::BoxImageFilter< itk::Image< #i_scalars# , #pdims# >, itk::Image< #o_scalars#, #pdims# > >
+
+** eof - $RCSfile$
\ No newline at end of file
index bac4523e85c425143527caeb25076672cf6952ec..96255b7f2b18c2e78777cc056931b0b34059484e 100644 (file)
@@ -23,6 +23,25 @@ namespace cpPlugins
       itkTypeMacro( Functor, ProcessObject );
       cpPlugins_Id_Macro( Functor, Object );
 
+    public:
+      virtual void Instantiate( itk::LightObject* filter ) = 0;
+
+      template< class _TFunctor >
+      _TFunctor* GetFunctor( )
+        {
+          return(
+            dynamic_cast< _TFunctor* >( this->m_Functor.GetPointer( ) )
+            );
+        }
+
+      template< class _TFunctor >
+      const _TFunctor* GetFunctor( ) const
+        {
+          return(
+            dynamic_cast< const _TFunctor* >( this->m_Functor.GetPointer( ) )
+            );
+        }
+
     protected:
       Functor( );
       virtual ~Functor( );
@@ -31,6 +50,9 @@ namespace cpPlugins
       // Purposely not implemented
       Functor( const Self& );
       Self& operator=( const Self& );
+
+    protected:
+      itk::LightObject::Pointer m_Functor;
     };
 
   } // ecapseman
diff --git a/plugins/cpExtensions/DiscontinuityMapImageFilter.cxx b/plugins/cpExtensions/DiscontinuityMapImageFilter.cxx
new file mode 100644 (file)
index 0000000..6d64365
--- /dev/null
@@ -0,0 +1,49 @@
+#include "DiscontinuityMapImageFilter.h"
+#include <cpInstances/DataObjects/Image.h>
+
+#include <cpExtensions/Algorithms/DiscontinuityMapImageFilter.h>
+
+// -------------------------------------------------------------------------
+cpPluginscpExtensions::DiscontinuityMapImageFilter::
+DiscontinuityMapImageFilter( )
+  : Superclass( )
+{
+  typedef cpInstances::DataObjects::Image _TImage;
+
+  this->_ConfigureInput< _TImage >( "Input", true, false );
+  this->_ConfigureOutput< _TImage >( "Output" );
+  this->m_Parameters.ConfigureAsInt( "Radius", 1 );
+  this->m_Parameters.ConfigureAsBool( "UseSquareRoot", false );
+}
+
+// -------------------------------------------------------------------------
+cpPluginscpExtensions::DiscontinuityMapImageFilter::
+~DiscontinuityMapImageFilter( )
+{
+}
+
+// -------------------------------------------------------------------------
+void cpPluginscpExtensions::DiscontinuityMapImageFilter::
+_GenerateData( )
+{
+  auto o = this->GetInputData( "Input" );
+  cpPlugins_Demangle_Image_ScalarPixels_AllDims_1( o, _GD0 )
+    this->_Error( "Invalid input image pixel type." );
+}
+
+// -------------------------------------------------------------------------
+template< class _TInputImage >
+void cpPluginscpExtensions::DiscontinuityMapImageFilter::
+_GD0( _TInputImage* input )
+{
+  typedef cpExtensions::Algorithms::DiscontinuityMapImageFilter< _TInputImage > _TFilter;
+
+  auto filter = this->_CreateITK< _TFilter >( );
+  filter->SetInput( input );
+  filter->SetRadius( this->m_Parameters.GetInt( "Radius" ) );
+  filter->SetUseSquareRoot( this->m_Parameters.GetBool( "UseSquareRoot" ) );
+  filter->Update( );
+  this->GetOutput( "Output" )->SetITK( filter->GetOutput( ) );
+}
+
+// eof - $RCSfile$
diff --git a/plugins/cpExtensions/DiscontinuityMapImageFilter.h b/plugins/cpExtensions/DiscontinuityMapImageFilter.h
new file mode 100644 (file)
index 0000000..409854a
--- /dev/null
@@ -0,0 +1,30 @@
+#ifndef __cpPluginscpExtensions__DiscontinuityMapImageFilter__h__
+#define __cpPluginscpExtensions__DiscontinuityMapImageFilter__h__
+
+#include <cpPlugins_cpExtensions_Export.h>
+#include <cpPlugins/Pipeline/ProcessObject.h>
+
+namespace cpPluginscpExtensions
+{
+  /**
+   */
+  class cpPlugins_cpExtensions_EXPORT DiscontinuityMapImageFilter
+    : public cpPlugins::Pipeline::ProcessObject
+  {
+    cpPluginsObject(
+      DiscontinuityMapImageFilter,
+      cpPlugins::Pipeline::ProcessObject,
+      cpExtensions
+      );
+
+  protected:
+    template< class _TInputImage >
+    inline void _GD0( _TInputImage* input );
+  };
+
+} // ecapseman
+
+#endif // __cpPluginscpExtensions__DiscontinuityMapImageFilter__h__
+
+
+// eof - $RCSfile$
index 85cc4945c65fcc1f18eda41c7815cdf2d3cfee6e..bafd384157840616e598de7074195629b39a757b 100644 (file)
@@ -10,12 +10,14 @@ cpPluginscpExtensions::ImageFunctorFilter::
 ImageFunctorFilter( )
   : Superclass( )
 {
-  typedef cpPlugins::DataObject _TFunctor;
-  typedef cpInstances::DataObjects::Image _TImage;
+  /* TODO
+     typedef cpPlugins::DataObject _TFunctor;
+     typedef cpInstances::DataObjects::Image _TImage;
 
-  this->_ConfigureInput< _TImage >( "Input", true, false );
-  this->_ConfigureInput< _TFunctor >( "Functor", true, false );
-  this->_ConfigureOutput< _TImage >( "Output" );
+     this->_ConfigureInput< _TImage >( "Input", true, false );
+     this->_ConfigureInput< _TFunctor >( "Functor", true, false );
+     this->_ConfigureOutput< _TImage >( "Output" );
+  */
 }
 
 // -------------------------------------------------------------------------
@@ -28,12 +30,14 @@ cpPluginscpExtensions::ImageFunctorFilter::
 void cpPluginscpExtensions::ImageFunctorFilter::
 _GenerateData( )
 {
-  auto o = this->GetInputData( "Input" );
-  cpPlugins_Demangle_Image_ScalarPixels_AllDims_1( o, _GD0 )
-    cpPlugins_Demangle_Image_ComplexPixels_AllDims_1( o, _GD0 )
-    cpPlugins_Demangle_Image_ColorPixels_AllDims_1( o, _GD0 )
-    cpPlugins_Demangle_Image_VectorPixels_AllDims_1( o, _GD0 )
-    this->_Error( "Invalid input image pixel type." );
+  /* TODO
+     auto o = this->GetInputData( "Input" );
+     cpPlugins_Demangle_Image_ScalarPixels_AllDims_1( o, _GD0 )
+     cpPlugins_Demangle_Image_ComplexPixels_AllDims_1( o, _GD0 )
+     cpPlugins_Demangle_Image_ColorPixels_AllDims_1( o, _GD0 )
+     cpPlugins_Demangle_Image_VectorPixels_AllDims_1( o, _GD0 )
+     this->_Error( "Invalid input image pixel type." );
+  */
 }
 
 // -------------------------------------------------------------------------
index 8ab94204f9792d11dc80e4d4bf43ac7b543a5d47..f0824670627d9548dcb1a9dd5be4c6cfae079fb2 100644 (file)
@@ -1,5 +1,6 @@
 header #define ITK_MANUAL_INSTANTIATION
 
+tinclude cpExtensions/Algorithms/DiscontinuityMapImageFilter:h|hxx
 tinclude cpExtensions/Algorithms/SkeletonToImageFilter:h|hxx
 tinclude cpExtensions/Algorithms/SkeletonReader:h|hxx
 tinclude cpExtensions/Algorithms/SkeletonWriter:h|hxx
@@ -7,6 +8,8 @@ tinclude cpExtensions/Algorithms/PolyLineParametricPathWriter:h|hxx
 cinclude cpExtensions/DataStructures/PolyLineParametricPath.h
 cinclude cpExtensions/DataStructures/Skeleton.h
 
+instances cpExtensions::Algorithms::DiscontinuityMapImageFilter< itk::Image< #scalar_types#, #pdims# > >
+
 instances cpExtensions::Algorithms::PolyLineParametricPathWriter< cpExtensions::DataStructures::PolyLineParametricPath< #pdims# >, itk::Image< #scalar_types#, #pdims# > >
 
 instances cpExtensions::Algorithms::SkeletonToImageFilter< cpExtensions::DataStructures::Skeleton< #pdims# >, itk::Image< unsigned char, #pdims# > >