FastJet  3.0.2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
fastjet_timing_plugins.cc
1 //STARTHEADER
2 // $Id: fastjet_timing_plugins.cc 2806 2011-12-01 17:21:00Z salam $
3 //
4 // Copyright (c) 2005-2011, 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 and are described in hep-ph/0512210. If you use
16 // FastJet as part of work towards a scientific publication, please
17 // include a citation to the FastJet paper.
18 //
19 // FastJet is distributed in the hope that it will be useful,
20 // but WITHOUT ANY WARRANTY; without even the implied warranty of
21 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
22 // GNU General Public License for more details.
23 //
24 // You should have received a copy of the GNU General Public License
25 // along with FastJet. If not, see <http://www.gnu.org/licenses/>.
26 //----------------------------------------------------------------------
27 //ENDHEADER
28 
29 
30 //----------------------------------------------------------------------
31 /// fastjet_timing.cc: Program to help time and test the fastjet package
32 ///
33 /// It reads files containing multiple events in the format
34 /// p1x p1y p1z E1
35 /// p2x p2y p2z E2
36 /// ...
37 /// #END
38 ///
39 /// An example input file containing 10 events is included as
40 /// data/Pythia-PtMin1000-LHC-10ev.dat
41 ///
42 /// Usage:
43 /// fastjet_timing [-strategy NUMBER] [-repeat nrepeats] [-massive] \
44 /// [-combine nevents] [-r Rparameter] [-incl ptmin] [...] \
45 /// < data_file
46 ///
47 /// where the clustering can be repeated to aid timing and multiple
48 /// events can be combined to get to larger multiplicities. Some options:
49 ///
50 /// Options for reading
51 /// -------------------
52 ///
53 /// -nev n number of events to run
54 ///
55 /// -combine n for combining multiple events from the data file in order
56 /// to get a single high-multipicity event to run.
57 ///
58 /// -massless read in only the 3-momenta and deduce energies assuming
59 /// that particles are massless
60 ///
61 /// -dense adds dense ghost coverage
62 ///
63 /// -repeat n repeats each event n times
64 ///
65 /// -nhardest n keep only the n hardest particles in the event
66 ///
67 /// -file name read from the corresponding file rather than stdin.
68 /// (The file will be reopened for each new jet alg.; in
69 /// constrast, if you use stdin, each new alg will take a
70 /// new event).
71 ///
72 /// Output Options
73 /// --------------
74 ///
75 /// -incl ptmin output of all inclusive jets with pt > ptmin is obtained
76 /// with the -incl option.
77 ///
78 /// -repeat-incl ptmin
79 /// same as -incl ptmin but do it for each repetition
80 /// of the clustering
81 ///
82 /// -excld dcut output of all exclusive jets as obtained in a clustering
83 /// with dcut
84 ///
85 /// -excly ycut output of all exclusive jets as obtained in a clustering
86 /// with ycut
87 ///
88 /// -excln n output of clustering to n exclusive jets
89 ///
90 /// -ee-print print things as px,py,pz,E
91 ///
92 /// -get-all-dij print out all dij values
93 /// -get-all-yij print out all yij values
94 ///
95 /// -const show jet constituents (works with excl jets)
96 ///
97 /// -write for writing out detailed clustering sequence (valuable
98 /// for testing purposes)
99 ///
100 /// -unique_write writes out the sequence of dij's according to the
101 /// "unique_history_order" (useful for verifying consistency
102 /// between different clustering strategies).
103 ///
104 /// -root file sends output to file that can be read in with the script in
105 /// root/ so as to show a lego-plot of the event
106 ///
107 /// -cones show extra info about internal steps for SISCone
108 ///
109 /// -area calculate areas. Additional options include
110 /// -area:active
111 /// -area:passive
112 /// -area:explicit
113 /// -area:voronoi Rfact
114 /// -area:repeat nrepeat
115 /// -ghost-area area
116 /// -ghost-maxrap maxrap
117 /// -area:fj2 place ghosts as in fj2
118 ///
119 /// -bkgd calculate the background density. Additional options include
120 /// -bkgd:csab use the old ClusterSequenceAreaBase methods
121 /// -bkgd:jetmedian use the new JetMedianBackgroundEstimator class
122 /// -bkgd:fj2 force jetmedian to calculate sigma as in fj2
123 /// -bkgd:gridmedian use GridMedianBackgroundEstimator with grid up to ghost_maxrap-ktR and grid spacing of 2ktR
124 ///
125 /// Algorithms
126 /// ----------
127 /// -all-algs runs all algorithms
128 ///
129 /// -kt switch to the longitudinally invariant kt algorithm
130 /// Note: this is the default one.
131 ///
132 /// -cam switch to the inclusive Cambridge/Aachen algorithm --
133 /// note that the option -excld dcut provides a clustering
134 /// up to the dcut which is the minimum squared
135 /// distance between any pair of jets.
136 ///
137 /// -antikt switch to the anti-kt clustering algorithm
138 ///
139 /// -genkt switch to the genkt algorithm
140 /// you can provide the parameter of the alg as an argument to
141 /// -genkt (1 by default)
142 ///
143 /// -eekt switch to the e+e- kt algorithm
144 ///
145 /// -eegenkt switch to the genkt algorithm
146 /// you can provide the parameter of the alg as an argument to
147 /// -ee_genkt (1 by default)
148 ///
149 /// plugins (don't delete this line)
150 ///
151 /// -pxcone switch to the PxCone jet algorithm
152 ///
153 /// -siscone switch to the SISCone jet algorithm (seedless cones)
154 /// -sisconespheri switch to the Spherical SISCone jet algorithm (seedless cones)
155 ///
156 /// -midpoint switch to CDF's midpoint code
157 /// -jetclu switch to CDF's jetclu code
158 ///
159 /// -d0runipre96cone switch to the D0RunIpre96Cone plugin
160 /// -d0runicone switch to the D0RunICone plugin
161 ///
162 /// -d0runiicone switch to D0's run II midpoint cone
163 ///
164 /// -trackjet switch to the TrackJet plugin
165 ///
166 /// -atlascone switch to the ATLASCone plugin
167 ///
168 /// -eecambridge switch to the EECambridge plugin
169 ///
170 /// -jade switch to the Jade plugin
171 ///
172 /// -cmsiterativecone switch to the CMSIterativeCone plugin
173 ///
174 /// -gridjet switch to the GridJet plugin
175 ///
176 /// end of plugins (don't delete this line)
177 ///
178 ///
179 /// Options for running algs
180 /// ------------------------
181 ///
182 /// -r sets the radius of the jet algorithm (default = 1.0)
183 ///
184 /// -overlap | -f sets the overlap fraction in cone algs with split-merge
185 ///
186 /// -seed sets the seed threshold
187 ///
188 /// -strategy N indicate stratgey from the enum fastjet::Strategy (see
189 /// fastjet/JetDefinition.hh).
190 ///
191 
192 
193 #include "fastjet/ClusterSequenceArea.hh"
194 #include "fastjet/tools/JetMedianBackgroundEstimator.hh"
195 #include "fastjet/tools/GridMedianBackgroundEstimator.hh"
196 #include "fastjet/Selector.hh"
197 #include<iostream>
198 #include<sstream>
199 #include<fstream>
200 #include<valarray>
201 #include<vector>
202 #include <cstdlib>
203 //#include<cstddef> // for size_t
204 #include "CmdLine.hh"
205 
206 // get info on how fastjet was configured
207 #include "fastjet/config.h"
208 
209 // include the installed plugins (don't delete this line)
210 #ifdef FASTJET_ENABLE_PLUGIN_SISCONE
211 #include "fastjet/SISConePlugin.hh"
212 #include "fastjet/SISConeSphericalPlugin.hh"
213 #endif
214 #ifdef FASTJET_ENABLE_PLUGIN_CDFCONES
215 #include "fastjet/CDFMidPointPlugin.hh"
216 #include "fastjet/CDFJetCluPlugin.hh"
217 #endif
218 #ifdef FASTJET_ENABLE_PLUGIN_PXCONE
219 #include "fastjet/PxConePlugin.hh"
220 #endif
221 #ifdef FASTJET_ENABLE_PLUGIN_D0RUNIICONE
222 #include "fastjet/D0RunIIConePlugin.hh"
223 #endif
224 #ifdef FASTJET_ENABLE_PLUGIN_TRACKJET
225 #include "fastjet/TrackJetPlugin.hh"
226 #endif
227 #ifdef FASTJET_ENABLE_PLUGIN_ATLASCONE
228 #include "fastjet/ATLASConePlugin.hh"
229 #endif
230 #ifdef FASTJET_ENABLE_PLUGIN_EECAMBRIDGE
231 #include "fastjet/EECambridgePlugin.hh"
232 #endif
233 #ifdef FASTJET_ENABLE_PLUGIN_JADE
234 #include "fastjet/JadePlugin.hh"
235 #endif
236 #ifdef FASTJET_ENABLE_PLUGIN_CMSITERATIVECONE
237 #include "fastjet/CMSIterativeConePlugin.hh"
238 #endif
239 #ifdef FASTJET_ENABLE_PLUGIN_D0RUNICONE
240 #include "fastjet/D0RunIpre96ConePlugin.hh"
241 #include "fastjet/D0RunIConePlugin.hh"
242 #endif
243 #ifdef FASTJET_ENABLE_PLUGIN_GRIDJET
244 #include "fastjet/GridJetPlugin.hh"
245 #endif
246 // end of installed plugins inclusion (don't delete this line)
247 
248 using namespace std;
249 
250 // to avoid excessive typing, define an abbreviation for the
251 // fastjet namespace
252 namespace fj = fastjet;
253 
254 inline double pow2(const double x) {return x*x;}
255 
256 // pretty print the jets and their subjets
257 void print_jets_and_sub (const vector<fj::PseudoJet> & jets, double dcut);
258 
259 string rootfile;
260 CmdLine * cmdline_p;
261 
262 bool do_areas;
263 
264 /// sort and pretty print jets, with exact behaviour depending on
265 /// whether ee_print is true or not
266 bool ee_print = false;
267 void print_jets(const vector<fj::PseudoJet> & jets, bool show_const = false);
268 
269 bool found_unavailable = false;
270 void is_unavailable(const string & algname) {
271  cerr << algname << " requested, but not available for this compilation" << endl;
272  found_unavailable = true;
273  //exit(0);
274 }
275 
276 
277 /// a program to test and time a range of algorithms as implemented or
278 /// wrapped in fastjet
279 int main (int argc, char ** argv) {
280 
282 
283  CmdLine cmdline(argc,argv);
284  cmdline_p = &cmdline;
285  // allow the use to specify the fj::Strategy either through the
286  // -clever or the -strategy options (both will take numerical
287  // values); the latter will override the former.
288  fj::Strategy strategy = fj::Strategy(cmdline.int_val("-strategy",
289  cmdline.int_val("-clever", fj::Best)));
290  int repeat = cmdline.int_val("-repeat",1);
291  int combine = cmdline.int_val("-combine",1);
292  bool write = cmdline.present("-write");
293  bool unique_write = cmdline.present("-unique_write");
294  bool hydjet = cmdline.present("-hydjet");
295  double ktR = cmdline.double_val("-r",1.0);
296  ktR = cmdline.double_val("-R",ktR); // allow -r and -R
297  double inclkt = cmdline.double_val("-incl",-1.0);
298  double repeatinclkt = cmdline.double_val("-repeat-incl",-1.0);
299  int excln = cmdline.int_val ("-excln",-1);
300  double excld = cmdline.double_val("-excld",-1.0);
301  double excly = cmdline.double_val("-excly",-1.0);
302  ee_print = cmdline.present("-ee-print");
303  bool get_all_dij = cmdline.present("-get-all-dij");
304  bool get_all_yij = cmdline.present("-get-all-yij");
305  double subdcut = cmdline.double_val("-subdcut",-1.0);
306  double etamax = cmdline.double_val("-etamax",1.0e305);
307  bool show_constituents = cmdline.present("-const");
308  bool massless = cmdline.present("-massless");
309  int nev = cmdline.int_val("-nev",1);
310  bool add_dense_coverage = cmdline.present("-dense");
311  double ghost_maxrap = cmdline.value("-ghost-maxrap",5.0);
312  bool all_algs = cmdline.present("-all-algs");
313 
314  fj::Selector particles_sel = (cmdline.present("-nhardest"))
315  ? fj::SelectorNHardest(cmdline.value<unsigned int>("-nhardest"))
316  : fj::SelectorIdentity();
317 
318  do_areas = cmdline.present("-area");
319  fj::AreaDefinition area_def;
320  if (do_areas) {
321  assert(!write); // it's incompatible
322  fj::GhostedAreaSpec ghost_spec(ghost_maxrap,
323  cmdline.value("-area:repeat", 1),
324  cmdline.value("-ghost-area", 0.01));
325  if (cmdline.present("-area:fj2")) ghost_spec.set_fj2_placement(true);
326  if (cmdline.present("-area:explicit")) {
327  area_def = fj::AreaDefinition(fj::active_area_explicit_ghosts, ghost_spec);
328  } else if (cmdline.present("-area:passive")) {
329  area_def = fj::AreaDefinition(fj::passive_area, ghost_spec);
330  } else if (cmdline.present("-area:voronoi")) {
331  double Rfact = cmdline.value<double>("-area:voronoi");
332  area_def = fj::AreaDefinition(fj::voronoi_area,
333  fj::VoronoiAreaSpec(Rfact));
334  } else {
335  cmdline.present("-area:active"); // allow, but do not require, arg
336  area_def = fj::AreaDefinition(fj::active_area, ghost_spec);
337  }
338  }
339  bool do_bkgd = cmdline.present("-bkgd"); // background estimation
340  bool do_bkgd_csab = false, do_bkgd_jetmedian = false, do_bkgd_fj2 = false;
341  bool do_bkgd_gridmedian = false;
342  fj::Selector bkgd_range;
343  if (do_bkgd) {
344  bkgd_range = fj::SelectorAbsRapMax(ghost_maxrap - ktR);
345  if (cmdline.present("-bkgd:csab")) {do_bkgd_csab = true;}
346  else if (cmdline.present("-bkgd:jetmedian")) {do_bkgd_jetmedian = true;
347  do_bkgd_fj2 = cmdline.present("-bkgd:fj2");
348  } else if (cmdline.present("-bkgd:gridmedian")) {do_bkgd_gridmedian = true;
349  } else {
350  throw fj::Error("with the -bkgd option, some particular background must be specified (csab or jetmedian)");
351  }
352  assert(do_areas || do_bkgd_gridmedian);
353  }
354 
355  bool show_cones = cmdline.present("-cones"); // only works for siscone
356 
357  // for cone algorithms
358  // allow -f and -overlap
359  double overlap_threshold = cmdline.double_val("-overlap",0.5);
360  overlap_threshold = cmdline.double_val("-f",overlap_threshold);
361  double seed_threshold = cmdline.double_val("-seed",1.0);
362 
363  // for ee algorithms, allow to specify ycut
364  double ycut = cmdline.double_val("-ycut",0.08);
365 
366  // for printing jets to a file for reading by root
367  rootfile = cmdline.value<string>("-root","");
368 
369  // out default scheme is the E_scheme
371 
372  // The following option causes the Cambridge algo to be used.
373  // Note that currently the only output that works sensibly here is
374  // "-incl 0"
375  vector<fj::JetDefinition> jet_defs;
376  if (all_algs || cmdline.present("-cam") || cmdline.present("-CA")) {
377  jet_defs.push_back( fj::JetDefinition(fj::cambridge_algorithm, ktR, scheme, strategy));
378  }
379  if (all_algs || cmdline.present("-antikt")) {
380  jet_defs.push_back( fj::JetDefinition(fj::antikt_algorithm, ktR, scheme, strategy));
381  }
382  if (all_algs || cmdline.present("-genkt")) {
383  double p;
384  if (cmdline.present("-genkt")) p = cmdline.value<double>("-genkt");
385  else p = -0.5;
386  jet_defs.push_back( fj::JetDefinition(fj::genkt_algorithm, ktR, p, scheme, strategy));
387  }
388  if (all_algs || cmdline.present("-eekt")) {
389  jet_defs.push_back( fj::JetDefinition(fj::ee_kt_algorithm));
390  }
391  if (all_algs || cmdline.present("-eegenkt")) {
392  double p;
393  if (cmdline.present("-eegenkt")) p = cmdline.value<double>("-eegenkt");
394  else p = -0.5;
395  jet_defs.push_back( fj::JetDefinition(fj::ee_genkt_algorithm, ktR, p, scheme, strategy));
396 
397 // checking if one asks to run a plugin (don't delete this line)
398  }
399  if (all_algs || cmdline.present("-midpoint")) {
400 #ifdef FASTJET_ENABLE_PLUGIN_CDFCONES
401  typedef fj::CDFMidPointPlugin MPPlug; // for brevity
402  double cone_area_fraction = 1.0;
403  int max_pair_size = 2;
404  int max_iterations = 100;
405  MPPlug::SplitMergeScale sm_scale = MPPlug::SM_pt;
406  if (cmdline.present("-sm-pttilde")) sm_scale = MPPlug::SM_pttilde;
407  if (cmdline.present("-sm-pt")) sm_scale = MPPlug::SM_pt; // default
408  if (cmdline.present("-sm-mt")) sm_scale = MPPlug::SM_mt;
409  if (cmdline.present("-sm-Et")) sm_scale = MPPlug::SM_Et;
410  jet_defs.push_back( fj::JetDefinition( new fj::CDFMidPointPlugin (
411  seed_threshold, ktR,
412  cone_area_fraction, max_pair_size,
413  max_iterations, overlap_threshold,
414  sm_scale)));
415 #else // FASTJET_ENABLE_PLUGIN_CDFCONES
416  is_unavailable("midpoint");
417 #endif // FASTJET_ENABLE_PLUGIN_CDFCONES
418  }
419  if (all_algs || cmdline.present("-pxcone")) {
420 #ifdef FASTJET_ENABLE_PLUGIN_PXCONE
421  double min_jet_energy = 5.0;
422  jet_defs.push_back( fj::JetDefinition( new fj::PxConePlugin (
423  ktR, min_jet_energy,
424  overlap_threshold)));
425 #else // FASTJET_ENABLE_PLUGIN_PXCONE
426  is_unavailable("pxcone");
427 #endif // FASTJET_ENABLE_PLUGIN_PXCONE
428  }
429  if (all_algs || cmdline.present("-jetclu")) {
430 #ifdef FASTJET_ENABLE_PLUGIN_CDFCONES
431  jet_defs.push_back( fj::JetDefinition( new fj::CDFJetCluPlugin (
432  ktR, overlap_threshold, seed_threshold)));
433 #else // FASTJET_ENABLE_PLUGIN_CDFCONES
434  is_unavailable("pxcone");
435 #endif // FASTJET_ENABLE_PLUGIN_CDFCONES
436  }
437  if (all_algs || cmdline.present("-siscone") || cmdline.present("-sisconespheri")) {
438 #ifdef FASTJET_ENABLE_PLUGIN_SISCONE
439  typedef fj::SISConePlugin SISPlug; // for brevity
440  int npass = cmdline.value("-npass",0);
441  if (all_algs || cmdline.present("-siscone")) {
442  double sisptmin = cmdline.value("-sisptmin",0.0);
443  SISPlug * plugin = new SISPlug (ktR, overlap_threshold,npass,sisptmin);
444  if (cmdline.present("-sm-pt")) plugin->set_split_merge_scale(SISPlug::SM_pt);
445  if (cmdline.present("-sm-mt")) plugin->set_split_merge_scale(SISPlug::SM_mt);
446  if (cmdline.present("-sm-Et")) plugin->set_split_merge_scale(SISPlug::SM_Et);
447  if (cmdline.present("-sm-pttilde")) plugin->set_split_merge_scale(SISPlug::SM_pttilde);
448  // cause it to use the jet-definition's own recombiner
449  plugin->set_use_jet_def_recombiner(true);
450  jet_defs.push_back( fj::JetDefinition(plugin));
451  }
452  if (all_algs || cmdline.present("-sisconespheri")) {
453  double sisEmin = cmdline.value("-sisEmin",0.0);
454  fj::SISConeSphericalPlugin * plugin =
455  new fj::SISConeSphericalPlugin(ktR, overlap_threshold,npass,sisEmin);
456  if (cmdline.present("-ghost-sep")) {
457  plugin->set_ghost_separation_scale(cmdline.value<double>("-ghost-sep"));
458  }
459  jet_defs.push_back( fj::JetDefinition(plugin));
460  }
461 #else // FASTJET_ENABLE_PLUGIN_SISCONE
462  is_unavailable("siscone");
463 #endif // FASTJET_ENABLE_PLUGIN_SISCONE
464  }
465  if (all_algs || cmdline.present("-d0runiicone")) {
466 #ifdef FASTJET_ENABLE_PLUGIN_D0RUNIICONE
467  double min_jet_Et = 6.0; // was 8 GeV in earlier work
468  jet_defs.push_back( fj::JetDefinition(new fj::D0RunIIConePlugin(ktR,min_jet_Et)));
469 #else // FASTJET_ENABLE_PLUGIN_D0RUNIICONE
470  is_unavailable("D0RunIICone");
471 #endif // FASTJET_ENABLE_PLUGIN_D0RUNIICONE
472  }
473  if (all_algs || cmdline.present("-trackjet")) {
474 #ifdef FASTJET_ENABLE_PLUGIN_TRACKJET
475  jet_defs.push_back( fj::JetDefinition(new fj::TrackJetPlugin(ktR)));
476 #else // FASTJET_ENABLE_PLUGIN_TRACKJET
477  is_unavailable("TrackJet");
478 #endif // FASTJET_ENABLE_PLUGIN_TRACKJET
479  }
480  if (all_algs || cmdline.present("-atlascone")) {
481 #ifdef FASTJET_ENABLE_PLUGIN_ATLASCONE
482  jet_defs.push_back( fj::JetDefinition(new fj::ATLASConePlugin(ktR)));
483 #else // FASTJET_ENABLE_PLUGIN_ATLASCONE
484  is_unavailable("ATLASCone");
485 #endif // FASTJET_ENABLE_PLUGIN_ATLASCONE
486  }
487  if (all_algs || cmdline.present("-eecambridge")) {
488 #ifdef FASTJET_ENABLE_PLUGIN_EECAMBRIDGE
489  jet_defs.push_back( fj::JetDefinition(new fj::EECambridgePlugin(ycut)));
490 #else // FASTJET_ENABLE_PLUGIN_EECAMBRIDGE
491  is_unavailable("EECambridge");
492 #endif // FASTJET_ENABLE_PLUGIN_EECAMBRIDGE
493  }
494  if (all_algs || cmdline.present("-jade")) {
495 #ifdef FASTJET_ENABLE_PLUGIN_JADE
496  jet_defs.push_back( fj::JetDefinition(new fj::JadePlugin()));
497 #else // FASTJET_ENABLE_PLUGIN_JADE
498  is_unavailable("Jade");
499 #endif // FASTJET_ENABLE_PLUGIN_JADE
500  }
501  if (all_algs || cmdline.present("-cmsiterativecone")) {
502 #ifdef FASTJET_ENABLE_PLUGIN_CMSITERATIVECONE
503  jet_defs.push_back( fj::JetDefinition(new fj::CMSIterativeConePlugin(ktR,seed_threshold)));
504 #else // FASTJET_ENABLE_PLUGIN_CMSITERATIVECONE
505  is_unavailable("CMSIterativeCone");
506 #endif // FASTJET_ENABLE_PLUGIN_CMSITERATIVECONE
507  }
508  if (all_algs || cmdline.present("-d0runipre96cone")) {
509 #ifdef FASTJET_ENABLE_PLUGIN_D0RUNICONE
510  jet_defs.push_back( fj::JetDefinition(new fj::D0RunIpre96ConePlugin(ktR, seed_threshold, overlap_threshold)));
511 #else // FASTJET_ENABLE_PLUGIN_D0RUNICONE
512  is_unavailable("D0RunICone");
513 #endif // FASTJET_ENABLE_PLUGIN_D0RUNICONE
514  }
515  if (all_algs || cmdline.present("-d0runicone")) {
516 #ifdef FASTJET_ENABLE_PLUGIN_D0RUNICONE
517  jet_defs.push_back( fj::JetDefinition(new fj::D0RunIConePlugin(ktR, seed_threshold, overlap_threshold)));
518 #else // FASTJET_ENABLE_PLUGIN_D0RUNICONE
519  is_unavailable("D0RunICone");
520 #endif // FASTJET_ENABLE_PLUGIN_D0RUNICONE
521  }
522  if (all_algs || cmdline.present("-gridjet")) {
523 #ifdef FASTJET_ENABLE_PLUGIN_GRIDJET
524  // we want a grid_ymax of 5.0, but when using R=0.4 (i.e. grid
525  // spacing of 0.8), this leads to 12.5 grid cells; depending on
526  // whether this is 12.499999999999 or 12.5000000....1 this gets
527  // converted either to 12 or 13, making the results sensitive to
528  // rounding errors.
529  //
530  // Instead we therefore take 4.9999999999, which avoids this problem.
531  double grid_ymax = 4.9999999999;
532  jet_defs.push_back( fj::JetDefinition(new fj::GridJetPlugin(grid_ymax, ktR*2.0)));
533 #else // FASTJET_ENABLE_PLUGIN_GRIDJET
534  is_unavailable("GridJet");
535 #endif // FASTJET_ENABLE_PLUGIN_GRIDJET
536 // end of checking if one asks to run a plugin (don't delete this line)
537  }
538  if (all_algs ||
539  cmdline.present("-kt") ||
540  (jet_defs.size() == 0 && !found_unavailable)) {
541  jet_defs.push_back( fj::JetDefinition(fj::kt_algorithm, ktR, strategy));
542  }
543 
544  string filename = cmdline.value<string>("-file", "");
545 
546 
547  if (!cmdline.all_options_used()) {cerr <<
548  "Error: some options were not recognized"<<endl;
549  exit(-1);}
550 
551  for (unsigned idef = 0; idef < jet_defs.size(); idef++) {
552  fj::JetDefinition & jet_def = jet_defs[idef];
553  istream * istr;
554  if (filename == "") istr = &cin;
555  else istr = new ifstream(filename.c_str());
556 
557  for (int iev = 0; iev < nev; iev++) {
558  vector<fj::PseudoJet> jets;
559  vector<fj::PseudoJet> particles;
560  string line;
561  int ndone = 0;
562  while (getline(*istr, line)) {
563  //cout << line<<endl;
564  istringstream linestream(line);
565  if (line == "#END") {
566  ndone += 1;
567  if (ndone == combine) {break;}
568  }
569  if (line.substr(0,1) == "#") {continue;}
570  valarray<double> fourvec(4);
571  if (hydjet) {
572  // special reading from hydjet.txt event record (though actually
573  // this is supposed to be a standard pythia event record, so
574  // being able to read from it is perhaps not so bad an idea...)
575  int ii, istat,id,m1,m2,d1,d2;
576  double mass;
577  linestream >> ii>> istat >> id >> m1 >> m2 >> d1 >> d2
578  >> fourvec[0] >> fourvec[1] >> fourvec[2] >> mass;
579  // current file contains mass of particle as 4th entry
580  if (istat == 1) {
581  fourvec[3] = sqrt(+pow2(fourvec[0])+pow2(fourvec[1])
582  +pow2(fourvec[2])+pow2(mass));
583  }
584  } else {
585  if (massless) {
586  linestream >> fourvec[0] >> fourvec[1] >> fourvec[2];
587  fourvec[3] = sqrt(pow2(fourvec[0])+pow2(fourvec[1])+pow2(fourvec[2]));}
588  else {
589  linestream >> fourvec[0] >> fourvec[1] >> fourvec[2] >> fourvec[3];
590  }
591  }
592  fj::PseudoJet psjet(fourvec);
593  if (abs(psjet.rap() < etamax)) {particles.push_back(psjet);}
594  }
595 
596  // add a fake underlying event which is very soft, uniformly distributed
597  // in eta,phi so as to allow one to reconstruct the area that is associated
598  // with each jet.
599  if (add_dense_coverage) {
600  fj::GhostedAreaSpec ghosted_area_spec(ghost_maxrap);
601  //fj::GhostedAreaSpec ghosted_area_spec(-2.0,4.0); // asymmetric range
602  // for plots, reduce the scatter default of 1, to avoid "holes"
603  // in the subsequent calorimeter view
604  ghosted_area_spec.set_grid_scatter(0.5);
605  ghosted_area_spec.add_ghosts(particles);
606  //----- old code ------------------
607  // srand(2);
608  // int nphi = 60;
609  // int neta = 100;
610  // double kt = 1e-1;
611  // for (int iphi = 0; iphi<nphi; iphi++) {
612  // for (int ieta = -neta; ieta<neta+1; ieta++) {
613  // double phi = (iphi+0.5) * (fj::twopi/nphi) + rand()*0.001/RAND_MAX;
614  // double eta = ieta * (10.0/neta) + rand()*0.001/RAND_MAX;
615  // kt = 1e-20*(1+rand()*0.1/RAND_MAX);
616  // double pminus = kt*exp(-eta);
617  // double pplus = kt*exp(+eta);
618  // double px = kt*sin(phi);
619  // double py = kt*cos(phi);
620  // //cout << kt<<" "<<eta<<" "<<phi<<"\n";
621  // fj::PseudoJet mom(px,py,0.5*(pplus-pminus),0.5*(pplus+pminus));
622  // particles.push_back(mom);
623  // }
624  // }
625  }
626 
627  // select the particles that pass the selection cut
628  particles = particles_sel(particles);
629 
630  for (int irepeat = 0; irepeat < repeat ; irepeat++) {
631  int nparticles = particles.size();
632  try {
633  auto_ptr<fj::ClusterSequence> clust_seq;
634  if (do_areas) {
635  clust_seq.reset(new fj::ClusterSequenceArea(particles,jet_def,area_def));
636  } else {
637  clust_seq.reset(new fj::ClusterSequence(particles,jet_def,write));
638  }
639 
640  // repetitive output
641  if (repeatinclkt >= 0.0) {
642  vector<fj::PseudoJet> jets_local = sorted_by_pt(clust_seq->inclusive_jets(repeatinclkt));
643  }
644 
645  if (irepeat != 0) {continue;}
646  cout << "iev "<<iev<< ": number of particles = "<< nparticles << endl;
647  cout << "strategy used = "<< clust_seq->strategy_string()<< endl;
648  if (iev == 0) cout << "Jet Definition: " << jet_def.description() << " (" << fj::fastjet_version_string() << ")" << endl;
649  if (do_areas && iev == 0) cout << "Area definition: " << area_def.description() << endl;
650 
651  // now provide some nice output...
652  if (inclkt >= 0.0) {
653  vector<fj::PseudoJet> jets_local = sorted_by_pt(clust_seq->inclusive_jets(inclkt));
654  print_jets(jets_local, show_constituents);
655 
656  }
657 
658  if (excln > 0) {
659  cout << "Printing "<<excln<<" exclusive jets\n";
660  print_jets(clust_seq->exclusive_jets(excln), show_constituents);
661  }
662 
663  if (excld > 0.0) {
664  cout << "Printing exclusive jets for d = "<<excld<<"\n";
665  print_jets(clust_seq->exclusive_jets(excld), show_constituents);
666  }
667 
668  if (excly > 0.0) {
669  cout << "Printing exclusive jets for ycut = "<<excly<<"\n";
670  print_jets(clust_seq->exclusive_jets_ycut(excly), show_constituents);
671  }
672 
673  if (get_all_dij) {
674  for (int i = nparticles-1; i >= 0; i--) {
675  printf("d for n = %4d -> %4d is %14.5e\n", i+1, i, clust_seq->exclusive_dmerge(i));
676  }
677  }
678  if (get_all_yij) {
679  for (int i = nparticles-1; i >= 0; i--) {
680  printf("y for n = %4d -> %4d is %14.5e\n", i+1, i, clust_seq->exclusive_ymerge(i));
681  }
682  }
683 
684  // have the option of printing out the subjets (at scale dcut) of
685  // each inclusive jet
686  if (subdcut >= 0.0) {
687  print_jets_and_sub(clust_seq->inclusive_jets(), subdcut);
688  }
689 
690  // useful for testing that recombination sequences are unique
691  if (unique_write) {
692  vector<int> unique_history = clust_seq->unique_history_order();
693  // construct the inverse of the above mapping
694  vector<int> inv_unique_history(clust_seq->history().size());
695  for (unsigned int i = 0; i < unique_history.size(); i++) {
696  inv_unique_history[unique_history[i]] = i;}
697 
698  for (unsigned int i = 0; i < unique_history.size(); i++) {
699  fj::ClusterSequence::history_element el =
700  clust_seq->history()[unique_history[i]];
701  int uhp1 = el.parent1>=0 ? inv_unique_history[el.parent1] : el.parent1;
702  int uhp2 = el.parent2>=0 ? inv_unique_history[el.parent2] : el.parent2;
703  printf("%7d u %15.8e %7d u %7d u\n",i,el.dij,uhp1, uhp2);
704  }
705  }
706 
707 
708 #ifdef FASTJET_ENABLE_PLUGIN_SISCONE
709  // provide some complementary information for SISCone
710  if (show_cones) {
711  const fj::SISConeExtras * extras =
712  dynamic_cast<const fj::SISConeExtras *>(clust_seq->extras());
713  cout << "most ambiguous split (difference in squared dist) = "
714  << extras->most_ambiguous_split() << endl;
715  vector<fastjet::PseudoJet> stable_cones(extras->stable_cones());
716  stable_cones = sorted_by_rapidity(stable_cones);
717  for (unsigned int i = 0; i < stable_cones.size(); i++) {
718  //if (stable_cones[i].phi() < 5.0 && stable_cones[i].phi() > 4.0) {
719  printf("%5u %15.8f %15.8f %15.8f\n",
720  i,stable_cones[i].rap(),stable_cones[i].phi(),
721  stable_cones[i].perp() );
722  //}
723  }
724 
725  // also show passes for jets
726  vector<fj::PseudoJet> sisjets = clust_seq->inclusive_jets();
727  printf("\n%15s %15s %15s %12s %8s %8s\n","rap","phi","pt","user-index","pass","nconst");
728  for (unsigned i = 0; i < sisjets.size(); i++) {
729  printf("%15.8f %15.8f %15.8f %12d %8d %8u\n",
730  sisjets[i].rap(), sisjets[i].phi(), sisjets[i].perp(),
731  sisjets[i].user_index(), extras->pass(sisjets[i]),
732  (unsigned int) clust_seq->constituents(sisjets[i]).size()
733  );
734 
735  }
736  }
737 #endif // FASTJET_ENABLE_PLUGIN_SISCONE
738 
739  if (do_bkgd) {
740  double rho, sigma, mean_area, empty_area, n_empty_jets;
742  dynamic_cast<fj::ClusterSequenceAreaBase *>(clust_seq.get());
743  if (do_bkgd_csab) {
744  csab->get_median_rho_and_sigma(bkgd_range, true, rho, sigma, mean_area);
745  empty_area = csab->empty_area(bkgd_range);
746  n_empty_jets = csab->n_empty_jets(bkgd_range);
747  } else if (do_bkgd_jetmedian) {
748  fj::JetMedianBackgroundEstimator bge(bkgd_range);
749  bge.set_provide_fj2_sigma(do_bkgd_fj2);
750  bge.set_cluster_sequence(*csab);
751  rho = bge.rho();
752  sigma = bge.sigma();
753  mean_area = bge.mean_area();
754  empty_area = bge.empty_area();
755  n_empty_jets = bge.n_empty_jets();
756  } else {
757  assert(do_bkgd_gridmedian);
758  double rapmin, rapmax;
759  bkgd_range.get_rapidity_extent(rapmin, rapmax);
760  fj::GridMedianBackgroundEstimator bge(rapmax, 2*ktR);
761  bge.set_particles(particles);
762  rho = bge.rho();
763  sigma = bge.sigma();
764  mean_area = bge.mean_area();
765  empty_area = 0;
766  n_empty_jets = 0;
767  }
768  cout << " rho = " << rho
769  << ", sigma = " << sigma
770  << ", mean_area = " << mean_area
771  << ", empty_area = " << empty_area
772  << ", n_empty_jets = " << n_empty_jets
773  << endl;
774  }
775  } // try
776  catch (fastjet::Error fjerr) {
777  cout << "Caught fastjet error, exiting gracefully" << endl;
778  exit(0);
779  }
780 
781  } // irepeat
782  } // iev
783  // if we've instantiated a plugin, delete it
784  if (jet_def.strategy()==fj::plugin_strategy){
785  delete jet_def.plugin();
786  }
787  // close any file that we've opened
788  if (istr != &cin) delete istr;
789  } // jet_defs
790 
791 }
792 
793 
794 
795 
796 //------ HELPER ROUTINES -----------------------------------------------
797 /// print a single jet
798 void print_jet (const fj::PseudoJet & jet) {
799  unsigned int n_constituents = jet.constituents().size();
800  printf("%15.8f %15.8f %15.8f %8u\n",
801  jet.rap(), jet.phi(), jet.perp(), n_constituents);
802 }
803 
804 
805 //----------------------------------------------------------------------
806 void print_jets(const vector<fj::PseudoJet> & jets_in, bool show_constituents) {
807  vector<fj::PseudoJet> jets;
808  if (ee_print) {
809  jets = sorted_by_E(jets_in);
810  for (unsigned int j = 0; j < jets.size(); j++) {
811  printf("%5u %15.8f %15.8f %15.8f %15.8f\n",
812  j,jets[j].px(),jets[j].py(),jets[j].pz(),jets[j].E());
813  if (show_constituents) {
814  vector<fj::PseudoJet> const_jets = jets[j].constituents();
815  for (unsigned int k = 0; k < const_jets.size(); k++) {
816  printf(" jet%03u %15.8f %15.8f %15.8f %15.8f\n",j,const_jets[k].px(),
817  const_jets[k].py(),const_jets[k].pz(),const_jets[k].E());
818  }
819  cout << "\n\n";
820  }
821 
822  }
823  } else {
824  jets = sorted_by_pt(jets_in);
825  for (unsigned int j = 0; j < jets.size(); j++) {
826  printf("%5u %15.8f %15.8f %15.8f",
827  j,jets[j].rap(),jets[j].phi(),jets[j].perp());
828  // also print out the scalar area and the perp component of the
829  // 4-vector (just enough to check a reasonable 4-vector?)
830  if (do_areas) printf(" %15.8f %15.8f", jets[j].area(),
831  jets[j].area_4vector().perp());
832  cout << "\n";
833 
834  if (show_constituents) {
835  vector<fj::PseudoJet> const_jets = jets[j].constituents();
836  for (unsigned int k = 0; k < const_jets.size(); k++) {
837  printf(" jet%03u %15.8f %15.8f %15.8f %5d\n",j,const_jets[k].rap(),
838  const_jets[k].phi(),sqrt(const_jets[k].kt2()), const_jets[k].cluster_hist_index());
839  }
840  cout << "\n\n";
841  }
842  }
843  }
844 
845  if (rootfile != "") {
846  ofstream ostr(rootfile.c_str());
847  ostr << "# " << cmdline_p->command_line() << endl;
848  ostr << "# output for root" << endl;
849  assert(jets.size() > 0);
850  jets[0].validated_cs()->print_jets_for_root(jets,ostr);
851  }
852 
853 }
854 
855 
856 //----- SUBJETS --------------------------------------------------------
857 /// a function that pretty prints a list of jets and the subjets for each
858 /// one
859 void print_jets_and_sub (const vector<fj::PseudoJet> & jets, double dcut) {
860 
861  // sort jets into increasing pt
862  vector<fj::PseudoJet> sorted_jets = sorted_by_pt(jets);
863 
864  // label the columns
865  printf("Printing jets and their subjets with subdcut = %10.5f\n",dcut);
866  printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity",
867  "phi", "pt", "n constituents");
868 
869  // have various kinds of subjet finding, to test consistency among them
870  enum SubType {internal, newclust_dcut, newclust_R};
871  SubType subtype = internal;
872  //SubType subtype = newclust_dcut;
873  //SubType subtype = newclust_R;
874 
875  // print out the details for each jet
876  //for (unsigned int i = 0; i < sorted_jets.size(); i++) {
877  for (vector<fj::PseudoJet>::const_iterator jet = sorted_jets.begin();
878  jet != sorted_jets.end(); jet++) {
879  const fj::JetDefinition & jet_def = jet->validated_cs()->jet_def();
880 
881  // if jet pt^2 < dcut with kt alg, then some methods of
882  // getting subjets will return nothing -- so skip the jet
883  if (jet_def.jet_algorithm() == fj::kt_algorithm
884  && jet->perp2() < dcut) continue;
885 
886 
887  printf("%5u ",(unsigned int) (jet - sorted_jets.begin()));
888  print_jet(*jet);
889  vector<fj::PseudoJet> subjets;
890  fj::ClusterSequence * cspoint;
891  if (subtype == internal) {
892  cspoint = 0;
893  subjets = jet->exclusive_subjets(dcut);
894  double ddnp1 = jet->exclusive_subdmerge_max(subjets.size());
895  double ddn = jet->exclusive_subdmerge_max(subjets.size()-1);
896  cout << " for " << ddnp1 << " < d < " << ddn << " one has " << endl;
897  } else if (subtype == newclust_dcut) {
898  cspoint = new fj::ClusterSequence(jet->constituents(), jet_def);
899  subjets = cspoint->exclusive_jets(dcut);
900  } else if (subtype == newclust_R) {
901  assert(jet_def.jet_algorithm() == fj::cambridge_algorithm);
902  fj::JetDefinition subjd(jet_def.jet_algorithm(), jet_def.R()*sqrt(dcut));
903  cspoint = new fj::ClusterSequence(jet->constituents(), subjd);
904  subjets = cspoint->inclusive_jets();
905  } else {
906  cerr << "unrecognized subtype for subjet finding" << endl;
907  exit(-1);
908  }
909 
910  subjets = sorted_by_pt(subjets);
911  for (unsigned int j = 0; j < subjets.size(); j++) {
912  printf(" -sub-%02u ",j);
913  print_jet(subjets[j]);
914  }
915 
916  if (cspoint != 0) delete cspoint;
917 
918  //fj::ClusterSequence subseq(clust_seq->constituents(sorted_jets[i]),
919  // fj::JetDefinition(fj::cambridge_algorithm, 0.4));
920  //vector<fj::PseudoJet> subjets = sorted_by_pt(subseq.inclusive_jets());
921  //for (unsigned int j = 0; j < subjets.size(); j++) {
922  // printf(" -sub-%02u ",j);
923  // print_jet(subseq, subjets[j]);
924  //}
925  }
926 
927 }
928