]> Creatis software - FrontAlgorithms.git/blob - lib/fpa/Base/ExtractEndPointsAndBifurcationsFromMinimumSpanningTree.hxx
3e8b6e591dc33fb916a1e22f849ddbeeb8af2bbf
[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 typename
71 fpa::Base::ExtractEndPointsAndBifurcationsFromMinimumSpanningTree< _TMST >::
72 TBranches*
73 fpa::Base::ExtractEndPointsAndBifurcationsFromMinimumSpanningTree< _TMST >::
74 GetBranches( )
75 {
76   return(
77     dynamic_cast< TBranches* >( this->itk::ProcessObject::GetOutput( 3 ) )
78     );
79 }
80
81 // -------------------------------------------------------------------------
82 template< class _TMST >
83 fpa::Base::ExtractEndPointsAndBifurcationsFromMinimumSpanningTree< _TMST >::
84 ExtractEndPointsAndBifurcationsFromMinimumSpanningTree( )
85   : Superclass( )
86 {
87   this->SetNumberOfRequiredInputs( 1 );
88   this->SetNumberOfRequiredOutputs( 4 );
89   typename TVertices::Pointer ep = TVertices::New( );
90   typename TVertices::Pointer bf = TVertices::New( );
91   typename TVertices::Pointer co = TVertices::New( );
92   typename TBranches::Pointer br = TBranches::New( );
93   this->itk::ProcessObject::SetNthOutput( 0, ep.GetPointer( ) );
94   this->itk::ProcessObject::SetNthOutput( 1, bf.GetPointer( ) );
95   this->itk::ProcessObject::SetNthOutput( 2, co.GetPointer( ) );
96   this->itk::ProcessObject::SetNthOutput( 3, br.GetPointer( ) );
97 }
98
99 // -------------------------------------------------------------------------
100 template< class _TMST >
101 fpa::Base::ExtractEndPointsAndBifurcationsFromMinimumSpanningTree< _TMST >::
102 ~ExtractEndPointsAndBifurcationsFromMinimumSpanningTree( )
103 {
104 }
105
106 // -------------------------------------------------------------------------
107 template< class _TMST >
108 void
109 fpa::Base::ExtractEndPointsAndBifurcationsFromMinimumSpanningTree< _TMST >::
110 GenerateData( )
111 {
112   typedef std::multimap< double, TVertex > _TQueue;
113   typedef std::pair< TVertex, TVertex >    _TBranch;
114   typedef std::vector< _TBranch >          _TBranches;
115   typedef typename _TQueue::value_type     _TQueueValue;
116
117   // 0. Useful objects and values
118   auto mst = this->GetMinimumSpanningTree( );
119   auto endpoints = this->GetEndPoints( );
120   auto bifurcations = this->GetBifurcations( );
121   auto collisions = this->GetCollisions( );
122   auto branches = this->GetBranches( );
123   endpoints->Get( ).clear( );
124   bifurcations->Get( ).clear( );
125   collisions->Get( ).clear( );
126   branches->Get( ).clear( );
127
128   // 1. Get priority queue
129   auto& q = mst->GetNodeQueue( );
130
131   // 2. Main loop: iterate in inverse order (max-priority)
132   unsigned long label = 0;
133   for( auto qIt = q.rbegin( ); qIt != q.rend( ); ++qIt )
134   {
135     auto gCost = qIt->first;
136     auto vertex = qIt->second;
137
138     // 2.1. Check if the vertex has already been visited
139     if( this->_Mark( vertex ) > 0 )
140       continue;
141
142     // 2.2. Get path from front seed
143     auto path = mst->GetPath( vertex );
144
145     // 2.3. Prepare new branch data and prepare new end-point
146     label++;
147     branches->Get( ).push_back( _TBranch( vertex, vertex ) );
148     endpoints->Get( ).push_back( vertex );
149
150     // 2.4. Backtracking
151     auto pIt = path.begin( );
152     TVertex last_collision;
153     bool inCollision = false;
154     while( pIt != path.end( ) )
155     {
156       auto mark = this->_SkeletonMark( *pIt );
157       if( mark > 0 && mark != label )
158       {
159         // 2.4.1. A bifurcation has been reached
160         bifurcations->Get( ).push_back( *pIt );
161
162         // Reorder labels
163         auto coll_branch = branches->Get( )[ mark ];
164         branches->Get( )[ mark  ] = _TBranch( coll_branch.first, *pIt );
165         branches->Get( )[ label - 1 ] = _TBranch( qIt->second, *pIt );
166         branches->Get( ).push_back( _TBranch( *pIt, coll_branch.second ) );
167
168         // Mark skeleton (b,coll_branch_second) with the new label
169         label++;
170         while( *pIt != coll_branch.second && pIt != path.end( ) )
171         {
172           this->_MarkSkeleton( *pIt, label );
173           pIt++;
174
175         } // elihw
176
177         // Force inner loop termination
178         pIt = path.end( );
179       }
180       else
181       {
182         if( !inCollision )
183         {
184           mark = this->_Mark( *pIt );
185           if( mark > 0 && mark != label )
186           {
187             // 2.4.2. A collision has been reached
188             last_collision = *pIt;
189             collisions->Get( ).push_back( *pIt );
190             inCollision = true;
191           }
192           else
193           {
194             // 2.4.3. Just mark sphere and skeleton
195             double r = this->_Radius( *pIt );
196             this->_MarkSkeleton( *pIt, label );
197             this->_MarkSphere( *pIt, r, label );
198             branches->Get( )[ label - 1 ].second = *pIt;
199
200             // 2.4.4. Is this a seed? -> add it to endpoints
201             auto sIt = pIt;
202             sIt++;
203             if( sIt == path.end( ) )
204               endpoints->Get( ).push_back( *pIt );
205
206           } // fi
207
208         } // fi
209         pIt++;
210
211       } // fi
212
213     } // elihw
214
215   } // rof
216
217   endpoints->SetReferenceImage( mst );
218   bifurcations->SetReferenceImage( mst );
219   collisions->SetReferenceImage( mst );
220 }
221
222 #endif // __FPA__BASE__EXTRACTENDPOINTSANDBIFURCATIONSFROMMINIMUMSPANNINGTREE__HXX__
223
224 // eof - $RCSfile$