07-subtraction.cc

Go to the documentation of this file.
00001 //----------------------------------------------------------------------
00002 /// \file
00003 /// \page Example07 07 - subtracting jet background contamination
00004 ///
00005 /// fastjet subtraction example program. 
00006 ///
00007 /// run it with    : ./07-subtraction < data/Pythia-Zp2jets-lhc-pileup-1ev.dat
00008 ///
00009 /// Source code: 07-subtraction.cc
00010 //----------------------------------------------------------------------
00011 
00012 //STARTHEADER
00013 // $Id: 07-subtraction.cc 1911 2011-01-28 18:18:24Z soyez $
00014 //
00015 // Copyright (c) 2005-2011, Matteo Cacciari, Gavin Salam and Gregory Soyez
00016 //
00017 //----------------------------------------------------------------------
00018 // This file is part of FastJet.
00019 //
00020 //  FastJet is free software; you can redistribute it and/or modify
00021 //  it under the terms of the GNU General Public License as published by
00022 //  the Free Software Foundation; either version 2 of the License, or
00023 //  (at your option) any later version.
00024 //
00025 //  The algorithms that underlie FastJet have required considerable
00026 //  development and are described in hep-ph/0512210. If you use
00027 //  FastJet as part of work towards a scientific publication, please
00028 //  include a citation to the FastJet paper.
00029 //
00030 //  FastJet is distributed in the hope that it will be useful,
00031 //  but WITHOUT ANY WARRANTY; without even the implied warranty of
00032 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00033 //  GNU General Public License for more details.
00034 //
00035 //  You should have received a copy of the GNU General Public License
00036 //  along with FastJet; if not, write to the Free Software
00037 //  Foundation, Inc.:
00038 //      59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00039 //----------------------------------------------------------------------
00040 //ENDHEADER
00041 
00042 #include "fastjet/PseudoJet.hh"
00043 #include "fastjet/ClusterSequenceArea.hh"
00044 #include <iostream> // needed for io
00045 
00046 using namespace std;
00047 
00048 int main (int argc, char ** argv) {
00049   
00050   // read in input particles
00051   //
00052   // since we use here simulated data we can split the hard event
00053   // from the full (i.e. with pileup added) one
00054   //----------------------------------------------------------
00055 
00056   vector<fastjet::PseudoJet> hard_event, full_event;
00057   
00058   // read in input particles. Keep the hard event generated by PYTHIA
00059   // separated from the full event, so as to be able to gauge the
00060   // "goodness" of the subtraction from the full event, which also
00061   // includes pileup
00062   double particle_maxrap = 5.0;
00063 
00064   string line;
00065   int  nsub  = 0; // counter to keep track of which sub-event we're reading
00066   while (getline(cin, line)) {
00067     istringstream linestream(line);
00068     // take substrings to avoid problems when there are extra "pollution"
00069     // characters (e.g. line-feed).
00070     if (line.substr(0,4) == "#END") {break;}
00071     if (line.substr(0,9) == "#SUBSTART") {
00072       // if more sub events follow, make copy of first one (the hard one) here
00073       if (nsub == 1) hard_event = full_event;
00074       nsub += 1;
00075     }
00076     if (line.substr(0,1) == "#") {continue;}
00077     double px,py,pz,E;
00078     linestream >> px >> py >> pz >> E;
00079     // you can construct 
00080     fastjet::PseudoJet particle(px,py,pz,E);
00081 
00082     // push event onto back of full_event vector
00083     if (abs(particle.rap()) <= particle_maxrap) full_event.push_back(particle);
00084   }
00085 
00086   // if we have read in only one event, copy it across here...
00087   if (nsub == 1) hard_event = full_event;
00088 
00089   // if there was nothing in the event 
00090   if (nsub == 0) {
00091     cerr << "Error: read empty event\n";
00092     exit(-1);
00093   }
00094   
00095   
00096   // create a jet definition for the clustering
00097   // We use the anti-kt algorithm with a radius of 0.5
00098   //----------------------------------------------------------
00099   double R = 0.5;
00100   fastjet::JetDefinition jet_def(fastjet::antikt_algorithm, R);
00101 
00102   // create an area definition for the clustering
00103   //----------------------------------------------------------
00104   double ghost_maxrap = 6.0;
00105   fastjet::ActiveAreaSpec area_spec(ghost_maxrap);
00106   fastjet::AreaDefinition area_def(fastjet::active_area, area_spec);
00107 
00108   // run the jet clustering with the above jet and area definitions
00109   // for both the hard and full event
00110   //
00111   // We retrieve the jets above 7 GeV in both case (note that the
00112   // 7-GeV cut we be applied again later on after we subtract the jets
00113   // from the full event)
00114   // ----------------------------------------------------------
00115   fastjet::ClusterSequenceArea clust_seq_hard(hard_event, jet_def, area_def);
00116   fastjet::ClusterSequenceArea clust_seq_full(full_event, jet_def, area_def);
00117 
00118   double ptmin = 7.0;
00119   vector<fastjet::PseudoJet> hard_jets = sorted_by_pt(clust_seq_hard.inclusive_jets(ptmin));
00120   vector<fastjet::PseudoJet> full_jets = sorted_by_pt(clust_seq_full.inclusive_jets(ptmin));
00121 
00122   // Now turn to the estimation of the background (for the full event)
00123   //
00124   // This also requires a ClusterSequenceArea.
00125   // In general, this ClusterSequenceArea does not need to be the same
00126   // as the one used (above) to cluster and extract the jets from the
00127   // event:
00128   //  - We strongly recommend using the kt or Cambridge/Aachen algorithm
00129   //    (a warning will be issued otherwise)
00130   //  - The choice of the radius is a bit more subtle. R=0.6 is a
00131   //    decent default, though for dense events, one can also
00132   //    decrease that value to 0.4-0.5.
00133   //  - For the area definition, we recommend the use of explicit
00134   //    ghosts (i.e. active_area_explicit_ghosts)
00135   //    As mentionned in the area example (06-area.cc), ghosts should
00136   //    extend sufficiently far in rapidity to cover the jets used in
00137   //    the computation of the background (see also the comment below)
00138   //
00139   // ----------------------------------------------------------
00140   fastjet::JetDefinition jet_def_bkgd(fastjet::kt_algorithm, 0.6);
00141   fastjet::GhostedAreaSpec area_spec_bkgd(ghost_maxrap);
00142   fastjet::AreaDefinition area_def_bkgd(fastjet::active_area_explicit_ghosts, area_spec_bkgd);
00143   fastjet::ClusterSequenceArea clust_seq_bkgd(full_event, jet_def_bkgd, area_def_bkgd);
00144 
00145   // Once you have the ClusterSequenceArea, you can compute the
00146   // background. This is estimated over a given range
00147   // (RangeDefinition) i.e. only jets within that range will be used
00148   // to estimate the background. You shold thus make sure the ghosts
00149   // extend far enough in rapidity to cover the range, a warning will
00150   // be issued otherwise.
00151   //
00152   // The simplest way to define a RangeDefinition is through its
00153   // maximal |y| extent but other options are possible e.g. through a
00154   // minimal and maximal rapidity and minimal and maximal azimuthal
00155   // angle. If needed, you can even define your own ranges (a few are
00156   // provided with FastJet)
00157   //
00158   // Finally, the estimation of the background properties rho (the
00159   // average density per unit area) and sigma (the average
00160   // fluctuations per unit area) is done using
00161   // ClusterSequenceArea::get_median_rho_and_sigma(). This takes
00162   // 2 main parameters: the range discussed above and a boolean
00163   // controlling the use of 4-vector or scalar areas (we suggest using
00164   // 4-vector areas)
00165   //
00166   // ----------------------------------------------------------
00167   double range_maxrap = 5.0;  // we have a ghost_maxrap of 6.0
00168   fastjet::RangeDefinition range(range_maxrap);
00169 
00170   bool use_4vector_area = true;
00171   
00172   double rho, sigma;
00173   clust_seq_bkgd.get_median_rho_and_sigma(range, use_4vector_area, rho, sigma);
00174 
00175   // show a summary of what was done so far
00176   //  - the description of the algorithms, areas and ranges used
00177   //  - the background properties
00178   //  - the jets in the hard event
00179   //----------------------------------------------------------
00180   cout << "Main clustering:" << endl;
00181   cout << "  Ran:   " << jet_def.description() << endl;
00182   cout << "  Area:  " << area_def.description() << endl;
00183   cout << "  Particles up to |y|=" << particle_maxrap << endl;
00184   cout << endl;
00185 
00186   cout << "Background estimation:" << endl;
00187   cout << "  Ran    " << jet_def_bkgd.description() << endl;
00188   cout << "  Area:  " << area_def_bkgd.description() << endl;
00189   cout << "  Range: " << range.description() << endl;
00190   cout << "  Giving, for the full event" << endl;
00191   cout << "    rho   = " << rho   << endl;
00192   cout << "    sigma = " << sigma << endl;
00193   cout << endl;
00194 
00195   cout << "Jets above " << ptmin << " GeV in the hard event (" << hard_event.size() << " particles)" << endl;
00196   cout << "---------------------------------------\n";
00197   printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt", "area");
00198    for (unsigned int i = 0; i < hard_jets.size(); i++) {
00199     printf("%5u %15.8f %15.8f %15.8f %15.8f\n", i,
00200            hard_jets[i].rap(), hard_jets[i].phi(), hard_jets[i].perp(),
00201            hard_jets[i].area());
00202   }
00203   cout << endl;
00204 
00205   // Once the background properties have been computed, subtraction
00206   // can be applied on the jets
00207   //
00208   // This uses ClusterSequenceArea::subtracted_jet(jet, rho), with the
00209   // ClusterSequence used to cluster the jet and the background
00210   // density we have just computed
00211   // 
00212   // (Note that when using scalar areas, subtracted_pt should be used
00213   // instead of subtracted_jet)
00214   // 
00215   // We output the jets before and after subtraction
00216   // ----------------------------------------------------------
00217   cout << "Jets above " << ptmin << " GeV in the full event (" << full_event.size() << " particles)" << endl;
00218   cout << "---------------------------------------\n";
00219   printf("%5s %15s %15s %15s %15s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt", "area", "rap_sub", "phi_sub", "pt_sub");
00220   unsigned int idx=0;
00221   for (unsigned int i=0; i<full_jets.size(); i++){
00222     // get the subtracted jet
00223     fastjet::PseudoJet subtracted_jet = clust_seq_full.subtracted_jet(full_jets[i], rho);
00224 
00225     // re-apply the pt cut
00226     if (subtracted_jet.perp2() >= ptmin*ptmin){
00227       printf("%5u %15.8f %15.8f %15.8f %15.8f %15.8f %15.8f %15.8f\n", idx,
00228              full_jets[i].rap(), full_jets[i].phi(), full_jets[i].perp(),
00229              full_jets[i].area(),
00230              subtracted_jet.rap(), subtracted_jet.phi(), 
00231              subtracted_jet.perp());
00232       idx++;
00233     }
00234   }
00235 
00236   return 0;
00237 }