18 #ifndef itkMRIBiasFieldCorrectionFilter_h
19 #define itkMRIBiasFieldCorrectionFilter_h
50 template <
typename TImage,
typename TImageMask,
typename TBiasField>
83 using typename Superclass::ParametersType;
91 static constexpr
unsigned int SpaceDimension = 3;
120 for (
unsigned int i = 0; i < SpaceDimension; ++i)
122 m_SamplingFactor[i] = factor[i];
132 return (*m_InternalEnergyFunction)(diff);
138 GetValue(
const ParametersType & parameters)
const override;
154 GetNumberOfParameters()
const override;
234 template <
typename TInputImage,
typename TOutputImage,
typename TMaskImage>
254 static constexpr
unsigned int ImageDimension = TOutputImage::ImageDimension;
318 #if defined(ITK_LEGACY_REMOVE)
325 m_BiasFieldMultiplicative = flag;
327 #endif // defined ( ITK_LEGACY_REMOVE )
332 itkSetMacro(BiasFieldMultiplicative,
bool);
333 itkGetConstMacro(BiasFieldMultiplicative,
bool);
334 itkBooleanMacro(BiasFieldMultiplicative);
337 #if defined(ITK_LEGACY_REMOVE)
342 return m_BiasFieldMultiplicative;
344 #endif // defined ( ITK_LEGACY_REMOVE )
350 itkSetMacro(UsingInterSliceIntensityCorrection,
bool);
351 itkGetConstMacro(UsingInterSliceIntensityCorrection,
bool);
352 itkBooleanMacro(UsingInterSliceIntensityCorrection);
360 itkSetMacro(UsingSlabIdentification,
bool);
361 itkGetConstMacro(UsingSlabIdentification,
bool);
362 itkBooleanMacro(UsingSlabIdentification);
365 itkSetMacro(SlabBackgroundMinimumThreshold, InputImagePixelType);
366 itkGetConstReferenceMacro(SlabBackgroundMinimumThreshold, InputImagePixelType);
368 itkSetMacro(SlabNumberOfSamples,
unsigned int);
369 itkGetConstReferenceMacro(SlabNumberOfSamples,
unsigned int);
371 itkSetMacro(SlabTolerance,
double);
372 itkGetConstReferenceMacro(SlabTolerance,
double);
379 itkSetMacro(UsingBiasFieldCorrection,
bool);
380 itkGetConstMacro(UsingBiasFieldCorrection,
bool);
381 itkBooleanMacro(UsingBiasFieldCorrection);
386 itkSetMacro(GeneratingOutput,
bool);
387 itkGetConstMacro(GeneratingOutput,
bool);
388 itkBooleanMacro(GeneratingOutput);
393 itkSetMacro(SlicingDirection,
int);
394 itkGetConstMacro(SlicingDirection,
int);
398 itkSetMacro(BiasFieldDegree,
int);
399 itkGetConstMacro(BiasFieldDegree,
int);
408 m_BiasFieldCoefficients = coefficients;
415 BiasFieldType::CoefficientArrayType
418 return m_EstimatedBiasFieldCoefficients;
429 itkSetMacro(VolumeCorrectionMaximumIteration,
int);
430 itkGetConstMacro(VolumeCorrectionMaximumIteration,
int);
435 itkSetMacro(InterSliceCorrectionMaximumIteration,
int);
436 itkGetConstMacro(InterSliceCorrectionMaximumIteration,
int);
440 itkSetMacro(OptimizerInitialRadius,
double);
441 itkGetConstMacro(OptimizerInitialRadius,
double);
445 itkSetMacro(OptimizerGrowthFactor,
double);
446 itkGetConstMacro(OptimizerGrowthFactor,
double);
451 itkSetMacro(OptimizerShrinkFactor,
double);
452 itkGetConstMacro(OptimizerShrinkFactor,
double);
461 SetNumberOfLevels(
unsigned int num);
464 itkGetConstMacro(NumberOfLevels,
unsigned int);
473 SetSchedule(
const ScheduleType & schedule);
476 itkGetConstReferenceMacro(Schedule, ScheduleType);
483 SetStartingShrinkFactors(
unsigned int factor);
486 SetStartingShrinkFactors(
const unsigned int * factors);
490 GetStartingShrinkFactors()
const;
496 IsScheduleDownwardDivisible(
const ScheduleType & schedule);
510 EstimateBiasField(InputImageRegionType region,
unsigned int degree,
int maximumIteration);
516 CorrectImage(BiasFieldType & bias, InputImageRegionType region);
521 CorrectInterSliceIntensityInhomogeneity(InputImageRegionType region);
527 PrintSelf(std::ostream & os,
Indent indent)
const override;
532 CheckMaskImage(ImageMaskType * mask);
539 Log1PImage(InternalImageType * source, InternalImageType * target);
544 ExpImage(InternalImageType * source, InternalImageType * target);
548 template <
typename TSource,
typename TTarget>
554 using TargetPixelType =
typename TTarget::PixelType;
556 SourceIterator s_iter(source, requestedRegion);
557 TargetIterator t_iter(target, requestedRegion);
561 while (!s_iter.IsAtEnd())
563 t_iter.Set(static_cast<TargetPixelType>(s_iter.Get()));
574 GetBiasFieldSize(InputImageRegionType region, BiasFieldType::DomainSizeType & biasSize);
580 AdjustSlabRegions(SlabRegionVectorType & slabs, OutputImageRegionType requestedRegion);
583 GenerateData()
override;
605 int m_SlicingDirection{};
608 bool m_BiasFieldMultiplicative{
true };
611 bool m_UsingInterSliceIntensityCorrection{};
612 bool m_UsingSlabIdentification{
false };
613 bool m_UsingBiasFieldCorrection{
true };
614 bool m_GeneratingOutput{
true };
616 unsigned int m_SlabNumberOfSamples{ 200 };
618 double m_SlabTolerance{ 0.0 };
621 int m_BiasFieldDegree{ 3 };
624 unsigned int m_NumberOfLevels{ 0 };
637 int m_VolumeCorrectionMaximumIteration{ 2000 };
639 int m_InterSliceCorrectionMaximumIteration{ 4000 };
642 double m_OptimizerInitialRadius{ 1.01 };
645 double m_OptimizerGrowthFactor{ 1.05 };
648 double m_OptimizerShrinkFactor{};
659 #ifndef ITK_MANUAL_INSTANTIATION
660 # include "itkMRIBiasFieldCorrectionFilter.hxx"