Rivet Analyses Reference

ATLAS_2014_I1310835

H(125) -> 4l at 8 TeV
Experiment: ATLAS (LHC)
Inspire ID: 1310835
Status: VALIDATED
Authors:
  • Jonathan Stahlman
  • Christian Gutschow
References:Beams: p+ p+
Beam energies: (4000.0, 4000.0) GeV
Run details:
  • p + p -> H (-> 4 l) + X at 8 TeV

Measurements of fiducial and differential cross sections of Higgs boson production in the $H\to ZZ^\ast\to 4\ell$ decay channel are presented. The cross sections are determined within a fiducial phase space and corrected for detection efficiency and resolution effects. They are based on 20.3 fb$^{-1}$ of $pp$ collision data, produced at $\sqrt{s}=8$ TeV centre-of-mass energy at the LHC and recorded by the ATLAS detector. The differential measurements are performed in bins of transverse momentum and rapidity of the four-lepton system, the invariant mass of the subleading lepton pair and the decay angle of the leading lepton pair with respect to the beam line in the four-lepton rest frame, as well as the number of jets and the transverse momentum of the leading jet. The measured cross sections are compared to selected theoretical calculations of the Standard Model expectations. No significant deviation from any of the tested predictions is found.

Source code: ATLAS_2014_I1310835.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
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/PromptFinalState.hh"
#include "Rivet/Projections/DressedLeptons.hh"

namespace Rivet {

  /// @brief H(125)->ZZ->4l at 8 TeV
  class ATLAS_2014_I1310835 : public Analysis {
  public:

    /// Default constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2014_I1310835);

    void init() {
      const FinalState fs(Cuts::abseta < 5.0);

      PromptFinalState photons(Cuts::abspid == PID::PHOTON);

      PromptFinalState bare_el(Cuts::abspid == PID::ELECTRON);

      PromptFinalState bare_mu(Cuts::abspid == PID::MUON);

      // Selection: lepton selection
      Cut etaranges_el = Cuts::abseta < 2.47 && Cuts::pT > 7*GeV; 
      DressedLeptons electron_sel4l(photons, bare_el, 0.1, etaranges_el, false);
      declare(electron_sel4l, "electrons");
 
      Cut etaranges_mu = Cuts::abseta < 2.7 && Cuts::pT > 6*GeV;
      DressedLeptons muon_sel4l(photons, bare_mu, 0.1, etaranges_mu, false);
      declare(muon_sel4l, "muons");

      FastJets jetpro(fs, FastJets::ANTIKT, 0.4, JetAlg::Muons::NONE, JetAlg::Invisibles::NONE);
      declare(jetpro, "jet");

      // Book histos
      book(_h_pt          , 1, 1, 1);
      book(_h_rapidity    , 2, 1, 1);
      book(_h_m34         , 3, 1, 1);
      book(_h_costheta    , 4, 1, 1);
      book(_h_njets       , 5, 1, 1);
      book(_h_leadingjetpt, 6, 1, 1);

    }



    /// Do the analysis
    void analyze(const Event& e) {
      
      ////////////////////////////////////////////////////////////////////
      // preselection of leptons for ZZ-> llll final state
      ////////////////////////////////////////////////////////////////////

      const vector<DressedLepton>& mu_sel4l = applyProjection<DressedLeptons>(e, "muons").dressedLeptons();
      const vector<DressedLepton>& el_sel4l = applyProjection<DressedLeptons>(e, "electrons").dressedLeptons();

      vector<DressedLepton> leptonsFS_sel4l;
      leptonsFS_sel4l.insert( leptonsFS_sel4l.end(), mu_sel4l.begin(), mu_sel4l.end() );
      leptonsFS_sel4l.insert( leptonsFS_sel4l.end(), el_sel4l.begin(), el_sel4l.end() );

      /////////////////////////////////////////////////////////////////////////////
      /// H->ZZ->4l pairing
      /////////////////////////////////////////////////////////////////////////////
 
      size_t el_p = 0;
      size_t el_n = 0;
      size_t mu_p = 0; 
      size_t mu_n = 0;
      
      for (const Particle& l : leptonsFS_sel4l) {
        if (l.abspid() == PID::ELECTRON) {
          if (l.pid() < 0)  ++el_n;
          if (l.pid() > 0)  ++el_p;
        }
        else if (l.abspid() == PID::MUON) {
          if (l.pid() < 0)  ++mu_n;
          if (l.pid() > 0)  ++mu_p;
        }
      }
            
      bool pass_sfos = ( (el_p >=2 && el_n >=2) || (mu_p >=2 && mu_n >=2) || (el_p >=1 && el_n >=1 && mu_p >=1 && mu_n >=1) );
      
      if (!pass_sfos)  vetoEvent;

      Zstate Z1, Z2, Zcand;
      size_t n_parts = leptonsFS_sel4l.size();
      size_t l1_index = 0;
      size_t l2_index = 0;

      // determine Z1 first
      double min_mass_diff = -1;
      for (size_t i = 0; i < n_parts; ++i) {
        for (size_t j = 0; j < n_parts; ++j) {
          if (i >= j)  continue;

          if (leptonsFS_sel4l[i].pid() != -1*leptonsFS_sel4l[j].pid())  continue; //only pair SFOS leptons

          Zcand = Zstate( ParticlePair(leptonsFS_sel4l[i], leptonsFS_sel4l[j]) );
          double mass_diff = fabs( Zcand.mom().mass() - 91.1876 );
         
          if (min_mass_diff == -1 || mass_diff < min_mass_diff) {
            min_mass_diff = mass_diff;
            Z1 = Zcand;
            l1_index = i;
            l2_index = j;
          }
        }
      }

      //determine Z2 second
      min_mass_diff = -1;
      for (size_t i = 0; i < n_parts; ++i) {
        if (i == l1_index || i == l2_index)  continue;
        for (size_t j = 0; j < n_parts; ++j) {
          if (j == l1_index || j == l2_index || i >= j)  continue;

          if (leptonsFS_sel4l[i].pid() != -1*leptonsFS_sel4l[j].pid())  continue; // only pair SFOS leptons

          Zcand = Zstate( ParticlePair(leptonsFS_sel4l[i], leptonsFS_sel4l[j]) );
          double mass_diff = fabs( Zcand.mom().mass() - 91.1876 );

          if (min_mass_diff == -1 || mass_diff < min_mass_diff) {
            min_mass_diff = mass_diff;
            Z2 = Zcand;
          }
        }
      }

      Particles leptons_sel4l;
      leptons_sel4l.push_back(Z1.first);
      leptons_sel4l.push_back(Z1.second);
      leptons_sel4l.push_back(Z2.first);
      leptons_sel4l.push_back(Z2.second);

      ////////////////////////////////////////////////////////////////////////////
      // Kinematic Requirements
      ///////////////////////////////////////////////////////////////////////////
      
      //leading lepton pT requirement
      std::vector<double> lepton_pt;
      for (const Particle& i : leptons_sel4l) lepton_pt.push_back(i.pT() / GeV);
      std::sort(lepton_pt.begin(), lepton_pt.end(), [](const double pT1, const double pT2) { return pT1 > pT2; });
      
      if (!(lepton_pt[0] > 20*GeV && lepton_pt[1] > 15*GeV && lepton_pt[2] > 10*GeV))  vetoEvent;
      
      //invariant mass requirements
      if (!(inRange(Z1.mom().mass(), 50*GeV, 106*GeV) && inRange(Z2.mom().mass(), 12*GeV, 115*GeV)))  vetoEvent;
      
      //lepton separation requirements
      for (unsigned int i = 0; i < 4; ++i) {
        for (unsigned int j = 0; j < 4; ++j) {
          if (i >= j) continue;
          double dR = deltaR(leptons_sel4l[i], leptons_sel4l[j]);
          bool sameflavor = leptons_sel4l[i].abspid() == leptons_sel4l[j].abspid();

          if ( sameflavor && dR < 0.1)  vetoEvent;
          if (!sameflavor && dR < 0.2)  vetoEvent;
        }
      }

      // J/Psi veto requirement
      for (unsigned int i = 0; i < 4; ++i) {
        for (unsigned int j = 0; j < 4; ++j) {
          if (i >= j) continue;
          if ( leptons_sel4l[i].pid() != -1*leptons_sel4l[j].pid() )  continue;
          if ((leptons_sel4l[i].momentum() + leptons_sel4l[j].momentum()).mass() <= 5*GeV)  vetoEvent;
        }
      }
 
      // 4-lepton invariant mass requirement
      double m4l = (Z1.mom() + Z2.mom()).mass();
      if (!(inRange(m4l, 118*GeV, 129*GeV)))  vetoEvent;
  
  
      ////////////////////////////////////////////////////////////////////////////
      // Higgs observables
      ///////////////////////////////////////////////////////////////////////////
      FourMomentum Higgs = Z1.mom() + Z2.mom();

      double H4l_pt       = Higgs.pt()/GeV; 
      double H4l_rapidity = Higgs.absrap(); 
      LorentzTransform HRF_boost;
      //HRF_boost.mkFrameTransformFromBeta(Higgs.betaVec());
      HRF_boost.setBetaVec(- Higgs.betaVec());
      FourMomentum Z1_in_HRF = HRF_boost.transform( Z1.mom() );
      double H4l_costheta = fabs(cos( Z1_in_HRF.theta())); 
      double H4l_m34      = Z2.mom().mass()/GeV;
      
      ////////////////////////////////////////////////////////////////////////////
      // Jet observables
      ///////////////////////////////////////////////////////////////////////////

      Jets jets;
      for (const Jet& jet : applyProjection<FastJets>(e, "jet").jetsByPt(Cuts::pT > 30*GeV && Cuts::absrap < 4.4)) {
        bool overlaps = false;
        for (const Particle& lep : leptonsFS_sel4l) {
          if (lep.abspid() != PID::ELECTRON)  continue;
          const double dR = deltaR(lep, jet);
          if (dR < 0.2) { overlaps = true; break; }
        }
        if (!overlaps) jets += jet;
      }
      size_t n_jets = jets.size();
      if (n_jets > 3)  n_jets = 3;

      std::vector<double> jet_pt;
      for (const Jet& i : jets) jet_pt.push_back(i.pT()/GeV);

      double leading_jet_pt = n_jets? jet_pt[0] : 0.;

      ////////////////////////////////////////////////////////////////////////////
      // End of H->ZZ->llll selection: now fill histograms
      ////////////////////////////////////////////////////////////////////////////


      _h_pt->fill(H4l_pt);
      _h_rapidity->fill(H4l_rapidity);
      _h_costheta->fill(H4l_costheta);
      _h_m34->fill(H4l_m34);
      _h_njets->fill(n_jets + 1);
      _h_leadingjetpt->fill(leading_jet_pt);


    }


    /// Generic Z candidate
    struct Zstate : public ParticlePair {
      Zstate() { }
      Zstate(ParticlePair _particlepair) : ParticlePair(_particlepair) { }
      FourMomentum mom() const { return first.momentum() + second.momentum(); }
      operator FourMomentum() const { return mom(); }
    };

    /// Finalize
    void finalize() {

      const double norm = crossSection()/sumOfWeights()/femtobarn;
      scale(_h_pt, norm);
      scale(_h_rapidity, norm);
      scale(_h_costheta, norm);
      scale(_h_m34, norm);
      scale(_h_njets, norm);
      scale(_h_leadingjetpt, norm);
    }


  private:

    Histo1DPtr _h_pt, _h_rapidity, _h_costheta;
    Histo1DPtr _h_m34, _h_njets, _h_leadingjetpt;

  };


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

}