Rivet Analyses Reference

CDF_1996_S3108457

Properties of high-mass multijet events
Experiment: CDF (Tevatron Run 1)
Inspire ID: 393345
Status: VALIDATED
Authors:
  • Frank Siegert
References:Beams: p- p+
Beam energies: (900.0, 900.0) GeV
Run details:
  • Pure QCD events without underlying event (the paper states that UE was corrected for). Several runs with different kinematic cuts might be needed to fill the 2,3,4,5 and 6-jet properly.

Properties of two-, three-, four-, five-, and six-jet events... Multijet-mass, leading jet angle, jet pT.

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

namespace Rivet {


  /// @brief CDF properties of high-mass multi-jet events
  class CDF_1996_S3108457 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(CDF_1996_S3108457);


    /// @name Analysis methods
    //@{

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

      /// Initialise and register projections here
      const FinalState fs(Cuts::abseta < 4.2);
      FastJets fj(fs, FastJets::CDFJETCLU, 0.7);

      // Smear energy and mass with the 10% uncertainty quoted in the paper
      SmearedJets sj_E(fj, [](const Jet& jet){ return P4_SMEAR_MASS_GAUSS(P4_SMEAR_E_GAUSS(jet, 0.1*jet.E()), 0.1*jet.mass()); });
      declare(sj_E, "SmearedJets_E");


      /// Book histograms here, e.g.:
      for (size_t i=0; i<5; ++i) {
        book(_h_m[i], 1+i, 1, 1);
        book(_h_costheta[i], 10+i, 1, 1);
        book(_h_pT[i], 15+i, 1, 1);
      }
      /// @todo Ratios of mass histograms left out: Binning doesn't work out
    }


    /// Perform the per-event analysis
    void analyze(const Event& event) {
      // Get the smeared jets
      const Jets SJets = apply<JetAlg>(event, "SmearedJets_E").jets(Cuts::Et > 20.0*GeV, cmpMomByEt);
      if (SJets.size() < 2 || SJets.size() > 6) vetoEvent;

      // Calculate Et, total jet 4 Momentum
      double sumEt(0), sumE(0);
      FourMomentum JS(0,0,0,0);

      for (const Jet& jet : SJets) {
        sumEt += jet.Et()/GeV;
        sumE  += jet.E()/GeV;
        JS+=jet.momentum();
      }

      if (sumEt < 420. || sumE > 2000.) vetoEvent;

      double mass = JS.mass()/GeV;

      LorentzTransform cms_boost = LorentzTransform::mkFrameTransformFromBeta(JS.betaVec());
      FourMomentum jet0boosted(cms_boost.transform(SJets[0].momentum()));
      double costheta0 = fabs(cos(jet0boosted.theta()));

      if (costheta0 < 2.0/3.0) _h_m[SJets.size()-2]->fill(mass);
      if (mass > 600.) _h_costheta[SJets.size()-2]->fill(costheta0);
      if (costheta0 < 2.0/3.0 && mass > 600.) {
        for (const Jet& jet : SJets) _h_pT[SJets.size()-2]->fill(jet.pT());
      }
    }


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

      /// Normalise, scale and otherwise manipulate histograms here
      for (size_t i=0; i<5; ++i) {
        normalize(_h_m[i], 40.0);
        normalize(_h_costheta[i], 2.0);
        normalize(_h_pT[i], 20.0);
      }

    }

    //@}


  private:

    /// @name Histograms
    //@{
    Histo1DPtr _h_m[5];
    Histo1DPtr _h_costheta[5];
    Histo1DPtr _h_pT[5];
    //@}

  };



  RIVET_DECLARE_ALIASED_PLUGIN(CDF_1996_S3108457, CDF_1996_I393345);

}