Logo ROOT   6.14/05
Reference Guide
HybridPlot.h
Go to the documentation of this file.
1 // @(#)root/roostats:$Id$
2 
3 /*************************************************************************
4  * Project: RooStats *
5  * Package: RooFit/RooStats *
6  * Authors: *
7  * Kyle Cranmer, Lorenzo Moneta, Gregory Schott, Wouter Verkerke *
8  *************************************************************************
9  * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
10  * All rights reserved. *
11  * *
12  * For the licensing terms see $ROOTSYS/LICENSE. *
13  * For the list of contributors see $ROOTSYS/README/CREDITS. *
14  *************************************************************************/
15 
16 #ifndef ROOSTATS_HybridPlot
17 #define ROOSTATS_HybridPlot
18 
19 #include <vector>
20 #include <iostream>
21 
22 #include "TNamed.h"
23 
24 // these should be maybe forward declared
25 // by moving implementations in source file
26 
27 #include "TH1.h"
28 
29 class TLine;
30 class TLegend;
31 class TH1;
32 class TVirtualPad;
33 
34 namespace RooStats {
35 
36  class HybridPlot : public TNamed {
37 
38  public:
39 
40  /// Constructor
41  HybridPlot(const char* name,
42  const char* title,
43  const std::vector<double> & sb_values,
44  const std::vector<double> & b_values,
45  double testStat_data,
46  int n_bins,
47  bool verbosity=true);
48 
49  /// Destructor
50  ~HybridPlot();
51 
52  /// Draw on current pad
53  void Draw (const char* options="");
54 
55  /// All the objects are written to rootfile
56  void DumpToFile (const char* RootFileName, const char* options);
57 
58  /// Get B histo mean
59  double GetBmean(){return fB_histo->GetMean();};
60 
61  /// Get B histo RMS
62  double GetBrms(){return fB_histo->GetRMS();};
63 
64  /// Get B histo
65  TH1F * GetBhisto(){return fB_histo;}
66 
67  /// Get B histo center
68  double GetBCenter(double n_sigmas=1, bool display=false)
69  {return GetHistoCenter(fB_histo,n_sigmas,display);};
70 
71  /// Get B histo integration extremes to obtain the requested area fraction
72  double* GetBIntExtremes(double frac)
73  {return GetHistoPvals(fB_histo,frac);};
74 
75  /// Get SB histo mean
76  double GetSBmean(){return fSb_histo->GetMean();};
77 
78  /// Get SB histo center
79  double GetSBCenter(double n_sigmas=1, bool display=false)
80  {return GetHistoCenter(fSb_histo,n_sigmas,display);};
81 
82  /// Get SB histo RMS
83  double GetSBrms(){return fSb_histo->GetRMS();};
84 
85  /// Get SB histo integration extremes to obtain the requested area fraction
86  double* GetSBIntExtremes(double frac)
87  {return GetHistoPvals(fSb_histo,frac);};
88 
89  /// Get B histo
91 
92  /// Get the pad (or canvas) where it has been drawn
93  TVirtualPad * GetCanvas() { return fPad; }
94 
95  /// Write an image on disk
96  void DumpToImage (const char* filename);
97 
98 
99  /// Get the center of the histo
100  double GetHistoCenter(TH1* histo, double n_rms=1,bool display_result=false);
101 
102  /// Get the "effective sigmas" of the histo
103  double* GetHistoPvals (TH1* histo, double percentage);
104 
105  /// Get the median of an histogram
106  double GetMedian(TH1* histo);
107 
108  private:
109 
110  TH1F* fSb_histo; // The sb Histo
111  TH1F* fSb_histo_shaded; // The sb Histo shaded
112  TH1F* fB_histo; // The b Histo
113  TH1F* fB_histo_shaded; // The b Histo shaded
114  TLine* fData_testStat_line; // The line for the data value of the test statistic
115  TLegend* fLegend; // The legend of the plot
116  TVirtualPad * fPad; // The pad where it has been drawn
117  bool fVerbose; // verbosity flag
118 
119  ClassDef(HybridPlot,1) // Provides the plots for an HybridResult
120  };
121 }
122 
123 #endif
double * GetBIntExtremes(double frac)
Get B histo integration extremes to obtain the requested area fraction.
Definition: HybridPlot.h:72
This class displays a legend box (TPaveText) containing several legend entries.
Definition: TLegend.h:23
double * GetHistoPvals(TH1 *histo, double percentage)
Get the "effective sigmas" of the histo.
Definition: HybridPlot.cxx:302
~HybridPlot()
Destructor.
Definition: HybridPlot.cxx:124
double GetBCenter(double n_sigmas=1, bool display=false)
Get B histo center.
Definition: HybridPlot.h:68
virtual Double_t GetMean(Int_t axis=1) const
For axis = 1,2 or 3 returns the mean value of the histogram along X,Y or Z axis.
Definition: TH1.cxx:6930
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:567
void Draw(const char *options="")
Draw on current pad.
Definition: HybridPlot.cxx:138
Double_t GetRMS(Int_t axis=1) const
Definition: TH1.h:310
void DumpToImage(const char *filename)
Write an image on disk.
Definition: HybridPlot.cxx:226
#define ClassDef(name, id)
Definition: Rtypes.h:320
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:29
This class provides the plots for the result of a study performed with the HybridCalculatorOriginal c...
Definition: HybridPlot.h:36
double GetSBCenter(double n_sigmas=1, bool display=false)
Get SB histo center.
Definition: HybridPlot.h:79
TVirtualPad is an abstract base class for the Pad and Canvas classes.
Definition: TVirtualPad.h:49
TVirtualPad * fPad
Definition: HybridPlot.h:116
double GetSBrms()
Get SB histo RMS.
Definition: HybridPlot.h:83
HybridPlot(const char *name, const char *title, const std::vector< double > &sb_values, const std::vector< double > &b_values, double testStat_data, int n_bins, bool verbosity=true)
Constructor.
Definition: HybridPlot.cxx:42
void DumpToFile(const char *RootFileName, const char *options)
All the objects are written to rootfile.
Definition: HybridPlot.cxx:198
A simple line.
Definition: TLine.h:23
TLine * fData_testStat_line
Definition: HybridPlot.h:114
double GetSBmean()
Get SB histo mean.
Definition: HybridPlot.h:76
Namespace for the RooStats classes.
Definition: Asimov.h:20
double * GetSBIntExtremes(double frac)
Get SB histo integration extremes to obtain the requested area fraction.
Definition: HybridPlot.h:86
double GetMedian(TH1 *histo)
Get the median of an histogram.
Definition: HybridPlot.cxx:352
The TH1 histogram class.
Definition: TH1.h:56
double GetHistoCenter(TH1 *histo, double n_rms=1, bool display_result=false)
Get the center of the histo.
Definition: HybridPlot.cxx:239
double GetBrms()
Get B histo RMS.
Definition: HybridPlot.h:62
double GetBmean()
Get B histo mean.
Definition: HybridPlot.h:59
TH1F * GetSBhisto()
Get B histo.
Definition: HybridPlot.h:90
TVirtualPad * GetCanvas()
Get the pad (or canvas) where it has been drawn.
Definition: HybridPlot.h:93
TH1F * GetBhisto()
Get B histo.
Definition: HybridPlot.h:65
char name[80]
Definition: TGX11.cxx:109