rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

MC_HHJETS

Monte Carlo validation observables for $HH$ (stable) + jets production
Experiment: ()
Status: VALIDATED
Authors:
  • Andreas Papaefstathiou
No references listed
Beams: * *
Beam energies: ANY
Run details:
  • $HH$ production with stable Higgses

The available observables are the Higgs boson pair invariant mass, the separation between the two Higgs bosons, pT of the di-Higgs (any, hardest, second hardest), pT and pseudorapidity of Higgs bosons (any, hardest, second hardest), pT of jets 1--4, jet multiplicity, $\Delta\eta(h, \text{jet1})$, $\Delta R(\text{jet2}, \text{jet3})$, differential jet rates 0->1, 1->2, 2->3, 3->4, and integrated 0--4 jet rates.

Source code: MC_HHJETS.cc
  1// -*- C++ -*-
  2#include "Rivet/Analyses/MC_JetAnalysis.hh"
  3#include "Rivet/Projections/ZFinder.hh"
  4#include "Rivet/Projections/VetoedFinalState.hh"
  5#include "Rivet/Projections/IdentifiedFinalState.hh"
  6#include "Rivet/Projections/FastJets.hh"
  7
  8namespace Rivet {
  9
 10  /// @brief MC validation analysis for higgs pairs events (stable Higgses)
 11  class MC_HHJETS : public MC_JetAnalysis {
 12  public:
 13
 14    /// Default constructor
 15    MC_HHJETS()
 16      : MC_JetAnalysis("MC_HHJETS", 4, "Jets")
 17    {    }
 18
 19
 20    /// @name Analysis methods
 21    //@{
 22
 23    /// Book histograms
 24    void init() {
 25      IdentifiedFinalState ifs(Cuts::abseta < 10.0 && Cuts::pT > 0*GeV);
 26      ifs.acceptId(25);
 27      declare(ifs,"IFS");
 28
 29      VetoedFinalState vfs;
 30      vfs.addVetoPairId(25);
 31      declare(FastJets(vfs, FastJets::ANTIKT, 0.4), "Jets");
 32
 33      book(_h_HH_mass ,"HH_mass", 250, 240, 4000.0);
 34      book(_h_HH_dR ,"HH_dR", 25, 0.5, 10.0);
 35      book(_h_HH_dPhi ,"HH_dPhi", 64, 0, 3.2);
 36      book(_h_HH_deta,"HH_deta", 50, -5, 5);
 37      book(_h_H_pT ,"H_pT", 50, 0, 2000.0);
 38      book(_h_HH_pT ,"HH_pT", 200, 0, 2000.0);
 39      book(_h_H_pT1 ,"H_pT1", 200, 0, 2000.0);
 40      book(_h_H_pT2 ,"H_pT2", 200, 0, 2000.0);
 41      book(_h_H_eta ,"H_eta", 50, -5.0, 5.0);
 42      book(_h_H_eta1 ,"H_eta1", 50, -5.0, 5.0);
 43      book(_h_H_eta2 ,"H_eta2", 50, -5.0, 5.0);
 44      book(_h_H_phi ,"H_phi", 25, 0.0, TWOPI);
 45      book(_h_H_jet1_deta ,"H_jet1_deta", 50, -5.0, 5.0);
 46      book(_h_H_jet1_dR ,"H_jet1_dR", 25, 0.5, 7.0);
 47
 48      MC_JetAnalysis::init();
 49    }
 50
 51
 52
 53    /// Do the analysis
 54    void analyze(const Event & e) {
 55
 56      const IdentifiedFinalState& ifs = apply<IdentifiedFinalState>(e, "IFS");
 57      Particles allp = ifs.particlesByPt();
 58      if (allp.empty()) vetoEvent;
 59
 60      const double weight = 1.0;
 61
 62      FourMomentum hmom = allp[0].momentum();
 63      if (allp.size() > 1) {
 64        FourMomentum hmom2(allp[1].momentum());
 65        _h_HH_dR->fill(deltaR(hmom, hmom2), weight);
 66        _h_HH_dPhi->fill(deltaPhi(hmom, hmom2), weight);
 67        _h_HH_deta->fill(hmom.eta()-hmom2.eta(), weight);
 68        _h_HH_pT->fill((hmom+hmom2).pT(), weight);
 69        _h_HH_mass->fill((hmom+hmom2).mass(), weight);
 70
 71        if (hmom.pT() > hmom2.pT()) {
 72          _h_H_pT1->fill(hmom.pT(), weight);
 73          _h_H_eta1->fill(hmom.eta(), weight);
 74          _h_H_pT2->fill(hmom2.pT(), weight);
 75          _h_H_eta2->fill(hmom2.eta(), weight);
 76        } else {
 77          _h_H_pT1->fill(hmom2.pT(), weight);
 78          _h_H_eta1->fill(hmom2.eta(), weight);
 79          _h_H_pT2->fill(hmom.pT(), weight);
 80          _h_H_eta2->fill(hmom.eta(), weight);
 81        }
 82      }
 83      _h_H_pT->fill(hmom.pT(), weight);
 84      _h_H_eta->fill(hmom.eta(), weight);
 85      _h_H_phi->fill(hmom.azimuthalAngle(), weight);
 86
 87
 88      // Get the jet candidates
 89      Jets jets = apply<FastJets>(e, "Jets").jetsByPt(20.0*GeV);
 90      if (!jets.empty()) {
 91        _h_H_jet1_deta->fill(deltaEta(hmom, jets[0]), weight);
 92        _h_H_jet1_dR->fill(deltaR(hmom, jets[0]), weight);
 93      }
 94
 95      MC_JetAnalysis::analyze(e);
 96    }
 97
 98
 99    /// Finalize
100    void finalize() {
101      scale(_h_HH_mass, crossSection()/sumOfWeights());
102      scale(_h_HH_dR, crossSection()/sumOfWeights());
103      scale(_h_HH_deta, crossSection()/sumOfWeights());
104      scale(_h_HH_dPhi, crossSection()/sumOfWeights());
105      scale(_h_H_pT, crossSection()/sumOfWeights());
106      scale(_h_H_pT1, crossSection()/sumOfWeights());
107      scale(_h_H_pT2, crossSection()/sumOfWeights());
108      scale(_h_HH_pT, crossSection()/sumOfWeights());
109      scale(_h_H_eta, crossSection()/sumOfWeights());
110      scale(_h_H_eta1, crossSection()/sumOfWeights());
111      scale(_h_H_eta2, crossSection()/sumOfWeights());
112      scale(_h_H_phi, crossSection()/sumOfWeights());
113      scale(_h_H_jet1_deta, crossSection()/sumOfWeights());
114      scale(_h_H_jet1_dR, crossSection()/sumOfWeights());
115
116      MC_JetAnalysis::finalize();
117    }
118
119    //@}
120
121
122  private:
123
124    /// @name Histograms
125    //@{
126    Histo1DPtr _h_HH_mass, _h_HH_pT, _h_HH_dR, _h_HH_deta, _h_HH_dPhi;
127    Histo1DPtr _h_H_pT, _h_H_pT1, _h_H_pT2, _h_H_eta, _h_H_eta1, _h_H_eta2, _h_H_phi;
128    Histo1DPtr _h_H_jet1_deta, _h_H_jet1_dR;
129    //@}
130
131  };
132
133
134  // The hook for the plugin system
135  RIVET_DECLARE_PLUGIN(MC_HHJETS);
136
137}