]> Creatis software - FrontAlgorithms.git/blob - lib/fpa/Base/ExtractEndPointsAndBifurcationsFromMinimumSpanningTree.hxx
...
[FrontAlgorithms.git] / lib / fpa / Base / ExtractEndPointsAndBifurcationsFromMinimumSpanningTree.hxx
1 #ifndef __FPA__BASE__EXTRACTENDPOINTSANDBIFURCATIONSFROMMINIMUMSPANNINGTREE__HXX__
2 #define __FPA__BASE__EXTRACTENDPOINTSANDBIFURCATIONSFROMMINIMUMSPANNINGTREE__HXX__
3
4 #include <cmath>
5 #include <map>
6
7 // -------------------------------------------------------------------------
8 template< class _TMST >
9 const typename
10 fpa::Base::ExtractEndPointsAndBifurcationsFromMinimumSpanningTree< _TMST >::
11 TMinimumSpanningTree*
12 fpa::Base::ExtractEndPointsAndBifurcationsFromMinimumSpanningTree< _TMST >::
13 GetMinimumSpanningTree( )
14 {
15   return(
16     dynamic_cast< const _TMST* >( this->itk::ProcessObject::GetInput( 0 ) )
17     );
18 }
19
20 // -------------------------------------------------------------------------
21 template< class _TMST >
22 void
23 fpa::Base::ExtractEndPointsAndBifurcationsFromMinimumSpanningTree< _TMST >::
24 SetMinimumSpanningTree( TMinimumSpanningTree* mst )
25 {
26   this->itk::ProcessObject::SetNthInput( 0, const_cast< _TMST* >( mst ) );
27 }
28
29 // -------------------------------------------------------------------------
30 template< class _TMST >
31 typename
32 fpa::Base::ExtractEndPointsAndBifurcationsFromMinimumSpanningTree< _TMST >::
33 TVertices*
34 fpa::Base::ExtractEndPointsAndBifurcationsFromMinimumSpanningTree< _TMST >::
35 GetEndPoints( )
36 {
37   return(
38     dynamic_cast< TVertices* >( this->itk::ProcessObject::GetOutput( 0 ) )
39     );
40 }
41
42 // -------------------------------------------------------------------------
43 template< class _TMST >
44 typename
45 fpa::Base::ExtractEndPointsAndBifurcationsFromMinimumSpanningTree< _TMST >::
46 TVertices*
47 fpa::Base::ExtractEndPointsAndBifurcationsFromMinimumSpanningTree< _TMST >::
48 GetBifurcations( )
49 {
50   return(
51     dynamic_cast< TVertices* >( this->itk::ProcessObject::GetOutput( 1 ) )
52     );
53 }
54
55 // -------------------------------------------------------------------------
56 template< class _TMST >
57 typename
58 fpa::Base::ExtractEndPointsAndBifurcationsFromMinimumSpanningTree< _TMST >::
59 TVertices*
60 fpa::Base::ExtractEndPointsAndBifurcationsFromMinimumSpanningTree< _TMST >::
61 GetCollisions( )
62 {
63   return(
64     dynamic_cast< TVertices* >( this->itk::ProcessObject::GetOutput( 2 ) )
65     );
66 }
67
68 // -------------------------------------------------------------------------
69 template< class _TMST >
70 fpa::Base::ExtractEndPointsAndBifurcationsFromMinimumSpanningTree< _TMST >::
71 ExtractEndPointsAndBifurcationsFromMinimumSpanningTree( )
72   : Superclass( )
73 {
74   this->SetNumberOfRequiredInputs( 1 );
75   this->SetNumberOfRequiredOutputs( 4 );
76   typename TVertices::Pointer ep = TVertices::New( );
77   typename TVertices::Pointer bf = TVertices::New( );
78   typename TVertices::Pointer co = TVertices::New( );
79   this->itk::ProcessObject::SetNthOutput( 0, ep.GetPointer( ) );
80   this->itk::ProcessObject::SetNthOutput( 1, bf.GetPointer( ) );
81   this->itk::ProcessObject::SetNthOutput( 2, co.GetPointer( ) );
82 }
83
84 // -------------------------------------------------------------------------
85 template< class _TMST >
86 fpa::Base::ExtractEndPointsAndBifurcationsFromMinimumSpanningTree< _TMST >::
87 ~ExtractEndPointsAndBifurcationsFromMinimumSpanningTree( )
88 {
89 }
90
91 // -------------------------------------------------------------------------
92 template< class _TMST >
93 void
94 fpa::Base::ExtractEndPointsAndBifurcationsFromMinimumSpanningTree< _TMST >::
95 GenerateData( )
96 {
97   typedef std::multimap< double, TVertex > _TQueue;
98   typedef std::pair< TVertex, TVertex >    _TBranch;
99   typedef std::vector< _TBranch >          _TBranches;
100   typedef typename _TQueue::value_type     _TQueueValue;
101
102   // 0. Useful objects and values
103   auto mst = this->GetMinimumSpanningTree( );
104   auto endpoints = this->GetEndPoints( );
105   auto bifurcations = this->GetBifurcations( );
106   auto collisions = this->GetCollisions( );
107   endpoints->Get( ).clear( );
108   bifurcations->Get( ).clear( );
109   collisions->Get( ).clear( );
110   this->m_Branches.clear( );
111
112   // 1. Get priority queue
113   auto& q = mst->GetNodeQueue( );
114
115   // 2. Main loop: iterate in inverse order (max-priority)
116   unsigned long label = 0;
117   for( auto qIt = q.rbegin( ); qIt != q.rend( ); ++qIt )
118   {
119     auto gCost = qIt->first;
120     auto vertex = qIt->second;
121
122     // 2.1. Check if the vertex has already been visited
123     if( this->_Mark( vertex ) > 0 || this->_SkeletonMark( vertex ) > 0 )
124       continue;
125
126     // 2.2. Get path from front seed
127     auto path = mst->GetPath( vertex );
128
129     // 2.3. Prepare new branch data and prepare new end-point
130     label++;
131     this->m_Branches.push_back( _TBranch( vertex, vertex ) );
132     endpoints->Get( ).push_back( vertex );
133
134     // 2.4. Backtracking
135     auto pIt = path.begin( );
136     TVertex last_collision;
137     bool inCollision = false;
138     while( pIt != path.end( ) )
139     {
140       auto mark = this->_SkeletonMark( *pIt );
141       if( mark > 0 && mark != label )
142       {
143         // 2.4.1. A bifurcation has been reached
144         bifurcations->Get( ).push_back( *pIt );
145
146         // Reorder labels
147         auto coll_branch = this->m_Branches[ mark - 1 ];
148         this->m_Branches[ mark - 1 ] = _TBranch( coll_branch.first, *pIt );
149         this->m_Branches[ label - 1 ] = _TBranch( qIt->second, *pIt );
150         this->m_Branches.push_back( _TBranch( *pIt, coll_branch.second ) );
151
152         // Mark skeleton (b,coll_branch_second) with the new label
153         label++;
154         while( *pIt != coll_branch.second && pIt != path.end( ) )
155         {
156           this->_MarkSkeleton( *pIt, label );
157           pIt++;
158
159         } // elihw
160         if( pIt != path.end( ) )
161           this->_MarkSkeleton( *pIt, label );
162
163         // Force inner loop termination
164         pIt = path.end( );
165       }
166       else
167       {
168         if( !inCollision )
169         {
170           mark = this->_Mark( *pIt );
171           if( mark > 0 && mark != label )
172           {
173             // 2.4.2. A collision has been reached
174             last_collision = *pIt;
175             collisions->Get( ).push_back( *pIt );
176             inCollision = true;
177           }
178           else
179           {
180             // 2.4.3. Just mark sphere and skeleton
181             double r = this->_Radius( *pIt );
182             this->_MarkSkeleton( *pIt, label );
183             this->_MarkSphere( *pIt, r, label );
184             this->m_Branches[ label - 1 ].second = *pIt;
185
186             // 2.4.4. Is this a seed? -> add it to endpoints
187             auto sIt = pIt;
188             sIt++;
189             if( sIt == path.end( ) )
190               endpoints->Get( ).push_back( *pIt );
191
192           } // fi
193
194         } // fi
195         pIt++;
196
197       } // fi
198
199     } // elihw
200
201   } // rof
202
203   // Remove illegal branches
204   auto dIt = this->m_Branches.begin( );
205   while( dIt != this->m_Branches.end( ) )
206   {
207     if( dIt->first == dIt->second )
208       dIt = this->m_Branches.erase( dIt );
209     else
210       dIt++;
211     
212   } // fi
213
214   // Finish
215   endpoints->SetReferenceImage( mst );
216   bifurcations->SetReferenceImage( mst );
217   collisions->SetReferenceImage( mst );
218 }
219
220 #endif // __FPA__BASE__EXTRACTENDPOINTSANDBIFURCATIONSFROMMINIMUMSPANNINGTREE__HXX__
221
222 // eof - $RCSfile$