Rivet Analyses Reference

ATLAS_2015_I1390114

$tt+b(b)$ at 8 TeV
Experiment: ATLAS (LHC)
Inspire ID: 1390114
Status: VALIDATED
Authors:
  • Matthias Danninger
  • Christian Gutschow
References:Beams: p+ p+
Beam energies: (4000.0, 4000.0) GeV
Run details:
  • Top-antitop production

Fiducial cross-sections for $t\bar{t}$ production with one or two additional $b$-jets are presented, using an integrated luminosity of 20.3 fb${}^{-1}$ of proton-proton collisions at a centre-of-mass energy of 8 TeV at the Large Hadron Collider, collected with the ATLAS detector. The cross-section times branching ratio for $t\bar{t}$ events with at least one additional $b$-jet is measured to be 950 $\pm$ 70 (stat.) ${}^{+240}_{-190}$ (syst.) fb in the lepton-plus-jets channel and 50 $\pm$ 10 (stat.) ${}^{+15}_{-10}$ (syst.) fb in the $e\mu$ channel. The cross-section times branching ratio for events with at least two additional $b$-jets is measured to be 19.3 $\pm$ 3.5 (stat.) $\pm$ 5.7 (syst.) fb in the dilepton channel ($e\mu$, $\mu\mu$, $ee$) using a method based on tight selection criteria, and 13.5 $\pm$ 3.3 (stat.) $\pm$ 3.6 (syst.) fb using a looser selection that allows the background normalisation to be extracted from data. The latter method also measures a value of 1.30 $\pm$ 0.33 (stat.) $\pm$ 0.28 (syst.)$\%$ for the ratio of $t\bar{t}$ production with two additional $b$-jets to $t\bar{t}$ production with any two additional jets. NB: When merging several output YODA files, the ratio (d02-x01-y01) needs to be reconstructed in a post-processing step.

Source code: ATLAS_2015_I1390114.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
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/VetoedFinalState.hh"
#include "Rivet/Projections/IdentifiedFinalState.hh"
#include "Rivet/Projections/PromptFinalState.hh"
#include "Rivet/Projections/DressedLeptons.hh"
#include "Rivet/Projections/FastJets.hh"

namespace Rivet {


  // ATLAS $tt+b(b)$ at 8~TeV
  class ATLAS_2015_I1390114 : public Analysis {
  public:

    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2015_I1390114);


    void init() {

      // Eta ranges
      Cut eta_full = Cuts::abseta < 5.0 && Cuts::pT >= 1.0*MeV;
      Cut eta_lep  = Cuts::abseta < 2.5;

      // All final state particles
      FinalState fs(eta_full);

      // Get photons to dress leptons
      IdentifiedFinalState photons(fs);
      photons.acceptIdPair(PID::PHOTON);

      // Projection to find the electrons
      IdentifiedFinalState el_id(fs);
      el_id.acceptIdPair(PID::ELECTRON);
      PromptFinalState electrons(el_id);
      electrons.acceptTauDecays(true);
      declare(electrons, "electrons");
      DressedLeptons dressedelectrons(photons, electrons, 0.1, eta_lep && Cuts::pT > 25*GeV, true);
      declare(dressedelectrons, "dressedelectrons");
      DressedLeptons ewdressedelectrons(photons, electrons, 0.1, eta_full, true);

      // Projection to find the muons
      IdentifiedFinalState mu_id(fs);
      mu_id.acceptIdPair(PID::MUON);
      PromptFinalState muons(mu_id);
      muons.acceptTauDecays(true);
      declare(muons, "muons");
      DressedLeptons dressedmuons(photons, muons, 0.1, eta_lep && Cuts::pT > 25*GeV, true);
      declare(dressedmuons, "dressedmuons");
      DressedLeptons ewdressedmuons(photons, muons, 0.1, eta_full, true);

      // Projection to find neutrinos and produce MET
      IdentifiedFinalState nu_id;
      nu_id.acceptNeutrinos();
      PromptFinalState neutrinos(nu_id);
      neutrinos.acceptTauDecays(true);

      // Jet clustering.
      VetoedFinalState vfs;
      vfs.addVetoOnThisFinalState(ewdressedelectrons);
      vfs.addVetoOnThisFinalState(ewdressedmuons);
      vfs.addVetoOnThisFinalState(neutrinos);
      FastJets jets(vfs,FastJets::ANTIKT, 0.4);
      jets.useInvisibles();
      declare(jets, "jets");

      book(_histo ,1,1,1);
      book(_ratio, 2,1,1, true);
      book(_aux   ,"_aux", 1, 0.5, 1.5);
    }


    void analyze(const Event& event) {

      // Get the selected objects, using the projections.
      vector<DressedLepton> electrons = apply<DressedLeptons>(event, "dressedelectrons").dressedLeptons();
      vector<DressedLepton> muons = apply<DressedLeptons>(event, "dressedmuons").dressedLeptons();
      if (electrons.empty() && muons.empty())  vetoEvent;

      Jets jets = apply<FastJets>(event, "jets").jets((Cuts::pT>20*GeV) && (Cuts::abseta < 2.5), cmpMomByPt);
      Jets bjets;

      // Check overlap of jets/leptons.
      for (Jet jet : jets) {
        // if dR(el,jet) < 0.4 skip the event
        for (DressedLepton el : electrons) {
          if (deltaR(jet, el) < 0.4)  vetoEvent;
        }
        // if dR(mu,jet) < 0.4 skip the event
        for (DressedLepton mu : muons) {
          if (deltaR(jet, mu) < 0.4)  vetoEvent;
        }
        // Count the number of b-tags
        // We have to check that the ghost-matched B hadrons have pT > 5 GeV
        // By default jet.bTags() returns all B hadrons without cuts
        bool btagged = jet.bTags(Cuts::pT >= 5*GeV).size();
        if (btagged)  bjets += jet;
      }

      // Evaluate basic event selection
      bool pass_1lep = (electrons.size() == 1 && muons.empty()) || (muons.size() == 1 && electrons.empty());
      if (pass_1lep && bjets.size() >= 3 && jets.size() >= 5)  _histo->fill(1);

      if (muons.size() == 1 && electrons.size() == 1 && bjets.size() >= 3) {
        if (muons[0].charge() * electrons[0].charge() < 0.0)  _histo->fill(2);
      }

      DressedLepton *lep1 = nullptr, *lep2 = nullptr;
      bool zveto = false;
      if (electrons.size() == 2 && muons.empty()) {
        lep1 = &electrons[0];
        lep2 = &electrons[1];
      }
      else if (muons.size() == 2 && electrons.empty()) {
        lep1 = &muons[0];
        lep2 = &muons[1];
      }
      else if (electrons.size() == 1 && muons.size() == 1) {
        lep1 = &electrons[0];
        lep2 = &muons[0];
        zveto = true;
      }
      if (lep1 && lep2 && lep1->charge() * lep2->charge() < 0.0 && bjets.size() >= 2) {
        double mass = (lep1->momentum() + lep2->momentum()).mass();
        bool pass_2lep = mass > 15*GeV && (zveto || !(mass > 81*GeV && mass < 101*GeV));
        if ( pass_2lep && bjets.size() >= 4) {
          _histo->fill(3);
          _histo->fill(4);
        }
        if ( pass_2lep && jets.size() >= 4) {
          _aux->fill(1);
        }
      }
    }


    void finalize() {
      const double sf(crossSection()/sumOfWeights()/femtobarn);
      scale(_histo, sf);
      scale(_aux,  sf);

      // construct ratio
      const double  n = _histo->bin(3).sumW();
      const double dN = _histo->bin(3).sumW2();
      const double  d = _aux->bin(0).sumW();
      const double dD = _aux->bin(0).sumW2();
      const double  r = safediv(n, d);
      double e = sqrt( safediv(r * (1 - r), d) );
      if ( _aux->effNumEntries() != _aux->numEntries() ) {
        // use F. James's approximation for weighted events:
        e = sqrt( safediv((1 - 2 * r) * dN + r * r * dD, d * d) );
      }
      _ratio->point(0).setY(100.0 * r, 100.0 * e); // convert into percentage
    }


  private:

    Histo1DPtr _histo, _aux;
    Scatter2DPtr _ratio;

  };


  // Hook for the plugin system
  RIVET_DECLARE_PLUGIN(ATLAS_2015_I1390114);

}