]> Creatis software - FrontAlgorithms.git/blob - lib/fpa/Image/SkeletonFilter.hxx
224829b93f5a37119a5ceea5f3421b26814b9f43
[FrontAlgorithms.git] / lib / fpa / Image / SkeletonFilter.hxx
1 // =========================================================================
2 // @author Leonardo Florez Valencia
3 // @email florez-l@javeriana.edu.co
4 // =========================================================================
5
6 #ifndef __fpa__Image__SkeletonFilter__hxx__
7 #define __fpa__Image__SkeletonFilter__hxx__
8
9 #include <itkImage.h>
10 #include <itkImageRegionIteratorWithIndex.h>
11 #include <itkMinimumMaximumImageCalculator.h>
12 #include <fpa/Image/Functors/Dijkstra/Invert.h>
13
14 // -------------------------------------------------------------------------
15 template< class _TInputImage, class _TDistanceMap >
16 fpa::Image::SkeletonFilter< _TInputImage, _TDistanceMap >::_TDijkstra::
17 _TDijkstra( )
18   : Superclass( )
19 {
20   // Prepare weight function
21   typedef fpa::Image::Functors::Dijkstra::Invert< TOutputImage, TScalar > _TWeight;
22   typename _TWeight::Pointer weight = _TWeight::New( );
23   weight->SetAlpha( 1 );
24   weight->SetBeta( 1 );
25   this->SetWeightFunction( weight );
26 }
27
28 // -------------------------------------------------------------------------
29 template< class _TInputImage, class _TDistanceMap >
30 fpa::Image::SkeletonFilter< _TInputImage, _TDistanceMap >::_TDijkstra::
31 ~_TDijkstra( )
32 {
33 }
34
35 // -------------------------------------------------------------------------
36 template< class _TInputImage, class _TDistanceMap >
37 void fpa::Image::SkeletonFilter< _TInputImage, _TDistanceMap >::_TDijkstra::
38 _BeforeGenerateData( )
39 {
40   this->Superclass::_BeforeGenerateData( );
41   this->m_SkeletonQueue.clear( );
42 }
43
44 // -------------------------------------------------------------------------
45 template< class _TInputImage, class _TDistanceMap >
46 void fpa::Image::SkeletonFilter< _TInputImage, _TDistanceMap >::_TDijkstra::
47 _UpdateOutputValue( const TNode& n )
48 {
49   typedef typename _TSkeletonQueue::value_type _TSkeletonQueueValue;
50
51   this->Superclass::_UpdateOutputValue( n );
52   double d = double( this->GetInput( )->GetPixel( n.Vertex ) );
53   if( d >= double( 0 ) )
54   {
55     // Update skeleton candidates
56     d += double( 1e-5 );
57     double v = double( n.Value ) / ( d * d );
58     this->m_SkeletonQueue.insert( _TSkeletonQueueValue( v, n.Vertex ) );
59
60   } // fi
61 }
62
63 // -------------------------------------------------------------------------
64 template< class _TInputImage, class _TDistanceMap >
65 itk::ModifiedTimeType
66 fpa::Image::SkeletonFilter< _TInputImage, _TDistanceMap >::
67 GetMTime( ) const
68 {
69   itk::ModifiedTimeType t = this->Superclass::GetMTime( );
70   itk::ModifiedTimeType q = this->m_DistanceMap->GetMTime( );
71   return( ( q < t )? q: t );
72 }
73
74 // -------------------------------------------------------------------------
75 template< class _TInputImage, class _TDistanceMap >
76 void fpa::Image::SkeletonFilter< _TInputImage, _TDistanceMap >::
77 SetInput( TInputImage* input )
78 {
79   this->Superclass::SetNthInput( 0, input );
80 }
81
82 // -------------------------------------------------------------------------
83 template< class _TInputImage, class _TDistanceMap >
84 typename fpa::Image::SkeletonFilter< _TInputImage, _TDistanceMap >::
85 TInputImage* fpa::Image::SkeletonFilter< _TInputImage, _TDistanceMap >::
86 GetInput( )
87 {
88   return( dynamic_cast< TInputImage* >( this->Superclass::GetInput( 0 ) ) );
89 }
90
91 // -------------------------------------------------------------------------
92 template< class _TInputImage, class _TDistanceMap >
93 const typename fpa::Image::SkeletonFilter< _TInputImage, _TDistanceMap >::
94 TInputImage* fpa::Image::SkeletonFilter< _TInputImage, _TDistanceMap >::
95 GetInput( ) const
96 {
97   return(
98     dynamic_cast< const TInputImage* >( this->Superclass::GetInput( 0 ) )
99     );
100 }
101
102 // -------------------------------------------------------------------------
103 template< class _TInputImage, class _TDistanceMap >
104 typename fpa::Image::SkeletonFilter< _TInputImage, _TDistanceMap >::
105 TSkeleton* fpa::Image::SkeletonFilter< _TInputImage, _TDistanceMap >::
106 GetOutput( )
107 {
108   return( dynamic_cast< TSkeleton* >( this->Superclass::GetOutput( 0 ) ) );
109 }
110
111 // -------------------------------------------------------------------------
112 template< class _TInputImage, class _TDistanceMap >
113 const typename fpa::Image::SkeletonFilter< _TInputImage, _TDistanceMap >::
114 TSkeleton* fpa::Image::SkeletonFilter< _TInputImage, _TDistanceMap >::
115 GetOutput( ) const
116 {
117   return(
118     dynamic_cast< const TSkeleton* >( this->Superclass::GetOutput( 0 ) )
119     );
120 }
121
122 // -------------------------------------------------------------------------
123 template< class _TInputImage, class _TDistanceMap >
124 fpa::Image::SkeletonFilter< _TInputImage, _TDistanceMap >::
125 SkeletonFilter( )
126   : Superclass( ),
127     m_SeedFromMaximumDistance( false )
128 {
129   this->SetNumberOfRequiredInputs( 1 );
130   this->SetNumberOfRequiredOutputs( 1 );
131   this->SetNthOutput( 0, TSkeleton::New( ) );
132
133   this->m_DistanceMap = TDistanceMap::New( );
134 }
135
136 // -------------------------------------------------------------------------
137 template< class _TInputImage, class _TDistanceMap >
138 fpa::Image::SkeletonFilter< _TInputImage, _TDistanceMap >::
139 ~SkeletonFilter( )
140 {
141 }
142
143 // -------------------------------------------------------------------------
144 template< class _TInputImage, class _TDistanceMap >
145 void fpa::Image::SkeletonFilter< _TInputImage, _TDistanceMap >::
146 GenerateData( )
147 {
148   // Update distance map
149   this->m_DistanceMap->SetInput( this->GetInput( ) );
150   this->m_DistanceMap->Update( );
151
152   // Correct seed
153   if( this->m_SeedFromMaximumDistance )
154   {
155     typedef itk::MinimumMaximumImageCalculator< TOutputImage > _TMinMax;
156     typename _TMinMax::Pointer minmax = _TMinMax::New( );
157     minmax->SetImage( this->m_DistanceMap->GetOutput( ) );
158     minmax->Compute( );
159     this->m_Seed = minmax->GetIndexOfMaximum( );
160
161   } // fi
162
163   // Compute MST
164   typename _TDijkstra::Pointer dijkstra = _TDijkstra::New( );
165   dijkstra->SetInput( this->m_DistanceMap->GetOutput( ) );
166   dijkstra->AddSeed( this->m_Seed );
167   dijkstra->Update( );
168
169   // Compute end-points
170   const _TMST* mst = dijkstra->GetMinimumSpanningTree( );
171   std::vector< TIndex > end_points;
172   this->_EndPoints(
173     end_points, this->m_DistanceMap->GetOutput( ),
174     mst, dijkstra->GetSkeletonQueue( )
175     );
176
177   // Compute symbolic branches
178   typedef std::map< TIndex, TIndex, typename TIndex::LexicographicCompare > _TTags;
179   _TTags tags, branches;
180   typename std::vector< TIndex >::const_iterator eIt = end_points.begin( );
181   for( ; eIt != end_points.end( ); ++eIt )
182   {
183     // Tag path
184     TIndex it = *eIt;
185     TIndex p = mst->GetParent( it );
186     typename _TTags::iterator bIt = tags.end( );
187     while( it != p && bIt == tags.end( ) )
188     {
189       typename _TTags::iterator tIt = tags.find( it );
190       if( tIt != tags.end( ) )
191       {
192         // Ok, a bifurcation point has been found
193         // branch1: tIt->second <-> it (ok)
194         branches[ tIt->second ] = it;
195
196         // branch2: *eit <-> it (ok)
197         branches[ *eIt ] = it;
198
199         // branch3: it <-> until next bifurcation
200         bIt = tIt;
201       }
202       else
203         tags[ it ] = *eIt;
204       it = p;
205       p = mst->GetParent( it );
206
207     } // elihw
208     if( bIt != tags.end( ) )
209     {
210       TIndex pTag = bIt->second;
211       TIndex nTag = bIt->first;
212       it = bIt->first;
213       p = it;
214       while( tags[ it ] == pTag )
215       {
216         tags[ it ] = nTag;
217         p = it;
218         it = mst->GetParent( it );
219
220       } // elihw
221       branches[ bIt->first ] = p;
222     }
223     else
224     {
225       tags[ it ] = *eIt;
226       branches[ *eIt ] = it;
227
228     } // fi
229
230   } // rof
231
232   // Fill full branches
233   typedef typename _TMST::TVertices _TVertices;
234   typedef typename TSkeleton::TPath _TPath;
235
236   TSkeleton* sk = this->GetOutput( );
237   typename _TTags::const_iterator bIt = branches.begin( );
238   for( ; bIt != branches.end( ); ++bIt )
239   {
240     _TVertices v = mst->GetPath( bIt->first, bIt->second );
241     typename _TPath::Pointer path = _TPath::New( );
242     path->SetReferenceImage( this->GetInput( ) );
243     typename _TVertices::const_reverse_iterator vIt = v.rbegin( );
244     for( ; vIt != v.rend( ); ++vIt )
245       path->AddVertex( *vIt );
246     sk->AddBranch( path );
247
248   } // rof
249 }
250
251 // -------------------------------------------------------------------------
252 template< class _TInputImage, class _TDistanceMap >
253 template< class _TMarksPointer >
254 void fpa::Image::SkeletonFilter< _TInputImage, _TDistanceMap >::
255 _MarkSphere(
256   _TMarksPointer& marks, const TOutputImage* dmap, const TIndex& center
257   )
258 {
259   typedef typename _TMarksPointer::ObjectType          _TMarks;
260   typedef itk::ImageRegionIteratorWithIndex< _TMarks > _TMarksIt;
261
262   static const double _eps = std::sqrt( double( Self::Dimension + 1 ) );
263   typename _TMarks::SpacingType spac = dmap->GetSpacing( );
264   typename _TMarks::RegionType region = dmap->GetRequestedRegion( );
265
266   typename _TMarks::PointType cnt;
267   dmap->TransformIndexToPhysicalPoint( center, cnt );
268   double r = double( dmap->GetPixel( center ) ) * _eps;
269
270   TIndex i0, i1;
271   for( unsigned int d = 0; d < Self::Dimension; ++d )
272   {
273     long off = long( std::ceil( r / double( spac[ d ] ) ) );
274     if( off < 3 )
275       off = 3;
276     i0[ d ] = center[ d ] - off;
277     i1[ d ] = center[ d ] + off;
278
279     if( i0[ d ] < region.GetIndex( )[ d ] )
280       i0[ d ] = region.GetIndex( )[ d ];
281
282     if( i1[ d ] >= region.GetIndex( )[ d ] + region.GetSize( )[ d ] )
283       i1[ d ] = region.GetIndex( )[ d ] + region.GetSize( )[ d ] - 1;
284
285   } // rof
286
287   typename _TMarks::SizeType size;
288   for( unsigned int d = 0; d < Self::Dimension; ++d )
289     size[ d ] = i1[ d ] - i0[ d ] + 1;
290
291   typename _TMarks::RegionType neighRegion;
292   neighRegion.SetIndex( i0 );
293   neighRegion.SetSize( size );
294
295   _TMarksIt mIt( marks, neighRegion );
296   for( mIt.GoToBegin( ); !mIt.IsAtEnd( ); ++mIt )
297     mIt.Set( true );
298 }
299
300 // -------------------------------------------------------------------------
301 template< class _TInputImage, class _TDistanceMap >
302 void fpa::Image::SkeletonFilter< _TInputImage, _TDistanceMap >::
303 _EndPoints(
304   std::vector< TIndex >& end_points,
305   const TOutputImage* dmap,
306   const _TMST* mst,
307   const _TSkeletonQueue& queue
308   )
309 {
310   typedef typename _TSkeletonQueue::value_type _TSkeletonQueueValue;
311
312   // Some values
313   typedef itk::Image< bool, Self::Dimension > _TMarks;
314   typename _TMarks::Pointer marks = _TMarks::New( );
315   marks->SetLargestPossibleRegion( mst->GetLargestPossibleRegion( ) );
316   marks->SetRequestedRegion( mst->GetRequestedRegion( ) );
317   marks->SetBufferedRegion( mst->GetBufferedRegion( ) );
318   marks->SetSpacing( mst->GetSpacing( ) );
319   marks->SetOrigin( mst->GetOrigin( ) );
320   marks->SetDirection( mst->GetDirection( ) );
321   marks->Allocate( );
322   marks->FillBuffer( false );
323
324   // Get candidates in maximum to minimum iteration
325   typename _TSkeletonQueue::const_reverse_iterator nIt = queue.rbegin( );
326   for( ; nIt != queue.rend( ); ++nIt )
327   {
328     // Mark it and update end-points
329     if( !( marks->GetPixel( nIt->second ) ) )
330     {
331       marks->SetPixel( nIt->second, true );
332       end_points.push_back( nIt->second );
333
334       // Mark path
335       TIndex it = nIt->second;
336       TIndex p = mst->GetParent( it );
337       while( it != p )
338       {
339         this->_MarkSphere( marks, dmap, it );
340         it = p;
341         p = mst->GetParent( it );
342
343       } // elihw
344       this->_MarkSphere( marks, dmap, it );
345
346     } // fi
347
348   } // rof
349 }
350
351 #endif // __fpa__Image__SkeletonFilter__hxx__
352
353 // eof - $RCSfile$