SimpleITK
1.2.4
|
Implementation of the N4 bias field correction algorithm. More...
#include <sitkN4BiasFieldCorrectionImageFilter.h>
Public Types | |
typedef RealPixelIDTypeList | PixelIDTypeList |
typedef N4BiasFieldCorrectionImageFilter | Self |
Public Types inherited from itk::simple::ImageFilter< 0 > | |
typedef ImageFilter | Self |
Public Types inherited from itk::simple::ProcessObject | |
typedef ProcessObject | Self |
Public Member Functions | |
Image | Execute (const Image &image, const Image &maskImage) |
Image | Execute (const Image &image) |
Image | Execute (const Image &image, const Image &maskImage, double convergenceThreshold, std::vector< uint32_t > maximumNumberOfIterations, double biasFieldFullWidthAtHalfMaximum, double wienerFilterNoise, uint32_t numberOfHistogramBins, const std::vector< uint32_t > &numberOfControlPoints, uint32_t splineOrder, bool useMaskLabel, uint8_t maskLabel) |
Image | Execute (const Image &image, double convergenceThreshold, std::vector< uint32_t > maximumNumberOfIterations, double biasFieldFullWidthAtHalfMaximum, double wienerFilterNoise, uint32_t numberOfHistogramBins, const std::vector< uint32_t > &numberOfControlPoints, uint32_t splineOrder, bool useMaskLabel, uint8_t maskLabel) |
double | GetBiasFieldFullWidthAtHalfMaximum () const |
double | GetConvergenceThreshold () const |
uint8_t | GetMaskLabel () const |
std::vector< uint32_t > | GetMaximumNumberOfIterations () const |
std::string | GetName () const |
std::vector< uint32_t > | GetNumberOfControlPoints () const |
uint32_t | GetNumberOfHistogramBins () const |
uint32_t | GetSplineOrder () const |
bool | GetUseMaskLabel () const |
double | GetWienerFilterNoise () const |
N4BiasFieldCorrectionImageFilter () | |
Self & | SetBiasFieldFullWidthAtHalfMaximum (double BiasFieldFullWidthAtHalfMaximum) |
Self & | SetConvergenceThreshold (double ConvergenceThreshold) |
Self & | SetMaskLabel (uint8_t MaskLabel) |
Self & | SetMaximumNumberOfIterations (std::vector< uint32_t > MaximumNumberOfIterations) |
Self & | SetNumberOfControlPoints (const std::vector< uint32_t > &NumberOfControlPoints) |
Self & | SetNumberOfControlPoints (uint32_t value) |
Self & | SetNumberOfHistogramBins (uint32_t NumberOfHistogramBins) |
Self & | SetSplineOrder (uint32_t SplineOrder) |
Self & | SetUseMaskLabel (bool UseMaskLabel) |
Self & | SetWienerFilterNoise (double WienerFilterNoise) |
std::string | ToString () const |
Self & | UseMaskLabelOff () |
Self & | UseMaskLabelOn () |
virtual | ~N4BiasFieldCorrectionImageFilter () |
Public Member Functions inherited from itk::simple::ImageFilter< 0 > | |
ImageFilter () | |
virtual | ~ImageFilter ()=0 |
Public Member Functions inherited from itk::simple::ProcessObject | |
virtual void | Abort () |
virtual int | AddCommand (itk::simple::EventEnum event, itk::simple::Command &cmd) |
Add a Command Object to observer the event. More... | |
virtual float | GetProgress () const |
An Active Measurement of the progress of execution. More... | |
virtual bool | HasCommand (itk::simple::EventEnum event) const |
Query of this object has any registered commands for event. More... | |
ProcessObject () | |
virtual void | RemoveAllCommands () |
Remove all registered commands. More... | |
virtual | ~ProcessObject () |
virtual void | DebugOn () |
virtual void | DebugOff () |
virtual bool | GetDebug () const |
virtual void | SetDebug (bool debugFlag) |
virtual void | SetNumberOfThreads (unsigned int n) |
virtual unsigned int | GetNumberOfThreads () const |
Private Types | |
typedef Image(Self::* | MemberFunctionType) (const Image *image, const Image *maskImage) |
Private Member Functions | |
template<class TImageType > | |
Image | ExecuteInternal (const Image *image, const Image *maskImage) |
Private Attributes | |
double | m_BiasFieldFullWidthAtHalfMaximum |
double | m_ConvergenceThreshold |
uint8_t | m_MaskLabel |
std::vector< uint32_t > | m_MaximumNumberOfIterations |
nsstd::auto_ptr< detail::MemberFunctionFactory< MemberFunctionType > > | m_MemberFactory |
std::vector< uint32_t > | m_NumberOfControlPoints |
uint32_t | m_NumberOfHistogramBins |
uint32_t | m_SplineOrder |
bool | m_UseMaskLabel |
double | m_WienerFilterNoise |
Friends | |
struct | detail::MemberFunctionAddressor< MemberFunctionType > |
Additional Inherited Members | |
Static Public Member Functions inherited from itk::simple::ProcessObject | |
static bool | GetGlobalDefaultDebug () |
static void | GlobalDefaultDebugOff () |
static void | GlobalDefaultDebugOn () |
static void | SetGlobalDefaultDebug (bool debugFlag) |
static void | GlobalWarningDisplayOn () |
static void | GlobalWarningDisplayOff () |
static void | SetGlobalWarningDisplay (bool flag) |
static bool | GetGlobalWarningDisplay () |
static void | SetGlobalDefaultNumberOfThreads (unsigned int n) |
static unsigned int | GetGlobalDefaultNumberOfThreads () |
static double | GetGlobalDefaultCoordinateTolerance () |
Access the global tolerance to determine congruent spaces. More... | |
static void | SetGlobalDefaultCoordinateTolerance (double) |
Access the global tolerance to determine congruent spaces. More... | |
static double | GetGlobalDefaultDirectionTolerance () |
Access the global tolerance to determine congruent spaces. More... | |
static void | SetGlobalDefaultDirectionTolerance (double) |
Access the global tolerance to determine congruent spaces. More... | |
Protected Member Functions inherited from itk::simple::ProcessObject | |
virtual unsigned long | AddITKObserver (const itk::EventObject &, itk::Command *) |
virtual itk::ProcessObject * | GetActiveProcess () |
virtual void | OnActiveProcessDelete () |
virtual void | onCommandDelete (const itk::simple::Command *cmd) SITK_NOEXCEPT |
virtual void | PreUpdate (itk::ProcessObject *p) |
virtual void | RemoveITKObserver (EventCommand &e) |
Protected Member Functions inherited from itk::simple::NonCopyable | |
NonCopyable () | |
Static Protected Member Functions inherited from itk::simple::ImageFilter< 0 > | |
static void | FixNonZeroIndex (TImageType *img) |
Static Protected Member Functions inherited from itk::simple::ProcessObject | |
template<class TImageType > | |
static TImageType::ConstPointer | CastImageToITK (const Image &img) |
template<class TImageType > | |
static Image | CastITKToImage (TImageType *img) |
template<class TPixelType , unsigned int VImageDimension, unsigned int VLength, template< typename, unsigned int > class TVector> | |
static Image | CastITKToImage (itk::Image< TVector< TPixelType, VLength >, VImageDimension > *img) |
static const itk::EventObject & | GetITKEventObject (EventEnum e) |
template<typename T > | |
static std::ostream & | ToStringHelper (std::ostream &os, const T &v) |
static std::ostream & | ToStringHelper (std::ostream &os, const char &v) |
static std::ostream & | ToStringHelper (std::ostream &os, const signed char &v) |
static std::ostream & | ToStringHelper (std::ostream &os, const unsigned char &v) |
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 71 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
private |
Setup for member function dispatching
Definition at line 204 of file sitkN4BiasFieldCorrectionImageFilter.h.
Define the pixels types supported by this filter
Definition at line 83 of file sitkN4BiasFieldCorrectionImageFilter.h.
Definition at line 73 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
virtual |
Destructor
itk::simple::N4BiasFieldCorrectionImageFilter::N4BiasFieldCorrectionImageFilter | ( | ) |
Default Constructor that takes no arguments and initializes default parameters
Image itk::simple::N4BiasFieldCorrectionImageFilter::Execute | ( | const Image & | image, |
const Image & | maskImage | ||
) |
Execute the filter on the input image
Image itk::simple::N4BiasFieldCorrectionImageFilter::Execute | ( | const Image & | image, |
const Image & | maskImage, | ||
double | convergenceThreshold, | ||
std::vector< uint32_t > | maximumNumberOfIterations, | ||
double | biasFieldFullWidthAtHalfMaximum, | ||
double | wienerFilterNoise, | ||
uint32_t | numberOfHistogramBins, | ||
const std::vector< uint32_t > & | numberOfControlPoints, | ||
uint32_t | splineOrder, | ||
bool | useMaskLabel, | ||
uint8_t | maskLabel | ||
) |
Execute the filter on the input image with the given parameters
Image itk::simple::N4BiasFieldCorrectionImageFilter::Execute | ( | const Image & | image, |
double | convergenceThreshold, | ||
std::vector< uint32_t > | maximumNumberOfIterations, | ||
double | biasFieldFullWidthAtHalfMaximum, | ||
double | wienerFilterNoise, | ||
uint32_t | numberOfHistogramBins, | ||
const std::vector< uint32_t > & | numberOfControlPoints, | ||
uint32_t | splineOrder, | ||
bool | useMaskLabel, | ||
uint8_t | maskLabel | ||
) |
|
private |
|
inline |
Get the full width at half maximum parameter characterizing the width of the Gaussian deconvolution. Default = 0.15.
Definition at line 115 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
inline |
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.
Definition at line 95 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
inline |
DeprecatedSet/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.
Definition at line 183 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
inline |
Get the maximum number of iterations specified at each fitting level. Default = 50.
Definition at line 105 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
inlinevirtual |
Name of this class
Implements itk::simple::ProcessObject.
Definition at line 185 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
inline |
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.
Definition at line 149 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
inline |
Get number of bins defining the log input intensity histogram. Default = 200.
Definition at line 135 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
inline |
Get the spline order defining the bias field estimate. Default = 3.
Definition at line 159 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
inline |
Use a mask label for identifying mask functionality. See SetMaskLabel. Defaults to true.
Definition at line 173 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
inline |
Get the noise estimate defining the Wiener filter. Default = 0.01.
Definition at line 125 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
inline |
Set the full width at half maximum parameter characterizing the width of the Gaussian deconvolution. Default = 0.15.
Definition at line 110 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
inline |
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.
Definition at line 90 of file sitkN4BiasFieldCorrectionImageFilter.h.
DeprecatedSet/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.
Definition at line 178 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
inline |
Set the maximum number of iterations specified at each fitting level. Default = 50.
Definition at line 100 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
inline |
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.
Definition at line 140 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
inline |
Set the values of the NumberOfControlPoints vector all to value
Definition at line 143 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
inline |
Set number of bins defining the log input intensity histogram. Default = 200.
Definition at line 130 of file sitkN4BiasFieldCorrectionImageFilter.h.
Set the spline order defining the bias field estimate. Default = 3.
Definition at line 154 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
inline |
Use a mask label for identifying mask functionality. See SetMaskLabel. Defaults to true.
Definition at line 164 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
inline |
Set the noise estimate defining the Wiener filter. Default = 0.01.
Definition at line 120 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
virtual |
Print ourselves out
Reimplemented from itk::simple::ProcessObject.
|
inline |
Definition at line 168 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
inline |
Set the value of UseMaskLabel to true or false respectfully.
Definition at line 167 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
friend |
Definition at line 208 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 218 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 214 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 228 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 216 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 210 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 224 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 222 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 226 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 227 of file sitkN4BiasFieldCorrectionImageFilter.h.
|
private |
Definition at line 220 of file sitkN4BiasFieldCorrectionImageFilter.h.