FastJet 3.0.6
CDFJetCluPlugin.cc
00001 //STARTHEADER
00002 // $Id: CDFJetCluPlugin.cc 2777 2011-11-25 15:01:56Z 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 #include "fastjet/CDFJetCluPlugin.hh"
00030 #include "fastjet/ClusterSequence.hh"
00031 #include <sstream>
00032 #include <cassert>
00033 
00034 // CDF stuff
00035 #include "JetCluAlgorithm.hh"
00036 #include "PhysicsTower.hh"
00037 #include "Cluster.hh"
00038 
00039 FASTJET_BEGIN_NAMESPACE      // defined in fastjet/internal/base.hh
00040 
00041 using namespace std;
00042 using namespace cdf;
00043 
00044 bool CDFJetCluPlugin::_first_time = true;
00045 
00046 string CDFJetCluPlugin::description () const {
00047   ostringstream desc;
00048   
00049   desc << "CDF JetClu jet algorithm with " 
00050        << "seed_threshold = "     << seed_threshold    () << ", "
00051        << "cone_radius = "        << cone_radius       () << ", "
00052        << "adjacency_cut = "      << adjacency_cut     () << ", " 
00053        << "max_iterations = "     << max_iterations    () << ", "
00054        << "iratch = "             << iratch            () << ", "
00055        << "overlap_threshold = "  << overlap_threshold () ;
00056 
00057   return desc.str();
00058 }
00059 
00060 
00061 void CDFJetCluPlugin::run_clustering(ClusterSequence & clust_seq) const {
00062   // print a banner if we run this for the first time
00063   _print_banner(clust_seq.fastjet_banner_stream());
00064  
00065   // create the physics towers needed by the CDF code
00066   vector<PhysicsTower> towers;
00067   towers.reserve(clust_seq.jets().size());
00068 
00069   // create a map to identify jets (actually just the input particles)...
00070   //map<double,int> jetmap;
00071 
00072   for (unsigned i = 0; i < clust_seq.jets().size(); i++) {
00073     PseudoJet particle(clust_seq.jets()[i]);
00074     //_insert_unique(particle, jetmap);
00075     LorentzVector fourvect(particle.px(), particle.py(),
00076                            particle.pz(), particle.E());
00077     PhysicsTower tower(fourvect);
00078     // add tracking information for later
00079     tower.fjindex = i;
00080     towers.push_back(tower);
00081   }
00082 
00083   // prepare the CDF algorithm
00084   JetCluAlgorithm j(seed_threshold(), cone_radius(), adjacency_cut(),
00085                     max_iterations(), iratch(), overlap_threshold());
00086     
00087   // run the CDF algorithm
00088   std::vector<Cluster> jets;
00089   j.run(towers,jets);
00090 
00091 
00092   // now transfer the jets back into our own structure -- we will
00093   // mimic the cone code with a sequential recombination sequence in
00094   // which the jets are built up by adding one particle at a time
00095 
00096   // NB: with g++-4.0, the reverse iterator code gave problems, so switch
00097   //     to indices instead
00098   //for(vector<Cluster>::const_reverse_iterator jetIter = jets.rbegin(); 
00099   //                                    jetIter != jets.rend(); jetIter++) {
00100   //  const vector<PhysicsTower> & tower_list = jetIter->towerList;
00101   //  int jet_k = jetmap[tower_list[0].fourVector.E];
00102   //
00103   //  int ntow = int(jetIter->towerList.size());
00104 
00105   for(int iCDFjets = jets.size()-1; iCDFjets >= 0; iCDFjets--) {
00106 
00107     const vector<PhysicsTower> & tower_list = jets[iCDFjets].towerList;
00108     int ntow = int(tower_list.size());
00109     
00110     // 2008-09-04: sort the towerList (according to fjindex) so as
00111     //             to have a consistent order for particles in jet
00112     //             (necessary because addition of ultra-soft particles
00113     //             sometimes often modifies the order, while maintaining
00114     //             the same overall set)
00115     vector<int>    jc_indices(ntow);
00116     vector<double> fj_indices(ntow); // use double: benefit from existing routine
00117     for (int itow = 0; itow < ntow; itow++) {
00118       jc_indices[itow] = itow;
00119       fj_indices[itow] = tower_list[itow].fjindex;
00120     }
00121     sort_indices(jc_indices, fj_indices);
00122 
00123     int jet_k = tower_list[jc_indices[0]].fjindex;
00124   
00125     for (int itow = 1; itow < ntow; itow++) {
00126       if (tower_list[jc_indices[itow]].Et() > 1e-50) {
00127       }
00128       int jet_i = jet_k;
00129       // retrieve our index for the jet
00130       int jet_j;
00131       jet_j = tower_list[jc_indices[itow]].fjindex;
00132 
00133       // safety check
00134       assert (jet_j >= 0 && jet_j < int(towers.size()));
00135 
00136       // do a fake recombination step with dij=0
00137       double dij = 0.0;
00138 
00139       // JetClu does E-scheme recombination so we can stick with the
00140       // simple option
00141       clust_seq.plugin_record_ij_recombination(jet_i, jet_j, dij, jet_k);
00142 
00143     }
00144   
00145     // NB: put a sensible looking d_iB just to be nice...
00146     double d_iB = clust_seq.jets()[jet_k].perp2();
00147     clust_seq.plugin_record_iB_recombination(jet_k, d_iB);
00148   }
00149 
00150 
00151   // following code is for testing only
00152   //cout << endl;
00153   //for(vector<Cluster>::const_iterator jetIter = jets.begin(); 
00154   //                                    jetIter != jets.end(); jetIter++) {
00155   //  cout << jetIter->fourVector.pt() << " " << jetIter->fourVector.y() << endl;
00156   //}
00157   //cout << "-----------------------------------------------------\n";
00158   //vector<PseudoJet> ourjets(clust_seq.inclusive_jets());
00159   //for (vector<PseudoJet>::const_iterator ourjet = ourjets.begin();
00160   //     ourjet != ourjets.end(); ourjet++) {
00161   //  cout << ourjet->perp() << " " << ourjet->rap() << endl;
00162   //}
00163   //cout << endl;
00164 }
00165 
00166 
00167 // print a banner for reference to the 3rd-party code
00168 void CDFJetCluPlugin::_print_banner(ostream *ostr) const{
00169   if (! _first_time) return;
00170   _first_time=false;
00171 
00172   // make sure the user has not set the banner stream to NULL
00173   if (!ostr) return;  
00174 
00175   (*ostr) << "#-------------------------------------------------------------------------" << endl;
00176   (*ostr) << "# You are running the CDF JetClu plugin for FastJet                       " << endl;
00177   (*ostr) << "# This is based on an implementation provided by Joey Huston.             " << endl;
00178   (*ostr) << "# If you use this plugin, please cite                                     " << endl;
00179   (*ostr) << "#   F. Abe et al. [CDF Collaboration], Phys. Rev. D 45 (1992) 1448.       " << endl;
00180   (*ostr) << "# in addition to the usual FastJet reference.                             " << endl;
00181   (*ostr) << "#-------------------------------------------------------------------------" << endl;
00182 
00183   // make sure we really have the output done.
00184   ostr->flush();
00185 }
00186 
00187 FASTJET_END_NAMESPACE      // defined in fastjet/internal/base.hh
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends