ATLAS_2011_I928289_W.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #include <cmath> 00003 00004 #include "Rivet/Analysis.hh" 00005 #include "Rivet/Projections/FinalState.hh" 00006 #include "Rivet/Projections/WFinder.hh" 00007 00008 00009 namespace Rivet { 00010 00011 using namespace Cuts; 00012 00013 00014 class ATLAS_2011_I928289_W : public Analysis { 00015 public: 00016 00017 /// Constructor 00018 ATLAS_2011_I928289_W() 00019 : Analysis("ATLAS_2011_I928289_W") 00020 { 00021 setNeedsCrossSection(true); 00022 } 00023 00024 00025 public: 00026 00027 /// @name Analysis methods 00028 //@{ 00029 00030 /// Book histograms and initialise projections before the run 00031 void init() { 00032 00033 ///Initialise and register projections here 00034 FinalState fs; 00035 00036 Cut cut = pT >= 20*GeV; 00037 00038 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); 00039 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); 00040 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); 00041 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); 00042 00043 addProjection(wfinder_el_bare , "WFinder_el_bare"); 00044 addProjection(wfinder_el_dressed, "WFinder_el_dressed"); 00045 addProjection(wfinder_mu_bare , "WFinder_mu_bare"); 00046 addProjection(wfinder_mu_dressed, "WFinder_mu_dressed"); 00047 00048 /// Book histograms here 00049 _h_Wminus_lepton_eta_el_bare = bookHisto1D(3, 1, 1); 00050 _h_Wminus_lepton_eta_el_dressed = bookHisto1D(3, 1, 2); 00051 _h_Wminus_lepton_eta_mu_bare = bookHisto1D(3, 1, 3); 00052 _h_Wminus_lepton_eta_mu_dressed = bookHisto1D(3, 1, 4); 00053 _h_Wplus_lepton_eta_el_bare = bookHisto1D(5, 1, 1); 00054 _h_Wplus_lepton_eta_el_dressed = bookHisto1D(5, 1, 2); 00055 _h_Wplus_lepton_eta_mu_bare = bookHisto1D(5, 1, 3); 00056 _h_Wplus_lepton_eta_mu_dressed = bookHisto1D(5, 1, 4); 00057 _h_W_asym_eta_el_bare = bookScatter2D(7, 1, 1); 00058 _h_W_asym_eta_el_dressed = bookScatter2D(7, 1, 2); 00059 _h_W_asym_eta_mu_bare = bookScatter2D(7, 1, 3); 00060 _h_W_asym_eta_mu_dressed = bookScatter2D(7, 1, 4); 00061 00062 } 00063 00064 00065 /// Perform the per-event analysis 00066 void analyze(const Event& event) { 00067 00068 const double weight = event.weight(); 00069 00070 ///Do the event by event analysis here 00071 const WFinder& wfinder_el_bare = applyProjection<WFinder>(event, "WFinder_el_bare"); 00072 const WFinder& wfinder_el_dressed = applyProjection<WFinder>(event, "WFinder_el_dressed"); 00073 const WFinder& wfinder_mu_bare = applyProjection<WFinder>(event, "WFinder_mu_bare"); 00074 const WFinder& wfinder_mu_dressed = applyProjection<WFinder>(event, "WFinder_mu_dressed"); 00075 00076 FillPlots1d(wfinder_el_bare , _h_Wplus_lepton_eta_el_bare , _h_Wminus_lepton_eta_el_bare , weight); 00077 FillPlots1d(wfinder_el_dressed, _h_Wplus_lepton_eta_el_dressed, _h_Wminus_lepton_eta_el_dressed, weight); 00078 FillPlots1d(wfinder_mu_bare , _h_Wplus_lepton_eta_mu_bare , _h_Wminus_lepton_eta_mu_bare , weight); 00079 FillPlots1d(wfinder_mu_dressed, _h_Wplus_lepton_eta_mu_dressed, _h_Wminus_lepton_eta_mu_dressed, weight); 00080 00081 00082 } 00083 00084 void FillPlots1d(const WFinder& wfinder, Histo1DPtr hist_plus, Histo1DPtr hist_minus, double weight) { 00085 00086 if (wfinder.bosons().size() != 1) return; 00087 00088 Particle l = wfinder.constituentLeptons()[0]; 00089 const FourMomentum& miss = wfinder.constituentNeutrinos()[0].momentum(); 00090 00091 if(l.momentum().pT() > 20*GeV && miss.Et() > 25*GeV && wfinder.mT() > 40*GeV) { 00092 int lepCharge = l.charge(); 00093 if (lepCharge > 0) hist_plus ->fill( fabs(l.eta()), weight); 00094 else if (lepCharge < 0) hist_minus->fill( fabs(l.eta()), weight); 00095 } 00096 00097 return; 00098 } 00099 00100 00101 /// Normalise histograms etc., after the run 00102 void finalize() { 00103 00104 // Construct asymmetry: (dsig+/deta - dsig-/deta) / (dsig+/deta + dsig-/deta) 00105 divide(*_h_Wplus_lepton_eta_el_bare - *_h_Wminus_lepton_eta_el_bare, 00106 *_h_Wplus_lepton_eta_el_bare + *_h_Wminus_lepton_eta_el_bare, 00107 _h_W_asym_eta_el_bare); 00108 divide(*_h_Wplus_lepton_eta_el_dressed - *_h_Wminus_lepton_eta_el_dressed, 00109 *_h_Wplus_lepton_eta_el_dressed + *_h_Wminus_lepton_eta_el_dressed, 00110 _h_W_asym_eta_el_dressed); 00111 divide(*_h_Wplus_lepton_eta_mu_bare - *_h_Wminus_lepton_eta_mu_bare, 00112 *_h_Wplus_lepton_eta_mu_bare + *_h_Wminus_lepton_eta_mu_bare, 00113 _h_W_asym_eta_mu_bare); 00114 divide(*_h_Wplus_lepton_eta_mu_dressed - *_h_Wminus_lepton_eta_mu_dressed, 00115 *_h_Wplus_lepton_eta_mu_dressed + *_h_Wminus_lepton_eta_mu_dressed, 00116 _h_W_asym_eta_mu_dressed); 00117 00118 // Print summary info 00119 const double xs_pb(crossSection() / picobarn); 00120 const double sumw(sumOfWeights()); 00121 MSG_INFO( "Cross-Section/pb : " << xs_pb ); 00122 MSG_INFO( "Sum of weights : " << sumw ); 00123 MSG_INFO( "nEvents : " << numEvents() ); 00124 00125 const double sf(0.5 * xs_pb / sumw); // 0.5 accounts for rapidity bin width 00126 00127 /// Normalise, scale and otherwise manipulate histograms here 00128 scale(_h_Wminus_lepton_eta_el_bare , sf); 00129 scale(_h_Wminus_lepton_eta_el_dressed, sf); 00130 scale(_h_Wminus_lepton_eta_mu_bare , sf); 00131 scale(_h_Wminus_lepton_eta_mu_dressed, sf); 00132 scale(_h_Wplus_lepton_eta_el_bare , sf); 00133 scale(_h_Wplus_lepton_eta_el_dressed , sf); 00134 scale(_h_Wplus_lepton_eta_mu_bare , sf); 00135 scale(_h_Wplus_lepton_eta_mu_dressed , sf); 00136 00137 } 00138 00139 //@} 00140 00141 00142 private: 00143 00144 // Data members like post-cuts event weight counters go here 00145 00146 00147 private: 00148 00149 /// @name Histograms 00150 //@{ 00151 Histo1DPtr _h_Wminus_lepton_eta_el_bare; 00152 Histo1DPtr _h_Wminus_lepton_eta_el_dressed; 00153 Histo1DPtr _h_Wminus_lepton_eta_mu_bare; 00154 Histo1DPtr _h_Wminus_lepton_eta_mu_dressed; 00155 Histo1DPtr _h_Wplus_lepton_eta_el_bare; 00156 Histo1DPtr _h_Wplus_lepton_eta_el_dressed; 00157 Histo1DPtr _h_Wplus_lepton_eta_mu_bare; 00158 Histo1DPtr _h_Wplus_lepton_eta_mu_dressed; 00159 Scatter2DPtr _h_W_asym_eta_el_bare; 00160 Scatter2DPtr _h_W_asym_eta_el_dressed; 00161 Scatter2DPtr _h_W_asym_eta_mu_bare; 00162 Scatter2DPtr _h_W_asym_eta_mu_dressed; 00163 00164 //@} 00165 00166 }; 00167 00168 00169 // The hook for the plugin system 00170 DECLARE_RIVET_PLUGIN(ATLAS_2011_I928289_W); 00171 00172 } Generated on Thu Mar 10 2016 08:29:46 for The Rivet MC analysis system by ![]() |