rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

MC_HJETS

Monte Carlo validation observables for $h[\tau^+ \, \tau^-]$ + jets production
Experiment: ()
Status: VALIDATED
Authors:
  • Frank Siegert
No references listed
Beams: * *
Beam energies: ANY
Run details:
  • $h [\to \tau^+ \tau^-]$ + jets.

The available observables are the Higgs mass, 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_HJETS.cc
 1// -*- C++ -*-
 2#include "Rivet/Analyses/MC_JetAnalysis.hh"
 3#include "Rivet/Projections/ZFinder.hh"
 4#include "Rivet/Projections/FastJets.hh"
 5
 6namespace Rivet {
 7
 8
 9  /// @brief MC validation analysis for Higgs [-> tau tau] + jets events
10  class MC_HJETS : public MC_JetAnalysis {
11  public:
12
13    /// Default constructor
14    MC_HJETS()
15      : MC_JetAnalysis("MC_HJETS", 4, "Jets")
16    {    }
17
18
19    /// @name Analysis methods
20    //@{
21
22    /// Book histograms
23    void init() {
24      Cut cut = Cuts::abseta < 3.5 && Cuts::pT > 25*GeV;
25      /// @todo Urk, abuse! Need explicit HiggsFinder (and TauFinder?)
26      ZFinder hfinder(FinalState(), cut, PID::TAU, 115*GeV, 135*GeV, 0.0, ZFinder::ClusterPhotons::NONE, ZFinder::AddPhotons::NO, 125*GeV);
27      declare(hfinder, "Hfinder");
28      FastJets jetpro(hfinder.remainingFinalState(), FastJets::ANTIKT, 0.4);
29      declare(jetpro, "Jets");
30
31      book(_h_H_jet1_deta ,"H_jet1_deta", 50, -5.0, 5.0);
32      book(_h_H_jet1_dR ,"H_jet1_dR", 25, 0.5, 7.0);
33
34      MC_JetAnalysis::init();
35    }
36
37
38
39    /// Do the analysis
40    void analyze(const Event & e) {
41      const ZFinder& hfinder = apply<ZFinder>(e, "Hfinder");
42      if (hfinder.bosons().size() != 1) vetoEvent;
43      const double weight = 1.0;
44
45      FourMomentum hmom(hfinder.bosons()[0].momentum());
46      const Jets& jets = apply<FastJets>(e, "Jets").jetsByPt(_jetptcut);
47      if (jets.size() > 0) {
48        _h_H_jet1_deta->fill(hmom.eta()-jets[0].eta(), weight);
49        _h_H_jet1_dR->fill(deltaR(hmom, jets[0].momentum()), weight);
50      }
51
52      MC_JetAnalysis::analyze(e);
53    }
54
55
56    /// Finalize
57    void finalize() {
58      normalize(_h_H_jet1_deta, crossSection()/picobarn);
59      normalize(_h_H_jet1_dR, crossSection()/picobarn);
60      MC_JetAnalysis::finalize();
61    }
62
63    //@}
64
65
66  private:
67
68    /// @name Histograms
69    //@{
70    Histo1DPtr _h_H_jet1_deta;
71    Histo1DPtr _h_H_jet1_dR;
72    //@}
73
74  };
75
76
77
78  // The hook for the plugin system
79  RIVET_DECLARE_PLUGIN(MC_HJETS);
80
81}