ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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);
90
96
98 itkNewMacro(Self);
99
101 itkOverrideGetNameOfClassMacro(WarpImageFilter);
102
104 using OutputImageRegionType = typename TOutputImage::RegionType;
105
107 using typename Superclass::InputImageType;
108 using typename Superclass::InputImagePointer;
109 using typename Superclass::OutputImageType;
110 using typename Superclass::OutputImagePointer;
112 using IndexType = typename OutputImageType::IndexType;
113 using IndexValueType = typename OutputImageType::IndexValueType;
114 using SizeType = typename OutputImageType::SizeType;
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;
129 using DisplacementFieldPointer = typename DisplacementFieldType::Pointer;
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
146 using DirectionType = typename TOutputImage::DirectionType;
147
148
151
154
156 itkSetObjectMacro(Interpolator, InterpolatorType);
157 itkGetModifiableObjectMacro(Interpolator, InterpolatorType);
159
161 itkSetMacro(OutputSpacing, SpacingType);
162 virtual void
163 SetOutputSpacing(const double * spacing);
165
167 itkGetConstReferenceMacro(OutputSpacing, SpacingType);
168
170 itkSetMacro(OutputOrigin, PointType);
171 virtual void
172 SetOutputOrigin(const double * origin);
174
176 itkGetConstReferenceMacro(OutputOrigin, PointType);
177
179 itkSetMacro(OutputDirection, DirectionType);
180 itkGetConstReferenceMacro(OutputDirection, DirectionType);
182
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
236 itkConceptMacro(DisplacementFieldHasNumericTraitsCheck,
238
239protected:
241 // ~WarpImageFilter() {} default implementation is ok
242
243 void
244 PrintSelf(std::ostream & os, Indent indent) const override;
245
249 void
250 DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
251
252
259 void
260 VerifyInputInformation() const override;
261
270 void
272
276 void
278 const DisplacementFieldType * fieldPtr,
279 DisplacementType & output);
280
282 // variables for deffield interpolator
284
285private:
290
292 SizeType m_OutputSize{}; // Size of the output image
293 IndexType m_OutputStartIndex{}; // output image start index
294};
295} // end namespace itk
296
297#ifndef ITK_MANUAL_INSTANTIATION
298# include "itkWarpImageFilter.hxx"
299#endif
300
301#endif
Base class for templated image classes.
TOutputImage OutputImageType
typename OutputImageType::Pointer OutputImagePointer
typename InputImageType::ConstPointer InputImageConstPointer
typename InputImageType::Pointer InputImagePointer
Control indentation during Print() invocation.
Definition itkIndent.h:50
Base class for all image interpolators.
Linearly interpolate an image at specified positions.
A templated class holding a geometric point in n-Dimensional space.
Definition itkPoint.h:54
Implements transparent reference counting.
void GenerateInputRequestedRegion() override
typename OutputImageType::SizeType SizeType
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
SmartPointer< Self > Pointer
ImageToImageFilter< TInputImage, TOutputImage > Superclass
void GenerateOutputInformation() override
void BeforeThreadedGenerateData() override
TDisplacementField DisplacementFieldType
typename InterpolatorType::Pointer InterpolatorPointer
itkGetInputMacro(DisplacementField, DisplacementFieldType)
typename OutputImageType::SpacingType SpacingType
ImageBase< Self::ImageDimension > ImageBaseType
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
Point< CoordinateType, Self::ImageDimension > PointType
typename DisplacementFieldType::PixelType DisplacementType
virtual void SetOutputOrigin(const double *origin)
InterpolateImageFunction< InputImageType, CoordinateType > InterpolatorType
void EvaluateDisplacementAtPhysicalPoint(const PointType &point, const DisplacementFieldType *fieldPtr, DisplacementType &output)
itkSetInputMacro(DisplacementField, DisplacementFieldType)
typename OutputImageType::InternalPixelType PixelComponentType
LinearInterpolateImageFunction< InputImageType, CoordinateType > DefaultInterpolatorType
typename OutputImageType::IndexValueType IndexValueType
typename OutputImageType::PixelType PixelType
void SetOutputParametersFromImage(const ImageBaseType *image)
SmartPointer< const Self > ConstPointer
void PrintSelf(std::ostream &os, Indent indent) const override
typename DisplacementFieldType::Pointer DisplacementFieldPointer
typename TOutputImage::RegionType OutputImageRegionType
void AfterThreadedGenerateData() override
#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