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