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() : ParticleBase(),
00022       _original(0), _id(0), _momentum()
00023     { }
00024 
00025     /// Constructor without GenParticle.
00026     Particle(PdgId pid, const FourMomentum& mom) : ParticleBase(),
00027       _original(0), _id(pid), _momentum(mom)
00028     { }
00029 
00030     /// Constructor from a HepMC GenParticle.
00031     Particle(const GenParticle& gp) : ParticleBase(),
00032     _original(&gp), _id(gp.pdg_id()),
00033     _momentum(gp.momentum())
00034     { }
00035 
00036 
00037   public:
00038 
00039     /// Get a const reference to the original GenParticle.
00040     const GenParticle& genParticle() const {
00041       assert(_original);
00042       return *_original;
00043     }
00044 
00045 
00046     /// Check if the particle corresponds to a GenParticle.
00047     bool hasGenParticle() const { 
00048       return bool(_original); 
00049     }
00050 
00051 
00052     /// The PDG ID code for this Particle.
00053     long pdgId() const { 
00054       return _id; 
00055     }
00056 
00057 
00058     /// The momentum of this Particle.
00059     const FourMomentum& momentum() const { 
00060       return _momentum; 
00061     }
00062 
00063 
00064     /// Set the momentum of this Particle.
00065     Particle& setMomentum(const FourMomentum& momentum) { 
00066       _momentum = momentum; 
00067       return *this; 
00068     }
00069 
00070 
00071     /// The mass of this Particle.
00072     double mass() const { 
00073       return momentum().mass(); 
00074     }
00075 
00076     // /// The charge of this Particle.
00077     // double charge() const { 
00078     //   return PID::charge(*this); 
00079     // }
00080 
00081     // /// Three times the charge of this Particle (i.e. integer multiple of smallest quark charge).
00082     // int threeCharge() const { 
00083     //   return PID::threeCharge(*this); 
00084     // }
00085 
00086 
00087     /// Check whether a given PID is found in the GenParticle's ancestor list
00088     bool hasAncestor(PdgId pdg_id) const;
00089 
00090 
00091   private:
00092 
00093     /// A pointer to the original GenParticle from which this Particle is projected.
00094     const GenParticle* _original;
00095 
00096     /// The PDG ID code for this Particle.
00097     long _id;
00098 
00099     /// The momentum of this projection of the Particle.
00100     FourMomentum _momentum;
00101   };
00102 
00103 
00104   /// @name String representation
00105   //@{
00106 
00107   /// Print a ParticlePair as a string.
00108   inline std::string toString(const ParticlePair& pair) {
00109     stringstream out;
00110     out << "[" 
00111         << toParticleName(pair.first.pdgId()) << " @ "
00112         << pair.first.momentum().E()/GeV << " GeV, "
00113         << toParticleName(pair.second.pdgId()) << " @ "
00114         << pair.second.momentum().E()/GeV << " GeV]";
00115     return out.str();
00116   }
00117 
00118   /// Allow ParticlePair to be passed to an ostream.
00119   inline std::ostream& operator<<(std::ostream& os, const ParticlePair& pp) {
00120     os << toString(pp);
00121     return os;
00122   }
00123 
00124   //@}
00125 
00126 
00127   /// @name Comparison functors
00128   //@{
00129   /// Sort by descending transverse momentum, \f$ p_\perp \f$
00130   inline bool cmpParticleByPt(const Particle& a, const Particle& b) {
00131     return a.momentum().pT() > b.momentum().pT();
00132   }
00133   /// Sort by ascending transverse momentum, \f$ p_\perp \f$
00134   inline bool cmpParticleByAscPt(const Particle& a, const Particle& b) {
00135     return a.momentum().pT() < b.momentum().pT();
00136   }
00137   /// Sort by descending transverse energy, \f$ E_\perp \f$
00138   inline bool cmpParticleByEt(const Particle& a, const Particle& b) {
00139     return a.momentum().Et() > b.momentum().Et();
00140   }
00141   /// Sort by ascending transverse energy, \f$ E_\perp \f$
00142   inline bool cmpParticleByAscEt(const Particle& a, const Particle& b) {
00143     return a.momentum().Et() < b.momentum().Et();
00144   }
00145   /// Sort by descending energy, \f$ E \f$
00146   inline bool cmpParticleByE(const Particle& a, const Particle& b) {
00147     return a.momentum().E() > b.momentum().E();
00148   }
00149   /// Sort by ascending energy, \f$ E \f$
00150   inline bool cmpParticleByAscE(const Particle& a, const Particle& b) {
00151     return a.momentum().E() < b.momentum().E();
00152   }
00153   /// Sort by descending pseudorapidity, \f$ \eta \f$
00154   inline bool cmpParticleByDescPseudorapidity(const Particle& a, const Particle& b) {
00155     return a.momentum().pseudorapidity() > b.momentum().pseudorapidity();
00156   }
00157   /// Sort by ascending pseudorapidity, \f$ \eta \f$
00158   inline bool cmpParticleByAscPseudorapidity(const Particle& a, const Particle& b) {
00159     return a.momentum().pseudorapidity() < b.momentum().pseudorapidity();
00160   }
00161   /// Sort by descending absolute pseudorapidity, \f$ |\eta| \f$
00162   inline bool cmpParticleByDescAbsPseudorapidity(const Particle& a, const Particle& b) {
00163     return fabs(a.momentum().pseudorapidity()) > fabs(b.momentum().pseudorapidity());
00164   }
00165   /// Sort by ascending absolute pseudorapidity, \f$ |\eta| \f$
00166   inline bool cmpParticleByAscAbsPseudorapidity(const Particle& a, const Particle& b) {
00167     return fabs(a.momentum().pseudorapidity()) < fabs(b.momentum().pseudorapidity());
00168   }
00169   /// Sort by descending rapidity, \f$ y \f$
00170   inline bool cmpParticleByDescRapidity(const Particle& a, const Particle& b) {
00171     return a.momentum().rapidity() > b.momentum().rapidity();
00172   }
00173   /// Sort by ascending rapidity, \f$ y \f$
00174   inline bool cmpParticleByAscRapidity(const Particle& a, const Particle& b) {
00175     return a.momentum().rapidity() < b.momentum().rapidity();
00176   }
00177   /// Sort by descending absolute rapidity, \f$ |y| \f$
00178   inline bool cmpParticleByDescAbsRapidity(const Particle& a, const Particle& b) {
00179     return fabs(a.momentum().rapidity()) > fabs(b.momentum().rapidity());
00180   }
00181   /// Sort by ascending absolute rapidity, \f$ |y| \f$
00182   inline bool cmpParticleByAscAbsRapidity(const Particle& a, const Particle& b) {
00183     return fabs(a.momentum().rapidity()) < fabs(b.momentum().rapidity());
00184   }
00185   //@}
00186 
00187 
00188 }
00189 
00190 #endif