rivet is hosted by Hepforge, IPPP Durham
MC_WINC.cc
Go to the documentation of this file.
00001 // -*- C++ -*-
00002 #include "Rivet/Analysis.hh"
00003 #include "Rivet/Projections/WFinder.hh"
00004 
00005 namespace Rivet {
00006 
00007   using namespace Cuts;
00008 
00009 
00010   /// @brief MC validation analysis for inclusive W events
00011   class MC_WINC : public Analysis {
00012   public:
00013 
00014     /// Default constructor
00015     MC_WINC()
00016       : Analysis("MC_WINC")
00017     {    }
00018 
00019 
00020     /// @name Analysis methods
00021     //@{
00022 
00023     /// Book histograms
00024     void init() {
00025       FinalState fs;
00026       WFinder wfinder(fs, etaIn(-3.5, 3.5) & (pT >= 25.0*GeV), PID::ELECTRON, 60.0*GeV, 100.0*GeV, 25.0*GeV, 0.2);
00027       addProjection(wfinder, "WFinder");
00028 
00029       _h_W_mass = bookHisto1D("W_mass", 50, 55.0, 105.0);
00030       _h_W_pT = bookHisto1D("W_pT", logspace(100, 1.0, 0.5*sqrtS()));
00031       _h_W_pT_peak = bookHisto1D("W_pT_peak", 25, 0.0, 125.0);
00032       _h_W_y = bookHisto1D("W_y", 40, -4.0, 4.0);
00033       _h_W_phi = bookHisto1D("W_phi", 25, 0.0, TWOPI);
00034       _h_Wplus_pT = bookHisto1D("Wplus_pT", logspace(25, 10.0, 0.5*sqrtS()));
00035       _h_Wminus_pT = bookHisto1D("Wminus_pT", logspace(25, 10.0, 0.5*sqrtS()));
00036       _h_lepton_pT = bookHisto1D("lepton_pT", logspace(100, 10.0, 0.25*sqrtS()));
00037       _h_lepton_eta = bookHisto1D("lepton_eta", 40, -4.0, 4.0);
00038       _htmp_dsigminus_deta = bookHisto1D("lepton_dsigminus_deta", 20, 0.0, 4.0);
00039       _htmp_dsigplus_deta  = bookHisto1D("lepton_dsigplus_deta", 20, 0.0, 4.0);
00040 
00041       _h_asym = bookScatter2D("W_chargeasymm_eta");
00042       _h_asym_pT = bookScatter2D("W_chargeasymm_pT");
00043     }
00044 
00045 
00046 
00047     /// Do the analysis
00048     void analyze(const Event & e) {
00049       const double weight = e.weight();
00050 
00051       const WFinder& wfinder = applyProjection<WFinder>(e, "WFinder");
00052       if (wfinder.bosons().size() != 1) {
00053         vetoEvent;
00054       }
00055 
00056       double charge3_x_eta = 0;
00057       int charge3 = 0;
00058       FourMomentum emom;
00059       FourMomentum wmom(wfinder.bosons().front().momentum());
00060       _h_W_mass->fill(wmom.mass(), weight);
00061       _h_W_pT->fill(wmom.pT(), weight);
00062       _h_W_pT_peak->fill(wmom.pT(), weight);
00063       _h_W_y->fill(wmom.rapidity(), weight);
00064       _h_W_phi->fill(wmom.phi(), weight);
00065       Particle l=wfinder.constituentLeptons()[0];
00066       _h_lepton_pT->fill(l.pT(), weight);
00067       _h_lepton_eta->fill(l.eta(), weight);
00068       if (PID::threeCharge(l.pid()) != 0) {
00069         emom = l.momentum();
00070         charge3_x_eta = PID::threeCharge(l.pid()) * emom.eta();
00071         charge3 = PID::threeCharge(l.pid());
00072       }
00073       assert(charge3_x_eta != 0);
00074       assert(charge3!=0);
00075       if (emom.Et() > 30/GeV) {
00076         if (charge3_x_eta < 0) {
00077           _htmp_dsigminus_deta->fill(emom.eta(), weight);
00078         } else {
00079           _htmp_dsigplus_deta->fill(emom.eta(), weight);
00080         }
00081       }
00082       if (charge3 < 0) {
00083         _h_Wminus_pT->fill(wmom.pT(), weight);
00084       } else {
00085         _h_Wplus_pT->fill(wmom.pT(), weight);
00086       }
00087     }
00088 
00089 
00090     /// Finalize
00091     void finalize() {
00092       scale(_h_W_mass, crossSection()/sumOfWeights());
00093       scale(_h_W_pT, crossSection()/sumOfWeights());
00094       scale(_h_W_pT_peak, crossSection()/sumOfWeights());
00095       scale(_h_W_y, crossSection()/sumOfWeights());
00096       scale(_h_W_phi, crossSection()/sumOfWeights());
00097       scale(_h_lepton_pT, crossSection()/sumOfWeights());
00098       scale(_h_lepton_eta, crossSection()/sumOfWeights());
00099 
00100       // Construct asymmetry: (dsig+/deta - dsig-/deta) / (dsig+/deta + dsig-/deta) for each Et region
00101       divide(*_htmp_dsigplus_deta - *_htmp_dsigminus_deta,
00102              *_htmp_dsigplus_deta + *_htmp_dsigminus_deta,
00103              _h_asym);
00104 
00105       // // W charge asymmetry vs. pTW: dsig+/dpT / dsig-/dpT
00106       divide(_h_Wplus_pT, _h_Wminus_pT,
00107              _h_asym_pT);
00108 
00109       scale(_h_Wplus_pT, crossSection()/sumOfWeights());
00110       scale(_h_Wminus_pT, crossSection()/sumOfWeights());
00111     }
00112 
00113     //@}
00114 
00115 
00116   private:
00117 
00118     /// @name Histograms
00119     //@{
00120     Histo1DPtr _h_W_mass;
00121     Histo1DPtr _h_W_pT;
00122     Histo1DPtr _h_W_pT_peak;
00123     Histo1DPtr _h_W_y;
00124     Histo1DPtr _h_W_phi;
00125     Histo1DPtr _h_Wplus_pT;
00126     Histo1DPtr _h_Wminus_pT;
00127     Histo1DPtr _h_lepton_pT;
00128     Histo1DPtr _h_lepton_eta;
00129 
00130     Histo1DPtr _htmp_dsigminus_deta;
00131     Histo1DPtr _htmp_dsigplus_deta;
00132 
00133     Scatter2DPtr _h_asym;
00134     Scatter2DPtr _h_asym_pT;
00135     //@}
00136 
00137   };
00138 
00139 
00140 
00141   // The hook for the plugin system
00142   DECLARE_RIVET_PLUGIN(MC_WINC);
00143 
00144 }