]> Creatis software - clitk.git/blob - itk/RelativePositionPropImageFilter.h
4b30db1a4fcc739f1796047fd12bd9fc0b3f13a7
[clitk.git] / itk / RelativePositionPropImageFilter.h
1 /*=========================================================================
2   Program:   vv                     http://www.creatis.insa-lyon.fr/rio/vv
3
4   Authors belong to:
5   - University of LYON              http://www.universite-lyon.fr/
6   - Léon Bérard cancer center       http://www.centreleonberard.fr
7   - CREATIS CNRS laboratory         http://www.creatis.insa-lyon.fr
8
9   This software is distributed WITHOUT ANY WARRANTY; without even
10   the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
11   PURPOSE.  See the copyright notices for more information.
12
13   It is distributed under dual licence
14
15   - BSD        See included LICENSE.txt file
16   - CeCILL-B   http://www.cecill.info/licences/Licence_CeCILL-B_V1-en.html
17 ===========================================================================*/
18 /*=========================================================================
19
20   Program:   Insight Segmentation & Registration Toolkit
21   Module:    $RCSfile: RelativePositionPropImageFilter.h,v $
22   Language:  C++
23   Date:      $Date: 2010/07/12 06:57:25 $
24   Version:   $Revision: 1.2 $
25
26   Copyright (c) Insight Software Consortium. All rights reserved.
27   See ITKCopyright.txt or http://www.itk.org/HTML/Copyright.htm for details.
28
29   This software is distributed WITHOUT ANY WARRANTY; without even 
30   the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR 
31   PURPOSE.  See the above copyright notices for more information.
32
33   =========================================================================*/
34
35 #ifndef __RelativePositionPropImageFilter_h
36 #define __RelativePositionPropImageFilter_h
37
38 #include "itkImageToImageFilter.h"
39 #include "itkImage.h"
40 #include "itkConceptChecking.h"
41 #include "itkPointSet.h"
42 #include "itkImageRegionConstIteratorWithIndex.h"
43 #include "itkImageRegionIteratorWithIndex.h"
44 #include "itkMinimumImageFilter.h"
45
46 namespace itk
47 {
48
49
50   /** \class RelativePositionPropImageFilter
51    * \brief Compute the fuzzy subset of an image which satisfies some directional relative position.
52    * \author Jamal Atif and Olivier Nempont
53    * 
54    * This filter computes a fuzzy subset of an image satisfying a particular directionnal relative position from an object (crisp or fuzzy).
55    *
56    * 
57    *    \par INPUT / OUTPUT
58    *    This filter takes a crisp or a fuzzy object as input. 
59    *  In fuzzy case, the values have to be defined between 0 and 1. 
60    *  
61    *  The result is a fuzzy subset which values are defined between
62    *  0 if the relation isn't fulfilled in this point to 1 is the relation is 
63    *  fully satisfied.
64    *  WARNING: the output image type as to be decimal.
65    *  
66    *    \par PARAMETERS
67    *  \par
68    *  The Alpha1 and Alpha2 parameters are used to specify the direction.
69    *  Alpha1 is the angle in 'xy' plane from 'x' unit vector.
70    *  Alpha2 is used in 3D to specify the angle with 'xy' plane
71    *  
72    *  \par
73    *    K is an opening parameter. Higher value enlarge the support of the result.
74    *  By default it is fixed at PI/2  
75    * 
76    *  \par REFERENCE
77    *   Fuzzy Relative Position Between Objects in Image Processing: A Morphological Approach
78    *     Isabelle Bloch
79    *   IEEE TRANSACTIONS ON PATTERN ANALYSIS AND MACHINE INTELLIGENCE, VOL. 21, NO. 7, JULY 1999
80    *   
81    *   This filter is implemented using the propagation algorithm
82    */
83
84   template <class TInputImage, class TOutputImage, class TtNorm=Functor::Minimum<
85                                                      typename TOutputImage::PixelType,
86                                                      typename TOutputImage::PixelType,
87                                                      typename TOutputImage::PixelType>  >
88   class ITK_EXPORT RelativePositionPropImageFilter :
89     public ImageToImageFilter< TInputImage, TOutputImage > 
90   {
91   public:
92     /** Standard class typedefs. */
93     typedef RelativePositionPropImageFilter Self;
94     typedef ImageToImageFilter< TInputImage, TOutputImage > Superclass;
95     typedef SmartPointer<Self> Pointer;
96     typedef SmartPointer<const Self>  ConstPointer;
97
98     /** Extract some information from the image types.  Dimensionality
99      * of the two images is assumed to be the same. */
100     typedef typename TOutputImage::PixelType OutputPixelType;
101     typedef typename TOutputImage::InternalPixelType OutputInternalPixelType;
102     typedef typename TInputImage::PixelType InputPixelType;
103     typedef typename TInputImage::InternalPixelType InputInternalPixelType;
104
105     /** Extract some information from the image types.  Dimensionality
106      * of the two images is assumed to be the same. */
107     itkStaticConstMacro(ImageDimension, unsigned int,
108                         TOutputImage::ImageDimension);
109   
110
111     typedef typename itk::Image<typename TInputImage::IndexType , ImageDimension>
112       CorrespondanceMapType;
113     typedef float TabulationPixelType;
114     typedef typename itk::Image<TabulationPixelType , ImageDimension> TabulationImageType;
115   
116   
117     /** Image typedef support. */
118     typedef TInputImage  InputImageType;
119     typedef TOutputImage OutputImageType;
120   
121     typedef TtNorm TNormType;
122   
123     typedef itk::Vector<double, ImageDimension> VectorType;
124
125
126     /** Method for creation through the object factory. */
127     itkNewMacro(Self);
128
129     /** Run-time type information (and related methods). */
130     itkTypeMacro(RelativePositionPropImageFilter, ImageToImageFilter);
131
132     /** The output pixel type must be signed. */
133     itkConceptMacro(SignedOutputPixelType, (Concept::Signed<OutputPixelType>));
134   
135     /** Standard get/set macros for filter parameters. */
136
137   
138     itkSetMacro(Alpha1, double);
139     itkGetMacro(Alpha1, double);
140     itkSetMacro(Alpha2, double);
141     itkGetMacro(Alpha2, double);
142   
143     itkSetMacro(K1, double);
144     itkGetMacro(K1, double);
145     //   itkSetMacro(K2, double);
146     //   itkGetMacro(K2, double);
147
148     itkSetMacro(Radius, double);
149     itkGetMacro(Radius, double);
150
151     itkSetMacro(VerboseProgress, bool);
152     itkGetMacro(VerboseProgress, bool);
153     itkBooleanMacro(VerboseProgress);
154
155     itkSetMacro(Fast, bool);
156     itkGetMacro(Fast, bool);
157     itkBooleanMacro(Fast);
158   
159     void computeDirection()
160     {
161       if (ImageDimension == 2) {
162           m_DirectionVector[0]=cos(m_Alpha1);
163           m_DirectionVector[1]=sin(m_Alpha1);           
164         }
165       else { // 3D
166         m_DirectionVector[0]=cos(m_Alpha1)*cos(m_Alpha2);
167         m_DirectionVector[1]=cos(m_Alpha2)*sin(m_Alpha1);
168         m_DirectionVector[2]=sin(m_Alpha2);
169       }
170     }
171   
172
173     virtual void GenerateInputRequestedRegion() ITK_OVERRIDE;
174     void EnlargeOutputRequestedRegion (DataObject * output) ITK_OVERRIDE;
175
176   protected:
177     RelativePositionPropImageFilter()
178       {
179         m_Alpha1 = 0;
180         m_Alpha2 = 0;
181         m_K1 = vcl_acos(-1.0)/2;
182         // m_K2 = 3.1417/2;
183         m_Radius = 2; // DS
184         m_Fast = true; // DS
185         m_VerboseProgress = false;
186       }
187     virtual ~RelativePositionPropImageFilter() {}
188     void PrintSelf(std::ostream& os, Indent indent) const ITK_OVERRIDE;
189
190     //void GenerateThreadedData(const typename TOutputImage::RegionType& outputRegionForThread, int threadId);
191     void GenerateData() ITK_OVERRIDE;
192
193   private:
194     RelativePositionPropImageFilter(const Self&); //purposely not implemented
195     void operator=(const Self&); //purposely not implemented
196
197
198     /** The angles*/
199     double m_Alpha1; 
200     double m_Alpha2; 
201     double m_K1;
202     // double m_K2;
203   
204     unsigned int m_Radius;
205     TNormType m_TNorm;
206     bool m_VerboseProgress;
207      
208     VectorType m_DirectionVector; 
209
210     /**
211      * 2 pass instead of 2^NDimension. Warning this may cause some artifacts 
212      */
213     bool m_Fast;
214
215     //allocation et initialisation de la carte de correspondance
216     typename CorrespondanceMapType::Pointer InitCorrespondanceMap();
217   
218     //compute the tabulation map
219     typename TabulationImageType::Pointer ComputeAngleTabulation();
220   
221   
222   };
223   
224 } // end namespace itk
225
226 #ifndef ITK_MANUAL_INSTANTIATION
227 #include "RelativePositionPropImageFilter.txx"
228 #endif
229
230 #endif