ITK  6.0.0
Insight Toolkit
itkVectorGradientMagnitudeImageFilter.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  * https://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 itkVectorGradientMagnitudeImageFilter_h
19 #define itkVectorGradientMagnitudeImageFilter_h
20 
22 #include "itkImageToImageFilter.h"
23 #include "itkImage.h"
24 #include "itkVector.h"
25 #include "vnl/vnl_matrix.h"
26 #include "vnl/vnl_vector_fixed.h"
27 #include "vnl/algo/vnl_symmetric_eigensystem.h"
28 #include "itkMath.h"
29 
30 namespace itk
31 {
105  * \par Constraints
106  * The filter requires an image with at least two dimensions and a vector
107  * length of at least 2. The theory supports extension to scalar images, but
108  * the implementation of the itk vector classes do not
109  *
110  * The template parameter TRealType must be floating point (float or double) or
111  * a user-defined "real" numerical type with arithmetic operations defined
112  * sufficient to compute derivatives.
113  *
114  * \par Performance
115  * This filter will automatically multithread if run with
116  * SetUsePrincipleComponents=Off or on 3D data in UsePrincipleComponents=On
117  * mode. Unfortunately the ND eigen solver used is not thread safe (a special
118  * 3D solver is), so it cannot multithread for data other than 3D in
120  *
121  * \par References
122  *
123  * [1] G. Sapiro and D. Ringach, "Anisotropic Diffusion of Multivalued Images
124  * with Application to Color Filtering," IEEE Transactions on Image Processing,
125  * Vol 5, No. 11 pp. 1582-1586, 1996
126 
127  * \ingroup GradientFilters
128  *
129  * \sa Image
130  * \sa Neighborhood
133  * \ingroup ITKImageGradient
134  */
135 template <typename TInputImage,
136  typename TRealType = float,
137  typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
138 class ITK_TEMPLATE_EXPORT VectorGradientMagnitudeImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
139 {
140 public:
141  ITK_DISALLOW_COPY_AND_MOVE(VectorGradientMagnitudeImageFilter);
142 
148 
150  itkNewMacro(Self);
151 
153  itkOverrideGetNameOfClassMacro(VectorGradientMagnitudeImageFilter);
154 
157  using OutputPixelType = typename TOutputImage::PixelType;
158  using InputPixelType = typename TInputImage::PixelType;
159 
161  using InputImageType = TInputImage;
162  using OutputImageType = TOutputImage;
165 
168  static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
169  static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
170 
172  static constexpr unsigned int VectorDimension = InputPixelType::Dimension;
173 
175  using RealType = TRealType;
178 
183 
185  using typename Superclass::OutputImageRegionType;
186 
195  void
196  GenerateInputRequestedRegion() override;
197 
203  void
204  SetUseImageSpacing(bool);
205  itkGetConstMacro(UseImageSpacing, bool);
206  itkBooleanMacro(UseImageSpacing);
209 #if !defined(ITK_FUTURE_LEGACY_REMOVE)
210 
215  void
217  {
218  this->SetUseImageSpacing(true);
219  }
220 
225  void
227  {
228  this->SetUseImageSpacing(false);
229  }
230 #endif
231 
234 #if !defined(ITK_LEGACY_REMOVE)
235  using WeightsType [[deprecated("Use DerivativeWeightsType or ComponentWeightsType instead.")]] = ComponentWeightsType;
236 #endif
237 
240  itkSetMacro(DerivativeWeights, DerivativeWeightsType);
241  itkGetConstReferenceMacro(DerivativeWeights, DerivativeWeightsType);
246  itkSetMacro(ComponentWeights, ComponentWeightsType);
247  itkGetConstReferenceMacro(ComponentWeights, ComponentWeightsType);
255  itkSetMacro(UsePrincipleComponents, bool);
256  itkGetConstMacro(UsePrincipleComponents, bool);
257  itkBooleanMacro(UsePrincipleComponents);
260 #if !defined(ITK_FUTURE_LEGACY_REMOVE)
261 
262  void
264  {
265  this->SetUsePrincipleComponents(true);
266  }
267 
269  void
271  {
272  this->SetUsePrincipleComponents(false);
273  }
274 #endif
275 
278  static int
279  CubicSolver(const double *, double *);
280 
281 #ifdef ITK_USE_CONCEPT_CHECKING
282  // Begin concept checking
285  itkConceptMacro(RealTypeHasNumericTraitsCheck, (Concept::HasNumericTraits<RealType>));
286  // End concept checking
287 #endif
288 
289 protected:
291  ~VectorGradientMagnitudeImageFilter() override = default;
292 
296  void
297  BeforeThreadedGenerateData() override;
298 
310  void
311  DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
312 
313 
314  void
315  PrintSelf(std::ostream & os, Indent indent) const override;
316 
318 
320  itkGetConstObjectMacro(RealValuedInputImage, RealVectorImageType);
321 
322  TRealType
324  {
325  unsigned int i;
326  unsigned int j;
327  TRealType dx;
328  TRealType sum;
329  auto accum = TRealType{};
330  for (i = 0; i < ImageDimension; ++i)
331  {
332  sum = TRealType{};
333  for (j = 0; j < VectorDimension; ++j)
334  {
335  dx = m_DerivativeWeights[i] * m_SqrtComponentWeights[j] * 0.5 * (it.GetNext(i)[j] - it.GetPrevious(i)[j]);
336  sum += dx * dx;
337  }
338  accum += sum;
339  }
340  return std::sqrt(accum);
341  }
342 
343  TRealType
345  {
346  // WARNING: ONLY CALL THIS METHOD WHEN PROCESSING A 3D IMAGE
347  unsigned int i;
348  unsigned int j;
349  double Lambda[3];
350  double CharEqn[3];
351  double ans;
352 
353  vnl_matrix<TRealType> g(ImageDimension, ImageDimension);
354  vnl_vector_fixed<TRealType, VectorDimension> d_phi_du[TInputImage::ImageDimension];
355 
356  // Calculate the directional derivatives for each vector component using
357  // central differences.
358  for (i = 0; i < ImageDimension; ++i)
359  {
360  for (j = 0; j < VectorDimension; ++j)
361  {
362  d_phi_du[i][j] =
363  m_DerivativeWeights[i] * m_SqrtComponentWeights[j] * 0.5 * (it.GetNext(i)[j] - it.GetPrevious(i)[j]);
364  }
365  }
366 
367  // Calculate the symmetric metric tensor g
368  for (i = 0; i < ImageDimension; ++i)
369  {
370  for (j = i; j < ImageDimension; ++j)
371  {
372  g[j][i] = g[i][j] = dot_product(d_phi_du[i], d_phi_du[j]);
373  }
374  }
375 
376  // Find the coefficients of the characteristic equation det(g - lambda I)=0
377  // CharEqn[3] = 1.0;
378 
379  CharEqn[2] = -(g[0][0] + g[1][1] + g[2][2]);
380 
381  CharEqn[1] = (g[0][0] * g[1][1] + g[0][0] * g[2][2] + g[1][1] * g[2][2]) -
382  (g[0][1] * g[1][0] + g[0][2] * g[2][0] + g[1][2] * g[2][1]);
383 
384  CharEqn[0] = g[0][0] * (g[1][2] * g[2][1] - g[1][1] * g[2][2]) + g[1][0] * (g[2][2] * g[0][1] - g[0][2] * g[2][1]) +
385  g[2][0] * (g[1][1] * g[0][2] - g[0][1] * g[1][2]);
386 
387  // Find the eigenvalues of g
388  const int numberOfDistinctRoots = this->CubicSolver(CharEqn, Lambda);
389 
390  // Define gradient magnitude as the difference between two largest
391  // eigenvalues. Other definitions may be appropriate here as well.
392  if (numberOfDistinctRoots == 3) // By far the most common case
393  {
394  if (Lambda[0] > Lambda[1])
395  {
396  if (Lambda[1] > Lambda[2]) // Most common, guaranteed?
397  {
398  ans = Lambda[0] - Lambda[1];
399  }
400  else
401  {
402  if (Lambda[0] > Lambda[2])
403  {
404  ans = Lambda[0] - Lambda[2];
405  }
406  else
407  {
408  ans = Lambda[2] - Lambda[0];
409  }
410  }
411  }
412  else
413  {
414  if (Lambda[0] > Lambda[2])
415  {
416  ans = Lambda[1] - Lambda[0];
417  }
418  else
419  {
420  if (Lambda[1] > Lambda[2])
421  {
422  ans = Lambda[1] - Lambda[2];
423  }
424  else
425  {
426  ans = Lambda[2] - Lambda[1];
427  }
428  }
429  }
430  }
431  else if (numberOfDistinctRoots == 2)
432  {
433  if (Lambda[0] > Lambda[1])
434  {
435  ans = Lambda[0] - Lambda[1];
436  }
437  else
438  {
439  ans = Lambda[1] - Lambda[0];
440  }
441  }
442  else if (numberOfDistinctRoots == 1)
443  {
444  ans = 0.0;
445  }
446  else
447  {
448  itkExceptionMacro("Undefined condition. Cubic root solver returned " << numberOfDistinctRoots
449  << " distinct roots.");
450  }
451 
452  return ans;
453  }
454 
455  // Function is defined here because the templating confuses gcc 2.96 when
456  // defined
457  // in .hxx file. jc 1/29/03
458  TRealType
460  {
461  unsigned int i;
462  unsigned int j;
463 
464  vnl_matrix<TRealType> g(ImageDimension, ImageDimension);
465  vnl_vector_fixed<TRealType, VectorDimension> d_phi_du[TInputImage::ImageDimension];
466 
467  // Calculate the directional derivatives for each vector component using
468  // central differences.
469  for (i = 0; i < ImageDimension; ++i)
470  {
471  for (j = 0; j < VectorDimension; ++j)
472  {
473  d_phi_du[i][j] =
474  m_DerivativeWeights[i] * m_SqrtComponentWeights[j] * 0.5 * (it.GetNext(i)[j] - it.GetPrevious(i)[j]);
475  }
476  }
477 
478  // Calculate the symmetric metric tensor g
479  for (i = 0; i < ImageDimension; ++i)
480  {
481  for (j = i; j < ImageDimension; ++j)
482  {
483  g[j][i] = g[i][j] = dot_product(d_phi_du[i], d_phi_du[j]);
484  }
485  }
486 
487  // Find the eigenvalues of g
488  const vnl_symmetric_eigensystem<TRealType> E(g);
489 
490  // Return the difference in length between the first two principle axes.
491  // Note that other edge strength metrics may be appropriate here instead..
492  return (E.get_eigenvalue(ImageDimension - 1) - E.get_eigenvalue(ImageDimension - 2));
493  }
494 
496  DerivativeWeightsType m_DerivativeWeights = DerivativeWeightsType::Filled(1);
497 
501  ComponentWeightsType m_ComponentWeights = ComponentWeightsType::Filled(1);
502  ComponentWeightsType m_SqrtComponentWeights = ComponentWeightsType::Filled(1);
505 private:
506  bool m_UseImageSpacing{ true };
507  bool m_UsePrincipleComponents{};
508 
509  ThreadIdType m_RequestedNumberOfWorkUnits{};
510 
511  typename RealVectorImageType::ConstPointer m_RealValuedInputImage{};
512 };
513 } // end namespace itk
514 
515 #ifndef ITK_MANUAL_INSTANTIATION
516 # include "itkVectorGradientMagnitudeImageFilter.hxx"
517 #endif
518 
519 #endif
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::VectorGradientMagnitudeImageFilter::RadiusType
typename ConstNeighborhoodIteratorType::RadiusType RadiusType
Definition: itkVectorGradientMagnitudeImageFilter.h:182
itk::Concept::HasNumericTraits
Definition: itkConceptChecking.h:717
itk::ImageSource::OutputImagePointer
typename OutputImageType::Pointer OutputImagePointer
Definition: itkImageSource.h:91
itk::UsePrincipleComponents
*par Constraints *The filter requires an image with at least two dimensions and a vector *length of at least The theory supports extension to scalar but *the implementation of the itk vector classes do not **The template parameter TRealType must be floating so it cannot multithread for data other than in * UsePrincipleComponents
Definition: itkVectorGradientMagnitudeImageFilter.h:119
itk::VectorGradientMagnitudeImageFilter
Computes a scalar, gradient magnitude image from a multiple channel (pixels are vectors) input.
Definition: itkVectorGradientMagnitudeImageFilter.h:138
itk::Size
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:69
itkNeighborhoodIterator.h
itk::NeighborhoodOperator
Virtual class that defines a common interface to all neighborhood operator subtypes.
Definition: itkNeighborhoodOperator.h:72
itk::Vector
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:62
itk::ConstNeighborhoodIterator::GetPrevious
PixelType GetPrevious(const unsigned int axis, NeighborIndexType i) const
Definition: itkConstNeighborhoodIterator.h:261
itk::Neighborhood
A light-weight container object for storing an N-dimensional neighborhood of values.
Definition: itkNeighborhood.h:54
itk::VectorGradientMagnitudeImageFilter::InputPixelType
typename TInputImage::PixelType InputPixelType
Definition: itkVectorGradientMagnitudeImageFilter.h:158
itkImage.h
itk::VectorGradientMagnitudeImageFilter::SetUsePrincipleComponentsOff
void SetUsePrincipleComponentsOff()
Definition: itkVectorGradientMagnitudeImageFilter.h:270
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::VectorGradientMagnitudeImageFilter::SetUseImageSpacingOn
void SetUseImageSpacingOn()
Definition: itkVectorGradientMagnitudeImageFilter.h:216
itk::Concept::SameDimension
Definition: itkConceptChecking.h:697
itk::ThreadIdType
unsigned int ThreadIdType
Definition: itkIntTypes.h:102
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
itk::VectorGradientMagnitudeImageFilter::EvaluateAtNeighborhood
TRealType EvaluateAtNeighborhood(const ConstNeighborhoodIteratorType &it) const
Definition: itkVectorGradientMagnitudeImageFilter.h:459
itk::ImageToImageFilter::InputImagePointer
typename InputImageType::Pointer InputImagePointer
Definition: itkImageToImageFilter.h:130
itk::point
*par Constraints *The filter requires an image with at least two dimensions and a vector *length of at least The theory supports extension to scalar but *the implementation of the itk vector classes do not **The template parameter TRealType must be floating point(float or double) or *a user-defined "real" numerical type with arithmetic operations defined *sufficient to compute derivatives. **\par Performance *This filter will automatically multithread if run with *SetUsePrincipleComponents
itk::VectorGradientMagnitudeImageFilter::SetUsePrincipleComponentsOn
void SetUsePrincipleComponentsOn()
Definition: itkVectorGradientMagnitudeImageFilter.h:263
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itkImageToImageFilter.h
itk::VectorGradientMagnitudeImageFilter::EvaluateAtNeighborhood3D
TRealType EvaluateAtNeighborhood3D(const ConstNeighborhoodIteratorType &it) const
Definition: itkVectorGradientMagnitudeImageFilter.h:344
itk::VectorGradientMagnitudeImageFilter::RealType
TRealType RealType
Definition: itkVectorGradientMagnitudeImageFilter.h:175
itk::FixedArray< TRealType, VectorDimension >
itk::Processing
*par Constraints *The filter requires an image with at least two dimensions and a vector *length of at least The theory supports extension to scalar but *the implementation of the itk vector classes do not **The template parameter TRealType must be floating so it cannot multithread for data other than in Anisotropic Diffusion of Multivalued Images *with Application to Color IEEE Transactions on Image Processing
Definition: itkVectorGradientMagnitudeImageFilter.h:119
itk::VectorGradientMagnitudeImageFilter::SetUseImageSpacingOff
void SetUseImageSpacingOff()
Definition: itkVectorGradientMagnitudeImageFilter.h:226
itk::VectorGradientMagnitudeImageFilter::ImageBaseType
typename InputImageType::Superclass ImageBaseType
Definition: itkVectorGradientMagnitudeImageFilter.h:317
itk::NeighborhoodIterator
Defines iteration of a local N-dimensional neighborhood of pixels across an itk::Image.
Definition: itkNeighborhoodIterator.h:212
itk::Vol
*par Constraints *The filter requires an image with at least two dimensions and a vector *length of at least The theory supports extension to scalar but *the implementation of the itk vector classes do not **The template parameter TRealType must be floating so it cannot multithread for data other than in Anisotropic Diffusion of Multivalued Images *with Application to Color IEEE Transactions on Image * Vol
Definition: itkVectorGradientMagnitudeImageFilter.h:119
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:65
itk::VectorGradientMagnitudeImageFilter::OutputPixelType
typename TOutputImage::PixelType OutputPixelType
Definition: itkVectorGradientMagnitudeImageFilter.h:157
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::ConstNeighborhoodIterator
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
Definition: itkConstNeighborhoodIterator.h:51
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:139
itkVector.h
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itk::VectorGradientMagnitudeImageFilter::NonPCEvaluateAtNeighborhood
TRealType NonPCEvaluateAtNeighborhood(const ConstNeighborhoodIteratorType &it) const
Definition: itkVectorGradientMagnitudeImageFilter.h:323
itk::ConstNeighborhoodIterator::GetNext
PixelType GetNext(const unsigned int axis, NeighborIndexType i) const
Definition: itkConstNeighborhoodIterator.h:243
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition: itkGTestTypedefsAndConstructors.h:44
Superclass
BinaryGeneratorImageFilter< TInputImage1, TInputImage2, TOutputImage > Superclass
Definition: itkAddImageFilter.h:90
itkMath.h
itk::images
*par Constraints *The filter requires an image with at least two dimensions and a vector *length of at least The theory supports extension to scalar images
Definition: itkVectorGradientMagnitudeImageFilter.h:107
itk::pp
*par Constraints *The filter requires an image with at least two dimensions and a vector *length of at least The theory supports extension to scalar but *the implementation of the itk vector classes do not **The template parameter TRealType must be floating so it cannot multithread for data other than in Anisotropic Diffusion of Multivalued Images *with Application to Color IEEE Transactions on Image No pp
Definition: itkVectorGradientMagnitudeImageFilter.h:119
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90