ITK  4.12.0
Insight Segmentation and Registration Toolkit
itkN4BiasFieldCorrectionImageFilter.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 itkN4BiasFieldCorrectionImageFilter_h
19 #define itkN4BiasFieldCorrectionImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 
23 #include "itkArray.h"
25 #include "itkPointSet.h"
26 #include "itkVector.h"
27 
28 #include "vnl/vnl_vector.h"
29 
30 namespace itk {
31 
92 template<typename TInputImage, typename TMaskImage =
93  Image<unsigned char, TInputImage::ImageDimension>,
94  class TOutputImage = TInputImage>
96  public ImageToImageFilter<TInputImage, TOutputImage>
97 {
98 public:
104 
107 
109  itkNewMacro( Self );
110 
112  itkStaticConstMacro( ImageDimension, unsigned int,
113  TInputImage::ImageDimension );
114 
116  typedef TInputImage InputImageType;
117  typedef TOutputImage OutputImageType;
118  typedef TMaskImage MaskImageType;
119  typedef typename MaskImageType::PixelType MaskPixelType;
120 
121  typedef float RealType;
125 
133 
139 
143  void SetInput1( const InputImageType *image ) { this->SetInput( image ); }
144 
145 
151  void SetMaskImage( const MaskImageType *mask )
152  {
153  this->SetNthInput( 1, const_cast<MaskImageType *>( mask ) );
154  }
155  void SetInput2( const MaskImageType *mask ) { this->SetMaskImage( mask ); }
157 
164  {
165  return static_cast<const MaskImageType*>( this->ProcessObject::GetInput( 1 ) );
166  }
167 
168 #if ! defined ( ITK_FUTURE_LEGACY_REMOVE )
169 
178  itkSetMacro( MaskLabel, MaskPixelType );
179  itkGetConstMacro( MaskLabel, MaskPixelType );
181 
185  itkSetMacro( UseMaskLabel, bool );
186  itkGetConstMacro( UseMaskLabel, bool );
187  itkBooleanMacro( UseMaskLabel );
188 #endif
189 
190 
201  void SetConfidenceImage( const RealImageType *image )
202  {
203  this->SetNthInput( 2, const_cast<RealImageType *>( image ) );
204  }
205  void SetInput3( const RealImageType *image ) { this->SetConfidenceImage( image ); }
207 
219  {
220  return static_cast<const RealImageType*>( this->ProcessObject::GetInput( 2 ) );
221  }
222 
223  // Sharpen histogram parameters: in estimating the bias field, the
224  // first step is to sharpen the intensity histogram by Wiener deconvolution
225  // with a 1-D Gaussian. The following parameters define this operation.
226  // These default values in N4 match the default values in N3.
227 
232  itkSetMacro( NumberOfHistogramBins, unsigned int );
233 
238  itkGetConstMacro( NumberOfHistogramBins, unsigned int );
239 
243  itkSetMacro( WienerFilterNoise, RealType );
244 
248  itkGetConstMacro( WienerFilterNoise, RealType );
249 
254  itkSetMacro( BiasFieldFullWidthAtHalfMaximum, RealType );
255 
260  itkGetConstMacro( BiasFieldFullWidthAtHalfMaximum, RealType );
261 
262  // B-spline parameters governing the fitting routine
263 
267  itkSetMacro( SplineOrder, unsigned int );
268 
272  itkGetConstMacro( SplineOrder, unsigned int );
273 
281  itkSetMacro( NumberOfControlPoints, ArrayType );
282 
290  itkGetConstMacro( NumberOfControlPoints, ArrayType );
291 
298  itkSetMacro( NumberOfFittingLevels, ArrayType );
299 
306  void SetNumberOfFittingLevels( unsigned int n )
307  {
308  ArrayType nlevels;
309 
310  nlevels.Fill( n );
311  this->SetNumberOfFittingLevels( nlevels );
312  }
313 
320  itkGetConstMacro( NumberOfFittingLevels, ArrayType );
321 
326  itkSetMacro( MaximumNumberOfIterations, VariableSizeArrayType );
327 
332  itkGetConstMacro( MaximumNumberOfIterations, VariableSizeArrayType );
333 
341  itkSetMacro( ConvergenceThreshold, RealType );
342 
350  itkGetConstMacro( ConvergenceThreshold, RealType );
351 
361  itkGetConstMacro( LogBiasFieldControlPointLattice,
363 
368  itkGetConstMacro( ElapsedIterations, unsigned int );
369 
374  itkGetConstMacro( CurrentConvergenceMeasurement, RealType );
375 
380  itkGetConstMacro( CurrentLevel, unsigned int );
381 
382 protected:
385  void PrintSelf( std::ostream& os, Indent indent ) const ITK_OVERRIDE;
386 
387  void GenerateData() ITK_OVERRIDE;
388 
389 private:
390  ITK_DISALLOW_COPY_AND_ASSIGN(N4BiasFieldCorrectionImageFilter);
391 
392  // N4 algorithm functions: The basic algorithm iterates between sharpening
393  // the intensity histogram of the corrected input image and spatially
394  // smoothing those results with a B-spline scalar field estimate of the
395  // bias field. The former is handled by the function SharpenImage()
396  // whereas the latter is handled by the function UpdateBiasFieldEstimate().
397  // Convergence is determined by the coefficient of variation of the difference
398  // image between the current bias field estimate and the previous estimate.
399 
406 
413 
418 
423  RealType CalculateConvergenceMeasurement( const RealImageType *, const RealImageType * ) const;
424 
425 #if ! defined ( ITK_FUTURE_LEGACY_REMOVE )
428 #endif
429 
430  // Parameters for deconvolution with Wiener filter
431 
435 
436  // Convergence parameters
437 
439  unsigned int m_ElapsedIterations;
442  unsigned int m_CurrentLevel;
443 
444  // B-spline fitting parameters
445 
446  typename
448 
449  unsigned int m_SplineOrder;
452 
453 };
454 
455 } // end namespace itk
456 
457 #ifndef ITK_MANUAL_INSTANTIATION
458 #include "itkN4BiasFieldCorrectionImageFilter.hxx"
459 #endif
460 
461 #endif
void PrintSelf(std::ostream &os, Indent indent) const override
PointSet< ScalarType, itkGetStaticConstMacro(ImageDimension)> PointSetType
SmartPointer< Self > Pointer
Definition: itkImage.h:81
Image< ScalarType, itkGetStaticConstMacro(ImageDimension)> ScalarImageType
Base class for all process objects that output image data.
void Fill(const ValueType &)
RealImagePointer SharpenImage(const RealImageType *) const
ImageToImageFilter< TInputImage, TOutputImage > Superclass
RealImagePointer UpdateBiasFieldEstimate(RealImageType *)
virtual void SetInput(const InputImageType *image)
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:62
BSplineScatteredDataPointSetToImageFilter< PointSetType, ScalarImageType > BSplineFilterType
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition: itkPointSet.h:84
virtual void SetNumberOfFittingLevels(ArrayType _arg)
MeshTraits::PointType PointType
Definition: itkPointSet.h:106
Image filter which provides a B-spline output approximation.
DataObject * GetInput(const DataObjectIdentifierType &key)
Return an input.
BSplineFilterType::PointDataImageType BiasFieldControlPointLatticeType
RealImagePointer ReconstructBiasField(BiasFieldControlPointLatticeType *)
Base class for filters that take an image as input and produce an image as output.
RealType CalculateConvergenceMeasurement(const RealImageType *, const RealImageType *) const
Control indentation during Print() invocation.
Definition: itkIndent.h:49
virtual void SetNthInput(DataObjectPointerArraySizeType num, DataObject *input)
BiasFieldControlPointLatticeType::Pointer m_LogBiasFieldControlPointLattice
Implementation of the N4 bias field correction algorithm.
Templated n-dimensional image class.
Definition: itkImage.h:75