ITK  6.0.0
Insight Toolkit
itkSliceBySliceImageFilter.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 itkSliceBySliceImageFilter_h
19 #define itkSliceBySliceImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 
23 namespace itk
24 {
73 template <typename TInputImage,
74  typename TOutputImage,
75  typename TInputFilter =
76  ImageToImageFilter<Image<typename TInputImage::PixelType, TInputImage::ImageDimension - 1>,
77  Image<typename TOutputImage::PixelType, TOutputImage::ImageDimension - 1>>,
78  class TOutputFilter = typename TInputFilter::Superclass,
79  class TInternalInputImage = typename TInputFilter::InputImageType,
80  class TInternalOutputImage = typename TOutputFilter::OutputImageType>
81 class ITK_TEMPLATE_EXPORT SliceBySliceImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
82 {
83 public:
84  ITK_DISALLOW_COPY_AND_MOVE(SliceBySliceImageFilter);
85 
91 
93  using typename Superclass::InputImagePointer;
94 
96  itkNewMacro(Self);
97 
99  itkOverrideGetNameOfClassMacro(SliceBySliceImageFilter);
100 
102  using InputImageType = TInputImage;
104  using SizeType = typename TInputImage::SizeType;
106  using PixelType = typename TInputImage::PixelType;
107  using OffsetType = typename TInputImage::OffsetType;
108 
109  using OutputImageType = TOutputImage;
110  using OutputPixelType = typename TOutputImage::PixelType;
111 
112  using InputFilterType = TInputFilter;
113  using OutputFilterType = TOutputFilter;
114 
115  using InternalInputImageType = TInternalInputImage;
119  using InternalOffsetType = typename InternalInputImageType::OffsetType;
120  using InternalInputPixelType = typename InternalInputImageType::PixelType;
121  using InternalSpacingType = typename InternalInputImageType::SpacingType;
123 
124  using InternalOutputImageType = TInternalOutputImage;
125  using InternalOutputPixelType = typename InternalOutputImageType::PixelType;
126 
128  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
129 
130  static constexpr unsigned int InternalImageDimension = InternalInputImageType::ImageDimension;
131 
132  itkSetMacro(Dimension, unsigned int);
133  itkGetConstMacro(Dimension, unsigned int);
134 
135  void
136  SetFilter(InputFilterType * filter);
137 
140  {
141  return this->m_InputFilter;
142  }
143 
144  const InputFilterType *
145  GetFilter() const
146  {
147  return this->m_InputFilter;
148  }
149 
150  void
151  SetInputFilter(InputFilterType * filter);
152  itkGetModifiableObjectMacro(InputFilter, InputFilterType);
153 
154  void
155  SetOutputFilter(OutputFilterType * filter);
156  itkGetModifiableObjectMacro(OutputFilter, OutputFilterType);
157 
162  itkGetConstMacro(SliceIndex, IndexValueType);
163 
164 protected:
166  ~SliceBySliceImageFilter() override = default;
167 
168  void
169  VerifyInputInformation() const override;
170 
171  void
172  GenerateData() override;
173 
174  void
175  PrintSelf(std::ostream & os, Indent indent) const override;
176 
177  void
178  GenerateInputRequestedRegion() override;
179 
180 private:
181  unsigned int m_Dimension{};
182 
183  typename InputFilterType::Pointer m_InputFilter{};
184 
185  typename OutputFilterType::Pointer m_OutputFilter{};
186 
187  IndexValueType m_SliceIndex{};
188 };
189 } // namespace itk
190 
191 #ifndef ITK_MANUAL_INSTANTIATION
192 # include "itkSliceBySliceImageFilter.hxx"
193 #endif
194 
195 #endif
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::SliceBySliceImageFilter::InternalPointType
typename InternalInputImageType::PointType InternalPointType
Definition: itkSliceBySliceImageFilter.h:122
itk::SliceBySliceImageFilter::RegionType
typename TInputImage::RegionType RegionType
Definition: itkSliceBySliceImageFilter.h:103
itk::SliceBySliceImageFilter::InternalInputPixelType
typename InternalInputImageType::PixelType InternalInputPixelType
Definition: itkSliceBySliceImageFilter.h:120
itk::SliceBySliceImageFilter::InternalOutputPixelType
typename InternalOutputImageType::PixelType InternalOutputPixelType
Definition: itkSliceBySliceImageFilter.h:125
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::SliceBySliceImageFilter::InternalIndexType
typename InternalInputImageType::IndexType InternalIndexType
Definition: itkSliceBySliceImageFilter.h:118
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::SliceBySliceImageFilter::InternalOffsetType
typename InternalInputImageType::OffsetType InternalOffsetType
Definition: itkSliceBySliceImageFilter.h:119
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::SliceBySliceImageFilter::IndexType
typename TInputImage::IndexType IndexType
Definition: itkSliceBySliceImageFilter.h:105
itk::SliceBySliceImageFilter::OutputFilterType
TOutputFilter OutputFilterType
Definition: itkSliceBySliceImageFilter.h:113
itk::IndexValueType
long IndexValueType
Definition: itkIntTypes.h:93
itk::SliceBySliceImageFilter::InternalRegionType
typename InternalInputImageType::RegionType InternalRegionType
Definition: itkSliceBySliceImageFilter.h:116
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
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::SliceBySliceImageFilter::InternalSizeType
typename InternalInputImageType::SizeType InternalSizeType
Definition: itkSliceBySliceImageFilter.h:117
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::SliceBySliceImageFilter::OutputPixelType
typename TOutputImage::PixelType OutputPixelType
Definition: itkSliceBySliceImageFilter.h:110
itk::SliceBySliceImageFilter::GetFilter
InputFilterType * GetFilter()
Definition: itkSliceBySliceImageFilter.h:139
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itkImageToImageFilter.h
itk::SliceBySliceImageFilter::InternalInputImageType
TInternalInputImage InternalInputImageType
Definition: itkSliceBySliceImageFilter.h:115
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::SliceBySliceImageFilter::InputFilterType
TInputFilter InputFilterType
Definition: itkSliceBySliceImageFilter.h:112
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:139
itk::SliceBySliceImageFilter::SizeType
typename TInputImage::SizeType SizeType
Definition: itkSliceBySliceImageFilter.h:104
itk::SliceBySliceImageFilter::InternalSpacingType
typename InternalInputImageType::SpacingType InternalSpacingType
Definition: itkSliceBySliceImageFilter.h:121
itk::SliceBySliceImageFilter::InternalOutputImageType
TInternalOutputImage InternalOutputImageType
Definition: itkSliceBySliceImageFilter.h:124
itk::SliceBySliceImageFilter::GetFilter
const InputFilterType * GetFilter() const
Definition: itkSliceBySliceImageFilter.h:145
itk::SliceBySliceImageFilter
Apply a filter or a pipeline slice by slice on an image.
Definition: itkSliceBySliceImageFilter.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
itk::SliceBySliceImageFilter::OffsetType
typename TInputImage::OffsetType OffsetType
Definition: itkSliceBySliceImageFilter.h:107
itk::SliceBySliceImageFilter::PixelType
typename TInputImage::PixelType PixelType
Definition: itkSliceBySliceImageFilter.h:106
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90