Rivet analyses referenceD0_2008_I791230Measurement of W charge asymmetry from D0 Run IIExperiment: D0 (Tevatron Run 2) Inspire ID: 791230 Status: VALIDATED Authors:
Beam energies: (980.0, 980.0) GeV Run details:
Measurement of the electron charge asymmetry in $p \bar p \to W + X \to e \nu_e + X$ events at a center of mass energy of 1.96 TeV. The asymmetry is measured as a function of the electron transverse momentum and pseudorapidity in the interval (-3.2, 3.2). This data is sensitive to proton parton distribution functions due to the valence asymmetry in the incoming quarks which produce the W. Initial state radiation should also affect the pT distribution. Source code: D0_2008_I791230.cc 1// -*- C++ -*-
2#include "Rivet/Analysis.hh"
3#include "Rivet/Projections/FinalState.hh"
4#include "Rivet/Projections/MissingMomentum.hh"
5#include "Rivet/Projections/LeptonFinder.hh"
6#include "Rivet/Projections/LeadingParticlesFinalState.hh"
7#include "Rivet/Projections/IdentifiedFinalState.hh"
8
9namespace Rivet {
10
11
12 /// @brief D0 Run II measurement of W charge asymmetry
13 ///
14 /// @author Andy Buckley
15 /// @author Gavin Hesketh
16 class D0_2008_I791230 : public Analysis {
17 public:
18
19 RIVET_DEFAULT_ANALYSIS_CTOR(D0_2008_I791230);
20
21
22 /// @name Analysis methods
23 /// @{
24
25 // Book histograms and set up projections
26 void init() {
27
28 // Projections
29 declare("MET", MissingMomentum());
30 LeptonFinder ef(0.2, Cuts::abseta < 5 && Cuts::pT > 25*GeV && Cuts::abspid == PID::ELECTRON);
31 declare(ef, "Elecs");
32
33 // Histograms (temporary +- charge histos and scatters to store the calculated asymmetries)
34 for (size_t pmindex = 0; pmindex < 2; ++pmindex) {
35 const string suffix = (pmindex == 0) ? "plus" : "minus";
36 book(_hs_dsigpm_deta_25_35[pmindex] ,"TMP/dsigpm_deta_25_35_" + suffix, refData(1, 1, 1));
37 book(_hs_dsigpm_deta_35[pmindex] ,"TMP/dsigpm_deta_35_" + suffix, refData(1, 1, 2));
38 book(_hs_dsigpm_deta_25[pmindex] ,"TMP/dsigpm_deta_25_" + suffix, refData(1, 1, 3));
39 }
40 book(_h_asym1, 1, 1, 1);
41 book(_h_asym2, 1, 1, 2);
42 book(_h_asym3, 1, 1, 3);
43 }
44
45
46 /// Do the analysis
47 void analyze(const Event & event) {
48
49 // W reco, starting with MET
50 const P4& pmiss = apply<MissingMom>(event, "MET").missingMom();
51 if (pmiss.pT() < 25*GeV) vetoEvent;
52
53 // Identify the closest-matching l+MET to m == mW
54 const Particles& es = apply<LeptonFinder>(event, "Elecs").particles();
55 const int ifound = closestMatchIndex(es, pmiss, Kin::mass, 80.4*GeV, 60*GeV, 100*GeV);
56 if (ifound < 0) {
57 MSG_DEBUG("No W candidates found: vetoing");
58 vetoEvent;
59 }
60
61 // Get the e+- momentum, and an effective charge including the eta sign
62 const Particle& e = es[ifound];
63 const int chg_e = sign(e.eta()) * sign(e.charge());
64 assert(chg_e == 1 || chg_e == -1);
65 MSG_TRACE("Charged lepton sign = " << chg_e);
66
67 // Fill histos with appropriate +- indexing
68 const size_t pmindex = (chg_e > 0) ? 0 : 1;
69 if (e.Et() < 35*GeV) _hs_dsigpm_deta_25_35[pmindex]->fill(e.abseta());
70 else _hs_dsigpm_deta_35[pmindex]->fill(e.abseta());
71 _hs_dsigpm_deta_25[pmindex]->fill(e.abseta());
72 }
73
74
75 /// @brief Finalize
76 ///
77 /// Construct asymmetry: (dsig+/deta - dsig-/deta) / (dsig+/deta + dsig-/deta) for each ET region
78 void finalize() {
79 asymm(_hs_dsigpm_deta_25_35[0], _hs_dsigpm_deta_25_35[1], _h_asym1);
80 asymm(_hs_dsigpm_deta_35[0], _hs_dsigpm_deta_35[1], _h_asym2);
81 asymm(_hs_dsigpm_deta_25[0], _hs_dsigpm_deta_25[1], _h_asym3);
82 _h_asym1->scale(100.);
83 _h_asym2->scale(100.);
84 _h_asym3->scale(100.);
85 }
86
87 /// @}
88
89
90 private:
91
92 /// @name Histograms
93 /// @{
94 Histo1DPtr _hs_dsigpm_deta_25_35[2], _hs_dsigpm_deta_35[2], _hs_dsigpm_deta_25[2];
95 Estimate1DPtr _h_asym1, _h_asym2, _h_asym3;
96 /// @}
97
98 };
99
100
101 RIVET_DECLARE_ALIASED_PLUGIN(D0_2008_I791230, D0_2008_S7837160);
102
103}
|