Rivet Analyses Reference

ATLAS_2015_I1364361

Total and differential Higgs cross sections at 8 TeV with ATLAS
Experiment: ATLAS (LHC)
Inspire ID: 1364361
Status: VALIDATED
Authors:
  • Michaela Queitsch-Maitland
References:Beams: p+ p+
Beam energies: (4000.0, 4000.0) GeV
Run details:
  • pp Higgs production at 8 TeV, include all processes (ggH, VH, VBF, ttH), stable Higgs

WARNING --- this analysis requires stable higgs bosons --- Measurements of the total and differential cross sections of Higgs boson production are performed using 20.3 $\text{fb}^{-1}$ of pp collisions produced by the Large Hadron Collider at a center-of-mass energy of 8 TeV and recorded by the ATLAS detector. Cross sections are obtained from measured Higgs to diphoton and ZZ* event yields, which are combined accounting for detector efficiencies, fiducial acceptances and branching fractions. Differential cross sections are reported as a function of Higgs boson transverse momentum, Higgs boson rapidity, number of jets in the event, and transverse momentum of the leading jet. The total production cross section is determined to be $33.0 \pm 5.3 (\text{stat}) \pm 1.6 (\text{sys})$ pb.

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

namespace Rivet {


  /// @brief Higgs differential cross-section combination between the ATLAS measurements in the yy and 4l channels
  ///
  /// Computes Higgs transverse momentum, rapidity, jet multiplicity and leading jet pT.
  ///
  /// @author Michaela Queitsch-Maitland <michaela.queitsch-maitland@cern.ch>
  /// @author Dag Gillberg <dag.gillberg@cern.ch>
  /// @author Florian Bernlochner <florian.bernlochner@cern.ch>
  /// @author Sarah Heim <sarah.heim@cern.ch>
  class ATLAS_2015_I1364361 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2015_I1364361);


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

      // All final state particles
      const FinalState fs;
      declare(fs, "FS");

      // Histograms with data bins
      book(_h_pTH_incl   ,1,1,1);
      book(_h_yH_incl    ,2,1,1);
      book(_h_Njets_incl ,3,1,1);
      book(_h_pTj1_incl  ,4,1,1);
    }


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

      // Get the final state particles ordered by pT
      const Particles& fs = apply<FinalState>(event, "FS").particlesByPt();

      // Find a stable Higgs (mandatory)
      const auto higgsiter = std::find_if(fs.begin(), fs.end(), [](const Particle& p){ return p.pid() == PID::HIGGSBOSON; });
      if (higgsiter == fs.end()) vetoEvent;
      const Particle& higgs = *higgsiter;

      // bool stable_higgs = false;
      // const Particle* higgs = 0;
      // for (const Particle& p : fs) {
      //   if (p.pid() == PID::HIGGSBOSON) {
      //     stable_higgs = true;
      //     higgs = &p;
      //     break;
      //   }
      // }
      // if (!stable_higgs) {
      //   MSG_WARNING("FATAL: No stable Higgs found in event record.\n");
      //   vetoEvent;
      // }


      // Loop over final state particles and fill various particle vectors
      Particles leptons, photons, jet_ptcls;
      for ( const Particle& ptcl : fs ) {
        // Do not include the Higgs in jet finding!
        if ( ptcl.pid() == PID::HIGGSBOSON ) continue;
        // Neutrinos not from hadronisation
        if ( ptcl.isNeutrino() && !ptcl.fromHadron() ) continue;
        // Electrons and muons not from hadronisation
        if ( ( ptcl.abspid() == PID::ELECTRON || ptcl.abspid() == PID::MUON ) && !ptcl.fromHadron() ) {
          leptons.push_back(ptcl);
          continue;
        }
        // Photons not from hadronisation
        if ( ptcl.abspid() == PID::PHOTON && !ptcl.fromHadron() ) {
          photons.push_back(ptcl);
          continue;
        }
        // Add particle to jet inputs
        jet_ptcls.push_back(ptcl);
      }

      // Match FS photons to leptons within cone R=0.1
      // If they are not 'dressing' photons, add to jet particle vector
      for ( const Particle& ph : photons ) {
        bool fsr_photon = false;
        for ( const Particle& lep : leptons ) {
          if ( deltaR(ph.momentum(),lep.momentum()) < 0.1 ){
            fsr_photon=true;
            continue;
          }
        }
        if ( !fsr_photon ) jet_ptcls.push_back(ph);
      }

      // Let's build the jets! By hand...
      const PseudoJets pjs_in = mkPseudoJets(jet_ptcls);
      const fastjet::JetDefinition jdef(fastjet::antikt_algorithm, 0.4);
      const Jets alljets = mkJets(fastjet::ClusterSequence(pjs_in, jdef).inclusive_jets());
      const Jets jets = sortByPt(filterBy(alljets, Cuts::pT > 30*GeV && Cuts::absrap < 4.4));
      // FastJets jet_pro(FastJets::ANTIKT, 0.4);
      // jet_pro.calc(jet_ptcls);
      // Jets jets = jet_pro.jetsByPt(Cuts::pT > 30*GeV && Cuts::absrap < 4.4);

      size_t njets = jets.size() > 3 ? 3 : jets.size();
      _h_pTH_incl->fill(higgs.pT());
      _h_yH_incl->fill(higgs.absrap());
      _h_Njets_incl->fill(njets + 1); // accounts for HEPData offset
      _h_pTj1_incl->fill(jets.empty() ? 0 : jets[0].pT());
    }


    /// Normalise histograms etc., after the run
    void finalize() {
      const double xs = crossSectionPerEvent();
      scale(_h_pTH_incl, xs);
      scale(_h_yH_incl, xs);
      scale(_h_Njets_incl, xs);
      scale(_h_pTj1_incl, xs);
    }


  private:

    Histo1DPtr _h_pTH_incl, _h_yH_incl, _h_Njets_incl, _h_pTj1_incl;

  };


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


}