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