ANIMA  4.0
animaBaseTractographyImageFilter.h
Go to the documentation of this file.
1 #pragma once
2 
3 #include <itkVectorImage.h>
4 #include <itkImage.h>
5 
6 #include <vtkPolyData.h>
7 #include <vtkSmartPointer.h>
8 #include <itkLinearInterpolateImageFunction.h>
9 #include <itkProcessObject.h>
10 #include <mutex>
11 #include <itkProgressReporter.h>
12 
13 #include "AnimaTractographyExport.h"
14 
15 #include <vector>
16 
17 namespace anima
18 {
19 
20 class ANIMATRACTOGRAPHY_EXPORT BaseTractographyImageFilter : public itk::ProcessObject
21 {
22 public:
25  typedef itk::ProcessObject Superclass;
26 
27  typedef itk::SmartPointer<Self> Pointer;
28  typedef itk::SmartPointer<const Self> ConstPointer;
29 
30  itkTypeMacro(BaseTractographyImageFilter,itk::ProcessObject)
31 
32  typedef enum {
35  Both
37 
38  typedef itk::VectorImage <double, 3> ModelImageType;
39  typedef ModelImageType::Pointer ModelImagePointer;
40  typedef ModelImageType::PixelType VectorType;
41  typedef ModelImageType::RegionType RegionType;
42  typedef itk::ContinuousIndex <double, 3> ContinuousIndexType;
43 
44  typedef itk::Image <unsigned short, 3> MaskImageType;
45  typedef MaskImageType::Pointer MaskImagePointer;
46  typedef MaskImageType::PointType PointType;
47  typedef MaskImageType::IndexType IndexType;
48 
49  typedef std::vector <PointType> FiberType;
50  typedef std::vector <std::pair < FiberProgressType, FiberType > > FiberProcessVectorType;
51 
52  typedef struct {
54  std::vector < std::vector <FiberType> > resultFibersFromThreads;
56 
57  virtual void SetInputImage(ModelImageType *input) {m_InputImage = input;}
58  ModelImageType *GetInputImage() {return m_InputImage;}
59 
60  void SetSeedingMask(MaskImageType *mask) {m_SeedingImage = mask;}
61  void SetFilteringMask(MaskImageType *mask) {m_FilteringImage = mask;}
62  void SetForbiddenMask(MaskImageType *mask) {m_ForbiddenMaskImage = mask;}
63  void SetCutMask(MaskImageType *mask) {m_CutMaskImage = mask;}
64 
65  void SetNumberOfFibersPerPixel(unsigned int num) {m_NumberOfFibersPerPixel = num;}
66 
67  void SetStepProgression(double num) {m_StepProgression = num;}
68 
69  void SetMaxFiberAngle(double num) {m_MaxFiberAngle = num;}
70  void SetMinLengthFiber(double num) {m_MinLengthFiber = num;}
71  void SetMaxLengthFiber(double num) {m_MaxLengthFiber = num;}
72 
73  void Update() ITK_OVERRIDE;
74 
75  void SetComputeLocalColors(bool flag) {m_ComputeLocalColors = flag;}
76  void createVTKOutput(std::vector < std::vector <PointType> > &filteredFibers);
77  vtkPolyData *GetOutput() {return m_Output;}
78 
79 protected:
81  virtual ~BaseTractographyImageFilter();
82 
83  static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION ThreadTracker(void *arg);
84  void ThreadTrack(unsigned int numThread, std::vector <FiberType> &resultFibers);
85  void ThreadedTrackComputer(unsigned int numThread, std::vector <FiberType> &resultFibers,
86  unsigned int startSeedIndex, unsigned int endSeedIndex);
87 
88  FiberProcessVectorType ComputeFiber(FiberType &fiber, FiberProgressType ways, itk::ThreadIdType threadId);
89 
90  virtual void PrepareTractography();
91  std::vector < FiberType > FilterOutputFibers(std::vector < FiberType > &fibers);
92 
93  virtual bool CheckModelCompatibility(VectorType &modelValue, itk::ThreadIdType threadId) = 0;
94 
95  typedef itk::ImageBase <ModelImageType::ImageDimension> ImageBaseType;
96  bool CheckIndexInImageBounds(IndexType &index, ImageBaseType *testImage);
97  virtual bool CheckIndexInImageBounds(ContinuousIndexType &index) = 0;
98 
100  virtual void GetModelValue(ContinuousIndexType &index, VectorType &modelValue) = 0;
101 
102  virtual PointType GetModelPrincipalDirection(VectorType &modelValue, bool is2d, itk::ThreadIdType threadId) = 0;
103  virtual PointType GetNextDirection(PointType &previousDirection, VectorType &modelValue, bool is2d, itk::ThreadIdType threadId) = 0;
104 
105  virtual void ComputeAdditionalScalarMaps() {}
106  bool isZero(VectorType &value);
107 
108 private:
109  ITK_DISALLOW_COPY_AND_ASSIGN(BaseTractographyImageFilter);
110 
111  unsigned int m_NumberOfFibersPerPixel;
112 
113  double m_StepProgression;
114  double m_MaxFiberAngle;
115 
116  double m_MinLengthFiber;
117  double m_MaxLengthFiber;
118 
119  ModelImagePointer m_InputImage;
120  MaskImagePointer m_SeedingImage, m_FilteringImage, m_ForbiddenMaskImage, m_CutMaskImage;
121 
122  FiberProcessVectorType m_PointsToProcess;
123  std::vector <unsigned int> m_FilteringValues;
124 
125  bool m_ComputeLocalColors;
126  vtkSmartPointer<vtkPolyData> m_Output;
127 
128  std::mutex m_LockHighestProcessedSeed;
129  int m_HighestProcessedSeed;
130  itk::ProgressReporter *m_ProgressReport;
131 };
132 
133 } // end of namespace anima
itk::ImageBase< ModelImageType::ImageDimension > ImageBaseType
virtual void SetInputImage(ModelImageType *input)
std::vector< std::pair< FiberProgressType, FiberType > > FiberProcessVectorType
itk::ContinuousIndex< double, 3 > ContinuousIndexType