ATLAS Offline Software
Loading...
Searching...
No Matches
Efficiency2D.h
Go to the documentation of this file.
1/* emacs: this is -*- c++ -*- */
10
11
12#ifndef TIDA_EFFICIENCY2D_H
13#define TIDA_EFFICIENCY2D_H
14
15#include "T_Efficiency.h"
16#include "TH2F.h"
17
18
19class Efficiency2D : public T_Efficiency<TH2F> {
20
21public:
22
23 Efficiency2D(TH2F* h, const std::string& n="") :
24 T_Efficiency<TH2F>( h, n ) { }
25
26 Efficiency2D(TH2F* hnum, TH2F* hden, const std::string& n, double scale=100) :
27 T_Efficiency<TH2F>( hnum, hden, n, scale ) {
28 getibinvec(true);
29 finalise( scale );
30 }
31
33
35
36 void Fill( double x, double y, double w=1) {
37 m_hnumer->Fill( float(x), float(y), float(w) );
38 m_hdenom->Fill( float(x), float(y), float(w) );
39 }
40
41 void FillDenom( double x, double y, float w=1) {
42 m_hmissed->Fill(float(x), float(y), float(w));
43 m_hdenom->Fill(float(x), float(y), float(w));
44 }
45
46
48
49 TGraphAsymmErrors* BayesX(int slice, double scale=100) {
50
51 if ( slice<0 || slice>=slicesX() ) return 0;
52
53 TH1D* hn = 0;
54 TH1D* hd = 0;
55
56 if ( m_hnumer ) {
57 // get some slice from the numerator histo
58 hn = m_hnumer->ProjectionY( slicename( name()+"_yslice", slice).c_str(), slice+1, slice+1, "e");
59 }
60
61 if ( m_hdenom ) {
62 // get some slice from the denominator histo
63 hd = m_hdenom->ProjectionY( slicename( name()+"_yslice", slice).c_str(), slice+1, slice+1, "e");
64 }
65
66 return BayesInternal( hn, hd, scale );
67
68 }
69
70
71
72 TGraphAsymmErrors* BayesY(int slice, double scale=100) {
73
74 if ( slice<0 || slice>=slicesY() ) return 0;
75
76 TH1D* hn = 0;
77 TH1D* hd = 0;
78
79 if ( m_hnumer ) {
80 // get some slice from the numerator histo
81 hn = m_hnumer->ProjectionX( slicename( name()+"_yslice", slice).c_str(), slice+1, slice+1, "e");
82 }
83
84 if ( m_hdenom ) {
85 // get some slice from the denominator histo
86 hd = m_hdenom->ProjectionX( slicename( name()+"_yslice", slice).c_str(), slice+1, slice+1, "e");
87 }
88
89 return BayesInternal( hn, hd, scale );
90
91 }
92
93
94 TH1D* sliceX( int i ) {
95 if ( i<0 || i>=slicesX() ) return 0;
96 return m_heff->ProjectionY( slicename( name()+"_xslice", i).c_str(), i+1, i+1, "e");
97 }
98
99 TH1D* sliceY( int i ) {
100 if ( i<0 || i>=slicesY() ) return 0;
101 return m_heff->ProjectionX( slicename( name()+"_yslice", i).c_str(), i+1, i+1, "e");
102 }
103
104
105 int slicesX() const {
106 if ( m_hdenom ) return m_hdenom->GetXaxis()->GetNbins();
107 return 0;
108 }
109
110
111 int slicesY() const {
112 if ( m_hdenom ) return m_hdenom->GetYaxis()->GetNbins();
113 return 0;
114 }
115
116protected:
117
118 virtual void getibinvec(bool force=false) {
119 if ( !force && !m_ibin.empty() ) return;
120 for ( int i=1 ; i<=m_hdenom->GetNbinsX() ; i++ ) {
121 for ( int j=1 ; j<=m_hdenom->GetNbinsY() ; j++ ) {
122 m_ibin.push_back( m_hdenom->GetBin(i,j) );
123 }
124 }
125 }
126
127
128 std::string slicename( const std::string& s, int i ) const {
129 if ( i<10 ) return s + "0" + std::to_string(i);
130 else return s + std::to_string(i);
131 }
132
133};
134
135
136
137
138
139inline std::ostream& operator<<( std::ostream& s, const Efficiency2D& ) {
140 return s;
141}
142
143
144#endif // TIDA_EFFICIENCY2D_H
145
146
147
148
149
150
151
152
153
154
std::ostream & operator<<(std::ostream &s, const Efficiency2D &)
#define y
#define x
Header file for AthHistogramAlgorithm.
int slicesX() const
TH1D * sliceX(int i)
void FillDenom(double x, double y, float w=1)
Efficiency2D(TH2F *h, const std::string &n="")
std::string slicename(const std::string &s, int i) const
int slicesY() const
virtual void getibinvec(bool force=false)
Efficiency2D(TH2F *hnum, TH2F *hden, const std::string &n, double scale=100)
TGraphAsymmErrors * BayesY(int slice, double scale=100)
TGraphAsymmErrors * BayesX(int slice, double scale=100)
evaluate the uncertainties correctly ...
TH1D * sliceY(int i)
void Fill(double x, double y, double w=1)
fill methods ...
const std::string & name() const
TGraphAsymmErrors * BayesInternal(TH1 *hn, TH1 *hd, double scale=100) const
std::vector< int > m_ibin
T_Efficiency(TH2F *h, const std::string &n)
void finalise(double scale=100)