ITK
6.0.0
Insight Toolkit
|
#include <itkSymmetricSecondRankTensor.h>
Represent a symmetric tensor of second rank.
This class implements a ND symmetric tensor of second rank.
Since SymmetricSecondRankTensor is a subclass of FixedArray, you can access its components as:
using TensorPixelType = itk::SymmetricSecondRankTensor< float >; TensorPixelType tensor;
tensor[0] = 1.233; tensor[1] = 1.456;
for convenience the indexed access is also available as
tensor(0,0) = 1.233; tensor(2,0) = 1.233;
The Tensor in principle represents a NxN matrix, but given that it is always symmetric the representation can be compacted into a N*(N+1)/2 elements array that derives from the itk::FixedArray<T>
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. The funding for creating this class was largely provided by NAMIC (National Alliance for Medical Image Computing) (https://www.na-mic.org). 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
Definition at line 75 of file itkSymmetricSecondRankTensor.h.
Public Member Functions | |
void | ComputeEigenAnalysis (EigenValuesArrayType &eigenValues, EigenVectorsMatrixType &eigenVectors) const |
void | ComputeEigenValues (EigenValuesArrayType &eigenValues) const |
ComponentType | GetNthComponent (int c) const |
AccumulateValueType | GetTrace () const |
ValueType & | operator() (unsigned int row, unsigned int col) |
const ValueType & | operator() (unsigned int row, unsigned int col) const |
Self | operator* (const RealValueType &r) const |
const Self & | operator*= (const RealValueType &r) |
Self | operator+ (const Self &r) const |
const Self & | operator+= (const Self &r) |
Self | operator- (const Self &r) const |
const Self & | operator-= (const Self &r) |
Self | operator/ (const RealValueType &r) const |
const Self & | operator/= (const RealValueType &r) |
Self & | operator= (const ComponentArrayType r) |
Self & | operator= (const ComponentType &r) |
MatrixType | PostMultiply (const MatrixType &m) const |
MatrixType | PreMultiply (const MatrixType &m) const |
void | SetIdentity () |
void | SetNthComponent (int c, const ComponentType &v) |
SymmetricSecondRankTensor (const ComponentArrayType r) | |
SymmetricSecondRankTensor (const ComponentType &r) | |
template<typename TCoordinateB > | |
SymmetricSecondRankTensor (const SymmetricSecondRankTensor< TCoordinateB, VDimension > &pa) | |
SymmetricSecondRankTensor ()=default | |
template<typename TCoordinateB > | |
Self & | operator= (const SymmetricSecondRankTensor< TCoordinateB, VDimension > &pa) |
template<typename TMatrixValueType > | |
Self | Rotate (const Matrix< TMatrixValueType, VDimension, VDimension > &m) const |
template<typename TMatrixValueType > | |
Self | Rotate (const vnl_matrix_fixed< TMatrixValueType, VDimension, VDimension > &m) const |
template<typename TMatrixValueType > | |
Self | Rotate (const vnl_matrix< TMatrixValueType > &m) const |
Public Member Functions inherited from itk::FixedArray< TComponent, VDimension *(VDimension+1)/2 > | |
Iterator | Begin () |
ConstIterator | Begin () const |
constexpr const_iterator | begin () const noexcept |
constexpr iterator | begin () noexcept |
constexpr const_iterator | cbegin () const noexcept |
constexpr const_iterator | cend () const noexcept |
const_reverse_iterator | crbegin () const |
const_reverse_iterator | crend () const |
ValueType * | data () |
const ValueType * | data () const |
Iterator | End () |
ConstIterator | End () const |
constexpr const_iterator | end () const noexcept |
constexpr iterator | end () noexcept |
void | Fill (const ValueType &) |
FixedArray ()=default | |
FixedArray (const FixedArray< TFixedArrayValueType, VLength > &r) | |
FixedArray (const std::array< ValueType, VLength > &stdArray) | |
FixedArray (const TScalarValue *r) | |
ValueType * | GetDataPointer () |
const ValueType * | GetDataPointer () const |
ITK_UNEQUAL_OPERATOR_MEMBER_FUNCTION (FixedArray) | |
itkLegacyMacro (ReverseIterator rBegin();) itkLegacyMacro(ConstReverseIterator rBegin() const | |
itkLegacyMacro (ReverseIterator rEnd();) itkLegacyMacro(ConstReverseIterator rEnd() const | |
FixedArray & | operator= (const FixedArray< TFixedArrayValueType, VLength > &r) |
FixedArray & | operator= (const ValueType r[VLength]) |
bool | operator== (const FixedArray &r) const |
reverse_iterator | rbegin () |
const_reverse_iterator | rbegin () const |
reverse_iterator | rend () |
const_reverse_iterator | rend () const |
SizeType | Size () const |
constexpr SizeType | size () const |
void | swap (FixedArray &other) noexcept |
FixedArray (const ValueType r[VLength]) | |
FixedArray (const ValueType &) | |
ITK_GCC_PRAGMA_PUSH constexpr ITK_GCC_SUPPRESS_Warray_bounds 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 |
Static Public Member Functions | |
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 Public Attributes | |
static constexpr unsigned int | Dimension = VDimension |
static constexpr unsigned int | InternalDimension = VDimension * (VDimension + 1) / 2 |
Static Public Attributes inherited from itk::FixedArray< TComponent, VDimension *(VDimension+1)/2 > | |
static constexpr unsigned int | Dimension |
static constexpr unsigned int | Length |
using itk::SymmetricSecondRankTensor< TComponent, VDimension >::AccumulateValueType = typename NumericTraits<ValueType>::RealType |
Definition at line 99 of file itkSymmetricSecondRankTensor.h.
using itk::SymmetricSecondRankTensor< TComponent, VDimension >::BaseArray = FixedArray<TComponent, Self::InternalDimension> |
Convenience type alias.
Definition at line 87 of file itkSymmetricSecondRankTensor.h.
using itk::SymmetricSecondRankTensor< TComponent, VDimension >::ComponentArrayType = ComponentType[Self::InternalDimension] |
Definition at line 122 of file itkSymmetricSecondRankTensor.h.
using itk::SymmetricSecondRankTensor< TComponent, VDimension >::ComponentType = TComponent |
Define the component type.
Definition at line 97 of file itkSymmetricSecondRankTensor.h.
using itk::SymmetricSecondRankTensor< TComponent, VDimension >::EigenValuesArrayType = FixedArray<TComponent, VDimension> |
Array of eigen-values.
Definition at line 90 of file itkSymmetricSecondRankTensor.h.
using itk::SymmetricSecondRankTensor< TComponent, VDimension >::EigenVectorsMatrixType = Matrix<TComponent, VDimension, VDimension> |
Definition at line 94 of file itkSymmetricSecondRankTensor.h.
using itk::SymmetricSecondRankTensor< TComponent, VDimension >::MatrixType = Matrix<TComponent, VDimension, VDimension> |
Matrix of eigen-vectors.
Definition at line 93 of file itkSymmetricSecondRankTensor.h.
using itk::SymmetricSecondRankTensor< TComponent, VDimension >::RealValueType = typename NumericTraits<ValueType>::RealType |
Definition at line 100 of file itkSymmetricSecondRankTensor.h.
using itk::SymmetricSecondRankTensor< TComponent, VDimension >::Self = SymmetricSecondRankTensor |
Standard class type aliases.
Definition at line 79 of file itkSymmetricSecondRankTensor.h.
using itk::SymmetricSecondRankTensor< TComponent, VDimension >::Superclass = FixedArray<TComponent, VDimension *(VDimension + 1) / 2> |
Definition at line 80 of file itkSymmetricSecondRankTensor.h.
using itk::SymmetricSecondRankTensor< TComponent, VDimension >::SymmetricEigenAnalysisType = SymmetricEigenAnalysisFixedDimension<Dimension, MatrixType, EigenValuesArrayType, EigenVectorsMatrixType> |
Definition at line 103 of file itkSymmetricSecondRankTensor.h.
|
default |
Default-constructor.
|
inline |
Definition at line 114 of file itkSymmetricSecondRankTensor.h.
|
inline |
Constructor to enable casting...
Definition at line 118 of file itkSymmetricSecondRankTensor.h.
|
inline |
Pass-through constructor for the Array base class.
Definition at line 125 of file itkSymmetricSecondRankTensor.h.
void itk::SymmetricSecondRankTensor< TComponent, VDimension >::ComputeEigenAnalysis | ( | EigenValuesArrayType & | eigenValues, |
EigenVectorsMatrixType & | eigenVectors | ||
) | const |
Return an array containing EigenValues, and a matrix containing Eigen vectors.
void itk::SymmetricSecondRankTensor< TComponent, VDimension >::ComputeEigenValues | ( | EigenValuesArrayType & | eigenValues | ) | const |
Return an array containing EigenValues.
|
inline |
Return the value for the Nth component.
Definition at line 182 of file itkSymmetricSecondRankTensor.h.
|
inlinestatic |
Return the number of components.
Definition at line 175 of file itkSymmetricSecondRankTensor.h.
AccumulateValueType itk::SymmetricSecondRankTensor< TComponent, VDimension >::GetTrace | ( | ) | const |
Get Trace value
ValueType& itk::SymmetricSecondRankTensor< TComponent, VDimension >::operator() | ( | unsigned int | row, |
unsigned int | col | ||
) |
Matrix notation, in const and non-const forms.
const ValueType& itk::SymmetricSecondRankTensor< TComponent, VDimension >::operator() | ( | unsigned int | row, |
unsigned int | col | ||
) | const |
Self itk::SymmetricSecondRankTensor< TComponent, VDimension >::operator* | ( | const RealValueType & | r | ) | const |
Arithmetic operations between tensors and scalars
const Self& itk::SymmetricSecondRankTensor< TComponent, VDimension >::operator*= | ( | const RealValueType & | r | ) |
Self itk::SymmetricSecondRankTensor< TComponent, VDimension >::operator+ | ( | const Self & | r | ) | const |
Arithmetic operations between pixels. Return a new SymmetricSecondRankTensor.
const Self& itk::SymmetricSecondRankTensor< TComponent, VDimension >::operator+= | ( | const Self & | r | ) |
Self itk::SymmetricSecondRankTensor< TComponent, VDimension >::operator- | ( | const Self & | r | ) | const |
const Self& itk::SymmetricSecondRankTensor< TComponent, VDimension >::operator-= | ( | const Self & | r | ) |
Self itk::SymmetricSecondRankTensor< TComponent, VDimension >::operator/ | ( | const RealValueType & | r | ) | const |
const Self& itk::SymmetricSecondRankTensor< TComponent, VDimension >::operator/= | ( | const RealValueType & | r | ) |
Self& itk::SymmetricSecondRankTensor< TComponent, VDimension >::operator= | ( | const ComponentArrayType | r | ) |
Self& itk::SymmetricSecondRankTensor< TComponent, VDimension >::operator= | ( | const ComponentType & | r | ) |
Pass-through assignment operator for the Array base class.
|
inline |
Templated Pass-through assignment for the Array base class.
Definition at line 132 of file itkSymmetricSecondRankTensor.h.
Referenced by itk::DiffusionTensor3D< TComponent >::operator=().
MatrixType itk::SymmetricSecondRankTensor< TComponent, VDimension >::PostMultiply | ( | const MatrixType & | m | ) | const |
MatrixType itk::SymmetricSecondRankTensor< TComponent, VDimension >::PreMultiply | ( | const MatrixType & | m | ) | const |
Self itk::SymmetricSecondRankTensor< TComponent, VDimension >::Rotate | ( | const Matrix< TMatrixValueType, VDimension, VDimension > & | m | ) | const |
Returns the tensor rotated by the provided matrix. ResultingTensor = Matrix * ThisTensor * Matrix.GetTranspose()
|
inline |
Returns the tensor rotated by the provided matrix. ResultingTensor = Matrix * ThisTensor * Matrix.GetTranspose()
Definition at line 233 of file itkSymmetricSecondRankTensor.h.
|
inline |
Returns the tensor rotated by the provided matrix. ResultingTensor = Matrix * ThisTensor * Matrix.GetTranspose()
Definition at line 227 of file itkSymmetricSecondRankTensor.h.
void itk::SymmetricSecondRankTensor< TComponent, VDimension >::SetIdentity | ( | ) |
Set the tensor to an identity tensor. This has 1 in its diagonal elements and zero elsewhere.
|
inline |
Set the Nth component to v.
Definition at line 189 of file itkSymmetricSecondRankTensor.h.
|
staticconstexpr |
Dimension of the vector space.
Definition at line 83 of file itkSymmetricSecondRankTensor.h.
|
staticconstexpr |
Definition at line 84 of file itkSymmetricSecondRankTensor.h.