FastJet 3.0.0
ATLASConePlugin.cc
00001 //STARTHEADER
00002 // $Id: ATLASConePlugin.cc 2577 2011-09-13 15:11:38Z salam $
00003 //
00004 // Copyright (c) 2007-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 // fastjet stuff
00030 #include "fastjet/ClusterSequence.hh"
00031 #include "fastjet/ATLASConePlugin.hh"
00032 
00033 // SpartyJet stuff
00034 #include "CommonUtils.hh"
00035 #include "JetConeFinderTool.hh"
00036 #include "JetSplitMergeTool.hh"
00037 
00038 // other stuff
00039 #include <vector>
00040 #include <sstream>
00041 
00042 FASTJET_BEGIN_NAMESPACE      // defined in fastjet/internal/base.hh
00043 
00044 using namespace std;
00045 
00046 string ATLASConePlugin::description () const {
00047   ostringstream desc;
00048   desc << "ATLASCone plugin with R = "<< _radius 
00049        << ", seed threshold = " << _seedPt
00050        << ", overlap threshold f = " << _f;
00051   return desc.str();
00052 }
00053 
00054 void ATLASConePlugin::run_clustering(ClusterSequence & clust_seq) const {
00055 
00056   // transfer the list of PseudoJet into a atlas::Jet::jet_list_t
00057   //  jet_list_t is a vector<Jet*>
00058   // We set the index of the 4-vect to trace the constituents at the end
00059   //------------------------------------------------------------------
00060   // cout << "ATLASConePlugin: transferring vectors from ClusterSequence" << endl;
00061   atlas::JetConeFinderTool::jetcollection_t jets_ptr;
00062   vector<atlas::Jet*> particles_ptr;
00063 
00064   for (unsigned int i=0 ; i<clust_seq.jets().size() ; i++) {
00065     const PseudoJet & mom = clust_seq.jets()[i];
00066     
00067     // first create the particle
00068     atlas::Jet *particle = new atlas::Jet(mom.px(), mom.py(), mom.pz(), mom.E(), i);
00069     particles_ptr.push_back(particle);
00070 
00071     // then add it to the list of particles we'll use for teh clustering
00072     atlas::Jet *jet = new atlas::Jet;
00073     jet->set_index(particle->index());
00074     jet->addConstituent(particle);
00075 
00076     // and finally add that one to the list of jets
00077     jets_ptr.push_back(jet);
00078   }
00079   // cout << "ATLASCone: " << jets_ptr.size() << " particles to cluster" << endl;
00080 
00081   // search the stable cones
00082   //------------------------------------------------------------------
00083   // cout << "ATLASConePlugin: searching for stable cones" << endl;
00084   atlas::JetConeFinderTool stable_cone_finder;
00085 
00086   // set the parameters
00087   stable_cone_finder.m_coneR  = _radius;
00088   stable_cone_finder.m_seedPt = _seedPt;
00089 
00090   // really do the search.
00091   // Note that the list of protocones is returned 
00092   // through the argument
00093   stable_cone_finder.execute(jets_ptr);
00094   // cout << "ATLASCone: " << jets_ptr.size() << " stable cones found" << endl;
00095 
00096   // perform the split-merge
00097   //------------------------------------------------------------------
00098   // cout << "ATLASConePlugin: running the split-merge" << endl;
00099   atlas::JetSplitMergeTool split_merge;
00100   
00101   // set the parameters
00102   split_merge.m_f = _f;
00103 
00104   // do the work
00105   // again, the list of jets is returned through the argument
00106   split_merge.execute(&jets_ptr);
00107   // cout << "ATLASCone: " << jets_ptr.size() << " jets after split--merge" << endl;
00108 
00109   // build the FastJet jets (a la SISConePlugin)
00110   //------------------------------------------------------------------
00111   // cout << "ATLASConePlugin: backporting jets to the ClusterSequence" << endl;
00112   for (atlas::Jet::jet_list_t::iterator jet_it = jets_ptr.begin() ;
00113          jet_it != jets_ptr.end(); jet_it++){
00114     // iterate over the constituents, starting from the first one
00115     // that we just take as a reference
00116     atlas::Jet::constit_vect_t::iterator constit_it = (*jet_it)->firstConstituent();
00117     // cout << " atlas: jet has " << (*jet_it)->getConstituentNum() << " constituents" << endl;
00118     int jet_k = (*constit_it)->index();
00119     constit_it++;
00120     
00121     // loop over the remaining particles
00122     while (constit_it != (*jet_it)->lastConstituent()){
00123       // take the last result of the merge
00124       int jet_i = jet_k;
00125       // and the next element of the jet
00126       int jet_j = (*constit_it)->index();
00127       // and merge them (with a fake dij)
00128       double dij = 0.0;
00129 
00130       // create the new jet by hand so that we can adjust its user index
00131       // Note again the use of the E-scheme recombination here!
00132       PseudoJet newjet = clust_seq.jets()[jet_i] + clust_seq.jets()[jet_j];
00133       clust_seq.plugin_record_ij_recombination(jet_i, jet_j, dij, newjet, jet_k);
00134 
00135       // jump to the next constituent
00136       constit_it++;
00137     }
00138 
00139     // we have merged all the jet's particles into a single object, so now
00140     // "declare" it to be a beam (inclusive) jet.
00141     // [NB: put a sensible looking d_iB just to be nice...]
00142     double d_iB = clust_seq.jets()[jet_k].perp2();
00143     clust_seq.plugin_record_iB_recombination(jet_k, d_iB);
00144   }
00145     
00146   // cout << "ATLASConePlugin: Bye" << endl;
00147   clear_list(particles_ptr);
00148 }
00149 
00150 FASTJET_END_NAMESPACE      // defined in fastjet/internal/base.hh
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends