ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkMultiResolutionPDEDeformableRegistration.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 itkMultiResolutionPDEDeformableRegistration_h
19 #define itkMultiResolutionPDEDeformableRegistration_h
20 
21 #include "itkImage.h"
25 #include "itkArray.h"
26 
27 namespace itk
28 {
80 template< typename TFixedImage, typename TMovingImage, typename TDisplacementField, typename TRealType = float >
81 class ITK_TEMPLATE_EXPORT MultiResolutionPDEDeformableRegistration:
82  public ImageToImageFilter< TDisplacementField, TDisplacementField >
83 {
84 public:
85  ITK_DISALLOW_COPY_AND_ASSIGN(MultiResolutionPDEDeformableRegistration);
86 
89  using Superclass =
93 
95  itkNewMacro(Self);
96 
100 
102  using FixedImageType = TFixedImage;
103  using FixedImagePointer = typename FixedImageType::Pointer;
104  using FixedImageConstPointer = typename FixedImageType::ConstPointer;
105 
107  using MovingImageType = TMovingImage;
108  using MovingImagePointer = typename MovingImageType::Pointer;
109  using MovingImageConstPointer = typename MovingImageType::ConstPointer;
110 
112  using DisplacementFieldType = TDisplacementField;
113  using DisplacementFieldPointer = typename DisplacementFieldType::Pointer;
114 
116  static constexpr unsigned int ImageDimension = FixedImageType::ImageDimension;
117 
120 
124 
127 
131 
135 
139 
141 
143  virtual void SetFixedImage(const FixedImageType *ptr);
144 
146  const FixedImageType * GetFixedImage() const;
147 
149  virtual void SetMovingImage(const MovingImageType *ptr);
150 
152  const MovingImageType * GetMovingImage() const;
153 
157  {
158  this->m_InitialDisplacementField = ptr;
159  }
160 
165  {
166  this->SetInput(ptr);
167  }
168 
171  { return this->GetOutput(); }
172 
179  std::vector< SmartPointer< DataObject > >::size_type GetNumberOfValidRequiredInputs() const override;
180 
182  itkSetObjectMacro(RegistrationFilter, RegistrationType);
183  itkGetModifiableObjectMacro(RegistrationFilter, RegistrationType);
185 
187  itkSetObjectMacro(FixedImagePyramid, FixedImagePyramidType);
188  itkGetModifiableObjectMacro(FixedImagePyramid, FixedImagePyramidType);
190 
192  itkSetObjectMacro(MovingImagePyramid, MovingImagePyramidType);
193  itkGetModifiableObjectMacro(MovingImagePyramid, MovingImagePyramidType);
195 
197  virtual void SetNumberOfLevels(unsigned int num);
198 
200  itkGetConstReferenceMacro(NumberOfLevels, unsigned int);
201 
203  itkGetConstReferenceMacro(CurrentLevel, unsigned int);
204 
206  itkSetObjectMacro(FieldExpander, FieldExpanderType);
207  itkGetModifiableObjectMacro(FieldExpander, FieldExpanderType);
209 
211  itkSetMacro(NumberOfIterations, NumberOfIterationsType);
212  itkSetVectorMacro(NumberOfIterations, unsigned int, m_NumberOfLevels);
214 
216  itkGetConstReferenceMacro(NumberOfIterations, NumberOfIterationsType);
217 
219  virtual void StopRegistration();
220 
221 protected:
223  ~MultiResolutionPDEDeformableRegistration() override = default;
224 
225  void PrintSelf(std::ostream & os, Indent indent) const override;
226 
229  void GenerateData() override;
230 
234  void GenerateInputRequestedRegion() override;
235 
242  void GenerateOutputInformation() override;
243 
247  void EnlargeOutputRequestedRegion(DataObject *ptr) override;
248 
251  virtual bool Halt();
252 
258  void VerifyInputInformation() ITKv5_CONST override {}
259 
260 private:
266 
267  unsigned int m_NumberOfLevels;
268  unsigned int m_CurrentLevel;
270 
273 };
274 } // end namespace itk
275 
276 #ifndef ITK_MANUAL_INSTANTIATION
277 #include "itkMultiResolutionPDEDeformableRegistration.hxx"
278 #endif
279 
280 #endif
Resample an image via a coordinate transform.
Deformably register two images using the demons algorithm.
Light weight base class for most itk classes.
Framework for performing multi-resolution PDE deformable registration.
virtual void SetArbitraryInitialDisplacementField(DisplacementFieldType *ptr)
Framework for creating images in a multi-resolution pyramid.
Deformably register two images using a PDE algorithm.
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75