OPAL_1993_S2692198.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #include "Rivet/Analysis.hh" 00003 #include "Rivet/RivetYODA.hh" 00004 #include "Rivet/Tools/ParticleIdUtils.hh" 00005 #include "Rivet/Projections/Beam.hh" 00006 #include "Rivet/Projections/FinalState.hh" 00007 #include "Rivet/Projections/ChargedFinalState.hh" 00008 #include "Rivet/Projections/FastJets.hh" 00009 #include "fastjet/JadePlugin.hh" 00010 00011 namespace Rivet { 00012 00013 00014 /// @brief OPAL photon production 00015 /// @author Peter Richardson 00016 class OPAL_1993_S2692198 : public Analysis { 00017 public: 00018 00019 /// Constructor 00020 OPAL_1993_S2692198() : Analysis("OPAL_1993_S2692198") 00021 { 00022 00023 } 00024 00025 00026 /// @name Analysis methods 00027 //@{ 00028 00029 void analyze(const Event& e) { 00030 // First, veto on leptonic events by requiring at least 4 charged FS particles 00031 // const ChargedFinalState& cfs = applyProjection<ChargedFinalState>(e, "CFS"); 00032 // const size_t numParticles = cfs.particles().size(); 00033 00034 // if (numParticles < 4) { 00035 // MSG_DEBUG("Failed ncharged cut"); 00036 // vetoEvent; 00037 // } 00038 // MSG_DEBUG("Passed ncharged cut"); 00039 00040 // Get event weight for histo filling 00041 const double weight = e.weight(); 00042 00043 // extract the photons 00044 ParticleVector photons; 00045 ParticleVector nonPhotons; 00046 FourMomentum ptotal; 00047 const FinalState& fs = applyProjection<FinalState>(e, "FS"); 00048 foreach (const Particle& p, fs.particles()) { 00049 ptotal+= p.momentum(); 00050 if(p.pdgId()==PHOTON) { 00051 photons.push_back(p); 00052 } 00053 else { 00054 nonPhotons.push_back(p); 00055 } 00056 } 00057 // no photon return but still count for normalisation 00058 if(photons.empty()) return; 00059 // definition of the Durham algorithm 00060 fastjet::JetDefinition durham_def(fastjet::ee_kt_algorithm,fastjet::E_scheme, 00061 fastjet::Best); 00062 // definition of the JADE algorithm 00063 fastjet::JadePlugin jade; 00064 fastjet::JetDefinition jade_def = fastjet::JetDefinition(&jade); 00065 // now for the weird jet algorithm 00066 double evis = ptotal.mass(); 00067 vector<fastjet::PseudoJet> input_particles; 00068 // pseudo jets from the non photons 00069 foreach (const Particle& p, nonPhotons) { 00070 input_particles.push_back(fastjet::PseudoJet(p.momentum().px(), 00071 p.momentum().py(), 00072 p.momentum().pz(), 00073 p.momentum().E())); 00074 } 00075 // pseudo jets from all bar the first photon 00076 for(unsigned int ix=1;ix<photons.size();++ix) { 00077 input_particles.push_back(fastjet::PseudoJet(photons[ix].momentum().px(), 00078 photons[ix].momentum().py(), 00079 photons[ix].momentum().pz(), 00080 photons[ix].momentum().E())); 00081 } 00082 // now loop over the photons 00083 for(unsigned int ix=0;ix<photons.size();++ix) { 00084 FourMomentum pgamma = photons[ix].momentum(); 00085 // run the jet clustering DURHAM 00086 fastjet::ClusterSequence clust_seq(input_particles, durham_def); 00087 // cluster the jets 00088 for (size_t j = 0; j < _nPhotonDurham->numBins(); ++j) { 00089 bool accept(true); 00090 double ycut = _nPhotonDurham->bin(j).midpoint(); 00091 double dcut = sqr(evis)*ycut; 00092 vector<fastjet::PseudoJet> exclusive_jets = 00093 sorted_by_E(clust_seq.exclusive_jets(dcut)); 00094 for(unsigned int iy=0;iy<exclusive_jets.size();++iy) { 00095 FourMomentum pjet(momentum(exclusive_jets[iy])); 00096 double cost = pjet.vector3().unit().dot(pgamma.vector3().unit()); 00097 double ygamma = 2.*min(sqr(pjet.E()/evis), 00098 sqr(pgamma.E()/evis))*(1.-cost); 00099 if(ygamma<ycut) { 00100 accept = false; 00101 break; 00102 } 00103 } 00104 if(!accept) continue; 00105 _nPhotonDurham->fill(ycut, weight*_nPhotonDurham->bin(j).width()); 00106 int njet = min(4,int(exclusive_jets.size())) - 1; 00107 if(j<_nPhotonJetDurham[njet]->numBins()) { 00108 _nPhotonJetDurham[njet]->fill(_nPhotonJetDurham[njet]->bin(j).midpoint(), weight*_nPhotonJetDurham[njet]->bin(j).width()); 00109 } 00110 } 00111 // run the jet clustering JADE 00112 fastjet::ClusterSequence clust_seq2(input_particles, jade_def); 00113 for (size_t j = 0; j < _nPhotonJade->numBins(); ++j) { 00114 bool accept(true); 00115 double ycut = _nPhotonJade->bin(j).midpoint(); 00116 double dcut = sqr(evis)*ycut; 00117 vector<fastjet::PseudoJet> exclusive_jets = 00118 sorted_by_E(clust_seq2.exclusive_jets(dcut)); 00119 for(unsigned int iy=0;iy<exclusive_jets.size();++iy) { 00120 FourMomentum pjet(momentum(exclusive_jets[iy])); 00121 double cost = pjet.vector3().unit().dot(pgamma.vector3().unit()); 00122 double ygamma = 2.*pjet.E()*pgamma.E()/sqr(evis)*(1.-cost); 00123 if(ygamma<ycut) { 00124 accept = false; 00125 break; 00126 } 00127 } 00128 if(!accept) continue; 00129 _nPhotonJade->fill(ycut, weight*_nPhotonJade->bin(j).width()); 00130 int njet = min(4,int(exclusive_jets.size())) - 1; 00131 if(j<_nPhotonJetJade[njet]->numBins()) { 00132 _nPhotonJetJade[njet]->fill(_nPhotonJetJade[njet]->bin(j).midpoint(), weight*_nPhotonJetJade[njet]->bin(j).width()); 00133 } 00134 } 00135 // add this photon back in for the next interation of the loop 00136 if(ix+1!=photons.size()) 00137 input_particles[nonPhotons.size()+ix] = 00138 fastjet::PseudoJet(photons[ix].momentum().px(),photons[ix].momentum().py(), 00139 photons[ix].momentum().pz(),photons[ix].momentum().E()); 00140 } 00141 } 00142 00143 00144 void init() { 00145 // Projections 00146 addProjection(FinalState(), "FS"); 00147 // addProjection(ChargedFinalState(), "CFS"); 00148 00149 // Book data sets 00150 _nPhotonJade = bookHisto1D(1, 1, 1); 00151 _nPhotonDurham = bookHisto1D(2, 1, 1); 00152 for(unsigned int ix=0;ix<4;++ix) { 00153 _nPhotonJetJade [ix] = bookHisto1D(3 , 1, 1+ix); 00154 _nPhotonJetDurham[ix] = bookHisto1D(4 , 1, 1+ix); 00155 } 00156 } 00157 00158 00159 /// Finalize 00160 void finalize() { 00161 const double fact = 1000./sumOfWeights(); 00162 scale(_nPhotonJade ,fact); 00163 scale(_nPhotonDurham,fact); 00164 for(unsigned int ix=0;ix<4;++ix) { 00165 scale(_nPhotonJetJade [ix],fact); 00166 scale(_nPhotonJetDurham[ix],fact); 00167 } 00168 } 00169 00170 //@} 00171 00172 private: 00173 00174 Histo1DPtr _nPhotonJade; 00175 Histo1DPtr _nPhotonDurham; 00176 Histo1DPtr _nPhotonJetJade [4]; 00177 Histo1DPtr _nPhotonJetDurham[4]; 00178 00179 }; 00180 00181 00182 00183 // The hook for the plugin system 00184 DECLARE_RIVET_PLUGIN(OPAL_1993_S2692198); 00185 00186 } Generated on Fri Dec 21 2012 15:03:41 for The Rivet MC analysis system by ![]() |