ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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);
94
100
101 using InputImageType = TInputImage;
102 using OutputImageType = TOutputImage;
103 using InputImagePointer = typename InputImageType::Pointer;
104 using InputImageConstPointer = typename InputImageType::ConstPointer;
105 using OutputImagePointer = typename OutputImageType::Pointer;
106 using InputImageRegionType = typename InputImageType::RegionType;
107
109 itkNewMacro(Self);
110
112 itkOverrideGetNameOfClassMacro(ResampleImageFilter);
113
115 static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
116
118 static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
119
120#if !defined(ITK_LEGACY_REMOVE)
121 static constexpr unsigned int ImageDimension = OutputImageDimension; // For backward compatibility
122#endif
123
126
134
135
139
141
143
145
148
152
155
157 using IndexType = typename TOutputImage::IndexType;
158
161 using OutputPointType = typename TOutputImage::PointType;
162
164 using PixelType = typename TOutputImage::PixelType;
165 using InputPixelType = typename TInputImage::PixelType;
166
168
170
173
175 using OutputImageRegionType = typename TOutputImage::RegionType;
176
178 using SpacingType = typename TOutputImage::SpacingType;
179 using OriginPointType = typename TOutputImage::PointType;
180 using DirectionType = typename TOutputImage::DirectionType;
181
184
185 /* See superclass for doxygen. This method adds the additional check
186 * that the output space is set */
187 void
188 VerifyPreconditions() const override;
189
197 itkSetGetDecoratedObjectInputMacro(Transform, TransformType);
198
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);
235 itkGetConstReferenceMacro(OutputSpacing, SpacingType);
236
239 itkSetMacro(OutputOrigin, OriginPointType);
240 virtual void
241 SetOutputOrigin(const double * origin);
244 itkGetConstReferenceMacro(OutputOrigin, OriginPointType);
245
248 itkSetMacro(OutputDirection, DirectionType);
249 itkGetConstReferenceMacro(OutputDirection, DirectionType);
252 void
254
257 itkSetMacro(OutputStartIndex, IndexType);
258
260 itkGetConstReferenceMacro(OutputStartIndex, IndexType);
261
270
273
277 itkSetMacro(UseReferenceImage, bool);
278 itkBooleanMacro(UseReferenceImage);
279 itkGetConstMacro(UseReferenceImage, bool);
282
283protected:
285 ~ResampleImageFilter() override = default;
286 void
287 PrintSelf(std::ostream & os, Indent indent) const override;
288
294 void
295 VerifyInputInformation() const override
296 {}
297
303 void
305
311 void
313
317 void
319
321 void
323
326 GetMTime() const override;
327
337 void
338 DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
339
340
343 virtual void
345
348 virtual void
349 LinearThreadedGenerateData(const OutputImageRegionType & outputRegionForThread);
350
351#if !defined(ITK_LEGACY_REMOVE)
353 itkLegacyMacro(virtual PixelType CastPixelWithBoundsChecking(const InterpolatorOutputType value,
354 const ComponentType minComponent,
355 const ComponentType maxComponent) const;)
356#endif
357
358private:
359 static PixelComponentType
361
362 template <typename TComponent>
363 static PixelComponentType
364 CastComponentWithBoundsChecking(const TComponent value);
365
366 static PixelType
368
369 template <typename TPixel>
370 static PixelType
371 CastPixelWithBoundsChecking(const TPixel value);
372
373 void
375
376 SizeType m_Size{}; // Size of the output image
378 // interpolation
380 // extrapolation
381 PixelType m_DefaultPixelValue{}; // default pixel value
382 // if the point is
383 // outside the image
384 SpacingType m_OutputSpacing{}; // output image spacing
385 OriginPointType m_OutputOrigin{}; // output image origin
386 DirectionType m_OutputDirection{}; // output image direction cosines
387 IndexType m_OutputStartIndex{}; // output image start index
388 bool m_UseReferenceImage{ false };
389};
390} // end namespace itk
391
392#ifndef ITK_MANUAL_INSTANTIATION
393# include "itkResampleImageFilter.hxx"
394#endif
395
396#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 InterpolatorOutputType::ComponentType ComponentType
Base class for all image extrapolaters.
Base class for templated image classes.
Control indentation during Print() invocation.
Definition itkIndent.h:50
Base class for all image interpolators.
Point< TInterpolatorPrecisionType, Self::ImageDimension > PointType
Linearly interpolate an image at specified positions.
typename InputImageType::ConstPointer InputImageConstPointer
DefaultConvertPixelTraits< PixelType > PixelConvertType
ModifiedTimeType GetMTime() const override
DefaultConvertPixelTraits< InterpolatorOutputType > InterpolatorConvertType
typename InterpolatorType::OutputType InterpolatorOutputType
typename DecoratedTransformType::Pointer DecoratedTransformPointer
virtual void NonlinearThreadedGenerateData(const OutputImageRegionType &outputRegionForThread)
ImageToImageFilter< TInputImage, TOutputImage > Superclass
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
typename InputImageType::RegionType InputImageRegionType
virtual void LinearThreadedGenerateData(const OutputImageRegionType &outputRegionForThread)
typename ExtrapolatorType::Pointer ExtrapolatorPointerType
static PixelType CastPixelWithBoundsChecking(const ComponentType value)
ImageBase< OutputImageDimension > ReferenceImageBaseType
void GenerateInputRequestedRegion() override
typename TOutputImage::PointType OutputPointType
SmartPointer< const Self > ConstPointer
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
static PixelComponentType CastComponentWithBoundsChecking(const TComponent value)
typename TOutputImage::RegionType OutputImageRegionType
typename TOutputImage::DirectionType DirectionType
typename TOutputImage::PointType OriginPointType
void PrintSelf(std::ostream &os, Indent indent) const override
void GenerateOutputInformation() override
~ResampleImageFilter() override=default
typename InputImageType::Pointer InputImagePointer
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
virtual void SetOutputOrigin(const double *origin)
static PixelComponentType CastComponentWithBoundsChecking(const PixelComponentType value)
DataObjectDecorator< TransformType > DecoratedTransformType
Transform< TTransformPrecisionType, Self::OutputImageDimension, Self::InputImageDimension > TransformType
typename PixelConvertType::ComponentType PixelComponentType
Size< Self::OutputImageDimension > SizeType
typename InterpolatorType::PointType InputPointType
typename OutputImageType::Pointer OutputImagePointer
typename TransformType::ConstPointer TransformPointerType
virtual void SetOutputSpacing(const double *spacing)
LinearInterpolateImageFunction< InputImageType, TInterpolatorPrecisionType > LinearInterpolatorType
ExtrapolateImageFunction< InputImageType, TInterpolatorPrecisionType > ExtrapolatorType
void VerifyInputInformation() const override
typename TOutputImage::IndexType IndexType
ContinuousIndex< TInterpolatorPrecisionType, InputImageDimension > ContinuousInputIndexType
ImageBase< Self::OutputImageDimension > ImageBaseType
typename TInputImage::PixelType InputPixelType
typename InterpolatorType::Pointer InterpolatorPointerType
InterpolateImageFunction< InputImageType, TInterpolatorPrecisionType > InterpolatorType
itkSetInputMacro(ReferenceImage, ReferenceImageBaseType)
Implements transparent reference counting.
Transform points and vectors from an input space to an output space.
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
SizeValueType ModifiedTimeType
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition itkSize.h:70