Medical Imaging Interaction Toolkit
2016.11.0
Medical Imaging Interaction Toolkit
|
Represents an ODF for Q-Ball imaging. More...
#include <itkOrientationDistributionFunction.h>
Public Types | |
enum | InterpolationMethods { ODF_NEAREST_NEIGHBOR_INTERP, ODF_TRILINEAR_BARYCENTRIC_INTERP, ODF_SPHERICAL_GAUSSIAN_BASIS_FUNCTIONS } |
typedef OrientationDistributionFunction | Self |
typedef FixedArray< TComponent, NOdfDirections > | Superclass |
typedef FixedArray< TComponent, itkGetStaticConstMacro(InternalDimension)> | BaseArray |
typedef TComponent | ComponentType |
typedef Superclass::ValueType | ValueType |
typedef NumericTraits< ValueType >::RealType | AccumulateValueType |
typedef NumericTraits< ValueType >::RealType | RealValueType |
typedef Matrix< TComponent, NOdfDirections, NOdfDirections > | MatrixType |
typedef vnl_matrix_fixed< double, 3, NOdfDirections > | DirectionsType |
typedef ComponentType | ComponentArrayType[itkGetStaticConstMacro(InternalDimension)] |
Public Member Functions | |
itkStaticConstMacro (InternalDimension, unsigned int, NOdfDirections) | |
OrientationDistributionFunction () | |
OrientationDistributionFunction (const ComponentType &r) | |
OrientationDistributionFunction (const Self &r) | |
OrientationDistributionFunction (const ComponentArrayType r) | |
Self & | operator= (const Self &r) |
Self & | operator= (const ComponentType &r) |
Self & | operator= (const ComponentArrayType r) |
Self | operator+ (const Self &vec) const |
Self | operator- (const Self &vec) const |
const Self & | operator+= (const Self &vec) |
const Self & | operator-= (const Self &vec) |
Self | operator* (const RealValueType &scalar) const |
Self | operator/ (const RealValueType &scalar) const |
const Self & | operator*= (const RealValueType &scalar) |
const Self & | operator/= (const RealValueType &scalar) |
ComponentType | GetNthComponent (int c) const |
ComponentType | GetInterpolatedComponent (vnl_vector_fixed< double, 3 > dir, InterpolationMethods method) const |
void | SetNthComponent (int c, const ComponentType &v) |
ValueType & | operator() (unsigned int row, unsigned int col) |
const ValueType & | operator() (unsigned int row, unsigned int col) const |
void | SetIsotropic () |
void | InitFromTensor (itk::DiffusionTensor3D< TComponent > tensor) |
void | InitFromEllipsoid (itk::DiffusionTensor3D< TComponent > tensor) |
Self | PreMultiply (const MatrixType &m) const |
Self | PostMultiply (const MatrixType &m) const |
void | Normalize () |
Self | MinMaxNormalize () const |
Self | MaxNormalize () const |
void | L2Normalize () |
int | GetPrincipleDiffusionDirection () const |
int | GetNthDiffusionDirection (int n, vnl_vector_fixed< double, 3 > rndVec) const |
TComponent | GetGeneralizedFractionalAnisotropy () const |
TComponent | GetGeneralizedGFA (int k, int p) const |
TComponent | GetNormalizedEntropy () const |
TComponent | GetNematicOrderParameter () const |
TComponent | GetStdDevByMaxValue () const |
ComponentType | GetMaxValue () const |
ComponentType | GetMinValue () const |
ComponentType | GetMeanValue () const |
TComponent | GetPrincipleCurvature (double alphaMinDegree, double alphaMaxDegree, int invert) const |
Static Public Member Functions | |
static DirectionsType * | GetDirections () |
static unsigned int | GetNumberOfComponents () |
static std::vector< int > | GetNeighbors (int idx) |
static vtkPolyData * | GetBaseMesh () |
static void | ComputeBaseMesh () |
static double | GetMaxChordLength () |
static vnl_vector_fixed< double, 3 > | GetDirection (int i) |
Represents an ODF for Q-Ball imaging.
Reference: David S. Tuch, Q-ball imaging, Magnetic Resonance in Medicine Volume 52 Issue 6, Pages 1358 - 1372
Definition at line 50 of file itkOrientationDistributionFunction.h.
typedef NumericTraits<ValueType>::RealType itk::OrientationDistributionFunction< TComponent, NOdfDirections >::AccumulateValueType |
Definition at line 75 of file itkOrientationDistributionFunction.h.
typedef FixedArray<TComponent, itkGetStaticConstMacro(InternalDimension)> itk::OrientationDistributionFunction< TComponent, NOdfDirections >::BaseArray |
Convenience typedefs.
Definition at line 70 of file itkOrientationDistributionFunction.h.
typedef ComponentType itk::OrientationDistributionFunction< TComponent, NOdfDirections >::ComponentArrayType[itkGetStaticConstMacro(InternalDimension)] |
Definition at line 87 of file itkOrientationDistributionFunction.h.
typedef TComponent itk::OrientationDistributionFunction< TComponent, NOdfDirections >::ComponentType |
Define the component type.
Definition at line 73 of file itkOrientationDistributionFunction.h.
typedef vnl_matrix_fixed<double, 3, NOdfDirections> itk::OrientationDistributionFunction< TComponent, NOdfDirections >::DirectionsType |
Definition at line 80 of file itkOrientationDistributionFunction.h.
typedef Matrix<TComponent, NOdfDirections, NOdfDirections> itk::OrientationDistributionFunction< TComponent, NOdfDirections >::MatrixType |
Definition at line 78 of file itkOrientationDistributionFunction.h.
typedef NumericTraits<ValueType>::RealType itk::OrientationDistributionFunction< TComponent, NOdfDirections >::RealValueType |
Definition at line 76 of file itkOrientationDistributionFunction.h.
typedef OrientationDistributionFunction itk::OrientationDistributionFunction< TComponent, NOdfDirections >::Self |
Standard class typedefs.
Definition at line 62 of file itkOrientationDistributionFunction.h.
typedef FixedArray<TComponent,NOdfDirections> itk::OrientationDistributionFunction< TComponent, NOdfDirections >::Superclass |
Definition at line 63 of file itkOrientationDistributionFunction.h.
typedef Superclass::ValueType itk::OrientationDistributionFunction< TComponent, NOdfDirections >::ValueType |
Definition at line 74 of file itkOrientationDistributionFunction.h.
enum itk::OrientationDistributionFunction::InterpolationMethods |
Enumerator | |
---|---|
ODF_NEAREST_NEIGHBOR_INTERP | |
ODF_TRILINEAR_BARYCENTRIC_INTERP | |
ODF_SPHERICAL_GAUSSIAN_BASIS_FUNCTIONS |
Definition at line 55 of file itkOrientationDistributionFunction.h.
|
inline |
Default constructor has nothing to do.
Definition at line 83 of file itkOrientationDistributionFunction.h.
|
inline |
Definition at line 85 of file itkOrientationDistributionFunction.h.
|
inline |
Pass-through constructor for the Array base class.
Definition at line 90 of file itkOrientationDistributionFunction.h.
|
inline |
Definition at line 91 of file itkOrientationDistributionFunction.h.
|
static |
|
inlinestatic |
Definition at line 183 of file itkOrientationDistributionFunction.h.
Referenced by QmitkODFRenderWidget::GenerateODF().
|
static |
Referenced by mitk::TrackingForestHandler< ShOrder, NumberOfSignalFeatures >::CalculateFeaturesForTraining(), itk::FiniteDiffOdfMaximaExtractionFilter< PixelType, ShOrder, NrOdfDirections >::FindCandidatePeaks(), itk::MLBSTrackingFilter< ShOrder, NumImageFeatures >::GetNewDirection(), and mitk::TrackingForestHandler< ShOrder, NumberOfSignalFeatures >::InitForTracking().
|
inlinestatic |
Return the number of components.
Definition at line 111 of file itkOrientationDistributionFunction.h.
Referenced by itk::FiniteDiffOdfMaximaExtractionFilter< PixelType, ShOrder, NrOdfDirections >::BeforeThreadedGenerateData(), and itk::ShCoefficientImageImporter< PixelType, ShOrder >::GetSphericalOdfDirections().
TComponent itk::OrientationDistributionFunction< TComponent, NOdfDirections >::GetGeneralizedFractionalAnisotropy | ( | ) | const |
TComponent itk::OrientationDistributionFunction< TComponent, NOdfDirections >::GetGeneralizedGFA | ( | int | k, |
int | p | ||
) | const |
ComponentType itk::OrientationDistributionFunction< TComponent, NOdfDirections >::GetInterpolatedComponent | ( | vnl_vector_fixed< double, 3 > | dir, |
InterpolationMethods | method | ||
) | const |
Return the value for the Nth component.
|
static |
ComponentType itk::OrientationDistributionFunction< TComponent, NOdfDirections >::GetMaxValue | ( | ) | const |
ComponentType itk::OrientationDistributionFunction< TComponent, NOdfDirections >::GetMeanValue | ( | ) | const |
ComponentType itk::OrientationDistributionFunction< TComponent, NOdfDirections >::GetMinValue | ( | ) | const |
|
static |
TComponent itk::OrientationDistributionFunction< TComponent, NOdfDirections >::GetNematicOrderParameter | ( | ) | const |
TComponent itk::OrientationDistributionFunction< TComponent, NOdfDirections >::GetNormalizedEntropy | ( | ) | const |
|
inline |
Return the value for the Nth component.
Definition at line 123 of file itkOrientationDistributionFunction.h.
int itk::OrientationDistributionFunction< TComponent, NOdfDirections >::GetNthDiffusionDirection | ( | int | n, |
vnl_vector_fixed< double, 3 > | rndVec | ||
) | const |
|
inlinestatic |
Return the number of components.
Definition at line 117 of file itkOrientationDistributionFunction.h.
TComponent itk::OrientationDistributionFunction< TComponent, NOdfDirections >::GetPrincipleCurvature | ( | double | alphaMinDegree, |
double | alphaMaxDegree, | ||
int | invert | ||
) | const |
int itk::OrientationDistributionFunction< TComponent, NOdfDirections >::GetPrincipleDiffusionDirection | ( | ) | const |
TComponent itk::OrientationDistributionFunction< TComponent, NOdfDirections >::GetStdDevByMaxValue | ( | ) | const |
void itk::OrientationDistributionFunction< TComponent, NOdfDirections >::InitFromEllipsoid | ( | itk::DiffusionTensor3D< TComponent > | tensor | ) |
Evaluate diffusion tensor as ellipsoid.
void itk::OrientationDistributionFunction< TComponent, NOdfDirections >::InitFromTensor | ( | itk::DiffusionTensor3D< TComponent > | tensor | ) |
itk::OrientationDistributionFunction< TComponent, NOdfDirections >::itkStaticConstMacro | ( | InternalDimension | , |
unsigned | int, | ||
NOdfDirections | |||
) |
Dimension of the vector space.
void itk::OrientationDistributionFunction< TComponent, NOdfDirections >::L2Normalize | ( | ) |
Self itk::OrientationDistributionFunction< TComponent, NOdfDirections >::MaxNormalize | ( | ) | const |
Referenced by vtkOdfSource::RequestData().
Self itk::OrientationDistributionFunction< TComponent, NOdfDirections >::MinMaxNormalize | ( | ) | const |
Referenced by vtkOdfSource::RequestData().
void itk::OrientationDistributionFunction< TComponent, NOdfDirections >::Normalize | ( | ) |
ValueType& itk::OrientationDistributionFunction< TComponent, NOdfDirections >::operator() | ( | unsigned int | row, |
unsigned int | col | ||
) |
Matrix notation, in const and non-const forms.
const ValueType& itk::OrientationDistributionFunction< TComponent, NOdfDirections >::operator() | ( | unsigned int | row, |
unsigned int | col | ||
) | const |
Self itk::OrientationDistributionFunction< TComponent, NOdfDirections >::operator* | ( | const RealValueType & | scalar | ) | const |
Arithmetic operations between ODFs and scalars
const Self& itk::OrientationDistributionFunction< TComponent, NOdfDirections >::operator*= | ( | const RealValueType & | scalar | ) |
Self itk::OrientationDistributionFunction< TComponent, NOdfDirections >::operator+ | ( | const Self & | vec | ) | const |
Aritmetic operations between pixels. Return a new OrientationDistributionFunction.
const Self& itk::OrientationDistributionFunction< TComponent, NOdfDirections >::operator+= | ( | const Self & | vec | ) |
Self itk::OrientationDistributionFunction< TComponent, NOdfDirections >::operator- | ( | const Self & | vec | ) | const |
const Self& itk::OrientationDistributionFunction< TComponent, NOdfDirections >::operator-= | ( | const Self & | vec | ) |
Self itk::OrientationDistributionFunction< TComponent, NOdfDirections >::operator/ | ( | const RealValueType & | scalar | ) | const |
const Self& itk::OrientationDistributionFunction< TComponent, NOdfDirections >::operator/= | ( | const RealValueType & | scalar | ) |
Self& itk::OrientationDistributionFunction< TComponent, NOdfDirections >::operator= | ( | const Self & | r | ) |
Pass-through assignment operator for the Array base class.
Self& itk::OrientationDistributionFunction< TComponent, NOdfDirections >::operator= | ( | const ComponentType & | r | ) |
Self& itk::OrientationDistributionFunction< TComponent, NOdfDirections >::operator= | ( | const ComponentArrayType | r | ) |
Self itk::OrientationDistributionFunction< TComponent, NOdfDirections >::PostMultiply | ( | const MatrixType & | m | ) | const |
Post-Multiply by a Matrix as ResultingTensor = ThisTensor * Matrix.
Self itk::OrientationDistributionFunction< TComponent, NOdfDirections >::PreMultiply | ( | const MatrixType & | m | ) | const |
Pre-Multiply by a Matrix as ResultingTensor = Matrix * ThisTensor.
void itk::OrientationDistributionFunction< TComponent, NOdfDirections >::SetIsotropic | ( | ) |
Set the distribution to isotropic.
|
inline |
Set the Nth component to v.
Definition at line 130 of file itkOrientationDistributionFunction.h.