ANIMA  4.0
animaFDRCorrectImageFilter.hxx
Go to the documentation of this file.
1 #pragma once
3 
4 #include <itkImageRegionConstIterator.h>
5 #include <itkImageRegionIterator.h>
6 
7 #include <animaFDRCorrection.h>
8 
9 namespace anima
10 {
11 
12 template <class PixelScalarType>
13 void
14 FDRCorrectImageFilter <PixelScalarType>
15 ::GenerateData()
16 {
17  this->AllocateOutputs();
18 
19  if (!m_MaskImage)
20  {
21  this->CreateFullMask();
22  this->GetOutput()->SetRequestedRegion(this->GetInput()->GetLargestPossibleRegion());
23  }
24 
25  typedef itk::ImageRegionConstIterator <TInputImage> InputImageIteratorType;
26  typedef itk::ImageRegionConstIterator <MaskImageType> MaskIteratorType;
27  typedef itk::ImageRegionIterator <TOutputImage> OutputImageIteratorType;
28 
29  MaskIteratorType maskItr(m_MaskImage, this->GetOutput()->GetRequestedRegion());
30  InputImageIteratorType inItr(this->GetInput(), this->GetOutput()->GetRequestedRegion());
31  OutputImageIteratorType outItr(this->GetOutput(), this->GetOutput()->GetRequestedRegion());
32 
33  std::vector <double> pvalues;
34 
35  while (!maskItr.IsAtEnd())
36  {
37  if (maskItr.Get() == 0)
38  {
39  ++maskItr;
40  ++inItr;
41  continue;
42  }
43 
44  pvalues.push_back(inItr.Get());
45 
46  ++inItr;
47  ++maskItr;
48  }
49 
50  if (m_BYCorrection)
51  anima::BYCorrection(pvalues, m_QValue);
52  else
53  anima::BHCorrection(pvalues, m_QValue);
54 
55  maskItr.GoToBegin();
56 
57  unsigned int pos = 0;
58  while (!maskItr.IsAtEnd())
59  {
60  if (maskItr.Get() == 0)
61  {
62  ++maskItr;
63  outItr.Set(0);
64  ++outItr;
65  continue;
66  }
67 
68  outItr.Set(pvalues[pos]);
69  ++pos;
70 
71  ++outItr;
72  ++maskItr;
73  }
74 }
75 
76 template <class PixelScalarType>
77 void
79 ::CreateFullMask()
80 {
81  if (m_MaskImage)
82  return;
83 
84  m_MaskImage = MaskImageType::New();
85  m_MaskImage->Initialize();
86  m_MaskImage->SetRegions(this->GetInput()->GetLargestPossibleRegion());
87  m_MaskImage->SetOrigin(this->GetInput()->GetOrigin());
88  m_MaskImage->SetDirection(this->GetInput()->GetDirection());
89  m_MaskImage->SetSpacing(this->GetInput()->GetSpacing());
90 
91  m_MaskImage->Allocate();
92 
93  m_MaskImage->FillBuffer(1);
94 }
95 
96 } // end namespace anima
void BHCorrection(std::vector< double > &pvalues, double qValue)
void BYCorrection(std::vector< double > &pvalues, double qValue)