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