]> Creatis software - FrontAlgorithms.git/blob - lib/fpa/Common/RandomWalker.h
...
[FrontAlgorithms.git] / lib / fpa / Common / RandomWalker.h
1 // =========================================================================
2 // @author Leonardo Florez Valencia
3 // @email florez-l@javeriana.edu.co
4 // =========================================================================
5 #ifndef __fpa__Common__RandomWalker__h__
6 #define __fpa__Common__RandomWalker__h__
7
8 #include <fpa/Config.h>
9 #include <map>
10 #include <vector>
11 #include <itkImage.h>
12 #include <itkImageToImageFilter.h>
13 #include <itkSimpleFastMutexLock.h>
14 #include <fpa/Functors/VertexFunction.h>
15
16 namespace fpa
17 {
18   namespace Common
19   {
20     /**
21      */
22     template< class _TImage, class _TLabels, class _TScalar = float >
23     class RandomWalker
24       : public itk::ImageToImageFilter< _TImage, _TLabels >
25     {
26     public:
27       typedef _TImage  TImage;
28       typedef _TLabels TLabels;
29       typedef _TScalar TScalar;
30
31       typedef RandomWalker                       Self;
32       typedef itk::ImageToImageFilter< TImage, TLabels > Superclass;
33       typedef itk::SmartPointer< Self >                  Pointer;
34       typedef itk::SmartPointer< const Self >            ConstPointer;
35
36       typedef itk::Image< TScalar, TImage::ImageDimension > TScalarImage;
37       typedef typename TImage::IndexType  TIndex;
38       typedef typename TImage::PixelType  TPixel;
39       typedef typename TImage::RegionType TRegion;
40       typedef typename TLabels::PixelType TLabel;
41
42       typedef fpa::Functors::VertexFunction< TIndex, TScalar > TEdgeFunction;
43
44     protected:
45       struct _TBoundaryThreadStruct
46       {
47         Pointer Filter;
48         void* Triplets;
49         std::map< TLabel, unsigned long >* Labels;
50       };
51
52       struct _TLaplacianThreadStruct
53       {
54         Pointer Filter;
55         void* A;
56         void* R;
57         const void* B;
58       };
59
60       struct _TOutputThreadStruct
61       {
62         Pointer Filter;
63         const void* X;
64         const void* S;
65         const std::vector< TLabel >* InvLabels;
66       };
67
68     public:
69       itkNewMacro( Self );
70       itkTypeMacro(
71         fpa::Common::RandomWalker, itk::ImageToImageFilter
72         );
73
74       itkGetConstObjectMacro( EdgeFunction, TEdgeFunction );
75       itkGetObjectMacro( EdgeFunction, TEdgeFunction );
76       itkSetObjectMacro( EdgeFunction, TEdgeFunction );
77
78       fpaFilterInputMacro( InputLabels, TLabels );
79       fpaFilterOutputMacro( OutputProbabilities, TScalarImage );
80
81     protected:
82       RandomWalker( );
83       virtual ~RandomWalker( );
84
85       virtual void GenerateData( ) override;
86
87       virtual TScalar _L( const TIndex& i, const TIndex& j );
88
89       // Boundary construction
90       template< class _TTriplets >
91       void _Boundary(
92         _TTriplets& B,
93         std::map< TLabel, unsigned long >& labels
94         );
95
96       template< class _TTriplets >
97       static ITK_THREAD_RETURN_TYPE _BoundaryCbk( void* arg );
98
99       template< class _TTriplets >
100       void _ThreadedBoundary(
101         const TRegion& region, const itk::ThreadIdType& id,
102         _TTriplets* B,
103         std::map< TLabel, unsigned long >* labels
104         );
105
106       // Laplacian construction
107       template< class _TTriplets >
108       void _Laplacian(
109         _TTriplets& A, _TTriplets& R, const _TTriplets& B
110         );
111
112       template< class _TTriplets >
113       static ITK_THREAD_RETURN_TYPE _LaplacianCbk( void* arg );
114
115       template< class _TTriplets >
116       void _ThreadedLaplacian(
117         const TRegion& region, const itk::ThreadIdType& id,
118         _TTriplets* A, _TTriplets* R, const _TTriplets* B
119         );
120
121       // Output filling
122       template< class _TMatrix, class _TTriplets >
123       void _Output(
124         const _TMatrix& X, const _TTriplets& S,
125         const std::vector< TLabel >& invLabels
126         );
127
128       template< class _TMatrix, class _TTriplets >
129       static ITK_THREAD_RETURN_TYPE _OutputCbk( void* arg );
130
131       template< class _TMatrix, class _TTriplets >
132       void _ThreadedOutput(
133         const TRegion& region, const itk::ThreadIdType& id,
134         const _TMatrix* X, const _TTriplets* S,
135         const std::vector< TLabel >* invLabels
136         );
137
138       // Array-Image conversion
139       static unsigned long _1D( const TIndex& idx, const TRegion& region );
140
141       template< class _TTriplets >
142       static unsigned long _SeedIndex(
143         const unsigned long& i, const _TTriplets& t
144         );
145
146       template< class _TTriplets >
147       static unsigned long _NearSeedIndex(
148         const unsigned long& i, const _TTriplets& t
149         );
150
151     private:
152       // Purposely not implemented
153       RandomWalker( const Self& other );
154       Self& operator=( const Self& other );
155
156     protected:
157       typename TEdgeFunction::Pointer m_EdgeFunction;
158       itk::SimpleFastMutexLock m_Mutex;
159     };
160
161   } // ecapseman
162
163 } // ecapseman
164
165 #ifndef ITK_MANUAL_INSTANTIATION
166 #  include <fpa/Common/RandomWalker.hxx>
167 #endif // ITK_MANUAL_INSTANTIATION
168 #endif // __fpa__Common__RandomWalker__h__
169 // eof - $RCSfile$