Public Types |
Public Member Functions |
Protected Member Functions |
Protected Attributes |
Private Member Functions |
Private Attributes |
Friends
FastJets Class Reference Project out jets found using the FastJet package jet algorithms. More...
Inheritance diagram for FastJets:
![]()
Collaboration diagram for FastJets:
![]()
Detailed DescriptionProject out jets found using the FastJet package jet algorithms. Definition at line 26 of file FastJets.hh. Member Typedef Documentation
Member Enumeration Documentation
Enum for the treatment of invisible particles: whether to include all, some, or none in jet-finding.
Wrapper enum for selected Fastjet jet algorithms.
Definition at line 30 of file FastJets.hh.
Enum for the treatment of muons: whether to include all, some, or none in jet-finding. Definition at line 19 of file JetAlg.hh. { NO_MUONS, DECAY_MUONS, ALL_MUONS }; Constructor & Destructor Documentation
"Wrapped" argument constructor using Rivet enums for most common jet alg choices (including some plugins). For the built-in algs, E-scheme recombination is used. For full control of FastJet built-in jet algs, use the native arg constructor. Definition at line 44 of file FastJets.hh. References FastJets::_init1(). Referenced by FastJets::clone().
Native argument constructor, using FastJet alg/scheme enums. Definition at line 51 of file FastJets.hh. References FastJets::_init2().
Explicitly pass in a JetDefinition. Definition at line 58 of file FastJets.hh. References FastJets::_init3().
Explicitly pass in an externally-constructed plugin (must be heap-allocated, Rivet will delete) Definition at line 64 of file FastJets.hh. References FastJets::_init4().
Explicitly pass in an externally-constructed plugin (must be heap-allocated, Rivet will delete)
Definition at line 71 of file FastJets.hh. References FastJets::_init4().
Same thing as above, but without an FS (for when we want to pass the particles directly to the calc method) Definition at line 78 of file FastJets.hh. References FastJets::_init1().
Same thing as above, but without an FS (for when we want to pass the particles directly to the calc method) Definition at line 81 of file FastJets.hh. References FastJets::_init2().
Same thing as above, but without an FS (for when we want to pass the particles directly to the calc method) Definition at line 84 of file FastJets.hh. References FastJets::_init3().
Same thing as above, but without an FS (for when we want to pass the particles directly to the calc method) Definition at line 87 of file FastJets.hh. References FastJets::_init3(). Member Function Documentation
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(). { if (!_allowProjReg) { cerr << "Trying to register projection '" << proj.name() << "' before init phase in '" << this->name() << "'." << endl; exit(2); } const Projection& reg = getProjHandler().registerProjection(*this, proj, name); return reg; }
Definition at line 17 of file FastJets.cc. References FastJets::_initBase(), FastJets::_jdef, FastJets::_plugin, FastJets::ANTIKT, FastJets::ATLASCONE, FastJets::CAM, FastJets::CDFJETCLU, FastJets::CDFMIDPOINT, FastJets::CMSCONE, FastJets::D0ILCONE, FastJets::DURHAM, FastJets::JADE, FastJets::KT, MSG_DEBUG, FastJets::SISCONE, and FastJets::TRACKJET. Referenced by FastJets::FastJets(). { _initBase(); MSG_DEBUG("JetAlg = " << alg); MSG_DEBUG("R parameter = " << rparameter); MSG_DEBUG("Seed threshold = " << seed_threshold); if (alg == KT) { _jdef = fastjet::JetDefinition(fastjet::kt_algorithm, rparameter, fastjet::E_scheme); } else if (alg == CAM) { _jdef = fastjet::JetDefinition(fastjet::cambridge_algorithm, rparameter, fastjet::E_scheme); } else if (alg == ANTIKT) { _jdef = fastjet::JetDefinition(fastjet::antikt_algorithm, rparameter, fastjet::E_scheme); } else if (alg == DURHAM) { _jdef = fastjet::JetDefinition(fastjet::ee_kt_algorithm, fastjet::E_scheme); } else { // Plugins: if (alg == SISCONE) { const double OVERLAP_THRESHOLD = 0.75; _plugin.reset(new fastjet::SISConePlugin(rparameter, OVERLAP_THRESHOLD)); // } else if (alg == PXCONE) { // string msg = "PxCone currently not supported, since FastJet doesn't install it by default. "; // msg += "Please notify the Rivet authors if this behaviour should be changed."; // throw Error(msg); //_plugin.reset(new fastjet::PxConePlugin(rparameter)); } else if (alg == ATLASCONE) { const double OVERLAP_THRESHOLD = 0.5; _plugin.reset(new fastjet::ATLASConePlugin(rparameter, seed_threshold, OVERLAP_THRESHOLD)); } else if (alg == CMSCONE) { _plugin.reset(new fastjet::CMSIterativeConePlugin(rparameter, seed_threshold)); } else if (alg == CDFJETCLU) { const double OVERLAP_THRESHOLD = 0.75; _plugin.reset(new fastjet::CDFJetCluPlugin(rparameter, OVERLAP_THRESHOLD, seed_threshold)); } else if (alg == CDFMIDPOINT) { const double OVERLAP_THRESHOLD = 0.5; _plugin.reset(new fastjet::CDFMidPointPlugin(rparameter, OVERLAP_THRESHOLD, seed_threshold)); } else if (alg == D0ILCONE) { const double min_jet_Et = 6.0; _plugin.reset(new fastjet::D0RunIIConePlugin(rparameter, min_jet_Et)); } else if (alg == JADE) { _plugin.reset(new fastjet::JadePlugin()); } else if (alg == TRACKJET) { _plugin.reset(new fastjet::TrackJetPlugin(rparameter)); } _jdef = fastjet::JetDefinition(_plugin.get()); } }
Definition at line 63 of file FastJets.cc. References FastJets::_initBase(), and FastJets::_jdef. Referenced by FastJets::FastJets().
Definition at line 69 of file FastJets.cc. References FastJets::_initBase(), and FastJets::_jdef. Referenced by FastJets::FastJets().
Definition at line 74 of file FastJets.cc. References FastJets::_initBase(), FastJets::_jdef, and FastJets::_plugin. Referenced by FastJets::FastJets().
Shared utility functions to implement constructor behaviour
Definition at line 11 of file FastJets.cc. References ProjectionApplier::addProjection(), TauFinder::HADRONIC, and Projection::setName(). Referenced by FastJets::_init1(), FastJets::_init2(), FastJets::_init3(), and FastJets::_init4(). { setName("FastJets"); addProjection(HeavyHadrons(), "HFHadrons"); addProjection(TauFinder(TauFinder::HADRONIC), "Taus"); } Get the jets (unordered).
Implements JetAlg. Definition at line 175 of file FastJets.cc. References FastJets::_particles, FastJets::clusterSeq(), and FastJets::pseudojets(). { Jets rtn; foreach (const fastjet::PseudoJet& pj, pseudojets()) { assert(clusterSeq() != NULL); const PseudoJets parts = clusterSeq()->constituents(pj); vector<Particle> constituents, tags; constituents.reserve(parts.size()); foreach (const fastjet::PseudoJet& p, parts) { map<int, Particle>::const_iterator found = _particles.find(p.user_index()); // assert(found != _particles.end()); if (found == _particles.end() && p.is_pure_ghost()) continue; //< Pure FJ ghosts are ok assert(found != _particles.end()); //< Anything else must be known assert(found->first != 0); //< All mapping IDs are pos-def (particles) or neg-def (tags) if (found->first > 0) constituents.push_back(found->second); else if (found->first < 0) tags.push_back(found->second); } Jet j(pj, constituents, tags); rtn.push_back(j); } /// @todo Cache? return rtn; }
Add a colliding beam pair. Definition at line 107 of file Projection.hh. References Projection::_beamPairs. Referenced by Projection::Projection(). { _beamPairs.insert(PdgIdPair(beam1, beam2)); return *this; }
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 Definition at line 116 of file ProjectionApplier.hh. References ProjectionApplier::_addProjection(). Referenced by FastJets::_initBase(), VetoedFinalState::addVetoOnThisFinalState(), BeamThrust::BeamThrust(), CDF_2009_S8057893::CDF_2009_S8057893::init(), CentralEtHCM::CentralEtHCM(), ChargedFinalState::ChargedFinalState(), ChargedLeptons::ChargedLeptons(), DISFinalState::DISFinalState(), DISKinematics::DISKinematics(), DISLepton::DISLepton(), DressedLeptons::DressedLeptons(), FinalState::FinalState(), FoxWolframMoments::FoxWolframMoments(), FParameter::FParameter(), HadronicFinalState::HadronicFinalState(), HeavyHadrons::HeavyHadrons(), Hemispheres::Hemispheres(), IdentifiedFinalState::IdentifiedFinalState(), CMS_2010_S8547297::init(), CMS_2010_S8656010::init(), CMS_2015_I1327224::init(), ATLAS_2010_S8894728::init(), ALICE_2012_I1181770::init(), ATLAS_2011_S8994773::init(), CMS_2011_S8950903::init(), CMS_2012_PAS_QCD_11_010::init(), CMS_2015_I1356998::init(), CMS_2011_S8968497::init(), CMS_2011_S8973270::init(), CMS_2012_I1184941::init(), CMS_2012_I1090423::init(), CMS_2012_I1193338::init(), ATLAS_2014_I1298811::init(), CMSTOTEM_2014_I1294140::init(), ALICE_2011_S8909580::init(), ATLAS_2011_I894867::init(), ALICE_2014_I1300380::init(), LHCB_2013_I1208105::init(), ATLAS_2010_CONF_2010_049::init(), CMS_2011_S8941262::init(), CDF_2007_S7057202::init(), TOTEM_2012_I1115294::init(), TOTEM_2014_I1328627::init(), CMS_2012_I1087342::init(), CMS_2011_S9086218::init(), CMS_2011_S9215166::init(), ATLAS_2014_I1282441::init(), ALICE_2011_S8945144::init(), ATLAS_2015_I1387176::init(), ALICE_2015_I1357424::init(), ATLAS_2012_I1091481::init(), CMS_2015_I1346843::init(), CMS_2015_I1384119::init(), D0_2011_I895662::init(), CMS_2011_S8957746::init(), LHCF_2012_I1115479::init(), MC_ELECTRONS::init(), MC_MUONS::init(), MC_TAUS::init(), ATLAS_2011_I930220::init(), ATLAS_2011_S9002537::init(), CMS_2011_I954992::init(), TOTEM_2012_002::init(), CMS_2011_S8978280::init(), ATLAS_2010_S8591806::init(), CMS_2011_S9088458::init(), CMS_2013_I1265659::init(), MC_JETTAGS::init(), ATLAS_2012_I1188891::init(), CDF_2012_NOTE10874::init(), CMS_2013_I1208923::init(), CMS_2013_I1273574::init(), D0_2010_S8570965::init(), ATLAS_2011_I925932::init(), CDF_1997_S3541940::init(), ATLAS_2012_I1124167::init(), STAR_2006_S6500200::init(), STAR_2008_S7993412::init(), UA5_1987_S1640666::init(), CMS_2013_I1256943::init(), CMS_QCD_10_024::init(), CDF_1993_S2742446::init(), ARGUS_1993_S2789213::init(), MC_HINC::init(), CDF_2000_S4155203::init(), MC_JETS::init(), MC_KTSPLITTINGS::init(), CDF_2006_S6450792::init(), CDF_2005_S6080774::init(), SFM_1984_S1178091::init(), CMS_2013_I1258128::init(), CMS_2013_I1261026::init(), CMS_2014_I1298810::init(), ATLAS_2011_S9128077::init(), CMS_2015_I1310737::init(), D0_2000_S4480767::init(), ATLAS_2014_I1319490::init(), BABAR_2007_S7266081::init(), BELLE_2008_I786560::init(), LHCB_2012_I1208102::init(), MC_HJETS::init(), ARGUS_1993_S2669951::init(), MC_WWINC::init(), CDF_2008_S7782535::init(), MC_ZZINC::init(), CDF_2008_S8093652::init(), STAR_2006_S6870392::init(), ATLAS_2012_I1204447::init(), ATLAS_2013_I1190187::init(), UA5_1989_S1926373::init(), CMS_2013_I1209721::init(), CMS_2015_I1385107::init(), D0_2008_S6879055::init(), UA5_1982_S875503::init(), CMS_2012_I1102908::init(), ATLAS_2013_I1219109::init(), ATLAS_2013_I1243871::init(), ATLAS_2014_I1268975::init(), CMS_2013_I1272853::init(), D0_2000_I499943::init(), ATLAS_2014_I1325553::init(), ATLAS_2014_I1327229::init(), ATLAS_2013_I1216670::init(), E735_1998_S3905616::init(), MC_DIPHOTON::init(), MC_DIJET::init(), ATLAS_2011_I944826::init(), MC_WWKTSPLITTINGS::init(), MC_ZZJETS::init(), ALEPH_1991_S2435284::init(), ATLAS_2012_I1183818::init(), PDG_TAUS::init(), ATLAS_2012_I1204784::init(), MC_ZZKTSPLITTINGS::init(), ATLAS_2011_S9035664::init(), CDF_1988_S1865951::init(), CDF_1990_S2089246::init(), MC_HKTSPLITTINGS::init(), MC_PHOTONINC::init(), ATLAS_2011_S8924791::init(), CDF_2005_S6217184::init(), MC_WKTSPLITTINGS::init(), MC_ZINC::init(), OPAL_1995_S3198391::init(), OPAL_1996_S3257789::init(), OPAL_1997_S3608263::init(), OPAL_1998_S3702294::init(), OPAL_2000_S4418603::init(), MC_ZKTSPLITTINGS::init(), ALEPH_1996_S3196992::init(), ATLAS_2010_S8817804::init(), OPAL_1998_S3749908::init(), ALEPH_2002_S4823664::init(), D0_1996_S3324664::init(), D0_2007_S7075677::init(), D0_2009_S8202443::init(), D0_2010_S8821313::init(), DELPHI_1999_S3960137::init(), EXAMPLE_CUTS::init(), UA5_1986_S1583476::init(), CDF_1994_S2952106::init(), MC_IDENTIFIED::init(), ATLAS_2011_I945498::init(), MC_PHOTONJETS::init(), MC_PHOTONKTSPLITTINGS::init(), ATLAS_2011_I954993::init(), MC_WINC::init(), CMS_2012_I1107658::init(), MC_ZJETS::init(), D0_2008_S7554427::init(), D0_2008_S7863608::init(), MC_WWJETS::init(), D0_2010_S8671338::init(), H1_1995_S3167097::init(), ATLAS_2011_S9131140::init(), MC_TTBAR::init(), LHCB_2011_I919315::init(), ATLAS_2012_I1180197::init(), MC_WJETS::init(), MC_LEADJETUE::init(), UA1_1990_S2044935::init(), ATLAS_2014_I1307756::init(), ATLAS_2014_I1306615::init(), CDF_2008_S7540469::init(), CDF_2008_S7828950::init(), ATLAS_2011_S9212353::init(), ATLAS_2012_I1119557::init(), ATLAS_2014_I1300647::init(), D0_2001_S4674421::init(), ATLAS_2012_CONF_2012_104::init(), MC_GENERIC::init(), D0_2008_S7837160::init(), ATLAS_2012_CONF_2012_105::init(), ZEUS_2001_S4815815::init(), ATLAS_2012_CONF_2012_103::init(), CMS_2014_I1303894::init(), ATLAS_2012_I1126136::init(), ATLAS_2011_CONF_2011_098::init(), ALICE_2010_S8624100::init(), STAR_2006_S6860818::init(), D0_1996_S3214044::init(), ALEPH_2004_S5765862::init(), STAR_2009_UE_HELEN::init(), ATLAS_2013_I1217863_Z::init(), ATLAS_2012_CONF_2012_109::init(), ATLAS_2011_I926145::init(), ATLAS_2013_I1244522::init(), DELPHI_1995_S3137023::init(), MC_QCD_PARTONS::init(), OPAL_1997_S3396100::init(), ATLAS_2012_I1186556::init(), ATLAS_2015_I1345452::init(), ATLAS_2013_I1217863_W::init(), ALEPH_1999_S4193598::init(), CMS_2012_I941555::init(), ATLAS_2012_I1095236::init(), CDF_2008_S8095620::init(), JADE_OPAL_2000_S4300807::init(), CDF_1996_S3108457::init(), EXAMPLE::init(), ATLAS_2011_S8983313::init(), JADE_1998_S3612880::init(), LHCB_2013_I1218996::init(), CDF_1996_S3349578::init(), MC_HFJETS::init(), CDF_2009_NOTE_9936::init(), ALICE_2010_S8625980::init(), ALICE_2010_S8706239::init(), ATLAS_2012_I1125961::init(), CMS_2013_I1218372::init(), ATLAS_2012_CONF_2012_001::init(), ATLAS_2012_I1125575::init(), ATLAS_2012_I1190891::init(), ATLAS_2012_I1112263::init(), ATLAS_2011_S9212183::init(), ATLAS_2014_I1288706::init(), ATLAS_2013_I1230812::init(), CDF_2000_S4266730::init(), CDF_2001_S4563131::init(), D0_2009_S8349509::init(), D0_2009_S8320160::init(), ATLAS_2013_I1263495::init(), ATLAS_2011_S8971293::init(), MC_SUSY::init(), CDF_1998_S3618439::init(), CDF_2009_S8436959::init(), CDF_2009_S8383952::init(), ATLAS_2014_I1298023::init(), ALEPH_2001_S4656318::init(), ATLAS_2015_CONF_2015_041::init(), DELPHI_2002_069_CONF_603::init(), ATLAS_2014_I1307243::init(), CDF_2001_S4517016::init(), D0_2006_S6438750::init(), D0_2010_S8566488::init(), SLD_2002_S4869273::init(), ATLAS_2012_I1082936::init(), ATLAS_2010_S8919674::init(), SLD_1996_S3398250::init(), LHCB_2011_I917009::init(), UA5_1988_S1867512::init(), LHCB_2014_I1281685::init(), LHCB_2012_I1119400::init(), CDF_1996_S3418421::init(), ATLAS_2014_I1326641::init(), MC_PHOTONS::init(), ATLAS_2011_I921594::init(), ATLAS_2014_I1306294::init(), ATLAS_2011_S9108483::init(), CDF_2006_S6653332::init(), D0_2008_S7662670::init(), MC_WPOL::init(), CDF_2008_S7541902::init(), ATLAS_2014_I1304688::init(), ATLAS_2012_I946427::init(), ATLAS_2012_I1082009::init(), ATLAS_2011_S9225137::init(), DELPHI_2000_S4328825::init(), ATLAS_2011_S9019561::init(), ATLAS_2012_I1199269::init(), ATLAS_2012_I943401::init(), ATLAS_2012_I1117704::init(), TASSO_1990_S2148048::init(), ATLAS_2011_CONF_2011_090::init(), OPAL_2002_S5361494::init(), ATLAS_2012_I1084540::init(), CMS_2013_I1224539_WJET::init(), ATLAS_2012_I1083318::init(), ATLAS_2012_CONF_2012_153::init(), CMS_2013_I1224539_DIJET::init(), ATLAS_2013_I1217867::init(), CDF_2009_S8233977::init(), ATLAS_2015_I1364361::init(), CMS_2013_I1224539_ZJET::init(), CDF_2010_S8591881_DY::init(), CDF_2010_S8591881_QCD::init(), D0_2004_S5992206::init(), ATLAS_2010_S8914702::init(), ATLAS_2014_I1312627::init(), CDF_2001_S4751469::init(), D0_2008_S7719523::init(), ATLAS_2011_S9120807::init(), ATLAS_2010_S8918562::init(), ATLAS_2011_S9041966::init(), ATLAS_2012_I1094568::init(), DELPHI_1996_S3430090::init(), ALEPH_1996_S3486095::init(), ATLAS_2012_I1093738::init(), OPAL_2004_S6132243::init(), OPAL_1994_S2927284::init(), BELLE_2013_I1216515::init(), STAR_2008_S7869363::init(), BABAR_2007_S6895344::init(), ATLAS_2011_S9126244::init(), MC_VH2BB::init(), BABAR_2005_S6181155::init(), BELLE_2001_S4598261::init(), ATLAS_2011_I919017::init(), OPAL_2001_S4553896::init(), BABAR_2013_I1238276::init(), ATLAS_2012_I1203852::init(), CDF_2004_S5839831::init(), BABAR_2003_I593379::init(), ATLAS_2012_I1094061::init(), ARGUS_1993_S2653028::init(), CLEO_2004_S5809304::init(), OPAL_1998_S3780481::init(), ATLAS_2012_I1093734::init(), OPAL_1993_S2692198::init(), ATLAS_2014_I1279489::init(), ATLAS_2014_I1282447::init(), ATLAS_2012_I1094564::init(), H1_1994_S2919893::init(), H1_2000_S4129130::init(), SLD_2004_S5693039::init(), SLD_1999_S3743934::init(), PDG_HADRON_MULTIPLICITIES::init(), PDG_HADRON_MULTIPLICITIES_RATIOS::init(), JetAlg::JetAlg(), JetShape::JetShape(), LeadingParticlesFinalState::LeadingParticlesFinalState(), LossyFinalState< ConstRandomFilter >::LossyFinalState(), MergedFinalState::MergedFinalState(), MissingMomentum::MissingMomentum(), NeutralFinalState::NeutralFinalState(), NonHadronicFinalState::NonHadronicFinalState(), ParisiTensor::ParisiTensor(), PrimaryHadrons::PrimaryHadrons(), PromptFinalState::PromptFinalState(), Sphericity::Sphericity(), Spherocity::Spherocity(), TauFinder::TauFinder(), Thrust::Thrust(), TriggerCDFRun0Run1::TriggerCDFRun0Run1(), TriggerCDFRun2::TriggerCDFRun2(), TriggerUA5::TriggerUA5(), VetoedFinalState::VetoedFinalState(), VisibleFinalState::VisibleFinalState(), WFinder::WFinder(), and ZFinder::ZFinder(). { const Projection& reg = _addProjection(proj, name); const PROJ& rtn = dynamic_cast<const PROJ&>(reg); return rtn; }
Apply the supplied projection on event. Definition at line 70 of file ProjectionApplier.hh. References ProjectionApplier::_applyProjection(). Referenced by DISFinalState::project(). { return pcast<PROJ>(_applyProjection(evt, proj)); }
Apply the supplied projection on event. Definition at line 77 of file ProjectionApplier.hh. References ProjectionApplier::_applyProjection(). { return pcast<PROJ>(_applyProjection(evt, proj)); }
Apply the named projection on event. Definition at line 84 of file ProjectionApplier.hh. References ProjectionApplier::_applyProjection(). { return pcast<PROJ>(_applyProjection(evt, name)); }
Return the area definition.. May be null. Definition at line 176 of file FastJets.hh. References FastJets::_adef. Referenced by FastJets::clusterSeqArea(). { return _adef; } Return the allowed beam pairs on which this projection can operate, 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 35 of file Projection.cc. References Projection::_beamPairs, Projection::beamPairs(), Projection::getLog(), ProjectionApplier::getProjections(), Rivet::intersection(), and Log::TRACE. Referenced by Projection::beamPairs(). { set<PdgIdPair> ret = _beamPairs; set<ConstProjectionPtr> projs = getProjections(); for (set<ConstProjectionPtr>::const_iterator ip = projs.begin(); ip != projs.end(); ++ip) { ConstProjectionPtr p = *ip; getLog() << Log::TRACE << "Proj addr = " << p << endl; if (p) ret = intersection(ret, p->beamPairs()); } return ret; }
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 24 of file Projection.cc. References Projection::compare(). Referenced by less< const Rivet::Projection * >::operator()(). { const std::type_info& thisid = typeid(*this); const std::type_info& otherid = typeid(p); if (thisid == otherid) { return compare(p) < 0; } else { return thisid.before(otherid); } } Do the calculation locally (no caching).
<
< Ghostify the momentum <
Implements JetAlg. Definition at line 125 of file FastJets.cc. References FastJets::_adef, FastJets::_cseq, FastJets::_jdef, FastJets::_particles, FourMomentum::E(), Particle::momentum(), MSG_DEBUG, FourMomentum::px(), FourMomentum::py(), and FourMomentum::pz(). Referenced by CDF_2008_S7540469::analyze(), ATLAS_2015_I1364361::analyze(), and FastJets::project(). { _particles.clear(); vector<fastjet::PseudoJet> pjs; /// @todo Use FastJet3's UserInfo system // Store 4 vector data about each particle into FastJet's PseudoJets int counter = 1; foreach (const Particle& p, fsparticles) { const FourMomentum fv = p.momentum(); fastjet::PseudoJet pj(fv.px(), fv.py(), fv.pz(), fv.E()); ///< @todo Eliminate with implicit cast? pj.set_user_index(counter); pjs.push_back(pj); _particles[counter] = p; counter += 1; } // And the same for ghost tagging particles (with negative user indices) counter = 1; foreach (const Particle& p, tagparticles) { const FourMomentum fv = 1e-20 * p.momentum(); ///< Ghostify the momentum fastjet::PseudoJet pj(fv.px(), fv.py(), fv.pz(), fv.E()); ///< @todo Eliminate with implicit cast? pj.set_user_index(-counter); pjs.push_back(pj); _particles[-counter] = p; counter += 1; } MSG_DEBUG("Running FastJet ClusterSequence construction"); // Choose cseq as basic or area-calculating if (_adef == NULL) { _cseq.reset(new fastjet::ClusterSequence(pjs, _jdef)); } else { _cseq.reset(new fastjet::ClusterSequenceArea(pjs, _jdef, *_adef)); } }
Clone on the heap. Implements JetAlg. Definition at line 92 of file FastJets.hh. References FastJets::FastJets(). { return new FastJets(*this); }
Return the cluster sequence. Definition at line 160 of file FastJets.hh. References FastJets::_cseq. Referenced by FastJets::_jets(), MC_JetSplittings::analyze(), ALEPH_1996_S3196992::analyze(), JADE_OPAL_2000_S4300807::analyze(), JADE_1998_S3612880::analyze(), OPAL_2001_S4553896::analyze(), OPAL_2004_S6132243::analyze(), ALEPH_1996_S3486095::analyze(), ALEPH_2004_S5765862::analyze(), DELPHI_1996_S3430090::analyze(), FastJets::pseudoJets(), and ATLAS_2012_I1094564::splitjet(). { return _cseq.get(); }
Return the area-enabled cluster sequence (if an area defn exists, otherwise returns a null ptr). Definition at line 165 of file FastJets.hh. References FastJets::_cseq, and FastJets::areaDef(). Referenced by ATLAS_2010_S8914702::analyze(), and ATLAS_2012_I1093738::analyze().
Compare projections. Implements JetAlg. Definition at line 82 of file FastJets.cc. References FastJets::_adef, FastJets::_jdef, JetAlg::_useInvisibles, JetAlg::_useMuons, Rivet::cmp(), and Projection::mkNamedPCmp(). { const FastJets& other = dynamic_cast<const FastJets&>(p); return \ cmp(_useMuons, other._useMuons) || cmp(_useInvisibles, other._useInvisibles) || mkNamedPCmp(other, "FS") || cmp(_jdef.jet_algorithm(), other._jdef.jet_algorithm()) || cmp(_jdef.recombination_scheme(), other._jdef.recombination_scheme()) || cmp(_jdef.plugin(), other._jdef.plugin()) || cmp(_jdef.R(), other._jdef.R()) || cmp(_adef, other._adef); }
Determine if the jet collection is empty. Definition at line 208 of file JetAlg.hh. References JetAlg::size(). { return size() != 0; }
Template-usable interface common to FinalState. Definition at line 217 of file JetAlg.hh. References JetAlg::jets(). { return jets(); } 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(), ChargedFinalState::project(), InitialQuarks::project(), PromptFinalState::project(), LossyFinalState< ConstRandomFilter >::project(), UnstableFinalState::project(), and VetoedFinalState::project(). { string logname = "Rivet.Projection." + name(); return Log::getLog(logname); }
Get the named projection, specifying return type via a template argument. Definition at line 52 of file ProjectionApplier.hh. References ProjectionHandler::getProjection(), and ProjectionApplier::getProjHandler(). Referenced by ProjectionApplier::_applyProjection(), Rivet::pcmp(), and Hemispheres::project(). { const Projection& p = getProjHandler().getProjection(*this, name); return pcast<PROJ>(p); }
Get the named projection (non-templated, so returns as a reference to a Projection base class). Definition at line 60 of file ProjectionApplier.hh. References ProjectionHandler::getProjection(), and ProjectionApplier::getProjHandler(). { return getProjHandler().getProjection(*this, name); }
Get the contained projections, including recursion. Definition at line 45 of file ProjectionApplier.hh. References ProjectionHandler::DEEP, ProjectionHandler::getChildProjections(), and ProjectionApplier::getProjHandler(). Referenced by Projection::beamPairs(). { return getProjHandler().getChildProjections(*this, ProjectionHandler::DEEP); }
Get a reference to the ProjectionHandler for this thread. Definition at line 97 of file ProjectionApplier.hh. References ProjectionApplier::_projhandler. Referenced by ProjectionApplier::_addProjection(), ProjectionApplier::getProjection(), ProjectionApplier::getProjections(), and ProjectionApplier::~ProjectionApplier(). { return _projhandler; }
Return the jet definition. Definition at line 171 of file FastJets.hh. References FastJets::_jdef. { return _jdef; }
Get jets in no guaranteed order, with optional cuts on
Definition at line 86 of file JetAlg.hh. References JetAlg::_jets(), Rivet::Cuts::open(), and JetAlg::size(). Referenced by CMS_2012_I1087342::analyze(), CMS_2011_S9086218::analyze(), CMS_2012_I1102908::analyze(), ZEUS_2001_S4815815::analyze(), CDF_2008_S7540469::analyze(), CDF_2008_S7782535::analyze(), D0_2008_S7662670::analyze(), CDF_2008_S7541902::analyze(), ATLAS_2013_I1217863_Z::analyze(), ATLAS_2013_I1217863_W::analyze(), JetAlg::entities(), JetAlg::jets(), JetAlg::jetsByE(), JetAlg::jetsByEt(), JetAlg::jetsByP(), and JetAlg::jetsByPt(). { const Jets rawjets = _jets(0.0); // arg means no pT cut // Just return a copy of rawjets if the cut is open if (c == Cuts::open()) return rawjets; // If there is a non-trivial cut... Jets rtn; rtn.reserve(size()); foreach (const Jet& j, rawjets) if (c->accept(j)) rtn.push_back(j); return rtn; }
Get the jets, ordered by supplied sorting function object, with optional cuts on
Definition at line 101 of file JetAlg.hh. References JetAlg::jets(), and Rivet::sortBy(). Get the jets, ordered by supplied sorting function object, with optional cuts on
Definition at line 109 of file JetAlg.hh. References JetAlg::jets(), and Rivet::sortBy().
Get jets in no guaranteed order, with optional cuts on
Definition at line 168 of file JetAlg.hh. References Rivet::Cuts::etaIn(), JetAlg::jets(), Rivet::Cuts::pT, Rivet::RAPIDITY, and Rivet::Cuts::rapIn(). { if (rapscheme == PSEUDORAPIDITY) { return jets((Cuts::pT >= ptmin) & (Cuts::pT < ptmax) & (Cuts::rapIn(rapmin, rapmax))); } else if (rapscheme == RAPIDITY) { return jets((Cuts::pT >= ptmin) & (Cuts::pT < ptmax) & (Cuts::etaIn(rapmin, rapmax))); } throw LogicError("Unknown rapidity scheme. This shouldn't be possible!"); }
Get the jets, ordered by
Definition at line 144 of file JetAlg.hh. References Rivet::cmpMomByE(), and JetAlg::jets().
Get the jets, ordered by
Definition at line 152 of file JetAlg.hh. References Rivet::cmpMomByEt(), and JetAlg::jets(). { return jets(c, cmpMomByEt); }
Get the jets, ordered by
Definition at line 136 of file JetAlg.hh. References Rivet::cmpMomByP(), and JetAlg::jets().
Get the jets, ordered by
This is a very common use-case, so is available as syntatic sugar for jets(c, cmpMomByPt).
Definition at line 121 of file JetAlg.hh. References Rivet::cmpMomByPt(), and JetAlg::jets(). Referenced by CDF_2008_S8093652::analyze(), STAR_2006_S6870392::analyze(), ATLAS_2010_CONF_2010_049::analyze(), CMS_2013_I1208923::analyze(), CMS_2011_S9215166::analyze(), CMS_2013_I1272853::analyze(), MC_DIJET::analyze(), MC_LEADJETUE::analyze(), D0_2008_S7863608::analyze(), ATLAS_2015_I1364361::analyze(), CMS_2013_I1209721::analyze(), CMS_2015_I1310737::analyze(), D0_2004_S5992206::analyze(), CMS_2013_I1218372::analyze(), CMS_2013_I1261026::analyze(), ATLAS_2012_I1183818::analyze(), MC_TTBAR::analyze(), CDF_2001_S4751469::analyze(), ATLAS_2012_I1083318::analyze(), ATLAS_2014_I1312627::analyze(), MC_SUSY::analyze(), ATLAS_2014_I1319490::analyze(), CDF_2009_S8057893::CDF_2009_S8057893::analyze(), CMS_2013_I1261026::eventDecomp(), CMS_2013_I1258128::makePhotonCut(), CMS_2013_I1258128::makeZCut(), and ATLAS_2011_I945498::selectJets(). { return jets(c, cmpMomByPt); } Get the jets, ordered by
This is a very common use-case, so is available as syntatic sugar for jets(Cuts::pT >= ptmin, cmpMomByPt).
Definition at line 186 of file JetAlg.hh. References Rivet::cmpMomByPt(), JetAlg::jets(), and Rivet::Cuts::pT. { return jets(Cuts::pT >= ptmin, cmpMomByPt); }
Shortcut to make a named Cmp<Projection> comparison with the Definition at line 47 of file Projection.cc. References Rivet::pcmp(). Referenced by BeamThrust::compare(), ChargedLeptons::compare(), FParameter::compare(), CentralEtHCM::compare(), MergedFinalState::compare(), ChargedFinalState::compare(), DISLepton::compare(), DISKinematics::compare(), DISFinalState::compare(), VisibleFinalState::compare(), TauFinder::compare(), PromptFinalState::compare(), NeutralFinalState::compare(), Spherocity::compare(), LeadingParticlesFinalState::compare(), ParisiTensor::compare(), FoxWolframMoments::compare(), LossyFinalState< ConstRandomFilter >::compare(), Thrust::compare(), InvMassFinalState::compare(), Sphericity::compare(), HeavyHadrons::compare(), MissingMomentum::compare(), Hemispheres::compare(), DressedLeptons::compare(), ZFinder::compare(), WFinder::compare(), IdentifiedFinalState::compare(), VetoedFinalState::compare(), JetShape::compare(), and FastJets::compare(). { return pcmp(*this, otherparent, pname); }
Shortcut to make a named Cmp<Projection> comparison with the Definition at line 53 of file Projection.cc. References Rivet::pcmp(). { return pcmp(*this, otherparent, pname); }
Get the name of the projection. Implements ProjectionApplier. Definition at line 101 of file Projection.hh. References Projection::_name. Referenced by ProjectionApplier::_addProjection(), ProjectionHandler::_checkDuplicate(), ProjectionHandler::_clone(), ProjectionHandler::_getEquiv(), VetoedFinalState::addVetoOnThisFinalState(), Projection::getLog(), ProjectionHandler::registerProjection(), and Projection::setName(). { return _name; }
Number of jets above the Definition at line 169 of file FastJets.cc. References FastJets::_cseq. Referenced by FastJets::size(). Perform the projection on the Event.
Implements JetAlg. Definition at line 103 of file FastJets.cc. References JetAlg::_useInvisibles, JetAlg::_useMuons, FastJets::calc(), JetAlg::DECAY_INVISIBLES, JetAlg::DECAY_MUONS, Rivet::PID::isMuon(), JetAlg::NO_INVISIBLES, JetAlg::NO_MUONS, and Rivet::particles(). { // Assemble final state particles const string fskey = (_useInvisibles == JetAlg::NO_INVISIBLES) ? "VFS" : "FS"; Particles fsparticles = applyProjection<FinalState>(e, fskey).particles(); // Remove prompt invisibles if needed (already done by VFS if using NO_INVISIBLES) if (_useInvisibles == JetAlg::DECAY_INVISIBLES) fsparticles.erase( std::remove_if(fsparticles.begin(), fsparticles.end(), isPromptInvisible), fsparticles.end() ); // Remove prompt/all muons if needed if (_useMuons == JetAlg::DECAY_MUONS) fsparticles.erase( std::remove_if(fsparticles.begin(), fsparticles.end(), isPromptMuon), fsparticles.end() ); else if (_useMuons == JetAlg::NO_MUONS) fsparticles.erase( std::remove_if(fsparticles.begin(), fsparticles.end(), isMuon), fsparticles.end() ); // Tagging particles /// @todo Allow the user to specify tag particle kinematic thresholds const Particles chadrons = applyProjection<HeavyHadrons>(e, "HFHadrons").cHadrons(); const Particles bhadrons = applyProjection<HeavyHadrons>(e, "HFHadrons").bHadrons(); const Particles taus = applyProjection<FinalState>(e, "Taus").particles(); calc(fsparticles, chadrons+bhadrons+taus); }
Get the pseudo jets (unordered). Definition at line 199 of file FastJets.cc. References FastJets::clusterSeq(). Referenced by ATLAS_2010_S8914702::analyze(), ATLAS_2012_I1093738::analyze(), FastJets::pseudojets(), FastJets::pseudoJetsByE(), FastJets::pseudoJetsByPt(), and FastJets::pseudoJetsByRapidity(). { return (clusterSeq() != NULL) ? clusterSeq()->inclusive_jets(ptmin) : PseudoJets(); }
Alias. Definition at line 130 of file FastJets.hh. References FastJets::pseudoJets(). Referenced by FastJets::_jets(). { return pseudoJets(ptmin); }
Get the pseudo jets, ordered by Definition at line 140 of file FastJets.hh. References FastJets::pseudoJets(). Referenced by FastJets::pseudojetsByE(). { return sorted_by_E(pseudoJets(ptmin)); }
Alias. Definition at line 144 of file FastJets.hh. References FastJets::pseudoJetsByE(). { return pseudoJetsByE(ptmin); }
Get the pseudo jets, ordered by Definition at line 133 of file FastJets.hh. References FastJets::pseudoJets(). Referenced by CDF_2008_S8095620::analyze(), CDF_2006_S6653332::analyze(), ATLAS_2012_I1094564::analyze(), and FastJets::pseudojetsByPt(). { return sorted_by_pt(pseudoJets(ptmin)); }
Alias. Definition at line 137 of file FastJets.hh. References FastJets::pseudoJetsByPt(). { return pseudoJetsByPt(ptmin); }
Get the pseudo jets, ordered by rapidity. Definition at line 147 of file FastJets.hh. References FastJets::pseudoJets(). Referenced by FastJets::pseudojetsByRapidity(). { return sorted_by_rapidity(pseudoJets(ptmin)); }
Alias. Definition at line 151 of file FastJets.hh. References FastJets::pseudoJetsByRapidity(). { return pseudoJetsByRapidity(ptmin); }
Reset the projection. Jet def, etc. are unchanged.
Implements JetAlg. Definition at line 162 of file FastJets.cc. References FastJets::_particles, and FastJets::_yscales. { _yscales.clear(); _particles.clear(); /// @todo _cseq = fastjet::ClusterSequence(); }
Used by derived classes to set their name. Definition at line 120 of file Projection.hh. References Projection::_name, and Projection::name(). Referenced by FastJets::_initBase(), Beam::Beam(), BeamThrust::BeamThrust(), CentralEtHCM::CentralEtHCM(), ChargedFinalState::ChargedFinalState(), ChargedLeptons::ChargedLeptons(), ConstLossyFinalState::ConstLossyFinalState(), DISFinalState::DISFinalState(), DISKinematics::DISKinematics(), DISLepton::DISLepton(), DressedLeptons::DressedLeptons(), FinalState::FinalState(), FoxWolframMoments::FoxWolframMoments(), FParameter::FParameter(), HadronicFinalState::HadronicFinalState(), HeavyHadrons::HeavyHadrons(), Hemispheres::Hemispheres(), IdentifiedFinalState::IdentifiedFinalState(), InitialQuarks::InitialQuarks(), JetAlg::JetAlg(), JetShape::JetShape(), LeadingParticlesFinalState::LeadingParticlesFinalState(), LossyFinalState< ConstRandomFilter >::LossyFinalState(), MergedFinalState::MergedFinalState(), MissingMomentum::MissingMomentum(), NeutralFinalState::NeutralFinalState(), NonHadronicFinalState::NonHadronicFinalState(), ParisiTensor::ParisiTensor(), PrimaryHadrons::PrimaryHadrons(), PromptFinalState::PromptFinalState(), Sphericity::Sphericity(), Spherocity::Spherocity(), TauFinder::TauFinder(), Thrust::Thrust(), TriggerCDFRun0Run1::TriggerCDFRun0Run1(), TriggerCDFRun2::TriggerCDFRun2(), TriggerUA5::TriggerUA5(), UnstableFinalState::UnstableFinalState(), VetoedFinalState::VetoedFinalState(), VisibleFinalState::VisibleFinalState(), WFinder::WFinder(), and ZFinder::ZFinder().
Number of jets. Implements JetAlg. Definition at line 120 of file FastJets.hh. References FastJets::numJets(). Referenced by CDF_2008_S8095620::analyze(), and CDF_2006_S6653332::analyze(). { return numJets(); }
Include (some) invisible particles in jet construction. The default behaviour is that jets are only constructed from visible particles. Some jet studies, including those from ATLAS, use a definition in which neutrinos from hadron decays are included via MC-based calibrations. Setting this flag to true avoids the automatic restriction to a VisibleFinalState. Definition at line 61 of file JetAlg.hh. References JetAlg::_useInvisibles. Referenced by CMS_2011_S8973270::init(), ATLAS_2015_I1387176::init(), ATLAS_2011_I930220::init(), D0_2011_I895662::init(), ATLAS_2011_S9128077::init(), ATLAS_2014_I1319490::init(), CDF_2008_S7782535::init(), ATLAS_2014_I1325553::init(), ATLAS_2013_I1219109::init(), ATLAS_2014_I1268975::init(), D0_2000_I499943::init(), ATLAS_2011_S8924791::init(), CDF_2005_S6217184::init(), ATLAS_2011_I945498::init(), ALEPH_2004_S5765862::init(), JADE_OPAL_2000_S4300807::init(), ATLAS_2013_I1217863_W::init(), ATLAS_2013_I1217863_Z::init(), ATLAS_2013_I1244522::init(), ATLAS_2015_I1345452::init(), CMS_2014_I1303894::init(), MC_HFJETS::init(), ATLAS_2013_I1230812::init(), ATLAS_2014_I1307243::init(), ATLAS_2015_CONF_2015_041::init(), ATLAS_2012_I1082936::init(), ATLAS_2014_I1326641::init(), ATLAS_2014_I1306294::init(), ATLAS_2014_I1304688::init(), ATLAS_2012_I1083318::init(), ATLAS_2013_I1217867::init(), ATLAS_2014_I1312627::init(), and ATLAS_2012_I1093738::init(). { _useInvisibles = useinvis; }
Include (some) invisible particles in jet construction. The default behaviour is that jets are only constructed from visible particles. Some jet studies, including those from ATLAS, use a definition in which neutrinos from hadron decays are included via MC-based calibrations. Setting this flag to true avoids the automatic restriction to a VisibleFinalState.
Definition at line 73 of file JetAlg.hh. References JetAlg::_useInvisibles, JetAlg::DECAY_INVISIBLES, and JetAlg::NO_INVISIBLES. { _useInvisibles = useinvis ? DECAY_INVISIBLES : NO_INVISIBLES; }
Use provided jet area definition.
Provide an adef null pointer to re-disable jet area calculation Definition at line 109 of file FastJets.hh. References FastJets::_adef. Referenced by ATLAS_2014_I1307756::init(), ATLAS_2013_I1244522::init(), ATLAS_2013_I1263495::init(), ATLAS_2011_I921594::init(), ATLAS_2012_I1199269::init(), ATLAS_2010_S8914702::init(), ATLAS_2011_S9120807::init(), and ATLAS_2012_I1093738::init(). { _adef = adef; }
Include (some) muons in jet construction. The default behaviour is that jets are only constructed from visible particles. Some jet studies, including those from ATLAS, use a definition in which neutrinos from hadron decays are included via MC-based calibrations. Setting this flag to true avoids the automatic restriction to a VisibleFinalState. Definition at line 51 of file JetAlg.hh. References JetAlg::_useMuons. { _useMuons = usemuons; } Friends And Related Function Documentation
The Cmp specialization for Projection is a friend. Definition at line 35 of file Projection.hh.
Event is a friend. Definition at line 32 of file Projection.hh. Member Data Documentation
Pointer to user-handled area definition. Definition at line 232 of file FastJets.hh. Referenced by FastJets::areaDef(), FastJets::calc(), FastJets::compare(), and FastJets::useJetArea().
Flag to forbid projection registration in analyses until the init phase. Definition at line 143 of file ProjectionApplier.hh. Referenced by ProjectionApplier::_addProjection(), and Analysis::Analysis().
Cluster sequence. Definition at line 235 of file FastJets.hh. Referenced by FastJets::calc(), FastJets::clusterSeq(), FastJets::clusterSeqArea(), and FastJets::numJets().
Jet definition. Definition at line 229 of file FastJets.hh. Referenced by FastJets::_init1(), FastJets::_init2(), FastJets::_init3(), FastJets::_init4(), FastJets::calc(), FastJets::compare(), and FastJets::jetDef().
set of particles sorted by their PT2 Definition at line 245 of file FastJets.hh. Referenced by FastJets::_jets(), FastJets::calc(), and FastJets::reset().
FastJet external plugin. Definition at line 238 of file FastJets.hh. Referenced by FastJets::_init1(), and FastJets::_init4().
Flag to determine whether or not to exclude (some) invisible particles from the would-be constituents. Definition at line 238 of file JetAlg.hh. Referenced by FastJets::compare(), FastJets::project(), and JetAlg::useInvisibles().
Flag to determine whether or not to exclude (some) muons from the would-be constituents. Definition at line 235 of file JetAlg.hh. Referenced by FastJets::compare(), FastJets::project(), and JetAlg::useMuons().
Map of vectors of y scales. This is mutable so we can use caching/lazy evaluation. Definition at line 241 of file FastJets.hh. Referenced by FastJets::reset(). The documentation for this class was generated from the following files: Generated on Wed Oct 7 2015 12:09:29 for The Rivet MC analysis system by ![]() |