ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkMRIBiasFieldCorrectionFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * http://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 #ifndef itkMRIBiasFieldCorrectionFilter_h
19 #define itkMRIBiasFieldCorrectionFilter_h
20 
21 #include <ctime>
22 
23 #include "itkImageToImageFilter.h"
24 #include "itkArray2D.h"
25 #include "itkMRASlabIdentifier.h"
30 #include "itkImageRegionIterator.h"
31 
32 namespace itk
33 {
46 template< typename TImage, typename TImageMask, typename TBiasField >
47 class ITK_TEMPLATE_EXPORT MRIBiasEnergyFunction : public SingleValuedCostFunction
48 {
49 public:
50  ITK_DISALLOW_COPY_AND_ASSIGN(MRIBiasEnergyFunction);
51 
57 
60 
62  itkNewMacro(Self);
63 
65  using ImageType = TImage;
66  using ImagePointer = typename ImageType::Pointer;
67  using ImageElementType = typename ImageType::PixelType;
70  using MaskType = TImageMask;
71  using MaskPointer = typename MaskType::Pointer;
72  using MaskElementType = typename MaskType::PixelType;
73 
75  using BiasFieldType = TBiasField;
76 
79  using ParametersType = typename Superclass::ParametersType;
80 
82  using DerivativeType = Superclass::DerivativeType;
83 
85  using MeasureType = Superclass::MeasureType;
86 
87  static constexpr unsigned int SpaceDimension = 3;
88 
90  using InternalEnergyFunction = CompositeValleyFunction;
91 
93  typedef unsigned int SamplingFactorType[SpaceDimension];
94 
96  itkSetObjectMacro(Image, ImageType);
97 
99  itkSetObjectMacro(Mask, MaskType);
100 
102  itkSetMacro(Region, ImageRegionType);
103 
106  { m_BiasField = bias; }
107 
110  void SetSamplingFactors(SamplingFactorType factor)
111  {
112  for ( unsigned int i = 0; i < SpaceDimension; i++ )
113  {
114  m_SamplingFactor[i] = factor[i];
115  }
116  }
118 
121  double GetEnergy0(double diff)
122  {
123  return ( *m_InternalEnergyFunction )( diff );
124  }
125 
128  MeasureType GetValue(const ParametersType & parameters) const override;
129 
132  void GetDerivative( const ParametersType & itkNotUsed(parameters),
133  DerivativeType & itkNotUsed(derivative) ) const override
134  {}
135 
140  void InitializeDistributions(Array< double > classMeans,
141  Array< double > classSigmas);
142 
143  unsigned int GetNumberOfParameters() const override;
144 
145 protected:
148 
150  ~MRIBiasEnergyFunction() override;
151 
152 private:
155 
158 
161 
164 
167 
169  SamplingFactorType m_SamplingFactor;
170 }; // end of class
171 
222 template< typename TInputImage, typename TOutputImage, typename TMaskImage >
223 class ITK_TEMPLATE_EXPORT MRIBiasFieldCorrectionFilter :
224  public ImageToImageFilter< TInputImage, TOutputImage >
225 {
226 public:
227  ITK_DISALLOW_COPY_AND_ASSIGN(MRIBiasFieldCorrectionFilter);
228 
234 
236  itkNewMacro(Self);
237 
240 
242  static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
243 
245  using OutputImageType = TOutputImage;
246  using OutputImagePointer = typename TOutputImage::Pointer;
248  using OutputImagePixelType = typename TOutputImage::PixelType;
251 
252  using InputImageType = TInputImage;
253  using InputImagePointer = typename TInputImage::Pointer;
255  using InputImagePixelType = typename TInputImage::PixelType;
258 
260  using ImageMaskType = TMaskImage;
261  using ImageMaskPointer = typename ImageMaskType::Pointer;
263 
269 
273  using SlabRegionVectorIteratorType = typename SlabRegionVectorType::iterator;
274 
277 
283 
286 
289 
292 
296  void SetInputMask(ImageMaskType *inputMask);
297  itkGetModifiableObjectMacro(InputMask, ImageMaskType);
299 
302  void SetOutputMask(ImageMaskType *outputMask);
303  itkGetModifiableObjectMacro(OutputMask, ImageMaskType);
305 
306 #if defined ( ITK_LEGACY_REMOVE )
307 
311  { m_BiasFieldMultiplicative = flag; }
312 #endif // defined ( ITK_LEGACY_REMOVE )
313 
317  itkSetMacro(BiasFieldMultiplicative, bool);
318  itkGetConstMacro(BiasFieldMultiplicative, bool);
319  itkBooleanMacro(BiasFieldMultiplicative);
321 
322 #if defined ( ITK_LEGACY_REMOVE )
323 
325  { return m_BiasFieldMultiplicative; }
326 #endif // defined ( ITK_LEGACY_REMOVE )
327 
332  itkSetMacro(UsingInterSliceIntensityCorrection, bool);
333  itkGetConstMacro(UsingInterSliceIntensityCorrection, bool);
334  itkBooleanMacro(UsingInterSliceIntensityCorrection);
336 
342  itkSetMacro(UsingSlabIdentification, bool);
343  itkGetConstMacro(UsingSlabIdentification, bool);
344  itkBooleanMacro(UsingSlabIdentification);
346 
347  itkSetMacro(SlabBackgroundMinimumThreshold, InputImagePixelType);
348  itkGetConstReferenceMacro(SlabBackgroundMinimumThreshold,
349  InputImagePixelType);
350 
351  itkSetMacro(SlabNumberOfSamples, unsigned int);
352  itkGetConstReferenceMacro(SlabNumberOfSamples, unsigned int);
353 
354  itkSetMacro(SlabTolerance, double);
355  itkGetConstReferenceMacro(SlabTolerance, double);
356 
362  itkSetMacro(UsingBiasFieldCorrection, bool);
363  itkGetConstMacro(UsingBiasFieldCorrection, bool);
364  itkBooleanMacro(UsingBiasFieldCorrection);
366 
369  itkSetMacro(GeneratingOutput, bool);
370  itkGetConstMacro(GeneratingOutput, bool);
371  itkBooleanMacro(GeneratingOutput);
373 
376  itkSetMacro(SlicingDirection, int);
377  itkGetConstMacro(SlicingDirection, int);
379 
381  itkSetMacro(BiasFieldDegree, int);
382  itkGetConstMacro(BiasFieldDegree, int);
384 
389  & coefficients)
390  { this->Modified(); m_BiasFieldCoefficients = coefficients; }
391 
396  { return m_EstimatedBiasFieldCoefficients; }
397 
401  void SetTissueClassStatistics(const Array< double > & means,
402  const Array< double > & sigmas);
403 
406  itkSetMacro(VolumeCorrectionMaximumIteration, int);
407  itkGetConstMacro(VolumeCorrectionMaximumIteration, int);
409 
412  itkSetMacro(InterSliceCorrectionMaximumIteration, int);
413  itkGetConstMacro(InterSliceCorrectionMaximumIteration, int);
415 
417  itkSetMacro(OptimizerInitialRadius, double);
418  itkGetConstMacro(OptimizerInitialRadius, double);
420 
422  itkSetMacro(OptimizerGrowthFactor, double);
423  itkGetConstMacro(OptimizerGrowthFactor, double);
425 
428  itkSetMacro(OptimizerShrinkFactor, double);
429  itkGetConstMacro(OptimizerShrinkFactor, double);
430 
437  void SetNumberOfLevels(unsigned int num);
438 
440  itkGetConstMacro(NumberOfLevels, unsigned int);
441 
448  void SetSchedule(const ScheduleType & schedule);
449 
451  itkGetConstReferenceMacro(Schedule, ScheduleType);
452 
457  void SetStartingShrinkFactors(unsigned int factor);
458 
459  void SetStartingShrinkFactors(unsigned int *factors);
460 
462  const unsigned int * GetStartingShrinkFactors() const;
463 
467  static bool IsScheduleDownwardDivisible(const ScheduleType & schedule);
468 
475  void Initialize();
476 
479  BiasFieldType EstimateBiasField(InputImageRegionType region,
480  unsigned int degree,
481  int maximumIteration);
482 
486  void CorrectImage(BiasFieldType & bias,
487  InputImageRegionType region);
488 
491  void CorrectInterSliceIntensityInhomogeneity(InputImageRegionType region);
492 
493 protected:
495  ~MRIBiasFieldCorrectionFilter() override = default;
496  void PrintSelf(std::ostream & os, Indent indent) const override;
497 
500  bool CheckMaskImage(ImageMaskType *mask);
501 
502 protected:
506  void Log1PImage(InternalImageType *source,
507  InternalImageType *target);
508 
511  void ExpImage(InternalImageType *source,
512  InternalImageType *target);
513 
516  template< typename TSource, typename TTarget >
517  void CopyAndConvertImage(const TSource *source,
518  TTarget *target,
519  typename TTarget::RegionType requestedRegion)
520  {
521  using SourceIterator = ImageRegionConstIterator< TSource >;
522  using TargetIterator = ImageRegionIterator< TTarget >;
523  using TargetPixelType = typename TTarget::PixelType;
524 
525  SourceIterator s_iter(source, requestedRegion);
526  TargetIterator t_iter(target, requestedRegion);
527 
528  s_iter.GoToBegin();
529  t_iter.GoToBegin();
530  while ( !s_iter.IsAtEnd() )
531  {
532  t_iter.Set( static_cast< TargetPixelType >( s_iter.Get() ) );
533  ++s_iter;
534  ++t_iter;
535  }
536  }
537 
542  void GetBiasFieldSize(InputImageRegionType region,
543  BiasFieldType::DomainSizeType & domainSize);
544 
548  void AdjustSlabRegions(SlabRegionVectorType & slabs,
549  OutputImageRegionType requestedRegion);
550 
551  void GenerateData() override;
552 
553 private:
556 
559 
562 
565 
568 
571 
574 
576  bool m_BiasFieldMultiplicative{ true };
577 
580  bool m_UsingSlabIdentification{ false };
581  bool m_UsingBiasFieldCorrection{ true };
582  bool m_GeneratingOutput{ true };
583 
584  unsigned int m_SlabNumberOfSamples{ 200 };
586  double m_SlabTolerance{ 0.0 };
587 
589  int m_BiasFieldDegree{ 3 };
590 
592  unsigned int m_NumberOfLevels{ 0 };
593 
596 
600 
604 
605  int m_VolumeCorrectionMaximumIteration{ 2000 };
606 
607  int m_InterSliceCorrectionMaximumIteration{ 4000 };
608 
610  double m_OptimizerInitialRadius{ 1.01 };
611 
613  double m_OptimizerGrowthFactor{ 1.05 };
614 
617 
620 
623 };
624 
625 } // end namespace itk
626 
627 #ifndef ITK_MANUAL_INSTANTIATION
628 #include "itkMRIBiasFieldCorrectionFilter.hxx"
629 #endif
630 
631 #endif
Array class with size defined at construction time.
Definition: itkArray.h:46
std::vector< ImageRegionType > SlabRegionVectorType
typename TOutputImage::IndexType OutputImageIndexType
TPixel PixelType
Definition: itkImage.h:95
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
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...
Base class for all process objects that output image data.
CompositeValleyFunction InternalEnergyFunction
typename ImageType::Pointer ImagePointer
typename OutputImageType::PixelType OutputImagePixelType
typename ImageType::RegionType ImageRegionType
typename InputImageType::PixelType InputImagePixelType
typename TInputImage::SizeType InputImageSizeType
BiasFieldType::CoefficientArrayType GetEstimatedBiasFieldCoefficients()
typename EnergyFunctionType::Pointer EnergyFunctionPointer
typename InputImageType::Pointer InputImagePointer
BiasFieldType::CoefficientArrayType m_BiasFieldCoefficients
typename MRASlabIdentifierType::SlabRegionVectorType SlabRegionVectorType
typename OutputImageType::RegionType OutputImageRegionType
typename MaskType::PixelType MaskElementType
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
A region represents some portion or piece of data.
Definition: itkRegion.h:64
typename ImageMaskType::Pointer ImageMaskPointer
identifies slab in MR images comparing minimum intensity averages
typename SlabRegionVectorType::iterator SlabRegionVectorIteratorType
void GetDerivative(const ParametersType &, DerivativeType &) const override
typename TInputImage::IndexType InputImageIndexType
typename InputImageType::RegionType InputImageRegionType
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.
Definition: itkIndent.h:49
typename ImageMaskType::RegionType ImageMaskRegionType
typename ImageType::IndexType ImageIndexType
BiasFieldType::CoefficientArrayType m_EstimatedBiasFieldCoefficients
void SetInitialBiasFieldCoefficients(const BiasFieldType::CoefficientArrayType &coefficients)
typename MaskType::Pointer MaskPointer
Templated n-dimensional image class.
Definition: itkImage.h:75
A multi-dimensional iterator templated over image type that walks a region of pixels.
typename TOutputImage::SizeType OutputImageSizeType