ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/UHHAnalysis/NtupleWriter/src/Qjets.cxx
Revision: 1.1
Committed: Wed Jun 19 13:22:06 2013 UTC (11 years, 10 months ago) by rkogler
Content type: text/plain
Branch: MAIN
CVS Tags: HEAD
Error occurred while calculating annotation data.
Log Message:
added substructure information and jet constituents

File Contents

# Content
1 #include "UHHAnalysis/NtupleWriter/interface/Qjets.h"
2
3 Qjets::Qjets(double zcut, double dcut_fctr, double exp_min, double exp_max, double rigidity)
4 : _zcut(zcut), _dcut_fctr(dcut_fctr), _exp_min(exp_min), _exp_max(exp_max), _rigidity(rigidity)
5 {
6 _dcut = -1.;
7 _rand_seed_set = false;
8 }
9
10 void Qjets::SetRandSeed(unsigned int seed){
11 _rand_seed_set = true;
12 _seed = seed;
13 }
14
15 bool Qjets::JetUnmerged(int num){
16 if(count(_merged_jets.begin(), _merged_jets.end(), num) > 0)
17 return false;
18 return true;
19 }
20
21 bool Qjets::JetsUnmerged(jet_distance& jd){
22 return JetUnmerged(jd.j1) && JetUnmerged(jd.j2);
23 }
24
25 void Qjets::ComputeNewDistanceMeasures(fastjet::ClusterSequence & cs, int new_jet){
26 // jet-jet distances
27 for(unsigned int i = 0; i < cs.jets().size(); i++)
28 if(JetUnmerged(i) && i != (unsigned int) new_jet){
29 jet_distance jd;
30 jd.j1 = new_jet;
31 jd.j2 = i;
32 jd.dij = d_ij(cs.jets()[jd.j1], cs.jets()[jd.j2]);
33 _distances.push_back(jd);
34 }
35 }
36
37 void Qjets::ComputeDCut(fastjet::ClusterSequence & cs){
38 // assume all jets in cs form a single jet. compute mass and pt
39 fastjet::PseudoJet sum(0.,0.,0.,0.);
40 for(vector<fastjet::PseudoJet>::const_iterator it = cs.jets().begin(); it != cs.jets().end(); it++)
41 sum += (*it);
42 _dcut = 2. * _dcut_fctr * sum.m()/sum.perp();
43 }
44
45 double Qjets::ComputeMinimumDistance(){
46 double dmin(-1.);
47 for(list<jet_distance>::iterator it = _distances.begin(); it != _distances.end(); )
48 if(JetsUnmerged(*it)){
49 if(dmin == -1. || (*it).dij < dmin)
50 dmin = (*it).dij;
51 it++;
52 } else
53 it = _distances.erase(it);
54 return dmin;
55 }
56
57 double Qjets::ComputeNormalization(double dmin){
58 double norm(0.);
59 for(list<jet_distance>::iterator it = _distances.begin(); it != _distances.end(); )
60 if(JetsUnmerged(*it)){
61 double inc = 0.;
62 /*
63 We need to be careful about dmin == 0 which happens for collinear jets.
64 (Thanks go to Nhan V Tran for finding this bug)
65 */
66 if(dmin == 0.){
67 if((*it).dij == 0.)
68 inc = 1.;
69 else
70 inc = 0.;
71 } else
72 inc = exp(-_rigidity*((*it).dij-dmin)/dmin);
73
74 assert(inc <= 1. && !isnan((float) inc));
75 norm += inc;
76 it++;
77 } else
78 it = _distances.erase(it);
79
80 return norm;
81 }
82
83 void Qjets::Cluster(fastjet::ClusterSequence & cs){
84 ComputeDCut(cs);
85 ComputeAllDistances(cs.jets());
86 while (!_distances.empty()){
87 double dmin = ComputeMinimumDistance();
88 double norm = ComputeNormalization(dmin);
89
90 // sometimes if the rigidity is too large the norm is zero. make sure this is not the case
91 if(_distances.size() == 0)
92 break;
93 assert(norm > 0.);
94
95 // Now compute a random number between 0 and 1 and find the corresponding measure
96 double rand = Rand();
97 double sum = 0.;
98
99 for(list<jet_distance>::iterator it = _distances.begin(); it != _distances.end(); it++){
100 /*
101 We need to be careful about dmin == 0 which happens for collinear jets.
102 (Thanks go to Nhan V Tran for finding this bug)
103 */
104 if(dmin == 0.){
105 if((*it).dij == 0.)
106 sum += 1./norm;
107 } else
108 sum += exp(-_rigidity*((*it).dij-dmin)/dmin)/norm;
109 assert(!isnan((float)sum));
110
111 if(sum > rand){
112 if(!Prune((*it),cs)){
113 _merged_jets.push_back((*it).j1);
114 _merged_jets.push_back((*it).j2);
115 int new_jet;
116 cs.plugin_record_ij_recombination((*it).j1, (*it).j2, 1., new_jet);
117 assert(JetUnmerged(new_jet));
118 ComputeNewDistanceMeasures(cs,new_jet);
119 } else {
120 double j1pt = cs.jets()[(*it).j1].perp();
121 double j2pt = cs.jets()[(*it).j2].perp();
122 if(j1pt>j2pt){
123 _merged_jets.push_back((*it).j2);
124 cs.plugin_record_iB_recombination((*it).j2, 1.);
125 } else {
126 _merged_jets.push_back((*it).j1);
127 cs.plugin_record_iB_recombination((*it).j1, 1.);
128 }
129 }
130 break;
131 }
132 }
133 }
134 // merge remaining jets with beam
135 int num_merged_final(0);
136 for(unsigned int i = 0 ; i < cs.jets().size(); i++)
137 if(JetUnmerged(i)){
138 num_merged_final++;
139 cs.plugin_record_iB_recombination(i,1.);
140 }
141 assert(num_merged_final < 2);
142 }
143
144 bool Qjets::Prune(jet_distance& jd,fastjet::ClusterSequence & cs){
145 double pt1 = cs.jets()[jd.j1].perp();
146 double pt2 = cs.jets()[jd.j2].perp();
147 fastjet::PseudoJet sum_jet = cs.jets()[jd.j1]+cs.jets()[jd.j2];
148 double sum_pt = sum_jet.perp();
149 double z = min(pt1,pt2)/sum_pt;
150 double d = sqrt(cs.jets()[jd.j1].plain_distance(cs.jets()[jd.j2]));
151 return (d > _dcut) && (z < _zcut);
152 }
153
154 void Qjets::ComputeAllDistances(const vector<fastjet::PseudoJet>& inp){
155 for(unsigned int i = 0 ; i < inp.size()-1; i++){
156 // jet-jet distances
157 for(unsigned int j = i+1 ; j < inp.size(); j++){
158 jet_distance jd;
159 jd.j1 = i;
160 jd.j2 = j;
161 if(jd.j1 != jd.j2){
162 jd.dij = d_ij(inp[i],inp[j]);
163 _distances.push_back(jd);
164 }
165 }
166 }
167 }
168
169 double Qjets::d_ij(const fastjet::PseudoJet& v1,const fastjet::PseudoJet& v2){
170 double p1 = v1.perp();
171 double p2 = v2.perp();
172 double ret = pow(min(p1,p2),_exp_min) * pow(max(p1,p2),_exp_max) * v1.squared_distance(v2);
173 assert(!isnan((float)ret));
174 return ret;
175 }
176
177 double Qjets::Rand(){
178 double ret = 0.;
179 if(_rand_seed_set)
180 ret = rand_r(&_seed)/(double)RAND_MAX;
181 else
182 ret = rand()/(double)RAND_MAX;
183 return ret;
184 }