2 #ifndef RIVET_FinalState_HH
3 #define RIVET_FinalState_HH
5 #include "Rivet/Projection.hh"
6 #include "Rivet/Particle.hh"
7 #include "Rivet/Event.hh"
23 double minpt = 0.0*GeV);
27 FinalState(
const vector<pair<double, double> >& etaRanges,
28 double minpt = 0.0*GeV);
39 virtual const ParticleVector&
particles()
const {
return _theParticles; }
44 std::sort(_theParticles.begin(), _theParticles.end(), sorter);
89 virtual size_t size()
const {
return _theParticles.size(); }
92 virtual bool empty()
const {
return _theParticles.empty(); }
94 virtual bool isEmpty()
const {
return _theParticles.empty(); }
97 virtual double ptMin()
const {
return _ptmin; }
103 typedef ParticleVector collection_type;
126 vector<pair<double,double> > _etaRanges;
132 mutable ParticleVector _theParticles;
Definition: MC_JetAnalysis.hh:9
bool cmpParticleByEt(const Particle &a, const Particle &b)
Sort by descending transverse energy, .
Definition: Particle.hh:170
const ParticleVector & particlesByEt() const
Get the final-state particles, ordered by decreasing .
Definition: FinalState.hh:64
const collection_type & entities() const
Template-usable interface common to JetAlg.
Definition: FinalState.hh:106
bool cmpParticleByAscAbsRapidity(const Particle &a, const Particle &b)
Sort by ascending absolute rapidity, .
Definition: Particle.hh:214
const ParticleVector & particlesByP() const
Get the final-state particles, ordered by decreasing .
Definition: FinalState.hh:54
virtual const Projection * clone() const
Clone on the heap.
Definition: FinalState.hh:31
bool cmpParticleByAscRapidity(const Particle &a, const Particle &b)
Sort by ascending rapidity, .
Definition: Particle.hh:206
const ParticleVector & particlesByModEta() const
Get the final-state particles, ordered by increasing .
Definition: FinalState.hh:74
virtual size_t size() const
Access the projected final-state particles.
Definition: FinalState.hh:89
Representation of particles from a HepMC::GenEvent.
Definition: Particle.hh:16
const ParticleVector & particlesByRapidity() const
Get the final-state particles, ordered by increasing .
Definition: FinalState.hh:79
const ParticleVector & particlesByPt() const
Get the final-state particles, ordered by decreasing .
Definition: FinalState.hh:49
virtual const ParticleVector & particles() const
Get the final-state particles.
Definition: FinalState.hh:39
const ParticleVector & particlesByModRapidity() const
Get the final-state particles, ordered by increasing .
Definition: FinalState.hh:84
virtual int compare(const Projection &p) const
Compare projections.
Definition: FinalState.cc:44
const ParticleVector & particles(F sorter) const
Get the final-state particles, ordered by supplied sorting function object.
Definition: FinalState.hh:43
bool cmpParticleByAscPseudorapidity(const Particle &a, const Particle &b)
Sort by ascending pseudorapidity, .
Definition: Particle.hh:190
virtual bool isEmpty() const
Definition: FinalState.hh:94
bool accept(const Particle &p) const
Decide if a particle is to be accepted or not.
Definition: FinalState.cc:95
Project out all final-state particles in an event. Probably the most important projection in Rivet! ...
Definition: FinalState.hh:14
bool cmpParticleByE(const Particle &a, const Particle &b)
Sort by descending energy, .
Definition: Particle.hh:178
FinalState(double mineta=-MAXRAPIDITY, double maxeta=MAXRAPIDITY, double minpt=0.0 *GeV)
Definition: FinalState.cc:9
static const double MAXRAPIDITY
A sensible default maximum value of rapidity for Rivet analyses to use.
Definition: Rivet.hh:24
virtual void project(const Event &e)
Apply the projection to the event.
Definition: FinalState.cc:63
const ParticleVector & particlesByEta() const
Get the final-state particles, ordered by increasing .
Definition: FinalState.hh:69
bool cmpParticleByP(const Particle &a, const Particle &b)
Sort by descending momentum, .
Definition: Particle.hh:162
virtual bool empty() const
Is this final state empty?
Definition: FinalState.hh:92
Base class for all Rivet projections.
Definition: Projection.hh:28
bool cmpParticleByPt(const Particle &a, const Particle &b)
Sort by descending transverse momentum, .
Definition: Particle.hh:154
virtual double ptMin() const
Minimum- requirement.
Definition: FinalState.hh:97
const ParticleVector & particlesByE() const
Get the final-state particles, ordered by decreasing .
Definition: FinalState.hh:59
bool cmpParticleByAscAbsPseudorapidity(const Particle &a, const Particle &b)
Sort by ascending absolute pseudorapidity, .
Definition: Particle.hh:198