ANIMA  4.0
animaPatientToGroupComparisonImageFilter.h
Go to the documentation of this file.
1 #pragma once
2 
3 #include <iostream>
5 #include <itkVectorImage.h>
6 #include <itkImage.h>
7 
8 #include <vector>
9 
10 namespace anima
11 {
12 
21 template <class PixelScalarType>
23  public anima::MaskedImageToImageFilter< itk::VectorImage <PixelScalarType, 3> , itk::Image <PixelScalarType, 3> >
24 {
25 public:
28  typedef itk::SmartPointer<Self> Pointer;
29  typedef itk::SmartPointer<const Self> ConstPointer;
30 
32  itkNewMacro(Self);
33 
36 
38  typedef itk::VectorImage <PixelScalarType, 3> InputImageType;
39  typedef itk::Image <PixelScalarType, 3> OutputImageType;
40 
41  typedef typename InputImageType::PixelType VectorType;
42 
43  enum TestType
44  {
47  };
48 
50  typedef typename InputImageType::Pointer InputImagePointer;
51  typedef typename InputImageType::IndexType InputImageIndexType;
52  typedef typename OutputImageType::Pointer OutputImagePointer;
53 
57 
58  typedef vnl_matrix<double> CovMatrixType;
59 
60  void AddDatabaseInput(InputImageType *tmpIm)
61  {
62  m_DatabaseImages.push_back(tmpIm);
63  }
64 
65  itkSetMacro(NumEigenValuesPCA, unsigned int);
66  itkSetMacro(ExplainedRatio, double);
67 
68  itkSetMacro(StatisticalTestType, TestType);
69 
70 protected:
72  : Superclass()
73  {
74  this->SetNumberOfRequiredOutputs(2);
75  this->SetNthOutput(0,this->MakeOutput(0));
76  this->SetNthOutput(1,this->MakeOutput(1));
77 
78  m_ExplainedRatio = 0.9;
79  m_NumEigenValuesPCA = 6;
80 
81  m_DatabaseImages.clear();
82  m_StatisticalTestType = FISHER;
83  }
84 
86 
87  virtual void BeforeThreadedGenerateData() ITK_OVERRIDE;
88  void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) ITK_OVERRIDE;
89 
90  virtual unsigned int SampleFromDiffusionModels(std::vector <VectorType> &databaseValues, VectorType &patientVectorValue)
91  {
92  // Reimplement in subclasses if they need to sample from the distribution before doing the projection and comparison
93  return patientVectorValue.GetSize();
94  }
95 
96 private:
97  ITK_DISALLOW_COPY_AND_ASSIGN(PatientToGroupComparisonImageFilter);
98 
99  bool isZero(const VectorType &vec);
100 
101  unsigned int GetPCAVectorsFromData(std::vector < itk::VariableLengthVector <double> > &databaseVectors,
102  itk::VariableLengthVector <double> &patientVectorValue);
103 
104  unsigned int m_NumEigenValuesPCA;
105  double m_ExplainedRatio;
106 
107  std::vector <InputImagePointer> m_DatabaseImages;
108  TestType m_StatisticalTestType;
109 };
110 
111 } // end namespace anima
112 
PatientToGroupComparisonImageFilter< PixelScalarType > Self
anima::MaskedImageToImageFilter< InputImageType, OutputImageType > Superclass
itk::Image< unsigned char, TInputImage::ImageDimension > MaskImageType
Superclass::OutputImageRegionType OutputImageRegionType
virtual unsigned int SampleFromDiffusionModels(std::vector< VectorType > &databaseValues, VectorType &patientVectorValue)
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) ITK_OVERRIDE
Implements patient to group comparison as in http://dx.doi.org/10.1007/978-3-540-85988-8_116.