ITK  6.0.0
Insight Toolkit
itkMultiResolutionPDEDeformableRegistration.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 itkMultiResolutionPDEDeformableRegistration_h
19#define itkMultiResolutionPDEDeformableRegistration_h
20
21#include "itkImage.h"
25#include "itkArray.h"
26
27namespace itk
28{
80template <
81 typename TFixedImage,
82 typename TMovingImage,
83 typename TDisplacementField,
84 typename TRealType = float,
85 typename TFloatImageType = Image<TRealType, TFixedImage::ImageDimension>,
86 typename TRegistrationType = PDEDeformableRegistrationFilter<TFloatImageType, TFloatImageType, TDisplacementField>,
87 typename TDefaultRegistrationType = DemonsRegistrationFilter<TFloatImageType, TFloatImageType, TDisplacementField>>
89 : public ImageToImageFilter<TDisplacementField, TDisplacementField>
90{
91public:
92 ITK_DISALLOW_COPY_AND_MOVE(MultiResolutionPDEDeformableRegistration);
93
99
101 itkNewMacro(Self);
102
104 itkOverrideGetNameOfClassMacro(MultiResolutionPDEDeformableRegistration);
105
107 using FixedImageType = TFixedImage;
110
112 using MovingImageType = TMovingImage;
115
117 using DisplacementFieldType = TDisplacementField;
119
121 static constexpr unsigned int ImageDimension = FixedImageType::ImageDimension;
122
124 using FloatImageType = TFloatImageType;
125
127 using RegistrationType = TRegistrationType;
129
131 using DefaultRegistrationType = TDefaultRegistrationType;
132
136
140
144
146
148 virtual void
150
152 const FixedImageType *
154
156 virtual void
158
160 const MovingImageType *
162
165 virtual void
167 {
168 this->m_InitialDisplacementField = ptr;
169 }
170
174 virtual void
176 {
177 this->SetInput(ptr);
178 }
179
181 const DisplacementFieldType *
183 {
184 return this->GetOutput();
185 }
186
193 std::vector<SmartPointer<DataObject>>::size_type
195
197 itkSetObjectMacro(RegistrationFilter, RegistrationType);
198 itkGetModifiableObjectMacro(RegistrationFilter, RegistrationType);
202 itkSetObjectMacro(FixedImagePyramid, FixedImagePyramidType);
203 itkGetModifiableObjectMacro(FixedImagePyramid, FixedImagePyramidType);
207 itkSetObjectMacro(MovingImagePyramid, MovingImagePyramidType);
208 itkGetModifiableObjectMacro(MovingImagePyramid, MovingImagePyramidType);
212 virtual void
213 SetNumberOfLevels(unsigned int num);
214
216 itkGetConstReferenceMacro(NumberOfLevels, unsigned int);
217
219 itkGetConstReferenceMacro(CurrentLevel, unsigned int);
220
222 itkSetObjectMacro(FieldExpander, FieldExpanderType);
223 itkGetModifiableObjectMacro(FieldExpander, FieldExpanderType);
227 itkSetMacro(NumberOfIterations, NumberOfIterationsType);
228 itkSetVectorMacro(NumberOfIterations, unsigned int, m_NumberOfLevels);
232 itkGetConstReferenceMacro(NumberOfIterations, NumberOfIterationsType);
233
235 virtual void
237
238protected:
241
242 void
243 PrintSelf(std::ostream & os, Indent indent) const override;
244
256 void
257 GenerateData() override;
258
262 void
264
271 void
273
277 void
279
282 virtual bool
284
290 void
291 VerifyInputInformation() const override
292 {}
293
294private:
295 RegistrationPointer m_RegistrationFilter{};
296 FixedImagePyramidPointer m_FixedImagePyramid{};
297 MovingImagePyramidPointer m_MovingImagePyramid{};
298 FieldExpanderPointer m_FieldExpander{};
299 DisplacementFieldPointer m_InitialDisplacementField{};
300
301 unsigned int m_NumberOfLevels{};
302 unsigned int m_CurrentLevel{};
303 NumberOfIterationsType m_NumberOfIterations{};
304
306 bool m_StopRegistrationFlag{};
307};
308} // end namespace itk
309
310#ifndef ITK_MANUAL_INSTANTIATION
311# include "itkMultiResolutionPDEDeformableRegistration.hxx"
312#endif
313
314#endif
Base class for all data objects in ITK.
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Light weight base class for most itk classes.
Framework for performing multi-resolution PDE deformable registration.
virtual void SetMovingImage(const MovingImageType *ptr)
void EnlargeOutputRequestedRegion(DataObject *ptr) override
void PrintSelf(std::ostream &os, Indent indent) const override
std::vector< SmartPointer< DataObject > >::size_type GetNumberOfValidRequiredInputs() const override
const MovingImageType * GetMovingImage() const
~MultiResolutionPDEDeformableRegistration() override=default
virtual void SetNumberOfLevels(unsigned int num)
virtual void SetArbitraryInitialDisplacementField(DisplacementFieldType *ptr)
const FixedImageType * GetFixedImage() const
virtual void SetFixedImage(const FixedImageType *ptr)
Framework for creating images in a multi-resolution pyramid.
Resample an image via a coordinate transform.
SmartPointer< const Self > ConstPointer
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....