]> Creatis software - clitk.git/blob - filters/clitkImageResampleGenericFilter.cxx
added the new headers
[clitk.git] / filters / clitkImageResampleGenericFilter.cxx
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://oncora1.lyon.fnclcc.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 #ifndef CLITKIMAGERESAMPLEGENERICFILTER2_CXX
19 #define CLITKIMAGERESAMPLEGENERICFILTER2_CXX
20 /**
21  -------------------------------------------------------------------
22  * @file   clitkImageResampleGenericFilter.cxx
23  * @author David Sarrut <David.Sarrut@creatis.insa-lyon.fr>
24  * @date   23 Feb 2008 08:37:53
25
26  * @brief  
27
28  -------------------------------------------------------------------*/
29
30 #include "clitkImageResampleGenericFilter.h"
31
32 // itk include
33 #include "itkImage.h"
34 #include "itkImageFileReader.h"
35 #include "itkImageSeriesReader.h"
36 #include "itkImageFileWriter.h"
37 #include "itkRecursiveGaussianImageFilter.h"
38 #include "itkResampleImageFilter.h"
39 #include "itkAffineTransform.h"
40 #include "itkNearestNeighborInterpolateImageFunction.h"
41 #include "itkLinearInterpolateImageFunction.h"
42 #include "itkBSplineInterpolateImageFunction.h"
43 #include "itkBSplineInterpolateImageFunctionWithLUT.h"
44 #include "itkCommand.h"
45
46 //--------------------------------------------------------------------
47 clitk::ImageResampleGenericFilter::ImageResampleGenericFilter():
48   ImageToImageGenericFilter<Self>("ImageResample") {
49   mApplyGaussianFilterBefore = false;
50   mDefaultPixelValue = 0.0;
51   mInterpolatorName = "NN";
52   mBSplineOrder=3;
53   InitializeImageTypeWithDim<2>();
54   InitializeImageTypeWithDim<3>();
55   InitializeImageTypeWithDim<4>();
56 }
57 //--------------------------------------------------------------------
58
59
60 //--------------------------------------------------------------------
61 template<unsigned int Dim>
62 void clitk::ImageResampleGenericFilter::InitializeImageTypeWithDim() {      
63   ADD_IMAGE_TYPE(Dim, char);
64   ADD_IMAGE_TYPE(Dim, uchar);
65   ADD_IMAGE_TYPE(Dim, short);
66   ADD_IMAGE_TYPE(Dim, int);
67   ADD_IMAGE_TYPE(Dim, float);
68 }
69 //--------------------------------------------------------------------
70
71
72 //--------------------------------------------------------------------
73 template<class InputImageType>
74 void clitk::ImageResampleGenericFilter::UpdateWithInputImageType() {
75
76   // Some typedefs
77   typedef typename InputImageType::SizeType    SizeType;
78   typedef typename InputImageType::SpacingType SpacingType;
79   typedef typename InputImageType::PointType   PointType;
80   typedef typename InputImageType::PixelType   PixelType;
81   static unsigned int dim = InputImageType::ImageDimension;
82
83   // Reading input
84   typename InputImageType::Pointer input = this->GetInput<InputImageType>(0);
85
86   // Warning
87   if (!std::numeric_limits<PixelType>::is_signed) {
88     if ((mInterpolatorName == "bspline") || (mInterpolatorName == "blut")) {
89       std::cerr << "Warning : input pixel type is not signed, use bspline interpolation at your own risk ..." << std::endl;
90     }
91   }
92
93   // Check options
94   if (mOutputSize.size() != dim) {
95     std::cerr << "Please set size with " << dim << " dimensions." << std::endl;
96     return;
97   }
98   if (mOutputSpacing.size() != dim) {
99     std::cerr << "Please set spacing with " << dim << " dimensions." << std::endl;
100     return;
101   }
102   mOutputOrigin.resize(dim);
103
104   if (mApplyGaussianFilterBefore && mSigma.size() != dim) {
105     std::cerr << "Please set sigma with " << dim << " dimensions." << std::endl;
106     return;
107   }
108
109   // Create Image Filter
110   typedef itk::ResampleImageFilter<InputImageType,InputImageType> FilterType;
111   typename FilterType::Pointer filter = FilterType::New();
112     
113   // Instance of the transform object to be passed to the resample
114   // filter. By default, identity transform is applied
115   typedef itk::AffineTransform<double, InputImageType::ImageDimension> TransformType;
116   typename TransformType::Pointer transform =  TransformType::New();
117   filter->SetTransform(transform);
118
119   // Set filter's parameters
120   SizeType outputSize;
121   SpacingType outputSpacing;
122   PointType outputOrigin;
123   for(unsigned int i=0; i<InputImageType::ImageDimension; i++) {
124     outputSize[i] = mOutputSize[i];
125     outputSpacing[i] = mOutputSpacing[i];
126     outputOrigin[i] = input->GetOrigin()[i];
127   }
128
129   filter->SetSize(outputSize);
130   filter->SetOutputSpacing(outputSpacing);
131   filter->SetOutputOrigin(outputOrigin);
132   filter->SetDefaultPixelValue(static_cast<PixelType>(mDefaultPixelValue));//DS TODO//JV comme ça?
133
134   // Select interpolator
135   if (mInterpolatorName == "nn") {
136     typedef itk::NearestNeighborInterpolateImageFunction<InputImageType, double> InterpolatorType;     
137     typename InterpolatorType::Pointer interpolator = InterpolatorType::New();
138     filter->SetInterpolator(interpolator);
139   }
140   else { 
141     if (mInterpolatorName == "linear") {
142       typedef itk::LinearInterpolateImageFunction<InputImageType, double> InterpolatorType;     
143       typename InterpolatorType::Pointer interpolator =  InterpolatorType::New();
144       filter->SetInterpolator(interpolator);
145     }
146     else {
147       if (mInterpolatorName == "bspline") {
148         typedef itk::BSplineInterpolateImageFunction<InputImageType, double> InterpolatorType; 
149         typename InterpolatorType::Pointer interpolator = InterpolatorType::New();
150         interpolator->SetSplineOrder(mBSplineOrder);
151         filter->SetInterpolator(interpolator);
152       }
153       else {
154         if (mInterpolatorName == "blut") {
155           typedef itk::BSplineInterpolateImageFunctionWithLUT<InputImageType, double> InterpolatorType; 
156           typename InterpolatorType::Pointer interpolator = InterpolatorType::New();
157           interpolator->SetSplineOrder(mBSplineOrder);
158           interpolator->SetLUTSamplingFactor(mSamplingFactors[0]);
159           filter->SetInterpolator(interpolator);
160         }
161         else {
162           std::cerr << "Sorry, I do not know the interpolator '" << mInterpolatorName 
163                     << "'. Known interpolators are :  nn, linear, bspline, blut" << std::endl;
164           exit(0);
165         }
166       }
167     }
168   }
169
170   // Build initial Gaussian bluring (if needed)
171   typedef itk::RecursiveGaussianImageFilter<InputImageType, InputImageType> GaussianFilterType;
172   std::vector<typename GaussianFilterType::Pointer> gaussianFilters;
173   if (mApplyGaussianFilterBefore) {
174     for(unsigned int i=0; i<InputImageType::ImageDimension; i++) {
175       // Create filter
176       gaussianFilters.push_back(GaussianFilterType::New());
177       // Set options
178       gaussianFilters[i]->SetDirection(i);
179       gaussianFilters[i]->SetOrder(GaussianFilterType::ZeroOrder);
180       gaussianFilters[i]->SetNormalizeAcrossScale(false);
181       gaussianFilters[i]->SetSigma(mSigma[i]); // in millimeter !
182       // Set input
183       if (i==0) gaussianFilters[i]->SetInput(input);
184       else gaussianFilters[i]->SetInput(gaussianFilters[i-1]->GetOutput());
185     }
186     filter->SetInput(gaussianFilters[InputImageType::ImageDimension-1]->GetOutput());
187   }
188   else {
189     filter->SetInput(input);
190   }
191
192   // Go !
193   try { 
194     filter->Update();
195   }
196   catch( itk::ExceptionObject & err ) {
197     std::cerr << "Error while filtering " << mInputFilenames[0].c_str() 
198               << " " << err << std::endl;
199     exit(0);
200   }
201
202   // Get result
203    typename InputImageType::Pointer outputImage = filter->GetOutput();
204
205   // Write/save results
206   this->SetNextOutput<InputImageType>(outputImage);
207
208 }
209 //--------------------------------------------------------------------
210
211
212 //--------------------------------------------------------------------
213 void clitk::ImageResampleGenericFilter::SetOutputSize(const std::vector<int> & size) {
214   mOutputSize.resize(size.size());
215   std::copy(size.begin(), size.end(), mOutputSize.begin());
216 }
217 //--------------------------------------------------------------------
218
219 //--------------------------------------------------------------------
220 void clitk::ImageResampleGenericFilter::SetOutputSpacing(const std::vector<double> & spacing) {
221   mOutputSpacing.resize(spacing.size());
222   std::copy(spacing.begin(), spacing.end(), mOutputSpacing.begin());
223 }
224 //--------------------------------------------------------------------
225
226 //--------------------------------------------------------------------
227 void clitk::ImageResampleGenericFilter::SetInterpolationName(const std::string & inter) {
228   mInterpolatorName = inter;
229 }
230 //--------------------------------------------------------------------
231
232 //--------------------------------------------------------------------
233 void clitk::ImageResampleGenericFilter::SetGaussianSigma(const std::vector<double> & sigma) {
234   mApplyGaussianFilterBefore = true;
235   mSigma.resize(sigma.size());
236   std::copy(sigma.begin(), sigma.end(), mSigma.begin());
237 }
238 //--------------------------------------------------------------------
239
240 #endif
241