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