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 __itkGradientMagnitudeRecursiveGaussianImageFilter_h 00019 #define __itkGradientMagnitudeRecursiveGaussianImageFilter_h 00020 00021 #include "itkNthElementImageAdaptor.h" 00022 #include "itkImage.h" 00023 #include "itkPixelTraits.h" 00024 #include "itkRecursiveGaussianImageFilter.h" 00025 #include "itkSqrtImageFilter.h" 00026 #include "itkBinaryFunctorImageFilter.h" 00027 00028 namespace itk 00029 { 00045 // NOTE that the typename macro has to be used here in lieu 00046 // of "typename" because VC++ doesn't like the typename keyword 00047 // on the defaults of template parameters 00048 template< typename TInputImage, 00049 typename TOutputImage = TInputImage > 00050 class ITK_EXPORT GradientMagnitudeRecursiveGaussianImageFilter: 00051 public InPlaceImageFilter< TInputImage, TOutputImage > 00052 { 00053 public: 00055 typedef GradientMagnitudeRecursiveGaussianImageFilter Self; 00056 typedef InPlaceImageFilter< TInputImage, TOutputImage > Superclass; 00057 typedef SmartPointer< Self > Pointer; 00058 typedef SmartPointer< const Self > ConstPointer; 00059 00061 typedef TInputImage InputImageType; 00062 typedef typename InputImageType::PixelType PixelType; 00063 00065 itkStaticConstMacro(ImageDimension, unsigned int, 00066 TInputImage::ImageDimension); 00067 00068 typedef typename NumericTraits< PixelType >::RealType RealType; 00069 00074 typedef float InternalRealType; 00075 typedef Image< InternalRealType, 00076 itkGetStaticConstMacro(ImageDimension) > RealImageType; 00077 00079 typedef RecursiveGaussianImageFilter< 00080 RealImageType, 00081 RealImageType 00082 > GaussianFilterType; 00083 00085 typedef RecursiveGaussianImageFilter< 00086 InputImageType, 00087 RealImageType 00088 > DerivativeFilterType; 00089 00091 typedef SqrtImageFilter< 00092 RealImageType, 00093 TOutputImage 00094 > SqrtFilterType; 00095 00097 typedef typename GaussianFilterType::Pointer GaussianFilterPointer; 00098 00100 typedef typename DerivativeFilterType::Pointer DerivativeFilterPointer; 00101 00102 typedef typename SqrtFilterType::Pointer SqrtFilterPointer; 00103 00105 typedef typename TOutputImage::Pointer OutputImagePointer; 00106 00108 typedef TOutputImage OutputImageType; 00109 typedef typename OutputImageType::PixelType OutputPixelType; 00110 00113 typedef Image< InternalRealType, 00114 itkGetStaticConstMacro(ImageDimension) > CumulativeImageType; 00115 typedef typename CumulativeImageType::Pointer CumulativeImagePointer; 00116 00118 itkNewMacro(Self); 00119 00121 itkTypeMacro(GradientMagnitudeRecursiveGaussianImageFilter, 00122 InPlaceImageFilter); 00123 00125 void SetSigma(RealType sigma); 00126 RealType GetSigma(); 00128 00132 void SetNormalizeAcrossScale(bool normalizeInScaleSpace); 00133 itkGetConstMacro(NormalizeAcrossScale, bool); 00135 00136 void SetNumberOfThreads(ThreadIdType nb); 00137 00138 #ifdef ITK_USE_CONCEPT_CHECKING 00139 00140 itkConceptMacro( InputHasNumericTraitsCheck, 00141 ( Concept::HasNumericTraits< PixelType > ) ); 00142 00144 #endif 00145 protected: 00146 GradientMagnitudeRecursiveGaussianImageFilter(); 00147 virtual ~GradientMagnitudeRecursiveGaussianImageFilter() {} 00148 void PrintSelf(std::ostream & os, Indent indent) const; 00150 00152 void GenerateData(void); 00153 00160 virtual void GenerateInputRequestedRegion() 00161 throw( InvalidRequestedRegionError ); 00162 00166 void EnlargeOutputRequestedRegion(DataObject *output); 00167 00168 private: 00169 GradientMagnitudeRecursiveGaussianImageFilter(const Self &); //purposely not 00170 // implemented 00171 void operator=(const Self &); //purposely not 00172 00173 // implemented 00174 00175 class SqrSpacing 00176 { 00177 public: 00178 SqrSpacing():m_Spacing(0) {} 00179 ~SqrSpacing() {} 00180 bool operator!=(const SqrSpacing & other) const 00181 { 00182 return !( *this == other ); 00183 } 00184 00185 bool operator==(const SqrSpacing & other) const 00186 { 00187 return other.m_Spacing == m_Spacing; 00188 } 00189 00190 inline InternalRealType operator()(const InternalRealType & a, const InternalRealType & b) 00191 { 00192 return a + vnl_math_sqr(b / m_Spacing); 00193 } 00194 00195 double m_Spacing; 00196 }; 00197 00198 typedef BinaryFunctorImageFilter< RealImageType, RealImageType, RealImageType, SqrSpacing > SqrSpacingFilterType; 00199 typedef typename SqrSpacingFilterType::Pointer SqrSpacingFilterPointer; 00200 00201 GaussianFilterPointer m_SmoothingFilters[ImageDimension - 1]; 00202 DerivativeFilterPointer m_DerivativeFilter; 00203 SqrSpacingFilterPointer m_SqrSpacingFilter; 00204 SqrtFilterPointer m_SqrtFilter; 00205 00207 bool m_NormalizeAcrossScale; 00208 }; 00209 } // end namespace itk 00210 00211 #ifndef ITK_MANUAL_INSTANTIATION 00212 #include "itkGradientMagnitudeRecursiveGaussianImageFilter.hxx" 00213 #endif 00214 00215 #endif 00216