18 #ifndef itkMRIBiasFieldCorrectionFilter_h
19 #define itkMRIBiasFieldCorrectionFilter_h
46 template<
typename TImage,
typename TImageMask,
typename TBiasField >
87 static constexpr
unsigned int SpaceDimension = 3;
93 typedef unsigned int SamplingFactorType[SpaceDimension];
106 { m_BiasField = bias; }
112 for (
unsigned int i = 0; i < SpaceDimension; i++ )
114 m_SamplingFactor[i] = factor[i];
123 return ( *m_InternalEnergyFunction )( diff );
128 MeasureType GetValue(
const ParametersType & parameters)
const override;
143 unsigned int GetNumberOfParameters()
const override;
222 template<
typename TInputImage,
typename TOutputImage,
typename TMaskImage >
242 static constexpr
unsigned int ImageDimension = TOutputImage::ImageDimension;
306 #if defined ( ITK_LEGACY_REMOVE )
311 { m_BiasFieldMultiplicative = flag; }
312 #endif // defined ( ITK_LEGACY_REMOVE )
317 itkSetMacro(BiasFieldMultiplicative,
bool);
318 itkGetConstMacro(BiasFieldMultiplicative,
bool);
319 itkBooleanMacro(BiasFieldMultiplicative);
322 #if defined ( ITK_LEGACY_REMOVE )
325 {
return m_BiasFieldMultiplicative; }
326 #endif // defined ( ITK_LEGACY_REMOVE )
332 itkSetMacro(UsingInterSliceIntensityCorrection,
bool);
333 itkGetConstMacro(UsingInterSliceIntensityCorrection,
bool);
334 itkBooleanMacro(UsingInterSliceIntensityCorrection);
342 itkSetMacro(UsingSlabIdentification,
bool);
343 itkGetConstMacro(UsingSlabIdentification,
bool);
344 itkBooleanMacro(UsingSlabIdentification);
347 itkSetMacro(SlabBackgroundMinimumThreshold, InputImagePixelType);
348 itkGetConstReferenceMacro(SlabBackgroundMinimumThreshold,
349 InputImagePixelType);
351 itkSetMacro(SlabNumberOfSamples,
unsigned int);
352 itkGetConstReferenceMacro(SlabNumberOfSamples,
unsigned int);
354 itkSetMacro(SlabTolerance,
double);
355 itkGetConstReferenceMacro(SlabTolerance,
double);
362 itkSetMacro(UsingBiasFieldCorrection,
bool);
363 itkGetConstMacro(UsingBiasFieldCorrection,
bool);
364 itkBooleanMacro(UsingBiasFieldCorrection);
369 itkSetMacro(GeneratingOutput,
bool);
370 itkGetConstMacro(GeneratingOutput,
bool);
371 itkBooleanMacro(GeneratingOutput);
376 itkSetMacro(SlicingDirection,
int);
377 itkGetConstMacro(SlicingDirection,
int);
381 itkSetMacro(BiasFieldDegree,
int);
382 itkGetConstMacro(BiasFieldDegree,
int);
390 { this->Modified(); m_BiasFieldCoefficients = coefficients; }
396 {
return m_EstimatedBiasFieldCoefficients; }
406 itkSetMacro(VolumeCorrectionMaximumIteration,
int);
407 itkGetConstMacro(VolumeCorrectionMaximumIteration,
int);
412 itkSetMacro(InterSliceCorrectionMaximumIteration,
int);
413 itkGetConstMacro(InterSliceCorrectionMaximumIteration,
int);
417 itkSetMacro(OptimizerInitialRadius,
double);
418 itkGetConstMacro(OptimizerInitialRadius,
double);
422 itkSetMacro(OptimizerGrowthFactor,
double);
423 itkGetConstMacro(OptimizerGrowthFactor,
double);
428 itkSetMacro(OptimizerShrinkFactor,
double);
429 itkGetConstMacro(OptimizerShrinkFactor,
double);
437 void SetNumberOfLevels(
unsigned int num);
440 itkGetConstMacro(NumberOfLevels,
unsigned int);
448 void SetSchedule(
const ScheduleType & schedule);
451 itkGetConstReferenceMacro(Schedule, ScheduleType);
457 void SetStartingShrinkFactors(
unsigned int factor);
459 void SetStartingShrinkFactors(
unsigned int *factors);
462 const unsigned int * GetStartingShrinkFactors()
const;
467 static bool IsScheduleDownwardDivisible(
const ScheduleType & schedule);
479 BiasFieldType EstimateBiasField(InputImageRegionType region,
481 int maximumIteration);
486 void CorrectImage(BiasFieldType & bias,
487 InputImageRegionType region);
491 void CorrectInterSliceIntensityInhomogeneity(InputImageRegionType region);
496 void PrintSelf(std::ostream & os,
Indent indent)
const override;
500 bool CheckMaskImage(ImageMaskType *mask);
506 void Log1PImage(InternalImageType *source,
507 InternalImageType *target);
511 void ExpImage(InternalImageType *source,
512 InternalImageType *target);
516 template<
typename TSource,
typename TTarget >
523 using TargetPixelType =
typename TTarget::PixelType;
525 SourceIterator s_iter(source, requestedRegion);
526 TargetIterator t_iter(target, requestedRegion);
530 while ( !s_iter.IsAtEnd() )
532 t_iter.Set( static_cast< TargetPixelType >( s_iter.Get() ) );
542 void GetBiasFieldSize(InputImageRegionType region,
543 BiasFieldType::DomainSizeType & domainSize);
548 void AdjustSlabRegions(SlabRegionVectorType & slabs,
549 OutputImageRegionType requestedRegion);
551 void GenerateData()
override;
576 bool m_BiasFieldMultiplicative{
true };
580 bool m_UsingSlabIdentification{
false };
581 bool m_UsingBiasFieldCorrection{
true };
582 bool m_GeneratingOutput{
true };
584 unsigned int m_SlabNumberOfSamples{ 200 };
586 double m_SlabTolerance{ 0.0 };
589 int m_BiasFieldDegree{ 3 };
592 unsigned int m_NumberOfLevels{ 0 };
605 int m_VolumeCorrectionMaximumIteration{ 2000 };
607 int m_InterSliceCorrectionMaximumIteration{ 4000 };
610 double m_OptimizerInitialRadius{ 1.01 };
613 double m_OptimizerGrowthFactor{ 1.05 };
627 #ifndef ITK_MANUAL_INSTANTIATION
628 #include "itkMRIBiasFieldCorrectionFilter.hxx"
Array class with size defined at construction time.
std::vector< ImageRegionType > SlabRegionVectorType
typename TOutputImage::IndexType OutputImageIndexType
This class is a base for the CostFunctions returning a single value.
typename OutputImageType::Pointer OutputImagePointer
void CopyAndConvertImage(const TSource *source, TTarget *target, typename TTarget::RegionType requestedRegion)
Light weight base class for most itk classes.
NormalVariateGeneratorType::Pointer m_NormalVariateGenerator
typename ImageType::PixelType ImageElementType
Array< double > m_TissueClassMeans
void IsBiasFieldMultiplicative(bool flag)
Array< double > m_TissueClassSigmas
typename InternalImageType::Pointer InternalImagePointer
typename InternalImageType::RegionType InternalImageRegionType
typename InternalImageType::PixelType InternalImagePixelType
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
ImageMaskPointer m_InputMask
double GetEnergy0(double diff)
Base class for all process objects that output image data.
CompositeValleyFunction InternalEnergyFunction
TInputImage InputImageType
typename ImageType::Pointer ImagePointer
InternalImagePointer m_InternalInput
typename OutputImageType::PixelType OutputImagePixelType
typename ImageType::RegionType ImageRegionType
typename InputImageType::PixelType InputImagePixelType
typename TInputImage::SizeType InputImageSizeType
BiasFieldType::CoefficientArrayType GetEstimatedBiasFieldCoefficients()
EnergyFunctionPointer m_EnergyFunction
InputImagePixelType m_SlabBackgroundMinimumThreshold
typename EnergyFunctionType::Pointer EnergyFunctionPointer
void SetBiasField(BiasFieldType *bias)
typename InputImageType::Pointer InputImagePointer
ImageBaseType::SizeType SizeType
BiasFieldType::CoefficientArrayType m_BiasFieldCoefficients
typename MRASlabIdentifierType::SlabRegionVectorType SlabRegionVectorType
double m_OptimizerShrinkFactor
typename OutputImageType::RegionType OutputImageRegionType
typename MaskType::PixelType MaskElementType
ImageBaseType::IndexType IndexType
A multi-dimensional iterator templated over image type that walks a region of pixels.
TOutputImage OutputImageType
InternalEnergyFunction * m_InternalEnergyFunction
Represents a cost function for MRI bias field correction optimization.
2D and 3D multivariate Legendre Polynomial
bool m_UsingInterSliceIntensityCorrection
A region represents some portion or piece of data.
typename ImageMaskType::Pointer ImageMaskPointer
identifies slab in MR images comparing minimum intensity averages
typename SlabRegionVectorType::iterator SlabRegionVectorIteratorType
1+1 evolutionary strategy optimizer
void GetDerivative(const ParametersType &, DerivativeType &) const override
bool IsBiasFieldMultiplicative()
typename TInputImage::IndexType InputImageIndexType
typename InputImageType::RegionType InputImageRegionType
ImageMaskPointer m_OutputMask
Base class for filters that take an image as input and produce an image as output.
void SetSamplingFactors(SamplingFactorType factor)
Control indentation during Print() invocation.
Corrects 3D MRI bias field.
DoubleArrayType CoefficientArrayType
SamplingFactorType m_SamplingFactor
typename ImageMaskType::RegionType ImageMaskRegionType
typename ImageType::IndexType ImageIndexType
SlabRegionVectorType m_Slabs
ImageBaseType::RegionType RegionType
BiasFieldType::CoefficientArrayType m_EstimatedBiasFieldCoefficients
void SetInitialBiasFieldCoefficients(const BiasFieldType::CoefficientArrayType &coefficients)
typename MaskType::Pointer MaskPointer
Templated n-dimensional image class.
A multi-dimensional iterator templated over image type that walks a region of pixels.
Normal random variate generator.
BiasFieldType * m_BiasField
typename TOutputImage::SizeType OutputImageSizeType