ANIMA  4.0
animaMultiCompartmentModel.h
Go to the documentation of this file.
1 #pragma once
2 
3 #include <animaBaseCompartment.h>
4 #include <AnimaMCMBaseExport.h>
5 
6 #include <itkLightObject.h>
7 #include <itkObjectFactory.h>
8 
9 namespace anima
10 {
11 
19 class ANIMAMCMBASE_EXPORT MultiCompartmentModel : public itk::LightObject
20 {
21 public:
22  // Useful typedefs
24  typedef itk::LightObject Superclass;
25  typedef itk::SmartPointer<Self> Pointer;
26  typedef itk::SmartPointer<const Self> ConstPointer;
27 
29  itkTypeMacro(MultiCompartmentModel, itk::LightObject)
30 
31  itkNewMacro(Self)
32 
37 
42  itk::LightObject::Pointer InternalClone() const ITK_OVERRIDE;
43 
44  void AddCompartment(double weight, BaseCompartment *compartment);
45 
46  unsigned int GetNumberOfCompartments() {return m_Compartments.size();}
47  itkGetMacro(NumberOfIsotropicCompartments, unsigned int)
48 
49  BaseCompartment *GetCompartment(unsigned int i);
50 
51  ListType &GetParametersAsVector();
52  void SetParametersFromVector(ListType &params);
53 
54  const ListType &GetCompartmentWeights() {return m_CompartmentWeights;}
55  double GetCompartmentWeight(unsigned int i);
56  void SetCompartmentWeights(const ListType &weights);
57 
58  ModelOutputVectorType &GetModelVector();
59 
60  void SetModelVector(const itk::VariableLengthVector <float> &mcmVec);
61  void SetModelVector(const ModelOutputVectorType &mcmVec);
62 
63  double GetPredictedSignal(double smallDelta, double bigDelta, double gradientStrength, const Vector3DType &gradient);
64  ListType &GetSignalJacobian(double smallDelta, double bigDelta, double gradientStrength, const Vector3DType &gradient);
65  double GetDiffusionProfile(Vector3DType &sample);
66 
67  ListType &GetParameterLowerBounds();
68  ListType &GetParameterUpperBounds();
69 
71  unsigned int GetNumberOfParameters();
72 
74  unsigned int GetNumberOfOptimizedWeights();
75 
77  unsigned int GetSize();
78 
80  void Reorient(vnl_matrix <double> &orientationMatrix, bool affineTransform);
81 
82  void SetOptimizeWeights(bool val) {m_OptimizeWeights = val;}
83  itkGetMacro(OptimizeWeights, bool)
84 
85  void SetCommonDiffusivityParameters(bool val) {m_CommonDiffusivityParameters = val;}
86  itkGetMacro(CommonDiffusivityParameters, bool)
87 
88  void SetCommonConcentrationParameters(bool val) {m_CommonConcentrationParameters = val;}
89  itkGetMacro(CommonConcentrationParameters, bool)
90 
91  void SetCommonExtraAxonalFractionParameters(bool val) {m_CommonExtraAxonalFractionParameters = val;}
92  itkGetMacro(CommonExtraAxonalFractionParameters, bool)
93 
94  void SetNegativeWeightBounds(bool val) {m_NegativeWeightBounds = val;}
95  itkGetMacro(NegativeWeightBounds, bool)
96 
97 protected:
100 
101 private:
102  std::vector <BaseCompartmentPointer> m_Compartments;
103  ListType m_CompartmentWeights;
104 
105  bool m_OptimizeWeights;
106  bool m_CommonDiffusivityParameters;
107  bool m_CommonConcentrationParameters;
108  bool m_CommonExtraAxonalFractionParameters;
109 
110  unsigned int m_NumberOfIsotropicCompartments;
111 
113  bool m_NegativeWeightBounds;
114 
116  ListType m_JacobianVector;
117 
119  ListType m_ParametersVector;
120 
122  ListType m_WorkVector;
123 
125  ListType m_ParametersLowerBoundsVector;
126 
128  ListType m_ParametersUpperBoundsVector;
129 
131  ModelOutputVectorType m_ModelVector;
132 };
133 
134 } // end namespace anima
BaseCompartment::Pointer BaseCompartmentPointer
BaseCompartment::ListType ListType
itk::SmartPointer< const Self > ConstPointer
BaseCompartment::Vector3DType Vector3DType
MultiCompartmentModel: holds several diffusion compartments, ordered by type It also handles weights ...
BaseCompartment::ModelOutputVectorType ModelOutputVectorType