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 #include "Rivet/Tools/ParticleIdUtils.hh" 00012 00013 namespace Rivet { 00014 00015 00016 /// Representation of particles from a HepMC::GenEvent. 00017 class Particle : public ParticleBase { 00018 public: 00019 00020 /// Default constructor. 00021 /// @deprecated A particle without info is useless. This only exists to keep STL containers happy. 00022 Particle() 00023 : ParticleBase(), 00024 _original(0), _id(0), _momentum() 00025 { } 00026 00027 /// Constructor without GenParticle. 00028 Particle(PdgId pid, const FourMomentum& mom) 00029 : ParticleBase(), 00030 _original(0), _id(pid), _momentum(mom) 00031 { } 00032 00033 /// Constructor from a HepMC GenParticle. 00034 Particle(const GenParticle& gp) 00035 : ParticleBase(), 00036 _original(&gp), _id(gp.pdg_id()), 00037 _momentum(gp.momentum()) 00038 { } 00039 00040 /// Constructor from a HepMC GenParticle pointer. 00041 Particle(const GenParticle* gp) 00042 : ParticleBase(), 00043 _original(gp), _id(gp->pdg_id()), 00044 _momentum(gp->momentum()) 00045 { } 00046 00047 00048 public: 00049 00050 /// @name Basic particle specific properties 00051 //@{ 00052 00053 /// Get a const reference to the original GenParticle. 00054 const GenParticle* genParticle() const { 00055 return _original; 00056 } 00057 00058 /// Cast operator for conversion to GenParticle* 00059 operator const GenParticle* () const { return genParticle(); } 00060 00061 /// The momentum. 00062 const FourMomentum& momentum() const { 00063 return _momentum; 00064 } 00065 /// Set the momentum. 00066 Particle& setMomentum(const FourMomentum& momentum) { 00067 _momentum = momentum; 00068 return *this; 00069 } 00070 00071 //@} 00072 00073 00074 /// @name Particle ID code accessors 00075 //@{ 00076 00077 /// This Particle's PDG ID code. 00078 PdgId pid() const { return _id; } 00079 /// Absolute value of the PDG ID code. 00080 PdgId abspid() const { return _id; } 00081 /// This Particle's PDG ID code (alias). 00082 /// @deprecatedThe pid/abspid form is nicer (don't need to remember 00083 /// lower/upper case, doesn't visually stick out in the code, etc. ...) 00084 PdgId pdgId() const { return _id; } 00085 00086 //@} 00087 00088 00089 /// @name Particle species properties 00090 //@{ 00091 00092 /// The charge of this Particle. 00093 double charge() const { 00094 return PID::charge(pdgId()); 00095 } 00096 /// Three times the charge of this Particle (i.e. integer multiple of smallest quark charge). 00097 int threeCharge() const { 00098 return PID::threeCharge(pdgId()); 00099 } 00100 00101 /// Is this a hadron? 00102 bool isHadron() const { return PID::isHadron(pdgId()); } 00103 00104 /// Is this a meson? 00105 bool isMeson() const { return PID::isMeson(pdgId()); } 00106 00107 /// Is this a baryon? 00108 bool isBaryon() const { return PID::isBaryon(pdgId()); } 00109 00110 /// Is this a lepton? 00111 bool isLepton() const { return PID::isLepton(pdgId()); } 00112 00113 /// Is this a neutrino? 00114 bool isNeutrino() const { return PID::isNeutrino(pdgId()); } 00115 00116 /// Does this (hadron) contain a b quark? 00117 bool hasBottom() const { return PID::hasBottom(pdgId()); } 00118 00119 /// Does this (hadron) contain a c quark? 00120 bool hasCharm() const { return PID::hasCharm(pdgId()); } 00121 00122 // /// Does this (hadron) contain an s quark? 00123 // bool hasStrange() const { return PID::hasStrange(pdgId()); } 00124 00125 //@} 00126 00127 00128 /// @name Ancestry properties 00129 //@{ 00130 00131 /// Check whether a given PID is found in the GenParticle's ancestor list 00132 /// 00133 /// @note This question is valid in MC, but may not be answerable 00134 /// experimentally -- use this function with care when replicating 00135 /// experimental analyses! 00136 bool hasAncestor(PdgId pdg_id) const; 00137 00138 /// @brief Determine whether the particle is from a hadron or tau decay 00139 /// 00140 /// Specifically, walk up the ancestor chain until a status 2 hadron or 00141 /// tau is found, if at all. 00142 /// 00143 /// @note This question is valid in MC, but may not be perfectly answerable 00144 /// experimentally -- use this function with care when replicating 00145 /// experimental analyses! 00146 bool fromDecay() const; 00147 00148 /// @brief Determine whether the particle is from a b-hadron decay 00149 /// 00150 /// @note This question is valid in MC, but may not be perfectly answerable 00151 /// experimentally -- use this function with care when replicating 00152 /// experimental analyses! 00153 bool fromBottom() const; 00154 00155 /// @brief Determine whether the particle is from a c-hadron decay 00156 /// 00157 /// @note If a hadron contains b and c quarks it is considered a bottom 00158 /// hadron and NOT a charm hadron. 00159 /// 00160 /// @note This question is valid in MC, but may not be perfectly answerable 00161 /// experimentally -- use this function with care when replicating 00162 /// experimental analyses! 00163 bool fromCharm() const; 00164 00165 // /// @brief Determine whether the particle is from a s-hadron decay 00166 // /// 00167 // /// @note If a hadron contains b or c quarks as well as strange it is 00168 // /// considered a b or c hadron, but NOT a strange hadron. 00169 // /// 00170 // /// @note This question is valid in MC, but may not be perfectly answerable 00171 // /// experimentally -- use this function with care when replicating 00172 // /// experimental analyses! 00173 // bool fromStrange() const; 00174 00175 /// @brief Determine whether the particle is from a tau decay 00176 /// 00177 /// @note This question is valid in MC, but may not be perfectly answerable 00178 /// experimentally -- use this function with care when replicating 00179 /// experimental analyses! 00180 bool fromTau() const; 00181 00182 //@} 00183 00184 00185 private: 00186 00187 /// A pointer to the original GenParticle from which this Particle is projected. 00188 const GenParticle* _original; 00189 00190 /// The PDG ID code for this Particle. 00191 PdgId _id; 00192 00193 /// The momentum of this projection of the Particle. 00194 FourMomentum _momentum; 00195 00196 /// @todo Also store production and decay positions and make them available. 00197 00198 }; 00199 00200 00201 /// @name String representation 00202 //@{ 00203 00204 /// Print a ParticlePair as a string. 00205 inline std::string toString(const ParticlePair& pair) { 00206 stringstream out; 00207 out << "[" 00208 << PID::toParticleName(pair.first.pdgId()) << " @ " 00209 << pair.first.momentum().E()/GeV << " GeV, " 00210 << PID::toParticleName(pair.second.pdgId()) << " @ " 00211 << pair.second.momentum().E()/GeV << " GeV]"; 00212 return out.str(); 00213 } 00214 00215 /// Allow ParticlePair to be passed to an ostream. 00216 inline std::ostream& operator<<(std::ostream& os, const ParticlePair& pp) { 00217 os << toString(pp); 00218 return os; 00219 } 00220 00221 //@} 00222 00223 00224 /// @name deltaR, deltaEta, deltaPhi functions specifically for Particle arguments 00225 //@{ 00226 00227 inline double deltaR(const Particle& p1, const Particle& p2, 00228 RapScheme scheme = PSEUDORAPIDITY) { 00229 return deltaR(p1.momentum(), p2.momentum(), scheme); 00230 } 00231 00232 inline double deltaR(const Particle& p, const FourMomentum& v, 00233 RapScheme scheme = PSEUDORAPIDITY) { 00234 return deltaR(p.momentum(), v, scheme); 00235 } 00236 00237 inline double deltaR(const Particle& p, const FourVector& v, 00238 RapScheme scheme = PSEUDORAPIDITY) { 00239 return deltaR(p.momentum(), v, scheme); 00240 } 00241 00242 inline double deltaR(const Particle& p, const Vector3& v) { 00243 return deltaR(p.momentum(), v); 00244 } 00245 00246 inline double deltaR(const Particle& p, double eta, double phi) { 00247 return deltaR(p.momentum(), eta, phi); 00248 } 00249 00250 inline double deltaR(const FourMomentum& v, const Particle& p, 00251 RapScheme scheme = PSEUDORAPIDITY) { 00252 return deltaR(v, p.momentum(), scheme); 00253 } 00254 00255 inline double deltaR(const FourVector& v, const Particle& p, 00256 RapScheme scheme = PSEUDORAPIDITY) { 00257 return deltaR(v, p.momentum(), scheme); 00258 } 00259 00260 inline double deltaR(const Vector3& v, const Particle& p) { 00261 return deltaR(v, p.momentum()); 00262 } 00263 00264 inline double deltaR(double eta, double phi, const Particle& p) { 00265 return deltaR(eta, phi, p.momentum()); 00266 } 00267 00268 00269 inline double deltaPhi(const Particle& p1, const Particle& p2) { 00270 return deltaPhi(p1.momentum(), p2.momentum()); 00271 } 00272 00273 inline double deltaPhi(const Particle& p, const FourMomentum& v) { 00274 return deltaPhi(p.momentum(), v); 00275 } 00276 00277 inline double deltaPhi(const Particle& p, const FourVector& v) { 00278 return deltaPhi(p.momentum(), v); 00279 } 00280 00281 inline double deltaPhi(const Particle& p, const Vector3& v) { 00282 return deltaPhi(p.momentum(), v); 00283 } 00284 00285 inline double deltaPhi(const Particle& p, double phi) { 00286 return deltaPhi(p.momentum(), phi); 00287 } 00288 00289 inline double deltaPhi(const FourMomentum& v, const Particle& p) { 00290 return deltaPhi(v, p.momentum()); 00291 } 00292 00293 inline double deltaPhi(const FourVector& v, const Particle& p) { 00294 return deltaPhi(v, p.momentum()); 00295 } 00296 00297 inline double deltaPhi(const Vector3& v, const Particle& p) { 00298 return deltaPhi(v, p.momentum()); 00299 } 00300 00301 inline double deltaPhi(double phi, const Particle& p) { 00302 return deltaPhi(phi, p.momentum()); 00303 } 00304 00305 00306 inline double deltaEta(const Particle& p1, const Particle& p2) { 00307 return deltaEta(p1.momentum(), p2.momentum()); 00308 } 00309 00310 inline double deltaEta(const Particle& p, const FourMomentum& v) { 00311 return deltaEta(p.momentum(), v); 00312 } 00313 00314 inline double deltaEta(const Particle& p, const FourVector& v) { 00315 return deltaEta(p.momentum(), v); 00316 } 00317 00318 inline double deltaEta(const Particle& p, const Vector3& v) { 00319 return deltaEta(p.momentum(), v); 00320 } 00321 00322 inline double deltaEta(const Particle& p, double eta) { 00323 return deltaEta(p.momentum(), eta); 00324 } 00325 00326 inline double deltaEta(const FourMomentum& v, const Particle& p) { 00327 return deltaEta(v, p.momentum()); 00328 } 00329 00330 inline double deltaEta(const FourVector& v, const Particle& p) { 00331 return deltaEta(v, p.momentum()); 00332 } 00333 00334 inline double deltaEta(const Vector3& v, const Particle& p) { 00335 return deltaEta(v, p.momentum()); 00336 } 00337 00338 inline double deltaEta(double eta, const Particle& p) { 00339 return deltaEta(eta, p.momentum()); 00340 } 00341 00342 //@} 00343 00344 00345 } 00346 00347 #endif Generated on Thu Feb 6 2014 17:38:46 for The Rivet MC analysis system by 1.7.6.1 |