ITK  4.6.0
Insight Segmentation and Registration Toolkit
itkGradientRecursiveGaussianImageFilter.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 __itkGradientRecursiveGaussianImageFilter_h
19 #define __itkGradientRecursiveGaussianImageFilter_h
20 
23 #include "itkImage.h"
24 #include "itkCovariantVector.h"
26 #include "itkProgressAccumulator.h"
27 #include "itkImageRegionIterator.h"
28 #include "itkVectorImage.h"
29 #include <vector>
30 
31 namespace itk
32 {
51 template< typename TInputImage,
52  typename TOutputImage = Image< CovariantVector<
53  typename NumericTraits< typename TInputImage::PixelType >::RealType,
54  TInputImage::ImageDimension >,
55  TInputImage::ImageDimension > >
57  public ImageToImageFilter< TInputImage, TOutputImage >
58 {
59 public:
65 
67  typedef TInputImage InputImageType;
68  typedef typename TInputImage::PixelType PixelType;
71 
77 
79  itkStaticConstMacro(ImageDimension, unsigned int,
80  TInputImage::ImageDimension);
81 
84 
88  typedef Image< InternalRealType,
89  itkGetStaticConstMacro(ImageDimension) > RealImageType;
90 
91 
96  typedef NthElementImageAdaptor< TOutputImage,
98 
100 
106 
112 
115 
118 
120  typedef typename TOutputImage::Pointer OutputImagePointer;
121 
123  typedef TOutputImage OutputImageType;
124  typedef typename OutputImageType::PixelType OutputPixelType;
128 
130  itkNewMacro(Self);
131 
135 
137  void SetSigma(ScalarRealType sigma);
138  RealType GetSigma() const;
140 
144  void SetNormalizeAcrossScale(bool normalizeInScaleSpace);
145  itkGetConstMacro(NormalizeAcrossScale, bool);
147 
153  virtual void GenerateInputRequestedRegion();
154 
165  itkSetMacro(UseImageDirection, bool);
166  itkGetConstMacro(UseImageDirection, bool);
167  itkBooleanMacro(UseImageDirection);
169 
170 #ifdef ITK_USE_CONCEPT_CHECKING
171  // Begin concept checking
172  itkConceptMacro( InputHasNumericTraitsCheck,
174  // End concept checking
175 #endif
176 
177 protected:
180  void PrintSelf(std::ostream & os, Indent indent) const;
181 
183  void GenerateData(void);
184 
185  // Override since the filter produces the entire dataset
187 
189 
190 private:
191 
192  template <typename TValue>
194  {
195  // To transform Variable length vector we need to convert to and
196  // fro the CovariantVectorType
197  const CovariantVectorType gradient( it.Get().GetDataPointer() );
198  CovariantVectorType physicalGradient;
199  it.GetImage()->TransformLocalVectorToPhysicalVector(gradient, physicalGradient );
200  it.Set( OutputPixelType( physicalGradient.GetDataPointer(), ImageDimension, false ) );
201  }
202 
203  template <typename T >
205  {
206  OutputPixelType correctedGradient;
207  const OutputPixelType & gradient = it.Get();
208 
209  const unsigned int nComponents = NumericTraits<OutputPixelType>::GetLength( gradient )/ImageDimension;
210 
211  for (unsigned int nc = 0; nc < nComponents; nc++ )
212  {
213  GradientVectorType componentGradient;
214  GradientVectorType correctedComponentGradient;
215  for (unsigned int dim = 0; dim < ImageDimension; dim++ )
216  {
217  componentGradient[dim] = DefaultConvertPixelTraits<OutputPixelType>::GetNthComponent( nc*ImageDimension+dim, gradient );
218  }
219  it.GetImage()->TransformLocalVectorToPhysicalVector(componentGradient, correctedComponentGradient );
220  for (unsigned int dim = 0; dim < ImageDimension; dim++ )
221  {
222  DefaultConvertPixelTraits<OutputPixelType>::SetNthComponent( nc*ImageDimension+dim, correctedGradient,
223  correctedComponentGradient[dim] );
224  }
225  }
226  it.Set(correctedGradient);
227  }
228 
229  template <template<typename, unsigned int> class P, class T, unsigned int N>
230  void TransformOutputPixel( ImageRegionIterator< Image< P<T,N>, N > > &it )
231  {
232  const OutputPixelType gradient = it.Get();
233  // This uses the more efficient set by reference method
234  it.GetImage()->TransformLocalVectorToPhysicalVector(gradient, it.Value() );
235  }
236 
237 
238  GradientRecursiveGaussianImageFilter(const Self &); //purposely not
239  // implemented
240  void operator=(const Self &); //purposely not
241 
242  // implemented
243 
244  std::vector< GaussianFilterPointer > m_SmoothingFilters;
247 
250 
253 };
254 } // end namespace itk
255 
256 #ifndef ITK_MANUAL_INSTANTIATION
257 #include "itkGradientRecursiveGaussianImageFilter.hxx"
258 #endif
259 
260 #endif
CovariantVector< OutputComponentType, ImageDimension > CovariantVectorType
RecursiveGaussianImageFilter< RealImageType, RealImageType > GaussianFilterType
static void SetNthComponent(int c, PixelType &pixel, const ComponentType &v)
RecursiveGaussianImageFilter< InputImageType, RealImageType > DerivativeFilterType
void Set(const PixelType &value) const
Base class for computing IIR convolution with an approximation of a Gaussian kernel.
Templated n-dimensional vector image class.
Computes the gradient of an image by convolution with the first derivative of a Gaussian.
Image< InternalRealType, itkGetStaticConstMacro(ImageDimension) > RealImageType
Base class for all process objects that output image data.
CovariantVector< ScalarRealType, ImageDimension > GradientVectorType
static ComponentType GetNthComponent(int c, const PixelType &pixel)
NumericTraits< PixelType >::ScalarRealType ScalarRealType
NumericTraits< OutputPixelType >::ValueType OutputComponentType
void EnlargeOutputRequestedRegion(DataObject *output)
void SetNormalizeAcrossScale(bool normalizeInScaleSpace)
static unsigned int GetLength()
void TransformOutputPixel(ImageRegionIterator< VectorImage< TValue, ImageDimension > > &it)
NumericTraits< InternalRealType >::ValueType InternalScalarRealType
void SetSigma(ScalarRealType sigma)
Presents an image as being composed of the N-th element of its pixels.
void TransformOutputPixel(ImageRegionIterator< Image< P< T, N >, N > > &it)
void PrintSelf(std::ostream &os, Indent indent) const
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
ImageToImageFilter< TInputImage, TOutputImage > Superclass
Define additional traits for native types such as int or float.
#define itkConceptMacro(name, concept)
PixelType Get(void) const
A templated class holding a n-Dimensional covariant vector.
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75
A multi-dimensional iterator templated over image type that walks a region of pixels.
NthElementImageAdaptor< TOutputImage, InternalScalarRealType > OutputImageAdaptorType
const ImageType * GetImage() const