ITK  4.1.0
Insight Segmentation and Registration Toolkit
itkGradientMagnitudeImageFilter.h
Go to the documentation of this file.
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 __itkGradientMagnitudeImageFilter_h
00019 #define __itkGradientMagnitudeImageFilter_h
00020 
00021 #include "itkImageToImageFilter.h"
00022 
00023 namespace itk
00024 {
00040 template< typename TInputImage, typename TOutputImage >
00041 class ITK_EXPORT GradientMagnitudeImageFilter:
00042   public ImageToImageFilter< TInputImage, TOutputImage >
00043 {
00044 public:
00046   typedef GradientMagnitudeImageFilter                    Self;
00047   typedef ImageToImageFilter< TInputImage, TOutputImage > Superclass;
00048   typedef SmartPointer< Self >                            Pointer;
00049   typedef SmartPointer< const Self >                      ConstPointer;
00050 
00052   itkNewMacro(Self);
00053 
00055   itkTypeMacro(GradientMagnitudeImageFilter, ImageToImageFilter);
00056 
00059   typedef typename TOutputImage::PixelType                   OutputPixelType;
00060   typedef typename TInputImage::PixelType                    InputPixelType;
00061   typedef typename NumericTraits< InputPixelType >::RealType RealType;
00062 
00065   itkStaticConstMacro(ImageDimension, unsigned int,
00066                       TOutputImage::ImageDimension);
00067 
00069   typedef TInputImage                       InputImageType;
00070   typedef TOutputImage                      OutputImageType;
00071   typedef typename InputImageType::Pointer  InputImagePointer;
00072   typedef typename OutputImageType::Pointer OutputImagePointer;
00073 
00075   typedef typename Superclass::OutputImageRegionType OutputImageRegionType;
00076 
00085   virtual void GenerateInputRequestedRegion()
00086   throw( InvalidRequestedRegionError );
00087 
00090   void SetUseImageSpacingOn()
00091   { this->SetUseImageSpacing(true); }
00092 
00095   void SetUseImageSpacingOff()
00096   { this->SetUseImageSpacing(false); }
00097 
00100   itkSetMacro(UseImageSpacing, bool);
00101   itkGetConstMacro(UseImageSpacing, bool);
00103 
00104 #ifdef ITK_USE_CONCEPT_CHECKING
00105 
00106   itkConceptMacro( InputHasNumericTraitsCheck,
00107                    ( Concept::HasNumericTraits< InputPixelType > ) );
00108 
00110 #endif
00111 protected:
00112   GradientMagnitudeImageFilter()
00113   {
00114     m_UseImageSpacing = true;
00115   }
00116 
00117   virtual ~GradientMagnitudeImageFilter() {}
00118 
00130   void ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread,
00131                             ThreadIdType threadId);
00132 
00133   void PrintSelf(std::ostream &, Indent) const;
00134 private:
00135   GradientMagnitudeImageFilter(const Self &); //purposely not implemented
00136   void operator=(const Self &);               //purposely not implemented
00137 
00138   bool m_UseImageSpacing;
00139 };
00140 } // end namespace itk
00141 
00142 #ifndef ITK_MANUAL_INSTANTIATION
00143 #include "itkGradientMagnitudeImageFilter.hxx"
00144 #endif
00145 
00146 #endif
00147