ITK  5.2.0
Insight Toolkit
itkGradientRecursiveGaussianImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
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 {
54 template <
55  typename TInputImage,
56  typename TOutputImage = Image<
57  CovariantVector<typename NumericTraits<typename TInputImage::PixelType>::RealType, TInputImage::ImageDimension>,
58  TInputImage::ImageDimension>>
59 class ITK_TEMPLATE_EXPORT GradientRecursiveGaussianImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
60 {
61 public:
62  ITK_DISALLOW_COPY_AND_MOVE(GradientRecursiveGaussianImageFilter);
63 
69 
71  using InputImageType = TInputImage;
72  using PixelType = typename TInputImage::PixelType;
75 
81 
83  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
84 
87 
92 
93 
99 
101 
104 
107 
110 
113 
116 
118  using OutputImagePointer = typename TOutputImage::Pointer;
119 
121  using OutputImageType = TOutputImage;
122  using OutputPixelType = typename OutputImageType::PixelType;
125 
127  itkNewMacro(Self);
128 
131 
133  void
134  SetSigmaArray(const SigmaArrayType & sigma);
135  void
136  SetSigma(ScalarRealType sigma);
138 
140  GetSigmaArray() const;
142  GetSigma() const;
143 
147  void
148  SetNormalizeAcrossScale(bool normalize);
149  itkGetConstMacro(NormalizeAcrossScale, bool);
151 
157  void
158  GenerateInputRequestedRegion() override;
159 
170  itkSetMacro(UseImageDirection, bool);
171  itkGetConstMacro(UseImageDirection, bool);
172  itkBooleanMacro(UseImageDirection);
174 
175 #ifdef ITK_USE_CONCEPT_CHECKING
176  // Begin concept checking
177  // Does not seem to work with wrappings, disabled
178  // itkConceptMacro( InputHasNumericTraitsCheck,
179  // ( Concept::HasNumericTraits< PixelType > ) );
180  // End concept checking
181 #endif
182 
183 protected:
185  ~GradientRecursiveGaussianImageFilter() override = default;
186  void
187  PrintSelf(std::ostream & os, Indent indent) const override;
188 
190  void
191  GenerateData() override;
192 
193  // Override since the filter produces the entire dataset
194  void
195  EnlargeOutputRequestedRegion(DataObject * output) override;
196 
197  void
198  GenerateOutputInformation() override;
199 
200 private:
201  template <typename TValue>
202  void
204  {
205  // To transform Variable length vector we need to convert to and
206  // fro the CovariantVectorType
207  const CovariantVectorType gradient(it.Get().GetDataPointer());
208  CovariantVectorType physicalGradient;
209  it.GetImage()->TransformLocalVectorToPhysicalVector(gradient, physicalGradient);
210  it.Set(OutputPixelType(physicalGradient.GetDataPointer(), ImageDimension, false));
211  }
212 
213  template <typename T>
214  void
216  {
217  OutputPixelType correctedGradient;
218  const OutputPixelType & gradient = it.Get();
219 
220  const unsigned int nComponents = NumericTraits<OutputPixelType>::GetLength(gradient) / ImageDimension;
221 
222  for (unsigned int nc = 0; nc < nComponents; nc++)
223  {
224  GradientVectorType componentGradient;
225  GradientVectorType correctedComponentGradient;
226  for (unsigned int dim = 0; dim < ImageDimension; dim++)
227  {
228  componentGradient[dim] =
229  DefaultConvertPixelTraits<OutputPixelType>::GetNthComponent(nc * ImageDimension + dim, gradient);
230  }
231  it.GetImage()->TransformLocalVectorToPhysicalVector(componentGradient, correctedComponentGradient);
232  for (unsigned int dim = 0; dim < ImageDimension; dim++)
233  {
235  nc * ImageDimension + dim, correctedGradient, correctedComponentGradient[dim]);
236  }
237  }
238  it.Set(correctedGradient);
239  }
240 
241  template <template <typename, unsigned int> class P, class T, unsigned int N>
242  void
244  {
245  const OutputPixelType gradient = it.Get();
246  // This uses the more efficient set by reference method
247  it.GetImage()->TransformLocalVectorToPhysicalVector(gradient, it.Value());
248  }
249 
250 
251  std::vector<GaussianFilterPointer> m_SmoothingFilters;
254 
257 
260 
263 };
264 } // end namespace itk
265 
266 #ifndef ITK_MANUAL_INSTANTIATION
267 # include "itkGradientRecursiveGaussianImageFilter.hxx"
268 #endif
269 
270 #endif
itk::NthElementImageAdaptor
Presents an image as being composed of the N-th element of its pixels.
Definition: itkNthElementImageAdaptor.h:58
itkRecursiveGaussianImageFilter.h
itkCovariantVector.h
itk::GradientRecursiveGaussianImageFilter::m_DerivativeFilter
DerivativeFilterPointer m_DerivativeFilter
Definition: itkGradientRecursiveGaussianImageFilter.h:252
itk::ImageSource::OutputImagePointer
typename OutputImageType::Pointer OutputImagePointer
Definition: itkImageSource.h:91
itk::NumericTraits::ValueType
T ValueType
Definition: itkNumericTraits.h:65
itkNthElementImageAdaptor.h
itk::GradientRecursiveGaussianImageFilter::InternalScalarRealType
typename NumericTraits< InternalRealType >::ValueType InternalScalarRealType
Definition: itkGradientRecursiveGaussianImageFilter.h:80
itk::VectorImage
Templated n-dimensional vector image class.
Definition: itkImageAlgorithm.h:29
itk::ImageConstIterator::GetImage
const ImageType * GetImage() const
Definition: itkImageConstIterator.h:336
itk::GradientRecursiveGaussianImageFilter::TransformOutputPixel
void TransformOutputPixel(ImageRegionIterator< Image< P< T, N >, N >> &it)
Definition: itkGradientRecursiveGaussianImageFilter.h:243
itkImage.h
itk::GradientRecursiveGaussianImageFilter::m_ImageAdaptor
OutputImageAdaptorPointer m_ImageAdaptor
Definition: itkGradientRecursiveGaussianImageFilter.h:253
itk::GradientRecursiveGaussianImageFilter::TransformOutputPixel
void TransformOutputPixel(ImageRegionIterator< T > &it)
Definition: itkGradientRecursiveGaussianImageFilter.h:215
itk::SmartPointer< Self >
itkImageRegionIterator.h
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::GradientRecursiveGaussianImageFilter::TransformOutputPixel
void TransformOutputPixel(ImageRegionIterator< VectorImage< TValue, ImageDimension >> &it)
Definition: itkGradientRecursiveGaussianImageFilter.h:203
itkVectorImage.h
itk::ImageRegionIterator
A multi-dimensional iterator templated over image type that walks a region of pixels.
Definition: itkImageRegionIterator.h:80
itk::GradientRecursiveGaussianImageFilter::m_NormalizeAcrossScale
bool m_NormalizeAcrossScale
Definition: itkGradientRecursiveGaussianImageFilter.h:256
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itkDefaultConvertPixelTraits.h
itk::GradientRecursiveGaussianImageFilter::OutputComponentType
typename NumericTraits< OutputPixelType >::ValueType OutputComponentType
Definition: itkGradientRecursiveGaussianImageFilter.h:123
itk::GradientRecursiveGaussianImageFilter::OutputPixelType
typename OutputImageType::PixelType OutputPixelType
Definition: itkGradientRecursiveGaussianImageFilter.h:122
itk::GradientRecursiveGaussianImageFilter::OutputImageAdaptorPointer
typename OutputImageAdaptorType::Pointer OutputImageAdaptorPointer
Definition: itkGradientRecursiveGaussianImageFilter.h:100
itk::GradientRecursiveGaussianImageFilter::m_SmoothingFilters
std::vector< GaussianFilterPointer > m_SmoothingFilters
Definition: itkGradientRecursiveGaussianImageFilter.h:251
itk::DefaultConvertPixelTraits::SetNthComponent
static void SetNthComponent(int c, PixelType &pixel, const ComponentType &v)
Definition: itkDefaultConvertPixelTraits.h:69
itk::GradientRecursiveGaussianImageFilter
Computes the gradient of an image by convolution with the first derivative of a Gaussian.
Definition: itkGradientRecursiveGaussianImageFilter.h:59
itk::GradientRecursiveGaussianImageFilter::ScalarRealType
typename NumericTraits< PixelType >::ScalarRealType ScalarRealType
Definition: itkGradientRecursiveGaussianImageFilter.h:74
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itk::GradientRecursiveGaussianImageFilter::m_UseImageDirection
bool m_UseImageDirection
Definition: itkGradientRecursiveGaussianImageFilter.h:259
itk::ImageRegionIterator::Set
void Set(const PixelType &value) const
Definition: itkImageRegionIterator.h:116
itk::FixedArray< ScalarRealType, Self::ImageDimension >
itk::GradientRecursiveGaussianImageFilter::GaussianFilterPointer
typename GaussianFilterType::Pointer GaussianFilterPointer
Definition: itkGradientRecursiveGaussianImageFilter.h:112
itk::NumericTraits
Define additional traits for native types such as int or float.
Definition: itkNumericTraits.h:58
itk::CovariantVector
A templated class holding a n-Dimensional covariant vector.
Definition: itkCovariantVector.h:70
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ImageConstIterator::Get
PixelType Get() const
Definition: itkImageConstIterator.h:343
itk::NumericTraits::GetLength
static unsigned int GetLength()
Definition: itkNumericTraits.h:208
itkProgressAccumulator.h
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:138
itk::DefaultConvertPixelTraits::GetNthComponent
static ComponentType GetNthComponent(int c, const PixelType &pixel)
Definition: itkDefaultConvertPixelTraits.h:62
itk::GradientRecursiveGaussianImageFilter::PixelType
typename TInputImage::PixelType PixelType
Definition: itkGradientRecursiveGaussianImageFilter.h:72
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86
itk::RecursiveGaussianImageFilter
Base class for computing IIR convolution with an approximation of a Gaussian kernel.
Definition: itkRecursiveGaussianImageFilter.h:100
itk::GradientRecursiveGaussianImageFilter::m_Sigma
SigmaArrayType m_Sigma
Definition: itkGradientRecursiveGaussianImageFilter.h:262
itk::GradientRecursiveGaussianImageFilter::RealType
typename NumericTraits< PixelType >::RealType RealType
Definition: itkGradientRecursiveGaussianImageFilter.h:73
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293
itk::GradientRecursiveGaussianImageFilter::InternalRealType
typename NumericTraits< RealType >::FloatType InternalRealType
Definition: itkGradientRecursiveGaussianImageFilter.h:79
itk::GradientRecursiveGaussianImageFilter::DerivativeFilterPointer
typename DerivativeFilterType::Pointer DerivativeFilterPointer
Definition: itkGradientRecursiveGaussianImageFilter.h:115