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 __itkSmoothingRecursiveGaussianImageFilter_h 00019 #define __itkSmoothingRecursiveGaussianImageFilter_h 00020 00021 #include "itkRecursiveGaussianImageFilter.h" 00022 #include "itkCastImageFilter.h" 00023 #include "itkImage.h" 00024 #include "itkPixelTraits.h" 00025 #include "itkCommand.h" 00026 00027 namespace itk 00028 { 00049 template< typename TInputImage, 00050 typename TOutputImage = TInputImage > 00051 class ITK_EXPORT SmoothingRecursiveGaussianImageFilter: 00052 public InPlaceImageFilter< TInputImage, TOutputImage > 00053 { 00054 public: 00056 typedef SmoothingRecursiveGaussianImageFilter Self; 00057 typedef InPlaceImageFilter< TInputImage, TOutputImage > Superclass; 00058 typedef SmartPointer< Self > Pointer; 00059 typedef SmartPointer< const Self > ConstPointer; 00060 00062 typedef TInputImage InputImageType; 00063 typedef TOutputImage OutputImageType; 00064 typedef typename TInputImage::PixelType PixelType; 00065 typedef typename NumericTraits< PixelType >::RealType RealType; 00066 typedef typename NumericTraits< PixelType >::ScalarRealType ScalarRealType; 00067 00069 itkTypeMacro(SmoothingRecursiveGaussianImageFilter, 00070 ImageToImageFilter); 00071 00073 itkStaticConstMacro(ImageDimension, unsigned int, 00074 TInputImage::ImageDimension); 00075 00077 typedef FixedArray< ScalarRealType, 00078 itkGetStaticConstMacro(ImageDimension) > SigmaArrayType; 00079 00084 typedef typename NumericTraits< PixelType >::FloatType InternalRealType; 00085 typedef typename InputImageType::template Rebind<InternalRealType>::Type RealImageType; 00086 00088 typedef RecursiveGaussianImageFilter< 00089 InputImageType, 00090 RealImageType 00091 > FirstGaussianFilterType; 00092 00094 typedef RecursiveGaussianImageFilter< 00095 RealImageType, 00096 RealImageType 00097 > InternalGaussianFilterType; 00098 00100 typedef CastImageFilter< 00101 RealImageType, 00102 OutputImageType 00103 > CastingFilterType; 00104 00106 typedef typename InternalGaussianFilterType::Pointer InternalGaussianFilterPointer; 00107 00109 typedef typename FirstGaussianFilterType::Pointer FirstGaussianFilterPointer; 00110 00112 typedef typename CastingFilterType::Pointer CastingFilterPointer; 00113 00115 typedef typename OutputImageType::Pointer OutputImagePointer; 00116 00118 itkNewMacro(Self); 00119 00124 void SetSigmaArray(const SigmaArrayType & sigmas); 00125 void SetSigma(ScalarRealType sigma); 00127 00128 SigmaArrayType GetSigmaArray() const; 00129 ScalarRealType GetSigma() const; 00130 00135 void SetNormalizeAcrossScale(bool normalizeInScaleSpace); 00136 itkGetConstMacro(NormalizeAcrossScale, bool); 00138 00139 // See super class for doxygen documentation 00140 // 00141 void SetNumberOfThreads(ThreadIdType nb); 00142 00143 // See super class for doxygen documentation 00144 // 00145 virtual bool CanRunInPlace( void ) const; 00146 00147 #ifdef ITK_USE_CONCEPT_CHECKING 00148 00149 // This concept does not work with variable length vector images 00150 //itkConceptMacro( InputHasNumericTraitsCheck, 00151 //( Concept::HasNumericTraits< PixelType > ) ); 00152 00154 #endif 00155 protected: 00156 SmoothingRecursiveGaussianImageFilter(); 00157 virtual ~SmoothingRecursiveGaussianImageFilter() {} 00158 void PrintSelf(std::ostream & os, Indent indent) const; 00160 00162 void GenerateData(void); 00163 00169 virtual void GenerateInputRequestedRegion() 00170 throw( InvalidRequestedRegionError ); 00171 00172 // Override since the filter produces the entire dataset 00173 void EnlargeOutputRequestedRegion(DataObject *output); 00174 00175 private: 00176 SmoothingRecursiveGaussianImageFilter(const Self &); //purposely not 00177 // implemented 00178 void operator=(const Self &); //purposely not 00179 00180 // implemented 00181 00182 InternalGaussianFilterPointer m_SmoothingFilters[ImageDimension - 1]; 00183 FirstGaussianFilterPointer m_FirstSmoothingFilter; 00184 CastingFilterPointer m_CastingFilter; 00185 00187 bool m_NormalizeAcrossScale; 00188 00190 SigmaArrayType m_Sigma; 00191 }; 00192 } // end namespace itk 00193 00194 #ifndef ITK_MANUAL_INSTANTIATION 00195 #include "itkSmoothingRecursiveGaussianImageFilter.hxx" 00196 #endif 00197 00198 #endif 00199