ANIMA  4.0
animaMultiTensorSmoothingCostFunction.h
Go to the documentation of this file.
1 #pragma once
2 
3 #include <itkSingleValuedCostFunction.h>
4 #include <AnimaMCMSimilarityExport.h>
6 
7 namespace anima
8 {
9 
10 class ANIMAMCMSIMILARITY_EXPORT MultiTensorSmoothingCostFunction : public itk::SingleValuedCostFunction
11 {
12 public:
15  typedef itk::SingleValuedCostFunction Superclass;
16  typedef itk::SmartPointer<Self> Pointer;
17  typedef itk::SmartPointer<const Self> ConstPointer;
18 
19  itkNewMacro(Self)
20 
21 
22  itkTypeMacro(MultiTensorSmoothingCostFunction, itk::SingleValuedCostFunction)
23 
24  typedef Superclass::MeasureType MeasureType;
25  typedef Superclass::ParametersType ParametersType;
26  typedef Superclass::DerivativeType DerivativeType;
27 
29  typedef MCModelType::Pointer MCMPointer;
30  typedef anima::BaseCompartment::Matrix3DType TensorType;
31 
32  virtual MeasureType GetValue(const ParametersType &parameters) const ITK_OVERRIDE;
33  virtual void GetDerivative(const ParametersType &parameters, DerivativeType &derivative) const ITK_OVERRIDE;
34 
35  void SetReferenceModels(const std::vector <MCMPointer> &refModels);
36  void SetMovingModels(const std::vector <MCMPointer> &movingModels);
37 
38  itkSetMacro(TensorsScale, double)
39  itkSetMacro(LowPassGaussianSigma, double)
40 
41  unsigned int GetNumberOfParameters() const ITK_OVERRIDE {return 1;}
42 
43 protected:
45  {
46  m_UpdatedReferenceData = false;
47  m_UpdatedMovingData = false;
48  m_RecomputeConstantTerm = true;
49 
50  m_ConstantTerm = 0;
51  m_TensorsScale = 1000;
52  m_LowPassGaussianSigma = 2000;
53  }
54 
56 
57  void UpdateDeterminants() const;
58  void ComputeMatrixTraces(const TensorType &matrix, double &matrixTrace, double &squaredMatrixTrace) const;
59 
60 private:
61  MultiTensorSmoothingCostFunction(const Self&); //purposely not implemented
62  void operator=(const Self&); //purposely not implemented
63 
64  std::vector < std::vector <TensorType> > m_ReferenceModels;
65  std::vector < std::vector <double> > m_ReferenceModelWeights;
66  std::vector <double> m_ReferenceNumberOfIsotropicCompartments;
67  std::vector < std::vector <TensorType> > m_MovingModels;
68  std::vector < std::vector <double> > m_MovingModelWeights;
69 
70  mutable std::vector < std::vector <double> > m_ReferenceReferenceDeterminants, m_ReferenceMovingDeterminants, m_MovingMovingDeterminants;
71  mutable std::vector < std::vector <double> > m_ReferenceReferenceTraces, m_ReferenceMovingTraces;
72  mutable std::vector < std::vector <double> > m_ReferenceReferenceTraceDifferences, m_ReferenceMovingTraceDifferences;
73 
74  mutable bool m_UpdatedReferenceData;
75  mutable bool m_UpdatedMovingData;
76  mutable bool m_RecomputeConstantTerm;
77 
78  mutable double m_ConstantTerm;
79  double m_TensorsScale;
80  double m_LowPassGaussianSigma;
81 };
82 
83 } // end namespace anima
STL namespace.
MultiCompartmentModel: holds several diffusion compartments, ordered by type It also handles weights ...