]> Creatis software - FrontAlgorithms.git/blob - lib/fpa/Filters/Image/ExtractSkeleton.h
...
[FrontAlgorithms.git] / lib / fpa / Filters / Image / ExtractSkeleton.h
1 // =========================================================================
2 // @author Leonardo Florez Valencia
3 // @email florez-l@javeriana.edu.co
4 // =========================================================================
5 #ifndef __fpa__Filters__Image__ExtractSkeleton__h__
6 #define __fpa__Filters__Image__ExtractSkeleton__h__
7
8 #include <map>
9 #include <itkProcessObject.h>
10 #include <itkSignedMaurerDistanceMapImageFilter.h>
11 #include <fpa/Filters/Image/Dijkstra.h>
12 #include <ivq/ITK/ImageSkeleton.h>
13
14 namespace fpa
15 {
16   namespace Filters
17   {
18     namespace Image
19     {
20       /**
21        */
22       template< class _TInputImage, class _TDistanceMap = itk::SignedMaurerDistanceMapImageFilter< _TInputImage, itk::Image< double, _TInputImage::ImageDimension > > >
23       class ExtractSkeleton
24         : public itk::ProcessObject
25       {
26       public:
27         typedef ExtractSkeleton                 Self;
28         typedef itk::ProcessObject              Superclass;
29         typedef itk::SmartPointer< Self >       Pointer;
30         typedef itk::SmartPointer< const Self > ConstPointer;
31
32         typedef _TInputImage  TInputImage;
33         typedef _TDistanceMap TDistanceMap;
34         itkStaticConstMacro(
35           Dimension,
36           unsigned int,
37           TInputImage::ImageDimension
38           );
39
40         typedef typename TDistanceMap::OutputImageType TOutputImage;
41         typedef typename TInputImage::IndexType        TIndex;
42         typedef typename TOutputImage::PixelType       TScalar;
43
44         typedef ivq::ITK::ImageSkeleton< Self::Dimension > TSkeleton;
45
46       protected:
47         typedef std::multimap< TScalar, TIndex > _TSkeletonQueue;
48
49         /**
50          */
51         class _TDijkstra
52           : public fpa::Filters::Image::Dijkstra< TOutputImage, TOutputImage >
53         {
54         public:
55           typedef fpa::Filters::Image::Dijkstra< TOutputImage, TOutputImage > Superclass;
56           typedef _TDijkstra                      Self;
57           typedef itk::SmartPointer< Self >       Pointer;
58           typedef itk::SmartPointer< const Self > ConstPointer;
59
60           typedef typename Superclass::TNode TNode;
61           typedef typename Superclass::TMST  TMST;
62
63         public:
64           itkNewMacro( Self );
65           itkTypeMacro( _TDijkstra, fpa::Image::Dijkstra );
66
67           itkGetConstReferenceMacro( SkeletonQueue, _TSkeletonQueue );
68
69         protected:
70           _TDijkstra( );
71           virtual ~_TDijkstra( );
72
73           virtual void _BeforeGenerateData( ) override;
74           virtual void _PostComputeOutputValue( TNode& n ) override;
75
76         private:
77           // Purposely not implemented
78           _TDijkstra( const Self& other );
79           Self& operator=( const Self& other );
80
81         protected:
82           _TSkeletonQueue m_SkeletonQueue;
83         };
84         typedef typename _TDijkstra::TMST _TMST;
85
86       public:
87         itkNewMacro( Self );
88         itkTypeMacro( fpa::Image::ExtractSkeleton, fpa::Image::Dijkstra );
89
90         itkBooleanMacro( SeedFromMaximumDistance );
91         itkGetConstMacro( SeedFromMaximumDistance, bool );
92         itkSetMacro( SeedFromMaximumDistance, bool );
93
94         itkGetConstObjectMacro( DistanceMap, TDistanceMap );
95         itkGetObjectMacro( DistanceMap, TDistanceMap );
96
97         itkGetConstMacro( Seed, TIndex );
98         itkSetMacro( Seed, TIndex );
99
100         itkGetConstMacro( EndPoints, std::vector< TIndex > );
101
102       public:
103         virtual itk::ModifiedTimeType GetMTime( ) const override;
104
105         void SetInput( TInputImage* input );
106         TInputImage* GetInput( );
107         const TInputImage* GetInput( ) const;
108
109         TSkeleton* GetOutput( );
110         const TSkeleton* GetOutput( ) const;
111
112       protected:
113         ExtractSkeleton( );
114         virtual ~ExtractSkeleton( );
115
116         virtual void GenerateData( ) override;
117
118         template< class _TMarksPointer >
119         void _MarkSphere(
120           _TMarksPointer& marks,
121           const TOutputImage* dmap,
122           const TIndex& center
123           );
124
125         void _EndPoints(
126           std::vector< TIndex >& end_points,
127           const TOutputImage* dmap,
128           const _TMST* mst,
129           const _TSkeletonQueue& queue
130           );
131
132       private:
133         // Purposely not implemented.
134         ExtractSkeleton( const Self& other );
135         Self& operator=( const Self& other );
136
137       protected:
138         typename TDistanceMap::Pointer m_DistanceMap;
139         bool   m_SeedFromMaximumDistance;
140         TIndex m_Seed;
141         std::vector< TIndex >  m_EndPoints;
142       };
143
144     } // ecapseman
145
146   } // ecapseman
147
148 } // ecapseman
149
150 #ifndef ITK_MANUAL_INSTANTIATION
151 #  include <fpa/Filters/Image/ExtractSkeleton.hxx>
152 #endif // ITK_MANUAL_INSTANTIATION
153
154 #endif // __fpa__Filters__Image__ExtractSkeleton__h__
155
156 // eof - $RCSfile$