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"
27#include "itkSize.h"
30
31
32namespace itk
33{
86template <typename TInputImage,
87 typename TOutputImage,
88 typename TInterpolatorPrecisionType = double,
89 typename TTransformPrecisionType = TInterpolatorPrecisionType>
90class ITK_TEMPLATE_EXPORT ResampleImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
91{
92public:
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
255
258 itkSetMacro(OutputStartIndex, IndexType);
259
261 itkGetConstReferenceMacro(OutputStartIndex, IndexType);
262
271
274
277 itkSetMacro(UseReferenceImage, bool);
278 itkBooleanMacro(UseReferenceImage);
279 itkGetConstMacro(UseReferenceImage, bool);
282#ifdef ITK_USE_CONCEPT_CHECKING
283 // Begin concept checking
285 // End concept checking
286#endif
287
288protected:
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
310
316 void
318
322 void
324
326 void
328
331 GetMTime() const override;
332
342 void
343 DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
344
345
348 virtual void
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
361private:
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
371
372 template <typename TPixel>
373 static PixelType
374 CastPixelWithBoundsChecking(const TPixel value);
375
376 void
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
A templated class holding a point in n-Dimensional image space.
Decorates any subclass of itkObject with a DataObject API.
Traits class used to by ConvertPixels to convert blocks of pixels.
typename PixelType::ComponentType ComponentType
Base class for all image extrapolaters.
Base class for templated image classes.
Definition: itkImageBase.h:115
Base class for all process objects that output image data.
typename OutputImageType::RegionType OutputImageRegionType
TOutputImage OutputImageType
typename OutputImageType::Pointer OutputImagePointer
Base class for filters that take an image as input and produce an image as output.
typename InputImageType::ConstPointer InputImageConstPointer
typename InputImageType::Pointer InputImagePointer
typename InputImageType::RegionType InputImageRegionType
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Base class for all image interpolators.
Linearly interpolate an image at specified positions.
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Resample an image via a coordinate transform.
ModifiedTimeType GetMTime() const override
typename InterpolatorType::OutputType InterpolatorOutputType
typename DecoratedTransformType::Pointer DecoratedTransformPointer
virtual void NonlinearThreadedGenerateData(const OutputImageRegionType &outputRegionForThread)
typename LinearInterpolatorType::Pointer LinearInterpolatorPointerType
typename InterpolatorConvertType::ComponentType ComponentType
typename TOutputImage::PixelType PixelType
void SetOutputParametersFromImage(const ImageBaseType *image)
static PixelType CastPixelWithBoundsChecking(const TPixel value)
void AfterThreadedGenerateData() override
virtual void LinearThreadedGenerateData(const OutputImageRegionType &outputRegionForThread)
typename ExtrapolatorType::Pointer ExtrapolatorPointerType
static PixelType CastPixelWithBoundsChecking(const ComponentType value)
void GenerateInputRequestedRegion() override
typename TOutputImage::PointType OutputPointType
void VerifyPreconditions() const override
Verifies that the process object has been configured correctly, that all required inputs are set,...
itkGetInputMacro(ReferenceImage, ReferenceImageBaseType)
void BeforeThreadedGenerateData() override
typename TOutputImage::SpacingType SpacingType
typename TOutputImage::DirectionType DirectionType
typename TOutputImage::PointType OriginPointType
void PrintSelf(std::ostream &os, Indent indent) const override
void GenerateOutputInformation() override
~ResampleImageFilter() override=default
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
virtual void SetOutputOrigin(const double *origin)
typename PixelConvertType::ComponentType PixelComponentType
typename InterpolatorType::PointType InputPointType
typename TransformType::ConstPointer TransformPointerType
virtual void SetOutputSpacing(const double *spacing)
void VerifyInputInformation() const override
static itkLegacyMacro(virtual PixelType CastPixelWithBoundsChecking(const InterpolatorOutputType value, const ComponentType minComponent, const ComponentType maxComponent) const ;) private PixelComponentType CastComponentWithBoundsChecking(const TComponent value)
typename TOutputImage::IndexType IndexType
typename TInputImage::PixelType InputPixelType
typename InterpolatorType::Pointer InterpolatorPointerType
itkSetInputMacro(ReferenceImage, ReferenceImageBaseType)
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:84
SmartPointer< const Self > ConstPointer
SmartPointer< Self > Pointer
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
SizeValueType ModifiedTimeType
Definition: itkIntTypes.h:105