rivet is hosted by Hepforge, IPPP Durham
FinalState.hh
Go to the documentation of this file.
00001 // -*- C++ -*-
00002 #ifndef RIVET_FinalState_HH
00003 #define RIVET_FinalState_HH
00004 
00005 #include "Rivet/Projection.hh"
00006 #include "Rivet/Particle.hh"
00007 #include "Rivet/Event.hh"
00008 
00009 namespace Rivet {
00010 
00011 
00012   /// @brief Project out all final-state particles in an event.
00013   /// Probably the most important projection in Rivet!
00014   class FinalState : public Projection {
00015   public:
00016 
00017     /// @name Standard constructors and destructors.
00018     //@{
00019     /// The default constructor. May specify the minimum and maximum
00020     /// pseudorapidity \f$ \eta \f$ and the min \f$ p_T \f$ (in GeV).
00021     FinalState(double mineta = -MAXDOUBLE,
00022                double maxeta =  MAXDOUBLE,
00023                double minpt  =  0.0*GeV);
00024 
00025     /// A constructor which allows to specify multiple eta ranges
00026     /// and the min \f$ p_T \f$ (in GeV).
00027     FinalState(const vector<pair<double, double> >& etaRanges,
00028                double minpt = 0.0*GeV);
00029 
00030     /// Clone on the heap.
00031     virtual const Projection* clone() const {
00032       return new FinalState(*this);
00033     }
00034 
00035     //@}
00036 
00037 
00038     /// Get the final-state particles.
00039     virtual const Particles& particles() const { return _theParticles; }
00040 
00041     /// Get the final-state particles, ordered by supplied sorting function object.
00042     template <typename F>
00043     const Particles& particles(F sorter) const {
00044       std::sort(_theParticles.begin(), _theParticles.end(), sorter);
00045       return _theParticles;
00046     }
00047 
00048     /// Get the final-state particles, ordered by decreasing \f$ p_T \f$.
00049     const Particles& particlesByPt() const {
00050       return particles(cmpMomByPt);
00051     }
00052 
00053     /// Get the final-state particles, ordered by decreasing \f$ p \f$.
00054     const Particles& particlesByP() const {
00055       return particles(cmpMomByP);
00056     }
00057 
00058     /// Get the final-state particles, ordered by decreasing \f$ E \f$.
00059     const Particles& particlesByE() const {
00060       return particles(cmpMomByE);
00061     }
00062 
00063     /// Get the final-state particles, ordered by decreasing \f$ E_T \f$.
00064     const Particles& particlesByEt() const {
00065       return particles(cmpMomByEt);
00066     }
00067 
00068     /// Get the final-state particles, ordered by increasing \f$ \eta \f$.
00069     const Particles& particlesByEta() const {
00070       return particles(cmpMomByAscPseudorapidity);
00071     }
00072 
00073     /// Get the final-state particles, ordered by increasing \f$ |\eta| \f$.
00074     const Particles& particlesByModEta() const {
00075       return particles(cmpMomByAscAbsPseudorapidity);
00076     }
00077 
00078     /// Get the final-state particles, ordered by increasing \f$ y \f$.
00079     const Particles& particlesByRapidity() const {
00080       return particles(cmpMomByAscRapidity);
00081     }
00082 
00083     /// Get the final-state particles, ordered by increasing \f$ |y| \f$.
00084     const Particles& particlesByModRapidity() const {
00085       return particles(cmpMomByAscAbsRapidity);
00086     }
00087 
00088     /// Access the projected final-state particles.
00089     virtual size_t size() const { return _theParticles.size(); }
00090 
00091     /// Is this final state empty?
00092     virtual bool empty() const { return _theParticles.empty(); }
00093     /// @deprecated Is this final state empty?
00094     virtual bool isEmpty() const { return _theParticles.empty(); }
00095 
00096     /// Minimum-\f$ p_\perp \f$ requirement.
00097     virtual double ptMin() const { return _ptmin; }
00098 
00099 
00100   public:
00101 
00102     typedef Particle entity_type;
00103     typedef Particles collection_type;
00104 
00105     /// Template-usable interface common to JetAlg.
00106     const collection_type& entities() const {
00107       return particles();
00108     }
00109 
00110 
00111   protected:
00112 
00113     /// Apply the projection to the event.
00114     virtual void project(const Event& e);
00115 
00116     /// Compare projections.
00117     virtual int compare(const Projection& p) const;
00118 
00119     /// Decide if a particle is to be accepted or not.
00120     bool accept(const Particle& p) const;
00121 
00122 
00123   protected:
00124 
00125     /// The ranges allowed for pseudorapidity.
00126     vector<pair<double,double> > _etaRanges;
00127 
00128     /// The minimum allowed transverse momentum.
00129     double _ptmin;
00130 
00131     /// The final-state particles.
00132     mutable Particles _theParticles;
00133 
00134   };
00135 
00136 
00137 }
00138 
00139 #endif