ATLAS Offline Software
Loading...
Searching...
No Matches
T_Efficiency.h
Go to the documentation of this file.
1/* emacs: this is -*- c++ -*- */
10
11
12#ifndef TIDA_T_EFFICIENCY_H
13#define TIDA_T_EFFICIENCY_H
14
15#include <iostream>
16#include <string>
17#include <cmath>
18
19#include "TPad.h"
20#include "TH1.h"
21#include "TGraphAsymmErrors.h"
22
23template<typename T>
25
26public:
27
28 T_Efficiency(T* h, const std::string& n) :
29 m_name(n.empty() ? std::string(h->GetName())+"_eff" : n) {
30
31 m_hnumer = (T*)h->Clone( (m_name+"_n").c_str() );
32 m_hdenom = (T*)h->Clone( (m_name+"_d").c_str() );
33 m_heff = (T*)h->Clone( m_name.c_str() );
34
35 m_hmissed = (T*)h->Clone( (m_name+"_missed").c_str() );
36
37 m_hnumer->Reset();
38 m_hdenom->Reset();
39 m_heff->Reset();
40 m_hmissed->Reset();
41
42 }
43
44
45 T_Efficiency(T* hnum, T* hden, const std::string& n, double ) :
46 m_name(n+"_eff") {
47
48 m_hnumer = (T*)hnum->Clone( (m_name+"_n").c_str() );
49 m_hdenom = (T*)hden->Clone( (m_name+"_d").c_str() );
50 m_heff = (T*)hnum->Clone( m_name.c_str() );
51
52 m_hmissed = (T*)hnum->Clone( (m_name+"_missed").c_str() );
53
54 m_heff->Reset();
55 m_hmissed->Reset();
56
60 }
61
62 virtual ~T_Efficiency() { }
63
64 const std::string& name() const { return m_name; }
65
66 T* Hist() { return m_heff; }
67
68
70
71 void finalise(double scale=100) {
72 getibinvec();
73 m_heff->Reset();
74 for ( size_t i=0 ; i<m_ibin.size() ; i++ ) {
75 double n = m_hnumer->GetBinContent(m_ibin[i]);
76 double d = m_hdenom->GetBinContent(m_ibin[i]);
77
78 double e = 0;
79 double ee = 0;
80 if ( d!=0 ) {
81 e = n/d;
82 ee = e*(1-e)/d;
83 }
84
85 // need proper error calculation...
86 m_heff->SetBinContent( m_ibin[i], scale*e );
87 m_heff->SetBinError( m_ibin[i], scale*std::sqrt(ee) );
88
89 }
90 }
91
92
95
96 // 1D fill
97 // virtual void Fill( double x, double w=1) { }
98
99 // virtual void FillDenom( double x, float w=1) { }
100
101 // 2D fill
102 // virtual void Fill( double x, double y, double w=1) { }
103
104 // virtual void FillDenom( double x, double y, float w=1) { }
105
106#if 0
107
110
111 void SetDenominator( T* h ) {
112 getibinvec();
113 for ( size_t i=0 ; i<m_ibin.size() ; i++ ) {
114 m_hdenom->SetBinContent( m_ibin[i], h->GetBinContent(m_ibin[i]) );
115 m_hdenom->SetBinError( m_ibin[i], h->GetBinError(m_ibin[i]) );
116 }
117 }
118
119 virtual void SetNumerator( T* h ) {
120 getibinvec();
121 for ( size_t i=0 ; i<m_ibin.size() ; i++ ) {
122 m_hnumer->SetBinContent( m_ibin[i], h->GetBinContent(m_ibin[i]) );
123 m_hnumer->SetBinError( m_ibin[i], h->GetBinError(m_ibin[i]) );
124 }
125 }
126
127#endif
128
130
132
133 getibinvec();
134
135 double n_tot = 0;
136 double d_tot = 0;
137
138 for ( size_t i=0 ; i<m_ibin.size() ; i++ ) {
139 double n = m_hnumer->GetBinContent(m_ibin[i]);
140 double d = m_hdenom->GetBinContent(m_ibin[i]);
141
142 n_tot += n;
143 d_tot += d;
144 }
145
146 if ( d_tot!=0 ) {
147 return n_tot / d_tot;
148 }
149
150 return 0.;
151 }
152
153 void Write() {
154 m_hnumer->Write();
155 m_hdenom->Write();
156 m_hmissed->Write();
157 m_heff->Write();
158 }
159
160protected:
161
162 virtual void getibinvec(bool force=false) = 0;
163
164protected:
165
166 TGraphAsymmErrors* BayesInternal( TH1* hn, TH1* hd, double scale=100) const {
167
176
177 for ( int i=1 ; i<=hd->GetNbinsX() ; i++ ) { // cppcheck-suppress [ctunullpointer, nullPointer]; false positive
178 double y = hd->GetBinContent(i);
179 if ( y==0 ) hd->SetBinContent(i, 1e-20);
180 }
181
182 TGraphAsymmErrors* tg = new TGraphAsymmErrors( hn, hd, "cl=0.683 b(1,1) mode" );
183
184
185 double* x = tg->GetX();
186 double* y = tg->GetY();
187
188 int n = tg->GetN();
189
190 for ( int i=0 ; i<n ; i++ ) {
191
192 y[i] *= scale;
193
194 double yeup = tg->GetErrorYhigh(i);
195 double yedown = tg->GetErrorYlow(i);
196
197 yeup *= scale;
198 yedown *= scale;
199
206
207 tg->SetPoint( i, x[i], y[i] );
208
209 tg->SetPointEYhigh( i, yeup );
210 tg->SetPointEYlow( i, yedown );
211
212 tg->SetPointEXhigh( i, 0 );
213 tg->SetPointEXlow( i, 0 );
214
215 }
216
217 return tg;
218
219 }
220
221
222protected:
223
224 std::string m_name;
225
228
230
232
233 std::vector<int> m_ibin;
234
235};
236
237
238
239
240#endif // TIDA_T_EFFICIENCY_H
241
242
243
244
245
246
247
248
249
250
#define y
#define x
static const Attributes_t empty
Header file for AthHistogramAlgorithm.
const std::string & name() const
TGraphAsymmErrors * BayesInternal(TH1 *hn, TH1 *hd, double scale=100) const
virtual void getibinvec(bool force=false)=0
virtual ~T_Efficiency()
std::vector< int > m_ibin
T_Efficiency(T *h, const std::string &n)
void finalise(double scale=100)
actually calculate the efficiencies
std::string m_name
double findTotalEfficiency()
these 1D and 2D Fill versionms should never be called directly in fact, are they even needed at all ?
T_Efficiency(T *hnum, T *hden, const std::string &n, double)
STL class.
STL namespace.