ITK  4.3.0
Insight Segmentation and Registration Toolkit
itkResampleImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * http://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 #ifndef __itkResampleImageFilter_h
19 #define __itkResampleImageFilter_h
20 
21 #include "itkFixedArray.h"
22 #include "itkTransform.h"
23 #include "itkImageRegionIterator.h"
24 #include "itkImageToImageFilter.h"
27 #include "itkSize.h"
29 
30 namespace itk
31 {
83 template< class TInputImage,
84  class TOutputImage,
85  class TInterpolatorPrecisionType = double >
86 class ITK_EXPORT ResampleImageFilter:
87  public ImageToImageFilter< TInputImage, TOutputImage >
88 {
89 public:
95 
96  typedef TInputImage InputImageType;
97  typedef TOutputImage OutputImageType;
98  typedef typename InputImageType::Pointer InputImagePointer;
99  typedef typename InputImageType::ConstPointer InputImageConstPointer;
100  typedef typename OutputImageType::Pointer OutputImagePointer;
101  typedef typename InputImageType::RegionType InputImageRegionType;
102 
104  itkNewMacro(Self);
105 
108 
110  itkStaticConstMacro(ImageDimension, unsigned int,
111  TOutputImage::ImageDimension);
112  itkStaticConstMacro(InputImageDimension, unsigned int,
113  TInputImage::ImageDimension);
115 
118 
122  typedef Transform< TInterpolatorPrecisionType,
123  itkGetStaticConstMacro(ImageDimension),
124  itkGetStaticConstMacro(ImageDimension) > TransformType;
127 
130  TInterpolatorPrecisionType > InterpolatorType;
132 
134 
136 
138 
140  TInterpolatorPrecisionType > LinearInterpolatorType;
141  typedef typename LinearInterpolatorType::Pointer
143 
146  TInterpolatorPrecisionType > ExtrapolatorType;
148 
151 
153  typedef typename TOutputImage::IndexType IndexType;
154 
157  //typedef typename TOutputImage::PointType PointType;
158 
160  typedef typename TOutputImage::PixelType PixelType;
161  typedef typename TInputImage::PixelType InputPixelType;
162 
164 
166 
170 
172  typedef typename TOutputImage::RegionType OutputImageRegionType;
173 
175  typedef typename TOutputImage::SpacingType SpacingType;
176  typedef typename TOutputImage::PointType OriginPointType;
177  typedef typename TOutputImage::DirectionType DirectionType;
178 
186  itkSetConstObjectMacro(Transform, TransformType);
187 
189  itkGetConstObjectMacro(Transform, TransformType);
190 
198  itkSetObjectMacro(Interpolator, InterpolatorType);
199 
201  itkGetConstObjectMacro(Interpolator, InterpolatorType);
202 
206  itkSetObjectMacro(Extrapolator, ExtrapolatorType);
207 
209  itkGetConstObjectMacro(Extrapolator, ExtrapolatorType);
210 
212  itkSetMacro(Size, SizeType);
213 
215  itkGetConstReferenceMacro(Size, SizeType);
216 
219  itkSetMacro(DefaultPixelValue, PixelType);
220 
222  itkGetConstReferenceMacro(DefaultPixelValue, PixelType);
223 
225  itkSetMacro(OutputSpacing, SpacingType);
226  virtual void SetOutputSpacing(const double *values);
228 
230  itkGetConstReferenceMacro(OutputSpacing, SpacingType);
231 
233  itkSetMacro(OutputOrigin, OriginPointType);
234  virtual void SetOutputOrigin(const double *values);
236 
238  itkGetConstReferenceMacro(OutputOrigin, OriginPointType);
239 
241  itkSetMacro(OutputDirection, DirectionType);
242  itkGetConstReferenceMacro(OutputDirection, DirectionType);
244 
246  void SetOutputParametersFromImage(const ImageBaseType *image);
247 
250  itkSetMacro(OutputStartIndex, IndexType);
251 
253  itkGetConstReferenceMacro(OutputStartIndex, IndexType);
254 
260  void SetReferenceImage(const TOutputImage *image);
261 
262  const TOutputImage * GetReferenceImage() const;
263 
264  itkSetMacro(UseReferenceImage, bool);
265  itkBooleanMacro(UseReferenceImage);
266  itkGetConstMacro(UseReferenceImage, bool);
267 
273  virtual void GenerateOutputInformation();
274 
280  virtual void GenerateInputRequestedRegion();
281 
284  virtual void BeforeThreadedGenerateData();
285 
288  virtual void AfterThreadedGenerateData();
289 
291  ModifiedTimeType GetMTime(void) const;
292 
293 #ifdef ITK_USE_CONCEPT_CHECKING
294 
295  itkConceptMacro( OutputHasNumericTraitsCheck,
297 
299 #endif
300 
301 protected:
304  void PrintSelf(std::ostream & os, Indent indent) const;
305 
306 
312  virtual void VerifyInputInformation() {}
313 
323  virtual void ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread,
324  ThreadIdType threadId);
325 
328  virtual void NonlinearThreadedGenerateData(const OutputImageRegionType &
329  outputRegionForThread,
330  ThreadIdType threadId);
331 
335  virtual void LinearThreadedGenerateData(const OutputImageRegionType &
336  outputRegionForThread,
337  ThreadIdType threadId);
338 
339  virtual PixelType CastPixelWithBoundsChecking( const InterpolatorOutputType value,
340  const ComponentType minComponent,
341  const ComponentType maxComponent) const;
342 
343 private:
344  ResampleImageFilter(const Self &); //purposely not implemented
345  void operator=(const Self &); //purposely not implemented
346 
347  SizeType m_Size; // Size of the output image
350  // interpolation
352  // extrapolation
353  PixelType m_DefaultPixelValue; // default pixel value
354  // if the point is
355  // outside the image
356  SpacingType m_OutputSpacing; // output image spacing
357  OriginPointType m_OutputOrigin; // output image origin
358  DirectionType m_OutputDirection; // output image direction cosines
359  IndexType m_OutputStartIndex; // output image start index
361 
362 };
363 } // end namespace itk
364 
365 #ifndef ITK_MANUAL_INSTANTIATION
366 #include "itkResampleImageFilter.hxx"
367 #endif
368 
369 #endif
370