Logo ROOT   6.08/07
Reference Guide
SamplingDistPlot.h
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 #ifndef ROOSTATS_SamplingDistPlot
13 #define ROOSTATS_SamplingDistPlot
14 
15 #include "RooList.h"
16 #include "RooPrintable.h"
17 #include "TNamed.h"
18 #include "TIterator.h"
19 #include "TH1F.h"
20 #include "TF1.h"
21 #include "TLegend.h"
22 
23 #ifndef ROOSTATS_SamplingDistribution
25 #endif
26 
27 #ifndef ROO_PLOT
28 #include "RooPlot.h"
29 #endif
30 
31 
32 namespace RooStats {
33 
34  /**
35 
36  \ingroup Roostats
37 
38 
39  This class provides simple and straightforward utilities to plot SamplingDistribution
40  objects.
41  */
42 
43 
44  class SamplingDistPlot : public TNamed, public RooPrintable {
45 
46  public:
47  /// Constructors for SamplingDistribution
50 // SamplingDistPlot(const char* name, const char* title, Int_t nbins, Double_t xmin, Double_t xmax);
51 
52  /// Destructor of SamplingDistribution
53  virtual ~SamplingDistPlot();
54 
55  /// adds the sampling distribution and returns the scale factor
56  Double_t AddSamplingDistribution(const SamplingDistribution *samplingDist, Option_t *drawOptions="NORMALIZE HIST");
57  /// Like AddSamplingDistribution, but also sets a shaded area in the
58  /// minShaded and maxShaded boundaries.
59  Double_t AddSamplingDistributionShaded(const SamplingDistribution *samplingDist, Double_t minShaded, Double_t maxShaded, Option_t *drawOptions="NORMALIZE HIST");
60 
61  /// add a line
62  void AddLine(Double_t x1, Double_t y1, Double_t x2, Double_t y2, const char* title = NULL);
63  /// add a TH1
64  void AddTH1(TH1* h, Option_t *drawOptions="");
65  /// add a TF1
66  void AddTF1(TF1* f, const char* title = NULL, Option_t *drawOptions="SAME");
67  /// set legend
68  void SetLegend(TLegend* l){ fLegend = l; }
69 
70  void Draw(Option_t *options=0);
71 
72  /// Applies a predefined style if fApplyStyle is kTRUE (default).
73  void ApplyDefaultStyle(void);
74 
75  void SetLineColor(Color_t color, const SamplingDistribution *samplDist = 0);
76  void SetLineWidth(Width_t lwidth, const SamplingDistribution *samplDist = 0);
77  void SetLineStyle(Style_t style, const SamplingDistribution *samplDist = 0);
78 
79  void SetMarkerColor(Color_t color, const SamplingDistribution *samplDist = 0);
80  void SetMarkerStyle(Style_t style, const SamplingDistribution *samplDist = 0);
81  void SetMarkerSize(Size_t size, const SamplingDistribution *samplDist = 0);
82 
83  void RebinDistribution(Int_t rebinFactor, const SamplingDistribution *samplDist = 0);
84 
85  void SetAxisTitle(char *varName) { fVarName = TString(varName); }
86 
87  /// If you do not want SamplingDistPlot to interfere with your style settings, call this
88  /// function with "false" before Draw().
89  void SetApplyStyle(Bool_t s) { fApplyStyle = s; }
90 
91  /// Returns the TH1F associated with the give SamplingDistribution.
92  /// Intended use: Access to member functions of TH1F like GetMean(),
93  /// GetRMS() etc.
94  /// The return objects is managed by SamplingDistPlot
95  TH1F* GetTH1F(const SamplingDistribution *samplDist = NULL);
96  TH1 * GetHistogram(const SamplingDistribution *samplDist = NULL) { return GetTH1F(samplDist); }
97 
98  /// return plotter class used to draw the sampling distribution histograms
99  /// object is managed by SamplingDistPlot
100  RooPlot * GetPlot() { return fRooPlot; }
101 
102  /// changes plot to log scale on x axis
103  void SetLogXaxis(Bool_t lx) { fLogXaxis = lx; }
104  /// changes plot to log scale on y axis
105  void SetLogYaxis(Bool_t ly) { fLogYaxis = ly; }
106 
107  /// change x range
108  void SetXRange( double mi, double ma ) { fXMin = mi; fXMax = ma; }
109  /// change y range
110  void SetYRange( double mi, double ma ) { fYMin = mi; fYMax = ma; }
111 
112  /// write to Root file
113  void DumpToFile(const char* RootFileName, Option_t *option="", const char *ftitle="", Int_t compress=1);
114 
115  private:
116  std::vector<Double_t> fSamplingDistr;
117  std::vector<Double_t> fSampleWeights;
118 
120 
124 
126 
127  protected:
128 
131 
132  RooList fItems; /// holds TH1Fs only
133  RooList fOtherItems; /// other objects to be drawn like TLine etc.
134  TIterator* fIterator; /// TODO remove class variable and instantiate locally as necessary
138 
139  double fXMin, fXMax, fYMin, fYMax;
140 
143 
144  void SetSampleWeights(const SamplingDistribution *samplingDist);
145 
146  void addObject(TObject *obj, Option_t *drawOptions=0); // for TH1Fs only
147  void addOtherObject(TObject *obj, Option_t *drawOptions=0);
148  void GetAbsoluteInterval(Double_t &theMin, Double_t &theMax, Double_t &theYMax) const;
149 
150  ClassDef(SamplingDistPlot,1) /// Class containing the results of the HybridCalculator
151  };
152 }
153 
154 #endif
RooPlot * fRooPlot
TODO remove class variable and instantiate locally as necessary.
void SetLogXaxis(Bool_t lx)
changes plot to log scale on x axis
This class displays a legend box (TPaveText) containing several legend entries.
Definition: TLegend.h:27
void ApplyDefaultStyle(void)
Applies a predefined style if fApplyStyle is kTRUE (default).
short Style_t
Definition: RtypesCore.h:76
RooPlot * GetPlot()
return plotter class used to draw the sampling distribution histograms object is managed by SamplingD...
void SetLineColor(Color_t color, const SamplingDistribution *samplDist=0)
Sets line color for given sampling distribution and fill color for the associated shaded TH1F...
float Size_t
Definition: RtypesCore.h:83
const char Option_t
Definition: RtypesCore.h:62
TH1 * h
Definition: legend2.C:5
Double_t AddSamplingDistribution(const SamplingDistribution *samplingDist, Option_t *drawOptions="NORMALIZE HIST")
adds the sampling distribution and returns the scale factor
Basic string class.
Definition: TString.h:137
tomato 1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:575
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
int nbins[3]
void SetLogYaxis(Bool_t ly)
changes plot to log scale on y axis
void SetSampleWeights(const SamplingDistribution *samplingDist)
Determine if the sampling distribution has weights and store them.
void GetAbsoluteInterval(Double_t &theMin, Double_t &theMax, Double_t &theYMax) const
Iterator abstract base class.
Definition: TIterator.h:32
static const double x2[5]
#define ClassDef(name, id)
Definition: Rtypes.h:254
RooPlotable is a &#39;mix-in&#39; base class that define the standard RooFit plotting and printing methods...
Definition: RooPrintable.h:26
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:33
void addObject(TObject *obj, Option_t *drawOptions=0)
void SetLineWidth(Width_t lwidth, const SamplingDistribution *samplDist=0)
void SetXRange(double mi, double ma)
change x range
void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
void AddTH1(TH1 *h, Option_t *drawOptions="")
add a TH1
void SetMarkerSize(Size_t size, const SamplingDistribution *samplDist=0)
short Color_t
Definition: RtypesCore.h:79
Double_t AddSamplingDistributionShaded(const SamplingDistribution *samplingDist, Double_t minShaded, Double_t maxShaded, Option_t *drawOptions="NORMALIZE HIST")
Like AddSamplingDistribution, but also sets a shaded area in the minShaded and maxShaded boundaries...
RooList fOtherItems
holds TH1Fs only
void DumpToFile(const char *RootFileName, Option_t *option="", const char *ftitle="", Int_t compress=1)
write to Root file
void SetMarkerStyle(Style_t style, const SamplingDistribution *samplDist=0)
void SetLegend(TLegend *l)
set legend
TIterator * fIterator
other objects to be drawn like TLine etc.
void SetAxisTitle(char *varName)
TLine * l
Definition: textangle.C:4
A RooList is a TList with extra support for working with options that are associated with each node...
Definition: RooList.h:21
TH1 * GetHistogram(const SamplingDistribution *samplDist=NULL)
void SetYRange(double mi, double ma)
change y range
This class simply holds a sampling distribution of some test statistic.
void RebinDistribution(Int_t rebinFactor, const SamplingDistribution *samplDist=0)
short Width_t
Definition: RtypesCore.h:78
A RooPlot is a plot frame and a container for graphics objects within that frame. ...
Definition: RooPlot.h:41
Namespace for the RooStats classes.
Definition: Asimov.h:20
void AddTF1(TF1 *f, const char *title=NULL, Option_t *drawOptions="SAME")
add a TF1
static const double x1[5]
TH1F * GetTH1F(const SamplingDistribution *samplDist=NULL)
Returns the TH1F associated with the give SamplingDistribution.
double f(double x)
double Double_t
Definition: RtypesCore.h:55
TCanvas * style()
Definition: style.C:1
The TH1 histogram class.
Definition: TH1.h:80
void AddLine(Double_t x1, Double_t y1, Double_t x2, Double_t y2, const char *title=NULL)
add a line
void SetMarkerColor(Color_t color, const SamplingDistribution *samplDist=0)
std::vector< Double_t > fSampleWeights
This class provides simple and straightforward utilities to plot SamplingDistribution objects...
Mother of all ROOT objects.
Definition: TObject.h:37
void SetLineStyle(Style_t style, const SamplingDistribution *samplDist=0)
1-Dim function class
Definition: TF1.h:149
#define NULL
Definition: Rtypes.h:82
void addOtherObject(TObject *obj, Option_t *drawOptions=0)
void SetApplyStyle(Bool_t s)
If you do not want SamplingDistPlot to interfere with your style settings, call this function with "f...
virtual ~SamplingDistPlot()
Destructor of SamplingDistribution.
SamplingDistPlot(Int_t nbins=100)
Constructors for SamplingDistribution.
std::vector< Double_t > fSamplingDistr