Rivet Analyses Reference

MC_WINC

Monte Carlo validation observables for inclusive $W[\ell \, \nu]$ production
Experiment: ()
Status: VALIDATED
Authors:
  • Frank Siegert
  • Christian Gutschow
No references listed
Beams: * *
Beam energies: ANY
Run details:
  • $\ell \, \nu$ + jets analysis.

Monte Carlo validation observables for inclusive $W[\ell \, \nu]$ production. 60 GeV $<m_W<$ 100 GeV cut applied.

Source code: MC_WINC.cc
  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/WFinder.hh"

namespace Rivet {



  /// @brief MC validation analysis for inclusive W events
  class MC_WINC : public Analysis {
  public:

    /// Default constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(MC_WINC);

    /// @name Analysis methods
    //@{

    /// Book histograms
    void init() {
		  _dR=0.2;
      if (getOption("SCHEME") == "BARE")  _dR = 0.0;
		  _lepton=PID::ELECTRON;
      if (getOption("LMODE") == "MU")  _lepton = PID::MUON;

      FinalState fs;
      WFinder wfinder(fs, Cuts::abseta < 3.5 && Cuts::pT > 25*GeV, _lepton, 60.0*GeV, 100.0*GeV, 25.0*GeV, _dR);
      declare(wfinder, "WFinder");

      double sqrts = sqrtS()>0. ? sqrtS() : 14000.;
      book(_h_W_mass ,"W_mass", 50, 55.0, 105.0);
      book(_h_W_mT ,"W_mT", 40, 60.0, 100.0);
      book(_h_W_pT ,"W_pT", logspace(100, 1.0, 0.5*sqrts));
      book(_h_W_pT_peak ,"W_pT_peak", 25, 0.0, 125.0);
      book(_h_W_y ,"W_y", 40, -4.0, 4.0);
      book(_h_W_phi ,"W_phi", 25, 0.0, TWOPI);
      book(_h_Wplus_pT ,"Wplus_pT", logspace(25, 10.0, 0.5*sqrts));
      book(_h_Wminus_pT ,"Wminus_pT", logspace(25, 10.0, 0.5*sqrts));
      book(_h_lepton_pT ,"lepton_pT", logspace(100, 10.0, 0.25*sqrts));
      book(_h_lepton_eta ,"lepton_eta", 40, -4.0, 4.0);
      book(_htmp_dsigminus_deta ,"lepton_dsigminus_deta", 20, 0.0, 4.0);
      book(_htmp_dsigplus_deta  ,"lepton_dsigplus_deta", 20, 0.0, 4.0);

      book(_h_asym, "W_chargeasymm_eta");
      book(_h_asym_pT, "W_chargeasymm_pT");
    }



    /// Do the analysis
    void analyze(const Event & e) {

      const WFinder& wfinder = apply<WFinder>(e, "WFinder");
      if (wfinder.bosons().size() != 1) {
        vetoEvent;
      }

      double charge3_x_eta = 0;
      int charge3 = 0;
      FourMomentum emom;
      FourMomentum wmom(wfinder.bosons().front().momentum());
      _h_W_mass->fill(wmom.mass()/GeV);
      _h_W_mT->fill(wfinder.mT()/GeV);
      _h_W_pT->fill(wmom.pT()/GeV);
      _h_W_pT_peak->fill(wmom.pT()/GeV);
      _h_W_y->fill(wmom.rapidity());
      _h_W_phi->fill(wmom.phi());
      Particle l=wfinder.constituentLeptons()[0];
      _h_lepton_pT->fill(l.pT()/GeV);
      _h_lepton_eta->fill(l.eta());
      if (PID::charge3(l.pid()) != 0) {
        emom = l.momentum();
        charge3_x_eta = PID::charge3(l.pid()) * emom.eta();
        charge3 = PID::charge3(l.pid());
      }
      assert(charge3_x_eta != 0);
      assert(charge3!=0);
      if (emom.Et() > 30/GeV) {
        if (charge3_x_eta < 0) {
          _htmp_dsigminus_deta->fill(emom.eta());
        } else {
          _htmp_dsigplus_deta->fill(emom.eta());
        }
      }
      if (charge3 < 0) {
        _h_Wminus_pT->fill(wmom.pT()/GeV);
      } else {
        _h_Wplus_pT->fill(wmom.pT()/GeV);
      }
    }


    /// Finalize
    void finalize() {
      scale(_h_W_mass, crossSection()/sumOfWeights());
      scale(_h_W_mT, crossSection()/sumOfWeights());
      scale(_h_W_pT, crossSection()/sumOfWeights());
      scale(_h_W_pT_peak, crossSection()/sumOfWeights());
      scale(_h_W_y, crossSection()/sumOfWeights());
      scale(_h_W_phi, crossSection()/sumOfWeights());
      scale(_h_lepton_pT, crossSection()/sumOfWeights());
      scale(_h_lepton_eta, crossSection()/sumOfWeights());

      // Construct asymmetry: (dsig+/deta - dsig-/deta) / (dsig+/deta + dsig-/deta) for each Et region
      divide(*_htmp_dsigplus_deta - *_htmp_dsigminus_deta,
             *_htmp_dsigplus_deta + *_htmp_dsigminus_deta,
             _h_asym);

      // // W charge asymmetry vs. pTW: dsig+/dpT / dsig-/dpT
      divide(_h_Wplus_pT, _h_Wminus_pT,
             _h_asym_pT);

      scale(_h_Wplus_pT, crossSection()/picobarn/sumOfWeights());
      scale(_h_Wminus_pT, crossSection()/picobarn/sumOfWeights());
    }

    //@}


  protected:

    /// @name Parameters for specialised e/mu and dressed/bare subclassing
    //@{
    double _dR;
    PdgId _lepton;
    //@}


  private:

    /// @name Histograms
    //@{
    Histo1DPtr _h_W_mass;
    Histo1DPtr _h_W_mT;
    Histo1DPtr _h_W_pT;
    Histo1DPtr _h_W_pT_peak;
    Histo1DPtr _h_W_y;
    Histo1DPtr _h_W_phi;
    Histo1DPtr _h_Wplus_pT;
    Histo1DPtr _h_Wminus_pT;
    Histo1DPtr _h_lepton_pT;
    Histo1DPtr _h_lepton_eta;

    Histo1DPtr _htmp_dsigminus_deta;
    Histo1DPtr _htmp_dsigplus_deta;

    Scatter2DPtr _h_asym;
    Scatter2DPtr _h_asym_pT;
    //@}

  };

  // The hooks for the plugin system
  RIVET_DECLARE_PLUGIN(MC_WINC);
}