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);
95
101
102 using InputImageType = TInputImage;
103 using OutputImageType = TOutputImage;
104 using InputImagePointer = typename InputImageType::Pointer;
105 using InputImageConstPointer = typename InputImageType::ConstPointer;
106 using OutputImagePointer = typename OutputImageType::Pointer;
107 using InputImageRegionType = typename InputImageType::RegionType;
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
158 using IndexType = typename TOutputImage::IndexType;
159
162 using OutputPointType = typename TOutputImage::PointType;
163
165 using PixelType = typename TOutputImage::PixelType;
166 using InputPixelType = typename TInputImage::PixelType;
167
169
171
174
176 using OutputImageRegionType = typename TOutputImage::RegionType;
177
179 using SpacingType = typename TOutputImage::SpacingType;
180 using OriginPointType = typename TOutputImage::PointType;
181 using DirectionType = typename TOutputImage::DirectionType;
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);
210
214 itkSetObjectMacro(Extrapolator, ExtrapolatorType);
215 itkGetModifiableObjectMacro(Extrapolator, ExtrapolatorType);
217
219 itkSetMacro(Size, SizeType);
220 itkGetConstReferenceMacro(Size, SizeType);
222
225 itkSetMacro(DefaultPixelValue, PixelType);
226 itkGetConstReferenceMacro(DefaultPixelValue, PixelType);
228
230 itkSetMacro(OutputSpacing, SpacingType);
231 virtual void
232 SetOutputSpacing(const double * spacing);
234
236 itkGetConstReferenceMacro(OutputSpacing, SpacingType);
237
239 itkSetMacro(OutputOrigin, OriginPointType);
240 virtual void
241 SetOutputOrigin(const double * origin);
243
245 itkGetConstReferenceMacro(OutputOrigin, OriginPointType);
246
248 itkSetMacro(OutputDirection, DirectionType);
249 itkGetConstReferenceMacro(OutputDirection, DirectionType);
251
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);
281
283
284protected:
286 ~ResampleImageFilter() override = default;
287 void
288 PrintSelf(std::ostream & os, Indent indent) const override;
289
295 void
296 VerifyInputInformation() const override
297 {}
298
304 void
306
312 void
314
318 void
320
322 void
324
327 GetMTime() const override;
328
338 void
339 DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
340
341
344 virtual void
346
349 virtual void
350 LinearThreadedGenerateData(const OutputImageRegionType & outputRegionForThread);
351
352#if !defined(ITK_LEGACY_REMOVE)
354 itkLegacyMacro(virtual PixelType CastPixelWithBoundsChecking(const InterpolatorOutputType value,
355 const ComponentType minComponent,
356 const ComponentType maxComponent) const;)
357#endif
358
359private:
360 static PixelComponentType
362
363 template <typename TComponent>
364 static PixelComponentType
365 CastComponentWithBoundsChecking(const TComponent value);
366
367 static PixelType
369
370 template <typename TPixel>
371 static PixelType
372 CastPixelWithBoundsChecking(const TPixel value);
373
374 void
376
377 SizeType m_Size{}; // Size of the output image
379 // interpolation
381 // extrapolation
382 PixelType m_DefaultPixelValue{}; // default pixel value
383 // if the point is
384 // outside the image
385 SpacingType m_OutputSpacing{}; // output image spacing
386 OriginPointType m_OutputOrigin{}; // output image origin
387 DirectionType m_OutputDirection{}; // output image direction cosines
388 IndexType m_OutputStartIndex{}; // output image start index
389 bool m_UseReferenceImage{ false };
390};
391} // end namespace itk
392
393#ifndef ITK_MANUAL_INSTANTIATION
394# include "itkResampleImageFilter.hxx"
395#endif
396
397#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