rivet is hosted by Hepforge, IPPP Durham
CMS_2012_I1107658.cc
Go to the documentation of this file.
00001 // -*- C++ -*-
00002 #include "Rivet/Analysis.hh"
00003 #include "Rivet/Projections/FinalState.hh"
00004 #include "Rivet/Projections/ZFinder.hh"
00005 #include "Rivet/Projections/FinalState.hh"
00006 #include "Rivet/Projections/ChargedFinalState.hh"
00007 #include "Rivet/Projections/VetoedFinalState.hh"
00008 #include "Rivet/ParticleName.hh"
00009 
00010 namespace Rivet {
00011 
00012 
00013   /// Underlying event activity in the Drell-Yan process at 7 TeV
00014   class CMS_2012_I1107658 : public Analysis {
00015   public:
00016 
00017     /// Constructor
00018     CMS_2012_I1107658()
00019       : Analysis("CMS_2012_I1107658")
00020     {   }
00021 
00022 
00023     /// Initialization
00024     void init() {
00025 
00026       /// @note Using a bare muon Z (but with a clustering radius!?)
00027       ZFinder zfinder(FinalState(), -2.4, 2.4, 20*GeV, PID::MUON, 4*GeV, 140*GeV, 0.2, ZFinder::NOCLUSTER);
00028       addProjection(zfinder, "ZFinder");
00029 
00030       ChargedFinalState cfs(-2, 2, 500*MeV);
00031       VetoedFinalState nonmuons(cfs);
00032       nonmuons.addVetoPairId(PID::MUON);
00033       addProjection(nonmuons, "nonmuons");
00034 
00035       _h_Nchg_towards_pTmumu                 = bookProfile1D(1, 1, 1);
00036       _h_Nchg_transverse_pTmumu              = bookProfile1D(2, 1, 1);
00037       _h_Nchg_away_pTmumu                    = bookProfile1D(3, 1, 1);
00038       _h_pTsum_towards_pTmumu                = bookProfile1D(4, 1, 1);
00039       _h_pTsum_transverse_pTmumu             = bookProfile1D(5, 1, 1);
00040       _h_pTsum_away_pTmumu                   = bookProfile1D(6, 1, 1);
00041       _h_avgpT_towards_pTmumu                = bookProfile1D(7, 1, 1);
00042       _h_avgpT_transverse_pTmumu             = bookProfile1D(8, 1, 1);
00043       _h_avgpT_away_pTmumu                   = bookProfile1D(9, 1, 1);
00044       _h_Nchg_towards_plus_transverse_Mmumu  = bookProfile1D(10, 1, 1);
00045       _h_pTsum_towards_plus_transverse_Mmumu = bookProfile1D(11, 1, 1);
00046       _h_avgpT_towards_plus_transverse_Mmumu = bookProfile1D(12, 1, 1);
00047       _h_Nchg_towards_zmass_81_101           = bookHisto1D(13, 1, 1);
00048       _h_Nchg_transverse_zmass_81_101        = bookHisto1D(14, 1, 1);
00049       _h_Nchg_away_zmass_81_101              = bookHisto1D(15, 1, 1);
00050       _h_pT_towards_zmass_81_101             = bookHisto1D(16, 1, 1);
00051       _h_pT_transverse_zmass_81_101          = bookHisto1D(17, 1, 1);
00052       _h_pT_away_zmass_81_101                = bookHisto1D(18, 1, 1);
00053       _h_Nchg_transverse_zpt_5               = bookHisto1D(19, 1, 1);
00054       _h_pT_transverse_zpt_5                 = bookHisto1D(20, 1, 1);
00055     }
00056 
00057 
00058     /// Perform the per-event analysis
00059     void analyze(const Event& event) {
00060       const double weight = event.weight();
00061       const ZFinder& zfinder = applyProjection<ZFinder>(event, "ZFinder");
00062 
00063       if (zfinder.bosons().size() != 1) vetoEvent;
00064 
00065       double Zpt = zfinder.bosons()[0].pT()/GeV;
00066       double Zphi = zfinder.bosons()[0].momentum().phi();
00067       double Zmass = zfinder.bosons()[0].momentum().mass()/GeV;
00068 
00069       Particles particles = applyProjection<VetoedFinalState>(event, "nonmuons").particles();
00070 
00071       int nTowards = 0;
00072       int nTransverse = 0;
00073       int nAway = 0;
00074       double ptSumTowards = 0;
00075       double ptSumTransverse = 0;
00076       double ptSumAway = 0;
00077 
00078       foreach (const Particle& p, particles) {
00079         double dphi = fabs(deltaPhi(Zphi, p.momentum().phi()));
00080         double pT = p.pT();
00081 
00082         if ( dphi < M_PI/3 ) {
00083           nTowards++;
00084           ptSumTowards += pT;
00085           if (Zmass > 81. && Zmass < 101.) _h_pT_towards_zmass_81_101->fill(pT, weight);
00086         } else if ( dphi < 2.*M_PI/3 ) {
00087           nTransverse++;
00088           ptSumTransverse += pT;
00089           if (Zmass > 81. && Zmass < 101.) _h_pT_transverse_zmass_81_101->fill(pT, weight);
00090           if (Zpt < 5.) _h_pT_transverse_zpt_5->fill(pT, weight);
00091         } else {
00092           nAway++;
00093           ptSumAway += pT;
00094           if (Zmass > 81. && Zmass < 101.) _h_pT_away_zmass_81_101->fill(pT, weight);
00095         }
00096 
00097       } // Loop over particles
00098 
00099 
00100       const double area = 8./3.*M_PI;
00101       if (Zmass > 81. && Zmass < 101.) {
00102         _h_Nchg_towards_pTmumu->         fill(Zpt, 1./area * nTowards, weight);
00103         _h_Nchg_transverse_pTmumu->      fill(Zpt, 1./area * nTransverse, weight);
00104         _h_Nchg_away_pTmumu->            fill(Zpt, 1./area * nAway, weight);
00105         _h_pTsum_towards_pTmumu->        fill(Zpt, 1./area * ptSumTowards, weight);
00106         _h_pTsum_transverse_pTmumu->     fill(Zpt, 1./area * ptSumTransverse, weight);
00107         _h_pTsum_away_pTmumu->           fill(Zpt, 1./area * ptSumAway, weight);
00108         if (nTowards > 0)    _h_avgpT_towards_pTmumu->    fill(Zpt, ptSumTowards/nTowards, weight);
00109         if (nTransverse > 0) _h_avgpT_transverse_pTmumu-> fill(Zpt, ptSumTransverse/nTransverse, weight);
00110         if (nAway > 0)       _h_avgpT_away_pTmumu->       fill(Zpt, ptSumAway/nAway, weight);
00111         _h_Nchg_towards_zmass_81_101->   fill(nTowards, weight);
00112         _h_Nchg_transverse_zmass_81_101->fill(nTransverse, weight);
00113         _h_Nchg_away_zmass_81_101->      fill(nAway, weight);
00114       }
00115 
00116       if (Zpt < 5.) {
00117         _h_Nchg_towards_plus_transverse_Mmumu->fill(Zmass, (nTowards + nTransverse)/(2.*area), weight);
00118         _h_pTsum_towards_plus_transverse_Mmumu->fill(Zmass, (ptSumTowards + ptSumTransverse)/(2.*area), weight);
00119         if ((nTowards + nTransverse) > 0) _h_avgpT_towards_plus_transverse_Mmumu->fill(Zmass, (ptSumTowards + ptSumTransverse)/(nTowards + nTransverse), weight);
00120         _h_Nchg_transverse_zpt_5->fill(nTransverse, weight);
00121       }
00122 
00123     }
00124 
00125 
00126     /// Normalise histograms etc., after the run
00127     void finalize() {
00128       if (integral(_h_Nchg_towards_zmass_81_101)    > 0) scale(_h_pT_towards_zmass_81_101,    1/integral(_h_Nchg_towards_zmass_81_101));
00129       if (integral(_h_Nchg_transverse_zmass_81_101) > 0) scale(_h_pT_transverse_zmass_81_101, 1/integral(_h_Nchg_transverse_zmass_81_101));
00130       if (integral(_h_Nchg_away_zmass_81_101)       > 0) scale(_h_pT_away_zmass_81_101,       1/integral(_h_Nchg_away_zmass_81_101));
00131       if (integral(_h_Nchg_transverse_zpt_5)        > 0) scale(_h_pT_transverse_zpt_5,        1/integral(_h_Nchg_transverse_zpt_5));
00132       normalize(_h_Nchg_towards_zmass_81_101);
00133       normalize(_h_Nchg_transverse_zmass_81_101);
00134       normalize(_h_Nchg_away_zmass_81_101);
00135       normalize(_h_Nchg_transverse_zpt_5);
00136     }
00137 
00138 
00139   private:
00140 
00141 
00142     /// @name Histogram objects
00143     //@{
00144     Profile1DPtr _h_Nchg_towards_pTmumu;
00145     Profile1DPtr _h_Nchg_transverse_pTmumu;
00146     Profile1DPtr _h_Nchg_away_pTmumu;
00147     Profile1DPtr _h_pTsum_towards_pTmumu;
00148     Profile1DPtr _h_pTsum_transverse_pTmumu;
00149     Profile1DPtr _h_pTsum_away_pTmumu;
00150     Profile1DPtr _h_avgpT_towards_pTmumu;
00151     Profile1DPtr _h_avgpT_transverse_pTmumu;
00152     Profile1DPtr _h_avgpT_away_pTmumu;
00153     Profile1DPtr _h_Nchg_towards_plus_transverse_Mmumu;
00154     Profile1DPtr _h_pTsum_towards_plus_transverse_Mmumu;
00155     Profile1DPtr _h_avgpT_towards_plus_transverse_Mmumu;
00156     Histo1DPtr _h_Nchg_towards_zmass_81_101;
00157     Histo1DPtr _h_Nchg_transverse_zmass_81_101;
00158     Histo1DPtr _h_Nchg_away_zmass_81_101;
00159     Histo1DPtr _h_pT_towards_zmass_81_101;
00160     Histo1DPtr _h_pT_transverse_zmass_81_101;
00161     Histo1DPtr _h_pT_away_zmass_81_101;
00162     Histo1DPtr _h_Nchg_transverse_zpt_5;
00163     Histo1DPtr _h_pT_transverse_zpt_5;
00164     //@}
00165 
00166   };
00167 
00168 
00169   // Hook for the plugin system
00170   DECLARE_RIVET_PLUGIN(CMS_2012_I1107658);
00171 
00172 }