ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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>,
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;
108 using FixedImagePointer = typename FixedImageType::Pointer;
109 using FixedImageConstPointer = typename FixedImageType::ConstPointer;
110
112 using MovingImageType = TMovingImage;
113 using MovingImagePointer = typename MovingImageType::Pointer;
114 using MovingImageConstPointer = typename MovingImageType::ConstPointer;
115
117 using DisplacementFieldType = TDisplacementField;
118 using DisplacementFieldPointer = typename DisplacementFieldType::Pointer;
119
121 static constexpr unsigned int ImageDimension = FixedImageType::ImageDimension;
122
124 using FloatImageType = TFloatImageType;
125
127 using RegistrationType = TRegistrationType;
128 using RegistrationPointer = typename RegistrationType::Pointer;
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
170
174 virtual void
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);
200
202 itkSetObjectMacro(FixedImagePyramid, FixedImagePyramidType);
203 itkGetModifiableObjectMacro(FixedImagePyramid, FixedImagePyramidType);
205
207 itkSetObjectMacro(MovingImagePyramid, MovingImagePyramidType);
208 itkGetModifiableObjectMacro(MovingImagePyramid, MovingImagePyramidType);
210
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);
225
227 virtual void
229 itkSetVectorMacro(NumberOfIterations, unsigned int, m_NumberOfLevels);
231
233 itkGetConstReferenceMacro(NumberOfIterations, NumberOfIterationsType);
234
236 virtual void
238
239protected:
242
243 void
244 PrintSelf(std::ostream & os, Indent indent) const override;
245
257 void
258 GenerateData() override;
259
263 void
265
272 void
274
278 void
280
283 virtual bool
285
291 void
292 VerifyInputInformation() const override
293 {}
294
295private:
301
302 unsigned int m_NumberOfLevels{};
303 unsigned int m_CurrentLevel{};
305
308};
309} // end namespace itk
310
311#ifndef ITK_MANUAL_INSTANTIATION
312# include "itkMultiResolutionPDEDeformableRegistration.hxx"
313#endif
314
315#endif
Array class with size defined at construction time.
Definition itkArray.h:48
Base class for all data objects in ITK.
Deformably register two images using the demons algorithm.
OutputImageType * GetOutput()
virtual void SetInput(const InputImageType *input)
Templated n-dimensional image class.
Definition itkImage.h:89
Control indentation during Print() invocation.
Definition itkIndent.h:50
MultiResolutionPyramidImageFilter< FixedImageType, FloatImageType > FixedImagePyramidType
virtual void SetMovingImage(const MovingImageType *ptr)
void EnlargeOutputRequestedRegion(DataObject *ptr) override
void PrintSelf(std::ostream &os, Indent indent) const override
MultiResolutionPyramidImageFilter< MovingImageType, FloatImageType > MovingImagePyramidType
std::vector< SmartPointer< DataObject > >::size_type GetNumberOfValidRequiredInputs() const override
const MovingImageType * GetMovingImage() const
~MultiResolutionPDEDeformableRegistration() override=default
virtual void SetNumberOfLevels(unsigned int num)
ResampleImageFilter< DisplacementFieldType, DisplacementFieldType > FieldExpanderType
virtual void SetArbitraryInitialDisplacementField(DisplacementFieldType *ptr)
virtual void SetNumberOfIterations(NumberOfIterationsType numberOfIterations)
ImageToImageFilter< TDisplacementField, TDisplacementField > Superclass
const FixedImageType * GetFixedImage() const
virtual void SetFixedImage(const FixedImageType *ptr)
Framework for creating images in a multi-resolution pyramid.
Deformably register two images using a PDE algorithm.
Resample an image via a coordinate transform.
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....