rivet is hosted by Hepforge, IPPP Durham
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/Particle.fhh"
00006 #include "Rivet/ParticleBase.hh"
00007 #include "Rivet/Config/RivetCommon.hh"
00008 #include "Rivet/Tools/Utils.hh"
00009 // NOTE: Rivet/Tools/ParticleUtils.hh included at the end
00010 #include "fastjet/PseudoJet.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     /// Constructor from a HepMC GenParticle pointer.
00040     Particle(const GenParticle* gp)
00041       : ParticleBase(),
00042         _original(gp), _id(gp->pdg_id()),
00043         _momentum(gp->momentum())
00044     { }
00045 
00046 
00047   public:
00048 
00049     /// Converter to FastJet3 PseudoJet
00050     virtual fastjet::PseudoJet pseudojet() const {
00051       return fastjet::PseudoJet(mom().px(), mom().py(), mom().pz(), mom().E());
00052     }
00053 
00054     /// Cast operator to FastJet3 PseudoJet
00055     operator PseudoJet () const { return pseudojet(); }
00056 
00057 
00058   public:
00059 
00060     /// @name Basic particle specific properties
00061     //@{
00062 
00063     /// Get a const reference to the original GenParticle.
00064     const GenParticle* genParticle() const {
00065       return _original;
00066     }
00067 
00068     /// Cast operator for conversion to GenParticle*
00069     operator const GenParticle* () const { return genParticle(); }
00070 
00071     /// The momentum.
00072     const FourMomentum& momentum() const {
00073       return _momentum;
00074     }
00075     /// Set the momentum.
00076     Particle& setMomentum(const FourMomentum& momentum) {
00077       _momentum = momentum;
00078       return *this;
00079     }
00080 
00081     //@}
00082 
00083 
00084     /// @name Particle ID code accessors
00085     //@{
00086 
00087     /// This Particle's PDG ID code.
00088     PdgId pid() const { return _id; }
00089     /// Absolute value of the PDG ID code.
00090     PdgId abspid() const { return std::abs(_id); }
00091     /// This Particle's PDG ID code (alias).
00092     /// @deprecatedThe pid/abspid form is nicer (don't need to remember
00093     ///    lower/upper case, doesn't visually stick out in the code, etc. ...)
00094     PdgId pdgId() const { return _id; }
00095 
00096     //@}
00097 
00098 
00099     /// @name Particle species properties
00100     //@{
00101 
00102     /// The charge of this Particle.
00103     double charge() const {
00104       return PID::charge(pid());
00105     }
00106     /// Three times the charge of this Particle (i.e. integer multiple of smallest quark charge).
00107     int threeCharge() const {
00108       return PID::threeCharge(pid());
00109     }
00110 
00111     /// Is this a hadron?
00112     bool isHadron() const { return PID::isHadron(pid()); }
00113 
00114     /// Is this a meson?
00115     bool isMeson() const { return PID::isMeson(pid()); }
00116 
00117     /// Is this a baryon?
00118     bool isBaryon() const { return PID::isBaryon(pid()); }
00119 
00120     /// Is this a lepton?
00121     bool isLepton() const { return PID::isLepton(pid()); }
00122 
00123     /// Is this a neutrino?
00124     bool isNeutrino() const { return PID::isNeutrino(pid()); }
00125 
00126     /// Does this (hadron) contain a b quark?
00127     bool hasBottom() const { return PID::hasBottom(pid()); }
00128 
00129     /// Does this (hadron) contain a c quark?
00130     bool hasCharm() const { return PID::hasCharm(pid()); }
00131 
00132     // /// Does this (hadron) contain an s quark?
00133     // bool hasStrange() const { return PID::hasStrange(pid()); }
00134 
00135     //@}
00136 
00137 
00138     /// @name Ancestry properties
00139     //@{
00140 
00141     /// Check whether a given PID is found in the GenParticle's ancestor list
00142     ///
00143     /// @note This question is valid in MC, but may not be answerable
00144     /// experimentally -- use this function with care when replicating
00145     /// experimental analyses!
00146     bool hasAncestor(PdgId pdg_id) const;
00147 
00148     /// @brief Determine whether the particle is from a hadron or tau decay
00149     ///
00150     /// Specifically, walk up the ancestor chain until a status 2 hadron or
00151     /// tau is found, if at all.
00152     ///
00153     /// @note This question is valid in MC, but may not be perfectly answerable
00154     /// experimentally -- use this function with care when replicating
00155     /// experimental analyses!
00156     bool fromDecay() const;
00157 
00158     /// @brief Determine whether the particle is from a b-hadron decay
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 fromBottom() const;
00164 
00165     /// @brief Determine whether the particle is from a c-hadron decay
00166     ///
00167     /// @note If a hadron contains b and c quarks it is considered a bottom
00168     /// hadron and NOT a charm 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 fromCharm() const;
00174 
00175     // /// @brief Determine whether the particle is from a s-hadron decay
00176     // ///
00177     // /// @note If a hadron contains b or c quarks as well as strange it is
00178     // /// considered a b or c hadron, but NOT a strange hadron.
00179     // ///
00180     // /// @note This question is valid in MC, but may not be perfectly answerable
00181     // /// experimentally -- use this function with care when replicating
00182     // /// experimental analyses!
00183     // bool fromStrange() const;
00184 
00185     /// @brief Determine whether the particle is from a tau decay
00186     ///
00187     /// @note This question is valid in MC, but may not be perfectly answerable
00188     /// experimentally -- use this function with care when replicating
00189     /// experimental analyses!
00190     bool fromTau() const;
00191 
00192     //@}
00193 
00194 
00195     /// @name Decay info
00196     //@{
00197 
00198     /// Whether this particle is stable according to the generator
00199     bool isStable() const {
00200       return genParticle() != NULL && genParticle()->status() == 1 && genParticle()->end_vertex() == NULL;
00201     }
00202 
00203     /// Get a list of the direct descendants from the current particle
00204     vector<Particle> children() const {
00205       vector<Particle> rtn;
00206       if (isStable()) return rtn;
00207       /// @todo Remove this const mess crap when HepMC doesn't suck
00208       HepMC::GenVertex* gv = const_cast<HepMC::GenVertex*>( genParticle()->end_vertex() );
00209       /// @todo Would like to do this, but the range objects are broken
00210       // foreach (const GenParticle* gp, gv->particles(HepMC::children))
00211       //   rtn += Particle(gp);
00212       for (GenVertex::particle_iterator it = gv->particles_begin(HepMC::children); it != gv->particles_end(HepMC::children); ++it)
00213         rtn += Particle(*it);
00214       return rtn;
00215     }
00216 
00217     /// Get a list of all the descendants (including duplication of parents and children) from the current particle
00218     /// @todo Use recursion through replica-avoiding MCUtils functions to avoid bookkeeping duplicates
00219     /// @todo Insist that the current particle is post-hadronization, otherwise throw an exception?
00220     vector<Particle> allDescendants() const {
00221       vector<Particle> rtn;
00222       if (isStable()) return rtn;
00223       /// @todo Remove this const mess crap when HepMC doesn't suck
00224       HepMC::GenVertex* gv = const_cast<HepMC::GenVertex*>( genParticle()->end_vertex() );
00225       /// @todo Would like to do this, but the range objects are broken
00226       // foreach (const GenParticle* gp, gv->particles(HepMC::descendants))
00227       //   rtn += Particle(gp);
00228       for (GenVertex::particle_iterator it = gv->particles_begin(HepMC::descendants); it != gv->particles_end(HepMC::descendants); ++it)
00229         rtn += Particle(*it);
00230       return rtn;
00231     }
00232 
00233     /// Get a list of all the stable descendants from the current particle
00234     /// @todo Use recursion through replica-avoiding MCUtils functions to avoid bookkeeping duplicates
00235     /// @todo Insist that the current particle is post-hadronization, otherwise throw an exception?
00236     vector<Particle> stableDescendants() const {
00237       vector<Particle> rtn;
00238       if (isStable()) return rtn;
00239       /// @todo Remove this const mess crap when HepMC doesn't suck
00240       HepMC::GenVertex* gv = const_cast<HepMC::GenVertex*>( genParticle()->end_vertex() );
00241       /// @todo Would like to do this, but the range objects are broken
00242       // foreach (const GenParticle* gp, gv->particles(HepMC::descendants))
00243       //   if (gp->status() == 1 && gp->end_vertex() == NULL)
00244       //     rtn += Particle(gp);
00245       for (GenVertex::particle_iterator it = gv->particles_begin(HepMC::descendants); it != gv->particles_end(HepMC::descendants); ++it)
00246         if ((*it)->status() == 1 && (*it)->end_vertex() == NULL)
00247           rtn += Particle(*it);
00248       return rtn;
00249     }
00250 
00251     /// Flight length (divide by mm or cm to get the appropriate units)
00252     double flightLength() const {
00253       if (isStable()) return -1;
00254       if (genParticle() == NULL) return 0;
00255       if (genParticle()->production_vertex() == NULL) return 0;
00256       const HepMC::FourVector v1 = genParticle()->production_vertex()->position();
00257       const HepMC::FourVector v2 = genParticle()->end_vertex()->position();
00258       return sqrt(sqr(v2.x()-v1.x()) + sqr(v2.y()-v1.y()) + sqr(v2.z()-v1.z()));
00259     }
00260 
00261     //@}
00262 
00263 
00264   private:
00265 
00266     /// A pointer to the original GenParticle from which this Particle is projected.
00267     const GenParticle* _original;
00268 
00269     /// The PDG ID code for this Particle.
00270     PdgId _id;
00271 
00272     /// The momentum of this projection of the Particle.
00273     FourMomentum _momentum;
00274 
00275   };
00276 
00277 
00278   /// @name String representation
00279   //@{
00280 
00281   /// Print a ParticlePair as a string.
00282   inline std::string to_str(const ParticlePair& pair) {
00283     stringstream out;
00284     out << "["
00285         << PID::toParticleName(pair.first.pid()) << " @ "
00286         << pair.first.momentum().E()/GeV << " GeV, "
00287         << PID::toParticleName(pair.second.pid()) << " @ "
00288         << pair.second.momentum().E()/GeV << " GeV]";
00289     return out.str();
00290   }
00291 
00292   /// Allow ParticlePair to be passed to an ostream.
00293   inline std::ostream& operator<<(std::ostream& os, const ParticlePair& pp) {
00294     os << to_str(pp);
00295     return os;
00296   }
00297 
00298   //@}
00299 
00300 
00301 }
00302 
00303 
00304 #include "Rivet/Tools/ParticleUtils.hh"
00305 
00306 #endif