ITK  4.8.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< class TFixedImage, class TMovingImage, class TDisplacementField, class TRealType = float>
82  public ImageToImageFilter< TDisplacementField, TDisplacementField >
83 {
84 public:
91 
93  itkNewMacro(Self);
94 
97 
99  typedef TFixedImage FixedImageType;
100  typedef typename FixedImageType::Pointer FixedImagePointer;
101  typedef typename FixedImageType::ConstPointer FixedImageConstPointer;
102 
104  typedef TMovingImage MovingImageType;
105  typedef typename MovingImageType::Pointer MovingImagePointer;
106  typedef typename MovingImageType::ConstPointer MovingImageConstPointer;
107 
109  typedef TDisplacementField DisplacementFieldType;
110  typedef typename DisplacementFieldType::Pointer DisplacementFieldPointer;
111 
113  itkStaticConstMacro(ImageDimension, unsigned int, FixedImageType::ImageDimension);
114 
116  typedef unsigned char MaskImagePixelType;
120 
123 
128 
132 
137 
142 
147 
152 
155 
157  virtual void SetFixedImage( const FixedImageType * ptr );
158 
160  const FixedImageType * GetFixedImage(void) const;
161 
163  virtual void SetMovingImage( const MovingImageType * ptr );
164 
166  const MovingImageType * GetMovingImage(void) const;
167 
169  void SetMaskImage( const MaskImageType * ptr );
170 
172  const MaskImageType * GetMaskImage(void) const;
173 
179  { this->SetInput( ptr ); }
180 
184  { return this->GetInput(); }
185 
189  { return this->GetOutput(); }
190 
192  itkGetObjectMacro( DisplacementField, DisplacementFieldType );
193 
199  virtual std::vector< SmartPointer<DataObject> >::size_type GetNumberOfValidRequiredInputs() const;
200 
202  itkSetObjectMacro( RegistrationFilter, RegistrationType );
203 
205  itkGetObjectMacro( RegistrationFilter, RegistrationType );
206 
208  itkSetObjectMacro( FixedImagePyramid, FixedImagePyramidType );
209 
211  itkGetObjectMacro( FixedImagePyramid, FixedImagePyramidType );
212 
214  itkSetObjectMacro( MovingImagePyramid, MovingImagePyramidType );
215 
217  itkGetObjectMacro( MovingImagePyramid, MovingImagePyramidType );
218 
220  itkSetObjectMacro( MaskImagePyramid, MaskImagePyramidType );
221 
223  itkGetObjectMacro( MaskImagePyramid, MaskImagePyramidType );
224 
226  virtual void SetNumberOfLevels( unsigned int num );
227 
229  itkGetConstReferenceMacro( NumberOfLevels, unsigned int );
230 
232  itkGetConstReferenceMacro( ElapsedLevels, unsigned int );
233 
235  itkSetMacro( NumberOfIterations, NumberOfIterationsType );
236 
239  itkSetVectorMacro( NumberOfIterations, unsigned int, m_NumberOfLevels );
240 
242  itkGetConstReferenceMacro( NumberOfIterations, NumberOfIterationsType );
243 
245  itkSetObjectMacro( FieldExpander, FieldExpanderType );
246 
248  itkGetObjectMacro( FieldExpander, FieldExpanderType );
249 
251  virtual void StopRegistration();
252 
253 protected:
256 
258  void PrintSelf( std::ostream& os, Indent indent ) const;
259 
262  virtual void GenerateData();
263 
267  virtual void GenerateInputRequestedRegion();
268 
275  virtual void GenerateOutputInformation();
276 
280  virtual void EnlargeOutputRequestedRegion( DataObject *ptr );
281 
284  virtual bool Halt();
285 
286 private:
287  VariationalRegistrationMultiResolutionFilter(const Self&); //purposely not implemented
288  void operator=( const Self& ); //purposely not implemented
289 
296 
297  unsigned int m_NumberOfLevels;
298  unsigned int m_ElapsedLevels;
300 
303 };
304 
305 } // end namespace itk
306 
307 #ifndef ITK_MANUAL_INSTANTIATION
308 #include "itkVariationalRegistrationMultiResolutionFilter.hxx"
309 #endif
310 
311 
312 #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...
MultiResolutionPyramidImageFilter< FixedImageType, FixedImageType > FixedImagePyramidType
MultiResolutionPyramidImageFilter< MovingImageType, MovingImageType > MovingImagePyramidType
Framework for creating images in a multi-resolution pyramid.
Image< TRealType, itkGetStaticConstMacro(ImageDimension) > FloatImageType
Framework for performing multi-resolution variational registration.
VariationalRegistrationFilter< FixedImageType, MovingImageType, DisplacementFieldType > RegistrationType
VariationalRegistrationFilter< FixedImageType, MovingImageType, DisplacementFieldType > DefaultRegistrationType
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
ImageToImageFilter< TDisplacementField, TDisplacementField > Superclass
MultiResolutionPyramidImageFilter< FloatImageType, FloatImageType > MaskImagePyramidType
VectorResampleImageFilter< DisplacementFieldType, DisplacementFieldType > FieldExpanderType
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75