ANIMA  4.0
animaMCMLinearInterpolateImageFunction.h
Go to the documentation of this file.
1 #pragma once
2 
3 #include <itkInterpolateImageFunction.h>
5 #include <mutex>
7 
8 namespace anima
9 {
10 
11 template <class TInputImage, class TCoordRep = double>
13 public itk::InterpolateImageFunction<TInputImage,TCoordRep>
14 {
15 public:
18  typedef itk::InterpolateImageFunction<TInputImage,TCoordRep> Superclass;
19  typedef itk::SmartPointer<Self> Pointer;
20  typedef itk::SmartPointer<const Self> ConstPointer;
21 
23  itkNewMacro(Self)
24 
25 
27  InterpolateImageFunction)
28 
29 
30  typedef typename Superclass::InputImageType InputImageType;
31  typedef typename TInputImage::PixelType PixelType;
32  typedef typename Superclass::RealType RealType;
33  typedef typename Superclass::SizeType SizeType;
34 
36  itkStaticConstMacro(ImageDimension, unsigned int,Superclass::ImageDimension);
37 
39  typedef typename Superclass::IndexType IndexType;
40  typedef typename Superclass::IndexValueType IndexValueType;
41 
43  typedef typename Superclass::ContinuousIndexType ContinuousIndexType;
44 
45  typedef typename Superclass::OutputType OutputType;
46 
47  // Multi-compartment models typedefs
49  typedef MCModelType::Pointer MCModelPointer;
50 
52  typedef AveragerType::Pointer AveragerPointer;
53 
62  virtual OutputType EvaluateAtContinuousIndex(const ContinuousIndexType & index ) const ITK_OVERRIDE;
63 
64  void SetReferenceOutputModel(MCModelPointer &model);
65 
66  virtual void SetInputImage(const InputImageType *ptr) ITK_OVERRIDE;
67 
68  std::vector <AveragerPointer> &GetAveragers() const {return m_MCMAveragers;}
69 
70  SizeType GetRadius() const override
71  {
72  return SizeType::Filled(1);
73  }
74 
75 protected:
78 
79  template <class T> bool isZero(const itk::VariableLengthVector <T> &value) const
80  {
81  for (unsigned int i = 0;i < value.GetNumberOfElements();++i)
82  {
83  if (value[i] != 0)
84  return false;
85  }
86 
87  return true;
88  }
89 
90  // Fake method for compilation purposes, should never go in there
91  template <class T> bool isZero(T &data) const
92  {
93  itkExceptionMacro("Access to unauthorized method");
94  return true;
95  }
96 
98  void TestModelsAdequation(MCModelPointer &inputModel, MCModelPointer &outputModel);
99 
101  virtual void ResetAveragePointers(MCModelPointer &model);
102 
104  virtual bool CheckModelCompatibility(MCModelPointer &model);
105 
107  virtual void SetSpecificAveragerParameters(unsigned int threadIndex) const {}
108 
109  unsigned int GetFreeWorkIndex() const;
110  void UnlockWorkIndex(unsigned int index) const;
111 
112 private:
113  MCMLinearInterpolateImageFunction(const Self&); //purposely not implemented
114  void operator=(const Self&); //purposely not implemented
115 
117  static const unsigned long m_Neighbors;
118  static const unsigned int m_SphereDimension = 3;
119 
120  mutable std::mutex m_LockUsedModels;
121  mutable std::vector <int> m_UsedModels;
122 
123  mutable std::vector < std::vector <MCModelPointer> > m_ReferenceInputModels;
124  mutable std::vector < std::vector <double> > m_ReferenceInputWeights;
125  mutable std::vector <AveragerPointer> m_MCMAveragers;
126 };
127 
128 } // end namespace anima
129 
itk::InterpolateImageFunction< TInputImage, TCoordRep > Superclass
bool isZero(const itk::VariableLengthVector< T > &value) const
Computes a weighted average of input multi-compartment models. The output model is at the same time g...
STL namespace.
void TestModelsAdequation(MCModelPointer &inputModel, MCModelPointer &outputModel)
Tests if input and output models of the interpolator are compatible.
virtual void SetInputImage(const InputImageType *ptr) ITK_OVERRIDE
virtual OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &index) const ITK_OVERRIDE
virtual void SetSpecificAveragerParameters(unsigned int threadIndex) const
Sets averager specific parameters if sub-classes are derived.
MultiCompartmentModel: holds several diffusion compartments, ordered by type It also handles weights ...
virtual void ResetAveragePointers(MCModelPointer &model)
Resets averager pointers, can be overloaded to handle more models.
virtual bool CheckModelCompatibility(MCModelPointer &model)
Check if model can actually be interpolated.
std::vector< AveragerPointer > & GetAveragers() const