ITK  6.0.0
Insight Toolkit
itkExpandImageFilter.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 itkExpandImageFilter_h
19#define itkExpandImageFilter_h
20
23
24namespace itk
25{
63template <typename TInputImage, typename TOutputImage>
64class ITK_TEMPLATE_EXPORT ExpandImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
65{
66public:
67 ITK_DISALLOW_COPY_AND_MOVE(ExpandImageFilter);
68
74
76 itkNewMacro(Self);
77
79 itkOverrideGetNameOfClassMacro(ExpandImageFilter);
80
83
85 static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
86
88 using typename Superclass::InputImageType;
89 using typename Superclass::OutputImageType;
90 using OutputPixelType = typename OutputImageType::PixelType;
93
95 using CoordinateType = double;
96#ifndef ITK_FUTURE_LEGACY_REMOVE
97 using CoordRepType ITK_FUTURE_DEPRECATED(
98 "ITK 6 discourages using `CoordRepType`. Please use `CoordinateType` instead!") = CoordinateType;
99#endif
103
105 itkSetObjectMacro(Interpolator, InterpolatorType);
106 itkGetModifiableObjectMacro(Interpolator, InterpolatorType);
111
114 itkSetMacro(ExpandFactors, ExpandFactorsType);
115 virtual void
116 SetExpandFactors(const unsigned int factor);
120 itkGetConstReferenceMacro(ExpandFactors, ExpandFactorsType);
121
122
129 void
131
137 void
139
140protected:
142 ~ExpandImageFilter() override = default;
143 void
144 PrintSelf(std::ostream & os, Indent indent) const override;
145
156 void
157 DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
158
159
162 void
164
165private:
166 ExpandFactorsType m_ExpandFactors{};
167 InterpolatorPointer m_Interpolator{};
168};
169} // end namespace itk
170
171#ifndef ITK_MANUAL_INSTANTIATION
172# include "itkExpandImageFilter.hxx"
173#endif
174
175#endif
Expand the size of an image by an integer factor in each dimension.
void GenerateOutputInformation() override
void PrintSelf(std::ostream &os, Indent indent) const override
void BeforeThreadedGenerateData() override
typename InterpolatorType::Pointer InterpolatorPointer
~ExpandImageFilter() override=default
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
void GenerateInputRequestedRegion() override
virtual void SetExpandFactors(const unsigned int factor)
typename OutputImageType::PixelType OutputPixelType
Base class for all process objects that output image data.
typename OutputImageType::RegionType OutputImageRegionType
typename OutputImageType::Pointer OutputImagePointer
Base class for filters that take an image as input and produce an image as output.
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.
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....