FastJet  3.4.0
fastjet_timing_plugins.cc
1 //STARTHEADER
2 // $Id$
3 //
4 // Copyright (c) 2005-2018, 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 /// -compare-strategy STRAT
126 /// compares the output of the default strategy (possibly as specified
127 /// with -strategy) with that from STRAT. Currently compares the history.
128 
129 /// Algorithms
130 /// ----------
131 /// -all-algs runs all algorithms
132 ///
133 /// -kt switch to the longitudinally invariant kt algorithm
134 /// Note: this is the default one.
135 ///
136 /// -cam switch to the inclusive Cambridge/Aachen algorithm --
137 /// note that the option -excld dcut provides a clustering
138 /// up to the dcut which is the minimum squared
139 /// distance between any pair of jets.
140 ///
141 /// -antikt switch to the anti-kt clustering algorithm
142 ///
143 /// -genkt switch to the genkt algorithm
144 /// you can provide the parameter of the alg as an argument to
145 /// -genkt (1 by default)
146 ///
147 /// -eekt switch to the e+e- kt algorithm
148 ///
149 /// -eegenkt switch to the genkt algorithm
150 /// you can provide the parameter of the alg as an argument to
151 /// -ee_genkt (1 by default)
152 ///
153 /// plugins (don't delete this line)
154 ///
155 /// -pxcone switch to the PxCone jet algorithm
156 ///
157 /// -siscone switch to the SISCone jet algorithm (seedless cones)
158 /// -sisconespheri switch to the Spherical SISCone jet algorithm (seedless cones)
159 ///
160 /// -midpoint switch to CDF's midpoint code
161 /// -jetclu switch to CDF's jetclu code
162 ///
163 /// -d0runipre96cone switch to the D0RunIpre96Cone plugin
164 /// -d0runicone switch to the D0RunICone plugin
165 ///
166 /// -d0runiicone switch to D0's run II midpoint cone
167 ///
168 /// -trackjet switch to the TrackJet plugin
169 ///
170 /// -atlascone switch to the ATLASCone plugin
171 ///
172 /// -eecambridge switch to the EECambridge plugin
173 ///
174 /// -jade switch to the Jade plugin
175 ///
176 /// -cmsiterativecone switch to the CMSIterativeCone plugin
177 ///
178 /// -gridjet switch to the GridJet plugin
179 ///
180 /// end of plugins (don't delete this line)
181 ///
182 ///
183 /// Options for running algs
184 /// ------------------------
185 ///
186 /// -r sets the radius of the jet algorithm (default = 1.0)
187 ///
188 /// -overlap | -f sets the overlap fraction in cone algs with split-merge
189 ///
190 /// -seed sets the seed threshold
191 ///
192 /// -strategy N indicate stratgey from the enum fastjet::Strategy (see
193 /// fastjet/JetDefinition.hh).
194 ///
195 
196 #ifndef __FJCORE__
197 #include "fastjet/ClusterSequenceArea.hh"
198 #include "fastjet/tools/JetMedianBackgroundEstimator.hh"
199 #include "fastjet/tools/GridMedianBackgroundEstimator.hh"
200 #include "fastjet/tools/Subtractor.hh"
201 #include "fastjet/Selector.hh"
202 #else
203 #include "fjcore.hh"
204 #endif
205 #include<iostream>
206 #include<iomanip>
207 #include<sstream>
208 #include<fstream>
209 #include<valarray>
210 #include<vector>
211 #include <cstdlib>
212 //#include<cstddef> // for size_t
213 #include "CmdLine.hh"
214 
215 #ifndef __FJCORE__
216 // get info on how fastjet was configured
217 #include "fastjet/config.h"
218 #endif
219 
220 // include the installed plugins (don't delete this line)
221 #ifdef FASTJET_ENABLE_PLUGIN_SISCONE
222 #include "fastjet/SISConePlugin.hh"
223 #include "fastjet/SISConeSphericalPlugin.hh"
224 #endif
225 #ifdef FASTJET_ENABLE_PLUGIN_CDFCONES
226 #include "fastjet/CDFMidPointPlugin.hh"
227 #include "fastjet/CDFJetCluPlugin.hh"
228 #endif
229 #ifdef FASTJET_ENABLE_PLUGIN_PXCONE
230 #include "fastjet/PxConePlugin.hh"
231 #endif
232 #ifdef FASTJET_ENABLE_PLUGIN_D0RUNIICONE
233 #include "fastjet/D0RunIIConePlugin.hh"
234 #endif
235 #ifdef FASTJET_ENABLE_PLUGIN_TRACKJET
236 #include "fastjet/TrackJetPlugin.hh"
237 #endif
238 #ifdef FASTJET_ENABLE_PLUGIN_ATLASCONE
239 #include "fastjet/ATLASConePlugin.hh"
240 #endif
241 #ifdef FASTJET_ENABLE_PLUGIN_EECAMBRIDGE
242 #include "fastjet/EECambridgePlugin.hh"
243 #endif
244 #ifdef FASTJET_ENABLE_PLUGIN_JADE
245 #include "fastjet/JadePlugin.hh"
246 #endif
247 #ifdef FASTJET_ENABLE_PLUGIN_CMSITERATIVECONE
248 #include "fastjet/CMSIterativeConePlugin.hh"
249 #endif
250 #ifdef FASTJET_ENABLE_PLUGIN_D0RUNICONE
251 #include "fastjet/D0RunIpre96ConePlugin.hh"
252 #include "fastjet/D0RunIConePlugin.hh"
253 #endif
254 #ifdef FASTJET_ENABLE_PLUGIN_GRIDJET
255 #include "fastjet/GridJetPlugin.hh"
256 #endif
257 // end of installed plugins inclusion (don't delete this line)
258 
259 using namespace std;
260 
261 // to avoid excessive typing, use the fastjet namespace
262 #ifndef __FJCORE__
263 using namespace fastjet;
264 #else
265 using namespace fjcore;
266 #endif
267 
268 inline double pow2(const double x) {return x*x;}
269 
270 // pretty print the jets and their subjets
271 void print_jets_and_sub (const vector<PseudoJet> & jets, double dcut);
272 
273 #ifndef __FJCORE__
274 void print_jets_bkgd(const vector<PseudoJet> &jets,
275  const vector<PseudoJet> &subtracted_jets,
276  BackgroundEstimatorBase * bge_ptr,
277  bool do_subtractor);
278 #endif // __FJCORE__
279 
280 // have various kinds of subjet finding, to test consistency among them
281 //
282 // this is needed in print_jets_and_sub and declaring it in the
283 // function scope results in errors with older intel compilers (due to
284 // the overloaded == operator in PseudoJet which results in the "a
285 // template argument may not reference a local type" error)
286 enum SubType {subtype_internal, subtype_newclust_dcut, subtype_newclust_R};
287 
288 void do_compare_strategy(int iev,
289  const vector<PseudoJet> & particles,
290  const JetDefinition & jet_def,
291  const ClusterSequence & cs,
292  int compare_strategy);
293 
294 
295 string rootfile;
296 CmdLine * cmdline_p;
297 
298 bool do_areas;
299 
300 /// sort and pretty print jets, with exact behaviour depending on
301 /// whether ee_print is true or not
302 bool ee_print = false;
303 void print_jets(const vector<PseudoJet> & jets, bool show_const = false);
304 
305 bool found_unavailable = false;
306 void is_unavailable(const string & algname) {
307  cerr << algname << " requested, but not available for this compilation" << endl;
308  found_unavailable = true;
309  //exit(0);
310 }
311 
312 
313 /// a program to test and time a range of algorithms as implemented or
314 /// wrapped in fastjet
315 int main (int argc, char ** argv) {
316 
317  CmdLine cmdline(argc,argv);
318  cout << "# " << cmdline.command_line() << endl;
319  ClusterSequence::print_banner();
320 
321  cmdline_p = &cmdline;
322  // allow the use to specify the Strategy either through the
323  // -clever or the -strategy options (both will take numerical
324  // values); the latter will override the former.
325  Strategy strategy = Strategy(cmdline.int_val("-strategy",
326  cmdline.int_val("-clever", Best)));
327  int repeat = cmdline.int_val("-repeat",1);
328  int combine = cmdline.int_val("-combine",1);
329  bool write = cmdline.present("-write");
330  bool unique_write = cmdline.present("-unique_write");
331  bool hydjet = cmdline.present("-hydjet");
332  double ktR = cmdline.double_val("-r",1.0);
333  ktR = cmdline.double_val("-R",ktR); // allow -r and -R
334  double inclkt = cmdline.double_val("-incl",-1.0);
335  double repeatinclkt = cmdline.double_val("-repeat-incl",-1.0);
336  int excln = cmdline.int_val ("-excln",-1);
337  double excld = cmdline.double_val("-excld",-1.0);
338  double excly = cmdline.double_val("-excly",-1.0);
339  ee_print = cmdline.present("-ee-print");
340  bool get_all_dij = cmdline.present("-get-all-dij");
341  bool get_all_yij = cmdline.present("-get-all-yij");
342  double subdcut = cmdline.double_val("-subdcut",-1.0);
343  double rapmax = cmdline.double_val("-rapmax",1.0e305);
344  if (cmdline.present("-etamax")) {
345  cerr << "WARNING: -etamax options actually sets maximum rapidity (and overrides -rapmax)\n";
346  rapmax = cmdline.double_val("-etamax");
347  }
348  bool show_constituents = cmdline.present("-const");
349  bool massless = cmdline.present("-massless");
350  int nev = cmdline.int_val("-nev",1);
351  int skip = cmdline.int_val("-skip",0);
352  bool add_dense_coverage = cmdline.present("-dense");
353  double ghost_maxrap = cmdline.value("-ghost-maxrap",5.0);
354  bool all_algs = cmdline.present("-all-algs");
355  // have the option of comparing the clustering results to those
356  // obtained with a different clustering strategy; the misuse the
357  // "plugin_strategy" to indicate that no comparison is needed.
358  // Does not currently support areas
359  int compare_strategy = cmdline.value<int>("-compare-strategy", plugin_strategy);
360 
361 
362  Selector particles_sel = (cmdline.present("-nhardest"))
363  ? SelectorNHardest(cmdline.value<unsigned int>("-nhardest"))
364  : SelectorIdentity();
365 
366  do_areas = cmdline.present("-area");
367 #ifndef __FJCORE__
368  AreaDefinition area_def;
369  if (do_areas) {
370  assert(!write); // it's incompatible
371  GhostedAreaSpec ghost_spec(ghost_maxrap,
372  cmdline.value("-area:repeat", 1),
373  cmdline.value("-ghost-area", 0.01));
374  if (cmdline.present("-area:fj2")) ghost_spec.set_fj2_placement(true);
375  if (cmdline.present("-area:explicit")) {
376  area_def = AreaDefinition(active_area_explicit_ghosts, ghost_spec);
377  } else if (cmdline.present("-area:passive")) {
378  area_def = AreaDefinition(passive_area, ghost_spec);
379  } else if (cmdline.present("-area:voronoi")) {
380  double Rfact = cmdline.value<double>("-area:voronoi");
381  area_def = AreaDefinition(voronoi_area,
382  VoronoiAreaSpec(Rfact));
383  } else {
384  cmdline.present("-area:active"); // allow, but do not require, arg
385  area_def = AreaDefinition(active_area, ghost_spec);
386  }
387  }
388 #else
389  do_areas=false;
390 #endif
391 
392  bool do_bkgd = cmdline.present("-bkgd"); // background estimation
393 #ifndef __FJCORE__
394  bool do_bkgd_csab = false, do_bkgd_jetmedian = false, do_bkgd_fj2 = false;
395  bool do_bkgd_gridmedian = false;
396  bool do_bkgd_localrange = false;
397  bool do_subtractor = false;
398  double bkgd_alt_ktR = -1.0;
399  BackgroundRescalingYPolynomial * bkgd_rescaling = 0;
400  Selector bkgd_range;
401  if (do_bkgd) {
402  bkgd_range = SelectorAbsRapMax(ghost_maxrap - ktR);
403  if (cmdline.present("-bkgd:csab")) {do_bkgd_csab = true;}
404  else if (cmdline.present("-bkgd:jetmedian")) {do_bkgd_jetmedian = true;
405  do_bkgd_fj2 = cmdline.present("-bkgd:fj2");
406  do_bkgd_localrange = cmdline.present("-bkgd:localrange");
407  bkgd_alt_ktR = cmdline.value("-bkgd:alt-ktR", bkgd_alt_ktR);
408  if (do_bkgd_localrange) bkgd_range = SelectorStrip(1.5);
409  } else if (cmdline.present("-bkgd:gridmedian")) {
410  do_bkgd_gridmedian = true;
411  } else {
412  throw Error("with the -bkgd option, some particular background must be specified (csab or jetmedian)");
413  }
414  if (cmdline.present("-bkgd:rescaling")) {
415  bkgd_rescaling = new BackgroundRescalingYPolynomial(1.157,0,-0.0266,0,0.000048);
416  }
417  assert(do_areas || do_bkgd_gridmedian);
418  do_subtractor = cmdline.present("-subtractor");
419  if (do_subtractor) assert(do_areas);
420  }
421 #else
422  do_bkgd = false;
423 #endif
424 
425  bool show_cones = cmdline.present("-cones"); // only works for siscone
426 
427  // for cone algorithms
428  // allow -f and -overlap
429  double overlap_threshold = cmdline.double_val("-overlap",0.5);
430  overlap_threshold = cmdline.double_val("-f",overlap_threshold);
431  double seed_threshold = cmdline.double_val("-seed",1.0);
432 
433 #ifdef __FJCORE__
434  show_cones = false;
435 #endif
436 
437  // for ee algorithms, allow to specify ycut
438  double ycut = cmdline.double_val("-ycut",0.08);
439 
440  // for printing jets to a file for reading by root
441  rootfile = cmdline.value<string>("-root","");
442 
443  // out default scheme is the E_scheme
444  RecombinationScheme scheme = E_scheme;
445 
446  // The following option causes the Cambridge algo to be used.
447  // Note that currently the only output that works sensibly here is
448  // "-incl 0"
449  vector<JetDefinition> jet_defs;
450  if (all_algs || cmdline.present("-cam") || cmdline.present("-CA")) {
451  jet_defs.push_back( JetDefinition(cambridge_algorithm, ktR, scheme, strategy));
452  }
453  if (all_algs || cmdline.present("-antikt")) {
454  jet_defs.push_back( JetDefinition(antikt_algorithm, ktR, scheme, strategy));
455  }
456  if (all_algs || cmdline.present("-genkt")) {
457  double p;
458  if (cmdline.present("-genkt")) p = cmdline.value<double>("-genkt");
459  else p = -0.5;
460  jet_defs.push_back( JetDefinition(genkt_algorithm, ktR, p, scheme, strategy));
461  }
462  if (all_algs || cmdline.present("-eekt")) {
463  jet_defs.push_back( JetDefinition(ee_kt_algorithm));
464  }
465  if (all_algs || cmdline.present("-eegenkt")) {
466  double p;
467  if (cmdline.present("-eegenkt")) p = cmdline.value<double>("-eegenkt");
468  else p = -0.5;
469  jet_defs.push_back( JetDefinition(ee_genkt_algorithm, ktR, p, scheme, strategy));
470 
471 // checking if one asks to run a plugin (don't delete this line)
472  }
473  if (all_algs || cmdline.present("-midpoint")) {
474 #ifdef FASTJET_ENABLE_PLUGIN_CDFCONES
475  typedef CDFMidPointPlugin MPPlug; // for brevity
476  double cone_area_fraction = 1.0;
477  int max_pair_size = 2;
478  int max_iterations = 100;
479  MPPlug::SplitMergeScale sm_scale = MPPlug::SM_pt;
480  if (cmdline.present("-sm-pttilde")) sm_scale = MPPlug::SM_pttilde;
481  if (cmdline.present("-sm-pt")) sm_scale = MPPlug::SM_pt; // default
482  if (cmdline.present("-sm-mt")) sm_scale = MPPlug::SM_mt;
483  if (cmdline.present("-sm-Et")) sm_scale = MPPlug::SM_Et;
484  jet_defs.push_back( JetDefinition( new CDFMidPointPlugin (
485  seed_threshold, ktR,
486  cone_area_fraction, max_pair_size,
487  max_iterations, overlap_threshold,
488  sm_scale)));
489 #else // FASTJET_ENABLE_PLUGIN_CDFCONES
490  is_unavailable("MidPoint");
491 #endif // FASTJET_ENABLE_PLUGIN_CDFCONES
492  }
493  if (all_algs || cmdline.present("-pxcone")) {
494 #ifdef FASTJET_ENABLE_PLUGIN_PXCONE
495  double min_jet_energy = 5.0;
496  // mode: 1=e+e-, 2=pp
497  int mode = cmdline.value("-pxcone-mode", 2);
498  bool E_scheme_jets = false;
499  jet_defs.push_back( JetDefinition( new PxConePlugin (
500  ktR, min_jet_energy,
501  overlap_threshold, E_scheme_jets, mode)));
502 #else // FASTJET_ENABLE_PLUGIN_PXCONE
503  is_unavailable("PxCone");
504 #endif // FASTJET_ENABLE_PLUGIN_PXCONE
505  }
506  if (all_algs || cmdline.present("-jetclu")) {
507 #ifdef FASTJET_ENABLE_PLUGIN_CDFCONES
508  jet_defs.push_back( JetDefinition( new CDFJetCluPlugin (
509  ktR, overlap_threshold, seed_threshold)));
510 #else // FASTJET_ENABLE_PLUGIN_CDFCONES
511  is_unavailable("JetClu");
512 #endif // FASTJET_ENABLE_PLUGIN_CDFCONES
513  }
514  if (all_algs || cmdline.present("-siscone") || cmdline.present("-sisconespheri")) {
515 #ifdef FASTJET_ENABLE_PLUGIN_SISCONE
516  typedef SISConePlugin SISPlug; // for brevity
517  int npass = cmdline.value("-npass",0);
518  if (all_algs || cmdline.present("-siscone")) {
519  double sisptmin = cmdline.value("-sisptmin",0.0);
520  bool cache = cmdline.present("-cache");
521  SISPlug * plugin = new SISPlug (ktR, overlap_threshold,npass,sisptmin,cache);
522  if (cmdline.present("-sm-pt")) plugin->set_split_merge_scale(SISPlug::SM_pt);
523  if (cmdline.present("-sm-mt")) plugin->set_split_merge_scale(SISPlug::SM_mt);
524  if (cmdline.present("-sm-Et")) plugin->set_split_merge_scale(SISPlug::SM_Et);
525  if (cmdline.present("-sm-pttilde")) plugin->set_split_merge_scale(SISPlug::SM_pttilde);
526  // cause it to use the jet-definition's own recombiner
527  plugin->set_use_jet_def_recombiner(true);
528  jet_defs.push_back( JetDefinition(plugin));
529  }
530  if (all_algs || cmdline.present("-sisconespheri")) {
531  double sisEmin = cmdline.value("-sisEmin",0.0);
532  SISConeSphericalPlugin * plugin =
533  new SISConeSphericalPlugin(ktR, overlap_threshold,npass,sisEmin);
534  if (cmdline.present("-ghost-sep")) {
535  plugin->set_ghost_separation_scale(cmdline.value<double>("-ghost-sep"));
536  }
537  jet_defs.push_back( JetDefinition(plugin));
538  }
539 #else // FASTJET_ENABLE_PLUGIN_SISCONE
540  is_unavailable("SISCone");
541 #endif // FASTJET_ENABLE_PLUGIN_SISCONE
542  }
543  if (all_algs || cmdline.present("-d0runiicone")) {
544 #ifdef FASTJET_ENABLE_PLUGIN_D0RUNIICONE
545  double min_jet_Et = 6.0; // was 8 GeV in earlier work
546  jet_defs.push_back( JetDefinition(new D0RunIIConePlugin(ktR,min_jet_Et)));
547 #else // FASTJET_ENABLE_PLUGIN_D0RUNIICONE
548  is_unavailable("D0RunIICone");
549 #endif // FASTJET_ENABLE_PLUGIN_D0RUNIICONE
550  }
551  if (all_algs || cmdline.present("-trackjet")) {
552 #ifdef FASTJET_ENABLE_PLUGIN_TRACKJET
553  jet_defs.push_back( JetDefinition(new TrackJetPlugin(ktR)));
554 #else // FASTJET_ENABLE_PLUGIN_TRACKJET
555  is_unavailable("TrackJet");
556 #endif // FASTJET_ENABLE_PLUGIN_TRACKJET
557  }
558  if (all_algs || cmdline.present("-atlascone")) {
559 #ifdef FASTJET_ENABLE_PLUGIN_ATLASCONE
560  jet_defs.push_back( JetDefinition(new ATLASConePlugin(ktR)));
561 #else // FASTJET_ENABLE_PLUGIN_ATLASCONE
562  is_unavailable("ATLASCone");
563 #endif // FASTJET_ENABLE_PLUGIN_ATLASCONE
564  }
565  if (all_algs || cmdline.present("-eecambridge")) {
566 #ifdef FASTJET_ENABLE_PLUGIN_EECAMBRIDGE
567  jet_defs.push_back( JetDefinition(new EECambridgePlugin(ycut)));
568 #else // FASTJET_ENABLE_PLUGIN_EECAMBRIDGE
569  is_unavailable("EECambridge");
570 #endif // FASTJET_ENABLE_PLUGIN_EECAMBRIDGE
571  }
572  if (all_algs || cmdline.present("-jade")) {
573 #ifdef FASTJET_ENABLE_PLUGIN_JADE
574  JadePlugin::Strategy jade_strategy =
575  JadePlugin::Strategy(cmdline.value<int>("-jade-strategy",
576  JadePlugin::strategy_NNFJN2Plain));
577  jet_defs.push_back( JetDefinition(new JadePlugin(jade_strategy)));
578 #else // FASTJET_ENABLE_PLUGIN_JADE
579  is_unavailable("Jade");
580 #endif // FASTJET_ENABLE_PLUGIN_JADE
581  }
582  if (all_algs || cmdline.present("-cmsiterativecone")) {
583 #ifdef FASTJET_ENABLE_PLUGIN_CMSITERATIVECONE
584  jet_defs.push_back( JetDefinition(new CMSIterativeConePlugin(ktR,seed_threshold)));
585 #else // FASTJET_ENABLE_PLUGIN_CMSITERATIVECONE
586  is_unavailable("CMSIterativeCone");
587 #endif // FASTJET_ENABLE_PLUGIN_CMSITERATIVECONE
588  }
589  if (all_algs || cmdline.present("-d0runipre96cone")) {
590 #ifdef FASTJET_ENABLE_PLUGIN_D0RUNICONE
591  jet_defs.push_back( JetDefinition(new D0RunIpre96ConePlugin(ktR, seed_threshold, overlap_threshold)));
592 #else // FASTJET_ENABLE_PLUGIN_D0RUNICONE
593  is_unavailable("D0RunIpre96Cone");
594 #endif // FASTJET_ENABLE_PLUGIN_D0RUNICONE
595  }
596  if (all_algs || cmdline.present("-d0runicone")) {
597 #ifdef FASTJET_ENABLE_PLUGIN_D0RUNICONE
598  jet_defs.push_back( JetDefinition(new D0RunIConePlugin(ktR, seed_threshold, overlap_threshold)));
599 #else // FASTJET_ENABLE_PLUGIN_D0RUNICONE
600  is_unavailable("D0RunICone");
601 #endif // FASTJET_ENABLE_PLUGIN_D0RUNICONE
602  }
603  if (all_algs || cmdline.present("-gridjet")) {
604 #ifdef FASTJET_ENABLE_PLUGIN_GRIDJET
605  // we want a grid_ymax of 5.0, but when using R=0.4 (i.e. grid
606  // spacing of 0.8), this leads to 12.5 grid cells; depending on
607  // whether this is 12.499999999999 or 12.5000000....1 this gets
608  // converted either to 12 or 13, making the results sensitive to
609  // rounding errors.
610  //
611  // Instead we therefore take 4.9999999999, which avoids this problem.
612  double grid_ymax = 4.9999999999;
613  jet_defs.push_back( JetDefinition(new GridJetPlugin(grid_ymax, ktR*2.0)));
614 #else // FASTJET_ENABLE_PLUGIN_GRIDJET
615  is_unavailable("GridJet");
616 #endif // FASTJET_ENABLE_PLUGIN_GRIDJET
617 // end of checking if one asks to run a plugin (don't delete this line)
618  }
619  if (all_algs ||
620  cmdline.present("-kt") ||
621  (jet_defs.size() == 0 && !found_unavailable)) {
622  jet_defs.push_back( JetDefinition(kt_algorithm, ktR, E_scheme, strategy));
623  }
624 
625  string filename = cmdline.value<string>("-file", "");
626 
627 
628  if (!cmdline.all_options_used()) {cerr <<
629  "Error: some options were not recognized"<<endl;
630  exit(-1);}
631 
632  for (unsigned idef = 0; idef < jet_defs.size(); idef++) {
633  JetDefinition & jet_def = jet_defs[idef];
634  istream * istr;
635  if (filename == "") istr = &cin;
636  else istr = new ifstream(filename.c_str());
637 
638  for (int iev = 0; iev < nev; iev++) {
639  vector<PseudoJet> jets;
640  vector<PseudoJet> particles;
641  string line;
642  int ndone = 0;
643  while (getline(*istr, line)) {
644  //cout << line<<endl;
645  istringstream linestream(line);
646  if (line == "#END") {
647  ndone += 1;
648  if (ndone == combine) {break;}
649  }
650  if (line.substr(0,1) == "#") {continue;}
651  valarray<double> fourvec(4);
652  if (hydjet) {
653  // special reading from hydjet.txt event record (though actually
654  // this is supposed to be a standard pythia event record, so
655  // being able to read from it is perhaps not so bad an idea...)
656  int ii, istat,id,m1,m2,d1,d2;
657  double mass;
658  linestream >> ii>> istat >> id >> m1 >> m2 >> d1 >> d2
659  >> fourvec[0] >> fourvec[1] >> fourvec[2] >> mass;
660  // current file contains mass of particle as 4th entry
661  if (istat == 1) {
662  fourvec[3] = sqrt(+pow2(fourvec[0])+pow2(fourvec[1])
663  +pow2(fourvec[2])+pow2(mass));
664  }
665  } else {
666  if (massless) {
667  linestream >> fourvec[0] >> fourvec[1] >> fourvec[2];
668  fourvec[3] = sqrt(pow2(fourvec[0])+pow2(fourvec[1])+pow2(fourvec[2]));}
669  else {
670  linestream >> fourvec[0] >> fourvec[1] >> fourvec[2] >> fourvec[3];
671  }
672  }
673  PseudoJet psjet(fourvec);
674  if (abs(psjet.rap()) < rapmax) {particles.push_back(psjet);}
675  }
676 
677 #ifndef __FJCORE__
678  // add a fake underlying event which is very soft, uniformly distributed
679  // in eta,phi so as to allow one to reconstruct the area that is associated
680  // with each jet.
681  if (add_dense_coverage) {
682  GhostedAreaSpec ghosted_area_spec(ghost_maxrap);
683  //GhostedAreaSpec ghosted_area_spec(-2.0,4.0); // asymmetric range
684  // for plots, reduce the scatter default of 1, to avoid "holes"
685  // in the subsequent calorimeter view
686  ghosted_area_spec.set_grid_scatter(0.5);
687  ghosted_area_spec.add_ghosts(particles);
688  //----- old code ------------------
689  // srand(2);
690  // int nphi = 60;
691  // int neta = 100;
692  // double kt = 1e-1;
693  // for (int iphi = 0; iphi<nphi; iphi++) {
694  // for (int ieta = -neta; ieta<neta+1; ieta++) {
695  // double phi = (iphi+0.5) * (twopi/nphi) + rand()*0.001/RAND_MAX;
696  // double eta = ieta * (10.0/neta) + rand()*0.001/RAND_MAX;
697  // kt = 1e-20*(1+rand()*0.1/RAND_MAX);
698  // double pminus = kt*exp(-eta);
699  // double pplus = kt*exp(+eta);
700  // double px = kt*sin(phi);
701  // double py = kt*cos(phi);
702  // //cout << kt<<" "<<eta<<" "<<phi<<"\n";
703  // PseudoJet mom(px,py,0.5*(pplus-pminus),0.5*(pplus+pminus));
704  // particles.push_back(mom);
705  // }
706  // }
707  }
708 #else
709  add_dense_coverage = false;
710 #endif
711 
712  // select the particles that pass the selection cut
713  particles = particles_sel(particles);
714 
715  // allow user to skip some number of events (e.g. for easier bug-chasing)
716  if (iev < skip) continue;
717 
718  for (int irepeat = 0; irepeat < repeat ; irepeat++) {
719  int nparticles = particles.size();
720  try {
721  // one could use a unique_ptr here, but SharedPtr is available independently of C++ standard
722  SharedPtr<ClusterSequence> clust_seq;
723  if (do_areas) {
724 #ifndef __FJCORE__
725  clust_seq.reset(new ClusterSequenceArea(particles,jet_def,area_def));
726 #endif
727  } else {
728  clust_seq.reset(new ClusterSequence(particles,jet_def,write));
729  }
730  if (compare_strategy != plugin_strategy) {
731  do_compare_strategy(iev, particles, jet_def, *clust_seq, compare_strategy);
732  }
733 
734  // repetitive output
735  if (repeatinclkt >= 0.0) {
736  vector<PseudoJet> jets_local = sorted_by_pt(clust_seq->inclusive_jets(repeatinclkt));
737  }
738 
739  if (irepeat != 0) {continue;}
740  cout << "iev "<<iev<< ": number of particles = "<< nparticles << endl;
741  cout << "strategy used = "<< clust_seq->strategy_string()<< endl;
742  if (iev == 0) cout << "Jet Definition: " << jet_def.description() << " (" << fastjet_version_string() << ")" << endl;
743 #ifndef __FJCORE__
744  if (do_areas && iev == 0) cout << "Area definition: " << area_def.description() << endl;
745 #endif
746 
747  // now provide some nice output...
748  if (inclkt >= 0.0) {
749  vector<PseudoJet> jets_local = sorted_by_pt(clust_seq->inclusive_jets(inclkt));
750  print_jets(jets_local, show_constituents);
751 
752  }
753 
754  if (excln > 0) {
755  cout << "Printing "<<excln<<" exclusive jets\n";
756  print_jets(clust_seq->exclusive_jets(excln), show_constituents);
757  }
758 
759  if (excld > 0.0) {
760  cout << "Printing exclusive jets for d = "<<excld<<"\n";
761  print_jets(clust_seq->exclusive_jets(excld), show_constituents);
762  }
763 
764  if (excly > 0.0) {
765  cout << "Printing exclusive jets for ycut = "<<excly<<"\n";
766  print_jets(clust_seq->exclusive_jets_ycut(excly), show_constituents);
767  }
768 
769  if (get_all_dij) {
770  for (int i = nparticles-1; i >= 0; i--) {
771  printf("d for n = %4d -> %4d is %14.5e\n", i+1, i, clust_seq->exclusive_dmerge(i));
772  }
773  }
774  if (get_all_yij) {
775  for (int i = nparticles-1; i >= 0; i--) {
776  printf("y for n = %4d -> %4d is %14.5e\n", i+1, i, clust_seq->exclusive_ymerge(i));
777  }
778  }
779 
780  // have the option of printing out the subjets (at scale dcut) of
781  // each inclusive jet
782  if (subdcut >= 0.0) {
783  print_jets_and_sub(clust_seq->inclusive_jets(), subdcut);
784  }
785 
786  // useful for testing that recombination sequences are unique
787  if (unique_write) {
788  vector<int> unique_history = clust_seq->unique_history_order();
789  // construct the inverse of the above mapping
790  vector<int> inv_unique_history(clust_seq->history().size());
791  for (unsigned int i = 0; i < unique_history.size(); i++) {
792  inv_unique_history[unique_history[i]] = i;}
793 
794  for (unsigned int i = 0; i < unique_history.size(); i++) {
796  clust_seq->history()[unique_history[i]];
797  int uhp1 = el.parent1>=0 ? inv_unique_history[el.parent1] : el.parent1;
798  int uhp2 = el.parent2>=0 ? inv_unique_history[el.parent2] : el.parent2;
799  printf("%7d u %15.8e %7d u %7d u\n",i,el.dij,uhp1, uhp2);
800  }
801  }
802 
803 
804 #ifdef FASTJET_ENABLE_PLUGIN_SISCONE
805  // provide some complementary information for SISCone
806  if (show_cones) {
807  const SISConeExtras * extras =
808  dynamic_cast<const SISConeExtras *>(clust_seq->extras());
809  if (extras == 0)
810  throw fastjet::Error("extras object for SISCone was null (this can happen with certain area types)");
811  cout << "most ambiguous split (difference in squared dist) = "
812  << extras->most_ambiguous_split() << endl;
813  vector<fastjet::PseudoJet> stable_cones(extras->stable_cones());
814  stable_cones = sorted_by_rapidity(stable_cones);
815  for (unsigned int i = 0; i < stable_cones.size(); i++) {
816  //if (stable_cones[i].phi() < 5.0 && stable_cones[i].phi() > 4.0) {
817  printf("%5u %15.8f %15.8f %15.8f\n",
818  i,stable_cones[i].rap(),stable_cones[i].phi(),
819  stable_cones[i].perp() );
820  //}
821  }
822 
823  // also show passes for jets
824  vector<PseudoJet> sisjets = clust_seq->inclusive_jets();
825  printf("\n%15s %15s %15s %12s %8s %8s\n","rap","phi","pt","user-index","pass","nconst");
826  for (unsigned i = 0; i < sisjets.size(); i++) {
827  printf("%15.8f %15.8f %15.8f %12d %8d %8u\n",
828  sisjets[i].rap(), sisjets[i].phi(), sisjets[i].perp(),
829  sisjets[i].user_index(), extras->pass(sisjets[i]),
830  (unsigned int) clust_seq->constituents(sisjets[i]).size()
831  );
832 
833  }
834  }
835 #endif // FASTJET_ENABLE_PLUGIN_SISCONE
836 
837 #ifndef __FJCORE__
838  if (do_bkgd) {
839  double rho, sigma, mean_area, empty_area, n_empty_jets;
840  ClusterSequenceAreaBase * csab =
841  dynamic_cast<ClusterSequenceAreaBase *>(clust_seq.get());
842  BackgroundEstimatorBase * bge_ptr = 0;
843  if (do_bkgd_csab) {
844  csab->get_median_rho_and_sigma(bkgd_range, true, rho, sigma, mean_area);
845  empty_area = csab->empty_area(bkgd_range);
846  n_empty_jets = csab->n_empty_jets(bkgd_range);
847  } else if (do_bkgd_jetmedian) {
849  bge_ptr = bge;
850  // may be null
851  bge->set_rescaling_class(bkgd_rescaling);
852  bge->set_provide_fj2_sigma(do_bkgd_fj2);
853  if (bkgd_alt_ktR > 0) {
854  ClusterSequenceAreaBase * clust_seq_bkgd =
855  new ClusterSequenceArea(particles, JetDefinition(kt_algorithm, bkgd_alt_ktR), area_def);
856  cout << "Alt JetDef for background-estimation CSAB: " << clust_seq_bkgd->jet_def().description() << endl;
857  bge->set_cluster_sequence(*clust_seq_bkgd);
858  clust_seq_bkgd->delete_self_when_unused();
859  } else {
860  bge->set_cluster_sequence(*csab);
861  }
862  if (!do_bkgd_localrange) {
863  rho = bge->rho();
864  sigma = bge->sigma();
865  mean_area = bge->mean_area();
866  empty_area = bge->empty_area();
867  n_empty_jets = bge->n_empty_jets();
868  }
869  } else {
870  assert(do_bkgd_gridmedian);
871  double grid_rapmin, grid_rapmax;
872  bkgd_range.get_rapidity_extent(grid_rapmin, grid_rapmax);
873  GridMedianBackgroundEstimator * bge = new GridMedianBackgroundEstimator(grid_rapmax, 2*ktR);
874  bge_ptr = bge;
875  bge->set_rescaling_class(bkgd_rescaling);
876  bge->set_particles(particles);
877  rho = bge->rho();
878  sigma = bge->sigma();
879  mean_area = bge->mean_area();
880  empty_area = 0;
881  n_empty_jets = 0;
882  }
883  if (do_bkgd_localrange || do_subtractor) {
884  assert(bge_ptr != 0);
885  cout << "Background estimator: " << bge_ptr->description() << endl;
886  //vector<PseudoJet>
887  jets = SelectorAbsRapMax(3.0)(sorted_by_pt(csab->inclusive_jets()));
888  vector<PseudoJet> subjets;
889  if (do_subtractor) {
890  Subtractor subtractor(bge_ptr);
891  subtractor.set_use_rho_m(true);
892  subtractor.set_safe_mass(true);
893  cout << "Subtractor: " << subtractor.description() << endl;
894  subjets = subtractor(jets);
895  }
896  print_jets_bkgd(jets, subjets, bge_ptr, do_subtractor);
897  // cout << "i pt rap phi m rho rho_m sigma sigma_m" << endl;
898  // if (do_subtractor) cout << "isub ptsub rapsub phisub msub area" << endl;
899  // for (unsigned i = 0; i < jets.size(); i++) {
900  // const PseudoJet & jet = jets[i];
901  // cout << i << " "
902  // << " " << jet.pt() << " " << jet.rap() << " " << jet.phi() << " " << jet.m()
903  // << " " << bge_ptr->rho(jet) << " " << bge_ptr->rho_m(jet)
904  // << " " << bge_ptr->sigma(jet) << " " << bge_ptr->sigma_m(jet) << endl;
905  // if (do_subtractor) {
906  // const PseudoJet & subjet = subjets[i];
907  // cout << i << "sub"
908  // << " " << subjet.pt() << " " << subjet.rap() << " " << subjet.phi() << " " << subjet.m()
909  // << " " << jet.area() << endl;
910  // }
911  // }
912  } else {
913  cout << " rho = " << rho
914  << ", sigma = " << sigma
915  << ", mean_area = " << mean_area
916  << ", empty_area = " << empty_area
917  << ", n_empty_jets = " << n_empty_jets
918  << endl;
919  }
920  if (bge_ptr != 0) delete bge_ptr;
921  }
922 #endif
923  } // try
924  catch (Error &fjerr) {
925  cout << "Caught fastjet error, exiting gracefully" << endl;
926  exit(0);
927  }
928 
929  } // irepeat
930  } // iev
931  // if we've instantiated a plugin, delete it
932  if (jet_def.strategy()==plugin_strategy){
933  delete jet_def.plugin();
934  }
935  // close any file that we've opened
936  if (istr != &cin) delete istr;
937  } //
938 
939 }
940 
941 
942 
943 
944 //------ HELPER ROUTINES -----------------------------------------------
945 /// print a single jet
946 void print_jet (const PseudoJet & jet) {
947  unsigned int n_constituents = jet.constituents().size();
948  printf("%15.8f %15.8f %15.8f %8u\n",
949  jet.rap(), jet.phi(), jet.perp(), n_constituents);
950 }
951 
952 
953 //----------------------------------------------------------------------
954 void print_jets(const vector<PseudoJet> & jets_in, bool show_constituents) {
955  vector<PseudoJet> jets;
956  if (ee_print) {
957  jets = sorted_by_E(jets_in);
958  for (unsigned int j = 0; j < jets.size(); j++) {
959  printf("%5u %15.8f %15.8f %15.8f %15.8f\n",
960  j,jets[j].px(),jets[j].py(),jets[j].pz(),jets[j].E());
961  if (show_constituents) {
962  vector<PseudoJet> const_jets = jets[j].constituents();
963  for (unsigned int k = 0; k < const_jets.size(); k++) {
964  printf(" jet%03u %15.8f %15.8f %15.8f %15.8f\n",j,const_jets[k].px(),
965  const_jets[k].py(),const_jets[k].pz(),const_jets[k].E());
966  }
967  cout << "\n\n";
968  }
969 
970  }
971  } else {
972  jets = sorted_by_pt(jets_in);
973  for (unsigned int j = 0; j < jets.size(); j++) {
974  printf("%5u %15.8f %15.8f %15.8f",
975  j,jets[j].rap(),jets[j].phi(),jets[j].perp());
976  // also print out the scalar area and the perp component of the
977  // 4-vector (just enough to check a reasonable 4-vector?)
978 #ifndef __FJCORE__
979  if (do_areas) printf(" %15.8f %15.8f", jets[j].area(),
980  jets[j].area_4vector().perp());
981  cout << "\n";
982 #endif
983 
984  if (show_constituents) {
985  vector<PseudoJet> const_jets = jets[j].constituents();
986  for (unsigned int k = 0; k < const_jets.size(); k++) {
987  printf(" jet%03u %15.8f %15.8f %15.8f %5d\n",j,const_jets[k].rap(),
988  const_jets[k].phi(),sqrt(const_jets[k].kt2()), const_jets[k].cluster_hist_index());
989  }
990  cout << "\n\n";
991  }
992  }
993  }
994 
995  if (rootfile != "") {
996  ofstream ostr(rootfile.c_str());
997  ostr << "# " << cmdline_p->command_line() << endl;
998  ostr << "# output for root" << endl;
999  assert(jets.size() > 0);
1000  jets[0].validated_cs()->print_jets_for_root(jets,ostr);
1001  }
1002 
1003 }
1004 
1005 
1006 //----- SUBJETS --------------------------------------------------------
1007 /// a function that pretty prints a list of jets and the subjets for each
1008 /// one
1009 void print_jets_and_sub (const vector<PseudoJet> & jets, double dcut) {
1010 
1011  // sort jets into increasing pt
1012  vector<PseudoJet> sorted_jets = sorted_by_pt(jets);
1013 
1014  // label the columns
1015  printf("Printing jets and their subjets with subdcut = %10.5f\n",dcut);
1016  printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity",
1017  "phi", "pt", "n constituents");
1018 
1019  // the kind of subjet finding used to test consistency among them
1020  SubType sub_type = subtype_internal;
1021  //SubType sub_type = subtype_newclust_dcut;
1022  //SubType sub_type = subtype_newclust_R;
1023 
1024  // print out the details for each jet
1025  //for (unsigned int i = 0; i < sorted_jets.size(); i++) {
1026  for (vector<PseudoJet>::const_iterator jet = sorted_jets.begin();
1027  jet != sorted_jets.end(); jet++) {
1028  const JetDefinition & jet_def = jet->validated_cs()->jet_def();
1029 
1030  // if jet pt^2 < dcut with kt alg, then some methods of
1031  // getting subjets will return nothing -- so skip the jet
1032  if (jet_def.jet_algorithm() == kt_algorithm
1033  && jet->perp2() < dcut) continue;
1034 
1035 
1036  printf("%5u ",(unsigned int) (jet - sorted_jets.begin()));
1037  print_jet(*jet);
1038  vector<PseudoJet> subjets;
1039  ClusterSequence * cspoint;
1040  if (sub_type == subtype_internal) {
1041  cspoint = 0;
1042  subjets = jet->exclusive_subjets(dcut);
1043  double ddnp1 = jet->exclusive_subdmerge_max(subjets.size());
1044  double ddn = jet->exclusive_subdmerge_max(subjets.size()-1);
1045  cout << " for " << ddnp1 << " < d < " << ddn << " one has " << endl;
1046  } else if (sub_type == subtype_newclust_dcut) {
1047  cspoint = new ClusterSequence(jet->constituents(), jet_def);
1048  subjets = cspoint->exclusive_jets(dcut);
1049  } else if (sub_type == subtype_newclust_R) {
1050  assert(jet_def.jet_algorithm() == cambridge_algorithm);
1051  JetDefinition subjd(jet_def.jet_algorithm(), jet_def.R()*sqrt(dcut));
1052  cspoint = new ClusterSequence(jet->constituents(), subjd);
1053  subjets = cspoint->inclusive_jets();
1054  } else {
1055  cerr << "unrecognized subtype for subjet finding" << endl;
1056  exit(-1);
1057  }
1058 
1059  subjets = sorted_by_pt(subjets);
1060  for (unsigned int j = 0; j < subjets.size(); j++) {
1061  printf(" -sub-%02u ",j);
1062  print_jet(subjets[j]);
1063  }
1064 
1065  if (cspoint != 0) delete cspoint;
1066 
1067  //ClusterSequence subseq(clust_seq->constituents(sorted_jets[i]),
1068  // JetDefinition(cambridge_algorithm, 0.4));
1069  //vector<PseudoJet> subjets = sorted_by_pt(subseq.inclusive_jets());
1070  //for (unsigned int j = 0; j < subjets.size(); j++) {
1071  // printf(" -sub-%02u ",j);
1072  // print_jet(subseq, subjets[j]);
1073  //}
1074  }
1075 
1076 }
1077 
1078 /// if abs(x)<precision/2, return 0
1079 double make_safe_zero_truncation(double x, double precision){
1080  return std::abs(x)<0.5*precision ? 0.0 : x;
1081 }
1082 
1083 #ifndef __FJCORE__
1084 void print_jets_bkgd(const vector<PseudoJet> &jets,
1085  const vector<PseudoJet> &subtracted_jets,
1086  BackgroundEstimatorBase * bge_ptr,
1087  bool do_subtractor){
1088  printf("Printing jets, background information");
1089  if (do_subtractor)
1090  printf(" and subtracted jets\n");
1091  printf("%5s %15s %15s %15s %15s %15s %15s %15s %15s\n","jet #",
1092  "rapidity", "phi", "pt", "pt^2+m^2",
1093  "rho", "rho_m", "sigma", "sigma_m");
1094  if (do_subtractor)
1095  printf("%5s %15s %15s %15s %15s %15ss\n","jet #",
1096  "rapidity", "phi", "pt", "sqrt(pt^2+m^2)", "area");
1097 
1098  for (unsigned i = 0; i < jets.size(); i++) {
1099  const PseudoJet & jet = jets[i];
1100  BackgroundEstimate estimate = bge_ptr->estimate(jet);
1101  // Note that the values of rho_m sometimes comes out as +- a very
1102  // small number and the format can produce either 0.00000000 or
1103  // -0.00000000. The call to "make_safe_zero_truncation" makes sure it is
1104  // printed wo the - sign in each case
1105  printf("%5u %15.8f %15.8f %15.8f %15.8f %15.8f %15.8f %15.8f %15.8f\n", i,
1106  jet.rap(), jet.phi(), jet.perp(), jet.mt(),
1107  estimate.rho(), make_safe_zero_truncation(estimate.rho_m(),1e-8),
1108  estimate.sigma(), estimate.sigma_m());
1109  if (do_subtractor) {
1110  const PseudoJet & subjet = subtracted_jets[i];
1111  printf("%5u %15.8f %15.8f %15.8f %15.8f %15.8f\n", i,
1112  subjet.rap(), subjet.phi(), subjet.perp(), subjet.mt(), jet.area());
1113  }
1114  }
1115 }
1116 #endif// __FJCORE__
1117 
1118 //----------------------------------------------------------------------
1119 void signal_failed_comparison(int iev,
1120  const string & message,
1121  const vector<PseudoJet> & particles) {
1122  cout << "# failed comparison, reason is " << message << endl;
1123  cout << "# iev = " << iev << endl;
1124  cout << setprecision(16);
1125  for (unsigned i = 0; i < particles.size(); i++) {
1126  const PseudoJet & p = particles[i];
1127  cout << p.px() << " "
1128  << p.py() << " "
1129  << p.pz() << " "
1130  << p.E () << endl;
1131  }
1132  cout << "#END" << endl;
1133 }
1134 
1135 //----------------------------------------------------------------------
1136 void do_compare_strategy(int iev,
1137  const vector<PseudoJet> & particles,
1138  const JetDefinition & jet_def,
1139  const ClusterSequence & cs,
1140  int compare_strategy) {
1141 
1142  // create a jet def with the reference comparison strategy
1143  JetDefinition jet_def_ref(jet_def.jet_algorithm(),
1144  jet_def.R(),
1145  jet_def.recombination_scheme(),
1146  Strategy(compare_strategy));
1147  // do the clustering
1148  ClusterSequence cs_ref(particles, jet_def_ref);
1149 
1150  // now compare the outputs. At this stage just based on the clustering
1151  // sequence - get more sophisticated later...
1152  const vector<ClusterSequence::history_element> & history_in = cs.history();
1153  const vector<ClusterSequence::history_element> & history_ref = cs_ref.history();
1154 
1155  if (history_in.size() != history_ref.size()) {
1156  signal_failed_comparison(iev, "history sizes do not match", particles);
1157  }
1158 
1159  // now run over each clustering step to do the comparison
1160  for (unsigned i = cs.n_particles(); i < history_in.size(); i++) {
1161  bool fail_parents = (history_in[i].parent1 != history_ref[i].parent1 ||
1162  history_in[i].parent2 != history_ref[i].parent2);
1163  bool fail_dij = (history_in[i].dij != history_ref[i].dij);
1164  if (fail_parents || fail_dij) {
1165  ostringstream ostr;
1166  ostr << "at step " << i << ", ";
1167  if (fail_parents) ostr << "parents do not match, ";
1168  if (fail_dij) ostr << "dij does not match, ";
1169  ostr << "history in (p1, p2, dij) = "
1170  << history_in[i].parent1 << " " << history_in[i].parent2 << " " << history_in[i].dij;
1171  ostr << ", history ref (p1, p2, dij) = "
1172  << history_ref[i].parent1 << " " << history_ref[i].parent2 << " " << history_ref[i].dij;
1173  signal_failed_comparison(iev, ostr.str(), particles);
1174  break;
1175  }
1176  }
1177 }
int main()
an example program showing how to use fastjet
Definition: 01-basic.cc:50
Implementation of the ATLAS Cone (plugin for fastjet v2.4 upwards)
class that holds a generic area definition
/// a class that holds the result of the calculation
double sigma_m() const
fluctuations in the purely longitudinal (particle-mass-induced) component of the background density p...
double rho() const
background density per unit area
double sigma() const
background fluctuations per unit square-root area must be multipled by sqrt(area) to get fluctuations...
double rho_m() const
purely longitudinal (particle-mass-induced) component of the background density per unit area
Abstract base class that provides the basic interface for classes that estimate levels of background ...
virtual std::string description() const =0
returns a textual description of the background estimator
virtual BackgroundEstimate estimate() const =0
get the full set of background properties
A background rescaling that is a simple polynomial in y.
Implementation of the JetClu algorithm from CDF (plugin for fastjet-v2.1 upwards)
Implementation of the MidPoint algorithm from CDF (plugin for fastjet-v2.1 upwards)
Implementation of the CMS Iterative Cone (plugin for fastjet v2.4 upwards)
base class that sets interface for extensions of ClusterSequence that provide information about the a...
virtual void get_median_rho_and_sigma(const Selector &selector, bool use_area_4vector, double &median, double &sigma, double &mean_area) const
using jets withing the selector range (and with 4-vector areas if use_area_4vector),...
virtual double n_empty_jets(const Selector &selector) const
return something similar to the number of pure ghost jets in the given selector's range in an active ...
virtual double empty_area(const Selector &selector) const
return the total area, corresponding to the given Selector, that is free of jets, in general based on...
General class for user to obtain ClusterSequence with additional area information.
deals with clustering
const std::vector< history_element > & history() const
allow the user to access the raw internal history.
const JetDefinition & jet_def() const
return a reference to the jet definition
void delete_self_when_unused()
by calling this routine you tell the ClusterSequence to delete itself when all the Pseudojets associa...
unsigned int n_particles() const
returns the number of particles that were provided to the clustering algorithm (helps the user find t...
std::vector< PseudoJet > inclusive_jets(const double ptmin=0.0) const
return a vector of all jets (in the sense of the inclusive algorithm) with pt >= ptmin.
std::vector< PseudoJet > exclusive_jets(const double dcut) const
return a vector of all jets (in the sense of the exclusive algorithm) that would be obtained when run...
A plugin for FastJet (v3.0 or later) that provides an interface to the D0 version of Run-I cone algor...
Implementation of the D0 Run II Cone (plugin for fastjet v2.1 upwards)
A plugin for FastJet (v3.0 or later) that provides an interface to the pre 1996 D0 version of Run-I c...
Implementation of the e+e- Cambridge algorithm (plugin for fastjet v2.4 upwards)
base class corresponding to errors that can be thrown by FastJet
Definition: Error.hh:52
Parameters to configure the computation of jet areas using ghosts.
plugin for fastjet (v3.0 upwards) that clusters particles such that all particles in a given cell of ...
Background Estimator based on the median pt/area of a set of grid cells.
double sigma() const override
returns sigma, the background fluctuations per unit area; must be multipled by sqrt(area) to get fluc...
virtual void set_rescaling_class(const FunctionOfPseudoJet< double > *rescaling_class) override
Set a pointer to a class that calculates the rescaling factor as a function of the jet (position).
void set_particles(const std::vector< PseudoJet > &particles) override
tell the background estimator that it has a new event, composed of the specified particles.
double mean_area() const
returns the area of the grid cells (all identical, but referred to as "mean" area for uniformity with...
double rho() const override
returns rho, the median background density per unit area
Implementation of the e+e- Jade algorithm (plugin for fastjet v2.4 upwards)
Definition: JadePlugin.hh:77
Strategy
enum that contains the two clustering strategy options; for higher multiplicities,...
Definition: JadePlugin.hh:82
class that is intended to hold a full definition of the jet clusterer
std::string description() const
return a textual description of the current jet definition
JetAlgorithm jet_algorithm() const
return information about the definition...
const Plugin * plugin() const
return a pointer to the plugin
Class to estimate the pt density of the background per unit area, using the median of the distributio...
void set_provide_fj2_sigma(bool provide_fj2_sigma=true)
The FastJet v2.X sigma calculation had a small spurious offset in the limit of a small number of jets...
double n_empty_jets() const
Returns the number of empty jets used when computing the background properties.
double empty_area() const
Returns the estimate of the area (within the range defined by the selector) that is not occupied by j...
double sigma() const override
get sigma, the background fluctuations per unit area
double mean_area() const
Returns the mean area of the jets used to actually compute the background properties in the last call...
double rho() const override
get rho, the median background density per unit area
virtual void set_rescaling_class(const FunctionOfPseudoJet< double > *rescaling_class_in) override
Set a pointer to a class that calculates the rescaling factor as a function of the jet (position).
void set_cluster_sequence(const ClusterSequenceAreaBase &csa)
(re)set the cluster sequence (with area support) to be used by future calls to rho() etc.
Class to contain pseudojets, including minimal information of use to jet-clustering routines.
Definition: PseudoJet.hh:68
virtual std::vector< PseudoJet > constituents() const
retrieve the constituents.
Definition: PseudoJet.cc:692
double mt() const
returns the transverse mass = sqrt(kt^2+m^2)
Definition: PseudoJet.hh:175
double rap() const
returns the rapidity or some large value when the rapidity is infinite
Definition: PseudoJet.hh:138
double phi() const
returns phi (in the range 0..2pi)
Definition: PseudoJet.hh:123
double perp() const
returns the scalar transverse momentum
Definition: PseudoJet.hh:158
double perp2() const
returns the squared transverse momentum
Definition: PseudoJet.hh:156
double exclusive_subdmerge_max(int nsub) const
Returns the maximum dij that occurred in the whole event at the stage that the nsub+1 -> nsub merge o...
Definition: PseudoJet.cc:774
virtual double area() const
return the jet (scalar) area.
Definition: PseudoJet.cc:829
std::vector< PseudoJet > exclusive_subjets(const double dcut) const
return a vector of all subjets of the current jet (in the sense of the exclusive algorithm) that woul...
Definition: PseudoJet.cc:715
const ClusterSequence * validated_cs() const
shorthand for validated_cluster_sequence()
Definition: PseudoJet.cc:565
Implementation of the PxCone algorithm (plugin for fastjet v2.1 upwards)
Definition: PxConePlugin.hh:69
int pass(const PseudoJet &jet) const
return the # of the pass at which a given jet was found; will return -1 if the pass is invalid
double most_ambiguous_split() const
return the smallest difference in squared distance encountered during splitting between a particle an...
const std::vector< PseudoJet > & stable_cones() const
returns a reference to the vector of stable cones (aka protocones)
Class that provides extra information about a SISCone clustering.
Implementation of the SISCone algorithm (plugin for fastjet v2.1 upwards)
Implementation of the spherical version of the SISCone algorithm (plugin for fastjet v2....
Class that encodes information about cuts and other selection criteria that can be applied to PseudoJ...
Definition: Selector.hh:149
void get_rapidity_extent(double &rapmin, double &rapmax) const
returns the rapidity range for which it may return "true"
Definition: Selector.hh:232
An implementation of shared pointers that is broadly similar to C++11 shared_ptr (https://en....
Definition: SharedPtr.hh:341
T * get() const
get the stored pointer
Definition: SharedPtr.hh:473
void reset()
reset the pointer to default value (NULL)
Definition: SharedPtr.hh:381
Class that helps perform jet background subtraction.
Definition: Subtractor.hh:62
Implementation of the TrackJet algorithm (plugin for fastjet v2.4 upwards)
Specification for the computation of the Voronoi jet area.
void print_jets(const vector< fastjet::PseudoJet > &)
a function that pretty prints a list of jets
Selector SelectorNHardest(unsigned int n)
select the n hardest objects
Definition: Selector.cc:1074
Selector SelectorStrip(const double half_width)
select objets within a rapidity distance 'half_width' from the location of the reference jet,...
Definition: Selector.cc:1267
Selector SelectorAbsRapMax(double absrapmax)
select objects with |rap| <= absrapmax
Definition: Selector.cc:880
the FastJet namespace
Strategy
the various options for the algorithmic strategy to adopt in clustering events with kt and cambridge ...
@ plugin_strategy
the plugin has been used...
vector< PseudoJet > sorted_by_pt(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing kt2
Definition: PseudoJet.cc:882
RecombinationScheme
The various recombination schemes.
@ E_scheme
summing the 4-momenta
vector< PseudoJet > sorted_by_rapidity(const vector< PseudoJet > &jets)
return a vector of jets sorted into increasing rapidity
Definition: PseudoJet.cc:890
string fastjet_version_string()
return a string containing information about the release
vector< PseudoJet > sorted_by_E(const vector< PseudoJet > &jets)
return a vector of jets sorted into decreasing energy
Definition: PseudoJet.cc:898
a single element in the clustering history
int parent2
index in _history where first parent of this jet was created (InexistentParent if this jet is an orig...
double dij
index in the _jets vector where we will find the