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