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