rivet is hosted by Hepforge, IPPP Durham
ATLAS_2011_CONF_2011_090.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 
00011 namespace Rivet {
00012 
00013 
00014   class ATLAS_2011_CONF_2011_090 : public Analysis {
00015   public:
00016 
00017     /// @name Constructors etc.
00018     //@{
00019 
00020     /// Constructor
00021 
00022     ATLAS_2011_CONF_2011_090()
00023       : Analysis("ATLAS_2011_CONF_2011_090")
00024     {    }
00025 
00026     //@}
00027 
00028 
00029   public:
00030 
00031     /// @name Analysis methods
00032     //@{
00033 
00034     /// Book histograms and initialize projections before the run
00035     void init() {
00036 
00037       // projection to find the electrons
00038       IdentifiedFinalState elecs(EtaIn(-2.47, 2.47) 
00039                  & (Cuts::pT >= 20.0*GeV));
00040       elecs.acceptIdPair(PID::ELECTRON);
00041       addProjection(elecs, "elecs");
00042 
00043 
00044       // veto region electrons (from 2010 arXiv:1102.2357v2)
00045       Cut vetocut = EtaIn(-1.52, -1.37) | EtaIn( 1.37,  1.52);
00046       IdentifiedFinalState veto_elecs(vetocut & (Cuts::pT >= 10.0*GeV));
00047       veto_elecs.acceptIdPair(PID::ELECTRON);
00048       addProjection(veto_elecs, "veto_elecs");
00049 
00050 
00051       // projection to find the muons
00052       IdentifiedFinalState muons(EtaIn(-2.4, 2.4) 
00053                  & (Cuts::pT >= 10.0*GeV));
00054       muons.acceptIdPair(PID::MUON);
00055       addProjection(muons, "muons");
00056 
00057 
00058       // Jet finder
00059       VetoedFinalState vfs;
00060       vfs.addVetoPairId(PID::MUON);
00061       addProjection(FastJets(vfs, FastJets::ANTIKT, 0.4),
00062                    "AntiKtJets04");
00063 
00064 
00065       // all tracks (to do deltaR with leptons)
00066       addProjection(ChargedFinalState(-3.0,3.0,0.5*GeV),"cfs");
00067 
00068 
00069       // for pTmiss
00070       addProjection(VisibleFinalState(-4.9,4.9),"vfs");
00071 
00072 
00073       /// Book histograms
00074       _count_mu_channel = bookHisto1D("count_muon_channel", 1, 0., 1.);
00075       _count_e_channel = bookHisto1D("count_electron_channel", 1, 0., 1.);
00076       _hist_eTmiss_e = bookHisto1D("Et_miss_e", 50, 0., 500.);
00077       _hist_eTmiss_mu = bookHisto1D("Et_miss_mu", 50, 0., 500.);
00078       _hist_m_eff_e = bookHisto1D("m_eff_e", 60, 0., 1500.);
00079       _hist_m_eff_mu = bookHisto1D("m_eff_mu", 60, 0., 1500.);
00080       _hist_m_eff_e_final = bookHisto1D("m_eff_e_final", 15, 0., 1500.);
00081       _hist_m_eff_mu_final = bookHisto1D("m_eff_mu_final", 15, 0., 1500.);
00082 
00083 
00084 
00085     }
00086 
00087 
00088 
00089     /// Perform the per-event analysis
00090     void analyze(const Event& event) {
00091 
00092       const double weight = event.weight();
00093 
00094 
00095       Particles veto_e
00096     = applyProjection<IdentifiedFinalState>(event, "veto_elecs").particles();
00097       if ( ! veto_e.empty() ) {
00098         MSG_DEBUG("electrons in veto region");
00099         vetoEvent;
00100       }
00101 
00102       Jets cand_jets;
00103       foreach ( const Jet& jet,
00104           applyProjection<FastJets>(event, "AntiKtJets04").jetsByPt(20.0*GeV) ) {
00105         if ( fabs( jet.eta() ) < 2.8 ) {
00106           cand_jets.push_back(jet);
00107         }
00108       }
00109 
00110       Particles candtemp_e =
00111     applyProjection<IdentifiedFinalState>(event, "elecs").particlesByPt();
00112       Particles candtemp_mu =
00113     applyProjection<IdentifiedFinalState>(event,"muons").particlesByPt();
00114       Particles chg_tracks =
00115     applyProjection<ChargedFinalState>(event, "cfs").particles();
00116       Particles cand_mu;
00117       Particles cand_e;
00118 
00119 
00120       // pTcone around muon track
00121       foreach ( const Particle & mu, candtemp_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 )
00128       cand_mu.push_back(mu);
00129       }
00130 
00131 
00132       // pTcone around electron
00133       foreach ( const Particle e, candtemp_e ) {
00134     double pTinCone = -e.pT();
00135     foreach ( const Particle & track, chg_tracks ) {
00136       if ( deltaR(e.momentum(),track.momentum()) < 0.2 )
00137         pTinCone += track.pT();
00138     }
00139     if ( pTinCone < 0.10 * e.pT() )
00140       cand_e.push_back(e);
00141       }
00142 
00143 
00144 
00145       // discard jets that overlap with electrons
00146       Jets cand_jets_2;
00147       foreach ( const Jet& jet, cand_jets ) {
00148       bool away_from_e = true;
00149       foreach ( const Particle & e, cand_e ) {
00150         if ( deltaR(e.momentum(),jet.momentum()) <= 0.2 ) {
00151           away_from_e = false;
00152           break;
00153         }
00154       }
00155       if ( away_from_e )
00156         cand_jets_2.push_back( jet );
00157       }
00158 
00159       // only consider leptons far from jet
00160       Particles recon_e, recon_mu;
00161       foreach ( const Particle & e, cand_e ) {
00162         bool e_near_jet = false;
00163     foreach ( const Jet& jet, cand_jets_2 ) {
00164           if ( deltaR(e.momentum(),jet.momentum()) < 0.4 &&
00165            deltaR(e.momentum(),jet.momentum()) > 0.2 )
00166         e_near_jet = true;
00167     }
00168         if ( ! e_near_jet )
00169           recon_e.push_back( e );
00170       }
00171 
00172       foreach ( const Particle & mu, cand_mu ) {
00173     bool mu_near_jet = false;
00174     foreach ( const Jet& jet, cand_jets_2 ) {
00175       if ( deltaR(mu.momentum(),jet.momentum()) < 0.4 )
00176         mu_near_jet = true;
00177     }
00178     if ( ! mu_near_jet )
00179       recon_mu.push_back( mu );
00180       }
00181 
00182       // pTmiss
00183       Particles vfs_particles
00184     = applyProjection<VisibleFinalState>(event, "vfs").particles();
00185       FourMomentum pTmiss;
00186       foreach ( const Particle & p, vfs_particles ) {
00187     pTmiss -= p.momentum();
00188       }
00189       double eTmiss = pTmiss.pT();
00190 
00191 
00192       // final jet filter
00193       Jets recon_jets;
00194       foreach ( const Jet& jet, cand_jets_2 ) {
00195       recon_jets.push_back( jet );
00196       }
00197 
00198 
00199 
00200 
00201       // ==================== observables ====================
00202 
00203 
00204       // Njets
00205 
00206       int Njets = 0;
00207       double pTmiss_phi = pTmiss.phi();
00208       foreach ( const Jet& jet, recon_jets ) {
00209     if ( fabs(jet.eta()) < 2.8 )
00210       Njets+=1;
00211       }
00212       if ( Njets < 3 ) {
00213     MSG_DEBUG("Only " << Njets << " jets w/ eta<2.8 left");
00214     vetoEvent;
00215       }
00216 
00217       if ( recon_jets[0].pT() <= 60.0 * GeV ) {
00218     MSG_DEBUG("No hard leading jet in " << recon_jets.size() << " jets");
00219     vetoEvent;
00220       }
00221       for ( int i = 1; i < 3; ++i ) {
00222     if ( recon_jets[i].pT() <= 25*GeV ) {
00223       vetoEvent;
00224     }
00225       }
00226 
00227       for ( int i = 0; i < 3; ++i ) {
00228     double dPhi = deltaPhi( pTmiss_phi, recon_jets[i].momentum().phi() );
00229     if ( dPhi <= 0.2 ) {
00230       MSG_DEBUG("dPhi too small");
00231       vetoEvent;
00232       break;
00233     }
00234       }
00235 
00236 
00237       Particles lepton;
00238       if ( recon_mu.empty() && recon_e.empty() ) {
00239     MSG_DEBUG("No leptons");
00240     vetoEvent;
00241       }
00242       else {
00243     foreach ( const Particle & mu, recon_mu )
00244         lepton.push_back(mu);
00245         foreach ( const Particle & e, recon_e )
00246         lepton.push_back(e);
00247       }
00248 
00249       std::sort(lepton.begin(), lepton.end(), cmpMomByPt);
00250 
00251       double e_id = 11;
00252       double mu_id = 13;
00253 
00254       // one hard leading lepton cut
00255       if ( fabs(lepton[0].pdgId()) == e_id &&
00256            lepton[0].pT() <= 25*GeV ) {
00257     vetoEvent;
00258       }
00259       else if ( fabs(lepton[0].pdgId()) == mu_id &&
00260                 lepton[0].pT() <= 20*GeV ) {
00261     vetoEvent;
00262       }
00263 
00264       // exactly one hard leading lepton cut
00265       if(lepton.size()>1) {
00266     if ( fabs(lepton[1].pdgId()) == e_id &&
00267          lepton[1].pT() > 20*GeV ) {
00268       vetoEvent;
00269     }
00270     else if ( fabs(lepton[1].pdgId()) == mu_id &&
00271           lepton[1].pT() > 10*GeV ) {
00272       vetoEvent;
00273     }
00274       }
00275 
00276 
00277       // ==================== FILL ====================
00278 
00279 
00280       FourMomentum pT_l = lepton[0].momentum();
00281 
00282 
00283       double dPhi = deltaPhi( pT_l.phi(), pTmiss_phi);
00284       double mT = sqrt( 2 * pT_l.pT() * eTmiss * (1 - cos(dPhi)) );
00285 
00286 
00287       // effective mass
00288       double m_eff = eTmiss + pT_l.pT()
00289     + recon_jets[0].pT()
00290     + recon_jets[1].pT()
00291     + recon_jets[2].pT();
00292 
00293 
00294       // Electron channel signal region
00295 
00296       if ( fabs( lepton[0].pdgId() ) == e_id ) {
00297 
00298         _hist_eTmiss_e->fill(eTmiss, weight);
00299         _hist_m_eff_e->fill(m_eff, weight);
00300 
00301         if ( mT > 100*GeV && eTmiss > 125*GeV ) {
00302       _hist_m_eff_e_final->fill(m_eff, weight);
00303       if ( m_eff > 500*GeV && eTmiss > 0.25*m_eff ) {
00304             _count_e_channel->fill(0.5,weight);
00305       }
00306         }
00307       }
00308 
00309       // Muon channel signal region
00310 
00311       else if ( fabs( lepton[0].pdgId() ) == mu_id ) {
00312 
00313         _hist_eTmiss_mu->fill(eTmiss, weight);
00314         _hist_m_eff_mu->fill(m_eff, weight);
00315 
00316         if ( mT > 100*GeV && eTmiss > 125*GeV ) {
00317           _hist_m_eff_mu_final->fill(m_eff, weight);
00318       if ( m_eff > 500*GeV && eTmiss > 0.25*m_eff ) {
00319             _count_mu_channel->fill(0.5,weight);
00320           }
00321         }
00322 
00323       }
00324 
00325 
00326     }
00327 
00328     //@}
00329 
00330 
00331     void finalize() {
00332     scale( _hist_eTmiss_e      , 10.  * 165. * crossSection()/picobarn/sumOfWeights() );
00333     scale( _hist_eTmiss_mu     , 10.  * 165. * crossSection()/picobarn/sumOfWeights() );
00334     scale( _hist_m_eff_e       , 25.  * 165. * crossSection()/picobarn/sumOfWeights() );
00335     scale( _hist_m_eff_mu      , 25.  * 165. * crossSection()/picobarn/sumOfWeights() );
00336     scale( _hist_m_eff_e_final , 100. * 165. * crossSection()/picobarn/sumOfWeights() );
00337     scale( _hist_m_eff_mu_final, 100. * 165. * crossSection()/picobarn/sumOfWeights() );
00338     }
00339 
00340   private:
00341 
00342     /// @name Histograms
00343     //@{
00344     Histo1DPtr _count_e_channel;
00345     Histo1DPtr _count_mu_channel;
00346 
00347     Histo1DPtr _hist_eTmiss_e;
00348     Histo1DPtr _hist_eTmiss_mu;
00349 
00350     Histo1DPtr _hist_m_eff_e;
00351     Histo1DPtr _hist_m_eff_mu;
00352     Histo1DPtr _hist_m_eff_e_final;
00353     Histo1DPtr _hist_m_eff_mu_final;
00354 
00355 
00356     //@}
00357 
00358 
00359   };
00360 
00361 
00362 
00363   // The hook for the plugin system
00364   DECLARE_RIVET_PLUGIN(ATLAS_2011_CONF_2011_090);
00365 
00366 }