Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
SamplingDistPlot.cxx
Go to the documentation of this file.
1 // @(#)root/roostats:$Id$
2 // Authors: Sven Kreiss June 2010
3 // Authors: Kyle Cranmer, Lorenzo Moneta, Gregory Schott, Wouter Verkerke
4 /*************************************************************************
5  * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 /** \class RooStats::SamplingDistPlot
13  \ingroup Roostats
14 
15 This class provides simple and straightforward utilities to plot SamplingDistribution
16 objects.
17 */
18 
20 
22 
23 #include "RooRealVar.h"
24 #include "TStyle.h"
25 #include "TLine.h"
26 #include "TFile.h"
27 #include "TVirtualPad.h" // for gPad
28 
29 #include <algorithm>
30 #include <iostream>
31 
32 
33 #include "RooMsgService.h"
34 
35 #include <limits>
36 #define NaN std::numeric_limits<float>::quiet_NaN()
37 #include "TMath.h"
38 #define IsNaN(a) TMath::IsNaN(a)
39 
40 ClassImp(RooStats::SamplingDistPlot);
41 
42 using namespace RooStats;
43 using namespace std;
44 
45 ////////////////////////////////////////////////////////////////////////////////
46 /// SamplingDistPlot default constructor with bin size
47 
48 SamplingDistPlot::SamplingDistPlot(Int_t nbins) :
49  fHist(0),
50  fLegend(NULL),
51  fItems(),
52  fOtherItems(),
53  fRooPlot(NULL),
54  fLogXaxis(kFALSE),
55  fLogYaxis(kFALSE),
56  fXMin(NaN), fXMax(NaN), fYMin(NaN), fYMax(NaN),
57  fApplyStyle(kTRUE),
58  fFillStyle(3004)
59 {
60  fIterator = fItems.MakeIterator();
61  fIsWeighted = kFALSE;
62  fBins = nbins;
63  fMarkerType = 20;
64  fColor = 1;
65 }
66 
67 ////////////////////////////////////////////////////////////////////////////////
68 /// SamplingDistPlot constructor with bin size, min and max
69 
70 SamplingDistPlot::SamplingDistPlot(Int_t nbins, Double_t min, Double_t max) :
71  fHist(0),
72  fLegend(NULL),
73  fItems(),
74  fOtherItems(),
75  fRooPlot(NULL),
76  fLogXaxis(kFALSE),
77  fLogYaxis(kFALSE),
78  fXMin(NaN), fXMax(NaN), fYMin(NaN), fYMax(NaN),
79  fApplyStyle(kTRUE),
80  fFillStyle(3004)
81 {
82  fIterator = fItems.MakeIterator();
83  fIsWeighted = kFALSE;
84  fBins = nbins;
85  fMarkerType = 20;
86  fColor = 1;
87 
88  SetXRange( min, max );
89 }
90 
91 ////////////////////////////////////////////////////////////////////////////////
92 /// destructors - delete objects contained in the list
93 
94 SamplingDistPlot::~SamplingDistPlot() {
95  fItems.Delete();
96  fOtherItems.Delete();
97  if (fRooPlot) delete fRooPlot;
98 }
99 
100 
101 ////////////////////////////////////////////////////////////////////////////////
102 /// adds sampling distribution (and normalizes if "NORMALIZE" is given as an option)
103 
104 Double_t SamplingDistPlot::AddSamplingDistribution(const SamplingDistribution *samplingDist, Option_t *drawOptions) {
105  fSamplingDistr = samplingDist->GetSamplingDistribution();
106  if( fSamplingDistr.empty() ) {
107  coutW(Plotting) << "Empty sampling distribution given to plot. Skipping." << endl;
108  return 0.0;
109  }
110  SetSampleWeights(samplingDist);
111 
112  TString options(drawOptions);
113  options.ToUpper();
114 
115  Double_t xmin(TMath::Infinity()), xmax(-TMath::Infinity());
116  // remove cases where xmin and xmax are +/- inf
117  for( unsigned int i=0; i < fSamplingDistr.size(); i++ ) {
118  if( fSamplingDistr[i] < xmin && fSamplingDistr[i] != -TMath::Infinity() ) {
119  xmin = fSamplingDistr[i];
120  }
121  if( fSamplingDistr[i] > xmax && fSamplingDistr[i] != TMath::Infinity() ) {
122  xmax = fSamplingDistr[i];
123  }
124  }
125  if( xmin >= xmax ) {
126  coutW(Plotting) << "Could not determine xmin and xmax of sampling distribution that was given to plot." << endl;
127  xmin = -1.0;
128  xmax = 1.0;
129  }
130 
131 
132  // add 1.5 bins left and right
133  assert(fBins > 1);
134  double binWidth = (xmax-xmin)/(fBins);
135  Double_t xlow = xmin - 1.5*binWidth;
136  Double_t xup = xmax + 1.5*binWidth;
137  if( !IsNaN(fXMin) ) xlow = fXMin;
138  if( !IsNaN(fXMax) ) xup = fXMax;
139 
140  fHist = new TH1F(samplingDist->GetName(), samplingDist->GetTitle(), fBins, xlow, xup);
141  fHist->SetDirectory(0); // make the object managed by this class
142 
143  if( fVarName.Length() == 0 ) fVarName = samplingDist->GetVarName();
144  fHist->GetXaxis()->SetTitle(fVarName.Data());
145 
146 
147  std::vector<Double_t>::iterator valuesIt = fSamplingDistr.begin();
148  for (int w_idx = 0; valuesIt != fSamplingDistr.end(); ++valuesIt, ++w_idx) {
149  if (fIsWeighted) fHist->Fill(*valuesIt, fSampleWeights[w_idx]);
150  else fHist->Fill(*valuesIt);
151  }
152 
153  // NORMALIZATION
154  fHist->Sumw2();
155  double weightSum = 1.0;
156  if(options.Contains("NORMALIZE")) {
157  weightSum = fHist->Integral("width");
158  fHist->Scale(1./weightSum);
159 
160  options.ReplaceAll("NORMALIZE", "");
161  options.Strip();
162  }
163 
164 
165  //some basic aesthetics
166  fHist->SetMarkerStyle(fMarkerType);
167  fHist->SetMarkerColor(fColor);
168  fHist->SetLineColor(fColor);
169 
170  fMarkerType++;
171  fColor++;
172 
173  fHist->SetStats(kFALSE);
174 
175  addObject(fHist, options.Data());
176 
177  TString title = samplingDist->GetTitle();
178  if(fLegend && title.Length() > 0) fLegend->AddEntry(fHist, title, "L");
179 
180  return 1./weightSum;
181 }
182 
183 ////////////////////////////////////////////////////////////////////////////////
184 
185 Double_t SamplingDistPlot::AddSamplingDistributionShaded(const SamplingDistribution *samplingDist, Double_t minShaded, Double_t maxShaded, Option_t *drawOptions) {
186  if( samplingDist->GetSamplingDistribution().empty() ) {
187  coutW(Plotting) << "Empty sampling distribution given to plot. Skipping." << endl;
188  return 0.0;
189  }
190  Double_t scaleFactor = AddSamplingDistribution(samplingDist, drawOptions);
191 
192  TH1F *shaded = (TH1F*)fHist->Clone((string(samplingDist->GetName())+string("_shaded")).c_str());
193  shaded->SetDirectory(0);
194  shaded->SetFillStyle(fFillStyle++);
195  shaded->SetLineWidth(1);
196 
197  for (int i=0; i<shaded->GetNbinsX(); ++i) {
198  if (shaded->GetBinCenter(i) < minShaded || shaded->GetBinCenter(i) > maxShaded){
199  shaded->SetBinContent(i,0);
200  }
201  }
202 
203  TString options(drawOptions);
204  options.ToUpper();
205  if(options.Contains("NORMALIZE")) {
206  options.ReplaceAll("NORMALIZE", "");
207  options.Strip();
208  }
209  addObject(shaded, options.Data());
210 
211  return scaleFactor;
212 }
213 
214 ////////////////////////////////////////////////////////////////////////////////
215 
216 void SamplingDistPlot::AddLine(Double_t x1, Double_t y1, Double_t x2, Double_t y2, const char* title) {
217  TLine *line = new TLine(x1, y1, x2, y2);
218  line->SetLineWidth(3);
219  line->SetLineColor(kBlack);
220 
221  if(fLegend && title) fLegend->AddEntry(line, title, "L");
222 
223  addOtherObject(line, ""); // no options
224 }
225 
226 ////////////////////////////////////////////////////////////////////////////////
227 /// add an histogram (it will be cloned);
228 
229 void SamplingDistPlot::AddTH1(TH1* h, Option_t *drawOptions) {
230  if(fLegend && h->GetTitle()) fLegend->AddEntry(h, h->GetTitle(), "L");
231  TH1 * hcopy = (TH1*) h->Clone();
232  hcopy->SetDirectory(0);
233  addObject(hcopy, drawOptions);
234 }
235 void SamplingDistPlot::AddTF1(TF1* f, const char* title, Option_t *drawOptions) {
236  if(fLegend && title) fLegend->AddEntry(f, title, "L");
237  addOtherObject(f->Clone(), drawOptions);
238 }
239 
240 ////////////////////////////////////////////////////////////////////////////////
241 ///Determine if the sampling distribution has weights and store them
242 
243 void SamplingDistPlot::SetSampleWeights(const SamplingDistribution* samplingDist)
244 {
245  fIsWeighted = kFALSE;
246 
247  if(samplingDist->GetSampleWeights().size() != 0){
248  fIsWeighted = kTRUE;
249  fSampleWeights = samplingDist->GetSampleWeights();
250  }
251 
252  return;
253 }
254 
255 ////////////////////////////////////////////////////////////////////////////////
256 /// Add a generic object to this plot. The specified options will be
257 /// used to Draw() this object later. The caller transfers ownership
258 /// of the object with this call, and the object will be deleted
259 /// when its containing plot object is destroyed.
260 
261 void SamplingDistPlot::addObject(TObject *obj, Option_t *drawOptions)
262 {
263 
264  if(0 == obj) {
265  std::cerr << fName << "::addObject: called with a null pointer" << std::endl;
266  return;
267  }
268 
269  fItems.Add(obj,drawOptions);
270 
271  return;
272 }
273 
274 ////////////////////////////////////////////////////////////////////////////////
275 /// Add a generic object to this plot. The specified options will be
276 /// used to Draw() this object later. The caller transfers ownership
277 /// of the object with this call, and the object will be deleted
278 /// when its containing plot object is destroyed.
279 
280 void SamplingDistPlot::addOtherObject(TObject *obj, Option_t *drawOptions)
281 {
282  if(0 == obj) {
283  oocoutE(this,InputArguments) << fName << "::addOtherObject: called with a null pointer" << std::endl;
284  return;
285  }
286 
287  fOtherItems.Add(obj,drawOptions);
288 
289  return;
290 }
291 
292 ////////////////////////////////////////////////////////////////////////////////
293 /// Draw this plot and all of the elements it contains. The specified options
294 /// only apply to the drawing of our frame. The options specified in our add...()
295 /// methods will be used to draw each object we contain.
296 
297 void SamplingDistPlot::Draw(Option_t * /*options */) {
298  ApplyDefaultStyle();
299 
300  Double_t theMin(0.), theMax(0.), theYMin(NaN), theYMax(0.);
301  GetAbsoluteInterval(theMin, theMax, theYMax);
302  if( !IsNaN(fXMin) ) theMin = fXMin;
303  if( !IsNaN(fXMax) ) theMax = fXMax;
304  if( !IsNaN(fYMin) ) theYMin = fYMin;
305  if( !IsNaN(fYMax) ) theYMax = fYMax;
306 
307  RooRealVar xaxis("xaxis", fVarName.Data(), theMin, theMax);
308 
309  //L.M. by drawing many times we create a memory leak ???
310  if (fRooPlot) delete fRooPlot;
311 
312  bool dirStatus = RooPlot::addDirectoryStatus();
313  // make the RooPlot managed by this class
314  if (dirStatus) RooPlot::setAddDirectoryStatus(false);
315  fRooPlot = xaxis.frame();
316  if (dirStatus) RooPlot::setAddDirectoryStatus(true);
317  if (!fRooPlot) {
318  oocoutE(this,InputArguments) << "invalid variable to plot" << std::endl;
319  return;
320  }
321  fRooPlot->SetTitle("");
322  if( !IsNaN(theYMax) ) {
323  //coutI(InputArguments) << "Setting maximum to " << theYMax << endl;
324  fRooPlot->SetMaximum(theYMax);
325  }
326  if( !IsNaN(theYMin) ) {
327  //coutI(InputArguments) << "Setting minimum to " << theYMin << endl;
328  fRooPlot->SetMinimum(theYMin);
329  }
330 
331  fIterator->Reset();
332  TH1F *obj = 0;
333  while ((obj = (TH1F*) fIterator->Next())) {
334  //obj->Draw(fIterator->GetOption());
335  // add cloned objects to avoid mem leaks
336  TH1 * cloneObj = (TH1*)obj->Clone();
337  if( !IsNaN(theYMax) ) {
338  //coutI(InputArguments) << "Setting maximum of TH1 to " << theYMax << endl;
339  cloneObj->SetMaximum(theYMax);
340  }
341  if( !IsNaN(theYMin) ) {
342  //coutI(InputArguments) << "Setting minimum of TH1 to " << theYMin << endl;
343  cloneObj->SetMinimum(theYMin);
344  }
345  cloneObj->SetDirectory(0); // transfer ownership of the object
346  fRooPlot->addTH1(cloneObj, fIterator->GetOption());
347  }
348 
349  TIterator *otherIt = fOtherItems.MakeIterator();
350  TObject *otherObj = NULL;
351  while ((otherObj = otherIt->Next())) {
352  TObject * cloneObj = otherObj->Clone();
353  fRooPlot->addObject(cloneObj, otherIt->GetOption());
354  }
355  delete otherIt;
356 
357 
358  if(fLegend) fRooPlot->addObject(fLegend);
359 
360  if(bool(gStyle->GetOptLogx()) != fLogXaxis) {
361  if(!fApplyStyle) coutW(Plotting) << "gStyle will be changed to adjust SetOptLogx(...)" << endl;
362  gStyle->SetOptLogx(fLogXaxis);
363  }
364  if(bool(gStyle->GetOptLogy()) != fLogYaxis) {
365  if(!fApplyStyle) coutW(Plotting) << "gStyle will be changed to adjust SetOptLogy(...)" << endl;
366  gStyle->SetOptLogy(fLogYaxis);
367  }
368  fRooPlot->Draw();
369 
370  // apply this since gStyle does not work for RooPlot
371  if (gPad) {
372  gPad->SetLogx(fLogXaxis);
373  gPad->SetLogy(fLogYaxis);
374  }
375 
376  return;
377 }
378 
379 ////////////////////////////////////////////////////////////////////////////////
380 
381 void SamplingDistPlot::ApplyDefaultStyle(void) {
382  if(fApplyStyle) {
383  // use plain black on white colors
384  Int_t icol = 0;
385  gStyle->SetFrameBorderMode( icol );
386  gStyle->SetCanvasBorderMode( icol );
387  gStyle->SetPadBorderMode( icol );
388  gStyle->SetPadColor( icol );
389  gStyle->SetCanvasColor( icol );
390  gStyle->SetStatColor( icol );
391  gStyle->SetFrameFillStyle( 0 );
392 
393  // set the paper & margin sizes
394  gStyle->SetPaperSize( 20, 26 );
395 
396  if(fLegend) {
397  fLegend->SetFillColor(0);
398  fLegend->SetBorderSize(1);
399  }
400  }
401 }
402 
403 ////////////////////////////////////////////////////////////////////////////////
404 
405 void SamplingDistPlot::GetAbsoluteInterval(Double_t &theMin, Double_t &theMax, Double_t &theYMax) const
406 {
407  Double_t tmpmin = TMath::Infinity();
408  Double_t tmpmax = -TMath::Infinity();
409  Double_t tmpYmax = -TMath::Infinity();
410 
411  fIterator->Reset();
412  TH1F *obj = 0;
413  while((obj = (TH1F*)fIterator->Next())) {
414  if(obj->GetXaxis()->GetXmin() < tmpmin) tmpmin = obj->GetXaxis()->GetXmin();
415  if(obj->GetXaxis()->GetXmax() > tmpmax) tmpmax = obj->GetXaxis()->GetXmax();
416  if(obj->GetMaximum() > tmpYmax) tmpYmax = obj->GetMaximum() + 0.1*obj->GetMaximum();
417  }
418 
419  theMin = tmpmin;
420  theMax = tmpmax;
421  theYMax = tmpYmax;
422 
423  return;
424 }
425 
426 ////////////////////////////////////////////////////////////////////////////////
427 /// Sets line color for given sampling distribution and
428 /// fill color for the associated shaded TH1F.
429 
430 void SamplingDistPlot::SetLineColor(Color_t color, const SamplingDistribution *samplDist) {
431  if (samplDist == 0) {
432  fHist->SetLineColor(color);
433 
434  fIterator->Reset();
435  TH1F *obj = 0;
436 
437  TString shadedName(fHist->GetName());
438  shadedName += "_shaded";
439 
440  while ((obj = (TH1F*) fIterator->Next())) {
441  if (!strcmp(obj->GetName(), shadedName.Data())) {
442  obj->SetLineColor(color);
443  obj->SetFillColor(color);
444  //break;
445  }
446  }
447  } else {
448  fIterator->Reset();
449  TH1F *obj = 0;
450 
451  TString shadedName(samplDist->GetName());
452  shadedName += "_shaded";
453 
454  while ((obj = (TH1F*) fIterator->Next())) {
455  if (!strcmp(obj->GetName(), samplDist->GetName())) {
456  obj->SetLineColor(color);
457  //break;
458  }
459  if (!strcmp(obj->GetName(), shadedName.Data())) {
460  obj->SetLineColor(color);
461  obj->SetFillColor(color);
462  //break;
463  }
464  }
465  }
466 
467  return;
468 }
469 
470 ////////////////////////////////////////////////////////////////////////////////
471 
472 void SamplingDistPlot::SetLineWidth(Width_t lwidth, const SamplingDistribution *samplDist)
473 {
474  if(samplDist == 0){
475  fHist->SetLineWidth(lwidth);
476  }
477  else{
478  fIterator->Reset();
479  TH1F *obj = 0;
480  while((obj = (TH1F*)fIterator->Next())) {
481  if(!strcmp(obj->GetName(),samplDist->GetName())){
482  obj->SetLineWidth(lwidth);
483  break;
484  }
485  }
486  }
487 
488  return;
489 }
490 
491 ////////////////////////////////////////////////////////////////////////////////
492 
493 void SamplingDistPlot::SetLineStyle(Style_t style, const SamplingDistribution *samplDist)
494 {
495  if(samplDist == 0){
496  fHist->SetLineStyle(style);
497  }
498  else{
499  fIterator->Reset();
500  TH1F *obj = 0;
501  while((obj = (TH1F*)fIterator->Next())) {
502  if(!strcmp(obj->GetName(),samplDist->GetName())){
503  obj->SetLineStyle(style);
504  break;
505  }
506  }
507  }
508 
509  return;
510 }
511 
512 ////////////////////////////////////////////////////////////////////////////////
513 
514 void SamplingDistPlot::SetMarkerStyle(Style_t style, const SamplingDistribution *samplDist)
515 {
516  if(samplDist == 0){
517  fHist->SetMarkerStyle(style);
518  }
519  else{
520  fIterator->Reset();
521  TH1F *obj = 0;
522  while((obj = (TH1F*)fIterator->Next())) {
523  if(!strcmp(obj->GetName(),samplDist->GetName())){
524  obj->SetMarkerStyle(style);
525  break;
526  }
527  }
528  }
529 
530  return;
531 }
532 
533 ////////////////////////////////////////////////////////////////////////////////
534 
535 void SamplingDistPlot::SetMarkerColor(Color_t color, const SamplingDistribution *samplDist)
536 {
537  if(samplDist == 0){
538  fHist->SetMarkerColor(color);
539  }
540  else{
541  fIterator->Reset();
542  TH1F *obj = 0;
543  while((obj = (TH1F*)fIterator->Next())) {
544  if(!strcmp(obj->GetName(),samplDist->GetName())){
545  obj->SetMarkerColor(color);
546  break;
547  }
548  }
549  }
550 
551  return;
552 }
553 
554 ////////////////////////////////////////////////////////////////////////////////
555 
556 void SamplingDistPlot::SetMarkerSize(Size_t size, const SamplingDistribution *samplDist)
557 {
558  if(samplDist == 0){
559  fHist->SetMarkerSize(size);
560  }
561  else{
562  fIterator->Reset();
563  TH1F *obj = 0;
564  while((obj = (TH1F*)fIterator->Next())) {
565  if(!strcmp(obj->GetName(),samplDist->GetName())){
566  obj->SetMarkerSize(size);
567  break;
568  }
569  }
570  }
571 
572  return;
573 }
574 
575 ////////////////////////////////////////////////////////////////////////////////
576 
577 TH1F* SamplingDistPlot::GetTH1F(const SamplingDistribution *samplDist)
578 {
579  if(samplDist == NULL){
580  return fHist;
581  }else{
582  fIterator->Reset();
583  TH1F *obj = 0;
584  while((obj = (TH1F*)fIterator->Next())) {
585  if(!strcmp(obj->GetName(),samplDist->GetName())){
586  return obj;
587  }
588  }
589  }
590 
591  return NULL;
592 }
593 
594 ////////////////////////////////////////////////////////////////////////////////
595 
596 void SamplingDistPlot::RebinDistribution(Int_t rebinFactor, const SamplingDistribution *samplDist)
597 {
598  if(samplDist == 0){
599  fHist->Rebin(rebinFactor);
600  }
601  else{
602  fIterator->Reset();
603  TH1F *obj = 0;
604  while((obj = (TH1F*)fIterator->Next())) {
605  if(!strcmp(obj->GetName(),samplDist->GetName())){
606  obj->Rebin(rebinFactor);
607  break;
608  }
609  }
610  }
611 
612  return;
613 }
614 
615 ////////////////////////////////////////////////////////////////////////////////
616 /// TODO test
617 
618 void SamplingDistPlot::DumpToFile(const char* RootFileName, Option_t *option, const char *ftitle, Int_t compress) {
619  // All the objects are written to rootfile
620 
621  if(!fRooPlot) {
622  cout << "Plot was not drawn yet. Dump can only be saved after it was drawn with Draw()." << endl;
623  return;
624  }
625 
626  TFile ofile(RootFileName, option, ftitle, compress);
627  ofile.cd();
628  fRooPlot->Write();
629  ofile.Close();
630 }