ITK  5.2.0
Insight Toolkit
itkN4BiasFieldCorrectionImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
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 itkN4BiasFieldCorrectionImageFilter_h
19 #define itkN4BiasFieldCorrectionImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 
23 #include "itkArray.h"
25 #include "itkPointSet.h"
26 #include "itkVector.h"
27 
28 #include "vnl/vnl_vector.h"
29 
30 namespace itk
31 {
32 
93 template <typename TInputImage,
94  typename TMaskImage = Image<unsigned char, TInputImage::ImageDimension>,
95  class TOutputImage = TInputImage>
96 class N4BiasFieldCorrectionImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
97 {
98 public:
99  ITK_DISALLOW_COPY_AND_MOVE(N4BiasFieldCorrectionImageFilter);
100 
106 
109 
111  itkNewMacro(Self);
112 
114  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
115 
117  using InputImageType = TInputImage;
118  using OutputImageType = TOutputImage;
119  using MaskImageType = TMaskImage;
120  using MaskPixelType = typename MaskImageType::PixelType;
121 
122  using RealType = float;
126 
133 
138 
140  void
142 
146  void
147  SetInput1(const InputImageType * image)
148  {
149  this->SetInput(image);
150  }
151 
157  itkSetInputMacro(MaskImage, MaskImageType);
158  void
159  SetInput2(const MaskImageType * mask)
160  {
161  this->SetMaskImage(mask);
162  }
164 
170  itkGetInputMacro(MaskImage, MaskImageType);
171 
180  itkSetMacro(MaskLabel, MaskPixelType);
181  itkGetConstMacro(MaskLabel, MaskPixelType);
183 
187  itkSetMacro(UseMaskLabel, bool);
188  itkGetConstMacro(UseMaskLabel, bool);
189  itkBooleanMacro(UseMaskLabel);
191 
202  itkSetInputMacro(ConfidenceImage, RealImageType);
203  void
204  SetInput3(const RealImageType * image)
205  {
206  this->SetConfidenceImage(image);
207  }
209 
220  itkGetInputMacro(ConfidenceImage, RealImageType);
221 
222  // Sharpen histogram parameters: in estimating the bias field, the
223  // first step is to sharpen the intensity histogram by Wiener deconvolution
224  // with a 1-D Gaussian. The following parameters define this operation.
225  // These default values in N4 match the default values in N3.
226 
231  itkSetMacro(NumberOfHistogramBins, unsigned int);
232 
237  itkGetConstMacro(NumberOfHistogramBins, unsigned int);
238 
242  itkSetMacro(WienerFilterNoise, RealType);
243 
247  itkGetConstMacro(WienerFilterNoise, RealType);
248 
253  itkSetMacro(BiasFieldFullWidthAtHalfMaximum, RealType);
254 
259  itkGetConstMacro(BiasFieldFullWidthAtHalfMaximum, RealType);
260 
261  // B-spline parameters governing the fitting routine
262 
266  itkSetMacro(SplineOrder, unsigned int);
267 
271  itkGetConstMacro(SplineOrder, unsigned int);
272 
280  itkSetMacro(NumberOfControlPoints, ArrayType);
281 
289  itkGetConstMacro(NumberOfControlPoints, ArrayType);
290 
297  itkSetMacro(NumberOfFittingLevels, ArrayType);
298 
305  void
306  SetNumberOfFittingLevels(unsigned int n)
307  {
308  ArrayType nlevels;
309 
310  nlevels.Fill(n);
311  this->SetNumberOfFittingLevels(nlevels);
312  }
313 
320  itkGetConstMacro(NumberOfFittingLevels, ArrayType);
321 
326  itkSetMacro(MaximumNumberOfIterations, VariableSizeArrayType);
327 
332  itkGetConstMacro(MaximumNumberOfIterations, VariableSizeArrayType);
333 
341  itkSetMacro(ConvergenceThreshold, RealType);
342 
350  itkGetConstMacro(ConvergenceThreshold, RealType);
351 
361  itkGetConstObjectMacro(LogBiasFieldControlPointLattice, BiasFieldControlPointLatticeType);
362 
367  itkGetConstMacro(ElapsedIterations, unsigned int);
368 
373  itkGetConstMacro(CurrentConvergenceMeasurement, RealType);
374 
379  itkGetConstMacro(CurrentLevel, unsigned int);
380 
386 
387 protected:
389  ~N4BiasFieldCorrectionImageFilter() override = default;
390  void
391  PrintSelf(std::ostream & os, Indent indent) const override;
392 
393  void
394  GenerateData() override;
395 
396 private:
397  // N4 algorithm functions: The basic algorithm iterates between sharpening
398  // the intensity histogram of the corrected input image and spatially
399  // smoothing those results with a B-spline scalar field estimate of the
400  // bias field. The former is handled by the function SharpenImage()
401  // whereas the latter is handled by the function UpdateBiasFieldEstimate().
402  // Convergence is determined by the coefficient of variation of the difference
403  // image between the current bias field estimate and the previous estimate.
404 
410  void
411  SharpenImage(const RealImageType * unsharpenedImage, RealImageType * sharpenedImage) const;
412 
419  UpdateBiasFieldEstimate(RealImageType *, std::size_t);
420 
425  RealType
427 
429  bool m_UseMaskLabel{ false };
430 
431  // Parameters for deconvolution with Wiener filter
432 
433  unsigned int m_NumberOfHistogramBins{ 200 };
434  RealType m_WienerFilterNoise{ static_cast<RealType>(0.01) };
435  RealType m_BiasFieldFullWidthAtHalfMaximum{ static_cast<RealType>(0.15) };
436 
437  // Convergence parameters
438 
440  unsigned int m_ElapsedIterations{ 0 };
441  RealType m_ConvergenceThreshold{ static_cast<RealType>(0.001) };
443  unsigned int m_CurrentLevel{ 0 };
444 
445  // B-spline fitting parameters
446 
447  typename BiasFieldControlPointLatticeType::Pointer m_LogBiasFieldControlPointLattice;
448 
449  unsigned int m_SplineOrder{ 3 };
452 };
453 
454 } // end namespace itk
455 
456 #ifndef ITK_MANUAL_INSTANTIATION
457 # include "itkN4BiasFieldCorrectionImageFilter.hxx"
458 #endif
459 
460 #endif
itk::N4BiasFieldCorrectionImageFilter
Implementation of the N4 bias field correction algorithm.
Definition: itkN4BiasFieldCorrectionImageFilter.h:96
itk::N4BiasFieldCorrectionImageFilter::SetConfidenceImage
virtual void SetConfidenceImage(const RealImageType *input)
itk::N4BiasFieldCorrectionImageFilter::PrintSelf
void PrintSelf(std::ostream &os, Indent indent) const override
itk::N4BiasFieldCorrectionImageFilter::RealImagePointer
typename RealImageType::Pointer RealImagePointer
Definition: itkN4BiasFieldCorrectionImageFilter.h:124
itk::PointSet
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition: itkPointSet.h:82
itk::N4BiasFieldCorrectionImageFilter::m_LogBiasFieldControlPointLattice
BiasFieldControlPointLatticeType::Pointer m_LogBiasFieldControlPointLattice
Definition: itkN4BiasFieldCorrectionImageFilter.h:447
itk::N4BiasFieldCorrectionImageFilter::m_MaximumNumberOfIterations
VariableSizeArrayType m_MaximumNumberOfIterations
Definition: itkN4BiasFieldCorrectionImageFilter.h:439
itk::N4BiasFieldCorrectionImageFilter::N4BiasFieldCorrectionImageFilter
N4BiasFieldCorrectionImageFilter()
itk::BSplineScatteredDataPointSetToImageFilter
Image filter which provides a B-spline output approximation.
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:132
itk::N4BiasFieldCorrectionImageFilter::MaskPixelType
typename MaskImageType::PixelType MaskPixelType
Definition: itkN4BiasFieldCorrectionImageFilter.h:120
itk::Vector
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:62
itk::N4BiasFieldCorrectionImageFilter::RealImageType
Image< RealType, ImageDimension > RealImageType
Definition: itkN4BiasFieldCorrectionImageFilter.h:123
itk::N4BiasFieldCorrectionImageFilter::ImageDimension
static constexpr unsigned int ImageDimension
Definition: itkN4BiasFieldCorrectionImageFilter.h:114
itk::N4BiasFieldCorrectionImageFilter::m_NumberOfControlPoints
ArrayType m_NumberOfControlPoints
Definition: itkN4BiasFieldCorrectionImageFilter.h:450
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::Image::Pointer
SmartPointer< Self > Pointer
Definition: itkImage.h:94
itk::N4BiasFieldCorrectionImageFilter::SetInput1
void SetInput1(const InputImageType *image)
Definition: itkN4BiasFieldCorrectionImageFilter.h:147
itk::N4BiasFieldCorrectionImageFilter::m_CurrentConvergenceMeasurement
RealType m_CurrentConvergenceMeasurement
Definition: itkN4BiasFieldCorrectionImageFilter.h:442
itk::N4BiasFieldCorrectionImageFilter::m_MaskLabel
MaskPixelType m_MaskLabel
Definition: itkN4BiasFieldCorrectionImageFilter.h:428
itk::BSplineScatteredDataPointSetToImageFilter::PointDataImageType
Image< PointDataType, Self::ImageDimension > PointDataImageType
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:174
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::N4BiasFieldCorrectionImageFilter::PointType
typename PointSetType::PointType PointType
Definition: itkN4BiasFieldCorrectionImageFilter.h:132
itk::N4BiasFieldCorrectionImageFilter::m_NumberOfFittingLevels
ArrayType m_NumberOfFittingLevels
Definition: itkN4BiasFieldCorrectionImageFilter.h:451
itk::N4BiasFieldCorrectionImageFilter::GenerateData
void GenerateData() override
itk::N4BiasFieldCorrectionImageFilter::~N4BiasFieldCorrectionImageFilter
~N4BiasFieldCorrectionImageFilter() override=default
itk::N4BiasFieldCorrectionImageFilter::m_UseMaskLabel
bool m_UseMaskLabel
Definition: itkN4BiasFieldCorrectionImageFilter.h:429
itk::N4BiasFieldCorrectionImageFilter::SetNumberOfFittingLevels
void SetNumberOfFittingLevels(unsigned int n)
Definition: itkN4BiasFieldCorrectionImageFilter.h:306
itk::N4BiasFieldCorrectionImageFilter::VariableSizeArrayType
Array< unsigned int > VariableSizeArrayType
Definition: itkN4BiasFieldCorrectionImageFilter.h:125
itk::N4BiasFieldCorrectionImageFilter::EnlargeOutputRequestedRegion
void EnlargeOutputRequestedRegion(DataObject *) override
itk::N4BiasFieldCorrectionImageFilter::m_ConvergenceThreshold
RealType m_ConvergenceThreshold
Definition: itkN4BiasFieldCorrectionImageFilter.h:441
itk::N4BiasFieldCorrectionImageFilter::BiasFieldControlPointLatticeType
typename BSplineFilterType::PointDataImageType BiasFieldControlPointLatticeType
Definition: itkN4BiasFieldCorrectionImageFilter.h:136
itk::BSplineScatteredDataPointSetToImageFilter::ArrayType
FixedArray< unsigned, Self::ImageDimension > ArrayType
Definition: itkBSplineScatteredDataPointSetToImageFilter.h:178
itk::PointSet::Pointer
SmartPointer< Self > Pointer
Definition: itkPointSet.h:90
itkBSplineScatteredDataPointSetToImageFilter.h
itk::N4BiasFieldCorrectionImageFilter::m_NumberOfHistogramBins
unsigned int m_NumberOfHistogramBins
Definition: itkN4BiasFieldCorrectionImageFilter.h:433
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itkImageToImageFilter.h
itk::N4BiasFieldCorrectionImageFilter::UpdateBiasFieldEstimate
RealImagePointer UpdateBiasFieldEstimate(RealImageType *, std::vcl_size_t)
itk::N4BiasFieldCorrectionImageFilter::SharpenImage
void SharpenImage(const RealImageType *unsharpenedImage, RealImageType *sharpenedImage) const
itk::ImageToImageFilter::SetInput
virtual void SetInput(const InputImageType *input)
itkArray.h
itk::N4BiasFieldCorrectionImageFilter::m_CurrentLevel
unsigned int m_CurrentLevel
Definition: itkN4BiasFieldCorrectionImageFilter.h:443
itk::N4BiasFieldCorrectionImageFilter::m_WienerFilterNoise
RealType m_WienerFilterNoise
Definition: itkN4BiasFieldCorrectionImageFilter.h:434
itk::N4BiasFieldCorrectionImageFilter::SetInput2
void SetInput2(const MaskImageType *mask)
Definition: itkN4BiasFieldCorrectionImageFilter.h:159
itk::N4BiasFieldCorrectionImageFilter::ReconstructBiasField
RealImagePointer ReconstructBiasField(const BiasFieldControlPointLatticeType *)
itk::N4BiasFieldCorrectionImageFilter::SetNumberOfFittingLevels
virtual void SetNumberOfFittingLevels(ArrayType _arg)
itk::N4BiasFieldCorrectionImageFilter::PointSetPointer
typename PointSetType::Pointer PointSetPointer
Definition: itkN4BiasFieldCorrectionImageFilter.h:131
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::N4BiasFieldCorrectionImageFilter::SetMaskImage
virtual void SetMaskImage(const MaskImageType *input)
itk::N4BiasFieldCorrectionImageFilter::m_SplineOrder
unsigned int m_SplineOrder
Definition: itkN4BiasFieldCorrectionImageFilter.h:449
itk::N4BiasFieldCorrectionImageFilter::RealType
float RealType
Definition: itkN4BiasFieldCorrectionImageFilter.h:122
itk::N4BiasFieldCorrectionImageFilter::m_BiasFieldFullWidthAtHalfMaximum
RealType m_BiasFieldFullWidthAtHalfMaximum
Definition: itkN4BiasFieldCorrectionImageFilter.h:435
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:138
itkVector.h
itk::N4BiasFieldCorrectionImageFilter::SetInput3
void SetInput3(const RealImageType *image)
Definition: itkN4BiasFieldCorrectionImageFilter.h:204
itk::N4BiasFieldCorrectionImageFilter::ArrayType
typename BSplineFilterType::ArrayType ArrayType
Definition: itkN4BiasFieldCorrectionImageFilter.h:137
itk::Array< unsigned int >
itk::N4BiasFieldCorrectionImageFilter::CalculateConvergenceMeasurement
RealType CalculateConvergenceMeasurement(const RealImageType *, const RealImageType *) const
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86
itk::N4BiasFieldCorrectionImageFilter::m_ElapsedIterations
unsigned int m_ElapsedIterations
Definition: itkN4BiasFieldCorrectionImageFilter.h:440
itk::PointSet::PointType
typename MeshTraits::PointType PointType
Definition: itkPointSet.h:106
itkPointSet.h
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293
itk::N4BiasFieldCorrectionImageFilter::MaskImageType
TMaskImage MaskImageType
Definition: itkN4BiasFieldCorrectionImageFilter.h:119