#include <valarray>
#include <vector>
#include <cassert>
#include <cmath>
#include <iostream>
#include "fastjet/internal/numconsts.hh"
Go to the source code of this file.
Classes | |
class | PseudoJet |
Class to contain pseudojets, including minimal information of use to to jet-clustering routines. More... | |
class | IndexedSortHelper |
a class that helps us carry out indexed sorting. More... | |
Functions | |
PseudoJet | operator+ (const PseudoJet &, const PseudoJet &) |
PseudoJet | operator- (const PseudoJet &, const PseudoJet &) |
PseudoJet | operator* (double, const PseudoJet &) |
PseudoJet | operator* (const PseudoJet &, double) |
PseudoJet | operator/ (const PseudoJet &, double) |
double | dot_product (const PseudoJet &a, const PseudoJet &b) |
bool | have_same_momentum (const PseudoJet &, const PseudoJet &) |
returns true if the momenta of the two input jets are identical | |
PseudoJet | PtYPhiM (double pt, double y, double phi, double m=0.0) |
return a pseudojet with the given pt, y, phi and mass | |
std::vector< PseudoJet > | sorted_by_pt (const std::vector< PseudoJet > &jets) |
return a vector of jets sorted into decreasing transverse momentum | |
std::vector< PseudoJet > | sorted_by_rapidity (const std::vector< PseudoJet > &jets) |
return a vector of jets sorted into increasing rapidity | |
std::vector< PseudoJet > | sorted_by_E (const std::vector< PseudoJet > &jets) |
return a vector of jets sorted into decreasing energy | |
std::vector< PseudoJet > | sorted_by_pz (const std::vector< PseudoJet > &jets) |
return a vector of jets sorted into increasing pz | |
void | sort_indices (std::vector< int > &indices, const std::vector< double > &values) |
sort the indices so that values[indices[0->n-1]] is sorted into increasing order | |
template<class T > | |
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't actually touch the values vector in the process). | |
Variables | |
FASTJET_BEGIN_NAMESPACE const double | MaxRap = 1e5 |
Used to protect against parton-level events where pt can be zero for some partons, giving rapidity=infinity. |
Definition at line 228 of file PseudoJet.hh.
References PseudoJet::E(), PseudoJet::px(), PseudoJet::py(), and PseudoJet::pz().
returns true if the momenta of the two input jets are identical
Definition at line 284 of file PseudoJet.cc.
References PseudoJet::E(), PseudoJet::px(), PseudoJet::py(), and PseudoJet::pz().
Referenced by fastjet::SISConeSphericalPlugin::run_clustering(), and fastjet::SISConePlugin::run_clustering().
std::vector<T> objects_sorted_by_values | ( | const std::vector< T > & | objects, | |
const std::vector< double > & | values | |||
) | [inline] |
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't actually touch the values vector in the process).
Definition at line 178 of file PseudoJet.cc.
Definition at line 168 of file PseudoJet.cc.
00168 { 00169 //return PseudoJet(coeff*jet.four_mom()); 00170 // the following code is hopefully more efficient 00171 PseudoJet coeff_times_jet(jet); 00172 coeff_times_jet *= coeff; 00173 return coeff_times_jet; 00174 }
Definition at line 148 of file PseudoJet.cc.
References PseudoJet::E(), PseudoJet::px(), PseudoJet::py(), and PseudoJet::pz().
00148 { 00149 //return PseudoJet(jet1.four_mom()+jet2.four_mom()); 00150 return PseudoJet(jet1.px()+jet2.px(), 00151 jet1.py()+jet2.py(), 00152 jet1.pz()+jet2.pz(), 00153 jet1.E() +jet2.E() ); 00154 }
Definition at line 158 of file PseudoJet.cc.
References PseudoJet::E(), PseudoJet::px(), PseudoJet::py(), and PseudoJet::pz().
00158 { 00159 //return PseudoJet(jet1.four_mom()-jet2.four_mom()); 00160 return PseudoJet(jet1.px()-jet2.px(), 00161 jet1.py()-jet2.py(), 00162 jet1.pz()-jet2.pz(), 00163 jet1.E() -jet2.E() ); 00164 }
Definition at line 184 of file PseudoJet.cc.
PseudoJet PtYPhiM | ( | double | pt, | |
double | y, | |||
double | phi, | |||
double | m = 0.0 | |||
) |
void sort_indices | ( | std::vector< int > & | indices, | |
const std::vector< double > & | values | |||
) |
sort the indices so that values[indices[0->n-1]] is sorted into increasing order
return a vector of jets sorted into decreasing energy
return a vector of jets sorted into decreasing transverse momentum
return a vector of jets sorted into increasing pz
return a vector of jets sorted into increasing rapidity
FASTJET_BEGIN_NAMESPACE const double MaxRap = 1e5 |
Used to protect against parton-level events where pt can be zero for some partons, giving rapidity=infinity.
KtJet fails in those cases.
Definition at line 48 of file PseudoJet.hh.
Referenced by PseudoJet::_finish_init(), and PseudoJet::pseudorapidity().