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 "RooPrintable.h"
17#include "TNamed.h"
18#include "TH1F.h"
19#include "TF1.h"
20#include "TLegend.h"
21#include "TList.h"
22
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 min, double max);
38
39 /// Destructor of SamplingDistribution
40 ~SamplingDistPlot() override;
41
42 /// adds the sampling distribution and returns the scale factor
43 double AddSamplingDistribution(const SamplingDistribution *samplingDist, Option_t *drawOptions="NORMALIZE HIST");
44 /// Like AddSamplingDistribution, but also sets a shaded area in the
45 /// minShaded and maxShaded boundaries.
46 double AddSamplingDistributionShaded(const SamplingDistribution *samplingDist, double minShaded, double maxShaded, Option_t *drawOptions="NORMALIZE HIST");
47
48 /// add a line
49 void AddLine(double x1, double y1, double x2, double y2, const char* title = nullptr);
50 /// add a TH1
51 void AddTH1(TH1* h, Option_t *drawOptions="");
52 /// add a TF1
53 void AddTF1(TF1* f, const char* title = nullptr, Option_t *drawOptions="SAME");
54 /// set legend
56
57 void Draw(Option_t *options=nullptr) override;
58
59 /// Applies a predefined style if fApplyStyle is true (default).
60 void ApplyDefaultStyle(void);
61
62 void SetLineColor(Color_t color, const SamplingDistribution *samplDist = nullptr);
63 void SetLineWidth(Width_t lwidth, const SamplingDistribution *samplDist = nullptr);
64 void SetLineStyle(Style_t style, const SamplingDistribution *samplDist = nullptr);
65
66 void SetMarkerColor(Color_t color, const SamplingDistribution *samplDist = nullptr);
67 void SetMarkerStyle(Style_t style, const SamplingDistribution *samplDist = nullptr);
68 void SetMarkerSize(Size_t size, const SamplingDistribution *samplDist = nullptr);
69
70 void RebinDistribution(Int_t rebinFactor, const SamplingDistribution *samplDist = nullptr);
71
72 void SetAxisTitle(char *varName) { fVarName = TString(varName); }
73
74 /// If you do not want SamplingDistPlot to interfere with your style settings, call this
75 /// function with "false" before Draw().
76 void SetApplyStyle(bool s) { fApplyStyle = s; }
77
78 /// Returns the TH1F associated with the give SamplingDistribution.
79 /// Intended use: Access to member functions of TH1F like GetMean(),
80 /// GetRMS() etc.
81 /// The return objects is managed by SamplingDistPlot
82 TH1F* GetTH1F(const SamplingDistribution *samplDist = nullptr);
83 TH1 * GetHistogram(const SamplingDistribution *samplDist = nullptr) { return GetTH1F(samplDist); }
84
85 /// return plotter class used to draw the sampling distribution histograms
86 /// object is managed by SamplingDistPlot
87 RooPlot * GetPlot() { return fRooPlot; }
88
89 /// changes plot to log scale on x axis
90 void SetLogXaxis(bool lx) { fLogXaxis = lx; }
91 /// changes plot to log scale on y axis
92 void SetLogYaxis(bool ly) { fLogYaxis = ly; }
93
94 /// change x range
95 void SetXRange( double mi, double ma ) { fXMin = mi; fXMax = ma; }
96 /// change y range
97 void SetYRange( double mi, double ma ) { fYMin = mi; fYMax = ma; }
98
99 /// write to Root file
100 void DumpToFile(const char* RootFileName, Option_t *option="", const char *ftitle="", Int_t compress = ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault);
101
102 private:
103 std::vector<double> fSamplingDistr;
104 std::vector<double> fSampleWeights;
105
107
111
113
114 protected:
115
118
119 TList fItems; ///< holds TH1Fs only
120 TList fOtherItems; ///< other objects to be drawn like TLine etc.
124
126
129
130 void SetSampleWeights(const SamplingDistribution *samplingDist);
131
132 void addObject(TObject *obj, Option_t *drawOptions=nullptr); // for TH1Fs only
133 void addOtherObject(TObject *obj, Option_t *drawOptions=nullptr);
134 void GetAbsoluteInterval(double &theMin, double &theMax, double &theYMax) const;
135
136 ClassDefOverride(SamplingDistPlot,2) /// Class containing the results of the HybridCalculator
137 };
138}
139
140#endif
#define f(i)
Definition RSha256.hxx:104
#define h(i)
Definition RSha256.hxx:106
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
short Style_t
Definition RtypesCore.h:89
short Color_t
Definition RtypesCore.h:92
float Size_t
Definition RtypesCore.h:96
short Width_t
Definition RtypesCore.h:91
const char Option_t
Definition RtypesCore.h:66
#define ClassDefOverride(name, id)
Definition Rtypes.h:341
Option_t Option_t option
Option_t Option_t SetLineWidth
Option_t Option_t SetLineColor
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t SetMarkerStyle
Option_t Option_t style
Option_t Option_t TPoint TPoint const char y1
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:43
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.
TList fItems
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 RebinDistribution(Int_t rebinFactor, const SamplingDistribution *samplDist=nullptr)
~SamplingDistPlot() override
Destructor of SamplingDistribution.
void SetLegend(TLegend *l)
set legend
void SetLogXaxis(bool lx)
changes plot to log scale on x axis
void GetAbsoluteInterval(double &theMin, double &theMax, double &theYMax) const
void SetLineStyle(Style_t style, const SamplingDistribution *samplDist=nullptr)
double AddSamplingDistributionShaded(const SamplingDistribution *samplingDist, double minShaded, double maxShaded, Option_t *drawOptions="NORMALIZE HIST")
Like AddSamplingDistribution, but also sets a shaded area in the minShaded and maxShaded boundaries.
void SetMarkerColor(Color_t color, const SamplingDistribution *samplDist=nullptr)
void SetYRange(double mi, double ma)
change y range
std::vector< double > fSampleWeights
void AddLine(double x1, double y1, double x2, double y2, const char *title=nullptr)
add a line
void ApplyDefaultStyle(void)
Applies a predefined style if fApplyStyle is true (default).
void AddTH1(TH1 *h, Option_t *drawOptions="")
add a TH1
RooPlot * GetPlot()
return plotter class used to draw the sampling distribution histograms object is managed by SamplingD...
void SetLogYaxis(bool ly)
changes plot to log scale on y axis
std::vector< double > fSamplingDistr
TH1 * GetHistogram(const SamplingDistribution *samplDist=nullptr)
void AddTF1(TF1 *f, const char *title=nullptr, Option_t *drawOptions="SAME")
add a TF1
double AddSamplingDistribution(const SamplingDistribution *samplingDist, Option_t *drawOptions="NORMALIZE HIST")
adds the sampling distribution and returns the scale factor
void addOtherObject(TObject *obj, Option_t *drawOptions=nullptr)
Add a generic object to this plot.
void addObject(TObject *obj, Option_t *drawOptions=nullptr)
Add a generic object to this plot.
void SetMarkerSize(Size_t size, const SamplingDistribution *samplDist=nullptr)
TH1F * GetTH1F(const SamplingDistribution *samplDist=nullptr)
Returns the TH1F associated with the give SamplingDistribution.
void SetAxisTitle(char *varName)
void SetXRange(double mi, double ma)
change x range
void SetApplyStyle(bool s)
If you do not want SamplingDistPlot to interfere with your style settings, call this function with "f...
TList fOtherItems
other objects to be drawn like TLine etc.
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:577
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:58
This class displays a legend box (TPaveText) containing several legend entries.
Definition TLegend.h:23
A doubly linked list.
Definition TList.h:38
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:139
Namespace for the RooStats classes.
Definition Asimov.h:19
@ kUseCompiledDefault
Use the compile-time default setting.
Definition Compression.h:50
th1 Draw()
TLine l
Definition textangle.C:4