]> 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 #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 OriginalRandomWalker
24       : public itk::ImageToImageFilter< _TImage, _TLabels >
25     {
26     public:
27       typedef _TImage  TImage;
28       typedef _TLabels TLabels;
29       typedef _TScalar TScalar;
30
31       typedef OriginalRandomWalker                       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::OriginalRandomWalker, 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       /* TODO
82          public:
83          void AddSeed( const TIndex& seed, const TLabel& label );
84       */
85
86     protected:
87       OriginalRandomWalker( );
88       virtual ~OriginalRandomWalker( );
89
90       virtual void GenerateData( ) override;
91
92       virtual TScalar _L( const TIndex& i, const TIndex& j );
93
94       // Boundary construction
95       template< class _TTriplets >
96       void _Boundary(
97         _TTriplets& B,
98         std::map< TLabel, unsigned long >& labels
99         );
100
101       template< class _TTriplets >
102       static ITK_THREAD_RETURN_TYPE _BoundaryCbk( void* arg );
103
104       template< class _TTriplets >
105       void _ThreadedBoundary(
106         const TRegion& region, const itk::ThreadIdType& id,
107         _TTriplets* B,
108         std::map< TLabel, unsigned long >* labels
109         );
110
111       // Laplacian construction
112       template< class _TTriplets >
113       void _Laplacian(
114         _TTriplets& A, _TTriplets& R, const _TTriplets& B
115         );
116
117       template< class _TTriplets >
118       static ITK_THREAD_RETURN_TYPE _LaplacianCbk( void* arg );
119
120       template< class _TTriplets >
121       void _ThreadedLaplacian(
122         const TRegion& region, const itk::ThreadIdType& id,
123         _TTriplets* A, _TTriplets* R, const _TTriplets* B
124         );
125
126       // Output filling
127       template< class _TMatrix, class _TTriplets >
128       void _Output(
129         const _TMatrix& X, const _TTriplets& S,
130         const std::vector< TLabel >& invLabels
131         );
132
133       template< class _TMatrix, class _TTriplets >
134       static ITK_THREAD_RETURN_TYPE _OutputCbk( void* arg );
135
136       template< class _TMatrix, class _TTriplets >
137       void _ThreadedOutput(
138         const TRegion& region, const itk::ThreadIdType& id,
139         const _TMatrix* X, const _TTriplets* S,
140         const std::vector< TLabel >* invLabels
141         );
142
143       // Array-Image conversion
144       static unsigned long _1D( const TIndex& idx, const TRegion& region );
145
146       template< class _TTriplets >
147       static unsigned long _SeedIndex(
148         const unsigned long& i, const _TTriplets& t
149         );
150
151       template< class _TTriplets >
152       static unsigned long _NearSeedIndex(
153         const unsigned long& i, const _TTriplets& t
154         );
155
156     private:
157       // Purposely not implemented
158       OriginalRandomWalker( const Self& other );
159       Self& operator=( const Self& other );
160
161     protected:
162       /* TODO
163          std::vector< TIndex > m_Seeds;
164          std::vector< TLabel > m_Labels;
165       */
166
167       typename TEdgeFunction::Pointer m_EdgeFunction;
168
169       itk::SimpleFastMutexLock m_Mutex;
170     };
171
172   } // ecapseman
173
174 } // ecapseman
175
176 #ifndef ITK_MANUAL_INSTANTIATION
177 #  include <fpa/Common/OriginalRandomWalker.hxx>
178 #endif // ITK_MANUAL_INSTANTIATION
179 #endif // __fpa__Common__OriginalRandomWalker__h__
180 // eof - $RCSfile$