#include <PseudoJet.hh>
Public Types | |
enum | { X = 0, Y = 1, Z = 2, T = 3, NUM_COORDINATES = 4, SIZE = NUM_COORDINATES } |
Public Member Functions | |
PseudoJet () | |
PseudoJet (const double px, const double py, const double pz, const double E) | |
construct a pseudojet from explicit components | |
template<class L> | |
PseudoJet (const L &some_four_vector) | |
constructor from any object that has px,py,pz,E = some_four_vector[0--3], | |
double | E () const |
double | e () const |
double | px () const |
double | py () const |
double | pz () const |
const double | phi () const |
returns phi (in the range 0..2pi) | |
const double | phi_std () const |
returns phi in the range -pi..pi | |
const double | phi_02pi () const |
returns phi in the range 0..2pi | |
double | rap () const |
returns the rapidity or some large value when the rapidity is infinite | |
double | rapidity () const |
the same as rap() | |
double | kt2 () const |
returns the squared transverse momentum | |
double | perp2 () const |
returns the squared transverse momentum | |
double | perp () const |
returns the scalar transverse momentum | |
double | m2 () const |
returns the squared invariant mass // like CLHEP | |
double | mperp2 () const |
returns the squared transverse mass = kt^2+m^2 | |
double | mperp () const |
returns the transverse mass = sqrt(kt^2+m^2) | |
double | m () const |
returns the invariant mass (If m2() is negative then -sqrt(-m2()) is returned, as in CLHEP) | |
double | operator() (int i) const |
returns component i, where X==0, Y==1, Z==2, E==3 | |
double | operator[] (int i) const |
returns component i, where X==0, Y==1, Z==2, E==3 | |
PseudoJet & | boost (const PseudoJet &prest) |
transform this jet (given in the rest frame of prest) into a jet in the lab frame; | |
PseudoJet & | unboost (const PseudoJet &prest) |
transform this jet (given in the rest frame of prest) into a jet in the lab frame; | |
const int & | cluster_hist_index () const |
return the cluster_hist_index, intended to be used by clustering routines. | |
void | set_cluster_hist_index (const int index) |
set the cluster_hist_index, intended to be used by clustering routines. | |
const int | cluster_sequence_history_index () const |
alternative name for cluster_hist_index() [perhaps more meaningful] | |
void | set_cluster_sequence_history_index (const int index) |
alternative name for set_cluster_hist_index(. | |
const int & | user_index () const |
return the user_index, intended to allow the user to "add" information | |
void | set_user_index (const int index) |
set the user_index, intended to allow the user to "add" information | |
std::valarray< double > | four_mom () const |
return a valarray containing the four-momentum (components 0-2 are 3-mom, component 3 is energy). | |
double | kt_distance (const PseudoJet &other) const |
returns kt distance (R=1) between this jet and another | |
double | plain_distance (const PseudoJet &other) const |
returns squared cylinder (rap-phi) distance between this jet and another | |
double | squared_distance (const PseudoJet &other) const |
returns squared cylinder (rap-phi) distance between this jet and another | |
double | beam_distance () const |
returns distance between this jet and the beam | |
void | operator *= (double) |
multiply the jet's momentum by the coefficient | |
void | operator/= (double) |
divide the jet's momentum by the coefficient | |
void | operator+= (const PseudoJet &) |
add the other jet's momentum to this jet | |
void | operator-= (const PseudoJet &) |
subtract the other jet's momentum from this jet | |
template<> | |
PseudoJet (const Cmomentum &four_vector) | |
shortcut for converting siscone Cmomentum into PseudoJet | |
Private Member Functions | |
void | _finish_init () |
do standard end of initialisation | |
Private Attributes | |
double | _px |
double | _py |
double | _pz |
double | _E |
double | _phi |
double | _rap |
double | _kt2 |
int | _cluster_hist_index |
int | _user_index |
Definition at line 54 of file PseudoJet.hh.
|
Definition at line 111 of file PseudoJet.hh. 00111 { X=0, Y=1, Z=2, T=3, NUM_COORDINATES=4, SIZE=NUM_COORDINATES };
|
|
Definition at line 56 of file PseudoJet.hh. Referenced by PseudoJet(). 00056 {};
|
|
construct a pseudojet from explicit components
Definition at line 46 of file PseudoJet.cc. References _E, _finish_init(), _px, _py, _pz, set_cluster_hist_index(), and set_user_index(). 00046 { 00047 00048 _E = E ; 00049 _px = px; 00050 _py = py; 00051 _pz = pz; 00052 00053 this->_finish_init(); 00054 00055 // some default values for these two indices 00056 set_cluster_hist_index(-1); 00057 set_user_index(-1); 00058 00059 }
|
|
constructor from any object that has px,py,pz,E = some_four_vector[0--3],
Definition at line 236 of file PseudoJet.hh. References _E, _finish_init(), _px, _py, _pz, set_cluster_hist_index(), and set_user_index(). 00236 { 00237 00238 _px = some_four_vector[0]; 00239 _py = some_four_vector[1]; 00240 _pz = some_four_vector[2]; 00241 _E = some_four_vector[3]; 00242 this->_finish_init(); 00243 // some default values for these two indices 00244 set_cluster_hist_index(-1); 00245 set_user_index(-1); 00246 }
|
|
shortcut for converting siscone Cmomentum into PseudoJet
Definition at line 51 of file SISConePlugin.cc. References PseudoJet(). 00051 { 00052 (*this) = PseudoJet(four_vector.px,four_vector.py,four_vector.pz, 00053 four_vector.E); 00054 }
|
|
do standard end of initialisation
Definition at line 64 of file PseudoJet.cc. References _E, _kt2, _phi, _pz, _rap, m2(), fastjet::MaxRap, px(), py(), and fastjet::twopi. Referenced by boost(), operator+=(), operator-=(), PseudoJet(), and unboost(). 00064 { 00065 _kt2 = this->px()*this->px() + this->py()*this->py(); 00066 if (_kt2 == 0.0) { 00067 _phi = 0.0; } 00068 else { 00069 _phi = atan2(this->py(),this->px()); 00070 } 00071 if (_phi < 0.0) {_phi += twopi;} 00072 if (_phi >= twopi) {_phi -= twopi;} // can happen if phi=-|eps<1e-15|? 00073 if (this->E() == abs(this->pz()) && _kt2 == 0) { 00074 // Point has infinite rapidity -- convert that into a very large 00075 // number, but in such a way that different 0-pt momenta will have 00076 // different rapidities (so as to lift the degeneracy between 00077 // them) [this can be relevant at parton-level] 00078 double MaxRapHere = MaxRap + abs(this->pz()); 00079 if (this->pz() >= 0.0) {_rap = MaxRapHere;} else {_rap = -MaxRapHere;} 00080 } else { 00081 // get the rapidity in a way that's modestly insensitive to roundoff 00082 // error when things pz,E are large (actually the best we can do without 00083 // explicit knowledge of mass) 00084 double effective_m2 = max(0.0,m2()); // force non tachyonic mass 00085 double E_plus_pz = _E + abs(_pz); // the safer of p+, p- 00086 // p+/p- = (p+ p-) / (p-)^2 = (kt^2+m^2)/(p-)^2 00087 _rap = 0.5*log((_kt2 + effective_m2)/(E_plus_pz*E_plus_pz)); 00088 if (_pz > 0) {_rap = - _rap;} 00089 } 00090 00092 //if (this->E() != abs(this->pz())) { 00093 // _rap = 0.5*log((this->E() + this->pz())/(this->E() - this->pz())); 00094 // } else { 00095 // // Overlapping points can give problems. Let's lift the degeneracy 00096 // // in case of multiple 0-pT points (can be found at parton-level) 00097 // double MaxRapHere = MaxRap + abs(this->pz()); 00098 // if (this->pz() >= 0.0) {_rap = MaxRapHere;} else {_rap = -MaxRapHere;} 00099 //} 00100 }
|
|
returns distance between this jet and the beam
Definition at line 161 of file PseudoJet.hh. 00161 {return _kt2;};
|
|
transform this jet (given in the rest frame of prest) into a jet in the lab frame;
Definition at line 222 of file PseudoJet.cc. References _E, _finish_init(), _px, _py, _pz, E(), m(), px(), py(), and pz(). 00222 { 00223 00224 if (prest.px() == 0.0 && prest.py() == 0.0 && prest.pz() == 0.0) 00225 return *this; 00226 00227 double m = prest.m(); 00228 assert(m != 0); 00229 00230 double pf4 = ( px()*prest.px() + py()*prest.py() 00231 + pz()*prest.pz() + E()*prest.E() )/m; 00232 double fn = (pf4 + E()) / (prest.E() + m); 00233 _px += fn*prest.px(); 00234 _py += fn*prest.py(); 00235 _pz += fn*prest.pz(); 00236 _E = pf4; 00237 00238 _finish_init(); // we need to recalculate phi,rap,kt2 00239 return *this; 00240 }
|
|
return the cluster_hist_index, intended to be used by clustering routines.
Definition at line 123 of file PseudoJet.hh. Referenced by fastjet::ClusterSequence::add_constituents(), fastjet::ClusterSequenceActiveArea::area(), fastjet::ClusterSequenceActiveAreaExplicitGhosts::area_4vector(), fastjet::ClusterSequenceActiveArea::area_4vector(), fastjet::ClusterSequenceActiveArea::area_error(), and fastjet::ClusterSequenceActiveAreaExplicitGhosts::is_pure_ghost(). 00123 {return _cluster_hist_index;};
|
|
alternative name for cluster_hist_index() [perhaps more meaningful]
Definition at line 128 of file PseudoJet.hh. 00128 { 00129 return cluster_hist_index();};
|
|
Definition at line 68 of file PseudoJet.hh. Referenced by operator()(). 00068 {return _E;}; // like CLHEP
|
|
Definition at line 67 of file PseudoJet.hh. Referenced by fastjet::ClusterSequenceActiveArea::_transfer_areas(), boost(), fastjet::have_same_momentum(), fastjet::operator+(), fastjet::operator-(), fastjet::JetDefinition::DefaultRecombiner::preprocess(), print_jets(), fastjet::SISConePlugin::run_clustering(), fastjet::CDFJetCluPlugin::run_clustering(), and unboost(). 00067 {return _E;};
|
|
return a valarray containing the four-momentum (components 0-2 are 3-mom, component 3 is energy).
Definition at line 105 of file PseudoJet.cc. References _E, _px, _py, and _pz. 00105 { 00106 valarray<double> mom(4); 00107 mom[0] = _px; 00108 mom[1] = _py; 00109 mom[2] = _pz; 00110 mom[3] = _E ; 00111 return mom; 00112 }
|
|
returns the squared transverse momentum
Definition at line 91 of file PseudoJet.hh. Referenced by fastjet::ClusterSequence::jet_scale_for_algorithm(). 00091 {return _kt2;};
|
|
returns kt distance (R=1) between this jet and another
Definition at line 282 of file PseudoJet.cc. References _kt2, _phi, _rap, fastjet::pi, and fastjet::twopi. 00282 { 00283 //double distance = min(this->kt2(), other.kt2()); 00284 double distance = min(_kt2, other._kt2); 00285 double dphi = abs(_phi - other._phi); 00286 if (dphi > pi) {dphi = twopi - dphi;} 00287 double drap = _rap - other._rap; 00288 distance = distance * (dphi*dphi + drap*drap); 00289 return distance; 00290 }
|
|
returns the invariant mass (If m2() is negative then -sqrt(-m2()) is returned, as in CLHEP)
Definition at line 258 of file PseudoJet.hh. References m2(). Referenced by boost(), and unboost(). 00258 { 00259 double mm = m2(); 00260 return mm < 0.0 ? -std::sqrt(-mm) : std::sqrt(mm); 00261 }
|
|
returns the squared invariant mass // like CLHEP
Definition at line 97 of file PseudoJet.hh. Referenced by _finish_init(), and m().
|
|
returns the transverse mass = sqrt(kt^2+m^2)
Definition at line 101 of file PseudoJet.hh. 00101 {return sqrt(std::abs(mperp2()));};
|
|
returns the squared transverse mass = kt^2+m^2
Definition at line 99 of file PseudoJet.hh.
|
|
multiply the jet's momentum by the coefficient
Definition at line 179 of file PseudoJet.cc. References _E, _kt2, _px, _py, and _pz. 00179 { 00180 _px *= coeff; 00181 _py *= coeff; 00182 _pz *= coeff; 00183 _E *= coeff; 00184 _kt2*= coeff*coeff; 00185 // phi and rap are unchanged 00186 }
|
|
returns component i, where X==0, Y==1, Z==2, E==3
Definition at line 117 of file PseudoJet.cc. References e(), px(), py(), pz(), T, X, Y, and Z. 00117 { 00118 switch(i) { 00119 case X: 00120 return px(); 00121 case Y: 00122 return py(); 00123 case Z: 00124 return pz(); 00125 case T: 00126 return e(); 00127 default: 00128 ostringstream err; 00129 err << "PseudoJet subscripting: bad index (" << i << ")"; 00130 throw Error(err.str()); 00131 } 00132 return 0.; 00133 }
|
|
add the other jet's momentum to this jet
Definition at line 197 of file PseudoJet.cc. References _E, _finish_init(), _px, _py, and _pz. 00197 { 00198 _px += other_jet._px; 00199 _py += other_jet._py; 00200 _pz += other_jet._pz; 00201 _E += other_jet._E ; 00202 _finish_init(); // we need to recalculate phi,rap,kt2 00203 }
|
|
subtract the other jet's momentum from this jet
Definition at line 208 of file PseudoJet.cc. References _E, _finish_init(), _px, _py, and _pz. 00208 { 00209 _px -= other_jet._px; 00210 _py -= other_jet._py; 00211 _pz -= other_jet._pz; 00212 _E -= other_jet._E ; 00213 _finish_init(); // we need to recalculate phi,rap,kt2 00214 }
|
|
divide the jet's momentum by the coefficient
Definition at line 190 of file PseudoJet.cc. 00190 { 00191 (*this) *= 1.0/coeff; 00192 }
|
|
returns component i, where X==0, Y==1, Z==2, E==3
Definition at line 108 of file PseudoJet.hh. 00108 { return (*this)(i); }; // this too
|
|
returns the scalar transverse momentum
Definition at line 95 of file PseudoJet.hh. Referenced by fastjet::ClusterSequenceActiveArea::_transfer_areas(), and fastjet::JetDefinition::DefaultRecombiner::recombine(). 00095 {return sqrt(_kt2);}; // like CLHEP
|
|
returns the squared transverse momentum
Definition at line 93 of file PseudoJet.hh. Referenced by fastjet::ClusterSequenceActiveArea::_transfer_areas(), fastjet::ClusterSequence::inclusive_jets(), fastjet::JetDefinition::DefaultRecombiner::preprocess(), print_jets(), and fastjet::JetDefinition::DefaultRecombiner::recombine(). 00093 {return _kt2;}; // like CLHEP
|
|
returns phi (in the range 0..2pi)
Definition at line 74 of file PseudoJet.hh. Referenced by fastjet::JetDefinition::DefaultRecombiner::recombine(). 00074 {return phi_02pi();};
|
|
returns phi in the range 0..2pi
Definition at line 81 of file PseudoJet.hh. 00081 {return _phi;};
|
|
returns phi in the range -pi..pi
Definition at line 77 of file PseudoJet.hh. References fastjet::pi, and twopi.
|
|
returns squared cylinder (rap-phi) distance between this jet and another
Definition at line 295 of file PseudoJet.cc. References _phi, _rap, fastjet::pi, and fastjet::twopi. 00295 { 00296 double dphi = abs(_phi - other._phi); 00297 if (dphi > pi) {dphi = twopi - dphi;} 00298 double drap = _rap - other._rap; 00299 return (dphi*dphi + drap*drap); 00300 }
|
|
Definition at line 69 of file PseudoJet.hh. Referenced by _finish_init(), fastjet::ClusterSequenceActiveArea::_transfer_areas(), boost(), fastjet::have_same_momentum(), operator()(), fastjet::operator+(), fastjet::operator-(), fastjet::JetDefinition::DefaultRecombiner::preprocess(), fastjet::SISConePlugin::run_clustering(), fastjet::CDFJetCluPlugin::run_clustering(), and unboost(). 00069 {return _px;};
|
|
Definition at line 70 of file PseudoJet.hh. Referenced by _finish_init(), fastjet::ClusterSequenceActiveArea::_transfer_areas(), boost(), fastjet::have_same_momentum(), operator()(), fastjet::operator+(), fastjet::operator-(), fastjet::JetDefinition::DefaultRecombiner::preprocess(), fastjet::SISConePlugin::run_clustering(), fastjet::CDFJetCluPlugin::run_clustering(), and unboost(). 00070 {return _py;};
|
|
Definition at line 71 of file PseudoJet.hh. Referenced by fastjet::ClusterSequenceActiveArea::_transfer_areas(), boost(), fastjet::have_same_momentum(), operator()(), fastjet::operator+(), fastjet::operator-(), fastjet::JetDefinition::DefaultRecombiner::preprocess(), fastjet::SISConePlugin::run_clustering(), fastjet::CDFJetCluPlugin::run_clustering(), and unboost(). 00071 {return _pz;};
|
|
returns the rapidity or some large value when the rapidity is infinite
Definition at line 85 of file PseudoJet.hh. Referenced by main(), and fastjet::JetDefinition::DefaultRecombiner::recombine(). 00085 {return _rap;};
|
|
the same as rap()
Definition at line 88 of file PseudoJet.hh. 00088 {return _rap;}; // like CLHEP
|
|
set the cluster_hist_index, intended to be used by clustering routines.
Definition at line 125 of file PseudoJet.hh. Referenced by PseudoJet(). 00125 {_cluster_hist_index = index;};
|
|
alternative name for set_cluster_hist_index(. ..) [perhaps more meaningful] Definition at line 132 of file PseudoJet.hh. 00132 { 00133 set_cluster_hist_index(index);};
|
|
set the user_index, intended to allow the user to "add" information
Definition at line 139 of file PseudoJet.hh. Referenced by main(), fastjet::JetDefinition::DefaultRecombiner::preprocess(), PseudoJet(), fastjet::JetDefinition::DefaultRecombiner::recombine(), and fastjet::SISConePlugin::run_clustering(). 00139 {_user_index = index;};
|
|
returns squared cylinder (rap-phi) distance between this jet and another
Definition at line 152 of file PseudoJet.hh. 00152 { 00153 return plain_distance(other);};
|
|
transform this jet (given in the rest frame of prest) into a jet in the lab frame;
Definition at line 249 of file PseudoJet.cc. References _E, _finish_init(), _px, _py, _pz, E(), m(), px(), py(), and pz(). 00249 { 00250 00251 if (prest.px() == 0.0 && prest.py() == 0.0 && prest.pz() == 0.0) 00252 return *this; 00253 00254 double m = prest.m(); 00255 assert(m != 0); 00256 00257 double pf4 = ( -px()*prest.px() - py()*prest.py() 00258 - pz()*prest.pz() + E()*prest.E() )/m; 00259 double fn = (pf4 + E()) / (prest.E() + m); 00260 _px -= fn*prest.px(); 00261 _py -= fn*prest.py(); 00262 _pz -= fn*prest.pz(); 00263 _E = pf4; 00264 00265 _finish_init(); // we need to recalculate phi,rap,kt2 00266 return *this; 00267 }
|
|
return the user_index, intended to allow the user to "add" information
Definition at line 137 of file PseudoJet.hh. Referenced by fastjet::JetDefinition::DefaultRecombiner::preprocess(). 00137 {return _user_index;};
|
|
Definition at line 173 of file PseudoJet.hh. |
|
Definition at line 171 of file PseudoJet.hh. Referenced by _finish_init(), boost(), four_mom(), operator *=(), operator+=(), operator-=(), PseudoJet(), and unboost(). |
|
Definition at line 172 of file PseudoJet.hh. Referenced by _finish_init(), kt_distance(), and operator *=(). |
|
Definition at line 172 of file PseudoJet.hh. Referenced by _finish_init(), kt_distance(), and plain_distance(). |
|
Definition at line 171 of file PseudoJet.hh. Referenced by boost(), four_mom(), operator *=(), operator+=(), operator-=(), PseudoJet(), and unboost(). |
|
Definition at line 171 of file PseudoJet.hh. Referenced by boost(), four_mom(), operator *=(), operator+=(), operator-=(), PseudoJet(), and unboost(). |
|
Definition at line 171 of file PseudoJet.hh. Referenced by _finish_init(), boost(), four_mom(), operator *=(), operator+=(), operator-=(), PseudoJet(), and unboost(). |
|
Definition at line 172 of file PseudoJet.hh. Referenced by _finish_init(), kt_distance(), and plain_distance(). |
|
Definition at line 173 of file PseudoJet.hh. |