ITK  4.6.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 
90 template<typename TInputImage, typename TMaskImage =
91  Image<unsigned char, TInputImage::ImageDimension>,
92  class TOutputImage = TInputImage>
94  public ImageToImageFilter<TInputImage, TOutputImage>
95 {
96 public:
102 
105 
107  itkNewMacro( Self );
108 
110  itkStaticConstMacro( ImageDimension, unsigned int,
111  TInputImage::ImageDimension );
112 
114  typedef TInputImage InputImageType;
115  typedef TOutputImage OutputImageType;
116  typedef TMaskImage MaskImageType;
117  typedef typename MaskImageType::PixelType MaskPixelType;
118 
119  typedef float RealType;
123 
131 
137 
141  void SetInput1( const InputImageType *image ) { this->SetInput( image ); }
142 
143 
149  void SetMaskImage( const MaskImageType *mask )
150  {
151  this->SetNthInput( 1, const_cast<MaskImageType *>( mask ) );
152  }
153  void SetInput2( const MaskImageType *mask ) { this->SetMaskImage( mask ); }
155 
162  {
163  return static_cast<const MaskImageType*>( this->ProcessObject::GetInput( 1 ) );
164  }
165 
176  void SetConfidenceImage( const RealImageType *image )
177  {
178  this->SetNthInput( 2, const_cast<RealImageType *>( image ) );
179  }
180  void SetInput3( const RealImageType *image ) { this->SetConfidenceImage( image ); }
182 
194  {
195  return static_cast<const RealImageType*>( this->ProcessObject::GetInput( 2 ) );
196  }
197 
203  itkSetMacro( MaskLabel, MaskPixelType );
204 
210  itkGetConstMacro( MaskLabel, MaskPixelType );
211 
212  // Sharpen histogram parameters: in estimating the bias field, the
213  // first step is to sharpen the intensity histogram by Wiener deconvolution
214  // with a 1-D Gaussian. The following parameters define this operation.
215  // These default values in N4 match the default values in N3.
216 
221  itkSetMacro( NumberOfHistogramBins, unsigned int );
222 
227  itkGetConstMacro( NumberOfHistogramBins, unsigned int );
228 
232  itkSetMacro( WienerFilterNoise, RealType );
233 
237  itkGetConstMacro( WienerFilterNoise, RealType );
238 
243  itkSetMacro( BiasFieldFullWidthAtHalfMaximum, RealType );
244 
249  itkGetConstMacro( BiasFieldFullWidthAtHalfMaximum, RealType );
250 
251  // B-spline parameters governing the fitting routine
252 
256  itkSetMacro( SplineOrder, unsigned int );
257 
261  itkGetConstMacro( SplineOrder, unsigned int );
262 
270  itkSetMacro( NumberOfControlPoints, ArrayType );
271 
279  itkGetConstMacro( NumberOfControlPoints, ArrayType );
280 
287  itkSetMacro( NumberOfFittingLevels, ArrayType );
288 
295  void SetNumberOfFittingLevels( unsigned int n )
296  {
297  ArrayType nlevels;
298 
299  nlevels.Fill( n );
300  this->SetNumberOfFittingLevels( nlevels );
301  }
302 
309  itkGetConstMacro( NumberOfFittingLevels, ArrayType );
310 
315  itkSetMacro( MaximumNumberOfIterations, VariableSizeArrayType );
316 
321  itkGetConstMacro( MaximumNumberOfIterations, VariableSizeArrayType );
322 
330  itkSetMacro( ConvergenceThreshold, RealType );
331 
339  itkGetConstMacro( ConvergenceThreshold, RealType );
340 
350  itkGetConstMacro( LogBiasFieldControlPointLattice,
352 
357  itkGetConstMacro( ElapsedIterations, unsigned int );
358 
363  itkGetConstMacro( CurrentConvergenceMeasurement, RealType );
364 
369  itkGetConstMacro( CurrentLevel, unsigned int );
370 
371 protected:
374  void PrintSelf( std::ostream& os, Indent indent ) const;
375 
376  void GenerateData();
377 
378 private:
379  N4BiasFieldCorrectionImageFilter( const Self& ); //purposely not
380  // implemented
381  void operator=( const Self& ); //purposely not
382  // implemented
383 
384  // N4 algorithm functions: The basic algorithm iterates between sharpening
385  // the intensity histogram of the corrected input image and spatially
386  // smoothing those results with a B-spline scalar field estimate of the
387  // bias field. The former is handled by the function SharpenImage()
388  // whereas the latter is handled by the function UpdateBiasFieldEstimate().
389  // Convergence is determined by the coefficient of variation of the difference
390  // image between the current bias field estimate and the previous estimate.
391 
398 
405 
411 
412 
414 
415  // Parameters for deconvolution with Wiener filter
416 
420 
421  // Convergence parameters
422 
424  unsigned int m_ElapsedIterations;
427  unsigned int m_CurrentLevel;
428 
429  // B-spline fitting parameters
430 
431  typename
433 
434  unsigned int m_SplineOrder;
437 
438 };
439 
440 } // end namespace itk
441 
442 #ifndef ITK_MANUAL_INSTANTIATION
443 #include "itkN4BiasFieldCorrectionImageFilter.hxx"
444 #endif
445 
446 #endif
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)
void PrintSelf(std::ostream &os, Indent indent) const
BSplineFilterType::PointDataImageType 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