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