ATLAS Offline Software
EventAnalysis.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2020 CERN for the benefit of the ATLAS collaboration
3 */
4 
5 
6 // Bit of an abstract class ( not in the C++ sense )
7 // designed to apply a bunch of cuts & chop up the data
8 
9 #ifndef IDPERFMON_EVENTANALYSIS_H
10 #define IDPERFMON_EVENTANALYSIS_H
11 
12 
13 // Crap that may stay or go
14 #include "CLHEP/Vector/LorentzVector.h"
15 
16 #include <map>
17 #include <vector>
18 
19 class TH1F;
20 class TH2F;
21 class TProfile;
22 class TProfile2D;
23 
24 namespace EAna
25 {
26  const float g_fMuonMass = 0.1056f;
27  const float g_fElecMass = 0.0005f;
28  const float CGeV = 1.0e-3;
29 }
30 
32 {
33  public:
34  EventAnalysis();
35  virtual ~EventAnalysis();
36 
37  // Overridden functions.
38  virtual void Init();
39  static constexpr float invalidAnswer{-999.9f};
40 
41  // Static Util. function declarations. Defined below class. Can use if no inheritance struct.
42  template<class T> static CLHEP::Hep3Vector calculateMomentum(const T * pP);
43  template<class T> static float EvalInvMass( const T* pxP1, const T* pxP2,
44  float fMass1, float fMass2 = invalidAnswer );
45  template<class T> static float EvalInvMass( const T* pxP1, const T* pxP2, const T* pxp3, const T* pxP4,
46  float fMass1, float fMass2 = -999.9, float fMass3 = -999.9, float fMass4 = invalidAnswer );
47  template<class T> static float EvalDiMuInvMass( const T* pxP1, const T* pxP2 );
48  template<class T> static float EvalFourMuInvMass( const T* pxP1, const T* pxP2, const T* pxP3, const T* pxP4);
49  template<class T> static float EvaluateAngle( const T* pxP1, const T* pxP2 );
50  template<class T> static float EvalPtDiff( const T* pxP1, const T* pxP2 );
51  template<class T> static float EvalPhiDiff( const T* pxP1, const T* pxP2 );
52  template<class T> static float EvalEtaDiff( const T* pxP1, const T* pxP2 );
53 
54  template<class T> static float EvalPt( const T* pxP1, const T* pxP2 );
55  template<class T> static float EvalPhi( const T* pxP1, const T* pxP2 );
56  template<class T> static float EvalEta( const T* pxP1, const T* pxP2 );
57  template<class T> static float EvalCharge( const T* pxP1, const T* pxP2 );
58 
59  template<class T> static float EvalTransverseMass( const T* pxP1,float fMETx, float fMETy,
60  float fMass1, float fMass2 = invalidAnswer);
61  template<class T> static float EvalTransverseMass( const T* pxP1, float fMETx, float fMETy );
62 
63  template<class T> static float EvalTransverseMass( const T* pxP1, const T* pxP2,
64  float fMETx, float fMETy,
65  float fMass1, float fMass2 = invalidAnswer);
66  template<class T> static float EvalTransverseMass( const T* pxP1, const T* pxP2,
67  float fMETx, float fMETy );
68 
69 
70 
71  protected:
72  virtual void BookHistograms();
73 
74  unsigned int m_uPassedEvents;
75  std::map<unsigned int, TH1F*> m_x1DHistograms;
76  std::map<unsigned int, TH2F*> m_x2DHistograms;
77  std::map<unsigned int, TProfile*> m_x1DProfHistograms;
78  std::map<unsigned int, TProfile2D*> m_x2DProfHistograms;
79 
80  std::string m_xSampleName;
81 
82  private:
83  void Register(); // Register the histograms.
84 };
85 
86 
87 //=============================================================================
88 // Useful static functions defined here
89 //=============================================================================
90 template<class T> CLHEP::Hep3Vector
92  const auto & p4(pP->p4());
93  return CLHEP::Hep3Vector(p4.Px() * EAna::CGeV , p4.Py() * EAna::CGeV, p4.Pz() * EAna::CGeV);
94 }
95 
96 // 2 Particle Invariant Mass
97 template<class T> float EventAnalysis::EvalDiMuInvMass( const T* pxP1, const T* pxP2 )
98 {
99  // Check integrity of inputs.
100  if ( !pxP1 || !pxP2 ) return invalidAnswer;
101 
102  // Evaluate Di-mu invariant mass.
103  return EvalInvMass( pxP1, pxP2, EAna::g_fMuonMass );
104 }
105 
106 template<class T> float EventAnalysis::EvalInvMass( const T* pxP1, const T* pxP2,
107  float fMass1, float fMass2 /* = -999.9 */ )
108 {
109  // Check integrity of inputs.No tachyons.
110  if ( fMass1 < 0.0f ) return invalidAnswer;
111  if ( !pxP1 || !pxP2 ) return invalidAnswer;
112  // Set masses equal if required by user
113  fMass2 = ( fMass2 < 0.0f ) ? fMass1 : fMass2;
114  // Evaluate invariant mass.
115  CLHEP::Hep3Vector xTmp1 = calculateMomentum(pxP1);
116  CLHEP::Hep3Vector xTmp2 = calculateMomentum(pxP2);
117  CLHEP::HepLorentzVector xLVec1; xLVec1.setVectM( xTmp1, fMass1 );
118  CLHEP::HepLorentzVector xLVec2; xLVec2.setVectM( xTmp2, fMass2 );
119  return static_cast<float>( xLVec1.invariantMass( xLVec2 ) );
120 }
121 
122 // 4 Particle Invariant Mass
123 template<class T> float EventAnalysis::EvalFourMuInvMass( const T* pxP1, const T* pxP2, const T* pxP3, const T* pxP4 )
124 {
125  // Check integrity of inputs.
126  if ( !pxP1 || !pxP2 || !pxP3 || !pxP4) return invalidAnswer;
127 
128  // Evaluate invariant mass.
129  return EvalInvMass( pxP1, pxP2, pxP3, pxP4, EAna::g_fMuonMass );
130 }
131 
132 template<class T> float EventAnalysis::EvalInvMass( const T* pxP1, const T* pxP2, const T* pxP3, const T* pxP4,
133  float fMass1, float fMass2, float fMass3, float fMass4 /* = -999.9 */ )
134 {
135  // Check integrity of inputs.No tachyons.
136  if ( fMass1 < 0.0f ) return invalidAnswer;
137  if ( !pxP1 || !pxP2 || !pxP3 || !pxP4) return invalidAnswer;
138 
139  // Set masses equal if required by user
140  fMass2 = ( fMass2 < 0.0f ) ? fMass1 : fMass2;
141  fMass3 = ( fMass3 < 0.0f ) ? fMass1 : fMass3;
142  fMass4 = ( fMass4 < 0.0f ) ? fMass1 : fMass4;
143 
144  // Evaluate invariant mass.
145  CLHEP::Hep3Vector xTmp1 = CLHEP::Hep3Vector( pxP1->p4().Px() * EAna::CGeV, pxP1->p4().Py() * EAna::CGeV, pxP1->p4().Pz() * EAna::CGeV );
146  CLHEP::Hep3Vector xTmp2 = CLHEP::Hep3Vector( pxP2->p4().Px() * EAna::CGeV, pxP2->p4().Py() * EAna::CGeV, pxP2->p4().Pz() * EAna::CGeV );
147  CLHEP::Hep3Vector xTmp3 = CLHEP::Hep3Vector( pxP3->p4().Px() * EAna::CGeV, pxP3->p4().Py() * EAna::CGeV, pxP3->p4().Pz() * EAna::CGeV );
148  CLHEP::Hep3Vector xTmp4 = CLHEP::Hep3Vector( pxP4->p4().Px() * EAna::CGeV, pxP4->p4().Py() * EAna::CGeV, pxP4->p4().Pz() * EAna::CGeV );
149 
150  CLHEP::HepLorentzVector xLVec; xLVec.setVectM ( xTmp1, fMass1 );
151  CLHEP::HepLorentzVector xLVec2; xLVec2.setVectM( xTmp2, fMass2 );
152  CLHEP::HepLorentzVector xLVec3; xLVec3.setVectM( xTmp3, fMass3 );
153  CLHEP::HepLorentzVector xLVec4; xLVec4.setVectM( xTmp4, fMass4 );
154 
155  xLVec += xLVec2;
156  xLVec += xLVec3;
157  xLVec += xLVec4;
158 
159  return static_cast<float>( xLVec.m() );
160 }
161 
162 // Angle Between two particles
163 template<class T> float EventAnalysis::EvaluateAngle( const T* pxP1, const T* pxP2 )
164 {
165  // Check integrity of inputs.
166  if ( !pxP1 || !pxP2 ) return invalidAnswer;
167  // Evaluate the angle.
168  CLHEP::Hep3Vector xTmp1 = calculateMomentum(pxP1);
169  CLHEP::Hep3Vector xTmp2 = calculateMomentum(pxP2);
170  return static_cast<float>( xTmp1.angle(xTmp2) );
171 }
172 
173 template<class T> float EventAnalysis::EvalPtDiff( const T* pxP1, const T* pxP2 )
174 {
175  // Check integrity of inputs.
176  if ( !pxP1 || !pxP2 ) return invalidAnswer;
177  // Evaluate the difference between the momenta. Signed using positive - negative if appropriate.
178  if ( pxP1->charge() > 0.5f )
179  {
180  return static_cast<float>( pxP1->pt() * EAna::CGeV - pxP2->pt() * EAna::CGeV );
181  }
182  else
183  {
184  return static_cast<float>( pxP2->pt() * EAna::CGeV - pxP1->pt() * EAna::CGeV );
185  }
186 }
187 
188 
189 
190 template<class T> float EventAnalysis::EvalPhiDiff( const T* pxP1, const T* pxP2 )
191 {
192  // Check integrity of inputs.
193  if ( !pxP1 || !pxP2 ) return invalidAnswer;
194  // Evaluate the angle.
195  CLHEP::Hep3Vector xTmp1 = calculateMomentum(pxP1);
196  CLHEP::Hep3Vector xTmp2 = calculateMomentum(pxP2);
197  return static_cast<float>( xTmp1.deltaPhi(xTmp2) );
198 }
199 
200 
201 template<class T> float EventAnalysis::EvalEtaDiff( const T* pxP1, const T* pxP2 )
202 {
203  // Check integrity of inputs.
204  if ( !pxP1 || !pxP2 ) return invalidAnswer;
205  // Evaluate the angle.
206  CLHEP::Hep3Vector xTmp1 = calculateMomentum(pxP1);
207  CLHEP::Hep3Vector xTmp2 = calculateMomentum(pxP2);
208  return static_cast<float>( xTmp1.polarAngle(xTmp2) );
209 }
210 
211 template<class T> float EventAnalysis::EvalPt( const T* pxP1, const T* pxP2 )
212 {
213  // Check integrity of inputs.
214  if ( !pxP1 || !pxP2 ) return invalidAnswer;
215  CLHEP::Hep3Vector xTmp1 = calculateMomentum(pxP1);
216  CLHEP::Hep3Vector xTmp2 = calculateMomentum(pxP2);
217  return static_cast<float>( (xTmp1 + xTmp2).perp() );
218 }
219 
220 template<class T> float EventAnalysis::EvalPhi( const T* pxP1, const T* pxP2 )
221 {
222  // Check integrity of inputs.
223  if ( !pxP1 || !pxP2 ) return invalidAnswer;
224  CLHEP::Hep3Vector xTmp1 = calculateMomentum(pxP1);
225  CLHEP::Hep3Vector xTmp2 = calculateMomentum(pxP2);
226  return static_cast<float>( (xTmp1 + xTmp2).phi() );
227 }
228 
229 template<class T> float EventAnalysis::EvalEta( const T* pxP1, const T* pxP2 )
230 {
231  // Check integrity of inputs.
232  if ( !pxP1 || !pxP2 ) return invalidAnswer;
233  CLHEP::Hep3Vector xTmp1 = calculateMomentum(pxP1);
234  CLHEP::Hep3Vector xTmp2 = calculateMomentum(pxP2);
235  return static_cast<float>( (xTmp1 + xTmp2).pseudoRapidity() );
236 }
237 
238 template<class T> float EventAnalysis::EvalCharge( const T* pxP1, const T* pxP2 )
239 {
240  // Check integrity of inputs.
241  if ( !pxP1 || !pxP2 ) return invalidAnswer;
242  return static_cast<float>( pxP1->charge() + pxP2->charge() );
243 }
244 
245 // Transverse Mass - 1 Particle + MET
246 template<class T> float EventAnalysis::EvalTransverseMass( const T* pxP1, float fMETx, float fMETy )
247 {
248  // Check integrity of inputs.
249  if ( !pxP1 ) return invalidAnswer;
250  // Evaluate Di-mu invariant mass.
251  return EvalInvMass( pxP1, fMETx, fMETy, EAna::g_fMuonMass );
252 }
253 
254 template<class T> float EventAnalysis::EvalTransverseMass( const T* pxP1, float fMETx, float fMETy,
255  float fMass1, float fMass2 /* = -999.9 */ )
256 {
257  // Check integrity of inputs.No tachyons.
258  if ( fMass1 < 0.0f ) return invalidAnswer;
259  if ( !pxP1 ) return invalidAnswer;
260  // Set masses equal if required by user.
261  fMass2 = ( fMass2 < 0.0f ) ? fMass1 : fMass2;
262  // Evaluate invariant mass.
263  CLHEP::Hep3Vector xTmp1 = CLHEP::Hep3Vector( pxP1->p4().Px() * EAna::CGeV, pxP1->p4().Py() * EAna::CGeV, 0.0f );
264  CLHEP::Hep3Vector xTmp2 = CLHEP::Hep3Vector( fMETx, fMETy, 0.0f );
265  CLHEP::HepLorentzVector xLVec1; xLVec1.setVectM( xTmp1, fMass1 );
266  CLHEP::HepLorentzVector xLVec2; xLVec2.setVectM( xTmp2, 0.0f );
267  return static_cast<float>( xLVec1.invariantMass( xLVec2 ) );
268 }
269 
270 // Transverse Mass - 2 Particles + MET
271 template<class T> float EventAnalysis::EvalTransverseMass( const T* pxP1, const T* pxP2, float fMETx, float fMETy )
272 {
273  // Check integrity of inputs.
274  if ( !pxP1 || !pxP2 ) return invalidAnswer;
275  // Evaluate Di-mu invariant mass.
276  return EvalTransverseMass( pxP1, pxP2, fMETx, fMETy, EAna::g_fMuonMass );
277 }
278 
279 template<class T> float EventAnalysis::EvalTransverseMass( const T* pxP1, const T* pxP2, float fMETx, float fMETy,
280  float fMass1, float fMass2 /* = -999.9 */ )
281 {
282  // Check integrity of inputs.No tachyons.
283  if ( fMass1 < 0.0f ) return invalidAnswer;
284  if ( !pxP1 || !pxP2 ) return invalidAnswer;
285  // Set masses equal if required by user.
286  fMass2 = ( fMass2 < 0.0f ) ? fMass1 : fMass2;
287  // Evaluate invariant mass.
288  CLHEP::Hep3Vector xTmp1 = CLHEP::Hep3Vector( pxP1->p4().Px() * EAna::CGeV, pxP1->p4().Py() * EAna::CGeV, 0.0f );
289  CLHEP::Hep3Vector xTmp2 = CLHEP::Hep3Vector( pxP2->p4().Px() * EAna::CGeV, pxP2->p4().Py() * EAna::CGeV, 0.0f );
290  CLHEP::Hep3Vector xTmp12 = xTmp1 + xTmp2;
291  CLHEP::Hep3Vector xTmp3 = CLHEP::Hep3Vector( fMETx, fMETy, 0.0f );
292  CLHEP::HepLorentzVector xLVec1; xLVec1.setVectM( xTmp12, fMass1 );
293  CLHEP::HepLorentzVector xLVec2; xLVec2.setVectM( xTmp3, 0.0f );
294  return static_cast<float>( xLVec1.invariantMass( xLVec2 ) );
295 }
296 
297 //=============================================================================
298 
299 #endif
EventAnalysis::EvalInvMass
static float EvalInvMass(const T *pxP1, const T *pxP2, float fMass1, float fMass2=invalidAnswer)
Definition: EventAnalysis.h:106
EAna::g_fElecMass
const float g_fElecMass
Definition: EventAnalysis.h:27
EventAnalysis::EvalEtaDiff
static float EvalEtaDiff(const T *pxP1, const T *pxP2)
Definition: EventAnalysis.h:201
EventAnalysis::EvaluateAngle
static float EvaluateAngle(const T *pxP1, const T *pxP2)
Definition: EventAnalysis.h:163
phi
Scalar phi() const
phi method
Definition: AmgMatrixBasePlugin.h:67
perp
Scalar perp() const
perp method - perpenticular length
Definition: AmgMatrixBasePlugin.h:44
EventAnalysis::~EventAnalysis
virtual ~EventAnalysis()
Definition: EventAnalysis.cxx:37
EventAnalysis::EvalDiMuInvMass
static float EvalDiMuInvMass(const T *pxP1, const T *pxP2)
Definition: EventAnalysis.h:97
EventAnalysis::EvalFourMuInvMass
static float EvalFourMuInvMass(const T *pxP1, const T *pxP2, const T *pxP3, const T *pxP4)
Definition: EventAnalysis.h:123
EventAnalysis
Definition: EventAnalysis.h:32
EventAnalysis::EvalCharge
static float EvalCharge(const T *pxP1, const T *pxP2)
Definition: EventAnalysis.h:238
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
EventAnalysis::m_x1DHistograms
std::map< unsigned int, TH1F * > m_x1DHistograms
Definition: EventAnalysis.h:75
EventAnalysis::EventAnalysis
EventAnalysis()
Definition: EventAnalysis.cxx:32
EventAnalysis::m_uPassedEvents
unsigned int m_uPassedEvents
Definition: EventAnalysis.h:74
EventAnalysis::m_x2DHistograms
std::map< unsigned int, TH2F * > m_x2DHistograms
Definition: EventAnalysis.h:76
EventAnalysis::EvalTransverseMass
static float EvalTransverseMass(const T *pxP1, float fMETx, float fMETy, float fMass1, float fMass2=invalidAnswer)
Definition: EventAnalysis.h:254
EventAnalysis::calculateMomentum
static CLHEP::Hep3Vector calculateMomentum(const T *pP)
Definition: EventAnalysis.h:91
EventAnalysis::EvalPtDiff
static float EvalPtDiff(const T *pxP1, const T *pxP2)
Definition: EventAnalysis.h:173
python.TrigEgammaMonitorHelper.TProfile
def TProfile(*args, **kwargs)
Definition: TrigEgammaMonitorHelper.py:81
EventAnalysis::BookHistograms
virtual void BookHistograms()
Definition: EventAnalysis.cxx:55
hist_file_dump.f
f
Definition: hist_file_dump.py:135
EventAnalysis::m_x1DProfHistograms
std::map< unsigned int, TProfile * > m_x1DProfHistograms
Definition: EventAnalysis.h:77
EventAnalysis::EvalEta
static float EvalEta(const T *pxP1, const T *pxP2)
Definition: EventAnalysis.h:229
EventAnalysis::Register
void Register()
Definition: EventAnalysis.cxx:63
EventAnalysis::EvalPt
static float EvalPt(const T *pxP1, const T *pxP2)
Definition: EventAnalysis.h:211
EventAnalysis::m_x2DProfHistograms
std::map< unsigned int, TProfile2D * > m_x2DProfHistograms
Definition: EventAnalysis.h:78
EAna::g_fMuonMass
const float g_fMuonMass
Definition: EventAnalysis.h:26
EventAnalysis::EvalPhiDiff
static float EvalPhiDiff(const T *pxP1, const T *pxP2)
Definition: EventAnalysis.h:190
EventAnalysis::EvalPhi
static float EvalPhi(const T *pxP1, const T *pxP2)
Definition: EventAnalysis.h:220
EAna
Definition: EventAnalysis.h:25
python.TrigEgammaMonitorHelper.TH1F
def TH1F(name, title, nxbins, bins_par2, bins_par3=None, path='', **kwargs)
Definition: TrigEgammaMonitorHelper.py:24
EventAnalysis::Init
virtual void Init()
Definition: EventAnalysis.cxx:44
EventAnalysis::m_xSampleName
std::string m_xSampleName
Definition: EventAnalysis.h:80
EAna::CGeV
const float CGeV
Definition: EventAnalysis.h:28
EventAnalysis::invalidAnswer
static constexpr float invalidAnswer
Definition: EventAnalysis.h:39
TSU::T
unsigned long long T
Definition: L1TopoDataTypes.h:35