ATLAS Offline Software
Loading...
Searching...
No Matches
LundVariablesTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4#include "fastjet/PseudoJet.hh"
6#include "fastjet/tools/Recluster.hh"
7#include "fastjet/contrib/SoftDrop.hh"
8#include "fastjet/JetDefinition.hh"
9#include <iostream>
10#include <queue>
12
15
17 fastjet::PseudoJet jet;
18 bool decorate = SetupDecoration(jet, injet);
19 int n_splits = 0;
20 std::vector<float> lund_all_lnR;
21 std::vector<float> lund_all_lnkT;
22 std::vector<float> lund_all_z;
23 std::vector<float> lund_all_kt;
24 std::vector<float> lund_all_deltaR;
25 std::vector<int> lund_all_idp1;
26 std::vector<int> lund_all_idp2;
27
28 if (decorate) {
29
30 std::vector<std::vector<LundVariablesTool::Declustering>> v_LundValues;
31 auto & constit_links = injet.constituentLinks();
32 std::vector<fastjet::PseudoJet> v_pj_constituents;
33
34 v_pj_constituents.clear();
35
36 for(size_t ij = 0; ij < injet.numConstituents(); ij++) {
37 const xAOD::FlowElement* constit = dynamic_cast<const xAOD::FlowElement*>(*constit_links[ij]);
38 if (constit) {
39 TLorentzVector tlv_const;
40 tlv_const.SetPtEtaPhiE(constit->pt() / 1.e3, constit->eta(), constit->phi(), constit->e() / 1.e3);
41
42 v_pj_constituents.push_back(fastjet::PseudoJet(tlv_const.Px(), tlv_const.Py(), tlv_const.Pz(), tlv_const.E()));
43
44 }
45 else {
46 ANA_MSG_WARNING("Failed to cast truth particle ------");
47 }
48 }
49 v_LundValues.push_back(LundVariablesTool::getLundVar(std::move(v_pj_constituents)));
50
51
52 std::vector<double> v_jj_pt;
53 std::vector<double> v_j1_pt;
54 std::vector<double> v_j2_pt;
55 std::vector<double> v_jj_kt;
56
57
58 for (const auto & declust : v_LundValues.at(0)){
59 v_jj_pt.push_back( declust.jj.pt() );
60 v_j1_pt.push_back( declust.j1.pt() );
61 v_j2_pt.push_back( declust.j2.pt() );
62 v_jj_kt.push_back( declust.kt );
63 }
64
65 for (int id=0; id < int(v_jj_pt.size()); id++) {
66
67 int idp1 = -1;
68 int idp2 = -1;
69
70 auto itp1 = find(v_jj_pt.begin(), v_jj_pt.end(), v_j1_pt.at(id));
71 auto itp2 = find(v_jj_pt.begin(), v_jj_pt.end(), v_j2_pt.at(id));
72
73
74 if (itp1 != v_jj_pt.end()) {
75 idp1 = itp1 - v_jj_pt.begin();
76
77 if (idp1 == id){ // in case second parent has a very very small pT
78 idp1 = -1;
79 }
80
81
82 }
83 if (itp2 != v_jj_pt.end()) {
84 idp2 = itp2 - v_jj_pt.begin();
85 }
86
87 v_LundValues.at(0).at(id).idp1 = idp1;
88 v_LundValues.at(0).at(id).idp2 = idp2;
89 }
90
91
92 n_splits = 0;
93 for (const auto& split : v_LundValues[0]) {
94 if (split.delta_R > 1e-8 && split.kt > 1e-8) {
95 // save declustering data
96 lund_all_lnR.push_back(std::log(1.0 / split.delta_R));
97 lund_all_lnkT.push_back(std::log(split.kt));
98 lund_all_z.push_back(split.z);
99 lund_all_kt.push_back(split.kt);
100 lund_all_deltaR.push_back(split.delta_R);
101 lund_all_idp1.push_back(split.idp1);
102 lund_all_idp2.push_back(split.idp2);
103 n_splits++;
104 }
105 }
106 }
107
108
109 // Decore the jet with the ljp variables
110 injet.setAttribute(m_prefix + "LundAllLnR", lund_all_lnR);
111 injet.setAttribute(m_prefix + "LundAllLnKT", lund_all_lnkT);
112 injet.setAttribute(m_prefix + "LundAllZ", lund_all_z);
113 injet.setAttribute(m_prefix + "LundAllKt", lund_all_kt);
114 injet.setAttribute(m_prefix + "LundAllDeltaR", lund_all_deltaR);
115 injet.setAttribute(m_prefix + "LundAllIDP1", lund_all_idp1);
116 injet.setAttribute(m_prefix + "LundAllIDP2", lund_all_idp2);
117 injet.setAttribute(m_prefix + "nSplits", n_splits);
118
119 return 0;
120}
121
124 ATH_MSG_INFO("Calculating Lund Plane Variables");
125}
126
127std::vector<LundVariablesTool::Declustering> LundVariablesTool::getLundVar( std::vector<fastjet::PseudoJet> v_jcs)
128{
129 fastjet::JetDefinition jd(fastjet::cambridge_algorithm, 1.0);
130 std::vector<LundVariablesTool::Declustering> result;
131 fastjet::ClusterSequence cs(v_jcs, jd);
132 std::vector<fastjet::PseudoJet> v_pj = sorted_by_pt(cs.inclusive_jets());
133
134 fastjet::PseudoJet j = v_pj.at(0);
135
136 std::queue< fastjet::PseudoJet > jetStack;
137 jetStack.push(j);
138
139 while ( jetStack.size() > 0 ) {
140 fastjet::PseudoJet thisJ, pJLeft, pJRight;
141 thisJ = jetStack.front();
142 jetStack.pop();
143
144 bool thisJHasParents = thisJ.has_parents(pJLeft, pJRight);
145
146 if (!thisJHasParents){
147 continue;
148 }
149
150 if (pJLeft.pt2() < pJRight.pt2()) {
151 // std::cout << pJLeft.pt2() << "\t" << pJRight.pt2() << std::endl;
152 fastjet::PseudoJet jTemp;
153 jTemp = pJLeft;
154 pJLeft = pJRight;
155 pJRight = jTemp;
156 }
157
158 jetStack.push(pJLeft);
159 jetStack.push(pJRight);
160
162
163 declust.jj = thisJ;
164 declust.j1 = pJLeft;
165 declust.j2 = pJRight;
166
167 // get info about the jet
168 declust.pt = thisJ.pt();
169 declust.m = thisJ.m();
170
171 // collect info about the declustering
172 declust.pt1 = pJLeft.pt();
173 declust.pt2 = pJRight.pt();
174 declust.eta = pJRight.eta();
175 declust.E = pJRight.E();
176 declust.delta_R = pJLeft.delta_R(pJRight);
177 declust.z = declust.pt2 / (declust.pt1 + declust.pt2);
178 declust.kt = pJRight.pt() * declust.delta_R;
179
180
181 declust.varphi = atan2(pJLeft.rap() - pJRight.rap(), pJLeft.delta_phi_to(pJRight));
182 result.push_back(std::move(declust));
183 }
184
185
186
187 return result;
188
189}
#define ATH_MSG_INFO(x)
#define ANA_MSG_WARNING(xmsg)
Macro printing warning messages.
JetSubStructureMomentToolsBase(const std::string &name)
bool SetupDecoration(fastjet::PseudoJet &pseudojet, const xAOD::Jet &jet, bool requireJetStructure=false) const
virtual void print() const
Print the state of the tool.
static std::vector< Declustering > getLundVar(std::vector< fastjet::PseudoJet > v_jcs)
void print() const override
Print the state of the tool.
int modifyJet(xAOD::Jet &injet) const override
Modify a single jet. This is obsolete and set to be removed.
LundVariablesTool(const std::string &name)
virtual double pt() const override
virtual double phi() const override
The azimuthal angle ( ) of the particle.
virtual double eta() const override
The pseudorapidity ( ) of the particle.
virtual double e() const override
The total energy of the particle.
void setAttribute(const std::string &name, const T &v)
size_t numConstituents() const
Number of constituents in this jets (this is valid even when reading a file where the constituents ha...
Definition Jet_v1.cxx:153
const std::vector< ElementLink< IParticleContainer > > & constituentLinks() const
Direct access to constituents. WARNING expert use only.
Definition Jet_v1.cxx:162
std::string find(const std::string &s)
return a remapped string
Definition hcg.cxx:138
std::vector< std::string > split(const std::string &s, const std::string &t=":")
Definition hcg.cxx:177
Jet_v1 Jet
Definition of the current "jet version".
FlowElement_v1 FlowElement
Definition of the current "pfo version".
Definition FlowElement.h:16