ANIMA  4.0
animaDWISimulatorFromDTIImageFilter.h
Go to the documentation of this file.
1 #pragma once
2 
3 #include <itkImageToImageFilter.h>
4 #include <itkVectorImage.h>
5 #include <itkImage.h>
6 
7 namespace anima
8 {
9 
10 template <class PixelScalarType>
12  public itk::ImageToImageFilter < itk::VectorImage<PixelScalarType,3>, itk::VectorImage<PixelScalarType,3> >
13 {
14 public:
17  typedef itk::VectorImage<PixelScalarType,3> TInputImage;
18  typedef TInputImage OutputB0ImageType;
19  typedef itk::VectorImage<PixelScalarType,3> DTIImageType;
20  typedef itk::Image<PixelScalarType,4> Image4DType;
21  typedef itk::Image<PixelScalarType,3> S0ImageType;
22  typedef typename S0ImageType::Pointer S0ImagePointer;
23  typedef itk::VectorImage<PixelScalarType,3> TOutputImage;
24  typedef itk::ImageToImageFilter< TInputImage, TOutputImage > Superclass;
25  typedef itk::SmartPointer<Self> Pointer;
26  typedef itk::SmartPointer<const Self> ConstPointer;
27 
29  itkNewMacro(Self)
30 
31 
32  itkTypeMacro(DWISimulatorFromDTIImageFilter, ImageToImageFilter)
33 
34 
35  typedef TInputImage InputImageType;
36  typedef typename InputImageType::PixelType InputPixelType;
37  typedef TOutputImage OutputImageType;
38  typedef typename InputImageType::Pointer InputImagePointer;
39  typedef typename OutputImageType::Pointer OutputImagePointer;
40 
42  typedef typename Superclass::InputImageRegionType InputImageRegionType;
43  typedef typename Superclass::OutputImageRegionType OutputImageRegionType;
44 
45  void SetBValuesList (std::vector <double> bValuesList) {m_BValuesList = bValuesList;}
46 
47  void AddGradientDirection(unsigned int i, std::vector <double> &grad);
48 
49  itkGetMacro(S0Value, double)
50  itkSetMacro(S0Value, double)
51 
52  itkSetMacro(S0Image, S0ImagePointer)
53 
54  Image4DType *GetOutputAs4DImage();
55 
56 protected:
58  : Superclass()
59  {
60  m_BValuesList.clear();
61  m_GradientDirections.clear();
62 
63  m_S0Value = 200;
64  m_S0Image = 0;
65  }
66 
68 
69  void GenerateOutputInformation() ITK_OVERRIDE;
70  void BeforeThreadedGenerateData() ITK_OVERRIDE;
71  void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) ITK_OVERRIDE;
72 
73  bool isZero(InputPixelType &vec);
74 
75 private:
76  ITK_DISALLOW_COPY_AND_ASSIGN(DWISimulatorFromDTIImageFilter);
77 
78  std::vector <double> m_BValuesList;
79  std::vector < std::vector <double> > m_GradientDirections;
80 
81  double m_S0Value;
82  S0ImagePointer m_S0Image;
83 
84  static const unsigned int m_NumberOfComponents = 6;
85 
86  typename Image4DType::Pointer m_Output4D;
87 };
88 
89 } // end namespace anima
90 
Superclass::OutputImageRegionType OutputImageRegionType
itk::VectorImage< PixelScalarType, 3 > TOutputImage
void SetBValuesList(std::vector< double > bValuesList)
STL namespace.
void AddGradientDirection(unsigned int i, std::vector< double > &grad)
DWISimulatorFromDTIImageFilter< PixelScalarType > Self
itk::VectorImage< PixelScalarType, 3 > TInputImage
itk::VectorImage< PixelScalarType, 3 > DTIImageType
itk::ImageToImageFilter< TInputImage, TOutputImage > Superclass
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) ITK_OVERRIDE