ITK  4.13.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 > >
56 class ITK_TEMPLATE_EXPORT GradientRecursiveGaussianImageFilter:
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 
102  typedef FixedArray< ScalarRealType,
103  itkGetStaticConstMacro(ImageDimension) > SigmaArrayType;
104 
110 
116 
119 
122 
124  typedef typename TOutputImage::Pointer OutputImagePointer;
125 
127  typedef TOutputImage OutputImageType;
128  typedef typename OutputImageType::PixelType OutputPixelType;
132 
134  itkNewMacro(Self);
135 
139 
141  void SetSigmaArray(const SigmaArrayType & sigmas);
142  void SetSigma(ScalarRealType sigma);
144 
145  SigmaArrayType GetSigmaArray() const;
146  ScalarRealType GetSigma() const;
147 
151  void SetNormalizeAcrossScale(bool normalizeInScaleSpace);
152  itkGetConstMacro(NormalizeAcrossScale, bool);
154 
160  virtual void GenerateInputRequestedRegion() ITK_OVERRIDE;
161 
172  itkSetMacro(UseImageDirection, bool);
173  itkGetConstMacro(UseImageDirection, bool);
174  itkBooleanMacro(UseImageDirection);
176 
177 #ifdef ITK_USE_CONCEPT_CHECKING
178  // Begin concept checking
179  // Does not seem to work with wrappings, disabled
180  // itkConceptMacro( InputHasNumericTraitsCheck,
181  // ( Concept::HasNumericTraits< PixelType > ) );
182  // End concept checking
183 #endif
184 
185 protected:
187  virtual ~GradientRecursiveGaussianImageFilter() ITK_OVERRIDE {}
188  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
189 
191  void GenerateData(void) ITK_OVERRIDE;
192 
193  // Override since the filter produces the entire dataset
194  void EnlargeOutputRequestedRegion(DataObject *output) ITK_OVERRIDE;
195 
196  void GenerateOutputInformation() ITK_OVERRIDE;
197 
198 private:
199 
200  template <typename TValue>
201  void TransformOutputPixel( ImageRegionIterator< VectorImage<TValue, ImageDimension> > &it )
202  {
203  // To transform Variable length vector we need to convert to and
204  // fro the CovariantVectorType
205  const CovariantVectorType gradient( it.Get().GetDataPointer() );
206  CovariantVectorType physicalGradient;
207  it.GetImage()->TransformLocalVectorToPhysicalVector(gradient, physicalGradient );
208  it.Set( OutputPixelType( physicalGradient.GetDataPointer(), ImageDimension, false ) );
209  }
210 
211  template <typename T >
213  {
214  OutputPixelType correctedGradient;
215  const OutputPixelType & gradient = it.Get();
216 
217  const unsigned int nComponents = NumericTraits<OutputPixelType>::GetLength( gradient )/ImageDimension;
218 
219  for (unsigned int nc = 0; nc < nComponents; nc++ )
220  {
221  GradientVectorType componentGradient;
222  GradientVectorType correctedComponentGradient;
223  for (unsigned int dim = 0; dim < ImageDimension; dim++ )
224  {
225  componentGradient[dim] = DefaultConvertPixelTraits<OutputPixelType>::GetNthComponent( nc*ImageDimension+dim, gradient );
226  }
227  it.GetImage()->TransformLocalVectorToPhysicalVector(componentGradient, correctedComponentGradient );
228  for (unsigned int dim = 0; dim < ImageDimension; dim++ )
229  {
230  DefaultConvertPixelTraits<OutputPixelType>::SetNthComponent( nc*ImageDimension+dim, correctedGradient,
231  correctedComponentGradient[dim] );
232  }
233  }
234  it.Set(correctedGradient);
235  }
236 
237  template <template<typename, unsigned int> class P, class T, unsigned int N>
238  void TransformOutputPixel( ImageRegionIterator< Image< P<T,N>, N > > &it )
239  {
240  const OutputPixelType gradient = it.Get();
241  // This uses the more efficient set by reference method
242  it.GetImage()->TransformLocalVectorToPhysicalVector(gradient, it.Value() );
243  }
244 
245 
246  ITK_DISALLOW_COPY_AND_ASSIGN(GradientRecursiveGaussianImageFilter);
247 
248  std::vector< GaussianFilterPointer > m_SmoothingFilters;
251 
254 
257 
260 };
261 } // end namespace itk
262 
263 #ifndef ITK_MANUAL_INSTANTIATION
264 #include "itkGradientRecursiveGaussianImageFilter.hxx"
265 #endif
266 
267 #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
FixedArray< ScalarRealType, itkGetStaticConstMacro(ImageDimension) > SigmaArrayType
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
Simulate a standard C array with copy semnatics.
Definition: itkFixedArray.h:50
static ComponentType GetNthComponent(int c, const PixelType &pixel)
NumericTraits< PixelType >::ScalarRealType ScalarRealType
NumericTraits< OutputPixelType >::ValueType OutputComponentType
static unsigned int GetLength()
NumericTraits< InternalRealType >::ValueType InternalScalarRealType
Presents an image as being composed of the N-th element of its pixels.
void TransformOutputPixel(ImageRegionIterator< Image< P< T, N >, N > > &it)
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.
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