ITK  4.8.0
Insight Segmentation and Registration Toolkit
itkVectorGradientMagnitudeImageFilter.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 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 "vnl/vnl_math.h"
29 
30 namespace itk
31 {
134 template< typename TInputImage,
135  typename TRealType = float,
136  typename TOutputImage = Image< TRealType,
137  TInputImage::ImageDimension >
138  >
140  public ImageToImageFilter< TInputImage, TOutputImage >
141 {
142 public:
148 
150  itkNewMacro(Self);
151 
154 
157  typedef typename TOutputImage::PixelType OutputPixelType;
158  typedef typename TInputImage::PixelType InputPixelType;
159 
161  typedef TInputImage InputImageType;
162  typedef TOutputImage OutputImageType;
163  typedef typename InputImageType::Pointer InputImagePointer;
164  typedef typename OutputImageType::Pointer OutputImagePointer;
165 
167  itkStaticConstMacro(ImageDimension, unsigned int,
168  TOutputImage::ImageDimension);
169 
171  itkStaticConstMacro(VectorDimension, unsigned int,
172  InputPixelType::Dimension);
173 
175  typedef TRealType RealType;
180 
185 
188 
197  virtual void GenerateInputRequestedRegion() ITK_OVERRIDE;
198 
204  { this->SetUseImageSpacing(true); }
205 
210  { this->SetUseImageSpacing(false); }
211 
214  void SetUseImageSpacing(bool);
215 
216  itkGetConstMacro(UseImageSpacing, bool);
217 
219 
222  itkSetMacro(DerivativeWeights, WeightsType);
223  itkGetConstReferenceMacro(DerivativeWeights, WeightsType);
225 
228  itkSetMacro(ComponentWeights, WeightsType);
229  itkGetConstReferenceMacro(ComponentWeights, WeightsType);
231 
233  itkGetConstReferenceMacro(NeighborhoodRadius, RadiusType);
234  itkSetMacro(NeighborhoodRadius, RadiusType);
236 
242  itkSetMacro(UsePrincipleComponents, bool);
243  itkGetConstMacro(UsePrincipleComponents, bool);
245  {
246  this->SetUsePrincipleComponents(true);
247  }
249 
251  {
252  this->SetUsePrincipleComponents(false);
253  }
254 
257  static int CubicSolver(double *, double *);
258 
259 #ifdef ITK_USE_CONCEPT_CHECKING
260  // Begin concept checking
261  itkConceptMacro( InputHasNumericTraitsCheck,
263  itkConceptMacro( RealTypeHasNumericTraitsCheck,
265  // End concept checking
266 #endif
267 
268 protected:
271 
275  void BeforeThreadedGenerateData() ITK_OVERRIDE;
276 
288  void ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread,
289  ThreadIdType threadId) ITK_OVERRIDE;
290 
291  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
292 
294 
296  itkGetConstObjectMacro(RealValuedInputImage, ImageBaseType);
297 
299  {
300  unsigned i, j;
301  TRealType dx, sum, accum;
302 
304  for ( i = 0; i < ImageDimension; ++i )
305  {
307  for ( j = 0; j < VectorDimension; ++j )
308  {
310  * 0.5 * ( it.GetNext(i)[j] - it.GetPrevious(i)[j] );
311  sum += dx * dx;
312  }
313  accum += sum;
314  }
315  return std::sqrt(accum);
316  }
317 
319  {
320  // WARNING: ONLY CALL THIS METHOD WHEN PROCESSING A 3D IMAGE
321  unsigned int i, j;
322  double Lambda[3];
323  double CharEqn[3];
324  double ans;
325 
326  vnl_matrix< TRealType > g(ImageDimension, ImageDimension);
327  vnl_vector_fixed< TRealType, VectorDimension >
328  d_phi_du[TInputImage::ImageDimension];
329 
330  // Calculate the directional derivatives for each vector component using
331  // central differences.
332  for ( i = 0; i < ImageDimension; i++ )
333  {
334  for ( j = 0; j < VectorDimension; j++ )
335  {
336  d_phi_du[i][j] = m_DerivativeWeights[i] * m_SqrtComponentWeights[j]
337  * 0.5 * ( it.GetNext(i)[j] - it.GetPrevious(i)[j] );
338  }
339  }
340 
341  // Calculate the symmetric metric tensor g
342  for ( i = 0; i < ImageDimension; i++ )
343  {
344  for ( j = i; j < ImageDimension; j++ )
345  {
346  g[j][i] = g[i][j] = dot_product(d_phi_du[i], d_phi_du[j]);
347  }
348  }
349 
350  // Find the coefficients of the characteristic equation det(g - lambda I)=0
351  // CharEqn[3] = 1.0;
352 
353  CharEqn[2] = -( g[0][0] + g[1][1] + g[2][2] );
354 
355  CharEqn[1] = ( g[0][0] * g[1][1] + g[0][0] * g[2][2] + g[1][1] * g[2][2] )
356  - ( g[0][1] * g[1][0] + g[0][2] * g[2][0] + g[1][2] * g[2][1] );
357 
358  CharEqn[0] = g[0][0] * ( g[1][2] * g[2][1] - g[1][1] * g[2][2] )
359  + g[1][0] * ( g[2][2] * g[0][1] - g[0][2] * g[2][1] )
360  + g[2][0] * ( g[1][1] * g[0][2] - g[0][1] * g[1][2] );
361 
362  // Find the eigenvalues of g
363  int numberOfDistinctRoots = this->CubicSolver(CharEqn, Lambda);
364 
365  // Define gradient magnitude as the difference between two largest
366  // eigenvalues. Other definitions may be appropriate here as well.
367  if ( numberOfDistinctRoots == 3 ) // By far the most common case
368  {
369  if ( Lambda[0] > Lambda[1] )
370  {
371  if ( Lambda[1] > Lambda[2] ) // Most common, guaranteed?
372  {
373  ans = Lambda[0] - Lambda[1];
374  }
375  else
376  {
377  if ( Lambda[0] > Lambda[2] )
378  {
379  ans = Lambda[0] - Lambda[2];
380  }
381  else
382  {
383  ans = Lambda[2] - Lambda[0];
384  }
385  }
386  }
387  else
388  {
389  if ( Lambda[0] > Lambda[2] )
390  {
391  ans = Lambda[1] - Lambda[0];
392  }
393  else
394  {
395  if ( Lambda[1] > Lambda[2] )
396  {
397  ans = Lambda[1] - Lambda[2];
398  }
399  else
400  {
401  ans = Lambda[2] - Lambda[1];
402  }
403  }
404  }
405  }
406  else if ( numberOfDistinctRoots == 2 )
407  {
408  if ( Lambda[0] > Lambda[1] )
409  {
410  ans = Lambda[0] - Lambda[1];
411  }
412  else
413  {
414  ans = Lambda[1] - Lambda[0];
415  }
416  }
417  else if ( numberOfDistinctRoots == 1 )
418  {
419  ans = 0.0;
420  }
421  else
422  {
423  itkExceptionMacro(<< "Undefined condition. Cubic root solver returned "
424  << numberOfDistinctRoots << " distinct roots.");
425  }
426 
427  return ans;
428  }
429 
430  // Function is defined here because the templating confuses gcc 2.96 when
431  // defined
432  // in .hxx file. jc 1/29/03
434  {
435  unsigned int i, j;
436 
437  vnl_matrix< TRealType > g(ImageDimension, ImageDimension);
438  vnl_vector_fixed< TRealType, VectorDimension >
439  d_phi_du[TInputImage::ImageDimension];
440 
441  // Calculate the directional derivatives for each vector component using
442  // central differences.
443  for ( i = 0; i < ImageDimension; i++ )
444  {
445  for ( j = 0; j < VectorDimension; j++ )
446  {
447  d_phi_du[i][j] = m_DerivativeWeights[i] * m_SqrtComponentWeights[j]
448  * 0.5 * ( it.GetNext(i)[j] - it.GetPrevious(i)[j] );
449  }
450  }
451 
452  // Calculate the symmetric metric tensor g
453  for ( i = 0; i < ImageDimension; i++ )
454  {
455  for ( j = i; j < ImageDimension; j++ )
456  {
457  g[j][i] = g[i][j] = dot_product(d_phi_du[i], d_phi_du[j]);
458  }
459  }
460 
461  // Find the eigenvalues of g
462  vnl_symmetric_eigensystem< TRealType > E(g);
463 
464  // Return the difference in length between the first two principle axes.
465  // Note that other edge strength metrics may be appropriate here instead..
466  return ( E.get_eigenvalue(ImageDimension - 1) - E.get_eigenvalue(ImageDimension - 2) );
467  }
468 
471 
477 
478 private:
481 
483 
484  typename ImageBaseType::ConstPointer m_RealValuedInputImage;
485 
486  VectorGradientMagnitudeImageFilter(const Self &); //purposely not implemented
487  void operator=(const Self &); //purposely not implemented
488 
490 };
491 } // end namespace itk
492 
493 #ifndef ITK_MANUAL_INSTANTIATION
494 #include "itkVectorGradientMagnitudeImageFilter.hxx"
495 #endif
496 
497 #endif
ConstNeighborhoodIterator< RealVectorImageType > ConstNeighborhoodIteratorType
void PrintSelf(std::ostream &os, Indent indent) const override
FixedArray< TRealType, VectorDimension > WeightsType
virtual void GenerateInputRequestedRegion() override
TRealType EvaluateAtNeighborhood3D(const ConstNeighborhoodIteratorType &it) const
virtual void SetUsePrincipleComponents(bool _arg)
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
ImageToImageFilter< TInputImage, TOutputImage > Superclass
Base class for all process objects that output image data.
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
ConstNeighborhoodIteratorType::RadiusType RadiusType
virtual PixelType GetPrevious(const unsigned axis, NeighborIndexType i) const
Computes a scalar, gradient magnitude image from a multiple channel (pixels are vectors) input...
virtual PixelType GetNext(const unsigned axis, NeighborIndexType i) const
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId) override
typedef(Concept::HasNumericTraits< typename InputPixelType::ValueType >) InputHasNumericTraitsCheck
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:62
static int CubicSolver(double *, double *)
Vector< TRealType, InputPixelType::Dimension > RealVectorType
TRealType NonPCEvaluateAtNeighborhood(const ConstNeighborhoodIteratorType &it) const
unsigned int ThreadIdType
Definition: itkIntTypes.h:159
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
Superclass::OutputImageRegionType OutputImageRegionType
#define itkConceptMacro(name, concept)
Image< RealVectorType, TInputImage::ImageDimension > RealVectorImageType
Templated n-dimensional image class.
Definition: itkImage.h:75
TRealType EvaluateAtNeighborhood(const ConstNeighborhoodIteratorType &it) const