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