rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

DELPHI_2000_I522656

Events shapes at MZ as function of thrust direction
Experiment: DELPHI (LEP)
Inspire ID: 522656
Status: VALIDATED
Authors:
  • Peter Richardson
References:
  • Eur.Phys.J.C 14 (2000) 557-584, 2000
Beams: e+ e-
Beam energies: (45.6, 45.6) GeV
Run details:
  • e+ e- to hadrons

Measurement of event shapes for different ranges of the angle between the thrust axis and the beam. Jet rate using a mass measure or the Geneva algorithm are not implemented.

Source code: DELPHI_2000_I522656.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/Beam.hh"
  4#include "Rivet/Projections/FinalState.hh"
  5#include "Rivet/Projections/FastJets.hh"
  6#include "Rivet/Projections/Thrust.hh"
  7#include "Rivet/Projections/Sphericity.hh"
  8#include "Rivet/Projections/Hemispheres.hh"
  9#include "Rivet/Projections/ParisiTensor.hh"
 10#include "fastjet/EECambridgePlugin.hh"
 11
 12namespace Rivet {
 13
 14
 15  /// @brief event shapes vs thrust direction
 16  class DELPHI_2000_I522656 : public Analysis {
 17  public:
 18
 19    /// Constructor
 20    RIVET_DEFAULT_ANALYSIS_CTOR(DELPHI_2000_I522656);
 21
 22
 23    /// @name Analysis methods
 24    ///@{
 25
 26    /// Book histograms and initialise projections before the run
 27    void init() {
 28      // Initialise and register projections.
 29      declare(Beam(), "Beams");
 30      const FinalState fs;
 31      declare(fs, "FS");
 32      const Thrust thrust(fs);
 33      declare(thrust, "Thrust");
 34      declare(Sphericity(fs), "Sphericity");
 35      declare(ParisiTensor(fs), "Parisi");
 36      declare(Hemispheres(thrust), "Hemispheres");
 37      declare(FastJets(fs, JetAlg::DURHAM, 0.7), "DurhamJets");
 38      declare(FastJets(fs, JetAlg::JADE  , 0.7), "JadeJets"  );
 39
 40      // book histograms
 41      vector<double> bins={0.00,0.12,0.24,0.36,0.48,0.60,0.72,0.84,0.96};
 42      book(_h_EEC,  bins);
 43      book(_h_AEEC, bins);
 44      book(_h_cone, bins);
 45      book(_h_thrust, {0.12, 0.24, 0.36});
 46      // thrust angle binned
 47      size_t iy=1, ioff=0;
 48      for (size_t ix = 0; ix < _h_EEC->numBins(); ++ix) {
 49        book(_h_EEC->bin(ix+1),  21+ioff, 1, iy);
 50        book(_h_AEEC->bin(ix+1), 25+ioff, 1, iy);
 51        book(_h_cone->bin(ix+1), 29+ioff, 1, iy);
 52        if (ioff==0) {
 53          book(_h_thrust->bin(iy), 33, 1, iy);
 54        }
 55        ++iy;
 56        if (iy==3) {
 57          ++ioff;
 58          iy=1;
 59        }
 60      }
 61      // total values
 62      book(_h_EEC_all   , 3,1,1);
 63      book(_h_AEEC_all  , 4,1,1);
 64      book(_h_cone_all  , 5,1,1);
 65      book(_h_thrust_all, 6,1,1);
 66      book(_h_Oblateness, 7,1,1);
 67      book(_h_C         , 8,1,1);
 68      book(_h_heavy     , 9,1,1);
 69      book(_h_sum       ,10,1,1);
 70      book(_h_diff      ,11,1,1);
 71      book(_h_wide      ,12,1,1);
 72      book(_h_total     ,13,1,1);
 73      book(_h_jade      ,17,1,1);
 74      book(_h_dur       ,18,1,1);
 75      book(_h_cam       ,20,1,1);
 76      book(_h_bin,"/TMP/hbin",bins);
 77    }
 78
 79
 80    /// Perform the per-event analysis
 81    void analyze(const Event& event) {
 82      const FinalState& fs = apply<FinalState>(event, "FS");
 83      if ( fs.particles().size() < 2) vetoEvent;
 84      // Get beams and average beam momentum
 85      const ParticlePair& beams = apply<Beam>(event, "Beams").beams();
 86      // thrust
 87      const Thrust& thrust = apply<Thrust>(event, "Thrust");
 88      // angle bettwen thrust and beam
 89      double cosThrust = abs(beams.first.p3().unit().dot(thrust.thrustAxis()));
 90      _h_bin->fill(cosThrust);
 91      // thrust and related
 92      _h_thrust_all->fill(           1.-thrust.thrust());
 93      _h_thrust->fill(cosThrust, 1.-thrust.thrust());
 94      _h_Oblateness->fill(thrust.oblateness() );
 95
 96      // visible energy and make pseudojets
 97      double Evis = 0.0;
 98      PseudoJets pjs;
 99      for (const Particle& p : fs.particles()) {
100        Evis += p.E();
101        fastjet::PseudoJet pj = p;
102        pjs.push_back(pj);
103      }
104      double Evis2 = sqr(Evis);
105      // (A)EEC
106      // Need iterators since second loop starts at current outer loop iterator, i.e. no "foreach" here!
107      for (Particles::const_iterator p_i = fs.particles().begin(); p_i != fs.particles().end(); ++p_i) {
108        for (Particles::const_iterator p_j = p_i; p_j != fs.particles().end(); ++p_j) {
109          if (p_i == p_j) continue;
110          const Vector3 mom3_i = p_i->momentum().p3();
111          const Vector3 mom3_j = p_j->momentum().p3();
112          const double energy_i = p_i->momentum().E();
113          const double energy_j = p_j->momentum().E();
114          const double thetaij = 180.*mom3_i.unit().angle(mom3_j.unit())/M_PI;
115          double eec = (energy_i*energy_j) / Evis2;
116          eec *= 2.;
117          _h_EEC_all->fill(thetaij, eec);
118          _h_EEC->fill(cosThrust, thetaij, eec);
119          if (thetaij <90.) {
120            _h_AEEC_all->fill(thetaij, -eec);
121            _h_AEEC->fill(cosThrust, thetaij, -eec);
122          }
123          else {
124            _h_AEEC_all->fill(180.-thetaij, eec);
125            _h_AEEC->fill(cosThrust,180.-thetaij, eec);
126          }
127        }
128      }
129      // hemisphere related
130      const Hemispheres& hemi = apply<Hemispheres>(event, "Hemispheres");
131      _h_heavy->fill(hemi.scaledM2high());
132      _h_diff ->fill(hemi.scaledM2diff());
133      _h_sum  ->fill(hemi.scaledM2low()+hemi.scaledM2high());
134      _h_wide ->fill(hemi.Bmax() );
135      _h_total->fill(hemi.Bsum() );
136      // C-parameter
137      const ParisiTensor& parisi = apply<ParisiTensor>(event, "Parisi");
138      _h_C->fill(parisi.C());
139      // jets
140      const FastJets&  durjet = apply<FastJets>(event, "DurhamJets");
141      const FastJets& jadejet = apply<FastJets>(event, "JadeJets");
142      if (durjet .clusterSeq())  _h_dur ->fill( durjet.clusterSeq()->exclusive_ymerge_max(2));
143      if (jadejet.clusterSeq())  _h_jade->fill(jadejet.clusterSeq()->exclusive_ymerge_max(2));
144      // Cambridge is more complicated, inclusive defn
145      for (size_t i = 0; i < _h_cam->numBins(); ++i) {
146        double ycut = _h_cam->bin(i).xMax();
147        // double width = _h_y_2_Cambridge->bin(i).xWidth();
148        fastjet::EECambridgePlugin plugin(ycut);
149        fastjet::JetDefinition jdef(&plugin);
150        fastjet::ClusterSequence cseq(pjs, jdef);
151        unsigned int njet = cseq.inclusive_jets().size();
152        if (njet==2) {
153          _h_cam->fill(_h_cam->bin(i).xMid());
154          break;
155        }
156      }
157      // jet cone
158      Vector3 jetAxis=thrust.thrustAxis();
159      if (hemi.highMassDirection()) jetAxis *=-1.;
160      for (const Particle& p : fs.particles()) {
161        const double thetaij = 180.*jetAxis.angle(p.p3().unit())/M_PI;
162        double jcef = p.E()/ Evis;
163        _h_cone_all->fill(thetaij,jcef);
164        _h_cone->fill(cosThrust,thetaij,jcef);
165      }
166    }
167
168
169    /// Normalise histograms etc., after the run
170    void finalize() {
171      for (size_t ix = 0; ix < _h_EEC->numBins(); ++ix) {
172        if (ix<2) scale(_h_thrust->bin(ix+1), 1./_h_bin->bin(ix).sumW());
173        scale(_h_EEC->bin(ix+1),  180./M_PI/_h_bin->bin(ix).sumW());
174        scale(_h_AEEC->bin(ix+1), 180./M_PI/_h_bin->bin(ix).sumW());
175        scale(_h_cone->bin(ix+1), 180./M_PI/_h_bin->bin(ix).sumW());
176      }
177
178      scale(_h_thrust_all, 1./sumOfWeights());
179      scale(_h_EEC_all, 180./M_PI/sumOfWeights());
180      scale(_h_AEEC_all, 180./M_PI/sumOfWeights());
181      scale(_h_cone_all, 180./M_PI/sumOfWeights());
182      scale(_h_Oblateness, 1./sumOfWeights());
183      scale(_h_C         , 1./sumOfWeights());
184      scale(_h_heavy     , 1./sumOfWeights());
185      scale(_h_sum       , 1./sumOfWeights());
186      scale(_h_diff      , 1./sumOfWeights());
187      scale(_h_wide      , 1./sumOfWeights());
188      scale(_h_total     , 1./sumOfWeights());
189      scale(_h_dur       , 1./sumOfWeights());
190      scale(_h_jade      , 1./sumOfWeights());
191      scale(_h_cam       , 1./sumOfWeights());
192    }
193
194    ///@}
195
196
197    /// @name Histograms
198    ///@{
199    Histo1DPtr _h_thrust_all,_h_EEC_all,_h_AEEC_all,_h_cone_all;
200    Histo1DPtr _h_Oblateness,_h_C,_h_heavy,_h_sum,_h_diff,_h_wide,_h_total;
201    Histo1DPtr _h_jade,_h_dur,_h_cam;
202    Histo1DGroupPtr _h_thrust, _h_EEC, _h_AEEC, _h_cone;
203    Histo1DPtr _h_bin;
204    ///@}
205
206
207  };
208
209
210  RIVET_DECLARE_PLUGIN(DELPHI_2000_I522656);
211
212}