Package org.itk.simple
Class N4BiasFieldCorrectionImageFilter
java.lang.Object
org.itk.simple.ProcessObject
org.itk.simple.ImageFilter
org.itk.simple.N4BiasFieldCorrectionImageFilter
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 ( https://www.bic.mni.mcgill.ca/ServicesSoftwareAdvancedImageProcessingT
ools/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:
Since much of the image manipulation is done in the log space of the
intensities, input images with negative and small values (< 1) can
produce poor results.
The original authors recommend performing the bias field correction on
a downsampled version of the original image.
A binary mask or a weighted image can be supplied. If a binary mask is
specified, those voxels in the input image which correspond to the
voxels in the mask image are used to estimate the bias field. If a
UseMaskLabel value is set to false (the default), all non-zero voxels
in the MaskImage will be masked; otherwise only voxels in the
MaskImage that match the MaskLabel will be used. If a confidence image
is specified, the input voxels are weighted in the b-spline fitting
routine according to the confidence voxel values.
The filter returns the corrected image. If the bias field is wanted,
one can reconstruct it using the class
itkBSplineControlPointImageFilter. See the IJ article and the test
file for an example.
The 'Z' parameter in Sled's 1998 paper is the square root of the class
variable 'm_WienerFilterNoise'.
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.
Nicholas J. Tustison
Contributed by Nicholas J. Tustison, James C. Gee in the Insight
Journal paper: https://www.insight-journal.org/browse/publication/640
REFERENCE
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.
See:
itk::simple::N4BiasFieldCorrection for the procedural interface
itk::N4BiasFieldCorrectionImageFilter for the Doxygen on the original ITK class.
C++ includes: sitkN4BiasFieldCorrectionImageFilter.h
-
Field Summary
Fields inherited from class org.itk.simple.ProcessObject
swigCMemOwn
-
Constructor Summary
ConstructorsModifierConstructorDescriptionitk::simple::N4BiasFieldCorrectionImageFilter::N4BiasFieldCorrectionImageFilter() Default Constructor that takes no arguments and initializes default parametersprotected
N4BiasFieldCorrectionImageFilter
(long cPtr, boolean cMemoryOwn) -
Method Summary
Modifier and TypeMethodDescriptionvoid
delete()
virtual itk::simple::N4BiasFieldCorrectionImageFilter::~N4BiasFieldCorrectionImageFilter() DestructorImage itk::simple::N4BiasFieldCorrectionImageFilter::Execute(const Image &image, const Image &maskImage) Execute the filter on the input imageImage itk::simple::N4BiasFieldCorrectionImageFilter::Execute(const Image &image, const Image &maskImage) Execute the filter on the input imageprotected void
finalize()
double
double itk::simple::N4BiasFieldCorrectionImageFilter::GetBiasFieldFullWidthAtHalfMaximum() const Get the full width at half maximum parameter characterizing the width of the Gaussian deconvolution.double
double itk::simple::N4BiasFieldCorrectionImageFilter::GetConvergenceThreshold() const Get the convergence threshold.protected static long
double
double itk::simple::N4BiasFieldCorrectionImageFilter::GetCurrentConvergenceMeasurement() const Get the current convergence measurement.long
uint32_t itk::simple::N4BiasFieldCorrectionImageFilter::GetCurrentLevel() const Get the current fitting level.long
uint32_t itk::simple::N4BiasFieldCorrectionImageFilter::GetElapsedIterations() const Get the number of elapsed iterations.getLogBiasFieldAsImage
(Image referenceImage) Image itk::simple::N4BiasFieldCorrectionImageFilter::GetLogBiasFieldAsImage(Image referenceImage) const The computed log bias field correction.short
uint8_t itk::simple::N4BiasFieldCorrectionImageFilter::GetMaskLabel() const Set/Get mask label value.std::vector<uint32_t> itk::simple::N4BiasFieldCorrectionImageFilter::GetMaximumNumberOfIterations() const Get the maximum number of iterations specified at each fitting level.getName()
std::string itk::simple::N4BiasFieldCorrectionImageFilter::GetName() const Name of this classstd::vector<uint32_t> itk::simple::N4BiasFieldCorrectionImageFilter::GetNumberOfControlPoints() const Get the control point grid size defining the B-spline estimate of the scalar bias field.long
uint32_t itk::simple::N4BiasFieldCorrectionImageFilter::GetNumberOfHistogramBins() const Get number of bins defining the log input intensity histogram.long
uint32_t itk::simple::N4BiasFieldCorrectionImageFilter::GetSplineOrder() const Get the spline order defining the bias field estimate.boolean
bool itk::simple::N4BiasFieldCorrectionImageFilter::GetUseMaskLabel() const Use a mask label for identifying mask functionality.double
double itk::simple::N4BiasFieldCorrectionImageFilter::GetWienerFilterNoise() const Get the noise estimate defining the Wiener filter.void
setBiasFieldFullWidthAtHalfMaximum
(double BiasFieldFullWidthAtHalfMaximum) Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetBiasFieldFullWidthAtHalfMaximum(double BiasFieldFullWidthAtHalfMaximum) Set the full width at half maximum parameter characterizing the width of the Gaussian deconvolution.void
setConvergenceThreshold
(double ConvergenceThreshold) Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetConvergenceThreshold(double ConvergenceThreshold) Set the convergence threshold.void
setMaskLabel
(short MaskLabel) Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetMaskLabel(uint8_t MaskLabel) Set/Get mask label value.void
setMaximumNumberOfIterations
(VectorUInt32 MaximumNumberOfIterations) Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetMaximumNumberOfIterations(std::vector< uint32_t > MaximumNumberOfIterations) Set the maximum number of iterations specified at each fitting level.void
setNumberOfControlPoints
(long value) Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetNumberOfControlPoints(uint32_t value) Set the values of the NumberOfControlPoints vector all to valuevoid
setNumberOfControlPoints
(VectorUInt32 NumberOfControlPoints) Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetNumberOfControlPoints(uint32_t value) Set the values of the NumberOfControlPoints vector all to valuevoid
setNumberOfHistogramBins
(long NumberOfHistogramBins) Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetNumberOfHistogramBins(uint32_t NumberOfHistogramBins) Set number of bins defining the log input intensity histogram.void
setSplineOrder
(long SplineOrder) Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetSplineOrder(uint32_t SplineOrder) Set the spline order defining the bias field estimate.void
setUseMaskLabel
(boolean UseMaskLabel) Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetUseMaskLabel(bool UseMaskLabel) Use a mask label for identifying mask functionality.void
setWienerFilterNoise
(double WienerFilterNoise) Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetWienerFilterNoise(double WienerFilterNoise) Set the noise estimate defining the Wiener filter.protected static long
toString()
std::string itk::simple::N4BiasFieldCorrectionImageFilter::ToString() const Print ourselves outvoid
Self& itk::simple::N4BiasFieldCorrectionImageFilter::UseMaskLabelOff()void
Self& itk::simple::N4BiasFieldCorrectionImageFilter::UseMaskLabelOn() Set the value of UseMaskLabel to true or false respectfully.Methods inherited from class org.itk.simple.ImageFilter
getCPtr, swigRelease
Methods inherited from class org.itk.simple.ProcessObject
abort, addCommand, debugOff, debugOn, getCPtr, getDebug, getGlobalDefaultCoordinateTolerance, getGlobalDefaultDebug, getGlobalDefaultDirectionTolerance, getGlobalDefaultNumberOfThreads, getGlobalDefaultThreader, getGlobalWarningDisplay, getNumberOfThreads, getNumberOfWorkUnits, getProgress, globalDefaultDebugOff, globalDefaultDebugOn, globalWarningDisplayOff, globalWarningDisplayOn, hasCommand, removeAllCommands, setDebug, setGlobalDefaultCoordinateTolerance, setGlobalDefaultDebug, setGlobalDefaultDirectionTolerance, setGlobalDefaultNumberOfThreads, setGlobalDefaultThreader, setGlobalWarningDisplay, setNumberOfThreads, setNumberOfWorkUnits, swigRelease
-
Constructor Details
-
N4BiasFieldCorrectionImageFilter
protected N4BiasFieldCorrectionImageFilter(long cPtr, boolean cMemoryOwn) -
N4BiasFieldCorrectionImageFilter
public N4BiasFieldCorrectionImageFilter()itk::simple::N4BiasFieldCorrectionImageFilter::N4BiasFieldCorrectionImageFilter() Default Constructor that takes no arguments and initializes default parameters
-
-
Method Details
-
getCPtr
-
swigRelease
-
finalize
protected void finalize()- Overrides:
finalize
in classImageFilter
-
delete
public void delete()virtual itk::simple::N4BiasFieldCorrectionImageFilter::~N4BiasFieldCorrectionImageFilter() Destructor- Overrides:
delete
in classImageFilter
-
setConvergenceThreshold
public void setConvergenceThreshold(double ConvergenceThreshold) Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetConvergenceThreshold(double ConvergenceThreshold) 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. -
getConvergenceThreshold
public double getConvergenceThreshold()double itk::simple::N4BiasFieldCorrectionImageFilter::GetConvergenceThreshold() const 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. -
setMaximumNumberOfIterations
Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetMaximumNumberOfIterations(std::vector< uint32_t > MaximumNumberOfIterations) Set the maximum number of iterations specified at each fitting level. Default = 50. -
getMaximumNumberOfIterations
std::vector<uint32_t> itk::simple::N4BiasFieldCorrectionImageFilter::GetMaximumNumberOfIterations() const Get the maximum number of iterations specified at each fitting level. Default = 50. -
setBiasFieldFullWidthAtHalfMaximum
public void setBiasFieldFullWidthAtHalfMaximum(double BiasFieldFullWidthAtHalfMaximum) Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetBiasFieldFullWidthAtHalfMaximum(double BiasFieldFullWidthAtHalfMaximum) Set the full width at half maximum parameter characterizing the width of the Gaussian deconvolution. Default = 0.15. -
getBiasFieldFullWidthAtHalfMaximum
public double getBiasFieldFullWidthAtHalfMaximum()double itk::simple::N4BiasFieldCorrectionImageFilter::GetBiasFieldFullWidthAtHalfMaximum() const Get the full width at half maximum parameter characterizing the width of the Gaussian deconvolution. Default = 0.15. -
setWienerFilterNoise
public void setWienerFilterNoise(double WienerFilterNoise) Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetWienerFilterNoise(double WienerFilterNoise) Set the noise estimate defining the Wiener filter. Default = 0.01. -
getWienerFilterNoise
public double getWienerFilterNoise()double itk::simple::N4BiasFieldCorrectionImageFilter::GetWienerFilterNoise() const Get the noise estimate defining the Wiener filter. Default = 0.01. -
setNumberOfHistogramBins
public void setNumberOfHistogramBins(long NumberOfHistogramBins) Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetNumberOfHistogramBins(uint32_t NumberOfHistogramBins) Set number of bins defining the log input intensity histogram. Default = 200. -
getNumberOfHistogramBins
public long getNumberOfHistogramBins()uint32_t itk::simple::N4BiasFieldCorrectionImageFilter::GetNumberOfHistogramBins() const Get number of bins defining the log input intensity histogram. Default = 200. -
setNumberOfControlPoints
Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetNumberOfControlPoints(uint32_t value) Set the values of the NumberOfControlPoints vector all to value -
setNumberOfControlPoints
public void setNumberOfControlPoints(long value) Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetNumberOfControlPoints(uint32_t value) Set the values of the NumberOfControlPoints vector all to value -
getNumberOfControlPoints
std::vector<uint32_t> itk::simple::N4BiasFieldCorrectionImageFilter::GetNumberOfControlPoints() const 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. -
setSplineOrder
public void setSplineOrder(long SplineOrder) Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetSplineOrder(uint32_t SplineOrder) Set the spline order defining the bias field estimate. Default = 3. -
getSplineOrder
public long getSplineOrder()uint32_t itk::simple::N4BiasFieldCorrectionImageFilter::GetSplineOrder() const Get the spline order defining the bias field estimate. Default = 3. -
setUseMaskLabel
public void setUseMaskLabel(boolean UseMaskLabel) Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetUseMaskLabel(bool UseMaskLabel) Use a mask label for identifying mask functionality. See SetMaskLabel. Defaults to true. -
useMaskLabelOn
public void useMaskLabelOn()Self& itk::simple::N4BiasFieldCorrectionImageFilter::UseMaskLabelOn() Set the value of UseMaskLabel to true or false respectfully. -
useMaskLabelOff
public void useMaskLabelOff()Self& itk::simple::N4BiasFieldCorrectionImageFilter::UseMaskLabelOff() -
getUseMaskLabel
public boolean getUseMaskLabel()bool itk::simple::N4BiasFieldCorrectionImageFilter::GetUseMaskLabel() const Use a mask label for identifying mask functionality. See SetMaskLabel. Defaults to true. -
setMaskLabel
public void setMaskLabel(short MaskLabel) Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetMaskLabel(uint8_t MaskLabel) 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. -
getMaskLabel
public short getMaskLabel()uint8_t itk::simple::N4BiasFieldCorrectionImageFilter::GetMaskLabel() const 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. -
getCurrentLevel
public long getCurrentLevel()uint32_t itk::simple::N4BiasFieldCorrectionImageFilter::GetCurrentLevel() const Get the current fitting level. This is a helper function for reporting observations. This is an active measurement. It may be accessed while the filter is being executing in command call-backs and can be accessed after execution. -
getElapsedIterations
public long getElapsedIterations()uint32_t itk::simple::N4BiasFieldCorrectionImageFilter::GetElapsedIterations() const Get the number of elapsed iterations. This is a helper function for reporting observations. This is an active measurement. It may be accessed while the filter is being executing in command call-backs and can be accessed after execution. -
getCurrentConvergenceMeasurement
public double getCurrentConvergenceMeasurement()double itk::simple::N4BiasFieldCorrectionImageFilter::GetCurrentConvergenceMeasurement() const Get the current convergence measurement. This is a helper function for reporting observations. This is an active measurement. It may be accessed while the filter is being executing in command call-backs and can be accessed after execution. -
getLogBiasFieldAsImage
Image itk::simple::N4BiasFieldCorrectionImageFilter::GetLogBiasFieldAsImage(Image referenceImage) const The computed log bias field correction. 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. Returns the b-spline log bias field reconstructioned onto the space of the referenceImage parameter. An input image can be corrected by: input/exp(bias_field). This is an active measurement. It may be accessed while the filter is being executing in command call-backs and can be accessed after execution. -
getName
std::string itk::simple::N4BiasFieldCorrectionImageFilter::GetName() const Name of this class- Overrides:
getName
in classProcessObject
-
toString
std::string itk::simple::N4BiasFieldCorrectionImageFilter::ToString() const Print ourselves out- Overrides:
toString
in classProcessObject
-
execute
Image itk::simple::N4BiasFieldCorrectionImageFilter::Execute(const Image &image, const Image &maskImage) Execute the filter on the input image -
execute
Image itk::simple::N4BiasFieldCorrectionImageFilter::Execute(const Image &image, const Image &maskImage) Execute the filter on the input image
-