ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkExpandImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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  * http://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 
21 #include "itkImageToImageFilter.h"
23 
24 namespace itk
25 {
66 template< typename TInputImage, typename TOutputImage >
67 class ITK_TEMPLATE_EXPORT ExpandImageFilter:
68  public ImageToImageFilter< TInputImage, TOutputImage >
69 {
70 public:
71  ITK_DISALLOW_COPY_AND_ASSIGN(ExpandImageFilter);
72 
78 
80  itkNewMacro(Self);
81 
84 
87 
89  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
90 
92  using InputImageType = typename Superclass::InputImageType;
93  using OutputImageType = typename Superclass::OutputImageType;
94  using OutputPixelType = typename OutputImageType::PixelType;
95  using InputImagePointer = typename InputImageType::Pointer;
96  using OutputImagePointer = typename OutputImageType::Pointer;
97 
99  using CoordRepType = double;
104 
106  itkSetObjectMacro(Interpolator, InterpolatorType);
107  itkGetModifiableObjectMacro(Interpolator, InterpolatorType);
109 
112 
115  itkSetMacro(ExpandFactors, ExpandFactorsType);
116  virtual void SetExpandFactors(const unsigned int factor);
118 
120  itkGetConstReferenceMacro(ExpandFactors, ExpandFactorsType);
121 
122 
129  void GenerateOutputInformation() override;
130 
136  void GenerateInputRequestedRegion() override;
137 
138 #ifdef ITK_USE_CONCEPT_CHECKING
139  // Begin concept checking
140  itkConceptMacro( InputHasNumericTraitsCheck,
142  itkConceptMacro( OutputHasNumericTraitsCheck,
144  // End concept checking
145 #endif
146 
147 protected:
149  ~ExpandImageFilter() override = default;
150  void PrintSelf(std::ostream & os, Indent indent) const override;
151 
162  void DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
163 
164 
167  void BeforeThreadedGenerateData() override;
168 
169 private:
172 };
173 } // end namespace itk
174 
175 #ifndef ITK_MANUAL_INSTANTIATION
176 #include "itkExpandImageFilter.hxx"
177 #endif
178 
179 #endif
typename OutputImageType::Pointer OutputImagePointer
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Base class for all process objects that output image data.
typename OutputImageType::PixelType OutputPixelType
typename InputImageType::Pointer InputImagePointer
typename OutputImageType::RegionType OutputImageRegionType
TOutputImage OutputImageType
Linearly interpolate an image at specified positions.
Base class for all image interpolaters.
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
ExpandFactorsType m_ExpandFactors
typename InterpolatorType::Pointer InterpolatorPointer
#define itkConceptMacro(name, concept)
Expand the size of an image by an integer factor in each dimension.
InterpolatorPointer m_Interpolator