11 return _merged_jets.find(num) == _merged_jets.end();
15 return JetUnmerged(jd.
j1) && JetUnmerged(jd.
j2);
19 vector< pair<JetDistance, double> > popped_distances;
28 while(!_distances.empty()){
31 if(JetsUnmerged(dist)){
37 popped_distances.push_back(make_pair(dist,weight));
39 if(weight/norm < _truncation_fctr)
44 double rand(Rand()), tot_weight(0.);
45 for(vector<pair<JetDistance, double> >::iterator it = popped_distances.begin(); it != popped_distances.end(); it++){
46 tot_weight += (*it).second/norm;
47 if(tot_weight >=
rand){
54 for(vector<pair<JetDistance, double> >::reverse_iterator it = popped_distances.rbegin(); it != popped_distances.rend(); it++)
55 if(JetsUnmerged((*it).first))
56 _distances.push((*it).first);
65 ComputeAllDistances(cs.jets());
68 while(!_distances.empty() && jd.dij != -1.){
72 _merged_jets[jd.j1] =
true;
73 _merged_jets[jd.j2] =
true;
76 cs.plugin_record_ij_recombination(jd.j1, jd.j2, 1., new_jet);
77 if(!JetUnmerged(new_jet))
78 edm::LogError(
"QJets Clustering") <<
"Problem with FastJet::plugin_record_ij_recombination";
79 ComputeNewDistanceMeasures(cs,new_jet);
81 double j1pt = cs.jets()[jd.j1].perp();
82 double j2pt = cs.jets()[jd.j2].perp();
85 _merged_jets[jd.j2] =
true;
86 cs.plugin_record_iB_recombination(jd.j2, 1.);
89 _merged_jets[jd.j1] =
true;
90 cs.plugin_record_iB_recombination(jd.j1, 1.);
93 jd = GetNextDistance();
97 int num_merged_final(0);
98 for(
unsigned int i = 0 ;
i < cs.jets().size();
i++)
101 cs.plugin_record_iB_recombination(
i,1.);
104 if(!(num_merged_final < 2))
105 edm::LogError(
"QJets Clustering") <<
"More than 1 jet remains after reclustering";
110 fastjet::PseudoJet sum(0.,0.,0.,0.);
111 for(vector<fastjet::PseudoJet>::const_iterator it = cs.jets().begin(); it != cs.jets().end(); it++)
113 _dcut = 2. * _dcut_fctr * sum.m()/sum.perp();
117 double pt1 = cs.jets()[jd.
j1].perp();
118 double pt2 = cs.jets()[jd.
j2].perp();
119 fastjet::PseudoJet sum_jet = cs.jets()[jd.
j1]+cs.jets()[jd.
j2];
120 double sum_pt = sum_jet.perp();
121 double z =
min(pt1,pt2)/sum_pt;
122 double d2 = cs.jets()[jd.
j1].plain_distance(cs.jets()[jd.
j2]);
123 return (d2 > _dcut*_dcut) && (z < _zcut);
127 for(
unsigned int i = 0 ;
i < inp.size()-1;
i++){
129 for(
unsigned int j =
i+1 ;
j < inp.size();
j++){
134 jd.
dij = d_ij(inp[
i],inp[
j]);
143 for(
unsigned int i = 0;
i < cs.jets().size();
i++)
144 if(JetUnmerged(
i) &&
i != new_jet){
148 jd.
dij = d_ij(cs.jets()[jd.
j1], cs.jets()[jd.
j2]);
153 double Qjets::d_ij(
const fastjet::PseudoJet& v1,
const fastjet::PseudoJet& v2)
const{
154 double p1 = v1.perp();
155 double p2 = v2.perp();
156 double ret =
pow(
min(p1,p2),_exp_min) *
pow(
max(p1,p2),_exp_max) * v1.squared_distance(v2);
168 ret = _rnEngine->flat();
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)
bool Prune(JetDistance &jd, fastjet::ClusterSequence &cs)
void computeDCut(fastjet::ClusterSequence &cs)
JetDistance GetNextDistance()
double d_ij(const fastjet::PseudoJet &v1, const fastjet::PseudoJet &v2) const
void ComputeAllDistances(const std::vector< fastjet::PseudoJet > &inp)
Power< A, B >::type pow(const A &a, const B &b)
bool JetsUnmerged(const JetDistance &jd) const