File indexing completed on 2025-05-12 09:05:02
0001
0002 #ifndef RIVET_FastJets_HH
0003 #define RIVET_FastJets_HH
0004
0005 #include "Rivet/Jet.hh"
0006 #include "Rivet/Particle.hh"
0007 #include "Rivet/Projection.hh"
0008 #include "Rivet/Projections/JetFinder.hh"
0009 #include "Rivet/Projections/FinalState.hh"
0010 #include "Rivet/Tools/RivetFastJet.hh"
0011
0012 #include "fastjet/SISConePlugin.hh"
0013 #include "fastjet/ATLASConePlugin.hh"
0014 #include "fastjet/CMSIterativeConePlugin.hh"
0015 #include "fastjet/CDFJetCluPlugin.hh"
0016 #include "fastjet/CDFMidPointPlugin.hh"
0017 #include "fastjet/D0RunIIConePlugin.hh"
0018 #include "fastjet/TrackJetPlugin.hh"
0019 #include "fastjet/JadePlugin.hh"
0020
0021 #include "Rivet/Projections/PxConePlugin.hh"
0022 #include "Rivet/Tools/TypeTraits.hh"
0023
0024 namespace Rivet {
0025
0026
0027
0028 class FastJets : public JetFinder {
0029 public:
0030
0031 using JetFinder::operator =;
0032
0033
0034
0035
0036
0037
0038
0039
0040 FastJets(const FinalState& fsp,
0041 const fastjet::JetDefinition& jdef,
0042 JetMuons usemuons=JetMuons::ALL,
0043 JetInvisibles useinvis=JetInvisibles::NONE,
0044 fastjet::AreaDefinition* adef=nullptr)
0045 : JetFinder(fsp, usemuons, useinvis), _jdef(jdef),
0046 _adef(adef), _cuts(Cuts::OPEN)
0047 {
0048 _initBase();
0049 }
0050
0051
0052
0053
0054 FastJets(const FinalState& fsp,
0055 const fastjet::JetDefinition& jdef,
0056 const Cut& c,
0057 JetMuons usemuons=JetMuons::ALL,
0058 JetInvisibles useinvis=JetInvisibles::NONE,
0059 fastjet::AreaDefinition* adef=nullptr)
0060 : JetFinder(fsp, usemuons, useinvis), _jdef(jdef),
0061 _adef(adef), _cuts(c)
0062 {
0063 _initBase();
0064 }
0065
0066
0067
0068
0069 FastJets(const FinalState& fsp,
0070 const fastjet::JetDefinition& jdef,
0071 fastjet::AreaDefinition* adef,
0072 JetMuons usemuons=JetMuons::ALL,
0073 JetInvisibles useinvis=JetInvisibles::NONE)
0074 : FastJets(fsp, jdef, usemuons, useinvis, adef)
0075 { }
0076
0077
0078
0079
0080
0081 FastJets(const FinalState& fsp,
0082 const fastjet::JetDefinition& jdef,
0083 fastjet::AreaDefinition* adef,
0084 const Cut& c,
0085 JetMuons usemuons=JetMuons::ALL,
0086 JetInvisibles useinvis=JetInvisibles::NONE)
0087 : FastJets(fsp, jdef, c, usemuons, useinvis, adef)
0088 { }
0089
0090
0091
0092
0093
0094 FastJets(const FinalState& fsp,
0095 fastjet::JetAlgorithm type,
0096 fastjet::RecombinationScheme recom, double rparameter,
0097 JetMuons usemuons=JetMuons::ALL,
0098 JetInvisibles useinvis=JetInvisibles::NONE,
0099 fastjet::AreaDefinition* adef=nullptr)
0100 : FastJets(fsp, fastjet::JetDefinition(type, rparameter, recom), usemuons, useinvis, adef)
0101 { }
0102
0103
0104
0105
0106 FastJets(const FinalState& fsp,
0107 fastjet::JetAlgorithm type,
0108 fastjet::RecombinationScheme recom, double rparameter,
0109 const Cut& c,
0110 JetMuons usemuons=JetMuons::ALL,
0111 JetInvisibles useinvis=JetInvisibles::NONE,
0112 fastjet::AreaDefinition* adef=nullptr)
0113 : FastJets(fsp, fastjet::JetDefinition(type, rparameter, recom), c, usemuons, useinvis, adef)
0114 { }
0115
0116
0117
0118
0119 FastJets(const FinalState& fsp,
0120 fastjet::JetAlgorithm type,
0121 fastjet::RecombinationScheme recom, double rparameter,
0122 fastjet::AreaDefinition* adef,
0123 JetMuons usemuons=JetMuons::ALL,
0124 JetInvisibles useinvis=JetInvisibles::NONE)
0125 : FastJets(fsp, type, recom, rparameter, usemuons, useinvis, adef)
0126 { }
0127
0128
0129
0130
0131 FastJets(const FinalState& fsp,
0132 fastjet::JetAlgorithm type,
0133 fastjet::RecombinationScheme recom, double rparameter,
0134 fastjet::AreaDefinition* adef,
0135 const Cut& c,
0136 JetMuons usemuons=JetMuons::ALL,
0137 JetInvisibles useinvis=JetInvisibles::NONE)
0138 : FastJets(fsp, type, recom, rparameter, c, usemuons, useinvis, adef)
0139 { }
0140
0141
0142
0143
0144 FastJets(const FinalState& fsp,
0145 fastjet::JetDefinition::Plugin* plugin,
0146 JetMuons usemuons=JetMuons::ALL,
0147 JetInvisibles useinvis=JetInvisibles::NONE,
0148 fastjet::AreaDefinition* adef=nullptr)
0149 : FastJets(fsp, fastjet::JetDefinition(plugin), usemuons, useinvis, adef)
0150 {
0151 _plugin.reset(plugin);
0152 }
0153
0154
0155
0156
0157 FastJets(const FinalState& fsp,
0158 fastjet::JetDefinition::Plugin* plugin,
0159 const Cut& c,
0160 JetMuons usemuons=JetMuons::ALL,
0161 JetInvisibles useinvis=JetInvisibles::NONE,
0162 fastjet::AreaDefinition* adef=nullptr)
0163 : FastJets(fsp, fastjet::JetDefinition(plugin), c, usemuons, useinvis, adef)
0164 {
0165 _plugin.reset(plugin);
0166 }
0167
0168
0169
0170
0171 FastJets(const FinalState& fsp,
0172 fastjet::JetDefinition::Plugin* plugin,
0173 fastjet::AreaDefinition* adef,
0174 JetMuons usemuons=JetMuons::ALL,
0175 JetInvisibles useinvis=JetInvisibles::NONE)
0176 : FastJets(fsp, plugin, usemuons, useinvis, adef)
0177 { }
0178
0179
0180
0181
0182 FastJets(const FinalState& fsp,
0183 fastjet::JetDefinition::Plugin* plugin,
0184 fastjet::AreaDefinition* adef,
0185 const Cut& c,
0186 JetMuons usemuons=JetMuons::ALL,
0187 JetInvisibles useinvis=JetInvisibles::NONE)
0188 : FastJets(fsp, plugin, c, usemuons, useinvis, adef)
0189 { }
0190
0191
0192
0193
0194
0195
0196
0197
0198 FastJets(const FinalState& fsp,
0199 JetAlg alg, double rparameter,
0200 JetMuons usemuons=JetMuons::ALL,
0201 JetInvisibles useinvis=JetInvisibles::NONE,
0202 fastjet::AreaDefinition* adef=nullptr,
0203 double seed_threshold=1.0)
0204 : JetFinder(fsp, usemuons, useinvis),
0205 _adef(adef), _cuts(Cuts::OPEN)
0206 {
0207 _initBase();
0208 _initJdef(alg, rparameter, seed_threshold);
0209 }
0210
0211
0212
0213
0214
0215
0216
0217
0218 FastJets(const FinalState& fsp,
0219 JetAlg alg, double rparameter,
0220 const Cut& c,
0221 JetMuons usemuons=JetMuons::ALL,
0222 JetInvisibles useinvis=JetInvisibles::NONE,
0223 fastjet::AreaDefinition* adef=nullptr,
0224 double seed_threshold=1.0)
0225 : JetFinder(fsp, usemuons, useinvis), _adef(adef), _cuts(c)
0226 {
0227 _initBase();
0228 _initJdef(alg, rparameter, seed_threshold);
0229 }
0230
0231
0232
0233
0234 RIVET_DEFAULT_PROJ_CLONE(FastJets);
0235
0236
0237
0238
0239
0240 using Projection::operator =;
0241
0242
0243
0244
0245
0246
0247 static PseudoJets mkClusterInputs(const Particles& fsparticles, const Particles& tagparticles=Particles());
0248
0249 static Jet mkJet(const PseudoJet& pj, const Particles& fsparticles, const Particles& tagparticles=Particles());
0250
0251 static Jets mkJets(const PseudoJets& pjs, const Particles& fsparticles, const Particles& tagparticles=Particles());
0252
0253
0254
0255
0256
0257 void reset();
0258
0259
0260
0261
0262
0263
0264
0265
0266
0267 void useJetArea(fastjet::AreaDefinition* adef) {
0268 _adef.reset(adef);
0269 }
0270
0271
0272 void clearJetArea() {
0273 _adef.reset();
0274 }
0275
0276
0277
0278
0279
0280
0281
0282
0283
0284
0285
0286 void addTrf(fastjet::Transformer* trf) {
0287 _trfs.push_back(shared_ptr<fastjet::Transformer>(trf));
0288 }
0289
0290
0291
0292
0293
0294 template<typename TRFS, typename TRF=typename TRFS::value_type>
0295 typename std::enable_if<Derefable<TRF>::value, void>::type
0296 addTrfs(const TRFS& trfs) {
0297 for (auto& trf : trfs) addTrf(trf);
0298 }
0299
0300
0301 void clearTrfs() {
0302 _trfs.clear();
0303 }
0304
0305
0306
0307
0308
0309
0310
0311
0312 Jets _jets() const;
0313
0314
0315
0316 PseudoJets pseudojets(double ptmin=0.0) const;
0317
0318
0319 PseudoJets pseudojetsByPt(double ptmin=0.0) const {
0320 return sorted_by_pt(pseudojets(ptmin));
0321 }
0322
0323
0324 PseudoJets pseudojetsByE(double ptmin=0.0) const {
0325 return sorted_by_E(pseudojets(ptmin));
0326 }
0327
0328
0329 PseudoJets pseudojetsByRapidity(double ptmin=0.0) const {
0330 return sorted_by_rapidity(pseudojets(ptmin));
0331 }
0332
0333
0334
0335
0336
0337
0338
0339
0340
0341 const shared_ptr<fastjet::ClusterSequence> clusterSeq() const {
0342 return _cseq;
0343 }
0344
0345
0346
0347 const shared_ptr<fastjet::ClusterSequenceArea> clusterSeqArea() const {
0348 return areaDef() ? dynamic_pointer_cast<fastjet::ClusterSequenceArea>(_cseq) : nullptr;
0349 }
0350
0351
0352 const fastjet::JetDefinition& jetDef() const {
0353 return _jdef;
0354 }
0355
0356
0357
0358
0359
0360 const shared_ptr<fastjet::AreaDefinition> areaDef() const {
0361 return _adef;
0362 }
0363
0364
0365
0366
0367 protected:
0368
0369
0370 void _initBase();
0371 void _initJdef(JetAlg alg, double rparameter, double seed_threshold);
0372
0373 protected:
0374
0375
0376 void project(const Event& e);
0377
0378
0379 CmpState compare(const Projection& p) const;
0380
0381 public:
0382
0383
0384 void calc(const Particles& fsparticles, const Particles& tagparticles=Particles());
0385
0386
0387 protected:
0388
0389
0390 fastjet::JetDefinition _jdef;
0391
0392
0393 std::shared_ptr<fastjet::AreaDefinition> _adef;
0394
0395
0396 std::shared_ptr<fastjet::ClusterSequence> _cseq;
0397
0398
0399 Cut _cuts;
0400
0401
0402 std::shared_ptr<fastjet::JetDefinition::Plugin> _plugin;
0403
0404
0405 std::vector< std::shared_ptr<fastjet::Transformer> > _trfs;
0406
0407
0408 mutable std::map<int, vector<double> > _yscales;
0409
0410
0411 Particles _fsparticles, _tagparticles;
0412
0413 };
0414
0415 }
0416
0417 #endif