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 TRealType accum{};
317 for (unsigned int i = 0; i < ImageDimension; ++i)
318 {
319 TRealType sum = TRealType{};
320 for (unsigned int j = 0; j < VectorDimension; ++j)
321 {
322 TRealType dx =
323 m_DerivativeWeights[i] * m_SqrtComponentWeights[j] * 0.5 * (it.GetNext(i)[j] - it.GetPrevious(i)[j]);
324 sum += dx * dx;
325 }
326 accum += sum;
327 }
328 return std::sqrt(accum);
329 }
330
331 TRealType
333 {
334 // WARNING: ONLY CALL THIS METHOD WHEN PROCESSING A 3D IMAGE
335 double Lambda[3];
336 double CharEqn[3];
337 double ans = NAN;
338
339 vnl_matrix<TRealType> g(ImageDimension, ImageDimension);
340 vnl_vector_fixed<TRealType, VectorDimension> d_phi_du[TInputImage::ImageDimension];
341
342 // Calculate the directional derivatives for each vector component using
343 // central differences.
344 for (unsigned int i = 0; i < ImageDimension; ++i)
345 {
346 for (unsigned int j = 0; j < VectorDimension; ++j)
347 {
348 d_phi_du[i][j] =
349 m_DerivativeWeights[i] * m_SqrtComponentWeights[j] * 0.5 * (it.GetNext(i)[j] - it.GetPrevious(i)[j]);
350 }
351 }
352
353 // Calculate the symmetric metric tensor g
354 for (unsigned int i = 0; i < ImageDimension; ++i)
355 {
356 for (unsigned int j = i; j < ImageDimension; ++j)
357 {
358 g[j][i] = g[i][j] = dot_product(d_phi_du[i], d_phi_du[j]);
359 }
360 }
361
362 // Find the coefficients of the characteristic equation det(g - lambda I)=0
363 // CharEqn[3] = 1.0;
364
365 CharEqn[2] = -(g[0][0] + g[1][1] + g[2][2]);
366
367 CharEqn[1] = (g[0][0] * g[1][1] + g[0][0] * g[2][2] + g[1][1] * g[2][2]) -
368 (g[0][1] * g[1][0] + g[0][2] * g[2][0] + g[1][2] * g[2][1]);
369
370 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]) +
371 g[2][0] * (g[1][1] * g[0][2] - g[0][1] * g[1][2]);
372
373 // Find the eigenvalues of g
374 const int numberOfDistinctRoots = this->CubicSolver(CharEqn, Lambda);
375
376 // Define gradient magnitude as the difference between two largest
377 // eigenvalues. Other definitions may be appropriate here as well.
378 if (numberOfDistinctRoots == 3) // By far the most common case
379 {
380 if (Lambda[0] > Lambda[1])
381 {
382 if (Lambda[1] > Lambda[2]) // Most common, guaranteed?
383 {
384 ans = Lambda[0] - Lambda[1];
385 }
386 else
387 {
388 if (Lambda[0] > Lambda[2])
389 {
390 ans = Lambda[0] - Lambda[2];
391 }
392 else
393 {
394 ans = Lambda[2] - Lambda[0];
395 }
396 }
397 }
398 else
399 {
400 if (Lambda[0] > Lambda[2])
401 {
402 ans = Lambda[1] - Lambda[0];
403 }
404 else
405 {
406 if (Lambda[1] > Lambda[2])
407 {
408 ans = Lambda[1] - Lambda[2];
409 }
410 else
411 {
412 ans = Lambda[2] - Lambda[1];
413 }
414 }
415 }
416 }
417 else if (numberOfDistinctRoots == 2)
418 {
419 if (Lambda[0] > Lambda[1])
420 {
421 ans = Lambda[0] - Lambda[1];
422 }
423 else
424 {
425 ans = Lambda[1] - Lambda[0];
426 }
427 }
428 else if (numberOfDistinctRoots == 1)
429 {
430 ans = 0.0;
431 }
432 else
433 {
434 itkExceptionMacro("Undefined condition. Cubic root solver returned " << numberOfDistinctRoots
435 << " distinct roots.");
436 }
437
438 return ans;
439 }
440
441 // Function is defined here because the templating confuses gcc 2.96 when
442 // defined
443 // in .hxx file. jc 1/29/03
444 TRealType
446 {
447 vnl_matrix<TRealType> g(ImageDimension, ImageDimension);
448 vnl_vector_fixed<TRealType, VectorDimension> d_phi_du[TInputImage::ImageDimension];
449
450 // Calculate the directional derivatives for each vector component using
451 // central differences.
452 for (unsigned int i = 0; i < ImageDimension; ++i)
453 {
454 for (unsigned int j = 0; j < VectorDimension; ++j)
455 {
456 d_phi_du[i][j] =
457 m_DerivativeWeights[i] * m_SqrtComponentWeights[j] * 0.5 * (it.GetNext(i)[j] - it.GetPrevious(i)[j]);
458 }
459 }
460
461 // Calculate the symmetric metric tensor g
462 for (unsigned int i = 0; i < ImageDimension; ++i)
463 {
464 for (unsigned int j = i; j < ImageDimension; ++j)
465 {
466 g[j][i] = g[i][j] = dot_product(d_phi_du[i], d_phi_du[j]);
467 }
468 }
469
470 // Find the eigenvalues of g
471 const vnl_symmetric_eigensystem<TRealType> E(g);
472
473 // Return the difference in length between the first two principle axes.
474 // Note that other edge strength metrics may be appropriate here instead..
475 return (E.get_eigenvalue(ImageDimension - 1) - E.get_eigenvalue(ImageDimension - 2));
476 }
477
480
488private:
489 bool m_UseImageSpacing{ true };
491
493
495};
496} // end namespace itk
497
498#ifndef ITK_MANUAL_INSTANTIATION
499# include "itkVectorGradientMagnitudeImageFilter.hxx"
500#endif
501
502#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