fastjet_timing_plugins.cc

00001 //STARTHEADER
00002 // $Id: fastjet_timing_plugins.cc 1909 2011-01-28 16:27:25Z soyez $
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 /// Algorithms
00101 /// ----------
00102 ///   -kt           switch to the longitudinally invariant kt algorithm
00103 ///                 Note: this is the default one.
00104 ///
00105 ///   -cam          switch to the inclusive Cambridge/Aachen algorithm --
00106 ///                 note that the option -excld dcut provides a clustering
00107 ///                 up to the dcut which is the minimum squared
00108 ///                 distance between any pair of jets.
00109 ///
00110 ///   -antikt       switch to the anti-kt clustering algorithm
00111 ///
00112 ///   -genkt        switch to the genkt algorithm
00113 ///                 you can provide the parameter of the alg as an argument to 
00114 ///                 -genkt (1 by default)
00115 ///                 
00116 ///   -eekt         switch to the e+e- kt algorithm
00117 ///
00118 ///   -eegenkt      switch to the genkt algorithm
00119 ///                 you can provide the parameter of the alg as an argument to 
00120 ///                 -ee_genkt (1 by default)
00121 ///                 
00122 /// plugins (don't delete this line)
00123 ///
00124 ///   -pxcone             switch to the PxCone jet algorithm
00125 /// 
00126 ///   -siscone            switch to the SISCone jet algorithm (seedless cones)
00127 ///   -sisconespheri      switch to the Spherical SISCone jet algorithm (seedless cones)
00128 ///
00129 ///   -midpoint           switch to CDF's midpoint code
00130 ///   -jetclu             switch to CDF's jetclu code
00131 ///
00132 ///   -d0runipre96cone    switch to the D0RunIpre96Cone plugin
00133 ///   -d0runicone         switch to the D0RunICone plugin
00134 ///
00135 ///   -d0runiicone        switch to D0's run II midpoint cone
00136 ///
00137 ///   -trackjet           switch to the TrackJet plugin
00138 ///
00139 ///   -atlascone          switch to the ATLASCone plugin
00140 ///
00141 ///   -eecambridge        switch to the EECambridge plugin
00142 ///
00143 ///   -jade               switch to the Jade plugin
00144 ///
00145 ///   -cmsiterativecone   switch to the CMSIterativeCone plugin
00146 ///
00147 ///  end of plugins (don't delete this line)
00148 ///
00149 ///
00150 /// Options for running algs
00151 /// ------------------------
00152 ///
00153 ///   -r            sets the radius of the jet algorithm (default = 1.0)
00154 ///
00155 ///   -overlap | -f sets the overlap fraction in cone algs with split-merge
00156 ///
00157 ///   -seed         sets the seed threshold
00158 ///
00159 ///   -strategy N   indicate stratgey from the enum fastjet::Strategy (see
00160 ///                 fastjet/JetDefinition.hh).
00161 ///
00162 
00163 
00164 #include "fastjet/PseudoJet.hh"
00165 #include "fastjet/ClusterSequence.hh"
00166 #include "fastjet/GhostedAreaSpec.hh"
00167 #include<iostream>
00168 #include<sstream>
00169 #include<fstream>
00170 #include<valarray>
00171 #include<vector>
00172 #include <cstdlib>
00173 #include<cstddef> // for size_t
00174 #include "CmdLine.hh"
00175 
00176 // get info on how fastjet was configured
00177 #include "fastjet/config.h"
00178 
00179 // include the installed plugins (don't delete this line)
00180 #ifdef FASTJET_ENABLE_PLUGIN_SISCONE
00181 #include "fastjet/SISConePlugin.hh"
00182 #include "fastjet/SISConeSphericalPlugin.hh"
00183 #endif
00184 #ifdef FASTJET_ENABLE_PLUGIN_CDFCONES
00185 #include "fastjet/CDFMidPointPlugin.hh"
00186 #include "fastjet/CDFJetCluPlugin.hh"
00187 #endif
00188 #ifdef FASTJET_ENABLE_PLUGIN_PXCONE
00189 #include "fastjet/PxConePlugin.hh"
00190 #endif
00191 #ifdef FASTJET_ENABLE_PLUGIN_D0RUNIICONE
00192 #include "fastjet/D0RunIIConePlugin.hh"
00193 #endif 
00194 #ifdef FASTJET_ENABLE_PLUGIN_TRACKJET
00195 #include "fastjet/TrackJetPlugin.hh"
00196 #endif
00197 #ifdef FASTJET_ENABLE_PLUGIN_ATLASCONE
00198 #include "fastjet/ATLASConePlugin.hh"
00199 #endif
00200 #ifdef FASTJET_ENABLE_PLUGIN_EECAMBRIDGE
00201 #include "fastjet/EECambridgePlugin.hh"
00202 #endif
00203 #ifdef FASTJET_ENABLE_PLUGIN_JADE
00204 #include "fastjet/JadePlugin.hh"
00205 #endif
00206 #ifdef FASTJET_ENABLE_PLUGIN_CMSITERATIVECONE
00207 #include "fastjet/CMSIterativeConePlugin.hh"
00208 #endif
00209 #ifdef FASTJET_ENABLE_PLUGIN_D0RUNICONE
00210 #include "fastjet/D0RunIpre96ConePlugin.hh"
00211 #include "fastjet/D0RunIConePlugin.hh"
00212 #endif
00213 // end of installed plugins inclusion (don't delete this line)
00214 
00215 using namespace std;
00216 
00217 // to avoid excessive typing, define an abbreviation for the 
00218 // fastjet namespace
00219 namespace fj = fastjet;
00220 
00221 inline double pow2(const double x) {return x*x;}
00222 
00223 // pretty print the jets and their subjets
00224 void print_jets_and_sub (const vector<fj::PseudoJet> & jets, double dcut);
00225 
00226 string rootfile;
00227 CmdLine * cmdline_p;
00228 
00229 /// sort and pretty print jets, with exact behaviour depending on 
00230 /// whether ee_print is true or not
00231 bool ee_print = false;
00232 void print_jets(const vector<fj::PseudoJet> & jets, bool show_const = false);
00233 
00234 void is_unavailable(const string & algname) {
00235   cerr << algname << " requested, but not available for this compilation";
00236   exit(0);
00237 }
00238 
00239 
00240 /// a program to test and time a range of algorithms as implemented or
00241 /// wrapped in fastjet
00242 int main (int argc, char ** argv) {
00243 
00244   CmdLine cmdline(argc,argv);
00245   cmdline_p = &cmdline;
00246   // allow the use to specify the fj::Strategy either through the
00247   // -clever or the -strategy options (both will take numerical
00248   // values); the latter will override the former.
00249   fj::Strategy  strategy  = fj::Strategy(cmdline.int_val("-strategy",
00250                                         cmdline.int_val("-clever", fj::Best)));
00251   int  repeat  = cmdline.int_val("-repeat",1);
00252   int  combine = cmdline.int_val("-combine",1);
00253   bool write   = cmdline.present("-write");
00254   bool unique_write = cmdline.present("-unique_write");
00255   bool hydjet  = cmdline.present("-hydjet");
00256   double ktR   = cmdline.double_val("-r",1.0);
00257   ktR   = cmdline.double_val("-R",ktR); // allow -r and -R
00258   double inclkt = cmdline.double_val("-incl",-1.0);
00259   int    excln  = cmdline.int_val   ("-excln",-1);
00260   double excld  = cmdline.double_val("-excld",-1.0);
00261   double excly  = cmdline.double_val("-excly",-1.0);
00262   ee_print = cmdline.present("-ee-print");
00263   bool   get_all_dij   = cmdline.present("-get-all-dij");
00264   bool   get_all_yij   = cmdline.present("-get-all-yij");
00265   double subdcut = cmdline.double_val("-subdcut",-1.0);
00266   double etamax = cmdline.double_val("-etamax",1.0e305);
00267   bool   show_constituents = cmdline.present("-const");
00268   bool   massless = cmdline.present("-massless");
00269   int    nev     = cmdline.int_val("-nev",1);
00270   bool   add_dense_coverage = cmdline.present("-dense");
00271   double ghost_maxrap = cmdline.value("-ghost-maxrap",5.0);
00272 
00273   bool show_cones = cmdline.present("-cones"); // only works for siscone
00274 
00275   // for cone algorithms
00276   // allow -f and -overlap
00277   double overlap_threshold = cmdline.double_val("-overlap",0.5);
00278   overlap_threshold = cmdline.double_val("-f",overlap_threshold); 
00279   double seed_threshold = cmdline.double_val("-seed",1.0);
00280 
00281   // for ee algorithms, allow to specify ycut
00282   double ycut = cmdline.double_val("-ycut",0.08);
00283 
00284   // for printing jets to a file for reading by root
00285   rootfile = cmdline.value<string>("-root","");
00286 
00287   // out default scheme is the E_scheme
00288   fj::RecombinationScheme scheme = fj::E_scheme;
00289 
00290   // The following option causes the Cambridge algo to be used.
00291   // Note that currently the only output that works sensibly here is
00292   // "-incl 0"
00293   fj::JetDefinition jet_def;
00294   if (cmdline.present("-cam") || cmdline.present("-CA")) {
00295     jet_def = fj::JetDefinition(fj::cambridge_algorithm, ktR, scheme, strategy);
00296   } else if (cmdline.present("-antikt")) {
00297     jet_def = fj::JetDefinition(fj::antikt_algorithm, ktR, scheme, strategy);
00298   } else if (cmdline.present("-genkt")) {
00299     double p = cmdline.value<double>("-genkt");
00300     jet_def = fj::JetDefinition(fj::genkt_algorithm, ktR, p, scheme, strategy);
00301   } else if (cmdline.present("-eekt")) {
00302     jet_def = fj::JetDefinition(fj::ee_kt_algorithm);
00303   } else if (cmdline.present("-eegenkt")) {
00304     double p = cmdline.value<double>("-eegenkt");
00305     jet_def = fj::JetDefinition(fj::ee_genkt_algorithm, ktR, p, scheme, strategy);
00306 
00307 // checking if one asks to run a plugin (don't delete this line)
00308   } else if (cmdline.present("-midpoint")) {
00309 #ifdef FASTJET_ENABLE_PLUGIN_CDFCONES
00310     typedef fj::CDFMidPointPlugin MPPlug; // for brevity
00311     double cone_area_fraction = 1.0;
00312     int    max_pair_size = 2;
00313     int    max_iterations = 100;
00314     MPPlug::SplitMergeScale sm_scale = MPPlug::SM_pt;
00315     if (cmdline.present("-sm-pttilde")) sm_scale = MPPlug::SM_pttilde;
00316     if (cmdline.present("-sm-pt")) sm_scale = MPPlug::SM_pt; // default
00317     if (cmdline.present("-sm-mt")) sm_scale = MPPlug::SM_mt;
00318     if (cmdline.present("-sm-Et")) sm_scale = MPPlug::SM_Et;
00319     jet_def = fj::JetDefinition( new fj::CDFMidPointPlugin (
00320                                       seed_threshold, ktR, 
00321                                       cone_area_fraction, max_pair_size,
00322                                       max_iterations, overlap_threshold,
00323                                       sm_scale));
00324 #else  // FASTJET_ENABLE_PLUGIN_CDFCONES
00325     is_unavailable("midpoint");
00326 #endif // FASTJET_ENABLE_PLUGIN_CDFCONES
00327   } else if (cmdline.present("-pxcone")) {
00328 #ifdef FASTJET_ENABLE_PLUGIN_PXCONE
00329     double min_jet_energy = 5.0;
00330     jet_def = fj::JetDefinition( new fj::PxConePlugin (
00331                                       ktR, min_jet_energy,
00332                                       overlap_threshold));
00333 #else  // FASTJET_ENABLE_PLUGIN_PXCONE
00334     is_unavailable("pxcone");
00335 #endif // FASTJET_ENABLE_PLUGIN_PXCONE
00336   } else if (cmdline.present("-jetclu")) {
00337 #ifdef FASTJET_ENABLE_PLUGIN_CDFCONES
00338     jet_def = fj::JetDefinition( new fj::CDFJetCluPlugin (
00339                                       ktR, overlap_threshold, seed_threshold));
00340 #else  // FASTJET_ENABLE_PLUGIN_CDFCONES
00341     is_unavailable("pxcone");
00342 #endif // FASTJET_ENABLE_PLUGIN_CDFCONES
00343   } else if (cmdline.present("-siscone") || cmdline.present("-sisconespheri")) {
00344 #ifdef FASTJET_ENABLE_PLUGIN_SISCONE
00345     typedef fj::SISConePlugin SISPlug; // for brevity
00346     int npass = cmdline.value("-npass",0);
00347     if (cmdline.present("-siscone")) {
00348       double sisptmin = cmdline.value("-sisptmin",0.0);
00349       SISPlug * plugin = new SISPlug (ktR, overlap_threshold,npass,sisptmin);
00350       if (cmdline.present("-sm-pt")) plugin->set_split_merge_scale(SISPlug::SM_pt);
00351       if (cmdline.present("-sm-mt")) plugin->set_split_merge_scale(SISPlug::SM_mt);
00352       if (cmdline.present("-sm-Et")) plugin->set_split_merge_scale(SISPlug::SM_Et);
00353       if (cmdline.present("-sm-pttilde")) plugin->set_split_merge_scale(SISPlug::SM_pttilde);
00354       // cause it to use the jet-definition's own recombiner
00355       plugin->set_use_jet_def_recombiner(true);
00356       jet_def = fj::JetDefinition(plugin);
00357     } else {
00358       double sisEmin = cmdline.value("-sisEmin",0.0);
00359       fj::SISConeSphericalPlugin * plugin = 
00360         new fj::SISConeSphericalPlugin(ktR, overlap_threshold,npass,sisEmin);
00361       if (cmdline.present("-ghost-sep")) {
00362         plugin->set_ghost_separation_scale(cmdline.value<double>("-ghost-sep"));
00363       }
00364       jet_def = fj::JetDefinition(plugin);
00365     }
00366 #else  // FASTJET_ENABLE_PLUGIN_SISCONE
00367     is_unavailable("siscone");
00368 #endif // FASTJET_ENABLE_PLUGIN_SISCONE
00369   } else if (cmdline.present("-d0runiicone")) {
00370 #ifdef FASTJET_ENABLE_PLUGIN_D0RUNIICONE
00371     double min_jet_Et = 6.0; // was 8 GeV in earlier work
00372     jet_def = fj::JetDefinition(new fj::D0RunIIConePlugin(ktR,min_jet_Et));
00373 #else  // FASTJET_ENABLE_PLUGIN_D0RUNIICONE
00374     is_unavailable("D0RunIICone");
00375 #endif // FASTJET_ENABLE_PLUGIN_D0RUNIICONE
00376   } else if (cmdline.present("-trackjet")) {
00377 #ifdef FASTJET_ENABLE_PLUGIN_TRACKJET
00378     jet_def = fj::JetDefinition(new fj::TrackJetPlugin(ktR));
00379 #else  // FASTJET_ENABLE_PLUGIN_TRACKJET
00380     is_unavailable("TrackJet");
00381 #endif // FASTJET_ENABLE_PLUGIN_TRACKJET
00382   } else if (cmdline.present("-atlascone")) {
00383 #ifdef FASTJET_ENABLE_PLUGIN_ATLASCONE
00384     jet_def = fj::JetDefinition(new fj::ATLASConePlugin(ktR));
00385 #else  // FASTJET_ENABLE_PLUGIN_ATLASCONE
00386     is_unavailable("ATLASCone");
00387 #endif // FASTJET_ENABLE_PLUGIN_ATLASCONE
00388   } else if (cmdline.present("-eecambridge")) {
00389 #ifdef FASTJET_ENABLE_PLUGIN_EECAMBRIDGE
00390     jet_def = fj::JetDefinition(new fj::EECambridgePlugin(ycut));
00391 #else  // FASTJET_ENABLE_PLUGIN_EECAMBRIDGE
00392     is_unavailable("EECambridge");
00393 #endif // FASTJET_ENABLE_PLUGIN_EECAMBRIDGE
00394   } else if (cmdline.present("-jade")) {
00395 #ifdef FASTJET_ENABLE_PLUGIN_JADE
00396     jet_def = fj::JetDefinition(new fj::JadePlugin());
00397 #else  // FASTJET_ENABLE_PLUGIN_JADE
00398     is_unavailable("Jade");
00399 #endif // FASTJET_ENABLE_PLUGIN_JADE
00400   } else if (cmdline.present("-cmsiterativecone")) {
00401 #ifdef FASTJET_ENABLE_PLUGIN_CMSITERATIVECONE
00402     jet_def = fj::JetDefinition(new fj::CMSIterativeConePlugin(ktR,seed_threshold));
00403 #else  // FASTJET_ENABLE_PLUGIN_CMSITERATIVECONE
00404     is_unavailable("CMSIterativeCone");
00405 #endif // FASTJET_ENABLE_PLUGIN_CMSITERATIVECONE
00406   } else if (cmdline.present("-d0runipre96cone")) {
00407 #ifdef FASTJET_ENABLE_PLUGIN_D0RUNICONE
00408     jet_def = fj::JetDefinition(new fj::D0RunIpre96ConePlugin(ktR, seed_threshold, overlap_threshold));
00409 #else  // FASTJET_ENABLE_PLUGIN_D0RUNICONE
00410     is_unavailable("D0RunICone");
00411 #endif // FASTJET_ENABLE_PLUGIN_D0RUNICONE
00412   } else if (cmdline.present("-d0runicone")) {
00413 #ifdef FASTJET_ENABLE_PLUGIN_D0RUNICONE
00414     jet_def = fj::JetDefinition(new fj::D0RunIConePlugin(ktR, seed_threshold, overlap_threshold));
00415 #else  // FASTJET_ENABLE_PLUGIN_D0RUNICONE
00416     is_unavailable("D0RunICone");
00417 #endif // FASTJET_ENABLE_PLUGIN_D0RUNICONE
00418 // end of checking if one asks to run a plugin (don't delete this line)
00419   } else {
00420     cmdline.present("-kt"); // kt is default, but allow user to specify it too [and ignore return value!]
00421     jet_def = fj::JetDefinition(fj::kt_algorithm, ktR, strategy);
00422   }
00423 
00424 
00425 
00426   if (!cmdline.all_options_used()) {cerr << 
00427       "Error: some options were not recognized"<<endl; 
00428     exit(-1);}
00429 
00430 
00431   for (int iev = 0; iev < nev; iev++) {
00432   vector<fj::PseudoJet> jets;
00433   string line;
00434   int  ndone = 0;
00435   while (getline(cin, line)) {
00436       //cout << line<<endl;
00437     istringstream linestream(line);
00438     if (line == "#END") {
00439       ndone += 1;
00440       if (ndone == combine) {break;}
00441     }
00442     if (line.substr(0,1) == "#") {continue;}
00443     valarray<double> fourvec(4);
00444     if (hydjet) {
00445       // special reading from hydjet.txt event record (though actually
00446       // this is supposed to be a standard pythia event record, so
00447       // being able to read from it is perhaps not so bad an idea...)
00448       int ii, istat,id,m1,m2,d1,d2;
00449       double mass;
00450       linestream >> ii>> istat >> id >> m1 >> m2 >> d1 >> d2
00451                  >> fourvec[0] >> fourvec[1] >> fourvec[2] >> mass;
00452       // current file contains mass of particle as 4th entry
00453       if (istat == 1) {
00454         fourvec[3] = sqrt(+pow2(fourvec[0])+pow2(fourvec[1])
00455                           +pow2(fourvec[2])+pow2(mass));
00456       }
00457     } else {
00458       if (massless) {
00459         linestream >> fourvec[0] >> fourvec[1] >> fourvec[2];
00460         fourvec[3] = sqrt(pow2(fourvec[0])+pow2(fourvec[1])+pow2(fourvec[2]));}
00461       else {
00462         linestream >> fourvec[0] >> fourvec[1] >> fourvec[2] >> fourvec[3];
00463       }
00464     }
00465     fj::PseudoJet psjet(fourvec);
00466     if (abs(psjet.rap() < etamax)) {jets.push_back(psjet);}
00467   }
00468 
00469   // add a fake underlying event which is very soft, uniformly distributed
00470   // in eta,phi so as to allow one to reconstruct the area that is associated
00471   // with each jet.
00472   if (add_dense_coverage) {
00473     fj::GhostedAreaSpec ghosted_area_spec(ghost_maxrap);
00474     //fj::GhostedAreaSpec ghosted_area_spec(-2.0,4.0); // asymmetric range
00475     // for plots, reduce the scatter default of 1, to avoid "holes"
00476     // in the subsequent calorimeter view
00477     ghosted_area_spec.set_grid_scatter(0.5); 
00478     ghosted_area_spec.add_ghosts(jets);
00479     //----- old code ------------------
00480     // srand(2);
00481     // int nphi = 60;
00482     // int neta = 100;
00483     // double kt = 1e-1;
00484     // for (int iphi = 0; iphi<nphi; iphi++) {
00485     //   for (int ieta = -neta; ieta<neta+1; ieta++) {
00486     //  double phi = (iphi+0.5) * (fj::twopi/nphi) + rand()*0.001/RAND_MAX;
00487     //  double eta = ieta * (10.0/neta)  + rand()*0.001/RAND_MAX;
00488     //  kt = 1e-20*(1+rand()*0.1/RAND_MAX);
00489     //  double pminus = kt*exp(-eta);
00490     //  double pplus  = kt*exp(+eta);
00491     //  double px = kt*sin(phi);
00492     //  double py = kt*cos(phi);
00493     //  //cout << kt<<" "<<eta<<" "<<phi<<"\n";
00494     //  fj::PseudoJet mom(px,py,0.5*(pplus-pminus),0.5*(pplus+pminus));
00495     //  jets.push_back(mom);
00496     //   }
00497     // }
00498   }
00499   
00500   for (int irepeat = 0; irepeat < repeat ; irepeat++) {
00501     int nparticles = jets.size();
00502     try {
00503     fj::ClusterSequence clust_seq(jets,jet_def,write);
00504     if (irepeat != 0) {continue;}
00505     cout << "iev "<<iev<< ": number of particles = "<< nparticles << endl;
00506     cout << "strategy used =  "<< clust_seq.strategy_string()<< endl;
00507     cout << "Algorithm: " << jet_def.description() << " (" << fj::fastjet_version_string() << ")" << endl;
00508 
00509     // now provide some nice output...
00510     if (inclkt >= 0.0) {
00511       vector<fj::PseudoJet> jets = sorted_by_pt(clust_seq.inclusive_jets(inclkt));
00512       print_jets(jets, show_constituents);
00513 
00514     }
00515 
00516     if (excln > 0) {
00517       cout << "Printing "<<excln<<" exclusive jets\n";
00518       print_jets(clust_seq.exclusive_jets(excln), show_constituents);
00519     }
00520 
00521     if (excld > 0.0) {
00522       cout << "Printing exclusive jets for d = "<<excld<<"\n";
00523       print_jets(clust_seq.exclusive_jets(excld), show_constituents);
00524     }
00525 
00526     if (excly > 0.0) {
00527       cout << "Printing exclusive jets for ycut = "<<excly<<"\n";
00528       print_jets(clust_seq.exclusive_jets_ycut(excly), show_constituents);
00529     }
00530 
00531     if (get_all_dij) {
00532       for (int i = nparticles-1; i >= 0; i--) {
00533         printf("d for n = %4d -> %4d is %14.5e\n", i+1, i, clust_seq.exclusive_dmerge(i));
00534       }
00535     }
00536     if (get_all_yij) {
00537       for (int i = nparticles-1; i >= 0; i--) {
00538         printf("y for n = %4d -> %4d is %14.5e\n", i+1, i, clust_seq.exclusive_ymerge(i));
00539       }
00540     }
00541 
00542     // have the option of printing out the subjets (at scale dcut) of
00543     // each inclusive jet
00544     if (subdcut >= 0.0) {
00545       print_jets_and_sub(clust_seq.inclusive_jets(), subdcut);
00546     }
00547     
00548     // useful for testing that recombination sequences are unique
00549     if (unique_write) {
00550       vector<int> unique_history = clust_seq.unique_history_order();
00551       // construct the inverse of the above mapping
00552       vector<int> inv_unique_history(clust_seq.history().size());
00553       for (unsigned int i = 0; i < unique_history.size(); i++) {
00554         inv_unique_history[unique_history[i]] = i;}
00555 
00556       for (unsigned int i = 0; i < unique_history.size(); i++) {
00557         fj::ClusterSequence::history_element el = 
00558           clust_seq.history()[unique_history[i]];
00559         int uhp1 = el.parent1>=0 ? inv_unique_history[el.parent1] : el.parent1;
00560         int uhp2 = el.parent2>=0 ? inv_unique_history[el.parent2] : el.parent2;
00561         printf("%7d u %15.8e %7d u %7d u\n",i,el.dij,uhp1, uhp2);
00562       }
00563     }
00564 
00565 
00566 #ifdef FASTJET_ENABLE_PLUGIN_SISCONE
00567     // provide some complementary information for SISCone 
00568     if (show_cones) {
00569       const fj::SISConeExtras * extras = 
00570         dynamic_cast<const fj::SISConeExtras *>(clust_seq.extras());
00571       cout << "most ambiguous split (difference in squared dist) = "
00572            << extras->most_ambiguous_split() << endl;
00573       vector<fastjet::PseudoJet> stable_cones(extras->stable_cones()); 
00574       stable_cones = sorted_by_rapidity(stable_cones);
00575       for (unsigned int i = 0; i < stable_cones.size(); i++) {
00576       //if (stable_cones[i].phi() < 5.0 && stable_cones[i].phi() > 4.0) {
00577         printf("%5u %15.8f %15.8f %15.8f\n",
00578                i,stable_cones[i].rap(),stable_cones[i].phi(),
00579                stable_cones[i].perp() );
00580       //}
00581       }
00582       
00583       // also show passes for jets
00584       vector<fj::PseudoJet> sisjets = clust_seq.inclusive_jets();
00585       printf("\n%15s %15s %15s %12s %8s %8s\n","rap","phi","pt","user-index","pass","nconst");
00586       for (unsigned i = 0; i < sisjets.size(); i++) {
00587         printf("%15.8f %15.8f %15.8f %12d %8d %8d\n",
00588                sisjets[i].rap(), sisjets[i].phi(), sisjets[i].perp(), 
00589                sisjets[i].user_index(), extras->pass(sisjets[i]),
00590                clust_seq.constituents(sisjets[i]).size()
00591                );
00592         
00593       }
00594     }
00595 #endif // FASTJET_ENABLE_PLUGIN_SISCONE
00596   } // try
00597   catch (fastjet::Error fjerr) {
00598     cout << "Caught fastjet error, exiting gracefully" << endl;
00599     exit(0);
00600   }
00601 
00602   } // irepeat
00603   } // iev
00604 
00605   // if we've instantiated a plugin, delete it
00606   if (jet_def.strategy()==fj::plugin_strategy){
00607     delete jet_def.plugin();
00608   }
00609 }
00610 
00611 
00612 
00613 
00614 //------ HELPER ROUTINES -----------------------------------------------
00615 /// print a single jet
00616 void print_jet (const fj::PseudoJet & jet) {
00617   int n_constituents = jet.constituents().size();
00618   printf("%15.8f %15.8f %15.8f %8u\n",
00619          jet.rap(), jet.phi(), jet.perp(), n_constituents);
00620 }
00621 
00622 
00623 //----------------------------------------------------------------------
00624 void print_jets(const vector<fj::PseudoJet> & jets_in, bool show_constituents) {
00625   vector<fj::PseudoJet> jets;
00626   if (ee_print) {
00627     jets = sorted_by_E(jets_in);
00628     for (size_t j = 0; j < jets.size(); j++) {
00629       printf("%5u %15.8f %15.8f %15.8f %15.8f\n",
00630              j,jets[j].px(),jets[j].py(),jets[j].pz(),jets[j].E());
00631       if (show_constituents) {
00632         vector<fj::PseudoJet> const_jets = jets[j].constituents();
00633         for (size_t k = 0; k < const_jets.size(); k++) {
00634           printf("        jet%03u %15.8f %15.8f %15.8f %15.8f\n",j,const_jets[k].px(),
00635                  const_jets[k].py(),const_jets[k].pz(),const_jets[k].E());
00636         }
00637         cout << "\n\n";
00638     }
00639 
00640     }
00641   } else {
00642     jets = sorted_by_pt(jets_in);
00643     for (size_t j = 0; j < jets.size(); j++) {
00644       printf("%5u %15.8f %15.8f %15.8f\n",
00645              j,jets[j].rap(),jets[j].phi(),jets[j].perp());
00646 
00647       if (show_constituents) {
00648         vector<fj::PseudoJet> const_jets = jets[j].constituents();
00649         for (size_t k = 0; k < const_jets.size(); k++) {
00650           printf("        jet%03u %15.8f %15.8f %15.8f %5d\n",j,const_jets[k].rap(),
00651                  const_jets[k].phi(),sqrt(const_jets[k].kt2()), const_jets[k].cluster_hist_index());
00652         }
00653         cout << "\n\n";
00654       }
00655     }
00656   }
00657 
00658   if (rootfile != "") {
00659     ofstream ostr(rootfile.c_str());
00660     ostr << "# " << cmdline_p->command_line() << endl;
00661     ostr << "# output for root" << endl;
00662     assert(jets.size() > 0);
00663     jets[0].validated_cs()->print_jets_for_root(jets,ostr);
00664   }
00665 
00666 }
00667 
00668 
00669 //----- SUBJETS --------------------------------------------------------
00670 /// a function that pretty prints a list of jets and the subjets for each
00671 /// one
00672 void print_jets_and_sub (const vector<fj::PseudoJet> & jets, double dcut) {
00673 
00674   // sort jets into increasing pt
00675   vector<fj::PseudoJet> sorted_jets = sorted_by_pt(jets);  
00676 
00677   // label the columns
00678   printf("Printing jets and their subjets with subdcut = %10.5f\n",dcut);
00679   printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity", 
00680          "phi", "pt", "n constituents");
00681 
00682   // have various kinds of subjet finding, to test consistency among them
00683   enum SubType {internal, newclust_dcut, newclust_R};
00684   SubType subtype = internal;
00685   //SubType subtype = newclust_dcut;
00686   //SubType subtype = newclust_R;
00687 
00688   // print out the details for each jet
00689   //for (unsigned int i = 0; i < sorted_jets.size(); i++) {
00690   for (vector<fj::PseudoJet>::const_iterator jet = sorted_jets.begin(); 
00691        jet != sorted_jets.end(); jet++) {
00692     const fj::JetDefinition & jet_def = jet->validated_cs()->jet_def();
00693 
00694     // if jet pt^2 < dcut with kt alg, then some methods of
00695     // getting subjets will return nothing -- so skip the jet
00696     if (jet_def.jet_algorithm() == fj::kt_algorithm 
00697         && jet->perp2() < dcut) continue;
00698 
00699 
00700     printf("%5u       ",jet - sorted_jets.begin());
00701     print_jet(*jet);
00702     vector<fj::PseudoJet> subjets;
00703     fj::ClusterSequence * cspoint;
00704     if (subtype == internal) {
00705       cspoint = 0;
00706       subjets = jet->exclusive_subjets(dcut);
00707       double ddnp1 = jet->exclusive_subdmerge_max(subjets.size());
00708       double ddn   = jet->exclusive_subdmerge_max(subjets.size()-1);
00709       cout << "     for " << ddnp1 << " < d < " << ddn << " one has " << endl;
00710     } else if (subtype == newclust_dcut) {
00711       cspoint = new fj::ClusterSequence(jet->constituents(), jet_def);
00712       subjets = cspoint->exclusive_jets(dcut);
00713     } else if (subtype == newclust_R) {
00714       assert(jet_def.jet_algorithm() == fj::cambridge_algorithm);
00715       fj::JetDefinition subjd(jet_def.jet_algorithm(), jet_def.R()*sqrt(dcut));
00716       cspoint = new fj::ClusterSequence(jet->constituents(), subjd);
00717       subjets = cspoint->inclusive_jets();
00718     } else {
00719       cerr << "unrecognized subtype for subjet finding" << endl;
00720       exit(-1);
00721     }
00722 
00723     subjets = sorted_by_pt(subjets);
00724     for (unsigned int j = 0; j < subjets.size(); j++) {
00725       printf("    -sub-%02u ",j);
00726       print_jet(subjets[j]);
00727     }
00728 
00729     if (cspoint != 0) delete cspoint;
00730 
00731     //fj::ClusterSequence subseq(clust_seq.constituents(sorted_jets[i]),
00732     //                          fj::JetDefinition(fj::cambridge_algorithm, 0.4));
00733     //vector<fj::PseudoJet> subjets = sorted_by_pt(subseq.inclusive_jets());
00734     //for (unsigned int j = 0; j < subjets.size(); j++) {
00735     //  printf("    -sub-%02u ",j);
00736     //  print_jet(subseq, subjets[j]);
00737     //}
00738   }
00739 
00740 }
00741