ITK  4.9.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> >
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 
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 
259  itkGetConstMacro( OptimizerWeights, OptimizerWeightsType );
261 
263  itkSetObjectMacro( Metric, MetricType );
264  itkGetModifiableObjectMacro( Metric, MetricType );
266 
268  itkSetMacro( MetricSamplingStrategy, MetricSamplingStrategyType );
269  itkGetConstMacro( MetricSamplingStrategy, MetricSamplingStrategyType );
271 
274 
276  itkSetMacro( MetricSamplingPercentagePerLevel, MetricSamplingPercentageArrayType );
277  itkGetConstMacro( MetricSamplingPercentagePerLevel, MetricSamplingPercentageArrayType );
279 
281  itkSetGetDecoratedObjectInputMacro( FixedInitialTransform, InitialTransformType );
282 
284  itkSetGetDecoratedObjectInputMacro( MovingInitialTransform, InitialTransformType );
285 
301  itkSetGetDecoratedObjectInputMacro(InitialTransform, InitialTransformType);
302 
307 
315  void SetNumberOfLevels( const SizeValueType );
316  itkGetConstMacro( NumberOfLevels, SizeValueType );
318 
327  {
328  for( unsigned int level = 0; level < factors.Size(); ++level )
329  {
331  shrinkFactors.Fill( factors[level] );
332  this->SetShrinkFactorsPerDimension( level, shrinkFactors );
333  }
334  }
336 
341  {
342  if( level >= this->m_ShrinkFactorsPerLevel.size() )
343  {
344  itkExceptionMacro( "Requesting level greater than the number of levels." );
345  }
346  return this->m_ShrinkFactorsPerLevel[level];
347  }
349 
354  {
355  if( level >= this->m_ShrinkFactorsPerLevel.size() )
356  {
357  this->m_ShrinkFactorsPerLevel.resize( level + 1 );
358  }
359  this->m_ShrinkFactorsPerLevel[level] = factors;
360  this->Modified();
361  }
363 
369  itkSetMacro( SmoothingSigmasPerLevel, SmoothingSigmasArrayType );
370  itkGetConstMacro( SmoothingSigmasPerLevel, SmoothingSigmasArrayType );
372 
377  itkSetMacro( SmoothingSigmasAreSpecifiedInPhysicalUnits, bool );
378  itkGetConstMacro( SmoothingSigmasAreSpecifiedInPhysicalUnits, bool );
379  itkBooleanMacro( SmoothingSigmasAreSpecifiedInPhysicalUnits );
381 
386 
389  virtual const DecoratedOutputTransformType * GetOutput() const;
391 
393  virtual const DecoratedOutputTransformType * GetTransformOutput() const { return this->GetOutput(); }
394 
396  virtual const OutputTransformType * GetTransform() const;
397 
399  itkGetConstMacro( CurrentLevel, SizeValueType );
400 
402  itkGetConstReferenceMacro( CurrentIteration, SizeValueType );
403 
404  /* Get the current metric value. This is a helper function for reporting observations. */
405  itkGetConstReferenceMacro( CurrentMetricValue, RealType );
406 
408  itkGetConstReferenceMacro( CurrentConvergenceValue, RealType );
409 
411  itkGetConstReferenceMacro( IsConverged, bool );
412 
416  itkSetMacro( InPlace, bool );
417  itkGetConstMacro( InPlace, bool );
418  itkBooleanMacro( InPlace );
420 
426  itkBooleanMacro( InitializeCenterOfLinearOutputTransform );
427  itkSetMacro( InitializeCenterOfLinearOutputTransform, bool );
428  itkGetConstMacro( InitializeCenterOfLinearOutputTransform, bool );
430 
444 
445 #ifdef ITKV3_COMPATIBILITY
446 
461  void StartRegistration(void) { this->Update(); }
462 #endif
463 
464 protected:
466  virtual ~ImageRegistrationMethodv4();
467  virtual void PrintSelf( std::ostream & os, Indent indent ) const ITK_OVERRIDE;
468 
470  virtual void GenerateData() ITK_OVERRIDE;
471 
472  virtual void AllocateOutputs();
473 
475  virtual void InitializeRegistrationAtEachLevel( const SizeValueType );
476 
479 
481  virtual void SetMetricSamplePoints();
482 
483  SizeValueType m_CurrentLevel;
484  SizeValueType m_NumberOfLevels;
485  SizeValueType m_CurrentIteration;
489 
497  SizeValueType m_NumberOfFixedObjects;
498  SizeValueType m_NumberOfMovingObjects;
499 
503 
507  SizeValueType m_NumberOfMetrics;
512 
514 
516 
517  //TODO: m_OutputTransform should be removed and replaced with a named input parameter for
518  // the pipeline
520 
521 
522 private:
523  ImageRegistrationMethodv4( const Self & ) ITK_DELETE_FUNCTION;
524  void operator=( const Self & ) ITK_DELETE_FUNCTION;
525 
526  bool m_InPlace;
527 
529 
530  // helper function to create the right kind of concrete transform
531  template<typename TTransform>
532  static void MakeOutputTransform(SmartPointer<TTransform> &ptr)
533  {
534  ptr = TTransform::New();
535  }
536 
538  {
540  }
541 
542 };
543 } // end namespace itk
544 
545 #ifndef ITK_MANUAL_INSTANTIATION
546 #include "itkImageRegistrationMethodv4.hxx"
547 #endif
548 
549 #endif
virtual void Update()
Bring this filter up-to-date.
virtual OutputTransformType * GetModifiableTransform()
FixedImageMaskType::ConstPointer FixedImageMaskConstPointer
Transform< RealType, ImageDimension, ImageDimension > InitialTransformType
ImageMetricType::FixedSampledPointSetType MetricSamplePointSetType
Light weight base class for most itk classes.
virtual DecoratedOutputTransformType * GetTransformOutput()
ProcessObject::DataObjectPointerArraySizeType DataObjectPointerArraySizeType
virtual void SetMovingImage(const MovingImageType *image)
virtual void PrintSelf(std::ostream &os, Indent indent) const override
virtual void SetMovingPointSet(const PointSetType *pointSet)
MovingImageMasksContainerType m_MovingImageMasks
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 &)
std::vector< ShrinkFactorsPerDimensionContainerType > m_ShrinkFactorsPerLevel
void SetShrinkFactorsPerLevel(ShrinkFactorsArrayType factors)
static Pointer New()
MetricSamplingPercentageArrayType m_MetricSamplingPercentagePerLevel
std::vector< MovingImageMaskConstPointer > MovingImageMasksContainerType
ImageToImageMetricv4< FixedImageType, MovingImageType, VirtualImageType, RealType > ImageMetricType
virtual DataObjectPointer MakeOutput(DataObjectPointerArraySizeType) override
Make a DataObject of the correct type to used as the specified output.
void SetMetricSamplingPercentage(const RealType)
virtual void SetFixedPointSet(const PointSetType *pointSet)
virtual void GenerateData() override
MovingImagesContainerType m_MovingSmoothImages
virtual void SetMetricSamplePoints()
static void MakeOutputTransform(SmartPointer< TTransform > &ptr)
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 void InitializeRegistrationAtEachLevel(const SizeValueType)
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
void SetOptimizerWeights(OptimizerWeightsType &)
virtual DataObjectPointer MakeOutput(DataObjectPointerArraySizeType idx)
Make a DataObject of the correct type to used as the specified output.
ObjectToObjectOptimizerBaseTemplate< RealType > OptimizerType
void SetTransformParametersAdaptorsPerLevel(TransformParametersAdaptorsContainerType &)
OutputTransformType::DerivativeType DerivativeType
virtual void Modified() const
TransformParametersAdaptorsContainerType m_TransformParametersAdaptorsPerLevel
CompositeTransformPointer m_CompositeTransform
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
MetricSamplingStrategyType m_MetricSamplingStrategy
Base class for templated image classes.
Definition: itkImageBase.h:112
const TransformParametersAdaptorsContainerType & GetTransformParametersAdaptorsPerLevel() const
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.
SmoothingSigmasArrayType m_SmoothingSigmasPerLevel
SizeValueType Size(void) const
Definition: itkArray.h:124
virtual const PointSetType * GetMovingPointSet() const
std::vector< FixedImagePointer > FixedImagesContainerType
VirtualImageBaseType::ConstPointer VirtualImageBaseConstPointer
FixedImageMasksContainerType m_FixedImageMasks
Interface method for the current registration framework.
virtual DecoratedOutputTransformType * GetOutput()
PointSetToPointSetMetricv4< PointSetType, PointSetType, RealType > PointSetMetricType
Abstract base for object-to-object optimizers.
CompositeTransform< RealType, ImageDimension > CompositeTransformType
void SetNumberOfLevels(const SizeValueType)
OutputTransformType::ScalarType RealType
std::vector< PointSetConstPointer > PointSetsContainerType
OutputTransformType::Pointer OutputTransformPointer
ObjectToObjectMultiMetricv4< ImageDimension, ImageDimension, VirtualImageType, RealType > MultiMetricType
ImageBase< ImageDimension > VirtualImageBaseType
virtual const MovingImageType * GetMovingImage() const
virtual const OutputTransformType * GetTransform() const
VirtualImageType::Pointer VirtualImagePointer
virtual VirtualImageBaseConstPointer GetCurrentLevelVirtualDomainImage()