ITK  5.0.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 >
70 class ITK_TEMPLATE_EXPORT ImageRegistrationMethod:public ProcessObject
71 {
72 public:
73  ITK_DISALLOW_COPY_AND_ASSIGN(ImageRegistrationMethod);
74 
80 
82  itkNewMacro(Self);
83 
86 
88  using FixedImageType = TFixedImage;
89  using FixedImageConstPointer = typename FixedImageType::ConstPointer;
90 
92  using MovingImageType = TMovingImage;
93  using MovingImageConstPointer = typename MovingImageType::ConstPointer;
94 
99 
102  using TransformPointer = typename TransformType::Pointer;
103 
109 
112  using InterpolatorPointer = typename InterpolatorType::Pointer;
113 
116 
120 
123 
125  void SetFixedImage(const FixedImageType *fixedImage);
126  itkGetConstObjectMacro(FixedImage, FixedImageType);
128 
130  void SetMovingImage(const MovingImageType *movingImage);
131  itkGetConstObjectMacro(MovingImage, MovingImageType);
133 
135  itkSetObjectMacro(Optimizer, OptimizerType);
136  itkGetModifiableObjectMacro(Optimizer, OptimizerType);
138 
140  itkSetObjectMacro(Metric, MetricType);
141  itkGetModifiableObjectMacro(Metric, MetricType);
143 
145  itkSetObjectMacro(Transform, TransformType);
146  itkGetModifiableObjectMacro(Transform, TransformType);
148 
150  itkSetObjectMacro(Interpolator, InterpolatorType);
151  itkGetModifiableObjectMacro(Interpolator, InterpolatorType);
153 
155  virtual void SetInitialTransformParameters(const ParametersType & param);
156 
157  itkGetConstReferenceMacro(InitialTransformParameters, ParametersType);
158 
161  itkGetConstReferenceMacro(LastTransformParameters, ParametersType);
162 
170  void SetFixedImageRegion(const FixedImageRegionType & region);
171 
176  itkGetConstReferenceMacro(FixedImageRegion, FixedImageRegionType);
177 
180  itkGetConstMacro(FixedImageRegionDefined, bool);
181 
186  itkSetMacro(FixedImageRegionDefined, bool);
187 
189  virtual void Initialize();
190 
192  const TransformOutputType * GetOutput() const;
193 
197  using Superclass::MakeOutput;
198  DataObjectPointer MakeOutput(DataObjectPointerArraySizeType idx) override;
199 
202  ModifiedTimeType GetMTime() const override;
203 
204 protected:
206  ~ImageRegistrationMethod() override = default;
207  void PrintSelf(std::ostream & os, Indent indent) const override;
210  void GenerateData() override;
211 
213  itkSetMacro(LastTransformParameters, ParametersType);
214 
215  /* Start the Optimization */
216  void StartOptimization();
217 
218 private:
221 
224 
227 
230 
233 };
234 } // end namespace itk
235 
236 #ifndef ITK_MANUAL_INSTANTIATION
237 #include "itkImageRegistrationMethod.hxx"
238 #endif
239 
240 #endif
typename FixedImageType::RegionType FixedImageRegionType
Light weight base class for most itk classes.
typename TransformOutputType::ConstPointer TransformOutputConstPointer
typename MetricType::TransformParametersType ParametersType
DataObjectPointerArray::size_type DataObjectPointerArraySizeType
typename MetricType::TransformType TransformType
typename TransformType::Pointer TransformPointer
This class is a base for the Optimization methods that optimize a single valued function.
typename MetricType::Pointer MetricPointer
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
typename TransformOutputType::Pointer TransformOutputPointer
class ITK_FORWARD_EXPORT ProcessObject
Definition: itkDataObject.h:40
Base class for Image Registration Methods.
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:83
typename MovingImageType::ConstPointer MovingImageConstPointer
Decorates any subclass of itkObject with a DataObject API.
typename TransformType::ParametersType TransformParametersType
typename FixedImageType::ConstPointer FixedImageConstPointer
Generic representation for an optimization method.
Definition: itkOptimizer.h:38
unsigned long ModifiedTimeType
Definition: itkIntTypes.h:104
Base class for all image interpolaters.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
typename MetricType::InterpolatorType InterpolatorType
typename InterpolatorType::Pointer InterpolatorPointer
Computes similarity between regions of two images.
SmartPointer< Self > Pointer
typename MetricType::FixedImageRegionType FixedImageRegionType