FastJet 3.5.0
Loading...
Searching...
No Matches
fastjet_timing_plugins.cc
1//FJSTARTHEADER
2// $Id$
3//
4// Copyright (c) 2005-2025, 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;
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 int prec = cmdline.int_val("-prec",-1);
358 if (prec > 0) cout.precision(prec);
359 // have the option of comparing the clustering results to those
360 // obtained with a different clustering strategy; the misuse the
361 // "plugin_strategy" to indicate that no comparison is needed.
362 // Does not currently support areas
363 int compare_strategy = cmdline.value<int>("-compare-strategy", plugin_strategy);
364
365
366 Selector particles_sel = (cmdline.present("-nhardest"))
367 ? SelectorNHardest(cmdline.value<unsigned int>("-nhardest"))
368 : SelectorIdentity();
369
370 do_areas = cmdline.present("-area");
371#ifndef __FJCORE__
372 AreaDefinition area_def;
373 if (do_areas) {
374 assert(!write); // it's incompatible
375 GhostedAreaSpec ghost_spec(ghost_maxrap,
376 cmdline.value("-area:repeat", 1),
377 cmdline.value("-ghost-area", 0.01));
378 if (cmdline.present("-area:fj2")) ghost_spec.set_fj2_placement(true);
379 if (cmdline.present("-area:explicit")) {
380 area_def = AreaDefinition(active_area_explicit_ghosts, ghost_spec);
381 } else if (cmdline.present("-area:passive")) {
382 area_def = AreaDefinition(passive_area, ghost_spec);
383 } else if (cmdline.present("-area:voronoi")) {
384 double Rfact = cmdline.value<double>("-area:voronoi");
385 area_def = AreaDefinition(voronoi_area,
386 VoronoiAreaSpec(Rfact));
387 } else {
388 cmdline.present("-area:active"); // allow, but do not require, arg
389 area_def = AreaDefinition(active_area, ghost_spec);
390 }
391 }
392#else
393 do_areas=false;
394#endif
395
396 bool do_bkgd = cmdline.present("-bkgd"); // background estimation
397#ifndef __FJCORE__
398 bool do_bkgd_csab = false, do_bkgd_jetmedian = false, do_bkgd_fj2 = false;
399 bool do_bkgd_gridmedian = false;
400 bool do_bkgd_localrange = false;
401 bool do_subtractor = false;
402 double bkgd_alt_ktR = -1.0;
403 BackgroundRescalingYPolynomial * bkgd_rescaling = 0;
404 Selector bkgd_range;
405 if (do_bkgd) {
406 bkgd_range = SelectorAbsRapMax(ghost_maxrap - ktR);
407 if (cmdline.present("-bkgd:csab")) {do_bkgd_csab = true;}
408 else if (cmdline.present("-bkgd:jetmedian")) {do_bkgd_jetmedian = true;
409 do_bkgd_fj2 = cmdline.present("-bkgd:fj2");
410 do_bkgd_localrange = cmdline.present("-bkgd:localrange");
411 bkgd_alt_ktR = cmdline.value("-bkgd:alt-ktR", bkgd_alt_ktR);
412 if (do_bkgd_localrange) bkgd_range = SelectorStrip(1.5);
413 } else if (cmdline.present("-bkgd:gridmedian")) {
414 do_bkgd_gridmedian = true;
415 } else {
416 throw Error("with the -bkgd option, some particular background must be specified (csab or jetmedian)");
417 }
418 if (cmdline.present("-bkgd:rescaling")) {
419 bkgd_rescaling = new BackgroundRescalingYPolynomial(1.157,0,-0.0266,0,0.000048);
420 }
421 assert(do_areas || do_bkgd_gridmedian);
422 do_subtractor = cmdline.present("-subtractor");
423 if (do_subtractor) assert(do_areas);
424 }
425#else
426 do_bkgd = false;
427#endif
428
429 bool show_cones = cmdline.present("-cones"); // only works for siscone
430
431 // for cone algorithms
432 // allow -f and -overlap
433 double overlap_threshold = cmdline.double_val("-overlap",0.5);
434 overlap_threshold = cmdline.double_val("-f",overlap_threshold);
435 double seed_threshold = cmdline.double_val("-seed",1.0);
436
437#ifdef __FJCORE__
438 show_cones = false;
439#endif
440
441 // for ee algorithms, allow to specify ycut
442 double ycut = cmdline.double_val("-ycut",0.08);
443
444 // for printing jets to a file for reading by root
445 rootfile = cmdline.value<string>("-root","");
446
447 // out default scheme is the E_scheme
449
450 // The following option causes the Cambridge algo to be used.
451 // Note that currently the only output that works sensibly here is
452 // "-incl 0"
453 vector<JetDefinition> jet_defs;
454 if (all_algs || cmdline.present("-cam") || cmdline.present("-CA")) {
455 jet_defs.push_back( JetDefinition(cambridge_algorithm, ktR, scheme, strategy));
456 }
457 if (all_algs || cmdline.present("-antikt")) {
458 jet_defs.push_back( JetDefinition(antikt_algorithm, ktR, scheme, strategy));
459 }
460 if (all_algs || cmdline.present("-genkt")) {
461 double p;
462 if (cmdline.present("-genkt")) p = cmdline.value<double>("-genkt");
463 else p = -0.5;
464 jet_defs.push_back( JetDefinition(genkt_algorithm, ktR, p, scheme, strategy));
465 }
466 if (all_algs || cmdline.present("-eekt")) {
467 jet_defs.push_back( JetDefinition(ee_kt_algorithm));
468 }
469 if (all_algs || cmdline.present("-eegenkt")) {
470 double p;
471 if (cmdline.present("-eegenkt")) p = cmdline.value<double>("-eegenkt");
472 else p = -0.5;
473 jet_defs.push_back( JetDefinition(ee_genkt_algorithm, ktR, p, scheme, strategy));
474
475// checking if one asks to run a plugin (don't delete this line)
476 }
477 if (all_algs || cmdline.present("-midpoint")) {
478#ifdef FASTJET_ENABLE_PLUGIN_CDFCONES
479 typedef CDFMidPointPlugin MPPlug; // for brevity
480 double cone_area_fraction = 1.0;
481 int max_pair_size = 2;
482 int max_iterations = 100;
483 MPPlug::SplitMergeScale sm_scale = MPPlug::SM_pt;
484 if (cmdline.present("-sm-pttilde")) sm_scale = MPPlug::SM_pttilde;
485 if (cmdline.present("-sm-pt")) sm_scale = MPPlug::SM_pt; // default
486 if (cmdline.present("-sm-mt")) sm_scale = MPPlug::SM_mt;
487 if (cmdline.present("-sm-Et")) sm_scale = MPPlug::SM_Et;
488 jet_defs.push_back( JetDefinition( new CDFMidPointPlugin (
489 seed_threshold, ktR,
490 cone_area_fraction, max_pair_size,
491 max_iterations, overlap_threshold,
492 sm_scale)));
493#else // FASTJET_ENABLE_PLUGIN_CDFCONES
494 is_unavailable("MidPoint");
495#endif // FASTJET_ENABLE_PLUGIN_CDFCONES
496 }
497 if (all_algs || cmdline.present("-pxcone")) {
498#ifdef FASTJET_ENABLE_PLUGIN_PXCONE
499 double min_jet_energy = 5.0;
500 // mode: 1=e+e-, 2=pp
501 int mode = cmdline.value("-pxcone-mode", 2);
502 bool E_scheme_jets = false;
503 jet_defs.push_back( JetDefinition( new PxConePlugin (
504 ktR, min_jet_energy,
505 overlap_threshold, E_scheme_jets, mode)));
506#else // FASTJET_ENABLE_PLUGIN_PXCONE
507 is_unavailable("PxCone");
508#endif // FASTJET_ENABLE_PLUGIN_PXCONE
509 }
510 if (all_algs || cmdline.present("-jetclu")) {
511#ifdef FASTJET_ENABLE_PLUGIN_CDFCONES
512 jet_defs.push_back( JetDefinition( new CDFJetCluPlugin (
513 ktR, overlap_threshold, seed_threshold)));
514#else // FASTJET_ENABLE_PLUGIN_CDFCONES
515 is_unavailable("JetClu");
516#endif // FASTJET_ENABLE_PLUGIN_CDFCONES
517 }
518 if (all_algs || cmdline.present("-siscone") || cmdline.present("-sisconespheri")) {
519#ifdef FASTJET_ENABLE_PLUGIN_SISCONE
520 typedef SISConePlugin SISPlug; // for brevity
521 int npass = cmdline.value("-npass",0);
522 if (all_algs || cmdline.present("-siscone")) {
523 double sisptmin = cmdline.value("-sisptmin",0.0);
524 bool cache = cmdline.present("-cache");
525 SISPlug * plugin = new SISPlug (ktR, overlap_threshold,npass,sisptmin,cache);
526 if (cmdline.present("-sm-pt")) plugin->set_split_merge_scale(SISPlug::SM_pt);
527 if (cmdline.present("-sm-mt")) plugin->set_split_merge_scale(SISPlug::SM_mt);
528 if (cmdline.present("-sm-Et")) plugin->set_split_merge_scale(SISPlug::SM_Et);
529 if (cmdline.present("-sm-pttilde")) plugin->set_split_merge_scale(SISPlug::SM_pttilde);
530 // cause it to use the jet-definition's own recombiner
531 plugin->set_use_jet_def_recombiner(true);
532 jet_defs.push_back( JetDefinition(plugin));
533 }
534 if (all_algs || cmdline.present("-sisconespheri")) {
535 double sisEmin = cmdline.value("-sisEmin",0.0);
536 SISConeSphericalPlugin * plugin =
537 new SISConeSphericalPlugin(ktR, overlap_threshold,npass,sisEmin);
538 if (cmdline.present("-ghost-sep")) {
539 plugin->set_ghost_separation_scale(cmdline.value<double>("-ghost-sep"));
540 }
541 jet_defs.push_back( JetDefinition(plugin));
542 }
543#else // FASTJET_ENABLE_PLUGIN_SISCONE
544 is_unavailable("SISCone");
545#endif // FASTJET_ENABLE_PLUGIN_SISCONE
546 }
547 if (all_algs || cmdline.present("-d0runiicone")) {
548#ifdef FASTJET_ENABLE_PLUGIN_D0RUNIICONE
549 double min_jet_Et = 6.0; // was 8 GeV in earlier work
550 jet_defs.push_back( JetDefinition(new D0RunIIConePlugin(ktR,min_jet_Et)));
551#else // FASTJET_ENABLE_PLUGIN_D0RUNIICONE
552 is_unavailable("D0RunIICone");
553#endif // FASTJET_ENABLE_PLUGIN_D0RUNIICONE
554 }
555 if (all_algs || cmdline.present("-trackjet")) {
556#ifdef FASTJET_ENABLE_PLUGIN_TRACKJET
557 jet_defs.push_back( JetDefinition(new TrackJetPlugin(ktR)));
558#else // FASTJET_ENABLE_PLUGIN_TRACKJET
559 is_unavailable("TrackJet");
560#endif // FASTJET_ENABLE_PLUGIN_TRACKJET
561 }
562 if (all_algs || cmdline.present("-atlascone")) {
563#ifdef FASTJET_ENABLE_PLUGIN_ATLASCONE
564 jet_defs.push_back( JetDefinition(new ATLASConePlugin(ktR)));
565#else // FASTJET_ENABLE_PLUGIN_ATLASCONE
566 is_unavailable("ATLASCone");
567#endif // FASTJET_ENABLE_PLUGIN_ATLASCONE
568 }
569 if (all_algs || cmdline.present("-eecambridge")) {
570#ifdef FASTJET_ENABLE_PLUGIN_EECAMBRIDGE
571 jet_defs.push_back( JetDefinition(new EECambridgePlugin(ycut)));
572#else // FASTJET_ENABLE_PLUGIN_EECAMBRIDGE
573 is_unavailable("EECambridge");
574#endif // FASTJET_ENABLE_PLUGIN_EECAMBRIDGE
575 }
576 if (all_algs || cmdline.present("-jade")) {
577#ifdef FASTJET_ENABLE_PLUGIN_JADE
578 JadePlugin::Strategy jade_strategy =
579 JadePlugin::Strategy(cmdline.value<int>("-jade-strategy",
580 JadePlugin::strategy_NNFJN2Plain));
581 jet_defs.push_back( JetDefinition(new JadePlugin(jade_strategy)));
582#else // FASTJET_ENABLE_PLUGIN_JADE
583 is_unavailable("Jade");
584#endif // FASTJET_ENABLE_PLUGIN_JADE
585 }
586 if (all_algs || cmdline.present("-cmsiterativecone")) {
587#ifdef FASTJET_ENABLE_PLUGIN_CMSITERATIVECONE
588 jet_defs.push_back( JetDefinition(new CMSIterativeConePlugin(ktR,seed_threshold)));
589#else // FASTJET_ENABLE_PLUGIN_CMSITERATIVECONE
590 is_unavailable("CMSIterativeCone");
591#endif // FASTJET_ENABLE_PLUGIN_CMSITERATIVECONE
592 }
593 if (all_algs || cmdline.present("-d0runipre96cone")) {
594#ifdef FASTJET_ENABLE_PLUGIN_D0RUNICONE
595 jet_defs.push_back( JetDefinition(new D0RunIpre96ConePlugin(ktR, seed_threshold, overlap_threshold)));
596#else // FASTJET_ENABLE_PLUGIN_D0RUNICONE
597 is_unavailable("D0RunIpre96Cone");
598#endif // FASTJET_ENABLE_PLUGIN_D0RUNICONE
599 }
600 if (all_algs || cmdline.present("-d0runicone")) {
601#ifdef FASTJET_ENABLE_PLUGIN_D0RUNICONE
602 jet_defs.push_back( JetDefinition(new D0RunIConePlugin(ktR, seed_threshold, overlap_threshold)));
603#else // FASTJET_ENABLE_PLUGIN_D0RUNICONE
604 is_unavailable("D0RunICone");
605#endif // FASTJET_ENABLE_PLUGIN_D0RUNICONE
606 }
607 if (all_algs || cmdline.present("-gridjet")) {
608#ifdef FASTJET_ENABLE_PLUGIN_GRIDJET
609 // we want a grid_ymax of 5.0, but when using R=0.4 (i.e. grid
610 // spacing of 0.8), this leads to 12.5 grid cells; depending on
611 // whether this is 12.499999999999 or 12.5000000....1 this gets
612 // converted either to 12 or 13, making the results sensitive to
613 // rounding errors.
614 //
615 // Instead we therefore take 4.9999999999, which avoids this problem.
616 double grid_ymax = 4.9999999999;
617 jet_defs.push_back( JetDefinition(new GridJetPlugin(grid_ymax, ktR*2.0)));
618#else // FASTJET_ENABLE_PLUGIN_GRIDJET
619 is_unavailable("GridJet");
620#endif // FASTJET_ENABLE_PLUGIN_GRIDJET
621// end of checking if one asks to run a plugin (don't delete this line)
622 }
623 if (all_algs ||
624 cmdline.present("-kt") ||
625 (jet_defs.size() == 0 && !found_unavailable)) {
626 jet_defs.push_back( JetDefinition(kt_algorithm, ktR, E_scheme, strategy));
627 }
628
629 string filename = cmdline.value<string>("-file", "");
630
631
632 if (!cmdline.all_options_used()) {cerr <<
633 "Error: some options were not recognized"<<endl;
634 exit(-1);}
635
636 for (unsigned idef = 0; idef < jet_defs.size(); idef++) {
637 JetDefinition & jet_def = jet_defs[idef];
638 istream * istr;
639 if (filename == "") istr = &cin;
640 else istr = new ifstream(filename.c_str());
641
642 for (int iev = 0; iev < nev; iev++) {
643 vector<PseudoJet> jets;
644 vector<PseudoJet> particles;
645 string line;
646 int ndone = 0;
647 while (getline(*istr, line)) {
648 //cout << line<<endl;
649 istringstream linestream(line);
650 if (line == "#END") {
651 ndone += 1;
652 if (ndone == combine) {break;}
653 }
654 if (line.substr(0,1) == "#") {continue;}
655 valarray<double> fourvec(4);
656 if (hydjet) {
657 // special reading from hydjet.txt event record (though actually
658 // this is supposed to be a standard pythia event record, so
659 // being able to read from it is perhaps not so bad an idea...)
660 int ii, istat,id,m1,m2,d1,d2;
661 double mass;
662 linestream >> ii>> istat >> id >> m1 >> m2 >> d1 >> d2
663 >> fourvec[0] >> fourvec[1] >> fourvec[2] >> mass;
664 // current file contains mass of particle as 4th entry
665 if (istat == 1) {
666 fourvec[3] = sqrt(+pow2(fourvec[0])+pow2(fourvec[1])
667 +pow2(fourvec[2])+pow2(mass));
668 }
669 } else {
670 if (massless) {
671 linestream >> fourvec[0] >> fourvec[1] >> fourvec[2];
672 fourvec[3] = sqrt(pow2(fourvec[0])+pow2(fourvec[1])+pow2(fourvec[2]));}
673 else {
674 linestream >> fourvec[0] >> fourvec[1] >> fourvec[2] >> fourvec[3];
675 }
676 }
677 PseudoJet psjet(fourvec);
678 if (abs(psjet.rap()) < rapmax) {particles.push_back(psjet);}
679 }
680
681#ifndef __FJCORE__
682 // add a fake underlying event which is very soft, uniformly distributed
683 // in eta,phi so as to allow one to reconstruct the area that is associated
684 // with each jet.
685 if (add_dense_coverage) {
686 GhostedAreaSpec ghosted_area_spec(ghost_maxrap);
687 //GhostedAreaSpec ghosted_area_spec(-2.0,4.0); // asymmetric range
688 // for plots, reduce the scatter default of 1, to avoid "holes"
689 // in the subsequent calorimeter view
690 ghosted_area_spec.set_grid_scatter(0.5);
691 ghosted_area_spec.add_ghosts(particles);
692 //----- old code ------------------
693 // srand(2);
694 // int nphi = 60;
695 // int neta = 100;
696 // double kt = 1e-1;
697 // for (int iphi = 0; iphi<nphi; iphi++) {
698 // for (int ieta = -neta; ieta<neta+1; ieta++) {
699 // double phi = (iphi+0.5) * (twopi/nphi) + rand()*0.001/RAND_MAX;
700 // double eta = ieta * (10.0/neta) + rand()*0.001/RAND_MAX;
701 // kt = 1e-20*(1+rand()*0.1/RAND_MAX);
702 // double pminus = kt*exp(-eta);
703 // double pplus = kt*exp(+eta);
704 // double px = kt*sin(phi);
705 // double py = kt*cos(phi);
706 // //cout << kt<<" "<<eta<<" "<<phi<<"\n";
707 // PseudoJet mom(px,py,0.5*(pplus-pminus),0.5*(pplus+pminus));
708 // particles.push_back(mom);
709 // }
710 // }
711 }
712#else
713 add_dense_coverage = false;
714#endif
715
716 // select the particles that pass the selection cut
717 particles = particles_sel(particles);
718
719 // allow user to skip some number of events (e.g. for easier bug-chasing)
720 if (iev < skip) continue;
721
722 for (int irepeat = 0; irepeat < repeat ; irepeat++) {
723 int nparticles = particles.size();
724 try {
725 // one could use a unique_ptr here, but SharedPtr is available independently of C++ standard
727 if (do_areas) {
728#ifndef __FJCORE__
729 clust_seq.reset(new ClusterSequenceArea(particles,jet_def,area_def));
730#endif
731 } else {
732 clust_seq.reset(new ClusterSequence(particles,jet_def,write));
733 }
734 if (compare_strategy != plugin_strategy) {
735 do_compare_strategy(iev, particles, jet_def, *clust_seq, compare_strategy);
736 }
737
738 // repetitive output
739 if (repeatinclkt >= 0.0) {
740 vector<PseudoJet> jets_local = sorted_by_pt(clust_seq->inclusive_jets(repeatinclkt));
741 }
742
743 if (irepeat != 0) {continue;}
744 cout << "iev "<<iev<< ": number of particles = "<< nparticles << endl;
745 cout << "strategy used = "<< clust_seq->strategy_string()<< endl;
746 if (iev == 0) cout << "Jet Definition: " << jet_def.description() << " (" << fastjet_version_string() << ")" << endl;
747#ifndef __FJCORE__
748 if (do_areas && iev == 0) cout << "Area definition: " << area_def.description() << endl;
749#endif
750
751 // now provide some nice output...
752 if (inclkt >= 0.0) {
753 vector<PseudoJet> jets_local = sorted_by_pt(clust_seq->inclusive_jets(inclkt));
754 print_jets(jets_local, show_constituents);
755
756 }
757
758 if (excln > 0) {
759 cout << "Printing "<<excln<<" exclusive jets\n";
760 print_jets(clust_seq->exclusive_jets(excln), show_constituents);
761 }
762
763 if (excld > 0.0) {
764 cout << "Printing exclusive jets for d = "<<excld<<"\n";
765 print_jets(clust_seq->exclusive_jets(excld), show_constituents);
766 }
767
768 if (excly > 0.0) {
769 cout << "Printing exclusive jets for ycut = "<<excly<<"\n";
770 print_jets(clust_seq->exclusive_jets_ycut(excly), show_constituents);
771 }
772
773 if (get_all_dij) {
774 for (int i = nparticles-1; i >= 0; i--) {
775 printf("d for n = %4d -> %4d is %14.5e\n", i+1, i, clust_seq->exclusive_dmerge(i));
776 }
777 }
778 if (get_all_yij) {
779 for (int i = nparticles-1; i >= 0; i--) {
780 printf("y for n = %4d -> %4d is %14.5e\n", i+1, i, clust_seq->exclusive_ymerge(i));
781 }
782 }
783
784 // have the option of printing out the subjets (at scale dcut) of
785 // each inclusive jet
786 if (subdcut >= 0.0) {
787 print_jets_and_sub(clust_seq->inclusive_jets(), subdcut);
788 }
789
790 // useful for testing that recombination sequences are unique
791 if (unique_write) {
792 vector<int> unique_history = clust_seq->unique_history_order();
793 // construct the inverse of the above mapping
794 vector<int> inv_unique_history(clust_seq->history().size());
795 for (unsigned int i = 0; i < unique_history.size(); i++) {
796 inv_unique_history[unique_history[i]] = i;}
797
798 for (unsigned int i = 0; i < unique_history.size(); i++) {
800 clust_seq->history()[unique_history[i]];
801 int uhp1 = el.parent1>=0 ? inv_unique_history[el.parent1] : el.parent1;
802 int uhp2 = el.parent2>=0 ? inv_unique_history[el.parent2] : el.parent2;
803 printf("%7d u %15.8e %7d u %7d u\n",i,el.dij,uhp1, uhp2);
804 }
805 }
806
807
808#ifdef FASTJET_ENABLE_PLUGIN_SISCONE
809 // provide some complementary information for SISCone
810 if (show_cones) {
811 const SISConeExtras * extras =
812 dynamic_cast<const SISConeExtras *>(clust_seq->extras());
813 if (extras == 0)
814 throw fastjet::Error("extras object for SISCone was null (this can happen with certain area types)");
815 cout << "most ambiguous split (difference in squared dist) = "
816 << extras->most_ambiguous_split() << endl;
817 vector<fastjet::PseudoJet> stable_cones(extras->stable_cones());
818 stable_cones = sorted_by_rapidity(stable_cones);
819 for (unsigned int i = 0; i < stable_cones.size(); i++) {
820 //if (stable_cones[i].phi() < 5.0 && stable_cones[i].phi() > 4.0) {
821 printf("%5u %15.8f %15.8f %15.8f\n",
822 i,stable_cones[i].rap(),stable_cones[i].phi(),
823 stable_cones[i].perp() );
824 //}
825 }
826
827 // also show passes for jets
828 vector<PseudoJet> sisjets = clust_seq->inclusive_jets();
829 printf("\n%15s %15s %15s %12s %8s %8s\n","rap","phi","pt","user-index","pass","nconst");
830 for (unsigned i = 0; i < sisjets.size(); i++) {
831 printf("%15.8f %15.8f %15.8f %12d %8d %8u\n",
832 sisjets[i].rap(), sisjets[i].phi(), sisjets[i].perp(),
833 sisjets[i].user_index(), extras->pass(sisjets[i]),
834 (unsigned int) clust_seq->constituents(sisjets[i]).size()
835 );
836
837 }
838 }
839#endif // FASTJET_ENABLE_PLUGIN_SISCONE
840
841#ifndef __FJCORE__
842 if (do_bkgd) {
843 double rho, sigma, mean_area, empty_area, n_empty_jets;
845 dynamic_cast<ClusterSequenceAreaBase *>(clust_seq.get());
846 BackgroundEstimatorBase * bge_ptr = 0;
847 if (do_bkgd_csab) {
848 csab->get_median_rho_and_sigma(bkgd_range, true, rho, sigma, mean_area);
849 empty_area = csab->empty_area(bkgd_range);
850 n_empty_jets = csab->n_empty_jets(bkgd_range);
851 } else if (do_bkgd_jetmedian) {
853 bge_ptr = bge;
854 // may be null
855 bge->set_rescaling_class(bkgd_rescaling);
856 bge->set_provide_fj2_sigma(do_bkgd_fj2);
857 if (bkgd_alt_ktR > 0) {
858 ClusterSequenceAreaBase * clust_seq_bkgd =
859 new ClusterSequenceArea(particles, JetDefinition(kt_algorithm, bkgd_alt_ktR), area_def);
860 cout << "Alt JetDef for background-estimation CSAB: " << clust_seq_bkgd->jet_def().description() << endl;
861 bge->set_cluster_sequence(*clust_seq_bkgd);
862 clust_seq_bkgd->delete_self_when_unused();
863 } else {
864 bge->set_cluster_sequence(*csab);
865 }
866 if (!do_bkgd_localrange) {
867 rho = bge->rho();
868 sigma = bge->sigma();
869 mean_area = bge->mean_area();
870 empty_area = bge->empty_area();
871 n_empty_jets = bge->n_empty_jets();
872 }
873 } else {
874 assert(do_bkgd_gridmedian);
875 double grid_rapmin, grid_rapmax;
876 bkgd_range.get_rapidity_extent(grid_rapmin, grid_rapmax);
877 GridMedianBackgroundEstimator * bge = new GridMedianBackgroundEstimator(grid_rapmax, 2*ktR);
878 bge_ptr = bge;
879 bge->set_rescaling_class(bkgd_rescaling);
880 bge->set_particles(particles);
881 rho = bge->rho();
882 sigma = bge->sigma();
883 mean_area = bge->mean_area();
884 empty_area = 0;
885 n_empty_jets = 0;
886 }
887 if (do_bkgd_localrange || do_subtractor) {
888 assert(bge_ptr != 0);
889 cout << "Background estimator: " << bge_ptr->description() << endl;
890 //vector<PseudoJet>
891 jets = SelectorAbsRapMax(3.0)(sorted_by_pt(csab->inclusive_jets()));
892 vector<PseudoJet> subjets;
893 if (do_subtractor) {
894 Subtractor subtractor(bge_ptr);
895 subtractor.set_use_rho_m(true);
896 subtractor.set_safe_mass(true);
897 cout << "Subtractor: " << subtractor.description() << endl;
898 subjets = subtractor(jets);
899 }
900 print_jets_bkgd(jets, subjets, bge_ptr, do_subtractor);
901 // cout << "i pt rap phi m rho rho_m sigma sigma_m" << endl;
902 // if (do_subtractor) cout << "isub ptsub rapsub phisub msub area" << endl;
903 // for (unsigned i = 0; i < jets.size(); i++) {
904 // const PseudoJet & jet = jets[i];
905 // cout << i << " "
906 // << " " << jet.pt() << " " << jet.rap() << " " << jet.phi() << " " << jet.m()
907 // << " " << bge_ptr->rho(jet) << " " << bge_ptr->rho_m(jet)
908 // << " " << bge_ptr->sigma(jet) << " " << bge_ptr->sigma_m(jet) << endl;
909 // if (do_subtractor) {
910 // const PseudoJet & subjet = subjets[i];
911 // cout << i << "sub"
912 // << " " << subjet.pt() << " " << subjet.rap() << " " << subjet.phi() << " " << subjet.m()
913 // << " " << jet.area() << endl;
914 // }
915 // }
916 } else {
917 cout << " rho = " << rho
918 << ", sigma = " << sigma
919 << ", mean_area = " << mean_area
920 << ", empty_area = " << empty_area
921 << ", n_empty_jets = " << n_empty_jets
922 << endl;
923 }
924 if (bge_ptr != 0) delete bge_ptr;
925 }
926#endif
927 } // try
928 catch (Error &fjerr) {
929 cout << "Caught fastjet error, exiting gracefully" << endl;
930 exit(0);
931 }
932
933 } // irepeat
934 } // iev
935 // if we've instantiated a plugin, delete it
936 if (jet_def.strategy()==plugin_strategy){
937 delete jet_def.plugin();
938 }
939 // close any file that we've opened
940 if (istr != &cin) delete istr;
941 } //
942
943}
944
945
946
947
948//------ HELPER ROUTINES -----------------------------------------------
949/// print a single jet
950void print_jet (const PseudoJet & jet) {
951 unsigned int n_constituents = jet.constituents().size();
952 printf("%15.8f %15.8f %15.8f %8u\n",
953 jet.rap(), jet.phi(), jet.perp(), n_constituents);
954}
955
956
957//----------------------------------------------------------------------
958void print_jets(const vector<PseudoJet> & jets_in, bool show_constituents) {
959 vector<PseudoJet> jets;
960 if (ee_print) {
961 jets = sorted_by_E(jets_in);
962 for (unsigned int j = 0; j < jets.size(); j++) {
963 printf("%5u %15.8f %15.8f %15.8f %15.8f\n",
964 j,jets[j].px(),jets[j].py(),jets[j].pz(),jets[j].E());
965 if (show_constituents) {
966 vector<PseudoJet> const_jets = jets[j].constituents();
967 for (unsigned int k = 0; k < const_jets.size(); k++) {
968 printf(" jet%03u %15.8f %15.8f %15.8f %15.8f\n",j,const_jets[k].px(),
969 const_jets[k].py(),const_jets[k].pz(),const_jets[k].E());
970 }
971 cout << "\n\n";
972 }
973
974 }
975 } else {
976 jets = sorted_by_pt(jets_in);
977 for (unsigned int j = 0; j < jets.size(); j++) {
978 printf("%5u %15.8f %15.8f %15.8f",
979 j,jets[j].rap(),jets[j].phi(),jets[j].perp());
980 // also print out the scalar area and the perp component of the
981 // 4-vector (just enough to check a reasonable 4-vector?)
982#ifndef __FJCORE__
983 if (do_areas) printf(" %15.8f %15.8f", jets[j].area(),
984 jets[j].area_4vector().perp());
985 cout << "\n";
986#endif
987
988 if (show_constituents) {
989 vector<PseudoJet> const_jets = jets[j].constituents();
990 for (unsigned int k = 0; k < const_jets.size(); k++) {
991 printf(" jet%03u %15.8f %15.8f %15.8f %5d\n",j,const_jets[k].rap(),
992 const_jets[k].phi(),sqrt(const_jets[k].kt2()), const_jets[k].cluster_hist_index());
993 }
994 cout << "\n\n";
995 }
996 }
997 }
998
999 if (rootfile != "") {
1000 ofstream ostr(rootfile.c_str());
1001 ostr << "# " << cmdline_p->command_line() << endl;
1002 ostr << "# output for root" << endl;
1003 assert(jets.size() > 0);
1004 jets[0].validated_cs()->print_jets_for_root(jets,ostr);
1005 }
1006
1007}
1008
1009
1010//----- SUBJETS --------------------------------------------------------
1011/// a function that pretty prints a list of jets and the subjets for each
1012/// one
1013void print_jets_and_sub (const vector<PseudoJet> & jets, double dcut) {
1014
1015 // sort jets into increasing pt
1016 vector<PseudoJet> sorted_jets = sorted_by_pt(jets);
1017
1018 // label the columns
1019 printf("Printing jets and their subjets with subdcut = %10.5f\n",dcut);
1020 printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity",
1021 "phi", "pt", "n constituents");
1022
1023 // the kind of subjet finding used to test consistency among them
1024 SubType sub_type = subtype_internal;
1025 //SubType sub_type = subtype_newclust_dcut;
1026 //SubType sub_type = subtype_newclust_R;
1027
1028 // print out the details for each jet
1029 //for (unsigned int i = 0; i < sorted_jets.size(); i++) {
1030 for (vector<PseudoJet>::const_iterator jet = sorted_jets.begin();
1031 jet != sorted_jets.end(); jet++) {
1032 const JetDefinition & jet_def = jet->validated_cs()->jet_def();
1033
1034 // if jet pt^2 < dcut with kt alg, then some methods of
1035 // getting subjets will return nothing -- so skip the jet
1036 if (jet_def.jet_algorithm() == kt_algorithm
1037 && jet->perp2() < dcut) continue;
1038
1039
1040 printf("%5u ",(unsigned int) (jet - sorted_jets.begin()));
1041 print_jet(*jet);
1042 vector<PseudoJet> subjets;
1043 ClusterSequence * cspoint;
1044 if (sub_type == subtype_internal) {
1045 cspoint = 0;
1046 subjets = jet->exclusive_subjets(dcut);
1047 double ddnp1 = jet->exclusive_subdmerge_max(subjets.size());
1048 double ddn = jet->exclusive_subdmerge_max(subjets.size()-1);
1049 cout << " for " << ddnp1 << " < d < " << ddn << " one has " << endl;
1050 } else if (sub_type == subtype_newclust_dcut) {
1051 cspoint = new ClusterSequence(jet->constituents(), jet_def);
1052 subjets = cspoint->exclusive_jets(dcut);
1053 } else if (sub_type == subtype_newclust_R) {
1054 assert(jet_def.jet_algorithm() == cambridge_algorithm);
1055 JetDefinition subjd(jet_def.jet_algorithm(), jet_def.R()*sqrt(dcut));
1056 cspoint = new ClusterSequence(jet->constituents(), subjd);
1057 subjets = cspoint->inclusive_jets();
1058 } else {
1059 cerr << "unrecognized subtype for subjet finding" << endl;
1060 exit(-1);
1061 }
1062
1063 subjets = sorted_by_pt(subjets);
1064 for (unsigned int j = 0; j < subjets.size(); j++) {
1065 printf(" -sub-%02u ",j);
1066 print_jet(subjets[j]);
1067 }
1068
1069 if (cspoint != 0) delete cspoint;
1070
1071 //ClusterSequence subseq(clust_seq->constituents(sorted_jets[i]),
1072 // JetDefinition(cambridge_algorithm, 0.4));
1073 //vector<PseudoJet> subjets = sorted_by_pt(subseq.inclusive_jets());
1074 //for (unsigned int j = 0; j < subjets.size(); j++) {
1075 // printf(" -sub-%02u ",j);
1076 // print_jet(subseq, subjets[j]);
1077 //}
1078 }
1079
1080}
1081
1082/// if abs(x)<precision/2, return 0
1083double make_safe_zero_truncation(double x, double precision){
1084 return std::abs(x)<0.5*precision ? 0.0 : x;
1085}
1086
1087#ifndef __FJCORE__
1088void print_jets_bkgd(const vector<PseudoJet> &jets,
1089 const vector<PseudoJet> &subtracted_jets,
1090 BackgroundEstimatorBase * bge_ptr,
1091 bool do_subtractor){
1092 printf("Printing jets, background information");
1093 if (do_subtractor)
1094 printf(" and subtracted jets\n");
1095 printf("%5s %15s %15s %15s %15s %15s %15s %15s %15s\n","jet #",
1096 "rapidity", "phi", "pt", "pt^2+m^2",
1097 "rho", "rho_m", "sigma", "sigma_m");
1098 if (do_subtractor)
1099 printf("%5s %15s %15s %15s %15s %15ss\n","jet #",
1100 "rapidity", "phi", "pt", "sqrt(pt^2+m^2)", "area");
1101
1102 for (unsigned i = 0; i < jets.size(); i++) {
1103 const PseudoJet & jet = jets[i];
1104 BackgroundEstimate estimate = bge_ptr->estimate(jet);
1105 // Note that the values of rho_m sometimes comes out as +- a very
1106 // small number and the format can produce either 0.00000000 or
1107 // -0.00000000. The call to "make_safe_zero_truncation" makes sure it is
1108 // printed wo the - sign in each case
1109 printf("%5u %15.8f %15.8f %15.8f %15.8f %15.8f %15.8f %15.8f %15.8f\n", i,
1110 jet.rap(), jet.phi(), jet.perp(), jet.mt(),
1111 estimate.rho(), make_safe_zero_truncation(estimate.rho_m(),1e-8),
1112 estimate.sigma(), estimate.sigma_m());
1113 if (do_subtractor) {
1114 const PseudoJet & subjet = subtracted_jets[i];
1115 printf("%5u %15.8f %15.8f %15.8f %15.8f %15.8f\n", i,
1116 subjet.rap(), subjet.phi(), subjet.perp(), subjet.mt(), jet.area());
1117 }
1118 }
1119}
1120#endif// __FJCORE__
1121
1122//----------------------------------------------------------------------
1123void signal_failed_comparison(int iev,
1124 const string & message,
1125 const vector<PseudoJet> & particles) {
1126 cout << "# failed comparison, reason is " << message << endl;
1127 cout << "# iev = " << iev << endl;
1128 cout << setprecision(16);
1129 for (unsigned i = 0; i < particles.size(); i++) {
1130 const PseudoJet & p = particles[i];
1131 cout << p.px() << " "
1132 << p.py() << " "
1133 << p.pz() << " "
1134 << p.E () << endl;
1135 }
1136 cout << "#END" << endl;
1137}
1138
1139//----------------------------------------------------------------------
1140void do_compare_strategy(int iev,
1141 const vector<PseudoJet> & particles,
1142 const JetDefinition & jet_def,
1143 const ClusterSequence & cs,
1144 int compare_strategy) {
1145
1146 // create a jet def with the reference comparison strategy
1147 JetDefinition jet_def_ref(jet_def.jet_algorithm(),
1148 jet_def.R(),
1149 jet_def.recombination_scheme(),
1150 Strategy(compare_strategy));
1151 // do the clustering
1152 ClusterSequence cs_ref(particles, jet_def_ref);
1153
1154 // now compare the outputs. At this stage just based on the clustering
1155 // sequence - get more sophisticated later...
1156 const vector<ClusterSequence::history_element> & history_in = cs.history();
1157 const vector<ClusterSequence::history_element> & history_ref = cs_ref.history();
1158
1159 if (history_in.size() != history_ref.size()) {
1160 signal_failed_comparison(iev, "history sizes do not match", particles);
1161 }
1162
1163 // now run over each clustering step to do the comparison
1164 for (unsigned i = cs.n_particles(); i < history_in.size(); i++) {
1165 bool fail_parents = (history_in[i].parent1 != history_ref[i].parent1 ||
1166 history_in[i].parent2 != history_ref[i].parent2);
1167 bool fail_dij = (history_in[i].dij != history_ref[i].dij);
1168 if (fail_parents || fail_dij) {
1169 ostringstream ostr;
1170 ostr << "at step " << i << ", ";
1171 if (fail_parents) ostr << "parents do not match, ";
1172 if (fail_dij) ostr << "dij does not match, ";
1173 ostr << "history in (p1, p2, dij) = "
1174 << history_in[i].parent1 << " " << history_in[i].parent2 << " " << history_in[i].dij;
1175 ostr << ", history ref (p1, p2, dij) = "
1176 << history_ref[i].parent1 << " " << history_ref[i].parent2 << " " << history_ref[i].dij;
1177 signal_failed_comparison(iev, ostr.str(), particles);
1178 break;
1179 }
1180 }
1181}
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.
static void print_banner()
This is the function that is automatically called during clustering to print the FastJet banner.
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 rho() const
returns rho, the median background density per unit area
virtual void set_rescaling_class(const FunctionOfPseudoJet< double > *rescaling_class)
Set a pointer to a class that calculates the rescaling factor as a function of the jet (position).
double mean_area() const
returns the area of the grid cells (all identical, but referred to as "mean" area for uniformity with...
void set_particles(const std::vector< PseudoJet > &particles)
tell the background estimator that it has a new event, composed of the specified particles.
double sigma() const
returns sigma, the background fluctuations per unit area; must be multipled by sqrt(area) to get fluc...
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
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
get rho, the median background density per unit area
void set_cluster_sequence(const ClusterSequenceAreaBase &csa)
(re)set the cluster sequence (with area support) to be used by future calls to rho() etc.
virtual void set_rescaling_class(const FunctionOfPseudoJet< double > *rescaling_class_in)
Set a pointer to a class that calculates the rescaling factor as a function of the jet (position).
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)
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.
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 ...
@ Best
automatic selection of the best (based on N), including the LazyTiled strategies that are new to FJ3....
@ 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
@ ee_genkt_algorithm
the e+e- genkt algorithm (R > 2 and p=1 gives ee_kt)
@ genkt_algorithm
like the k_t but with distance measures dij = min(kti^{2p},ktj^{2p}) Delta R_{ij}^2 / R^2 diB = 1/kti...
@ ee_kt_algorithm
the e+e- kt algorithm
@ cambridge_algorithm
the longitudinally invariant variant of the cambridge algorithm (aka Aachen algoithm).
@ antikt_algorithm
like the k_t but with distance measures dij = min(1/kti^2,1/ktj^2) Delta R_{ij}^2 / R^2 diB = 1/kti^2
@ kt_algorithm
the longitudinally invariant kt algorithm
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.