FastJet 3.0beta1
GridJetPlugin.cc
00001 //STARTHEADER
00002 // $Id: GridJetPlugin.cc 2268 2011-06-20 15:12:26Z salam $
00003 //
00004 // Copyright (c) 2011, 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/GridJetPlugin.hh"
00034 
00035 // other stuff
00036 #include <vector>
00037 #include <sstream>
00038 
00039 FASTJET_BEGIN_NAMESPACE      // defined in fastjet/internal/base.hh
00040 
00041 using namespace std;
00042 
00043 //----------------------------------------------------------------------
00044 GridJetPlugin::GridJetPlugin (double ymax,
00045                               double requested_grid_spacing, 
00046                               const JetDefinition & post_jet_def) :
00047   _ymin(-ymax), _ymax(ymax), 
00048   _requested_grid_spacing(requested_grid_spacing) ,
00049   _post_jet_def(post_jet_def)
00050 {
00051   setup_grid();
00052 }
00053 
00054 void GridJetPlugin::setup_grid() {
00055   // since we've exchanged the arguments of the constructor,
00056   // there's a danger of calls with exchanged ymax,spacing arguments -- 
00057   // the following check should catch most such situations.
00058   assert(_ymax>0 && _ymax - _ymin >= _requested_grid_spacing);
00059 
00060   double ny_double = (_ymax-_ymin) / _requested_grid_spacing;
00061   _ny = int(ny_double+0.49999);
00062   _dy = (_ymax-_ymin) / _ny;
00063   
00064   _nphi = int (twopi / _requested_grid_spacing + 0.5);
00065   _dphi = twopi / _nphi;
00066 
00067   // some sanity checking (could throw a fastjet::Error)
00068   assert(_ny >= 1 && _nphi >= 1);
00069 
00070   _ntotal = _nphi * _ny;
00071 }
00072 
00073 
00074 //----------------------------------------------------------------------
00075 string GridJetPlugin::description () const {
00076   ostringstream desc;
00077   desc << "GridJetPlugin plugin with ymax = " << _ymax << ", dy = " << _dy << ", dphi = " << _dphi << " (requested grid spacing was " << _requested_grid_spacing << ")";
00078   if (_post_jet_def.jet_algorithm() != undefined_jet_algorithm) {
00079     desc << ", followed by " << _post_jet_def.description();
00080   }
00081   return desc.str();
00082 }
00083 
00084 
00085 //----------------------------------------------------------------------
00086 double GridJetPlugin::R() const {return sqrt(_dy*_dphi/pi);}
00087 
00088 
00089 //----------------------------------------------------------------------
00090 int GridJetPlugin::igrid(const PseudoJet & p) const {
00091   // directly taking int does not work for values between -1 and 0
00092   // so use floor instead
00093   // double iy_double = (p.rap() - _ymin) / _dy;
00094   // if (iy_double < 0.0) return -1;
00095   // int iy = int(iy_double);
00096   // if (iy >= _ny) return -1;
00097 
00098   // writing it as below gives a huge speed gain (factor two!). Even
00099   // though answers are identical and the routine here is not the
00100   // speed-critical step. It's not at all clear why.
00101   int iy = int(floor( (p.rap() - _ymin) / _dy ));
00102   if (iy < 0 || iy >= _ny) return -1;
00103 
00104   int iphi = int( p.phi()/_dphi );
00105   assert(iphi >= 0 && iphi <= _nphi);
00106   if (iphi == _nphi) iphi = 0; // just in case of rounding errors
00107 
00108   int igrid_res = iy*_nphi + iphi;
00109   assert (igrid_res >= 0 && igrid_res < _ny*_nphi);
00110   return igrid_res;
00111 }
00112 
00113 
00114 //----------------------------------------------------------------------
00115 void GridJetPlugin::run_clustering(ClusterSequence & cs) const {
00116   
00117   // we will create a grid; 
00118   //  * -1 will indicate there is no jet here currently
00119   //  * a number >= 0 will mean that particle indicated by the index
00120   //    is currently the jet on the grid
00121   vector<int> grid(_ntotal, -1);
00122  
00123   int nparticles = cs.jets().size();
00124   double dij_or_diB = 1.0;
00125 
00126   int ngrid_active = 0;
00127 
00128   // combine particles with whatever is in the grid
00129   for (int i = 0; i < nparticles; i++) {
00130     int igrd = igrid(cs.jets()[i]);
00131     //cout << i << " " << cs.jets()[i].rap() << " " << cs.jets()[i].phi() 
00132     //   << " " << igrd << " " << grid.size() << " " << _ntotal << endl;
00133     if (igrd < 0) continue;
00134     assert(igrd <= _ntotal);
00135     if (grid[igrd] == -1) {
00136       grid[igrd] = i; // jet index of initial particle i is i
00137       ngrid_active++;
00138     } else {
00139       int k;
00140       cs.plugin_record_ij_recombination(grid[igrd], i, dij_or_diB, k);
00141       grid[igrd] = k; // grid takes jet index of new particle
00142       //cout << "  res: " << cs.jets()[k].rap() << " " << cs.jets()[k].phi() << endl;
00143     }
00144   }
00145 
00146   if (_post_jet_def.jet_algorithm() == undefined_jet_algorithm) {
00147     // make the final jets via iB recombinations
00148     for (unsigned igrd = 0; igrd < grid.size(); igrd++) {
00149       if (grid[igrd] != -1) cs.plugin_record_iB_recombination(grid[igrd], 
00150                                                               dij_or_diB);
00151     }
00152   } else {
00153     // otherwise post-cluster the grid elements with a normal jet algorithm
00154     vector<PseudoJet> inputs;
00155     vector<int>       cs_indices;
00156     inputs.reserve(ngrid_active);
00157     cs_indices.reserve(2*ngrid_active);
00158     for (unsigned igrd = 0; igrd < grid.size(); igrd++) {
00159       if (grid[igrd] != -1) {
00160         inputs.push_back(cs.jets()[grid[igrd]]);
00161         cs_indices.push_back(grid[igrd]);
00162       }
00163     }
00164     ClusterSequence post_cs(inputs, _post_jet_def);
00165     const vector<ClusterSequence::history_element> & post_history = post_cs.history();
00166     const vector<PseudoJet>                        & post_jets = post_cs.jets();
00167     for (unsigned ihist = ngrid_active; ihist < post_history.size(); ihist++) {
00168       const ClusterSequence::history_element & hist = post_history[ihist];
00169       int post_ij1 = post_history[hist.parent1].jetp_index;
00170       int ij1 = cs_indices[post_ij1];
00171       if (hist.parent2 >= 0) {
00172         int post_ij2 = post_history[hist.parent2].jetp_index;
00173         int ij2 = cs_indices[post_ij2];
00174         int k;
00175         cs.plugin_record_ij_recombination(ij1, ij2, hist.dij, post_jets[hist.jetp_index], k); 
00176         assert(int(cs_indices.size()) == hist.jetp_index);
00177         cs_indices.push_back(k);
00178       } else {
00179         cs.plugin_record_iB_recombination(ij1, hist.dij);
00180       }
00181     }
00182     
00183   }
00184 }
00185 
00186 FASTJET_END_NAMESPACE      // defined in fastjet/internal/base.hh
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends