SimpleITK  1.0.1
sitkImageRegistrationMethod.h
Go to the documentation of this file.
1 /*=========================================================================
2 *
3 * Copyright Insight Software Consortium
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 
34 namespace itk
35 {
36 
37 #ifndef SWIG
38 template< typename TInternalComputationValueType> 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,
47  typename TMovingImageType,
48  typename TVirtualImageType,
49  typename TCoordRep>
50 class DefaultImageToImageMetricTraitsv4;
51 
52 template<typename TMetric> class RegistrationParameterScalesEstimator;
53 
54 template<unsigned int VDimension> class SpatialObject;
55 
56 class Command;
57 class EventObject;
58 #endif
59 
60 
61 
62 namespace simple
63 {
64 
89  : public ProcessObject
90  {
91  public:
92 
95 
97  virtual ~ImageRegistrationMethod();
98 
99  std::string GetName() const { return std::string("ImageRegistrationMethod"); }
100 
107  std::string ToString() const;
108 
114  { return this->m_Interpolator; }
115  SITK_RETURN_SELF_TYPE_HEADER SetInterpolator ( InterpolatorEnum Interpolator )
116  { this->m_Interpolator = Interpolator; return *this; }
135 #if !defined(SWIG) || defined(JAVASWIG) || defined(CSHARPSWIG)
136  // Only wrap this method if the wrapping language is strongly typed
137  SITK_RETURN_SELF_TYPE_HEADER SetInitialTransform ( const Transform &transform );
138 #endif
139  SITK_RETURN_SELF_TYPE_HEADER SetInitialTransform ( Transform &transform, bool inPlace=true );
141  { return this->m_InitialTransform; }
143  { return this->m_InitialTransformInPlace;}
159  SITK_RETURN_SELF_TYPE_HEADER SetMovingInitialTransform( const Transform &transform )
160  { this->m_MovingInitialTransform = transform; return *this; }
162  { return this->m_MovingInitialTransform; }
175  SITK_RETURN_SELF_TYPE_HEADER SetFixedInitialTransform( const Transform &transform )
176  { this->m_FixedInitialTransform = transform; return *this; }
178  { return this->m_FixedInitialTransform; }
186  SITK_RETURN_SELF_TYPE_HEADER SetVirtualDomain( const std::vector<uint32_t> &virtualSize,
187  const std::vector<double> &virtualOrigin,
188  const std::vector<double> &virtualSpacing,
189  const std::vector<double> &virtualDirection );
190  SITK_RETURN_SELF_TYPE_HEADER SetVirtualDomainFromImage( const Image &virtualImage );
199  SITK_RETURN_SELF_TYPE_HEADER SetMetricAsANTSNeighborhoodCorrelation( unsigned int radius );
200 
205  SITK_RETURN_SELF_TYPE_HEADER SetMetricAsCorrelation( );
206 
211  SITK_RETURN_SELF_TYPE_HEADER SetMetricAsDemons( double intensityDifferenceThreshold = 0.001 );
212 
217  SITK_RETURN_SELF_TYPE_HEADER SetMetricAsJointHistogramMutualInformation( unsigned int numberOfHistogramBins = 20,
218  double varianceForJointPDFSmoothing = 1.5);
219 
224  SITK_RETURN_SELF_TYPE_HEADER SetMetricAsMeanSquares( );
225 
231  SITK_RETURN_SELF_TYPE_HEADER SetMetricAsMattesMutualInformation( unsigned int numberOfHistogramBins = 50 );
232 
233 
236  EachIteration
237  };
238 
243  SITK_RETURN_SELF_TYPE_HEADER SetOptimizerAsConjugateGradientLineSearch( double learningRate,
244  unsigned int numberOfIterations,
245  double convergenceMinimumValue = 1e-6,
246  unsigned int convergenceWindowSize = 10,
247  double lineSearchLowerLimit = 0,
248  double lineSearchUpperLimit = 5.0,
249  double lineSearchEpsilon = 0.01,
250  unsigned int lineSearchMaximumIterations = 20,
251  EstimateLearningRateType estimateLearningRate = Once,
252  double maximumStepSizeInPhysicalUnits = 0.0 );
253 
258  SITK_RETURN_SELF_TYPE_HEADER SetOptimizerAsRegularStepGradientDescent( double learningRate,
259  double minStep,
260  unsigned int numberOfIterations,
261  double relaxationFactor=0.5,
262  double gradientMagnitudeTolerance = 1e-4,
263  EstimateLearningRateType estimateLearningRate = Never,
264  double maximumStepSizeInPhysicalUnits = 0.0);
265 
270  SITK_RETURN_SELF_TYPE_HEADER SetOptimizerAsGradientDescent( double learningRate,
271  unsigned int numberOfIterations,
272  double convergenceMinimumValue = 1e-6,
273  unsigned int convergenceWindowSize = 10,
274  EstimateLearningRateType estimateLearningRate = Once,
275  double maximumStepSizeInPhysicalUnits = 0.0);
276 
281  SITK_RETURN_SELF_TYPE_HEADER SetOptimizerAsGradientDescentLineSearch( double learningRate,
282  unsigned int numberOfIterations,
283  double convergenceMinimumValue = 1e-6,
284  unsigned int convergenceWindowSize = 10,
285  double lineSearchLowerLimit = 0,
286  double lineSearchUpperLimit = 5.0,
287  double lineSearchEpsilon = 0.01,
288  unsigned int lineSearchMaximumIterations = 20,
289  EstimateLearningRateType estimateLearningRate = Once,
290  double maximumStepSizeInPhysicalUnits = 0.0 );
291 
298  SITK_RETURN_SELF_TYPE_HEADER SetOptimizerAsLBFGSB(double gradientConvergenceTolerance = 1e-5,
299  unsigned int numberOfIterations = 500,
300  unsigned int maximumNumberOfCorrections = 5,
301  unsigned int maximumNumberOfFunctionEvaluations = 2000,
302  double costFunctionConvergenceFactor = 1e+7,
303  double lowerBound = std::numeric_limits<double>::min(),
304  double upperBound = std::numeric_limits<double>::max(),
305  bool trace = false );
306 
323  SITK_RETURN_SELF_TYPE_HEADER SetOptimizerAsExhaustive( const std::vector<unsigned int> &numberOfSteps,
324  double stepLength = 1.0 );
325 
330  SITK_RETURN_SELF_TYPE_HEADER SetOptimizerAsAmoeba(double simplexDelta,
331  unsigned int numberOfIterations,
332  double parametersConvergenceTolerance=1e-8,
333  double functionConvergenceTolerance=1e-4,
334  bool withRestarts = false );
335 
347  SITK_RETURN_SELF_TYPE_HEADER SetOptimizerWeights(const std::vector<double> &weights);
348  std::vector<double> GetOptimizerWeights( ) const;
355  SITK_RETURN_SELF_TYPE_HEADER SetOptimizerAsPowell(unsigned int numberOfIterations = 100,
356  unsigned int maximumLineIterations = 100,
357  double stepLength = 1,
358  double stepTolerance = 1e-6,
359  double valueTolerance = 1e-6 );
360 
369  SITK_RETURN_SELF_TYPE_HEADER SetOptimizerAsOnePlusOneEvolutionary(unsigned int numberOfIterations = 100,
370  double epsilon = 1.5e-4,
371  double initialRadius = 1.01,
372  double growthFactor = -1.0,
373  double shrinkFactor = -1.0,
374  unsigned int seed = sitkWallClock);
375 
376 
377 
379  SITK_RETURN_SELF_TYPE_HEADER SetOptimizerScales( const std::vector<double> &scales );
380 
387  SITK_RETURN_SELF_TYPE_HEADER SetOptimizerScalesFromJacobian( unsigned int centralRegionRadius = 5 );
388 
394  SITK_RETURN_SELF_TYPE_HEADER SetOptimizerScalesFromIndexShift( unsigned int centralRegionRadius = 5,
395  double smallParameterVariation = 0.01 );
396 
403  SITK_RETURN_SELF_TYPE_HEADER SetOptimizerScalesFromPhysicalShift( unsigned int centralRegionRadius = 5,
404  double smallParameterVariation = 0.01 );
405 
406 
416  SITK_RETURN_SELF_TYPE_HEADER SetMetricFixedMask( const Image &binaryMask );
417 
427  SITK_RETURN_SELF_TYPE_HEADER SetMetricMovingMask( const Image &binaryMask );
428 
439  SITK_RETURN_SELF_TYPE_HEADER SetMetricSamplingPercentage(double percentage, unsigned int seed = sitkWallClock);
440  SITK_RETURN_SELF_TYPE_HEADER SetMetricSamplingPercentagePerLevel(const std::vector<double> &percentage, unsigned int seed = sitkWallClock);
446  RANDOM
447  };
448 
453  SITK_RETURN_SELF_TYPE_HEADER SetMetricSamplingStrategy( MetricSamplingStrategyType strategy);
454 
464  SITK_RETURN_SELF_TYPE_HEADER SetMetricUseFixedImageGradientFilter(bool);
465  SITK_RETURN_SELF_TYPE_HEADER MetricUseFixedImageGradientFilterOn() {return this->SetMetricUseFixedImageGradientFilter(true);}
466  SITK_RETURN_SELF_TYPE_HEADER MetricUseFixedImageGradientFilterOff() {return this->SetMetricUseFixedImageGradientFilter(false);}
479  SITK_RETURN_SELF_TYPE_HEADER SetMetricUseMovingImageGradientFilter(bool);
480  SITK_RETURN_SELF_TYPE_HEADER MetricUseMovingImageGradientFilterOn() {return this->SetMetricUseMovingImageGradientFilter(true);}
481  SITK_RETURN_SELF_TYPE_HEADER MetricUseMovingImageGradientFilterOff() {return this->SetMetricUseMovingImageGradientFilter(false);}
490  SITK_RETURN_SELF_TYPE_HEADER SetShrinkFactorsPerLevel( const std::vector<unsigned int> &shrinkFactors );
491 
497  SITK_RETURN_SELF_TYPE_HEADER SetSmoothingSigmasPerLevel( const std::vector<double> &smoothingSigmas );
498 
505  SITK_RETURN_SELF_TYPE_HEADER SetSmoothingSigmasAreSpecifiedInPhysicalUnits(bool arg);
506  SITK_RETURN_SELF_TYPE_HEADER SmoothingSigmasAreSpecifiedInPhysicalUnitsOn() { this->SetSmoothingSigmasAreSpecifiedInPhysicalUnits(true); return *this;}
507  SITK_RETURN_SELF_TYPE_HEADER SmoothingSigmasAreSpecifiedInPhysicalUnitsOff() { this->SetSmoothingSigmasAreSpecifiedInPhysicalUnits(false); return *this;}
512  Transform Execute ( const Image &fixed, const Image & moving );
513 
514 
524  double MetricEvaluate( const Image &fixed, const Image & moving );
525 
526 
533  unsigned int GetOptimizerIteration() const;
534  std::vector<double> GetOptimizerPosition() const;
535  double GetOptimizerLearningRate() const;
536  double GetOptimizerConvergenceValue() const;
537  double GetMetricValue() const;
538 
539 
540  unsigned int GetCurrentLevel() const;
541 
549  std::vector<double> GetOptimizerScales() const;
550 
553  std::string GetOptimizerStopConditionDescription() const;
554 
555  protected:
556 
557  template<class TImage>
558  Transform ExecuteInternal ( const Image &fixed, const Image &moving );
559 
560  template<class TImage>
561  double EvaluateInternal ( const Image &fixed, const Image &moving );
562 
563 
564  itk::ObjectToObjectOptimizerBaseTemplate<double> *CreateOptimizer( unsigned int numberOfTransformParameters );
565 
566  template<unsigned int VDimension>
567  itk::SpatialObject<VDimension> *CreateSpatialObjectMask(const Image &mask);
568 
569  template <class TImageType>
570  itk::ImageToImageMetricv4<TImageType,
571  TImageType,
572  TImageType,
573  double,
575  >* CreateMetric( );
576 
577  template <class TImageType>
578  void SetupMetric(
579  itk::ImageToImageMetricv4<TImageType,
580  TImageType,
581  TImageType,
582  double,
584  >*, const TImageType*, const TImageType* );
585 
586  template <typename TMetric>
588 
589  template <typename TTransformAdaptorPointer, typename TRegistrationMethod >
590  std::vector< TTransformAdaptorPointer >
591  CreateTransformParametersAdaptor(
592  TRegistrationMethod* method);
593 
594  virtual void PreUpdate( itk::ProcessObject *p );
595  virtual void OnActiveProcessDelete( ) throw();
596  virtual unsigned long AddITKObserver(const itk::EventObject &, itk::Command *);
597  virtual void RemoveITKObserver( EventCommand &e );
598 
599  private:
600 
601  nsstd::function<unsigned int()> m_pfGetOptimizerIteration;
602  nsstd::function<std::vector<double>()> m_pfGetOptimizerPosition;
603  nsstd::function<double()> m_pfGetOptimizerLearningRate;
604  nsstd::function<double()> m_pfGetOptimizerConvergenceValue;
605  nsstd::function<double()> m_pfGetMetricValue;
606  nsstd::function<std::vector<double>()> m_pfGetOptimizerScales;
607  nsstd::function<std::string()> m_pfGetOptimizerStopConditionDescription;
608 
609 
610  nsstd::function<unsigned int()> m_pfGetCurrentLevel;
611 
612  nsstd::function<void (itk::TransformBase *outTransform)> m_pfUpdateWithBestValue;
613 
614  template < class TMemberFunctionPointer >
616  {
617  typedef typename ::detail::FunctionTraits<TMemberFunctionPointer>::ClassType ObjectType;
618 
619  template< typename TImageType >
620  TMemberFunctionPointer operator() ( void ) const
621  {
622  return &ObjectType::template EvaluateInternal< TImageType >;
623  }
624  };
625 
626  typedef Transform (ImageRegistrationMethod::*MemberFunctionType)( const Image &fixed, const Image &moving );
627  typedef double (ImageRegistrationMethod::*EvaluateMemberFunctionType)( const Image &fixed, const Image &moving );
628  friend struct detail::MemberFunctionAddressor<MemberFunctionType>;
629  nsstd::auto_ptr<detail::MemberFunctionFactory<MemberFunctionType> > m_MemberFactory;
630  nsstd::auto_ptr<detail::MemberFunctionFactory<EvaluateMemberFunctionType> > m_EvaluateMemberFactory;
631 
637 
638  std::vector<uint32_t> m_VirtualDomainSize;
639  std::vector<double> m_VirtualDomainOrigin;
640  std::vector<double> m_VirtualDomainSpacing;
641  std::vector<double> m_VirtualDomainDirection;
642 
643  // optimizer
644  enum OptimizerType { ConjugateGradientLineSearch,
652  OnePlusOneEvolutionary
653  };
675  std::vector<unsigned int> m_OptimizerNumberOfSteps;
688  unsigned int m_OptimizerSeed;
689 
690  std::vector<double> m_OptimizerWeights;
691 
696  PhysicalShift
697  };
699  std::vector<double> m_OptimizerScales;
702 
703  // metric
704  enum MetricType { ANTSNeighborhoodCorrelation,
709  MattesMutualInformation
710  };
712  unsigned int m_MetricRadius;
716 
719 
720  std::vector<double> m_MetricSamplingPercentage;
722  unsigned int m_MetricSamplingSeed;
723 
726 
727  std::vector<unsigned int> m_ShrinkFactorsPerLevel;
728  std::vector<double> m_SmoothingSigmasPerLevel;
730 
733  unsigned int m_Iteration;
734 
736  };
737 
738 }
739 }
740 
741 #endif // sitkImageRegistrationMethod_h
std::vector< unsigned int > m_OptimizerNumberOfSteps
Never run estimation, use provided value.
Self & SetFixedInitialTransform(const Transform &transform)
Set transform mapping to the fixed domain.
STL namespace.
std::vector< unsigned int > m_ShrinkFactorsPerLevel
Transform GetInitialTransform()
Set the initial transform and parameters to optimize.
itk::ObjectToObjectOptimizerBaseTemplate< double > * m_ActiveOptimizer
#define SITKRegistration_EXPORT
::detail::FunctionTraits< TMemberFunctionPointer >::ClassType ObjectType
Self & SetInterpolator(InterpolatorEnum Interpolator)
Set and get the interpolator to use.
An interface method to the modular ITKv4 registration framework.
nsstd::auto_ptr< detail::MemberFunctionFactory< MemberFunctionType > > m_MemberFactory
Transform GetMovingInitialTransform() const
Set a fixed transform component towards moving domain.
An implementation of the Command design pattern for callback.
Definition: sitkCommand.h:44
A simplified wrapper around a variety of ITK transforms.
Definition: sitkTransform.h:83
Self & SmoothingSigmasAreSpecifiedInPhysicalUnitsOff()
Enable the smoothing sigmas for each level in physical units (default) or in terms of voxels...
Self & MetricUseFixedImageGradientFilterOff()
Enable image gradient computation by a filter.
Self & MetricUseMovingImageGradientFilterOff()
Enable image gradient computation by a filter.
The main Image class for SimpleITK.
Definition: sitkImage.h:54
Self & MetricUseFixedImageGradientFilterOn()
Enable image gradient computation by a filter.
Self & MetricUseMovingImageGradientFilterOn()
Enable image gradient computation by a filter.
bool GetInitialTransformInPlace() const
Set the initial transform and parameters to optimize.
Self & SetMovingInitialTransform(const Transform &transform)
Set a fixed transform component towards moving domain.
Estimate learning once each level, ignore provided values.
Self & SmoothingSigmasAreSpecifiedInPhysicalUnitsOn()
Enable the smoothing sigmas for each level in physical units (default) or in terms of voxels...
Base class for SimpleITK classes based on ProcessObject.
Transform GetFixedInitialTransform() const
Set transform mapping to the fixed domain.
class ITK_FORWARD_EXPORT Command
nsstd::auto_ptr< detail::MemberFunctionFactory< EvaluateMemberFunctionType > > m_EvaluateMemberFactory
InterpolatorEnum GetInterpolator()
Set and get the interpolator to use.
EstimateLearningRateType m_OptimizerEstimateLearningRate