SimpleITK  
Public Types | Public Member Functions | Private Types | Private Member Functions | Private Attributes | Friends | List of all members
itk::simple::N4BiasFieldCorrectionImageFilter Class Reference

Implementation of the N4 bias field correction algorithm. More...

#include <sitkN4BiasFieldCorrectionImageFilter.h>

+ Inheritance diagram for itk::simple::N4BiasFieldCorrectionImageFilter:
+ Collaboration diagram for itk::simple::N4BiasFieldCorrectionImageFilter:

Public Types

using PixelIDTypeList = RealPixelIDTypeList
 
using Self = N4BiasFieldCorrectionImageFilter
 
- Public Types inherited from itk::simple::ImageFilter
using Self = ImageFilter
 
- Public Types inherited from itk::simple::ProcessObject
using Self = ProcessObject
 

Public Member Functions

Image Execute (const Image &image)
 
Image Execute (const Image &image, const Image &maskImage)
 
double GetBiasFieldFullWidthAtHalfMaximum () const
 
double GetConvergenceThreshold () const
 
double GetCurrentConvergenceMeasurement () const
 
uint32_t GetCurrentLevel () const
 
uint32_t GetElapsedIterations () const
 
Image 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). More...
 
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 ()
 
SelfSetBiasFieldFullWidthAtHalfMaximum (double BiasFieldFullWidthAtHalfMaximum)
 
SelfSetConvergenceThreshold (double ConvergenceThreshold)
 
SelfSetMaskLabel (uint8_t MaskLabel)
 
SelfSetMaximumNumberOfIterations (std::vector< uint32_t > MaximumNumberOfIterations)
 
SelfSetNumberOfControlPoints (std::vector< uint32_t > NumberOfControlPoints)
 
SelfSetNumberOfControlPoints (uint32_t value)
 
SelfSetNumberOfHistogramBins (uint32_t NumberOfHistogramBins)
 
SelfSetSplineOrder (uint32_t SplineOrder)
 
SelfSetUseMaskLabel (bool UseMaskLabel)
 
SelfSetWienerFilterNoise (double WienerFilterNoise)
 
std::string ToString () const
 
SelfUseMaskLabelOff ()
 
SelfUseMaskLabelOn ()
 
virtual ~N4BiasFieldCorrectionImageFilter ()
 
- Public Member Functions inherited from itk::simple::ImageFilter
 ImageFilter ()
 
virtual ~ImageFilter ()=0
 
- Public Member Functions inherited from itk::simple::ProcessObject
virtual void Abort ()
 
virtual int AddCommand (itk::simple::EventEnum event, const std::function< void()> &func)
 Directly add a callback to observe an event. More...
 
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
 
virtual void SetNumberOfWorkUnits (unsigned int n)
 
virtual unsigned int GetNumberOfWorkUnits () const
 

Private Types

using MemberFunctionType = Image(Self::*)(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 {0.15}
 
double m_ConvergenceThreshold {0.001}
 
itk::ProcessObjectm_Filter {nullptr}
 
uint8_t m_MaskLabel {1}
 
std::vector< uint32_t > m_MaximumNumberOfIterations {std::vector<uint32_t>(4,50)}
 
std::unique_ptr< detail::MemberFunctionFactory< MemberFunctionType > > m_MemberFactory
 
std::vector< uint32_t > m_NumberOfControlPoints {std::vector<uint32_t>(3, 4)}
 
uint32_t m_NumberOfHistogramBins {200u}
 
std::function< double()> m_pfGetCurrentConvergenceMeasurement
 
std::function< uint32_t()> m_pfGetCurrentLevel
 
std::function< uint32_t()> m_pfGetElapsedIterations
 
std::function< Image(Image)> m_pfGetLogBiasFieldAsImage
 
uint32_t m_SplineOrder {3u}
 
bool m_UseMaskLabel {true}
 
double m_WienerFilterNoise {0.01}
 

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 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...
 
static bool SetGlobalDefaultThreader (const std::string &threader)
 Set/Get the default threader used for process objects. More...
 
static std::string GetGlobalDefaultThreader ()
 Set/Get the default threader used for process objects. More...
 
static void SetGlobalDefaultNumberOfThreads (unsigned int n)
 
static unsigned int GetGlobalDefaultNumberOfThreads ()
 Set/Get the default threader used for process objects. More...
 
- Protected Member Functions inherited from itk::simple::ImageFilter
void CheckImageMatchingDimension (const Image &image1, const Image &image2, const std::string &image2Name)
 
void CheckImageMatchingPixelType (const Image &image1, const Image &image2, const std::string &image2Name)
 
void CheckImageMatchingSize (const Image &image1, const Image &image2, const std::string &image2Name)
 
- Protected Member Functions inherited from itk::simple::ProcessObject
virtual unsigned long AddITKObserver (const itk::EventObject &, itk::Command *)
 
virtual itk::ProcessObjectGetActiveProcess ()
 
virtual void OnActiveProcessDelete ()
 
virtual void onCommandDelete (const itk::simple::Command *cmd) noexcept
 
virtual void PreUpdate (itk::ProcessObject *p)
 
virtual void RemoveITKObserver (EventCommand &e)
 
- Protected Member Functions inherited from itk::simple::NonCopyable
 NonCopyable ()=default
 
 NonCopyable (const NonCopyable &)=delete
 
NonCopyableoperator= (const NonCopyable &)=delete
 
- Static Protected Member Functions inherited from itk::simple::ImageFilter
template<class TImageType >
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 TPixelType , unsigned int VImageDimension, unsigned int VLength, template< typename, unsigned int > class TVector>
static Image CastITKToImage (itk::Image< TVector< TPixelType, VLength >, VImageDimension > *img)
 
template<unsigned int VImageDimension, unsigned int VLength, template< unsigned int > class TVector>
static Image CastITKToImage (itk::Image< TVector< VLength >, VImageDimension > *img)
 
template<class TImageType >
static Image CastITKToImage (TImageType *img)
 
static const itk::EventObjectGetITKEventObject (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)
 

Detailed Description

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/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.

Author
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 also
itk::simple::N4BiasFieldCorrection for the procedural interface
itk::N4BiasFieldCorrectionImageFilter for the Doxygen on the original ITK class.
Examples
N4BiasFieldCorrection/N4BiasFieldCorrection.cxx.

Definition at line 73 of file sitkN4BiasFieldCorrectionImageFilter.h.

Member Typedef Documentation

◆ MemberFunctionType

using itk::simple::N4BiasFieldCorrectionImageFilter::MemberFunctionType = Image (Self::*)( const Image * image, const Image * maskImage )
private

Setup for member function dispatching

Definition at line 238 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ PixelIDTypeList

Define the pixels types supported by this filter

Definition at line 85 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ Self

Definition at line 75 of file sitkN4BiasFieldCorrectionImageFilter.h.

Constructor & Destructor Documentation

◆ ~N4BiasFieldCorrectionImageFilter()

virtual itk::simple::N4BiasFieldCorrectionImageFilter::~N4BiasFieldCorrectionImageFilter ( )
virtual

Destructor

◆ N4BiasFieldCorrectionImageFilter()

itk::simple::N4BiasFieldCorrectionImageFilter::N4BiasFieldCorrectionImageFilter ( )

Default Constructor that takes no arguments and initializes default parameters

Member Function Documentation

◆ Execute() [1/2]

Image itk::simple::N4BiasFieldCorrectionImageFilter::Execute ( const Image image)

◆ Execute() [2/2]

Image itk::simple::N4BiasFieldCorrectionImageFilter::Execute ( const Image image,
const Image maskImage 
)

Execute the filter on the input image

Examples
N4BiasFieldCorrection/N4BiasFieldCorrection.cxx.

◆ ExecuteInternal()

template<class TImageType >
Image itk::simple::N4BiasFieldCorrectionImageFilter::ExecuteInternal ( const Image image,
const Image maskImage 
)
private

◆ GetBiasFieldFullWidthAtHalfMaximum()

double itk::simple::N4BiasFieldCorrectionImageFilter::GetBiasFieldFullWidthAtHalfMaximum ( ) const
inline

Get the full width at half maximum parameter characterizing the width of the Gaussian deconvolution. Default = 0.15.

Definition at line 116 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ GetConvergenceThreshold()

double itk::simple::N4BiasFieldCorrectionImageFilter::GetConvergenceThreshold ( ) const
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 96 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ GetCurrentConvergenceMeasurement()

double itk::simple::N4BiasFieldCorrectionImageFilter::GetCurrentConvergenceMeasurement ( ) const
inline

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.

Definition at line 209 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ GetCurrentLevel()

uint32_t itk::simple::N4BiasFieldCorrectionImageFilter::GetCurrentLevel ( ) const
inline

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.

Definition at line 191 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ GetElapsedIterations()

uint32_t itk::simple::N4BiasFieldCorrectionImageFilter::GetElapsedIterations ( ) const
inline

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.

Definition at line 200 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ GetLogBiasFieldAsImage()

Image itk::simple::N4BiasFieldCorrectionImageFilter::GetLogBiasFieldAsImage ( Image  referenceImage) const
inline

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.
Examples
N4BiasFieldCorrection/N4BiasFieldCorrection.cxx.

Definition at line 219 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ GetMaskLabel()

uint8_t itk::simple::N4BiasFieldCorrectionImageFilter::GetMaskLabel ( ) const
inline

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.

Definition at line 183 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ GetMaximumNumberOfIterations()

std::vector<uint32_t> itk::simple::N4BiasFieldCorrectionImageFilter::GetMaximumNumberOfIterations ( ) const
inline

Get the maximum number of iterations specified at each fitting level. Default = 50.

Definition at line 106 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ GetName()

std::string itk::simple::N4BiasFieldCorrectionImageFilter::GetName ( ) const
inlinevirtual

Name of this class

Implements itk::simple::ProcessObject.

Definition at line 223 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ GetNumberOfControlPoints()

std::vector<uint32_t> itk::simple::N4BiasFieldCorrectionImageFilter::GetNumberOfControlPoints ( ) const
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.

◆ GetNumberOfHistogramBins()

uint32_t itk::simple::N4BiasFieldCorrectionImageFilter::GetNumberOfHistogramBins ( ) const
inline

Get number of bins defining the log input intensity histogram. Default = 200.

Definition at line 136 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ GetSplineOrder()

uint32_t itk::simple::N4BiasFieldCorrectionImageFilter::GetSplineOrder ( ) const
inline

Get the spline order defining the bias field estimate. Default = 3.

Definition at line 159 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ GetUseMaskLabel()

bool itk::simple::N4BiasFieldCorrectionImageFilter::GetUseMaskLabel ( ) const
inline

Use a mask label for identifying mask functionality. See SetMaskLabel. Defaults to true.

Definition at line 173 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ GetWienerFilterNoise()

double itk::simple::N4BiasFieldCorrectionImageFilter::GetWienerFilterNoise ( ) const
inline

Get the noise estimate defining the Wiener filter. Default = 0.01.

Definition at line 126 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ SetBiasFieldFullWidthAtHalfMaximum()

Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetBiasFieldFullWidthAtHalfMaximum ( double  BiasFieldFullWidthAtHalfMaximum)
inline

Set the full width at half maximum parameter characterizing the width of the Gaussian deconvolution. Default = 0.15.

Definition at line 111 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ SetConvergenceThreshold()

Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetConvergenceThreshold ( double  ConvergenceThreshold)
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 91 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ SetMaskLabel()

Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetMaskLabel ( uint8_t  MaskLabel)
inline

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.

Definition at line 178 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ SetMaximumNumberOfIterations()

Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetMaximumNumberOfIterations ( std::vector< uint32_t >  MaximumNumberOfIterations)
inline

Set the maximum number of iterations specified at each fitting level. Default = 50.

Examples
N4BiasFieldCorrection/N4BiasFieldCorrection.cxx.

Definition at line 101 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ SetNumberOfControlPoints() [1/2]

Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetNumberOfControlPoints ( std::vector< uint32_t >  NumberOfControlPoints)
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 141 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ SetNumberOfControlPoints() [2/2]

Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetNumberOfControlPoints ( uint32_t  value)
inline

Set the values of the NumberOfControlPoints vector all to value

Definition at line 144 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ SetNumberOfHistogramBins()

Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetNumberOfHistogramBins ( uint32_t  NumberOfHistogramBins)
inline

Set number of bins defining the log input intensity histogram. Default = 200.

Definition at line 131 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ SetSplineOrder()

Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetSplineOrder ( uint32_t  SplineOrder)
inline

Set the spline order defining the bias field estimate. Default = 3.

Definition at line 154 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ SetUseMaskLabel()

Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetUseMaskLabel ( bool  UseMaskLabel)
inline

Use a mask label for identifying mask functionality. See SetMaskLabel. Defaults to true.

Definition at line 164 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ SetWienerFilterNoise()

Self& itk::simple::N4BiasFieldCorrectionImageFilter::SetWienerFilterNoise ( double  WienerFilterNoise)
inline

Set the noise estimate defining the Wiener filter. Default = 0.01.

Definition at line 121 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ ToString()

std::string itk::simple::N4BiasFieldCorrectionImageFilter::ToString ( ) const
virtual

Print ourselves out

Reimplemented from itk::simple::ProcessObject.

◆ UseMaskLabelOff()

Self& itk::simple::N4BiasFieldCorrectionImageFilter::UseMaskLabelOff ( )
inline

Definition at line 168 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ UseMaskLabelOn()

Self& itk::simple::N4BiasFieldCorrectionImageFilter::UseMaskLabelOn ( )
inline

Set the value of UseMaskLabel to true or false respectfully.

Definition at line 167 of file sitkN4BiasFieldCorrectionImageFilter.h.

Friends And Related Function Documentation

◆ detail::MemberFunctionAddressor< MemberFunctionType >

Definition at line 242 of file sitkN4BiasFieldCorrectionImageFilter.h.

Member Data Documentation

◆ m_BiasFieldFullWidthAtHalfMaximum

double itk::simple::N4BiasFieldCorrectionImageFilter::m_BiasFieldFullWidthAtHalfMaximum {0.15}
private

Definition at line 254 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ m_ConvergenceThreshold

double itk::simple::N4BiasFieldCorrectionImageFilter::m_ConvergenceThreshold {0.001}
private

Definition at line 248 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ m_Filter

itk::ProcessObject* itk::simple::N4BiasFieldCorrectionImageFilter::m_Filter {nullptr}
private

Definition at line 282 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ m_MaskLabel

uint8_t itk::simple::N4BiasFieldCorrectionImageFilter::m_MaskLabel {1}
private

Definition at line 270 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ m_MaximumNumberOfIterations

std::vector<uint32_t> itk::simple::N4BiasFieldCorrectionImageFilter::m_MaximumNumberOfIterations {std::vector<uint32_t>(4,50)}
private

Definition at line 251 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ m_MemberFactory

std::unique_ptr<detail::MemberFunctionFactory<MemberFunctionType> > itk::simple::N4BiasFieldCorrectionImageFilter::m_MemberFactory
private

Definition at line 244 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ m_NumberOfControlPoints

std::vector<uint32_t> itk::simple::N4BiasFieldCorrectionImageFilter::m_NumberOfControlPoints {std::vector<uint32_t>(3, 4)}
private

Definition at line 263 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ m_NumberOfHistogramBins

uint32_t itk::simple::N4BiasFieldCorrectionImageFilter::m_NumberOfHistogramBins {200u}
private

Definition at line 260 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ m_pfGetCurrentConvergenceMeasurement

std::function<double()> itk::simple::N4BiasFieldCorrectionImageFilter::m_pfGetCurrentConvergenceMeasurement
private

Definition at line 277 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ m_pfGetCurrentLevel

std::function<uint32_t()> itk::simple::N4BiasFieldCorrectionImageFilter::m_pfGetCurrentLevel
private

Definition at line 273 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ m_pfGetElapsedIterations

std::function<uint32_t()> itk::simple::N4BiasFieldCorrectionImageFilter::m_pfGetElapsedIterations
private

Definition at line 275 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ m_pfGetLogBiasFieldAsImage

std::function<Image(Image)> itk::simple::N4BiasFieldCorrectionImageFilter::m_pfGetLogBiasFieldAsImage
private

Definition at line 279 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ m_SplineOrder

uint32_t itk::simple::N4BiasFieldCorrectionImageFilter::m_SplineOrder {3u}
private

Definition at line 266 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ m_UseMaskLabel

bool itk::simple::N4BiasFieldCorrectionImageFilter::m_UseMaskLabel {true}
private

Definition at line 268 of file sitkN4BiasFieldCorrectionImageFilter.h.

◆ m_WienerFilterNoise

double itk::simple::N4BiasFieldCorrectionImageFilter::m_WienerFilterNoise {0.01}
private

Definition at line 257 of file sitkN4BiasFieldCorrectionImageFilter.h.


The documentation for this class was generated from the following file: