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