ATLAS Offline Software
Loading...
Searching...
No Matches
EGTransverseMassTool.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// EGTransverseMassTool.cxx
7// Author: Giovanni Marchiori (giovanni.marchiori@cern.ch)
8//
9// The tool computes the transverse mass between the MET and the
10// particles belonging to some collection.
11// The MET met and phi values can either be stored in SG under
12// keys specified with METPtBranchName and METPhiBranchName,
13// otherwise the met() and phi() of the MET object stored in
14// SG under the METContainerName (default = MET_LocHadTopo)
15// will be used.
16// A MET>METmin cut can be applied.
17// For the particles, only those passing the ObjectRequirements
18// selection are retained. Their pt and phi can be stored
19// in SG under keys specified with ObjectPtBranchName and
20// ObjectPhiBranchName, otherwise the Pt() and Phi() methods
21// of the IParticles object of the containred stored in SG
22// with the key specified by ObjectContainerName will be used
24
26
27#include <cmath>
28using std::abs;
29using std::sqrt;
30
31#include "TLorentzVector.h"
32
33namespace DerivationFramework {
34
37{
38 if (m_sgName.key().empty()) {
40 "No SG name provided for the output of the transverse mass tool!");
41 return StatusCode::FAILURE;
42 }
43 ATH_CHECK(m_sgName.initialize());
44
45 if (!m_container1Name.key().empty()) {
46 ATH_CHECK(m_container1Name.initialize());
47 }
48 if (!m_container2Name.key().empty()) {
49 ATH_CHECK(m_container2Name.initialize());
50 }
51 if (!m_pt1BranchName.key().empty()) {
52 ATH_CHECK(m_pt1BranchName.initialize());
53 }
54 if (!m_phi1BranchName.key().empty()) {
55 ATH_CHECK(m_phi1BranchName.initialize());
56 }
57 if (!m_pt2BranchName.key().empty()) {
58 ATH_CHECK(m_pt2BranchName.initialize());
59 }
60 if (!m_phi2BranchName.key().empty()) {
61 ATH_CHECK(m_phi2BranchName.initialize());
62 }
63
64 ATH_CHECK(initializeParser(m_expression1));
65
66 return StatusCode::SUCCESS;
67}
68
69StatusCode
70EGTransverseMassTool::addBranches(const EventContext& ctx) const
71{
73
74 // create the vector which will hold the values invariant masses
75 auto masses = std::make_unique<std::vector<float>>();
76 // compute the invariant mass values
77 ATH_CHECK(getTransverseMasses(ctx, *masses));
78
79 ATH_CHECK(writeHandle.record(std::move(masses)));
80
81 return StatusCode::SUCCESS;
82}
83
84StatusCode
86 std::vector<float>& masses) const
87{
88 // Get optional payload
89 const std::vector<float>* pt1 = nullptr;
90 if (!m_pt1BranchName.key().empty()) {
92 pt1 = readHandle.ptr();
93 }
94 const std::vector<float>* pt2 = nullptr;
95 if (!m_pt2BranchName.key().empty()) {
97 pt2 = readHandle.ptr();
98 }
99
100 const std::vector<float>* phi1 = nullptr;
101 if (!m_phi1BranchName.key().empty()) {
103 phi1 = readHandle.ptr();
104 }
105 const std::vector<float>* phi2 = nullptr;
106 if (!m_phi2BranchName.key().empty()) {
108 phi2 = readHandle.ptr();
109 }
110
111 // Get the input particle and MET
113 ctx };
115 ctx };
116 const xAOD::IParticleContainer* particles1 = inputParticles1.ptr();
117 const xAOD::MissingETContainer* particles2 = inputParticles2.ptr();
118
119 // compute MET
120 if (particles2->empty()) {
121 if (!pt2) {
122 ATH_MSG_WARNING("No MET info found");
123 return StatusCode::SUCCESS;
124 }
125 if (pt2->empty()) {
126 ATH_MSG_WARNING("No MET info found");
127 return StatusCode::SUCCESS;
128 }
129 }
130
131 float MET = pt2 ? (*pt2)[0] : particles2->at(0)->met();
132 float MET_phi = phi2 ? (*phi2)[0] : particles2->at(0)->phi();
133
134 // apply MET requirement
135 if (MET < m_METmin)
136 return StatusCode::SUCCESS;
137
138 // get the positions of the particles which pass the requirement
139 std::vector<int> entries1 = m_parser->evaluateAsVector();
140 unsigned int nEntries1 = entries1.size();
141
142 // if there are no particles in one of the two lists to combine, just leave
143 // function
144 if (nEntries1 == 0)
145 return StatusCode::SUCCESS;
146
147 // check that the sizes are compatible
148 if (particles1->size() != nEntries1) {
149 ATH_MSG_ERROR("Branch sizes incompatible - returning zero");
150 return StatusCode::FAILURE;
151 }
152 if ((pt1 && pt1->size() != nEntries1) ||
153 (phi1 && phi1->size() != nEntries1)) {
154 ATH_MSG_ERROR("Branch sizes incompatible - returning zero");
155 return StatusCode::FAILURE;
156 }
157
158 // Loop over the objects and calculate the transverse mass
159 unsigned int iter;
160 for (iter = 0; iter < nEntries1; ++iter) {
161 if (entries1[iter] != 1)
162 continue;
163 float apt1 = pt1 ? (*pt1)[iter] : ((*particles1)[iter])->p4().Pt();
164 float aphi1 = phi1 ? (*phi1)[iter] : ((*particles1)[iter])->p4().Phi();
165
166 TLorentzVector v1, v2;
167 v1.SetPtEtaPhiM(apt1, 0., aphi1, m_mass1Hypothesis);
168 v2.SetPtEtaPhiM(MET, 0., MET_phi, 0.);
169 float mass = (v1 + v2).M();
170 masses.push_back(mass);
171 }
172 return StatusCode::SUCCESS;
173}
174}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_WARNING(x)
const T * at(size_type n) const
Access an element, as an rvalue.
size_type size() const noexcept
Returns the number of elements in the collection.
bool empty() const noexcept
Returns true if the collection is empty.
SG::ReadHandleKey< std::vector< float > > m_phi1BranchName
SG::ReadHandleKey< std::vector< float > > m_phi2BranchName
virtual StatusCode addBranches(const EventContext &ctx) const override final
SG::ReadHandleKey< std::vector< float > > m_pt1BranchName
SG::ReadHandleKey< std::vector< float > > m_pt2BranchName
SG::ReadHandleKey< xAOD::IParticleContainer > m_container1Name
Gaudi::Property< std::string > m_expression1
SG::WriteHandleKey< std::vector< float > > m_sgName
StatusCode getTransverseMasses(const EventContext &ctx, std::vector< float > &) const
SG::ReadHandleKey< xAOD::MissingETContainer > m_container2Name
virtual StatusCode initialize() override final
const_pointer_type ptr()
Dereference the pointer.
StatusCode record(std::unique_ptr< T > data)
Record a const object to the store.
THE reconstruction tool.
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition MET.py:1
DataVector< IParticle > IParticleContainer
Simple convenience declaration of IParticleContainer.