rivet is hosted by Hepforge, IPPP Durham
ATLAS_2011_I928289_W.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/WFinder.hh"
00005 
00006 namespace Rivet {
00007 
00008 
00009   class ATLAS_2011_I928289_W : public Analysis {
00010   public:
00011 
00012     /// Constructor
00013     ATLAS_2011_I928289_W()
00014       : Analysis("ATLAS_2011_I928289_W")
00015     {
00016       setNeedsCrossSection(true);
00017     }
00018 
00019 
00020     /// @name Analysis methods
00021     //@{
00022 
00023     /// Book histograms and initialise projections before the run
00024     void init() {
00025 
00026       ///Initialise and register projections here
00027       FinalState fs;
00028 
00029       Cut cut = (Cuts::pT >= 20*GeV);
00030 
00031       WFinder wfinder_el_bare(   fs, cut, PID::ELECTRON, 40.0*GeV, 7000.0*GeV, 25.0*GeV, 0.0, WFinder::CLUSTERNODECAY, WFinder::NOTRACK, WFinder::TRANSMASS);
00032       WFinder wfinder_el_dressed(fs, cut, PID::ELECTRON, 40.0*GeV, 7000.0*GeV, 25.0*GeV, 0.1, WFinder::CLUSTERNODECAY, WFinder::NOTRACK, WFinder::TRANSMASS);
00033       WFinder wfinder_mu_bare   (fs, cut, PID::MUON    , 40.0*GeV, 7000.0*GeV, 25.0*GeV, 0.0, WFinder::CLUSTERNODECAY, WFinder::NOTRACK, WFinder::TRANSMASS);
00034       WFinder wfinder_mu_dressed(fs, cut, PID::MUON    , 40.0*GeV, 7000.0*GeV, 25.0*GeV, 0.1, WFinder::CLUSTERNODECAY, WFinder::NOTRACK, WFinder::TRANSMASS);
00035 
00036       declare(wfinder_el_bare   , "WFinder_el_bare");
00037       declare(wfinder_el_dressed, "WFinder_el_dressed");
00038       declare(wfinder_mu_bare   , "WFinder_mu_bare");
00039       declare(wfinder_mu_dressed, "WFinder_mu_dressed");
00040 
00041       /// Book histograms here
00042       _h_Wminus_lepton_eta_el_bare       = bookHisto1D(3, 1, 1);
00043       _h_Wminus_lepton_eta_el_dressed    = bookHisto1D(3, 1, 2);
00044       _h_Wminus_lepton_eta_mu_bare       = bookHisto1D(3, 1, 3);
00045       _h_Wminus_lepton_eta_mu_dressed    = bookHisto1D(3, 1, 4);
00046       _h_Wplus_lepton_eta_el_bare        = bookHisto1D(5, 1, 1);
00047       _h_Wplus_lepton_eta_el_dressed     = bookHisto1D(5, 1, 2);
00048       _h_Wplus_lepton_eta_mu_bare        = bookHisto1D(5, 1, 3);
00049       _h_Wplus_lepton_eta_mu_dressed     = bookHisto1D(5, 1, 4);
00050       _h_W_asym_eta_el_bare              = bookScatter2D(7, 1, 1);
00051       _h_W_asym_eta_el_dressed           = bookScatter2D(7, 1, 2);
00052       _h_W_asym_eta_mu_bare              = bookScatter2D(7, 1, 3);
00053       _h_W_asym_eta_mu_dressed           = bookScatter2D(7, 1, 4);
00054 
00055     }
00056 
00057 
00058     /// Perform the per-event analysis
00059     void analyze(const Event& event) {
00060 
00061       const WFinder& wfinder_el_bare     = apply<WFinder>(event, "WFinder_el_bare");
00062       const WFinder& wfinder_el_dressed  = apply<WFinder>(event, "WFinder_el_dressed");
00063       const WFinder& wfinder_mu_bare     = apply<WFinder>(event, "WFinder_mu_bare");
00064       const WFinder& wfinder_mu_dressed  = apply<WFinder>(event, "WFinder_mu_dressed");
00065 
00066       const double weight = event.weight();
00067       fillPlots1D(wfinder_el_bare   , _h_Wplus_lepton_eta_el_bare   , _h_Wminus_lepton_eta_el_bare   , weight);
00068       fillPlots1D(wfinder_el_dressed, _h_Wplus_lepton_eta_el_dressed, _h_Wminus_lepton_eta_el_dressed, weight);
00069       fillPlots1D(wfinder_mu_bare   , _h_Wplus_lepton_eta_mu_bare   , _h_Wminus_lepton_eta_mu_bare   , weight);
00070       fillPlots1D(wfinder_mu_dressed, _h_Wplus_lepton_eta_mu_dressed, _h_Wminus_lepton_eta_mu_dressed, weight);
00071     }
00072 
00073 
00074     void fillPlots1D(const WFinder& wfinder, Histo1DPtr hist_plus, Histo1DPtr hist_minus, double weight) {
00075       if (wfinder.bosons().size() != 1) return;
00076       const Particle l = wfinder.constituentLeptons()[0];
00077       const FourMomentum& miss = wfinder.constituentNeutrinos()[0].momentum();
00078       if (l.pT() > 20*GeV && miss.Et() > 25*GeV && wfinder.mT() > 40*GeV)
00079         (l.charge3() > 0 ? hist_plus : hist_minus)->fill(l.abseta(), weight);
00080     }
00081 
00082 
00083     /// Normalise histograms etc., after the run
00084     void finalize() {
00085 
00086       // Construct asymmetry: (dsig+/deta - dsig-/deta) / (dsig+/deta + dsig-/deta)
00087       divide(*_h_Wplus_lepton_eta_el_bare - *_h_Wminus_lepton_eta_el_bare,
00088              *_h_Wplus_lepton_eta_el_bare + *_h_Wminus_lepton_eta_el_bare,
00089              _h_W_asym_eta_el_bare);
00090       divide(*_h_Wplus_lepton_eta_el_dressed - *_h_Wminus_lepton_eta_el_dressed,
00091              *_h_Wplus_lepton_eta_el_dressed + *_h_Wminus_lepton_eta_el_dressed,
00092              _h_W_asym_eta_el_dressed);
00093       divide(*_h_Wplus_lepton_eta_mu_bare - *_h_Wminus_lepton_eta_mu_bare,
00094              *_h_Wplus_lepton_eta_mu_bare + *_h_Wminus_lepton_eta_mu_bare,
00095              _h_W_asym_eta_mu_bare);
00096       divide(*_h_Wplus_lepton_eta_mu_dressed - *_h_Wminus_lepton_eta_mu_dressed,
00097              *_h_Wplus_lepton_eta_mu_dressed + *_h_Wminus_lepton_eta_mu_dressed,
00098              _h_W_asym_eta_mu_dressed);
00099 
00100       // Print summary info
00101       const double xs_pb(crossSection() / picobarn);
00102       const double sumw(sumOfWeights());
00103       MSG_DEBUG( "Cross-section/pb     : " << xs_pb       );
00104       MSG_DEBUG( "Sum of weights       : " << sumw        );
00105       MSG_DEBUG( "nEvents              : " << numEvents() );
00106 
00107       ///  Normalise, scale and otherwise manipulate histograms here
00108       const double sf = 0.5 * xs_pb / sumw; // 0.5 accounts for rapidity bin width
00109       scale(_h_Wminus_lepton_eta_el_bare   , sf);
00110       scale(_h_Wminus_lepton_eta_el_dressed, sf);
00111       scale(_h_Wminus_lepton_eta_mu_bare   , sf);
00112       scale(_h_Wminus_lepton_eta_mu_dressed, sf);
00113       scale(_h_Wplus_lepton_eta_el_bare    , sf);
00114       scale(_h_Wplus_lepton_eta_el_dressed , sf);
00115       scale(_h_Wplus_lepton_eta_mu_bare    , sf);
00116       scale(_h_Wplus_lepton_eta_mu_dressed , sf);
00117 
00118     }
00119 
00120     //@}
00121 
00122 
00123   private:
00124 
00125     /// @name Histograms
00126     //@{
00127     Histo1DPtr _h_Wminus_lepton_eta_el_bare;
00128     Histo1DPtr _h_Wminus_lepton_eta_el_dressed;
00129     Histo1DPtr _h_Wminus_lepton_eta_mu_bare;
00130     Histo1DPtr _h_Wminus_lepton_eta_mu_dressed;
00131     Histo1DPtr _h_Wplus_lepton_eta_el_bare;
00132     Histo1DPtr _h_Wplus_lepton_eta_el_dressed;
00133     Histo1DPtr _h_Wplus_lepton_eta_mu_bare;
00134     Histo1DPtr _h_Wplus_lepton_eta_mu_dressed;
00135     Scatter2DPtr _h_W_asym_eta_el_bare;
00136     Scatter2DPtr _h_W_asym_eta_el_dressed;
00137     Scatter2DPtr _h_W_asym_eta_mu_bare;
00138     Scatter2DPtr _h_W_asym_eta_mu_dressed;
00139 
00140     //@}
00141 
00142   };
00143 
00144 
00145   // The hook for the plugin system
00146   DECLARE_RIVET_PLUGIN(ATLAS_2011_I928289_W);
00147 
00148 }