]> Creatis software - cpPlugins.git/blob - lib/cpExtensions/Algorithms/GaussianDensityImageFilter.h
308b3c5c084f8ccbb2ab232bf21c4828d5b6d490
[cpPlugins.git] / lib / cpExtensions / Algorithms / GaussianDensityImageFilter.h
1 // -------------------------------------------------------------------------
2 // @author Leonardo Florez-Valencia (florez-l@javeriana.edu.co)
3 // -------------------------------------------------------------------------
4
5 #ifndef __CPEXTENSIONS__ALGORITHMS__GAUSSIANDENSITYIMAGEFILTER__H__
6 #define __CPEXTENSIONS__ALGORITHMS__GAUSSIANDENSITYIMAGEFILTER__H__
7
8 #include <itkImageToImageFilter.h>
9 #include <itkImageScanlineConstIterator.h>
10 #include <itkImageScanlineIterator.h>
11 #include <itkProgressReporter.h>
12
13 namespace cpExtensions
14 {
15   namespace Algorithms
16   {
17     /**
18      */
19     template< class _TInputImage, class _TEstimator >
20     class GaussianDensityImageFilter
21       : public itk::ImageToImageFilter< _TInputImage, itk::Image< typename _TEstimator::TScalar, _TInputImage::ImageDimension > >
22     {
23     public:
24       typedef typename _TEstimator::TScalar                         TScalar;
25       typedef itk::Image< TScalar, _TInputImage::ImageDimension >   TOutputImage;
26       typedef GaussianDensityImageFilter                            Self;
27       typedef itk::ImageToImageFilter< _TInputImage, TOutputImage > Superclass;
28       typedef itk::SmartPointer< Self >                             Pointer;
29       typedef itk::SmartPointer< const Self >                       ConstPointer;
30
31       typedef _TInputImage TInputImage;
32       typedef _TEstimator  TEstimator;
33       typedef typename TOutputImage::RegionType TRegion;
34
35     public:
36       itkNewMacro( Self );
37       itkTypeMacro( GaussianDensityImageFilter, itkImageToImageFilter );
38
39       itkGetConstObjectMacro( Estimator, TEstimator );
40       itkSetConstObjectMacro( Estimator, TEstimator );
41
42     protected:
43       GaussianDensityImageFilter( )
44         : Superclass( )
45         {
46         }
47       virtual ~GaussianDensityImageFilter( )
48         {
49         }
50
51       virtual void ThreadedGenerateData(
52         const TRegion& region,
53         itk::ThreadIdType threadId
54         )
55         {
56           const typename TRegion::SizeType& regionSize = region.GetSize( );
57           if( regionSize[ 0 ] == 0 )
58             return;
59           const TInputImage* in = this->GetInput( );
60           TOutputImage* out = this->GetOutput( 0 );
61
62           const size_t nLines = region.GetNumberOfPixels( ) / regionSize[ 0 ];
63           itk::ProgressReporter progress( this, threadId, nLines );
64
65           // Define the iterators
66           itk::ImageScanlineConstIterator< TInputImage > inIt( in, region );
67           itk::ImageScanlineIterator< TOutputImage > outIt( out, region );
68
69           inIt.GoToBegin( );
70           outIt.GoToBegin( );
71           while( !inIt.IsAtEnd( ) )
72           {
73             while( !inIt.IsAtEndOfLine( ) )
74             {
75               if( this->m_Estimator.IsNotNull( ) )
76                 outIt.Set( this->m_Estimator->Density( inIt.Get( ) ) );
77               ++inIt;
78               ++outIt;
79
80             } // elihw
81             inIt.NextLine( );
82             outIt.NextLine( );
83             progress.CompletedPixel( );
84
85           } // elihw
86         }
87
88     private:
89       // Purposely not implemented.
90       GaussianDensityImageFilter( const Self& );
91       void operator=( const Self& );
92
93     protected:
94       typename TEstimator::ConstPointer m_Estimator;
95     };
96
97   } // ecapseman
98
99 } // ecapseman
100
101 #endif // __CPEXTENSIONS__ALGORITHMS__GAUSSIANDENSITYIMAGEFILTER__H__
102
103 // eof - $RCSfile$