CMS 3D CMS Logo

Plot1D.h
Go to the documentation of this file.
1 #ifndef PLOT_1D__H
2 #define PLOT_1D__H
3 
4 #include "PlotCompareUtility.h"
5 #include "PlotTypes.h"
6 
7 #include <TCanvas.h>
8 #include <TH1F.h>
9 #include <TLegend.h>
10 #include <TString.h>
11 #include <TText.h>
12 
13 #include <fstream>
14 #include <iostream>
15 #include <string>
16 
17 template <>
18 bool PlotCompareUtility::compare<Plot1D>(HistoData *HD) {
19  // get the reference and comparison histograms
20  TH1F *href = (TH1F *)HD->getRefHisto();
21  TH1F *hnew = (TH1F *)HD->getNewHisto();
22 
23  // do not run comparisons if either histogram is empty/broken
24  if (hnew == nullptr || href == nullptr || hnew->GetEntries() <= 1 || href->GetEntries() <= 1) {
25  // std::cerr << HD->getName() << " error: unable to retrieve histogram (or
26  // no entries)\n";
27  HD->setIsEmpty(true);
28  return false;
29  }
30 
31  // rebin histograms and center on common range
32  if (HD->getDoAllow1DRebinning())
33  centerRebin(href, hnew);
34 
35  // run statistical comparisons
36  double ksScore = hnew->KolmogorovTest(href, "D");
37  double chi2Score = hnew->Chi2Test(href, "uup");
38 
39  // renormalize histograms for common display
40  renormalize(href, hnew);
41 
42  // set ks/chi2 and determine high/low scores
43  HD->setKSScore(ksScore);
44  HD->setChi2Score(chi2Score);
45  if (ksThreshold > 0 && chi2Threshold > 0) {
46  HD->setLowScore(ksScore < chi2Score ? ksScore : chi2Score);
47  HD->setHighScore(ksScore > chi2Score ? ksScore : chi2Score);
48  } else if (ksThreshold > 0) {
49  HD->setLowScore(ksScore);
50  HD->setHighScore(ksScore);
51  } else if (chi2Threshold > 0) {
52  HD->setLowScore(chi2Score);
53  HD->setHighScore(chi2Score);
54  } else
55  std::cerr << "error: no test performed? chi2Threshold and ksThreshold <= 0\n";
56 
57  // check overall result
58  bool passed = (ksScore >= ksThreshold && chi2Score >= chi2Threshold);
59  HD->setResult(passed);
60  // if (!passed) std::cout << "NOTPASSED!!!!"; else std::cout << "YESPASSED!";
61  // returns true on test passed and false on test failed
62  HD->setIsEmpty(false);
63  return passed;
64 }
65 
66 template <>
67 void PlotCompareUtility::makePlots<Plot1D>(HistoData *HD) {
68  std::cerr << HD->getName() << "makePlots<Plot1D>\n";
69 
70  // do not make any new plot if empty
71  if (HD->getIsEmpty()) {
72  HD->setResultImage("NoData_Results.gif");
73  HD->setResultTarget("NoData_Results.gif");
74  return;
75  }
76 
77  // get the reference and comparison histograms
78  TH1F *href = (TH1F *)HD->getRefHisto();
79  TH1F *hnew = (TH1F *)HD->getNewHisto();
80 
81  // set drawing options on the reference histogram
82  href->SetStats(false);
83  href->SetLineWidth(2);
84  href->SetLineColor(14);
85  href->SetMarkerColor(14);
86  href->SetFillColor(18);
87 
88  // set drawing options on the new histogram
89  hnew->SetStats(false);
90  hnew->SetLineWidth(2);
91  hnew->SetLineColor(HD->getShadedLineColor());
92  hnew->SetFillStyle(HD->getShadedFillStyle());
93  hnew->SetFillColor(HD->getShadedFillColor());
94 
95  // place the test results as the title
96  TString title = HD->getName();
97  if (ksThreshold > 0) {
98  title += " KS Score = ";
99  title += HD->getKSScore();
100  }
101  if (chi2Threshold > 0) {
102  title += " Chi^2 Score = ";
103  title += HD->getChi2Score();
104  }
105 
106  // the canvas is rescaled during gif conversion, so add padding to Canvas
107  // dimensions
108  int plotsCanvasWidth = plotsWidth + 4;
109  int plotsCanvasHeight = plotsHeight + 28;
110 
111  // setup canvas for displaying the compared histograms
112  TCanvas hCanvas("hCanvas", "hCanvas", plotsCanvasWidth, plotsCanvasHeight);
113  hCanvas.SetTopMargin(float(plotsTopMargin) / plotsHeight);
114  hCanvas.SetLeftMargin(float(plotsLeftMargin) / plotsWidth);
115  hCanvas.SetRightMargin(float(plotsRightMargin) / plotsWidth);
116  hCanvas.SetBottomMargin(float(plotsBottomMargin) / plotsHeight);
117  hCanvas.SetFrameFillColor(10);
118  hCanvas.SetGrid();
119  hCanvas.SetLogy(1);
120  hCanvas.Draw();
121 
122  TText canvasTitle(0.1, 0.97, title.Data());
123  canvasTitle.Draw("SAME");
124 
125  // draw the histograms
126  href->Draw();
127  hnew->Draw("SAME");
128  if (HD->getDoDrawErrorBars())
129  hnew->Draw("E1SAME");
130 
131  // draw a legend
132  TLegend legend(0.15, 0.01, 0.3, 0.08);
133  legend.AddEntry(hnew, "New", "lF");
134  legend.AddEntry(href, "Reference", "lF");
135  legend.SetFillColor(kNone);
136  legend.Draw("SAME");
137 
138  // create the plots overlay image
139  std::string gifName = HD->getName() + "_Results.gif";
140  if (HD->getResultImage().empty())
141  HD->setResultImage(gifName);
142  if (HD->getResultTarget().empty())
143  HD->setResultTarget(gifName);
144  std::cerr << "About to print" << gifName << "\n";
145  hCanvas.Print(gifName.c_str());
146 }
147 
148 template <>
149 void PlotCompareUtility::makeHTML<Plot1D>(HistoData *HD) {
150  /* HTML is not presently required for 1D histograms -- do nothing
151 
152  // create HTML support code for this HistoData
153  std::string Name = hd->getName();
154  std::string gifName = Name + "_Results.gif";
155  std::string html = Name + "_Results.html";
156  ofstream fout(html.c_str());
157 
158  // simply link to the appropriate image overlay
159  fout << "<html><body><img src=\"" << gifName << "\"></body></html>" << endl;
160 
161  // close the file
162  fout.close();
163 
164  */
165 }
166 
167 #endif // PLOT_1D__H
HistoData::getResultTarget
std::string getResultTarget() const
Definition: HistoData.h:24
HistoData::setResultImage
void setResultImage(std::string Image)
Definition: HistoData.h:59
runGCPTkAlMap.title
string title
Definition: runGCPTkAlMap.py:94
HistoData::getKSScore
float getKSScore() const
Definition: HistoData.h:38
HistoData::getRefHisto
TH1 * getRefHisto() const
Definition: HistoData.h:22
HistoData::getIsEmpty
bool getIsEmpty() const
Definition: HistoData.h:43
PlotCompareUtility::ksThreshold
double ksThreshold
Definition: PlotCompareUtility.h:136
HistoData::setHighScore
void setHighScore(float Score)
Definition: HistoData.h:80
PlotCompareUtility::centerRebin
void centerRebin(TH1 *, TH1 *)
HistoData::getShadedFillStyle
int getShadedFillStyle() const
Definition: HistoData.h:52
PlotCompareUtility::plotsHeight
int plotsHeight
Definition: PlotCompareUtility.h:169
HistoData::getChi2Score
float getChi2Score() const
Definition: HistoData.h:39
PlotCompareUtility::plotsLeftMargin
int plotsLeftMargin
Definition: PlotCompareUtility.h:171
PlotCompareUtility::plotsTopMargin
int plotsTopMargin
Definition: PlotCompareUtility.h:170
HistoData::setLowScore
void setLowScore(float Score)
Definition: HistoData.h:79
PlotCompareUtility::plotsBottomMargin
int plotsBottomMargin
Definition: PlotCompareUtility.h:173
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
HistoData::getNewHisto
TH1 * getNewHisto() const
Definition: HistoData.h:21
HistoData::getName
std::string getName() const
Definition: HistoData.h:17
PlotCompareUtility::plotsRightMargin
int plotsRightMargin
Definition: PlotCompareUtility.h:172
HistoData
Definition: HistoData.h:9
HistoData::setResultTarget
void setResultTarget(std::string Target)
Definition: HistoData.h:60
listHistos.legend
legend
Definition: listHistos.py:41
PlotCompareUtility::chi2Threshold
double chi2Threshold
Definition: PlotCompareUtility.h:137
HistoData::getDoAllow1DRebinning
bool getDoAllow1DRebinning() const
Definition: HistoData.h:28
WDecay::kNone
Definition: TopGenEvent.h:27
PlotTypes.h
PlotCompareUtility::plotsWidth
int plotsWidth
Definition: PlotCompareUtility.h:168
postprocess-scan-build.href
href
Definition: postprocess-scan-build.py:26
TriggerAnalyzer.passed
passed
Definition: TriggerAnalyzer.py:62
PlotCompareUtility.h
HistoData::getShadedFillColor
int getShadedFillColor() const
Definition: HistoData.h:51
HistoData::setChi2Score
void setChi2Score(float Score)
Definition: HistoData.h:78
HistoData::setIsEmpty
void setIsEmpty(bool Toggle)
Definition: HistoData.h:82
HistoData::setKSScore
void setKSScore(float Score)
Definition: HistoData.h:77
HistoData::getShadedLineColor
int getShadedLineColor() const
Definition: HistoData.h:50
HistoData::setResult
void setResult(bool Result)
EcnaPython_AdcPeg12_S1_10_R170298_1_0_150_Dee0.cerr
cerr
Definition: EcnaPython_AdcPeg12_S1_10_R170298_1_0_150_Dee0.py:8
PlotCompareUtility::renormalize
void renormalize(TH1 *, TH1 *)
HistoData::getResultImage
std::string getResultImage() const
Definition: HistoData.h:23
HistoData::getDoDrawErrorBars
bool getDoDrawErrorBars() const
Definition: HistoData.h:27