ITK  4.9.0
Insight Segmentation and Registration Toolkit
itkImageRegistrationMethod.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 itkImageRegistrationMethod_h
19 #define itkImageRegistrationMethod_h
20 
21 #include "itkProcessObject.h"
22 #include "itkImage.h"
23 #include "itkImageToImageMetric.h"
25 #include "itkDataObjectDecorator.h"
26 
27 namespace itk
28 {
69 template< typename TFixedImage, typename TMovingImage >
71 {
72 public:
78 
80  itkNewMacro(Self);
81 
84 
86  typedef TFixedImage FixedImageType;
87  typedef typename FixedImageType::ConstPointer FixedImageConstPointer;
88 
90  typedef TMovingImage MovingImageType;
91  typedef typename MovingImageType::ConstPointer MovingImageConstPointer;
92 
97 
101 
107 
111 
114 
118 
121 
123  void SetFixedImage(const FixedImageType *fixedImage);
124  itkGetConstObjectMacro(FixedImage, FixedImageType);
126 
128  void SetMovingImage(const MovingImageType *movingImage);
129  itkGetConstObjectMacro(MovingImage, MovingImageType);
131 
133  itkSetObjectMacro(Optimizer, OptimizerType);
134  itkGetModifiableObjectMacro(Optimizer, OptimizerType);
136 
138  itkSetObjectMacro(Metric, MetricType);
139  itkGetModifiableObjectMacro(Metric, MetricType);
141 
143  itkSetObjectMacro(Transform, TransformType);
144  itkGetModifiableObjectMacro(Transform, TransformType);
146 
148  itkSetObjectMacro(Interpolator, InterpolatorType);
149  itkGetModifiableObjectMacro(Interpolator, InterpolatorType);
151 
153  virtual void SetInitialTransformParameters(const ParametersType & param);
154 
155  itkGetConstReferenceMacro(InitialTransformParameters, ParametersType);
156 
159  itkGetConstReferenceMacro(LastTransformParameters, ParametersType);
160 
168  void SetFixedImageRegion(const FixedImageRegionType & region);
169 
174  itkGetConstReferenceMacro(FixedImageRegion, FixedImageRegionType);
175 
178  itkGetConstMacro(FixedImageRegionDefined, bool);
179 
184  itkSetMacro(FixedImageRegionDefined, bool);
185 
187  virtual void Initialize()
188  throw ( ExceptionObject );
189 
191  const TransformOutputType * GetOutput() const;
192 
196  using Superclass::MakeOutput;
197  virtual DataObjectPointer MakeOutput(DataObjectPointerArraySizeType idx) ITK_OVERRIDE;
198 
201  virtual ModifiedTimeType GetMTime() const ITK_OVERRIDE;
202 
203 #ifdef ITKV3_COMPATIBILITY
204 
219  void StartRegistration(void) { this->Update(); }
220 #endif
221 
222 protected:
225  virtual void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
228  virtual void GenerateData() ITK_OVERRIDE;
229 
231  itkSetMacro(LastTransformParameters, ParametersType);
232 
233  /* Start the Optimization */
234  void StartOptimization();
235 
236 private:
237  ImageRegistrationMethod(const Self &) ITK_DELETE_FUNCTION;
238  void operator=(const Self &) ITK_DELETE_FUNCTION;
239 
242 
245 
248 
251 
254 };
255 } // end namespace itk
256 
257 #ifndef ITK_MANUAL_INSTANTIATION
258 #include "itkImageRegistrationMethod.hxx"
259 #endif
260 
261 #endif
virtual void Update()
Bring this filter up-to-date.
void SetMovingImage(const MovingImageType *movingImage)
FixedImageType::ConstPointer FixedImageConstPointer
virtual ModifiedTimeType GetMTime() const override
Light weight base class for most itk classes.
MetricType::FixedImageRegionType FixedImageRegionType
unsigned long ModifiedTimeType
Definition: itkIntTypes.h:164
virtual void PrintSelf(std::ostream &os, Indent indent) const override
This class is a base for the Optimization methods that optimize a single valued function.
TransformOutputType::ConstPointer TransformOutputConstPointer
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
DataObjectDecorator< TransformType > TransformOutputType
virtual void GenerateData() override
const TransformOutputType * GetOutput() const
ImageToImageMetric< FixedImageType, MovingImageType > MetricType
TransformOutputType::Pointer TransformOutputPointer
MetricType::TransformType TransformType
Base class for Image Registration Methods.
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:82
DataObjectPointerArray::size_type DataObjectPointerArraySizeType
Class to hold and manage different parameter types used during optimization.
virtual void SetInitialTransformParameters(const ParametersType &param)
Decorates any subclass of itkObject with a DataObject API.
MovingImageType::ConstPointer MovingImageConstPointer
void SetFixedImage(const FixedImageType *fixedImage)
Standard exception handling object.
Generic representation for an optimization method.
Definition: itkOptimizer.h:38
virtual DataObjectPointer MakeOutput(DataObjectPointerArraySizeType idx) override
Make a DataObject of the correct type to used as the specified output.
void SetFixedImageRegion(const FixedImageRegionType &region)
Base class for all image interpolaters.
SingleValuedNonLinearOptimizer OptimizerType
MetricType::TransformParametersType ParametersType
MetricType::InterpolatorType InterpolatorType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
InterpolatorType::Pointer InterpolatorPointer
Computes similarity between regions of two images.
SmartPointer< const Self > ConstPointer
FixedImageType::RegionType FixedImageRegionType