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"
21 #include "itkImageToImageFilter.h"
23 
24 namespace itk
25 {
84 template <typename TInputImage, typename TOutputImage, typename TDisplacementField>
85 class ITK_TEMPLATE_EXPORT WarpImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
86 {
87 public:
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 
150  itkSetInputMacro(DisplacementField, DisplacementFieldType);
151 
153  itkGetInputMacro(DisplacementField, DisplacementFieldType);
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
185  SetOutputParametersFromImage(const ImageBaseType * image);
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
212  GenerateOutputInformation() override;
213 
220  void
221  GenerateInputRequestedRegion() override;
222 
225  void
226  BeforeThreadedGenerateData() override;
227 
230  void
231  AfterThreadedGenerateData() override;
232 
233 #ifdef ITK_USE_CONCEPT_CHECKING
234  // Begin concept checking
238  itkConceptMacro(DisplacementFieldHasNumericTraitsCheck,
240  // End concept checking
241 #endif
242 
243 protected:
244  WarpImageFilter();
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
275  EvaluateDisplacementAtPhysicalPoint(const PointType & point, DisplacementType & output);
276 
280  void
281  EvaluateDisplacementAtPhysicalPoint(const PointType & point,
282  const DisplacementFieldType * fieldPtr,
283  DisplacementType & output);
284 
285  bool m_DefFieldSameInformation{};
286  // variables for deffield interpolator
287  IndexType m_StartIndex, m_EndIndex{};
288 
289 private:
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
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::GTest::TypedefsAndConstructors::Dimension2::DirectionType
ImageBaseType::DirectionType DirectionType
Definition: itkGTestTypedefsAndConstructors.h:52
itk::Concept::HasNumericTraits
Definition: itkConceptChecking.h:717
itkLinearInterpolateImageFunction.h
itk::WarpImageFilter::m_StartIndex
IndexType m_StartIndex
Definition: itkWarpImageFilter.h:287
itk::WarpImageFilter::DirectionType
typename TOutputImage::DirectionType DirectionType
Definition: itkWarpImageFilter.h:146
itk::ImageBase
Base class for templated image classes.
Definition: itkImageBase.h:114
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::WarpImageFilter::IndexValueType
typename OutputImageType::IndexValueType IndexValueType
Definition: itkWarpImageFilter.h:113
itk::WarpImageFilter::PixelType
typename OutputImageType::PixelType PixelType
Definition: itkWarpImageFilter.h:115
itk::IndexValueType
long IndexValueType
Definition: itkIntTypes.h:93
itk::Concept::SameDimension
Definition: itkConceptChecking.h:697
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::WarpImageFilter::PixelComponentType
typename OutputImageType::InternalPixelType PixelComponentType
Definition: itkWarpImageFilter.h:116
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::WarpImageFilter::DisplacementType
typename DisplacementFieldType::PixelType DisplacementType
Definition: itkWarpImageFilter.h:130
itk::LinearInterpolateImageFunction
Linearly interpolate an image at specified positions.
Definition: itkLinearInterpolateImageFunction.h:51
itk::WarpImageFilter::DisplacementFieldType
TDisplacementField DisplacementFieldType
Definition: itkWarpImageFilter.h:128
itk::point
*par Constraints *The filter requires an 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
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::WarpImageFilter::CoordinateType
double CoordinateType
Definition: itkWarpImageFilter.h:133
itk::WarpImageFilter::SizeType
typename OutputImageType::SizeType SizeType
Definition: itkWarpImageFilter.h:114
itkImageToImageFilter.h
itk::WarpImageFilter::DisplacementFieldPointer
typename DisplacementFieldType::Pointer DisplacementFieldPointer
Definition: itkWarpImageFilter.h:129
itk::ImageSource::OutputImageRegionType
typename OutputImageType::RegionType OutputImageRegionType
Definition: itkImageSource.h:92
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:65
itk::WarpImageFilter
Warps an image using an input displacement field.
Definition: itkWarpImageFilter.h:85
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::WarpImageFilter::InterpolatorPointer
typename InterpolatorType::Pointer InterpolatorPointer
Definition: itkWarpImageFilter.h:139
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:139
itk::WarpImageFilter::SpacingType
typename OutputImageType::SpacingType SpacingType
Definition: itkWarpImageFilter.h:117
itk::Point< CoordinateType, Self::ImageDimension >
itk::TransformCategory::DisplacementField
itk::InterpolateImageFunction
Base class for all image interpolators.
Definition: itkInterpolateImageFunction.h:45
itk::WarpImageFilter::IndexType
typename OutputImageType::IndexType IndexType
Definition: itkWarpImageFilter.h:112
itkImageBase.h