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
22
23namespace itk
24{
73template <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>
81class ITK_TEMPLATE_EXPORT SliceBySliceImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
82{
83public:
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;
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
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
152 itkGetModifiableObjectMacro(InputFilter, InputFilterType);
153
154 void
156 itkGetModifiableObjectMacro(OutputFilter, OutputFilterType);
157
162 itkGetConstMacro(SliceIndex, IndexValueType);
163
164protected:
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
179
180private:
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
Base class for all process objects that output image data.
TOutputImage OutputImageType
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:50
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Apply a filter or a pipeline slice by slice on an image.
typename TInputImage::PixelType PixelType
~SliceBySliceImageFilter() override=default
typename InternalInputImageType::IndexType InternalIndexType
typename InternalInputImageType::OffsetType InternalOffsetType
typename InternalInputImageType::SizeType InternalSizeType
typename TInputImage::SizeType SizeType
typename TInputImage::OffsetType OffsetType
void GenerateData() override
void VerifyInputInformation() const override
Verifies that the input images occupy the same physical space and the each index is at the same physi...
typename InternalInputImageType::PointType InternalPointType
void GenerateInputRequestedRegion() override
typename InternalInputImageType::RegionType InternalRegionType
void PrintSelf(std::ostream &os, Indent indent) const override
typename InternalInputImageType::PixelType InternalInputPixelType
const InputFilterType * GetFilter() const
typename TInputImage::RegionType RegionType
typename InternalOutputImageType::PixelType InternalOutputPixelType
typename InternalInputImageType::SpacingType InternalSpacingType
void SetFilter(InputFilterType *filter)
typename TInputImage::IndexType IndexType
typename TOutputImage::PixelType OutputPixelType
void SetOutputFilter(OutputFilterType *filter)
void SetInputFilter(InputFilterType *filter)
BinaryGeneratorImageFilter< TInputImage1, TInputImage2, TOutputImage > Superclass
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
long IndexValueType
Definition: itkIntTypes.h:93