FastJet  3.4.0-beta.1
PseudoJet.hh
1 //FJSTARTHEADER
2 // $Id$
3 //
4 // Copyright (c) 2005-2021, 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 #ifndef __FASTJET_PSEUDOJET_HH__
33 #define __FASTJET_PSEUDOJET_HH__
34 
35 #include<valarray>
36 #include<vector>
37 #include<cassert>
38 #include<cmath>
39 #include<iostream>
40 #include "fastjet/config.h"
41 #include "fastjet/internal/numconsts.hh"
42 #include "fastjet/internal/IsBase.hh"
43 #include "fastjet/SharedPtr.hh"
44 #include "fastjet/Error.hh"
45 #include "fastjet/PseudoJetStructureBase.hh"
46 
47 FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
48 
49 //using namespace std;
50 
51 /// Used to protect against parton-level events where pt can be zero
52 /// for some partons, giving rapidity=infinity. KtJet fails in those cases.
53 const double MaxRap = 1e5;
54 
55 /// default value for phi, meaning it (and rapidity) have yet to be calculated)
56 const double pseudojet_invalid_phi = -100.0;
57 const double pseudojet_invalid_rap = -1e200;
58 
59 #ifndef __FJCORE__
60 // forward definition
62 #endif // __FJCORE__
63 
64 /// @ingroup basic_classes
65 /// \class PseudoJet
66 /// Class to contain pseudojets, including minimal information of use to
67 /// jet-clustering routines.
68 class PseudoJet {
69 
70  public:
71  //----------------------------------------------------------------------
72  /// @name Constructors and destructor
73  //\{
74  /// default constructor, which as of FJ3.0 provides an object for
75  /// which all operations are now valid and which has zero momentum
76  ///
77  // (cf. this is actually OK from a timing point of view and in some
78  // cases better than just having the default constructor for the
79  // internal shared pointer: see PJtiming.cc and the notes therein)
80  //
81  // note: no reset of shared pointers needed
82  PseudoJet() : _px(0), _py(0), _pz(0), _E(0) {_finish_init(); _reset_indices();}
83  //PseudoJet() : _px(0), _py(0), _pz(0), _E(0), _phi(pseudojet_invalid_phi), _rap(pseudojet_invalid_rap), _kt2(0) {_reset_indices();}
84  /// construct a pseudojet from explicit components
85  PseudoJet(const double px, const double py, const double pz, const double E);
86 
87  /// constructor from any object that has px,py,pz,E = some_four_vector[0--3],
88  #ifndef SWIG
89  template <class L> PseudoJet(const L & some_four_vector);
90  #endif
91 
92  // Constructor that performs minimal initialisation (only that of
93  // the shared pointers), of use in certain speed-critical contexts
94  //
95  // NB: "dummy" is commented to avoid unused-variable compiler warnings
96  PseudoJet(bool /* dummy */) {}
97 
98 #ifdef FASTJET_HAVE_THREAD_SAFETY
99  PseudoJet(const PseudoJet &other){ (*this)=other; }
100  PseudoJet& operator=(const PseudoJet& other);
101 #endif
102 
103  /// default (virtual) destructor
104  virtual ~PseudoJet(){}
105 //std::shared_ptr: #ifdef FASTJET_HAVE_THREAD_SAFETY
106 //std::shared_ptr: ;
107 //std::shared_ptr: #else
108 //std::shared_ptr: {}
109 //std::shared_ptr: #endif
110  //\} ---- end of constructors and destructors --------------------------
111 
112  //----------------------------------------------------------------------
113  /// @name Kinematic access functions
114  //\{
115  //----------------------------------------------------------------------
116  inline double E() const {return _E;}
117  inline double e() const {return _E;} // like CLHEP
118  inline double px() const {return _px;}
119  inline double py() const {return _py;}
120  inline double pz() const {return _pz;}
121 
122  /// returns phi (in the range 0..2pi)
123  inline double phi() const {return phi_02pi();}
124 
125  /// returns phi in the range -pi..pi
126  inline double phi_std() const {
127  _ensure_valid_rap_phi();
128  return _phi > pi ? _phi-twopi : _phi;}
129 
130  /// returns phi in the range 0..2pi
131  inline double phi_02pi() const {
132  _ensure_valid_rap_phi();
133  return _phi;
134  }
135 
136  /// returns the rapidity or some large value when the rapidity
137  /// is infinite
138  inline double rap() const {
139  _ensure_valid_rap_phi();
140  return _rap;
141  }
142 
143  /// the same as rap()
144  inline double rapidity() const {return rap();} // like CLHEP
145 
146  /// returns the pseudo-rapidity or some large value when the
147  /// rapidity is infinite
148  double pseudorapidity() const;
149  double eta() const {return pseudorapidity();}
150 
151  /// returns the squared transverse momentum
152  inline double pt2() const {return _kt2;}
153  /// returns the scalar transverse momentum
154  inline double pt() const {return sqrt(_kt2);}
155  /// returns the squared transverse momentum
156  inline double perp2() const {return _kt2;} // like CLHEP
157  /// returns the scalar transverse momentum
158  inline double perp() const {return sqrt(_kt2);} // like CLHEP
159  /// returns the squared transverse momentum
160  inline double kt2() const {return _kt2;} // for bkwds compatibility
161 
162  /// returns the squared invariant mass // like CLHEP
163  inline double m2() const {return (_E+_pz)*(_E-_pz)-_kt2;}
164  /// returns the invariant mass
165  /// (If m2() is negative then -sqrt(-m2()) is returned, as in CLHEP)
166  inline double m() const;
167 
168  /// returns the squared transverse mass = kt^2+m^2
169  inline double mperp2() const {return (_E+_pz)*(_E-_pz);}
170  /// returns the transverse mass = sqrt(kt^2+m^2)
171  inline double mperp() const {return sqrt(std::abs(mperp2()));}
172  /// returns the squared transverse mass = kt^2+m^2
173  inline double mt2() const {return (_E+_pz)*(_E-_pz);}
174  /// returns the transverse mass = sqrt(kt^2+m^2)
175  inline double mt() const {return sqrt(std::abs(mperp2()));}
176 
177  /// return the squared 3-vector modulus = px^2+py^2+pz^2
178  inline double modp2() const {return _kt2+_pz*_pz;}
179  /// return the 3-vector modulus = sqrt(px^2+py^2+pz^2)
180  inline double modp() const {return sqrt(_kt2+_pz*_pz);}
181 
182  /// return the transverse energy
183  inline double Et() const {return (_kt2==0) ? 0.0 : _E/sqrt(1.0+_pz*_pz/_kt2);}
184  /// return the transverse energy squared
185  inline double Et2() const {return (_kt2==0) ? 0.0 : _E*_E/(1.0+_pz*_pz/_kt2);}
186 
187  /// cos of the polar angle
188  /// should we have: min(1.0,max(-1.0,_pz/sqrt(modp2())));
189  inline double cos_theta() const {
190  return std::min(1.0, std::max(-1.0, _pz/sqrt(modp2())));
191  }
192  /// polar angle
193  inline double theta() const { return acos(cos_theta()); }
194 
195  /// returns component i, where X==0, Y==1, Z==2, E==3
196  double operator () (int i) const ;
197  /// returns component i, where X==0, Y==1, Z==2, E==3
198  inline double operator [] (int i) const { return (*this)(i); }; // this too
199 
200 
201 
202  /// returns kt distance (R=1) between this jet and another
203  double kt_distance(const PseudoJet & other) const;
204 
205  /// returns squared cylinder (rap-phi) distance between this jet and another
206  double plain_distance(const PseudoJet & other) const;
207  /// returns squared cylinder (rap-phi) distance between this jet and
208  /// another
209  inline double squared_distance(const PseudoJet & other) const {
210  return plain_distance(other);}
211 
212  /// return the cylinder (rap-phi) distance between this jet and another,
213  /// \f$\Delta_R = \sqrt{\Delta y^2 + \Delta \phi^2}\f$.
214  inline double delta_R(const PseudoJet & other) const {
215  return sqrt(squared_distance(other));
216  }
217 
218  /// returns other.phi() - this.phi(), constrained to be in
219  /// range -pi .. pi
220  double delta_phi_to(const PseudoJet & other) const;
221 
222  //// this seemed to compile except if it was used
223  //friend inline double
224  // kt_distance(const PseudoJet & jet1, const PseudoJet & jet2) {
225  // return jet1.kt_distance(jet2);}
226 
227  /// returns distance between this jet and the beam
228  inline double beam_distance() const {return _kt2;}
229 
230  /// return a valarray containing the four-momentum (components 0-2
231  /// are 3-mom, component 3 is energy).
232  std::valarray<double> four_mom() const;
233 
234  //\} ------- end of kinematic access functions
235 
236  // taken from CLHEP
237  enum { X=0, Y=1, Z=2, T=3, NUM_COORDINATES=4, SIZE=NUM_COORDINATES };
238 
239 
240  //----------------------------------------------------------------------
241  /// @name Kinematic modification functions
242  //\{
243  //----------------------------------------------------------------------
244  /// transform this jet (given in the rest frame of prest) into a jet
245  /// in the lab frame
246  PseudoJet & boost(const PseudoJet & prest);
247  /// transform this jet (given in lab) into a jet in the rest
248  /// frame of prest
249  PseudoJet & unboost(const PseudoJet & prest);
250 
251  PseudoJet & operator*=(double);
252  PseudoJet & operator/=(double);
253  PseudoJet & operator+=(const PseudoJet &);
254  PseudoJet & operator-=(const PseudoJet &);
255 
256 //std::shared_ptr: #ifdef FASTJET_HAVE_THREAD_SAFETY
257 //std::shared_ptr: /// overload the assignment through the = operator
258 //std::shared_ptr: ///
259 //std::shared_ptr: /// this is needed to make sure that release_from_cs is called
260 //std::shared_ptr: /// before copying the structure shared pointer!
261 //std::shared_ptr: void operator=(const PseudoJet &other){
262 //std::shared_ptr: _release_jet_from_cs();
263 //std::shared_ptr: _structure = other._structure;
264 //std::shared_ptr: _user_info = other._user_info;
265 //std::shared_ptr:
266 //std::shared_ptr: _px = other._px;
267 //std::shared_ptr: _py = other._py;
268 //std::shared_ptr: _pz = other._pz;
269 //std::shared_ptr: _E = other._E;
270 //std::shared_ptr:
271 //std::shared_ptr: _phi = other._phi;
272 //std::shared_ptr: _rap = other._rap;
273 //std::shared_ptr: _kt2 = other._kt2;
274 //std::shared_ptr:
275 //std::shared_ptr: _cluster_hist_index = other._cluster_hist_index;
276 //std::shared_ptr: _user_index = other._user_index;
277 //std::shared_ptr: }
278 //std::shared_ptr:
279 //std::shared_ptr: /// force resetting the structure pointer to an empty structure
280 //std::shared_ptr: ///
281 //std::shared_ptr: /// THIS IS STRICTLY MEANT FOR INTERNAL USAGE IN SELF-DELETING
282 //std::shared_ptr: /// ClusterSequences. IT SHOULD NOT BE USED BY END-USERS.
283 //std::shared_ptr: void force_reset_structure(){
284 //std::shared_ptr: _structure.reset();
285 //std::shared_ptr: }
286 //std::shared_ptr: #endif
287 
288  /// reset the 4-momentum according to the supplied components and
289  /// put the user and history indices back to their default values
290  inline void reset(double px, double py, double pz, double E);
291 
292  /// reset the PseudoJet to be equal to psjet (including its
293  /// indices); NB if the argument is derived from a PseudoJet then
294  /// the "reset" used will be the templated version
295  ///
296  /// Note: this is included on top of the templated version because
297  /// PseudoJet is not "derived" from PseudoJet, so the templated
298  /// reset would not handle this case properly.
299  inline void reset(const PseudoJet & psjet) {
300  (*this) = psjet;
301  }
302 
303  /// reset the 4-momentum according to the supplied generic 4-vector
304  /// (accessible via indexing, [0]==px,...[3]==E) and put the user
305  /// and history indices back to their default values.
306 #ifndef SWIG
307  template <class L> inline void reset(const L & some_four_vector) {
308  // check if some_four_vector can be cast to a PseudoJet
309  //
310  // Note that a regular dynamic_cast would not work here because
311  // there is no guarantee that L is polymorphic. We use a more
312  // complex construct here that works also in such a case. As for
313  // dynamic_cast, NULL is returned if L is not derived from
314  // PseudoJet
315  //
316  // Note the explicit request for fastjet::cast_if_derived; when
317  // combining fastjet and fjcore, this avoids ambiguity in which of
318  // the two cast_if_derived calls to use.
319  const PseudoJet * pj = fastjet::cast_if_derived<const PseudoJet>(&some_four_vector);
320 
321  if (pj){
322  (*this) = *pj;
323  } else {
324  reset(some_four_vector[0], some_four_vector[1],
325  some_four_vector[2], some_four_vector[3]);
326  }
327  }
328 #endif // SWIG
329 
330  /// reset the PseudoJet according to the specified pt, rapidity,
331  /// azimuth and mass (also resetting indices, etc.)
332  /// (phi should satisfy -2pi<phi<4pi)
333  inline void reset_PtYPhiM(double pt_in, double y_in, double phi_in, double m_in=0.0) {
334  reset_momentum_PtYPhiM(pt_in, y_in, phi_in, m_in);
335  _reset_indices();
336  //std::shared_ptr: _reset_shared_pointers();
337  }
338 
339  /// reset the 4-momentum according to the supplied components
340  /// but leave all other information (indices, user info, etc.)
341  /// untouched
342  inline void reset_momentum(double px, double py, double pz, double E);
343 
344  /// reset the 4-momentum according to the components of the supplied
345  /// PseudoJet, including cached components; note that the template
346  /// version (below) will be called for classes derived from PJ.
347  inline void reset_momentum(const PseudoJet & pj);
348 
349  /// reset the 4-momentum according to the specified pt, rapidity,
350  /// azimuth and mass (phi should satisfy -2pi<phi<4pi)
351  void reset_momentum_PtYPhiM(double pt, double y, double phi, double m=0.0);
352 
353  /// reset the 4-momentum according to the supplied generic 4-vector
354  /// (accessible via indexing, [0]==px,...[3]==E), but leave all
355  /// other information (indices, user info, etc.) untouched
356  template <class L> inline void reset_momentum(const L & some_four_vector) {
357  reset_momentum(some_four_vector[0], some_four_vector[1],
358  some_four_vector[2], some_four_vector[3]);
359  }
360 
361  /// in some cases when setting a 4-momentum, the user/program knows
362  /// what rapidity and azimuth are associated with that 4-momentum;
363  /// by calling this routine the user can provide the information
364  /// directly to the PseudoJet and avoid expensive rap-phi
365  /// recalculations.
366  ///
367  /// - \param rap rapidity
368  /// - \param phi (in range -twopi...4*pi)
369  ///
370  /// USE WITH CAUTION: there are no checks that the rapidity and
371  /// azimuth supplied are sensible, nor does this reset the
372  /// 4-momentum components if things don't match.
373  void set_cached_rap_phi(double rap, double phi);
374 
375 
376  //\} --- end of kin mod functions ------------------------------------
377 
378  //----------------------------------------------------------------------
379  /// @name User index functions
380  ///
381  /// To allow the user to set and access an integer index which can
382  /// be exploited by the user to associate extra information with a
383  /// particle/jet (for example pdg id, or an indication of a
384  /// particle's origin within the user's analysis)
385  //
386  //\{
387 
388  /// return the user_index,
389  inline int user_index() const {return _user_index;}
390  /// set the user_index, intended to allow the user to add simple
391  /// identifying information to a particle/jet
392  inline void set_user_index(const int index) {_user_index = index;}
393 
394  //\} ----- end of use index functions ---------------------------------
395 
396  //----------------------------------------------------------------------
397  /// @name User information types and functions
398  ///
399  /// Allows PseudoJet to carry extra user info (as an object derived from
400  /// UserInfoBase).
401  //\{
402 
403  /// @ingroup user_info
404  /// \class UserInfoBase
405  /// a base class to hold extra user information in a PseudoJet
406  ///
407  /// This is a base class to help associate extra user information
408  /// with a jet. The user should store their information in a class
409  /// derived from this. This allows information of arbitrary
410  /// complexity to be easily associated with a PseudoJet (in contrast
411  /// to the user index). For example, in a Monte Carlo simulation,
412  /// the user information might include the PDG ID, and the position
413  /// of the production vertex for the particle.
414  ///
415  /// The PseudoJet is able to store a shared pointer to any object
416  /// derived from UserInfo. The use of a shared pointer frees the
417  /// user of the need to handle the memory management associated with
418  /// the information.
419  ///
420  /// Having the user information derive from a common base class also
421  /// facilitates dynamic casting, etc.
422  ///
424  public:
425  // dummy ctor
426  UserInfoBase(){};
427 
428  // dummy virtual dtor
429  // makes it polymorphic to allow for dynamic_cast
430  virtual ~UserInfoBase(){};
431  };
432 
433  /// error class to be thrown if accessing user info when it doesn't
434  /// exist
435  class InexistentUserInfo : public Error {
436  public:
438  };
439 
440  /// sets the internal shared pointer to the user information.
441  ///
442  /// Note that the PseudoJet will now _own_ the pointer, and delete
443  /// the corresponding object when it (the jet, and any copies of the jet)
444  /// goes out of scope.
445  void set_user_info(UserInfoBase * user_info_in) {
446  _user_info.reset(user_info_in);
447  }
448 
449  /// returns a reference to the dynamic cast conversion of user_info
450  /// to type L.
451  ///
452  /// Usage: suppose you have previously set the user info with a pointer
453  /// to an object of type MyInfo,
454  ///
455  /// class MyInfo: public PseudoJet::UserInfoBase {
456  /// MyInfo(int id) : _pdg_id(id);
457  /// int pdg_id() const {return _pdg_id;}
458  /// int _pdg_id;
459  /// };
460  ///
461  /// PseudoJet particle(...);
462  /// particle.set_user_info(new MyInfo(its_pdg_id));
463  ///
464  /// Then you would access that pdg_id() as
465  ///
466  /// particle.user_info<MyInfo>().pdg_id();
467  ///
468  /// It's overkill for just a single integer, but scales easily to
469  /// more extensive information.
470  ///
471  /// Note that user_info() throws an InexistentUserInfo() error if
472  /// there is no user info; throws a std::bad_cast if the conversion
473  /// doesn't work
474  ///
475  /// If this behaviour does not fit your needs, use instead the the
476  /// user_info_ptr() or user_info_shared_ptr() member functions.
477  template<class L>
478  const L & user_info() const{
479  if (_user_info.get() == 0) throw InexistentUserInfo();
480  return dynamic_cast<const L &>(* _user_info.get());
481  }
482 
483  /// returns true if the PseudoJet has user information
484  bool has_user_info() const{
485  return _user_info.get();
486  }
487 
488  /// returns true if the PseudoJet has user information than can be
489  /// cast to the template argument type.
490  template<class L>
491  bool has_user_info() const{
492  return _user_info.get() && dynamic_cast<const L *>(_user_info.get());
493  }
494 
495  /// retrieve a pointer to the (const) user information
496  const UserInfoBase * user_info_ptr() const{
497  // the line below is not needed since the next line would anyway
498  // return NULL in that case
499  //if (!_user_info) return NULL;
500  return _user_info.get();
501  }
502 
503 
504  /// retrieve a (const) shared pointer to the user information
506  return _user_info;
507  }
508 
509  /// retrieve a (non-const) shared pointer to the user information;
510  /// you can use this, for example, to set the shared pointer, eg
511  ///
512  /// \code
513  /// p2.user_info_shared_ptr() = p1.user_info_shared_ptr();
514  /// \endcode
515  ///
516  /// or
517  ///
518  /// \code
519  /// SharedPtr<PseudoJet::UserInfoBase> info_shared(new MyInfo(...));
520  /// p2.user_info_shared_ptr() = info_shared;
521  /// \endcode
523  return _user_info;
524  }
525 
526  // \} --- end of extra info functions ---------------------------------
527 
528  //----------------------------------------------------------------------
529  /// @name Description
530  ///
531  /// Since a PseudoJet can have a structure that contains a variety
532  /// of information, we provide a description that allows one to check
533  /// exactly what kind of PseudoJet we are dealing with
534  //
535  //\{
536 
537  /// return a string describing what kind of PseudoJet we are dealing with
538  std::string description() const;
539 
540  //\} ----- end of description functions ---------------------------------
541 
542  //-------------------------------------------------------------
543  /// @name Access to the associated ClusterSequence object.
544  ///
545  /// In addition to having kinematic information, jets may contain a
546  /// reference to an associated ClusterSequence (this is the case,
547  /// for example, if the jet has been returned by a ClusterSequence
548  /// member function).
549  //\{
550  //-------------------------------------------------------------
551  /// returns true if this PseudoJet has an associated ClusterSequence.
552  bool has_associated_cluster_sequence() const;
553  /// shorthand for has_associated_cluster_sequence()
554  bool has_associated_cs() const {return has_associated_cluster_sequence();}
555 
556  /// returns true if this PseudoJet has an associated and still
557  /// valid(ated) ClusterSequence.
558  bool has_valid_cluster_sequence() const;
559  /// shorthand for has_valid_cluster_sequence()
560  bool has_valid_cs() const {return has_valid_cluster_sequence();}
561 
562  /// get a (const) pointer to the parent ClusterSequence (NULL if
563  /// inexistent)
564  const ClusterSequence* associated_cluster_sequence() const;
565  // shorthand for associated_cluster_sequence()
566  const ClusterSequence* associated_cs() const {return associated_cluster_sequence();}
567 
568  /// if the jet has a valid associated cluster sequence then return a
569  /// pointer to it; otherwise throw an error
571  return validated_cs();
572  }
573  /// shorthand for validated_cluster_sequence()
574  const ClusterSequence * validated_cs() const;
575 
576 #ifndef __FJCORE__
577  /// if the jet has valid area information then return a pointer to
578  /// the associated ClusterSequenceAreaBase object; otherwise throw an error
580  return validated_csab();
581  }
582 
583  /// shorthand for validated_cluster_sequence_area_base()
584  const ClusterSequenceAreaBase * validated_csab() const;
585 #endif // __FJCORE__
586 
587  //\}
588 
589  //-------------------------------------------------------------
590  /// @name Access to the associated PseudoJetStructureBase object.
591  ///
592  /// In addition to having kinematic information, jets may contain a
593  /// reference to an associated ClusterSequence (this is the case,
594  /// for example, if the jet has been returned by a ClusterSequence
595  /// member function).
596  //\{
597  //-------------------------------------------------------------
598 
599  /// set the associated structure
600  void set_structure_shared_ptr(const SharedPtr<PseudoJetStructureBase> &structure_in);
601 
602  /// return true if there is some structure associated with this PseudoJet
603  bool has_structure() const;
604 
605  /// return a pointer to the structure (of type
606  /// PseudoJetStructureBase*) associated with this PseudoJet.
607  ///
608  /// return NULL if there is no associated structure
609  const PseudoJetStructureBase* structure_ptr() const;
610 
611  /// return a non-const pointer to the structure (of type
612  /// PseudoJetStructureBase*) associated with this PseudoJet.
613  ///
614  /// return NULL if there is no associated structure
615  ///
616  /// Only use this if you know what you are doing. In any case,
617  /// prefer the 'structure_ptr()' (the const version) to this method,
618  /// unless you really need a write access to the PseudoJet's
619  /// underlying structure.
620  PseudoJetStructureBase* structure_non_const_ptr();
621 
622  /// return a pointer to the structure (of type
623  /// PseudoJetStructureBase*) associated with this PseudoJet.
624  ///
625  /// throw an error if there is no associated structure
626  const PseudoJetStructureBase* validated_structure_ptr() const;
627 
628  /// return a reference to the shared pointer to the
629  /// PseudoJetStructureBase associated with this PseudoJet
630  const SharedPtr<PseudoJetStructureBase> & structure_shared_ptr() const;
631 
632  /// returns a reference to the structure casted to the requested
633  /// structure type
634  ///
635  /// If there is no structure associated, an Error is thrown.
636  /// If the type is not met, a std::bad_cast error is thrown.
637  template<typename StructureType>
638  const StructureType & structure() const;
639 
640  /// check if the PseudoJet has the structure resulting from a Transformer
641  /// (that is, its structure is compatible with a Transformer::StructureType).
642  /// If there is no structure, false is returned.
643  template<typename TransformerType>
644  bool has_structure_of() const;
645 
646  /// this is a helper to access any structure created by a Transformer
647  /// (that is, of type Transformer::StructureType).
648  ///
649  /// If there is no structure, or if the structure is not compatible
650  /// with TransformerType, an error is thrown.
651  template<typename TransformerType>
652  const typename TransformerType::StructureType & structure_of() const;
653 
654  //\}
655 
656  //-------------------------------------------------------------
657  /// @name Methods for access to information about jet structure
658  ///
659  /// These allow access to jet constituents, and other jet
660  /// subtructure information. They only work if the jet is associated
661  /// with a ClusterSequence.
662  //-------------------------------------------------------------
663  //\{
664 
665  /// check if it has been recombined with another PseudoJet in which
666  /// case, return its partner through the argument. Otherwise,
667  /// 'partner' is set to 0.
668  ///
669  /// an Error is thrown if this PseudoJet has no currently valid
670  /// associated ClusterSequence
671  virtual bool has_partner(PseudoJet &partner) const;
672 
673  /// check if it has been recombined with another PseudoJet in which
674  /// case, return its child through the argument. Otherwise, 'child'
675  /// is set to 0.
676  ///
677  /// an Error is thrown if this PseudoJet has no currently valid
678  /// associated ClusterSequence
679  virtual bool has_child(PseudoJet &child) const;
680 
681  /// check if it is the product of a recombination, in which case
682  /// return the 2 parents through the 'parent1' and 'parent2'
683  /// arguments. Otherwise, set these to 0.
684  ///
685  /// an Error is thrown if this PseudoJet has no currently valid
686  /// associated ClusterSequence
687  virtual bool has_parents(PseudoJet &parent1, PseudoJet &parent2) const;
688 
689  /// check if the current PseudoJet contains the one passed as
690  /// argument.
691  ///
692  /// an Error is thrown if this PseudoJet has no currently valid
693  /// associated ClusterSequence
694  virtual bool contains(const PseudoJet &constituent) const;
695 
696  /// check if the current PseudoJet is contained the one passed as
697  /// argument.
698  ///
699  /// an Error is thrown if this PseudoJet has no currently valid
700  /// associated ClusterSequence
701  virtual bool is_inside(const PseudoJet &jet) const;
702 
703 
704  /// returns true if the PseudoJet has constituents
705  virtual bool has_constituents() const;
706 
707  /// retrieve the constituents.
708  ///
709  /// an Error is thrown if this PseudoJet has no currently valid
710  /// associated ClusterSequence or other substructure information
711  virtual std::vector<PseudoJet> constituents() const;
712 
713 
714  /// returns true if the PseudoJet has support for exclusive subjets
715  virtual bool has_exclusive_subjets() const;
716 
717  /// return a vector of all subjets of the current jet (in the sense
718  /// of the exclusive algorithm) that would be obtained when running
719  /// the algorithm with the given dcut.
720  ///
721  /// Time taken is O(m ln m), where m is the number of subjets that
722  /// are found. If m gets to be of order of the total number of
723  /// constituents in the jet, this could be substantially slower than
724  /// just getting that list of constituents.
725  ///
726  /// an Error is thrown if this PseudoJet has no currently valid
727  /// associated ClusterSequence
728  std::vector<PseudoJet> exclusive_subjets (const double dcut) const;
729 
730  /// return the size of exclusive_subjets(...); still n ln n with same
731  /// coefficient, but marginally more efficient than manually taking
732  /// exclusive_subjets.size()
733  ///
734  /// an Error is thrown if this PseudoJet has no currently valid
735  /// associated ClusterSequence
736  int n_exclusive_subjets(const double dcut) const;
737 
738  /// return the list of subjets obtained by unclustering the supplied
739  /// jet down to nsub subjets. Throws an error if there are fewer than
740  /// nsub particles in the jet.
741  ///
742  /// For ClusterSequence type jets, requires nsub ln nsub time
743  ///
744  /// An Error is thrown if this PseudoJet has no currently valid
745  /// associated ClusterSequence
746  std::vector<PseudoJet> exclusive_subjets (int nsub) const;
747 
748  /// return the list of subjets obtained by unclustering the supplied
749  /// jet down to nsub subjets (or all constituents if there are fewer
750  /// than nsub).
751  ///
752  /// For ClusterSequence type jets, requires nsub ln nsub time
753  ///
754  /// An Error is thrown if this PseudoJet has no currently valid
755  /// associated ClusterSequence
756  std::vector<PseudoJet> exclusive_subjets_up_to (int nsub) const;
757 
758  /// Returns the dij that was present in the merging nsub+1 -> nsub
759  /// subjets inside this jet.
760  ///
761  /// Returns 0 if there were nsub or fewer constituents in the jet.
762  ///
763  /// an Error is thrown if this PseudoJet has no currently valid
764  /// associated ClusterSequence
765  double exclusive_subdmerge(int nsub) const;
766 
767  /// Returns the maximum dij that occurred in the whole event at the
768  /// stage that the nsub+1 -> nsub merge of subjets occurred inside
769  /// this jet.
770  ///
771  /// Returns 0 if there were nsub or fewer constituents in the jet.
772  ///
773  /// an Error is thrown if this PseudoJet has no currently valid
774  /// associated ClusterSequence
775  double exclusive_subdmerge_max(int nsub) const;
776 
777 
778  /// returns true if a jet has pieces
779  ///
780  /// By default a single particle or a jet coming from a
781  /// ClusterSequence have no pieces and this methos will return false.
782  ///
783  /// In practice, this is equivalent to have an structure of type
784  /// CompositeJetStructure.
785  virtual bool has_pieces() const;
786 
787 
788  /// retrieve the pieces that make up the jet.
789  ///
790  /// If the jet does not support pieces, an error is throw
791  virtual std::vector<PseudoJet> pieces() const;
792 
793 
794  // the following ones require a computation of the area in the
795  // parent ClusterSequence (See ClusterSequenceAreaBase for details)
796  //------------------------------------------------------------------
797 #ifndef __FJCORE__
798 
799  /// check if it has a defined area
800  virtual bool has_area() const;
801 
802  /// return the jet (scalar) area.
803  /// throws an Error if there is no support for area in the parent CS
804  virtual double area() const;
805 
806  /// return the error (uncertainty) associated with the determination
807  /// of the area of this jet.
808  /// throws an Error if there is no support for area in the parent CS
809  virtual double area_error() const;
810 
811  /// return the jet 4-vector area.
812  /// throws an Error if there is no support for area in the parent CS
813  virtual PseudoJet area_4vector() const;
814 
815  /// true if this jet is made exclusively of ghosts.
816  /// throws an Error if there is no support for area in the parent CS
817  virtual bool is_pure_ghost() const;
818 
819 #endif // __FJCORE__
820  //\} --- end of jet structure -------------------------------------
821 
822 
823 
824  //----------------------------------------------------------------------
825  /// @name Members mainly intended for internal use
826  //----------------------------------------------------------------------
827  //\{
828  /// return the cluster_hist_index, intended to be used by clustering
829  /// routines.
830  inline int cluster_hist_index() const {return _cluster_hist_index;}
831  /// set the cluster_hist_index, intended to be used by clustering routines.
832  inline void set_cluster_hist_index(const int index) {_cluster_hist_index = index;}
833 
834  /// alternative name for cluster_hist_index() [perhaps more meaningful]
835  inline int cluster_sequence_history_index() const {
836  return cluster_hist_index();}
837  /// alternative name for set_cluster_hist_index(...) [perhaps more
838  /// meaningful]
839  inline void set_cluster_sequence_history_index(const int index) {
840  set_cluster_hist_index(index);}
841 
842  //\} ---- end of internal use functions ---------------------------
843 
844  protected:
845 
847  SharedPtr<UserInfoBase> _user_info;
848 
849 
850  private:
851  // NB: following order must be kept for things to behave sensibly...
852  double _px,_py,_pz,_E;
853  mutable double _phi, _rap;
854  double _kt2;
855  int _cluster_hist_index, _user_index;
856 
857 #ifdef FASTJET_HAVE_THREAD_SAFETY
858  enum {
859  Init_Done=1,
860  Init_NotDone=0,
861  Init_InProgress=-1
862  };
863 
864  mutable std::atomic<int> _init_status;
865 #endif
866 
867  /// calculate phi, rap, kt2 based on the 4-momentum components
868  void _finish_init();
869  /// set the indices to default values
870  void _reset_indices();
871 
872 //std::shared_ptr: /// reset the shared pointers to empty ones
873 //std::shared_ptr: void _reset_shared_pointers();
874 //std::shared_ptr:
875 //std::shared_ptr: #ifdef FASTJET_HAVE_THREAD_SAFETY
876 //std::shared_ptr: /// For jets associated with a ClusterSequence, this will "free" the
877 //std::shared_ptr: /// jet from the ClusterSequence. This means that, for self-deleting
878 //std::shared_ptr: /// cluster sequences, it will reset the structure pointer and check
879 //std::shared_ptr: /// if the CS needs to be deleted.
880 //std::shared_ptr: ///
881 //std::shared_ptr: /// This is the replacement mechanism for self-deleting cs (with
882 //std::shared_ptr: /// set_count not working with std shared_ptr) and...
883 //std::shared_ptr: ///
884 //std::shared_ptr: /// IT HAS TO BE CALLED BEFORE ANY CHANGE OF THE JET STRUCTURE
885 //std::shared_ptr: /// POINTER
886 //std::shared_ptr: void _release_jet_from_cs();
887 //std::shared_ptr: #endif //FASTJET_HAVE_THREAD_SAFETY
888 
889  /// ensure that the internal values for rapidity and phi
890  /// correspond to 4-momentum structure
891 #ifdef FASTJET_HAVE_THREAD_SAFETY
892  void _ensure_valid_rap_phi() const;
893 #else
894  inline void _ensure_valid_rap_phi() const{
895  if (_phi == pseudojet_invalid_phi) _set_rap_phi();
896  }
897 #endif
898 
899  /// set cached rapidity and phi values
900  void _set_rap_phi() const;
901 
902  // needed for operator* to have access to _ensure_valid_rap_phi()
903  friend PseudoJet operator*(double, const PseudoJet &);
904 };
905 
906 
907 //----------------------------------------------------------------------
908 // routines for basic binary operations
909 
910 PseudoJet operator+(const PseudoJet &, const PseudoJet &);
911 PseudoJet operator-(const PseudoJet &, const PseudoJet &);
912 PseudoJet operator*(double, const PseudoJet &);
913 PseudoJet operator*(const PseudoJet &, double);
914 PseudoJet operator/(const PseudoJet &, double);
915 
916 /// returns true if the 4 momentum components of the two PseudoJets
917 /// are identical and all the internal indices (user, cluster_history)
918 /// + structure and user-info shared pointers are too
919 bool operator==(const PseudoJet &, const PseudoJet &);
920 
921 /// inequality test which is exact opposite of operator==
922 inline bool operator!=(const PseudoJet & a, const PseudoJet & b) {return !(a==b);}
923 
924 /// Can only be used with val=0 and tests whether all four
925 /// momentum components are equal to val (=0.0)
926 bool operator==(const PseudoJet & jet, const double val);
927 inline bool operator==(const double val, const PseudoJet & jet) {return jet == val;}
928 
929 /// Can only be used with val=0 and tests whether at least one of the
930 /// four momentum components is different from val (=0.0)
931 inline bool operator!=(const PseudoJet & a, const double val) {return !(a==val);}
932 inline bool operator!=( const double val, const PseudoJet & a) {return !(a==val);}
933 
934 /// returns the 4-vector dot product of a and b
935 inline double dot_product(const PseudoJet & a, const PseudoJet & b) {
936  return a.E()*b.E() - a.px()*b.px() - a.py()*b.py() - a.pz()*b.pz();
937 }
938 
939 /// returns the cosine of the angle between a and b
940 inline double cos_theta(const PseudoJet & a, const PseudoJet & b) {
941  double dot_3d = a.px()*b.px() + a.py()*b.py() + a.pz()*b.pz();
942  return std::min(1.0, std::max(-1.0, dot_3d/sqrt(a.modp2()*b.modp2())));
943 }
944 
945 /// returns the angle between a and b
946 inline double theta(const PseudoJet & a, const PseudoJet & b) {
947  return acos(cos_theta(a,b));
948 }
949 
950 /// returns true if the momenta of the two input jets are identical
951 bool have_same_momentum(const PseudoJet &, const PseudoJet &);
952 
953 /// return a pseudojet with the given pt, y, phi and mass
954 /// (phi should satisfy -2pi<phi<4pi)
955 PseudoJet PtYPhiM(double pt, double y, double phi, double m = 0.0);
956 
957 //----------------------------------------------------------------------
958 // Routines to do with providing sorted arrays of vectors.
959 
960 /// return a vector of jets sorted into decreasing transverse momentum
961 std::vector<PseudoJet> sorted_by_pt(const std::vector<PseudoJet> & jets);
962 
963 /// return a vector of jets sorted into increasing rapidity
964 std::vector<PseudoJet> sorted_by_rapidity(const std::vector<PseudoJet> & jets);
965 
966 /// return a vector of jets sorted into decreasing energy
967 std::vector<PseudoJet> sorted_by_E(const std::vector<PseudoJet> & jets);
968 
969 /// return a vector of jets sorted into increasing pz
970 std::vector<PseudoJet> sorted_by_pz(const std::vector<PseudoJet> & jets);
971 
972 //----------------------------------------------------------------------
973 // some code to help sorting
974 
975 /// sort the indices so that values[indices[0->n-1]] is sorted
976 /// into increasing order
977 void sort_indices(std::vector<int> & indices,
978  const std::vector<double> & values);
979 
980 /// given a vector of values with a one-to-one correspondence with the
981 /// vector of objects, sort objects into an order such that the
982 /// associated values would be in increasing order (but don't actually
983 /// touch the values vector in the process).
984 template<class T> std::vector<T> objects_sorted_by_values(const std::vector<T> & objects,
985  const std::vector<double> & values) {
986  //assert(objects.size() == values.size());
987  if (objects.size() != values.size()){
988  throw Error("fastjet::objects_sorted_by_values(...): the size of the 'objects' vector must match the size of the 'values' vector");
989  }
990 
991  // get a vector of indices
992  std::vector<int> indices(values.size());
993  for (size_t i = 0; i < indices.size(); i++) {indices[i] = i;}
994 
995  // sort the indices
996  sort_indices(indices, values);
997 
998  // copy the objects
999  std::vector<T> objects_sorted(objects.size());
1000 
1001  // place the objects in the correct order
1002  for (size_t i = 0; i < indices.size(); i++) {
1003  objects_sorted[i] = objects[indices[i]];
1004  }
1005 
1006  return objects_sorted;
1007 }
1008 
1009 /// \if internal_doc
1010 /// @ingroup internal
1011 /// \class IndexedSortHelper
1012 /// a class that helps us carry out indexed sorting.
1013 /// \endif
1014 class IndexedSortHelper {
1015 public:
1016  inline IndexedSortHelper (const std::vector<double> * reference_values) {
1017  _ref_values = reference_values;
1018  };
1019  inline int operator() (const int i1, const int i2) const {
1020  return (*_ref_values)[i1] < (*_ref_values)[i2];
1021  };
1022 private:
1023  const std::vector<double> * _ref_values;
1024 };
1025 
1026 
1027 //----------------------------------------------------------------------
1028 /// constructor from any object that has px,py,pz,E = some_four_vector[0--3],
1029 // NB: do not know if it really needs to be inline, but when it wasn't
1030 // linking failed with g++ (who knows what was wrong...)
1031 #ifndef SWIG
1032 template <class L> inline PseudoJet::PseudoJet(const L & some_four_vector) {
1033  reset(some_four_vector);
1034 }
1035 #endif
1036 
1037 //----------------------------------------------------------------------
1038 inline void PseudoJet::_reset_indices() {
1039  set_cluster_hist_index(-1);
1040  set_user_index(-1);
1041 
1042  _structure.reset();
1043  _user_info.reset();
1044 }
1045 
1046 //std::shared_ptr: inline void PseudoJet::_reset_shared_pointers() {
1047 //std::shared_ptr: #ifdef FASTJET_HAVE_THREAD_SAFETY
1048 //std::shared_ptr: // if the jet currently belongs to a cs, we need to release it before any chenge
1049 //std::shared_ptr: _release_jet_from_cs();
1050 //std::shared_ptr: #endif // FASTJET_HAVE_THREAD_SAFETY
1051 //std::shared_ptr: _structure.reset();
1052 //std::shared_ptr: _user_info.reset();
1053 //std::shared_ptr: // and do not forget to remove it in reset_indices above
1054 //std::shared_ptr: }
1055 
1056 
1057 
1058 
1059 // taken literally from CLHEP
1060 inline double PseudoJet::m() const {
1061  double mm = m2();
1062  return mm < 0.0 ? -std::sqrt(-mm) : std::sqrt(mm);
1063 }
1064 
1065 
1066 inline void PseudoJet::reset(double px_in, double py_in, double pz_in, double E_in) {
1067  _px = px_in;
1068  _py = py_in;
1069  _pz = pz_in;
1070  _E = E_in;
1071  _finish_init();
1072  _reset_indices();
1073  //std::shared_ptr: _reset_shared_pointers();
1074 }
1075 
1076 inline void PseudoJet::reset_momentum(double px_in, double py_in, double pz_in, double E_in) {
1077  _px = px_in;
1078  _py = py_in;
1079  _pz = pz_in;
1080  _E = E_in;
1081  _finish_init();
1082 }
1083 
1084 inline void PseudoJet::reset_momentum(const PseudoJet & pj) {
1085  _px = pj._px ;
1086  _py = pj._py ;
1087  _pz = pj._pz ;
1088  _E = pj._E ;
1089  _phi = pj._phi;
1090  _rap = pj._rap;
1091  _kt2 = pj._kt2;
1092 }
1093 
1094 //-------------------------------------------------------------------------------
1095 // implementation of the templated accesses to the underlying structyre
1096 //-------------------------------------------------------------------------------
1097 
1098 // returns a reference to the structure casted to the requested
1099 // structure type
1100 //
1101 // If there is no sructure associated, an Error is thrown.
1102 // If the type is not met, a std::bad_cast error is thrown.
1103 template<typename StructureType>
1104 const StructureType & PseudoJet::structure() const{
1105  return dynamic_cast<const StructureType &>(* validated_structure_ptr());
1106 
1107 }
1108 
1109 // check if the PseudoJet has the structure resulting from a Transformer
1110 // (that is, its structure is compatible with a Transformer::StructureType)
1111 template<typename TransformerType>
1112 bool PseudoJet::has_structure_of() const{
1113  if (!_structure) return false;
1114 
1115  return dynamic_cast<const typename TransformerType::StructureType *>(_structure.get()) != 0;
1116 }
1117 
1118 // this is a helper to access a structure created by a Transformer
1119 // (that is, of type Transformer::StructureType)
1120 // NULL is returned if the corresponding type is not met
1121 template<typename TransformerType>
1122 const typename TransformerType::StructureType & PseudoJet::structure_of() const{
1123  if (!_structure)
1124  throw Error("Trying to access the structure of a PseudoJet without an associated structure");
1125 
1126  return dynamic_cast<const typename TransformerType::StructureType &>(*_structure);
1127 }
1128 
1129 
1130 
1131 //-------------------------------------------------------------------------------
1132 // helper functions to build a jet made of pieces
1133 //
1134 // Note that there are more complete versions of these functions, with
1135 // an additional argument for a recombination scheme, in
1136 // JetDefinition.hh
1137 // -------------------------------------------------------------------------------
1138 
1139 /// build a "CompositeJet" from the vector of its pieces
1140 ///
1141 /// In this case, E-scheme recombination is assumed to compute the
1142 /// total momentum
1143 PseudoJet join(const std::vector<PseudoJet> & pieces);
1144 
1145 /// build a MergedJet from a single PseudoJet
1146 PseudoJet join(const PseudoJet & j1);
1147 
1148 /// build a MergedJet from 2 PseudoJet
1149 PseudoJet join(const PseudoJet & j1, const PseudoJet & j2);
1150 
1151 /// build a MergedJet from 3 PseudoJet
1152 PseudoJet join(const PseudoJet & j1, const PseudoJet & j2, const PseudoJet & j3);
1153 
1154 /// build a MergedJet from 4 PseudoJet
1155 PseudoJet join(const PseudoJet & j1, const PseudoJet & j2, const PseudoJet & j3, const PseudoJet & j4);
1156 
1157 
1158 
1159 FASTJET_END_NAMESPACE
1160 
1161 #endif // __FASTJET_PSEUDOJET_HH__
fastjet::PseudoJet::modp
double modp() const
return the 3-vector modulus = sqrt(px^2+py^2+pz^2)
Definition: PseudoJet.hh:180
fastjet::PseudoJet::set_user_info
void set_user_info(UserInfoBase *user_info_in)
sets the internal shared pointer to the user information.
Definition: PseudoJet.hh:445
fastjet::SharedPtr::reset
void reset()
reset the pointer to default value (NULL)
Definition: SharedPtr.hh:381
fastjet::objects_sorted_by_values
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:984
fastjet::MaxRap
const double MaxRap
Used to protect against parton-level events where pt can be zero for some partons,...
Definition: PseudoJet.hh:53
fastjet::PseudoJet::phi_02pi
double phi_02pi() const
returns phi in the range 0..2pi
Definition: PseudoJet.hh:131
fastjet::PseudoJet::mperp2
double mperp2() const
returns the squared transverse mass = kt^2+m^2
Definition: PseudoJet.hh:169
fastjet::PseudoJet::user_info_ptr
const UserInfoBase * user_info_ptr() const
retrieve a pointer to the (const) user information
Definition: PseudoJet.hh:496
fastjet::PseudoJet::modp2
double modp2() const
return the squared 3-vector modulus = px^2+py^2+pz^2
Definition: PseudoJet.hh:178
fastjet::PseudoJet::Et
double Et() const
return the transverse energy
Definition: PseudoJet.hh:183
fastjet::PseudoJet::delta_R
double delta_R(const PseudoJet &other) const
return the cylinder (rap-phi) distance between this jet and another, .
Definition: PseudoJet.hh:214
fastjet::sorted_by_rapidity
vector< PseudoJet > sorted_by_rapidity(const vector< PseudoJet > &jets)
return a vector of jets sorted into increasing rapidity
Definition: PseudoJet.cc:890
fastjet::PseudoJet::has_valid_cs
bool has_valid_cs() const
shorthand for has_valid_cluster_sequence()
Definition: PseudoJet.hh:560
fastjet::sorted_by_pz
vector< PseudoJet > sorted_by_pz(const vector< PseudoJet > &jets)
return a vector of jets sorted into increasing pz
Definition: PseudoJet.cc:906
fastjet::SharedPtr::get
T * get() const
get the stored pointer
Definition: SharedPtr.hh:473
fastjet::PseudoJet::~PseudoJet
virtual ~PseudoJet()
default (virtual) destructor
Definition: PseudoJet.hh:104
fastjet::PseudoJet::reset_momentum
void reset_momentum(const L &some_four_vector)
reset the 4-momentum according to the supplied generic 4-vector (accessible via indexing,...
Definition: PseudoJet.hh:356
fastjet::ClusterSequence
deals with clustering
Definition: ClusterSequence.hh:63
fastjet::PseudoJet::mt2
double mt2() const
returns the squared transverse mass = kt^2+m^2
Definition: PseudoJet.hh:173
fastjet::PseudoJet::user_info_shared_ptr
SharedPtr< UserInfoBase > & user_info_shared_ptr()
retrieve a (non-const) shared pointer to the user information; you can use this, for example,...
Definition: PseudoJet.hh:522
fastjet::PseudoJet::squared_distance
double squared_distance(const PseudoJet &other) const
returns squared cylinder (rap-phi) distance between this jet and another
Definition: PseudoJet.hh:209
fastjet::operator*
Selector operator*(const Selector &s1, const Selector &s2)
successive application of 2 selectors
Definition: Selector.cc:559
fastjet::dot_product
double dot_product(const PseudoJet &a, const PseudoJet &b)
returns the 4-vector dot product of a and b
Definition: PseudoJet.hh:935
fastjet::PseudoJet::pt
double pt() const
returns the scalar transverse momentum
Definition: PseudoJet.hh:154
fastjet::PseudoJet::rapidity
double rapidity() const
the same as rap()
Definition: PseudoJet.hh:144
fastjet::PseudoJet::m2
double m2() const
returns the squared invariant mass // like CLHEP
Definition: PseudoJet.hh:163
fastjet::sorted_by_E
vector< PseudoJet > sorted_by_E(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing energy
Definition: PseudoJet.cc:898
fastjet::SharedPtr
An implementation of shared pointers that is broadly similar to C++11 shared_ptr (https://en....
Definition: SharedPtr.hh:341
fastjet::PseudoJet::perp
double perp() const
returns the scalar transverse momentum
Definition: PseudoJet.hh:158
fastjet::PseudoJet::beam_distance
double beam_distance() const
returns distance between this jet and the beam
Definition: PseudoJet.hh:228
fastjet::PseudoJetStructureBase
Contains any information related to the clustering that should be directly accessible to PseudoJet.
Definition: PseudoJetStructureBase.hh:59
fastjet::PseudoJet::cos_theta
double cos_theta() const
cos of the polar angle should we have: min(1.0,max(-1.0,_pz/sqrt(modp2())));
Definition: PseudoJet.hh:189
fastjet::PseudoJet::set_user_index
void set_user_index(const int index)
set the user_index, intended to allow the user to add simple identifying information to a particle/je...
Definition: PseudoJet.hh:392
fastjet::PseudoJet::user_index
int user_index() const
return the user_index,
Definition: PseudoJet.hh:389
fastjet::operator==
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:356
fastjet::PseudoJet::mt
double mt() const
returns the transverse mass = sqrt(kt^2+m^2)
Definition: PseudoJet.hh:175
fastjet::PseudoJet::mperp
double mperp() const
returns the transverse mass = sqrt(kt^2+m^2)
Definition: PseudoJet.hh:171
fastjet::PseudoJet::PseudoJet
PseudoJet()
default constructor, which as of FJ3.0 provides an object for which all operations are now valid and ...
Definition: PseudoJet.hh:82
fastjet::PseudoJet::UserInfoBase
a base class to hold extra user information in a PseudoJet
Definition: PseudoJet.hh:423
fastjet::PseudoJet::kt2
double kt2() const
returns the squared transverse momentum
Definition: PseudoJet.hh:160
fastjet::PseudoJet::theta
double theta() const
polar angle
Definition: PseudoJet.hh:193
fastjet::PseudoJet::cluster_hist_index
int cluster_hist_index() const
return the cluster_hist_index, intended to be used by clustering routines.
Definition: PseudoJet.hh:830
fastjet::PseudoJet::InexistentUserInfo
error class to be thrown if accessing user info when it doesn't exist
Definition: PseudoJet.hh:435
fastjet::PseudoJet::has_associated_cs
bool has_associated_cs() const
shorthand for has_associated_cluster_sequence()
Definition: PseudoJet.hh:554
fastjet::PseudoJet::has_user_info
bool has_user_info() const
returns true if the PseudoJet has user information than can be cast to the template argument type.
Definition: PseudoJet.hh:491
fastjet::PseudoJet::validated_cluster_sequence
const ClusterSequence * validated_cluster_sequence() const
if the jet has a valid associated cluster sequence then return a pointer to it; otherwise throw an er...
Definition: PseudoJet.hh:570
fastjet::PseudoJet
Class to contain pseudojets, including minimal information of use to jet-clustering routines.
Definition: PseudoJet.hh:68
fastjet::ClusterSequenceAreaBase
base class that sets interface for extensions of ClusterSequence that provide information about the a...
Definition: ClusterSequenceAreaBase.hh:48
fastjet::PseudoJet::rap
double rap() const
returns the rapidity or some large value when the rapidity is infinite
Definition: PseudoJet.hh:138
fastjet::PseudoJet::user_info_shared_ptr
const SharedPtr< UserInfoBase > & user_info_shared_ptr() const
retrieve a (const) shared pointer to the user information
Definition: PseudoJet.hh:505
fastjet::PseudoJet::cluster_sequence_history_index
int cluster_sequence_history_index() const
alternative name for cluster_hist_index() [perhaps more meaningful]
Definition: PseudoJet.hh:835
fastjet::operator!=
bool operator!=(const PseudoJet &a, const PseudoJet &b)
inequality test which is exact opposite of operator==
Definition: PseudoJet.hh:922
fastjet::PseudoJet::set_cluster_sequence_history_index
void set_cluster_sequence_history_index(const int index)
alternative name for set_cluster_hist_index(...) [perhaps more meaningful]
Definition: PseudoJet.hh:839
fastjet::pseudojet_invalid_phi
const double pseudojet_invalid_phi
default value for phi, meaning it (and rapidity) have yet to be calculated)
Definition: PseudoJet.hh:56
fastjet::cos_theta
double cos_theta(const PseudoJet &a, const PseudoJet &b)
returns the cosine of the angle between a and b
Definition: PseudoJet.hh:940
fastjet::PseudoJet::has_user_info
bool has_user_info() const
returns true if the PseudoJet has user information
Definition: PseudoJet.hh:484
fastjet::PseudoJet::set_cluster_hist_index
void set_cluster_hist_index(const int index)
set the cluster_hist_index, intended to be used by clustering routines.
Definition: PseudoJet.hh:832
fastjet::PseudoJet::reset_PtYPhiM
void reset_PtYPhiM(double pt_in, double y_in, double phi_in, double m_in=0.0)
reset the PseudoJet according to the specified pt, rapidity, azimuth and mass (also resetting indices...
Definition: PseudoJet.hh:333
fastjet::PseudoJet::phi
double phi() const
returns phi (in the range 0..2pi)
Definition: PseudoJet.hh:123
fastjet::PseudoJet::validated_cluster_sequence_area_base
const ClusterSequenceAreaBase * validated_cluster_sequence_area_base() const
if the jet has valid area information then return a pointer to the associated ClusterSequenceAreaBase...
Definition: PseudoJet.hh:579
fastjet::PseudoJet::reset
void reset(const L &some_four_vector)
reset the 4-momentum according to the supplied generic 4-vector (accessible via indexing,...
Definition: PseudoJet.hh:307
fastjet::PseudoJet::perp2
double perp2() const
returns the squared transverse momentum
Definition: PseudoJet.hh:156
fastjet::PseudoJet::phi_std
double phi_std() const
returns phi in the range -pi..pi
Definition: PseudoJet.hh:126
fastjet::theta
double theta(const PseudoJet &a, const PseudoJet &b)
returns the angle between a and b
Definition: PseudoJet.hh:946
fastjet::PseudoJet::pt2
double pt2() const
returns the squared transverse momentum
Definition: PseudoJet.hh:152
fastjet::PseudoJet::Et2
double Et2() const
return the transverse energy squared
Definition: PseudoJet.hh:185
fastjet::sorted_by_pt
vector< PseudoJet > sorted_by_pt(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing kt2
Definition: PseudoJet.cc:882
fastjet::PseudoJet::reset
void reset(const PseudoJet &psjet)
reset the PseudoJet to be equal to psjet (including its indices); NB if the argument is derived from ...
Definition: PseudoJet.hh:299
fastjet::Error
base class corresponding to errors that can be thrown by FastJet
Definition: Error.hh:51
fastjet::PseudoJet::user_info
const L & user_info() const
returns a reference to the dynamic cast conversion of user_info to type L.
Definition: PseudoJet.hh:478