Particle.hh
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #ifndef RIVET_Particle_HH 00003 #define RIVET_Particle_HH 00004 00005 #include "Rivet/Rivet.hh" 00006 #include "Rivet/Particle.fhh" 00007 #include "Rivet/ParticleBase.hh" 00008 #include "Rivet/ParticleName.hh" 00009 #include "Rivet/Math/Vectors.hh" 00010 #include "Rivet/Tools/Logging.hh" 00011 00012 namespace Rivet { 00013 00014 00015 /// Representation of particles from a HepMC::GenEvent. 00016 class Particle : public ParticleBase { 00017 public: 00018 00019 /// Default constructor. 00020 /// @deprecated A particle without info is useless. This only exists to keep STL containers happy. 00021 Particle() 00022 : ParticleBase(), 00023 _original(0), _id(0), _momentum() 00024 { } 00025 00026 /// Constructor without GenParticle. 00027 Particle(PdgId pid, const FourMomentum& mom) 00028 : ParticleBase(), 00029 _original(0), _id(pid), _momentum(mom) 00030 { } 00031 00032 /// Constructor from a HepMC GenParticle. 00033 Particle(const GenParticle& gp) 00034 : ParticleBase(), 00035 _original(&gp), _id(gp.pdg_id()), 00036 _momentum(gp.momentum()) 00037 { } 00038 00039 00040 public: 00041 00042 /// Get a const reference to the original GenParticle. 00043 const GenParticle& genParticle() const { 00044 assert(_original); 00045 return *_original; 00046 } 00047 00048 00049 /// Check if the particle corresponds to a GenParticle. 00050 bool hasGenParticle() const { 00051 return bool(_original); 00052 } 00053 00054 00055 /// The PDG ID code for this Particle. 00056 PdgId pdgId() const { 00057 return _id; 00058 } 00059 00060 00061 /// Set the momentum of this Particle. 00062 Particle& setMomentum(const FourMomentum& momentum) { 00063 _momentum = momentum; 00064 return *this; 00065 } 00066 00067 /// The momentum of this Particle. 00068 const FourMomentum& momentum() const { 00069 return _momentum; 00070 } 00071 00072 /// The energy of this Particle. 00073 double energy() const { 00074 return momentum().E(); 00075 } 00076 00077 /// The mass of this Particle. 00078 double mass() const { 00079 return momentum().mass(); 00080 } 00081 00082 // /// The charge of this Particle. 00083 // double charge() const { 00084 // return PID::charge(*this); 00085 // } 00086 00087 // /// Three times the charge of this Particle (i.e. integer multiple of smallest quark charge). 00088 // int threeCharge() const { 00089 // return PID::threeCharge(*this); 00090 // } 00091 00092 00093 /// Check whether a given PID is found in the GenParticle's ancestor list 00094 bool hasAncestor(PdgId pdg_id) const; 00095 00096 00097 private: 00098 00099 /// A pointer to the original GenParticle from which this Particle is projected. 00100 const GenParticle* _original; 00101 00102 /// The PDG ID code for this Particle. 00103 PdgId _id; 00104 00105 /// The momentum of this projection of the Particle. 00106 FourMomentum _momentum; 00107 }; 00108 00109 00110 /// @name String representation 00111 //@{ 00112 00113 /// Print a ParticlePair as a string. 00114 inline std::string toString(const ParticlePair& pair) { 00115 stringstream out; 00116 out << "[" 00117 << toParticleName(pair.first.pdgId()) << " @ " 00118 << pair.first.momentum().E()/GeV << " GeV, " 00119 << toParticleName(pair.second.pdgId()) << " @ " 00120 << pair.second.momentum().E()/GeV << " GeV]"; 00121 return out.str(); 00122 } 00123 00124 /// Allow ParticlePair to be passed to an ostream. 00125 inline std::ostream& operator<<(std::ostream& os, const ParticlePair& pp) { 00126 os << toString(pp); 00127 return os; 00128 } 00129 00130 //@} 00131 00132 00133 /// @name Comparison functors 00134 //@{ 00135 /// Sort by descending transverse momentum, \f$ p_\perp \f$ 00136 inline bool cmpParticleByPt(const Particle& a, const Particle& b) { 00137 return a.momentum().pT() > b.momentum().pT(); 00138 } 00139 /// Sort by ascending transverse momentum, \f$ p_\perp \f$ 00140 inline bool cmpParticleByAscPt(const Particle& a, const Particle& b) { 00141 return a.momentum().pT() < b.momentum().pT(); 00142 } 00143 /// Sort by descending momentum, \f$ p \f$ 00144 inline bool cmpParticleByP(const Particle& a, const Particle& b) { 00145 return a.momentum().vector3().mod() > b.momentum().vector3().mod(); 00146 } 00147 /// Sort by ascending momentum, \f$ p \f$ 00148 inline bool cmpParticleByAscP(const Particle& a, const Particle& b) { 00149 return a.momentum().vector3().mod() < b.momentum().vector3().mod(); 00150 } 00151 /// Sort by descending transverse energy, \f$ E_\perp \f$ 00152 inline bool cmpParticleByEt(const Particle& a, const Particle& b) { 00153 return a.momentum().Et() > b.momentum().Et(); 00154 } 00155 /// Sort by ascending transverse energy, \f$ E_\perp \f$ 00156 inline bool cmpParticleByAscEt(const Particle& a, const Particle& b) { 00157 return a.momentum().Et() < b.momentum().Et(); 00158 } 00159 /// Sort by descending energy, \f$ E \f$ 00160 inline bool cmpParticleByE(const Particle& a, const Particle& b) { 00161 return a.momentum().E() > b.momentum().E(); 00162 } 00163 /// Sort by ascending energy, \f$ E \f$ 00164 inline bool cmpParticleByAscE(const Particle& a, const Particle& b) { 00165 return a.momentum().E() < b.momentum().E(); 00166 } 00167 /// Sort by descending pseudorapidity, \f$ \eta \f$ 00168 inline bool cmpParticleByDescPseudorapidity(const Particle& a, const Particle& b) { 00169 return a.momentum().pseudorapidity() > b.momentum().pseudorapidity(); 00170 } 00171 /// Sort by ascending pseudorapidity, \f$ \eta \f$ 00172 inline bool cmpParticleByAscPseudorapidity(const Particle& a, const Particle& b) { 00173 return a.momentum().pseudorapidity() < b.momentum().pseudorapidity(); 00174 } 00175 /// Sort by descending absolute pseudorapidity, \f$ |\eta| \f$ 00176 inline bool cmpParticleByDescAbsPseudorapidity(const Particle& a, const Particle& b) { 00177 return fabs(a.momentum().pseudorapidity()) > fabs(b.momentum().pseudorapidity()); 00178 } 00179 /// Sort by ascending absolute pseudorapidity, \f$ |\eta| \f$ 00180 inline bool cmpParticleByAscAbsPseudorapidity(const Particle& a, const Particle& b) { 00181 return fabs(a.momentum().pseudorapidity()) < fabs(b.momentum().pseudorapidity()); 00182 } 00183 /// Sort by descending rapidity, \f$ y \f$ 00184 inline bool cmpParticleByDescRapidity(const Particle& a, const Particle& b) { 00185 return a.momentum().rapidity() > b.momentum().rapidity(); 00186 } 00187 /// Sort by ascending rapidity, \f$ y \f$ 00188 inline bool cmpParticleByAscRapidity(const Particle& a, const Particle& b) { 00189 return a.momentum().rapidity() < b.momentum().rapidity(); 00190 } 00191 /// Sort by descending absolute rapidity, \f$ |y| \f$ 00192 inline bool cmpParticleByDescAbsRapidity(const Particle& a, const Particle& b) { 00193 return fabs(a.momentum().rapidity()) > fabs(b.momentum().rapidity()); 00194 } 00195 /// Sort by ascending absolute rapidity, \f$ |y| \f$ 00196 inline bool cmpParticleByAscAbsRapidity(const Particle& a, const Particle& b) { 00197 return fabs(a.momentum().rapidity()) < fabs(b.momentum().rapidity()); 00198 } 00199 //@} 00200 00201 inline double deltaR(const Particle& p1, const Particle& p2, 00202 RapScheme scheme = PSEUDORAPIDITY) { 00203 return deltaR(p1.momentum(), p2.momentum(), scheme); 00204 } 00205 00206 inline double deltaR(const Particle& p, const FourMomentum& v, 00207 RapScheme scheme = PSEUDORAPIDITY) { 00208 return deltaR(p.momentum(), v, scheme); 00209 } 00210 00211 inline double deltaR(const Particle& p, const FourVector& v, 00212 RapScheme scheme = PSEUDORAPIDITY) { 00213 return deltaR(p.momentum(), v, scheme); 00214 } 00215 00216 inline double deltaR(const Particle& p, const Vector3& v) { 00217 return deltaR(p.momentum(), v); 00218 } 00219 00220 inline double deltaR(const Particle& p, double eta, double phi) { 00221 return deltaR(p.momentum(), eta, phi); 00222 } 00223 00224 inline double deltaR(const FourMomentum& v, const Particle& p, 00225 RapScheme scheme = PSEUDORAPIDITY) { 00226 return deltaR(v, p.momentum(), scheme); 00227 } 00228 00229 inline double deltaR(const FourVector& v, const Particle& p, 00230 RapScheme scheme = PSEUDORAPIDITY) { 00231 return deltaR(v, p.momentum(), scheme); 00232 } 00233 00234 inline double deltaR(const Vector3& v, const Particle& p) { 00235 return deltaR(v, p.momentum()); 00236 } 00237 00238 inline double deltaR(double eta, double phi, const Particle& p) { 00239 return deltaR(eta, phi, p.momentum()); 00240 } 00241 00242 00243 inline double deltaPhi(const Particle& p1, const Particle& p2) { 00244 return deltaPhi(p1.momentum(), p2.momentum()); 00245 } 00246 00247 inline double deltaPhi(const Particle& p, const FourMomentum& v) { 00248 return deltaPhi(p.momentum(), v); 00249 } 00250 00251 inline double deltaPhi(const Particle& p, const FourVector& v) { 00252 return deltaPhi(p.momentum(), v); 00253 } 00254 00255 inline double deltaPhi(const Particle& p, const Vector3& v) { 00256 return deltaPhi(p.momentum(), v); 00257 } 00258 00259 inline double deltaPhi(const Particle& p, double phi) { 00260 return deltaPhi(p.momentum(), phi); 00261 } 00262 00263 inline double deltaPhi(const FourMomentum& v, const Particle& p) { 00264 return deltaPhi(v, p.momentum()); 00265 } 00266 00267 inline double deltaPhi(const FourVector& v, const Particle& p) { 00268 return deltaPhi(v, p.momentum()); 00269 } 00270 00271 inline double deltaPhi(const Vector3& v, const Particle& p) { 00272 return deltaPhi(v, p.momentum()); 00273 } 00274 00275 inline double deltaPhi(double phi, const Particle& p) { 00276 return deltaPhi(phi, p.momentum()); 00277 } 00278 00279 00280 inline double deltaEta(const Particle& p1, const Particle& p2) { 00281 return deltaEta(p1.momentum(), p2.momentum()); 00282 } 00283 00284 inline double deltaEta(const Particle& p, const FourMomentum& v) { 00285 return deltaEta(p.momentum(), v); 00286 } 00287 00288 inline double deltaEta(const Particle& p, const FourVector& v) { 00289 return deltaEta(p.momentum(), v); 00290 } 00291 00292 inline double deltaEta(const Particle& p, const Vector3& v) { 00293 return deltaEta(p.momentum(), v); 00294 } 00295 00296 inline double deltaEta(const Particle& p, double eta) { 00297 return deltaEta(p.momentum(), eta); 00298 } 00299 00300 inline double deltaEta(const FourMomentum& v, const Particle& p) { 00301 return deltaEta(v, p.momentum()); 00302 } 00303 00304 inline double deltaEta(const FourVector& v, const Particle& p) { 00305 return deltaEta(v, p.momentum()); 00306 } 00307 00308 inline double deltaEta(const Vector3& v, const Particle& p) { 00309 return deltaEta(v, p.momentum()); 00310 } 00311 00312 inline double deltaEta(double eta, const Particle& p) { 00313 return deltaEta(eta, p.momentum()); 00314 } 00315 00316 } 00317 00318 #endif Generated on Fri Dec 21 2012 15:03:41 for The Rivet MC analysis system by ![]() |