ITK  6.0.0
Insight Toolkit
itkVectorExpandImageFilter.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 itkVectorExpandImageFilter_h
19 #define itkVectorExpandImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
23 
24 #if !defined(ITK_LEGACY_REMOVE)
25 namespace itk
26 {
73 template <typename TInputImage, typename TOutputImage>
74 class ITK_TEMPLATE_EXPORT VectorExpandImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
75 {
76 public:
77  ITK_DISALLOW_COPY_AND_MOVE(VectorExpandImageFilter);
78 
80  using Self = VectorExpandImageFilter;
81  using Superclass = ImageToImageFilter<TInputImage, TOutputImage>;
82  using Pointer = SmartPointer<Self>;
83  using ConstPointer = SmartPointer<const Self>;
84 
86  itkNewMacro(Self);
87 
89  using InputImagePointer = typename TInputImage::Pointer;
90  using OutputImagePointer = typename TOutputImage::Pointer;
91  using OutputImageRegionType = typename TOutputImage::RegionType;
92 
94  itkOverrideGetNameOfClassMacro(VectorExpandImageFilter);
95 
97  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
98 
100  using typename Superclass::InputImageType;
101  using typename Superclass::OutputImageType;
102 
104  using OutputPixelType = typename OutputImageType::PixelType;
105  using OutputValueType = typename OutputPixelType::ValueType;
106  using InputPixelType = typename InputImageType::PixelType;
107  using InputValueType = typename InputPixelType::ValueType;
108 
110  enum
111  {
112  VectorDimension = InputPixelType::Dimension
113  };
114 
116  using ExpandFactorsType = float;
117  using ExpandFactorsArrayType = FixedArray<ExpandFactorsType, ImageDimension>;
118 
120  using CoordRepType = double;
121  using InterpolatorType = VectorInterpolateImageFunction<InputImageType, CoordRepType>;
122  using InterpolatorPointer = typename InterpolatorType::Pointer;
123  using DefaultInterpolatorType = VectorLinearInterpolateImageFunction<InputImageType, CoordRepType>;
124 
126  itkSetObjectMacro(Interpolator, InterpolatorType);
127  itkGetModifiableObjectMacro(Interpolator, InterpolatorType);
132  itkSetMacro(ExpandFactors, ExpandFactorsArrayType);
133  virtual void
134  SetExpandFactors(const float factor);
135  itkSetVectorMacro(ExpandFactors, const unsigned int, ImageDimension);
139  itkGetConstReferenceMacro(ExpandFactors, ExpandFactorsArrayType);
140 
147  void
148  GenerateOutputInformation() override;
149 
155  void
156  GenerateInputRequestedRegion() override;
157 
158 # ifdef ITK_USE_CONCEPT_CHECKING
159  // Begin concept checking
160  itkConceptMacro(InputHasNumericTraitsCheck, (Concept::HasNumericTraits<InputValueType>));
161  itkConceptMacro(OutputHasNumericTraitsCheck, (Concept::HasNumericTraits<OutputValueType>));
162  // End concept checking
163 # endif
164 
165 protected:
166  VectorExpandImageFilter();
167  ~VectorExpandImageFilter() override = default;
168  void
169  PrintSelf(std::ostream & os, Indent indent) const override;
170 
179  void
180  DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
181 
182 
184  void
185  BeforeThreadedGenerateData() override;
186 
187 private:
188  ExpandFactorsArrayType m_ExpandFactors{};
189  InterpolatorPointer m_Interpolator{};
190 };
191 } // end namespace itk
192 
193 # ifndef ITK_MANUAL_INSTANTIATION
194 # include "itkVectorExpandImageFilter.hxx"
195 # endif
196 
197 #endif // ITK_LEGACY_REMOVE
198 #endif
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itkVectorLinearInterpolateImageFunction.h
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itkImageToImageFilter.h
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:65
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
AddImageFilter
Definition: itkAddImageFilter.h:81
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition: itkGTestTypedefsAndConstructors.h:44
Superclass
BinaryGeneratorImageFilter< TInputImage1, TInputImage2, TOutputImage > Superclass
Definition: itkAddImageFilter.h:90