FastJet 3.0beta1
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 2472 2011-07-26 11:25:42Z cacciari $
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 "fastjet/Selector.hh"
00045 #include "fastjet/tools/JetMedianBackgroundEstimator.hh"
00046 #include "fastjet/tools/Subtractor.hh" 
00047 #include <iostream> // needed for io
00048 
00049 using namespace std;
00050 using namespace fastjet;
00051 
00052 int main (int argc, char ** argv) {
00053   
00054   // read in input particles
00055   //
00056   // since we use here simulated data we can split the hard event
00057   // from the full (i.e. with pileup added) one
00058   //----------------------------------------------------------
00059 
00060   vector<PseudoJet> hard_event, full_event;
00061   
00062   // read in input particles. Keep the hard event generated by PYTHIA
00063   // separated from the full event, so as to be able to gauge the
00064   // "goodness" of the subtraction from the full event, which also
00065   // includes pileup
00066   double particle_maxrap = 5.0;
00067 
00068   string line;
00069   int  nsub  = 0; // counter to keep track of which sub-event we're reading
00070   while (getline(cin, line)) {
00071     istringstream linestream(line);
00072     // take substrings to avoid problems when there are extra "pollution"
00073     // characters (e.g. line-feed).
00074     if (line.substr(0,4) == "#END") {break;}
00075     if (line.substr(0,9) == "#SUBSTART") {
00076       // if more sub events follow, make copy of first one (the hard one) here
00077       if (nsub == 1) hard_event = full_event;
00078       nsub += 1;
00079     }
00080     if (line.substr(0,1) == "#") {continue;}
00081     double px,py,pz,E;
00082     linestream >> px >> py >> pz >> E;
00083     // you can construct 
00084     PseudoJet particle(px,py,pz,E);
00085 
00086     // push event onto back of full_event vector
00087     if (abs(particle.rap()) <= particle_maxrap) full_event.push_back(particle);
00088   }
00089 
00090   // if we have read in only one event, copy it across here...
00091   if (nsub == 1) hard_event = full_event;
00092 
00093   // if there was nothing in the event 
00094   if (nsub == 0) {
00095     cerr << "Error: read empty event\n";
00096     exit(-1);
00097   }
00098   
00099   
00100   // create a jet definition for the clustering
00101   // We use the anti-kt algorithm with a radius of 0.5
00102   //----------------------------------------------------------
00103   double R = 0.5;
00104   JetDefinition jet_def(antikt_algorithm, R);
00105 
00106   // create an area definition for the clustering
00107   //----------------------------------------------------------
00108   // ghosts should go up to the acceptance of the detector or
00109   // (with infinite acceptance) at least 2R beyond the region
00110   // where you plan to investigate jets.
00111   double ghost_maxrap = 6.0;
00112   GhostedAreaSpec area_spec(ghost_maxrap);
00113   AreaDefinition area_def(active_area, area_spec);
00114 
00115   // run the jet clustering with the above jet and area definitions
00116   // for both the hard and full event
00117   //
00118   // We retrieve the jets above 7 GeV in both case (note that the
00119   // 7-GeV cut we be applied again later on after we subtract the jets
00120   // from the full event)
00121   // ----------------------------------------------------------
00122   ClusterSequenceArea clust_seq_hard(hard_event, jet_def, area_def);
00123   ClusterSequenceArea clust_seq_full(full_event, jet_def, area_def);
00124 
00125   double ptmin = 7.0;
00126   vector<PseudoJet> hard_jets = sorted_by_pt(clust_seq_hard.inclusive_jets(ptmin));
00127   vector<PseudoJet> full_jets = sorted_by_pt(clust_seq_full.inclusive_jets(ptmin));
00128 
00129   // Now turn to the estimation of the background (for the full event)
00130   //
00131   // There are different ways to do that. In general, this also
00132   // requires clustering the particles that will be handled internally
00133   // in FastJet. 
00134   //
00135   // The suggested way to proceed is to use a BackgroundEstimator
00136   // constructed from the following 3 arguments:
00137   //  - a jet definition used to cluster the particles.
00138   //    . We strongly recommend using the kt or Cambridge/Aachen
00139   //      algorithm (a warning will be issued otherwise)
00140   //    . The choice of the radius is a bit more subtle. R=0.4 has
00141   //      been chosen to limit the impact of hard jets; in samples of
00142   //      dominantly sparse events it may cause the UE/pileup to be
00143   //      underestimated a little, a slightly larger value (0.5 or
00144   //      0.6) may be better.
00145   //  - An area definition for which we recommend the use of explicit
00146   //    ghosts (i.e. active_area_explicit_ghosts)
00147   //    As mentionned in the area example (06-area.cc), ghosts should
00148   //    extend sufficiently far in rapidity to cover the jets used in
00149   //    the computation of the background (see also the comment below)
00150   //  - A Selector specifying the range over which we will keep the
00151   //    jets entering the estimation of the background (you should
00152   //    thus make sure the ghosts extend far enough in rapidity to
00153   //    cover the range, a warning will be issued otherwise).
00154   //    In this particular example, the two hardest jets in the event
00155   //    are removed from the background estimation
00156   // ----------------------------------------------------------
00157   JetDefinition jet_def_bkgd(kt_algorithm, 0.4);
00158   AreaDefinition area_def_bkgd(active_area_explicit_ghosts, 
00159                                GhostedAreaSpec(ghost_maxrap));
00160   Selector selector = SelectorAbsRapMax(4.5) * (!SelectorNHardest(2));
00161   JetMedianBackgroundEstimator bkgd_estimator(selector, jet_def_bkgd, area_def_bkgd);
00162 
00163   // To help manipulate the background estimator, we also provide a
00164   // transformer that allows to apply directly the background
00165   // subtraction on the jets. This will use the background estimator
00166   // to compute rho for the jets to be subtracted.
00167   // ----------------------------------------------------------
00168   Subtractor subtractor(&bkgd_estimator);
00169 
00170   // Finally, once we have an event, we can just tell the background
00171   // estimator to use that list of particles
00172   // This could be done directly when declaring the background
00173   // estimator but the usage below can more easily be accomodated to a
00174   // loop over a set of events.
00175   // ----------------------------------------------------------
00176   bkgd_estimator.set_particles(full_event);
00177 
00178   // show a summary of what was done so far
00179   //  - the description of the algorithms, areas and ranges used
00180   //  - the background properties
00181   //  - the jets in the hard event
00182   //----------------------------------------------------------
00183   cout << "Main clustering:" << endl;
00184   cout << "  Ran:   " << jet_def.description() << endl;
00185   cout << "  Area:  " << area_def.description() << endl;
00186   cout << "  Particles up to |y|=" << particle_maxrap << endl;
00187   cout << endl;
00188 
00189   cout << "Background estimation:" << endl;
00190   cout << "  " << bkgd_estimator.description() << endl << endl;;
00191   cout << "  Giving, for the full event" << endl;
00192   cout << "    rho   = " << bkgd_estimator.rho()   << endl;
00193   cout << "    sigma = " << bkgd_estimator.sigma() << endl;
00194   cout << endl;
00195 
00196   cout << "Jets above " << ptmin << " GeV in the hard event (" << hard_event.size() << " particles)" << endl;
00197   cout << "---------------------------------------\n";
00198   printf("%5s %15s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt", "area");
00199    for (unsigned int i = 0; i < hard_jets.size(); i++) {
00200     printf("%5u %15.8f %15.8f %15.8f %15.8f\n", i,
00201            hard_jets[i].rap(), hard_jets[i].phi(), hard_jets[i].perp(),
00202            hard_jets[i].area());
00203   }
00204   cout << endl;
00205 
00206   // Once the background properties have been computed, subtraction
00207   // can be applied on the jets. Subtraction is performed on the
00208   // full 4-vector
00209   //
00210   // We output the jets before and after subtraction
00211   // ----------------------------------------------------------
00212   cout << "Jets above " << ptmin << " GeV in the full event (" << full_event.size() << " particles)" << endl;
00213   cout << "---------------------------------------\n";
00214   printf("%5s %15s %15s %15s %15s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt", "area", "rap_sub", "phi_sub", "pt_sub");
00215   unsigned int idx=0;
00216 
00217   // get the subtracted jets
00218   vector<PseudoJet> subtracted_jets = subtractor(full_jets);
00219 
00220   for (unsigned int i=0; i<full_jets.size(); i++){
00221     // re-apply the pt cut
00222     if (subtracted_jets[i].perp2() >= ptmin*ptmin){
00223       printf("%5u %15.8f %15.8f %15.8f %15.8f %15.8f %15.8f %15.8f\n", idx,
00224              full_jets[i].rap(), full_jets[i].phi(), full_jets[i].perp(),
00225              full_jets[i].area(),
00226              subtracted_jets[i].rap(), subtracted_jets[i].phi(), 
00227              subtracted_jets[i].perp());
00228       idx++;
00229     }
00230   }
00231 
00232   return 0;
00233 }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends