ITK  4.9.0
Insight Segmentation and Registration Toolkit
itkSliceBySliceImageFilter.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 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 = ImageToImageFilter<
76  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 >
82  public ImageToImageFilter< TInputImage, TOutputImage >
83 {
84 public:
90 
93 
95  itkNewMacro(Self);
96 
99 
101  typedef TInputImage InputImageType;
102  typedef typename TInputImage::RegionType RegionType;
103  typedef typename TInputImage::SizeType SizeType;
104  typedef typename TInputImage::IndexType IndexType;
105  typedef typename TInputImage::PixelType PixelType;
106  typedef typename TInputImage::OffsetType OffsetType;
107 
108  typedef TOutputImage OutputImageType;
109  typedef typename TOutputImage::PixelType OutputPixelType;
110 
111  typedef TInputFilter InputFilterType;
112  typedef TOutputFilter OutputFilterType;
113 
114  typedef TInternalInputImage InternalInputImageType;
115  typedef typename InternalInputImageType::RegionType InternalRegionType;
116  typedef typename InternalInputImageType::SizeType InternalSizeType;
117  typedef typename InternalInputImageType::IndexType InternalIndexType;
118  typedef typename InternalInputImageType::OffsetType InternalOffsetType;
119  typedef typename InternalInputImageType::PixelType InternalInputPixelType;
120  typedef typename InternalInputImageType::SpacingType InternalSpacingType;
121  typedef typename InternalInputImageType::PointType InternalPointType;
122 
123  typedef TInternalOutputImage InternalOutputImageType;
124  typedef typename InternalOutputImageType::PixelType InternalOutputPixelType;
125 
127  itkStaticConstMacro(ImageDimension, unsigned int,
128  TInputImage::ImageDimension);
129 
130  itkStaticConstMacro(InternalImageDimension, unsigned int,
131  InternalInputImageType::ImageDimension);
132 
133  itkSetMacro(Dimension, unsigned int);
134  itkGetConstMacro(Dimension, unsigned int);
135 
136  void SetFilter(InputFilterType *filter);
137 
139  {
140  return this->m_InputFilter;
141  }
142 
143  const InputFilterType * GetFilter() const
144  {
145  return this->m_InputFilter;
146  }
147 
148  void SetInputFilter(InputFilterType *filter);
149  itkGetModifiableObjectMacro(InputFilter, InputFilterType);
150 
151  void SetOutputFilter(OutputFilterType *filter);
152  itkGetModifiableObjectMacro(OutputFilter, OutputFilterType);
153 
158  itkGetConstMacro(SliceIndex, IndexValueType);
159 
160 protected:
163 
164  void VerifyInputInformation() ITK_OVERRIDE;
165 
166  void GenerateData() ITK_OVERRIDE;
167 
168  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
169 
170  virtual void GenerateInputRequestedRegion() ITK_OVERRIDE;
171 
172 private:
173  SliceBySliceImageFilter(const Self &) ITK_DELETE_FUNCTION;
174  void operator=(const Self &) ITK_DELETE_FUNCTION;
175 
176  unsigned int m_Dimension;
177 
179 
181 
183 };
184 }
185 
186 #ifndef ITK_MANUAL_INSTANTIATION
187 #include "itkSliceBySliceImageFilter.hxx"
188 #endif
189 
190 #endif
static const unsigned int ImageDimension
void PrintSelf(std::ostream &os, Indent indent) const override
void GenerateData() override
void VerifyInputInformation() override
Verifies that the input images occupy the same physical space and the each index is at the same physi...
InternalOutputImageType::PixelType InternalOutputPixelType
InputImageType::Pointer InputImagePointer
signed long IndexValueType
Definition: itkIntTypes.h:150
void SetInputFilter(InputFilterType *filter)
void SetOutputFilter(OutputFilterType *filter)
static const unsigned int InternalImageDimension
Apply a filter or a pipeline slice by slice on an image.
InternalInputImageType::RegionType InternalRegionType
Base class for all process objects that output image data.
Superclass::InputImagePointer InputImagePointer
ImageToImageFilter< TInputImage, TOutputImage > Superclass
InternalInputImageType::IndexType InternalIndexType
void SetFilter(InputFilterType *filter)
const InputFilterType * GetFilter() const
InputFilterType::Pointer m_InputFilter
InternalInputImageType::SpacingType InternalSpacingType
InternalInputImageType::OffsetType InternalOffsetType
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
OutputFilterType::Pointer m_OutputFilter
virtual void GenerateInputRequestedRegion() override
InternalInputImageType::SizeType InternalSizeType
InternalInputImageType::PointType InternalPointType
SmartPointer< const Self > ConstPointer
InternalInputImageType::PixelType InternalInputPixelType