rivet is hosted by Hepforge, IPPP Durham
ATLAS_2012_I1186556.cc
Go to the documentation of this file.
00001 // -*- C++ -*-
00002 #include "Rivet/Analysis.hh"
00003 #include "Rivet/Tools/BinnedHistogram.hh"
00004 #include "Rivet/Projections/FinalState.hh"
00005 #include "Rivet/Projections/ChargedFinalState.hh"
00006 #include "Rivet/Projections/VisibleFinalState.hh"
00007 #include "Rivet/Projections/IdentifiedFinalState.hh"
00008 #include "Rivet/Projections/VetoedFinalState.hh"
00009 #include "Rivet/Projections/FastJets.hh"
00010 #include "Rivet/Tools/RivetMT2.hh"
00011 
00012 namespace Rivet {
00013 
00014   using namespace Cuts;
00015 
00016 
00017   class ATLAS_2012_I1186556 : public Analysis {
00018   public:
00019 
00020     /// @name Constructors etc.
00021     //@{
00022 
00023     /// Constructor
00024 
00025     ATLAS_2012_I1186556()
00026       : Analysis("ATLAS_2012_I1186556")
00027     {    }
00028 
00029     //@}
00030 
00031 
00032   public:
00033 
00034     /// @name Analysis methods
00035     //@{
00036 
00037     /// Book histograms and initialize projections before the run
00038     void init() {
00039 
00040       // projection to find the electrons
00041       IdentifiedFinalState elecs(etaIn(-2.47, 2.47) 
00042                  & (pT >= 20.0*GeV));
00043       elecs.acceptIdPair(PID::ELECTRON);
00044       addProjection(elecs, "elecs");
00045 
00046       // projection to find the muons
00047       IdentifiedFinalState muons(etaIn(-2.4, 2.4) 
00048                  & (pT >= 10.0*GeV));
00049       muons.acceptIdPair(PID::MUON);
00050       addProjection(muons, "muons");
00051 
00052       // Jet finder
00053       VetoedFinalState vfs;
00054       vfs.addVetoPairId(PID::MUON);
00055       addProjection(FastJets(vfs, FastJets::ANTIKT, 0.4), "AntiKtJets04");
00056 
00057       // all tracks (to do deltaR with leptons)
00058       addProjection(ChargedFinalState(-3.0,3.0,1.*GeV),"cfs");
00059 
00060       // for pTmiss
00061       addProjection(VisibleFinalState(-4.9,4.9),"vfs");
00062 
00063       // Book histograms
00064       _count_SR_SF     = bookHisto1D("count_SR_SF"    , 1, 0., 1.);
00065       _count_SR_OF     = bookHisto1D("count_SR_OF"    , 1, 0., 1.);
00066 
00067       _hist_mT2_SF_exp = bookHisto1D("hist_mT2_SF_exp", 40 , 0., 200. );
00068       _hist_mT2_OF_exp = bookHisto1D("hist_mT2_OF_exp", 40 , 0., 200. );
00069       _hist_mT2_SF_MC  = bookHisto1D("hist_mT2_SF_MC" , 500, 0., 1000.);
00070       _hist_mT2_OF_MC  = bookHisto1D("hist_mT2_OF_MC" , 500, 0., 1000.);
00071 
00072     }
00073 
00074     /// Perform the per-event analysis
00075     void analyze(const Event& event) {
00076       const double weight = event.weight();
00077 
00078       // get the candiate jets
00079       Jets cand_jets;
00080       foreach ( const Jet& jet,
00081                 applyProjection<FastJets>(event, "AntiKtJets04").jetsByPt(20.0*GeV) ) {
00082         if ( fabs( jet.eta() ) < 4.5 ) {
00083           cand_jets.push_back(jet);
00084         }
00085       }
00086       // charged tracks for isolation
00087       Particles chg_tracks =
00088         applyProjection<ChargedFinalState>(event, "cfs").particles();
00089       // find the electrons
00090       Particles cand_e;
00091       foreach( const Particle & e,
00092                applyProjection<IdentifiedFinalState>(event, "elecs").particlesByPt()) {
00093         // remove any leptons within 0.4 of any candidate jets
00094         bool e_near_jet = false;
00095         foreach ( const Jet& jet, cand_jets ) {
00096           double dR = deltaR(e.momentum(),jet.momentum());
00097           if ( dR < 0.4 && dR > 0.2 ) {
00098             e_near_jet = true;
00099             break;
00100           }
00101         }
00102     if ( e_near_jet ) continue;
00103     cand_e.push_back(e);
00104       }
00105       Particles cand_mu;
00106       foreach( const Particle & mu,
00107                applyProjection<IdentifiedFinalState>(event, "muons").particlesByPt()) {
00108         // remove any leptons within 0.4 of any candidate jets
00109         bool mu_near_jet = false;
00110         foreach ( const Jet& jet, cand_jets ) {
00111           if ( deltaR(mu.momentum(),jet.momentum()) < 0.4 ) {
00112             mu_near_jet = true;
00113             break;
00114           }
00115         }
00116         if ( mu_near_jet ) continue;
00117     cand_mu.push_back(mu);
00118       }
00119       // pTcone around muon track
00120       Particles recon_mu;
00121       foreach ( const Particle & mu, cand_mu ) {
00122         double pTinCone = -mu.pT();
00123         foreach ( const Particle & track, chg_tracks ) {
00124           if ( deltaR(mu.momentum(),track.momentum()) < 0.2 )
00125             pTinCone += track.pT();
00126         }
00127         if ( pTinCone < 1.8*GeV ) recon_mu.push_back(mu);
00128       }
00129       // pTcone around electron track
00130       Particles recon_e;
00131       foreach ( const Particle & e, cand_e ) {
00132         double pTinCone = -e.pT();
00133         foreach ( const Particle & track, chg_tracks ) {
00134           if ( deltaR(e.momentum(),track.momentum()) < 0.2 )
00135             pTinCone += track.pT();
00136         }
00137         if ( pTinCone < 0.1 * e.pT() ) recon_e.push_back(e);
00138       }
00139 
00140       // pTmiss
00141       FourMomentum pTmiss;
00142       foreach ( const Particle & p,
00143                 applyProjection<VisibleFinalState>(event, "vfs").particles() ) {
00144         pTmiss -= p.momentum();
00145       }
00146 
00147       // discard jets that overlap with electrons
00148       Jets recon_jets;
00149       foreach ( const Jet& jet, cand_jets ) {
00150         if(jet.abseta()>2.5||
00151            jet.perp()<20.) continue;
00152     bool away_from_e = true;
00153     foreach ( const Particle & e, cand_e ) {
00154       if ( deltaR(e.momentum(),jet.momentum()) < 0.2 ) {
00155         away_from_e = false;
00156         break;
00157       }
00158     }
00159     if ( away_from_e ) recon_jets.push_back( jet );
00160       }
00161 
00162       // put leptons into 1 vector and order by pT
00163       Particles leptons(recon_e.begin(),recon_e.end());
00164       leptons.insert(leptons.begin(),recon_mu.begin(),recon_mu.end());
00165       sort(leptons.begin(),leptons.end(),cmpMomByPt);
00166 
00167       // exactly two leptons
00168       if(leptons.size() !=2) vetoEvent;
00169 
00170       // hardest lepton pT greater the 25 (20) e(mu)
00171       if( (leptons[0].abspid()==PID::ELECTRON && leptons[0].perp()<25.) ||
00172       (leptons[0].abspid()==PID::ELECTRON && leptons[0].perp()<20.))
00173     vetoEvent;
00174 
00175       // require opposite sign
00176       if(leptons[0].pid()*leptons[1].pid()>0) vetoEvent;
00177 
00178       // and invariant mass > 20
00179       double mll = (leptons[0].momentum()+leptons[1].momentum()).mass();
00180       if(mll<20.) vetoEvent;
00181 
00182       // two jets 1st pT > 50 and second pT> 25
00183       if(recon_jets.size()<2 || recon_jets[0].perp()<50. ||
00184      recon_jets[1].perp()<25.) vetoEvent;
00185 
00186       // calculate mT2
00187       double m_T2 = mT2::mT2( leptons[0].momentum(),leptons[1].momentum(),
00188                   pTmiss,0.0 ); // zero mass invisibles
00189 
00190       // same flavour region
00191       if(leptons[0].pid()==-leptons[1].pid()) {
00192     // remove Z region
00193     if(mll>71.&&mll<111.) vetoEvent;
00194     // require at least 1 b jet
00195     unsigned int n_b=0;
00196     for(unsigned int ix=0;ix<recon_jets.size();++ix) {
00197        if(recon_jets[ix].bTagged() && rand()/static_cast<double>(RAND_MAX)<=0.60)
00198          ++n_b;
00199     }
00200     if(n_b==0) vetoEvent;
00201     _hist_mT2_SF_exp->fill(m_T2,weight);
00202     _hist_mT2_SF_MC ->fill(m_T2,weight);
00203     if(m_T2>120.) _count_SR_SF->fill(0.5,weight);
00204       }
00205       // opposite flavour region
00206       else {
00207     _hist_mT2_OF_exp->fill(m_T2,weight);
00208     _hist_mT2_OF_MC ->fill(m_T2,weight);
00209     if(m_T2>120.) _count_SR_OF->fill(0.5,weight);
00210       }
00211     }
00212     //@}
00213 
00214 
00215     void finalize() {
00216 
00217       double norm = 4.7* crossSection()/sumOfWeights()/femtobarn;
00218       scale(_count_SR_SF    ,   norm);
00219       scale(_count_SR_OF    ,   norm);
00220       scale(_hist_mT2_SF_exp,5.*norm);
00221       scale(_hist_mT2_OF_exp,5.*norm);
00222       scale(_hist_mT2_SF_MC ,   norm/4.7);
00223       scale(_hist_mT2_OF_MC ,   norm/4.7);
00224 
00225     }
00226 
00227   private:
00228 
00229     /// @name Histograms
00230     //@{
00231     Histo1DPtr _count_SR_SF;
00232     Histo1DPtr _count_SR_OF;
00233 
00234     Histo1DPtr _hist_mT2_SF_exp;
00235     Histo1DPtr _hist_mT2_OF_exp;
00236     Histo1DPtr _hist_mT2_SF_MC;
00237     Histo1DPtr _hist_mT2_OF_MC;
00238     //@}
00239 
00240   };
00241 
00242   // The hook for the plugin system
00243   DECLARE_RIVET_PLUGIN(ATLAS_2012_I1186556);
00244 
00245 }