]> Creatis software - cpPlugins.git/blob - lib/ivq/ITK/RasterContourFilter.h
...
[cpPlugins.git] / lib / ivq / ITK / RasterContourFilter.h
1 // =======================================================================
2 // @author: Leonardo Florez-Valencia
3 // @email: florez-l@javeriana.edu.co
4 // =======================================================================
5 #ifndef __ivq__ITK__RasterContourFilter__h__
6 #define __ivq__ITK__RasterContourFilter__h__
7
8 #include <deque>
9 #include <itkImageSource.h>
10
11 namespace ivq
12 {
13   namespace ITK
14   {
15     /**
16      */
17     template< class _TImage >
18     class RasterContourFilter
19       : public itk::ImageSource< _TImage >
20     {
21     public:
22       // Basic types
23       typedef RasterContourFilter             Self;
24       typedef itk::ImageSource< _TImage >     Superclass;
25       typedef itk::SmartPointer< Self >       Pointer;
26       typedef itk::SmartPointer< const Self > ConstPointer;
27
28       typedef _TImage TImage;
29       typedef typename _TImage::IndexType  TIndex;
30       typedef typename _TImage::PixelType  TPixel;
31       typedef typename _TImage::PointType  TPoint;
32       typedef typename _TImage::RegionType TRegion;
33       typedef itk::ImageBase< 2 > TImageBase;
34
35     public:
36       itkNewMacro( Self );
37       itkTypeMacro( ivq::ITK::RasterContourFilter, itk::ImageSource );
38
39       itkGetConstObjectMacro( Template, TImageBase );
40       itkGetConstMacro( InsideValue, TPixel );
41       itkGetConstMacro( OutsideValue, TPixel );
42
43       itkSetConstObjectMacro( Template, TImageBase );
44       itkSetMacro( InsideValue, TPixel );
45       itkSetMacro( OutsideValue, TPixel );
46
47     public:
48       void AddPoint( double x, double y );
49       void AddPoint( double p[ 2 ] );
50
51       template< class _TPoint >
52       inline void AddPoint( const _TPoint& p );
53
54       void ClearPoints( );
55
56     protected:
57       RasterContourFilter( );
58       virtual ~RasterContourFilter( );
59
60       virtual void AllocateOutputs( ) override;
61       virtual void BeforeThreadedGenerateData( ) override;
62       virtual void AfterThreadedGenerateData( ) override;
63       virtual void ThreadedGenerateData(
64         const TRegion& region, itk::ThreadIdType id
65         ) override;
66
67     private:
68       // Purposely not implemented
69       RasterContourFilter( const Self& );
70       void operator=( const Self& );
71
72     protected:
73       std::deque< TPoint > m_Contour;
74       std::deque< TIndex > m_Polygon;
75       TRegion m_ROI;
76       typename TImageBase::ConstPointer m_Template;
77       TPixel m_InsideValue;
78       TPixel m_OutsideValue;
79     };
80
81   } // ecapseman
82
83 } // ecapseman
84
85 // -------------------------------------------------------------------------
86 template< class _TImage >
87 template< class _TPoint >
88 void ivq::ITK::RasterContourFilter< _TImage >::
89 AddPoint( const _TPoint& p )
90 {
91   TPoint pnt;
92   pnt[ 0 ] = p[ 0 ];
93   pnt[ 1 ] = p[ 2 ];
94   this->m_Contour.push_back( pnt );
95   this->Modified( );
96 }
97
98 // -------------------------------------------------------------------------
99 #ifndef ITK_MANUAL_INSTANTIATION
100 #  include <ivq/ITK/RasterContourFilter.hxx>
101 #endif // ITK_MANUAL_INSTANTIATION
102
103 #endif // __ivq__ITK__RasterContourFilter__h__
104
105 // eof - $RCSfile$