ANIMA  4.0
animaPyramidalDistortionCorrectionBlockMatchingBridge.h
Go to the documentation of this file.
1 #pragma once
2 
3 #include <itkImage.h>
6 #include <itkAffineTransform.h>
8 #include <rpiDisplacementFieldTransform.h>
9 
10 enum Metric
11 {
15 };
16 
18 {
19  Direction = 0,
22 };
23 
25 {
26  Baloo = 0,
28 };
29 
30 namespace anima
31 {
32 
33 template <unsigned int ImageDimension = 3>
34 class PyramidalDistortionCorrectionBlockMatchingBridge : public itk::ProcessObject
35 {
36 public:
37  typedef itk::Image <double,ImageDimension> InputImageType;
38  typedef typename InputImageType::IOPixelType InputPixelType;
39  typedef typename InputImageType::Pointer InputImagePointer;
40  typedef typename InputImageType::ConstPointer InputImageConstPointer;
41 
45 
46  typedef itk::AffineTransform<typename BaseAgregatorType::InternalScalarType,ImageDimension> AffineTransformType;
47  typedef typename AffineTransformType::Pointer AffineTransformPointer;
48 
49  typedef rpi::DisplacementFieldTransform<typename BaseAgregatorType::ScalarType, ImageDimension> DisplacementFieldTransformType;
50  typedef typename DisplacementFieldTransformType::Pointer DisplacementFieldTransformPointer;
51  typedef typename DisplacementFieldTransformType::VectorFieldType VectorFieldType;
52 
55 
58 
59  void Update();
60 
61  void WriteOutputs();
62 
66  void SetBackwardImage(InputImageConstPointer backwardImage) {m_BackwardImage = backwardImage;}
67  void SetForwardImage(InputImageConstPointer forwardImage) {m_ForwardImage = forwardImage;}
68 
69  InputImagePointer GetOutputImage() {return m_OutputImage;}
70 
74  DisplacementFieldTransformPointer GetOutputTransform() {return m_OutputTransform;}
75 
80  std::string GetResultFile() {return m_resultFile;}
81  void SetResultFile(std::string resultFile) {m_resultFile=resultFile;}
82 
83  std::string GetOutputTransformFile() {return m_outputTransformFile;}
84  void SetOutputTransformFile(std::string outputTransformFile) {m_outputTransformFile=outputTransformFile;}
85 
86  unsigned int GetBlockSize() {return m_BlockSize;}
87  void SetBlockSize(int blockSize) {m_BlockSize=blockSize;}
88 
89  unsigned int GetBlockSpacing() {return m_BlockSpacing;}
90  void SetBlockSpacing(unsigned int blockSpacing) {m_BlockSpacing=blockSpacing;}
91 
92  double GetStDevThreshold() {return m_StDevThreshold;}
93  void SetStDevThreshold(double StDevThreshold) {m_StDevThreshold=StDevThreshold;}
94 
95  unsigned int GetTransformDirection() {return m_TransformDirection;}
96  void SetTransformDirection(unsigned int TransformDirection) {m_TransformDirection = TransformDirection;}
97 
98  unsigned int GetOptimizerMaximumIterations() {return m_OptimizerMaximumIterations;}
99  void SetOptimizerMaximumIterations(unsigned int OptimizerMaximumIterations) {m_OptimizerMaximumIterations=OptimizerMaximumIterations;}
100 
101  unsigned int GetMaximumIterations() {return m_MaximumIterations;}
102  void SetMaximumIterations(unsigned int MaximumIterations) {m_MaximumIterations=MaximumIterations;}
103 
104  double GetSearchRadius() {return m_SearchRadius;}
105  void SetSearchRadius(double SearchRadius) {m_SearchRadius=SearchRadius;}
106 
107  double GetSearchScaleRadius() {return m_SearchScaleRadius;}
108  void SetSearchScaleRadius(double SearchScaleRadius) {m_SearchScaleRadius=SearchScaleRadius;}
109 
110  double GetSearchSkewRadius() {return m_SearchSkewRadius;}
111  void SetSearchSkewRadius(double SearchSkewRadius) {m_SearchSkewRadius=SearchSkewRadius;}
112 
113  double GetFinalRadius() {return m_FinalRadius;}
114  void SetFinalRadius(double FinalRadius) {m_FinalRadius=FinalRadius;}
115 
116  double GetTranlateUpperBound() {return m_TranlateUpperBound;}
117  void SetTranlateUpperBound(double TranlateUpperBound) {m_TranlateUpperBound=TranlateUpperBound;}
118 
119  double GetScaleUpperBound() {return m_ScaleUpperBound;}
120  void SetScaleUpperBound(double ScaleUpperBound) {m_ScaleUpperBound=ScaleUpperBound;}
121 
122  double GetSkewUpperBound() {return m_SkewUpperBound;}
123  void SetSkewUpperBound(double SkewUpperBound) {m_SkewUpperBound=SkewUpperBound;}
124 
125  Metric GetMetric() {return m_Metric;}
126  void SetMetric(Metric metric) {m_Metric=metric;}
127 
128  TransformKind GetTransformKind() {return m_TransformKind;}
129  void SetTransformKind(TransformKind tr) {m_TransformKind=tr;}
130 
131  Agregator GetAgregator() {return m_Agregator;}
132  void SetAgregator(Agregator agregator) {m_Agregator=agregator;}
133 
134  double GetExtrapolationSigma() {return m_ExtrapolationSigma;}
135  void SetExtrapolationSigma(double extrapolationSigma) {m_ExtrapolationSigma = extrapolationSigma;}
136 
137  double GetElasticSigma() {return m_ElasticSigma;}
138  void SetElasticSigma(double elasticSigma) {m_ElasticSigma = elasticSigma;}
139 
140  double GetOutlierSigma() {return m_OutlierSigma;}
141  void SetOutlierSigma(double outlierSigma) {m_OutlierSigma = outlierSigma;}
142 
143  double GetMEstimateConvergenceThreshold() {return m_MEstimateConvergenceThreshold;}
144  void SetMEstimateConvergenceThreshold(double mEstimateConvergenceThreshold) {m_MEstimateConvergenceThreshold = mEstimateConvergenceThreshold;}
145 
146  double GetNeighborhoodApproximation() {return m_NeighborhoodApproximation;}
147  void SetNeighborhoodApproximation(double neighborhoodApproximation) {m_NeighborhoodApproximation = neighborhoodApproximation;}
148 
149  unsigned int GetExponentiationOrder() {return m_ExponentiationOrder;}
150  void SetExponentiationOrder(unsigned int order) {m_ExponentiationOrder = order;}
151 
152  bool GetWeightedAgregation() {return m_WeightedAgregation;}
153  void SetWeightedAgregation(bool WeightedAgregation) {m_WeightedAgregation=WeightedAgregation;}
154 
155  unsigned int GetNumberOfPyramidLevels() {return m_NumberOfPyramidLevels;}
156  void SetNumberOfPyramidLevels(unsigned int NumberOfPyramidLevels) {m_NumberOfPyramidLevels=NumberOfPyramidLevels;}
157 
158  unsigned int GetLastPyramidLevel() {return m_LastPyramidLevel;}
159  void SetLastPyramidLevel(unsigned int LastPyramidLevel) {m_LastPyramidLevel=LastPyramidLevel;}
160 
161  double GetPercentageKept() {return m_PercentageKept;}
162  void SetPercentageKept(double PercentageKept) {m_PercentageKept=PercentageKept;}
163 
164  void SetInitialTransformField(VectorFieldType *field);
165 
166 private:
167  ITK_DISALLOW_COPY_AND_ASSIGN(PyramidalDistortionCorrectionBlockMatchingBridge);
168 
169  void SetupPyramids();
170 
171  DisplacementFieldTransformPointer m_InitialTransform;
172  DisplacementFieldTransformPointer m_OutputTransform;
173  InputImagePointer m_OutputImage;
174 
175  InputImageConstPointer m_BackwardImage, m_ForwardImage;
176  PyramidPointer m_BackwardPyramid, m_ForwardPyramid;
177 
178  std::string m_outputTransformFile;
179  std::string m_resultFile;
180 
181  unsigned int m_TransformDirection;
182  unsigned int m_BlockSize;
183  unsigned int m_BlockSpacing;
184  double m_StDevThreshold;
185  unsigned int m_MaximumIterations;
186  unsigned int m_OptimizerMaximumIterations;
187  unsigned int m_HistogramSize;
188  double m_SearchRadius;
189  double m_SearchScaleRadius;
190  double m_SearchSkewRadius;
191  double m_FinalRadius;
192  double m_TranlateUpperBound;
193  double m_ScaleUpperBound;
194  double m_SkewUpperBound;
195 
196  Metric m_Metric;
197  TransformKind m_TransformKind;
198  Agregator m_Agregator;
199 
200  bool m_WeightedAgregation;
201  double m_ExtrapolationSigma;
202  double m_ElasticSigma;
203  double m_OutlierSigma;
204  double m_MEstimateConvergenceThreshold;
205  double m_NeighborhoodApproximation;
206  unsigned int m_ExponentiationOrder;
207 
208  unsigned int m_NumberOfPyramidLevels;
209  unsigned int m_LastPyramidLevel;
210  double m_PercentageKept;
211 };
212 
213 } //end namespace anima
214 
Computes a pyramid of images using the provided resampler to perform resampling.
itk::AffineTransform< typename BaseAgregatorType::InternalScalarType, ImageDimension > AffineTransformType
itk::SmartPointer< Self > Pointer
rpi::DisplacementFieldTransform< typename BaseAgregatorType::ScalarType, ImageDimension > DisplacementFieldTransformType