ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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
23
24namespace itk
25{
44template <typename TInputImage, typename TOutputImage>
45class ITK_TEMPLATE_EXPORT InterpolateImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
46{
47public:
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;
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
81 {
82 this->SetInput(image);
83 }
84 const InputImageType *
86 {
87 return this->GetInput();
88 }
89
90
92 void
93 SetInput2(const InputImageType * image);
94
95 const InputImageType *
97
100 itkSetClampMacro(Distance, double, 0.0, 1.0);
101 itkGetConstMacro(Distance, double);
103
105 itkSetObjectMacro(Interpolator, InterpolatorType);
106 itkGetModifiableObjectMacro(Interpolator, InterpolatorType);
108
114 void
116
118 void
120
121 itkConceptMacro(InputHasNumericTraitsCheck, (Concept::HasNumericTraits<InputPixelType>));
122
123protected:
125 ~InterpolateImageFilter() override = default;
126 void
127 PrintSelf(std::ostream & os, Indent indent) const override;
128
130 void
131 DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
132
133
134private:
136
138
139 double m_Distance{};
140};
141} // end namespace itk
142
143#ifndef ITK_MANUAL_INSTANTIATION
144# include "itkInterpolateImageFilter.hxx"
145#endif
146
147#endif
TOutputImage OutputImageType
typename OutputImageType::Pointer OutputImagePointer
typename OutputImageType::RegionType OutputImageRegionType
typename InputImageType::Pointer InputImagePointer
virtual void SetInput(const InputImageType *input)
const InputImageType * GetInput() const
Templated n-dimensional image class.
Definition itkImage.h:89
Control indentation during Print() invocation.
Definition itkIndent.h:50
static constexpr unsigned int ImageDimension
const InputImageType * GetInput2()
Image< InputPixelType, Self::IntermediateImageDimension > IntermediateImageType
SmartPointer< const Self > ConstPointer
void BeforeThreadedGenerateData() override
typename OutputImageType::RegionType OutputImageRegionType
static constexpr unsigned int IntermediateImageDimension
ImageToImageFilter< TInputImage, TOutputImage > Superclass
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
InterpolatorType::Pointer m_Interpolator
InterpolateImageFunction< IntermediateImageType > InterpolatorType
typename TInputImage::PixelType InputPixelType
~InterpolateImageFilter() override=default
void PrintSelf(std::ostream &os, Indent indent) const override
void AfterThreadedGenerateData() override
IntermediateImageType::Pointer m_IntermediateImage
void SetInput1(const InputImageType *image)
void SetInput2(const InputImageType *image)
Base class for all image interpolators.
Implements transparent reference counting.
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....