fastjet 2.4.5
TrackJetPlugin.cc
Go to the documentation of this file.
00001 //STARTHEADER
00002 // $Id: TrackJetPlugin.cc 1558 2009-05-10 14:19:05Z salam $
00003 //
00004 // Copyright (c) 2007-2008, 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 // fastjet stuff
00032 #include "fastjet/ClusterSequence.hh"
00033 #include "fastjet/TrackJetPlugin.hh"
00034 
00035 // other stuff
00036 #include <list>
00037 #include <memory>
00038 #include <cmath>
00039 #include <vector>
00040 #include <sstream>
00041 
00042 FASTJET_BEGIN_NAMESPACE      // defined in fastjet/internal/base.hh
00043 
00044 using namespace std;
00045 
00046 class TrackJetParticlePtr{
00047 public:
00048   TrackJetParticlePtr(int i_index, double i_perp2)
00049     :  index(i_index), perp2(i_perp2){}
00050 
00051   int index;
00052   double perp2;
00053 
00054   bool operator <(const TrackJetParticlePtr &other) const { 
00055     return perp2>other.perp2;
00056   }
00057 };
00058 
00059 string TrackJetPlugin::description () const {
00060   ostringstream desc;
00061   desc << "TrackJet algorithm with R = " << R();
00062   return desc.str();
00063 }
00064 
00065 void TrackJetPlugin::run_clustering(ClusterSequence & clust_seq) const {
00066   // we first need to sort the particles in pt
00067   vector<TrackJetParticlePtr> particle_list;
00068 
00069   const vector<PseudoJet> & jets = clust_seq.jets();  
00070   int index=0;
00071   for (vector<PseudoJet>::const_iterator mom_it = jets.begin(); mom_it != jets.end(); mom_it++){
00072     particle_list.push_back(TrackJetParticlePtr(index, mom_it->perp2()));
00073     index++;
00074   }
00075 
00076   // sort the particles into decreasing pt
00077   sort(particle_list.begin(), particle_list.end());
00078 
00079 
00080   // if we're using a recombination scheme different from the E scheme,
00081   // we first need to update the particles' energy so that they
00082   // are massless (and rapidity = pseudorapidity)
00083   vector<PseudoJet> tuned_particles = clust_seq.jets();
00084   vector<PseudoJet> tuned_tracks = clust_seq.jets();
00085   for (vector<PseudoJet>::iterator pit = tuned_particles.begin();
00086        pit != tuned_particles.end(); pit++)
00087     _jet_recombiner.preprocess(*pit);
00088   for (vector<PseudoJet>::iterator pit = tuned_tracks.begin();
00089        pit != tuned_tracks.end(); pit++)
00090     _track_recombiner.preprocess(*pit);
00091 
00092 
00093   // we'll just need the particle indices for what follows
00094   list<int> sorted_pt_index;
00095   for (vector<TrackJetParticlePtr>::iterator mom_it = particle_list.begin();
00096        mom_it != particle_list.end(); mom_it++)
00097     sorted_pt_index.push_back(mom_it->index);
00098   
00099   // now start building the jets
00100   while (sorted_pt_index.size()){
00101     // note that here 'track' refers to the direction we're using to test if a particle belongs to the jet
00102     // 'jet' refers to the momentum of the jet
00103     // the difference between the two is in the recombination scheme used to compute the sum of 4-vectors
00104     int current_jet_index = sorted_pt_index.front();
00105     PseudoJet current_jet   = tuned_particles[current_jet_index];
00106     PseudoJet current_track = tuned_tracks[current_jet_index];
00107 
00108     // remove the first particle from the available ones    
00109     list<int>::iterator index_it = sorted_pt_index.begin();
00110     sorted_pt_index.erase(index_it);
00111 
00112     // start browsing the remaining ones
00113     index_it = sorted_pt_index.begin();
00114     while (index_it != sorted_pt_index.end()){
00115       const PseudoJet & current_particle = tuned_particles[*index_it];
00116       const PseudoJet & current_particle_track = tuned_tracks[*index_it];
00117 
00118       // check if the particle is within a distance R of the jet
00119       double distance2 = current_track.plain_distance(current_particle_track);
00120       if (distance2 <= _radius2){
00121         // add the particle to the jet
00122         PseudoJet new_track;
00123         PseudoJet new_jet;
00124         _jet_recombiner.recombine(current_jet, current_particle, new_jet);
00125         _track_recombiner.recombine(current_track, current_particle_track, new_track);
00126 
00127         int new_jet_index;
00128         clust_seq.plugin_record_ij_recombination(current_jet_index, *index_it, distance2, new_jet, new_jet_index);
00129 
00130         current_jet = new_jet;
00131         current_track = new_track;
00132         current_jet_index = new_jet_index;
00133 
00134         // particle has been clustered so remove it from the list
00135         sorted_pt_index.erase(index_it);
00136 
00137         // and don't forget to start again from the beginning
00138         //  as the jet axis may have changed
00139         index_it = sorted_pt_index.begin();
00140       } else {
00141         index_it++;
00142       }
00143     }
00144 
00145     // now we have a final jet, so cluster it with the beam
00146     clust_seq.plugin_record_iB_recombination(current_jet_index, _radius2);
00147   }
00148     
00149 }
00150 
00151 FASTJET_END_NAMESPACE      // defined in fastjet/internal/base.hh
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines