ITK  6.0.0
Insight Toolkit
itkImageRegistrationMethod.h
Go to the documentation of this file.
1/*=========================================================================
2 *
3 * Copyright NumFOCUS
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 * https://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"
26
27namespace itk
28{
69template <typename TFixedImage, typename TMovingImage>
70class ITK_TEMPLATE_EXPORT ImageRegistrationMethod : public ProcessObject
71{
72public:
73 ITK_DISALLOW_COPY_AND_MOVE(ImageRegistrationMethod);
74
80
82 itkNewMacro(Self);
83
85 itkOverrideGetNameOfClassMacro(ImageRegistrationMethod);
86
88 using FixedImageType = TFixedImage;
90
92 using MovingImageType = TMovingImage;
94
99
103
109
113
116
120
123
125 void
126 SetFixedImage(const FixedImageType * fixedImage);
127 itkGetConstObjectMacro(FixedImage, FixedImageType);
131 void
132 SetMovingImage(const MovingImageType * movingImage);
133 itkGetConstObjectMacro(MovingImage, MovingImageType);
137 itkSetObjectMacro(Optimizer, OptimizerType);
138 itkGetModifiableObjectMacro(Optimizer, OptimizerType);
142 itkSetObjectMacro(Metric, MetricType);
143 itkGetModifiableObjectMacro(Metric, MetricType);
147 itkSetObjectMacro(Transform, TransformType);
148 itkGetModifiableObjectMacro(Transform, TransformType);
152 itkSetObjectMacro(Interpolator, InterpolatorType);
153 itkGetModifiableObjectMacro(Interpolator, InterpolatorType);
157 virtual void
159
160 itkGetConstReferenceMacro(InitialTransformParameters, ParametersType);
161
164 itkGetConstReferenceMacro(LastTransformParameters, ParametersType);
165
173 void
175
180 itkGetConstReferenceMacro(FixedImageRegion, FixedImageRegionType);
181
184 itkGetConstMacro(FixedImageRegionDefined, bool);
185
190 itkSetMacro(FixedImageRegionDefined, bool);
191
193 virtual void
195
197 const TransformOutputType *
198 GetOutput() const;
199
203 using Superclass::MakeOutput;
206
210 GetMTime() const override;
211
212protected:
214 ~ImageRegistrationMethod() override = default;
215 void
216 PrintSelf(std::ostream & os, Indent indent) const override;
219 void
220 GenerateData() override;
221
223 itkSetMacro(LastTransformParameters, ParametersType);
224
226 void
228
229private:
230 MetricPointer m_Metric{};
232
233 MovingImageConstPointer m_MovingImage{};
235
236 TransformPointer m_Transform{};
237 InterpolatorPointer m_Interpolator{};
238
239 ParametersType m_InitialTransformParameters{};
240 ParametersType m_LastTransformParameters{};
241
242 bool m_FixedImageRegionDefined{};
243 FixedImageRegionType m_FixedImageRegion{};
244};
245} // end namespace itk
246
247#ifndef ITK_MANUAL_INSTANTIATION
248# include "itkImageRegistrationMethod.hxx"
249#endif
250
251#endif
Decorates any subclass of itkObject with a DataObject API.
SmartPointer< Self > Pointer
Base class for Image Registration Methods.
typename MetricType::FixedImageRegionType FixedImageRegionType
~ImageRegistrationMethod() override=default
void SetFixedImage(const FixedImageType *fixedImage)
typename InterpolatorType::Pointer InterpolatorPointer
typename TransformType::Pointer TransformPointer
typename TransformOutputType::Pointer TransformOutputPointer
typename FixedImageType::ConstPointer FixedImageConstPointer
virtual void SetInitialTransformParameters(const ParametersType &param)
void SetMovingImage(const MovingImageType *movingImage)
typename MetricType::TransformParametersType ParametersType
typename MovingImageType::ConstPointer MovingImageConstPointer
ModifiedTimeType GetMTime() const override
typename MetricType::TransformType TransformType
void GenerateData() override
void PrintSelf(std::ostream &os, Indent indent) const override
typename MetricType::InterpolatorType InterpolatorType
void SetFixedImageRegion(const FixedImageRegionType &region)
typename MetricType::Pointer MetricPointer
const TransformOutputType * GetOutput() const
typename TransformOutputType::ConstPointer TransformOutputConstPointer
DataObjectPointer MakeOutput(DataObjectPointerArraySizeType output) override
Make a DataObject of the correct type to used as the specified output.
Computes similarity between regions of two images.
typename FixedImageType::RegionType FixedImageRegionType
typename TransformType::ParametersType TransformParametersType
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Base class for all image interpolators.
Light weight base class for most itk classes.
Generic representation for an optimization method.
Definition: itkOptimizer.h:39
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
DataObjectPointerArray::size_type DataObjectPointerArraySizeType
This class is a base for the Optimization methods that optimize a single valued function.
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:84
SmartPointer< const Self > ConstPointer
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
class ITK_FORWARD_EXPORT ProcessObject
Definition: itkDataObject.h:41
SizeValueType ModifiedTimeType
Definition: itkIntTypes.h:105