ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkPatchBasedDenoisingImageFilter.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 itkPatchBasedDenoisingImageFilter_h
19 #define itkPatchBasedDenoisingImageFilter_h
20 
22 #include "itkImageRegionIterator.h"
24 #include "itkVector.h"
25 #include "itkVectorImage.h"
26 #include "itkRGBPixel.h"
27 #include "itkRGBAPixel.h"
28 #include "itkDiffusionTensor3D.h"
29 #include "itkFixedArray.h"
30 #include "itkMatrix.h"
32 #include "itkEnableIf.h"
33 #include "itkIsSame.h"
34 
35 #include <vector>
36 
37 namespace itk
38 {
60 template <typename TInputImage, typename TOutputImage>
61 class ITK_TEMPLATE_EXPORT PatchBasedDenoisingImageFilter :
62  public PatchBasedDenoisingBaseImageFilter<TInputImage, TOutputImage>
63 {
64 public:
70  typedef typename Superclass::OutputImagePointer OutputImagePointer;
71 
73  itkNewMacro(Self);
74 
76  itkTypeMacro(PatchBasedDenoisingImageFilter,
78 
82 
84  itkStaticConstMacro(ImageDimension, unsigned int,
85  Superclass::ImageDimension);
86 
88  typedef typename InputImageType::RegionType InputImageRegionType;
89 
93 
96  typedef typename Superclass::PixelType PixelType;
97  typedef typename Superclass::PixelValueType PixelValueType;
98 
104 
107  typedef typename Superclass::PatchRadiusType PatchRadiusType;
108  typedef typename Superclass::InputImagePatchIterator InputImagePatchIterator;
110  typedef typename Superclass::PatchWeightsType PatchWeightsType;
111 
117 
125  typedef std::vector<EigenValuesArrayType> EigenValuesCacheType;
126  typedef std::vector<EigenVectorsMatrixType> EigenVectorsCacheType;
127 
129  {
139  };
140 
145  itkSetMacro(UseSmoothDiscPatchWeights, bool);
146  itkBooleanMacro(UseSmoothDiscPatchWeights);
147  itkGetConstMacro(UseSmoothDiscPatchWeights, bool);
149 
154  void SetKernelBandwidthSigma(const RealArrayType& kernelSigma);
155  itkGetConstMacro(KernelBandwidthSigma, RealArrayType);
157 
162  itkSetClampMacro(KernelBandwidthFractionPixelsForEstimation, double, 0.01, 1.0);
163  itkGetConstReferenceMacro(KernelBandwidthFractionPixelsForEstimation, double);
165 
168  itkSetMacro(ComputeConditionalDerivatives, bool);
169  itkBooleanMacro(ComputeConditionalDerivatives);
170  itkGetConstMacro(ComputeConditionalDerivatives, bool);
172 
186  itkSetMacro(UseFastTensorComputations, bool);
187  itkBooleanMacro(UseFastTensorComputations);
188  itkGetConstMacro(UseFastTensorComputations, bool);
190 
192  itkStaticConstMacro(MaxSigmaUpdateIterations, unsigned int,
193  20);
194 
201  itkSetClampMacro(KernelBandwidthMultiplicationFactor, double, 0.01, 100);
202  itkGetConstReferenceMacro(KernelBandwidthMultiplicationFactor, double);
204 
208  void SetNoiseSigma(const RealType& sigma);
209 
210  itkGetConstMacro(NoiseSigma, RealType);
211 
213  itkSetObjectMacro(Sampler, BaseSamplerType);
214  itkGetModifiableObjectMacro(Sampler, BaseSamplerType);
216 
218  itkGetConstMacro(NumIndependentComponents, unsigned int);
219 
220 protected:
222  ~PatchBasedDenoisingImageFilter() ITK_OVERRIDE;
223  virtual void PrintSelf(std::ostream& os, Indent indent) const ITK_OVERRIDE;
224 
226  virtual void EmptyCaches();
227 
229  virtual void AllocateUpdateBuffer() ITK_OVERRIDE;
230 
231  virtual void CopyInputToOutput() ITK_OVERRIDE;
232 
233  virtual void GenerateInputRequestedRegion() ITK_OVERRIDE;
234 
241  template< typename T>
242  typename EnableIfC<
243  IsSame<T, typename NumericTraits<T>::ValueType>::Value,
244  T >::Type
245  GetComponent(const T pix,
246  unsigned int itkNotUsed( idx ) ) const
247  {
248  // The enable if idiom is used to overload this method for both
249  // scalars and multi-component types. By exploiting that
250  // NumericTraits' ValueType typedef (defines the per-element type
251  // for multi-component types ) is different then the parameterize
252  // type, the bracket operator is used only for multi-component
253  // types.
254  return pix;
255  }
256 
257  template< typename T>
258  typename DisableIfC<
259  IsSame<T, typename NumericTraits<T>::ValueType>::Value,
260  typename NumericTraits<T>::ValueType>::Type
261  GetComponent(const T& pix,
262  unsigned int idx ) const
263  {
264  return pix[idx];
265  }
266 
268  template< typename T >
269  void
270  SetComponent( T &pix,
271  unsigned int itkNotUsed( idx ),
272  typename EnableIfC< IsSame<T,
273  typename NumericTraits<T>::ValueType>::Value, RealValueType>::Type val) const
274  {
275  pix = val;
276  }
277 
278  template< typename T >
279  void
280  SetComponent( T &pix,
281  unsigned int idx,
282  typename DisableIfC< IsSame<T,
283  typename NumericTraits<T>::ValueType>::Value,
284  RealValueType>::Type val) const
285  {
286  pix[idx] = val;
287  }
288 
292  void ComputeMinMax(const Image< DiffusionTensor3D<PixelValueType> , ImageDimension>* img)
293  {
294  if( this->GetComponentSpace() == Superclass::RIEMANNIAN )
295  {
296  DispatchedRiemannianMinMax(img);
297  }
298  else
299  {
300  DispatchedArrayMinMax(img);
301  }
302  }
303 
304  template< typename TImageType>
305  typename EnableIfC<
306  IsSame<typename TImageType::PixelType, typename NumericTraits<typename TImageType::PixelType>::ValueType>::Value
307  >::Type
308  ComputeMinMax(const TImageType* img)
309  {
310  DispatchedMinMax(img);
311  }
312 
313  template< typename TImageType>
314  typename DisableIfC<
315  IsSame<typename TImageType::PixelType, typename NumericTraits<typename TImageType::PixelType>::ValueType>::Value
316  >::Type
317  ComputeMinMax(const TImageType* img)
318  {
319  DispatchedArrayMinMax(img);
320  }
321 
332  const RealArrayType& weight,
333  bool useCachedComputations,
334  SizeValueType cacheIndex,
335  EigenValuesCacheType& eigenValsCache,
336  EigenVectorsCacheType& eigenVecsCache,
337  RealType& diff, RealArrayType& norm)
338  {
339  if( this->GetComponentSpace() == Superclass::RIEMANNIAN )
340  {
341  ComputeLogMapAndWeightedSquaredGeodesicDifference(a, b, weight,
342  useCachedComputations, cacheIndex,
343  eigenValsCache, eigenVecsCache,
344  diff, norm);
345  }
346  else
347  {
348  ComputeSignedEuclideanDifferenceAndWeightedSquaredNorm(a, b, weight,
349  useCachedComputations, cacheIndex,
350  eigenValsCache, eigenVecsCache,
351  diff, norm);
352  }
353  }
355 
356  template <typename PixelT>
358  const PixelT& b,
359  const RealArrayType& weight,
360  bool useCachedComputations,
361  SizeValueType cacheIndex,
362  EigenValuesCacheType& eigenValsCache,
363  EigenVectorsCacheType& eigenVecsCache,
364  RealType& diff, RealArrayType& norm)
365  {
366  ComputeSignedEuclideanDifferenceAndWeightedSquaredNorm(a, b, weight,
367  useCachedComputations, cacheIndex,
368  eigenValsCache, eigenVecsCache,
369  diff, norm);
370  }
371 
376  const RealType& b)
377  {
378  if( this->GetComponentSpace() == Superclass::RIEMANNIAN )
379  {
380  return this->AddExponentialMapUpdate(a, b);
381  }
382  else
383  {
384  return this->AddEuclideanUpdate(a, b);
385  }
386  }
388 
389  template <typename RealT>
390  RealType AddUpdate(const RealT& a,
391  const RealType& b)
392  {
393  return this->AddEuclideanUpdate(a, b);
394  }
395 
396  virtual void EnforceConstraints();
397 
398  virtual void Initialize() ITK_OVERRIDE;
399 
400  virtual void InitializeKernelSigma();
401 
402  virtual void InitializePatchWeights() ITK_OVERRIDE;
403 
404  virtual void InitializePatchWeightsSmoothDisc();
405 
406  virtual void InitializeIteration() ITK_OVERRIDE;
407 
408  virtual void ComputeKernelBandwidthUpdate() ITK_OVERRIDE; // derived from base class;
409 
410  // define here
411 
412  virtual ThreadDataStruct ThreadedComputeSigmaUpdate(const InputImageRegionType& regionToProcess,
413  const int itkNotUsed(threadId),
414  ThreadDataStruct threadData);
415 
416  virtual RealArrayType ResolveSigmaUpdate();
417 
418  virtual void ComputeImageUpdate() ITK_OVERRIDE;
419 
420  virtual ThreadDataStruct ThreadedComputeImageUpdate(const InputImageRegionType& regionToProcess,
421  const int threadId,
422  ThreadDataStruct threadData);
423 
424  virtual RealType ComputeGradientJointEntropy(InstanceIdentifier id,
425  typename ListAdaptorType::Pointer& inList,
426  BaseSamplerPointer& sampler,
427  ThreadDataStruct& threadData);
428 
429  virtual void ApplyUpdate() ITK_OVERRIDE;
430 
431  virtual void ThreadedApplyUpdate(const InputImageRegionType& regionToProcess,
432  const int itkNotUsed(threadId) );
433 
434  virtual void PostProcessOutput() ITK_OVERRIDE;
435 
436  virtual void SetThreadData(int threadId, const ThreadDataStruct& data);
437 
438  virtual ThreadDataStruct GetThreadData(int threadId);
439 
440 private:
441  ITK_DISALLOW_COPY_AND_ASSIGN(PatchBasedDenoisingImageFilter);
442 
445  static ITK_THREAD_RETURN_TYPE ComputeSigmaUpdateThreaderCallback( void *arg );
446 
449  static ITK_THREAD_RETURN_TYPE ComputeImageUpdateThreaderCallback( void *arg );
450 
453  static ITK_THREAD_RETURN_TYPE ApplyUpdateThreaderCallback( void *arg );
454 
455  template <typename TInputImageType>
456  void DispatchedMinMax(const TInputImageType* img);
457 
458  template <typename TInputImageType>
459  void DispatchedArrayMinMax(const TInputImageType* img);
460 
461  template <typename TInputImageType>
462  void DispatchedVectorMinMax(const TInputImageType* img);
463 
464  template <typename TInputImageType>
465  void DispatchedRiemannianMinMax(const TInputImageType* img);
466 
469  static ITK_THREAD_RETURN_TYPE RiemannianMinMaxThreaderCallback( void *arg );
470 
471  ThreadDataStruct ThreadedRiemannianMinMax(const InputImageRegionType& regionToProcess,
472  const int itkNotUsed(threadId),
473  const InputImageType* img,
474  ThreadDataStruct threadData);
475 
476  virtual void ResolveRiemannianMinMax();
477 
478  void ComputeSignedEuclideanDifferenceAndWeightedSquaredNorm(const PixelType& a, const PixelType& b,
479  const RealArrayType& weight,
480  bool useCachedComputations,
481  SizeValueType cacheIndex,
482  EigenValuesCacheType& eigenValsCache,
483  EigenVectorsCacheType& eigenVecsCache,
484  RealType& diff, RealArrayType& norm);
485 
487  void ComputeLogMapAndWeightedSquaredGeodesicDifference(const DiffusionTensor3D<PixelValueType>& spdMatrixA,
488  const DiffusionTensor3D<PixelValueType>& spdMatrixB,
489  const RealArrayType& weight,
490  bool useCachedComputations,
491  SizeValueType cacheIndex,
492  EigenValuesCacheType& eigenValsCache,
493  EigenVectorsCacheType& eigenVecsCache,
494  RealType& symMatrixLogMap, RealArrayType& geodesicDist);
495 
496  template <typename TensorValueT>
497  void Compute3x3EigenAnalysis(const DiffusionTensor3D<TensorValueT>& spdMatrix,
498  FixedArray< TensorValueT, 3 >& eigenVals,
499  Matrix< TensorValueT, 3, 3 >& eigenVecs);
500 
501  RealType AddEuclideanUpdate(const RealType& a, const RealType& b);
502 
504  RealType AddExponentialMapUpdate(const DiffusionTensor3D<RealValueType>& spdMatrix,
505  const DiffusionTensor3D<RealValueType>& symMatrix);
506 
508  {
509  PatchBasedDenoisingImageFilter *Filter;
511  };
512 
513  std::vector<ThreadDataStruct> m_ThreadData;
514 
516  typename OutputImageType::Pointer m_UpdateBuffer;
517 
518  unsigned int m_NumPixelComponents;
520  unsigned int m_TotalNumberPixels;
521 
523 
525 
534  double m_MinSigma;
540 
544 
546  typename ListAdaptorType::Pointer m_SearchSpaceList;
547 };
548 } // end namespace itk
549 
550 #ifndef ITK_MANUAL_INSTANTIATION
551 #include "itkPatchBasedDenoisingImageFilter.hxx"
552 #endif
553 
554 #endif
A templated class holding a M x N size Matrix.
Definition: itkMatrix.h:53
Derived class implementing a specific patch-based denoising algorithm, as detailed below...
BaseSamplerType::InstanceIdentifier InstanceIdentifier
void ComputeDifferenceAndWeightedSquaredNorm(const PixelT &a, const PixelT &b, const RealArrayType &weight, bool useCachedComputations, SizeValueType cacheIndex, EigenValuesCacheType &eigenValsCache, EigenVectorsCacheType &eigenVecsCache, RealType &diff, RealArrayType &norm)
Superclass::InputImagePatchIterator InputImagePatchIterator
Base class for patch-based denoising algorithms.
DisableIfC< IsSame< T, typename NumericTraits< T >::ValueType >::Value, typename NumericTraits< T >::ValueType >::Type GetComponent(const T &pix, unsigned int idx) const
std::vector< EigenVectorsMatrixType > EigenVectorsCacheType
Base class for all process objects that output image data.
unsigned long SizeValueType
Definition: itkIntTypes.h:143
#define ITK_THREAD_RETURN_TYPE
Simulate a standard C array with copy semnatics.
Definition: itkFixedArray.h:50
std::vector< EigenValuesArrayType > EigenValuesCacheType
PatchBasedDenoisingBaseImageFilter< TInputImage, TOutputImage > Superclass
RealType AddUpdate(const DiffusionTensor3D< RealValueType > &a, const RealType &b)
typename::itk::Statistics::ImageToNeighborhoodSampleAdaptor< OutputImageType, BoundaryConditionType > ListAdaptorType
void ComputeMinMax(const Image< DiffusionTensor3D< PixelValueType >, ImageDimension > *img)
InputImageType::RegionType InputImageRegionType
This an abstract subsampler that constrains subsamples to be contained within a given image region...
A multi-dimensional iterator templated over image type that walks a region of pixels.
void SetComponent(T &pix, unsigned int, typename EnableIfC< IsSame< T, typename NumericTraits< T >::ValueType >::Value, RealValueType >::Type val) const
A method to generically set a component.
FixedArray< PixelValueType, 3 > EigenValuesArrayType
void SetComponent(T &pix, unsigned int idx, typename DisableIfC< IsSame< T, typename NumericTraits< T >::ValueType >::Value, RealValueType >::Type val) const
NumericTraits< PixelType >::RealType RealType
ImageRegionIterator< OutputImageType > OutputImageRegionIteratorType
NumericTraits< PixelValueType >::RealType RealValueType
ImageRegionConstIterator< InputImageType > InputImageRegionConstIteratorType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Define additional traits for native types such as int or float.
EnableIfC< IsSame< typename TImageType::PixelType, typename NumericTraits< typename TImageType::PixelType >::ValueType >::Value >::Type ComputeMinMax(const TImageType *img)
DisableIfC< IsSame< typename TImageType::PixelType, typename NumericTraits< typename TImageType::PixelType >::ValueType >::Value >::Type ComputeMinMax(const TImageType *img)
void ComputeDifferenceAndWeightedSquaredNorm(const DiffusionTensor3D< PixelValueType > &a, const DiffusionTensor3D< PixelValueType > &b, const RealArrayType &weight, bool useCachedComputations, SizeValueType cacheIndex, EigenValuesCacheType &eigenValsCache, EigenVectorsCacheType &eigenVecsCache, RealType &diff, RealArrayType &norm)
Represent a diffusion tensor as used in DTI images.
An implementation of the negation of the enable if idiom.
Definition: itkEnableIf.h:71
Templated n-dimensional image class.
Definition: itkImage.h:75
A multi-dimensional iterator templated over image type that walks a region of pixels.
itk::Statistics::RegionConstrainedSubsampler< PatchSampleType, InputImageRegionType > BaseSamplerType
RealType AddUpdate(const RealT &a, const RealType &b)