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     /// @note 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     /// @deprecated Prefer the pid/abspid form
00093     PdgId pdgId() const { return _id; }
00094 
00095     //@}
00096 
00097 
00098     /// @name Particle species properties
00099     //@{
00100 
00101     /// The charge of this Particle.
00102     double charge() const {
00103       return PID::charge(pid());
00104     }
00105     /// Three times the charge of this Particle (i.e. integer multiple of smallest quark charge).
00106     int threeCharge() const {
00107       return PID::threeCharge(pid());
00108     }
00109 
00110     /// Is this a hadron?
00111     bool isHadron() const { return PID::isHadron(pid()); }
00112 
00113     /// Is this a meson?
00114     bool isMeson() const { return PID::isMeson(pid()); }
00115 
00116     /// Is this a baryon?
00117     bool isBaryon() const { return PID::isBaryon(pid()); }
00118 
00119     /// Is this a lepton?
00120     bool isLepton() const { return PID::isLepton(pid()); }
00121 
00122     /// Is this a neutrino?
00123     bool isNeutrino() const { return PID::isNeutrino(pid()); }
00124 
00125     /// Does this (hadron) contain a b quark?
00126     bool hasBottom() const { return PID::hasBottom(pid()); }
00127 
00128     /// Does this (hadron) contain a c quark?
00129     bool hasCharm() const { return PID::hasCharm(pid()); }
00130 
00131     // /// Does this (hadron) contain an s quark?
00132     // bool hasStrange() const { return PID::hasStrange(pid()); }
00133 
00134     //@}
00135 
00136 
00137     /// @name Ancestry properties
00138     //@{
00139 
00140     /// Check whether a given PID is found in the GenParticle's ancestor list
00141     ///
00142     /// @note This question is valid in MC, but may not be answerable
00143     /// experimentally -- use this function with care when replicating
00144     /// experimental analyses!
00145     bool hasAncestor(PdgId pdg_id) const;
00146 
00147     /// @brief Determine whether the particle is from a b-hadron decay
00148     ///
00149     /// @note This question is valid in MC, but may not be perfectly answerable
00150     /// experimentally -- use this function with care when replicating
00151     /// experimental analyses!
00152     bool fromBottom() const;
00153 
00154     /// @brief Determine whether the particle is from a c-hadron decay
00155     ///
00156     /// @note If a hadron contains b and c quarks it is considered a bottom
00157     /// hadron and NOT a charm hadron.
00158     ///
00159     /// @note This question is valid in MC, but may not be perfectly answerable
00160     /// experimentally -- use this function with care when replicating
00161     /// experimental analyses!
00162     bool fromCharm() const;
00163 
00164     // /// @brief Determine whether the particle is from a s-hadron decay
00165     // ///
00166     // /// @note If a hadron contains b or c quarks as well as strange it is
00167     // /// considered a b or c hadron, but NOT a strange hadron.
00168     // ///
00169     // /// @note This question is valid in MC, but may not be perfectly answerable
00170     // /// experimentally -- use this function with care when replicating
00171     // /// experimental analyses!
00172     // bool fromStrange() const;
00173 
00174     /// @brief Determine whether the particle is from a hadron decay
00175     ///
00176     /// @note This question is valid in MC, but may not be perfectly answerable
00177     /// experimentally -- use this function with care when replicating
00178     /// experimental analyses!
00179     bool fromHadron() const;
00180 
00181     /// @brief Determine whether the particle is from a tau decay
00182     ///
00183     /// @note This question is valid in MC, but may not be perfectly answerable
00184     /// experimentally -- use this function with care when replicating
00185     /// experimental analyses!
00186     bool fromTau(bool prompt_taus_only=false) const;
00187 
00188     /// @brief Determine whether the particle is from a prompt tau decay
00189     ///
00190     /// @note This question is valid in MC, but may not be perfectly answerable
00191     /// experimentally -- use this function with care when replicating
00192     /// experimental analyses!
00193     bool fromPromptTau() const { return fromTau(true); }
00194 
00195     /// @brief Determine whether the particle is from a hadron or tau decay
00196     ///
00197     /// Specifically, walk up the ancestor chain until a status 2 hadron or
00198     /// tau is found, if at all.
00199     ///
00200     /// @note This question is valid in MC, but may not be perfectly answerable
00201     /// experimentally -- use this function with care when replicating
00202     /// experimental analyses!
00203     bool fromDecay() const { return fromHadron() || fromPromptTau(); }
00204 
00205     //@}
00206 
00207 
00208     /// @name Decay info
00209     //@{
00210 
00211     /// Whether this particle is stable according to the generator
00212     bool isStable() const {
00213       return genParticle() != NULL && genParticle()->status() == 1 && genParticle()->end_vertex() == NULL;
00214     }
00215 
00216     /// Get a list of the direct descendants from the current particle
00217     vector<Particle> children() const {
00218       vector<Particle> rtn;
00219       if (isStable()) return rtn;
00220       /// @todo Remove this const mess crap when HepMC doesn't suck
00221       HepMC::GenVertex* gv = const_cast<HepMC::GenVertex*>( genParticle()->end_vertex() );
00222       /// @todo Would like to do this, but the range objects are broken
00223       // foreach (const GenParticle* gp, gv->particles(HepMC::children))
00224       //   rtn += Particle(gp);
00225       for (GenVertex::particle_iterator it = gv->particles_begin(HepMC::children); it != gv->particles_end(HepMC::children); ++it)
00226         rtn += Particle(*it);
00227       return rtn;
00228     }
00229 
00230     /// Get a list of all the descendants (including duplication of parents and children) from the current particle
00231     /// @todo Use recursion through replica-avoiding MCUtils functions to avoid bookkeeping duplicates
00232     /// @todo Insist that the current particle is post-hadronization, otherwise throw an exception?
00233     vector<Particle> allDescendants() const {
00234       vector<Particle> rtn;
00235       if (isStable()) return rtn;
00236       /// @todo Remove this const mess crap when HepMC doesn't suck
00237       HepMC::GenVertex* gv = const_cast<HepMC::GenVertex*>( genParticle()->end_vertex() );
00238       /// @todo Would like to do this, but the range objects are broken
00239       // foreach (const GenParticle* gp, gv->particles(HepMC::descendants))
00240       //   rtn += Particle(gp);
00241       for (GenVertex::particle_iterator it = gv->particles_begin(HepMC::descendants); it != gv->particles_end(HepMC::descendants); ++it)
00242         rtn += Particle(*it);
00243       return rtn;
00244     }
00245 
00246     /// Get a list of all the stable descendants from the current particle
00247     /// @todo Use recursion through replica-avoiding MCUtils functions to avoid bookkeeping duplicates
00248     /// @todo Insist that the current particle is post-hadronization, otherwise throw an exception?
00249     vector<Particle> stableDescendants() const {
00250       vector<Particle> rtn;
00251       if (isStable()) return rtn;
00252       /// @todo Remove this const mess crap when HepMC doesn't suck
00253       HepMC::GenVertex* gv = const_cast<HepMC::GenVertex*>( genParticle()->end_vertex() );
00254       /// @todo Would like to do this, but the range objects are broken
00255       // foreach (const GenParticle* gp, gv->particles(HepMC::descendants))
00256       //   if (gp->status() == 1 && gp->end_vertex() == NULL)
00257       //     rtn += Particle(gp);
00258       for (GenVertex::particle_iterator it = gv->particles_begin(HepMC::descendants); it != gv->particles_end(HepMC::descendants); ++it)
00259         if ((*it)->status() == 1 && (*it)->end_vertex() == NULL)
00260           rtn += Particle(*it);
00261       return rtn;
00262     }
00263 
00264     /// Flight length (divide by mm or cm to get the appropriate units)
00265     double flightLength() const {
00266       if (isStable()) return -1;
00267       if (genParticle() == NULL) return 0;
00268       if (genParticle()->production_vertex() == NULL) return 0;
00269       const HepMC::FourVector v1 = genParticle()->production_vertex()->position();
00270       const HepMC::FourVector v2 = genParticle()->end_vertex()->position();
00271       return sqrt(sqr(v2.x()-v1.x()) + sqr(v2.y()-v1.y()) + sqr(v2.z()-v1.z()));
00272     }
00273 
00274     //@}
00275 
00276 
00277   private:
00278 
00279     /// A pointer to the original GenParticle from which this Particle is projected.
00280     const GenParticle* _original;
00281 
00282     /// The PDG ID code for this Particle.
00283     PdgId _id;
00284 
00285     /// The momentum of this projection of the Particle.
00286     FourMomentum _momentum;
00287 
00288   };
00289 
00290 
00291   /// @name String representation
00292   //@{
00293 
00294   /// Print a ParticlePair as a string.
00295   inline std::string to_str(const ParticlePair& pair) {
00296     stringstream out;
00297     out << "["
00298         << PID::toParticleName(pair.first.pid()) << " @ "
00299         << pair.first.momentum().E()/GeV << " GeV, "
00300         << PID::toParticleName(pair.second.pid()) << " @ "
00301         << pair.second.momentum().E()/GeV << " GeV]";
00302     return out.str();
00303   }
00304 
00305   /// Allow ParticlePair to be passed to an ostream.
00306   inline std::ostream& operator<<(std::ostream& os, const ParticlePair& pp) {
00307     os << to_str(pp);
00308     return os;
00309   }
00310 
00311   //@}
00312 
00313 
00314 }
00315 
00316 
00317 #include "Rivet/Tools/ParticleUtils.hh"
00318 
00319 #endif