ATLAS Offline Software
Loading...
Searching...
No Matches
PhotonVertexSelectionWrapper.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
7#include <vector>
8
9namespace DerivationFramework {
10
13{
14 ATH_MSG_DEBUG("Initializing " << name() << "...");
16 ATH_MSG_DEBUG("Retrieved tool " << m_photonPointingTool);
17
18 ATH_CHECK(m_photonContainer.initialize());
19 ATH_CHECK(m_vertexContainer.initialize());
20
21 ATH_CHECK(m_vtxPt.initialize());
22 ATH_CHECK(m_vtxEta.initialize());
23 ATH_CHECK(m_vtxPhi.initialize());
24 ATH_CHECK(m_vtxSumPt.initialize());
25 ATH_CHECK(m_vtxSumPt2.initialize());
26
27 ATH_MSG_DEBUG("Initialization successful");
28
29 return StatusCode::SUCCESS;
30}
31
32StatusCode
33PhotonVertexSelectionWrapper::addBranches(const EventContext& ctx) const
34{
35 // retrieve the input containers
38
39 // Update calo pointing auxdata for photons
40 if (m_photonPointingTool->updatePointingAuxdata(*photons).isFailure()) {
41 ATH_MSG_ERROR("Couldn't update photon calo pointing auxdata");
42 return StatusCode::FAILURE;
43 }
44
45 // create the decorators
51 bool isMomentum_available = vtxPt.isAvailable();
52 bool isSumPt_available = vtxSumPt.isAvailable();
53 bool isSumPt2_available = vtxSumPt2.isAvailable();
54 bool found_PV = false;
55
56 // Loop over vertices and update auxdata
57 for (const auto *vertex : *vertices) {
58
59 float pt = -999.;
60 float eta = -999.;
61 float phi = -999.;
62 float sumPt = -999.;
63 float sumPt2 = -999.;
64
65 if(vertex->vertexType() == xAOD::VxType::VertexType::PriVtx or
66 vertex->vertexType() == xAOD::VxType::VertexType::PileUp){
67
68 found_PV = true;
69 // Get momentum vector of vertex and add it as a decoration
70 TLorentzVector vmom = xAOD::PVHelpers::getVertexMomentum(vertex, isMomentum_available);
71 pt = sqrt(vmom.Px() * vmom.Px() + vmom.Py() * vmom.Py());
72 if(pt>0.){
73 eta = asinh(vmom.Pz() / pt);
74 phi = acos(vmom.Px() / pt); // in [0,Pi]
75 }
76 // Calculate additional quantities
77 sumPt = xAOD::PVHelpers::getVertexSumPt(vertex, 1, isSumPt_available);
78 sumPt2 = xAOD::PVHelpers::getVertexSumPt(vertex, 2, isSumPt2_available);
79
80 }
81
82 // write decorations
83 if(!isMomentum_available){
84 vtxPt(*vertex) = pt;
85 vtxEta(*vertex) = eta;
86 vtxPhi(*vertex) = phi;
87 }
88 if(!isSumPt_available) vtxSumPt(*vertex) = sumPt;
89 // For events where no PV is reconstructed, beamspot is saved and isSumPt2_available = False
90 // Need to also check found_PV to avoid modifying locked store
91 if(!isSumPt2_available && found_PV) vtxSumPt2(*vertex) = sumPt2;
92
93 } // end loop o vertices
94
95 return StatusCode::SUCCESS;
96}
97
98}
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_DEBUG(x)
SG::WriteDecorHandleKey< xAOD::VertexContainer > m_vtxSumPt
SG::WriteDecorHandleKey< xAOD::VertexContainer > m_vtxPhi
ToolHandle< CP::IPhotonPointingTool > m_photonPointingTool
PhotonPointingTool.
SG::WriteDecorHandleKey< xAOD::VertexContainer > m_vtxSumPt2
SG::ReadHandleKey< xAOD::PhotonContainer > m_photonContainer
Input photon container.
SG::ReadHandleKey< xAOD::VertexContainer > m_vertexContainer
Input primary vertex container.
SG::WriteDecorHandleKey< xAOD::VertexContainer > m_vtxEta
SG::WriteDecorHandleKey< xAOD::VertexContainer > m_vtxPt
virtual StatusCode addBranches(const EventContext &ctx) const override final
Handle class for adding a decoration to an object.
bool isAvailable()
Test to see if this variable exists in the store, for the referenced object.
THE reconstruction tool.
static const SG::AuxElement::Decorator< float > sumPt2("sumPt2")
::StatusCode StatusCode
StatusCode definition for legacy code.
float getVertexSumPt(const xAOD::Vertex *vertex, int power=1, bool useAux=true)
Loop over track particles associated with vertex and return scalar sum of pT^power in GeV (from auxda...
TLorentzVector getVertexMomentum(const xAOD::Vertex *vertex, bool useAux=true, const std::string &derivationPrefix="")
Return vector sum of tracks associated with vertex (from auxdata if available and useAux = true)
@ PileUp
Pile-up vertex.
@ PriVtx
Primary vertex.