ITK  6.0.0
Insight Toolkit
itkResampleImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
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  * https://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 {
86 template <typename TInputImage,
87  typename TOutputImage,
88  typename TInterpolatorPrecisionType = double,
89  typename TTransformPrecisionType = TInterpolatorPrecisionType>
90 class ITK_TEMPLATE_EXPORT ResampleImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
91 {
92 public:
93  ITK_DISALLOW_COPY_AND_MOVE(ResampleImageFilter);
101 
102  using InputImageType = TInputImage;
103  using OutputImageType = TOutputImage;
108 
110  itkNewMacro(Self);
111 
113  itkOverrideGetNameOfClassMacro(ResampleImageFilter);
114 
116  static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
117 
119  static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
120 
121 #if !defined(ITK_LEGACY_REMOVE)
122  static constexpr unsigned int ImageDimension = OutputImageDimension; // For backward compatibility
123 #endif
124 
127 
135 
136 
140 
142 
144 
146 
149 
153 
156 
159 
163 
165  using PixelType = typename TOutputImage::PixelType;
166  using InputPixelType = typename TInputImage::PixelType;
167 
169 
171 
174 
177 
179  using SpacingType = typename TOutputImage::SpacingType;
182 
185 
186  /* See superclass for doxygen. This method adds the additional check
187  * that the output space is set */
188  void
189  VerifyPreconditions() const override;
190 
198  itkSetGetDecoratedObjectInputMacro(Transform, TransformType);
199 
207  itkSetObjectMacro(Interpolator, InterpolatorType);
208  itkGetModifiableObjectMacro(Interpolator, InterpolatorType);
214  itkSetObjectMacro(Extrapolator, ExtrapolatorType);
215  itkGetModifiableObjectMacro(Extrapolator, ExtrapolatorType);
219  itkSetMacro(Size, SizeType);
220  itkGetConstReferenceMacro(Size, SizeType);
225  itkSetMacro(DefaultPixelValue, PixelType);
226  itkGetConstReferenceMacro(DefaultPixelValue, PixelType);
230  itkSetMacro(OutputSpacing, SpacingType);
231  virtual void
232  SetOutputSpacing(const double * spacing);
236  itkGetConstReferenceMacro(OutputSpacing, SpacingType);
237 
239  itkSetMacro(OutputOrigin, OriginPointType);
240  virtual void
241  SetOutputOrigin(const double * origin);
245  itkGetConstReferenceMacro(OutputOrigin, OriginPointType);
246 
248  itkSetMacro(OutputDirection, DirectionType);
249  itkGetConstReferenceMacro(OutputDirection, DirectionType);
253  void
254  SetOutputParametersFromImage(const ImageBaseType * image);
255 
258  itkSetMacro(OutputStartIndex, IndexType);
259 
261  itkGetConstReferenceMacro(OutputStartIndex, IndexType);
262 
270  itkSetInputMacro(ReferenceImage, ReferenceImageBaseType);
271 
273  itkGetInputMacro(ReferenceImage, ReferenceImageBaseType);
274 
277  itkSetMacro(UseReferenceImage, bool);
278  itkBooleanMacro(UseReferenceImage);
279  itkGetConstMacro(UseReferenceImage, bool);
282 #ifdef ITK_USE_CONCEPT_CHECKING
283  // Begin concept checking
284  itkConceptMacro(OutputHasNumericTraitsCheck, (Concept::HasNumericTraits<PixelComponentType>));
285  // End concept checking
286 #endif
287 
288 protected:
290  ~ResampleImageFilter() override = default;
291  void
292  PrintSelf(std::ostream & os, Indent indent) const override;
293 
299  void
300  VerifyInputInformation() const override
301  {}
302 
308  void
309  GenerateOutputInformation() override;
310 
316  void
317  GenerateInputRequestedRegion() override;
318 
322  void
323  BeforeThreadedGenerateData() override;
324 
326  void
327  AfterThreadedGenerateData() override;
328 
331  GetMTime() const override;
332 
342  void
343  DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
344 
345 
348  virtual void
349  NonlinearThreadedGenerateData(const OutputImageRegionType & outputRegionForThread);
350 
353  virtual void
354  LinearThreadedGenerateData(const OutputImageRegionType & outputRegionForThread);
355 
357  itkLegacyMacro(virtual PixelType CastPixelWithBoundsChecking(const InterpolatorOutputType value,
358  const ComponentType minComponent,
359  const ComponentType maxComponent) const;)
360 
361 private:
362  static PixelComponentType
363  CastComponentWithBoundsChecking(const PixelComponentType value);
364 
365  template <typename TComponent>
366  static PixelComponentType
367  CastComponentWithBoundsChecking(const TComponent value);
368 
369  static PixelType
370  CastPixelWithBoundsChecking(const ComponentType value);
371 
372  template <typename TPixel>
373  static PixelType
374  CastPixelWithBoundsChecking(const TPixel value);
375 
376  void
377  InitializeTransform();
378 
379  SizeType m_Size{}; // Size of the output image
380  InterpolatorPointerType m_Interpolator{}; // Image function for
381  // interpolation
382  ExtrapolatorPointerType m_Extrapolator{}; // Image function for
383  // extrapolation
384  PixelType m_DefaultPixelValue{}; // default pixel value
385  // if the point is
386  // outside the image
387  SpacingType m_OutputSpacing{}; // output image spacing
388  OriginPointType m_OutputOrigin{}; // output image origin
389  DirectionType m_OutputDirection{}; // output image direction cosines
390  IndexType m_OutputStartIndex{}; // output image start index
391  bool m_UseReferenceImage{ false };
392 };
393 } // end namespace itk
394 
395 #ifndef ITK_MANUAL_INSTANTIATION
396 # include "itkResampleImageFilter.hxx"
397 #endif
398 
399 #endif
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::ExtrapolateImageFunction
Base class for all image extrapolaters.
Definition: itkExtrapolateImageFunction.h:44
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itk::GTest::TypedefsAndConstructors::Dimension2::DirectionType
ImageBaseType::DirectionType DirectionType
Definition: itkGTestTypedefsAndConstructors.h:52
itk::Concept::HasNumericTraits
Definition: itkConceptChecking.h:716
itk::ImageSource::OutputImagePointer
typename OutputImageType::Pointer OutputImagePointer
Definition: itkImageSource.h:91
itk::ModifiedTimeType
SizeValueType ModifiedTimeType
Definition: itkIntTypes.h:105
itk::ResampleImageFilter::InterpolatorOutputType
typename InterpolatorType::OutputType InterpolatorOutputType
Definition: itkResampleImageFilter.h:141
itk::Size< Self::OutputImageDimension >
itkLinearInterpolateImageFunction.h
itk::ResampleImageFilter::TransformPointerType
typename TransformType::ConstPointer TransformPointerType
Definition: itkResampleImageFilter.h:132
itk::ImageBase
Base class for templated image classes.
Definition: itkImageBase.h:114
itk::ResampleImageFilter::PixelType
typename TOutputImage::PixelType PixelType
Definition: itkResampleImageFilter.h:165
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::SmartPointer< Self >
itkImageRegionIterator.h
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::FunctionBase< Point< TCoordRep, TInputImage::ImageDimension >, NumericTraits< TInputImage::PixelType >::RealType >::OutputType
NumericTraits< TInputImage::PixelType >::RealType OutputType
Definition: itkFunctionBase.h:62
itk::DataObjectDecorator
Decorates any subclass of itkObject with a DataObject API.
Definition: itkDataObjectDecorator.h:66
itk::ResampleImageFilter::PixelComponentType
typename PixelConvertType::ComponentType PixelComponentType
Definition: itkResampleImageFilter.h:170
itk::ResampleImageFilter::IndexType
typename TOutputImage::IndexType IndexType
Definition: itkResampleImageFilter.h:158
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itkDefaultConvertPixelTraits.h
itk::LinearInterpolateImageFunction
Linearly interpolate an image at specified positions.
Definition: itkLinearInterpolateImageFunction.h:51
itk::ResampleImageFilter::DecoratedTransformPointer
typename DecoratedTransformType::Pointer DecoratedTransformPointer
Definition: itkResampleImageFilter.h:134
itk::ImageToImageFilter::InputImagePointer
typename InputImageType::Pointer InputImagePointer
Definition: itkImageToImageFilter.h:130
itk::ResampleImageFilter::InputPixelType
typename TInputImage::PixelType InputPixelType
Definition: itkResampleImageFilter.h:166
itk::DefaultConvertPixelTraits
Traits class used to by ConvertPixels to convert blocks of pixels.
Definition: itkDefaultConvertPixelTraits.h:41
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itkDataObjectDecorator.h
itkFixedArray.h
itkExtrapolateImageFunction.h
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itkImageToImageFilter.h
itk::ResampleImageFilter::VerifyInputInformation
void VerifyInputInformation() const override
Definition: itkResampleImageFilter.h:300
itk::ResampleImageFilter::OriginPointType
typename TOutputImage::PointType OriginPointType
Definition: itkResampleImageFilter.h:180
itk::ResampleImageFilter::ExtrapolatorPointerType
typename ExtrapolatorType::Pointer ExtrapolatorPointerType
Definition: itkResampleImageFilter.h:152
itk::ResampleImageFilter::DirectionType
typename TOutputImage::DirectionType DirectionType
Definition: itkResampleImageFilter.h:181
itk::ImageSource::OutputImageRegionType
typename OutputImageType::RegionType OutputImageRegionType
Definition: itkImageSource.h:92
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:65
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::ContinuousIndex
A templated class holding a point in n-Dimensional image space.
Definition: itkContinuousIndex.h:46
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:139
itk::ResampleImageFilter::ComponentType
typename InterpolatorConvertType::ComponentType ComponentType
Definition: itkResampleImageFilter.h:145
itk::ResampleImageFilter
Resample an image via a coordinate transform.
Definition: itkResampleImageFilter.h:90
itk::Transform
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:83
itkTransform.h
itk::ResampleImageFilter::OutputPointType
typename TOutputImage::PointType OutputPointType
Definition: itkResampleImageFilter.h:162
itk::ResampleImageFilter::SpacingType
typename TOutputImage::SpacingType SpacingType
Definition: itkResampleImageFilter.h:179
itk::ResampleImageFilter::InputPointType
typename InterpolatorType::PointType InputPointType
Definition: itkResampleImageFilter.h:161
itk::ImageToImageFilter::InputImageRegionType
typename InputImageType::RegionType InputImageRegionType
Definition: itkImageToImageFilter.h:132
itk::InterpolateImageFunction
Base class for all image interpolators.
Definition: itkInterpolateImageFunction.h:45
itk::DefaultConvertPixelTraits::ComponentType
typename PixelType::ComponentType ComponentType
Definition: itkDefaultConvertPixelTraits.h:45
itk::ImageToImageFilter::InputImageConstPointer
typename InputImageType::ConstPointer InputImageConstPointer
Definition: itkImageToImageFilter.h:131
itk::ResampleImageFilter::InterpolatorPointerType
typename InterpolatorType::Pointer InterpolatorPointerType
Definition: itkResampleImageFilter.h:139
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90
itkSize.h
itk::ResampleImageFilter::LinearInterpolatorPointerType
typename LinearInterpolatorType::Pointer LinearInterpolatorPointerType
Definition: itkResampleImageFilter.h:148