ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkVariationalRegistrationMultiResolutionFilter.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 itkVariationalRegistrationMultiResolutionFilter_h
19 #define itkVariationalRegistrationMultiResolutionFilter_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  ITK_DISALLOW_COPY_AND_ASSIGN(VariationalRegistrationMultiResolutionFilter);
86 
92 
94  itkNewMacro(Self);
95 
98 
100  using FixedImageType = TFixedImage;
101  using FixedImagePointer = typename FixedImageType::Pointer;
102  using FixedImageConstPointer = typename FixedImageType::ConstPointer;
103 
105  using MovingImageType = TMovingImage;
106  using MovingImagePointer = typename MovingImageType::Pointer;
107  using MovingImageConstPointer = typename MovingImageType::ConstPointer;
108 
110  using DisplacementFieldType = TDisplacementField;
111  using DisplacementFieldPointer = typename DisplacementFieldType::Pointer;
112 
114  static constexpr unsigned int ImageDimension = FixedImageType::ImageDimension;
115 
117  using MaskImagePixelType = unsigned char;
121 
124 
128 
131 
135 
139 
143 
147 
150 
152  virtual void SetFixedImage( const FixedImageType * ptr );
153 
155  const FixedImageType * GetFixedImage(void) const;
156 
158  virtual void SetMovingImage( const MovingImageType * ptr );
159 
161  const MovingImageType * GetMovingImage(void) const;
162 
164  void SetMaskImage( const MaskImageType * ptr );
165 
167  const MaskImageType * GetMaskImage(void) const;
168 
174  { this->SetInput( ptr ); }
175 
179  { return this->GetInput(); }
180 
184  { return this->GetOutput(); }
185 
187  itkGetModifiableObjectMacro( DisplacementField, DisplacementFieldType );
188 
194  std::vector< SmartPointer<DataObject> >::size_type GetNumberOfValidRequiredInputs() const override;
195 
197  itkSetObjectMacro( RegistrationFilter, RegistrationType );
198 
200  itkGetConstObjectMacro( RegistrationFilter, RegistrationType );
201 
203  itkSetObjectMacro( FixedImagePyramid, FixedImagePyramidType );
204 
206  itkGetConstObjectMacro( FixedImagePyramid, FixedImagePyramidType );
207 
209  itkSetObjectMacro( MovingImagePyramid, MovingImagePyramidType );
210 
212  itkGetConstObjectMacro( MovingImagePyramid, MovingImagePyramidType );
213 
215  itkSetObjectMacro( MaskImagePyramid, MaskImagePyramidType );
216 
218  itkGetConstObjectMacro( MaskImagePyramid, MaskImagePyramidType );
219 
221  virtual void SetNumberOfLevels( unsigned int num );
222 
224  itkGetConstReferenceMacro( NumberOfLevels, unsigned int );
225 
227  itkGetConstReferenceMacro( ElapsedLevels, unsigned int );
228 
230  itkSetMacro( NumberOfIterations, NumberOfIterationsType );
231 
234  itkSetVectorMacro( NumberOfIterations, unsigned int, m_NumberOfLevels );
235 
237  itkGetConstReferenceMacro( NumberOfIterations, NumberOfIterationsType );
238 
240  itkSetObjectMacro( FieldExpander, FieldExpanderType );
241 
243  itkGetConstObjectMacro( FieldExpander, FieldExpanderType );
244 
246  virtual void StopRegistration();
247 
248 protected:
251 
253  void PrintSelf( std::ostream& os, Indent indent ) const override;
254 
257  void GenerateData() override;
258 
262  void GenerateInputRequestedRegion() override;
263 
270  void GenerateOutputInformation() override;
271 
275  void EnlargeOutputRequestedRegion( DataObject *ptr ) override;
276 
279  virtual bool Halt();
280 
281 private:
288 
289  unsigned int m_NumberOfLevels;
290  unsigned int m_ElapsedLevels;
292 
295 };
296 
297 } // end namespace itk
298 
299 #ifndef ITK_MANUAL_INSTANTIATION
300 #include "itkVariationalRegistrationMultiResolutionFilter.hxx"
301 #endif
302 
303 
304 #endif
Resample an image via a coordinate transform.
Light weight base class for most itk classes.
Flexible framework for deformable registration of two images using PDE-based variational registration...
Framework for creating images in a multi-resolution pyramid.
Framework for performing multi-resolution variational registration.
MultiResolutionPyramidImageFilter< FloatImageType, FloatImageType > MaskImagePyramidType
const FixedImageType * GetFixedImage(void) const
virtual void SetInput(const InputImageType *image)
MultiResolutionPyramidImageFilter< MovingImageType, MovingImageType > MovingImagePyramidType
std::vector< SmartPointer< DataObject > >::size_type GetNumberOfValidRequiredInputs() const override
VariationalRegistrationFilter< FixedImageType, MovingImageType, DisplacementFieldType > RegistrationType
virtual void SetFixedImage(const FixedImageType *ptr)
MultiResolutionPyramidImageFilter< FixedImageType, FixedImageType > FixedImagePyramidType
const MaskImageType * GetMaskImage(void) const
void PrintSelf(std::ostream &os, Indent indent) const override
VectorResampleImageFilter< DisplacementFieldType, DisplacementFieldType > FieldExpanderType
SmartPointer< const Self > ConstPointer
Definition: itkImage.h:84
virtual void SetNumberOfLevels(unsigned int num)
virtual void SetMovingImage(const MovingImageType *ptr)
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
SmartPointer< Self > Pointer
Definition: itkImage.h:83
const MovingImageType * GetMovingImage(void) const
void SetMaskImage(const MaskImageType *ptr)
void EnlargeOutputRequestedRegion(DataObject *ptr) override
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75