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
33namespace itk
34{
35
36#ifndef SWIG
37template <typename TInternalComputationValueType>
39template <typename TFixedImage,
40 typename TMovingImage,
41 typename TVirtualImage,
42 typename TInternalComputationValueType,
43 typename TMetricTraits>
45
46template <typename TFixedImageType, typename TMovingImageType, typename TVirtualImageType, typename TCoordRep>
48
49template <typename TMetric>
51
52template <unsigned int VDimension>
53class SpatialObject;
54
55class Command;
56class EventObject;
57#endif
58
59
60namespace simple
61{
63
88{
89public:
92
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 }
127
128
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);
154 {
155 return this->m_InitialTransform;
156 }
157 bool
159 {
160 return this->m_InitialTransformInPlace;
161 }
162
163
179 void
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 }
205 {
206 return this->m_MovingInitialTransform;
207 }
208
209
220 SITK_RETURN_SELF_TYPE_HEADER
222 {
223 this->m_FixedInitialTransform = transform;
224 return *this;
225 }
228 {
229 return this->m_FixedInitialTransform;
230 }
231
232
233
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);
246
253 SITK_RETURN_SELF_TYPE_HEADER
255
260 SITK_RETURN_SELF_TYPE_HEADER
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
284
290 SITK_RETURN_SELF_TYPE_HEADER
291 SetMetricAsMattesMutualInformation(unsigned int numberOfHistogramBins = 50);
292
293
300
305 SITK_RETURN_SELF_TYPE_HEADER
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
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
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>
449
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);
556
558 const std::vector<double> &
560
577
582 SITK_RETURN_SELF_TYPE_HEADER
584
594 SITK_RETURN_SELF_TYPE_HEADER
596 SITK_RETURN_SELF_TYPE_HEADER
598 SITK_RETURN_SELF_TYPE_HEADER
601
602
613 SITK_RETURN_SELF_TYPE_HEADER
615 SITK_RETURN_SELF_TYPE_HEADER
617 SITK_RETURN_SELF_TYPE_HEADER
620
621
629 SITK_RETURN_SELF_TYPE_HEADER
630 SetShrinkFactorsPerLevel(const std::vector<unsigned int> & shrinkFactors);
631
640 SITK_RETURN_SELF_TYPE_HEADER
641 SetSmoothingSigmasPerLevel(const std::vector<double> & smoothingSigmas);
642
649 SITK_RETURN_SELF_TYPE_HEADER
651 SITK_RETURN_SELF_TYPE_HEADER
657 SITK_RETURN_SELF_TYPE_HEADER
663
664
665
668 Execute(const Image & fixed, const Image & moving);
669
670
680 double
681 MetricEvaluate(const Image & fixed, const Image & moving);
682
683
690 unsigned int
692 std::vector<double>
694 double
696 double
698 double
700
709 uint64_t
711
712
713 unsigned int
715
723 std::vector<double>
725
728 std::string
730
731
740 bool
742
743protected:
744 template <class TImage>
746 ExecuteInternal(const Image & fixed, const Image & moving);
747
748 template <class TImage>
749 double
750 EvaluateInternal(const Image & fixed, const Image & moving);
751
752
754 CreateOptimizer(unsigned int numberOfTransformParameters);
755
756 template <unsigned int VDimension>
759
760 template <class TImageType>
761 itk::ImageToImageMetricv4<TImageType,
762 TImageType,
763 TImageType,
764 double,
767
768 template <class TImageType>
769 void
771 itk::ImageToImageMetricv4<TImageType,
772 TImageType,
773 TImageType,
774 double,
776 const TImageType *,
777 const TImageType *);
778
779 template <typename TMetric>
782
783 template <typename TTransformAdaptorPointer, typename TRegistrationMethod>
784 std::vector<TTransformAdaptorPointer>
785 CreateTransformParametersAdaptor(TRegistrationMethod * method);
786
787 void
789 void
790 OnActiveProcessDelete() noexcept override;
791 unsigned long
792 AddITKObserver(const itk::EventObject &, itk::Command *) override;
793 void
795
796private:
797 std::function<unsigned int()> m_pfGetOptimizerIteration;
798 std::function<std::vector<double>()> m_pfGetOptimizerPosition;
801 std::function<double()> m_pfGetMetricValue;
803 std::function<std::vector<double>()> m_pfGetOptimizerScales;
806
807
808 std::function<unsigned int()> m_pfGetCurrentLevel;
809
810 std::function<void(itk::TransformBase * outTransform)> m_pfUpdateWithBestValue;
811
812 template <class TMemberFunctionPointer>
814 {
815 using ObjectType = typename ::detail::FunctionTraits<TMemberFunctionPointer>::ClassType;
816
817 template <typename TImageType>
818 TMemberFunctionPointer
820 {
821 return &ObjectType::template EvaluateInternal<TImageType>;
822 }
823 };
824
825 typedef Transform (ImageRegistrationMethod::*MemberFunctionType)(const Image & fixed, const Image & moving);
826 typedef double (ImageRegistrationMethod::*EvaluateMemberFunctionType)(const Image & fixed, const Image & moving);
828 std::unique_ptr<detail::MemberFunctionFactory<MemberFunctionType>> m_MemberFactory;
829 std::unique_ptr<detail::MemberFunctionFactory<EvaluateMemberFunctionType>> m_EvaluateMemberFactory;
830
836
837 std::vector<uint32_t> m_VirtualDomainSize;
838 std::vector<double> m_VirtualDomainOrigin;
839 std::vector<double> m_VirtualDomainSpacing;
840 std::vector<double> m_VirtualDomainDirection;
841
842 // optimizer
877 std::vector<unsigned int> m_OptimizerNumberOfSteps;
890 unsigned int m_OptimizerSeed;
899
900
901 std::vector<unsigned int> m_TransformBSplineScaleFactors;
902
903 std::vector<double> m_OptimizerWeights;
904
913 std::vector<double> m_OptimizerScales;
916
917 // metric
928 unsigned int m_MetricRadius;
932
935
936 std::vector<double> m_MetricSamplingPercentage;
939
942
943 std::vector<unsigned int> m_ShrinkFactorsPerLevel;
944 std::vector<double> m_SmoothingSigmasPerLevel;
946
949 unsigned int m_Iteration;
951
953};
954
955} // namespace simple
956} // namespace itk
957
958#endif // sitkImageRegistrationMethod_h
A deformable transform over a bounded spatial domain using a BSpline representation for a 2D or 3D co...
Self & SetInitialTransform(Transform &transform, bool inPlace=true)
Set the initial transform and parameters to optimize.
double EvaluateInternal(const Image &fixed, const Image &moving)
Self & SetMetricSamplingStrategy(MetricSamplingStrategyType strategy)
Set sampling strategy for sample generation.
std::function< std::vector< double >()> m_pfGetOptimizerPosition
Self & SetInterpolator(InterpolatorEnum Interpolator)
Set and get the interpolator to use.
Self & SetVirtualDomain(const std::vector< uint32_t > &virtualSize, const std::vector< double > &virtualOrigin, const std::vector< double > &virtualSpacing, const std::vector< double > &virtualDirection)
Set the virtual domain used for sampling.
Self & SetSmoothingSigmasAreSpecifiedInPhysicalUnits(bool arg)
Enable the smoothing sigmas for each level in physical units (default) or in terms of voxels.
Transform(ImageRegistrationMethod::* MemberFunctionType)(const Image &fixed, const Image &moving)
Self & SetMetricUseFixedImageGradientFilter(bool)
Enable image gradient computation by a filter.
double(ImageRegistrationMethod::* EvaluateMemberFunctionType)(const Image &fixed, const Image &moving)
std::vector< double > GetOptimizerScales() const
Get the OptimizerScales.
Self & SetMetricMovingMask(const Image &binaryMask)
Set an image mask in order to restrict the sampled points for the metric in the moving image space.
uint64_t GetMetricNumberOfValidPoints() const
Self & SetMetricAsDemons(double intensityDifferenceThreshold=0.001)
Use demons image metric.
itk::SpatialObject< VDimension > * CreateSpatialObjectMask(const Image &mask)
const std::vector< double > & GetMetricSamplingPercentagePerLevel() const
Get the percentage of pixels used for metric evaluation.
Self & MetricUseFixedImageGradientFilterOn()
Enable image gradient computation by a filter.
itk::ImageToImageMetricv4< TImageType, TImageType, TImageType, double, itk::DefaultImageToImageMetricTraitsv4< TImageType, TImageType, TImageType, double > > * CreateMetric()
Self & SetOptimizerAsAmoeba(double simplexDelta, unsigned int numberOfIterations, double parametersConvergenceTolerance=1e-8, double functionConvergenceTolerance=1e-4, bool withRestarts=false)
Set optimizer to Nelder-Mead downhill simplex algorithm.
bool GetInitialTransformInPlace() const
Set the initial transform and parameters to optimize.
Transform Execute(const Image &fixed, const Image &moving)
Optimize the configured registration problem.
std::vector< double > GetOptimizerPosition() const
Self & SetOptimizerWeights(const std::vector< double > &weights)
A per parameter weighting array for the optimizer.
Self & SetOptimizerAsExhaustive(const std::vector< unsigned int > &numberOfSteps, double stepLength=1.0)
Set the optimizer to sample the metric at regular steps.
Self & SetMetricAsMattesMutualInformation(unsigned int numberOfHistogramBins=50)
Use the mutual information between two images to be registered using the method of Mattes et al.
Self & SetOptimizerAsOnePlusOneEvolutionary(unsigned int numberOfIterations=100, double epsilon=1.5e-4, double initialRadius=1.01, double growthFactor=-1.0, double shrinkFactor=-1.0, unsigned int seed=sitkWallClock)
1+1 evolutionary optimizer strategy.
void PreUpdate(itk::ProcessObject *p) override
Self & SetMetricSamplingPercentagePerLevel(const std::vector< double > &percentage, unsigned int seed=sitkWallClock)
Set percentage of pixels sampled for metric evaluation.
Self & SetFixedInitialTransform(const Transform &transform)
Set transform mapping to the fixed domain.
Self & SetOptimizerAsRegularStepGradientDescent(double learningRate, double minStep, unsigned int numberOfIterations, double relaxationFactor=0.5, double gradientMagnitudeTolerance=1e-4, EstimateLearningRateType estimateLearningRate=Never, double maximumStepSizeInPhysicalUnits=0.0)
Regular Step Gradient descent optimizer.
Transform GetFixedInitialTransform() const
Set transform mapping to the fixed domain.
itk::ObjectToObjectOptimizerBaseTemplate< double > * m_ActiveOptimizer
Self & SetShrinkFactorsPerLevel(const std::vector< unsigned int > &shrinkFactors)
Set the isotropic shrink factors for each level.
Self & SmoothingSigmasAreSpecifiedInPhysicalUnitsOff()
Enable the smoothing sigmas for each level in physical units (default) or in terms of voxels.
Self & SetMetricAsCorrelation()
Use negative normalized cross correlation image metric.
std::string ToString() const override
Print the information about the object to a string.
Self & SetMetricUseMovingImageGradientFilter(bool)
Enable image gradient computation by a filter.
Self & SetSmoothingSigmasPerLevel(const std::vector< double > &smoothingSigmas)
Set the sigmas of Gaussian used for smoothing.
Transform GetInitialTransform()
Set the initial transform and parameters to optimize.
Self & SetOptimizerAsPowell(unsigned int numberOfIterations=100, unsigned int maximumLineIterations=100, double stepLength=1, double stepTolerance=1e-6, double valueTolerance=1e-6)
Powell optimization using Brent line search.
Self & MetricUseMovingImageGradientFilterOff()
Enable image gradient computation by a filter.
std::string GetOptimizerStopConditionDescription() const
Self & SetMetricAsANTSNeighborhoodCorrelation(unsigned int radius)
Use normalized cross correlation using a small neighborhood for each voxel between two images,...
itk::RegistrationParameterScalesEstimator< TMetric > * CreateScalesEstimator()
std::function< unsigned int()> m_pfGetOptimizerIteration
MetricSamplingStrategyType
Sampling strategies for obtaining points.
Self & SetMetricSamplingPercentage(double percentage, unsigned int seed=sitkWallClock)
Set percentage of pixels sampled for metric evaluation.
std::function< void(itk::TransformBase *outTransform)> m_pfUpdateWithBestValue
itk::ObjectToObjectOptimizerBaseTemplate< double > * CreateOptimizer(unsigned int numberOfTransformParameters)
Self & SetVirtualDomainFromImage(const Image &virtualImage)
Set the virtual domain used for sampling.
void OnActiveProcessDelete() noexcept override
Transform ExecuteInternal(const Image &fixed, const Image &moving)
unsigned long AddITKObserver(const itk::EventObject &, itk::Command *) override
double MetricEvaluate(const Image &fixed, const Image &moving)
Get the value of the metric given the state of the method.
Self & SetMovingInitialTransform(const Transform &transform)
Set a fixed transform component towards moving domain.
std::vector< unsigned int > m_TransformBSplineScaleFactors
Self & SmoothingSigmasAreSpecifiedInPhysicalUnitsOn()
Enable the smoothing sigmas for each level in physical units (default) or in terms of voxels.
std::vector< double > GetOptimizerWeights() const
A per parameter weighting array for the optimizer.
void SetupMetric(itk::ImageToImageMetricv4< TImageType, TImageType, TImageType, double, itk::DefaultImageToImageMetricTraitsv4< TImageType, TImageType, TImageType, double > > *, const TImageType *, const TImageType *)
Self & SetOptimizerAsLBFGSB(double gradientConvergenceTolerance=1e-5, unsigned int numberOfIterations=500, unsigned int maximumNumberOfCorrections=5, unsigned int maximumNumberOfFunctionEvaluations=2000, double costFunctionConvergenceFactor=1e+7, double lowerBound=std::numeric_limits< double >::min(), double upperBound=std::numeric_limits< double >::max(), bool trace=false)
Limited memory Broyden Fletcher Goldfarb Shannon minimization with simple bounds.
void RemoveITKObserver(EventCommand &e) override
Self & SetMetricAsMeanSquares()
Use negative means squares image metric.
InterpolatorEnum GetInterpolator()
Set and get the interpolator to use.
Self & SetOptimizerScalesFromPhysicalShift(unsigned int centralRegionRadius=5, double smallParameterVariation=0.01)
Estimating scales of transform parameters a step sizes, from the maximum voxel shift in physical spac...
std::function< uint64_t()> m_pfGetMetricNumberOfValidPoints
Self & MetricUseMovingImageGradientFilterOn()
Enable image gradient computation by a filter.
Self & SetMetricAsJointHistogramMutualInformation(unsigned int numberOfHistogramBins=20, double varianceForJointPDFSmoothing=1.5)
Use mutual information between two images.
unsigned int GetCurrentLevel() const
Self & SetMetricFixedMask(const Image &binaryMask)
Set an image mask in order to restrict the sampled points for the metric.
@ Once
Estimate learning once each level, ignore provided values.
@ EachIteration
Estimate learning rate at each iteration.
@ Never
Never run estimation, use provided value.
std::function< unsigned int()> m_pfGetCurrentLevel
Self & SetOptimizerAsConjugateGradientLineSearch(double learningRate, unsigned int numberOfIterations, double convergenceMinimumValue=1e-6, unsigned int convergenceWindowSize=10, double lineSearchLowerLimit=0, double lineSearchUpperLimit=5.0, double lineSearchEpsilon=0.01, unsigned int lineSearchMaximumIterations=20, EstimateLearningRateType estimateLearningRate=Once, double maximumStepSizeInPhysicalUnits=0.0)
Conjugate gradient descent optimizer with a golden section line search for nonlinear optimization.
unsigned int GetOptimizerIteration() const
std::function< std::string()> m_pfGetOptimizerStopConditionDescription
Self & SetOptimizerAsGradientDescent(double learningRate, unsigned int numberOfIterations, double convergenceMinimumValue=1e-6, unsigned int convergenceWindowSize=10, EstimateLearningRateType estimateLearningRate=Once, double maximumStepSizeInPhysicalUnits=0.0)
Gradient descent optimizer.
std::unique_ptr< detail::MemberFunctionFactory< EvaluateMemberFunctionType > > m_EvaluateMemberFactory
void SetInitialTransformAsBSpline(BSplineTransform &transform, bool inPlace=true, const std::vector< unsigned int > &scaleFactors=std::vector< unsigned int >())
Set an initial BSpline transform to optimize.
Self & SetOptimizerScalesFromJacobian(unsigned int centralRegionRadius=5)
Estimate scales from Jacobian norms.
Self & SetOptimizerScalesFromIndexShift(unsigned int centralRegionRadius=5, double smallParameterVariation=0.01)
Estimate scales from maximum voxel shift in index space cause by parameter change.
std::vector< TTransformAdaptorPointer > CreateTransformParametersAdaptor(TRegistrationMethod *method)
Self & SetOptimizerScales(const std::vector< double > &scales)
Manually set per parameter weighting for the transform parameters.
std::function< std::vector< double >()> m_pfGetOptimizerScales
std::unique_ptr< detail::MemberFunctionFactory< MemberFunctionType > > m_MemberFactory
Transform GetMovingInitialTransform() const
Set a fixed transform component towards moving domain.
Self & SetOptimizerAsLBFGS2(double solutionAccuracy=1e-5, unsigned int numberOfIterations=0, unsigned int hessianApproximateAccuracy=6, unsigned int deltaConvergenceDistance=0, double deltaConvergenceTolerance=1e-5, unsigned int lineSearchMaximumEvaluations=40, double lineSearchMinimumStep=1e-20, double lineSearchMaximumStep=1e20, double lineSearchAccuracy=1e-4)
Limited memory Broyden Fletcher Goldfarb Shannon minimization without bounds.
Self & SetInitialTransform(const Transform &transform)
Set the initial transform and parameters to optimize.
Self & MetricUseFixedImageGradientFilterOff()
Enable image gradient computation by a filter.
Self & SetOptimizerAsGradientDescentLineSearch(double learningRate, unsigned int numberOfIterations, double convergenceMinimumValue=1e-6, unsigned int convergenceWindowSize=10, double lineSearchLowerLimit=0, double lineSearchUpperLimit=5.0, double lineSearchEpsilon=0.01, unsigned int lineSearchMaximumIterations=20, EstimateLearningRateType estimateLearningRate=Once, double maximumStepSizeInPhysicalUnits=0.0)
Gradient descent optimizer with a golden section line search.
The Image class for SimpleITK.
Definition sitkImage.h:77
friend class itk::simple::Command
A simplified wrapper around a variety of ITK transforms.
TransformBaseTemplate< double > TransformBase
STL namespace.
#define SITKRegistration_EXPORT
typename ::detail::FunctionTraits< TMemberFunctionPointer >::ClassType ObjectType