FastJet 3.0beta1
|
00001 //STARTHEADER 00002 // $Id: fastjet_timing_plugins.cc 2429 2011-07-20 17:50:33Z salam $ 00003 // 00004 // Copyright (c) 2005-2009, Matteo Cacciari, Gavin Salam and Gregory Soyez 00005 // 00006 //---------------------------------------------------------------------- 00007 // This file is part of FastJet. 00008 // 00009 // FastJet is free software; you can redistribute it and/or modify 00010 // it under the terms of the GNU General Public License as published by 00011 // the Free Software Foundation; either version 2 of the License, or 00012 // (at your option) any later version. 00013 // 00014 // The algorithms that underlie FastJet have required considerable 00015 // development and are described in hep-ph/0512210. If you use 00016 // FastJet as part of work towards a scientific publication, please 00017 // include a citation to the FastJet paper. 00018 // 00019 // FastJet is distributed in the hope that it will be useful, 00020 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00021 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00022 // GNU General Public License for more details. 00023 // 00024 // You should have received a copy of the GNU General Public License 00025 // along with FastJet; if not, write to the Free Software 00026 // Foundation, Inc.: 00027 // 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00028 //---------------------------------------------------------------------- 00029 //ENDHEADER 00030 00031 00032 //---------------------------------------------------------------------- 00033 /// fastjet_timing.cc: Program to help time and test the fastjet package 00034 /// 00035 /// It reads files containing multiple events in the format 00036 /// p1x p1y p1z E1 00037 /// p2x p2y p2z E2 00038 /// ... 00039 /// #END 00040 /// 00041 /// An example input file containing 10 events is included as 00042 /// data/Pythia-PtMin1000-LHC-10ev.dat 00043 /// 00044 /// Usage: 00045 /// fastjet_timing [-strategy NUMBER] [-repeat nrepeats] [-massive] \ 00046 /// [-combine nevents] [-r Rparameter] [-incl ptmin] [...] \ 00047 /// < data_file 00048 /// 00049 /// where the clustering can be repeated to aid timing and multiple 00050 /// events can be combined to get to larger multiplicities. Some options: 00051 /// 00052 /// Options for reading 00053 /// ------------------- 00054 /// 00055 /// -nev n number of events to run 00056 /// 00057 /// -combine n for combining multiple events from the data file in order 00058 /// to get a single high-multipicity event to run. 00059 /// 00060 /// -massless read in only the 3-momenta and deduce energies assuming 00061 /// that particles are massless 00062 /// 00063 /// -dense adds dense ghost coverage 00064 /// 00065 /// -repeat n repeats each event n times 00066 /// 00067 /// Output Options 00068 /// -------------- 00069 /// 00070 /// -incl ptmin output of all inclusive jets with pt > ptmin is obtained 00071 /// with the -incl option. 00072 /// 00073 /// -excld dcut output of all exclusive jets as obtained in a clustering 00074 /// with dcut 00075 /// 00076 /// -excly ycut output of all exclusive jets as obtained in a clustering 00077 /// with ycut 00078 /// 00079 /// -excln n output of clustering to n exclusive jets 00080 /// 00081 /// -ee-print print things as px,py,pz,E 00082 /// 00083 /// -get-all-dij print out all dij values 00084 /// -get-all-yij print out all yij values 00085 /// 00086 /// -const show jet constituents (works with excl jets) 00087 /// 00088 /// -write for writing out detailed clustering sequence (valuable 00089 /// for testing purposes) 00090 /// 00091 /// -unique_write writes out the sequence of dij's according to the 00092 /// "unique_history_order" (useful for verifying consistency 00093 /// between different clustering strategies). 00094 /// 00095 /// -root file sends output to file that can be read in with the script in 00096 /// root/ so as to show a lego-plot of the event 00097 /// 00098 /// -cones show extra info about internal steps for SISCone 00099 /// 00100 /// -area calculate areas. Additional options include 00101 /// -area:active 00102 /// -area:passive 00103 /// -area:explicit 00104 /// -area:voronoi Rfact 00105 /// -area:repeat nrepeat 00106 /// -ghost-area area 00107 /// -ghost-maxrap maxrap 00108 /// -area:fj2 place ghosts as in fj2 00109 /// 00110 /// -bkgd calculate the background density. Additional options include 00111 /// -bkgd:csab use the old ClusterSequenceAreaBase methods 00112 /// -bkgd:jetmedian use the new JetMedianBackgroundEstimator class 00113 /// -bkgd:fj2 force jetmedian to calculate sigma as in fj2 00114 /// -bkgd:gridmedian use GridMedianBackgroundEstimator with grid up to ghost_maxrap-ktR and grid spacing of 2ktR 00115 /// 00116 /// Algorithms 00117 /// ---------- 00118 /// -kt switch to the longitudinally invariant kt algorithm 00119 /// Note: this is the default one. 00120 /// 00121 /// -cam switch to the inclusive Cambridge/Aachen algorithm -- 00122 /// note that the option -excld dcut provides a clustering 00123 /// up to the dcut which is the minimum squared 00124 /// distance between any pair of jets. 00125 /// 00126 /// -antikt switch to the anti-kt clustering algorithm 00127 /// 00128 /// -genkt switch to the genkt algorithm 00129 /// you can provide the parameter of the alg as an argument to 00130 /// -genkt (1 by default) 00131 /// 00132 /// -eekt switch to the e+e- kt algorithm 00133 /// 00134 /// -eegenkt switch to the genkt algorithm 00135 /// you can provide the parameter of the alg as an argument to 00136 /// -ee_genkt (1 by default) 00137 /// 00138 /// plugins (don't delete this line) 00139 /// 00140 /// -pxcone switch to the PxCone jet algorithm 00141 /// 00142 /// -siscone switch to the SISCone jet algorithm (seedless cones) 00143 /// -sisconespheri switch to the Spherical SISCone jet algorithm (seedless cones) 00144 /// 00145 /// -midpoint switch to CDF's midpoint code 00146 /// -jetclu switch to CDF's jetclu code 00147 /// 00148 /// -d0runipre96cone switch to the D0RunIpre96Cone plugin 00149 /// -d0runicone switch to the D0RunICone plugin 00150 /// 00151 /// -d0runiicone switch to D0's run II midpoint cone 00152 /// 00153 /// -trackjet switch to the TrackJet plugin 00154 /// 00155 /// -atlascone switch to the ATLASCone plugin 00156 /// 00157 /// -eecambridge switch to the EECambridge plugin 00158 /// 00159 /// -jade switch to the Jade plugin 00160 /// 00161 /// -cmsiterativecone switch to the CMSIterativeCone plugin 00162 /// 00163 /// -gridjet switch to the GridJet plugin 00164 /// 00165 /// end of plugins (don't delete this line) 00166 /// 00167 /// 00168 /// Options for running algs 00169 /// ------------------------ 00170 /// 00171 /// -r sets the radius of the jet algorithm (default = 1.0) 00172 /// 00173 /// -overlap | -f sets the overlap fraction in cone algs with split-merge 00174 /// 00175 /// -seed sets the seed threshold 00176 /// 00177 /// -strategy N indicate stratgey from the enum fastjet::Strategy (see 00178 /// fastjet/JetDefinition.hh). 00179 /// 00180 00181 00182 #include "fastjet/ClusterSequenceArea.hh" 00183 #include "fastjet/tools/JetMedianBackgroundEstimator.hh" 00184 #include "fastjet/tools/GridMedianBackgroundEstimator.hh" 00185 #include<iostream> 00186 #include<sstream> 00187 #include<fstream> 00188 #include<valarray> 00189 #include<vector> 00190 #include <cstdlib> 00191 #include<cstddef> // for size_t 00192 #include "CmdLine.hh" 00193 00194 // get info on how fastjet was configured 00195 #include "fastjet/config.h" 00196 00197 // include the installed plugins (don't delete this line) 00198 #ifdef FASTJET_ENABLE_PLUGIN_SISCONE 00199 #include "fastjet/SISConePlugin.hh" 00200 #include "fastjet/SISConeSphericalPlugin.hh" 00201 #endif 00202 #ifdef FASTJET_ENABLE_PLUGIN_CDFCONES 00203 #include "fastjet/CDFMidPointPlugin.hh" 00204 #include "fastjet/CDFJetCluPlugin.hh" 00205 #endif 00206 #ifdef FASTJET_ENABLE_PLUGIN_PXCONE 00207 #include "fastjet/PxConePlugin.hh" 00208 #endif 00209 #ifdef FASTJET_ENABLE_PLUGIN_D0RUNIICONE 00210 #include "fastjet/D0RunIIConePlugin.hh" 00211 #endif 00212 #ifdef FASTJET_ENABLE_PLUGIN_TRACKJET 00213 #include "fastjet/TrackJetPlugin.hh" 00214 #endif 00215 #ifdef FASTJET_ENABLE_PLUGIN_ATLASCONE 00216 #include "fastjet/ATLASConePlugin.hh" 00217 #endif 00218 #ifdef FASTJET_ENABLE_PLUGIN_EECAMBRIDGE 00219 #include "fastjet/EECambridgePlugin.hh" 00220 #endif 00221 #ifdef FASTJET_ENABLE_PLUGIN_JADE 00222 #include "fastjet/JadePlugin.hh" 00223 #endif 00224 #ifdef FASTJET_ENABLE_PLUGIN_CMSITERATIVECONE 00225 #include "fastjet/CMSIterativeConePlugin.hh" 00226 #endif 00227 #ifdef FASTJET_ENABLE_PLUGIN_D0RUNICONE 00228 #include "fastjet/D0RunIpre96ConePlugin.hh" 00229 #include "fastjet/D0RunIConePlugin.hh" 00230 #endif 00231 #ifdef FASTJET_ENABLE_PLUGIN_GRIDJET 00232 #include "fastjet/GridJetPlugin.hh" 00233 #endif 00234 // end of installed plugins inclusion (don't delete this line) 00235 00236 using namespace std; 00237 00238 // to avoid excessive typing, define an abbreviation for the 00239 // fastjet namespace 00240 namespace fj = fastjet; 00241 00242 inline double pow2(const double x) {return x*x;} 00243 00244 // pretty print the jets and their subjets 00245 void print_jets_and_sub (const vector<fj::PseudoJet> & jets, double dcut); 00246 00247 string rootfile; 00248 CmdLine * cmdline_p; 00249 00250 bool do_areas; 00251 00252 /// sort and pretty print jets, with exact behaviour depending on 00253 /// whether ee_print is true or not 00254 bool ee_print = false; 00255 void print_jets(const vector<fj::PseudoJet> & jets, bool show_const = false); 00256 00257 void is_unavailable(const string & algname) { 00258 cerr << algname << " requested, but not available for this compilation" << endl; 00259 exit(0); 00260 } 00261 00262 00263 /// a program to test and time a range of algorithms as implemented or 00264 /// wrapped in fastjet 00265 int main (int argc, char ** argv) { 00266 00267 CmdLine cmdline(argc,argv); 00268 cmdline_p = &cmdline; 00269 // allow the use to specify the fj::Strategy either through the 00270 // -clever or the -strategy options (both will take numerical 00271 // values); the latter will override the former. 00272 fj::Strategy strategy = fj::Strategy(cmdline.int_val("-strategy", 00273 cmdline.int_val("-clever", fj::Best))); 00274 int repeat = cmdline.int_val("-repeat",1); 00275 int combine = cmdline.int_val("-combine",1); 00276 bool write = cmdline.present("-write"); 00277 bool unique_write = cmdline.present("-unique_write"); 00278 bool hydjet = cmdline.present("-hydjet"); 00279 double ktR = cmdline.double_val("-r",1.0); 00280 ktR = cmdline.double_val("-R",ktR); // allow -r and -R 00281 double inclkt = cmdline.double_val("-incl",-1.0); 00282 int excln = cmdline.int_val ("-excln",-1); 00283 double excld = cmdline.double_val("-excld",-1.0); 00284 double excly = cmdline.double_val("-excly",-1.0); 00285 ee_print = cmdline.present("-ee-print"); 00286 bool get_all_dij = cmdline.present("-get-all-dij"); 00287 bool get_all_yij = cmdline.present("-get-all-yij"); 00288 double subdcut = cmdline.double_val("-subdcut",-1.0); 00289 double etamax = cmdline.double_val("-etamax",1.0e305); 00290 bool show_constituents = cmdline.present("-const"); 00291 bool massless = cmdline.present("-massless"); 00292 int nev = cmdline.int_val("-nev",1); 00293 bool add_dense_coverage = cmdline.present("-dense"); 00294 double ghost_maxrap = cmdline.value("-ghost-maxrap",5.0); 00295 do_areas = cmdline.present("-area"); 00296 fj::AreaDefinition area_def; 00297 if (do_areas) { 00298 assert(!write); // it's incompatible 00299 fj::GhostedAreaSpec ghost_spec(ghost_maxrap, 00300 cmdline.value("-area:repeat", 1), 00301 cmdline.value("-ghost-area", 0.01)); 00302 if (cmdline.present("-area:fj2")) ghost_spec.set_fj2_placement(true); 00303 if (cmdline.present("-area:explicit")) { 00304 area_def = fj::AreaDefinition(fj::active_area_explicit_ghosts, ghost_spec); 00305 } else if (cmdline.present("-area:passive")) { 00306 area_def = fj::AreaDefinition(fj::passive_area, ghost_spec); 00307 } else if (cmdline.present("-area:voronoi")) { 00308 double Rfact = cmdline.value<double>("-area:voronoi"); 00309 area_def = fj::AreaDefinition(fj::voronoi_area, 00310 fj::VoronoiAreaSpec(Rfact)); 00311 } else { 00312 cmdline.present("-area:active"); // allow, but do not require, arg 00313 area_def = fj::AreaDefinition(fj::active_area, ghost_spec); 00314 } 00315 } 00316 bool do_bkgd = cmdline.present("-bkgd"); // background estimation 00317 bool do_bkgd_csab = false, do_bkgd_jetmedian = false, do_bkgd_fj2 = false; 00318 bool do_bkgd_gridmedian = false; 00319 fj::Selector bkgd_range; 00320 if (do_bkgd) { 00321 bkgd_range = fj::SelectorAbsRapMax(ghost_maxrap - ktR); 00322 if (cmdline.present("-bkgd:csab")) {do_bkgd_csab = true;} 00323 else if (cmdline.present("-bkgd:jetmedian")) {do_bkgd_jetmedian = true; 00324 do_bkgd_fj2 = cmdline.present("-bkgd:fj2"); 00325 } else if (cmdline.present("-bkgd:gridmedian")) {do_bkgd_gridmedian = true; 00326 } else { 00327 throw fj::Error("with the -bkgd option, some particular background must be specified (csab or jetmedian)"); 00328 } 00329 assert(do_areas || do_bkgd_gridmedian); 00330 } 00331 00332 bool show_cones = cmdline.present("-cones"); // only works for siscone 00333 00334 // for cone algorithms 00335 // allow -f and -overlap 00336 double overlap_threshold = cmdline.double_val("-overlap",0.5); 00337 overlap_threshold = cmdline.double_val("-f",overlap_threshold); 00338 double seed_threshold = cmdline.double_val("-seed",1.0); 00339 00340 // for ee algorithms, allow to specify ycut 00341 double ycut = cmdline.double_val("-ycut",0.08); 00342 00343 // for printing jets to a file for reading by root 00344 rootfile = cmdline.value<string>("-root",""); 00345 00346 // out default scheme is the E_scheme 00347 fj::RecombinationScheme scheme = fj::E_scheme; 00348 00349 // The following option causes the Cambridge algo to be used. 00350 // Note that currently the only output that works sensibly here is 00351 // "-incl 0" 00352 fj::JetDefinition jet_def; 00353 if (cmdline.present("-cam") || cmdline.present("-CA")) { 00354 jet_def = fj::JetDefinition(fj::cambridge_algorithm, ktR, scheme, strategy); 00355 } else if (cmdline.present("-antikt")) { 00356 jet_def = fj::JetDefinition(fj::antikt_algorithm, ktR, scheme, strategy); 00357 } else if (cmdline.present("-genkt")) { 00358 double p = cmdline.value<double>("-genkt"); 00359 jet_def = fj::JetDefinition(fj::genkt_algorithm, ktR, p, scheme, strategy); 00360 } else if (cmdline.present("-eekt")) { 00361 jet_def = fj::JetDefinition(fj::ee_kt_algorithm); 00362 } else if (cmdline.present("-eegenkt")) { 00363 double p = cmdline.value<double>("-eegenkt"); 00364 jet_def = fj::JetDefinition(fj::ee_genkt_algorithm, ktR, p, scheme, strategy); 00365 00366 // checking if one asks to run a plugin (don't delete this line) 00367 } else if (cmdline.present("-midpoint")) { 00368 #ifdef FASTJET_ENABLE_PLUGIN_CDFCONES 00369 typedef fj::CDFMidPointPlugin MPPlug; // for brevity 00370 double cone_area_fraction = 1.0; 00371 int max_pair_size = 2; 00372 int max_iterations = 100; 00373 MPPlug::SplitMergeScale sm_scale = MPPlug::SM_pt; 00374 if (cmdline.present("-sm-pttilde")) sm_scale = MPPlug::SM_pttilde; 00375 if (cmdline.present("-sm-pt")) sm_scale = MPPlug::SM_pt; // default 00376 if (cmdline.present("-sm-mt")) sm_scale = MPPlug::SM_mt; 00377 if (cmdline.present("-sm-Et")) sm_scale = MPPlug::SM_Et; 00378 jet_def = fj::JetDefinition( new fj::CDFMidPointPlugin ( 00379 seed_threshold, ktR, 00380 cone_area_fraction, max_pair_size, 00381 max_iterations, overlap_threshold, 00382 sm_scale)); 00383 #else // FASTJET_ENABLE_PLUGIN_CDFCONES 00384 is_unavailable("midpoint"); 00385 #endif // FASTJET_ENABLE_PLUGIN_CDFCONES 00386 } else if (cmdline.present("-pxcone")) { 00387 #ifdef FASTJET_ENABLE_PLUGIN_PXCONE 00388 double min_jet_energy = 5.0; 00389 jet_def = fj::JetDefinition( new fj::PxConePlugin ( 00390 ktR, min_jet_energy, 00391 overlap_threshold)); 00392 #else // FASTJET_ENABLE_PLUGIN_PXCONE 00393 is_unavailable("pxcone"); 00394 #endif // FASTJET_ENABLE_PLUGIN_PXCONE 00395 } else if (cmdline.present("-jetclu")) { 00396 #ifdef FASTJET_ENABLE_PLUGIN_CDFCONES 00397 jet_def = fj::JetDefinition( new fj::CDFJetCluPlugin ( 00398 ktR, overlap_threshold, seed_threshold)); 00399 #else // FASTJET_ENABLE_PLUGIN_CDFCONES 00400 is_unavailable("pxcone"); 00401 #endif // FASTJET_ENABLE_PLUGIN_CDFCONES 00402 } else if (cmdline.present("-siscone") || cmdline.present("-sisconespheri")) { 00403 #ifdef FASTJET_ENABLE_PLUGIN_SISCONE 00404 typedef fj::SISConePlugin SISPlug; // for brevity 00405 int npass = cmdline.value("-npass",0); 00406 if (cmdline.present("-siscone")) { 00407 double sisptmin = cmdline.value("-sisptmin",0.0); 00408 SISPlug * plugin = new SISPlug (ktR, overlap_threshold,npass,sisptmin); 00409 if (cmdline.present("-sm-pt")) plugin->set_split_merge_scale(SISPlug::SM_pt); 00410 if (cmdline.present("-sm-mt")) plugin->set_split_merge_scale(SISPlug::SM_mt); 00411 if (cmdline.present("-sm-Et")) plugin->set_split_merge_scale(SISPlug::SM_Et); 00412 if (cmdline.present("-sm-pttilde")) plugin->set_split_merge_scale(SISPlug::SM_pttilde); 00413 // cause it to use the jet-definition's own recombiner 00414 plugin->set_use_jet_def_recombiner(true); 00415 jet_def = fj::JetDefinition(plugin); 00416 } else { 00417 double sisEmin = cmdline.value("-sisEmin",0.0); 00418 fj::SISConeSphericalPlugin * plugin = 00419 new fj::SISConeSphericalPlugin(ktR, overlap_threshold,npass,sisEmin); 00420 if (cmdline.present("-ghost-sep")) { 00421 plugin->set_ghost_separation_scale(cmdline.value<double>("-ghost-sep")); 00422 } 00423 jet_def = fj::JetDefinition(plugin); 00424 } 00425 #else // FASTJET_ENABLE_PLUGIN_SISCONE 00426 is_unavailable("siscone"); 00427 #endif // FASTJET_ENABLE_PLUGIN_SISCONE 00428 } else if (cmdline.present("-d0runiicone")) { 00429 #ifdef FASTJET_ENABLE_PLUGIN_D0RUNIICONE 00430 double min_jet_Et = 6.0; // was 8 GeV in earlier work 00431 jet_def = fj::JetDefinition(new fj::D0RunIIConePlugin(ktR,min_jet_Et)); 00432 #else // FASTJET_ENABLE_PLUGIN_D0RUNIICONE 00433 is_unavailable("D0RunIICone"); 00434 #endif // FASTJET_ENABLE_PLUGIN_D0RUNIICONE 00435 } else if (cmdline.present("-trackjet")) { 00436 #ifdef FASTJET_ENABLE_PLUGIN_TRACKJET 00437 jet_def = fj::JetDefinition(new fj::TrackJetPlugin(ktR)); 00438 #else // FASTJET_ENABLE_PLUGIN_TRACKJET 00439 is_unavailable("TrackJet"); 00440 #endif // FASTJET_ENABLE_PLUGIN_TRACKJET 00441 } else if (cmdline.present("-atlascone")) { 00442 #ifdef FASTJET_ENABLE_PLUGIN_ATLASCONE 00443 jet_def = fj::JetDefinition(new fj::ATLASConePlugin(ktR)); 00444 #else // FASTJET_ENABLE_PLUGIN_ATLASCONE 00445 is_unavailable("ATLASCone"); 00446 #endif // FASTJET_ENABLE_PLUGIN_ATLASCONE 00447 } else if (cmdline.present("-eecambridge")) { 00448 #ifdef FASTJET_ENABLE_PLUGIN_EECAMBRIDGE 00449 jet_def = fj::JetDefinition(new fj::EECambridgePlugin(ycut)); 00450 #else // FASTJET_ENABLE_PLUGIN_EECAMBRIDGE 00451 is_unavailable("EECambridge"); 00452 #endif // FASTJET_ENABLE_PLUGIN_EECAMBRIDGE 00453 } else if (cmdline.present("-jade")) { 00454 #ifdef FASTJET_ENABLE_PLUGIN_JADE 00455 jet_def = fj::JetDefinition(new fj::JadePlugin()); 00456 #else // FASTJET_ENABLE_PLUGIN_JADE 00457 is_unavailable("Jade"); 00458 #endif // FASTJET_ENABLE_PLUGIN_JADE 00459 } else if (cmdline.present("-cmsiterativecone")) { 00460 #ifdef FASTJET_ENABLE_PLUGIN_CMSITERATIVECONE 00461 jet_def = fj::JetDefinition(new fj::CMSIterativeConePlugin(ktR,seed_threshold)); 00462 #else // FASTJET_ENABLE_PLUGIN_CMSITERATIVECONE 00463 is_unavailable("CMSIterativeCone"); 00464 #endif // FASTJET_ENABLE_PLUGIN_CMSITERATIVECONE 00465 } else if (cmdline.present("-d0runipre96cone")) { 00466 #ifdef FASTJET_ENABLE_PLUGIN_D0RUNICONE 00467 jet_def = fj::JetDefinition(new fj::D0RunIpre96ConePlugin(ktR, seed_threshold, overlap_threshold)); 00468 #else // FASTJET_ENABLE_PLUGIN_D0RUNICONE 00469 is_unavailable("D0RunICone"); 00470 #endif // FASTJET_ENABLE_PLUGIN_D0RUNICONE 00471 } else if (cmdline.present("-d0runicone")) { 00472 #ifdef FASTJET_ENABLE_PLUGIN_D0RUNICONE 00473 jet_def = fj::JetDefinition(new fj::D0RunIConePlugin(ktR, seed_threshold, overlap_threshold)); 00474 #else // FASTJET_ENABLE_PLUGIN_D0RUNICONE 00475 is_unavailable("D0RunICone"); 00476 #endif // FASTJET_ENABLE_PLUGIN_D0RUNICONE 00477 } else if (cmdline.present("-gridjet")) { 00478 #ifdef FASTJET_ENABLE_PLUGIN_GRIDJET 00479 // we want a grid_ymax of 5.0, but when using R=0.4 (i.e. grid 00480 // spacing of 0.8), this leads to 12.5 grid cells; depending on 00481 // whether this is 12.499999999999 or 12.5000000....1 this gets 00482 // converted either to 12 or 13, making the results sensitive to 00483 // rounding errors. 00484 // 00485 // Instead we therefore take 4.9999999999, which avoids this problem. 00486 double grid_ymax = 4.9999999999; 00487 jet_def = fj::JetDefinition(new fj::GridJetPlugin(grid_ymax, ktR*2.0)); 00488 #else // FASTJET_ENABLE_PLUGIN_GRIDJET 00489 is_unavailable("GridJet"); 00490 #endif // FASTJET_ENABLE_PLUGIN_GRIDJET 00491 // end of checking if one asks to run a plugin (don't delete this line) 00492 } else { 00493 cmdline.present("-kt"); // kt is default, but allow user to specify it too [and ignore return value!] 00494 jet_def = fj::JetDefinition(fj::kt_algorithm, ktR, strategy); 00495 } 00496 00497 00498 00499 if (!cmdline.all_options_used()) {cerr << 00500 "Error: some options were not recognized"<<endl; 00501 exit(-1);} 00502 00503 00504 for (int iev = 0; iev < nev; iev++) { 00505 vector<fj::PseudoJet> jets; 00506 vector<fj::PseudoJet> particles; 00507 string line; 00508 int ndone = 0; 00509 while (getline(cin, line)) { 00510 //cout << line<<endl; 00511 istringstream linestream(line); 00512 if (line == "#END") { 00513 ndone += 1; 00514 if (ndone == combine) {break;} 00515 } 00516 if (line.substr(0,1) == "#") {continue;} 00517 valarray<double> fourvec(4); 00518 if (hydjet) { 00519 // special reading from hydjet.txt event record (though actually 00520 // this is supposed to be a standard pythia event record, so 00521 // being able to read from it is perhaps not so bad an idea...) 00522 int ii, istat,id,m1,m2,d1,d2; 00523 double mass; 00524 linestream >> ii>> istat >> id >> m1 >> m2 >> d1 >> d2 00525 >> fourvec[0] >> fourvec[1] >> fourvec[2] >> mass; 00526 // current file contains mass of particle as 4th entry 00527 if (istat == 1) { 00528 fourvec[3] = sqrt(+pow2(fourvec[0])+pow2(fourvec[1]) 00529 +pow2(fourvec[2])+pow2(mass)); 00530 } 00531 } else { 00532 if (massless) { 00533 linestream >> fourvec[0] >> fourvec[1] >> fourvec[2]; 00534 fourvec[3] = sqrt(pow2(fourvec[0])+pow2(fourvec[1])+pow2(fourvec[2]));} 00535 else { 00536 linestream >> fourvec[0] >> fourvec[1] >> fourvec[2] >> fourvec[3]; 00537 } 00538 } 00539 fj::PseudoJet psjet(fourvec); 00540 if (abs(psjet.rap() < etamax)) {particles.push_back(psjet);} 00541 } 00542 00543 // add a fake underlying event which is very soft, uniformly distributed 00544 // in eta,phi so as to allow one to reconstruct the area that is associated 00545 // with each jet. 00546 if (add_dense_coverage) { 00547 fj::GhostedAreaSpec ghosted_area_spec(ghost_maxrap); 00548 //fj::GhostedAreaSpec ghosted_area_spec(-2.0,4.0); // asymmetric range 00549 // for plots, reduce the scatter default of 1, to avoid "holes" 00550 // in the subsequent calorimeter view 00551 ghosted_area_spec.set_grid_scatter(0.5); 00552 ghosted_area_spec.add_ghosts(particles); 00553 //----- old code ------------------ 00554 // srand(2); 00555 // int nphi = 60; 00556 // int neta = 100; 00557 // double kt = 1e-1; 00558 // for (int iphi = 0; iphi<nphi; iphi++) { 00559 // for (int ieta = -neta; ieta<neta+1; ieta++) { 00560 // double phi = (iphi+0.5) * (fj::twopi/nphi) + rand()*0.001/RAND_MAX; 00561 // double eta = ieta * (10.0/neta) + rand()*0.001/RAND_MAX; 00562 // kt = 1e-20*(1+rand()*0.1/RAND_MAX); 00563 // double pminus = kt*exp(-eta); 00564 // double pplus = kt*exp(+eta); 00565 // double px = kt*sin(phi); 00566 // double py = kt*cos(phi); 00567 // //cout << kt<<" "<<eta<<" "<<phi<<"\n"; 00568 // fj::PseudoJet mom(px,py,0.5*(pplus-pminus),0.5*(pplus+pminus)); 00569 // particles.push_back(mom); 00570 // } 00571 // } 00572 } 00573 00574 for (int irepeat = 0; irepeat < repeat ; irepeat++) { 00575 int nparticles = particles.size(); 00576 try { 00577 auto_ptr<fj::ClusterSequence> clust_seq; 00578 if (do_areas) { 00579 clust_seq.reset(new fj::ClusterSequenceArea(particles,jet_def,area_def)); 00580 } else { 00581 clust_seq.reset(new fj::ClusterSequence(particles,jet_def,write)); 00582 } 00583 if (irepeat != 0) {continue;} 00584 cout << "iev "<<iev<< ": number of particles = "<< nparticles << endl; 00585 cout << "strategy used = "<< clust_seq->strategy_string()<< endl; 00586 if (iev == 0) cout << "Jet Definition: " << jet_def.description() << " (" << fj::fastjet_version_string() << ")" << endl; 00587 if (do_areas && iev == 0) cout << "Area definition: " << area_def.description() << endl; 00588 00589 // now provide some nice output... 00590 if (inclkt >= 0.0) { 00591 vector<fj::PseudoJet> jets = sorted_by_pt(clust_seq->inclusive_jets(inclkt)); 00592 print_jets(jets, show_constituents); 00593 00594 } 00595 00596 if (excln > 0) { 00597 cout << "Printing "<<excln<<" exclusive jets\n"; 00598 print_jets(clust_seq->exclusive_jets(excln), show_constituents); 00599 } 00600 00601 if (excld > 0.0) { 00602 cout << "Printing exclusive jets for d = "<<excld<<"\n"; 00603 print_jets(clust_seq->exclusive_jets(excld), show_constituents); 00604 } 00605 00606 if (excly > 0.0) { 00607 cout << "Printing exclusive jets for ycut = "<<excly<<"\n"; 00608 print_jets(clust_seq->exclusive_jets_ycut(excly), show_constituents); 00609 } 00610 00611 if (get_all_dij) { 00612 for (int i = nparticles-1; i >= 0; i--) { 00613 printf("d for n = %4d -> %4d is %14.5e\n", i+1, i, clust_seq->exclusive_dmerge(i)); 00614 } 00615 } 00616 if (get_all_yij) { 00617 for (int i = nparticles-1; i >= 0; i--) { 00618 printf("y for n = %4d -> %4d is %14.5e\n", i+1, i, clust_seq->exclusive_ymerge(i)); 00619 } 00620 } 00621 00622 // have the option of printing out the subjets (at scale dcut) of 00623 // each inclusive jet 00624 if (subdcut >= 0.0) { 00625 print_jets_and_sub(clust_seq->inclusive_jets(), subdcut); 00626 } 00627 00628 // useful for testing that recombination sequences are unique 00629 if (unique_write) { 00630 vector<int> unique_history = clust_seq->unique_history_order(); 00631 // construct the inverse of the above mapping 00632 vector<int> inv_unique_history(clust_seq->history().size()); 00633 for (unsigned int i = 0; i < unique_history.size(); i++) { 00634 inv_unique_history[unique_history[i]] = i;} 00635 00636 for (unsigned int i = 0; i < unique_history.size(); i++) { 00637 fj::ClusterSequence::history_element el = 00638 clust_seq->history()[unique_history[i]]; 00639 int uhp1 = el.parent1>=0 ? inv_unique_history[el.parent1] : el.parent1; 00640 int uhp2 = el.parent2>=0 ? inv_unique_history[el.parent2] : el.parent2; 00641 printf("%7d u %15.8e %7d u %7d u\n",i,el.dij,uhp1, uhp2); 00642 } 00643 } 00644 00645 00646 #ifdef FASTJET_ENABLE_PLUGIN_SISCONE 00647 // provide some complementary information for SISCone 00648 if (show_cones) { 00649 const fj::SISConeExtras * extras = 00650 dynamic_cast<const fj::SISConeExtras *>(clust_seq->extras()); 00651 cout << "most ambiguous split (difference in squared dist) = " 00652 << extras->most_ambiguous_split() << endl; 00653 vector<fastjet::PseudoJet> stable_cones(extras->stable_cones()); 00654 stable_cones = sorted_by_rapidity(stable_cones); 00655 for (unsigned int i = 0; i < stable_cones.size(); i++) { 00656 //if (stable_cones[i].phi() < 5.0 && stable_cones[i].phi() > 4.0) { 00657 printf("%5u %15.8f %15.8f %15.8f\n", 00658 i,stable_cones[i].rap(),stable_cones[i].phi(), 00659 stable_cones[i].perp() ); 00660 //} 00661 } 00662 00663 // also show passes for jets 00664 vector<fj::PseudoJet> sisjets = clust_seq->inclusive_jets(); 00665 printf("\n%15s %15s %15s %12s %8s %8s\n","rap","phi","pt","user-index","pass","nconst"); 00666 for (unsigned i = 0; i < sisjets.size(); i++) { 00667 printf("%15.8f %15.8f %15.8f %12d %8d %8d\n", 00668 sisjets[i].rap(), sisjets[i].phi(), sisjets[i].perp(), 00669 sisjets[i].user_index(), extras->pass(sisjets[i]), 00670 clust_seq->constituents(sisjets[i]).size() 00671 ); 00672 00673 } 00674 } 00675 #endif // FASTJET_ENABLE_PLUGIN_SISCONE 00676 00677 if (do_bkgd) { 00678 double rho, sigma, mean_area, empty_area, n_empty_jets; 00679 fj::ClusterSequenceAreaBase * csab = 00680 dynamic_cast<fj::ClusterSequenceAreaBase *>(clust_seq.get()); 00681 if (do_bkgd_csab) { 00682 csab->get_median_rho_and_sigma(bkgd_range, true, rho, sigma, mean_area); 00683 empty_area = csab->empty_area(bkgd_range); 00684 n_empty_jets = csab->n_empty_jets(bkgd_range); 00685 } else if (do_bkgd_jetmedian) { 00686 fj::JetMedianBackgroundEstimator bge(bkgd_range); 00687 bge.set_provide_fj2_sigma(do_bkgd_fj2); 00688 bge.set_cluster_sequence(*csab); 00689 rho = bge.rho(); 00690 sigma = bge.sigma(); 00691 mean_area = bge.mean_area(); 00692 empty_area = bge.empty_area(); 00693 n_empty_jets = bge.n_empty_jets(); 00694 } else { 00695 assert(do_bkgd_gridmedian); 00696 double rapmin, rapmax; 00697 bkgd_range.get_rapidity_extent(rapmin, rapmax); 00698 fj::GridMedianBackgroundEstimator bge(rapmax, 2*ktR); 00699 bge.set_particles(particles); 00700 rho = bge.rho(); 00701 sigma = bge.sigma(); 00702 mean_area = bge.mean_area(); 00703 empty_area = 0; 00704 n_empty_jets = 0; 00705 } 00706 cout << " rho = " << rho 00707 << ", sigma = " << sigma 00708 << ", mean_area = " << mean_area 00709 << ", empty_area = " << empty_area 00710 << ", n_empty_jets = " << n_empty_jets 00711 << endl; 00712 } 00713 } // try 00714 catch (fastjet::Error fjerr) { 00715 cout << "Caught fastjet error, exiting gracefully" << endl; 00716 exit(0); 00717 } 00718 00719 } // irepeat 00720 } // iev 00721 00722 // if we've instantiated a plugin, delete it 00723 if (jet_def.strategy()==fj::plugin_strategy){ 00724 delete jet_def.plugin(); 00725 } 00726 } 00727 00728 00729 00730 00731 //------ HELPER ROUTINES ----------------------------------------------- 00732 /// print a single jet 00733 void print_jet (const fj::PseudoJet & jet) { 00734 int n_constituents = jet.constituents().size(); 00735 printf("%15.8f %15.8f %15.8f %8u\n", 00736 jet.rap(), jet.phi(), jet.perp(), n_constituents); 00737 } 00738 00739 00740 //---------------------------------------------------------------------- 00741 void print_jets(const vector<fj::PseudoJet> & jets_in, bool show_constituents) { 00742 vector<fj::PseudoJet> jets; 00743 if (ee_print) { 00744 jets = sorted_by_E(jets_in); 00745 for (size_t j = 0; j < jets.size(); j++) { 00746 printf("%5u %15.8f %15.8f %15.8f %15.8f\n", 00747 j,jets[j].px(),jets[j].py(),jets[j].pz(),jets[j].E()); 00748 if (show_constituents) { 00749 vector<fj::PseudoJet> const_jets = jets[j].constituents(); 00750 for (size_t k = 0; k < const_jets.size(); k++) { 00751 printf(" jet%03u %15.8f %15.8f %15.8f %15.8f\n",j,const_jets[k].px(), 00752 const_jets[k].py(),const_jets[k].pz(),const_jets[k].E()); 00753 } 00754 cout << "\n\n"; 00755 } 00756 00757 } 00758 } else { 00759 jets = sorted_by_pt(jets_in); 00760 for (size_t j = 0; j < jets.size(); j++) { 00761 printf("%5u %15.8f %15.8f %15.8f", 00762 j,jets[j].rap(),jets[j].phi(),jets[j].perp()); 00763 // also print out the scalar area and the perp component of the 00764 // 4-vector (just enough to check a reasonable 4-vector?) 00765 if (do_areas) printf(" %15.8f %15.8f", jets[j].area(), 00766 jets[j].area_4vector().perp()); 00767 cout << "\n"; 00768 00769 if (show_constituents) { 00770 vector<fj::PseudoJet> const_jets = jets[j].constituents(); 00771 for (size_t k = 0; k < const_jets.size(); k++) { 00772 printf(" jet%03u %15.8f %15.8f %15.8f %5d\n",j,const_jets[k].rap(), 00773 const_jets[k].phi(),sqrt(const_jets[k].kt2()), const_jets[k].cluster_hist_index()); 00774 } 00775 cout << "\n\n"; 00776 } 00777 } 00778 } 00779 00780 if (rootfile != "") { 00781 ofstream ostr(rootfile.c_str()); 00782 ostr << "# " << cmdline_p->command_line() << endl; 00783 ostr << "# output for root" << endl; 00784 assert(jets.size() > 0); 00785 jets[0].validated_cs()->print_jets_for_root(jets,ostr); 00786 } 00787 00788 } 00789 00790 00791 //----- SUBJETS -------------------------------------------------------- 00792 /// a function that pretty prints a list of jets and the subjets for each 00793 /// one 00794 void print_jets_and_sub (const vector<fj::PseudoJet> & jets, double dcut) { 00795 00796 // sort jets into increasing pt 00797 vector<fj::PseudoJet> sorted_jets = sorted_by_pt(jets); 00798 00799 // label the columns 00800 printf("Printing jets and their subjets with subdcut = %10.5f\n",dcut); 00801 printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity", 00802 "phi", "pt", "n constituents"); 00803 00804 // have various kinds of subjet finding, to test consistency among them 00805 enum SubType {internal, newclust_dcut, newclust_R}; 00806 SubType subtype = internal; 00807 //SubType subtype = newclust_dcut; 00808 //SubType subtype = newclust_R; 00809 00810 // print out the details for each jet 00811 //for (unsigned int i = 0; i < sorted_jets.size(); i++) { 00812 for (vector<fj::PseudoJet>::const_iterator jet = sorted_jets.begin(); 00813 jet != sorted_jets.end(); jet++) { 00814 const fj::JetDefinition & jet_def = jet->validated_cs()->jet_def(); 00815 00816 // if jet pt^2 < dcut with kt alg, then some methods of 00817 // getting subjets will return nothing -- so skip the jet 00818 if (jet_def.jet_algorithm() == fj::kt_algorithm 00819 && jet->perp2() < dcut) continue; 00820 00821 00822 printf("%5u ",jet - sorted_jets.begin()); 00823 print_jet(*jet); 00824 vector<fj::PseudoJet> subjets; 00825 fj::ClusterSequence * cspoint; 00826 if (subtype == internal) { 00827 cspoint = 0; 00828 subjets = jet->exclusive_subjets(dcut); 00829 double ddnp1 = jet->exclusive_subdmerge_max(subjets.size()); 00830 double ddn = jet->exclusive_subdmerge_max(subjets.size()-1); 00831 cout << " for " << ddnp1 << " < d < " << ddn << " one has " << endl; 00832 } else if (subtype == newclust_dcut) { 00833 cspoint = new fj::ClusterSequence(jet->constituents(), jet_def); 00834 subjets = cspoint->exclusive_jets(dcut); 00835 } else if (subtype == newclust_R) { 00836 assert(jet_def.jet_algorithm() == fj::cambridge_algorithm); 00837 fj::JetDefinition subjd(jet_def.jet_algorithm(), jet_def.R()*sqrt(dcut)); 00838 cspoint = new fj::ClusterSequence(jet->constituents(), subjd); 00839 subjets = cspoint->inclusive_jets(); 00840 } else { 00841 cerr << "unrecognized subtype for subjet finding" << endl; 00842 exit(-1); 00843 } 00844 00845 subjets = sorted_by_pt(subjets); 00846 for (unsigned int j = 0; j < subjets.size(); j++) { 00847 printf(" -sub-%02u ",j); 00848 print_jet(subjets[j]); 00849 } 00850 00851 if (cspoint != 0) delete cspoint; 00852 00853 //fj::ClusterSequence subseq(clust_seq->constituents(sorted_jets[i]), 00854 // fj::JetDefinition(fj::cambridge_algorithm, 0.4)); 00855 //vector<fj::PseudoJet> subjets = sorted_by_pt(subseq.inclusive_jets()); 00856 //for (unsigned int j = 0; j < subjets.size(); j++) { 00857 // printf(" -sub-%02u ",j); 00858 // print_jet(subseq, subjets[j]); 00859 //} 00860 } 00861 00862 } 00863