1 #ifndef RecoJets_JetAlgorithms_QJets_h
2 #define RecoJets_JetAlgorithms_QJets_h
5 #include "fastjet/JetDefinition.hh"
6 #include "fastjet/PseudoJet.hh"
7 #include "fastjet/ClusterSequence.hh"
12 #include "CLHEP/Random/RandomEngine.h"
23 bool operator() (
const JetDistance& lhs,
const JetDistance&rhs)
const {
return lhs.dij > rhs.dij;};
35 double d_ij(
const fastjet::PseudoJet& v1,
const fastjet::PseudoJet& v2)
const;
39 bool Prune(JetDistance& jd,fastjet::ClusterSequence &
cs);
47 bool Same(
const JetDistance& lhs,
const JetDistance& rhs);
50 Qjets(
double zcut,
double dcut_fctr,
double exp_min,
double exp_max,
double rigidity,
double truncation_fctr, CLHEP::HepRandomEngine* rnEngine) :
_rand_seed_set(
false),
61 void Cluster(fastjet::ClusterSequence &
cs);
void ComputeNewDistanceMeasures(fastjet::ClusterSequence &cs, unsigned int new_jet)
void Cluster(fastjet::ClusterSequence &cs)
auto_ptr< ClusterSequence > cs
bool JetUnmerged(int num) const
void SetRandSeed(unsigned int seed)
double ComputeNormalization(double dmin)
bool Prune(JetDistance &jd, fastjet::ClusterSequence &cs)
bool operator()(const JetDistance &lhs, const JetDistance &rhs) const
void computeDCut(fastjet::ClusterSequence &cs)
JetDistance GetNextDistance()
Qjets(double zcut, double dcut_fctr, double exp_min, double exp_max, double rigidity, double truncation_fctr, CLHEP::HepRandomEngine *rnEngine)
CLHEP::HepRandomEngine * _rnEngine
double d_ij(const fastjet::PseudoJet &v1, const fastjet::PseudoJet &v2) const
double ComputeMinimumDistance()
void ComputeAllDistances(const std::vector< fastjet::PseudoJet > &inp)
std::priority_queue< JetDistance, std::vector< JetDistance >, JetDistanceCompare > _distances
std::map< int, bool > _merged_jets
volatile std::atomic< bool > shutdown_flag false
bool Same(const JetDistance &lhs, const JetDistance &rhs)
bool JetsUnmerged(const JetDistance &jd) const