ITK  6.0.0
Insight Toolkit
itkWarpImageFilter.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 itkWarpImageFilter_h
19#define itkWarpImageFilter_h
20#include "itkImageBase.h"
23
24namespace itk
25{
84template <typename TInputImage, typename TOutputImage, typename TDisplacementField>
85class ITK_TEMPLATE_EXPORT WarpImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
86{
87public:
88 ITK_DISALLOW_COPY_AND_MOVE(WarpImageFilter);
96
98 itkNewMacro(Self);
99
101 itkOverrideGetNameOfClassMacro(WarpImageFilter);
102
105
107 using typename Superclass::InputImageType;
108 using typename Superclass::InputImagePointer;
109 using typename Superclass::OutputImageType;
110 using typename Superclass::OutputImagePointer;
111 using typename Superclass::InputImageConstPointer;
115 using PixelType = typename OutputImageType::PixelType;
116 using PixelComponentType = typename OutputImageType::InternalPixelType;
117 using SpacingType = typename OutputImageType::SpacingType;
118
120 static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
121 static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
122 static constexpr unsigned int DisplacementFieldDimension = TDisplacementField::ImageDimension;
123
126
128 using DisplacementFieldType = TDisplacementField;
130 using DisplacementType = typename DisplacementFieldType::PixelType;
131
133 using CoordinateType = double;
134#ifndef ITK_FUTURE_LEGACY_REMOVE
135 using CoordRepType ITK_FUTURE_DEPRECATED(
136 "ITK 6 discourages using `CoordRepType`. Please use `CoordinateType` instead!") = CoordinateType;
137#endif
141
144
147
148
151
154
156 itkSetObjectMacro(Interpolator, InterpolatorType);
157 itkGetModifiableObjectMacro(Interpolator, InterpolatorType);
161 itkSetMacro(OutputSpacing, SpacingType);
162 virtual void
163 SetOutputSpacing(const double * spacing);
167 itkGetConstReferenceMacro(OutputSpacing, SpacingType);
168
170 itkSetMacro(OutputOrigin, PointType);
171 virtual void
172 SetOutputOrigin(const double * origin);
176 itkGetConstReferenceMacro(OutputOrigin, PointType);
177
179 itkSetMacro(OutputDirection, DirectionType);
180 itkGetConstReferenceMacro(OutputDirection, DirectionType);
184 void
186
189 itkSetMacro(OutputStartIndex, IndexType);
190
192 itkGetConstReferenceMacro(OutputStartIndex, IndexType);
193
195 itkSetMacro(OutputSize, SizeType);
196
198 itkGetConstReferenceMacro(OutputSize, SizeType);
199
201 itkSetMacro(EdgePaddingValue, PixelType);
202
204 itkGetConstMacro(EdgePaddingValue, PixelType);
205
211 void
213
220 void
222
225 void
227
230 void
232
233#ifdef ITK_USE_CONCEPT_CHECKING
234 // Begin concept checking
238 itkConceptMacro(DisplacementFieldHasNumericTraitsCheck,
240 // End concept checking
241#endif
242
243protected:
245 // ~WarpImageFilter() {} default implementation is ok
246
247 void
248 PrintSelf(std::ostream & os, Indent indent) const override;
249
253 void
254 DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
255
256
263 void
264 VerifyInputInformation() const override;
265
274 void
276
280 void
282 const DisplacementFieldType * fieldPtr,
283 DisplacementType & output);
284
285 bool m_DefFieldSameInformation{};
286 // variables for deffield interpolator
287 IndexType m_StartIndex, m_EndIndex{};
288
289private:
290 PixelType m_EdgePaddingValue{};
291 SpacingType m_OutputSpacing{};
292 PointType m_OutputOrigin{};
293 DirectionType m_OutputDirection{};
294
295 InterpolatorPointer m_Interpolator{};
296 SizeType m_OutputSize{}; // Size of the output image
297 IndexType m_OutputStartIndex{}; // output image start index
298};
299} // end namespace itk
300
301#ifndef ITK_MANUAL_INSTANTIATION
302# include "itkWarpImageFilter.hxx"
303#endif
304
305#endif
Base class for templated image classes.
Definition: itkImageBase.h:115
Base class for all process objects that output image data.
typename OutputImageType::RegionType OutputImageRegionType
Base class for filters that take an image as input and produce an image as output.
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...
Warps an image using an input displacement field.
void GenerateInputRequestedRegion() override
typename OutputImageType::SizeType SizeType
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
void GenerateOutputInformation() override
void BeforeThreadedGenerateData() override
TDisplacementField DisplacementFieldType
typename InterpolatorType::Pointer InterpolatorPointer
itkGetInputMacro(DisplacementField, DisplacementFieldType)
typename OutputImageType::SpacingType SpacingType
void EvaluateDisplacementAtPhysicalPoint(const PointType &point, DisplacementType &output)
void VerifyInputInformation() const override
typename OutputImageType::IndexType IndexType
virtual void SetOutputSpacing(const double *spacing)
typename TOutputImage::DirectionType DirectionType
typename DisplacementFieldType::PixelType DisplacementType
virtual void SetOutputOrigin(const double *origin)
void EvaluateDisplacementAtPhysicalPoint(const PointType &point, const DisplacementFieldType *fieldPtr, DisplacementType &output)
itkSetInputMacro(DisplacementField, DisplacementFieldType)
typename OutputImageType::InternalPixelType PixelComponentType
typename OutputImageType::IndexValueType IndexValueType
typename OutputImageType::PixelType PixelType
void SetOutputParametersFromImage(const ImageBaseType *image)
void PrintSelf(std::ostream &os, Indent indent) const override
typename DisplacementFieldType::Pointer DisplacementFieldPointer
void AfterThreadedGenerateData() override
SmartPointer< Self > Pointer
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
*par Constraints *The filter image with at least two dimensions and a vector *length of at least The theory supports extension to scalar but *the implementation of the itk vector classes do not **The template parameter TRealType must be floating point(float or double) or *a user-defined "real" numerical type with arithmetic operations defined *sufficient to compute derivatives. **\par Performance *This filter will automatically multithread if run with *SetUsePrincipleComponents
long IndexValueType
Definition: itkIntTypes.h:93