ANIMA  4.0
animaTensorCompartment.h
Go to the documentation of this file.
1 #pragma once
2 
3 #include <animaBaseCompartment.h>
5 #include <animaBaseTensorTools.h>
6 #include <AnimaMCMExport.h>
7 
8 namespace anima
9 {
10 
11 class ANIMAMCM_EXPORT TensorCompartment : public BaseCompartment
12 {
13 public:
14  // Useful typedefs
18  typedef itk::SmartPointer<Self> Pointer;
19  typedef itk::SmartPointer<const Self> ConstPointer;
25 
26  // New macro
27  itkNewMacro(Self)
28 
29 
30  itkTypeMacro(TensorCompartment, BaseCompartment)
31 
32  DiffusionModelCompartmentType GetCompartmentType() ITK_OVERRIDE {return Tensor;}
33 
34  virtual double GetFourierTransformedDiffusionProfile(double smallDelta, double bigDelta, double gradientStrength, const Vector3DType &gradient) ITK_OVERRIDE;
35  virtual ListType &GetSignalAttenuationJacobian(double smallDelta, double bigDelta, double gradientStrength, const Vector3DType &gradient) ITK_OVERRIDE;
36  virtual double GetLogDiffusionProfile(const Vector3DType &sample) ITK_OVERRIDE;
37 
38  virtual void SetParametersFromVector(const ListType &params) ITK_OVERRIDE;
39  virtual ListType &GetParametersAsVector() ITK_OVERRIDE;
40 
41  virtual ListType &GetParameterLowerBounds() ITK_OVERRIDE;
42  virtual ListType &GetParameterUpperBounds() ITK_OVERRIDE;
43 
44  // Set constraints
45  void SetEstimateDiffusivities(bool arg);
46 
47  void SetCompartmentVector(ModelOutputVectorType &compartmentVector) ITK_OVERRIDE;
48  void UpdateParametersFromTensor();
49 
50  unsigned int GetCompartmentSize() ITK_OVERRIDE;
51  unsigned int GetNumberOfParameters() ITK_OVERRIDE;
52  ModelOutputVectorType &GetCompartmentVector() ITK_OVERRIDE;
53 
54  void Reorient(vnl_matrix <double> &orientationMatrix, bool affineTransform) ITK_OVERRIDE;
55 
56  //Reimplement for handling modification flags
57  void SetOrientationTheta(double num) ITK_OVERRIDE;
58  void SetOrientationPhi(double num) ITK_OVERRIDE;
59  void SetPerpendicularAngle(double num) ITK_OVERRIDE;
60  void SetAxialDiffusivity(double num) ITK_OVERRIDE;
61  void SetRadialDiffusivity1(double num) ITK_OVERRIDE;
62  void SetRadialDiffusivity2(double num) ITK_OVERRIDE;
63  double GetOrientationTheta() ITK_OVERRIDE;
64  double GetOrientationPhi() ITK_OVERRIDE;
65  double GetPerpendicularAngle() ITK_OVERRIDE;
66  double GetAxialDiffusivity() ITK_OVERRIDE;
67  double GetRadialDiffusivity1() ITK_OVERRIDE;
68  double GetRadialDiffusivity2() ITK_OVERRIDE;
69 
70  const Matrix3DType &GetDiffusionTensor() ITK_OVERRIDE;
71  double GetApparentFractionalAnisotropy() ITK_OVERRIDE;
72  double GetApparentMeanDiffusivity() ITK_OVERRIDE;
73  double GetApparentPerpendicularDiffusivity() ITK_OVERRIDE;
74  double GetApparentParallelDiffusivity() ITK_OVERRIDE;
75 
76 protected:
77  TensorCompartment() : Superclass()
78  {
79  m_EstimateDiffusivities = true;
80  m_ChangedConstraints = true;
81 
82  m_ModifiedTensor = true;
83  m_UpdateInverseTensor = true;
84  m_ModifiedAngles = true;
85  m_UpdatedCompartment = false;
86  m_TensorDeterminant = 0;
87 
88  m_leCalculator = LEcalculatorType::New();
89  }
90 
91  virtual ~TensorCompartment() {}
92 
94  void UpdateDiffusionTensor();
95 
97  void UpdateInverseDiffusionTensor();
98 
100  void UpdateAngleConfiguration();
101 
102 private:
103  bool m_EstimateDiffusivities;
104  bool m_ChangedConstraints;
105  unsigned int m_NumberOfParameters;
106 
107  // Internal work variables for faster processing
108 
110  bool m_ModifiedTensor;
111 
113  bool m_UpdateInverseTensor;
114 
116  bool m_ModifiedAngles;
117 
119  bool m_UpdatedCompartment;
120 
121  vnl_matrix <double> m_WorkVnlMatrix1, m_WorkVnlMatrix2;
122  Matrix3DType m_DiffusionTensor;
123  Matrix3DType m_InverseDiffusionTensor;
124  Vector3DType m_EigenVector1, m_EigenVector2;
125  double m_SinTheta, m_CosTheta, m_SinPhi, m_CosPhi, m_SinAlpha, m_CosAlpha;
126  double m_TensorDeterminant;
127 
128  LEcalculatorPointer m_leCalculator;
129 };
130 
131 } //end namespace anima
132 
std::vector< double > ListType
Superclass::ModelOutputVectorType ModelOutputVectorType
Superclass::Vector3DType Vector3DType
itk::VariableLengthVector< double > ModelOutputVectorType
itk::Matrix< double, 3, 3 > Matrix3DType
vnl_vector_fixed< double, 3 > Vector3DType
LEcalculatorType::Pointer LEcalculatorPointer
DiffusionModelCompartmentType
itk::SmartPointer< Self > Pointer
itk::SmartPointer< const Self > ConstPointer
Superclass::Pointer BasePointer
Superclass::Matrix3DType Matrix3DType
itk::SmartPointer< Self > Pointer