ATLAS Offline Software
Loading...
Searching...
No Matches
JetMSVAugmentation.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
6
12
13
14namespace DerivationFramework {
15
16
17JetMSVAugmentation::JetMSVAugmentation(const std::string& t, const std::string& n, const IInterface* p):
18 base_class(t,n,p)
19{
20}
21
22
23
25
26
27
29 ATH_MSG_DEBUG("Initialize " );
30 ATH_CHECK(m_jetCollectionName.initialize());
31
32 return StatusCode::SUCCESS;
33
34}
35
36
37
38StatusCode JetMSVAugmentation::addBranches(const EventContext& ctx) const{
39
41 if ( !jets.isValid() ) {
42 ATH_MSG_ERROR ("Couldn't retrieve jets with key: " << m_jetCollectionName );
43 return StatusCode::FAILURE;
44 }
45
46
52
56
59
60
61 for (auto jet : *jets) {
63 if (!bjet) {
64 ATH_MSG_WARNING("btagging information not available" );
65 continue;
66 }
67
68 std::vector< ElementLink< xAOD::VertexContainer > > msvVertices;
69 bjet->variable<std::vector<ElementLink<xAOD::VertexContainer> > >(m_vtxAlgName, "vertices", msvVertices);
70
71 std::vector<float> vtx_mass;
72 std::vector<float> vtx_pt;
73 std::vector<float> vtx_eta;
74 std::vector<float> vtx_phi;
75 std::vector<float> vtx_efrac;
76 std::vector<float> vtx_x;
77 std::vector<float> vtx_y;
78 std::vector<float> vtx_z;
79 std::vector<int> vtx_ntrk;
80 std::vector<float> vtx_dls;
81
82
83 for (auto vtx : msvVertices) {//loop in vertices
84
85 int ntrk = xAOD::SecVtxHelper::VtxNtrk(*vtx);
86 float mass = xAOD::SecVtxHelper::VertexMass(*vtx);
87 float efrc = xAOD::SecVtxHelper::EnergyFraction(*vtx);
88 float pt = xAOD::SecVtxHelper::Vtxpt(*vtx);
89 float eta = xAOD::SecVtxHelper::Vtxeta(*vtx);
90 float phi = xAOD::SecVtxHelper::Vtxphi(*vtx);
91 float dls = xAOD::SecVtxHelper::VtxnormDist(*vtx);
92 float xp = (*vtx)->x();
93 float yp = (*vtx)->y();
94 float zp = (*vtx)->z();
95 // float chi = (*vtx)->chiSquared();
96 // float ndf = (*vtx)->numberDoF();
97
98 TLorentzVector p;
99 p.SetPtEtaPhiM(pt,eta,phi,mass);
100
101 vtx_mass.push_back(mass);
102 vtx_pt.push_back(pt);
103 vtx_eta.push_back(eta);
104 vtx_phi.push_back(phi);
105 vtx_efrac.push_back(efrc);
106 vtx_x.push_back(xp);
107 vtx_y.push_back(yp);
108 vtx_z.push_back(zp);
109 vtx_ntrk.push_back(ntrk);
110 vtx_dls.push_back(dls);
111 }
112
113 dec_vtxmass(*bjet)=vtx_mass;
114 dec_vtxpt(*bjet)=vtx_pt;
115 dec_vtxeta(*bjet)=vtx_eta;
116 dec_vtxphi(*bjet)=vtx_phi;
117 dec_vtxefrac(*bjet)=vtx_efrac;
118 dec_vtxx(*bjet)=vtx_x;
119 dec_vtxy(*bjet)=vtx_y;
120 dec_vtxz(*bjet)=vtx_z;
121 dec_vtxdls(*bjet)=vtx_dls;
122 dec_vtxntrk(*bjet)=vtx_ntrk;
123
124 }
125
126 return StatusCode::SUCCESS;
127
128}
129
130
131
132}
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
tool to add some MSV variables to jets
Handle class for adding a decoration to an object.
JetMSVAugmentation(const std::string &t, const std::string &n, const IInterface *p)
virtual StatusCode addBranches(const EventContext &ctx) const override final
SG::WriteDecorHandleKey< xAOD::JetContainer > m_dec_vtxpt
virtual StatusCode initialize() override final
SG::ReadHandleKey< xAOD::JetContainer > m_jetCollectionName
Gaudi::Property< std::string > m_vtxAlgName
SG::WriteDecorHandleKey< xAOD::JetContainer > m_dec_vtxmass
SG::WriteDecorHandleKey< xAOD::JetContainer > m_dec_vtxx
SG::WriteDecorHandleKey< xAOD::JetContainer > m_dec_vtxntrk
SG::WriteDecorHandleKey< xAOD::JetContainer > m_dec_vtxz
SG::WriteDecorHandleKey< xAOD::JetContainer > m_dec_vtxy
SG::WriteDecorHandleKey< xAOD::JetContainer > m_dec_vtxeta
SG::WriteDecorHandleKey< xAOD::JetContainer > m_dec_vtxphi
SG::WriteDecorHandleKey< xAOD::JetContainer > m_dec_vtxefrac
SG::WriteDecorHandleKey< xAOD::JetContainer > m_dec_vtxdls
Handle class for adding a decoration to an object.
THE reconstruction tool.
const BTagging * getBTagging(const SG::AuxElement &part)
Access the default xAOD::BTagging object associated to an object.
int VtxNtrk(const xAOD::Vertex *)
float Vtxpt(const xAOD::Vertex *)
float Vtxeta(const xAOD::Vertex *)
float EnergyFraction(const xAOD::Vertex *)
float VertexMass(const xAOD::Vertex *)
float Vtxphi(const xAOD::Vertex *)
float VtxnormDist(const xAOD::Vertex *)
BTagging_v1 BTagging
Definition of the current "BTagging version".
Definition BTagging.h:17