ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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
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
30namespace itk
31{
130template <typename TInputImage,
131 typename TRealType = float,
132 typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
133class ITK_TEMPLATE_EXPORT VectorGradientMagnitudeImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
134{
135public:
136 ITK_DISALLOW_COPY_AND_MOVE(VectorGradientMagnitudeImageFilter);
137
143
145 itkNewMacro(Self);
146
148 itkOverrideGetNameOfClassMacro(VectorGradientMagnitudeImageFilter);
149
152 using OutputPixelType = typename TOutputImage::PixelType;
153 using InputPixelType = typename TInputImage::PixelType;
154
156 using InputImageType = TInputImage;
157 using OutputImageType = TOutputImage;
158 using InputImagePointer = typename InputImageType::Pointer;
159 using OutputImagePointer = typename OutputImageType::Pointer;
160
163 static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
164 static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
165
167 static constexpr unsigned int VectorDimension = InputPixelType::Dimension;
168
170 using RealType = TRealType;
173
178
181
190 void
192
199 void
201 itkGetConstMacro(UseImageSpacing, bool);
202 itkBooleanMacro(UseImageSpacing);
204#if !defined(ITK_FUTURE_LEGACY_REMOVE)
210 void
211 SetUseImageSpacingOn()
212 {
213 this->SetUseImageSpacing(true);
214 }
215
220 void
221 SetUseImageSpacingOff()
222 {
223 this->SetUseImageSpacing(false);
224 }
225#endif
226
229#if !defined(ITK_LEGACY_REMOVE)
230 using WeightsType [[deprecated("Use DerivativeWeightsType or ComponentWeightsType instead.")]] = ComponentWeightsType;
231#endif
232
236 itkSetMacro(DerivativeWeights, DerivativeWeightsType);
237 itkGetConstReferenceMacro(DerivativeWeights, DerivativeWeightsType);
242 itkSetMacro(ComponentWeights, ComponentWeightsType);
243 itkGetConstReferenceMacro(ComponentWeights, ComponentWeightsType);
251 itkSetMacro(UsePrincipleComponents, bool);
252 itkGetConstMacro(UsePrincipleComponents, bool);
253 itkBooleanMacro(UsePrincipleComponents);
255#if !defined(ITK_FUTURE_LEGACY_REMOVE)
257 void
258 SetUsePrincipleComponentsOn()
259 {
260 this->SetUsePrincipleComponents(true);
261 }
262
264 void
265 SetUsePrincipleComponentsOff()
266 {
267 this->SetUsePrincipleComponents(false);
268 }
269#endif
270
273 static int
274 CubicSolver(const double *, double *);
275
278 itkConceptMacro(RealTypeHasNumericTraitsCheck, (Concept::HasNumericTraits<RealType>));
279
280protected:
283
287 void
289
301 void
302 DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
303
304
305 void
306 PrintSelf(std::ostream & os, Indent indent) const override;
307
308 using ImageBaseType = typename InputImageType::Superclass;
309
311 itkGetConstObjectMacro(RealValuedInputImage, RealVectorImageType);
312
313 TRealType
315 {
316 unsigned int i;
317 unsigned int j;
318 TRealType dx;
319 TRealType sum;
320 auto accum = TRealType{};
321 for (i = 0; i < ImageDimension; ++i)
322 {
323 sum = TRealType{};
324 for (j = 0; j < VectorDimension; ++j)
325 {
326 dx = m_DerivativeWeights[i] * m_SqrtComponentWeights[j] * 0.5 * (it.GetNext(i)[j] - it.GetPrevious(i)[j]);
327 sum += dx * dx;
328 }
329 accum += sum;
330 }
331 return std::sqrt(accum);
332 }
333
334 TRealType
336 {
337 // WARNING: ONLY CALL THIS METHOD WHEN PROCESSING A 3D IMAGE
338 unsigned int i;
339 unsigned int j;
340 double Lambda[3];
341 double CharEqn[3];
342 double ans;
343
344 vnl_matrix<TRealType> g(ImageDimension, ImageDimension);
345 vnl_vector_fixed<TRealType, VectorDimension> d_phi_du[TInputImage::ImageDimension];
346
347 // Calculate the directional derivatives for each vector component using
348 // central differences.
349 for (i = 0; i < ImageDimension; ++i)
350 {
351 for (j = 0; j < VectorDimension; ++j)
352 {
353 d_phi_du[i][j] =
354 m_DerivativeWeights[i] * m_SqrtComponentWeights[j] * 0.5 * (it.GetNext(i)[j] - it.GetPrevious(i)[j]);
355 }
356 }
357
358 // Calculate the symmetric metric tensor g
359 for (i = 0; i < ImageDimension; ++i)
360 {
361 for (j = i; j < ImageDimension; ++j)
362 {
363 g[j][i] = g[i][j] = dot_product(d_phi_du[i], d_phi_du[j]);
364 }
365 }
366
367 // Find the coefficients of the characteristic equation det(g - lambda I)=0
368 // CharEqn[3] = 1.0;
369
370 CharEqn[2] = -(g[0][0] + g[1][1] + g[2][2]);
371
372 CharEqn[1] = (g[0][0] * g[1][1] + g[0][0] * g[2][2] + g[1][1] * g[2][2]) -
373 (g[0][1] * g[1][0] + g[0][2] * g[2][0] + g[1][2] * g[2][1]);
374
375 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]) +
376 g[2][0] * (g[1][1] * g[0][2] - g[0][1] * g[1][2]);
377
378 // Find the eigenvalues of g
379 const int numberOfDistinctRoots = this->CubicSolver(CharEqn, Lambda);
380
381 // Define gradient magnitude as the difference between two largest
382 // eigenvalues. Other definitions may be appropriate here as well.
383 if (numberOfDistinctRoots == 3) // By far the most common case
384 {
385 if (Lambda[0] > Lambda[1])
386 {
387 if (Lambda[1] > Lambda[2]) // Most common, guaranteed?
388 {
389 ans = Lambda[0] - Lambda[1];
390 }
391 else
392 {
393 if (Lambda[0] > Lambda[2])
394 {
395 ans = Lambda[0] - Lambda[2];
396 }
397 else
398 {
399 ans = Lambda[2] - Lambda[0];
400 }
401 }
402 }
403 else
404 {
405 if (Lambda[0] > Lambda[2])
406 {
407 ans = Lambda[1] - Lambda[0];
408 }
409 else
410 {
411 if (Lambda[1] > Lambda[2])
412 {
413 ans = Lambda[1] - Lambda[2];
414 }
415 else
416 {
417 ans = Lambda[2] - Lambda[1];
418 }
419 }
420 }
421 }
422 else if (numberOfDistinctRoots == 2)
423 {
424 if (Lambda[0] > Lambda[1])
425 {
426 ans = Lambda[0] - Lambda[1];
427 }
428 else
429 {
430 ans = Lambda[1] - Lambda[0];
431 }
432 }
433 else if (numberOfDistinctRoots == 1)
434 {
435 ans = 0.0;
436 }
437 else
438 {
439 itkExceptionMacro("Undefined condition. Cubic root solver returned " << numberOfDistinctRoots
440 << " distinct roots.");
441 }
442
443 return ans;
444 }
445
446 // Function is defined here because the templating confuses gcc 2.96 when
447 // defined
448 // in .hxx file. jc 1/29/03
449 TRealType
451 {
452 unsigned int i;
453 unsigned int j;
454
455 vnl_matrix<TRealType> g(ImageDimension, ImageDimension);
456 vnl_vector_fixed<TRealType, VectorDimension> d_phi_du[TInputImage::ImageDimension];
457
458 // Calculate the directional derivatives for each vector component using
459 // central differences.
460 for (i = 0; i < ImageDimension; ++i)
461 {
462 for (j = 0; j < VectorDimension; ++j)
463 {
464 d_phi_du[i][j] =
465 m_DerivativeWeights[i] * m_SqrtComponentWeights[j] * 0.5 * (it.GetNext(i)[j] - it.GetPrevious(i)[j]);
466 }
467 }
468
469 // Calculate the symmetric metric tensor g
470 for (i = 0; i < ImageDimension; ++i)
471 {
472 for (j = i; j < ImageDimension; ++j)
473 {
474 g[j][i] = g[i][j] = dot_product(d_phi_du[i], d_phi_du[j]);
475 }
476 }
477
478 // Find the eigenvalues of g
479 const vnl_symmetric_eigensystem<TRealType> E(g);
480
481 // Return the difference in length between the first two principle axes.
482 // Note that other edge strength metrics may be appropriate here instead..
483 return (E.get_eigenvalue(ImageDimension - 1) - E.get_eigenvalue(ImageDimension - 2));
484 }
485
488
496private:
497 bool m_UseImageSpacing{ true };
499
501
503};
504} // end namespace itk
505
506#ifndef ITK_MANUAL_INSTANTIATION
507# include "itkVectorGradientMagnitudeImageFilter.hxx"
508#endif
509
510#endif
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
PixelType GetPrevious(const unsigned int axis, NeighborIndexType i) const
PixelType GetNext(const unsigned int axis, NeighborIndexType i) const
Simulate a standard C array with copy semantics.
static constexpr FixedArray Filled(const ValueType &value)
typename OutputImageType::RegionType OutputImageRegionType
Templated n-dimensional image class.
Definition itkImage.h:89
Control indentation during Print() invocation.
Definition itkIndent.h:50
Implements transparent reference counting.
~VectorGradientMagnitudeImageFilter() override=default
void PrintSelf(std::ostream &os, Indent indent) const override
Image< RealVectorType, TInputImage::ImageDimension > RealVectorImageType
TRealType EvaluateAtNeighborhood3D(const ConstNeighborhoodIteratorType &it) const
typename OutputImageType::RegionType OutputImageRegionType
ConstNeighborhoodIterator< RealVectorImageType > ConstNeighborhoodIteratorType
virtual void SetUsePrincipleComponents(bool _arg)
static int CubicSolver(const double *, double *)
FixedArray< TRealType, ImageDimension > DerivativeWeightsType
typename ConstNeighborhoodIteratorType::RadiusType RadiusType
TRealType EvaluateAtNeighborhood(const ConstNeighborhoodIteratorType &it) const
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
TRealType NonPCEvaluateAtNeighborhood(const ConstNeighborhoodIteratorType &it) const
Vector< TRealType, InputPixelType::Dimension > RealVectorType
FixedArray< TRealType, VectorDimension > ComponentWeightsType
ImageToImageFilter< TInputImage, TOutputImage > Superclass
A templated class holding a n-Dimensional vector.
Definition itkVector.h:63
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned int ThreadIdType