rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

L3_1992_I334954

Event Shapes at LEPI
Experiment: L3 (LEP)
Inspire ID: 334954
Status: VALIDATED
Authors:
  • Peter Richardson
References:
  • Z.Phys.C 55 (1992) 39-62, 1992
Beams: e+ e-
Beam energies: (45.6, 45.6) GeV
Run details:
  • E+e- to hadrons

Measurement of a range of event shapes by L3 at LEP1. One analysis specific observable and Fox-Wolfram moments not implemented.

Source code: L3_1992_I334954.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FastJets.hh"
  4#include "Rivet/Projections/FinalState.hh"
  5#include "Rivet/Projections/ChargedFinalState.hh"
  6#include "Rivet/Projections/Sphericity.hh"
  7#include "Rivet/Projections/Thrust.hh"
  8#include "Rivet/Projections/Hemispheres.hh"
  9#include "Rivet/Projections/ParisiTensor.hh"
 10
 11namespace Rivet {
 12
 13
 14  /// @brief Event shapes at MZ
 15  class L3_1992_I334954 : public Analysis {
 16  public:
 17
 18    /// Constructor
 19    RIVET_DEFAULT_ANALYSIS_CTOR(L3_1992_I334954);
 20
 21
 22    /// @name Analysis methods
 23    ///@{
 24
 25    /// Book histograms and initialise projections before the run
 26    void init() {
 27      // projections
 28      const FinalState fs;
 29      declare(fs, "FS");
 30      const ChargedFinalState cfs;
 31      declare(cfs, "CFS");
 32      declare(Sphericity(fs), "Sphericity");
 33      Thrust thrust(fs);
 34      declare(thrust, "Thrust"    );
 35      declare(Hemispheres(thrust), "Hemispheres");
 36      declare(ParisiTensor(fs), "Parisi");
 37      declare(FastJets(cfs, JetAlg::DURHAM, 0.7), "DurhamJets");
 38      declare(FastJets(cfs, JetAlg::JADE  , 0.7), "JadeJets"  );
 39      // histograms
 40      book(_histThrust    ,  1, 1, 1);
 41      book(_histMajor     ,  2, 1, 1);
 42      book(_histMinor     ,  3, 1, 1);
 43      book(_histOblateness,  4, 1, 1);
 44      book(_histJade      ,  6, 1, 1);
 45      book(_histDurham    ,  7, 1, 1);
 46      book(_histSphericity, 10, 1, 1);
 47      book(_histAplanarity, 11, 1, 1);
 48      book(_histC         , 12, 1, 1);
 49      book(_histD         , 13, 1, 1);
 50      book(_histMJetHeavy , 14, 1, 1);
 51      book(_histMJetLight , 15, 1, 1);
 52      book(_histMult      , 16, 1, 1);
 53    }
 54
 55
 56    /// Perform the per-event analysis
 57    void analyze(const Event& event) {
 58      // 5 charged particles
 59      const FinalState& cfs = apply<FinalState>(event, "CFS");
 60      if(cfs.particles().size()<5) vetoEvent;
 61      // charged particle mult
 62      const string multi_edge = std::to_string(cfs.particles().size()) + ".0";
 63      _histMult->fill(multi_edge);
 64      // Sphericity related
 65      const Sphericity& sphericity = apply<Sphericity>(event, "Sphericity");
 66      _histSphericity->fill(sphericity.sphericity());
 67      _histAplanarity->fill(sphericity.aplanarity());
 68      // thrust related
 69      const Thrust& thrust = apply<Thrust>(event, "Thrust");
 70      _histThrust    ->fill(thrust.thrust());
 71      _histMajor     ->fill(thrust.thrustMajor());
 72      _histMinor     ->fill(thrust.thrustMinor());
 73      _histOblateness->fill(thrust.oblateness());
 74      // hemisphere related
 75      const Hemispheres& hemi = apply<Hemispheres>(event, "Hemispheres");
 76      _histMJetHeavy->fill(hemi.scaledM2high());
 77      _histMJetLight->fill(hemi.scaledM2low());
 78      // C and D
 79      const ParisiTensor& parisi = apply<ParisiTensor>(event, "Parisi");
 80      _histC->fill(parisi.C());
 81      _histD->fill(parisi.D());
 82      // jet rates
 83      const FastJets& durjet = apply<FastJets>(event, "DurhamJets");
 84      double y23 = durjet.clusterSeq()->exclusive_ymerge_max(2);
 85      _histDurham->fill(y23);
 86      const FastJets& jadejet = apply<FastJets>(event, "JadeJets");
 87      y23 = jadejet.clusterSeq()->exclusive_ymerge_max(2);
 88      _histJade->fill(y23);
 89    }
 90
 91
 92    /// Normalise histograms etc., after the run
 93    void finalize() {
 94      scale(_histThrust    , 1./sumOfWeights());
 95      scale(_histMajor     , 1./sumOfWeights());
 96      scale(_histMinor     , 1./sumOfWeights());
 97      scale(_histOblateness, 1./sumOfWeights());
 98      scale(_histJade      , 1./sumOfWeights());
 99      scale(_histDurham    , 1./sumOfWeights());
100      scale(_histSphericity, 1./sumOfWeights());
101      scale(_histAplanarity, 1./sumOfWeights());
102      scale(_histC         , 1./sumOfWeights());
103      scale(_histD         , 1./sumOfWeights());
104      scale(_histMJetHeavy , 1./sumOfWeights());
105      scale(_histMJetLight , 1./sumOfWeights());
106      // percentage and bin width
107      scale(_histMult, 100./sumOfWeights());
108    }
109
110    ///@}
111
112
113    /// @name Histograms
114    ///@{
115    Histo1DPtr _histThrust, _histMajor, _histMinor, _histOblateness;
116    Histo1DPtr _histJade,_histDurham;
117    Histo1DPtr _histSphericity, _histAplanarity;
118    Histo1DPtr _histC, _histD;
119    Histo1DPtr _histMJetHeavy, _histMJetLight;
120    BinnedHistoPtr<string> _histMult;
121    ///@}
122
123
124  };
125
126
127  RIVET_DECLARE_PLUGIN(L3_1992_I334954);
128
129}