BayesianAnalysisToolkit  0.9.3
Static Private Member Functions | Private Attributes | Static Private Attributes | List of all members
BCH2D Class Reference

A class for handling 2D distributions. More...

#include <BCH2D.h>

Public Member Functions

Constructors and destructors
 BCH2D (TH2D *h=0)
 
 ~BCH2D ()
 
Member functions (get)
TH2D * GetHistogram ()
 
int GetColor (int index)
 
Member functions (set)
void SetColorScheme (int scheme)
 
void SetHistogram (TH2D *hist)
 
void SetGlobalMode (double mode[2])
 
Member functions (miscellaneous methods)
void Print (const char *filename, std::string options="BTfB3CS1meangmode", std::vector< double > intervals=std::vector< double >(0), int ww=0, int wh=0)
 
void Print (const char *filename, std::string options, double interval, int ww=0, int wh=0)
 
void Draw (std::string options="BTfB3CS1meangmodelmode", std::vector< double > intervals=std::vector< double >(0))
 
void Draw (std::string options, double interval)
 
void CalculateIntegratedHistogram ()
 
void PrintIntegratedHistogram (const char *filename)
 
double GetLevel (double p)
 
double GetArea (double p)
 
std::vector< int > GetNIntervalsY (TH2D *h, int &nfoundmax)
 
TGraph * CalculateProfileGraph (int axis, std::string options="mode")
 
TGraph * DrawProfile (int axis, std::string options, std::string drawoptions="blacksolid")
 
TGraph * DrawProfileX (std::string options, std::string drawoptions)
 
TGraph * DrawProfileY (std::string options, std::string drawoptions)
 

Static Private Member Functions

static unsigned int getNextIndex ()
 

Private Attributes

TH2D * fHistogram
 
TH1D * fIntegratedHistogram
 
double fMode [2]
 
int fModeFlag
 
std::vector< int > fColors
 
std::vector< TObject * > fROOTObjects
 

Static Private Attributes

static unsigned int fHCounter =0
 

Detailed Description

A class for handling 2D distributions.

Author
Daniel Kollar
Kevin Kröninger
Version
1.0
Date
08.2008 This class contains a TH2D histogram and some additional functions. It is used for marginalized distributions.

Definition at line 37 of file BCH2D.h.

Constructor & Destructor Documentation

BCH2D::BCH2D ( TH2D *  h = 0)

The complete constructor.

Definition at line 36 of file BCH2D.cxx.

37  : fHistogram(0)
39  , fModeFlag(0)
40 {
41  if (h)
42  SetHistogram(h);
43 }
BCH2D::~BCH2D ( )

The default destructor.

Definition at line 46 of file BCH2D.cxx.

47 {
49  delete fIntegratedHistogram;
50 
51  for (unsigned i = 0; i < fROOTObjects.size(); ++i)
52  delete fROOTObjects[i];
53 }

Member Function Documentation

void BCH2D::CalculateIntegratedHistogram ( )

Calculates the integral of the distribution as a function of the height.

Definition at line 564 of file BCH2D.cxx.

565 {
566  int nz = 100;
567 
568  double zmin = fHistogram->GetMinimum();
569  double zmax = fHistogram->GetMaximum();
570  double dz = (zmax-zmin);
571 
572  double nx = fHistogram->GetNbinsX();
573  double ny = fHistogram->GetNbinsY();
574 
575  // create histogram
577  delete fIntegratedHistogram;
578 
579  fIntegratedHistogram = new TH1D(
580  TString::Format("%s_int_prob_%d",fHistogram->GetName(),BCLog::GetHIndex()), "", nz, zmin, zmax+0.05*dz);
581  fIntegratedHistogram->SetXTitle("z");
582  fIntegratedHistogram->SetYTitle("Integrated probability");
583  fIntegratedHistogram->SetStats(kFALSE);
584 
585  // loop over histogram
586  for (int ix = 1; ix <= nx; ix++) {
587  for (int iy = 1; iy <= ny; iy++) {
588  double p = fHistogram->GetBinContent(ix, iy);
589  fIntegratedHistogram->Fill(p, p);
590  }
591  }
592 
593  // normalize histogram
594  fIntegratedHistogram->Scale(1.0/fIntegratedHistogram->GetEntries());
595 
596 }
TGraph * BCH2D::CalculateProfileGraph ( int  axis,
std::string  options = "mode" 
)

Return a graph of the profile along x or y. The profile is calculated by scanning through the one axis, e.g., y, and finding the mode, mean, median, etc. with respect to the other axis, e.g. x.

Parameters
axisx-axis (0) or y-axis (1)
optionsOptions:
mode : calculate the profile using the mode [default]
mean : calculate the profile using the mean
median : calculate the profile using the median

Definition at line 691 of file BCH2D.cxx.

692 {
693  // option flags
694  bool flag_mode = true;
695  bool flag_mean = false;
696  bool flag_median = false;
697 
698  // check content of options string
699  if (options.find("mode") < options.size()) {
700  flag_mode = true;
701  }
702 
703  if (options.find("mean") < options.size()) {
704  flag_mean = true;
705  flag_mode = false;
706  flag_median = false;
707  }
708 
709  if (options.find("median") < options.size()) {
710  flag_median = true;
711  flag_mode = false;
712  flag_mean = false;
713  }
714 
715  // define profile graph
716  TGraph* graph_profile = new TGraph();
717 
718  // define limits of running
719  int nx = fHistogram->GetNbinsX();
720  int ny = fHistogram->GetNbinsY();
721 
722  double xmin = fHistogram->GetXaxis()->GetBinLowEdge(1);
723  double xmax = fHistogram->GetXaxis()->GetBinLowEdge(nx+1);
724 
725  double ymin = fHistogram->GetYaxis()->GetBinLowEdge(1);
726  double ymax = fHistogram->GetYaxis()->GetBinLowEdge(nx+1);
727 
728  int nbins_outer = 0;
729  int nbins_inner = 0;
730  double axis_min = 0;
731  double axis_max = 0;
732 
733  if (axis==0) {
734  nbins_outer = nx;
735  nbins_inner = ny;
736  axis_min = ymin;
737  axis_max = ymax;
738  }
739  else {
740  nbins_outer = ny;
741  nbins_inner = nx;
742  axis_min = xmin;
743  axis_max = xmax;
744  }
745 
746  // loop over outer axis of choice
747  for (int ibin_outer = 1; ibin_outer <= nbins_outer ; ibin_outer++) {
748 
749  // copy slice at a fixed value into a 1D histogram
750  TH1D* hist_temp = new TH1D("", "", nbins_inner, axis_min, axis_max);
751 
752  for (int ibin_inner = 1; ibin_inner <= nbins_inner; ibin_inner++) {
753  int ix = 0;
754  int iy = 0;
755  if (axis == 0) {
756  ix = ibin_outer;
757  iy = ibin_inner;
758  }
759  else {
760  ix = ibin_inner;
761  iy = ibin_outer;
762  }
763  double content = fHistogram->GetBinContent(ix, iy);
764  hist_temp->SetBinContent(ibin_inner, content);
765  }
766  // normalize to unity
767  hist_temp->Scale(1.0/hist_temp->Integral());
768 
769  // create BAT histogram
770  BCH1D* bchist_temp = new BCH1D(hist_temp);
771 
772  // calculate (x,y) of new point
773  double x = fHistogram->GetXaxis()->GetBinCenter(ibin_outer);
774  double y = fHistogram->GetYaxis()->GetBinCenter(ibin_outer);
775 
776  double temp = 0;
777 
778  if (flag_mode)
779  temp = bchist_temp->GetMode();
780  else if (flag_mean)
781  temp = bchist_temp->GetMean();
782  else if (flag_median)
783  temp = bchist_temp->GetMedian();
784 
785  if (axis == 0)
786  y = temp;
787  else
788  x = temp;
789 
790  // add new point to graph
791  graph_profile->SetPoint(ibin_outer-1, x, y);
792 
793  // clean up
794  delete bchist_temp;
795  }
796 
797  // return the graph
798  return graph_profile;
799 }
void BCH2D::Draw ( std::string  options = "BTfB3CS1meangmodelmode",
std::vector< double >  intervals = std::vector<double>(0) 
)

Draw distribution into the active canvas.

Parameters
optionsDrawing options:
BTf : band type a filled area [default]
BTc : band type is a contour
B1 : draw one band corresponding to an integrated probability specified in intervals
B2 : draw two bands corresponding to an integrated probability specified in intervals
B3 : draw three bands corresponding to an integrated probability specified in intervals [default]
CS0 : choose color scheme 0 (B&W)
CS1 : choose color scheme 1 (green/yellow/red) [default]
CS2 : choose color scheme 2 (blueish colors)
CS3 : choose color scheme 3 (redish colors)
smooth1 : use ROOT smoothing algorithm once
smooth3 : use ROOT smoothing algorithm three times
smooth5 : use ROOT smoothing algorithm five times
smooth10 : use ROOT smoothing algorithm ten times
mean : draw mean value and standard deviation [default]
gmode : draw global mode [default]
lmode : draw global mode [default]
profilex : draw the profile line vs. x using the mode
profiley : draw the profile line vs. y using the mode
nL : remove legend
intervalsthe intervals for the bands

Definition at line 171 of file BCH2D.cxx.

172 {
173  // option flags
174  bool flag_legend = true;
175  bool flag_mode_global = false;
176  bool flag_mode_local = false;
177  bool flag_mean = false;
178  bool flag_smooth1 = false;
179  bool flag_smooth3 = false;
180  bool flag_smooth5 = false;
181  bool flag_smooth10 = false;
182  bool flag_profilex = false;
183  bool flag_profiley = false;
184 
185  // band type
186  int bandtype = 0;
187 
188  // number of bands
189  int nbands = 1; // number of shaded bands
190  intervals.push_back(0.6827);
191 
192  // define draw options called in TH1D::Draw(...)
193  std::string draw_options = "COLZ";
194 
195  // check content of options string
196  if (options.find("smooth1") < options.size()) {
197  flag_smooth1 = true;
198  }
199 
200  if (options.find("smooth3") < options.size()) {
201  flag_smooth3 = true;
202  }
203 
204  if (options.find("smooth5") < options.size()) {
205  flag_smooth5 = true;
206  }
207 
208  if (options.find("smooth10") < options.size()) {
209  flag_smooth10 = true;
210  }
211 
212  if (options.find("nL") < options.size()) {
213  flag_legend = false;
214  }
215 
216  if (options.find("BTf") < options.size()) {
217  bandtype = 0;
218  }
219  else if (options.find("BTc") < options.size()) {
220  bandtype = 1;
221  }
222  else {
223  bandtype = 0;
224  }
225 
226  if (options.find("profilex") < options.size()) {
227  flag_profilex = true;
228  }
229 
230  if (options.find("profiley") < options.size()) {
231  flag_profiley = true;
232  }
233 
234  if (options.find("gmode") < options.size()) {
235  if (fModeFlag)
236  flag_mode_global = true;
237  }
238 
239  if (options.find("lmode") < options.size()) {
240  flag_mode_local = true;
241  }
242 
243  if (options.find("mean") < options.size()) {
244  flag_mean = true;
245  }
246 
247  if (options.find("B1") < options.size()) {
248  nbands = 1;
249  if (intervals.size() != 1) {
250  intervals.clear();
251  intervals.push_back(0.6827);
252  }
253  }
254 
255  if (options.find("B2") < options.size()) {
256  nbands = 2;
257  if (intervals.size() != 2) {
258  intervals.clear();
259  intervals.push_back(0.6827);
260  intervals.push_back(0.9545);
261  }
262  }
263 
264  if (options.find("B3") < options.size()) {
265  nbands = 3;
266  if (intervals.size() != 3) {
267  intervals.clear();
268  intervals.push_back(0.6827);
269  intervals.push_back(0.9545);
270  intervals.push_back(0.9973);
271  }
272  }
273 
274  if (options.find("CS0") < options.size()) {
275  SetColorScheme(0);
276  }
277  else if (options.find("CS1") < options.size()) {
278  SetColorScheme(1);
279  }
280  else if (options.find("CS2") < options.size()) {
281  SetColorScheme(2);
282  }
283  else if (options.find("CS3") < options.size()) {
284  SetColorScheme(3);
285  }
286  else {
287  SetColorScheme(1);
288  }
289 
290  // prepare size of histogram
291  double xmin = fHistogram->GetXaxis()->GetXmin();
292  double xmax = fHistogram->GetXaxis()->GetXmax();
293  double ymin = fHistogram->GetYaxis()->GetXmin();
294  double ymaxhist = fHistogram->GetYaxis()->GetXmax();
295  double ymax = ymaxhist;
296 
297  // prepare legend
298  TLegend* legend = new TLegend();
299  legend->SetBorderSize(0);
300  legend->SetFillColor(kWhite);
301  legend->SetTextAlign(12);
302  legend->SetTextFont(62);
303  legend->SetTextSize(0.03);
304 
305  // add legend to list of objects
306  fROOTObjects.push_back(legend);
307 
308  // copy histograms for bands
309  TH2D* hist_band = new TH2D(*fHistogram);
310 
311  // add hist_band to list of ROOT objects
312  fROOTObjects.push_back(hist_band);
313 
314  // calculate integrated histogram
316 
317  for (int ix = 1; ix <= fHistogram->GetNbinsX(); ++ix) {
318  for (int iy = 1; iy <= fHistogram->GetNbinsY(); ++iy) {
319  double p = fHistogram->GetBinContent(ix, iy);
320  hist_band->SetBinContent(ix, iy, p);
321  }
322  }
323 
324  // define levels and colors
325  std::vector<double> levels(nbands+2);
326  levels[0] = 0.;
327 
328  std::vector<int> colors(nbands+1);
329  colors[0] = kWhite;
330 
331  for (int i = 1; i <= nbands; ++i) {
332  levels[i] = GetLevel((1.-intervals[nbands-i]));
333  colors[i] = GetColor(nbands-i);
334  }
335  levels[nbands+1] = fIntegratedHistogram->GetXaxis()->GetXmax();
336 
337  // set contour
338  hist_band->SetContour(nbands+2, &levels[0]);
339 
340  // set colors
341  gStyle->SetPalette(nbands+1, &colors[0]);
342 
343  // add hist_band to legend
344  for (int i = 0; i < nbands; ++i) {
345  if (bandtype == 0) {
346  TLegendEntry* le = legend->AddEntry((TObject*)0, Form("smallest %.1f%% interval(s)", intervals[nbands-1-i]*100), "F");
347  le->SetFillColor(GetColor(nbands-1-i));
348  le->SetFillStyle(1001);
349  le->SetLineColor(GetColor(nbands-1-i));
350  le->SetTextAlign(12);
351  le->SetTextFont(62);
352  le->SetTextSize(0.03);
353  }
354  else if (bandtype == 1) {
355  TLegendEntry* le = legend->AddEntry((TObject*)0, Form("smallest %.1f%% interval(s)", intervals[nbands-1-i]*100), "F");
356  fROOTObjects.push_back(le);
357  le->SetLineColor(GetColor(nbands-1-i));
358  le->SetTextAlign(12);
359  le->SetTextFont(62);
360  le->SetTextSize(0.03);
361  }
362  }
363 
364  // mean, mode, median
365  TMarker* marker_mode_global = new TMarker(fMode[0], fMode[1], 24);
366  marker_mode_global->SetMarkerColor(GetColor(4));
367  marker_mode_global->SetMarkerSize(1.5);
368 
369  int binx, biny, binz;
370  fHistogram->GetBinXYZ( fHistogram->GetMaximumBin(), binx, biny, binz);
371  TMarker* marker_mode_local = new TMarker(fHistogram->GetXaxis()->GetBinCenter(binx), fHistogram->GetYaxis()->GetBinCenter(biny), 25);
372  marker_mode_local->SetMarkerColor(GetColor(4));
373  marker_mode_local->SetMarkerSize(1.5);
374 
375  double xmean = fHistogram->GetMean(1);
376  double ymean = fHistogram->GetMean(2);
377  double xrms = fHistogram->GetRMS(1);
378  double yrms = fHistogram->GetRMS(2);
379 
380  TMarker* marker_mean = new TMarker(xmean, ymean, 32);
381  marker_mean->SetMarkerColor(GetColor(4));
382  marker_mean->SetMarkerSize(1.5);
383 
384  // standard deviation
385  TArrow* arrow_std1 = new TArrow(xmean-xrms, ymean,
386  xmean+xrms, ymean,
387  0.02, "<|>");
388  arrow_std1->SetLineColor(GetColor(4));
389  arrow_std1->SetFillColor(GetColor(4));
390 
391  TArrow* arrow_std2 = new TArrow(xmean, ymean-yrms,
392  xmean, ymean+yrms,
393  0.02, "<|>");
394  arrow_std2->SetLineColor(GetColor(4));
395  arrow_std2->SetFillColor(GetColor(4));
396 
397  // add marker_mean and arrow_std to list of ROOT objects
398  fROOTObjects.push_back(marker_mean);
399  fROOTObjects.push_back(marker_mode_global);
400  fROOTObjects.push_back(marker_mode_local);
401  fROOTObjects.push_back(arrow_std1);
402  fROOTObjects.push_back(arrow_std2);
403 
404  if (flag_mode_global) {
405  TLegendEntry* le = legend->AddEntry(marker_mode_global, "global mode", "P");
406  le->SetMarkerStyle(24);
407  le->SetMarkerSize(1.5);
408  le->SetMarkerColor(GetColor(4));
409  }
410 
411  if (flag_mode_local) {
412  TLegendEntry* le = legend->AddEntry(marker_mode_local, "local mode", "P");
413  le->SetMarkerStyle(25);
414  le->SetMarkerSize(1.5);
415  le->SetMarkerColor(GetColor(4));
416  }
417 
418  if (flag_mean) {
419  TLegendEntry* le = legend->AddEntry(arrow_std1, "mean and standard deviation", "PL");
420  le->SetLineColor(GetColor(4));
421  le->SetMarkerStyle(32);
422  le->SetMarkerSize(1.5);
423  le->SetMarkerColor(GetColor(4));
424  }
425 
426  TGraph* graph_profilex = 0;
427  TGraph* graph_profiley = 0;
428 
429  if (flag_profilex) {
430  TLegendEntry* le = legend->AddEntry(graph_profilex, "profile x", "L");
431  le->SetLineStyle(2);
432  }
433 
434  if (flag_profiley) {
435  TLegendEntry* le = legend->AddEntry(graph_profiley, "profile y", "L");
436  le->SetLineStyle(3);
437  }
438 
439  // calculate legend height in NDC coordinates
440  double height = 0.03*legend->GetNRows();
441 
442  // make room for legend
443  if (flag_legend)
444  ymax+=(ymax-ymin)*(0.1+height);
445 
446  double deltax = 0.0015*(xmax - xmin);
447  double deltay = 0.0015*(ymaxhist - ymin);
448  TH2D* hist_axes = new TH2D("", "", 1, xmin-deltax, xmax+deltax, 1, ymin-deltay, ymaxhist+deltay);
449  hist_axes->SetXTitle(fHistogram->GetXaxis()->GetTitle());
450  hist_axes->SetYTitle(fHistogram->GetYaxis()->GetTitle());
451  hist_axes->SetLineWidth(fHistogram->GetLineWidth());
452  hist_axes->SetStats(kFALSE);
453  fROOTObjects.push_back(hist_axes);
454 
455  // draw axes
456  hist_axes->Draw("COL");
457 
458  // smooth
459  if (flag_smooth1) {
460  fHistogram->Smooth(1);
461  fHistogram->Scale(1.0/fHistogram->Integral("width"));
462  }
463  if (flag_smooth3) {
464  fHistogram->Smooth(3);
465  fHistogram->Scale(1.0/fHistogram->Integral("width"));
466  }
467  if (flag_smooth5) {
468  fHistogram->Smooth(5);
469  fHistogram->Scale(1.0/fHistogram->Integral("width"));
470  }
471  if (flag_smooth10) {
472  fHistogram->Smooth(10);
473  fHistogram->Scale(1.0/fHistogram->Integral("width"));
474  }
475 
476  // draw histogram
477  if (bandtype == 0)
478  hist_band->Draw("COL SAME");
479  else if (bandtype == 1)
480  hist_band->Draw("CONT1 SAME");
481 
482  // draw profiles
483  if (flag_profilex) {
484  graph_profilex = DrawProfileX("mode", "dashed");
485 
486  // add graph to list of objects
487  fROOTObjects.push_back(graph_profilex);
488  }
489 
490  if (flag_profiley) {
491  graph_profiley = DrawProfileY("mode", "dotted");
492 
493  // add graph to list of objects
494  fROOTObjects.push_back(graph_profiley);
495  }
496 
497  // mean, mode, median
498  if (flag_mode_global) {
499  marker_mode_global->Draw();
500  }
501 
502  if (flag_mode_local) {
503  marker_mode_local->Draw();
504  }
505 
506  if (flag_mean) {
507  arrow_std1->Draw();
508  arrow_std2->Draw();
509  marker_mean->Draw();
510  }
511 
512  if (flag_legend)
513  gPad->SetTopMargin(0.02);
514 
515  double xlegend1 = gPad->GetLeftMargin();
516  double xlegend2 = 1.0-gPad->GetRightMargin();
517  double ylegend1 = 1.-gPad->GetTopMargin()-height;
518  double ylegend2 = 1.-gPad->GetTopMargin();
519 
520  // place legend on top of histogram
521  legend->SetX1NDC(xlegend1);
522  legend->SetX2NDC(xlegend2);
523  legend->SetY1NDC(ylegend1);
524  legend->SetY2NDC(ylegend2);
525 
526  // draw legend
527  if (flag_legend) {
528  legend->Draw();
529  }
530 
531  // draw axes again
532  hist_axes->Draw("SAME");
533 
534  // rescale
535  gPad->SetTopMargin(1.-ylegend1+0.01);
536 
537  gPad->RedrawAxis();
538 
539  return;
540 }
void BCH2D::Draw ( std::string  options,
double  interval 
)

Draw distribution into the active canvas.

Parameters
optionsDrawing options,
See Also
Draw(std::string options="BTfB3CS1meangmodelmode", std::vector<double> intervals=std::vector<double>(0))
Parameters
intervalan upper or lower limit
See Also
Draw(std::string options="BTfB3CS1meangmodelmode", std::vector<double> intervals=std::vector<double>(0))

Definition at line 543 of file BCH2D.cxx.

544 {
545  std::vector<double> tempvec;
546  tempvec.push_back(interval);
547  Draw(options, tempvec);
548 }
TGraph * BCH2D::DrawProfile ( int  axis,
std::string  options,
std::string  drawoptions = "blacksolid" 
)

Draw the profile along x or y

Parameters
optionsThe options of the calculation of the profile
drawoptionsThe drawing options:
black : line color is black [default]
red : line color is red
solid : line is solid [default]
dashed : line is dashed
dotted : line is dotted
Returns
The profile graph

Definition at line 802 of file BCH2D.cxx.

803 {
804  // option flags
805  bool flag_black = false;
806  bool flag_red = false;
807  bool flag_solid = false;
808  bool flag_dashed = false;
809  bool flag_dotted = false;
810 
811  // check content of options string
812  if (drawoptions.find("black") < options.size()) {
813  flag_black = true;
814  }
815 
816  else if (drawoptions.find("red") < options.size()) {
817  flag_red = true;
818  }
819 
820  else {
821  flag_black = true;
822  }
823 
824  if (drawoptions.find("solid") < options.size()) {
825  flag_solid = true;
826  }
827 
828  else if (drawoptions.find("dashed") < options.size()) {
829  flag_dashed = true;
830  }
831 
832  else if (drawoptions.find("dotted") < options.size()) {
833  flag_dotted = true;
834  }
835 
836  else {
837  flag_solid = true;
838  }
839 
840  // get the profile graph
841  TGraph* graph_profile = CalculateProfileGraph(axis, options);
842 
843  // change drawing options
844  if (flag_black)
845  graph_profile->SetLineColor(kBlack);
846 
847  if (flag_red)
848  graph_profile->SetLineColor(kRed);
849 
850  if (flag_solid)
851  graph_profile->SetLineStyle(1);
852 
853  if (flag_dashed)
854  graph_profile->SetLineStyle(2);
855 
856  if (flag_dotted)
857  graph_profile->SetLineStyle(3);
858 
859  // draw
860  graph_profile->Draw("L");
861 
862  // return graph
863  return graph_profile;
864 }
TGraph* BCH2D::DrawProfileX ( std::string  options,
std::string  drawoptions 
)
inline

Draw the profile along x

Parameters
optionsThe options defined in DrawProfile
drawoptionsThe drawoptions defined DrawProfile
Returns
The profile graph

Definition at line 218 of file BCH2D.h.

219  { return DrawProfile(0, options, drawoptions); };
TGraph* BCH2D::DrawProfileY ( std::string  options,
std::string  drawoptions 
)
inline

Draw the profile along y

Parameters
optionsThe options defined in DrawProfile
drawoptionsThe drawoptions defined DrawProfile
Returns
The profile graph

Definition at line 226 of file BCH2D.h.

227  { return DrawProfile(1, options, drawoptions); };
double BCH2D::GetArea ( double  p)

Calculate the smallest area over which the integral of the function is p

Parameters
pThe integrated probability in the region below the height to be estimated.

Definition at line 611 of file BCH2D.cxx.

612 {
613  // copy histograms for bands
614  TH2D hist_temp(*fHistogram);
615 
616  // calculate total number of bins
617  int nbins = hist_temp.GetNbinsX() * hist_temp.GetNbinsY();
618 
619  // area
620  double area = 0;
621 
622  // integrated probability
623  double intp = 0;
624 
625  // a counter
626  int counter = 0;
627 
628  // loop over bins
629  while ( (intp < p) && (counter < nbins) ) {
630 
631  // find maximum bin
632  int binx;
633  int biny;
634  int binz;
635  hist_temp.GetBinXYZ(hist_temp.GetMaximumBin(), binx, biny, binz);
636 
637  // increase probability
638  double dp = hist_temp.GetBinContent(binx, biny);
639  intp += dp * hist_temp.GetXaxis()->GetBinWidth(binx) * hist_temp.GetYaxis()->GetBinWidth(biny);
640 
641  // reset maximum bin
642  hist_temp.SetBinContent(binx, biny, 0.);
643 
644  // increase area
645  area += hist_temp.GetXaxis()->GetBinWidth(binx) * hist_temp.GetYaxis()->GetBinWidth(biny);
646 
647  // increase counter
648  counter++;
649  }
650 
651  // return area
652  return area;
653 }
int BCH2D::GetColor ( int  index)
inline

Returns a color of the current color scheme.

Parameters
indexthe color index
Returns
the color number.

Definition at line 67 of file BCH2D.h.

68  { return fColors.at(index); };
TH2D* BCH2D::GetHistogram ( )
inline

Return the TH2D histogram

Returns
The TH2D histogram.

Definition at line 60 of file BCH2D.h.

61  { return fHistogram; };
double BCH2D::GetLevel ( double  p)

Calculates the height below which the integrated probability has a certain value.

Parameters
pThe integrated probability in the region below the height to be estimated.

Definition at line 599 of file BCH2D.cxx.

600 {
601  double quantiles[1];
602  double probsum[1];
603  probsum[0] = p;
604 
605  fIntegratedHistogram->GetQuantiles( 1, quantiles, probsum);
606 
607  return quantiles[0];
608 }
static unsigned int BCH2D::getNextIndex ( )
inlinestaticprivate

Helper method to get an unique number to be used in histogram naming

Definition at line 258 of file BCH2D.h.

259  { return ++fHCounter; }
std::vector< int > BCH2D::GetNIntervalsY ( TH2D *  h,
int &  nfoundmax 
)

Returns the number of intervals as a function of x

Parameters
hThe histogram.
nfoundmaxThe maximum number of intervals.
Returns
A vector containing the number of intervals for all bins in x.

Definition at line 656 of file BCH2D.cxx.

657 {
658  std::vector<int> nint;
659 
660  int nx = h->GetNbinsX();
661  int ny = h->GetNbinsY();
662 
663  nfoundmax=0;
664 
665  // loop over histogram bins in x
666  for (int ix=1; ix<=nx; ix++)
667  {
668  int nfound=0;
669 
670  // loop over histogram bins in y
671  // count nonzero intervals in y
672  for (int iy=1; iy<=ny; iy++)
673  if(h->GetBinContent(ix,iy)>0.)
674  {
675  while(h->GetBinContent(ix,++iy)>0.)
676  ;
677  nfound++;
678  }
679 
680  // store maximum number of nonzero intervals for the histogram
681  if(nfound>nfoundmax)
682  nfoundmax=nfound;
683 
684  nint.push_back(nfound);
685  }
686 
687  return nint;
688 }
void BCH2D::Print ( const char *  filename,
std::string  options = "BTfB3CS1meangmode",
std::vector< double >  intervals = std::vector<double>(0),
int  ww = 0,
int  wh = 0 
)

Print distribution into a PostScript file.

Parameters
filenameOutput filename
optionthe draw options (see myDraw()) plus
intervalsthe intervals for the bands logz : draw z-axis in log-scale
R : rescale canvas to have a squared histogram
wwcanvas size in pixels along X
wwcanvas size in pixels along Y If ww and wh are set to 0, default ROOT canvas size is used.

Definition at line 107 of file BCH2D.cxx.

108 {
109  // option flags
110  bool flag_logz = false;
111  bool flag_rescale = false;
112 
113  // check content of options string
114  if (options.find("logz") < options.size()) {
115  flag_logz = true;
116  }
117 
118  if (options.find("R") < options.size()) {
119  flag_rescale = true;
120  }
121 
122  // create temporary canvas
123  TCanvas * c;
124  unsigned int cindex = getNextIndex();
125  if(ww > 0 && wh > 0)
126  c = new TCanvas(TString::Format("c_bch2d_%d",cindex), TString::Format("c_bch2d_%d",cindex), ww, wh);
127  else
128  c = new TCanvas(TString::Format("c_bch2d_%d",cindex));
129 
130  // add c to list of objects
131  fROOTObjects.push_back(c);
132 
133  // set log axis
134  if (flag_logz) {
135  c->SetLogz();
136  }
137 
138  // draw histogram
139  Draw(options, intervals);
140 
141  if (flag_rescale) {
142  double top = gPad->GetTopMargin();
143  double bottom = gPad->GetBottomMargin();
144  double left = gPad->GetLeftMargin();
145  double right = gPad->GetRightMargin();
146 
147  double dx = 1.-right - left;
148  double dy = 1.-top-bottom;
149  double ratio = dy/dx;
150  double ynew = c->GetWindowWidth()/ratio;
151  c->SetCanvasSize(c->GetWindowWidth(), (int) ynew);
152  gPad->RedrawAxis();
153 
154  c->Modified();
155  c->Update();
156  }
157 
158  // print to file
159  c->Print(filename);
160 }
void BCH2D::Print ( const char *  filename,
std::string  options,
double  interval,
int  ww = 0,
int  wh = 0 
)

Print distribution into a PostScript file.

Parameters
filenameOutput filename
optionthe draw options,
See Also
Print(const char* filename, std::string options="BTfB3CS1meangmode", std::vector<double> intervals=std::vector<double>(0), int ww=0, int wh=0)
Parameters
intervalan upper or lower limit
wwcanvas size in pixels along X
wwcanvas size in pixels along Y
See Also
Print(const char* filename, std::string options="BTfB3CS1meangmode", std::vector<double> intervals=std::vector<double>(0), int ww=0, int wh=0)

Definition at line 163 of file BCH2D.cxx.

164 {
165  std::vector<double> tempvec;
166  tempvec.push_back(interval);
167  Print(filename, options, tempvec, ww, wh);
168 }
void BCH2D::PrintIntegratedHistogram ( const char *  filename)

Print the integrated histogram.

Parameters
filenamethe name of the file.

Definition at line 551 of file BCH2D.cxx.

552 {
553  TCanvas c;
554  c.Flush();
555 
556  c.cd();
557  c.SetLogy();
559  fIntegratedHistogram->Draw();
560  c.Print(filename);
561 }
void BCH2D::SetColorScheme ( int  scheme)

Sets the color scheme.

Parameters
schemethe scheme index
0 : black and white 1 : yellow-green-red 2 : blueish colors 2 : redish colors 2 : blueish colors

Definition at line 56 of file BCH2D.cxx.

57 {
58  fColors.clear();
59 
60  // color numbering
61  // 0,1,2 : intervals
62  // 3 : quantile lines
63  // 4 : mean, mode, median
64 
65  if (scheme == 0) {
66  fColors.push_back(18);
67  fColors.push_back(16);
68  fColors.push_back(14);
69  fColors.push_back(kBlack);
70  fColors.push_back(kBlack);
71  }
72  else if (scheme == 1) {
73  fColors.push_back(kGreen);
74  fColors.push_back(kYellow);
75  fColors.push_back(kRed);
76  fColors.push_back(kBlack);
77  fColors.push_back(kBlack);
78  }
79  else if (scheme == 2) {
80  fColors.push_back(kBlue+4);
81  fColors.push_back(kBlue+2);
82  fColors.push_back(kBlue);
83  fColors.push_back(kOrange);
84  fColors.push_back(kOrange);
85  }
86  else if (scheme == 3) {
87  fColors.push_back(kRed+4);
88  fColors.push_back(kRed+2);
89  fColors.push_back(kRed);
90  fColors.push_back(kGreen);
91  fColors.push_back(kGreen);
92  }
93  else {
94  SetColorScheme(1);
95  }
96 }
void BCH2D::SetGlobalMode ( double  mode[2])
inline

Set global mode.

Parameters
Theglobal mode.

Definition at line 92 of file BCH2D.h.

93  { fMode[0] = mode[0];
94  fMode[1] = mode[1];
95  fModeFlag =1; };
void BCH2D::SetHistogram ( TH2D *  hist)

Set the 2D histogram.

Definition at line 99 of file BCH2D.cxx.

100 {
101  fHistogram = hist;
102  if (fHistogram and fHistogram->Integral() > 0)
103  fHistogram->Scale(1.0/fHistogram->Integral("width"));
104 }

Member Data Documentation

std::vector<int> BCH2D::fColors
private

The colors of the color scheme.

Definition at line 251 of file BCH2D.h.

unsigned int BCH2D::fHCounter =0
staticprivate

helper variable to get an unique number to be used in histogram naming

Definition at line 262 of file BCH2D.h.

TH2D* BCH2D::fHistogram
private

The 2D histogram

Definition at line 227 of file BCH2D.h.

TH1D* BCH2D::fIntegratedHistogram
private

The integrated 2D histogram

Definition at line 239 of file BCH2D.h.

double BCH2D::fMode[2]
private

Global mode

Definition at line 243 of file BCH2D.h.

int BCH2D::fModeFlag
private

"Is there a global mode?" flag

Definition at line 247 of file BCH2D.h.

std::vector<TObject*> BCH2D::fROOTObjects
private

The colors of the color scheme.

Definition at line 255 of file BCH2D.h.


The documentation for this class was generated from the following files: