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