Rivet Analyses Reference

ATLAS_2016_I1444991

Higgs-to-WW differential cross sections at 8 TeV
Experiment: ATLAS (LHC)
Inspire ID: 1444991
Status: VALIDATED
Authors:
  • Kathrin Becker
References:Beams: p+ p+
Beam energies: (4000.0, 4000.0) GeV
Run details:
  • gg -> H -> W W* -> enu munu production at 8 TeV

This paper describes a measurement of fiducial and differential cross sections of gluon-fusion Higgs boson production in the $H\rightarrow W W^\ast \rightarrow e\nu\mu\nu$ channel, using 20.3 fb$^{-1}$ of proton-proton collision data. The data were produced at a centre-of-mass energy of $\sqrt{s} = 8$ TeV at the CERN Large Hadron Collider and recorded by the ATLAS detector in 2012. Cross sections are measured from the observed $H\rightarrow W W^\ast \rightarrow e\nu\mu\nu$ signal yield in categories distinguished by the number of associated jets. The total cross section is measured in a fiducial region defined by the kinematic properties of the charged leptons and neutrinos. Differential cross sections are reported as a function of the number of jets, the Higgs boson transverse momentum, the dilepton rapidity, and the transverse momentum of the leading jet. The jet-veto efficiency, or fraction of events with no jets above a given transverse momentum threshold, is also reported. All measurements are compared to QCD predictions from Monte Carlo generators and fixed-order calculations, and are in agreement with the Standard Model predictions.

Source code: ATLAS_2016_I1444991.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
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/DressedLeptons.hh"
#include "Rivet/Projections/IdentifiedFinalState.hh"
#include "Rivet/Projections/PromptFinalState.hh"
#include "Rivet/Projections/VetoedFinalState.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Projections/VisibleFinalState.hh"

namespace Rivet {


  class ATLAS_2016_I1444991 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2016_I1444991);


  public:

    /// Book histograms and initialise projections before the run
    void init() {

      // All particles within |eta| < 5.0
      const FinalState FS(Cuts::abseta < 5.0);

      // Project photons for dressing
      IdentifiedFinalState photon_id(FS);
      photon_id.acceptIdPair(PID::PHOTON);

      // Project dressed electrons with pT > 15 GeV and |eta| < 2.47
      IdentifiedFinalState el_id(FS);
      el_id.acceptIdPair(PID::ELECTRON);
      PromptFinalState el_bare(el_id);
      Cut cuts = (Cuts::abseta < 2.47) && ( (Cuts::abseta <= 1.37) || (Cuts::abseta >= 1.52) ) && (Cuts::pT > 15*GeV);
      DressedLeptons el_dressed_FS(photon_id, el_bare, 0.1, cuts, true);
      declare(el_dressed_FS,"EL_DRESSED_FS");

      // Project dressed muons with pT > 15 GeV and |eta| < 2.5
      IdentifiedFinalState mu_id(FS);
      mu_id.acceptIdPair(PID::MUON);
      PromptFinalState mu_bare(mu_id);
      DressedLeptons mu_dressed_FS(photon_id, mu_bare, 0.1, Cuts::abseta < 2.5 && Cuts::pT > 15*GeV, true);
      declare(mu_dressed_FS,"MU_DRESSED_FS");

      // get MET from generic invisibles
      VetoedFinalState inv_fs(FS);
      inv_fs.addVetoOnThisFinalState(VisibleFinalState(FS));
      declare(inv_fs, "InvisibleFS");

      // Project jets
      FastJets jets(FS, FastJets::ANTIKT, 0.4);
      jets.useInvisibles(JetAlg::Invisibles::NONE);
      jets.useMuons(JetAlg::Muons::NONE);
      declare(jets, "jets");

      // Book histograms
      book(_h_Njets        , 2,1,1);
      book(_h_PtllMET      , 3,1,1);
      book(_h_Yll          , 4,1,1);
      book(_h_PtLead       , 5,1,1);
      book(_h_Njets_norm   , 6,1,1);
      book(_h_PtllMET_norm , 7,1,1);
      book(_h_Yll_norm     , 8,1,1);
      book(_h_PtLead_norm  , 9,1,1);
      book(_h_JetVeto      , 10, 1, 1, true);

      //histos for jetveto
      std::vector<double> ptlead25_bins = { 0., 25., 300. };
      std::vector<double> ptlead40_bins = { 0., 40., 300. };
      book(_h_pTj1_sel25 , "pTj1_sel25", ptlead25_bins);
      book(_h_pTj1_sel40 , "pTj1_sel40", ptlead40_bins);
    }


    /// Perform the per-event analysis
    void analyze(const Event& event) {

      // Get final state particles
      const FinalState& ifs = applyProjection<FinalState>(event, "InvisibleFS");
      const vector<DressedLepton>& good_mu = applyProjection<DressedLeptons>(event, "MU_DRESSED_FS").dressedLeptons();
      const vector<DressedLepton>& el_dressed = applyProjection<DressedLeptons>(event, "EL_DRESSED_FS").dressedLeptons();
      const Jets& jets = applyProjection<FastJets>(event, "jets").jetsByPt(Cuts::pT>25*GeV && Cuts::abseta < 4.5);

      //find good electrons
      vector<DressedLepton> good_el;
      for (const DressedLepton& el : el_dressed){
        bool keep = true;
        for (const DressedLepton& mu : good_mu) {
          keep &= deltaR(el, mu) >= 0.1;
        }
        if (keep)  good_el += el;
      }

      // select only emu events
      if ((good_el.size() != 1) || good_mu.size() != 1)  vetoEvent;

      //built dilepton
      FourMomentum dilep = good_el[0].momentum() + good_mu[0].momentum();
      double Mll = dilep.mass();
      double Yll = dilep.rapidity();
      double DPhill = fabs(deltaPhi(good_el[0], good_mu[0]));
      double pTl1 = (good_el[0].pT() > good_mu[0].pT())? good_el[0].pT() : good_mu[0].pT();

      //get MET
      FourMomentum met;
      for (const Particle& p : ifs.particles())  met += p.momentum();

      // do a few cuts before looking at jets
      if (pTl1 <= 22. || DPhill >= 1.8 || met.pT() <= 20.)  vetoEvent;
      if (Mll <= 10. || Mll >= 55.)  vetoEvent;

      Jets jets_selected;
      for (const Jet &j : jets) {
        if( j.abseta() > 2.4 && j.pT()<=30*GeV ) continue;
        bool keep = true;
        for(DressedLepton el : good_el) {
          keep &= deltaR(j, el) >= 0.3;
        }
        if (keep)  jets_selected += j;
      }

      double PtllMET = (met + good_el[0].momentum() + good_mu[0].momentum()).pT();

      double Njets = jets_selected.size() > 2 ? 2 : jets_selected.size();
      double pTj1 = jets_selected.size()? jets_selected[0].pT() : 0.1;

      // Fill histograms
      _h_Njets->fill(Njets);
      _h_PtllMET->fill(PtllMET);
      _h_Yll->fill(fabs(Yll));
      _h_PtLead->fill(pTj1);
      _h_Njets_norm->fill(Njets);
      _h_PtllMET_norm->fill(PtllMET);
      _h_Yll_norm->fill(fabs(Yll));
      _h_PtLead_norm->fill(pTj1);
      _h_pTj1_sel25->fill(pTj1);
      _h_pTj1_sel40->fill(pTj1);
    }


    /// Normalise histograms etc., after the run
    void finalize() {

      const double xs = crossSectionPerEvent()/femtobarn;

      /// @todo Normalise, scale and otherwise manipulate histograms here
      scale(_h_Njets, xs);
      scale(_h_PtllMET, xs);
      scale(_h_Yll, xs);
      scale(_h_PtLead, xs);
      normalize(_h_Njets_norm);
      normalize(_h_PtllMET_norm);
      normalize(_h_Yll_norm);
      normalize(_h_PtLead_norm);
      scale(_h_pTj1_sel25, xs);
      scale(_h_pTj1_sel40, xs);
      normalize(_h_pTj1_sel25);
      normalize(_h_pTj1_sel40);
      // fill jet veto efficiency histogram
      _h_JetVeto->point(0).setY(_h_pTj1_sel25->bin(0).sumW(), sqrt(_h_pTj1_sel25->bin(0).sumW2()));
      _h_JetVeto->point(1).setY(_h_PtLead_norm->bin(0).sumW(), sqrt(_h_PtLead_norm->bin(0).sumW2()));
      _h_JetVeto->point(2).setY(_h_pTj1_sel40->bin(0).sumW(), sqrt(_h_pTj1_sel25->bin(0).sumW2()));

      scale(_h_PtLead_norm , 1000.); // curveball unit change in HepData, just for this one
      scale(_h_PtllMET_norm, 1000.); // curveball unit change in HepData, and this one
    }

  private:

    /// @name Histograms
    //@{
    Histo1DPtr _h_Njets;
    Histo1DPtr _h_PtllMET;
    Histo1DPtr _h_Yll;
    Histo1DPtr _h_PtLead;
    Histo1DPtr _h_Njets_norm;
    Histo1DPtr _h_PtllMET_norm;
    Histo1DPtr _h_Yll_norm;
    Histo1DPtr _h_PtLead_norm;

    Scatter2DPtr _h_JetVeto;

    Histo1DPtr _h_pTj1_sel25;
    Histo1DPtr _h_pTj1_sel40;

  };


  // The hook for the plugin system
  RIVET_DECLARE_PLUGIN(ATLAS_2016_I1444991);

}