CMS 3D CMS Logo

Qjets.cc
Go to the documentation of this file.
2 
3 using namespace std;
4 
5 void Qjets::SetRandSeed(unsigned int seed) {
6  _rand_seed_set = true;
7  _seed = seed;
8 }
9 
10 bool Qjets::JetUnmerged(int num) const { return _merged_jets.find(num) == _merged_jets.end(); }
11 
12 bool Qjets::JetsUnmerged(const JetDistance& jd) const { return JetUnmerged(jd.j1) && JetUnmerged(jd.j2); }
13 
15  vector<pair<JetDistance, double> > popped_distances;
16  double norm(0.);
18  ret.j1 = -1;
19  ret.j2 = -1;
20  ret.dij = -1.;
21  bool dmin_set(false);
22  double dmin(0.);
23 
24  while (!_distances.empty()) {
25  JetDistance dist = _distances.top();
26  _distances.pop();
27  if (JetsUnmerged(dist)) {
28  if (!dmin_set) {
29  dmin = dist.dij;
30  dmin_set = true;
31  }
32  double weight = exp(-_rigidity * (dist.dij - dmin) / dmin);
33  popped_distances.push_back(make_pair(dist, weight));
34  norm += weight;
35  if (weight / norm < _truncation_fctr)
36  break;
37  }
38  }
39 
40  double rand(Rand()), tot_weight(0.);
41  for (vector<pair<JetDistance, double> >::iterator it = popped_distances.begin(); it != popped_distances.end(); it++) {
42  tot_weight += (*it).second / norm;
43  if (tot_weight >= rand) {
44  ret = (*it).first;
45  omega *= ((*it).second);
46  break;
47  }
48  }
49 
50  // repopulate in reverse (maybe quicker?)
51  for (vector<pair<JetDistance, double> >::reverse_iterator it = popped_distances.rbegin();
52  it != popped_distances.rend();
53  it++)
54  if (JetsUnmerged((*it).first))
55  _distances.push((*it).first);
56 
57  return ret;
58 }
59 
60 void Qjets::Cluster(fastjet::ClusterSequence& cs) {
61  omega = 1.;
62  QjetsBaseExtras* extras = new QjetsBaseExtras();
63  computeDCut(cs);
64 
65  // Populate all the distances
66  ComputeAllDistances(cs.jets());
67  JetDistance jd = GetNextDistance();
68 
69  while (!_distances.empty() && jd.dij != -1.) {
70  if (!Prune(jd, cs)) {
71  // _merged_jets.push_back(jd.j1);
72  // _merged_jets.push_back(jd.j2);
73  _merged_jets[jd.j1] = true;
74  _merged_jets[jd.j2] = true;
75 
76  int new_jet = -1;
77  cs.plugin_record_ij_recombination(jd.j1, jd.j2, 1., new_jet);
78  if (!JetUnmerged(new_jet))
79  edm::LogError("QJets Clustering") << "Problem with FastJet::plugin_record_ij_recombination";
80  ComputeNewDistanceMeasures(cs, new_jet);
81  } else {
82  double j1pt = cs.jets()[jd.j1].perp();
83  double j2pt = cs.jets()[jd.j2].perp();
84  if (j1pt > j2pt) {
85  // _merged_jets.push_back(jd.j2);
86  _merged_jets[jd.j2] = true;
87  cs.plugin_record_iB_recombination(jd.j2, 1.);
88  } else {
89  // _merged_jets.push_back(jd.j1);
90  _merged_jets[jd.j1] = true;
91  cs.plugin_record_iB_recombination(jd.j1, 1.);
92  }
93  }
94  jd = GetNextDistance();
95  }
96 
97  extras->_wij = omega;
98  cs.plugin_associate_extras(extras);
99 
100  // merge remaining jets with beam
101  int num_merged_final(0);
102  for (unsigned int i = 0; i < cs.jets().size(); i++)
103  if (JetUnmerged(i)) {
104  num_merged_final++;
105  cs.plugin_record_iB_recombination(i, 1.);
106  }
107 
108  if (!(num_merged_final < 2))
109  edm::LogError("QJets Clustering") << "More than 1 jet remains after reclustering";
110 }
111 
112 void Qjets::computeDCut(fastjet::ClusterSequence& cs) {
113  // assume all jets in cs form a single jet. compute mass and pt
114  fastjet::PseudoJet sum(0., 0., 0., 0.);
115  for (vector<fastjet::PseudoJet>::const_iterator it = cs.jets().begin(); it != cs.jets().end(); it++)
116  sum += (*it);
117  _dcut = 2. * _dcut_fctr * sum.m() / sum.perp();
118 }
119 
120 bool Qjets::Prune(JetDistance& jd, fastjet::ClusterSequence& cs) {
121  double pt1 = cs.jets()[jd.j1].perp();
122  double pt2 = cs.jets()[jd.j2].perp();
123  fastjet::PseudoJet sum_jet = cs.jets()[jd.j1] + cs.jets()[jd.j2];
124  double sum_pt = sum_jet.perp();
125  double z = min(pt1, pt2) / sum_pt;
126  double d2 = cs.jets()[jd.j1].plain_distance(cs.jets()[jd.j2]);
127  return (d2 > _dcut * _dcut) && (z < _zcut);
128 }
129 
130 void Qjets::ComputeAllDistances(const vector<fastjet::PseudoJet>& inp) {
131  for (unsigned int i = 0; i < inp.size() - 1; i++) {
132  // jet-jet distances
133  for (unsigned int j = i + 1; j < inp.size(); j++) {
134  JetDistance jd;
135  jd.j1 = i;
136  jd.j2 = j;
137  if (jd.j1 != jd.j2) {
138  jd.dij = d_ij(inp[i], inp[j]);
139  _distances.push(jd);
140  }
141  }
142  }
143 }
144 
145 void Qjets::ComputeNewDistanceMeasures(fastjet::ClusterSequence& cs, unsigned int new_jet) {
146  // jet-jet distances
147  for (unsigned int i = 0; i < cs.jets().size(); i++)
148  if (JetUnmerged(i) && i != new_jet) {
149  JetDistance jd;
150  jd.j1 = new_jet;
151  jd.j2 = i;
152  jd.dij = d_ij(cs.jets()[jd.j1], cs.jets()[jd.j2]);
153  _distances.push(jd);
154  }
155 }
156 
157 double Qjets::d_ij(const fastjet::PseudoJet& v1, const fastjet::PseudoJet& v2) const {
158  double p1 = v1.perp();
159  double p2 = v2.perp();
160  double ret = pow(min(p1, p2), _exp_min) * pow(max(p1, p2), _exp_max) * v1.squared_distance(v2);
161  if (edm::isNotFinite(ret))
162  edm::LogError("QJets Clustering") << "d_ij is not finite";
163  return ret;
164 }
165 
166 double Qjets::Rand() {
167  double ret = 0.;
168  //if(_rand_seed_set)
169  // ret = rand_r(&_seed)/(double)RAND_MAX;
170  //else
171  //ret = rand()/(double)RAND_MAX;
172  ret = _rnEngine->flat();
173  return ret;
174 }
runTheMatrix.ret
ret
prodAgent to be discontinued
Definition: runTheMatrix.py:543
QjetsBaseExtras::_wij
double _wij
Definition: Qjets.h:82
Qjets.h
mps_fire.i
i
Definition: mps_fire.py:428
Qjets::SetRandSeed
void SetRandSeed(unsigned int seed)
Definition: Qjets.cc:5
Qjets::Prune
bool Prune(JetDistance &jd, fastjet::ClusterSequence &cs)
Definition: Qjets.cc:120
edm::isNotFinite
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
fwrapper::cs
unique_ptr< ClusterSequence > cs
Definition: fastjetfortran_madfks.cc:47
min
T min(T a, T b)
Definition: MathUtil.h:58
Qjets::JetsUnmerged
bool JetsUnmerged(const JetDistance &jd) const
Definition: Qjets.cc:12
mps_merge.weight
weight
Definition: mps_merge.py:88
Qjets::computeDCut
void computeDCut(fastjet::ClusterSequence &cs)
Definition: Qjets.cc:112
fileCollector.seed
seed
Definition: fileCollector.py:127
JetDistance::j2
int j2
Definition: Qjets.h:19
p2
double p2[4]
Definition: TauolaWrapper.h:90
JetDistance::dij
double dij
Definition: Qjets.h:17
HLT_FULL_cff.pt1
pt1
Definition: HLT_FULL_cff.py:9893
Qjets::GetNextDistance
JetDistance GetNextDistance()
Definition: Qjets.cc:14
Qjets::d_ij
double d_ij(const fastjet::PseudoJet &v1, const fastjet::PseudoJet &v2) const
Definition: Qjets.cc:157
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
QjetsBaseExtras
Definition: Qjets.h:74
JetDistance::j1
int j1
Definition: Qjets.h:18
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
p1
double p1[4]
Definition: TauolaWrapper.h:89
JetDistance
Definition: Qjets.h:16
HLT_FULL_cff.pt2
pt2
Definition: HLT_FULL_cff.py:9895
edm::LogError
Log< level::Error, false > LogError
Definition: MessageLogger.h:123
EgammaValidation_cff.num
num
Definition: EgammaValidation_cff.py:33
std
Definition: JetResolutionObject.h:76
Qjets::JetUnmerged
bool JetUnmerged(int num) const
Definition: Qjets.cc:10
Qjets::Rand
double Rand()
Definition: Qjets.cc:166
Qjets::ComputeNewDistanceMeasures
void ComputeNewDistanceMeasures(fastjet::ClusterSequence &cs, unsigned int new_jet)
Definition: Qjets.cc:145
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
Qjets::ComputeAllDistances
void ComputeAllDistances(const std::vector< fastjet::PseudoJet > &inp)
Definition: Qjets.cc:130
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
JetChargeProducer_cfi.exp
exp
Definition: JetChargeProducer_cfi.py:6
Qjets::Cluster
void Cluster(fastjet::ClusterSequence &cs)
Definition: Qjets.cc:60
weight
Definition: weight.py:1