ITK  5.2.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  * 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 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 {
104 
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 
154 
157  using OutputPixelType = typename TOutputImage::PixelType;
158  using InputPixelType = typename TInputImage::PixelType;
159 
161  using InputImageType = TInputImage;
162  using OutputImageType = TOutputImage;
163  using InputImagePointer = typename InputImageType::Pointer;
164  using OutputImagePointer = typename OutputImageType::Pointer;
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 OutputImageRegionType = typename Superclass::OutputImageRegionType;
186 
195  void
196  GenerateInputRequestedRegion() override;
197 
203  void
204  SetUseImageSpacing(bool);
205  itkGetConstMacro(UseImageSpacing, bool);
206  itkBooleanMacro(UseImageSpacing);
208 
209 #if !defined(ITK_FUTURE_LEGACY_REMOVE)
210 
215  void
216  SetUseImageSpacingOn()
217  {
218  this->SetUseImageSpacing(true);
219  }
220 
225  void
226  SetUseImageSpacingOff()
227  {
228  this->SetUseImageSpacing(false);
229  }
230 #endif
231 
234 #if !defined(ITK_LEGACY_REMOVE)
235  using WeightsType ITK_DEPRECATED_MSG("Use DerivativeWeightsType or ComponentWeightsType instead.") =
237 #endif
238 
241  itkSetMacro(DerivativeWeights, DerivativeWeightsType);
242  itkGetConstReferenceMacro(DerivativeWeights, DerivativeWeightsType);
244 
247  itkSetMacro(ComponentWeights, ComponentWeightsType);
248  itkGetConstReferenceMacro(ComponentWeights, ComponentWeightsType);
250 
256  itkSetMacro(UsePrincipleComponents, bool);
257  itkGetConstMacro(UsePrincipleComponents, bool);
258  itkBooleanMacro(UsePrincipleComponents);
260 
261 #if !defined(ITK_FUTURE_LEGACY_REMOVE)
262 
263  void
264  SetUsePrincipleComponentsOn()
265  {
266  this->SetUsePrincipleComponents(true);
267  }
268 
270  void
271  SetUsePrincipleComponentsOff()
272  {
273  this->SetUsePrincipleComponents(false);
274  }
275 #endif
276 
279  static int
280  CubicSolver(const double *, double *);
281 
282 #ifdef ITK_USE_CONCEPT_CHECKING
283  // Begin concept checking
284  itkConceptMacro(SameDimensionCheck, (Concept::SameDimension<InputImageDimension, ImageDimension>));
285  itkConceptMacro(InputHasNumericTraitsCheck, (Concept::HasNumericTraits<typename InputPixelType::ValueType>));
286  itkConceptMacro(RealTypeHasNumericTraitsCheck, (Concept::HasNumericTraits<RealType>));
287  // End concept checking
288 #endif
289 
290 protected:
291  VectorGradientMagnitudeImageFilter();
292  ~VectorGradientMagnitudeImageFilter() override = default;
293 
297  void
298  BeforeThreadedGenerateData() override;
299 
311  void
312  DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
313 
314 
315  void
316  PrintSelf(std::ostream & os, Indent indent) const override;
317 
318  using ImageBaseType = typename InputImageType::Superclass;
319 
321  itkGetConstObjectMacro(RealValuedInputImage, RealVectorImageType);
322 
323  TRealType
325  {
326  unsigned i, j;
327  TRealType dx, sum, accum;
328 
330  for (i = 0; i < ImageDimension; ++i)
331  {
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, j;
348  double Lambda[3];
349  double CharEqn[3];
350  double ans;
351 
352  vnl_matrix<TRealType> g(ImageDimension, ImageDimension);
353  vnl_vector_fixed<TRealType, VectorDimension> d_phi_du[TInputImage::ImageDimension];
354 
355  // Calculate the directional derivatives for each vector component using
356  // central differences.
357  for (i = 0; i < ImageDimension; i++)
358  {
359  for (j = 0; j < VectorDimension; j++)
360  {
361  d_phi_du[i][j] =
362  m_DerivativeWeights[i] * m_SqrtComponentWeights[j] * 0.5 * (it.GetNext(i)[j] - it.GetPrevious(i)[j]);
363  }
364  }
365 
366  // Calculate the symmetric metric tensor g
367  for (i = 0; i < ImageDimension; i++)
368  {
369  for (j = i; j < ImageDimension; j++)
370  {
371  g[j][i] = g[i][j] = dot_product(d_phi_du[i], d_phi_du[j]);
372  }
373  }
374 
375  // Find the coefficients of the characteristic equation det(g - lambda I)=0
376  // CharEqn[3] = 1.0;
377 
378  CharEqn[2] = -(g[0][0] + g[1][1] + g[2][2]);
379 
380  CharEqn[1] = (g[0][0] * g[1][1] + g[0][0] * g[2][2] + g[1][1] * g[2][2]) -
381  (g[0][1] * g[1][0] + g[0][2] * g[2][0] + g[1][2] * g[2][1]);
382 
383  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]) +
384  g[2][0] * (g[1][1] * g[0][2] - g[0][1] * g[1][2]);
385 
386  // Find the eigenvalues of g
387  int numberOfDistinctRoots = this->CubicSolver(CharEqn, Lambda);
388 
389  // Define gradient magnitude as the difference between two largest
390  // eigenvalues. Other definitions may be appropriate here as well.
391  if (numberOfDistinctRoots == 3) // By far the most common case
392  {
393  if (Lambda[0] > Lambda[1])
394  {
395  if (Lambda[1] > Lambda[2]) // Most common, guaranteed?
396  {
397  ans = Lambda[0] - Lambda[1];
398  }
399  else
400  {
401  if (Lambda[0] > Lambda[2])
402  {
403  ans = Lambda[0] - Lambda[2];
404  }
405  else
406  {
407  ans = Lambda[2] - Lambda[0];
408  }
409  }
410  }
411  else
412  {
413  if (Lambda[0] > Lambda[2])
414  {
415  ans = Lambda[1] - Lambda[0];
416  }
417  else
418  {
419  if (Lambda[1] > Lambda[2])
420  {
421  ans = Lambda[1] - Lambda[2];
422  }
423  else
424  {
425  ans = Lambda[2] - Lambda[1];
426  }
427  }
428  }
429  }
430  else if (numberOfDistinctRoots == 2)
431  {
432  if (Lambda[0] > Lambda[1])
433  {
434  ans = Lambda[0] - Lambda[1];
435  }
436  else
437  {
438  ans = Lambda[1] - Lambda[0];
439  }
440  }
441  else if (numberOfDistinctRoots == 1)
442  {
443  ans = 0.0;
444  }
445  else
446  {
447  itkExceptionMacro(<< "Undefined condition. Cubic root solver returned " << numberOfDistinctRoots
448  << " distinct roots.");
449  }
450 
451  return ans;
452  }
453 
454  // Function is defined here because the templating confuses gcc 2.96 when
455  // defined
456  // in .hxx file. jc 1/29/03
457  TRealType
459  {
460  unsigned int i, j;
461 
462  vnl_matrix<TRealType> g(ImageDimension, ImageDimension);
463  vnl_vector_fixed<TRealType, VectorDimension> d_phi_du[TInputImage::ImageDimension];
464 
465  // Calculate the directional derivatives for each vector component using
466  // central differences.
467  for (i = 0; i < ImageDimension; i++)
468  {
469  for (j = 0; j < VectorDimension; j++)
470  {
471  d_phi_du[i][j] =
472  m_DerivativeWeights[i] * m_SqrtComponentWeights[j] * 0.5 * (it.GetNext(i)[j] - it.GetPrevious(i)[j]);
473  }
474  }
475 
476  // Calculate the symmetric metric tensor g
477  for (i = 0; i < ImageDimension; i++)
478  {
479  for (j = i; j < ImageDimension; j++)
480  {
481  g[j][i] = g[i][j] = dot_product(d_phi_du[i], d_phi_du[j]);
482  }
483  }
484 
485  // Find the eigenvalues of g
486  vnl_symmetric_eigensystem<TRealType> E(g);
487 
488  // Return the difference in length between the first two principle axes.
489  // Note that other edge strength metrics may be appropriate here instead..
490  return (E.get_eigenvalue(ImageDimension - 1) - E.get_eigenvalue(ImageDimension - 2));
491  }
492 
494  DerivativeWeightsType m_DerivativeWeights = DerivativeWeightsType::Filled(1);
495 
499  ComponentWeightsType m_ComponentWeights = ComponentWeightsType::Filled(1);
500  ComponentWeightsType m_SqrtComponentWeights = ComponentWeightsType::Filled(1);
502 
503 private:
506 
508 
510 };
511 } // end namespace itk
512 
513 #ifndef ITK_MANUAL_INSTANTIATION
514 # include "itkVectorGradientMagnitudeImageFilter.hxx"
515 #endif
516 
517 #endif
itk::VectorGradientMagnitudeImageFilter::RadiusType
typename ConstNeighborhoodIteratorType::RadiusType RadiusType
Definition: itkVectorGradientMagnitudeImageFilter.h:182
itk::VectorGradientMagnitudeImageFilter::m_RequestedNumberOfThreads
ThreadIdType m_RequestedNumberOfThreads
Definition: itkVectorGradientMagnitudeImageFilter.h:507
itk::ConstNeighborhoodIterator::RadiusType
typename Superclass::RadiusType RadiusType
Definition: itkConstNeighborhoodIterator.h:71
itk::ImageSource::OutputImagePointer
typename OutputImageType::Pointer OutputImagePointer
Definition: itkImageSource.h:91
itk::VectorGradientMagnitudeImageFilter
Computes a scalar, gradient magnitude image from a multiple channel (pixels are vectors) input.
Definition: itkVectorGradientMagnitudeImageFilter.h:138
itkNeighborhoodIterator.h
itk::Vector
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:62
itk::VectorGradientMagnitudeImageFilter::InputPixelType
typename TInputImage::PixelType InputPixelType
Definition: itkVectorGradientMagnitudeImageFilter.h:158
itkImage.h
itk::SmartPointer< Self >
itk::ConstNeighborhoodIterator::GetPrevious
ITK_ITERATOR_VIRTUAL PixelType GetPrevious(const unsigned axis, NeighborIndexType i) const ITK_ITERATOR_FINAL
Definition: itkConstNeighborhoodIterator.h:261
itk::VectorGradientMagnitudeImageFilter::m_UsePrincipleComponents
bool m_UsePrincipleComponents
Definition: itkVectorGradientMagnitudeImageFilter.h:505
itk::ThreadIdType
unsigned int ThreadIdType
Definition: itkIntTypes.h:99
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:458
itk::ImageToImageFilter::InputImagePointer
typename InputImageType::Pointer InputImagePointer
Definition: itkImageToImageFilter.h:130
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::VectorGradientMagnitudeImageFilter::ImageBaseType
typename InputImageType::Superclass ImageBaseType
Definition: itkVectorGradientMagnitudeImageFilter.h:318
itk::VectorGradientMagnitudeImageFilter::m_UseImageSpacing
bool m_UseImageSpacing
Definition: itkVectorGradientMagnitudeImageFilter.h:504
itk::ImageSource::OutputImageRegionType
typename OutputImageType::RegionType OutputImageRegionType
Definition: itkImageSource.h:92
itk::NumericTraits::ZeroValue
static T ZeroValue()
Definition: itkNumericTraits.h:148
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: itkAnnulusOperator.h:24
itk::ConstNeighborhoodIterator
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
Definition: itkConstNeighborhoodIterator.h:51
itk::VectorGradientMagnitudeImageFilter::m_RealValuedInputImage
RealVectorImageType::ConstPointer m_RealValuedInputImage
Definition: itkVectorGradientMagnitudeImageFilter.h:509
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:138
itkVector.h
itk::ConstNeighborhoodIterator::GetNext
ITK_ITERATOR_VIRTUAL PixelType GetNext(const unsigned axis, NeighborIndexType i) const ITK_ITERATOR_FINAL
Definition: itkConstNeighborhoodIterator.h:243
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86
itk::VectorGradientMagnitudeImageFilter::NonPCEvaluateAtNeighborhood
TRealType NonPCEvaluateAtNeighborhood(const ConstNeighborhoodIteratorType &it) const
Definition: itkVectorGradientMagnitudeImageFilter.h:324
itk::GTest::TypedefsAndConstructors::Dimension2::Dimension
constexpr unsigned int Dimension
Definition: itkGTestTypedefsAndConstructors.h:44
itkMath.h
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90