FastJet 3.4.1
PseudoJet.cc
1//FJSTARTHEADER
2// $Id$
3//
4// Copyright (c) 2005-2023, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
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. They are described in the original FastJet paper,
16// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
17// FastJet as part of work towards a scientific publication, please
18// quote the version you use and include a citation to the manual and
19// optionally also to hep-ph/0512210.
20//
21// FastJet is distributed in the hope that it will be useful,
22// but WITHOUT ANY WARRANTY; without even the implied warranty of
23// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24// GNU General Public License for more details.
25//
26// You should have received a copy of the GNU General Public License
27// along with FastJet. If not, see <http://www.gnu.org/licenses/>.
28//----------------------------------------------------------------------
29//FJENDHEADER
30
31
32#include "fastjet/Error.hh"
33#include "fastjet/PseudoJet.hh"
34#include "fastjet/ClusterSequence.hh"
35#ifndef __FJCORE__
36#include "fastjet/ClusterSequenceAreaBase.hh"
37#endif // __FJCORE__
38#include "fastjet/CompositeJetStructure.hh"
39#include<valarray>
40#include<iostream>
41#include<sstream>
42#include<cmath>
43#include<algorithm>
44#include <cstdarg>
45
46FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
47
48using namespace std;
49
50
51//----------------------------------------------------------------------
52// another constructor...
53PseudoJet::PseudoJet(const double px_in, const double py_in, const double pz_in, const double E_in) {
54
55 _E = E_in ;
56 _px = px_in;
57 _py = py_in;
58 _pz = pz_in;
59
60 this->_finish_init();
61
62 // some default values for the history and user indices
63 // note: no reset of shared pointers needed
64 _reset_indices();
65
66}
67
68#ifdef FASTJET_HAVE_THREAD_SAFETY
69/// copy-assignmemt
70///
71/// this has to be explicitly specified since atomic does not support it.
72PseudoJet & PseudoJet::operator=(const PseudoJet & other_pj){
73 _structure = other_pj._structure;
74 _user_info = other_pj._user_info;
75
76 _cluster_hist_index = other_pj._cluster_hist_index;
77 _user_index = other_pj._user_index;
78
79 // copy the remaining information through the reset_momentum call
80 // which properly handles the thread safety, notably related
81 // to caching of the phi and rap values
82 reset_momentum(other_pj);
83
84 return *this;
85}
86#endif // FASTJET_HAVE_THREAD_SAFETY
87
88//std::shared_ptr: //----------------------------------------------------------------------
89//std::shared_ptr: #ifdef FASTJET_HAVE_THREAD_SAFETY
90//std::shared_ptr: PseudoJet::~PseudoJet(){
91//std::shared_ptr: _release_jet_from_cs();
92//std::shared_ptr: }
93//std::shared_ptr:
94//std::shared_ptr: // this has to be called everytime one tries to alter the jet
95//std::shared_ptr: // structural info
96//std::shared_ptr: void PseudoJet::_release_jet_from_cs(){
97//std::shared_ptr: // check if the jet has the structure of type CSstruct in which case
98//std::shared_ptr: // we have to check if there is a need for self-deletion of the CS
99//std::shared_ptr: ClusterSequenceStructure * assoc_css = dynamic_cast<ClusterSequenceStructure *>(_structure.get());
100//std::shared_ptr: if (assoc_css) {
101//std::shared_ptr: const ClusterSequence * assoc_cs = assoc_css->associated_cluster_sequence();
102//std::shared_ptr: if (assoc_cs) assoc_cs->release_pseudojet(*this);
103//std::shared_ptr: }
104//std::shared_ptr:
105//std::shared_ptr: // slightly less efficient version
106//std::shared_ptr: // if ((has_structure_of<ClusterSequence>()) && (has_valid_cluster_sequence())){
107//std::shared_ptr: // associated_cs()->release_pseudojet(*this);
108//std::shared_ptr: // }
109//std::shared_ptr: }
110//std::shared_ptr:
111//std::shared_ptr: #endif // FASTJET_HAVE_THREAD_SAFETY
112
113//----------------------------------------------------------------------
114/// do standard end of initialisation
115void PseudoJet::_finish_init () {
116 _kt2 = this->px()*this->px() + this->py()*this->py();
117 _phi = pseudojet_invalid_phi;
118 // strictly speaking, _rap does not need initialising, because
119 // it's never used as long as _phi == pseudojet_invalid_phi
120 // (and gets set when _phi is requested). However ATLAS
121 // 2013-03-28 complained that they sometimes have NaN's in
122 // _rap and this interferes with some of their internal validation.
123 // So we initialise it; penalty is about 0.3ns per PseudoJet out of about
124 // 10ns total initialisation time (on a intel Core i7 2.7GHz)
125 _rap = pseudojet_invalid_rap;
126
127#ifdef FASTJET_HAVE_THREAD_SAFETY
128 _init_status = Init_NotDone;
129#endif
130}
131
132//----------------------------------------------------------------------
133#ifdef FASTJET_HAVE_THREAD_SAFETY
134void PseudoJet::_ensure_valid_rap_phi() const{
135 //TODO: for _init_status we can use memory_order_release when
136 //writing and memory_order_acquire when reading. Check if that has
137 //an impact on timing
138 if (_init_status!=Init_Done){
139 int expected = Init_NotDone;
140 // if it's 0 (not done), set thing to "operation in progress"
141 // (-1) and do the init
142 //
143 // Note:
144 // - this has to be the strong version so we can make a single
145 // test
146 // - we need to make sure that the -1 is loaded properly
147 // (success), hence a std::memory_order_seq_cst model
148 // - the failure model can be anything since we're not going
149 // to use the value anyway
150 if (_init_status.compare_exchange_strong(expected, Init_InProgress,
151 std::memory_order_seq_cst,
152 std::memory_order_relaxed)){
153 _set_rap_phi();
154 // Now we can safely set the status flag to Init_Done (after all
155 // physics variables are set). The memory ordering guarantees
156 // synchronisation across threads.
157 _init_status.store(Init_Done, memory_order_release);
158 } else {
159 // wait until done
160 //
161 // Here we use memory_order_acquire to guarantee that the
162 // rap/phi values set in the thread which actually does the
163 // calculation are synchronised.
164 // Reading https://gcc.gnu.org/wiki/Atomic/GCCMM/AtomicSync in
165 // interesting in this context.
166 while (_init_status.load(memory_order_acquire) != Init_Done);
167 }
168
169 }
170}
171#endif
172
173//----------------------------------------------------------------------
174void PseudoJet::_set_rap_phi() const {
175
176 if (_kt2 == 0.0) {
177 _phi = 0.0; }
178 else {
179 _phi = atan2(this->py(),this->px());
180 }
181 if (_phi < 0.0) {_phi += twopi;}
182 if (_phi >= twopi) {_phi -= twopi;} // can happen if phi=-|eps<1e-15|?
183 if (this->E() == abs(this->pz()) && _kt2 == 0) {
184 // Point has infinite rapidity -- convert that into a very large
185 // number, but in such a way that different 0-pt momenta will have
186 // different rapidities (so as to lift the degeneracy between
187 // them) [this can be relevant at parton-level]
188 double MaxRapHere = MaxRap + abs(this->pz());
189 if (this->pz() >= 0.0) {_rap = MaxRapHere;} else {_rap = -MaxRapHere;}
190 } else {
191 // get the rapidity in a way that's modestly insensitive to roundoff
192 // error when things pz,E are large (actually the best we can do without
193 // explicit knowledge of mass)
194 double effective_m2 = max(0.0,m2()); // force non tachyonic mass
195 double E_plus_pz = _E + abs(_pz); // the safer of p+, p-
196 // p+/p- = (p+ p-) / (p-)^2 = (kt^2+m^2)/(p-)^2
197 _rap = 0.5*log((_kt2 + effective_m2)/(E_plus_pz*E_plus_pz));
198 if (_pz > 0) {_rap = - _rap;}
199 }
200
201}
202
203
204//----------------------------------------------------------------------
205// return a valarray four-momentum
206valarray<double> PseudoJet::four_mom() const {
207 valarray<double> mom(4);
208 mom[0] = _px;
209 mom[1] = _py;
210 mom[2] = _pz;
211 mom[3] = _E ;
212 return mom;
213}
214
215//----------------------------------------------------------------------
216// Return the component corresponding to the specified index.
217// taken from CLHEP
218double PseudoJet::operator () (int i) const {
219 switch(i) {
220 case X:
221 return px();
222 case Y:
223 return py();
224 case Z:
225 return pz();
226 case T:
227 return e();
228 default:
229 ostringstream err;
230 err << "PseudoJet subscripting: bad index (" << i << ")";
231 throw Error(err.str());
232 }
233 return 0.;
234}
235
236//----------------------------------------------------------------------
237// return the pseudorapidity
238double PseudoJet::pseudorapidity() const {
239 if (px() == 0.0 && py() ==0.0) return MaxRap;
240 if (pz() == 0.0) return 0.0;
241
242 double theta = atan(perp()/pz());
243 if (theta < 0) theta += pi;
244 return -log(tan(theta/2));
245}
246
247//----------------------------------------------------------------------
248// return "sum" of two pseudojets
249PseudoJet operator+ (const PseudoJet & jet1, const PseudoJet & jet2) {
250 //return PseudoJet(jet1.four_mom()+jet2.four_mom());
251 return PseudoJet(jet1.px()+jet2.px(),
252 jet1.py()+jet2.py(),
253 jet1.pz()+jet2.pz(),
254 jet1.E() +jet2.E() );
255}
256
257//----------------------------------------------------------------------
258// return difference of two pseudojets
259PseudoJet operator- (const PseudoJet & jet1, const PseudoJet & jet2) {
260 //return PseudoJet(jet1.four_mom()-jet2.four_mom());
261 return PseudoJet(jet1.px()-jet2.px(),
262 jet1.py()-jet2.py(),
263 jet1.pz()-jet2.pz(),
264 jet1.E() -jet2.E() );
265}
266
267//----------------------------------------------------------------------
268// return the product, coeff * jet
269PseudoJet operator* (double coeff, const PseudoJet & jet) {
270 // see the comment in operator*= about ensuring valid rap phi
271 // before a multiplication to handle case of multiplication by
272 // zero, while maintaining rapidity and phi
273 jet._ensure_valid_rap_phi();
274 //return PseudoJet(coeff*jet.four_mom());
275 // the following code is hopefully more efficient
276 PseudoJet coeff_times_jet(jet);
277 coeff_times_jet *= coeff;
278 return coeff_times_jet;
279}
280
281//----------------------------------------------------------------------
282// return the product, coeff * jet
283PseudoJet operator* (const PseudoJet & jet, double coeff) {
284 return coeff*jet;
285}
286
287//----------------------------------------------------------------------
288// return the ratio, jet / coeff
289PseudoJet operator/ (const PseudoJet & jet, double coeff) {
290 return (1.0/coeff)*jet;
291}
292
293//----------------------------------------------------------------------
294/// multiply the jet's momentum by the coefficient
295PseudoJet & PseudoJet::operator*=(double coeff) {
296 // operator*= aims to maintain the rapidity and azimuth
297 // for the PseudoJet; if they have already been evaluated
298 // this is fine, but if they haven't and coeff is sufficiently
299 // small as to cause a zero or underflow result, then a subsequent
300 // invocation of rap or phi will lead to a non-sensical result.
301 // So, here, we preemptively ensure that rapidity and phi
302 // are correctly cached
303 _ensure_valid_rap_phi();
304 _px *= coeff;
305 _py *= coeff;
306 _pz *= coeff;
307 _E *= coeff;
308 _kt2*= coeff*coeff;
309 // phi and rap are unchanged
310 return *this;
311}
312
313//----------------------------------------------------------------------
314/// divide the jet's momentum by the coefficient
315PseudoJet & PseudoJet::operator/=(double coeff) {
316 (*this) *= 1.0/coeff;
317 return *this;
318}
319
320
321//----------------------------------------------------------------------
322/// add the other jet's momentum to this jet
323PseudoJet & PseudoJet::operator+=(const PseudoJet & other_jet) {
324 _px += other_jet._px;
325 _py += other_jet._py;
326 _pz += other_jet._pz;
327 _E += other_jet._E ;
328 _finish_init(); // we need to recalculate phi,rap,kt2
329 return *this;
330}
331
332
333//----------------------------------------------------------------------
334/// subtract the other jet's momentum from this jet
335PseudoJet & PseudoJet::operator-=(const PseudoJet & other_jet) {
336 _px -= other_jet._px;
337 _py -= other_jet._py;
338 _pz -= other_jet._pz;
339 _E -= other_jet._E ;
340 _finish_init(); // we need to recalculate phi,rap,kt2
341 return *this;
342}
343
344//----------------------------------------------------------------------
345bool operator==(const PseudoJet & a, const PseudoJet & b) {
346 if (a.px() != b.px()) return false;
347 if (a.py() != b.py()) return false;
348 if (a.pz() != b.pz()) return false;
349 if (a.E () != b.E ()) return false;
350
351 if (a.user_index() != b.user_index()) return false;
352 if (a.cluster_hist_index() != b.cluster_hist_index()) return false;
353 if (a.user_info_ptr() != b.user_info_ptr()) return false;
354 if (a.structure_ptr() != b.structure_ptr()) return false;
355
356 return true;
357}
358
359//----------------------------------------------------------------------
360// check if the jet has zero momentum
361bool operator==(const PseudoJet & jet, const double val) {
362 if (val != 0)
363 throw Error("comparing a PseudoJet with a non-zero constant (double) is not allowed.");
364 return (jet.px() == 0 && jet.py() == 0 &&
365 jet.pz() == 0 && jet.E() == 0);
366}
367
368
369
370//----------------------------------------------------------------------
371/// transform this jet (given in the rest frame of prest) into a jet
372/// in the lab frame
373//
374// NB: code adapted from that in herwig f77 (checked how it worked
375// long ago)
376PseudoJet & PseudoJet::boost(const PseudoJet & prest) {
377
378 if (prest.px() == 0.0 && prest.py() == 0.0 && prest.pz() == 0.0)
379 return *this;
380
381 double m_local = prest.m();
382 assert(m_local != 0);
383
384 double pf4 = ( px()*prest.px() + py()*prest.py()
385 + pz()*prest.pz() + E()*prest.E() )/m_local;
386 double fn = (pf4 + E()) / (prest.E() + m_local);
387 _px += fn*prest.px();
388 _py += fn*prest.py();
389 _pz += fn*prest.pz();
390 _E = pf4;
391
392 _finish_init(); // we need to recalculate phi,rap,kt2
393 return *this;
394}
395
396
397//----------------------------------------------------------------------
398/// transform this jet (given in lab) into a jet in the rest
399/// frame of prest
400//
401// NB: code adapted from that in herwig f77 (checked how it worked
402// long ago)
403PseudoJet & PseudoJet::unboost(const PseudoJet & prest) {
404
405 if (prest.px() == 0.0 && prest.py() == 0.0 && prest.pz() == 0.0)
406 return *this;
407
408 double m_local = prest.m();
409 assert(m_local != 0);
410
411 double pf4 = ( -px()*prest.px() - py()*prest.py()
412 - pz()*prest.pz() + E()*prest.E() )/m_local;
413 double fn = (pf4 + E()) / (prest.E() + m_local);
414 _px -= fn*prest.px();
415 _py -= fn*prest.py();
416 _pz -= fn*prest.pz();
417 _E = pf4;
418
419 _finish_init(); // we need to recalculate phi,rap,kt2
420 return *this;
421}
422
423
424//----------------------------------------------------------------------
425/// returns true if the momenta of the two input jets are identical
426bool have_same_momentum(const PseudoJet & jeta, const PseudoJet & jetb) {
427 return jeta.px() == jetb.px()
428 && jeta.py() == jetb.py()
429 && jeta.pz() == jetb.pz()
430 && jeta.E() == jetb.E();
431}
432
433//----------------------------------------------------------------------
434void PseudoJet::set_cached_rap_phi(double rap_in, double phi_in) {
435 _rap = rap_in; _phi = phi_in;
436 if (_phi >= twopi) _phi -= twopi;
437 if (_phi < 0) _phi += twopi;
438#ifdef FASTJET_HAVE_THREAD_SAFETY
439 _init_status = Init_Done;
440#endif
441}
442
443//----------------------------------------------------------------------
444void PseudoJet::reset_momentum_PtYPhiM(double pt_in, double y_in, double phi_in, double m_in) {
445 assert(phi_in < 2*twopi && phi_in > -twopi);
446 double ptm = (m_in == 0) ? pt_in : sqrt(pt_in*pt_in+m_in*m_in);
447 double exprap = exp(y_in);
448 double pminus = ptm/exprap;
449 double pplus = ptm*exprap;
450 double px_local = pt_in*cos(phi_in);
451 double py_local = pt_in*sin(phi_in);
452 reset_momentum(px_local,py_local,0.5*(pplus-pminus),0.5*(pplus+pminus));
453 set_cached_rap_phi(y_in,phi_in);
454}
455
456//----------------------------------------------------------------------
457/// return a pseudojet with the given pt, y, phi and mass
458PseudoJet PtYPhiM(double pt, double y, double phi, double m) {
459 assert(phi < 2*twopi && phi > -twopi);
460 double ptm = (m == 0) ? pt : sqrt(pt*pt+m*m);
461 double exprap = exp(y);
462 double pminus = ptm/exprap;
463 double pplus = ptm*exprap;
464 double px = pt*cos(phi);
465 double py = pt*sin(phi);
466 PseudoJet mom(px,py,0.5*(pplus-pminus),0.5*(pplus+pminus));
467 mom.set_cached_rap_phi(y,phi);
468 return mom;
469 //return PseudoJet(pt*cos(phi), pt*sin(phi), ptm*sinh(y), ptm*cosh(y));
470}
471
472
473//----------------------------------------------------------------------
474// return kt-distance between this jet and another one
475double PseudoJet::kt_distance(const PseudoJet & other) const {
476 //double distance = min(this->kt2(), other.kt2());
477 double distance = min(_kt2, other._kt2);
478 double dphi = abs(phi() - other.phi());
479 if (dphi > pi) {dphi = twopi - dphi;}
480 double drap = rap() - other.rap();
481 distance = distance * (dphi*dphi + drap*drap);
482 return distance;
483}
484
485
486//----------------------------------------------------------------------
487// return squared cylinder (eta-phi) distance between this jet and another one
488double PseudoJet::plain_distance(const PseudoJet & other) const {
489 double dphi = abs(phi() - other.phi());
490 if (dphi > pi) {dphi = twopi - dphi;}
491 double drap = rap() - other.rap();
492 return (dphi*dphi + drap*drap);
493}
494
495//----------------------------------------------------------------------
496/// returns other.phi() - this.phi(), i.e. the phi distance to
497/// other, constrained to be in range -pi .. pi
498double PseudoJet::delta_phi_to(const PseudoJet & other) const {
499 double dphi = other.phi() - phi();
500 if (dphi > pi) dphi -= twopi;
501 if (dphi < -pi) dphi += twopi;
502 return dphi;
503}
504
505
506string PseudoJet::description() const{
507 // the "default" case of a PJ which does not belong to any cluster sequence
508 if (!_structure)
509 return "standard PseudoJet (with no associated clustering information)";
510
511 // for all the other cases, the description comes from the structure
512 return _structure->description();
513}
514
515
516
517//----------------------------------------------------------------------
518//
519// The following methods access the associated jet structure (if any)
520//
521//----------------------------------------------------------------------
522
523
524//----------------------------------------------------------------------
525// check whether this PseudoJet has an associated parent
526// ClusterSequence
527bool PseudoJet::has_associated_cluster_sequence() const{
528 return (_structure) && (_structure->has_associated_cluster_sequence());
529}
530
531//----------------------------------------------------------------------
532// get a (const) pointer to the associated ClusterSequence (NULL if
533// inexistent)
534const ClusterSequence* PseudoJet::associated_cluster_sequence() const{
535 if (! has_associated_cluster_sequence()) return NULL;
536
537 return _structure->associated_cluster_sequence();
538}
539
540
541//----------------------------------------------------------------------
542// check whether this PseudoJet has an associated parent
543// ClusterSequence that is still valid
544bool PseudoJet::has_valid_cluster_sequence() const{
545 return (_structure) && (_structure->has_valid_cluster_sequence());
546}
547
548//----------------------------------------------------------------------
549// If there is a valid cluster sequence associated with this jet,
550// returns a pointer to it; otherwise throws an Error.
551//
552// Open question: should these errors be upgraded to classes of their
553// own so that they can be caught? [Maybe, but later]
554const ClusterSequence * PseudoJet::validated_cs() const {
555 return validated_structure_ptr()->validated_cs();
556}
557
558
559//----------------------------------------------------------------------
560// set the associated structure
561void PseudoJet::set_structure_shared_ptr(const SharedPtr<PseudoJetStructureBase> &structure_in){
562//std::shared_ptr: #ifdef FASTJET_HAVE_THREAD_SAFETY
563//std::shared_ptr: // if the jet currently belongs to a cs, we need to release it before any chenge
564//std::shared_ptr: _release_jet_from_cs();
565//std::shared_ptr: #endif //FASTJET_HAVE_THREAD_SAFETY
566 _structure = structure_in;
567}
568
569//----------------------------------------------------------------------
570// return true if there is some structure associated with this PseudoJet
571bool PseudoJet::has_structure() const{
572 return (bool) _structure; // cast to bool has to be made explicit
573}
574
575//----------------------------------------------------------------------
576// return a pointer to the structure (of type
577// PseudoJetStructureBase*) associated with this PseudoJet.
578//
579// return NULL if there is no associated structure
580const PseudoJetStructureBase* PseudoJet::structure_ptr() const {
581 //if (!_structure) return NULL;
582 return _structure.get();
583}
584
585//----------------------------------------------------------------------
586// return a non-const pointer to the structure (of type
587// PseudoJetStructureBase*) associated with this PseudoJet.
588//
589// return NULL if there is no associated structure
590//
591// Only use this if you know what you are doing. In any case,
592// prefer the 'structure_ptr()' (the const version) to this method,
593// unless you really need a write access to the PseudoJet's
594// underlying structure.
595PseudoJetStructureBase* PseudoJet::structure_non_const_ptr(){
596 //if (!_structure) return NULL;
597 return _structure.get();
598}
599
600//----------------------------------------------------------------------
601// return a pointer to the structure (of type
602// PseudoJetStructureBase*) associated with this PseudoJet.
603//
604// throw an error if there is no associated structure
605const PseudoJetStructureBase* PseudoJet::validated_structure_ptr() const {
606 if (!_structure)
607 throw Error("Trying to access the structure of a PseudoJet which has no associated structure");
608 return _structure.get();
609}
610
611//----------------------------------------------------------------------
612// return a reference to the shared pointer to the
613// PseudoJetStructureBase associated with this PseudoJet
614const SharedPtr<PseudoJetStructureBase> & PseudoJet::structure_shared_ptr() const {
615 return _structure;
616}
617
618
619//----------------------------------------------------------------------
620// check if it has been recombined with another PseudoJet in which
621// case, return its partner through the argument. Otherwise,
622// 'partner' is set to 0.
623//
624// false is also returned if this PseudoJet has no associated
625// ClusterSequence
626bool PseudoJet::has_partner(PseudoJet &partner) const{
627 return validated_structure_ptr()->has_partner(*this, partner);
628}
629
630//----------------------------------------------------------------------
631// check if it has been recombined with another PseudoJet in which
632// case, return its child through the argument. Otherwise, 'child'
633// is set to 0.
634//
635// false is also returned if this PseudoJet has no associated
636// ClusterSequence, with the child set to 0
637bool PseudoJet::has_child(PseudoJet &child) const{
638 return validated_structure_ptr()->has_child(*this, child);
639}
640
641//----------------------------------------------------------------------
642// check if it is the product of a recombination, in which case
643// return the 2 parents through the 'parent1' and 'parent2'
644// arguments. Otherwise, set these to 0.
645//
646// false is also returned if this PseudoJet has no parent
647// ClusterSequence
648bool PseudoJet::has_parents(PseudoJet &parent1, PseudoJet &parent2) const{
649 return validated_structure_ptr()->has_parents(*this, parent1, parent2);
650}
651
652//----------------------------------------------------------------------
653// check if the current PseudoJet contains the one passed as
654// argument
655//
656// false is also returned if this PseudoJet has no associated
657// ClusterSequence.
658bool PseudoJet::contains(const PseudoJet &constituent) const{
659 return validated_structure_ptr()->object_in_jet(constituent, *this);
660}
661
662//----------------------------------------------------------------------
663// check if the current PseudoJet is contained the one passed as
664// argument
665//
666// false is also returned if this PseudoJet has no associated
667// ClusterSequence
668bool PseudoJet::is_inside(const PseudoJet &jet) const{
669 return validated_structure_ptr()->object_in_jet(*this, jet);
670}
671
672
673//----------------------------------------------------------------------
674// returns true if the PseudoJet has constituents
675bool PseudoJet::has_constituents() const{
676 return (_structure) && (_structure->has_constituents());
677}
678
679//----------------------------------------------------------------------
680// retrieve the constituents.
681vector<PseudoJet> PseudoJet::constituents() const{
682 return validated_structure_ptr()->constituents(*this);
683}
684
685
686//----------------------------------------------------------------------
687// returns true if the PseudoJet has support for exclusive subjets
688bool PseudoJet::has_exclusive_subjets() const{
689 return (_structure) && (_structure->has_exclusive_subjets());
690}
691
692//----------------------------------------------------------------------
693// return a vector of all subjets of the current jet (in the sense
694// of the exclusive algorithm) that would be obtained when running
695// the algorithm with the given dcut.
696//
697// Time taken is O(m ln m), where m is the number of subjets that
698// are found. If m gets to be of order of the total number of
699// constituents in the jet, this could be substantially slower than
700// just getting that list of constituents.
701//
702// an Error is thrown if this PseudoJet has no currently valid
703// associated ClusterSequence
704std::vector<PseudoJet> PseudoJet::exclusive_subjets (const double dcut) const {
705 return validated_structure_ptr()->exclusive_subjets(*this, dcut);
706}
707
708//----------------------------------------------------------------------
709// return the size of exclusive_subjets(...); still n ln n with same
710// coefficient, but marginally more efficient than manually taking
711// exclusive_subjets.size()
712//
713// an Error is thrown if this PseudoJet has no currently valid
714// associated ClusterSequence
715int PseudoJet::n_exclusive_subjets(const double dcut) const {
716 return validated_structure_ptr()->n_exclusive_subjets(*this, dcut);
717}
718
719//----------------------------------------------------------------------
720// return the list of subjets obtained by unclustering the supplied
721// jet down to n subjets (or all constituents if there are fewer
722// than n).
723//
724// requires n ln n time
725//
726// an Error is thrown if this PseudoJet has no currently valid
727// associated ClusterSequence
728std::vector<PseudoJet> PseudoJet::exclusive_subjets_up_to (int nsub) const {
729 return validated_structure_ptr()->exclusive_subjets_up_to(*this, nsub);
730}
731
732//----------------------------------------------------------------------
733// Same as exclusive_subjets_up_to but throws an error if there are
734// fewer than nsub particles in the jet
735std::vector<PseudoJet> PseudoJet::exclusive_subjets (int nsub) const {
736 vector<PseudoJet> subjets = exclusive_subjets_up_to(nsub);
737 if (int(subjets.size()) < nsub) {
738 ostringstream err;
739 err << "Requested " << nsub << " exclusive subjets, but there were only "
740 << subjets.size() << " particles in the jet";
741 throw Error(err.str());
742 }
743 return subjets;
744}
745
746//----------------------------------------------------------------------
747// return the dij that was present in the merging nsub+1 -> nsub
748// subjets inside this jet.
749//
750// an Error is thrown if this PseudoJet has no currently valid
751// associated ClusterSequence
752double PseudoJet::exclusive_subdmerge(int nsub) const {
753 return validated_structure_ptr()->exclusive_subdmerge(*this, nsub);
754}
755
756//----------------------------------------------------------------------
757// return the maximum dij that occurred in the whole event at the
758// stage that the nsub+1 -> nsub merge of subjets occurred inside
759// this jet.
760//
761// an Error is thrown if this PseudoJet has no currently valid
762// associated ClusterSequence
763double PseudoJet::exclusive_subdmerge_max(int nsub) const {
764 return validated_structure_ptr()->exclusive_subdmerge_max(*this, nsub);
765}
766
767
768// returns true if a jet has pieces
769//
770// By default a single particle or a jet coming from a
771// ClusterSequence have no pieces and this methos will return false.
772bool PseudoJet::has_pieces() const{
773 return ((_structure) && (_structure->has_pieces(*this)));
774}
775
776// retrieve the pieces that make up the jet.
777//
778// By default a jet does not have pieces.
779// If the underlying interface supports "pieces" retrieve the
780// pieces from there.
781std::vector<PseudoJet> PseudoJet::pieces() const{
782 return validated_structure_ptr()->pieces(*this);
783 // if (!has_pieces())
784 // throw Error("Trying to retrieve the pieces of a PseudoJet that has no support for pieces.");
785 //
786 // return _structure->pieces(*this);
787}
788
789
790//----------------------------------------------------------------------
791// the following ones require a computation of the area in the
792// associated ClusterSequence (See ClusterSequenceAreaBase for details)
793//----------------------------------------------------------------------
794
795#ifndef __FJCORE__
796
797//----------------------------------------------------------------------
798// if possible, return a valid ClusterSequenceAreaBase pointer; otherwise
799// throw an error
800const ClusterSequenceAreaBase * PseudoJet::validated_csab() const {
801 const ClusterSequenceAreaBase *csab = dynamic_cast<const ClusterSequenceAreaBase*>(validated_cs());
802 if (csab == NULL) throw Error("you requested jet-area related information, but the PseudoJet does not have associated area information.");
803 return csab;
804}
805
806
807//----------------------------------------------------------------------
808// check if it has a defined area
809bool PseudoJet::has_area() const{
810 //if (! has_associated_cluster_sequence()) return false;
811 if (! has_structure()) return false;
812 return (validated_structure_ptr()->has_area() != 0);
813}
814
815//----------------------------------------------------------------------
816// return the jet (scalar) area.
817// throw an Error if there is no support for area in the associated CS
818double PseudoJet::area() const{
819 return validated_structure_ptr()->area(*this);
820}
821
822//----------------------------------------------------------------------
823// return the error (uncertainty) associated with the determination
824// of the area of this jet.
825// throws an Error if there is no support for area in the associated CS
826double PseudoJet::area_error() const{
827 return validated_structure_ptr()->area_error(*this);
828}
829
830//----------------------------------------------------------------------
831// return the jet 4-vector area
832// throws an Error if there is no support for area in the associated CS
833PseudoJet PseudoJet::area_4vector() const{
834 return validated_structure_ptr()->area_4vector(*this);
835}
836
837//----------------------------------------------------------------------
838// true if this jet is made exclusively of ghosts
839// throws an Error if there is no support for area in the associated CS
840bool PseudoJet::is_pure_ghost() const{
841 return validated_structure_ptr()->is_pure_ghost(*this);
842}
843
844#endif // __FJCORE__
845
846//----------------------------------------------------------------------
847//
848// end of the methods accessing the information in the associated
849// Cluster Sequence
850//
851//----------------------------------------------------------------------
852
853//----------------------------------------------------------------------
854/// provide a meaningful error message for InexistentUserInfo
855PseudoJet::InexistentUserInfo::InexistentUserInfo() : Error("you attempted to perform a dynamic cast of a PseudoJet's extra info, but the extra info pointer was null")
856{}
857
858
859//----------------------------------------------------------------------
860// sort the indices so that values[indices[0..n-1]] is sorted
861// into increasing order
862void sort_indices(vector<int> & indices,
863 const vector<double> & values) {
864 IndexedSortHelper index_sort_helper(&values);
865 sort(indices.begin(), indices.end(), index_sort_helper);
866}
867
868
869//----------------------------------------------------------------------
870/// return a vector of jets sorted into decreasing kt2
871vector<PseudoJet> sorted_by_pt(const vector<PseudoJet> & jets) {
872 vector<double> minus_kt2(jets.size());
873 for (size_t i = 0; i < jets.size(); i++) {minus_kt2[i] = -jets[i].kt2();}
874 return objects_sorted_by_values(jets, minus_kt2);
875}
876
877//----------------------------------------------------------------------
878/// return a vector of jets sorted into increasing rapidity
879vector<PseudoJet> sorted_by_rapidity(const vector<PseudoJet> & jets) {
880 vector<double> rapidities(jets.size());
881 for (size_t i = 0; i < jets.size(); i++) {rapidities[i] = jets[i].rap();}
882 return objects_sorted_by_values(jets, rapidities);
883}
884
885//----------------------------------------------------------------------
886/// return a vector of jets sorted into decreasing energy
887vector<PseudoJet> sorted_by_E(const vector<PseudoJet> & jets) {
888 vector<double> energies(jets.size());
889 for (size_t i = 0; i < jets.size(); i++) {energies[i] = -jets[i].E();}
890 return objects_sorted_by_values(jets, energies);
891}
892
893//----------------------------------------------------------------------
894/// return a vector of jets sorted into increasing pz
895vector<PseudoJet> sorted_by_pz(const vector<PseudoJet> & jets) {
896 vector<double> pz(jets.size());
897 for (size_t i = 0; i < jets.size(); i++) {pz[i] = jets[i].pz();}
898 return objects_sorted_by_values(jets, pz);
899}
900
901
902
903//-------------------------------------------------------------------------------
904// helper functions to build a jet made of pieces
905//-------------------------------------------------------------------------------
906
907// build a "CompositeJet" from the vector of its pieces
908//
909// In this case, E-scheme recombination is assumed to compute the
910// total momentum
911PseudoJet join(const vector<PseudoJet> & pieces){
912 // compute the total momentum
913 //--------------------------------------------------
914 PseudoJet result; // automatically initialised to 0
915 for (unsigned int i=0; i<pieces.size(); i++)
916 result += pieces[i];
917
918 // attach a CompositeJetStructure to the result
919 //--------------------------------------------------
920 CompositeJetStructure *cj_struct = new CompositeJetStructure(pieces);
921
922 result.set_structure_shared_ptr(SharedPtr<PseudoJetStructureBase>(cj_struct));
923
924 return result;
925}
926
927// build a "CompositeJet" from a single PseudoJet
928PseudoJet join(const PseudoJet & j1){
929 return join(vector<PseudoJet>(1,j1));
930}
931
932// build a "CompositeJet" from two PseudoJet
933PseudoJet join(const PseudoJet & j1, const PseudoJet & j2){
934 vector<PseudoJet> pieces;
935 pieces.reserve(2);
936 pieces.push_back(j1);
937 pieces.push_back(j2);
938 return join(pieces);
939}
940
941// build a "CompositeJet" from 3 PseudoJet
942PseudoJet join(const PseudoJet & j1, const PseudoJet & j2, const PseudoJet & j3){
943 vector<PseudoJet> pieces;
944 pieces.reserve(3);
945 pieces.push_back(j1);
946 pieces.push_back(j2);
947 pieces.push_back(j3);
948 return join(pieces);
949}
950
951// build a "CompositeJet" from 4 PseudoJet
952PseudoJet join(const PseudoJet & j1, const PseudoJet & j2, const PseudoJet & j3, const PseudoJet & j4){
953 vector<PseudoJet> pieces;
954 pieces.reserve(4);
955 pieces.push_back(j1);
956 pieces.push_back(j2);
957 pieces.push_back(j3);
958 pieces.push_back(j4);
959 return join(pieces);
960}
961
962
963
964
965FASTJET_END_NAMESPACE
966
base class that sets interface for extensions of ClusterSequence that provide information about the a...
virtual double area(const PseudoJet &) const
return the area associated with the given jet; this base class returns 0.
deals with clustering
base class corresponding to errors that can be thrown by FastJet
Definition: Error.hh:52
Contains any information related to the clustering that should be directly accessible to PseudoJet.
Class to contain pseudojets, including minimal information of use to jet-clustering routines.
Definition: PseudoJet.hh:68
PseudoJet()
default constructor, which as of FJ3.0 provides an object for which all operations are now valid and ...
Definition: PseudoJet.hh:82
void set_cached_rap_phi(double rap, double phi)
in some cases when setting a 4-momentum, the user/program knows what rapidity and azimuth are associa...
Definition: PseudoJet.cc:434
double rap() const
returns the rapidity or some large value when the rapidity is infinite
Definition: PseudoJet.hh:138
const UserInfoBase * user_info_ptr() const
retrieve a pointer to the (const) user information
Definition: PseudoJet.hh:496
virtual bool is_pure_ghost() const
true if this jet is made exclusively of ghosts.
Definition: PseudoJet.cc:840
double phi() const
returns phi (in the range 0..2pi)
Definition: PseudoJet.hh:123
const PseudoJetStructureBase * structure_ptr() const
return a pointer to the structure (of type PseudoJetStructureBase*) associated with this PseudoJet.
Definition: PseudoJet.cc:580
int cluster_hist_index() const
return the cluster_hist_index, intended to be used by clustering routines.
Definition: PseudoJet.hh:834
double m() const
returns the invariant mass (If m2() is negative then -sqrt(-m2()) is returned, as in CLHEP)
Definition: PseudoJet.hh:1064
int user_index() const
return the user_index,
Definition: PseudoJet.hh:389
virtual std::vector< PseudoJet > pieces() const
retrieve the pieces that make up the jet.
Definition: PseudoJet.cc:781
virtual PseudoJet area_4vector() const
return the jet 4-vector area.
Definition: PseudoJet.cc:833
An implementation of shared pointers that is broadly similar to C++11 shared_ptr (https://en....
Definition: SharedPtr.hh:341
Selector operator*(const Selector &s1, const Selector &s2)
successive application of 2 selectors
Definition: Selector.cc:559
vector< PseudoJet > sorted_by_rapidity(const vector< PseudoJet > &jets)
return a vector of jets sorted into increasing rapidity
Definition: PseudoJet.cc:879
const double MaxRap
Used to protect against parton-level events where pt can be zero for some partons,...
Definition: PseudoJet.hh:53
vector< PseudoJet > sorted_by_E(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing energy
Definition: PseudoJet.cc:887
bool have_same_momentum(const PseudoJet &jeta, const PseudoJet &jetb)
returns true if the momenta of the two input jets are identical
Definition: PseudoJet.cc:426
bool operator==(const PseudoJet &a, const PseudoJet &b)
returns true if the 4 momentum components of the two PseudoJets are identical and all the internal in...
Definition: PseudoJet.cc:345
std::vector< T > objects_sorted_by_values(const std::vector< T > &objects, const std::vector< double > &values)
given a vector of values with a one-to-one correspondence with the vector of objects,...
Definition: PseudoJet.hh:988
vector< PseudoJet > sorted_by_pt(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing kt2
Definition: PseudoJet.cc:871
double theta(const PseudoJet &a, const PseudoJet &b)
returns the angle between a and b
Definition: PseudoJet.hh:950
PseudoJet PtYPhiM(double pt, double y, double phi, double m)
return a pseudojet with the given pt, y, phi and mass
Definition: PseudoJet.cc:458
vector< PseudoJet > sorted_by_pz(const vector< PseudoJet > &jets)
return a vector of jets sorted into increasing pz
Definition: PseudoJet.cc:895