FastJet 3.0beta1
fastjet_subtraction.cc
00001 
00002 //STARTHEADER
00003 // $Id: fastjet_subtraction.cc 2444 2011-07-21 19:07:31Z salam $
00004 //
00005 // Copyright (c) 2005-2006, Matteo Cacciari and Gavin Salam
00006 //
00007 //----------------------------------------------------------------------
00008 // This file is part of FastJet.
00009 //
00010 //  FastJet is free software; you can redistribute it and/or modify
00011 //  it under the terms of the GNU General Public License as published by
00012 //  the Free Software Foundation; either version 2 of the License, or
00013 //  (at your option) any later version.
00014 //
00015 //  The algorithms that underlie FastJet have required considerable
00016 //  development and are described in hep-ph/0512210. If you use
00017 //  FastJet as part of work towards a scientific publication, please
00018 //  include a citation to the FastJet paper.
00019 //
00020 //  FastJet is distributed in the hope that it will be useful,
00021 //  but WITHOUT ANY WARRANTY; without even the implied warranty of
00022 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00023 //  GNU General Public License for more details.
00024 //
00025 //  You should have received a copy of the GNU General Public License
00026 //  along with FastJet; if not, write to the Free Software
00027 //  Foundation, Inc.:
00028 //      59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00029 //----------------------------------------------------------------------
00030 //ENDHEADER
00031 
00032 
00033 
00034 
00035 
00036 
00037 
00038 //----------------------------------------------------------------------
00039 // fastjet subtraction example program. 
00040 //
00041 // Compile it with: make fastjet_subtraction
00042 // run it with    : ./fastjet_subtraction < data/Pythia-dijet-ptmin100-lhc-pileup-1ev.dat
00043 //
00044 //----------------------------------------------------------------------
00045 #include "fastjet/PseudoJet.hh"
00046 #include "fastjet/ClusterSequenceArea.hh"
00047 #include<iostream> // needed for io
00048 #include<sstream>  // needed for internal io
00049 #include<vector> 
00050 
00051 using namespace std;
00052 
00053 // A declaration of a function that pretty prints a list of jets
00054 // The subtraction is also performed inside this function
00055 void print_jets (const fastjet::ClusterSequenceAreaBase &, 
00056                  const vector<fastjet::PseudoJet> &);
00057 
00058 /// an example program showing how to use fastjet
00059 int main (int argc, char ** argv) {
00060   
00061   // since we use here simulated data we can split the hard event
00062   // from the full (i.e. with pileup added) one
00063   vector<fastjet::PseudoJet> hard_event, full_event;
00064   
00065   // maximum rapidity that we accept for the input particles
00066   double etamax = 6.0;
00067   
00068   // read in input particles. Keep the hard event generated by PYTHIA
00069   // separated from the full event, so as to be able to gauge the
00070   // "goodness" of the subtraction from the full event which also
00071   // includes pileup
00072   string line;
00073   int  nsub  = 0;
00074   while (getline(cin, line)) {
00075     istringstream linestream(line);
00076     // take substrings to avoid problems when there are extra "pollution"
00077     // characters (e.g. line-feed).
00078     if (line.substr(0,4) == "#END") {break;}
00079     if (line.substr(0,9) == "#SUBSTART") {
00080       // if more sub events follow, make copy of hard one here
00081       if (nsub == 1) hard_event = full_event;
00082       nsub += 1;
00083     }
00084     if (line.substr(0,1) == "#") {continue;}
00085     valarray<double> fourvec(4);
00086     linestream >> fourvec[0] >> fourvec[1] >> fourvec[2] >> fourvec[3];
00087     fastjet::PseudoJet psjet(fourvec);
00088     psjet.set_user_index(0);
00089     // push event onto back of full_event vector
00090     if (abs(psjet.rap()) < etamax) {full_event.push_back(psjet);}
00091   }
00092 
00093   // if we have read in only one event, copy it across here...
00094   if (nsub == 1) hard_event = full_event;
00095 
00096   // if there was nothing in the event 
00097   if (nsub == 0) {
00098     cerr << "Error: read empty event\n";
00099     exit(-1);
00100   }
00101   
00102   
00103   // create an object that represents your choice of jet algorithm and 
00104   // the associated parameters
00105   double R = 0.7;
00106   fastjet::Strategy strategy = fastjet::Best;
00107   fastjet::JetDefinition jet_def(fastjet::kt_algorithm, R, strategy);
00108   //fastjet::JetDefinition jet_def(fastjet::cambridge_algorithm, R, strategy);
00109   //fastjet::JetDefinition jet_def(fastjet::antikt_algorithm, R, strategy);
00110 
00111   // create an object that specifies how we to define the (active) area
00112   double ghost_etamax = 6.0;
00113   int    active_area_repeats = 1;
00114   double ghost_area    = 0.01;
00115   fastjet::GhostedAreaSpec ghost_spec(ghost_etamax, active_area_repeats, 
00116                                     ghost_area);
00117   fastjet::AreaDefinition area_def(ghost_spec);
00118   //fastjet::AreaDefinition area_def(fastjet::VoronoiAreaSpec(1.0));
00119 
00120   // run the jet clustering with the above jet definition. hard event first
00121   fastjet::ClusterSequenceArea clust_seq(hard_event, jet_def, area_def);
00122   //fastjet::ClusterSequencePassiveArea clust_seq(hard_event, jet_def);
00123 
00124 
00125   // extract the inclusive jets with pt > 5 GeV, sorted by pt
00126   double ptmin = 5.0;
00127   vector<fastjet::PseudoJet> inclusive_jets = clust_seq.inclusive_jets(ptmin);
00128 
00129   // print them out
00130   cout << "" <<endl;
00131   cout << "Hard event only"<<endl;
00132   cout << "Number of input particles: "<<hard_event.size()<<endl;
00133   cout << "Strategy used: "<<clust_seq.strategy_string()<<endl;
00134   cout << "Printing inclusive jets with pt > "<< ptmin<<" GeV\n";
00135   cout << "---------------------------------------\n";
00136   print_jets(clust_seq, inclusive_jets);
00137   cout << endl;
00138 
00139 
00140   // repeat everything on the full event
00141 
00142   // run the jet clustering with the above jet definition
00143   fastjet::ClusterSequenceArea clust_seq_full(full_event, jet_def, area_def);
00144 
00145   // extract the inclusive jets with pt > 20 GeV, sorted by pt
00146   ptmin = 20.0;
00147   inclusive_jets = clust_seq_full.inclusive_jets(ptmin);
00148 
00149   // print them out
00150   cout << "Full event, with pileup, and its subtraction"<<endl;
00151   cout << "Number of input particles: "<<full_event.size()<<endl;
00152   cout << "Strategy used: "<<clust_seq_full.strategy_string()<<endl;
00153   cout << "Printing inclusive jets with pt > "<< ptmin<<" GeV (before subtraction)\n";
00154   cout << "---------------------------------------\n";
00155   print_jets(clust_seq_full, inclusive_jets);
00156   cout << endl;
00157 
00158 
00159 }
00160 
00161 
00162 //----------------------------------------------------------------------
00163 /// a function that pretty prints a list of jets, and performs the subtraction
00164 /// in two different ways, using a generic ClusterSequenceAreaBase
00165 /// type object.
00166 void print_jets (const fastjet::ClusterSequenceAreaBase & clust_seq, 
00167                  const vector<fastjet::PseudoJet> & unsorted_jets ) {
00168 
00169   // sort jets into increasing pt
00170   vector<fastjet::PseudoJet> jets = sorted_by_pt(unsorted_jets);  
00171 
00172   // the corrected jets will go in here
00173   vector<fastjet::PseudoJet> corrected_jets(jets.size());
00174   
00175   // get median pt per unit area -- you must decide the rapidity range
00176   // in which you measure the activity. (If doing a ghost-based active
00177   // area, make sure that your ghosts go up at least to \sim range+R).
00178   double range = 5.0;
00179   double median_pt_per_area = clust_seq.median_pt_per_unit_area(range);
00180   double median_pt_per_area4vector = clust_seq.median_pt_per_unit_area_4vector(range);
00181 
00182   printf(" ijet     rap     phi        Pt    area  Pt corr  (rap corr phi corr Pt corr)ext\n");
00183   for (size_t j = 0; j < jets.size(); j++) {
00184 
00185     // get area of each jet
00186     double area     = jets[j].area();
00187 
00188     // "standard" correction. Subtract only the Pt
00189     double pt_corr  = jets[j].perp() - area*median_pt_per_area;
00190 
00191     // "extended" correction
00192     fastjet::PseudoJet sub_4vect = 
00193                        median_pt_per_area4vector*jets[j].area_4vector();
00194     if (sub_4vect.perp2() >= jets[j].perp2() || 
00195         sub_4vect.E()     >= jets[j].E()) {
00196       // if the correction is too large, set the jet to zero
00197       corrected_jets[j] =  0.0 * jets[j];
00198     } else {
00199       // otherwise do an E-scheme subtraction
00200       corrected_jets[j] = jets[j] - sub_4vect;
00201     }
00202     // NB We could also leave out the above "if": 
00203     // corrected_jets[j] = jets[j] - sub_4vect;
00204     // but the result would be different, since we would not avoid
00205     // jets with negative Pt or energy
00206     
00207     printf("%5u %7.3f %7.3f %9.3f %7.3f %9.3f %7.3f %7.3f %9.3f\n",
00208      j,jets[j].rap(),jets[j].phi(),jets[j].perp(), area, pt_corr,
00209      corrected_jets[j].rap(),corrected_jets[j].phi(), corrected_jets[j].perp());
00210   }
00211 
00212   cout << endl;
00213   cout << "median pt_over_area = " << median_pt_per_area << endl;
00214   cout << "median pt_over_area4vector = " << median_pt_per_area4vector << endl << endl;
00215 
00216 
00217 }
00218 
00219 
00220 
00221 
00222 
00223 
00224 
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends