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