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 51 of file FastJets.hh. Member Typedef Documentation
Member Enumeration Documentation
Wrapper enum for selected Fastjet jet algorithms.
Definition at line 55 of file FastJets.hh. 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 69 of file FastJets.hh. References FastJets::_init1(). Referenced by FastJets::clone().
Native argument constructor, using FastJet alg/scheme enums. Definition at line 74 of file FastJets.hh. References FastJets::_init2().
Explicitly pass in an externally-constructed plugin (must be heap-allocated, Rivet will delete) Definition at line 79 of file FastJets.hh. References FastJets::_init3().
Explicitly pass in an externally-constructed plugin (must be heap-allocated, Rivet will delete) Definition at line 82 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::_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 90 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 93 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 96 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; }
Shared utility functions to implement constructor behaviour. Definition at line 9 of file FastJets.cc. References 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::PXCONE, Projection::setName(), FastJets::SISCONE, and FastJets::TRACKJET. Referenced by FastJets::FastJets(). { setName("FastJets"); 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 55 of file FastJets.cc. References FastJets::_jdef, and Projection::setName(). Referenced by FastJets::FastJets().
Definition at line 61 of file FastJets.cc. References FastJets::_jdef, FastJets::_plugin, and Projection::setName(). Referenced by FastJets::FastJets(). Get the jets (unordered). Implements JetAlg. Definition at line 156 of file FastJets.cc. References FastJets::_pseudojetsToJets(), and FastJets::pseudoJets(). { Jets rtn = _pseudojetsToJets(pseudoJets(ptmin)); return rtn; }
Definition at line 119 of file FastJets.cc. References FastJets::_particles, and FastJets::clusterSeq(). Referenced by FastJets::_jets(). { Jets rtn; foreach (const fastjet::PseudoJet& pj, pjets) { assert(clusterSeq()); const PseudoJets parts = clusterSeq()->constituents(pj); vector<Particle> constituents; 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()); constituents.push_back(found->second); } FourMomentum pjet(pj.E(), pj.px(), pj.py(), pj.pz()); Jet j(constituents, pjet); rtn.push_back(j); } /// @todo Cache? return rtn; }
Add a colliding beam pair. Definition at line 109 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 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_S8656010::init(), CMS_2010_S8547297::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_2011_S8941262::init(), CMS_2011_S8968497::init(), CMS_2011_S8973270::init(), CMS_2011_S8978280::init(), ATLAS_2010_CONF_2010_049::init(), ATLAS_2011_I894867::init(), LHCB_2013_I1208105::init(), ALICE_2011_S8909580::init(), CMS_2011_S9086218::init(), CMS_2011_S9215166::init(), TOTEM_2012_I1115294::init(), CMS_2012_I1087342::init(), CMS_2012_I1184941::init(), CMS_2012_I1193338::init(), ALICE_2011_S8945144::init(), CDF_2007_S7057202::init(), TOTEM_2012_002::init(), CMS_2011_I954992::init(), CMS_2011_S8957746::init(), ATLAS_2010_S8591806::init(), ATLAS_2011_S9002537::init(), D0_2011_I895662::init(), CMS_2011_S9088458::init(), LHCF_2012_I1115479::init(), ATLAS_2011_I930220::init(), D0_2010_S8570965::init(), MC_JETS::init(), MC_KTSPLITTINGS::init(), CDF_2012_NOTE10874::init(), STAR_2006_S6500200::init(), STAR_2008_S7993412::init(), UA5_1987_S1640666::init(), CMS_2013_I1209721::init(), CMS_QCD_10_024::init(), CDF_1997_S3541940::init(), SFM_1984_S1178091::init(), ARGUS_1993_S2789213::init(), ATLAS_2012_I1204784::init(), CMS_2013_I1258128::init(), BABAR_2007_S7266081::init(), ATLAS_2011_S9128077::init(), CDF_1993_S2742446::init(), CDF_2000_S4155203::init(), CDF_2005_S6080774::init(), MC_DIJET::init(), MC_DIPHOTON::init(), CDF_2006_S6450792::init(), MC_HINC::init(), MC_WINC::init(), MC_ZINC::init(), MC_ZZINC::init(), MC_ZKTSPLITTINGS::init(), STAR_2006_S6870392::init(), UA5_1982_S875503::init(), UA5_1989_S1926373::init(), D0_2008_S6879055::init(), D0_2010_S8821313::init(), MC_HJETS::init(), MC_HKTSPLITTINGS::init(), MC_PHOTONINC::init(), CDF_2008_S8093652::init(), MC_WKTSPLITTINGS::init(), MC_WWINC::init(), MC_ZJETS::init(), CDF_2008_S7782535::init(), ARGUS_1993_S2669951::init(), MC_ZZJETS::init(), MC_ZZKTSPLITTINGS::init(), ATLAS_2012_I1091481::init(), ALEPH_1991_S2435284::init(), CMS_2012_I1102908::init(), MC_WWKTSPLITTINGS::init(), ATLAS_2012_I1183818::init(), ATLAS_2013_I1243871::init(), D0_2008_S7863608::init(), ATLAS_2011_S9035664::init(), D0_2010_S8671338::init(), CMS_2012_I1107658::init(), E735_1998_S3905616::init(), MC_PHOTONJETS::init(), MC_PHOTONKTSPLITTINGS::init(), MC_WJETS::init(), MC_WWJETS::init(), ATLAS_2011_I944826::init(), MC_TTBAR::init(), MC_PHOTONJETUE::init(), OPAL_1995_S3198391::init(), OPAL_1996_S3257789::init(), OPAL_1997_S3608263::init(), OPAL_1998_S3702294::init(), OPAL_1998_S3749908::init(), OPAL_2000_S4418603::init(), ATLAS_2011_S8924791::init(), ATLAS_2010_S8817804::init(), D0_1996_S3324664::init(), D0_2007_S7075677::init(), ALEPH_1996_S3196992::init(), D0_2009_S8202443::init(), DELPHI_1999_S3960137::init(), CDF_1990_S2089246::init(), UA5_1986_S1583476::init(), CDF_1988_S1865951::init(), EXAMPLE_CUTS::init(), ALEPH_2002_S4823664::init(), MC_GENERIC::init(), CDF_2005_S6217184::init(), MC_IDENTIFIED::init(), MC_LEADJETUE::init(), ATLAS_2011_I945498::init(), ATLAS_2011_I954993::init(), ATLAS_2012_I1119557::init(), D0_2001_S4674421::init(), D0_2008_S7554427::init(), CDF_1994_S2952106::init(), H1_1995_S3167097::init(), MC_QCD_PARTONS::init(), UA1_1990_S2044935::init(), ZEUS_2001_S4815815::init(), ATLAS_2012_I1188891::init(), D0_2008_S7837160::init(), ATLAS_2011_S9131140::init(), LHCB_2011_I919315::init(), CDF_2008_S7540469::init(), CDF_2008_S7828950::init(), CMS_2012_I941555::init(), CDF_2008_S8095620::init(), STAR_2006_S6860818::init(), JADE_1998_S3612880::init(), D0_1996_S3214044::init(), JADE_OPAL_2000_S4300807::init(), DELPHI_1995_S3137023::init(), STAR_2009_UE_HELEN::init(), ALICE_2010_S8624100::init(), OPAL_1997_S3396100::init(), EXAMPLE::init(), MC_PHOTONS::init(), MC_SUSY::init(), ATLAS_2011_I925932::init(), ALEPH_2004_S5765862::init(), ATLAS_2012_I1125575::init(), CMS_2013_I1218372::init(), ALICE_2010_S8625980::init(), D0_2000_S4480767::init(), ALICE_2010_S8706239::init(), LHCB_2013_I1218996::init(), CDF_2009_NOTE_9936::init(), ATLAS_2011_S8971293::init(), D0_2009_S8349509::init(), CDF_1996_S3108457::init(), CDF_2001_S4563131::init(), CDF_1998_S3618439::init(), CDF_2009_S8383952::init(), CDF_2000_S4266730::init(), ATLAS_2011_I926145::init(), D0_2009_S8320160::init(), CDF_2009_S8436959::init(), SLD_2002_S4869273::init(), CDF_1996_S3349578::init(), MC_WPOL::init(), D0_2006_S6438750::init(), ALEPH_2001_S4656318::init(), DELPHI_2002_069_CONF_603::init(), ATLAS_2010_S8919674::init(), CDF_2001_S4517016::init(), ATLAS_2012_I1082936::init(), D0_2010_S8566488::init(), SLD_1996_S3398250::init(), CDF_1996_S3418421::init(), UA5_1988_S1867512::init(), D0_2008_S7662670::init(), ATLAS_2011_S9108483::init(), CDF_2006_S6653332::init(), CDF_2008_S7541902::init(), ATLAS_2012_I1125961::init(), ATLAS_2013_I1230812::init(), ATLAS_2012_I1095236::init(), ATLAS_2012_CONF_2012_109::init(), ATLAS_2012_I1083318::init(), ATLAS_2011_S9212183::init(), ATLAS_2012_I1180197::init(), DELPHI_2000_S4328825::init(), ATLAS_2012_I943401::init(), ATLAS_2012_I1190891::init(), ATLAS_2011_S9225137::init(), ATLAS_2013_I1217867::init(), ATLAS_2012_I1112263::init(), ATLAS_2012_CONF_2012_103::init(), ATLAS_2011_S9212353::init(), ATLAS_2012_CONF_2012_105::init(), TASSO_1990_S2148048::init(), ATLAS_2011_S9019561::init(), ATLAS_2012_CONF_2012_001::init(), ATLAS_2012_I1082009::init(), ATLAS_2012_I1117704::init(), ATLAS_2011_CONF_2011_090::init(), OPAL_2002_S5361494::init(), ATLAS_2012_CONF_2012_104::init(), ATLAS_2012_I946427::init(), ATLAS_2011_S8983313::init(), ATLAS_2012_I1186556::init(), ATLAS_2011_CONF_2011_098::init(), ATLAS_2012_I1084540::init(), ATLAS_2012_I1126136::init(), CMS_2013_I1224539_ZJET::init(), CMS_2013_I1224539_WJET::init(), CMS_2013_I1224539_DIJET::init(), CDF_2009_S8233977::init(), CDF_2010_S8591881_DY::init(), ATLAS_2012_CONF_2012_153::init(), CDF_2010_S8591881_QCD::init(), D0_2004_S5992206::init(), ATLAS_2010_S8914702::init(), LHCB_2011_I917009::init(), CDF_2001_S4751469::init(), ATLAS_2011_S9120807::init(), D0_2008_S7719523::init(), ATLAS_2010_S8918562::init(), LHCB_2012_I1119400::init(), ATLAS_2012_I1094568::init(), ATLAS_2011_S9041966::init(), DELPHI_1996_S3430090::init(), ALEPH_1996_S3486095::init(), ALEPH_1999_S4193598::init(), ATLAS_2012_I1093738::init(), OPAL_2004_S6132243::init(), OPAL_1994_S2927284::init(), STAR_2008_S7869363::init(), ATLAS_2011_S9126244::init(), BABAR_2007_S6895344::init(), MC_VH2BB::init(), BABAR_2005_S6181155::init(), BELLE_2001_S4598261::init(), ATLAS_2011_I919017::init(), OPAL_2001_S4553896::init(), DELPHI_2003_WUD_03_11::init(), CDF_2004_S5839831::init(), BABAR_2003_I593379::init(), ARGUS_1993_S2653028::init(), CLEO_2004_S5809304::init(), OPAL_1998_S3780481::init(), OPAL_1993_S2692198::init(), ATLAS_2012_I1093734::init(), ATLAS_2012_I1094564::init(), H1_1994_S2919893::init(), H1_2000_S4129130::init(), SLD_2004_S5693039::init(), BELLE_2006_S6265367::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(), Sphericity::Sphericity(), Spherocity::Spherocity(), 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 (FastJet-specific). May be null. Definition at line 167 of file FastJets.hh. References FastJets::_adef. { 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 33 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 22 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). Implements JetAlg. Definition at line 95 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(), and FastJets::project(). { _particles.clear(); vector<fastjet::PseudoJet> vecs; // Store 4 vector data about each particle into vecs int counter = 1; foreach (const Particle& p, ps) { const FourMomentum fv = p.momentum(); fastjet::PseudoJet pJet(fv.px(), fv.py(), fv.pz(), fv.E()); pJet.set_user_index(counter); vecs.push_back(pJet); _particles[counter] = p; ++counter; } MSG_DEBUG("Running FastJet ClusterSequence construction"); // Choose CSeq as basic or area-calculating depending on whether _adef pointer is non-null. if (_adef == 0) { _cseq.reset(new fastjet::ClusterSequence(vecs, _jdef)); } else { _cseq.reset(new fastjet::ClusterSequenceArea(vecs, _jdef, *_adef)); } }
Clone on the heap. Implements JetAlg. Definition at line 101 of file FastJets.hh. References FastJets::FastJets(). { return new FastJets(*this); }
Return the cluster sequence (FastJet-specific). Definition at line 150 of file FastJets.hh. References FastJets::_cseq. Referenced by FastJets::_pseudojetsToJets(), MC_JetSplittings::analyze(), ALEPH_1996_S3196992::analyze(), JADE_OPAL_2000_S4300807::analyze(), JADE_1998_S3612880::analyze(), OPAL_2001_S4553896::analyze(), OPAL_2004_S6132243::analyze(), DELPHI_2003_WUD_03_11::analyze(), ALEPH_1996_S3486095::analyze(), ALEPH_2004_S5765862::analyze(), DELPHI_1996_S3430090::analyze(), FastJets::filterJet(), ATLAS_2012_I1094564::splitjet(), FastJets::splitJet(), and FastJets::ySubJet(). { return _cseq.get(); }
Return the cluster sequence (FastJet-specific).
Definition at line 155 of file FastJets.hh. References FastJets::_adef, and FastJets::_cseq.
Compare projections. Implements JetAlg. Definition at line 70 of file FastJets.cc. References FastJets::_adef, FastJets::_jdef, JetAlg::_useInvisibles, Rivet::cmp(), and Projection::mkNamedPCmp(). { const FastJets& other = dynamic_cast<const FastJets&>(p); return \ cmp(_useInvisibles, other._useInvisibles) || (_useInvisibles ? mkNamedPCmp(other, "FS") : mkNamedPCmp(other, "VFS")) || 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 117 of file JetAlg.hh. References JetAlg::size(). { return size() != 0; }
Template-usable interface common to FinalState. Definition at line 126 of file JetAlg.hh. References JetAlg::jets(). { return jets(); }
Filter a jet a la PRL100,242001(2008). Based on code from G.Salam, A.Davison. Definition at line 241 of file FastJets.cc. References FastJets::clusterSeq(). { assert(clusterSeq()); if (jet.E() <= 0.0 || clusterSeq()->constituents(jet).size() == 0) { return jet; } if (stingy_R == 0.0) { stingy_R = def_R; } stingy_R = def_R < stingy_R ? def_R : stingy_R; fastjet::JetDefinition stingy_jet_def(fastjet::cambridge_algorithm, stingy_R); //FlavourRecombiner recom; //stingy_jet_def.set_recombiner(&recom); fastjet::ClusterSequence scs(clusterSeq()->constituents(jet), stingy_jet_def); std::vector<fastjet::PseudoJet> stingy_jets = sorted_by_pt(scs.inclusive_jets()); fastjet::PseudoJet reconst_jet(0.0, 0.0, 0.0, 0.0); for (unsigned isj = 0; isj < std::min(3U, (unsigned int) stingy_jets.size()); ++isj) { reconst_jet += stingy_jets[isj]; } return reconst_jet; } Get a Log object based on the getName() property of the calling projection object. Reimplemented from ProjectionApplier. Definition at line 116 of file Projection.hh. References Projection::name(). Referenced by Projection::beamPairs(), InvMassFinalState::calc(), ChargedFinalState::project(), InitialQuarks::project(), UnstableFinalState::project(), LossyFinalState< ConstRandomFilter >::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 (FastJet-specific). Definition at line 162 of file FastJets.hh. References FastJets::_jdef. { return _jdef; }
Get jets in no guaranteed order, with optional cuts on and rapidity.
Definition at line 41 of file JetAlg.hh. References JetAlg::_jets(), FourVector::eta(), Rivet::GeV, Rivet::inRange(), Jet::momentum(), MSG_DEBUG, Rivet::PSEUDORAPIDITY, FourMomentum::pT(), Rivet::RAPIDITY, and FourMomentum::rapidity(). Referenced by CMS_2012_I1087342::analyze(), CMS_2011_S9086218::analyze(), CMS_2012_I1102908::analyze(), ZEUS_2001_S4815815::analyze(), CDF_2008_S7540469::analyze(), D0_2008_S7662670::analyze(), JetAlg::entities(), JetAlg::jets(), JetAlg::jetsByE(), JetAlg::jetsByEt(), JetAlg::jetsByP(), and JetAlg::jetsByPt(). { const Jets rawjets = _jets(ptmin); Jets rtn; MSG_DEBUG("Raw jet size (with pTmin cut = " << ptmin/GeV << " GeV) = " << rawjets.size()); foreach (const Jet& j, rawjets) { const FourMomentum pj = j.momentum(); if (!inRange(pj.pT(), ptmin, ptmax)) continue; if (rapscheme == PSEUDORAPIDITY && !inRange(pj.eta(), rapmin, rapmax)) continue; if (rapscheme == RAPIDITY && !inRange(pj.rapidity(), rapmin, rapmax)) continue; rtn += j; } return rtn; }
Get the jets, ordered by supplied sorting function object, with optional cuts on and rapidity.
Definition at line 60 of file JetAlg.hh. References JetAlg::jets().
Get the jets, ordered by , with optional cuts on and rapidity.
Definition at line 88 of file JetAlg.hh. References Rivet::cmpMomByE(), and JetAlg::jets().
Get the jets, ordered by , with optional cuts on and rapidity.
Definition at line 96 of file JetAlg.hh. References Rivet::cmpMomByEt(), and JetAlg::jets(). Referenced by CDF_2008_S7541902::analyze(). { return jets(cmpMomByEt, ptmin, ptmax, rapmin, rapmax, rapscheme); }
Get the jets, ordered by , with optional cuts on and rapidity.
Definition at line 80 of file JetAlg.hh. References Rivet::cmpMomByP(), and JetAlg::jets().
Get the jets, ordered by , with optional cuts on and rapidity.
Definition at line 72 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_2011_S9215166::analyze(), MC_LEADJETUE::analyze(), MC_DIJET::analyze(), D0_2008_S7863608::analyze(), CMS_2013_I1209721::analyze(), D0_2004_S5992206::analyze(), CMS_2013_I1218372::analyze(), ATLAS_2012_I1183818::analyze(), MC_TTBAR::analyze(), CDF_2001_S4751469::analyze(), ATLAS_2012_I1083318::analyze(), MC_SUSY::analyze(), CDF_2009_S8057893::CDF_2009_S8057893::analyze(), CMS_2013_I1258128::makePhotonCut(), CMS_2013_I1258128::makeZCut(), and ATLAS_2011_I945498::selectJets(). { return jets(cmpMomByPt, ptmin, ptmax, rapmin, rapmax, rapscheme); }
Shortcut to make a named Cmp<Projection> comparison with the Definition at line 45 of file Projection.cc. References Rivet::pcmp(). Referenced by BeamThrust::compare(), ChargedLeptons::compare(), FParameter::compare(), CentralEtHCM::compare(), MergedFinalState::compare(), DISLepton::compare(), DISKinematics::compare(), ChargedFinalState::compare(), DISFinalState::compare(), NeutralFinalState::compare(), VisibleFinalState::compare(), LeadingParticlesFinalState::compare(), ParisiTensor::compare(), Spherocity::compare(), DressedLeptons::compare(), FoxWolframMoments::compare(), LossyFinalState< ConstRandomFilter >::compare(), MissingMomentum::compare(), Thrust::compare(), InvMassFinalState::compare(), Sphericity::compare(), Hemispheres::compare(), ZFinder::compare(), IdentifiedFinalState::compare(), WFinder::compare(), JetShape::compare(), VetoedFinalState::compare(), and FastJets::compare(). { return pcmp(*this, otherparent, pname); }
Shortcut to make a named Cmp<Projection> comparison with the Definition at line 51 of file Projection.cc. References Rivet::pcmp(). { return pcmp(*this, otherparent, pname); }
Get the name of the projection. Implements ProjectionApplier. Definition at line 103 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 cut. Definition at line 147 of file FastJets.cc. References FastJets::_cseq. Referenced by FastJets::size(). Perform the projection on the Event. Implements JetAlg. Definition at line 84 of file FastJets.cc. References JetAlg::_useInvisibles, FastJets::calc(), and Rivet::particles().
Get the pseudo jets (unordered). Definition at line 177 of file FastJets.cc. References FastJets::_cseq. Referenced by FastJets::_jets(), FastJets::pseudoJetsByE(), FastJets::pseudoJetsByPt(), and FastJets::pseudoJetsByRapidity(). { if (_cseq.get() != 0) { return _cseq->inclusive_jets(ptmin); } else { return PseudoJets(); } }
Get the pseudo jets, ordered by . Definition at line 140 of file FastJets.hh. References FastJets::pseudoJets(). { return sorted_by_E(pseudoJets(ptmin)); }
Get the pseudo jets, ordered by . Definition at line 135 of file FastJets.hh. References FastJets::pseudoJets(). Referenced by CDF_2008_S8095620::analyze(), CDF_2006_S6653332::analyze(), and ATLAS_2012_I1094564::analyze(). { return sorted_by_pt(pseudoJets(ptmin)); }
Get the pseudo jets, ordered by rapidity. Definition at line 145 of file FastJets.hh. References FastJets::pseudoJets(). { return sorted_by_rapidity(pseudoJets(ptmin)); }
Reset the projection. Jet def, etc. are unchanged.
Implements JetAlg. Definition at line 140 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 122 of file Projection.hh. References Projection::_name, and Projection::name(). Referenced by FastJets::_init1(), FastJets::_init2(), FastJets::_init3(), 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(), Sphericity::Sphericity(), Spherocity::Spherocity(), 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 124 of file FastJets.hh. References FastJets::numJets(). Referenced by CDF_2008_S8095620::analyze(), and CDF_2006_S6653332::analyze(). { return numJets(); }
Split a jet a la PRL100,242001(2008). Based on code from G.Salam, A.Davison. Definition at line 202 of file FastJets.cc. References FastJets::_cseq, FastJets::_jdef, FastJets::clusterSeq(), and MSG_DEBUG. { // Sanity cuts if (jet.E() <= 0 || _cseq->constituents(jet).size() <= 1) { return jet; } // Build a new cluster sequence just using the consituents of this jet. assert(clusterSeq()); fastjet::ClusterSequence cs(clusterSeq()->constituents(jet), _jdef); // Get the jet back again fastjet::PseudoJet remadeJet = cs.inclusive_jets()[0]; MSG_DEBUG("Jet2:" << remadeJet.m() << "," << remadeJet.e()); fastjet::PseudoJet parent1, parent2; fastjet::PseudoJet split(0.0, 0.0, 0.0, 0.0); while (cs.has_parents(remadeJet, parent1, parent2)) { MSG_DEBUG("Parents:" << parent1.m() << "," << parent2.m()); if (parent1.m2() < parent2.m2()) { fastjet::PseudoJet tmp; tmp = parent1; parent1 = parent2; parent2 = tmp; } double ktdist = parent1.kt_distance(parent2); double rtycut2 = 0.3*0.3; if (parent1.m() < ((2.0*remadeJet.m())/3.0) && ktdist > rtycut2*remadeJet.m2()) { break; } else { remadeJet = parent1; } } last_R = 0.5 * sqrt(parent1.squared_distance(parent2)); split.reset(remadeJet.px(), remadeJet.py(), remadeJet.pz(), remadeJet.E()); return split; }
Include invisible particles in jet construction. The default behaviour is that jets are only constructed from visible (i.e. charged under an SM gauge group) particles. Some jet studies, including those from ATLAS, use a definition in which neutrinos from hadron decays are included (via MC correction) in the experimental jet definition. Setting this flag to true avoids the automatic restriction to a VisibleFinalState. Definition at line 35 of file JetAlg.hh. References JetAlg::_useInvisibles. Referenced by CMS_2011_S8973270::init(), ATLAS_2011_I930220::init(), D0_2011_I895662::init(), ATLAS_2011_S9128077::init(), CDF_2008_S7782535::init(), ATLAS_2011_S8924791::init(), CDF_2005_S6217184::init(), ATLAS_2011_I945498::init(), ALEPH_2004_S5765862::init(), JADE_OPAL_2000_S4300807::init(), ATLAS_2012_I1082936::init(), ATLAS_2012_I1083318::init(), ATLAS_2013_I1230812::init(), ATLAS_2013_I1217867::init(), and ATLAS_2012_I1093738::init(). { _useInvisibles = useinvis; }
Use provided jet area definition.
Definition at line 116 of file FastJets.hh. References FastJets::_adef. Referenced by ATLAS_2010_S8914702::init(), ATLAS_2011_S9120807::init(), and ATLAS_2012_I1093738::init(). { _adef = adef; }
Get the subjet splitting variables for the given jet. Definition at line 186 of file FastJets.cc. References FastJets::_jdef, FastJets::_yscales, and FastJets::clusterSeq(). { assert(clusterSeq()); fastjet::ClusterSequence subjet_cseq(clusterSeq()->constituents(jet), _jdef); vector<double> yMergeVals; for (int i = 1; i < 4; ++i) { // Multiply the dmerge value by R^2 so that it corresponds to a // relative k_T (fastjet has 1/R^2 in the d_ij distance by default) const double ktmerge = subjet_cseq.exclusive_dmerge(i) * _jdef.R()*_jdef.R(); yMergeVals.push_back(ktmerge/jet.perp2()); } _yscales.insert(make_pair( jet.cluster_hist_index(), yMergeVals )); return yMergeVals; } Friends And Related Function Documentation
The Cmp specialization for Projection is a friend. Definition at line 37 of file Projection.hh.
Event is a friend. Definition at line 34 of file Projection.hh. Member Data Documentation
Pointer to user-handled area definition. Definition at line 212 of file FastJets.hh. Referenced by FastJets::areaDef(), FastJets::calc(), FastJets::clusterSeqArea(), 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 215 of file FastJets.hh. Referenced by FastJets::calc(), FastJets::clusterSeq(), FastJets::clusterSeqArea(), FastJets::numJets(), FastJets::pseudoJets(), and FastJets::splitJet().
Jet definition. Definition at line 209 of file FastJets.hh. Referenced by FastJets::_init1(), FastJets::_init2(), FastJets::_init3(), FastJets::calc(), FastJets::compare(), FastJets::jetDef(), FastJets::splitJet(), and FastJets::ySubJet().
set of particles sorted by their PT2 Definition at line 225 of file FastJets.hh. Referenced by FastJets::_pseudojetsToJets(), FastJets::calc(), and FastJets::reset().
FastJet external plugin. Definition at line 218 of file FastJets.hh. Referenced by FastJets::_init1(), and FastJets::_init3().
Flag to determine whether or not the VFS wrapper is to be used. Definition at line 144 of file JetAlg.hh. Referenced by FastJets::compare(), FastJets::project(), and JetAlg::useInvisibles().
Map of vectors of y scales. This is mutable so we can use caching/lazy evaluation. Definition at line 221 of file FastJets.hh. Referenced by FastJets::reset(), and FastJets::ySubJet(). The documentation for this class was generated from the following files: Generated on Thu Feb 6 2014 17:38:59 for The Rivet MC analysis system by 1.7.6.1 |