]> Creatis software - FrontAlgorithms.git/blob - lib/fpa/Image/SkeletonFilter.hxx
0f33388f80d04b8202fa3791a787beb66531e4ba
[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   std::vector< TIndex > end_points;
171   this->_EndPoints(
172     end_points,
173     this->m_DistanceMap->GetOutput( ),
174     dijkstra->GetMinimumSpanningTree( ),
175     dijkstra->GetSkeletonQueue( )
176     );
177 }
178
179 // -------------------------------------------------------------------------
180 template< class _TInputImage, class _TDistanceMap >
181 template< class _TMarksPointer >
182 void fpa::Image::SkeletonFilter< _TInputImage, _TDistanceMap >::
183 _MarkSphere(
184   _TMarksPointer& marks, const TOutputImage* dmap, const TIndex& center
185   )
186 {
187   typedef typename _TMarksPointer::ObjectType          _TMarks;
188   typedef itk::ImageRegionIteratorWithIndex< _TMarks > _TMarksIt;
189
190   static const double _eps = std::sqrt( double( Self::Dimension + 1 ) );
191   typename _TMarks::SpacingType spac = dmap->GetSpacing( );
192   typename _TMarks::RegionType region = dmap->GetRequestedRegion( );
193
194   typename _TMarks::PointType cnt;
195   dmap->TransformIndexToPhysicalPoint( center, cnt );
196   double r = double( dmap->GetPixel( center ) ) * _eps;
197
198   TIndex i0, i1;
199   for( unsigned int d = 0; d < Self::Dimension; ++d )
200   {
201     long off = long( std::ceil( r / double( spac[ d ] ) ) );
202     if( off < 3 )
203       off = 3;
204     i0[ d ] = center[ d ] - off;
205     i1[ d ] = center[ d ] + off;
206
207     if( i0[ d ] < region.GetIndex( )[ d ] )
208       i0[ d ] = region.GetIndex( )[ d ];
209
210     if( i1[ d ] >= region.GetIndex( )[ d ] + region.GetSize( )[ d ] )
211       i1[ d ] = region.GetIndex( )[ d ] + region.GetSize( )[ d ] - 1;
212
213   } // rof
214
215   typename _TMarks::SizeType size;
216   for( unsigned int d = 0; d < Self::Dimension; ++d )
217     size[ d ] = i1[ d ] - i0[ d ] + 1;
218
219   typename _TMarks::RegionType neighRegion;
220   neighRegion.SetIndex( i0 );
221   neighRegion.SetSize( size );
222
223   _TMarksIt mIt( marks, neighRegion );
224   for( mIt.GoToBegin( ); !mIt.IsAtEnd( ); ++mIt )
225     mIt.Set( true );
226 }
227
228 // -------------------------------------------------------------------------
229 template< class _TInputImage, class _TDistanceMap >
230 void fpa::Image::SkeletonFilter< _TInputImage, _TDistanceMap >::
231 _EndPoints(
232   std::vector< TIndex >& end_points,
233   const TOutputImage* dmap,
234   const _TMST* mst,
235   const _TSkeletonQueue& queue
236   )
237 {
238   typedef typename _TSkeletonQueue::value_type _TSkeletonQueueValue;
239
240   // Some values
241   typedef itk::Image< bool, Self::Dimension > _TMarks;
242   typename _TMarks::Pointer marks = _TMarks::New( );
243   marks->SetLargestPossibleRegion( mst->GetLargestPossibleRegion( ) );
244   marks->SetRequestedRegion( mst->GetRequestedRegion( ) );
245   marks->SetBufferedRegion( mst->GetBufferedRegion( ) );
246   marks->SetSpacing( mst->GetSpacing( ) );
247   marks->SetOrigin( mst->GetOrigin( ) );
248   marks->SetDirection( mst->GetDirection( ) );
249   marks->Allocate( );
250   marks->FillBuffer( false );
251
252   // Get candidates in maximum to minimum iteration
253   typename _TSkeletonQueue::const_reverse_iterator nIt = queue.rbegin( );
254   for( ; nIt != queue.rend( ); ++nIt )
255   {
256     // Mark it and update end-points
257     if( !( marks->GetPixel( nIt->second ) ) )
258     {
259       marks->SetPixel( nIt->second, true );
260       end_points.push_back( nIt->second );
261
262       // Mark path
263       TIndex it = nIt->second;
264       TIndex p = mst->GetParent( it );
265       while( it != p )
266       {
267         this->_MarkSphere( marks, dmap, it );
268         it = p;
269         p = mst->GetParent( it );
270
271       } // elihw
272       this->_MarkSphere( marks, dmap, it );
273
274     } // fi
275
276   } // rof
277 }
278
279
280
281 /* TODO
282    this->m_DistanceMap = TDistanceMap::New( );
283    this->m_DistanceMap->InsideIsPositiveOn( );
284    this->m_DistanceMap->SquaredDistanceOff( );
285    this->m_DistanceMap->UseImageSpacingOn( );
286
287    template< class _TImage, class _TScalar >
288    void fpa::Image::SkeletonFilter< _TImage, _TScalar >::
289    _Skeleton( const std::vector< TVertex >& end_points, _TAdjacencies& A )
290    {
291    typedef typename TSkeleton::TPath _TPath;
292    typedef itk::Image< unsigned long, TImage::ImageDimension > _TTagsImage;
293
294    TMST* mst = this->GetMinimumSpanningTree( );
295    TSkeleton* skeleton = this->GetSkeleton( );
296
297    // Tag branches
298    typename _TTagsImage::Pointer tags = _TTagsImage::New( );
299    tags->SetLargestPossibleRegion( mst->GetLargestPossibleRegion( ) );
300    tags->SetRequestedRegion( mst->GetRequestedRegion( ) );
301    tags->SetBufferedRegion( mst->GetBufferedRegion( ) );
302    tags->Allocate( );
303    tags->FillBuffer( 0 );
304    typename std::vector< TVertex >::const_iterator eIt = end_points.begin( );
305    for( ; eIt != end_points.end( ); ++eIt )
306    {
307    TVertex it = *eIt;
308    TVertex p = mst->GetParent( it );
309    while( it != p )
310    {
311    tags->SetPixel( it, tags->GetPixel( it ) + 1 );
312    it = p;
313    p = mst->GetParent( it );
314    
315    } // elihw
316    tags->SetPixel( it, tags->GetPixel( it ) + 1 );
317
318    } // rof
319
320    // Build paths (branches)
321    eIt = end_points.begin( );
322    for( ; eIt != end_points.end( ); ++eIt )
323    {
324    TVertex it = *eIt;
325    TVertex p = mst->GetParent( it );
326    TVertex sIdx = it;
327    typename _TPath::Pointer path = _TPath::New( );
328    path->SetReferenceImage( mst );
329    while( it != p )
330    {
331    if( tags->GetPixel( sIdx ) != tags->GetPixel( it ) )
332    {
333    // Ok, a new branch has been added
334    path->AddVertex( it );
335    skeleton->AddBranch( path );
336    
337    // Update a new starting index
338    path = _TPath::New( );
339    path->SetReferenceImage( mst );
340    sIdx = it;
341    }
342    else
343    path->AddVertex( it );
344    it = p;
345    p = mst->GetParent( it );
346
347    } // elihw
348    
349    // Finally, add last branch
350    path->AddVertex( it );
351    skeleton->AddBranch( path );
352    
353    } // rof
354    }
355 */
356
357 #endif // __fpa__Image__SkeletonFilter__hxx__
358
359 // eof - $RCSfile$