OPAL_1993_S2692198.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #include "Rivet/Analysis.hh" 00003 #include "Rivet/Projections/Beam.hh" 00004 #include "Rivet/Projections/FinalState.hh" 00005 #include "Rivet/Projections/ChargedFinalState.hh" 00006 #include "Rivet/Projections/FastJets.hh" 00007 #include "fastjet/JadePlugin.hh" 00008 00009 namespace Rivet { 00010 00011 00012 /// @brief OPAL photon production 00013 /// @author Peter Richardson 00014 class OPAL_1993_S2692198 : public Analysis { 00015 public: 00016 00017 /// Constructor 00018 OPAL_1993_S2692198() 00019 : Analysis("OPAL_1993_S2692198") 00020 { } 00021 00022 00023 /// @name Analysis methods 00024 //@{ 00025 00026 void analyze(const Event& e) { 00027 // Get event weight for histo filling 00028 const double weight = e.weight(); 00029 00030 // Extract the photons 00031 Particles photons; 00032 Particles nonPhotons; 00033 FourMomentum ptotal; 00034 const FinalState& fs = applyProjection<FinalState>(e, "FS"); 00035 foreach (const Particle& p, fs.particles()) { 00036 ptotal+= p.momentum(); 00037 if (p.pid() == PID::PHOTON) { 00038 photons.push_back(p); 00039 } else { 00040 nonPhotons.push_back(p); 00041 } 00042 } 00043 // No photon return but still count for normalisation 00044 if (photons.empty()) return; 00045 // Definition of the Durham algorithm 00046 fastjet::JetDefinition durham_def(fastjet::ee_kt_algorithm, fastjet::E_scheme, fastjet::Best); 00047 // Definition of the JADE algorithm 00048 fastjet::JadePlugin jade; 00049 fastjet::JetDefinition jade_def = fastjet::JetDefinition(&jade); 00050 // Now for the weird jet algorithm 00051 double evis = ptotal.mass(); 00052 vector<fastjet::PseudoJet> input_particles; 00053 // Pseudo-jets from the non photons 00054 foreach (const Particle& p, nonPhotons) { 00055 const FourMomentum p4 = p.momentum(); 00056 input_particles.push_back(fastjet::PseudoJet(p4.px(), p4.py(), p4.pz(), p4.E())); 00057 } 00058 // Pseudo-jets from all bar the first photon 00059 for (size_t ix = 1; ix < photons.size(); ++ix) { 00060 const FourMomentum p4 = photons[ix].momentum(); 00061 input_particles.push_back(fastjet::PseudoJet(p4.px(), p4.py(), p4.pz(), p4.E())); 00062 } 00063 // Now loop over the photons 00064 for (size_t ix = 0; ix < photons.size(); ++ix) { 00065 FourMomentum pgamma = photons[ix].momentum(); 00066 // Run the jet clustering DURHAM 00067 fastjet::ClusterSequence clust_seq(input_particles, durham_def); 00068 // Cluster the jets 00069 for (size_t j = 0; j < _nPhotonDurham->numBins(); ++j) { 00070 bool accept(true); 00071 double ycut = _nPhotonDurham->bin(j).xMid(); ///< @todo Should this be xMin? 00072 double dcut = sqr(evis)*ycut; 00073 vector<fastjet::PseudoJet> exclusive_jets = sorted_by_E(clust_seq.exclusive_jets(dcut)); 00074 for (size_t iy = 0; iy < exclusive_jets.size(); ++iy) { 00075 FourMomentum pjet(momentum(exclusive_jets[iy])); 00076 double cost = pjet.p3().unit().dot(pgamma.p3().unit()); 00077 double ygamma = 2 * min(sqr(pjet.E()/evis), sqr(pgamma.E()/evis)) * (1 - cost); 00078 if (ygamma < ycut) { 00079 accept = false; 00080 break; 00081 } 00082 } 00083 if (!accept) continue; 00084 _nPhotonDurham->fill(ycut, weight*_nPhotonDurham->bin(j).xWidth()); 00085 size_t njet = min(size_t(4), exclusive_jets.size()) - 1; 00086 if (j < _nPhotonJetDurham[njet]->numBins()) { 00087 _nPhotonJetDurham[njet]->fillBin(j, weight*_nPhotonJetDurham[njet]->bin(j).xWidth()); 00088 } 00089 } 00090 // Run the jet clustering JADE 00091 fastjet::ClusterSequence clust_seq2(input_particles, jade_def); 00092 for (size_t j = 0; j < _nPhotonJade->numBins(); ++j) { 00093 bool accept(true); 00094 double ycut = _nPhotonJade->bin(j).xMid(); ///< @todo Should this be xMin? 00095 double dcut = sqr(evis)*ycut; 00096 vector<fastjet::PseudoJet> exclusive_jets = sorted_by_E(clust_seq2.exclusive_jets(dcut)); 00097 for (size_t iy = 0; iy < exclusive_jets.size(); ++iy) { 00098 FourMomentum pjet(momentum(exclusive_jets[iy])); 00099 double cost = pjet.p3().unit().dot(pgamma.p3().unit()); 00100 double ygamma = 2.*pjet.E()*pgamma.E()/sqr(evis)*(1.-cost); 00101 if (ygamma < ycut) { 00102 accept = false; 00103 break; 00104 } 00105 } 00106 if (!accept) continue; 00107 /// @todo Really want to use a "bar graph" here (i.e. ignore bin width) 00108 _nPhotonJade->fill(ycut, weight*_nPhotonJade->bin(j).xWidth()); 00109 size_t njet = min(size_t(4), exclusive_jets.size()) - 1; 00110 if (j < _nPhotonJetJade[njet]->numBins()) { 00111 _nPhotonJetJade[njet]->fillBin(j, weight*_nPhotonJetJade[njet]->bin(j).xWidth()); 00112 } 00113 } 00114 // Add this photon back in for the next iteration of the loop 00115 if (ix+1 != photons.size()) { 00116 input_particles[nonPhotons.size()+ix] = fastjet::PseudoJet(pgamma.px(), pgamma.py(), pgamma.pz(), pgamma.E()); 00117 } 00118 } 00119 } 00120 00121 00122 void init() { 00123 // Projections 00124 addProjection(FinalState(), "FS"); 00125 00126 // Book datasets 00127 _nPhotonJade = bookHisto1D(1, 1, 1); 00128 _nPhotonDurham = bookHisto1D(2, 1, 1); 00129 for (size_t ix = 0; ix < 4; ++ix) { 00130 _nPhotonJetJade [ix] = bookHisto1D(3 , 1, 1+ix); 00131 _nPhotonJetDurham[ix] = bookHisto1D(4 , 1, 1+ix); 00132 } 00133 } 00134 00135 00136 /// Finalize 00137 void finalize() { 00138 const double fact = 1000/sumOfWeights(); 00139 scale(_nPhotonJade, fact); 00140 scale(_nPhotonDurham, fact); 00141 for (size_t ix = 0; ix < 4; ++ix) { 00142 scale(_nPhotonJetJade [ix],fact); 00143 scale(_nPhotonJetDurham[ix],fact); 00144 } 00145 } 00146 00147 //@} 00148 00149 private: 00150 00151 Histo1DPtr _nPhotonJade; 00152 Histo1DPtr _nPhotonDurham; 00153 Histo1DPtr _nPhotonJetJade [4]; 00154 Histo1DPtr _nPhotonJetDurham[4]; 00155 00156 }; 00157 00158 00159 00160 // The hook for the plugin system 00161 DECLARE_RIVET_PLUGIN(OPAL_1993_S2692198); 00162 00163 } Generated on Thu Mar 10 2016 08:29:52 for The Rivet MC analysis system by ![]() |