#include <itkTensorDerivedMeasurementsFilter.h>
Public Types | |
| enum | Measure { FA, RA, AD, RD, CA, L2, L3 } |
| typedef itk::DiffusionTensor3D < TPixel > | TensorType |
| typedef itk::Image< TensorType, 3 > | TensorImageType |
| typedef itk::Image< TPixel, 3 > | OutputImageType |
| typedef TensorDerivedMeasurementsFilter | Self |
| typedef SmartPointer< Self > | Pointer |
| typedef SmartPointer< const Self > | ConstPointer |
| typedef ImageToImageFilter < TensorType, OutputImageType > | SuperClass |
Public Member Functions | |
| virtual const char * | GetClassName () const |
| virtual void | SetMeasure (Measure _arg) |
Static Public Member Functions | |
| static Pointer | New () |
Protected Member Functions | |
| TensorDerivedMeasurementsFilter () | |
| ~TensorDerivedMeasurementsFilter () | |
| void | GenerateData () |
Protected Attributes | |
| Measure | m_Measure |
Definition at line 30 of file itkTensorDerivedMeasurementsFilter.h.
| typedef SmartPointer<const Self> itk::TensorDerivedMeasurementsFilter< TPixel >::ConstPointer |
Definition at line 54 of file itkTensorDerivedMeasurementsFilter.h.
| typedef itk::Image< TPixel, 3 > itk::TensorDerivedMeasurementsFilter< TPixel >::OutputImageType |
Definition at line 50 of file itkTensorDerivedMeasurementsFilter.h.
| typedef SmartPointer<Self> itk::TensorDerivedMeasurementsFilter< TPixel >::Pointer |
Definition at line 53 of file itkTensorDerivedMeasurementsFilter.h.
| typedef TensorDerivedMeasurementsFilter itk::TensorDerivedMeasurementsFilter< TPixel >::Self |
Definition at line 52 of file itkTensorDerivedMeasurementsFilter.h.
| typedef ImageToImageFilter< TensorType, OutputImageType > itk::TensorDerivedMeasurementsFilter< TPixel >::SuperClass |
Definition at line 57 of file itkTensorDerivedMeasurementsFilter.h.
| typedef itk::Image< TensorType, 3 > itk::TensorDerivedMeasurementsFilter< TPixel >::TensorImageType |
Definition at line 49 of file itkTensorDerivedMeasurementsFilter.h.
| typedef itk::DiffusionTensor3D<TPixel> itk::TensorDerivedMeasurementsFilter< TPixel >::TensorType |
Definition at line 48 of file itkTensorDerivedMeasurementsFilter.h.
| enum itk::TensorDerivedMeasurementsFilter::Measure |
| itk::TensorDerivedMeasurementsFilter< TPixel >::TensorDerivedMeasurementsFilter | ( | ) | [protected] |
Definition at line 12 of file itkTensorDerivedMeasurementsFilter.txx.
| itk::TensorDerivedMeasurementsFilter< TPixel >::~TensorDerivedMeasurementsFilter | ( | ) | [inline, protected] |
Definition at line 70 of file itkTensorDerivedMeasurementsFilter.h.
{};
| void itk::TensorDerivedMeasurementsFilter< TPixel >::GenerateData | ( | ) | [protected] |
Definition at line 17 of file itkTensorDerivedMeasurementsFilter.txx.
{
typename TensorImageType::Pointer tensorImage = static_cast< TensorImageType * >( this->ProcessObject::GetInput(0) );
typedef ImageRegionConstIterator< TensorImageType > TensorImageIteratorType;
typedef ImageRegionIterator< OutputImageType > OutputImageIteratorType;
typename OutputImageType::Pointer outputImage =
static_cast< OutputImageType * >(this->ProcessObject::GetOutput(0));
typename TensorImageType::RegionType region = tensorImage->GetLargestPossibleRegion();
outputImage->SetSpacing( tensorImage->GetSpacing() ); // Set the image spacing
outputImage->SetOrigin( tensorImage->GetOrigin() ); // Set the image origin
outputImage->SetDirection( tensorImage->GetDirection() ); // Set the image direction
outputImage->SetRegions( tensorImage->GetLargestPossibleRegion());
outputImage->Allocate();
TensorImageIteratorType tensorIt(tensorImage, tensorImage->GetLargestPossibleRegion());
OutputImageIteratorType outputIt(outputImage, outputImage->GetLargestPossibleRegion());
tensorIt.GoToBegin();
outputIt.GoToBegin();
while(!tensorIt.IsAtEnd() && !outputIt.IsAtEnd()){
TensorType tensor = tensorIt.Get();
switch(m_Measure)
{
case FA:
{
TPixel diffusionIndex = tensor.GetFractionalAnisotropy();
outputIt.Set(diffusionIndex);
break;
}
case RA:
{
TPixel diffusionIndex = tensor.GetRelativeAnisotropy();
outputIt.Set(diffusionIndex);
break;
}
case AD:
{
// eigenvalues are sorted in ascending order by default because the
// itk::SymmetricEigenAnalysis defaults are not touched in the tensor implementation
typename TensorType::EigenValuesArrayType evs;
tensor.ComputeEigenValues(evs);
outputIt.Set(evs[2]);
break;
}
case RD:
{
// eigenvalues are sorted in ascending order by default because the
// itk::SymmetricEigenAnalysis defaults are not touched in the tensor implementation
typename TensorType::EigenValuesArrayType evs;
tensor.ComputeEigenValues(evs);
outputIt.Set((evs[0]+evs[1])/2);
break;
}
case CA:
{
// eigenvalues are sorted in ascending order by default because the
// itk::SymmetricEigenAnalysis defaults are not touched in the tensor implementation
typename TensorType::EigenValuesArrayType evs;
tensor.ComputeEigenValues(evs);
if (evs[2] == 0)
{
outputIt.Set(0);
break;
}
outputIt.Set(1-(evs[0]+evs[1])/(2*evs[2]));
break;
}
case L2:
{
// eigenvalues are sorted in ascending order by default because the
// itk::SymmetricEigenAnalysis defaults are not touched in the tensor implementation
typename TensorType::EigenValuesArrayType evs;
tensor.ComputeEigenValues(evs);
outputIt.Set(evs[1]);
break;
}
case L3:
{
// eigenvalues are sorted in ascending order by default because the
// itk::SymmetricEigenAnalysis defaults are not touched in the tensor implementation
typename TensorType::EigenValuesArrayType evs;
tensor.ComputeEigenValues(evs);
outputIt.Set(evs[0]);
break;
}
}
++tensorIt;
++outputIt;
}
}
| virtual const char* itk::TensorDerivedMeasurementsFilter< TPixel >::GetClassName | ( | ) | const [virtual] |
Runtime information support.
| static Pointer itk::TensorDerivedMeasurementsFilter< TPixel >::New | ( | ) | [static] |
Method for creation through the object factory.
| virtual void itk::TensorDerivedMeasurementsFilter< TPixel >::SetMeasure | ( | Measure | _arg ) | [virtual] |
Measure itk::TensorDerivedMeasurementsFilter< TPixel >::m_Measure [protected] |
Definition at line 70 of file itkTensorDerivedMeasurementsFilter.h.
1.7.2