SimpleITK  
sitkImageRegistrationMethod.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * http://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 #ifndef sitkImageRegistrationMethod_h
19 #define sitkImageRegistrationMethod_h
20 
21 #include "sitkRegistration.h"
22 
23 #include "sitkDetail.h"
24 #include "sitkImage.h"
25 #include "sitkPixelIDTokens.h"
27 #include "sitkProcessObject.h"
28 
29 #include "sitkRandomSeed.h"
30 #include "sitkInterpolator.h"
31 #include "sitkTransform.h"
32 
33 namespace itk
34 {
35 
36 #ifndef SWIG
37 template <typename TInternalComputationValueType>
38 class ObjectToObjectOptimizerBaseTemplate;
39 template <typename TFixedImage,
40  typename TMovingImage,
41  typename TVirtualImage,
42  typename TInternalComputationValueType,
43  typename TMetricTraits>
44 class ImageToImageMetricv4;
45 
46 template <typename TFixedImageType, typename TMovingImageType, typename TVirtualImageType, typename TCoordRep>
47 class DefaultImageToImageMetricTraitsv4;
48 
49 template <typename TMetric>
50 class RegistrationParameterScalesEstimator;
51 
52 template <unsigned int VDimension>
53 class SpatialObject;
54 
55 class Command;
56 class EventObject;
57 #endif
58 
59 
60 namespace simple
61 {
62 class BSplineTransform;
63 
88 {
89 public:
92 
93  ~ImageRegistrationMethod() override;
94 
96 
97  std::string
98  GetName() const override
99  {
100  return std::string("ImageRegistrationMethod");
101  }
102 
109  std::string
110  ToString() const override;
111 
118  {
119  return this->m_Interpolator;
120  }
121  SITK_RETURN_SELF_TYPE_HEADER
123  {
124  this->m_Interpolator = Interpolator;
125  return *this;
126  }
145 #if !defined(SWIG) || defined(JAVASWIG) || defined(CSHARPSWIG)
146  // Only wrap this method if the wrapping language is strongly typed
147  SITK_RETURN_SELF_TYPE_HEADER
148  SetInitialTransform(const Transform & transform);
149 #endif
150  SITK_RETURN_SELF_TYPE_HEADER
151  SetInitialTransform(Transform & transform, bool inPlace = true);
152  Transform
154  {
155  return this->m_InitialTransform;
156  }
157  bool
159  {
160  return this->m_InitialTransformInPlace;
161  }
179  void
180  SetInitialTransformAsBSpline(BSplineTransform & transform,
181  bool inPlace = true,
182  const std::vector<unsigned int> & scaleFactors = std::vector<unsigned int>());
183 
197  SITK_RETURN_SELF_TYPE_HEADER
199  {
200  this->m_MovingInitialTransform = transform;
201  return *this;
202  }
203  Transform
205  {
206  return this->m_MovingInitialTransform;
207  }
220  SITK_RETURN_SELF_TYPE_HEADER
222  {
223  this->m_FixedInitialTransform = transform;
224  return *this;
225  }
226  Transform
228  {
229  return this->m_FixedInitialTransform;
230  }
238  SITK_RETURN_SELF_TYPE_HEADER
239  SetVirtualDomain(const std::vector<uint32_t> & virtualSize,
240  const std::vector<double> & virtualOrigin,
241  const std::vector<double> & virtualSpacing,
242  const std::vector<double> & virtualDirection);
243  SITK_RETURN_SELF_TYPE_HEADER
244  SetVirtualDomainFromImage(const Image & virtualImage);
253  SITK_RETURN_SELF_TYPE_HEADER
254  SetMetricAsANTSNeighborhoodCorrelation(unsigned int radius);
255 
260  SITK_RETURN_SELF_TYPE_HEADER
261  SetMetricAsCorrelation();
262 
267  SITK_RETURN_SELF_TYPE_HEADER
268  SetMetricAsDemons(double intensityDifferenceThreshold = 0.001);
269 
274  SITK_RETURN_SELF_TYPE_HEADER
275  SetMetricAsJointHistogramMutualInformation(unsigned int numberOfHistogramBins = 20,
276  double varianceForJointPDFSmoothing = 1.5);
277 
282  SITK_RETURN_SELF_TYPE_HEADER
283  SetMetricAsMeanSquares();
284 
290  SITK_RETURN_SELF_TYPE_HEADER
291  SetMetricAsMattesMutualInformation(unsigned int numberOfHistogramBins = 50);
292 
293 
295  {
298  EachIteration
299  };
300 
305  SITK_RETURN_SELF_TYPE_HEADER
306  SetOptimizerAsConjugateGradientLineSearch(double learningRate,
307  unsigned int numberOfIterations,
308  double convergenceMinimumValue = 1e-6,
309  unsigned int convergenceWindowSize = 10,
310  double lineSearchLowerLimit = 0,
311  double lineSearchUpperLimit = 5.0,
312  double lineSearchEpsilon = 0.01,
313  unsigned int lineSearchMaximumIterations = 20,
314  EstimateLearningRateType estimateLearningRate = Once,
315  double maximumStepSizeInPhysicalUnits = 0.0);
316 
321  SITK_RETURN_SELF_TYPE_HEADER
322  SetOptimizerAsRegularStepGradientDescent(double learningRate,
323  double minStep,
324  unsigned int numberOfIterations,
325  double relaxationFactor = 0.5,
326  double gradientMagnitudeTolerance = 1e-4,
327  EstimateLearningRateType estimateLearningRate = Never,
328  double maximumStepSizeInPhysicalUnits = 0.0);
329 
334  SITK_RETURN_SELF_TYPE_HEADER
335  SetOptimizerAsGradientDescent(double learningRate,
336  unsigned int numberOfIterations,
337  double convergenceMinimumValue = 1e-6,
338  unsigned int convergenceWindowSize = 10,
339  EstimateLearningRateType estimateLearningRate = Once,
340  double maximumStepSizeInPhysicalUnits = 0.0);
341 
346  SITK_RETURN_SELF_TYPE_HEADER
347  SetOptimizerAsGradientDescentLineSearch(double learningRate,
348  unsigned int numberOfIterations,
349  double convergenceMinimumValue = 1e-6,
350  unsigned int convergenceWindowSize = 10,
351  double lineSearchLowerLimit = 0,
352  double lineSearchUpperLimit = 5.0,
353  double lineSearchEpsilon = 0.01,
354  unsigned int lineSearchMaximumIterations = 20,
355  EstimateLearningRateType estimateLearningRate = Once,
356  double maximumStepSizeInPhysicalUnits = 0.0);
357 
364  SITK_RETURN_SELF_TYPE_HEADER
365  SetOptimizerAsLBFGSB(double gradientConvergenceTolerance = 1e-5,
366  unsigned int numberOfIterations = 500,
367  unsigned int maximumNumberOfCorrections = 5,
368  unsigned int maximumNumberOfFunctionEvaluations = 2000,
369  double costFunctionConvergenceFactor = 1e+7,
370  double lowerBound = std::numeric_limits<double>::min(),
371  double upperBound = std::numeric_limits<double>::max(),
372  bool trace = false);
373 
386  SITK_RETURN_SELF_TYPE_HEADER
387  SetOptimizerAsLBFGS2(double solutionAccuracy = 1e-5,
388  unsigned int numberOfIterations = 0,
389  unsigned int hessianApproximateAccuracy = 6,
390  unsigned int deltaConvergenceDistance = 0,
391  double deltaConvergenceTolerance = 1e-5,
392  unsigned int lineSearchMaximumEvaluations = 40,
393  double lineSearchMinimumStep = 1e-20,
394  double lineSearchMaximumStep = 1e20,
395  double lineSearchAccuracy = 1e-4);
396 
397 
414  SITK_RETURN_SELF_TYPE_HEADER
415  SetOptimizerAsExhaustive(const std::vector<unsigned int> & numberOfSteps, double stepLength = 1.0);
416 
426  SITK_RETURN_SELF_TYPE_HEADER
427  SetOptimizerAsAmoeba(double simplexDelta,
428  unsigned int numberOfIterations,
429  double parametersConvergenceTolerance = 1e-8,
430  double functionConvergenceTolerance = 1e-4,
431  bool withRestarts = false);
432 
444  SITK_RETURN_SELF_TYPE_HEADER
445  SetOptimizerWeights(const std::vector<double> & weights);
446  std::vector<double>
447  GetOptimizerWeights() const;
454  SITK_RETURN_SELF_TYPE_HEADER
455  SetOptimizerAsPowell(unsigned int numberOfIterations = 100,
456  unsigned int maximumLineIterations = 100,
457  double stepLength = 1,
458  double stepTolerance = 1e-6,
459  double valueTolerance = 1e-6);
460 
469  SITK_RETURN_SELF_TYPE_HEADER
470  SetOptimizerAsOnePlusOneEvolutionary(unsigned int numberOfIterations = 100,
471  double epsilon = 1.5e-4,
472  double initialRadius = 1.01,
473  double growthFactor = -1.0,
474  double shrinkFactor = -1.0,
475  unsigned int seed = sitkWallClock);
476 
477 
479  SITK_RETURN_SELF_TYPE_HEADER
480  SetOptimizerScales(const std::vector<double> & scales);
481 
488  SITK_RETURN_SELF_TYPE_HEADER
489  SetOptimizerScalesFromJacobian(unsigned int centralRegionRadius = 5);
490 
496  SITK_RETURN_SELF_TYPE_HEADER
497  SetOptimizerScalesFromIndexShift(unsigned int centralRegionRadius = 5, double smallParameterVariation = 0.01);
498 
505  SITK_RETURN_SELF_TYPE_HEADER
506  SetOptimizerScalesFromPhysicalShift(unsigned int centralRegionRadius = 5, double smallParameterVariation = 0.01);
507 
508 
518  SITK_RETURN_SELF_TYPE_HEADER
519  SetMetricFixedMask(const Image & binaryMask);
520 
530  SITK_RETURN_SELF_TYPE_HEADER
531  SetMetricMovingMask(const Image & binaryMask);
532 
551  SITK_RETURN_SELF_TYPE_HEADER
552  SetMetricSamplingPercentage(double percentage, unsigned int seed = sitkWallClock);
553  SITK_RETURN_SELF_TYPE_HEADER
554  SetMetricSamplingPercentagePerLevel(const std::vector<double> & percentage, unsigned int seed = sitkWallClock);
558  const std::vector<double> &
559  GetMetricSamplingPercentagePerLevel() const;
560 
572  {
575  RANDOM
576  };
577 
582  SITK_RETURN_SELF_TYPE_HEADER
583  SetMetricSamplingStrategy(MetricSamplingStrategyType strategy);
584 
594  SITK_RETURN_SELF_TYPE_HEADER
595  SetMetricUseFixedImageGradientFilter(bool);
596  SITK_RETURN_SELF_TYPE_HEADER
598  {
599  return this->SetMetricUseFixedImageGradientFilter(true);
600  }
601  SITK_RETURN_SELF_TYPE_HEADER
603  {
604  return this->SetMetricUseFixedImageGradientFilter(false);
605  }
619  SITK_RETURN_SELF_TYPE_HEADER
620  SetMetricUseMovingImageGradientFilter(bool);
621  SITK_RETURN_SELF_TYPE_HEADER
623  {
624  return this->SetMetricUseMovingImageGradientFilter(true);
625  }
626  SITK_RETURN_SELF_TYPE_HEADER
628  {
629  return this->SetMetricUseMovingImageGradientFilter(false);
630  }
641  SITK_RETURN_SELF_TYPE_HEADER
642  SetShrinkFactorsPerLevel(const std::vector<unsigned int> & shrinkFactors);
643 
652  SITK_RETURN_SELF_TYPE_HEADER
653  SetSmoothingSigmasPerLevel(const std::vector<double> & smoothingSigmas);
654 
661  SITK_RETURN_SELF_TYPE_HEADER
662  SetSmoothingSigmasAreSpecifiedInPhysicalUnits(bool arg);
663  SITK_RETURN_SELF_TYPE_HEADER
665  {
666  this->SetSmoothingSigmasAreSpecifiedInPhysicalUnits(true);
667  return *this;
668  }
669  SITK_RETURN_SELF_TYPE_HEADER
671  {
672  this->SetSmoothingSigmasAreSpecifiedInPhysicalUnits(false);
673  return *this;
674  }
679  Transform
680  Execute(const Image & fixed, const Image & moving);
681 
682 
692  double
693  MetricEvaluate(const Image & fixed, const Image & moving);
694 
695 
702  unsigned int
703  GetOptimizerIteration() const;
704  std::vector<double>
705  GetOptimizerPosition() const;
706  double
707  GetOptimizerLearningRate() const;
708  double
709  GetOptimizerConvergenceValue() const;
710  double
711  GetMetricValue() const;
712 
721  uint64_t
722  GetMetricNumberOfValidPoints() const;
723 
724 
725  unsigned int
726  GetCurrentLevel() const;
727 
735  std::vector<double>
736  GetOptimizerScales() const;
737 
740  std::string
741  GetOptimizerStopConditionDescription() const;
742 
743 
752  bool
753  StopRegistration();
754 
755 protected:
756  template <class TImage>
757  Transform
758  ExecuteInternal(const Image & fixed, const Image & moving);
759 
760  template <class TImage>
761  double
762  EvaluateInternal(const Image & fixed, const Image & moving);
763 
764 
766  CreateOptimizer(unsigned int numberOfTransformParameters);
767 
768  template <unsigned int VDimension>
770  CreateSpatialObjectMask(const Image & mask);
771 
772  template <class TImageType>
773  itk::ImageToImageMetricv4<TImageType,
774  TImageType,
775  TImageType,
776  double,
778  CreateMetric();
779 
780  template <class TImageType>
781  void
782  SetupMetric(
783  itk::ImageToImageMetricv4<TImageType,
784  TImageType,
785  TImageType,
786  double,
788  const TImageType *,
789  const TImageType *);
790 
791  template <typename TMetric>
793  CreateScalesEstimator();
794 
795  template <typename TTransformAdaptorPointer, typename TRegistrationMethod>
796  std::vector<TTransformAdaptorPointer>
797  CreateTransformParametersAdaptor(TRegistrationMethod * method);
798 
799  void
800  PreUpdate(itk::ProcessObject * p) override;
801  void
802  OnActiveProcessDelete() noexcept override;
803  unsigned long
804  AddITKObserver(const itk::EventObject &, itk::Command *) override;
805  void
806  RemoveITKObserver(EventCommand & e) override;
807 
808 private:
809  std::function<unsigned int()> m_pfGetOptimizerIteration;
810  std::function<std::vector<double>()> m_pfGetOptimizerPosition;
811  std::function<double()> m_pfGetOptimizerLearningRate;
812  std::function<double()> m_pfGetOptimizerConvergenceValue;
813  std::function<double()> m_pfGetMetricValue;
814  std::function<uint64_t()> m_pfGetMetricNumberOfValidPoints;
815  std::function<std::vector<double>()> m_pfGetOptimizerScales;
816  std::function<std::string()> m_pfGetOptimizerStopConditionDescription;
817  std::function<void()> m_pfOptimizerStopRegistration;
818 
819 
820  std::function<unsigned int()> m_pfGetCurrentLevel;
821 
822  std::function<void(itk::TransformBase * outTransform)> m_pfUpdateWithBestValue;
823 
824  template <class TMemberFunctionPointer>
826  {
827  using ObjectType = typename ::detail::FunctionTraits<TMemberFunctionPointer>::ClassType;
828 
829  template <typename TImageType>
830  TMemberFunctionPointer
831  operator()() const
832  {
833  return &ObjectType::template EvaluateInternal<TImageType>;
834  }
835  };
836 
837  typedef Transform (ImageRegistrationMethod::*MemberFunctionType)(const Image & fixed, const Image & moving);
838  typedef double (ImageRegistrationMethod::*EvaluateMemberFunctionType)(const Image & fixed, const Image & moving);
839  friend struct detail::MemberFunctionAddressor<MemberFunctionType>;
840  std::unique_ptr<detail::MemberFunctionFactory<MemberFunctionType>> m_MemberFactory;
841  std::unique_ptr<detail::MemberFunctionFactory<EvaluateMemberFunctionType>> m_EvaluateMemberFactory;
842 
848 
849  std::vector<uint32_t> m_VirtualDomainSize;
850  std::vector<double> m_VirtualDomainOrigin;
851  std::vector<double> m_VirtualDomainSpacing;
852  std::vector<double> m_VirtualDomainDirection;
853 
854  // optimizer
856  {
866  LBFGS2
867  };
889  std::vector<unsigned int> m_OptimizerNumberOfSteps;
902  unsigned int m_OptimizerSeed;
911 
912 
913  std::vector<unsigned int> m_TransformBSplineScaleFactors;
914 
915  std::vector<double> m_OptimizerWeights;
916 
918  {
922  PhysicalShift
923  };
925  std::vector<double> m_OptimizerScales;
928 
929  // metric
931  {
937  MattesMutualInformation
938  };
940  unsigned int m_MetricRadius;
944 
947 
948  std::vector<double> m_MetricSamplingPercentage;
950  unsigned int m_MetricSamplingSeed;
951 
954 
955  std::vector<unsigned int> m_ShrinkFactorsPerLevel;
956  std::vector<double> m_SmoothingSigmasPerLevel;
958 
961  unsigned int m_Iteration;
963 
965 };
966 
967 } // namespace simple
968 } // namespace itk
969 
970 #endif // sitkImageRegistrationMethod_h
itk::simple::Image
The Image class for SimpleITK.
Definition: sitkImage.h:76
itk::simple::ImageRegistrationMethod::m_InitialTransformInPlace
bool m_InitialTransformInPlace
Definition: sitkImageRegistrationMethod.h:845
itk::simple::ImageRegistrationMethod::m_Iteration
unsigned int m_Iteration
Definition: sitkImageRegistrationMethod.h:961
itk::simple::ImageRegistrationMethod::SmoothingSigmasAreSpecifiedInPhysicalUnitsOn
Self & SmoothingSigmasAreSpecifiedInPhysicalUnitsOn()
Enable the smoothing sigmas for each level in physical units (default) or in terms of voxels.
Definition: sitkImageRegistrationMethod.h:664
itk::simple::ImageRegistrationMethod::m_MetricNumberOfHistogramBins
unsigned int m_MetricNumberOfHistogramBins
Definition: sitkImageRegistrationMethod.h:942
itk::simple::ImageRegistrationMethod::EvaluateMemberFunctionAddressor
Definition: sitkImageRegistrationMethod.h:825
sitkRandomSeed.h
itk::simple::ImageRegistrationMethod::SmoothingSigmasAreSpecifiedInPhysicalUnitsOff
Self & SmoothingSigmasAreSpecifiedInPhysicalUnitsOff()
Enable the smoothing sigmas for each level in physical units (default) or in terms of voxels.
Definition: sitkImageRegistrationMethod.h:670
itk::simple::ImageRegistrationMethod::m_OptimizerMaximumNumberOfFunctionEvaluations
unsigned int m_OptimizerMaximumNumberOfFunctionEvaluations
Definition: sitkImageRegistrationMethod.h:884
itk::simple::ImageRegistrationMethod::MetricUseMovingImageGradientFilterOff
Self & MetricUseMovingImageGradientFilterOff()
Enable image gradient computation by a filter.
Definition: sitkImageRegistrationMethod.h:627
itk::simple::InterpolatorEnum
InterpolatorEnum
Definition: sitkInterpolator.h:28
itk::simple::ImageRegistrationMethod::m_MetricValue
double m_MetricValue
Definition: sitkImageRegistrationMethod.h:960
itk::simple::ImageRegistrationMethod::m_OptimizerInitialRadius
double m_OptimizerInitialRadius
Definition: sitkImageRegistrationMethod.h:899
itk::simple::ImageRegistrationMethod::m_MetricSamplingStrategy
MetricSamplingStrategyType m_MetricSamplingStrategy
Definition: sitkImageRegistrationMethod.h:949
itk::simple::ImageRegistrationMethod::OptimizerType
OptimizerType
Definition: sitkImageRegistrationMethod.h:855
itk::simple::ImageRegistrationMethod::m_OptimizerMinimumStepLength
double m_OptimizerMinimumStepLength
Definition: sitkImageRegistrationMethod.h:870
itk::simple::ImageRegistrationMethod::m_OptimizerWeights
std::vector< double > m_OptimizerWeights
Definition: sitkImageRegistrationMethod.h:915
itk::simple::ImageRegistrationMethod::MetricType
MetricType
Definition: sitkImageRegistrationMethod.h:930
itk::simple::ImageRegistrationMethod::m_OptimizerConvergenceMinimumValue
double m_OptimizerConvergenceMinimumValue
Definition: sitkImageRegistrationMethod.h:880
itk::simple::Transform
A simplified wrapper around a variety of ITK transforms.
Definition: sitkTransform.h:86
itk::simple::detail::MemberFunctionAddressor
Definition: sitkDetail.h:29
itk::simple::ImageRegistrationMethod::m_MetricUseFixedImageGradientFilter
bool m_MetricUseFixedImageGradientFilter
Definition: sitkImageRegistrationMethod.h:952
itk::simple::ImageRegistrationMethod::Jacobian
@ Jacobian
Definition: sitkImageRegistrationMethod.h:920
itk::simple::ImageRegistrationMethod::m_OptimizerLowerBound
double m_OptimizerLowerBound
Definition: sitkImageRegistrationMethod.h:886
itk::simple::ImageRegistrationMethod::Exhaustive
@ Exhaustive
Definition: sitkImageRegistrationMethod.h:862
itk::simple::ImageRegistrationMethod::m_OptimizerLineSearchMaximumEvaluations
unsigned int m_OptimizerLineSearchMaximumEvaluations
Definition: sitkImageRegistrationMethod.h:907
itk::simple::ImageRegistrationMethod::m_OptimizerNumberOfIterations
unsigned int m_OptimizerNumberOfIterations
Definition: sitkImageRegistrationMethod.h:871
itk::simple::ImageRegistrationMethod::m_OptimizerSimplexDelta
double m_OptimizerSimplexDelta
Definition: sitkImageRegistrationMethod.h:891
itk::simple::ImageRegistrationMethod::RegularStepGradientDescent
@ RegularStepGradientDescent
Definition: sitkImageRegistrationMethod.h:858
itk::simple::ImageRegistrationMethod::m_OptimizerHessianApproximationAccuracy
unsigned int m_OptimizerHessianApproximationAccuracy
Definition: sitkImageRegistrationMethod.h:904
itk::simple::ImageRegistrationMethod::MetricUseFixedImageGradientFilterOn
Self & MetricUseFixedImageGradientFilterOn()
Enable image gradient computation by a filter.
Definition: sitkImageRegistrationMethod.h:597
SITKRegistration_EXPORT
#define SITKRegistration_EXPORT
Definition: sitkRegistration.h:32
sitkImage.h
itk::simple::ImageRegistrationMethod::OptimizerScalesType
OptimizerScalesType
Definition: sitkImageRegistrationMethod.h:917
itk::simple::ImageRegistrationMethod::EvaluateMemberFunctionAddressor::ObjectType
typename ::detail::FunctionTraits< TMemberFunctionPointer >::ClassType ObjectType
Definition: sitkImageRegistrationMethod.h:827
itk::RegistrationParameterScalesEstimator
Definition: sitkImageRegistrationMethod.h:50
itk::simple::ImageRegistrationMethod::m_OptimizerLineSearchMaximumStep
double m_OptimizerLineSearchMaximumStep
Definition: sitkImageRegistrationMethod.h:909
itk::simple::ImageRegistrationMethod::m_VirtualDomainDirection
std::vector< double > m_VirtualDomainDirection
Definition: sitkImageRegistrationMethod.h:852
itk::simple::ImageRegistrationMethod::m_OptimizerDeltaConvergenceTolerance
double m_OptimizerDeltaConvergenceTolerance
Definition: sitkImageRegistrationMethod.h:906
itk::simple::ImageRegistrationMethod::m_OptimizerConvergenceWindowSize
unsigned int m_OptimizerConvergenceWindowSize
Definition: sitkImageRegistrationMethod.h:881
itk::simple::ImageRegistrationMethod::m_OptimizerLineSearchMaximumIterations
unsigned int m_OptimizerLineSearchMaximumIterations
Definition: sitkImageRegistrationMethod.h:875
itk::simple::ImageRegistrationMethod::GetInterpolator
InterpolatorEnum GetInterpolator()
Set and get the interpolator to use.
Definition: sitkImageRegistrationMethod.h:117
itk::simple::ImageRegistrationMethod::m_MetricVarianceForJointPDFSmoothing
double m_MetricVarianceForJointPDFSmoothing
Definition: sitkImageRegistrationMethod.h:943
itk::simple::ImageRegistrationMethod::m_OptimizerTrace
bool m_OptimizerTrace
Definition: sitkImageRegistrationMethod.h:888
itk::simple::ImageRegistrationMethod::m_OptimizerSolutionAccuracy
double m_OptimizerSolutionAccuracy
Definition: sitkImageRegistrationMethod.h:903
itk::simple::ImageRegistrationMethod::m_MetricIntensityDifferenceThreshold
double m_MetricIntensityDifferenceThreshold
Definition: sitkImageRegistrationMethod.h:941
itk::simple::ImageRegistrationMethod::m_OptimizerLearningRate
double m_OptimizerLearningRate
Definition: sitkImageRegistrationMethod.h:869
itk::simple::ImageRegistrationMethod::MetricUseFixedImageGradientFilterOff
Self & MetricUseFixedImageGradientFilterOff()
Enable image gradient computation by a filter.
Definition: sitkImageRegistrationMethod.h:602
itk::simple::ImageRegistrationMethod::Correlation
@ Correlation
Definition: sitkImageRegistrationMethod.h:933
itk::simple::ImageRegistrationMethod::MetricUseMovingImageGradientFilterOn
Self & MetricUseMovingImageGradientFilterOn()
Enable image gradient computation by a filter.
Definition: sitkImageRegistrationMethod.h:622
itk::simple::Command
An implementation of the Command design pattern for callback.
Definition: sitkCommand.h:44
itk::simple::ImageRegistrationMethod::m_MetricMovingMaskImage
Image m_MetricMovingMaskImage
Definition: sitkImageRegistrationMethod.h:946
itk::simple::ImageRegistrationMethod::m_OptimizerMaximumNumberOfCorrections
unsigned int m_OptimizerMaximumNumberOfCorrections
Definition: sitkImageRegistrationMethod.h:883
itk::simple::ImageRegistrationMethod::m_MetricSamplingSeed
unsigned int m_MetricSamplingSeed
Definition: sitkImageRegistrationMethod.h:950
itk::simple::ImageRegistrationMethod::ConjugateGradientLineSearch
@ ConjugateGradientLineSearch
Definition: sitkImageRegistrationMethod.h:857
itk::simple::ImageRegistrationMethod::LBFGSB
@ LBFGSB
Definition: sitkImageRegistrationMethod.h:861
itk::simple::ImageRegistrationMethod::SetMovingInitialTransform
Self & SetMovingInitialTransform(const Transform &transform)
Set a fixed transform component towards moving domain.
Definition: sitkImageRegistrationMethod.h:198
sitkMemberFunctionFactory.h
itk::simple::ImageRegistrationMethod::REGULAR
@ REGULAR
Definition: sitkImageRegistrationMethod.h:574
itk::simple::ImageRegistrationMethod::m_SmoothingSigmasPerLevel
std::vector< double > m_SmoothingSigmasPerLevel
Definition: sitkImageRegistrationMethod.h:956
itk::simple::ImageRegistrationMethod::NONE
@ NONE
Definition: sitkImageRegistrationMethod.h:573
itk::simple::ImageRegistrationMethod::m_OptimizerType
OptimizerType m_OptimizerType
Definition: sitkImageRegistrationMethod.h:868
itk::simple::ImageRegistrationMethod::m_OptimizerWithRestarts
bool m_OptimizerWithRestarts
Definition: sitkImageRegistrationMethod.h:894
sitkDetail.h
itk::simple::ImageRegistrationMethod::m_OptimizerLineSearchLowerLimit
double m_OptimizerLineSearchLowerLimit
Definition: sitkImageRegistrationMethod.h:872
itk::simple::ImageRegistrationMethod::m_OptimizerFunctionConvergenceTolerance
double m_OptimizerFunctionConvergenceTolerance
Definition: sitkImageRegistrationMethod.h:893
itk::simple::ImageRegistrationMethod::EvaluateMemberFunctionAddressor::operator()
TMemberFunctionPointer operator()() const
Definition: sitkImageRegistrationMethod.h:831
itk::simple::ImageRegistrationMethod::m_OptimizerEstimateLearningRate
EstimateLearningRateType m_OptimizerEstimateLearningRate
Definition: sitkImageRegistrationMethod.h:876
itk::simple::ImageRegistrationMethod::OnePlusOneEvolutionary
@ OnePlusOneEvolutionary
Definition: sitkImageRegistrationMethod.h:865
sitkProcessObject.h
sitkRegistration.h
itk::simple::ImageRegistrationMethod::m_OptimizerScalesType
OptimizerScalesType m_OptimizerScalesType
Definition: sitkImageRegistrationMethod.h:924
itk::simple::ImageRegistrationMethod::Demons
@ Demons
Definition: sitkImageRegistrationMethod.h:934
itk::simple::ImageRegistrationMethod::m_OptimizerSeed
unsigned int m_OptimizerSeed
Definition: sitkImageRegistrationMethod.h:902
itk::simple::ImageRegistrationMethod::GetFixedInitialTransform
Transform GetFixedInitialTransform() const
Set transform mapping to the fixed domain.
Definition: sitkImageRegistrationMethod.h:227
itk::simple::ImageRegistrationMethod::m_VirtualDomainSpacing
std::vector< double > m_VirtualDomainSpacing
Definition: sitkImageRegistrationMethod.h:851
itk::simple::ImageRegistrationMethod::GradientDescentLineSearch
@ GradientDescentLineSearch
Definition: sitkImageRegistrationMethod.h:860
itk::Command
class ITK_FORWARD_EXPORT Command
itk::simple::ImageRegistrationMethod::m_OptimizerEpsilon
double m_OptimizerEpsilon
Definition: sitkImageRegistrationMethod.h:898
itk::simple::ImageRegistrationMethod::GetInitialTransform
Transform GetInitialTransform()
Set the initial transform and parameters to optimize.
Definition: sitkImageRegistrationMethod.h:153
itk::simple::ImageRegistrationMethod::m_MetricSamplingPercentage
std::vector< double > m_MetricSamplingPercentage
Definition: sitkImageRegistrationMethod.h:948
itk::simple::ImageRegistrationMethod::m_EvaluateMemberFactory
std::unique_ptr< detail::MemberFunctionFactory< EvaluateMemberFunctionType > > m_EvaluateMemberFactory
Definition: sitkImageRegistrationMethod.h:841
itk::simple::ImageRegistrationMethod::m_ActiveOptimizer
itk::ObjectToObjectOptimizerBaseTemplate< double > * m_ActiveOptimizer
Definition: sitkImageRegistrationMethod.h:964
itk::simple::ImageRegistrationMethod::m_OptimizerUpperBound
double m_OptimizerUpperBound
Definition: sitkImageRegistrationMethod.h:887
itk::DefaultImageToImageMetricTraitsv4
Definition: sitkImageRegistrationMethod.h:47
itk::simple::ImageRegistrationMethod::m_OptimizerNumberOfSteps
std::vector< unsigned int > m_OptimizerNumberOfSteps
Definition: sitkImageRegistrationMethod.h:889
itk::simple::ImageRegistrationMethod::m_VirtualDomainOrigin
std::vector< double > m_VirtualDomainOrigin
Definition: sitkImageRegistrationMethod.h:850
itk::SpatialObject
Definition: sitkImageRegistrationMethod.h:53
itk::simple::ImageRegistrationMethod::IndexShift
@ IndexShift
Definition: sitkImageRegistrationMethod.h:921
itk::simple::ImageRegistrationMethod::ANTSNeighborhoodCorrelation
@ ANTSNeighborhoodCorrelation
Definition: sitkImageRegistrationMethod.h:932
itk::simple::ImageRegistrationMethod::m_MemberFactory
std::unique_ptr< detail::MemberFunctionFactory< MemberFunctionType > > m_MemberFactory
Definition: sitkImageRegistrationMethod.h:840
itk::simple::ImageRegistrationMethod::m_ShrinkFactorsPerLevel
std::vector< unsigned int > m_ShrinkFactorsPerLevel
Definition: sitkImageRegistrationMethod.h:955
itk::simple::ImageRegistrationMethod::m_OptimizerMaximumLineIterations
unsigned int m_OptimizerMaximumLineIterations
Definition: sitkImageRegistrationMethod.h:895
itk::simple::ImageRegistrationMethod::m_OptimizerGrowthFactor
double m_OptimizerGrowthFactor
Definition: sitkImageRegistrationMethod.h:900
itk::simple::ImageRegistrationMethod::m_InitialTransform
Transform m_InitialTransform
Definition: sitkImageRegistrationMethod.h:844
itk::simple::ImageRegistrationMethod::m_MovingInitialTransform
Transform m_MovingInitialTransform
Definition: sitkImageRegistrationMethod.h:846
sitkInterpolator.h
itk::simple::ImageRegistrationMethod::m_OptimizerGradientMagnitudeTolerance
double m_OptimizerGradientMagnitudeTolerance
Definition: sitkImageRegistrationMethod.h:879
itk::simple::ImageRegistrationMethod::m_OptimizerLineSearchEpsilon
double m_OptimizerLineSearchEpsilon
Definition: sitkImageRegistrationMethod.h:874
itk::simple::ImageRegistrationMethod::Amoeba
@ Amoeba
Definition: sitkImageRegistrationMethod.h:863
itk::simple::ImageRegistrationMethod::m_MetricRadius
unsigned int m_MetricRadius
Definition: sitkImageRegistrationMethod.h:940
itk::simple::ImageRegistrationMethod::m_OptimizerRelaxationFactor
double m_OptimizerRelaxationFactor
Definition: sitkImageRegistrationMethod.h:878
itk::simple::ImageRegistrationMethod::m_OptimizerDeltaConvergenceDistance
unsigned int m_OptimizerDeltaConvergenceDistance
Definition: sitkImageRegistrationMethod.h:905
itk::simple::ImageRegistrationMethod::m_SmoothingSigmasAreSpecifiedInPhysicalUnits
bool m_SmoothingSigmasAreSpecifiedInPhysicalUnits
Definition: sitkImageRegistrationMethod.h:957
itk::simple::ImageRegistrationMethod::m_OptimizerShrinkFactor
double m_OptimizerShrinkFactor
Definition: sitkImageRegistrationMethod.h:901
itk::simple::ImageRegistrationMethod::m_OptimizerGradientConvergenceTolerance
double m_OptimizerGradientConvergenceTolerance
Definition: sitkImageRegistrationMethod.h:882
itk::simple::ImageRegistrationMethod::m_VirtualDomainSize
std::vector< uint32_t > m_VirtualDomainSize
Definition: sitkImageRegistrationMethod.h:849
sitkTransform.h
itk::simple::ImageRegistrationMethod::GetMovingInitialTransform
Transform GetMovingInitialTransform() const
Set a fixed transform component towards moving domain.
Definition: sitkImageRegistrationMethod.h:204
itk::simple::ImageRegistrationMethod::m_OptimizerCostFunctionConvergenceFactor
double m_OptimizerCostFunctionConvergenceFactor
Definition: sitkImageRegistrationMethod.h:885
itk::simple::ImageRegistrationMethod::Never
@ Never
Never run estimation, use provided value.
Definition: sitkImageRegistrationMethod.h:296
sitkPixelIDTokens.h
itk::simple::ImageRegistrationMethod
An interface method to the modular ITKv4 registration framework.
Definition: sitkImageRegistrationMethod.h:87
itk::simple::sitkWallClock
@ sitkWallClock
Definition: sitkRandomSeed.h:30
itk::simple::ImageRegistrationMethod::Powell
@ Powell
Definition: sitkImageRegistrationMethod.h:864
itk
itk::simple::ImageRegistrationMethod::GetInitialTransformInPlace
bool GetInitialTransformInPlace() const
Set the initial transform and parameters to optimize.
Definition: sitkImageRegistrationMethod.h:158
itk::simple::ImageRegistrationMethod::m_OptimizerLineSearchMinimumStep
double m_OptimizerLineSearchMinimumStep
Definition: sitkImageRegistrationMethod.h:908
itk::simple::ImageRegistrationMethod::m_OptimizerLineSearchUpperLimit
double m_OptimizerLineSearchUpperLimit
Definition: sitkImageRegistrationMethod.h:873
itk::simple::ImageRegistrationMethod::m_OptimizerLineSearchAccuracy
double m_OptimizerLineSearchAccuracy
Definition: sitkImageRegistrationMethod.h:910
itk::ProcessObject
itk::simple::ImageRegistrationMethod::m_OptimizerMaximumStepSizeInPhysicalUnits
double m_OptimizerMaximumStepSizeInPhysicalUnits
Definition: sitkImageRegistrationMethod.h:877
itk::simple::ImageRegistrationMethod::m_OptimizerScalesCentralRegionRadius
unsigned int m_OptimizerScalesCentralRegionRadius
Definition: sitkImageRegistrationMethod.h:926
itk::ProcessObject
class ITK_FORWARD_EXPORT ProcessObject
itk::ObjectToObjectOptimizerBaseTemplate< double >
itk::simple::ProcessObject
Base class for SimpleITK classes based on ProcessObject.
Definition: sitkProcessObject.h:54
itk::simple::ImageRegistrationMethod::m_OptimizerScalesSmallParameterVariation
double m_OptimizerScalesSmallParameterVariation
Definition: sitkImageRegistrationMethod.h:927
itk::ImageToImageMetricv4
Definition: sitkImageRegistrationMethod.h:44
itk::simple::ImageRegistrationMethod::m_MetricUseMovingImageGradientFilter
bool m_MetricUseMovingImageGradientFilter
Definition: sitkImageRegistrationMethod.h:953
itk::simple::ImageRegistrationMethod::SetInterpolator
Self & SetInterpolator(InterpolatorEnum Interpolator)
Set and get the interpolator to use.
Definition: sitkImageRegistrationMethod.h:122
itk::simple::ImageRegistrationMethod::Manual
@ Manual
Definition: sitkImageRegistrationMethod.h:919
itk::simple::ImageRegistrationMethod::m_MetricFixedMaskImage
Image m_MetricFixedMaskImage
Definition: sitkImageRegistrationMethod.h:945
itk::EventObject
itk::simple::ImageRegistrationMethod::m_OptimizerValueTolerance
double m_OptimizerValueTolerance
Definition: sitkImageRegistrationMethod.h:897
itk::simple::ImageRegistrationMethod::Once
@ Once
Estimate learning once each level, ignore provided values.
Definition: sitkImageRegistrationMethod.h:297
itk::simple::ImageRegistrationMethod::EstimateLearningRateType
EstimateLearningRateType
Definition: sitkImageRegistrationMethod.h:294
itk::simple::ImageRegistrationMethod::m_FixedInitialTransform
Transform m_FixedInitialTransform
Definition: sitkImageRegistrationMethod.h:847
itk::simple::BSplineTransform
A deformable transform over a bounded spatial domain using a BSpline representation for a 2D or 3D co...
Definition: sitkBSplineTransform.h:33
itk::simple::ImageRegistrationMethod::GradientDescent
@ GradientDescent
Definition: sitkImageRegistrationMethod.h:859
itk::simple::ImageRegistrationMethod::m_OptimizerScales
std::vector< double > m_OptimizerScales
Definition: sitkImageRegistrationMethod.h:925
itk::simple::ImageRegistrationMethod::GetName
std::string GetName() const override
Definition: sitkImageRegistrationMethod.h:98
itk::simple::ImageRegistrationMethod::MetricSamplingStrategyType
MetricSamplingStrategyType
Sampling strategies for obtaining points.
Definition: sitkImageRegistrationMethod.h:571
itk::simple::ImageRegistrationMethod::m_TransformBSplineScaleFactors
std::vector< unsigned int > m_TransformBSplineScaleFactors
Definition: sitkImageRegistrationMethod.h:913
itk::simple::ImageRegistrationMethod::m_OptimizerStepLength
double m_OptimizerStepLength
Definition: sitkImageRegistrationMethod.h:890
itk::simple::ImageRegistrationMethod::m_StopConditionDescription
std::string m_StopConditionDescription
Definition: sitkImageRegistrationMethod.h:959
itk::simple::ImageRegistrationMethod::m_Interpolator
InterpolatorEnum m_Interpolator
Definition: sitkImageRegistrationMethod.h:843
itk::TransformBaseTemplate
Definition: sitkTransform.h:32
itk::simple::ImageRegistrationMethod::m_MetricType
MetricType m_MetricType
Definition: sitkImageRegistrationMethod.h:939
itk::simple::ImageRegistrationMethod::m_NumberOfValidPoints
uint64_t m_NumberOfValidPoints
Definition: sitkImageRegistrationMethod.h:962
itk::simple::ImageRegistrationMethod::MeanSquares
@ MeanSquares
Definition: sitkImageRegistrationMethod.h:936
itk::simple::ImageRegistrationMethod::SetFixedInitialTransform
Self & SetFixedInitialTransform(const Transform &transform)
Set transform mapping to the fixed domain.
Definition: sitkImageRegistrationMethod.h:221
itk::simple::ImageRegistrationMethod::m_OptimizerParametersConvergenceTolerance
double m_OptimizerParametersConvergenceTolerance
Definition: sitkImageRegistrationMethod.h:892
itk::simple::ImageRegistrationMethod::m_OptimizerStepTolerance
double m_OptimizerStepTolerance
Definition: sitkImageRegistrationMethod.h:896
itk::simple::ImageRegistrationMethod::JointHistogramMutualInformation
@ JointHistogramMutualInformation
Definition: sitkImageRegistrationMethod.h:935