22 #ifndef __itkStreamlineTrackingFilter_h_
23 #define __itkStreamlineTrackingFilter_h_
25 #include <MitkFiberTrackingExports.h>
26 #include <itkImageToImageFilter.h>
27 #include <itkVectorContainer.h>
28 #include <itkVectorImage.h>
29 #include <itkDiffusionTensor3D.h>
30 #include <vtkSmartPointer.h>
31 #include <vtkPolyData.h>
32 #include <vtkCellArray.h>
33 #include <vtkPoints.h>
34 #include <vtkPolyLine.h>
41 template<
class TTensorPixelType,
class TPDPixelType=
double>
43 public ImageToImageFilter< Image< DiffusionTensor3D<TTensorPixelType>, 3 >,
44 Image< Vector< TPDPixelType, 3 >, 3 > >
52 typedef ImageToImageFilter< Image< DiffusionTensor3D<TTensorPixelType>, 3 >, Image< Vector< TPDPixelType, 3 >, 3 > >
Superclass;
55 itkFactorylessNewMacro(Self)
72 itkGetMacro( FiberPolyData, FiberPolyDataType )
73 itkSetMacro(
SeedImage, ItkUcharImgType::Pointer)
74 itkSetMacro( MaskImage, ItkUcharImgType::Pointer)
75 itkSetMacro( FaImage, ItkFloatImgType::Pointer)
76 itkSetMacro( SeedsPerVoxel,
int)
77 itkSetMacro( FaThreshold,
double)
78 itkSetMacro( StepSize,
double)
79 itkSetMacro( F,
double )
80 itkSetMacro( G,
double )
81 itkSetMacro( Interpolate,
bool )
82 itkSetMacro( MinTractLength,
double )
83 itkGetMacro( MinTractLength,
double )
84 itkSetMacro( MinCurvatureRadius,
double )
85 itkGetMacro( MinCurvatureRadius,
double )
86 itkSetMacro( ResampleFibers,
bool )
89 StreamlineTrackingFilter();
90 ~StreamlineTrackingFilter() {}
91 void PrintSelf(std::ostream& os, Indent indent)
const;
93 void CalculateNewPosition(itk::ContinuousIndex<double, 3>& pos, vnl_vector_fixed<double,3>& dir,
typename InputImageType::IndexType& index);
94 double FollowStreamline(itk::ContinuousIndex<double, 3> pos,
int dirSign, vtkPoints* points, std::vector< vtkIdType >& ids,
int imageIdx);
95 bool IsValidPosition(itk::ContinuousIndex<double, 3>& pos,
typename InputImageType::IndexType& index, vnl_vector_fixed< double, 8 >& interpWeights,
int imageIdx);
138 #ifndef ITK_MANUAL_INSTANTIATION
142 #endif //__itkStreamlineTrackingFilter_h_
std::vector< typename InputImageType::Pointer > m_InputImage
Input tensor images. For multi tensor tracking provide multiple tensor images.
itk::SmartPointer< Self > Pointer
Performes deterministic streamline tracking on the input tensor image.
itk::VectorContainer< int, FiberPolyDataType >::Pointer m_PolyDataContainer
void CalculateNewPosition(itk::ContinuousIndex< double, 3 > &pos, vnl_vector_fixed< double, 3 > &dir, typename InputImageType::IndexType &index)
Calculate next integration step.
itk::Image< unsigned char, 3 > SeedImage
double FollowStreamline(itk::ContinuousIndex< double, 3 > pos, int dirSign, vtkPoints *points, std::vector< vtkIdType > &ids, int imageIdx)
Start streamline in one direction.
vtkSmartPointer< vtkCellArray > m_Cells
double RoundToNearest(double num)
itk::Image< double, 3 > ItkDoubleImgType
std::vector< ItkPDImgType::Pointer > m_PdImage
Stores principal direction of each tensor in each voxel.
StreamlineTrackingFilter Self
TPDPixelType DirectionPixelType
ItkUcharImgType::Pointer m_MaskImage
vtkSmartPointer< vtkPoints > m_Points
TTensorPixelType TensorComponentType
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId)
ImageToImageFilter< Image< DiffusionTensor3D< TTensorPixelType >, 3 >, Image< Vector< TPDPixelType, 3 >, 3 > > Superclass
FiberPolyDataType AddPolyData(FiberPolyDataType poly1, FiberPolyDataType poly2)
Combine tracking results generated by the individual threads.
void PrintSelf(std::ostream &os, Indent indent) const
ItkFloatImgType::Pointer m_FaImage
FA image used to determine streamline termination.
SmartPointer< Self > Pointer
ItkUcharImgType::Pointer m_SeedImage
Superclass::OutputImageType OutputImageType
Superclass::InputImageType InputImageType
itk::Image< unsigned char, 3 > ItkUcharImgType
void BeforeThreadedGenerateData()
SmartPointer< const Self > ConstPointer
std::vector< int > m_ImageSize
vtkSmartPointer< vtkPolyData > FiberPolyDataType
void AfterThreadedGenerateData()
itk::Image< vnl_vector_fixed< double, 3 >, 3 > ItkPDImgType
itk::Image< float, 3 > ItkFloatImgType
double m_MinCurvatureRadius
bool IsValidPosition(itk::ContinuousIndex< double, 3 > &pos, typename InputImageType::IndexType &index, vnl_vector_fixed< double, 8 > &interpWeights, int imageIdx)
Are we outside of the mask image? Is the FA too low?
std::vector< double > m_ImageSpacing
std::vector< ItkDoubleImgType::Pointer > m_EmaxImage
Stores largest eigenvalues per voxel (one for each tensor)
FiberPolyDataType m_FiberPolyData
Superclass::OutputImageRegionType OutputImageRegionType