InvMassFinalState Class Reference

Identify particles which can be paired to fit within a given invariant mass window. More...

#include <InvMassFinalState.hh>

Inheritance diagram for InvMassFinalState:
Inheritance graph
[legend]
Collaboration diagram for InvMassFinalState:
Collaboration graph
[legend]

List of all members.

Public Types

typedef Particle entity_type
typedef ParticleVector collection_type

Public Member Functions

 InvMassFinalState (const FinalState &fsp, const std::pair< PdgId, PdgId > &idpair, double minmass, double maxmass, double masstarget=-1.0)
 Constructor for a single inv-mass pair.
 InvMassFinalState (const FinalState &fsp, const std::vector< std::pair< PdgId, PdgId > > &idpairs, double minmass, double maxmass, double masstarget=-1.0)
 Constructor for multiple inv-mass pairs.
 InvMassFinalState (const std::pair< PdgId, PdgId > &idpair, double minmass, double maxmass, double masstarget=-1.0)
 Same thing as above, but we want to pass the particles directly to the calc method.
 InvMassFinalState (const std::vector< std::pair< PdgId, PdgId > > &idpairs, double minmass, double maxmass, double masstarget=-1.0)
virtual const Projectionclone () const
 Clone on the heap.
const std::vector< std::pair
< Particle, Particle > > & 
particlePairs () const
 Constituent pairs.
void useTransverseMass (bool usetrans=true)
 Choose whether to use the full inv mass or just the transverse mass.
void calc (const ParticleVector &inparticles)
 Operate on a given particle vector directly instead of through project (no caching).
virtual const ParticleVectorparticles () const
 Get the final-state particles.
template<typename F >
const ParticleVectorparticles (F sorter) const
 Get the final-state particles, ordered by supplied sorting function object.
const ParticleVectorparticlesByPt () const
 Get the final-state particles, ordered by decreasing $ p_T $.
const ParticleVectorparticlesByP () const
 Get the final-state particles, ordered by decreasing $ p $.
const ParticleVectorparticlesByE () const
 Get the final-state particles, ordered by decreasing $ E $.
const ParticleVectorparticlesByEt () const
 Get the final-state particles, ordered by decreasing $ E_T $.
const ParticleVectorparticlesByEta () const
 Get the final-state particles, ordered by increasing $ \eta $.
const ParticleVectorparticlesByModEta () const
 Get the final-state particles, ordered by increasing $ |\eta| $.
virtual size_t size () const
 Access the projected final-state particles.
virtual bool empty () const
 Is this final state empty?
virtual bool isEmpty () const
virtual double ptMin () const
 Minimum-$ p_\perp $ requirement.
const collection_typeentities () const
 Template-usable interface common to JetAlg.
bool before (const Projection &p) const
virtual const std::set< PdgIdPairbeamPairs () const
virtual std::string name () const
 Get the name of the projection.
ProjectionaddPdgIdPair (PdgId beam1, PdgId beam2)
 Add a colliding beam pair.
LoggetLog () const
 Get a Log object based on the getName() property of the calling projection object.
void setName (const std::string &name)
 Used by derived classes to set their name.
Projection "getting" functions

std::set< ConstProjectionPtrgetProjections () const
 Get the contained projections, including recursion.
template<typename PROJ >
const PROJ & getProjection (const std::string &name) const
 Get the named projection, specifying return type via a template argument.
const ProjectiongetProjection (const std::string &name) const
Projection applying functions

template<typename PROJ >
const PROJ & applyProjection (const Event &evt, const PROJ &proj) const
 Apply the supplied projection on event.
template<typename PROJ >
const PROJ & applyProjection (const Event &evt, const Projection &proj) const
 Apply the supplied projection on event.
template<typename PROJ >
const PROJ & applyProjection (const Event &evt, const std::string &name) const
 Apply the named projection on event.

Protected Member Functions

void project (const Event &e)
 Apply the projection on the supplied event.
int compare (const Projection &p) const
 Compare projections.
bool accept (const Particle &p) const
 Decide if a particle is to be accepted or not.
Cmp< ProjectionmkNamedPCmp (const Projection &otherparent, const std::string &pname) const
Cmp< ProjectionmkPCmp (const Projection &otherparent, const std::string &pname) const
ProjectionHandlergetProjHandler () const
 Get a reference to the ProjectionHandler for this thread.
Projection registration functions

template<typename PROJ >
const PROJ & addProjection (const PROJ &proj, const std::string &name)
const Projection_addProjection (const Projection &proj, const std::string &name)
 Untemplated function to do the work...

Protected Attributes

vector< pair< double, double > > _etaRanges
 The ranges allowed for pseudorapidity.
double _ptmin
 The minimum allowed transverse momentum.
ParticleVector _theParticles
 The final-state particles.
bool _allowProjReg
 Flag to forbid projection registration in analyses until the init phase.

Private Attributes

std::vector< PdgIdPair_decayids
 IDs of the decay products.
std::vector< std::pair
< Particle, Particle > > 
_particlePairs
 Constituent pairs.
double _minmass
 Min inv mass.
double _maxmass
 Max inv mass.
double _masstarget
 Target mass if only one pair should be returned.
bool _useTransverseMass
 Flag to decide whether to use the full inv mass or just the transverse mass.

Friends

class Event
 Event is a friend.
class Cmp< Projection >
 The Cmp specialization for Projection is a friend.

Detailed Description

Identify particles which can be paired to fit within a given invariant mass window.

Definition at line 11 of file InvMassFinalState.hh.


Member Typedef Documentation

typedef ParticleVector collection_type [inherited]

Definition at line 93 of file FinalState.hh.

typedef Particle entity_type [inherited]

Definition at line 92 of file FinalState.hh.


Constructor & Destructor Documentation

InvMassFinalState ( const FinalState fsp,
const std::pair< PdgId, PdgId > &  idpair,
double  minmass,
double  maxmass,
double  masstarget = -1.0 
)

Constructor for a single inv-mass pair.

Referenced by InvMassFinalState::clone().

InvMassFinalState ( const FinalState fsp,
const std::vector< std::pair< PdgId, PdgId > > &  idpairs,
double  minmass,
double  maxmass,
double  masstarget = -1.0 
)

Constructor for multiple inv-mass pairs.

InvMassFinalState ( const std::pair< PdgId, PdgId > &  idpair,
double  minmass,
double  maxmass,
double  masstarget = -1.0 
)

Same thing as above, but we want to pass the particles directly to the calc method.

InvMassFinalState ( const std::vector< std::pair< PdgId, PdgId > > &  idpairs,
double  minmass,
double  maxmass,
double  masstarget = -1.0 
)

Member Function Documentation

const Projection & _addProjection ( const Projection proj,
const std::string &  name 
) [protected, inherited]

Untemplated function to do the work...

Definition at line 33 of file ProjectionApplier.cc.

References ProjectionApplier::_allowProjReg, ProjectionApplier::getProjHandler(), ProjectionApplier::name(), Projection::name(), and ProjectionHandler::registerProjection().

Referenced by ProjectionApplier::addProjection().

00034                                                                              {
00035     if (!_allowProjReg) {
00036       cerr << "Trying to register projection '"
00037            << proj.name() << "' before init phase in '" << this->name() << "'." << endl;
00038       exit(2);
00039     }
00040     const Projection& reg = getProjHandler().registerProjection(*this, proj, name);
00041     return reg;
00042   }

bool accept ( const Particle p  )  const [protected, inherited]

Decide if a particle is to be accepted or not.

Definition at line 95 of file FinalState.cc.

References FinalState::_etaRanges, FinalState::_ptmin, FourVector::eta(), Rivet::eta(), Particle::genParticle(), Particle::hasGenParticle(), Particle::momentum(), FourMomentum::pT(), and Rivet::pT().

Referenced by InvMassFinalState::calc(), LeptonClusters::project(), LeadingParticlesFinalState::project(), and FinalState::project().

00095                                                  {
00096     // Not having s.c. == 1 should never happen!
00097     assert(!p.hasGenParticle() || p.genParticle().status() == 1);
00098 
00099     // Check pT cut
00100     if (_ptmin > 0.0) {
00101       const double pT = p.momentum().pT();
00102       if (pT < _ptmin) return false;
00103     }
00104 
00105     // Check eta cuts
00106     if (!_etaRanges.empty()) {
00107       bool eta_pass = false;
00108       const double eta = p.momentum().eta();
00109       typedef pair<double,double> EtaPair;
00110       foreach (const EtaPair& etacuts, _etaRanges) {
00111         if (eta > etacuts.first && eta < etacuts.second) {
00112           eta_pass = true;
00113           break;
00114         }
00115       }
00116       if (!eta_pass) return false;
00117     }
00118 
00119     return true;
00120   }

Projection& addPdgIdPair ( PdgId  beam1,
PdgId  beam2 
) [inline, inherited]

Add a colliding beam pair.

Definition at line 107 of file Projection.hh.

References Projection::_beamPairs.

Referenced by Projection::Projection().

00107                                                        {
00108       _beamPairs.insert(PdgIdPair(beam1, beam2));
00109       return *this;
00110     }

const PROJ& addProjection ( const PROJ &  proj,
const std::string &  name 
) [inline, protected, inherited]

Register a contained projection. The type of the argument is used to instantiate a new projection internally: this new object is applied to events rather than the argument object. Hence you are advised to only use locally-scoped Projection objects in your Projection and Analysis constructors, and to avoid polymorphism (e.g. handling ConcreteProjection via a pointer or reference to type Projection) since this will screw up the internal type management.

Definition at line 113 of file ProjectionApplier.hh.

References ProjectionApplier::_addProjection().

Referenced by ZFinder::_init(), WFinder::_init(), VetoedFinalState::addVetoOnThisFinalState(), BeamThrust::BeamThrust(), CDF_2009_S8057893::CDF_2009_S8057893::init(), CentralEtHCM::CentralEtHCM(), ChargedFinalState::ChargedFinalState(), ChargedLeptons::ChargedLeptons(), ClusteredPhotons::ClusteredPhotons(), DISFinalState::DISFinalState(), DISKinematics::DISKinematics(), DISLepton::DISLepton(), FinalState::FinalState(), FoxWolframMoments::FoxWolframMoments(), FParameter::FParameter(), HadronicFinalState::HadronicFinalState(), Hemispheres::Hemispheres(), IdentifiedFinalState::IdentifiedFinalState(), ZEUS_2001_S4815815::init(), UA5_1989_S1926373::init(), UA5_1988_S1867512::init(), UA5_1987_S1640666::init(), UA5_1986_S1583476::init(), UA5_1982_S875503::init(), UA1_1990_S2044935::init(), TASSO_1990_S2148048::init(), STAR_2009_UE_HELEN::init(), STAR_2008_S7993412::init(), STAR_2008_S7869363::init(), STAR_2006_S6870392::init(), STAR_2006_S6860818::init(), STAR_2006_S6500200::init(), SFM_1984_S1178091::init(), PDG_HADRON_MULTIPLICITIES_RATIOS::init(), PDG_HADRON_MULTIPLICITIES::init(), OPAL_2004_S6132243::init(), OPAL_2001_S4553896::init(), OPAL_1998_S3780481::init(), OPAL_1993_S2692198::init(), MC_ZZJETS::init(), MC_ZJETS::init(), MC_WWJETS::init(), MC_WPOL::init(), MC_WJETS::init(), MC_VH2BB::init(), MC_TTBAR::init(), MC_SUSY::init(), MC_PHOTONJETUE::init(), MC_PHOTONJETS::init(), MC_LEADINGJETS::init(), MC_JETS::init(), MC_HJETS::init(), MC_GENERIC::init(), MC_DIPHOTON::init(), MC_DIJET::init(), LHCB_2010_S8758301::init(), JADE_OPAL_2000_S4300807::init(), JADE_1998_S3612880::init(), H1_2000_S4129130::init(), H1_1995_S3167097::init(), H1_1994_S2919893::init(), ExampleAnalysis::init(), E735_1998_S3905616::init(), DELPHI_2003_WUD_03_11::init(), DELPHI_2002_069_CONF_603::init(), DELPHI_1996_S3430090::init(), DELPHI_1995_S3137023::init(), D0_2010_S8821313::init(), D0_2010_S8671338::init(), D0_2010_S8570965::init(), D0_2010_S8566488::init(), D0_2009_S8349509::init(), D0_2009_S8320160::init(), D0_2009_S8202443::init(), D0_2008_S7863608::init(), D0_2008_S7837160::init(), D0_2008_S7719523::init(), D0_2008_S7662670::init(), D0_2008_S7554427::init(), D0_2008_S6879055::init(), D0_2007_S7075677::init(), D0_2006_S6438750::init(), D0_2004_S5992206::init(), D0_2001_S4674421::init(), D0_2000_S4480767::init(), D0_1996_S3324664::init(), D0_1996_S3214044::init(), CMS_2011_S9120041::init(), CMS_2011_S9088458::init(), CMS_2011_S9086218::init(), CMS_2011_S8978280::init(), CMS_2011_S8968497::init(), CMS_2011_S8957746::init(), CMS_2011_S8950903::init(), CMS_2011_S8884919::init(), CMS_2010_S8656010::init(), CMS_2010_S8547297::init(), CDF_2010_S8591881_QCD::init(), CDF_2010_S8591881_DY::init(), CDF_2009_S8436959::init(), CDF_2009_S8383952::init(), CDF_2009_S8233977::init(), CDF_2009_NOTE_9936::init(), CDF_2008_S8095620::init(), CDF_2008_S8093652::init(), CDF_2008_S7828950::init(), CDF_2008_S7782535::init(), CDF_2008_S7541902::init(), CDF_2008_S7540469::init(), CDF_2008_NOTE_9351::init(), CDF_2008_LEADINGJETS::init(), CDF_2007_S7057202::init(), CDF_2006_S6653332::init(), CDF_2006_S6450792::init(), CDF_2005_S6217184::init(), CDF_2005_S6080774::init(), CDF_2004_S5839831::init(), CDF_2001_S4751469::init(), CDF_2001_S4563131::init(), CDF_2001_S4517016::init(), CDF_2000_S4266730::init(), CDF_2000_S4155203::init(), CDF_1998_S3618439::init(), CDF_1997_S3541940::init(), CDF_1996_S3418421::init(), CDF_1996_S3349578::init(), CDF_1996_S3108457::init(), CDF_1994_S2952106::init(), CDF_1993_S2742446::init(), CDF_1990_S2089246::init(), CDF_1988_S1865951::init(), BELLE_2006_S6265367::init(), ATLAS_2011_S9225137::init(), ATLAS_2011_S9212183::init(), ATLAS_2011_S9131140::init(), ATLAS_2011_S9128077::init(), ATLAS_2011_S9126244::init(), ATLAS_2011_S9120807::init(), ATLAS_2011_S9108483::init(), ATLAS_2011_S9041966::init(), ATLAS_2011_S9019561::init(), ATLAS_2011_S9002537::init(), ATLAS_2011_S8994773::init(), ATLAS_2011_S8983313::init(), ATLAS_2011_S8971293::init(), ATLAS_2011_S8924791::init(), ATLAS_2011_I925932::init(), ATLAS_2011_I919017::init(), ATLAS_2011_CONF_2011_098::init(), ATLAS_2011_CONF_2011_090::init(), ATLAS_2010_S8919674::init(), ATLAS_2010_S8918562::init(), ATLAS_2010_S8914702::init(), ATLAS_2010_S8894728::init(), ATLAS_2010_S8817804::init(), ATLAS_2010_S8591806::init(), ATLAS_2010_CONF_2010_049::init(), ALICE_2011_S8945144::init(), ALICE_2011_S8909580::init(), ALICE_2010_S8706239::init(), ALICE_2010_S8625980::init(), ALICE_2010_S8624100::init(), ALEPH_2004_S5765862::init(), ALEPH_1996_S3486095::init(), ALEPH_1996_S3196992::init(), ALEPH_1991_S2435284::init(), IsolationProjection< PROJ1, PROJ2, EST >::IsolationProjection(), JetAlg::JetAlg(), JetShape::JetShape(), LeadingParticlesFinalState::LeadingParticlesFinalState(), LeptonClusters::LeptonClusters(), LossyFinalState< ConstRandomFilter >::LossyFinalState(), MergedFinalState::MergedFinalState(), MissingMomentum::MissingMomentum(), Multiplicity::Multiplicity(), NeutralFinalState::NeutralFinalState(), NonHadronicFinalState::NonHadronicFinalState(), ParisiTensor::ParisiTensor(), Sphericity::Sphericity(), Spherocity::Spherocity(), SVertex::SVertex(), Thrust::Thrust(), TotalVisibleMomentum::TotalVisibleMomentum(), TriggerCDFRun0Run1::TriggerCDFRun0Run1(), TriggerCDFRun2::TriggerCDFRun2(), TriggerUA5::TriggerUA5(), VetoedFinalState::VetoedFinalState(), and VisibleFinalState::VisibleFinalState().

00113                                                                        {
00114       const Projection& reg = _addProjection(proj, name);
00115       const PROJ& rtn = dynamic_cast<const PROJ&>(reg);
00116       return rtn;
00117     }

const PROJ& applyProjection ( const Event evt,
const std::string &  name 
) const [inline, inherited]

Apply the named projection on event.

Definition at line 81 of file ProjectionApplier.hh.

References ProjectionApplier::_applyProjection().

00081                                                                                {
00082       return pcast<PROJ>(_applyProjection(evt, name));
00083     }

const PROJ& applyProjection ( const Event evt,
const Projection proj 
) const [inline, inherited]

Apply the supplied projection on event.

Definition at line 74 of file ProjectionApplier.hh.

References ProjectionApplier::_applyProjection().

00074                                                                                 {
00075       return pcast<PROJ>(_applyProjection(evt, proj));
00076     }

const PROJ& applyProjection ( const Event evt,
const PROJ &  proj 
) const [inline, inherited]

Apply the supplied projection on event.

Definition at line 67 of file ProjectionApplier.hh.

References ProjectionApplier::_applyProjection().

Referenced by DISFinalState::project().

00067                                                                           {
00068       return pcast<PROJ>(_applyProjection(evt, proj));
00069     }

const set< PdgIdPair > beamPairs (  )  const [virtual, inherited]

Return the BeamConstraints for this projection, not including recursion. Derived classes should ensure that all contained projections are registered in the _projections set for the beam constraint chaining to work.

Definition at line 39 of file Projection.cc.

References Projection::_beamPairs, Projection::beamPairs(), Projection::getLog(), ProjectionApplier::getProjections(), Rivet::intersection(), and Log::TRACE.

Referenced by Projection::beamPairs().

00039                                                    {
00040     set<PdgIdPair> ret = _beamPairs;
00041     set<ConstProjectionPtr> projs = getProjections();
00042     for (set<ConstProjectionPtr>::const_iterator ip = projs.begin(); ip != projs.end(); ++ip) {
00043       ConstProjectionPtr p = *ip;
00044       getLog() << Log::TRACE << "Proj addr = " << p << endl;
00045       if (p) ret = intersection(ret, p->beamPairs());
00046     }
00047     return ret;
00048   }

bool before ( const Projection p  )  const [inherited]

Determine whether this object should be ordered before the object p given as argument. If p is of a different class than this, the before() function of the corresponding type_info objects is used. Otherwise, if the objects are of the same class, the virtual compare(const Projection &) will be returned.

Definition at line 28 of file Projection.cc.

References Projection::compare().

Referenced by less< const Rivet::Projection * >::operator()().

00028                                                    {
00029     const std::type_info& thisid = typeid(*this);
00030     const std::type_info& otherid = typeid(p);
00031     if (thisid == otherid) {
00032       return compare(p) < 0;
00033     } else {
00034       return thisid.before(otherid);
00035     }
00036   }

void calc ( const ParticleVector inparticles  ) 

Operate on a given particle vector directly instead of through project (no caching).

Definition at line 88 of file InvMassFinalState.cc.

References InvMassFinalState::_decayids, InvMassFinalState::_masstarget, InvMassFinalState::_maxmass, InvMassFinalState::_minmass, InvMassFinalState::_particlePairs, FinalState::_theParticles, InvMassFinalState::_useTransverseMass, FinalState::accept(), Particle::genParticle(), Projection::getLog(), Rivet::GeV, Rivet::inRange(), Log::isActive(), FourMomentum::mass(), FourMomentum::mass2(), FourMomentum::massT(), Particle::momentum(), MSG_DEBUG, MSG_TRACE, Particle::pdgId(), and Log::TRACE.

Referenced by WFinder::project(), and InvMassFinalState::project().

00088                                                                 {
00089     _theParticles.clear();
00090     _particlePairs.clear();
00091 
00092     // Containers for the particles of type specified in the pair
00093     vector<const Particle*> type1;
00094     vector<const Particle*> type2;
00095     // Get all the particles of the type specified in the pair from the particle list
00096     foreach (const Particle& ipart, inparticles) {
00097       // Loop around possible particle pairs (typedef needed to keep BOOST_FOREACH happy)
00098       foreach (const PdgIdPair& ipair, _decayids) {
00099         if (ipart.pdgId() == ipair.first) {
00100           if (accept(ipart)) {
00101             type1 += &ipart;
00102           }
00103         } else if (ipart.pdgId() == ipair.second) {
00104           if (accept(ipart)) {
00105             type2 += &ipart;
00106           }
00107         }
00108       }
00109     }
00110     if (type1.empty() || type2.empty()) return;
00111 
00112     // Temporary container of selected particles iterators
00113     // Useful to compare iterators and avoid double occurrences of the same
00114     // particle in case it matches with more than another particle
00115     vector<const Particle*> tmp;
00116 
00117     // Now calculate the inv mass
00118     pair<double, pair<Particle, Particle> > closestPair;
00119     closestPair.first = 1e30;
00120     foreach (const Particle* i1, type1) {
00121       foreach (const Particle* i2, type2) {
00122         // Check this is actually a pair
00123         // (if more than one pair in vector particles can be unrelated)
00124         bool found = false;
00125         foreach (const PdgIdPair& ipair, _decayids) {
00126           if (i1->pdgId() == ipair.first && i2->pdgId() == ipair.second) {
00127             found = true;
00128             break;
00129           }
00130         }
00131         if (!found) continue;
00132 
00133         FourMomentum v4 = i1->momentum() + i2->momentum();
00134         if (v4.mass2() < 0) {
00135           MSG_DEBUG("Constructed negative inv mass2: skipping!");
00136           continue;
00137         }
00138         bool passedMassCut = false;
00139         if (_useTransverseMass) {
00140           passedMassCut = inRange(v4.massT(), _minmass, _maxmass);
00141         } else {
00142           passedMassCut = inRange(v4.mass(), _minmass, _maxmass);
00143         }
00144 
00145         if (passedMassCut) {
00146           MSG_DEBUG("Selecting particles with IDs " << i1->pdgId() << " & " << i2->pdgId()
00147                     << " and mass = " << v4.mass()/GeV << " GeV");
00148           // Store accepted particles, avoiding duplicates
00149           if (find(tmp.begin(), tmp.end(), i1) == tmp.end()) {
00150             tmp.push_back(i1);
00151             _theParticles += *i1;
00152           }
00153           if (find(tmp.begin(), tmp.end(), i2) == tmp.end()) {
00154             tmp.push_back(i2);
00155             _theParticles += *i2;
00156           }
00157           // Store accepted particle pairs
00158           _particlePairs += make_pair(*i1, *i2);
00159           if (_masstarget>0.0) {
00160             double diff=fabs(v4.mass()-_masstarget);
00161             if (diff<closestPair.first) {
00162               closestPair.first=diff;
00163               closestPair.second=make_pair(*i1, *i2);
00164             }
00165           }
00166         }
00167       }
00168     }
00169     if (_masstarget > 0.0 && closestPair.first < 1e30) {
00170       _theParticles.clear();
00171       _particlePairs.clear();
00172       _theParticles += closestPair.second.first;
00173       _theParticles += closestPair.second.second;
00174       _particlePairs += closestPair.second;
00175     }
00176 
00177     MSG_DEBUG("Selected " << _theParticles.size() << " particles " << "(" << _particlePairs.size() << " pairs)");
00178     if (getLog().isActive(Log::TRACE)) {
00179       foreach (const Particle& p, _theParticles) {
00180         MSG_TRACE("ID: " << p.pdgId() << ", barcode: " << p.genParticle().barcode());
00181       }
00182     }
00183   }

virtual const Projection* clone (  )  const [inline, virtual]

Clone on the heap.

Reimplemented from FinalState.

Definition at line 42 of file InvMassFinalState.hh.

References InvMassFinalState::InvMassFinalState().

00042                                             {
00043         return new InvMassFinalState(*this);
00044     }

int compare ( const Projection p  )  const [protected, virtual]

Compare projections.

Reimplemented from FinalState.

Definition at line 55 of file InvMassFinalState.cc.

References InvMassFinalState::_decayids, InvMassFinalState::_maxmass, InvMassFinalState::_minmass, InvMassFinalState::_useTransverseMass, Rivet::cmp(), Rivet::EQUIVALENT, and Projection::mkNamedPCmp().

00055                                                           {
00056     // First compare the final states we are running on
00057     int fscmp = mkNamedPCmp(p, "FS");
00058     if (fscmp != EQUIVALENT) return fscmp;
00059 
00060     // Then compare the two as final states
00061     const InvMassFinalState& other = dynamic_cast<const InvMassFinalState&>(p);
00062     fscmp = FinalState::compare(other);
00063     if (fscmp != EQUIVALENT) return fscmp;
00064 
00065     // Compare the mass limits
00066     int masstypecmp = cmp(_useTransverseMass, other._useTransverseMass);
00067     if (masstypecmp != EQUIVALENT) return masstypecmp;
00068     int massllimcmp = cmp(_minmass, other._minmass);
00069     if (massllimcmp != EQUIVALENT) return massllimcmp;
00070     int masshlimcmp = cmp(_maxmass, other._maxmass);
00071     if (masshlimcmp != EQUIVALENT) return masshlimcmp;
00072 
00073     // Compare the decay species
00074     int decaycmp = cmp(_decayids, other._decayids);
00075     if (decaycmp != EQUIVALENT) return decaycmp;
00076 
00077     // Finally compare them as final states
00078     return FinalState::compare(other);
00079   }

virtual bool empty (  )  const [inline, virtual, inherited]
const collection_type& entities (  )  const [inline, inherited]

Template-usable interface common to JetAlg.

Definition at line 96 of file FinalState.hh.

References FinalState::particles().

00096                                             {
00097       return particles();
00098     }

Log& getLog (  )  const [inline, inherited]

Get a Log object based on the getName() property of the calling projection object.

Reimplemented from ProjectionApplier.

Definition at line 114 of file Projection.hh.

References Projection::name().

Referenced by Projection::beamPairs(), InvMassFinalState::calc(), VetoedFinalState::project(), UnstableFinalState::project(), LossyFinalState< ConstRandomFilter >::project(), IsolationProjection< PROJ1, PROJ2, EST >::project(), InitialQuarks::project(), ChargedFinalState::project(), and TotalVisibleMomentum::TotalVisibleMomentum().

00114                         {
00115       string logname = "Rivet.Projection." + name();
00116       return Log::getLog(logname);
00117     }

const Projection& getProjection ( const std::string &  name  )  const [inline, inherited]

Get the named projection (non-templated, so returns as a reference to a Projection base class).

Definition at line 57 of file ProjectionApplier.hh.

References ProjectionHandler::getProjection(), and ProjectionApplier::getProjHandler().

00057                                                                  {
00058       return getProjHandler().getProjection(*this, name);
00059     }

const PROJ& getProjection ( const std::string &  name  )  const [inline, inherited]

Get the named projection, specifying return type via a template argument.

Definition at line 49 of file ProjectionApplier.hh.

References ProjectionHandler::getProjection(), and ProjectionApplier::getProjHandler().

Referenced by ProjectionApplier::_applyProjection(), Rivet::pcmp(), and Hemispheres::project().

00049                                                            {
00050       const Projection& p = getProjHandler().getProjection(*this, name);
00051       return pcast<PROJ>(p);
00052     }

std::set<ConstProjectionPtr> getProjections (  )  const [inline, inherited]

Get the contained projections, including recursion.

Definition at line 42 of file ProjectionApplier.hh.

References ProjectionHandler::DEEP, ProjectionHandler::getChildProjections(), and ProjectionApplier::getProjHandler().

Referenced by Projection::beamPairs().

00042                                                       {
00043       return getProjHandler().getChildProjections(*this, ProjectionHandler::DEEP);
00044     }

ProjectionHandler& getProjHandler (  )  const [inline, protected, inherited]
virtual bool isEmpty (  )  const [inline, virtual, inherited]
Deprecated:
Is this final state empty?

Reimplemented in UnstableFinalState.

Definition at line 84 of file FinalState.hh.

References FinalState::_theParticles.

00084 { return _theParticles.empty(); }

Cmp< Projection > mkNamedPCmp ( const Projection otherparent,
const std::string &  pname 
) const [protected, inherited]
Cmp< Projection > mkPCmp ( const Projection otherparent,
const std::string &  pname 
) const [protected, inherited]

Shortcut to make a named Cmp<Projection> comparison with the *this object automatically passed as one of the parent projections.

Definition at line 57 of file Projection.cc.

References Rivet::pcmp().

00058                                                                 {
00059     return pcmp(*this, otherparent, pname);
00060   }

virtual std::string name (  )  const [inline, virtual, inherited]
const std::vector< std::pair< Particle, Particle > > & particlePairs (  )  const

Constituent pairs.

Definition at line 187 of file InvMassFinalState.cc.

References InvMassFinalState::_particlePairs.

Referenced by WFinder::project().

00187                                                                                         {
00188     return _particlePairs;
00189   }

const ParticleVector& particles ( sorter  )  const [inline, inherited]

Get the final-state particles, ordered by supplied sorting function object.

Definition at line 43 of file FinalState.hh.

References FinalState::_theParticles.

00043                                                     {
00044       std::sort(_theParticles.begin(), _theParticles.end(), sorter);
00045       return _theParticles;
00046     }

virtual const ParticleVector& particles (  )  const [inline, virtual, inherited]

Get the final-state particles.

Reimplemented in UnstableFinalState.

Definition at line 39 of file FinalState.hh.

References FinalState::_theParticles.

Referenced by UA5_1986_S1583476::analyze(), UA5_1982_S875503::analyze(), UA1_1990_S2044935::analyze(), TASSO_1990_S2148048::analyze(), STAR_2009_UE_HELEN::analyze(), STAR_2008_S7993412::analyze(), STAR_2008_S7869363::analyze(), STAR_2006_S6500200::analyze(), SFM_1984_S1178091::analyze(), PDG_HADRON_MULTIPLICITIES_RATIOS::analyze(), PDG_HADRON_MULTIPLICITIES::analyze(), OPAL_1998_S3780481::analyze(), OPAL_1993_S2692198::analyze(), MC_VH2BB::analyze(), MC_SUSY::analyze(), MC_PHOTONJETUE::analyze(), MC_PHOTONJETS::analyze(), MC_LEADINGJETS::analyze(), MC_GENERIC::analyze(), MC_DIJET::analyze(), JADE_1998_S3612880::analyze(), H1_2000_S4129130::analyze(), H1_1995_S3167097::analyze(), H1_1994_S2919893::analyze(), E735_1998_S3905616::analyze(), DELPHI_2003_WUD_03_11::analyze(), DELPHI_2002_069_CONF_603::analyze(), DELPHI_1996_S3430090::analyze(), DELPHI_1995_S3137023::analyze(), D0_2008_S7719523::analyze(), D0_2006_S6438750::analyze(), D0_2001_S4674421::analyze(), CMS_2011_S8884919::analyze(), CMS_2010_S8656010::analyze(), CMS_2010_S8547297::analyze(), CDF_2010_S8591881_QCD::analyze(), CDF_2010_S8591881_DY::analyze(), CDF_2009_S8233977::analyze(), CDF_2008_S8095620::analyze(), CDF_2008_S7541902::analyze(), CDF_2008_S7540469::analyze(), CDF_2008_NOTE_9351::analyze(), CDF_2008_LEADINGJETS::analyze(), CDF_2006_S6653332::analyze(), CDF_1990_S2089246::analyze(), CDF_1988_S1865951::analyze(), ATLAS_2011_S9002537::analyze(), ATLAS_2011_I925932::analyze(), ATLAS_2010_S8919674::analyze(), ATLAS_2010_S8894728::analyze(), ATLAS_2010_S8591806::analyze(), ALICE_2011_S8945144::analyze(), ALICE_2010_S8706239::analyze(), ALICE_2010_S8625980::analyze(), ALEPH_2004_S5765862::analyze(), ALEPH_1996_S3486095::analyze(), Thrust::calc(), Spherocity::calc(), Sphericity::calc(), FParameter::calc(), BeamThrust::calc(), FinalState::entities(), ATLAS_2010_S8918562::fillPtEtaNch(), FinalState::particlesByE(), FinalState::particlesByEt(), FinalState::particlesByEta(), FinalState::particlesByModEta(), FinalState::particlesByP(), FinalState::particlesByPt(), VisibleFinalState::project(), VetoedFinalState::project(), TriggerUA5::project(), TriggerCDFRun2::project(), TriggerCDFRun0Run1::project(), TotalVisibleMomentum::project(), SVertex::project(), NonHadronicFinalState::project(), NeutralFinalState::project(), Multiplicity::project(), MissingMomentum::project(), MergedFinalState::project(), LossyFinalState< ConstRandomFilter >::project(), LeptonClusters::project(), LeadingParticlesFinalState::project(), InvMassFinalState::project(), IdentifiedFinalState::project(), Hemispheres::project(), HadronicFinalState::project(), FoxWolframMoments::project(), FinalState::project(), DISLepton::project(), ClusteredPhotons::project(), ChargedLeptons::project(), ChargedFinalState::project(), and CentralEtHCM::project().

00039 { return _theParticles; }

const ParticleVector& particlesByE (  )  const [inline, inherited]

Get the final-state particles, ordered by decreasing $ E $.

Definition at line 59 of file FinalState.hh.

References Rivet::cmpParticleByE(), and FinalState::particles().

00059                                                {
00060       return particles(cmpParticleByE);
00061     }

const ParticleVector& particlesByEt (  )  const [inline, inherited]

Get the final-state particles, ordered by decreasing $ E_T $.

Definition at line 64 of file FinalState.hh.

References Rivet::cmpParticleByEt(), and FinalState::particles().

00064                                                 {
00065       return particles(cmpParticleByEt);
00066     }

const ParticleVector& particlesByEta (  )  const [inline, inherited]

Get the final-state particles, ordered by increasing $ \eta $.

Definition at line 69 of file FinalState.hh.

References Rivet::cmpParticleByAscPseudorapidity(), and FinalState::particles().

00069                                                  {
00070       return particles(cmpParticleByAscPseudorapidity);
00071     }

const ParticleVector& particlesByModEta (  )  const [inline, inherited]

Get the final-state particles, ordered by increasing $ |\eta| $.

Definition at line 74 of file FinalState.hh.

References Rivet::cmpParticleByAscAbsPseudorapidity(), and FinalState::particles().

00074                                                     {
00075       return particles(cmpParticleByAscAbsPseudorapidity);
00076     }

const ParticleVector& particlesByP (  )  const [inline, inherited]

Get the final-state particles, ordered by decreasing $ p $.

Definition at line 54 of file FinalState.hh.

References Rivet::cmpParticleByP(), and FinalState::particles().

00054                                                {
00055       return particles(cmpParticleByP);
00056     }

const ParticleVector& particlesByPt (  )  const [inline, inherited]

Get the final-state particles, ordered by decreasing $ p_T $.

Definition at line 49 of file FinalState.hh.

References Rivet::cmpParticleByPt(), and FinalState::particles().

Referenced by MC_PHOTONJETUE::analyze(), ATLAS_2011_S8994773::analyze(), and ATLAS_2010_S8894728::analyze().

00049                                                 {
00050       return particles(cmpParticleByPt);
00051     }

void project ( const Event e  )  [protected, virtual]

Apply the projection on the supplied event.

Reimplemented from FinalState.

Definition at line 82 of file InvMassFinalState.cc.

References InvMassFinalState::calc(), and FinalState::particles().

00082                                                 {
00083     const FinalState& fs = applyProjection<FinalState>(e, "FS");
00084     calc(fs.particles());
00085   }

virtual double ptMin (  )  const [inline, virtual, inherited]

Minimum-$ p_\perp $ requirement.

Definition at line 87 of file FinalState.hh.

References FinalState::_ptmin.

00087 { return _ptmin; }

void setName ( const std::string &  name  )  [inline, inherited]

Used by derived classes to set their name.

Definition at line 120 of file Projection.hh.

References Projection::_name.

Referenced by ZFinder::_init(), WFinder::_init(), FastJets::_init1(), FastJets::_init2(), FastJets::_init3(), Beam::Beam(), BeamThrust::BeamThrust(), CentralEtHCM::CentralEtHCM(), ChargedFinalState::ChargedFinalState(), ChargedLeptons::ChargedLeptons(), ClusteredPhotons::ClusteredPhotons(), ConstLossyFinalState::ConstLossyFinalState(), DISFinalState::DISFinalState(), DISKinematics::DISKinematics(), DISLepton::DISLepton(), FinalState::FinalState(), FoxWolframMoments::FoxWolframMoments(), FParameter::FParameter(), HadronicFinalState::HadronicFinalState(), Hemispheres::Hemispheres(), IdentifiedFinalState::IdentifiedFinalState(), InitialQuarks::InitialQuarks(), IsolationProjection< PROJ1, PROJ2, EST >::IsolationProjection(), JetAlg::JetAlg(), JetShape::JetShape(), LeadingParticlesFinalState::LeadingParticlesFinalState(), LeptonClusters::LeptonClusters(), LossyFinalState< ConstRandomFilter >::LossyFinalState(), MergedFinalState::MergedFinalState(), MissingMomentum::MissingMomentum(), Multiplicity::Multiplicity(), NeutralFinalState::NeutralFinalState(), NonHadronicFinalState::NonHadronicFinalState(), ParisiTensor::ParisiTensor(), PVertex::PVertex(), Sphericity::Sphericity(), Spherocity::Spherocity(), SVertex::SVertex(), Thrust::Thrust(), TotalVisibleMomentum::TotalVisibleMomentum(), TriggerCDFRun0Run1::TriggerCDFRun0Run1(), TriggerCDFRun2::TriggerCDFRun2(), TriggerUA5::TriggerUA5(), UnstableFinalState::UnstableFinalState(), VetoedFinalState::VetoedFinalState(), and VisibleFinalState::VisibleFinalState().

00120                                         {
00121       _name = name;
00122     }

virtual size_t size (  )  const [inline, virtual, inherited]
void useTransverseMass ( bool  usetrans = true  )  [inline]

Choose whether to use the full inv mass or just the transverse mass.

Definition at line 54 of file InvMassFinalState.hh.

References InvMassFinalState::_useTransverseMass.

Referenced by WFinder::project().

00054                                                {
00055       _useTransverseMass = usetrans;
00056     }


Friends And Related Function Documentation

friend class Cmp< Projection > [friend, inherited]

The Cmp specialization for Projection is a friend.

Definition at line 36 of file Projection.hh.

friend class Event [friend, inherited]

Event is a friend.

Definition at line 33 of file Projection.hh.


Member Data Documentation

bool _allowProjReg [protected, inherited]

Flag to forbid projection registration in analyses until the init phase.

Definition at line 140 of file ProjectionApplier.hh.

Referenced by ProjectionApplier::_addProjection(), and Analysis::Analysis().

std::vector<PdgIdPair> _decayids [private]

IDs of the decay products.

Definition at line 74 of file InvMassFinalState.hh.

Referenced by InvMassFinalState::calc(), and InvMassFinalState::compare().

vector<pair<double,double> > _etaRanges [protected, inherited]

The ranges allowed for pseudorapidity.

Definition at line 116 of file FinalState.hh.

Referenced by FinalState::accept(), FinalState::compare(), FinalState::FinalState(), and FinalState::project().

double _masstarget [private]

Target mass if only one pair should be returned.

Definition at line 86 of file InvMassFinalState.hh.

Referenced by InvMassFinalState::calc().

double _maxmass [private]

Max inv mass.

Definition at line 83 of file InvMassFinalState.hh.

Referenced by InvMassFinalState::calc(), and InvMassFinalState::compare().

double _minmass [private]

Min inv mass.

Definition at line 80 of file InvMassFinalState.hh.

Referenced by InvMassFinalState::calc(), and InvMassFinalState::compare().

std::vector<std::pair<Particle, Particle> > _particlePairs [private]

Constituent pairs.

Definition at line 77 of file InvMassFinalState.hh.

Referenced by InvMassFinalState::calc(), and InvMassFinalState::particlePairs().

double _ptmin [protected, inherited]

The minimum allowed transverse momentum.

Reimplemented in UnstableFinalState.

Definition at line 119 of file FinalState.hh.

Referenced by FinalState::accept(), FinalState::compare(), FinalState::project(), and FinalState::ptMin().

ParticleVector _theParticles [mutable, protected, inherited]
bool _useTransverseMass [private]

Flag to decide whether to use the full inv mass or just the transverse mass.

Definition at line 89 of file InvMassFinalState.hh.

Referenced by InvMassFinalState::calc(), InvMassFinalState::compare(), and InvMassFinalState::useTransverseMass().


The documentation for this class was generated from the following files: