11-boosted_higgs.cc

Go to the documentation of this file.
00001 //----------------------------------------------------------------------
00002 /// \file
00003 /// \page Example11 11 - boosted Higgs tagging
00004 ///
00005 /// fastjet example program, illustration of carrying out boosted
00006 /// Higgs subjet ID analysis
00007 ///
00008 /// It illustrates two kinds of functionality: 
00009 ///
00010 ///  - following the decomposition of a jet into pieces
00011 ///  - following information on a b-tag through the jet
00012 ///
00013 /// This kind of functionality was used in arXiv:0802.2470
00014 /// (Butterworth, Davison, Rubin & Salam) for boosted Higgs searches,
00015 /// and related functionality was used in arXiv:0806.0848 (Kaplan,
00016 /// Rehermann, Schwartz & Tweedie) in searching for boosted tops
00017 /// (without b-tag assumptions).
00018 ///
00019 /// run it with    : ./11-boosted_higgs < data/HZ-event-Hmass115.dat
00020 ///
00021 /// Source code: 11-boosted_higgs.cc
00022 //----------------------------------------------------------------------
00023 
00024 
00025 //STARTHEADER
00026 // $Id: 11-boosted_higgs.cc 1911 2011-01-28 18:18:24Z soyez $
00027 //
00028 // Copyright (c) 2005-2011, Matteo Cacciari, Gavin Salam and Gregory Soyez
00029 //
00030 //----------------------------------------------------------------------
00031 // This file is part of FastJet.
00032 //
00033 //  FastJet is free software; you can redistribute it and/or modify
00034 //  it under the terms of the GNU General Public License as published by
00035 //  the Free Software Foundation; either version 2 of the License, or
00036 //  (at your option) any later version.
00037 //
00038 //  The algorithms that underlie FastJet have required considerable
00039 //  development and are described in hep-ph/0512210. If you use
00040 //  FastJet as part of work towards a scientific publication, please
00041 //  include a citation to the FastJet paper.
00042 //
00043 //  FastJet is distributed in the hope that it will be useful,
00044 //  but WITHOUT ANY WARRANTY; without even the implied warranty of
00045 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00046 //  GNU General Public License for more details.
00047 //
00048 //  You should have received a copy of the GNU General Public License
00049 //  along with FastJet; if not, write to the Free Software
00050 //  Foundation, Inc.:
00051 //      59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00052 //----------------------------------------------------------------------
00053 //ENDHEADER
00054 
00055 #include "fastjet/ClusterSequence.hh"
00056 #include <iostream> // needed for io
00057 #include <sstream>  // needed for internal io
00058 #include <iomanip>  
00059 #include <cmath>
00060 
00061 using namespace std;
00062 using namespace fastjet;
00063 
00064 
00065 //----------------------------------------------------------------------
00066 // set up a class to give standard (by default E-scheme)
00067 // recombination, with additional tracking of flavour information in
00068 // the user_index. 
00069 //
00070 // If you use this, you must explicitly set the user index to 0 for
00071 // non-flavoured particles (the default value is -1);
00072 //
00073 // This will work for native algorithms, but not for all plugins
00074 //----------------------------------------------------------------------
00075 typedef JetDefinition::DefaultRecombiner DefRecomb;
00076 
00077 class FlavourRecombiner : public  DefRecomb {
00078 public:
00079   FlavourRecombiner(RecombinationScheme recomb_scheme = E_scheme) : 
00080     DefRecomb(recomb_scheme) {};
00081 
00082   virtual std::string description() const {
00083     return DefRecomb::description()+" (with user index addition)";}
00084 
00085   /// recombine pa and pb and put result into pab
00086   virtual void recombine(const PseudoJet & pa, const PseudoJet & pb, 
00087                          PseudoJet & pab) const {
00088     DefRecomb::recombine(pa,pb,pab);
00089     pab.set_user_index(pa.user_index() + pb.user_index());
00090   }
00091 };
00092 
00093 
00094 //----------------------------------------------------------------------
00095 // forward declaration for printing out info about a jet
00096 //----------------------------------------------------------------------
00097 ostream & operator<<(ostream &, PseudoJet &);
00098 
00099 
00100 //----------------------------------------------------------------------
00101 // core of the program
00102 //----------------------------------------------------------------------
00103 int main (int argc, char ** argv) {
00104 
00105   vector<PseudoJet> particles;
00106 
00107   // read in data in format px py pz E b-tag [last of these is optional]
00108   // lines starting with "#" are considered as comments and discarded
00109   //----------------------------------------------------------
00110 
00111   string line;
00112   while (getline(cin,line)) {
00113     if (line.substr(0,1) == "#") {continue;}
00114     istringstream linestream(line);
00115     double px,py,pz,E;
00116     linestream >> px >> py >> pz >> E;
00117 
00118     // optionally read in btag information
00119     int    btag;
00120     if (! (linestream >> btag)) btag = 0;
00121 
00122     // construct the particle
00123     PseudoJet particle(px,py,pz,E);
00124     particle.set_user_index(btag); // btag info goes in user index, for flavour tracking
00125     particles.push_back(particle);
00126   }
00127 
00128 
00129   // set up the jet finding
00130   //
00131   // This also shows how to use the "FlavourRecombiner" user-defined
00132   // recombiner 
00133   // ----------------------------------------------------------
00134   double R = 1.2;
00135   FlavourRecombiner flav_recombiner; // for tracking flavour
00136   JetDefinition jet_def(cambridge_algorithm, R, &flav_recombiner);
00137 
00138   
00139   // run the jet finding; find the hardest jet
00140   ClusterSequence cs(particles, jet_def);
00141   vector<PseudoJet> jets = sorted_by_pt(cs.inclusive_jets());
00142 
00143   cout << "Ran: " << jet_def.description() << endl << endl;
00144   cout << "Hardest jet: " << jets[0] << endl << endl;
00145 
00146   // now do the subjet decomposition
00147   //----------------------------------------------------------
00148   //
00149   // when unpeeling a C/A jet, often only a very soft piece may break off;
00150   // the mass_drop_threshold indicates how much "lighter" the heavier of the two
00151   // resulting pieces must be in order for us to consider that we've really
00152   // seen some form of substructure
00153   double mass_drop_threshold = 0.667; 
00154   // QCD backgrounds that give larger jet masses have a component
00155   // where a quite soft gluon is emitted; to eliminate part of this
00156   // one can place a cut on the asymmetry of the branching; 
00157   //
00158   // Here the cut is expressed in terms of y, the kt-distance scaled
00159   // to the squared jet mass; an easier way to see it is in terms of
00160   // a requirement on the momentum fraction in the splitting: z/(1-z)
00161   // and (1-z)/z > rtycut^2 [the correspondence holds only at LO]
00162   double rtycut              = 0.3;
00163 
00164   PseudoJet this_jet = jets[0];
00165   PseudoJet parent1, parent2;
00166   bool had_parents;
00167 
00168   while ((had_parents = this_jet.has_parents(parent1,parent2))) {
00169     // make parent1 the more massive jet
00170     if (parent1.m() < parent2.m()) swap(parent1,parent2);
00171 
00172     // if we pass the conditions on the mass drop and its degree of
00173     // asymmetry (z/(1-z) \sim kt_dist/m^2 > rtycut), then we've found
00174     // something interesting, so exit the loop
00175     if (parent1.m() < mass_drop_threshold * this_jet.m() &&
00176         parent1.kt_distance(parent2) > pow(rtycut,2) * this_jet.m2()) {
00177       break;
00178     } else {
00179       // otherwise try a futher decomposition on the more massive jet
00180       this_jet = parent1;
00181     }
00182   }
00183 
00184   // look to see what we found
00185   if (!had_parents) {
00186     cout << "Did not find suitable hard substructure in this event." << endl;
00187     return 0;
00188   }
00189 
00190   cout << "Found suitable pair of subjets: " << endl;
00191   cout << " " << parent1 << endl;
00192   cout << " " << parent2 << endl;
00193   cout << "Total = " << endl;
00194   cout << " " << this_jet << endl << endl;
00195 
00196   // next we "filter" it, to remove UE & pileup contamination
00197   //----------------------------------------------------------
00198   //
00199   // [there are two ways of doing this; here we directly use the
00200   // exsiting cluster sequence and find the exclusive subjets of
00201   // this_jet (i.e. work backwards within the cs starting from
00202   // this_jet); alternatively one can recluster just the
00203   // constituents of the jet]
00204   //
00205   // first get separation between the subjets (called Rbb -- assuming it's a Higgs!)
00206   double   Rbb = sqrt(parent1.squared_distance(parent2));
00207   double   Rfilt = min(Rbb/2, 0.3); // somewhat arbitrary choice
00208   unsigned nfilt = 3;               // number of pieces we'll take
00209   cout << "Subjet separation (Rbb) = " << Rbb << ", Rfilt = " << Rfilt << endl;
00210 
00211   double   dcut  = pow(Rfilt/R,2);  // for C/A get a view at Rfilt by
00212   // using a dcut=(Rfilt/R)^2
00213   vector<PseudoJet> filt_subjets = sorted_by_pt(this_jet.exclusive_subjets(dcut));
00214 
00215   // now print out the filtered jets and reconstruct total 
00216   // at the same time
00217   cout << "Filtered pieces are " << endl;
00218   cout << " " << filt_subjets[0] << endl;
00219   PseudoJet filtered_total = filt_subjets[0];
00220   for (unsigned i = 1; i < nfilt && i < filt_subjets.size(); i++) {
00221     cout << " " << filt_subjets[i] << endl;
00222     flav_recombiner.plus_equal(filtered_total, filt_subjets[i]);
00223   }
00224   cout << "Filtered total is " << endl;
00225   cout << " " << filtered_total << endl;
00226 
00227 }
00228 
00229 
00230 //----------------------------------------------------------------------
00231 // does the actual work for printing out a jet
00232 //----------------------------------------------------------------------
00233 ostream & operator<<(ostream & ostr, PseudoJet & jet) {
00234   ostr << "pt, y, phi =" 
00235        << " " << setw(10) << jet.perp() 
00236        << " " << setw(6) <<  jet.rap()  
00237        << " " << setw(6) <<  jet.phi()  
00238        << ", mass = " << setw(10) << jet.m()
00239        << ", btag = " << jet.user_index();
00240   return ostr;
00241 }