ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches

#include <itkDiffusionTensor3D.h>

Detailed Description

template<typename TComponent>
class itk::DiffusionTensor3D< TComponent >

Represent a diffusion tensor as used in DTI images.

This class implements a 3D symmetric tensor as it is used for representing diffusion of water molecules in Diffusion Tensor Images.

This class derives from the SymmetricSecondRankTensor, inheriting most of the Tensor-related behavior. At this level we add the methods that are specific to 3D and that are closely related to the concept of diffusion.

Author
Jeffrey Duda from School of Engineering at University of Pennsylvania
Torsten Rohlfing from SRI International Neuroscience Program.

This class was mostly based on files that Jeffrey Duda, Torsten Rohlfing and Martin Styner contributed to the ITK users list during a discussion on support for DiffusionTensorImages. A discussion on the design of this class can be found in the WIKI pages of NAMIC:

https://www.na-mic.org/Wiki/index.php/NAMIC_Wiki:DTI:ITK-DiffusionTensorPixelType

Note
This work is part of the National Alliance for Medical Image Computing (NAMIC), funded by the National Institutes of Health through the NIH Roadmap for Medical Research, Grant U54 EB005149. Information on the National Centers for Biomedical Computing can be obtained from http://commonfund.nih.gov/bioinformatics.
Contributions by Torsten Rohlfing were funded by the following NIH grants

Alcohol, HIV and the Brain, NIAAA AA12999, PI: A. Pfefferbaum

Normal Aging of Brain Structure and Function NIA AG 17919, PI: E.V. Sullivan.

For algorithmic details see [78].

See also
SymmetricSecondRankTensor

Definition at line 76 of file itkDiffusionTensor3D.h.

+ Inheritance diagram for itk::DiffusionTensor3D< TComponent >:
+ Collaboration diagram for itk::DiffusionTensor3D< TComponent >:

Public Types

using AccumulateValueType
 
using ComponentArrayType
 
using ComponentType
 
using EigenValuesArrayType
 
using EigenVectorsMatrixType
 
using RealValueType
 
using Self = DiffusionTensor3D
 
using Superclass = SymmetricSecondRankTensor<TComponent, 3>
 
- Public Types inherited from itk::SymmetricSecondRankTensor< TComponent, 3 >
using AccumulateValueType
 
using AccumulateValueType
 
using BaseArray
 
using BaseArray
 
using ComponentArrayType
 
using ComponentArrayType
 
using ComponentType
 
using ComponentType
 
using EigenValuesArrayType
 
using EigenValuesArrayType
 
using EigenVectorsMatrixType
 
using EigenVectorsMatrixType
 
using MatrixType
 
using MatrixType
 
using RealValueType
 
using RealValueType
 
using Self
 
using Self
 
using Superclass
 
using Superclass
 
using SymmetricEigenAnalysisType
 
using SymmetricEigenAnalysisType
 
- Public Types inherited from itk::FixedArray< TComponent, VDimension *(VDimension+1)/2 >
using CArray
 
using const_iterator
 
using const_pointer
 
using const_reference
 
using const_reverse_iterator
 
using ConstIterator
 
using Iterator
 
using iterator
 
using pointer
 
using reference
 
using reverse_iterator
 
using SizeType
 
using value_type
 
using ValueType
 

Public Member Functions

 DiffusionTensor3D ()=default
 
template<typename TCoordinateB>
 DiffusionTensor3D (const DiffusionTensor3D< TCoordinateB > &pa)
 
RealValueType GetFractionalAnisotropy () const
 
RealValueType GetInnerScalarProduct () const
 
RealValueType GetRelativeAnisotropy () const
 
AccumulateValueType GetTrace () const
 
Selfoperator= (const ComponentArrayType r)
 
Selfoperator= (const ComponentType &r)
 
Selfoperator= (const Superclass &r)
 
 DiffusionTensor3D (const Superclass &r)
 
 DiffusionTensor3D (const ComponentType &r)
 
 DiffusionTensor3D (const ComponentArrayType r)
 
template<typename TCoordinateB>
Selfoperator= (const DiffusionTensor3D< TCoordinateB > &pa)
 
- Public Member Functions inherited from itk::SymmetricSecondRankTensor< TComponent, 3 >
void ComputeEigenAnalysis (EigenValuesArrayType &eigenValues, EigenVectorsMatrixType &eigenVectors) const
 
void ComputeEigenAnalysis (EigenValuesArrayType &eigenValues, EigenVectorsMatrixType &eigenVectors) const
 
void ComputeEigenValues (EigenValuesArrayType &eigenValues) const
 
void ComputeEigenValues (EigenValuesArrayType &eigenValues) const
 
ComponentType GetNthComponent (int c) const
 
ComponentType GetNthComponent (int c) const
 
AccumulateValueType GetTrace () const
 
AccumulateValueType GetTrace () const
 
ValueTypeoperator() (unsigned int row, unsigned int col)
 
ValueTypeoperator() (unsigned int row, unsigned int col)
 
const ValueTypeoperator() (unsigned int row, unsigned int col) const
 
const ValueTypeoperator() (unsigned int row, unsigned int col) const
 
Self operator* (const RealValueType &r) const
 
Self operator* (const RealValueType &r) const
 
const Selfoperator*= (const RealValueType &r)
 
const Selfoperator*= (const RealValueType &r)
 
Self operator+ (const Self &r) const
 
Self operator+ (const Self &r) const
 
const Selfoperator+= (const Self &r)
 
const Selfoperator+= (const Self &r)
 
Self operator- (const Self &r) const
 
Self operator- (const Self &r) const
 
const Selfoperator-= (const Self &r)
 
const Selfoperator-= (const Self &r)
 
Self operator/ (const RealValueType &r) const
 
Self operator/ (const RealValueType &r) const
 
const Selfoperator/= (const RealValueType &r)
 
const Selfoperator/= (const RealValueType &r)
 
Selfoperator= (const ComponentArrayType r)
 
Selfoperator= (const ComponentArrayType r)
 
Selfoperator= (const ComponentType &r)
 
Selfoperator= (const ComponentType &r)
 
MatrixType PostMultiply (const MatrixType &m) const
 
MatrixType PostMultiply (const MatrixType &m) const
 
MatrixType PreMultiply (const MatrixType &m) const
 
MatrixType PreMultiply (const MatrixType &m) const
 
void SetIdentity ()
 
void SetIdentity ()
 
void SetNthComponent (int c, const ComponentType &v)
 
void SetNthComponent (int c, const ComponentType &v)
 
 SymmetricSecondRankTensor (const ComponentArrayType r)
 
 SymmetricSecondRankTensor (const ComponentArrayType r)
 
 SymmetricSecondRankTensor (const ComponentType &r)
 
 SymmetricSecondRankTensor (const ComponentType &r)
 
 SymmetricSecondRankTensor (const SymmetricSecondRankTensor< TCoordinateB, VDimension > &pa)
 
 SymmetricSecondRankTensor (const SymmetricSecondRankTensor< TCoordinateB, VDimension > &pa)
 
 SymmetricSecondRankTensor ()=default
 
 SymmetricSecondRankTensor ()=default
 
Selfoperator= (const SymmetricSecondRankTensor< TCoordinateB, VDimension > &pa)
 
Selfoperator= (const SymmetricSecondRankTensor< TCoordinateB, VDimension > &pa)
 
Self Rotate (const Matrix< TMatrixValueType, VDimension, VDimension > &m) const
 
Self Rotate (const vnl_matrix_fixed< TMatrixValueType, VDimension, VDimension > &m) const
 
Self Rotate (const vnl_matrix< TMatrixValueType > &m) const
 
Self Rotate (const Matrix< TMatrixValueType, VDimension, VDimension > &m) const
 
Self Rotate (const vnl_matrix_fixed< TMatrixValueType, VDimension, VDimension > &m) const
 
Self Rotate (const vnl_matrix< TMatrixValueType > &m) const
 
- Public Member Functions inherited from itk::FixedArray< TComponent, VDimension *(VDimension+1)/2 >
Iterator Begin ()
 
Iterator Begin ()
 
ConstIterator Begin () const
 
ConstIterator Begin () const
 
constexpr const_iterator begin () const noexcept
 
constexpr const_iterator begin () const noexcept
 
constexpr iterator begin () noexcept
 
constexpr iterator begin () noexcept
 
constexpr const_iterator cbegin () const noexcept
 
constexpr const_iterator cbegin () const noexcept
 
constexpr const_iterator cend () const noexcept
 
constexpr const_iterator cend () const noexcept
 
const_reverse_iterator crbegin () const
 
const_reverse_iterator crbegin () const
 
const_reverse_iterator crend () const
 
const_reverse_iterator crend () const
 
ValueTypedata ()
 
ValueTypedata ()
 
const ValueTypedata () const
 
const ValueTypedata () const
 
Iterator End ()
 
Iterator End ()
 
ConstIterator End () const
 
ConstIterator End () const
 
constexpr const_iterator end () const noexcept
 
constexpr const_iterator end () const noexcept
 
constexpr iterator end () noexcept
 
constexpr iterator end () noexcept
 
void Fill (const ValueType &)
 
void Fill (const ValueType &)
 
 FixedArray ()=default
 
 FixedArray ()=default
 
 FixedArray (const FixedArray< TFixedArrayValueType, VLength > &r)
 
 FixedArray (const FixedArray< TFixedArrayValueType, VLength > &r)
 
 FixedArray (const std::array< ValueType, VLength > &stdArray)
 
 FixedArray (const std::array< ValueType, VLength > &stdArray)
 
 FixedArray (const TScalarValue *r)
 
 FixedArray (const TScalarValue *r)
 
ValueTypeGetDataPointer ()
 
ValueTypeGetDataPointer ()
 
const ValueTypeGetDataPointer () const
 
const ValueTypeGetDataPointer () const
 
 ITK_UNEQUAL_OPERATOR_MEMBER_FUNCTION (FixedArray)
 
 ITK_UNEQUAL_OPERATOR_MEMBER_FUNCTION (FixedArray)
 
FixedArrayoperator= (const FixedArray< TFixedArrayValueType, VLength > &r)
 
FixedArrayoperator= (const FixedArray< TFixedArrayValueType, VLength > &r)
 
FixedArrayoperator= (const ValueType r[VLength])
 
FixedArrayoperator= (const ValueType r[VLength])
 
bool operator== (const FixedArray &r) const
 
bool operator== (const FixedArray &r) const
 
ReverseIterator rBegin ()
 
ReverseIterator rBegin ()
 
ConstReverseIterator rBegin () const
 
ConstReverseIterator rBegin () const
 
reverse_iterator rbegin ()
 
reverse_iterator rbegin ()
 
const_reverse_iterator rbegin () const
 
const_reverse_iterator rbegin () const
 
ReverseIterator rEnd ()
 
ReverseIterator rEnd ()
 
ConstReverseIterator rEnd () const
 
ConstReverseIterator rEnd () const
 
reverse_iterator rend ()
 
reverse_iterator rend ()
 
const_reverse_iterator rend () const
 
const_reverse_iterator rend () const
 
SizeType Size () const
 
SizeType Size () const
 
constexpr SizeType size () const
 
constexpr SizeType size () const
 
void swap (FixedArray &other) noexcept
 
void swap (FixedArray &other) noexcept
 
 FixedArray (const ValueType r[VLength])
 
 FixedArray (const ValueType &)
 
 FixedArray (const ValueType r[VLength])
 
 FixedArray (const ValueType &)
 
ITK_GCC_PRAGMA_PUSH ITK_GCC_SUPPRESS_Warray_bounds constexpr reference operator[] (unsigned int index)
 
constexpr const_reference operator[] (unsigned int index) const
 
ITK_GCC_PRAGMA_PUSH ITK_GCC_SUPPRESS_Warray_bounds constexpr reference operator[] (unsigned int index)
 
constexpr const_reference operator[] (unsigned int index) const
 
ITK_GCC_PRAGMA_POP void SetElement (unsigned int index, const_reference value)
 
const_reference GetElement (unsigned int index) const
 
ITK_GCC_PRAGMA_POP void SetElement (unsigned int index, const_reference value)
 
const_reference GetElement (unsigned int index) const
 

Additional Inherited Members

- Static Public Member Functions inherited from itk::SymmetricSecondRankTensor< TComponent, 3 >
static unsigned int GetNumberOfComponents ()
 
static unsigned int GetNumberOfComponents ()
 
- Static Public Member Functions inherited from itk::FixedArray< TComponent, VDimension *(VDimension+1)/2 >
static constexpr FixedArray Filled (const ValueType &value)
 
static constexpr FixedArray Filled (const ValueType &value)
 
- Static Public Attributes inherited from itk::SymmetricSecondRankTensor< TComponent, 3 >
static constexpr unsigned int Dimension
 
static constexpr unsigned int Dimension
 
static constexpr unsigned int InternalDimension
 
static constexpr unsigned int InternalDimension
 
- Static Public Attributes inherited from itk::FixedArray< TComponent, VDimension *(VDimension+1)/2 >
static constexpr unsigned int Dimension
 
static constexpr unsigned int Dimension
 
static constexpr unsigned int Length
 
static constexpr unsigned int Length
 

Member Typedef Documentation

◆ AccumulateValueType

template<typename TComponent>
using itk::SymmetricSecondRankTensor< TComponent, VDimension >::AccumulateValueType

Definition at line 99 of file itkSymmetricSecondRankTensor.h.

◆ ComponentArrayType

template<typename TComponent>
using itk::SymmetricSecondRankTensor< TComponent, VDimension >::ComponentArrayType

Definition at line 122 of file itkSymmetricSecondRankTensor.h.

◆ ComponentType

template<typename TComponent>
using itk::SymmetricSecondRankTensor< TComponent, VDimension >::ComponentType

Define the component type.

Definition at line 97 of file itkSymmetricSecondRankTensor.h.

◆ EigenValuesArrayType

template<typename TComponent>
using itk::SymmetricSecondRankTensor< TComponent, VDimension >::EigenValuesArrayType

Array of eigen-values.

Definition at line 90 of file itkSymmetricSecondRankTensor.h.

◆ EigenVectorsMatrixType

template<typename TComponent>
using itk::SymmetricSecondRankTensor< TComponent, VDimension >::EigenVectorsMatrixType

Definition at line 94 of file itkSymmetricSecondRankTensor.h.

◆ RealValueType

template<typename TComponent>
using itk::SymmetricSecondRankTensor< TComponent, VDimension >::RealValueType

Definition at line 100 of file itkSymmetricSecondRankTensor.h.

◆ Self

template<typename TComponent>
using itk::DiffusionTensor3D< TComponent >::Self = DiffusionTensor3D

Standard class type aliases.

Definition at line 80 of file itkDiffusionTensor3D.h.

◆ Superclass

template<typename TComponent>
using itk::DiffusionTensor3D< TComponent >::Superclass = SymmetricSecondRankTensor<TComponent, 3>

Definition at line 81 of file itkDiffusionTensor3D.h.

Constructor & Destructor Documentation

◆ DiffusionTensor3D() [1/5]

template<typename TComponent>
itk::DiffusionTensor3D< TComponent >::DiffusionTensor3D ( )
default

Default Constructor.

◆ DiffusionTensor3D() [2/5]

template<typename TComponent>
itk::DiffusionTensor3D< TComponent >::DiffusionTensor3D ( const Superclass & r)

Constructor with initialization.

◆ DiffusionTensor3D() [3/5]

template<typename TComponent>
itk::DiffusionTensor3D< TComponent >::DiffusionTensor3D ( const ComponentType & r)

Constructor with initialization.

◆ DiffusionTensor3D() [4/5]

template<typename TComponent>
itk::DiffusionTensor3D< TComponent >::DiffusionTensor3D ( const ComponentArrayType r)

Constructor with initialization.

◆ DiffusionTensor3D() [5/5]

template<typename TComponent>
template<typename TCoordinateB>
itk::DiffusionTensor3D< TComponent >::DiffusionTensor3D ( const DiffusionTensor3D< TCoordinateB > & pa)
inline

Constructor to enable casting...

Definition at line 105 of file itkDiffusionTensor3D.h.

Member Function Documentation

◆ GetFractionalAnisotropy()

template<typename TComponent>
RealValueType itk::DiffusionTensor3D< TComponent >::GetFractionalAnisotropy ( ) const

Get the value of Fractional Anisotropy from the Tensor.

◆ GetInnerScalarProduct()

template<typename TComponent>
RealValueType itk::DiffusionTensor3D< TComponent >::GetInnerScalarProduct ( ) const

Get the inner scalar product from the Tensor.

◆ GetRelativeAnisotropy()

template<typename TComponent>
RealValueType itk::DiffusionTensor3D< TComponent >::GetRelativeAnisotropy ( ) const

Get the value of Relative Anisotropy from the Tensor.

◆ GetTrace()

template<typename TComponent>
AccumulateValueType itk::DiffusionTensor3D< TComponent >::GetTrace ( ) const

Get the trace value.

Note that the indices are related to the fact that we store only the upper-right triangle of the matrix. Like

  | 0  1  2  |
  | X  3  4  |
  | X  X  5  |

The trace is therefore the sum of the components M[0], M[3] and M[5].

◆ operator=() [1/4]

template<typename TComponent>
Self & itk::DiffusionTensor3D< TComponent >::operator= ( const ComponentArrayType r)

◆ operator=() [2/4]

template<typename TComponent>
Self & itk::DiffusionTensor3D< TComponent >::operator= ( const ComponentType & r)

◆ operator=() [3/4]

template<typename TComponent>
template<typename TCoordinateB>
Self & itk::DiffusionTensor3D< TComponent >::operator= ( const DiffusionTensor3D< TCoordinateB > & pa)
inline

Templated Pass-through assignment for the Array base class.

Definition at line 122 of file itkDiffusionTensor3D.h.

◆ operator=() [4/4]

template<typename TComponent>
Self & itk::DiffusionTensor3D< TComponent >::operator= ( const Superclass & r)

Pass-through assignment operator for the Array base class.


The documentation for this class was generated from the following file: