ITK  4.8.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 >
82  public ImageToImageFilter< TDisplacementField, TDisplacementField >
83 {
84 public:
85 
92 
94  itkNewMacro(Self);
95 
99 
101  typedef TFixedImage FixedImageType;
102  typedef typename FixedImageType::Pointer FixedImagePointer;
103  typedef typename FixedImageType::ConstPointer FixedImageConstPointer;
104 
106  typedef TMovingImage MovingImageType;
107  typedef typename MovingImageType::Pointer MovingImagePointer;
108  typedef typename MovingImageType::ConstPointer MovingImageConstPointer;
109 
111  typedef TDisplacementField DisplacementFieldType;
112  typedef typename DisplacementFieldType::Pointer DisplacementFieldPointer;
113 #ifdef ITKV3_COMPATIBILITY
114  typedef TDisplacementField DeformationFieldType;
115  typedef typename DeformationFieldType::Pointer DeformationFieldPointer;
116 #endif
117 
119  itkStaticConstMacro(ImageDimension, unsigned int, FixedImageType::ImageDimension);
120 
123 
127 
130 
134 
138 
142 
144 
146  virtual void SetFixedImage(const FixedImageType *ptr);
147 
149  const FixedImageType * GetFixedImage() const;
150 
152  virtual void SetMovingImage(const MovingImageType *ptr);
153 
155  const MovingImageType * GetMovingImage() const;
156 
160  {
161  this->m_InitialDisplacementField = ptr;
162  }
163 
168  {
169  this->SetInput(ptr);
170  }
171 
174  { return this->GetOutput(); }
175 
176 #ifdef ITKV3_COMPATIBILITY
177  virtual void SetInitialDeformationField(DisplacementFieldType *ptr)
178  {
179  this->SetInitialDisplacementField(ptr);
180  }
181 
182  virtual void SetArbitraryInitialDeformationField(DisplacementFieldType *ptr)
183  {
185  }
186 
188  const DeformationFieldType * GetDeformationField(void)
189  {
190  return itkDynamicCastInDebugMode<DeformationFieldType *> (this->GetDisplacementField());
191  }
192 #endif
193 
200  virtual std::vector< SmartPointer< DataObject > >::size_type GetNumberOfValidRequiredInputs() const ITK_OVERRIDE;
201 
203  itkSetObjectMacro(RegistrationFilter, RegistrationType);
204  itkGetModifiableObjectMacro(RegistrationFilter, RegistrationType);
206 
208  itkSetObjectMacro(FixedImagePyramid, FixedImagePyramidType);
209  itkGetModifiableObjectMacro(FixedImagePyramid, FixedImagePyramidType);
211 
213  itkSetObjectMacro(MovingImagePyramid, MovingImagePyramidType);
214  itkGetModifiableObjectMacro(MovingImagePyramid, MovingImagePyramidType);
216 
218  virtual void SetNumberOfLevels(unsigned int num);
219 
221  itkGetConstReferenceMacro(NumberOfLevels, unsigned int);
222 
224  itkGetConstReferenceMacro(CurrentLevel, unsigned int);
225 
227  itkSetObjectMacro(FieldExpander, FieldExpanderType);
228  itkGetModifiableObjectMacro(FieldExpander, FieldExpanderType);
230 
232  itkSetMacro(NumberOfIterations, NumberOfIterationsType);
233  itkSetVectorMacro(NumberOfIterations, unsigned int, m_NumberOfLevels);
235 
237  itkGetConstReferenceMacro(NumberOfIterations, NumberOfIterationsType);
238 
240  virtual void StopRegistration();
241 
242 protected:
244  // ~MultiResolutionPDEDeformableRegistration() {} default implementation ok
245 
246  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
247 
250  virtual void GenerateData() ITK_OVERRIDE;
251 
255  virtual void GenerateInputRequestedRegion() ITK_OVERRIDE;
256 
263  virtual void GenerateOutputInformation() ITK_OVERRIDE;
264 
268  virtual void EnlargeOutputRequestedRegion(DataObject *ptr) ITK_OVERRIDE;
269 
272  virtual bool Halt();
273 
279  virtual void VerifyInputInformation() ITK_OVERRIDE {}
280 
281 private:
282  MultiResolutionPDEDeformableRegistration(const Self &); //purposely not
283  // implemented
284  void operator=(const Self &); //purposely not
285 
286  // implemented
287 
293 
294  unsigned int m_NumberOfLevels;
295  unsigned int m_CurrentLevel;
297 
300 };
301 } // end namespace itk
302 
303 #ifndef ITK_MANUAL_INSTANTIATION
304 #include "itkMultiResolutionPDEDeformableRegistration.hxx"
305 #endif
306 
307 #endif
Resample an image via a coordinate transform.
Deformably register two images using the demons algorithm.
virtual std::vector< SmartPointer< DataObject > >::size_type GetNumberOfValidRequiredInputs() const override
VectorResampleImageFilter< DisplacementFieldType, DisplacementFieldType > FieldExpanderType
Light weight base class for most itk classes.
Framework for performing multi-resolution PDE deformable registration.
virtual void SetMovingImage(const MovingImageType *ptr)
void PrintSelf(std::ostream &os, Indent indent) const override
virtual void SetArbitraryInitialDisplacementField(DisplacementFieldType *ptr)
DemonsRegistrationFilter< FloatImageType, FloatImageType, DisplacementFieldType > DefaultRegistrationType
MultiResolutionPyramidImageFilter< MovingImageType, FloatImageType > MovingImagePyramidType
virtual void SetFixedImage(const FixedImageType *ptr)
const MovingImageType * GetMovingImage() const
Framework for creating images in a multi-resolution pyramid.
friend class DataObject
ImageToImageFilter< TDisplacementField, TDisplacementField > Superclass
virtual void SetInput(const InputImageType *image)
virtual void GenerateOutputInformation() override
const FixedImageType * GetFixedImage() const
MultiResolutionPyramidImageFilter< FixedImageType, FloatImageType > FixedImagePyramidType
virtual void GenerateInputRequestedRegion() override
virtual void SetNumberOfLevels(unsigned int num)
Image< TRealType, itkGetStaticConstMacro(ImageDimension) > FloatImageType
Deformably register two images using a PDE algorithm.
virtual void EnlargeOutputRequestedRegion(DataObject *ptr) override
Base class for filters that take an image as input and produce an image as output.
PDEDeformableRegistrationFilter< FloatImageType, FloatImageType, DisplacementFieldType > RegistrationType
Templated n-dimensional image class.
Definition: itkImage.h:75