ANIMA  4.0
animaSegPerfCAnalyzer.h
Go to the documentation of this file.
1 #pragma once
2 
3 #include <math.h>
4 #include <vector>
5 #include <algorithm>
6 #include <itkImageFileReader.h>
7 #include <itkImageRegionConstIterator.h>
9 #include <itkLabelContourImageFilter.h>
10 #include <itkBinaryContourImageFilter.h>
11 #include <itkSimpleFilterWatcher.h>
12 #include <itkHausdorffDistanceImageFilter.h>
13 #include <itkDirectedHausdorffDistanceImageFilter.h>
14 #include <itkContourMeanDistanceImageFilter.h>
15 #include <itkFlipImageFilter.h>
16 #include <itkImageDuplicator.h>
17 
18 namespace anima
19 {
20 
26 {
27 public:
28  SegPerfCAnalyzer(std::string &pi_pchImageTestName, std::string &pi_pchImageRefName, bool advancedEvaluation);
30 
32 
33  void setNbThreads(int pi_iNbThreads);
34  void selectCluster(unsigned int);
35 
36  void setDetectionThresholdAlpha(double pi_fVal)
37  {
38  m_dfDetectionThresholdAlpha = pi_fVal;
39  }
40 
41  void setDetectionThresholdBeta(double pi_fVal)
42  {
43  m_dfDetectionThresholdBeta = pi_fVal;
44  }
45 
46  void setDetectionThresholdGamma(double pi_fVal)
47  {
48  m_dfDetectionThresholdGamma = pi_fVal;
49  }
50 
51  void setMinLesionVolumeDetection(double pi_fVal)
52  {
53  m_dfMinLesionVolumeDetection = pi_fVal;
54  }
55 
56  double computeHausdorffDist();
57  double computeMeanDist();
59  void computeITKMeasures();
60 
61  //getters
62  double getUnionOverlap();
63  double getMeanOverlap();
64  double getSensitivity();
65  double getSpecificity();
66  double getPPV();
67  double getNPV();
68  double getDiceCoefficient();
69  double getJaccardCoefficient();
70  double getRelativeVolumeError();
71  bool getDetectionMarks(double&po_fPPVL, double&po_fSensL, double&po_fF1);
72 
73  int getNumberOfClusters();
74 
75 protected:
76  void formatLabels();
77  void contourDectection();
78  void checkNumberOfLabels(int, int);
79 
80  int getTruePositiveLesions(int pi_iNbLabelsRef, int pi_iNbLabelsTest, int **pi_ppiOverlapTab);
81  bool falsePositiveRatioTester(int pi_iLessionReference, int pi_iNbLabelsTest, int pi_iNbLabelsRef, int **pi_ppiOverlapTab, int *pi_piTPRowSumTab, int *pi_piColumnSumTab, double pi_dfBeta, double pi_dfGamma);
82  void getOverlapTab(int&po_iNbLabelsRef, int&po_iNbLabelsTest, int **&po_ppiOverlapTab);
83 
84 private:
85  SegPerfCAnalyzer() {}
86  void transposer(int pi_iNbLabelsRef, int pi_iNbLabelsTest, int * *pi_ppiOverlapTab, int* *&po_rppiOverlapTabTransposed);
87  void removeOverlapTab(int **pi_ppiOverlapTab, int pi_iNbLabelsRef);
88 
89  unsigned int m_uiNbLabels;
90  bool m_bValuesComputed;
91  bool m_bContourDetected;
93  double m_dfDetectionThresholdAlpha;
94  double m_dfDetectionThresholdBeta;
95  double m_dfDetectionThresholdGamma;
96  double m_dfMinLesionVolumeDetection;
97 
98 
99  typedef itk::Image <unsigned short, 3> ImageType;
100  typedef itk::ImageFileReader <ImageType> ImageReaderType;
101  typedef itk::ImageRegionConstIterator <ImageType> ImageIteratorType;
103 
104  ImageType::Pointer m_imageTest;
105  ImageType::Pointer m_imageRef;
106  ImageType::Pointer m_imageTestContour;
107  ImageType::Pointer m_imageRefContour;
108  ImageType::Pointer m_imageRefDuplicated;
109  ImageType::Pointer m_imageTestDuplicated;
110 
111  FilterType::Pointer m_oFilter;
112  itk::ThreadIdType m_ThreadNb;
113 };
114 
115 } // end namespace anima
double getPPV()
Getter of Positive predictive value.
double getDiceCoefficient()
Getter of Dice coefficient.
double computeAverageSurfaceDistance()
Compute average surface distance.
void setDetectionThresholdBeta(double pi_fVal)
bool falsePositiveRatioTester(int pi_iLessionReference, int pi_iNbLabelsTest, int pi_iNbLabelsRef, int **pi_ppiOverlapTab, int *pi_piTPRowSumTab, int *pi_piColumnSumTab, double pi_dfBeta, double pi_dfGamma)
Compute if a lesion is detected or not with beta and gamma thresholds.
double getJaccardCoefficient()
Getter of Jaccard coefficient.
void computeITKMeasures()
Compute different measures with ITK to evaluate segmentation.
double computeMeanDist()
Compute mean distance.
double getNPV()
Getter of Negative predictive value.
void selectCluster(unsigned int)
Select the cluster we want to use to compute evaluation results.
double getMeanOverlap()
Getter of Mean overlap.
bool checkImagesMatrixAndVolumes()
Check if the 2 inputs images are compatible.
void contourDectection()
Extract the contour of the image to evaluate and the ground truth.
void getOverlapTab(int &po_iNbLabelsRef, int &po_iNbLabelsTest, int **&po_ppiOverlapTab)
Compute a table of lesion overlap between 2 images.
int getNumberOfClusters()
Return the number of clusters.
void formatLabels()
Format labels values for the image to evaluate and the ground truth to obtain : background pixels val...
double getRelativeVolumeError()
Getter of Relative volume error.
void setNbThreads(int pi_iNbThreads)
Set the number of threads to use for the computation of ITK.
bool getDetectionMarks(double &po_fPPVL, double &po_fSensL, double &po_fF1)
Compute useful variables to get detection scores.
Class to compute various metrics to evaluate segmentation results.
void setDetectionThresholdAlpha(double pi_fVal)
double getSensitivity()
Getter of Sensibility.
int getTruePositiveLesions(int pi_iNbLabelsRef, int pi_iNbLabelsTest, int **pi_ppiOverlapTab)
Getter of True positive lesions.
double getSpecificity()
Getter of Specificity.
double getUnionOverlap()
Getter of Union overlap.
void setMinLesionVolumeDetection(double pi_fVal)
double computeHausdorffDist()
Compute Haussdorf distance.
void checkNumberOfLabels(int, int)
Check if the number of labels is the same for both input images.
void setDetectionThresholdGamma(double pi_fVal)