ATLAS Offline Software
Loading...
Searching...
No Matches
Qjets.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5#include "Qjets.h"
6
7using namespace std;
8using namespace JetSubStructureUtils;
9
10Qjets::Qjets(double zcut, double dcut_fctr, double exp_min, double exp_max, double rigidity, double truncation_fctr)
11: m_rand(),
12 m_zcut(zcut),
13 m_dcut(-1.),
14 m_dcut_fctr(dcut_fctr),
15 m_exp_min(exp_min),
16 m_exp_max(exp_max),
17 m_rigidity(rigidity),
18 m_truncation_fctr(truncation_fctr)
19{
20}
21
22void Qjets::SetRandSeed(unsigned int seed){
23 m_rand.SetSeed(seed);
24}
25
26bool Qjets::JetUnmerged(int num) const{
27 return m_merged_jets.find(num) == m_merged_jets.end();
28}
29
30bool Qjets::JetsUnmerged(const jet_distance& jd) const{
31 return JetUnmerged(jd.j1) && JetUnmerged(jd.j2);
32}
33
35 vector< pair<jet_distance, double> > popped_distances;
36 double norm(0.);
37 jet_distance ret{};
38 ret.j1 = -1;
39 ret.j2 = -1;
40 ret.dij = -1.;
41 bool dmin_set(false);
42 double dmin(0.);
43
44 while(!m_distances.empty()){
45 jet_distance dist = m_distances.top();
46 m_distances.pop();
47 if(JetsUnmerged(dist)){
48 if(!dmin_set){
49 dmin = dist.dij;
50 dmin_set = true;
51 }
52 double weight = exp(-m_rigidity * (dist.dij-dmin) /dmin);
53 popped_distances.emplace_back(dist,weight);
54 norm += weight;
55 if(weight/norm < m_truncation_fctr)
56 break;
57 }
58 }
59
60 double rand(Rand()), tot_weight(0.);
61 if (!popped_distances.empty()) {
62 const double inv_norm = norm == 0 ? 1 : 1. / norm;
63 for(vector<pair<jet_distance, double> >::iterator it = popped_distances.begin(); it != popped_distances.end(); ++it){
64 tot_weight += (*it).second * inv_norm;
65 if(tot_weight >= rand){
66 ret = (*it).first;
67 break;
68 }
69 }
70 }
71
72 // repopulate in reverse (maybe quicker?)
73 for(vector<pair<jet_distance, double> >::reverse_iterator it = popped_distances.rbegin(); it != popped_distances.rend(); ++it)
74 if(JetsUnmerged((*it).first))
75 m_distances.push((*it).first);
76
77 return ret;
78}
79
80void Qjets::Cluster(fastjet::ClusterSequence & cs){
81 ComputeDCut(cs);
82
83 // Populate all the distances
84 ComputeAllDistances(cs.jets());
86
87 while(!m_distances.empty() && jd.dij != -1.){
88 if(!Prune(jd,cs)){
89 // m_merged_jets.push_back(jd.j1);
90 // m_merged_jets.push_back(jd.j2);
91 m_merged_jets[jd.j1] = true;
92 m_merged_jets[jd.j2] = true;
93
94 int new_jet;
95 cs.plugin_record_ij_recombination(jd.j1, jd.j2, 1., new_jet);
96 assert(JetUnmerged(new_jet));
98 } else {
99 double j1pt = cs.jets()[jd.j1].perp();
100 double j2pt = cs.jets()[jd.j2].perp();
101 if(j1pt>j2pt){
102 // m_merged_jets.push_back(jd.j2);
103 m_merged_jets[jd.j2] = true;
104 cs.plugin_record_iB_recombination(jd.j2, 1.);
105 } else {
106 // m_merged_jets.push_back(jd.j1);
107 m_merged_jets[jd.j1] = true;
108 cs.plugin_record_iB_recombination(jd.j1, 1.);
109 }
110 }
111 jd = GetNextDistance();
112 }
113
114 // merge remaining jets with beam
115 [[maybe_unused]] int num_merged_final(0);
116
117 for(unsigned int i = 0 ; i < cs.jets().size(); i++)
118 if(JetUnmerged(i)){
119 num_merged_final++;
120 cs.plugin_record_iB_recombination(i,1.);
121 }
122
123 assert(num_merged_final < 2);
124}
125
126void Qjets::ComputeDCut(fastjet::ClusterSequence & cs){
127 // assume all jets in cs form a single jet. compute mass and pt
128 fastjet::PseudoJet sum(0.,0.,0.,0.);
129 for(vector<fastjet::PseudoJet>::const_iterator it = cs.jets().begin(); it != cs.jets().end(); ++it)
130 sum += (*it);
131 m_dcut = 2. * m_dcut_fctr * sum.m()/sum.perp();
132}
133
134bool Qjets::Prune(jet_distance& jd,fastjet::ClusterSequence & cs){
135 double pt1 = cs.jets()[jd.j1].perp();
136 double pt2 = cs.jets()[jd.j2].perp();
137 fastjet::PseudoJet sum_jet = cs.jets()[jd.j1]+cs.jets()[jd.j2];
138 double sum_pt = sum_jet.perp();
139 double z = min(pt1,pt2)/sum_pt;
140 double d = sqrt(cs.jets()[jd.j1].plain_distance(cs.jets()[jd.j2]));
141 return (d > m_dcut) && (z < m_zcut);
142}
143
144void Qjets::ComputeAllDistances(const vector<fastjet::PseudoJet>& inp){
145 for(unsigned int i = 0 ; i < inp.size()-1; i++){
146 // jet-jet distances
147 for(unsigned int j = i+1 ; j < inp.size(); j++){
148 jet_distance jd{};
149 jd.j1 = i;
150 jd.j2 = j;
151 if(jd.j1 != jd.j2){
152 jd.dij = d_ij(inp[i],inp[j]);
153 m_distances.push(jd);
154 }
155 }
156 }
157}
158
159void Qjets::ComputeNewDistanceMeasures(fastjet::ClusterSequence & cs, int new_jet){
160 // jet-jet distances
161 for(unsigned int i = 0; i < cs.jets().size(); i++)
162 if(JetUnmerged(i) && ((int)i) != new_jet){
163 jet_distance jd{};
164 jd.j1 = new_jet;
165 jd.j2 = i;
166 jd.dij = d_ij(cs.jets()[jd.j1], cs.jets()[jd.j2]);
167 m_distances.push(jd);
168 }
169}
170
171double Qjets::d_ij(const fastjet::PseudoJet& v1,const fastjet::PseudoJet& v2) const{
172 double p1 = v1.perp();
173 double p2 = v2.perp();
174 double ret = pow(min(p1,p2),m_exp_min) * pow(max(p1,p2),m_exp_max) * v1.squared_distance(v2);
175 assert(!std::isnan(ret));
176 return ret;
177}
178
179double Qjets::Rand(){
180 return m_rand.Rndm();
181}
#define z
constexpr int pow(int base, int exp) noexcept
#define min(a, b)
Definition cfImp.cxx:40
#define max(a, b)
Definition cfImp.cxx:41
void ComputeNewDistanceMeasures(fastjet::ClusterSequence &cs, int new_jet)
Definition Qjets.cxx:159
double d_ij(const fastjet::PseudoJet &v1, const fastjet::PseudoJet &v2) const
Definition Qjets.cxx:171
void SetRandSeed(unsigned int seed)
Definition Qjets.cxx:22
void ComputeAllDistances(const std::vector< fastjet::PseudoJet > &inp)
Definition Qjets.cxx:144
bool JetsUnmerged(const jet_distance &jd) const
Definition Qjets.cxx:30
std::priority_queue< jet_distance, std::vector< jet_distance >, JetDistanceCompare > m_distances
Definition Qjets.h:36
bool JetUnmerged(int num) const
Definition Qjets.cxx:26
jet_distance GetNextDistance()
Definition Qjets.cxx:34
std::map< int, bool > m_merged_jets
Definition Qjets.h:35
void ComputeDCut(fastjet::ClusterSequence &cs)
Definition Qjets.cxx:126
void Cluster(fastjet::ClusterSequence &cs)
Definition Qjets.cxx:80
bool Prune(jet_distance &jd, fastjet::ClusterSequence &cs)
Definition Qjets.cxx:134
Qjets(double zcut, double dcut_fctr, double exp_min, double exp_max, double rigidity, double truncation_fctr)
Definition Qjets.cxx:10
STL namespace.