ITK
5.0.0
Insight Segmentation and Registration Toolkit
|
#include <itkN4BiasFieldCorrectionImageFilter.h>
Implementation of the N4 bias field correction algorithm.
The nonparametric nonuniform intensity normalization (N3) algorithm, as introduced by Sled et al. in 1998 is a method for correcting nonuniformity associated with MR images. The algorithm assumes a simple parametric model (Gaussian) for the bias field and does not require tissue class segmentation. In addition, there are only a couple of parameters to tune with the default values performing quite well. N3 has been publicly available as a set of perl scripts (http://www.bic.mni.mcgill.ca/ServicesSoftwareAdvancedImageProcessingTools/HomePage)
The N4 algorithm, encapsulated with this class, is a variation of the original N3 algorithm with the additional benefits of an improved B-spline fitting routine which allows for multiple resolutions to be used during the correction process. We also modify the iterative update component of algorithm such that the residual bias field is continually updated
Notes for the user:
The basic algorithm iterates between sharpening the intensity histogram of the corrected input image and spatially smoothing those results with a B-spline scalar field estimate of the bias field.
Contributed by Nicholas J. Tustison, James C. Gee in the Insight Journal paper: https://hdl.handle.net/10380/3053
J.G. Sled, A.P. Zijdenbos and A.C. Evans. "A Nonparametric Method for Automatic Correction of Intensity Nonuniformity in Data" IEEE Transactions on Medical Imaging, Vol 17, No 1. Feb 1998.
N.J. Tustison, B.B. Avants, P.A. Cook, Y. Zheng, A. Egan, P.A. Yushkevich, and J.C. Gee. "N4ITK: Improved N3 Bias Correction" IEEE Transactions on Medical Imaging, 29(6):1310-1320, June 2010.
Definition at line 95 of file itkN4BiasFieldCorrectionImageFilter.h.
Static Public Member Functions | |
static Pointer | New () |
Static Public Member Functions inherited from itk::Object | |
static bool | GetGlobalWarningDisplay () |
static void | GlobalWarningDisplayOff () |
static void | GlobalWarningDisplayOn () |
static Pointer | New () |
static void | SetGlobalWarningDisplay (bool flag) |
Static Public Member Functions inherited from itk::LightObject | |
static void | BreakOnError () |
static Pointer | New () |
Static Public Attributes | |
static constexpr unsigned int | ImageDimension = TInputImage::ImageDimension |
Static Public Attributes inherited from itk::ImageToImageFilter< TInputImage, TOutputImage > | |
static constexpr unsigned int | InputImageDimension = TInputImage::ImageDimension |
static constexpr unsigned int | OutputImageDimension = TOutputImage::ImageDimension |
Static Public Attributes inherited from itk::ImageSource< TOutputImage > | |
static constexpr unsigned int | OutputImageDimension = TOutputImage::ImageDimension |
Private Member Functions | |
RealType | CalculateConvergenceMeasurement (const RealImageType *, const RealImageType *) const |
RealImagePointer | ReconstructBiasField (BiasFieldControlPointLatticeType *) |
RealImagePointer | SharpenImage (const RealImageType *) const |
RealImagePointer | UpdateBiasFieldEstimate (RealImageType *, std::vcl_size_t) |
Private Attributes | |
RealType | m_BiasFieldFullWidthAtHalfMaximum { static_cast< RealType >( 0.15 ) } |
RealType | m_ConvergenceThreshold { static_cast< RealType >( 0.001 ) } |
RealType | m_CurrentConvergenceMeasurement |
unsigned int | m_CurrentLevel { 0 } |
unsigned int | m_ElapsedIterations { 0 } |
BiasFieldControlPointLatticeType::Pointer | m_LogBiasFieldControlPointLattice |
MaskPixelType | m_MaskLabel |
VariableSizeArrayType | m_MaximumNumberOfIterations |
ArrayType | m_NumberOfControlPoints |
ArrayType | m_NumberOfFittingLevels |
unsigned int | m_NumberOfHistogramBins { 200 } |
unsigned int | m_SplineOrder { 3 } |
bool | m_UseMaskLabel { false } |
RealType | m_WienerFilterNoise { static_cast< RealType >( 0.01 ) } |
Additional Inherited Members | |
Protected Types inherited from itk::ImageToImageFilter< TInputImage, TOutputImage > | |
using | InputToOutputRegionCopierType = ImageToImageFilterDetail::ImageRegionCopier< Self::OutputImageDimension, Self::InputImageDimension > |
using | OutputToInputRegionCopierType = ImageToImageFilterDetail::ImageRegionCopier< Self::InputImageDimension, Self::OutputImageDimension > |
Static Protected Member Functions inherited from itk::ImageSource< TOutputImage > | |
static const ImageRegionSplitterBase * | GetGlobalDefaultSplitter () |
static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION | ThreaderCallback (void *arg) |
Protected Attributes inherited from itk::ImageSource< TOutputImage > | |
bool | m_DynamicMultiThreading |
Protected Attributes inherited from itk::ProcessObject | |
TimeStamp | m_OutputInformationMTime |
bool | m_Updating |
Protected Attributes inherited from itk::LightObject | |
std::atomic< int > | m_ReferenceCount |
using itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::ArrayType = typename BSplineFilterType::ArrayType |
Definition at line 138 of file itkN4BiasFieldCorrectionImageFilter.h.
using itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::BiasFieldControlPointLatticeType = typename BSplineFilterType::PointDataImageType |
Definition at line 137 of file itkN4BiasFieldCorrectionImageFilter.h.
using itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::BSplineFilterType = BSplineScatteredDataPointSetToImageFilter< PointSetType, ScalarImageType> |
B-sline filter type alias
Definition at line 136 of file itkN4BiasFieldCorrectionImageFilter.h.
using itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::ConstPointer = SmartPointer<const Self> |
Definition at line 105 of file itkN4BiasFieldCorrectionImageFilter.h.
using itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::InputImageType = TInputImage |
Some convenient type alias.
Definition at line 117 of file itkN4BiasFieldCorrectionImageFilter.h.
using itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::MaskImageType = TMaskImage |
Definition at line 119 of file itkN4BiasFieldCorrectionImageFilter.h.
using itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::MaskPixelType = typename MaskImageType::PixelType |
Definition at line 120 of file itkN4BiasFieldCorrectionImageFilter.h.
using itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::OutputImageType = TOutputImage |
Definition at line 118 of file itkN4BiasFieldCorrectionImageFilter.h.
using itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::Pointer = SmartPointer<Self> |
Definition at line 104 of file itkN4BiasFieldCorrectionImageFilter.h.
using itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::PointSetPointer = typename PointSetType::Pointer |
Definition at line 131 of file itkN4BiasFieldCorrectionImageFilter.h.
using itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::PointSetType = PointSet<ScalarType, Self::ImageDimension > |
Definition at line 129 of file itkN4BiasFieldCorrectionImageFilter.h.
using itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::PointType = typename PointSetType::PointType |
Definition at line 132 of file itkN4BiasFieldCorrectionImageFilter.h.
using itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::RealImagePointer = typename RealImageType::Pointer |
Definition at line 124 of file itkN4BiasFieldCorrectionImageFilter.h.
using itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::RealImageType = Image<RealType, ImageDimension> |
Definition at line 123 of file itkN4BiasFieldCorrectionImageFilter.h.
using itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::RealType = float |
Definition at line 122 of file itkN4BiasFieldCorrectionImageFilter.h.
using itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::ScalarImageType = Image<ScalarType, Self::ImageDimension > |
Definition at line 130 of file itkN4BiasFieldCorrectionImageFilter.h.
using itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::ScalarType = Vector<RealType, 1> |
B-spline smoothing filter argument type alias
Definition at line 128 of file itkN4BiasFieldCorrectionImageFilter.h.
using itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::Self = N4BiasFieldCorrectionImageFilter |
Standard class type aliases.
Definition at line 102 of file itkN4BiasFieldCorrectionImageFilter.h.
using itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::Superclass = ImageToImageFilter<TInputImage, TOutputImage> |
Definition at line 103 of file itkN4BiasFieldCorrectionImageFilter.h.
using itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::VariableSizeArrayType = Array<unsigned int> |
Definition at line 125 of file itkN4BiasFieldCorrectionImageFilter.h.
|
protected |
|
overrideprotecteddefault |
|
private |
Convergence is determined by the coefficient of variation of the difference image between the current bias field estimate and the previous estimate.
|
virtual |
Create an object from an instance, potentially deferring to a factory. This method allows you to create an instance of an object that is exactly the same type as the referring object. This is useful in cases where an object has been cast back to a base class.
Reimplemented from itk::Object.
|
overridevirtual |
Ensures that this filter can compute the entire output at once.
Reimplemented from itk::ProcessObject.
|
overrideprotectedvirtual |
A version of GenerateData() specific for image processing filters. This implementation will split the processing across multiple threads. The buffer is allocated by this method. Then the BeforeThreadedGenerateData() method is called (if provided). Then, a series of threads are spawned each calling DynamicThreadedGenerateData(). After all the threads have completed processing, the AfterThreadedGenerateData() method is called (if provided). If an image processing filter cannot be threaded, the filter should provide an implementation of GenerateData(). That implementation is responsible for allocating the output buffer. If a filter can be threaded, it should NOT provide a GenerateData() method but should provide a DynamicThreadedGenerateData() instead.
Reimplemented from itk::ImageSource< TOutputImage >.
|
virtual |
Get the full width at half maximum parameter characterizing the width of the Gaussian deconvolution. Default = 0.15.
|
virtual |
Get confidence image function. If a confidence image is specified, estimation of the bias field weights the contribution of each voxel according the value of the corresponding voxel in the confidence image. For example, when estimating the bias field using brain , one can use a soft segmentation of the white matter as the confidence image instead of using a hard segmentation of the white matter as the mask image (as has been done in the literature) as an alternative strategy to estimating the bias field.
|
virtual |
Get the convergence threshold. Convergence is determined by the coefficient of variation of the difference image between the current bias field estimate and the previous estimate. If this value is less than the specified threshold, the algorithm proceeds to the next fitting level or terminates if it is at the last level.
|
virtual |
Get the current convergence measurement. This is a helper function for reporting observations.
|
virtual |
Get the current fitting level. This is a helper function for reporting observations.
|
virtual |
Get the number of elapsed iterations. This is a helper function for reporting observations.
|
virtual |
Typically, a reduced size image is used as input to the N4 filter using something like itkShrinkImageFilter. Since the output is a corrected version of the input, the user will probably want to apply the bias field correction to the full resolution image. This can be done by using the LogBiasFieldControlPointLattice to reconstruct the bias field at the full image resolution (using the class BSplineControlPointImageFilter).
|
virtual |
Get mask image function. If a binary mask image is specified, only those input image voxels inside the mask image values are used in estimating the bias field.
|
virtual |
Set/Get mask label value. If a binary mask image is specified and if UseMaskValue is true, only those input image voxels corresponding with mask image values equal to MaskLabel are used in estimating the bias field. If a MaskImage is specified and UseMaskLabel is false, all input image voxels corresponding to non-zero voxels in the MaskImage are used in estimating the bias field. Default = 1.
|
virtual |
Get the maximum number of iterations specified at each fitting level. Default = 50.
|
virtual |
Runtime information support.
Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.
|
virtual |
Get the control point grid size defining the B-spline estimate of the scalar bias field. In each dimension, the B-spline mesh size is equal to the number of control points in that dimension minus the spline order. Default = 4 control points in each dimension for a mesh size of 1 in each dimension.
|
virtual |
Get the number of fitting levels. One of the contributions of N4 is the introduction of a multi-scale approach to fitting. This allows one to specify a B-spline mesh size for initial fitting followed by a doubling of the mesh resolution for each subsequent fitting level. Default = 1 level.
|
virtual |
Get number of bins defining the log input intensity histogram. Default = 200.
|
virtual |
Get the spline order defining the bias field estimate. Default = 3.
|
virtual |
Use a mask label for identifying mask functionality. See SetMaskLabel. Defaults to false.
|
virtual |
Get the noise estimate defining the Wiener filter. Default = 0.01.
|
static |
Standard New method.
|
overrideprotectedvirtual |
Methods invoked by Print() to print information about the object including superclasses. Typically not called by the user (use Print() instead) but used in the hierarchical print process to combine the output of several classes.
Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.
|
private |
Reconstruct bias field given the control point lattice.
|
virtual |
Set the full width at half maximum parameter characterizing the width of the Gaussian deconvolution. Default = 0.15.
|
virtual |
Set confidence image function. If a confidence image is specified, estimation of the bias field weights the contribution of each voxel according the value of the corresponding voxel in the confidence image. For example, when estimating the bias field using brain , one can use a soft segmentation of the white matter as the confidence image instead of using a hard segmentation of the white matter as the mask image (as has been done in the literature) as an alternative strategy to estimating the bias field.
Referenced by itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::SetInput3().
|
virtual |
Set the convergence threshold. Convergence is determined by the coefficient of variation of the difference image between the current bias field estimate and the previous estimate. If this value is less than the specified threshold, the algorithm proceeds to the next fitting level or terminates if it is at the last level.
|
inline |
The image expected for input for bias correction.
Definition at line 146 of file itkN4BiasFieldCorrectionImageFilter.h.
References itk::ImageToImageFilter< TInputImage, TOutputImage >::SetInput().
|
inline |
Set mask image function. If a binary mask image is specified, only those input image voxels inside the mask image values are used in estimating the bias field.
Definition at line 154 of file itkN4BiasFieldCorrectionImageFilter.h.
References itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::SetMaskImage().
|
inline |
Set confidence image function. If a confidence image is specified, estimation of the bias field weights the contribution of each voxel according the value of the corresponding voxel in the confidence image. For example, when estimating the bias field using brain , one can use a soft segmentation of the white matter as the confidence image instead of using a hard segmentation of the white matter as the mask image (as has been done in the literature) as an alternative strategy to estimating the bias field.
Definition at line 195 of file itkN4BiasFieldCorrectionImageFilter.h.
References itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::SetConfidenceImage().
|
virtual |
Set mask image function. If a binary mask image is specified, only those input image voxels inside the mask image values are used in estimating the bias field.
Referenced by itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::SetInput2().
|
virtual |
Set/Get mask label value. If a binary mask image is specified and if UseMaskValue is true, only those input image voxels corresponding with mask image values equal to MaskLabel are used in estimating the bias field. If a MaskImage is specified and UseMaskLabel is false, all input image voxels corresponding to non-zero voxels in the MaskImage are used in estimating the bias field. Default = 1.
|
virtual |
Set the maximum number of iterations specified at each fitting level. Default = 50.
|
virtual |
Set the control point grid size defining the B-spline estimate of the scalar bias field. In each dimension, the B-spline mesh size is equal to the number of control points in that dimension minus the spline order. Default = 4 control points in each dimension for a mesh size of 1 in each dimension.
|
virtual |
Set the number of fitting levels. One of the contributions of N4 is the introduction of a multi-scale approach to fitting. This allows one to specify a B-spline mesh size for initial fitting followed by a doubling of the mesh resolution for each subsequent fitting level. Default = 1 level.
Referenced by itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::SetNumberOfFittingLevels().
|
inline |
Set the number of fitting levels. One of the contributions of N4 is the introduction of a multi-scale approach to fitting. This allows one to specify a B-spline mesh size for initial fitting followed by a doubling of the mesh resolution for each subsequent fitting level. Default = 1 level.
Definition at line 293 of file itkN4BiasFieldCorrectionImageFilter.h.
References itk::N4BiasFieldCorrectionImageFilter< TInputImage, TMaskImage, TOutputImage >::SetNumberOfFittingLevels().
|
virtual |
Set number of bins defining the log input intensity histogram. Default = 200.
|
virtual |
Set the spline order defining the bias field estimate. Default = 3.
|
virtual |
Use a mask label for identifying mask functionality. See SetMaskLabel. Defaults to false.
|
virtual |
Set the noise estimate defining the Wiener filter. Default = 0.01.
|
private |
Sharpen the intensity histogram of the current estimate of the corrected image and map those results to a new estimate of the unsmoothed corrected image.
|
private |
Given the unsmoothed estimate of the bias field, this function smooths the estimate and adds the resulting control point values to the total bias field estimate.
|
virtual |
Use a mask label for identifying mask functionality. See SetMaskLabel. Defaults to false.
|
virtual |
Use a mask label for identifying mask functionality. See SetMaskLabel. Defaults to false.
|
static |
ImageDimension constants
Definition at line 114 of file itkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 416 of file itkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 422 of file itkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 423 of file itkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 424 of file itkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 421 of file itkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 429 of file itkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 409 of file itkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 420 of file itkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 432 of file itkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 433 of file itkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 414 of file itkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 431 of file itkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 410 of file itkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 415 of file itkN4BiasFieldCorrectionImageFilter.h.