Rivet Analyses Reference

ATLAS_2014_I1306294

Measurement of Z boson in association with b-jets at 7 TeV in ATLAS (electron channel)
Experiment: ATLAS (LHC)
Inspire ID: 1306294
Status: VALIDATED
Authors:
  • Gavin Hesketh
References:Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • Z+b(b) production in $pp$ collisions at $\sqrt{s} = 7$ TeV, electronic Z-decays

Measurements of differential production cross-sections of a $Z$ boson in association with $b$-jets in $pp$ collisions at $\sqrt{s}=7$ TeV are reported. The data analysed correspond to an integrated luminosity of 4.6 fb$^{-1}$ recorded with the ATLAS detector at the Large Hadron Collider. Particle-level cross-sections are determined for events with a $Z$ boson decaying into an electron or muon pair, and containing $b$-jets. For events with at least one $b$-jet, the cross-section is presented as a function of the $Z$ boson transverse momentum and rapidity, together with the inclusive $b$-jet cross-section as a function of $b$-jet transverse momentum, rapidity and angular separations between the $b$-jet and the $Z$ boson. For events with at least two $b$-jets, the cross-section is determined as a function of the invariant mass and angular separation of the two highest transverse momentum $b$-jets, and as a function of the $Z$ boson transverse momentum and rapidity. Results are compared to leading-order and next-to-leading-order perturbative QCD calculations. This Rivet module implements the event selection for Z decaying into electrons. If you want to use muonic events, please refer to ATLAS\_2014\_I1306294\_MU

Source code: ATLAS_2014_I1306294.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
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/ZFinder.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Projections/HeavyHadrons.hh"
#include "Rivet/Projections/VetoedFinalState.hh"

namespace Rivet {



  /// Electroweak Wjj production at 8 TeV
  class ATLAS_2014_I1306294 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2014_I1306294);


    /// @name Analysis methods
    //@{

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

      // Get options from the new option system
      _mode = 1;
      if ( getOption("LMODE") == "EL" ) _mode = 1;
      if ( getOption("LMODE") == "MU" ) _mode = 2;

      FinalState fs;
      Cut cuts = Cuts::abseta < 2.5 && Cuts::pT > 20*GeV;

      ZFinder zfinder(fs, cuts, _mode==1? PID::ELECTRON : PID::MUON, 76.0*GeV, 106.0*GeV, 0.1, 
                      ZFinder::ChargedLeptons::ALL, ZFinder::ClusterPhotons::NODECAY, ZFinder::AddPhotons::NO);
      declare(zfinder, "ZFinder");

      VetoedFinalState jet_fs(fs);
      jet_fs.addVetoOnThisFinalState(getProjection<ZFinder>("ZFinder"));
      FastJets jetpro1(jet_fs, FastJets::ANTIKT, 0.4, JetAlg::Muons::ALL, JetAlg::Invisibles::ALL);
      declare(jetpro1, "AntiKtJets04");
      declare(HeavyHadrons(), "BHadrons");

      // Histograms with data binning
      book(_h_bjet_Pt      , 3, 1, 1);
      book(_h_bjet_Y       , 5, 1, 1);
      book(_h_bjet_Yboost  , 7, 1, 1);
      book(_h_bjet_DY20    , 9, 1, 1);
      book(_h_bjet_ZdPhi20 ,11, 1, 1);
      book(_h_bjet_ZdR20   ,13, 1, 1);
      book(_h_bjet_ZPt     ,15, 1, 1);
      book(_h_bjet_ZY      ,17, 1, 1);
      book(_h_2bjet_dR     ,21, 1, 1);
      book(_h_2bjet_Mbb    ,23, 1, 1);
      book(_h_2bjet_ZPt    ,25, 1, 1);
      book(_h_2bjet_ZY     ,27, 1, 1);
    }


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

      // Check we have a Z:
      const ZFinder& zfinder = apply<ZFinder>(e, "ZFinder");
      if (zfinder.bosons().size() != 1) vetoEvent;
      
      const Particles boson_s =  zfinder.bosons();
      const Particle boson_f =  boson_s[0];
      const Particles zleps   =  zfinder.constituents();

      // Stop processing the event if no true b-partons or hadrons are found
      const Particles allBs = apply<HeavyHadrons>(e, "BHadrons").bHadrons(5.0*GeV);
      Particles stableBs = filter_select(allBs, Cuts::abseta < 2.5);
      if (stableBs.empty()) vetoEvent;

      // Get the b-jets
      const Jets& jets = apply<JetAlg>(e, "AntiKtJets04").jetsByPt(Cuts::pT >20.0*GeV && Cuts::abseta <2.4);
      Jets b_jets;
      for (const Jet& jet : jets) {
        //veto overlaps with Z leptons:
        bool veto = false;
        for (const Particle& zlep : zleps) {
          if (deltaR(jet, zlep) < 0.5) veto = true;
        }
        if (veto) continue;

        for (const Particle& bhadron : stableBs) {
          if (deltaR(jet, bhadron) <= 0.3) {
            b_jets.push_back(jet);
            break; // match
          }
	}
      }

      // Make sure we have at least 1
      if (b_jets.empty()) vetoEvent;

      // Fill the plots
      const double ZpT = boson_f.pT()/GeV;
      const double ZY  = boson_f.absrap();

      _h_bjet_ZPt->fill(ZpT);
      _h_bjet_ZY ->fill(ZY);

      for (const Jet& jet : b_jets) {
        _h_bjet_Pt->fill(jet.pT()/GeV);
        _h_bjet_Y ->fill(jet.absrap());

        const double Yboost = 0.5 * fabs(boson_f.rapidity() + jet.rapidity());

        _h_bjet_Yboost->fill(Yboost);

        if(ZpT > 20.) {

          const double ZBDY   = fabs( boson_f.rapidity() - jet.rapidity() );
          const double ZBDPHI = fabs( deltaPhi(jet.phi(), boson_f.phi()) );
          const double ZBDR   = deltaR(jet, boson_f, RAPIDITY);
          _h_bjet_DY20->fill(   ZBDY);
          _h_bjet_ZdPhi20->fill(ZBDPHI);
          _h_bjet_ZdR20->fill(  ZBDR);
        }

      } //loop over b-jets

      if (b_jets.size() < 2) return;

      _h_2bjet_ZPt->fill(ZpT);
      _h_2bjet_ZY ->fill(ZY);

      const double BBDR = deltaR(b_jets[0], b_jets[1], RAPIDITY);
      const double Mbb  = (b_jets[0].momentum() + b_jets[1].momentum()).mass();

      _h_2bjet_dR ->fill(BBDR);
      _h_2bjet_Mbb->fill(Mbb);

    } // end of analysis loop


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

      const double normfac = crossSection() / sumOfWeights();

      scale( _h_bjet_Pt,      normfac);
      scale( _h_bjet_Y,       normfac);
      scale( _h_bjet_Yboost,  normfac);
      scale( _h_bjet_DY20,    normfac);
      scale( _h_bjet_ZdPhi20, normfac);
      scale( _h_bjet_ZdR20,   normfac);
      scale( _h_bjet_ZPt,     normfac);
      scale( _h_bjet_ZY,      normfac);
      scale( _h_2bjet_dR,     normfac);
      scale( _h_2bjet_Mbb,    normfac);
      scale( _h_2bjet_ZPt,    normfac);
      scale( _h_2bjet_ZY,     normfac);
    }

    //@}


  protected:

    // Data members like post-cuts event weight counters go here
    size_t _mode;


  private:

    Histo1DPtr _h_bjet_Pt;
    Histo1DPtr _h_bjet_Y;
    Histo1DPtr _h_bjet_Yboost;
    Histo1DPtr _h_bjet_DY20;
    Histo1DPtr _h_bjet_ZdPhi20;
    Histo1DPtr _h_bjet_ZdR20;
    Histo1DPtr _h_bjet_ZPt;
    Histo1DPtr _h_bjet_ZY;
    Histo1DPtr _h_2bjet_dR;
    Histo1DPtr _h_2bjet_Mbb;
    Histo1DPtr _h_2bjet_ZPt;
    Histo1DPtr _h_2bjet_ZY;

  };


  RIVET_DECLARE_PLUGIN(ATLAS_2014_I1306294);

}