ATLAS Offline Software
Loading...
Searching...
No Matches
BPhysHelper.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3*/
4
35 * @code
36 * void myFunction(xAOD::Vertex* vtx) {
37 * // Let "vtx" be some xAOD::Vertex created by the b-physics software
38 * // We gain access to augmentations through the helper class:
39 * xAOD::BPhysHelper bVtx(vtx);
40 *
41 * std::cout << "Refitted pT of the 1st track: ";
42 * if(bVtx.nRefTrks()>0) {
43 * std::cout << bVtx.refTrk(0).Pt();
44 * }
45 * std::cout << std::endl;
46 *
47 * }
48 * @endcode
49 *
50 */
51
52#ifndef XAOD_BPHYSHELPER_H
53#define XAOD_BPHYSHELPER_H
54
60
61
62#include "TVector3.h"
63#include "TLorentzVector.h"
64#include "TMatrixTSym.h"
65
66#include <assert.h>
67
68
70namespace xAOD {
72
73
74 public:
75
76 /************************************************************************/
84
86 m_b(b),
87 m_covCached(false),
88 m_cachedCov(0),
89 m_refTracksCached(false),
91 m_muonsCached(false),
92 m_electronsCached(false),
99 {
100 assert(m_b!=0); // sanity check: m_b must not be null
101 }
102
103 /************************************************************************/
107
108 const xAOD::Vertex* vtx() const { return m_b; }
109
110 /************************************************************************/
122
123 const TMatrixTSym<double>& covariance();
124
125 /************************************************************************/
130
134
135 int nRefTrks();
136
141
142 TVector3 refTrk(const size_t index);
143
149
150 const std::vector<TVector3>& refTrks();
151
157
158 TLorentzVector refTrk(const size_t index, const float mass);
159
185
186 const xAOD::IParticle* refTrkOrigin(const size_t index) const;
187
197
198 TVector3 refTrkOriginP(const size_t index) const;
199
215
216 TLorentzVector refTrkOriginP(const size_t index, const float mass) const;
217
227
228 float refTrkCharge(const size_t index) const;
229
234
235 bool setRefTrks(std::vector<float> px,
236 std::vector<float> py,
237 std::vector<float> pz);
238
243
244 bool setRefTrks(const std::vector<TVector3>& refTrks);
245
246#ifndef XAOD_ANALYSIS
247
262
263 bool setRefTrks();
264
265
266#endif // not XAOD_ANALYSIS
267
269
270 /************************************************************************/
274
285
286 TVector3 totalP();
287
294
295 TLorentzVector totalP(std::span<const double> masses);
296
301
302 float ptErr();
303
309
310 bool setPtErr(const float val);
311
313
314 /************************************************************************/
322
326
327 int nMuons();
328 int nElectrons();
329
334
335 const xAOD::Muon* muon(const size_t index);
336 const xAOD::Electron* electron(const size_t index);
337
341
342 const std::vector<const xAOD::Muon*>& muons();
343 const std::vector<const xAOD::Electron*>& electrons();
344
350
351 bool setMuons(const std::vector<const xAOD::Muon*>& muons,
353
354 bool setElectrons(const std::vector<const xAOD::Electron*>& electrons,
357
358 /************************************************************************/
370
375
376 int nPrecedingVertices();
377
383
384 const xAOD::Vertex* precedingVertex(const size_t index);
385
390
391 const std::vector<const xAOD::Vertex*>& precedingVertices();
392
398
399 bool setPrecedingVertices(const std::vector<const xAOD::Vertex*>& vertices,
400 const xAOD::VertexContainer* vertexContainer);
401
402
404
405 /************************************************************************/
414
419
420 int nCascadeVertices();
421
427
428 const xAOD::Vertex* cascadeVertex(const size_t index);
429
434
435 const std::vector<const xAOD::Vertex*>& cascadeVertices();
436
442
443 bool setCascadeVertices(const std::vector<const xAOD::Vertex*>& vertices,
444 const xAOD::VertexContainer* vertexContainer);
445
451
452 int nRefTrksCascade();
453
455
456 /************************************************************************/
474
476
477 /************************************************************************/
484
490
493
501 bool setOrigPv(const xAOD::Vertex* pv,
502 const xAOD::VertexContainer* vertexContainer,
504
512
513 bool setPv(const xAOD::Vertex* pv,
514 const xAOD::VertexContainer* vertexContainer,
516
523
530
532
533 /************************************************************************/
537
544
547
555
556 bool setLxy (const float val, const pv_type vertexType = BPhysHelper::PV_MIN_A0);
557 bool setLxyErr(const float val, const pv_type vertexType = BPhysHelper::PV_MIN_A0);
558
561 bool setLxyz (const float val, const pv_type vertexType = BPhysHelper::PV_MIN_A0);
562 bool setLxyzErr(const float val, const pv_type vertexType = BPhysHelper::PV_MIN_A0);
563
565
566 /************************************************************************/
570
576
583
590
591 float setA0 (const float val, const pv_type vertexType = BPhysHelper::PV_MIN_A0);
592 float setA0Err (const float val, const pv_type vertexType = BPhysHelper::PV_MIN_A0);
593 float setA0xy (const float val, const pv_type vertexType = BPhysHelper::PV_MIN_A0);
594 float setA0xyErr(const float val, const pv_type vertexType = BPhysHelper::PV_MIN_A0);
595 float setZ0 (const float val, const pv_type vertexType = BPhysHelper::PV_MIN_A0);
596 float setZ0Err (const float val, const pv_type vertexType = BPhysHelper::PV_MIN_A0);
597
599
600 /************************************************************************/
604
605 static const std::string pv_type_str[];
607
608 /************************************************************************/
612
613 static const unsigned int n_pv_types;
615
616 /************************************************************************/
617 /************************************************************************/
618 protected:
619
627
628 bool cacheCov();
629
637
638 bool cacheRefTracks();
639
647
648 bool cacheMuons();
649 bool cacheElectrons();
650
651
659
661
669
671
672 /************************************************************************/
675
676 /************************************************************************/
680
682 TMatrixTSym<double> m_cachedCov;
683
685
686 /************************************************************************/
690
692 std::vector<TVector3> m_cachedRefTracks;
693
695
696 /************************************************************************/
700
703 std::vector<const xAOD::Muon*> m_cachedMuons;
704 std::vector<const xAOD::Electron*> m_cachedElectrons;
705
707
708 /************************************************************************/
712
714 std::vector<const xAOD::Vertex*> m_cachedPrecedingVertices;
716 std::vector<const xAOD::Vertex*> m_cachedCascadeVertices;
717
719
720 /************************************************************************/
724
725 static const std::vector<TVector3> s_emptyVectorOfTVector3;
726 static const std::vector<const xAOD::Muon*> s_emptyVectorOfMuons;
727 static const std::vector<const xAOD::Electron*> s_emptyVectorOfElectrons;
728 static const TMatrixTSym<double> s_emptyMatrix;
729 static const std::vector<const xAOD::Vertex*> s_emptyVectorOfVertices;
731
732 }; // BPhysHelper
733
734} // namespace xAOD
735
738#define BPHYS_CHECK( EXP ) { if( ! EXP ) ATH_MSG_WARNING ( "Call of \"" << #EXP << "\" failed" ); }
739
740#endif // XAOD_BPHYSHELPER_H
xAOD::ElectronContainer * electronContainer
xAOD::MuonContainer * muonContainer
static const std::vector< const xAOD::Electron * > s_emptyVectorOfElectrons
bool setPtErr(const float val)
Set pT error.
int nMuons()
: Methods providing access to the linked muons
const std::vector< const xAOD::Muon * > & muons()
Returns linked muons.
bool cacheCascadeVertices()
: Cache cascade vertices
bool setElectrons(const std::vector< const xAOD::Electron * > &electrons, const xAOD::ElectronContainer *electronContainer)
float lxyzErr(const pv_type vertexType=BPhysHelper::PV_MIN_A0)
its error
float refTrkCharge(const size_t index) const
Returns charge of the i-th track.
float setZ0(const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
longitudinal impact parameter
int nCascadeVertices()
: Links to cascade vertices
const xAOD::Vertex * origPv(const pv_type vertexType=BPhysHelper::PV_MIN_A0)
original PV
const xAOD::Vertex * cascadeVertex(const size_t index)
Returns pointer to a cascade vertex.
std::vector< const xAOD::Vertex * > m_cachedCascadeVertices
float z0(const pv_type vertexType=BPhysHelper::PV_MIN_A0)
longitudinal impact parameter
const xAOD::Vertex * m_b
Cached B decay xAOD vertex.
float lxyz(const pv_type vertexType=BPhysHelper::PV_MIN_A0)
decay distance
float a0xyErr(const pv_type vertexType=BPhysHelper::PV_MIN_A0)
rtansverse impact parameter error
bool cacheMuons()
: Cache linked muons
bool setLxyErr(const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
its error
const TMatrixTSym< double > & covariance()
: Returns full covariance matrix
const std::vector< const xAOD::Vertex * > & precedingVertices()
Returns vector of pointers to preceding vertices.
bool setRefTrks()
: Sets refitted track momenta
bool setLxy(const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
Set the transverse decay distance and its error measured between the refitted primary vertex of type ...
float lxy(const pv_type vertexType=BPhysHelper::PV_MIN_A0)
Get the transverse decay distance and its error measured between the refitted primary vertex of type ...
static const std::vector< const xAOD::Muon * > s_emptyVectorOfMuons
static const unsigned int n_pv_types
TMatrixTSym< double > m_cachedCov
bool cacheRefTracks()
: Cache refitted tracks
static const std::vector< const xAOD::Vertex * > s_emptyVectorOfVertices
bool setMuons(const std::vector< const xAOD::Muon * > &muons, const xAOD::MuonContainer *muonContainer)
Set links to muons.
const xAOD::Vertex * precedingVertex(const size_t index)
Returns pointer to a preceding vertex.
std::vector< const xAOD::Electron * > m_cachedElectrons
float setZ0Err(const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
longitudinal impact parameter error
const xAOD::Vertex * vtx() const
Getter method for the cached xAOD::Vertex.
const std::vector< const xAOD::Electron * > & electrons()
bool cachePrecedingVertices()
: Cache preceding vertices
static const std::string pv_type_str[]
Definition BPhysHelper.h:35
static const TMatrixTSym< double > s_emptyMatrix
bool setPrecedingVertices(const std::vector< const xAOD::Vertex * > &vertices, const xAOD::VertexContainer *vertexContainer)
Sets links to preceding vertices.
float setA0(const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
Set the 3D and transverse impact parameters and their error.
const std::vector< TVector3 > & refTrks()
Returns refitted track momenta.
bool setPv(const xAOD::Vertex *pv, const xAOD::VertexContainer *vertexContainer, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
Set the refitted collision vertex of type pv_type.
float z0Err(const pv_type vertexType=BPhysHelper::PV_MIN_A0)
longitudinal impact parameter error
const xAOD::Vertex * pv(const pv_type vertexType=BPhysHelper::PV_MIN_A0)
Get the refitted collision vertex of type pv_type.
TVector3 refTrk(const size_t index)
Returns i-th refitted track 3-momentum.
bool m_precedingVerticesCached
pv_type
: Enum type of the PV
float a0Err(const pv_type vertexType=BPhysHelper::PV_MIN_A0)
3D impact parameter error
BPhysHelper(const xAOD::Vertex *b)
: Main constructor
Definition BPhysHelper.h:85
float setA0Err(const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
3D impact parameter error
bool setCascadeVertices(const std::vector< const xAOD::Vertex * > &vertices, const xAOD::VertexContainer *vertexContainer)
Sets links to cascade vertices.
int nPrecedingVertices()
: Links to preceding vertices
float ptErr()
Returns pT error.
bool cacheCov()
: Cache covariance matrix
bool setRefitPVStatus(int code, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
Set the exitCode of the refitter for vertex of type pv_type.
std::vector< TVector3 > m_cachedRefTracks
const std::vector< const xAOD::Vertex * > & cascadeVertices()
Returns vector of pointers to cascade vertices.
float a0xy(const pv_type vertexType=BPhysHelper::PV_MIN_A0)
transverse impact parameter
TVector3 refTrkOriginP(const size_t index) const
Returns perigee 3-momentum of the original track corresponding i-th refitted track.
std::vector< const xAOD::Vertex * > m_cachedPrecedingVertices
float setA0xyErr(const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
transverse impact parameter error
int nRefTrks()
Returns number of stored refitted track momenta.
std::vector< const xAOD::Muon * > m_cachedMuons
float setA0xy(const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
transverse impact parameter
TVector3 totalP()
: Returns total 3-momentum calculated from the refitted tracks
int RefitPVStatus(const pv_type vertexType=BPhysHelper::PV_MIN_A0)
Get the exitCode of the refitter for vertex of type pv_type.
const xAOD::IParticle * refTrkOrigin(const size_t index) const
: Returns the original track (charged or neutral) corresponding to the i-th refitted track
int nRefTrksCascade()
Returns number of stored refitted tracks INCLUDING those from the linked cascade vertices.
bool setOrigPv(const xAOD::Vertex *pv, const xAOD::VertexContainer *vertexContainer, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
Set the original collision vertex of type pv_type.
float lxyErr(const pv_type vertexType=BPhysHelper::PV_MIN_A0)
its error
bool setLxyz(const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
decay distance
static const std::vector< TVector3 > s_emptyVectorOfTVector3
bool setLxyzErr(const float val, const pv_type vertexType=BPhysHelper::PV_MIN_A0)
its error
Class providing the definition of the 4-vector interface.
double a0
Definition globals.cxx:27
Definition index.py:1
ICaloAffectedTool is abstract interface for tools checking if 4 mom is in calo affected region.
ElectronContainer_v1 ElectronContainer
Definition of the current "electron container version".
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
Vertex_v1 Vertex
Define the latest version of the vertex class.
Muon_v1 Muon
Reference the current persistent version:
MuonContainer_v1 MuonContainer
Definition of the current "Muon container version".
Electron_v1 Electron
Definition of the current "egamma version".