ANIMA  4.0
animaBaseBlockMatcher.h
Go to the documentation of this file.
1 #pragma once
3 
4 #include <itkSingleValuedNonLinearOptimizer.h>
5 #include <itkSingleValuedCostFunction.h>
6 #include <mutex>
7 
8 namespace anima
9 {
10 
11 template <typename TInputImageType>
13 {
14 public:
17  virtual ~BaseBlockMatcher() {}
18 
20  {
23  };
24 
25  typedef TInputImageType InputImageType;
26  typedef typename InputImageType::Pointer InputImagePointer;
27  typedef typename InputImageType::RegionType ImageRegionType;
28  typedef typename InputImageType::PointType PointType;
29 
33  typedef typename BaseInputTransformType::Pointer BaseInputTransformPointer;
34 
35  typedef itk::Image <unsigned char, TInputImageType::ImageDimension> MaskImageType;
36  typedef typename MaskImageType::Pointer MaskImagePointer;
37 
39  typedef itk::SingleValuedNonLinearOptimizer OptimizerType;
40  typedef typename OptimizerType::Pointer OptimizerPointer;
41 
43  typedef itk::SingleValuedCostFunction MetricType;
44  typedef typename MetricType::Pointer MetricPointer;
45 
46  void SetReferenceImage(InputImageType *image) {m_ReferenceImage = image;}
47  void SetMovingImage(InputImageType *image) {m_MovingImage = image;}
48  void SetBlockGenerationMask(MaskImageType *image) {m_BlockGenerationMask = image;}
49 
51 
52  void SetForceComputeBlocks(bool val) {m_ForceComputeBlocks = val;}
53  void SetNumberOfWorkUnits(unsigned int val) {m_NumberOfThreads = val;}
54  unsigned int GetNumberOfWorkUnits() {return m_NumberOfThreads;}
55 
56  void SetBlockVarianceThreshold(double val) {m_BlockVarianceThreshold = val;}
57  double GetBlockVarianceThreshold() {return m_BlockVarianceThreshold;}
58 
59  void SetBlockPercentageKept(double val) {m_BlockPercentageKept = val;}
60  double GetBlockPercentageKept() {return m_BlockPercentageKept;}
61 
62  void SetBlockSize(unsigned int val) {m_BlockSize = val;}
63  unsigned int GetBlockSize() {return m_BlockSize;}
64 
65  void SetBlockSpacing(unsigned int val) {m_BlockSpacing = val;}
66  unsigned int GetBlockSpacing() {return m_BlockSpacing;}
67 
68  void SetSearchRadius(double val) {m_SearchRadius = val;}
69  double GetSearchRadius() {return m_SearchRadius;}
70  void SetFinalRadius(double val) {m_FinalRadius = val;}
71  void SetStepSize (double val) {m_StepSize = val;}
72 
73  void SetOptimizerMaximumIterations (unsigned int val) {m_OptimizerMaximumIterations = val;}
74 
75  void Update();
76 
77  std::vector <PointType> &GetBlockPositions() {return m_BlockPositions;}
78  std::vector <ImageRegionType> &GetBlockRegions() {return m_BlockRegions;}
79  ImageRegionType &GetBlockRegion(unsigned int i) {return m_BlockRegions[i];}
80 
81  std::vector <BaseInputTransformPointer> &GetBlockTransformPointers() {return m_BlockTransformPointers;}
82  BaseInputTransformPointer &GetBlockTransformPointer(unsigned int i) {return m_BlockTransformPointers[i];}
83 
84  const std::vector <double> &GetBlockWeights() {return m_BlockWeights;}
85 
86  void SetOptimizerType(OptimizerDefinition val) {m_OptimizerType = val;}
87  OptimizerDefinition GetOptimizerType() {return m_OptimizerType;}
88 
89  InputImagePointer &GetReferenceImage() {return m_ReferenceImage;}
90  InputImagePointer &GetMovingImage() {return m_MovingImage;}
91  MaskImagePointer &GetBlockGenerationMask() {return m_BlockGenerationMask;}
92 
93  virtual bool GetMaximizedMetric() = 0;
94 
95  void SetVerbose(bool value) {m_Verbose = value;}
96  bool GetVerbose() {return m_Verbose;}
97 
98 protected:
100  {
101  Self *BlockMatch;
102  };
103 
105  static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION ThreadedMatching(void *arg);
106 
107  void ProcessBlockMatch();
108  void BlockMatch(unsigned int startIndex, unsigned int endIndex);
109 
110  virtual void InitializeBlocks();
111 
112  virtual MetricPointer SetupMetric() = 0;
113  virtual double ComputeBlockWeight(double val, unsigned int block) = 0;
114  virtual BaseInputTransformPointer GetNewBlockTransform(PointType &blockCenter) = 0;
115 
116  // May be overloaded but in practice, much easier if this superclass implementation is always called
117  virtual OptimizerPointer SetupOptimizer();
118 
119  virtual void BlockMatchingSetup(MetricPointer &metric, unsigned int block) = 0;
120  virtual void TransformDependantOptimizerSetup(OptimizerPointer &optimizer) = 0;
121 
122  // Internal setters for re-implementations of block initialization
123  void SetBlockWeights(std::vector <double> &val) {m_BlockWeights = val;}
124  void SetBlockRegions(std::vector <ImageRegionType> &val) {m_BlockRegions = val;}
125  void SetBlockPositions(std::vector <PointType> &val) {m_BlockPositions = val;}
126 
127 private:
128  InputImagePointer m_ReferenceImage;
129  InputImagePointer m_MovingImage;
130  MaskImagePointer m_BlockGenerationMask;
131 
132  bool m_ForceComputeBlocks;
133  unsigned int m_NumberOfThreads;
134 
135  // The origins of the blocks
136  std::vector <PointType> m_BlockPositions;
137  std::vector <ImageRegionType> m_BlockRegions;
138 
139  // The weights associated to blocks
140  std::vector <double> m_BlockWeights;
141 
142  // Parameters fo block creation
143  double m_BlockVarianceThreshold;
144  double m_BlockPercentageKept;
145  unsigned int m_BlockSize;
146  unsigned int m_BlockSpacing;
147 
148  bool m_Verbose;
149 
150  // The matched transform for each of the BlockRegions
151  std::vector <BaseInputTransformPointer> m_BlockTransformPointers;
152 
153  // Optimizer parameters
154  OptimizerDefinition m_OptimizerType;
155  double m_SearchRadius;
156  double m_FinalRadius;
157  unsigned int m_OptimizerMaximumIterations;
158  double m_StepSize;
159 
160  std::mutex m_LockHighestProcessedBlock;
161  int m_HighestProcessedBlock;
162 };
163 
164 } // end namespace anima
165 
166 #include "animaBaseBlockMatcher.hxx"
MetricType::Pointer MetricPointer
void SetForceComputeBlocks(bool val)
void SetBlockPositions(std::vector< PointType > &val)
void SetReferenceImage(InputImageType *image)
InputImagePointer & GetMovingImage()
void SetBlockSize(unsigned int val)
std::vector< BaseInputTransformPointer > & GetBlockTransformPointers()
ImageRegionType & GetBlockRegion(unsigned int i)
BaseInputTransformType::Pointer BaseInputTransformPointer
itk::Image< unsigned char, TInputImageType::ImageDimension > MaskImageType
virtual void TransformDependantOptimizerSetup(OptimizerPointer &optimizer)=0
void SetNumberOfWorkUnits(unsigned int val)
BaseInputTransformPointer & GetBlockTransformPointer(unsigned int i)
itk::SingleValuedCostFunction MetricType
virtual double ComputeBlockWeight(double val, unsigned int block)=0
void SetOptimizerType(OptimizerDefinition val)
MaskImagePointer & GetBlockGenerationMask()
std::vector< ImageRegionType > & GetBlockRegions()
anima::BaseTransformAgregator< TInputImageType::ImageDimension > AgregatorType
void SetBlockWeights(std::vector< double > &val)
std::vector< PointType > & GetBlockPositions()
InputImageType::PointType PointType
void SetBlockRegions(std::vector< ImageRegionType > &val)
void SetBlockSpacing(unsigned int val)
MaskImageType::Pointer MaskImagePointer
virtual AgregatorType::TRANSFORM_TYPE GetAgregatorInputTransformType()=0
OptimizerType::Pointer OptimizerPointer
AgregatorType::BaseInputTransformType BaseInputTransformType
virtual MetricPointer SetupMetric()=0
itk::SingleValuedNonLinearOptimizer OptimizerType
InputImageType::RegionType ImageRegionType
void SetBlockPercentageKept(double val)
AgregatorType::InternalScalarType InternalScalarType
void SetOptimizerMaximumIterations(unsigned int val)
const std::vector< double > & GetBlockWeights()
virtual OptimizerPointer SetupOptimizer()
void SetBlockGenerationMask(MaskImageType *image)
virtual BaseInputTransformPointer GetNewBlockTransform(PointType &blockCenter)=0
InputImageType::Pointer InputImagePointer
void SetMovingImage(InputImageType *image)
virtual void BlockMatchingSetup(MetricPointer &metric, unsigned int block)=0
OptimizerDefinition GetOptimizerType()
itk::Transform< InternalScalarType, NDimensions, NDimensions > BaseInputTransformType
virtual bool GetMaximizedMetric()=0
void SetBlockVarianceThreshold(double val)
static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION ThreadedMatching(void *arg)
InputImagePointer & GetReferenceImage()
void SetSearchRadius(double val)