ITK  4.5.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< typename TInputImage,
84  typename TOutputImage,
85  typename TInterpolatorPrecisionType = double,
86  typename TTransformPrecisionType = TInterpolatorPrecisionType>
88  public ImageToImageFilter< TInputImage, TOutputImage >
89 {
90 public:
96 
97  typedef TInputImage InputImageType;
98  typedef TOutputImage OutputImageType;
99  typedef typename InputImageType::Pointer InputImagePointer;
100  typedef typename InputImageType::ConstPointer InputImageConstPointer;
101  typedef typename OutputImageType::Pointer OutputImagePointer;
102  typedef typename InputImageType::RegionType InputImageRegionType;
103 
105  itkNewMacro(Self);
106 
109 
111  itkStaticConstMacro(ImageDimension, unsigned int,
112  TOutputImage::ImageDimension);
113  itkStaticConstMacro(InputImageDimension, unsigned int,
114  TInputImage::ImageDimension);
116 
119 
123  typedef Transform< TTransformPrecisionType,
124  itkGetStaticConstMacro(ImageDimension),
125  itkGetStaticConstMacro(ImageDimension) > TransformType;
128 
131  TInterpolatorPrecisionType > InterpolatorType;
133 
135 
137 
139 
141  TInterpolatorPrecisionType > LinearInterpolatorType;
142  typedef typename LinearInterpolatorType::Pointer
144 
147  TInterpolatorPrecisionType > ExtrapolatorType;
149 
152 
154  typedef typename TOutputImage::IndexType IndexType;
155 
158  //typedef typename TOutputImage::PointType PointType;
159 
161  typedef typename TOutputImage::PixelType PixelType;
162  typedef typename TInputImage::PixelType InputPixelType;
163 
165 
167 
171 
173  typedef typename TOutputImage::RegionType OutputImageRegionType;
174 
176  typedef typename TOutputImage::SpacingType SpacingType;
177  typedef typename TOutputImage::PointType OriginPointType;
178  typedef typename TOutputImage::DirectionType DirectionType;
179 
182 
190  itkSetConstObjectMacro(Transform, TransformType);
191  itkGetConstObjectMacro(Transform, TransformType);
193 
201  itkSetObjectMacro(Interpolator, InterpolatorType);
202  itkGetModifiableObjectMacro(Interpolator, InterpolatorType);
204 
208  itkSetObjectMacro(Extrapolator, ExtrapolatorType);
209  itkGetModifiableObjectMacro(Extrapolator, ExtrapolatorType);
211 
213  itkSetMacro(Size, SizeType);
214  itkGetConstReferenceMacro(Size, SizeType);
216 
219  itkSetMacro(DefaultPixelValue, PixelType);
220  itkGetConstReferenceMacro(DefaultPixelValue, PixelType);
222 
224  itkSetMacro(OutputSpacing, SpacingType);
225  virtual void SetOutputSpacing(const double *values);
227 
229  itkGetConstReferenceMacro(OutputSpacing, SpacingType);
230 
232  itkSetMacro(OutputOrigin, OriginPointType);
233  virtual void SetOutputOrigin(const double *values);
235 
237  itkGetConstReferenceMacro(OutputOrigin, OriginPointType);
238 
240  itkSetMacro(OutputDirection, DirectionType);
241  itkGetConstReferenceMacro(OutputDirection, DirectionType);
243 
245  void SetOutputParametersFromImage(const ImageBaseType *image);
246 
249  itkSetMacro(OutputStartIndex, IndexType);
250 
252  itkGetConstReferenceMacro(OutputStartIndex, IndexType);
253 
260  void SetReferenceImage(const ReferenceImageBaseType *image);
261 
264 
267  itkSetMacro(UseReferenceImage, bool);
268  itkBooleanMacro(UseReferenceImage);
269  itkGetConstMacro(UseReferenceImage, bool);
271 
277  virtual void GenerateOutputInformation();
278 
284  virtual void GenerateInputRequestedRegion();
285 
288  virtual void BeforeThreadedGenerateData();
289 
292  virtual void AfterThreadedGenerateData();
293 
295  ModifiedTimeType GetMTime(void) const;
296 
297 #ifdef ITK_USE_CONCEPT_CHECKING
298  // Begin concept checking
299  itkConceptMacro( OutputHasNumericTraitsCheck,
301  // End concept checking
302 #endif
303 
304 protected:
307  }
308  void PrintSelf(std::ostream & os, Indent indent) const;
309 
315  virtual void VerifyInputInformation() {
316  }
317 
327  virtual void ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread,
328  ThreadIdType threadId);
329 
333  outputRegionForThread,
334  ThreadIdType threadId);
335 
340  outputRegionForThread,
341  ThreadIdType threadId);
342 
344  const ComponentType minComponent,
345  const ComponentType maxComponent) const;
346 
347 private:
348  ResampleImageFilter(const Self &); //purposely not implemented
349  void operator=(const Self &); //purposely not implemented
350 
351  SizeType m_Size; // Size of the output image
354  // interpolation
356  // extrapolation
357  PixelType m_DefaultPixelValue; // default pixel value
358  // if the point is
359  // outside the image
360  SpacingType m_OutputSpacing; // output image spacing
361  OriginPointType m_OutputOrigin; // output image origin
362  DirectionType m_OutputDirection; // output image direction cosines
363  IndexType m_OutputStartIndex; // output image start index
365 
366 };
367 } // end namespace itk
368 
369 #ifndef ITK_MANUAL_INSTANTIATION
370 #include "itkResampleImageFilter.hxx"
371 #endif
372 
373 #endif
ImageBase< ImageDimension > ReferenceImageBaseType
virtual void SetOutputSpacing(SpacingType _arg)
TOutputImage::DirectionType DirectionType
TOutputImage::IndexType IndexType
ImageToImageFilter< TInputImage, TOutputImage > Superclass
InterpolatorPointerType m_Interpolator
virtual void GenerateOutputInformation()
virtual void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId)
DefaultConvertPixelTraits< PixelType > PixelConvertType
ExtrapolateImageFunction< InputImageType, TInterpolatorPrecisionType > ExtrapolatorType
const ReferenceImageBaseType * GetReferenceImage() const
virtual void GenerateInputRequestedRegion()
InterpolatorConvertType::ComponentType ComponentType
TransformPointerType m_Transform
Size< itkGetStaticConstMacro(ImageDimension) > SizeType
Resample an image via a coordinate transform.
LinearInterpolateImageFunction< InputImageType, TInterpolatorPrecisionType > LinearInterpolatorType
unsigned long ModifiedTimeType
Definition: itkIntTypes.h:164
InputImageType::RegionType InputImageRegionType
Traits class used to by ConvertPixels to convert blocks of pixels.
SmartPointer< Self > Pointer
TOutputImage::RegionType OutputImageRegionType
ContinuousIndex< TTransformPrecisionType, ImageDimension > ContinuousInputIndexType
Base class for all process objects that output image data.
ExtrapolatorType::Pointer ExtrapolatorPointerType
void operator=(const Self &)
InterpolatorType::Pointer InterpolatorPointerType
Transform< TTransformPrecisionType, itkGetStaticConstMacro(ImageDimension), itkGetStaticConstMacro(ImageDimension) > TransformType
void SetOutputParametersFromImage(const ImageBaseType *image)
virtual PixelType CastPixelWithBoundsChecking(const InterpolatorOutputType value, const ComponentType minComponent, const ComponentType maxComponent) const
ImageBase< itkGetStaticConstMacro(ImageDimension) > ImageBaseType
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:82
PixelConvertType::ComponentType PixelComponentType
DefaultConvertPixelTraits< InterpolatorOutputType > InterpolatorConvertType
LinearInterpolatorType::Pointer LinearInterpolatorPointerType
static const unsigned int ImageDimension
static const unsigned int InputImageDimension
virtual void SetOutputOrigin(OriginPointType _arg)
TransformType::ConstPointer TransformPointerType
void PrintSelf(std::ostream &os, Indent indent) const
TOutputImage::SpacingType SpacingType
virtual void BeforeThreadedGenerateData()
Linearly interpolate an image at specified positions.
void SetReferenceImage(const ReferenceImageBaseType *image)
SmartPointer< const Self > ConstPointer
Base class for all image interpolaters.
virtual void LinearThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId)
TOutputImage::PixelType PixelType
InputImageType::Pointer InputImagePointer
Base class for templated image classes.
Definition: itkImageBase.h:112
virtual void NonlinearThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId)
A templated class holding a point in n-Dimensional image space.
Base class for all image extrapolaters.
InputImageType::ConstPointer InputImageConstPointer
InterpolatorType::PointType PointType
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
ExtrapolatorPointerType m_Extrapolator
ModifiedTimeType GetMTime(void) const
OutputImageType::Pointer OutputImagePointer
#define itkConceptMacro(name, concept)
A templated class holding a geometric point in n-Dimensional space.
Definition: itkPoint.h:51
TOutputImage::PointType OriginPointType
InterpolateImageFunction< InputImageType, TInterpolatorPrecisionType > InterpolatorType
virtual void AfterThreadedGenerateData()
TInputImage::PixelType InputPixelType
InterpolatorType::OutputType InterpolatorOutputType
unsigned int ThreadIdType
Definition: itkIntTypes.h:159