ATLAS Offline Software
Efficiency2D.h
Go to the documentation of this file.
1 /* emacs: this is -*- c++ -*- */
12 #ifndef TIDA_EFFICIENCY2D_H
13 #define TIDA_EFFICIENCY2D_H
14 
15 #include "T_Efficiency.h"
16 #include "TH2F.h"
17 
18 
19 class Efficiency2D : public T_Efficiency<TH2F> {
20 
21 public:
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 
116 protected:
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 
139 inline 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 
python.SystemOfUnits.s
int s
Definition: SystemOfUnits.py:131
calibdata.force
bool force
Definition: calibdata.py:19
Efficiency2D::sliceY
TH1D * sliceY(int i)
Definition: Efficiency2D.h:99
T_Efficiency< TH2F >::finalise
void finalise(double scale=100)
actually calculate the efficiencies
Definition: T_Efficiency.h:80
yodamerge_tmp.scale
scale
Definition: yodamerge_tmp.py:138
Efficiency2D::slicesY
int slicesY() const
Definition: Efficiency2D.h:111
python.TrigEgammaMonitorHelper.TH2F
def TH2F(name, title, nxbins, bins_par2, bins_par3, bins_par4, bins_par5=None, bins_par6=None, path='', **kwargs)
Definition: TrigEgammaMonitorHelper.py:45
x
#define x
T_Efficiency
Definition: T_Efficiency.h:24
operator<<
std::ostream & operator<<(std::ostream &s, const Efficiency2D &)
Definition: Efficiency2D.h:139
T_Efficiency< TH2F >::m_hmissed
TH2F * m_hmissed
Definition: T_Efficiency.h:238
Efficiency2D::Fill
void Fill(double x, double y, double w=1)
fill methods ...
Definition: Efficiency2D.h:36
Efficiency2D::sliceX
TH1D * sliceX(int i)
Definition: Efficiency2D.h:94
Efficiency2D::~Efficiency2D
~Efficiency2D()
Definition: Efficiency2D.h:32
perfmonmt-refit.slice
slice
Definition: perfmonmt-refit.py:52
lumiFormat.i
int i
Definition: lumiFormat.py:85
beamspotman.n
n
Definition: beamspotman.py:731
Efficiency2D::FillDenom
void FillDenom(double x, double y, float w=1)
Definition: Efficiency2D.h:41
Efficiency2D::Efficiency2D
Efficiency2D(TH2F *hnum, TH2F *hden, const std::string &n, double scale=100)
Definition: Efficiency2D.h:26
ActsTrk::to_string
std::string to_string(const DetectorType &type)
Definition: GeometryDefs.h:34
T_Efficiency.h
T_Efficiency< TH2F >::m_heff
TH2F * m_heff
Definition: T_Efficiency.h:240
Efficiency2D::BayesY
TGraphAsymmErrors * BayesY(int slice, double scale=100)
Definition: Efficiency2D.h:72
y
#define y
h
Efficiency2D::slicesX
int slicesX() const
Definition: Efficiency2D.h:105
T_Efficiency< TH2F >::m_ibin
std::vector< int > m_ibin
Definition: T_Efficiency.h:242
Efficiency2D::getibinvec
virtual void getibinvec(bool force=false)
Definition: Efficiency2D.h:118
T_Efficiency< TH2F >::name
std::string name() const
Definition: T_Efficiency.h:73
Efficiency2D
Definition: Efficiency2D.h:19
python.IoTestsLib.w
def w
Definition: IoTestsLib.py:200
Efficiency2D::slicename
std::string slicename(const std::string &s, int i) const
Definition: Efficiency2D.h:128
Efficiency2D::Efficiency2D
Efficiency2D(TH2F *h, const std::string &n="")
Definition: Efficiency2D.h:23
Efficiency2D::BayesX
TGraphAsymmErrors * BayesX(int slice, double scale=100)
evaluate the uncertainties correctly ...
Definition: Efficiency2D.h:49
T_Efficiency< TH2F >::m_hdenom
TH2F * m_hdenom
Definition: T_Efficiency.h:236
T_Efficiency< TH2F >::BayesInternal
TGraphAsymmErrors * BayesInternal(TH1 *hn, TH1 *hd, double scale=100) const
Definition: T_Efficiency.h:175
T_Efficiency< TH2F >::m_hnumer
TH2F * m_hnumer
Definition: T_Efficiency.h:235