ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkMRIBiasFieldCorrectionFilter.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 * https://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
24#include "itkArray2D.h"
31
32#include <memory> // For unique_ptr.
33
34
35namespace itk
36{
50template <typename TImage, typename TImageMask, typename TBiasField>
51class ITK_TEMPLATE_EXPORT MRIBiasEnergyFunction : public SingleValuedCostFunction
52{
53public:
54 ITK_DISALLOW_COPY_AND_MOVE(MRIBiasEnergyFunction);
55
61
63 itkOverrideGetNameOfClassMacro(MRIBiasEnergyFunction);
64
66 itkNewMacro(Self);
67
69 using ImageType = TImage;
70 using ImagePointer = typename ImageType::Pointer;
71 using ImageElementType = typename ImageType::PixelType;
72 using ImageIndexType = typename ImageType::IndexType;
73 using ImageRegionType = typename ImageType::RegionType;
74 using MaskType = TImageMask;
75 using MaskPointer = typename MaskType::Pointer;
76 using MaskElementType = typename MaskType::PixelType;
77
79 using BiasFieldType = TBiasField;
80
83 using typename Superclass::ParametersType;
84
87
90
91 static constexpr unsigned int SpaceDimension = 3;
92
95
97 using SamplingFactorType = unsigned int[SpaceDimension];
98
100 itkSetObjectMacro(Image, ImageType);
101
103 itkSetObjectMacro(Mask, MaskType);
104
106 itkSetMacro(Region, ImageRegionType);
107
109 void
111 {
112 m_BiasField = bias;
113 }
114
117 void
119 {
120 for (unsigned int i = 0; i < SpaceDimension; ++i)
121 {
122 m_SamplingFactor[i] = factor[i];
123 }
124 }
125
128 double
129 GetEnergy0(double diff)
130 {
131 return (*m_InternalEnergyFunction)(diff);
132 }
133
136 MeasureType
137 GetValue(const ParametersType & parameters) const override;
138
141 void
142 GetDerivative(const ParametersType & itkNotUsed(parameters), DerivativeType & itkNotUsed(derivative)) const override
143 {}
144
149 void
151
152 unsigned int
153 GetNumberOfParameters() const override;
154
155protected:
158
160 ~MRIBiasEnergyFunction() override = default;
161
162private:
165
168
171
174
176 std::unique_ptr<InternalEnergyFunction> m_InternalEnergyFunction;
177
180}; // end of class
181
225template <typename TInputImage, typename TOutputImage, typename TMaskImage>
226class ITK_TEMPLATE_EXPORT MRIBiasFieldCorrectionFilter : public ImageToImageFilter<TInputImage, TOutputImage>
227{
228public:
229 ITK_DISALLOW_COPY_AND_MOVE(MRIBiasFieldCorrectionFilter);
230
236
238 itkNewMacro(Self);
239
241 itkOverrideGetNameOfClassMacro(MRIBiasFieldCorrectionFilter);
242
244 static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
245
247 using OutputImageType = TOutputImage;
248 using OutputImagePointer = typename TOutputImage::Pointer;
249 using OutputImageIndexType = typename TOutputImage::IndexType;
250 using OutputImagePixelType = typename TOutputImage::PixelType;
251 using OutputImageSizeType = typename TOutputImage::SizeType;
252 using OutputImageRegionType = typename TOutputImage::RegionType;
253
254 using InputImageType = TInputImage;
255 using InputImagePointer = typename TInputImage::Pointer;
256 using InputImageIndexType = typename TInputImage::IndexType;
257 using InputImagePixelType = typename TInputImage::PixelType;
258 using InputImageSizeType = typename TInputImage::SizeType;
259 using InputImageRegionType = typename TInputImage::RegionType;
260
262 using ImageMaskType = TMaskImage;
263 using ImageMaskPointer = typename ImageMaskType::Pointer;
264 using ImageMaskRegionType = typename ImageMaskType::RegionType;
265
271
275 using SlabRegionVectorIteratorType = typename SlabRegionVectorType::iterator;
276
279
283
286
289
292
297 void
299 itkGetModifiableObjectMacro(InputMask, ImageMaskType);
301
305 void
307 itkGetModifiableObjectMacro(OutputMask, ImageMaskType);
309
310#if defined(ITK_LEGACY_REMOVE)
314 void
316 {
318 }
319#endif // defined ( ITK_LEGACY_REMOVE )
320
325 itkSetMacro(BiasFieldMultiplicative, bool);
326 itkGetConstMacro(BiasFieldMultiplicative, bool);
327 itkBooleanMacro(BiasFieldMultiplicative);
329
330#if defined(ITK_LEGACY_REMOVE)
332 bool
337#endif // defined ( ITK_LEGACY_REMOVE )
338
344 itkSetMacro(UsingInterSliceIntensityCorrection, bool);
345 itkGetConstMacro(UsingInterSliceIntensityCorrection, bool);
346 itkBooleanMacro(UsingInterSliceIntensityCorrection);
348
355 itkSetMacro(UsingSlabIdentification, bool);
356 itkGetConstMacro(UsingSlabIdentification, bool);
357 itkBooleanMacro(UsingSlabIdentification);
359
360 itkSetMacro(SlabBackgroundMinimumThreshold, InputImagePixelType);
361 itkGetConstReferenceMacro(SlabBackgroundMinimumThreshold, InputImagePixelType);
362
363 itkSetMacro(SlabNumberOfSamples, unsigned int);
364 itkGetConstReferenceMacro(SlabNumberOfSamples, unsigned int);
365
366 itkSetMacro(SlabTolerance, double);
367 itkGetConstReferenceMacro(SlabTolerance, double);
368
375 itkSetMacro(UsingBiasFieldCorrection, bool);
376 itkGetConstMacro(UsingBiasFieldCorrection, bool);
377 itkBooleanMacro(UsingBiasFieldCorrection);
379
383 itkSetMacro(GeneratingOutput, bool);
384 itkGetConstMacro(GeneratingOutput, bool);
385 itkBooleanMacro(GeneratingOutput);
387
391 itkSetMacro(SlicingDirection, int);
392 itkGetConstMacro(SlicingDirection, int);
394
397 itkSetMacro(BiasFieldDegree, int);
398 itkGetConstMacro(BiasFieldDegree, int);
400
403 void
405 {
406 this->Modified();
407 m_BiasFieldCoefficients = coefficients;
408 }
409
413 BiasFieldType::CoefficientArrayType
418
422 void
424
428 itkSetMacro(VolumeCorrectionMaximumIteration, int);
429 itkGetConstMacro(VolumeCorrectionMaximumIteration, int);
431
435 itkSetMacro(InterSliceCorrectionMaximumIteration, int);
436 itkGetConstMacro(InterSliceCorrectionMaximumIteration, int);
438
441 itkSetMacro(OptimizerInitialRadius, double);
442 itkGetConstMacro(OptimizerInitialRadius, double);
444
447 itkSetMacro(OptimizerGrowthFactor, double);
448 itkGetConstMacro(OptimizerGrowthFactor, double);
450
453 itkSetMacro(OptimizerShrinkFactor, double);
454 itkGetConstMacro(OptimizerShrinkFactor, double);
456
463 void
464 SetNumberOfLevels(unsigned int num);
465
467 itkGetConstMacro(NumberOfLevels, unsigned int);
468
475 void
476 SetSchedule(const ScheduleType & schedule);
477
479 itkGetConstReferenceMacro(Schedule, ScheduleType);
480
485 void
486 SetStartingShrinkFactors(unsigned int factor);
487
488 void
489 SetStartingShrinkFactors(const unsigned int * factors);
490
492 const unsigned int *
494
498 static bool
500
507 void
509
513 EstimateBiasField(InputImageRegionType region, unsigned int degree, int maximumIteration);
514
518 void
520
523 void
525
526protected:
528 ~MRIBiasFieldCorrectionFilter() override = default;
529 void
530 PrintSelf(std::ostream & os, Indent indent) const override;
531
534 bool
536
537protected:
541 void
543
546 void
548
551 template <typename TSource, typename TTarget>
552 void
553 CopyAndConvertImage(const TSource * source, TTarget * target, typename TTarget::RegionType requestedRegion)
554 {
555 using SourceIterator = ImageRegionConstIterator<TSource>;
556 using TargetIterator = ImageRegionIterator<TTarget>;
557 using TargetPixelType = typename TTarget::PixelType;
558
559 SourceIterator s_iter(source, requestedRegion);
560 TargetIterator t_iter(target, requestedRegion);
561
562 s_iter.GoToBegin();
563 t_iter.GoToBegin();
564 while (!s_iter.IsAtEnd())
565 {
566 t_iter.Set(static_cast<TargetPixelType>(s_iter.Get()));
567 ++s_iter;
568 ++t_iter;
569 }
570 }
571
576 void
578
582 void
584
585 void
586 GenerateData() override;
587
588private:
591
594
597
600
603
606
609
612
617 bool m_GeneratingOutput{ true };
618
619 unsigned int m_SlabNumberOfSamples{ 200 };
621 double m_SlabTolerance{ 0.0 };
622
625
627 unsigned int m_NumberOfLevels{ 0 };
628
631
635
639
641
643
646
649
652
655
658};
659
660} // end namespace itk
661
662#ifndef ITK_MANUAL_INSTANTIATION
663# include "itkMRIBiasFieldCorrectionFilter.hxx"
664#endif
665
666#endif
Array2D class representing a 2D array.
Definition itkArray2D.h:43
Array class with size defined at construction time.
Definition itkArray.h:48
Multiple valley shaped curve function.
A multi-dimensional iterator templated over image type that walks a region of pixels.
A multi-dimensional iterator templated over image type that walks a region of pixels.
Templated n-dimensional image class.
Definition itkImage.h:89
ImageRegion< VImageDimension > RegionType
Control indentation during Print() invocation.
Definition itkIndent.h:50
Identifies slabs in MR images comparing minimum intensity averages.
std::vector< ImageRegionType > SlabRegionVectorType
Represents a cost function for MRI bias field correction optimization.
typename ImageType::IndexType ImageIndexType
void SetSamplingFactors(const SamplingFactorType factor)
Superclass::DerivativeType DerivativeType
void GetDerivative(const ParametersType &parameters, DerivativeType &derivative) const override
unsigned int[SpaceDimension] SamplingFactorType
~MRIBiasEnergyFunction() override=default
unsigned int GetNumberOfParameters() const override
void InitializeDistributions(Array< double > classMeans, Array< double > classSigmas)
typename ImageType::PixelType ImageElementType
typename ImageType::RegionType ImageRegionType
typename MaskType::PixelType MaskElementType
MeasureType GetValue(const ParametersType &parameters) const override
typename ImageType::Pointer ImagePointer
void SetNumberOfLevels(unsigned int num)
NormalVariateGeneratorType::Pointer m_NormalVariateGenerator
ImageToImageFilter< TInputImage, TOutputImage > Superclass
typename TOutputImage::IndexType OutputImageIndexType
void SetStartingShrinkFactors(unsigned int factor)
Statistics::NormalVariateGenerator NormalVariateGeneratorType
MRASlabIdentifier< InputImageType > MRASlabIdentifierType
typename TInputImage::IndexType InputImageIndexType
void SetInitialBiasFieldCoefficients(const BiasFieldType::CoefficientArrayType &coefficients)
void CorrectImage(BiasFieldType &bias, InputImageRegionType region)
void GetBiasFieldSize(InputImageRegionType region, BiasFieldType::DomainSizeType &biasSize)
static bool IsScheduleDownwardDivisible(const ScheduleType &schedule)
const unsigned int * GetStartingShrinkFactors() const
bool CheckMaskImage(ImageMaskType *mask)
void SetOutputMask(ImageMaskType *outputMask)
typename InternalImageType::PixelType InternalImagePixelType
Image< float, Self::ImageDimension > InternalImageType
BiasFieldType EstimateBiasField(InputImageRegionType region, unsigned int degree, int maximumIteration)
typename ImageMaskType::Pointer ImageMaskPointer
MRIBiasEnergyFunction< InternalImageType, ImageMaskType, BiasFieldType > EnergyFunctionType
void PrintSelf(std::ostream &os, Indent indent) const override
typename TOutputImage::PixelType OutputImagePixelType
void SetSchedule(const ScheduleType &schedule)
BiasFieldType::CoefficientArrayType GetEstimatedBiasFieldCoefficients()
typename MRASlabIdentifierType::SlabRegionVectorType SlabRegionVectorType
typename ImageMaskType::RegionType ImageMaskRegionType
typename TInputImage::SizeType InputImageSizeType
void CorrectInterSliceIntensityInhomogeneity(InputImageRegionType region)
typename TOutputImage::SizeType OutputImageSizeType
void ExpImage(InternalImageType *source, InternalImageType *target)
typename InternalImageType::Pointer InternalImagePointer
typename TInputImage::PixelType InputImagePixelType
typename InternalImageType::RegionType InternalImageRegionType
void AdjustSlabRegions(SlabRegionVectorType &slabs, OutputImageRegionType requestedRegion)
typename TOutputImage::Pointer OutputImagePointer
typename SlabRegionVectorType::iterator SlabRegionVectorIteratorType
typename TOutputImage::RegionType OutputImageRegionType
BiasFieldType::CoefficientArrayType m_EstimatedBiasFieldCoefficients
void SetStartingShrinkFactors(const unsigned int *factors)
BiasFieldType::CoefficientArrayType m_BiasFieldCoefficients
void CopyAndConvertImage(const TSource *source, TTarget *target, typename TTarget::RegionType requestedRegion)
~MRIBiasFieldCorrectionFilter() override=default
typename EnergyFunctionType::Pointer EnergyFunctionPointer
void Log1PImage(InternalImageType *source, InternalImageType *target)
void SetInputMask(ImageMaskType *inputMask)
void SetTissueClassStatistics(const Array< double > &means, const Array< double > &sigmas)
typename TInputImage::RegionType InputImageRegionType
2D and 3D multivariate Legendre Polynomial
virtual void Modified() const
A region represents some portion or piece of data.
Definition itkRegion.h:66
Array< ParametersValueType > DerivativeType
Superclass::ParametersType ParametersType
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....