ITK  4.6.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 #include "itkDataObjectDecorator.h"
30 
31 
32 namespace itk
33 {
85 template< typename TInputImage,
86  typename TOutputImage,
87  typename TInterpolatorPrecisionType = double,
88  typename TTransformPrecisionType = TInterpolatorPrecisionType>
90  public ImageToImageFilter< TInputImage, TOutputImage >
91 {
92 public:
98 
99  typedef TInputImage InputImageType;
100  typedef TOutputImage OutputImageType;
101  typedef typename InputImageType::Pointer InputImagePointer;
102  typedef typename InputImageType::ConstPointer InputImageConstPointer;
103  typedef typename OutputImageType::Pointer OutputImagePointer;
104  typedef typename InputImageType::RegionType InputImageRegionType;
105 
107  itkNewMacro(Self);
108 
111 
113  itkStaticConstMacro(ImageDimension, unsigned int,
114  TOutputImage::ImageDimension);
115  itkStaticConstMacro(InputImageDimension, unsigned int,
116  TInputImage::ImageDimension);
118 
121 
125  typedef Transform< TTransformPrecisionType,
126  itkGetStaticConstMacro(ImageDimension),
127  itkGetStaticConstMacro(ImageDimension) > TransformType;
132 
133 
136  TInterpolatorPrecisionType > InterpolatorType;
138 
140 
142 
144 
146  TInterpolatorPrecisionType > LinearInterpolatorType;
147  typedef typename LinearInterpolatorType::Pointer
149 
152  TInterpolatorPrecisionType > ExtrapolatorType;
154 
157 
159  typedef typename TOutputImage::IndexType IndexType;
160 
163  //typedef typename TOutputImage::PointType PointType;
164 
166  typedef typename TOutputImage::PixelType PixelType;
167  typedef typename TInputImage::PixelType InputPixelType;
168 
170 
172 
176 
178  typedef typename TOutputImage::RegionType OutputImageRegionType;
179 
181  typedef typename TOutputImage::SpacingType SpacingType;
182  typedef typename TOutputImage::PointType OriginPointType;
183  typedef typename TOutputImage::DirectionType DirectionType;
184 
187 
196 
204  itkSetObjectMacro(Interpolator, InterpolatorType);
205  itkGetModifiableObjectMacro(Interpolator, InterpolatorType);
207 
211  itkSetObjectMacro(Extrapolator, ExtrapolatorType);
212  itkGetModifiableObjectMacro(Extrapolator, ExtrapolatorType);
214 
216  itkSetMacro(Size, SizeType);
217  itkGetConstReferenceMacro(Size, SizeType);
219 
222  itkSetMacro(DefaultPixelValue, PixelType);
223  itkGetConstReferenceMacro(DefaultPixelValue, PixelType);
225 
227  itkSetMacro(OutputSpacing, SpacingType);
228  virtual void SetOutputSpacing(const double *values);
230 
232  itkGetConstReferenceMacro(OutputSpacing, SpacingType);
233 
235  itkSetMacro(OutputOrigin, OriginPointType);
236  virtual void SetOutputOrigin(const double *values);
238 
240  itkGetConstReferenceMacro(OutputOrigin, OriginPointType);
241 
243  itkSetMacro(OutputDirection, DirectionType);
244  itkGetConstReferenceMacro(OutputDirection, DirectionType);
246 
248  void SetOutputParametersFromImage(const ImageBaseType *image);
249 
252  itkSetMacro(OutputStartIndex, IndexType);
253 
255  itkGetConstReferenceMacro(OutputStartIndex, IndexType);
256 
263  itkSetInputMacro(ReferenceImage, ReferenceImageBaseType);
264 
266  itkGetInputMacro(ReferenceImage, ReferenceImageBaseType);
267 
270  itkSetMacro(UseReferenceImage, bool);
271  itkBooleanMacro(UseReferenceImage);
272  itkGetConstMacro(UseReferenceImage, bool);
274 
280  virtual void GenerateOutputInformation();
281 
287  virtual void GenerateInputRequestedRegion();
288 
291  virtual void BeforeThreadedGenerateData();
292 
295  virtual void AfterThreadedGenerateData();
296 
298  ModifiedTimeType GetMTime(void) const;
299 
300 #ifdef ITK_USE_CONCEPT_CHECKING
301  // Begin concept checking
302  itkConceptMacro( OutputHasNumericTraitsCheck,
304  // End concept checking
305 #endif
306 
307 protected:
310  }
311  void PrintSelf(std::ostream & os, Indent indent) const;
312 
318  virtual void VerifyInputInformation() {
319  }
320 
330  virtual void ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread,
331  ThreadIdType threadId);
332 
336  outputRegionForThread,
337  ThreadIdType threadId);
338 
343  outputRegionForThread,
344  ThreadIdType threadId);
345 
347  const ComponentType minComponent,
348  const ComponentType maxComponent) const;
349 
350 private:
351  ResampleImageFilter(const Self &); //purposely not implemented
352  void operator=(const Self &); //purposely not implemented
353 
354  SizeType m_Size; // Size of the output image
356  // interpolation
358  // extrapolation
359  PixelType m_DefaultPixelValue; // default pixel value
360  // if the point is
361  // outside the image
362  SpacingType m_OutputSpacing; // output image spacing
363  OriginPointType m_OutputOrigin; // output image origin
364  DirectionType m_OutputDirection; // output image direction cosines
365  IndexType m_OutputStartIndex; // output image start index
367 
368 };
369 } // end namespace itk
370 
371 #ifndef ITK_MANUAL_INSTANTIATION
372 #include "itkResampleImageFilter.hxx"
373 #endif
374 
375 #endif
ImageBase< ImageDimension > ReferenceImageBaseType
virtual void SetOutputSpacing(SpacingType _arg)
TOutputImage::DirectionType DirectionType
TOutputImage::IndexType IndexType
ImageToImageFilter< TInputImage, TOutputImage > Superclass
InterpolatorPointerType m_Interpolator
DecoratedTransformType::Pointer DecoratedTransformPointer
virtual void GenerateOutputInformation()
virtual void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId)
DefaultConvertPixelTraits< PixelType > PixelConvertType
ExtrapolateImageFunction< InputImageType, TInterpolatorPrecisionType > ExtrapolatorType
virtual void GenerateInputRequestedRegion()
InterpolatorConvertType::ComponentType ComponentType
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
DataObjectDecorator< TransformType > DecoratedTransformType
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
Decorates any subclass of itkObject with a DataObject API.
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.
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
itkSetGetDecoratedObjectInputMacro(Transform, TransformType)
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