File indexing completed on 2024-09-07 04:37:32
0001 #ifndef RecoJets_JetAlgorithms_QJETSPLUGIN_h
0002 #define RecoJets_JetAlgorithms_QJETSPLUGIN_h
0003 #include "fastjet/JetDefinition.hh"
0004 #include "fastjet/PseudoJet.hh"
0005 #include "fastjet/ClusterSequence.hh"
0006 #include "RecoJets/JetAlgorithms/interface/Qjets.h"
0007
0008 class QjetsPlugin : public fastjet::JetDefinition::Plugin {
0009 private:
0010 bool _rand_seed_set;
0011 unsigned int _seed;
0012 int _truncated_length;
0013 double _zcut, _dcut_fctr, _exp_min, _exp_max, _rigidity, _truncation_fctr;
0014 CLHEP::HepRandomEngine* _rnEngine;
0015
0016 public:
0017 QjetsPlugin(
0018 double zcut, double dcut_fctr, double exp_min, double exp_max, double rigidity, double truncation_fctr = 0.)
0019 : _rand_seed_set(false),
0020 _zcut(zcut),
0021 _dcut_fctr(dcut_fctr),
0022 _exp_min(exp_min),
0023 _exp_max(exp_max),
0024 _rigidity(rigidity),
0025 _truncation_fctr(truncation_fctr),
0026 _rnEngine(nullptr) {}
0027 void SetRandSeed(unsigned int seed);
0028 void SetRNEngine(CLHEP::HepRandomEngine* rnEngine) { _rnEngine = rnEngine; };
0029 double R() const override;
0030 std::string description() const override;
0031 void run_clustering(fastjet::ClusterSequence& cs) const override;
0032 };
0033 #endif