00001 #ifndef jetsubstructureutils_qjets_header
00002 #define jetsubstructureutils_qjets_header
00003
00004 #include <queue>
00005 #include <vector>
00006 #include <list>
00007 #include <algorithm>
00008 #include "fastjet/JetDefinition.hh"
00009 #include "fastjet/PseudoJet.hh"
00010 #include "fastjet/ClusterSequence.hh"
00011
00012 using namespace std;
00013
00014 namespace JetSubStructureUtils {
00015 struct jet_distance{
00016 double dij;
00017 int j1;
00018 int j2;
00019 };
00020
00021 class JetDistanceCompare {
00022 public:
00023 JetDistanceCompare() {};
00024 bool operator() (const jet_distance& lhs, const jet_distance &rhs) const {return lhs.dij > rhs.dij;};
00025 };
00026
00027 class Qjets {
00028 private:
00029 bool _rand_seed_set;
00030 unsigned int _seed;
00031 double _zcut, _dcut, _dcut_fctr, _exp_min, _exp_max, _rigidity, _truncation_fctr;
00032 std::map<int,bool> _merged_jets;
00033 std::priority_queue <jet_distance, std::vector<jet_distance>, JetDistanceCompare> _distances;
00034
00035 double d_ij(const fastjet::PseudoJet& v1, const fastjet::PseudoJet& v2) const;
00036 void ComputeDCut(fastjet::ClusterSequence & cs);
00037
00038 double Rand();
00039 bool Prune(jet_distance& jd,fastjet::ClusterSequence & cs);
00040 bool JetsUnmerged(const jet_distance& jd) const;
00041 bool JetUnmerged(int num) const;
00042 void ComputeNewDistanceMeasures(fastjet::ClusterSequence & cs, int new_jet);
00043 void ComputeAllDistances(const std::vector<fastjet::PseudoJet>& inp);
00044 double ComputeMinimumDistance();
00045 double ComputeNormalization(double dmin);
00046 jet_distance GetNextDistance();
00047 bool Same(const jet_distance& lhs, const jet_distance& rhs);
00048
00049 public:
00050 Qjets(double zcut, double dcut_fctr, double exp_min, double exp_max, double rigidity, double truncation_fctr);
00051 void Cluster(fastjet::ClusterSequence & cs);
00052 void SetRandSeed(unsigned int seed);
00053 };
00054 }
00055 #endif