ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkImageRegistrationMethodv4.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 itkImageRegistrationMethodv4_h
19 #define itkImageRegistrationMethodv4_h
20 
21 #include "itkProcessObject.h"
22 
23 #include "itkCompositeTransform.h"
24 #include "itkDataObjectDecorator.h"
30 #include "itkShrinkImageFilter.h"
31 #include "itkIdentityTransform.h"
33 
34 #include <vector>
35 
36 namespace itk
37 {
38 
89 template<typename TFixedImage,
90  typename TMovingImage,
91  typename TOutputTransform = Transform<double, TFixedImage::ImageDimension, TFixedImage::ImageDimension>,
92  typename TVirtualImage = TFixedImage,
93  typename TPointSet = PointSet<unsigned int, TFixedImage::ImageDimension> >
94 class ITK_TEMPLATE_EXPORT ImageRegistrationMethodv4
95 :public ProcessObject
96 {
97 public:
103 
105  itkNewMacro( Self );
106 
108  itkStaticConstMacro( ImageDimension, unsigned int, TFixedImage::ImageDimension );
109 
112 
114  typedef TFixedImage FixedImageType;
115  typedef typename FixedImageType::Pointer FixedImagePointer;
116  typedef std::vector<FixedImagePointer> FixedImagesContainerType;
117  typedef TMovingImage MovingImageType;
118  typedef typename MovingImageType::Pointer MovingImagePointer;
119  typedef std::vector<MovingImagePointer> MovingImagesContainerType;
120 
121  typedef TPointSet PointSetType;
122  typedef typename PointSetType::ConstPointer PointSetConstPointer;
123  typedef std::vector<PointSetConstPointer> PointSetsContainerType;
124 
126  typedef TOutputTransform OutputTransformType;
127  typedef typename OutputTransformType::Pointer OutputTransformPointer;
128  typedef typename OutputTransformType::ScalarType RealType;
129  typedef typename OutputTransformType::DerivativeType DerivativeType;
130  typedef typename DerivativeType::ValueType DerivativeValueType;
131 
134 
137 
140 
142 
143  typedef TVirtualImage VirtualImageType;
144  typedef typename VirtualImageType::Pointer VirtualImagePointer;
147 
151 
154  typedef std::vector<FixedImageMaskConstPointer> FixedImageMasksContainerType;
157  typedef std::vector<MovingImageMaskConstPointer> MovingImageMasksContainerType;
158 
167 
170 
172 
175 
179  typedef std::vector<TransformParametersAdaptorPointer> TransformParametersAdaptorsContainerType;
180 
184 
187 
189  enum MetricSamplingStrategyType { NONE, REGULAR, RANDOM };
190 
192 
194  virtual void SetFixedImage( const FixedImageType *image )
195  {
196  this->SetFixedImage( 0, image );
197  }
198  virtual const FixedImageType * GetFixedImage() const
199  {
200  return this->GetFixedImage( 0 );
201  }
202  virtual void SetFixedImage( SizeValueType, const FixedImageType * );
203  virtual const FixedImageType * GetFixedImage( SizeValueType ) const;
205 
207  virtual void SetMovingImage( const MovingImageType *image )
208  {
209  this->SetMovingImage( 0, image );
210  }
211  virtual const MovingImageType * GetMovingImage() const
212  {
213  return this->GetMovingImage( 0 );
214  }
215  virtual void SetMovingImage( SizeValueType, const MovingImageType * );
216  virtual const MovingImageType * GetMovingImage( SizeValueType ) const;
218 
220  virtual void SetFixedPointSet( const PointSetType *pointSet )
221  {
222  this->SetFixedPointSet( 0, pointSet );
223  }
224  virtual const PointSetType * GetFixedPointSet() const
225  {
226  return this->GetFixedPointSet( 0 );
227  }
228  virtual void SetFixedPointSet( SizeValueType, const PointSetType * );
229  virtual const PointSetType * GetFixedPointSet( SizeValueType ) const;
231 
233  virtual void SetMovingPointSet( const PointSetType *pointSet )
234  {
235  this->SetMovingPointSet( 0, pointSet );
236  }
237  virtual const PointSetType * GetMovingPointSet() const
238  {
239  return this->GetMovingPointSet( 0 );
240  }
241  virtual void SetMovingPointSet( SizeValueType, const PointSetType * );
242  virtual const PointSetType * GetMovingPointSet( SizeValueType ) const;
244 
246  itkSetObjectMacro( Optimizer, OptimizerType );
247  itkGetModifiableObjectMacro( Optimizer, OptimizerType );
249 
258  void SetOptimizerWeights( OptimizerWeightsType & );
259  itkGetConstMacro( OptimizerWeights, OptimizerWeightsType );
261 
263  itkSetObjectMacro( Metric, MetricType );
264  itkGetModifiableObjectMacro( Metric, MetricType );
266 
268  itkSetMacro( MetricSamplingStrategy, MetricSamplingStrategyType );
269  itkGetConstMacro( MetricSamplingStrategy, MetricSamplingStrategyType );
271 
283  void MetricSamplingReinitializeSeed();
284  void MetricSamplingReinitializeSeed(int seed);
286 
288  void SetMetricSamplingPercentage( const RealType );
289 
291  virtual void SetMetricSamplingPercentagePerLevel( const MetricSamplingPercentageArrayType &samplingPercentages );
292  itkGetConstMacro( MetricSamplingPercentagePerLevel, MetricSamplingPercentageArrayType );
294 
296  itkSetGetDecoratedObjectInputMacro( FixedInitialTransform, InitialTransformType );
297 
299  itkSetGetDecoratedObjectInputMacro( MovingInitialTransform, InitialTransformType );
300 
316  itkSetGetDecoratedObjectInputMacro(InitialTransform, InitialTransformType);
317 
319  void SetTransformParametersAdaptorsPerLevel( TransformParametersAdaptorsContainerType & );
320  const TransformParametersAdaptorsContainerType & GetTransformParametersAdaptorsPerLevel() const;
322 
330  void SetNumberOfLevels( const SizeValueType );
331  itkGetConstMacro( NumberOfLevels, SizeValueType );
333 
342  {
343  for( unsigned int level = 0; level < factors.Size(); ++level )
344  {
346  shrinkFactors.Fill( factors[level] );
347  this->SetShrinkFactorsPerDimension( level, shrinkFactors );
348  }
349  }
351 
356  {
357  if( level >= this->m_ShrinkFactorsPerLevel.size() )
358  {
359  itkExceptionMacro( "Requesting level greater than the number of levels." );
360  }
361  return this->m_ShrinkFactorsPerLevel[level];
362  }
364 
369  {
370  if( level >= this->m_ShrinkFactorsPerLevel.size() )
371  {
372  this->m_ShrinkFactorsPerLevel.resize( level + 1 );
373  }
374  this->m_ShrinkFactorsPerLevel[level] = factors;
375  this->Modified();
376  }
378 
384  itkSetMacro( SmoothingSigmasPerLevel, SmoothingSigmasArrayType );
385  itkGetConstMacro( SmoothingSigmasPerLevel, SmoothingSigmasArrayType );
387 
392  itkSetMacro( SmoothingSigmasAreSpecifiedInPhysicalUnits, bool );
393  itkGetConstMacro( SmoothingSigmasAreSpecifiedInPhysicalUnits, bool );
394  itkBooleanMacro( SmoothingSigmasAreSpecifiedInPhysicalUnits );
396 
399  using Superclass::MakeOutput;
400  virtual DataObjectPointer MakeOutput( DataObjectPointerArraySizeType ) ITK_OVERRIDE;
401 
403  virtual DecoratedOutputTransformType * GetOutput();
404  virtual const DecoratedOutputTransformType * GetOutput() const;
406 
407  virtual DecoratedOutputTransformType * GetTransformOutput() { return this->GetOutput(); }
408  virtual const DecoratedOutputTransformType * GetTransformOutput() const { return this->GetOutput(); }
409 
410  virtual OutputTransformType * GetModifiableTransform();
411  virtual const OutputTransformType * GetTransform() const;
412 
414  itkGetConstMacro( CurrentLevel, SizeValueType );
415 
417  itkGetConstReferenceMacro( CurrentIteration, SizeValueType );
418 
419  /* Get the current metric value. This is a helper function for reporting observations. */
420  itkGetConstReferenceMacro( CurrentMetricValue, RealType );
421 
423  itkGetConstReferenceMacro( CurrentConvergenceValue, RealType );
424 
426  itkGetConstReferenceMacro( IsConverged, bool );
427 
431  itkSetMacro( InPlace, bool );
432  itkGetConstMacro( InPlace, bool );
433  itkBooleanMacro( InPlace );
435 
441  itkBooleanMacro( InitializeCenterOfLinearOutputTransform );
442  itkSetMacro( InitializeCenterOfLinearOutputTransform, bool );
443  itkGetConstMacro( InitializeCenterOfLinearOutputTransform, bool );
445 
458  void InitializeCenterOfLinearOutputTransform();
459 
460 #ifdef ITKV3_COMPATIBILITY
461 
476  void StartRegistration(void) { this->Update(); }
477 #endif
478 
479 protected:
480  ImageRegistrationMethodv4();
481  virtual ~ImageRegistrationMethodv4() ITK_OVERRIDE;
482  virtual void PrintSelf( std::ostream & os, Indent indent ) const ITK_OVERRIDE;
483 
485  virtual void GenerateData() ITK_OVERRIDE;
486 
487  virtual void AllocateOutputs();
488 
490  virtual void InitializeRegistrationAtEachLevel( const SizeValueType );
491 
493  virtual VirtualImageBaseConstPointer GetCurrentLevelVirtualDomainImage();
494 
496  virtual void SetMetricSamplePoints();
497 
498  SizeValueType m_CurrentLevel;
499  SizeValueType m_NumberOfLevels;
500  SizeValueType m_CurrentIteration;
501  RealType m_CurrentMetricValue;
502  RealType m_CurrentConvergenceValue;
503  bool m_IsConverged;
504 
505  FixedImagesContainerType m_FixedSmoothImages;
506  MovingImagesContainerType m_MovingSmoothImages;
507  FixedImageMasksContainerType m_FixedImageMasks;
508  MovingImageMasksContainerType m_MovingImageMasks;
509  VirtualImagePointer m_VirtualDomainImage;
510  PointSetsContainerType m_FixedPointSets;
511  PointSetsContainerType m_MovingPointSets;
512  SizeValueType m_NumberOfFixedObjects;
513  SizeValueType m_NumberOfMovingObjects;
514 
515  OptimizerPointer m_Optimizer;
516  OptimizerWeightsType m_OptimizerWeights;
517  bool m_OptimizerWeightsAreIdentity;
518 
519  MetricPointer m_Metric;
520  MetricSamplingStrategyType m_MetricSamplingStrategy;
521  MetricSamplingPercentageArrayType m_MetricSamplingPercentagePerLevel;
522  SizeValueType m_NumberOfMetrics;
523  int m_FirstImageMetricIndex;
524  std::vector<ShrinkFactorsPerDimensionContainerType> m_ShrinkFactorsPerLevel;
525  SmoothingSigmasArrayType m_SmoothingSigmasPerLevel;
526  bool m_SmoothingSigmasAreSpecifiedInPhysicalUnits;
527 
528  bool m_ReseedIterator;
529  int m_RandomSeed;
530  int m_CurrentRandomSeed;
531 
532 
533  TransformParametersAdaptorsContainerType m_TransformParametersAdaptorsPerLevel;
534 
535  CompositeTransformPointer m_CompositeTransform;
536 
537  //TODO: m_OutputTransform should be removed and replaced with a named input parameter for
538  // the pipeline
539  OutputTransformPointer m_OutputTransform;
540 
541 
542 private:
543  ITK_DISALLOW_COPY_AND_ASSIGN(ImageRegistrationMethodv4);
544 
545  bool m_InPlace;
546 
547  bool m_InitializeCenterOfLinearOutputTransform;
548 
549  // helper function to create the right kind of concrete transform
550  template<typename TTransform>
551  static void MakeOutputTransform(SmartPointer<TTransform> &ptr)
552  {
553  ptr = TTransform::New();
554  }
555 
557  {
559  }
560 
561 };
562 } // end namespace itk
563 
564 #ifndef ITK_MANUAL_INSTANTIATION
565 #include "itkImageRegistrationMethodv4.hxx"
566 #endif
567 
568 #endif
FixedImageMaskType::ConstPointer FixedImageMaskConstPointer
Transform< RealType, ImageDimension, ImageDimension > InitialTransformType
ImageMetricType::FixedSampledPointSetType MetricSamplePointSetType
Light weight base class for most itk classes.
ProcessObject::DataObjectPointerArraySizeType DataObjectPointerArraySizeType
virtual void SetMovingImage(const MovingImageType *image)
virtual void SetMovingPointSet(const PointSetType *pointSet)
DerivativeType::ValueType DerivativeValueType
This class takes one ore more ObjectToObject metrics and assigns weights to their derivatives to comp...
virtual const PointSetType * GetFixedPointSet() const
ShrinkFilterType::ShrinkFactorsType ShrinkFactorsPerDimensionContainerType
Computes similarity between two point sets.
TransformParametersAdaptorType::Pointer TransformParametersAdaptorPointer
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
MovingImageMaskType::ConstPointer MovingImageMaskConstPointer
DecoratedInitialTransformType::Pointer DecoratedInitialTransformPointer
CompositeTransformType::Pointer CompositeTransformPointer
ImageMetricType::FixedImageMaskType FixedImageMaskType
ObjectType * GetPointer() const
Base helper class intended for multi-resolution image registration.
Vector< RealType, ImageDimension > VectorType
Base class for all object-to-object similarlity metrics added in ITKv4.
unsigned long SizeValueType
Definition: itkIntTypes.h:143
void Fill(const ValueType &)
void SetShrinkFactorsPerLevel(ShrinkFactorsArrayType factors)
static Pointer New()
std::vector< MovingImageMaskConstPointer > MovingImageMasksContainerType
ImageToImageMetricv4< FixedImageType, MovingImageType, VirtualImageType, RealType > ImageMetricType
virtual void SetFixedPointSet(const PointSetType *pointSet)
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:82
ShrinkFactorsPerDimensionContainerType GetShrinkFactorsPerDimension(const unsigned int level) const
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:62
std::vector< FixedImageMaskConstPointer > FixedImageMasksContainerType
DataObjectPointerArray::size_type DataObjectPointerArraySizeType
virtual const FixedImageType * GetFixedImage() const
std::vector< TransformParametersAdaptorPointer > TransformParametersAdaptorsContainerType
This class contains a list of transforms and concatenates them by composition.
ImageMetricType::MovingImageMaskType MovingImageMaskType
Decorates any subclass of itkObject with a DataObject API.
static void MakeOutputTransform(SmartPointer< InitialTransformType > &ptr)
Implementation of the composite pattern.
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition: itkPointSet.h:84
InitialTransformType::Pointer InitialTransformPointer
Generic representation for an optimization method.
Definition: itkOptimizer.h:38
OptimizerType::ScalesType OptimizerWeightsType
DecoratedOutputTransformType::Pointer DecoratedOutputTransformPointer
ObjectToObjectOptimizerBaseTemplate< RealType > OptimizerType
OutputTransformType::DerivativeType DerivativeType
virtual void SetFixedImage(const FixedImageType *image)
void SetShrinkFactorsPerDimension(unsigned int level, ShrinkFactorsPerDimensionContainerType factors)
PointSetType::ConstPointer PointSetConstPointer
TransformParametersAdaptorBase< InitialTransformType > TransformParametersAdaptorType
std::vector< MovingImagePointer > MovingImagesContainerType
ObjectToObjectMetricBaseTemplate< RealType > MetricType
Base class for templated image classes.
Definition: itkImageBase.h:114
DataObjectDecorator< InitialTransformType > DecoratedInitialTransformType
virtual const DecoratedOutputTransformType * GetTransformOutput() const
DataObjectDecorator< OutputTransformType > DecoratedOutputTransformType
ShrinkImageFilter< FixedImageType, VirtualImageType > ShrinkFilterType
Reduce the size of an image by an integer factor in each dimension.
SizeValueType Size(void) const
Definition: itkArray.h:124
virtual const PointSetType * GetMovingPointSet() const
std::vector< FixedImagePointer > FixedImagesContainerType
VirtualImageBaseType::ConstPointer VirtualImageBaseConstPointer
Interface method for the current registration framework.
PointSetToPointSetMetricv4< PointSetType, PointSetType, RealType > PointSetMetricType
Abstract base for object-to-object optimizers.
CompositeTransform< RealType, ImageDimension > CompositeTransformType
OutputTransformType::ScalarType RealType
std::vector< PointSetConstPointer > PointSetsContainerType
OutputTransformType::Pointer OutputTransformPointer
ObjectToObjectMultiMetricv4< ImageDimension, ImageDimension, VirtualImageType, RealType > MultiMetricType
ImageBase< ImageDimension > VirtualImageBaseType
virtual const MovingImageType * GetMovingImage() const
VirtualImageType::Pointer VirtualImagePointer