ITK  6.0.0
Insight Toolkit
itkInterpolateImageFilter.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 itkInterpolateImageFilter_h
19 #define itkInterpolateImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
23 
24 namespace itk
25 {
44 template <typename TInputImage, typename TOutputImage>
45 class ITK_TEMPLATE_EXPORT InterpolateImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
46 {
47 public:
48  ITK_DISALLOW_COPY_AND_MOVE(InterpolateImageFilter);
49 
55 
57  itkNewMacro(Self);
58 
60  itkOverrideGetNameOfClassMacro(InterpolateImageFilter);
61 
63  using typename Superclass::InputImageType;
64  using typename Superclass::InputImagePointer;
65  using typename Superclass::OutputImageType;
66  using typename Superclass::OutputImagePointer;
67  using typename Superclass::OutputImageRegionType;
68 
70  static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
71  static constexpr unsigned int IntermediateImageDimension = TOutputImage::ImageDimension + 1;
72 
74  using InputPixelType = typename TInputImage::PixelType;
77 
79  void
80  SetInput1(const InputImageType * image)
81  {
82  this->SetInput(image);
83  }
84  const InputImageType *
86  {
87  return this->GetInput();
88  }
92  void
93  SetInput2(const InputImageType * image);
94 
95  const InputImageType *
96  GetInput2();
97 
100  itkSetClampMacro(Distance, double, 0.0, 1.0);
101  itkGetConstMacro(Distance, double);
105  itkSetObjectMacro(Interpolator, InterpolatorType);
106  itkGetModifiableObjectMacro(Interpolator, InterpolatorType);
114  void
115  BeforeThreadedGenerateData() override;
116 
118  void
119  AfterThreadedGenerateData() override;
120 
121 #ifdef ITK_USE_CONCEPT_CHECKING
122  // Begin concept checking
123  itkConceptMacro(InputHasNumericTraitsCheck, (Concept::HasNumericTraits<InputPixelType>));
124  // End concept checking
125 #endif
126 
127 protected:
129  ~InterpolateImageFilter() override = default;
130  void
131  PrintSelf(std::ostream & os, Indent indent) const override;
132 
134  void
135  DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
136 
137 
138 private:
139  typename InterpolatorType::Pointer m_Interpolator{};
140 
141  typename IntermediateImageType::Pointer m_IntermediateImage{};
142 
143  double m_Distance{};
144 };
145 } // end namespace itk
146 
147 #ifndef ITK_MANUAL_INSTANTIATION
148 # include "itkInterpolateImageFilter.hxx"
149 #endif
150 
151 #endif
itk::Concept::HasNumericTraits
Definition: itkConceptChecking.h:716
itkLinearInterpolateImageFunction.h
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::InterpolateImageFilter::GetInput1
const InputImageType * GetInput1()
Definition: itkInterpolateImageFilter.h:85
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::InterpolateImageFilter::InputPixelType
typename TInputImage::PixelType InputPixelType
Definition: itkInterpolateImageFilter.h:74
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itkImageToImageFilter.h
itk::InterpolateImageFilter
Interpolate an image from two N-D images.
Definition: itkInterpolateImageFilter.h:45
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:65
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:139
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itk::InterpolateImageFilter::SetInput1
void SetInput1(const InputImageType *image)
Definition: itkInterpolateImageFilter.h:80
itk::InterpolateImageFunction
Base class for all image interpolators.
Definition: itkInterpolateImageFunction.h:45