]> Creatis software - clitk.git/blob - itk/clitkResampleImageWithOptionsFilter.h
Ensure compatibility with VTK6 for Image2DicomRTStruct tool
[clitk.git] / itk / clitkResampleImageWithOptionsFilter.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 #ifndef CLITKRESAMPLEIMAGEWITHOPTIONSFILTER_H
20 #define CLITKRESAMPLEIMAGEWITHOPTIONSFILTER_H
21
22 #include "itkImageToImageFilter.h"
23 #include "itkAffineTransform.h"
24
25 namespace clitk {
26   
27   //--------------------------------------------------------------------
28   /*
29     Image resampling with several interpolations and Gaussian filtering included.
30   */
31   //-------------------------------------------------------------------- 
32   template <class InputImageType, class OutputImageType=InputImageType>
33     class ITK_EXPORT ResampleImageWithOptionsFilter: 
34   public itk::ImageToImageFilter<InputImageType, OutputImageType> {
35
36   public:
37   /** Standard class typedefs. */
38   typedef ResampleImageWithOptionsFilter                     Self;
39   typedef itk::ImageToImageFilter<InputImageType, OutputImageType> Superclass;
40   typedef itk::SmartPointer<Self>                            Pointer;
41   typedef itk::SmartPointer<const Self>                      ConstPointer;
42     
43   /** Method for creation through the object factory. */
44   itkNewMacro(Self);
45     
46   /** Run-time type information (and related methods). */
47   itkTypeMacro(ResampleImageWithOptionsFilter, ImageToImageFilter);
48
49   /** Some convenient typedefs. */
50   typedef typename InputImageType::ConstPointer InputImageConstPointer;
51   typedef typename InputImageType::Pointer      InputImagePointer;
52   typedef typename InputImageType::RegionType   InputImageRegionType; 
53   typedef typename InputImageType::PixelType    InputImagePixelType;
54   typedef typename InputImageType::SpacingType  InputImageSpacingType;
55   typedef typename InputImageType::SizeType     InputImageSizeType;
56     
57   typedef typename OutputImageType::ConstPointer OutputImageConstPointer;
58   typedef typename OutputImageType::Pointer      OutputImagePointer;
59   typedef typename OutputImageType::RegionType   OutputImageRegionType; 
60   typedef typename OutputImageType::PixelType    OutputImagePixelType;
61   typedef typename OutputImageType::SpacingType  OutputImageSpacingType;
62   typedef typename OutputImageType::SizeType     OutputImageSizeType;
63   typedef typename OutputImageType::PointType    OutputImageOriginType;
64   typedef typename OutputImageType::DirectionType        OutputImageDirectionType;
65     
66   typedef itk::AffineTransform<double, InputImageType::ImageDimension> TransformType;
67   typedef typename InputImageType::SpacingType                         GaussianSigmaType;
68
69   /** Interpolation types */
70   typedef enum {
71     NearestNeighbor = 0,
72     Linear = 1, 
73     BSpline = 2, 
74     B_LUT = 3,
75     WSINC = 4
76   } InterpolationTypeEnumeration;
77
78   /** ImageDimension constants */
79   itkStaticConstMacro(InputImageDimension, unsigned int,
80                       InputImageType::ImageDimension);
81   itkStaticConstMacro(OutputImageDimension, unsigned int,
82                       OutputImageType::ImageDimension);
83   itkConceptMacro(SameDimensionCheck,
84                   (itk::Concept::SameDimension<InputImageDimension, OutputImageDimension>));
85
86   /** Input : image to resample */
87   void SetInput(const InputImageType * image);
88     
89   /** ImageDimension constants */
90   itkStaticConstMacro(ImageDimension, unsigned int, InputImageType::ImageDimension);
91
92   // Options
93   itkGetMacro(LastDimensionIsTime, bool);
94   itkSetMacro(LastDimensionIsTime, bool);
95   itkSetMacro(OutputIsoSpacing, double);
96   itkGetMacro(OutputIsoSpacing, double);
97   itkSetMacro(OutputSpacing, OutputImageSpacingType);
98   itkGetMacro(OutputSpacing, OutputImageSpacingType);
99   itkSetMacro(OutputSize, OutputImageSizeType);
100   itkGetMacro(OutputSize, OutputImageSizeType);
101   itkSetMacro(OutputOrigin, OutputImageOriginType);
102   itkGetMacro(OutputOrigin, OutputImageOriginType);
103   itkSetMacro(OutputDirection, OutputImageDirectionType);
104   itkGetMacro(OutputDirection, OutputImageDirectionType);
105   itkGetMacro(InterpolationType, InterpolationTypeEnumeration);
106   itkSetMacro(InterpolationType, InterpolationTypeEnumeration);    
107   itkGetMacro(GaussianFilteringEnabled, bool);
108   itkSetMacro(GaussianFilteringEnabled, bool);
109   itkGetMacro(BSplineOrder, int);
110   itkSetMacro(BSplineOrder, int);
111   itkGetMacro(BLUTSamplingFactor, int);
112   itkSetMacro(BLUTSamplingFactor, int);
113   itkGetMacro(Transform, typename TransformType::Pointer);
114   itkSetMacro(Transform, typename TransformType::Pointer);
115   itkGetMacro(GaussianSigma, GaussianSigmaType);
116   itkSetMacro(GaussianSigma, GaussianSigmaType);
117   itkGetMacro(DefaultPixelValue, OutputImagePixelType);
118   itkSetMacro(DefaultPixelValue, OutputImagePixelType);
119   itkGetMacro(VerboseOptions, bool);
120   itkSetMacro(VerboseOptions, bool);
121     
122   protected:
123   ResampleImageWithOptionsFilter();
124   virtual ~ResampleImageWithOptionsFilter() {}
125     
126   bool m_LastDimensionIsTime;
127   double m_OutputIsoSpacing;
128   InterpolationTypeEnumeration m_InterpolationType;
129   bool m_GaussianFilteringEnabled;
130   int m_BSplineOrder;
131   int m_BLUTSamplingFactor;    
132   OutputImageSizeType m_OutputSize;
133   OutputImageSpacingType m_OutputSpacing;
134   OutputImageOriginType    m_OutputOrigin;
135   OutputImageDirectionType m_OutputDirection;
136   typename TransformType::Pointer m_Transform;
137   GaussianSigmaType m_GaussianSigma;
138   OutputImagePixelType m_DefaultPixelValue;
139   bool m_VerboseOptions;
140   OutputImageRegionType m_OutputRegion;
141
142   virtual void GenerateInputRequestedRegion();
143   virtual void GenerateOutputInformation();
144   virtual void GenerateData();
145     
146   private:
147   ResampleImageWithOptionsFilter(const Self&); //purposely not implemented
148   void operator=(const Self&); //purposely not implemented
149     
150   }; // end class
151   //--------------------------------------------------------------------
152
153   // Convenient function 
154   template<class InputImageType>
155     typename InputImageType::Pointer ResampleImageSpacing(typename InputImageType::Pointer input, 
156                                                           typename InputImageType::SpacingType spacing, 
157                                                           int interpolationType=0);
158
159 } // end namespace clitk
160 //--------------------------------------------------------------------
161
162 #ifndef ITK_MANUAL_INSTANTIATION
163 #include "clitkResampleImageWithOptionsFilter.txx"
164 #endif
165
166 #endif