ATLAS Offline Software
Loading...
Searching...
No Matches
EventAnalysis Class Reference

#include <EventAnalysis.h>

Inheritance diagram for EventAnalysis:
Collaboration diagram for EventAnalysis:

Public Member Functions

 EventAnalysis ()
virtual ~EventAnalysis ()
virtual void Init ()

Static Public Member Functions

template<class T>
static CLHEP::Hep3Vector calculateMomentum (const T *pP)
template<class T>
static float EvalInvMass (const T *pxP1, const T *pxP2, float fMass1, float fMass2=invalidAnswer)
template<class T>
static float EvalInvMass (const T *pxP1, const T *pxP2, const T *pxp3, const T *pxP4, float fMass1, float fMass2=-999.9, float fMass3=-999.9, float fMass4=invalidAnswer)
template<class T>
static float EvalDiMuInvMass (const T *pxP1, const T *pxP2)
template<class T>
static float EvalFourMuInvMass (const T *pxP1, const T *pxP2, const T *pxP3, const T *pxP4)
template<class T>
static float EvaluateAngle (const T *pxP1, const T *pxP2)
template<class T>
static float EvalPtDiff (const T *pxP1, const T *pxP2)
template<class T>
static float EvalPhiDiff (const T *pxP1, const T *pxP2)
template<class T>
static float EvalEtaDiff (const T *pxP1, const T *pxP2)
template<class T>
static float EvalPt (const T *pxP1, const T *pxP2)
template<class T>
static float EvalPhi (const T *pxP1, const T *pxP2)
template<class T>
static float EvalEta (const T *pxP1, const T *pxP2)
template<class T>
static float EvalCharge (const T *pxP1, const T *pxP2)
template<class T>
static float EvalTransverseMass (const T *pxP1, float fMETx, float fMETy, float fMass1, float fMass2=invalidAnswer)
template<class T>
static float EvalTransverseMass (const T *pxP1, float fMETx, float fMETy)
template<class T>
static float EvalTransverseMass (const T *pxP1, const T *pxP2, float fMETx, float fMETy, float fMass1, float fMass2=invalidAnswer)
template<class T>
static float EvalTransverseMass (const T *pxP1, const T *pxP2, float fMETx, float fMETy)

Static Public Attributes

static constexpr float invalidAnswer {-999.9f}

Protected Member Functions

virtual void BookHistograms ()

Protected Attributes

unsigned int m_uPassedEvents
std::map< unsigned int, TH1F * > m_x1DHistograms
std::map< unsigned int, TH2F * > m_x2DHistograms
std::map< unsigned int, TProfile * > m_x1DProfHistograms
std::map< unsigned int, TProfile2D * > m_x2DProfHistograms
std::string m_xSampleName

Private Member Functions

void Register ()

Detailed Description

Definition at line 32 of file EventAnalysis.h.

Constructor & Destructor Documentation

◆ EventAnalysis()

EventAnalysis::EventAnalysis ( )

Definition at line 32 of file EventAnalysis.cxx.

33{
35}
unsigned int m_uPassedEvents

◆ ~EventAnalysis()

EventAnalysis::~EventAnalysis ( )
virtual

Definition at line 37 of file EventAnalysis.cxx.

38{
39}

Member Function Documentation

◆ BookHistograms()

void EventAnalysis::BookHistograms ( )
protectedvirtual

Reimplemented in FourMuonEvent, MuonSelector, and ZmumuEvent.

Definition at line 55 of file EventAnalysis.cxx.

56{
57 // This must be overriden by an inheriting class.
58}

◆ calculateMomentum()

template<class T>
CLHEP::Hep3Vector EventAnalysis::calculateMomentum ( const T * pP)
static

Definition at line 92 of file EventAnalysis.h.

92 {
93 const auto & p4(pP->p4());
94 return CLHEP::Hep3Vector(p4.Px() * EAna::CGeV , p4.Py() * EAna::CGeV, p4.Pz() * EAna::CGeV);
95}
const float CGeV

◆ EvalCharge()

template<class T>
float EventAnalysis::EvalCharge ( const T * pxP1,
const T * pxP2 )
static

Definition at line 239 of file EventAnalysis.h.

240{
241 // Check integrity of inputs.
242 if ( !pxP1 || !pxP2 ) return invalidAnswer;
243 return static_cast<float>( pxP1->charge() + pxP2->charge() );
244}
static constexpr float invalidAnswer

◆ EvalDiMuInvMass()

template<class T>
float EventAnalysis::EvalDiMuInvMass ( const T * pxP1,
const T * pxP2 )
static

Definition at line 98 of file EventAnalysis.h.

99{
100 // Check integrity of inputs.
101 if ( !pxP1 || !pxP2 ) return invalidAnswer;
102
103 // Evaluate Di-mu invariant mass.
104 return EvalInvMass( pxP1, pxP2, EAna::g_fMuonMass );
105}
static float EvalInvMass(const T *pxP1, const T *pxP2, float fMass1, float fMass2=invalidAnswer)
const float g_fMuonMass

◆ EvalEta()

template<class T>
float EventAnalysis::EvalEta ( const T * pxP1,
const T * pxP2 )
static

Definition at line 230 of file EventAnalysis.h.

231{
232 // Check integrity of inputs.
233 if ( !pxP1 || !pxP2 ) return invalidAnswer;
234 CLHEP::Hep3Vector xTmp1 = calculateMomentum(pxP1);
235 CLHEP::Hep3Vector xTmp2 = calculateMomentum(pxP2);
236 return static_cast<float>( (xTmp1 + xTmp2).pseudoRapidity() );
237}
static CLHEP::Hep3Vector calculateMomentum(const T *pP)

◆ EvalEtaDiff()

template<class T>
float EventAnalysis::EvalEtaDiff ( const T * pxP1,
const T * pxP2 )
static

Definition at line 202 of file EventAnalysis.h.

203{
204 // Check integrity of inputs.
205 if ( !pxP1 || !pxP2 ) return invalidAnswer;
206 // Evaluate the angle.
207 CLHEP::Hep3Vector xTmp1 = calculateMomentum(pxP1);
208 CLHEP::Hep3Vector xTmp2 = calculateMomentum(pxP2);
209 return static_cast<float>( xTmp1.polarAngle(xTmp2) );
210}

◆ EvalFourMuInvMass()

template<class T>
float EventAnalysis::EvalFourMuInvMass ( const T * pxP1,
const T * pxP2,
const T * pxP3,
const T * pxP4 )
static

Definition at line 124 of file EventAnalysis.h.

125{
126 // Check integrity of inputs.
127 if ( !pxP1 || !pxP2 || !pxP3 || !pxP4) return invalidAnswer;
128
129 // Evaluate invariant mass.
130 return EvalInvMass( pxP1, pxP2, pxP3, pxP4, EAna::g_fMuonMass );
131}

◆ EvalInvMass() [1/2]

template<class T>
float EventAnalysis::EvalInvMass ( const T * pxP1,
const T * pxP2,
const T * pxp3,
const T * pxP4,
float fMass1,
float fMass2 = -999.9,
float fMass3 = -999.9,
float fMass4 = invalidAnswer )
static

Definition at line 133 of file EventAnalysis.h.

135{
136 // Check integrity of inputs.No tachyons.
137 if ( fMass1 < 0.0f ) return invalidAnswer;
138 if ( !pxP1 || !pxP2 || !pxP3 || !pxP4) return invalidAnswer;
139
140 // Set masses equal if required by user
141 fMass2 = ( fMass2 < 0.0f ) ? fMass1 : fMass2;
142 fMass3 = ( fMass3 < 0.0f ) ? fMass1 : fMass3;
143 fMass4 = ( fMass4 < 0.0f ) ? fMass1 : fMass4;
144
145 // Evaluate invariant mass.
146 CLHEP::Hep3Vector xTmp1 = CLHEP::Hep3Vector( pxP1->p4().Px() * EAna::CGeV, pxP1->p4().Py() * EAna::CGeV, pxP1->p4().Pz() * EAna::CGeV );
147 CLHEP::Hep3Vector xTmp2 = CLHEP::Hep3Vector( pxP2->p4().Px() * EAna::CGeV, pxP2->p4().Py() * EAna::CGeV, pxP2->p4().Pz() * EAna::CGeV );
148 CLHEP::Hep3Vector xTmp3 = CLHEP::Hep3Vector( pxP3->p4().Px() * EAna::CGeV, pxP3->p4().Py() * EAna::CGeV, pxP3->p4().Pz() * EAna::CGeV );
149 CLHEP::Hep3Vector xTmp4 = CLHEP::Hep3Vector( pxP4->p4().Px() * EAna::CGeV, pxP4->p4().Py() * EAna::CGeV, pxP4->p4().Pz() * EAna::CGeV );
150
151 CLHEP::HepLorentzVector xLVec; xLVec.setVectM ( xTmp1, fMass1 );
152 CLHEP::HepLorentzVector xLVec2; xLVec2.setVectM( xTmp2, fMass2 );
153 CLHEP::HepLorentzVector xLVec3; xLVec3.setVectM( xTmp3, fMass3 );
154 CLHEP::HepLorentzVector xLVec4; xLVec4.setVectM( xTmp4, fMass4 );
155
156 xLVec += xLVec2;
157 xLVec += xLVec3;
158 xLVec += xLVec4;
159
160 return static_cast<float>( xLVec.m() );
161}

◆ EvalInvMass() [2/2]

template<class T>
float EventAnalysis::EvalInvMass ( const T * pxP1,
const T * pxP2,
float fMass1,
float fMass2 = invalidAnswer )
static

Definition at line 107 of file EventAnalysis.h.

109{
110 // Check integrity of inputs.No tachyons.
111 if ( fMass1 < 0.0f ) return invalidAnswer;
112 if ( !pxP1 || !pxP2 ) return invalidAnswer;
113 // Set masses equal if required by user
114 fMass2 = ( fMass2 < 0.0f ) ? fMass1 : fMass2;
115 // Evaluate invariant mass.
116 CLHEP::Hep3Vector xTmp1 = calculateMomentum(pxP1);
117 CLHEP::Hep3Vector xTmp2 = calculateMomentum(pxP2);
118 CLHEP::HepLorentzVector xLVec1; xLVec1.setVectM( xTmp1, fMass1 );
119 CLHEP::HepLorentzVector xLVec2; xLVec2.setVectM( xTmp2, fMass2 );
120 return static_cast<float>( xLVec1.invariantMass( xLVec2 ) );
121}

◆ EvalPhi()

template<class T>
float EventAnalysis::EvalPhi ( const T * pxP1,
const T * pxP2 )
static

Definition at line 221 of file EventAnalysis.h.

222{
223 // Check integrity of inputs.
224 if ( !pxP1 || !pxP2 ) return invalidAnswer;
225 CLHEP::Hep3Vector xTmp1 = calculateMomentum(pxP1);
226 CLHEP::Hep3Vector xTmp2 = calculateMomentum(pxP2);
227 return static_cast<float>( (xTmp1 + xTmp2).phi() );
228}

◆ EvalPhiDiff()

template<class T>
float EventAnalysis::EvalPhiDiff ( const T * pxP1,
const T * pxP2 )
static

Definition at line 191 of file EventAnalysis.h.

192{
193 // Check integrity of inputs.
194 if ( !pxP1 || !pxP2 ) return invalidAnswer;
195 // Evaluate the angle.
196 CLHEP::Hep3Vector xTmp1 = calculateMomentum(pxP1);
197 CLHEP::Hep3Vector xTmp2 = calculateMomentum(pxP2);
198 return static_cast<float>( xTmp1.deltaPhi(xTmp2) );
199}

◆ EvalPt()

template<class T>
float EventAnalysis::EvalPt ( const T * pxP1,
const T * pxP2 )
static

Definition at line 212 of file EventAnalysis.h.

213{
214 // Check integrity of inputs.
215 if ( !pxP1 || !pxP2 ) return invalidAnswer;
216 CLHEP::Hep3Vector xTmp1 = calculateMomentum(pxP1);
217 CLHEP::Hep3Vector xTmp2 = calculateMomentum(pxP2);
218 return static_cast<float>( (xTmp1 + xTmp2).perp() );
219}

◆ EvalPtDiff()

template<class T>
float EventAnalysis::EvalPtDiff ( const T * pxP1,
const T * pxP2 )
static

Definition at line 174 of file EventAnalysis.h.

175{
176 // Check integrity of inputs.
177 if ( !pxP1 || !pxP2 ) return invalidAnswer;
178 // Evaluate the difference between the momenta. Signed using positive - negative if appropriate.
179 if ( pxP1->charge() > 0.5f )
180 {
181 return static_cast<float>( pxP1->pt() * EAna::CGeV - pxP2->pt() * EAna::CGeV );
182 }
183 else
184 {
185 return static_cast<float>( pxP2->pt() * EAna::CGeV - pxP1->pt() * EAna::CGeV );
186 }
187}

◆ EvalTransverseMass() [1/4]

template<class T>
float EventAnalysis::EvalTransverseMass ( const T * pxP1,
const T * pxP2,
float fMETx,
float fMETy )
static

Definition at line 272 of file EventAnalysis.h.

273{
274 // Check integrity of inputs.
275 if ( !pxP1 || !pxP2 ) return invalidAnswer;
276 // Evaluate Di-mu invariant mass.
277 return EvalTransverseMass( pxP1, pxP2, fMETx, fMETy, EAna::g_fMuonMass );
278}
static float EvalTransverseMass(const T *pxP1, float fMETx, float fMETy, float fMass1, float fMass2=invalidAnswer)

◆ EvalTransverseMass() [2/4]

template<class T>
float EventAnalysis::EvalTransverseMass ( const T * pxP1,
const T * pxP2,
float fMETx,
float fMETy,
float fMass1,
float fMass2 = invalidAnswer )
static

Definition at line 280 of file EventAnalysis.h.

282{
283 // Check integrity of inputs.No tachyons.
284 if ( fMass1 < 0.0f ) return invalidAnswer;
285 if ( !pxP1 || !pxP2 ) return invalidAnswer;
286 // Set masses equal if required by user.
287 fMass2 = ( fMass2 < 0.0f ) ? fMass1 : fMass2;
288 // Evaluate invariant mass.
289 CLHEP::Hep3Vector xTmp1 = CLHEP::Hep3Vector( pxP1->p4().Px() * EAna::CGeV, pxP1->p4().Py() * EAna::CGeV, 0.0f );
290 CLHEP::Hep3Vector xTmp2 = CLHEP::Hep3Vector( pxP2->p4().Px() * EAna::CGeV, pxP2->p4().Py() * EAna::CGeV, 0.0f );
291 CLHEP::Hep3Vector xTmp12 = xTmp1 + xTmp2;
292 CLHEP::Hep3Vector xTmp3 = CLHEP::Hep3Vector( fMETx, fMETy, 0.0f );
293 CLHEP::HepLorentzVector xLVec1; xLVec1.setVectM( xTmp12, fMass1 );
294 CLHEP::HepLorentzVector xLVec2; xLVec2.setVectM( xTmp3, 0.0f );
295 return static_cast<float>( xLVec1.invariantMass( xLVec2 ) );
296}

◆ EvalTransverseMass() [3/4]

template<class T>
float EventAnalysis::EvalTransverseMass ( const T * pxP1,
float fMETx,
float fMETy )
static

Definition at line 247 of file EventAnalysis.h.

248{
249 // Check integrity of inputs.
250 if ( !pxP1 ) return invalidAnswer;
251 // Evaluate Di-mu invariant mass.
252 return EvalInvMass( pxP1, fMETx, fMETy, EAna::g_fMuonMass );
253}

◆ EvalTransverseMass() [4/4]

template<class T>
float EventAnalysis::EvalTransverseMass ( const T * pxP1,
float fMETx,
float fMETy,
float fMass1,
float fMass2 = invalidAnswer )
static

Definition at line 255 of file EventAnalysis.h.

257{
258 // Check integrity of inputs.No tachyons.
259 if ( fMass1 < 0.0f ) return invalidAnswer;
260 if ( !pxP1 ) return invalidAnswer;
261 // Set masses equal if required by user.
262 fMass2 = ( fMass2 < 0.0f ) ? fMass1 : fMass2;
263 // Evaluate invariant mass.
264 CLHEP::Hep3Vector xTmp1 = CLHEP::Hep3Vector( pxP1->p4().Px() * EAna::CGeV, pxP1->p4().Py() * EAna::CGeV, 0.0f );
265 CLHEP::Hep3Vector xTmp2 = CLHEP::Hep3Vector( fMETx, fMETy, 0.0f );
266 CLHEP::HepLorentzVector xLVec1; xLVec1.setVectM( xTmp1, fMass1 );
267 CLHEP::HepLorentzVector xLVec2; xLVec2.setVectM( xTmp2, 0.0f );
268 return static_cast<float>( xLVec1.invariantMass( xLVec2 ) );
269}

◆ EvaluateAngle()

template<class T>
float EventAnalysis::EvaluateAngle ( const T * pxP1,
const T * pxP2 )
static

Definition at line 164 of file EventAnalysis.h.

165{
166 // Check integrity of inputs.
167 if ( !pxP1 || !pxP2 ) return invalidAnswer;
168 // Evaluate the angle.
169 CLHEP::Hep3Vector xTmp1 = calculateMomentum(pxP1);
170 CLHEP::Hep3Vector xTmp2 = calculateMomentum(pxP2);
171 return static_cast<float>( xTmp1.angle(xTmp2) );
172}

◆ Init()

void EventAnalysis::Init ( )
virtual

Reimplemented in ElectronSelector, FourMuonEvent, MuonSelector, and ZmumuEvent.

Definition at line 44 of file EventAnalysis.cxx.

45{
46 // This must be called by an inheriting class in order to book
47 // & register histograms.
49 Register();
50}
virtual void BookHistograms()

◆ Register()

void EventAnalysis::Register ( )
private

Definition at line 63 of file EventAnalysis.cxx.

64{
65 ServiceHandle<ITHistSvc> histSvc ("THistSvc", "EventAnalysis");
66
67 // Register histograms in monitoring tool
68 registerHistogramType(*histSvc, m_x1DHistograms, m_xSampleName, "/1dhisto_");
69 registerHistogramType(*histSvc, m_x2DHistograms, m_xSampleName, "/2dhisto_");
70
71 registerHistogramType(*histSvc, m_x1DProfHistograms, m_xSampleName, "/1dprof_");
72 registerHistogramType(*histSvc, m_x2DProfHistograms, m_xSampleName, "/2dprof_");
73}
std::map< unsigned int, TProfile2D * > m_x2DProfHistograms
std::map< unsigned int, TH2F * > m_x2DHistograms
std::map< unsigned int, TH1F * > m_x1DHistograms
std::map< unsigned int, TProfile * > m_x1DProfHistograms
std::string m_xSampleName

Member Data Documentation

◆ invalidAnswer

float EventAnalysis::invalidAnswer {-999.9f}
staticconstexpr

Definition at line 40 of file EventAnalysis.h.

40{-999.9f};

◆ m_uPassedEvents

unsigned int EventAnalysis::m_uPassedEvents
protected

Definition at line 75 of file EventAnalysis.h.

◆ m_x1DHistograms

std::map<unsigned int, TH1F*> EventAnalysis::m_x1DHistograms
protected

Definition at line 76 of file EventAnalysis.h.

◆ m_x1DProfHistograms

std::map<unsigned int, TProfile*> EventAnalysis::m_x1DProfHistograms
protected

Definition at line 78 of file EventAnalysis.h.

◆ m_x2DHistograms

std::map<unsigned int, TH2F*> EventAnalysis::m_x2DHistograms
protected

Definition at line 77 of file EventAnalysis.h.

◆ m_x2DProfHistograms

std::map<unsigned int, TProfile2D*> EventAnalysis::m_x2DProfHistograms
protected

Definition at line 79 of file EventAnalysis.h.

◆ m_xSampleName

std::string EventAnalysis::m_xSampleName
protected

Definition at line 81 of file EventAnalysis.h.


The documentation for this class was generated from the following files: