ITK  6.0.0
Insight Toolkit
itkMultiResolutionImageRegistrationMethod.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 itkMultiResolutionImageRegistrationMethod_h
19#define itkMultiResolutionImageRegistrationMethod_h
20
21#include "itkProcessObject.h"
25#include "itkNumericTraits.h"
27
28namespace itk
29{
71template <typename TFixedImage, typename TMovingImage>
72class ITK_TEMPLATE_EXPORT MultiResolutionImageRegistrationMethod : public ProcessObject
73{
74public:
75 ITK_DISALLOW_COPY_AND_MOVE(MultiResolutionImageRegistrationMethod);
76
82
84 itkNewMacro(Self);
85
87 itkOverrideGetNameOfClassMacro(MultiResolutionImageRegistrationMethod);
88
90 using FixedImageType = TFixedImage;
93
95 using MovingImageType = TMovingImage;
97
101
105
111
115
118
122
125
129
133
136
138 void
140
142 itkSetConstObjectMacro(FixedImage, FixedImageType);
143 itkGetConstObjectMacro(FixedImage, FixedImageType);
147 itkSetConstObjectMacro(MovingImage, MovingImageType);
148 itkGetConstObjectMacro(MovingImage, MovingImageType);
152 itkSetObjectMacro(Optimizer, OptimizerType);
153 itkGetModifiableObjectMacro(Optimizer, OptimizerType);
157 itkSetObjectMacro(Metric, MetricType);
158 itkGetModifiableObjectMacro(Metric, MetricType);
162 itkSetMacro(FixedImageRegion, FixedImageRegionType);
163 itkGetConstReferenceMacro(FixedImageRegion, FixedImageRegionType);
167 itkSetObjectMacro(Transform, TransformType);
168 itkGetModifiableObjectMacro(Transform, TransformType);
172 itkSetObjectMacro(Interpolator, InterpolatorType);
173 itkGetModifiableObjectMacro(Interpolator, InterpolatorType);
177 itkSetObjectMacro(FixedImagePyramid, FixedImagePyramidType);
178 itkGetModifiableObjectMacro(FixedImagePyramid, FixedImagePyramidType);
182 itkSetObjectMacro(MovingImagePyramid, MovingImagePyramidType);
183 itkGetModifiableObjectMacro(MovingImagePyramid, MovingImagePyramidType);
187 void
188 SetSchedules(const ScheduleType & fixedImagePyramidSchedule, const ScheduleType & movingImagePyramidSchedule);
189
190 itkGetConstMacro(FixedImagePyramidSchedule, ScheduleType);
191 itkGetConstMacro(MovingImagePyramidSchedule, ScheduleType);
192
194 void
196
197 itkGetConstMacro(NumberOfLevels, SizeValueType);
198
200 itkGetConstMacro(CurrentLevel, SizeValueType);
201
203 itkSetMacro(InitialTransformParameters, ParametersType);
204 itkGetConstReferenceMacro(InitialTransformParameters, ParametersType);
210 itkSetMacro(InitialTransformParametersOfNextLevel, ParametersType);
211 itkGetConstReferenceMacro(InitialTransformParametersOfNextLevel, ParametersType);
216 itkGetConstReferenceMacro(LastTransformParameters, ParametersType);
217
219 const TransformOutputType *
220 GetOutput() const;
221
225 using Superclass::MakeOutput;
228
232 GetMTime() const override;
233
234protected:
237 void
238 PrintSelf(std::ostream & os, Indent indent) const override;
239
242 void
243 GenerateData() override;
244
249 void
251
253 void
255
257 itkSetMacro(CurrentLevel, SizeValueType);
258
259private:
260 MetricPointer m_Metric{};
262
263 MovingImageConstPointer m_MovingImage{};
265
266 TransformPointer m_Transform{};
267 InterpolatorPointer m_Interpolator{};
268
269 MovingImagePyramidPointer m_MovingImagePyramid{};
270 FixedImagePyramidPointer m_FixedImagePyramid{};
271
272 ParametersType m_InitialTransformParameters{};
273 ParametersType m_InitialTransformParametersOfNextLevel{};
274 ParametersType m_LastTransformParameters{};
275
276 FixedImageRegionType m_FixedImageRegion{};
277 std::vector<FixedImageRegionType> m_FixedImageRegionPyramid{};
278
279 SizeValueType m_NumberOfLevels{};
280 SizeValueType m_CurrentLevel{};
281
282 bool m_Stop{};
283
284 ScheduleType m_FixedImagePyramidSchedule{};
285 ScheduleType m_MovingImagePyramidSchedule{};
286
287 bool m_ScheduleSpecified{};
288 bool m_NumberOfLevelsSpecified{};
289};
290} // end namespace itk
291
292#ifndef ITK_MANUAL_INSTANTIATION
293# include "itkMultiResolutionImageRegistrationMethod.hxx"
294#endif
295
296#endif
Decorates any subclass of itkObject with a DataObject API.
SmartPointer< Self > Pointer
Computes similarity between regions of two images.
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.
Base class for multi-resolution image registration methods.
void SetNumberOfLevels(SizeValueType numberOfLevels)
typename MovingImagePyramidType::Pointer MovingImagePyramidPointer
typename MetricType::TransformParametersType ParametersType
typename TransformOutputType::ConstPointer TransformOutputConstPointer
const TransformOutputType * GetOutput() const
void PrintSelf(std::ostream &os, Indent indent) const override
~MultiResolutionImageRegistrationMethod() override=default
void SetSchedules(const ScheduleType &fixedImagePyramidSchedule, const ScheduleType &movingImagePyramidSchedule)
DataObjectPointer MakeOutput(DataObjectPointerArraySizeType output) override
Make a DataObject of the correct type to used as the specified output.
ModifiedTimeType GetMTime() const override
Framework for creating images in a multi-resolution pyramid.
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
unsigned long SizeValueType
Definition: itkIntTypes.h:86