D0_2008_S7837160.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 #include "Rivet/Projections/LeadingParticlesFinalState.hh" 00006 #include "Rivet/Projections/IdentifiedFinalState.hh" 00007 00008 namespace Rivet { 00009 00010 00011 /// @brief D0 Run II measurement of W charge asymmetry 00012 /// @author Andy Buckley 00013 /// @author Gavin Hesketh 00014 class D0_2008_S7837160 : public Analysis { 00015 public: 00016 00017 /// Default constructor. 00018 D0_2008_S7837160() 00019 : Analysis("D0_2008_S7837160") 00020 { } 00021 00022 00023 /// @name Analysis methods 00024 //@{ 00025 00026 // Book histograms and set up projections 00027 void init() { 00028 // Projections 00029 FinalState fs; 00030 /// @todo Use separate pT and ETmiss cuts in WFinder 00031 const WFinder wfe(fs, Cuts::abseta < 5 && Cuts::pT > 25*GeV, PID::ELECTRON, 60.0*GeV, 100.0*GeV, 25.0*GeV, 0.2); 00032 addProjection(wfe, "WFe"); 00033 00034 // Histograms (temporary +- charge histos and scatters to store the calculated asymmetries) 00035 for (size_t pmindex = 0; pmindex < 2; ++pmindex) { 00036 const string suffix = (pmindex == 0) ? "plus" : "minus"; 00037 _hs_dsigpm_deta_25_35[pmindex] = bookHisto1D("TMP/dsigpm_deta_25_35_" + suffix, refData(1, 1, 1)); 00038 _hs_dsigpm_deta_35[pmindex] = bookHisto1D("TMP/dsigpm_deta_35_" + suffix, refData(1, 1, 2)); 00039 _hs_dsigpm_deta_25[pmindex] = bookHisto1D("TMP/dsigpm_deta_25_" + suffix, refData(1, 1, 3)); 00040 } 00041 _h_asym1 = bookScatter2D(1, 1, 1); 00042 _h_asym2 = bookScatter2D(1, 1, 2); 00043 _h_asym3 = bookScatter2D(1, 1, 3); 00044 } 00045 00046 00047 /// Do the analysis 00048 void analyze(const Event & event) { 00049 const WFinder& wf = applyProjection<WFinder>(event, "WFe"); 00050 if (wf.bosons().size() == 0) { 00051 MSG_DEBUG("No W candidates found: vetoing"); 00052 vetoEvent; 00053 } 00054 00055 // Get the e+- momentum, and an effective charge including the eta sign 00056 /// @todo Is it correct to multiply the eta sign into the charge to "fold" the plot? 00057 const FourMomentum p_e = wf.constituentLeptons()[0].momentum(); 00058 const int chg_e = sign(p_e.eta()) * sign(charge(wf.constituentLeptons()[0])); 00059 assert(chg_e == 1 || chg_e == -1); 00060 MSG_TRACE("Charged lepton sign = " << chg_e); 00061 00062 // Fill histos with appropriate +- indexing 00063 const double weight = event.weight(); 00064 const size_t pmindex = (chg_e > 0) ? 0 : 1; 00065 if (p_e.Et() < 35*GeV) _hs_dsigpm_deta_25_35[pmindex]->fill(fabs(p_e.eta()), weight); 00066 else _hs_dsigpm_deta_35[pmindex]->fill(fabs(p_e.eta()), weight); 00067 _hs_dsigpm_deta_25[pmindex]->fill(fabs(p_e.eta()), weight); 00068 } 00069 00070 00071 /// @name Helper functions for constructing asymmetry histograms in finalize() 00072 //@{ 00073 void calc_asymm(const Histo1DPtr plus, const Histo1DPtr minus, Scatter2DPtr target) { 00074 divide(*plus - *minus, *plus + *minus, target); 00075 } 00076 void calc_asymm(const Histo1DPtr histos[2], Scatter2DPtr target) { 00077 calc_asymm(histos[0], histos[1], target); 00078 } 00079 //@} 00080 00081 00082 /// @brief Finalize 00083 /// 00084 /// Construct asymmetry: (dsig+/deta - dsig-/deta) / (dsig+/deta + dsig-/deta) for each ET region 00085 void finalize() { 00086 calc_asymm(_hs_dsigpm_deta_25_35, _h_asym1); 00087 calc_asymm(_hs_dsigpm_deta_35, _h_asym2); 00088 calc_asymm(_hs_dsigpm_deta_25, _h_asym3); 00089 _h_asym1->scale(1.,100.); 00090 _h_asym2->scale(1.,100.); 00091 _h_asym3->scale(1.,100.); 00092 } 00093 00094 //@} 00095 00096 00097 private: 00098 00099 /// @name Histograms 00100 //@{ 00101 Histo1DPtr _hs_dsigpm_deta_25_35[2], _hs_dsigpm_deta_35[2], _hs_dsigpm_deta_25[2]; 00102 Scatter2DPtr _h_asym1, _h_asym2, _h_asym3; 00103 //@} 00104 00105 }; 00106 00107 00108 00109 // The hook for the plugin system 00110 DECLARE_RIVET_PLUGIN(D0_2008_S7837160); 00111 00112 } Generated on Thu Mar 10 2016 08:29:49 for The Rivet MC analysis system by ![]() |