]> Creatis software - FrontAlgorithms.git/blob - lib/fpa/Image/SkeletonFilter.hxx
...
[FrontAlgorithms.git] / lib / fpa / Image / SkeletonFilter.hxx
1 #ifndef __fpa__Image__SkeletonFilter__hxx__
2 #define __fpa__Image__SkeletonFilter__hxx__
3
4 #include <fpa/Base/Functors/Inverse.h>
5 #include <fpa/Image/Functors/SimpleNeighborhood.h>
6 #include <itkImageRegionIteratorWithIndex.h>
7
8 // -------------------------------------------------------------------------
9 #define fpa_Image_SkeletonFilter_OutputMacro( o_n, o_t )                \
10   template< class _TImage >                                             \
11   typename fpa::Image::SkeletonFilter< _TImage >::                      \
12   o_t* fpa::Image::SkeletonFilter< _TImage >::                          \
13   Get##o_n( )                                                           \
14   {                                                                     \
15     return(                                                             \
16       dynamic_cast< o_t* >(                                             \
17         this->itk::ProcessObject::GetOutput( this->m_##o_n##Idx )       \
18         )                                                               \
19       );                                                                \
20   }
21
22 fpa_Image_SkeletonFilter_OutputMacro( Skeleton, TSkeleton );
23 fpa_Image_SkeletonFilter_OutputMacro( Marks, TMarks );
24
25 // -------------------------------------------------------------------------
26 template< class _TImage >
27 fpa::Image::SkeletonFilter< _TImage >::
28 SkeletonFilter( )
29   : Superclass( )
30 {
31   typedef typename _TImage::PixelType _TPixel;
32   typedef typename _TImage::Superclass _TImageBase;
33   typedef fpa::Image::Functors::SimpleNeighborhood< _TImageBase > _TNeighFunc;
34   typedef fpa::Base::Functors::Inverse< _TPixel, _TPixel > _TInvFunc;
35
36   unsigned int nOutputs = this->GetNumberOfRequiredOutputs( );
37   this->SetNumberOfRequiredOutputs( nOutputs + 2 );
38   this->m_SkeletonIdx = nOutputs;
39   this->m_MarksIdx = nOutputs + 1;
40
41   typename TSkeleton::Pointer skeleton = TSkeleton::New( );
42   typename TMarks::Pointer marks = TMarks::New( );
43   this->SetNthOutput( this->m_SkeletonIdx, skeleton.GetPointer( ) );
44   this->SetNthOutput( this->m_MarksIdx, marks.GetPointer( ) );
45
46   typename _TNeighFunc::Pointer nfunc = _TNeighFunc::New( );
47   nfunc->SetOrder( 2 );
48   this->SetNeighborhoodFunction( nfunc );
49
50   typename _TInvFunc::Pointer ifunc = _TInvFunc::New( );
51   this->SetCostConversionFunction( ifunc );
52 }
53
54 // -------------------------------------------------------------------------
55 template< class _TImage >
56 fpa::Image::SkeletonFilter< _TImage >::
57 ~SkeletonFilter( )
58 {
59 }
60
61 // -------------------------------------------------------------------------
62 template< class _TImage >
63 void fpa::Image::SkeletonFilter< _TImage >::
64 _BeforeGenerateData( )
65 {
66   this->Superclass::_BeforeGenerateData( );
67   this->m_SkeletonQueue.clear( );
68 }
69
70 // -------------------------------------------------------------------------
71 template< class _TImage >
72 void fpa::Image::SkeletonFilter< _TImage >::
73 _UpdateResult( const _TQueueNode& n )
74 {
75   typedef typename _TSkeletonQueue::value_type _TSkeletonQueueValue;
76
77   this->Superclass::_UpdateResult( n );
78   double d = double( this->GetInput( )->GetPixel( n.Vertex ) );
79   if( d >= double( 0 ) )
80   {
81     // Update skeleton candidates
82     d += double( 1e-5 );
83     double v = double( n.Result ) / ( d * d );
84     this->m_SkeletonQueue.insert( _TSkeletonQueueValue( v, n.Vertex ) );
85
86   } // fi
87 }
88
89 // -------------------------------------------------------------------------
90 template< class _TImage >
91 void fpa::Image::SkeletonFilter< _TImage >::
92 _AfterGenerateData( )
93 {
94   this->Superclass::_AfterGenerateData( );
95
96   std::vector< TIndex > end_points;
97   this->_EndPoints( end_points );
98   this->_Skeleton( end_points );
99 }
100
101 // -------------------------------------------------------------------------
102 template< class _TImage >
103 void fpa::Image::SkeletonFilter< _TImage >::
104 _EndPoints( std::vector< TIndex >& end_points )
105 {
106   typedef itk::ImageRegionIteratorWithIndex< TMarks > _TMarksIt;
107
108   static const double _eps = std::sqrt( double( Self::Dimension + 1 ) );
109   auto input = this->GetInput( );
110   auto mst = this->GetMinimumSpanningTree( );
111   auto marks = this->GetMarks( );
112
113   // Some values
114   marks->SetLargestPossibleRegion( input->GetLargestPossibleRegion( ) );
115   marks->SetRequestedRegion( input->GetRequestedRegion( ) );
116   marks->SetBufferedRegion( input->GetBufferedRegion( ) );
117   marks->SetSpacing( input->GetSpacing( ) );
118   marks->SetOrigin( input->GetOrigin( ) );
119   marks->SetDirection( input->GetDirection( ) );
120   marks->Allocate( );
121   marks->FillBuffer( 0 );
122   auto spac = marks->GetSpacing( );
123
124   // BFS from maximum queue
125   auto region = input->GetRequestedRegion( );
126   while( this->m_SkeletonQueue.size( ) > 0 )
127   {
128     // Get node
129     auto nIt = this->m_SkeletonQueue.begin( );
130     auto n = *nIt;
131     this->m_SkeletonQueue.erase( nIt );
132
133     // Mark it and update end-points
134     unsigned char m = marks->GetPixel( n.second );
135     if( m != 0 )
136       continue;
137     marks->SetPixel( n.second, 1 );
138     end_points.push_back( n.second );
139
140     // Mark path
141     TIndex it = n.second;
142     TIndex p = mst->GetParent( it );
143     while( it != p )
144     {
145       // TODO: tags->SetPixel( it, tags->GetPixel( it ) + 1 );
146       it = p;
147       p = mst->GetParent( it );
148
149     } // elihw
150     // TODO: tags->SetPixel( it, tags->GetPixel( it ) + 1 );
151
152     /* TODO: use mst directly rather than computing paths
153        typename TMST::TPath::Pointer path;
154        mst->GetPath( path, n.second );
155        for( unsigned long i = 0; i < path->GetSize( ); ++i )
156        {
157        TIndex idx = path->GetVertex( i );
158        typename _TImage::PointType cnt;
159        input->TransformIndexToPhysicalPoint( idx, cnt );
160        double r = double( input->GetPixel( idx ) ) * _eps;
161
162        TIndex i0, i1;
163        for( unsigned int d = 0; d < Self::Dimension; ++d )
164        {
165        long off = long( std::ceil( r / double( spac[ d ] ) ) );
166        if( off < 3 )
167        off = 3;
168        i0[ d ] = idx[ d ] - off;
169        i1[ d ] = idx[ d ] + off;
170
171        if( i0[ d ] < region.GetIndex( )[ d ] )
172        i0[ d ] = region.GetIndex( )[ d ];
173
174        if( i1[ d ] >= region.GetIndex( )[ d ] + region.GetSize( )[ d ] )
175        i1[ d ] = region.GetIndex( )[ d ] + region.GetSize( )[ d ] - 1;
176
177        } // rof
178
179        typename _TImage::SizeType size;
180        for( unsigned int d = 0; d < Self::Dimension; ++d )
181        size[ d ] = i1[ d ] - i0[ d ] + 1;
182
183        typename _TImage::RegionType neighRegion;
184        neighRegion.SetIndex( i0 );
185        neighRegion.SetSize( size );
186
187        _TMarksIt mIt( marks, neighRegion );
188        for( mIt.GoToBegin( ); !mIt.IsAtEnd( ); ++mIt )
189        {
190        TIndex w = mIt.GetIndex( );
191        typename _TImage::PointType p;
192        marks->TransformIndexToPhysicalPoint( w, p );
193        mIt.Set( 1 );
194
195        } // rof
196
197        } // rof
198     */
199
200   } // elihw
201 }
202
203 // -------------------------------------------------------------------------
204 template< class _TImage >
205 void fpa::Image::SkeletonFilter< _TImage >::
206 _Skeleton( const std::vector< TIndex >& end_points )
207 {
208   typedef itk::Image< unsigned long, Self::Dimension > _TTagsImage;
209   typedef typename TMST::TPath _TPath;
210
211   auto mst = this->GetMinimumSpanningTree( );
212   auto skeleton = this->GetSkeleton( );
213
214   // Tag branches
215   typename _TTagsImage::Pointer tags = _TTagsImage::New( );
216   tags->SetLargestPossibleRegion( mst->GetLargestPossibleRegion( ) );
217   tags->SetRequestedRegion( mst->GetRequestedRegion( ) );
218   tags->SetBufferedRegion( mst->GetBufferedRegion( ) );
219   tags->Allocate( );
220   tags->FillBuffer( 0 );
221   for( auto eIt = end_points.begin( ); eIt != end_points.end( ); ++eIt )
222   {
223     TIndex it = *eIt;
224     TIndex p = mst->GetParent( it );
225     while( it != p )
226     {
227       tags->SetPixel( it, tags->GetPixel( it ) + 1 );
228       it = p;
229       p = mst->GetParent( it );
230
231     } // elihw
232     tags->SetPixel( it, tags->GetPixel( it ) + 1 );
233
234   } // rof
235
236   // Build paths (branches)
237   for( auto eIt = end_points.begin( ); eIt != end_points.end( ); ++eIt )
238   {
239     TIndex it = *eIt;
240     TIndex p = mst->GetParent( it );
241     TIndex sIdx = *eIt;
242     typename _TPath::Pointer path = _TPath::New( );
243     path->SetReferenceImage( mst );
244     while( it != p )
245     {
246       if( tags->GetPixel( sIdx ) != tags->GetPixel( it ) )
247       {
248         // Ok, a new branch has been added
249         path->AddVertex( it );
250         skeleton->AddBranch( path );
251
252         // Update a new starting index
253         path = _TPath::New( );
254         path->SetReferenceImage( mst );
255         sIdx = it;
256       }
257       else
258         path->AddVertex( it );
259       it = p;
260       p = mst->GetParent( it );
261
262     } // elihw
263
264     // Finally, add last branch
265     path->AddVertex( it );
266     skeleton->AddBranch( path );
267
268   } // rof
269 }
270
271 #endif // __fpa__Image__SkeletonFilter__hxx__
272
273 // eof - $RCSfile$