Fork me on GitHub

source: svn/trunk/Utilities/Fastjet/src/PseudoJet.cc@ 634

Last change on this file since 634 was 11, checked in by severine ovyn, 16 years ago

Fastjet added; CDFCones directory has been changed

File size: 12.4 KB
Line 
1//STARTHEADER
2// $Id: PseudoJet.cc,v 1.1 2008-11-06 14:32:15 ovyn Exp $
3//
4// Copyright (c) 2005-2006, Matteo Cacciari and Gavin Salam
5//
6//----------------------------------------------------------------------
7// This file is part of FastJet.
8//
9// FastJet is free software; you can redistribute it and/or modify
10// it under the terms of the GNU General Public License as published by
11// the Free Software Foundation; either version 2 of the License, or
12// (at your option) any later version.
13//
14// The algorithms that underlie FastJet have required considerable
15// development and are described in hep-ph/0512210. If you use
16// FastJet as part of work towards a scientific publication, please
17// include a citation to the FastJet paper.
18//
19// FastJet is distributed in the hope that it will be useful,
20// but WITHOUT ANY WARRANTY; without even the implied warranty of
21// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
22// GNU General Public License for more details.
23//
24// You should have received a copy of the GNU General Public License
25// along with FastJet; if not, write to the Free Software
26// Foundation, Inc.:
27// 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
28//----------------------------------------------------------------------
29//ENDHEADER
30
31
32#include "../include/fastjet/Error.hh"
33#include "../include/fastjet/PseudoJet.hh"
34#include<valarray>
35#include<iostream>
36#include<sstream>
37#include<cmath>
38#include<algorithm>
39
40FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
41
42using namespace std;
43
44
45//----------------------------------------------------------------------
46// another constructor...
47PseudoJet::PseudoJet(const double px, const double py, const double pz, const double E) {
48
49 _E = E ;
50 _px = px;
51 _py = py;
52 _pz = pz;
53
54 this->_finish_init();
55
56 // some default values for the history and user indices
57 _reset_indices();
58
59}
60
61
62//----------------------------------------------------------------------
63/// do standard end of initialisation
64void PseudoJet::_finish_init () {
65 _kt2 = this->px()*this->px() + this->py()*this->py();
66 if (_kt2 == 0.0) {
67 _phi = 0.0; }
68 else {
69 _phi = atan2(this->py(),this->px());
70 }
71 if (_phi < 0.0) {_phi += twopi;}
72 if (_phi >= twopi) {_phi -= twopi;} // can happen if phi=-|eps<1e-15|?
73 if (this->E() == abs(this->pz()) && _kt2 == 0) {
74 // Point has infinite rapidity -- convert that into a very large
75 // number, but in such a way that different 0-pt momenta will have
76 // different rapidities (so as to lift the degeneracy between
77 // them) [this can be relevant at parton-level]
78 double MaxRapHere = MaxRap + abs(this->pz());
79 if (this->pz() >= 0.0) {_rap = MaxRapHere;} else {_rap = -MaxRapHere;}
80 } else {
81 // get the rapidity in a way that's modestly insensitive to roundoff
82 // error when things pz,E are large (actually the best we can do without
83 // explicit knowledge of mass)
84 double effective_m2 = max(0.0,m2()); // force non tachyonic mass
85 double E_plus_pz = _E + abs(_pz); // the safer of p+, p-
86 // p+/p- = (p+ p-) / (p-)^2 = (kt^2+m^2)/(p-)^2
87 _rap = 0.5*log((_kt2 + effective_m2)/(E_plus_pz*E_plus_pz));
88 if (_pz > 0) {_rap = - _rap;}
89 }
90
91 //// original determination
92 //if (this->E() != abs(this->pz())) {
93 // _rap = 0.5*log((this->E() + this->pz())/(this->E() - this->pz()));
94 // } else {
95 // // Overlapping points can give problems. Let's lift the degeneracy
96 // // in case of multiple 0-pT points (can be found at parton-level)
97 // double MaxRapHere = MaxRap + abs(this->pz());
98 // if (this->pz() >= 0.0) {_rap = MaxRapHere;} else {_rap = -MaxRapHere;}
99 //}
100}
101
102
103//----------------------------------------------------------------------
104// return a valarray four-momentum
105valarray<double> PseudoJet::four_mom() const {
106 valarray<double> mom(4);
107 mom[0] = _px;
108 mom[1] = _py;
109 mom[2] = _pz;
110 mom[3] = _E ;
111 return mom;
112}
113
114//----------------------------------------------------------------------
115// Return the component corresponding to the specified index.
116// taken from CLHEP
117double PseudoJet::operator () (int i) const {
118 switch(i) {
119 case X:
120 return px();
121 case Y:
122 return py();
123 case Z:
124 return pz();
125 case T:
126 return e();
127 default:
128 ostringstream err;
129 err << "PseudoJet subscripting: bad index (" << i << ")";
130 throw Error(err.str());
131 }
132 return 0.;
133}
134
135//----------------------------------------------------------------------
136// return the pseudorapidity
137double PseudoJet::pseudorapidity() const {
138 if (px() == 0.0 && py() ==0.0) return MaxRap;
139 if (pz() == 0.0) return 0.0;
140
141 double theta = atan(perp()/pz());
142 if (theta < 0) theta += pi;
143 return -log(tan(theta/2));
144}
145
146//----------------------------------------------------------------------
147// return "sum" of two pseudojets
148PseudoJet operator+ (const PseudoJet & jet1, const PseudoJet & jet2) {
149 //return PseudoJet(jet1.four_mom()+jet2.four_mom());
150 return PseudoJet(jet1.px()+jet2.px(),
151 jet1.py()+jet2.py(),
152 jet1.pz()+jet2.pz(),
153 jet1.E() +jet2.E() );
154}
155
156//----------------------------------------------------------------------
157// return difference of two pseudojets
158PseudoJet operator- (const PseudoJet & jet1, const PseudoJet & jet2) {
159 //return PseudoJet(jet1.four_mom()-jet2.four_mom());
160 return PseudoJet(jet1.px()-jet2.px(),
161 jet1.py()-jet2.py(),
162 jet1.pz()-jet2.pz(),
163 jet1.E() -jet2.E() );
164}
165
166//----------------------------------------------------------------------
167// return the product, coeff * jet
168PseudoJet operator* (double coeff, const PseudoJet & jet) {
169 //return PseudoJet(coeff*jet.four_mom());
170 // the following code is hopefully more efficient
171 PseudoJet coeff_times_jet(jet);
172 coeff_times_jet *= coeff;
173 return coeff_times_jet;
174}
175
176//----------------------------------------------------------------------
177// return the product, coeff * jet
178PseudoJet operator* (const PseudoJet & jet, double coeff) {
179 return coeff*jet;
180}
181
182//----------------------------------------------------------------------
183// return the ratio, jet / coeff
184PseudoJet operator/ (const PseudoJet & jet, double coeff) {
185 return (1.0/coeff)*jet;
186}
187
188//----------------------------------------------------------------------
189/// multiply the jet's momentum by the coefficient
190void PseudoJet::operator*=(double coeff) {
191 _px *= coeff;
192 _py *= coeff;
193 _pz *= coeff;
194 _E *= coeff;
195 _kt2*= coeff*coeff;
196 // phi and rap are unchanged
197}
198
199//----------------------------------------------------------------------
200/// divide the jet's momentum by the coefficient
201void PseudoJet::operator/=(double coeff) {
202 (*this) *= 1.0/coeff;
203}
204
205
206//----------------------------------------------------------------------
207/// add the other jet's momentum to this jet
208void PseudoJet::operator+=(const PseudoJet & other_jet) {
209 _px += other_jet._px;
210 _py += other_jet._py;
211 _pz += other_jet._pz;
212 _E += other_jet._E ;
213 _finish_init(); // we need to recalculate phi,rap,kt2
214}
215
216
217//----------------------------------------------------------------------
218/// subtract the other jet's momentum from this jet
219void PseudoJet::operator-=(const PseudoJet & other_jet) {
220 _px -= other_jet._px;
221 _py -= other_jet._py;
222 _pz -= other_jet._pz;
223 _E -= other_jet._E ;
224 _finish_init(); // we need to recalculate phi,rap,kt2
225}
226
227
228//----------------------------------------------------------------------
229/// transform this jet (given in the rest frame of prest) into a jet
230/// in the lab frame;
231//
232// NB: code adapted from that in herwig f77 (checked how it worked
233// long ago)
234PseudoJet & PseudoJet::boost(const PseudoJet & prest) {
235
236 if (prest.px() == 0.0 && prest.py() == 0.0 && prest.pz() == 0.0)
237 return *this;
238
239 double m = prest.m();
240 assert(m != 0);
241
242 double pf4 = ( px()*prest.px() + py()*prest.py()
243 + pz()*prest.pz() + E()*prest.E() )/m;
244 double fn = (pf4 + E()) / (prest.E() + m);
245 _px += fn*prest.px();
246 _py += fn*prest.py();
247 _pz += fn*prest.pz();
248 _E = pf4;
249
250 _finish_init(); // we need to recalculate phi,rap,kt2
251 return *this;
252}
253
254
255//----------------------------------------------------------------------
256/// transform this jet (given in the rest frame of prest) into a jet
257/// in the lab frame;
258//
259// NB: code adapted from that in herwig f77 (checked how it worked
260// long ago)
261PseudoJet & PseudoJet::unboost(const PseudoJet & prest) {
262
263 if (prest.px() == 0.0 && prest.py() == 0.0 && prest.pz() == 0.0)
264 return *this;
265
266 double m = prest.m();
267 assert(m != 0);
268
269 double pf4 = ( -px()*prest.px() - py()*prest.py()
270 - pz()*prest.pz() + E()*prest.E() )/m;
271 double fn = (pf4 + E()) / (prest.E() + m);
272 _px -= fn*prest.px();
273 _py -= fn*prest.py();
274 _pz -= fn*prest.pz();
275 _E = pf4;
276
277 _finish_init(); // we need to recalculate phi,rap,kt2
278 return *this;
279}
280
281
282//----------------------------------------------------------------------
283/// returns true if the momenta of the two input jets are identical
284bool have_same_momentum(const PseudoJet & jeta, const PseudoJet & jetb) {
285 return jeta.px() == jetb.px()
286 && jeta.py() == jetb.py()
287 && jeta.pz() == jetb.pz()
288 && jeta.E() == jetb.E();
289}
290
291
292//----------------------------------------------------------------------
293/// return a pseudojet with the given pt, y, phi and mass
294PseudoJet PtYPhiM(double pt, double y, double phi, double m) {
295 double ptm = sqrt(pt*pt+m*m);
296 return PseudoJet(pt*cos(phi), pt*sin(phi), ptm*sinh(y), ptm*cosh(y));
297}
298
299
300//----------------------------------------------------------------------
301// return kt-distance between this jet and another one
302double PseudoJet::kt_distance(const PseudoJet & other) const {
303 //double distance = min(this->kt2(), other.kt2());
304 double distance = min(_kt2, other._kt2);
305 double dphi = abs(_phi - other._phi);
306 if (dphi > pi) {dphi = twopi - dphi;}
307 double drap = _rap - other._rap;
308 distance = distance * (dphi*dphi + drap*drap);
309 return distance;
310}
311
312
313//----------------------------------------------------------------------
314// return squared cylinder (eta-phi) distance between this jet and another one
315double PseudoJet::plain_distance(const PseudoJet & other) const {
316 double dphi = abs(_phi - other._phi);
317 if (dphi > pi) {dphi = twopi - dphi;}
318 double drap = _rap - other._rap;
319 return (dphi*dphi + drap*drap);
320}
321
322//----------------------------------------------------------------------
323// sort the indices so that values[indices[0..n-1]] is sorted
324// into increasing order
325void sort_indices(vector<int> & indices,
326 const vector<double> & values) {
327 IndexedSortHelper index_sort_helper(&values);
328 sort(indices.begin(), indices.end(), index_sort_helper);
329}
330
331//----------------------------------------------------------------------
332/// given a vector of values with a one-to-one correspondence with the
333/// vector of objects, sort objects into an order such that the
334/// associated values would be in increasing order
335template<class T> vector<T> objects_sorted_by_values(
336 const vector<T> & objects,
337 const vector<double> & values) {
338
339 assert(objects.size() == values.size());
340
341 // get a vector of indices
342 vector<int> indices(values.size());
343 for (size_t i = 0; i < indices.size(); i++) {indices[i] = i;}
344
345 // sort the indices
346 sort_indices(indices, values);
347
348 // copy the objects
349 vector<T> objects_sorted(objects.size());
350
351 // place the objects in the correct order
352 for (size_t i = 0; i < indices.size(); i++) {
353 objects_sorted[i] = objects[indices[i]];
354 }
355
356 return objects_sorted;
357}
358
359//----------------------------------------------------------------------
360/// return a vector of jets sorted into decreasing kt2
361vector<PseudoJet> sorted_by_pt(const vector<PseudoJet> & jets) {
362 vector<double> minus_kt2(jets.size());
363 for (size_t i = 0; i < jets.size(); i++) {minus_kt2[i] = -jets[i].kt2();}
364 return objects_sorted_by_values(jets, minus_kt2);
365}
366
367//----------------------------------------------------------------------
368/// return a vector of jets sorted into increasing rapidity
369vector<PseudoJet> sorted_by_rapidity(const vector<PseudoJet> & jets) {
370 vector<double> rapidities(jets.size());
371 for (size_t i = 0; i < jets.size(); i++) {rapidities[i] = jets[i].rap();}
372 return objects_sorted_by_values(jets, rapidities);
373}
374
375//----------------------------------------------------------------------
376/// return a vector of jets sorted into decreasing energy
377vector<PseudoJet> sorted_by_E(const vector<PseudoJet> & jets) {
378 vector<double> energies(jets.size());
379 for (size_t i = 0; i < jets.size(); i++) {energies[i] = -jets[i].E();}
380 return objects_sorted_by_values(jets, energies);
381}
382
383
384FASTJET_END_NAMESPACE
385
Note: See TracBrowser for help on using the repository browser.