CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Plot2D.h
Go to the documentation of this file.
1 #ifndef PLOT_2D__H
2 #define PLOT_2D__H
3 
4 #include "PlotCompareUtility.h"
5 #include "Plot1D.h"
6 #include "PlotTypes.h"
7 
8 #include <TProfile.h>
9 #include <TH1F.h>
10 #include <TH2F.h>
11 #include <TGaxis.h>
12 #include <TLine.h>
13 
14 #include <iostream>
15 #include <fstream>
16 #include <string>
17 #include <math.h>
18 #include <stdio.h>
19 
20 template < >
21 bool PlotCompareUtility::compare<Plot2D>(HistoData *HD) {
22 
23  // get the reference and comparison histograms
24  TH2F *href2d = (TH2F *)HD->getRefHisto();
25  TH2F *hnew2d = (TH2F *)HD->getNewHisto();
26 
27  // do not run comparisons if either histogram is empty/broken
28  if (hnew2d == NULL || href2d == NULL || hnew2d->GetEntries() <= 1 || href2d->GetEntries() <= 1) {
29  //std::cerr << HD->getName() << " error: unable to retrieve histogram (or no entries)\n";
30  HD->setIsEmpty(true); return false;
31  }
32 
33  // prepare an overall result
34  bool projectionsPassed = true;
35 
36  // loop over axes (projections on one or both may be requested)
37  for (int axis = axisX; axis <= axisY; ++axis) {
38 
39  // for X: verify projections requested and proper Y binning of href2d and hnew2d
40  if (axis == axisX && !HD->getDoProjectionsX()) continue;
41  if (axis == axisX && href2d->GetNbinsY() != hnew2d->GetNbinsY()) {
42  std::cerr << HD->getName() << " error: incorrect number of bins for X projection tests\n";
43  projectionsPassed = false; continue;
44  }
45 
46  // for Y: verify projections requested and proper X binning of href2d and hnew2d
47  if (axis == axisY && !HD->getDoProjectionsY()) continue;
48  if (axis == axisY && href2d->GetNbinsX() != hnew2d->GetNbinsX()) {
49  std::cerr << HD->getName() << " error: incorrect number of bins for Y projection tests\n";
50  projectionsPassed = false; continue;
51  }
52 
53  // setup the rebinning variables
54  int nBins = (axis == axisX) ? href2d->GetNbinsY() : href2d->GetNbinsX();
55  int nProjections = (axis == axisX) ? HD->getMaxProjectionsX() : HD->getMaxProjectionsY();
56  int nGroups = (int)ceil(float(nBins) / nProjections);
57  bool rebinned = false;
58 
59  // for X projections: if required rebin a clone of href2d and hnew2d
60  if (axis == axisX && HD->getDoAllow2DRebinningX() && nGroups > 1) {
61  href2d = (TH2F *)(((TH2F *)(href2d->Clone()))->RebinY(nGroups));
62  hnew2d = (TH2F *)(((TH2F *)(hnew2d->Clone()))->RebinY(nGroups));
63  nBins = href2d->GetNbinsY();
64  rebinned = true;
65  }
66 
67  // for Y projections: if required rebin a clone of href2d and hnew2d
68  if (axis == axisY && HD->getDoAllow2DRebinningY() && nGroups > 1) {
69  href2d = (TH2F *)(((TH2F *)(href2d->Clone()))->RebinX(nGroups));
70  hnew2d = (TH2F *)(((TH2F *)(hnew2d->Clone()))->RebinX(nGroups));
71  nBins = href2d->GetNbinsX();
72  rebinned = true;
73  }
74 
75  // loop over bins in histograms (go backwords to keep in order)
76  //for (int bin = nBins; bin >= 1; --bin) {
77  for (int bin = 1; bin <= nBins; ++bin) {
78 
79  std::cout << "bin " << bin << " of " << nBins << std::endl;
80  // create some unique identifiers for the histogram names
81  TString projName = HD->getName() + (axis == axisX ? "_px" : "_py"); projName += bin;
82  TString newProjName = "new_"; newProjName += projName;
83  TString refProjName = "ref_"; refProjName += projName;
84 
85  // get the 1d projections for this bin out of the histogram
86  TH1D *hnew = (axis == axisX) ? hnew2d->ProjectionX(newProjName.Data(),bin,bin) : hnew2d->ProjectionY(newProjName.Data(),bin,bin);
87  TH1D *href = (axis == axisX) ? href2d->ProjectionX(refProjName.Data(),bin,bin) : href2d->ProjectionY(refProjName.Data(),bin,bin);
88 
89  // set histogram axis labels
90  hnew->GetXaxis()->SetTitle((axis == axisX ? hnew2d->GetXaxis()->GetTitle() : hnew2d->GetYaxis()->GetTitle()));
91  href->GetXaxis()->SetTitle((axis == axisX ? href2d->GetXaxis()->GetTitle() : href2d->GetYaxis()->GetTitle()));
92 
93  // allow Root to delete these histograms after display
94  hnew->SetBit(kCanDelete);
95  href->SetBit(kCanDelete);
96 
97  // create a new HistoData based on this projection
98  HistoData *proj = (axis == axisX) ? addProjectionXData(HD,projName.Data(),Plot1D,bin,hnew,href) \
99  : addProjectionYData(HD,projName.Data(),Plot1D,bin,hnew,href);
100 
101  // ignore empty bins
102  //if (hnew->Integral() == 0 || href->Integral() == 0) continue;
103  if (hnew->GetEntries() <= 1 || href->GetEntries() <= 1||hnew->Integral() == 0 || href->Integral() == 0) continue;
104 
105  // run this new HistoData through compare<Plot1D>
106  projectionsPassed &= compare<Plot1D>(proj);
107 
108  // get the high and low scores from this comparison
109  float lowScore = proj->getLowScore(); float highScore = proj->getHighScore();
110  if (lowScore < HD->getLowScore()) HD->setLowScore(lowScore);
111  if (highScore > HD->getHighScore()) HD->setHighScore(highScore);
112 
113  }
114 
115  // if 2d histograms were rebinned, delete the clone and re-get the original
116  if (rebinned) {
117  delete href2d; href2d = (TH2F *)HD->getRefHisto();
118  delete hnew2d; hnew2d = (TH2F *)HD->getNewHisto();
119  }
120 
121  }
122 
123  // check overall result
124  HD->setResult(projectionsPassed);
125  HD->setIsEmpty(false);
126 
127  // returns true on test passed and false on test failed
128  return projectionsPassed;
129 
130 }
131 
132 template < >
133 void PlotCompareUtility::makePlots<Plot2D>(HistoData *HD) {
134 
135  // do not make any new plot if empty
136  if (HD->getIsEmpty()) {
137  HD->setResultImage("NoData_Results.gif");
138  HD->setResultTarget("NoData_Results.gif");
139  return;
140  }
141 
142  // loop over the projections to make 1D plots
143  std::vector<HistoData>::iterator hd;
144  for (hd = projectionsX[HD].begin(); hd != projectionsX[HD].end(); hd++) makePlots<Plot1D>(&(*hd));
145  for (hd = projectionsY[HD].begin(); hd != projectionsY[HD].end(); hd++) makePlots<Plot1D>(&(*hd));
146 
147  // make projection summaries
148  for (int axis = axisX; axis <= axisY; ++axis) {
149 
150  // get the list of projections associated with this HistoData
151  std::vector<HistoData> *proj = (axis == axisX) ? &projectionsX[HD] : &projectionsY[HD];
152  if (proj == NULL || proj->size() == 0) continue;
153 
154  // get the 2d histograms
155  TH2F *hnew2d = (TH2F *)HD->getNewHisto();
156  TH2F *href2d = (TH2F *)HD->getRefHisto();
157 
158  // generate a reasonable width for the projections summary
159  int numHistos = proj->size();
160  int bodyWidth = int(float(numHistos * projectionsBarsThickness) * 1.5);
161  projectionsWidth = projectionsLeftMargin + projectionsRightMargin + bodyWidth;
162 
163  // the canvas is rescaled during gif conversion, so add padding to Canvas dimensions
164  int projectionsCanvasWidth = projectionsWidth + 4;
165  int projectionsCanvasHeight = projectionsHeight + 28;
166 
167  // create and setup projections canvas
168  TCanvas projectionsCanvas("projectionsCanvas","projectionsCanvas",projectionsCanvasWidth,projectionsCanvasHeight);
169  projectionsCanvas.SetFrameFillColor(10);
170  projectionsCanvas.SetLogy(1);
171  projectionsCanvas.SetTopMargin(float(projectionsTopMargin) / projectionsHeight);
172  projectionsCanvas.SetLeftMargin(float(projectionsLeftMargin) / projectionsWidth);
173  projectionsCanvas.SetRightMargin(float(projectionsRightMargin) / projectionsWidth);
174  projectionsCanvas.SetBottomMargin(float(projectionsBottomMargin) / projectionsHeight);
175  projectionsCanvas.Draw();
176 
177  // create and setup the summary histogram
178  TH1F projectionsSummary("projectionsSummary","Compatibility with Reference Histograms",numHistos,1,numHistos+1);
179  projectionsSummary.GetYaxis()->SetRangeUser(getThreshold()/10,2);
180  projectionsSummary.SetStats(0);
181 
182  // display histogram (take axis from original histogram)
183  projectionsSummary.Draw("AH");
184 
185  // draw X axis
186  float xMin = hnew2d->GetXaxis()->GetXmin();
187  float xMax = hnew2d->GetXaxis()->GetXmax();
188  int ticksNDiv = numHistos * 20 + bodyWidth / 50;//formerly *20
189  TGaxis *xAxis = new TGaxis(1,0,numHistos + 1,0,xMin,xMax,ticksNDiv,"");
190  if (axis == axisX) xAxis->SetTitle(hnew2d->GetYaxis()->GetTitle());
191  if (axis == axisY) xAxis->SetTitle(hnew2d->GetXaxis()->GetTitle());
192  xAxis->Draw();
193 
194  // draw Y axis
195  float yMin = getThreshold()/10;
196  float yMax = 2;
197  TGaxis *yAxis = new TGaxis(1,yMin,1,yMax,yMin,yMax,510,"G");
198  yAxis->SetTitle("Compatibility");
199  yAxis->Draw();
200 
201  // loop over projections and draw result
202  std::vector<HistoData>::iterator pd;
203  for (pd = proj->begin(); pd != proj->end(); pd++)
204  pd->drawResult(&projectionsSummary,true,false);
205 
206  // draw the pass/fail cutoff line
207  TLine passLine(1, getThreshold(), numHistos + 1, getThreshold());
208  passLine.SetLineColor(kRed);
209  passLine.SetLineWidth(2);
210  passLine.SetLineStyle(2);
211  passLine.Draw("SAME");
212 
213  // create the summary image
214  std::string gifName = HD->getName() + (axis == axisX ? "_Results_px.gif" : "_Results_py.gif");
215  projectionsCanvas.Print(gifName.c_str());
216 
217  // make overall projection plot of the original 2d histogram
218  std::string projName = HD->getName() + (axis == axisX ? "_py" : "_px");
219  std::string newBinsProj = projName + "_new";
220  std::string refBinsProj = projName + "_ref";
221  TH1D *href = (axis == axisX) ? href2d->ProjectionY(refBinsProj.c_str()) : href2d->ProjectionX(refBinsProj.c_str());
222  TH1D *hnew = (axis == axisX) ? hnew2d->ProjectionY(newBinsProj.c_str()) : hnew2d->ProjectionX(newBinsProj.c_str());
223 
224  // allow Root to delete these histograms after display
225  href->SetBit(kCanDelete);
226  hnew->SetBit(kCanDelete);
227 
228  // create a new HistoData based on this projection and plot it
229  HistoData allBins(projName,Plot1D,0,hnew,href);
230  allBins.setIsEmpty(false);
231  allBins.setShadedFillColor(HD->getShadedFillColor());
232  allBins.setShadedFillStyle(HD->getShadedFillStyle());
233  allBins.setShadedLineColor(HD->getShadedLineColor());
234  makePlots<Plot1D>(&allBins);
235 
236  // set the default image (axisY takes priority by default)
237  if (HD->getResultImage() == "" || axis == axisY) HD->setResultImage(projName + "_Results.gif");
238 
239  // set the default target (in case additional HTML code is/was not produced)
240  std::string currentTarget = HD->getResultTarget();
241  std::string xImgTarget = HD->getName() + "_px_Results.gif";
242  if (currentTarget == "" || (axis == axisY && currentTarget == xImgTarget)) HD->setResultTarget(projName + "_Results.gif");
243 
244  }
245 
246  /*
247  // make overall profile plot of the original 2d histogram
248  for (int axis = axisX; axis <= axisY; ++axis) {
249 
250  // make profile plots out of original 2D histograms
251  TProfile *pref = (axis == axisX) ? ((TH2F *)HD->getRefHisto())->ProfileY() : ((TH2F *)HD->getRefHisto())->ProfileX();
252  TProfile *pnew = (axis == axisX) ? ((TH2F *)HD->getNewHisto())->ProfileY() : ((TH2F *)HD->getNewHisto())->ProfileX();
253 
254  // renormalize results for display
255  renormalize(pref,pnew);
256 
257  // do not allow Root to deallocate this memory after drawing (tries to free twice?)
258  pref->SetBit(kCanDelete);
259  pnew->SetBit(kCanDelete);
260 
261  // set drawing options on the reference histogram
262  pref->SetStats(0);
263  pref->SetLineColor(kBlack);
264  pref->SetMarkerColor(kBlack);
265 
266  // set drawing options on the new histogram
267  pnew->SetStats(0);
268  pnew->SetLineColor(HD->getSolidLineColor());
269 
270  // place the test results as the title
271  TString title = HD->getName();
272 
273  // the canvas is rescaled during gif conversion, so add padding to Canvas dimensions
274  int plotsCanvasWidth = plotsWidth + 4;
275  int plotsCanvasHeight = plotsHeight + 28;
276 
277  // setup canvas for displaying the compared histograms
278  TCanvas hCanvas("hCanvas",title.Data(),plotsCanvasWidth,plotsCanvasHeight);
279  hCanvas.SetTopMargin(float(plotsTopMargin) / plotsHeight);
280  hCanvas.SetLeftMargin(float(plotsLeftMargin) / plotsWidth);
281  hCanvas.SetRightMargin(float(plotsRightMargin) / plotsWidth);
282  hCanvas.SetBottomMargin(float(plotsBottomMargin) / plotsHeight);
283  hCanvas.SetFrameFillColor(10);
284  hCanvas.SetGrid();
285  //hCanvas.SetLogy(1);
286  hCanvas.Draw();
287 
288  // draw the profiles
289  pref->Draw();
290  pnew->Draw("SAME");
291  if (HD->getDoDrawErrorBars()) pnew->Draw("E1SAME");
292 
293  // draw a legend
294  TLegend legend(0.15,0.01,0.3, 0.08);
295  legend.AddEntry(pnew,"New","lF");
296  legend.AddEntry(pref,"Reference","lF");
297  legend.SetFillColor(kNone);
298  legend.Draw("SAME");
299 
300  // create the plots overlay image
301  std::string gifName = HD->getName() + (axis == axisX ? "_pfx.gif" : "_pfy.gif");
302  hCanvas.Print(gifName.c_str());
303 
304  // set the default image (axisY takes priority by default)
305  if (HD->getResultImage() == "" || axis == axisY) HD->setResultImage(gifName);
306 
307  // set the default target (in case additional HTML code is/was not produced)
308  std::string currentTarget = HD->getResultTarget();
309  std::string xImgTarget = HD->getName() + "_pfx.gif";
310  if (currentTarget == "" || (axis == axisY && currentTarget == xImgTarget)) HD->setResultTarget(gifName);
311 
312  }
313  */
314 
315 
316 }
317 
318 template < >
319 void PlotCompareUtility::makeHTML<Plot2D>(HistoData *HD) {
320 
321  /* at present, makeHTML<Plot1D> does nothing, so don't waste the CPU cycles
322  // loop over projections and produce HTML
323  std::vector<HistoData>::iterator hd;
324  for (hd = projectionsX[HD].begin(); hd != projectionsX[HD].end(); hd++) makePlots<Plot1D>(&(*hd));
325  for (hd = projectionsY[HD].begin(); hd != projectionsY[HD].end(); hd++) makePlots<Plot1D>(&(*hd));
326  */
327 
328  // get the HistoData name for later reuse
329  std::string Name = HD->getName();
330 
331  // loop over the axes to see if projections were produced
332  bool pfDone[2] = { false, false };
333  for (int axis = axisX; axis <= axisY; axis++) {
334 
335  // get the list of projections associated with this HistoData
336  std::vector<HistoData> *proj = (axis == axisX) ? &projectionsX[HD] : &projectionsY[HD];
337  if (proj == NULL || proj->size() == 0) continue; else pfDone[axis] = true;
338 
339  // setup some names, etc. for insertion into the HTML
340  std::string gifNameProjections = Name + (axis == axisX ? "_Results_px.gif" : "_Results_py.gif");
341  std::string gifNameAllProj = Name + (axis == axisX ? "_py_Results.gif" : "_px_Results.gif");
342  std::string gifNameProfile = Name + (axis == axisX ? "_pfx.gif" : "_pfy.gif");
343  std::string gifBinPrefix = Name + (axis == axisX ? "_px" : "_py");
344 
345  // setup some locations to put thumbnails, etc.
346  int offset = 10;
347  int thumbWidth = plotsWidth / 4;
348  int thumbHeight = plotsHeight / 4;
349  int bodyWidth = projectionsWidth - projectionsLeftMargin - projectionsRightMargin;
350  int leftThumbPos = offset + projectionsLeftMargin + bodyWidth / 4 - thumbWidth / 2;
351  int rightThumbPos = leftThumbPos + bodyWidth / 2;
352  int thumbsLoc = projectionsTopMargin + thumbHeight / 2;
353 
354  // create the profile's HTML document
355  std::string htmlNameProfile = Name + (axis == axisX ? "_Results_px.html" : "_Results_py.html");
356  std::ofstream fout(htmlNameProfile.c_str());
357 
358  // print top portion of document
359  fout << "<!DOCTYPE gif PUBLIC '-//W3C//DTD HTML 4.01 Transitional//EN'>" << std::endl
360  << "<html>" << std::endl
361  << " <head>" << std::endl
362  << " <title>Compatibility of Projections for " << HD->getRefHisto()->GetTitle() << "</title>" << std::endl
363  << " <script type='text/javascript'>" << std::endl << std::endl
364  << " function tn(target,image,class) {" << std::endl
365  << " clear()" << std::endl
366  << " document.getElementById('thumb_div').setAttribute('class',class)" << std::endl
367  << " document.getElementById('thumb_div').setAttribute('className',class)" << std::endl
368  << " document.getElementById('thumb_link').href = target" << std::endl
369  << " document.getElementById('thumb_img').src = image" << std::endl
370  << " document.getElementById('thumb_img').width = '" << thumbWidth << "'" << std::endl
371  << " document.getElementById('thumb_img').height = '" << thumbHeight << "'" << std::endl
372  << " document.getElementById('thumb_img').border = '1'" << std::endl
373  << " }" << std::endl << std::endl
374  << " function clear() {" << std::endl
375  << " document.getElementById('thumb_link').href = '#'" << std::endl
376  << " document.getElementById('thumb_img').src = ''" << std::endl
377  << " document.getElementById('thumb_img').width = '0'" << std::endl
378  << " document.getElementById('thumb_img').height = '0'" << std::endl
379  << " document.getElementById('thumb_img').border = '0'" << std::endl
380  << " }" << std::endl << std::endl
381  << " </script>" << std::endl
382  << " </head>" << std::endl
383  << " <body onClick=\"window.location.href='index.html'\">" << std::endl
384  // << "<a href='index.html'>"
385  << " <style type='text/css'>" << std::endl
386  << " #thumb_div {}" << std::endl
387  << " div.thumb_left {position: absolute; left: " << leftThumbPos << "px; top: " << thumbsLoc << "px;}" << std::endl
388  << " div.thumb_right {position: absolute; left: " << rightThumbPos << "px; top: " << thumbsLoc << "px;}" << std::endl
389  << " #main_d {position: absolute; left: " << offset << "px;}" << std::endl
390  << " a:link {color: #000000}" << std::endl
391  << " a:visited {color: #000000}" << std::endl
392  << " a:hover {color: #000000}" << std::endl
393  << " a:active {color: #000000}" << std::endl
394  << " </style>" << std::endl
395  << " <div id='main_d'>" << std::endl
396  // << " <p>" << HD->getRefHisto()->GetTitle() << "</p>" //include the Title of the Plot as a title of the page
397  << " <img src='" << gifNameProjections << "' usemap='#results' alt=''"
398  << " height=" << projectionsHeight << " width=" << projectionsWidth << " border=0>" << std::endl
399  << " <map id='#results' name='results' onMouseOut=\"clear()\">" << std::endl;
400 
401  // loop over projections
402  std::vector<HistoData>::iterator pd;
403  for (pd = proj->begin(); pd != proj->end(); pd++) {
404 
405  // determine map coordinates for this bin (1 pixel offset due to borders?)
406  int bin = pd->getBin();
407  int x1 = projectionsLeftMargin + int(float(bin * 1.5 - 1.25) * projectionsBarsThickness);
408  int x2 = x1 + projectionsBarsThickness;
409  int y1 = projectionsTopMargin + 1;
410  int y2 = projectionsHeight - projectionsBottomMargin;
411  std::string image = pd->getResultImage();
412  std::string target = pd->getResultTarget();
413 
414  // add coordinates area to image map
415  std::string tnClass = (bin - 1 >= float(proj->size()) / 2 ? "thumb_left" : "thumb_right");
416  fout << " <area shape='rect' alt='' coords='" << x1 << "," << y1 << "," << x2 << "," << y2 << "'"
417  << " href='" << target << "' onMouseOver=\"tn('" << target << "','" << image << "','" << tnClass
418  << "')\" "
419  << "onMouseDown=\"window.location.href='" << target << "'\">" << std::endl;
420 
421  }
422 
423  fout << " <area shape='default' nohref='nohref' onMouseDown='window.location.reload()' alt=''>" << std::endl
424  << " </map>" << std::endl
425  << " <br><img src=\"" << gifNameAllProj << "\">" << std::endl
426  << " </div>" << std::endl
427  << " <div id='thumb_div'><a href='#' id='thumb_link'><img src='' id='thumb_img' width=0 height=0 border=0></a></div>" << std::endl
428  // << " </a>"
429  << " </body>" << std::endl
430  << "</html>" << std::endl;
431 
432  // close the file
433  HD->setResultTarget(htmlNameProfile);
434  fout.close();
435 
436  }
437 
438  // if both profile dimensions were filled, we need an additional HTML document
439  if (pfDone[axisX] && pfDone[axisY]) {
440 
441  // create HTML support code for this HistoData
442  std::string html = Name + "_Results_Profiles.html";
443  std::ofstream fout(html.c_str());
444 
445  // make a simple frames portal to show both profile
446  fout << "<html>" << std::endl
447  << " <frameset rows=\"50%,50%\">" << std::endl
448  << " <frame src=\"" << Name << "_Results_py.html\">" << std::endl
449  << " <frame src=\"" << Name << "_Results_px.html\">" << std::endl
450  << " <noframes><body>" << std::endl
451  << " unable to display frames -- click to select desired page" << std::endl
452  << " <br><a href =\"" << Name << "_Results_py.html\">Y Projections</a>" << std::endl
453  << " <br><a href =\"" << Name << "_Results_px.html\">X Projections</a>" << std::endl
454  << " </body></noframes>" << std::endl
455  << " </frameset>" << std::endl
456  << "</html>" << std::endl;
457 
458  // close the file
459  HD->setResultTarget(html);
460  fout.close();
461 
462  }
463 
464 }
465 
466 #endif // PLOT_2D__H
#define NULL
Definition: scimark2.h:8
float getLowScore() const
Definition: HistoData.h:41
unsigned int offset(bool)
float getHighScore() const
Definition: HistoData.h:42
#define begin
Definition: vmac.h:30
tuple cout
Definition: gather_cfg.py:121
TH1 * getNewHisto() const
Definition: HistoData.h:22