Rivet Analyses Reference

JADE_OPAL_2000_S4300807

Jet rates in $e^+e^-$ at JADE [35--44 GeV] and OPAL [91--189 GeV].
Experiment: JADE_OPAL (PETRA and LEP)
Inspire ID: 513337
Status: VALIDATED
Authors:
  • Frank Siegert
  • Andy Buckley
References:Beams: e+ e-
Beam energies: (17.5, 17.5); (22.0, 22.0); (45.6, 45.6); (66.5, 66.5); (80.5, 80.5); (86.0, 86.0); (91.5, 91.5); (94.5, 94.5) GeV
Run details:
  • $e^+ e^- \to$ jet jet (+ jets)

Differential and integrated jet rates for Durham and JADE jet algorithms. The integration cut value used for the integrated rate observables is not well-defined in the paper: the midpoint of the differential bin has been used thanks to information from Stefan Kluth and Christoph Pahl. We anyway recommend that the differential plots be preferred over the integrated ones for MC generator validation and tuning, to minimise correlations.

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

namespace Rivet {


  /// @brief Jet rates in \f$ e^+ e^- \f$ at OPAL and JADE
  ///
  /// @author Frank Siegert
  class JADE_OPAL_2000_S4300807 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(JADE_OPAL_2000_S4300807)


    /// @name Analysis methods
    //@{

    void init() {
      // Projections
      const FinalState fs;
      declare(fs, "FS");
      FastJets jadeJets = FastJets(fs, FastJets::JADE, 0.7, JetAlg::Muons::ALL, JetAlg::Invisibles::DECAY);
      FastJets durhamJets = FastJets(fs, FastJets::DURHAM, 0.7, JetAlg::Muons::ALL, JetAlg::Invisibles::DECAY);
      declare(jadeJets, "JadeJets");
      declare(durhamJets, "DurhamJets");

      // Histos
      int offset = 0;
      switch (int(sqrtS()/GeV + 0.5)) {
      case 35: offset = 7; break;
      case 44: offset = 8; break;
      case 91: offset = 9; break;
      case 133: offset = 10; break;
      case 161: offset = 11; break;
      case 172: offset = 12; break;
      case 183: offset = 13; break;
      case 189: offset = 14; break;
      default: break;
      }
      for (size_t i = 0; i < 5; ++i) {
        book(_h_R_Jade[i] ,offset, 1, i+1);
        book(_h_R_Durham[i] ,offset+9, 1, i+1);
        if (i < 4) book(_h_y_Durham[i] ,offset+17, 1, i+1);
      }
    }



    void analyze(const Event& e) {
      MSG_DEBUG("Num particles = " << apply<FinalState>(e, "FS").particles().size());

      const FastJets& jadejet = apply<FastJets>(e, "JadeJets");
      if (jadejet.clusterSeq()) {
        const double y_23 = jadejet.clusterSeq()->exclusive_ymerge_max(2);
        const double y_34 = jadejet.clusterSeq()->exclusive_ymerge_max(3);
        const double y_45 = jadejet.clusterSeq()->exclusive_ymerge_max(4);
        const double y_56 = jadejet.clusterSeq()->exclusive_ymerge_max(5);

        for (size_t i = 0; i < _h_R_Jade[0]->numBins(); ++i) {
          double ycut = _h_R_Jade[0]->bin(i).xMid();
          double width = _h_R_Jade[0]->bin(i).xWidth();
          if (y_23 < ycut) {
            _h_R_Jade[0]->fillBin(i, width);
          }
        }
        for (size_t i = 0; i < _h_R_Jade[1]->numBins(); ++i) {
          double ycut = _h_R_Jade[1]->bin(i).xMid();
          double width = _h_R_Jade[1]->bin(i).xWidth();
          if (y_34 < ycut && y_23 > ycut) {
            _h_R_Jade[1]->fillBin(i, width);
          }
        }
        for (size_t i = 0; i < _h_R_Jade[2]->numBins(); ++i) {
          double ycut = _h_R_Jade[2]->bin(i).xMid();
          double width = _h_R_Jade[2]->bin(i).xWidth();
          if (y_45 < ycut && y_34 > ycut) {
            _h_R_Jade[2]->fillBin(i, width);
          }
        }
        for (size_t i = 0; i < _h_R_Jade[3]->numBins(); ++i) {
          double ycut = _h_R_Jade[3]->bin(i).xMid();
          double width = _h_R_Jade[3]->bin(i).xWidth();
          if (y_56 < ycut && y_45 > ycut) {
            _h_R_Jade[3]->fillBin(i, width);
          }
        }
        for (size_t i = 0; i < _h_R_Jade[4]->numBins(); ++i) {
          double ycut = _h_R_Jade[4]->bin(i).xMid();
          double width = _h_R_Jade[4]->bin(i).xWidth();
          if (y_56 > ycut) {
            _h_R_Jade[4]->fillBin(i, width);
          }
        }
      }

      const FastJets& durjet = apply<FastJets>(e, "DurhamJets");
      if (durjet.clusterSeq()) {
        const double y_23 = durjet.clusterSeq()->exclusive_ymerge_max(2);
        const double y_34 = durjet.clusterSeq()->exclusive_ymerge_max(3);
        const double y_45 = durjet.clusterSeq()->exclusive_ymerge_max(4);
        const double y_56 = durjet.clusterSeq()->exclusive_ymerge_max(5);

        _h_y_Durham[0]->fill(y_23);
        _h_y_Durham[1]->fill(y_34);
        _h_y_Durham[2]->fill(y_45);
        _h_y_Durham[3]->fill(y_56);

        for (size_t i = 0; i < _h_R_Durham[0]->numBins(); ++i) {
          double ycut = _h_R_Durham[0]->bin(i).xMid();
          double width = _h_R_Durham[0]->bin(i).xWidth();
          if (y_23 < ycut) {
            _h_R_Durham[0]->fillBin(i, width);
          }
        }
        for (size_t i = 0; i < _h_R_Durham[1]->numBins(); ++i) {
          double ycut = _h_R_Durham[1]->bin(i).xMid();
          double width = _h_R_Durham[1]->bin(i).xWidth();
          if (y_34 < ycut && y_23 > ycut) {
            _h_R_Durham[1]->fillBin(i, width);
          }
        }
        for (size_t i = 0; i < _h_R_Durham[2]->numBins(); ++i) {
          double ycut = _h_R_Durham[2]->bin(i).xMid();
          double width = _h_R_Durham[2]->bin(i).xWidth();
          if (y_45 < ycut && y_34 > ycut) {
            _h_R_Durham[2]->fillBin(i, width);
          }
        }
        for (size_t i = 0; i < _h_R_Durham[3]->numBins(); ++i) {
          double ycut = _h_R_Durham[3]->bin(i).xMid();
          double width = _h_R_Durham[3]->bin(i).xWidth();
          if (y_56 < ycut && y_45 > ycut) {
            _h_R_Durham[3]->fillBin(i, width);
          }
        }
        for (size_t i = 0; i < _h_R_Durham[4]->numBins(); ++i) {
          double ycut = _h_R_Durham[4]->bin(i).xMid();
          double width = _h_R_Durham[4]->bin(i).xWidth();
          if (y_56 > ycut) {
            _h_R_Durham[4]->fillBin(i, width);
          }
        }
      }
    }



    /// Finalize
    void finalize() {
      for (size_t n = 0; n < 4; ++n) normalize(_h_y_Durham[n]);
      for (size_t n = 0; n < 5; ++n) scale(_h_R_Jade[n], 100/sumOfWeights());
      for (size_t n = 0; n < 5; ++n) scale(_h_R_Durham[n], 100/sumOfWeights());
    }

    //@}


  private:

    /// @name Histograms
    //@{
    Histo1DPtr _h_R_Jade[5];
    Histo1DPtr _h_R_Durham[5];
    Histo1DPtr _h_y_Durham[4];
    //@}

  };



  RIVET_DECLARE_ALIASED_PLUGIN(JADE_OPAL_2000_S4300807, JADE_OPAL_2000_I513337);

}