ITK
4.1.0
Insight Segmentation and Registration Toolkit
|
00001 /*========================================================================= 00002 * 00003 * Copyright Insight Software Consortium 00004 * 00005 * Licensed under the Apache License, Version 2.0 (the "License"); 00006 * you may not use this file except in compliance with the License. 00007 * You may obtain a copy of the License at 00008 * 00009 * http://www.apache.org/licenses/LICENSE-2.0.txt 00010 * 00011 * Unless required by applicable law or agreed to in writing, software 00012 * distributed under the License is distributed on an "AS IS" BASIS, 00013 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. 00014 * See the License for the specific language governing permissions and 00015 * limitations under the License. 00016 * 00017 *=========================================================================*/ 00018 #ifndef __itkSliceBySliceImageFilter_h 00019 #define __itkSliceBySliceImageFilter_h 00020 00021 #include "itkImageToImageFilter.h" 00022 00023 namespace itk 00024 { 00072 template< class TInputImage, 00073 class TOutputImage, 00074 class TInputFilter = ImageToImageFilter< 00075 Image< typename TInputImage::PixelType, ::itk::GetImageDimension< TInputImage >::ImageDimension - 1 >, 00076 Image< typename TOutputImage::PixelType, ::itk::GetImageDimension< TOutputImage >::ImageDimension - 1 > >, 00077 class TOutputFilter = typename TInputFilter::Superclass, 00078 class TInternalInputImage = typename TInputFilter::InputImageType, 00079 class TInternalOutputImage = typename TOutputFilter::OutputImageType > 00080 class ITK_EXPORT SliceBySliceImageFilter: 00081 public ImageToImageFilter< TInputImage, TOutputImage > 00082 { 00083 public: 00085 typedef SliceBySliceImageFilter Self; 00086 typedef ImageToImageFilter< TInputImage, TOutputImage > Superclass; 00087 typedef SmartPointer< Self > Pointer; 00088 typedef SmartPointer< const Self > ConstPointer; 00089 00091 typedef typename Superclass::InputImagePointer InputImagePointer; 00092 00094 itkNewMacro(Self); 00095 00097 itkTypeMacro(SliceBySliceImageFilter, ImageToImageFilter); 00098 00100 typedef TInputImage InputImageType; 00101 typedef typename TInputImage::RegionType RegionType; 00102 typedef typename TInputImage::SizeType SizeType; 00103 typedef typename TInputImage::IndexType IndexType; 00104 typedef typename TInputImage::PixelType PixelType; 00105 typedef typename TInputImage::OffsetType OffsetType; 00106 00107 typedef TOutputImage OutputImageType; 00108 typedef typename TOutputImage::PixelType OutputPixelType; 00109 00110 typedef TInputFilter InputFilterType; 00111 typedef TOutputFilter OutputFilterType; 00112 00113 typedef TInternalInputImage InternalInputImageType; 00114 typedef typename InternalInputImageType::RegionType InternalRegionType; 00115 typedef typename InternalInputImageType::SizeType InternalSizeType; 00116 typedef typename InternalInputImageType::IndexType InternalIndexType; 00117 typedef typename InternalInputImageType::OffsetType InternalOffsetType; 00118 typedef typename InternalInputImageType::PixelType InternalInputPixelType; 00119 00120 typedef TInternalOutputImage InternalOutputImageType; 00121 typedef typename InternalOutputImageType::PixelType InternalOutputPixelType; 00122 00124 itkStaticConstMacro(ImageDimension, unsigned int, 00125 TInputImage::ImageDimension); 00126 00127 itkStaticConstMacro(InternalImageDimension, unsigned int, 00128 InternalInputImageType::ImageDimension); 00129 00130 itkSetMacro(Dimension, unsigned int); 00131 itkGetConstMacro(Dimension, unsigned int); 00132 00133 void SetFilter(InputFilterType *filter); 00134 00135 InputFilterType * GetFilter() 00136 { 00137 return this->m_InputFilter; 00138 } 00139 00140 const InputFilterType * GetFilter() const 00141 { 00142 return this->m_InputFilter; 00143 } 00144 00145 void SetInputFilter(InputFilterType *filter); 00146 00147 itkGetObjectMacro(InputFilter, InputFilterType); 00148 00149 void SetOutputFilter(OutputFilterType *filter); 00150 00151 itkGetObjectMacro(OutputFilter, OutputFilterType); 00152 00157 itkGetConstMacro(SliceIndex, IndexValueType); 00158 protected: 00159 SliceBySliceImageFilter(); 00160 ~SliceBySliceImageFilter() {} 00162 00163 void GenerateData(); 00164 00165 void PrintSelf(std::ostream & os, Indent indent) const; 00166 00167 void GenerateInputRequestedRegion(); 00168 00169 void EnlargeOutputRequestedRegion( DataObject *itkNotUsed(output) ); 00170 private: 00171 SliceBySliceImageFilter(const Self &); //purposely not implemented 00172 void operator=(const Self &); //purposely not implemented 00173 00174 unsigned int m_Dimension; 00175 00176 typename InputFilterType::Pointer m_InputFilter; 00177 00178 typename OutputFilterType::Pointer m_OutputFilter; 00179 00180 IndexValueType m_SliceIndex; 00181 }; 00182 } 00183 00184 #ifndef ITK_MANUAL_INSTANTIATION 00185 #include "itkSliceBySliceImageFilter.hxx" 00186 #endif 00187 00188 #endif 00189