FastJet  3.1.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Modules Pages
CDFJetCluPlugin.cc
1 //FJSTARTHEADER
2 // $Id: CDFJetCluPlugin.cc 3433 2014-07-23 08:17:03Z salam $
3 //
4 // Copyright (c) 2005-2014, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
5 //
6 //----------------------------------------------------------------------
7 // This file is part of FastJet.
8 //
9 // FastJet is free software; you can redistribute it and/or modify
10 // it under the terms of the GNU General Public License as published by
11 // the Free Software Foundation; either version 2 of the License, or
12 // (at your option) any later version.
13 //
14 // The algorithms that underlie FastJet have required considerable
15 // development. They are described in the original FastJet paper,
16 // hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
17 // FastJet as part of work towards a scientific publication, please
18 // quote the version you use and include a citation to the manual and
19 // optionally also to hep-ph/0512210.
20 //
21 // FastJet is distributed in the hope that it will be useful,
22 // but WITHOUT ANY WARRANTY; without even the implied warranty of
23 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24 // GNU General Public License for more details.
25 //
26 // You should have received a copy of the GNU General Public License
27 // along with FastJet. If not, see <http://www.gnu.org/licenses/>.
28 //----------------------------------------------------------------------
29 //FJENDHEADER
30 
31 #include "fastjet/CDFJetCluPlugin.hh"
32 #include "fastjet/ClusterSequence.hh"
33 #include <sstream>
34 #include <cassert>
35 
36 // CDF stuff
37 #include "JetCluAlgorithm.hh"
38 #include "PhysicsTower.hh"
39 #include "Cluster.hh"
40 
41 FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
42 
43 using namespace std;
44 using namespace cdf;
45 
46 bool CDFJetCluPlugin::_first_time = true;
47 
48 string CDFJetCluPlugin::description () const {
49  ostringstream desc;
50 
51  desc << "CDF JetClu jet algorithm with "
52  << "seed_threshold = " << seed_threshold () << ", "
53  << "cone_radius = " << cone_radius () << ", "
54  << "adjacency_cut = " << adjacency_cut () << ", "
55  << "max_iterations = " << max_iterations () << ", "
56  << "iratch = " << iratch () << ", "
57  << "overlap_threshold = " << overlap_threshold () ;
58 
59  return desc.str();
60 }
61 
62 
63 void CDFJetCluPlugin::run_clustering(ClusterSequence & clust_seq) const {
64  // print a banner if we run this for the first time
65  _print_banner(clust_seq.fastjet_banner_stream());
66 
67  // create the physics towers needed by the CDF code
68  vector<PhysicsTower> towers;
69  towers.reserve(clust_seq.jets().size());
70 
71  // create a map to identify jets (actually just the input particles)...
72  //map<double,int> jetmap;
73 
74  for (unsigned i = 0; i < clust_seq.jets().size(); i++) {
75  PseudoJet particle(clust_seq.jets()[i]);
76  //_insert_unique(particle, jetmap);
77  LorentzVector fourvect(particle.px(), particle.py(),
78  particle.pz(), particle.E());
79  PhysicsTower tower(fourvect);
80  // add tracking information for later
81  tower.fjindex = i;
82  towers.push_back(tower);
83  }
84 
85  // prepare the CDF algorithm
86  JetCluAlgorithm j(seed_threshold(), cone_radius(), adjacency_cut(),
87  max_iterations(), iratch(), overlap_threshold());
88 
89  // run the CDF algorithm
90  std::vector<Cluster> jets;
91  j.run(towers,jets);
92 
93 
94  // now transfer the jets back into our own structure -- we will
95  // mimic the cone code with a sequential recombination sequence in
96  // which the jets are built up by adding one particle at a time
97 
98  // NB: with g++-4.0, the reverse iterator code gave problems, so switch
99  // to indices instead
100  //for(vector<Cluster>::const_reverse_iterator jetIter = jets.rbegin();
101  // jetIter != jets.rend(); jetIter++) {
102  // const vector<PhysicsTower> & tower_list = jetIter->towerList;
103  // int jet_k = jetmap[tower_list[0].fourVector.E];
104  //
105  // int ntow = int(jetIter->towerList.size());
106 
107  for(int iCDFjets = jets.size()-1; iCDFjets >= 0; iCDFjets--) {
108 
109  const vector<PhysicsTower> & tower_list = jets[iCDFjets].towerList;
110  int ntow = int(tower_list.size());
111 
112  // 2008-09-04: sort the towerList (according to fjindex) so as
113  // to have a consistent order for particles in jet
114  // (necessary because addition of ultra-soft particles
115  // sometimes often modifies the order, while maintaining
116  // the same overall set)
117  vector<int> jc_indices(ntow);
118  vector<double> fj_indices(ntow); // use double: benefit from existing routine
119  for (int itow = 0; itow < ntow; itow++) {
120  jc_indices[itow] = itow;
121  fj_indices[itow] = tower_list[itow].fjindex;
122  }
123  sort_indices(jc_indices, fj_indices);
124 
125  int jet_k = tower_list[jc_indices[0]].fjindex;
126 
127  for (int itow = 1; itow < ntow; itow++) {
128  if (tower_list[jc_indices[itow]].Et() > 1e-50) {
129  }
130  int jet_i = jet_k;
131  // retrieve our index for the jet
132  int jet_j;
133  jet_j = tower_list[jc_indices[itow]].fjindex;
134 
135  // safety check
136  assert (jet_j >= 0 && jet_j < int(towers.size()));
137 
138  // do a fake recombination step with dij=0
139  double dij = 0.0;
140 
141  // JetClu does E-scheme recombination so we can stick with the
142  // simple option
143  clust_seq.plugin_record_ij_recombination(jet_i, jet_j, dij, jet_k);
144 
145  }
146 
147  // NB: put a sensible looking d_iB just to be nice...
148  double d_iB = clust_seq.jets()[jet_k].perp2();
149  clust_seq.plugin_record_iB_recombination(jet_k, d_iB);
150  }
151 
152 
153  // following code is for testing only
154  //cout << endl;
155  //for(vector<Cluster>::const_iterator jetIter = jets.begin();
156  // jetIter != jets.end(); jetIter++) {
157  // cout << jetIter->fourVector.pt() << " " << jetIter->fourVector.y() << endl;
158  //}
159  //cout << "-----------------------------------------------------\n";
160  //vector<PseudoJet> ourjets(clust_seq.inclusive_jets());
161  //for (vector<PseudoJet>::const_iterator ourjet = ourjets.begin();
162  // ourjet != ourjets.end(); ourjet++) {
163  // cout << ourjet->perp() << " " << ourjet->rap() << endl;
164  //}
165  //cout << endl;
166 }
167 
168 
169 // print a banner for reference to the 3rd-party code
170 void CDFJetCluPlugin::_print_banner(ostream *ostr) const{
171  if (! _first_time) return;
172  _first_time=false;
173 
174  // make sure the user has not set the banner stream to NULL
175  if (!ostr) return;
176 
177  (*ostr) << "#-------------------------------------------------------------------------" << endl;
178  (*ostr) << "# You are running the CDF JetClu plugin for FastJet " << endl;
179  (*ostr) << "# This is based on an implementation provided by Joey Huston. " << endl;
180  (*ostr) << "# If you use this plugin, please cite " << endl;
181  (*ostr) << "# F. Abe et al. [CDF Collaboration], Phys. Rev. D 45 (1992) 1448. " << endl;
182  (*ostr) << "# in addition to the usual FastJet reference. " << endl;
183  (*ostr) << "#-------------------------------------------------------------------------" << endl;
184 
185  // make sure we really have the output done.
186  ostr->flush();
187 }
188 
189 FASTJET_END_NAMESPACE // defined in fastjet/internal/base.hh
const std::vector< PseudoJet > & jets() const
allow the user to access the internally stored _jets() array, which contains both the initial particl...
deals with clustering
void plugin_record_ij_recombination(int jet_i, int jet_j, double dij, int &newjet_k)
record the fact that there has been a recombination between jets()[jet_i] and jets()[jet_k], with the specified dij, and return the index (newjet_k) allocated to the new jet, whose momentum is assumed to be the 4-vector sum of that of jet_i and jet_j
void plugin_record_iB_recombination(int jet_i, double diB)
record the fact that there has been a recombination between jets()[jet_i] and the beam...
static std::ostream * fastjet_banner_stream()
returns a pointer to the stream to be used to print banners (cout by default).
Class to contain pseudojets, including minimal information of use to jet-clustering routines...
Definition: PseudoJet.hh:67