rivet is hosted by Hepforge, IPPP Durham
DELPHI_2003_WUD_03_11.cc
Go to the documentation of this file.
00001 // -*- C++ -*-
00002 #include "Rivet/Analysis.hh"
00003 #include "Rivet/Projections/FastJets.hh"
00004 #include "Rivet/Projections/FinalState.hh"
00005 #include "Rivet/Projections/ChargedFinalState.hh"
00006 
00007 namespace Rivet {
00008 
00009 
00010   namespace {
00011 
00012     /// @name Jet angle calculator functions
00013     //@{
00014     /// @todo Move to utils?
00015 
00016     /// @todo Use Jet or FourMomentum interface rather than PseudoJet
00017     /// @todo Move to utils?
00018     double calc_BZ(const vector<fastjet::PseudoJet>& jets) {
00019       assert(jets.size() == 4);
00020       Vector3 p12 = cross( momentum3(jets[0]), momentum3(jets[1]));
00021       Vector3 p34 = cross( momentum3(jets[2]), momentum3(jets[3]));
00022       return dot(p12,p34) / (p12.mod()*p34.mod());
00023     }
00024 
00025 
00026     /// @todo Use Jet or FourMomentum interface rather than PseudoJet
00027     /// @todo Move to utils?
00028     double calc_KSW(const vector<fastjet::PseudoJet>& jets) {
00029       assert(jets.size() == 4);
00030       Vector3 p13 = cross( momentum3(jets[0]), momentum3(jets[2]));
00031       Vector3 p24 = cross( momentum3(jets[1]), momentum3(jets[3]));
00032       Vector3 p14 = cross( momentum3(jets[0]), momentum3(jets[3]));
00033       Vector3 p23 = cross( momentum3(jets[1]), momentum3(jets[2]));
00034       return cos (0.5*( acos (dot(p14,p23) / (p14.mod()*p23.mod())) +
00035                         acos (dot(p13,p24) / (p13.mod()*p24.mod())) ));
00036     }
00037 
00038 
00039     /// @todo Use Jet or FourMomentum interface rather than PseudoJet
00040     /// @todo Move to utils?
00041     double calc_NR(const vector<fastjet::PseudoJet>& jets) {
00042       assert(jets.size() == 4);
00043       Vector3 p12 = momentum3(jets[0]) - momentum3(jets[1]);
00044       Vector3 p34 = momentum3(jets[2]) - momentum3(jets[3]);
00045       return dot(p12,p34) / (p12.mod()*p34.mod());
00046     }
00047 
00048     /// @todo Use Jet or FourMomentum interface rather than PseudoJet
00049     /// @todo Move to utils?
00050     double calc_ALPHA34(const vector<fastjet::PseudoJet>& jets) {
00051       assert(jets.size() == 4);
00052       Vector3 p3 = momentum3(jets[2]);
00053       Vector3 p4 = momentum3(jets[3]);
00054       return dot(p3,p4) / (p3.mod()*p4.mod());
00055     }
00056 
00057     //@}
00058 
00059   }
00060 
00061 
00062   /**
00063    * @brief DELPHI 4-jet angular distributions
00064    * @author Hendrik Hoeth
00065    *
00066    * This is Hendrik Hoeth's Diploma thesis, measuring the 4-jet angular
00067    * distributions at LEP-1.
00068    *
00069    * @par Run conditions
00070    *
00071    * @arg LEP1 beam energy: \f$ \sqrt{s} = \f$ 91.2 GeV
00072    * @arg Run with generic QCD events.
00073    * @arg No \f$ p_\perp^\text{min} \f$ cutoff is required
00074    */
00075   class DELPHI_2003_WUD_03_11 : public Analysis {
00076   public:
00077 
00078     /// Constructor
00079     DELPHI_2003_WUD_03_11()
00080       : Analysis("DELPHI_2003_WUD_03_11")
00081     {
00082       _numdurjets = 0;
00083       _numjadejets = 0;
00084     }
00085 
00086 
00087     /// @name Analysis methods
00088     //@{
00089 
00090     void init() {
00091       const ChargedFinalState cfs;
00092       addProjection(cfs, "FS");
00093       addProjection(FastJets(cfs, FastJets::JADE, 0.7), "JadeJets");
00094       addProjection(FastJets(cfs, FastJets::DURHAM, 0.7), "DurhamJets");
00095 
00096       _histDurhamBZ      = bookHisto1D(1, 1, 1);
00097       _histDurhamKSW     = bookHisto1D(2, 1, 1);
00098       _histDurhamNR      = bookHisto1D(3, 1, 1);
00099       _histDurhamALPHA34 = bookHisto1D(4, 1, 1);
00100       _histJadeBZ        = bookHisto1D(1, 2, 1);
00101       _histJadeKSW       = bookHisto1D(2, 2, 1);
00102       _histJadeNR        = bookHisto1D(3, 2, 1);
00103       _histJadeALPHA34   = bookHisto1D(4, 2, 1);
00104     }
00105 
00106 
00107     void analyze(const Event& e) {
00108       // First, veto on leptonic events by requiring at least 4 charged FS particles
00109       const FinalState& fs = applyProjection<FinalState>(e, "FS");
00110       const size_t numParticles = fs.particles().size();
00111 
00112       // Even if we only generate hadronic events, we still need a cut on numCharged >= 2.
00113       if (numParticles < 2) {
00114         MSG_DEBUG("Failed multiplicity cut");
00115         vetoEvent;
00116       }
00117       MSG_DEBUG("Passed multiplicity cut");
00118 
00119       // Get event weight for histo filling
00120       const double weight = e.weight();
00121 
00122       // Jets
00123       const FastJets& durjet = applyProjection<FastJets>(e, "DurhamJets");
00124       vector<fastjet::PseudoJet> jets_durham;
00125       if (durjet.clusterSeq()) {
00126         jets_durham = fastjet::sorted_by_E(durjet.clusterSeq()->exclusive_jets_ycut(0.008));
00127         if (jets_durham.size() == 4) {
00128           _histDurhamBZ->fill(fabs(calc_BZ(jets_durham)), weight);
00129           _histDurhamKSW->fill(calc_KSW(jets_durham), weight);
00130           _histDurhamNR->fill(fabs(calc_NR(jets_durham)), weight);
00131           _histDurhamALPHA34->fill(calc_ALPHA34(jets_durham), weight);
00132         }
00133         if (durjet.clusterSeq()->exclusive_ymerge_max(3) > 0.008 &&
00134             durjet.clusterSeq()->exclusive_ymerge_max(4) < 0.008) {
00135           _numdurjets++;
00136         }
00137       }
00138 
00139       const FastJets& jadejet = applyProjection<FastJets>(e, "JadeJets");
00140       vector<fastjet::PseudoJet> jets_jade;
00141       if (jadejet.clusterSeq()) {
00142         jets_jade = fastjet::sorted_by_E(jadejet.clusterSeq()->exclusive_jets_ycut(0.015));
00143         if (jets_jade.size() == 4) {
00144           _histJadeBZ->fill(fabs(calc_BZ(jets_jade)), weight);
00145           _histJadeKSW->fill(calc_KSW(jets_jade), weight);
00146           _histJadeNR->fill(fabs(calc_NR(jets_jade)), weight);
00147           _histJadeALPHA34->fill(calc_ALPHA34(jets_jade), weight);
00148         }
00149         if (jadejet.clusterSeq()->exclusive_ymerge_max(3) > 0.015 &&
00150             jadejet.clusterSeq()->exclusive_ymerge_max(4) < 0.015) {
00151           _numjadejets++;
00152         }
00153       }
00154 
00155     }
00156 
00157 
00158     // Finalize
00159     void finalize() {
00160       // Normalize inclusive single particle distributions to the average number
00161       // of charged particles per event.
00162 
00163       MSG_INFO("Number of Durham jets = " << _numdurjets);
00164       MSG_INFO("Number of Jade jets   = " << _numjadejets);
00165 
00166       /// @todo Scale rather than normalize?
00167       normalize(_histDurhamBZ      , 0.0785);
00168       normalize(_histDurhamKSW     , 0.0785);
00169       normalize(_histDurhamNR      , 0.0785);
00170       normalize(_histDurhamALPHA34 , 0.0785);
00171       normalize(_histJadeBZ        , 0.0277);
00172       normalize(_histJadeKSW       , 0.0277);
00173       normalize(_histJadeNR        , 0.0277);
00174       normalize(_histJadeALPHA34   , 0.0277);
00175     }
00176 
00177     //@}
00178 
00179 
00180   private:
00181 
00182     int _numdurjets;
00183     int _numjadejets;
00184 
00185     /// @name Histograms
00186     //@{
00187     Histo1DPtr _histDurhamBZ;
00188     Histo1DPtr _histDurhamKSW;
00189     Histo1DPtr _histDurhamNR;
00190     Histo1DPtr _histDurhamALPHA34;
00191     Histo1DPtr _histJadeBZ;
00192     Histo1DPtr _histJadeKSW;
00193     Histo1DPtr _histJadeNR;
00194     Histo1DPtr _histJadeALPHA34;
00195     //@}
00196 
00197   };
00198 
00199 
00200 
00201   // The hook for the plugin system
00202   DECLARE_RIVET_PLUGIN(DELPHI_2003_WUD_03_11);
00203 
00204 }