]> Creatis software - clitk.git/blob - itk/RelativePositionPropImageFilter.h
Fix deprecated warning (Begin -> GoToBegin)
[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 #if ITK_VERSION_MAJOR == 4
85   template <class TInputImage, class TOutputImage, class TtNorm=Functor::Minimum<
86                                                      typename TOutputImage::PixelType,
87                                                      typename TOutputImage::PixelType,
88                                                      typename TOutputImage::PixelType>  >
89 #else
90   template <class TInputImage, class TOutputImage, class TtNorm=Function::Minimum<
91                                                      typename TOutputImage::PixelType,
92                                                      typename TOutputImage::PixelType,
93                                                      typename TOutputImage::PixelType>  >
94 #endif
95   class ITK_EXPORT RelativePositionPropImageFilter :
96     public ImageToImageFilter< TInputImage, TOutputImage > 
97   {
98   public:
99     /** Standard class typedefs. */
100     typedef RelativePositionPropImageFilter Self;
101     typedef ImageToImageFilter< TInputImage, TOutputImage > Superclass;
102     typedef SmartPointer<Self> Pointer;
103     typedef SmartPointer<const Self>  ConstPointer;
104
105     /** Extract some information from the image types.  Dimensionality
106      * of the two images is assumed to be the same. */
107     typedef typename TOutputImage::PixelType OutputPixelType;
108     typedef typename TOutputImage::InternalPixelType OutputInternalPixelType;
109     typedef typename TInputImage::PixelType InputPixelType;
110     typedef typename TInputImage::InternalPixelType InputInternalPixelType;
111
112     /** Extract some information from the image types.  Dimensionality
113      * of the two images is assumed to be the same. */
114     itkStaticConstMacro(ImageDimension, unsigned int,
115                         TOutputImage::ImageDimension);
116   
117
118     typedef typename itk::Image<typename TInputImage::IndexType , ImageDimension>
119       CorrespondanceMapType;
120     typedef float TabulationPixelType;
121     typedef typename itk::Image<TabulationPixelType , ImageDimension> TabulationImageType;
122   
123   
124     /** Image typedef support. */
125     typedef TInputImage  InputImageType;
126     typedef TOutputImage OutputImageType;
127   
128     typedef TtNorm TNormType;
129   
130     typedef itk::Vector<double, ImageDimension> VectorType;
131
132
133     /** Method for creation through the object factory. */
134     itkNewMacro(Self);
135
136     /** Run-time type information (and related methods). */
137     itkTypeMacro(RelativePositionPropImageFilter, ImageToImageFilter);
138
139     /** The output pixel type must be signed. */
140     itkConceptMacro(SignedOutputPixelType, (Concept::Signed<OutputPixelType>));
141   
142     /** Standard get/set macros for filter parameters. */
143
144   
145     itkSetMacro(Alpha1, double);
146     itkGetMacro(Alpha1, double);
147     itkSetMacro(Alpha2, double);
148     itkGetMacro(Alpha2, double);
149   
150     itkSetMacro(K1, double);
151     itkGetMacro(K1, double);
152     //   itkSetMacro(K2, double);
153     //   itkGetMacro(K2, double);
154
155     itkSetMacro(Radius, double);
156     itkGetMacro(Radius, double);
157
158     itkSetMacro(VerboseProgress, bool);
159     itkGetMacro(VerboseProgress, bool);
160     itkBooleanMacro(VerboseProgress);
161
162     itkSetMacro(Fast, bool);
163     itkGetMacro(Fast, bool);
164     itkBooleanMacro(Fast);
165   
166     void computeDirection()
167     {
168       if (ImageDimension == 2) {
169           m_DirectionVector[0]=cos(m_Alpha1);
170           m_DirectionVector[1]=sin(m_Alpha1);           
171         }
172       else { // 3D
173         m_DirectionVector[0]=cos(m_Alpha1)*cos(m_Alpha2);
174         m_DirectionVector[1]=cos(m_Alpha2)*sin(m_Alpha1);
175         m_DirectionVector[2]=sin(m_Alpha2);
176       }
177     }
178   
179
180     virtual void GenerateInputRequestedRegion() throw(InvalidRequestedRegionError);
181     void EnlargeOutputRequestedRegion (DataObject * output);
182
183   protected:
184     RelativePositionPropImageFilter()
185       {
186         m_Alpha1 = 0;
187         m_Alpha2 = 0;
188         m_K1 = vcl_acos(-1.0)/2;
189         // m_K2 = 3.1417/2;
190         m_Radius = 2; // DS
191         m_Fast = true; // DS
192         m_VerboseProgress = false;
193       }
194     virtual ~RelativePositionPropImageFilter() {}
195     void PrintSelf(std::ostream& os, Indent indent) const;
196
197     //void GenerateThreadedData(const typename TOutputImage::RegionType& outputRegionForThread, int threadId);
198     void GenerateData();
199
200   private:
201     RelativePositionPropImageFilter(const Self&); //purposely not implemented
202     void operator=(const Self&); //purposely not implemented
203
204
205     /** The angles*/
206     double m_Alpha1; 
207     double m_Alpha2; 
208     double m_K1;
209     // double m_K2;
210   
211     unsigned int m_Radius;
212     TNormType m_TNorm;
213     bool m_VerboseProgress;
214      
215     VectorType m_DirectionVector; 
216
217     /**
218      * 2 pass instead of 2^NDimension. Warning this may cause some artifacts 
219      */
220     bool m_Fast;
221
222     //allocation et initialisation de la carte de correspondance
223     typename CorrespondanceMapType::Pointer InitCorrespondanceMap();
224   
225     //compute the tabulation map
226     typename TabulationImageType::Pointer ComputeAngleTabulation();
227   
228   
229   };
230   
231 } // end namespace itk
232
233 #ifndef ITK_MANUAL_INSTANTIATION
234 #include "RelativePositionPropImageFilter.txx"
235 #endif
236
237 #endif